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