#include <CascadeMinimizer.h>
Classes | |
struct | Algo |
compact information about an algorithm More... | |
Public Types | |
enum | Mode { Constrained, Unconstrained } |
Public Member Functions | |
CascadeMinimizer (RooAbsReal &nll, Mode mode, RooRealVar *poi=0, int initialStrategy=0) | |
bool | improve (int verbose=0, bool cascade=true) |
bool | minimize (int verbose=0, bool cascade=true) |
RooMinimizerOpt & | minimizer () |
RooFitResult * | save () |
void | setErrorLevel (float errorLevel) |
void | setStrategy (int strategy) |
void | trivialMinimize (const RooAbsReal &nll, RooRealVar &r, int points=100) const |
Static Public Member Functions | |
static void | applyOptions (const boost::program_options::variables_map &vm) |
static void | initOptions () |
static const boost::program_options::options_description & | options () |
Private Member Functions | |
bool | improveOnce (int verbose) |
Private Attributes | |
RooMinimizerOpt | minimizer_ |
Mode | mode_ |
RooAbsReal & | nll_ |
RooRealVar * | poi_ |
int | strategy_ |
Static Private Attributes | |
static std::vector< Algo > | fallbacks_ |
list of algorithms to run if the default one fails | |
static bool | oldFallback_ = true |
don't do old fallback using robustMinimize | |
static boost::program_options::options_description | options_ |
options configured from command line | |
static bool | poiOnlyFit_ |
do first a fit of only the POI | |
static bool | preScan_ |
do a pre-scan | |
static bool | setZeroPoint_ = true |
do first a fit of only the POI | |
static bool | singleNuisFit_ |
do first a minimization of each nuisance individually |
Definition at line 11 of file CascadeMinimizer.h.
Definition at line 14 of file CascadeMinimizer.h.
{ return minimizer_; }
CascadeMinimizer::CascadeMinimizer | ( | RooAbsReal & | nll, |
Mode | mode, | ||
RooRealVar * | poi = 0 , |
||
int | initialStrategy = 0 |
||
) |
Definition at line 20 of file CascadeMinimizer.cc.
void CascadeMinimizer::applyOptions | ( | const boost::program_options::variables_map & | vm | ) | [static] |
Definition at line 118 of file CascadeMinimizer.cc.
References algo, gather_cfg::cout, CascadeMinimizer::Algo::default_strategy(), CascadeMinimizer::Algo::default_tolerance(), fallbacks_, min, poiOnlyFit_, preScan_, setZeroPoint_, and singleNuisFit_.
{ using namespace std; preScan_ = vm.count("cminPreScan"); poiOnlyFit_ = vm.count("cminPoiOnlyFit"); singleNuisFit_ = vm.count("cminSingleNuisFit"); setZeroPoint_ = vm.count("cminSetZeroPoint"); if (vm.count("cminFallbackAlgo")) { vector<string> falls(vm["cminFallbackAlgo"].as<vector<string> >()); for (vector<string>::const_iterator it = falls.begin(), ed = falls.end(); it != ed; ++it) { std::string algo = *it; float tolerance = Algo::default_tolerance(); int strategy = Algo::default_strategy(); string::size_type idx = std::min(algo.find(";"), algo.find(":")); if (idx != string::npos && idx < algo.length()) { tolerance = atof(algo.substr(idx+1).c_str()); algo = algo.substr(0,idx); // DON'T SWAP THESE TWO LINES } idx = algo.find(","); if (idx != string::npos && idx < algo.length()) { // if after the comma there's a number, then it's a strategy if ( '0' <= algo[idx+1] && algo[idx+1] <= '9' ) { strategy = atoi(algo.substr(idx+1).c_str()); algo = algo.substr(0,idx); // DON'T SWAP THESE TWO LINES } else { // otherwise, it could be Name,subname,strategy idx = algo.find(",",idx+1); if (idx != string::npos && idx < algo.length()) { strategy = atoi(algo.substr(idx+1).c_str()); algo = algo.substr(0,idx); // DON'T SWAP THESE TWO LINES } } } fallbacks_.push_back(Algo(algo, tolerance, strategy)); std::cout << "Configured fallback algorithm " << fallbacks_.back().algo << ", strategy " << fallbacks_.back().strategy << ", tolerance " << fallbacks_.back().tolerance << std::endl; } } //if (vm.count("cminDefaultIntegratorEpsAbs")) RooAbsReal::defaultIntegratorConfig()->setEpsAbs(vm["cminDefaultIntegratorEpsAbs"].as<double>()); //if (vm.count("cminDefaultIntegratorEpsRel")) RooAbsReal::defaultIntegratorConfig()->setEpsRel(vm["cminDefaultIntegratorEpsRel"].as<double>()); //if (vm.count("cminDefaultIntegrator1D")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->method1D(), vm["cminDefaultIntegrator1D"].as<std::string>()); //if (vm.count("cminDefaultIntegrator1DOpen")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->method1DOpen(), vm["cminDefaultIntegrator1DOpen"].as<std::string>()); //if (vm.count("cminDefaultIntegrator2D")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->method2D(), vm["cminDefaultIntegrator2D"].as<std::string>()); //if (vm.count("cminDefaultIntegrator2DOpen")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->method2DOpen(), vm["cminDefaultIntegrator2DOpen"].as<std::string>()); //if (vm.count("cminDefaultIntegratorND")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->methodND(), vm["cminDefaultIntegratorND"].as<std::string>()); //if (vm.count("cminDefaultIntegratorNDOpen")) setDefaultIntegrator(RooAbsReal::defaultIntegratorConfig()->methodNDOpen(), vm["cminDefaultIntegratorNDOpen"].as<std::string>()); }
bool CascadeMinimizer::improve | ( | int | verbose = 0 , |
bool | cascade = true |
||
) |
Definition at line 29 of file CascadeMinimizer.cc.
References dtNoiseDBValidation_cfg::cerr, cacheutils::CachingSimNLL::clearZeroPoint(), CascadeMinimizer::Algo::default_strategy(), CascadeMinimizer::Algo::default_tolerance(), fallbacks_, improveOnce(), minimizer_, nll_, setZeroPoint_, and strategy_.
Referenced by FitterAlgoBase::doFit(), FitterAlgoBase::findCrossing(), Asymptotic::findExpectedLimitFromCrossing(), Asymptotic::getCLs(), minimize(), ProfileLikelihood::significanceBruteForce(), and ProfileLikelihood::significanceFromScan().
{ minimizer_.setPrintLevel(verbose-2); minimizer_.setStrategy(strategy_); bool outcome = improveOnce(verbose-2); if (cascade && !outcome && !fallbacks_.empty()) { std::string nominalType(ROOT::Math::MinimizerOptions::DefaultMinimizerType()); std::string nominalAlgo(ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo()); float nominalTol(ROOT::Math::MinimizerOptions::DefaultTolerance()); int nominalStrat(strategy_); if (verbose > 0) std::cerr << "Failed minimization with " << nominalType << "," << nominalAlgo << " and tolerance " << nominalTol << std::endl; for (std::vector<Algo>::const_iterator it = fallbacks_.begin(), ed = fallbacks_.end(); it != ed; ++it) { ProfileLikelihood::MinimizerSentry minimizerConfig(it->algo, it->tolerance != Algo::default_tolerance() ? it->tolerance : nominalTol); int myStrategy = it->strategy; if (myStrategy == Algo::default_strategy()) myStrategy = nominalStrat; if (nominalType != ROOT::Math::MinimizerOptions::DefaultMinimizerType() || nominalAlgo != ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo() || nominalTol != ROOT::Math::MinimizerOptions::DefaultTolerance() || myStrategy != nominalStrat) { if (verbose > 0) std::cerr << "Will fallback to minimization using " << it->algo << ", strategy " << myStrategy << " and tolerance " << it->tolerance << std::endl; minimizer_.setStrategy(myStrategy); outcome = improveOnce(verbose-2); if (outcome) break; } } } if (setZeroPoint_) { cacheutils::CachingSimNLL *simnll = dynamic_cast<cacheutils::CachingSimNLL *>(&nll_); if (simnll) simnll->clearZeroPoint(); } return outcome; }
bool CascadeMinimizer::improveOnce | ( | int | verbose | ) | [private] |
Definition at line 61 of file CascadeMinimizer.cc.
References cacheutils::CachingSimNLL::clearZeroPoint(), minimizer_, nll_, oldFallback_, nllutils::robustMinimize(), cacheutils::CachingSimNLL::setZeroPoint(), setZeroPoint_, and ntuplemaker::status.
Referenced by improve().
{ std::string myType(ROOT::Math::MinimizerOptions::DefaultMinimizerType()); std::string myAlgo(ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo()); if (setZeroPoint_) { cacheutils::CachingSimNLL *simnll = dynamic_cast<cacheutils::CachingSimNLL *>(&nll_); if (simnll) { simnll->setZeroPoint(); } } bool outcome = false; if (oldFallback_){ outcome = nllutils::robustMinimize(nll_, minimizer_, verbose); } else { int status = minimizer_.minimize(myType.c_str(), myAlgo.c_str()); outcome = (status == 0); } if (setZeroPoint_) { cacheutils::CachingSimNLL *simnll = dynamic_cast<cacheutils::CachingSimNLL *>(&nll_); if (simnll) simnll->clearZeroPoint(); } return outcome; }
void CascadeMinimizer::initOptions | ( | ) | [static] |
Definition at line 98 of file CascadeMinimizer.cc.
References oldFallback_, options_, setZeroPoint_, and relativeConstraints::value.
{ options_.add_options() ("cminPoiOnlyFit", "Do first a fit floating only the parameter of interest") ("cminPreScan", "Do a scan before first minimization") ("cminSingleNuisFit", "Do first a minimization of each nuisance parameter individually") ("cminFallbackAlgo", boost::program_options::value<std::vector<std::string> >(), "Fallback algorithms if the default minimizer fails (can use multiple ones). Syntax is algo[,subalgo][,strategy][:tolerance]") ("cminSetZeroPoint", boost::program_options::value<bool>(&setZeroPoint_)->default_value(setZeroPoint_), "Change the reference point of the NLL to be zero during minimization") ("cminOldRobustMinimize", boost::program_options::value<bool>(&oldFallback_)->default_value(oldFallback_), "Use the old 'robustMinimize' logic in addition to the cascade") //("cminDefaultIntegratorEpsAbs", boost::program_options::value<double>(), "RooAbsReal::defaultIntegratorConfig()->setEpsAbs(x)") //("cminDefaultIntegratorEpsRel", boost::program_options::value<double>(), "RooAbsReal::defaultIntegratorConfig()->setEpsRel(x)") //("cminDefaultIntegrator1D", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->method1D().setLabel(x)") //("cminDefaultIntegrator1DOpen", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->method1DOpen().setLabel(x)") //("cminDefaultIntegrator2D", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->method2D().setLabel(x)") //("cminDefaultIntegrator2DOpen", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->method2DOpen().setLabel(x)") //("cminDefaultIntegratorND", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->methodND().setLabel(x)") //("cminDefaultIntegratorNDOpen", boost::program_options::value<std::string>(), "RooAbsReal::defaultIntegratorConfig()->methodNDOpen().setLabel(x)") ; }
bool CascadeMinimizer::minimize | ( | int | verbose = 0 , |
bool | cascade = true |
||
) |
Definition at line 85 of file CascadeMinimizer.cc.
References improve(), minimizer_, mode_, nll_, poi_, poiOnlyFit_, preScan_, strategy_, trivialMinimize(), and Unconstrained.
Referenced by MultiDimFit::doContour2D(), FitterAlgoBase::doFit(), MultiDimFit::doGrid(), MultiDimFit::doRandomPoints(), Asymptotic::findExpectedLimitFromCrossing(), ProfiledLikelihoodTestStatOpt::minNLL(), BestFitSigmaTestStat::minNLL(), Asymptotic::runLimit(), ProfileLikelihood::significanceBruteForce(), and ProfileLikelihood::significanceFromScan().
{ minimizer_.setPrintLevel(verbose-2); minimizer_.setStrategy(strategy_); if (preScan_) minimizer_.minimize("Minuit2","Scan"); // FIXME can be made smarter than this if (mode_ == Unconstrained && poiOnlyFit_) { trivialMinimize(nll_, *poi_, 200); } return improve(verbose, cascade); }
RooMinimizerOpt& CascadeMinimizer::minimizer | ( | ) | [inline] |
Definition at line 20 of file CascadeMinimizer.h.
References minimizer().
Referenced by FitterAlgoBase::doFit(), Asymptotic::findExpectedLimitFromCrossing(), and minimizer().
{ return minimizer().save(); }
static const boost::program_options::options_description& CascadeMinimizer::options | ( | ) | [inline, static] |
Definition at line 26 of file CascadeMinimizer.h.
: RooAbsReal & nll_;
RooFitResult* CascadeMinimizer::save | ( | ) | [inline] |
Definition at line 21 of file CascadeMinimizer.h.
References strategy_.
Referenced by FitterAlgoBase::doFit(), FitterAlgoBase::findCrossing(), ProfileLikelihood::significanceBruteForce(), and ProfileLikelihood::significanceFromScan().
{ strategy_ = strategy; }
void CascadeMinimizer::setErrorLevel | ( | float | errorLevel | ) | [inline] |
Definition at line 23 of file CascadeMinimizer.h.
Referenced by FitterAlgoBase::doFit(), and Asymptotic::findExpectedLimitFromCrossing().
void CascadeMinimizer::setStrategy | ( | int | strategy | ) | [inline] |
Definition at line 22 of file CascadeMinimizer.h.
References minimizer_.
Referenced by MultiDimFit::doBox(), MultiDimFit::doContour2D(), FitterAlgoBase::doFit(), MultiDimFit::doGrid(), MultiDimFit::doRandomPoints(), Asymptotic::findExpectedLimitFromCrossing(), Asymptotic::getCLs(), Asymptotic::runLimit(), ProfileLikelihood::significanceBruteForce(), and ProfileLikelihood::significanceFromScan().
{ minimizer_.setErrorLevel(errorLevel); }
void CascadeMinimizer::trivialMinimize | ( | const RooAbsReal & | nll, |
RooRealVar & | r, | ||
int | points = 100 |
||
) | const |
Definition at line 183 of file CascadeMinimizer.cc.
References i, x, and detailsBasic3DVector::y.
Referenced by minimize().
{ double rMin = r.getMin(), rMax = r.getMax(), rStep = (rMax-rMin)/(points-1); int iMin = -1; double minnll = 0; for (int i = 0; i < points; ++i) { double x = rMin + (i+0.5)*rStep; r.setVal(x); double y = nll.getVal(); if (iMin == -1 || y < minnll) { minnll = y; iMin = i; } } r.setVal( rMin + (iMin+0.5)*rStep ); }
std::vector< CascadeMinimizer::Algo > CascadeMinimizer::fallbacks_ [static, private] |
list of algorithms to run if the default one fails
Definition at line 48 of file CascadeMinimizer.h.
Referenced by applyOptions(), and improve().
RooMinimizerOpt CascadeMinimizer::minimizer_ [private] |
Definition at line 30 of file CascadeMinimizer.h.
Referenced by improve(), improveOnce(), minimize(), and setStrategy().
Mode CascadeMinimizer::mode_ [private] |
Definition at line 31 of file CascadeMinimizer.h.
Referenced by minimize().
RooAbsReal& CascadeMinimizer::nll_ [private] |
Definition at line 29 of file CascadeMinimizer.h.
Referenced by improve(), improveOnce(), and minimize().
bool CascadeMinimizer::oldFallback_ = true [static, private] |
don't do old fallback using robustMinimize
Definition at line 58 of file CascadeMinimizer.h.
Referenced by improveOnce(), and initOptions().
boost::program_options::options_description CascadeMinimizer::options_ [static, private] |
options configured from command line
Definition at line 38 of file CascadeMinimizer.h.
Referenced by initOptions().
RooRealVar* CascadeMinimizer::poi_ [private] |
Definition at line 33 of file CascadeMinimizer.h.
Referenced by minimize().
bool CascadeMinimizer::poiOnlyFit_ [static, private] |
do first a fit of only the POI
Definition at line 52 of file CascadeMinimizer.h.
Referenced by applyOptions(), and minimize().
bool CascadeMinimizer::preScan_ [static, private] |
do a pre-scan
Definition at line 50 of file CascadeMinimizer.h.
Referenced by applyOptions(), and minimize().
bool CascadeMinimizer::setZeroPoint_ = true [static, private] |
do first a fit of only the POI
Definition at line 56 of file CascadeMinimizer.h.
Referenced by applyOptions(), improve(), improveOnce(), and initOptions().
bool CascadeMinimizer::singleNuisFit_ [static, private] |
do first a minimization of each nuisance individually
Definition at line 54 of file CascadeMinimizer.h.
Referenced by applyOptions().
int CascadeMinimizer::strategy_ [private] |
Definition at line 32 of file CascadeMinimizer.h.
Referenced by improve(), minimize(), and save().