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