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