Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RooGenContext.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 RooGenContext.cxx
19\class RooGenContext
20\ingroup Roofitcore
21
22Implements a universal generator context for all
23RooAbsPdf classes that do not have or need a specialized generator
24context. This generator context queries the input p.d.f which observables
25it can generate internally and delegates generation of those observables
26to the p.d.f if it deems that safe. The other observables are generated
27use a RooAcceptReject sampling technique.
28**/
29
30#include "RooMsgService.h"
31#include "Riostream.h"
32
33#include "RooGenContext.h"
34#include "RooAbsPdf.h"
35#include "RooDataSet.h"
36#include "RooRealIntegral.h"
37#include "RooAcceptReject.h"
38#include "RooRealVar.h"
39#include "RooDataHist.h"
40#include "RooErrorHandler.h"
41#include "RooNumGenConfig.h"
42#include "RooNumGenFactory.h"
43
44#include "TString.h"
45
46
47using std::endl, std::string, std::ostream;
48
49
50
51////////////////////////////////////////////////////////////////////////////////
52/// Initialize a new context for generating events with the specified
53/// variables, using the specified PDF model. A prototype dataset (if provided)
54/// is not cloned and still belongs to the caller. The contents and shape
55/// of this dataset can be changed between calls to generate() as long as the
56/// expected columns to be copied to the generated dataset are present.
57/// Any argument supplied in the forceDirect RooArgSet are always offered
58/// for internal generation to the p.d.f., even if this is deemed unsafe by
59/// the logic of RooGenContext.
60
62 const RooDataSet *prototype, const RooArgSet* auxProto,
63 bool verbose, const RooArgSet* forceDirect) :
64 RooAbsGenContext(model,vars,prototype,auxProto,verbose),
66{
67 cxcoutI(Generation) << "RooGenContext::ctor() setting up event generator context for p.d.f. " << model.GetName()
68 << " for generation of observable(s) " << vars ;
69 if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
70 if (auxProto && !auxProto->empty()) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
71 if (forceDirect && !forceDirect->empty()) ccxcoutI(Generation) << " with internal generation forced for observables " << *forceDirect ;
72 ccxcoutI(Generation) << std::endl ;
73
74
75 // Clone the model and all nodes that it depends on so that this context
76 // is independent of any existing objects.
77 RooArgSet nodes(model,model.GetName());
78 if (nodes.snapshot(_cloneSet, true)) {
79 coutE(Generation) << "RooGenContext::RooGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << std::endl ;
81 }
82
83 // Find the clone in the snapshot list
84 _pdfClone = static_cast<RooAbsPdf*>(_cloneSet.find(model.GetName()));
85
86 // Optionally fix RooAddPdf normalizations
87 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
88 RooArgSet fullNormSet(vars) ;
89 fullNormSet.add(*prototype->get(),true) ;
90 _pdfClone->fixAddCoefNormalization(fullNormSet) ;
91 }
92
93 // Analyze the list of variables to generate...
94 _isValid= true;
95 const RooAbsArg *arg = nullptr;
96 for(RooAbsArg * tmp : vars) {
97 if(!_isValid) break;
98
99 // is this argument derived?
100 if(!tmp->isFundamental()) {
101 coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << std::endl;
102 _isValid= false;
103 continue;
104 }
105 // lookup this argument in the cloned set of PDF dependents
106 arg= static_cast<RooAbsArg const*>(_cloneSet.find(tmp->GetName()));
107 if(nullptr == arg) {
108 coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
109 << "\" which will have uniform distribution" << std::endl;
110 _uniformVars.add(*tmp);
111 }
112 else {
113
114 // does the model depend on this variable directly, ie, like "x" in
115 // f(x) or f(x,g(x,y)) or even f(x,x) ?
116 const RooAbsArg *direct= arg ;
117 if (forceDirect==nullptr || !forceDirect->find(direct->GetName())) {
118 if (!_pdfClone->isDirectGenSafe(*arg)) {
119 cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << std::endl;
120 direct=nullptr ;
121 }
122 }
123
124 // does the model depend indirectly on this variable through an lvalue chain?
125 // otherwise, this variable will have to be generated with accept/reject
126 if(direct) {
127 _directVars.add(*arg);
128 } else {
129 _otherVars.add(*arg);
130 }
131 }
132 }
133
134 if(!isValid()) {
135 coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << std::endl;
136 return;
137 }
138
139 // At this point directVars are all variables that are safe to be generated directly
140 // otherVars are all variables that are _not_ safe to be generated directly
141 if (!_directVars.empty()) {
142 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << std::endl ;
143 }
144 if (!_otherVars.empty()) {
145 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << std::endl ;
146 }
147
148 // If PDF depends on prototype data, direct generator cannot use static initialization
149 // in initGenerator()
150 bool staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
151 if (!staticInitOK) {
152 cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << std::endl ;
153 }
154
155 // Can the model generate any of the direct variables itself?
156 RooArgSet generatedVars;
157 if (!_directVars.empty()) {
158 _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
159
160 cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
161 << generatedVars << " with configuration identifier " << _code << std::endl ;
162 } else {
163 _code = 0 ;
164 }
165
166 // Move variables which cannot be generated into the list to be generated with accept/reject
167 _directVars.remove(generatedVars);
169
170 // Update _directVars to only include variables that will actually be directly generated
171 _directVars.removeAll();
172 _directVars.add(generatedVars);
173
174 cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
175 if (!_directVars.empty()) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
176 if (!_directVars.empty() && !_otherVars.empty()) ccxcoutI(Generation) << " and" ;
177 if (!_otherVars.empty()) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
178 if (!_uniformVars.empty()) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
179 ccxcoutI(Generation) << std::endl ;
180
181 // initialize the accept-reject generator
182 RooArgSet depList;
183 _pdfClone->getObservables(&_theEvent, depList);
184 depList.remove(_otherVars);
185
186 TString nname(_pdfClone->GetName()) ;
187 nname.Append("_AccRej") ;
188 TString ntitle(_pdfClone->GetTitle()) ;
189 ntitle.Append(" (Accept/Reject)") ;
190
191
192 RooArgSet protoDeps;
193 model.getObservables(&_protoVars, protoDeps);
194
195
196 if (_protoVars.empty()) {
197
198 // No prototype variables
199
200 if(depList.empty()) {
201 // All variable are generated with accept-reject
202
203 // Check if PDF supports maximum finding
204 Int_t maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
205 if (maxFindCode != 0) {
206 coutI(Generation) << "RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
207 << "model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << std::endl ;
208 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&_theEvent) ;
209 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",maxVal) ;
210 cxcoutD(Generation) << "RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << std::endl ;
211 }
212 }
213
214 if (!_otherVars.empty()) {
215 _pdfClone->getVal(&vars) ; // WVE debug
216 _acceptRejectFunc= std::make_unique<RooRealIntegral>(nname,ntitle,*_pdfClone,depList,&vars);
217 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << std::endl ;
218 } else {
219 _acceptRejectFunc = nullptr;
220 }
221
222 } else {
223
224 // Generation _with_ prototype variable
225 depList.remove(_protoVars,true,true);
226 _acceptRejectFunc= std::unique_ptr<RooAbsReal>{_pdfClone->createIntegral(depList,vars)};
227 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << std::endl ;
228
229 // Check if PDF supports maximum finding for the entire phase space
230 RooArgSet allVars(_otherVars);
231 allVars.add(_directVars);
232 Int_t maxFindCode = _pdfClone->getMaxVal(allVars) ;
233 if (maxFindCode != 0) {
234 // Special case: PDF supports max-finding in otherVars, no need to scan other+proto space for maximum
235 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
236 << "model supports analytical maximum finding in the full phase space: "
237 << "can provide analytical pdf maximum to numeric generator" << std::endl ;
238 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&allVars);
239 _maxVar = std::make_unique<RooRealVar>("funcMax", "function maximum", maxVal);
240 } else {
241 maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
242 if (maxFindCode != 0) {
243 _updateFMaxPerEvent = maxFindCode ;
244 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
245 << "model supports analytical maximum finding in the variables which are not"
246 << " internally generated. Can provide analytical pdf maximum to numeric "
247 << "generator" << std::endl;
248 cxcoutD(Generation) << "RooGenContext::ctor() maximum value must be reevaluated for each "
249 << "event with configuration code " << maxFindCode << std::endl ;
250 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",1) ;
251 }
252 }
253
254 if (!_maxVar) {
255
256 // Regular case: First find maximum in other+proto space
257 RooArgSet otherAndProto(_otherVars) ;
258
259 otherAndProto.add(protoDeps) ;
260
261 if (!_otherVars.empty()) {
262
263 cxcoutD(Generation) << "RooGenContext::ctor() prototype data provided, observables are generated numerically no "
264 << "analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
265 << "of accept/reject observables plus prototype observables: " << otherAndProto << std::endl ;
266
267 // Calculate maximum in other+proto space if there are any accept/reject generated observables
268 std::unique_ptr<RooAbsNumGenerator> maxFinder{RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,otherAndProto,RooArgSet(_protoVars),
269 *model.getGeneratorConfig(),_verbose)};
270// RooAcceptReject maxFinder(*_acceptRejectFunc,otherAndProto,RooNumGenConfig::defaultConfig(),_verbose) ;
271 double max = maxFinder->getFuncMax() ;
272 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",max) ;
273
274 if (max==0) {
275 oocoutE(nullptr, Generation) << "RooGenContext::ctor(" << model.GetName()
276 << ") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
277 << "but chosen numeric generator (" << maxFinder->generatorName() << ") does not support maximum finding" << std::endl ;
278 throw string("RooGenContext::ctor()") ;
279 }
280
281 cxcoutD(Generation) << "RooGenContext::ctor() maximum function value found through initial sampling is " << max << std::endl ;
282 }
283 }
284
285 }
286
287 if (_acceptRejectFunc && !_otherVars.empty()) {
288 _generator = std::unique_ptr<RooAbsNumGenerator>{RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,_otherVars,RooArgSet(_protoVars),*model.getGeneratorConfig(),_verbose,_maxVar.get())};
289 cxcoutD(Generation) << "RooGenContext::ctor() creating MC sampling generator " << _generator->generatorName() << " from function for observables " << _otherVars << std::endl ;
290 //_generator= new RooAcceptReject(*_acceptRejectFunc,_otherVars,RooNumGenConfig::defaultConfig(),_verbose,_maxVar);
291 } else {
292 _generator = nullptr ;
293 }
294
296}
297
299
300////////////////////////////////////////////////////////////////////////////////
301/// Attach the cloned model to the event buffer we will be filling.
302
304{
305 _pdfClone->recursiveRedirectServers(args,false);
306 if (_acceptRejectFunc) {
307 _acceptRejectFunc->recursiveRedirectServers(args,false) ; // WVE DEBUG
308 }
309
310 // Attach the RooAcceptReject generator the event buffer
311 if (_generator) {
312 _generator->attachParameters(args) ;
313 }
314
315}
316
317
318
319////////////////////////////////////////////////////////////////////////////////
320/// Perform one-time initialization of the generator context
321
323{
324 attach(theEvent) ;
325
326 // Reset the cloned model's error counters.
327 _pdfClone->resetErrorCounters();
328
329 // Initialize the PDFs internal generator
330 if (!_directVars.empty()) {
331 cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << std::endl ;
332 _pdfClone->initGenerator(_code) ;
333 }
334}
335
336
337////////////////////////////////////////////////////////////////////////////////
338/// Generate one event. The 'remaining' integer is not used other than
339/// for printing messages
340
342{
343 if(!_otherVars.empty()) {
344 // call the accept-reject generator to generate its variables
345
346 if (_updateFMaxPerEvent!=0) {
347 double max = _pdfClone->maxVal(_updateFMaxPerEvent)/_pdfClone->getNorm(_otherVars) ;
348 cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << std::endl ;
349 _maxVar->setVal(max) ;
350 }
351
352 if (_generator) {
353 double resampleRatio(1) ;
354 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
355 if (resampleRatio<1) {
356 coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
357 << resampleRatio << " due to increased maximum weight" << std::endl ;
358 resampleData(resampleRatio) ;
359 }
360 if(nullptr == subEvent) {
361 coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << std::endl ;
362 return;
363 }
364 theEvent.assignValueOnly(*subEvent) ;
365 //theEvent= *subEvent;
366
367 }
368 }
369
370 // Use the model's optimized generator, if one is available.
371 // The generator writes directly into our local 'event' since we attached it above.
372 if(!_directVars.empty()) {
373 _pdfClone->generateEvent(_code);
374 }
375
376 // Generate uniform variables (non-dependents)
377 for(auto * uniVar : _uniformVars) {
378 auto * arglv = dynamic_cast<RooAbsLValue*>(uniVar);
379 if (!arglv) {
380 coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << std::endl ;
382 }
383 arglv->randomize() ;
384 }
385 theEvent.assign(_uniformVars) ;
386
387}
388
389
390////////////////////////////////////////////////////////////////////////////////
391/// Printing interface
392
393void RooGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
394{
395 RooAbsGenContext::printMultiline(os,content,verbose,indent);
396 os << indent << " --- RooGenContext --- " << std::endl ;
397 os << indent << "Using PDF ";
399
400 if(verbose) {
401 os << indent << "Use PDF generator for " << _directVars << std::endl ;
402 os << indent << "Use MC sampling generator " << (_generator ? _generator->generatorName() : "<none>") << " for " << _otherVars << std::endl ;
403 if (!_protoVars.empty()) {
404 os << indent << "Prototype observables are " << _protoVars << std::endl ;
405 }
406 }
407}
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define oocoutE(o, 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)
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, bool forceIfSizeOne=false)
Sets the value of any argument in our set that also appears in the other set.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
RooAbsArg * find(const char *name) const
Find object with given name in 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.
RooArgSet _protoVars
Prototype observables.
bool _verbose
Verbose messaging?
bool _isValid
Is context in valid state?
bool isValid() const
void resampleData(double &ratio)
Rescale existing output buffer with given ratio.
Abstract base class for objects that are lvalues, i.e.
virtual void randomize(const char *rangeName=nullptr)=0
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
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'.
static void softAbort()
Soft abort function that interrupts macro execution but doesn't kill ROOT.
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Printing interface.
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
Generate one event.
RooArgSet _otherVars
List of observables generated internally, randomly, and by accept/reject sampling.
Int_t _code
Internal generation code.
RooArgSet _uniformVars
void initGenerator(const RooArgSet &theEvent) override
Perform one-time initialization of the generator context.
~RooGenContext() override
std::unique_ptr< RooAbsNumGenerator > _generator
MC sampling generation engine.
RooArgSet _directVars
RooAbsPdf * _pdfClone
Clone of input p.d.f.
RooArgSet _cloneSet
Clone of all nodes of input p.d.f.
Int_t _updateFMaxPerEvent
If true, maximum p.d.f value needs to be recalculated for each event.
std::unique_ptr< RooRealVar > _maxVar
Variable holding maximum value of p.d.f.
std::unique_ptr< RooAbsReal > _acceptRejectFunc
Projection function to be passed to accept/reject sampler.
void attach(const RooArgSet &params) override
Attach the cloned model to the event buffer we will be filling.
RooGenContext(const RooAbsPdf &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false, const RooArgSet *forceDirect=nullptr)
Initialize a new context for generating events with the specified variables, using the specified PDF ...
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Basic string class.
Definition TString.h:138
TString & Append(const char *cs)
Definition TString.h:581