// @(#)root/hist:$Name: $:$Id: TF1.cxx,v 1.108 2005/09/02 19:18:11 brun Exp $ // Author: Rene Brun 18/08/95 /************************************************************************* * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. * * All rights reserved. * * * * For the licensing terms see $ROOTSYS/LICENSE. * * For the list of contributors see $ROOTSYS/README/CREDITS. * *************************************************************************/ #include "Riostream.h" #include "TROOT.h" #include "TMath.h" #include "TF1.h" #include "TH1.h" #include "TVirtualPad.h" #include "TStyle.h" #include "TRandom.h" #include "Api.h" #include "TPluginManager.h" #include "TVirtualUtilPad.h" #include "TBrowser.h" #include "TColor.h" Bool_t TF1::fgAbsValue = kFALSE; Bool_t TF1::fgRejectPoint = kFALSE; static TF1 *gHelper = 0; static Double_t gErrorTF1 = 0; ClassImp(TF1) //______________________________________________________________________________ // // a TF1 object is a 1-Dim function defined between a lower and upper limit. // The function may be a simple function (see TFormula) or a precompiled // user function. // The function may have associated parameters. // TF1 graphics function is via the TH1/TGraph drawing functions. // // The following types of functions can be created: // A- Expression using variable x and no parameters // B- Expression using variable x with parameters // C- A general C function with parameters // // +++++++++++++++++++++++++++++++++++ // ===> + Example of a function of type A + // +++++++++++++++++++++++++++++++++++ // // Case A1 (inline expression using standard C++ functions/operators) // ------------------------------------------------------------------ // TF1 *fa1 = new TF1("fa1","sin(x)/x",0,10); // fa1->Draw(); ///* */ // // // Case A2 (inline expression using TMath functions without parameters) // -------------------------------------------------------------------- // TF1 *fa2 = new TF1("fa2","TMath::DiLog(x)",0,10); // fa2->Draw(); // // Case A3 (inline expression using a CINT function by name // -------------------------------------------------------- // Double_t myFunc(x) { // return x+sin(x); // } // TF1 *fa3 = new TF1("fa4","myFunc(x)",-3,5); // fa3->Draw(); // // // +++++++++++++++++++++++++++++++++++ // ===> + Example of a function of type B+ // +++++++++++++++++++++++++++++++++++ // // Case B1 (inline expression using standard C++ functions/operators) // ------------------------------------------------------------------ // TF1 *f1 = new TF1("f1","[0]*x*sin([1]*x)",-3,3); // This creates a function of variable x with 2 parameters. // The parameters must be initialized via: // f1->SetParameter(0,value_first_parameter); // f1->SetParameter(1,value_second_parameter); // Parameters may be given a name: // f1->SetParName(0,"Constant"); // // Case B2 (inline expression using TMath functions with parameters) // -------------------------------------------------------------------- // TF1 *fb2 = new TF1("fa3","TMath::Landau(x,[0],[1],0)",-5,10); // fb2->SetParameters(0.2,1.3); // fb2->Draw(); // // // +++++++++++++++++++++++++++++++++++ // ===> + Example of a function of type C+ // +++++++++++++++++++++++++++++++++++ // // Consider the macro myfunc.C below //-------------macro myfunc.C----------------------------- //Double_t myfunction(Double_t *x, Double_t *par) //{ // Float_t xx =x[0]; // Double_t f = TMath::Abs(par[0]*sin(par[1]*xx)/xx); // return f; //} //void myfunc() //{ // TF1 *f1 = new TF1("myfunc",myfunction,0,10,2); // f1->SetParameters(2,1); // f1->SetParNames("constant","coefficient"); // f1->Draw(); //} //void myfit() //{ // TH1F *h1=new TH1F("h1","test",100,0,10); // h1->FillRandom("myfunc",20000); // TF1 *f1=gROOT->GetFunction("myfunc"); // f1->SetParameters(800,1); // h1.Fit("myfunc"); //} //--------end of macro myfunc.C--------------------------------- // // In an interactive session you can do: // Root > .L myfunc.C // Root > myfunc(); // Root > myfit(); // // // TF1 objects can reference other TF1 objects (thanks John Odonnell) // of type A or B defined above.This excludes CINT interpreted functions // and compiled functions. // However, there is a restriction. A function cannot reference a basic // function if the basic function is a polynomial polN. // Example: //{ // TF1 *fcos = new TF1 ("fcos", "[0]*cos(x)", 0., 10.); // fcos->SetParNames( "cos"); // fcos->SetParameter( 0, 1.1); // // TF1 *fsin = new TF1 ("fsin", "[0]*sin(x)", 0., 10.); // fsin->SetParNames( "sin"); // fsin->SetParameter( 0, 2.1); // // TF1 *fsincos = new TF1 ("fsc", "fcos+fsin"); // // TF1 *fs2 = new TF1 ("fs2", "fsc+fsc"); //} // // WHY TF1 CANNOT ACCEPT A CLASS MEMBER FUNCTION ? // =============================================== // This is a frequently asked question. // C++ is a strongly typed language. There is no way for TF1 (without // recompiling this class) to know about all possible user defined data types. // This also apply to the case of a static class function. // //------------------------------------------------------------------------ TF1 *TF1::fgCurrent = 0; //______________________________________________________________________________ TF1::TF1(): TFormula(), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*-*-*-*-*F1 default constructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ====================== fXmin = 0; fXmax = 0; fNpx = 100; fType = 0; fNpfits = 0; fNDF = 0; fNsave = 0; fChisquare = 0; fIntegral = 0; fFunction = 0; fParErrors = 0; fParMin = 0; fParMax = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fParent = 0; fSave = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fMethodCall = 0; SetFillStyle(0); } //______________________________________________________________________________ TF1::TF1(const char *name,const char *formula, Double_t xmin, Double_t xmax) :TFormula(name,formula), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*F1 constructor using a formula definition*-*-*-*-*-*-*-*-*-*-* //*-* ========================================= //*-* //*-* See TFormula constructor for explanation of the formula syntax. //*-* //*-* See tutorials: fillrandom, first, fit1, formula1, multifit //*-* for real examples. //*-* //*-* Creates a function of type A or B between xmin and xmax //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* if (xmin < xmax ) { fXmin = xmin; fXmax = xmax; } else { fXmin = xmax; //when called from TF2,TF3 fXmax = xmin; } fNpx = 100; fType = 0; fFunction = 0; if (fNpar) { fParErrors = new Double_t[fNpar]; fParMin = new Double_t[fNpar]; fParMax = new Double_t[fNpar]; for (int i = 0; i < fNpar; i++) { fParErrors[i] = 0; fParMin[i] = 0; fParMax[i] = 0; } } else { fParErrors = 0; fParMin = 0; fParMax = 0; } fChisquare = 0; fIntegral = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fParent = 0; fNpfits = 0; fNDF = 0; fNsave = 0; fSave = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fMethodCall = 0; if (fNdim != 1 && xmin < xmax) { Error("TF1","function: %s/%s has %d parameters instead of 1",name,formula,fNdim); MakeZombie(); } if (!gStyle) return; SetLineColor(gStyle->GetFuncColor()); SetLineWidth(gStyle->GetFuncWidth()); SetLineStyle(gStyle->GetFuncStyle()); SetFillStyle(0); } //______________________________________________________________________________ TF1::TF1(const char *name, Double_t xmin, Double_t xmax, Int_t npar) :TFormula(), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*F1 constructor using name of an interpreted function*-*-*-* //*-* ======================================================= //*-* //*-* Creates a function of type C between xmin and xmax. //*-* name is the name of an interpreted CINT cunction. //*-* The function is defined with npar parameters //*-* fcn must be a function of type: //*-* Double_t fcn(Double_t *x, Double_t *params) //*-* //*-* This constructor is called for functions of type C by CINT. //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* fXmin = xmin; fXmax = xmax; fNpx = 100; fType = 2; fFunction = 0; if (npar > 0 ) fNpar = npar; if (fNpar) { fNames = new TString[fNpar]; fParams = new Double_t[fNpar]; fParErrors = new Double_t[fNpar]; fParMin = new Double_t[fNpar]; fParMax = new Double_t[fNpar]; for (int i = 0; i < fNpar; i++) { fParams[i] = 0; fParErrors[i] = 0; fParMin[i] = 0; fParMax[i] = 0; } } else { fParErrors = 0; fParMin = 0; fParMax = 0; } fChisquare = 0; fIntegral = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fParent = 0; fNpfits = 0; fNDF = 0; fNsave = 0; fSave = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fMethodCall = 0; fNdim = 1; TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name); if (f1old) delete f1old; SetName(name); if (gStyle) { SetLineColor(gStyle->GetFuncColor()); SetLineWidth(gStyle->GetFuncWidth()); SetLineStyle(gStyle->GetFuncStyle()); } SetFillStyle(0); SetTitle(name); if (name) { if (*name == '*') return; //case happens via SavePrimitive fMethodCall = new TMethodCall(); fMethodCall->InitWithPrototype(name,"Double_t*,Double_t*"); fNumber = -1; gROOT->GetListOfFunctions()->Add(this); if (! fMethodCall->IsValid() ) { Error("TF1","No function found with the signature %s(Double_t*,Double_t*)",name); } } else { Error("TF1","requires a proper function name!"); } } //______________________________________________________________________________ TF1::TF1(const char *name,void *fcn, Double_t xmin, Double_t xmax, Int_t npar) :TFormula(), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*F1 constructor using pointer to an interpreted function*-*-*-* //*-* ======================================================= //*-* //*-* See TFormula constructor for explanation of the formula syntax. //*-* //*-* Creates a function of type C between xmin and xmax. //*-* The function is defined with npar parameters //*-* fcn must be a function of type: //*-* Double_t fcn(Double_t *x, Double_t *params) //*-* //*-* see tutorial; myfit for an example of use //*-* also test/stress.cxx (see function stress1) //*-* //*-* //*-* This constructor is called for functions of type C by CINT. //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* fXmin = xmin; fXmax = xmax; fNpx = 100; fType = 2; fFunction = 0; if (npar > 0 ) fNpar = npar; if (fNpar) { fNames = new TString[fNpar]; fParams = new Double_t[fNpar]; fParErrors = new Double_t[fNpar]; fParMin = new Double_t[fNpar]; fParMax = new Double_t[fNpar]; for (int i = 0; i < fNpar; i++) { fParams[i] = 0; fParErrors[i] = 0; fParMin[i] = 0; fParMax[i] = 0; } } else { fParErrors = 0; fParMin = 0; fParMax = 0; } fChisquare = 0; fIntegral = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fParent = 0; fNpfits = 0; fNDF = 0; fNsave = 0; fSave = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fMethodCall = 0; fNdim = 1; TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name); if (f1old) delete f1old; SetName(name); if (gStyle) { SetLineColor(gStyle->GetFuncColor()); SetLineWidth(gStyle->GetFuncWidth()); SetLineStyle(gStyle->GetFuncStyle()); } SetFillStyle(0); if (!fcn) return; char *funcname = G__p2f2funcname(fcn); SetTitle(funcname); if (funcname) { fMethodCall = new TMethodCall(); fMethodCall->InitWithPrototype(funcname,"Double_t*,Double_t*"); fNumber = -1; gROOT->GetListOfFunctions()->Add(this); if (! fMethodCall->IsValid() ) { Error("TF1","No function found with the signature %s(Double_t*,Double_t*)",funcname); } } else { Error("TF1","can not find any function at the address 0x%x. This function requested for %s",fcn,name); } } //______________________________________________________________________________ TF1::TF1(const char *name,Double_t (*fcn)(Double_t *, Double_t *), Double_t xmin, Double_t xmax, Int_t npar) :TFormula(), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*F1 constructor using a pointer to real function*-*-*-*-*-*-*-* //*-* =============================================== //*-* //*-* npar is the number of free parameters used by the function //*-* //*-* This constructor creates a function of type C when invoked //*-* with the normal C++ compiler. //*-* //*-* see test program test/stress.cxx (function stress1) for an example. //*-* note the interface with an intermediate pointer. //*-* //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* fXmin = xmin; fXmax = xmax; fNpx = 100; fType = 1; fMethodCall = 0; fFunction = fcn; if (npar > 0 ) fNpar = npar; if (fNpar) { fNames = new TString[fNpar]; fParams = new Double_t[fNpar]; fParErrors = new Double_t[fNpar]; fParMin = new Double_t[fNpar]; fParMax = new Double_t[fNpar]; for (int i = 0; i < fNpar; i++) { fParams[i] = 0; fParErrors[i] = 0; fParMin[i] = 0; fParMax[i] = 0; } } else { fParErrors = 0; fParMin = 0; fParMax = 0; } fChisquare = 0; fIntegral = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fNsave = 0; fSave = 0; fParent = 0; fNpfits = 0; fNDF = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fNdim = 1; //*-*- Store formula in linked list of formula in ROOT TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name); if (f1old) delete f1old; SetName(name); gROOT->GetListOfFunctions()->Add(this); if (!gStyle) return; SetLineColor(gStyle->GetFuncColor()); SetLineWidth(gStyle->GetFuncWidth()); SetLineStyle(gStyle->GetFuncStyle()); SetFillStyle(0); } //______________________________________________________________________________ TF1::TF1(const char *name,Double_t (*fcn)(const Double_t *, const Double_t *), Double_t xmin, Double_t xmax, Int_t npar) :TFormula(), TAttLine(), TAttFill(), TAttMarker() { //*-*-*-*-*-*-*F1 constructor using a pointer to real function*-*-*-*-*-*-*-* //*-* =============================================== //*-* //*-* npar is the number of free parameters used by the function //*-* //*-* This constructor creates a function of type C when invoked //*-* with the normal C++ compiler. //*-* //*-* see test program test/stress.cxx (function stress1) for an example. //*-* note the interface with an intermediate pointer. //*-* //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* fXmin = xmin; fXmax = xmax; fNpx = 100; fType = 1; fMethodCall = 0; typedef Double_t (*Function_t) (Double_t *, Double_t *); fFunction = (Function_t)fcn; if (npar > 0 ) fNpar = npar; if (fNpar) { fNames = new TString[fNpar]; fParams = new Double_t[fNpar]; fParErrors = new Double_t[fNpar]; fParMin = new Double_t[fNpar]; fParMax = new Double_t[fNpar]; for (int i = 0; i < fNpar; i++) { fParams[i] = 0; fParErrors[i] = 0; fParMin[i] = 0; fParMax[i] = 0; } } else { fParErrors = 0; fParMin = 0; fParMax = 0; } fChisquare = 0; fIntegral = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fNsave = 0; fSave = 0; fParent = 0; fNpfits = 0; fNDF = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fNdim = 1; //*-*- Store formula in linked list of formula in ROOT TF1 *f1old = (TF1*)gROOT->GetListOfFunctions()->FindObject(name); if (f1old) delete f1old; SetName(name); gROOT->GetListOfFunctions()->Add(this); if (!gStyle) return; SetLineColor(gStyle->GetFuncColor()); SetLineWidth(gStyle->GetFuncWidth()); SetLineStyle(gStyle->GetFuncStyle()); SetFillStyle(0); } //______________________________________________________________________________ TF1& TF1::operator=(const TF1 &rhs) { if (this != &rhs) { rhs.Copy(*this); } return *this; } //______________________________________________________________________________ TF1::~TF1() { //*-*-*-*-*-*-*-*-*-*-*F1 default destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ===================== if (fParMin) delete [] fParMin; if (fParMax) delete [] fParMax; if (fParErrors) delete [] fParErrors; if (fIntegral) delete [] fIntegral; if (fAlpha) delete [] fAlpha; if (fBeta) delete [] fBeta; if (fGamma) delete [] fGamma; if (fSave) delete [] fSave; delete fHistogram; delete fMethodCall; if (fParent) { if (fParent->InheritsFrom(TH1::Class())) { ((TH1*)fParent)->GetListOfFunctions()->Remove(this); return; } //parent may be a graph, or ?? //The pad utility manager is required (a plugin) TVirtualUtilPad *util = (TVirtualUtilPad*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilPad"); if (!util) { TPluginHandler *h; if ((h = gROOT->GetPluginManager()->FindHandler("TVirtualUtilPad"))) { if (h->LoadPlugin() == -1) return; h->ExecPlugin(0); util = (TVirtualUtilPad*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilPad"); } } util->RemoveObject(fParent,this); fParent = 0; } } //______________________________________________________________________________ TF1::TF1(const TF1 &f1) : TFormula(), TAttLine(f1), TAttFill(f1), TAttMarker(f1) { fXmin = 0; fXmax = 0; fNpx = 100; fType = 0; fNpfits = 0; fNDF = 0; fNsave = 0; fChisquare = 0; fIntegral = 0; fFunction = 0; fParErrors = 0; fParMin = 0; fParMax = 0; fAlpha = 0; fBeta = 0; fGamma = 0; fParent = 0; fSave = 0; fHistogram = 0; fMinimum = -1111; fMaximum = -1111; fMethodCall = 0; SetFillStyle(0); ((TF1&)f1).Copy(*this); } //______________________________________________________________________________ void TF1::AbsValue(Bool_t flag) { // static function: set the fgAbsValue flag. // By default TF1::Integral uses the original function value to compute the integral // However, TF1::Moment, CentralMoment require to compute the integral // using the absolute value of the function. fgAbsValue = flag; } //______________________________________________________________________________ void TF1::Browse(TBrowser *b) { Draw(b ? b->GetDrawOption() : ""); gPad->Update(); } //______________________________________________________________________________ void TF1::Copy(TObject &obj) const { //*-*-*-*-*-*-*-*-*-*-*Copy this F1 to a new F1*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ======================== if (((TF1&)obj).fParMin) delete [] ((TF1&)obj).fParMin; if (((TF1&)obj).fParMax) delete [] ((TF1&)obj).fParMax; if (((TF1&)obj).fParErrors) delete [] ((TF1&)obj).fParErrors; if (((TF1&)obj).fIntegral) delete [] ((TF1&)obj).fIntegral; if (((TF1&)obj).fAlpha) delete [] ((TF1&)obj).fAlpha; if (((TF1&)obj).fBeta) delete [] ((TF1&)obj).fBeta; if (((TF1&)obj).fGamma) delete [] ((TF1&)obj).fGamma; if (((TF1&)obj).fSave) delete [] ((TF1&)obj).fSave; delete ((TF1&)obj).fHistogram; delete ((TF1&)obj).fMethodCall; TFormula::Copy(obj); TAttLine::Copy((TF1&)obj); TAttFill::Copy((TF1&)obj); TAttMarker::Copy((TF1&)obj); ((TF1&)obj).fXmin = fXmin; ((TF1&)obj).fXmax = fXmax; ((TF1&)obj).fNpx = fNpx; ((TF1&)obj).fType = fType; ((TF1&)obj).fFunction = fFunction; ((TF1&)obj).fChisquare = fChisquare; ((TF1&)obj).fNpfits = fNpfits; ((TF1&)obj).fNDF = fNDF; ((TF1&)obj).fMinimum = fMinimum; ((TF1&)obj).fMaximum = fMaximum; ((TF1&)obj).fParErrors = 0; ((TF1&)obj).fParMin = 0; ((TF1&)obj).fParMax = 0; ((TF1&)obj).fIntegral = 0; ((TF1&)obj).fAlpha = 0; ((TF1&)obj).fBeta = 0; ((TF1&)obj).fGamma = 0; ((TF1&)obj).fParent = fParent; ((TF1&)obj).fNsave = fNsave; ((TF1&)obj).fSave = 0; ((TF1&)obj).fHistogram = 0; ((TF1&)obj).fMethodCall = 0; if (fNsave) { ((TF1&)obj).fSave = new Double_t[fNsave]; for (Int_t j=0;j<fNsave;j++) ((TF1&)obj).fSave[j] = fSave[j]; } if (fNpar) { ((TF1&)obj).fParErrors = new Double_t[fNpar]; ((TF1&)obj).fParMin = new Double_t[fNpar]; ((TF1&)obj).fParMax = new Double_t[fNpar]; Int_t i; for (i=0;i<fNpar;i++) ((TF1&)obj).fParErrors[i] = fParErrors[i]; for (i=0;i<fNpar;i++) ((TF1&)obj).fParMin[i] = fParMin[i]; for (i=0;i<fNpar;i++) ((TF1&)obj).fParMax[i] = fParMax[i]; } if (fMethodCall) { TMethodCall *m = new TMethodCall(); m->InitWithPrototype(fMethodCall->GetMethodName(),fMethodCall->GetProto()); ((TF1&)obj).fMethodCall = m; } } //______________________________________________________________________________ Double_t TF1::Derivative(Double_t x, Double_t *params, Double_t eps) const { // returns the first derivative of the function at point x, // computed by Richardson's extrapolation method (use 2 derivative estimates // to compute a third, more accurate estimation) // first, derivatives with steps h and h/2 are computed by central difference formulas // D(h) = (f(x+h) - f(x-h))/2h // the final estimate D = (4*D(h/2) - D(h))/3 // "Numerical Methods for Scientists and Engineers", H.M.Antia, 2nd edition" // // if the argument params is null, the current function parameters are used, // otherwise the parameters in params are used. // // the argument eps may be specified to control the step size (precision). // the step size is taken as eps*(xmax-xmin). // the default value (0.001) should be good enough for the vast majority // of functions. Give a smaller value if your function has many changes // of the second derivative in the function range. // // Getting the error via TF1::DerivativeError // ----------------- // (total error = roundoff error + interpolation error) // the estimate of the roundoff error is taken as follows: // err = k*Sqrt(f(x)*f(x) + x*x*deriv*deriv)*Sqrt(Sum(ai)*(ai)), // where k is the double precision, ai are coefficients used in // central difference formulas // interpolation error is decreased by making the step size h smaller. // // Author: Anna Kreshuk const Double_t kC1 = 1e-15; if(eps< 1e-10 || eps > 1e-2) { Warning("Derivative","parameter esp=%g out of allowed range[1e-10,1e-2], reset to 0.001",eps); eps = 0.001; } Double_t xmin, xmax; GetRange(xmin, xmax); Double_t h = eps*(xmax-xmin); Double_t xx[1]; TF1 *func = (TF1*)this; func->InitArgs(xx, params); xx[0] = x+h; Double_t f1 = func->EvalPar(xx, params); xx[0] = x; Double_t fx = func->EvalPar(xx, params); xx[0] = x-h; Double_t f2 = func->EvalPar(xx, params); xx[0] = x+h/2; Double_t g1 = func->EvalPar(xx, params); xx[0] = x-h/2; Double_t g2 = func->EvalPar(xx, params); //compute the central differences Double_t h2 = 1/(2.*h); Double_t d0 = f1 - f2; Double_t d2 = 2*(g1 - g2); gErrorTF1 = kC1*h2*fx; //compute the error Double_t deriv = h2*(4*d2 - d0)/3.; return deriv; } //______________________________________________________________________________ Double_t TF1::Derivative2(Double_t x, Double_t *params, Double_t eps) const { // returns the first derivative of the function at point x, // computed by Richardson's extrapolation method (use 2 derivative estimates // to compute a third, more accurate estimation) // first, derivatives with steps h and h/2 are computed by central difference formulas // D(h) = (f(x+h) - 2*f(x) + f(x-h))/(h*h) // the final estimate D = (4*D(h/2) - D(h))/3 // "Numerical Methods for Scientists and Engineers", H.M.Antia, 2nd edition" // // if the argument params is null, the current function parameters are used, // otherwise the parameters in params are used. // // the argument eps may be specified to control the step size (precision). // the step size is taken as eps*(xmax-xmin). // the default value (0.001) should be good enough for the vast majority // of functions. Give a smaller value if your function has many changes // of the second derivative in the function range. // // Getting the error via TF1::DerivativeError // ----------------- // (total error = roundoff error + interpolation error) // the estimate of the roundoff error is taken as follows: // err = k*Sqrt(f(x)*f(x) + x*x*deriv*deriv)*Sqrt(Sum(ai)*(ai)), // where k is the double precision, ai are coefficients used in // central difference formulas // interpolation error is decreased by making the step size h smaller. // // Author: Anna Kreshuk const Double_t kC1 = 2*1e-15; if(eps< 1e-6 || eps > 1e-2) { Warning("Derivative2","parameter esp=%g out of allowed range[1e-6,1e-2], reset to 0.001",eps); eps = 0.001; } Double_t xmin, xmax; GetRange(xmin, xmax); Double_t h = eps*(xmax-xmin); Double_t xx[1]; TF1 *func = (TF1*)this; func->InitArgs(xx, params); xx[0] = x+h; Double_t f1 = func->EvalPar(xx, params); xx[0] = x; Double_t f2 = func->EvalPar(xx, params); xx[0] = x-h; Double_t f3 = func->EvalPar(xx, params); xx[0] = x+h/2; Double_t g1 = func->EvalPar(xx, params); xx[0] = x-h/2; Double_t g3 = func->EvalPar(xx, params); //compute the central differences Double_t hh = 1/(h*h); Double_t d0 = f3 - 2*f2 + f1; Double_t d2 = 4*g3 - 8*f2 +4*g1; gErrorTF1 = kC1*hh*f2; //compute the error Double_t deriv = hh*(4*d2 - d0)/3.; return deriv; } //______________________________________________________________________________ Double_t TF1::Derivative3(Double_t x, Double_t *params, Double_t eps) const { // returns the first derivative of the function at point x, // computed by Richardson's extrapolation method (use 2 derivative estimates // to compute a third, more accurate estimation) // first, derivatives with steps h and h/2 are computed by central difference formulas // D(h) = (f(x+2h) - 2*f(x+h) + 2*f(x-h) - f(x-2h))/(2*h*h*h) // the final estimate D = (4*D(h/2) - D(h))/3 // "Numerical Methods for Scientists and Engineers", H.M.Antia, 2nd edition" // // if the argument params is null, the current function parameters are used, // otherwise the parameters in params are used. // // the argument eps may be specified to control the step size (precision). // the step size is taken as eps*(xmax-xmin). // the default value (0.001) should be good enough for the vast majority // of functions. Give a smaller value if your function has many changes // of the second derivative in the function range. // // Getting the error via TF1::DerivativeError // ----------------- // (total error = roundoff error + interpolation error) // the estimate of the roundoff error is taken as follows: // err = k*Sqrt(f(x)*f(x) + x*x*deriv*deriv)*Sqrt(Sum(ai)*(ai)), // where k is the double precision, ai are coefficients used in // central difference formulas // interpolation error is decreased by making the step size h smaller. // // Author: Anna Kreshuk //const Double_t C1 = (1e-16)*TMath::Sqrt(5./2.)*TMath::Sqrt(16*64 + 1.)/3; const Double_t kC1 = 1e-15; if(eps< 1e-4 || eps > 1e-2) { Warning("Derivative3","parameter esp=%g out of allowed range[1e-4,1e-2], reset to 0.001",eps); eps = 0.001; } Double_t xmin, xmax; GetRange(xmin, xmax); Double_t h = eps*(xmax-xmin); Double_t xx[1]; TF1 *func = (TF1*)this; func->InitArgs(xx, params); xx[0] = x+2*h; Double_t f1 = func->EvalPar(xx, params); xx[0] = x+h; Double_t f2 = func->EvalPar(xx, params); xx[0] = x-h; Double_t f3 = func->EvalPar(xx, params); xx[0] = x-2*h; Double_t f4 = func->EvalPar(xx, params); xx[0] = x; Double_t fx = func->EvalPar(xx, params); xx[0] = x+h/2; Double_t g2 = func->EvalPar(xx, params); xx[0] = x-h/2; Double_t g3 = func->EvalPar(xx, params); //compute the central differences Double_t hhh = 1/(h*h*h); Double_t d0 = 0.5*f1 - f2 +f3 - 0.5*f4; Double_t d2 = 4*f2 - 8*g2 +8*g3 - 4*f3; gErrorTF1 = kC1*hhh*fx; //compute the error Double_t deriv = hhh*(4*d2 - d0)/3.; return deriv; } //______________________________________________________________________________ Double_t TF1::DerivativeError() { //static function returning the error of the last call to the Derivative functions return gErrorTF1; } //______________________________________________________________________________ Int_t TF1::DistancetoPrimitive(Int_t px, Int_t py) { //*-*-*-*-*-*-*-*-*-*-*Compute distance from point px,py to a function*-*-*-*-* //*-* =============================================== //*-* Compute the closest distance of approach from point px,py to this function. //*-* The distance is computed in pixels units. //*-* //*-* Note that px is called with a negative value when the TF1 is in //*-* TGraph or TH1 list of functions. In this case there is no point //*-* looking at the histogram axis. //*-* //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* if (!fHistogram) return 9999; Int_t distance = 9999; if (px >= 0) { distance = fHistogram->DistancetoPrimitive(px,py); if (distance <= 1) return distance; } else { px = -px; } Double_t xx[1]; Double_t x = gPad->AbsPixeltoX(px); xx[0] = gPad->PadtoX(x); if (xx[0] < fXmin || xx[0] > fXmax) return distance; Double_t fval = Eval(xx[0]); Double_t y = gPad->YtoPad(fval); Int_t pybin = gPad->YtoAbsPixel(y); return TMath::Abs(py - pybin); } //______________________________________________________________________________ void TF1::Draw(Option_t *option) { //*-*-*-*-*-*-*-*-*-*-*Draw this function with its current attributes*-*-*-*-* //*-* ============================================== //*-* //*-* Possible option values are: //*-* "SAME" superimpose on top of existing picture //*-* "L" connect all computed points with a straight line //*-* "C" connect all computed points with a smooth curve. //*-* "FC" draw a fill area below a smooth curve //*-* //*-* Note that the default value is "L". Therefore to draw on top //*-* of an existing picture, specify option "LSAME" //*-* //*-* NB. You must use DrawCopy if you want to draw several times the same //*-* function in the current canvas. //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* TString opt = option; opt.ToLower(); if (gPad && !opt.Contains("same")) gPad->Clear(); AppendPad(option); } //______________________________________________________________________________ TF1 *TF1::DrawCopy(Option_t *option) const { //*-*-*-*-*-*-*-*Draw a copy of this function with its current attributes*-*-* //*-* ======================================================== //*-* //*-* This function MUST be used instead of Draw when you want to draw //*-* the same function with different parameters settings in the same canvas. //*-* //*-* Possible option values are: //*-* "SAME" superimpose on top of existing picture //*-* "L" connect all computed points with a straight line //*-* "C" connect all computed points with a smooth curve. //*-* "FC" draw a fill area below a smooth curve //*-* //*-* Note that the default value is "L". Therefore to draw on top //*-* of an existing picture, specify option "LSAME" //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* TF1 *newf1 = new TF1(); Copy(*newf1); newf1->AppendPad(option); newf1->SetBit(kCanDelete); return newf1; } //______________________________________________________________________________ void TF1::DrawDerivative(Option_t *option) { // Draw derivative of this function // // An intermediate TGraph object is built and drawn with option. // // The resulting graph will be drawn into the current pad. // If this function is used via the context menu, it recommended // to create a new canvas/pad before invoking this function. TVirtualPad *pad = gROOT->GetSelectedPad(); TVirtualPad *padsav = gPad; if (pad) pad->cd(); char cmd[512]; sprintf(cmd,"{TGraph *R__%s_Derivative = new TGraph((TF1*)0x%lx,\"d\");R__%s_Derivative->Draw(\"%s\");}",GetName(),(Long_t)this,GetName(),option); gROOT->ProcessLine(cmd); if (padsav) padsav->cd(); } //______________________________________________________________________________ void TF1::DrawIntegral(Option_t *option) { // Draw integral of this function // // An intermediate TGraph object is built and drawn with option. // // The resulting graph will be drawn into the current pad. // If this function is used via the context menu, it recommended // to create a new canvas/pad before invoking this function. TVirtualPad *pad = gROOT->GetSelectedPad(); TVirtualPad *padsav = gPad; if (pad) pad->cd(); char cmd[512]; sprintf(cmd,"{TGraph *R__%s_Integral = new TGraph((TF1*)0x%lx,\"i\");R__%s_Integral->Draw(\"%s\");}",GetName(),(Long_t)this,GetName(),option); gROOT->ProcessLine(cmd); if (padsav) padsav->cd(); } //______________________________________________________________________________ void TF1::DrawF1(const char *formula, Double_t xmin, Double_t xmax, Option_t *option) { //*-*-*-*-*-*-*-*-*-*Draw formula between xmin and xmax*-*-*-*-*-*-*-*-*-*-*-* //*-* ================================== //*-* if (Compile(formula)) return; SetRange(xmin, xmax); Draw(option); } //______________________________________________________________________________ void TF1::DrawPanel() { //*-*-*-*-*-*-*Display a panel with all function drawing options*-*-*-*-*-* //*-* ================================================= //*-* //*-* See class TDrawPanelHist for example //The pad utility manager is required (a plugin) TVirtualUtilPad *util = (TVirtualUtilPad*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilPad"); if (!util) { TPluginHandler *h; if ((h = gROOT->GetPluginManager()->FindHandler("TVirtualUtilPad"))) { if (h->LoadPlugin() == -1) return; h->ExecPlugin(0); util = (TVirtualUtilPad*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilPad"); } } util->DrawPanel(gPad,this); } //______________________________________________________________________________ Double_t TF1::Eval(Double_t x, Double_t y, Double_t z, Double_t t) const { //*-*-*-*-*-*-*-*-*-*-*Evaluate this formula*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ===================== //*-* //*-* Computes the value of this function (general case for a 3-d function) //*-* at point x,y,z. //*-* For a 1-d function give y=0 and z=0 //*-* The current value of variables x,y,z is passed through x, y and z. //*-* The parameters used will be the ones in the array params if params is given //*-* otherwise parameters will be taken from the stored data members fParams //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* Double_t xx[4]; xx[0] = x; xx[1] = y; xx[2] = z; xx[3] = t; ((TF1*)this)->InitArgs(xx,fParams); return ((TF1*)this)->EvalPar(xx,fParams); } //______________________________________________________________________________ Double_t TF1::EvalPar(const Double_t *x, const Double_t *params) { //*-*-*-*-*-*Evaluate function with given coordinates and parameters*-*-*-*-*-* //*-* ======================================================= //*-* // Compute the value of this function at point defined by array x // and current values of parameters in array params. // If argument params is omitted or equal 0, the internal values // of parameters (array fParams) will be used instead. // For a 1-D function only x[0] must be given. // In case of a multi-dimemsional function, the arrays x must be // filled with the corresponding number of dimensions. // // WARNING. In case of an interpreted function (fType=2), it is the // user's responsability to initialize the parameters via InitArgs // before calling this function. // InitArgs should be called at least once to specify the addresses // of the arguments x and params. // InitArgs should be called everytime these addresses change. // fgCurrent = this; if (fType == 0) return TFormula::EvalPar(x,params); Double_t result = 0; if (fType == 1) { if (fFunction) { if (params) result = (*fFunction)((Double_t*)x,(Double_t*)params); else result = (*fFunction)((Double_t*)x,fParams); }else result = GetSave(x); } if (fType == 2) { if (fMethodCall) fMethodCall->Execute(result); else result = GetSave(x); } return result; } //______________________________________________________________________________ void TF1::ExecuteEvent(Int_t event, Int_t px, Int_t py) { //*-*-*-*-*-*-*-*-*-*-*Execute action corresponding to one event*-*-*-* //*-* ========================================= //*-* This member function is called when a F1 is clicked with the locator //*-* //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* fHistogram->ExecuteEvent(event,px,py); if (!gPad->GetView()) { if (event == kMouseMotion) gPad->SetCursor(kHand); } } //______________________________________________________________________________ void TF1::FixParameter(Int_t ipar, Double_t value) { // Fix the value of a parameter // The specified value will be used in a fit operation if (ipar < 0 || ipar > fNpar-1) return; SetParameter(ipar,value); if (value != 0) SetParLimits(ipar,value,value); else SetParLimits(ipar,1,1); } //______________________________________________________________________________ TF1 *TF1::GetCurrent() { // static function returning the current function being processed return fgCurrent; } //______________________________________________________________________________ TH1 *TF1::GetHistogram() const { // return a pointer to the histogram used to vusualize the function if (fHistogram) return fHistogram; // may be function has not yet be painted. force a pad update //gPad->Modified(); //gPad->Update(); ((TF1*)this)->Paint(); return fHistogram; } //______________________________________________________________________________ Double_t TF1::GetMaximum(Double_t xmin, Double_t xmax) const { // return the maximum value of the function // Method: // First, the grid search is used to bracket the maximum // with the step size = (xmax-xmin)/fNpx. This way, the step size // can be controlled via the SetNpx() function. If the function is // unimodal or if its extrema are far apart, setting the fNpx to // a small value speeds the algorithm up many times. // Then, Brent's method is applied on the bracketed interval if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;} Double_t x; Int_t niter=0; x = MinimStep(3, xmin, xmax, 0); Bool_t ok = kTRUE; x = MinimBrent(3, xmin, xmax, x, 0, ok); while (!ok){ if (niter>10){ Error("GetMaximum", "maximum search didn't converge"); break; } x=MinimStep(3, xmin, xmax,0); x = MinimBrent(3, xmin, xmax, x, 0, ok); niter++; } return x; } //______________________________________________________________________________ Double_t TF1::GetMaximumX(Double_t xmin, Double_t xmax) const { // return the X value corresponding to the maximum value of the function // Method: // First, the grid search is used to bracket the maximum // with the step size = (xmax-xmin)/fNpx. This way, the step size // can be controlled via the SetNpx() function. If the function is // unimodal or if its extrema are far apart, setting the fNpx to // a small value speeds the algorithm up many times. // Then, Brent's method is applied on the bracketed interval if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;} Double_t x; Int_t niter=0; x = MinimStep(2, xmin, xmax, 0); Bool_t ok = kTRUE; x = MinimBrent(2, xmin, xmax, x, 0, ok); while (!ok){ if (niter>10){ Error("GetMaximumX", "maximum search didn't converge"); break; } x=MinimStep(2, xmin, xmax, 0); x = MinimBrent(2, xmin, xmax, x, 0, ok); niter++; } return x; } //______________________________________________________________________________ Double_t TF1::GetMinimum(Double_t xmin, Double_t xmax) const { // Returns the minimum value of the function on the (xmin, xmax) interval // Method: // First, the grid search is used to bracket the maximum // with the step size = (xmax-xmin)/fNpx. This way, the step size // can be controlled via the SetNpx() function. If the function is // unimodal or if its extrema are far apart, setting the fNpx to // a small value speeds the algorithm up many times. // Then, Brent's method is applied on the bracketed interval if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;} Double_t x; Int_t niter=0; x = MinimStep(1, xmin, xmax, 0); Bool_t ok = kTRUE; x = MinimBrent(1, xmin, xmax, x, 0, ok); while (!ok){ if (niter>10){ Error("GetMinimum", "minimum search didn't converge"); break; } x=MinimStep(1, xmin, xmax,0); x = MinimBrent(1, xmin, xmax, x, 0, ok); niter++; } return x; } //______________________________________________________________________________ Double_t TF1::GetMinimumX(Double_t xmin, Double_t xmax) const { // Returns the X value corresponding to the minimum value of the function on the // (xmin, xmax) interval // Method: // First, the grid search is used to bracket the maximum // with the step size = (xmax-xmin)/fNpx. This way, the step size // can be controlled via the SetNpx() function. If the function is // unimodal or if its extrema are far apart, setting the fNpx to // a small value speeds the algorithm up many times. // Then, Brent's method is applied on the bracketed interval if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;} Int_t niter=0; Double_t x; x = MinimStep(0, xmin, xmax, 0); Bool_t ok = kTRUE; x = MinimBrent(0, xmin, xmax, x, 0, ok); while (!ok){ if (niter>10){ Error("GetMinimumX", "minimum search didn't converge"); break; } x=MinimStep(0, xmin, xmax,0); x = MinimBrent(0, xmin, xmax, x, 0, ok); niter++; } return x; } //______________________________________________________________________________ Double_t TF1::GetX(Double_t fy, Double_t xmin, Double_t xmax) const { // Returns the X value corresponding to the function value fy for (xmin<x<xmax). // Method: // First, the grid search is used to bracket the maximum // with the step size = (xmax-xmin)/fNpx. This way, the step size // can be controlled via the SetNpx() function. If the function is // unimodal or if its extrema are far apart, setting the fNpx to // a small value speeds the algorithm up many times. // Then, Brent's method is applied on the bracketed interval if (xmin >= xmax) {xmin = fXmin; xmax = fXmax;} Int_t niter=0; Double_t x; x = MinimStep(4, xmin, xmax, fy); Bool_t ok = kTRUE; x = MinimBrent(4, xmin, xmax, x, fy, ok); while (!ok){ if (niter>10){ Error("GetX", "Search didn't converge"); break; } x=MinimStep(4, xmin, xmax, fy); x = MinimBrent(4, xmin, xmax, x, fy, ok); niter++; } return x; } //______________________________________________________________________________ Double_t TF1::MinimStep(Int_t type, Double_t &xmin, Double_t &xmax, Double_t fy) const { // Grid search implementation, used to bracket the minimum and later // use Brent's method with the bracketed interval // The step of the search is set to (xmax-xmin)/fNpx // type: 0-returns MinimumX // 1-returns Minimum // 2-returns MaximumX // 3-returns Maximum // 4-returns X corresponding to fy Double_t x,y, dx; dx = (xmax-xmin)/(fNpx-1); Double_t xxmin = xmin; Double_t yymin; if (type < 2) yymin = Eval(xmin); else if (type < 4) yymin = -Eval(xmin); else yymin = TMath::Abs(Eval(xmin)-fy); for (Int_t i=1; i<=fNpx-1; i++) { x = xmin + i*dx; if (type < 2) y = Eval(x); else if (type < 4) y = -Eval(x); else y = TMath::Abs(Eval(x)-fy); if (y < yymin) {xxmin = x; yymin = y;} } xmin = TMath::Max(xmin,xxmin-dx); xmax = TMath::Min(xmax,xxmin+dx); return TMath::Min(xxmin, xmax); } //______________________________________________________________________________ Double_t TF1::MinimBrent(Int_t type, Double_t &xmin, Double_t &xmax, Double_t xmiddle, Double_t fy, Bool_t &ok) const { //Finds a minimum of a function, if the function is unimodal between xmin and xmax //This method uses a combination of golden section search and parabolic interpolation //Details about convergence and properties of this algorithm can be //found in the book by R.P.Brent "Algorithms for Minimization Without Derivatives" //or in the "Numerical Recipes", chapter 10.2 // //type: 0-returns MinimumX // 1-returns Minimum // 2-returns MaximumX // 3-returns Maximum // 4-returns X corresponding to fy //if ok=true the method has converged Double_t eps = 1e-10; Double_t t = 1e-8; Int_t itermax = 100; Double_t c = (3.-TMath::Sqrt(5.))/2.; //comes from golden section Double_t u, v, w, x, fv, fu, fw, fx, e, p, q, r, t2, d=0, m, tol; v = w = x = xmiddle; e=0; Double_t a=xmin; Double_t b=xmax; if (type < 2) fv = fw = fx = Eval(x); else if (type < 4) fv = fw = fx = -Eval(x); else fv = fw = fx = TMath::Abs(Eval(x)-fy); for (Int_t i=0; i<itermax; i++){ m=0.5*(a + b); tol = eps*(TMath::Abs(x))+t; t2 = 2*tol; if (TMath::Abs(x-m) <= (t2-0.5*(b-a))) { //converged, return x ok=kTRUE; if (type==1) return fx; else if (type==3) return -fx; else return x; } if (TMath::Abs(e)>tol){ //fit parabola r = (x-w)*(fx-fv); q = (x-v)*(fx-fw); p = (x-v)*q - (x-w)*r; q = 2*(q-r); if (q>0) p=-p; else q=-q; r=e; e=d; if (TMath::Abs(p) < TMath::Abs(0.5*q*r) || p < q*(a-x) || p < q*(b-x)) { //a parabolic interpolation step d = p/q; u = x+d; if (u-a < t2 || b-u < t2) d=TMath::Sign(tol, m-x); } else { e=(x>=m ? a-x : b-x); d = c*e; } } else { e=(x>=m ? a-x : b-x); d = c*e; } u = (TMath::Abs(d)>=tol ? x+d : x+TMath::Sign(tol, d)); if (type < 2) fu = Eval(u); else if (type < 4) fu = -Eval(u); else fu = TMath::Abs(Eval(u)-fy); //update a, b, v, w and x if (fu<=fx){ if (u<x) b=x; else a=x; v=w; fv=fw; w=x; fw=fx; x=u; fx=fu; } else { if (u<x) a=u; else b=u; if (fu<=fw || w==x){ v=w; fv=fw; w=u; fw=fu; } else if (fu<=fv || v==x || v==w){ v=u; fv=fu; } } } //didn't converge ok = kFALSE; xmin = a; xmax = b; return x; } //______________________________________________________________________________ Int_t TF1::GetNDF() const { // return the number of degrees of freedom in the fit // the fNDF parameter has been previously computed during a fit. // The number of degrees of freedom corresponds to the number of points // used in the fit minus the number of free parameters. if (fNDF == 0) return fNpfits-fNpar; return fNDF; } //______________________________________________________________________________ Int_t TF1::GetNumberFreeParameters() const { // return the number of free parameters Int_t nfree = fNpar; Double_t al,bl; for (Int_t i=0;i<fNpar;i++) { ((TF1*)this)->GetParLimits(i,al,bl); if (al*bl != 0 && al >= bl) nfree--; } return nfree; } //______________________________________________________________________________ char *TF1::GetObjectInfo(Int_t px, Int_t /* py */) const { // Redefines TObject::GetObjectInfo. // Displays the function info (x, function value // corresponding to cursor position px,py // static char info[64]; Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px)); sprintf(info,"(x=%g, f=%g)",x,((TF1*)this)->Eval(x)); return info; } //______________________________________________________________________________ Double_t TF1::GetParError(Int_t ipar) const { //return value of parameter number ipar if (ipar < 0 || ipar > fNpar-1) return 0; return fParErrors[ipar]; } //______________________________________________________________________________ void TF1::GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const { //*-*-*-*-*-*Return limits for parameter ipar*-*-*-* //*-* ================================ parmin = 0; parmax = 0; if (ipar < 0 || ipar > fNpar-1) return; if (fParMin) parmin = fParMin[ipar]; if (fParMax) parmax = fParMax[ipar]; } //______________________________________________________________________________ Double_t TF1::GetProb() const { // return the fit probability if (fNDF <= 0) return 0; return TMath::Prob(fChisquare,fNDF); } //______________________________________________________________________________ Int_t TF1::GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum) { // Compute Quantiles for density distribution of this function // Quantile x_q of a probability distribution Function F is defined as // // F(x_q) = Integral_{xmin}^(x_q) f dx = q with 0 <= q <= 1. // // For instance the median x_0.5 of a distribution is defined as that value // of the random variable for which the distribution function equals 0.5: // // F(x_0.5) = Probability(x < x_0.5) = 0.5 // // code from Eddy Offermann, Renaissance // // input parameters // - this TF1 function // - nprobSum maximum size of array q and size of array probSum // - probSum array of positions where quantiles will be computed. // It is assumed to contain at least nprobSum values. // output // - return value nq (<=nprobSum) with the number of quantiles computed // - array q filled with nq quantiles // // Getting quantiles from two histograms and storing results in a TGraph, // a so-called QQ-plot // // TGraph *gr = new TGraph(nprob); // f1->GetQuantiles(nprob,gr->GetX()); // f2->GetQuantiles(nprob,gr->GetY()); // gr->Draw("alp"); const Int_t npx = TMath::Min(250,TMath::Max(50,2*nprobSum)); const Double_t xMin = GetXmin(); const Double_t xMax = GetXmax(); const Double_t dx = (xMax-xMin)/npx; TArrayD integral(npx+1); TArrayD alpha(npx); TArrayD beta(npx); TArrayD gamma(npx); integral[0] = 0; Int_t intNegative = 0; Int_t i; for (i = 0; i < npx; i++) { const Double_t *params = 0; Double_t integ = Integral(Double_t(xMin+i*dx),Double_t(xMin+i*dx+dx),params); if (integ < 0) {intNegative++; integ = -integ;} integral[i+1] = integral[i] + integ; } if (intNegative > 0) Warning("GetQuantiles","function:%s has %d negative values: abs assumed", GetName(),intNegative); if (integral[npx] == 0) { Error("GetQuantiles","Integral of function is zero"); return 0; } const Double_t total = integral[npx]; for (i = 1; i <= npx; i++) integral[i] /= total; //the integral r for each bin is approximated by a parabola // x = alpha + beta*r +gamma*r**2 // compute the coefficients alpha, beta, gamma for each bin for (i = 0; i < npx; i++) { const Double_t x0 = xMin+dx*i; const Double_t r2 = integral[i+1]-integral[i]; const Double_t r1 = Integral(x0,x0+0.5*dx)/total; gamma[i] = (2*r2-4*r1)/(dx*dx); beta[i] = r2/dx-gamma[i]*dx; alpha[i] = x0; gamma[i] *= 2; } // Be careful because of finite precision in the integral; Use the fact that the integral // is monotone increasing for (i = 0; i < nprobSum; i++) { const Double_t r = probSum[i]; Int_t bin = TMath::Max(TMath::BinarySearch(npx+1,integral.GetArray(),r)-1,(Long64_t)0); while (bin < npx-1 && integral[bin+1] == r) { if (integral[bin+2] == r) bin++; else break; } const Double_t rr = r-integral[bin]; if (rr != 0.0) { Double_t xx; if (gamma[bin] && beta[bin]*beta[bin]+2*gamma[bin]*rr >= 0.0) xx = (-beta[bin]+TMath::Sqrt(beta[bin]*beta[bin]+2*gamma[bin]*rr))/gamma[bin]; else xx = rr/beta[bin]; q[i] = alpha[bin]+xx; } else { q[i] = alpha[bin]; if (integral[bin+1] == r) q[i] += dx; } } return nprobSum; } //______________________________________________________________________________ Double_t TF1::GetRandom() { // Return a random number following this function shape //*-* //*-* The distribution contained in the function fname (TF1) is integrated //*-* over the channel contents. //*-* It is normalized to 1. //*-* For each bin the integral is approximated by a parabola. //*-* The parabola coefficients are stored as non persistent data members //*-* Getting one random number implies: //*-* - Generating a random number between 0 and 1 (say r1) //*-* - Look in which bin in the normalized integral r1 corresponds to //*-* - Evaluate the parabolic curve in the selected bin to find //*-* the corresponding X value. //*-* The parabolic approximation is very good as soon as the number //*-* of bins is greater than 50. //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-* // Check if integral array must be build if (fIntegral == 0) { Double_t dx = (fXmax-fXmin)/fNpx; fIntegral = new Double_t[fNpx+1]; fAlpha = new Double_t[fNpx]; fBeta = new Double_t[fNpx]; fGamma = new Double_t[fNpx]; fIntegral[0] = 0; Double_t integ; Int_t intNegative = 0; Int_t i; for (i=0;i<fNpx;i++) { integ = Integral(Double_t(fXmin+i*dx), Double_t(fXmin+i*dx+dx)); if (integ < 0) {intNegative++; integ = -integ;} fIntegral[i+1] = fIntegral[i] + integ; } if (intNegative > 0) { Warning("GetRandom","function:%s has %d negative values: abs assumed",GetName(),intNegative); } if (fIntegral[fNpx] == 0) { Error("GetRandom","Integral of function is zero"); return 0; } Double_t total = fIntegral[fNpx]; for (i=1;i<=fNpx;i++) { // normalize integral to 1 fIntegral[i] /= total; } //the integral r for each bin is approximated by a parabola // x = alpha + beta*r +gamma*r**2 // compute the coefficients alpha, beta, gamma for each bin Double_t x0,r1,r2,r3; for (i=0;i<fNpx;i++) { x0 = fXmin+i*dx; r2 = fIntegral[i+1] - fIntegral[i]; r1 = Integral(x0,x0+0.5*dx)/total; r3 = 2*r2 - 4*r1; if (TMath::Abs(r3) > 1e-8) fGamma[i] = r3/(dx*dx); else fGamma[i] = 0; fBeta[i] = r2/dx - fGamma[i]*dx; fAlpha[i] = x0; fGamma[i] *= 2; } } // return random number Double_t r = gRandom->Rndm(); Int_t bin = TMath::BinarySearch(fNpx,fIntegral,r); Double_t rr = r - fIntegral[bin]; Double_t xx; if(fGamma[bin] != 0) xx = (-fBeta[bin] + TMath::Sqrt(fBeta[bin]*fBeta[bin]+2*fGamma[bin]*rr))/fGamma[bin]; else xx = rr/fBeta[bin]; Double_t x = fAlpha[bin] + xx; return x; } //______________________________________________________________________________ Double_t TF1::GetRandom(Double_t xmin, Double_t xmax) { // Return a random number following this function shape in [xmin,xmax] //*-* //*-* The distribution contained in the function fname (TF1) is integrated //*-* over the channel contents. //*-* It is normalized to 1. //*-* For each bin the integral is approximated by a parabola. //*-* The parabola coefficients are stored as non persistent data members //*-* Getting one random number implies: //*-* - Generating a random number between 0 and 1 (say r1) //*-* - Look in which bin in the normalized integral r1 corresponds to //*-* - Evaluate the parabolic curve in the selected bin to find //*-* the corresponding X value. //*-* The parabolic approximation is very good as soon as the number //*-* of bins is greater than 50. //*-* //*-* IMPORTANT NOTE //*-* The integral of the function is computed at fNpx points. If the function //*-* has sharp peaks, you should increase the number of points (SetNpx) //*-* such that the peak is correctly tabulated at several points. //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-* // Check if integral array must be build if (fIntegral == 0) { Double_t dx = (fXmax-fXmin)/fNpx; fIntegral = new Double_t[fNpx+1]; fAlpha = new Double_t[fNpx]; fBeta = new Double_t[fNpx]; fGamma = new Double_t[fNpx]; fIntegral[0] = 0; Double_t integ; Int_t intNegative = 0; Int_t i; for (i=0;i<fNpx;i++) { integ = Integral(Double_t(fXmin+i*dx), Double_t(fXmin+i*dx+dx)); if (integ < 0) {intNegative++; integ = -integ;} fIntegral[i+1] = fIntegral[i] + integ; } if (intNegative > 0) { Warning("GetRandom","function:%s has %d negative values: abs assumed",GetName(),intNegative); } if (fIntegral[fNpx] == 0) { Error("GetRandom","Integral of function is zero"); return 0; } Double_t total = fIntegral[fNpx]; for (i=1;i<=fNpx;i++) { // normalize integral to 1 fIntegral[i] /= total; } //the integral r for each bin is approximated by a parabola // x = alpha + beta*r +gamma*r**2 // compute the coefficients alpha, beta, gamma for each bin Double_t x0,r1,r2,r3; for (i=0;i<fNpx;i++) { x0 = fXmin+i*dx; r2 = fIntegral[i+1] - fIntegral[i]; r1 = Integral(x0,x0+0.5*dx)/total; r3 = 2*r2 - 4*r1; if (TMath::Abs(r3) > 1e-8) fGamma[i] = r3/(dx*dx); else fGamma[i] = 0; fBeta[i] = r2/dx - fGamma[i]*dx; fAlpha[i] = x0; fGamma[i] *= 2; } } // return random number Double_t dx = (fXmax-fXmin)/fNpx; Int_t nbinmin = (Int_t)((xmin-fXmin)/dx); Int_t nbinmax = (Int_t)((xmax-fXmin)/dx)+2; if(nbinmax>fNpx) nbinmax=fNpx; Double_t pmin=fIntegral[nbinmin]; Double_t pmax=fIntegral[nbinmax]; Double_t r,x,xx,rr; do { r = gRandom->Uniform(pmin,pmax); Int_t bin = TMath::BinarySearch(fNpx,fIntegral,r); rr = r - fIntegral[bin]; if(fGamma[bin] != 0) xx = (-fBeta[bin] + TMath::Sqrt(fBeta[bin]*fBeta[bin]+2*fGamma[bin]*rr))/fGamma[bin]; else xx = rr/fBeta[bin]; x = fAlpha[bin] + xx; } while(x<xmin || x>xmax); return x; } //______________________________________________________________________________ void TF1::GetRange(Double_t &xmin, Double_t &xmax) const { //*-*-*-*-*-*-*-*-*-*-*Return range of a 1-D function*-*-*-*-*-*-*-*-*-*-*-* //*-* ============================== xmin = fXmin; xmax = fXmax; } //______________________________________________________________________________ void TF1::GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const { //*-*-*-*-*-*-*-*-*-*-*Return range of a 2-D function*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ============================== xmin = fXmin; xmax = fXmax; ymin = 0; ymax = 0; } //______________________________________________________________________________ void TF1::GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const { //*-*-*-*-*-*-*-*-*-*-*Return range of function*-*-*-*-*-*-*-*-*-*-*-*-*-*-* //*-* ======================== xmin = fXmin; xmax = fXmax; ymin = 0; ymax = 0; zmin = 0; zmax = 0; } //______________________________________________________________________________ Double_t TF1::GetSave(const Double_t *xx) { // Get value corresponding to X in array of fSave values if (fNsave <= 0) return 0; if (fSave == 0) return 0; Int_t np = fNsave - 3; Double_t xmin = Double_t(fSave[np+1]); Double_t xmax = Double_t(fSave[np+2]); Double_t x = Double_t(xx[0]); Double_t dx = (xmax-xmin)/np; if (x < xmin || x > xmax) return 0; if (dx <= 0) return 0; Int_t bin = Int_t((x-xmin)/dx); Double_t xlow = xmin + bin*dx; Double_t xup = xlow + dx; Double_t ylow = fSave[bin]; Double_t yup = fSave[bin+1]; Double_t y = ((xup*ylow-xlow*yup) + x*(yup-ylow))/dx; return y; } //______________________________________________________________________________ TAxis *TF1::GetXaxis() const { // Get x axis of the function. //if (!gPad) return 0; TH1 *h = GetHistogram(); if (!h) return 0; return h->GetXaxis(); } //______________________________________________________________________________ TAxis *TF1::GetYaxis() const { // Get y axis of the function. //if (!gPad) return 0; TH1 *h = GetHistogram(); if (!h) return 0; return h->GetYaxis(); } //______________________________________________________________________________ TAxis *TF1::GetZaxis() const { // Get z axis of the function. (In case this object is a TF2 or TF3) //if (!gPad) return 0; TH1 *h = GetHistogram(); if (!h) return 0; return h->GetZaxis(); } //______________________________________________________________________________ void TF1::InitArgs(const Double_t *x, const Double_t *params) { //*-*-*-*-*-*-*-*-*-*-*Initialize parameters addresses*-*-*-*-*-*-*-*-*-*-*-* //*-* =============================== if (fMethodCall) { Long_t args[2]; args[0] = (Long_t)x; if (params) args[1] = (Long_t)params; else args[1] = (Long_t)fParams; fMethodCall->SetParamPtrs(args); } } //______________________________________________________________________________ void TF1::InitStandardFunctions() { // Create the basic function objects TF1 *f1; if (!gROOT->GetListOfFunctions()->FindObject("gaus")) { f1 = new TF1("gaus","gaus",-1,1); f1->SetParameters(1,0,1); f1 = new TF1("gausn","gausn",-1,1); f1->SetParameters(1,0,1); f1 = new TF1("landau","landau",-1,1); f1->SetParameters(1,0,1); f1 = new TF1("landaun","landaun",-1,1); f1->SetParameters(1,0,1); f1 = new TF1("expo","expo",-1,1); f1->SetParameters(1,1); for (Int_t i=0;i<10;i++) { f1 = new TF1(Form("pol%d",i),Form("pol%d",i),-1,1); f1->SetParameters(1,1,1,1,1,1,1,1,1,1); } } } //______________________________________________________________________________ Double_t TF1::Integral(Double_t a, Double_t b, const Double_t *params, Double_t epsilon) { //*-*-*-*-*-*-*-*-*Return Integral of function between a and b*-*-*-*-*-*-*-* // // based on original CERNLIB routine DGAUSS by Sigfried Kolbig // converted to C++ by Rene Brun // // /*
This function computes, to an attempted specified accuracy, the value of the integral
Usage:
In any arithmetic expression, this function has the approximate value of the integral I.
Method:
For any interval [a,b] we define and to be the 8-point and 16-point Gaussian quadrature approximations to
and define
Then,
where, starting with and finishing with , the subdivision points are given by
with equal to the first member of the sequence for which . If, at any stage in the process of subdivision, the ratio
is so small that 1+0.005q is indistinguishable from 1 to machine accuracy, an error exit occurs with the function value set equal to zero.
Accuracy:
Unless there is severe cancellation of positive and negative values of f(x) over the interval [A,B], the argument EPS may be considered as specifying a bound on the relative error of I in the case |I|>1, and a bound on the absolute error in the case |I|<1. More precisely, if k is the number of sub-intervals contributing to the approximation (see Method), and if
then the relation
will nearly always be true, provided the routine terminates without printing an error message. For functions f having no singularities in the closed interval [A,B] the accuracy will usually be much higher than this.
Error handling:
The requested accuracy cannot be obtained (see Method). The function value is set equal to zero.
Notes:
Values of the function f(x) at the interval end-points
A and B are not required. The subprogram may therefore
be used when these values are undefined.