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 "RooObjCacheManager.h"
21#include "RooCmdArg.h"
22
23class RooDataSet;
24class RooDataHist ;
25class RooArgSet ;
26class RooAbsGenContext ;
27class RooFitResult ;
28class RooExtendPdf ;
29class RooCategory ;
30class TPaveText;
31class TH1F;
32class TH2F;
33class TList ;
34class RooLinkedList ;
35class RooMinimizer ;
36class RooNumGenConfig ;
37class RooRealIntegral ;
38namespace RooBatchCompute {
39struct RunContext;
40}
41
42class RooAbsPdf : public RooAbsReal {
43public:
44
45 // Constructors, assignment etc
46 RooAbsPdf() ;
47 RooAbsPdf(const char *name, const char *title=0) ;
48 RooAbsPdf(const char *name, const char *title, Double_t minVal, Double_t maxVal) ;
49 // RooAbsPdf(const RooAbsPdf& other, const char* name=0);
50 virtual ~RooAbsPdf();
51
52 // Toy MC generation
53
54 ////////////////////////////////////////////////////////////////////////////////
55 /// See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&)
56 /// \param[in] nEvents How many events to generate
57 RooDataSet *generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg& arg1,
58 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
59 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) {
60 return generate(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5) ;
61 }
62 RooDataSet *generate(const RooArgSet &whatVars,
63 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
64 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
65 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
66 RooDataSet *generate(const RooArgSet &whatVars, Double_t nEvents = 0, Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE,
67 const char* binnedTag="", Bool_t expectedData=kFALSE, Bool_t extended = kFALSE) const;
68 RooDataSet *generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents= 0,
69 Bool_t verbose=kFALSE, Bool_t randProtoOrder=kFALSE, Bool_t resampleProto=kFALSE) const;
70
71
72 class GenSpec {
73 public:
74 virtual ~GenSpec() ;
76 private:
77 GenSpec(RooAbsGenContext* context, const RooArgSet& whatVars, RooDataSet* protoData, Int_t nGen, Bool_t extended,
78 Bool_t randProto, Bool_t resampleProto, TString dsetName, Bool_t init=kFALSE) ;
79 GenSpec(const GenSpec& other) ;
80
81 friend class RooAbsPdf ;
91 ClassDef(GenSpec,0) // Generation specification
92 } ;
93
94 ///Prepare GenSpec configuration object for efficient generation of multiple datasets from identical specification.
95 GenSpec* prepareMultiGen(const RooArgSet &whatVars,
96 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
97 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
98 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
99 ///Generate according to GenSpec obtained from prepareMultiGen().
100 RooDataSet* generate(GenSpec&) const ;
101
102
103 ////////////////////////////////////////////////////////////////////////////////
104 /// As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&)
105 /// \param[in] nEvents How many events to generate
106 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg& arg1,
107 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
108 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) const {
109 return generateBinned(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5);
110 }
111 virtual RooDataHist *generateBinned(const RooArgSet &whatVars,
112 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
113 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
114 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) const;
115 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, Bool_t expectedData=kFALSE, Bool_t extended=kFALSE) const;
116
117 virtual RooDataSet* generateSimGlobal(const RooArgSet& whatVars, Int_t nEvents) ;
118
119 ///Helper calling plotOn(RooPlot*, RooLinkedList&) const
120 virtual RooPlot* plotOn(RooPlot* frame,
121 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
122 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
123 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
124 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none(),
125 const RooCmdArg& arg9=RooCmdArg::none(), const RooCmdArg& arg10=RooCmdArg::none()
126 ) const {
127 return RooAbsReal::plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9,arg10) ;
128 }
129 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
130
131 /// Add a box with parameter values (and errors) to the specified frame
132 virtual RooPlot* paramOn(RooPlot* frame,
133 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
134 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
135 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
136 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
137
138 virtual RooPlot* paramOn(RooPlot* frame, const RooAbsData* data, const char *label= "", Int_t sigDigits = 2,
139 Option_t *options = "NELU", Double_t xmin=0.65,
140 Double_t xmax = 0.9, Double_t ymax = 0.9) ;
141
142 // Built-in generator support
143 virtual Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
144 virtual void initGenerator(Int_t code) ;
145 virtual void generateEvent(Int_t code);
146 virtual Bool_t isDirectGenSafe(const RooAbsArg& arg) const ;
147
148 // Configuration of MC generators used for this pdf
149 const RooNumGenConfig* getGeneratorConfig() const ;
153 void setGeneratorConfig() ;
154 void setGeneratorConfig(const RooNumGenConfig& config) ;
155
156 // -log(L) fits to binned and unbinned data
157 virtual RooFitResult* fitTo(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
158 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
159 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
160 virtual RooFitResult* fitTo(RooAbsData& data, const RooLinkedList& cmdList) ;
161
162 virtual RooAbsReal* createNLL(RooAbsData& data, const RooLinkedList& cmdList) ;
163 virtual RooAbsReal* createNLL(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
164 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
165 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
166
167 // Chi^2 fits to histograms
170 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
171 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
172 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
173 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
174
175 // Chi^2 fits to X-Y datasets
176 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
177
178
179
180
181
182 // Constraint management
183 virtual RooArgSet* getConstraints(const RooArgSet& /*observables*/, RooArgSet& /*constrainedParams*/, Bool_t /*stripDisconnected*/) const {
184 // Interface to retrieve constraint terms on this pdf. Default implementation returns null
185 return 0 ;
186 }
187 virtual RooArgSet* getAllConstraints(const RooArgSet& observables, RooArgSet& constrainedParams, Bool_t stripDisconnected=kTRUE) const ;
188
189 // Project p.d.f into lower dimensional p.d.f
190 virtual RooAbsPdf* createProjection(const RooArgSet& iset) ;
191
192 // Create cumulative density function from p.d.f
193 RooAbsReal* createCdf(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
194 RooAbsReal* createCdf(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
195 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
196 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
197 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
198 RooAbsReal* createScanCdf(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
199
200 // Function evaluation support
201 virtual Double_t getValV(const RooArgSet* set=0) const ;
202 virtual Double_t getLogVal(const RooArgSet* set=0) const ;
203
205 RooSpan<const double> getLogValBatch(std::size_t begin, std::size_t batchSize,
206 const RooArgSet* normSet = nullptr) const;
207 RooSpan<const double> getLogProbabilities(RooBatchCompute::RunContext& evalData, const RooArgSet* normSet = nullptr) const;
208
209 /// \copydoc getNorm(const RooArgSet*) const
210 Double_t getNorm(const RooArgSet& nset) const {
211 return getNorm(&nset) ;
212 }
213 virtual Double_t getNorm(const RooArgSet* set=0) const ;
214
215 virtual void resetErrorCounters(Int_t resetValue=10) ;
216 void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE) ;
217private:
218 Bool_t traceEvalPdf(Double_t value) const;
219
220public:
221
222 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
223
224 /// Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
225 /// Always returns false, unless a PDF overrides this function.
226 virtual Bool_t selfNormalized() const {
227 return kFALSE ;
228 }
229
230 // Support for extended maximum likelihood, switched off by default
232 /// Returns ability of PDF to provide extended likelihood terms. Possible
233 /// answers are in the enumerator RooAbsPdf::ExtendMode.
234 /// This default implementation always returns CanNotBeExtended.
235 virtual ExtendMode extendMode() const { return CanNotBeExtended; }
236 /// If true, PDF can provide extended likelihood term.
237 inline Bool_t canBeExtended() const {
238 return (extendMode() != CanNotBeExtended) ;
239 }
240 /// If true PDF must provide extended likelihood term.
241 inline Bool_t mustBeExtended() const {
242 return (extendMode() == MustBeExtended) ;
243 }
244 /// Return expected number of events to be used in calculation of extended
245 /// likelihood.
246 virtual Double_t expectedEvents(const RooArgSet* nset) const ;
247 /// Return expected number of events to be used in calculation of extended
248 /// likelihood. This function should not be overridden, as it just redirects
249 /// to the actual virtual function but takes a RooArgSet reference instead of
250 /// pointer (\see expectedEvents(const RooArgSet*) const).
251 double expectedEvents(const RooArgSet& nset) const {
252 return expectedEvents(&nset) ;
253 }
254
255 // Printing interface (human readable)
256 virtual void printValue(std::ostream& os) const ;
257 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
258
259 static void verboseEval(Int_t stat) ;
260 static int verboseEval() ;
261
262 virtual Double_t extendedTerm(Double_t observedEvents, const RooArgSet* nset=0) const ;
263
264 void setNormRange(const char* rangeName) ;
265 const char* normRange() const {
266 return _normRange.Length()>0 ? _normRange.Data() : 0 ;
267 }
268 void setNormRangeOverride(const char* rangeName) ;
269
270 const RooAbsReal* getNormIntegral(const RooArgSet& nset) const { return getNormObj(0,&nset,0) ; }
271
272protected:
273
274public:
275 virtual const RooAbsReal* getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName=0) const ;
276
277 virtual RooAbsGenContext* binnedGenContext(const RooArgSet &vars, Bool_t verbose= kFALSE) const ;
278
279 virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=0,
280 const RooArgSet* auxProto=0, Bool_t verbose= kFALSE) const ;
281
282 virtual RooAbsGenContext* autoGenContext(const RooArgSet &vars, const RooDataSet* prototype=0, const RooArgSet* auxProto=0,
283 Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE, const char* binnedTag="") const ;
284
285private:
286
287 RooDataSet *generate(RooAbsGenContext& context, const RooArgSet& whatVars, const RooDataSet* prototype,
288 Double_t nEvents, Bool_t verbose, Bool_t randProtoOrder, Bool_t resampleProto, Bool_t skipInit=kFALSE,
289 Bool_t extended=kFALSE) const ;
290
291 // Implementation version
292 virtual RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants=kFALSE,
293 const char *label= "", Int_t sigDigits = 2, Option_t *options = "NELU", Double_t xmin=0.65,
294 Double_t xmax= 0.99,Double_t ymax=0.95, const RooCmdArg* formatCmd=0) ;
295
296 void logBatchComputationErrors(RooSpan<const double>& outputs, std::size_t begin) const;
297
298protected:
299 virtual RooPlot *plotOn(RooPlot *frame, PlotOpt o) const;
300
301 friend class RooEffGenContext ;
302 friend class RooAddGenContext ;
303 friend class RooProdGenContext ;
304 friend class RooSimGenContext ;
306 friend class RooConvGenContext ;
307 friend class RooSimultaneous ;
308 friend class RooAddGenContextOrig ;
309 friend class RooProdPdf ;
310 friend class RooMCStudy ;
311
312 Int_t* randomizeProtoOrder(Int_t nProto,Int_t nGen,Bool_t resample=kFALSE) const ;
313
314 friend class RooExtendPdf ;
315 // This also forces the definition of a copy ctor in derived classes
316 RooAbsPdf(const RooAbsPdf& other, const char* name = 0);
317
318 friend class RooRealIntegral ;
320
321 virtual Bool_t syncNormalization(const RooArgSet* dset, Bool_t adjustProxies=kTRUE) const ;
322
323 friend class RooAbsAnaConvPdf ;
325 mutable RooAbsReal* _norm ; //! Normalization integral (owned by _normMgr)
326 mutable RooArgSet* _normSet ; //! Normalization set with for above integral
327
329 public:
330 CacheElem(RooAbsReal& norm) : _norm(&norm) {} ;
331 virtual ~CacheElem() ;
334 } ;
335 mutable RooObjCacheManager _normMgr ; // The cache manager
336
337 friend class CacheElem ; // Cache needs to be able to clear _norm pointer
338
340 // Hook function intercepting redirectServer calls. Discard current normalization
341 // object if any server is redirected
342
343 // Object is own by _normCacheManager that will delete object as soon as cache
344 // is sterilized by server redirect
345 _norm = 0 ;
346 return kFALSE ;
347 } ;
348
349
350 mutable Int_t _errorCount ; // Number of errors remaining to print
351 mutable Int_t _traceCount ; // Number of traces remaining to print
352 mutable Int_t _negCount ; // Number of negative probablities remaining to print
353
354 Bool_t _selectComp ; // Component selection flag for RooAbsPdf::plotCompOn
355
356 RooNumGenConfig* _specGeneratorConfig ; //! MC generator configuration specific for this object
357
358 TString _normRange ; // Normalization range
360
361private:
363 int calculateSumW2CorrectedCovMatrix(RooMinimizer& minimizer, RooAbsReal const& nll) const;
364
365 ClassDef(RooAbsPdf,4) // Abstract PDF with normalization support
366};
367
368
369#endif
int Int_t
Definition: CPyCppyy.h:43
const Bool_t kFALSE
Definition: RtypesCore.h:101
bool Bool_t
Definition: RtypesCore.h:63
double Double_t
Definition: RtypesCore.h:59
const Bool_t kTRUE
Definition: RtypesCore.h:100
const char Option_t
Definition: RtypesCore.h:66
#define ClassDef(name, id)
Definition: Rtypes.h:325
static void indent(ostringstream &buf, int indent_level)
float xmin
Definition: THbookFile.cxx:95
float xmax
Definition: THbookFile.cxx:95
float ymax
Definition: THbookFile.cxx:95
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 and a "shape" in RooFi...
Definition: RooAbsArg.h:72
friend class RooArgSet
Definition: RooAbsArg.h:600
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:67
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
Normalization set with for above integral.
Definition: RooAbsPdf.h:328
virtual ~CacheElem()
Destructor of normalization cache element.
Definition: RooAbsPdf.cxx:3251
virtual RooArgList containedArgs(Action)
Definition: RooAbsPdf.h:332
CacheElem(RooAbsReal &norm)
Definition: RooAbsPdf.h:330
RooAbsReal * _norm
Definition: RooAbsPdf.h:333
RooArgSet _whatVars
Definition: RooAbsPdf.h:83
RooAbsGenContext * _genContext
Definition: RooAbsPdf.h:82
Bool_t _resampleProto
Definition: RooAbsPdf.h:88
GenSpec(const GenSpec &other)
TString _dsetName
Definition: RooAbsPdf.h:89
RooDataSet * _protoData
Definition: RooAbsPdf.h:84
const RooAbsReal * getNormIntegral(const RooArgSet &nset) const
Definition: RooAbsPdf.h:270
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:2061
RooObjCacheManager _normMgr
Definition: RooAbsPdf.h:335
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:120
int calculateSumW2CorrectedCovMatrix(RooMinimizer &minimizer, RooAbsReal const &nll) const
Definition: RooAbsPdf.cxx:1218
Double_t getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition: RooAbsPdf.h:210
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition: RooAbsPdf.cxx:2319
virtual ~RooAbsPdf()
Destructor.
Definition: RooAbsPdf.cxx:264
double expectedEvents(const RooArgSet &nset) const
Return expected number of events to be used in calculation of extended likelihood.
Definition: RooAbsPdf.h:251
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:717
RooSpan< const double > getLogProbabilities(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
Compute the log-likelihoods for all events in the requested batch.
Definition: RooAbsPdf.cxx:740
RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute batch of values for given input data, and normalise by integrating over the observables in no...
Definition: RooAbsPdf.cxx:354
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
Definition: RooAbsPdf.cxx:3485
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:1738
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooLinkedList &cmdList)
Calls RooAbsPdf::createChi2(RooDataSet& data, const RooLinkedList& cmdList) and returns fit result.
Definition: RooAbsPdf.cxx:1692
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:3239
Bool_t traceEvalPdf(Double_t value) const
Check that passed value is positive and not 'not-a-number'.
Definition: RooAbsPdf.cxx:432
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Definition: RooAbsPdf.cxx:977
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:1889
RooAbsReal * createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
Definition: RooAbsPdf.cxx:3381
TString _normRange
MC generator configuration specific for this object.
Definition: RooAbsPdf.h:358
void setNormRange(const char *rangeName)
Definition: RooAbsPdf.cxx:3516
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:2297
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Definition: RooAbsPdf.cxx:3430
Bool_t canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition: RooAbsPdf.h:237
virtual Bool_t selfNormalized() const
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
Definition: RooAbsPdf.h:226
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:1862
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:412
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:1898
Int_t _traceCount
Definition: RooAbsPdf.h:351
RooAbsReal * _norm
Definition: RooAbsPdf.h:325
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:783
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:1843
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:3398
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:1410
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:526
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:282
Int_t _errorCount
Definition: RooAbsPdf.h:347
Bool_t _selectComp
Definition: RooAbsPdf.h:354
const char * normRange() const
Definition: RooAbsPdf.h:265
@ CanBeExtended
Definition: RooAbsPdf.h:231
@ MustBeExtended
Definition: RooAbsPdf.h:231
@ CanNotBeExtended
Definition: RooAbsPdf.h:231
int calculateAsymptoticCorrectedCovMatrix(RooMinimizer &minimizer, RooAbsData const &data)
Definition: RooAbsPdf.cxx:1151
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:3067
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
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:3304
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:491
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:2253
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
Definition: RooAbsPdf.cxx:3219
Int_t _negCount
Definition: RooAbsPdf.h:352
Bool_t mustBeExtended() const
If true PDF must provide extended likelihood term.
Definition: RooAbsPdf.h:241
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet &, Bool_t) const
Definition: RooAbsPdf.h:183
void setNormRangeOverride(const char *rangeName)
Definition: RooAbsPdf.cxx:3533
virtual RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
Definition: RooAbsPdf.cxx:2576
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:2332
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
Definition: RooAbsPdf.cxx:3458
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
Definition: RooAbsPdf.cxx:2307
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:106
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
Definition: RooAbsPdf.h:235
RooAbsPdf()
Default constructor.
Definition: RooAbsPdf.cxx:206
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, Bool_t verbose=kFALSE) const
Return a binned generator context.
Definition: RooAbsPdf.cxx:1879
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
Definition: RooAbsPdf.cxx:3420
RooNumGenConfig * _specGeneratorConfig
Definition: RooAbsPdf.h:356
RooSpan< const double > getLogValBatch(std::size_t begin, std::size_t batchSize, const RooArgSet *normSet=nullptr) const
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:326
friend class RooAddGenContextOrig
Definition: RooAbsPdf.h:308
Double_t _rawValue
Definition: RooAbsPdf.h:324
static TString _normRangeOverride
Definition: RooAbsPdf.h:359
static Int_t _verboseEval
Definition: RooAbsPdf.h:319
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Definition: RooAbsPdf.h:339
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:3269
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:57
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:61
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:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:35
RooCategory is an object to represent discrete states.
Definition: RooCategory.h:27
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: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:39
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:36
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:37
RooMCStudy is a helper class to facilitate Monte Carlo studies such as 'goodness-of-fit' studies,...
Definition: RooMCStudy.h:32
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Definition: RooMinimizer.h:49
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:37
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:34
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
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:136
Ssiz_t Length() const
Definition: TString.h:410
const char * Data() const
Definition: TString.h:369
RooCmdArg NumEvents(Int_t numEvents)
static const std::string name("name")
Namespace for dispatching RooFit computations to various backends.
This struct enables passing computation data around between elements of a computation graph.
Definition: RunContext.h:31