MINT2
DalitzEventBase.h
Go to the documentation of this file.
1 #ifndef DALITZ_EVENT_BASE_HH
2 #define DALITZ_EVENT_BASE_HH
3 // author: Jonas Rademacker (Jonas.Rademacker@bristol.ac.uk)
4 // status: Mon 9 Feb 2009 19:18:00 GMT
5 
6 #include <iostream>
7 #include <string>
8 #include <vector>
9 
10 #include "TLorentzVector.h"
11 #include "TVector3.h"
12 
13 #include "Mint/IDalitzEvent.h"
16 
17 #include "Mint/RememberAnything.h"
18 
19 #include "Mint/Permutator.h"
20 
21 class TRandom;
22 class TNtupleD;
23 
24 class DalitzEventBase : virtual public IDalitzEvent{
25  protected:
26 
27  static const char ntpNameChars[];
28  static const char prtNameChars[];
29  static std::string prtToNtpName(const std::string& s_in);
30  static std::string ntpToPrtName(const std::string& s_in);
31 
32  static int singleParticleNtpArraySize();
33 
34  // begin all data members:
36  std::vector<TLorentzVector> _p;
37 
38  mutable double _rememberPhaseSpace;
39 
40  RememberAnything<std::complex<double> > _rememberAmps;
41 
42  double _aValue;
43  double _weight;
45 
46  mutable std::vector< std::vector<double> > _s;
47  mutable std::vector< std::vector<double> > _t;
48 
51 
52  // end all data members;
53 
54 
55  const TLorentzVector& p_intern(unsigned int i) const; // 4-vectors
56  double s_intern(unsigned int i, unsigned int j) const; // sij
57  double t_intern(unsigned int i, unsigned int j) const; // tij
58 
59  bool shoutAndKill();
60  double BDet() const;
61  double Delta4() const;
62  double phaseSpace3(double epsilon = 1.e-9) const;
63  double phaseSpace4() const;
64  bool setMomenta(const std::vector<TLorentzVector>& mumAndDgtr_p4);
65  bool setMomenta(const std::vector<TVector3>& mumAndDgtr_p3);
66  void setP(unsigned int i, const TLorentzVector& p4);
67 public:
68  DalitzEvent();
70  , const std::vector<TLorentzVector> mumAndDgtr_p4
71  );
73  , double t01_in
74  , double s12_in, double s23_in, double s34_in
75  , double t40_in
76  );
78  , const std::vector<TVector3> mumAndDgtr_p4
79  );
81  );
83  , TRandom* rnd
84  );
85  DalitzEvent(const DalitzEvent& other);
86 
87  virtual double getWeight() const;
88  virtual void setWeight(double w);
89 
90  virtual void setGeneratorPdfRelativeToPhaseSpace(double gpdf);
91  virtual double getGeneratorPdfRelativeToPhaseSpace()const;
92 
93  bool resetST();
94 
96  return _perm[_permutationIndex];
97  }
98 
99  void setPermutationIndex(int i);
100  int numPermutations() const;
101  int permutationIndex() const{
102  return _permutationIndex;
103  }
104 
105 
106  DalitzEvent(TNtupleD* ntp);
107 
108  bool fromNtuple(TNtupleD* ntp);
109 
110  // required
111  virtual void setAValue(double a){_aValue=a;}
112  virtual double getAValue()const{return _aValue;}
113 
114  virtual const DalitzEventPattern& eventPattern() const{return _pat;}
115  virtual const TLorentzVector& p(unsigned int i) const; // 4-vectors
116  virtual double s(unsigned int i, unsigned int j) const; // sij
117  virtual double t(unsigned int i, unsigned int j) const; // tij
118  virtual double sij(const std::vector<int>& indices) const;
119 
120  virtual double phaseSpace() const;
121 
122  virtual bool retrieveComplex(void* key, std::complex<double>& value){
123  return _rememberAmps.find(key, value);
124  }
125 
126  virtual void setComplex(void* key, const std::complex<double>& value){
127  _rememberAmps.set(key, value);
128  }
129 
130  // helpful
131  double sijMin(const std::vector<int>& indices) const;
132  double sijMax(const std::vector<int>& indices) const;
133 
134  int nDgtr() const{
135  return _pat.size()-1;
136  }
137  double m(unsigned int i) const; // nominal masses (from _pat)
138  double m2(unsigned int i) const; // nominal masses squared
139  TLorentzVector p_dgtr_sum() const;
140  bool kinematicallyAllowed(double epsilon = 1.e-9) const;
141 
142  void print(std::ostream& os = std::cout) const;
143  // generate yourself:
144  void generateThisToPhaseSpace(TRandom* rnd = 0);
145  void generateThisToPhaseSpace(double p_of_D, TRandom* rnd = 0);
146 
147  void generateThisToPhaseSpace(const TVector3& pVec_of_D
148  , TRandom* rnd = 0);
149 
150 
151  std::string makeNtupleVarnames()const;
152  bool fillNtupleVarArray(Double_t* array, unsigned int arraySize) const;
153  unsigned int ntupleVarArraySize() const;
154 
155  // mainly for debugging and x-check with prev version
157 
158 
159 };
160 
161 std::ostream& operator<<(std::ostream& os, const DalitzEvent& de);
162 
163 #endif
164 //
void print(std::ostream &os=std::cout) const
double phaseSpace4() const
int numPermutations() const
virtual void setGeneratorPdfRelativeToPhaseSpace(double gpdf)
std::string makeNtupleVarnames() const
int permutationIndex() const
RememberAnything< std::complex< double > > _rememberAmps
virtual double s(unsigned int i, unsigned int j) const
virtual const DalitzEventPattern & eventPattern() const
DalitzEventPattern _pat
std::vector< TLorentzVector > _p
bool kinematicallyAllowed(double epsilon=1.e-9) const
static const char ntpNameChars[]
Permutator _perm
double sijMin(const std::vector< int > &indices) const
virtual double getAValue() const
void setPermutationIndex(int i)
int nDgtr() const
static std::string prtToNtpName(const std::string &s_in)
double Delta4() const
const TLorentzVector & p_intern(unsigned int i) const
std::vector< std::vector< double > > _s
double _rememberPhaseSpace
virtual double phaseSpace() const
virtual const TLorentzVector & p(unsigned int i) const
double BDet() const
const Permutation & currentPermutation() const
std::vector< std::vector< double > > _t
double s_intern(unsigned int i, unsigned int j) const
Calculate4BodyProps makeCalculate4BodyProps() const
bool fillNtupleVarArray(Double_t *array, unsigned int arraySize) const
TLorentzVector p_dgtr_sum() const
static int singleParticleNtpArraySize()
virtual double getWeight() const
bool setMomenta(const std::vector< TLorentzVector > &mumAndDgtr_p4)
unsigned int ntupleVarArraySize() const
static const char prtNameChars[]
double _generatorPdfRelativeToPhaseSpace
bool fromNtuple(TNtupleD *ntp)
static std::string ntpToPrtName(const std::string &s_in)
double t_intern(unsigned int i, unsigned int j) const
unsigned int size() const
void setP(unsigned int i, const TLorentzVector &p4)
double m(unsigned int i) const
std::ostream & operator<<(std::ostream &os, const DalitzEvent &de)
double phaseSpace3(double epsilon=1.e-9) const
double m2(unsigned int i) const
virtual void setAValue(double a)
virtual double getGeneratorPdfRelativeToPhaseSpace() const
virtual double t(unsigned int i, unsigned int j) const
virtual void setWeight(double w)
virtual bool retrieveComplex(void *key, std::complex< double > &value)
virtual double sij(const std::vector< int > &indices) const
double sijMax(const std::vector< int > &indices) const
virtual void setComplex(void *key, const std::complex< double > &value)
void generateThisToPhaseSpace(TRandom *rnd=0)