13 #ifndef ROOT_Fit_FitResult 14 #define ROOT_Fit_FitResult 88 void FillResult(
const std::shared_ptr<ROOT::Math::Minimizer> & min,
const FitConfig & fconfig,
const std::shared_ptr<IModelFunction> & f,
98 bool Update(
const std::shared_ptr<ROOT::Math::Minimizer> & min,
bool isValid,
unsigned int ncalls = 0 );
118 bool IsEmpty()
const {
return (fParams.size() == 0); }
124 unsigned int NCalls()
const {
return fNCalls; }
127 double Edm()
const {
return fEdm; }
132 unsigned int NPar()
const {
return NTotalParameters(); }
150 return fFitFunc.get();
156 const BinData * FittedBinData()
const;
161 double Chi2()
const {
return fChi2; }
164 unsigned int Ndf()
const {
return fNdf; }
170 const std::vector<double> &
Errors()
const {
return fErrors; }
172 const double *
GetErrors()
const {
return (fErrors.empty()) ? 0 : &fErrors.front(); }
175 const std::vector<double> &
Parameters()
const {
return fParams; }
177 const double *
GetParams()
const {
return &fParams.front(); }
180 double Value(
unsigned int i)
const {
return fParams[i]; }
182 double Parameter(
unsigned int i)
const {
return fParams[i]; }
187 double Error(
unsigned int i)
const {
188 return (i < fErrors.size() ) ? fErrors[i] : 0;
192 return (i < fErrors.size() ) ? fErrors[i] : 0;
196 std::string ParName(
unsigned int i)
const;
199 void SetMinosError(
unsigned int i,
double elow,
double eup);
202 bool HasMinosError(
unsigned int i)
const;
205 double LowerError(
unsigned int i)
const;
208 double UpperError(
unsigned int i)
const;
212 return (i < fGlobalCC.size() ) ? fGlobalCC[i] : -1;
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;
221 return fCovMatrix[j + i* (i+1) / 2];
223 return fCovMatrix[i + j* (j+1) / 2];
228 if ( i >= fErrors.size() || j >= fErrors.size() )
return 0;
229 if (fCovMatrix.size() == 0)
return 0;
230 double tmp = CovMatrix(i,i)*CovMatrix(j,j);
231 return ( tmp > 0) ? CovMatrix(i,j)/
std::sqrt(tmp) : 0;
236 template<
class Matrix>
238 unsigned int npar = fErrors.size();
239 if (fCovMatrix.size() != npar*(npar+1)/2 )
return;
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);
250 template<
class Matrix>
252 unsigned int npar = fErrors.size();
253 if (fCovMatrix.size() != npar*(npar+1)/2)
return;
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);
275 void GetConfidenceIntervals(
unsigned int n,
unsigned int stride1,
unsigned int stride2,
const double *
x,
double * ci,
double cl=0.95,
bool norm =
true )
const;
293 int Index(
const std::string &
name)
const;
297 void NormalizeErrors();
303 void Print(std::ostream & os,
bool covmat =
false)
const;
306 void PrintCovMatrix(std::ostream & os)
const;
309 bool IsParameterBound(
unsigned int ipar)
const;
312 bool IsParameterFixed(
unsigned int ipar)
const;
315 bool ParameterBounds(
unsigned int ipar,
double &lower,
double &upper)
const;
320 return ParName(ipar);
347 std::shared_ptr<ROOT::Math::IMultiGenFunction>
fObjFunc;
std::shared_ptr< ROOT::Math::IMultiGenFunction > fObjFunc
minimizer object used for fitting
Namespace for new ROOT classes and functions.
double CovMatrix(unsigned int i, unsigned int j) const
retrieve covariance matrix element
double Error(unsigned int i) const
parameter error by index
const std::vector< double > & Errors() const
parameter errors (return st::vector)
double Edm() const
Expected distance from minimum.
unsigned int NPar() const
total number of parameters (abbreviation)
double Value(unsigned int i) const
parameter value by index
unsigned int Ndf() const
Number of degree of freedom.
double Parameter(unsigned int i) const
parameter value by index
double MinFcnValue() const
Return value of the objective function (chi2 or likelihood) used in the fit.
const IModelFunction * FittedFunction() const
fitting quantities
void SetModelFunction(const std::shared_ptr< IModelFunction > &func)
unsigned int NTotalParameters() const
get total number of parameters
Double_t Prob(Double_t chi2, Int_t ndf)
Computation of the probability for a certain Chi-squared (chi2) and number of degrees of freedom (ndf...
bool GetConfidenceIntervals(const TH1 *h1, const ROOT::Fit::FitResult &r, TGraphErrors *gr, double cl=0.95)
compute confidence intervals at level cl for a fitted histogram h1 in a TGraphErrors gr ...
std::map< unsigned int, bool > fFixedParams
data set used in the fit
std::vector< double > fErrors
unsigned int NCalls() const
Number of function calls to find minimum.
std::shared_ptr< FitData > fFitData
model function resulting from the fit.
std::string GetParameterName(unsigned int ipar) const
get name of parameter (deprecated)
bool NormalizedErrors() const
flag to chek if errors are normalized
std::shared_ptr< ROOT::Math::Minimizer > fMinimizer
Documentation for the abstract class IBaseFunctionMultiDim.
std::shared_ptr< IModelFunction > ModelFunction()
Return pointer non const pointer to model (fit) function with fitted parameter values.
std::shared_ptr< IModelFunction > fFitFunc
objective function used for fitting
double Correlation(unsigned int i, unsigned int j) const
retrieve correlation elements
std::map< unsigned int, std::pair< double, double > > fMinosErrors
unsigned int NFreeParameters() const
get total number of free parameters
RooCmdArg Minimizer(const char *type, const char *alg=0)
Fitter class, entry point for performing all type of fits.
const double * GetParams() const
parameter values (return const pointer)
Class describing the binned data sets : vectors of x coordinates, y values and optionally error on y ...
int Status() const
minimizer status code
int CovMatrixStatus() const
covariance matrix status code using Minuit convention : =0 not calculated, =1 approximated, =2 made pos def , =3 accurate
RooCmdArg Index(RooCategory &icat)
const std::string & MinimizerType() const
minimization quantities
class containg the result of the fit and all the related information (fitted parameter values...
void GetCovarianceMatrix(Matrix &mat) const
fill covariance matrix elements using a generic matrix class implementing operator(i,j) the matrix must be previously allocates with right size (npar * npar)
void Print(std::ostream &os, const OptionType &opt)
bool IsValid() const
True if fit successful, otherwise false.
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
const std::vector< double > & Parameters() const
parameter values (return std::vector)
Namespace for new Math classes and functions.
std::vector< std::pair< double, double > > fParamBounds
Binding & operator=(OUT(*fun)(void))
bool IsEmpty() const
True if a fit result does not exist (even invalid) with parameter values.
void GetCorrelationMatrix(Matrix &mat) const
fill a correlation matrix elements using a generic symmetric matrix class implementing operator(i...
std::vector< std::string > fParNames
std::vector< double > fGlobalCC
std::vector< double > fCovMatrix
ROOT::Math::IParamMultiFunction IModelFunction
std::map< unsigned int, unsigned int > fBoundParams
std::vector< double > fParams
double GlobalCC(unsigned int i) const
parameter global correlation coefficient
double Chi2() const
Chi2 fit value in case of likelihood must be computed ?
double ParError(unsigned int i) const
parameter error by index
Class describing the configuration of the fit, options and parameter settings using the ROOT::Fit::Pa...
const double * GetErrors() const
parameter errors (return const pointer)