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"
27
28#include <ROOT/RSpan.hxx>
29
30class RooDataSet ;
31class RooPlot;
32class RooRealVar;
33class RooAbsFunc;
35class RooLinkedList ;
36class RooNumIntConfig ;
37class RooDataHist ;
38class RooFunctor ;
39class RooFitResult ;
40class RooAbsMoment ;
41class RooDerivative ;
43struct TreeReadBuffer; /// A space to attach TBranches
44namespace RooBatchCompute {
45struct RunContext;
46}
47
48class TH1;
49class TH1F;
50class TH2F;
51class TH3F;
52
53#include <iostream>
54#include <list>
55#include <map>
56#include <string>
57#include <sstream>
58
59class RooAbsReal : public RooAbsArg {
60public:
62
63 /// A RooAbsReal::Ref can be constructed from a `RooAbsReal&` or a `double`
64 /// that will be implicitly converted to a RooConstVar&. The RooAbsReal::Ref
65 /// can be used as a replacement for `RooAbsReal&`. With this type
66 /// definition, you can write RooFit interfaces that accept both RooAbsReal,
67 /// or simply a number that will be implicitly converted to a RooConstVar&.
68 class Ref {
69 public:
70 inline Ref(RooAbsReal &ref) : _ref{ref} {}
71 Ref(double val);
72 inline operator RooAbsReal &() const { return _ref; }
73
74 private:
76 };
77
78 // Constructors, assignment etc
79 RooAbsReal() ;
80 RooAbsReal(const char *name, const char *title, const char *unit= "") ;
81 RooAbsReal(const char *name, const char *title, double minVal, double maxVal,
82 const char *unit= "") ;
83 RooAbsReal(const RooAbsReal& other, const char* name=nullptr);
84 ~RooAbsReal() override;
85
86
87
88
89 //////////////////////////////////////////////////////////////////////////////////
90 /// Evaluate object. Returns either cached value or triggers a recalculation.
91 /// The recalculation happens by calling getValV(), which in the end calls the
92 /// virtual evaluate() functions of the respective PDFs.
93 /// \param[in] normalisationSet getValV() reacts differently depending on the value of the normalisation set.
94 /// If the set is `nullptr`, an unnormalised value is returned.
95 /// \note The normalisation is arbitrary, because it is up to the implementation
96 /// of the PDF to e.g. leave out normalisation constants for speed reasons. The range
97 /// of the variables is also ignored.
98 ///
99 /// To normalise the result properly, a RooArgSet has to be passed, which contains
100 /// the variables to normalise over.
101 /// These are integrated over their current ranges to compute the normalisation constant,
102 /// and the unnormalised result is divided by this value.
103 inline double getVal(const RooArgSet* normalisationSet = nullptr) const {
104 // Sometimes, the calling code uses an empty RooArgSet to request evaluation
105 // without normalization set instead of following the `nullptr` convention.
106 // To remove this ambiguity which might not always be correctly handled in
107 // downstream code, we set `normalisationSet` to nullptr if it is pointing
108 // to an empty set.
109 if(normalisationSet && normalisationSet->empty()) {
110 normalisationSet = nullptr;
111 }
112#ifdef ROOFIT_CHECK_CACHED_VALUES
113 return _DEBUG_getVal(normalisationSet);
114#else
115
116#ifndef _WIN32
117 return (_fast && !_inhibitDirty) ? _value : getValV(normalisationSet) ;
118#else
119 return (_fast && !inhibitDirty()) ? _value : getValV(normalisationSet) ;
120#endif
121
122#endif
123 }
124
125 /// Like getVal(const RooArgSet*), but always requires an argument for normalisation.
126 inline double getVal(const RooArgSet& normalisationSet) const {
127 // Sometimes, the calling code uses an empty RooArgSet to request evaluation
128 // without normalization set instead of following the `nullptr` convention.
129 // To remove this ambiguity which might not always be correctly handled in
130 // downstream code, we set `normalisationSet` to nullptr if it is an empty set.
131 return _fast ? _value : getValV(normalisationSet.empty() ? nullptr : &normalisationSet) ;
132 }
133
134 virtual double getValV(const RooArgSet* normalisationSet = nullptr) const ;
135
136 double getPropagatedError(const RooFitResult &fr, const RooArgSet &nset = {}) const;
137
138 bool operator==(double value) const ;
139 bool operator==(const RooAbsArg& other) const override;
140 bool isIdentical(const RooAbsArg& other, bool assumeSameType=false) const override;
141
142
143 inline const Text_t *getUnit() const {
144 // Return string with unit description
145 return _unit.Data();
146 }
147 inline void setUnit(const char *unit) {
148 // Set unit description to given string
149 _unit= unit;
150 }
151 TString getTitle(bool appendUnit= false) const;
152
153 // Lightweight interface adaptors (caller takes ownership)
154 RooFit::OwningPtr<RooAbsFunc> bindVars(const RooArgSet &vars, const RooArgSet* nset=nullptr, bool clipInvalid=false) const;
155
156 // Create a fundamental-type object that can hold our value.
157 RooFit::OwningPtr<RooAbsArg> createFundamental(const char* newname=nullptr) const override;
158
159 // Analytical integration support
160 virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=nullptr) const ;
161 virtual double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=nullptr) const ;
162 virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const ;
163 virtual double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const ;
164 virtual bool forceAnalyticalInt(const RooAbsArg& /*dep*/) const {
165 // Interface to force RooRealIntegral to offer given observable for internal integration
166 // even if this is deemed unsafe. This default implementation returns always false
167 return false ;
168 }
169 virtual void forceNumInt(bool flag=true) {
170 // If flag is true, all advertised analytical integrals will be ignored
171 // and all integrals are calculated numerically
172 _forceNumInt = flag ;
173 }
174 bool getForceNumInt() const { return _forceNumInt ; }
175
176 // Chi^2 fits to histograms
177 virtual RooFit::OwningPtr<RooFitResult> chi2FitTo(RooDataHist& data, const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
178 const RooCmdArg& arg3={}, const RooCmdArg& arg4={}, const RooCmdArg& arg5={},
179 const RooCmdArg& arg6={}, const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
181
183 virtual RooFit::OwningPtr<RooAbsReal> createChi2(RooDataHist& data, const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
184 const RooCmdArg& arg3={}, const RooCmdArg& arg4={}, const RooCmdArg& arg5={},
185 const RooCmdArg& arg6={}, const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
186
187 // Chi^2 fits to X-Y datasets
188 virtual RooFit::OwningPtr<RooFitResult> chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
189 const RooCmdArg& arg3={}, const RooCmdArg& arg4={}, const RooCmdArg& arg5={},
190 const RooCmdArg& arg6={}, const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
191 virtual RooFit::OwningPtr<RooFitResult> chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;
192
194 virtual RooFit::OwningPtr<RooAbsReal> createChi2(RooDataSet& data, const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
195 const RooCmdArg& arg3={}, const RooCmdArg& arg4={}, const RooCmdArg& arg5={},
196 const RooCmdArg& arg6={}, const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
197
198 virtual RooFit::OwningPtr<RooAbsReal> createProfile(const RooArgSet& paramsOfInterest) ;
199
200
201 RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2={},
202 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
203 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
204 const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) const ;
205
206 /// Create integral over observables in iset in range named rangeName.
207 RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const char* rangeName) const {
208 return createIntegral(iset,nullptr,nullptr,rangeName) ;
209 }
210 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
211 RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=nullptr) const {
212 return createIntegral(iset,&nset,nullptr,rangeName) ;
213 }
214 /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
215 /// using specified configuration for any numeric integration.
216 RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=nullptr) const {
217 return createIntegral(iset,&nset,&cfg,rangeName) ;
218 }
219 /// Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration.
220 RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=nullptr) const {
221 return createIntegral(iset,nullptr,&cfg,rangeName) ;
222 }
223 virtual RooFit::OwningPtr<RooAbsReal> createIntegral(const RooArgSet& iset, const RooArgSet* nset=nullptr, const RooNumIntConfig* cfg=nullptr, const char* rangeName=nullptr) const ;
224
225
226 void setParameterizeIntegral(const RooArgSet& paramVars) ;
227
228 // Create running integrals
231 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
232 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
233 const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) ;
234 RooFit::OwningPtr<RooAbsReal> createIntRI(const RooArgSet& iset, const RooArgSet& nset={}) ;
235 RooFit::OwningPtr<RooAbsReal> createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
236
237
238 // Optimized accept/reject generator support
239 virtual Int_t getMaxVal(const RooArgSet& vars) const ;
240 virtual double maxVal(Int_t code) const ;
241 virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }
242
243
244 // Plotting options
245 void setPlotLabel(const char *label);
246 const char *getPlotLabel() const;
247
248 virtual double defaultErrorLevel() const {
249 // Return default level for MINUIT error analysis
250 return 1.0 ;
251 }
252
253 const RooNumIntConfig* getIntegratorConfig() const ;
257 RooNumIntConfig* specialIntegratorConfig(bool createOnTheFly) ;
258 void setIntegratorConfig() ;
259 void setIntegratorConfig(const RooNumIntConfig& config) ;
260
261 virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),bool force=true) ;
262 virtual void fixAddCoefRange(const char* rangeName=nullptr,bool force=true) ;
263
264 virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;
265
266 // User entry point for plotting
267 virtual RooPlot* plotOn(RooPlot* frame,
268 const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
269 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
270 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
271 const RooCmdArg& arg7={}, const RooCmdArg& arg8={},
272 const RooCmdArg& arg9={}, const RooCmdArg& arg10={}
273 ) const ;
274
275
277
278 // Forwarder function for backward compatibility
279 virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L",
280 double scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=nullptr) const;
281
282 // Fill an existing histogram
283 TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
284 double scaleFactor= 1, const RooArgSet *projectedVars= nullptr, bool scaling=true,
285 const RooArgSet* condObs=nullptr, bool setError=true) const;
286
287 // Create 1,2, and 3D histograms from and fill it
288 TH1 *createHistogram(RooStringView varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
289 TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
290 TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
291 const RooCmdArg& arg1={}, const RooCmdArg& arg2={},
292 const RooCmdArg& arg3={}, const RooCmdArg& arg4={},
293 const RooCmdArg& arg5={}, const RooCmdArg& arg6={},
294 const RooCmdArg& arg7={}, const RooCmdArg& arg8={}) const ;
295
296 // Fill a RooDataHist
297 RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, double scaleFactor,
298 bool correctForBinVolume=false, bool showProgress=false) const ;
299
300 // I/O streaming interface (machine readable)
301 bool readFromStream(std::istream& is, bool compact, bool verbose=false) override ;
302 void writeToStream(std::ostream& os, bool compact) const override ;
303
304 // Printing interface (human readable)
305 void printValue(std::ostream& os) const override ;
306 void printMultiline(std::ostream& os, Int_t contents, bool verbose=false, TString indent="") const override ;
307
308 inline void setCachedValue(double value, bool notifyClients = true) final;
309
310 // Evaluation error logging
311 class EvalError {
312 public:
314 EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
315 void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
316 void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
317 std::string _msg;
318 std::string _srvval;
319 } ;
320
322
323 /// Context to temporarily change the error logging mode as long as the context is alive.
325 public:
327
332
334 private:
336 };
337
340 void logEvalError(const char* message, const char* serverValueString=nullptr) const ;
341 static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=nullptr) ;
342 static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
343 static Int_t numEvalErrors() ;
344 static Int_t numEvalErrorItems() ;
345
346
347 typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ;
349
350 static void clearEvalErrorLog() ;
351
352 /// Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
353 virtual bool isBinnedDistribution(const RooArgSet& /*obs*/) const { return false ; }
354 virtual std::list<double>* binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const;
355 virtual std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const;
356
357 RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
358 TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
359
360 RooDerivative* derivative(RooRealVar& obs, Int_t order=1, double eps=0.001) ;
361 RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, double eps=0.001) ;
362
363 RooAbsMoment* moment(RooRealVar& obs, Int_t order, bool central, bool takeRoot) ;
364 RooAbsMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, bool central, bool takeRoot, bool intNormObs) ;
365
366 RooAbsMoment* mean(RooRealVar& obs) { return moment(obs,1,false,false) ; }
367 RooAbsMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,false,false,true) ; }
368 RooAbsMoment* sigma(RooRealVar& obs) { return moment(obs,2,true,true) ; }
369 RooAbsMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,true,true,true) ; }
370
371 double findRoot(RooRealVar& x, double xmin, double xmax, double yval) ;
372
373
374 virtual bool setData(RooAbsData& /*data*/, bool /*cloneData*/=true) { return true ; }
375
376 virtual void enableOffsetting(bool);
377 virtual bool isOffsetting() const { return false ; }
378 virtual double offset() const { return 0 ; }
379
380 static void setHideOffset(bool flag);
381 static bool hideOffset() ;
382
383 bool isSelectedComp() const ;
384 void selectComp(bool flag) {
385 // If flag is true, only selected component will be included in evaluates of RooAddPdf components
386 _selectComp = flag ;
387 }
388
389 const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
390 const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
391 RooArgSet *&cloneSet, const char* rangeName=nullptr, const RooArgSet* condObs=nullptr) const;
392 virtual void computeBatch(double* output, size_t size, RooFit::Detail::DataMap const&) const;
393
394 virtual bool hasGradient() const { return false; }
395 virtual void gradient(double *) const {
396 if(!hasGradient()) throw std::runtime_error("RooAbsReal::gradient(double *) not implemented by this class!");
397 }
398
399 virtual std::string
400 buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const;
401
402 // PlotOn with command list
403 virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
404
405protected:
407 friend class RooVectorDataStore;
408 friend class RooRealBinding;
409 friend class RooRealSumPdf;
410 friend class RooRealSumFunc;
411 friend class RooAddHelpers;
412 friend class RooAddPdf;
413 friend class RooAddModel;
414 friend class AddCacheElem;
416
417 // Hook for objects with normalization-dependent parameters interpretation
418 virtual void selectNormalization(const RooArgSet* depSet=nullptr, bool force=false) ;
419 virtual void selectNormalizationRange(const char* rangeName=nullptr, bool force=false) ;
420
421 // Helper functions for plotting
422 bool plotSanityChecks(RooPlot* frame) const ;
423 void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
424 RooArgSet& projectedVars, bool silent) const ;
425
426 TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=nullptr, const char* rangeName=nullptr, bool omitEmpty=false) const ;
427
428 void plotOnCompSelect(RooArgSet* selNodes) const ;
429 RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, double Z, const RooArgSet* params, const RooLinkedList& argList, bool method1) const ;
430
431 // Support interface for subclasses to advertise their analytic integration
432 // and generator capabilities in their analyticalIntegral() and generateEvent()
433 // implementations.
434 bool matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
435 const RooArgProxy& a) const ;
436 bool matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
437 const RooArgProxy& a, const RooArgProxy& b) const ;
438 bool matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
439 const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
440 bool matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
441 const RooArgProxy& a, const RooArgProxy& b,
442 const RooArgProxy& c, const RooArgProxy& d) const ;
443
444 bool matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
445 const RooArgSet& set) const ;
446
447 RooFit::OwningPtr<RooAbsReal> createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
448 void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;
449
450 // Internal consistency checking (needed by RooDataSet)
451 /// Check if current value is valid.
452 bool isValid() const override { return isValidReal(_value); }
453 /// Interface function to check if given value is a valid value for this object. Returns true unless overridden.
454 virtual bool isValidReal(double /*value*/, bool printError = false) const { (void)printError; return true; }
455
456 // Function evaluation and error tracing
457 double traceEval(const RooArgSet* set) const ;
458
459 /// Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
460 virtual double evaluate() const = 0;
461
462 // Hooks for RooDataSet interface
463 void syncCache(const RooArgSet* set=nullptr) override { getVal(set) ; }
464 void copyCache(const RooAbsArg* source, bool valueOnly=false, bool setValDirty=true) override ;
465 void attachToTree(TTree& t, Int_t bufSize=32000) override ;
466 void attachToVStore(RooVectorDataStore& vstore) override ;
467 void setTreeBranchStatus(TTree& t, bool active) override ;
468 void fillTreeBranch(TTree& t) override ;
469
470 struct PlotOpt {
472 double scaleFactor = 1.0;
474 const RooAbsData *projData = nullptr;
475 bool binProjData = false;
476 const RooArgSet *projSet = nullptr;
477 double precision = 1e-3;
478 bool shiftToZero = false;
479 const RooArgSet *projDataSet = nullptr;
480 const char *normRangeName = nullptr;
481 double rangeLo = 0.0;
482 double rangeHi = 0.0;
483 bool postRangeFracScale = false;
485 const char *projectionRangeName = nullptr;
486 bool curveInvisible = false;
487 const char *curveName = nullptr;
488 const char *addToCurveName = nullptr;
489 double addToWgtSelf = 1.0;
490 double addToWgtOther = 1.0;
493 const char *curveNameSuffix = "";
495 double eeval = 0.0;
496 bool doeeval = false;
497 bool progress = false;
498 const RooFitResult *errorFR = nullptr;
499 };
500
501 // Plot implementation functions
502 virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;
503
504 virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;
505
506 bool matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;
507
508 bool redirectServersHook(const RooAbsCollection & newServerList, bool mustReplaceAll,
509 bool nameChange, bool isRecursiveStep) override;
510
511 static void globalSelectComp(bool flag) ;
512
513 // This struct can be used to flip the global switch to select components.
514 // Doing this with RAII prevents forgetting to reset the state.
520 }
521
525 }
526
528 };
529
530
531private:
532
533 /// Debug version of getVal(), which is slow and does error checking.
534 double _DEBUG_getVal(const RooArgSet* normalisationSet) const;
535
536 //--------------------------------------------------------------------
537
538 protected:
539
540 double _plotMin = 0.0; ///< Minimum of plot range
541 double _plotMax = 0.0; ///< Maximum of plot range
542 Int_t _plotBins = 100; ///< Number of plot bins
543 mutable double _value = 0.0; ///< Cache for current value of object
544 TString _unit; ///< Unit for objects value
545 TString _label; ///< Plot label for objects value
546 bool _forceNumInt = false; ///< Force numerical integration if flag set
547 std::unique_ptr<RooNumIntConfig> _specIntegratorConfig; // Numeric integrator configuration specific for this object
548 std::unique_ptr<TreeReadBuffer> _treeReadBuffer; //! A buffer for reading values from trees
549 bool _selectComp = true; //! Component selection flag for RooAbsPdf::plotCompOn
551
553 static std::map<const RooAbsArg *, std::pair<std::string, std::list<EvalError>>> _evalErrorList;
555 static bool _globalSelectComp; // Global activation switch for component selection
556 static bool _hideOffset; ///< Offset hiding flag
557
558 ClassDefOverride(RooAbsReal,3); // Abstract real-valued variable
559};
560
561
562////////////////////////////////////////////////////////////////////////////////
563/// Overwrite the value stored in this object's cache.
564/// This can be used to fake a computation that resulted in `value`.
565/// \param[in] value Value to write.
566/// \param[in] notifyClients If true, notify users of this object that its value changed.
567/// This is the default.
568void RooAbsReal::setCachedValue(double value, bool notifyClients) {
569 _value = value;
570
571 if (notifyClients) {
573 _valueDirty = false;
574 }
575}
576
577
578#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
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
bool _fast
Definition RooAbsArg.h:717
bool _valueDirty
Definition RooAbsArg.h:713
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:490
bool inhibitDirty() const
Delete watch flag.
static bool _inhibitDirty
Definition RooAbsArg.h:696
Abstract base class for objects that represent a discrete value that can be set from the outside,...
Abstract container object that can hold multiple RooAbsArg objects.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
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:324
EvalErrorContext(ErrorLoggingMode m)
Definition RooAbsReal.h:326
EvalErrorContext & operator=(EvalErrorContext const &)=delete
EvalErrorContext & operator=(EvalErrorContext &&)=delete
EvalErrorContext(EvalErrorContext const &)=delete
EvalErrorContext(EvalErrorContext &&)=delete
void setServerValues(const char *tmp)
Definition RooAbsReal.h:316
void setMessage(const char *tmp)
Definition RooAbsReal.h:315
EvalError(const EvalError &other)
Definition RooAbsReal.h:314
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:68
Ref(RooAbsReal &ref)
Definition RooAbsReal.h:70
RooAbsReal & _ref
Definition RooAbsReal.h:75
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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:548
virtual std::list< double > * binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
void selectComp(bool flag)
Definition RooAbsReal.h:384
TString _label
Plot label for objects value.
Definition RooAbsReal.h:545
bool _selectComp
A buffer for reading values from trees.
Definition RooAbsReal.h:549
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
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:366
virtual bool isOffsetting() const
Definition RooAbsReal.h:377
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()
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...
static std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > > _evalErrorList
Definition RooAbsReal.h:553
virtual bool hasGradient() const
Definition RooAbsReal.h:394
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 ...
bool isValid() const override
Check if current value is valid.
Definition RooAbsReal.h:452
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:546
~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()
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:207
virtual double defaultErrorLevel() const
Definition RooAbsReal.h:248
double _DEBUG_getVal(const RooArgSet *normalisationSet) const
Debug version of getVal(), which is slow and does error checking.
RooFit::OwningPtr< RooAbsArg > createFundamental(const char *newname=nullptr) const override
Create a RooRealVar fundamental object with our properties.
bool plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
RooFit::OwningPtr< 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).
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()
virtual RooFit::OwningPtr< RooFitResult > chi2FitTo(RooDataHist &data, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={})
Perform a fit to given histogram.
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 bool setData(RooAbsData &, bool=true)
Definition RooAbsReal.h:374
TString _unit
Unit for objects value.
Definition RooAbsReal.h:544
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
virtual void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
double getVal(const RooArgSet &normalisationSet) const
Like getVal(const RooArgSet*), but always requires an argument for normalisation.
Definition RooAbsReal.h:126
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.
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:211
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
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
Function that is called at the end of redirectServers().
virtual Int_t minTrialSamples(const RooArgSet &) const
Definition RooAbsReal.h:241
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:454
void setIntegratorConfig()
Remove the specialized numeric integration configuration associated with this object.
void syncCache(const RooArgSet *set=nullptr) override
Definition RooAbsReal.h:463
void printValue(std::ostream &os) const override
Print object value.
virtual bool forceAnalyticalInt(const RooAbsArg &) const
Definition RooAbsReal.h:164
bool isIdentical(const RooAbsArg &other, bool assumeSameType=false) const override
void setUnit(const char *unit)
Definition RooAbsReal.h:147
bool getForceNumInt() const
Definition RooAbsReal.h:174
virtual RooFit::OwningPtr< RooAbsReal > createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
static bool _hideOffset
Offset hiding flag.
Definition RooAbsReal.h:556
static ErrorLoggingMode _evalErrorMode
Definition RooAbsReal.h:552
void attachToVStore(RooVectorDataStore &vstore) override
static Int_t _evalErrorCount
Definition RooAbsReal.h:554
void copyCache(const RooAbsArg *source, bool valueOnly=false, bool setValDirty=true) override
Copy the cached value of another RooAbsArg to our cache.
virtual void gradient(double *) const
Definition RooAbsReal.h:395
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:543
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:347
friend class BatchInterfaceAccessor
Definition RooAbsReal.h:406
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.
double getPropagatedError(const RooFitResult &fr, const RooArgSet &nset={}) const
Propagates parameter uncertainties to an uncertainty estimate for this RooAbsReal.
static void setHideOffset(bool flag)
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:216
RooAbsMoment * moment(RooRealVar &obs, Int_t order, bool central, bool takeRoot)
Return function representing moment of function of given order.
virtual std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const
This function defines the analytical integral translation for the class.
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 ...
RooFit::UniqueId< RooArgSet >::Value_t _lastNormSetId
Component selection flag for RooAbsPdf::plotCompOn.
Definition RooAbsReal.h:550
const char * getPlotLabel() const
Get the label associated with the variable.
virtual void forceNumInt(bool flag=true)
Definition RooAbsReal.h:169
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:220
std::unique_ptr< RooNumIntConfig > _specIntegratorConfig
Definition RooAbsReal.h:547
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:541
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 string with unique suffix name to give to integral object that encodes integrated observabl...
virtual double evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
TString getTitle(bool appendUnit=false) const
Return this variable's title 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:143
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:353
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:369
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:568
Int_t _plotBins
Number of plot bins.
Definition RooAbsReal.h:542
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.
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const
Plot (project) PDF on specified frame.
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:540
RooAbsMoment * mean(RooRealVar &obs, const RooArgSet &nset)
Definition RooAbsReal.h:367
virtual RooFit::OwningPtr< RooAbsReal > createChi2(RooDataHist &data, const RooLinkedList &cmdList)
virtual double offset() const
Definition RooAbsReal.h:378
RooAbsMoment * sigma(RooRealVar &obs)
Definition RooAbsReal.h:368
static bool _globalSelectComp
Definition RooAbsReal.h:555
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:27
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
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
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
@ Extended
Definition RooCurve.h:39
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...
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
A class to maintain the context for squashing of RooFit models into code.
Lightweight interface adaptor that exports a RooAbsPdf as a functor.
Definition RooFunctor.h:25
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
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...
Implements a PDF constructed from a sum of functions:
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
The RooStringView is a wrapper around a C-style string that can also be constructed from a std::strin...
RooVectorDataStore uses std::vectors to store data columns.
1-Dim function class
Definition TF1.h:214
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:258
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
Namespace for dispatching RooFit computations to various backends.
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:484
const char * normRangeName
Definition RooAbsReal.h:480
RooFit::MPSplit interleave
Definition RooAbsReal.h:492
const char * projectionRangeName
Definition RooAbsReal.h:485
const RooArgSet * projDataSet
Definition RooAbsReal.h:479
const char * curveNameSuffix
Definition RooAbsReal.h:493
const char * addToCurveName
Definition RooAbsReal.h:488
const RooFitResult * errorFR
Definition RooAbsReal.h:498
const RooArgSet * projSet
Definition RooAbsReal.h:476
const char * curveName
Definition RooAbsReal.h:487
const RooAbsData * projData
Definition RooAbsReal.h:474
Option_t * drawOptions
Definition RooAbsReal.h:471
A UniqueId can be added as a class member to enhance any class with a unique identifier for each inst...
Definition UniqueId.h:39
unsigned long Value_t
Definition UniqueId.h:41
TMarker m
Definition textangle.C:8
static void output()