meerkat is hosted by Hepforge, IPPP Durham
Meerkat  v1r3
Multidimensional kernel density estimation package
BinnedKernelDensity.hh
Go to the documentation of this file.
1 #ifndef BINNED_KERNEL_DENSITY
2 #define BINNED_KERNEL_DENSITY
3 
4 #include "AbsDensity.hh"
5 
6 #include "TMath.h"
7 
8 #include <vector>
9 
11 
12  public:
13 
15 
27  BinnedKernelDensity(const char* pdfName,
28  AbsPhaseSpace* thePhaseSpace,
29  TTree* tree,
30  const char *vars1,
31  UInt_t bins1,
32  Double_t width1,
33  AbsDensity* approx = 0,
34  UInt_t toyEvents = 0,
35  UInt_t maxEvents = 0,
36  UInt_t skipEvents = 0
37  );
38 
40 
53  BinnedKernelDensity(const char* pdfName,
54  AbsPhaseSpace* thePhaseSpace,
55  TTree* tree,
56  const char *vars1,
57  const char *weight,
58  UInt_t bins1,
59  Double_t width1,
60  AbsDensity* approx = 0,
61  UInt_t toyEvents = 0,
62  UInt_t maxEvents = 0,
63  UInt_t skipEvents = 0
64  );
65 
67 
82  BinnedKernelDensity(const char* pdfName,
83  AbsPhaseSpace* thePhaseSpace,
84  TTree* tree,
85  const char *vars1, const char *vars2,
86  UInt_t bins1, UInt_t bins2,
87  Double_t width1, Double_t width2,
88  AbsDensity* approx = 0,
89  UInt_t toyEvents = 0,
90  UInt_t maxEvents = 0,
91  UInt_t skipEvents = 0
92  );
93 
95 
111  BinnedKernelDensity(const char* pdfName,
112  AbsPhaseSpace* thePhaseSpace,
113  TTree* tree,
114  const char *vars1, const char *vars2,
115  const char* weight,
116  UInt_t bins1, UInt_t bins2,
117  Double_t width1, Double_t width2,
118  AbsDensity* approx = 0,
119  UInt_t toyEvents = 0,
120  UInt_t maxEvents = 0,
121  UInt_t skipEvents = 0
122  );
123 
125 
144  BinnedKernelDensity(const char* pdfName,
145  AbsPhaseSpace* thePhaseSpace,
146  TTree* tree,
147  const char *vars1, const char *vars2, const char *vars3,
148  UInt_t bins1, UInt_t bins2, UInt_t bins3,
149  Double_t width1, Double_t width2, Double_t width3,
150  AbsDensity* approx = 0,
151  UInt_t toyEvents = 0,
152  UInt_t maxEvents = 0,
153  UInt_t skipEvents = 0
154  );
155 
157 
177  BinnedKernelDensity(const char* pdfName,
178  AbsPhaseSpace* thePhaseSpace,
179  TTree* tree,
180  const char *vars1, const char *vars2, const char *vars3,
181  const char *weight,
182  UInt_t bins1, UInt_t bins2, UInt_t bins3,
183  Double_t width1, Double_t width2, Double_t width3,
184  AbsDensity* approx = 0,
185  UInt_t toyEvents = 0,
186  UInt_t maxEvents = 0,
187  UInt_t skipEvents = 0
188  );
189 
191 
212  BinnedKernelDensity(const char* pdfName,
213  AbsPhaseSpace* thePhaseSpace,
214  TTree* tree,
215  const char *vars1, const char *vars2, const char *vars3, const char *vars4,
216  UInt_t bins1, UInt_t bins2, UInt_t bins3, UInt_t bins4,
217  Double_t width1, Double_t width2, Double_t width3, Double_t width4,
218  AbsDensity* approx = 0,
219  UInt_t toyEvents = 0,
220  UInt_t maxEvents = 0,
221  UInt_t skipEvents = 0
222  );
223 
225 
247  BinnedKernelDensity(const char* pdfName,
248  AbsPhaseSpace* thePhaseSpace,
249  TTree* tree,
250  const char *vars1, const char *vars2, const char *vars3, const char *vars4,
251  const char *weight,
252  UInt_t bins1, UInt_t bins2, UInt_t bins3, UInt_t bins4,
253  Double_t width1, Double_t width2, Double_t width3, Double_t width4,
254  AbsDensity* approx = 0,
255  UInt_t toyEvents = 0,
256  UInt_t maxEvents = 0,
257  UInt_t skipEvents = 0
258  );
259 
261 
285  BinnedKernelDensity(const char* pdfName,
286  AbsPhaseSpace* thePhaseSpace,
287  TTree* tree,
288  const char *vars1, const char *vars2, const char *vars3, const char *vars4, const char *vars5,
289  UInt_t bins1, UInt_t bins2, UInt_t bins3, UInt_t bins4, UInt_t bins5,
290  Double_t width1, Double_t width2, Double_t width3, Double_t width4, Double_t width5,
291  AbsDensity* approx = 0,
292  UInt_t toyEvents = 0,
293  UInt_t maxEvents = 0,
294  UInt_t skipEvents = 0
295  );
296 
298 
323  BinnedKernelDensity(const char* pdfName,
324  AbsPhaseSpace* thePhaseSpace,
325  TTree* tree,
326  const char *vars1, const char *vars2, const char *vars3, const char *vars4, const char *vars5,
327  const char *weight,
328  UInt_t bins1, UInt_t bins2, UInt_t bins3, UInt_t bins4, UInt_t bins5,
329  Double_t width1, Double_t width2, Double_t width3, Double_t width4, Double_t width5,
330  AbsDensity* approx = 0,
331  UInt_t toyEvents = 0,
332  UInt_t maxEvents = 0,
333  UInt_t skipEvents = 0
334  );
335 
337 
349  BinnedKernelDensity(const char* pdfName,
350  AbsPhaseSpace* thePhaseSpace,
351  TTree* tree,
352  std::vector<TString> &vars,
353  std::vector<UInt_t> &binning,
354  std::vector<Double_t> &width,
355  AbsDensity* approx = 0,
356  UInt_t toyEvents = 0,
357  UInt_t maxEvents = 0,
358  UInt_t skipEvents = 0
359  );
360 
362  virtual ~BinnedKernelDensity();
363 
366 
369  void writeToFile(const char* fileName);
370 
372 
375  void writeToTextFile(const char* fileName);
376 
378 
381  void writeToRootFile(const char* fileName);
382 
384 
388  Double_t density(std::vector<Double_t> &x);
389 
391 
395 
397 
400  void setFractionalMode(Bool_t mode = true) { m_fractionalMode = mode; }
401 
403  void normalise(void);
404 
405  private:
406 
408 
419  void init(AbsPhaseSpace* thePhaseSpace,
420  TTree* tree,
421  std::vector<TString> &vars,
422  std::vector<UInt_t> &binning,
423  std::vector<Double_t> &width,
424  AbsDensity* approx = 0,
425  UInt_t toyEvents = 0,
426  UInt_t maxEvents = 0,
427  UInt_t skipEvents = 0
428  );
429 
431 
435  UInt_t iterToIndex( std::vector<UInt_t> &iter );
436 
438 
443  void addToMap(std::vector<Double_t> &map, std::vector<Double_t> &point, Double_t weight = 1.);
444 
445  void fillMapFromTree( TTree* tree, std::vector<TString> &vars, UInt_t maxEvents = 0, UInt_t skipEvents = 0);
446 
447  void fillMapFromDensity(AbsDensity* density, UInt_t toyEvents = 0);
448 
450 
454  Double_t mapDensity(std::vector<Double_t> &map, std::vector<Double_t> &x);
455 
457  std::vector<Double_t> m_map;
458 
460  std::vector<Double_t> m_approxMap;
461 
463  std::vector<UInt_t> m_binning;
464 
467 
470 
472  std::vector<Double_t> m_width;
473 
475  UInt_t m_dim;
476 
479 
480 };
481 
482 #endif
void writeToRootFile(const char *fileName)
Write the binned density into a ROOT file.
AbsPhaseSpace * m_phaseSpace
Reference to phase space.
std::vector< Double_t > m_width
Vector of kernel widths.
Abstract class which defines probability density interface.
Definition: AbsDensity.hh:16
std::vector< Double_t > m_map
Bin map of estimated PDF.
BinnedKernelDensity(const char *pdfName, AbsPhaseSpace *thePhaseSpace, TTree *tree, const char *vars1, UInt_t bins1, Double_t width1, AbsDensity *approx=0, UInt_t toyEvents=0, UInt_t maxEvents=0, UInt_t skipEvents=0)
Constructor for 1-dimensional kernel PDF with binned interpolation from the sample of points in an NT...
void normalise(void)
Normalise the PDF such that the average PDF value over the allowed phase space equals to 1...
Bool_t m_fractionalMode
Fractional mode flag.
std::vector< UInt_t > m_binning
Vector of bin numbers in each variable.
UInt_t iterToIndex(std::vector< UInt_t > &iter)
Convert an N-dimensional iterator vector into a linear bin index in the bin map.
Double_t density(std::vector< Double_t > &x)
Calculate PDF density at the point. Multilinear binned interpolation is used which is calculated in t...
Abstract class which defines phase space interface.
void fillMapFromTree(TTree *tree, std::vector< TString > &vars, UInt_t maxEvents=0, UInt_t skipEvents=0)
void addToMap(std::vector< Double_t > &map, std::vector< Double_t > &point, Double_t weight=1.)
Add a kernel density of a single data point to the binned map.
virtual ~BinnedKernelDensity()
Destructor.
void init(AbsPhaseSpace *thePhaseSpace, TTree *tree, std::vector< TString > &vars, std::vector< UInt_t > &binning, std::vector< Double_t > &width, AbsDensity *approx=0, UInt_t toyEvents=0, UInt_t maxEvents=0, UInt_t skipEvents=0)
Common initialise method used by all constructors.
Double_t mapDensity(std::vector< Double_t > &map, std::vector< Double_t > &x)
Calculate the raw density using the binned map (estimated or approximation) at a given point...
void writeToTextFile(const char *fileName)
Write the binned PDF to a text file.
void writeToFile(const char *fileName)
Write density map to a file, ROOT or text depending on extension.
void fillMapFromDensity(AbsDensity *density, UInt_t toyEvents=0)
AbsPhaseSpace * phaseSpace()
Return the phase space.
UInt_t m_dim
Cached dimensionality of the phase space.
std::vector< Double_t > m_approxMap
Bin map of approximation PDF convolved with the kernel.
void setFractionalMode(Bool_t mode=true)
Set fractional mode of operation.
AbsDensity * m_approxDensity
Reference to approximation PDF.