Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <RooDataHist.h>
21#include <RooDataSet.h>
22#include <RooFit/UniqueId.h>
23#include <RooGlobalFunc.h>
24#include <RooObjCacheManager.h>
25
26class RooArgSet ;
27class RooAbsGenContext ;
28class RooFitResult ;
29class RooExtendPdf ;
30class RooCategory ;
31class TPaveText;
32class TH1F;
33class TH2F;
34class TList ;
35class RooMinimizer ;
36class RooNumGenConfig ;
37class RooRealIntegral ;
38
39
40class RooAbsPdf : public RooAbsReal {
41public:
42
43 // Constructors, assignment etc
44 RooAbsPdf() ;
45 RooAbsPdf(const char *name, const char *title=nullptr) ;
46 RooAbsPdf(const char *name, const char *title, double minVal, double maxVal) ;
47 // RooAbsPdf(const RooAbsPdf& other, const char* name=nullptr);
48 ~RooAbsPdf() override;
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] whatVars Set of observables to generate for each event according to this model.
55 /// \param[in] nEvents How many events to generate
56 /// \param arg1,arg2,arg3,arg4,arg5 Optional command arguments.
57 RooFit::OwningPtr<RooDataSet> generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg& arg1,
58 const RooCmdArg& arg2={}, const RooCmdArg& arg3={},
59 const RooCmdArg& arg4={}, const RooCmdArg& arg5={}) {
60 return generate(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5) ;
61 }
63 const RooCmdArg& arg1={},const RooCmdArg& arg2={},
64 const RooCmdArg& arg3={},const RooCmdArg& arg4={},
65 const RooCmdArg& arg5={},const RooCmdArg& arg6={}) ;
66 RooFit::OwningPtr<RooDataSet> generate(const RooArgSet &whatVars, double nEvents = 0, bool verbose=false, bool autoBinned=true,
67 const char* binnedTag="", bool expectedData=false, bool extended = false) const;
68 RooFit::OwningPtr<RooDataSet> generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents= 0,
69 bool verbose=false, bool randProtoOrder=false, bool resampleProto=false) const;
70
71
72 class GenSpec {
73 public:
74 virtual ~GenSpec() ;
75 GenSpec() = default;
76
77 private:
78 GenSpec(RooAbsGenContext* context, const RooArgSet& whatVars, RooDataSet* protoData, Int_t nGen, bool extended,
79 bool randProto, bool resampleProto, TString dsetName, bool init=false) ;
80 GenSpec(const GenSpec& other) ;
81
82 friend class RooAbsPdf ;
83 std::unique_ptr<RooAbsGenContext> _genContext;
87 bool _extended = false;
88 bool _randProto = false;
89 bool _resampleProto = false;
91 bool _init = false;
92
93 ClassDef(GenSpec,0) // Generation specification
94 } ;
95
96 ///Prepare GenSpec configuration object for efficient generation of multiple datasets from identical specification.
97 GenSpec* prepareMultiGen(const RooArgSet &whatVars,
98 const RooCmdArg& arg1={},const RooCmdArg& arg2={},
99 const RooCmdArg& arg3={},const RooCmdArg& arg4={},
100 const RooCmdArg& arg5={},const RooCmdArg& arg6={}) ;
101 ///Generate according to GenSpec obtained from prepareMultiGen().
103
104
105 ////////////////////////////////////////////////////////////////////////////////
106 /// As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&) const.
107 /// \param[in] whatVars set
108 /// \param[in] nEvents How many events to generate
109 /// \param arg1,arg2,arg3,arg4,arg5 ordered arguments
110 virtual RooFit::OwningPtr<RooDataHist> generateBinned(const RooArgSet &whatVars, double nEvents, const RooCmdArg& arg1,
111 const RooCmdArg& arg2={}, const RooCmdArg& arg3={},
112 const RooCmdArg& arg4={}, const RooCmdArg& arg5={}) const {
113 return generateBinned(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5);
114 }
116 const RooCmdArg& arg1={},const RooCmdArg& arg2={},
117 const RooCmdArg& arg3={},const RooCmdArg& arg4={},
118 const RooCmdArg& arg5={},const RooCmdArg& arg6={}) const;
119 virtual RooFit::OwningPtr<RooDataHist> generateBinned(const RooArgSet &whatVars, double nEvents, bool expectedData=false, bool extended=false) const;
120
121 virtual RooFit::OwningPtr<RooDataSet> generateSimGlobal(const RooArgSet& whatVars, Int_t nEvents) ;
122
123 ///Helper calling plotOn(RooPlot*, RooLinkedList&) const
125 const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
126 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
127 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
128 const RooCmdArg& arg7={}, const RooCmdArg& arg8={},
129 const RooCmdArg& arg9={}, const RooCmdArg& arg10={}
130 ) const override {
131 return RooAbsReal::plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9,arg10) ;
132 }
133 RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const override ;
134
135 /// Add a box with parameter values (and errors) to the specified frame
136 virtual RooPlot* paramOn(RooPlot* frame,
137 const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
138 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
139 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
140 const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
141
142 // Built-in generator support
143 virtual Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK=true) const;
144 virtual void initGenerator(Int_t code) ;
145 virtual void generateEvent(Int_t code);
146 virtual bool isDirectGenSafe(const RooAbsArg& arg) const ;
147
148 // Configuration of MC generators used for this pdf
149 const RooNumGenConfig* getGeneratorConfig() const ;
152 RooNumGenConfig* specialGeneratorConfig(bool createOnTheFly) ;
153 void setGeneratorConfig() ;
154 void setGeneratorConfig(const RooNumGenConfig& config) ;
155
156 template <typename... CmdArgs_t>
158 {
160 }
161
162 template <typename... CmdArgs_t>
164 {
166 }
167
168 // Constraint management
169 virtual RooArgSet* getConstraints(const RooArgSet& /*observables*/, RooArgSet const& /*constrainedParams*/, RooArgSet& /*pdfParams*/) const
170 {
171 // Interface to retrieve constraint terms on this pdf. Default implementation returns null
172 return nullptr ;
173 }
174 RooArgSet* getAllConstraints(const RooArgSet& observables, RooArgSet& constrainedParams,
175 bool stripDisconnected=true) const ;
176
177 // Project p.d.f into lower dimensional p.d.f
178 virtual RooAbsPdf* createProjection(const RooArgSet& iset) ;
179
180 // Create cumulative density function from p.d.f
182 RooFit::OwningPtr<RooAbsReal> createCdf(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2={},
183 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
184 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
185 const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
186 RooFit::OwningPtr<RooAbsReal> createScanCdf(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
187
188 // Function evaluation support
189 double getValV(const RooArgSet* set=nullptr) const override ;
190 virtual double getLogVal(const RooArgSet* set=nullptr) const ;
191
192 void getLogProbabilities(std::span<const double> pdfValues, double * output) const;
193
194 /// \copydoc getNorm(const RooArgSet*) const
195 double getNorm(const RooArgSet& nset) const {
196 return getNorm(&nset) ;
197 }
198 virtual double getNorm(const RooArgSet* set=nullptr) const ;
199
200 virtual void resetErrorCounters(Int_t resetValue=10) ;
201 void setTraceCounter(Int_t value, bool allNodes=false) ;
202
203 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=nullptr) const override ;
204
205 /// Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
206 /// Always returns false, unless a PDF overrides this function.
207 virtual bool selfNormalized() const {
208 return false ;
209 }
210
211 // Support for extended maximum likelihood, switched off by default
213 /// Returns ability of PDF to provide extended likelihood terms. Possible
214 /// answers are in the enumerator RooAbsPdf::ExtendMode.
215 /// This default implementation always returns CanNotBeExtended.
216 virtual ExtendMode extendMode() const { return CanNotBeExtended; }
217 /// If true, PDF can provide extended likelihood term.
218 inline bool canBeExtended() const {
219 return (extendMode() != CanNotBeExtended) ;
220 }
221 /// If true PDF must provide extended likelihood term.
222 inline bool mustBeExtended() const {
223 return (extendMode() == MustBeExtended) ;
224 }
225 /// Return expected number of events to be used in calculation of extended
226 /// likelihood.
227 virtual double expectedEvents(const RooArgSet* nset) const ;
228 /// Return expected number of events to be used in calculation of extended
229 /// likelihood. This function should not be overridden, as it just redirects
230 /// to the actual virtual function but takes a RooArgSet reference instead of
231 /// pointer. \see expectedEvents(const RooArgSet*) const
232 double expectedEvents(const RooArgSet& nset) const {
233 return expectedEvents(&nset) ;
234 }
235
236 virtual std::unique_ptr<RooAbsReal> createExpectedEventsFunc(const RooArgSet* nset) const;
237
238 // Printing interface (human readable)
239 void printValue(std::ostream& os) const override ;
240 void printMultiline(std::ostream& os, Int_t contents, bool verbose=false, TString indent="") const override ;
241
242 static void verboseEval(Int_t stat) ;
243 static int verboseEval() ;
244
245 double extendedTerm(double sumEntries, double expected, double sumEntriesW2=0.0, bool doOffset=false) const;
246 double extendedTerm(double sumEntries, RooArgSet const* nset, double sumEntriesW2=0.0, bool doOffset=false) const;
247 double extendedTerm(RooAbsData const& data, bool weightSquared, bool doOffset=false) const;
248
249 void setNormRange(const char* rangeName) ;
250 const char* normRange() const {
251 return _normRange.Length()>0 ? _normRange.Data() : nullptr ;
252 }
253 void setNormRangeOverride(const char* rangeName) ;
254
255 const RooAbsReal* getNormIntegral(const RooArgSet& nset) const { return getNormObj(nullptr,&nset,nullptr) ; }
256
257 virtual const RooAbsReal* getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName=nullptr) const ;
258
259 virtual RooAbsGenContext* binnedGenContext(const RooArgSet &vars, bool verbose= false) const ;
260
261 virtual RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr,
262 const RooArgSet* auxProto=nullptr, bool verbose= false) const ;
263
264 virtual RooAbsGenContext* autoGenContext(const RooArgSet &vars, const RooDataSet* prototype=nullptr, const RooArgSet* auxProto=nullptr,
265 bool verbose=false, bool autoBinned=true, const char* binnedTag="") const ;
266
267 std::unique_ptr<RooAbsArg> compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext & ctx) const override;
268
269private:
270
271 std::unique_ptr<RooDataSet> generate(RooAbsGenContext& context, const RooArgSet& whatVars, const RooDataSet* prototype,
272 double nEvents, bool verbose, bool randProtoOrder, bool resampleProto, bool skipInit=false,
273 bool extended=false) const ;
274
275 // Implementation version
276 virtual RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, bool showConstants=false,
277 const char *label= "", double xmin=0.65,
278 double xmax= 0.99,double ymax=0.95, const RooCmdArg* formatCmd=nullptr) ;
279
280 void logBatchComputationErrors(std::span<const double>& outputs, std::size_t begin) const;
281 bool traceEvalPdf(double value) const;
282
283 /// Setter for the _normSet member, which should never be set directly.
284 inline void setActiveNormSet(RooArgSet const* normSet) const {
285 _normSet = normSet;
286 // Also store the unique ID of the _normSet. This makes it possible to
287 // detect if the pointer was invalidated.
289 }
290
291protected:
292
293 virtual std::unique_ptr<RooAbsReal> createNLLImpl(RooAbsData& data, const RooLinkedList& cmdList);
294 virtual std::unique_ptr<RooFitResult> fitToImpl(RooAbsData& data, const RooLinkedList& cmdList);
295
296 /// Checks if `normSet` is the currently active normalization set of this
297 /// PDF, meaning is exactly the same object as the one the `_normSet` member
298 /// points to (or both are `nullptr`).
299 inline bool isActiveNormSet(RooArgSet const* normSet) const {
300 return RooFit::getUniqueId(normSet).value() == _normSetId;
301 }
302
303 double normalizeWithNaNPacking(double rawVal, double normVal) const;
304
305 RooPlot *plotOn(RooPlot *frame, PlotOpt o) const override;
306
307 friend class RooMCStudy ;
308
309 Int_t* randomizeProtoOrder(Int_t nProto,Int_t nGen,bool resample=false) const ;
310
311 // This also forces the definition of a copy ctor in derived classes
312 RooAbsPdf(const RooAbsPdf& other, const char* name = nullptr);
313
315
316 virtual bool syncNormalization(const RooArgSet* dset, bool adjustProxies=true) const ;
317
318 mutable double _rawValue = 0;
319 mutable RooAbsReal* _norm = nullptr; //! Normalization integral (owned by _normMgr)
320 mutable RooArgSet const* _normSet = nullptr; //! Normalization set with for above integral
321
323 public:
324 CacheElem(RooAbsReal& norm) : _norm(&norm) {} ;
325 ~CacheElem() override ;
327 std::unique_ptr<RooAbsReal> _norm;
328 } ;
329 mutable RooObjCacheManager _normMgr ; //! The cache manager
330
331 bool redirectServersHook(const RooAbsCollection & newServerList, bool mustReplaceAll,
332 bool nameChange, bool isRecursiveStep) override;
333
334 mutable Int_t _errorCount = 0; ///< Number of errors remaining to print
335 mutable Int_t _traceCount = 0; ///< Number of traces remaining to print
336 mutable Int_t _negCount = 0; ///< Number of negative probabilities remaining to print
337
338 bool _selectComp = false; ///< Component selection flag for RooAbsPdf::plotCompOn
339
340 std::unique_ptr<RooNumGenConfig> _specGeneratorConfig ; ///<! MC generator configuration specific for this object
341
342 TString _normRange ; ///< Normalization range
344
345private:
346 mutable RooFit::UniqueId<RooArgSet>::Value_t _normSetId = RooFit::UniqueId<RooArgSet>::nullval; ///<! Unique ID of the currently-active normalization set
347
348 friend class RooAbsReal;
349 friend class RooChi2Var;
350
351 ClassDefOverride(RooAbsPdf,5) // Abstract PDF with normalization support
352};
353
354
355
356
357#endif
#define ClassDef(name, id)
Definition Rtypes.h:342
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
static void indent(ostringstream &buf, int indent_level)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
float ymax
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
Abstract container object that can hold multiple RooAbsArg objects.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract base class for generator contexts of RooAbsPdf objects.
Normalization set with for above integral.
Definition RooAbsPdf.h:322
std::unique_ptr< RooAbsReal > _norm
Definition RooAbsPdf.h:327
~CacheElem() override
Destructor of normalization cache element.
CacheElem(RooAbsReal &norm)
Definition RooAbsPdf.h:324
RooArgList containedArgs(Action) override
Definition RooAbsPdf.h:326
std::unique_ptr< RooAbsGenContext > _genContext
Definition RooAbsPdf.h:83
RooArgSet _whatVars
Definition RooAbsPdf.h:84
GenSpec(const GenSpec &other)
RooDataSet * _protoData
Definition RooAbsPdf.h:85
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual bool syncNormalization(const RooArgSet *dset, bool adjustProxies=true) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:195
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
RooObjCacheManager _normMgr
Definition RooAbsPdf.h:329
RooFit::UniqueId< RooArgSet >::Value_t _normSetId
! Unique ID of the currently-active normalization set
Definition RooAbsPdf.h:346
std::unique_ptr< RooNumGenConfig > _specGeneratorConfig
! MC generator configuration specific for this object
Definition RooAbsPdf.h:340
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
Definition RooAbsPdf.h:163
double getValV(const RooArgSet *set=nullptr) const override
Return current value, normalized by integrating over the observables in nset.
virtual std::unique_ptr< RooFitResult > fitToImpl(RooAbsData &data, const RooLinkedList &cmdList)
Protected implementation of the likelihood fitting routine.
bool _selectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition RooAbsPdf.h:338
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
RooFit::OwningPtr< RooAbsReal > createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
double expectedEvents(const RooArgSet &nset) const
Return expected number of events to be used in calculation of extended likelihood.
Definition RooAbsPdf.h:232
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet const &, RooArgSet &) const
Definition RooAbsPdf.h:169
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
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...
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
RooFit::OwningPtr< 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.
bool isActiveNormSet(RooArgSet const *normSet) const
Checks if normSet is the currently active normalization set of this PDF, meaning is exactly the same ...
Definition RooAbsPdf.h:299
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, bool verbose=false) const
Return a binned generator context.
TString _normRange
Normalization range.
Definition RooAbsPdf.h:342
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Int_t * randomizeProtoOrder(Int_t nProto, Int_t nGen, bool resample=false) const
Return lookup table with randomized order for nProto prototype events.
void setNormRange(const char *rangeName)
const RooAbsReal * getNormIntegral(const RooArgSet &nset) const
Definition RooAbsPdf.h:255
~RooAbsPdf() override
Destructor.
bool mustBeExtended() const
If true PDF must provide extended likelihood term.
Definition RooAbsPdf.h:222
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:320
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:124
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
virtual bool selfNormalized() const
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
Definition RooAbsPdf.h:207
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print multi line detailed information of this RooAbsPdf.
Int_t _traceCount
Number of traces remaining to print.
Definition RooAbsPdf.h:335
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:218
RooAbsReal * _norm
Definition RooAbsPdf.h:319
void setTraceCounter(Int_t value, bool allNodes=false)
Reset trace counter to given value, limiting the number of future trace messages for this pdf to 'val...
GenSpec * prepareMultiGen(const RooArgSet &whatVars, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={})
Prepare GenSpec configuration object for efficient generation of multiple datasets from identical spe...
Int_t _errorCount
Number of errors remaining to print.
Definition RooAbsPdf.h:334
@ CanBeExtended
Definition RooAbsPdf.h:212
@ MustBeExtended
Definition RooAbsPdf.h:212
@ CanNotBeExtended
Definition RooAbsPdf.h:212
double _rawValue
Definition RooAbsPdf.h:318
const char * normRange() const
Definition RooAbsPdf.h:250
virtual std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const
Returns an object that represents the expected number of events for a given normalization set,...
virtual RooPlot * paramOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={})
Add a box with parameter values (and errors) to the specified frame.
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:157
Int_t _negCount
Number of negative probabilities remaining to print.
Definition RooAbsPdf.h:336
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:57
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=nullptr) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
void setActiveNormSet(RooArgSet const *normSet) const
Setter for the _normSet member, which should never be set directly.
Definition RooAbsPdf.h:284
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
void setNormRangeOverride(const char *rangeName)
virtual RooFit::OwningPtr< RooDataSet > generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
double normalizeWithNaNPacking(double rawVal, double normVal) const
virtual RooAbsGenContext * autoGenContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false, bool autoBinned=true, const char *binnedTag="") const
RooArgSet * getAllConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, bool stripDisconnected=true) const
This helper function finds and collects all constraints terms of all component p.d....
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
Definition RooAbsPdf.h:216
RooAbsPdf()
Default constructor.
virtual RooFit::OwningPtr< RooDataHist > generateBinned(const RooArgSet &whatVars, double nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}) const
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition RooAbsPdf.h:110
bool traceEvalPdf(double value) const
Check that passed value is positive and not 'not-a-number'.
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
void printValue(std::ostream &os) const override
Print value of p.d.f, also print normalization integral that was last used, if any.
virtual std::unique_ptr< RooAbsReal > createNLLImpl(RooAbsData &data, const RooLinkedList &cmdList)
Protected implementation of the NLL creation routine.
void logBatchComputationErrors(std::span< const double > &outputs, std::size_t begin) const
Scan through outputs and fix+log all nans and negative values.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const
Interface function to create a generator context from a p.d.f.
void getLogProbabilities(std::span< const double > pdfValues, double *output) const
static TString _normRangeOverride
Definition RooAbsPdf.h:343
static Int_t _verboseEval
Definition RooAbsPdf.h:314
double extendedTerm(double sumEntries, double expected, double sumEntriesW2=0.0, bool doOffset=false) const
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
virtual double getLogVal(const RooArgSet *set=nullptr) const
Return the log of the current value with given normalization An error message is printed if the argum...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
virtual double maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const
Plot (project) PDF on specified frame.
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:24
Object to represent discrete states.
Definition RooCategory.h:28
Simple calculation from a binned dataset and a PDF.
Definition RooChi2Var.h:19
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Container class to hold unbinned data.
Definition RooDataSet.h:33
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
Helper class to facilitate Monte Carlo studies such as 'goodness-of-fit' studies, that involve fittin...
Definition RooMCStudy.h:32
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:622
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
A doubly linked list.
Definition TList.h:38
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:139
Ssiz_t Length() const
Definition TString.h:417
const char * Data() const
Definition TString.h:376
RooCmdArg NumEvents(Int_t numEvents)
std::unique_ptr< RooLinkedList > createCmdList()
UniqueId_t const & getUniqueId(Class const *ptr)
A helper function to replace pointer comparisons with UniqueId comparisons.
Definition UniqueId.h:89
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35
A UniqueId can be added as a class member to enhance any class with a unique identifier for each inst...
Definition UniqueId.h:39
unsigned long Value_t
Definition UniqueId.h:41
static void output()