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