Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TF1 Class Reference

1-Dim function class

TF1: 1-Dim function class

A TF1 object is a 1-Dim function defined between a lower and upper limit. The function may be a simple function based on a TFormula expression or a precompiled user function. The function may have associated parameters. TF1 graphics function is via the TH1 and TGraph drawing functions.

The following types of functions can be created:

  1. Expression using variable x and no parameters
  2. Expression using variable x with parameters
  3. Lambda Expression with variable x and parameters
  4. A general C function with parameters
  5. A general C++ function object (functor) with parameters
  6. A member function with parameters of a general C++ class

1 - Expression using variable x and no parameters

Case 1: inline expression using standard C++ functions/operators

{
auto fa1 = new TF1("fa1","sin(x)/x",0,10);
fa1->Draw();
}
TF1()
TF1 default constructor.
Definition TF1.cxx:489

Case 2: inline expression using a ROOT function (e.g. from TMath) without parameters

{
auto fa2 = new TF1("fa2","TMath::DiLog(x)",0,10);
fa2->Draw();
}

Case 3: inline expression using a user defined CLING function by name

Double_t myFunc(double x) { return x+sin(x); }
....
auto fa3 = new TF1("fa3","myFunc(x)",-3,5);
fa3->Draw();
Double_t myFunc(Double_t x)
Definition ROOTR.C:4
Double_t x[n]
Definition legend1.C:17

2 - Expression using variable x with parameters

Case 1: inline expression using standard C++ functions/operators

Example a:

auto fa = new TF1("fa","[0]*x*sin([1]*x)",-3,3);

This creates a function of variable x with 2 parameters. The parameters must be initialized via:

fa->SetParameter(0,value_first_parameter);
fa->SetParameter(1,value_second_parameter);

Parameters may be given a name:

fa->SetParName(0,"Constant");

Example b:

auto fb = new TF1("fb","gaus(0)*expo(3)",0,10);

gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2) and (0) means start numbering parameters at 0. expo(3) is a substitute for exp([3]+[4]*x).

Case 2: inline expression using TMath functions with parameters

{
auto fb2 = new TF1("fa3","TMath::Landau(x,[0],[1],0)",-5,10);
fb2->SetParameters(0.2,1.3);
fb2->Draw();
}

3 - A lambda expression with variables and parameters

Since
6.00/00: TF1 supports using lambda expressions in the formula. This allows, by using a full C++ syntax the full power of lambda functions and still maintain the capability of storing the function in a file which cannot be done with function pointer or lambda written not as expression, but as code (see items below).

Example on how using lambda to define a sum of two functions. Note that is necessary to provide the number of parameters

TF1 f1("f1","sin(x)",0,10);
TF1 f2("f2","cos(x)",0,10);
TF1 fsum("f1","[&](double *x, double *p){ return p[0]*f1(x) + p[1]*f2(x); }",0,10,2);
1-Dim function class
Definition TF1.h:214
TF1 * f1
Definition legend1.C:11

4 - A general C function with parameters

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()
{
auto f1 = new TF1("myfunc",myfunction,0,10,2);
f1->SetParNames("constant","coefficient");
f1->Draw();
}
void myfit()
{
auto h1 = new TH1F("h1","test",100,0,10);
h1->FillRandom("myfunc",20000);
TF1 *f1 = (TF1 *)gROOT->GetFunction("myfunc");
f1->SetParameters(800,1);
h1->Fit("myfunc");
}
#define f(i)
Definition RSha256.hxx:104
float Float_t
Definition RtypesCore.h:57
#define gROOT
Definition TROOT.h:407
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1335
virtual void SetParameters(const Double_t *params)
Definition TF1.h:650
virtual void SetParNames(const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
Set up to 10 parameter names.
Definition TF1.cxx:3462
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
virtual void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Fill histogram following distribution in function fname.
Definition TH1.cxx:3520
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition TH1.cxx:3901
TH1F * h1
Definition legend1.C:5
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123

In an interactive session you can do:

Root > .L myfunc.C
Root > myfunc();
Root > myfit();

TF1 objects can reference other TF1 objects of type A or B defined above. This excludes CLing or compiled functions. However, there is a restriction. A function cannot reference a basic function if the basic function is a polynomial polN.

Example:

{
auto fcos = new TF1 ("fcos", "[0]*cos(x)", 0., 10.);
fcos->SetParNames( "cos");
fcos->SetParameter( 0, 1.1);
auto fsin = new TF1 ("fsin", "[0]*sin(x)", 0., 10.);
fsin->SetParNames( "sin");
fsin->SetParameter( 0, 2.1);
auto fsincos = new TF1 ("fsc", "fcos+fsin");
auto fs2 = new TF1 ("fs2", "fsc+fsc");
}

5 - A general C++ function object (functor) with parameters

A TF1 can be created from any C++ class implementing the operator()(double *x, double *p). The advantage of the function object is that he can have a state and reference therefore what-ever other object. In this way the user can customize his function.

Example:

class MyFunctionObject {
public:
// use constructor to customize your function object
double operator() (double *x, double *p) {
// function implementation using class data members
}
};
{
....
MyFunctionObject fobj;
auto f = new TF1("f",fobj,0,1,npar); // create TF1 class.
.....
}
winID h TVirtualViewer3D TVirtualGLPainter p
TRObject operator()(const T1 &t1) const

Using a lambda function as a general C++ functor object

From C++11 we can use both std::function or even better lambda functions to create the TF1. As above the lambda must have the right signature but can capture whatever we want. For example we can make a TF1 from the TGraph::Eval function as shown below where we use as function parameter the graph normalization.

auto g = new TGraph(npointx, xvec, yvec);
auto f = new TF1("f",[&](double*x, double *p){ return p[0]*g->Eval(x[0]); }, xmin, xmax, 1);
#define g(i)
Definition RSha256.hxx:105
float xmin
float xmax
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41

6 - A member function with parameters of a general C++ class

A TF1 can be created in this case from any member function of a class which has the signature of (double * , double *) and returning a double.

Example:

class MyFunction {
public:
...
double Evaluate() (double *x, double *p) {
// function implementation
}
};
{
....
MyFunction *fptr = new MyFunction(....); // create the user function class
auto f = new TF1("f",fptr,&MyFunction::Evaluate,0,1,npar,"MyFunction","Evaluate"); // create TF1 class.
.....
}

See also the tutorial math/exampleFunctor.C for a running example.

Definition at line 214 of file TF1.h.

Classes

struct  TF1FunctorPointer
 
struct  TF1FunctorPointerImpl
 

Public Types

enum class  EAddToList { kDefault , kAdd , kNo }
 Add to list behavior. More...
 
enum  EStatusBits { kNotGlobal = (1ULL << ( 10 )) , kNotDraw = (1ULL << ( 9 )) }
 
- Public Types inherited from TObject
enum  {
  kIsOnHeap = 0x01000000 , kNotDeleted = 0x02000000 , kZombie = 0x04000000 , kInconsistent = 0x08000000 ,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = (1ULL << ( 0 )) , kOverwrite = (1ULL << ( 1 )) , kWriteDelete = (1ULL << ( 2 )) }
 
enum  EDeprecatedStatusBits { kObjInCanvas = (1ULL << ( 3 )) }
 
enum  EStatusBits {
  kCanDelete = (1ULL << ( 0 )) , kMustCleanup = (1ULL << ( 3 )) , kIsReferenced = (1ULL << ( 4 )) , kHasUUID = (1ULL << ( 5 )) ,
  kCannotPick = (1ULL << ( 6 )) , kNoContextMenu = (1ULL << ( 8 )) , kInvalidObject = (1ULL << ( 13 ))
}
 

Public Member Functions

 TF1 ()
 TF1 default constructor.
 
 TF1 (const char *name, const char *formula, Double_t xmin, Double_t xmax, Option_t *option)
 Same constructor as above (for TFormula based function) but passing an option strings available options VEC - vectorize the formula expressions (not possible for lambda based expressions) NL - function is not stores in the global list of functions GL - function will be always stored in the global list of functions , independently of the global setting of TF1::DefaultAddToGlobalList.
 
 TF1 (const char *name, const char *formula, Double_t xmin=0, Double_t xmax=1, EAddToList addToGlobList=EAddToList::kDefault, bool vectorize=false)
 F1 constructor using a formula definition.
 
template<class PtrObj , typename MemFn >
 TF1 (const char *name, const PtrObj &p, MemFn memFn, Double_t xmin, Double_t xmax, Int_t npar, const char *, const char *, EAddToList addToGlobList=EAddToList::kDefault)
 
template<class PtrObj , typename MemFn >
 TF1 (const char *name, const PtrObj &p, MemFn memFn, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 
 TF1 (const char *name, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 F1 constructor using name of an interpreted function.
 
 TF1 (const char *name, Double_t(*fcn)(const Double_t *, const Double_t *), Double_t xmin=0, Double_t xmax=1, Int_t npar=0, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 Constructor using a pointer to (const) real function.
 
 TF1 (const char *name, Double_t(*fcn)(Double_t *, Double_t *), Double_t xmin=0, Double_t xmax=1, Int_t npar=0, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 Constructor using a pointer to a real function.
 
template<typename Func >
 TF1 (const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *, EAddToList addToGlobList=EAddToList::kDefault)
 
template<typename Func >
 TF1 (const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 
 TF1 (const char *name, ROOT::Math::ParamFunctor f, Double_t xmin=0, Double_t xmax=1, Int_t npar=0, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 Constructor using the Functor class.
 
template<class T >
 TF1 (const char *name, std::function< T(const T *data, const Double_t *param)> &fcn, Double_t xmin=0, Double_t xmax=1, Int_t npar=0, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 
template<class T >
 TF1 (const char *name, T(*fcn)(const T *, const Double_t *), Double_t xmin=0, Double_t xmax=1, Int_t npar=0, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
 Constructor using a pointer to function.
 
 TF1 (const TF1 &f1)
 
 ~TF1 () override
 TF1 default destructor.
 
virtual void AddParameter (const TString &name, Double_t value)
 
virtual Bool_t AddToGlobalList (Bool_t on=kTRUE)
 Add to global list of functions (gROOT->GetListOfFunctions() ) return previous status (true if the function was already in the list false if not)
 
void Browse (TBrowser *b) override
 Browse.
 
virtual Double_t CentralMoment (Double_t n, Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
 Return nth central moment of function between a and b (i.e the n-th moment around the mean value)
 
TObjectClone (const char *newname=nullptr) const override
 Make a complete copy of the underlying object.
 
void Copy (TObject &f1) const override
 Copy this F1 to a new F1.
 
virtual TH1CreateHistogram ()
 
virtual Double_t Derivative (Double_t x, Double_t *params=nullptr, Double_t epsilon=0.001) 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.
 
virtual Double_t Derivative2 (Double_t x, Double_t *params=nullptr, Double_t epsilon=0.001) const
 Returns the second 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.
 
virtual Double_t Derivative3 (Double_t x, Double_t *params=nullptr, Double_t epsilon=0.001) const
 Returns the third 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.
 
Int_t DistancetoPrimitive (Int_t px, Int_t py) override
 Compute distance from point px,py to a function.
 
void Draw (Option_t *option="") override
 Draw this function with its current attributes.
 
virtual TF1DrawCopy (Option_t *option="") const
 Draw a copy of this function with its current attributes.
 
virtual TObjectDrawDerivative (Option_t *option="al")
 Draw derivative of this function.
 
virtual void DrawF1 (Double_t xmin, Double_t xmax, Option_t *option="")
 Draw function between xmin and xmax.
 
virtual TObjectDrawIntegral (Option_t *option="al")
 Draw integral of this function.
 
virtual Double_t Eval (Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
 Evaluate this function.
 
virtual Double_t EvalPar (const Double_t *x, const Double_t *params=nullptr)
 Evaluate function with given coordinates and parameters.
 
template<class T >
EvalPar (const T *x, const Double_t *params=nullptr)
 EvalPar for vectorized.
 
void ExecuteEvent (Int_t event, Int_t px, Int_t py) override
 Execute action corresponding to one event.
 
virtual void FixParameter (Int_t ipar, Double_t value)
 Fix the value of a parameter for a fit operation The specified value will be used in the fit and the parameter will be constant (nor varying) during fitting Note that when using pre-defined functions (e.g gaus), one needs to use the fit option 'B' to have the fix of the paramter effective.
 
Double_t GetChisquare () const
 
virtual TString GetExpFormula (Option_t *option="") const
 
virtual TFormulaGetFormula ()
 
virtual const TFormulaGetFormula () const
 
virtual TH1GetHistogram () const
 Return a pointer to the histogram used to visualise the function Note that this histogram is managed by the function and in same case it is automatically deleted when some TF1 functions are called such as TF1::SetParameters, TF1::SetNpx, TF1::SetRange It is then reccomended either to clone the return object or calling again teh GetHistogram function whenever is needed.
 
virtual const TObjectGetLinearPart (Int_t i) const
 
virtual Double_t GetMaximum (Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
 Returns the maximum value of the function.
 
virtual Double_t GetMaximumStored () const
 
virtual Double_t GetMaximumX (Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
 Returns the X value corresponding to the maximum value of the function.
 
TMethodCallGetMethodCall () const
 
virtual Double_t GetMinimum (Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
 Returns the minimum value of the function on the (xmin, xmax) interval.
 
virtual Double_t GetMinimumStored () const
 
virtual Double_t GetMinimumX (Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
 Returns the X value corresponding to the minimum value of the function on the (xmin, xmax) interval.
 
virtual Int_t GetNDF () const
 Return the number of degrees of freedom in the fit the fNDF parameter has been previously computed during a fit.
 
virtual Int_t GetNdim () const
 
virtual Int_t GetNpar () const
 
virtual Int_t GetNpx () const
 
virtual Int_t GetNumber () const
 
virtual Int_t GetNumberFitPoints () const
 
virtual Int_t GetNumberFreeParameters () const
 Return the number of free parameters.
 
char * GetObjectInfo (Int_t px, Int_t py) const override
 Redefines TObject::GetObjectInfo.
 
virtual Double_t GetParameter (const TString &name) const
 
virtual Double_t GetParameter (Int_t ipar) const
 
virtual Double_tGetParameters () const
 
virtual void GetParameters (Double_t *params)
 
TObjectGetParent () const
 
virtual Double_t GetParError (Int_t ipar) const
 Return value of parameter number ipar.
 
virtual const Double_tGetParErrors () const
 
virtual void GetParLimits (Int_t ipar, Double_t &parmin, Double_t &parmax) const
 Return limits for parameter ipar.
 
virtual const char * GetParName (Int_t ipar) const
 
virtual Int_t GetParNumber (const char *name) const
 
virtual Double_t GetProb () const
 Return the fit probability.
 
virtual Int_t GetQuantiles (Int_t nprobSum, Double_t *q, const Double_t *probSum)
 Compute Quantiles for density distribution of this function.
 
virtual Double_t GetRandom (Double_t xmin, Double_t xmax, TRandom *rng=nullptr, Option_t *opt=nullptr)
 Return a random number following this function shape in [xmin,xmax].
 
virtual Double_t GetRandom (TRandom *rng=nullptr, Option_t *opt=nullptr)
 Return a random number following this function shape.
 
virtual void GetRange (Double_t &xmin, Double_t &xmax) const
 Return range of a 1-D function.
 
virtual void GetRange (Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
 Return range of a 2-D function.
 
virtual void GetRange (Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const
 Return range of function.
 
virtual Double_t GetSave (const Double_t *x)
 Get value corresponding to X in array of fSave values.
 
virtual Double_t GetVariable (const TString &name)
 
virtual Double_t GetX (Double_t y, Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
 Returns the X value corresponding to the function value fy for (xmin<x<xmax).
 
TAxisGetXaxis () const
 Get x axis of the function.
 
virtual Double_t GetXmax () const
 
virtual Double_t GetXmin () const
 
TAxisGetYaxis () const
 Get y axis of the function.
 
TAxisGetZaxis () const
 Get z axis of the function. (In case this object is a TF2 or TF3)
 
virtual void GradientPar (const Double_t *x, Double_t *grad, Double_t eps=0.01)
 Compute the gradient wrt parameters If the TF1 object is based on a formula expression (TFormula) and TFormula::GenerateGradientPar() has been successfully called automatic differentiation using CLAD is used instead of the default numerical differentiation.
 
template<class T >
void GradientPar (const T *x, T *grad, Double_t eps=0.01)
 
virtual Double_t GradientPar (Int_t ipar, const Double_t *x, Double_t eps=0.01)
 Compute the gradient (derivative) wrt a parameter ipar.
 
template<class T >
GradientPar (Int_t ipar, const T *x, Double_t eps=0.01)
 
template<class T >
void GradientParTempl (const T *x, T *grad, Double_t eps=0.01)
 
template<class T >
GradientParTempl (Int_t ipar, const T *x, Double_t eps=0.01)
 
virtual void InitArgs (const Double_t *x, const Double_t *params)
 Initialize parameters addresses.
 
virtual Double_t Integral (Double_t a, Double_t b, Double_t epsrel=1.e-12)
 IntegralOneDim or analytical integral.
 
virtual Double_t IntegralError (Double_t a, Double_t b, const Double_t *params=nullptr, const Double_t *covmat=nullptr, Double_t epsilon=1.E-2)
 Return Error on Integral of a parametric function between a and b due to the parameter uncertainties and their covariance matrix from the fit.
 
virtual Double_t IntegralError (Int_t n, const Double_t *a, const Double_t *b, const Double_t *params=nullptr, const Double_t *covmat=nullptr, Double_t epsilon=1.E-2)
 Return Error on Integral of a parametric function with dimension larger than one between a[] and b[] due to the parameters uncertainties.
 
virtual Double_t IntegralFast (Int_t num, Double_t *x, Double_t *w, Double_t a, Double_t b, Double_t *params=nullptr, Double_t epsilon=1e-12)
 Gauss-Legendre integral, see CalcGaussLegendreSamplingPoints.
 
virtual Double_t IntegralMultiple (Int_t n, const Double_t *a, const Double_t *b, Double_t epsrel, Double_t &relerr)
 See more general prototype below.
 
virtual Double_t IntegralMultiple (Int_t n, const Double_t *a, const Double_t *b, Int_t maxpts, Double_t epsrel, Double_t epsabs, Double_t &relerr, Int_t &nfnevl, Int_t &ifail)
 This function computes, to an attempted specified accuracy, the value of the integral.
 
virtual Double_t IntegralMultiple (Int_t n, const Double_t *a, const Double_t *b, Int_t, Int_t maxpts, Double_t epsrel, Double_t &relerr, Int_t &nfnevl, Int_t &ifail)
 
virtual Double_t IntegralOneDim (Double_t a, Double_t b, Double_t epsrel, Double_t epsabs, Double_t &err)
 Return Integral of function between a and b using the given parameter values and relative and absolute tolerance.
 
TClassIsA () const override
 
virtual Bool_t IsEvalNormalized () const
 
virtual Bool_t IsInside (const Double_t *x) const
 return kTRUE if the point is inside the function range
 
virtual Bool_t IsLinear () const
 
virtual Bool_t IsValid () const
 Return kTRUE if the function is valid.
 
bool IsVectorized ()
 
virtual Double_t Mean (Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
 
virtual Double_t Moment (Double_t n, Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
 Return nth moment of function between a and b.
 
template<class T >
operator() (const T *x, const Double_t *params=nullptr)
 
virtual Double_t operator() (Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
 
TF1operator= (const TF1 &rhs)
 Operator =.
 
void Paint (Option_t *option="") override
 Paint this function with its current attributes.
 
void Print (Option_t *option="") const override
 This method must be overridden when a class wants to print itself.
 
virtual void ReleaseParameter (Int_t ipar)
 Release parameter number ipar during a fit operation.
 
virtual void Save (Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax)
 Save values of function in array fSave.
 
void SavePrimitive (std::ostream &out, Option_t *option="") override
 Save primitive as a C++ statement(s) on output stream out.
 
virtual void SetChisquare (Double_t chi2)
 
virtual void SetFitResult (const ROOT::Fit::FitResult &result, const Int_t *indpar=nullptr)
 Set the result from the fit parameter values, errors, chi2, etc... Optionally a pointer to a vector (with size fNpar) of the parameter indices in the FitResult can be passed This is useful in the case of a combined fit with different functions, and the FitResult contains the global result By default it is assume that indpar = {0,1,2,....,fNpar-1}.
 
template<typename Func >
void SetFunction (Func f)
 
template<class PtrObj , typename MemFn >
void SetFunction (PtrObj &p, MemFn memFn)
 
virtual void SetMaximum (Double_t maximum=-1111)
 Set the maximum value along Y for this function In case the function is already drawn, set also the maximum in the helper histogram.
 
virtual void SetMinimum (Double_t minimum=-1111)
 Set the minimum value along Y for this function In case the function is already drawn, set also the minimum in the helper histogram.
 
virtual void SetNDF (Int_t ndf)
 Set the number of degrees of freedom ndf should be the number of points used in a fit - the number of free parameters.
 
virtual void SetNormalized (Bool_t flag)
 
virtual void SetNpx (Int_t npx=100)
 Set the number of points used to draw the function.
 
virtual void SetNumberFitPoints (Int_t npfits)
 
virtual void SetParameter (const TString &name, Double_t value)
 
virtual void SetParameter (Int_t param, Double_t value)
 
virtual void SetParameters (const Double_t *params)
 
virtual void SetParameters (double p0, double p1=0.0, double p2=0.0, double p3=0.0, double p4=0.0, double p5=0.0, double p6=0.0, double p7=0.0, double p8=0.0, double p9=0.0, double p10=0.0)
 
virtual void SetParent (TObject *p=nullptr)
 
virtual void SetParError (Int_t ipar, Double_t error)
 Set error for parameter number ipar.
 
virtual void SetParErrors (const Double_t *errors)
 Set errors for all active parameters when calling this function, the array errors must have at least fNpar values.
 
virtual void SetParLimits (Int_t ipar, Double_t parmin, Double_t parmax)
 Set lower and upper limits for parameter ipar.
 
virtual void SetParName (Int_t ipar, const char *name)
 Set name of parameter number ipar.
 
virtual void SetParNames (const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
 Set up to 10 parameter names.
 
virtual void SetRange (Double_t xmin, Double_t xmax)
 Initialize the upper and lower bounds to draw the function.
 
virtual void SetRange (Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax)
 
virtual void SetRange (Double_t xmin, Double_t ymin, Double_t zmin, Double_t xmax, Double_t ymax, Double_t zmax)
 
virtual void SetSavedPoint (Int_t point, Double_t value)
 Restore value of function saved at point.
 
void SetTitle (const char *title="") override
 Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title is "fffffff" and "xxxx" and "yyyy" are the titles for the X and Y axis respectively.
 
virtual void SetVectorized (Bool_t vectorized)
 
void Streamer (TBuffer &) override
 Stream a class object.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
virtual void Update ()
 Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated histogram or Integral.
 
virtual Double_t Variance (Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
 
- Public Member Functions inherited from TNamed
 TNamed ()
 
 TNamed (const char *name, const char *title)
 
 TNamed (const TNamed &named)
 TNamed copy ctor.
 
 TNamed (const TString &name, const TString &title)
 
virtual ~TNamed ()
 TNamed destructor.
 
void Clear (Option_t *option="") override
 Set name and title to empty strings ("").
 
TObjectClone (const char *newname="") const override
 Make a clone of an object using the Streamer facility.
 
Int_t Compare (const TObject *obj) const override
 Compare two TNamed objects.
 
void Copy (TObject &named) const override
 Copy this to obj.
 
virtual void FillBuffer (char *&buffer)
 Encode TNamed into output buffer.
 
const char * GetName () const override
 Returns name of object.
 
const char * GetTitle () const override
 Returns title of object.
 
ULong_t Hash () const override
 Return hash value for this object.
 
TClassIsA () const override
 
Bool_t IsSortable () const override
 
void ls (Option_t *option="") const override
 List TNamed name and title.
 
TNamedoperator= (const TNamed &rhs)
 TNamed assignment operator.
 
void Print (Option_t *option="") const override
 Print TNamed name and title.
 
virtual void SetName (const char *name)
 Set the name of the TNamed.
 
virtual void SetNameTitle (const char *name, const char *title)
 Set all the TNamed parameters (name and title).
 
virtual Int_t Sizeof () const
 Return size of the TNamed part of the TObject.
 
void Streamer (TBuffer &) override
 Stream an object of class TObject.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TObject
 TObject ()
 TObject constructor.
 
 TObject (const TObject &object)
 TObject copy ctor.
 
virtual ~TObject ()
 TObject destructor.
 
void AbstractMethod (const char *method) const
 Use this method to implement an "abstract" method that you don't want to leave purely abstract.
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad.
 
ULong_t CheckedHash ()
 Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object.
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs.
 
virtual void Delete (Option_t *option="")
 Delete this object.
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs.
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad with: gROOT->SetSelectedPad(c1).
 
virtual void Dump () const
 Dump contents of object on stdout.
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message.
 
virtual void Execute (const char *method, const char *params, Int_t *error=nullptr)
 Execute method on this object with the given parameter string, e.g.
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=nullptr)
 Execute method on this object with parameters stored in the TObjArray.
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message.
 
virtual TObjectFindObject (const char *name) const
 Must be redefined in derived classes.
 
virtual TObjectFindObject (const TObject *obj) const
 Must be redefined in derived classes.
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object.
 
virtual const char * GetIconName () const
 Returns mime type name of object.
 
virtual Option_tGetOption () const
 
virtual UInt_t GetUniqueID () const
 Return the unique object id.
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out.
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e.
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message.
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname".
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl.
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas.
 
void InvertBit (UInt_t f)
 
Bool_t IsDestructed () const
 IsDestructed.
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory).
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary).
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification (the base implementation is no-op).
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete.
 
void operator delete (void *ptr)
 Operator delete.
 
void operator delete (void *ptr, void *vp)
 Only called by placement new when throwing an exception.
 
void operator delete[] (void *ptr)
 Operator delete [].
 
void operator delete[] (void *ptr, void *vp)
 Only called by placement new[] when throwing an exception.
 
void * operator new (size_t sz)
 
void * operator new (size_t sz, void *vp)
 
void * operator new[] (size_t sz)
 
void * operator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 TObject assignment operator.
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list.
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory.
 
virtual void RecursiveRemove (TObject *obj)
 Recursively remove this object from a list.
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename.
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f.
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object.
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message.
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void UseCurrentStyle ()
 Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked.
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory.
 
- Public Member Functions inherited from TAttLine
 TAttLine ()
 AttLine default constructor.
 
 TAttLine (Color_t lcolor, Style_t lstyle, Width_t lwidth)
 AttLine normal constructor.
 
virtual ~TAttLine ()
 AttLine destructor.
 
void Copy (TAttLine &attline) const
 Copy this line attributes to a new TAttLine.
 
Int_t DistancetoLine (Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
 Compute distance from point px,py to a line.
 
virtual Color_t GetLineColor () const
 Return the line color.
 
virtual Style_t GetLineStyle () const
 Return the line style.
 
virtual Width_t GetLineWidth () const
 Return the line width.
 
virtual void Modify ()
 Change current line attributes if necessary.
 
virtual void ResetAttLine (Option_t *option="")
 Reset this line attributes to default values.
 
virtual void SaveLineAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
 Save line attributes as C++ statement(s) on output stream out.
 
virtual void SetLineAttributes ()
 Invoke the DialogCanvas Line attributes.
 
virtual void SetLineColor (Color_t lcolor)
 Set the line color.
 
virtual void SetLineColorAlpha (Color_t lcolor, Float_t lalpha)
 Set a transparent line color.
 
virtual void SetLineStyle (Style_t lstyle)
 Set the line style.
 
virtual void SetLineWidth (Width_t lwidth)
 Set the line width.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TAttFill
 TAttFill ()
 AttFill default constructor.
 
 TAttFill (Color_t fcolor, Style_t fstyle)
 AttFill normal constructor.
 
virtual ~TAttFill ()
 AttFill destructor.
 
void Copy (TAttFill &attfill) const
 Copy this fill attributes to a new TAttFill.
 
virtual Color_t GetFillColor () const
 Return the fill area color.
 
virtual Style_t GetFillStyle () const
 Return the fill area style.
 
virtual Bool_t IsTransparent () const
 
virtual void Modify ()
 Change current fill area attributes if necessary.
 
virtual void ResetAttFill (Option_t *option="")
 Reset this fill attributes to default values.
 
virtual void SaveFillAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
 Save fill attributes as C++ statement(s) on output stream out.
 
virtual void SetFillAttributes ()
 Invoke the DialogCanvas Fill attributes.
 
virtual void SetFillColor (Color_t fcolor)
 Set the fill area color.
 
virtual void SetFillColorAlpha (Color_t fcolor, Float_t falpha)
 Set a transparent fill color.
 
virtual void SetFillStyle (Style_t fstyle)
 Set the fill area style.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TAttMarker
 TAttMarker ()
 TAttMarker default constructor.
 
 TAttMarker (Color_t color, Style_t style, Size_t msize)
 TAttMarker normal constructor.
 
virtual ~TAttMarker ()
 TAttMarker destructor.
 
void Copy (TAttMarker &attmarker) const
 Copy this marker attributes to a new TAttMarker.
 
virtual Color_t GetMarkerColor () const
 Return the marker color.
 
virtual Size_t GetMarkerSize () const
 Return the marker size.
 
virtual Style_t GetMarkerStyle () const
 Return the marker style.
 
virtual void Modify ()
 Change current marker attributes if necessary.
 
virtual void ResetAttMarker (Option_t *toption="")
 Reset this marker attributes to the default values.
 
virtual void SaveMarkerAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
 Save line attributes as C++ statement(s) on output stream out.
 
virtual void SetMarkerAttributes ()
 Invoke the DialogCanvas Marker attributes.
 
virtual void SetMarkerColor (Color_t mcolor=1)
 Set the marker color.
 
virtual void SetMarkerColorAlpha (Color_t mcolor, Float_t malpha)
 Set a transparent marker color.
 
virtual void SetMarkerSize (Size_t msize=1)
 Set the marker size.
 
virtual void SetMarkerStyle (Style_t mstyle=1)
 Set the marker style.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 

Static Public Member Functions

static void AbsValue (Bool_t reject=kTRUE)
 Static function: set the fgAbsValue flag.
 
static void CalcGaussLegendreSamplingPoints (Int_t num, Double_t *x, Double_t *w, Double_t eps=3.0e-11)
 Type safe interface (static method) The number of sampling points are taken from the TGraph.
 
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Bool_t DefaultAddToGlobalList (Bool_t on=kTRUE)
 Static method to add/avoid to add automatically functions to the global list (gROOT->GetListOfFunctions() ) After having called this static method, all the functions created afterwards will follow the desired behaviour.
 
static Double_t DerivativeError ()
 Static function returning the error of the last call to the of Derivative's functions.
 
static TF1GetCurrent ()
 Static function returning the current function being processed.
 
static void InitStandardFunctions ()
 Create the basic function objects.
 
static Bool_t RejectedPoint ()
 See TF1::RejectPoint above.
 
static void RejectPoint (Bool_t reject=kTRUE)
 Static function to set the global flag to reject points the fgRejectPoint global flag is tested by all fit functions if TRUE the point is not included in the fit.
 
static void SetCurrent (TF1 *f1)
 Static function setting the current function.
 
- Static Public Member Functions inherited from TNamed
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TObject
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Longptr_t GetDtorOnly ()
 Return destructor only flag.
 
static Bool_t GetObjectStat ()
 Get status of object stat flag.
 
static void SetDtorOnly (void *obj)
 Set destructor only flag.
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable.
 
- Static Public Member Functions inherited from TAttLine
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TAttFill
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TAttMarker
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Width_t GetMarkerLineWidth (Style_t style)
 Internal helper function that returns the line width of the given marker style (0 = filled marker)
 
static Style_t GetMarkerStyleBase (Style_t style)
 Internal helper function that returns the corresponding marker style with line width 1 for the given style.
 

Protected Types

enum  EFType {
  kFormula = 0 , kPtrScalarFreeFcn , kInterpreted , kTemplVec ,
  kTemplScalar , kCompositionFcn
}
 
- Protected Types inherited from TObject
enum  { kOnlyPrepStep = (1ULL << ( 3 )) }
 

Protected Member Functions

 TF1 (EFType functionType, const char *name, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim, EAddToList addToGlobList, TF1Parameters *params=nullptr, TF1FunctorPointer *functor=nullptr)
 General constructor for TF1. Most of the other constructors delegate on it.
 
Bool_t ComputeCdfTable (Option_t *opt)
 Compute the cumulative function at fNpx points between fXmin and fXmax.
 
virtual TH1DoCreateHistogram (Double_t xmin, Double_t xmax, Bool_t recreate=kFALSE)
 Create histogram with bin content equal to function value computed at the bin center This histogram will be used to paint the function A re-creation is forced and a new histogram is done if recreate=true.
 
void DoInitialize (EAddToList addToGlobList)
 Common initialization of the TF1.
 
virtual Double_t GetMinMaxNDim (Double_t *x, Bool_t findmax, Double_t epsilon=0, Int_t maxiter=0) const
 Find the minimum of a function of whatever dimension.
 
virtual void GetRange (Double_t *xmin, Double_t *xmax) const
 Return range of a generic N-D function.
 
void IntegrateForNormalization ()
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected).
 
void MakeZombie ()
 

Protected Attributes

std::vector< Double_tfAlpha
 ! Array alpha. for each bin in x the deconvolution r of fIntegral
 
std::vector< Double_tfBeta
 ! Array beta. is approximated by x = alpha +beta*r *gamma*r**2
 
Double_t fChisquare {}
 Function fit chisquare.
 
std::unique_ptr< TF1AbsCompositionfComposition
 Pointer to composition (NSUM or CONV)
 
std::unique_ptr< TFormulafFormula
 Pointer to TFormula in case when user define formula.
 
std::unique_ptr< TF1FunctorPointerfFunctor
 ! Functor object to wrap any C++ callable object
 
std::vector< Double_tfGamma
 ! Array gamma.
 
TH1fHistogram {nullptr}
 ! Pointer to histogram used for visualisation
 
std::vector< Double_tfIntegral
 ! Integral of function binned on fNpx bins
 
Double_t fMaximum {-1111}
 Maximum value for plotting.
 
std::unique_ptr< TMethodCallfMethodCall
 ! Pointer to MethodCall in case of interpreted function
 
Double_t fMinimum {-1111}
 Minimum value for plotting.
 
Int_t fNDF {}
 Number of degrees of freedom in the fit.
 
Int_t fNdim {}
 Function dimension.
 
Bool_t fNormalized {false}
 Normalization option (false by default)
 
Double_t fNormIntegral {}
 Integral of the function before being normalized.
 
Int_t fNpar {}
 Number of parameters.
 
Int_t fNpfits {}
 Number of points used in the fit.
 
Int_t fNpx {100}
 Number of points used for the graphical representation.
 
std::unique_ptr< TF1ParametersfParams
 Pointer to Function parameters object (exists only for not-formula functions)
 
TObjectfParent {nullptr}
 ! Parent object hooking this function (if one)
 
std::vector< Double_tfParErrors
 Array of errors of the fNpar parameters.
 
std::vector< Double_tfParMax
 Array of upper limits of the fNpar parameters.
 
std::vector< Double_tfParMin
 Array of lower limits of the fNpar parameters.
 
std::vector< Double_tfSave
 Array of fNsave function values.
 
EFType fType {EFType::kTemplScalar}
 
Double_t fXmax {-1111}
 Upper bounds for the range.
 
Double_t fXmin {-1111}
 Lower bounds for the range.
 
- Protected Attributes inherited from TNamed
TString fName
 
TString fTitle
 
- Protected Attributes inherited from TAttLine
Color_t fLineColor
 Line color.
 
Style_t fLineStyle
 Line style.
 
Width_t fLineWidth
 Line width.
 
- Protected Attributes inherited from TAttFill
Color_t fFillColor
 Fill area color.
 
Style_t fFillStyle
 Fill area style.
 
- Protected Attributes inherited from TAttMarker
Color_t fMarkerColor
 Marker color.
 
Size_t fMarkerSize
 Marker size.
 
Style_t fMarkerStyle
 Marker style.
 

Static Protected Attributes

static std::atomic< Bool_tfgAbsValue
 
static std::atomic< Bool_tfgAddToGlobList
 
static TF1fgCurrent = nullptr
 
static Bool_t fgRejectPoint = kFALSE
 

Private Member Functions

void DefineNSUMTerm (TObjArray *newFuncs, TObjArray *coeffNames, TString &fullFormula, TString &formula, int termStart, int termEnd, Double_t xmin, Double_t xmax)
 Helper functions for NSUM parsing.
 
template<class T >
EvalParTempl (const T *data, const Double_t *params=nullptr)
 Eval for vectorized functions.
 
int TermCoeffLength (TString &term)
 

Friends

template<class Func >
struct ROOT::Internal::TF1Builder
 

#include <TF1.h>

Inheritance diagram for TF1:
[legend]

Member Enumeration Documentation

◆ EAddToList

enum class TF1::EAddToList
strong

Add to list behavior.

Enumerator
kDefault 
kAdd 
kNo 

Definition at line 221 of file TF1.h.

◆ EFType

enum TF1::EFType
protected
Enumerator
kFormula 

Formula functions which can be stored,.

kPtrScalarFreeFcn 

Pointer to scalar free function,.

kInterpreted 

Interpreted functions constructed by name,.

kTemplVec 

Vectorized free functions or TemplScalar functors evaluating on vectorized parameters,.

kTemplScalar 

TemplScalar functors evaluating on scalar parameters.

kCompositionFcn 

Definition at line 235 of file TF1.h.

◆ EStatusBits

Enumerator
kNotGlobal 
kNotDraw 

Definition at line 325 of file TF1.h.

Constructor & Destructor Documentation

◆ TF1() [1/15]

TF1::TF1 ( EFType  functionType,
const char *  name,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
Int_t  ndim,
EAddToList  addToGlobList,
TF1Parameters params = nullptr,
TF1FunctorPointer functor = nullptr 
)
inlineprotected

General constructor for TF1. Most of the other constructors delegate on it.

Definition at line 274 of file TF1.h.

◆ TF1() [2/15]

TF1::TF1 ( )

TF1 default constructor.

Definition at line 489 of file TF1.cxx.

◆ TF1() [3/15]

TF1::TF1 ( const char *  name,
const char *  formula,
Double_t  xmin = 0,
Double_t  xmax = 1,
EAddToList  addToGlobList = EAddToList::kDefault,
bool  vectorize = false 
)

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 formula has the form "fffffff;xxxx;yyyy", it is assumed that the formula string is "fffffff" and "xxxx" and "yyyy" are the titles for the X and Y axis respectively.

Definition at line 511 of file TF1.cxx.

◆ TF1() [4/15]

TF1::TF1 ( const char *  name,
const char *  formula,
Double_t  xmin,
Double_t  xmax,
Option_t opt 
)

Same constructor as above (for TFormula based function) but passing an option strings available options VEC - vectorize the formula expressions (not possible for lambda based expressions) NL - function is not stores in the global list of functions GL - function will be always stored in the global list of functions , independently of the global setting of TF1::DefaultAddToGlobalList.

Definition at line 690 of file TF1.cxx.

◆ TF1() [5/15]

TF1::TF1 ( const char *  name,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)

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 C++ function. 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 the C++ interpreter.

Warning
A function created with this constructor cannot be Cloned.

Definition at line 716 of file TF1.cxx.

◆ TF1() [6/15]

TF1::TF1 ( const char *  name,
Double_t(*)(Double_t *, Double_t *)  fcn,
Double_t  xmin = 0,
Double_t  xmax = 1,
Int_t  npar = 0,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)

Constructor using a pointer to a real function.

Parameters
[in]nameobject name
[in]fcnpointer to function
[in]xmin,xmaxx axis limits
[in]nparis the number of free parameters used by the function
[in]ndimnumber of dimensions
[in]addToGlobListboolean marking if it should be added to global list

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.

Warning
A function created with this constructor cannot be Cloned.

Definition at line 755 of file TF1.cxx.

◆ TF1() [7/15]

TF1::TF1 ( const char *  name,
Double_t(*)(const Double_t *, const Double_t *)  fcn,
Double_t  xmin = 0,
Double_t  xmax = 1,
Int_t  npar = 0,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)

Constructor using a pointer to (const) real function.

Parameters
[in]nameobject name
[in]fcnpointer to function
[in]xmin,xmaxx axis limits
[in]nparis the number of free parameters used by the function
[in]ndimnumber of dimensions
[in]addToGlobListboolean marking if it should be added to global list

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.

Warning
A function created with this constructor cannot be Cloned.

Definition at line 777 of file TF1.cxx.

◆ TF1() [8/15]

template<class T >
TF1::TF1 ( const char *  name,
std::function< T(const T *data, const Double_t *param)> &  fcn,
Double_t  xmin = 0,
Double_t  xmax = 1,
Int_t  npar = 0,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Definition at line 338 of file TF1.h.

◆ TF1() [9/15]

template<class T >
TF1::TF1 ( const char *  name,
T(*)(const T *, const Double_t *)  fcn,
Double_t  xmin = 0,
Double_t  xmax = 1,
Int_t  npar = 0,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Constructor using a pointer to function.

Parameters
[in]nameobject name
[in]fcnpointer to function
[in]xmin,xmaxx axis limits
[in]nparis the number of free parameters used by the function
[in]ndimnumber of dimensions
[in]addToGlobListboolean marking if it should be added to global list

This constructor creates a function of type C when invoked with the normal C++ compiler.

Warning
A function created with this constructor cannot be Cloned

Definition at line 362 of file TF1.h.

◆ TF1() [10/15]

TF1::TF1 ( const char *  name,
ROOT::Math::ParamFunctor  f,
Double_t  xmin = 0,
Double_t  xmax = 1,
Int_t  npar = 0,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)

Constructor using the Functor class.

Parameters
[in]nameobject name
fparameterized functor
xminand
xmaxdefine the plotting range of the function
[in]nparis the number of free parameters used by the function
[in]ndimnumber of dimensions
[in]addToGlobListboolean marking if it should be added to global list

This constructor can be used only in compiled code

WARNING! A function created with this constructor cannot be Cloned.

Definition at line 796 of file TF1.cxx.

◆ TF1() [11/15]

template<typename Func >
TF1::TF1 ( const char *  name,
Func  f,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Definition at line 378 of file TF1.h.

◆ TF1() [12/15]

template<typename Func >
TF1::TF1 ( const char *  name,
Func  f,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
const char *  ,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Definition at line 387 of file TF1.h.

◆ TF1() [13/15]

template<class PtrObj , typename MemFn >
TF1::TF1 ( const char *  name,
const PtrObj &  p,
MemFn  memFn,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
Int_t  ndim = 1,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Definition at line 403 of file TF1.h.

◆ TF1() [14/15]

template<class PtrObj , typename MemFn >
TF1::TF1 ( const char *  name,
const PtrObj &  p,
MemFn  memFn,
Double_t  xmin,
Double_t  xmax,
Int_t  npar,
const char *  ,
const char *  ,
EAddToList  addToGlobList = EAddToList::kDefault 
)
inline

Definition at line 409 of file TF1.h.

◆ TF1() [15/15]

TF1::TF1 ( const TF1 f1)

Definition at line 972 of file TF1.cxx.

◆ ~TF1()

TF1::~TF1 ( )
override

TF1 default destructor.

Definition at line 955 of file TF1.cxx.

Member Function Documentation

◆ AbsValue()

void TF1::AbsValue ( Bool_t  flag = kTRUE)
static

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.

Definition at line 986 of file TF1.cxx.

◆ AddParameter()

virtual void TF1::AddParameter ( const TString name,
Double_t  value 
)
inlinevirtual

Definition at line 416 of file TF1.h.

◆ AddToGlobalList()

Bool_t TF1::AddToGlobalList ( Bool_t  on = kTRUE)
virtual

Add to global list of functions (gROOT->GetListOfFunctions() ) return previous status (true if the function was already in the list false if not)

Definition at line 849 of file TF1.cxx.

◆ Browse()

void TF1::Browse ( TBrowser b)
overridevirtual

Browse.

Reimplemented from TObject.

Definition at line 995 of file TF1.cxx.

◆ CalcGaussLegendreSamplingPoints()

void TF1::CalcGaussLegendreSamplingPoints ( Int_t  num,
Double_t x,
Double_t w,
Double_t  eps = 3.0e-11 
)
static

Type safe interface (static method) The number of sampling points are taken from the TGraph.

Type: unsafe but fast interface filling the arrays x and w (static method)

Given the number of sampling points this routine fills the arrays x and w of length num, containing the abscissa and weight of the Gauss-Legendre n-point quadrature formula.

Gauss-Legendre:

\[ W(x)=1 -1<x<1 \\ (j+1)P_{j+1} = (2j+1)xP_j-jP_{j-1} \]

num is the number of sampling points (>0) x and w are arrays of size num eps is the relative precision

If num<=0 or eps<=0 no action is done.

Reference: Numerical Recipes in C, Second Edition

Definition at line 3791 of file TF1.cxx.

◆ CentralMoment()

Double_t TF1::CentralMoment ( Double_t  n,
Double_t  a,
Double_t  b,
const Double_t params = nullptr,
Double_t  epsilon = 0.000001 
)
virtual

Return nth central moment of function between a and b (i.e the n-th moment around the mean value)

See TF1::Integral() for parameter definitions

Author
Gene Van Buren gene@.nosp@m.bnl..nosp@m.gov

Definition at line 3704 of file TF1.cxx.

◆ Class()

static TClass * TF1::Class ( )
static
Returns
TClass describing this class

◆ Class_Name()

static const char * TF1::Class_Name ( )
static
Returns
Name of this class

◆ Class_Version()

static constexpr Version_t TF1::Class_Version ( )
inlinestaticconstexpr
Returns
Version of this class

Definition at line 720 of file TF1.h.

◆ Clone()

TObject * TF1::Clone ( const char *  newname = nullptr) const
overridevirtual

Make a complete copy of the underlying object.

If 'newname' is set, the copy's name will be set to that name.

Reimplemented from TObject.

Definition at line 1066 of file TF1.cxx.

◆ ComputeCdfTable()

Bool_t TF1::ComputeCdfTable ( Option_t option)
protected

Compute the cumulative function at fNpx points between fXmin and fXmax.

Option can be used to force a log scale (option = "log"), linear (option = "lin") or automatic if empty.

Definition at line 2081 of file TF1.cxx.

◆ Copy()

void TF1::Copy ( TObject obj) const
overridevirtual

Copy this F1 to a new F1.

Note that the cached integral with its related arrays are not copied (they are also set as transient data members)

Reimplemented from TObject.

Reimplemented in TF12, TF2, and TF3.

Definition at line 1007 of file TF1.cxx.

◆ CreateHistogram()

virtual TH1 * TF1::CreateHistogram ( )
inlinevirtual

Reimplemented in TF2, and TF3.

Definition at line 455 of file TF1.h.

◆ DeclFileName()

static const char * TF1::DeclFileName ( )
inlinestatic
Returns
Name of the file containing the class declaration

Definition at line 720 of file TF1.h.

◆ DefaultAddToGlobalList()

Bool_t TF1::DefaultAddToGlobalList ( Bool_t  on = kTRUE)
static

Static method to add/avoid to add automatically functions to the global list (gROOT->GetListOfFunctions() ) After having called this static method, all the functions created afterwards will follow the desired behaviour.

By default the functions are added automatically It returns the previous status (true if the functions are added automatically)

Definition at line 840 of file TF1.cxx.

◆ DefineNSUMTerm()

void TF1::DefineNSUMTerm ( TObjArray newFuncs,
TObjArray coeffNames,
TString fullFormula,
TString formula,
int  termStart,
int  termEnd,
Double_t  xmin,
Double_t  xmax 
)
private

Helper functions for NSUM parsing.

Definition at line 885 of file TF1.cxx.

◆ Derivative()

Double_t TF1::Derivative ( Double_t  x,
Double_t params = nullptr,
Double_t  eps = 0.001 
) const
virtual

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) = \frac{f(x+h) - f(x-h)}{2h} \]

the final estimate

\[ D = \frac{4D(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)^{2} + x^{2}deriv^{2}}\sqrt{\sum ai^{2}}, \]

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

Definition at line 1115 of file TF1.cxx.

◆ Derivative2()

Double_t TF1::Derivative2 ( Double_t  x,
Double_t params = nullptr,
Double_t  eps = 0.001 
) const
virtual

Returns the second 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) = \frac{f(x+h) - 2f(x) + f(x-h)}{h^{2}} \]

the final estimate

\[ D = \frac{4D(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)^{2} + x^{2}deriv^{2}}\sqrt{\sum ai^{2}}, \]

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

Definition at line 1180 of file TF1.cxx.

◆ Derivative3()

Double_t TF1::Derivative3 ( Double_t  x,
Double_t params = nullptr,
Double_t  eps = 0.001 
) const
virtual

Returns the third 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) = \frac{f(x+2h) - 2f(x+h) + 2f(x-h) - f(x-2h)}{2h^{3}} \]

the final estimate

\[ D = \frac{4D(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)^{2} + x^{2}deriv^{2}}\sqrt{\sum ai^{2}}, \]

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

Definition at line 1245 of file TF1.cxx.

◆ DerivativeError()

Double_t TF1::DerivativeError ( )
static

Static function returning the error of the last call to the of Derivative's functions.

Definition at line 1279 of file TF1.cxx.

◆ DistancetoPrimitive()

Int_t TF1::DistancetoPrimitive ( Int_t  px,
Int_t  py 
)
overridevirtual

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.

Reimplemented from TObject.

Reimplemented in TF2, and TF3.

Definition at line 1295 of file TF1.cxx.

◆ DoCreateHistogram()

TH1 * TF1::DoCreateHistogram ( Double_t  xmin,
Double_t  xmax,
Bool_t  recreate = kFALSE 
)
protectedvirtual

Create histogram with bin content equal to function value computed at the bin center This histogram will be used to paint the function A re-creation is forced and a new histogram is done if recreate=true.

Definition at line 3048 of file TF1.cxx.

◆ DoInitialize()

void TF1::DoInitialize ( EAddToList  addToGlobalList)
protected

Common initialization of the TF1.

Add to the global list and set the default style

Definition at line 804 of file TF1.cxx.

◆ Draw()

void TF1::Draw ( Option_t option = "")
overridevirtual

Draw this function with its current attributes.

Possible option values are:

option description
"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.

Reimplemented from TObject.

Reimplemented in TF2, and TF3.

Definition at line 1335 of file TF1.cxx.

◆ DrawCopy()

TF1 * TF1::DrawCopy ( Option_t option = "") const
virtual

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:

option description
"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"

Reimplemented in TF12, and TF2.

Definition at line 1365 of file TF1.cxx.

◆ DrawDerivative()

TObject * TF1::DrawDerivative ( Option_t option = "al")
virtual

Draw derivative of this function.

An intermediate TGraph object is built and drawn with option. The function returns a pointer to the TGraph object. Do:

TGraph *g = (TGraph*)myfunc.DrawDerivative(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.

Reimplemented in TF2, and TF3.

Definition at line 1387 of file TF1.cxx.

◆ DrawF1()

void TF1::DrawF1 ( Double_t  xmin,
Double_t  xmax,
Option_t option = "" 
)
virtual

Draw function between xmin and xmax.

Definition at line 1422 of file TF1.cxx.

◆ DrawIntegral()

TObject * TF1::DrawIntegral ( Option_t option = "al")
virtual

Draw integral of this function.

An intermediate TGraph object is built and drawn with option. The function returns a pointer to the TGraph object. Do:

TGraph *g = (TGraph*)myfunc.DrawIntegral(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.

Reimplemented in TF2, and TF3.

Definition at line 1409 of file TF1.cxx.

◆ Eval()

Double_t TF1::Eval ( Double_t  x,
Double_t  y = 0,
Double_t  z = 0,
Double_t  t = 0 
) const
virtual

Evaluate this function.

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

Reimplemented in TF12.

Definition at line 1441 of file TF1.cxx.

◆ EvalPar() [1/2]

Double_t TF1::EvalPar ( const Double_t x,
const Double_t params = nullptr 
)
virtual

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-dimensional 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 responsibility 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 every time these addresses change.

Reimplemented in TF12.

Definition at line 1470 of file TF1.cxx.

◆ EvalPar() [2/2]

template<class T >
T TF1::EvalPar ( const T *  x,
const Double_t params = nullptr 
)

EvalPar for vectorized.

Definition at line 776 of file TF1.h.

◆ EvalParTempl()

template<class T >
T TF1::EvalParTempl ( const T *  data,
const Double_t params = nullptr 
)
inlineprivate

Eval for vectorized functions.

Definition at line 801 of file TF1.h.

◆ ExecuteEvent()

void TF1::ExecuteEvent ( Int_t  event,
Int_t  px,
Int_t  py 
)
overridevirtual

Execute action corresponding to one event.

This member function is called when a F1 is clicked with the locator

Reimplemented from TObject.

Reimplemented in TF2, and TF3.

Definition at line 1538 of file TF1.cxx.

◆ FixParameter()

void TF1::FixParameter ( Int_t  ipar,
Double_t  value 
)
virtual

Fix the value of a parameter for a fit operation The specified value will be used in the fit and the parameter will be constant (nor varying) during fitting Note that when using pre-defined functions (e.g gaus), one needs to use the fit option 'B' to have the fix of the paramter effective.

See TH1::Fit(TF1*, Option_t *, Option_t *, Double_t, Double_t) for the fitting documentation and the fitting options.

Definition at line 1559 of file TF1.cxx.

◆ GetChisquare()

Double_t TF1::GetChisquare ( ) const
inline

Definition at line 450 of file TF1.h.

◆ GetCurrent()

TF1 * TF1::GetCurrent ( )
static

Static function returning the current function being processed.

Definition at line 1571 of file TF1.cxx.

◆ GetExpFormula()

virtual TString TF1::GetExpFormula ( Option_t option = "") const
inlinevirtual

Definition at line 467 of file TF1.h.

◆ GetFormula() [1/2]

virtual TFormula * TF1::GetFormula ( )
inlinevirtual

Definition at line 459 of file TF1.h.

◆ GetFormula() [2/2]

virtual const TFormula * TF1::GetFormula ( ) const
inlinevirtual

Definition at line 463 of file TF1.h.

◆ GetHistogram()

TH1 * TF1::GetHistogram ( ) const
virtual

Return a pointer to the histogram used to visualise the function Note that this histogram is managed by the function and in same case it is automatically deleted when some TF1 functions are called such as TF1::SetParameters, TF1::SetNpx, TF1::SetRange It is then reccomended either to clone the return object or calling again teh GetHistogram function whenever is needed.

Definition at line 1586 of file TF1.cxx.

◆ GetLinearPart()

virtual const TObject * TF1::GetLinearPart ( Int_t  i) const
inlinevirtual

Definition at line 471 of file TF1.h.

◆ GetMaximum()

Double_t TF1::GetMaximum ( Double_t  xmin = 0,
Double_t  xmax = 0,
Double_t  epsilon = 1.E-10,
Int_t  maxiter = 100,
Bool_t  logx = false 
) const
virtual

Returns 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 epsilon (default = 1.E-10) controls the relative accuracy (if |x| > 1 ) and absolute (if |x| < 1) and maxiter (default = 100) controls the maximum number of iteration of the Brent algorithm If the flag logx is set the grid search is done in log step size This is done automatically if the log scale is set in the current Pad

NOTE: see also TF1::GetMaximumX and TF1::GetX

Reimplemented in TF2.

Definition at line 1616 of file TF1.cxx.

◆ GetMaximumStored()

virtual Double_t TF1::GetMaximumStored ( ) const
inlinevirtual

Definition at line 479 of file TF1.h.

◆ GetMaximumX()

Double_t TF1::GetMaximumX ( Double_t  xmin = 0,
Double_t  xmax = 0,
Double_t  epsilon = 1.E-10,
Int_t  maxiter = 100,
Bool_t  logx = false 
) const
virtual

Returns 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 epsilon (default = 1.E-10) controls the relative accuracy (if |x| > 1 ) and absolute (if |x| < 1) and maxiter (default = 100) controls the maximum number of iteration of the Brent algorithm If the flag logx is set the grid search is done in log step size This is done automatically if the log scale is set in the current Pad

NOTE: see also TF1::GetX

Definition at line 1657 of file TF1.cxx.

◆ GetMethodCall()

TMethodCall * TF1::GetMethodCall ( ) const
inline

Definition at line 500 of file TF1.h.

◆ GetMinimum()

Double_t TF1::GetMinimum ( Double_t  xmin = 0,
Double_t  xmax = 0,
Double_t  epsilon = 1.E-10,
Int_t  maxiter = 100,
Bool_t  logx = false 
) const
virtual

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 epsilon (default = 1.E-10) controls the relative accuracy (if |x| > 1 ) and absolute (if |x| < 1) and maxiter (default = 100) controls the maximum number of iteration of the Brent algorithm If the flag logx is set the grid search is done in log step size This is done automatically if the log scale is set in the current Pad

NOTE: see also TF1::GetMaximumX and TF1::GetX

Reimplemented in TF2.

Definition at line 1698 of file TF1.cxx.

◆ GetMinimumStored()

virtual Double_t TF1::GetMinimumStored ( ) const
inlinevirtual

Definition at line 483 of file TF1.h.

◆ GetMinimumX()

Double_t TF1::GetMinimumX ( Double_t  xmin = 0,
Double_t  xmax = 0,
Double_t  epsilon = 1.E-10,
Int_t  maxiter = 100,
Bool_t  logx = false 
) const
virtual

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 epsilon (default = 1.E-10) controls the relative accuracy (if |x| > 1 ) and absolute (if |x| < 1) and maxiter (default = 100) controls the maximum number of iteration of the Brent algorithm If the flag logx is set the grid search is done in log step size This is done automatically if the log scale is set in the current Pad

NOTE: see also TF1::GetX

Definition at line 1825 of file TF1.cxx.

◆ GetMinMaxNDim()

Double_t TF1::GetMinMaxNDim ( Double_t x,
Bool_t  findmax,
Double_t  epsilon = 0,
Int_t  maxiter = 0 
) const
protectedvirtual

Find the minimum of a function of whatever dimension.

While GetMinimum works only for 1D function , GetMinimumNDim works for all dimensions since it uses the minimizer interface vector x at beginning will contained the initial point, on exit will contain the result

Definition at line 1725 of file TF1.cxx.

◆ GetNDF()

Int_t TF1::GetNDF ( ) const
virtual

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.

Definition at line 1891 of file TF1.cxx.

◆ GetNdim()

virtual Int_t TF1::GetNdim ( ) const
inlinevirtual

Definition at line 491 of file TF1.h.

◆ GetNpar()

virtual Int_t TF1::GetNpar ( ) const
inlinevirtual

Definition at line 487 of file TF1.h.

◆ GetNpx()

virtual Int_t TF1::GetNpx ( ) const
inlinevirtual

Definition at line 496 of file TF1.h.

◆ GetNumber()

virtual Int_t TF1::GetNumber ( ) const
inlinevirtual

Definition at line 504 of file TF1.h.

◆ GetNumberFitPoints()

virtual Int_t TF1::GetNumberFitPoints ( ) const
inlinevirtual

Definition at line 509 of file TF1.h.

◆ GetNumberFreeParameters()

Int_t TF1::GetNumberFreeParameters ( ) const
virtual

Return the number of free parameters.

Definition at line 1902 of file TF1.cxx.

◆ GetObjectInfo()

char * TF1::GetObjectInfo ( Int_t  px,
Int_t  py 
) const
overridevirtual

Redefines TObject::GetObjectInfo.

Displays the function info (x, function value) corresponding to cursor position px,py

Reimplemented from TObject.

Reimplemented in TF2.

Definition at line 1920 of file TF1.cxx.

◆ GetParameter() [1/2]

virtual Double_t TF1::GetParameter ( const TString name) const
inlinevirtual

Definition at line 522 of file TF1.h.

◆ GetParameter() [2/2]

virtual Double_t TF1::GetParameter ( Int_t  ipar) const
inlinevirtual

Definition at line 518 of file TF1.h.

◆ GetParameters() [1/2]

virtual Double_t * TF1::GetParameters ( ) const
inlinevirtual

Definition at line 526 of file TF1.h.

◆ GetParameters() [2/2]

virtual void TF1::GetParameters ( Double_t params)
inlinevirtual

Definition at line 530 of file TF1.h.

◆ GetParent()

TObject * TF1::GetParent ( ) const
inline

Definition at line 514 of file TF1.h.

◆ GetParError()

Double_t TF1::GetParError ( Int_t  ipar) const
virtual

Return value of parameter number ipar.

Definition at line 1932 of file TF1.cxx.

◆ GetParErrors()

virtual const Double_t * TF1::GetParErrors ( ) const
inlinevirtual

Definition at line 544 of file TF1.h.

◆ GetParLimits()

void TF1::GetParLimits ( Int_t  ipar,
Double_t parmin,
Double_t parmax 
) const
virtual

Return limits for parameter ipar.

Definition at line 1942 of file TF1.cxx.

◆ GetParName()

virtual const char * TF1::GetParName ( Int_t  ipar) const
inlinevirtual

Definition at line 535 of file TF1.h.

◆ GetParNumber()

virtual Int_t TF1::GetParNumber ( const char *  name) const
inlinevirtual

Definition at line 539 of file TF1.h.

◆ GetProb()

Double_t TF1::GetProb ( ) const
virtual

Return the fit probability.

Definition at line 1957 of file TF1.cxx.

◆ GetQuantiles()

Int_t TF1::GetQuantiles ( Int_t  nprobSum,
Double_t q,
const Double_t probSum 
)
virtual

Compute Quantiles for density distribution of this function.

Quantile x_q of a probability distribution Function F is defined as

\[ F(x_{q}) = \int_{xmin}^{x_{q}} f dx = q with 0 <= q <= 1. \]

For instance the median \( x_{\frac{1}{2}} \) of a distribution is defined as that value of the random variable for which the distribution function equals 0.5:

\[ F(x_{\frac{1}{2}}) = \prod(x < x_{\frac{1}{2}}) = \frac{1}{2} \]

Parameters
[in]nprobSummaximum size of array q and size of array probSum
[out]qarray filled with nq quantiles
[in]probSumarray of positions where quantiles will be computed. It is assumed to contain at least nprobSum values.
Returns
value nq (<=nprobSum) with the number of quantiles computed

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");
Author
Eddy Offermann

Definition at line 1994 of file TF1.cxx.

◆ GetRandom() [1/2]

Double_t TF1::GetRandom ( Double_t  xmin,
Double_t  xmax,
TRandom rng = nullptr,
Option_t option = nullptr 
)
virtual

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.

Parameters
xminminimum value for generated random numbers
xmaxmaximum value for generated random numbers
rng(optional) random number generator pointer
option(optional) : LOG or LIN to force the usage of a log or linear scale for computing the cumulative integral table

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.

Reimplemented in TF2.

Definition at line 2245 of file TF1.cxx.

◆ GetRandom() [2/2]

Double_t TF1::GetRandom ( TRandom rng = nullptr,
Option_t option = nullptr 
)
virtual

Return a random number following this function shape.

Parameters
rngRandom number generator. By default (or when passing a nullptr) the global gRandom is used
optionOption string which controls the binning used to compute the integral. Default mode is automatic depending of xmax, xmin and Npx (function points). Possible values are:
  • "LOG" to force usage of log scale for tabulating the integral
  • "LIN" to force usage of linear scale when tabulating the integral

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 user can provide as optional parameter a Random number generator. By default gRandom is used

If the ratio fXmax/fXmin > fNpx the integral is tabulated in log scale in x A log scale for the intergral is also always used if a user specifies the "LOG" option Instead if a user requestes a "LIN" option the integral binning is never done in log scale whatever the fXmax/fXmin ratio is

Note that the parabolic approximation is very good as soon as the number of bins is greater than 50.

Reimplemented in TF2.

Definition at line 2192 of file TF1.cxx.

◆ GetRange() [1/4]

void TF1::GetRange ( Double_t xmin,
Double_t xmax 
) const
virtual

Return range of a 1-D function.

Reimplemented in TF2, TF3, and TF3.

Definition at line 2308 of file TF1.cxx.

◆ GetRange() [2/4]

void TF1::GetRange ( Double_t xmin,
Double_t ymin,
Double_t xmax,
Double_t ymax 
) const
virtual

Return range of a 2-D function.

Reimplemented in TF2, TF3, TF2, and TF3.

Definition at line 2318 of file TF1.cxx.

◆ GetRange() [3/4]

void TF1::GetRange ( Double_t xmin,
Double_t ymin,
Double_t zmin,
Double_t xmax,
Double_t ymax,
Double_t zmax 
) const
virtual

Return range of function.

Reimplemented in TF2, TF3, TF2, and TF3.

Definition at line 2330 of file TF1.cxx.

◆ GetRange() [4/4]

void TF1::GetRange ( Double_t xmin,
Double_t xmax 
) const
protectedvirtual

Return range of a generic N-D function.

Reimplemented in TF2, and TF3.

Definition at line 2281 of file TF1.cxx.

◆ GetSave()

Double_t TF1::GetSave ( const Double_t x)
virtual

Get value corresponding to X in array of fSave values.

Reimplemented in TF2, and TF3.

Definition at line 2344 of file TF1.cxx.

◆ GetVariable()

virtual Double_t TF1::GetVariable ( const TString name)
inlinevirtual

Definition at line 569 of file TF1.h.

◆ GetX()

Double_t TF1::GetX ( Double_t  fy,
Double_t  xmin = 0,
Double_t  xmax = 0,
Double_t  epsilon = 1.E-10,
Int_t  maxiter = 100,
Bool_t  logx = false 
) const
virtual

Returns the X value corresponding to the function value fy for (xmin<x<xmax).

in other words it can find the roots of the function when fy=0 and successive calls by changing the next call to [xmin+eps,xmax] where xmin is the previous root.

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 epsilon (default = 1.E-10) controls the relative accuracy (if |x| > 1 ) and absolute (if |x| < 1) and maxiter (default = 100) controls the maximum number of iteration of the Brent algorithm If the flag logx is set the grid search is done in log step size This is done automatically if the log scale is set in the current Pad

NOTE: see also TF1::GetMaximumX, TF1::GetMinimumX

Definition at line 1865 of file TF1.cxx.

◆ GetXaxis()

TAxis * TF1::GetXaxis ( ) const

Get x axis of the function.

Definition at line 2400 of file TF1.cxx.

◆ GetXmax()

virtual Double_t TF1::GetXmax ( ) const
inlinevirtual

Definition at line 562 of file TF1.h.

◆ GetXmin()

virtual Double_t TF1::GetXmin ( ) const
inlinevirtual

Definition at line 558 of file TF1.h.

◆ GetYaxis()

TAxis * TF1::GetYaxis ( ) const

Get y axis of the function.

Definition at line 2411 of file TF1.cxx.

◆ GetZaxis()

TAxis * TF1::GetZaxis ( ) const

Get z axis of the function. (In case this object is a TF2 or TF3)

Definition at line 2422 of file TF1.cxx.

◆ GradientPar() [1/4]

void TF1::GradientPar ( const Double_t x,
Double_t grad,
Double_t  eps = 0.01 
)
virtual

Compute the gradient wrt parameters If the TF1 object is based on a formula expression (TFormula) and TFormula::GenerateGradientPar() has been successfully called automatic differentiation using CLAD is used instead of the default numerical differentiation.

Parameters
xpoint, were the gradient is computed
gradused to return the computed gradient, assumed to be of at least fNpar size
epsif the errors of parameters have been computed, the step used in numerical differentiation is eps*parameter_error.

if the errors have not been computed, step=eps is used default value of eps = 0.01 Method is the same as in Derivative() function

If a parameter is fixed, the gradient on this parameter = 0

Definition at line 2468 of file TF1.cxx.

◆ GradientPar() [2/4]

template<class T >
void TF1::GradientPar ( const T *  x,
T *  grad,
Double_t  eps = 0.01 
)
inline

Definition at line 927 of file TF1.h.

◆ GradientPar() [3/4]

Double_t TF1::GradientPar ( Int_t  ipar,
const Double_t x,
Double_t  eps = 0.01 
)
virtual

Compute the gradient (derivative) wrt a parameter ipar.

Parameters
iparindex of parameter for which the derivative is computed
xpoint, where the derivative is computed
eps- if the errors of parameters have been computed, the step used in numerical differentiation is eps*parameter_error.

if the errors have not been computed, step=eps is used default value of eps = 0.01 Method is the same as in Derivative() function

If a parameter is fixed, the gradient on this parameter = 0

Definition at line 2445 of file TF1.cxx.

◆ GradientPar() [4/4]

template<class T >
T TF1::GradientPar ( Int_t  ipar,
const T *  x,
Double_t  eps = 0.01 
)
inline

Definition at line 860 of file TF1.h.

◆ GradientParTempl() [1/2]

template<class T >
void TF1::GradientParTempl ( const T *  x,
T *  grad,
Double_t  eps = 0.01 
)
inline

Definition at line 936 of file TF1.h.

◆ GradientParTempl() [2/2]

template<class T >
T TF1::GradientParTempl ( Int_t  ipar,
const T *  x,
Double_t  eps = 0.01 
)
inline

Definition at line 869 of file TF1.h.

◆ InitArgs()

void TF1::InitArgs ( const Double_t x,
const Double_t params 
)
virtual

Initialize parameters addresses.

Definition at line 2482 of file TF1.cxx.

◆ InitStandardFunctions()

void TF1::InitStandardFunctions ( )
static

Create the basic function objects.

Definition at line 2497 of file TF1.cxx.

◆ Integral()

Double_t TF1::Integral ( Double_t  a,
Double_t  b,
Double_t  epsrel = 1.e-12 
)
virtual

IntegralOneDim or analytical integral.

Reimplemented in TF2, and TF3.

Definition at line 2531 of file TF1.cxx.

◆ IntegralError() [1/2]

Double_t TF1::IntegralError ( Double_t  a,
Double_t  b,
const Double_t params = nullptr,
const Double_t covmat = nullptr,
Double_t  epsilon = 1.E-2 
)
virtual

Return Error on Integral of a parametric function between a and b due to the parameter uncertainties and their covariance matrix from the fit.

In addition to the integral limits, this method takes as input a pointer to the fitted parameter values and a pointer the covariance matrix from the fit. These pointers should be retrieved from the previously performed fit using the TFitResult class. Note that to get the TFitResult, te fit should be done using the fit option S. Example:

TFitResultPtr r = histo->Fit(func, "S");
func->IntegralError(x1,x2,r->GetParams(), r->GetCovarianceMatrix()->GetMatrixArray() );
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...

IMPORTANT NOTE1:

A null pointer to the parameter values vector and to the covariance matrix can be passed. In this case, when the parameter values pointer is null, the parameter values stored in this TF1 function object are used in the integral error computation. When the poassed pointer to the covariance matrix is null, a covariance matrix from the last fit is retrieved from a global fitter instance when it exists. Note that the global fitter instance esists only when ROOT is not running with multi-threading enabled (ROOT::IsImplicitMTEnabled() == True). When the ovariance matrix from the last fit cannot be retrieved, an error message is printed and a zero value is returned.

IMPORTANT NOTE2:

When no covariance matrix is passed and in the meantime a fit is done using another function, the routine will signal an error and it will return zero only when the number of fit parameter is different than the values stored in TF1 (TF1::GetNpar() ). In the case that npar is the same, an incorrect result is returned.

IMPORTANT NOTE3:

The user must pass a pointer to the elements of the full covariance matrix dimensioned with the right size (npar*npar), where npar is the total number of parameters (TF1::GetNpar()), including also the fixed parameters. The covariance matrix must be retrieved from the TFitResult class as shown above and not from TVirtualFitter::GetCovarianceMatrix() function.

Definition at line 2708 of file TF1.cxx.

◆ IntegralError() [2/2]

Double_t TF1::IntegralError ( Int_t  n,
const Double_t a,
const Double_t b,
const Double_t params = nullptr,
const Double_t covmat = nullptr,
Double_t  epsilon = 1.E-2 
)
virtual

Return Error on Integral of a parametric function with dimension larger than one between a[] and b[] due to the parameters uncertainties.

For a TF1 with dimension larger than 1 (for example a TF2 or TF3) TF1::IntegralMultiple is used for the integral calculation

In addition to the integral limits, this method takes as input a pointer to the fitted parameter values and a pointer the covariance matrix from the fit. These pointers should be retrieved from the previously performed fit using the TFitResult class. Note that to get the TFitResult, te fit should be done using the fit option S. Example:

TFitResultPtr r = histo2d->Fit(func2, "S");
func2->IntegralError(a,b,r->GetParams(), r->GetCovarianceMatrix()->GetMatrixArray() );
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99

IMPORTANT NOTE1:

A null pointer to the parameter values vector and to the covariance matrix can be passed. In this case, when the parameter values pointer is null, the parameter values stored in this TF1 function object are used in the integral error computation. When the poassed pointer to the covariance matrix is null, a covariance matrix from the last fit is retrieved from a global fitter instance when it exists. Note that the global fitter instance esists only when ROOT is not running with multi-threading enabled (ROOT::IsImplicitMTEnabled() == True). When the ovariance matrix from the last fit cannot be retrieved, an error message is printed and a zero value is returned.

IMPORTANT NOTE2:

When no covariance matrix is passed and in the meantime a fit is done using another function, the routine will signal an error and it will return zero only when the number of fit parameter is different than the values stored in TF1 (TF1::GetNpar() ). In the case that npar is the same, an incorrect result is returned.

IMPORTANT NOTE3:

The user must pass a pointer to the elements of the full covariance matrix dimensioned with the right size (npar*npar), where npar is the total number of parameters (TF1::GetNpar()), including also the fixed parameters. The covariance matrix must be retrieved from the TFitResult class as shown above and not from TVirtualFitter::GetCovarianceMatrix() function.

Definition at line 2758 of file TF1.cxx.

◆ IntegralFast()

Double_t TF1::IntegralFast ( Int_t  num,
Double_t x,
Double_t w,
Double_t  a,
Double_t  b,
Double_t params = nullptr,
Double_t  epsilon = 1e-12 
)
virtual

Gauss-Legendre integral, see CalcGaussLegendreSamplingPoints.

Definition at line 2778 of file TF1.cxx.

◆ IntegralMultiple() [1/3]

Double_t TF1::IntegralMultiple ( Int_t  n,
const Double_t a,
const Double_t b,
Double_t  epsrel,
Double_t relerr 
)
virtual

See more general prototype below.

This interface kept for back compatibility It is recommended to use the other interface where one can specify also epsabs and the maximum number of points

Definition at line 2798 of file TF1.cxx.

◆ IntegralMultiple() [2/3]

Double_t TF1::IntegralMultiple ( Int_t  n,
const Double_t a,
const Double_t b,
Int_t  maxpts,
Double_t  epsrel,
Double_t  epsabs,
Double_t relerr,
Int_t nfnevl,
Int_t ifail 
)
virtual

This function computes, to an attempted specified accuracy, the value of the integral.

Parameters
[in]nNumber of dimensions [2,15]
[in]a,bOne-dimensional arrays of length >= N . On entry A[i], and B[i], contain the lower and upper limits of integration, respectively.
[in]maxptsMaximum number of function evaluations to be allowed. maxpts >= 2^n +2*n*(n+1) +1 if maxpts<minpts, maxpts is set to 10*minpts
[in]epsrelSpecified relative accuracy.
[in]epsabsSpecified absolute accuracy. The integration algorithm will attempt to reach either the relative or the absolute accuracy. In case the maximum function called is reached the algorithm will stop earlier without having reached the desired accuracy
[out]relerrContains, on exit, an estimation of the relative accuracy of the result.
[out]nfnevlnumber of function evaluations performed.
[out]ifail

0 Normal exit. At least minpts and at most maxpts calls to the function were performed.

1 maxpts is too small for the specified accuracy eps. The result and relerr contain the values obtainable for the specified value of maxpts.

3 n<2 or n>15

Method:

The default method used is the Genz-Mallik adaptive multidimensional algorithm using the class ROOT::Math::AdaptiveIntegratorMultiDim (see the reference documentation of the class)

Other methods can be used by setting ROOT::Math::IntegratorMultiDimOptions::SetDefaultIntegrator() to different integrators. Other possible integrators are MC integrators based on the ROOT::Math::GSLMCIntegrator class Possible methods are : Vegas, Miser or Plain IN case of MC integration the accuracy is determined by the number of function calls, one should be careful not to use a too large value of maxpts

Definition at line 2851 of file TF1.cxx.

◆ IntegralMultiple() [3/3]

virtual Double_t TF1::IntegralMultiple ( Int_t  n,
const Double_t a,
const Double_t b,
Int_t  ,
Int_t  maxpts,
Double_t  epsrel,
Double_t relerr,
Int_t nfnevl,
Int_t ifail 
)
inlinevirtual

Definition at line 594 of file TF1.h.

◆ IntegralOneDim()

Double_t TF1::IntegralOneDim ( Double_t  a,
Double_t  b,
Double_t  epsrel,
Double_t  epsabs,
Double_t error 
)
virtual

Return Integral of function between a and b using the given parameter values and relative and absolute tolerance.

The default integrator defined in ROOT::Math::IntegratorOneDimOptions::DefaultIntegrator() is used If ROOT contains the MathMore library the default integrator is set to be the adaptive ROOT::Math::GSLIntegrator (based on QUADPACK) or otherwise the ROOT::Math::GaussIntegrator is used See the reference documentation of these classes for more information about the integration algorithms To change integration algorithm just do : ROOT::Math::IntegratorOneDimOptions::SetDefaultIntegrator(IntegratorName); Valid integrator names are:

In order to use the GSL integrators one needs to have the MathMore library installed

Note 1:

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.

Note 2:

Instead of TF1::Integral, you may want to use the combination of TF1::CalcGaussLegendreSamplingPoints and TF1::IntegralFast. See an example with the following script:

void gint() {
TF1 *g = new TF1("g","gaus",-5,5);
g->SetParameters(1,0,1);
//default gaus integration method uses 6 points
//not suitable to integrate on a large domain
double r1 = g->Integral(0,5);
double r2 = g->Integral(0,1000);
//try with user directives computing more points
Int_t np = 1000;
double *x=new double[np];
double *w=new double[np];
g->CalcGaussLegendreSamplingPoints(np,x,w,1e-15);
double r3 = g->IntegralFast(np,x,w,0,5);
double r4 = g->IntegralFast(np,x,w,0,1000);
double r5 = g->IntegralFast(np,x,w,0,10000);
double r6 = g->IntegralFast(np,x,w,0,100000);
printf("g->Integral(0,5) = %g\n",r1);
printf("g->Integral(0,1000) = %g\n",r2);
printf("g->IntegralFast(n,x,w,0,5) = %g\n",r3);
printf("g->IntegralFast(n,x,w,0,1000) = %g\n",r4);
printf("g->IntegralFast(n,x,w,0,10000) = %g\n",r5);
printf("g->IntegralFast(n,x,w,0,100000)= %g\n",r6);
delete [] x;
delete [] w;
}
#define e(i)
Definition RSha256.hxx:103
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np

This example produces the following results:

g->Integral(0,5) = 1.25331
g->Integral(0,1000) = 1.25319
g->IntegralFast(n,x,w,0,5) = 1.25331
g->IntegralFast(n,x,w,0,1000) = 1.25331
g->IntegralFast(n,x,w,0,10000) = 1.25331
g->IntegralFast(n,x,w,0,100000)= 1.253
const Int_t n
Definition legend1.C:16

Definition at line 2621 of file TF1.cxx.

◆ IntegrateForNormalization()

void TF1::IntegrateForNormalization ( )
protected

◆ IsA()

TClass * TF1::IsA ( ) const
inlineoverridevirtual
Returns
TClass describing current object

Reimplemented from TObject.

Reimplemented in TF12, TF2, and TF3.

Definition at line 720 of file TF1.h.

◆ IsEvalNormalized()

virtual Bool_t TF1::IsEvalNormalized ( ) const
inlinevirtual

Definition at line 599 of file TF1.h.

◆ IsInside()

virtual Bool_t TF1::IsInside ( const Double_t x) const
inlinevirtual

return kTRUE if the point is inside the function range

Reimplemented in TF2, and TF3.

Definition at line 604 of file TF1.h.

◆ IsLinear()

virtual Bool_t TF1::IsLinear ( ) const
inlinevirtual

Definition at line 608 of file TF1.h.

◆ IsValid()

Bool_t TF1::IsValid ( ) const
virtual

Return kTRUE if the function is valid.

Definition at line 2882 of file TF1.cxx.

◆ IsVectorized()

bool TF1::IsVectorized ( )
inline

Definition at line 446 of file TF1.h.

◆ Mean()

virtual Double_t TF1::Mean ( Double_t  a,
Double_t  b,
const Double_t params = nullptr,
Double_t  epsilon = 0.000001 
)
inlinevirtual

Definition at line 698 of file TF1.h.

◆ Moment()

Double_t TF1::Moment ( Double_t  n,
Double_t  a,
Double_t  b,
const Double_t params = nullptr,
Double_t  epsilon = 0.000001 
)
virtual

Return nth moment of function between a and b.

See TF1::Integral() for parameter definitions

Definition at line 3667 of file TF1.cxx.

◆ operator()() [1/2]

template<class T >
T TF1::operator() ( const T *  x,
const Double_t params = nullptr 
)
inline

Definition at line 768 of file TF1.h.

◆ operator()() [2/2]

Double_t TF1::operator() ( Double_t  x,
Double_t  y = 0,
Double_t  z = 0,
Double_t  t = 0 
) const
inlinevirtual

Definition at line 762 of file TF1.h.

◆ operator=()

TF1 & TF1::operator= ( const TF1 rhs)

Operator =.

Definition at line 944 of file TF1.cxx.

◆ Paint()

void TF1::Paint ( Option_t option = "")
overridevirtual

Paint this function with its current attributes.

The function is going to be converted in an histogram and the corresponding histogram is painted. The painted histogram can be retrieved calling afterwards the method TF1::GetHistogram()

Reimplemented from TObject.

Reimplemented in TF2, and TF3.

Definition at line 2953 of file TF1.cxx.

◆ Print()

void TF1::Print ( Option_t option = "") const
overridevirtual

This method must be overridden when a class wants to print itself.

Reimplemented from TObject.

Definition at line 2897 of file TF1.cxx.

◆ RejectedPoint()

Bool_t TF1::RejectedPoint ( )
static

See TF1::RejectPoint above.

Definition at line 3657 of file TF1.cxx.

◆ RejectPoint()

void TF1::RejectPoint ( Bool_t  reject = kTRUE)
static

Static function to set the global flag to reject points the fgRejectPoint global flag is tested by all fit functions if TRUE the point is not included in the fit.

This flag can be set by a user in a fitting function. The fgRejectPoint flag is reset by the TH1 and TGraph fitting functions.

Definition at line 3648 of file TF1.cxx.

◆ ReleaseParameter()

void TF1::ReleaseParameter ( Int_t  ipar)
virtual

Release parameter number ipar during a fit operation.

After releasing it, the parameter can vary freely in the fit. The parameter limits are reset to 0,0.

Definition at line 3151 of file TF1.cxx.

◆ Save()

void TF1::Save ( Double_t  xmin,
Double_t  xmax,
Double_t  ymin,
Double_t  ymax,
Double_t  zmin,
Double_t  zmax 
)
virtual

Save values of function in array fSave.

Reimplemented in TF2, and TF3.

Definition at line 3161 of file TF1.cxx.

◆ SavePrimitive()

void TF1::SavePrimitive ( std::ostream &  out,
Option_t option = "" 
)
overridevirtual

Save primitive as a C++ statement(s) on output stream out.

Reimplemented from TObject.

Reimplemented in TF12, TF2, and TF3.

Definition at line 3218 of file TF1.cxx.

◆ SetChisquare()

virtual void TF1::SetChisquare ( Double_t  chi2)
inlinevirtual

Definition at line 618 of file TF1.h.

◆ SetCurrent()

void TF1::SetCurrent ( TF1 f1)
static

Static function setting the current function.

the current function may be accessed in static C-like functions when fitting or painting a function.

Definition at line 3343 of file TF1.cxx.

◆ SetFitResult()

void TF1::SetFitResult ( const ROOT::Fit::FitResult result,
const Int_t indpar = nullptr 
)
virtual

Set the result from the fit parameter values, errors, chi2, etc... Optionally a pointer to a vector (with size fNpar) of the parameter indices in the FitResult can be passed This is useful in the case of a combined fit with different functions, and the FitResult contains the global result By default it is assume that indpar = {0,1,2,....,fNpar-1}.

Definition at line 3355 of file TF1.cxx.

◆ SetFunction() [1/2]

template<typename Func >
void TF1::SetFunction ( Func  f)

Definition at line 845 of file TF1.h.

◆ SetFunction() [2/2]

template<class PtrObj , typename MemFn >
void TF1::SetFunction ( PtrObj &  p,
MemFn  memFn 
)

Definition at line 852 of file TF1.h.

◆ SetMaximum()

void TF1::SetMaximum ( Double_t  maximum = -1111)
virtual

Set the maximum value along Y for this function In case the function is already drawn, set also the maximum in the helper histogram.

Definition at line 3394 of file TF1.cxx.

◆ SetMinimum()

void TF1::SetMinimum ( Double_t  minimum = -1111)
virtual

Set the minimum value along Y for this function In case the function is already drawn, set also the minimum in the helper histogram.

Definition at line 3407 of file TF1.cxx.

◆ SetNDF()

void TF1::SetNDF ( Int_t  ndf)
virtual

Set the number of degrees of freedom ndf should be the number of points used in a fit - the number of free parameters.

Definition at line 3419 of file TF1.cxx.

◆ SetNormalized()

virtual void TF1::SetNormalized ( Bool_t  flag)
inlinevirtual

Definition at line 634 of file TF1.h.

◆ SetNpx()

void TF1::SetNpx ( Int_t  npx = 100)
virtual

Set the number of points used to draw the function.

The default number of points along x is 100 for 1-d functions and 30 for 2-d/3-d functions You can increase this value to get a better resolution when drawing pictures with sharp peaks or to get a better result when using TF1::GetRandom the minimum number of points is 4, the maximum is 10000000 for 1-d and 10000 for 2-d/3-d functions

Definition at line 3433 of file TF1.cxx.

◆ SetNumberFitPoints()

virtual void TF1::SetNumberFitPoints ( Int_t  npfits)
inlinevirtual

Definition at line 630 of file TF1.h.

◆ SetParameter() [1/2]

virtual void TF1::SetParameter ( const TString name,
Double_t  value 
)
inlinevirtual

Definition at line 645 of file TF1.h.

◆ SetParameter() [2/2]

virtual void TF1::SetParameter ( Int_t  param,
Double_t  value 
)
inlinevirtual

Definition at line 640 of file TF1.h.

◆ SetParameters() [1/2]

virtual void TF1::SetParameters ( const Double_t params)
inlinevirtual

Definition at line 650 of file TF1.h.

◆ SetParameters() [2/2]

virtual void TF1::SetParameters ( double  p0,
double  p1 = 0.0,
double  p2 = 0.0,
double  p3 = 0.0,
double  p4 = 0.0,
double  p5 = 0.0,
double  p6 = 0.0,
double  p7 = 0.0,
double  p8 = 0.0,
double  p9 = 0.0,
double  p10 = 0.0 
)
inlinevirtual

Definition at line 655 of file TF1.h.

◆ SetParent()

virtual void TF1::SetParent ( TObject p = nullptr)
inlinevirtual

Definition at line 671 of file TF1.h.

◆ SetParError()

void TF1::SetParError ( Int_t  ipar,
Double_t  error 
)
virtual

Set error for parameter number ipar.

Definition at line 3473 of file TF1.cxx.

◆ SetParErrors()

void TF1::SetParErrors ( const Double_t errors)
virtual

Set errors for all active parameters when calling this function, the array errors must have at least fNpar values.

Definition at line 3484 of file TF1.cxx.

◆ SetParLimits()

void TF1::SetParLimits ( Int_t  ipar,
Double_t  parmin,
Double_t  parmax 
)
virtual

Set lower and upper limits for parameter ipar.

The specified limits will be used in a fit operation. Note that when this function is a pre-defined function (e.g. gaus) one needs to use the fit option "B" to have the limits used in the fit. See TH1::Fit(TF1*, Option_t *, Option_t *, Double_t, Double_t) for the fitting documentation and the fitting options

To fix a parameter, use TF1::FixParameter

Definition at line 3501 of file TF1.cxx.

◆ SetParName()

void TF1::SetParName ( Int_t  ipar,
const char *  name 
)
virtual

Set name of parameter number ipar.

Definition at line 3450 of file TF1.cxx.

◆ SetParNames()

void TF1::SetParNames ( const char *  name0 = "p0",
const char *  name1 = "p1",
const char *  name2 = "p2",
const char *  name3 = "p3",
const char *  name4 = "p4",
const char *  name5 = "p5",
const char *  name6 = "p6",
const char *  name7 = "p7",
const char *  name8 = "p8",
const char *  name9 = "p9",
const char *  name10 = "p10" 
)
virtual

Set up to 10 parameter names.

Definition at line 3462 of file TF1.cxx.

◆ SetRange() [1/3]

void TF1::SetRange ( Double_t  xmin,
Double_t  xmax 
)
virtual

Initialize the upper and lower bounds to draw the function.

The function range is also used in an histogram fit operation when the option "R" is specified.

Reimplemented in TF2, and TF3.

Definition at line 3522 of file TF1.cxx.

◆ SetRange() [2/3]

void TF1::SetRange ( Double_t  xmin,
Double_t  ymin,
Double_t  xmax,
Double_t  ymax 
)
inlinevirtual

Reimplemented in TF2, and TF3.

Definition at line 835 of file TF1.h.

◆ SetRange() [3/3]

void TF1::SetRange ( Double_t  xmin,
Double_t  ymin,
Double_t  zmin,
Double_t  xmax,
Double_t  ymax,
Double_t  zmax 
)
inlinevirtual

Reimplemented in TF2, and TF3.

Definition at line 839 of file TF1.h.

◆ SetSavedPoint()

void TF1::SetSavedPoint ( Int_t  point,
Double_t  value 
)
virtual

Restore value of function saved at point.

Definition at line 3536 of file TF1.cxx.

◆ SetTitle()

void TF1::SetTitle ( const char *  title = "")
overridevirtual

Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title is "fffffff" and "xxxx" and "yyyy" are the titles for the X and Y axis respectively.

Reimplemented from TNamed.

Definition at line 3552 of file TF1.cxx.

◆ SetVectorized()

virtual void TF1::SetVectorized ( Bool_t  vectorized)
inlinevirtual

Definition at line 680 of file TF1.h.

◆ Streamer()

void TF1::Streamer ( TBuffer b)
overridevirtual

Stream a class object.

Reimplemented from TObject.

Reimplemented in TF12, TF2, and TF3.

Definition at line 3565 of file TF1.cxx.

◆ StreamerNVirtual()

void TF1::StreamerNVirtual ( TBuffer ClassDef_StreamerNVirtual_b)
inline

Definition at line 720 of file TF1.h.

◆ TermCoeffLength()

int TF1::TermCoeffLength ( TString term)
private

Definition at line 926 of file TF1.cxx.

◆ Update()

void TF1::Update ( )
virtual

Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated histogram or Integral.

Definition at line 3613 of file TF1.cxx.

◆ Variance()

virtual Double_t TF1::Variance ( Double_t  a,
Double_t  b,
const Double_t params = nullptr,
Double_t  epsilon = 0.000001 
)
inlinevirtual

Definition at line 702 of file TF1.h.

Friends And Related Symbol Documentation

◆ ROOT::Internal::TF1Builder

template<class Func >
friend struct ROOT::Internal::TF1Builder
friend

Definition at line 217 of file TF1.h.

Member Data Documentation

◆ fAlpha

std::vector<Double_t> TF1::fAlpha
protected

! Array alpha. for each bin in x the deconvolution r of fIntegral

Definition at line 260 of file TF1.h.

◆ fBeta

std::vector<Double_t> TF1::fBeta
protected

! Array beta. is approximated by x = alpha +beta*r *gamma*r**2

Definition at line 261 of file TF1.h.

◆ fChisquare

Double_t TF1::fChisquare {}
protected

Function fit chisquare.

Definition at line 252 of file TF1.h.

◆ fComposition

std::unique_ptr<TF1AbsComposition> TF1::fComposition
protected

Pointer to composition (NSUM or CONV)

Definition at line 271 of file TF1.h.

◆ fFormula

std::unique_ptr<TFormula> TF1::fFormula
protected

Pointer to TFormula in case when user define formula.

Definition at line 269 of file TF1.h.

◆ fFunctor

std::unique_ptr<TF1FunctorPointer> TF1::fFunctor
protected

! Functor object to wrap any C++ callable object

Definition at line 268 of file TF1.h.

◆ fgAbsValue

std::atomic< Bool_t > TF1::fgAbsValue
staticprotected

Definition at line 305 of file TF1.h.

◆ fgAddToGlobList

std::atomic< Bool_t > TF1::fgAddToGlobList
staticprotected

Definition at line 307 of file TF1.h.

◆ fGamma

std::vector<Double_t> TF1::fGamma
protected

! Array gamma.

Definition at line 262 of file TF1.h.

◆ fgCurrent

TF1 * TF1::fgCurrent = nullptr
staticprotected

Definition at line 308 of file TF1.h.

◆ fgRejectPoint

Bool_t TF1::fgRejectPoint = kFALSE
staticprotected

Definition at line 306 of file TF1.h.

◆ fHistogram

TH1* TF1::fHistogram {nullptr}
protected

! Pointer to histogram used for visualisation

Definition at line 264 of file TF1.h.

◆ fIntegral

std::vector<Double_t> TF1::fIntegral
protected

! Integral of function binned on fNpx bins

Definition at line 259 of file TF1.h.

◆ fMaximum

Double_t TF1::fMaximum {-1111}
protected

Maximum value for plotting.

Definition at line 254 of file TF1.h.

◆ fMethodCall

std::unique_ptr<TMethodCall> TF1::fMethodCall
protected

! Pointer to MethodCall in case of interpreted function

Definition at line 265 of file TF1.h.

◆ fMinimum

Double_t TF1::fMinimum {-1111}
protected

Minimum value for plotting.

Definition at line 253 of file TF1.h.

◆ fNDF

Int_t TF1::fNDF {}
protected

Number of degrees of freedom in the fit.

Definition at line 251 of file TF1.h.

◆ fNdim

Int_t TF1::fNdim {}
protected

Function dimension.

Definition at line 247 of file TF1.h.

◆ fNormalized

Bool_t TF1::fNormalized {false}
protected

Normalization option (false by default)

Definition at line 266 of file TF1.h.

◆ fNormIntegral

Double_t TF1::fNormIntegral {}
protected

Integral of the function before being normalized.

Definition at line 267 of file TF1.h.

◆ fNpar

Int_t TF1::fNpar {}
protected

Number of parameters.

Definition at line 246 of file TF1.h.

◆ fNpfits

Int_t TF1::fNpfits {}
protected

Number of points used in the fit.

Definition at line 250 of file TF1.h.

◆ fNpx

Int_t TF1::fNpx {100}
protected

Number of points used for the graphical representation.

Definition at line 248 of file TF1.h.

◆ fParams

std::unique_ptr<TF1Parameters> TF1::fParams
protected

Pointer to Function parameters object (exists only for not-formula functions)

Definition at line 270 of file TF1.h.

◆ fParent

TObject* TF1::fParent {nullptr}
protected

! Parent object hooking this function (if one)

Definition at line 263 of file TF1.h.

◆ fParErrors

std::vector<Double_t> TF1::fParErrors
protected

Array of errors of the fNpar parameters.

Definition at line 255 of file TF1.h.

◆ fParMax

std::vector<Double_t> TF1::fParMax
protected

Array of upper limits of the fNpar parameters.

Definition at line 257 of file TF1.h.

◆ fParMin

std::vector<Double_t> TF1::fParMin
protected

Array of lower limits of the fNpar parameters.

Definition at line 256 of file TF1.h.

◆ fSave

std::vector<Double_t> TF1::fSave
protected

Array of fNsave function values.

Definition at line 258 of file TF1.h.

◆ fType

EFType TF1::fType {EFType::kTemplScalar}
protected

Definition at line 249 of file TF1.h.

◆ fXmax

Double_t TF1::fXmax {-1111}
protected

Upper bounds for the range.

Definition at line 245 of file TF1.h.

◆ fXmin

Double_t TF1::fXmin {-1111}
protected

Lower bounds for the range.

Definition at line 244 of file TF1.h.

Libraries for TF1:

The documentation for this class was generated from the following files: