ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
26 class RooArgList ;
27 class RooDataSet ;
28 class RooPlot;
29 class RooRealVar;
30 class RooAbsFunc;
32 class RooCategory ;
33 class RooLinkedList ;
34 class RooNumIntConfig ;
35 class RooDataHist ;
36 class RooFunctor ;
37 class RooGenFunction ;
38 class RooMultiGenFunction ;
39 class RooFitResult ;
40 class RooAbsMoment ;
41 class RooDerivative ;
42 class RooVectorDataStore ;
43 
44 class TH1;
45 class TH1F;
46 class TH2F;
47 class TH3F;
48 
49 #include <list>
50 #include <string>
51 #include <iostream>
52 
53 class RooAbsReal : public RooAbsArg {
54 public:
55  // Constructors, assignment etc
56  RooAbsReal() ;
57  RooAbsReal(const char *name, const char *title, const char *unit= "") ;
58  RooAbsReal(const char *name, const char *title, Double_t minVal, Double_t maxVal,
59  const char *unit= "") ;
60  RooAbsReal(const RooAbsReal& other, const char* name=0);
61  virtual ~RooAbsReal();
62 
63  // Return value and unit accessors
64  inline Double_t getVal(const RooArgSet* set=0) const {
65 /* if (_fast && !_inhibitDirty && std::string("RooHistFunc")==IsA()->GetName()) std::cout << "RooAbsReal::getVal(" << GetName() << ") CLEAN value = " << _value << std::endl ; */
66 #ifndef _WIN32
67  return (_fast && !_inhibitDirty) ? _value : getValV(set) ;
68 #else
69  return (_fast && !inhibitDirty()) ? _value : getValV(set) ;
70 #endif
71  }
72  inline Double_t getVal(const RooArgSet& set) const { return _fast ? _value : getValV(&set) ; }
73 
74  virtual Double_t getValV(const RooArgSet* set=0) const ;
75 
77 
79  virtual Bool_t operator==(const RooAbsArg& other) ;
80  virtual Bool_t isIdentical(const RooAbsArg& other, Bool_t assumeSameType=kFALSE) ;
81 
82 
83  inline const Text_t *getUnit() const {
84  // Return string with unit description
85  return _unit.Data();
86  }
87  inline void setUnit(const char *unit) {
88  // Set unit description to given string
89  _unit= unit;
90  }
91  TString getTitle(Bool_t appendUnit= kFALSE) const;
92 
93  // Lightweight interface adaptors (caller takes ownership)
94  RooAbsFunc *bindVars(const RooArgSet &vars, const RooArgSet* nset=0, Bool_t clipInvalid=kFALSE) const;
95 
96  // Create a fundamental-type object that can hold our value.
97  RooAbsArg *createFundamental(const char* newname=0) const;
98 
99  // Analytical integration support
100  virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=0) const ;
101  virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
102  virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
103  virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
104  virtual Bool_t forceAnalyticalInt(const RooAbsArg& /*dep*/) const {
105  // Interface to force RooRealIntegral to offer given observable for internal integration
106  // even if this is deemed unsafe. This default implementation returns always flase
107  return kFALSE ;
108  }
109  virtual void forceNumInt(Bool_t flag=kTRUE) {
110  // If flag is true, all advertised analytical integrals will be ignored
111  // and all integrals are calculated numerically
112  _forceNumInt = flag ;
113  }
114  Bool_t getForceNumInt() const { return _forceNumInt ; }
115 
116  // Chi^2 fits to histograms
117  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
118  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
119  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
120  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;
121 
122  virtual RooAbsReal* createChi2(RooDataHist& data, const RooLinkedList& cmdList) ;
123  virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
124  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
125  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
126 
127  // Chi^2 fits to X-Y datasets
128  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
129  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
130  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
131  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;
132 
133  virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
134  virtual RooAbsReal* createChi2(RooDataSet& data, const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
135  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),
136  const RooCmdArg& arg6=RooCmdArg::none(), const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
137 
138 
139  virtual RooAbsReal* createProfile(const RooArgSet& paramsOfInterest) ;
140 
141 
142  RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
143  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
144  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
145  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
146 
147  RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const {
148  // Create integral over observables in iset in range named rangeName
149  return createIntegral(iset,0,0,rangeName) ;
150  }
151  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=0) const {
152  // Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
153  return createIntegral(iset,&nset,0,rangeName) ;
154  }
155  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
156  // Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
157  // using specified configuration for any numeric integration
158  return createIntegral(iset,&nset,&cfg,rangeName) ;
159  }
160  RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=0) const {
161  // Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration
162  return createIntegral(iset,0,&cfg,rangeName) ;
163  }
164  virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset=0, const RooNumIntConfig* cfg=0, const char* rangeName=0) const ;
165 
166 
167  void setParameterizeIntegral(const RooArgSet& paramVars) ;
168 
169  // Create running integrals
170  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
171  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
172  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
173  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
174  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
175  RooAbsReal* createIntRI(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
176  RooAbsReal* createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;
177 
178 
179  // Optimized accept/reject generator support
180  virtual Int_t getMaxVal(const RooArgSet& vars) const ;
181  virtual Double_t maxVal(Int_t code) const ;
182  virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }
183 
184 
185  // Plotting options
186  void setPlotLabel(const char *label);
187  const char *getPlotLabel() const;
188 
189  virtual Double_t defaultErrorLevel() const {
190  // Return default level for MINUIT error analysis
191  return 1.0 ;
192  }
193 
194  const RooNumIntConfig* getIntegratorConfig() const ;
199  void setIntegratorConfig() ;
200  void setIntegratorConfig(const RooNumIntConfig& config) ;
201 
202  virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),Bool_t force=kTRUE) ;
203  virtual void fixAddCoefRange(const char* rangeName=0,Bool_t force=kTRUE) ;
204 
205  virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;
206 
207  // User entry point for plotting
208  virtual RooPlot* plotOn(RooPlot* frame,
209  const RooCmdArg& arg1=RooCmdArg(), const RooCmdArg& arg2=RooCmdArg(),
210  const RooCmdArg& arg3=RooCmdArg(), const RooCmdArg& arg4=RooCmdArg(),
211  const RooCmdArg& arg5=RooCmdArg(), const RooCmdArg& arg6=RooCmdArg(),
212  const RooCmdArg& arg7=RooCmdArg(), const RooCmdArg& arg8=RooCmdArg(),
213  const RooCmdArg& arg9=RooCmdArg(), const RooCmdArg& arg10=RooCmdArg()
214  ) const ;
215 
216 
218 
219  // Forwarder function for backward compatibility
220  virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L",
221  Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=0) const;
222 
223  // Fill an existing histogram
224  TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
225  Double_t scaleFactor= 1, const RooArgSet *projectedVars= 0, Bool_t scaling=kTRUE,
226  const RooArgSet* condObs=0, Bool_t setError=kTRUE) const;
227 
228  // Create 1,2, and 3D histograms from and fill it
229  TH1 *createHistogram(const char* varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
230  TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
231  TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
232  const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(),
233  const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(),
234  const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(),
235  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;
236 
237  // Fill a RooDataHist
238  RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, Double_t scaleFactor,
239  Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const ;
240 
241  // I/O streaming interface (machine readable)
242  virtual Bool_t readFromStream(std::istream& is, Bool_t compact, Bool_t verbose=kFALSE) ;
243  virtual void writeToStream(std::ostream& os, Bool_t compact) const ;
244 
245  // Printing interface (human readable)
246  virtual void printValue(std::ostream& os) const ;
247  virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;
248 
249  static void setCacheCheck(Bool_t flag) ;
250 
251  // Evaluation error logging
252  class EvalError {
253  public:
254  EvalError() { }
255  EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
256  void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
257  void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
258  std::string _msg;
259  std::string _srvval;
260  } ;
261 
265  void logEvalError(const char* message, const char* serverValueString=0) const ;
266  static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=0) ;
267  static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
268  static Int_t numEvalErrors() ;
269  static Int_t numEvalErrorItems() ;
270 
271 
272  typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ;
273  static EvalErrorIter evalErrorIter() ;
274 
275  static void clearEvalErrorLog() ;
276 
277  virtual Bool_t isBinnedDistribution(const RooArgSet& /*obs*/) const { return kFALSE ; }
278  virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const { return 0 ; }
279  virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const {
280  // Interface for returning an optional hint for initial sampling points when constructing a curve
281  // projected on observable.
282  return 0 ;
283  }
284 
286  RooMultiGenFunction* iGenFunction(const RooArgSet& observables, const RooArgSet& nset=RooArgSet()) ;
287 
288  RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
289  TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
290 
291  RooDerivative* derivative(RooRealVar& obs, Int_t order=1, Double_t eps=0.001) ;
292  RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, Double_t eps=0.001) ;
293 
294  RooAbsMoment* moment(RooRealVar& obs, Int_t order, Bool_t central, Bool_t takeRoot) ;
295  RooAbsMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, Bool_t central, Bool_t takeRoot, Bool_t intNormObs) ;
296 
297  RooAbsMoment* mean(RooRealVar& obs) { return moment(obs,1,kFALSE,kFALSE) ; }
298  RooAbsMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,kFALSE,kFALSE,kTRUE) ; }
299  RooAbsMoment* sigma(RooRealVar& obs) { return moment(obs,2,kTRUE,kTRUE) ; }
300  RooAbsMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,kTRUE,kTRUE,kTRUE) ; }
301 
303 
304 
305  virtual Bool_t setData(RooAbsData& /*data*/, Bool_t /*cloneData*/=kTRUE) { return kTRUE ; }
306 
307  virtual void enableOffsetting(Bool_t) {} ;
308  virtual Bool_t isOffsetting() const { return kFALSE ; }
309  virtual Double_t offset() const { return 0 ; }
310 
311  static void setHideOffset(Bool_t flag);
312  static Bool_t hideOffset() ;
313 
314 protected:
315  // Hook for objects with normalization-dependent parameters interperetation
316  virtual void selectNormalization(const RooArgSet* depSet=0, Bool_t force=kFALSE) ;
317  virtual void selectNormalizationRange(const char* rangeName=0, Bool_t force=kFALSE) ;
318 
319  // Helper functions for plotting
320  Bool_t plotSanityChecks(RooPlot* frame) const ;
321  void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
322  RooArgSet& projectedVars, Bool_t silent) const ;
323 
324  TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=0, const char* rangeName=0, Bool_t omitEmpty=kFALSE) const ;
325 
326 
327  Bool_t isSelectedComp() const ;
328 
329 
330  public:
331  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars) const ;
332  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
333  const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
334  RooArgSet *&cloneSet, const char* rangeName=0, const RooArgSet* condObs=0) const;
335  protected:
336 
338 
339  void plotOnCompSelect(RooArgSet* selNodes) const ;
340  RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, Double_t Z, const RooArgSet* params, const RooLinkedList& argList, Bool_t method1) const ;
341 
342  // Support interface for subclasses to advertise their analytic integration
343  // and generator capabilities in their analticalIntegral() and generateEvent()
344  // implementations.
345  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
346  const RooArgProxy& a) const ;
347  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
348  const RooArgProxy& a, const RooArgProxy& b) const ;
349  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
350  const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
351  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
352  const RooArgProxy& a, const RooArgProxy& b,
353  const RooArgProxy& c, const RooArgProxy& d) const ;
354 
355  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps,
356  const RooArgSet& set) const ;
357 
358 
359  RooAbsReal* createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
360  void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;
361 
362 
363  // Internal consistency checking (needed by RooDataSet)
364  virtual Bool_t isValid() const ;
365  virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const ;
366 
367  // Function evaluation and error tracing
368  Double_t traceEval(const RooArgSet* set) const ;
369  virtual Bool_t traceEvalHook(Double_t /*value*/) const {
370  // Hook function to add functionality to evaluation tracing in derived classes
371  return kFALSE ;
372  }
373  virtual Double_t evaluate() const = 0 ;
374 
375  // Hooks for RooDataSet interface
376  friend class RooRealIntegral ;
377  friend class RooVectorDataStore ;
378  virtual void syncCache(const RooArgSet* set=0) { getVal(set) ; }
379  virtual void copyCache(const RooAbsArg* source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE) ;
380  virtual void attachToTree(TTree& t, Int_t bufSize=32000) ;
381  virtual void attachToVStore(RooVectorDataStore& vstore) ;
382  virtual void setTreeBranchStatus(TTree& t, Bool_t active) ;
383  virtual void fillTreeBranch(TTree& t) ;
384 
385  friend class RooRealBinding ;
386  Double_t _plotMin ; // Minimum of plot range
387  Double_t _plotMax ; // Maximum of plot range
388  Int_t _plotBins ; // Number of plot bins
389  mutable Double_t _value ; // Cache for current value of object
390  TString _unit ; // Unit for objects value
391  TString _label ; // Plot label for objects value
392  Bool_t _forceNumInt ; // Force numerical integration if flag set
393 
394  mutable Float_t _floatValue ; //! Transient cache for floating point values from tree branches
395  mutable Int_t _intValue ; //! Transient cache for integer values from tree branches
396  mutable Bool_t _boolValue ; //! Transient cache for bool values from tree branches
397  mutable UChar_t _byteValue ; //! Transient cache for byte values from tree branches
398  mutable Char_t _sbyteValue ; //! Transient cache for signed byte values from tree branches
399  mutable UInt_t _uintValue ; //! Transient cache for unsigned integer values from tree branches
400 
401  friend class RooAbsPdf ;
402  friend class RooAbsAnaConvPdf ;
403  friend class RooRealProxy ;
404 
405  RooNumIntConfig* _specIntegratorConfig ; // Numeric integrator configuration specific for this object
406 
407  Bool_t _treeVar ; // !do not persist
408 
409  static Bool_t _cacheCheck ; // If true, always validate contents of clean which outcome of evaluate()
410 
411  friend class RooDataProjBinding ;
412  friend class RooAbsOptGoodnessOfFit ;
413 
414  struct PlotOpt {
428  const char* normRangeName ;
433  const char* projectionRangeName ;
435  const char* curveName ;
436  const char* addToCurveName ;
441  const char* curveNameSuffix ;
446  } ;
447 
448  // Plot implementation functions
449  virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;
450 
451 public:
452  // PlotOn with command list
453  virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;
454 
455  protected:
456  virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;
457 
458 
459 private:
460 
462  static std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > > _evalErrorList ;
464 
465  Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;
466 
467 protected:
468 
469 
470  friend class RooRealSumPdf ;
471  friend class RooAddPdf ;
472  friend class RooAddModel ;
473  void selectComp(Bool_t flag) {
474  // If flag is true, only selected component will be included in evaluates of RooAddPdf components
475  _selectComp = flag ;
476  }
477  static void globalSelectComp(Bool_t flag) ;
478  Bool_t _selectComp ; //! Component selection flag for RooAbsPdf::plotCompOn
479  static Bool_t _globalSelectComp ; // Global activation switch for component selection
480 
481  mutable RooArgSet* _lastNSet ; //!
482  static Bool_t _hideOffset ; // Offset hiding flag
483 
484  ClassDef(RooAbsReal,2) // Abstract real-valued variable
485 };
486 
487 #endif
virtual Bool_t isValid() const
Check if current value is valid.
Definition: RooAbsReal.cxx:439
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.
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
OBSOLETE – RETAINED FOR BACKWARD COMPATIBILITY. Use the plotOn(frame,Slice(...)) instead.
Bool_t postRangeFracScale
Definition: RooAbsReal.h:431
static void globalSelectComp(Bool_t flag)
Global switch controlling the activation of the selectComp() functionality.
RooAbsMoment * sigma(RooRealVar &obs, const RooArgSet &nset)
Definition: RooAbsReal.h:300
float xmin
Definition: THbookFile.cxx:93
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:278
void setServerValues(const char *tmp)
Definition: RooAbsReal.h:257
Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const
Check if allArgs contains matching elements for each name in nameList.
A RooCurve is a one-dimensional graphical representation of a real-valued function.
Definition: RooCurve.h:32
static EvalErrorIter evalErrorIter()
Definition: RooAbsReal.cxx:279
Bool_t _boolValue
Transient cache for integer values from tree branches.
Definition: RooAbsReal.h:396
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
Definition: Ifit.C:26
const char * addToCurveName
Definition: RooAbsReal.h:436
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Double_t scaleFactor
Definition: RooAbsReal.h:420
adaptor that projects a real function via summation of states provided in a dataset.
float Float_t
Definition: RtypesCore.h:53
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Definition: RooAbsReal.h:277
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:385
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:834
3-D histogram with a float per channel (see TH1 documentation)}
Definition: TH3.h:272
virtual Double_t evaluate() const =0
return c
const char Option_t
Definition: RtypesCore.h:62
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Definition: RooAbsReal.cxx:424
Bool_t isSelectedComp() const
If true, the current pdf is a selected component (for use in plotting)
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.
static std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > > _evalErrorList
Definition: RooAbsReal.h:462
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:337
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:757
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
virtual RooAbsReal * createChi2(RooDataHist &data, const RooLinkedList &cmdList)
Internal back-end function to create a chi2.
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:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
static Int_t numEvalErrorItems()
Definition: RooAbsReal.cxx:271
void selectComp(Bool_t flag)
Definition: RooAbsReal.h:473
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
RooAbsReal * createIntegral(const RooArgSet &iset, const char *rangeName) const
Definition: RooAbsReal.h:147
const Bool_t kFALSE
Definition: Rtypes.h:92
RooAbsReal * createIntObj(const RooArgSet &iset, const RooArgSet *nset, const RooNumIntConfig *cfg, const char *rangeName) const
Utility function for createIntegral that creates the actual integreal object.
Definition: RooAbsReal.cxx:578
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...
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:481
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
Definition: RooAbsMoment.h:27
Bool_t _fast
Definition: RooAbsArg.h:571
RooCmdArg Extended(Bool_t flag=kTRUE)
static void setHideOffset(Bool_t flag)
Definition: RooAbsReal.cxx:118
std::string _srvval
Definition: RooAbsReal.h:259
friend class RooAbsOptGoodnessOfFit
Definition: RooAbsReal.h:412
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 ~RooAbsReal()
Destructor.
Definition: RooAbsReal.cxx:186
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:378
RooAbsReal * createIntegral(const RooArgSet &iset, const RooArgSet &nset, const RooNumIntConfig &cfg, const char *rangeName=0) const
Definition: RooAbsReal.h:155
RooAbsReal * createRunningIntegral(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Create a running integral over this function, i.e.
RooAbsArg * createFundamental(const char *newname=0) const
Create a RooRealVar fundamental object with our properties.
const char * curveNameSuffix
Definition: RooAbsReal.h:441
Double_t traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
Definition: RooAbsReal.cxx:288
Bool_t getForceNumInt() const
Definition: RooAbsReal.h:114
const char * Data() const
Definition: TString.h:349
Int_t _plotBins
Definition: RooAbsReal.h:388
static const RooCmdArg & none()
Return reference to null argument.
Definition: RooCmdArg.cxx:50
RooAbsReal * createIntRI(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Utility function for createRunningIntegral that construct an object implementing the standard (analyt...
Bool_t _treeVar
Definition: RooAbsReal.h:407
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
const char * normRangeName
Definition: RooAbsReal.h:428
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects The class perfor...
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:395
#define ClassDef(name, id)
Definition: Rtypes.h:254
const char Int_t const char TProof Int_t stype
Definition: TXSlave.cxx:46
Int_t _intValue
Transient cache for floating point values from tree branches.
Definition: RooAbsReal.h:395
virtual void preferredObservableScanOrder(const RooArgSet &obs, RooArgSet &orderedObs) const
Interface method for function objects to indicate their prefferred order of observables for scanning ...
virtual void writeToStream(std::ostream &os, Bool_t compact) const
Write object contents to stream (dummy for now)
Definition: RooAbsReal.cxx:405
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:349
RooVectorDataStore is the abstract base class for data collection that use a TTree as internal storag...
const Text_t * getUnit() const
Definition: RooAbsReal.h:83
int d
Definition: tornado.py:11
std::map< const RooAbsArg *, std::pair< std::string, std::list< EvalError > > >::const_iterator EvalErrorIter
Definition: RooAbsReal.h:272
RooAbsMoment * sigma(RooRealVar &obs)
Definition: RooAbsReal.h:299
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...
void setIntegratorConfig()
Remove the specialized numeric integration configuration associated with this object.
friend class RooArgSet
Definition: RooAbsArg.h:469
static Bool_t _inhibitDirty
Definition: RooAbsArg.h:553
virtual Bool_t isIdentical(const RooAbsArg &other, Bool_t assumeSameType=kFALSE)
Definition: RooAbsReal.cxx:216
virtual Bool_t forceAnalyticalInt(const RooAbsArg &) const
Definition: RooAbsReal.h:104
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
void setParameterizeIntegral(const RooArgSet &paramVars)
virtual Double_t getValV(const RooArgSet *set=0) const
Return value of object.
Definition: RooAbsReal.cxx:249
Double_t getVal(const RooArgSet &set) const
Definition: RooAbsReal.h:72
virtual Int_t minTrialSamples(const RooArgSet &) const
Definition: RooAbsReal.h:182
Double_t getPropagatedError(const RooFitResult &fr)
Calculate error on self by propagated errors on parameters with correlations as given by fit result T...
RooAbsReal * createIntegral(const RooArgSet &iset, const RooNumIntConfig &cfg, const char *rangeName=0) const
Definition: RooAbsReal.h:160
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
const char * projectionRangeName
Definition: RooAbsReal.h:433
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:375
virtual void enableOffsetting(Bool_t)
Definition: RooAbsReal.h:307
virtual Bool_t setData(RooAbsData &, Bool_t=kTRUE)
Definition: RooAbsReal.h:305
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
A doubly linked list.
Definition: TList.h:47
virtual void forceNumInt(Bool_t flag=kTRUE)
Definition: RooAbsReal.h:109
static Int_t _evalErrorCount
Definition: RooAbsReal.h:463
Double_t _plotMin
Definition: RooAbsReal.h:386
Double_t addToWgtSelf
Definition: RooAbsReal.h:437
const char * curveName
Definition: RooAbsReal.h:435
virtual Bool_t isOffsetting() const
Definition: RooAbsReal.h:308
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:320
TThread * t[5]
Definition: threadsh1.C:13
Bool_t operator==(Double_t value) const
Equality operator comparing to a Double_t.
Definition: RooAbsReal.cxx:196
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Definition: RooAbsReal.cxx:363
Bool_t _forceNumInt
Definition: RooAbsReal.h:392
virtual Double_t defaultErrorLevel() const
Definition: RooAbsReal.h:189
RooNumIntConfig * _specIntegratorConfig
Definition: RooAbsReal.h:405
RooFit::MPSplit interleave
Definition: RooAbsReal.h:440
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
RooPlot * plotOnWithErrorBand(RooPlot *frame, const RooFitResult &fr, Double_t Z, const RooArgSet *params, const RooLinkedList &argList, Bool_t method1) const
Plot function or p.d.f.
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
Bool_t _selectComp
Definition: RooAbsReal.h:478
Lightweight interface adaptor that exports a RooAbsReal as a ROOT::Math::IGenFunction.
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
bool verbose
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:424
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:279
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
float xmax
Definition: THbookFile.cxx:93
RooCurve::WingMode wmode
Definition: RooAbsReal.h:432
const RooAbsData * projData
Definition: RooAbsReal.h:422
static void indent(ostringstream &buf, int indent_level)
Float_t _floatValue
Definition: RooAbsReal.h:394
virtual Bool_t traceEvalHook(Double_t) const
Definition: RooAbsReal.h:369
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
virtual void attachToTree(TTree &t, Int_t bufSize=32000)
Attach object to a branch of given TTree.
Bool_t plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
void setMessage(const char *tmp)
Definition: RooAbsReal.h:256
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:25
UChar_t _byteValue
Transient cache for bool values from tree branches.
Definition: RooAbsReal.h:397
static Bool_t hideOffset()
Definition: RooAbsReal.cxx:119
TH1 * fillHistogram(TH1 *hist, const RooArgList &plotVars, Double_t scaleFactor=1, const RooArgSet *projectedVars=0, Bool_t scaling=kTRUE, const RooArgSet *condObs=0, Bool_t setError=kTRUE) const
Fill the ROOT histogram 'hist' with values sampled from this function at the bin centers.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:35
virtual Double_t offset() const
Definition: RooAbsReal.h:309
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:53
RooAbsReal()
coverity[UNINIT_CTOR] Default constructor
Definition: RooAbsReal.cxx:130
virtual void attachToVStore(RooVectorDataStore &vstore)
static ErrorLoggingMode evalErrorLoggingMode()
Return current evaluation error logging mode.
The TH1 histogram class.
Definition: TH1.h:80
RooDataHist * fillDataHist(RooDataHist *hist, const RooArgSet *nset, Double_t scaleFactor, Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const
Fill a RooDataHist with values sampled from this function at the bin centers.
Double_t addToWgtOther
Definition: RooAbsReal.h:438
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 wi...
static Bool_t _cacheCheck
Definition: RooAbsReal.h:409
Double_t _value
Definition: RooAbsReal.h:389
static Bool_t _globalSelectComp
Component selection flag for RooAbsPdf::plotCompOn.
Definition: RooAbsReal.h:479
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:699
RooArgProxy is the abstact interface for RooAbsArg proxy classes.
Definition: RooArgProxy.h:24
virtual void printValue(std::ostream &os) const
Print object value.
Definition: RooAbsReal.cxx:414
static Bool_t _hideOffset
Definition: RooAbsReal.h:482
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
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...
EvalError(const EvalError &other)
Definition: RooAbsReal.h:255
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooFitResult * chi2FitDriver(RooAbsReal &fcn, RooLinkedList &cmdList)
Internal driver function for chi2 fits.
TString getTitle(Bool_t appendUnit=kFALSE) const
Return this variable's title string.
Definition: RooAbsReal.cxx:231
1-Dim function class
Definition: TF1.h:149
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function...
Definition: RooAbsReal.cxx:463
virtual void setTreeBranchStatus(TTree &t, Bool_t active)
(De)Activate associated tree branch
static void setCacheCheck(Bool_t flag)
Activate cache validation mode.
Bool_t inhibitDirty() const
Delete watch flag.
Definition: RooAbsArg.cxx:84
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
char Text_t
Definition: RtypesCore.h:58
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
RooAbsReal * createIntegral(const RooArgSet &iset, const RooArgSet &nset, const char *rangeName=0) const
Definition: RooAbsReal.h:151
A TTree object has a header with a name and a title.
Definition: TTree.h:98
Lightweight interface adaptor that exports a RooAbsPdf as a functor.
Definition: RooFunctor.h:25
void setUnit(const char *unit)
Definition: RooAbsReal.h:87
UInt_t _uintValue
Transient cache for signed byte values from tree branches.
Definition: RooAbsReal.h:399
unsigned char UChar_t
Definition: RtypesCore.h:34
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:450
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:418
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
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...
string message
Definition: ROOT.py:94
void plotOnCompSelect(RooArgSet *selNodes) const
Helper function for plotting of composite p.d.fs.
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
std::complex< float_v > Z
Definition: main.cpp:120
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...
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())
const Bool_t kTRUE
Definition: Rtypes.h:91
void logEvalError(const char *message, const char *serverValueString=0) const
Log evaluation error message.
static ErrorLoggingMode _evalErrorMode
Definition: RooAbsReal.h:461
Char_t _sbyteValue
Transient cache for byte values from tree branches.
Definition: RooAbsReal.h:398
float value
Definition: math.cpp:443
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:23
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)...
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:297
Double_t _plotMax
Definition: RooAbsReal.h:387
const RooArgSet * projDataSet
Definition: RooAbsReal.h:427
RooAbsMoment * mean(RooRealVar &obs, const RooArgSet &nset)
Definition: RooAbsReal.h:298
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
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 chi^2 fit to given histogram By default the fit is executed through the MINUIT commands MIG...
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:503
virtual void selectNormalizationRange(const char *rangeName=0, Bool_t force=kFALSE)
Interface function to force use of a given normalization range to interpret function value...
TString _unit
Definition: RooAbsReal.h:390
TString _label
Definition: RooAbsReal.h:391
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:27