1 #include "../interface/ProfileLikelihood.h"
2 #include "RooRealVar.h"
5 #include "RooDataSet.h"
6 #include "RooFitResult.h"
7 #include "../interface/RooMinimizerOpt.h"
10 #include "TFitResult.h"
12 #include "RooStats/ProfileLikelihoodCalculator.h"
13 #include "RooStats/LikelihoodInterval.h"
14 #include "RooStats/LikelihoodIntervalPlot.h"
15 #include "RooStats/HypoTestResult.h"
16 #include "RooStats/RooStatsUtils.h"
17 #include "../interface/Combine.h"
18 #include "../interface/CloseCoutSentry.h"
19 #include "../interface/utils.h"
20 #include "../interface/ProfiledLikelihoodRatioTestStatExt.h"
21 #include "../interface/CascadeMinimizer.h"
24 #include <Math/MinimizerOptions.h>
26 using namespace RooStats;
48 LimitAlgo(
"Profile Likelihood specific options")
51 (
"minimizerAlgo", boost::program_options::value<std::string>(&
minimizerAlgo_)->default_value(
minimizerAlgo_),
"Choice of minimizer (Minuit vs Minuit2)")
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")
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")
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")
74 if (vm.count(
"usePLC"))
useMinos_ =
false;
75 else if (vm.count(
"useMinos"))
useMinos_ =
true;
79 mass_ = vm[
"mass"].as<
float>();
83 minimizerTypeBackup(
ROOT::Math::MinimizerOptions::DefaultMinimizerType()),
84 minimizerAlgoBackup(
ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo()),
85 minimizerTollBackup(
ROOT::Math::MinimizerOptions::DefaultTolerance())
87 ROOT::Math::MinimizerOptions::SetDefaultTolerance(tolerance);
88 if (minimizerAlgo.find(
",") != std::string::npos) {
89 size_t idx = minimizerAlgo.find(
",");
90 std::string
type = minimizerAlgo.substr(0,idx),
algo = minimizerAlgo.substr(idx+1);
91 if (
verbose > 1)
std::cout <<
"Set default minimizer to " << type <<
", algorithm " <<
algo << std::endl;
92 ROOT::Math::MinimizerOptions::SetDefaultMinimizer(type.c_str(),
algo.c_str());
94 if (
verbose > 1)
std::cout <<
"Set default minimizer to " << minimizerAlgo << std::endl;
95 ROOT::Math::MinimizerOptions::SetDefaultMinimizer(minimizerAlgo.c_str());
101 ROOT::Math::MinimizerOptions::SetDefaultTolerance(minimizerTollBackup);
102 ROOT::Math::MinimizerOptions::SetDefaultMinimizer(minimizerTypeBackup.c_str(),minimizerAlgoBackup.empty() ? 0 : minimizerAlgoBackup.c_str());
105 bool ProfileLikelihood::run(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
109 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(mc_s->GetParametersOfInterest()->first());
111 std::vector<double> limits;
double rMax = r->getMax();
112 std::auto_ptr<RooAbsPdf> nuisancePdf(0);
114 w->loadSnapshot(
"clean");
116 r->setMax(rMax*(0.5+RooRandom::uniform()));
117 r->setVal((0.1+0.5*RooRandom::uniform())*r->getMax());
120 RooArgSet
set(*mc_s->GetNuisanceParameters());
121 RooDataSet *randoms = nuisancePdf->generate(
set, 1);
122 set = *randoms->get(0);
124 std::cout <<
"Starting minimization from point " << std::endl;
133 RooFitResult *res = mc_s->GetPdf()->fitTo(data, RooFit::Save(1), RooFit::Minimizer(
"Minuit2"));
135 if (
verbose > 1)
std::cout <<
"Fit failed (covQual " << (res ? res->covQual() : -1) <<
", edm " << (res ? res->edm() : 0) <<
")" << std::endl;
140 std::cout <<
"Covariance quality: " << res->covQual() <<
", Edm = " << res->edm() << std::endl;
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;
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;
156 spreadIn =
max(spreadIn, diff);
159 spreadOut =
max(spreadOut, diff);
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;
181 RooArgSet poi(*mc_s->GetParametersOfInterest());
182 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(poi.first());
183 double rMax = r->getMax();
188 ProfileLikelihoodCalculator plcB(data, *mc_s, 1.0-
cl);
189 std::auto_ptr<LikelihoodInterval> plInterval;
199 limitErr = le.second;
202 plInterval.reset(plcB.GetInterval());
203 if (plInterval.get() == 0)
break;
204 limit =
lowerLimit_ ? plInterval->LowerLimit(*r) : plInterval->UpperLimit(*r);
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);
212 if (limit == r->getMin()) {
213 std::cerr <<
"ProfileLikelihoodCalculator failed (returned upper limit equal to the lower bound)" << std::endl;
217 if (!
plot_.empty()) {
218 TCanvas *
c1 =
new TCanvas(
"c1",
"c1");
219 LikelihoodIntervalPlot
plot(&*plInterval);
221 c1->Print(
plot_.c_str());
228 std::cout <<
"\n -- Profile Likelihood -- " <<
"\n";
230 std::cout <<
"Limit: " << r->GetName() << (
lowerLimit_ ?
" > " :
" < ") << limit <<
" +/- " << limitErr <<
" @ " <<
cl * 100 <<
"% CL" << std::endl;
232 std::cout <<
"Limit: " << r->GetName() << (
lowerLimit_ ?
" > " :
" < ") << limit <<
" @ " <<
cl * 100 <<
"% CL" << std::endl;
240 RooArgSet poi(*mc_s->GetParametersOfInterest());
241 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(poi.first());
243 ProfileLikelihoodCalculator plcS(data, *mc_s, 1.0-
cl);
244 RooArgSet nullParamValues;
246 nullParamValues.addClone(*r);
247 plcS.SetNullParameters(nullParamValues);
255 if (q0 == -1)
return false;
256 limit = (q0 > 0 ?
sqrt(2*q0) : 0);
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;
263 std::auto_ptr<HypoTestResult>
result(plcS.GetHypoTest());
264 if (
result.get() == 0)
return false;
265 limit =
result->Significance();
268 if (limit == 0 && std::signbit(limit)) {
273 if (
reportPVal_) limit = RooStats::SignificanceToPValue(limit);
275 std::cout <<
"\n -- Profile Likelihood -- " <<
"\n";
276 std::cout << (
reportPVal_ ?
"p-value of background: " :
"Significance: ") << limit << std::endl;
278 if (
reportPVal_)
std::cout <<
" (Significance = " << RooStats::PValueToSignificance(limit) <<
")" << std::endl;
279 else std::cout <<
" (p-value = " << RooStats::SignificanceToPValue(limit) <<
")" << std::endl;
286 std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
288 minim.setStrategy(0);
289 minim.setPrintLevel(
verbose-1);
290 minim.setErrorLevel(0.5*TMath::ChisquareQuantile(cl,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();
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;
302 return std::numeric_limits<double>::quiet_NaN();
304 if (
verbose > 1) res->Print(
"V");
309 poi.setConstant(
false);
310 std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
312 minim0.setStrategy(0);
313 minim0.setPrintLevel(-1);
315 poi.setConstant(
true);
317 minim.setPrintLevel(-1);
319 std::cerr <<
"Initial minimization failed. Aborting." << std::endl;
320 return std::pair<double,double>(0, -1);
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);
331 printf(
" %-6s delta(NLL)\n",poi.GetName());
332 printf(
"%8.5f %8.5f\n", rval, 0.);
337 minim.setStrategy(0);
339 if (success ==
false) {
340 std::cerr <<
"Minimization failed at " << poi.getVal() <<
". exiting the bisection loop" << std::endl;
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) {
350 rval = 0.5*(rval + rhigh);
353 rval = 0.5*(rval + rlow);
355 }
while (fabs(rhigh - rlow) > tolerance);
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;
364 minim.setStrategy(0);
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);
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) {
380 }
while (rval < poi.getMax() && rval > poi.getMin());
381 return std::pair<double,double>(poi.getMax(), 0);
383 return std::pair<double,double>(poi.getVal(), fabs(rhigh - rlow)*0.5);
389 poi.setConstant(
false);
391 poi.setVal(0.05*poi.getMax());
392 std::auto_ptr<RooAbsReal> nll(pdf.createNLL(data, RooFit::Constrain(*nuisances)));
394 minim0.setStrategy(0);
396 if (poi.getVal() < 0) {
397 printf(
"Minimum found at %s = %8.5f < 0: significance will be zero\n", poi.GetName(), poi.getVal());
400 poi.setConstant(
true);
403 std::cerr <<
"Initial minimization failed. Aborting." << std::endl;
406 printf(
"Minimum found at %s = %8.5f\n", poi.GetName(), poi.getVal());
409 std::auto_ptr<RooFitResult>
start(minim.
save());
410 double minnll = nll->getVal(), thisnll = minnll, lastnll = thisnll;
411 double rbest = poi.getVal(),
rval = rbest;
414 printf(
" %-6s delta(NLL)\n", poi.GetName());
415 printf(
"%8.5f %8.5f\n",
rval, 0.);
417 points =
new TGraph(1);
418 points->SetName(Form(
"nll_scan_%g",
mass_));
419 points->SetPoint(0,
rval, 0);
421 while (
rval >= tolerance * poi.getMax()) {
427 thisnll = nll->getVal();
428 if (success ==
false) {
429 std::cerr <<
"Minimization failed at " << poi.getVal() <<
". exiting the loop" << std::endl;
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);
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;
446 if (lastnll > thisnll) {
447 std::cout <<
"Secondary minimum found, back to original minimization" << std::endl;
450 poi.setConstant(
false);
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();
457 poi.setConstant(
true);
463 return (thisnll - minnll);
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());
477 if (poi.getVal() < 0) {
478 printf(
"Minimum found at %s = %8.5f < 0: significance will be zero\n", poi.GetName(), poi.getVal());
481 maxScan = poi.getVal()*1.4;
485 poi.setConstant(
true);
488 std::cerr <<
"Initial minimization failed. Aborting." << std::endl;
491 printf(
"Minimum found at %s = %8.5f\n", poi.GetName(), poi.getVal());
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);
502 printf(
" %-6s delta(NLL)\n", poi.GetName());
503 printf(
"%8.5f %8.5f\n",
rval, 0.);
507 rval = (maxScan * (stepDown ?
i : steps-
i-1))/
steps;
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; }
519 printf(
"\nPlay it again, sam.\n");
520 printf(
" %-6s delta(NLL)\n", poi.GetName());
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;
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; }
536 fit->SetParameters(1, rbest, 2.0, minnll-refnll);
538 double ret = (maxnll - minnll);
542 res = points->Fit(fit,
"S0");
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;
551 std::cout <<
"Using first and last value" << std::endl;
bool runLimit(RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
static int maxOutliers_
Stop trying after finding N outliers.
static float signalForSignificance_
bool runSignificance(RooWorkspace *w, RooStats::ModelConfig *mc, RooAbsData &data, double &limit, double &limitErr)
virtual Double_t Evaluate(RooAbsData &data, RooArgSet &nullPOI)
static std::string static float minimizerToleranceForBF_
bool robustMinimize(RooAbsReal &nll, RooMinimizerOpt &minimizer, int verbosity=0)
static float maxOutlierFraction_
Ignore up to this fraction of results if they're too far from the median.
bool minimize(int verbose=0, bool cascade=true)
double significanceFromScan(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, int npoints) const
virtual void applyOptions(const boost::program_options::variables_map &vm)
RooAbsPdf * makeNuisancePdf(RooStats::ModelConfig &model, const char *name="nuisancePdf")
bool improve(int verbose=0, bool cascade=true)
static int maxTries_
trying up to M times from different points
static bool preFit_
Try first a plain fit.
virtual bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
double upperLimitWithMinos(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
const T & max(const T &a, const T &b)
std::pair< double, double > upperLimitBruteForce(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance, double cl) const
unsigned long long int rval
static float maxRelDeviation_
maximum relative deviation of the different points from the median to accept
MinimizerSentry(const std::string &algo, double tolerance)
double significanceBruteForce(RooAbsPdf &pdf, RooAbsData &data, RooRealVar &poi, const RooArgSet *nuisances, double tolerance) const
Setup Minimizer configuration on creation, reset the previous one on destruction. ...
static std::string bfAlgo_
static std::string minimizerAlgoForBF_
static std::string static float minimizerTolerance_
static std::string minimizerAlgo_
char data[epos_bytes_allocation]
static int tries_
compute the limit N times
boost::program_options::options_description options_
void setStrategy(int strategy)
void set(const std::string &name, int value)
set the flag, with a run-time name
static bool reportPVal_
Report p-value instead of significance.