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 "TMatrixDSymfwd.h"
38#include "Math/Types.h"
39#include "Math/ParamFunctor.h"
40
41#include <string>
42
43class TF1;
44class TH1;
45class TAxis;
46class TRandom;
47
48namespace ROOT {
49 namespace Fit {
50 class FitResult;
51 }
52}
53
55public:
56 TF1Parameters() {} // needed for the I/O
58 fParameters(std::vector<Double_t>(npar)),
59 fParNames(std::vector<std::string>(npar))
60 {
61 for (int i = 0; i < npar; ++i) {
62 fParNames[i] = std::string(TString::Format("p%d", i).Data());
63 }
64 }
65 // copy constructor
70 // assignment
72 {
73 if (&rhs == this) return *this;
74 fParameters = rhs.fParameters;
75 fParNames = rhs.fParNames;
76 return *this;
77 }
78 virtual ~TF1Parameters() {}
79
80 // getter methods
82 {
83 return (CheckIndex(iparam)) ? fParameters[iparam] : 0;
84 }
85 Double_t GetParameter(const char *name) const
86 {
88 }
89 const Double_t *GetParameters() const
90 {
91 return fParameters.data();
92 }
93 const std::vector<double> &ParamsVec() const
94 {
95 return fParameters;
96 }
97
98 Int_t GetParNumber(const char *name) const;
99
100 const char *GetParName(Int_t iparam) const
101 {
102 return (CheckIndex(iparam)) ? fParNames[iparam].c_str() : "";
103 }
104
105
106 // setter methods
108 {
109 if (!CheckIndex(iparam)) return;
111 }
112 void SetParameters(const Double_t *params)
113 {
114 std::copy(params, params + fParameters.size(), fParameters.begin());
115 }
116 template <typename... Args>
117 void SetParameters(Double_t arg1, Args &&... args);
118
119 void SetParameter(const char *name, Double_t value)
120 {
122 }
123 void SetParName(Int_t iparam, const char *name)
124 {
125 if (!CheckIndex(iparam)) return;
126 fParNames[iparam] = std::string(name);
127 }
128 template <typename... Args>
129 void SetParNames(Args &&... args);
130
131 ClassDef(TF1Parameters, 1) // The Parameters of a parameteric function
132private:
133
134 bool CheckIndex(Int_t i) const
135 {
136 return (i >= 0 && i < int(fParameters.size()));
137 }
138
139 std::vector<Double_t> fParameters; // parameter values
140 std::vector<std::string> fParNames; // parameter names
141};
142
143/// Set parameter values.
144/// NaN values will be skipped, meaning that the corresponding parameters will not be changed.
145template <typename... Args>
147{
148 int i = 0;
149 for (double val : {arg1, static_cast<Double_t>(args)...}) {
150 if(!TMath::IsNaN(val)) SetParameter(i++, val);
151 }
152}
153
154/// Set parameter names.
155/// Empty strings will be skipped, meaning that the corresponding name will not be changed.
156template <typename... Args>
157void TF1Parameters::SetParNames(Args &&...args)
158{
159 int i = 0;
160 for (auto name : {static_cast<std::string const&>(args)...}) {
161 if(!name.empty()) SetParName(i++, name.c_str());
162 }
163}
164
165namespace ROOT {
166 namespace Internal {
167 /// %Internal class used by TF1 for defining
168 /// template specialization for different TF1 constructors
169 template<class Func>
170 struct TF1Builder {
171 static void Build(TF1 *f, Func func);
172 };
173
174 template<class Func>
175 struct TF1Builder<Func *> {
176 static void Build(TF1 *f, Func *func);
177 };
178 }
179}
180
181
182class TF1 : public TNamed, public TAttLine, public TAttFill, public TAttMarker {
183
184 template<class Func>
186
187public:
188 /// Add to list behavior
189 enum class EAddToList {
190 kDefault,
191 kAdd,
192 kNo
193 };
194
195
198 virtual TF1FunctorPointer * Clone() const = 0;
199 };
200
201protected:
202
203 enum EFType {
204 kFormula = 0, ///< Formula functions which can be stored,
205 kPtrScalarFreeFcn, ///< Pointer to scalar free function,
206 kInterpreted, ///< Interpreted functions constructed by name,
207 kTemplVec, ///< Vectorized free functions or TemplScalar functors evaluating on vectorized parameters,
208 kTemplScalar, ///< TemplScalar functors evaluating on scalar parameters
210 }; // formula based on composition class (e.g. NSUM, CONV)
211
212 Double_t fXmin{-1111}; ///< Lower bounds for the range
213 Double_t fXmax{-1111}; ///< Upper bounds for the range
214 Int_t fNpar{}; ///< Number of parameters
215 Int_t fNdim{}; ///< Function dimension
216 Int_t fNpx{100}; ///< Number of points used for the graphical representation
218 Int_t fNpfits{}; ///< Number of points used in the fit
219 Int_t fNDF{}; ///< Number of degrees of freedom in the fit
220 Double_t fChisquare{}; ///< Function fit chisquare
221 Double_t fMinimum{-1111}; ///< Minimum value for plotting
222 Double_t fMaximum{-1111}; ///< Maximum value for plotting
223 std::vector<Double_t> fParErrors; ///< Array of errors of the fNpar parameters
224 std::vector<Double_t> fParMin; ///< Array of lower limits of the fNpar parameters
225 std::vector<Double_t> fParMax; ///< Array of upper limits of the fNpar parameters
226 std::vector<Double_t> fSave; ///< Array of fNsave function values
227 std::vector<Double_t> fIntegral; ///<! Integral of function binned on fNpx bins
228 std::vector<Double_t> fAlpha; ///<! Array alpha. for each bin in x the deconvolution r of fIntegral
229 std::vector<Double_t> fBeta; ///<! Array beta. is approximated by x = alpha +beta*r *gamma*r**2
230 std::vector<Double_t> fGamma; ///<! Array gamma.
231 TObject *fParent{nullptr}; ///<! Parent object hooking this function (if one)
232 TH1 *fHistogram{nullptr}; ///<! Pointer to histogram used for visualisation
233 std::unique_ptr<TMethodCall> fMethodCall; ///<! Pointer to MethodCall in case of interpreted function
234 Bool_t fNormalized{false}; ///< Normalization option (false by default)
235 Double_t fNormIntegral{}; ///< Integral of the function before being normalized
236 std::unique_ptr<TF1FunctorPointer> fFunctor; ///<! Functor object to wrap any C++ callable object
237 std::unique_ptr<TFormula> fFormula; ///< Pointer to TFormula in case when user define formula
238 std::unique_ptr<TF1Parameters> fParams; ///< Pointer to Function parameters object (exists only for not-formula functions)
239 std::unique_ptr<TF1AbsComposition> fComposition; ///< Pointer to composition (NSUM or CONV)
240
241 /// General constructor for TF1. Most of the other constructors delegate on it
242 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):
245 {
246 fParams.reset(params);
247 fFunctor.reset(functor);
249 }
250
251private:
252 // NSUM parsing helper functions
255 TString &formula, int termStart, int termEnd,
258
259protected:
260
261 template <class T>
269
270
271
272
273 static std::atomic<Bool_t> fgAbsValue; //use absolute value of function when computing integral
274 static Bool_t fgRejectPoint; //True if point must be rejected in a fit
275 static std::atomic<Bool_t> fgAddToGlobList; //True if we want to register the function in the global list
276 static TF1 *fgCurrent; //pointer to current function being processed
277
278
279 //void CreateFromFunctor(const char *name, Int_t npar, Int_t ndim = 1);
281
283 // tabulate the cumulative function integral at fNpx points. Used by GetRandom
285
286 virtual Double_t GetMinMaxNDim(Double_t *x , Bool_t findmax, Double_t epsilon = 0, Int_t maxiter = 0) const;
287 virtual void GetRange(Double_t *xmin, Double_t *xmax) const;
289
291
292public:
293
294 // TF1 status bits
296 kNotGlobal = BIT(10), // don't register in global list of functions
297 kNotDraw = BIT(9) // don't draw the function when in a TH1
298 };
299
300 TF1();
301 TF1(const char *name, const char *formula, Double_t xmin = 0, Double_t xmax = 1, EAddToList addToGlobList = EAddToList::kDefault, bool vectorize = false);
302 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
305 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);
306
307 template <class T>
308 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):
310 {
311 fType = std::is_same<T, double>::value ? TF1::EFType::kTemplScalar : TF1::EFType::kTemplVec;
312 }
313
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Constructor using a pointer to function.
316 ///
317 /// \param[in] name object name
318 /// \param[in] fcn pointer to function
319 /// \param[in] xmin,xmax x axis limits
320 /// \param[in] npar is the number of free parameters used by the function
321 /// \param[in] ndim number of dimensions
322 /// \param[in] addToGlobList boolean marking if it should be added to global list
323 ///
324 /// This constructor creates a function of type C when invoked
325 /// with the normal C++ compiler.
326 ///
327 ///
328 /// \warning A function created with this constructor cannot be Cloned
329
330
331 template <class T>
335
336 // Constructors using functors (compiled mode only)
338
339 /// Template constructors from any C++ callable object, defining the operator() (double * , double *)
340 /// and returning a double.
341 /// The class name is not needed when using compile code, while it is required when using
342 /// interpreted code via the specialized constructor with void *.
343 /// An instance of the C++ function class or its pointer can both be used. The former is recommended when using
344 /// C++ compiled code, but if CINT-compatibility is needed, then a pointer to the function class must be used.
345 /// xmin and xmax specify the plotting range, npar is the number of parameters.
346 /// See the tutorial math/exampleFunctor.C for an example of using this constructor
347 template <typename Func>
350 {
351 //actual fType set in TF1Builder
353 }
354
355 /// Template constructors from a pointer to any C++ class of type PtrObj with a specific member function of type
356 /// MemFn.
357 /// The member function must have the signature of (double * , double *) and returning a double.
358 /// The class name and the method name are not needed when using compile code
359 /// (the member function pointer is used in this case), while they are required when using interpreted
360 /// code via the specialized constructor with void *.
361 /// xmin and xmax specify the plotting range, npar is the number of parameters.
362 /// See the tutorial math/exampleFunctor.C for an example of using this constructor
363 template <class PtrObj, typename MemFn>
367
368 TF1(const TF1 &f1);
369 TF1 &operator=(const TF1 &rhs);
370 ~TF1() override;
371 virtual void AddParameter(const TString &name, Double_t value)
372 {
373 if (fFormula) fFormula->AddParameter(name, value);
374 }
375 // virtual void AddParameters(const pair<TString,Double_t> *pairs, Int_t size) { fFormula->AddParameters(pairs,size); }
376 // virtual void AddVariable(const TString &name, Double_t value = 0) { if (fFormula) fFormula->AddVariable(name,value); }
377 // virtual void AddVariables(const TString *vars, Int_t size) { if (fFormula) fFormula->AddVariables(vars,size); }
380 void Browse(TBrowser *b) override;
381 void Copy(TObject &f1) const override;
382 TObject *Clone(const char *newname = nullptr) const override;
383 virtual Double_t Derivative(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
384 virtual Double_t Derivative2(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
385 virtual Double_t Derivative3(Double_t x, Double_t *params = nullptr, Double_t epsilon = 0.001) const;
386 static Double_t DerivativeError();
387 Int_t DistancetoPrimitive(Int_t px, Int_t py) override;
388 void Draw(Option_t *option = "") override;
389 virtual TF1 *DrawCopy(Option_t *option = "") const;
390 virtual TObject *DrawDerivative(Option_t *option = "al"); // *MENU*
391 virtual TObject *DrawIntegral(Option_t *option = "al"); // *MENU*
392 virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option = "");
393 virtual Double_t Eval(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
394 //template <class T> T Eval(T x, T y = 0, T z = 0, T t = 0) const;
395 virtual Double_t EvalPar(const Double_t *x, const Double_t *params = nullptr);
396 template <class T> T EvalPar(const T *x, const Double_t *params = nullptr);
397 virtual Double_t operator()(Double_t x, Double_t y = 0, Double_t z = 0, Double_t t = 0) const;
398 template <class T> T operator()(const T *x, const Double_t *params = nullptr);
399 Double_t EvalUncertainty(Double_t x, const TMatrixDSym *covMatrix = nullptr) const;
400 void ExecuteEvent(Int_t event, Int_t px, Int_t py) override;
401 virtual void FixParameter(Int_t ipar, Double_t value);
402 /// Return true if function has data in fSave buffer
403 Bool_t HasSave() const { return !fSave.empty(); }
405 {
406 return (fType == EFType::kTemplVec) || (fType == EFType::kFormula && fFormula && fFormula->IsVectorized());
407 }
408 /// Return the Chisquare after fitting. See ROOT::Fit::FitResult::Chi2()
410 {
411 return fChisquare;
412 }
413 virtual TH1 *GetHistogram() const;
415 {
417 }
419 {
420 return fFormula.get();
421 }
422 virtual const TFormula *GetFormula() const
423 {
424 return fFormula.get();
425 }
427 {
428 return (fFormula) ? fFormula->GetExpFormula(option) : TString();
429 }
430 virtual const TObject *GetLinearPart(Int_t i) const
431 {
432 return (fFormula) ? fFormula->GetLinearPart(i) : nullptr;
433 }
434 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;
435 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;
436 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;
437 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;
439 {
440 return fMaximum;
441 }
443 {
444 return fMinimum;
445 }
446 virtual Int_t GetNpar() const
447 {
448 return fNpar;
449 }
450 virtual Int_t GetNdim() const
451 {
452 return fNdim;
453 }
454 virtual Int_t GetNDF() const;
455 virtual Int_t GetNpx() const
456 {
457 return fNpx;
458 }
460 {
461 return fMethodCall.get();
462 }
463 virtual Int_t GetNumber() const
464 {
465 return (fFormula) ? fFormula->GetNumber() : 0;
466 }
467 virtual Int_t GetNumberFreeParameters() const;
468 virtual Int_t GetNumberFitPoints() const
469 {
470 return fNpfits;
471 }
472 char *GetObjectInfo(Int_t px, Int_t py) const override;
474 {
475 return fParent;
476 }
477 virtual Double_t GetParameter(Int_t ipar) const
478 {
479 return (fFormula) ? fFormula->GetParameter(ipar) : fParams->GetParameter(ipar);
480 }
481 virtual Double_t GetParameter(const TString &name) const
482 {
483 return (fFormula) ? fFormula->GetParameter(name) : fParams->GetParameter(name);
484 }
485 virtual Double_t *GetParameters() const
486 {
487 return (fFormula) ? fFormula->GetParameters() : const_cast<Double_t *>(fParams->GetParameters());
488 }
489 virtual void GetParameters(Double_t *params)
490 {
491 if (fFormula) fFormula->GetParameters(params);
492 else std::copy(fParams->ParamsVec().begin(), fParams->ParamsVec().end(), params);
493 }
494 virtual const char *GetParName(Int_t ipar) const
495 {
496 return (fFormula) ? fFormula->GetParName(ipar) : fParams->GetParName(ipar);
497 }
498 virtual Int_t GetParNumber(const char *name) const
499 {
500 return (fFormula) ? fFormula->GetParNumber(name) : fParams->GetParNumber(name);
501 }
502 virtual Double_t GetParError(Int_t ipar) const;
503 virtual Double_t GetParError(const char *name) const
504 {
506 }
507 virtual const Double_t *GetParErrors() const
508 {
509 return fParErrors.data();
510 }
511 virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const;
512 virtual Double_t GetProb() const;
513 virtual Int_t GetQuantiles(Int_t n, Double_t *xp, const Double_t *p);
514 virtual Double_t GetRandom(TRandom * rng = nullptr, Option_t * opt = nullptr);
515 virtual Double_t GetRandom(Double_t xmin, Double_t xmax, TRandom * rng = nullptr, Option_t * opt = nullptr);
516 virtual void GetRange(Double_t &xmin, Double_t &xmax) const;
517 virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const;
518 virtual void GetRange(Double_t &xmin, Double_t &ymin, Double_t &zmin, Double_t &xmax, Double_t &ymax, Double_t &zmax) const;
519 virtual Double_t GetSave(const Double_t *x);
520 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;
521 virtual Double_t GetXmin() const
522 {
523 return fXmin;
524 }
525 virtual Double_t GetXmax() const
526 {
527 return fXmax;
528 }
529 TAxis *GetXaxis() const ;
530 TAxis *GetYaxis() const ;
531 TAxis *GetZaxis() const ;
533 {
534 return (fFormula) ? fFormula->GetVariable(name) : 0;
535 }
536 virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps = 0.01) const;
537 template <class T>
538 T GradientPar(Int_t ipar, const T *x, Double_t eps = 0.01) const;
539 template <class T>
540 T GradientParTempl(Int_t ipar, const T *x, Double_t eps = 0.01) const;
541
542 virtual void GradientPar(const Double_t *x, Double_t *grad, Double_t eps = 0.01) const;
543 template <class T>
544 void GradientPar(const T *x, T *grad, Double_t eps = 0.01) const;
545 template <class T>
546 void GradientParTempl(const T *x, T *grad, Double_t eps = 0.01) const;
547
548 virtual void InitArgs(const Double_t *x, const Double_t *params);
549 static void InitStandardFunctions();
550 virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel = 1.e-12);
552 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);
553 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);
554 // virtual Double_t IntegralFast(const TGraph *g, Double_t a, Double_t b, Double_t *params = nullptr);
555 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);
562 virtual Bool_t IsEvalNormalized() const
563 {
564 return fNormalized;
565 }
566 /// return kTRUE if the point is inside the function range
567 virtual Bool_t IsInside(const Double_t *x) const
568 {
569 return !((x[0] < fXmin) || (x[0] > fXmax));
570 }
571 virtual Bool_t IsLinear() const
572 {
573 return (fFormula) ? fFormula->IsLinear() : false;
574 }
575 virtual Bool_t IsValid() const;
576 void Print(Option_t *option = "") const override;
577 void Paint(Option_t *option = "") override;
578 virtual void ReleaseParameter(Int_t ipar);
580 void SavePrimitive(std::ostream &out, Option_t *option = "") override;
582 {
584 }
585 virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar = nullptr);
586 template <class PtrObj, typename MemFn>
588 template <typename Func>
589 void SetFunction(Func f);
590 virtual void SetMaximum(Double_t maximum = -1111); // *MENU*
591 virtual void SetMinimum(Double_t minimum = -1111); // *MENU*
592 virtual void SetNDF(Int_t ndf);
594 {
595 fNpfits = npfits;
596 }
598 {
600 Update();
601 }
602 inline void SetNdim(Int_t ndim)
603 {
604 fNdim = ndim;
605 Update();
606 }
607 virtual void SetNpx(Int_t npx = 100); // *MENU*
608 virtual void SetParameter(Int_t param, Double_t value)
609 {
610 (fFormula) ? fFormula->SetParameter(param, value) : fParams->SetParameter(param, value);
611 Update();
612 }
613 virtual void SetParameter(const TString &name, Double_t value)
614 {
615 (fFormula) ? fFormula->SetParameter(name, value) : fParams->SetParameter(name, value);
616 Update();
617 }
618 virtual void SetParameters(const Double_t *params)
619 {
620 (fFormula) ? fFormula->SetParameters(params) : fParams->SetParameters(params);
621 Update();
622 }
623 /// Set parameter values.
624 /// NaN values will be skipped, meaning that the corresponding parameters will not be changed.
625 virtual void SetParameters(double p0, double p1 = TMath::QuietNaN(), double p2 = TMath::QuietNaN(),
626 double p3 = TMath::QuietNaN(), double p4 = TMath::QuietNaN(), double p5 = TMath::QuietNaN(),
627 double p6 = TMath::QuietNaN(), double p7 = TMath::QuietNaN(), double p8 = TMath::QuietNaN(),
628 double p9 = TMath::QuietNaN(), double p10 = TMath::QuietNaN())
629 {
630 // Note: this is not made a variadic template method because it would
631 // presumably break the context menu in the TBrowser. Also, probably this
632 // method should not be virtual, because if the user wants to change
633 // parameter setting behavior, the SetParameter() method can be
634 // overridden.
635 if (fFormula) fFormula->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
636 else fParams->SetParameters(p0, p1, p2, p3, p4, p5, p6, p7, p8, p9, p10);
637 Update();
638 } // *MENU*
639 virtual void SetParName(Int_t ipar, const char *name);
640 virtual void SetParNames(const char *name0 = "", const char *name1 = "", const char *name2 = "",
641 const char *name3 = "", const char *name4 = "", const char *name5 = "",
642 const char *name6 = "", const char *name7 = "", const char *name8 = "",
643 const char *name9 = "", const char *name10 = ""); // *MENU*
644 virtual void SetParError(Int_t ipar, Double_t error);
645 virtual void SetParErrors(const Double_t *errors);
646 virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax);
647 virtual void SetParent(TObject *p = nullptr)
648 {
649 fParent = p;
650 }
651 virtual void SetRange(Double_t xmin, Double_t xmax); // *MENU*
654 virtual void SetSavedPoint(Int_t point, Double_t value);
655 void SetTitle(const char *title = "") override; // *MENU*
657 {
659 fFormula->SetVectorized(vectorized);
660 else
661 Warning("SetVectorized", "Can only set vectorized flag on formula-based TF1");
662 }
663 virtual void Update();
664
665 static TF1 *GetCurrent();
666 static void AbsValue(Bool_t reject = kTRUE);
667 static void RejectPoint(Bool_t reject = kTRUE);
668 static Bool_t RejectedPoint();
669 static void SetCurrent(TF1 *f1);
670
671 //Moments
672 virtual Double_t Moment(Double_t n, Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001);
673 virtual Double_t CentralMoment(Double_t n, Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001);
674 virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001)
675 {
676 return Moment(1, a, b, params, epsilon);
677 }
678 virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params = nullptr, Double_t epsilon = 0.000001)
679 {
680 return CentralMoment(2, a, b, params, epsilon);
681 }
682
683 //some useful static utility functions to compute sampling points for Integral
684 //static void CalcGaussLegendreSamplingPoints(TGraph *g, Double_t eps=3.0e-11);
685 //static TGraph *CalcGaussLegendreSamplingPoints(Int_t num=21, Double_t eps=3.0e-11);
686 static void CalcGaussLegendreSamplingPoints(Int_t num, Double_t *x, Double_t *w, Double_t eps = 3.0e-11);
687
688private:
689 template <class T>
690 T EvalParTempl(const T *data, const Double_t *params = nullptr);
691
692 double EvalParVec(const Double_t *data, const Double_t *params);
693
694 ClassDefOverride(TF1, 12) // The Parametric 1-D function
695};
696
697namespace ROOT {
698 namespace Internal {
699
700 template<class Func>
701 void TF1Builder<Func>::Build(TF1 *f, Func func)
702 {
703 // check if vector interface is supported by Func
704#ifdef R__HAS_STD_EXPERIMENTAL_SIMD
705 if constexpr(std::is_invocable_r_v<Double_v, Func, Double_v*, double *>) {
706 f->fType = TF1::EFType::kTemplVec;
708 } else
709#endif
710 {
713 }
714
715 f->fParams = std::make_unique<TF1Parameters>(f->fNpar);
716 }
717
718 template<class Func>
719 void TF1Builder<Func *>::Build(TF1 *f, Func *func)
720 {
721 // check if vector interface is supported by Func
722#ifdef R__HAS_STD_EXPERIMENTAL_SIMD
723 if constexpr(std::is_invocable_r_v<Double_v, Func, Double_v*, double *>) {
724 f->fType = TF1::EFType::kTemplVec;
726 } else
727#endif
728 {
731 }
732
733 f->fParams = std::make_unique<TF1Parameters>(f->fNpar);
734 }
735
736 /// TF1 building from a string
737 /// used to build a TFormula based on a lambda function
738 template<>
739 struct TF1Builder<const char *> {
740 static void Build(TF1 *f, const char *formula)
741 {
742 f->fType = TF1::EFType::kFormula;
743 f->fFormula = std::make_unique<TFormula>("tf1lambda", formula, f->fNdim, f->fNpar, false);
744 TString formulaExpression(formula);
745 Ssiz_t first = formulaExpression.Index("return") + 7;
746 Ssiz_t last = formulaExpression.Last(';');
747 TString title = formulaExpression(first, last - first);
748 f->SetTitle(title);
749 }
750 };
751
752 inline void
753 EvalParMultiDim(TF1 *func, Double_t *out, const Double_t *x, std::size_t size, std::size_t rows, Double_t *params)
754 {
755 for (size_t i = 0; i < rows; i++) {
756 out[i] = func->EvalPar(x + i * size, params);
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
821
822template <typename Func>
824{
825 // set function from a generic C++ callable object
827 fFunctor = std::make_unique<TF1::TF1FunctorPointerImpl<double>>(ROOT::Math::ParamFunctor(f));
828}
829template <class PtrObj, typename MemFn>
831{
832 // set from a pointer to a member function
834 fFunctor = std::make_unique<TF1::TF1FunctorPointerImpl<double>>(ROOT::Math::ParamFunctor(p, memFn));
835}
836
837template <class T>
838inline T TF1::GradientPar(Int_t ipar, const T *x, Double_t eps) const
839{
841 return GradientParTempl<T>(ipar, x, eps);
842 } else
843 return GradientParTempl<Double_t>(ipar, (const Double_t *)x, eps);
844}
845
846template <class T>
847inline T TF1::GradientParTempl(Int_t ipar, const T *x, Double_t eps) const
848{
849 if (GetNpar() == 0)
850 return 0;
851
852 if (eps < 1e-10 || eps > 1) {
853 Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
854 eps = 0.01;
855 }
856 Double_t h;
857 TF1 *func = (TF1 *)this;
858 Double_t *parameters = GetParameters();
859
860 // Copy parameters for thread safety
861 std::vector<Double_t> parametersCopy(parameters, parameters + GetNpar());
862 parameters = parametersCopy.data();
863
864 Double_t al, bl, h2;
865 T f1, f2, g1, g2, d0, d2;
866
867 ((TF1 *)this)->GetParLimits(ipar, al, bl);
868 if (al * bl != 0 && al >= bl) {
869 // this parameter is fixed
870 return 0;
871 }
872
873 // check if error has been computer (is not zero)
874 if (func->GetParError(ipar) != 0)
875 h = eps * func->GetParError(ipar);
876 else
877 h = eps;
878
879 // save original parameters
880 Double_t par0 = parameters[ipar];
881
882 parameters[ipar] = par0 + h;
883 f1 = func->EvalPar(x, parameters);
884 parameters[ipar] = par0 - h;
885 f2 = func->EvalPar(x, parameters);
886 parameters[ipar] = par0 + h / 2;
887 g1 = func->EvalPar(x, parameters);
888 parameters[ipar] = par0 - h / 2;
889 g2 = func->EvalPar(x, parameters);
890
891 // compute the central differences
892 h2 = 1 / (2. * h);
893 d0 = f1 - f2;
894 d2 = 2 * (g1 - g2);
895
896 T grad = h2 * (4 * d2 - d0) / 3.;
897
898 // restore original value
899 parameters[ipar] = par0;
900
901 return grad;
902}
903
904template <class T>
905inline void TF1::GradientPar(const T *x, T *grad, Double_t eps) const
906{
908 GradientParTempl<T>(x, grad, eps);
909 } else
910 GradientParTempl<Double_t>((const Double_t *)x, (Double_t *)grad, eps);
911}
912
913template <class T>
914inline void TF1::GradientParTempl(const T *x, T *grad, Double_t eps) const
915{
916 if (eps < 1e-10 || eps > 1) {
917 Warning("Derivative", "parameter esp=%g out of allowed range[1e-10,1], reset to 0.01", eps);
918 eps = 0.01;
919 }
920
921 for (Int_t ipar = 0; ipar < GetNpar(); ipar++) {
922 grad[ipar] = GradientParTempl<T>(ipar, x, eps);
923 }
924}
925
926#endif
#define 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
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
#define ClassDef(name, id)
Definition Rtypes.h:344
#define BIT(n)
Definition Rtypes.h:91
#define ClassDefOverride(name, id)
Definition Rtypes.h:348
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
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:44
Param Functor class for Multidimensional functions.
Fill Area Attributes class.
Definition TAttFill.h:20
Line Attributes class.
Definition TAttLine.h:20
Marker Attributes class.
Definition TAttMarker.h:20
Class to manage histogram axis.
Definition TAxis.h:32
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
TF1 Parameters class.
Definition TF1.h:54
std::vector< Double_t > fParameters
Definition TF1.h:139
const char * GetParName(Int_t iparam) const
Definition TF1.h:100
const std::vector< double > & ParamsVec() const
Definition TF1.h:93
virtual ~TF1Parameters()
Definition TF1.h:78
Double_t GetParameter(const char *name) const
Definition TF1.h:85
Double_t GetParameter(Int_t iparam) const
Definition TF1.h:81
void SetParName(Int_t iparam, const char *name)
Definition TF1.h:123
TF1Parameters & operator=(const TF1Parameters &rhs)
Definition TF1.h:71
void SetParameter(Int_t iparam, Double_t value)
Definition TF1.h:107
void SetParameter(const char *name, Double_t value)
Definition TF1.h:119
TF1Parameters(Int_t npar)
Definition TF1.h:57
std::vector< std::string > fParNames
Definition TF1.h:140
TF1Parameters()
Definition TF1.h:56
bool CheckIndex(Int_t i) const
Definition TF1.h:134
void SetParNames(Args &&... args)
Set parameter names.
Definition TF1.h:157
void SetParameters(const Double_t *params)
Definition TF1.h:112
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:3886
TF1Parameters(const TF1Parameters &rhs)
Definition TF1.h:66
const Double_t * GetParameters() const
Definition TF1.h:89
1-Dim function class
Definition TF1.h:182
std::unique_ptr< TF1FunctorPointer > fFunctor
! Functor object to wrap any C++ callable object
Definition TF1.h:236
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:1873
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:1746
virtual Double_t GetXmax() const
Definition TF1.h:525
virtual void ReleaseParameter(Int_t ipar)
Release parameter number ipar during a fit operation.
Definition TF1.cxx:3201
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition TF1.cxx:3534
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:3723
EAddToList
Add to list behavior.
Definition TF1.h:189
virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
Definition TF1.h:674
T GradientParTempl(Int_t ipar, const T *x, Double_t eps=0.01) const
Definition TF1.h:847
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:1120
virtual const Double_t * GetParErrors() const
Definition TF1.h:507
virtual Int_t GetNumber() const
Definition TF1.h:463
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:1939
std::vector< Double_t > fParErrors
Array of errors of the fNpar parameters.
Definition TF1.h:223
Int_t fNdim
Function dimension.
Definition TF1.h:215
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:3866
static void AbsValue(Bool_t reject=kTRUE)
Static function: set the fgAbsValue flag.
Definition TF1.cxx:985
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition TF1.cxx:1634
virtual const TFormula * GetFormula() const
Definition TF1.h:422
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1990
Int_t fNpar
Number of parameters.
Definition TF1.h:214
virtual Double_t GetParameter(const TString &name) const
Definition TF1.h:481
TAxis * GetYaxis() const
Get y axis of the function.
Definition TF1.cxx:2459
virtual void GetParameters(Double_t *params)
Definition TF1.h:489
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:3474
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1980
static std::atomic< Bool_t > fgAddToGlobList
Definition TF1.h:275
virtual Bool_t IsEvalNormalized() const
Definition TF1.h:562
virtual Double_t GetVariable(const TString &name)
Definition TF1.h:532
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:2757
virtual void SetChisquare(Double_t chi2)
Definition TF1.h:581
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:2827
Double_t fNormIntegral
Integral of the function before being normalized.
Definition TF1.h:235
Double_t GetChisquare() const
Return the Chisquare after fitting. See ROOT::Fit::FitResult::Chi2()
Definition TF1.h:409
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:3449
void Print(Option_t *option="") const override
This method must be overridden when a class wants to print itself.
Definition TF1.cxx:2946
virtual TH1 * CreateHistogram()
Definition TF1.h:414
Double_t fXmin
Lower bounds for the range.
Definition TF1.h:212
std::unique_ptr< TMethodCall > fMethodCall
! Pointer to MethodCall in case of interpreted function
Definition TF1.h:233
virtual void Update()
Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated h...
Definition TF1.cxx:3674
virtual Double_t GetProb() const
Return the fit probability.
Definition TF1.cxx:2005
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:308
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:2042
TAxis * GetZaxis() const
Get z axis of the function. (In case this object is a TF2 or TF3)
Definition TF1.cxx:2470
virtual Double_t GetRandom(TRandom *rng=nullptr, Option_t *opt=nullptr)
Return a random number following this function shape.
Definition TF1.cxx:2240
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF1.cxx:3583
virtual Double_t GetMaximumStored() const
Definition TF1.h:438
virtual Int_t GetNpar() const
Definition TF1.h:446
virtual TString GetExpFormula(Option_t *option="") const
Definition TF1.h:426
std::vector< Double_t > fBeta
! Array beta. is approximated by x = alpha +beta*r *gamma*r**2
Definition TF1.h:229
Double_t EvalUncertainty(Double_t x, const TMatrixDSym *covMatrix=nullptr) const
Evaluate the uncertainty of the function at location x due to the parameter uncertainties.
Definition TF1.cxx:1567
TString ProvideSaveName(Option_t *option)
Provide variable name for function for saving as primitive When TH1 or TGraph stores list of function...
Definition TF1.cxx:3269
Int_t fNDF
Number of degrees of freedom in the fit.
Definition TF1.h:219
TH1 * fHistogram
! Pointer to histogram used for visualisation
Definition TF1.h:232
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:557
std::unique_ptr< TF1AbsComposition > fComposition
Pointer to composition (NSUM or CONV)
Definition TF1.h:239
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:3545
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:3097
TObject * GetParent() const
Definition TF1.h:473
Int_t fNpfits
Number of points used in the fit.
Definition TF1.h:218
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:1185
static void SetCurrent(TF1 *f1)
Static function setting the current function.
Definition TF1.cxx:3398
void SetFunction(PtrObj &p, MemFn memFn)
Definition TF1.h:830
virtual void SetParent(TObject *p=nullptr)
Definition TF1.h:647
std::vector< Double_t > fAlpha
! Array alpha. for each bin in x the deconvolution r of fIntegral
Definition TF1.h:228
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition TF1.cxx:2580
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:3613
virtual Int_t GetNumberFitPoints() const
Definition TF1.h:468
std::unique_ptr< TFormula > fFormula
Pointer to TFormula in case when user define formula.
Definition TF1.h:237
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:3518
static Double_t DerivativeError()
Static function returning the error of the last call to the of Derivative's functions.
Definition TF1.cxx:1284
std::vector< Double_t > fParMin
Array of lower limits of the fNpar parameters.
Definition TF1.h:224
static void InitStandardFunctions()
Create the basic function objects.
Definition TF1.cxx:2545
Double_t fMaximum
Maximum value for plotting.
Definition TF1.h:222
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition TF1.cxx:3488
virtual Double_t * GetParameters() const
Definition TF1.h:485
Double_t fMinimum
Minimum value for plotting.
Definition TF1.h:221
int TermCoeffLength(TString &term)
Definition TF1.cxx:925
static Bool_t fgRejectPoint
Definition TF1.h:274
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1006
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:593
double EvalParVec(const Double_t *data, const Double_t *params)
virtual Double_t GetMinimumStored() const
Definition TF1.h:442
virtual void SetNormalized(Bool_t flag)
Definition TF1.h:597
void Paint(Option_t *option="") override
Paint this function with its current attributes.
Definition TF1.cxx:3002
TF1 & operator=(const TF1 &rhs)
Operator =.
Definition TF1.cxx:943
virtual Int_t GetNumberFreeParameters() const
Return the number of free parameters.
Definition TF1.cxx:1950
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:3742
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:3779
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)
Template constructors from a pointer to any C++ class of type PtrObj with a specific member function ...
Definition TF1.h:364
Double_t fChisquare
Function fit chisquare.
Definition TF1.h:220
EStatusBits
Definition TF1.h:295
@ kNotGlobal
Definition TF1.h:296
@ kNotDraw
Definition TF1.h:297
virtual TFormula * GetFormula()
Definition TF1.h:418
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition TF1.cxx:2530
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:2900
TMethodCall * GetMethodCall() const
Definition TF1.h:459
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a function.
Definition TF1.cxx:1300
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:217
Bool_t fNormalized
Normalization option (false by default)
Definition TF1.h:234
virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
Definition TF1.h:678
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:3462
virtual void AddParameter(const TString &name, Double_t value)
Definition TF1.h:371
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2329
void SetNdim(Int_t ndim)
Definition TF1.h:602
void Browse(TBrowser *b) override
Browse.
Definition TF1.cxx:994
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:494
~TF1() override
TF1 default destructor.
Definition TF1.cxx:954
static TF1 * fgCurrent
Definition TF1.h:276
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1498
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:242
Int_t fNpx
Number of points used for the graphical representation.
Definition TF1.h:216
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set lower and upper limits for parameter ipar.
Definition TF1.cxx:3562
void DoInitialize(EAddToList addToGlobList)
Common initialization of the TF1.
Definition TF1.cxx:803
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:1913
static TF1 * GetCurrent()
Static function returning the current function being processed.
Definition TF1.cxx:1619
virtual void SetParName(Int_t ipar, const char *name)
Set name of parameter number ipar.
Definition TF1.cxx:3505
char * GetObjectInfo(Int_t px, Int_t py) const override
Redefines TObject::GetObjectInfo.
Definition TF1.cxx:1968
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TF1.cxx:1586
virtual Double_t GetSave(const Double_t *x)
Get value corresponding to X in array of fSave values.
Definition TF1.cxx:2392
static std::atomic< Bool_t > fgAbsValue
Definition TF1.h:273
virtual Bool_t IsLinear() const
Definition TF1.h:571
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:625
TF1()
TF1 default constructor.
Definition TF1.cxx:490
virtual TF1 * DrawCopy(Option_t *option="") const
Draw a copy of this function with its current attributes.
Definition TF1.cxx:1370
std::vector< Double_t > fParMax
Array of upper limits of the fNpar parameters.
Definition TF1.h:225
bool IsVectorized()
Definition TF1.h:404
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TF1.cxx:3286
virtual Bool_t IsValid() const
Return kTRUE if the function is valid.
Definition TF1.cxx:2931
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:839
std::vector< Double_t > fSave
Array of fNsave function values.
Definition TF1.h:226
virtual Double_t GetParError(const char *name) const
Definition TF1.h:503
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition TF1.cxx:3732
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:884
std::vector< Double_t > fGamma
! Array gamma.
Definition TF1.h:230
TObject * fParent
! Parent object hooking this function (if one)
Definition TF1.h:231
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:1773
virtual void DrawF1(Double_t xmin, Double_t xmax, Option_t *option="")
Draw function between xmin and xmax.
Definition TF1.cxx:1427
Bool_t ComputeCdfTable(Option_t *opt)
Compute the cumulative function at fNpx points between fXmin and fXmax.
Definition TF1.cxx:2129
virtual void SetParameters(const Double_t *params)
Definition TF1.h:618
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:1414
std::vector< Double_t > fIntegral
! Integral of function binned on fNpx bins
Definition TF1.h:227
virtual TObject * DrawDerivative(Option_t *option="al")
Draw derivative of this function.
Definition TF1.cxx:1392
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:1446
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:1664
std::unique_ptr< TF1Parameters > fParams
Pointer to Function parameters object (exists only for not-formula functions)
Definition TF1.h:238
virtual void SetParameter(const TString &name, Double_t value)
Definition TF1.h:613
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:608
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:1250
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:3211
TObject * Clone(const char *newname=nullptr) const override
Make a complete copy of the underlying object.
Definition TF1.cxx:1071
EFType
Definition TF1.h:203
@ kCompositionFcn
Definition TF1.h:209
@ kFormula
Formula functions which can be stored,.
Definition TF1.h:204
@ kPtrScalarFreeFcn
Pointer to scalar free function,.
Definition TF1.h:205
@ kTemplScalar
TemplScalar functors evaluating on scalar parameters.
Definition TF1.h:208
@ kTemplVec
Vectorized free functions or TemplScalar functors evaluating on vectorized parameters,...
Definition TF1.h:207
@ kInterpreted
Interpreted functions constructed by name,.
Definition TF1.h:206
virtual Double_t GradientPar(Int_t ipar, const Double_t *x, Double_t eps=0.01) const
Compute the gradient (derivative) wrt a parameter ipar.
Definition TF1.cxx:2493
virtual void SetSavedPoint(Int_t point, Double_t value)
Restore value of function saved at point.
Definition TF1.cxx:3597
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:1607
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition TF1.h:567
virtual Int_t GetNpx() const
Definition TF1.h:455
Double_t fXmax
Upper bounds for the range.
Definition TF1.h:213
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:1705
virtual Int_t GetNdim() const
Definition TF1.h:450
virtual Double_t GetXmin() const
Definition TF1.h:521
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:848
virtual const TObject * GetLinearPart(Int_t i) const
Definition TF1.h:430
virtual void SetVectorized(Bool_t vectorized)
Definition TF1.h:656
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:2670
TF1(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, Int_t ndim=1, EAddToList addToGlobList=EAddToList::kDefault)
Template constructors from any C++ callable object, defining the operator() (double * ,...
Definition TF1.h:348
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:477
virtual Int_t GetParNumber(const char *name) const
Definition TF1.h:498
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:332
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:3410
TAxis * GetXaxis() const
Get x axis of the function.
Definition TF1.cxx:2448
Bool_t HasSave() const
Return true if function has data in fSave buffer.
Definition TF1.h:403
The Formula class.
Definition TFormula.h:89
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
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:1074
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
Basic string class.
Definition TString.h:138
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:2384
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
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:753
ParamFunctorTempl< double > ParamFunctor
Bool_t IsNaN(Double_t x)
Definition TMath.h:903
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:913
Double_t SignalingNaN()
Returns a signaling NaN as defined by IEEE 754](http://en.wikipedia.org/wiki/NaN#Signaling_NaN).
Definition TMath.h:921
static void Build(TF1 *f, const char *formula)
Definition TF1.h:740
Internal class used by TF1 for defining template specialization for different TF1 constructors
Definition TF1.h:170
static void Build(TF1 *f, Func func)
Definition TF1.h:701
TF1FunctorPointer * Clone() const override
Definition TF1.h:266
ROOT::Math::ParamFunctorTempl< T > fImpl
Definition TF1.h:267
TF1FunctorPointerImpl(const std::function< T(const T *f, const Double_t *param)> &func)
Definition TF1.h:264
~TF1FunctorPointerImpl() override
Definition TF1.h:265
TF1FunctorPointerImpl(const ROOT::Math::ParamFunctorTempl< T > &func)
Definition TF1.h:263
virtual ~TF1FunctorPointer()
Definition TF1.h:197
virtual TF1FunctorPointer * Clone() const =0
th1 Draw()