library: libRooFit
#include "RooAbsReal.h"

RooAbsReal


class description - source file - inheritance tree (.pdf)

class RooAbsReal : public RooAbsArg

Inheritance Chart:
TObject
<-
TNamed
RooPrintable
<-
RooAbsArg
<-
RooAbsReal
<-
RooAbsGoodnessOfFit
<-
RooAbsOptGoodnessOfFit
<-
RooChi2Var
RooNLLVar
RooAbsHiddenReal
<-
RooUnblindCPAsymVar
RooUnblindOffset
RooUnblindPrecision
RooUnblindUniform
RooAbsPdf
 [more...]
 
    This is an abstract class, constructors will not be documented.
    Look at the header to check for available constructors.

    private:
Bool_t matchArgsByName(const RooArgSet& allArgs, RooArgSet& matchedArgs, const TList& nameList) const protected:
Bool_t allClientsCached(RooAbsArg* var, RooArgSet& cacheList) virtual void attachToTree(TTree& t, Int_t bufSize = 32000) virtual void copyCache(const RooAbsArg* source) const RooAbsReal* createProjection(const RooArgSet& dependentVars, const RooArgSet* projectedVars, RooArgSet*& cloneSet, const char* rangeName = "0") const void doConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose) virtual Double_t evaluate() const virtual void fillTreeBranch(TTree& t) Bool_t findCacheableBranches(RooAbsArg* arg, RooAbsData* dset, RooArgSet& cacheList, const RooArgSet* normSet, Bool_t verbose) void findRedundantCacheServers(RooAbsData* dset, RooArgSet& cacheList, RooArgSet& pruneList, Bool_t verbose) void findUnusedDataVariables(RooAbsData* dset, RooArgSet& pruneList, Bool_t verbose) TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset = 0, const char* rangeName = "0") const virtual Bool_t isValid() const virtual Bool_t isValidReal(Double_t value, Bool_t printError = kFALSE) const void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars, RooArgSet& projectedVars, Bool_t silent) const Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a) const Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b) const Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c, const RooArgProxy& d) const Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgSet& set) const void optimizeDirty(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose) virtual RooPlot* plotAsymOn(RooPlot* frame, const RooAbsCategoryLValue& asymCat, RooAbsReal::PlotOpt o) const virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const virtual RooPlot* plotOn(RooPlot* frame, RooAbsReal::PlotOpt o) const Bool_t plotSanityChecks(RooPlot* frame) const virtual void selectNormalization(const RooArgSet* depSet = 0, Bool_t force = kFALSE) virtual void selectNormalizationRange(const char* rangeName = "0", Bool_t force = kFALSE) virtual void setTreeBranchStatus(TTree& t, Bool_t active) virtual void syncCache(const RooArgSet* set = 0) Double_t traceEval(const RooArgSet* set) const virtual Bool_t traceEvalHook(Double_t) const void undoConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose) public:
virtual ~RooAbsReal() virtual Double_t analyticalIntegral(Int_t code, const char* rangeName = "0") const virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName = "0") const RooAbsFunc* bindVars(const RooArgSet& vars, const RooArgSet* nset = 0, Bool_t clipInvalid = kFALSE) const static TClass* Class() virtual RooAbsArg* createFundamental(const char* newname = "0") const TH1* createHistogram(const char* name, const RooAbsRealLValue& xvar, 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) const 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 RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName = "0") const RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, RooNumIntConfig& cfg, const char* rangeName = "0") const RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName = "0") const virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset = 0, const RooNumIntConfig* cfg = 0, const char* rangeName = "0") const const RooAbsReal* createProjection(const RooArgSet& depVars, const RooArgSet& projVars) const const RooAbsReal* createProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const virtual Double_t defaultErrorLevel() const static RooNumIntConfig* defaultIntegratorConfig() TH1* fillHistogram(TH1* hist, const RooArgList& plotVars, Double_t scaleFactor = 1, const RooArgSet* projectedVars = 0) const virtual Bool_t forceAnalyticalInt(const RooAbsArg&) const virtual void forceNumInt(Bool_t flag = kTRUE) virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName = "0") const virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName = "0") const const RooNumIntConfig* getIntegratorConfig() const virtual Int_t getMaxVal(const RooArgSet& vars) const virtual Int_t getPlotBins() const const char* getPlotLabel() const Double_t getPlotMax() const Double_t getPlotMin() const TString getTitle(Bool_t appendUnit = kFALSE) const const Text_t* getUnit() const virtual Double_t getVal(const RooArgSet* set = 0) const Double_t getVal(const RooArgSet& set) const virtual Bool_t inPlotRange(Double_t value) const virtual TClass* IsA() const virtual Double_t maxVal(Int_t code) Bool_t operator==(Double_t value) const virtual Bool_t operator==(const RooAbsArg& other) 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 virtual RooPlot* plotSliceOn(RooPlot* frame, const RooArgSet& sliceSet, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0) const virtual void printToStream(ostream& stream, RooPrintable::PrintOption opt = Standard, TString indent = ) const virtual Bool_t readFromStream(istream& is, Bool_t compact, Bool_t verbose = kFALSE) static void setCacheCheck(Bool_t flag) void setIntegratorConfig() void setIntegratorConfig(const RooNumIntConfig& config) void setPlotBins(Int_t value) void setPlotLabel(const char* label) void setPlotMax(Double_t value) void setPlotMin(Double_t value) void setPlotRange(Double_t min, Double_t max) void setUnit(const char* unit) virtual void ShowMembers(TMemberInspector& insp, char* parent) RooNumIntConfig* specialIntegratorConfig() const virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) virtual void writeToStream(ostream& os, Bool_t compact) const

Data Members


    protected:
Double_t _plotMin Minimum of plot range Double_t _plotMax Maximum of plot range Int_t _plotBins Number of plot bins Double_t _value Cache for current value of object TString _unit Unit for objects value TString _label Plot label for objects value Bool_t _forceNumInt Force numerical integration if flag set RooNumIntConfig* _specIntegratorConfig ! Numeric integrator configuration specific for this object static Bool_t _cacheCheck public:
static const RooAbsReal::ScaleType Raw static const RooAbsReal::ScaleType Relative static const RooAbsReal::ScaleType NumEvent static const RooAbsReal::ScaleType RelativeExpected

Class Description

 This class holds in addition a unit and label string, as well
 as a plot range and number of plot bins and plot creation methods.

~RooAbsReal()

TString getTitle(Bool_t appendUnit) const
 Return this variable's title string. If appendUnit is true and
 this variable has units, also append a string " (<unit>)".

Double_t getVal(const RooArgSet* set) const
 Return value of object. Calculated if dirty, otherwise cached value is returned.

Double_t traceEval(const RooArgSet* /*nset*/) const
 Calculate current value of object, with error tracing wrapper

Int_t getAnalyticalIntegralWN(RooArgSet& allDeps, RooArgSet& analDeps, const RooArgSet* /*normSet*/, const char* rangeName) const
 Default implementation of getAnalyticalIntegralWN for real valued objects defers to
 normalization invariant getAnalyticalIntegral()

Int_t getAnalyticalIntegral(RooArgSet& /*allDeps*/, RooArgSet& /*analDeps*/, const char* /*rangeName*/) const
 By default we do not supply any analytical integrals

Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
 Default implementation of analyticalIntegralWN handles only the pass-through
 scenario (code =0). All other codes are deferred to to the normalization
 invariant analyticalIntegral()

Double_t analyticalIntegral(Int_t code, const char* /*rangeName*/) const
 By default no analytical integrals are implemented

const char* getPlotLabel() const
 Get the label associated with the variable

void setPlotLabel(const char *label)
 Set the label associated with this variable

Bool_t readFromStream(istream& /*is*/, Bool_t /*compact*/, Bool_t /*verbose*/)
Read object contents from stream (dummy for now)

void writeToStream(ostream& /*os*/, Bool_t /*compact*/) const
Write object contents to stream (dummy for now)

void printToStream(ostream& os, PrintOption opt, TString indent) const
 Print info about this object to the specified stream. In addition to the info
 from RooAbsArg::printToStream() we add:

     Shape : value, units, plot range
   Verbose : default binning and print label

void setPlotMin(Double_t value)
 Set minimum value of output associated with this object

void setPlotMax(Double_t value)
 Set maximum value of output associated with this object

void setPlotRange(Double_t, Double_t)
 Set a new plot range

void setPlotBins(Int_t /*value*/)
 Set number of histogram bins

Bool_t inPlotRange(Double_t value) const
 Check if given value is in the min-max range for this object

Bool_t isValid() const
 Check if current value is valid

Bool_t isValidReal(Double_t /*value*/, Bool_t /*printError*/) const
 Check if given value is valid

RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg arg1, const RooCmdArg arg2, const RooCmdArg arg3, const RooCmdArg arg4, const RooCmdArg arg5, const RooCmdArg arg6, const RooCmdArg arg7, const RooCmdArg arg8) const
 Create an object that represents the integral of the function over one or more observables listed in iset
 The actual integration calculation is only performed when the return object is evaluated. The name
 of the integral object is automatically constructed from the name of the input function, the variables
 it integrates and the range integrates over

 The following named arguments are accepted

 NormSet(const RooArgSet&)            -- Specify normalization set, mostly useful when working with PDFS
 NumIntConfig(const RooNumIntConfig&) -- Use given configuration for any numeric integration, if necessary
 Range(const char* name)              -- Integrate only over given range. Multiple ranges may be specified
                                         by passing multiple Range() arguments

RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const

TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset, const char* rangeName) const

const RooAbsReal* createProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const

const RooAbsReal* createProjection(const RooArgSet& depVars, const RooArgSet& projVars) const

const RooAbsReal* createProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars, RooArgSet *&cloneSet, const char* rangeName) const
 Create a new object G that represents the normalized projection:

             Integral [ F[x,y,p] , { y } ]
  G[x,p] = ---------------------------------
            Integral [ F[x,y,p] , { x,y } ]

 where F[x,y,p] is the function we represent, "x" are the
 specified dependentVars, "y" are the specified projectedVars, and
 "p" are our remaining variables ("parameters"). Return a
 pointer to the newly created object, or else zero in case of an
 error.  The caller is responsible for deleting the contents of
 cloneSet (which includes the returned projection object) whatever
 the return value. Note that you should normally call getVal()
 on the returned object, without providing any set of normalization
 variables. Otherwise you are requesting an additional normalization
 beyond what is already specified in the equation above.

TH1* fillHistogram(TH1 *hist, const RooArgList &plotVars, Double_t scaleFactor, const RooArgSet *projectedVars) const
 Loop over the bins of the input histogram and add an amount equal to our value evaluated
 at the bin center to each one. Our value is calculated by first integrating out any variables
 in projectedVars and then scaling the result by scaleFactor. Returns a pointer to the
 input histogram, or zero in case of an error. The input histogram can be any TH1 subclass, and
 therefore of arbitrary dimension. Variables are matched with the (x,y,...) dimensions of the input
 histogram according to the order in which they appear in the input plotVars list.

TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
 Create and fill a ROOT histogram TH1,TH2 or TH3 with the values of this function.

 This function accepts the following arguments

 name -- Name of the ROOT histogram
 xvar -- Observable to be mapped on x axis of ROOT histogram

 Binning(const char* name)                    -- Apply binning with given name to x axis of histogram
 Binning(RooAbsBinning& binning)              -- Apply specified binning to x axis of histogram
 Binning(double lo, double hi, int nbins)     -- Apply specified binning to x axis of histogram
 ConditionalObservables(const RooArgSet& set) -- Do not normalized PDF over following observables when projecting PDF into histogram

 YVar(const RooAbsRealLValue& var,...)    -- Observable to be mapped on y axis of ROOT histogram
 ZVar(const RooAbsRealLValue& var,...)    -- Observable to be mapped on z axis of ROOT histogram

 The YVar() and ZVar() arguments can be supplied with optional Binning() arguments to control the binning of the Y and Z axes, e.g.
 createHistogram("histo",x,Binning(-1,1,20), YVar(y,Binning(-1,1,30)), ZVar(z,Binning("zbinning")))

 The caller takes ownership of the returned histogram

RooPlot* plotOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8, const RooCmdArg& arg9, const RooCmdArg& arg10) const
 Plot (project) PDF on specified frame. If a PDF is plotted in an empty frame, it
 will show a unit normalized curve in the frame variable, taken at the present value
 of other observables defined for this PDF

 If a PDF is plotted in a frame in which a dataset has already been plotted, it will
 show a projected curve integrated over all variables that were present in the shown
 dataset except for the one on the x-axis. The normalization of the curve will also
 be adjusted to the event count of the plotted dataset. An informational message
 will be printed for each projection step that is performed

 This function takes the following named arguments

 Projection control
 ------------------
 Slice(const RooArgSet& set)     -- Override default projection behaviour by omittting observables listed
                                    in set from the projection, resulting a 'slice' plot. Slicing is usually
                                    only sensible in discrete observables
 Project(const RooArgSet& set)   -- Override default projection behaviour by projecting over observables
                                    given in set and complete ignoring the default projection behavior. Advanced use only.
 ProjWData(const RooAbsData& d)  -- Override default projection _technique_ (integration). For observables present in given dataset
                                    projection of PDF is achieved by constructing an average over all observable values in given set.
                                    Consult RooFit plotting tutorial for further explanation of meaning & use of this technique
 ProjWData(const RooArgSet& s,   -- As above but only consider subset 's' of observables in dataset 'd' for projection through data averaging
           const RooAbsData& d)
 ProjectionRange(const char* rn) -- Override default range of projection integrals to a different range speficied by given range name.
                                    This technique allows you to project a finite width slice in a real-valued observable

 Misc content control
 --------------------
 Normalization(Double_t scale,   -- Adjust normalization by given scale factor. Interpretation of number depends on code: Relative:
                ScaleType code)     relative adjustment factor, NumEvent: scale to match given number of events.
 Name(const chat* name)          -- Give curve specified name in frame. Useful if curve is to be referenced later
 Asymmetry(const RooCategory& c) -- Show the asymmetry of the PDF in given two-state category [F(+)-F(-)] / [F(+)+F(-)] rather than
                                    the PDF projection. Category must have two states with indices -1 and +1 or three states with
                                    indeces -1,0 and +1.
 ShiftToZero(Bool_t flag)        -- Shift entire curve such that lowest visible point is at exactly zero. Mostly useful when
                                    plotting -log(L) or chi^2 distributions
 AddTo(const char* name,         -- Add constructed projection to already existing curve with given name and relative weight factors
       double_t wgtSelf, double_t wgtOther)

 Plotting control
 ----------------
 DrawOption(const char* opt)     -- Select ROOT draw option for resulting TGraph object
 LineStyle(Int_t style)          -- Select line style by ROOT line style code, default is solid
 LineColor(Int_t color)          -- Select line color by ROOT color code, default is blue
 LineWidth(Int_t width)          -- Select line with in pixels, default is 3
 FillStyle(Int_t style)          -- Select fill style, default is not filled. If a filled style is selected, also use VLines()
                                    to add vertical downward lines at end of curve to ensure proper closure
 FillColor(Int_t color)          -- Select fill color by ROOT color code
 Range(const char* name)         -- Only draw curve in range defined by given name
 Range(double lo, double hi)     -- Only draw curve in specified range
 VLines()                        -- Add vertical lines to y=0 at end points of curve
 Precision(Double_t eps)         -- Control precision of drawn curve w.r.t to scale of plot, default is 1e-3. Higher precision
                                    will result in more and more densely spaced curve points
 Invisble(Bool_t flag)           -- Add curve to frame, but do not display. Useful in combination AddTo()

RooPlot* plotOn(RooPlot* frame, RooLinkedList& argList) const

RooPlot* plotOn(RooPlot *frame, PlotOpt o) const
 Plot ourselves on given frame. If frame contains a histogram, all dimensions of the plotted
 function that occur in the previously plotted dataset are projected via partial integration,
 otherwise no projections are performed. Optionally, certain projections can be performed
 by summing over the values present in a provided dataset ('projData'), to correctly
 project out data dependents that are not properly described by the PDF (e.g. per-event errors).

 The functions value can be multiplied with an optional scale factor. The interpretation
 of the scale factor is unique for generic real functions, for PDFs there are various interpretations
 possible, which can be selection with 'stype' (see RooAbsPdf::plotOn() for details).

 The default projection behaviour can be overriden by supplying an optional set of dependents
 to project. For most cases, plotSliceOn() and plotProjOn() provide a more intuitive interface
 to modify the default projection behavour.

RooPlot* plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData) const
 Plot ourselves on given frame, as done in plotOn(), except that the variables
 listed in 'sliceSet' are taken out from the default list of projected dimensions created
 by plotOn().

RooPlot* plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const
 Plot asymmetry of ourselves, defined as

   asym = f(asymCat=-1) - f(asymCat=+1) / ( f(asymCat=-1) + f(asymCat=+1) )

 on frame. If frame contains a histogram, all dimensions of the plotted
 asymmetry function that occur in the previously plotted dataset are projected via partial integration.
 Otherwise no projections are performed,

 The asymmetry function can be multiplied with an optional scale factor. The default projection
 behaviour can be overriden by supplying an optional set of dependents to project.

Bool_t plotSanityChecks(RooPlot* frame) const
 Perform general sanity check on frame to ensure safe plotting operations

void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars, RooArgSet& projectedVars, Bool_t silent) const
 Construct the set of dependents to project when plotting ourselves as function
 of 'plotVar'. 'allVars' is the list of variables that must be projected, but
 may contain variables that we do not depend on. If 'silent' is cleared,
 warnings about inconsistent input parameters will be printed.

RooAbsFunc* bindVars(const RooArgSet &vars, const RooArgSet* nset, Bool_t clipInvalid) const
 Create an interface adaptor f(vars) that binds us to the specified variables
 (in arbitrary order). For example, calling bindVars({x1,x3}) on an object
 F(x1,x2,x3,x4) returns an object f(x1,x3) that is evaluated using the
 current values of x2 and x4. The caller takes ownership of the returned adaptor.

void copyCache(const RooAbsArg* source)
 Copy the cached value of another RooAbsArg to our cache

void attachToTree(TTree& t, Int_t bufSize)
 Attach object to a branch of given TTree

void fillTreeBranch(TTree& t)
 Attach object to a branch of given TTree

void setTreeBranchStatus(TTree& t, Bool_t active)
 (De)Activate associate tree branch

RooAbsArg* createFundamental(const char* newname) const
 Create a RooRealVar fundamental object with our properties. The new
 object will be created without any fit limits.

Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps, const RooArgProxy& a) const
 Wrapper function for matchArgsByName()

Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps, const RooArgProxy& a, const RooArgProxy& b) const
 Wrapper function for matchArgsByName()

Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const
 Wrapper function for matchArgsByName()

Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c, const RooArgProxy& d) const
 Wrapper function for matchArgsByName()

Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps, const RooArgSet& set) const
 Wrapper function for matchArgsByName()

Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const
 Check if allArgs contains matching elements for each name in nameList. If it does,
 add the corresponding args from allArgs to matchedArgs and return kTRUE. Otherwise
 return kFALSE and do not change matchedArgs.

RooNumIntConfig* defaultIntegratorConfig()

RooNumIntConfig* specialIntegratorConfig() const

const RooNumIntConfig* getIntegratorConfig() const

void setIntegratorConfig(const RooNumIntConfig& config)

void setIntegratorConfig()

void optimizeDirty(RooAbsData& dataset, const RooArgSet* normSet, Bool_t /*verbose*/)

void doConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose)
 optimizeDirty must have been run first!

void undoConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose)
 Delete the cache

Bool_t findCacheableBranches(RooAbsArg* arg, RooAbsData* dset, RooArgSet& cacheList, const RooArgSet* normSet, Bool_t verbose)
 Find branch PDFs with all-constant parameters, and add them
 to the dataset cache list

void findUnusedDataVariables(RooAbsData* dset,RooArgSet& pruneList, Bool_t /*verbose*/)

void findRedundantCacheServers(RooAbsData* dset,RooArgSet& cacheList, RooArgSet& pruneList, Bool_t /*verbose*/)

Bool_t allClientsCached(RooAbsArg* var, RooArgSet& cacheList)

void selectNormalization(const RooArgSet*, Bool_t)

void selectNormalizationRange(const char*, Bool_t)

Int_t getMaxVal(const RooArgSet& /*vars*/) const // Advertise capability to determine maximum value of function for given set of // observables. If no direct generator method is provided, this information // will assist the accept/reject generator to operate more efficiently as // it can skip the initial trial sampling phase to empirically find the function // maximum

Double_t maxVal(Int_t /*code*/)
 Return maximum value for set of observables identified by code assigned
 in getMaxVal



Inline Functions


             Double_t getVal(const RooArgSet& set) const
               Bool_t operator==(Double_t value) const
               Bool_t operator==(const RooAbsArg& other)
        const Text_t* getUnit() const
                 void setUnit(const char* unit)
               Bool_t forceAnalyticalInt(const RooAbsArg&) const
                 void forceNumInt(Bool_t flag = kTRUE)
          RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName = "0") const
          RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, RooNumIntConfig& cfg, const char* rangeName = "0") const
          RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName = "0") const
          RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset = 0, const RooNumIntConfig* cfg = 0, const char* rangeName = "0") const
             Double_t getPlotMin() const
             Double_t getPlotMax() const
                Int_t getPlotBins() const
             Double_t defaultErrorLevel() const
                 void setCacheCheck(Bool_t flag)
               Bool_t traceEvalHook(Double_t) const
             Double_t evaluate() const
                 void syncCache(const RooArgSet* set = 0)
              TClass* Class()
              TClass* IsA() const
                 void ShowMembers(TMemberInspector& insp, char* parent)
                 void Streamer(TBuffer& b)
                 void StreamerNVirtual(TBuffer& b)


Last update: Tue Jun 28 18:10:11 2005
Copyright (c) 2000-2005, Regents of the University of California *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.