47using std::endl, std::string, std::ostream;
63 bool verbose,
const RooArgSet* forceDirect) :
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 ;
79 coutE(Generation) <<
"RooGenContext::RooGenContext(" <<
GetName() <<
") Couldn't deep-clone PDF, abort," << std::endl ;
87 if (prototype&&
_pdfClone->dependsOn(*prototype->
get())) {
89 fullNormSet.
add(*prototype->
get(),
true) ;
90 _pdfClone->fixAddCoefNormalization(fullNormSet) ;
100 if(!tmp->isFundamental()) {
101 coutE(Generation) <<
"RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() <<
"\"" << std::endl;
108 coutI(Generation) <<
"RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
109 <<
"\" which will have uniform distribution" << std::endl;
117 if (forceDirect==
nullptr || !forceDirect->
find(direct->
GetName())) {
119 cxcoutD(Generation) <<
"RooGenContext::ctor() observable " << arg->
GetName() <<
" has been determined to be unsafe for internal generation" << std::endl;
135 coutE(Generation) <<
"RooGenContext::ctor() constructor failed with errors" << std::endl;
142 cxcoutD(Generation) <<
"RooGenContext::ctor() observables " <<
_directVars <<
" are safe for internal generator (if supported by p.d.f)" << std::endl ;
145 cxcoutD(Generation) <<
"RooGenContext::ctor() observables " <<
_otherVars <<
" are NOT safe for internal generator (if supported by p.d.f)" << std::endl ;
152 cxcoutD(Generation) <<
"RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << std::endl ;
160 cxcoutD(Generation) <<
"RooGenContext::ctor() Model indicates that it can internally generate observables "
161 << generatedVars <<
" with configuration identifier " <<
_code << std::endl ;
174 cxcoutI(Generation) <<
"RooGenContext::ctor() Context will" ;
189 ntitle.
Append(
" (Accept/Reject)") ;
200 if(depList.
empty()) {
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 ;
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 ;
217 cxcoutI(Generation) <<
"RooGenContext::ctor() accept/reject sampling function is " <<
_acceptRejectFunc->GetName() << std::endl ;
227 cxcoutI(Generation) <<
"RooGenContext::ctor() accept/reject sampling function is " <<
_acceptRejectFunc->GetName() << std::endl ;
233 if (maxFindCode != 0) {
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 ;
239 _maxVar = std::make_unique<RooRealVar>(
"funcMax",
"function maximum", maxVal);
242 if (maxFindCode != 0) {
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) ;
259 otherAndProto.
add(protoDeps) ;
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 ;
269 *model.getGeneratorConfig(),
_verbose)};
271 double max = maxFinder->getFuncMax() ;
272 _maxVar = std::make_unique<RooRealVar>(
"funcMax",
"function maximum",max) ;
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()") ;
281 cxcoutD(Generation) <<
"RooGenContext::ctor() maximum function value found through initial sampling is " << max << std::endl ;
289 cxcoutD(Generation) <<
"RooGenContext::ctor() creating MC sampling generator " <<
_generator->generatorName() <<
" from function for observables " <<
_otherVars << std::endl ;
305 _pdfClone->recursiveRedirectServers(args,
false);
331 cxcoutD(Generation) <<
"RooGenContext::initGenerator() initializing internal generator of model with code " <<
_code << std::endl ;
348 cxcoutD(Generation) <<
"RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << std::endl ;
353 double resampleRatio(1) ;
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 ;
360 if(
nullptr == subEvent) {
361 coutE(Generation) <<
"RooGenContext::generateEvent ERROR accept/reject generator failed" << std::endl ;
380 coutE(Generation) <<
"RooGenContext::generateEvent(" <<
GetName() <<
") ERROR: uniform variable " << uniVar->GetName() <<
" is not an lvalue" << std::endl ;
396 os <<
indent <<
" --- RooGenContext --- " << std::endl ;
397 os <<
indent <<
"Using PDF ";
int Int_t
Signed integer 4 bytes (int).
static void indent(ostringstream &buf, int indent_level)
Common abstract base class for objects that represent a value and a "shape" in RooFit.
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?
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.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Container class to hold unbinned data.
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.
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.
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 ¶ms) 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.
TString & Append(const char *cs)