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 Bool_t R__DEPRECATED(6,22,"Call traceEvalPdf() instead.") traceEvalHook(Double_t value) const ;
200 virtual Double_t getValV(const RooArgSet* set=0) const ;
201 virtual Double_t getLogVal(const RooArgSet* set=0) const ;
202
203 virtual RooSpan<const double> getValBatch(std::size_t begin, std::size_t batchSize,
204 const RooArgSet* normSet = nullptr) const;
205 RooSpan<const double> getLogValBatch(std::size_t begin, std::size_t batchSize,
206 const RooArgSet* normSet = nullptr) const;
207
208 Double_t getNorm(const RooArgSet& nset) const {
209 // Get p.d.f normalization term needed for observables 'nset'
210 return getNorm(&nset) ;
211 }
212 virtual Double_t getNorm(const RooArgSet* set=0) const ;
213
214 virtual void resetErrorCounters(Int_t resetValue=10) ;
215 void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE) ;
216private:
217 Bool_t traceEvalPdf(Double_t value) const;
218
219public:
220
221 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
222
223 virtual Bool_t selfNormalized() const {
224 // If true, p.d.f is taken as self-normalized and no attempt is made to add a normalization term
225 // This default implementation return false
226 return kFALSE ;
227 }
228
229 // Support for extended maximum likelihood, switched off by default
231 virtual ExtendMode extendMode() const {
232 // Returns ability of p.d.f to provided extended likelihood terms. Possible
233 // answers are CanNotBeExtended, CanBeExtended or MustBeExtended. This
234 // default implementation always return CanNotBeExtended
235 return CanNotBeExtended ;
236 }
237 inline Bool_t canBeExtended() const {
238 // If true p.d.f can provide extended likelihood term
239 return (extendMode() != CanNotBeExtended) ;
240 }
241 inline Bool_t mustBeExtended() const {
242 // If true p.d.f must extended likelihood term
243 return (extendMode() == MustBeExtended) ;
244 }
245 virtual Double_t expectedEvents(const RooArgSet* nset) const ;
246 virtual Double_t expectedEvents(const RooArgSet& nset) const {
247 // Return expecteded number of p.d.fs to be used in calculated of extended likelihood
248 return expectedEvents(&nset) ;
249 }
250
251 // Printing interface (human readable)
252 virtual void printValue(std::ostream& os) const ;
253 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
254
255 static void verboseEval(Int_t stat) ;
256 static int verboseEval() ;
257
258 virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet* nset=0) const ;
259
260 void setNormRange(const char* rangeName) ;
261 const char* normRange() const {
262 return _normRange.Length()>0 ? _normRange.Data() : 0 ;
263 }
264 void setNormRangeOverride(const char* rangeName) ;
265
266 const RooAbsReal* getNormIntegral(const RooArgSet& nset) const { return getNormObj(0,&nset,0) ; }
267
268protected:
269
270public:
271 virtual const RooAbsReal* getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName=0) const ;
272
273 virtual RooAbsGenContext* binnedGenContext(const RooArgSet &vars, Bool_t verbose= kFALSE) const ;
274
275 virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=0,
276 const RooArgSet* auxProto=0, Bool_t verbose= kFALSE) const ;
277
278 virtual RooAbsGenContext* autoGenContext(const RooArgSet &vars, const RooDataSet* prototype=0, const RooArgSet* auxProto=0,
279 Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE, const char* binnedTag="") const ;
280
281private:
282
283 RooDataSet *generate(RooAbsGenContext& context, const RooArgSet& whatVars, const RooDataSet* prototype,
284 Double_t nEvents, Bool_t verbose, Bool_t randProtoOrder, Bool_t resampleProto, Bool_t skipInit=kFALSE,
285 Bool_t extended=kFALSE) const ;
286
287 // Implementation version
288 virtual RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants=kFALSE,
289 const char *label= "", Int_t sigDigits = 2, Option_t *options = "NELU", Double_t xmin=0.65,
290 Double_t xmax= 0.99,Double_t ymax=0.95, const RooCmdArg* formatCmd=0) ;
291
292 void logBatchComputationErrors(RooSpan<const double>& outputs, std::size_t begin) const;
293
294protected:
295 virtual RooPlot *plotOn(RooPlot *frame, PlotOpt o) const;
296
297 friend class RooEffGenContext ;
298 friend class RooAddGenContext ;
299 friend class RooProdGenContext ;
300 friend class RooSimGenContext ;
302 friend class RooConvGenContext ;
303 friend class RooSimultaneous ;
304 friend class RooAddGenContextOrig ;
305 friend class RooProdPdf ;
306 friend class RooMCStudy ;
307
308 Int_t* randomizeProtoOrder(Int_t nProto,Int_t nGen,Bool_t resample=kFALSE) const ;
309
310 friend class RooExtendPdf ;
311 // This also forces the definition of a copy ctor in derived classes
312 RooAbsPdf(const RooAbsPdf& other, const char* name = 0);
313
314 friend class RooRealIntegral ;
316
317 virtual Bool_t syncNormalization(const RooArgSet* dset, Bool_t adjustProxies=kTRUE) const ;
318
319 friend class RooAbsAnaConvPdf ;
321 mutable RooAbsReal* _norm ; //! Normalization integral (owned by _normMgr)
322 mutable RooArgSet* _normSet ; //! Normalization set with for above integral
323
325 public:
326 CacheElem(RooAbsReal& norm) : _norm(&norm) {} ;
328 virtual ~CacheElem() ;
331 } ;
332 mutable RooObjCacheManager _normMgr ; // The cache manager
333
334 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
335
337 // Hook function intercepting redirectServer calls. Discard current normalization
338 // object if any server is redirected
339
340 // Object is own by _normCacheManager that will delete object as soon as cache
341 // is sterilized by server redirect
342 _norm = 0 ;
343 return kFALSE ;
344 } ;
345
346
347 mutable Int_t _errorCount ; // Number of errors remaining to print
348 mutable Int_t _traceCount ; // Number of traces remaining to print
349 mutable Int_t _negCount ; // Number of negative probablities remaining to print
350
351 Bool_t _selectComp ; // Component selection flag for RooAbsPdf::plotCompOn
352
353 RooNumGenConfig* _specGeneratorConfig ; //! MC generator configuration specific for this object
354
355 TString _normRange ; // Normalization range
357
358 ClassDef(RooAbsPdf,4) // Abstract PDF with normalization support
359};
360
361
362#endif
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:71
friend class RooArgSet
Definition: RooAbsArg.h:551
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:39
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Normalization set with for above integral.
Definition: RooAbsPdf.h:324
virtual ~CacheElem()
Destructor of normalization cache element.
Definition: RooAbsPdf.cxx:3344
virtual RooArgList containedArgs(Action)
Definition: RooAbsPdf.h:329
CacheElem(RooAbsReal &norm)
Definition: RooAbsPdf.h:326
void operModeHook(RooAbsArg::OperMode)
Dummy implementation.
Definition: RooAbsPdf.cxx:3333
RooAbsReal * _norm
Definition: RooAbsPdf.h:330
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:2195
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:781
RooObjCacheManager _normMgr
Definition: RooAbsPdf.h:332
Double_t getNorm(const RooArgSet &nset) const
Definition: RooAbsPdf.h:208
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition: RooAbsPdf.cxx:2453
virtual ~RooAbsPdf()
Destructor.
Definition: RooAbsPdf.cxx:263
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:758
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
Definition: RooAbsPdf.cxx:3582
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:662
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:1882
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooLinkedList &cmdList)
Calls RooAbsPdf::createChi2(RooDataSet& data, const RooLinkedList& cmdList) and returns fit result.
Definition: RooAbsPdf.cxx:1837
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:3323
Bool_t traceEvalPdf(Double_t value) const
Check that passed value is positive and not 'not-a-number'.
Definition: RooAbsPdf.cxx:437
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:957
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:2033
RooAbsReal * createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
Definition: RooAbsPdf.cxx:3478
TString _normRange
MC generator configuration specific for this object.
Definition: RooAbsPdf.h:355
void setNormRange(const char *rangeName)
Definition: RooAbsPdf.cxx:3613
const RooAbsReal * getNormIntegral(const RooArgSet &nset) const
Definition: RooAbsPdf.h:266
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:2431
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Definition: RooAbsPdf.cxx:3527
virtual RooSpan< const double > getValBatch(std::size_t begin, std::size_t batchSize, const RooArgSet *normSet=nullptr) const
Compute batch of values for given range, and normalise by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:348
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:237
virtual Bool_t selfNormalized() const
Definition: RooAbsPdf.h:223
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:2006
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:417
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:2042
Int_t _traceCount
Definition: RooAbsPdf.h:348
virtual Double_t expectedEvents(const RooArgSet &nset) const
Definition: RooAbsPdf.h:246
RooAbsReal * _norm
Definition: RooAbsPdf.h:321
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:817
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:1990
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:3495
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:1286
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:530
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:281
Int_t _errorCount
Definition: RooAbsPdf.h:344
Bool_t _selectComp
Definition: RooAbsPdf.h:351
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet &, Bool_t) const
Definition: RooAbsPdf.h:181
@ CanBeExtended
Definition: RooAbsPdf.h:230
@ MustBeExtended
Definition: RooAbsPdf.h:230
@ CanNotBeExtended
Definition: RooAbsPdf.h:230
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:3152
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:674
const char * normRange() const
Definition: RooAbsPdf.h:261
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:3400
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:495
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:700
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:2387
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:3303
Int_t _negCount
Definition: RooAbsPdf.h:349
virtual Bool_t R__DEPRECATED(6, 22,"Call traceEvalPdf() instead.") traceEvalHook(Double_t value) const
Bool_t mustBeExtended() const
Definition: RooAbsPdf.h:241
void setNormRangeOverride(const char *rangeName)
Definition: RooAbsPdf.cxx:3630
virtual RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
Definition: RooAbsPdf.cxx:2689
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:2466
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
Definition: RooAbsPdf.cxx:3555
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2441
virtual ExtendMode extendMode() const
Definition: RooAbsPdf.h:231
RooAbsPdf()
Default constructor.
Definition: RooAbsPdf.cxx:205
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, Bool_t verbose=kFALSE) const
Return a binned generator context.
Definition: RooAbsPdf.cxx:2023
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
Definition: RooAbsPdf.cxx:3517
RooNumGenConfig * _specGeneratorConfig
Definition: RooAbsPdf.h:353
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:322
friend class RooAddGenContextOrig
Definition: RooAbsPdf.h:304
Double_t _rawValue
Definition: RooAbsPdf.h:320
static TString _normRangeOverride
Definition: RooAbsPdf.h:356
static Int_t _verboseEval
Definition: RooAbsPdf.h:315
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Definition: RooAbsPdf.h:336
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:3362
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:59
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: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: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: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
EvaluateInfo init(std::vector< RooRealProxy > parameters, std::vector< ArrayWrapper * > wrappers, std::vector< double * > arrays, size_t begin, size_t batchSize)
RooCmdArg NumEvents(Int_t numEvents)