meerkat is hosted by Hepforge, IPPP Durham
Meerkat  v1r3
Multidimensional kernel density estimation package
AbsDensity.hh
Go to the documentation of this file.
1 #ifndef ABS_DENSITY
2 #define ABS_DENSITY
3 
4 #include "AbsPhaseSpace.hh"
5 
6 #include "TMath.h"
7 #include "TH1F.h"
8 #include "TH2F.h"
9 #include "TNtuple.h"
10 #include "TRandom3.h"
11 
12 #include <vector>
13 
15 
16 class AbsDensity {
17 
18  public:
19 
21 
24  AbsDensity(const char* pdfName);
25 
27  virtual ~AbsDensity();
28 
30 
34  virtual Double_t density(std::vector<Double_t> &x) = 0;
35 
37 
40  virtual AbsPhaseSpace* phaseSpace() = 0;
41 
43 
48  void slice(std::vector<Double_t> &x, UInt_t num, TH1F* hist);
49 
51 
58  void slice(std::vector<Double_t> &x, UInt_t numx, UInt_t numy, TH2F* hist, Bool_t inPhaseSpace = true);
59 
60  double transform(TH1F* hist1, TH1F* hist2, double x);
61 
63 
66  void project(TH1F* hist);
67 
69 
73  void project(TH2F* hist, Bool_t inPhaseSpace = true);
74 
76 
79  void setMajorant(Double_t majorant) { m_majorant = majorant; }
80 
82 
85  void setMaxTries(UInt_t maxTries) { m_maxTries = maxTries; }
86 
88 
92  Double_t generate(std::vector<Double_t> &x);
93 
95 
100  void generate(TNtuple* tree, UInt_t numEvents);
101 
103 
106  const char* name(void) { return m_name; }
107 
109 
113  void setSeed(UInt_t seed = 0) { m_rnd.SetSeed(seed); }
114 
115  protected:
116 
118  char m_name[256];
119 
121  Double_t m_majorant;
122 
124  UInt_t m_maxTries;
125 
127  TRandom3 m_rnd;
128 };
129 
130 #endif
const char * name(void)
Return the name of the PDF.
Definition: AbsDensity.hh:106
void project(TH1F *hist)
Calculate projection of the 1D PDF.
Definition: AbsDensity.cpp:76
double transform(TH1F *hist1, TH1F *hist2, double x)
Definition: AbsDensity.cpp:179
void setSeed(UInt_t seed=0)
Set random seed.
Definition: AbsDensity.hh:113
Abstract class which defines probability density interface.
Definition: AbsDensity.hh:16
virtual AbsPhaseSpace * phaseSpace()=0
Return phase space definition for this PDF.
virtual ~AbsDensity()
Destructor.
Definition: AbsDensity.cpp:21
char m_name[256]
PDF name.
Definition: AbsDensity.hh:118
TRandom3 m_rnd
Random number generator.
Definition: AbsDensity.hh:127
Abstract class which defines phase space interface.
void setMajorant(Double_t majorant)
Set majorant for accept-reject method.
Definition: AbsDensity.hh:79
UInt_t m_maxTries
Maximum number of tries for accept-reject method.
Definition: AbsDensity.hh:124
AbsDensity(const char *pdfName)
Constructor.
Definition: AbsDensity.cpp:14
Double_t generate(std::vector< Double_t > &x)
Generate a single point within the phase space according to the PDF using accept-reject method...
Definition: AbsDensity.cpp:99
void setMaxTries(UInt_t maxTries)
Set maximum number of tries for accept-reject method.
Definition: AbsDensity.hh:85
virtual Double_t density(std::vector< Double_t > &x)=0
Calculate PDF value at the given point.
Double_t m_majorant
PDF majorant (maximum PDF value needed for accept-reject).
Definition: AbsDensity.hh:121
void slice(std::vector< Double_t > &x, UInt_t num, TH1F *hist)
Calculate 1D slice of the PDF.
Definition: AbsDensity.cpp:25