Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsReal.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $Id: RooAbsReal.h,v 1.75 2007/07/13 21:50:24 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_REAL
17#define ROO_ABS_REAL
18
19#include "RooAbsArg.h"
20#include "RooCmdArg.h"
21#include "RooCurve.h"
22#include "RooArgSet.h"
23#include "RooArgList.h"
24#include "RooGlobalFunc.h"
25#include "RooSpan.h"
28
29class RooArgList ;
30class RooDataSet ;
31class RooPlot;
32class RooRealVar;
33class RooAbsFunc;
35class RooCategory ;
36class RooLinkedList ;
37class RooNumIntConfig ;
38class RooDataHist ;
39class RooFunctor ;
40class RooGenFunction ;
42class RooFitResult ;
43class RooAbsMoment ;
44class RooDerivative ;
46struct TreeReadBuffer; /// A space to attach TBranches
47namespace ROOT {
48namespace Experimental {
49class RooFitDriver ;
50}
51}
52
53class TH1;
54class TH1F;
55class TH2F;
56class TH3F;
57
58#include <iostream>
59#include <list>
60#include <map>
61#include <string>
62#include <sstream>
63
64class RooAbsReal : public RooAbsArg {
65public:
67
68 // Constructors, assignment etc
69 RooAbsReal() ;
70 RooAbsReal(const char *name, const char *title, const char *unit= "") ;
71 RooAbsReal(const char *name, const char *title, Double_t minVal, Double_t maxVal,
72 const char *unit= "") ;
73 RooAbsReal(const RooAbsReal& other, const char* name=0);
74 RooAbsReal& operator=(const RooAbsReal& other);
75 virtual ~RooAbsReal();
76
77
78
79
80 //////////////////////////////////////////////////////////////////////////////////
81 /// Evaluate object. Returns either cached value or triggers a recalculation.
82 /// The recalculation happens by calling getValV(), which in the end calls the
83 /// virtual evaluate() functions of the respective PDFs.
84 /// \param[in] normalisationSet getValV() reacts differently depending on the value of the normalisation set.
85 /// If the set is `nullptr`, an unnormalised value is returned.
86 /// \note The normalisation is arbitrary, because it is up to the implementation
87 /// of the PDF to e.g. leave out normalisation constants for speed reasons. The range
88 /// of the variables is also ignored.
89 ///
90 /// To normalise the result properly, a RooArgSet has to be passed, which contains
91 /// the variables to normalise over.
92 /// These are integrated over their current ranges to compute the normalisation constant,
93 /// and the unnormalised result is divided by this value.
94 inline Double_t getVal(const RooArgSet* normalisationSet = nullptr) const {
95 // Sometimes, the calling code uses an empty RooArgSet to request evaluation
96 // without normalization set instead of following the `nullptr` convention.
97 // To remove this ambiguity which might not always be correctly handled in
98 // downstream code, we set `normalisationSet` to nullptr if it is pointing
99 // to an empty set.
100 if(normalisationSet && normalisationSet->empty()) {
101 normalisationSet = nullptr;
102 }
103#ifdef ROOFIT_CHECK_CACHED_VALUES
104 return _DEBUG_getVal(normalisationSet);
105#else
106
107#ifndef _WIN32
108 return (_fast && !_inhibitDirty) ? _value : getValV(normalisationSet) ;
109#else
110 return (_fast && !inhibitDirty()) ? _value : getValV(normalisationSet) ;
111#endif
112
113#endif
114 }
115
116 /// Like getVal(const RooArgSet*), but always requires an argument for normalisation.
117 inline Double_t getVal(const RooArgSet& normalisationSet) const {
118 // Sometimes, the calling code uses an empty RooArgSet to request evaluation
119 // without normalization set instead of following the `nullptr` convention.
120 // To remove this ambiguity which might not always be correctly handled in
121 // downstream code, we set `normalisationSet` to nullptr if it is an empty set.
122 return _fast ? _value : getValV(normalisationSet.empty() ? nullptr : &normalisationSet) ;
123 }
124
125 virtual Double_t getValV(const RooArgSet* normalisationSet = nullptr) const ;
126
127 /// \deprecated getValBatch() has been removed in favour of the faster getValues(). If your code is affected
128 /// by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this transition.
129 /// https://root.cern/doc/v624/release-notes.html
130#ifndef R__MACOSX
131 virtual RooSpan<const double> getValBatch(std::size_t /*begin*/, std::size_t /*maxSize*/, const RooArgSet* /*normSet*/ = nullptr) = delete;
132#else
133 //AppleClang in MacOS10.14 has a linker bug and fails to link programs that create objects of classes containing virtual deleted methods.
134 //This can be safely deleted when MacOS10.14 is no longer supported by ROOT. See https://reviews.llvm.org/D37830
135 virtual RooSpan<const double> getValBatch(std::size_t /*begin*/, std::size_t /*maxSize*/, const RooArgSet* /*normSet*/ = nullptr) final {
136 throw std::logic_error("Deprecated getValBatch() has been removed in favour of the faster getValues(). If your code is affected by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this transition. https://root.cern/doc/v624/release-notes.html");
137 }
138#endif
139 /// \copydoc getValBatch(std::size_t, std::size_t, const RooArgSet*)
140 virtual RooSpan<const double> getValues(RooBatchCompute::RunContext& evalData, const RooArgSet* normSet = nullptr) const;
141 std::vector<double> getValues(RooAbsData const& data, RooFit::BatchModeOption batchMode=RooFit::BatchModeOption::Cpu) const;
142
143 Double_t getPropagatedError(const RooFitResult &fr, const RooArgSet &nset = RooArgSet()) const;
144
145 Bool_t operator==(Double_t value) const ;
146 virtual Bool_t operator==(const RooAbsArg& other) const;
147 virtual Bool_t isIdentical(const RooAbsArg& other, Bool_t assumeSameType=kFALSE) const;
148
149
150 inline const Text_t *getUnit() const {
151 // Return string with unit description
152 return _unit.Data();
153 }
154 inline void setUnit(const char *unit) {
155 // Set unit description to given string
156 _unit= unit;
157 }
158 TString getTitle(Bool_t appendUnit= kFALSE) const;
159
160 // Lightweight interface adaptors (caller takes ownership)
161 RooAbsFunc *bindVars(const RooArgSet &vars, const RooArgSet* nset=0, Bool_t clipInvalid=kFALSE) const;
162
163 // Create a fundamental-type object that can hold our value.
164 RooAbsArg *createFundamental(const char* newname=0) const;
165
166 // Analytical integration support
167 virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=0) const ;
168 virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
169 virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
170 virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
171 virtual Bool_t forceAnalyticalInt(const RooAbsArg& /*dep*/) const {
172 // Interface to force RooRealIntegral to offer given observable for internal integration
173 // even if this is deemed unsafe. This default implementation returns always flase
174 return kFALSE ;
175 }
176 virtual void forceNumInt(Bool_t flag=kTRUE) {
177 // If flag is true, all advertised analytical integrals will be ignored
178 // and all integrals are calculated numerically
179 _forceNumInt = flag ;
180 }
181 Bool_t getForceNumInt() const { return _forceNumInt ; }
182
183 // Chi^2 fits to histograms
184 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
185 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
186 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
187 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
188
189 virtual RooAbsReal* createChi2(RooDataHist& data, const RooLinkedList& cmdList) ;
190 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
191 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
192 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
193
194 // Chi^2 fits to X-Y datasets
195 virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
196 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
197 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
198 virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;
199
200 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
201 virtual RooAbsReal* createChi2(RooDataSet& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
202 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
203 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
204
205
206 virtual RooAbsReal* createProfile(const RooArgSet& paramsOfInterest) ;
207
208
209 RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
210 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
211 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
212 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
213
214 /// Create integral over observables in iset in range named rangeName.
215 RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const {
216 return createIntegral(iset,0,0,rangeName) ;
217 }
218 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
219 RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=0) const {
220 return createIntegral(iset,&nset,0,rangeName) ;
221 }
222 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
223 /// using specified configuration for any numeric integration.
224 RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
225 return createIntegral(iset,&nset,&cfg,rangeName) ;
226 }
227 /// Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration.
228 RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
229 return createIntegral(iset,0,&cfg,rangeName) ;
230 }
231 virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset=0, const RooNumIntConfig* cfg=0, const char* rangeName=0) const ;
232
233
234 void setParameterizeIntegral(const RooArgSet& paramVars) ;
235
236 // Create running integrals
237 RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
238 RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
239 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
240 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
241 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
242 RooAbsReal* createIntRI(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
243 RooAbsReal* createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
244
245
246 // Optimized accept/reject generator support
247 virtual Int_t getMaxVal(const RooArgSet& vars) const ;
248 virtual Double_t maxVal(Int_t code) const ;
249 virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }
250
251
252 // Plotting options
253 void setPlotLabel(const char *label);
254 const char *getPlotLabel() const;
255
256 virtual Double_t defaultErrorLevel() const {
257 // Return default level for MINUIT error analysis
258 return 1.0 ;
259 }
260
261 const RooNumIntConfig* getIntegratorConfig() const ;
266 void setIntegratorConfig() ;
267 void setIntegratorConfig(const RooNumIntConfig& config) ;
268
269 virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),Bool_t force=kTRUE) ;
270 virtual void fixAddCoefRange(const char* rangeName=0,Bool_t force=kTRUE) ;
271
272 virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;
273
274 // User entry point for plotting
275 virtual RooPlot* plotOn(RooPlot* frame,
276 const RooCmdArg& arg1=RooCmdArg(), const RooCmdArg& arg2=RooCmdArg(),
277 const RooCmdArg& arg3=RooCmdArg(), const RooCmdArg& arg4=RooCmdArg(),
278 const RooCmdArg& arg5=RooCmdArg(), const RooCmdArg& arg6=RooCmdArg(),
279 const RooCmdArg& arg7=RooCmdArg(), const RooCmdArg& arg8=RooCmdArg(),
280 const RooCmdArg& arg9=RooCmdArg(), const RooCmdArg& arg10=RooCmdArg()
281 ) const ;
282
283
285
286 // Forwarder function for backward compatibility
287 virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L",
288 Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=0) const;
289
290 // Fill an existing histogram
291 TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
292 Double_t scaleFactor= 1, const RooArgSet *projectedVars= 0, Bool_t scaling=kTRUE,
293 const RooArgSet* condObs=0, Bool_t setError=kTRUE) const;
294
295 // Create 1,2, and 3D histograms from and fill it
296 TH1 *createHistogram(const char* varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
297 TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
298 TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
299 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
300 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
301 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
302 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
303
304 // Fill a RooDataHist
305 RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, Double_t scaleFactor,
306 Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const ;
307
308 // I/O streaming interface (machine readable)
309 virtual Bool_t readFromStream(std::istream& is, Bool_t compact, Bool_t verbose=kFALSE) ;
310 virtual void writeToStream(std::ostream& os, Bool_t compact) const ;
311
312 // Printing interface (human readable)
313 virtual void printValue(std::ostream& os) const ;
314 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
315
316 inline void setCachedValue(double value, bool notifyClients = true) final;
317
318 // Evaluation error logging
319 class EvalError {
320 public:
322 EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
323 void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
324 void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
325 std::string _msg;
326 std::string _srvval;
327 } ;
328
332 void logEvalError(const char* message, const char* serverValueString=0) const ;
333 static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=0) ;
334 static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
335 static Int_t numEvalErrors() ;
336 static Int_t numEvalErrorItems() ;
337
338
339 typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ;
341
342 static void clearEvalErrorLog() ;
343
344 /// Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
345 virtual Bool_t isBinnedDistribution(const RooArgSet& /*obs*/) const { return kFALSE ; }
346 virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const;
347 virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const;
348
350 RooMultiGenFunction* iGenFunction(const RooArgSet& observables, const RooArgSet& nset=RooArgSet()) ;
351
352 RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
353 TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
354
355 RooDerivative* derivative(RooRealVar& obs, Int_t order=1, Double_t eps=0.001) ;
356 RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, Double_t eps=0.001) ;
357
358 RooAbsMoment* moment(RooRealVar& obs, Int_t order, Bool_t central, Bool_t takeRoot) ;
359 RooAbsMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, Bool_t central, Bool_t takeRoot, Bool_t intNormObs) ;
360
361 RooAbsMoment* mean(RooRealVar& obs) { return moment(obs,1,kFALSE,kFALSE) ; }
362 RooAbsMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,kFALSE,kFALSE,kTRUE) ; }
363 RooAbsMoment* sigma(RooRealVar& obs) { return moment(obs,2,kTRUE,kTRUE) ; }
364 RooAbsMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,kTRUE,kTRUE,kTRUE) ; }
365
367
368
369 virtual Bool_t setData(RooAbsData& /*data*/, Bool_t /*cloneData*/=kTRUE) { return kTRUE ; }
370
371 virtual void enableOffsetting(Bool_t) {} ;
372 virtual Bool_t isOffsetting() const { return kFALSE ; }
373 virtual Double_t offset() const { return 0 ; }
374
375 static void setHideOffset(Bool_t flag);
376 static Bool_t hideOffset() ;
377
378protected:
379 // Hook for objects with normalization-dependent parameters interperetation
380 virtual void selectNormalization(const RooArgSet* depSet=0, Bool_t force=kFALSE) ;
381 virtual void selectNormalizationRange(const char* rangeName=0, Bool_t force=kFALSE) ;
382
383 // Helper functions for plotting
384 Bool_t plotSanityChecks(RooPlot* frame) const ;
385 void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
386 RooArgSet& projectedVars, Bool_t silent) const ;
387
388 TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=0, const char* rangeName=0, Bool_t omitEmpty=kFALSE) const ;
389
390
391 Bool_t isSelectedComp() const ;
392
393
394 public:
395 const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
396 const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
397 RooArgSet *&cloneSet, const char* rangeName=0, const RooArgSet* condObs=0) const;
398 virtual void computeBatch(cudaStream_t*, double* output, size_t size, RooFit::Detail::DataMap const&) const;
399
400 protected:
401
403
404 void plotOnCompSelect(RooArgSet* selNodes) const ;
405 RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, Double_t Z, const RooArgSet* params, const RooLinkedList& argList, Bool_t method1) const ;
406
407 // Support interface for subclasses to advertise their analytic integration
408 // and generator capabilities in their analyticalIntegral() and generateEvent()
409 // implementations.
410 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
411 const RooArgProxy& a) const ;
412 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
413 const RooArgProxy& a, const RooArgProxy& b) const ;
414 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
415 const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
416 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
417 const RooArgProxy& a, const RooArgProxy& b,
418 const RooArgProxy& c, const RooArgProxy& d) const ;
419
420 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
421 const RooArgSet& set) const ;
422
423
424 RooAbsReal* createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
425 void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;
426
427
428 // Internal consistency checking (needed by RooDataSet)
429 /// Check if current value is valid.
430 virtual bool isValid() const { return isValidReal(_value); }
431 /// Interface function to check if given value is a valid value for this object. Returns true unless overridden.
432 virtual bool isValidReal(double /*value*/, bool printError = false) const { (void)printError; return true; }
433
434
435 // Function evaluation and error tracing
436 Double_t traceEval(const RooArgSet* set) const ;
437
438 /// Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
439 virtual Double_t evaluate() const = 0;
440
441 /// \deprecated evaluateBatch() has been removed in favour of the faster evaluateSpan(). If your code is affected
442 /// by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this transition.
443 /// https://root.cern/doc/v624/release-notes.html
444#ifndef R__MACOSX
445 virtual RooSpan<double> evaluateBatch(std::size_t /*begin*/, std::size_t /*maxSize*/) = delete;
446#else
447 //AppleClang in MacOS10.14 has a linker bug and fails to link programs that create objects of classes containing virtual deleted methods.
448 //This can be safely deleted when MacOS10.14 is no longer supported by ROOT. See https://reviews.llvm.org/D37830
449 virtual RooSpan<double> evaluateBatch(std::size_t /*begin*/, std::size_t /*maxSize*/) final {
450 throw std::logic_error("Deprecated evaluatedBatch() has been removed in favour of the faster evaluateSpan(). If your code is affected by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this transition. https://root.cern/doc/v624/release-notes.html");
451 }
452#endif
453
454 virtual RooSpan<double> evaluateSpan(RooBatchCompute::RunContext& evalData, const RooArgSet* normSet) const;
455
456
457 //---------- Interface to access batch data ---------------------------
458 //
460
461 private:
462 void checkBatchComputation(const RooBatchCompute::RunContext& evalData, std::size_t evtNo, const RooArgSet* normSet = nullptr, double relAccuracy = 1.E-13) const;
463
464 /// Debug version of getVal(), which is slow and does error checking.
465 Double_t _DEBUG_getVal(const RooArgSet* normalisationSet) const;
466
467 //--------------------------------------------------------------------
468
469 protected:
470 // Hooks for RooDataSet interface
471 friend class RooRealIntegral ;
472 friend class RooVectorDataStore ;
473 virtual void syncCache(const RooArgSet* set=0) { getVal(set) ; }
474 virtual void copyCache(const RooAbsArg* source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE) ;
475 virtual void attachToTree(TTree& t, Int_t bufSize=32000) ;
476 virtual void attachToVStore(RooVectorDataStore& vstore) ;
477 virtual void setTreeBranchStatus(TTree& t, Bool_t active) ;
478 virtual void fillTreeBranch(TTree& t) ;
479
480 friend class RooRealBinding ;
481 Double_t _plotMin ; // Minimum of plot range
482 Double_t _plotMax ; // Maximum of plot range
483 Int_t _plotBins ; // Number of plot bins
484 mutable Double_t _value ; // Cache for current value of object
485 TString _unit ; // Unit for objects value
486 TString _label ; // Plot label for objects value
487 Bool_t _forceNumInt ; // Force numerical integration if flag set
488
489 friend class RooAbsPdf ;
490 friend class RooAbsAnaConvPdf ;
491
492 RooNumIntConfig* _specIntegratorConfig ; // Numeric integrator configuration specific for this object
493
494 friend class RooDataProjBinding ;
496
497 struct PlotOpt {
501 numCPU(1),interleave(RooFit::Interleave),curveNameSuffix(""), numee(10), eeval(0), doeeval(kFALSE), progress(kFALSE), errorFR(0) {} ;
511 const char* normRangeName ;
518 const char* curveName ;
519 const char* addToCurveName ;
524 const char* curveNameSuffix ;
530 } ;
531
532 // Plot implementation functions
533 virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;
534
535public:
536 // PlotOn with command list
537 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
538
539 protected:
540 virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;
541
542
543private:
544
546 static std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > > _evalErrorList ;
548
549 Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;
550
551 std::unique_ptr<TreeReadBuffer> _treeReadBuffer; //! A buffer for reading values from trees
552
553protected:
554
555
556 friend class RooRealSumPdf ;
557 friend class RooRealSumFunc;
558 friend class RooAddPdf ;
559 friend class RooAddModel ;
560 void selectComp(Bool_t flag) {
561 // If flag is true, only selected component will be included in evaluates of RooAddPdf components
562 _selectComp = flag ;
563 }
564 static void globalSelectComp(Bool_t flag) ;
565 Bool_t _selectComp ; //! Component selection flag for RooAbsPdf::plotCompOn
566 static Bool_t _globalSelectComp ; // Global activation switch for component selection
567 // This struct can be used to flip the global switch to select components.
568 // Doing this with RAII prevents forgetting to reset the state.
574 }
575
579 }
580
582 };
583
584
585 mutable RooArgSet* _lastNSet ; //!
586 static Bool_t _hideOffset ; // Offset hiding flag
587
588 ClassDef(RooAbsReal,2) // Abstract real-valued variable
589};
590
591
592/// Helper class to access a batch-related part of RooAbsReal's interface, which should not leak to the outside world.
594 public:
595 static void checkBatchComputation(const RooAbsReal& theReal, const RooBatchCompute::RunContext& evalData, std::size_t evtNo,
596 const RooArgSet* normSet = nullptr, double relAccuracy = 1.E-13) {
597 theReal.checkBatchComputation(evalData, evtNo, normSet, relAccuracy);
598 }
599};
600
601
602////////////////////////////////////////////////////////////////////////////////
603/// Overwrite the value stored in this object's cache.
604/// This can be used to fake a computation that resulted in `value`.
605/// \param[in] value Value to write.
606/// \param[in] setValDirty If true, notify users of this object that its value changed.
607/// This is the default.
608void RooAbsReal::setCachedValue(double value, bool notifyClients) {
609 _value = value;
610
611 if (notifyClients) {
613 _valueDirty = false;
614 }
615}
616
617
618#endif
double
typedef void(GLAPIENTRYP _GLUfuncptr)(void)
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
char Text_t
Definition RtypesCore.h:62
const Bool_t kFALSE
Definition RtypesCore.h:101
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)
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Helper class to access a batch-related part of RooAbsReal's interface, which should not leak to the o...
Definition RooAbsReal.h:593
static void checkBatchComputation(const RooAbsReal &theReal, const RooBatchCompute::RunContext &evalData, std::size_t evtNo, const RooArgSet *normSet=nullptr, double relAccuracy=1.E-13)
Definition RooAbsReal.h:595
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:69
Bool_t _fast
Definition RooAbsArg.h:735
friend class RooArgSet
Definition RooAbsArg.h:642
Bool_t inhibitDirty() const
Delete watch flag.
Bool_t _valueDirty
Definition RooAbsArg.h:731
static Bool_t _inhibitDirty
Definition RooAbsArg.h:714
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:505
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:82
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition RooAbsFunc.h:27
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
void setServerValues(const char *tmp)
Definition RooAbsReal.h:324
void setMessage(const char *tmp)
Definition RooAbsReal.h:323
EvalError(const EvalError &other)
Definition RooAbsReal.h:322
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
virtual RooSpan< const double > getValBatch(std::size_t, std::size_t, const RooArgSet *=nullptr)=delete
void plotOnCompSelect(RooArgSet *selNodes) const
Helper function for plotting of composite p.d.fs.
std::unique_ptr< TreeReadBuffer > _treeReadBuffer
Definition RooAbsReal.h:551
RooGenFunction * iGenFunction(RooRealVar &x, const RooArgSet &nset=RooArgSet())
virtual RooPlot * plotSliceOn(RooPlot *frame, const RooArgSet &sliceSet, Option_t *drawOptions="L", Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData *projData=0) const
virtual Bool_t readFromStream(std::istream &is, Bool_t compact, Bool_t verbose=kFALSE)
Read object contents from stream (dummy for now)
TString _label
Definition RooAbsReal.h:486
RooAbsReal * createIntegral(const RooArgSet &iset, const RooNumIntConfig &cfg, const char *rangeName=0) const
Create integral over observables in iset in range named rangeName using specified configuration for a...
Definition RooAbsReal.h:228
RooAbsMoment * mean(RooRealVar &obs)
Definition RooAbsReal.h:361
TH1 * fillHistogram(TH1 *hist, const RooArgList &plotVars, Double_t scaleFactor=1, const RooArgSet *projectedVars=0, Bool_t scaling=kTRUE, const RooArgSet *condObs=0, Bool_t setError=kTRUE) const
Fill the ROOT histogram 'hist' with values sampled from this function at the bin centers.
virtual void selectNormalizationRange(const char *rangeName=0, Bool_t force=kFALSE)
Interface function to force use of a given normalization range to interpret function value.
TString getTitle(Bool_t appendUnit=kFALSE) const
Return this variable's title string.
virtual Bool_t isOffsetting() const
Definition RooAbsReal.h:372
static Int_t numEvalErrorItems()
static void setHideOffset(Bool_t flag)
RooAbsReal()
coverity[UNINIT_CTOR] Default constructor
virtual Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
Bool_t _forceNumInt
Definition RooAbsReal.h:487
virtual RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Evaluate this object for a batch/span of data points.
RooAbsReal * createScanRI(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
Utility function for createRunningIntegral that construct an object implementing the numeric scanning...
virtual Double_t defaultErrorLevel() const
Definition RooAbsReal.h:256
virtual Bool_t isIdentical(const RooAbsArg &other, Bool_t assumeSameType=kFALSE) const
Bool_t isSelectedComp() const
If true, the current pdf is a selected component (for use in plotting)
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, 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
Create an object that represents the integral of the function over one or more observables listed in ...
virtual RooAbsReal * createChi2(RooDataHist &data, const RooLinkedList &cmdList)
Create a variable from a histogram and this function.
RooArgSet * _lastNSet
Definition RooAbsReal.h:585
RooDerivative * derivative(RooRealVar &obs, Int_t order=1, Double_t eps=0.001)
Return function representing first, second or third order derivative of this function.
double value_type
Definition RooAbsReal.h:66
void setParameterizeIntegral(const RooArgSet &paramVars)
static Bool_t hideOffset()
virtual Double_t evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
virtual void setTreeBranchStatus(TTree &t, Bool_t active)
(De)Activate associated tree branch
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
RooAbsReal * createIntegral(const RooArgSet &iset, const RooArgSet &nset, const RooNumIntConfig &cfg, const char *rangeName=0) const
Create integral over observables in iset in range named rangeName with integrand normalized over obse...
Definition RooAbsReal.h:224
TString integralNameSuffix(const RooArgSet &iset, const RooArgSet *nset=0, const char *rangeName=0, Bool_t omitEmpty=kFALSE) const
Construct string with unique suffix name to give to integral object that encodes integrated observabl...
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
Double_t _DEBUG_getVal(const RooArgSet *normalisationSet) const
Debug version of getVal(), which is slow and does error checking.
Double_t traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
Double_t getVal(const RooArgSet &normalisationSet) const
Like getVal(const RooArgSet*), but always requires an argument for normalisation.
Definition RooAbsReal.h:117
virtual void copyCache(const RooAbsArg *source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE)
Copy the cached value of another RooAbsArg to our cache.
virtual void fixAddCoefRange(const char *rangeName=0, Bool_t force=kTRUE)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
TH1 * createHistogram(const char *varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const
Create and fill a ROOT histogram TH1, TH2 or TH3 with the values of this function for the variables w...
static EvalErrorIter evalErrorIter()
RooAbsReal & operator=(const RooAbsReal &other)
Assign values, name and configs from another RooAbsReal.
Double_t findRoot(RooRealVar &x, Double_t xmin, Double_t xmax, Double_t yval)
Return value of x (in range xmin,xmax) at which function equals yval.
RooAbsReal * createIntRI(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Utility function for createRunningIntegral.
TF1 * asTF(const RooArgList &obs, const RooArgList &pars=RooArgList(), const RooArgSet &nset=RooArgSet()) const
Return a ROOT TF1,2,3 object bound to this RooAbsReal with given definition of observables and parame...
virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
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 void attachToVStore(RooVectorDataStore &vstore)
TString _unit
Definition RooAbsReal.h:485
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
RooFitResult * chi2FitDriver(RooAbsReal &fcn, RooLinkedList &cmdList)
Internal driver function for chi2 fits.
friend class RooAbsOptGoodnessOfFit
Definition RooAbsReal.h:495
Double_t _value
Definition RooAbsReal.h:484
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
virtual void attachToTree(TTree &t, Int_t bufSize=32000)
Attach object to a branch of given TTree.
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooNumIntConfig * _specIntegratorConfig
Definition RooAbsReal.h:492
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
virtual RooSpan< double > evaluateBatch(std::size_t, std::size_t)=delete
static ErrorLoggingMode evalErrorLoggingMode()
Return current evaluation error logging mode.
virtual Bool_t forceAnalyticalInt(const RooAbsArg &) const
Definition RooAbsReal.h:171
virtual Int_t minTrialSamples(const RooArgSet &) const
Definition RooAbsReal.h:249
virtual bool isValidReal(double, bool printError=false) const
Interface function to check if given value is a valid value for this object. Returns true unless over...
Definition RooAbsReal.h:432
void setIntegratorConfig()
Remove the specialized numeric integration configuration associated with this object.
virtual void selectNormalization(const RooArgSet *depSet=0, Bool_t force=kFALSE)
Interface function to force use of a given set of observables to interpret function value.
void setUnit(const char *unit)
Definition RooAbsReal.h:154
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
virtual void syncCache(const RooArgSet *set=0)
Definition RooAbsReal.h:473
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
RooAbsMoment * moment(RooRealVar &obs, Int_t order, Bool_t central, Bool_t takeRoot)
Return function representing moment of function of given order.
static ErrorLoggingMode _evalErrorMode
Definition RooAbsReal.h:545
virtual void enableOffsetting(Bool_t)
Definition RooAbsReal.h:371
Double_t getPropagatedError(const RooFitResult &fr, const RooArgSet &nset=RooArgSet()) const
Calculate error on self by linearly propagating errors on parameters using the covariance matrix from...
static Int_t _evalErrorCount
Definition RooAbsReal.h:547
static void globalSelectComp(Bool_t flag)
Global switch controlling the activation of the selectComp() functionality.
RooAbsArg * createFundamental(const char *newname=0) const
Create a RooRealVar fundamental object with our properties.
virtual Bool_t setData(RooAbsData &, Bool_t=kTRUE)
Definition RooAbsReal.h:369
std::map< constRooAbsArg *, std::pair< std::string, std::list< EvalError > > >::const_iterator EvalErrorIter
Definition RooAbsReal.h:339
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
static Bool_t _globalSelectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition RooAbsReal.h:566
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
RooAbsReal * createIntObj(const RooArgSet &iset, const RooArgSet *nset, const RooNumIntConfig *cfg, const char *rangeName) const
Internal utility function for createIntegral() that creates the actual integral object.
RooAbsReal * createIntegral(const RooArgSet &iset, const RooArgSet &nset, const char *rangeName=0) const
Create integral over observables in iset in range named rangeName with integrand normalized over obse...
Definition RooAbsReal.h:219
Bool_t getForceNumInt() const
Definition RooAbsReal.h:181
void checkBatchComputation(const RooBatchCompute::RunContext &evalData, std::size_t evtNo, const RooArgSet *normSet=nullptr, double relAccuracy=1.E-13) const
Walk through expression tree headed by the this object, and check a batch computation.
static std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > > _evalErrorList
Definition RooAbsReal.h:546
virtual void forceNumInt(Bool_t flag=kTRUE)
Definition RooAbsReal.h:176
Bool_t plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
Double_t _plotMin
Definition RooAbsReal.h:481
const char * getPlotLabel() const
Get the label associated with the variable.
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.
virtual void fillTreeBranch(TTree &t)
Fill the tree branch that associated with this object with its current value.
Bool_t operator==(Double_t value) const
Equality operator comparing to a Double_t.
RooAbsFunc * bindVars(const RooArgSet &vars, const RooArgSet *nset=0, Bool_t clipInvalid=kFALSE) const
Create an interface adaptor f(vars) that binds us to the specified variables (in arbitrary order).
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const
Check if allArgs contains matching elements for each name in nameList.
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
virtual void preferredObservableScanOrder(const RooArgSet &obs, RooArgSet &orderedObs) const
Interface method for function objects to indicate their preferred order of observables for scanning t...
void findInnerMostIntegration(const RooArgSet &allObs, RooArgSet &innerObs, const char *rangeName) const
Utility function for createIntObj() that aids in the construct of recursive integrals over functions ...
virtual void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
void selectComp(Bool_t flag)
Definition RooAbsReal.h:560
RooAbsReal * createIntegral(const RooArgSet &iset, const char *rangeName) const
Create integral over observables in iset in range named rangeName.
Definition RooAbsReal.h:215
static Bool_t _hideOffset
Definition RooAbsReal.h:586
const Text_t * getUnit() const
Definition RooAbsReal.h:150
virtual void printValue(std::ostream &os) const
Print object value.
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:345
virtual bool isValid() const
Check if current value is valid.
Definition RooAbsReal.h:430
RooPlot * plotOnWithErrorBand(RooPlot *frame, const RooFitResult &fr, Double_t Z, const RooArgSet *params, const RooLinkedList &argList, Bool_t method1) const
Plot function or PDF on frame with support for visualization of the uncertainty encoded in the given ...
virtual Double_t getValV(const RooArgSet *normalisationSet=nullptr) const
Return value of object.
static void printEvalErrors(std::ostream &os=std::cout, Int_t maxPerNode=10000000)
Print all outstanding logged evaluation error on the given ostream.
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), Bool_t force=kTRUE)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
void logEvalError(const char *message, const char *serverValueString=0) const
Log evaluation error message.
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
RooAbsMoment * sigma(RooRealVar &obs, const RooArgSet &nset)
Definition RooAbsReal.h:364
Bool_t _selectComp
Definition RooAbsReal.h:565
RooFunctor * functor(const RooArgList &obs, const RooArgList &pars=RooArgList(), const RooArgSet &nset=RooArgSet()) const
Return a RooFunctor object bound to this RooAbsReal with given definition of observables and paramete...
void setCachedValue(double value, bool notifyClients=true) final
Overwrite the value stored in this object's cache.
Definition RooAbsReal.h:608
void makeProjectionSet(const RooAbsArg *plotVar, const RooArgSet *allVars, RooArgSet &projectedVars, Bool_t silent) const
Utility function for plotOn() that constructs the set of observables to project when plotting ourselv...
Int_t _plotBins
Definition RooAbsReal.h:483
const RooAbsReal * createPlotProjection(const RooArgSet &depVars, const RooArgSet &projVars, RooArgSet *&cloneSet) const
Utility function for plotOn() that creates a projection of a function or p.d.f to be plotted on a Roo...
void setPlotLabel(const char *label)
Set the label associated with this variable.
virtual ~RooAbsReal()
Destructor.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
RooAbsMoment * mean(RooRealVar &obs, const RooArgSet &nset)
Definition RooAbsReal.h:362
virtual void writeToStream(std::ostream &os, Bool_t compact) const
Write object contents to stream (dummy for now)
RooDataHist * fillDataHist(RooDataHist *hist, const RooArgSet *nset, Double_t scaleFactor, Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const
Fill a RooDataHist with values sampled from this function at the bin centers.
Double_t _plotMax
Definition RooAbsReal.h:482
virtual Double_t offset() const
Definition RooAbsReal.h:373
RooAbsMoment * sigma(RooRealVar &obs)
Definition RooAbsReal.h:363
RooAbsReal * createRunningIntegral(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Calls createRunningIntegral(const RooArgSet&, const RooCmdArg&, const RooCmdArg&, const RooCmdArg&,...
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:26
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgProxy is the abstract interface for RooAbsArg proxy classes.
Definition RooArgProxy.h:24
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
A RooCurve is a one-dimensional graphical representation of a real-valued function.
Definition RooCurve.h:32
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:45
adaptor that projects a real function via summation of states provided in a dataset.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:36
RooDerivative represents the first, second, or third order derivative of any RooAbsReal as calculated...
This class can evaluate a RooAbsReal object in other ways than recursive graph traversal.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Lightweight interface adaptor that exports a RooAbsPdf as a functor.
Definition RooFunctor.h:25
Lightweight interface adaptor that exports a RooAbsReal as a ROOT::Math::IGenFunction.
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Lightweight interface adaptor that exports a RooAbsReal as a ROOT::Math::IMultiGenFunction.
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
Lightweight interface adaptor that binds a RooAbsReal object to a subset of its servers and present i...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
A simple container to hold a batch of data values.
Definition RooSpan.h:34
RooVectorDataStore uses std::vectors to store data columns.
1-Dim function class
Definition TF1.h:213
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:251
3-D histogram with a float per channel (see TH1 documentation)}
Definition TH3.h:268
A doubly linked list.
Definition TList.h:38
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
A TTree represents a columnar dataset.
Definition TTree.h:79
Double_t x[n]
Definition legend1.C:17
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
BatchModeOption
For setting the batch mode flag with the BatchMode() command argument to RooAbsPdf::fitTo();.
RooCurve::WingMode wmode
Definition RooAbsReal.h:515
const char * normRangeName
Definition RooAbsReal.h:511
RooFit::MPSplit interleave
Definition RooAbsReal.h:523
const char * projectionRangeName
Definition RooAbsReal.h:516
const RooArgSet * projDataSet
Definition RooAbsReal.h:510
const char * curveNameSuffix
Definition RooAbsReal.h:524
const char * addToCurveName
Definition RooAbsReal.h:519
const RooFitResult * errorFR
Definition RooAbsReal.h:529
const RooArgSet * projSet
Definition RooAbsReal.h:507
const char * curveName
Definition RooAbsReal.h:518
const RooAbsData * projData
Definition RooAbsReal.h:505
Option_t * drawOptions
Definition RooAbsReal.h:502
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:32
auto * m
Definition textangle.C:8
static void output(int code)
Definition gifencode.c:226