13 #ifndef ROOT_Fit_FitResult 14 #define ROOT_Fit_FitResult 16 #ifndef ROOT_Fit_IFunctionfwd 19 #ifndef ROOT_Fit_IParamFunctionfwd 92 void FillResult(
const std::shared_ptr<ROOT::Math::Minimizer> & min,
const FitConfig & fconfig,
const std::shared_ptr<IModelFunction> &
f,
102 bool Update(
const std::shared_ptr<ROOT::Math::Minimizer> & min,
bool isValid,
unsigned int ncalls = 0 );
122 bool IsEmpty()
const {
return (fParams.size() == 0); }
128 unsigned int NCalls()
const {
return fNCalls; }
131 double Edm()
const {
return fEdm; }
136 unsigned int NPar()
const {
return NTotalParameters(); }
154 return fFitFunc.get();
160 const BinData * FittedBinData()
const;
165 double Chi2()
const {
return fChi2; }
168 unsigned int Ndf()
const {
return fNdf; }
174 const std::vector<double> &
Errors()
const {
return fErrors; }
176 const double *
GetErrors()
const {
return (fErrors.empty()) ? 0 : &fErrors.front(); }
179 const std::vector<double> &
Parameters()
const {
return fParams; }
181 const double *
GetParams()
const {
return &fParams.front(); }
184 double Value(
unsigned int i)
const {
return fParams[i]; }
186 double Parameter(
unsigned int i)
const {
return fParams[i]; }
191 double Error(
unsigned int i)
const {
192 return (i < fErrors.size() ) ? fErrors[i] : 0;
196 return (i < fErrors.size() ) ? fErrors[i] : 0;
200 std::string ParName(
unsigned int i)
const;
203 void SetMinosError(
unsigned int i,
double elow,
double eup);
206 bool HasMinosError(
unsigned int i)
const;
209 double LowerError(
unsigned int i)
const;
212 double UpperError(
unsigned int i)
const;
216 return (i < fGlobalCC.size() ) ? fGlobalCC[i] : -1;
221 double CovMatrix (
unsigned int i,
unsigned int j)
const {
222 if ( i >= fErrors.size() || j >= fErrors.size() )
return 0;
223 if (fCovMatrix.size() == 0)
return 0;
225 return fCovMatrix[j + i* (i+1) / 2];
227 return fCovMatrix[i + j* (j+1) / 2];
232 if ( i >= fErrors.size() || j >= fErrors.size() )
return 0;
233 if (fCovMatrix.size() == 0)
return 0;
234 double tmp = CovMatrix(i,i)*CovMatrix(j,j);
235 return ( tmp > 0) ? CovMatrix(i,j)/
std::sqrt(tmp) : 0;
240 template<
class Matrix>
242 unsigned int npar = fErrors.size();
243 if (fCovMatrix.size() != npar*(npar+1)/2 )
return;
244 for (
unsigned int i = 0; i< npar; ++i) {
245 for (
unsigned int j = 0; j<=i; ++j) {
246 mat(i,j) = fCovMatrix[j + i*(i+1)/2 ];
247 if (i != j) mat(j,i) = mat(i,j);
254 template<
class Matrix>
256 unsigned int npar = fErrors.size();
257 if (fCovMatrix.size() != npar*(npar+1)/2)
return;
258 for (
unsigned int i = 0; i< npar; ++i) {
259 for (
unsigned int j = 0; j<=i; ++j) {
260 double tmp = fCovMatrix[i * (i +3)/2 ] * fCovMatrix[ j * (j+3)/2 ];
261 mat(i,j) = (tmp > 0) ? fCovMatrix[j + i*(i+1)/2 ] /
std::sqrt(tmp) : 0;
262 if (i != j) mat(j,i) = mat(i,j);
279 void GetConfidenceIntervals(
unsigned int n,
unsigned int stride1,
unsigned int stride2,
const double *
x,
double * ci,
double cl=0.95,
bool norm =
true )
const;
297 int Index(
const std::string &
name)
const;
301 void NormalizeErrors();
307 void Print(std::ostream & os,
bool covmat =
false)
const;
310 void PrintCovMatrix(std::ostream & os)
const;
313 bool IsParameterBound(
unsigned int ipar)
const;
316 bool IsParameterFixed(
unsigned int ipar)
const;
319 bool ParameterBounds(
unsigned int ipar,
double &lower,
double &upper)
const;
324 return ParName(ipar);
351 std::shared_ptr<ROOT::Math::IMultiGenFunction>
fObjFunc;
std::shared_ptr< ROOT::Math::IMultiGenFunction > fObjFunc
minimizer object used for fitting
This namespace contains pre-defined functions to be used in conjuction with TExecutor::Map and TExecu...
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
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
IParamFunction interface (abstract class) describing multi-dimensional parameteric functions It is a ...
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)
double func(double *x, double *p)
Namespace for new Math classes and functions.
std::vector< std::pair< double, double > > fParamBounds
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
Documentation for the abstract class IBaseFunctionMultiDim.
double GlobalCC(unsigned int i) const
parameter global correlation coefficient
double norm(double *x, double *p)
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)