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 _cloneSet(0), _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());
88 if (!_cloneSet) {
89 coutE(Generation) << "RooGenContext::RooGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
91 }
92
93 // Find the clone in the snapshot list
96
97 // Optionally fix RooAddPdf normalizations
98 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
99 RooArgSet fullNormSet(vars) ;
100 fullNormSet.add(*prototype->get()) ;
101 _pdfClone->fixAddCoefNormalization(fullNormSet) ;
102 }
103
104 // Analyze the list of variables to generate...
106 TIterator *iterator= vars.createIterator();
108 const RooAbsArg *tmp = 0;
109 const RooAbsArg *arg = 0;
110 while((_isValid && (tmp= (const RooAbsArg*)iterator->Next()))) {
111 // is this argument derived?
112 if(!tmp->isFundamental()) {
113 coutE(Generation) << "RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() << "\"" << endl;
115 continue;
116 }
117 // lookup this argument in the cloned set of PDF dependents
118 arg= (const RooAbsArg*)_cloneSet->find(tmp->GetName());
119 if(0 == arg) {
120 coutI(Generation) << "RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
121 << "\" which will have uniform distribution" << endl;
122 _uniformVars.add(*tmp);
123 }
124 else {
125
126 // does the model depend on this variable directly, ie, like "x" in
127 // f(x) or f(x,g(x,y)) or even f(x,x) ?
128 const RooAbsArg *direct= arg ;
129 if (forceDirect==0 || !forceDirect->find(direct->GetName())) {
130 if (!_pdfClone->isDirectGenSafe(*arg)) {
131 cxcoutD(Generation) << "RooGenContext::ctor() observable " << arg->GetName() << " has been determined to be unsafe for internal generation" << endl;
132 direct=0 ;
133 }
134 }
135
136 // does the model depend indirectly on this variable through an lvalue chain?
137 // otherwise, this variable will have to be generated with accept/reject
138 if(direct) {
139 _directVars.add(*arg);
140 } else {
141 _otherVars.add(*arg);
142 }
143 }
144 }
145 delete servers;
146 delete iterator;
147 if(!isValid()) {
148 coutE(Generation) << "RooGenContext::ctor() constructor failed with errors" << endl;
149 return;
150 }
151
152 // At this point directVars are all variables that are safe to be generated directly
153 // otherVars are all variables that are _not_ safe to be generated directly
154 if (_directVars.getSize()>0) {
155 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _directVars << " are safe for internal generator (if supported by p.d.f)" << endl ;
156 }
157 if (_otherVars.getSize()>0) {
158 cxcoutD(Generation) << "RooGenContext::ctor() observables " << _otherVars << " are NOT safe for internal generator (if supported by p.d.f)" << endl ;
159 }
160
161 // If PDF depends on prototype data, direct generator cannot use static initialization
162 // in initGenerator()
163 Bool_t staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
164 if (!staticInitOK) {
165 cxcoutD(Generation) << "RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << endl ;
166 }
167
168 // Can the model generate any of the direct variables itself?
169 RooArgSet generatedVars;
170 if (_directVars.getSize()>0) {
171 _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
172
173 cxcoutD(Generation) << "RooGenContext::ctor() Model indicates that it can internally generate observables "
174 << generatedVars << " with configuration identifier " << _code << endl ;
175 } else {
176 _code = 0 ;
177 }
178
179 // Move variables which cannot be generated into the list to be generated with accept/reject
180 _directVars.remove(generatedVars);
182
183 // Update _directVars to only include variables that will actually be directly generated
185 _directVars.add(generatedVars);
186
187 cxcoutI(Generation) << "RooGenContext::ctor() Context will" ;
188 if (_directVars.getSize()>0) ccxcoutI(Generation) << " let model internally generate observables " << _directVars ;
189 if (_directVars.getSize()>0 && _otherVars.getSize()>0) ccxcoutI(Generation) << " and" ;
190 if (_otherVars.getSize()>0) ccxcoutI(Generation) << " generate variables " << _otherVars << " with accept/reject sampling" ;
191 if (_uniformVars.getSize()>0) ccxcoutI(Generation) << ", non-observable variables " << _uniformVars << " will be generated with flat distribution" ;
192 ccxcoutI(Generation) << endl ;
193
194 // initialize the accept-reject generator
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 = model.getObservables(_protoVars) ;
205
206
207 if (_protoVars.getSize()==0) {
208
209 // No prototype variables
210
211 if(depList->getSize()==0) {
212 // All variable are generated with accept-reject
213
214 // Check if PDF supports maximum finding
215 Int_t maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
216 if (maxFindCode != 0) {
217 coutI(Generation) << "RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
218 << "model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << endl ;
219 Double_t maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(_theEvent) ;
220 _maxVar = new RooRealVar("funcMax","function maximum",maxVal) ;
221 cxcoutD(Generation) << "RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << endl ;
222 }
223 }
224
225 if (_otherVars.getSize()>0) {
226 _pdfClone->getVal(&vars) ; // WVE debug
227 _acceptRejectFunc= new RooRealIntegral(nname,ntitle,*_pdfClone,*depList,&vars);
228 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
229 } else {
231 }
232
233 } else {
234
235 // Generation _with_ prototype variable
236 depList->remove(_protoVars,kTRUE,kTRUE) ;
238 cxcoutI(Generation) << "RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
239
240 // Check if PDF supports maximum finding for the entire phase space
241 RooArgSet allVars(_otherVars);
242 allVars.add(_directVars);
243 Int_t maxFindCode = _pdfClone->getMaxVal(allVars) ;
244 if (maxFindCode != 0) {
245 // Special case: PDF supports max-finding in otherVars, no need to scan other+proto space for maximum
246 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
247 << "model supports analytical maximum finding in the full phase space: "
248 << "can provide analytical pdf maximum to numeric generator" << endl ;
249 _maxVar = new RooRealVar("funcMax","function maximum",_pdfClone->maxVal(maxFindCode)) ;
250 } else {
251 maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
252 if (maxFindCode != 0) {
253 _updateFMaxPerEvent = maxFindCode ;
254 coutI(Generation) << "RooGenContext::ctor() prototype data provided, and "
255 << "model supports analytical maximum finding in the variables which are not"
256 << " internally generated. Can provide analytical pdf maximum to numeric "
257 << "generator" << endl;
258 cxcoutD(Generation) << "RooGenContext::ctor() maximum value must be reevaluated for each "
259 << "event with configuration code " << maxFindCode << endl ;
260 _maxVar = new RooRealVar("funcMax","function maximum",1) ;
261 }
262 }
263
264 if (!_maxVar) {
265
266 // Regular case: First find maximum in other+proto space
267 RooArgSet otherAndProto(_otherVars) ;
268
269 otherAndProto.add(*protoDeps) ;
270
271 if (_otherVars.getSize()>0) {
272
273 cxcoutD(Generation) << "RooGenContext::ctor() prototype data provided, observables are generated numericaly no "
274 << "analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
275 << "of accept/reject observables plus prototype observables: " << otherAndProto << endl ;
276
277 // Calculate maximum in other+proto space if there are any accept/reject generated observables
279 *model.getGeneratorConfig(),_verbose) ;
280// RooAcceptReject maxFinder(*_acceptRejectFunc,otherAndProto,RooNumGenConfig::defaultConfig(),_verbose) ;
281 Double_t max = maxFinder->getFuncMax() ;
282 _maxVar = new RooRealVar("funcMax","function maximum",max) ;
283
284 if (max==0) {
285 coutE(Generation) << "RooGenContext::ctor(" << model.GetName()
286 << ") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
287 << "but chosen numeric generator (" << maxFinder->IsA()->GetName() << ") does not support maximum finding" << endl ;
288 delete maxFinder ;
289 throw string("RooGenContext::ctor()") ;
290 }
291 delete maxFinder ;
292
293 cxcoutD(Generation) << "RooGenContext::ctor() maximum function value found through initial sampling is " << max << endl ;
294 }
295 }
296
297 }
298
301 cxcoutD(Generation) << "RooGenContext::ctor() creating MC sampling generator " << _generator->IsA()->GetName() << " from function for observables " << _otherVars << endl ;
302 //_generator= new RooAcceptReject(*_acceptRejectFunc,_otherVars,RooNumGenConfig::defaultConfig(),_verbose,_maxVar);
303 } else {
304 _generator = 0 ;
305 }
306
307 delete protoDeps ;
308 delete depList;
310}
311
312
313////////////////////////////////////////////////////////////////////////////////
314/// Destructor.
315
317{
318 // Clean up the cloned objects used in this context.
319 delete _cloneSet;
320
321 // Clean up our accept/reject generator
322 if (_generator) delete _generator;
324 if (_maxVar) delete _maxVar ;
325 if (_uniIter) delete _uniIter ;
326}
327
328
329
330////////////////////////////////////////////////////////////////////////////////
331/// Attach the cloned model to the event buffer we will be filling.
332
334{
336 if (_acceptRejectFunc) {
338 }
339
340 // Attach the RooAcceptReject generator the event buffer
341 if (_generator) {
343 }
344
345}
346
347
348
349////////////////////////////////////////////////////////////////////////////////
350/// Perform one-time initialization of the generator context
351
353{
354 RooFIter iter = theEvent.fwdIterator() ;
355 RooAbsArg* arg ;
356 while((arg=iter.next())) {
358 }
359
360 attach(theEvent) ;
361
362 // Reset the cloned model's error counters.
364
365 // Initialize the PDFs internal generator
366 if (_directVars.getSize()>0) {
367 cxcoutD(Generation) << "RooGenContext::initGenerator() initializing internal generator of model with code " << _code << endl ;
369 }
370
371 // Create iterator for uniform vars (if any)
372 if (_uniformVars.getSize()>0) {
374 }
375}
376
377
378////////////////////////////////////////////////////////////////////////////////
379/// Generate one event. The 'remaining' integer is not used other than
380/// for printing messages
381
383{
384 if(_otherVars.getSize() > 0) {
385 // call the accept-reject generator to generate its variables
386
387 if (_updateFMaxPerEvent!=0) {
389 cxcoutD(Generation) << "RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << endl ;
390 _maxVar->setVal(max) ;
391 }
392
393 if (_generator) {
394 Double_t resampleRatio(1) ;
395 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
396 if (resampleRatio<1) {
397 coutI(Generation) << "RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
398 << resampleRatio << " due to increased maximum weight" << endl ;
399 resampleData(resampleRatio) ;
400 }
401 if(0 == subEvent) {
402 coutE(Generation) << "RooGenContext::generateEvent ERROR accept/reject generator failed" << endl ;
403 return;
404 }
405 theEvent.assignValueOnly(*subEvent) ;
406 //theEvent= *subEvent;
407
408 }
409 }
410
411 // Use the model's optimized generator, if one is available.
412 // The generator writes directly into our local 'event' since we attached it above.
413 if(_directVars.getSize() > 0) {
415 }
416
417 // Generate uniform variables (non-dependents)
418 if (_uniIter) {
419 _uniIter->Reset() ;
420 RooAbsArg* uniVar ;
421 while((uniVar=(RooAbsArg*)_uniIter->Next())) {
422 RooAbsLValue* arglv = dynamic_cast<RooAbsLValue*>(uniVar) ;
423 if (!arglv) {
424 coutE(Generation) << "RooGenContext::generateEvent(" << GetName() << ") ERROR: uniform variable " << uniVar->GetName() << " is not an lvalue" << endl ;
426 }
427 arglv->randomize() ;
428 }
429 theEvent = _uniformVars ;
430 }
431
432}
433
434
435////////////////////////////////////////////////////////////////////////////////
436/// Printing interface
437
438void RooGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
439{
440 RooAbsGenContext::printMultiline(os,content,verbose,indent);
441 os << indent << " --- RooGenContext --- " << endl ;
442 os << indent << "Using PDF ";
444
445 if(verbose) {
446 os << indent << "Use PDF generator for " << _directVars << endl ;
447 os << indent << "Use MC sampling generator " << (_generator ? _generator->IsA()->GetName() : "<none>") << " for " << _otherVars << endl ;
448 if (_protoVars.getSize()>0) {
449 os << indent << "Prototype observables are " << _protoVars << endl ;
450 }
451 }
452}
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define coutE(a)
#define ccxcoutI(a)
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#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:72
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:313
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:243
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
TIterator * serverIterator() const
Definition RooAbsArg.h:140
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_t oneSafe=kFALSE)
The assignment operator 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 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.
RooArgSet * _theEvent
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:211
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.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
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 ...
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:29
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:118
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
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.
RooArgSet * _cloneSet
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
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