CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BayesianToyMC.cc
Go to the documentation of this file.
1 #include <stdexcept>
2 #include <cmath>
3 #include "../interface/BayesianToyMC.h"
4 #include "RooRealVar.h"
5 #include "RooArgSet.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"
13 
14 using namespace RooStats;
15 
16 int BayesianToyMC::numIters_ = 1000;
17 std::string BayesianToyMC::integrationType_ = "toymc";
18 unsigned int BayesianToyMC::tries_ = 1;
20 
22  LimitAlgo("BayesianToyMC specific options")
23 {
24  options_.add_options()
25  ("integrationType", boost::program_options::value<std::string>(&integrationType_)->default_value(integrationType_), "Integration algorithm to use")
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)")
28  ("hintSafetyFactor",
29  boost::program_options::value<float>(&hintSafetyFactor_)->default_value(hintSafetyFactor_),
30  "set range of integration equal to this number of times the hinted limit")
31  ;
32 }
33 
34 void BayesianToyMC::applyOptions(const boost::program_options::variables_map &vm) {
35  if (!withSystematics) {
36  std::cout << "BayesianToyMC: when running with no systematics, BayesianToyMC is identical to BayesianSimple." << std::endl;
37  tries_ = 1;
38  numIters_ = 0;
39  integrationType_ = "";
40  }
41 }
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());
44 
45  RooRealVar *r = dynamic_cast<RooRealVar *>(poi.first());
46  if ((hint != 0) && (*hint > r->getMin())) {
47  r->setMax(hintSafetyFactor_*(*hint));
48  }
49  double rMax = r->getMax();
50 
51  std::auto_ptr<RooStats::ModelConfig> mc_noNuis(0);
52  if (!withSystematics && mc_s->GetNuisanceParameters() != 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();
59  }
60  std::auto_ptr<RooAbsPdf> nuisancePdf;
61  for (;;) {
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);
67  if (!integrationType_.empty()) bcalc.SetIntegrationType(integrationType_.c_str());
68  if (numIters_) bcalc.SetNumIters(numIters_);
69  if (integrationType_ == "toymc") {
70  if (nuisancePdf.get() == 0) nuisancePdf.reset(utils::makeNuisancePdf(*mc_s));
71  bcalc.ForceNuisancePdf(*nuisancePdf);
72  }
73  // get the interval
74  std::auto_ptr<SimpleInterval> bcInterval(bcalc.GetInterval());
75  if (bcInterval.get() == 0) return false;
76  double lim = bcInterval->UpperLimit();
77  // check against bound
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);
82  rerun = true; break;
83  }
84  // add to running sum(x) and sum(x2)
85  limit += lim;
86  limitErr += lim*lim;
87  if (tries_ > 1 && verbose > 1) std::cout << " - limit from try " << i << ": " << lim << std::endl;
88  }
89  if (rerun) continue;
90  limit /= tries_;
91  limitErr = (tries_ > 1 ? std::sqrt((limitErr/tries_ - limit*limit)/((tries_-1)*tries_)) : 0);
92  if (verbose > -1) {
93  std::cout << "\n -- BayesianToyMC -- " << "\n";
94  if (limitErr > 0) {
95  std::cout << "Limit: " << r->GetName() << " < " << limit << " +/- " << limitErr << " @ " << cl * 100 << "% CL" << std::endl;
96  } else {
97  std::cout << "Limit: " << r->GetName() << " < " << limit << " @ " << cl * 100 << "% CL" << std::endl;
98  }
99  }
100  break;
101  }
102  return true;
103 }
int i
Definition: DBlmapReader.cc:9
bool withSystematics
Definition: Combine.cc:68
static int numIters_
number of iterations for each toy mc computation
Definition: BayesianToyMC.h:26
RooAbsPdf * makeNuisancePdf(RooStats::ModelConfig &model, const char *name="nuisancePdf")
Definition: utils.cc:200
static unsigned int tries_
number of toy mc computations to run
Definition: BayesianToyMC.h:28
T sqrt(T t)
Definition: SSEVec.h:46
float cl
Definition: Combine.cc:71
virtual void applyOptions(const boost::program_options::variables_map &vm)
static std::string integrationType_
numerical integration algorithm
Definition: BayesianToyMC.h:24
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
boost::program_options::options_description options_
Definition: LimitAlgo.h:31
virtual bool run(RooWorkspace *w, RooStats::ModelConfig *mc_s, RooStats::ModelConfig *mc_b, RooAbsData &data, double &limit, double &limitErr, const double *hint)
T w() const
static float hintSafetyFactor_
Safety factor for hint (integrate up to this number of times the hinted limit)
Definition: BayesianToyMC.h:30