61 if (prevErrorIgnoreLevel < 1001) {
63 return prevErrorIgnoreLevel;
86 : Minimizer(), fDim(0), fMinimizer(0), fMinuitFCN(0), fMinimum(0)
96 std::string algoname(
type);
98 std::transform(algoname.begin(), algoname.end(), algoname.begin(), (
int (*)(
int))tolower);
101 if (algoname ==
"simplex")
103 if (algoname ==
"minimize")
105 if (algoname ==
"scan")
107 if (algoname ==
"fumili")
109 if (algoname ==
"bfgs")
192 print.
Info(
"Parameter",
name,
"has zero or invalid step size - consider it as constant");
198 if (minuit2Index != ivar) {
199 print.
Warn(
"Wrong index", minuit2Index,
"used for the variable",
name);
229 double lower,
double upper)
243 double step = (val != 0) ? 0.1 * std::abs(val) : 0.1;
255 return std::string();
281 for (
unsigned int ivar = 0; ivar <
n; ++ivar)
348 print.
Error(
"Wrong variable index");
359 print.
Error(
"Wrong variable index");
391 print.
Error(
"Wrong Fit method function for Fumili");
411 print.
Error(
"Wrong Fit method function for Fumili");
426 print.
Error(
"FCN function has not been set");
439 const int strategyLevel =
Strategy();
445 int maxfcn_used = maxfcn;
446 if (maxfcn_used == 0) {
448 maxfcn_used = 200 + 100 * nvar + 5 * nvar * nvar;
450 std::cout <<
"Minuit2Minimizer: Minimize with max-calls " << maxfcn_used <<
" convergence for edm < " << tol
451 <<
" strategy " << strategyLevel << std::endl;
479 minuit2Opt->
GetValue(
"GradientNCycles", nGradCycles);
480 minuit2Opt->
GetValue(
"HessianNCycles", nHessCycles);
481 minuit2Opt->
GetValue(
"HessianGradientNCycles", nHessGradCycles);
483 minuit2Opt->
GetValue(
"GradientTolerance", gradTol);
484 minuit2Opt->
GetValue(
"GradientStepTolerance", gradStepTol);
485 minuit2Opt->
GetValue(
"HessianStepTolerance", hessStepTol);
486 minuit2Opt->
GetValue(
"HessianG2Tolerance", hessG2Tol);
497 int storageLevel = 1;
498 bool ret = minuit2Opt->
GetValue(
"StorageLevel", storageLevel);
502 if (printLevel > 0) {
503 std::cout <<
"Minuit2Minimizer::Minuit - Changing default options" << std::endl;
512 if (printLevel == 10 &&
gROOT) {
520 if (printLevel == 20 || printLevel == 30 || printLevel == 40 || (printLevel >= 20000 && printLevel < 30000)) {
521 int parNumber = printLevel - 20000;
522 if (printLevel == 20)
524 if (printLevel == 30)
526 if (printLevel == 40)
531 if (printLevel == 100 || (printLevel >= 10000 && printLevel < 20000)) {
532 int parNumber = printLevel - 10000;
581 if (debugLevel >= 3) {
583 const std::vector<ROOT::Minuit2::MinimumState> &iterationStates = min.
States();
584 std::cout <<
"Number of iterations " << iterationStates.size() << std::endl;
585 for (
unsigned int i = 0; i < iterationStates.size(); ++i) {
588 std::cout <<
"----------> Iteration " << i << std::endl;
589 int pr = std::cout.precision(12);
590 std::cout <<
" FVAL = " << st.
Fval() <<
" Edm = " << st.
Edm() <<
" Nfcn = " << st.
NFcn()
592 std::cout.precision(pr);
594 std::cout <<
" Error matrix change = " << st.
Error().
Dcovar() << std::endl;
596 std::cout <<
" Parameters : ";
598 for (
int j = 0; j < st.
size(); ++j)
600 std::cout << std::endl;
610 txt =
"Covar is not pos def";
614 txt =
"Covar was made pos def";
618 txt =
"Hesse is not valid";
622 txt =
"Edm is above max";
626 txt =
"Reached call limit";
630 MnPrint print(
"Minuit2Minimizer::Minimize", debugLevel);
631 bool validMinimum = min.
IsValid();
634 if (
fStatus != 0 && debugLevel > 0)
640 txt =
"unknown failure";
643 print.
Warn(
"Minimization did NOT converge,", txt);
651 if (paramsObj.size() == 0)
653 assert(
fDim == paramsObj.size());
657 for (
unsigned int i = 0; i <
fDim; ++i) {
658 fValues[i] = paramsObj[i].Value();
671 std::cout <<
"Minuit2Minimizer : Valid minimum - status = " <<
fStatus << std::endl;
672 int pr = std::cout.precision(18);
673 std::cout <<
"FVAL = " <<
fState.
Fval() << std::endl;
674 std::cout <<
"Edm = " <<
fState.
Edm() << std::endl;
675 std::cout.precision(pr);
676 std::cout <<
"Nfcn = " <<
fState.
NFcn() << std::endl;
679 std::cout << par.
Name() <<
"\t = " << par.
Value() <<
"\t ";
681 std::cout <<
"(fixed)" << std::endl;
683 std::cout <<
"(const)" << std::endl;
685 std::cout <<
"+/- " << par.
Error() <<
"\t(limited)" << std::endl;
687 std::cout <<
"+/- " << par.
Error() << std::endl;
690 std::cout <<
"Minuit2Minimizer : Invalid Minimum - status = " <<
fStatus << std::endl;
691 std::cout <<
"FVAL = " <<
fState.
Fval() << std::endl;
692 std::cout <<
"Edm = " <<
fState.
Edm() << std::endl;
693 std::cout <<
"Nfcn = " <<
fState.
NFcn() << std::endl;
701 if (paramsObj.size() == 0)
703 assert(
fDim == paramsObj.size());
708 for (
unsigned int i = 0; i <
fDim; ++i) {
740 for (
unsigned int i = 0; i <
fDim; ++i) {
742 for (
unsigned int j = 0; j <
fDim; ++j) {
743 cov[i *
fDim + j] = 0;
747 for (
unsigned int j = 0; j <
fDim; ++j) {
749 int k = i *
fDim + j;
770 for (
unsigned int i = 0; i <
fDim; ++i) {
772 for (
unsigned int j = 0; j <
fDim; ++j) {
773 hess[i *
fDim + j] = 0;
777 for (
unsigned int j = 0; j <
fDim; ++j) {
779 int k = i *
fDim + j;
857 print.
Error(
"Failed - no function minimum existing");
862 print.
Error(
"Failed - invalid function minimum");
875 if ((mstatus & 8) != 0) {
876 print.
Info([&](std::ostream &os) {
877 os <<
"Found a new minimum: run again the Minimization starting from the new point";
880 os <<
'\n' << par.Name() <<
"\t = " << par.Value();
889 print.
Info(
"Run now again Minos from the new found Minimum");
899 bool isValid = ((mstatus & 1) == 0) && ((mstatus & 2) == 0);
906 bool runLower = runopt != 2;
907 bool runUpper = runopt != 1;
931 tol = std::max(tol, 0.01);
934 int maxfcn_used = maxfcn;
935 if (maxfcn_used == 0) {
937 maxfcn_used = 2 * (nvar + 1) * (200 + 100 * nvar + 5 * nvar * nvar);
941 if (debugLevel >= 1) {
942 std::cout <<
"************************************************************************************************"
944 std::cout <<
"Minuit2Minimizer::GetMinosError - Run MINOS LOWER error for parameter #" << i <<
" : "
945 << par_name <<
" using max-calls " << maxfcn_used <<
", tolerance " << tol << std::endl;
947 low = minos.
Loval(i, maxfcn, tol);
950 if (debugLevel >= 1) {
951 std::cout <<
"************************************************************************************************"
953 std::cout <<
"Minuit2Minimizer::GetMinosError - Run MINOS UPPER error for parameter #" << i <<
" : "
954 << par_name <<
" using max-calls " << maxfcn_used <<
", tolerance " << tol << std::endl;
956 up = minos.
Upval(i, maxfcn, tol);
973 if (debugLevel > 0) {
976 std::cout <<
"Minos: Invalid lower error for parameter " << par_name << std::endl;
978 std::cout <<
"Minos: Parameter : " << par_name <<
" is at Lower limit; error is " << me.
Lower()
981 std::cout <<
"Minos: Maximum number of function calls exceeded when running for lower error for parameter "
982 << par_name << std::endl;
984 std::cout <<
"Minos: New Minimum found while running Minos for lower error for parameter " << par_name
988 std::cout <<
"Minos: Lower error for parameter " << par_name <<
" : " << me.
Lower() << std::endl;
992 std::cout <<
"Minos: Invalid upper error for parameter " << par_name << std::endl;
994 std::cout <<
"Minos: Parameter " << par_name <<
" is at Upper limit; error is " << me.
Upper() << std::endl;
996 std::cout <<
"Minos: Maximum number of function calls exceeded when running for upper error for parameter "
997 << par_name << std::endl;
999 std::cout <<
"Minos: New Minimum found while running Minos for upper error for parameter " << par_name
1003 std::cout <<
"Minos: Upper error for parameter " << par_name <<
" : " << me.
Upper() << std::endl;
1008 bool lowerInvalid = (runLower && !me.
LowerValid());
1009 bool upperInvalid = (runUpper && !me.
UpperValid());
1026 if (lowerInvalid || upperInvalid) {
1052 errLow = me.
Lower();
1078 print.
Error(
"Function must be set before using Scan");
1083 print.
Error(
"Invalid number; minimizer variables must be set before using Scan");
1096 double amin = scan.
Fval();
1099 std::vector<std::pair<double, double>> result = scan(ipar, nstep - 1,
xmin,
xmax);
1102 if (prev_level > -2)
1106 if (result.size() != nstep) {
1107 print.
Error(
"Invalid result from MnParameterScan");
1111 std::sort(result.begin(), result.end());
1113 for (
unsigned int i = 0; i < nstep; ++i) {
1114 x[i] = result[i].first;
1115 y[i] = result[i].second;
1120 if (scan.
Fval() < amin) {
1121 print.
Info(
"A new minimum has been found");
1136 print.
Error(
"No function minimum existing; must minimize function before");
1141 print.
Error(
"Invalid function minimum");
1166 if (prev_level > -2)
1171 std::vector<std::pair<double, double>> result = contour(ipar, jpar, npoints);
1172 if (result.size() != npoints) {
1173 print.
Error(
"Invalid result from MnContours");
1176 for (
unsigned int i = 0; i < npoints; ++i) {
1177 x[i] = result[i].first;
1178 y[i] = result[i].second;
1194 print.
Error(
"FCN function has not been set");
1200 print.
Info(
"Using max-calls", maxfcn);
1233 if (prev_level > -2)
1238 std::cout <<
"Minuit2Minimizer::Hesse - State returned from Hesse " << std::endl;
1239 std::cout <<
fState << std::endl;
1243 std::string covStatusType =
"not valid";
1245 covStatusType =
"approximate";
1247 covStatusType =
"full but made positive defined";
1249 covStatusType =
"accurate";
1265 print.
Warn(
"Hesse failed - matrix is", covStatusType);
1266 print.
Warn(hstatus);
1272 print.
Info(
"Hesse is valid - matrix is", covStatusType);
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
R__EXTERN Int_t gErrorIgnoreLevel
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
void Set(const std::string &name, double value, double step)
set value and name (unlimited parameter)
void SetLimits(double low, double up)
set a double side limit, if low == up the parameter is fixed if low > up the limits are removed The c...
void SetUpperLimit(double up)
set a single upper limit
void Fix()
fix the parameter
void SetLowerLimit(double low)
set a single lower limit
FitMethodFunction class Interface for objective functions (like chi2 and likelihood used in the fit) ...
Documentation for the abstract class IBaseFunctionMultiDim.
virtual unsigned int NDim() const =0
Retrieve the dimension of the function.
Interface (abstract class) for multi-dimensional functions providing a gradient calculation.
virtual unsigned int NDim() const=0
Retrieve the dimension of the function.
Generic interface for defining configuration options of a numerical algorithm.
virtual void Print(std::ostream &=std::cout) const
print options
bool GetValue(const char *name, T &t) const
static ROOT::Math::IOptions * FindDefault(const char *name)
double Tolerance() const
absolute tolerance
unsigned int MaxFunctionCalls() const
max number of function calls
double Precision() const
precision of minimizer in the evaluation of the objective function ( a value <=0 corresponds to the l...
int Strategy() const
strategy
double ErrorDef() const
return the statistical scale used for calculate the error is typically 1 for Chi2 and 0....
bool IsValidError() const
return true if Minimizer has performed a detailed error validation (e.g. run Hesse for Minuit)
int PrintLevel() const
minimizer configuration parameters
Combined minimizer: combination of Migrad and Simplex.
template wrapped class for adapting to FCNBase signature
virtual void SetErrorDef(double)
add interface to set dynamically a new error definition Re-implement this function if needed.
template wrapped class for adapting to FCNBase signature a IGradFunction
template wrapped class for adapting to FumiliFCNBase signature
Instantiates the seed generator and Minimum builder for the Fumili minimization method.
class holding the full result of the minimization; both internal and external (MnUserParameterState) ...
const std::vector< MinimumState > & States() const
const MinimumError & Error() const
bool HasAccurateCovar() const
bool HasReachedCallLimit() const
const MnUserParameterState & UserState() const
const MinimumState & State() const
bool HasMadePosDefCovar() const
bool IsAboveMaxEdm() const
bool HasValidCovariance() const
bool HasPosDefCovar() const
bool HasCovariance() const
void SetErrorDef(double up)
virtual void SetStorageLevel(int level)
virtual void SetPrintLevel(int level)
virtual void SetTraceObject(MnTraceObject &obj)
bool InvertFailed() const
MinimumState keeps the information (position, Gradient, 2nd deriv, etc) after one minimization step (...
bool HasParameters() const
const MinimumError & Error() const
const MnAlgebraicVector & Vec() const
bool HasCovariance() const
Class holding the result of Minos (lower and upper values) for a specific parameter.
bool AtUpperLimit() const
bool AtLowerMaxFcn() const
bool AtUpperMaxFcn() const
bool AtLowerLimit() const
Minuit2Minimizer class implementing the ROOT::Math::Minimizer interface for Minuit2 minimization algo...
bool ExamineMinimum(const ROOT::Minuit2::FunctionMinimum &min)
examine the minimum result
Minuit2Minimizer(ROOT::Minuit2::EMinimizerType type=ROOT::Minuit2::kMigrad)
Default constructor.
void SetStorageLevel(int level)
set storage level = 1 : store all iteration states (default) = 0 : store only first and last state to...
Minuit2Minimizer & operator=(const Minuit2Minimizer &rhs)
Assignment operator.
std::vector< double > fValues
virtual bool SetVariableUpperLimit(unsigned int ivar, double upper)
set the upper-limit of an already existing variable
virtual double GlobalCC(unsigned int i) const
get global correlation coefficient for the variable i.
virtual int VariableIndex(const std::string &name) const
get index of variable given a variable given a name return -1 if variable is not found
virtual bool SetVariable(unsigned int ivar, const std::string &name, double val, double step)
set free variable
virtual bool SetFixedVariable(unsigned int, const std::string &, double)
set fixed variable (override if minimizer supports them )
virtual void SetFunction(const ROOT::Math::IMultiGenFunction &func)
set the function to minimize
virtual const ROOT::Minuit2::FCNBase * GetFCN() const
virtual bool SetLowerLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double lower)
set lower limit variable (override if minimizer supports them )
virtual bool SetLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double, double)
set upper/lower limited variable (override if minimizer supports them )
virtual bool SetVariableLimits(unsigned int ivar, double lower, double upper)
set the limits of an already existing variable
virtual bool SetVariableValues(const double *val)
set the values of all existing variables (array must be dimensioned to the size of the existing param...
virtual int CovMatrixStatus() const
return the status of the covariance matrix status = -1 : not available (inversion failed or Hesse fai...
virtual void Clear()
reset for consecutive minimizations - implement if needed
virtual double CovMatrix(unsigned int i, unsigned int j) const
return covariance matrix elements if the variable is fixed or const the value is zero The ordering of...
virtual ~Minuit2Minimizer()
Destructor (no operations)
virtual bool ReleaseVariable(unsigned int ivar)
release an existing variable
std::vector< double > fErrors
ROOT::Minuit2::ModularFunctionMinimizer * fMinimizer
virtual bool Scan(unsigned int i, unsigned int &nstep, double *x, double *y, double xmin=0, double xmax=0)
scan a parameter i around the minimum.
virtual std::string VariableName(unsigned int ivar) const
get name of variables (override if minimizer support storing of variable names)
int RunMinosError(unsigned int i, double &errLow, double &errUp, int runopt)
virtual bool GetCovMatrix(double *cov) const
Fill the passed array with the covariance matrix elements if the variable is fixed or const the value...
virtual bool SetUpperLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double upper)
set upper limit variable (override if minimizer supports them )
virtual bool Minimize()
method to perform the minimization.
void SetTraceObject(MnTraceObject &obj)
set an object to trace operation for each iteration The object must be a (or inherit from) ROOT::Minu...
virtual bool SetVariableLowerLimit(unsigned int ivar, double lower)
set the lower-limit of an already existing variable
virtual const ROOT::Minuit2::ModularFunctionMinimizer * GetMinimizer() const
void SetMinimizerType(ROOT::Minuit2::EMinimizerType type)
virtual void PrintResults()
return reference to the objective function virtual const ROOT::Math::IGenFunction & Function() const;
virtual bool GetHessianMatrix(double *h) const
Fill the passed array with the Hessian matrix elements The Hessian matrix is the matrix of the second...
ROOT::Minuit2::MnUserParameterState fState
virtual bool GetMinosError(unsigned int i, double &errLow, double &errUp, int=0)
get the minos error for parameter i, return false if Minos failed A minimizaiton must be performed be...
virtual bool IsFixedVariable(unsigned int ivar) const
query if an existing variable is fixed (i.e.
virtual void SetMinimizer(ROOT::Minuit2::ModularFunctionMinimizer *m)
virtual const double * Errors() const
return errors at the minimum
virtual bool GetVariableSettings(unsigned int ivar, ROOT::Fit::ParameterSettings &varObj) const
get variable settings in a variable object (like ROOT::Fit::ParamsSettings)
ROOT::Minuit2::FunctionMinimum * fMinimum
virtual bool SetVariableValue(unsigned int ivar, double val)
set variable
ROOT::Minuit2::FCNBase * fMinuitFCN
virtual bool Contour(unsigned int i, unsigned int j, unsigned int &npoints, double *xi, double *xj)
find the contour points (xi,xj) of the function for parameter i and j around the minimum The contour ...
virtual double Correlation(unsigned int i, unsigned int j) const
return correlation coefficient between variable i and j.
virtual bool SetVariableStepSize(unsigned int ivar, double step)
set the step size of an already existing variable
virtual bool FixVariable(unsigned int ivar)
fix an existing variable
virtual bool Hesse()
perform a full calculation of the Hessian matrix for error calculation If a valid minimum exists the ...
class for the individual Minuit Parameter with Name and number; contains the input numbers for the mi...
double LowerLimit() const
double UpperLimit() const
const char * Name() const
bool HasLowerLimit() const
bool HasUpperLimit() const
API class for Contours Error analysis (2-dim errors); minimization has to be done before and Minimum ...
const MnUserParameterState & State() const
const std::vector< double > & GlobalCC() const
API class for calculating the numerical covariance matrix (== 2x Inverse Hessian == 2x Inverse 2nd de...
API class for Minos Error analysis (asymmetric errors); minimization has to be done before and Minimu...
MnCross Loval(unsigned int, unsigned int maxcalls=0, double toler=0.1) const
MnCross Upval(unsigned int, unsigned int maxcalls=0, double toler=0.1) const
Scans the values of FCN as a function of one Parameter and retains the best function and Parameter va...
const MnUserParameters & Parameters() const
void Error(const Ts &... args)
void Info(const Ts &... args)
static int SetGlobalLevel(int level)
void Warn(const Ts &... args)
API class for defining three levels of strategies: low (0), medium (1), high (>=2); acts on: Migrad (...
double HessianG2Tolerance() const
unsigned int HessianGradientNCycles() const
double GradientStepTolerance() const
void SetHessianNCycles(unsigned int n)
void SetHessianStepTolerance(double stp)
double GradientTolerance() const
void SetGradientTolerance(double toler)
double HessianStepTolerance() const
unsigned int HessianNCycles() const
unsigned int GradientNCycles() const
void SetGradientNCycles(unsigned int n)
void SetGradientStepTolerance(double stp)
void SetHessianGradientNCycles(unsigned int n)
void SetHessianG2Tolerance(double toler)
virtual void Init(const MnUserParameterState &state)
class which holds the external user and/or internal Minuit representation of the parameters and error...
void SetLimits(unsigned int, double, double)
const MnUserParameters & Parameters() const
double Value(unsigned int) const
void RemoveLimits(unsigned int)
unsigned int NFcn() const
unsigned int Index(const std::string &) const
const std::string & GetName(unsigned int) const
double Int2ext(unsigned int, double) const
const MnGlobalCorrelationCoeff & GlobalCC() const
void Release(unsigned int)
unsigned int VariableParameters() const
const MinuitParameter & Parameter(unsigned int i) const
void SetValue(unsigned int, double)
void Add(const std::string &name, double val, double err)
const char * Name(unsigned int) const
MnUserCovariance Hessian() const
const std::vector< ROOT::Minuit2::MinuitParameter > & MinuitParameters() const
facade: forward interface of MnUserParameters and MnUserTransformation
void SetPrecision(double eps)
unsigned int IntOfExt(unsigned int) const
const MnUserTransformation & Trafo() const
void SetUpperLimit(unsigned int, double)
const MnUserCovariance & IntCovariance() const
const MnUserCovariance & Covariance() const
int CovarianceStatus() const
void SetError(unsigned int, double)
void SetLowerLimit(unsigned int, double)
bool HasCovariance() const
double Value(unsigned int) const
virtual const MinimumBuilder & Builder() const =0
virtual FunctionMinimum Minimize(const FCNBase &, const std::vector< double > &, const std::vector< double > &, unsigned int stra=1, unsigned int maxfcn=0, double toler=0.1) const
Class implementing the required methods for a minimization using SCAN API is provided in the upper RO...
Class implementing the required methods for a minimization using Simplex.
Instantiates the SeedGenerator and MinimumBuilder for Variable Metric Minimization method.
Mother of all ROOT objects.
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Namespace for new Math classes and functions.
void RestoreGlobalPrintLevel(int)
int TurnOffPrintInfoLevel()
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...