Logo ROOT  
Reference Guide
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"
26#include "BatchData.h"
27
28class RooArgList ;
29class RooDataSet ;
30class RooPlot;
31class RooRealVar;
32class RooAbsFunc;
34class RooCategory ;
35class RooLinkedList ;
36class RooNumIntConfig ;
37class RooDataHist ;
38class RooFunctor ;
39class RooGenFunction ;
41class RooFitResult ;
42class RooAbsMoment ;
43class RooDerivative ;
45namespace RooHelpers {
47}
48
49class TH1;
50class TH1F;
51class TH2F;
52class TH3F;
53
54#include <list>
55#include <string>
56#include <iostream>
57#include <sstream>
58
59class RooAbsReal : public RooAbsArg {
60public:
61 // Constructors, assignment etc
62 RooAbsReal() ;
63 RooAbsReal(const char *name, const char *title, const char *unit= "") ;
64 RooAbsReal(const char *name, const char *title, Double_t minVal, Double_t maxVal,
65 const char *unit= "") ;
66 RooAbsReal(const RooAbsReal& other, const char* name=0);
67 RooAbsReal& operator=(const RooAbsReal& other);
68 virtual ~RooAbsReal();
69
70
71
72
73 //////////////////////////////////////////////////////////////////////////////////
74 /// Evaluate object. Returns either cached value or triggers a recalculation.
75 /// The recalculation happens by calling getValV(), which in the end calls the
76 /// virtual evaluate() functions of the respective PDFs.
77 /// \param[in] normalisationSet getValV() reacts differently depending on the value of the normalisation set.
78 /// If the set is `nullptr`, an unnormalised value is returned.
79 /// \note The normalisation is arbitrary, because it is up to the implementation
80 /// of the PDF to e.g. leave out normalisation constants for speed reasons. The range
81 /// of the variables is also ignored.
82 ///
83 /// To normalise the result properly, a RooArgSet has to be passed, which contains
84 /// the variables to normalise over.
85 /// These are integrated over their current ranges to compute the normalisation constant,
86 /// and the unnormalised result is divided by this value.
87 inline Double_t getVal(const RooArgSet* normalisationSet = nullptr) const {
88#ifdef ROOFIT_CHECK_CACHED_VALUES
89 return _DEBUG_getVal(normalisationSet);
90#else
91
92#ifndef _WIN32
93 return (_fast && !_inhibitDirty) ? _value : getValV(normalisationSet) ;
94#else
95 return (_fast && !inhibitDirty()) ? _value : getValV(normalisationSet) ;
96#endif
97
98#endif
99 }
100
101 /// Like getVal(const RooArgSet*), but always requires an argument for normalisation.
102 inline Double_t getVal(const RooArgSet& normalisationSet) const { return _fast ? _value : getValV(&normalisationSet) ; }
103
104 virtual Double_t getValV(const RooArgSet* normalisationSet = nullptr) const ;
105
106 virtual RooSpan<const double> getValBatch(std::size_t begin, std::size_t maxSize, const RooArgSet* normSet = nullptr) const;
107
108 Double_t getPropagatedError(const RooFitResult &fr, const RooArgSet &nset = RooArgSet()) const;
109
110 Bool_t operator==(Double_t value) const ;
111 virtual Bool_t operator==(const RooAbsArg& other) ;
112 virtual Bool_t isIdentical(const RooAbsArg& other, Bool_t assumeSameType=kFALSE) ;
113
114
115 inline const Text_t *getUnit() const {
116 // Return string with unit description
117 return _unit.Data();
118 }
119 inline void setUnit(const char *unit) {
120 // Set unit description to given string
121 _unit= unit;
122 }
123 TString getTitle(Bool_t appendUnit= kFALSE) const;
124
125 // Lightweight interface adaptors (caller takes ownership)
126 RooAbsFunc *bindVars(const RooArgSet &vars, const RooArgSet* nset=0, Bool_t clipInvalid=kFALSE) const;
127
128 // Create a fundamental-type object that can hold our value.
129 RooAbsArg *createFundamental(const char* newname=0) const;
130
131 // Analytical integration support
132 virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=0) const ;
133 virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
134 virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
135 virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
136 virtual Bool_t forceAnalyticalInt(const RooAbsArg& /*dep*/) const {
137 // Interface to force RooRealIntegral to offer given observable for internal integration
138 // even if this is deemed unsafe. This default implementation returns always flase
139 return kFALSE ;
140 }
141 virtual void forceNumInt(Bool_t flag=kTRUE) {
142 // If flag is true, all advertised analytical integrals will be ignored
143 // and all integrals are calculated numerically
144 _forceNumInt = flag ;
145 }
146 Bool_t getForceNumInt() const { return _forceNumInt ; }
147
148 // Chi^2 fits to histograms
149 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
150 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
151 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
152 virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
153
154 virtual RooAbsReal* createChi2(RooDataHist& data, const RooLinkedList& cmdList) ;
155 virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
156 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
157 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
158
159 // Chi^2 fits to X-Y datasets
160 virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
161 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
162 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
163 virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;
164
165 virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
166 virtual RooAbsReal* createChi2(RooDataSet& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
167 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
168 const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
169
170
171 virtual RooAbsReal* createProfile(const RooArgSet& paramsOfInterest) ;
172
173
174 RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
175 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
176 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
177 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
178
179 /// Create integral over observables in iset in range named rangeName.
180 RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const {
181 return createIntegral(iset,0,0,rangeName) ;
182 }
183 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
184 RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=0) const {
185 return createIntegral(iset,&nset,0,rangeName) ;
186 }
187 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
188 /// using specified configuration for any numeric integration.
189 RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
190 return createIntegral(iset,&nset,&cfg,rangeName) ;
191 }
192 /// Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration.
193 RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
194 return createIntegral(iset,0,&cfg,rangeName) ;
195 }
196 virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset=0, const RooNumIntConfig* cfg=0, const char* rangeName=0) const ;
197
198
199 void setParameterizeIntegral(const RooArgSet& paramVars) ;
200
201 // Create running integrals
202 RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
203 RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
204 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
205 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
206 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
207 RooAbsReal* createIntRI(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
208 RooAbsReal* createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
209
210
211 // Optimized accept/reject generator support
212 virtual Int_t getMaxVal(const RooArgSet& vars) const ;
213 virtual Double_t maxVal(Int_t code) const ;
214 virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }
215
216
217 // Plotting options
218 void setPlotLabel(const char *label);
219 const char *getPlotLabel() const;
220
221 virtual Double_t defaultErrorLevel() const {
222 // Return default level for MINUIT error analysis
223 return 1.0 ;
224 }
225
226 const RooNumIntConfig* getIntegratorConfig() const ;
231 void setIntegratorConfig() ;
232 void setIntegratorConfig(const RooNumIntConfig& config) ;
233
234 virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),Bool_t force=kTRUE) ;
235 virtual void fixAddCoefRange(const char* rangeName=0,Bool_t force=kTRUE) ;
236
237 virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;
238
239 // User entry point for plotting
240 virtual RooPlot* plotOn(RooPlot* frame,
241 const RooCmdArg& arg1=RooCmdArg(), const RooCmdArg& arg2=RooCmdArg(),
242 const RooCmdArg& arg3=RooCmdArg(), const RooCmdArg& arg4=RooCmdArg(),
243 const RooCmdArg& arg5=RooCmdArg(), const RooCmdArg& arg6=RooCmdArg(),
244 const RooCmdArg& arg7=RooCmdArg(), const RooCmdArg& arg8=RooCmdArg(),
245 const RooCmdArg& arg9=RooCmdArg(), const RooCmdArg& arg10=RooCmdArg()
246 ) const ;
247
248
250
251 // Forwarder function for backward compatibility
252 virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L",
253 Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=0) const;
254
255 // Fill an existing histogram
256 TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
257 Double_t scaleFactor= 1, const RooArgSet *projectedVars= 0, Bool_t scaling=kTRUE,
258 const RooArgSet* condObs=0, Bool_t setError=kTRUE) const;
259
260 // Create 1,2, and 3D histograms from and fill it
261 TH1 *createHistogram(const char* varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
262 TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
263 TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
264 const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
265 const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
266 const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
267 const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
268
269 // Fill a RooDataHist
270 RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, Double_t scaleFactor,
271 Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const ;
272
273 // I/O streaming interface (machine readable)
274 virtual Bool_t readFromStream(std::istream& is, Bool_t compact, Bool_t verbose=kFALSE) ;
275 virtual void writeToStream(std::ostream& os, Bool_t compact) const ;
276
277 // Printing interface (human readable)
278 virtual void printValue(std::ostream& os) const ;
279 virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
280
281 static void setCacheCheck(Bool_t flag) ;
282
283 // Evaluation error logging
284 class EvalError {
285 public:
287 EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
288 void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
289 void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
290 std::string _msg;
291 std::string _srvval;
292 } ;
293
297 void logEvalError(const char* message, const char* serverValueString=0) const ;
298 static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=0) ;
299 static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
300 static Int_t numEvalErrors() ;
301 static Int_t numEvalErrorItems() ;
302
303
304 typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ;
306
307 static void clearEvalErrorLog() ;
308
309 virtual Bool_t isBinnedDistribution(const RooArgSet& /*obs*/) const { return kFALSE ; }
310 virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const { return 0 ; }
311 virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const {
312 // Interface for returning an optional hint for initial sampling points when constructing a curve
313 // projected on observable.
314 return 0 ;
315 }
316
318 RooMultiGenFunction* iGenFunction(const RooArgSet& observables, const RooArgSet& nset=RooArgSet()) ;
319
320 RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
321 TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
322
323 RooDerivative* derivative(RooRealVar& obs, Int_t order=1, Double_t eps=0.001) ;
324 RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, Double_t eps=0.001) ;
325
326 RooAbsMoment* moment(RooRealVar& obs, Int_t order, Bool_t central, Bool_t takeRoot) ;
327 RooAbsMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, Bool_t central, Bool_t takeRoot, Bool_t intNormObs) ;
328
329 RooAbsMoment* mean(RooRealVar& obs) { return moment(obs,1,kFALSE,kFALSE) ; }
330 RooAbsMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,kFALSE,kFALSE,kTRUE) ; }
331 RooAbsMoment* sigma(RooRealVar& obs) { return moment(obs,2,kTRUE,kTRUE) ; }
332 RooAbsMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,kTRUE,kTRUE,kTRUE) ; }
333
335
336
337 virtual Bool_t setData(RooAbsData& /*data*/, Bool_t /*cloneData*/=kTRUE) { return kTRUE ; }
338
339 virtual void enableOffsetting(Bool_t) {} ;
340 virtual Bool_t isOffsetting() const { return kFALSE ; }
341 virtual Double_t offset() const { return 0 ; }
342
343 static void setHideOffset(Bool_t flag);
344 static Bool_t hideOffset() ;
345
346protected:
347 // Hook for objects with normalization-dependent parameters interperetation
348 virtual void selectNormalization(const RooArgSet* depSet=0, Bool_t force=kFALSE) ;
349 virtual void selectNormalizationRange(const char* rangeName=0, Bool_t force=kFALSE) ;
350
351 // Helper functions for plotting
352 Bool_t plotSanityChecks(RooPlot* frame) const ;
353 void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
354 RooArgSet& projectedVars, Bool_t silent) const ;
355
356 TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=0, const char* rangeName=0, Bool_t omitEmpty=kFALSE) const ;
357
358
359 Bool_t isSelectedComp() const ;
360
361
362 public:
363 const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars) const ;
364 const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
365 const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
366 RooArgSet *&cloneSet, const char* rangeName=0, const RooArgSet* condObs=0) const;
367 protected:
368
370
371 void plotOnCompSelect(RooArgSet* selNodes) const ;
372 RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, Double_t Z, const RooArgSet* params, const RooLinkedList& argList, Bool_t method1) const ;
373
374 // Support interface for subclasses to advertise their analytic integration
375 // and generator capabilities in their analticalIntegral() and generateEvent()
376 // implementations.
377 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
378 const RooArgProxy& a) const ;
379 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
380 const RooArgProxy& a, const RooArgProxy& b) const ;
381 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
382 const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
383 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
384 const RooArgProxy& a, const RooArgProxy& b,
385 const RooArgProxy& c, const RooArgProxy& d) const ;
386
387 Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
388 const RooArgSet& set) const ;
389
390
391 RooAbsReal* createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
392 void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;
393
394
395 // Internal consistency checking (needed by RooDataSet)
396 virtual Bool_t isValid() const ;
397 virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const ;
398
399 // Function evaluation and error tracing
400 Double_t traceEval(const RooArgSet* set) const ;
401
402 /// Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
403 virtual Double_t evaluate() const = 0;
404 virtual RooSpan<double> evaluateBatch(std::size_t begin, std::size_t maxSize) const;
405
406 //---------- Interface to access batch data ---------------------------
407 //
411 for (auto arg : _serverList) {
412 //TODO get rid of this cast?
413 auto absReal = dynamic_cast<RooAbsReal*>(arg);
414 if (absReal)
415 absReal->clearBatchMemory();
416 }
417 }
418
419 private:
420 void checkBatchComputation(std::size_t evtNo, const RooArgSet* normSet = nullptr, double relAccuracy = 1.E-13) const;
421
423 return _batchData;
424 }
425
426 /// Debug version of getVal(), which is slow and does error checking.
427 Double_t _DEBUG_getVal(const RooArgSet* normalisationSet) const;
428
429 //--------------------------------------------------------------------
430
431 protected:
432 // Hooks for RooDataSet interface
433 friend class RooRealIntegral ;
434 friend class RooVectorDataStore ;
435 virtual void syncCache(const RooArgSet* set=0) { getVal(set) ; }
436 virtual void copyCache(const RooAbsArg* source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE) ;
437 virtual void attachToTree(TTree& t, Int_t bufSize=32000) ;
438 virtual void attachToVStore(RooVectorDataStore& vstore) ;
439 virtual void setTreeBranchStatus(TTree& t, Bool_t active) ;
440 virtual void fillTreeBranch(TTree& t) ;
441
442 friend class RooRealBinding ;
443 Double_t _plotMin ; // Minimum of plot range
444 Double_t _plotMax ; // Maximum of plot range
445 Int_t _plotBins ; // Number of plot bins
446 mutable Double_t _value ; // Cache for current value of object
447 mutable BatchHelpers::BatchData _batchData; //! Value storage for batches of events
448 TString _unit ; // Unit for objects value
449 TString _label ; // Plot label for objects value
450 Bool_t _forceNumInt ; // Force numerical integration if flag set
451
452 mutable Float_t _floatValue{0.}; //! Transient cache for floating point values from tree branches
453 mutable Int_t _intValue{0}; //! Transient cache for integer values from tree branches
454 mutable Bool_t _boolValue{false}; //! Transient cache for bool values from tree branches
455 mutable UChar_t _byteValue{'\0'}; //! Transient cache for byte values from tree branches
456 mutable Char_t _sbyteValue{'\0'}; //! Transient cache for signed byte values from tree branches
457 mutable UInt_t _uintValue{0u}; //! Transient cache for unsigned integer values from tree branches
458
459 friend class RooAbsPdf ;
460 friend class RooAbsAnaConvPdf ;
461
462 RooNumIntConfig* _specIntegratorConfig ; // Numeric integrator configuration specific for this object
463
464 Bool_t _treeVar ; // !do not persist
465
466 friend class RooDataProjBinding ;
468
469 struct PlotOpt {
483 const char* normRangeName ;
490 const char* curveName ;
491 const char* addToCurveName ;
496 const char* curveNameSuffix ;
502 } ;
503
504 // Plot implementation functions
505 virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;
506
507public:
508 // PlotOn with command list
509 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
510
511 protected:
512 virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;
513
514
515private:
516
518 static std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > > _evalErrorList ;
520
521 Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;
522
523protected:
524
525
526 friend class RooRealSumPdf ;
527 friend class RooRealSumFunc;
528 friend class RooAddPdf ;
529 friend class RooAddModel ;
530 void selectComp(Bool_t flag) {
531 // If flag is true, only selected component will be included in evaluates of RooAddPdf components
532 _selectComp = flag ;
533 }
534 static void globalSelectComp(Bool_t flag) ;
535 Bool_t _selectComp ; //! Component selection flag for RooAbsPdf::plotCompOn
536 static Bool_t _globalSelectComp ; // Global activation switch for component selection
537 // This struct can be used to flip the global switch to select components.
538 // Doing this with RAII prevents forgetting to reset the state.
544 }
545
549 }
550
552 };
553
554
555 mutable RooArgSet* _lastNSet ; //!
556 static Bool_t _hideOffset ; // Offset hiding flag
557
558 ClassDef(RooAbsReal,2) // Abstract real-valued variable
559};
560
561#endif
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define c(i)
Definition: RSha256.hxx:101
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
char Text_t
Definition: RtypesCore.h:58
unsigned char UChar_t
Definition: RtypesCore.h:34
char Char_t
Definition: RtypesCore.h:29
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassDef(name, id)
Definition: Rtypes.h:326
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
void clear()
Discard all storage.
Definition: BatchData.h:40
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 (of arbitrary type) an...
Definition: RooAbsArg.h:71
Bool_t _fast
Definition: RooAbsArg.h:636
friend class RooArgSet
Definition: RooAbsArg.h:551
Bool_t inhibitDirty() const
Delete watch flag.
Definition: RooAbsArg.cxx:85
static Bool_t _inhibitDirty
Definition: RooAbsArg.h:620
RefCountList_t _serverList
Definition: RooAbsArg.h:559
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:39
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:23
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
Definition: RooAbsMoment.h:27
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:289
std::string _srvval
Definition: RooAbsReal.h:291
void setMessage(const char *tmp)
Definition: RooAbsReal.h:288
EvalError(const EvalError &other)
Definition: RooAbsReal.h:287
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
void plotOnCompSelect(RooArgSet *selNodes) const
Helper function for plotting of composite p.d.fs.
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)
Definition: RooAbsReal.cxx:448
TString _label
Definition: RooAbsReal.h:449
@ RelativeExpected
Definition: RooAbsReal.h:249
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:193
virtual Bool_t isIdentical(const RooAbsArg &other, Bool_t assumeSameType=kFALSE)
Definition: RooAbsReal.cxx:239
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:329
void clearBatchMemory()
Definition: RooAbsReal.h:409
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.
Definition: RooAbsReal.cxx:254
virtual Bool_t isOffsetting() const
Definition: RooAbsReal.h:340
static Int_t numEvalErrorItems()
Definition: RooAbsReal.cxx:324
static void setHideOffset(Bool_t flag)
Definition: RooAbsReal.cxx:117
static void setCacheCheck(Bool_t flag)
RooAbsReal()
coverity[UNINIT_CTOR] Default constructor
Definition: RooAbsReal.cxx:129
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 ...
Definition: RooAbsReal.cxx:373
Bool_t _forceNumInt
Definition: RooAbsReal.h:450
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:221
UChar_t _byteValue
Transient cache for bool values from tree branches.
Definition: RooAbsReal.h:455
Bool_t isSelectedComp() const
If true, the current pdf is a selected component (for use in plotting)
virtual RooAbsReal * createChi2(RooDataHist &data, const RooLinkedList &cmdList)
Create a variable from a histogram and this function.
RooArgSet * _lastNSet
Definition: RooAbsReal.h:555
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.
Float_t _floatValue
Definition: RooAbsReal.h:452
virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const
Interface function to check if given value is a valid value for this object.
Definition: RooAbsReal.cxx:504
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:311
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.
void setParameterizeIntegral(const RooArgSet &paramVars)
static Bool_t hideOffset()
Definition: RooAbsReal.cxx:118
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.
Definition: RooAbsReal.cxx:517
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:189
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...
Definition: RooAbsReal.cxx:809
virtual RooSpan< double > evaluateBatch(std::size_t begin, std::size_t maxSize) const
Evaluate function for a batch of input data points.
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.
Definition: RooAbsReal.cxx:341
Int_t _intValue
Transient cache for floating point values from tree branches.
Definition: RooAbsReal.h:453
Double_t getVal(const RooArgSet &normalisationSet) const
Like getVal(const RooArgSet*), but always requires an argument for normalisation.
Definition: RooAbsReal.h:102
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()
Definition: RooAbsReal.cxx:332
RooAbsReal & operator=(const RooAbsReal &other)
Assign values, name and configs from another RooAbsReal.
Definition: RooAbsReal.cxx:181
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.
Definition: RooAbsReal.cxx:402
virtual RooSpan< const double > getValBatch(std::size_t begin, std::size_t maxSize, const RooArgSet *normSet=nullptr) const
Return value of object for all data events in the batch.
Definition: RooAbsReal.cxx:296
virtual void attachToVStore(RooVectorDataStore &vstore)
TString _unit
Value storage for batches of events.
Definition: RooAbsReal.h:448
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:467
Double_t _value
Definition: RooAbsReal.h:446
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:310
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:462
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
static ErrorLoggingMode evalErrorLoggingMode()
Return current evaluation error logging mode.
virtual Bool_t forceAnalyticalInt(const RooAbsArg &) const
Definition: RooAbsReal.h:136
virtual Int_t minTrialSamples(const RooArgSet &) const
Definition: RooAbsReal.h:214
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:119
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Definition: RooAbsReal.cxx:477
virtual void syncCache(const RooArgSet *set=0)
Definition: RooAbsReal.h:435
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
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:517
virtual void enableOffsetting(Bool_t)
Definition: RooAbsReal.h:339
Double_t getPropagatedError(const RooFitResult &fr, const RooArgSet &nset=RooArgSet()) const
Calculate error on self by propagated errors on parameters with correlations as given by fit result T...
static Int_t _evalErrorCount
Definition: RooAbsReal.h:519
virtual Bool_t isValid() const
Check if current value is valid.
Definition: RooAbsReal.cxx:493
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:337
std::map< constRooAbsArg *, std::pair< std::string, std::list< EvalError > > >::const_iterator EvalErrorIter
Definition: RooAbsReal.h:304
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:536
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 ...
Definition: RooAbsReal.cxx:562
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.
Definition: RooAbsReal.cxx:633
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:184
BatchHelpers::BatchData _batchData
Definition: RooAbsReal.h:447
Bool_t getForceNumInt() const
Definition: RooAbsReal.h:146
static std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > > _evalErrorList
Definition: RooAbsReal.h:518
virtual void forceNumInt(Bool_t flag=kTRUE)
Definition: RooAbsReal.h:141
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:443
Char_t _sbyteValue
Transient cache for byte values from tree branches.
Definition: RooAbsReal.h:456
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:428
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.
Definition: RooAbsReal.cxx:219
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).
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.
Definition: RooAbsReal.cxx:416
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition: RooAbsReal.cxx:390
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 ...
Definition: RooAbsReal.cxx:754
const BatchHelpers::BatchData & batchData() const
Definition: RooAbsReal.h:422
void selectComp(Bool_t flag)
Definition: RooAbsReal.h:530
RooAbsReal * createIntegral(const RooArgSet &iset, const char *rangeName) const
Create integral over observables in iset in range named rangeName.
Definition: RooAbsReal.h:180
static Bool_t _hideOffset
Definition: RooAbsReal.h:556
const Text_t * getUnit() const
Definition: RooAbsReal.h:115
virtual void printValue(std::ostream &os) const
Print object value.
Definition: RooAbsReal.cxx:467
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
Bool_t _boolValue
Transient cache for integer values from tree branches.
Definition: RooAbsReal.h:454
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:309
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.
Definition: RooAbsReal.cxx:272
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:332
void checkBatchComputation(std::size_t evtNo, const RooArgSet *normSet=nullptr, double relAccuracy=1.E-13) const
Bool_t _selectComp
Definition: RooAbsReal.h:535
Bool_t _treeVar
Definition: RooAbsReal.h:464
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...
UInt_t _uintValue
Transient cache for signed byte values from tree branches.
Definition: RooAbsReal.h:457
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:445
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:438
virtual ~RooAbsReal()
Destructor.
Definition: RooAbsReal.cxx:209
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:330
virtual void writeToStream(std::ostream &os, Bool_t compact) const
Write object contents to stream (dummy for now)
Definition: RooAbsReal.cxx:458
const RooAbsReal * createPlotProjection(const RooArgSet &depVars, const RooArgSet &projVars) const
Utility function for plotOn() that creates a projection of a function or p.d.f to be plotted on a Roo...
Definition: RooAbsReal.cxx:888
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:444
virtual Double_t offset() const
Definition: RooAbsReal.h:341
RooAbsMoment * sigma(RooRealVar &obs)
Definition: RooAbsReal.h:331
RooAbsReal * createRunningIntegral(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Calls createRunningIntegral(const RooArgSet&, const RooCmdArg&, const RooCmdArg&, const RooCmdArg&,...
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgProxy is the abstact interface for RooAbsArg proxy classes.
Definition: RooArgProxy.h:24
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:28
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:40
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:31
RooDerivative represents the first, second, or third order derivative of any RooAbsReal as calculated...
Definition: RooDerivative.h:31
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
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.
Helper class to access a batch-related part of RooAbsReal's interface, which should not leak to the o...
Definition: RooHelpers.h:159
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:36
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:
Definition: RooRealSumPdf.h:24
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
A simple container to hold a batch of data values.
Definition: RooSpan.h:32
RooVectorDataStore is the abstract base class for data collection that use a TTree as internal storag...
1-Dim function class
Definition: TF1.h:211
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
The TH1 histogram class.
Definition: TH1.h:56
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:267
A doubly linked list.
Definition: TList.h:44
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
A TTree represents a columnar dataset.
Definition: TTree.h:72
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg Extended(Bool_t flag=kTRUE)
@ Interleave
Definition: RooGlobalFunc.h:70
static constexpr double s
Bool_t postRangeFracScale
Definition: RooAbsReal.h:486
RooCurve::WingMode wmode
Definition: RooAbsReal.h:487
const char * normRangeName
Definition: RooAbsReal.h:483
RooFit::MPSplit interleave
Definition: RooAbsReal.h:495
const char * projectionRangeName
Definition: RooAbsReal.h:488
Double_t scaleFactor
Definition: RooAbsReal.h:475
const RooArgSet * projDataSet
Definition: RooAbsReal.h:482
Double_t addToWgtSelf
Definition: RooAbsReal.h:492
const char * curveNameSuffix
Definition: RooAbsReal.h:496
Double_t addToWgtOther
Definition: RooAbsReal.h:493
const char * addToCurveName
Definition: RooAbsReal.h:491
const RooFitResult * errorFR
Definition: RooAbsReal.h:501
const RooArgSet * projSet
Definition: RooAbsReal.h:479
const char * curveName
Definition: RooAbsReal.h:490
const RooAbsData * projData
Definition: RooAbsReal.h:477
Option_t * drawOptions
Definition: RooAbsReal.h:473
auto * m
Definition: textangle.C:8
auto * a
Definition: textangle.C:12