46 static ROOT::Math::IMultiGenFunction *&GetGlobalFuncPtr() {
47 TTHREAD_TLS(ROOT::Math::IMultiGenFunction *) fgFunc =
nullptr;
50 TMinuit * TMinuitMinimizer::fgMinuit = 0;
51 bool TMinuitMinimizer::fgUsed = false;
52 bool TMinuitMinimizer::fgUseStaticMinuit = true;
54 ClassImp(TMinuitMinimizer);
57 TMinuitMinimizer::TMinuitMinimizer(ROOT::Minuit::EMinimizerType type,
unsigned int ndim ) :
70 if (fDim > 0) InitTMinuit(fDim);
73 TMinuitMinimizer::TMinuitMinimizer(
const char * type,
unsigned int ndim ) :
84 std::string algoname(type);
85 std::transform(algoname.begin(), algoname.end(), algoname.begin(), (int(*)(int)) tolower );
87 ROOT::Minuit::EMinimizerType algoType = ROOT::Minuit::kMigrad;
88 if (algoname ==
"simplex") algoType = ROOT::Minuit::kSimplex;
89 if (algoname ==
"minimize" ) algoType = ROOT::Minuit::kCombined;
90 if (algoname ==
"migradimproved" ) algoType = ROOT::Minuit::kMigradImproved;
91 if (algoname ==
"scan" ) algoType = ROOT::Minuit::kScan;
92 if (algoname ==
"seek" ) algoType = ROOT::Minuit::kSeek;
97 if (fDim > 0) InitTMinuit(fDim);
101 TMinuitMinimizer::~TMinuitMinimizer()
104 if (fMinuit && !fgUseStaticMinuit) {
110 TMinuitMinimizer::TMinuitMinimizer(
const TMinuitMinimizer &) :
116 TMinuitMinimizer & TMinuitMinimizer::operator = (
const TMinuitMinimizer &rhs)
119 if (
this == &rhs)
return *
this;
123 bool TMinuitMinimizer::UseStaticMinuit(
bool on ) {
125 bool prev = fgUseStaticMinuit;
126 fgUseStaticMinuit = on;
130 void TMinuitMinimizer::InitTMinuit(
int dim) {
134 if (fMinuit ==0 || dim > fMinuit->fMaxpar) {
137 if (fgUseStaticMinuit) {
145 if (fgMinuit != gMinuit) {
148 if (gROOT->GetListOfSpecials()->FindObject(fgMinuit) == 0) {
170 fgMinuit =
new TMinuit(dim);
172 else if (fgMinuit->GetNumPars() != int(dim) ) {
175 fgMinuit =
new TMinuit(dim);
183 if (fMinuit)
delete fMinuit;
184 fMinuit =
new TMinuit(dim);
199 arglist[0] = PrintLevel() - 1;
200 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
201 if (PrintLevel() <= 0) SuppressMinuitWarnings();
205 void TMinuitMinimizer::SetFunction(
const ROOT::Math::IMultiGenFunction & func) {
218 GetGlobalFuncPtr() =
const_cast<ROOT::Math::IMultiGenFunction *
>(&func);
219 fMinuit->SetFCN(&TMinuitMinimizer::Fcn);
224 fMinuit->mnexcm(
"SET NOGrad",arglist,0,ierr);
227 void TMinuitMinimizer::SetFunction(
const ROOT::Math::IMultiGradFunction & func) {
238 GetGlobalFuncPtr() =
const_cast<ROOT::Math::IMultiGradFunction *
>(&func);
239 fMinuit->SetFCN(&TMinuitMinimizer::FcnGrad);
247 fMinuit->mnexcm(
"SET GRAD",arglist,1,ierr);
250 void TMinuitMinimizer::Fcn(
int &,
double * ,
double & f,
double * x ,
int ) {
253 f = GetGlobalFuncPtr()->operator()(x);
256 void TMinuitMinimizer::FcnGrad(
int &,
double * g,
double & f,
double * x ,
int iflag ) {
260 ROOT::Math::IMultiGradFunction * gFunc =
dynamic_cast<ROOT::Math::IMultiGradFunction *
> ( GetGlobalFuncPtr());
263 f = gFunc->operator()(x);
266 if (iflag == 2) gFunc->Gradient(x,g);
269 bool TMinuitMinimizer::SetVariable(
unsigned int ivar,
const std::string & name,
double val,
double step) {
271 if (!CheckMinuitInstance())
return false;
273 if (fgUseStaticMinuit) fUsed = fgUsed;
276 if (fUsed) DoClear();
279 DoReleaseFixParameter(ivar);
281 int iret = fMinuit->DefineParameter(ivar , name.c_str(), val, step, 0., 0. );
285 bool TMinuitMinimizer::SetLimitedVariable(
unsigned int ivar,
const std::string & name,
double val,
double step,
double lower,
double upper) {
287 if (!CheckMinuitInstance())
return false;
289 if (fgUseStaticMinuit) fUsed = fgUsed;
292 if (fUsed) DoClear();
295 DoReleaseFixParameter(ivar);
297 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, lower, upper );
301 bool TMinuitMinimizer::SetLowerLimitedVariable(
unsigned int ivar ,
const std::string & name ,
double val ,
double step ,
double lower ) {
304 Warning(
"TMinuitMinimizer::SetLowerLimitedVariable",
"not implemented - use as upper limit 1.E7 instead of +inf");
305 return SetLimitedVariable(ivar, name, val , step, lower, lower+ 1.E7);
308 bool TMinuitMinimizer::SetUpperLimitedVariable(
unsigned int ivar ,
const std::string & name ,
double val ,
double step ,
double upper ) {
311 Warning(
"TMinuitMinimizer::SetUpperLimitedVariable",
"not implemented - use as lower limit -1.E7 instead of -inf");
312 return SetLimitedVariable(ivar, name, val , step, upper -1.E7, upper);
316 bool TMinuitMinimizer::CheckMinuitInstance()
const {
319 Error(
"TMinuitMinimizer::CheckMinuitInstance",
"Invalid TMinuit pointer. Need to call first SetFunction");
325 bool TMinuitMinimizer::CheckVarIndex(
unsigned int ivar)
const {
327 if ((
int) ivar >= fMinuit->fNu ) {
328 Error(
"TMinuitMinimizer::CheckVarIndex",
"Invalid parameter index");
335 bool TMinuitMinimizer::SetFixedVariable(
unsigned int ivar,
const std::string & name,
double val) {
337 if (!CheckMinuitInstance())
return false;
340 if (fgUseStaticMinuit) fUsed = fgUsed;
343 if (fUsed) DoClear();
348 double step = ( val != 0) ? 0.1 * std::abs(val) : 0.1;
349 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, 0., 0. );
350 if (iret == 0) iret = fMinuit->FixParameter(ivar);
354 bool TMinuitMinimizer::SetVariableValue(
unsigned int ivar,
double val) {
357 if (!CheckMinuitInstance())
return false;
364 fMinuit->mnexcm(
"SET PAR",arglist,2,ierr);
368 bool TMinuitMinimizer::SetVariableStepSize(
unsigned int ivar,
double step) {
371 if (!CheckMinuitInstance())
return false;
374 double curval,err, lowlim, uplim;
377 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
378 if (iuint == -1)
return false;
379 int iret = fMinuit->DefineParameter(ivar, name, curval, step, lowlim, uplim );
384 bool TMinuitMinimizer::SetVariableLowerLimit(
unsigned int ivar,
double lower ) {
387 Warning(
"TMinuitMinimizer::SetVariableLowerLimit",
"not implemented - use as upper limit 1.E30 instead of +inf");
388 return SetVariableLimits(ivar, lower, 1.E30);
390 bool TMinuitMinimizer::SetVariableUpperLimit(
unsigned int ivar,
double upper ) {
393 Warning(
"TMinuitMinimizer::SetVariableUpperLimit",
"not implemented - - use as lower limit -1.E30 instead of +inf");
394 return SetVariableLimits(ivar, -1.E30, upper);
397 bool TMinuitMinimizer::SetVariableLimits(
unsigned int ivar,
double lower,
double upper) {
401 if (!CheckMinuitInstance())
return false;
404 double curval,err, lowlim, uplim;
407 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
408 if (iuint == -1)
return false;
409 int iret = fMinuit->DefineParameter(ivar, name, curval, err, lower, upper );
414 bool TMinuitMinimizer::FixVariable(
unsigned int ivar) {
416 if (!CheckMinuitInstance())
return false;
417 if (!CheckVarIndex(ivar))
return false;
418 int iret = fMinuit->FixParameter(ivar);
422 bool TMinuitMinimizer::ReleaseVariable(
unsigned int ivar) {
424 if (!CheckMinuitInstance())
return false;
425 if (!CheckVarIndex(ivar))
return false;
426 int iret = fMinuit->Release(ivar);
430 bool TMinuitMinimizer::IsFixedVariable(
unsigned int ivar)
const {
432 if (!CheckMinuitInstance())
return false;
433 if (!CheckVarIndex(ivar))
return false;
434 return (fMinuit->fNiofex[ivar] == 0 );
437 bool TMinuitMinimizer::GetVariableSettings(
unsigned int ivar, ROOT::Fit::ParameterSettings & var)
const {
439 if (!CheckMinuitInstance())
return false;
440 if (!CheckVarIndex(ivar))
return false;
441 double curval,err, lowlim, uplim;
444 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
445 if (iuint == -1)
return false;
446 var.Set(name.Data(), curval, err, lowlim, uplim);
447 if (IsFixedVariable(ivar)) var.Fix();
453 std::string TMinuitMinimizer::VariableName(
unsigned int ivar)
const {
455 if (!CheckMinuitInstance())
return std::string();
456 if (!CheckVarIndex(ivar))
return std::string();
457 return fMinuit->fCpnam[ivar].Data();
460 int TMinuitMinimizer::VariableIndex(
const std::string & )
const {
462 Error(
"TMinuitMinimizer::VariableIndex",
" find index of a variable from its name is not implemented in TMinuit");
466 bool TMinuitMinimizer::Minimize() {
476 Error(
"TMinuitMinimizer::Minimize",
"invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
482 if (fMinuit->fNu < static_cast<int>(fDim) ) {
483 Error(
"TMinuitMinimizer::Minimize",
"The total number of defined parameters is different than the function dimension, npar = %d, dim = %d",fMinuit->fNu, fDim);
487 int printlevel = PrintLevel();
490 if (fMinuit->fNpar <= 0) {
493 fMinuit->fAmin = (*GetGlobalFuncPtr())(&fParams.front());
494 if (printlevel > 0) Info(
"TMinuitMinimizer::Minimize",
"There are no free parameter - just compute the function value");
504 arglist[0] = ErrorDef();
505 fMinuit->mnexcm(
"SET Err",arglist,1,ierr);
507 arglist[0] = printlevel - 1;
508 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
511 if (printlevel == 0) fMinuit->mnexcm(
"SET NOW",arglist,0,ierr);
514 if (Precision() > 0) {
515 arglist[0] = Precision();
516 fMinuit->mnexcm(
"SET EPS",arglist,1,ierr);
520 int strategy = Strategy();
521 if (strategy >=0 && strategy <=2 ) {
522 arglist[0] = strategy;
523 fMinuit->mnexcm(
"SET STR",arglist,1,ierr);
526 arglist[0] = MaxFunctionCalls();
527 arglist[1] = Tolerance();
532 case ROOT::Minuit::kMigrad:
534 fMinuit->mnexcm(
"MIGRAD",arglist,nargs,ierr);
536 case ROOT::Minuit::kCombined:
538 fMinuit->mnexcm(
"MINIMIZE",arglist,nargs,ierr);
540 case ROOT::Minuit::kSimplex:
542 fMinuit->mnexcm(
"SIMPLEX",arglist,nargs,ierr);
544 case ROOT::Minuit::kScan:
547 fMinuit->mnexcm(
"SCAN",arglist,nargs,ierr);
549 case ROOT::Minuit::kSeek:
553 if (arglist[1] >= 1.) nargs = 2;
554 fMinuit->mnexcm(
"SEEK",arglist,nargs,ierr);
558 fMinuit->mnexcm(
"MIGRAD",arglist,nargs,ierr);
566 int minErrStatus = ierr;
568 if (printlevel>2) Info(
"TMinuitMinimizer::Minimize",
"Finished to run MIGRAD - status %d",ierr);
571 if (ierr == 0 && fType == ROOT::Minuit::kMigradImproved) {
572 fMinuit->mnexcm(
"IMPROVE",arglist,1,ierr);
573 fStatus += 1000*ierr;
574 if (printlevel>2) Info(
"TMinuitMinimizer::Minimize",
"Finished to run IMPROVE - status %d",ierr);
581 if (minErrStatus == 0 && (IsValidError() || ( strategy >=1 && CovMatrixStatus() < 3) ) ) {
582 fMinuit->mnexcm(
"HESSE",arglist,1,ierr);
584 if (printlevel>2) Info(
"TMinuitMinimizer::Minimize",
"Finished to run HESSE - status %d",ierr);
590 if (minErrStatus == 0) {
594 RetrieveErrorMatrix();
606 void TMinuitMinimizer::RetrieveParams() {
610 assert(fMinuit != 0);
613 if (fParams.size() != fDim) fParams.resize( fDim);
614 if (fErrors.size() != fDim) fErrors.resize( fDim);
615 for (
unsigned int i = 0; i < fDim; ++i) {
616 fMinuit->GetParameter( i, fParams[i], fErrors[i]);
620 void TMinuitMinimizer::RetrieveErrorMatrix() {
624 assert(fMinuit != 0);
626 unsigned int nfree = NFree();
628 unsigned int ndim2 = fDim*fDim;
629 if (fCovar.size() != ndim2 ) fCovar.resize(fDim*fDim);
631 fMinuit->mnemat(&fCovar.front(), fDim);
635 std::vector<double> tmpMat(nfree*nfree);
636 fMinuit->mnemat(&tmpMat.front(), nfree);
639 for (
unsigned int i = 0; i < fDim; ++i) {
641 if ( fMinuit->fNiofex[i] > 0 ) {
643 for (
unsigned int j = 0; j <= i; ++j) {
644 if ( fMinuit->fNiofex[j] > 0 ) {
645 fCovar[i*fDim + j] = tmpMat[l*nfree + m];
646 fCovar[j*fDim + i] = fCovar[i*fDim + j];
657 unsigned int TMinuitMinimizer::NCalls()
const {
659 if (fMinuit == 0)
return 0;
660 return fMinuit->fNfcn;
663 double TMinuitMinimizer::MinValue()
const {
667 if (!fMinuit)
return 0;
668 double minval = fMinuit->fAmin;
669 if (minval == fMinuit->fUndefi)
return 0;
673 double TMinuitMinimizer::Edm()
const {
677 if (!fMinuit)
return -1;
678 if (fMinuit->fAmin == fMinuit->fUndefi || fMinuit->fEDM == fMinuit->fBigedm)
return fMinuit->fUp;
679 return fMinuit->fEDM;
682 unsigned int TMinuitMinimizer::NFree()
const {
684 if (!fMinuit)
return 0;
685 if (fMinuit->fNpar < 0)
return 0;
686 return fMinuit->fNpar;
689 bool TMinuitMinimizer::GetCovMatrix(
double * cov)
const {
691 int covStatus = CovMatrixStatus();
692 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
693 Error(
"TMinuitMinimizer::GetHessianMatrix",
"Hessian matrix has not been computed - status %d",covStatus);
696 std::copy(fCovar.begin(), fCovar.end(), cov);
697 TMatrixDSym cmat(fDim,cov);
701 bool TMinuitMinimizer::GetHessianMatrix(
double * hes)
const {
705 int covStatus = CovMatrixStatus();
706 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
707 Error(
"TMinuitMinimizer::GetHessianMatrix",
"Hessian matrix has not been computed - status %d",covStatus);
711 unsigned int nfree = NFree();
712 TMatrixDSym mat(nfree);
713 fMinuit->mnemat(mat.GetMatrixArray(), nfree);
719 for (
unsigned int i = 0; i < fDim; ++i) {
720 if ( fMinuit->fNiofex[i] > 0 ) {
722 for (
unsigned int j = 0; j <= i; ++j) {
723 if ( fMinuit->fNiofex[j] > 0 ) {
724 hes[i*fDim + j] = mat(l,m);
725 hes[j*fDim + i] = hes[i*fDim + j];
740 int TMinuitMinimizer::CovMatrixStatus()
const {
748 if (!fMinuit)
return 0;
749 if (fMinuit->fAmin == fMinuit->fUndefi)
return 0;
750 return fMinuit->fISW[1];
753 double TMinuitMinimizer::GlobalCC(
unsigned int i)
const {
755 if (!fMinuit)
return 0;
756 if (!fMinuit->fGlobcc)
return 0;
757 if (
int(i) >= fMinuit->fNu)
return 0;
759 int iin = fMinuit->fNiofex[i];
761 if (iin < 1)
return 0;
762 return fMinuit->fGlobcc[iin-1];
765 bool TMinuitMinimizer::GetMinosError(
unsigned int i,
double & errLow,
double & errUp,
int ) {
769 Error(
"TMinuitMinimizer::GetMinosError",
"invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
774 if (fMinuit->fNiofex[i] == 0 ) {
775 if (PrintLevel() > 0) Info(
"TMinuitMinimizer::GetMinosError",
"Parameter %s is fixed. There are no Minos error to calculate. Ignored.",VariableName(i).c_str());
776 errLow = 0; errUp = 0;
785 if (fMinuit->fUp != ErrorDef() ) {
786 arglist[0] = ErrorDef();
787 fMinuit->mnexcm(
"SET Err",arglist,1,ierr);
790 if (fMinuit->fISW[4] != (PrintLevel()-1) ) {
791 arglist[0] = PrintLevel()-1;
792 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
794 if (PrintLevel() == 0) fMinuit->mnexcm(
"SET NOW",arglist,0,ierr);
796 if (fMinuit->fIstrat != Strategy() ) {
797 arglist[0] = Strategy();
798 fMinuit->mnexcm(
"SET STR",arglist,1,ierr);
801 if (Precision() > 0 && fMinuit->fEpsma2 != Precision() ) {
802 arglist[0] = Precision();
803 fMinuit->mnexcm(
"SET EPS",arglist,1,ierr);
809 arglist[0] = MaxFunctionCalls();
813 fMinuit->mnexcm(
"MINOS",arglist,nargs,ierr);
814 bool isValid = (ierr == 0);
816 if (isValid && fMinuit->fCstatu !=
"SUCCESSFUL") {
817 if (fMinuit->fCstatu ==
"FAILURE" ) {
822 if (fMinuit->fCstatu ==
"PROBLEMS") ierr = 6;
833 fMinuit->mnerrs(i,errUp,errLow, errParab, gcor);
840 void TMinuitMinimizer::DoClear() {
848 fMinuit->mnrn15(val,inseed);
855 void TMinuitMinimizer::DoReleaseFixParameter(
int ivar) {
859 if (fMinuit == 0)
return;
860 if (fMinuit->GetNumFixedPars() == 0)
return;
862 if (
int(ivar) >= fMinuit->GetNumPars() )
return;
865 for (
int i = 0; i < fMinuit->fNpfix; ++i) {
866 if (fMinuit->fIpfix[i] == ivar+1 ) {
868 fMinuit->Release(ivar);
876 void TMinuitMinimizer::PrintResults() {
878 if (fMinuit == 0)
return;
881 if (PrintLevel() > 2)
882 fMinuit->mnprin(4,fMinuit->fAmin);
884 fMinuit->mnprin(3,fMinuit->fAmin);
887 void TMinuitMinimizer::SuppressMinuitWarnings(
bool nowarn) {
892 fMinuit->mnexcm(
"SET NOW",&arglist,0,ierr);
894 fMinuit->mnexcm(
"SET WAR",&arglist,0,ierr);
898 bool TMinuitMinimizer::Contour(
unsigned int ipar,
unsigned int jpar,
unsigned int &npoints,
double * x,
double * y) {
902 Error(
"TMinuitMinimizer::Contour",
" invalid TMinuit instance");
909 arglist[0] = ErrorDef();
910 fMinuit->mnexcm(
"SET Err",arglist,1,ierr);
912 arglist[0] = PrintLevel()-1;
913 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
916 if (PrintLevel() == 0) fMinuit->mnexcm(
"SET NOW",arglist,0,ierr);
919 if (Precision() > 0) {
920 arglist[0] = Precision();
921 fMinuit->mnexcm(
"SET EPS",arglist,1,ierr);
926 Error(
"TMinuitMinimizer::Contour",
"Cannot make contour with so few points");
931 fMinuit->mncont( ipar,jpar,npoints, x, y,npfound);
934 Error(
"TMinuitMinimizer::Contour",
"Cannot find more than 4 points");
937 if (npfound!=(
int)npoints) {
939 Warning(
"TMinuitMinimizer::Contour",
"Returning only %d points ",npfound);
946 bool TMinuitMinimizer::Scan(
unsigned int ipar,
unsigned int & nstep,
double * x,
double * y,
double xmin,
double xmax) {
956 Error(
"TMinuitMinimizer::Scan",
" invalid TMinuit instance");
961 if (xmin >= xmax && (
int) ipar < fMinuit->GetNumPars() ) {
962 double val = 0;
double err = 0;
964 double xlow = 0;
double xup = 0 ;
967 fMinuit->mnpout( ipar, name, val, err, xlow, xup, iuint);
969 if (iuint > 0 && err > 0) {
971 xmax = val + 2 * err;
978 arglist[0] = PrintLevel()-1;
979 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
981 if (PrintLevel() == 0) fMinuit->mnexcm(
"SET NOW",arglist,0,ierr);
984 if (Precision() > 0) {
985 arglist[0] = Precision();
986 fMinuit->mnexcm(
"SET EPS",arglist,1,ierr);
989 if (nstep == 0)
return false;
998 fMinuit->mnexcm(
"SCAN",arglist,nargs,ierr);
1000 Error(
"TMinuitMinimizer::Scan",
" Error executing command SCAN");
1004 TGraph * gr =
dynamic_cast<TGraph *
>(fMinuit->GetPlot() );
1006 Error(
"TMinuitMinimizer::Scan",
" Error in returned graph object");
1009 nstep = std::min(gr->GetN(), (int) nstep);
1012 std::copy(gr->GetX(), gr->GetX()+nstep, x);
1013 std::copy(gr->GetY(), gr->GetY()+nstep, y);
1018 bool TMinuitMinimizer::Hesse() {
1022 Error(
"TMinuitMinimizer::Hesse",
"invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1031 arglist[0] = ErrorDef();
1032 fMinuit->mnexcm(
"SET ERR",arglist,1,ierr);
1034 int printlevel = PrintLevel();
1035 arglist[0] = printlevel - 1;
1036 fMinuit->mnexcm(
"SET PRINT",arglist,1,ierr);
1039 if (printlevel == 0) fMinuit->mnexcm(
"SET NOW",arglist,0,ierr);
1042 if (Precision() > 0) {
1043 arglist[0] = Precision();
1044 fMinuit->mnexcm(
"SET EPS",arglist,1,ierr);
1047 arglist[0] = MaxFunctionCalls();
1049 fMinuit->mnexcm(
"HESSE",arglist,1,ierr);
1050 fStatus += 100*ierr;
1052 if (ierr != 0)
return false;
1058 RetrieveErrorMatrix();
1063 bool TMinuitMinimizer::SetDebug(
bool on) {
1067 Error(
"TMinuitMinimizer::SetDebug",
"invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1074 fMinuit->mnexcm(
"SET DEBUG",arglist,1,ierr);
1076 fMinuit->mnexcm(
"SET NODEBUG",arglist,1,ierr);