3 #include "../interface/Asymptotic.h"
4 #include <RooRealVar.h>
7 #include <RooFitResult.h>
9 #include <RooStats/ModelConfig.h>
10 #include <Math/DistFuncMathCore.h>
11 #include "../interface/Combine.h"
12 #include "../interface/CloseCoutSentry.h"
13 #include "../interface/RooFitGlobalKillSentry.h"
14 #include "../interface/ProfiledLikelihoodRatioTestStatExt.h"
15 #include "../interface/ToyMCSamplerOpt.h"
16 #include "../interface/ProfileLikelihood.h"
17 #include "../interface/CascadeMinimizer.h"
18 #include "../interface/utils.h"
19 #include "../interface/AsimovUtils.h"
21 using namespace RooStats;
38 LimitAlgo(
"Asymptotic specific options") {
40 (
"rAbsAcc", boost::program_options::value<double>(&
rAbsAccuracy_)->default_value(
rAbsAccuracy_),
"Absolute accuracy on r to reach to terminate the scan")
41 (
"rRelAcc", boost::program_options::value<double>(&
rRelAccuracy_)->default_value(
rRelAccuracy_),
"Relative accuracy on r to reach to terminate the scan")
42 (
"run", boost::program_options::value<std::string>(&
what_)->default_value(
what_),
"What to run: both (default), observed, expected, blind.")
43 (
"singlePoint", boost::program_options::value<double>(&
rValue_),
"Just compute CLs for the given value of r")
44 (
"minimizerAlgo", boost::program_options::value<std::string>(&
minimizerAlgo_)->default_value(
minimizerAlgo_),
"Choice of minimizer used for profiling (Minuit vs Minuit2)")
47 (
"qtilde", boost::program_options::value<bool>(&
qtilde_)->default_value(
qtilde_),
"Allow only non-negative signal strengths (default is true).")
48 (
"picky",
"Abort on fit failures")
49 (
"noFitAsimov",
"Use the pre-fit asimov dataset")
50 (
"newExpected",
"Use the new formula for expected limits")
51 (
"minosAlgo", boost::program_options::value<std::string>(&
minosAlgo_)->default_value(
minosAlgo_),
"Algorithm to use to get the median expected limit: 'minos' (fastest), 'bisection', 'stepping' (default, most robust)")
56 if (vm.count(
"singlePoint") && !vm[
"singlePoint"].defaulted()) {
57 if (!vm[
"run"].defaulted())
throw std::invalid_argument(
"Asymptotic: when using --singlePoint you can't use --run (at least for now)");
58 what_ =
"singlePoint";
61 throw std::invalid_argument(
"Asymptotic: option 'run' can only be 'observed', 'expected' or 'both' (the default) or 'blind' (a-priori expected)");
63 picky_ = vm.count(
"picky");
67 if (
noFitAsimov_)
std::cout <<
"Will use a-priori expected background instead of a-posteriori one." << std::endl;
74 bool Asymptotic::run(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
81 std::vector<std::pair<float,float> > expected;
83 if (
what_ !=
"expected") ret =
runLimit(w, mc_s, mc_b, data, limit, limitErr, hint);
86 const char *rname = mc_s->GetParametersOfInterest()->first()->GetName();
87 std::cout <<
"\n -- Asymptotic -- " <<
"\n";
88 if (
what_ ==
"singlePoint") {
89 printf(
"Observed CLs for %s = %.1f: %6.4f \n", rname,
rValue_, limit);
90 }
else if (ret &&
what_ !=
"expected") {
91 printf(
"Observed Limit: %s < %6.4f\n", rname, limit);
93 for (std::vector<std::pair<float,float> >::const_iterator it = expected.begin(), ed = expected.end(); it != ed; ++it) {
94 printf(
"Expected %4.1f%%: %s < %6.4f\n", it->first*100, rname, it->second);
103 bool Asymptotic::runLimit(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
104 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(mc_s->GetParametersOfInterest()->first());
105 w->loadSnapshot(
"clean");
108 r->setConstant(
false);
109 r->setVal(0.1*r->getMax());
110 r->setMin(
qtilde_ ? 0 : -r->getMax());
112 if (
params_.get() == 0)
params_.reset(mc_s->GetPdf()->getParameters(data));
115 std::auto_ptr<TIterator> itparam(
params_->createIterator());
116 for (RooAbsArg *
a = (RooAbsArg *) itparam->Next();
a != 0;
a = (RooAbsArg *) itparam->Next()) {
117 RooRealVar *rrv =
dynamic_cast<RooRealVar *
>(
a);
123 nllD_.reset(mc_s->GetPdf()->createNLL(data, RooFit::Constrain(constraints)));
124 nllA_.reset(mc_s->GetPdf()->createNLL(asimov, RooFit::Constrain(constraints)));
126 if (
verbose > 0)
std::cout << (
qtilde_ ?
"Restricting" :
"Not restricting") <<
" " << r->GetName() <<
" to positive values." << std::endl;
129 if (
verbose > 0)
std::cout <<
"\nMake global fit of real data" << std::endl;
139 if (
verbose > 0)
std::cout <<
"NLL at global minimum of data: " <<
minNllD_ <<
" (" << r->GetName() <<
" = " << r->getVal() <<
")" << std::endl;
140 double rErr = std::max<double>(r->getError(), 0.02 * (r->getMax() - r->getMin()));
145 if (
verbose > 0)
std::cout <<
"\nMake global fit of asimov data" << std::endl;
155 if (
verbose > 0)
std::cout <<
"NLL at global minimum of asimov: " <<
minNllA_ <<
" (" << r->GetName() <<
" = " << r->getVal() <<
")" << std::endl;
159 r->setConstant(
true);
161 if (
what_ ==
"singlePoint") {
166 double clsTarget = 1-
cl;
167 double rMin = std::max<double>(0, r->getVal()), rMax = rMin + 3 * rErr;
168 double clsMax = 1, clsMin = 0;
169 for (
int tries = 0; tries < 5; ++tries) {
170 double cls =
getCLs(*r, rMax);
171 if (cls == -999) {
std::cerr <<
"Minimization failed in an unrecoverable way" << std::endl;
break; }
172 if (cls < clsTarget) { clsMin = cls;
break; }
177 if (clsMax < 3*clsTarget && clsMin > 0.3*clsTarget) {
178 double rCross = rMin + (rMax-rMin)*
log(clsMax/clsTarget)/
log(clsMax/clsMin);
179 if ((rCross-rMin) < (rMax - rCross)) {
180 limit = 0.8*rCross + 0.2*rMax;
182 limit = 0.8*rCross + 0.2*rMin;
184 limitErr = 0.5*(rMax - rMin);
186 limit = 0.5*(rMin + rMax);
187 limitErr = 0.5*(rMax - rMin);
189 double cls =
getCLs(*r, limit);
190 if (cls == -999) {
std::cerr <<
"Minimization failed in an unrecoverable way" << std::endl;
break; }
191 if (cls > clsTarget) {
204 r.setMax(1.1 * rVal);
221 double qmu = 2*(
nllD_->getVal() -
minNllD_);
if (qmu < 0) qmu = 0;
235 double qA = 2*(
nllA_->getVal() -
minNllA_);
if (qA < 0) qA = 0;
237 double CLsb = ROOT::Math::normal_cdf_c(
sqrt(qmu));
238 double CLb = ROOT::Math::normal_cdf(
sqrt(qA)-
sqrt(qmu));
241 double mos =
sqrt(qA);
242 CLsb = ROOT::Math::normal_cdf_c( (qmu + qA)/(2*mos) );
243 CLb = ROOT::Math::normal_cdf_c( (qmu - qA)/(2*mos) );
245 double CLs = (CLb == 0 ? 0 : CLsb/CLb);
247 if (
verbose > 0) printf(
"At %s = %f:\tq_mu = %.5f\tq_A = %.5f\tCLsb = %7.5f\tCLb = %7.5f\tCLs = %7.5f\n", r.GetName(), rVal, qmu, qA, CLsb, CLb, CLs);
249 if (getAlsoExpected) {
250 const double quantiles[5] = { 0.025, 0.16, 0.50, 0.84, 0.975 };
251 for (
int iq = 0; iq < 5; ++iq) {
252 double N = ROOT::Math::normal_quantile(quantiles[iq], 1.0);
253 double clb = quantiles[iq];
254 double clsplusb = ROOT::Math::normal_cdf_c(
sqrt(qA) - N, 1.);
255 *limit = (clb != 0 ? clsplusb/clb : 0); *limitErr = 0;
257 if (
verbose > 0) printf(
"Expected %4.1f%%: CLsb = %.5f CLb = %.5f CLs = %.5f\n", quantiles[iq]*100, clsplusb, clb, clsplusb/clb);
263 std::vector<std::pair<float,float> >
Asymptotic::runLimitExpected(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
276 RooArgSet poi(*mc_s->GetParametersOfInterest());
277 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(poi.first());
283 if (
params_.get() == 0)
params_.reset(mc_s->GetPdf()->getParameters(data));
287 r->setConstant(
false);
289 r->setMin(
qtilde_ ? 0 : -r->getMax());
290 r->setVal(0.01*r->getMax());
291 r->setError(0.1*r->getMax());
294 std::auto_ptr<RooAbsReal> nll(mc_s->GetPdf()->createNLL(*asimov, RooFit::Constrain(*mc_s->GetNuisanceParameters())));
297 minim.setErrorLevel(0.5*
pow(ROOT::Math::normal_quantile(1-0.5*(1-
cl),1.0), 2));
303 std::cout <<
"Fit to asimov dataset:" << std::endl;
304 std::auto_ptr<RooFitResult> res(minim.save());
307 if (r->getVal()/r->getMax() > 1
e-3) {
308 if (
verbose) printf(
"WARNING: Best fit of asimov dataset is at %s = %f (%f times %sMax), while it should be at zero\n",
309 r->GetName(), r->getVal(), r->getVal()/r->getMax(), r->GetName());
314 double nll0 = nll->getVal();
316 double sigma = median / ROOT::Math::normal_quantile(1-0.5*(1-
cl),1.0);
319 std::cout <<
"Median for expected limits: " << median << std::endl;
320 std::cout <<
"Sigma for expected limits: " << sigma << std::endl;
323 std::vector<std::pair<float,float> > expected;
324 const double quantiles[5] = { 0.025, 0.16, 0.50, 0.84, 0.975 };
325 for (
int iq = 0; iq < 5; ++iq) {
326 double N = ROOT::Math::normal_quantile(quantiles[iq], 1.0);
337 if (
std::isnan(limit)) { expected.clear();
break; }
339 limit = sigma*(ROOT::Math::normal_quantile(1 - alpha * quantiles[iq], 1.0) +
N);
343 expected.push_back(std::pair<float,float>(quantiles[iq], limit));
356 double N = ROOT::Math::normal_quantile(clb, 1.0);
357 double errorlevel = 0.5 *
pow(N+ROOT::Math::normal_quantile_c(clb*(1-
cl),1.0), 2);
360 double rMax0 = r->getMax();
368 for (
int tries = 0; tries < 3; ++tries) {
369 minosStat = minim.
minimizer().minos(RooArgSet(*r));
370 if (minosStat != -1) {
371 while ((minosStat != -1) && (r->getVal()+r->getAsymErrorHi())/r->getMax() > 0.9) {
372 if (r->getMax() >= 100*rMax0) { minosStat = -1;
break; }
373 r->setMax(2*r->getMax());
378 minosStat = minim2.
minimizer().minos(RooArgSet(*r));
385 minim.
minimizer().minimize(
"Minuit",
"minimize");
387 minim.
minimizer().minimize(
"Minuit2",
"minmize");
391 if (minosStat != -1)
return r->getVal()+r->getAsymErrorHi();
394 double rCross = 0.5*(rMin+rMax), rErr = 0.5*(rMax-rMin);
395 r->setVal(rCross); r->setConstant(
true);
399 if (
verbose > 1) printf(
"Will search for NLL crossing by bisection\n");
401 if (rCross >= r->getMax()) r->setMax(rCross*1.1);
408 if (!ok &&
picky_)
break;
else minosStat = 0;
409 double here = nll.getVal();
410 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll) = %.5f\n", r->GetName(), rCross, here-nll0);
412 if (here < threshold) rMin = rCross;
else rMax = rCross;
413 rCross = 0.5*(rMin+rMax); rErr = 0.5*(rMax-rMin);
416 if (
verbose > 1) printf(
"Will search for NLL crossing by stepping.\n");
417 rCross = 0.05 * rMax; rErr = rMax;
418 double stride = rCross;
bool overstepped =
false;
420 if (rCross >= r->getMax()) r->setMax(rCross*1.1);
421 double there = nll.getVal();
428 if (!ok &&
picky_)
break;
else minosStat = 0;
429 double here = nll.getVal();
430 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll) = %.5f\n", r->GetName(), rCross, here-nll0);
432 if (here < threshold) {
433 if ((threshold-here) < 0.5*fabs(threshold-there)) stride *= 0.5;
436 stride *= 0.5; overstepped =
true;
439 if (overstepped) rErr = stride;
442 if (
verbose > 1) printf(
"Will search for NLL crossing with new algorithm.\n");
458 double nll_0 = nll.getVal();
459 double rMax0 = rMax*100;
462 double rStep = 0.05 * (rMax - r_0);
463 double r_1 = r_0, nll_1 = nll_0;
466 if (r_1 >= r->getMax()) r->setMax(r_1*1.1);
468 nll_1 = nll.getVal();
470 bool binNLLchange = (nll_1 < threshold && nll_1 - nll_0 > 0.5);
471 bool aboveThresh = (nll_1 > threshold + kappa*
std::pow(r_1-r_0,2));
472 if (binNLLchange || aboveThresh) {
473 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll unprof) = %.5f\t \tkappa=%.5f\n", r->GetName(), r_1, nll_1-nll0, kappa);
477 if (!ok &&
picky_)
return std::numeric_limits<float>::quiet_NaN();
479 double nll_1_prof = nll.getVal();
480 kappa = (nll_1 - nll_1_prof) /
std::pow(r_1 - r_0,2);
481 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll unprof) = %.5f\tdelta(nll prof) = %.5f\tkappa=%.5f\n", r->GetName(), r_1, nll_1-nll0, nll.getVal()-nll0, kappa);
482 if (nll_1_prof > threshold) {
488 if (aboveThresh) rStep *= 2;
491 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll unprof) = %.5f\t \tkappa=%.5f\n", r->GetName(), r_1, nll_1-nll0, kappa);
493 if (r_1 > rMax0)
return std::numeric_limits<float>::quiet_NaN();
496 if (
verbose > 1) printf(
"At %s = %f:\t \tdelta(nll prof) = %.5f\tkappa=%.5f\n", r->GetName(), r_0, nll_0-nll0, kappa);
497 if (
verbose > 1) printf(
"At %s = %f:\t \tdelta(nll prof) = %.5f\tkappa=%.5f\n", r->GetName(), r_1, nll_1-nll0, kappa);
506 while (r_hi - r_lo > rEps) {
507 double r_2 = 0.5*(r_hi+r_lo);
509 double y0 = nll.getVal(),
y = y0 - kappa*
std::pow(r_2-r_1,2);
510 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll unprof) = %.5f\tdelta(nll appr) = %.5f\tkappa=%.5f\n", r->GetName(), r_2, y0-nll0,
y-nll0, kappa);
511 if (
y < threshold) { r_lo = r_2; }
else { r_hi = r_2; }
514 rCross = r->getVal();
515 double nll_unprof = nll.getVal();
521 if (!ok &&
picky_)
return std::numeric_limits<float>::quiet_NaN();
522 double nll_prof = nll.getVal();
523 if (
verbose > 1) printf(
"At %s = %f:\tdelta(nll unprof) = %.5f\tdelta(nll prof) = %.5f\tdelta(nll appr) = %.5f\n", r->GetName(), rCross, nll_unprof-nll0, nll_prof-nll0, nll_unprof-nll0 - kappa*
std::pow(rCross-r_1,2));
526 kappa = (nll_unprof - nll_prof)/
std::pow(rCross-r_1,2);
528 if ((nll_prof < threshold) == (nll_0 <
threshold)) {
534 r_1 = rCross; nll_1 = nll_prof;
537 if (minosStat != -1)
return rCross;
539 return std::numeric_limits<float>::quiet_NaN();
544 if (w->data(
"_Asymptotic_asimovDataset_") != 0) {
545 return w->data(
"_Asymptotic_asimovDataset_");
550 gobs.add(*mc_s->GetGlobalObservables());
558 asimovData->SetName(
"_Asymptotic_asimovDataset_");
559 w->import(*asimovData);
561 return w->data(
"_Asymptotic_asimovDataset_");
std::auto_ptr< RooArgSet > params_
static std::string minimizerAlgo_
static std::string minosAlgo_
bool setAllConstant(const RooAbsCollection &coll, bool constant=true)
set all RooRealVars to constants. return true if at least one changed status
RooArgSet snapGlobalObsAsimov
static int minimizerStrategy_
float findExpectedLimitFromCrossing(RooAbsReal &nll, RooRealVar *r, double rMin, double rMax, double nll0, double quantile)
virtual bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
bool minimize(int verbose=0, bool cascade=true)
static float minimizerTolerance_
void writeTo(const RooAbsCollection &) const
bool improve(int verbose=0, bool cascade=true)
utils::CheapValueSnapshot fitFixD_
virtual void applyOptions(const boost::program_options::variables_map &vm)
RooAbsData * asimovDatasetWithFit(RooStats::ModelConfig *mc, RooAbsData &realdata, RooAbsCollection &snapshot, double poiValue=0.0, int verbose=0)
Generate asimov dataset from best fit value of nuisance parameters, and fill in snapshot of correspon...
utils::CheapValueSnapshot fitFixA_
static void commitPoint(bool expected, float quantile)
Save a point into the output tree. Usually if expected = false, quantile should be set to -1 (except ...
std::auto_ptr< RooAbsReal > nllA_
const T & max(const T &a, const T &b)
static double rRelAccuracy_
RooAbsData * asimovDataset(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data)
static double rAbsAccuracy_
virtual void applyDefaultOptions()
void readFrom(const RooAbsCollection &)
utils::CheapValueSnapshot fitFreeA_
Setup Minimizer configuration on creation, reset the previous one on destruction. ...
RooAbsData * asimovDatasetNominal(RooStats::ModelConfig *mc, double poiValue=0.0, int verbose=0)
Generate asimov dataset from nominal value of nuisance parameters.
std::auto_ptr< RooAbsReal > nllD_
RooMinimizerOpt & minimizer()
char data[epos_bytes_allocation]
virtual bool runLimit(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
utils::CheapValueSnapshot fitFreeD_
void Print(const char *fmt) const
boost::program_options::options_description options_
void setErrorLevel(float errorLevel)
RooArgSet snapGlobalObsData
Power< A, B >::type pow(const A &a, const B &b)
void setStrategy(int strategy)
std::vector< std::pair< float, float > > runLimitExpected(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
double getCLs(RooRealVar &r, double rVal, bool getAlsoExpected=false, double *limit=0, double *limitErr=0)