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
40#include <string>
41
42class TF1;
43class TH1;
44class TAxis;
45class TRandom;
46
47namespace ROOT {
48 namespace Fit {
49 class FitResult;
50 }
51}
52
54public:
55 TF1Parameters() {} // needed for the I/O
57 fParameters(std::vector<Double_t>(npar)),
58 fParNames(std::vector<std::string>(npar))
59 {
60 for (int i = 0; i < npar; ++i) {
61 fParNames[i] = std::string(TString::Format("p%d", i).Data());
62 }
63 }
64 // copy constructor
68 {}
69 // assignment
71 {
72 if (&rhs == this) return *this;
74 fParNames = rhs.fParNames;
75 return *this;
76 }
77 virtual ~TF1Parameters() {}
78
79 // getter methods
81 {
82 return (CheckIndex(iparam)) ? fParameters[iparam] : 0;
83 }
84 Double_t GetParameter(const char *name) const
85 {
87 }
88 const Double_t *GetParameters() const
89 {
90 return fParameters.data();
91 }
92 const std::vector<double> &ParamsVec() const
93 {
94 return fParameters;
95 }
96
97 Int_t GetParNumber(const char *name) const;
98
99 const char *GetParName(Int_t iparam) const
100 {
101 return (CheckIndex(iparam)) ? fParNames[iparam].c_str() : "";
102 }
103
104
105 // setter methods
107 {
108 if (!CheckIndex(iparam)) return;
109 fParameters[iparam] = value;
110 }
111 void SetParameters(const Double_t *params)
112 {
113 std::copy(params, params + fParameters.size(), fParameters.begin());
114 }
115 template <typename... Args>
116 void SetParameters(Double_t arg1, Args &&... args);
117
118 void SetParameter(const char *name, Double_t value)
119 {
121 }
122 void SetParName(Int_t iparam, const char *name)
123 {
124 if (!CheckIndex(iparam)) return;
125 fParNames[iparam] = std::string(name);
126 }
127 template <typename... Args>
128 void SetParNames(Args &&... args);
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
142/// Set parameter values.
143/// NaN values will be skipped, meaning that the corresponding parameters will not be changed.
144template <typename... Args>
145void TF1Parameters::SetParameters(Double_t arg1, Args &&...args)
146{
147 int i = 0;
148 for (double val : {arg1, static_cast<Double_t>(args)...}) {
149 if(!TMath::IsNaN(val)) SetParameter(i++, val);
150 }
151}
152
153/// Set parameter names.
154/// Empty strings will be skipped, meaning that the corresponding name will not be changed.
155template <typename... Args>
156void TF1Parameters::SetParNames(Args &&...args)
157{
158 int i = 0;
159 for (auto name : {static_cast<std::string const&>(args)...}) {
160 if(!name.empty()) SetParName(i++, name.c_str());
161 }
162}
163
164namespace ROOT {
165 namespace Internal {
166 /// %Internal class used by TF1 for defining
167 /// template specialization for different TF1 constructors
168 template<class Func>
169 struct TF1Builder {
170 static void Build(TF1 *f, Func func);
171 };
172
173 template<class Func>
174 struct TF1Builder<Func *> {
175 static void Build(TF1 *f, Func *func);
176 };
177
178 /// %Internal class used by TF1 for obtaining the type from a functor
179 /// out of the set of valid operator() signatures.
180 template<typename T>
182 };
183
184 template<typename F, typename T>
185 struct GetFunctorType<T(F::*)(const T *, const double *)> {
186 using type = T;
187 };
188
189 template<typename F, typename T>
190 struct GetFunctorType<T(F::*)(const T *, const double *) const> {
191 using type = T;
192 };
193
194 template<typename F, typename T>
195 struct GetFunctorType<T(F::*)(T *, double *)> {
196 using type = T;
197 };
198
199 template<typename F, typename T>
200 struct GetFunctorType<T(F::*)(T *, double *) const> {
201 using type = T;
202 };
203
204 /// %Internal class used by TF1 to get the right operator() signature
205 /// from a Functor with several ones.
206 template<typename T, typename F>
207 auto GetTheRightOp(T(F::*opPtr)(const T *, const double *)) -> decltype(opPtr)
208 {
209 return opPtr;
210 }
211
212 template<typename T, typename F>
213 auto GetTheRightOp(T(F::*opPtr)(const T *, const double *) const) -> decltype(opPtr)
214 {
215 return opPtr;
216 }
217
218 template<typename T, typename F>
219 auto GetTheRightOp(T(F::*opPtr)(T *, double *)) -> decltype(opPtr)
220 {
221 return opPtr;
222 }
223
224 template<typename T, typename F>
225 auto GetTheRightOp(T(F::*opPtr)(T *, double *) const) -> decltype(opPtr)
226 {
227 return opPtr;
228 }
229 }
230}
231
232
233class TF1 : public TNamed, public TAttLine, public TAttFill, public TAttMarker {
234
235 template<class Func>
237
238public:
239 /// Add to list behavior
240 enum class EAddToList {
241 kDefault,
242 kAdd,
243 kNo
244 };
245
246
249 virtual TF1FunctorPointer * Clone() const = 0;
250 };
251
252protected:
253
254 enum EFType {
255 kFormula = 0, ///< Formula functions which can be stored,
256 kPtrScalarFreeFcn, ///< Pointer to scalar free function,
257 kInterpreted, ///< Interpreted functions constructed by name,
258 kTemplVec, ///< Vectorized free functions or TemplScalar functors evaluating on vectorized parameters,
259 kTemplScalar, ///< TemplScalar functors evaluating on scalar parameters
261 }; // formula based on composition class (e.g. NSUM, CONV)
262
263 Double_t fXmin{-1111}; ///< Lower bounds for the range
264 Double_t fXmax{-1111}; ///< Upper bounds for the range
265 Int_t fNpar{}; ///< Number of parameters
266 Int_t fNdim{}; ///< Function dimension
267 Int_t fNpx{100}; ///< Number of points used for the graphical representation
269 Int_t fNpfits{}; ///< Number of points used in the fit
270 Int_t fNDF{}; ///< Number of degrees of freedom in the fit
271 Double_t fChisquare{}; ///< Function fit chisquare
272 Double_t fMinimum{-1111}; ///< Minimum value for plotting
273 Double_t fMaximum{-1111}; ///< Maximum value for plotting
274 std::vector<Double_t> fParErrors; ///< Array of errors of the fNpar parameters
275 std::vector<Double_t> fParMin; ///< Array of lower limits of the fNpar parameters
276 std::vector<Double_t> fParMax; ///< Array of upper limits of the fNpar parameters
277 std::vector<Double_t> fSave; ///< Array of fNsave function values
278 std::vector<Double_t> fIntegral; ///<! Integral of function binned on fNpx bins
279 std::vector<Double_t> fAlpha; ///<! Array alpha. for each bin in x the deconvolution r of fIntegral
280 std::vector<Double_t> fBeta; ///<! Array beta. is approximated by x = alpha +beta*r *gamma*r**2
281 std::vector<Double_t> fGamma; ///<! Array gamma.
282 TObject *fParent{nullptr}; ///<! Parent object hooking this function (if one)
283 TH1 *fHistogram{nullptr}; ///<! Pointer to histogram used for visualisation
284 std::unique_ptr<TMethodCall> fMethodCall; ///<! Pointer to MethodCall in case of interpreted function
285 Bool_t fNormalized{false}; ///< Normalization option (false by default)
286 Double_t fNormIntegral{}; ///< Integral of the function before being normalized
287 std::unique_ptr<TF1FunctorPointer> fFunctor; ///<! Functor object to wrap any C++ callable object
288 std::unique_ptr<TFormula> fFormula; ///< Pointer to TFormula in case when user define formula
289 std::unique_ptr<TF1Parameters> fParams; ///< Pointer to Function parameters object (exists only for not-formula functions)
290 std::unique_ptr<TF1AbsComposition> fComposition; ///< Pointer to composition (NSUM or CONV)
291
292 /// General constructor for TF1. Most of the other constructors delegate on it
293 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):
294 TNamed(name, name), TAttLine(), TAttFill(), TAttMarker(), fXmin(xmin), fXmax(xmax), fNpar(npar), fNdim(ndim),
295 fType(functionType), fParErrors(npar), fParMin(npar), fParMax(npar)
296 {
297 fParams.reset(params);
298 fFunctor.reset(functor);
299 DoInitialize(addToGlobList);
300 }
301
302private:
303 // NSUM parsing helper functions
304 void DefineNSUMTerm(TObjArray *newFuncs, TObjArray *coeffNames,
305 TString &fullFormula,
306 TString &formula, int termStart, int termEnd,
308 int TermCoeffLength(TString &term);
309
310protected:
311
312 template <class T>
315 TF1FunctorPointerImpl(const std::function<T(const T *f, const Double_t *param)> &func) : fImpl(func){};
317 TF1FunctorPointer * Clone() const override { return new TF1FunctorPointerImpl<T>(fImpl); }
319 };
320
321
322
323
324 static std::atomic<Bool_t> fgAbsValue; //use absolute value of function when computing integral
325 static Bool_t fgRejectPoint; //True if point must be rejected in a fit
326 static std::atomic<Bool_t> fgAddToGlobList; //True if we want to register the function in the global list
327 static TF1 *fgCurrent; //pointer to current function being processed
328
329
330 //void CreateFromFunctor(const char *name, Int_t npar, Int_t ndim = 1);
331 void DoInitialize(EAddToList addToGlobList);
332
334 // tabulate the cumulative function integral at fNpx points. Used by GetRandom
336
337 virtual Double_t GetMinMaxNDim(Double_t *x , Bool_t findmax, Double_t epsilon = 0, Int_t maxiter = 0) const;
338 virtual void GetRange(Double_t *xmin, Double_t *xmax) const;
340
341public:
342
343 // TF1 status bits
345 kNotGlobal = BIT(10), // don't register in global list of functions
346 kNotDraw = BIT(9) // don't draw the function when in a TH1
347 };
348
349 TF1();
350 TF1(const char *name, const char *formula, Double_t xmin = 0, Double_t xmax = 1, EAddToList addToGlobList = EAddToList::kDefault, bool vectorize = false);
351 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
352 TF1(const char *name, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault);
353 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);
354 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);
355
356 template <class T>
357 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):
358 TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
359 {
360 fType = std::is_same<T, double>::value ? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
361 }
362
363 ////////////////////////////////////////////////////////////////////////////////
364 /// Constructor using a pointer to function.
365 ///
366 /// \param[in] name object name
367 /// \param[in] fcn pointer to function
368 /// \param[in] xmin,xmax x axis limits
369 /// \param[in] npar is the number of free parameters used by the function
370 /// \param[in] ndim number of dimensions
371 /// \param[in] addToGlobList boolean marking if it should be added to global list
372 ///
373 /// This constructor creates a function of type C when invoked
374 /// with the normal C++ compiler.
375 ///
376 ///
377 /// \warning A function created with this constructor cannot be Cloned
378
379
380 template <class T>
381 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):
382 TF1(EFType::kTemplVec, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<T>(fcn))
383 {}
384
385 // Constructors using functors (compiled mode only)
386 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);
387
388 // Template constructors from any C++ callable object, defining the operator() (double * , double *)
389 // and returning a double.
390 // The class name is not needed when using compile code, while it is required when using
391 // interpreted code via the specialized constructor with void *.
392 // An instance of the C++ function class or its pointer can both be used. The former is reccomended when using
393 // C++ compiled code, but if CINT compatibility is needed, then a pointer to the function class must be used.
394 // xmin and xmax specify the plotting range, npar is the number of parameters.
395 // See the tutorial math/exampleFunctor.C for an example of using this constructor
396 template <typename Func>
397 TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim = 1, EAddToList addToGlobList = EAddToList::kDefault) :
398 TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList)
399 {
400 //actual fType set in TF1Builder
402 }
403
404 // backward compatible interface
405 template <typename Func>
406 TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *, EAddToList addToGlobList = EAddToList::kDefault) :
407 TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar))
408 {
410 }
411
412
413 // Template constructors from a pointer to any C++ class of type PtrObj with a specific member function of type
414 // MemFn.
415 // The member function must have the signature of (double * , double *) and returning a double.
416 // The class name and the method name are not needed when using compile code
417 // (the member function pointer is used in this case), while they are required when using interpreted
418 // code via the specialized constructor with void *.
419 // xmin and xmax specify the plotting range, npar is the number of parameters.
420 // See the tutorial math/exampleFunctor.C for an example of using this constructor
421 template <class PtrObj, typename MemFn>
422 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) :
423 TF1(EFType::kTemplScalar, name, xmin, xmax, npar, ndim, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
424 {}
425
426 // backward compatible interface
427 template <class PtrObj, typename MemFn>
428 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) :
429 TF1(EFType::kTemplScalar, name, xmin, xmax, npar, 1, addToGlobList, new TF1Parameters(npar), new TF1FunctorPointerImpl<double>(ROOT::Math::ParamFunctor(p, memFn)))
430 {}
431
432 TF1(const TF1 &f1);
433 TF1 &operator=(const TF1 &rhs);
434 ~TF1() override;
435 virtual void AddParameter(const TString &name, Double_t value)
436 {
437 if (fFormula) fFormula->AddParameter(name, value);
438 }
439 // virtual void AddParameters(const pair<TString,Double_t> *pairs, Int_t size) { fFormula->AddParameters(pairs,size); }
440 // virtual void AddVariable(const TString &name, Double_t value = 0) { if (fFormula) fFormula->AddVariable(name,value); }
441 // virtual void AddVariables(const TString *vars, Int_t size) { if (fFormula) fFormula->AddVariables(vars,size); }
444 void Browse(TBrowser *b) override;
445 void Copy(TObject &f1) const override;
446 TObject *Clone(const char *newname = nullptr) const override;
447 virtual Double_t Derivative(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
448 virtual Double_t Derivative2(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
449 virtual Double_t Derivative3(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
450 static Double_t DerivativeError();
451 Int_t DistancetoPrimitive(Int_t px, Int_t py) override;
452 void Draw(Option_t *option = "") override;
453 virtual TF1 *DrawCopy(Option_t *option = "") const;
454 virtual TObject *DrawDerivative(Option_t *option = "al"); // *MENU*
455 virtual TObject *DrawIntegral(Option_t *option = "al"); // *MENU*
456 virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option = "");
457 virtual Double_t Eval(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
458 //template <class T> T Eval(T x, T y = 0, T z = 0, T t = 0) const;
459 virtual Double_t EvalPar(const Double_t *x, const Double_t *params = nullptr);
460 template <class T> T EvalPar(const T *x, const Double_t *params = nullptr);
461 virtual Double_t operator()(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
462 template <class T> T operator()(const T *x, const Double_t *params = nullptr);
463 void ExecuteEvent(Int_t event, Int_t px, Int_t py) override;
464 virtual void FixParameter(Int_t ipar, Double_t value);
465 /// Return true if function has data in fSave buffer
466 Bool_t HasSave() const { return !fSave.empty(); }
468 {
469 return (fType == EFType::kTemplVec) || (fType == EFType::kFormula && fFormula && fFormula->IsVectorized());
470 }
471 /// Return the Chisquare after fitting. See ROOT::Fit::FitResult::Chi2()
473 {
474 return fChisquare;
475 }
476 virtual TH1 *GetHistogram() const;
478 {
480 }
482 {
483 return fFormula.get();
484 }
485 virtual const TFormula *GetFormula() const
486 {
487 return fFormula.get();
488 }
490 {
491 return (fFormula) ? fFormula->GetExpFormula(option) : TString();
492 }
493 virtual const TObject *GetLinearPart(Int_t i) const
494 {
495 return (fFormula) ? fFormula->GetLinearPart(i) : nullptr;
496 }
497 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;
498 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;
499 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;
500 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;
502 {
503 return fMaximum;
504 }
506 {
507 return fMinimum;
508 }
509 virtual Int_t GetNpar() const
510 {
511 return fNpar;
512 }
513 virtual Int_t GetNdim() const
514 {
515 return fNdim;
516 }
517 virtual Int_t GetNDF() const;
518 virtual Int_t GetNpx() const
519 {
520 return fNpx;
521 }
523 {
524 return fMethodCall.get();
525 }
526 virtual Int_t GetNumber() const
527 {
528 return (fFormula) ? fFormula->GetNumber() : 0;
529 }
530 virtual Int_t GetNumberFreeParameters() const;
531 virtual Int_t GetNumberFitPoints() const
532 {
533 return fNpfits;
534 }
535 char *GetObjectInfo(Int_t px, Int_t py) const override;
537 {
538 return fParent;
539 }
540 virtual Double_t GetParameter(Int_t ipar) const
541 {
542 return (fFormula) ? fFormula->GetParameter(ipar) : fParams->GetParameter(ipar);
543 }
544 virtual Double_t GetParameter(const TString &name) const
545 {
546 return (fFormula) ? fFormula->GetParameter(name) : fParams->GetParameter(name);
547 }
548 virtual Double_t *GetParameters() const
549 {
550 return (fFormula) ? fFormula->GetParameters() : const_cast<Double_t *>(fParams->GetParameters());
551 }
552 virtual void GetParameters(Double_t *params)
553 {
554 if (fFormula) fFormula->GetParameters(params);
555 else std::copy(fParams->ParamsVec().begin(), fParams->ParamsVec().end(), params);
556 }
557 virtual const char *GetParName(Int_t ipar) const
558 {
559 return (fFormula) ? fFormula->GetParName(ipar) : fParams->GetParName(ipar);
560 }
561 virtual Int_t GetParNumber(const char *name) const
562 {
563 return (fFormula) ? fFormula->GetParNumber(name) : fParams->GetParNumber(name);
564 }
565 virtual Double_t GetParError(Int_t ipar) const;
566 virtual const Double_t *GetParErrors() const
567 {
568 return fParErrors.data();
569 }
570 virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const;
571 virtual Double_t GetProb() const;
572 virtual Int_t GetQuantiles(Int_t n, Double_t *xp, const Double_t *p);
573 virtual Double_t GetRandom(TRandom * rng = nullptr, Option_t * opt = nullptr);
574 virtual Double_t GetRandom(Double_t xmin, Double_t xmax, TRandom * rng = nullptr, Option_t * opt = nullptr);
575 virtual void GetRange(Double_t &xmin, Double_t &xmax) const;
576 virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const;
577 virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const;
578 virtual Double_t GetSave(const Double_t *x);
579 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;
580 virtual Double_t GetXmin() const
581 {
582 return fXmin;
583 }
584 virtual Double_t GetXmax() const
585 {
586 return fXmax;
587 }
588 TAxis *GetXaxis() const ;
589 TAxis *GetYaxis() const ;
590 TAxis *GetZaxis() const ;
592 {
593 return (fFormula) ? fFormula->GetVariable(name) : 0;
594 }
595 virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps = 0.01);
596 template <class T>
597 T GradientPar(Int_t ipar, const T *x, Double_t eps = 0.01);
598 template <class T>
599 T GradientParTempl(Int_t ipar, const T *x, Double_t eps = 0.01);
600
601 virtual void GradientPar(const Double_t *x, Double_t *grad, Double_t eps = 0.01);
602 template <class T>
603 void GradientPar(const T *x, T *grad, Double_t eps = 0.01);
604 template <class T>
605 void GradientParTempl(const T *x, T *grad, Double_t eps = 0.01);
606
607 virtual void InitArgs(const Double_t *x, const Double_t *params);
608 static void InitStandardFunctions();
609 virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel = 1.e-12);
610 virtual Double_t IntegralOneDim(Double_t a, Double_t b, Double_t epsrel, Double_t epsabs, Double_t &err);
611 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);
612 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);
613 // virtual Double_t IntegralFast(const TGraph *g, Double_t a, Double_t b, Double_t *params = nullptr);
614 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);
615 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);
616 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)
617 {
618 return IntegralMultiple(n, a, b, maxpts, epsrel, epsrel, relerr, nfnevl, ifail);
619 }
620 virtual Double_t IntegralMultiple(Int_t n, const Double_t *a, const Double_t *b, Double_t epsrel, Double_t &relerr);
621 virtual Bool_t IsEvalNormalized() const
622 {
623 return fNormalized;
624 }
625 /// return kTRUE if the point is inside the function range
626 virtual Bool_t IsInside(const Double_t *x) const
627 {
628 return !((x[0] < fXmin) || (x[0] > fXmax));
629 }
630 virtual Bool_t IsLinear() const
631 {
632 return (fFormula) ? fFormula->IsLinear() : false;
633 }
634 virtual Bool_t IsValid() const;
635 void Print(Option_t *option = "") const override;
636 void Paint(Option_t *option = "") override;
637 virtual void ReleaseParameter(Int_t ipar);
639 void SavePrimitive(std::ostream &out, Option_t *option = "") override;
640 virtual void SetChisquare(Double_t chi2)
641 {
642 fChisquare = chi2;
643 }
644 virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar = nullptr);
645 template <class PtrObj, typename MemFn>
646 void SetFunction(PtrObj &p, MemFn memFn);
647 template <typename Func>
648 void SetFunction(Func f);
649 virtual void SetMaximum(Double_t maximum = -1111); // *MENU*
650 virtual void SetMinimum(Double_t minimum = -1111); // *MENU*
651 virtual void SetNDF(Int_t ndf);
652 virtual void SetNumberFitPoints(Int_t npfits)
653 {
654 fNpfits = npfits;
655 }
656 virtual void SetNormalized(Bool_t flag)
657 {
658 fNormalized = flag;
659 Update();
660 }
661 inline void SetNdim(Int_t ndim)
662 {
663 fNdim = ndim;
664 Update();
665 }
666 virtual void SetNpx(Int_t npx = 100); // *MENU*
667 virtual void SetParameter(Int_t param, Double_t value)
668 {
669 (fFormula) ? fFormula->SetParameter(param, value) : fParams->SetParameter(param, value);
670 Update();
671 }
672 virtual void SetParameter(const TString &name, Double_t value)
673 {
674 (fFormula) ? fFormula->SetParameter(name, value) : fParams->SetParameter(name, value);
675 Update();
676 }
677 virtual void SetParameters(const Double_t *params)
678 {
679 (fFormula) ? fFormula->SetParameters(params) : fParams->SetParameters(params);
680 Update();
681 }
682 /// Set parameter values.
683 /// NaN values will be skipped, meaning that the corresponding parameters will not be changed.
684 virtual void SetParameters(double p0, double p1 = TMath::QuietNaN(), double p2 = TMath::QuietNaN(),
685 double p3 = TMath::QuietNaN(), double p4 = TMath::QuietNaN(), double p5 = TMath::QuietNaN(),
686 double p6 = TMath::QuietNaN(), double p7 = TMath::QuietNaN(), double p8 = TMath::QuietNaN(),
687 double p9 = TMath::QuietNaN(), double p10 = TMath::QuietNaN())
688 {
689 // Note: this is not made a variadic template method because it would
690 // presumably break the context menu in the TBrowser. Also, probably this
691 // method should not be virtual, because if the user wants to change
692 // parameter setting behavior, the SetParameter() method can be
693 // overridden.
694 if (fFormula) fFormula->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
695 else fParams->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
696 Update();
697 } // *MENU*
698 virtual void SetParName(Int_t ipar, const char *name);
699 virtual void SetParNames(const char *name0 = "", const char *name1 = "", const char *name2 = "",
700 const char *name3 = "", const char *name4 = "", const char *name5 = "",
701 const char *name6 = "", const char *name7 = "", const char *name8 = "",
702 const char *name9 = "", const char *name10 = ""); // *MENU*
703 virtual void SetParError(Int_t ipar, Double_t error);
704 virtual void SetParErrors(const Double_t *errors);
705 virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax);
706 virtual void SetParent(TObject *p = nullptr)
707 {
708 fParent = p;
709 }
710 virtual void SetRange(Double_t xmin, Double_t xmax); // *MENU*
713 virtual void SetSavedPoint(Int_t point, Double_t value);
714 void SetTitle(const char *title = "") override; // *MENU*
715 virtual void SetVectorized(Bool_t vectorized)
716 {
718 fFormula->SetVectorized(vectorized);
719 else
720 Warning("SetVectorized", "Can only set vectorized flag on formula-based TF1");
721 }
722 virtual void Update();
723
724 static TF1 *GetCurrent();
725 static void AbsValue(Bool_t reject = kTRUE);
726 static void RejectPoint(Bool_t reject = kTRUE);
727 static Bool_t RejectedPoint();
728 static void SetCurrent(TF1 *f1);
729
730 //Moments
731 virtual Double_t Moment(Double_t n, Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001);
732 virtual Double_t CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001);
733 virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001)
734 {
735 return Moment(1, a, b, params, epsilon);
736 }
737 virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001)
738 {
739 return CentralMoment(2, a, b, params, epsilon);
740 }
741
742 //some useful static utility functions to compute sampling points for Integral
743 //static void CalcGaussLegendreSamplingPoints(TGraph *g, Double_t eps=3.0e-11);
744 //static TGraph *CalcGaussLegendreSamplingPoints(Int_t num=21, Double_t eps=3.0e-11);
745 static void CalcGaussLegendreSamplingPoints(Int_t num, Double_t *x, Double_t *w, Double_t eps = 3.0e-11);
746
747private:
748 template <class T>
749 T EvalParTempl(const T *data, const Double_t *params = nullptr);
750
751#ifdef R__HAS_VECCORE
752 inline double EvalParVec(const Double_t *data, const Double_t *params);
753#endif
754
755 ClassDefOverride(TF1, 12) // The Parametric 1-D function
756};
757
758namespace ROOT {
759 namespace Internal {
760
761 template<class Func>
762 void TF1Builder<Func>::Build(TF1 *f, Func func)
763 {
764 using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
765 f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
767 f->fParams = std::make_unique<TF1Parameters>(f->fNpar);
768 }
769
770 template<class Func>
771 void TF1Builder<Func *>::Build(TF1 *f, Func *func)
772 {
773 using Fnc_t = typename ROOT::Internal::GetFunctorType<decltype(ROOT::Internal::GetTheRightOp(&Func::operator()))>::type;
774 f->fType = std::is_same<Fnc_t, double>::value? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
776 f->fParams = std::make_unique<TF1Parameters>(f->fNpar);
777 }
778
779 /// TF1 building from a string
780 /// used to build a TFormula based on a lambda function
781 template<>
782 struct TF1Builder<const char *> {
783 static void Build(TF1 *f, const char *formula)
784 {
785 f->fType = TF1::EFType::kFormula;
786 f->fFormula = std::make_unique<TFormula>("tf1lambda", formula, f->fNdim, f->fNpar, false);
787 TString formulaExpression(formula);
788 Ssiz_t first = formulaExpression.Index("return") + 7;
789 Ssiz_t last = formulaExpression.Last(';');
790 TString title = formulaExpression(first, last - first);
791 f->SetTitle(title);
792 }
793 };
794
795 inline void
796 EvalParMultiDim(TF1 *func, Double_t *out, const Double_t *x, std::size_t size, std::size_t rows, Double_t *params)
797 {
798 for (size_t i = 0; i < rows; i++) {
799 out[i] = func->EvalPar(x + i * size, params);
800 }
801 }
802 }
803}
804
806{
807 return Eval(x, y, z, t);
808}
809
810template <class T>
811inline T TF1::operator()(const T *x, const Double_t *params)
812{
813 return EvalPar(x, params);
814}
815
816////////////////////////////////////////////////////////////////////////////////
817/// EvalPar for vectorized
818template <class T>
819T TF1::EvalPar(const T *x, const Double_t *params)
820{
822 return EvalParTempl(x, params);
823 } else if (fType == EFType::kFormula) {
824 return fFormula->EvalPar(x, params);
825 } else
826 return TF1::EvalPar((double *)x, params);
827}
828
829////////////////////////////////////////////////////////////////////////////////
830/// Eval for vectorized functions
831// template <class T>
832// T TF1::Eval(T x, T y, T z, T t) const
833// {
834// if (fType == EFType::kFormula)
835// return fFormula->Eval(x, y, z, t);
836
837// T xx[] = {x, y, z, t};
838// Double_t *pp = (Double_t *)fParams->GetParameters();
839// return ((TF1 *)this)->EvalPar(xx, pp);
840// }
841
842// Internal to TF1. Evaluates Templated interfaces
843template <class T>
844inline T TF1::EvalParTempl(const T *data, const Double_t *params)
845{
847 if (!params) params = (Double_t *)fParams->GetParameters();
848 if (fFunctor)
849 return ((TF1FunctorPointerImpl<T> *)fFunctor.get())->fImpl(data, params);
850
851 // this should throw an error
852 // we nned to implement a vectorized GetSave(x)
853 return TMath::SignalingNaN();
854}
855
856#ifdef R__HAS_VECCORE
857// Internal to TF1. Evaluates Vectorized TF1 on data of type Double_v
858inline double TF1::EvalParVec(const Double_t *data, const Double_t *params)
859{
860 assert(fType == EFType::kTemplVec);
861 std::vector<ROOT::Double_v> d(fNdim);
862 ROOT::Double_v res;
863
864 for(auto i=0; i<fNdim; i++) {
865 d[i] = ROOT::Double_v(data[i]);
866 }
867
868 if (fFunctor) {
869 res = ((TF1FunctorPointerImpl<ROOT::Double_v> *) fFunctor.get())->fImpl(d.data(), params);
870 } else {
871 // res = GetSave(x);
872 return TMath::SignalingNaN();
873 }
874 return vecCore::Get<ROOT::Double_v>(res, 0);
875}
876#endif
877
879{
881}
883{
885}
886
887template <typename Func>
889{
890 // set function from a generic C++ callable object
892 fFunctor = std::make_unique<TF1::TF1FunctorPointerImpl<double>>(ROOT::Math::ParamFunctor(f));
893}
894template <class PtrObj, typename MemFn>
895void TF1::SetFunction(PtrObj &p, MemFn memFn)
896{
897 // set from a pointer to a member function
899 fFunctor = std::make_unique<TF1::TF1FunctorPointerImpl<double>>(ROOT::Math::ParamFunctor(p, memFn));
900}
901
902template <class T>
903inline T TF1::GradientPar(Int_t ipar, const T *x, Double_t eps)
904{
906 return GradientParTempl<T>(ipar, x, eps);
907 } else
908 return GradientParTempl<Double_t>(ipar, (const Double_t *)x, eps);
909}
910
911template <class T>
912inline T TF1::GradientParTempl(Int_t ipar, const T *x, Double_t eps)
913{
914 if (GetNpar() == 0)
915 return 0;
916
917 if (eps < 1e-10 || eps > 1) {
918 Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
919 eps = 0.01;
920 }
921 Double_t h;
922 TF1 *func = (TF1 *)this;
923 Double_t *parameters = GetParameters();
924
925 // Copy parameters for thread safety
926 std::vector<Double_t> parametersCopy(parameters, parameters + GetNpar());
927 parameters = parametersCopy.data();
928
929 Double_t al, bl, h2;
930 T f1, f2, g1, g2, d0, d2;
931
932 ((TF1 *)this)->GetParLimits(ipar, al, bl);
933 if (al * bl != 0 && al >= bl) {
934 // this parameter is fixed
935 return 0;
936 }
937
938 // check if error has been computer (is not zero)
939 if (func->GetParError(ipar) != 0)
940 h = eps * func->GetParError(ipar);
941 else
942 h = eps;
943
944 // save original parameters
945 Double_t par0 = parameters[ipar];
946
947 parameters[ipar] = par0 + h;
948 f1 = func->EvalPar(x, parameters);
949 parameters[ipar] = par0 - h;
950 f2 = func->EvalPar(x, parameters);
951 parameters[ipar] = par0 + h / 2;
952 g1 = func->EvalPar(x, parameters);
953 parameters[ipar] = par0 - h / 2;
954 g2 = func->EvalPar(x, parameters);
955
956 // compute the central differences
957 h2 = 1 / (2. * h);
958 d0 = f1 - f2;
959 d2 = 2 * (g1 - g2);
960
961 T grad = h2 * (4 * d2 - d0) / 3.;
962
963 // restore original value
964 parameters[ipar] = par0;
965
966 return grad;
967}
968
969template <class T>
970inline void TF1::GradientPar(const T *x, T *grad, Double_t eps)
971{
973 GradientParTempl<T>(x, grad, eps);
974 } else
975 GradientParTempl<Double_t>((const Double_t *)x, (Double_t *)grad, eps);
976}
977
978template <class T>
979inline void TF1::GradientParTempl(const T *x, T *grad, Double_t eps)
980{
981 if (eps < 1e-10 || eps > 1) {
982 Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
983 eps = 0.01;
984 }
985
986 for (Int_t ipar = 0; ipar < GetNpar(); ipar++) {
987 grad[ipar] = GradientParTempl<T>(ipar, x, eps);
988 }
989}
990
991#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
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassDef(name, id)
Definition Rtypes.h:342
#define BIT(n)
Definition Rtypes.h:90
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
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 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:53
std::vector< Double_t > fParameters
Definition TF1.h:138
const char * GetParName(Int_t iparam) const
Definition TF1.h:99
const std::vector< double > & ParamsVec() const
Definition TF1.h:92
virtual ~TF1Parameters()
Definition TF1.h:77
Double_t GetParameter(const char *name) const
Definition TF1.h:84
Double_t GetParameter(Int_t iparam) const
Definition TF1.h:80
void SetParName(Int_t iparam, const char *name)
Definition TF1.h:122
TF1Parameters & operator=(const TF1Parameters &rhs)
Definition TF1.h:70
void SetParameter(Int_t iparam, Double_t value)
Definition TF1.h:106
void SetParameter(const char *name, Double_t value)
Definition TF1.h:118
TF1Parameters(Int_t npar)
Definition TF1.h:56
std::vector< std::string > fParNames
Definition TF1.h:139
TF1Parameters()
Definition TF1.h:55
bool CheckIndex(Int_t i) const
Definition TF1.h:133
void SetParNames(Args &&... args)
Set parameter names.
Definition TF1.h:156
void SetParameters(const Double_t *params)
Definition TF1.h:111
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:3846
TF1Parameters(const TF1Parameters &rhs)
Definition TF1.h:65
const Double_t * GetParameters() const
Definition TF1.h:88
1-Dim function class
Definition TF1.h:233
std::unique_ptr< TF1FunctorPointer > fFunctor
! Functor object to wrap any C++ callable object
Definition TF1.h:287
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:1823
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:1696
virtual Double_t GetXmax() const
Definition TF1.h:584
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:3479
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:3683
EAddToList
Add to list behavior.
Definition TF1.h:240
virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
Definition TF1.h:733
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:1113
virtual const Double_t * GetParErrors() const
Definition TF1.h:566
virtual Int_t GetNumber() const
Definition TF1.h:526
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:1889
std::vector< Double_t > fParErrors
Array of errors of the fNpar parameters.
Definition TF1.h:274
Int_t fNdim
Function dimension.
Definition TF1.h:266
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:3826
static void AbsValue(Bool_t reject=kTRUE)
Static function: set the fgAbsValue flag.
Definition TF1.cxx:984
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition TF1.cxx:1584
virtual const TFormula * GetFormula() const
Definition TF1.h:485
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1940
Int_t fNpar
Number of parameters.
Definition TF1.h:265
virtual Double_t GetParameter(const TString &name) const
Definition TF1.h:544
TAxis * GetYaxis() const
Get y axis of the function.
Definition TF1.cxx:2411
virtual void GetParameters(Double_t *params)
Definition TF1.h:552
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:1930
static std::atomic< Bool_t > fgAddToGlobList
Definition TF1.h:326
virtual Bool_t IsEvalNormalized() const
Definition TF1.h:621
virtual Double_t GetVariable(const TString &name)
Definition TF1.h:591
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:640
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:286
Double_t GetChisquare() const
Return the Chisquare after fitting. See ROOT::Fit::FitResult::Chi2()
Definition TF1.h:472
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:912
virtual TH1 * CreateHistogram()
Definition TF1.h:477
Double_t fXmin
Lower bounds for the range.
Definition TF1.h:263
std::unique_ptr< TMethodCall > fMethodCall
! Pointer to MethodCall in case of interpreted function
Definition TF1.h:284
virtual void Update()
Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated h...
Definition TF1.cxx:3619
virtual Double_t GetProb() const
Return the fit probability.
Definition TF1.cxx:1955
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:357
virtual Int_t GetQuantiles(Int_t n, Double_t *xp, const Double_t *p)
Compute Quantiles for density distribution of this function.
Definition TF1.cxx:1994
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:3528
virtual Double_t GetMaximumStored() const
Definition TF1.h:501
virtual Int_t GetNpar() const
Definition TF1.h:509
virtual TString GetExpFormula(Option_t *option="") const
Definition TF1.h:489
std::vector< Double_t > fBeta
! Array beta. is approximated by x = alpha +beta*r *gamma*r**2
Definition TF1.h:280
Int_t fNDF
Number of degrees of freedom in the fit.
Definition TF1.h:270
TH1 * fHistogram
! Pointer to histogram used for visualisation
Definition TF1.h:283
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:616
std::unique_ptr< TF1AbsComposition > fComposition
Pointer to composition (NSUM or CONV)
Definition TF1.h:290
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:3490
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:536
Int_t fNpfits
Number of points used in the fit.
Definition TF1.h:269
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:1178
static void SetCurrent(TF1 *f1)
Static function setting the current function.
Definition TF1.cxx:3343
void SetFunction(PtrObj &p, MemFn memFn)
Definition TF1.h:895
virtual void SetParent(TObject *p=nullptr)
Definition TF1.h:706
std::vector< Double_t > fAlpha
! Array alpha. for each bin in x the deconvolution r of fIntegral
Definition TF1.h:279
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:428
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:3558
virtual Int_t GetNumberFitPoints() const
Definition TF1.h:531
std::unique_ptr< TFormula > fFormula
Pointer to TFormula in case when user define formula.
Definition TF1.h:288
virtual void SetParNames(const char *name0="", const char *name1="", const char *name2="", const char *name3="", const char *name4="", const char *name5="", const char *name6="", const char *name7="", const char *name8="", const char *name9="", const char *name10="")
Set up to 10 parameter names.
Definition TF1.cxx:3463
static Double_t DerivativeError()
Static function returning the error of the last call to the of Derivative's functions.
Definition TF1.cxx:1277
std::vector< Double_t > fParMin
Array of lower limits of the fNpar parameters.
Definition TF1.h:275
static void InitStandardFunctions()
Create the basic function objects.
Definition TF1.cxx:2497
Double_t fMaximum
Maximum value for plotting.
Definition TF1.h:273
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:548
Double_t fMinimum
Minimum value for plotting.
Definition TF1.h:272
int TermCoeffLength(TString &term)
Definition TF1.cxx:924
static Bool_t fgRejectPoint
Definition TF1.h:325
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1005
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:652
virtual Double_t GetMinimumStored() const
Definition TF1.h:505
virtual void SetNormalized(Bool_t flag)
Definition TF1.h:656
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:942
virtual Int_t GetNumberFreeParameters() const
Return the number of free parameters.
Definition TF1.cxx:1900
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:3702
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:3739
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:422
Double_t fChisquare
Function fit chisquare.
Definition TF1.h:271
EStatusBits
Definition TF1.h:344
@ kNotGlobal
Definition TF1.h:345
@ kNotDraw
Definition TF1.h:346
virtual TFormula * GetFormula()
Definition TF1.h:481
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:406
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:522
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a function.
Definition TF1.cxx:1293
virtual Double_t operator()(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Definition TF1.h:805
EFType fType
Definition TF1.h:268
Bool_t fNormalized
Normalization option (false by default)
Definition TF1.h:285
virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
Definition TF1.h:737
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:435
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2281
void SetNdim(Int_t ndim)
Definition TF1.h:661
void Browse(TBrowser *b) override
Browse.
Definition TF1.cxx:993
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:557
~TF1() override
TF1 default destructor.
Definition TF1.cxx:953
static TF1 * fgCurrent
Definition TF1.h:327
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1468
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:293
Int_t fNpx
Number of points used for the graphical representation.
Definition TF1.h:267
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set lower and upper limits for parameter ipar.
Definition TF1.cxx:3507
void DoInitialize(EAddToList addToGlobList)
Common initialization of the TF1.
Definition TF1.cxx:802
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:1863
static TF1 * GetCurrent()
Static function returning the current function being processed.
Definition TF1.cxx:1569
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:1918
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TF1.cxx:1536
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:324
virtual Bool_t IsLinear() const
Definition TF1.h:630
virtual void SetParameters(double p0, double p1=TMath::QuietNaN(), double p2=TMath::QuietNaN(), double p3=TMath::QuietNaN(), double p4=TMath::QuietNaN(), double p5=TMath::QuietNaN(), double p6=TMath::QuietNaN(), double p7=TMath::QuietNaN(), double p8=TMath::QuietNaN(), double p9=TMath::QuietNaN(), double p10=TMath::QuietNaN())
Set parameter values.
Definition TF1.h:684
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:1363
std::vector< Double_t > fParMax
Array of upper limits of the fNpar parameters.
Definition TF1.h:276
bool IsVectorized()
Definition TF1.h:467
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:838
std::vector< Double_t > fSave
Array of fNsave function values.
Definition TF1.h:277
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition TF1.cxx:3692
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:883
std::vector< Double_t > fGamma
! Array gamma.
Definition TF1.h:281
TObject * fParent
! Parent object hooking this function (if one)
Definition TF1.h:282
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:1723
virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option="")
Draw function between xmin and xmax.
Definition TF1.cxx:1420
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:677
T EvalParTempl(const T *data, const Double_t *params=nullptr)
Eval for vectorized functions.
Definition TF1.h:844
virtual TObject * DrawIntegral(Option_t *option="al")
Draw integral of this function.
Definition TF1.cxx:1407
std::vector< Double_t > fIntegral
! Integral of function binned on fNpx bins
Definition TF1.h:278
virtual TObject * DrawDerivative(Option_t *option="al")
Draw derivative of this function.
Definition TF1.cxx:1385
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:1439
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:1614
std::unique_ptr< TF1Parameters > fParams
Pointer to Function parameters object (exists only for not-formula functions)
Definition TF1.h:289
virtual void SetParameter(const TString &name, Double_t value)
Definition TF1.h:672
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:667
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:1243
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:1064
EFType
Definition TF1.h:254
@ kCompositionFcn
Definition TF1.h:260
@ kFormula
Formula functions which can be stored,.
Definition TF1.h:255
@ kPtrScalarFreeFcn
Pointer to scalar free function,.
Definition TF1.h:256
@ kTemplScalar
TemplScalar functors evaluating on scalar parameters.
Definition TF1.h:259
@ kTemplVec
Vectorized free functions or TemplScalar functors evaluating on vectorized parameters,...
Definition TF1.h:258
@ kInterpreted
Interpreted functions constructed by name,.
Definition TF1.h:257
virtual void SetSavedPoint(Int_t point, Double_t value)
Restore value of function saved at point.
Definition TF1.cxx:3542
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:1557
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition TF1.h:626
virtual Int_t GetNpx() const
Definition TF1.h:518
Double_t fXmax
Upper bounds for the range.
Definition TF1.h:264
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:1655
virtual Int_t GetNdim() const
Definition TF1.h:513
virtual Double_t GetXmin() const
Definition TF1.h:580
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:847
virtual const TObject * GetLinearPart(Int_t i) const
Definition TF1.h:493
virtual void SetVectorized(Bool_t vectorized)
Definition TF1.h:715
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:397
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:540
virtual Int_t GetParNumber(const char *name) const
Definition TF1.h:561
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:381
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
Bool_t HasSave() const
Return true if function has data in fSave buffer.
Definition TF1.h:466
The Formula class.
Definition TFormula.h:89
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
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:991
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:931
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:2378
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.
void EvalParMultiDim(TF1 *func, Double_t *out, const Double_t *x, std::size_t size, std::size_t rows, Double_t *params)
Definition TF1.h:796
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:207
ParamFunctorTempl< double > ParamFunctor
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Double_t Double_v
Definition Types.h:55
Bool_t IsNaN(Double_t x)
Definition TMath.h:896
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:906
Double_t SignalingNaN()
Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
Definition TMath.h:914
Internal class used by TF1 for obtaining the type from a functor out of the set of valid operator() s...
Definition TF1.h:181
static void Build(TF1 *f, const char *formula)
Definition TF1.h:783
Internal class used by TF1 for defining template specialization for different TF1 constructors
Definition TF1.h:169
static void Build(TF1 *f, Func func)
Definition TF1.h:762
TF1FunctorPointer * Clone() const override
Definition TF1.h:317
ROOT::Math::ParamFunctorTempl< T > fImpl
Definition TF1.h:318
TF1FunctorPointerImpl(const std::function< T(const T *f, const Double_t *param)> &func)
Definition TF1.h:315
~TF1FunctorPointerImpl() override
Definition TF1.h:316
TF1FunctorPointerImpl(const ROOT::Math::ParamFunctorTempl< T > &func)
Definition TF1.h:314
virtual ~TF1FunctorPointer()
Definition TF1.h:248
virtual TF1FunctorPointer * Clone() const =0
th1 Draw()