Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RooBinnedGenContext.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooBinnedGenContext.cxx
19\class RooBinnedGenContext
20\ingroup Roofitcore
21
22Efficient implementation of the generator context specific for binned pdfs.
23**/
24
25#include "Riostream.h"
26
27
28#include "RooMsgService.h"
29#include "RooBinnedGenContext.h"
30#include "RooAbsPdf.h"
31#include "RooRealVar.h"
32#include "RooDataHist.h"
33#include "RooDataSet.h"
34#include "RooRandom.h"
35
36using std::endl, std::vector, std::ostream;
37
38
39
40////////////////////////////////////////////////////////////////////////////////
41/// Constructor
42
44 const RooDataSet *prototype, const RooArgSet* auxProto,
45 bool verbose) :
46 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
47{
48 cxcoutI(Generation) << "RooBinnedGenContext::ctor() setting up event special generator context for sum p.d.f. " << model.GetName()
49 << " for generation of observable(s) " << vars ;
50 if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
51 if (auxProto && !auxProto->empty()) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
52 ccxcoutI(Generation) << std::endl ;
53
54 // Constructor. Build an array of generator contexts for each product component PDF
55 RooArgSet(model).snapshot(_pdfSet, true);
56 _pdf = static_cast<RooAbsPdf*>(_pdfSet.find(model.GetName())) ;
57 _pdf->setOperMode(RooAbsArg::ADirty,true) ;
58
59 // Fix normalization set of this RooAddPdf
60 if (prototype)
61 {
62 RooArgSet coefNSet(vars) ;
63 coefNSet.add(*prototype->get(),true) ;
64 _pdf->fixAddCoefNormalization(coefNSet) ;
65 }
66
67 _pdf->recursiveRedirectServers(_theEvent) ;
68 _vars = std::unique_ptr<RooArgSet>{_pdf->getObservables(vars)};
69
70 // Create empty RooDataHist
71 _hist = std::make_unique<RooDataHist>("genData","genData",*_vars);
72
73 _expectedData = false ;
74}
75
76
78
79
80////////////////////////////////////////////////////////////////////////////////
81/// Attach given set of variables to internal p.d.f. clone
82
84{
85 _pdf->recursiveRedirectServers(args) ;
86}
87
88
89
90////////////////////////////////////////////////////////////////////////////////
91/// One-time initialization of generator context. Attach theEvent
92/// to internal p.d.f clone and forward initialization call to
93/// the component generators.
94
96{
97 _pdf->recursiveRedirectServers(theEvent) ;
98
99}
100
101
102////////////////////////////////////////////////////////////////////////////////
103
105{
106 _expectedData = flag ;
107}
108
109
110////////////////////////////////////////////////////////////////////////////////
111
112RooDataSet *RooBinnedGenContext::generate(double nEvt, bool /*skipInit*/, bool extended)
113{
114 // Scale to number of events and introduce Poisson fluctuations
115 _hist->reset() ;
116
117 double nEvents = nEvt ;
118
119 if (nEvents<=0) {
120 if (!_pdf->canBeExtended()) {
121 coutE(InputArguments) << "RooAbsPdf::generateBinned(" << GetName()
122 << ") ERROR: No event count provided and p.d.f does not provide expected number of events" << std::endl ;
123 return nullptr ;
124 } else {
125 // Don't round in expectedData mode
126 if (_expectedData || extended) {
127 nEvents = _pdf->expectedEvents(_vars.get());
128 } else {
129 nEvents = Int_t(_pdf->expectedEvents(_vars.get())+0.5) ;
130 }
131 }
132 }
133
134 // Sample p.d.f. distribution
135 _pdf->fillDataHist(_hist.get(),_vars.get(),1,true) ;
136
137 // Output container
138 RooDataSet* wudata = new RooDataSet("wu","wu",*_vars,RooFit::WeightVar()) ;
139
140 vector<int> histOut(_hist->numEntries()) ;
141 double histMax(-1) ;
142 Int_t histOutSum(0) ;
143 for (int i=0 ; i<_hist->numEntries() ; i++) {
144 const RooArgSet* coords = _hist->get(i) ;
145 const double wi = _hist->weight(i) ;
146 if (_expectedData) {
147
148 // Expected data, multiply p.d.f by nEvents
149 double w=wi*nEvents ;
150 wudata->add(*coords,w) ;
151
152 } else if (extended) {
153
154 // Extended mode, set contents to Poisson(pdf*nEvents)
155 double w = RooRandom::randomGenerator()->Poisson(wi*nEvents) ;
156 wudata->add(*coords,w) ;
157
158 } else {
159
160 // Regular mode, fill array of weights with Poisson(pdf*nEvents), but to not fill
161 // histogram yet.
162 if (wi>histMax) {
163 histMax = wi ;
164 }
165 histOut[i] = RooRandom::randomGenerator()->Poisson(wi*nEvents) ;
166 histOutSum += histOut[i] ;
167 }
168 }
169
170
171 if (!_expectedData && !extended) {
172
173 // Second pass for regular mode - Trim/Extend dataset to exact number of entries
174
175 // Calculate difference between what is generated so far and what is requested
176 Int_t nEvtExtra = std::abs(Int_t(nEvents)-histOutSum) ;
177 Int_t wgt = (histOutSum>nEvents) ? -1 : 1 ;
178
179 // Perform simple binned accept/reject procedure to get to exact event count
180 while(nEvtExtra>0) {
181
182 Int_t ibinRand = RooRandom::randomGenerator()->Integer(_hist->numEntries()) ;
183 _hist->get(ibinRand) ;
184 double ranY = RooRandom::randomGenerator()->Uniform(histMax) ;
185
186 if (ranY<_hist->weight(ibinRand)) {
187 if (wgt==1) {
188 histOut[ibinRand]++ ;
189 } else {
190 // If weight is negative, prior bin content must be at least 1
191 if (histOut[ibinRand]>0) {
192 histOut[ibinRand]-- ;
193 } else {
194 continue ;
195 }
196 }
197 nEvtExtra-- ;
198 }
199 }
200
201 // Transfer working array to histogram
202 for (int i=0 ; i<_hist->numEntries() ; i++) {
203 _hist->get(i) ;
204 wudata->add(*_hist->get(),histOut[i]) ;
205 }
206
207 }
208
209 return wudata ;
210
211}
212
213
214
215////////////////////////////////////////////////////////////////////////////////
216/// this method is not implemented for this context
217
219{
220 assert(0) ;
221}
222
223
224
225////////////////////////////////////////////////////////////////////////////////
226/// Print the details of the context
227
228void RooBinnedGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
229{
230 RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
231 os << indent << "--- RooBinnedGenContext ---" << std::endl ;
232 os << indent << "Using PDF ";
233 _pdf->printStream(os,kName|kArgs|kClassName,kSingleLine,indent);
234}
#define cxcoutI(a)
#define coutE(a)
#define ccxcoutI(a)
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
static void indent(ostringstream &buf, int indent_level)
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooArgSet _theEvent
Pointer to observable event being generated.
RooAbsGenContext(const RooAbsPdf &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool _verbose=false)
Constructor.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Interface for multi-line printing.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
RooAbsPdf * _pdf
Pointer to cloned p.d.f.
~RooBinnedGenContext() override
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
this method is not implemented for this context
std::unique_ptr< RooDataHist > _hist
Histogram.
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Print the details of the context.
RooArgSet _pdfSet
Set owned all nodes of internal clone of p.d.f.
RooBinnedGenContext(const RooAbsPdf &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool _verbose=false)
Constructor.
void setExpectedData(bool) override
std::unique_ptr< RooArgSet > _vars
void initGenerator(const RooArgSet &theEvent) override
One-time initialization of generator context.
RooDataSet * generate(double nEvents=0.0, bool skipInit=false, bool extendedMode=false) override
Generate the specified number of events with nEvents>0 and and return a dataset containing the genera...
void attach(const RooArgSet &params) override
Attach given set of variables to internal p.d.f. clone.
Container class to hold unbinned data.
Definition RooDataSet.h:32
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual ULong64_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition TRandom.cxx:403
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:681
virtual UInt_t Integer(UInt_t imax)
Returns a random integer uniformly distributed on the interval [ 0, imax-1 ].
Definition TRandom.cxx:360
Basic string class.
Definition TString.h:138
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)