Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TBackCompFitter.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Lorenzo Moneta
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2012, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 ////////////////////////////////////////////////////////////////////////////////
13 /** \class TBackCompFitter
14  \ingroup Hist
15  \brief Backward compatible implementation of TVirtualFitter
16 
17 Backward compatible implementation of TVirtualFitter using the
18 class ROOT::Fit::Fitter. This class is created after fitting an
19 histogram (TH1), TGraph or TTree and provides in addition to the
20 methods of the TVirtualFitter hooks to access the fit result class
21 (ROOT::Fit::FitResult), the fit configuration
22 (ROOT::Fit::FitConfig) or the fit data (ROOT::Fit::FitData) using
23 
24 ~~~~~~~~{.cpp}
25  TBackCompFitter * fitter = (TBackCompFitter *) TVirtualFitter::GetFitter();
26  ROOT::Fit::FitResult & result = fitter->GetFitResult();
27  result.Print(std::cout);
28 ~~~~~~~~
29 
30 Methods for getting the confidence level or contours are also
31 provided. Note that after a new calls to TH1::Fit (or similar) the
32 class will be deleted and all reference to the FitResult, FitConfig
33 or minimizer will be invalid. One could eventually copying the
34 class before issuing a new fit to avoid deleting this information.
35 *///////////////////////////////////////////////////////////////////////////////
36 
37 #include "TROOT.h"
38 #include "TBackCompFitter.h"
39 
40 #include "Math/Util.h"
41 
42 #include <iostream>
43 #include <cassert>
44 
45 //needed by GetCondifenceLevel
46 #include "Math/IParamFunction.h"
47 #include "TH1.h"
48 #include "TH2.h"
49 #include "TH3.h"
50 #include "TMath.h"
51 #include "TGraph.h"
52 #include "TGraphErrors.h"
53 #include "TGraph2DErrors.h"
54 #include "TMultiGraph.h"
55 #include "HFitInterface.h"
56 #include "Math/Minimizer.h"
57 #include "Fit/BinData.h"
58 #include "Fit/FitConfig.h"
59 #include "Fit/UnBinData.h"
61 #include "Fit/LogLikelihoodFCN.h"
62 #include "Fit/Chi2FCN.h"
63 #include "Fit/FcnAdapter.h"
64 #include "TFitResult.h"
65 
66 //#define DEBUG 1
67 
68 
69 ClassImp(TBackCompFitter);
70 
71 
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// Constructor needed by TVirtualFitter interface. Same behavior as default constructor.
75 /// initialize setting name and the global pointer
76 
77 TBackCompFitter::TBackCompFitter( ) :
78  fMinimizer(0),
79  fObjFunc(0),
80  fModelFunc(0)
81 {
82  SetName("BCFitter");
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Constructor used after having fit using directly ROOT::Fit::Fitter
87 /// will create a dummy fitter copying configuration and parameter settings
88 
89 TBackCompFitter::TBackCompFitter(const std::shared_ptr<ROOT::Fit::Fitter> & fitter, const std::shared_ptr<ROOT::Fit::FitData> & data) :
90  fFitData(data),
91  fFitter(fitter),
92  fMinimizer(0),
93  fObjFunc(0),
94  fModelFunc(0)
95 {
96  SetName("LastFitter");
97 }
98 
99 ////////////////////////////////////////////////////////////////////////////////
100 /// Destructor - delete the managed objects
101 
102 TBackCompFitter::~TBackCompFitter() {
103  if (fMinimizer) delete fMinimizer;
104  if (fObjFunc) delete fObjFunc;
105  if (fModelFunc) delete fModelFunc;
106 }
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 /// Do chisquare calculations in case of likelihood fits
110 /// Do evaluation a the minimum only
111 
112 Double_t TBackCompFitter::Chisquare(Int_t npar, Double_t *params) const {
113  const std::vector<double> & minpar = fFitter->Result().Parameters();
114  assert (npar == (int) minpar.size() );
115  double diff = 0;
116  double s = 0;
117  for (int i =0; i < npar; ++i) {
118  diff += std::abs( params[i] - minpar[i] );
119  s += minpar[i];
120  }
121 
122  if (diff > s * 1.E-12 ) Warning("Chisquare","given parameter values are not at minimum - chi2 at minimum is returned");
123  return fFitter->Result().Chi2();
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Clear resources for consecutive fits
128 
129 void TBackCompFitter::Clear(Option_t*) {
130  // need to do something here ??? to be seen
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// Execute the command (Fortran Minuit compatible interface)
135 
136 Int_t TBackCompFitter::ExecuteCommand(const char *command, Double_t *args, Int_t nargs) {
137 #ifdef DEBUG
138  std::cout<<"Execute command= "<<command<<std::endl;
139 #endif
140 
141  // set also number of parameters in obj function
142  DoSetDimension();
143 
144  TString scommand(command);
145  scommand.ToUpper();
146 
147  // MIGRAD
148  if (scommand.Contains("MIG")) {
149  if (!fObjFunc) {
150  Error("ExecuteCommand","FCN must set before executing this command");
151  return -1;
152  }
153 
154  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad");
155  bool ret = fFitter->FitFCN(*fObjFunc);
156  return (ret) ? 0 : -1;
157  }
158 
159  //Minimize
160  if (scommand.Contains("MINI")) {
161 
162  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Minimize");
163  if (!fObjFunc) {
164  Error("ExecuteCommand","FCN must set before executing this command");
165  return -1;
166  }
167  bool ret = fFitter->FitFCN(*fObjFunc);
168  return (ret) ? 0 : -1;
169  }
170  //Simplex
171  if (scommand.Contains("SIM")) {
172 
173  if (!fObjFunc) {
174  Error("ExecuteCommand","FCN must set before executing this command");
175  return -1;
176  }
177 
178  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Simplex");
179  bool ret = fFitter->FitFCN(*fObjFunc);
180  return (ret) ? 0 : -1;
181  }
182  //SCan
183  if (scommand.Contains("SCA")) {
184 
185  if (!fObjFunc) {
186  Error("ExecuteCommand","FCN must set before executing this command");
187  return -1;
188  }
189 
190  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Scan");
191  bool ret = fFitter->FitFCN(*fObjFunc);
192  return (ret) ? 0 : -1;
193  }
194  // MINOS
195  else if (scommand.Contains("MINO")) {
196 
197  if (fFitter->Config().MinosErrors() ) return 0;
198 
199  if (!fObjFunc) {
200  Error("ExecuteCommand","FCN must set before executing this command");
201  return -1;
202  }
203  // do only MINOS. need access to minimizer. For the moment re-run fitting with minos options
204  fFitter->Config().SetMinosErrors(true);
205  // set new parameter values
206 
207  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad"); // redo -minimization with Minos
208  bool ret = fFitter->FitFCN(*fObjFunc);
209  return (ret) ? 0 : -1;
210 
211  }
212  //HESSE
213  else if (scommand.Contains("HES")) {
214 
215  if (fFitter->Config().ParabErrors() ) return 0;
216 
217  if (!fObjFunc) {
218  Error("ExecuteCommand","FCN must set before executing this command");
219  return -1;
220  }
221 
222  // do only HESSE. need access to minimizer. For the moment re-run fitting with hesse options
223  fFitter->Config().SetParabErrors(true);
224  fFitter->Config().SetMinimizer(GetDefaultFitter(), "Migrad"); // redo -minimization with Minos
225  bool ret = fFitter->FitFCN(*fObjFunc);
226  return (ret) ? 0 : -1;
227  }
228 
229  // FIX
230  else if (scommand.Contains("FIX")) {
231  for(int i = 0; i < nargs; i++) {
232  FixParameter(int(args[i])-1);
233  }
234  return 0;
235  }
236  // SET LIMIT (upper and lower)
237  else if (scommand.Contains("SET LIM")) {
238  if (nargs < 3) {
239  Error("ExecuteCommand","Invalid parameters given in SET LIMIT");
240  return -1;
241  }
242  int ipar = int(args[0]);
243  if (!ValidParameterIndex(ipar) ) return -1;
244  double low = args[1];
245  double up = args[2];
246  fFitter->Config().ParSettings(ipar).SetLimits(low,up);
247  return 0;
248  }
249  // SET PRINT
250  else if (scommand.Contains("SET PRIN")) {
251  if (nargs < 1) return -1;
252  fFitter->Config().MinimizerOptions().SetPrintLevel(int(args[0]) );
253  return 0;
254  }
255  // SET ERR
256  else if (scommand.Contains("SET ERR")) {
257  if (nargs < 1) return -1;
258  fFitter->Config().MinimizerOptions().SetPrintLevel(int( args[0]) );
259  return 0;
260  }
261  // SET STRATEGY
262  else if (scommand.Contains("SET STR")) {
263  if (nargs < 1) return -1;
264  fFitter->Config().MinimizerOptions().SetStrategy(int(args[0]) );
265  return 0;
266  }
267  //SET GRAD (not impl.)
268  else if (scommand.Contains("SET GRA")) {
269  // not yet available
270  // fGradient = true;
271  return -1;
272  }
273  //SET NOW (not impl.)
274  else if (scommand.Contains("SET NOW")) {
275  // no warning (works only for TMinuit)
276  // fGradient = true;
277  return -1;
278  }
279  // CALL FCN
280  else if (scommand.Contains("CALL FCN")) {
281  // call fcn function (global pointer to free function)
282 
283  if (nargs < 1 || fFCN == 0 ) return -1;
284  int npar = fObjFunc->NDim();
285  // use values in fit result if existing otherwise in ParameterSettings
286  std::vector<double> params(npar);
287  for (int i = 0; i < npar; ++i)
288  params[i] = GetParameter(i);
289 
290  double fval = 0;
291  (*fFCN)(npar, 0, fval, &params[0],int(args[0]) ) ;
292  return 0;
293  }
294  else {
295  // other commands passed
296  Error("ExecuteCommand","Invalid or not supported command given %s",command);
297  return -1;
298  }
299 }
300 
301 ////////////////////////////////////////////////////////////////////////////////
302 /// Check if ipar is a valid parameter index
303 
304 bool TBackCompFitter::ValidParameterIndex(int ipar) const {
305  int nps = fFitter->Config().ParamsSettings().size();
306  if (ipar < 0 || ipar >= nps ) {
307  std::string msg = ROOT::Math::Util::ToString(ipar) + " is an invalid Parameter index";
308  Error("ValidParameterIndex","%s",msg.c_str());
309  return false;
310  }
311  return true;
312 }
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Fix the parameter
316 
317 void TBackCompFitter::FixParameter(Int_t ipar) {
318  if (ValidParameterIndex(ipar) )
319  fFitter->Config().ParSettings(ipar).Fix();
320 }
321 
322 ////////////////////////////////////////////////////////////////////////////////
323 /// Computes point-by-point confidence intervals for the fitted function.
324 /// \param n number of points
325 /// \param ndim dimensions of points
326 /// \param x points, at which to compute the intervals, for ndim > 1
327 /// should be in order: (x0,y0, x1, y1, ... xn, yn)
328 /// \param ci computed intervals are returned in this array
329 /// \param cl confidence level, default=0.95
330 ///
331 /// NOTE, that the intervals are approximate for nonlinear(in parameters) models
332 
333 void TBackCompFitter::GetConfidenceIntervals(Int_t n, Int_t ndim, const Double_t *x, Double_t *ci, Double_t cl)
334 {
335  if (!fFitter->Result().IsValid()) {
336  Error("GetConfidenceIntervals","Cannot compute confidence intervals with an invalide fit result");
337  return;
338  }
339 
340  fFitter->Result().GetConfidenceIntervals(n,ndim,1,x,ci,cl,false);
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// Computes confidence intervals at level cl. Default is 0.95
345 /// The TObject parameter can be a TGraphErrors, a TGraph2DErrors or a TH1,2,3.
346 /// For Graphs, confidence intervals are computed for each point,
347 /// the value of the graph at that point is set to the function value at that
348 /// point, and the graph y-errors (or z-errors) are set to the value of
349 /// the confidence interval at that point.
350 /// For Histograms, confidence intervals are computed for each bin center
351 /// The bin content of this bin is then set to the function value at the bin
352 /// center, and the bin error is set to the confidence interval value.
353 //
354 /// NOTE: confidence intervals are approximate for nonlinear models!
355 ///
356 /// Allowed combinations:
357 ///
358 /// Fitted object | Passed object
359 /// --------------------------|------------------
360 /// TGraph | TGraphErrors, TH1
361 /// TGraphErrors, AsymmErrors | TGraphErrors, TH1
362 /// TH1 | TGraphErrors, TH1
363 /// TGraph2D | TGraph2DErrors, TH2
364 /// TGraph2DErrors | TGraph2DErrors, TH2
365 /// TH2 | TGraph2DErrors, TH2
366 /// TH3 | TH3
367 
368 void TBackCompFitter::GetConfidenceIntervals(TObject *obj, Double_t cl)
369 {
370  if (!fFitter->Result().IsValid() ) {
371  Error("GetConfidenceIntervals","Cannot compute confidence intervals with an invalide fit result");
372  return;
373  }
374 
375  // get data dimension from fit object
376  int datadim = 1;
377  TObject * fitobj = GetObjectFit();
378  if (!fitobj) {
379  Error("GetConfidenceIntervals","Cannot compute confidence intervals without a fitting object");
380  return;
381  }
382 
383  if (fitobj->InheritsFrom(TGraph2D::Class())) datadim = 2;
384  if (fitobj->InheritsFrom(TH1::Class())) {
385  TH1 * h1 = dynamic_cast<TH1*>(fitobj);
386  assert(h1 != 0);
387  datadim = h1->GetDimension();
388  }
389 
390  if (datadim == 1) {
391  if (!obj->InheritsFrom(TGraphErrors::Class()) && !obj->InheritsFrom(TH1::Class() ) ) {
392  Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TGraphErrors or a TH1");
393  return;
394  }
395  }
396  if (datadim == 2) {
397  if (!obj->InheritsFrom(TGraph2DErrors::Class()) && !obj->InheritsFrom(TH2::Class() ) ) {
398  Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TGraph2DErrors or a TH2");
399  return;
400  }
401  }
402  if (datadim == 3) {
403  if (!obj->InheritsFrom(TH3::Class() ) ) {
404  Error("GetConfidenceIntervals", "Invalid object passed for storing confidence level data, must be a TH3");
405  return;
406  }
407  }
408 
409  // fill bin data (for the moment use all ranges) according to object passed
410  ROOT::Fit::BinData data;
411  data.Opt().fUseEmpty = true; // need to use all bins of given histograms
412  // call appropriate function according to type of object
413  if (obj->InheritsFrom(TGraph::Class()) )
414  ROOT::Fit::FillData(data, dynamic_cast<TGraph *>(obj) );
415  else if (obj->InheritsFrom(TGraph2D::Class()) )
416  ROOT::Fit::FillData(data, dynamic_cast<TGraph2D *>(obj) );
417 // else if (obj->InheritsFrom(TMultiGraph::Class()) )
418 // ROOT::Fit::FillData(data, dynamic_cast<TMultiGraph *>(obj) );
419  else if (obj->InheritsFrom(TH1::Class()) )
420  ROOT::Fit::FillData(data, dynamic_cast<TH1 *>(obj) );
421 
422 
423  unsigned int n = data.Size();
424 
425  std::vector<double> ci( n );
426 
427  fFitter->Result().GetConfidenceIntervals(data,&ci[0],cl,false);
428 
429  const ROOT::Math::IParamMultiFunction * func = fFitter->Result().FittedFunction();
430  assert(func != 0);
431 
432  // fill now the object with cl data
433  for (unsigned int i = 0; i < n; ++i) {
434  const double * x = data.Coords(i);
435  double y = (*func)( x ); // function is evaluated using its parameters
436 
437  if (obj->InheritsFrom(TGraphErrors::Class()) ) {
438  TGraphErrors * gr = dynamic_cast<TGraphErrors *> (obj);
439  assert(gr != 0);
440  gr->SetPoint(i, *x, y);
441  gr->SetPointError(i, 0, ci[i]);
442  }
443  if (obj->InheritsFrom(TGraph2DErrors::Class()) ) {
444  TGraph2DErrors * gr = dynamic_cast<TGraph2DErrors *> (obj);
445  assert(gr != 0);
446  gr->SetPoint(i, x[0], x[1], y);
447  gr->SetPointError(i, 0, 0, ci[i]);
448  }
449  if (obj->InheritsFrom(TH1::Class()) ) {
450  TH1 * h1 = dynamic_cast<TH1 *> (obj);
451  assert(h1 != 0);
452  int ibin = 0;
453  if (datadim == 1) ibin = h1->FindBin(*x);
454  if (datadim == 2) ibin = h1->FindBin(x[0],x[1]);
455  if (datadim == 3) ibin = h1->FindBin(x[0],x[1],x[2]);
456  h1->SetBinContent(ibin, y);
457  h1->SetBinError(ibin, ci[i]);
458  }
459  }
460 
461 }
462 
463 ////////////////////////////////////////////////////////////////////////////////
464 /// Get the error matrix in a pointer to a NxN array.
465 /// excluding the fixed parameters
466 
467 Double_t* TBackCompFitter::GetCovarianceMatrix() const {
468  unsigned int nfreepar = GetNumberFreeParameters();
469  unsigned int ntotpar = GetNumberTotalParameters();
470 
471  if (fCovar.size() != nfreepar*nfreepar )
472  fCovar.resize(nfreepar*nfreepar);
473 
474  if (!fFitter->Result().IsValid() ) {
475  Warning("GetCovarianceMatrix","Invalid fit result");
476  return 0;
477  }
478 
479  unsigned int l = 0;
480  for (unsigned int i = 0; i < ntotpar; ++i) {
481  if (fFitter->Config().ParSettings(i).IsFixed() ) continue;
482  unsigned int m = 0;
483  for (unsigned int j = 0; j < ntotpar; ++j) {
484  if (fFitter->Config().ParSettings(j).IsFixed() ) continue;
485  unsigned int index = nfreepar*l + m;
486  assert(index < fCovar.size() );
487  fCovar[index] = fFitter->Result().CovMatrix(i,j);
488  m++;
489  }
490  l++;
491  }
492  return &(fCovar.front());
493 }
494 
495 ////////////////////////////////////////////////////////////////////////////////
496 /// Get error matrix element (return all zero if matrix is not available)
497 
498 Double_t TBackCompFitter::GetCovarianceMatrixElement(Int_t i, Int_t j) const {
499  unsigned int np2 = fCovar.size();
500  unsigned int npar = GetNumberFreeParameters();
501  if ( np2 == 0 || np2 != npar *npar ) {
502  double * c = GetCovarianceMatrix();
503  if (c == 0) return 0;
504  }
505  return fCovar[i*npar + j];
506 }
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Get fit errors
510 
511 Int_t TBackCompFitter::GetErrors(Int_t ipar,Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &globcc) const {
512  if (!ValidParameterIndex(ipar) ) return -1;
513 
514  const ROOT::Fit::FitResult & result = fFitter->Result();
515  if (!result.IsValid() ) {
516  Warning("GetErrors","Invalid fit result");
517  return -1;
518  }
519 
520  eparab = result.Error(ipar);
521  eplus = result.UpperError(ipar);
522  eminus = result.LowerError(ipar);
523  globcc = result.GlobalCC(ipar);
524  return 0;
525 }
526 
527 ////////////////////////////////////////////////////////////////////////////////
528 /// Number of total parameters
529 
530 Int_t TBackCompFitter::GetNumberTotalParameters() const {
531  return fFitter->Result().NTotalParameters();
532 }
533 Int_t TBackCompFitter::GetNumberFreeParameters() const {
534  // number of variable parameters
535  return fFitter->Result().NFreeParameters();
536 }
537 
538 ////////////////////////////////////////////////////////////////////////////////
539 /// Parameter error
540 
541 Double_t TBackCompFitter::GetParError(Int_t ipar) const {
542  if (fFitter->Result().IsEmpty() ) {
543  if (ValidParameterIndex(ipar) ) return fFitter->Config().ParSettings(ipar).StepSize();
544  else return 0;
545  }
546  return fFitter->Result().Error(ipar);
547 }
548 
549 ////////////////////////////////////////////////////////////////////////////////
550 /// Parameter value
551 
552 Double_t TBackCompFitter::GetParameter(Int_t ipar) const {
553  if (fFitter->Result().IsEmpty() ) {
554  if (ValidParameterIndex(ipar) ) return fFitter->Config().ParSettings(ipar).Value();
555  else return 0;
556  }
557  return fFitter->Result().Value(ipar);
558 }
559 
560 ////////////////////////////////////////////////////////////////////////////////
561 /// Get all parameter info (name, value, errors)
562 
563 Int_t TBackCompFitter::GetParameter(Int_t ipar,char *name,Double_t &value,Double_t &verr,Double_t &vlow, Double_t &vhigh) const {
564  if (!ValidParameterIndex(ipar) ) {
565  return -1;
566  }
567  const std::string & pname = fFitter->Config().ParSettings(ipar).Name();
568  const char * c = pname.c_str();
569  std::copy(c,c + pname.size(),name);
570 
571  if (fFitter->Result().IsEmpty() ) {
572  value = fFitter->Config().ParSettings(ipar).Value();
573  verr = fFitter->Config().ParSettings(ipar).Value(); // error is step size in this case
574  vlow = fFitter->Config().ParSettings(ipar).LowerLimit(); // vlow is lower limit in this case
575  vhigh = fFitter->Config().ParSettings(ipar).UpperLimit(); // vlow is lower limit in this case
576  return 1;
577  }
578  else {
579  value = fFitter->Result().Value(ipar);
580  verr = fFitter->Result().Error(ipar);
581  vlow = fFitter->Result().LowerError(ipar);
582  vhigh = fFitter->Result().UpperError(ipar);
583  }
584  return 0;
585 }
586 
587 ////////////////////////////////////////////////////////////////////////////////
588 /// Return name of parameter ipar
589 
590 const char *TBackCompFitter::GetParName(Int_t ipar) const {
591  if (!ValidParameterIndex(ipar) ) {
592  return 0;
593  }
594  return fFitter->Config().ParSettings(ipar).Name().c_str();
595 }
596 
597 ////////////////////////////////////////////////////////////////////////////////
598 /// Get fit statistical information
599 
600 Int_t TBackCompFitter::GetStats(Double_t &amin, Double_t &edm, Double_t &errdef, Int_t &nvpar, Int_t &nparx) const {
601  const ROOT::Fit::FitResult & result = fFitter->Result();
602  amin = result.MinFcnValue();
603  edm = result.Edm();
604  errdef = fFitter->Config().MinimizerOptions().ErrorDef();
605  nvpar = result.NFreeParameters();
606  nparx = result.NTotalParameters();
607  return 0;
608 }
609 
610 ////////////////////////////////////////////////////////////////////////////////
611 /// Sum of log (un-needed)
612 
613 Double_t TBackCompFitter::GetSumLog(Int_t) {
614  Warning("GetSumLog","Dummy method - returned 0");
615  return 0.;
616 }
617 
618 ////////////////////////////////////////////////////////////////////////////////
619 /// Query if parameter ipar is fixed
620 
621 Bool_t TBackCompFitter::IsFixed(Int_t ipar) const {
622  if (!ValidParameterIndex(ipar) ) {
623  return false;
624  }
625  return fFitter->Config().ParSettings(ipar).IsFixed();
626 }
627 
628 ////////////////////////////////////////////////////////////////////////////////
629 /// Print the fit result.
630 /// Use PrintResults function in case of Minuit for old -style printing
631 
632 void TBackCompFitter::PrintResults(Int_t level, Double_t ) const {
633  if (fFitter->GetMinimizer() && fFitter->Config().MinimizerType() == "Minuit")
634  fFitter->GetMinimizer()->PrintResults();
635  else {
636  if (level > 0) fFitter->Result().Print(std::cout);
637  if (level > 1) fFitter->Result().PrintCovMatrix(std::cout);
638  }
639  // need to print minos errors and globalCC + other info
640 }
641 
642 ////////////////////////////////////////////////////////////////////////////////
643 /// Release a fit parameter
644 
645 void TBackCompFitter::ReleaseParameter(Int_t ipar) {
646  if (ValidParameterIndex(ipar) )
647  fFitter->Config().ParSettings(ipar).Release();
648 }
649 
650 ////////////////////////////////////////////////////////////////////////////////
651 /// Set fit method (chi2 or likelihood).
652 /// According to the method the appropriate FCN function will be created
653 
654 void TBackCompFitter::SetFitMethod(const char *) {
655  Info("SetFitMethod","non supported method");
656 }
657 
658 ////////////////////////////////////////////////////////////////////////////////
659 /// Set (add) a new fit parameter passing initial value, step size (verr) and parameter limits
660 /// if vlow > vhigh the parameter is unbounded
661 /// if the stepsize (verr) == 0 the parameter is treated as fixed
662 
663 Int_t TBackCompFitter::SetParameter(Int_t ipar,const char *parname,Double_t value,Double_t verr,Double_t vlow, Double_t vhigh) {
664  std::vector<ROOT::Fit::ParameterSettings> & parlist = fFitter->Config().ParamsSettings();
665  if ( ipar >= (int) parlist.size() ) parlist.resize(ipar+1);
666  ROOT::Fit::ParameterSettings ps(parname, value, verr);
667  if (verr == 0) ps.Fix();
668  if (vlow < vhigh) ps.SetLimits(vlow, vhigh);
669  parlist[ipar] = ps;
670  return 0;
671 }
672 
673 //______________________________________________________________________________
674 // static method evaluating FCN
675 // void TBackCompFitter::FCN( int &, double * , double & f, double * x , int /* iflag */) {
676 // // get static instance of fitter
677 // TBackCompFitter * fitter = dynamic_cast<TBackCompFitter *>(TVirtualFitter::GetFitter());
678 // assert(fitter);
679 // if (fitter->fObjFunc == 0) fitter->RecreateFCN();
680 // assert(fitter->fObjFunc);
681 // f = (*(fitter.fObjFunc) )(x);
682 // }
683 
684 ////////////////////////////////////////////////////////////////////////////////
685 /// Recreate a minimizer instance using the function and data
686 /// set objective function in minimizers function to re-create FCN from stored data object and fit options
687 
688 void TBackCompFitter::ReCreateMinimizer() {
689  assert(fFitData.get());
690 
691  // case of standard fits (not made fia Fitter::FitFCN)
692  if (fFitter->Result().FittedFunction() != 0) {
693 
694  if (fModelFunc) delete fModelFunc;
695  fModelFunc = dynamic_cast<ROOT::Math::IParamMultiFunction *>((fFitter->Result().FittedFunction())->Clone());
696  assert(fModelFunc);
697 
698  // create fcn functions, should consider also gradient case
699  const ROOT::Fit::BinData * bindata = dynamic_cast<const ROOT::Fit::BinData *>(fFitData.get());
700  if (bindata) {
701  if (GetFitOption().Like )
702  fObjFunc = new ROOT::Fit::PoissonLikelihoodFCN<ROOT::Math::IMultiGenFunction>(*bindata, *fModelFunc);
703  else
704  fObjFunc = new ROOT::Fit::Chi2FCN<ROOT::Math::IMultiGenFunction>(*bindata, *fModelFunc);
705  }
706  else {
707  const ROOT::Fit::UnBinData * unbindata = dynamic_cast<const ROOT::Fit::UnBinData *>(fFitData.get());
708  assert(unbindata);
709  fObjFunc = new ROOT::Fit::LogLikelihoodFCN<ROOT::Math::IMultiGenFunction>(*unbindata, *fModelFunc);
710  }
711  }
712 
713  // recreate the minimizer
714  fMinimizer = fFitter->Config().CreateMinimizer();
715  if (fMinimizer == 0) {
716  Error("SetMinimizerFunction","cannot create minimizer %s",fFitter->Config().MinimizerType().c_str() );
717  }
718  else {
719  if (!fObjFunc) {
720  Error("SetMinimizerFunction","Object Function pointer is NULL");
721  }
722  else
723  fMinimizer->SetFunction(*fObjFunc);
724  }
725 
726 }
727 
728 ////////////////////////////////////////////////////////////////////////////////
729 /// Override setFCN to use the Adapter to Minuit2 FCN interface
730 /// To set the address of the minimization function
731 
732 void TBackCompFitter::SetFCN(void (*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
733 {
734  fFCN = fcn;
735  if (fObjFunc) delete fObjFunc;
736  fObjFunc = new ROOT::Fit::FcnAdapter(fFCN);
737  DoSetDimension();
738 }
739 
740 
741 ////////////////////////////////////////////////////////////////////////////////
742 /// Set the objective function for fitting
743 /// Needed if fitting directly using TBackCompFitter class
744 /// The class clones a copy of the function and manages it
745 
746 void TBackCompFitter::SetObjFunction(ROOT::Math::IMultiGenFunction * fcn) {
747  if (fObjFunc) delete fObjFunc;
748  fObjFunc = fcn->Clone();
749 }
750 
751 ////////////////////////////////////////////////////////////////////////////////
752 /// Private method to set dimension in objective function
753 
754 void TBackCompFitter::DoSetDimension() {
755  if (!fObjFunc) return;
756  ROOT::Fit::FcnAdapter * fobj = dynamic_cast<ROOT::Fit::FcnAdapter*>(fObjFunc);
757  assert(fobj != 0);
758  int ndim = fFitter->Config().ParamsSettings().size();
759  if (ndim != 0) fobj->SetDimension(ndim);
760 }
761 
762 ////////////////////////////////////////////////////////////////////////////////
763 /// Return a pointer to the objective function (FCN)
764 /// If fitting directly using TBackCompFitter the pointer is managed by the class,
765 /// which has been set previously when calling SetObjFunction or SetFCN
766 /// Otherwise if the class is used in the backward compatible mode (e.g. after having fitted a TH1)
767 /// the return pointer will be valid after fitting and as long a new fit will not be done.
768 
769 ROOT::Math::IMultiGenFunction * TBackCompFitter::GetObjFunction( ) const {
770  if (fObjFunc) return fObjFunc;
771  return fFitter->GetFCN();
772 }
773 
774 ////////////////////////////////////////////////////////////////////////////////
775 /// Return a pointer to the minimizer.
776 /// the return pointer will be valid after fitting and as long a new fit will not be done.
777 /// For keeping a minimizer pointer the method ReCreateMinimizer() could eventually be used
778 
779 ROOT::Math::Minimizer * TBackCompFitter::GetMinimizer( ) const {
780  if (fMinimizer) return fMinimizer;
781  return fFitter->GetMinimizer();
782 }
783 
784 ////////////////////////////////////////////////////////////////////////////////
785 /// Return a new copy of the TFitResult object which needs to be deleted later by the user
786 
787 TFitResult * TBackCompFitter::GetTFitResult( ) const {
788  if (!fFitter.get() ) return 0;
789  return new TFitResult( fFitter->Result() );
790 }
791 
792 ////////////////////////////////////////////////////////////////////////////////
793 /// Scan parameter ipar between value of xmin and xmax
794 /// A graph must be given which will be on return filled with the scan resul
795 /// If the graph size is zero, a default size n = 40 will be used
796 
797 bool TBackCompFitter::Scan(unsigned int ipar, TGraph * gr, double xmin, double xmax )
798 {
799 
800  if (!gr) return false;
801  ROOT::Math::Minimizer * minimizer = fFitter->GetMinimizer();
802  if (!minimizer) {
803  Error("Scan","Minimizer is not available - cannot scan before fitting");
804  return false;
805  }
806 
807  unsigned int npoints = gr->GetN();
808  if (npoints == 0) {
809  npoints = 40;
810  gr->Set(npoints);
811  }
812  bool ret = minimizer->Scan( ipar, npoints, gr->GetX(), gr->GetY(), xmin, xmax);
813  if ((int) npoints < gr->GetN() ) gr->Set(npoints);
814  return ret;
815 }
816 
817 //______________________________________________________________________________
818 // bool TBackCompFitter::Scan2D(unsigned int ipar, unsigned int jpar, TGraph2D * gr,
819 // double xmin = 0, double xmax = 0, double ymin = 0, double ymax = 0) {
820 // // scan the parameters ipar between values of [xmin,xmax] and
821 // // jpar between values of [ymin,ymax] and
822 // // a graph2D must be given which will be on return filled with the scan resul
823 // // If the graph size is zero, a default size n = 20x20 will be used
824 // //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
825 
826 // if (!gr) return false;
827 // if (!fMinimizer) {
828 // Error("Scan","Minimizer is not available - cannot scan before fitting");
829 // return false;
830 // }
831 // unsigned int npoints = gr->GetN();
832 // if (npoints == 0) {
833 // npoints = 40;
834 // gr->Set(npoints);
835 // }
836 // // to be implemented
837 // for (unsigned int ix = 0; ix < npoints; ++ix) {
838 // return fMinimizer->Scan( ipar, npoints, gr->GetX(), gr->GetY(), xmin, xmax);
839 
840 // }
841 
842 ////////////////////////////////////////////////////////////////////////////////
843 /// Create a 2D contour around the minimum for the parameter ipar and jpar
844 /// if a minimum does not exist or is invalid it will return false
845 /// on exit a TGraph is filled with the contour points
846 /// the number of contour points is determined by the size of the TGraph.
847 /// if the size is zero a default number of points = 20 is used
848 /// pass optionally the confidence level, default is 0.683
849 /// it is assumed that ErrorDef() defines the right error definition
850 /// (i.e 1 sigma error for one parameter). If not the confidence level are scaled to new level
851 
852 bool TBackCompFitter::Contour(unsigned int ipar, unsigned int jpar, TGraph * gr, double confLevel) {
853  if (!gr) return false;
854  ROOT::Math::Minimizer * minimizer = fFitter->GetMinimizer();
855  if (!minimizer) {
856  Error("Scan","Minimizer is not available - cannot scan before fitting");
857  return false;
858  }
859 
860  // get error level used for fitting
861  double upScale = fFitter->Config().MinimizerOptions().ErrorDef();
862 
863  double upVal = TMath::ChisquareQuantile( confLevel, 2); // 2 is number of parameter we do the contour
864 
865  // set required error definition in minimizer
866  minimizer->SetErrorDef (upScale * upVal);
867 
868  unsigned int npoints = gr->GetN();
869  if (npoints == 0) {
870  npoints = 40;
871  gr->Set(npoints);
872  }
873  bool ret = minimizer->Contour( ipar, jpar, npoints, gr->GetX(), gr->GetY());
874  if ((int) npoints < gr->GetN() ) gr->Set(npoints);
875 
876  // restore the error level used for fitting
877  minimizer->SetErrorDef ( upScale);
878 
879  return ret;
880 }
881 
882