CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Protected Member Functions | Static Protected Attributes
ProfileLikelihood Class Reference

#include <ProfileLikelihood.h>

Inheritance diagram for ProfileLikelihood:
LimitAlgo

Classes

class  MinimizerSentry
 Setup Minimizer configuration on creation, reset the previous one on destruction. More...
 

Public Member Functions

virtual void applyOptions (const boost::program_options::variables_map &vm)
 
virtual const std::string & name () const
 
 ProfileLikelihood ()
 
virtual bool run (RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
 
- Public Member Functions inherited from LimitAlgo
virtual void applyDefaultOptions ()
 
 LimitAlgo ()
 
 LimitAlgo (const char *desc)
 
const
boost::program_options::options_description & 
options () const
 
virtual void setNToys (const int)
 
virtual void setToyNumber (const int)
 
virtual ~LimitAlgo ()
 

Protected Member Functions

bool runLimit (RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
 
bool runSignificance (RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
 
double significanceBruteForce (RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance) const
 
double significanceFromScan (RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, int npoints) const
 
std::pair< double, double > upperLimitBruteForce (RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
 
double upperLimitWithMinos (RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
 

Static Protected Attributes

static std::string bfAlgo_ = "scale"
 
static bool bruteForce_ = false
 
static float mass_
 
static float maxOutlierFraction_ = 0.25
 Ignore up to this fraction of results if they're too far from the median. More...
 
static int maxOutliers_ = 3
 Stop trying after finding N outliers. More...
 
static float maxRelDeviation_ = 0.05
 maximum relative deviation of the different points from the median to accept More...
 
static int maxTries_ = 1
 trying up to M times from different points More...
 
static std::string minimizerAlgo_ = "Minuit2"
 
static std::string minimizerAlgoForBF_ = "Minuit2,simplex"
 
static std::string static float minimizerTolerance_ = 1e-2
 
static std::string static float minimizerToleranceForBF_ = 1e-4
 
static std::string plot_ = ""
 
static int points_ = 20
 
static bool preFit_ = false
 Try first a plain fit. More...
 
static bool reportPVal_ = false
 Report p-value instead of significance. More...
 
static float signalForSignificance_ = 0
 
static int tries_ = 1
 compute the limit N times More...
 
static bool useMinos_ = true
 

Additional Inherited Members

- Protected Attributes inherited from LimitAlgo
boost::program_options::options_description options_
 

Detailed Description

Class for computing limits and significances from naive Profile Likelihood asymptotics (i.e. no Asimov dataset)

Author
Giovanni Petrucciani (UCSD)

Definition at line 15 of file ProfileLikelihood.h.

Constructor & Destructor Documentation

ProfileLikelihood::ProfileLikelihood ( )

Definition at line 47 of file ProfileLikelihood.cc.

References bfAlgo_, maxOutlierFraction_, maxOutliers_, maxRelDeviation_, maxTries_, minimizerAlgo_, minimizerAlgoForBF_, minimizerTolerance_, minimizerToleranceForBF_, LimitAlgo::options_, plot_, points_, signalForSignificance_, and tries_.

47  :
48  LimitAlgo("Profile Likelihood specific options")
49 {
50  options_.add_options()
51  ("minimizerAlgo", boost::program_options::value<std::string>(&minimizerAlgo_)->default_value(minimizerAlgo_), "Choice of minimizer (Minuit vs Minuit2)")
52  ("minimizerTolerance", boost::program_options::value<float>(&minimizerTolerance_)->default_value(minimizerTolerance_), "Tolerance for minimizer")
53  ("tries", boost::program_options::value<int>(&tries_)->default_value(tries_), "Compute PL limit N times, to check for numerical instabilities")
54  ("maxTries", boost::program_options::value<int>(&maxTries_)->default_value(maxTries_), "Stop trying after N attempts per point")
55  ("maxRelDeviation", boost::program_options::value<float>(&maxRelDeviation_)->default_value(maxOutlierFraction_), "Max absolute deviation of the results from the median")
56  ("maxOutlierFraction", boost::program_options::value<float>(&maxOutlierFraction_)->default_value(maxOutlierFraction_), "Ignore up to this fraction of results if they're too far from the median")
57  ("signalForSignificance", boost::program_options::value<float>(&signalForSignificance_)->default_value(signalForSignificance_), "Signal strength used when computing significances (default is zero, just background)")
58  ("maxOutliers", boost::program_options::value<int>(&maxOutliers_)->default_value(maxOutliers_), "Stop trying after finding N outliers")
59  ("plot", boost::program_options::value<std::string>(&plot_)->default_value(plot_), "Save a plot of the negative log of the profiled likelihood into the specified file")
60  ("pvalue", "Report p-value instead of significance (when running with --significance)")
61  ("preFit", "Attept a fit before running the ProfileLikelihood calculator")
62  ("usePLC", "Compute PL limit using the ProfileLikelihoodCalculator (not default)")
63  ("useMinos", "Compute PL limit using Minos directly, bypassing the ProfileLikelihoodCalculator (default)")
64  ("bruteForce", "Compute PL limit by brute force, bypassing the ProfileLikelihoodCalculator and Minos")
65  ("bfAlgo", boost::program_options::value<std::string>(&bfAlgo_)->default_value(bfAlgo_), "NLL scan algorithm used for --bruteForce. Supported values are 'scale' (default), 'stepUp[Twice]', 'stepDown[Twice]'")
66  ("scanPoints", boost::program_options::value<int>(&points_)->default_value(points_), "Points for the scan")
67  ("minimizerAlgoForBF", boost::program_options::value<std::string>(&minimizerAlgoForBF_)->default_value(minimizerAlgoForBF_), "Choice of minimizer for brute-force search")
68  ("minimizerToleranceForBF", boost::program_options::value<float>(&minimizerToleranceForBF_)->default_value(minimizerToleranceForBF_), "Tolerance for minimizer when doing brute-force search")
69  ;
70 }
static int maxOutliers_
Stop trying after finding N outliers.
static float signalForSignificance_
static std::string static float minimizerToleranceForBF_
static float maxOutlierFraction_
Ignore up to this fraction of results if they&#39;re too far from the median.
static int maxTries_
trying up to M times from different points
static std::string plot_
LimitAlgo()
Definition: LimitAlgo.h:19
static float maxRelDeviation_
maximum relative deviation of the different points from the median to accept
static std::string bfAlgo_
static std::string minimizerAlgoForBF_
static std::string static float minimizerTolerance_
static std::string minimizerAlgo_
static int tries_
compute the limit N times
boost::program_options::options_description options_
Definition: LimitAlgo.h:32

Member Function Documentation

void ProfileLikelihood::applyOptions ( const boost::program_options::variables_map &  vm)
virtual

Reimplemented from LimitAlgo.

Definition at line 72 of file ProfileLikelihood.cc.

References bruteForce_, mass_, reportPVal_, and useMinos_.

73 {
74  if (vm.count("usePLC")) useMinos_ = false;
75  else if (vm.count("useMinos")) useMinos_ = true;
76  else useMinos_ = true;
77  bruteForce_ = vm.count("bruteForce");
78  reportPVal_ = vm.count("pvalue");
79  mass_ = vm["mass"].as<float>();
80 }
static bool bruteForce_
static bool reportPVal_
Report p-value instead of significance.
virtual const std::string& ProfileLikelihood::name ( void  ) const
inlinevirtual

Implements LimitAlgo.

Definition at line 19 of file ProfileLikelihood.h.

19  {
20  static const std::string name("ProfileLikelihood");
21  return name;
22  }
virtual const std::string & name() const
bool ProfileLikelihood::run ( RooWorkspace *  w,
RooStats::ModelConfig *  mc_s,
RooStats::ModelConfig *  mc_b,
RooAbsData &  data,
double &  limit,
double &  limitErr,
const double *  hint 
)
virtual

Implements LimitAlgo.

Definition at line 105 of file ProfileLikelihood.cc.

References gather_cfg::cout, diffTreeTool::diff, doSignificance_, i, j, utils::makeNuisancePdf(), max(), maxOutlierFraction_, maxOutliers_, maxRelDeviation_, maxTries_, minimizerAlgo_, minimizerTolerance_, preFit_, alignCSCRings::r, runLimit(), runSignificance(), runtimedef::set(), python.multivaluedict::sort(), summarizeEdmComparisonLogfiles::success, tries_, and withSystematics.

105  {
106  MinimizerSentry minimizerConfig(minimizerAlgo_, minimizerTolerance_);
107  CloseCoutSentry sentry(verbose < 0);
108 
109  RooRealVar *r = dynamic_cast<RooRealVar *>(mc_s->GetParametersOfInterest()->first());
110  bool success = false;
111  std::vector<double> limits; double rMax = r->getMax();
112  std::auto_ptr<RooAbsPdf> nuisancePdf(0);
113  for (int i = 0; i < maxTries_; ++i) {
114  w->loadSnapshot("clean");
115  if (i > 0) { // randomize starting point
116  r->setMax(rMax*(0.5+RooRandom::uniform()));
117  r->setVal((0.1+0.5*RooRandom::uniform())*r->getMax());
118  if (withSystematics) {
119  if (nuisancePdf.get() == 0) nuisancePdf.reset(utils::makeNuisancePdf(*mc_s));
120  RooArgSet set(*mc_s->GetNuisanceParameters());
121  RooDataSet *randoms = nuisancePdf->generate(set, 1);
122  set = *randoms->get(0);
123  if (verbose > 2) {
124  std::cout << "Starting minimization from point " << std::endl;
125  r->Print("V");
126  set.Print("V");
127  }
128  delete randoms;
129  }
130  }
131  if (preFit_) {
132  CloseCoutSentry sentry(verbose < 2);
133  RooFitResult *res = mc_s->GetPdf()->fitTo(data, RooFit::Save(1), RooFit::Minimizer("Minuit2"));
134  if (res == 0 || res->covQual() != 3 || res->edm() > minimizerTolerance_) {
135  if (verbose > 1) std::cout << "Fit failed (covQual " << (res ? res->covQual() : -1) << ", edm " << (res ? res->edm() : 0) << ")" << std::endl;
136  continue;
137  }
138  if (verbose > 1) {
139  res->Print("V");
140  std::cout << "Covariance quality: " << res->covQual() << ", Edm = " << res->edm() << std::endl;
141  }
142  delete res;
143  }
144  bool thisTry = (doSignificance_ ? runSignificance(w,mc_s,data,limit,limitErr) : runLimit(w,mc_s,data,limit,limitErr));
145  if (!thisTry) continue;
146  if (tries_ == 1) { success = true; break; }
147  limits.push_back(limit);
148  int nresults = limits.size();
149  if (nresults < tries_) continue;
150  std::sort(limits.begin(), limits.end());
151  double median = (nresults % 2 ? limits[nresults/2] : 0.5*(limits[nresults/2] + limits[nresults/2+1]));
152  int noutlier = 0; double spreadIn = 0, spreadOut = 0;
153  for (int j = 0; j < nresults; ++j) {
154  double diff = fabs(limits[j]-median)/median;
155  if (diff < maxRelDeviation_) {
156  spreadIn = max(spreadIn, diff);
157  } else {
158  noutlier++;
159  spreadOut = max(spreadOut, diff);
160  }
161  }
162  if (verbose > 0) {
163  std::cout << "Numer of tries: " << i << " Number of successes: " << nresults
164  << ", Outliers: " << noutlier << " (frac = " << noutlier/double(nresults) << ")"
165  << ", Spread of non-outliers: " << spreadIn <<" / of outliers: " << spreadOut << std::endl;
166  }
167  if (noutlier <= maxOutlierFraction_*nresults) {
168  if (verbose > 0) std::cout << " \\--> success! " << std::endl;
169  success = true;
170  limit = median;
171  break;
172  } else if (noutlier > maxOutliers_) {
173  if (verbose > 0) std::cout << " \\--> failure! " << std::endl;
174  break;
175  }
176  }
177  return success;
178 }
bool runLimit(RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
int i
Definition: DBlmapReader.cc:9
static int maxOutliers_
Stop trying after finding N outliers.
bool runSignificance(RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
bool doSignificance_
Definition: Combine.cc:69
static float maxOutlierFraction_
Ignore up to this fraction of results if they&#39;re too far from the median.
bool withSystematics
Definition: Combine.cc:68
RooAbsPdf * makeNuisancePdf(RooStats::ModelConfig &model, const char *name="nuisancePdf")
Definition: utils.cc:200
static int maxTries_
trying up to M times from different points
static bool preFit_
Try first a plain fit.
const T & max(const T &a, const T &b)
int j
Definition: DBlmapReader.cc:9
static float maxRelDeviation_
maximum relative deviation of the different points from the median to accept
static std::string static float minimizerTolerance_
static std::string minimizerAlgo_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
static int tries_
compute the limit N times
T w() const
void set(const std::string &name, int value)
set the flag, with a run-time name
bool ProfileLikelihood::runLimit ( RooWorkspace *  w,
RooStats::ModelConfig *  mc,
RooAbsData &  data,
double &  limit,
double &  limitErr 
)
protected

Definition at line 180 of file ProfileLikelihood.cc.

References bruteForce_, alignmentValidation::c1, dtNoiseDBValidation_cfg::cerr, cl, CloseCoutSentry::clear(), gather_cfg::cout, data, alignCSCRings::e, edm::detail::isnan(), vdt::le, lowerLimit_, minimizerTolerance_, plotResiduals::plot(), plot_, alignCSCRings::r, summarizeEdmComparisonLogfiles::success, upperLimitBruteForce(), upperLimitWithMinos(), and useMinos_.

Referenced by run().

180  {
181  RooArgSet poi(*mc_s->GetParametersOfInterest());
182  RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
183  double rMax = r->getMax();
184  bool success = false;
185  CloseCoutSentry coutSentry(verbose <= 1); // close standard output and error, so that we don't flood them with minuit messages
186 
187  while (!success) {
188  ProfileLikelihoodCalculator plcB(data, *mc_s, 1.0-cl);
189  std::auto_ptr<LikelihoodInterval> plInterval;
190  if (useMinos_ || bruteForce_) {
191  // try first with Minos, unless brute force requested
192  if (!bruteForce_) {
193  limit = upperLimitWithMinos(*mc_s->GetPdf(), data, *r, mc_s->GetNuisanceParameters(), minimizerTolerance_, cl);
194  }
195  // if brute force forced, or minos failed, go with the next one
196  if (std::isnan(limit) || bruteForce_) {
197  std::pair<double,double> le = upperLimitBruteForce(*mc_s->GetPdf(), data, *r, mc_s->GetNuisanceParameters(), 1e-3*minimizerTolerance_, cl);
198  limit = le.first;
199  limitErr = le.second;
200  }
201  } else {
202  plInterval.reset(plcB.GetInterval());
203  if (plInterval.get() == 0) break;
204  limit = lowerLimit_ ? plInterval->LowerLimit(*r) : plInterval->UpperLimit(*r);
205  }
206  if (limit >= 0.75*r->getMax()) {
207  std::cout << "Limit " << r->GetName() << " < " << limit << "; " << r->GetName() << " max < " << r->getMax() << std::endl;
208  if (r->getMax()/rMax > 20) break;
209  r->setMax(r->getMax()*2);
210  continue;
211  }
212  if (limit == r->getMin()) {
213  std::cerr << "ProfileLikelihoodCalculator failed (returned upper limit equal to the lower bound)" << std::endl;
214  break;
215  }
216  success = true;
217  if (!plot_.empty()) {
218  TCanvas *c1 = new TCanvas("c1","c1");
219  LikelihoodIntervalPlot plot(&*plInterval);
220  plot.Draw();
221  c1->Print(plot_.c_str());
222  delete c1;
223  }
224  }
225  coutSentry.clear();
226  if (verbose >= 0) {
227  if (success) {
228  std::cout << "\n -- Profile Likelihood -- " << "\n";
229  if (limitErr) {
230  std::cout << "Limit: " << r->GetName() << (lowerLimit_ ? " > " : " < ") << limit << " +/- " << limitErr << " @ " << cl * 100 << "% CL" << std::endl;
231  } else {
232  std::cout << "Limit: " << r->GetName() << (lowerLimit_ ? " > " : " < ") << limit << " @ " << cl * 100 << "% CL" << std::endl;
233  }
234  }
235  }
236  return success;
237 }
bool lowerLimit_
Definition: Combine.cc:70
static bool bruteForce_
static std::string plot_
double upperLimitWithMinos(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
bool isnan(float x)
Definition: math.h:13
float cl
Definition: Combine.cc:71
std::pair< double, double > upperLimitBruteForce(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
static std::string static float minimizerTolerance_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
unsigned long long le
Definition: VDTMath.h:202
bool ProfileLikelihood::runSignificance ( RooWorkspace *  w,
RooStats::ModelConfig *  mc,
RooAbsData &  data,
double &  limit,
double &  limitErr 
)
protected

Definition at line 239 of file ProfileLikelihood.cc.

References bfAlgo_, bruteForce_, dtNoiseDBValidation_cfg::cerr, cl, CloseCoutSentry::clear(), gather_cfg::cout, data, ProfiledLikelihoodTestStatOpt::Evaluate(), minimizerTolerance_, points_, alignCSCRings::r, reportPVal_, query::result, signalForSignificance_, significanceBruteForce(), significanceFromScan(), mathSSE::sqrt(), and useMinos_.

Referenced by run().

239  {
240  RooArgSet poi(*mc_s->GetParametersOfInterest());
241  RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
242 
243  ProfileLikelihoodCalculator plcS(data, *mc_s, 1.0-cl);
244  RooArgSet nullParamValues;
245  r->setVal(signalForSignificance_);
246  nullParamValues.addClone(*r);
247  plcS.SetNullParameters(nullParamValues);
248 
249  CloseCoutSentry coutSentry(verbose <= 1); // close standard output and error, so that we don't flood them with minuit messages
250 
251  if (bruteForce_) {
252  double q0 = -1;
253  if (bfAlgo_ == "scale") q0 = significanceBruteForce(*mc_s->GetPdf(), data, *r, mc_s->GetNuisanceParameters(), 0.1*minimizerTolerance_);
254  else q0 = significanceFromScan(*mc_s->GetPdf(), data, *r, mc_s->GetNuisanceParameters(), 0.1*minimizerTolerance_, points_);
255  if (q0 == -1) return false;
256  limit = (q0 > 0 ? sqrt(2*q0) : 0);
257  } else if (useMinos_) {
258  ProfiledLikelihoodTestStatOpt testStat(*mc_s->GetObservables(), *mc_s->GetPdf(), mc_s->GetNuisanceParameters(),
259  nullParamValues, *mc_s->GetParametersOfInterest(), RooArgList(), RooArgList(), verbose-1);
260  Double_t q0 = testStat.Evaluate(data, nullParamValues);
261  limit = q0 > 0 ? sqrt(2*q0) : 0;
262  } else {
263  std::auto_ptr<HypoTestResult> result(plcS.GetHypoTest());
264  if (result.get() == 0) return false;
265  limit = result->Significance();
266  }
267  coutSentry.clear();
268  if (limit == 0 && std::signbit(limit)) {
269  //..... This is not an error, it just means we have a deficit of events.....
270  std::cerr << "The minimum of the likelihood is for r <= " << signalForSignificance_ << ", so the significance is zero" << std::endl;
271  limit = 0;
272  }
273  if (reportPVal_) limit = RooStats::SignificanceToPValue(limit);
274 
275  std::cout << "\n -- Profile Likelihood -- " << "\n";
276  std::cout << (reportPVal_ ? "p-value of background: " : "Significance: ") << limit << std::endl;
277  if (verbose > 0) {
278  if (reportPVal_) std::cout << " (Significance = " << RooStats::PValueToSignificance(limit) << ")" << std::endl;
279  else std::cout << " (p-value = " << RooStats::SignificanceToPValue(limit) << ")" << std::endl;
280  }
281  return true;
282 }
static float signalForSignificance_
virtual Double_t Evaluate(RooAbsData &data, RooArgSet &nullPOI)
double significanceFromScan(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, int npoints) const
static bool bruteForce_
T sqrt(T t)
Definition: SSEVec.h:46
tuple result
Definition: query.py:137
float cl
Definition: Combine.cc:71
double significanceBruteForce(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance) const
static std::string bfAlgo_
static std::string static float minimizerTolerance_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
static bool reportPVal_
Report p-value instead of significance.
double ProfileLikelihood::significanceBruteForce ( RooAbsPdf &  pdf,
RooAbsData &  data,
RooRealVar &  poi,
const RooArgSet *  nuisances,
double  tolerance 
) const
protected

Definition at line 388 of file ProfileLikelihood.cc.

References dtNoiseDBValidation_cfg::cerr, CascadeMinimizer::Constrained, gather_cfg::cout, CascadeMinimizer::improve(), mass_, CascadeMinimizer::minimize(), minimizerAlgo_, minimizerAlgoForBF_, minimizerTolerance_, minimizerToleranceForBF_, download_sqlite_cfg::outputFile, CascadeMinimizer::save(), CascadeMinimizer::setStrategy(), errorMatrix2Lands_multiChannel::start, summarizeEdmComparisonLogfiles::success, and CascadeMinimizer::Unconstrained.

Referenced by runSignificance().

388  {
389  poi.setConstant(false);
390  //poi.setMin(0);
391  poi.setVal(0.05*poi.getMax());
392  std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
394  minim0.setStrategy(0);
395  minim0.minimize(verbose-2);
396  if (poi.getVal() < 0) {
397  printf("Minimum found at %s = %8.5f < 0: significance will be zero\n", poi.GetName(), poi.getVal());
398  return 0;
399  }
400  poi.setConstant(true);
402  if (!minim.minimize(verbose-2)) {
403  std::cerr << "Initial minimization failed. Aborting." << std::endl;
404  return -1;
405  } else if (verbose > 0) {
406  printf("Minimum found at %s = %8.5f\n", poi.GetName(), poi.getVal());
407  }
408  MinimizerSentry minimizerConfig(minimizerAlgoForBF_, minimizerToleranceForBF_);
409  std::auto_ptr<RooFitResult> start(minim.save());
410  double minnll = nll->getVal(), thisnll = minnll, lastnll = thisnll;
411  double rbest = poi.getVal(), rval = rbest;
412  TGraph *points = 0;
413  if (verbose) {
414  printf(" %-6s delta(NLL)\n", poi.GetName());
415  printf("%8.5f %8.5f\n", rval, 0.);
416  fflush(stdout);
417  points = new TGraph(1);
418  points->SetName(Form("nll_scan_%g", mass_));
419  points->SetPoint(0, rval, 0);
420  }
421  while (rval >= tolerance * poi.getMax()) {
422  rval *= 0.8;
423  poi.setVal(rval);
424  minim.setStrategy(0);
425  bool success = minim.improve(verbose-2, /*cascade=*/false);
426  lastnll = thisnll;
427  thisnll = nll->getVal();
428  if (success == false) {
429  std::cerr << "Minimization failed at " << poi.getVal() <<". exiting the loop" << std::endl;
430  return -1;
431  }
432  if (verbose) {
433  printf("%8.5f %8.5f\n", rval, thisnll-minnll); fflush(stdout);
434  points->Set(points->GetN()+1);
435  points->SetPoint(points->GetN()-1, rval, thisnll - minnll);
436  }
437  if (fabs(lastnll - thisnll) < 7*minimizerToleranceForBF_) {
438  std::cout << "This is enough." << std::endl;
439  if (thisnll < lastnll) {
440  std::cout << "Linear extrapolation from " << thisnll << " to " << ( thisnll - (lastnll-thisnll)*rval )<< std::endl;
441  thisnll -= (lastnll-thisnll)*rval;
442  }
443  break;
444  }
445 #if 0
446  if (lastnll > thisnll) {
447  std::cout << "Secondary minimum found, back to original minimization" << std::endl;
448  {
449  MinimizerSentry minimizerConfig(minimizerAlgo_, minimizerTolerance_);
450  poi.setConstant(false);
451  poi.setVal(rbest);
452  minim.improve(verbose - 1);
453  printf("New global minimum at %8.5f (was %8.5f), the shift in nll is %8.5f\n", poi.getVal(), rbest, nll->getVal()-minnll);
454  minnll = nll->getVal();
455  rbest = poi.getVal();
456  rval = rbest;
457  poi.setConstant(true);
458  }
459  }
460 #endif
461  }
462  if (points) outputFile->WriteTObject(points);
463  return (thisnll - minnll);
464 }
static std::string static float minimizerToleranceForBF_
unsigned long long int rval
Definition: vlib.h:23
static std::string minimizerAlgoForBF_
static std::string static float minimizerTolerance_
static std::string minimizerAlgo_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
double ProfileLikelihood::significanceFromScan ( RooAbsPdf &  pdf,
RooAbsData &  data,
RooRealVar &  poi,
const RooArgSet *  nuisances,
double  tolerance,
int  npoints 
) const
protected

Definition at line 466 of file ProfileLikelihood.cc.

References bfAlgo_, dtNoiseDBValidation_cfg::cerr, CascadeMinimizer::Constrained, gather_cfg::cout, i, CascadeMinimizer::improve(), mass_, CascadeMinimizer::minimize(), minimizerAlgo_, minimizerAlgoForBF_, minimizerTolerance_, minimizerToleranceForBF_, download_sqlite_cfg::outputFile, run_regression::ret, CascadeMinimizer::save(), CascadeMinimizer::setStrategy(), errorMatrix2Lands_multiChannel::start, relval_steps::steps, summarizeEdmComparisonLogfiles::success, and CascadeMinimizer::Unconstrained.

Referenced by runSignificance().

466  {
467  std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
468  double maxScan = poi.getMax()*0.7;
469  bool stepDown = (bfAlgo_.find("stepDown") != std::string::npos);
470  bool twice = (bfAlgo_.find("Twice") != std::string::npos);
471  poi.setConstant(false);
472  poi.setVal(0.05*poi.getMax());
474  minim0.setStrategy(0);
475  minim0.minimize(verbose-2);
476  if (!stepDown) {
477  if (poi.getVal() < 0) {
478  printf("Minimum found at %s = %8.5f < 0: significance will be zero\n", poi.GetName(), poi.getVal());
479  return 0;
480  }
481  maxScan = poi.getVal()*1.4;
482  } else {
483  poi.setVal(0);
484  }
485  poi.setConstant(true);
487  if (!minim.minimize(verbose-2)) {
488  std::cerr << "Initial minimization failed. Aborting." << std::endl;
489  return -1;
490  } else if (verbose > 0) {
491  printf("Minimum found at %s = %8.5f\n", poi.GetName(), poi.getVal());
492  }
493  MinimizerSentry minimizerConfig(minimizerAlgoForBF_, minimizerToleranceForBF_);
494  std::auto_ptr<RooFitResult> start(minim.save());
495  double minnll = nll->getVal(), thisnll = minnll, refnll = thisnll, maxnll = thisnll;
496  double rbest = poi.getVal(), rval = rbest;
497  TGraph *points = new TGraph(steps+1); points->SetName(Form("nll_scan_%g", mass_));
498  TF1 *fit = new TF1("fit","[0]*pow(abs(x-[1]), [2])+[3]",0,poi.getMax());
499  fit->SetParNames("norm","bestfit","power","offset");
500  points->SetPoint(0, rval, 0);
501  if (verbose) {
502  printf(" %-6s delta(NLL)\n", poi.GetName());
503  printf("%8.5f %8.5f\n", rval, 0.);
504  fflush(stdout);
505  }
506  for (int i = 1; i < steps; ++i) {
507  rval = (maxScan * (stepDown ? i : steps-i-1))/steps;
508  poi.setVal(rval);
509  bool success = minim.improve(verbose-2, /*cascade=*/false);
510  thisnll = nll->getVal();
511  if (success == false) std::cerr << "Minimization failed at " << poi.getVal() <<"." << std::endl;
512  if (verbose) { printf("%8.5f %8.5f\n", rval, thisnll-refnll); fflush(stdout); }
513  points->SetPoint(i, rval, thisnll-refnll);
514  if (thisnll < minnll) { minnll = thisnll; rbest = rval; }
515  if (rval <= rbest && thisnll > maxnll) { maxnll = thisnll; }
516  }
517  if (twice) {
518  if (verbose) {
519  printf("\nPlay it again, sam.\n");
520  printf(" %-6s delta(NLL)\n", poi.GetName());
521  fflush(stdout);
522  }
523  for (int i = steps-1; i >= 0; --i) {
524  rval = (maxScan * (stepDown ? i : steps-i-1))/steps;
525  if (i == 0 && !stepDown) rval = rbest;
526  poi.setVal(rval);
527  bool success = minim.improve(verbose-2, /*cascade=*/false);
528  thisnll = nll->getVal();
529  if (success == false) std::cerr << "Minimization failed at " << poi.getVal() <<"." << std::endl;
530  if (verbose) { printf("%8.5f %8.5f\n", rval, thisnll-refnll); fflush(stdout); }
531  points->SetPoint(i, rval, thisnll-refnll);
532  if (thisnll < minnll) { minnll = thisnll; rbest = rval; }
533  if (rval <= rbest && thisnll > maxnll) { maxnll = thisnll; }
534  }
535  }
536  fit->SetParameters(1, rbest, 2.0, minnll-refnll);
537  points->Sort();
538  double ret = (maxnll - minnll);
539  TFitResultPtr res;
540  {
541  MinimizerSentry minimizerConfig(minimizerAlgo_, minimizerTolerance_);
542  res = points->Fit(fit,"S0");
543  }
544  outputFile->WriteTObject(points);
545  outputFile->WriteTObject(fit);
546  if (res.Get()->Status() == 0) {
547  std::cout << "Using first and last value, the result would be " << ret << std::endl;
548  ret = fit->Eval(0) - fit->Eval(fit->GetParameter(1)) ;
549  std::cout << "Using output of fit, the result is " << ret << std::endl;
550  } else {
551  std::cout << "Using first and last value" << std::endl;
552  }
553  return ret;
554 }
int i
Definition: DBlmapReader.cc:9
static std::string static float minimizerToleranceForBF_
unsigned long long int rval
Definition: vlib.h:23
static std::string bfAlgo_
static std::string minimizerAlgoForBF_
static std::string static float minimizerTolerance_
static std::string minimizerAlgo_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
std::pair< double, double > ProfileLikelihood::upperLimitBruteForce ( RooAbsPdf &  pdf,
RooAbsData &  data,
RooRealVar &  poi,
const RooArgSet *  nuisances,
double  tolerance,
double  cl 
) const
protected

Definition at line 308 of file ProfileLikelihood.cc.

References dtNoiseDBValidation_cfg::cerr, cmsPerfPublish::fail(), lowerLimit_, minimizerAlgoForBF_, minimizerToleranceForBF_, nllutils::robustMinimize(), errorMatrix2Lands_multiChannel::start, summarizeEdmComparisonLogfiles::success, and filterCSVwithJSON::target.

Referenced by runLimit().

308  {
309  poi.setConstant(false);
310  std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
311  RooMinimizerOpt minim0(*nll);
312  minim0.setStrategy(0);
313  minim0.setPrintLevel(-1);
314  nllutils::robustMinimize(*nll, minim0, verbose-2);
315  poi.setConstant(true);
316  RooMinimizerOpt minim(*nll);
317  minim.setPrintLevel(-1);
318  if (!nllutils::robustMinimize(*nll, minim, verbose-2)) {
319  std::cerr << "Initial minimization failed. Aborting." << std::endl;
320  return std::pair<double,double>(0, -1);
321  }
322  std::auto_ptr<RooFitResult> start(minim.save());
323  double minnll = nll->getVal();
324  double rval = poi.getVal() + (lowerLimit_ ? -3 : +3)*poi.getError(), rlow = poi.getVal(), rhigh = lowerLimit_ ? poi.getMin() : poi.getMax();
325  if (rval >= rhigh || rval <= rlow) rval = 0.5*(rlow + rhigh);
326  double target = minnll + 0.5*TMath::ChisquareQuantile(cl,1);
327  //minim.setPrintLevel(verbose-2);
328  MinimizerSentry minimizerConfig(minimizerAlgoForBF_, minimizerToleranceForBF_);
329  bool fail = false;
330  if (verbose) {
331  printf(" %-6s delta(NLL)\n",poi.GetName());
332  printf("%8.5f %8.5f\n", rval, 0.);
333  fflush(stdout);
334  }
335  do {
336  poi.setVal(rval);
337  minim.setStrategy(0);
338  bool success = nllutils::robustMinimize(*nll, minim, verbose-2);
339  if (success == false) {
340  std::cerr << "Minimization failed at " << poi.getVal() <<". exiting the bisection loop" << std::endl;
341  fail = true;
342  break;
343  }
344  double nllthis = nll->getVal();
345  if (verbose) { printf("%8.5f %8.5f\n", rval, nllthis-minnll); fflush(stdout); }
346  if (fabs(nllthis - target) < tolerance) {
347  return std::pair<double,double>(rval, (rhigh - rlow)*0.5);
348  } else if (nllthis < target) {
349  (lowerLimit_ ? rhigh : rlow) = rval;
350  rval = 0.5*(rval + rhigh);
351  } else {
352  (lowerLimit_ ? rlow : rhigh) = rval;
353  rval = 0.5*(rval + rlow);
354  }
355  } while (fabs(rhigh - rlow) > tolerance);
356  if (fail) {
357  // try do do it in small steps instead
358  std::auto_ptr<RooArgSet> pars(nll->getParameters((const RooArgSet *)0));
359  double dx = (lowerLimit_ ? -0.05 : +0.05)*poi.getError();
360  *pars = start->floatParsFinal();
361  rval = poi.getVal() + dx;
362  do {
363  poi.setVal(rval);
364  minim.setStrategy(0);
365  bool success = nllutils::robustMinimize(*nll, minim, verbose-2);
366  if (success == false) {
367  std::cerr << "Minimization failed at " << poi.getVal() <<". exiting the stepping loop" << std::endl;
368  return std::pair<double,double>(poi.getVal(), fabs(rhigh - rlow)*0.5);
369  }
370  double nllthis = nll->getVal();
371  if (verbose) { printf("%8.5f %8.5f\n", rval, nllthis-minnll); fflush(stdout); }
372  if (fabs(nllthis - target) < tolerance) {
373  return std::pair<double,double>(rval, fabs(dx));
374  } else if (nllthis < target) {
375  rval += dx;
376  } else {
377  dx *= 0.5;
378  rval -= dx;
379  }
380  } while (rval < poi.getMax() && rval > poi.getMin());
381  return std::pair<double,double>(poi.getMax(), 0);
382  } else {
383  return std::pair<double,double>(poi.getVal(), fabs(rhigh - rlow)*0.5);
384  }
385 }
bool lowerLimit_
Definition: Combine.cc:70
static std::string static float minimizerToleranceForBF_
bool robustMinimize(RooAbsReal &nll, RooMinimizerOpt &minimizer, int verbosity=0)
float cl
Definition: Combine.cc:71
unsigned long long int rval
Definition: vlib.h:23
static std::string minimizerAlgoForBF_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double ProfileLikelihood::upperLimitWithMinos ( RooAbsPdf &  pdf,
RooAbsData &  data,
RooRealVar &  poi,
const RooArgSet *  nuisances,
double  tolerance,
double  cl 
) const
protected

Definition at line 285 of file ProfileLikelihood.cc.

References dtNoiseDBValidation_cfg::cerr, MessageLogger_cff::limit, lowerLimit_, and nllutils::robustMinimize().

Referenced by runLimit().

285  {
286  std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
287  RooMinimizerOpt minim(*nll);
288  minim.setStrategy(0);
289  minim.setPrintLevel(verbose-1);
290  minim.setErrorLevel(0.5*TMath::ChisquareQuantile(cl,1));
291  nllutils::robustMinimize(*nll, minim, verbose-1);
292  int minosStat = minim.minos(RooArgSet(poi));
293  if (minosStat == -1) return std::numeric_limits<double>::quiet_NaN();
294  std::auto_ptr<RooFitResult> res(minim.save());
295  double muhat = poi.getVal(), limit = poi.getVal() + (lowerLimit_ ? poi.getAsymErrorLo() : poi.getAsymErrorHi());
296  double nll0 = nll->getVal();
297  poi.setVal(limit);
298  double nll2 = nll->getVal();
299  if (nll2 < nll0 + 0.75 * 0.5*TMath::ChisquareQuantile(cl,1)) {
300  std::cerr << "ERROR: unprofiled likelihood gives better result than profiled one. deltaNLL = " << (nll2-nll0) << ". will try brute force." << std::endl;
301  poi.setVal(muhat);
302  return std::numeric_limits<double>::quiet_NaN();
303  }
304  if (verbose > 1) res->Print("V");
305  return limit;
306 }
bool lowerLimit_
Definition: Combine.cc:70
bool robustMinimize(RooAbsReal &nll, RooMinimizerOpt &minimizer, int verbosity=0)
float cl
Definition: Combine.cc:71
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82

Member Data Documentation

std::string ProfileLikelihood::bfAlgo_ = "scale"
staticprotected

Definition at line 40 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), runSignificance(), and significanceFromScan().

bool ProfileLikelihood::bruteForce_ = false
staticprotected

Definition at line 39 of file ProfileLikelihood.h.

Referenced by applyOptions(), runLimit(), and runSignificance().

float ProfileLikelihood::mass_
staticprotected

Definition at line 61 of file ProfileLikelihood.h.

Referenced by applyOptions(), significanceBruteForce(), and significanceFromScan().

float ProfileLikelihood::maxOutlierFraction_ = 0.25
staticprotected

Ignore up to this fraction of results if they're too far from the median.

Definition at line 51 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and run().

int ProfileLikelihood::maxOutliers_ = 3
staticprotected

Stop trying after finding N outliers.

Definition at line 53 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and run().

float ProfileLikelihood::maxRelDeviation_ = 0.05
staticprotected

maximum relative deviation of the different points from the median to accept

Definition at line 49 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and run().

int ProfileLikelihood::maxTries_ = 1
staticprotected

trying up to M times from different points

Definition at line 47 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and run().

std::string ProfileLikelihood::minimizerAlgo_ = "Minuit2"
staticprotected
std::string ProfileLikelihood::minimizerAlgoForBF_ = "Minuit2,simplex"
staticprotected
float ProfileLikelihood::minimizerTolerance_ = 1e-2
staticprotected
float ProfileLikelihood::minimizerToleranceForBF_ = 1e-4
staticprotected
std::string ProfileLikelihood::plot_ = ""
staticprotected

Definition at line 63 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and runLimit().

int ProfileLikelihood::points_ = 20
staticprotected

Definition at line 41 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and runSignificance().

bool ProfileLikelihood::preFit_ = false
staticprotected

Try first a plain fit.

Definition at line 55 of file ProfileLikelihood.h.

Referenced by run().

bool ProfileLikelihood::reportPVal_ = false
staticprotected

Report p-value instead of significance.

Definition at line 58 of file ProfileLikelihood.h.

Referenced by applyOptions(), and runSignificance().

float ProfileLikelihood::signalForSignificance_ = 0
staticprotected

Definition at line 60 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and runSignificance().

int ProfileLikelihood::tries_ = 1
staticprotected

compute the limit N times

Definition at line 45 of file ProfileLikelihood.h.

Referenced by ProfileLikelihood(), and run().

bool ProfileLikelihood::useMinos_ = true
staticprotected

Definition at line 39 of file ProfileLikelihood.h.

Referenced by applyOptions(), runLimit(), and runSignificance().