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 "TMethodCall.h"
31 #include "TAttLine.h"
32 #include "TAttFill.h"
33 #include "TAttMarker.h"
34 #include "TF1AbsComposition.h"
35 #include "TMath.h"
36 #include "Math/Types.h"
37 #include "Math/ParamFunctor.h"
38 
39 class TF1;
40 class TH1;
41 class TAxis;
42 class TRandom;
43 
44 namespace ROOT {
45  namespace Fit {
46  class FitResult;
47  }
48 }
49 
51 public:
52  TF1Parameters() {} // needed for the I/O
54  fParameters(std::vector<Double_t>(npar)),
55  fParNames(std::vector<std::string>(npar))
56  {
57  for (int i = 0; i < npar; ++i) {
58  fParNames[i] = std::string(TString::Format("p%d", i).Data());
59  }
60  }
61  // copy constructor
64  fParNames(rhs.fParNames)
65  {}
66  // assignment
68  {
69  if (&rhs == this) return *this;
71  fParNames = rhs.fParNames;
72  return *this;
73  }
74  virtual ~TF1Parameters() {}
75 
76  // getter methods
77  Double_t GetParameter(Int_t iparam) const
78  {
79  return (CheckIndex(iparam)) ? fParameters[iparam] : 0;
80  }
81  Double_t GetParameter(const char *name) const
82  {
84  }
85  const Double_t *GetParameters() const
86  {
87  return fParameters.data();
88  }
89  const std::vector<double> &ParamsVec() const
90  {
91  return fParameters;
92  }
93 
94  Int_t GetParNumber(const char *name) const;
95 
96  const char *GetParName(Int_t iparam) const
97  {
98  return (CheckIndex(iparam)) ? fParNames[iparam].c_str() : "";
99  }
100 
101 
102  // setter methods
103  void SetParameter(Int_t iparam, Double_t value)
104  {
105  if (!CheckIndex(iparam)) return;
106  fParameters[iparam] = value;
107  }
108  void SetParameters(const Double_t *params)
109  {
110  std::copy(params, params + fParameters.size(), fParameters.begin());
111  }
112  void SetParameters(Double_t p0, Double_t p1, Double_t p2 = 0, Double_t p3 = 0, Double_t p4 = 0,
113  Double_t p5 = 0, Double_t p6 = 0, Double_t p7 = 0, Double_t p8 = 0,
114  Double_t p9 = 0, Double_t p10 = 0);
115 
116  void SetParameter(const char *name, Double_t value)
117  {
118  SetParameter(GetParNumber(name), value);
119  }
120  void SetParName(Int_t iparam, const char *name)
121  {
122  if (!CheckIndex(iparam)) return;
123  fParNames[iparam] = std::string(name);
124  }
125  void SetParNames(const char *name0 = "p0", const char *name1 = "p1", const char *name2 = "p2",
126  const char *name3 = "p3", const char *name4 = "p4", const char *name5 = "p5",
127  const char *name6 = "p6", const char *name7 = "p7", const char *name8 = "p8",
128  const char *name9 = "p9", const char *name10 = "p10");
129 
130 
131 
132  ClassDef(TF1Parameters, 1) // The Parameters of a parameteric function
133 private:
134 
135  bool CheckIndex(Int_t i) const
136  {
137  return (i >= 0 && i < int(fParameters.size()));
138  }
139 
140  std::vector<Double_t> fParameters; // parameter values
141  std::vector<std::string> fParNames; // parameter names
142 };
143 
144 namespace ROOT {
145  namespace Internal {
146  /// Internal class used by TF1 for defining
147  /// template specialization for different TF1 constructors
148  template<class Func>
149  struct TF1Builder {
150  static void Build(TF1 *f, Func func);
151  };
152 
153  template<class Func>
154  struct TF1Builder<Func *> {
155  static void Build(TF1 *f, Func *func);
156  };
157 
158  // Internal class used by TF1 for obtaining the type from a functor
159  // out of the set of valid operator() signatures.
160  template<typename T>
161  struct GetFunctorType {
162  };
163 
164  template<typename F, typename T>
165  struct GetFunctorType<T(F::*)(const T *, const double *)> {
166  using type = T;
167  };
168 
169  template<typename F, typename T>
170  struct GetFunctorType<T(F::*)(const T *, const double *) const> {
171  using type = T;
172  };
173 
174  template<typename F, typename T>
175  struct GetFunctorType<T(F::*)(T *, double *)> {
176  using type = T;
177  };
178 
179  template<typename F, typename T>
180  struct GetFunctorType<T(F::*)(T *, double *) const> {
181  using type = T;
182  };
183 
184  // Internal class used by TF1 to get the right operator() signature
185  // from a Functor with several ones.
186  template<typename T, typename F>
187  auto GetTheRightOp(T(F::*opPtr)(const T *, const double *)) -> decltype(opPtr)
188  {
189  return opPtr;
190  }
191 
192  template<typename T, typename F>
193  auto GetTheRightOp(T(F::*opPtr)(const T *, const double *) const) -> decltype(opPtr)
194  {
195  return opPtr;
196  }
197 
198  template<typename T, typename F>
199  auto GetTheRightOp(T(F::*opPtr)(T *, double *)) -> decltype(opPtr)
200  {
201  return opPtr;
202  }
203 
204  template<typename T, typename F>
205  auto GetTheRightOp(T(F::*opPtr)(T *, double *) const) -> decltype(opPtr)
206  {
207  return opPtr;
208  }
209  }
210 }
211 
212 
213 class TF1 : public TNamed, public TAttLine, public TAttFill, public TAttMarker {
214 
215  template<class Func>
217 
218 public:
219  // Add to list behavior
220  enum class EAddToList {
221  kDefault,
222  kAdd,
223  kNo
224  };
225 
226 
228  virtual ~TF1FunctorPointer() {}
229  virtual TF1FunctorPointer * Clone() const = 0;
230  };
231 
232 protected:
233 
234  enum EFType {
235  kFormula = 0, // formula functions which can be stored,
236  kPtrScalarFreeFcn, // pointer to scalar free function,
237  kInterpreted, // interpreted functions constructed by name,
238  kTemplVec, // vectorized free functions or TemplScalar functors evaluating on vectorized parameters,
239  kTemplScalar, // TemplScalar functors evaluating on scalar parameters
241  }; // formula based on composition class (e.g. NSUM, CONV)
242 
243  Double_t fXmin{-1111}; //Lower bounds for the range
244  Double_t fXmax{-1111}; //Upper bounds for the range
245  Int_t fNpar{}; //Number of parameters
246  Int_t fNdim{}; //Function dimension
247  Int_t fNpx{100}; //Number of points used for the graphical representation
248  EFType fType{EFType::kTemplScalar};
249  Int_t fNpfits{}; //Number of points used in the fit
250  Int_t fNDF{}; //Number of degrees of freedom in the fit
251  Double_t fChisquare{}; //Function fit chisquare
252  Double_t fMinimum{-1111}; //Minimum value for plotting
253  Double_t fMaximum{-1111}; //Maximum value for plotting
254  std::vector<Double_t> fParErrors; //Array of errors of the fNpar parameters
255  std::vector<Double_t> fParMin; //Array of lower limits of the fNpar parameters
256  std::vector<Double_t> fParMax; //Array of upper limits of the fNpar parameters
257  std::vector<Double_t> fSave; //Array of fNsave function values
258  std::vector<Double_t> fIntegral; //!Integral of function binned on fNpx bins
259  std::vector<Double_t> fAlpha; //!Array alpha. for each bin in x the deconvolution r of fIntegral
260  std::vector<Double_t> fBeta; //!Array beta. is approximated by x = alpha +beta*r *gamma*r**2
261  std::vector<Double_t> fGamma; //!Array gamma.
262  TObject *fParent{nullptr}; //!Parent object hooking this function (if one)
263  TH1 *fHistogram{nullptr}; //!Pointer to histogram used for visualisation
264  std::unique_ptr<TMethodCall> fMethodCall; //!Pointer to MethodCall in case of interpreted function
265  Bool_t fNormalized{false}; //Normalization option (false by default)
266  Double_t fNormIntegral{}; //Integral of the function before being normalized
267  std::unique_ptr<TF1FunctorPointer> fFunctor; //! Functor object to wrap any C++ callable object
268  std::unique_ptr<TFormula> fFormula; //Pointer to TFormula in case when user define formula
269  std::unique_ptr<TF1Parameters> fParams; //Pointer to Function parameters object (exists only for not-formula functions)
270  std::unique_ptr<TF1AbsComposition> fComposition; //Pointer to composition (NSUM or CONV)
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)
276  {
277  fParams.reset(params);
278  fFunctor.reset(functor);
279  DoInitialize(addToGlobList);
280  }
281 
282 private:
283  // NSUM parsing helper functions
284  void DefineNSUMTerm(TObjArray *newFuncs, TObjArray *coeffNames,
285  TString &fullFormula,
286  TString &formula, int termStart, int termEnd,
288  int TermCoeffLength(TString &term);
289 
290 protected:
291 
292  template <class T>
295  TF1FunctorPointerImpl(const std::function<T(const T *f, const Double_t *param)> &func) : fImpl(func){};
297  virtual TF1FunctorPointer * Clone() const { return new TF1FunctorPointerImpl<T>(fImpl); }
299  };
300 
301 
302 
303 
304  static std::atomic<Bool_t> fgAbsValue; //use absolute value of function when computing integral
305  static Bool_t fgRejectPoint; //True if point must be rejected in a fit
306  static std::atomic<Bool_t> fgAddToGlobList; //True if we want to register the function in the global list
307  static TF1 *fgCurrent; //pointer to current function being processed
308 
309 
310  //void CreateFromFunctor(const char *name, Int_t npar, Int_t ndim = 1);
311  void DoInitialize(EAddToList addToGlobList);
312 
314  // tabulate the cumulative function integral at fNpx points. Used by GetRandom
316 
317  virtual Double_t GetMinMaxNDim(Double_t *x , Bool_t findmax, Double_t epsilon = 0, Int_t maxiter = 0) const;
318  virtual void GetRange(Double_t *xmin, Double_t *xmax) const;
319  virtual TH1 *DoCreateHistogram(Double_t xmin, Double_t xmax, Bool_t recreate = kFALSE);
320 
321 public:
322 
323  // TF1 status bits
324  enum EStatusBits {
325  kNotGlobal = BIT(10), // don't register in global list of functions
326  kNotDraw = BIT(9) // don't draw the function when in a TH1
327  };
328 
329  TF1();
330  TF1(const char *name, const char *formula, Double_t xmin = 0, Double_t xmax = 1, EAddToList addToGlobList = EAddToList::kDefault, bool vectorize = false);
331  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
332  TF1(const char *name, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault);
333  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);
334  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);
335 
336  template <class T>
337  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):
338  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
339  {
340  fType = std::is_same<T, double>::value ? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
341  }
342 
343  ////////////////////////////////////////////////////////////////////////////////
344  /// Constructor using a pointer to function.
345  ///
346  /// \param npar is the number of free parameters used by the function
347  ///
348  /// This constructor creates a function of type C when invoked
349  /// with the normal C++ compiler.
350  ///
351  ///
352  /// WARNING! A function created with this constructor cannot be Cloned
353 
354 
355  template <class T>
356  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):
357  TF1(EFType::kTemplVec, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
358  {}
359 
360  // Constructors using functors (compiled mode only)
361  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);
362 
363  // Template constructors from any C++ callable object, defining the operator() (double * , double *)
364  // and returning a double.
365  // The class name is not needed when using compile code, while it is required when using
366  // interpreted code via the specialized constructor with void *.
367  // An instance of the C++ function class or its pointer can both be used. The former is reccomended when using
368  // C++ compiled code, but if CINT compatibility is needed, then a pointer to the function class must be used.
369  // xmin and xmax specify the plotting range, npar is the number of parameters.
370  // See the tutorial math/exampleFunctor.C for an example of using this constructor
371  template <typename Func>
372  TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault) :
373  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList)
374  {
375  //actual fType set in TF1Builder
377  }
378 
379  // backward compatible interface
380  template <typename Func>
381  TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *, EAddToList addToGlobList = EAddToList::kDefault) :
382  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar))
383  {
385  }
386 
387 
388  // Template constructors from a pointer to any C++ class of type PtrObj with a specific member function of type
389  // MemFn.
390  // The member function must have the signature of (double * , double *) and returning a double.
391  // The class name and the method name are not needed when using compile code
392  // (the member function pointer is used in this case), while they are required when using interpreted
393  // code via the specialized constructor with void *.
394  // xmin and xmax specify the plotting range, npar is the number of parameters.
395  // See the tutorial math/exampleFunctor.C for an example of using this constructor
396  template <class PtrObj, typename MemFn>
397  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) :
398  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
399  {}
400 
401  // backward compatible interface
402  template <class PtrObj, typename MemFn>
403  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) :
404  TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
405  {}
406 
407  TF1(const TF1 &f1);
408  TF1 &operator=(const TF1 &rhs);
409  virtual ~TF1();
410  virtual void AddParameter(const TString &name, Double_t value)
411  {
412  if (fFormula) fFormula->AddParameter(name, value);
413  }
414  // virtual void AddParameters(const pair<TString,Double_t> *pairs, Int_t size) { fFormula->AddParameters(pairs,size); }
415  // virtual void AddVariable(const TString &name, Double_t value = 0) { if (fFormula) fFormula->AddVariable(name,value); }
416  // virtual void AddVariables(const TString *vars, Int_t size) { if (fFormula) fFormula->AddVariables(vars,size); }
417  virtual Bool_t AddToGlobalList(Bool_t on = kTRUE);
419  virtual void Browse(TBrowser *b);
420  virtual void Copy(TObject &f1) const;
421  TObject* Clone(const char* newname=0) const;
422  virtual Double_t Derivative(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
423  virtual Double_t Derivative2(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
424  virtual Double_t Derivative3(Double_t x, Double_t *params = 0, Double_t epsilon = 0.001) const;
425  static Double_t DerivativeError();
426  virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
427  virtual void Draw(Option_t *option = "");
428  virtual TF1 *DrawCopy(Option_t *option = "") const;
429  virtual TObject *DrawDerivative(Option_t *option = "al"); // *MENU*
430  virtual TObject *DrawIntegral(Option_t *option = "al"); // *MENU*
431  virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option = "");
432  virtual Double_t Eval(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
433  //template <class T> T Eval(T x, T y = 0, T z = 0, T t = 0) const;
434  virtual Double_t EvalPar(const Double_t *x, const Double_t *params = 0);
435  template <class T> T EvalPar(const T *x, const Double_t *params = 0);
436  virtual Double_t operator()(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
437  template <class T> T operator()(const T *x, const Double_t *params = nullptr);
438  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
439  virtual void FixParameter(Int_t ipar, Double_t value);
441  {
442  return (fType == EFType::kTemplVec) || (fType == EFType::kFormula && fFormula && fFormula->IsVectorized());
443  }
445  {
446  return fChisquare;
447  }
448  virtual TH1 *GetHistogram() const;
449  virtual TH1 *CreateHistogram()
450  {
451  return DoCreateHistogram(fXmin, fXmax);
452  }
453  virtual TFormula *GetFormula()
454  {
455  return fFormula.get();
456  }
457  virtual const TFormula *GetFormula() const
458  {
459  return fFormula.get();
460  }
461  virtual TString GetExpFormula(Option_t *option = "") const
462  {
463  return (fFormula) ? fFormula->GetExpFormula(option) : "";
464  }
465  virtual const TObject *GetLinearPart(Int_t i) const
466  {
467  return (fFormula) ? fFormula->GetLinearPart(i) : nullptr;
468  }
469  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;
470  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;
471  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;
472  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;
473  virtual Double_t GetMaximumStored() const
474  {
475  return fMaximum;
476  }
477  virtual Double_t GetMinimumStored() const
478  {
479  return fMinimum;
480  }
481  virtual Int_t GetNpar() const
482  {
483  return fNpar;
484  }
485  virtual Int_t GetNdim() const
486  {
487  return fNdim;
488  }
489  virtual Int_t GetNDF() const;
490  virtual Int_t GetNpx() const
491  {
492  return fNpx;
493  }
495  {
496  return fMethodCall.get();
497  }
498  virtual Int_t GetNumber() const
499  {
500  return (fFormula) ? fFormula->GetNumber() : 0;
501  }
502  virtual Int_t GetNumberFreeParameters() const;
503  virtual Int_t GetNumberFitPoints() const
504  {
505  return fNpfits;
506  }
507  virtual char *GetObjectInfo(Int_t px, Int_t py) const;
509  {
510  return fParent;
511  }
512  virtual Double_t GetParameter(Int_t ipar) const
513  {
514  return (fFormula) ? fFormula->GetParameter(ipar) : fParams->GetParameter(ipar);
515  }
516  virtual Double_t GetParameter(const TString &name) const
517  {
518  return (fFormula) ? fFormula->GetParameter(name) : fParams->GetParameter(name);
519  }
520  virtual Double_t *GetParameters() const
521  {
522  return (fFormula) ? fFormula->GetParameters() : const_cast<Double_t *>(fParams->GetParameters());
523  }
524  virtual void GetParameters(Double_t *params)
525  {
526  if (fFormula) fFormula->GetParameters(params);
527  else std::copy(fParams->ParamsVec().begin(), fParams->ParamsVec().end(), params);
528  }
529  virtual const char *GetParName(Int_t ipar) const
530  {
531  return (fFormula) ? fFormula->GetParName(ipar) : fParams->GetParName(ipar);
532  }
533  virtual Int_t GetParNumber(const char *name) const
534  {
535  return (fFormula) ? fFormula->GetParNumber(name) : fParams->GetParNumber(name);
536  }
537  virtual Double_t GetParError(Int_t ipar) const;
538  virtual const Double_t *GetParErrors() const
539  {
540  return fParErrors.data();
541  }
542  virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const;
543  virtual Double_t GetProb() const;
544  virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum);
545  virtual Double_t GetRandom(TRandom * rng = nullptr, Option_t * opt = nullptr);
546  virtual Double_t GetRandom(Double_t xmin, Double_t xmax, TRandom * rng = nullptr, Option_t * opt = nullptr);
547  virtual void GetRange(Double_t &xmin, Double_t &xmax) const;
548  virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const;
549  virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const;
550  virtual Double_t GetSave(const Double_t *x);
551  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;
552  virtual Double_t GetXmin() const
553  {
554  return fXmin;
555  }
556  virtual Double_t GetXmax() const
557  {
558  return fXmax;
559  }
560  TAxis *GetXaxis() const ;
561  TAxis *GetYaxis() const ;
562  TAxis *GetZaxis() const ;
564  {
565  return (fFormula) ? fFormula->GetVariable(name) : 0;
566  }
567  virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps = 0.01);
568  template <class T>
569  T GradientPar(Int_t ipar, const T *x, Double_t eps = 0.01);
570  template <class T>
571  T GradientParTempl(Int_t ipar, const T *x, Double_t eps = 0.01);
572 
573  virtual void GradientPar(const Double_t *x, Double_t *grad, Double_t eps = 0.01);
574  template <class T>
575  void GradientPar(const T *x, T *grad, Double_t eps = 0.01);
576  template <class T>
577  void GradientParTempl(const T *x, T *grad, Double_t eps = 0.01);
578 
579  virtual void InitArgs(const Double_t *x, const Double_t *params);
580  static void InitStandardFunctions();
581  virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel = 1.e-12);
582  virtual Double_t IntegralOneDim(Double_t a, Double_t b, Double_t epsrel, Double_t epsabs, Double_t &err);
583  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);
584  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);
585  // virtual Double_t IntegralFast(const TGraph *g, Double_t a, Double_t b, Double_t *params=0);
586  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);
587  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);
588  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)
589  {
590  return IntegralMultiple(n, a, b, maxpts, epsrel, epsrel, relerr, nfnevl, ifail);
591  }
592  virtual Double_t IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Double_t epsrel, Double_t &relerr);
593  virtual Bool_t IsEvalNormalized() const
594  {
595  return fNormalized;
596  }
597  /// return kTRUE if the point is inside the function range
598  virtual Bool_t IsInside(const Double_t *x) const
599  {
600  return !((x[0] < fXmin) || (x[0] > fXmax));
601  }
602  virtual Bool_t IsLinear() const
603  {
604  return (fFormula) ? fFormula->IsLinear() : false;
605  }
606  virtual Bool_t IsValid() const;
607  virtual void Print(Option_t *option = "") const;
608  virtual void Paint(Option_t *option = "");
609  virtual void ReleaseParameter(Int_t ipar);
610  virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax);
611  virtual void SavePrimitive(std::ostream &out, Option_t *option = "");
612  virtual void SetChisquare(Double_t chi2)
613  {
614  fChisquare = chi2;
615  }
616  virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar = 0);
617  template <class PtrObj, typename MemFn>
618  void SetFunction(PtrObj &p, MemFn memFn);
619  template <typename Func>
620  void SetFunction(Func f);
621  virtual void SetMaximum(Double_t maximum = -1111); // *MENU*
622  virtual void SetMinimum(Double_t minimum = -1111); // *MENU*
623  virtual void SetNDF(Int_t ndf);
624  virtual void SetNumberFitPoints(Int_t npfits)
625  {
626  fNpfits = npfits;
627  }
628  virtual void SetNormalized(Bool_t flag)
629  {
630  fNormalized = flag;
631  Update();
632  }
633  virtual void SetNpx(Int_t npx = 100); // *MENU*
634  virtual void SetParameter(Int_t param, Double_t value)
635  {
636  (fFormula) ? fFormula->SetParameter(param, value) : fParams->SetParameter(param, value);
637  Update();
638  }
639  virtual void SetParameter(const TString &name, Double_t value)
640  {
641  (fFormula) ? fFormula->SetParameter(name, value) : fParams->SetParameter(name, value);
642  Update();
643  }
644  virtual void SetParameters(const Double_t *params)
645  {
646  (fFormula) ? fFormula->SetParameters(params) : fParams->SetParameters(params);
647  Update();
648  }
649  virtual void SetParameters(Double_t p0, Double_t p1, Double_t p2 = 0, Double_t p3 = 0, Double_t p4 = 0,
650  Double_t p5 = 0, Double_t p6 = 0, Double_t p7 = 0, Double_t p8 = 0,
651  Double_t p9 = 0, Double_t p10 = 0)
652  {
653  if (fFormula) fFormula->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
654  else fParams->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
655  Update();
656  } // *MENU*
657  virtual void SetParName(Int_t ipar, const char *name);
658  virtual void SetParNames(const char *name0 = "p0", const char *name1 = "p1", const char *name2 = "p2",
659  const char *name3 = "p3", const char *name4 = "p4", const char *name5 = "p5",
660  const char *name6 = "p6", const char *name7 = "p7", const char *name8 = "p8",
661  const char *name9 = "p9", const char *name10 = "p10"); // *MENU*
662  virtual void SetParError(Int_t ipar, Double_t error);
663  virtual void SetParErrors(const Double_t *errors);
664  virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax);
665  virtual void SetParent(TObject *p = 0)
666  {
667  fParent = p;
668  }
669  virtual void SetRange(Double_t xmin, Double_t xmax); // *MENU*
672  virtual void SetSavedPoint(Int_t point, Double_t value);
673  virtual void SetTitle(const char *title = ""); // *MENU*
674  virtual void SetVectorized(Bool_t vectorized)
675  {
676  if (fType == EFType::kFormula && fFormula)
677  fFormula->SetVectorized(vectorized);
678  else
679  Warning("SetVectorized", "Can only set vectorized flag on formula-based TF1");
680  }
681  virtual void Update();
682 
683  static TF1 *GetCurrent();
684  static void AbsValue(Bool_t reject = kTRUE);
685  static void RejectPoint(Bool_t reject = kTRUE);
686  static Bool_t RejectedPoint();
687  static void SetCurrent(TF1 *f1);
688 
689  //Moments
690  virtual Double_t Moment(Double_t n, Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001);
691  virtual Double_t CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001);
692  virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001)
693  {
694  return Moment(1, a, b, params, epsilon);
695  }
696  virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params = 0, Double_t epsilon = 0.000001)
697  {
698  return CentralMoment(2, a, b, params, epsilon);
699  }
700 
701  //some useful static utility functions to compute sampling points for Integral
702  //static void CalcGaussLegendreSamplingPoints(TGraph *g, Double_t eps=3.0e-11);
703  //static TGraph *CalcGaussLegendreSamplingPoints(Int_t num=21, Double_t eps=3.0e-11);
704  static void CalcGaussLegendreSamplingPoints(Int_t num, Double_t *x, Double_t *w, Double_t eps = 3.0e-11);
705 
706 private:
707  template <class T>
708  T EvalParTempl(const T *data, const Double_t *params = 0);
709 
710 #ifdef R__HAS_VECCORE
711  inline double EvalParVec(const Double_t *data, const Double_t *params);
712 #endif
713 
714  ClassDef(TF1, 12) // The Parametric 1-D function
715 };
716 
717 namespace ROOT {
718  namespace Internal {
719 
720  template<class Func>
721  void TF1Builder<Func>::Build(TF1 *f, Func func)
722  {
723  using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
724  f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
726  f->fParams.reset(new TF1Parameters(f->fNpar));
727  }
728 
729  template<class Func>
730  void TF1Builder<Func *>::Build(TF1 *f, Func *func)
731  {
732  using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
733  f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
735  f->fParams.reset(new TF1Parameters(f->fNpar));
736  }
737 
738  /// TF1 building from a string
739  /// used to build a TFormula based on a lambda function
740  template<>
741  struct TF1Builder<const char *> {
742  static void Build(TF1 *f, const char *formula)
743  {
744  f->fType = TF1::EFType::kFormula;
745  f->fFormula.reset(new TFormula("tf1lambda", formula, f->fNdim, f->fNpar, false));
746  TString formulaExpression(formula);
747  Ssiz_t first = formulaExpression.Index("return") + 7;
748  Ssiz_t last = formulaExpression.Last(';');
749  TString title = formulaExpression(first, last - first);
750  f->SetTitle(title);
751  }
752  };
753  }
754 }
755 
757 {
758  return Eval(x, y, z, t);
759 }
760 
761 template <class T>
762 inline T TF1::operator()(const T *x, const Double_t *params)
763 {
764  return EvalPar(x, params);
765 }
766 
767 ////////////////////////////////////////////////////////////////////////////////
768 /// EvalPar for vectorized
769 template <class T>
770 T TF1::EvalPar(const T *x, const Double_t *params)
771 {
772  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
773  return EvalParTempl(x, params);
774  } else if (fType == EFType::kFormula) {
775  return fFormula->EvalPar(x, params);
776  } else
777  return TF1::EvalPar((double *)x, params);
778 }
779 
780 ////////////////////////////////////////////////////////////////////////////////
781 /// Eval for vectorized functions
782 // template <class T>
783 // T TF1::Eval(T x, T y, T z, T t) const
784 // {
785 // if (fType == EFType::kFormula)
786 // return fFormula->Eval(x, y, z, t);
787 
788 // T xx[] = {x, y, z, t};
789 // Double_t *pp = (Double_t *)fParams->GetParameters();
790 // return ((TF1 *)this)->EvalPar(xx, pp);
791 // }
792 
793 // Internal to TF1. Evaluates Templated interfaces
794 template <class T>
795 inline T TF1::EvalParTempl(const T *data, const Double_t *params)
796 {
797  assert(fType == EFType::kTemplScalar || fType == EFType::kTemplVec);
798  if (!params) params = (Double_t *)fParams->GetParameters();
799  if (fFunctor)
800  return ((TF1FunctorPointerImpl<T> *)fFunctor.get())->fImpl(data, params);
801 
802  // this should throw an error
803  // we nned to implement a vectorized GetSave(x)
804  return TMath::SignalingNaN();
805 }
806 
807 #ifdef R__HAS_VECCORE
808 // Internal to TF1. Evaluates Vectorized TF1 on data of type Double_v
809 inline double TF1::EvalParVec(const Double_t *data, const Double_t *params)
810 {
811  assert(fType == EFType::kTemplVec);
812  std::vector<ROOT::Double_v> d(fNdim);
813  ROOT::Double_v res;
814 
815  for(auto i=0; i<fNdim; i++) {
816  d[i] = ROOT::Double_v(data[i]);
817  }
818 
819  if (fFunctor) {
820  res = ((TF1FunctorPointerImpl<ROOT::Double_v> *) fFunctor.get())->fImpl(d.data(), params);
821  } else {
822  // res = GetSave(x);
823  return TMath::SignalingNaN();
824  }
825  return vecCore::Get<ROOT::Double_v>(res, 0);
826 }
827 #endif
828 
830 {
832 }
834 {
836 }
837 
838 template <typename Func>
839 void TF1::SetFunction(Func f)
840 {
841  // set function from a generic C++ callable object
842  fType = EFType::kPtrScalarFreeFcn;
844 }
845 template <class PtrObj, typename MemFn>
846 void TF1::SetFunction(PtrObj &p, MemFn memFn)
847 {
848  // set from a pointer to a member function
849  fType = EFType::kPtrScalarFreeFcn;
851 }
852 
853 template <class T>
854 inline T TF1::GradientPar(Int_t ipar, const T *x, Double_t eps)
855 {
856  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
857  return GradientParTempl<T>(ipar, x, eps);
858  } else
859  return GradientParTempl<Double_t>(ipar, (const Double_t *)x, eps);
860 }
861 
862 template <class T>
863 inline T TF1::GradientParTempl(Int_t ipar, const T *x, Double_t eps)
864 {
865  if (GetNpar() == 0)
866  return 0;
867 
868  if (eps < 1e-10 || eps > 1) {
869  Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
870  eps = 0.01;
871  }
872  Double_t h;
873  TF1 *func = (TF1 *)this;
874  Double_t *parameters = GetParameters();
875 
876  // Copy parameters for thread safety
877  std::vector<Double_t> parametersCopy(parameters, parameters + GetNpar());
878  parameters = parametersCopy.data();
879 
880  Double_t al, bl, h2;
881  T f1, f2, g1, g2, d0, d2;
882 
883  ((TF1 *)this)->GetParLimits(ipar, al, bl);
884  if (al * bl != 0 && al >= bl) {
885  // this parameter is fixed
886  return 0;
887  }
888 
889  // check if error has been computer (is not zero)
890  if (func->GetParError(ipar) != 0)
891  h = eps * func->GetParError(ipar);
892  else
893  h = eps;
894 
895  // save original parameters
896  Double_t par0 = parameters[ipar];
897 
898  parameters[ipar] = par0 + h;
899  f1 = func->EvalPar(x, parameters);
900  parameters[ipar] = par0 - h;
901  f2 = func->EvalPar(x, parameters);
902  parameters[ipar] = par0 + h / 2;
903  g1 = func->EvalPar(x, parameters);
904  parameters[ipar] = par0 - h / 2;
905  g2 = func->EvalPar(x, parameters);
906 
907  // compute the central differences
908  h2 = 1 / (2. * h);
909  d0 = f1 - f2;
910  d2 = 2 * (g1 - g2);
911 
912  T grad = h2 * (4 * d2 - d0) / 3.;
913 
914  // restore original value
915  parameters[ipar] = par0;
916 
917  return grad;
918 }
919 
920 template <class T>
921 inline void TF1::GradientPar(const T *x, T *grad, Double_t eps)
922 {
923  if (fType == EFType::kTemplVec || fType == EFType::kTemplScalar) {
924  GradientParTempl<T>(x, grad, eps);
925  } else
926  GradientParTempl<Double_t>((const Double_t *)x, (Double_t *)grad, eps);
927 }
928 
929 template <class T>
930 inline void TF1::GradientParTempl(const T *x, T *grad, Double_t eps)
931 {
932  if (eps < 1e-10 || eps > 1) {
933  Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
934  eps = 0.01;
935  }
936 
937  for (Int_t ipar = 0; ipar < GetNpar(); ipar++) {
938  grad[ipar] = GradientParTempl<T>(ipar, x, eps);
939  }
940 }
941 
942 #endif
TF1::TF1FunctorPointerImpl::fImpl
ROOT::Math::ParamFunctorTempl< T > fImpl
Definition: TF1.h:298
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:1813
TF1::fNDF
Int_t fNDF
Definition: TF1.h:250
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:372
TF1::GetHistogram
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition: TF1.cxx:1574
TF1::DrawF1
virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option="")
Draw function between xmin and xmax.
Definition: TF1.cxx:1415
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:1559
TF1::fgAbsValue
static std::atomic< Bool_t > fgAbsValue
Definition: TF1.h:304
first
Definition: first.py:1
TF1Parameters::SetParName
void SetParName(Int_t iparam, const char *name)
Definition: TF1.h:120
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:403
TBrowser
Using a TBrowser one can browse all ROOT objects.
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:588
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:3658
TF1::ExecuteEvent
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TF1.cxx:1531
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:91
e
#define e(i)
Definition: RSha256.hxx:103
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:1908
TF1::kFormula
@ kFormula
Definition: TF1.h:235
TF1::fSave
std::vector< Double_t > fSave
Definition: TF1.h:257
TObjArray
An array of TObjects.
Definition: TObjArray.h:37
TF1::fHistogram
TH1 * fHistogram
Parent object hooking this function (if one)
Definition: TF1.h:263
TF1::GetParErrors
virtual const Double_t * GetParErrors() const
Definition: TF1.h:538
f
#define f(i)
Definition: RSha256.hxx:104
TF1Parameters::SetParameters
void SetParameters(const Double_t *params)
Definition: TF1.h:108
TF1::fNpar
Int_t fNpar
Definition: TF1.h:245
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:1604
Option_t
const char Option_t
Definition: RtypesCore.h:66
TF1Parameters::GetParameter
Double_t GetParameter(const char *name) const
Definition: TF1.h:81
TF1::TF1FunctorPointerImpl::TF1FunctorPointerImpl
TF1FunctorPointerImpl(const std::function< T(const T *f, const Double_t *param)> &func)
Definition: TF1.h:295
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:696
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:1853
TFormula
The Formula class.
Definition: TFormula.h:87
TF1::kPtrScalarFreeFcn
@ kPtrScalarFreeFcn
Definition: TF1.h:236
ROOT::Internal::GetFunctorType< T(F::*)(const T *, const double *) const >::type
T type
Definition: TF1.h:171
TF1::fNormalized
Bool_t fNormalized
Pointer to MethodCall in case of interpreted function.
Definition: TF1.h:265
TF1Parameters::~TF1Parameters
virtual ~TF1Parameters()
Definition: TF1.h:74
TF1::fgCurrent
static TF1 * fgCurrent
Definition: TF1.h:307
TF1::SetNpx
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3446
TF1Parameters::ParamsVec
const std::vector< double > & ParamsVec() const
Definition: TF1.h:89
TF1::GetMaximumStored
virtual Double_t GetMaximumStored() const
Definition: TF1.h:473
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:2515
TF1::GetParameters
virtual Double_t * GetParameters() const
Definition: TF1.h:520
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:649
ROOT::Internal::GetFunctorType< T(F::*)(T *, double *)>::type
T type
Definition: TF1.h:176
TF1::EvalParTempl
T EvalParTempl(const T *data, const Double_t *params=0)
Eval for vectorized functions.
Definition: TF1.h:795
TF1::SetParameter
virtual void SetParameter(const TString &name, Double_t value)
Definition: TF1.h:639
F
#define F(x, y, z)
TF1::GetParent
TObject * GetParent() const
Definition: TF1.h:508
TF1::fNormIntegral
Double_t fNormIntegral
Definition: TF1.h:266
TF1::GetNdim
virtual Int_t GetNdim() const
Definition: TF1.h:485
TF1::CreateHistogram
virtual TH1 * CreateHistogram()
Definition: TF1.h:449
TF1::GetNpar
virtual Int_t GetNpar() const
Definition: TF1.h:481
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:3475
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:3407
TF1Parameters::GetParameters
const Double_t * GetParameters() const
Definition: TF1.h:85
xmax
float xmax
Definition: THbookFile.cxx:95
TF1::~TF1
virtual ~TF1()
TF1 default destructor.
Definition: TF1.cxx:942
TF1::GetProb
virtual Double_t GetProb() const
Return the fit probability.
Definition: TF1.cxx:1945
TF1::kTemplVec
@ kTemplVec
Definition: TF1.h:238
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:1232
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:3677
TF1::fType
EFType fType
Definition: TF1.h:248
TF1::fIntegral
std::vector< Double_t > fIntegral
Definition: TF1.h:258
TAttMarker.h
TF1::SetCurrent
static void SetCurrent(TF1 *f1)
Static function setting the current function.
Definition: TF1.cxx:3356
TF1::kInterpreted
@ kInterpreted
Definition: TF1.h:237
TF1::GetNumberFitPoints
virtual Int_t GetNumberFitPoints() const
Definition: TF1.h:503
Int_t
int Int_t
Definition: RtypesCore.h:45
TF1::GetXmax
virtual Double_t GetXmax() const
Definition: TF1.h:556
TF1::GetFormula
virtual const TFormula * GetFormula() const
Definition: TF1.h:457
TF1::DistancetoPrimitive
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a function.
Definition: TF1.cxx:1282
TF1Parameters::GetParName
const char * GetParName(Int_t iparam) const
Definition: TF1.h:96
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:3200
x
Double_t x[n]
Definition: legend1.C:17
TF1::fXmin
Double_t fXmin
Definition: TF1.h:243
TF1::AddParameter
virtual void AddParameter(const TString &name, Double_t value)
Definition: TF1.h:410
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:1713
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:3497
TF1Parameters::TF1Parameters
TF1Parameters(const TF1Parameters &rhs)
Definition: TF1.h:62
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:598
TF1::kNotDraw
@ kNotDraw
Definition: TF1.h:326
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::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:3801
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:2331
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:3420
TF1::kNotGlobal
@ kNotGlobal
Definition: TF1.h:325
TString
Basic string class.
Definition: TString.h:136
TF1::IsVectorized
bool IsVectorized()
Definition: TF1.h:440
TF1::DerivativeError
static Double_t DerivativeError()
Static function returning the error of the last call to the of Derivative's functions.
Definition: TF1.cxx:1266
TF1Parameters::TF1Parameters
TF1Parameters(Int_t npar)
Definition: TF1.h:53
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:2937
b
#define b(i)
Definition: RSha256.hxx:100
TF1::ComputeCdfTable
Bool_t ComputeCdfTable(Option_t *opt)
Compute the cumulative function at fNpx points between fXmin and fXmax.
Definition: TF1.cxx:2069
TF1Parameters::operator=
TF1Parameters & operator=(const TF1Parameters &rhs)
Definition: TF1.h:67
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:756
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:2605
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:356
TF1::GetParLimits
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1930
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:3432
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:337
q
float * q
Definition: THbookFile.cxx:89
TF1::fMethodCall
std::unique_ptr< TMethodCall > fMethodCall
Pointer to histogram used for visualisation.
Definition: TF1.h:264
TString::Last
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:912
TF1::GetRange
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition: TF1.cxx:2269
TF1::GetNumber
virtual Int_t GetNumber() const
Definition: TF1.h:498
TF1::Clone
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TF1.cxx:1053
TF1::GetRandom
virtual Double_t GetRandom(TRandom *rng=nullptr, Option_t *opt=nullptr)
Return a random number following this function shape.
Definition: TF1.cxx:2180
TF1::SetParameter
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:634
TF1::SetParameters
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:644
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:2433
TF1::kTemplScalar
@ kTemplScalar
Definition: TF1.h:239
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:2762
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:3368
TAttLine.h
TF1::IsValid
virtual Bool_t IsValid() const
Return kTRUE if the function is valid.
Definition: TF1.cxx:2866
TF1Parameters::GetParameter
Double_t GetParameter(Int_t iparam) const
Definition: TF1.h:77
TF1::GetParNumber
virtual Int_t GetParNumber(const char *name) const
Definition: TF1.h:533
TF1::fNpfits
Int_t fNpfits
Definition: TF1.h:249
TF1::fNdim
Int_t fNdim
Definition: TF1.h:246
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:3821
TF1::SetParError
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition: TF1.cxx:3486
TF1::GetYaxis
TAxis * GetYaxis() const
Get y axis of the function.
Definition: TF1.cxx:2399
TF1::fXmax
Double_t fXmax
Definition: TF1.h:244
TF1::fAlpha
std::vector< Double_t > fAlpha
Integral of function binned on fNpx bins.
Definition: TF1.h:259
TF1::GetExpFormula
virtual TString GetExpFormula(Option_t *option="") const
Definition: TF1.h:461
TF1::TF1FunctorPointer::~TF1FunctorPointer
virtual ~TF1FunctorPointer()
Definition: TF1.h:228
TAttFill.h
TRandom
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
TAttLine
Line Attributes class.
Definition: TAttLine.h:18
ROOT::Internal::TF1Builder::Build
static void Build(TF1 *f, Func func)
Definition: TF1.h:721
TF1::InitStandardFunctions
static void InitStandardFunctions()
Create the basic function objects.
Definition: TF1.cxx:2482
TF1::SetFunction
void SetFunction(PtrObj &p, MemFn memFn)
Definition: TF1.h:846
TF1::GetMethodCall
TMethodCall * GetMethodCall() const
Definition: TF1.h:494
TF1::DrawDerivative
virtual TObject * DrawDerivative(Option_t *option="al")
Draw derivative of this function.
Definition: TF1.cxx:1374
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:3623
h
#define h(i)
Definition: RSha256.hxx:106
TF1::SetNormalized
virtual void SetNormalized(Bool_t flag)
Definition: TF1.h:628
TF1::GetXaxis
TAxis * GetXaxis() const
Get x axis of the function.
Definition: TF1.cxx:2388
TF1::InitArgs
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition: TF1.cxx:2467
epsilon
REAL epsilon
Definition: triangle.c:617
ROOT::Internal::GetFunctorType
Definition: TF1.h:161
a
auto * a
Definition: textangle.C:12
TF1::fMinimum
Double_t fMinimum
Definition: TF1.h:252
TNamed
The TNamed class is the base class for all named ROOT classes.
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:1167
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:92
TF1::GetParName
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:529
TF1::GetNpx
virtual Int_t GetNpx() const
Definition: TF1.h:490
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:397
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:3714
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:3036
TF1::Browse
virtual void Browse(TBrowser *b)
Browse.
Definition: TF1.cxx:982
ROOT::Fit::FitResult
class containg the result of the fit and all the related information (fitted parameter values,...
Definition: FitResult.h:47
TF1::SetChisquare
virtual void SetChisquare(Double_t chi2)
Definition: TF1.h:612
TF1::Copy
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition: TF1.cxx:994
TF1::EFType
EFType
Definition: TF1.h:234
TAttMarker
Marker Attributes class.
Definition: TAttMarker.h:19
TF1::fgRejectPoint
static Bool_t fgRejectPoint
Definition: TF1.h:305
ROOT::Internal::TF1Builder< const char * >::Build
static void Build(TF1 *f, const char *formula)
Definition: TF1.h:742
ROOT::Internal::GetFunctorType< T(F::*)(const T *, const double *)>::type
T type
Definition: TF1.h:166
BIT
#define BIT(n)
Definition: Rtypes.h:85
TF1::IsEvalNormalized
virtual Bool_t IsEvalNormalized() const
Definition: TF1.h:593
double
double
Definition: Converters.cxx:921
y
Double_t y[n]
Definition: legend1.C:17
ROOT::Double_v
Double_t Double_v
Definition: Types.h:51
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:381
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:2410
TF1::fFunctor
std::unique_ptr< TF1FunctorPointer > fFunctor
Definition: TF1.h:267
TF1::GetLinearPart
virtual const TObject * GetLinearPart(Int_t i) const
Definition: TF1.h:465
TF1::kCompositionFcn
@ kCompositionFcn
Definition: TF1.h:240
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:908
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:254
TF1::fMaximum
Double_t fMaximum
Definition: TF1.h:253
TF1::SetNumberFitPoints
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:624
TF1::TF1FunctorPointerImpl::Clone
virtual TF1FunctorPointer * Clone() const
Definition: TF1.h:297
TF1::fBeta
std::vector< Double_t > fBeta
Array alpha. for each bin in x the deconvolution r of fIntegral.
Definition: TF1.h:260
TF1::SetRange
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF1.cxx:3532
TF1::GetSave
virtual Double_t GetSave(const Double_t *x)
Get value corresponding to X in array of fSave values.
Definition: TF1.cxx:2332
ROOT::Internal::GetFunctorType< T(F::*)(T *, double *) const >::type
T type
Definition: TF1.h:181
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:1463
TF1Parameters::CheckIndex
bool CheckIndex(Int_t i) const
Definition: TF1.h:135
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:3853
TF1::fChisquare
Double_t fChisquare
Definition: TF1.h:251
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:255
TF1::EAddToList
EAddToList
Definition: TF1.h:220
TF1Parameters::fParNames
std::vector< std::string > fParNames
Definition: TF1.h:141
TF1::fParams
std::unique_ptr< TF1Parameters > fParams
Definition: TF1.h:269
TF1::GetMinimumStored
virtual Double_t GetMinimumStored() const
Definition: TF1.h:477
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:692
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:3546
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:1920
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:3562
TF1::TF1FunctorPointerImpl
Definition: TF1.h:293
TF1::DrawCopy
virtual TF1 * DrawCopy(Option_t *option="") const
Draw a copy of this function with its current attributes.
Definition: TF1.cxx:1352
TF1::TF1
TF1()
TF1 default constructor.
Definition: TF1.cxx:494
TF1::RejectedPoint
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3667
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:1434
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:1645
TF1::Draw
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1322
TF1::TF1FunctorPointerImpl::TF1FunctorPointerImpl
TF1FunctorPointerImpl(const ROOT::Math::ParamFunctorTempl< T > &func)
Definition: TF1.h:294
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:1686
TObject
Mother of all ROOT objects.
Definition: TObject.h:37
TF1::fFormula
std::unique_ptr< TFormula > fFormula
Functor object to wrap any C++ callable object.
Definition: TF1.h:268
TF1::fgAddToGlobList
static std::atomic< Bool_t > fgAddToGlobList
Definition: TF1.h:306
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:973
TH1
TH1 is the base class of all histogram classes in ROOT.
Definition: TH1.h:58
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:1102
TF1Parameters
TF1 Parameters class.
Definition: TF1.h:50
TF1::fParent
TObject * fParent
Array gamma.
Definition: TF1.h:262
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:3463
ROOT::Internal::GetTheRightOp
auto GetTheRightOp(T(F::*opPtr)(const T *, const double *)) -> decltype(opPtr)
Definition: TF1.h:187
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:34
TF1::operator=
TF1 & operator=(const TF1 &rhs)
Operator =.
Definition: TF1.cxx:930
d
#define d(i)
Definition: RSha256.hxx:102
TF1Parameters::TF1Parameters
TF1Parameters()
Definition: TF1.h:52
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:1879
ROOT::Math::ParamFunctor
ParamFunctorTempl< double > ParamFunctor
Definition: ParamFunctor.h:387
TF1Parameters::fParameters
std::vector< Double_t > fParameters
Definition: TF1.h:140
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::fGamma
std::vector< Double_t > fGamma
Array beta. is approximated by x = alpha +beta*r *gamma*r**2.
Definition: TF1.h:261
TF1::GetNumberFreeParameters
virtual Int_t GetNumberFreeParameters() const
Return the number of free parameters.
Definition: TF1.cxx:1890
TF1::GetVariable
virtual Double_t GetVariable(const TString &name)
Definition: TF1.h:563
TMethodCall
Method or function calling interface.
Definition: TMethodCall.h:37
TAttFill
Fill Area Attributes class.
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:3511
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:2835
TF1::GetParameter
virtual Double_t GetParameter(const TString &name) const
Definition: TF1.h:516
TF1::SetParent
virtual void SetParent(TObject *p=0)
Definition: TF1.h:665
ROOT::Internal::TF1Builder
Internal class used by TF1 for defining template specialization for different TF1 constructors.
Definition: TF1.h:149
TF1::IntegrateForNormalization
void IntegrateForNormalization()
type
int type
Definition: TGX11.cxx:121
TF1
1-Dim function class
Definition: TF1.h:213
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:1982
TF1::DrawIntegral
virtual TObject * DrawIntegral(Option_t *option="al")
Draw integral of this function.
Definition: TF1.cxx:1399
TF1Parameters::SetParameter
void SetParameter(const char *name, Double_t value)
Definition: TF1.h:116
ROOT::Math::ParamFunctorTempl
Param Functor class for Multidimensional functions.
Definition: ParamFunctor.h:273
TF1::TF1FunctorPointerImpl::~TF1FunctorPointerImpl
virtual ~TF1FunctorPointerImpl()
Definition: TF1.h:296
TF1::GetParameter
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:512
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:3147
TF1::GetParameters
virtual void GetParameters(Double_t *params)
Definition: TF1.h:524
TF1::fNpx
Int_t fNpx
Definition: TF1.h:247
TF1::fComposition
std::unique_ptr< TF1AbsComposition > fComposition
Definition: TF1.h:270
TF1::SetVectorized
virtual void SetVectorized(Bool_t vectorized)
Definition: TF1.h:674
TF1::GetChisquare
Double_t GetChisquare() const
Definition: TF1.h:444
ROOT
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Definition: EExecutionPolicy.hxx:4
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:1547
TF1Parameters::SetParameter
void SetParameter(Int_t iparam, Double_t value)
Definition: TF1.h:103
TF1::IsLinear
virtual Bool_t IsLinear() const
Definition: TF1.h:602
Math
Namespace for new Math classes and functions.
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:2692
TMath.h
int
TF1::Print
virtual void Print(Option_t *option="") const
Print TNamed name and title.
Definition: TF1.cxx:2881
TF1::GradientParTempl
T GradientParTempl(Int_t ipar, const T *x, Double_t eps=0.01)
Definition: TF1.h:863
TF1::GetFormula
virtual TFormula * GetFormula()
Definition: TF1.h:453
TF1::fParMax
std::vector< Double_t > fParMax
Definition: TF1.h:256
TF1::GetXmin
virtual Double_t GetXmin() const
Definition: TF1.h:552
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:3137
TMethodCall.h