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