Logo ROOT  
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 RooAbsGenContext ;
29class RooFitResult ;
30class RooExtendPdf ;
31class RooCategory ;
32class TPaveText;
33class TH1F;
34class TH2F;
35class TList ;
36class RooLinkedList ;
37class RooNumGenConfig ;
38class RooRealIntegral ;
39
40class RooAbsPdf : public RooAbsReal {
41public:
42
43 // Constructors, assignment etc
44 RooAbsPdf() ;
45 RooAbsPdf(const char *name, const char *title=0) ;
46 RooAbsPdf(const char *name, const char *title, Double_t minVal, Double_t maxVal) ;
47 // RooAbsPdf(const RooAbsPdf& other, const char* name=0);
48 virtual ~RooAbsPdf();
49
50 // Toy MC generation
51
52 ////////////////////////////////////////////////////////////////////////////////
53 /// See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&)
54 /// \param[in] nEvents How many events to generate
55 RooDataSet *generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg& arg1,
56 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
57 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) {
58 return generate(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5) ;
59 }
60 RooDataSet *generate(const RooArgSet &whatVars,
61 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
62 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
63 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
64 RooDataSet *generate(const RooArgSet &whatVars, Double_t nEvents = 0, Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE,
65 const char* binnedTag="", Bool_t expectedData=kFALSE, Bool_t extended = kFALSE) const;
66 RooDataSet *generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents= 0,
67 Bool_t verbose=kFALSE, Bool_t randProtoOrder=kFALSE, Bool_t resampleProto=kFALSE) const;
68
69
70 class GenSpec {
71 public:
72 virtual ~GenSpec() ;
74 private:
75 GenSpec(RooAbsGenContext* context, const RooArgSet& whatVars, RooDataSet* protoData, Int_t nGen, Bool_t extended,
76 Bool_t randProto, Bool_t resampleProto, TString dsetName, Bool_t init=kFALSE) ;
77 GenSpec(const GenSpec& other) ;
78
79 friend class RooAbsPdf ;
89 ClassDef(GenSpec,0) // Generation specification
90 } ;
91
92 ///Prepare GenSpec configuration object for efficient generation of multiple datasets from identical specification.
93 GenSpec* prepareMultiGen(const RooArgSet &whatVars,
94 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
95 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
96 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
97 ///Generate according to GenSpec obtained from prepareMultiGen().
98 RooDataSet* generate(GenSpec&) const ;
99
100
101 ////////////////////////////////////////////////////////////////////////////////
102 /// As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&)
103 /// \param[in] nEvents How many events to generate
104 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg& arg1,
105 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
106 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) const {
107 return generateBinned(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5);
108 }
109 virtual RooDataHist *generateBinned(const RooArgSet &whatVars,
110 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
111 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
112 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) const;
113 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, Bool_t expectedData=kFALSE, Bool_t extended=kFALSE) const;
114
115 virtual RooDataSet* generateSimGlobal(const RooArgSet& whatVars, Int_t nEvents) ;
116
117 ///Helper calling plotOn(RooPlot*, RooLinkedList&) const
118 virtual RooPlot* plotOn(RooPlot* frame,
119 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
120 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
121 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
122 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none(),
123 const RooCmdArg& arg9=RooCmdArg::none(), const RooCmdArg& arg10=RooCmdArg::none()
124 ) const {
125 return RooAbsReal::plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9,arg10) ;
126 }
127 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
128
129 /// Add a box with parameter values (and errors) to the specified frame
130 virtual RooPlot* paramOn(RooPlot* frame,
131 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
132 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
133 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
134 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
135
136 virtual RooPlot* paramOn(RooPlot* frame, const RooAbsData* data, const char *label= "", Int_t sigDigits = 2,
137 Option_t *options = "NELU", Double_t xmin=0.50,
138 Double_t xmax= 0.99,Double_t ymax=0.95) ;
139
140 // Built-in generator support
141 virtual Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
142 virtual void initGenerator(Int_t code) ;
143 virtual void generateEvent(Int_t code);
144 virtual Bool_t isDirectGenSafe(const RooAbsArg& arg) const ;
145
146 // Configuration of MC generators used for this pdf
147 const RooNumGenConfig* getGeneratorConfig() const ;
151 void setGeneratorConfig() ;
152 void setGeneratorConfig(const RooNumGenConfig& config) ;
153
154 // -log(L) fits to binned and unbinned data
155 virtual RooFitResult* fitTo(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
156 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
157 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
158 virtual RooFitResult* fitTo(RooAbsData& data, const RooLinkedList& cmdList) ;
159
160 virtual RooAbsReal* createNLL(RooAbsData& data, const RooLinkedList& cmdList) ;
161 virtual RooAbsReal* createNLL(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
162 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
163 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
164
165 // Chi^2 fits to histograms
168 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
169 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
170 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
171 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
172
173 // Chi^2 fits to X-Y datasets
174 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
175
176
177
178
179
180 // Constraint management
181 virtual RooArgSet* getConstraints(const RooArgSet& /*observables*/, RooArgSet& /*constrainedParams*/, Bool_t /*stripDisconnected*/) const {
182 // Interface to retrieve constraint terms on this pdf. Default implementation returns null
183 return 0 ;
184 }
185 virtual RooArgSet* getAllConstraints(const RooArgSet& observables, RooArgSet& constrainedParams, Bool_t stripDisconnected=kTRUE) const ;
186
187 // Project p.d.f into lower dimensional p.d.f
188 virtual RooAbsPdf* createProjection(const RooArgSet& iset) ;
189
190 // Create cumulative density function from p.d.f
191 RooAbsReal* createCdf(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
192 RooAbsReal* createCdf(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
193 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
194 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
195 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
196 RooAbsReal* createScanCdf(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
197
198 // Function evaluation support
199 virtual Double_t getValV(const RooArgSet* set=0) const ;
200 virtual Double_t getLogVal(const RooArgSet* set=0) const ;
201
202 RooSpan<const double> getValBatch(std::size_t begin, std::size_t batchSize,
203 const RooArgSet* normSet = nullptr) const final;
204 RooSpan<const double> getLogValBatch(std::size_t begin, std::size_t batchSize,
205 const RooArgSet* normSet = nullptr) const;
206
207 Double_t getNorm(const RooArgSet& nset) const {
208 // Get p.d.f normalization term needed for observables 'nset'
209 return getNorm(&nset) ;
210 }
211 virtual Double_t getNorm(const RooArgSet* set=0) const ;
212
213 virtual void resetErrorCounters(Int_t resetValue=10) ;
214 void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE) ;
215private:
216 Bool_t traceEvalPdf(Double_t value) const;
217
218public:
219
220 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
221
222 virtual Bool_t selfNormalized() const {
223 // If true, p.d.f is taken as self-normalized and no attempt is made to add a normalization term
224 // This default implementation return false
225 return kFALSE ;
226 }
227
228 // Support for extended maximum likelihood, switched off by default
230 virtual ExtendMode extendMode() const {
231 // Returns ability of p.d.f to provided extended likelihood terms. Possible
232 // answers are CanNotBeExtended, CanBeExtended or MustBeExtended. This
233 // default implementation always return CanNotBeExtended
234 return CanNotBeExtended ;
235 }
236 inline Bool_t canBeExtended() const {
237 // If true p.d.f can provide extended likelihood term
238 return (extendMode() != CanNotBeExtended) ;
239 }
240 inline Bool_t mustBeExtended() const {
241 // If true p.d.f must extended likelihood term
242 return (extendMode() == MustBeExtended) ;
243 }
244 virtual Double_t expectedEvents(const RooArgSet* nset) const ;
245 virtual Double_t expectedEvents(const RooArgSet& nset) const {
246 // Return expecteded number of p.d.fs to be used in calculated of extended likelihood
247 return expectedEvents(&nset) ;
248 }
249
250 // Printing interface (human readable)
251 virtual void printValue(std::ostream& os) const ;
252 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
253
254 static void verboseEval(Int_t stat) ;
255 static int verboseEval() ;
256
257 virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet* nset=0) const ;
258
259 void setNormRange(const char* rangeName) ;
260 const char* normRange() const {
261 return _normRange.Length()>0 ? _normRange.Data() : 0 ;
262 }
263 void setNormRangeOverride(const char* rangeName) ;
264
265 const RooAbsReal* getNormIntegral(const RooArgSet& nset) const { return getNormObj(0,&nset,0) ; }
266
267protected:
268
269public:
270 virtual const RooAbsReal* getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName=0) const ;
271
272 virtual RooAbsGenContext* binnedGenContext(const RooArgSet &vars, Bool_t verbose= kFALSE) const ;
273
274 virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=0,
275 const RooArgSet* auxProto=0, Bool_t verbose= kFALSE) const ;
276
277 virtual RooAbsGenContext* autoGenContext(const RooArgSet &vars, const RooDataSet* prototype=0, const RooArgSet* auxProto=0,
278 Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE, const char* binnedTag="") const ;
279
280private:
281
282 RooDataSet *generate(RooAbsGenContext& context, const RooArgSet& whatVars, const RooDataSet* prototype,
283 Double_t nEvents, Bool_t verbose, Bool_t randProtoOrder, Bool_t resampleProto, Bool_t skipInit=kFALSE,
284 Bool_t extended=kFALSE) const ;
285
286 // Implementation version
287 virtual RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants=kFALSE,
288 const char *label= "", Int_t sigDigits = 2, Option_t *options = "NELU", Double_t xmin=0.65,
289 Double_t xmax= 0.99,Double_t ymax=0.95, const RooCmdArg* formatCmd=0) ;
290
291 void logBatchComputationErrors(RooSpan<const double>& outputs, std::size_t begin) const;
292
293protected:
294 virtual RooPlot *plotOn(RooPlot *frame, PlotOpt o) const;
295
296 friend class RooEffGenContext ;
297 friend class RooAddGenContext ;
298 friend class RooProdGenContext ;
299 friend class RooSimGenContext ;
301 friend class RooConvGenContext ;
302 friend class RooSimultaneous ;
303 friend class RooAddGenContextOrig ;
304 friend class RooProdPdf ;
305 friend class RooMCStudy ;
306
307 Int_t* randomizeProtoOrder(Int_t nProto,Int_t nGen,Bool_t resample=kFALSE) const ;
308
309 friend class RooExtendPdf ;
310 // This also forces the definition of a copy ctor in derived classes
311 RooAbsPdf(const RooAbsPdf& other, const char* name = 0);
312
313 friend class RooRealIntegral ;
315
316 virtual Bool_t syncNormalization(const RooArgSet* dset, Bool_t adjustProxies=kTRUE) const ;
317
318 friend class RooAbsAnaConvPdf ;
320 mutable RooAbsReal* _norm ; //! Normalization integral (owned by _normMgr)
321 mutable RooArgSet* _normSet ; //! Normalization set with for above integral
322
324 public:
325 CacheElem(RooAbsReal& norm) : _norm(&norm) {} ;
327 virtual ~CacheElem() ;
330 } ;
331 mutable RooObjCacheManager _normMgr ; // The cache manager
332
333 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
334
336 // Hook function intercepting redirectServer calls. Discard current normalization
337 // object if any server is redirected
338
339 // Object is own by _normCacheManager that will delete object as soon as cache
340 // is sterilized by server redirect
341 _norm = 0 ;
342 return kFALSE ;
343 } ;
344
345
346 mutable Int_t _errorCount ; // Number of errors remaining to print
347 mutable Int_t _traceCount ; // Number of traces remaining to print
348 mutable Int_t _negCount ; // Number of negative probablities remaining to print
349
350 Bool_t _selectComp ; // Component selection flag for RooAbsPdf::plotCompOn
351
352 RooNumGenConfig* _specGeneratorConfig ; //! MC generator configuration specific for this object
353
354 TString _normRange ; // Normalization range
356
357 ClassDef(RooAbsPdf,4) // Abstract PDF with normalization support
358};
359
360
361#endif
const Bool_t kFALSE
Definition: RtypesCore.h:90
double Double_t
Definition: RtypesCore.h:57
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassDef(name, id)
Definition: Rtypes.h:322
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:73
friend class RooArgSet
Definition: RooAbsArg.h:572
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:44
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Normalization set with for above integral.
Definition: RooAbsPdf.h:323
virtual ~CacheElem()
Destructor of normalization cache element.
Definition: RooAbsPdf.cxx:3324
virtual RooArgList containedArgs(Action)
Definition: RooAbsPdf.h:328
CacheElem(RooAbsReal &norm)
Definition: RooAbsPdf.h:325
void operModeHook(RooAbsArg::OperMode)
Dummy implementation.
Definition: RooAbsPdf.cxx:3313
RooAbsReal * _norm
Definition: RooAbsPdf.h:329
RooArgSet _whatVars
Definition: RooAbsPdf.h:81
RooAbsGenContext * _genContext
Definition: RooAbsPdf.h:80
Bool_t _resampleProto
Definition: RooAbsPdf.h:86
GenSpec(const GenSpec &other)
TString _dsetName
Definition: RooAbsPdf.h:87
RooDataSet * _protoData
Definition: RooAbsPdf.h:82
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:2160
RooSpan< const double > getLogValBatch(std::size_t begin, std::size_t batchSize, const RooArgSet *normSet=nullptr) const
Compute the log-likelihoods for all events in the requested batch.
Definition: RooAbsPdf.cxx:741
RooObjCacheManager _normMgr
Definition: RooAbsPdf.h:331
Double_t getNorm(const RooArgSet &nset) const
Definition: RooAbsPdf.h:207
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition: RooAbsPdf.cxx:2418
virtual ~RooAbsPdf()
Destructor.
Definition: RooAbsPdf.cxx:262
void logBatchComputationErrors(RooSpan< const double > &outputs, std::size_t begin) const
Scan through outputs and fix+log all nans and negative values.
Definition: RooAbsPdf.cxx:718
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
Definition: RooAbsPdf.cxx:3562
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:622
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:1850
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooLinkedList &cmdList)
Calls RooAbsPdf::createChi2(RooDataSet& data, const RooLinkedList& cmdList) and returns fit result.
Definition: RooAbsPdf.cxx:1805
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:3303
Bool_t traceEvalPdf(Double_t value) const
Check that passed value is positive and not 'not-a-number'.
Definition: RooAbsPdf.cxx:430
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:918
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:1998
RooAbsReal * createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
Definition: RooAbsPdf.cxx:3458
TString _normRange
MC generator configuration specific for this object.
Definition: RooAbsPdf.h:354
void setNormRange(const char *rangeName)
Definition: RooAbsPdf.cxx:3593
const RooAbsReal * getNormIntegral(const RooArgSet &nset) const
Definition: RooAbsPdf.h:265
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()) const
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition: RooAbsPdf.h:104
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:2396
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Definition: RooAbsPdf.cxx:3507
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:236
virtual Bool_t selfNormalized() const
Definition: RooAbsPdf.h:222
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:1971
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:410
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:55
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:2007
Int_t _traceCount
Definition: RooAbsPdf.h:347
virtual Double_t expectedEvents(const RooArgSet &nset) const
Definition: RooAbsPdf.h:245
RooAbsReal * _norm
Definition: RooAbsPdf.h:320
virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet *nset=0) const
Return the extended likelihood term ( ) of this PDF for the given number of observed events.
Definition: RooAbsPdf.cxx:778
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:1955
virtual RooArgSet * getAllConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, Bool_t stripDisconnected=kTRUE) const
This helper function finds and collects all constraints terms of all component p.d....
Definition: RooAbsPdf.cxx:3475
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:1254
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:523
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:280
Int_t _errorCount
Definition: RooAbsPdf.h:343
Bool_t _selectComp
Definition: RooAbsPdf.h:350
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet &, Bool_t) const
Definition: RooAbsPdf.h:181
RooSpan< const double > getValBatch(std::size_t begin, std::size_t batchSize, const RooArgSet *normSet=nullptr) const final
Compute batch of values for given range, and normalise by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:347
@ CanBeExtended
Definition: RooAbsPdf.h:229
@ MustBeExtended
Definition: RooAbsPdf.h:229
@ CanNotBeExtended
Definition: RooAbsPdf.h:229
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:3132
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:634
const char * normRange() const
Definition: RooAbsPdf.h:260
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:3380
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:488
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:660
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:2352
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:3283
Int_t _negCount
Definition: RooAbsPdf.h:348
Bool_t mustBeExtended() const
Definition: RooAbsPdf.h:240
void setNormRangeOverride(const char *rangeName)
Definition: RooAbsPdf.cxx:3610
virtual RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
Definition: RooAbsPdf.cxx:2654
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:2431
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
Definition: RooAbsPdf.cxx:3535
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2406
virtual ExtendMode extendMode() const
Definition: RooAbsPdf.h:230
RooAbsPdf()
Default constructor.
Definition: RooAbsPdf.cxx:204
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, Bool_t verbose=kFALSE) const
Return a binned generator context.
Definition: RooAbsPdf.cxx:1988
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
Definition: RooAbsPdf.cxx:3497
RooNumGenConfig * _specGeneratorConfig
Definition: RooAbsPdf.h:352
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:321
friend class RooAddGenContextOrig
Definition: RooAbsPdf.h:303
Double_t _rawValue
Definition: RooAbsPdf.h:319
static TString _normRangeOverride
Definition: RooAbsPdf.h:355
static Int_t _verboseEval
Definition: RooAbsPdf.h:314
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Definition: RooAbsPdf.h:335
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:3342
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:118
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
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 is an object to represent discrete states.
Definition: RooCategory.h:23
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:28
static const RooCmdArg & none()
Return reference to null argument.
Definition: RooCmdArg.cxx:52
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:33
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:35
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:44
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.
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.
A simple container to hold a batch of data values.
Definition: RooSpan.h:32
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:251
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)
EvaluateInfo init(std::vector< RooRealProxy > parameters, std::vector< ArrayWrapper * > wrappers, std::vector< double * > arrays, size_t begin, size_t batchSize)