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