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