1 #include "../interface/ToyMCSamplerOpt.h"
2 #include "../interface/utils.h"
8 #include <RooSimultaneous.h>
9 #include <RooRealVar.h>
10 #include <RooProdPdf.h>
11 #include <RooPoisson.h>
12 #include <RooDataHist.h>
13 #include <RooDataSet.h>
14 #include <RooRandom.h>
15 #include <../interface/ProfilingTools.h>
18 ToyMCSampler(ts, ntoys),
19 globalObsPdf_(globalObsPdf),
20 globalObsValues_(0), globalObsIndex_(-1),
21 nuisValues_(0), nuisIndex_(-1),
24 if (!generateNuisances) fPriorNuisance = 0;
31 globalObsValues_(0), globalObsIndex_(-1),
39 globalObsValues_(0), globalObsIndex_(-1),
47 for (std::map<RooAbsPdf *, toymcoptutils::SimPdfGenInfo *>::iterator it =
genCache_.begin(), ed =
genCache_.end(); it != ed; ++it) {
51 delete _allVars; _allVars = 0;
57 mode_(pdf.canBeExtended() ? (preferBinned ? Binned : Unbinned) : Counting),
59 spec_(0),histoSpec_(0),keepHistoSpec_(0),weightVar_(0)
61 if (pdf.canBeExtended()) {
62 if (pdf.getAttribute(
"forceGenBinned"))
mode_ =
Binned;
63 else if (pdf.getAttribute(
"forceGenPoisson"))
mode_ =
Poisson;
64 else if (pdf.getAttribute(
"forceGenUnbinned"))
mode_ =
Unbinned;
67 RooArgSet *obs = pdf.getObservables(observables);
93 assert(forceEvents == 0 &&
"SinglePdfGenInfo: forceEvents must be zero at least for now");
97 if (spec_ == 0) spec_ = protoData ? pdf_->prepareMultiGen(observables_, RooFit::Extended(), RooFit::ProtoData(*protoData,
true,
true))
98 : pdf_->prepareMultiGen(observables_, RooFit::Extended());
99 if (spec_) ret = pdf_->generate(*spec_);
100 else ret = pdf_->generate(observables_, RooFit::Extended());
104 RooDataSet *
data = pdf_->generate(observables_, RooFit::Extended());
105 ret =
new RooDataHist(data->GetName(),
"", *data->get(), *
data);
109 case BinnedNoWorkaround:
110 ret = protoData ? pdf_->generateBinned(observables_, RooFit::Extended(), RooFit::ProtoData(*protoData,
true,
true))
111 : pdf_->generateBinned(observables_, RooFit::Extended());
114 ret = generateWithHisto(weightVar_,
false);
117 ret = pdf_->generate(observables_, 1);
120 throw std::logic_error(
"Mode not foreseen in SinglePdfGenInfo::generate");
131 if (nPA)
return generatePseudoAsimov(weightVar, nPA);
132 return generateWithHisto(weightVar,
true);
138 if (mode_ == Unbinned) {
139 double expEvents = pdf_->expectedEvents(observables_);
140 std::auto_ptr<RooDataSet>
data(pdf_->generate(observables_, nPoints));
141 if (weightVar == 0) weightVar =
new RooRealVar(
"_weight_",
"",1.0);
142 RooArgSet obsPlusW(observables_); obsPlusW.add(*weightVar);
143 RooDataSet *rds =
new RooDataSet(
data->GetName(),
"", obsPlusW, weightVar->GetName());
144 for (
int i = 0;
i < nPoints; ++
i) {
145 observables_ = *
data->get(
i);
146 rds->add(observables_, expEvents/nPoints);
150 return generateWithHisto(weightVar,
true);
158 if (mode_ == Counting)
return generateCountingAsimov();
159 if (observables_.getSize() > 3)
throw std::invalid_argument(std::string(
"ERROR in SinglePdfGenInfo::generateWithHisto for ") + pdf_->GetName() +
", more than 3 observable");
160 RooArgList obs(observables_);
161 RooRealVar *
x = (RooRealVar*)obs.at(0);
162 RooRealVar *
y = obs.getSize() > 1 ? (RooRealVar*)obs.at(1) : 0;
163 RooRealVar *
z = obs.getSize() > 2 ? (RooRealVar*)obs.at(2) : 0;
164 if (weightVar == 0) weightVar =
new RooRealVar(
"_weight_",
"",1.0);
169 if (histoSpec_ == 0) {
170 histoSpec_ = pdf_->createHistogram(
"htemp", *x, ay, az);
171 histoSpec_->SetDirectory(0);
174 double expectedEvents = pdf_->expectedEvents(observables_);
175 histoSpec_->Scale(expectedEvents/ histoSpec_->Integral());
176 RooArgSet obsPlusW(obs); obsPlusW.add(*weightVar);
177 RooDataSet *
data =
new RooDataSet(TString::Format(
"%sData", pdf_->GetName()),
"", obsPlusW, weightVar->GetName());
178 RooAbsArg::setDirtyInhibit(
true);
179 switch (obs.getSize()) {
181 for (
int i = 1,
n = histoSpec_->GetNbinsX();
i <=
n; ++
i) {
182 x->setVal(histoSpec_->GetXaxis()->GetBinCenter(
i));
183 data->add(observables_, asimov ? histoSpec_->GetBinContent(
i) : RooRandom::randomGenerator()->Poisson(histoSpec_->GetBinContent(
i)) );
188 TH2& h2 =
dynamic_cast<TH2&
>(*histoSpec_);
189 for (
int ix = 1, nx = h2.GetNbinsX(); ix <= nx; ++ix) {
190 for (
int iy = 1, ny = h2.GetNbinsY(); iy <= ny; ++iy) {
191 x->setVal(h2.GetXaxis()->GetBinCenter(ix));
192 y->setVal(h2.GetYaxis()->GetBinCenter(iy));
193 data->add(observables_, asimov ? h2.GetBinContent(ix,iy) : RooRandom::randomGenerator()->Poisson(h2.GetBinContent(ix,iy)) );
199 TH3& h3 =
dynamic_cast<TH3&
>(*histoSpec_);
200 for (
int ix = 1, nx = h3.GetNbinsX(); ix <= nx; ++ix) {
201 for (
int iy = 1, ny = h3.GetNbinsY(); iy <= ny; ++iy) {
202 for (
int iz = 1, nz = h3.GetNbinsZ(); iz <= nz; ++iz) {
203 x->setVal(h3.GetXaxis()->GetBinCenter(ix));
204 y->setVal(h3.GetYaxis()->GetBinCenter(iy));
205 z->setVal(h3.GetYaxis()->GetBinCenter(iz));
206 data->add(observables_, asimov ? h3.GetBinContent(ix,iy,iz) : RooRandom::randomGenerator()->Poisson(h3.GetBinContent(ix,iy,iz)) );
210 RooAbsArg::setDirtyInhibit(
false);
211 if (!keepHistoSpec_) {
delete histoSpec_; histoSpec_ = 0; }
221 RooArgSet obs(observables_);
222 RooProdPdf *
prod =
dynamic_cast<RooProdPdf *
>(pdf_);
223 RooPoisson *pois = 0;
225 setToExpected(*prod, observables_);
226 }
else if ((pois = dynamic_cast<RooPoisson *>(pdf_)) != 0) {
227 setToExpected(*pois, observables_);
228 }
else throw std::logic_error(
"A counting model pdf must be either a RooProdPdf or a RooPoisson");
229 RooDataSet *
ret =
new RooDataSet(TString::Format(
"%sData", pdf_->GetName()),
"", obs);
237 std::auto_ptr<TIterator> iter(prod.pdfList().createIterator());
238 for (RooAbsArg *
a = (RooAbsArg *) iter->Next();
a != 0;
a = (RooAbsArg *) iter->Next()) {
239 if (!
a->dependsOn(obs))
continue;
240 RooPoisson *pois = 0;
241 if ((pois = dynamic_cast<RooPoisson *>(
a)) != 0) {
242 setToExpected(*pois, obs);
244 RooProdPdf *subprod =
dynamic_cast<RooProdPdf *
>(
a);
245 if (subprod) setToExpected(*subprod, obs);
246 else throw std::logic_error(
"Illegal term in counting model: depends on observables, but not Poisson or Product");
254 RooRealVar *myobs = 0;
255 RooAbsReal *myexp = 0;
256 std::auto_ptr<TIterator> iter(pois.serverIterator());
257 for (RooAbsArg *
a = (RooAbsArg *) iter->Next();
a != 0;
a = (RooAbsArg *) iter->Next()) {
258 if (obs.contains(*
a)) {
259 assert(myobs == 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): Two observables??");
260 myobs =
dynamic_cast<RooRealVar *
>(
a);
261 assert(myobs != 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): Observables is not a RooRealVar??");
263 assert(myexp == 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): Two expecteds??");
264 myexp =
dynamic_cast<RooAbsReal *
>(
a);
265 assert(myexp != 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): Expectedis not a RooAbsReal??");
268 assert(myobs != 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): No observable?");
269 assert(myexp != 0 &&
"SinglePdfGenInfo::setToExpected(RooPoisson): No expected?");
270 myobs->setVal(myexp->getVal());
276 observables_(observables),
279 assert(forceEvents == 0 &&
"SimPdfGenInfo: forceEvents must be zero at least for now");
280 RooSimultaneous *simPdf =
dynamic_cast<RooSimultaneous *
>(&pdf);
282 cat_ =
const_cast<RooAbsCategoryLValue *
>(&simPdf->indexCat());
283 int nbins =
cat_->numBins((
const char *)0);
284 pdfs_.resize(nbins, 0);
286 for (
int ic = 0; ic <
nbins; ++ic) {
288 RooAbsPdf *pdfi = simPdf->getPdf(
cat_->getLabel());
289 if (pdfi == 0)
throw std::logic_error(std::string(
"Unmapped category state: ") +
cat_->getLabel());
292 if (newpdf != 0 && newpdf != pdfi) {
303 for (std::vector<SinglePdfGenInfo *>::iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
308 for (std::map<std::string,RooAbsData*>::iterator it = datasetPieces_.begin(), ed = datasetPieces_.end(); it != ed; ++it) {
311 datasetPieces_.clear();
319 TString retName = TString::Format(
"%sData", pdf_->GetName());
322 for (
int i = 0,
n = cat_->numBins((
const char *)0);
i <
n; ++
i) {
323 if (pdfs_[
i] == 0)
continue;
325 RooAbsData *&
data = datasetPieces_[cat_->getLabel()];
delete data;
326 assert(protoData == 0);
327 data = pdfs_[
i]->generate(protoData);
328 if (data->isWeighted()) {
329 if (weightVar == 0) weightVar =
new RooRealVar(
"_weight_",
"",1.0);
330 RooArgSet obs(*data->get());
332 RooDataSet *wdata =
new RooDataSet(data->GetName(),
"", obs,
"_weight_");
333 for (
int i = 0, n = data->numEntries();
i <
n; ++
i) {
335 if (data->weight()) wdata->add(obs, data->weight());
346 RooArgSet vars(observables_), varsPlusWeight(observables_);
347 if (weightVar) varsPlusWeight.add(*weightVar);
348 ret =
new RooDataSet(retName,
"", varsPlusWeight, (weightVar ? weightVar->GetName() : 0));
349 for (std::map<std::string,RooAbsData*>::iterator it = datasetPieces_.begin(), ed = datasetPieces_.end(); it != ed; ++it) {
350 cat_->setLabel(it->first.c_str());
351 for (
unsigned int i = 0, n = it->second->numEntries();
i <
n; ++
i) {
352 vars = *it->second->get(
i);
353 ret->add(vars, it->second->weight());
360 ret =
new RooDataSet(retName,
"", observables_, RooFit::Index((RooCategory&)*cat_), RooFit::Link(datasetPieces_) );
362 }
else ret = pdfs_[0]->generate(protoData, forceEvents);
371 TString retName = TString::Format(
"%sData", pdf_->GetName());
374 for (
int i = 0,
n = cat_->numBins((
const char *)0);
i <
n; ++
i) {
375 if (pdfs_[
i] == 0)
continue;
377 RooAbsData *&
data = datasetPieces_[cat_->getLabel()];
delete data;
378 data = pdfs_[
i]->generateAsimov(weightVar);
381 RooArgSet vars(observables_), varsPlusWeight(observables_); varsPlusWeight.add(*weightVar);
382 ret =
new RooDataSet(retName,
"", varsPlusWeight, (weightVar ? weightVar->GetName() : 0));
383 for (std::map<std::string,RooAbsData*>::iterator it = datasetPieces_.begin(), ed = datasetPieces_.end(); it != ed; ++it) {
384 cat_->setLabel(it->first.c_str());
385 for (
unsigned int i = 0, n = it->second->numEntries();
i <
n; ++
i) {
386 vars = *it->second->get(
i);
387 ret->add(vars, it->second->weight());
394 ret =
new RooDataSet(retName,
"", observables_, RooFit::Index((RooCategory&)*cat_), RooFit::Link(datasetPieces_) );
396 }
else ret = pdfs_[0]->generateAsimov(weightVar);
405 for (std::vector<SinglePdfGenInfo *>::const_iterator it = pdfs_.begin(), ed = pdfs_.end(); it != ed; ++it) {
406 if (*it) (*it)->setCacheTemplates(cache);
413 ToyMCSampler::SetPdf(pdf);
414 delete _allVars; _allVars = 0;
415 delete globalObsValues_; globalObsValues_ = 0; globalObsIndex_ = -1;
416 delete nuisValues_; nuisValues_ = 0; nuisIndex_ = -1;
424 if (fObservables ==
NULL) {
425 ooccoutE((TObject*)
NULL,InputArguments) <<
"Observables not set." << endl;
431 if(fPriorNuisance && fNuisancePars && fNuisancePars->getSize() > 0) {
432 if (nuisValues_ == 0 || nuisIndex_ == nuisValues_->numEntries()) {
434 nuisValues_ = fPriorNuisance->generate(*fNuisancePars, fNToys);
437 fNuisancePars->snapshot(saveNuis);
438 const RooArgSet *
values = nuisValues_->get(nuisIndex_++);
439 RooArgSet pars(*fNuisancePars); pars = *
values;
442 RooArgSet observables(*fObservables);
443 if(fGlobalObservables && fGlobalObservables->getSize()) {
444 observables.remove(*fGlobalObservables);
447 assert(globalObsPdf_);
448 if (globalObsValues_ == 0 || globalObsIndex_ == globalObsValues_->numEntries()) {
449 delete globalObsValues_;
450 globalObsValues_ = (globalObsPdf_ ? globalObsPdf_ : fPdf)->
generate(*fGlobalObservables, fNToys);
453 const RooArgSet *
values = globalObsValues_->get(globalObsIndex_++);
454 if (!_allVars) _allVars = fPdf->getObservables(*fGlobalObservables);
460 if(!fImportanceDensity) {
462 data = Generate(*fPdf, observables);
464 throw std::runtime_error(
"No importance sampling yet");
466 RooArgSet* allVars = fPdf->getVariables();
467 RooArgSet* allVars2 = fImportanceDensity->getVariables();
468 allVars->add(*allVars2);
469 const RooArgSet* saveVars = (
const RooArgSet*)allVars->snapshot();
476 forceEvents = (int)fPdf->expectedEvents(observables);
477 forceEvents = RooRandom::randomGenerator()->Poisson(forceEvents);
482 if(fImportanceSnapshot) *allVars = *fImportanceSnapshot;
487 data = Generate(*fImportanceDensity, observables,
NULL, forceEvents);
489 *allVars = *saveVars;
495 if (saveNuis.getSize()) { RooArgSet pars(*fNuisancePars); pars = saveNuis; }
504 protoData = fProtoData;
505 forceEvents = protoData->numEntries();
508 if (events == 0) events = fNEvents;
511 assert(protoData == 0);
512 RooAbsData *
ret = pdf.generate(observables, events);
521 return info->
generate(weightVar_, protoData, forceEvents);
RooAbsCategoryLValue * cat_
ToyMCSamplerOpt(RooStats::TestStatistic &ts, Int_t ntoys, RooAbsPdf *globalObsPdf=0, bool generateNuisances=false)
RooAbsData * Generate(RooAbsPdf &pdf, RooArgSet &observables, const RooDataSet *protoData=NULL, int forceEvents=0) const
int get(const char *name)
RooAbsData * generateAsimov(RooRealVar *&weightVar)
RooAbsData * generate(const RooDataSet *protoData=NULL, int forceEvents=0)
RooDataSet * globalObsValues_
virtual RooAbsData * GenerateToyData(RooArgSet &, double &weight) const
RooDataSet * generatePseudoAsimov(RooRealVar *&weightVar, int nPoints)
void setToExpected(RooProdPdf &prod, RooArgSet &obs)
std::map< RooAbsPdf *, toymcoptutils::SimPdfGenInfo * > genCache_
RooAbsPdf * factorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &constraints)
RooDataSet * generateCountingAsimov()
void setCacheTemplates(bool cache)
RooAbsData * generate(RooRealVar *&weightVar, const RooDataSet *protoData=NULL, int forceEvents=0)
virtual void SetPdf(RooAbsPdf &pdf)
RooDataSet * generateWithHisto(RooRealVar *&weightVar, bool asimov)
void setCopyData(bool copyData)
RooDataSet * generateAsimov(RooRealVar *&weightVar)
std::vector< SinglePdfGenInfo * > pdfs_
SimPdfGenInfo(RooAbsPdf &pdf, const RooArgSet &observables, bool preferBinned, const RooDataSet *protoData=NULL, int forceEvents=0)
char data[epos_bytes_allocation]
SinglePdfGenInfo(RooAbsPdf &pdf, const RooArgSet &observables, bool preferBinned, const RooDataSet *protoData=NULL, int forceEvents=0)