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

Class that describes the "divide" density: the ratio of two densities in the same phase spaces. More...

#include <DivideDensity.hh>

Inheritance diagram for DivideDensity:
AbsDensity

Public Member Functions

 DivideDensity (const char *pdfName, AbsPhaseSpace *thePhaseSpace, AbsDensity *d1, AbsDensity *d2)
 Constructor of Divide density. More...
 
virtual ~DivideDensity ()
 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 Attributes

AbsPhaseSpacem_phaseSpace
 Reference to phase space. More...
 
AbsDensitym_density1
 Density component 1. More...
 
AbsDensitym_density2
 Density component 2. More...
 
UInt_t m_dim
 Cached dimensionality of the phase space. 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 "divide" density: the ratio of two densities in the same phase spaces.

Definition at line 13 of file DivideDensity.hh.

Constructor & Destructor Documentation

DivideDensity::DivideDensity ( const char *  pdfName,
AbsPhaseSpace thePhaseSpace,
AbsDensity d1,
AbsDensity d2 
)

Constructor of Divide density.

Parameters
[in]pdfNamePDF name
[in]thePhaseSpacephase space. Dimensionality of the phase space should be equal to the ones of the components.
[in]d1numerator density component.
[in]d2denominator density component.

Definition at line 13 of file DivideDensity.cpp.

DivideDensity::~DivideDensity ( )
virtual

Destructor.

Definition at line 37 of file DivideDensity.cpp.

Member Function Documentation

Double_t DivideDensity::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

Implements AbsDensity.

Definition at line 41 of file DivideDensity.cpp.

AbsPhaseSpace* DivideDensity::phaseSpace ( )
inlinevirtual

Return phase space definition for this PDF.

Returns
PDF phase space

Implements AbsDensity.

Definition at line 43 of file DivideDensity.hh.

Member Data Documentation

AbsDensity* DivideDensity::m_density1
private

Density component 1.

Definition at line 51 of file DivideDensity.hh.

AbsDensity* DivideDensity::m_density2
private

Density component 2.

Definition at line 54 of file DivideDensity.hh.

UInt_t DivideDensity::m_dim
private

Cached dimensionality of the phase space.

Definition at line 57 of file DivideDensity.hh.

AbsPhaseSpace* DivideDensity::m_phaseSpace
private

Reference to phase space.

Definition at line 48 of file DivideDensity.hh.


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