59 int TurnOffPrintInfoLevel() {
61 int prevErrorIgnoreLevel = gErrorIgnoreLevel;
62 if (prevErrorIgnoreLevel < 1001) {
63 gErrorIgnoreLevel = 1001;
64 return prevErrorIgnoreLevel;
69 void RestoreGlobalPrintLevel(
int value) {
70 gErrorIgnoreLevel = value;
74 int TurnOffPrintInfoLevel() {
return -1; }
75 int ControlPrintLevel( ) {
return -1;}
76 void RestoreGlobalPrintLevel(
int ) {}
82 Minuit2Minimizer::Minuit2Minimizer(ROOT::Minuit2::EMinimizerType type ) :
90 SetMinimizerType(type);
93 Minuit2Minimizer::Minuit2Minimizer(
const char * type ) :
102 std::string algoname(type);
104 std::transform(algoname.begin(), algoname.end(), algoname.begin(), (int(*)(int)) tolower );
106 EMinimizerType algoType = kMigrad;
107 if (algoname ==
"simplex") algoType = kSimplex;
108 if (algoname ==
"minimize" ) algoType = kCombined;
109 if (algoname ==
"scan" ) algoType = kScan;
110 if (algoname ==
"fumili" ) algoType = kFumili;
111 if (algoname ==
"bfgs" ) algoType = kMigradBFGS;
113 SetMinimizerType(algoType);
116 void Minuit2Minimizer::SetMinimizerType(ROOT::Minuit2::EMinimizerType type) {
120 case ROOT::Minuit2::kMigrad:
122 SetMinimizer(
new ROOT::Minuit2::VariableMetricMinimizer() );
124 case ROOT::Minuit2::kMigradBFGS:
126 SetMinimizer(
new ROOT::Minuit2::VariableMetricMinimizer(VariableMetricMinimizer::BFGSType()) );
128 case ROOT::Minuit2::kSimplex:
130 SetMinimizer(
new ROOT::Minuit2::SimplexMinimizer() );
132 case ROOT::Minuit2::kCombined:
133 SetMinimizer(
new ROOT::Minuit2::CombinedMinimizer() );
135 case ROOT::Minuit2::kScan:
136 SetMinimizer(
new ROOT::Minuit2::ScanMinimizer() );
138 case ROOT::Minuit2::kFumili:
139 SetMinimizer(
new ROOT::Minuit2::FumiliMinimizer() );
144 SetMinimizer(
new ROOT::Minuit2::VariableMetricMinimizer() );
150 Minuit2Minimizer::~Minuit2Minimizer()
153 if (fMinimizer)
delete fMinimizer;
154 if (fMinuitFCN)
delete fMinuitFCN;
155 if (fMinimum)
delete fMinimum;
158 Minuit2Minimizer::Minuit2Minimizer(
const Minuit2Minimizer &) :
159 ROOT::Math::Minimizer()
164 Minuit2Minimizer & Minuit2Minimizer::operator = (
const Minuit2Minimizer &rhs)
167 if (
this == &rhs)
return *
this;
172 void Minuit2Minimizer::Clear() {
174 fState = MnUserParameterState();
176 if (fMinimum)
delete fMinimum;
183 bool Minuit2Minimizer::SetVariable(
unsigned int ivar,
const std::string & name,
double val,
double step) {
193 std::string txtmsg =
"Parameter " + name +
" has zero or invalid step size - consider it as constant ";
194 MN_INFO_MSG2(
"Minuit2Minimizer::SetVariable",txtmsg);
195 fState.Add(name.c_str(), val);
198 fState.Add(name.c_str(), val, step);
200 unsigned int minuit2Index = fState.Index(name.c_str() );
201 if ( minuit2Index != ivar) {
202 std::string txtmsg(
"Wrong index used for the variable " + name);
203 MN_INFO_MSG2(
"Minuit2Minimizer::SetVariable",txtmsg);
204 MN_INFO_VAL2(
"Minuit2Minimizer::SetVariable",minuit2Index);
208 fState.RemoveLimits(ivar);
213 bool Minuit2Minimizer::SetLowerLimitedVariable(
unsigned int ivar ,
const std::string & name ,
double val ,
double step ,
double lower ) {
215 if (!SetVariable(ivar, name, val, step) )
return false;
216 fState.SetLowerLimit(ivar, lower);
220 bool Minuit2Minimizer::SetUpperLimitedVariable(
unsigned int ivar ,
const std::string & name ,
double val ,
double step ,
double upper ) {
222 if (!SetVariable(ivar, name, val, step) )
return false;
223 fState.SetUpperLimit(ivar, upper);
229 bool Minuit2Minimizer::SetLimitedVariable(
unsigned int ivar ,
const std::string & name ,
double val ,
double step ,
double lower ,
double upper) {
231 if (!SetVariable(ivar, name, val, step) )
return false;
232 fState.SetLimits(ivar, lower, upper);
236 bool Minuit2Minimizer::SetFixedVariable(
unsigned int ivar ,
const std::string & name ,
double val ) {
240 double step = ( val != 0) ? 0.1 * std::abs(val) : 0.1;
241 if (!SetVariable(ivar, name, val, step ) ) {
242 ivar = fState.Index(name.c_str() );
248 std::string Minuit2Minimizer::VariableName(
unsigned int ivar)
const {
250 if (ivar >= fState.MinuitParameters().size() )
return std::string();
251 return fState.GetName(ivar);
255 int Minuit2Minimizer::VariableIndex(
const std::string & name)
const {
258 return fState.Trafo().FindIndex(name);
262 bool Minuit2Minimizer::SetVariableValue(
unsigned int ivar,
double val) {
264 if (ivar >= fState.MinuitParameters().size() )
return false;
265 fState.SetValue(ivar, val);
269 bool Minuit2Minimizer::SetVariableValues(
const double * x) {
271 unsigned int n = fState.MinuitParameters().size();
272 if (n== 0)
return false;
273 for (
unsigned int ivar = 0; ivar < n; ++ivar)
274 fState.SetValue(ivar, x[ivar]);
278 bool Minuit2Minimizer::SetVariableStepSize(
unsigned int ivar,
double step) {
281 if (ivar >= fState.MinuitParameters().size() )
return false;
282 fState.SetError(ivar, step);
286 bool Minuit2Minimizer::SetVariableLowerLimit(
unsigned int ivar,
double lower) {
289 if (ivar >= fState.MinuitParameters().size() )
return false;
290 fState.SetLowerLimit(ivar, lower);
293 bool Minuit2Minimizer::SetVariableUpperLimit(
unsigned int ivar,
double upper ) {
296 if (ivar >= fState.MinuitParameters().size() )
return false;
297 fState.SetUpperLimit(ivar, upper);
301 bool Minuit2Minimizer::SetVariableLimits(
unsigned int ivar,
double lower,
double upper) {
304 if (ivar >= fState.MinuitParameters().size() )
return false;
305 fState.SetLimits(ivar, lower,upper);
309 bool Minuit2Minimizer::FixVariable(
unsigned int ivar) {
311 if (ivar >= fState.MinuitParameters().size() )
return false;
316 bool Minuit2Minimizer::ReleaseVariable(
unsigned int ivar) {
318 if (ivar >= fState.MinuitParameters().size() )
return false;
319 fState.Release(ivar);
323 bool Minuit2Minimizer::IsFixedVariable(
unsigned int ivar)
const {
325 if (ivar >= fState.MinuitParameters().size() ) {
326 MN_ERROR_MSG2(
"Minuit2Minimizer",
"wrong variable index");
329 return (fState.Parameter(ivar).IsFixed() || fState.Parameter(ivar).IsConst() );
332 bool Minuit2Minimizer::GetVariableSettings(
unsigned int ivar, ROOT::Fit::ParameterSettings & varObj)
const {
334 if (ivar >= fState.MinuitParameters().size() ) {
335 MN_ERROR_MSG2(
"Minuit2Minimizer",
"wrong variable index");
338 const MinuitParameter & par = fState.Parameter(ivar);
339 varObj.Set( par.Name(), par.Value(), par.Error() );
340 if (par.HasLowerLimit() ) {
341 if (par.HasUpperLimit() ) {
342 varObj.SetLimits(par.LowerLimit(), par.UpperLimit() );
344 varObj.SetLowerLimit(par.LowerLimit() );
346 }
else if (par.HasUpperLimit() ) {
347 varObj.SetUpperLimit(par.UpperLimit() );
349 if (par.IsConst() || par.IsFixed() ) varObj.Fix();
355 void Minuit2Minimizer::SetFunction(
const ROOT::Math::IMultiGenFunction & func) {
357 if (fMinuitFCN)
delete fMinuitFCN;
360 fMinuitFCN =
new ROOT::Minuit2::FCNAdapter<ROOT::Math::IMultiGenFunction> (func, ErrorDef() );
364 const ROOT::Math::FitMethodFunction * fcnfunc =
dynamic_cast<const ROOT::Math::FitMethodFunction *
>(&func);
366 MN_ERROR_MSG(
"Minuit2Minimizer: Wrong Fit method function for Fumili");
369 fMinuitFCN =
new ROOT::Minuit2::FumiliFCNAdapter<ROOT::Math::FitMethodFunction> (*fcnfunc, fDim, ErrorDef() );
373 void Minuit2Minimizer::SetFunction(
const ROOT::Math::IMultiGradFunction & func) {
376 if (fMinuitFCN)
delete fMinuitFCN;
378 fMinuitFCN =
new ROOT::Minuit2::FCNGradAdapter<ROOT::Math::IMultiGradFunction> (func, ErrorDef() );
382 const ROOT::Math::FitMethodGradFunction * fcnfunc =
dynamic_cast<const ROOT::Math::FitMethodGradFunction*
>(&func);
384 MN_ERROR_MSG(
"Minuit2Minimizer: Wrong Fit method function for Fumili");
387 fMinuitFCN =
new ROOT::Minuit2::FumiliFCNAdapter<ROOT::Math::FitMethodGradFunction> (*fcnfunc, fDim, ErrorDef() );
391 bool Minuit2Minimizer::Minimize() {
395 MN_ERROR_MSG2(
"Minuit2Minimizer::Minimize",
"FCN function has not been set");
399 assert(GetMinimizer() != 0 );
402 if (fMinimum)
delete fMinimum;
406 int maxfcn = MaxFunctionCalls();
407 double tol = Tolerance();
408 int strategyLevel = Strategy();
409 fMinuitFCN->SetErrorDef(ErrorDef() );
411 int printLevel = PrintLevel();
412 if (printLevel >=1) {
414 int maxfcn_used = maxfcn;
415 if (maxfcn_used == 0) {
416 int nvar = fState.VariableParameters();
417 maxfcn_used = 200 + 100*nvar + 5*nvar*nvar;
419 std::cout <<
"Minuit2Minimizer: Minimize with max-calls " << maxfcn_used
420 <<
" convergence for edm < " << tol <<
" strategy "
421 << strategyLevel << std::endl;
425 MnPrint::SetLevel(printLevel );
426 fMinimizer->Builder().SetPrintLevel(printLevel);
429 int prev_level = (printLevel <= 0 ) ? TurnOffPrintInfoLevel() : -2;
432 if (Precision() > 0) fState.SetPrecision(Precision());
435 ROOT::Minuit2::MnStrategy strategy(strategyLevel);
436 ROOT::Math::IOptions * minuit2Opt = ROOT::Math::MinimizerOptions::FindDefault(
"Minuit2");
439 int nGradCycles = strategy.GradientNCycles();
440 int nHessCycles = strategy.HessianNCycles();
441 int nHessGradCycles = strategy.HessianGradientNCycles();
443 double gradTol = strategy.GradientTolerance();
444 double gradStepTol = strategy.GradientStepTolerance();
445 double hessStepTol = strategy.HessianStepTolerance();
446 double hessG2Tol = strategy.HessianG2Tolerance();
448 minuit2Opt->GetValue(
"GradientNCycles",nGradCycles);
449 minuit2Opt->GetValue(
"HessianNCycles",nHessCycles);
450 minuit2Opt->GetValue(
"HessianGradientNCycles",nHessGradCycles);
452 minuit2Opt->GetValue(
"GradientTolerance",gradTol);
453 minuit2Opt->GetValue(
"GradientStepTolerance",gradStepTol);
454 minuit2Opt->GetValue(
"HessianStepTolerance",hessStepTol);
455 minuit2Opt->GetValue(
"HessianG2Tolerance",hessG2Tol);
457 strategy.SetGradientNCycles(nGradCycles);
458 strategy.SetHessianNCycles(nHessCycles);
459 strategy.SetHessianGradientNCycles(nHessGradCycles);
461 strategy.SetGradientTolerance(gradTol);
462 strategy.SetGradientStepTolerance(gradStepTol);
463 strategy.SetHessianStepTolerance(hessStepTol);
464 strategy.SetHessianG2Tolerance(hessStepTol);
466 int storageLevel = 1;
467 bool ret = minuit2Opt->GetValue(
"StorageLevel",storageLevel);
468 if (ret) SetStorageLevel(storageLevel);
470 if (printLevel > 0) {
471 std::cout <<
"Minuit2Minimizer::Minuit - Changing default options" << std::endl;
480 MnTraceObject * traceObj = 0;
481 #ifdef USE_ROOT_ERROR
482 if (printLevel == 10 && gROOT) {
483 TObject * obj = gROOT->FindObject(
"Minuit2TraceObject");
484 traceObj =
dynamic_cast<ROOT::Minuit2::MnTraceObject*
>(obj);
490 if (printLevel == 20 || printLevel == 30 || printLevel == 40 || (printLevel >= 20000 && printLevel < 30000) ) {
491 int parNumber = printLevel-20000;
492 if (printLevel == 20) parNumber = -1;
493 if (printLevel == 30) parNumber = -2;
494 if (printLevel == 40) parNumber = 0;
495 traceObj =
new TMinuit2TraceObject(parNumber);
498 if (printLevel == 100 || (printLevel >= 10000 && printLevel < 20000)) {
499 int parNumber = printLevel-10000;
500 traceObj =
new MnTraceObject(parNumber);
503 traceObj->Init(fState);
504 SetTraceObject(*traceObj);
507 const ROOT::Minuit2::FCNGradientBase * gradFCN =
dynamic_cast<const ROOT::Minuit2::FCNGradientBase *
>( fMinuitFCN );
511 ROOT::Minuit2::FunctionMinimum min = GetMinimizer()->Minimize(*gradFCN, fState, strategy, maxfcn, tol);
512 fMinimum =
new ROOT::Minuit2::FunctionMinimum (min);
515 ROOT::Minuit2::FunctionMinimum min = GetMinimizer()->Minimize(*GetFCN(), fState, strategy, maxfcn, tol);
516 fMinimum =
new ROOT::Minuit2::FunctionMinimum (min);
520 if (fMinimum->IsValid() && IsValidError() && fMinimum->State().Error().Dcovar() != 0 ) {
522 ROOT::Minuit2::MnHesse hesse(strategy );
523 hesse( *fMinuitFCN, *fMinimum, maxfcn);
527 if (prev_level > -2) RestoreGlobalPrintLevel(prev_level);
529 fState = fMinimum->UserState();
530 bool ok = ExamineMinimum(*fMinimum);
534 if (traceObj) {
delete traceObj; }
538 bool Minuit2Minimizer::ExamineMinimum(
const ROOT::Minuit2::FunctionMinimum & min) {
542 int debugLevel = PrintLevel();
543 if (debugLevel >= 3) {
545 const std::vector<ROOT::Minuit2::MinimumState>& iterationStates = min.States();
546 std::cout <<
"Number of iterations " << iterationStates.size() << std::endl;
547 for (
unsigned int i = 0; i < iterationStates.size(); ++i) {
549 const ROOT::Minuit2::MinimumState & st = iterationStates[i];
550 std::cout <<
"----------> Iteration " << i << std::endl;
551 int pr = std::cout.precision(12);
552 std::cout <<
" FVAL = " << st.Fval() <<
" Edm = " << st.Edm() <<
" Nfcn = " << st.NFcn() << std::endl;
553 std::cout.precision(pr);
554 if (st.HasCovariance() )
555 std::cout <<
" Error matrix change = " << st.Error().Dcovar() << std::endl;
556 if (st.HasParameters() ) {
557 std::cout <<
" Parameters : ";
559 for (
int j = 0; j < st.size() ; ++j) std::cout <<
" p" << j <<
" = " << fState.Int2ext( j, st.Vec()(j) );
560 std::cout << std::endl;
567 if (!min.HasPosDefCovar() ) {
570 txt =
"Covar is not pos def";
573 if (min.HasMadePosDefCovar() ) {
574 txt =
"Covar was made pos def";
577 if (min.HesseFailed() ) {
578 txt =
"Hesse is not valid";
581 if (min.IsAboveMaxEdm() ) {
582 txt =
"Edm is above max";
585 if (min.HasReachedCallLimit() ) {
586 txt =
"Reached call limit";
591 bool validMinimum = min.IsValid();
594 if (fStatus != 0 && debugLevel > 0) MN_INFO_MSG2(
"Minuit2Minimizer::Minimize",txt);
600 txt =
"unknown failure";
603 std::string msg =
"Minimization did NOT converge, " + txt;
604 MN_INFO_MSG2(
"Minuit2Minimizer::Minimize",msg);
607 if (debugLevel >= 1) PrintResults();
612 void Minuit2Minimizer::PrintResults() {
614 if (!fMinimum)
return;
615 if (fMinimum->IsValid() ) {
617 std::cout <<
"Minuit2Minimizer : Valid minimum - status = " << fStatus << std::endl;
618 int pr = std::cout.precision(18);
619 std::cout <<
"FVAL = " << fState.Fval() << std::endl;
620 std::cout <<
"Edm = " << fState.Edm() << std::endl;
621 std::cout.precision(pr);
622 std::cout <<
"Nfcn = " << fState.NFcn() << std::endl;
623 for (
unsigned int i = 0; i < fState.MinuitParameters().size(); ++i) {
624 const MinuitParameter & par = fState.Parameter(i);
625 std::cout << par.Name() <<
"\t = " << par.Value() <<
"\t ";
626 if (par.IsFixed() ) std::cout <<
"(fixed)" << std::endl;
627 else if (par.IsConst() ) std::cout <<
"(const)" << std::endl;
628 else if (par.HasLimits() )
629 std::cout <<
"+/- " << par.Error() <<
"\t(limited)"<< std::endl;
631 std::cout <<
"+/- " << par.Error() << std::endl;
635 std::cout <<
"Minuit2Minimizer : Invalid Minimum - status = " << fStatus << std::endl;
636 std::cout <<
"FVAL = " << fState.Fval() << std::endl;
637 std::cout <<
"Edm = " << fState.Edm() << std::endl;
638 std::cout <<
"Nfcn = " << fState.NFcn() << std::endl;
642 const double * Minuit2Minimizer::X()
const {
644 const std::vector<MinuitParameter> & paramsObj = fState.MinuitParameters();
645 if (paramsObj.size() == 0)
return 0;
646 assert(fDim == paramsObj.size());
649 if (fValues.size() != fDim) fValues.resize(fDim);
650 for (
unsigned int i = 0; i < fDim; ++i) {
651 fValues[i] = paramsObj[i].Value();
654 return &fValues.front();
658 const double * Minuit2Minimizer::Errors()
const {
660 const std::vector<MinuitParameter> & paramsObj = fState.MinuitParameters();
661 if (paramsObj.size() == 0)
return 0;
662 assert(fDim == paramsObj.size());
665 if (fErrors.size() != fDim) fErrors.resize( fDim );
666 for (
unsigned int i = 0; i < fDim; ++i) {
667 const MinuitParameter & par = paramsObj[i];
668 if (par.IsFixed() || par.IsConst() )
671 fErrors[i] = par.Error();
674 return &fErrors.front();
678 double Minuit2Minimizer::CovMatrix(
unsigned int i,
unsigned int j)
const {
680 if ( i >= fDim || j >= fDim)
return 0;
681 if ( !fState.HasCovariance() )
return 0;
682 if (fState.Parameter(i).IsFixed() || fState.Parameter(i).IsConst() )
return 0;
683 if (fState.Parameter(j).IsFixed() || fState.Parameter(j).IsConst() )
return 0;
684 unsigned int k = fState.IntOfExt(i);
685 unsigned int l = fState.IntOfExt(j);
686 return fState.Covariance()(k,l);
689 bool Minuit2Minimizer::GetCovMatrix(
double * cov)
const {
691 if ( !fState.HasCovariance() )
return false;
692 for (
unsigned int i = 0; i < fDim; ++i) {
693 if (fState.Parameter(i).IsFixed() || fState.Parameter(i).IsConst() ) {
694 for (
unsigned int j = 0; j < fDim; ++j) { cov[i*fDim + j] = 0; }
698 unsigned int l = fState.IntOfExt(i);
699 for (
unsigned int j = 0; j < fDim; ++j) {
702 if (fState.Parameter(j).IsFixed() || fState.Parameter(j).IsConst() )
707 unsigned int m = fState.IntOfExt(j);
708 cov[k] = fState.Covariance()(l,m);
716 bool Minuit2Minimizer::GetHessianMatrix(
double * hess)
const {
719 if ( !fState.HasCovariance() )
return false;
720 for (
unsigned int i = 0; i < fDim; ++i) {
721 if (fState.Parameter(i).IsFixed() || fState.Parameter(i).IsConst() ) {
722 for (
unsigned int j = 0; j < fDim; ++j) { hess[i*fDim + j] = 0; }
725 unsigned int l = fState.IntOfExt(i);
726 for (
unsigned int j = 0; j < fDim; ++j) {
729 if (fState.Parameter(j).IsFixed() || fState.Parameter(j).IsConst() )
734 unsigned int m = fState.IntOfExt(j);
735 hess[k] = fState.Hessian()(l,m);
745 double Minuit2Minimizer::Correlation(
unsigned int i,
unsigned int j)
const {
747 if ( i >= fDim || j >= fDim)
return 0;
748 if ( !fState.HasCovariance() )
return 0;
749 if (fState.Parameter(i).IsFixed() || fState.Parameter(i).IsConst() )
return 0;
750 if (fState.Parameter(j).IsFixed() || fState.Parameter(j).IsConst() )
return 0;
751 unsigned int k = fState.IntOfExt(i);
752 unsigned int l = fState.IntOfExt(j);
753 double cij = fState.IntCovariance()(k,l);
754 double tmp = std::sqrt( std::abs ( fState.IntCovariance()(k,k) * fState.IntCovariance()(l,l) ) );
755 if (tmp > 0 )
return cij/tmp;
759 double Minuit2Minimizer::GlobalCC(
unsigned int i)
const {
764 if ( i >= fDim )
return 0;
766 if ( !fState.HasGlobalCC() )
return 0;
767 if (fState.Parameter(i).IsFixed() || fState.Parameter(i).IsConst() )
return 0;
768 unsigned int k = fState.IntOfExt(i);
769 return fState.GlobalCC().GlobalCC()[k];
773 bool Minuit2Minimizer::GetMinosError(
unsigned int i,
double & errLow,
double & errUp,
int runopt) {
778 errLow = 0; errUp = 0;
779 bool runLower = runopt != 2;
780 bool runUpper = runopt != 1;
782 assert( fMinuitFCN );
785 if ( fState.Parameter(i).IsConst() || fState.Parameter(i).IsFixed() ) {
789 int debugLevel = PrintLevel();
791 MnPrint::SetLevel( debugLevel );
799 MN_ERROR_MSG(
"Minuit2Minimizer::GetMinosErrors: failed - no function minimum existing");
803 if (!fMinimum->IsValid() ) {
804 MN_ERROR_MSG(
"Minuit2Minimizer::MINOS failed due to invalid function minimum");
808 fMinuitFCN->SetErrorDef(ErrorDef() );
810 if (ErrorDef() != fMinimum->Up() )
811 fMinimum->SetErrorDef(ErrorDef() );
814 int prev_level = (PrintLevel() <= 0 ) ? TurnOffPrintInfoLevel() : -2;
817 if (Precision() > 0) fState.SetPrecision(Precision());
820 ROOT::Minuit2::MnMinos minos( *fMinuitFCN, *fMinimum);
825 int maxfcn = MaxFunctionCalls();
826 double tol = Tolerance();
828 const char * par_name = fState.Name(i);
833 tol = std::max(tol, 0.01);
835 if (PrintLevel() >=1) {
837 int maxfcn_used = maxfcn;
838 if (maxfcn_used == 0) {
839 int nvar = fState.VariableParameters();
840 maxfcn_used = 2*(nvar+1)*(200 + 100*nvar + 5*nvar*nvar);
842 std::cout <<
"Minuit2Minimizer::GetMinosError for parameter " << i <<
" " << par_name
843 <<
" using max-calls " << maxfcn_used <<
", tolerance " << tol << std::endl;
847 if (runLower) low = minos.Loval(i,maxfcn,tol);
848 if (runUpper) up = minos.Upval(i,maxfcn,tol);
850 ROOT::Minuit2::MinosError me(i, fMinimum->UserState().Value(i),low, up);
852 if (prev_level > -2) RestoreGlobalPrintLevel(prev_level);
861 if (debugLevel >= 1) {
863 if (!me.LowerValid() )
864 std::cout <<
"Minos: Invalid lower error for parameter " << par_name << std::endl;
865 if(me.AtLowerLimit())
866 std::cout <<
"Minos: Parameter : " << par_name <<
" is at Lower limit."<<std::endl;
867 if(me.AtLowerMaxFcn())
868 std::cout <<
"Minos: Maximum number of function calls exceeded when running for lower error" <<std::endl;
869 if(me.LowerNewMin() )
870 std::cout <<
"Minos: New Minimum found while running Minos for lower error" <<std::endl;
872 if (debugLevel > 1) std::cout <<
"Minos: Lower error for parameter " << par_name <<
" : " << me.Lower() << std::endl;
876 if (!me.UpperValid() )
877 std::cout <<
"Minos: Invalid upper error for parameter " << par_name << std::endl;
878 if(me.AtUpperLimit())
879 std::cout <<
"Minos: Parameter " << par_name <<
" is at Upper limit."<<std::endl;
880 if(me.AtUpperMaxFcn())
881 std::cout <<
"Minos: Maximum number of function calls exceeded when running for upper error" <<std::endl;
882 if(me.UpperNewMin() )
883 std::cout <<
"Minos: New Minimum found while running Minos for upper error" <<std::endl;
885 if (debugLevel > 1) std::cout <<
"Minos: Upper error for parameter " << par_name <<
" : " << me.Upper() << std::endl;
890 bool lowerInvalid = (runLower && !me.LowerValid() );
891 bool upperInvalid = (runUpper && !me.UpperValid() );
893 if (lowerInvalid || upperInvalid ) {
901 if (me.AtLowerMaxFcn() ) mstatus |= 4;
902 if (me.LowerNewMin() ) mstatus |= 8;
906 if (me.AtUpperMaxFcn() ) mstatus |= 4;
907 if (me.UpperNewMin() ) mstatus |= 8;
910 fStatus += 10*mstatus;
916 bool isValid = (runLower && me.LowerValid() ) || (runUpper && me.UpperValid() );
920 bool Minuit2Minimizer::Scan(
unsigned int ipar,
unsigned int & nstep,
double * x,
double * y,
double xmin,
double xmax) {
927 MN_ERROR_MSG2(
"Minuit2Minimizer::Scan",
" Function must be set before using Scan");
931 if ( ipar > fState.MinuitParameters().size() ) {
932 MN_ERROR_MSG2(
"Minuit2Minimizer::Scan",
" Invalid number. Minimizer variables must be set before using Scan");
937 int prev_level = (PrintLevel() <= 0 ) ? TurnOffPrintInfoLevel() : -2;
939 MnPrint::SetLevel( PrintLevel() );
943 if (Precision() > 0) fState.SetPrecision(Precision());
945 MnParameterScan scan( *fMinuitFCN, fState.Parameters() );
946 double amin = scan.Fval();
949 std::vector<std::pair<double, double> > result = scan(ipar, nstep-1, xmin, xmax);
951 if (prev_level > -2) RestoreGlobalPrintLevel(prev_level);
953 if (result.size() != nstep) {
954 MN_ERROR_MSG2(
"Minuit2Minimizer::Scan",
" Invalid result from MnParameterScan");
958 std::sort(result.begin(), result.end() );
961 for (
unsigned int i = 0; i < nstep; ++i ) {
962 x[i] = result[i].first;
963 y[i] = result[i].second;
968 if (scan.Fval() < amin ) {
969 MN_INFO_MSG2(
"Minuit2Minimizer::Scan",
"A new minimum has been found");
970 fState.SetValue(ipar, scan.Parameters().Value(ipar) );
978 bool Minuit2Minimizer::Contour(
unsigned int ipar,
unsigned int jpar,
unsigned int & npoints,
double * x,
double * y) {
982 MN_ERROR_MSG2(
"Minuit2Minimizer::Contour",
" no function minimum existing. Must minimize function before");
986 if (!fMinimum->IsValid() ) {
987 MN_ERROR_MSG2(
"Minuit2Minimizer::Contour",
"Invalid function minimum");
992 fMinuitFCN->SetErrorDef(ErrorDef() );
994 if (ErrorDef() != fMinimum->Up() ) {
995 fMinimum->SetErrorDef(ErrorDef() );
998 if ( PrintLevel() >= 1 )
999 MN_INFO_VAL2(
"Minuit2Minimizer::Contour - computing contours - ",ErrorDef());
1002 int prev_level = (PrintLevel() <= 1 ) ? TurnOffPrintInfoLevel() : -2;
1005 MnPrint::SetLevel( PrintLevel() -1 );
1008 if (Precision() > 0) fState.SetPrecision(Precision());
1011 MnContours contour(*fMinuitFCN, *fMinimum, Strategy() );
1013 if (prev_level > -2) RestoreGlobalPrintLevel(prev_level);
1016 std::vector<std::pair<double,double> > result = contour(ipar,jpar, npoints);
1017 if (result.size() != npoints) {
1018 MN_ERROR_MSG2(
"Minuit2Minimizer::Contour",
" Invalid result from MnContours");
1021 for (
unsigned int i = 0; i < npoints; ++i ) {
1022 x[i] = result[i].first;
1023 y[i] = result[i].second;
1027 MnPrint::SetLevel( PrintLevel() );
1035 bool Minuit2Minimizer::Hesse( ) {
1042 MN_ERROR_MSG2(
"Minuit2Minimizer::Hesse",
"FCN function has not been set");
1046 int strategy = Strategy();
1047 int maxfcn = MaxFunctionCalls();
1050 int prev_level = (PrintLevel() <= 0 ) ? TurnOffPrintInfoLevel() : -2;
1052 MnPrint::SetLevel( PrintLevel() );
1055 if (Precision() > 0) fState.SetPrecision(Precision());
1057 ROOT::Minuit2::MnHesse hesse( strategy );
1059 if (PrintLevel() >= 1)
1060 std::cout <<
"Minuit2Minimizer::Hesse using max-calls " << maxfcn << std::endl;
1071 hesse( *fMinuitFCN, *fMinimum, maxfcn );
1073 fState = fMinimum->UserState();
1079 fState = hesse( *fMinuitFCN, fState, maxfcn);
1082 if (prev_level > -2) RestoreGlobalPrintLevel(prev_level);
1084 if (PrintLevel() >= 3) {
1085 std::cout <<
"Minuit2Minimizer::Hesse - State returned from Hesse " << std::endl;
1086 std::cout << fState << std::endl;
1089 int covStatus = fState.CovarianceStatus();
1090 std::string covStatusType =
"not valid";
1091 if (covStatus == 1) covStatusType =
"approximate";
1092 if (covStatus == 2) covStatusType =
"full but made positive defined";
1093 if (covStatus == 3) covStatusType =
"accurate";
1095 if (!fState.HasCovariance() ) {
1101 if (fMinimum->Error().HesseFailed() ) hstatus = 1;
1102 if (fMinimum->Error().InvertFailed() ) hstatus = 2;
1103 else if (!(fMinimum->Error().IsPosDef()) ) hstatus = 3;
1105 if (PrintLevel() > 0) {
1106 std::string msg =
"Hesse failed - matrix is " + covStatusType;
1107 MN_INFO_MSG2(
"Minuit2Minimizer::Hesse",msg);
1108 MN_INFO_VAL2(
"MInuit2Minimizer::Hesse",hstatus);
1110 fStatus += 100*hstatus;
1113 if (PrintLevel() > 0) {
1114 std::string msg =
"Hesse is valid - matrix is " + covStatusType;
1115 MN_INFO_MSG2(
"Minuit2Minimizer::Hesse",msg);
1121 int Minuit2Minimizer::CovMatrixStatus()
const {
1131 if (fMinimum->HasAccurateCovar() )
return 3;
1132 else if (fMinimum->HasMadePosDefCovar() )
return 2;
1133 else if (fMinimum->HasValidCovariance() )
return 1;
1134 else if (fMinimum->HasCovariance() )
return 0;
1139 return fState.CovarianceStatus();
1144 void Minuit2Minimizer::SetTraceObject(MnTraceObject & obj) {
1146 if (!fMinimizer)
return;
1147 fMinimizer->Builder().SetTraceObject(obj);
1150 void Minuit2Minimizer::SetStorageLevel(
int level) {
1152 if (!fMinimizer)
return;
1153 fMinimizer->Builder().SetStorageLevel(level);