BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
FitResult.h
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 // Header file for class FitResult
12 
13 #ifndef ROOT_Fit_FitResult
14 #define ROOT_Fit_FitResult
15 
16 #include "Math/IFunctionfwd.h"
17 #include "Math/IParamFunctionfwd.h"
18 
19 #include <vector>
20 #include <map>
21 #include <string>
22 #include <cmath>
23 #include <cassert>
24 #include <memory>
25 
26 namespace ROOT {
27 
28  namespace Math {
29  class Minimizer;
30  }
31 
32 
33  namespace Fit {
34 
35  class FitConfig;
36  class FitData;
37  class BinData;
38 
39 //___________________________________________________________________________________
40 /**
41  class containg the result of the fit and all the related information
42  (fitted parameter values, error, covariance matrix and minimizer result information)
43  Contains a pointer also to the fitted (model) function, modified with the fit parameter values.
44  When the fit is valid, it is constructed from a Minimizer and a model function pointer
45 
46  @ingroup FitMain
47 */
48 class FitResult {
49 
50 public:
51 
53 
54  /**
55  Default constructor for an empty (non valid) fit result
56  */
58 
59  /**
60  Constructor from a fit-config for a dummy fit
61  (e.g. when only one fcn evaluation is done)
62  */
63  FitResult (const FitConfig & fconfig);
64 
65 
66  /**
67  Copy constructor.
68  */
69  FitResult(const FitResult & rhs);
70 
71  /**
72  Assignment operator
73  */
75 
76  /**
77  Destructor
78  */
79  virtual ~FitResult ();
80 
81 
82 public:
83 
84  /**
85  Fill the fit result from a Minimizer instance after fitting
86  Run also Minos if requested from the configuration
87  */
88  void FillResult(const std::shared_ptr<ROOT::Math::Minimizer> & min, const FitConfig & fconfig, const std::shared_ptr<IModelFunction> & f,
89  bool isValid, unsigned int sizeOfData = 0, bool binFit = true, const ROOT::Math::IMultiGenFunction * chi2func = 0, unsigned int ncalls = 0);
90 
91 
92  /**
93  Update the fit result with a new minimization status
94  To be run only if same fit is performed with same configuration
95  Note that in this case MINOS is not re-run. If one wants to run also MINOS
96  a new result must be created
97  */
98  bool Update(const std::shared_ptr<ROOT::Math::Minimizer> & min, bool isValid, unsigned int ncalls = 0 );
99 
100  /** minimization quantities **/
101 
102  /// minimizer type
103  const std::string & MinimizerType() const { return fMinimType; }
104 
105  /**
106  True if fit successful, otherwise false.
107  A fit is considered successful if the minimizer succeded in finding the
108  minimum. It could happen that subsequent operations like error analysis (e.g. Minos)
109  failed. In that case the status can be still true if the original minimization algorithm
110  succeeded in finding the minimum.
111  One can query in that case the minimizer return status using Status().
112  It is responability to the Minimizer class to tag a found minimum as valid or not
113  and to produce also a status code.
114  */
115  bool IsValid() const { return fValid; }
116 
117  /// True if a fit result does not exist (even invalid) with parameter values
118  bool IsEmpty() const { return (fParams.size() == 0); }
119 
120  /// Return value of the objective function (chi2 or likelihood) used in the fit
121  double MinFcnValue() const { return fVal; }
122 
123  ///Number of function calls to find minimum
124  unsigned int NCalls() const { return fNCalls; }
125 
126  ///Expected distance from minimum
127  double Edm() const { return fEdm; }
128 
129  /// get total number of parameters
130  unsigned int NTotalParameters() const { return fParams.size(); }
131  /// total number of parameters (abbreviation)
132  unsigned int NPar() const { return NTotalParameters(); }
133 
134  /// get total number of free parameters
135  unsigned int NFreeParameters() const { return fNFree; }
136 
137  /// minimizer status code
138  int Status() const { return fStatus; }
139 
140  ///covariance matrix status code
141  /// using Minuit convention : =0 not calculated, =1 approximated, =2 made pos def , =3 accurate
142 
143  int CovMatrixStatus() const { return fCovStatus; }
144 
145  /** fitting quantities **/
146 
147  /// Return pointer to model (fit) function with fitted parameter values.
148  /// Pointer is managed internally. I must not be deleted
149  const IModelFunction * FittedFunction() const {
150  return fFitFunc.get();
151  }
152 
153  /// return BinData used in the fit (return a nullptr in case a different fit is done
154  /// or the data are not available
155  /// Pointer is managed internally, it must not be deleted
156  const BinData * FittedBinData() const;
157 
158 
159  /// Chi2 fit value
160  /// in case of likelihood must be computed ?
161  double Chi2() const { return fChi2; }
162 
163  /// Number of degree of freedom
164  unsigned int Ndf() const { return fNdf; }
165 
166  /// p value of the fit (chi2 probability)
167  double Prob() const;
168 
169  /// parameter errors (return st::vector)
170  const std::vector<double> & Errors() const { return fErrors; }
171  /// parameter errors (return const pointer)
172  const double * GetErrors() const { return (fErrors.empty()) ? 0 : &fErrors.front(); }
173 
174  /// parameter values (return std::vector)
175  const std::vector<double> & Parameters() const { return fParams; }
176  /// parameter values (return const pointer)
177  const double * GetParams() const { return &fParams.front(); }
178 
179  /// parameter value by index
180  double Value(unsigned int i) const { return fParams[i]; }
181  /// parameter value by index
182  double Parameter(unsigned int i) const { return fParams[i]; }
183 
184  /// parameter error by index
185  // (NOTE: this due to conflict with TObject::Error cannot used in derived class which
186  // inherits from TObject. Use instead ParError (or Errors()[i] )
187  double Error(unsigned int i) const {
188  return (i < fErrors.size() ) ? fErrors[i] : 0;
189  }
190  /// parameter error by index
191  double ParError(unsigned int i) const {
192  return (i < fErrors.size() ) ? fErrors[i] : 0;
193  }
194 
195  /// name of the parameter
196  std::string ParName(unsigned int i) const;
197 
198  /// set the Minos errors for parameter i (called by the Fitter class when running Minos)
199  void SetMinosError(unsigned int i, double elow, double eup);
200 
201  /// query if parameter i has the Minos error
202  bool HasMinosError(unsigned int i) const;
203 
204  /// lower Minos error. If Minos has not run for parameter i return the parabolic error
205  double LowerError(unsigned int i) const;
206 
207  /// upper Minos error. If Minos has not run for parameter i return the parabolic error
208  double UpperError(unsigned int i) const;
209 
210  /// parameter global correlation coefficient
211  double GlobalCC(unsigned int i) const {
212  return (i < fGlobalCC.size() ) ? fGlobalCC[i] : -1;
213  }
214 
215 
216  /// retrieve covariance matrix element
217  double CovMatrix (unsigned int i, unsigned int j) const {
218  if ( i >= fErrors.size() || j >= fErrors.size() ) return 0;
219  if (fCovMatrix.size() == 0) return 0; // no matrix is available in case of non-valid fits
220  if ( j < i )
221  return fCovMatrix[j + i* (i+1) / 2];
222  else
223  return fCovMatrix[i + j* (j+1) / 2];
224  }
225 
226  /// retrieve correlation elements
227  double Correlation(unsigned int i, unsigned int j ) const {
228  if ( i >= fErrors.size() || j >= fErrors.size() ) return 0;
229  if (fCovMatrix.size() == 0) return 0; // no matrix is available in case of non-valid fits
230  double tmp = CovMatrix(i,i)*CovMatrix(j,j);
231  return ( tmp > 0) ? CovMatrix(i,j)/ std::sqrt(tmp) : 0;
232  }
233 
234  /// fill covariance matrix elements using a generic matrix class implementing operator(i,j)
235  /// the matrix must be previously allocates with right size (npar * npar)
236  template<class Matrix>
237  void GetCovarianceMatrix(Matrix & mat) const {
238  unsigned int npar = fErrors.size();
239  if (fCovMatrix.size() != npar*(npar+1)/2 ) return; // do nothing
240  for (unsigned int i = 0; i< npar; ++i) {
241  for (unsigned int j = 0; j<=i; ++j) {
242  mat(i,j) = fCovMatrix[j + i*(i+1)/2 ];
243  if (i != j) mat(j,i) = mat(i,j);
244  }
245  }
246  }
247 
248  /// fill a correlation matrix elements using a generic symmetric matrix class implementing operator(i,j)
249  /// the matrix must be previously allocates with right size (npar * npar)
250  template<class Matrix>
251  void GetCorrelationMatrix(Matrix & mat) const {
252  unsigned int npar = fErrors.size();
253  if (fCovMatrix.size() != npar*(npar+1)/2) return; // do nothing
254  for (unsigned int i = 0; i< npar; ++i) {
255  for (unsigned int j = 0; j<=i; ++j) {
256  double tmp = fCovMatrix[i * (i +3)/2 ] * fCovMatrix[ j * (j+3)/2 ];
257  mat(i,j) = (tmp > 0) ? fCovMatrix[j + i*(i+1)/2 ] / std::sqrt(tmp) : 0;
258  if (i != j) mat(j,i) = mat(i,j);
259  }
260  }
261  }
262 
263  /**
264  get confidence intervals for an array of n points x.
265  stride1 indicates the stride in the coordinate space while stride2 the stride in dimension space.
266  For 1-dim points : stride1=1, stride2=1
267  for multi-dim points arranged as (x0,x1,...,xN,y0,....yN) stride1=1 stride2=n
268  for multi-dim points arraged as (x0,y0,..,x1,y1,...,xN,yN,..) stride1=ndim, stride2=1
269 
270  the confidence interval are returned in the array ci
271  cl is the desired confidedence interval value
272  norm is a flag to control if the intervals need to be normalized to the chi2/ndf value
273  The intervals can be corrected optionally using the chi2/ndf value of the fit if a chi2 fit is performed.
274  This has changed since ROOT 6.14, before the interval were corrected by default.
275  */
276  void GetConfidenceIntervals(unsigned int n, unsigned int stride1, unsigned int stride2, const double * x, double * ci, double cl=0.95, bool norm = false ) const;
277 
278  /**
279  evaluate confidence interval for the point specified in the passed data sets
280  the confidence interval are returned in the array ci
281  cl is the desired confidence interval value.
282  This method is mantained for backward compatibility and will be deprecated
283  */
284  void GetConfidenceIntervals(const BinData & data, double * ci, double cl=0.95, bool norm = false ) const;
285 
286  /**
287  evaluate confidence interval for the data set used in the last fit
288  the confidence interval are returned as a vector of data points
289  */
290  std::vector<double> GetConfidenceIntervals(double cl=0.95, bool norm = false ) const;
291 
292  /**
293  scan likelihood value of parameter and fill the given graph.
294  */
295  bool Scan(unsigned int ipar, unsigned int &npoints, double *pntsx, double *pntsy, double xmin = 0, double xmax = 0);
296 
297  /**
298  create contour of two parameters around the minimum
299  pass as option confidence level: default is a value of 0.683
300  */
301  bool Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double *pntsx, double *pntsy, double confLevel = 0.683);
302 
303  /// get index for parameter name (return -1 if not found)
304  int Index(const std::string & name) const;
305 
306  ///normalize errors using chi2/ndf for chi2 fits
308 
309  /// flag to chek if errors are normalized
310  bool NormalizedErrors() const { return fNormalized; }
311 
312  /// print the result and optionaly covariance matrix and correlations
313  void Print(std::ostream & os, bool covmat = false) const;
314 
315  ///print error matrix and correlations
316  void PrintCovMatrix(std::ostream & os) const;
317 
318  /// query if a parameter is bound
319  bool IsParameterBound(unsigned int ipar) const;
320 
321  /// query if a parameter is fixed
322  bool IsParameterFixed(unsigned int ipar) const;
323 
324  /// retrieve parameter bounds - return false if parameter is not bound
325  bool ParameterBounds(unsigned int ipar, double &lower, double &upper) const;
326 
327 
328  /// get name of parameter (deprecated)
329  std::string GetParameterName(unsigned int ipar) const {
330  return ParName(ipar);
331  }
332 
333 
334 protected:
335 
336 
337  /// Return pointer non const pointer to model (fit) function with fitted parameter values.
338  /// used by Fitter class
339  std::shared_ptr<IModelFunction> ModelFunction() { return fFitFunc; }
340  void SetModelFunction(const std::shared_ptr<IModelFunction> & func) { fFitFunc = func; }
341 
342 
343  friend class Fitter;
344 
345 
346  bool fValid; // flag for indicating valid fit
347  bool fNormalized; // flag for indicating is errors are normalized
348  unsigned int fNFree; // number of fit free parameters (total parameters are in size of parameter vector)
349  unsigned int fNdf; // number of degree of freedom
350  unsigned int fNCalls; // number of function calls
351  int fStatus; // minimizer status code
352  int fCovStatus; // covariance matrix status code
353  double fVal; // minimum function value
354  double fEdm; // expected distance from mimimum
355  double fChi2; // fit chi2 value (different than fval in case of chi2 fits)
356  std::shared_ptr<ROOT::Math::Minimizer> fMinimizer; //! minimizer object used for fitting
357  std::shared_ptr<ROOT::Math::IMultiGenFunction> fObjFunc; //! objective function used for fitting
358  std::shared_ptr<IModelFunction> fFitFunc; //! model function resulting from the fit.
359  std::shared_ptr<FitData> fFitData; //! data set used in the fit
360  std::map<unsigned int, bool> fFixedParams; // list of fixed parameters
361  std::map<unsigned int, unsigned int> fBoundParams; // list of limited parameters
362  std::vector<std::pair<double,double> > fParamBounds; // parameter bounds
363  std::vector<double> fParams; // parameter values. Size is total number of parameters
364  std::vector<double> fErrors; // errors
365  std::vector<double> fCovMatrix; // covariance matrix (size is npar*(npar+1)/2) where npar is total parameters
366  std::vector<double> fGlobalCC; // global Correlation coefficient
367  std::map<unsigned int, std::pair<double,double> > fMinosErrors; // map contains the two Minos errors
368  std::string fMinimType; // string indicating type of minimizer
369  std::vector<std::string> fParNames; // parameter names (only with FCN only fits, when fFitFunc=0)
370 
371 };
372 
373 
374  } // end namespace Fit
375 
376 } // end namespace ROOT
377 
378 
379 
380 
381 
382 #endif /* ROOT_Fit_FitResult */
std::vector< double > fGlobalCC
Definition: FitResult.h:366
bool IsValid() const
Definition: FitResult.h:115
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
Definition: FitResult.h:118
unsigned int fNFree
Definition: FitResult.h:348
std::map< unsigned int, unsigned int > fBoundParams
Definition: FitResult.h:361
double ParError(unsigned int i) const
parameter error by index
Definition: FitResult.h:191
FitResult & operator=(const FitResult &rhs)
const std::string & MinimizerType() const
minimizer type
Definition: FitResult.h:103
double UpperError(unsigned int i) const
upper Minos error. If Minos has not run for parameter i return the parabolic error
void FillResult(const std::shared_ptr< ROOT::Math::Minimizer > &min, const FitConfig &fconfig, const std::shared_ptr< IModelFunction > &f, bool isValid, unsigned int sizeOfData=0, bool binFit=true, const ROOT::Math::IMultiGenFunction *chi2func=0, unsigned int ncalls=0)
std::vector< double > fErrors
Definition: FitResult.h:364
bool NormalizedErrors() const
flag to chek if errors are normalized
Definition: FitResult.h:310
std::shared_ptr< ROOT::Math::Minimizer > fMinimizer
Definition: FitResult.h:356
bool IsParameterFixed(unsigned int ipar) const
query if a parameter is fixed
void GetConfidenceIntervals(const BinData &data, double *ci, double cl=0.95, bool norm=false) const
unsigned int fNdf
Definition: FitResult.h:349
ROOT::Math::IParamMultiFunction IModelFunction
Definition: FitResult.h:52
const BinData * FittedBinData() const
return BinData used in the fit (return a nullptr in case a different fit is done or the data are not ...
double Error(unsigned int i) const
parameter error by index
Definition: FitResult.h:187
double CovMatrix(unsigned int i, unsigned int j) const
retrieve covariance matrix element
Definition: FitResult.h:217
void GetConfidenceIntervals(unsigned int n, unsigned int stride1, unsigned int stride2, const double *x, double *ci, double cl=0.95, bool norm=false) const
double Value(unsigned int i) const
parameter value by index
Definition: FitResult.h:180
bool Scan(unsigned int ipar, unsigned int &npoints, double *pntsx, double *pntsy, double xmin=0, double xmax=0)
FitResult(const FitResult &rhs)
std::shared_ptr< FitData > fFitData
model function resulting from the fit.
Definition: FitResult.h:359
std::string GetParameterName(unsigned int ipar) const
get name of parameter (deprecated)
Definition: FitResult.h:329
bool ParameterBounds(unsigned int ipar, double &lower, double &upper) const
retrieve parameter bounds - return false if parameter is not bound
unsigned int Ndf() const
Number of degree of freedom.
Definition: FitResult.h:164
double Chi2() const
Chi2 fit value in case of likelihood must be computed ?
Definition: FitResult.h:161
std::vector< double > fParams
Definition: FitResult.h:363
std::vector< double > fCovMatrix
Definition: FitResult.h:365
std::vector< double > GetConfidenceIntervals(double cl=0.95, bool norm=false) const
void SetMinosError(unsigned int i, double elow, double eup)
set the Minos errors for parameter i (called by the Fitter class when running Minos)
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
double LowerError(unsigned int i) const
lower Minos error. If Minos has not run for parameter i return the parabolic error
void PrintCovMatrix(std::ostream &os) const
print error matrix and correlations
unsigned int fNCalls
Definition: FitResult.h:350
const double * GetErrors() const
parameter errors (return const pointer)
Definition: FitResult.h:172
bool Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double *pntsx, double *pntsy, double confLevel=0.683)
void GetCorrelationMatrix(Matrix &mat) const
fill a correlation matrix elements using a generic symmetric matrix class implementing operator(i,...
Definition: FitResult.h:251
bool HasMinosError(unsigned int i) const
query if parameter i has the Minos error
std::vector< std::pair< double, double > > fParamBounds
Definition: FitResult.h:362
double MinFcnValue() const
Return value of the objective function (chi2 or likelihood) used in the fit.
Definition: FitResult.h:121
unsigned int NPar() const
total number of parameters (abbreviation)
Definition: FitResult.h:132
FitResult(const FitConfig &fconfig)
int CovMatrixStatus() const
covariance matrix status code using Minuit convention : =0 not calculated, =1 approximated,...
Definition: FitResult.h:143
double Edm() const
Expected distance from minimum.
Definition: FitResult.h:127
std::shared_ptr< IModelFunction > fFitFunc
objective function used for fitting
Definition: FitResult.h:358
std::map< unsigned int, bool > fFixedParams
data set used in the fit
Definition: FitResult.h:360
std::string fMinimType
Definition: FitResult.h:368
const std::vector< double > & Errors() const
parameter errors (return st::vector)
Definition: FitResult.h:170
unsigned int NCalls() const
Number of function calls to find minimum.
Definition: FitResult.h:124
double Correlation(unsigned int i, unsigned int j) const
retrieve correlation elements
Definition: FitResult.h:227
bool Update(const std::shared_ptr< ROOT::Math::Minimizer > &min, bool isValid, unsigned int ncalls=0)
std::shared_ptr< ROOT::Math::IMultiGenFunction > fObjFunc
minimizer object used for fitting
Definition: FitResult.h:357
unsigned int NTotalParameters() const
get total number of parameters
Definition: FitResult.h:130
int Index(const std::string &name) const
get index for parameter name (return -1 if not found)
const IModelFunction * FittedFunction() const
Return pointer to model (fit) function with fitted parameter values. Pointer is managed internally....
Definition: FitResult.h:149
double Prob() const
p value of the fit (chi2 probability)
std::string ParName(unsigned int i) const
name of the parameter
void GetCovarianceMatrix(Matrix &mat) const
fill covariance matrix elements using a generic matrix class implementing operator(i,...
Definition: FitResult.h:237
void NormalizeErrors()
normalize errors using chi2/ndf for chi2 fits
void SetModelFunction(const std::shared_ptr< IModelFunction > &func)
Definition: FitResult.h:340
unsigned int NFreeParameters() const
get total number of free parameters
Definition: FitResult.h:135
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Definition: FitResult.h:175
bool IsParameterBound(unsigned int ipar) const
query if a parameter is bound
const double * GetParams() const
parameter values (return const pointer)
Definition: FitResult.h:177
double Parameter(unsigned int i) const
parameter value by index
Definition: FitResult.h:182
int Status() const
minimizer status code
Definition: FitResult.h:138
std::vector< std::string > fParNames
Definition: FitResult.h:369
std::map< unsigned int, std::pair< double, double > > fMinosErrors
Definition: FitResult.h:367
std::shared_ptr< IModelFunction > ModelFunction()
Return pointer non const pointer to model (fit) function with fitted parameter values....
Definition: FitResult.h:339
double GlobalCC(unsigned int i) const
parameter global correlation coefficient
Definition: FitResult.h:211
Various mathematical functions.
Definition: Bessel.h:26
Definition: TUUID.h:7