Logo ROOT   6.18/05
Reference Guide
RooAbsPdf.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $Id: RooAbsPdf.h,v 1.90 2007/07/21 21:32:52 wouter Exp $
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#ifndef ROO_ABS_PDF
17#define ROO_ABS_PDF
18
19#include "RooAbsReal.h"
20//#include "RooRealIntegral.h"
21#include "RooNameSet.h"
22#include "RooObjCacheManager.h"
23#include "RooCmdArg.h"
24
25class RooDataSet;
26class RooDataHist ;
27class RooArgSet ;
28class RooRealProxy ;
29class RooAbsGenContext ;
30class RooFitResult ;
31class RooExtendPdf ;
32class RooCategory ;
33class TPaveText;
34class TH1F;
35class TH2F;
36class TList ;
37class RooLinkedList ;
38class RooNumGenConfig ;
39class RooRealIntegral ;
40
41class RooAbsPdf : public RooAbsReal {
42public:
43
44 // Constructors, assignment etc
45 RooAbsPdf() ;
46 RooAbsPdf(const char *name, const char *title=0) ;
47 RooAbsPdf(const char *name, const char *title, Double_t minVal, Double_t maxVal) ;
48 // RooAbsPdf(const RooAbsPdf& other, const char* name=0);
49 virtual ~RooAbsPdf();
50
51 // Toy MC generation
52
53 ////////////////////////////////////////////////////////////////////////////////
54 /// See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&)
55 /// \param[in] nEvents How many events to generate
56 RooDataSet *generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg& arg1,
57 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
58 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) {
59 return generate(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5) ;
60 }
61 RooDataSet *generate(const RooArgSet &whatVars,
62 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
63 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
64 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
65 RooDataSet *generate(const RooArgSet &whatVars, Double_t nEvents = 0, Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE,
66 const char* binnedTag="", Bool_t expectedData=kFALSE, Bool_t extended = kFALSE) const;
67 RooDataSet *generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents= 0,
68 Bool_t verbose=kFALSE, Bool_t randProtoOrder=kFALSE, Bool_t resampleProto=kFALSE) const;
69
70
71 class GenSpec {
72 public:
73 virtual ~GenSpec() ;
75 private:
76 GenSpec(RooAbsGenContext* context, const RooArgSet& whatVars, RooDataSet* protoData, Int_t nGen, Bool_t extended,
77 Bool_t randProto, Bool_t resampleProto, TString dsetName, Bool_t init=kFALSE) ;
78 GenSpec(const GenSpec& other) ;
79
80 friend class RooAbsPdf ;
90 ClassDef(GenSpec,0) // Generation specification
91 } ;
92
93 ///Prepare GenSpec configuration object for efficient generation of multiple datasets from identical specification.
94 GenSpec* prepareMultiGen(const RooArgSet &whatVars,
95 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
96 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
97 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
98 ///Generate according to GenSpec obtained from prepareMultiGen().
99 RooDataSet* generate(GenSpec&) const ;
100
101
102 ////////////////////////////////////////////////////////////////////////////////
103 /// As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&)
104 /// \param[in] nEvents How many events to generate
105 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg& arg1,
106 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
107 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) {
108 return generateBinned(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5);
109 }
110 virtual RooDataHist *generateBinned(const RooArgSet &whatVars,
111 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
112 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
113 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
114 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, Bool_t expectedData=kFALSE, Bool_t extended=kFALSE) const;
115
116 virtual RooDataSet* generateSimGlobal(const RooArgSet& whatVars, Int_t nEvents) ;
117
118 ///Helper calling plotOn(RooPlot*, RooLinkedList&) const
119 virtual RooPlot* plotOn(RooPlot* frame,
120 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
121 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
122 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
123 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none(),
124 const RooCmdArg& arg9=RooCmdArg::none(), const RooCmdArg& arg10=RooCmdArg::none()
125 ) const {
126 return RooAbsReal::plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9,arg10) ;
127 }
128 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
129
130 /// Add a box with parameter values (and errors) to the specified frame
131 virtual RooPlot* paramOn(RooPlot* frame,
132 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
133 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
134 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
135 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
136
137 virtual RooPlot* paramOn(RooPlot* frame, const RooAbsData* data, const char *label= "", Int_t sigDigits = 2,
138 Option_t *options = "NELU", Double_t xmin=0.50,
139 Double_t xmax= 0.99,Double_t ymax=0.95) ;
140
141 // Built-in generator support
142 virtual Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
143 virtual void initGenerator(Int_t code) ;
144 virtual void generateEvent(Int_t code);
145 virtual Bool_t isDirectGenSafe(const RooAbsArg& arg) const ;
146
147 // Configuration of MC generators used for this pdf
148 const RooNumGenConfig* getGeneratorConfig() const ;
152 void setGeneratorConfig() ;
153 void setGeneratorConfig(const RooNumGenConfig& config) ;
154
155 // -log(L) fits to binned and unbinned data
156 virtual RooFitResult* fitTo(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
157 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
158 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
159 virtual RooFitResult* fitTo(RooAbsData& data, const RooLinkedList& cmdList) ;
160
161 virtual RooAbsReal* createNLL(RooAbsData& data, const RooLinkedList& cmdList) ;
162 virtual RooAbsReal* createNLL(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
163 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
164 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
165
166 // Chi^2 fits to histograms
169 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
170 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
171 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
172 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
173
174 // Chi^2 fits to X-Y datasets
175 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
176
177
178
179
180
181 // Constraint management
182 virtual RooArgSet* getConstraints(const RooArgSet& /*observables*/, RooArgSet& /*constrainedParams*/, Bool_t /*stripDisconnected*/) const {
183 // Interface to retrieve constraint terms on this pdf. Default implementation returns null
184 return 0 ;
185 }
186 virtual RooArgSet* getAllConstraints(const RooArgSet& observables, RooArgSet& constrainedParams, Bool_t stripDisconnected=kTRUE) const ;
187
188 // Project p.d.f into lower dimensional p.d.f
189 virtual RooAbsPdf* createProjection(const RooArgSet& iset) ;
190
191 // Create cumulative density function from p.d.f
192 RooAbsReal* createCdf(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
193 RooAbsReal* createCdf(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
194 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
195 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
196 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
197 RooAbsReal* createScanCdf(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
198
199 // Function evaluation support
200 virtual Bool_t traceEvalHook(Double_t value) const ;
201 virtual Double_t getValV(const RooArgSet* set=0) const ;
202 virtual Double_t getLogVal(const RooArgSet* set=0) const ;
203
204 Double_t getNorm(const RooArgSet& nset) const {
205 // Get p.d.f normalization term needed for observables 'nset'
206 return getNorm(&nset) ;
207 }
208 virtual Double_t getNorm(const RooArgSet* set=0) const ;
209
210 virtual void resetErrorCounters(Int_t resetValue=10) ;
211 void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE) ;
212 Bool_t traceEvalPdf(Double_t value) const ;
213
214 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
215
216 virtual Bool_t selfNormalized() const {
217 // If true, p.d.f is taken as self-normalized and no attempt is made to add a normalization term
218 // This default implementation return false
219 return kFALSE ;
220 }
221
222 // Support for extended maximum likelihood, switched off by default
224 virtual ExtendMode extendMode() const {
225 // Returns ability of p.d.f to provided extended likelihood terms. Possible
226 // answers are CanNotBeExtended, CanBeExtended or MustBeExtended. This
227 // default implementation always return CanNotBeExtended
228 return CanNotBeExtended ;
229 }
230 inline Bool_t canBeExtended() const {
231 // If true p.d.f can provide extended likelihood term
232 return (extendMode() != CanNotBeExtended) ;
233 }
234 inline Bool_t mustBeExtended() const {
235 // If true p.d.f must extended likelihood term
236 return (extendMode() == MustBeExtended) ;
237 }
238 virtual Double_t expectedEvents(const RooArgSet* nset) const ;
239 virtual Double_t expectedEvents(const RooArgSet& nset) const {
240 // Return expecteded number of p.d.fs to be used in calculated of extended likelihood
241 return expectedEvents(&nset) ;
242 }
243
244 // Printing interface (human readable)
245 virtual void printValue(std::ostream& os) const ;
246 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
247
248 static void verboseEval(Int_t stat) ;
249 static int verboseEval() ;
250
251 virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet* nset=0) const ;
252
253 static void clearEvalError() ;
254 static Bool_t evalError() ;
255
256 void setNormRange(const char* rangeName) ;
257 const char* normRange() const {
258 return _normRange.Length()>0 ? _normRange.Data() : 0 ;
259 }
260 void setNormRangeOverride(const char* rangeName) ;
261
262 const RooAbsReal* getNormIntegral(const RooArgSet& nset) const { return getNormObj(0,&nset,0) ; }
263
264protected:
265
266public:
267 virtual const RooAbsReal* getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName=0) const ;
268
269 virtual RooAbsGenContext* binnedGenContext(const RooArgSet &vars, Bool_t verbose= kFALSE) const ;
270
271 virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=0,
272 const RooArgSet* auxProto=0, Bool_t verbose= kFALSE) const ;
273
274 virtual RooAbsGenContext* autoGenContext(const RooArgSet &vars, const RooDataSet* prototype=0, const RooArgSet* auxProto=0,
275 Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE, const char* binnedTag="") const ;
276
277protected:
278
279 RooDataSet *generate(RooAbsGenContext& context, const RooArgSet& whatVars, const RooDataSet* prototype,
280 Double_t nEvents, Bool_t verbose, Bool_t randProtoOrder, Bool_t resampleProto, Bool_t skipInit=kFALSE,
281 Bool_t extended=kFALSE) const ;
282
283 // Implementation version
284 virtual RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants=kFALSE,
285 const char *label= "", Int_t sigDigits = 2, Option_t *options = "NELU", Double_t xmin=0.65,
286 Double_t xmax= 0.99,Double_t ymax=0.95, const RooCmdArg* formatCmd=0) ;
287
288
289 virtual RooPlot *plotOn(RooPlot *frame, PlotOpt o) const;
290
291 friend class RooEffGenContext ;
292 friend class RooAddGenContext ;
293 friend class RooProdGenContext ;
294 friend class RooSimGenContext ;
296 friend class RooConvGenContext ;
297 friend class RooSimultaneous ;
298 friend class RooAddGenContextOrig ;
299 friend class RooProdPdf ;
300 friend class RooMCStudy ;
301
302 Int_t* randomizeProtoOrder(Int_t nProto,Int_t nGen,Bool_t resample=kFALSE) const ;
303
304 friend class RooExtendPdf ;
305 // This also forces the definition of a copy ctor in derived classes
306 RooAbsPdf(const RooAbsPdf& other, const char* name = 0);
307
308 friend class RooRealIntegral ;
310
311 virtual Bool_t syncNormalization(const RooArgSet* dset, Bool_t adjustProxies=kTRUE) const ;
312
313 friend class RooAbsAnaConvPdf ;
315 mutable RooAbsReal* _norm ; //! Normalization integral (owned by _normMgr)
316 mutable RooArgSet* _normSet ; //! Normalization set with for above integral
317
319 public:
320 CacheElem(RooAbsReal& norm) : _norm(&norm) {} ;
322 virtual ~CacheElem() ;
325 } ;
326 mutable RooObjCacheManager _normMgr ; // The cache manager
327
328 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
329
331 // Hook function intercepting redirectServer calls. Discard current normalization
332 // object if any server is redirected
333
334 // Object is own by _normCacheManager that will delete object as soon as cache
335 // is sterilized by server redirect
336 _norm = 0 ;
337 return kFALSE ;
338 } ;
339
340
341 mutable Int_t _errorCount ; // Number of errors remaining to print
342 mutable Int_t _traceCount ; // Number of traces remaining to print
343 mutable Int_t _negCount ; // Number of negative probablities remaining to print
344
345 Bool_t _selectComp ; // Component selection flag for RooAbsPdf::plotCompOn
346
347 static void raiseEvalError() ;
348
350
351 RooNumGenConfig* _specGeneratorConfig ; //! MC generator configuration specific for this object
352
353 TString _normRange ; // Normalization range
355
356 ClassDef(RooAbsPdf,4) // Abstract PDF with normalization support
357};
358
359
360#endif
static Int_t init()
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassDef(name, id)
Definition: Rtypes.h:326
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float ymax
Definition: THbookFile.cxx:93
RooAbsAnaConvPdf is the base class for PDFs that represent a physics model that can be analytically c...
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:70
friend class RooArgSet
Definition: RooAbsArg.h:516
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Normalization set with for above integral.
Definition: RooAbsPdf.h:318
virtual ~CacheElem()
Destructor of normalization cache element.
Definition: RooAbsPdf.cxx:2961
virtual RooArgList containedArgs(Action)
Definition: RooAbsPdf.h:323
CacheElem(RooAbsReal &norm)
Definition: RooAbsPdf.h:320
void operModeHook(RooAbsArg::OperMode)
Dummy implementation.
Definition: RooAbsPdf.cxx:2950
RooAbsReal * _norm
Definition: RooAbsPdf.h:324
RooArgSet _whatVars
Definition: RooAbsPdf.h:82
RooAbsGenContext * _genContext
Definition: RooAbsPdf.h:81
Bool_t _resampleProto
Definition: RooAbsPdf.h:87
GenSpec(const GenSpec &other)
TString _dsetName
Definition: RooAbsPdf.h:88
RooDataSet * _protoData
Definition: RooAbsPdf.h:83
GenSpec * prepareMultiGen(const RooArgSet &whatVars, const RooCmdArg &arg1=RooCmdArg::none(), 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())
Prepare GenSpec configuration object for efficient generation of multiple datasets from identical spe...
Definition: RooAbsPdf.cxx:1871
RooObjCacheManager _normMgr
Definition: RooAbsPdf.h:326
Double_t getNorm(const RooArgSet &nset) const
Definition: RooAbsPdf.h:204
static Bool_t evalError()
Return the evaluation error flag.
Definition: RooAbsPdf.cxx:3149
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition: RooAbsPdf.cxx:2129
virtual ~RooAbsPdf()
Destructor.
Definition: RooAbsPdf.cxx:249
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
Definition: RooAbsPdf.cxx:3234
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...
Definition: RooAbsPdf.cxx:583
virtual RooAbsReal * createChi2(RooDataHist &data, const RooCmdArg &arg1=RooCmdArg::none(), 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())
Create a from a histogram and this function.
Definition: RooAbsPdf.cxx:1558
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooLinkedList &cmdList)
Calls RooAbsPdf::createChi2(RooDataSet& data, const RooLinkedList& cmdList) and returns fit result.
Definition: RooAbsPdf.cxx:1513
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:2940
Bool_t traceEvalPdf(Double_t value) const
Check that passed value is positive and not 'not-a-number'.
Definition: RooAbsPdf.cxx:356
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:800
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Interface function to create a generator context from a p.d.f.
Definition: RooAbsPdf.cxx:1709
RooAbsReal * createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
Definition: RooAbsPdf.cxx:3095
TString _normRange
MC generator configuration specific for this object.
Definition: RooAbsPdf.h:353
void setNormRange(const char *rangeName)
Definition: RooAbsPdf.cxx:3265
const RooAbsReal * getNormIntegral(const RooArgSet &nset) const
Definition: RooAbsPdf.h:262
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:2107
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Definition: RooAbsPdf.cxx:3179
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
virtual Bool_t selfNormalized() const
Definition: RooAbsPdf.h:216
virtual Bool_t traceEvalHook(Double_t value) const
WVE 08/21/01 Probably obsolete now.
Definition: RooAbsPdf.cxx:549
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print multi line detailed information of this RooAbsPdf.
Definition: RooAbsPdf.cxx:1682
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
Definition: RooAbsPdf.cxx:336
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition: RooAbsPdf.h:56
virtual RooAbsGenContext * autoGenContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE, const char *binnedTag="") const
Definition: RooAbsPdf.cxx:1718
Int_t _traceCount
Definition: RooAbsPdf.h:342
virtual Double_t expectedEvents(const RooArgSet &nset) const
Definition: RooAbsPdf.h:239
RooAbsReal * _norm
Definition: RooAbsPdf.h:315
virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet *nset=0) const
Returned the extended likelihood term (Nexpect - Nobserved*log(NExpected) of this PDF for the given n...
Definition: RooAbsPdf.cxx:661
virtual void printValue(std::ostream &os) const
Print value of p.d.f, also print normalization integral that was last used, if any.
Definition: RooAbsPdf.cxx:1666
virtual RooArgSet * getAllConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, Bool_t stripDisconnected=kTRUE) const
This helper function finds and collects all constraints terms of all coponent p.d....
Definition: RooAbsPdf.cxx:3112
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), 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())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1119
static void raiseEvalError()
Raise the evaluation error flag.
Definition: RooAbsPdf.cxx:3159
virtual Bool_t syncNormalization(const RooArgSet *dset, Bool_t adjustProxies=kTRUE) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
Definition: RooAbsPdf.cxx:450
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:267
Int_t _errorCount
Definition: RooAbsPdf.h:338
Bool_t _selectComp
Definition: RooAbsPdf.h:345
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet &, Bool_t) const
Definition: RooAbsPdf.h:182
@ CanBeExtended
Definition: RooAbsPdf.h:223
@ MustBeExtended
Definition: RooAbsPdf.h:223
@ CanNotBeExtended
Definition: RooAbsPdf.h:223
virtual RooDataHist * generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition: RooAbsPdf.h:105
virtual RooPlot * paramOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), 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())
Add a box with parameter values (and errors) to the specified frame.
Definition: RooAbsPdf.cxx:2769
void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE)
Reset trace counter to given value, limiting the number of future trace messages for this pdf to 'val...
Definition: RooAbsPdf.cxx:595
const char * normRange() const
Definition: RooAbsPdf.h:257
RooAbsReal * createCdf(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Create a cumulative distribution function of this p.d.f in terms of the observables listed in iset.
Definition: RooAbsPdf.cxx:3017
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=0) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
Definition: RooAbsPdf.cxx:415
virtual Double_t getLogVal(const RooArgSet *set=0) const
Return the log of the current value with given normalization An error message is printed if the argum...
Definition: RooAbsPdf.cxx:621
Int_t * randomizeProtoOrder(Int_t nProto, Int_t nGen, Bool_t resample=kFALSE) const
Return lookup table with randomized access order for prototype events, given nProto prototype data ev...
Definition: RooAbsPdf.cxx:2063
static Bool_t _evalError
Definition: RooAbsPdf.h:349
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Definition: RooAbsPdf.cxx:2920
Int_t _negCount
Definition: RooAbsPdf.h:343
Bool_t mustBeExtended() const
Definition: RooAbsPdf.h:234
void setNormRangeOverride(const char *rangeName)
Definition: RooAbsPdf.cxx:3282
virtual RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
Definition: RooAbsPdf.cxx:2368
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:2142
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
Definition: RooAbsPdf.cxx:3207
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2117
virtual ExtendMode extendMode() const
Definition: RooAbsPdf.h:224
RooAbsPdf()
Default constructor.
Definition: RooAbsPdf.cxx:191
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, Bool_t verbose=kFALSE) const
Return a binned generator context.
Definition: RooAbsPdf.cxx:1699
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
Definition: RooAbsPdf.cxx:3169
RooNumGenConfig * _specGeneratorConfig
Definition: RooAbsPdf.h:351
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:316
friend class RooAddGenContextOrig
Definition: RooAbsPdf.h:298
Double_t _rawValue
Definition: RooAbsPdf.h:314
static TString _normRangeOverride
Definition: RooAbsPdf.h:354
static Int_t _verboseEval
Definition: RooAbsPdf.h:309
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Definition: RooAbsPdf.h:330
static void clearEvalError()
Clear the evaluation error flag.
Definition: RooAbsPdf.cxx:3139
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:2979
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), 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 RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition: RooAbsPdf.h:119
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
virtual RooAbsReal * createChi2(RooDataHist &data, const RooLinkedList &cmdList)
Create a variable from a histogram and this function.
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooCmdArg &arg1=RooCmdArg::none(), 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())
Perform a fit to given histogram.
RooAddGenContext is an efficient implementation of the generator context specific for RooAddPdf PDFs.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:27
static const RooCmdArg & none()
Return reference to null argument.
Definition: RooCmdArg.cxx:50
RooConvGenContext is an efficient implementation of the generator context specific for RooAbsAnaConvP...
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
RooEffGenContext is a specialized generator context for p.d.fs represented by class RooEffProd,...
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Definition: RooExtendPdf.h:22
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:36
RooMCStudy is a helper class to facilitate Monte Carlo studies such as 'goodness-of-fit' studies,...
Definition: RooMCStudy.h:32
RooNumGenConfig holds the configuration parameters of the various numeric integrators used by RooReal...
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
RooProdGenContext is an efficient implementation of the generator context specific for RooProdPdf PDF...
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
RooSimGenContext is an efficient implementation of the generator context specific for RooSimultaneous...
RooSimSplitGenContext is an efficient implementation of the generator context specific for RooSimulta...
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:248
A doubly linked list.
Definition: TList.h:44
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
A Pave (see TPave) with text, lines or/and boxes inside.
Definition: TPaveText.h:21
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
const char * Data() const
Definition: TString.h:364
RooCmdArg NumEvents(Int_t numEvents)