3 #include "../interface/BayesianToyMC.h"
4 #include "RooRealVar.h"
6 #include "RooUniform.h"
7 #include "RooWorkspace.h"
8 #include "RooStats/BayesianCalculator.h"
9 #include "RooStats/SimpleInterval.h"
10 #include "RooStats/ModelConfig.h"
11 #include "../interface/Combine.h"
12 #include "../interface/utils.h"
14 using namespace RooStats;
22 LimitAlgo(
"BayesianToyMC specific options")
26 (
"tries", boost::program_options::value<unsigned int>(&
tries_)->default_value(
tries_),
"Number of times to run the ToyMC on the same data")
27 (
"numIters,i", boost::program_options::value<int>(&
numIters_)->default_value(
numIters_),
"Number of iterations or calls used within iteration (0=ROOT Default)")
30 "set range of integration equal to this number of times the hinted limit")
36 std::cout <<
"BayesianToyMC: when running with no systematics, BayesianToyMC is identical to BayesianSimple." << std::endl;
42 bool BayesianToyMC::run(RooWorkspace *
w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &
data,
double &
limit,
double &limitErr,
const double *hint) {
43 RooArgSet poi(*mc_s->GetParametersOfInterest());
45 RooRealVar *
r =
dynamic_cast<RooRealVar *
>(poi.first());
46 if ((hint != 0) && (*hint > r->getMin())) {
49 double rMax = r->getMax();
51 std::auto_ptr<RooStats::ModelConfig> mc_noNuis(0);
53 mc_noNuis.reset(
new RooStats::ModelConfig(w));
54 mc_noNuis->SetPdf(*mc_s->GetPdf());
55 mc_noNuis->SetObservables(*mc_s->GetObservables());
56 mc_noNuis->SetParametersOfInterest(*mc_s->GetParametersOfInterest());
57 mc_noNuis->SetPriorPdf(*mc_s->GetPriorPdf());
58 mc_s = mc_noNuis.get();
60 std::auto_ptr<RooAbsPdf> nuisancePdf;
62 limit = 0; limitErr = 0;
bool rerun =
false;
63 for (
unsigned int i = 0;
i <
tries_; ++
i) {
64 BayesianCalculator bcalc(data, *mc_s);
65 bcalc.SetLeftSideTailFraction(0);
66 bcalc.SetConfidenceLevel(
cl);
71 bcalc.ForceNuisancePdf(*nuisancePdf);
74 std::auto_ptr<SimpleInterval> bcInterval(bcalc.GetInterval());
75 if (bcInterval.get() == 0)
return false;
76 double lim = bcInterval->UpperLimit();
78 if (lim >= 0.5*r->getMax()) {
79 std::cout <<
"Limit " << r->GetName() <<
" < " << lim <<
"; " << r->GetName() <<
" max < " << r->getMax() << std::endl;
80 if (r->getMax()/rMax > 20)
return false;
81 r->setMax(r->getMax()*2);
87 if (tries_ > 1 &&
verbose > 1)
std::cout <<
" - limit from try " <<
i <<
": " << lim << std::endl;
91 limitErr = (tries_ > 1 ?
std::sqrt((limitErr/tries_ - limit*limit)/((tries_-1)*tries_)) : 0);
93 std::cout <<
"\n -- BayesianToyMC -- " <<
"\n";
95 std::cout <<
"Limit: " << r->GetName() <<
" < " << limit <<
" +/- " << limitErr <<
" @ " <<
cl * 100 <<
"% CL" << std::endl;
97 std::cout <<
"Limit: " << r->GetName() <<
" < " << limit <<
" @ " <<
cl * 100 <<
"% CL" << std::endl;
static int numIters_
number of iterations for each toy mc computation
RooAbsPdf * makeNuisancePdf(RooStats::ModelConfig &model, const char *name="nuisancePdf")
static unsigned int tries_
number of toy mc computations to run
virtual void applyOptions(const boost::program_options::variables_map &vm)
static std::string integrationType_
numerical integration algorithm
char data[epos_bytes_allocation]
boost::program_options::options_description options_
virtual bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
static float hintSafetyFactor_
Safety factor for hint (integrate up to this number of times the hinted limit)