|
Meerkat
v1r3
Multidimensional kernel density estimation package
|
Class that describes the unbinned kernel density.
More...
#include <KernelDensity.hh>
|
| KernelDensity (const char *pdfName, AbsPhaseSpace *thePhaseSpace, std::vector< Double_t > &width, UInt_t approxSize, AbsDensity *approxDensity=0) |
|
| KernelDensity (const char *pdfName, AbsPhaseSpace *thePhaseSpace, UInt_t approxSize, Double_t width1, Double_t width2=0, Double_t width3=0, Double_t width4=0, Double_t width5=0) |
|
| KernelDensity (const char *pdfName, AbsPhaseSpace *thePhaseSpace, UInt_t approxSize, AbsDensity *approxDensity, Double_t width1, Double_t width2=0, Double_t width3=0, Double_t width4=0, Double_t width5=0) |
|
virtual | ~KernelDensity () |
|
void | setWidth (std::vector< Double_t > &width) |
| Set kernel width. More...
|
|
Bool_t | generateApproximation (UInt_t approxSize) |
| Create normalisation vector. More...
|
|
Bool_t | readTuple (TTree *tree, std::vector< TString > &vars, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, const char *var2, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, const char *var2, const char *var3, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, const char *var2, const char *var3, const char *var4, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, const char *var2, const char *var3, const char *var4, const char *var5, UInt_t maxEvents=0) |
|
Bool_t | readTuple (TTree *tree, const char *var1, const char *var2, const char *var3, const char *var4, const char *var5, const char *var6, UInt_t maxEvents=0) |
|
Double_t | density (std::vector< Double_t > &x) |
| Calculate PDF value at the given point. More...
|
|
AbsPhaseSpace * | phaseSpace () |
| Return phase space definition for this PDF. More...
|
|
| 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...
|
|
Class that describes the unbinned kernel density.
Definition at line 19 of file KernelDensity.hh.
KernelDensity::KernelDensity |
( |
const char * |
pdfName, |
|
|
AbsPhaseSpace * |
thePhaseSpace, |
|
|
std::vector< Double_t > & |
width, |
|
|
UInt_t |
approxSize, |
|
|
AbsDensity * |
approxDensity = 0 |
|
) |
| |
KernelDensity::KernelDensity |
( |
const char * |
pdfName, |
|
|
AbsPhaseSpace * |
thePhaseSpace, |
|
|
UInt_t |
approxSize, |
|
|
Double_t |
width1, |
|
|
Double_t |
width2 = 0 , |
|
|
Double_t |
width3 = 0 , |
|
|
Double_t |
width4 = 0 , |
|
|
Double_t |
width5 = 0 |
|
) |
| |
KernelDensity::KernelDensity |
( |
const char * |
pdfName, |
|
|
AbsPhaseSpace * |
thePhaseSpace, |
|
|
UInt_t |
approxSize, |
|
|
AbsDensity * |
approxDensity, |
|
|
Double_t |
width1, |
|
|
Double_t |
width2 = 0 , |
|
|
Double_t |
width3 = 0 , |
|
|
Double_t |
width4 = 0 , |
|
|
Double_t |
width5 = 0 |
|
) |
| |
KernelDensity::~KernelDensity |
( |
| ) |
|
|
virtual |
Int_t KernelDensity::cellIndex |
( |
std::vector< Double_t > & |
x | ) |
|
|
private |
Double_t KernelDensity::density |
( |
std::vector< Double_t > & |
x | ) |
|
|
virtual |
Calculate PDF value at the given point.
- Parameters
-
[in] | x | the point at which to calculate the PDF |
- Returns
- PDF value
Implements AbsDensity.
Definition at line 485 of file KernelDensity.cpp.
Bool_t KernelDensity::generateApproximation |
( |
UInt_t |
approxSize | ) |
|
UInt_t KernelDensity::numCells |
( |
void |
| ) |
|
|
private |
Double_t KernelDensity::rawDensity |
( |
std::vector< Double_t > & |
x, |
|
|
std::vector< TCell > & |
vector |
|
) |
| |
|
private |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
std::vector< TString > & |
vars, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
const char * |
var2, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
const char * |
var2, |
|
|
const char * |
var3, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
const char * |
var2, |
|
|
const char * |
var3, |
|
|
const char * |
var4, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
const char * |
var2, |
|
|
const char * |
var3, |
|
|
const char * |
var4, |
|
|
const char * |
var5, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
Bool_t KernelDensity::readTuple |
( |
TTree * |
tree, |
|
|
const char * |
var1, |
|
|
const char * |
var2, |
|
|
const char * |
var3, |
|
|
const char * |
var4, |
|
|
const char * |
var5, |
|
|
const char * |
var6, |
|
|
UInt_t |
maxEvents = 0 |
|
) |
| |
void KernelDensity::setWidth |
( |
std::vector< Double_t > & |
width | ) |
|
std::vector<TCell> KernelDensity::m_apprVector |
|
private |
std::vector<TCell> KernelDensity::m_dataVector |
|
private |
The documentation for this class was generated from the following files:
|