ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TBinomialEfficiencyFitter.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Frank Filthaut, Rene Brun 30/05/2007
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2007, 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 
12 ////////////////////////////////////////////////////////////////////////////////
13 /** \class TBinomialEfficiencyFitter
14  \ingroup Hist
15  \brief Binomial fitter for the division of two histograms.
16 
17 Use when you need to calculate a selection's efficiency from two histograms,
18 one containing all entries, and one containing the subset of these entries
19 that pass the selection, and when you have a parametrization available for
20 the efficiency as a function of the variable(s) under consideration.
21 
22 A very common problem when estimating efficiencies is that of error estimation:
23 when no other information is available than the total number of events N and
24 the selected number n, the best estimate for the selection efficiency p is n/N.
25 Standard binomial statistics dictates that the uncertainty (this presupposes
26 sufficiently high statistics that an approximation by a normal distribution is
27 reasonable) on p, given N, is
28 \f[
29  \sqrt{\frac{p(1-p)}{N}}
30 \f]
31 However, when p is estimated as n/N, fluctuations from the true p to its
32 estimate become important, especially for low numbers of events, and giving
33 rise to biased results.
34 
35 When fitting a parametrized efficiency, these problems can largely be overcome,
36 as a hypothesized true efficiency is available by construction. Even so, simply
37 using the corresponding uncertainty still presupposes that Gaussian errors
38 yields a reasonable approximation. When using, instead of binned efficiency
39 histograms, the original numerator and denominator histograms, a binned maximum
40 likelihood can be constructed as the product of bin-by-bin binomial probabilities
41 to select n out of N events. Assuming that a correct parametrization of the
42 efficiency is provided, this construction in general yields less biased results
43 (and is much less sensitive to binning details).
44 
45 A generic use of this method is given below (note that the method works for 2D
46 and 3D histograms as well):
47 
48 ~~~~~~~~~~~~~~~{.cpp}
49  {
50  TH1* denominator; // denominator histogram
51  TH1* numerator; // corresponding numerator histogram
52  TF1* eff; // efficiency parametrization
53  .... // set step sizes and initial parameter
54  .... // values for the fit function
55  .... // possibly also set ranges, see TF1::SetRange()
56  TBinomialEfficiencyFitter* f = new TBinomialEfficiencyFitter(
57  numerator, denominator);
58  Int_t status = f->Fit(eff, "I");
59  if (status == 0) {
60  // if the fit was successful, display bin-by-bin efficiencies
61  // as well as the result of the fit
62  numerator->Sumw2();
63  TH1* hEff = dynamic_cast<TH1*>(numerator->Clone("heff"));
64  hEff->Divide(hEff, denominator, 1.0, 1.0, "B");
65  hEff->Draw("E");
66  eff->Draw("same");
67  }
68  }
69 ~~~~~~~~~~~~~~~
70 
71 Note that this method cannot be expected to yield reliable results when using
72 weighted histograms (because the likelihood computation will be incorrect).
73 
74 *///////////////////////////////////////////////////////////////////////////////////
75 
77 
78 #include "TMath.h"
79 #include "TPluginManager.h"
80 #include "TROOT.h"
81 #include "TH1.h"
82 #include "TF1.h"
83 #include "TF2.h"
84 #include "TF3.h"
85 #include "Fit/FitConfig.h"
86 #include "Fit/Fitter.h"
87 #include "TFitResult.h"
88 #include "Math/Functor.h"
89 
90 #include <limits>
91 
92 
94 
96 
97 
98 ////////////////////////////////////////////////////////////////////////////////
99 /// default constructor
100 
102  fNumerator = 0;
103  fDenominator = 0;
104  fFunction = 0;
105  fFitDone = kFALSE;
106  fAverage = kFALSE;
107  fRange = kFALSE;
108  fEpsilon = kDefaultEpsilon;
109  fFitter = 0;
110 }
111 
112 ////////////////////////////////////////////////////////////////////////////////
113 /// Constructor.
114 ///
115 /// Note that no objects are copied, so it is up to the user to ensure that the
116 /// histogram pointers remain valid.
117 ///
118 /// Both histograms need to be "consistent". This is not checked here, but in
119 /// TBinomialEfficiencyFitter::Fit().
120 
121 TBinomialEfficiencyFitter::TBinomialEfficiencyFitter(const TH1 *numerator, const TH1 *denominator) {
123  fFunction = 0;
124  fFitter = 0;
125  Set(numerator,denominator);
126 }
127 
128 ////////////////////////////////////////////////////////////////////////////////
129 /// destructor
130 
132  if (fFitter) delete fFitter;
133  fFitter = 0;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Initialize with a new set of inputs.
138 
139 void TBinomialEfficiencyFitter::Set(const TH1 *numerator, const TH1 *denominator)
140 {
141  fNumerator = (TH1*)numerator;
142  fDenominator = (TH1*)denominator;
143 
144  fFitDone = kFALSE;
145  fAverage = kFALSE;
146  fRange = kFALSE;
147 }
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 /// Set the required integration precision, see TF1::Integral()
151 
153 {
154  fEpsilon = epsilon;
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Provide access to the underlying fitter object.
159 /// This may be useful e.g. for the retrieval of additional information (such
160 /// as the output covariance matrix of the fit).
161 
163 {
164  if (!fFitter) fFitter = new ROOT::Fit::Fitter();
165  return fFitter;
166 
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 /// Carry out the fit of the given function to the given histograms.
171 ///
172 /// If option "I" is used, the fit function will be averaged over the
173 /// bin (the default is to evaluate it simply at the bin center).
174 ///
175 /// If option "R" is used, the fit range will be taken from the fit
176 /// function (the default is to use the entire histogram).
177 ///
178 /// If option "S" a TFitResult object is returned and it can be used to obtain
179 /// additional fit information, like covariance or correlation matrix.
180 ///
181 /// Note that all parameter values, limits, and step sizes are copied
182 /// from the input fit function f1 (so they should be set before calling
183 /// this method. This is particularly relevant for the step sizes, taken
184 /// to be the "error" set on input, as a null step size usually fixes the
185 /// corresponding parameter. That is protected against, but in such cases
186 /// an arbitrary starting step size will be used, and the reliability of
187 /// the fit should be questioned). If parameters are to be fixed, this
188 /// should be done by specifying non-null parameter limits, with lower
189 /// limits larger than upper limits.
190 ///
191 /// On output, f1 contains the fitted parameters and errors, as well as
192 /// the number of degrees of freedom, and the goodness-of-fit estimator
193 /// as given by S. Baker and R. Cousins, Nucl. Instr. Meth. A221 (1984) 437.
194 
196 {
197  TString opt = option;
198  opt.ToUpper();
199  fAverage = opt.Contains("I");
200  fRange = opt.Contains("R");
201  Bool_t verbose = opt.Contains("V");
202  Bool_t quiet = opt.Contains("Q");
203  Bool_t saveResult = opt.Contains("S");
204  if (!f1) return -1;
205  fFunction = (TF1*)f1;
206  Int_t i, npar;
207  npar = f1->GetNpar();
208  if (npar <= 0) {
209  Error("Fit", "function %s has illegal number of parameters = %d",
210  f1->GetName(), npar);
211  return -3;
212  }
213 
214  // Check that function has same dimension as histogram
215  if (!fNumerator || !fDenominator) {
216  Error("Fit","No numerator or denominator histograms set");
217  return -5;
218  }
219  if (f1->GetNdim() != fNumerator->GetDimension()) {
220  Error("Fit","function %s dimension, %d, does not match histogram dimension, %d",
221  f1->GetName(), f1->GetNdim(), fNumerator->GetDimension());
222  return -4;
223  }
224  // Check that the numbers of bins for the histograms match
226  (f1->GetNdim() > 1 && fNumerator->GetNbinsY() != fDenominator->GetNbinsY()) ||
227  (f1->GetNdim() > 2 && fNumerator->GetNbinsZ() != fDenominator->GetNbinsZ())) {
228  Error("Fit", "numerator and denominator histograms do not have identical numbers of bins");
229  return -6;
230  }
231 
232  // initialize the fitter
233 
234  if (!fFitter) {
235  fFitter = new ROOT::Fit::Fitter();
236  }
237 
238 
239  std::vector<ROOT::Fit::ParameterSettings> & parameters = fFitter->Config().ParamsSettings();
240  parameters.reserve(npar);
241  for (i = 0; i < npar; i++) {
242 
243  // assign an ARBITRARY starting error to ensure the parameter won't be fixed!
244  Double_t we = f1->GetParError(i);
245  if (we <= 0) we = 0.3*TMath::Abs(f1->GetParameter(i));
246  if (we == 0) we = 0.01;
247 
248  parameters.push_back(ROOT::Fit::ParameterSettings(f1->GetParName(i), f1->GetParameter(i), we) );
249 
250  Double_t plow, pup;
251  f1->GetParLimits(i,plow,pup);
252  if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
253  parameters.back().Fix();
254  }
255  else if (plow < pup ) {
256  parameters.back().SetLimits(plow,pup);
257  }
258  }
259 
260  // fcn must be set after setting the parameters
262  fFitter->SetFCN(static_cast<ROOT::Math::IMultiGenFunction&>(fcnFunction));
263 
264 
265  // in case default value of 1.0 is used
266  if (fFitter->Config().MinimizerOptions().ErrorDef() == 1.0 ) {
268  }
269 
270  if (verbose) {
272  }
273  else if (quiet) {
275  }
276 
277 
278 
279  // perform the actual fit
280 
281  fFitDone = kTRUE;
283  if ( !status && !quiet)
284  Warning("Fit","Abnormal termination of minimization.");
285 
286 
287  //Store fit results in fitFunction
288  const ROOT::Fit::FitResult & fitResult = fFitter->Result();
289  if (!fitResult.IsEmpty() ) {
290  // set in f1 the result of the fit
291  f1->SetNDF(fitResult.Ndf() );
292 
293  //f1->SetNumberFitPoints(...); // this is set in ComputeFCN
294 
295  f1->SetParameters( &(fitResult.Parameters().front()) );
296  if ( int( fitResult.Errors().size()) >= f1->GetNpar() )
297  f1->SetParErrors( &(fitResult.Errors().front()) );
298 
299  f1->SetChisquare(2.*fitResult.MinFcnValue()); // store goodness of fit (Baker&Cousins)
300  f1->SetNDF(f1->GetNumberFitPoints()- fitResult.NFreeParameters());
301  Info("result"," chi2 %f ndf %d ",2.*fitResult.MinFcnValue(), fitResult.Ndf() );
302 
303  }
304  // create a new result class if needed
305  if (saveResult) {
306  TFitResult* fr = new TFitResult(fitResult);
307  TString name = TString::Format("TBinomialEfficiencyFitter_result_of_%s",f1->GetName() );
308  fr->SetName(name); fr->SetTitle(name);
309  return TFitResultPtr(fr);
310  }
311  else {
312  return TFitResultPtr(fitResult.Status() );
313  }
314 
315 }
316 
317 ////////////////////////////////////////////////////////////////////////////////
318 /// Compute the likelihood.
319 
321 {
322  int nDim = fDenominator->GetDimension();
323 
324  int xlowbin = fDenominator->GetXaxis()->GetFirst();
325  int xhighbin = fDenominator->GetXaxis()->GetLast();
326  int ylowbin = 0, yhighbin = 0, zlowbin = 0, zhighbin = 0;
327  if (nDim > 1) {
328  ylowbin = fDenominator->GetYaxis()->GetFirst();
329  yhighbin = fDenominator->GetYaxis()->GetLast();
330  if (nDim > 2) {
331  zlowbin = fDenominator->GetZaxis()->GetFirst();
332  zhighbin = fDenominator->GetZaxis()->GetLast();
333  }
334  }
335 
336  fFunction->SetParameters(par);
337 
338  if (fRange) {
339  double xmin, xmax, ymin, ymax, zmin, zmax;
340 
341  // This way to ensure that a minimum range chosen exactly at a
342  // bin boundary is far from elegant, but is hopefully adequate.
343 
344  if (nDim == 1) {
345  fFunction->GetRange(xmin, xmax);
346  xlowbin = fDenominator->GetXaxis()->FindBin(xmin);
347  xhighbin = fDenominator->GetXaxis()->FindBin(xmax);
348  } else if (nDim == 2) {
349  fFunction->GetRange(xmin, ymin, xmax, ymax);
350  xlowbin = fDenominator->GetXaxis()->FindBin(xmin);
351  xhighbin = fDenominator->GetXaxis()->FindBin(xmax);
352  ylowbin = fDenominator->GetYaxis()->FindBin(ymin);
353  yhighbin = fDenominator->GetYaxis()->FindBin(ymax);
354  } else if (nDim == 3) {
355  fFunction->GetRange(xmin, ymin, zmin, xmax, ymax, zmax);
356  xlowbin = fDenominator->GetXaxis()->FindBin(xmin);
357  xhighbin = fDenominator->GetXaxis()->FindBin(xmax);
358  ylowbin = fDenominator->GetYaxis()->FindBin(ymin);
359  yhighbin = fDenominator->GetYaxis()->FindBin(ymax);
360  zlowbin = fDenominator->GetZaxis()->FindBin(zmin);
361  zhighbin = fDenominator->GetZaxis()->FindBin(zmax);
362  }
363  }
364 
365  // The coding below is perhaps somewhat awkward -- but it is done
366  // so that 1D, 2D, and 3D cases can be covered in the same loops.
367 
368  f = 0.;
369 
370  Int_t npoints = 0;
371  Double_t nmax = 0;
372  for (int xbin = xlowbin; xbin <= xhighbin; ++xbin) {
373 
374  // compute the bin edges
375  Double_t xlow = fDenominator->GetXaxis()->GetBinLowEdge(xbin);
376  Double_t xup = fDenominator->GetXaxis()->GetBinLowEdge(xbin+1);
377 
378  for (int ybin = ylowbin; ybin <= yhighbin; ++ybin) {
379 
380  // compute the bin edges (if applicable)
381  Double_t ylow = (nDim > 1) ? fDenominator->GetYaxis()->GetBinLowEdge(ybin) : 0;
382  Double_t yup = (nDim > 1) ? fDenominator->GetYaxis()->GetBinLowEdge(ybin+1) : 0;
383 
384  for (int zbin = zlowbin; zbin <= zhighbin; ++zbin) {
385 
386  // compute the bin edges (if applicable)
387  Double_t zlow = (nDim > 2) ? fDenominator->GetZaxis()->GetBinLowEdge(zbin) : 0;
388  Double_t zup = (nDim > 2) ? fDenominator->GetZaxis()->GetBinLowEdge(zbin+1) : 0;
389 
390  int bin = fDenominator->GetBin(xbin, ybin, zbin);
391  Double_t nDen = fDenominator->GetBinContent(bin);
392  Double_t nNum = fNumerator->GetBinContent(bin);
393 
394  // count maximum value to use in the likelihood for inf
395  // i.e. a number much larger than the other terms
396  if (nDen> nmax) nmax = nDen;
397  if (nDen <= 0.) continue;
398  npoints++;
399 
400  // mu is the average of the function over the bin OR
401  // the function evaluated at the bin centre
402  // As yet, there is nothing to prevent mu from being
403  // outside the range <0,1> !!
404 
405  Double_t mu = 0;
406  switch (nDim) {
407  case 1:
408  mu = (fAverage) ?
409  fFunction->Integral(xlow, xup, fEpsilon)
410  / (xup-xlow) :
412  break;
413  case 2:
414  {
415  mu = (fAverage) ?
416  ((TF2*)fFunction)->Integral(xlow, xup, ylow, yup, fEpsilon)
417  / ((xup-xlow)*(yup-ylow)) :
420  }
421  break;
422  case 3:
423  {
424  mu = (fAverage) ?
425  ((TF3*)fFunction)->Integral(xlow, xup, ylow, yup, zlow, zup, fEpsilon)
426  / ((xup-xlow)*(yup-ylow)*(zup-zlow)) :
430  }
431  }
432 
433  // binomial formula (forgetting about the factorials)
434  if (nNum != 0.) {
435  if (mu > 0.)
436  f -= nNum * TMath::Log(mu*nDen/nNum);
437  else
438  f -= nmax * -1E30; // crossing our fingers
439  }
440  if (nDen - nNum != 0.) {
441  if (1. - mu > 0.)
442  f -= (nDen - nNum) * TMath::Log((1. - mu)*nDen/(nDen-nNum));
443  else
444  f -= nmax * -1E30; // crossing our fingers
445  }
446  }
447  }
448  }
449 
450  fFunction->SetNumberFitPoints(npoints);
451 }
452 
Double_t EvaluateFCN(const Double_t *par)
const FitConfig & Config() const
access to the fit configuration (const method)
Definition: Fitter.h:382
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:429
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
float xmin
Definition: THbookFile.cxx:93
void SetPrintLevel(int level)
set print level
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4629
Double_t Log(Double_t x)
Definition: TMath.h:526
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition: FitResult.h:174
const char Option_t
Definition: RtypesCore.h:62
float ymin
Definition: THbookFile.cxx:93
Documentation for class Functor class.
Definition: Functor.h:394
virtual Int_t GetNumberFitPoints() const
Definition: TF1.h:349
ClassImp(TBinomialEfficiencyFitter) TBinomialEfficiencyFitter
default constructor
Binomial fitter for the division of two histograms.
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:363
virtual ~TBinomialEfficiencyFitter()
destructor
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
virtual Int_t GetDimension() const
Definition: TH1.h:283
unsigned int Ndf() const
Number of degree of freedom.
Definition: FitResult.h:168
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
bool FitFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Fit using the a generic FCN function as a C++ callable object implementing double () (const double *)...
Definition: Fitter.h:538
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1088
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:489
Basic string class.
Definition: TString.h:137
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:421
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
double ErrorDef() const
error definition
bool Bool_t
Definition: RtypesCore.h:59
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1608
const Bool_t kFALSE
Definition: Rtypes.h:92
void ComputeFCN(Double_t &f, const Double_t *par)
Compute the likelihood.
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
ROOT::Math::MinimizerOptions & MinimizerOptions()
access to the minimizer control parameter (non const method)
Definition: FitConfig.h:138
void SetErrorDef(double err)
set error def
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TFile * f
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Definition: FitResult.h:179
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:2321
void Set(const TH1 *numerator, const TH1 *denominator)
Initialize with a new set of inputs.
Extends the ROOT::Fit::Result class with a TNamed inheritance providing easy possibility for I/O...
Definition: TFitResult.h:36
const FitResult & Result() const
get fit result
Definition: Fitter.h:354
const std::vector< ROOT::Fit::ParameterSettings > & ParamsSettings() const
get the vector of parameter settings (const method)
Definition: FitConfig.h:90
virtual Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition: TH1.cxx:4535
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition: TF1.cxx:1953
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1618
unsigned int NFreeParameters() const
get total number of free parameters
Definition: FitResult.h:139
void SetPrecision(Double_t epsilon)
Set the required integration precision, see TF1::Integral()
virtual void SetChisquare(Double_t chi2)
Definition: TF1.h:412
virtual Double_t GetBinCenter(Int_t bin) const
return bin center for 1D historam Better to use h1.GetXaxis().GetBinCenter(bin)
Definition: TH1.cxx:8470
float ymax
Definition: THbookFile.cxx:93
A 3-Dim function with parameters.
Definition: TF3.h:30
int Status() const
minimizer status code
Definition: FitResult.h:142
Fitter class, entry point for performing all type of fits.
Definition: Fitter.h:94
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:33
bool verbose
virtual Int_t GetNdim() const
Definition: TF1.h:343
Double_t E()
Definition: TMath.h:54
virtual Int_t GetNbinsZ() const
Definition: TH1.h:298
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
TAxis * GetYaxis()
Definition: TH1.h:320
float xmax
Definition: THbookFile.cxx:93
A 2-Dim function with parameters.
Definition: TF2.h:33
REAL epsilon
Definition: triangle.c:617
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition: TAxis.cxx:264
class containg the result of the fit and all the related information (fitted parameter values...
Definition: FitResult.h:52
double Double_t
Definition: RtypesCore.h:55
The TH1 histogram class.
Definition: TH1.h:80
#define name(a, b)
Definition: linkTestLib0.cpp:5
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:440
TAxis * GetZaxis()
Definition: TH1.h:321
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:352
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
double MinFcnValue() const
Return value of the objective function (chi2 or likelihood) used in the fit.
Definition: FitResult.h:125
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
bool SetFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Set a generic FCN function as a C++ callable object implementing double () (const double *) Note that...
Definition: Fitter.h:543
TFitResultPtr Fit(TF1 *f1, Option_t *option="")
Carry out the fit of the given function to the given histograms.
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:3169
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
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:1162
TF1 * f1
Definition: legend1.C:11
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:3102
ROOT::Fit::Fitter * GetFitter()
Provide access to the underlying fitter object.
const Bool_t kTRUE
Definition: Rtypes.h:91
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
Definition: FitResult.h:122
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual Int_t GetNpar() const
Definition: TF1.h:342
TAxis * GetXaxis()
Definition: TH1.h:319
const Double_t kDefaultEpsilon
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904