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