meerkat is hosted by Hepforge, IPPP Durham
Meerkat  v1r3
Multidimensional kernel density estimation package

#include <PolynomialDensity.hh>

Inheritance diagram for PolynomialDensity:
AbsDensity

Public Member Functions

 PolynomialDensity (const char *pdfName, OneDimPhaseSpace *thePhaseSpace, UInt_t maxPower, TTree *tree, const char *var, UInt_t maxEvents=0)
 
 PolynomialDensity (const char *pdfName, AbsPhaseSpace *thePhaseSpace, UInt_t maxPower, TTree *tree, const char *var1, const char *var2, UInt_t integEvents, UInt_t maxEvents=0)
 
virtual ~PolynomialDensity ()
 Destructor. More...
 
Double_t density (std::vector< Double_t > &x)
 Calculate PDF value at the given point. More...
 
AbsPhaseSpacephaseSpace ()
 Return phase space definition for this PDF. More...
 
- Public Member Functions inherited from AbsDensity
 AbsDensity (const char *pdfName)
 Constructor. More...
 
virtual ~AbsDensity ()
 Destructor. More...
 
void slice (std::vector< Double_t > &x, UInt_t num, TH1F *hist)
 Calculate 1D slice of the PDF. More...
 
void slice (std::vector< Double_t > &x, UInt_t numx, UInt_t numy, TH2F *hist, Bool_t inPhaseSpace=true)
 Calculate 2D slice of the PDF. More...
 
double transform (TH1F *hist1, TH1F *hist2, double x)
 
void project (TH1F *hist)
 Calculate projection of the 1D PDF. More...
 
void project (TH2F *hist, Bool_t inPhaseSpace=true)
 Calculate projection of the 2D PDF. More...
 
void setMajorant (Double_t majorant)
 Set majorant for accept-reject method. More...
 
void setMaxTries (UInt_t maxTries)
 Set maximum number of tries for accept-reject method. More...
 
Double_t generate (std::vector< Double_t > &x)
 Generate a single point within the phase space according to the PDF using accept-reject method. More...
 
void generate (TNtuple *tree, UInt_t numEvents)
 Generate a sample of points within the phase space according to the PDF using accept-reject method. More...
 
const char * name (void)
 Return the name of the PDF. More...
 
void setSeed (UInt_t seed=0)
 Set random seed. More...
 

Private Member Functions

void readTuple (TTree *tree, std::vector< TString > &vars, UInt_t maxEvents=0)
 Read Ntuple contaning data points to be used for the fit. More...
 

Private Attributes

AbsPhaseSpacem_phaseSpace
 Reference to phase space. More...
 
UInt_t m_dim
 Cached dimensionality of the phase space. More...
 
Double_t m_par [100]
 Vector of fitted parameters. More...
 
UInt_t m_power
 Power of the polynomial. More...
 

Additional Inherited Members

- Protected Attributes inherited from AbsDensity
char m_name [256]
 PDF name. More...
 
Double_t m_majorant
 PDF majorant (maximum PDF value needed for accept-reject). More...
 
UInt_t m_maxTries
 Maximum number of tries for accept-reject method. More...
 
TRandom3 m_rnd
 Random number generator. More...
 

Detailed Description

Class that describes the polynomial density which can be fitted to the data sample from NTuple The phase space dimensionality must be either 1 or 2

Definition at line 14 of file PolynomialDensity.hh.

Constructor & Destructor Documentation

PolynomialDensity::PolynomialDensity ( const char *  pdfName,
OneDimPhaseSpace thePhaseSpace,
UInt_t  maxPower,
TTree *  tree,
const char *  var,
UInt_t  maxEvents = 0 
)

Constructor for 1D polynomial density. Integration over the 1D phase space is performed analytically

Parameters
[in]pdfNamePDF name
[in]thePhaseSpacephase space definition
[in]maxPowermaximum power of the polynomial
[in]treeROOT tree containing data
[in]varname of the tree variable
[in]maxEventsmaximum number of events to read from tree

Definition at line 156 of file PolynomialDensity.cpp.

PolynomialDensity::PolynomialDensity ( const char *  pdfName,
AbsPhaseSpace thePhaseSpace,
UInt_t  maxPower,
TTree *  tree,
const char *  var1,
const char *  var2,
UInt_t  integEvents,
UInt_t  maxEvents = 0 
)

Constructor for 2D polynomial density. Integration over the 2D phase space is performed numerically (using MC integration).

Parameters
[in]pdfNamePDF name
[in]thePhaseSpacephase space definition
[in]maxPowermaximum power of the polynomial
[in]treeROOT tree containing data
[in]var1name of the 1st tree variable
[in]var2name of the 2nd tree variable
[in]integEventsnumber of events for MC integration
[in]maxEventsmaximum number of events to read from tree

Definition at line 210 of file PolynomialDensity.cpp.

PolynomialDensity::~PolynomialDensity ( )
virtual

Destructor.

Definition at line 295 of file PolynomialDensity.cpp.

Member Function Documentation

Double_t PolynomialDensity::density ( std::vector< Double_t > &  x)
virtual

Calculate PDF value at the given point.

Parameters
[in]xthe point at which to calculate the PDF
Returns
PDF value (equals 1 inside phase space, 0 outside it)

Implements AbsDensity.

Definition at line 366 of file PolynomialDensity.cpp.

AbsPhaseSpace* PolynomialDensity::phaseSpace ( )
inlinevirtual

Return phase space definition for this PDF.

Returns
PDF phase space

Implements AbsDensity.

Definition at line 70 of file PolynomialDensity.hh.

void PolynomialDensity::readTuple ( TTree *  tree,
std::vector< TString > &  vars,
UInt_t  maxEvents = 0 
)
private

Read Ntuple contaning data points to be used for the fit.

Definition at line 299 of file PolynomialDensity.cpp.

Member Data Documentation

UInt_t PolynomialDensity::m_dim
private

Cached dimensionality of the phase space.

Definition at line 81 of file PolynomialDensity.hh.

Double_t PolynomialDensity::m_par[100]
private

Vector of fitted parameters.

Definition at line 84 of file PolynomialDensity.hh.

AbsPhaseSpace* PolynomialDensity::m_phaseSpace
private

Reference to phase space.

Definition at line 78 of file PolynomialDensity.hh.

UInt_t PolynomialDensity::m_power
private

Power of the polynomial.

Definition at line 87 of file PolynomialDensity.hh.


The documentation for this class was generated from the following files: