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