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
22Class RooGenContext implement 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
31#include "RooFit.h"
32#include "RooMsgService.h"
33#include "Riostream.h"
34
35#include "RooGenContext.h"
36#include "RooAbsPdf.h"
37#include "RooDataSet.h"
38#include "RooRealIntegral.h"
39#include "RooAcceptReject.h"
40#include "RooRealVar.h"
41#include "RooDataHist.h"
42#include "RooErrorHandler.h"
43#include "RooNumGenConfig.h"
44#include "RooNumGenFactory.h"
45
46#include "TString.h"
47#include "TIterator.h"
48#include "TClass.h"
49
50
51
52using namespace std;
53
55 ;
56
57
58
59////////////////////////////////////////////////////////////////////////////////
60/// Initialize a new context for generating events with the specified
61/// variables, using the specified PDF model. A prototype dataset (if provided)
62/// is not cloned and still belongs to the caller. The contents and shape
63/// of this dataset can be changed between calls to generate() as long as the
64/// expected columns to be copied to the generated dataset are present.
65/// Any argument supplied in the forceDirect RooArgSet are always offered
66/// for internal generation to the p.d.f., even if this is deemed unsafe by
67/// the logic of RooGenContext.
68
70 const RooDataSet *prototype, const RooArgSet* auxProto,
71 Bool_t verbose, const RooArgSet* forceDirect) :
72 RooAbsGenContext(model,vars,prototype,auxProto,verbose),
73 _pdfClone(0), _acceptRejectFunc(0), _generator(0),
74 _maxVar(0), _uniIter(0), _updateFMaxPerEvent(0)
75{
76 cxcoutI(Generation) << "RooGenContext::ctor() setting up event generator context for p.d.f. " << model.GetName()
77 << " for generation of observable(s) " << vars ;
78 if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
79 if (auxProto && auxProto->getSize()>0) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
80 if (forceDirect && forceDirect->getSize()>0) ccxcoutI(Generation) << " with internal generation forced for observables " << *forceDirect ;
81 ccxcoutI(Generation) << endl ;
82
83
84 // Clone the model and all nodes that it depends on so that this context
85 // is independent of any existing objects.
86 RooArgSet nodes(model,model.GetName());
87 if (nodes.snapshot(_cloneSet, true)) {
88 coutE(Generation) << "RooGenContext::RooGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
90 }
91
92 // Find the clone in the snapshot list
93 _pdfClone = static_cast<RooAbsPdf*>(_cloneSet.find(model.GetName()));
95
96 // Optionally fix RooAddPdf normalizations
97 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
98 RooArgSet fullNormSet(vars) ;
99 fullNormSet.add(*prototype->get()) ;
100 _pdfClone->fixAddCoefNormalization(fullNormSet) ;
101 }
102
103 // Analyze the list of variables to generate...
105 TIterator *iterator= vars.createIterator();
107 const RooAbsArg *tmp = 0;
108 const RooAbsArg *arg = 0;
109 while((_isValid && (tmp= (const RooAbsArg*)iterator->Next()))) {
110 // is this argument derived?
111 if(!tmp->isFundamental()) {
112 coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << endl;
114 continue;
115 }
116 // lookup this argument in the cloned set of PDF dependents
117 arg= static_cast<RooAbsArg const*>(_cloneSet.find(tmp->GetName()));
118 if(0 == arg) {
119 coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
120 << "\" which will have uniform distribution" << endl;
121 _uniformVars.add(*tmp);
122 }
123 else {
124
125 // does the model depend on this variable directly, ie, like "x" in
126 // f(x) or f(x,g(x,y)) or even f(x,x) ?
127 const RooAbsArg *direct= arg ;
128 if (forceDirect==0 || !forceDirect->find(direct->GetName())) {
129 if (!_pdfClone->isDirectGenSafe(*arg)) {
130 cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << endl;
131 direct=0 ;
132 }
133 }
134
135 // does the model depend indirectly on this variable through an lvalue chain?
136 // otherwise, this variable will have to be generated with accept/reject
137 if(direct) {
138 _directVars.add(*arg);
139 } else {
140 _otherVars.add(*arg);
141 }
142 }
143 }
144 delete servers;
145 delete iterator;
146 if(!isValid()) {
147 coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << endl;
148 return;
149 }
150
151 // At this point directVars are all variables that are safe to be generated directly
152 // otherVars are all variables that are _not_ safe to be generated directly
153 if (_directVars.getSize()>0) {
154 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << endl ;
155 }
156 if (_otherVars.getSize()>0) {
157 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << endl ;
158 }
159
160 // If PDF depends on prototype data, direct generator cannot use static initialization
161 // in initGenerator()
162 Bool_t staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
163 if (!staticInitOK) {
164 cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << endl ;
165 }
166
167 // Can the model generate any of the direct variables itself?
168 RooArgSet generatedVars;
169 if (_directVars.getSize()>0) {
170 _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
171
172 cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
173 << generatedVars << " with configuration identifier " << _code << endl ;
174 } else {
175 _code = 0 ;
176 }
177
178 // Move variables which cannot be generated into the list to be generated with accept/reject
179 _directVars.remove(generatedVars);
181
182 // Update _directVars to only include variables that will actually be directly generated
184 _directVars.add(generatedVars);
185
186 cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
187 if (_directVars.getSize()>0) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
188 if (_directVars.getSize()>0 && _otherVars.getSize()>0) ccxcoutI(Generation) << " and" ;
189 if (_otherVars.getSize()>0) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
190 if (_uniformVars.getSize()>0) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
191 ccxcoutI(Generation) << endl ;
192
193 // initialize the accept-reject generator
194 RooArgSet depList;
196 depList.remove(_otherVars);
197
198 TString nname(_pdfClone->GetName()) ;
199 nname.Append("_AccRej") ;
200 TString ntitle(_pdfClone->GetTitle()) ;
201 ntitle.Append(" (Accept/Reject)") ;
202
203
204 RooArgSet protoDeps;
205 model.getObservables(&_protoVars, protoDeps);
206
207
208 if (_protoVars.getSize()==0) {
209
210 // No prototype variables
211
212 if(depList.empty()) {
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 {
232 }
233
234 } else {
235
236 // Generation _with_ prototype variable
237 depList.remove(_protoVars,true,true);
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
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
309}
310
311
312////////////////////////////////////////////////////////////////////////////////
313/// Destructor.
314
316{
317 // Clean up our accept/reject generator
318 if (_generator) delete _generator;
320 if (_maxVar) delete _maxVar ;
321 if (_uniIter) delete _uniIter ;
322}
323
324
325
326////////////////////////////////////////////////////////////////////////////////
327/// Attach the cloned model to the event buffer we will be filling.
328
330{
332 if (_acceptRejectFunc) {
334 }
335
336 // Attach the RooAcceptReject generator the event buffer
337 if (_generator) {
339 }
340
341}
342
343
344
345////////////////////////////////////////////////////////////////////////////////
346/// Perform one-time initialization of the generator context
347
349{
350 RooFIter iter = theEvent.fwdIterator() ;
351 RooAbsArg* arg ;
352 while((arg=iter.next())) {
354 }
355
356 attach(theEvent) ;
357
358 // Reset the cloned model's error counters.
360
361 // Initialize the PDFs internal generator
362 if (_directVars.getSize()>0) {
363 cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << endl ;
365 }
366
367 // Create iterator for uniform vars (if any)
368 if (_uniformVars.getSize()>0) {
369 delete _uniIter;
371 }
372}
373
374
375////////////////////////////////////////////////////////////////////////////////
376/// Generate one event. The 'remaining' integer is not used other than
377/// for printing messages
378
380{
381 if(_otherVars.getSize() > 0) {
382 // call the accept-reject generator to generate its variables
383
384 if (_updateFMaxPerEvent!=0) {
386 cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << endl ;
387 _maxVar->setVal(max) ;
388 }
389
390 if (_generator) {
391 Double_t resampleRatio(1) ;
392 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
393 if (resampleRatio<1) {
394 coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
395 << resampleRatio << " due to increased maximum weight" << endl ;
396 resampleData(resampleRatio) ;
397 }
398 if(0 == subEvent) {
399 coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << endl ;
400 return;
401 }
402 theEvent.assignValueOnly(*subEvent) ;
403 //theEvent= *subEvent;
404
405 }
406 }
407
408 // Use the model's optimized generator, if one is available.
409 // The generator writes directly into our local 'event' since we attached it above.
410 if(_directVars.getSize() > 0) {
412 }
413
414 // Generate uniform variables (non-dependents)
415 if (_uniIter) {
416 _uniIter->Reset() ;
417 RooAbsArg* uniVar ;
418 while((uniVar=(RooAbsArg*)_uniIter->Next())) {
419 RooAbsLValue* arglv = dynamic_cast<RooAbsLValue*>(uniVar) ;
420 if (!arglv) {
421 coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << endl ;
423 }
424 arglv->randomize() ;
425 }
426 theEvent.assign(_uniformVars) ;
427 }
428
429}
430
431
432////////////////////////////////////////////////////////////////////////////////
433/// Printing interface
434
435void RooGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
436{
437 RooAbsGenContext::printMultiline(os,content,verbose,indent);
438 os << indent << " --- RooGenContext --- " << endl ;
439 os << indent << "Using PDF ";
441
442 if(verbose) {
443 os << indent << "Use PDF generator for " << _directVars << endl ;
444 os << indent << "Use MC sampling generator " << (_generator ? _generator->IsA()->GetName() : "<none>") << " for " << _otherVars << endl ;
445 if (_protoVars.getSize()>0) {
446 os << indent << "Prototype observables are " << _protoVars << endl ;
447 }
448 }
449}
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define coutE(a)
#define ccxcoutI(a)
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
static void indent(ostringstream &buf, int indent_level)
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
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.
virtual Bool_t isFundamental() const
Is this object a fundamental type that can be added to a dataset? Fundamental-type subclasses overrid...
Definition RooAbsArg.h:239
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
TIterator * serverIterator() const
Definition RooAbsArg.h:137
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Recursively replace all servers with the new servers in newSet.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, bool forceIfSizeOne=false)
Sets the value of any argument in our set that also appears in the other set.
Int_t getSize() const
RooFIter fwdIterator() const
One-time forward iterator.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
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...
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Bool_t isValid() const
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for multi-line printing.
void resampleData(Double_t &ratio)
Rescale existing output buffer with given ratio.
Abstract base class for objects that are lvalues, i.e.
virtual void randomize(const char *rangeName=0)=0
Class RooAbsNumGenerator is the abstract base class for MC event generator implementations like RooAc...
virtual Double_t getFuncMax()
void attachParameters(const RooArgSet &vars)
Reattach original parameters to function clone.
virtual const RooArgSet * generateEvent(UInt_t remaining, Double_t &resampleRatio)=0
Double_t getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:239
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 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...
virtual Bool_t 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.
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 ...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
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...
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:35
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:158
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:36
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
static void softAbort()
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
RooArgSet _otherVars
virtual void attach(const RooArgSet &params)
Attach the cloned model to the event buffer we will be filling.
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Printing interface.
RooArgSet _uniformVars
RooArgSet _directVars
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)
Generate one event.
TIterator * _uniIter
virtual ~RooGenContext()
Destructor.
RooAbsPdf * _pdfClone
RooArgSet _cloneSet
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 ...
Int_t _updateFMaxPerEvent
RooRealIntegral * _acceptRejectFunc
RooAbsNumGenerator * _generator
RooRealVar * _maxVar
virtual void initGenerator(const RooArgSet &theEvent)
Perform one-time initialization of the generator context.
static RooNumGenFactory & instance()
Static method returning reference to singleton instance of factory.
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 'func' and is configured with 'conf...
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,...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Iterator abstract base class.
Definition TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
virtual const char * GetTitle() const
Returns title of object.
Definition TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:136
TString & Append(const char *cs)
Definition TString.h:564