Logo ROOT  
Reference Guide
TF1.h
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Rene Brun 18/08/95
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 // ---------------------------------- F1.h
12 
13 #ifndef ROOT_TF1
14 #define ROOT_TF1
15 
16 //////////////////////////////////////////////////////////////////////////
17 // //
18 // TF1 //
19 // //
20 // The Parametric 1-D function //
21 // //
22 //////////////////////////////////////////////////////////////////////////
23 
24 #include "RConfigure.h"
25 #include <functional>
26 #include <cassert>
27 #include <string>
28 #include <vector>
29 #include "TFormula.h"
30 #include "TAttLine.h"
31 #include "TAttFill.h"
32 #include "TAttMarker.h"
33 #include "TF1AbsComposition.h"
34 #include "TMath.h"
35 #include "Math/Types.h"
36 #include "Math/ParamFunctor.h"
37 
38 class TF1;
39 class TH1;
40 class TAxis;
41 class TMethodCall;
42 
43 namespace ROOT {
44  namespace Fit {
45  class FitResult;
46  }
47 }
48 
50 public:
51  TF1Parameters() {} // needed for the I/O
53  fParameters(std::vector<Double_t>(npar)),
54  fParNames(std::vector<std::string>(npar))
55  {
56  for (int i = 0; i < npar; ++i) {
57  fParNames[i] = std::string(TString::Format("p%d", i).Data());
58  }
59  }
60  // copy constructor
63  fParNames(rhs.fParNames)
64  {}
65  // assignment
67  {
68  if (&rhs == this) return *this;
70  fParNames = rhs.fParNames;
71  return *this;
72  }
73  virtual ~TF1Parameters() {}
74 
75  // getter methods
76  Double_t GetParameter(Int_t iparam) const
77  {
78  return (CheckIndex(iparam)) ? fParameters[iparam] : 0;
79  }
80  Double_t GetParameter(const char *name) const
81  {
83  }
84  const Double_t *GetParameters() const
85  {
86  return fParameters.data();
87  }
88  const std::vector<double> &ParamsVec() const
89  {
90  return fParameters;
91  }
92 
93  Int_t GetParNumber(const char *name) const;
94 
95  const char *GetParName(Int_t iparam) const
96  {
97  return (CheckIndex(iparam)) ? fParNames[iparam].c_str() : "";
98  }
99 
100 
101  // setter methods
102  void SetParameter(Int_t iparam, Double_t value)
103  {
104  if (!CheckIndex(iparam)) return;
105  fParameters[iparam] = value;
106  }
107  void SetParameters(const Double_t *params)
108  {
109  std::copy(params, params + fParameters.size(), fParameters.begin());
110  }
111  void SetParameters(Double_t p0, Double_t p1, Double_t p2 = 0, Double_t p3 = 0, Double_t p4 = 0,
112  Double_t p5 = 0, Double_t p6 = 0, Double_t p7 = 0, Double_t p8 = 0,
113  Double_t p9 = 0, Double_t p10 = 0);
114 
115  void SetParameter(const char *name, Double_t value)
116  {
117  SetParameter(GetParNumber(name), value);
118  }
119  void SetParName(Int_t iparam, const char *name)
120  {
121  if (!CheckIndex(iparam)) return;
122  fParNames[iparam] = std::string(name);
123  }
124  void SetParNames(const char *name0 = "p0", const char *name1 = "p1", const char *name2 = "p2",
125  const char *name3 = "p3", const char *name4 = "p4", const char *name5 = "p5",
126  const char *name6 = "p6", const char *name7 = "p7", const char *name8 = "p8",
127  const char *name9 = "p9", const char *name10 = "p10");
128 
129 
130 
131  ClassDef(TF1Parameters, 1) // The Parameters of a parameteric function
132 private:
133 
134  bool CheckIndex(Int_t i) const
135  {
136  return (i >= 0 && i < int(fParameters.size()));
137  }
138 
139  std::vector<Double_t> fParameters; // parameter values
140  std::vector<std::string> fParNames; // parameter names
141 };
142 
143 namespace ROOT {
144  namespace Internal {
145  /// Internal class used by TF1 for defining
146  /// template specialization for different TF1 constructors
147  template<class Func>
148  struct TF1Builder {
149  static void Build(TF1 *f, Func func);
150  };
151 
152  template<class Func>
153  struct TF1Builder<Func *> {
154  static void Build(TF1 *f, Func *func);
155  };
156 
157  // Internal class used by TF1 for obtaining the type from a functor
158  // out of the set of valid operator() signatures.
159  template<typename T>
160  struct GetFunctorType {
161  };
162 
163  template<typename F, typename T>
164  struct GetFunctorType<T(F::*)(const T *, const double *)> {
165  using type = T;
166  };
167 
168  template<typename F, typename T>
169  struct GetFunctorType<T(F::*)(const T *, const double *) const> {
170  using type = T;
171  };
172 
173  template<typename F, typename T>
174  struct GetFunctorType<T(F::*)(T *, double *)> {
175  using type = T;
176  };
177 
178  template<typename F, typename T>
179  struct GetFunctorType<T(F::*)(T *, double *) const> {
180  using type = T;
181  };
182 
183  // Internal class used by TF1 to get the right operator() signature
184  // from a Functor with several ones.
185  template<typename T, typename F>
186  auto GetTheRightOp(T(F::*opPtr)(const T *, const double *)) -> decltype(opPtr)
187  {
188  return opPtr;
189  }
190 
191  template<typename T, typename F>
192  auto GetTheRightOp(T(F::*opPtr)(const T *, const double *) const) -> decltype(opPtr)
193  {
194  return opPtr;
195  }
196 
197  template<typename T, typename F>
198  auto GetTheRightOp(T(F::*opPtr)(T *, double *)) -> decltype(opPtr)
199  {
200  return opPtr;
201  }
202 
203  template<typename T, typename F>
204  auto GetTheRightOp(T(F::*opPtr)(T *, double *) const) -> decltype(opPtr)
205  {
206  return opPtr;
207  }
208  }
209 }
210 
211 
212 class TF1 : public TNamed, public TAttLine, public TAttFill, public TAttMarker {
213 
214  template<class Func>
216 
217 public:
218  // Add to list behavior
219  enum class EAddToList {
220  kDefault,
221  kAdd,
222  kNo
223  };
224 
225 protected:
226 
228  virtual ~TF1FunctorPointer() {}
229  virtual TF1FunctorPointer * Clone() const = 0;
230  };
231 
232 
233  enum EFType {
234  kFormula = 0, // formula functions which can be stored,
235  kPtrScalarFreeFcn, // pointer to scalar free function,
236  kInterpreted, // interpreted functions constructed by name,
237  kTemplVec, // vectorized free functions or TemplScalar functors evaluating on vectorized parameters,
238  kTemplScalar, // TemplScalar functors evaluating on scalar parameters
240  }; // formula based on composition class (e.g. NSUM, CONV)
241 
242  Double_t fXmin{-1111}; //Lower bounds for the range
243  Double_t fXmax{-1111}; //Upper bounds for the range
244  Int_t fNpar{}; //Number of parameters
245  Int_t fNdim{}; //Function dimension
246  Int_t fNpx{100}; //Number of points used for the graphical representation
247  EFType fType{EFType::kTemplScalar};
248  Int_t fNpfits{}; //Number of points used in the fit
249  Int_t fNDF{}; //Number of degrees of freedom in the fit
250  Double_t fChisquare{}; //Function fit chisquare
251  Double_t fMinimum{-1111}; //Minimum value for plotting
252  Double_t fMaximum{-1111}; //Maximum value for plotting
253  std::vector<Double_t> fParErrors; //Array of errors of the fNpar parameters
254  std::vector<Double_t> fParMin; //Array of lower limits of the fNpar parameters
255  std::vector<Double_t> fParMax; //Array of upper limits of the fNpar parameters
256  std::vector<Double_t> fSave; //Array of fNsave function values
257  std::vector<Double_t> fIntegral; //!Integral of function binned on fNpx bins
258  std::vector<Double_t> fAlpha; //!Array alpha. for each bin in x the deconvolution r of fIntegral
259  std::vector<Double_t> fBeta; //!Array beta. is approximated by x = alpha +beta*r *gamma*r**2
260  std::vector<Double_t> fGamma; //!Array gamma.
261  TObject *fParent{nullptr}; //!Parent object hooking this function (if one)
262  TH1 *fHistogram{nullptr}; //!Pointer to histogram used for visualisation
263  TMethodCall *fMethodCall{nullptr}; //!Pointer to MethodCall in case of interpreted function
264  Bool_t fNormalized{false}; //Normalization option (false by default)
265  Double_t fNormIntegral{}; //Integral of the function before being normalized
266  TF1FunctorPointer *fFunctor{nullptr}; //! Functor object to wrap any C++ callable object
267  TFormula *fFormula{nullptr}; //Pointer to TFormula in case when user define formula
268  TF1Parameters *fParams{nullptr}; //Pointer to Function parameters object (exists only for not-formula functions)
269  std::unique_ptr<TF1AbsComposition> fComposition; //! Pointer to composition (NSUM or CONV)
270  TF1AbsComposition *fComposition_ptr{nullptr}; // saved pointer (unique_ptr is transient)
271 
272  /// General constructor for TF1. Most of the other constructors delegate on it
273  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):
274  TNamed(name, name), TAttLine(), TAttFill(), TAttMarker(), fXmin(xmin), fXmax(xmax), fNpar(npar), fNdim(ndim),
275  fType(functionType), fParErrors(npar), fParMin(npar), fParMax(npar), fFunctor(functor), fParams(params)
276  {
277  DoInitialize(addToGlobList);
278  };
279 
280 private:
281  // NSUM parsing helper functions
282  void DefineNSUMTerm(TObjArray *newFuncs, TObjArray *coeffNames,
283  TString &fullFormula,
284  TString &formula, int termStart, int termEnd,
286  int TermCoeffLength(TString &term);
287 
288 protected:
289 
290  template <class T>
293  TF1FunctorPointerImpl(const std::function<T(const T *f, const Double_t *param)> &func) : fImpl(func){};
295  virtual TF1FunctorPointer * Clone() const { return new TF1FunctorPointerImpl<T>(fImpl); }
297  };
298 
299 
300 
301 
302  static std::atomic<Bool_t> fgAbsValue; //use absolute value of function when computing integral
303  static Bool_t fgRejectPoint; //True if point must be rejected in a fit
304  static std::atomic<Bool_t> fgAddToGlobList; //True if we want to register the function in the global list
305  static TF1 *fgCurrent; //pointer to current function being processed
306 
307 
308  //void CreateFromFunctor(const char *name, Int_t npar, Int_t ndim = 1);
309  void DoInitialize(EAddToList addToGlobList);
310 
312 
313  virtual Double_t GetMinMaxNDim(Double_t *x , Bool_t findmax, Double_t epsilon = 0, Int_t maxiter = 0) const;
314  virtual void GetRange(Double_t *xmin, Double_t *xmax) const;
315  virtual TH1 *DoCreateHistogram(Double_t xmin, Double_t xmax, Bool_t recreate = kFALSE);
316 
317 public:
318 
319  // TF1 status bits
320  enum EStatusBits {
321  kNotGlobal = BIT(10), // don't register in global list of functions
322  kNotDraw = BIT(9) // don't draw the function when in a TH1
323  };
324 
325  TF1();
326  TF1(const char *name, const char *formula, Double_t xmin = 0, Double_t xmax = 1, EAddToList addToGlobList = EAddToList::kDefault, bool vectorize = false);
327  TF1(const char *name, const char *formula, Double_t xmin, Double_t xmax, Option_t * option); // same as above but using a string for option
328  TF1(const char *name, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault);
329  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);
330  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);
331 
332  template <class T>
333  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):
334  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
335  {
336  fType = std::is_same<T, double>::value ? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
337  }
338 
339  ////////////////////////////////////////////////////////////////////////////////
340  /// Constructor using a pointer to function.
341  ///
342  /// \param npar is the number of free parameters used by the function
343  ///
344  /// This constructor creates a function of type C when invoked
345  /// with the normal C++ compiler.
346  ///
347  ///
348  /// WARNING! A function created with this constructor cannot be Cloned
349 
350 
351  template <class T>
352  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):
353  TF1(EFType::kTemplVec, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
354  {}
355 
356  // Constructors using functors (compiled mode only)
357  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);
358 
359  // Template constructors from any C++ callable object, defining the operator() (double * , double *)
360  // and returning a double.
361  // The class name is not needed when using compile code, while it is required when using
362  // interpreted code via the specialized constructor with void *.
363  // An instance of the C++ function class or its pointer can both be used. The former is reccomended when using
364  // C++ compiled code, but if CINT compatibility is needed, then a pointer to the function class must be used.
365  // xmin and xmax specify the plotting range, npar is the number of parameters.
366  // See the tutorial math/exampleFunctor.C for an example of using this constructor
367  template <typename Func>
368  TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault) :
369  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList)
370  {
371  //actual fType set in TF1Builder
373  }
374 
375  // backward compatible interface
376  template <typename Func>
377  TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *, EAddToList addToGlobList = EAddToList::kDefault) :
378  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar))
379  {
381  }
382 
383 
384  // Template constructors from a pointer to any C++ class of type PtrObj with a specific member function of type
385  // MemFn.
386  // The member function must have the signature of (double * , double *) and returning a double.
387  // The class name and the method name are not needed when using compile code
388  // (the member function pointer is used in this case), while they are required when using interpreted
389  // code via the specialized constructor with void *.
390  // xmin and xmax specify the plotting range, npar is the number of parameters.
391  // See the tutorial math/exampleFunctor.C for an example of using this constructor
392  template <class PtrObj, typename MemFn>
393  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) :
394  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
395  {}
396 
397  // backward compatible interface
398  template <class PtrObj, typename MemFn>
399  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) :
400  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
401  {}
402 
403  TF1(const TF1 &f1);
404  TF1 &operator=(const TF1 &rhs);
405  virtual ~TF1();
406  virtual void AddParameter(const TString &name, Double_t value)
407  {
408  if (fFormula) fFormula->AddParameter(name, value);
409  }
410  // virtual void AddParameters(const pair<TString,Double_t> *pairs, Int_t size) { fFormula->AddParameters(pairs,size); }
411  // virtual void AddVariable(const TString &name, Double_t value = 0) { if (fFormula) fFormula->AddVariable(name,value); }
412  // virtual void AddVariables(const TString *vars, Int_t size) { if (fFormula) fFormula->AddVariables(vars,size); }
413  virtual Bool_t AddToGlobalList(Bool_t on = kTRUE);
415  virtual void Browse(TBrowser *b);
416  virtual void Copy(TObject &f1) const;
417  TObject* Clone(const char* newname=0) const;
418  virtual Double_t Derivative(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
419  virtual Double_t Derivative2(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
420  virtual Double_t Derivative3(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
421  static Double_t DerivativeError();
422  virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
423  virtual void Draw(Option_t *option = "");
424  virtual TF1 *DrawCopy(Option_t *option = "") const;
425  virtual TObject *DrawDerivative(Option_t *option = "al"); // *MENU*
426  virtual TObject *DrawIntegral(Option_t *option = "al"); // *MENU*
427  virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option = "");
428  virtual Double_t Eval(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
429  //template <class T> T Eval(T x, T y = 0, T z = 0, T t = 0) const;
430  virtual Double_t EvalPar(const Double_t *x, const Double_t *params = 0);
431  template <class T> T EvalPar(const T *x, const Double_t *params = 0);
432  virtual Double_t operator()(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
433  template <class T> T operator()(const T *x, const Double_t *params = nullptr);
434  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
435  virtual void FixParameter(Int_t ipar, Double_t value);
437  {
438  return (fType == EFType::kTemplVec) || (fType == EFType::kFormula && fFormula && fFormula->IsVectorized());
439  }
441  {
442  return fChisquare;
443  }
444  virtual TH1 *GetHistogram() const;
445  virtual TH1 *CreateHistogram()
446  {
447  return DoCreateHistogram(fXmin, fXmax);
448  }
449  virtual TFormula *GetFormula()
450  {
451  return fFormula;
452  }
453  virtual const TFormula *GetFormula() const
454  {
455  return fFormula;
456  }
457  virtual TString GetExpFormula(Option_t *option = "") const
458  {
459  return (fFormula) ? fFormula->GetExpFormula(option) : "";
460  }
461  virtual const TObject *GetLinearPart(Int_t i) const
462  {
463  return (fFormula) ? fFormula->GetLinearPart(i) : nullptr;
464  }
465  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;
466  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;
467  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;
468  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;
469  virtual Double_t GetMaximumStored() const
470  {
471  return fMaximum;
472  }
473  virtual Double_t GetMinimumStored() const
474  {
475  return fMinimum;
476  }
477  virtual Int_t GetNpar() const
478  {
479  return fNpar;
480  }
481  virtual Int_t GetNdim() const
482  {
483  return fNdim;
484  }
485  virtual Int_t GetNDF() const;
486  virtual Int_t GetNpx() const
487  {
488  return fNpx;
489  }
491  {
492  return fMethodCall;
493  }
494  virtual Int_t GetNumber() const
495  {
496  return (fFormula) ? fFormula->GetNumber() : 0;
497  }
498  virtual Int_t GetNumberFreeParameters() const;
499  virtual Int_t GetNumberFitPoints() const
500  {
501  return fNpfits;
502  }
503  virtual char *GetObjectInfo(Int_t px, Int_t py) const;
505  {
506  return fParent;
507  }
508  virtual Double_t GetParameter(Int_t ipar) const
509  {
510  return (fFormula) ? fFormula->GetParameter(ipar) : fParams->GetParameter(ipar);
511  }
512  virtual Double_t GetParameter(const TString &name) const
513  {
515  }
516  virtual Double_t *GetParameters() const
517  {
518  return (fFormula) ? fFormula->GetParameters() : const_cast<Double_t *>(fParams->GetParameters());
519  }
520  virtual void GetParameters(Double_t *params)
521  {
522  if (fFormula) fFormula->GetParameters(params);
523  else std::copy(fParams->ParamsVec().begin(), fParams->ParamsVec().end(), params);
524  }
525  virtual const char *GetParName(Int_t ipar) const
526  {
527  return (fFormula) ? fFormula->GetParName(ipar) : fParams->GetParName(ipar);
528  }
529  virtual Int_t GetParNumber(const char *name) const
530  {
532  }
533  virtual Double_t GetParError(Int_t ipar) const;
534  virtual const Double_t *GetParErrors() const
535  {
536  return fParErrors.data();
537  }
538  virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const;
539  virtual Double_t GetProb() const;
540  virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum);
541  virtual Double_t GetRandom();
543  virtual void GetRange(Double_t &xmin, Double_t &xmax) const;
544  virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const;
545  virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const;
546  virtual Double_t GetSave(const Double_t *x);
547  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;
548  virtual Double_t GetXmin() const
549  {
550  return fXmin;
551  }
552  virtual Double_t GetXmax() const
553  {
554  return fXmax;
555  }
556  TAxis *GetXaxis() const ;
557  TAxis *GetYaxis() const ;
558  TAxis *GetZaxis() const ;
560  {
561  return (fFormula) ? fFormula->GetVariable(name) : 0;
562  }
563  virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps = 0.01);
564  template <class T>
565  T GradientPar(Int_t ipar, const T *x, Double_t eps = 0.01);
566  template <class T>
567  T GradientParTempl(Int_t ipar, const T *x, Double_t eps = 0.01);
568 
569  virtual void GradientPar(const Double_t *x, Double_t *grad, Double_t eps = 0.01);
570  template <class T>
571  void GradientPar(const T *x, T *grad, Double_t eps = 0.01);
572  template <class T>
573  void GradientParTempl(const T *x, T *grad, Double_t eps = 0.01);
574 
575  virtual void InitArgs(const Double_t *x, const Double_t *params);
576  static void InitStandardFunctions();
577  virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel = 1.e-12);
578  virtual Double_t IntegralOneDim(Double_t a, Double_t b, Double_t epsrel, Double_t epsabs, Double_t &err);
579  virtual Double_t IntegralError(Double_t a, Double_t b, const Double_t *params = 0, const Double_t *covmat = 0, Double_t epsilon = 1.E-2);
580  virtual Double_t IntegralError(Int_t n, const Double_t *a, const Double_t *b, const Double_t *params = 0, const Double_t *covmat = 0, Double_t epsilon = 1.E-2);
581  // virtual Double_t IntegralFast(const TGraph *g, Double_t a, Double_t b, Double_t *params=0);
582  virtual Double_t IntegralFast(Int_t num, Double_t *x, Double_t *w, Double_t a, Double_t b, Double_t *params = 0, Double_t epsilon = 1e-12);
583  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);
584  virtual Double_t IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Int_t /*minpts*/, Int_t maxpts, Double_t epsrel, Double_t &relerr, Int_t &nfnevl, Int_t &ifail)
585  {
586  return IntegralMultiple(n, a, b, maxpts, epsrel, epsrel, relerr, nfnevl, ifail);
587  }
588  virtual Double_t IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Double_t epsrel, Double_t &relerr);
589  virtual Bool_t IsEvalNormalized() const
590  {
591  return fNormalized;
592  }
593  /// return kTRUE if the point is inside the function range
594  virtual Bool_t IsInside(const Double_t *x) const
595  {
596  return !((x[0] < fXmin) || (x[0] > fXmax));
597  }
598  virtual Bool_t IsLinear() const
599  {
600  return (fFormula) ? fFormula->IsLinear() : false;
601  }
602  virtual Bool_t IsValid() const;
603  virtual void Print(Option_t *option = "") const;
604  virtual void Paint(Option_t *option = "");
605  virtual void ReleaseParameter(Int_t ipar);
606  virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax);
607  virtual void SavePrimitive(std::ostream &out, Option_t *option = "");
608  virtual void SetChisquare(Double_t chi2)
609  {
610  fChisquare = chi2;
611  }
612  virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar = 0);
613  template <class PtrObj, typename MemFn>
614  void SetFunction(PtrObj &p, MemFn memFn);
615  template <typename Func>
616  void SetFunction(Func f);
617  virtual void SetMaximum(Double_t maximum = -1111); // *MENU*
618  virtual void SetMinimum(Double_t minimum = -1111); // *MENU*
619  virtual void SetNDF(Int_t ndf);
620  virtual void SetNumberFitPoints(Int_t npfits)
621  {
622  fNpfits = npfits;
623  }
624  virtual void SetNormalized(Bool_t flag)
625  {
626  fNormalized = flag;
627  Update();
628  }
629  virtual void SetNpx(Int_t npx = 100); // *MENU*
630  virtual void SetParameter(Int_t param, Double_t value)
631  {
632  (fFormula) ? fFormula->SetParameter(param, value) : fParams->SetParameter(param, value);
633  Update();
634  }
635  virtual void SetParameter(const TString &name, Double_t value)
636  {
637  (fFormula) ? fFormula->SetParameter(name, value) : fParams->SetParameter(name, value);
638  Update();
639  }
640  virtual void SetParameters(const Double_t *params)
641  {
642  (fFormula) ? fFormula->SetParameters(params) : fParams->SetParameters(params);
643  Update();
644  }
645  virtual void SetParameters(Double_t p0, Double_t p1, Double_t p2 = 0, Double_t p3 = 0, Double_t p4 = 0,
646  Double_t p5 = 0, Double_t p6 = 0, Double_t p7 = 0, Double_t p8 = 0,
647  Double_t p9 = 0, Double_t p10 = 0)
648  {
649  if (fFormula) fFormula->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
650  else fParams->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
651  Update();
652  } // *MENU*
653  virtual void SetParName(Int_t ipar, const char *name);
654  virtual void SetParNames(const char *name0 = "p0", const char *name1 = "p1", const char *name2 = "p2",
655  const char *name3 = "p3", const char *name4 = "p4", const char *name5 = "p5",
656  const char *name6 = "p6", const char *name7 = "p7", const char *name8 = "p8",
657  const char *name9 = "p9", const char *name10 = "p10"); // *MENU*
658  virtual void SetParError(Int_t ipar, Double_t error);
659  virtual void SetParErrors(const Double_t *errors);
660  virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax);
661  virtual void SetParent(TObject *p = 0)
662  {
663  fParent = p;
664  }
665  virtual void SetRange(Double_t xmin, Double_t xmax); // *MENU*
668  virtual void SetSavedPoint(Int_t point, Double_t value);
669  virtual void SetTitle(const char *title = ""); // *MENU*
670  virtual void SetVectorized(Bool_t vectorized)
671  {
672  if (fType == EFType::kFormula && fFormula)
673  fFormula->SetVectorized(vectorized);
674  else
675  Warning("SetVectorized", "Can only set vectorized flag on formula-based TF1");
676  }
677  virtual void Update();
678 
679  static TF1 *GetCurrent();
680  static void AbsValue(Bool_t reject = kTRUE);
681  static void RejectPoint(Bool_t reject = kTRUE);
682  static Bool_t RejectedPoint();
683  static void SetCurrent(TF1 *f1);
684 
685  //Moments
686  virtual Double_t Moment(Double_t n, Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001);
687  virtual Double_t CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001);
688  virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001)
689  {
690  return Moment(1, a, b, params, epsilon);
691  }
692  virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001)
693  {
694  return CentralMoment(2, a, b, params, epsilon);
695  }
696 
697  //some useful static utility functions to compute sampling points for Integral
698  //static void CalcGaussLegendreSamplingPoints(TGraph *g, Double_t eps=3.0e-11);
699  //static TGraph *CalcGaussLegendreSamplingPoints(Int_t num=21, Double_t eps=3.0e-11);
700  static void CalcGaussLegendreSamplingPoints(Int_t num, Double_t *x, Double_t *w, Double_t eps = 3.0e-11);
701 
702 private:
703  template <class T>
704  T EvalParTempl(const T *data, const Double_t *params = 0);
705 
706 #ifdef R__HAS_VECCORE
707  inline double EvalParVec(const Double_t *data, const Double_t *params);
708 #endif
709 
710  ClassDef(TF1, 10) // The Parametric 1-D function
711 };
712 
713 namespace ROOT {
714  namespace Internal {
715 
716  template<class Func>
717  void TF1Builder<Func>::Build(TF1 *f, Func func)
718  {
719  using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
720  f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
722  f->fParams = new TF1Parameters(f->fNpar);
723  }
724 
725  template<class Func>
726  void TF1Builder<Func *>::Build(TF1 *f, Func *func)
727  {
728  using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
729  f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
731  f->fParams = new TF1Parameters(f->fNpar);
732  }
733 
734  /// TF1 building from a string
735  /// used to build a TFormula based on a lambda function
736  template<>
737  struct TF1Builder<const char *> {
738  static void Build(TF1 *f, const char *formula)
739  {
740  f->fType = TF1::EFType::kFormula;
741  f->fFormula = new TFormula("tf1lambda", formula, f->fNdim, f->fNpar, false);
742  TString formulaExpression(formula);
743  Ssiz_t first = formulaExpression.Index("return") + 7;
744  Ssiz_t last = formulaExpression.Last(';');
745  TString title = formulaExpression(first, last - first);
746  f->SetTitle(title);
747  }
748  };
749  }
750 }
751 
753 {
754  return Eval(x, y, z, t);
755 }
756 
757 template <class T>
758 inline T TF1::operator()(const T *x, const Double_t *params)
759 {
760  return EvalPar(x, params);
761 }
762 
763 ////////////////////////////////////////////////////////////////////////////////
764 /// EvalPar for vectorized
765 template <class T>
766 T TF1::EvalPar(const T *x, const Double_t *params)
767 {
768  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
769  return EvalParTempl(x, params);
770  } else if (fType == EFType::kFormula) {
771  return fFormula->EvalPar(x, params);
772  } else
773  return TF1::EvalPar((double *)x, params);
774 }
775 
776 ////////////////////////////////////////////////////////////////////////////////
777 /// Eval for vectorized functions
778 // template <class T>
779 // T TF1::Eval(T x, T y, T z, T t) const
780 // {
781 // if (fType == EFType::kFormula)
782 // return fFormula->Eval(x, y, z, t);
783 
784 // T xx[] = {x, y, z, t};
785 // Double_t *pp = (Double_t *)fParams->GetParameters();
786 // return ((TF1 *)this)->EvalPar(xx, pp);
787 // }
788 
789 // Internal to TF1. Evaluates Templated interfaces
790 template <class T>
791 inline T TF1::EvalParTempl(const T *data, const Double_t *params)
792 {
793  assert(fType == EFType::kTemplScalar || fType == EFType::kTemplVec);
794  if (!params) params = (Double_t *)fParams->GetParameters();
795  if (fFunctor)
796  return ((TF1FunctorPointerImpl<T> *)fFunctor)->fImpl(data, params);
797 
798  // this should throw an error
799  // we nned to implement a vectorized GetSave(x)
800  return TMath::SignalingNaN();
801 }
802 
803 #ifdef R__HAS_VECCORE
804 // Internal to TF1. Evaluates Vectorized TF1 on data of type Double_v
805 inline double TF1::EvalParVec(const Double_t *data, const Double_t *params)
806 {
807  assert(fType == EFType::kTemplVec);
808  std::vector<ROOT::Double_v> d(fNdim);
809  ROOT::Double_v res;
810 
811  for(auto i=0; i<fNdim; i++) {
812  d[i] = ROOT::Double_v(data[i]);
813  }
814 
815  if (fFunctor) {
816  res = ((TF1FunctorPointerImpl<ROOT::Double_v> *) fFunctor)->fImpl(d.data(), params);
817  } else {
818  // res = GetSave(x);
819  return TMath::SignalingNaN();
820  }
821  return vecCore::Get<ROOT::Double_v>(res, 0);
822 }
823 #endif
824 
826 {
828 }
830 {
832 }
833 
834 template <typename Func>
835 void TF1::SetFunction(Func f)
836 {
837  // set function from a generic C++ callable object
838  fType = EFType::kPtrScalarFreeFcn;
840 }
841 template <class PtrObj, typename MemFn>
842 void TF1::SetFunction(PtrObj &p, MemFn memFn)
843 {
844  // set from a pointer to a member function
845  fType = EFType::kPtrScalarFreeFcn;
847 }
848 
849 template <class T>
850 inline T TF1::GradientPar(Int_t ipar, const T *x, Double_t eps)
851 {
852  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
853  return GradientParTempl<T>(ipar, x, eps);
854  } else
855  return GradientParTempl<Double_t>(ipar, (const Double_t *)x, eps);
856 }
857 
858 template <class T>
859 inline T TF1::GradientParTempl(Int_t ipar, const T *x, Double_t eps)
860 {
861  if (GetNpar() == 0)
862  return 0;
863 
864  if (eps < 1e-10 || eps > 1) {
865  Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
866  eps = 0.01;
867  }
868  Double_t h;
869  TF1 *func = (TF1 *)this;
870  Double_t *parameters = GetParameters();
871 
872  // Copy parameters for thread safety
873  std::vector<Double_t> parametersCopy(parameters, parameters + GetNpar());
874  parameters = parametersCopy.data();
875 
876  Double_t al, bl, h2;
877  T f1, f2, g1, g2, d0, d2;
878 
879  ((TF1 *)this)->GetParLimits(ipar, al, bl);
880  if (al * bl != 0 && al >= bl) {
881  // this parameter is fixed
882  return 0;
883  }
884 
885  // check if error has been computer (is not zero)
886  if (func->GetParError(ipar) != 0)
887  h = eps * func->GetParError(ipar);
888  else
889  h = eps;
890 
891  // save original parameters
892  Double_t par0 = parameters[ipar];
893 
894  parameters[ipar] = par0 + h;
895  f1 = func->EvalPar(x, parameters);
896  parameters[ipar] = par0 - h;
897  f2 = func->EvalPar(x, parameters);
898  parameters[ipar] = par0 + h / 2;
899  g1 = func->EvalPar(x, parameters);
900  parameters[ipar] = par0 - h / 2;
901  g2 = func->EvalPar(x, parameters);
902 
903  // compute the central differences
904  h2 = 1 / (2. * h);
905  d0 = f1 - f2;
906  d2 = 2 * (g1 - g2);
907 
908  T grad = h2 * (4 * d2 - d0) / 3.;
909 
910  // restore original value
911  parameters[ipar] = par0;
912 
913  return grad;
914 }
915 
916 template <class T>
917 inline void TF1::GradientPar(const T *x, T *grad, Double_t eps)
918 {
919  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
920  GradientParTempl<T>(x, grad, eps);
921  } else
922  GradientParTempl<Double_t>((const Double_t *)x, (Double_t *)grad, eps);
923 }
924 
925 template <class T>
926 inline void TF1::GradientParTempl(const T *x, T *grad, Double_t eps)
927 {
928  if (eps < 1e-10 || eps > 1) {
929  Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
930  eps = 0.01;
931  }
932 
933  for (Int_t ipar = 0; ipar < GetNpar(); ipar++) {
934  grad[ipar] = GradientParTempl<T>(ipar, x, eps);
935  }
936 }
937 
938 #endif
TF1::TF1FunctorPointerImpl::fImpl
ROOT::Math::ParamFunctorTempl< T > fImpl
Definition: TF1.h:296
TF1::GetMinimumX
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,...
Definition: TF1.cxx:1824
TF1::fNDF
Int_t fNDF
Definition: TF1.h:249
TF1::TF1
TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
Definition: TF1.h:368
TF1::GetHistogram
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function.
Definition: TF1.cxx:1585
TF1::DrawF1
virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option="")
Draw function between xmin and xmax.
Definition: TF1.cxx:1431
TF1::EAddToList::kDefault
@ kDefault
n
const Int_t n
Definition: legend1.C:16
TF1::EAddToList::kAdd
@ kAdd
TAxis
Class to manage histogram axis.
Definition: TAxis.h:30
TF1::GetCurrent
static TF1 * GetCurrent()
Static function returning the current function being processed.
Definition: TF1.cxx:1575
TF1::fgAbsValue
static std::atomic< Bool_t > fgAbsValue
Definition: TF1.h:302
first
Definition: first.py:1
TF1Parameters::SetParName
void SetParName(Int_t iparam, const char *name)
Definition: TF1.h:119
TF1::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)
Definition: TF1.h:399
TBrowser
Definition: TBrowser.h:37
TF1::IntegralMultiple
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)
Definition: TF1.h:584
ymax
float ymax
Definition: THbookFile.cxx:95
TF1::RejectPoint
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition: TF1.cxx:3672
TF1::ExecuteEvent
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TF1.cxx:1547
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:91
TFormula::IsLinear
Bool_t IsLinear() const
Definition: TFormula.h:237
e
#define e(i)
Definition: RSha256.hxx:121
HFit::Fit
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:133
TF1::GetObjectInfo
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Redefines TObject::GetObjectInfo.
Definition: TF1.cxx:1919
TF1::kFormula
@ kFormula
Definition: TF1.h:234
TFormula::IsVectorized
Bool_t IsVectorized() const
Definition: TFormula.h:236
TF1::fSave
std::vector< Double_t > fSave
Definition: TF1.h:256
TObjArray
Definition: TObjArray.h:37
TF1::fHistogram
TH1 * fHistogram
Parent object hooking this function (if one)
Definition: TF1.h:262
TF1::GetParErrors
virtual const Double_t * GetParErrors() const
Definition: TF1.h:534
f
#define f(i)
Definition: RSha256.hxx:122
TF1Parameters::SetParameters
void SetParameters(const Double_t *params)
Definition: TF1.h:107
TF1::fNpar
Int_t fNpar
Definition: TF1.h:244
TF1::GetMaximum
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.
Definition: TF1.cxx:1615
TFormula::SetVectorized
void SetVectorized(Bool_t vectorized)
TF1Parameters::GetParameter
Double_t GetParameter(const char *name) const
Definition: TF1.h:80
TF1::TF1FunctorPointerImpl::TF1FunctorPointerImpl
TF1FunctorPointerImpl(const std::function< T(const T *f, const Double_t *param)> &func)
Definition: TF1.h:293
Types.h
TF1::Variance
virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params=0, Double_t epsilon=0.000001)
Definition: TF1.h:692
TF1::GetX
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).
Definition: TF1.cxx:1864
TFormula
Definition: TFormula.h:85
TF1::kPtrScalarFreeFcn
@ kPtrScalarFreeFcn
Definition: TF1.h:235
ROOT::Internal::GetFunctorType< T(F::*)(const T *, const double *) const >::type
T type
Definition: TF1.h:170
TF1::fNormalized
Bool_t fNormalized
Pointer to MethodCall in case of interpreted function.
Definition: TF1.h:264
TF1Parameters::~TF1Parameters
virtual ~TF1Parameters()
Definition: TF1.h:73
TF1::fgCurrent
static TF1 * fgCurrent
Definition: TF1.h:305
TF1::SetNpx
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3454
TF1Parameters::ParamsVec
const std::vector< double > & ParamsVec() const
Definition: TF1.h:88
TF1::GetMaximumStored
virtual Double_t GetMaximumStored() const
Definition: TF1.h:469
TF1::Integral
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2523
TF1::GetParameters
virtual Double_t * GetParameters() const
Definition: TF1.h:516
TF1::SetParameters
virtual void SetParameters(Double_t p0, Double_t p1, Double_t p2=0, Double_t p3=0, Double_t p4=0, Double_t p5=0, Double_t p6=0, Double_t p7=0, Double_t p8=0, Double_t p9=0, Double_t p10=0)
Definition: TF1.h:645
ROOT::Internal::GetFunctorType< T(F::*)(T *, double *)>::type
T type
Definition: TF1.h:175
TF1::EvalParTempl
T EvalParTempl(const T *data, const Double_t *params=0)
Eval for vectorized functions.
Definition: TF1.h:791
TF1::SetParameter
virtual void SetParameter(const TString &name, Double_t value)
Definition: TF1.h:635
F
#define F(x, y, z)
TF1::GetParent
TObject * GetParent() const
Definition: TF1.h:504
TF1::fNormIntegral
Double_t fNormIntegral
Definition: TF1.h:265
TFormula::GetParameters
Double_t * GetParameters() const
TF1::GetNdim
virtual Int_t GetNdim() const
Definition: TF1.h:481
TF1::CreateHistogram
virtual TH1 * CreateHistogram()
Definition: TF1.h:445
TF1::GetNpar
virtual Int_t GetNpar() const
Definition: TF1.h:477
TF1::SetParNames
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:3483
TF1::SetMaximum
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum value along Y for this function In case the function is already drawn,...
Definition: TF1.cxx:3415
TF1Parameters::GetParameters
const Double_t * GetParameters() const
Definition: TF1.h:84
xmax
float xmax
Definition: THbookFile.cxx:95
TF1::~TF1
virtual ~TF1()
TF1 default destructor.
Definition: TF1.cxx:942
TF1AbsComposition
Definition: TF1AbsComposition.h:16
TF1::GetProb
virtual Double_t GetProb() const
Return the fit probability.
Definition: TF1.cxx:1956
TF1::kTemplVec
@ kTemplVec
Definition: TF1.h:237
TF1::Derivative3
virtual Double_t Derivative3(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the third derivative of the function at point x, computed by Richardson's extrapolation metho...
Definition: TF1.cxx:1248
TF1::Moment
virtual Double_t Moment(Double_t n, Double_t a, Double_t b, const Double_t *params=0, Double_t epsilon=0.000001)
Return nth moment of function between a and b.
Definition: TF1.cxx:3691
TF1::fType
EFType fType
Definition: TF1.h:247
TF1::fIntegral
std::vector< Double_t > fIntegral
Definition: TF1.h:257
TAttMarker.h
TF1::SetCurrent
static void SetCurrent(TF1 *f1)
Static function setting the current function.
Definition: TF1.cxx:3364
TF1::kInterpreted
@ kInterpreted
Definition: TF1.h:236
TF1::GetNumberFitPoints
virtual Int_t GetNumberFitPoints() const
Definition: TF1.h:499
Int_t
int Int_t
Definition: RtypesCore.h:45
TF1::GetXmax
virtual Double_t GetXmax() const
Definition: TF1.h:552
TF1::GetFormula
virtual const TFormula * GetFormula() const
Definition: TF1.h:453
TF1::DistancetoPrimitive
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a function.
Definition: TF1.cxx:1298
TF1Parameters::GetParName
const char * GetParName(Int_t iparam) const
Definition: TF1.h:95
TF1::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)
General constructor for TF1. Most of the other constructors delegate on it.
Definition: TF1.h:273
TF1::TF1FunctorPointer
Definition: TF1.h:227
TF1::SavePrimitive
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TF1.cxx:3208
x
Double_t x[n]
Definition: legend1.C:17
TF1::fXmin
Double_t fXmin
Definition: TF1.h:242
TF1::AddParameter
virtual void AddParameter(const TString &name, Double_t value)
Definition: TF1.h:406
TF1::GetMinMaxNDim
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.
Definition: TF1.cxx:1724
TF1::SetParErrors
virtual void SetParErrors(const Double_t *errors)
Set errors for all active parameters when calling this function, the array errors must have at least ...
Definition: TF1.cxx:3505
TFormula::GetNumber
Int_t GetNumber() const
Definition: TFormula.h:225
TF1Parameters::TF1Parameters
TF1Parameters(const TF1Parameters &rhs)
Definition: TF1.h:61
TF1::EAddToList::kNo
@ kNo
TF1::IsInside
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition: TF1.h:594
TF1::kNotDraw
@ kNotDraw
Definition: TF1.h:322
TF1::DefineNSUMTerm
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.
Definition: TF1.cxx:871
TF1::fMethodCall
TMethodCall * fMethodCall
Pointer to histogram used for visualisation.
Definition: TF1.h:263
TF1::CalcGaussLegendreSamplingPoints
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.
Definition: TF1.cxx:3815
TString::Format
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition: TString.cxx:2311
TF1::SetMinimum
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum value along Y for this function In case the function is already drawn,...
Definition: TF1.cxx:3428
TF1::kNotGlobal
@ kNotGlobal
Definition: TF1.h:321
TString
Definition: TString.h:136
TF1::IsVectorized
bool IsVectorized()
Definition: TF1.h:436
TF1::DerivativeError
static Double_t DerivativeError()
Static function returning the error of the last call to the of Derivative's functions.
Definition: TF1.cxx:1282
TF1Parameters::TF1Parameters
TF1Parameters(Int_t npar)
Definition: TF1.h:52
Bool_t
bool Bool_t
Definition: RtypesCore.h:63
TF1::Paint
virtual void Paint(Option_t *option="")
Paint this function with its current attributes.
Definition: TF1.cxx:2945
b
#define b(i)
Definition: RSha256.hxx:118
TF1Parameters::operator=
TF1Parameters & operator=(const TF1Parameters &rhs)
Definition: TF1.h:66
TF1::operator()
virtual Double_t operator()(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Definition: TF1.h:752
bool
TF1::IntegralOneDim
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 absolut...
Definition: TF1.cxx:2613
TF1::TF1
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.
Definition: TF1.h:352
TF1::GetParLimits
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1941
TF1::SetNDF
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...
Definition: TF1.cxx:3440
TF1::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)
Definition: TF1.h:333
q
float * q
Definition: THbookFile.cxx:89
TString::Last
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:892
TF1::GetRange
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition: TF1.cxx:2284
TF1::GetNumber
virtual Int_t GetNumber() const
Definition: TF1.h:494
TF1::Clone
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TF1.cxx:1069
TF1::SetParameter
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:630
TF1::SetParameters
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:640
TF1::GradientPar
virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps=0.01)
Compute the gradient (derivative) wrt a parameter ipar.
Definition: TF1.cxx:2448
TF1::kTemplScalar
@ kTemplScalar
Definition: TF1.h:238
TF1::fFormula
TFormula * fFormula
Functor object to wrap any C++ callable object.
Definition: TF1.h:267
TF1::IntegralFast
virtual Double_t IntegralFast(Int_t num, Double_t *x, Double_t *w, Double_t a, Double_t b, Double_t *params=0, Double_t epsilon=1e-12)
Gauss-Legendre integral, see CalcGaussLegendreSamplingPoints.
Definition: TF1.cxx:2770
TF1::SetFitResult
virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar=0)
Set the result from the fit parameter values, errors, chi2, etc...
Definition: TF1.cxx:3376
TAttLine.h
TF1::IsValid
virtual Bool_t IsValid() const
Return kTRUE if the function is valid.
Definition: TF1.cxx:2874
TF1Parameters::GetParameter
Double_t GetParameter(Int_t iparam) const
Definition: TF1.h:76
TF1::GetParNumber
virtual Int_t GetParNumber(const char *name) const
Definition: TF1.h:529
TF1::fNpfits
Int_t fNpfits
Definition: TF1.h:248
TF1::fNdim
Int_t fNdim
Definition: TF1.h:245
TFormula::EvalPar
Double_t EvalPar(const Double_t *x, const Double_t *params=0) const
TF1Parameters::GetParNumber
Int_t GetParNumber(const char *name) const
Returns the parameter number given a name not very efficient but list of parameters is typically smal...
Definition: TF1.cxx:3835
TF1::SetParError
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition: TF1.cxx:3494
TF1::GetYaxis
TAxis * GetYaxis() const
Get y axis of the function.
Definition: TF1.cxx:2414
TF1::fXmax
Double_t fXmax
Definition: TF1.h:243
TF1::fAlpha
std::vector< Double_t > fAlpha
Integral of function binned on fNpx bins.
Definition: TF1.h:258
TF1::GetExpFormula
virtual TString GetExpFormula(Option_t *option="") const
Definition: TF1.h:457
Option_t
const typedef char Option_t
Definition: RtypesCore.h:66
TFormula::GetExpFormula
TString GetExpFormula(Option_t *option="") const
TFormula::GetParName
const char * GetParName(Int_t ipar) const
TF1::TF1FunctorPointer::~TF1FunctorPointer
virtual ~TF1FunctorPointer()
Definition: TF1.h:228
TAttFill.h
TAttLine
Definition: TAttLine.h:18
ROOT::Internal::TF1Builder::Build
static void Build(TF1 *f, Func func)
Definition: TF1.h:717
TF1::InitStandardFunctions
static void InitStandardFunctions()
Create the basic function objects.
Definition: TF1.cxx:2490
TF1::SetFunction
void SetFunction(PtrObj &p, MemFn memFn)
Definition: TF1.h:842
TF1::GetMethodCall
TMethodCall * GetMethodCall() const
Definition: TF1.h:490
TF1::DrawDerivative
virtual TObject * DrawDerivative(Option_t *option="al")
Draw derivative of this function.
Definition: TF1.cxx:1390
xmin
float xmin
Definition: THbookFile.cxx:95
TF1::Update
virtual void Update()
Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated h...
Definition: TF1.cxx:3637
h
#define h(i)
Definition: RSha256.hxx:124
TF1::SetNormalized
virtual void SetNormalized(Bool_t flag)
Definition: TF1.h:624
TF1::GetXaxis
TAxis * GetXaxis() const
Get x axis of the function.
Definition: TF1.cxx:2403
TF1::InitArgs
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition: TF1.cxx:2475
epsilon
REAL epsilon
Definition: triangle.c:617
ROOT::Internal::GetFunctorType
Definition: TF1.h:160
a
auto * a
Definition: textangle.C:12
TF1::fMinimum
Double_t fMinimum
Definition: TF1.h:251
TNamed
Definition: TNamed.h:29
TF1::Derivative2
virtual Double_t Derivative2(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the second derivative of the function at point x, computed by Richardson's extrapolation meth...
Definition: TF1.cxx:1183
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:92
TF1::GetParName
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:525
TF1::GetNpx
virtual Int_t GetNpx() const
Definition: TF1.h:486
ROOT::R::function
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
TF1::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)
Definition: TF1.h:393
TF1::TF1FunctorPointer::Clone
virtual TF1FunctorPointer * Clone() const =0
TF1::CentralMoment
virtual Double_t CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params=0, 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)
Definition: TF1.cxx:3728
TF1::DoCreateHistogram
virtual TH1 * DoCreateHistogram(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 w...
Definition: TF1.cxx:3044
TF1::Browse
virtual void Browse(TBrowser *b)
Browse.
Definition: TF1.cxx:986
ROOT::Fit::FitResult
class containg the result of the fit and all the related information (fitted parameter values,...
Definition: FitResult.h:52
TF1::SetChisquare
virtual void SetChisquare(Double_t chi2)
Definition: TF1.h:608
TF1::Copy
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition: TF1.cxx:998
TF1::EFType
EFType
Definition: TF1.h:233
TAttMarker
Definition: TAttMarker.h:19
TF1::fgRejectPoint
static Bool_t fgRejectPoint
Definition: TF1.h:303
ROOT::Internal::TF1Builder< const char * >::Build
static void Build(TF1 *f, const char *formula)
Definition: TF1.h:738
ROOT::Internal::GetFunctorType< T(F::*)(const T *, const double *)>::type
T type
Definition: TF1.h:165
BIT
#define BIT(n)
Definition: Rtypes.h:85
TF1::IsEvalNormalized
virtual Bool_t IsEvalNormalized() const
Definition: TF1.h:589
double
double
Definition: Converters.cxx:921
y
Double_t y[n]
Definition: legend1.C:17
TFormula::GetVariable
Double_t GetVariable(const char *name) const
ROOT::Double_v
Double_t Double_v
Definition: Types.h:51
TF1::fFunctor
TF1FunctorPointer * fFunctor
Definition: TF1.h:266
TF1::TF1
TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *, EAddToList addToGlobList=EAddToList::kDefault)
Definition: TF1.h:377
TF1AbsComposition.h
TF1::GetZaxis
TAxis * GetZaxis() const
Get z axis of the function. (In case this object is a TF2 or TF3)
Definition: TF1.cxx:2425
TF1::GetLinearPart
virtual const TObject * GetLinearPart(Int_t i) const
Definition: TF1.h:461
TF1::kCompositionFcn
@ kCompositionFcn
Definition: TF1.h:239
TMath::SignalingNaN
Double_t SignalingNaN()
Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN)
Definition: TMath.h:897
TObject::Warning
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:876
TF1::fParErrors
std::vector< Double_t > fParErrors
Definition: TF1.h:253
TF1::fMaximum
Double_t fMaximum
Definition: TF1.h:252
TF1::SetNumberFitPoints
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:620
TF1::TF1FunctorPointerImpl::Clone
virtual TF1FunctorPointer * Clone() const
Definition: TF1.h:295
TF1::fBeta
std::vector< Double_t > fBeta
Array alpha. for each bin in x the deconvolution r of fIntegral.
Definition: TF1.h:259
TF1::SetRange
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF1.cxx:3540
TF1::GetSave
virtual Double_t GetSave(const Double_t *x)
Get value corresponding to X in array of fSave values.
Definition: TF1.cxx:2347
ROOT::Internal::GetFunctorType< T(F::*)(T *, double *) const >::type
T type
Definition: TF1.h:180
TF1::EvalPar
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1479
TF1Parameters::CheckIndex
bool CheckIndex(Int_t i) const
Definition: TF1.h:134
TF1Parameters::SetParNames
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 parameter names.
Definition: TF1.cxx:3867
TF1::fChisquare
Double_t fChisquare
Definition: TF1.h:250
TFormula::GetLinearPart
const TObject * GetLinearPart(Int_t i) const
ymin
float ymin
Definition: THbookFile.cxx:95
TString::Index
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:639
TF1::fParMin
std::vector< Double_t > fParMin
Definition: TF1.h:254
TF1::EAddToList
EAddToList
Definition: TF1.h:219
TF1Parameters::fParNames
std::vector< std::string > fParNames
Definition: TF1.h:140
TF1::GetMinimumStored
virtual Double_t GetMinimumStored() const
Definition: TF1.h:473
TF1::AddToGlobalList
virtual Bool_t AddToGlobalList(Bool_t on=kTRUE)
Add to global list of functions (gROOT->GetListOfFunctions() ) return previous status (true if the fu...
Definition: TF1.cxx:835
TF1::Mean
virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params=0, Double_t epsilon=0.000001)
Definition: TF1.h:688
f1
TF1 * f1
Definition: legend1.C:11
TF1::SetSavedPoint
virtual void SetSavedPoint(Int_t point, Double_t value)
Restore value of function saved at point.
Definition: TF1.cxx:3554
Double_t
double Double_t
Definition: RtypesCore.h:59
TF1::GetParError
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1931
TF1::SetTitle
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition: TF1.cxx:3570
TF1::TF1FunctorPointerImpl
Definition: TF1.h:291
TFormula::GetParNumber
Int_t GetParNumber(const char *name) const
TF1::DrawCopy
virtual TF1 * DrawCopy(Option_t *option="") const
Draw a copy of this function with its current attributes.
Definition: TF1.cxx:1368
TF1::TF1
TF1()
TF1 default constructor.
Definition: TF1.cxx:494
TF1::RejectedPoint
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3681
TF1::Eval
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1450
TFormula::AddParameter
void AddParameter(const TString &name, Double_t value=0)
Definition: TFormula.h:186
TF1::GetMaximumX
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.
Definition: TF1.cxx:1656
TFormula::SetParameters
void SetParameters(const Double_t *params)
TF1::Draw
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1338
TF1::TF1FunctorPointerImpl::TF1FunctorPointerImpl
TF1FunctorPointerImpl(const ROOT::Math::ParamFunctorTempl< T > &func)
Definition: TF1.h:292
TF1::GetMinimum
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.
Definition: TF1.cxx:1697
TObject
Definition: TObject.h:37
TF1::fgAddToGlobList
static std::atomic< Bool_t > fgAddToGlobList
Definition: TF1.h:304
ClassDef
#define ClassDef(name, id)
Definition: Rtypes.h:325
TF1::AbsValue
static void AbsValue(Bool_t reject=kTRUE)
Static function: set the fgAbsValue flag.
Definition: TF1.cxx:977
TH1
Definition: TH1.h:57
TF1::Derivative
virtual Double_t Derivative(Double_t x, Double_t *params=0, Double_t epsilon=0.001) const
Returns the first derivative of the function at point x, computed by Richardson's extrapolation metho...
Definition: TF1.cxx:1118
TF1Parameters
Definition: TF1.h:49
TF1::fParent
TObject * fParent
Array gamma.
Definition: TF1.h:261
TFormula::SetParameter
void SetParameter(const char *name, Double_t value)
name
char name[80]
Definition: TGX11.cxx:110
TF1::SetParName
virtual void SetParName(Int_t ipar, const char *name)
Set name of parameter number ipar.
Definition: TF1.cxx:3471
ROOT::Internal::GetTheRightOp
auto GetTheRightOp(T(F::*opPtr)(const T *, const double *)) -> decltype(opPtr)
Definition: TF1.h:186
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:52
TF1::operator=
TF1 & operator=(const TF1 &rhs)
Operator =.
Definition: TF1.cxx:930
d
#define d(i)
Definition: RSha256.hxx:120
TF1Parameters::TF1Parameters
TF1Parameters()
Definition: TF1.h:51
TF1::GetNDF
virtual Int_t GetNDF() const
Return the number of degrees of freedom in the fit the fNDF parameter has been previously computed du...
Definition: TF1.cxx:1890
ROOT::Math::ParamFunctor
ParamFunctorTempl< double > ParamFunctor
Definition: ParamFunctor.h:397
TF1Parameters::fParameters
std::vector< Double_t > fParameters
Definition: TF1.h:139
TF1::DoInitialize
void DoInitialize(EAddToList addToGlobList)
Common initialization of the TF1.
Definition: TF1.cxx:790
TF1::DefaultAddToGlobalList
static Bool_t DefaultAddToGlobalList(Bool_t on=kTRUE)
Static method to add/avoid to add automatically functions to the global list (gROOT->GetListOfFunctio...
Definition: TF1.cxx:826
TF1::fComposition_ptr
TF1AbsComposition * fComposition_ptr
Pointer to composition (NSUM or CONV)
Definition: TF1.h:270
TF1::fGamma
std::vector< Double_t > fGamma
Array beta. is approximated by x = alpha +beta*r *gamma*r**2.
Definition: TF1.h:260
TF1::GetNumberFreeParameters
virtual Int_t GetNumberFreeParameters() const
Return the number of free parameters.
Definition: TF1.cxx:1901
TF1::GetVariable
virtual Double_t GetVariable(const TString &name)
Definition: TF1.h:559
TMethodCall
Definition: TMethodCall.h:37
TAttFill
Definition: TAttFill.h:19
TF1::SetParLimits
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set limits for parameter ipar.
Definition: TF1.cxx:3519
TF1::IntegralMultiple
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.
Definition: TF1.cxx:2843
TF1::GetParameter
virtual Double_t GetParameter(const TString &name) const
Definition: TF1.h:512
TF1::SetParent
virtual void SetParent(TObject *p=0)
Definition: TF1.h:661
ROOT::Internal::TF1Builder
Internal class used by TF1 for defining template specialization for different TF1 constructors.
Definition: TF1.h:148
TF1::IntegrateForNormalization
void IntegrateForNormalization()
type
int type
Definition: TGX11.cxx:121
TF1
1-Dim function class
Definition: TF1.h:212
ParamFunctor.h
TF1::TermCoeffLength
int TermCoeffLength(TString &term)
Definition: TF1.cxx:912
TF1::GetQuantiles
virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum)
Compute Quantiles for density distribution of this function.
Definition: TF1.cxx:1994
TF1::DrawIntegral
virtual TObject * DrawIntegral(Option_t *option="al")
Draw integral of this function.
Definition: TF1.cxx:1415
TF1Parameters::SetParameter
void SetParameter(const char *name, Double_t value)
Definition: TF1.h:115
ROOT::Math::ParamFunctorTempl
Param Functor class for Multidimensional functions.
Definition: ParamFunctor.h:283
TF1::TF1FunctorPointerImpl::~TF1FunctorPointerImpl
virtual ~TF1FunctorPointerImpl()
Definition: TF1.h:294
TF1::GetParameter
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:508
TF1::Save
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.
Definition: TF1.cxx:3155
TF1::GetParameters
virtual void GetParameters(Double_t *params)
Definition: TF1.h:520
TF1::GetRandom
virtual Double_t GetRandom()
Return a random number following this function shape.
Definition: TF1.cxx:2095
TF1::fNpx
Int_t fNpx
Definition: TF1.h:246
TF1::fComposition
std::unique_ptr< TF1AbsComposition > fComposition
Definition: TF1.h:269
TF1::SetVectorized
virtual void SetVectorized(Bool_t vectorized)
Definition: TF1.h:670
TF1::GetChisquare
Double_t GetChisquare() const
Definition: TF1.h:440
ROOT
VSD Structures.
Definition: StringConv.hxx:21
TMath::E
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:102
TF1::FixParameter
virtual void FixParameter(Int_t ipar, Double_t value)
Fix the value of a parameter The specified value will be used in a fit operation.
Definition: TF1.cxx:1563
TF1Parameters::SetParameter
void SetParameter(Int_t iparam, Double_t value)
Definition: TF1.h:102
TF1::IsLinear
virtual Bool_t IsLinear() const
Definition: TF1.h:598
Math
TObject::EStatusBits
EStatusBits
Definition: TObject.h:57
TF1::IntegralError
virtual Double_t IntegralError(Double_t a, Double_t b, const Double_t *params=0, const Double_t *covmat=0, Double_t epsilon=1.E-2)
Return Error on Integral of a parametric function between a and b due to the parameter uncertainties ...
Definition: TF1.cxx:2700
TMath.h
TF1::fParams
TF1Parameters * fParams
Definition: TF1.h:268
int
TF1::Print
virtual void Print(Option_t *option="") const
Print TNamed name and title.
Definition: TF1.cxx:2889
TF1::GradientParTempl
T GradientParTempl(Int_t ipar, const T *x, Double_t eps=0.01)
Definition: TF1.h:859
TF1::GetFormula
virtual TFormula * GetFormula()
Definition: TF1.h:449
TF1::fParMax
std::vector< Double_t > fParMax
Definition: TF1.h:255
TF1::GetXmin
virtual Double_t GetXmin() const
Definition: TF1.h:548
TF1::ReleaseParameter
virtual void ReleaseParameter(Int_t ipar)
Release parameter number ipar If used in a fit, the parameter can vary freely.
Definition: TF1.cxx:3145
TFormula::GetParameter
Double_t GetParameter(const char *name) const