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

Abstract class which defines probability density interface. More...

#include <AbsDensity.hh>

Inheritance diagram for AbsDensity:
AdaptiveKernelDensity BinnedDensity BinnedKernelDensity DivideDensity FactorisedDensity FormulaDensity HistogramDensity KernelDensity PolynomialDensity ProductDensity SumDensity TransposedFactorisedDensity UniformDensity

Public Member Functions

 AbsDensity (const char *pdfName)
 Constructor. More...
 
virtual ~AbsDensity ()
 Destructor. More...
 
virtual Double_t density (std::vector< Double_t > &x)=0
 Calculate PDF value at the given point. More...
 
virtual AbsPhaseSpacephaseSpace ()=0
 Return phase space definition for this PDF. 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...
 

Protected Attributes

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

Abstract class which defines probability density interface.

Definition at line 16 of file AbsDensity.hh.

Constructor & Destructor Documentation

AbsDensity::AbsDensity ( const char *  pdfName)

Constructor.

Parameters
[in]pdfNamethe PDF name

Definition at line 14 of file AbsDensity.cpp.

AbsDensity::~AbsDensity ( )
virtual

Destructor.

Definition at line 21 of file AbsDensity.cpp.

Member Function Documentation

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

Calculate PDF value at the given point.

Parameters
[in]xthe point at which to calculate the PDF
Returns
PDF value

Implemented in AdaptiveKernelDensity, BinnedKernelDensity, BinnedDensity, SumDensity, KernelDensity, PolynomialDensity, FactorisedDensity, ProductDensity, FormulaDensity, TransposedFactorisedDensity, DivideDensity, HistogramDensity, and UniformDensity.

Double_t AbsDensity::generate ( std::vector< Double_t > &  x)

Generate a single point within the phase space according to the PDF using accept-reject method.

Parameters
[out]xthe generated point
Returns
the PDF value at this point.

Definition at line 99 of file AbsDensity.cpp.

void AbsDensity::generate ( TNtuple *  tree,
UInt_t  numEvents 
)

Generate a sample of points within the phase space according to the PDF using accept-reject method.

Store them in the ROOT NTuple.

Parameters
[out]treeROOT NTuple to fill
[in]numEventsnumber of events to generate

Definition at line 150 of file AbsDensity.cpp.

const char* AbsDensity::name ( void  )
inline

Return the name of the PDF.

Returns
PDF name

Definition at line 106 of file AbsDensity.hh.

virtual AbsPhaseSpace* AbsDensity::phaseSpace ( )
pure virtual
void AbsDensity::project ( TH1F *  hist)

Calculate projection of the 1D PDF.

Parameters
[out]histoutput 1D ROOT histogram

Definition at line 76 of file AbsDensity.cpp.

void AbsDensity::project ( TH2F *  hist,
Bool_t  inPhaseSpace = true 
)

Calculate projection of the 2D PDF.

Parameters
[out]histoutput 2D ROOT histogram
[in]inPhaseSpaceIf False, estimates the PDF somewhat outside of the phase space. Otherwise the bins outside of phase space are assigned zero.

Definition at line 87 of file AbsDensity.cpp.

void AbsDensity::setMajorant ( Double_t  majorant)
inline

Set majorant for accept-reject method.

Parameters
[in]majorantthe majorant

Definition at line 79 of file AbsDensity.hh.

void AbsDensity::setMaxTries ( UInt_t  maxTries)
inline

Set maximum number of tries for accept-reject method.

Parameters
[in]maxTriesmaximum number of tries

Definition at line 85 of file AbsDensity.hh.

void AbsDensity::setSeed ( UInt_t  seed = 0)
inline

Set random seed.

Set random seed

Parameters
[in]seedSeed

Definition at line 113 of file AbsDensity.hh.

void AbsDensity::slice ( std::vector< Double_t > &  x,
UInt_t  num,
TH1F *  hist 
)

Calculate 1D slice of the PDF.

Parameters
[in]xthe point at which to calculate the slice
[in]numnumber of the phase space variable to scan (value of x for this variable will be ignored)
[out]histoutput 1D ROOT histogram

Definition at line 25 of file AbsDensity.cpp.

void AbsDensity::slice ( std::vector< Double_t > &  x,
UInt_t  numx,
UInt_t  numy,
TH2F *  hist,
Bool_t  inPhaseSpace = true 
)

Calculate 2D slice of the PDF.

Parameters
[in]xthe point at which to calculate the slice
[in]numxnumber of the 1st phase space variable to scan (value of x for this variable will be ignored)
[in]numynumber of the 2nd phase space variable to scan (value of x for this variable will be ignored)
[out]histoutput 2D ROOT histogram
[in]inPhaseSpaceIf False, estimates the PDF somewhat outside of the phase space. Otherwise the bins outside of phase space are assigned zero.

Definition at line 48 of file AbsDensity.cpp.

double AbsDensity::transform ( TH1F *  hist1,
TH1F *  hist2,
double  x 
)

Definition at line 179 of file AbsDensity.cpp.

Member Data Documentation

Double_t AbsDensity::m_majorant
protected

PDF majorant (maximum PDF value needed for accept-reject).

Definition at line 121 of file AbsDensity.hh.

UInt_t AbsDensity::m_maxTries
protected

Maximum number of tries for accept-reject method.

Definition at line 124 of file AbsDensity.hh.

char AbsDensity::m_name[256]
protected

PDF name.

Definition at line 118 of file AbsDensity.hh.

TRandom3 AbsDensity::m_rnd
protected

Random number generator.

Definition at line 127 of file AbsDensity.hh.


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