Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
FitResult.cxx
Go to the documentation of this file.
1 // @(#)root/mathcore:$Id$
2 // Author: L. Moneta Wed Aug 30 11:05:34 2006
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
7  * *
8  * *
9  **********************************************************************/
10 
11 // Implementation file for class FitResult
12 
13 #include "Fit/FitResult.h"
14 
15 #include "Fit/FitConfig.h"
16 
17 #include "Fit/BinData.h"
18 
19 //#include "Fit/Chi2FCN.h"
20 
21 #include "Math/Minimizer.h"
22 
23 #include "Math/IParamFunction.h"
25 
26 #include "Math/ProbFuncMathCore.h"
27 #include "Math/QuantFuncMathCore.h"
28 
29 #include "TMath.h"
31 #include "Math/Error.h"
32 
33 #include <cassert>
34 #include <cmath>
35 #include <iostream>
36 #include <iomanip>
37 
38 namespace ROOT {
39 
40  namespace Fit {
41 
42 
43 const int gInitialResultStatus = -99; // use this special convention to flag it when printing result
44 
45 FitResult::FitResult() :
46  fValid(false), fNormalized(false), fNFree(0), fNdf(0), fNCalls(0),
47  fStatus(-1), fCovStatus(0), fVal(0), fEdm(-1), fChi2(-1)
48 {
49  // Default constructor implementation.
50 }
51 
52 FitResult::FitResult(const FitConfig & fconfig) :
53  fValid(false),
54  fNormalized(false),
55  fNFree(0),
56  fNdf(0),
57  fNCalls(0),
58  fStatus(gInitialResultStatus),
59  fCovStatus(0),
60  fVal(0),
61  fEdm(-1),
62  fChi2(-1),
63  fFitFunc(0),
64  fParams(std::vector<double>( fconfig.NPar() ) ),
65  fErrors(std::vector<double>( fconfig.NPar() ) ),
66  fParNames(std::vector<std::string> ( fconfig.NPar() ) )
67 {
68  // create a Fit result from a fit config (i.e. with initial parameter values
69  // and errors equal to step values
70  // The model function is NULL in this case
71 
72  // set minimizer type and algorithm
73  fMinimType = fconfig.MinimizerType();
74  // append algorithm name for minimizer that support it
75  if ( (fMinimType.find("Fumili") == std::string::npos) &&
76  (fMinimType.find("GSLMultiFit") == std::string::npos)
77  ) {
78  if (fconfig.MinimizerAlgoType() != "") fMinimType += " / " + fconfig.MinimizerAlgoType();
79  }
80 
81  // get parameter values and errors (step sizes)
82  unsigned int npar = fconfig.NPar();
83  for (unsigned int i = 0; i < npar; ++i ) {
84  const ParameterSettings & par = fconfig.ParSettings(i);
85  fParams[i] = par.Value();
86  fErrors[i] = par.StepSize();
87  fParNames[i] = par.Name();
88  if (par.IsFixed() ) fFixedParams[i] = true;
89  else fNFree++;
90  if (par.IsBound() ) {
91  double lower = (par.HasLowerLimit()) ? par.LowerLimit() : - std::numeric_limits<double>::infinity() ;
92  double upper = (par.HasUpperLimit()) ? par.UpperLimit() : std::numeric_limits<double>::infinity() ;
93  fBoundParams[i] = fParamBounds.size();
94  fParamBounds.push_back(std::make_pair(lower,upper));
95  }
96  }
97  std::cout << "create fit result from config - nfree " << fNFree << std::endl;
98 }
99 
100 void FitResult::FillResult(const std::shared_ptr<ROOT::Math::Minimizer> & min, const FitConfig & fconfig, const std::shared_ptr<IModelFunction> & func,
101  bool isValid, unsigned int sizeOfData, bool binnedFit, const ROOT::Math::IMultiGenFunction * chi2func, unsigned int ncalls )
102 {
103  // Fill the FitResult after minimization using result from Minimizers
104 
105  // minimizer must exist
106  assert(min);
107 
108  fValid = isValid;
109  fNFree= min->NFree();
110  fNCalls = min->NCalls();
111  fStatus = min->Status();
112  fCovStatus= min->CovMatrixStatus();
113  fVal = min->MinValue();
114  fEdm = min->Edm();
115 
116  fMinimizer= min;
117  fFitFunc = func;
118 
119  fMinimType = fconfig.MinimizerName();
120 
121  // replace ncalls if minimizer does not support it (they are taken then from the FitMethodFunction)
122  if (fNCalls == 0) fNCalls = ncalls;
123 
124  const unsigned int npar = min->NDim();
125  if (npar == 0) return;
126 
127  if (min->X() )
128  fParams = std::vector<double>(min->X(), min->X() + npar);
129  else {
130  // case minimizer does not provide minimum values (it failed) take from configuration
131  fParams.resize(npar);
132  for (unsigned int i = 0; i < npar; ++i ) {
133  fParams[i] = ( fconfig.ParSettings(i).Value() );
134  }
135  }
136 
137  if (sizeOfData > min->NFree() ) fNdf = sizeOfData - min->NFree();
138 
139 
140  // set right parameters in function (in case minimizer did not do before)
141  // do also when fit is not valid
142  if (func ) {
143  // I think we can avoid cloning the model function
144  //fFitFunc = dynamic_cast<IModelFunction *>( func->Clone() );
145  //assert(fFitFunc);
146  fFitFunc->SetParameters(&fParams.front());
147  }
148  else {
149  // when no fFitFunc is present take parameters from FitConfig
150  fParNames.reserve( npar );
151  for (unsigned int i = 0; i < npar; ++i ) {
152  fParNames.push_back( fconfig.ParSettings(i).Name() );
153  }
154  }
155 
156 
157  // check for fixed or limited parameters
158  unsigned int nfree = 0;
159  for (unsigned int ipar = 0; ipar < npar; ++ipar) {
160  const ParameterSettings & par = fconfig.ParSettings(ipar);
161  if (par.IsFixed() ) fFixedParams[ipar] = true;
162  else nfree++;
163  if (par.IsBound() ) {
164  double lower = (par.HasLowerLimit()) ? par.LowerLimit() : - std::numeric_limits<double>::infinity() ;
165  double upper = (par.HasUpperLimit()) ? par.UpperLimit() : std::numeric_limits<double>::infinity() ;
166  fBoundParams[ipar] = fParamBounds.size();
167  fParamBounds.push_back(std::make_pair(lower,upper));
168  }
169  }
170  // check if nfree (from FitConfig) and fNFree (from minimizer) are consistent
171  if (nfree != fNFree ) {
172  MATH_ERROR_MSG("FitResult","FitConfiguration and Minimizer result are not consistent");
173  std::cout << "Number of free parameters from FitConfig = " << nfree << std::endl;
174  std::cout << "Number of free parameters from Minimizer = " << fNFree << std::endl;
175  }
176 
177  // if flag is binned compute a chi2 when a chi2 function is given
178  if (binnedFit) {
179  if (chi2func == 0)
180  fChi2 = fVal;
181  else {
182  // compute chi2 equivalent for likelihood fits
183  // NB: empty bins are considered
184  fChi2 = (*chi2func)(&fParams[0]);
185  }
186  }
187 
188  // fill error matrix
189  // if minimizer provides error provides also error matrix
190  if (min->Errors() != 0) {
191 
192  fErrors = std::vector<double>(min->Errors(), min->Errors() + npar ) ;
193 
194  if (fCovStatus != 0) {
195  unsigned int r = npar * ( npar + 1 )/2;
196  fCovMatrix.reserve(r);
197  for (unsigned int i = 0; i < npar; ++i)
198  for (unsigned int j = 0; j <= i; ++j)
199  fCovMatrix.push_back(min->CovMatrix(i,j) );
200  }
201 
202  // minos errors
203  if (fValid && fconfig.MinosErrors()) {
204  const std::vector<unsigned int> & ipars = fconfig.MinosParams();
205  unsigned int n = (ipars.size() > 0) ? ipars.size() : npar;
206  for (unsigned int i = 0; i < n; ++i) {
207  double elow, eup;
208  unsigned int index = (ipars.size() > 0) ? ipars[i] : i;
209  bool ret = min->GetMinosError(index, elow, eup);
210  if (ret) SetMinosError(index, elow, eup);
211  }
212  }
213 
214  // globalCC
215  fGlobalCC.reserve(npar);
216  for (unsigned int i = 0; i < npar; ++i) {
217  double globcc = min->GlobalCC(i);
218  if (globcc < 0) break; // it is not supported by that minimizer
219  fGlobalCC.push_back(globcc);
220  }
221 
222  }
223 
224 }
225 
226 FitResult::~FitResult() {
227  // destructor. FitResult manages the fit Function pointer
228  //if (fFitFunc) delete fFitFunc;
229 }
230 
231 FitResult::FitResult(const FitResult &rhs) :
232  fFitFunc(0)
233 {
234  // Implementation of copy constructor
235  (*this) = rhs;
236 }
237 
238 FitResult & FitResult::operator = (const FitResult &rhs) {
239  // Implementation of assignment operator.
240  if (this == &rhs) return *this; // time saving self-test
241 
242  // Manages the fitted function
243  // if (fFitFunc) delete fFitFunc;
244  // fFitFunc = 0;
245  // if (rhs.fFitFunc != 0 ) {
246  // fFitFunc = dynamic_cast<IModelFunction *>( (rhs.fFitFunc)->Clone() );
247  // assert(fFitFunc != 0);
248  // }
249 
250  // copy all other data members
251  fValid = rhs.fValid;
252  fNormalized = rhs.fNormalized;
253  fNFree = rhs.fNFree;
254  fNdf = rhs.fNdf;
255  fNCalls = rhs.fNCalls;
256  fCovStatus = rhs.fCovStatus;
257  fStatus = rhs.fStatus;
258  fVal = rhs.fVal;
259  fEdm = rhs.fEdm;
260  fChi2 = rhs.fChi2;
261  fMinimizer = rhs.fMinimizer;
262  fObjFunc = rhs.fObjFunc;
263  fFitFunc = rhs.fFitFunc;
264 
265  fFixedParams = rhs.fFixedParams;
266  fBoundParams = rhs.fBoundParams;
267  fParamBounds = rhs.fParamBounds;
268  fParams = rhs.fParams;
269  fErrors = rhs.fErrors;
270  fCovMatrix = rhs.fCovMatrix;
271  fGlobalCC = rhs.fGlobalCC;
272  fMinosErrors = rhs.fMinosErrors;
273 
274  fMinimType = rhs.fMinimType;
275  fParNames = rhs.fParNames;
276 
277  return *this;
278 
279 }
280 
281 bool FitResult::Update(const std::shared_ptr<ROOT::Math::Minimizer> & min, const ROOT::Fit::FitConfig & fconfig, bool isValid, unsigned int ncalls) {
282  // update fit result with new status from minimizer
283  // ncalls if it is not zero is used instead of value from minimizer
284 
285  fMinimizer = min;
286 
287  // in case minimizer changes
288  fMinimType = fconfig.MinimizerName();
289 
290  const unsigned int npar = fParams.size();
291  if (min->NDim() != npar ) {
292  MATH_ERROR_MSG("FitResult::Update","Wrong minimizer status ");
293  return false;
294  }
295  if (min->X() == 0 ) {
296  MATH_ERROR_MSG("FitResult::Update","Invalid minimizer status ");
297  return false;
298  }
299  //fNFree = min->NFree();
300  if (fNFree != min->NFree() ) {
301  MATH_ERROR_MSG("FitResult::Update","Configuration has changed ");
302  return false;
303  }
304 
305  fValid = isValid;
306  // update minimum value
307  fVal = min->MinValue();
308  fEdm = min->Edm();
309  fStatus = min->Status();
310  fCovStatus = min->CovMatrixStatus();
311 
312  // update number of function calls
313  if ( min->NCalls() > 0) fNCalls = min->NCalls();
314  else fNCalls = ncalls;
315 
316  // copy parameter value and errors
317  std::copy(min->X(), min->X() + npar, fParams.begin());
318 
319 
320  // set parameters in fit model function
321  if (fFitFunc) fFitFunc->SetParameters(&fParams.front());
322 
323  if (min->Errors() != 0) {
324 
325  if (fErrors.size() != npar) fErrors.resize(npar);
326 
327  std::copy(min->Errors(), min->Errors() + npar, fErrors.begin() ) ;
328 
329  if (fCovStatus != 0) {
330 
331  // update error matrix
332  unsigned int r = npar * ( npar + 1 )/2;
333  if (fCovMatrix.size() != r) fCovMatrix.resize(r);
334  unsigned int l = 0;
335  for (unsigned int i = 0; i < npar; ++i) {
336  for (unsigned int j = 0; j <= i; ++j)
337  fCovMatrix[l++] = min->CovMatrix(i,j);
338  }
339  }
340 
341  // update global CC
342  if (fGlobalCC.size() != npar) fGlobalCC.resize(npar);
343  for (unsigned int i = 0; i < npar; ++i) {
344  double globcc = min->GlobalCC(i);
345  if (globcc < 0) {
346  fGlobalCC.clear();
347  break; // it is not supported by that minimizer
348  }
349  fGlobalCC[i] = globcc;
350  }
351 
352  }
353  return true;
354 }
355 
356 void FitResult::NormalizeErrors() {
357  // normalize errors and covariance matrix according to chi2 value
358  if (fNdf == 0 || fChi2 <= 0) return;
359  double s2 = fChi2/fNdf;
360  double s = std::sqrt(fChi2/fNdf);
361  for (unsigned int i = 0; i < fErrors.size() ; ++i)
362  fErrors[i] *= s;
363  for (unsigned int i = 0; i < fCovMatrix.size() ; ++i)
364  fCovMatrix[i] *= s2;
365 
366  fNormalized = true;
367 }
368 
369 
370 double FitResult::Prob() const {
371  // fit probability
372  return ROOT::Math::chisquared_cdf_c(fChi2, static_cast<double>(fNdf) );
373 }
374 
375 bool FitResult::HasMinosError(unsigned int i) const {
376  // query if the parameter i has the Minos error
377  std::map<unsigned int, std::pair<double,double> >::const_iterator itr = fMinosErrors.find(i);
378  return (itr != fMinosErrors.end() );
379 }
380 
381 
382 double FitResult::LowerError(unsigned int i) const {
383  // return lower Minos error for parameter i
384  // return the parabolic error if Minos error has not been calculated for the parameter i
385  std::map<unsigned int, std::pair<double,double> >::const_iterator itr = fMinosErrors.find(i);
386  return ( itr != fMinosErrors.end() ) ? itr->second.first : Error(i) ;
387 }
388 
389 double FitResult::UpperError(unsigned int i) const {
390  // return upper Minos error for parameter i
391  // return the parabolic error if Minos error has not been calculated for the parameter i
392  std::map<unsigned int, std::pair<double,double> >::const_iterator itr = fMinosErrors.find(i);
393  return ( itr != fMinosErrors.end() ) ? itr->second.second : Error(i) ;
394 }
395 
396 void FitResult::SetMinosError(unsigned int i, double elow, double eup) {
397  // set the Minos error for parameter i
398  fMinosErrors[i] = std::make_pair(elow,eup);
399 }
400 
401 int FitResult::Index(const std::string & name) const {
402  // find index for given parameter name
403  if (! fFitFunc) return -1;
404  unsigned int npar = fParams.size();
405  for (unsigned int i = 0; i < npar; ++i)
406  if ( fFitFunc->ParameterName(i) == name) return i;
407 
408  return -1; // case name is not found
409 }
410 
411 bool FitResult::IsParameterBound(unsigned int ipar) const {
412  return fBoundParams.find(ipar) != fBoundParams.end();
413 }
414 
415 bool FitResult::IsParameterFixed(unsigned int ipar) const {
416  return fFixedParams.find(ipar) != fFixedParams.end();
417 }
418 
419 bool FitResult::ParameterBounds(unsigned int ipar, double & lower, double & upper) const {
420  std::map<unsigned int, unsigned int>::const_iterator itr = fBoundParams.find(ipar);
421  if (itr == fBoundParams.end() ) {
422  lower = -std::numeric_limits<Double_t>::infinity();
423  upper = std::numeric_limits<Double_t>::infinity();
424  return false;
425  }
426  assert(itr->second < fParamBounds.size() );
427  lower = fParamBounds[itr->second].first;
428  upper = fParamBounds[itr->second].second;
429  return false;
430 }
431 
432 std::string FitResult::ParName(unsigned int ipar) const {
433  // return parameter name
434  if (fFitFunc) return fFitFunc->ParameterName(ipar);
435  else if (ipar < fParNames.size() ) return fParNames[ipar];
436  return "param_" + ROOT::Math::Util::ToString(ipar);
437 }
438 
439 void FitResult::Print(std::ostream & os, bool doCovMatrix) const {
440  // print the result in the given stream
441  // need to add also minos errors , globalCC, etc..
442  unsigned int npar = fParams.size();
443  if (npar == 0) {
444  os << "<Empty FitResult>\n";
445  return;
446  }
447  os << "\n****************************************\n";
448  if (!fValid) {
449  if (fStatus != gInitialResultStatus) {
450  os << " Invalid FitResult";
451  os << " (status = " << fStatus << " )";
452  }
453  else {
454  os << " FitResult before fitting";
455  }
456  os << "\n****************************************\n";
457  }
458 
459  //os << " FitResult \n\n";
460  os << "Minimizer is " << fMinimType << std::endl;
461  const unsigned int nw = 25; // spacing for text
462  const unsigned int nn = 12; // spacing for numbers
463  const std::ios_base::fmtflags prFmt = os.setf(std::ios::left,std::ios::adjustfield); // set left alignment
464 
465  if (fVal != fChi2 || fChi2 < 0)
466  os << std::left << std::setw(nw) << "MinFCN" << " = " << std::right << std::setw(nn) << fVal << std::endl;
467  if (fChi2 >= 0)
468  os << std::left << std::setw(nw) << "Chi2" << " = " << std::right << std::setw(nn) << fChi2 << std::endl;
469  os << std::left << std::setw(nw) << "NDf" << " = " << std::right << std::setw(nn) << fNdf << std::endl;
470  if (fMinimType.find("Linear") == std::string::npos) { // no need to print this for linear fits
471  if (fEdm >=0) os << std::left << std::setw(nw) << "Edm" << " = " << std::right << std::setw(nn) << fEdm << std::endl;
472  os << std::left << std::setw(nw) << "NCalls" << " = " << std::right << std::setw(nn) << fNCalls << std::endl;
473  }
474  for (unsigned int i = 0; i < npar; ++i) {
475  os << std::left << std::setw(nw) << GetParameterName(i);
476  os << " = " << std::right << std::setw(nn) << fParams[i];
477  if (IsParameterFixed(i) )
478  os << std::setw(9) << " " << std::setw(nn) << " " << " \t (fixed)";
479  else {
480  if (fErrors.size() != 0)
481  os << " +/- " << std::left << std::setw(nn) << fErrors[i] << std::right;
482  if (IsParameterBound(i) )
483  os << " \t (limited)";
484  }
485  os << std::endl;
486  }
487 
488  // restore stremam adjustfield
489  if (prFmt != os.flags() ) os.setf(prFmt, std::ios::adjustfield);
490 
491  if (doCovMatrix) PrintCovMatrix(os);
492 }
493 
494 void FitResult::PrintCovMatrix(std::ostream &os) const {
495  // print the covariance and correlation matrix
496  if (!fValid) return;
497  if (fCovMatrix.size() == 0) return;
498 // os << "****************************************\n";
499  os << "\nCovariance Matrix:\n\n";
500  unsigned int npar = fParams.size();
501  const int kPrec = 5;
502  const int kWidth = 8;
503  const int parw = 12;
504  const int matw = kWidth+4;
505 
506  // query previous precision and format flags
507  int prevPrec = os.precision(kPrec);
508  const std::ios_base::fmtflags prevFmt = os.flags();
509 
510  os << std::setw(parw) << " " << "\t";
511  for (unsigned int i = 0; i < npar; ++i) {
512  if (!IsParameterFixed(i) ) {
513  os << std::right << std::setw(matw) << GetParameterName(i) ;
514  }
515  }
516  os << std::endl;
517  for (unsigned int i = 0; i < npar; ++i) {
518  if (!IsParameterFixed(i) ) {
519  os << std::left << std::setw(parw) << GetParameterName(i) << "\t";
520  for (unsigned int j = 0; j < npar; ++j) {
521  if (!IsParameterFixed(j) ) {
522  os.precision(kPrec); os.width(kWidth); os << std::right << std::setw(matw) << CovMatrix(i,j);
523  }
524  }
525  os << std::endl;
526  }
527  }
528 // os << "****************************************\n";
529  os << "\nCorrelation Matrix:\n\n";
530  os << std::setw(parw) << " " << "\t";
531  for (unsigned int i = 0; i < npar; ++i) {
532  if (!IsParameterFixed(i) ) {
533  os << std::right << std::setw(matw) << GetParameterName(i) ;
534  }
535  }
536  os << std::endl;
537  for (unsigned int i = 0; i < npar; ++i) {
538  if (!IsParameterFixed(i) ) {
539  os << std::left << std::setw(parw) << std::left << GetParameterName(i) << "\t";
540  for (unsigned int j = 0; j < npar; ++j) {
541  if (!IsParameterFixed(j) ) {
542  os.precision(kPrec); os.width(kWidth); os << std::right << std::setw(matw) << Correlation(i,j);
543  }
544  }
545  os << std::endl;
546  }
547  }
548  // restore alignment and precision
549  os.setf(prevFmt, std::ios::adjustfield);
550  os.precision(prevPrec);
551 }
552 
553 void FitResult::GetConfidenceIntervals(unsigned int n, unsigned int stride1, unsigned int stride2, const double * x, double * ci, double cl, bool norm ) const {
554  // stride1 stride in coordinate stride2 stride in dimension space
555  // i.e. i-th point in k-dimension is x[ stride1 * i + stride2 * k]
556  // compute the confidence interval of the fit on the given data points
557  // the dimension of the data points must match the dimension of the fit function
558  // confidence intervals are returned in array ci
559 
560  if (!fFitFunc) {
561  // check if model function exists
562  MATH_ERROR_MSG("FitResult::GetConfidenceIntervals","Cannot compute Confidence Intervals without fit model function");
563  return;
564  }
565  assert(fFitFunc);
566 
567  // use student quantile in case of normalized errors
568  double corrFactor = 1;
569  if (fChi2 <= 0 || fNdf == 0) norm = false;
570  if (norm)
571  corrFactor = TMath::StudentQuantile(0.5 + cl/2, fNdf) * std::sqrt( fChi2/fNdf );
572  else
573  // value to go up in chi2 (1: 1 sigma error(CL=0.683) , 4: 2 sigma errors
574  corrFactor = ROOT::Math::chisquared_quantile(cl, 1);
575 
576 
577 
578  unsigned int ndim = fFitFunc->NDim();
579  unsigned int npar = fFitFunc->NPar();
580 
581  std::vector<double> xpoint(ndim);
582  std::vector<double> grad(npar);
583  std::vector<double> vsum(npar);
584 
585  // loop on the points
586  for (unsigned int ipoint = 0; ipoint < n; ++ipoint) {
587 
588  for (unsigned int kdim = 0; kdim < ndim; ++kdim) {
589  unsigned int i = ipoint * stride1 + kdim * stride2;
590  assert(i < ndim*n);
591  xpoint[kdim] = x[i];
592  }
593 
594  // calculate gradient of fitted function w.r.t the parameters
595  ROOT::Math::RichardsonDerivator d;
596  for (unsigned int ipar = 0; ipar < npar; ++ipar) {
597  if (!IsParameterFixed(ipar)) {
598  ROOT::Math::OneDimParamFunctionAdapter<const ROOT::Math::IParamMultiFunction &> fadapter(*fFitFunc,&xpoint.front(),&fParams.front(),ipar);
599  d.SetFunction(fadapter);
600  // compute step size as a small fraction of the error
601  // (see numerical recipes in C 5.7.8) 1.E-5 is ~ (eps)^1/3
602  if ( fErrors[ipar] > 0 )
603  d.SetStepSize( std::max( fErrors[ipar]*1.E-5, 1.E-15) );
604  else
605  d.SetStepSize( std::min(std::max(fParams[ipar]*1.E-5, 1.E-15), 0.0001 ) );
606 
607  grad[ipar] = d(fParams[ipar] ); // evaluate df/dp
608  }
609  else
610  grad[ipar] = 0.; // for fixed parameters
611  }
612 
613  // multiply covariance matrix with gradient
614  vsum.assign(npar,0.0);
615  for (unsigned int ipar = 0; ipar < npar; ++ipar) {
616  for (unsigned int jpar = 0; jpar < npar; ++jpar) {
617  vsum[ipar] += CovMatrix(ipar,jpar) * grad[jpar];
618  }
619  }
620  // multiply gradient by vsum
621  double r2 = 0;
622  for (unsigned int ipar = 0; ipar < npar; ++ipar) {
623  r2 += grad[ipar] * vsum[ipar];
624  }
625  double r = std::sqrt(r2);
626  ci[ipoint] = r * corrFactor;
627  }
628 }
629 
630 void FitResult::GetConfidenceIntervals(const BinData & data, double * ci, double cl, bool norm ) const {
631  // implement confidence intervals from a given bin data sets
632  // currently copy the data from Bindata.
633  // could implement otherwise directly
634  unsigned int ndim = data.NDim();
635  unsigned int np = data.NPoints();
636  std::vector<double> xdata( ndim * np );
637  for (unsigned int i = 0; i < np ; ++i) {
638  const double * x = data.Coords(i);
639  std::vector<double>::iterator itr = xdata.begin()+ ndim * i;
640  std::copy(x,x+ndim,itr);
641  }
642  // points are arraned as x0,y0,z0, ....xN,yN,zN (stride1=ndim, stride2=1)
643  GetConfidenceIntervals(np,ndim,1,&xdata.front(),ci,cl,norm);
644 }
645 
646 std::vector<double> FitResult::GetConfidenceIntervals(double cl, bool norm ) const {
647  // implement confidence intervals using stored data sets (if can be retrieved from objective function)
648  // it works only in case of chi2 or binned likelihood fits
649  const BinData * data = FittedBinData();
650  std::vector<double> result;
651  if (data) {
652  result.resize(data->NPoints() );
653  GetConfidenceIntervals(*data, result.data(), cl, norm);
654  }
655  else {
656  MATH_ERROR_MSG("FitResult::GetConfidenceIntervals","Cannot compute Confidence Intervals without the fit bin data");
657  }
658  return result;
659 }
660 
661 // const BinData * GetFitBinData() const {
662 // // return a pointer to the binned data used in the fit
663 // // works only for chi2 or binned likelihood fits
664 // // thus when the objective function stored is a Chi2Func or a PoissonLikelihood
665 // ROOT::Math::IMultiGenFunction * f = fObjFunc->get();
666 // Chi2Function * chi2func = dynamic_cast<Chi2Function*>(f);
667 // if (chi2func) return &(chi2func->Data());
668 // PoissonLLFunction * pllfunc = dynamic_cast<PoissonLLFunction*>(f);
669 // if (pllfunc) return &(pllfunc->Data());
670 // Chi2GradFunction * chi2gradfunc = dynamic_cast<Chi2GradFunction*>(f);
671 // if (chi2gradfunc) return &(chi2gradfunc->Data());
672 // PoissonLLGradFunction * pllgradfunc = dynamic_cast<PoissonLLFunction*>(f);
673 // if (pllgradfunc) return &(pllgradfunc->Data());
674 // MATH_WARN_MSG("FitResult::GetFitBinData","Cannot retrun fit bin data set if objective function is not of a known type");
675 // return nullptr;
676 // }
677 
678 const BinData * FitResult::FittedBinData() const {
679  return dynamic_cast<const BinData*> ( fFitData.get() );
680 }
681 
682 ////////////////////////////////////////////////////////////////////////////////
683 /// Scan parameter ipar between value of xmin and xmax
684 /// A array for x and y points should be provided
685 
686 bool FitResult::Scan(unsigned int ipar, unsigned int &npoints, double *pntsx, double *pntsy, double xmin, double xmax)
687 {
688  if (!pntsx || !pntsy || !npoints)
689  return false;
690 
691  if (!fMinimizer) {
692  MATH_ERROR_MSG("FitResult::Scan", "Minimizer is not available - cannot Scan");
693  return false;
694  }
695 
696  return fMinimizer->Scan(ipar, npoints, pntsx, pntsy, xmin, xmax);
697 }
698 
699 ////////////////////////////////////////////////////////////////////////////////
700 /// Create a 2D contour around the minimum for the parameter ipar and jpar
701 /// if a minimum does not exist or is invalid it will return false
702 /// A array for x and y points should be provided
703 /// Pass optionally the confidence level, default is 0.683
704 /// it is assumed that ErrorDef() defines the right error definition
705 /// (i.e 1 sigma error for one parameter). If not the confidence level are scaled to new level
706 
707 bool FitResult::Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double *pntsx, double *pntsy, double confLevel)
708 {
709  if (!pntsx || !pntsy || !npoints)
710  return false;
711 
712  if (!fMinimizer) {
713  MATH_ERROR_MSG("FitResult::Contour", "Minimizer is not available - cannot produce Contour");
714  return false;
715  }
716 
717  // get error level used for fitting
718  double upScale = fMinimizer->ErrorDef();
719 
720  double upVal = TMath::ChisquareQuantile(confLevel, 2); // 2 is number of parameter we do the contour
721 
722  // set required error definition in minimizer
723  fMinimizer->SetErrorDef(upScale * upVal);
724 
725  bool ret = fMinimizer->Contour(ipar, jpar, npoints, pntsx, pntsy);
726 
727  // restore the error level used for fitting
728  fMinimizer->SetErrorDef(upScale);
729 
730  return ret;
731 }
732 
733  } // end namespace Fit
734 
735 } // end namespace ROOT