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 "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 ;
39namespace RooBatchCompute {
40struct RunContext;
41}
42
43class RooAbsPdf : public RooAbsReal {
44public:
45
46 // Constructors, assignment etc
47 RooAbsPdf() ;
48 RooAbsPdf(const char *name, const char *title=0) ;
49 RooAbsPdf(const char *name, const char *title, Double_t minVal, Double_t maxVal) ;
50 // RooAbsPdf(const RooAbsPdf& other, const char* name=0);
51 virtual ~RooAbsPdf();
52
53 // Toy MC generation
54
55 ////////////////////////////////////////////////////////////////////////////////
56 /// See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&)
57 /// \param[in] nEvents How many events to generate
58 RooDataSet *generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg& arg1,
59 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
60 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) {
61 return generate(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5) ;
62 }
63 RooDataSet *generate(const RooArgSet &whatVars,
64 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
65 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
66 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
67 RooDataSet *generate(const RooArgSet &whatVars, Double_t nEvents = 0, Bool_t verbose=kFALSE, Bool_t autoBinned=kTRUE,
68 const char* binnedTag="", Bool_t expectedData=kFALSE, Bool_t extended = kFALSE) const;
69 RooDataSet *generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents= 0,
70 Bool_t verbose=kFALSE, Bool_t randProtoOrder=kFALSE, Bool_t resampleProto=kFALSE) const;
71
72
73 class GenSpec {
74 public:
75 virtual ~GenSpec() ;
77 private:
78 GenSpec(RooAbsGenContext* context, const RooArgSet& whatVars, RooDataSet* protoData, Int_t nGen, Bool_t extended,
79 Bool_t randProto, Bool_t resampleProto, TString dsetName, Bool_t init=kFALSE) ;
80 GenSpec(const GenSpec& other) ;
81
82 friend class RooAbsPdf ;
92 ClassDef(GenSpec,0) // Generation specification
93 } ;
94
95 ///Prepare GenSpec configuration object for efficient generation of multiple datasets from identical specification.
96 GenSpec* prepareMultiGen(const RooArgSet &whatVars,
97 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
98 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
99 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) ;
100 ///Generate according to GenSpec obtained from prepareMultiGen().
101 RooDataSet* generate(GenSpec&) const ;
102
103
104 ////////////////////////////////////////////////////////////////////////////////
105 /// As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,const RooCmdArg&)
106 /// \param[in] nEvents How many events to generate
107 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg& arg1,
108 const RooCmdArg& arg2=RooCmdArg::none(), const RooCmdArg& arg3=RooCmdArg::none(),
109 const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none()) const {
110 return generateBinned(whatVars,RooFit::NumEvents(nEvents),arg1,arg2,arg3,arg4,arg5);
111 }
112 virtual RooDataHist *generateBinned(const RooArgSet &whatVars,
113 const RooCmdArg& arg1=RooCmdArg::none(),const RooCmdArg& arg2=RooCmdArg::none(),
114 const RooCmdArg& arg3=RooCmdArg::none(),const RooCmdArg& arg4=RooCmdArg::none(),
115 const RooCmdArg& arg5=RooCmdArg::none(),const RooCmdArg& arg6=RooCmdArg::none()) const;
116 virtual RooDataHist *generateBinned(const RooArgSet &whatVars, Double_t nEvents, Bool_t expectedData=kFALSE, Bool_t extended=kFALSE) const;
117
118 virtual RooDataSet* generateSimGlobal(const RooArgSet& whatVars, Int_t nEvents) ;
119
120 ///Helper calling plotOn(RooPlot*, RooLinkedList&) const
121 virtual RooPlot* plotOn(RooPlot* frame,
122 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
123 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
124 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
125 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none(),
126 const RooCmdArg& arg9=RooCmdArg::none(), const RooCmdArg& arg10=RooCmdArg::none()
127 ) const {
128 return RooAbsReal::plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8,arg9,arg10) ;
129 }
130 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
131
132 /// Add a box with parameter values (and errors) to the specified frame
133 virtual RooPlot* paramOn(RooPlot* frame,
134 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
135 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
136 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
137 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
138
139 virtual RooPlot* paramOn(RooPlot* frame, const RooAbsData* data, const char *label= "", Int_t sigDigits = 2,
140 Option_t *options = "NELU", Double_t xmin=0.50,
141 Double_t xmax= 0.99,Double_t ymax=0.95) ;
142
143 // Built-in generator support
144 virtual Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
145 virtual void initGenerator(Int_t code) ;
146 virtual void generateEvent(Int_t code);
147 virtual Bool_t isDirectGenSafe(const RooAbsArg& arg) const ;
148
149 // Configuration of MC generators used for this pdf
150 const RooNumGenConfig* getGeneratorConfig() const ;
154 void setGeneratorConfig() ;
155 void setGeneratorConfig(const RooNumGenConfig& config) ;
156
157 // -log(L) fits to binned and unbinned data
158 virtual RooFitResult* fitTo(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
159 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
160 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
161 virtual RooFitResult* fitTo(RooAbsData& data, const RooLinkedList& cmdList) ;
162
163 virtual RooAbsReal* createNLL(RooAbsData& data, const RooLinkedList& cmdList) ;
164 virtual RooAbsReal* createNLL(RooAbsData& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
165 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
166 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
167
168 // Chi^2 fits to histograms
171 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
172 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
173 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
174 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
175
176 // Chi^2 fits to X-Y datasets
177 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
178
179
180
181
182
183 // Constraint management
184 virtual RooArgSet* getConstraints(const RooArgSet& /*observables*/, RooArgSet& /*constrainedParams*/, Bool_t /*stripDisconnected*/) const {
185 // Interface to retrieve constraint terms on this pdf. Default implementation returns null
186 return 0 ;
187 }
188 virtual RooArgSet* getAllConstraints(const RooArgSet& observables, RooArgSet& constrainedParams, Bool_t stripDisconnected=kTRUE) const ;
189
190 // Project p.d.f into lower dimensional p.d.f
191 virtual RooAbsPdf* createProjection(const RooArgSet& iset) ;
192
193 // Create cumulative density function from p.d.f
194 RooAbsReal* createCdf(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
195 RooAbsReal* createCdf(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
196 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
197 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
198 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
199 RooAbsReal* createScanCdf(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
200
201 // Function evaluation support
202 virtual Double_t getValV(const RooArgSet* set=0) const ;
203 virtual Double_t getLogVal(const RooArgSet* set=0) const ;
204
206 RooSpan<const double> getLogValBatch(std::size_t begin, std::size_t batchSize,
207 const RooArgSet* normSet = nullptr) const;
208 RooSpan<const double> getLogProbabilities(RooBatchCompute::RunContext& evalData, const RooArgSet* normSet = nullptr) const;
209
210 /// \copydoc getNorm(const RooArgSet*) const
211 Double_t getNorm(const RooArgSet& nset) const {
212 return getNorm(&nset) ;
213 }
214 virtual Double_t getNorm(const RooArgSet* set=0) const ;
215
216 virtual void resetErrorCounters(Int_t resetValue=10) ;
217 void setTraceCounter(Int_t value, Bool_t allNodes=kFALSE) ;
218private:
219 Bool_t traceEvalPdf(Double_t value) const;
220
221public:
222
223 Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
224
225 /// Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
226 /// Always returns false, unless a PDF overrides this function.
227 virtual Bool_t selfNormalized() const {
228 return kFALSE ;
229 }
230
231 // Support for extended maximum likelihood, switched off by default
233 /// Returns ability of PDF to provide extended likelihood terms. Possible
234 /// answers are in the enumerator RooAbsPdf::ExtendMode.
235 /// This default implementation always returns CanNotBeExtended.
236 virtual ExtendMode extendMode() const { return CanNotBeExtended; }
237 /// If true, PDF can provide extended likelihood term.
238 inline Bool_t canBeExtended() const {
239 return (extendMode() != CanNotBeExtended) ;
240 }
241 /// If true PDF must provide extended likelihood term.
242 inline Bool_t mustBeExtended() const {
243 return (extendMode() == MustBeExtended) ;
244 }
245 virtual Double_t expectedEvents(const RooArgSet* nset) const ;
246 /// Return expected number of events to be used in calculation of extended likelihood.
247 virtual Double_t expectedEvents(const RooArgSet& nset) const {
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 = nullptr; //! Normalization integral (owned by _normMgr)
322 mutable RooArgSet const* _normSet = nullptr; //! 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 = nullptr ;
343
344 // Similar to the situation with the normalization integral above: if a
345 // server is redirected, the cached normalization set might not point to
346 // the right observables anymore. We need to reset it.
347 _normSet = nullptr ;
348 return false ;
349 } ;
350
351
352 mutable Int_t _errorCount ; // Number of errors remaining to print
353 mutable Int_t _traceCount ; // Number of traces remaining to print
354 mutable Int_t _negCount ; // Number of negative probablities remaining to print
355
356 Bool_t _selectComp ; // Component selection flag for RooAbsPdf::plotCompOn
357
358 RooNumGenConfig* _specGeneratorConfig ; //! MC generator configuration specific for this object
359
360 TString _normRange ; // Normalization range
362
363private:
364 template<class Minimizer>
365 int calculateAsymptoticCorrectedCovMatrix(Minimizer& minimizer, RooAbsData const& data);
366
367 template<class Minimizer>
368 int calculateSumW2CorrectedCovMatrix(Minimizer& minimizer, RooAbsReal const& nll) const;
369
370 ClassDef(RooAbsPdf,4) // Abstract PDF with normalization support
371};
372
373
374#endif
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:91
const char Option_t
Definition RtypesCore.h:66
#define ClassDef(name, id)
Definition Rtypes.h:325
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
float ymax
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:606
virtual void operModeHook()
Definition RooAbsArg.h:567
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:49
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.
virtual RooArgList containedArgs(Action)
Definition RooAbsPdf.h:329
CacheElem(RooAbsReal &norm)
Definition RooAbsPdf.h:326
RooAbsReal * _norm
Definition RooAbsPdf.h:330
RooArgSet _whatVars
Definition RooAbsPdf.h:84
RooAbsGenContext * _genContext
Definition RooAbsPdf.h:83
Bool_t _resampleProto
Definition RooAbsPdf.h:89
GenSpec(const GenSpec &other)
RooDataSet * _protoData
Definition RooAbsPdf.h:85
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...
RooSpan< const double > getLogValBatch(std::size_t begin, std::size_t batchSize, const RooArgSet *normSet=nullptr) const
RooObjCacheManager _normMgr
Definition RooAbsPdf.h:332
Double_t getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:211
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
virtual ~RooAbsPdf()
Destructor.
void logBatchComputationErrors(RooSpan< const double > &outputs, std::size_t begin) const
Scan through outputs and fix+log all nans and negative values.
RooSpan< const double > getLogProbabilities(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
Compute the log-likelihoods for all events in the requested batch.
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...
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...
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.
virtual RooFitResult * chi2FitTo(RooDataHist &data, const RooLinkedList &cmdList)
Calls RooAbsPdf::createChi2(RooDataSet& data, const RooLinkedList& cmdList) and returns fit result.
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Bool_t traceEvalPdf(Double_t value) const
Check that passed value is positive and not 'not-a-number'.
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
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.
RooAbsReal * createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
TString _normRange
MC generator configuration specific for this object.
Definition RooAbsPdf.h:360
void setNormRange(const char *rangeName)
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:107
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:322
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...
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Bool_t canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:238
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:227
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print multi line detailed information of this RooAbsPdf.
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
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:58
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
Int_t _traceCount
Definition RooAbsPdf.h:353
virtual Double_t expectedEvents(const RooArgSet &nset) const
Return expected number of events to be used in calculation of extended likelihood.
Definition RooAbsPdf.h:247
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.
virtual void printValue(std::ostream &os) const
Print value of p.d.f, also print normalization integral that was last used, if any.
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....
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.
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...
int calculateAsymptoticCorrectedCovMatrix(Minimizer &minimizer, RooAbsData const &data)
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Int_t _errorCount
Definition RooAbsPdf.h:352
Bool_t _selectComp
Definition RooAbsPdf.h:356
virtual RooArgSet * getConstraints(const RooArgSet &, RooArgSet &, Bool_t) const
Definition RooAbsPdf.h:184
@ CanBeExtended
Definition RooAbsPdf.h:232
@ MustBeExtended
Definition RooAbsPdf.h:232
@ CanNotBeExtended
Definition RooAbsPdf.h:232
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.
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...
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.
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 ...
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...
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...
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Int_t _negCount
Definition RooAbsPdf.h:354
Bool_t mustBeExtended() const
If true PDF must provide extended likelihood term.
Definition RooAbsPdf.h:242
void setNormRangeOverride(const char *rangeName)
virtual RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
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:236
int calculateSumW2CorrectedCovMatrix(Minimizer &minimizer, RooAbsReal const &nll) const
RooAbsPdf()
Default constructor.
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, Bool_t verbose=kFALSE) const
Return a binned generator context.
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
RooNumGenConfig * _specGeneratorConfig
Definition RooAbsPdf.h:358
friend class RooAddGenContextOrig
Definition RooAbsPdf.h:304
Double_t _rawValue
Definition RooAbsPdf.h:320
static TString _normRangeOverride
Definition RooAbsPdf.h:361
static Int_t _verboseEval
Definition RooAbsPdf.h:315
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Function that is called at the end of redirectServers().
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.
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:121
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:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
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:37
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
RooEffGenContext is a specialized generator context for p.d.fs represented by class RooEffProd,...
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
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: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)
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