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