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
63 bool verbose, const RooArgSet* forceDirect) :
64 RooAbsGenContext(model,vars,prototype,auxProto,verbose),
65 _updateFMaxPerEvent(0)
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()));
86
87 // Optionally fix RooAddPdf normalizations
88 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
90 fullNormSet.add(*prototype->get()) ;
92 }
93
94 // Analyze the list of variables to generate...
95 _isValid= true;
96 const RooAbsArg *arg = nullptr;
97 for(RooAbsArg * tmp : vars) {
98 if(!_isValid) break;
99
100 // is this argument derived?
101 if(!tmp->isFundamental()) {
102 coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << std::endl;
103 _isValid= false;
104 continue;
105 }
106 // lookup this argument in the cloned set of PDF dependents
107 arg= static_cast<RooAbsArg const*>(_cloneSet.find(tmp->GetName()));
108 if(nullptr == arg) {
109 coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
110 << "\" which will have uniform distribution" << std::endl;
112 }
113 else {
114
115 // does the model depend on this variable directly, ie, like "x" in
116 // f(x) or f(x,g(x,y)) or even f(x,x) ?
117 const RooAbsArg *direct= arg ;
118 if (forceDirect==nullptr || !forceDirect->find(direct->GetName())) {
119 if (!_pdfClone->isDirectGenSafe(*arg)) {
120 cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << std::endl;
121 direct=nullptr ;
122 }
123 }
124
125 // does the model depend indirectly on this variable through an lvalue chain?
126 // otherwise, this variable will have to be generated with accept/reject
127 if(direct) {
128 _directVars.add(*arg);
129 } else {
130 _otherVars.add(*arg);
131 }
132 }
133 }
134
135 if(!isValid()) {
136 coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << std::endl;
137 return;
138 }
139
140 // At this point directVars are all variables that are safe to be generated directly
141 // otherVars are all variables that are _not_ safe to be generated directly
142 if (!_directVars.empty()) {
143 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << std::endl ;
144 }
145 if (!_otherVars.empty()) {
146 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << std::endl ;
147 }
148
149 // If PDF depends on prototype data, direct generator cannot use static initialization
150 // in initGenerator()
152 if (!staticInitOK) {
153 cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << std::endl ;
154 }
155
156 // Can the model generate any of the direct variables itself?
158 if (!_directVars.empty()) {
160
161 cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
162 << generatedVars << " with configuration identifier " << _code << std::endl ;
163 } else {
164 _code = 0 ;
165 }
166
167 // Move variables which cannot be generated into the list to be generated with accept/reject
170
171 // Update _directVars to only include variables that will actually be directly generated
174
175 cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
176 if (!_directVars.empty()) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
177 if (!_directVars.empty() && !_otherVars.empty()) ccxcoutI(Generation) << " and" ;
178 if (!_otherVars.empty()) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
179 if (!_uniformVars.empty()) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
180 ccxcoutI(Generation) << std::endl ;
181
182 // initialize the accept-reject generator
185 depList.remove(_otherVars);
186
188 nname.Append("_AccRej") ;
190 ntitle.Append(" (Accept/Reject)") ;
191
192
195
196
197 if (_protoVars.empty()) {
198
199 // No prototype variables
200
201 if(depList.empty()) {
202 // All variable are generated with accept-reject
203
204 // Check if PDF supports maximum finding
206 if (maxFindCode != 0) {
207 coutI(Generation) << "RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
208 << "model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << std::endl ;
209 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&_theEvent) ;
210 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",maxVal) ;
211 cxcoutD(Generation) << "RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << std::endl ;
212 }
213 }
214
215 if (!_otherVars.empty()) {
216 _pdfClone->getVal(&vars) ; // WVE debug
217 _acceptRejectFunc= std::make_unique<RooRealIntegral>(nname,ntitle,*_pdfClone,depList,&vars);
218 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << std::endl ;
219 } else {
220 _acceptRejectFunc = nullptr;
221 }
222
223 } else {
224
225 // Generation _with_ prototype variable
226 depList.remove(_protoVars,true,true);
227 _acceptRejectFunc= std::unique_ptr<RooAbsReal>{_pdfClone->createIntegral(depList,vars)};
228 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << std::endl ;
229
230 // Check if PDF supports maximum finding for the entire phase space
231 RooArgSet allVars(_otherVars);
232 allVars.add(_directVars);
234 if (maxFindCode != 0) {
235 // Special case: PDF supports max-finding in otherVars, no need to scan other+proto space for maximum
236 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
237 << "model supports analytical maximum finding in the full phase space: "
238 << "can provide analytical pdf maximum to numeric generator" << std::endl ;
239 double maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(&allVars);
240 _maxVar = std::make_unique<RooRealVar>("funcMax", "function maximum", maxVal);
241 } else {
243 if (maxFindCode != 0) {
245 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
246 << "model supports analytical maximum finding in the variables which are not"
247 << " internally generated. Can provide analytical pdf maximum to numeric "
248 << "generator" << std::endl;
249 cxcoutD(Generation) << "RooGenContext::ctor() maximum value must be reevaluated for each "
250 << "event with configuration code " << maxFindCode << std::endl ;
251 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",1) ;
252 }
253 }
254
255 if (!_maxVar) {
256
257 // Regular case: First find maximum in other+proto space
259
261
262 if (!_otherVars.empty()) {
263
264 cxcoutD(Generation) << "RooGenContext::ctor() prototype data provided, observables are generated numerically no "
265 << "analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
266 << "of accept/reject observables plus prototype observables: " << otherAndProto << std::endl ;
267
268 // Calculate maximum in other+proto space if there are any accept/reject generated observables
269 std::unique_ptr<RooAbsNumGenerator> maxFinder{RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,otherAndProto,RooArgSet(_protoVars),
270 *model.getGeneratorConfig(),_verbose)};
271// RooAcceptReject maxFinder(*_acceptRejectFunc,otherAndProto,RooNumGenConfig::defaultConfig(),_verbose) ;
272 double max = maxFinder->getFuncMax() ;
273 _maxVar = std::make_unique<RooRealVar>("funcMax","function maximum",max) ;
274
275 if (max==0) {
276 oocoutE(nullptr, Generation) << "RooGenContext::ctor(" << model.GetName()
277 << ") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
278 << "but chosen numeric generator (" << maxFinder->generatorName() << ") does not support maximum finding" << std::endl ;
279 throw string("RooGenContext::ctor()") ;
280 }
281
282 cxcoutD(Generation) << "RooGenContext::ctor() maximum function value found through initial sampling is " << max << std::endl ;
283 }
284 }
285
286 }
287
289 _generator = std::unique_ptr<RooAbsNumGenerator>{RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,_otherVars,RooArgSet(_protoVars),*model.getGeneratorConfig(),_verbose,_maxVar.get())};
290 cxcoutD(Generation) << "RooGenContext::ctor() creating MC sampling generator " << _generator->generatorName() << " from function for observables " << _otherVars << std::endl ;
291 //_generator= new RooAcceptReject(*_acceptRejectFunc,_otherVars,RooNumGenConfig::defaultConfig(),_verbose,_maxVar);
292 } else {
293 _generator = nullptr ;
294 }
295
297}
298
300
301////////////////////////////////////////////////////////////////////////////////
302/// Attach the cloned model to the event buffer we will be filling.
303
305{
307 if (_acceptRejectFunc) {
308 _acceptRejectFunc->recursiveRedirectServers(args,false) ; // WVE DEBUG
309 }
310
311 // Attach the RooAcceptReject generator the event buffer
312 if (_generator) {
313 _generator->attachParameters(args) ;
314 }
315
316}
317
318
319
320////////////////////////////////////////////////////////////////////////////////
321/// Perform one-time initialization of the generator context
322
324{
325 for (auto* arg : theEvent) {
326 arg->setOperMode(RooAbsArg::ADirty) ;
327 }
328
330
331 // Reset the cloned model's error counters.
333
334 // Initialize the PDFs internal generator
335 if (!_directVars.empty()) {
336 cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << std::endl ;
338 }
339}
340
341
342////////////////////////////////////////////////////////////////////////////////
343/// Generate one event. The 'remaining' integer is not used other than
344/// for printing messages
345
347{
348 if(!_otherVars.empty()) {
349 // call the accept-reject generator to generate its variables
350
351 if (_updateFMaxPerEvent!=0) {
353 cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << std::endl ;
354 _maxVar->setVal(max) ;
355 }
356
357 if (_generator) {
358 double resampleRatio(1) ;
359 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
360 if (resampleRatio<1) {
361 coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
362 << resampleRatio << " due to increased maximum weight" << std::endl ;
364 }
365 if(nullptr == subEvent) {
366 coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << std::endl ;
367 return;
368 }
369 theEvent.assignValueOnly(*subEvent) ;
370 //theEvent= *subEvent;
371
372 }
373 }
374
375 // Use the model's optimized generator, if one is available.
376 // The generator writes directly into our local 'event' since we attached it above.
377 if(!_directVars.empty()) {
379 }
380
381 // Generate uniform variables (non-dependents)
382 for(auto * uniVar : _uniformVars) {
383 auto * arglv = dynamic_cast<RooAbsLValue*>(uniVar);
384 if (!arglv) {
385 coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << std::endl ;
387 }
388 arglv->randomize() ;
389 }
390 theEvent.assign(_uniformVars) ;
391
392}
393
394
395////////////////////////////////////////////////////////////////////////////////
396/// Printing interface
397
398void RooGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
399{
401 os << indent << " --- RooGenContext --- " << std::endl ;
402 os << indent << "Using PDF ";
404
405 if(verbose) {
406 os << indent << "Use PDF generator for " << _directVars << std::endl ;
407 os << indent << "Use MC sampling generator " << (_generator ? _generator->generatorName() : "<none>") << " for " << _otherVars << std::endl ;
408 if (!_protoVars.empty()) {
409 os << indent << "Prototype observables are " << _protoVars << std::endl ;
410 }
411 }
412}
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define oocoutE(o, a)
#define coutE(a)
#define ccxcoutI(a)
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h direct
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
bool recursiveRedirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool recurseInNewSet=true)
Recursively replace all servers with the new servers in newSet.
void setOperMode(OperMode mode, bool recurseADirty=true)
Set the operation mode of this node.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for generator contexts of RooAbsPdf objects.
RooArgSet _theEvent
Pointer to observable event being generated.
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.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:195
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), bool force=true)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
virtual double maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
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:34
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 ...
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
Basic string class.
Definition TString.h:139