CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
nllutils Namespace Reference

Functions

bool robustMinimize (RooAbsReal &nll, RooMinimizerOpt &minimizer, int verbosity=0)
 

Function Documentation

bool nllutils::robustMinimize ( RooAbsReal &  nll,
RooMinimizerOpt minimizer,
int  verbosity = 0 
)

Definition at line 383 of file ProfiledLikelihoodRatioTestStatExt.cc.

References COUNT_ONE, DBG, RooMinimizerOpt::edm(), runtimedef::get(), run_regression::ret, and ntuplemaker::status.

Referenced by CascadeMinimizer::improveOnce(), ProfileLikelihood::upperLimitBruteForce(), and ProfileLikelihood::upperLimitWithMinos().

384 {
385  static bool do_debug = (runtimedef::get("DEBUG_MINIM") || runtimedef::get("DEBUG_PLTSO") > 1);
386  double initialNll = nll.getVal();
387  std::auto_ptr<RooArgSet> pars;
388  bool ret = false;
389  for (int tries = 0, maxtries = 4; tries <= maxtries; ++tries) {
390  int status = minim.minimize(ROOT::Math::MinimizerOptions::DefaultMinimizerType().c_str(), ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str());
391  //if (verbosity > 1) res->Print("V");
392  if (status == 0 && nll.getVal() > initialNll + 0.02) {
393  std::auto_ptr<RooFitResult> res(minim.save());
394  //PerfCounter::add("Minimizer.save() called for false minimum");
395  DBG(DBG_PLTestStat_main, (printf("\n --> false minimum, status %d, cov. quality %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->covQual(), res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal())))
396  if (pars.get() == 0) pars.reset(nll.getParameters((const RooArgSet*)0));
397  *pars = res->floatParsInit();
398  if (tries == 0) {
399  COUNT_ONE("nllutils::robustMinimize: false minimum (first try)")
400  DBG(DBG_PLTestStat_main, (printf(" ----> Doing a re-scan and re-trying\n")))
401  minim.minimize("Minuit2","Scan");
402  } else if (tries == 1) {
403  COUNT_ONE("nllutils::robustMinimize: false minimum (second try)")
404  DBG(DBG_PLTestStat_main, (printf(" ----> Re-trying with strategy = 1\n")))
405  minim.setStrategy(1);
406  } else if (tries == 2) {
407  COUNT_ONE("nllutils::robustMinimize: false minimum (third try)")
408  DBG(DBG_PLTestStat_main, (printf(" ----> Re-trying with strategy = 2\n")))
409  minim.setStrategy(2);
410  } else {
411  COUNT_ONE("nllutils::robustMinimize: false minimum (third try)")
412  DBG(DBG_PLTestStat_main, (printf(" ----> Last attempt: simplex method \n")))
413  status = minim.minimize("Minuit2","Simplex");
414  if (nll.getVal() < initialNll + 0.02) {
415  DBG(DBG_PLTestStat_main, (printf("\n --> success: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal())))
416  if (do_debug) printf("\n --> success: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal());
417  COUNT_ONE("nllutils::robustMinimize: final success")
418  ret = true;
419  break;
420  } else {
421  COUNT_ONE("nllutils::robustMinimize: final fail")
422  DBG(DBG_PLTestStat_main, (printf("\n --> final fail: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal())))
423  if (do_debug) printf("\n --> final fail: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal());
424  return false;
425  }
426  }
427  } else if (status == 0) {
428  DBG(DBG_PLTestStat_main, (printf("\n --> success: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal())))
429  if (do_debug) printf("\n --> success: status %d, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, initialNll, nll.getVal(), initialNll - nll.getVal());
430  COUNT_ONE("nllutils::robustMinimize: final success")
431  ret = true;
432  break;
433  } else if (tries != maxtries) {
434  std::auto_ptr<RooFitResult> res(do_debug ? minim.save() : 0);
435  //PerfCounter::add("Minimizer.save() called for failed minimization");
436  if (tries > 0 && minim.edm() < 0.05*ROOT::Math::MinimizerOptions::DefaultTolerance()) {
437  DBG(DBG_PLTestStat_main, (printf("\n --> acceptable: status %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal())))
438  if (do_debug) printf("\n --> acceptable: status %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal());
439  COUNT_ONE("nllutils::robustMinimize: accepting fit with bad status but good EDM")
440  COUNT_ONE("nllutils::robustMinimize: final success")
441  ret = true;
442  break;
443  }
444  DBG(DBG_PLTestStat_main, (printf("\n --> partial fail: status %d, cov. quality %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->covQual(), res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal())))
445  if (tries == 1) {
446  COUNT_ONE("nllutils::robustMinimize: failed first attempt")
447  DBG(DBG_PLTestStat_main, (printf(" ----> Doing a re-scan first, and switching to strategy 1\n")))
448  minim.minimize("Minuit2","Scan");
449  minim.setStrategy(1);
450  }
451  if (tries == 2) {
452  COUNT_ONE("nllutils::robustMinimize: failed second attempt")
453  DBG(DBG_PLTestStat_main, (printf(" ----> trying with strategy = 2\n")))
454  minim.minimize("Minuit2","Scan");
455  minim.setStrategy(2);
456  }
457  } else {
458  std::auto_ptr<RooFitResult> res(do_debug ? minim.save() : 0);
459  DBG(DBG_PLTestStat_main, (printf("\n --> final fail: status %d, cov. quality %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->covQual(), res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal())))
460  if (do_debug) printf("\n --> final fail: status %d, cov. quality %d, edm %10.7f, nll initial % 10.4f, nll final % 10.4f, change %10.5f\n", status, res->covQual(), res->edm(), initialNll, nll.getVal(), initialNll - nll.getVal());
461  COUNT_ONE("nllutils::robustMinimize: final fail")
462  }
463  }
464  return ret;
465 }
bool robustMinimize(RooAbsReal &nll, RooMinimizerOpt &minimizer, int verbosity=0)
int get(const char *name)
int bad(Items const &cont)
double f[11][100]
bool first
Definition: L1TdeRCT.cc:94
double a
Definition: hdecay.h:121
perl if(1 lt scalar(@::datatypes))
Definition: edlooper.cc:31
#define DBG(X, Z)
tuple status
Definition: ntuplemaker.py:245