Logo ROOT   master
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 #include "TROOT.h"
32 
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);
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  // error normalization also in case of W1 option (weights = 1)
248  if (fitdata->Opt().fErrors1) fitConfig.SetNormErrors(true);
249  // normalize errors also in case you are fitting a Ndim histo with a N-1 function
250  if (int(fitdata->NDim()) == hdim -1 ) fitConfig.SetNormErrors(true);
251 
252 
253  // here need to get some static extra information (like max iterations, error def, etc...)
254 
255 
256  // parameter settings and transfer the parameters values, names and limits from the functions
257  // is done automatically in the Fitter.cxx
258  for (int i = 0; i < npar; ++i) {
259  ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
260 
261  // check limits
262  double plow,pup;
263  f1->GetParLimits(i,plow,pup);
264  if (plow*pup != 0 && plow >= pup) { // this is a limitation - cannot fix a parameter to zero value
265  parSettings.Fix();
266  }
267  else if (plow < pup ) {
268  if (!TMath::Finite(pup) && TMath::Finite(plow) )
269  parSettings.SetLowerLimit(plow);
270  else if (!TMath::Finite(plow) && TMath::Finite(pup) )
271  parSettings.SetUpperLimit(pup);
272  else
273  parSettings.SetLimits(plow,pup);
274  }
275 
276  // set the parameter step size (by default are set to 0.3 of value)
277  // if function provides meaningful error values
278  double err = f1->GetParError(i);
279  if ( err > 0)
280  parSettings.SetStepSize(err);
281  else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
282  double step = 0.1 * (pup - plow);
283  // check if value is not too close to limit otherwise trim value
284  if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
285  step = (pup - parSettings.Value() ) / 2;
286  else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
287  step = (parSettings.Value() - plow ) / 2;
288 
289  parSettings.SetStepSize(step);
290  }
291 
292 
293  }
294 
295  // needed for setting precision ?
296  // - Compute sum of squares of errors in the bin range
297  // should maybe use stat[1] ??
298  // Double_t ey, sumw2=0;
299 // for (i=hxfirst;i<=hxlast;i++) {
300 // ey = GetBinError(i);
301 // sumw2 += ey*ey;
302 // }
303 
304 
305  // set all default minimizer options (tolerance, max iterations, etc..)
306  fitConfig.SetMinimizerOptions(minOption);
307 
308  // specific print level options
309  if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
310  if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
311 
312  // specific minimizer options depending on minimizer
313  if (linear) {
314  if (fitOption.Robust ) {
315  // robust fitting
316  std::string type = "Robust";
317  // if an h is specified print out the value adding to the type
318  if (fitOption.hRobust > 0 && fitOption.hRobust < 1.)
319  type += " (h=" + ROOT::Math::Util::ToString(fitOption.hRobust) + ")";
320  fitConfig.SetMinimizer("Linear",type.c_str());
321  fitConfig.MinimizerOptions().SetTolerance(fitOption.hRobust); // use tolerance for passing robust parameter
322  }
323  else
324  fitConfig.SetMinimizer("Linear","");
325  }
326  else {
327  if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
328  }
329 
330 
331  // check if Error option (run Hesse and Minos) then
332  if (fitOption.Errors) {
333  // run Hesse and Minos
334  fitConfig.SetParabErrors(true);
335  fitConfig.SetMinosErrors(true);
336  }
337 
338 
339  // do fitting
340 
341 #ifdef DEBUG
342  if (fitOption.Like) printf("do likelihood fit...\n");
343  if (linear) printf("do linear fit...\n");
344  printf("do now fit...\n");
345 #endif
346 
347  bool fitok = false;
348 
349 
350  // check if can use option user
351  //typedef void (* MinuitFCN_t )(int &npar, double *gin, double &f, double *u, int flag);
352  TVirtualFitter::FCNFunc_t userFcn = 0;
353  if (fitOption.User && TVirtualFitter::GetFitter() ) {
354  userFcn = (TVirtualFitter::GetFitter())->GetFCN();
355  (TVirtualFitter::GetFitter())->SetUserFunc(f1);
356  }
357 
358 
359  if (fitOption.User && userFcn) // user provided fit objective function
360  fitok = fitter->FitFCN( userFcn );
361  else if (fitOption.Like) {// likelihood fit
362  // perform a weighted likelihood fit by applying weight correction to errors
363  bool weight = ((fitOption.Like & 2) == 2);
364  fitConfig.SetWeightCorrection(weight);
365  bool extended = ((fitOption.Like & 4 ) != 4 );
366  //if (!extended) Info("HFitImpl","Do a not -extended binned fit");
367 
368  // pass fitdata as a shared pointer so ownership is shared with Fitter and FitResult class
369  fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
370  }
371  else{ // standard least square fit
372  fitok = fitter->Fit(fitdata, fitOption.ExecPolicy);
373  }
374  if ( !fitok && !fitOption.Quiet )
375  Warning("Fit","Abnormal termination of minimization.");
376  iret |= !fitok;
377 
378 
379  const ROOT::Fit::FitResult & fitResult = fitter->Result();
380  // one could set directly the fit result in TF1
381  iret = fitResult.Status();
382  if (!fitResult.IsEmpty() ) {
383  // set in f1 the result of the fit
384  f1->SetChisquare(fitResult.Chi2() );
385  f1->SetNDF(fitResult.Ndf() );
386  f1->SetNumberFitPoints(fitdata->Size() );
387 
388  assert((Int_t)fitResult.Parameters().size() >= f1->GetNpar() );
389  f1->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
390  if ( int( fitResult.Errors().size()) >= f1->GetNpar() )
391  f1->SetParErrors( &(fitResult.Errors().front()) );
392 
393 
394  }
395 
396 // - Store fitted function in histogram functions list and draw
397  if (!fitOption.Nostore) {
398  HFit::GetDrawingRange(h1, range);
399  HFit::StoreAndDrawFitFunction(h1, f1, range, !fitOption.Plus, !fitOption.Nograph, goption);
400  }
401 
402  // print the result
403  // if using Fitter class must be done here
404  // use old style Minuit for TMinuit and if no corrections have been applied
405  if (!fitOption.Quiet) {
406  if (fitter->GetMinimizer() && fitConfig.MinimizerType() == "Minuit" &&
407  !fitConfig.NormalizeErrors() && fitOption.Like <= 1) {
408  fitter->GetMinimizer()->PrintResults(); // use old style Minuit
409  }
410  else {
411  // print using FitResult class
412  if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
413  fitResult.Print(std::cout);
414  }
415  }
416 
417 
418  // store result in the backward compatible VirtualFitter
419  // in case multi-thread is not enabled
420  if (!gGlobalMutex) {
421  TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
422  TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
423  bcfitter->SetFitOption(fitOption);
424  bcfitter->SetObjectFit(h1);
425  bcfitter->SetUserFunc(f1);
427  if (userFcn) {
428  bcfitter->SetFCN(userFcn);
429  // for interpreted FCN functions
430  if (lastFitter->GetMethodCall() ) bcfitter->SetMethodCall(lastFitter->GetMethodCall() );
431  }
432 
433  // delete last fitter if it has been created here before
434  if (lastFitter) {
435  TBackCompFitter * lastBCFitter = dynamic_cast<TBackCompFitter *> (lastFitter);
436  if (lastBCFitter && lastBCFitter->TestBit(TBackCompFitter::kCanDeleteLast) )
437  delete lastBCFitter;
438  }
439  //N.B= this might create a memory leak if user does not delete the fitter they create
440  TVirtualFitter::SetFitter( bcfitter );
441  }
442 
443  // use old-style for printing the results
444  // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
445  // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
446 
447  if (fitOption.StoreResult)
448  {
449  TString name = "TFitResult-";
450  name = name + h1->GetName() + "-" + f1->GetName();
451  TString title = "TFitResult-";
452  title += h1->GetTitle();
453  tfr->SetName(name);
454  tfr->SetTitle(title);
455  return TFitResultPtr(tfr);
456  }
457  else
458  return TFitResultPtr(iret);
459 }
460 
461 
463  // get range from histogram and update the DataRange class
464  // if a ranges already exist in that dimension use that one
465 
466  Int_t ndim = GetDimension(h1);
467 
468  double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
469  if (range.Size(0) == 0) {
470  TAxis & xaxis = *(h1->GetXaxis());
471  Int_t hxfirst = xaxis.GetFirst();
472  Int_t hxlast = xaxis.GetLast();
473  Double_t binwidx = xaxis.GetBinWidth(hxlast);
474  xmin = xaxis.GetBinLowEdge(hxfirst);
475  xmax = xaxis.GetBinLowEdge(hxlast) +binwidx;
476  range.AddRange(xmin,xmax);
477  }
478 
479  if (ndim > 1) {
480  if (range.Size(1) == 0) {
481  TAxis & yaxis = *(h1->GetYaxis());
482  Int_t hyfirst = yaxis.GetFirst();
483  Int_t hylast = yaxis.GetLast();
484  Double_t binwidy = yaxis.GetBinWidth(hylast);
485  ymin = yaxis.GetBinLowEdge(hyfirst);
486  ymax = yaxis.GetBinLowEdge(hylast) +binwidy;
487  range.AddRange(1,ymin,ymax);
488  }
489  }
490  if (ndim > 2) {
491  if (range.Size(2) == 0) {
492  TAxis & zaxis = *(h1->GetZaxis());
493  Int_t hzfirst = zaxis.GetFirst();
494  Int_t hzlast = zaxis.GetLast();
495  Double_t binwidz = zaxis.GetBinWidth(hzlast);
496  zmin = zaxis.GetBinLowEdge(hzfirst);
497  zmax = zaxis.GetBinLowEdge(hzlast) +binwidz;
498  range.AddRange(2,zmin,zmax);
499  }
500  }
501 #ifdef DEBUG
502  std::cout << "xmin,xmax" << xmin << " " << xmax << std::endl;
503 #endif
504 
505 }
506 
508  // get range for graph (used sub-set histogram)
509  // N.B. : this is different than in previous implementation of TGraph::Fit where range used was from xmin to xmax.
510  TH1 * h1 = gr->GetHistogram();
511  // an histogram is normally always returned for a TGraph
512  if (h1) HFit::GetDrawingRange(h1, range);
513 }
515  // get range for multi-graph (used sub-set histogram)
516  // N.B. : this is different than in previous implementation of TMultiGraph::Fit where range used was from data xmin to xmax.
517  TH1 * h1 = mg->GetHistogram();
518  if (h1) {
519  HFit::GetDrawingRange(h1, range);
520  }
521  else if (range.Size(0) == 0) {
522  // compute range from all the TGraph's belonging to the MultiGraph
523  double xmin = std::numeric_limits<double>::infinity();
524  double xmax = -std::numeric_limits<double>::infinity();
525  TIter next(mg->GetListOfGraphs() );
526  TGraph * g = 0;
527  while ( (g = (TGraph*) next() ) ) {
528  double x1 = 0, x2 = 0, y1 = 0, y2 = 0;
529  g->ComputeRange(x1,y1,x2,y2);
530  if (x1 < xmin) xmin = x1;
531  if (x2 > xmax) xmax = x2;
532  }
533  range.AddRange(xmin,xmax);
534  }
535 }
537  // get range for graph2D (used sub-set histogram)
538  // N.B. : this is different than in previous implementation of TGraph2D::Fit. There range used was always(0,0)
539  // cannot use TGraph2D::GetHistogram which makes an interpolation
540  //TH1 * h1 = gr->GetHistogram();
541  //if (h1) HFit::GetDrawingRange(h1, range);
542  // not very efficient (t.b.i.)
543  if (range.Size(0) == 0) {
544  double xmin = gr->GetXmin();
545  double xmax = gr->GetXmax();
546  range.AddRange(0,xmin,xmax);
547  }
548  if (range.Size(1) == 0) {
549  double ymin = gr->GetYmin();
550  double ymax = gr->GetYmax();
551  range.AddRange(1,ymin,ymax);
552  }
553 }
554 
556  // get range from histogram and update the DataRange class
557  // if a ranges already exist in that dimension use that one
558 
559  Int_t ndim = GetDimension(s1);
560 
561  for ( int i = 0; i < ndim; ++i ) {
562  if ( range.Size(i) == 0 ) {
563  TAxis *axis = s1->GetAxis(i);
564  range.AddRange(i, axis->GetXmin(), axis->GetXmax());
565  }
566  }
567 }
568 
569 template<class FitObject>
570 void HFit::StoreAndDrawFitFunction(FitObject * h1, TF1 * f1, const ROOT::Fit::DataRange & range, bool delOldFunction, bool drawFunction, const char *goption) {
571 // - Store fitted function in histogram functions list and draw
572 // should have separate functions for 1,2,3d ? t.b.d in case
573 
574 #ifdef DEBUG
575  std::cout <<"draw and store fit function " << f1->GetName() << std::endl;
576 #endif
577 
578 
579  Int_t ndim = GetDimension(h1);
580  double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
581  if (range.Size(0) ) range.GetRange(0,xmin,xmax);
582  if (range.Size(1) ) range.GetRange(1,ymin,ymax);
583  if (range.Size(2) ) range.GetRange(2,zmin,zmax);
584 
585 
586 #ifdef DEBUG
587  std::cout <<"draw and store fit function " << f1->GetName()
588  << " Range in x = [ " << xmin << " , " << xmax << " ]" << std::endl;
589 #endif
590 
591  TList * funcList = h1->GetListOfFunctions();
592  if (funcList == 0){
593  Error("StoreAndDrawFitFunction","Function list has not been created - cannot store the fitted function");
594  return;
595  }
596 
597  // delete the function in the list only if
598  // the function we are fitting is not in that list
599  // If this is the case we re-use that function object and
600  // we do not create a new one (if delOldFunction is true)
601  bool reuseOldFunction = false;
602  if (delOldFunction) {
603  TIter next(funcList, kIterBackward);
604  TObject *obj;
605  while ((obj = next())) {
606  if (obj->InheritsFrom(TF1::Class())) {
607  if (obj != f1) {
608  funcList->Remove(obj);
609  delete obj;
610  }
611  else {
612  reuseOldFunction = true;
613  }
614  }
615  }
616  }
617 
618  TF1 *fnew1 = 0;
619  TF2 *fnew2 = 0;
620  TF3 *fnew3 = 0;
621 
622  // copy TF1 using TClass to avoid slicing in case of derived classes
623  if (ndim < 2) {
624  if (!reuseOldFunction) {
625  fnew1 = (TF1*)f1->IsA()->New();
626  R__ASSERT(fnew1);
627  f1->Copy(*fnew1);
628  funcList->Add(fnew1);
629  }
630  else {
631  fnew1 = f1;
632  }
633  fnew1->SetParent( h1 );
634  fnew1->SetRange(xmin,xmax);
635  fnew1->Save(xmin,xmax,0,0,0,0);
636  if (!drawFunction) fnew1->SetBit(TF1::kNotDraw);
637  fnew1->AddToGlobalList(false);
638  } else if (ndim < 3) {
639  if (!reuseOldFunction) {
640  fnew2 = (TF2*)f1->IsA()->New();
641  R__ASSERT(fnew2);
642  f1->Copy(*fnew2);
643  funcList->Add(fnew2);
644  }
645  else {
646  fnew2 = dynamic_cast<TF2*>(f1);
647  R__ASSERT(fnew2);
648  }
649  fnew2->SetRange(xmin,ymin,xmax,ymax);
650  fnew2->SetParent( h1 );
651  fnew2->Save(xmin,xmax,ymin,ymax,0,0);
652  if (!drawFunction) fnew2->SetBit(TF1::kNotDraw);
653  fnew2->AddToGlobalList(false);
654  } else {
655  if (!reuseOldFunction) {
656  fnew3 = (TF3*)f1->IsA()->New();
657  R__ASSERT(fnew3);
658  f1->Copy(*fnew3);
659  funcList->Add(fnew3);
660  }
661  else {
662  fnew2 = dynamic_cast<TF3*>(f1);
663  R__ASSERT(fnew3);
664  }
665  fnew3->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
666  fnew3->SetParent( h1 );
667  fnew3->Save(xmin,xmax,ymin,ymax,zmin,zmax);
668  if (!drawFunction) fnew3->SetBit(TF1::kNotDraw);
669  fnew3->AddToGlobalList(false);
670  }
671  if (h1->TestBit(kCanDelete)) return;
672  // draw only in case of histograms
673  if (drawFunction && ndim < 3 && h1->InheritsFrom(TH1::Class() ) ) {
674  // no need to re-draw the histogram if the histogram is already in the pad
675  // in that case the function will be just drawn (if option N is not set)
676  if (!gPad || (gPad && gPad->GetListOfPrimitives()->FindObject(h1) == NULL ) )
677  h1->Draw(goption);
678  }
679  if (gPad) gPad->Modified(); // this is not in TH1 code (needed ??)
680 
681  return;
682 }
683 
684 
685 void ROOT::Fit::FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption) {
686  // - Decode list of options into fitOption (used by both TGraph and TH1)
687  // works for both histograms and graph depending on the enum FitObjectType defined in HFit
690  }
691 
692  if (option == 0) return;
693  if (!option[0]) return;
694 
695  TString opt = option;
696  opt.ToUpper();
697 
698  // parse firt the specific options
699  if (type == kHistogram) {
700 
701  if (opt.Contains("WIDTH")) {
702  fitOption.BinVolume = 1; // scale content by the bin width
703  if (opt.Contains("NORMWIDTH")) {
704  // for variable bins: scale content by the bin width normalized by a reference value (typically the minimum bin)
705  // this option is for variable bin widths
706  fitOption.BinVolume = 2;
707  opt.ReplaceAll("NORMWIDTH","");
708  }
709  else
710  opt.ReplaceAll("WIDTH","");
711  }
712 
713  // if (opt.Contains("MULTIPROC")) {
714  // fitOption.ExecPolicy = ROOT::Fit::kMultiprocess;
715  // opt.ReplaceAll("MULTIPROC","");
716  // }
717 
718  if (opt.Contains("SERIAL")) {
720  opt.ReplaceAll("SERIAL","");
721  }
722 
723  if (opt.Contains("MULTITHREAD")) {
725  opt.ReplaceAll("MULTITHREAD","");
726  }
727 
728  if (opt.Contains("I")) fitOption.Integral= 1; // integral of function in the bin (no sense for graph)
729  if (opt.Contains("WW")) fitOption.W1 = 2; //all bins have weight=1, even empty bins
730  }
731 
732  // specific Graph options (need to be parsed before)
733  else if (type == kGraph) {
734  opt.ReplaceAll("ROB", "H");
735  opt.ReplaceAll("EX0", "T");
736 
737  //for robust fitting, see if # of good points is defined
738  // decode parameters for robust fitting
739  Double_t h=0;
740  if (opt.Contains("H=0.")) {
741  int start = opt.Index("H=0.");
742  int numpos = start + strlen("H=0.");
743  int numlen = 0;
744  int len = opt.Length();
745  while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
746  TString num = opt(numpos,numlen);
747  opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
748  h = atof(num.Data());
749  h*=TMath::Power(10, -numlen);
750  }
751 
752  if (opt.Contains("H")) { fitOption.Robust = 1; fitOption.hRobust = h; }
753  if (opt.Contains("T")) fitOption.NoErrX = 1; // no error in X
754 
755  }
756 
757  if (opt.Contains("U")) fitOption.User = 1;
758  if (opt.Contains("Q")) fitOption.Quiet = 1;
759  if (opt.Contains("V")) {fitOption.Verbose = 1; fitOption.Quiet = 0;}
760  if (opt.Contains("L")) fitOption.Like = 1;
761  if (opt.Contains("X")) fitOption.Chi2 = 1;
762  if (opt.Contains("P")) fitOption.PChi2 = 1;
763 
764 
765  // likelihood fit options
766  if (fitOption.Like == 1) {
767  //if (opt.Contains("LL")) fitOption.Like = 2;
768  if (opt.Contains("W")){ fitOption.Like = 2; fitOption.W1=0;}// (weighted likelihood)
769  if (opt.Contains("MULTI")) {
770  if (fitOption.Like == 2) fitOption.Like = 6; // weighted multinomial
771  else fitOption.Like = 4; // multinomial likelihood fit instead of Poisson
772  opt.ReplaceAll("MULTI","");
773  }
774  // in case of histogram give precedence for likelihood options
775  if (type == kHistogram) {
776  if (fitOption.Chi2 == 1 || fitOption.PChi2 == 1)
777  Warning("Fit","Cannot use P or X option in combination of L. Ignore the chi2 option and perform a likelihood fit");
778  }
779 
780  } else {
781  if (opt.Contains("W")) fitOption.W1 = 1; // all non-empty bins have weight =1 (for chi2 fit)
782  }
783 
784 
785  if (opt.Contains("E")) fitOption.Errors = 1;
786  if (opt.Contains("R")) fitOption.Range = 1;
787  if (opt.Contains("G")) fitOption.Gradient= 1;
788  if (opt.Contains("M")) fitOption.More = 1;
789  if (opt.Contains("N")) fitOption.Nostore = 1;
790  if (opt.Contains("0")) fitOption.Nograph = 1;
791  if (opt.Contains("+")) fitOption.Plus = 1;
792  if (opt.Contains("B")) fitOption.Bound = 1;
793  if (opt.Contains("C")) fitOption.Nochisq = 1;
794  if (opt.Contains("F")) fitOption.Minuit = 1;
795  if (opt.Contains("S")) fitOption.StoreResult = 1;
796 
797 }
798 
800  if (foption.Like) {
801  Info("CheckGraphFitOptions","L (Log Likelihood fit) is an invalid option when fitting a graph. It is ignored");
802  foption.Like = 0;
803  }
804  if (foption.Integral) {
805  Info("CheckGraphFitOptions","I (use function integral) is an invalid option when fitting a graph. It is ignored");
806  foption.Integral = 0;
807  }
808  return;
809 }
810 
811 // implementation of unbin fit function (defined in HFitInterface)
812 
814  // do unbin fit, ownership of fitdata is passed later to the TBackFitter class
815 
816  // create a shared pointer to the fit data to managed it
817  std::shared_ptr<ROOT::Fit::UnBinData> fitdata(data);
818 
819 #ifdef DEBUG
820  printf("tree data size is %d \n",fitdata->Size());
821  for (unsigned int i = 0; i < fitdata->Size(); ++i) {
822  if (fitdata->NDim() == 1) printf(" x[%d] = %f \n", i,*(fitdata->Coords(i) ) );
823  }
824 #endif
825  if (fitdata->Size() == 0 ) {
826  Warning("Fit","Fit data is empty ");
827  return -1;
828  }
829 
830  // create an empty TFitResult
831  std::shared_ptr<TFitResult> tfr(new TFitResult() );
832  // create the fitter
833  std::shared_ptr<ROOT::Fit::Fitter> fitter(new ROOT::Fit::Fitter(tfr) );
834  ROOT::Fit::FitConfig & fitConfig = fitter->Config();
835 
836  // dimension is given by data because TF1 pointer can have wrong one
837  unsigned int dim = fitdata->NDim();
838 
839  // set the fit function
840  // if option grad is specified use gradient
841  // need to create a wrapper for an automatic normalized TF1 ???
842  if ( fitOption.Gradient ) {
843  assert ( (int) dim == fitfunc->GetNdim() );
844  fitter->SetFunction(ROOT::Math::WrappedMultiTF1(*fitfunc) );
845  }
846  else
847  fitter->SetFunction(static_cast<const ROOT::Math::IParamMultiFunction &>(ROOT::Math::WrappedMultiTF1(*fitfunc, dim) ) );
848 
849  // parameter setting is done automaticaly in the Fitter class
850  // need only to set limits
851  int npar = fitfunc->GetNpar();
852  for (int i = 0; i < npar; ++i) {
853  ROOT::Fit::ParameterSettings & parSettings = fitConfig.ParSettings(i);
854  double plow,pup;
855  fitfunc->GetParLimits(i,plow,pup);
856  // this is a limitation of TF1 interface - cannot fix a parameter to zero value
857  if (plow*pup != 0 && plow >= pup) {
858  parSettings.Fix();
859  }
860  else if (plow < pup ) {
861  if (!TMath::Finite(pup) && TMath::Finite(plow) )
862  parSettings.SetLowerLimit(plow);
863  else if (!TMath::Finite(plow) && TMath::Finite(pup) )
864  parSettings.SetUpperLimit(pup);
865  else
866  parSettings.SetLimits(plow,pup);
867  }
868 
869  // set the parameter step size (by default are set to 0.3 of value)
870  // if function provides meaningful error values
871  double err = fitfunc->GetParError(i);
872  if ( err > 0)
873  parSettings.SetStepSize(err);
874  else if (plow < pup && TMath::Finite(plow) && TMath::Finite(pup) ) { // in case of limits improve step sizes
875  double step = 0.1 * (pup - plow);
876  // check if value is not too close to limit otherwise trim value
877  if ( parSettings.Value() < pup && pup - parSettings.Value() < 2 * step )
878  step = (pup - parSettings.Value() ) / 2;
879  else if ( parSettings.Value() > plow && parSettings.Value() - plow < 2 * step )
880  step = (parSettings.Value() - plow ) / 2;
881 
882  parSettings.SetStepSize(step);
883  }
884 
885  }
886 
887  fitConfig.SetMinimizerOptions(minOption);
888 
889  if (fitOption.Verbose) fitConfig.MinimizerOptions().SetPrintLevel(3);
890  if (fitOption.Quiet) fitConfig.MinimizerOptions().SetPrintLevel(0);
891 
892  // more
893  if (fitOption.More) fitConfig.SetMinimizer("Minuit","MigradImproved");
894 
895  // chech if Minos or more options
896  if (fitOption.Errors) {
897  // run Hesse and Minos
898  fitConfig.SetParabErrors(true);
899  fitConfig.SetMinosErrors(true);
900  }
901  // use weight correction
902  if ( (fitOption.Like & 2) == 2)
903  fitConfig.SetWeightCorrection(true);
904 
905  bool extended = (fitOption.Like & 1) == 1;
906 
907  bool fitok = false;
908  fitok = fitter->LikelihoodFit(fitdata, extended, fitOption.ExecPolicy);
909  if ( !fitok && !fitOption.Quiet )
910  Warning("UnBinFit","Abnormal termination of minimization.");
911 
912  const ROOT::Fit::FitResult & fitResult = fitter->Result();
913  // one could set directly the fit result in TF1
914  int iret = fitResult.Status();
915  if (!fitResult.IsEmpty() ) {
916  // set in fitfunc the result of the fit
917  fitfunc->SetNDF(fitResult.Ndf() );
918  fitfunc->SetNumberFitPoints(fitdata->Size() );
919 
920  assert( (Int_t)fitResult.Parameters().size() >= fitfunc->GetNpar() );
921  fitfunc->SetParameters( const_cast<double*>(&(fitResult.Parameters().front())));
922  if ( int( fitResult.Errors().size()) >= fitfunc->GetNpar() )
923  fitfunc->SetParErrors( &(fitResult.Errors().front()) );
924 
925  }
926 
927  // store result in the backward compatible VirtualFitter
928  // in case not running in a multi-thread enabled mode
929  if (gGlobalMutex) {
930  TVirtualFitter * lastFitter = TVirtualFitter::GetFitter();
931  // pass ownership of Fitter and Fitdata to TBackCompFitter (fitter pointer cannot be used afterwards)
932  TBackCompFitter * bcfitter = new TBackCompFitter(fitter, fitdata);
933  // cannot use anymore now fitdata (given away ownership)
934  fitdata = 0;
935  bcfitter->SetFitOption(fitOption);
936  //bcfitter->SetObjectFit(fTree);
937  bcfitter->SetUserFunc(fitfunc);
938 
939  if (lastFitter) delete lastFitter;
940  TVirtualFitter::SetFitter( bcfitter );
941 
942  // use old-style for printing the results
943  // if (fitOption.Verbose) bcfitter->PrintResults(2,0.);
944  // else if (!fitOption.Quiet) bcfitter->PrintResults(1,0.);
945 
946  }
947  // print results
948  if (fitOption.Verbose) fitResult.PrintCovMatrix(std::cout);
949  else if (!fitOption.Quiet) fitResult.Print(std::cout);
950 
951  if (fitOption.StoreResult)
952  {
953  TString name = "TFitResult-";
954  name = name + "UnBinData-" + fitfunc->GetName();
955  TString title = "TFitResult-";
956  title += name;
957  tfr->SetName(name);
958  tfr->SetTitle(title);
959  return TFitResultPtr(tfr);
960  }
961  else
962  return TFitResultPtr(iret);
963 }
964 
965 
966 // implementations of ROOT::Fit::FitObject functions (defined in HFitInterface) in terms of the template HFit::Fit
967 
969 moption, const char *goption, ROOT::Fit::DataRange & range) {
970  // check fit options
971  // check if have weights in case of weighted likelihood
972  if ( ((foption.Like & 2) == 2) && h1->GetSumw2N() == 0) {
973  Warning("HFit::FitObject","A weighted likelihood fit is requested but histogram is not weighted - do a standard Likelihood fit");
974  foption.Like = 1;
975  }
976  // histogram fitting
977  return HFit::Fit(h1,f1,foption,moption,goption,range);
978 }
979 
980 TFitResultPtr ROOT::Fit::FitObject(TGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
981  // exclude options not valid for graphs
983  // TGraph fitting
984  return HFit::Fit(gr,f1,foption,moption,goption,range);
985 }
986 
987 TFitResultPtr ROOT::Fit::FitObject(TMultiGraph * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
988  // exclude options not valid for graphs
990  // TMultiGraph fitting
991  return HFit::Fit(gr,f1,foption,moption,goption,range);
992 }
993 
994 TFitResultPtr ROOT::Fit::FitObject(TGraph2D * gr, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
995  // exclude options not valid for graphs
997  // TGraph2D fitting
998  return HFit::Fit(gr,f1,foption,moption,goption,range);
999 }
1000 
1001 TFitResultPtr ROOT::Fit::FitObject(THnBase * s1, TF1 *f1 , Foption_t & foption , const ROOT::Math::MinimizerOptions & moption, const char *goption, ROOT::Fit::DataRange & range) {
1002  // sparse histogram fitting
1003  return HFit::Fit(s1,f1,foption,moption,goption,range);
1004 }
1005 
1006 
1007 
1008 // Int_t TGraph2D::DoFit(TF2 *f2 ,Option_t *option ,Option_t *goption) {
1009 // // internal graph2D fitting methods
1010 // Foption_t fitOption;
1011 // ROOT::Fit::FitOptionsMake(option,fitOption);
1012 
1013 // // create range and minimizer options with default values
1014 // ROOT::Fit::DataRange range(2);
1015 // ROOT::Math::MinimizerOptions minOption;
1016 // return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
1017 // }
1018 
1019 
1020 // function to compute the simple chi2 for graphs and histograms
1021 
1022 double ROOT::Fit::Chisquare(const TH1 & h1, TF1 & f1, bool useRange, bool usePL) {
1023  return HFit::ComputeChi2(h1,f1,useRange, usePL);
1024 }
1025 
1026 double ROOT::Fit::Chisquare(const TGraph & g, TF1 & f1, bool useRange) {
1027  return HFit::ComputeChi2(g,f1, useRange, false);
1028 }
1029 
1030 template<class FitObject>
1031 double HFit::ComputeChi2(const FitObject & obj, TF1 & f1, bool useRange, bool usePL ) {
1032 
1033  // implement using the fitting classes
1035  if (usePL) opt.fUseEmpty=true;
1036  ROOT::Fit::DataRange range;
1037  // get range of function
1038  if (useRange) HFit::GetFunctionRange(f1,range);
1039  // fill the data set
1040  ROOT::Fit::BinData data(opt,range);
1041  ROOT::Fit::FillData(data, &obj, &f1);
1042  if (data.Size() == 0 ) {
1043  Warning("Chisquare","data set is empty - return -1");
1044  return -1;
1045  }
1047  if (usePL) {
1048  // use the poisson log-lokelihood (Baker-Cousins chi2)
1049  ROOT::Fit::PoissonLLFunction nll(data, wf1);
1050  return 2.* nll( f1.GetParameters() ) ;
1051  }
1052  ROOT::Fit::Chi2Function chi2(data, wf1);
1053  return chi2(f1.GetParameters() );
1054 
1055 }
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:638
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:455
bool IsVectorized()
Definition: TF1.h:434
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:1022
void SetTolerance(double tol)
set the tolerance
virtual void Copy(TObject &f1) const
Copy this F1 to a new F1.
Definition: TF1.cxx:994
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:529
virtual void SetMethodCall(TMethodCall *m)
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:515
float ymin
Definition: THbookFile.cxx:93
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition: FitResult.h:169
#define g(i)
Definition: RSha256.hxx:105
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
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:203
int Verbose
Definition: Foption.h:30
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:36
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF1.cxx:3535
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1138
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:3150
#define R__ASSERT(e)
Definition: TError.h:96
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:610
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:634
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:163
Basic string class.
Definition: TString.h:131
virtual void SetNumberFitPoints(Int_t npfits)
Definition: TF1.h:618
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.
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1927
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:725
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:693
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:813
double hRobust
Definition: Foption.h:51
virtual Int_t GetDimension() const
Definition: TH1.h:278
virtual void SetParent(TObject *p=0)
Definition: TF1.h:659
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:32
class evaluating the log likelihood for binned Poisson likelihood fits it is template to distinguish ...
int Chi2
Definition: Foption.h:32
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:462
void Info(const char *location, const char *msgfmt,...)
void Fix()
fix the parameter
Int_t Finite(Double_t x)
Check if it is finite with a mask in order to be consistent in presence of fast math.
Definition: TMath.h:761
virtual Int_t GetNdim() const
Definition: TF1.h:479
static constexpr double s
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 GetDimension(const THnBase *s1)
Definition: HFitImpl.cxx:55
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:230
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:606
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.
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:466
Class to manage histogram axis.
Definition: TAxis.h:30
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2980
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:819
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:442
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
static constexpr double mg
void FitOptionsMake(const char *option, Foption_t &fitOption)
Ssiz_t Length() const
Definition: TString.h:405
virtual Bool_t IsLinear() const
Definition: TF1.h:596
static TVirtualFitter * GetFitter()
static: return the current Fitter
int More
Definition: Foption.h:38
#define s1(x)
Definition: RSha256.hxx:91
int Like
Definition: Foption.h:34
TAxis * GetYaxis()
Definition: TH1.h:317
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:137
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
#define h(i)
Definition: RSha256.hxx:106
void SetMinimizerOptions(const ROOT::Math::MinimizerOptions &minopt)
set all the minimizer options using class MinimizerOptions
Definition: FitConfig.cxx:256
const Bool_t kFALSE
Definition: RtypesCore.h:88
TString & Remove(Ssiz_t pos)
Definition: TString.h:668
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:310
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:47
void PrintCovMatrix(std::ostream &os) const
print error matrix and correlations
Definition: FitResult.cxx:486
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:1473
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...
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:174
void SetWeightCorrection(bool on=true)
apply the weight correction for error matric computation
Definition: FitConfig.h:233
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:619
The TH1 histogram class.
Definition: TH1.h:56
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:685
int Nochisq
Definition: Foption.h:45
Abstract Base Class for Fitting.
TAxis * GetZaxis()
Definition: TH1.h:318
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
Definition: FitResult.h:117
Mother of all ROOT objects.
Definition: TObject.h:37
int NoErrX
Definition: Foption.h:47
virtual Int_t GetNpar() const
Definition: TF1.h:475
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition: TF1.cxx:1937
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition: TAxis.cxx:537
void StoreAndDrawFitFunction(FitObject *h1, TF1 *f1, const ROOT::Fit::DataRange &range, bool, bool, const char *goption)
Definition: HFitImpl.cxx:570
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
Definition: FitResult.cxx:428
std::string ToString(const T &val)
Utility function for conversion to strings.
Definition: Util.h:50
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition: TF3.h:146
virtual void Add(TObject *obj)
Definition: TList.h:87
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:3500
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition: TF1.cxx:2280
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:831
1-Dim function class
Definition: TF1.h:210
void SetNormErrors(bool on=true)
set the option to normalize the error on the result according to chi2/ndf
Definition: FitConfig.h:224
void SetUpperLimit(double up)
set a single upper limit
A TGraph is an 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:1031
TF1 * f1
Definition: legend1.C:11
#define gPad
Definition: TVirtualPad.h:287
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:3435
virtual Double_t * GetParameters() const
Definition: TF1.h:514
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:968
const Bool_t kIterBackward
Definition: TCollection.h:41
void SetParabErrors(bool on=true)
set parabolic erros
Definition: FitConfig.h:227
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
TList * GetListOfFunctions() const
Definition: TH1.h:239
const Bool_t kTRUE
Definition: RtypesCore.h:87
Double_t GetXmax() const
Definition: TAxis.h:134
virtual Int_t GetNumber() const
Definition: TF1.h:492
double Chi2() const
Chi2 fit value in case of likelihood must be computed ?
Definition: FitResult.h:160
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:316
void CheckGraphFitOptions(Foption_t &fitOption)
Definition: HFitImpl.cxx:799
Class describing the configuration of the fit, options and parameter settings using the ROOT::Fit::Pa...
Definition: FitConfig.h:46
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
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:775
const char * Data() const
Definition: TString.h:364