#include "TROOT.h"
#include "TBackCompFitter.h"
#include "TMethodCall.h"
#include "TInterpreter.h"
#include "Math/Util.h"
#include <iostream>
#include <cassert>
#include "Math/IParamFunction.h"
#include "TH1.h"
#include "TH2.h"
#include "TH3.h"
#include "TMath.h"
#include "TGraph.h"
#include "TGraphErrors.h"
#include "TGraph2DErrors.h"
#include "TMultiGraph.h"
#include "HFitInterface.h"
#include "Math/Minimizer.h"
#include "Fit/BinData.h"
#include "Fit/UnBinData.h"
#include "Fit/PoissonLikelihoodFCN.h"
#include "Fit/LogLikelihoodFCN.h"
#include "Fit/Chi2FCN.h"
#include "Fit/FcnAdapter.h"
#include "TFitResult.h"
ClassImp(TBackCompFitter);
TBackCompFitter::TBackCompFitter( ) :
fMinimizer(0),
fObjFunc(0),
fModelFunc(0)
{
SetName("BCFitter");
}
TBackCompFitter::TBackCompFitter(std::auto_ptr<ROOT::Fit::Fitter> fitter, std::auto_ptr<ROOT::Fit::FitData> data) :
fFitData(data),
fFitter(fitter),
fMinimizer(0),
fObjFunc(0),
fModelFunc(0)
{
SetName("LastFitter");
}
TBackCompFitter::~TBackCompFitter() {
if (fMinimizer) delete fMinimizer;
if (fObjFunc) delete fObjFunc;
if (fModelFunc) delete fModelFunc;
}
Double_t TBackCompFitter::Chisquare(Int_t npar, Double_t *params) const {
const std::vector<double> & minpar = fFitter->Result().Parameters();
assert (npar == (int) minpar.size() );
double diff = 0;
double s = 0;
for (int i =0; i < npar; ++i) {
diff += std::abs( params[i] - minpar[i] );
s += minpar[i];
}
if (diff > s * 1.E-12 ) Warning("Chisquare","given parameter values are not at minimum - chi2 at minimum is returned");
return fFitter->Result().Chi2();
}
void TBackCompFitter::Clear(Option_t*) {
}
Int_t TBackCompFitter::ExecuteCommand(const char *command, Double_t *args, Int_t nargs){
#ifdef DEBUG
std::cout<<"Execute command= "<<command<<std::endl;
#endif
DoSetDimension();
TString scommand(command);
scommand.ToUpper();
if (scommand.Contains("MIG")) {
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad");
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
if (scommand.Contains("MINI")) {
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Minimize");
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
if (scommand.Contains("SIM")) {
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Simplex");
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
if (scommand.Contains("SCA")) {
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Scan");
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
else if (scommand.Contains("MINO")) {
if (fFitter->Config().MinosErrors() ) return 0;
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
fFitter->Config().SetMinosErrors(true);
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad");
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
else if (scommand.Contains("HES")) {
if (fFitter->Config().ParabErrors() ) return 0;
if (!fObjFunc) {
Error("ExecuteCommand","FCN must set before executing this command");
return -1;
}
fFitter->Config().SetParabErrors(true);
fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad");
bool ret = fFitter->FitFCN(*fObjFunc);
return (ret) ? 0 : -1;
}
else if (scommand.Contains("FIX")) {
for(int i = 0; i < nargs; i++) {
FixParameter(int(args[i])-1);
}
return 0;
}
else if (scommand.Contains("SET LIM")) {
if (nargs < 3) {
Error("ExecuteCommand","Invalid parameters given in SET LIMIT");
return -1;
}
int ipar = int(args[0]);
if (!ValidParameterIndex(ipar) ) return -1;
double low = args[1];
double up = args[2];
fFitter->Config().ParSettings(ipar).SetLimits(low,up);
return 0;
}
else if (scommand.Contains("SET PRIN")) {
if (nargs < 1) return -1;
fFitter->Config().MinimizerOptions().SetPrintLevel(int(args[0]) );
return 0;
}
else if (scommand.Contains("SET ERR")) {
if (nargs < 1) return -1;
fFitter->Config().MinimizerOptions().SetPrintLevel(int( args[0]) );
return 0;
}
else if (scommand.Contains("SET STR")) {
if (nargs < 1) return -1;
fFitter->Config().MinimizerOptions().SetStrategy(int(args[0]) );
return 0;
}
else if (scommand.Contains("SET GRA")) {
return -1;
}
else if (scommand.Contains("SET NOW")) {
return -1;
}
else if (scommand.Contains("CALL FCN")) {
if (nargs < 1 || fFCN == 0 ) return -1;
int npar = fObjFunc->NDim();
std::vector<double> params(npar);
for (int i = 0; i < npar; ++i)
params[i] = GetParameter(i);
double fval = 0;
(*fFCN)(npar, 0, fval, ¶ms[0],int(args[0]) ) ;
return 0;
}
else {
Error("ExecuteCommand","Invalid or not supported command given %s",command);
return -1;
}
}
bool TBackCompFitter::ValidParameterIndex(int ipar) const {
int nps = fFitter->Config().ParamsSettings().size();
if (ipar < 0 || ipar >= nps ) {
std::string msg = ROOT::Math::Util::ToString(ipar) + " is an invalid Parameter index";
Error("ValidParameterIndex","%s",msg.c_str());
return false;
}
return true;
}
void TBackCompFitter::FixParameter(Int_t ipar) {
if (ValidParameterIndex(ipar) )
fFitter->Config().ParSettings(ipar).Fix();
}
void TBackCompFitter::GetConfidenceIntervals(Int_t n, Int_t ndim, const Double_t *x, Double_t *ci, Double_t cl)
{
if (!fFitter->Result().IsValid()) {
Error("GetConfidenceIntervals","Cannot compute confidence intervals with an invalide fit result");
return;
}
fFitter->Result().GetConfidenceIntervals(n,ndim,1,x,ci,cl);
}
void TBackCompFitter::GetConfidenceIntervals(TObject *obj, Double_t cl)
{
if (!fFitter->Result().IsValid() ) {
Error("GetConfidenceIntervals","Cannot compute confidence intervals with an invalide fit result");
return;
}
int datadim = 1;
TObject * fitobj = GetObjectFit();
if (!fitobj) {
Error("GetConfidenceIntervals","Cannot compute confidence intervals without a fitting object");
return;
}
if (fitobj->InheritsFrom(TGraph2D::Class())) datadim = 2;
if (fitobj->InheritsFrom(TH1::Class())) {
TH1 * h1 = dynamic_cast<TH1*>(fitobj);
assert(h1 != 0);
datadim = h1->GetDimension();
}
if (datadim == 1) {
if (!obj->InheritsFrom(TGraphErrors::Class()) && !obj->InheritsFrom(TH1::Class() ) ) {
Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TGraphErrors or a TH1");
return;
}
}
if (datadim == 2) {
if (!obj->InheritsFrom(TGraph2DErrors::Class()) && !obj->InheritsFrom(TH2::Class() ) ) {
Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TGraph2DErrors or a TH2");
return;
}
}
if (datadim == 3) {
if (!obj->InheritsFrom(TH3::Class() ) ) {
Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TH3");
return;
}
}
ROOT::Fit::BinData data;
data.Opt().fUseEmpty = true;
if (obj->InheritsFrom(TGraph::Class()) )
ROOT::Fit::FillData(data, dynamic_cast<TGraph *>(obj) );
else if (obj->InheritsFrom(TGraph2D::Class()) )
ROOT::Fit::FillData(data, dynamic_cast<TGraph2D *>(obj) );
else if (obj->InheritsFrom(TH1::Class()) )
ROOT::Fit::FillData(data, dynamic_cast<TH1 *>(obj) );
unsigned int n = data.Size();
std::vector<double> ci( n );
fFitter->Result().GetConfidenceIntervals(data,&ci[0],cl);
const ROOT::Math::IParamMultiFunction * func = fFitter->Result().FittedFunction();
assert(func != 0);
for (unsigned int i = 0; i < n; ++i) {
const double * x = data.Coords(i);
double y = (*func)( x );
if (obj->InheritsFrom(TGraphErrors::Class()) ) {
TGraphErrors * gr = dynamic_cast<TGraphErrors *> (obj);
assert(gr != 0);
gr->SetPoint(i, *x, y);
gr->SetPointError(i, 0, ci[i]);
}
if (obj->InheritsFrom(TGraph2DErrors::Class()) ) {
TGraph2DErrors * gr = dynamic_cast<TGraph2DErrors *> (obj);
assert(gr != 0);
gr->SetPoint(i, x[0], x[1], y);
gr->SetPointError(i, 0, 0, ci[i]);
}
if (obj->InheritsFrom(TH1::Class()) ) {
TH1 * h1 = dynamic_cast<TH1 *> (obj);
assert(h1 != 0);
int ibin = 0;
if (datadim == 1) ibin = h1->FindBin(*x);
if (datadim == 2) ibin = h1->FindBin(x[0],x[1]);
if (datadim == 3) ibin = h1->FindBin(x[0],x[1],x[2]);
h1->SetBinContent(ibin, y);
h1->SetBinError(ibin, ci[i]);
}
}
}
Double_t* TBackCompFitter::GetCovarianceMatrix() const {
unsigned int nfreepar = GetNumberFreeParameters();
unsigned int ntotpar = GetNumberTotalParameters();
if (fCovar.size() != nfreepar*nfreepar )
fCovar.resize(nfreepar*nfreepar);
if (!fFitter->Result().IsValid() ) {
Warning("GetCovarianceMatrix","Invalid fit result");
return 0;
}
unsigned int l = 0;
for (unsigned int i = 0; i < ntotpar; ++i) {
if (fFitter->Config().ParSettings(i).IsFixed() ) continue;
unsigned int m = 0;
for (unsigned int j = 0; j < ntotpar; ++j) {
if (fFitter->Config().ParSettings(j).IsFixed() ) continue;
unsigned int index = nfreepar*l + m;
assert(index < fCovar.size() );
fCovar[index] = fFitter->Result().CovMatrix(i,j);
m++;
}
l++;
}
return &(fCovar.front());
}
Double_t TBackCompFitter::GetCovarianceMatrixElement(Int_t i, Int_t j) const {
unsigned int np2 = fCovar.size();
unsigned int npar = GetNumberFreeParameters();
if ( np2 == 0 || np2 != npar *npar ) {
double * c = GetCovarianceMatrix();
if (c == 0) return 0;
}
return fCovar[i*npar + j];
}
Int_t TBackCompFitter::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const {
if (!ValidParameterIndex(ipar) ) return -1;
const ROOT::Fit::FitResult & result = fFitter->Result();
if (!result.IsValid() ) {
Warning("GetErrors","Invalid fit result");
return -1;
}
eparab = result.Error(ipar);
eplus = result.UpperError(ipar);
eminus = result.LowerError(ipar);
globcc = result.GlobalCC(ipar);
return 0;
}
Int_t TBackCompFitter::GetNumberTotalParameters() const {
return fFitter->Result().NTotalParameters();
}
Int_t TBackCompFitter::GetNumberFreeParameters() const {
return fFitter->Result().NFreeParameters();
}
Double_t TBackCompFitter::GetParError(Int_t ipar) const {
if (fFitter->Result().IsEmpty() ) {
if (ValidParameterIndex(ipar) ) return fFitter->Config().ParSettings(ipar).StepSize();
else return 0;
}
return fFitter->Result().Error(ipar);
}
Double_t TBackCompFitter::GetParameter(Int_t ipar) const {
if (fFitter->Result().IsEmpty() ) {
if (ValidParameterIndex(ipar) ) return fFitter->Config().ParSettings(ipar).Value();
else return 0;
}
return fFitter->Result().Value(ipar);
}
Int_t TBackCompFitter::GetParameter(Int_t ipar,char *name,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const {
if (!ValidParameterIndex(ipar) ) {
return -1;
}
const std::string & pname = fFitter->Config().ParSettings(ipar).Name();
const char * c = pname.c_str();
std::copy(c,c + pname.size(),name);
if (fFitter->Result().IsEmpty() ) {
value = fFitter->Config().ParSettings(ipar).Value();
verr = fFitter->Config().ParSettings(ipar).Value();
vlow = fFitter->Config().ParSettings(ipar).LowerLimit();
vhigh = fFitter->Config().ParSettings(ipar).UpperLimit();
return 1;
}
else {
value = fFitter->Result().Value(ipar);
verr = fFitter->Result().Error(ipar);
vlow = fFitter->Result().LowerError(ipar);
vhigh = fFitter->Result().UpperError(ipar);
}
return 0;
}
const char *TBackCompFitter::GetParName(Int_t ipar) const {
if (!ValidParameterIndex(ipar) ) {
return 0;
}
return fFitter->Config().ParSettings(ipar).Name().c_str();
}
Int_t TBackCompFitter::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const {
const ROOT::Fit::FitResult & result = fFitter->Result();
amin = result.MinFcnValue();
edm = result.Edm();
errdef = fFitter->Config().MinimizerOptions().ErrorDef();
nvpar = result.NFreeParameters();
nparx = result.NTotalParameters();
return 0;
}
Double_t TBackCompFitter::GetSumLog(Int_t) {
Warning("GetSumLog","Dummy method - returned 0");
return 0.;
}
Bool_t TBackCompFitter::IsFixed(Int_t ipar) const {
if (!ValidParameterIndex(ipar) ) {
return false;
}
return fFitter->Config().ParSettings(ipar).IsFixed();
}
void TBackCompFitter::PrintResults(Int_t level, Double_t ) const {
if (fFitter->GetMinimizer() && fFitter->Config().MinimizerType() == "Minuit")
fFitter->GetMinimizer()->PrintResults();
else {
if (level > 0) fFitter->Result().Print(std::cout);
if (level > 1) fFitter->Result().PrintCovMatrix(std::cout);
}
}
void TBackCompFitter::ReleaseParameter(Int_t ipar) {
if (ValidParameterIndex(ipar) )
fFitter->Config().ParSettings(ipar).Release();
}
void TBackCompFitter::SetFitMethod(const char *) {
Info("SetFitMethod","non supported method");
}
Int_t TBackCompFitter::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh) {
std::vector<ROOT::Fit::ParameterSettings> & parlist = fFitter->Config().ParamsSettings();
if ( ipar >= (int) parlist.size() ) parlist.resize(ipar+1);
ROOT::Fit::ParameterSettings ps(parname, value, verr);
if (verr == 0) ps.Fix();
if (vlow < vhigh) ps.SetLimits(vlow, vhigh);
parlist[ipar] = ps;
return 0;
}
void TBackCompFitter::ReCreateMinimizer() {
assert(fFitData.get());
if (fFitter->Result().FittedFunction() != 0) {
if (fModelFunc) delete fModelFunc;
fModelFunc = dynamic_cast<ROOT::Math::IParamMultiFunction *>((fFitter->Result().FittedFunction())->Clone());
assert(fModelFunc);
const ROOT::Fit::BinData * bindata = dynamic_cast<const ROOT::Fit::BinData *>(fFitData.get());
if (bindata) {
if (GetFitOption().Like )
fObjFunc = new ROOT::Fit::PoissonLikelihoodFCN<ROOT::Math::IMultiGenFunction>(*bindata, *fModelFunc);
else
fObjFunc = new ROOT::Fit::Chi2FCN<ROOT::Math::IMultiGenFunction>(*bindata, *fModelFunc);
}
else {
const ROOT::Fit::UnBinData * unbindata = dynamic_cast<const ROOT::Fit::UnBinData *>(fFitData.get());
assert(unbindata);
fObjFunc = new ROOT::Fit::LogLikelihoodFCN<ROOT::Math::IMultiGenFunction>(*unbindata, *fModelFunc);
}
}
fMinimizer = fFitter->Config().CreateMinimizer();
if (fMinimizer == 0) {
Error("SetMinimizerFunction","cannot create minimizer %s",fFitter->Config().MinimizerType().c_str() );
}
else {
if (!fObjFunc) {
Error("SetMinimizerFunction","Object Function pointer is NULL");
}
else
fMinimizer->SetFunction(*fObjFunc);
}
}
void TBackCompFitter::SetFCN(void (*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
{
fFCN = fcn;
if (fObjFunc) delete fObjFunc;
fObjFunc = new ROOT::Fit::FcnAdapter(fFCN);
DoSetDimension();
}
void InteractiveFCNm2(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
{
TMethodCall *m = (TVirtualFitter::GetFitter())->GetMethodCall();
if (!m) return;
Long_t args[5];
args[0] = (Long_t)∦
args[1] = (Long_t)gin;
args[2] = (Long_t)&f;
args[3] = (Long_t)u;
args[4] = (Long_t)flag;
m->SetParamPtrs(args);
Double_t result;
m->Execute(result);
}
void TBackCompFitter::SetFCN(void *fcn)
{
if (!fcn) return;
const char *funcname = gCint->Getp2f2funcname(fcn);
if (funcname) {
fMethodCall = new TMethodCall();
fMethodCall->InitWithPrototype(funcname,"Int_t&,Double_t*,Double_t&,Double_t*,Int_t");
}
fFCN = InteractiveFCNm2;
TVirtualFitter::SetFitter(this);
if (fObjFunc) delete fObjFunc;
fObjFunc = new ROOT::Fit::FcnAdapter(fFCN);
DoSetDimension();
}
void TBackCompFitter::SetObjFunction(ROOT::Math::IMultiGenFunction * fcn) {
if (fObjFunc) delete fObjFunc;
fObjFunc = fcn->Clone();
}
void TBackCompFitter::DoSetDimension() {
if (!fObjFunc) return;
ROOT::Fit::FcnAdapter * fobj = dynamic_cast<ROOT::Fit::FcnAdapter*>(fObjFunc);
assert(fobj != 0);
int ndim = fFitter->Config().ParamsSettings().size();
if (ndim != 0) fobj->SetDimension(ndim);
}
ROOT::Math::IMultiGenFunction * TBackCompFitter::GetObjFunction( ) const {
if (fObjFunc) return fObjFunc;
return fFitter->GetFCN();
}
ROOT::Math::Minimizer * TBackCompFitter::GetMinimizer( ) const {
if (fMinimizer) return fMinimizer;
return fFitter->GetMinimizer();
}
TFitResult * TBackCompFitter::GetTFitResult( ) const {
if (!fFitter.get() ) return 0;
return new TFitResult( fFitter->Result() );
}
bool TBackCompFitter::Scan(unsigned int ipar, TGraph * gr, double xmin, double xmax )
{
if (!gr) return false;
ROOT::Math::Minimizer * minimizer = fFitter->GetMinimizer();
if (!minimizer) {
Error("Scan","Minimizer is not available - cannot scan before fitting");
return false;
}
unsigned int npoints = gr->GetN();
if (npoints == 0) {
npoints = 40;
gr->Set(npoints);
}
bool ret = minimizer->Scan( ipar, npoints, gr->GetX(), gr->GetY(), xmin, xmax);
if ((int) npoints < gr->GetN() ) gr->Set(npoints);
return ret;
}
bool TBackCompFitter::Contour(unsigned int ipar, unsigned int jpar, TGraph * gr, double confLevel) {
if (!gr) return false;
ROOT::Math::Minimizer * minimizer = fFitter->GetMinimizer();
if (!minimizer) {
Error("Scan","Minimizer is not available - cannot scan before fitting");
return false;
}
double upScale = fFitter->Config().MinimizerOptions().ErrorDef();
double upVal = TMath::ChisquareQuantile( confLevel, 2);
minimizer->SetErrorDef (upScale * upVal);
unsigned int npoints = gr->GetN();
if (npoints == 0) {
npoints = 40;
gr->Set(npoints);
}
bool ret = minimizer->Contour( ipar, jpar, npoints, gr->GetX(), gr->GetY());
if ((int) npoints < gr->GetN() ) gr->Set(npoints);
minimizer->SetErrorDef ( upScale);
return ret;
}