Logo ROOT   6.12/07
Reference Guide
HFitImpl.cxx
Go to the documentation of this file.
1 // new HFit function
2 //______________________________________________________________________________
3 
4 
5 #include "TH1.h"
6 #include "TH2.h"
7 #include "TF1.h"
8 #include "TF2.h"
9 #include "TF3.h"
10 #include "TError.h"
11 #include "TGraph.h"
12 #include "TMultiGraph.h"
13 #include "TGraph2D.h"
14 #include "THnBase.h"
15 
16 #include "Fit/Fitter.h"
17 #include "Fit/FitConfig.h"
18 #include "Fit/BinData.h"
19 #include "Fit/UnBinData.h"
20 #include "Fit/Chi2FCN.h"
22 #include "HFitInterface.h"
23 #include "Math/MinimizerOptions.h"
24 #include "Math/Minimizer.h"
25 
26 #include "Math/WrappedTF1.h"
27 #include "Math/WrappedMultiTF1.h"
28 
29 #include "TList.h"
30 #include "TMath.h"
31 
32 #include "TClass.h"
33 #include "TVirtualPad.h" // for gPad
34 
35 #include "TBackCompFitter.h"
36 #include "TFitResultPtr.h"
37 #include "TFitResult.h"
38 
39 #include <stdlib.h>
40 #include <cmath>
41 #include <memory>
42 #include <limits>
43 
44 //#define DEBUG
45 
46 // utility functions used in TH1::Fit
47 
48 namespace HFit {
49 
50 
51  int GetDimension(const TH1 * h1) { return h1->GetDimension(); }
52  int GetDimension(const TGraph * ) { return 1; }
53  int GetDimension(const TMultiGraph * ) { return 1; }
54  int GetDimension(const TGraph2D * ) { return 2; }
55  int GetDimension(const THnBase * s1) { return s1->GetNdimensions(); }
56 
57  int CheckFitFunction(const TF1 * f1, int hdim);
58 
59 
60  void GetFunctionRange(const TF1 & f1, ROOT::Fit::DataRange & range);
61 
62  void FitOptionsMake(const char *option, Foption_t &fitOption);
63 
64  void CheckGraphFitOptions(Foption_t &fitOption);
65 
66 
67  void GetDrawingRange(TH1 * h1, ROOT::Fit::DataRange & range);
70  void GetDrawingRange(TGraph2D * gr, ROOT::Fit::DataRange & range);
72 
73 
74  template <class FitObject>
75  TFitResultPtr Fit(FitObject * h1, TF1 *f1 , Foption_t & option , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range);
76 
77  template <class FitObject>
78  void StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool, bool, const char *goption);
79 
80  template <class FitObject>
81  double ComputeChi2(const FitObject & h1, TF1 &f1, bool useRange, bool usePL );
82 
83 
84 
85 }
86 
87 int HFit::CheckFitFunction(const TF1 * f1, int dim) {
88  // Check validity of fitted function
89  if (!f1) {
90  Error("Fit", "function may not be null pointer");
91  return -1;
92  }
93  if (f1->IsZombie()) {
94  Error("Fit", "function is zombie");
95  return -2;
96  }
97 
98  int npar = f1->GetNpar();
99  if (npar <= 0) {
100  Error("Fit", "function %s has illegal number of parameters = %d", f1->GetName(), npar);
101  return -3;
102  }
103 
104  // Check that function has same dimension as histogram
105  if (f1->GetNdim() > dim) {
106  Error("Fit","function %s dimension, %d, is greater than fit object dimension, %d",
107  f1->GetName(), f1->GetNdim(), dim);
108  return -4;
109  }
110  if (f1->GetNdim() < dim-1) {
111  Error("Fit","function %s dimension, %d, is smaller than fit object dimension -1, %d",
112  f1->GetName(), f1->GetNdim(), dim);
113  return -5;
114  }
115 
116  return 0;
117 
118 }
119 
120 
122  // get the range form the function and fill and return the DataRange object
123  Double_t fxmin, fymin, fzmin, fxmax, fymax, fzmax;
124  f1.GetRange(fxmin, fymin, fzmin, fxmax, fymax, fzmax);
125  // support only one range - so add only if was not set before
126  if (range.Size(0) == 0) range.AddRange(0,fxmin,fxmax);
127  if (range.Size(1) == 0) range.AddRange(1,fymin,fymax);
128  if (range.Size(2) == 0) range.AddRange(2,fzmin,fzmax);
129  return;
130 }
131 
132 
133 template<class FitObject>
134 TFitResultPtr HFit::Fit(FitObject * h1, TF1 *f1 , Foption_t & fitOption , const ROOT::Math::MinimizerOptions & minOption, const char *goption, ROOT::Fit::DataRange & range)
135 {
136  // perform fit of histograms, or graphs using new fitting classes
137  // use same routines for fitting both graphs and histograms
138 
139 #ifdef DEBUG
140  printf("fit function %s\n",f1->GetName() );
141 #endif
142 
143  // replacement function using new fitter
144  int hdim = HFit::GetDimension(h1);
145  int iret = HFit::CheckFitFunction(f1, hdim);
146  if (iret != 0) return iret;
147 
148 
149 
150  // integral option is not supported in this case
151  if (f1->GetNdim() < hdim ) {
152  if (fitOption.Integral) Info("Fit","Ignore Integral option. Model function dimension is less than the data object dimension");
153  if (fitOption.Like) Info("Fit","Ignore Likelihood option. Model function dimension is less than the data object dimension");
154  fitOption.Integral = 0;
155  fitOption.Like = 0;
156  }
157 
158  Int_t special = f1->GetNumber();
159  Bool_t linear = f1->IsLinear();
160  Int_t npar = f1->GetNpar();
161  if (special==299+npar) linear = kTRUE; // for polynomial functions
162  // do not use linear fitter in these case
163  if (fitOption.Bound || fitOption.Like || fitOption.Errors || fitOption.Gradient || fitOption.More || fitOption.User|| fitOption.Integral || fitOption.Minuit)
164  linear = kFALSE;
165 
166  // create an empty TFitResult
167  std::shared_ptr<TFitResult> tfr(new TFitResult() );
168  // create the fitter from an empty fit result
169  std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(std::static_pointer_cast<ROOT::Fit::FitResult>(tfr) ) );
170  ROOT::Fit::FitConfig & fitConfig = fitter->Config();
171 
172  // create options
174  opt.fIntegral = fitOption.Integral;
175  opt.fUseRange = fitOption.Range;
176  opt.fExpErrors = fitOption.PChi2; // pearson chi2 with expected errors
177  if (fitOption.Like || fitOption.PChi2) opt.fUseEmpty = true; // use empty bins in log-likelihood fits
178  if (special==300) opt.fCoordErrors = false; // no need to use coordinate errors in a pol0 fit
179  if (fitOption.NoErrX) opt.fCoordErrors = false; // do not use coordinate errors when requested
180  if (fitOption.W1 ) opt.fErrors1 = true;
181  if (fitOption.W1 > 1) opt.fUseEmpty = true; // use empty bins with weight=1
182 
183  if (fitOption.BinVolume) {
184  opt.fBinVolume = true; // scale by bin volume
185  if (fitOption.BinVolume == 2) opt.fNormBinVolume = true; // scale by normalized bin volume
186  }
187 
188  if (opt.fUseRange) {
189 #ifdef DEBUG
190  printf("use range \n" );
191 #endif
192  HFit::GetFunctionRange(*f1,range);
193  }
194 #ifdef DEBUG
195  printf("range size %d\n", range.Size(0) );
196  if (range.Size(0)) {
197  double x1; double x2; range.GetRange(0,x1,x2);
198  printf(" range in x = [%f,%f] \n",x1,x2);
199  }
200 #endif
201 
202  // fill data
203  std::shared_ptr<ROOT::Fit::BinData> fitdata(new ROOT::Fit::BinData(opt,range) );
204  ROOT::Fit::FillData(*fitdata, h1, f1);
205  if (fitdata->Size() == 0 ) {
206  Warning("Fit","Fit data is empty ");
207  return -1;
208  }
209 
210 #ifdef DEBUG
211  printf("HFit:: data size is %d \n",fitdata->Size());
212  for (unsigned int i = 0; i < fitdata->Size(); ++i) {
213  if (fitdata->NDim() == 1) printf(" x[%d] = %f - value = %f \n", i,*(fitdata->Coords(i)),fitdata->Value(i) );
214  }
215 #endif
216 
217  // switch off linear fitting in case data has coordinate errors and the option is set
218  if (fitdata->GetErrorType() == ROOT::Fit::BinData::kCoordError && fitdata->Opt().fCoordErrors ) linear = false;
219  // linear fit cannot be done also in case of asymmetric errors
220  if (fitdata->GetErrorType() == ROOT::Fit::BinData::kAsymError && fitdata->Opt().fAsymErrors ) linear = false;
221 
222  // this functions use the TVirtualFitter
223  if (special != 0 && !fitOption.Bound && !linear) {
224  if (special == 100) ROOT::Fit::InitGaus (*fitdata,f1); // gaussian
225  else if (special == 110 || special == 112) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D gaussians ( xygaus or bigaus)
226  else if (special == 400) ROOT::Fit::InitGaus (*fitdata,f1); // landau (use the same)
227  else if (special == 410) ROOT::Fit::Init2DGaus(*fitdata,f1); // 2D landau (use the same)
228 
229  else if (special == 200) ROOT::Fit::InitExpo (*fitdata, f1); // exponential
230 
231  }
232 
233 
234  // set the fit function
235  // if option grad is specified use gradient
236  if ( (linear || fitOption.Gradient) )
237  fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*f1));
238 #ifdef R__HAS_VECCORE
239  else if(f1->IsVectorized())
241 #endif
242  else
243  fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*f1) ) );
244 
245  // error normalization in case of zero error in the data
246  if (fitdata->GetErrorType() == ROOT::Fit::BinData::kNoError) fitConfig.SetNormErrors(true);
247  // normalize errors also in case you are fitting a Ndim histo with a N-1 function
248  if (int(fitdata->NDim()) == hdim -1 ) fitConfig.SetNormErrors(true);
249 
250 
251  // here need to get some static extra information (like max iterations, error def, etc...)
252 
253 
254  // parameter settings and transfer the parameters values, names and limits from the functions
255  // is done automatically in the Fitter.cxx
256  for (int i = 0; i < npar; ++i) {
257  ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
258 
259  // check limits
260  double plow,pup;
261  f1->GetParLimits(i,plow,pup);
262  if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
263  parSettings.Fix();
264  }
265  else if (plow < pup ) {
266  if (!TMath::Finite(pup) && TMath::Finite(plow) )
267  parSettings.SetLowerLimit(plow);
268  else if (!TMath::Finite(plow) && TMath::Finite(pup) )
269  parSettings.SetUpperLimit(pup);
270  else
271  parSettings.SetLimits(plow,pup);
272  }
273 
274  // set the parameter step size (by default are set to 0.3 of value)
275  // if function provides meaningful error values
276  double err = f1->GetParError(i);
277  if ( err > 0)
278  parSettings.SetStepSize(err);
279  else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
280  double step = 0.1 * (pup - plow);
281  // check if value is not too close to limit otherwise trim value
282  if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
283  step = (pup - parSettings.Value() ) / 2;
284  else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
285  step = (parSettings.Value() - plow ) / 2;
286 
287  parSettings.SetStepSize(step);
288  }
289 
290 
291  }
292 
293  // needed for setting precision ?
294  // - Compute sum of squares of errors in the bin range
295  // should maybe use stat[1] ??
296  // Double_t ey, sumw2=0;
297 // for (i=hxfirst;i<=hxlast;i++) {
298 // ey = GetBinError(i);
299 // sumw2 += ey*ey;
300 // }
301 
302 
303  // set all default minimizer options (tolerance, max iterations, etc..)
304  fitConfig.SetMinimizerOptions(minOption);
305 
306  // specific print level options
307  if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
308  if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
309 
310  // specific minimizer options depending on minimizer
311  if (linear) {
312  if (fitOption.Robust ) {
313  // robust fitting
314  std::string type = "Robust";
315  // if an h is specified print out the value adding to the type
316  if (fitOption.hRobust > 0 && fitOption.hRobust < 1.)
317  type += " (h=" + ROOT::Math::Util::ToString(fitOption.hRobust) + ")";
318  fitConfig.SetMinimizer("Linear",type.c_str());
319  fitConfig.MinimizerOptions().SetTolerance(fitOption.hRobust); // use tolerance for passing robust parameter
320  }
321  else
322  fitConfig.SetMinimizer("Linear","");
323  }
324  else {
325  if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
326  }
327 
328 
329  // check if Error option (run Hesse and Minos) then
330  if (fitOption.Errors) {
331  // run Hesse and Minos
332  fitConfig.SetParabErrors(true);
333  fitConfig.SetMinosErrors(true);
334  }
335 
336 
337  // do fitting
338 
339 #ifdef DEBUG
340  if (fitOption.Like) printf("do likelihood fit...\n");
341  if (linear) printf("do linear fit...\n");
342  printf("do now fit...\n");
343 #endif
344 
345  bool fitok = false;
346 
347 
348  // check if can use option user
349  //typedef void (* MinuitFCN_t )(int &npar, double *gin, double &f, double *u, int flag);
350  TVirtualFitter::FCNFunc_t userFcn = 0;
351  if (fitOption.User && TVirtualFitter::GetFitter() ) {
352  userFcn = (TVirtualFitter::GetFitter())->GetFCN();
353  (TVirtualFitter::GetFitter())->SetUserFunc(f1);
354  }
355 
356 
357  if (fitOption.User && userFcn) // user provided fit objective function
358  fitok = fitter->FitFCN( userFcn );
359  else if (fitOption.Like) {// likelihood fit
360  // perform a weighted likelihood fit by applying weight correction to errors
361  bool weight = ((fitOption.Like & 2) == 2);
362  fitConfig.SetWeightCorrection(weight);
363  bool extended = ((fitOption.Like & 4 ) != 4 );
364  //if (!extended) Info("HFitImpl","Do a not -extended binned fit");
365  fitok = fitter->LikelihoodFit(*fitdata, extended, fitOption.ExecPolicy);
366  }
367  else{ // standard least square fit
368  fitok = fitter->Fit(*fitdata, fitOption.ExecPolicy);
369  }
370  if ( !fitok && !fitOption.Quiet )
371  Warning("Fit","Abnormal termination of minimization.");
372  iret |= !fitok;
373 
374 
375  const ROOT::Fit::FitResult & fitResult = fitter->Result();
376  // one could set directly the fit result in TF1
377  iret = fitResult.Status();
378  if (!fitResult.IsEmpty() ) {
379  // set in f1 the result of the fit
380  f1->SetChisquare(fitResult.Chi2() );
381  f1->SetNDF(fitResult.Ndf() );
382  f1->SetNumberFitPoints(fitdata->Size() );
383 
384  assert((Int_t)fitResult.Parameters().size() >= f1->GetNpar() );
385  f1->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
386  if ( int( fitResult.Errors().size()) >= f1->GetNpar() )
387  f1->SetParErrors( &(fitResult.Errors().front()) );
388 
389 
390  }
391 
392 // - Store fitted function in histogram functions list and draw
393  if (!fitOption.Nostore) {
394  HFit::GetDrawingRange(h1, range);
395  HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption);
396  }
397 
398  // print the result
399  // if using Fitter class must be done here
400  // use old style Minuit for TMinuit and if no corrections have been applied
401  if (!fitOption.Quiet) {
402  if (fitter->GetMinimizer() && fitConfig.MinimizerType() == "Minuit" &&
403  !fitConfig.NormalizeErrors() && fitOption.Like <= 1) {
404  fitter->GetMinimizer()->PrintResults(); // use old style Minuit
405  }
406  else {
407  // print using FitResult class
408  if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
409  fitResult.Print(std::cout);
410  }
411  }
412 
413 
414  // store result in the backward compatible VirtualFitter
415  // in case multi-thread is not enabled
416  if (!gGlobalMutex) {
417  TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
418  TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
419  bcfitter->SetFitOption(fitOption);
420  bcfitter->SetObjectFit(h1);
421  bcfitter->SetUserFunc(f1);
423  if (userFcn) {
424  bcfitter->SetFCN(userFcn);
425  // for interpreted FCN functions
426  if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
427  }
428 
429  // delete last fitter if it has been created here before
430  if (lastFitter) {
431  TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter);
432  if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) )
433  delete lastBCFitter;
434  }
435  //N.B= this might create a memory leak if user does not delete the fitter they create
436  TVirtualFitter::SetFitter( bcfitter );
437  }
438 
439  // use old-style for printing the results
440  // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
441  // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
442 
443  if (fitOption.StoreResult)
444  {
445  TString name = "TFitResult-";
446  name = name + h1->GetName() + "-" + f1->GetName();
447  TString title = "TFitResult-";
448  title += h1->GetTitle();
449  tfr->SetName(name);
450  tfr->SetTitle(title);
451  return TFitResultPtr(tfr);
452  }
453  else
454  return TFitResultPtr(iret);
455 }
456 
457 
459  // get range from histogram and update the DataRange class
460  // if a ranges already exist in that dimension use that one
461 
462  Int_t ndim = GetDimension(h1);
463 
464  double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
465  if (range.Size(0) == 0) {
466  TAxis & xaxis = *(h1->GetXaxis());
467  Int_t hxfirst = xaxis.GetFirst();
468  Int_t hxlast = xaxis.GetLast();
469  Double_t binwidx = xaxis.GetBinWidth(hxlast);
470  xmin = xaxis.GetBinLowEdge(hxfirst);
471  xmax = xaxis.GetBinLowEdge(hxlast) +binwidx;
472  range.AddRange(xmin,xmax);
473  }
474 
475  if (ndim > 1) {
476  if (range.Size(1) == 0) {
477  TAxis & yaxis = *(h1->GetYaxis());
478  Int_t hyfirst = yaxis.GetFirst();
479  Int_t hylast = yaxis.GetLast();
480  Double_t binwidy = yaxis.GetBinWidth(hylast);
481  ymin = yaxis.GetBinLowEdge(hyfirst);
482  ymax = yaxis.GetBinLowEdge(hylast) +binwidy;
483  range.AddRange(1,ymin,ymax);
484  }
485  }
486  if (ndim > 2) {
487  if (range.Size(2) == 0) {
488  TAxis & zaxis = *(h1->GetZaxis());
489  Int_t hzfirst = zaxis.GetFirst();
490  Int_t hzlast = zaxis.GetLast();
491  Double_t binwidz = zaxis.GetBinWidth(hzlast);
492  zmin = zaxis.GetBinLowEdge(hzfirst);
493  zmax = zaxis.GetBinLowEdge(hzlast) +binwidz;
494  range.AddRange(2,zmin,zmax);
495  }
496  }
497 #ifdef DEBUG
498  std::cout << "xmin,xmax" << xmin << " " << xmax << std::endl;
499 #endif
500 
501 }
502 
504  // get range for graph (used sub-set histogram)
505  // N.B. : this is different than in previous implementation of TGraph::Fit where range used was from xmin to xmax.
506  TH1 * h1 = gr->GetHistogram();
507  // an histogram is normally always returned for a TGraph
508  if (h1) HFit::GetDrawingRange(h1, range);
509 }
511  // get range for multi-graph (used sub-set histogram)
512  // N.B. : this is different than in previous implementation of TMultiGraph::Fit where range used was from data xmin to xmax.
513  TH1 * h1 = mg->GetHistogram();
514  if (h1) {
515  HFit::GetDrawingRange(h1, range);
516  }
517  else if (range.Size(0) == 0) {
518  // compute range from all the TGraph's belonging to the MultiGraph
519  double xmin = std::numeric_limits<double>::infinity();
520  double xmax = -std::numeric_limits<double>::infinity();
521  TIter next(mg->GetListOfGraphs() );
522  TGraph * g = 0;
523  while ( (g = (TGraph*) next() ) ) {
524  double x1 = 0, x2 = 0, y1 = 0, y2 = 0;
525  g->ComputeRange(x1,y1,x2,y2);
526  if (x1 < xmin) xmin = x1;
527  if (x2 > xmax) xmax = x2;
528  }
529  range.AddRange(xmin,xmax);
530  }
531 }
533  // get range for graph2D (used sub-set histogram)
534  // N.B. : this is different than in previous implementation of TGraph2D::Fit. There range used was always(0,0)
535  // cannot use TGraph2D::GetHistogram which makes an interpolation
536  //TH1 * h1 = gr->GetHistogram();
537  //if (h1) HFit::GetDrawingRange(h1, range);
538  // not very efficient (t.b.i.)
539  if (range.Size(0) == 0) {
540  double xmin = gr->GetXmin();
541  double xmax = gr->GetXmax();
542  range.AddRange(0,xmin,xmax);
543  }
544  if (range.Size(1) == 0) {
545  double ymin = gr->GetYmin();
546  double ymax = gr->GetYmax();
547  range.AddRange(1,ymin,ymax);
548  }
549 }
550 
552  // get range from histogram and update the DataRange class
553  // if a ranges already exist in that dimension use that one
554 
555  Int_t ndim = GetDimension(s1);
556 
557  for ( int i = 0; i < ndim; ++i ) {
558  if ( range.Size(i) == 0 ) {
559  TAxis *axis = s1->GetAxis(i);
560  range.AddRange(i, axis->GetXmin(), axis->GetXmax());
561  }
562  }
563 }
564 
565 template<class FitObject>
566 void HFit::StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool delOldFunction, bool drawFunction, const char *goption) {
567 // - Store fitted function in histogram functions list and draw
568 // should have separate functions for 1,2,3d ? t.b.d in case
569 
570 #ifdef DEBUG
571  std::cout <<"draw and store fit function " << f1->GetName() << std::endl;
572 #endif
573 
574 
575  Int_t ndim = GetDimension(h1);
576  double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
577  if (range.Size(0) ) range.GetRange(0,xmin,xmax);
578  if (range.Size(1) ) range.GetRange(1,ymin,ymax);
579  if (range.Size(2) ) range.GetRange(2,zmin,zmax);
580 
581 
582 #ifdef DEBUG
583  std::cout <<"draw and store fit function " << f1->GetName()
584  << " Range in x = [ " << xmin << " , " << xmax << " ]" << std::endl;
585 #endif
586 
587  TList * funcList = h1->GetListOfFunctions();
588  if (funcList == 0){
589  Error("StoreAndDrawFitFunction","Function list has not been created - cannot store the fitted function");
590  return;
591  }
592 
593  // delete the function in the list only if
594  // the function we are fitting is not in that list
595  // If this is the case we re-use that function object and
596  // we do not create a new one (if delOldFunction is true)
597  bool reuseOldFunction = false;
598  if (delOldFunction) {
599  TIter next(funcList, kIterBackward);
600  TObject *obj;
601  while ((obj = next())) {
602  if (obj->InheritsFrom(TF1::Class())) {
603  if (obj != f1) {
604  funcList->Remove(obj);
605  delete obj;
606  }
607  else {
608  reuseOldFunction = true;
609  }
610  }
611  }
612  }
613 
614  TF1 *fnew1 = 0;
615  TF2 *fnew2 = 0;
616  TF3 *fnew3 = 0;
617 
618  // copy TF1 using TClass to avoid slicing in case of derived classes
619  if (ndim < 2) {
620  if (!reuseOldFunction) {
621  fnew1 = (TF1*)f1->IsA()->New();
622  R__ASSERT(fnew1);
623  f1->Copy(*fnew1);
624  funcList->Add(fnew1);
625  }
626  else {
627  fnew1 = f1;
628  }
629  fnew1->SetParent( h1 );
630  fnew1->SetRange(xmin,xmax);
631  fnew1->Save(xmin,xmax,0,0,0,0);
632  if (!drawFunction) fnew1->SetBit(TF1::kNotDraw);
633  fnew1->AddToGlobalList(false);
634  } else if (ndim < 3) {
635  if (!reuseOldFunction) {
636  fnew2 = (TF2*)f1->IsA()->New();
637  R__ASSERT(fnew2);
638  f1->Copy(*fnew2);
639  funcList->Add(fnew2);
640  }
641  else {
642  fnew2 = dynamic_cast<TF2*>(f1);
643  R__ASSERT(fnew2);
644  }
645  fnew2->SetRange(xmin,ymin,xmax,ymax);
646  fnew2->SetParent( h1 );
647  fnew2->Save(xmin,xmax,ymin,ymax,0,0);
648  if (!drawFunction) fnew2->SetBit(TF1::kNotDraw);
649  fnew2->AddToGlobalList(false);
650  } else {
651  if (!reuseOldFunction) {
652  fnew3 = (TF3*)f1->IsA()->New();
653  R__ASSERT(fnew3);
654  f1->Copy(*fnew3);
655  funcList->Add(fnew3);
656  }
657  else {
658  fnew2 = dynamic_cast<TF3*>(f1);
659  R__ASSERT(fnew3);
660  }
661  fnew3->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
662  fnew3->SetParent( h1 );
663  fnew3->Save(xmin,xmax,ymin,ymax,zmin,zmax);
664  if (!drawFunction) fnew3->SetBit(TF1::kNotDraw);
665  fnew3->AddToGlobalList(false);
666  }
667  if (h1->TestBit(kCanDelete)) return;
668  // draw only in case of histograms
669  if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) {
670  // no need to re-draw the histogram if the histogram is already in the pad
671  // in that case the function will be just drawn (if option N is not set)
672  if (!gPad || (gPad && gPad->GetListOfPrimitives()->FindObject(h1) == NULL ) )
673  h1->Draw(goption);
674  }
675  if (gPad) gPad->Modified(); // this is not in TH1 code (needed ??)
676 
677  return;
678 }
679 
680 
681 void ROOT::Fit::FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption) {
682  // - Decode list of options into fitOption (used by both TGraph and TH1)
683  // works for both histograms and graph depending on the enum FitObjectType defined in HFit
686  }
687 
688  if (option == 0) return;
689  if (!option[0]) return;
690 
691  TString opt = option;
692  opt.ToUpper();
693 
694  // parse firt the specific options
695  if (type == kHistogram) {
696 
697  if (opt.Contains("WIDTH")) {
698  fitOption.BinVolume = 1; // scale content by the bin width
699  if (opt.Contains("NORMWIDTH")) {
700  // for variable bins: scale content by the bin width normalized by a reference value (typically the minimum bin)
701  // this option is for variable bin widths
702  fitOption.BinVolume = 2;
703  opt.ReplaceAll("NORMWIDTH","");
704  }
705  else
706  opt.ReplaceAll("WIDTH","");
707  }
708 
709  // if (opt.Contains("MULTIPROC")) {
710  // fitOption.ExecPolicy = ROOT::Fit::kMultiprocess;
711  // opt.ReplaceAll("MULTIPROC","");
712  // }
713 
714  if (opt.Contains("SERIAL")) {
716  opt.ReplaceAll("SERIAL","");
717  }
718 
719  if (opt.Contains("MULTITHREAD")) {
721  opt.ReplaceAll("MULTITHREAD","");
722  }
723 
724  if (opt.Contains("I")) fitOption.Integral= 1; // integral of function in the bin (no sense for graph)
725  if (opt.Contains("WW")) fitOption.W1 = 2; //all bins have weight=1, even empty bins
726  }
727 
728  // specific Graph options (need to be parsed before)
729  else if (type == kGraph) {
730  opt.ReplaceAll("ROB", "H");
731  opt.ReplaceAll("EX0", "T");
732 
733  //for robust fitting, see if # of good points is defined
734  // decode parameters for robust fitting
735  Double_t h=0;
736  if (opt.Contains("H=0.")) {
737  int start = opt.Index("H=0.");
738  int numpos = start + strlen("H=0.");
739  int numlen = 0;
740  int len = opt.Length();
741  while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
742  TString num = opt(numpos,numlen);
743  opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
744  h = atof(num.Data());
745  h*=TMath::Power(10, -numlen);
746  }
747 
748  if (opt.Contains("H")) { fitOption.Robust = 1; fitOption.hRobust = h; }
749  if (opt.Contains("T")) fitOption.NoErrX = 1; // no error in X
750 
751  }
752 
753  if (opt.Contains("U")) fitOption.User = 1;
754  if (opt.Contains("Q")) fitOption.Quiet = 1;
755  if (opt.Contains("V")) {fitOption.Verbose = 1; fitOption.Quiet = 0;}
756  if (opt.Contains("L")) fitOption.Like = 1;
757  if (opt.Contains("X")) fitOption.Chi2 = 1;
758  if (opt.Contains("P")) fitOption.PChi2 = 1;
759 
760 
761  // likelihood fit options
762  if (fitOption.Like == 1) {
763  //if (opt.Contains("LL")) fitOption.Like = 2;
764  if (opt.Contains("W")){ fitOption.Like = 2; fitOption.W1=0;}// (weighted likelihood)
765  if (opt.Contains("MULTI")) {
766  if (fitOption.Like == 2) fitOption.Like = 6; // weighted multinomial
767  else fitOption.Like = 4; // multinomial likelihood fit instead of Poisson
768  opt.ReplaceAll("MULTI","");
769  }
770  // in case of histogram give precedence for likelihood options
771  if (type == kHistogram) {
772  if (fitOption.Chi2 == 1 || fitOption.PChi2 == 1)
773  Warning("Fit","Cannot use P or X option in combination of L. Ignore the chi2 option and perform a likelihood fit");
774  }
775 
776  } else {
777  if (opt.Contains("W")) fitOption.W1 = 1; // all non-empty bins have weight =1 (for chi2 fit)
778  }
779 
780 
781  if (opt.Contains("E")) fitOption.Errors = 1;
782  if (opt.Contains("R")) fitOption.Range = 1;
783  if (opt.Contains("G")) fitOption.Gradient= 1;
784  if (opt.Contains("M")) fitOption.More = 1;
785  if (opt.Contains("N")) fitOption.Nostore = 1;
786  if (opt.Contains("0")) fitOption.Nograph = 1;
787  if (opt.Contains("+")) fitOption.Plus = 1;
788  if (opt.Contains("B")) fitOption.Bound = 1;
789  if (opt.Contains("C")) fitOption.Nochisq = 1;
790  if (opt.Contains("F")) fitOption.Minuit = 1;
791  if (opt.Contains("S")) fitOption.StoreResult = 1;
792 
793 }
794 
796  if (foption.Like) {
797  Info("CheckGraphFitOptions","L (Log Likelihood fit) is an invalid option when fitting a graph. It is ignored");
798  foption.Like = 0;
799  }
800  if (foption.Integral) {
801  Info("CheckGraphFitOptions","I (use function integral) is an invalid option when fitting a graph. It is ignored");
802  foption.Integral = 0;
803  }
804  return;
805 }
806 
807 // implementation of unbin fit function (defined in HFitInterface)
808 
810  // do unbin fit, ownership of fitdata is passed later to the TBackFitter class
811 
812  // create a shared pointer to the fit data to managed it
813  std::shared_ptr<ROOT::Fit::UnBinData> fitdata(data);
814 
815 #ifdef DEBUG
816  printf("tree data size is %d \n",fitdata->Size());
817  for (unsigned int i = 0; i < fitdata->Size(); ++i) {
818  if (fitdata->NDim() == 1) printf(" x[%d] = %f \n", i,*(fitdata->Coords(i) ) );
819  }
820 #endif
821  if (fitdata->Size() == 0 ) {
822  Warning("Fit","Fit data is empty ");
823  return -1;
824  }
825 
826  // create an empty TFitResult
827  std::shared_ptr<TFitResult> tfr(new TFitResult() );
828  // create the fitter
829  std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(tfr) );
830  ROOT::Fit::FitConfig & fitConfig = fitter->Config();
831 
832  // dimension is given by data because TF1 pointer can have wrong one
833  unsigned int dim = fitdata->NDim();
834 
835  // set the fit function
836  // if option grad is specified use gradient
837  // need to create a wrapper for an automatic normalized TF1 ???
838  if ( fitOption.Gradient ) {
839  assert ( (int) dim == fitfunc->GetNdim() );
840  fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*fitfunc) );
841  }
842  else
843  fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*fitfunc, dim) ) );
844 
845  // parameter setting is done automaticaly in the Fitter class
846  // need only to set limits
847  int npar = fitfunc->GetNpar();
848  for (int i = 0; i < npar; ++i) {
849  ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
850  double plow,pup;
851  fitfunc->GetParLimits(i,plow,pup);
852  // this is a limitation of TF1 interface - cannot fix a parameter to zero value
853  if (plow*pup != 0 && plow >= pup) {
854  parSettings.Fix();
855  }
856  else if (plow < pup ) {
857  if (!TMath::Finite(pup) && TMath::Finite(plow) )
858  parSettings.SetLowerLimit(plow);
859  else if (!TMath::Finite(plow) && TMath::Finite(pup) )
860  parSettings.SetUpperLimit(pup);
861  else
862  parSettings.SetLimits(plow,pup);
863  }
864 
865  // set the parameter step size (by default are set to 0.3 of value)
866  // if function provides meaningful error values
867  double err = fitfunc->GetParError(i);
868  if ( err > 0)
869  parSettings.SetStepSize(err);
870  else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
871  double step = 0.1 * (pup - plow);
872  // check if value is not too close to limit otherwise trim value
873  if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
874  step = (pup - parSettings.Value() ) / 2;
875  else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
876  step = (parSettings.Value() - plow ) / 2;
877 
878  parSettings.SetStepSize(step);
879  }
880 
881  }
882 
883  fitConfig.SetMinimizerOptions(minOption);
884 
885  if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
886  if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
887 
888  // more
889  if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
890 
891  // chech if Minos or more options
892  if (fitOption.Errors) {
893  // run Hesse and Minos
894  fitConfig.SetParabErrors(true);
895  fitConfig.SetMinosErrors(true);
896  }
897  // use weight correction
898  if ( (fitOption.Like & 2) == 2)
899  fitConfig.SetWeightCorrection(true);
900 
901  bool extended = (fitOption.Like & 1) == 1;
902 
903  bool fitok = false;
904  fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
905  if ( !fitok && !fitOption.Quiet )
906  Warning("UnBinFit","Abnormal termination of minimization.");
907 
908  const ROOT::Fit::FitResult & fitResult = fitter->Result();
909  // one could set directly the fit result in TF1
910  int iret = fitResult.Status();
911  if (!fitResult.IsEmpty() ) {
912  // set in fitfunc the result of the fit
913  fitfunc->SetNDF(fitResult.Ndf() );
914  fitfunc->SetNumberFitPoints(fitdata->Size() );
915 
916  assert( (Int_t)fitResult.Parameters().size() >= fitfunc->GetNpar() );
917  fitfunc->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
918  if ( int( fitResult.Errors().size()) >= fitfunc->GetNpar() )
919  fitfunc->SetParErrors( &(fitResult.Errors().front()) );
920 
921  }
922 
923  // store result in the backward compatible VirtualFitter
924  // in case not running in a multi-thread enabled mode
925  if (gGlobalMutex) {
926  TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
927  // pass ownership of Fitter and Fitdata to TBackCompFitter (fitter pointer cannot be used afterwards)
928  TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
929  // cannot use anymore now fitdata (given away ownership)
930  fitdata = 0;
931  bcfitter->SetFitOption(fitOption);
932  //bcfitter->SetObjectFit(fTree);
933  bcfitter->SetUserFunc(fitfunc);
934 
935  if (lastFitter) delete lastFitter;
936  TVirtualFitter::SetFitter( bcfitter );
937 
938  // use old-style for printing the results
939  // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
940  // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
941 
942  }
943  // print results
944  if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
945  else if (!fitOption.Quiet) fitResult.Print(std::cout);
946 
947  if (fitOption.StoreResult)
948  {
949  TString name = "TFitResult-";
950  name = name + "UnBinData-" + fitfunc->GetName();
951  TString title = "TFitResult-";
952  title += name;
953  tfr->SetName(name);
954  tfr->SetTitle(title);
955  return TFitResultPtr(tfr);
956  }
957  else
958  return TFitResultPtr(iret);
959 }
960 
961 
962 // implementations of ROOT::Fit::FitObject functions (defined in HFitInterface) in terms of the template HFit::Fit
963 
965 moption, const char *goption, ROOT::Fit::DataRange & range) {
966  // check fit options
967  // check if have weights in case of weighted likelihood
968  if ( ((foption.Like & 2) == 2) && h1->GetSumw2N() == 0) {
969  Warning("HFit::FitObject","A weighted likelihood fit is requested but histogram is not weighted - do a standard Likelihood fit");
970  foption.Like = 1;
971  }
972  // histogram fitting
973  return HFit::Fit(h1,f1,foption,moption,goption,range);
974 }
975 
976 TFitResultPtr ROOT::Fit::FitObject(TGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
977  // exclude options not valid for graphs
979  // TGraph fitting
980  return HFit::Fit(gr,f1,foption,moption,goption,range);
981 }
982 
983 TFitResultPtr ROOT::Fit::FitObject(TMultiGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
984  // exclude options not valid for graphs
986  // TMultiGraph fitting
987  return HFit::Fit(gr,f1,foption,moption,goption,range);
988 }
989 
990 TFitResultPtr ROOT::Fit::FitObject(TGraph2D * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
991  // exclude options not valid for graphs
993  // TGraph2D fitting
994  return HFit::Fit(gr,f1,foption,moption,goption,range);
995 }
996 
997 TFitResultPtr ROOT::Fit::FitObject(THnBase * s1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
998  // sparse histogram fitting
999  return HFit::Fit(s1,f1,foption,moption,goption,range);
1000 }
1001 
1002 
1003 
1004 // Int_t TGraph2D::DoFit(TF2 *f2 ,Option_t *option ,Option_t *goption) {
1005 // // internal graph2D fitting methods
1006 // Foption_t fitOption;
1007 // ROOT::Fit::FitOptionsMake(option,fitOption);
1008 
1009 // // create range and minimizer options with default values
1010 // ROOT::Fit::DataRange range(2);
1011 // ROOT::Math::MinimizerOptions minOption;
1012 // return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
1013 // }
1014 
1015 
1016 // function to compute the simple chi2 for graphs and histograms
1017 
1018 double ROOT::Fit::Chisquare(const TH1 & h1, TF1 & f1, bool useRange, bool usePL) {
1019  return HFit::ComputeChi2(h1,f1,useRange, usePL);
1020 }
1021 
1022 double ROOT::Fit::Chisquare(const TGraph & g, TF1 & f1, bool useRange) {
1023  return HFit::ComputeChi2(g,f1, useRange, false);
1024 }
1025 
1026 template<class FitObject>
1027 double HFit::ComputeChi2(const FitObject & obj, TF1 & f1, bool useRange, bool usePL ) {
1028 
1029  // implement using the fitting classes
1031  if (usePL) opt.fUseEmpty=true;
1032  ROOT::Fit::DataRange range;
1033  // get range of function
1034  if (useRange) HFit::GetFunctionRange(f1,range);
1035  // fill the data set
1036  ROOT::Fit::BinData data(opt,range);
1037  ROOT::Fit::FillData(data, &obj, &f1);
1038  if (data.Size() == 0 ) {
1039  Warning("Chisquare","data set is empty - return -1");
1040  return -1;
1041  }
1043  if (usePL) {
1044  // use the poisson log-lokelihood (Baker-Cousins chi2)
1045  ROOT::Fit::PoissonLLFunction nll(data, wf1);
1046  return 2.* nll( f1.GetParameters() ) ;
1047  }
1048  ROOT::Fit::Chi2Function chi2(data, wf1);
1049  return chi2(f1.GetParameters() );
1050 
1051 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
int Errors
Definition: Foption.h:37
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:628
float xmin
Definition: THbookFile.cxx:93
void SetPrintLevel(int level)
set print level
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
bool IsVectorized()
Definition: TF1.h:424
double Chisquare(const TH1 &h1, TF1 &f1, bool useRange, bool usePL=false)
compute the chi2 value for an histogram given a function (see TH1::Chisquare for the documentation) ...
Definition: HFitImpl.cxx:1018
void SetTolerance(double tol)
set the tolerance
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition: TF1.cxx:904
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: TF3.cxx:534
virtual void SetMethodCall(TMethodCall *m)
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
float ymin
Definition: THbookFile.cxx:93
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition: FitResult.h:170
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:638
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
virtual void SetUserFunc(TObject *userfunc)
bool NormalizeErrors() const
flag to check if resulting errors are be normalized according to chi2/ndf
Definition: FitConfig.h:199
int Verbose
Definition: Foption.h:30
TH1 * h
Definition: legend2.C:5
Double_t GetXmax() const
Returns the X maximum.
Definition: TGraph2D.cxx:1198
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:35
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF1.cxx:3406
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1112
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:3021
#define R__ASSERT(e)
Definition: TError.h:96
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:585
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
double Value() const
copy constructor and assignment operators (leave them to the compiler)
unsigned int Ndf() const
Number of degree of freedom.
Definition: FitResult.h:164
Basic string class.
Definition: TString.h:125
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:608
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
Override setFCN to use the Adapter to Minuit2 FCN interface To set the address of the minimization fu...
Class describing the unbinned data sets (just x coordinates values) of any dimensions.
Definition: UnBinData.h:42
Backward compatible implementation of TVirtualFitter.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1815
Double_t GetYmax() const
Returns the Y maximum.
Definition: TGraph2D.cxx:1220
ROOT::Math::MinimizerOptions & MinimizerOptions()
access to the minimizer control parameter (non const method)
Definition: FitConfig.h:166
int Nograph
Definition: Foption.h:42
int Minuit
Definition: Foption.h:46
const ParameterSettings & ParSettings(unsigned int i) const
get the parameter settings for the i-th parameter (const method)
Definition: FitConfig.h:75
int Nostore
Definition: Foption.h:41
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:627
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
static constexpr double mg
TMethodCall * GetMethodCall() const
TFitResultPtr UnBinFit(ROOT::Fit::UnBinData *data, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption)
fit an unbin data set (from tree or from histogram buffer) using a TF1 pointer and fit options...
Definition: HFitImpl.cxx:809
double hRobust
Definition: Foption.h:51
virtual Int_t GetDimension() const
Definition: TH1.h:277
virtual void SetParent(TObject *p=0)
Definition: TF1.h:649
Double_t GetXmin() const
Definition: TAxis.h:133
static const double x2[5]
void GetRange(unsigned int irange, unsigned int icoord, double &xmin, double &xmax) const
get the i-th range for given coordinate.
Definition: DataRange.h:103
void Class()
Definition: Class.C:29
int GetDimension(const TH1 *h1)
Definition: HFitImpl.cxx:51
Extends the ROOT::Fit::Result class with a TNamed inheritance providing easy possibility for I/O...
Definition: TFitResult.h:30
class evaluating the log likelihood for binned Poisson likelihood fits it is template to distinguish ...
int Chi2
Definition: Foption.h:32
Double_t GetXmin() const
Returns the X minimum.
Definition: TGraph2D.cxx:1209
Class to Wrap a ROOT Function class (like TF1) in a IParamMultiFunction interface of multi-dimensions...
int Plus
Definition: Foption.h:43
void GetDrawingRange(TH1 *h1, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:458
void Info(const char *location, const char *msgfmt,...)
void Fix()
fix the parameter
Int_t Finite(Double_t x)
Definition: TMath.h:654
virtual Int_t GetNdim() const
Definition: TF1.h:469
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
TH1F * h1
Definition: legend1.C:5
ROOT::Fit::ExecutionPolicy ExecPolicy
Definition: Foption.h:52
int BinVolume
Definition: Foption.h:50
unsigned int Size() const
return number of fit points
Definition: FitData.h:303
void Error(const char *location, const char *msgfmt,...)
R__ALWAYS_INLINE Bool_t IsZombie() const
Definition: TObject.h:134
A doubly linked list.
Definition: TList.h:44
void SetLowerLimit(double low)
set a single lower limit
Definition: Buttons.h:34
static void SetFitter(TVirtualFitter *fitter, Int_t maxpar=25)
Static function to set an alternative fitter.
int CheckFitFunction(const TF1 *f1, int hdim)
Definition: HFitImpl.cxx:87
R__EXTERN TVirtualMutex * gGlobalMutex
Definition: TVirtualMutex.h:29
void SetMinosErrors(bool on=true)
set Minos erros computation to be performed after fitting
Definition: FitConfig.h:226
Chi2FCN class for binnned fits using the least square methods.
Definition: Chi2FCN.h:49
virtual void SetObjectFit(TObject *obj)
void InitGaus(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for gaussian function given the fit data Set the sigma limits for zero top ...
virtual void SetChisquare(Double_t chi2)
Definition: TF1.h:596
float ymax
Definition: THbookFile.cxx:93
int User
Definition: Foption.h:35
void SetMinimizer(const char *type, const char *algo=0)
set minimizer type
Definition: FitConfig.h:180
void SetStepSize(double err)
set the step size
void FillData(BinData &dv, const TH1 *hist, TF1 *func=0)
fill the data vector from a TH1.
TAxis * GetAxis(Int_t dim) const
Definition: THnBase.h:125
DataOptions : simple structure holding the options on how the data are filled.
Definition: DataOptions.h:28
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
Class to manage histogram axis.
Definition: TAxis.h:30
int PChi2
Definition: Foption.h:33
A 3-Dim function with parameters.
Definition: TF3.h:28
int Gradient
Definition: Foption.h:40
Fitter class, entry point for performing all type of fits.
Definition: Fitter.h:77
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:818
virtual void SetFitOption(Foption_t option)
int StoreResult
Definition: Foption.h:49
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:31
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
unsigned int Size(unsigned int icoord=0) const
return range size for coordinate icoord (starts from zero) Size == 0 indicates no range is present [-...
Definition: DataRange.h:70
void FitOptionsMake(const char *option, Foption_t &fitOption)
Ssiz_t Length() const
Definition: TString.h:386
virtual Bool_t IsLinear() const
Definition: TF1.h:586
static TVirtualFitter * GetFitter()
static: return the current Fitter
int More
Definition: Foption.h:38
int Like
Definition: Foption.h:34
TAxis * GetYaxis()
Definition: TH1.h:316
int Integral
Definition: Foption.h:44
Class describing the binned data sets : vectors of x coordinates, y values and optionally error on y ...
Definition: BinData.h:53
float xmax
Definition: THbookFile.cxx:93
int Status() const
minimizer status code
Definition: FitResult.h:138
A 2-Dim function with parameters.
Definition: TF2.h:29
void InitExpo(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for an exponential function given the fit data Set the constant and slope a...
void Warning(const char *location, const char *msgfmt,...)
void(* FCNFunc_t)(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag)
TGraphErrors * gr
Definition: legend1.C:25
int W1
Definition: Foption.h:36
void SetMinimizerOptions(const ROOT::Math::MinimizerOptions &minopt)
set all the minimizer options using class MinimizerOptions
Definition: FitConfig.cxx:236
const Bool_t kFALSE
Definition: RtypesCore.h:88
TString & Remove(Ssiz_t pos)
Definition: TString.h:619
void AddRange(unsigned int icoord, double xmin, double xmax)
add a range [xmin,xmax] for the new coordinate icoord Adding a range does not delete existing one...
Definition: DataRange.cxx:94
IParamFunction interface (abstract class) describing multi-dimensional parameteric functions It is a ...
int Robust
Definition: Foption.h:48
virtual Int_t GetSumw2N() const
Definition: TH1.h:309
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF2.h:150
class containg the result of the fit and all the related information (fitted parameter values...
Definition: FitResult.h:48
void PrintCovMatrix(std::ostream &os) const
print error matrix and correlations
Definition: FitResult.cxx:501
int Range
Definition: Foption.h:39
void GetFunctionRange(const TF1 &f1, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:121
static const double x1[5]
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition: DataRange.h:34
TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases...
Definition: TGraph.cxx:1469
void Init2DGaus(const ROOT::Fit::BinData &data, TF1 *f1)
compute initial parameter for 2D gaussian function given the fit data Set the sigma limits for zero t...
double Double_t
Definition: RtypesCore.h:55
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:134
int type
Definition: TGX11.cxx:120
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Definition: FitResult.h:175
void SetWeightCorrection(bool on=true)
apply the weight correction for error matric computation
Definition: FitConfig.h:229
const std::string & MinimizerType() const
return type of minimizer package
Definition: FitConfig.h:188
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:570
The TH1 histogram class.
Definition: TH1.h:56
static constexpr double s
int Bound
Definition: Foption.h:31
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition: HFitImpl.cxx:681
int Nochisq
Definition: Foption.h:45
Int_t GetNdimensions() const
Definition: THnBase.h:135
Abstract Base Class for Fitting.
TAxis * GetZaxis()
Definition: TH1.h:317
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
Definition: FitResult.h:118
Mother of all ROOT objects.
Definition: TObject.h:37
int NoErrX
Definition: Foption.h:47
virtual Int_t GetNpar() const
Definition: TF1.h:465
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1825
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition: TAxis.cxx:526
TList * GetListOfGraphs() const
Definition: TMultiGraph.h:69
void StoreAndDrawFitFunction(FitObject *h1, TF1 *f1, const ROOT::Fit::DataRange &range, bool, bool, const char *goption)
Definition: HFitImpl.cxx:566
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
Definition: FitResult.cxx:446
std::string ToString(const T &val)
Utility function for conversion to strings.
Definition: Util.h:42
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF3.h:143
virtual void Add(TObject *obj)
Definition: TList.h:87
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:584
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:3371
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition: TF1.cxx:2163
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:742
1-Dim function class
Definition: TF1.h:211
void SetNormErrors(bool on=true)
set the option to normalize the error on the result according to chi2/ndf
Definition: FitConfig.h:220
void SetUpperLimit(double up)
set a single upper limit
Double_t GetYmin() const
Returns the Y minimum.
Definition: TGraph2D.cxx:1231
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
double ComputeChi2(const FitObject &h1, TF1 &f1, bool useRange, bool usePL)
Definition: HFitImpl.cxx:1027
TF1 * f1
Definition: legend1.C:11
#define gPad
Definition: TVirtualPad.h:285
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:3306
virtual Double_t * GetParameters() const
Definition: TF1.h:504
void SetLimits(double low, double up)
set a double side limit, if low == up the parameter is fixed if low > up the limits are removed The c...
Multidimensional histogram base.
Definition: THnBase.h:43
TFitResultPtr FitObject(TH1 *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
fitting function for a TH1 (called from TH1::Fit)
Definition: HFitImpl.cxx:964
const Bool_t kIterBackward
Definition: TCollection.h:41
void SetParabErrors(bool on=true)
set parabolic erros
Definition: FitConfig.h:223
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:40
int Quiet
Definition: Foption.h:29
const Bool_t kTRUE
Definition: RtypesCore.h:87
Double_t GetXmax() const
Definition: TAxis.h:134
virtual Int_t GetNumber() const
Definition: TF1.h:482
double Chi2() const
Chi2 fit value in case of likelihood must be computed ?
Definition: FitResult.h:161
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
void CheckGraphFitOptions(Foption_t &fitOption)
Definition: HFitImpl.cxx:795
Class describing the configuration of the fit, options and parameter settings using the ROOT::Fit::Pa...
Definition: FitConfig.h:46
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: TF2.cxx:752
const char * Data() const
Definition: TString.h:345
static constexpr double g