39 #ifndef __ROOFIT_NOROOMINIMIZER
73 #if (__GNUC__==3&&__GNUC_MINOR__==2&&__GNUC_PATCHLEVEL__==3)
74 char* operator+( streampos&,
char* );
79 ClassImp(RooMinimizer);
82 ROOT::Fit::Fitter *RooMinimizer::_theFitter = 0 ;
90 void RooMinimizer::cleanup()
112 RooMinimizer::RooMinimizer(RooAbsReal&
function)
114 RooSentinel::activate() ;
122 _profileStart = kFALSE ;
124 _minimizerType =
"Minuit";
126 if (_theFitter)
delete _theFitter ;
127 _theFitter =
new ROOT::Fit::Fitter;
128 _fcn =
new RooMinimizerFcn(_func,
this,_verbose);
129 _theFitter->Config().SetMinimizer(_minimizerType.c_str());
132 _theFitter->Config().MinimizerOptions().SetMaxIterations(500*_fcn->NDim());
133 _theFitter->Config().MinimizerOptions().SetMaxFunctionCalls(500*_fcn->NDim());
139 setErrorLevel(_func->defaultErrorLevel()) ;
142 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
143 _optConst,_verbose) ;
146 if (RooMsgService::instance().silentMode()) {
158 RooMinimizer::~RooMinimizer()
178 void RooMinimizer::setStrategy(Int_t istrat)
180 _theFitter->Config().MinimizerOptions().SetStrategy(istrat);
190 void RooMinimizer::setMaxIterations(Int_t n)
192 _theFitter->Config().MinimizerOptions().SetMaxIterations(n);
202 void RooMinimizer::setMaxFunctionCalls(Int_t n)
204 _theFitter->Config().MinimizerOptions().SetMaxFunctionCalls(n);
216 void RooMinimizer::setErrorLevel(Double_t level)
218 _theFitter->Config().MinimizerOptions().SetErrorDef(level);
227 void RooMinimizer::setEps(Double_t eps)
229 _theFitter->Config().MinimizerOptions().SetTolerance(eps);
237 void RooMinimizer::setOffsetting(Bool_t flag)
239 _func->enableOffsetting(flag) ;
248 void RooMinimizer::setMinimizerType(
const char* type)
250 _minimizerType = type;
259 ROOT::Fit::Fitter* RooMinimizer::fitter()
268 const ROOT::Fit::Fitter* RooMinimizer::fitter()
const
286 RooFitResult* RooMinimizer::fit(
const char* options)
288 TString opts(options) ;
292 if (opts.Contains(
"v")) setVerbose(1) ;
293 if (opts.Contains(
"t")) setProfile(1) ;
294 if (opts.Contains(
"l")) setLogFile(Form(
"%s.log",_func->GetName())) ;
295 if (opts.Contains(
"c")) optimizeConst(1) ;
298 if (opts.Contains(
"0")) setStrategy(0) ;
300 if (opts.Contains(
"0")) setStrategy(1) ;
301 if (opts.Contains(
"h")||!opts.Contains(
"m")) hesse() ;
302 if (!opts.Contains(
"m")) minos() ;
304 return (opts.Contains(
"r")) ? save() : 0 ;
312 Int_t RooMinimizer::minimize(
const char* type,
const char* alg)
314 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
315 _optConst,_verbose) ;
317 _theFitter->Config().SetMinimizer(type,alg);
320 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
321 RooAbsReal::clearEvalErrorLog() ;
323 bool ret = _theFitter->FitFCN(*_fcn);
324 _status = ((ret) ? _theFitter->Result().Status() : -1);
326 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
328 _fcn->BackProp(_theFitter->Result());
330 saveStatus(
"MINIMIZE",_status) ;
343 Int_t RooMinimizer::migrad()
345 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
346 _optConst,_verbose) ;
348 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
349 RooAbsReal::clearEvalErrorLog() ;
351 _theFitter->Config().SetMinimizer(_minimizerType.c_str(),
"migrad");
352 bool ret = _theFitter->FitFCN(*_fcn);
353 _status = ((ret) ? _theFitter->Result().Status() : -1);
355 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
357 _fcn->BackProp(_theFitter->Result());
359 saveStatus(
"MIGRAD",_status) ;
372 Int_t RooMinimizer::hesse()
374 if (_theFitter->GetMinimizer()==0) {
375 coutW(Minimization) <<
"RooMinimizer::hesse: Error, run Migrad before Hesse!"
381 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
382 _optConst,_verbose) ;
384 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
385 RooAbsReal::clearEvalErrorLog() ;
387 _theFitter->Config().SetMinimizer(_minimizerType.c_str());
388 bool ret = _theFitter->CalculateHessErrors();
389 _status = ((ret) ? _theFitter->Result().Status() : -1);
391 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
393 _fcn->BackProp(_theFitter->Result());
395 saveStatus(
"HESSE",_status) ;
409 Int_t RooMinimizer::minos()
411 if (_theFitter->GetMinimizer()==0) {
412 coutW(Minimization) <<
"RooMinimizer::minos: Error, run Migrad before Minos!"
418 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
419 _optConst,_verbose) ;
421 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
422 RooAbsReal::clearEvalErrorLog() ;
424 _theFitter->Config().SetMinimizer(_minimizerType.c_str());
425 bool ret = _theFitter->CalculateMinosErrors();
426 _status = ((ret) ? _theFitter->Result().Status() : -1);
428 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
430 _fcn->BackProp(_theFitter->Result());
432 saveStatus(
"MINOS",_status) ;
447 Int_t RooMinimizer::minos(
const RooArgSet& minosParamList)
449 if (_theFitter->GetMinimizer()==0) {
450 coutW(Minimization) <<
"RooMinimizer::minos: Error, run Migrad before Minos!"
454 else if (minosParamList.getSize()>0) {
456 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
457 _optConst,_verbose) ;
459 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
460 RooAbsReal::clearEvalErrorLog() ;
463 TIterator* aIter = minosParamList.createIterator() ;
465 std::vector<unsigned int> paramInd;
466 while((arg=(RooAbsArg*)aIter->Next())) {
467 RooAbsArg* par = _fcn->GetFloatParamList()->find(arg->GetName());
468 if (par && !par->isConstant()) {
469 Int_t index = _fcn->GetFloatParamList()->index(par);
470 paramInd.push_back(index);
475 if (paramInd.size()) {
477 _theFitter->Config().SetMinosErrors(paramInd);
479 _theFitter->Config().SetMinimizer(_minimizerType.c_str());
480 bool ret = _theFitter->CalculateMinosErrors();
481 _status = ((ret) ? _theFitter->Result().Status() : -1);
483 _theFitter->Config().SetMinosErrors(kFALSE);
487 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
489 _fcn->BackProp(_theFitter->Result());
491 saveStatus(
"MINOS",_status) ;
506 Int_t RooMinimizer::seek()
508 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
509 _optConst,_verbose) ;
511 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
512 RooAbsReal::clearEvalErrorLog() ;
514 _theFitter->Config().SetMinimizer(_minimizerType.c_str(),
"seek");
515 bool ret = _theFitter->FitFCN(*_fcn);
516 _status = ((ret) ? _theFitter->Result().Status() : -1);
518 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
520 _fcn->BackProp(_theFitter->Result());
522 saveStatus(
"SEEK",_status) ;
535 Int_t RooMinimizer::simplex()
537 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
538 _optConst,_verbose) ;
540 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
541 RooAbsReal::clearEvalErrorLog() ;
543 _theFitter->Config().SetMinimizer(_minimizerType.c_str(),
"simplex");
544 bool ret = _theFitter->FitFCN(*_fcn);
545 _status = ((ret) ? _theFitter->Result().Status() : -1);
547 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
549 _fcn->BackProp(_theFitter->Result());
551 saveStatus(
"SEEK",_status) ;
564 Int_t RooMinimizer::improve()
566 _fcn->Synchronize(_theFitter->Config().ParamsSettings(),
567 _optConst,_verbose) ;
569 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
570 RooAbsReal::clearEvalErrorLog() ;
572 _theFitter->Config().SetMinimizer(_minimizerType.c_str(),
"migradimproved");
573 bool ret = _theFitter->FitFCN(*_fcn);
574 _status = ((ret) ? _theFitter->Result().Status() : -1);
576 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
578 _fcn->BackProp(_theFitter->Result());
580 saveStatus(
"IMPROVE",_status) ;
590 Int_t RooMinimizer::setPrintLevel(Int_t newLevel)
592 Int_t ret = _printLevel ;
593 _theFitter->Config().MinimizerOptions().SetPrintLevel(newLevel+1);
594 _printLevel = newLevel+1 ;
602 void RooMinimizer::optimizeConst(Int_t flag)
604 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
606 if (_optConst && !flag){
607 if (_printLevel>-1) coutI(Minimization) <<
"RooMinimizer::optimizeConst: deactivating const optimization" << endl ;
608 _func->constOptimizeTestStatistic(RooAbsArg::DeActivate) ;
610 }
else if (!_optConst && flag) {
611 if (_printLevel>-1) coutI(Minimization) <<
"RooMinimizer::optimizeConst: activating const optimization" << endl ;
612 _func->constOptimizeTestStatistic(RooAbsArg::Activate,flag>1) ;
614 }
else if (_optConst && flag) {
615 if (_printLevel>-1) coutI(Minimization) <<
"RooMinimizer::optimizeConst: const optimization already active" << endl ;
617 if (_printLevel>-1) coutI(Minimization) <<
"RooMinimizer::optimizeConst: const optimization wasn't active" << endl ;
620 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
634 RooFitResult* RooMinimizer::save(
const char* userName,
const char* userTitle)
636 if (_theFitter->GetMinimizer()==0) {
637 coutW(Minimization) <<
"RooMinimizer::save: Error, run minimization before!"
643 name = userName ? userName : Form(
"%s", _func->GetName()) ;
644 title = userTitle ? userTitle : Form(
"%s", _func->GetTitle()) ;
645 RooFitResult* fitRes =
new RooFitResult(name,title) ;
649 RooArgList saveConstList(*(_fcn->GetConstParamList())) ;
650 RooArgList saveFloatInitList(*(_fcn->GetInitFloatParamList())) ;
651 RooArgList saveFloatFinalList(*(_fcn->GetFloatParamList())) ;
652 for (i=0 ; i<_fcn->GetFloatParamList()->getSize() ; i++) {
653 RooAbsArg* par = _fcn->GetFloatParamList()->at(i) ;
654 if (par->isConstant()) {
655 saveFloatInitList.remove(*saveFloatInitList.find(par->GetName()),kTRUE) ;
656 saveFloatFinalList.remove(*par) ;
657 saveConstList.add(*par) ;
660 saveConstList.sort() ;
662 fitRes->setConstParList(saveConstList) ;
663 fitRes->setInitParList(saveFloatInitList) ;
665 fitRes->setStatus(_status) ;
666 fitRes->setCovQual(_theFitter->GetMinimizer()->CovMatrixStatus()) ;
667 fitRes->setMinNLL(_theFitter->Result().MinFcnValue()) ;
668 fitRes->setNumInvalidNLL(_fcn->GetNumInvalidNLL()) ;
669 fitRes->setEDM(_theFitter->Result().Edm()) ;
670 fitRes->setFinalParList(saveFloatFinalList) ;
672 std::vector<double> globalCC;
673 TMatrixDSym corrs(_theFitter->Result().Parameters().size()) ;
674 TMatrixDSym covs(_theFitter->Result().Parameters().size()) ;
675 for (UInt_t ic=0; ic<_theFitter->Result().Parameters().size(); ic++) {
676 globalCC.push_back(_theFitter->Result().GlobalCC(ic));
677 for (UInt_t ii=0; ii<_theFitter->Result().Parameters().size(); ii++) {
678 corrs(ic,ii) = _theFitter->Result().Correlation(ic,ii);
679 covs(ic,ii) = _theFitter->Result().CovMatrix(ic,ii);
682 fitRes->fillCorrMatrix(globalCC,corrs,covs) ;
684 fitRes->setCovarianceMatrix(*_extV) ;
687 fitRes->setStatusHistory(_statusHistory) ;
709 RooPlot* RooMinimizer::contour(RooRealVar& var1, RooRealVar& var2,
710 Double_t n1, Double_t n2, Double_t n3,
711 Double_t n4, Double_t n5, Double_t n6,
unsigned int npoints)
715 RooArgList* params = _fcn->GetFloatParamList() ;
716 RooArgList* paramSave = (RooArgList*) params->snapshot() ;
719 Int_t index1= _fcn->GetFloatParamList()->index(&var1);
721 coutE(Minimization) <<
"RooMinimizer::contour(" << GetName()
722 <<
") ERROR: " << var1.GetName()
723 <<
" is not a floating parameter of "
724 << _func->GetName() << endl ;
728 Int_t index2= _fcn->GetFloatParamList()->index(&var2);
730 coutE(Minimization) <<
"RooMinimizer::contour(" << GetName()
731 <<
") ERROR: " << var2.GetName()
732 <<
" is not a floating parameter of PDF "
733 << _func->GetName() << endl ;
738 RooPlot* frame =
new RooPlot(var1,var2) ;
741 TMarker *point=
new TMarker(var1.getVal(), var2.getVal(), 8);
742 frame->addObject(point) ;
746 if (_theFitter->GetMinimizer()==0) {
747 coutW(Minimization) <<
"RooMinimizer::contour: Error, run Migrad before contours!"
754 Double_t errdef= _theFitter->GetMinimizer()->ErrorDef();
757 n[0] = n1 ; n[1] = n2 ; n[2] = n3 ; n[3] = n4 ; n[4] = n5 ; n[5] = n6 ;
759 for (Int_t ic = 0 ; ic<6 ; ic++) {
763 _theFitter->GetMinimizer()->SetErrorDef(n[ic]*n[ic]*errdef);
766 Double_t *xcoor =
new Double_t[npoints+1];
767 Double_t *ycoor =
new Double_t[npoints+1];
768 bool ret = _theFitter->GetMinimizer()->Contour(index1,index2,npoints,xcoor,ycoor);
771 coutE(Minimization) <<
"RooMinimizer::contour("
773 <<
") ERROR: MINUIT did not return a contour graph for n="
776 xcoor[npoints] = xcoor[0];
777 ycoor[npoints] = ycoor[0];
778 TGraph* graph =
new TGraph(npoints+1,xcoor,ycoor);
780 graph->SetName(Form(
"contour_%s_n%f",_func->GetName(),n[ic])) ;
781 graph->SetLineStyle(ic+1) ;
782 graph->SetLineWidth(2) ;
783 graph->SetLineColor(kBlue) ;
784 frame->addObject(graph,
"L") ;
794 _theFitter->Config().MinimizerOptions().SetErrorDef(errdef);
797 *params = *paramSave ;
808 void RooMinimizer::profileStart()
812 _cumulTimer.Start(_profileStart?kFALSE:kTRUE) ;
813 _profileStart = kTRUE ;
821 void RooMinimizer::profileStop()
826 coutI(Minimization) <<
"Command timer: " ; _timer.Print() ;
827 coutI(Minimization) <<
"Session timer: " ; _cumulTimer.Print() ;
840 void RooMinimizer::applyCovarianceMatrix(TMatrixDSym& V)
842 _extV = (TMatrixDSym*) V.Clone() ;
843 _fcn->ApplyCovarianceMatrix(*_extV);
849 RooFitResult* RooMinimizer::lastMinuitFit(
const RooArgList& varList)
854 if (_theFitter==0 || _theFitter->GetMinimizer()==0) {
855 oocoutE((TObject*)0,InputArguments) <<
"RooMinimizer::save: Error, run minimization before!"
861 if (varList.getSize()>0 && varList.getSize()!=Int_t(_theFitter->Result().NTotalParameters())) {
862 oocoutE((TObject*)0,InputArguments)
863 <<
"RooMinimizer::lastMinuitFit: ERROR: supplied variable list must be either empty " << endl
864 <<
" or match the number of variables of the last fit ("
865 << _theFitter->Result().NTotalParameters() <<
")" << endl ;
871 TIterator* iter = varList.createIterator() ;
873 while((arg=(RooAbsArg*)iter->Next())) {
874 if (!dynamic_cast<RooRealVar*>(arg)) {
875 oocoutE((TObject*)0,InputArguments) <<
"RooMinimizer::lastMinuitFit: ERROR: variable '"
876 << arg->GetName() <<
"' is not of type RooRealVar" << endl ;
882 RooFitResult* res =
new RooFitResult(
"lastMinuitFit",
"Last MINUIT fit") ;
886 RooArgList constPars(
"constPars") ;
887 RooArgList floatPars(
"floatPars") ;
890 for (i = 0; i < _theFitter->Result().NTotalParameters(); ++i) {
892 TString varName(_theFitter->Result().GetParameterName(i));
893 Bool_t isConst(_theFitter->Result().IsParameterFixed(i)) ;
895 Double_t xlo = _theFitter->Config().ParSettings(i).LowerLimit();
896 Double_t xhi = _theFitter->Config().ParSettings(i).UpperLimit();
897 Double_t xerr = _theFitter->Result().Error(i);
898 Double_t xval = _theFitter->Result().Value(i);
901 if (varList.getSize()==0) {
903 if ((xlo<xhi) && !isConst) {
904 var =
new RooRealVar(varName,varName,xval,xlo,xhi) ;
906 var =
new RooRealVar(varName,varName,xval) ;
908 var->setConstant(isConst) ;
911 var = (RooRealVar*) varList.at(i)->Clone() ;
912 var->setConstant(isConst) ;
915 var->setRange(xlo,xhi) ;
918 if (varName.CompareTo(var->GetName())) {
919 oocoutI((TObject*)0,Eval) <<
"RooMinimizer::lastMinuitFit: fit parameter '" << varName
920 <<
"' stored in variable '" << var->GetName() <<
"'" << endl ;
926 constPars.addOwned(*var) ;
928 var->setError(xerr) ;
929 floatPars.addOwned(*var) ;
933 res->setConstParList(constPars) ;
934 res->setInitParList(floatPars) ;
935 res->setFinalParList(floatPars) ;
936 res->setMinNLL(_theFitter->Result().MinFcnValue()) ;
937 res->setEDM(_theFitter->Result().Edm()) ;
938 res->setCovQual(_theFitter->GetMinimizer()->CovMatrixStatus()) ;
939 res->setStatus(_theFitter->Result().Status()) ;
940 std::vector<double> globalCC;
941 TMatrixDSym corrs(_theFitter->Result().Parameters().size()) ;
942 TMatrixDSym covs(_theFitter->Result().Parameters().size()) ;
943 for (UInt_t ic=0; ic<_theFitter->Result().Parameters().size(); ic++) {
944 globalCC.push_back(_theFitter->Result().GlobalCC(ic));
945 for (UInt_t ii=0; ii<_theFitter->Result().Parameters().size(); ii++) {
946 corrs(ic,ii) = _theFitter->Result().Correlation(ic,ii);
947 covs(ic,ii) = _theFitter->Result().CovMatrix(ic,ii);
950 res->fillCorrMatrix(globalCC,corrs,covs) ;