library: libRooFit #include "RooResolutionModel.h" |
RooResolutionModel
class description - source file - inheritance tree (.pdf)
This is an abstract class, constructors will not be documented.
Look at the header to check for available constructors.
protected:
virtual void changeBasis(RooFormulaVar* basis)
virtual Bool_t redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
virtual Bool_t syncNormalizationPreHook(RooAbsReal* norm, const RooArgSet* nset) const
virtual Bool_t traceEvalHook(Double_t value) const
public:
virtual ~RooResolutionModel()
const RooFormulaVar& basis() const
virtual Int_t basisCode(const char* name) const
const RooRealVar& basisConvVar() const
static TClass* Class()
virtual TObject* clone(const char* newname) const
virtual RooResolutionModel* convolution(RooFormulaVar* basis, RooAbsArg* owner) const
RooRealVar& convVar() const
virtual Double_t getNorm(const RooArgSet* nset = 0) const
virtual Double_t getVal(const RooArgSet* nset = 0) const
virtual TClass* IsA() const
Bool_t isBasisSupported(const char* name) const
virtual void normLeafServerList(RooArgSet& list) const
virtual void printToStream(ostream& stream, RooPrintable::PrintOption opt = Standard, TString indent = ) const
virtual void ShowMembers(TMemberInspector& insp, char* parent)
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
protected:
static RooFormulaVar* _identity Identity basis function pointer
RooRealProxy x Dependent/convolution variable
Int_t _basisCode Identifier code for selected basis function
RooFormulaVar* _basis Basis function convolved with this resolution model
Bool_t _ownBasis Flag indicating ownership of _basis
RooResolutionModel is the base class of for PDFs that represent a
resolution model that can be convoluted with physics a physics model of the form
Phys(x,a,b) = Sum_k coef_k(a) * basis_k(x,b)
where basis_k are a limited number of functions in terms of the variable
to be convoluted and coef_k are coefficients independent of the convolution
variable.
Classes derived from RooResolutionModel implement
_ _ _ _
R_k(x,b,c) = Int(dx') basis_k(x',b) * resModel(x-x',c)
which RooAbsAnaConvPdf uses to construct the pdf for [ Phys (x) R ] :
_ _ _ _ _ _
PDF(x,a,b,c) = Sum_k coef_k(a) * R_k(x,b,c)
A minimal implementation of a RooResolutionModel consists of a
Int_t basisCode(const char* name)
function indication which basis functions this resolution model supports, and
Double_t evaluate()
Implementing the resolution model, optionally convoluted with one of the
supported basis functions. RooResolutionModel objects can be used as regular
PDFs (They inherit from RooAbsPdf), or as resolution model convoluted with
a basis function. The implementation of evaluate() can identify the requested
from of use from the basisCode() function. If zero, the regular PDF value
should be calculated. If non-zero, the models value convoluted with the
basis function identified by the code should be calculated.
Optionally, analytical integrals can be advertised and implemented, in the
same way as done for regular PDFS (see RooAbsPdf for further details).
Also in getAnalyticalIntegral()/analyticalIntegral() the implementation
should use basisCode() to determine for which scenario the integral is
requested.
The choice of basis returned by basisCode() is guaranteed not to change
of the lifetime of a RooResolutionModel object.
~RooResolutionModel()
Destructor
RooResolutionModel* convolution(RooFormulaVar* basis, RooAbsArg* owner) const
Instantiate a clone of this resolution model representing a convolution with given
basis function. The owners object name is incorporated in the clones name
to avoid multiple convolution objects with the same name in complex PDF structures.
void changeBasis(RooFormulaVar* basis)
Change the basis function we convolute with.
For one-time use by convolution() only.
const RooRealVar& basisConvVar() const
Return the convolution variable of the selection basis function.
This is, by definition, the first parameter of the basis function
RooRealVar& convVar() const
Return the convolution variable of the resolution model
Double_t getVal(const RooArgSet* nset) const
Modified version of RooAbsPdf::getVal(). If used as regular PDF,
call RooAbsPdf::getVal(), otherwise return unnormalized value
regardless of specified normalization set
Bool_t redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t /*isRecursive*/)
Forward redirectServers call to our basis function, which is not connected to either resolution
model or the physics model.
Bool_t traceEvalHook(Double_t value) const
Floating point error checking and tracing for given float value
void normLeafServerList(RooArgSet& list) const
Return the list of servers used by our normalization integral
Double_t getNorm(const RooArgSet* nset) const
Return the integral of this PDF over all elements of 'nset'.
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
Bool_t syncNormalizationPreHook(RooAbsReal*,const RooArgSet*) const
Inline Functions
TObject* clone(const char* newname) const
Bool_t isBasisSupported(const char* name) const
Int_t basisCode(const char* name) const
const RooFormulaVar& basis() const
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:15:02 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.