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 #include <map>
28 
29 class RooArgList ;
30 class RooDataSet ;
31 class RooPlot;
32 class RooRealVar;
33 class RooAbsFunc;
35 class RooCategory ;
36 class RooLinkedList ;
37 class RooNumIntConfig ;
38 class RooDataHist ;
39 class RooFunctor ;
40 class RooGenFunction ;
41 class RooMultiGenFunction ;
42 class RooFitResult ;
43 class RooAbsMoment ;
44 class RooDerivative ;
45 class RooVectorDataStore ;
46 namespace RooHelpers {
48 }
49 struct TreeReadBuffer; /// A space to attach TBranches
50 
51 class TH1;
52 class TH1F;
53 class TH2F;
54 class TH3F;
55 
56 #include <list>
57 #include <string>
58 #include <iostream>
59 #include <sstream>
60 
61 class RooAbsReal : public RooAbsArg {
62 public:
63  using value_type = double;
64 
65  // Constructors, assignment etc
66  RooAbsReal() ;
67  RooAbsReal(const char *name, const char *title, const char *unit= "") ;
68  RooAbsReal(const char *name, const char *title, Double_t minVal, Double_t maxVal,
69  const char *unit= "") ;
70  RooAbsReal(const RooAbsReal& other, const char* name=0);
71  RooAbsReal& operator=(const RooAbsReal& other);
72  virtual ~RooAbsReal();
73 
74 
75 
76 
77  //////////////////////////////////////////////////////////////////////////////////
78  /// Evaluate object. Returns either cached value or triggers a recalculation.
79  /// The recalculation happens by calling getValV(), which in the end calls the
80  /// virtual evaluate() functions of the respective PDFs.
81  /// \param[in] normalisationSet getValV() reacts differently depending on the value of the normalisation set.
82  /// If the set is `nullptr`, an unnormalised value is returned.
83  /// \note The normalisation is arbitrary, because it is up to the implementation
84  /// of the PDF to e.g. leave out normalisation constants for speed reasons. The range
85  /// of the variables is also ignored.
86  ///
87  /// To normalise the result properly, a RooArgSet has to be passed, which contains
88  /// the variables to normalise over.
89  /// These are integrated over their current ranges to compute the normalisation constant,
90  /// and the unnormalised result is divided by this value.
91  inline Double_t getVal(const RooArgSet* normalisationSet = nullptr) const {
92 #ifdef ROOFIT_CHECK_CACHED_VALUES
93  return _DEBUG_getVal(normalisationSet);
94 #else
95 
96 #ifndef _WIN32
97  return (_fast && !_inhibitDirty) ? _value : getValV(normalisationSet) ;
98 #else
99  return (_fast && !inhibitDirty()) ? _value : getValV(normalisationSet) ;
100 #endif
101 
102 #endif
103  }
104 
105  /// Like getVal(const RooArgSet*), but always requires an argument for normalisation.
106  inline Double_t getVal(const RooArgSet& normalisationSet) const { return _fast ? _value : getValV(&normalisationSet) ; }
107 
108  virtual Double_t getValV(const RooArgSet* normalisationSet = nullptr) const ;
109 
110  virtual RooSpan<const double> getValBatch(std::size_t begin, std::size_t maxSize, const RooArgSet* normSet = nullptr) const;
111 
112  Double_t getPropagatedError(const RooFitResult &fr, const RooArgSet &nset = RooArgSet()) const;
113 
114  Bool_t operator==(Double_t value) const ;
115  virtual Bool_t operator==(const RooAbsArg& other) const;
116  virtual Bool_t isIdentical(const RooAbsArg& other, Bool_t assumeSameType=kFALSE) const;
117 
118 
119  inline const Text_t *getUnit() const {
120  // Return string with unit description
121  return _unit.Data();
122  }
123  inline void setUnit(const char *unit) {
124  // Set unit description to given string
125  _unit= unit;
126  }
127  TString getTitle(Bool_t appendUnit= kFALSE) const;
128 
129  // Lightweight interface adaptors (caller takes ownership)
130  RooAbsFunc *bindVars(const RooArgSet &vars, const RooArgSet* nset=0, Bool_t clipInvalid=kFALSE) const;
131 
132  // Create a fundamental-type object that can hold our value.
133  RooAbsArg *createFundamental(const char* newname=0) const;
134 
135  // Analytical integration support
136  virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=0) const ;
137  virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
138  virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
139  virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
140  virtual Bool_t forceAnalyticalInt(const RooAbsArg& /*dep*/) const {
141  // Interface to force RooRealIntegral to offer given observable for internal integration
142  // even if this is deemed unsafe. This default implementation returns always flase
143  return kFALSE ;
144  }
145  virtual void forceNumInt(Bool_t flag=kTRUE) {
146  // If flag is true, all advertised analytical integrals will be ignored
147  // and all integrals are calculated numerically
148  _forceNumInt = flag ;
149  }
150  Bool_t getForceNumInt() const { return _forceNumInt ; }
151 
152  // Chi^2 fits to histograms
153  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
154  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
155  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
156  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
157 
158  virtual RooAbsReal* createChi2(RooDataHist& data, const RooLinkedList& cmdList) ;
159  virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
160  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
161  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
162 
163  // Chi^2 fits to X-Y datasets
164  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
165  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
166  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
167  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;
168 
169  virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
170  virtual RooAbsReal* createChi2(RooDataSet& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
171  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
172  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
173 
174 
175  virtual RooAbsReal* createProfile(const RooArgSet& paramsOfInterest) ;
176 
177 
178  RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
179  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
180  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
181  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
182 
183  /// Create integral over observables in iset in range named rangeName.
184  RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const {
185  return createIntegral(iset,0,0,rangeName) ;
186  }
187  /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
188  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=0) const {
189  return createIntegral(iset,&nset,0,rangeName) ;
190  }
191  /// Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
192  /// using specified configuration for any numeric integration.
193  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
194  return createIntegral(iset,&nset,&cfg,rangeName) ;
195  }
196  /// Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration.
197  RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
198  return createIntegral(iset,0,&cfg,rangeName) ;
199  }
200  virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset=0, const RooNumIntConfig* cfg=0, const char* rangeName=0) const ;
201 
202 
203  void setParameterizeIntegral(const RooArgSet& paramVars) ;
204 
205  // Create running integrals
206  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
207  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
208  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
209  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
210  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
211  RooAbsReal* createIntRI(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
212  RooAbsReal* createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
213 
214 
215  // Optimized accept/reject generator support
216  virtual Int_t getMaxVal(const RooArgSet& vars) const ;
217  virtual Double_t maxVal(Int_t code) const ;
218  virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }
219 
220 
221  // Plotting options
222  void setPlotLabel(const char *label);
223  const char *getPlotLabel() const;
224 
225  virtual Double_t defaultErrorLevel() const {
226  // Return default level for MINUIT error analysis
227  return 1.0 ;
228  }
229 
230  const RooNumIntConfig* getIntegratorConfig() const ;
235  void setIntegratorConfig() ;
236  void setIntegratorConfig(const RooNumIntConfig& config) ;
237 
238  virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),Bool_t force=kTRUE) ;
239  virtual void fixAddCoefRange(const char* rangeName=0,Bool_t force=kTRUE) ;
240 
241  virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;
242 
243  // User entry point for plotting
244  virtual RooPlot* plotOn(RooPlot* frame,
245  const RooCmdArg& arg1=RooCmdArg(), const RooCmdArg& arg2=RooCmdArg(),
246  const RooCmdArg& arg3=RooCmdArg(), const RooCmdArg& arg4=RooCmdArg(),
247  const RooCmdArg& arg5=RooCmdArg(), const RooCmdArg& arg6=RooCmdArg(),
248  const RooCmdArg& arg7=RooCmdArg(), const RooCmdArg& arg8=RooCmdArg(),
249  const RooCmdArg& arg9=RooCmdArg(), const RooCmdArg& arg10=RooCmdArg()
250  ) const ;
251 
252 
254 
255  // Forwarder function for backward compatibility
256  virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L",
257  Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=0) const;
258 
259  // Fill an existing histogram
260  TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
261  Double_t scaleFactor= 1, const RooArgSet *projectedVars= 0, Bool_t scaling=kTRUE,
262  const RooArgSet* condObs=0, Bool_t setError=kTRUE) const;
263 
264  // Create 1,2, and 3D histograms from and fill it
265  TH1 *createHistogram(const char* varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
266  TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
267  TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
268  const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
269  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
270  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
271  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
272 
273  // Fill a RooDataHist
274  RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, Double_t scaleFactor,
275  Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const ;
276 
277  // I/O streaming interface (machine readable)
278  virtual Bool_t readFromStream(std::istream& is, Bool_t compact, Bool_t verbose=kFALSE) ;
279  virtual void writeToStream(std::ostream& os, Bool_t compact) const ;
280 
281  // Printing interface (human readable)
282  virtual void printValue(std::ostream& os) const ;
283  virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
284 
285  static void setCacheCheck(Bool_t flag) ;
286 
287  // Evaluation error logging
288  class EvalError {
289  public:
290  EvalError() { }
291  EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
292  void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
293  void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
294  std::string _msg;
295  std::string _srvval;
296  } ;
297 
301  void logEvalError(const char* message, const char* serverValueString=0) const ;
302  static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=0) ;
303  static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
304  static Int_t numEvalErrors() ;
305  static Int_t numEvalErrorItems() ;
306 
307 
308  typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ;
309  static EvalErrorIter evalErrorIter() ;
310 
311  static void clearEvalErrorLog() ;
312 
313  virtual Bool_t isBinnedDistribution(const RooArgSet& /*obs*/) const { return kFALSE ; }
314  virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const { return 0 ; }
315  virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const {
316  // Interface for returning an optional hint for initial sampling points when constructing a curve
317  // projected on observable.
318  return 0 ;
319  }
320 
322  RooMultiGenFunction* iGenFunction(const RooArgSet& observables, const RooArgSet& nset=RooArgSet()) ;
323 
324  RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
325  TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
326 
327  RooDerivative* derivative(RooRealVar& obs, Int_t order=1, Double_t eps=0.001) ;
328  RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, Double_t eps=0.001) ;
329 
330  RooAbsMoment* moment(RooRealVar& obs, Int_t order, Bool_t central, Bool_t takeRoot) ;
331  RooAbsMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, Bool_t central, Bool_t takeRoot, Bool_t intNormObs) ;
332 
333  RooAbsMoment* mean(RooRealVar& obs) { return moment(obs,1,kFALSE,kFALSE) ; }
334  RooAbsMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,kFALSE,kFALSE,kTRUE) ; }
335  RooAbsMoment* sigma(RooRealVar& obs) { return moment(obs,2,kTRUE,kTRUE) ; }
336  RooAbsMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,kTRUE,kTRUE,kTRUE) ; }
337 
339 
340 
341  virtual Bool_t setData(RooAbsData& /*data*/, Bool_t /*cloneData*/=kTRUE) { return kTRUE ; }
342 
343  virtual void enableOffsetting(Bool_t) {} ;
344  virtual Bool_t isOffsetting() const { return kFALSE ; }
345  virtual Double_t offset() const { return 0 ; }
346 
347  static void setHideOffset(Bool_t flag);
348  static Bool_t hideOffset() ;
349 
350 protected:
351  // Hook for objects with normalization-dependent parameters interperetation
352  virtual void selectNormalization(const RooArgSet* depSet=0, Bool_t force=kFALSE) ;
353  virtual void selectNormalizationRange(const char* rangeName=0, Bool_t force=kFALSE) ;
354 
355  // Helper functions for plotting
356  Bool_t plotSanityChecks(RooPlot* frame) const ;
357  void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
358  RooArgSet& projectedVars, Bool_t silent) const ;
359 
360  TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=0, const char* rangeName=0, Bool_t omitEmpty=kFALSE) const ;
361 
362 
363  Bool_t isSelectedComp() const ;
364 
365 
366  public:
367  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars) const ;
368  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
369  const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
370  RooArgSet *&cloneSet, const char* rangeName=0, const RooArgSet* condObs=0) const;
371  protected:
372 
374 
375  void plotOnCompSelect(RooArgSet* selNodes) const ;
376  RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, Double_t Z, const RooArgSet* params, const RooLinkedList& argList, Bool_t method1) const ;
377 
378  // Support interface for subclasses to advertise their analytic integration
379  // and generator capabilities in their analticalIntegral() and generateEvent()
380  // implementations.
381  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
382  const RooArgProxy& a) const ;
383  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
384  const RooArgProxy& a, const RooArgProxy& b) const ;
385  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
386  const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
387  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
388  const RooArgProxy& a, const RooArgProxy& b,
389  const RooArgProxy& c, const RooArgProxy& d) const ;
390 
391  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
392  const RooArgSet& set) const ;
393 
394 
395  RooAbsReal* createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
396  void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;
397 
398 
399  // Internal consistency checking (needed by RooDataSet)
400  virtual Bool_t isValid() const ;
401  virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const ;
402 
403  // Function evaluation and error tracing
404  Double_t traceEval(const RooArgSet* set) const ;
405 
406  /// Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
407  virtual Double_t evaluate() const = 0;
408  virtual RooSpan<double> evaluateBatch(std::size_t begin, std::size_t maxSize) const;
409 
410  //---------- Interface to access batch data ---------------------------
411  //
414  _batchData.clear();
415  for (auto arg : _serverList) {
416  //TODO get rid of this cast?
417  auto absReal = dynamic_cast<RooAbsReal*>(arg);
418  if (absReal)
419  absReal->clearBatchMemory();
420  }
421  }
422 
423  private:
424  void checkBatchComputation(std::size_t evtNo, const RooArgSet* normSet = nullptr, double relAccuracy = 1.E-13) const;
425 
427  return _batchData;
428  }
429 
430  /// Debug version of getVal(), which is slow and does error checking.
431  Double_t _DEBUG_getVal(const RooArgSet* normalisationSet) const;
432 
433  //--------------------------------------------------------------------
434 
435  protected:
436  // Hooks for RooDataSet interface
437  friend class RooRealIntegral ;
438  friend class RooVectorDataStore ;
439  virtual void syncCache(const RooArgSet* set=0) { getVal(set) ; }
440  virtual void copyCache(const RooAbsArg* source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE) ;
441  virtual void attachToTree(TTree& t, Int_t bufSize=32000) ;
442  virtual void attachToVStore(RooVectorDataStore& vstore) ;
443  virtual void setTreeBranchStatus(TTree& t, Bool_t active) ;
444  virtual void fillTreeBranch(TTree& t) ;
445 
446  friend class RooRealBinding ;
447  Double_t _plotMin ; // Minimum of plot range
448  Double_t _plotMax ; // Maximum of plot range
449  Int_t _plotBins ; // Number of plot bins
450  mutable Double_t _value ; // Cache for current value of object
451  mutable BatchHelpers::BatchData _batchData; //! Value storage for batches of events
452  TString _unit ; // Unit for objects value
453  TString _label ; // Plot label for objects value
454  Bool_t _forceNumInt ; // Force numerical integration if flag set
455 
456  friend class RooAbsPdf ;
457  friend class RooAbsAnaConvPdf ;
458 
459  RooNumIntConfig* _specIntegratorConfig ; // Numeric integrator configuration specific for this object
460 
461  friend class RooDataProjBinding ;
462  friend class RooAbsOptGoodnessOfFit ;
463 
464  struct PlotOpt {
478  const char* normRangeName ;
483  const char* projectionRangeName ;
485  const char* curveName ;
486  const char* addToCurveName ;
491  const char* curveNameSuffix ;
497  } ;
498 
499  // Plot implementation functions
500  virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;
501 
502 public:
503  // PlotOn with command list
504  virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
505 
506  protected:
507  virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;
508 
509 
510 private:
511 
513  static std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > > _evalErrorList ;
515 
516  Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;
517 
518  std::unique_ptr<TreeReadBuffer> _treeReadBuffer; //! A buffer for reading values from trees
519 
520 protected:
521 
522 
523  friend class RooRealSumPdf ;
524  friend class RooRealSumFunc;
525  friend class RooAddPdf ;
526  friend class RooAddModel ;
527  void selectComp(Bool_t flag) {
528  // If flag is true, only selected component will be included in evaluates of RooAddPdf components
529  _selectComp = flag ;
530  }
531  static void globalSelectComp(Bool_t flag) ;
532  Bool_t _selectComp ; //! Component selection flag for RooAbsPdf::plotCompOn
533  static Bool_t _globalSelectComp ; // Global activation switch for component selection
534  // This struct can be used to flip the global switch to select components.
535  // Doing this with RAII prevents forgetting to reset the state.
539  if (state != RooAbsReal::_globalSelectComp)
541  }
542 
546  }
547 
548  bool _oldState;
549  };
550 
551 
552  mutable RooArgSet* _lastNSet ; //!
553  static Bool_t _hideOffset ; // Offset hiding flag
554 
555  ClassDef(RooAbsReal,2) // Abstract real-valued variable
556 };
557 
558 #endif
RooAbsReal::analyticalIntegral
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition: RooAbsReal.cxx:415
c
#define c(i)
Definition: RSha256.hxx:119
RooAbsReal::EvalError
Definition: RooAbsReal.h:288
RooAbsReal::setData
virtual Bool_t setData(RooAbsData &, Bool_t=kTRUE)
Definition: RooAbsReal.h:341
m
auto * m
Definition: textangle.C:8
RooAbsReal::copyCache
virtual void copyCache(const RooAbsArg *source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE)
Copy the cached value of another RooAbsArg to our cache.
Definition: RooAbsReal.cxx:3226
RooAbsReal::EvalError::_srvval
std::string _srvval
Definition: RooAbsReal.h:295
RooAbsReal::NumEvent
@ NumEvent
Definition: RooAbsReal.h:253
RooAbsReal::~RooAbsReal
virtual ~RooAbsReal()
Destructor.
Definition: RooAbsReal.cxx:208
RooAbsReal::CountErrors
@ CountErrors
Definition: RooAbsReal.h:298
RooMultiGenFunction
Definition: RooMultiGenFunction.h:27
RooAbsReal::selectNormalization
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.
Definition: RooAbsReal.cxx:3625
RooCmdArg
Definition: RooCmdArg.h:27
RooAbsReal::plotSamplingHint
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:315
RooAbsArg::_inhibitDirty
static Bool_t _inhibitDirty
Definition: RooAbsArg.h:648
RooAbsReal::PlotOpt::curveNameSuffix
const char * curveNameSuffix
Definition: RooAbsReal.h:491
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooAbsReal::printEvalErrors
static void printEvalErrors(std::ostream &os=std::cout, Int_t maxPerNode=10000000)
Print all outstanding logged evaluation error on the given ostream.
Definition: RooAbsReal.cxx:3825
e
#define e(i)
Definition: RSha256.hxx:121
RooAbsReal::PlotOpt::doeeval
Bool_t doeeval
Definition: RooAbsReal.h:494
RooAbsReal::setEvalErrorLoggingMode
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
Definition: RooAbsReal.cxx:4838
RooAddPdf
Definition: RooAddPdf.h:32
RooAbsReal::getPropagatedError
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...
Definition: RooAbsReal.cxx:2703
RooAbsReal::createIntegral
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:188
RooAbsReal::logEvalError
void logEvalError(const char *message, const char *serverValueString=0) const
Log evaluation error message.
Definition: RooAbsReal.cxx:3728
RooAbsData
Definition: RooAbsData.h:46
RooGenFunction
Definition: RooGenFunction.h:22
RooAbsReal::_evalErrorMode
static ErrorLoggingMode _evalErrorMode
Definition: RooAbsReal.h:512
RooDerivative
Definition: RooDerivative.h:31
RooAbsReal::GlobalSelectComponentRAII::GlobalSelectComponentRAII
GlobalSelectComponentRAII(bool state)
Definition: RooAbsReal.h:537
RooAbsReal::functor
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...
Definition: RooAbsReal.cxx:4148
RooAbsReal::syncCache
virtual void syncCache(const RooArgSet *set=0)
Definition: RooAbsReal.h:439
TH2F
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
RooArgSet.h
RooAbsReal::isBinnedDistribution
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:313
RooAbsReal::writeToStream
virtual void writeToStream(std::ostream &os, Bool_t compact) const
Write object contents to stream (dummy for now)
Definition: RooAbsReal.cxx:457
TString::Data
const char * Data() const
Definition: TString.h:369
RooAbsReal::createIntegral
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:197
RooAbsReal::sigma
RooAbsMoment * sigma(RooRealVar &obs)
Definition: RooAbsReal.h:335
RooAbsReal::createIntegral
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:561
RooAbsReal::getValV
virtual Double_t getValV(const RooArgSet *normalisationSet=nullptr) const
Return value of object.
Definition: RooAbsReal.cxx:271
RooAbsReal::operator==
Bool_t operator==(Double_t value) const
Equality operator comparing to a Double_t.
Definition: RooAbsReal.cxx:218
xmax
float xmax
Definition: THbookFile.cxx:95
RooAbsReal::PlotOpt::projectionRangeName
const char * projectionRangeName
Definition: RooAbsReal.h:483
RooAbsArg::_serverList
RefCountList_t _serverList
Definition: RooAbsArg.h:587
RooAbsReal::CollectErrors
@ CollectErrors
Definition: RooAbsReal.h:298
RooAbsReal::getPlotLabel
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:427
RooArgList
Definition: RooArgList.h:21
RooAbsReal::attachToTree
virtual void attachToTree(TTree &t, Int_t bufSize=32000)
Attach object to a branch of given TTree.
Definition: RooAbsReal.cxx:3285
TTree
Definition: TTree.h:79
RooAbsReal::defaultErrorLevel
virtual Double_t defaultErrorLevel() const
Definition: RooAbsReal.h:225
RooAbsReal::EvalError::setServerValues
void setServerValues(const char *tmp)
Definition: RooAbsReal.h:293
RooAbsReal::PlotOpt::numCPU
Int_t numCPU
Definition: RooAbsReal.h:489
RooAbsReal::plotSanityChecks
Bool_t plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
Definition: RooAbsReal.cxx:3075
RooAbsReal::getVal
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:91
RooAbsReal::PlotOpt::projDataSet
const RooArgSet * projDataSet
Definition: RooAbsReal.h:477
TGeant4Unit::s
static constexpr double s
Definition: TGeant4SystemOfUnits.h:168
RooAbsArg.h
RooAbsReal::_unit
TString _unit
Value storage for batches of events.
Definition: RooAbsReal.h:452
RooFit::MPSplit
MPSplit
Definition: RooGlobalFunc.h:70
RooAbsReal::isSelectedComp
Bool_t isSelectedComp() const
If true, the current pdf is a selected component (for use in plotting)
Definition: RooAbsReal.cxx:3178
RooAbsReal::makeProjectionSet
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...
Definition: RooAbsReal.cxx:3118
RooAbsReal::getVal
Double_t getVal(const RooArgSet &normalisationSet) const
Like getVal(const RooArgSet*), but always requires an argument for normalisation.
Definition: RooAbsReal.h:106
RooCurve.h
x
Double_t x[n]
Definition: legend1.C:17
RooAbsReal::_forceNumInt
Bool_t _forceNumInt
Definition: RooAbsReal.h:454
RooAbsReal::setPlotLabel
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:437
RooAbsReal::_treeReadBuffer
std::unique_ptr< TreeReadBuffer > _treeReadBuffer
Definition: RooAbsReal.h:518
indent
static void indent(ostringstream &buf, int indent_level)
Definition: TClingCallFunc.cxx:87
RooHelpers::BatchInterfaceAccessor
Helper class to access a batch-related part of RooAbsReal's interface, which should not leak to the o...
Definition: RooHelpers.h:173
RooAbsReal
Definition: RooAbsReal.h:61
RooAbsReal::evalErrorLoggingMode
static ErrorLoggingMode evalErrorLoggingMode()
Return current evaluation error logging mode.
Definition: RooAbsReal.cxx:4824
RooAbsMoment
Definition: RooAbsMoment.h:27
RooAbsReal::PlotOpt::curveInvisible
Bool_t curveInvisible
Definition: RooAbsReal.h:484
RooAbsReal::isValidReal
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:503
RooAbsArg::inhibitDirty
Bool_t inhibitDirty() const
Delete watch flag.
Definition: RooAbsArg.cxx:83
RooAbsReal::isOffsetting
virtual Bool_t isOffsetting() const
Definition: RooAbsReal.h:344
RooAbsReal::plotOnWithErrorBand
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 ...
Definition: RooAbsReal.cxx:2828
RooAbsReal::selectNormalizationRange
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.
Definition: RooAbsReal.cxx:3638
RooAbsReal::mean
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:333
TString
Definition: TString.h:136
RooAbsReal::clearBatchMemory
void clearBatchMemory()
Definition: RooAbsReal.h:413
RooAbsReal::printValue
virtual void printValue(std::ostream &os) const
Print object value.
Definition: RooAbsReal.cxx:466
RooAbsArg::_fast
Bool_t _fast
Definition: RooAbsArg.h:664
RooAbsReal::_globalSelectComp
static Bool_t _globalSelectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition: RooAbsReal.h:533
RooCmdArg::none
static const RooCmdArg & none()
Return reference to null argument.
Definition: RooCmdArg.cxx:52
RooAbsReal::fixAddCoefNormalization
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...
Definition: RooAbsReal.cxx:3893
RooAbsReal::PlotOpt::projData
const RooAbsData * projData
Definition: RooAbsReal.h:472
b
#define b(i)
Definition: RSha256.hxx:118
RooAbsReal::clearEvalErrorLog
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
Definition: RooAbsReal.cxx:3805
RooAbsReal::RooAbsReal
RooAbsReal()
coverity[UNINIT_CTOR] Default constructor
Definition: RooAbsReal.cxx:129
bool
RooFitResult
Definition: RooFitResult.h:40
RooAbsReal::ScaleType
ScaleType
Definition: RooAbsReal.h:253
RooAbsReal::PlotOpt::wmode
RooCurve::WingMode wmode
Definition: RooAbsReal.h:482
RooAbsReal::PlotOpt::scaleFactor
Double_t scaleFactor
Definition: RooAbsReal.h:470
RooAbsFunc
Definition: RooAbsFunc.h:23
RooAbsReal::GlobalSelectComponentRAII
Definition: RooAbsReal.h:536
RooAbsReal::numEvalErrorItems
static Int_t numEvalErrorItems()
Definition: RooAbsReal.cxx:323
RooAbsReal::getValBatch
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:295
RooAbsReal::fillTreeBranch
virtual void fillTreeBranch(TTree &t)
Fill the tree branch that associated with this object with its current value.
Definition: RooAbsReal.cxx:3353
RooAbsReal::PlotOpt::addToWgtOther
Double_t addToWgtOther
Definition: RooAbsReal.h:488
RooAbsReal::createIntegral
RooAbsReal * createIntegral(const RooArgSet &iset, const char *rangeName) const
Create integral over observables in iset in range named rangeName.
Definition: RooAbsReal.h:184
RooAbsReal::EvalError::setMessage
void setMessage(const char *tmp)
Definition: RooAbsReal.h:292
RooRealIntegral
Definition: RooRealIntegral.h:34
RooAbsReal::_selectComp
Bool_t _selectComp
Definition: RooAbsReal.h:532
RooAbsReal::forceAnalyticalInt
virtual Bool_t forceAnalyticalInt(const RooAbsArg &) const
Definition: RooAbsReal.h:140
RooAbsReal::value_type
double value_type
Definition: RooAbsReal.h:63
RooAbsReal::_plotMax
Double_t _plotMax
Definition: RooAbsReal.h:448
RooAbsReal::PlotOpt::numee
Int_t numee
Definition: RooAbsReal.h:492
RooAbsReal::moment
RooAbsMoment * moment(RooRealVar &obs, Int_t order, Bool_t central, Bool_t takeRoot)
Return function representing moment of function of given order.
Definition: RooAbsReal.cxx:4280
RooAbsReal::createIntObj
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:632
RooAbsReal::_evalErrorList
static std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > > _evalErrorList
Definition: RooAbsReal.h:513
RooAbsReal::plotAsymOn
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
Definition: RooAbsReal.cxx:2412
RooAbsReal::PlotOpt::addToCurveName
const char * addToCurveName
Definition: RooAbsReal.h:486
RooAbsReal::_label
TString _label
Definition: RooAbsReal.h:453
RooAbsReal::PlotOpt
Definition: RooAbsReal.h:464
RooAbsReal::isValid
virtual Bool_t isValid() const
Check if current value is valid.
Definition: RooAbsReal.cxx:492
RooDataHist
Definition: RooDataHist.h:39
RooRealBinding
Definition: RooRealBinding.h:26
Option_t
const typedef char Option_t
Definition: RtypesCore.h:66
RooAbsReal::minTrialSamples
virtual Int_t minTrialSamples(const RooArgSet &) const
Definition: RooAbsReal.h:218
RooAbsReal::PlotOpt::PlotOpt
PlotOpt()
Definition: RooAbsReal.h:465
RooAbsReal::EvalError::_msg
std::string _msg
Definition: RooAbsReal.h:294
BatchHelpers::BatchData
A class to store batches of data points that can be accessed via RooSpan.
Definition: BatchData.h:32
RooAbsReal::ErrorLoggingMode
ErrorLoggingMode
Definition: RooAbsReal.h:298
RooAbsReal::selectComp
void selectComp(Bool_t flag)
Definition: RooAbsReal.h:527
RooAbsReal::GlobalSelectComponentRAII::~GlobalSelectComponentRAII
~GlobalSelectComponentRAII()
Definition: RooAbsReal.h:543
RooRealSumFunc
Definition: RooRealSumFunc.h:25
RooAbsReal::readFromStream
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:447
RooAbsReal::PlotOpt::rangeLo
Double_t rangeLo
Definition: RooAbsReal.h:479
xmin
float xmin
Definition: THbookFile.cxx:95
RooAbsReal::batchData
const BatchHelpers::BatchData & batchData() const
Definition: RooAbsReal.h:426
RooCurve::WingMode
WingMode
Definition: RooCurve.h:35
RooArgProxy
Definition: RooArgProxy.h:24
RooAbsReal::matchArgsByName
Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const
Check if allArgs contains matching elements for each name in nameList.
Definition: RooAbsReal.cxx:3498
RooFit::Interleave
@ Interleave
Definition: RooGlobalFunc.h:70
RooFit
Definition: RooCFunction1Binding.h:29
RooAbsReal::setTreeBranchStatus
virtual void setTreeBranchStatus(TTree &t, Bool_t active)
(De)Activate associated tree branch
Definition: RooAbsReal.cxx:3370
RooAbsReal::asTF
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...
Definition: RooAbsReal.cxx:4174
RooAbsReal::evalErrorIter
static EvalErrorIter evalErrorIter()
Definition: RooAbsReal.cxx:331
RooAbsReal::checkBatchComputation
void checkBatchComputation(std::size_t evtNo, const RooArgSet *normSet=nullptr, double relAccuracy=1.E-13) const
Definition: RooAbsReal.cxx:4979
RooAbsReal::binBoundaries
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:314
RooAbsReal::setHideOffset
static void setHideOffset(Bool_t flag)
Definition: RooAbsReal.cxx:117
RooAbsReal::integralNameSuffix
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:808
RooAbsReal::getUnit
const Text_t * getUnit() const
Definition: RooAbsReal.h:119
RooAbsReal::analyticalIntegralWN
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:401
a
auto * a
Definition: textangle.C:12
RooAbsReal::PlotOpt::drawOptions
Option_t * drawOptions
Definition: RooAbsReal.h:468
RooAbsReal::createChi2
virtual RooAbsReal * createChi2(RooDataHist &data, const RooLinkedList &cmdList)
Create a variable from a histogram and this function.
Definition: RooAbsReal.cxx:4449
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:92
RooAbsReal::derivative
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.
Definition: RooAbsReal.cxx:4252
RooAbsReal::bindVars
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).
Definition: RooAbsReal.cxx:3202
RooAbsReal::getTitle
TString getTitle(Bool_t appendUnit=kFALSE) const
Return this variable's title string.
Definition: RooAbsReal.cxx:253
RooLinkedList
Definition: RooLinkedList.h:35
RooAbsReal::EvalError::EvalError
EvalError(const EvalError &other)
Definition: RooAbsReal.h:291
RooAbsReal::evaluate
virtual Double_t evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsReal::fillDataHist
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.
Definition: RooAbsReal.cxx:1229
RooAbsReal::chi2FitDriver
RooFitResult * chi2FitDriver(RooAbsReal &fcn, RooLinkedList &cmdList)
Internal driver function for chi2 fits.
Definition: RooAbsReal.cxx:4609
RooAbsReal::PlotOpt::stype
ScaleType stype
Definition: RooAbsReal.h:471
RooAbsReal::Ignore
@ Ignore
Definition: RooAbsReal.h:298
RooAbsReal::createIntRI
RooAbsReal * createIntRI(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Utility function for createRunningIntegral.
Definition: RooAbsReal.cxx:4082
RooAbsReal::getMaxVal
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
Definition: RooAbsReal.cxx:3651
RooAbsReal::_plotBins
Int_t _plotBins
Definition: RooAbsReal.h:449
RooAbsReal::setParameterizeIntegral
void setParameterizeIntegral(const RooArgSet &paramVars)
Definition: RooAbsReal.cxx:4846
RooPlot
Definition: RooPlot.h:44
RooAbsArg::RooArgSet
friend class RooArgSet
Definition: RooAbsArg.h:579
RooAbsReal::mean
RooAbsMoment * mean(RooRealVar &obs, const RooArgSet &nset)
Definition: RooAbsReal.h:334
RooAbsReal::getIntegratorConfig
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
Definition: RooAbsReal.cxx:3568
RooAbsReal::getAnalyticalIntegral
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:389
RooCurve
Definition: RooCurve.h:32
double
double
Definition: Converters.cxx:921
RooAbsReal::chi2FitTo
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.
Definition: RooAbsReal.cxx:4381
RooAbsReal::_value
Double_t _value
Definition: RooAbsReal.h:450
RooAbsReal::enableOffsetting
virtual void enableOffsetting(Bool_t)
Definition: RooAbsReal.h:343
RooAbsReal::RooAbsOptGoodnessOfFit
friend class RooAbsOptGoodnessOfFit
Definition: RooAbsReal.h:462
RooAbsReal::plotOnCompSelect
void plotOnCompSelect(RooArgSet *selNodes) const
Helper function for plotting of composite p.d.fs.
Definition: RooAbsReal.cxx:1550
RooAbsReal::forceNumInt
virtual void forceNumInt(Bool_t flag=kTRUE)
Definition: RooAbsReal.h:145
RooAbsReal::offset
virtual Double_t offset() const
Definition: RooAbsReal.h:345
RooAbsReal::createRunningIntegral
RooAbsReal * createRunningIntegral(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Calls createRunningIntegral(const RooArgSet&, const RooCmdArg&, const RooCmdArg&, const RooCmdArg&,...
Definition: RooAbsReal.cxx:3959
RooAbsReal::findRoot
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.
Definition: RooAbsReal.cxx:4316
RooAbsReal::_plotMin
Double_t _plotMin
Definition: RooAbsReal.h:447
RooCmdArg.h
RooAbsReal::setIntegratorConfig
void setIntegratorConfig()
Remove the specialized numeric integration configuration associated with this object.
Definition: RooAbsReal.cxx:3608
RooAbsReal::createIntegral
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:193
RooAbsReal::_batchData
BatchHelpers::BatchData _batchData
Definition: RooAbsReal.h:451
RooAbsReal::plotOn
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.
Definition: RooAbsReal.cxx:1714
TH3F
3-D histogram with a float per channel (see TH1 documentation)}
Definition: TH3.h:267
RooGlobalFunc.h
RooAbsReal::createProfile
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
Definition: RooAbsReal.cxx:516
RooAbsReal::_specIntegratorConfig
RooNumIntConfig * _specIntegratorConfig
Definition: RooAbsReal.h:459
RooAbsReal::maxVal
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Definition: RooAbsReal.cxx:3662
RooAbsReal::Raw
@ Raw
Definition: RooAbsReal.h:253
RooAbsReal::PlotOpt::normRangeName
const char * normRangeName
Definition: RooAbsReal.h:478
RooAbsReal::findInnerMostIntegration
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:753
RooAbsReal::setCacheCheck
static void setCacheCheck(Bool_t flag)
RooAbsReal::traceEval
Double_t traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
Definition: RooAbsReal.cxx:340
RooAbsReal::defaultIntegratorConfig
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Definition: RooAbsReal.cxx:3531
RooAbsReal::createHistogram
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...
Definition: RooAbsReal.cxx:1291
RooAbsReal::getAnalyticalIntegralWN
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:372
Double_t
double Double_t
Definition: RtypesCore.h:59
RooFunctor
Definition: RooFunctor.h:25
RooAbsReal::specialIntegratorConfig
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Definition: RooAbsReal.cxx:3541
RooAbsReal::setUnit
void setUnit(const char *unit)
Definition: RooAbsReal.h:123
RooAbsReal::operator=
RooAbsReal & operator=(const RooAbsReal &other)
Assign values, name and configs from another RooAbsReal.
Definition: RooAbsReal.cxx:181
RooAbsReal::PlotOpt::curveName
const char * curveName
Definition: RooAbsReal.h:485
RooAbsReal::EvalErrorIter
std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > >::const_iterator EvalErrorIter
Definition: RooAbsReal.h:308
RooDataProjBinding
Definition: RooDataProjBinding.h:25
RooAbsReal::EvalError::EvalError
EvalError()
Definition: RooAbsReal.h:290
RooAbsReal::_evalErrorCount
static Int_t _evalErrorCount
Definition: RooAbsReal.h:514
RooAbsReal::evaluateBatch
virtual RooSpan< double > evaluateBatch(std::size_t begin, std::size_t maxSize) const
Evaluate function for a batch of input data points.
Definition: RooAbsReal.cxx:4872
RooAbsReal::iGenFunction
RooGenFunction * iGenFunction(RooRealVar &x, const RooArgSet &nset=RooArgSet())
Definition: RooAbsReal.cxx:4328
RooAbsReal::PlotOpt::postRangeFracScale
Bool_t postRangeFracScale
Definition: RooAbsReal.h:481
RooCategory
Definition: RooCategory.h:27
RooAbsReal::PlotOpt::precision
Double_t precision
Definition: RooAbsReal.h:475
TH1F
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:572
RooAddModel
Definition: RooAddModel.h:28
RooAbsReal::GlobalSelectComponentRAII::_oldState
bool _oldState
Definition: RooAbsReal.h:548
RooVectorDataStore
Definition: RooVectorDataStore.h:35
RooAbsReal::isIdentical
virtual Bool_t isIdentical(const RooAbsArg &other, Bool_t assumeSameType=kFALSE) const
Definition: RooAbsReal.cxx:238
RooAbsReal::createScanRI
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...
Definition: RooAbsReal.cxx:4065
ClassDef
#define ClassDef(name, id)
Definition: Rtypes.h:325
RooAbsReal::PlotOpt::errorFR
const RooFitResult * errorFR
Definition: RooAbsReal.h:496
RooAbsReal::getForceNumInt
Bool_t getForceNumInt() const
Definition: RooAbsReal.h:150
TH1
Definition: TH1.h:57
RooAbsReal::createFundamental
RooAbsArg * createFundamental(const char *newname=0) const
Create a RooRealVar fundamental object with our properties.
Definition: RooAbsReal.cxx:3384
RooAbsReal::plotSliceOn
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
Definition: RooAbsReal.cxx:2378
RooAbsReal::createPlotProjection
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:887
name
char name[80]
Definition: TGX11.cxx:110
RooAbsReal::hideOffset
static Bool_t hideOffset()
Definition: RooAbsReal.cxx:118
BatchHelpers::BatchData::clear
void clear()
Discard all storage.
Definition: BatchData.h:47
BatchData.h
d
#define d(i)
Definition: RSha256.hxx:120
RooAbsReal::printMultiline
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Definition: RooAbsReal.cxx:476
genreflex::verbose
bool verbose
Definition: rootcling_impl.cxx:133
RooDataSet
Definition: RooDataSet.h:33
RooAbsReal::preferredObservableScanOrder
virtual void preferredObservableScanOrder(const RooArgSet &obs, RooArgSet &orderedObs) const
Interface method for function objects to indicate their preferred order of observables for scanning t...
Definition: RooAbsReal.cxx:3947
RooAbsArg
Definition: RooAbsArg.h:73
RooAbsReal::PlotOpt::addToWgtSelf
Double_t addToWgtSelf
Definition: RooAbsReal.h:487
RooNumIntConfig
Definition: RooNumIntConfig.h:25
RooAbsReal::_DEBUG_getVal
Double_t _DEBUG_getVal(const RooArgSet *normalisationSet) const
Debug version of getVal(), which is slow and does error checking.
Definition: RooAbsReal.cxx:4941
RooAbsReal::matchArgs
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Definition: RooAbsReal.cxx:3400
RooAbsReal::globalSelectComp
static void globalSelectComp(Bool_t flag)
Global switch controlling the activation of the selectComp() functionality.
Definition: RooAbsReal.cxx:3188
RooAbsPdf
Definition: RooAbsPdf.h:40
RooAbsReal::PlotOpt::interleave
RooFit::MPSplit interleave
Definition: RooAbsReal.h:490
RooAbsReal::PlotOpt::eeval
Double_t eeval
Definition: RooAbsReal.h:493
RooAbsReal::fixAddCoefRange
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...
Definition: RooAbsReal.cxx:3923
RooSpan.h
RooAbsReal::sigma
RooAbsMoment * sigma(RooRealVar &obs, const RooArgSet &nset)
Definition: RooAbsReal.h:336
RooAbsReal::Relative
@ Relative
Definition: RooAbsReal.h:253
RooAbsReal::numEvalErrors
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
Definition: RooAbsReal.cxx:3867
Text_t
char Text_t
Definition: RtypesCore.h:62
RooHelpers
Definition: RooAbsReal.h:46
TF1
1-Dim function class
Definition: TF1.h:212
RooAbsReal::PrintErrors
@ PrintErrors
Definition: RooAbsReal.h:298
RooAbsReal::PlotOpt::binProjData
Bool_t binProjData
Definition: RooAbsReal.h:473
RooAbsCategoryLValue
Definition: RooAbsCategoryLValue.h:25
RooRealVar
Definition: RooRealVar.h:36
RooAbsRealLValue
Definition: RooAbsRealLValue.h:31
RooArgList.h
RooRealSumPdf
Definition: RooRealSumPdf.h:25
TMath::E
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:102
RooAbsReal::PlotOpt::projSet
const RooArgSet * projSet
Definition: RooAbsReal.h:474
RooAbsReal::fillHistogram
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.
Definition: RooAbsReal.cxx:1052
RooAbsReal::_hideOffset
static Bool_t _hideOffset
Definition: RooAbsReal.h:553
RooAbsAnaConvPdf
Definition: RooAbsAnaConvPdf.h:34
TList
Definition: TList.h:44
RooSpan
A simple container to hold a batch of data values.
Definition: RooSpan.h:33
RooAbsReal::PlotOpt::shiftToZero
Bool_t shiftToZero
Definition: RooAbsReal.h:476
RooAbsReal::PlotOpt::progress
Bool_t progress
Definition: RooAbsReal.h:495
RooArgSet
Definition: RooArgSet.h:28
RooAbsReal::_lastNSet
RooArgSet * _lastNSet
Definition: RooAbsReal.h:552
RooAbsReal::RelativeExpected
@ RelativeExpected
Definition: RooAbsReal.h:253
int
RooAbsReal::PlotOpt::rangeHi
Double_t rangeHi
Definition: RooAbsReal.h:480
RooAbsReal::attachToVStore
virtual void attachToVStore(RooVectorDataStore &vstore)
Definition: RooAbsReal.cxx:3242
RooFit::Extended
RooCmdArg Extended(Bool_t flag=kTRUE)
Definition: RooGlobalFunc.cxx:155