Logo ROOT   6.08/07
Reference Guide
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 
22 Class RooGenContext implement a universal generator context for all
23 RooAbsPdf classes that do not have or need a specialized generator
24 context. This generator context queries the input p.d.f which observables
25 it can generate internally and delegates generation of those observables
26 to the p.d.f if it deems that safe. The other observables are generated
27 use a RooAcceptReject sampling technique.
28 **/
29 
30 
31 #include "RooFit.h"
32 #include "RooMsgService.h"
33 #include "Riostream.h"
34 
35 #include "RooGenContext.h"
36 #include "RooGenContext.h"
37 #include "RooAbsPdf.h"
38 #include "RooDataSet.h"
39 #include "RooRealIntegral.h"
40 #include "RooAcceptReject.h"
41 #include "RooRealVar.h"
42 #include "RooDataHist.h"
43 #include "RooErrorHandler.h"
44 #include "RooNumGenConfig.h"
45 #include "RooNumGenFactory.h"
46 
47 #include "TString.h"
48 #include "TIterator.h"
49 #include "TClass.h"
50 
51 
52 
53 using namespace std;
54 
56  ;
57 
58 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Initialize a new context for generating events with the specified
62 /// variables, using the specified PDF model. A prototype dataset (if provided)
63 /// is not cloned and still belongs to the caller. The contents and shape
64 /// of this dataset can be changed between calls to generate() as long as the
65 /// expected columns to be copied to the generated dataset are present.
66 /// Any argument supplied in the forceDirect RooArgSet are always offered
67 /// for internal generation to the p.d.f., even if this is deemed unsafe by
68 /// the logic of RooGenContext.
69 
71  const RooDataSet *prototype, const RooArgSet* auxProto,
72  Bool_t verbose, const RooArgSet* forceDirect) :
73  RooAbsGenContext(model,vars,prototype,auxProto,verbose),
74  _cloneSet(0), _pdfClone(0), _acceptRejectFunc(0), _generator(0),
75  _maxVar(0), _uniIter(0), _updateFMaxPerEvent(0)
76 {
77  cxcoutI(Generation) << "RooGenContext::ctor() setting up event generator context for p.d.f. " << model.GetName()
78  << " for generation of observable(s) " << vars ;
79  if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
80  if (auxProto && auxProto->getSize()>0) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
81  if (forceDirect && forceDirect->getSize()>0) ccxcoutI(Generation) << " with internal generation forced for observables " << *forceDirect ;
82  ccxcoutI(Generation) << endl ;
83 
84 
85  // Clone the model and all nodes that it depends on so that this context
86  // is independent of any existing objects.
87  RooArgSet nodes(model,model.GetName());
88  _cloneSet= (RooArgSet*) nodes.snapshot(kTRUE);
89  if (!_cloneSet) {
90  coutE(Generation) << "RooGenContext::RooGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
92  }
93 
94  // Find the clone in the snapshot list
95  _pdfClone = (RooAbsPdf*)_cloneSet->find(model.GetName());
97 
98  // Optionally fix RooAddPdf normalizations
99  if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
100  RooArgSet fullNormSet(vars) ;
101  fullNormSet.add(*prototype->get()) ;
102  _pdfClone->fixAddCoefNormalization(fullNormSet) ;
103  }
104 
105  // Analyze the list of variables to generate...
106  _isValid= kTRUE;
107  TIterator *iterator= vars.createIterator();
108  TIterator *servers= _pdfClone->serverIterator();
109  const RooAbsArg *tmp = 0;
110  const RooAbsArg *arg = 0;
111  while((_isValid && (tmp= (const RooAbsArg*)iterator->Next()))) {
112  // is this argument derived?
113  if(!tmp->isFundamental()) {
114  coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << endl;
115  _isValid= kFALSE;
116  continue;
117  }
118  // lookup this argument in the cloned set of PDF dependents
119  arg= (const RooAbsArg*)_cloneSet->find(tmp->GetName());
120  if(0 == arg) {
121  coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
122  << "\" which will have uniform distribution" << endl;
123  _uniformVars.add(*tmp);
124  }
125  else {
126 
127  // does the model depend on this variable directly, ie, like "x" in
128  // f(x) or f(x,g(x,y)) or even f(x,x) ?
129  const RooAbsArg *direct= arg ;
130  if (forceDirect==0 || !forceDirect->find(direct->GetName())) {
131  if (!_pdfClone->isDirectGenSafe(*arg)) {
132  cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << endl;
133  direct=0 ;
134  }
135  }
136 
137  // does the model depend indirectly on this variable through an lvalue chain?
138  // otherwise, this variable will have to be generated with accept/reject
139  if(direct) {
140  _directVars.add(*arg);
141  } else {
142  _otherVars.add(*arg);
143  }
144  }
145  }
146  delete servers;
147  delete iterator;
148  if(!isValid()) {
149  coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << endl;
150  return;
151  }
152 
153  // At this point directVars are all variables that are safe to be generated directly
154  // otherVars are all variables that are _not_ safe to be generated directly
155  if (_directVars.getSize()>0) {
156  cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << endl ;
157  }
158  if (_otherVars.getSize()>0) {
159  cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << endl ;
160  }
161 
162  // If PDF depends on prototype data, direct generator cannot use static initialization
163  // in initGenerator()
164  Bool_t staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
165  if (!staticInitOK) {
166  cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << endl ;
167  }
168 
169  // Can the model generate any of the direct variables itself?
170  RooArgSet generatedVars;
171  if (_directVars.getSize()>0) {
172  _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
173 
174  cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
175  << generatedVars << " with configuration identifier " << _code << endl ;
176  } else {
177  _code = 0 ;
178  }
179 
180  // Move variables which cannot be generated into the list to be generated with accept/reject
181  _directVars.remove(generatedVars);
183 
184  // Update _directVars to only include variables that will actually be directly generated
186  _directVars.add(generatedVars);
187 
188  cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
189  if (_directVars.getSize()>0) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
190  if (_directVars.getSize()>0 && _otherVars.getSize()>0) ccxcoutI(Generation) << " and" ;
191  if (_otherVars.getSize()>0) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
192  if (_uniformVars.getSize()>0) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
193  ccxcoutI(Generation) << endl ;
194 
195  // initialize the accept-reject generator
197  depList->remove(_otherVars);
198 
199  TString nname(_pdfClone->GetName()) ;
200  nname.Append("_AccRej") ;
201  TString ntitle(_pdfClone->GetTitle()) ;
202  ntitle.Append(" (Accept/Reject)") ;
203 
204 
205  RooArgSet* protoDeps = model.getObservables(_protoVars) ;
206 
207 
208  if (_protoVars.getSize()==0) {
209 
210  // No prototype variables
211 
212  if(depList->getSize()==0) {
213  // All variable are generated with accept-reject
214 
215  // Check if PDF supports maximum finding
216  Int_t maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
217  if (maxFindCode != 0) {
218  coutI(Generation) << "RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
219  << "model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << endl ;
220  Double_t maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(_theEvent) ;
221  _maxVar = new RooRealVar("funcMax","function maximum",maxVal) ;
222  cxcoutD(Generation) << "RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << endl ;
223  }
224  }
225 
226  if (_otherVars.getSize()>0) {
227  _pdfClone->getVal(&vars) ; // WVE debug
228  _acceptRejectFunc= new RooRealIntegral(nname,ntitle,*_pdfClone,*depList,&vars);
229  cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
230  } else {
231  _acceptRejectFunc = 0 ;
232  }
233 
234  } else {
235 
236  // Generation _with_ prototype variable
237  depList->remove(_protoVars,kTRUE,kTRUE) ;
239  cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
240 
241  // Check if PDF supports maximum finding for the entire phase space
242  RooArgSet allVars(_otherVars);
243  allVars.add(_directVars);
244  Int_t maxFindCode = _pdfClone->getMaxVal(allVars) ;
245  if (maxFindCode != 0) {
246  // Special case: PDF supports max-finding in otherVars, no need to scan other+proto space for maximum
247  coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
248  << "model supports analytical maximum finding in the full phase space: "
249  << "can provide analytical pdf maximum to numeric generator" << endl ;
250  _maxVar = new RooRealVar("funcMax","function maximum",_pdfClone->maxVal(maxFindCode)) ;
251  } else {
252  maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
253  if (maxFindCode != 0) {
254  _updateFMaxPerEvent = maxFindCode ;
255  coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
256  << "model supports analytical maximum finding in the variables which are not"
257  << " internally generated. Can provide analytical pdf maximum to numeric "
258  << "generator" << endl;
259  cxcoutD(Generation) << "RooGenContext::ctor() maximum value must be reevaluated for each "
260  << "event with configuration code " << maxFindCode << endl ;
261  _maxVar = new RooRealVar("funcMax","function maximum",1) ;
262  }
263  }
264 
265  if (!_maxVar) {
266 
267  // Regular case: First find maximum in other+proto space
268  RooArgSet otherAndProto(_otherVars) ;
269 
270  otherAndProto.add(*protoDeps) ;
271 
272  if (_otherVars.getSize()>0) {
273 
274  cxcoutD(Generation) << "RooGenContext::ctor() prototype data provided, observables are generated numericaly no "
275  << "analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
276  << "of accept/reject observables plus prototype observables: " << otherAndProto << endl ;
277 
278  // Calculate maximum in other+proto space if there are any accept/reject generated observables
280  *model.getGeneratorConfig(),_verbose) ;
281 // RooAcceptReject maxFinder(*_acceptRejectFunc,otherAndProto,RooNumGenConfig::defaultConfig(),_verbose) ;
282  Double_t max = maxFinder->getFuncMax() ;
283  _maxVar = new RooRealVar("funcMax","function maximum",max) ;
284 
285  if (max==0) {
286  coutE(Generation) << "RooGenContext::ctor(" << model.GetName()
287  << ") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
288  << "but chosen numeric generator (" << maxFinder->IsA()->GetName() << ") does not support maximum finding" << endl ;
289  delete maxFinder ;
290  throw string("RooGenContext::ctor()") ;
291  }
292  delete maxFinder ;
293 
294  cxcoutD(Generation) << "RooGenContext::ctor() maximum function value found through initial sampling is " << max << endl ;
295  }
296  }
297 
298  }
299 
300  if (_acceptRejectFunc && _otherVars.getSize()>0) {
302  cxcoutD(Generation) << "RooGenContext::ctor() creating MC sampling generator " << _generator->IsA()->GetName() << " from function for observables " << _otherVars << endl ;
303  //_generator= new RooAcceptReject(*_acceptRejectFunc,_otherVars,RooNumGenConfig::defaultConfig(),_verbose,_maxVar);
304  } else {
305  _generator = 0 ;
306  }
307 
308  delete protoDeps ;
309  delete depList;
311 }
312 
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Destructor.
316 
318 {
319  // Clean up the cloned objects used in this context.
320  delete _cloneSet;
321 
322  // Clean up our accept/reject generator
323  if (_generator) delete _generator;
325  if (_maxVar) delete _maxVar ;
326  if (_uniIter) delete _uniIter ;
327 }
328 
329 
330 
331 ////////////////////////////////////////////////////////////////////////////////
332 /// Attach the cloned model to the event buffer we will be filling.
333 
335 {
337  if (_acceptRejectFunc) {
339  }
340 
341  // Attach the RooAcceptReject generator the event buffer
342  if (_generator) {
344  }
345 
346 }
347 
348 
349 
350 ////////////////////////////////////////////////////////////////////////////////
351 /// Perform one-time initialization of the generator context
352 
354 {
355  RooFIter iter = theEvent.fwdIterator() ;
356  RooAbsArg* arg ;
357  while((arg=iter.next())) {
359  }
360 
361  attach(theEvent) ;
362 
363  // Reset the cloned model's error counters.
365 
366  // Initialize the PDFs internal generator
367  if (_directVars.getSize()>0) {
368  cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << endl ;
370  }
371 
372  // Create iterator for uniform vars (if any)
373  if (_uniformVars.getSize()>0) {
375  }
376 }
377 
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 /// Generate one event. The 'remaining' integer is not used other than
381 /// for printing messages
382 
383 void RooGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
384 {
385  if(_otherVars.getSize() > 0) {
386  // call the accept-reject generator to generate its variables
387 
388  if (_updateFMaxPerEvent!=0) {
390  cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << endl ;
391  _maxVar->setVal(max) ;
392  }
393 
394  if (_generator) {
395  Double_t resampleRatio(1) ;
396  const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
397  if (resampleRatio<1) {
398  coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
399  << resampleRatio << " due to increased maximum weight" << endl ;
400  resampleData(resampleRatio) ;
401  }
402  if(0 == subEvent) {
403  coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << endl ;
404  return;
405  }
406  theEvent.assignValueOnly(*subEvent) ;
407  //theEvent= *subEvent;
408 
409  }
410  }
411 
412  // Use the model's optimized generator, if one is available.
413  // The generator writes directly into our local 'event' since we attached it above.
414  if(_directVars.getSize() > 0) {
416  }
417 
418  // Generate uniform variables (non-dependents)
419  if (_uniIter) {
420  _uniIter->Reset() ;
421  RooAbsArg* uniVar ;
422  while((uniVar=(RooAbsArg*)_uniIter->Next())) {
423  RooAbsLValue* arglv = dynamic_cast<RooAbsLValue*>(uniVar) ;
424  if (!arglv) {
425  coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << endl ;
427  }
428  arglv->randomize() ;
429  }
430  theEvent = _uniformVars ;
431  }
432 
433 }
434 
435 
436 ////////////////////////////////////////////////////////////////////////////////
437 /// Printing interface
438 
439 void RooGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
440 {
441  RooAbsGenContext::printMultiline(os,content,verbose,indent);
442  os << indent << " --- RooGenContext --- " << endl ;
443  os << indent << "Using PDF ";
445 
446  if(verbose) {
447  os << indent << "Use PDF generator for " << _directVars << endl ;
448  os << indent << "Use MC sampling generator " << (_generator ? _generator->IsA()->GetName() : "<none>") << " for " << _otherVars << endl ;
449  if (_protoVars.getSize()>0) {
450  os << indent << "Prototype observables are " << _protoVars << endl ;
451  }
452  }
453 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
RooAbsPdf * _pdfClone
Definition: RooGenContext.h:46
TIterator * createIterator(Bool_t dir=kIterForward) const
#define coutE(a)
Definition: RooMsgService.h:35
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, which is interpreted as an OR of &#39;enum ContentsOptions&#39; values and in the style given by &#39;enum StyleOption&#39;.
static void softAbort()
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:86
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
Definition: RooAbsArg.cxx:744
#define cxcoutI(a)
Definition: RooMsgService.h:84
virtual void Reset()=0
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
RooArgSet _uniformVars
Definition: RooGenContext.h:47
#define coutI(a)
Definition: RooMsgService.h:32
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
#define cxcoutD(a)
Definition: RooMsgService.h:80
Bool_t isValid() const
int Int_t
Definition: RtypesCore.h:41
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)
Generate one event.
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
RooAbsNumGenerator * createSampler(RooAbsReal &func, const RooArgSet &genVars, const RooArgSet &condVars, const RooNumGenConfig &config, Bool_t verbose=kFALSE, RooAbsReal *maxFuncVal=0)
Construct a numeric integrator instance that operates on function &#39;func&#39; and is configured with &#39;conf...
Iterator abstract base class.
Definition: TIterator.h:32
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects The class perfor...
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
void attachParameters(const RooArgSet &vars)
Reattach original parameters to function clone.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
RooArgSet * _cloneSet
Definition: RooGenContext.h:45
TIterator * _uniIter
Definition: RooGenContext.h:53
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual Bool_t isFundamental() const
Definition: RooAbsArg.h:157
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:205
Int_t getSize() const
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:503
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
virtual void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
Definition: RooAbsPdf.cxx:2123
virtual Double_t getFuncMax()
virtual const RooArgSet * generateEvent(UInt_t remaining, Double_t &resampleRatio)=0
bool verbose
Double_t getNorm(const RooArgSet &nset) const
Definition: RooAbsPdf.h:190
Class RooAbsNumGenerator is the abstract base class for MC event generator implementations like RooAc...
RooGenContext(const RooAbsPdf &model, const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE, const RooArgSet *forceDirect=0)
Initialize a new context for generating events with the specified variables, using the specified PDF ...
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to &#39;res...
Definition: RooAbsPdf.cxx:568
static void indent(ostringstream &buf, int indent_level)
RooAbsNumGenerator * _generator
Definition: RooGenContext.h:51
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
TIterator * serverIterator() const
Definition: RooAbsArg.h:112
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event &#39;index&#39;.
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
Definition: RooGenContext.h:30
RooAbsArg * next()
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Definition: RooAbsPdf.cxx:2136
virtual void randomize(const char *rangeName=0)=0
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
Definition: RooAbsPdf.cxx:3218
RooArgSet * _theEvent
#define ClassImp(name)
Definition: Rtypes.h:279
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, Bool_t oneSafe=kFALSE)
The assignment operator sets the value of any argument in our set that also appears in the other set...
double Double_t
Definition: RtypesCore.h:55
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition: RooAbsPdf.cxx:2101
RooFIter fwdIterator() const
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for multi-line printing.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
virtual ~RooGenContext()
Destructor.
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Printing interface.
Int_t _updateFMaxPerEvent
Definition: RooGenContext.h:54
void resampleData(Double_t &ratio)
Rescale existing output buffer with given ratio.
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
#define ccxcoutI(a)
Definition: RooMsgService.h:85
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2111
virtual TObject * Next()=0
static RooNumGenFactory & instance()
Static method returning reference to singleton instance of factory.
virtual void initGenerator(const RooArgSet &theEvent)
Perform one-time initialization of the generator context.
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Change cache operation mode to given mode.
Definition: RooAbsArg.cxx:1754
RooRealIntegral * _acceptRejectFunc
Definition: RooGenContext.h:50
RooArgSet _directVars
Definition: RooGenContext.h:47
RooArgSet _otherVars
Definition: RooGenContext.h:47
Abstract base class for objects that are lvalues, i.e.
Definition: RooAbsLValue.h:27
virtual void attach(const RooArgSet &params)
Attach the cloned model to the event buffer we will be filling.
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), Bool_t force=kTRUE)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Definition: RooAbsArg.cxx:1088
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
RooRealVar * _maxVar
Definition: RooGenContext.h:52
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52