|
Meerkat
v1r3
Multidimensional kernel density estimation package
|
Abstract class which defines probability density interface.
More...
#include <AbsDensity.hh>
|
| 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 AbsPhaseSpace * | phaseSpace ()=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...
|
|
Abstract class which defines probability density interface.
Definition at line 16 of file AbsDensity.hh.
AbsDensity::AbsDensity |
( |
const char * |
pdfName | ) |
|
AbsDensity::~AbsDensity |
( |
| ) |
|
|
virtual |
virtual Double_t AbsDensity::density |
( |
std::vector< Double_t > & |
x | ) |
|
|
pure virtual |
Calculate PDF value at the given point.
- Parameters
-
[in] | x | the 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] | x | the 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] | tree | ROOT NTuple to fill |
[in] | numEvents | number 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.
Return phase space definition for this PDF.
- Returns
- PDF phase space
Implemented in AdaptiveKernelDensity, BinnedKernelDensity, BinnedDensity, SumDensity, KernelDensity, PolynomialDensity, FactorisedDensity, ProductDensity, FormulaDensity, TransposedFactorisedDensity, DivideDensity, HistogramDensity, and UniformDensity.
void AbsDensity::project |
( |
TH1F * |
hist | ) |
|
Calculate projection of the 1D PDF.
- Parameters
-
[out] | hist | output 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] | hist | output 2D ROOT histogram |
[in] | inPhaseSpace | If 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
-
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] | maxTries | maximum 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
-
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] | x | the point at which to calculate the slice |
[in] | num | number of the phase space variable to scan (value of x for this variable will be ignored) |
[out] | hist | output 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] | x | the point at which to calculate the slice |
[in] | numx | number of the 1st phase space variable to scan (value of x for this variable will be ignored) |
[in] | numy | number of the 2nd phase space variable to scan (value of x for this variable will be ignored) |
[out] | hist | output 2D ROOT histogram |
[in] | inPhaseSpace | If 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 |
|
) |
| |
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 |
TRandom3 AbsDensity::m_rnd |
|
protected |
The documentation for this class was generated from the following files:
|