MINT2
Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
DalitzEvent Class Reference

#include <DalitzEvent.h>

Inheritance diagram for DalitzEvent:
IDalitzEvent MINT::IWeightedEvent TimeDependentGenerator::GenTimeEvent TimeDependentGeneratorOld::GenTimeEvent

Public Member Functions

 DalitzEvent ()
 
 DalitzEvent (const DalitzEventPattern &pat, const std::vector< TLorentzVector > &mumAndDgtr_p4)
 
 DalitzEvent (const DalitzEventPattern &pat, double t01_in, double s12_in, double s23_in, double s34_in, double t40_in)
 Constructor for a 4-body decay at a given point in phase space. More...
 
 DalitzEvent (const DalitzEventPattern &pat, const double s13, const double s23)
 Constructor for a 3-body decay at a given point in phase space. More...
 
 DalitzEvent (const DalitzEventPattern &pat, const std::vector< TVector3 > &mumAndDgtr_p3)
 
 DalitzEvent (const DalitzEventPattern &pat)
 
 DalitzEvent (const DalitzEventPattern &pat, TRandom *rnd)
 
 DalitzEvent (TNtupleD *ntp)
 
 DalitzEvent (const IDalitzEvent &ievt)
 
 DalitzEvent (const IDalitzEvent *ievt)
 
 DalitzEvent (const DalitzEvent *other)
 
 DalitzEvent (const DalitzEvent &other)
 
virtual ~DalitzEvent ()
 
IDalitzEventclone () const
 
virtual double getWeight () const
 
virtual void setWeight (double w)
 
virtual void setGeneratorPdfRelativeToPhaseSpace (double gpdf)
 
virtual double getGeneratorPdfRelativeToPhaseSpace () const
 
virtual const std::vector< double > & getVectorOfValues () const
 
virtual std::vector< double > & getVectorOfValues ()
 
virtual const std::vector< double > & getVectorOfWeights () const
 
virtual std::vector< double > & getVectorOfWeights ()
 
virtual void setValueInVector (unsigned int i, double value)
 
virtual void setWeightInVector (unsigned int i, double weight)
 
virtual double getValueFromVector (unsigned int i) const
 
virtual double getWeightFromVector (unsigned int i) const
 
virtual void P_conjugateYourself ()
 
virtual void C_conjugateYourself ()
 
virtual void CP_conjugateYourself ()
 
bool resetST ()
 
const PermutationcurrentPermutation () const
 
void setPermutationIndex (int i)
 
int numPermutations () const
 
int permutationIndex () const
 
bool fromTree (TTree *tree)
 
bool fromNtuple (TNtupleD *ntp)
 
bool fromParasTree (TTree *ntp)
 
bool fromParasTreeOld (TTree *ntp)
 
virtual void setAValue (double a)
 
virtual double getAValue () const
 
virtual const DalitzEventPatterneventPattern () const
 
virtual const TLorentzVector & p (unsigned int i) const
 
virtual TLorentzVector & p (unsigned int i)
 
virtual double s (unsigned int i, unsigned int j) const
 
virtual double t (unsigned int i, unsigned int j) const
 
virtual double sij (const std::vector< int > &indices) const
 
virtual double sij (const MINT::PolymorphVector< int > &indices) const
 
virtual void setMothers3Momentum (const TVector3 &mp3)
 
virtual double phaseSpace () const
 
virtual bool retrieveValue (int i, std::complex< double > &value, long int configNumber=1)
 
virtual void setValue (int i, const std::complex< double > &value, long int configNumber=1)
 
virtual bool retrieveValue (int i, double value, long int configNumber=1)
 
virtual void setValue (int i, double value, long int configNumber=1)
 
double sijMin (const std::vector< int > &indices) const
 
double sijMin (const MINT::PolymorphVector< int > &indices) const
 
double sijMax (const std::vector< int > &indices) const
 
double sijMax (const MINT::PolymorphVector< int > &indices) const
 
int nDgtr () const
 
double m (unsigned int i) const
 
double m2 (unsigned int i) const
 
TLorentzVector p_dgtr_sum () const
 
bool kinematicallyAllowed (double epsilon=1.e-9) const
 
void print (std::ostream &os=std::cout) const
 
void generateThisToPhaseSpace (TRandom *rnd=0)
 
void generateThisToPhaseSpace (double p_of_D, TRandom *rnd=0)
 
void generateThisToPhaseSpace (const TVector3 &pVec_of_D, TRandom *rnd=0)
 
std::string makeNtupleVarnames () const
 
bool fillNtupleVarArray (std::vector< Double_t > &array) const
 
unsigned int ntupleVarArraySize () const
 
Calculate4BodyProps makeCalculate4BodyProps () const
 
- Public Member Functions inherited from IDalitzEvent
virtual ~IDalitzEvent ()
 

Static Public Member Functions

static long int assignUniqueRememberNumber ()
 
static long int eventCounter ()
 

Protected Member Functions

const TLorentzVector & p_intern (unsigned int i) const
 
TLorentzVector & p_intern (unsigned int i)
 
double s_intern (unsigned int i, unsigned int j) const
 
double t_intern (unsigned int i, unsigned int j) const
 
double evalsij_intern (const std::vector< int > &intern_indices) const
 
double sijMap_intern (const std::vector< int > &intern_indices) const
 
double sijMap (const std::vector< int > &indices) const
 
bool shoutAndKill ()
 
double BDet () const
 
double Delta4 () const
 
double phaseSpace3 (double epsilon=1.e-9) const
 
double phaseSpace4 () const
 
bool setMomenta (const std::vector< TLorentzVector > &mumAndDgtr_p4)
 
bool setMomenta (const std::vector< TVector3 > &mumAndDgtr_p3)
 
void setP (unsigned int i, const TLorentzVector &p4)
 
- Protected Member Functions inherited from IDalitzEvent
 IDalitzEvent ()
 

Static Protected Member Functions

static std::string prtToNtpName (const std::string &s_in)
 
static std::string ntpToPrtName (const std::string &s_in)
 
static int singleParticleNtpArraySize ()
 
static bool parseNtpEntryName (const std::string &entry, std::string &part1, std::string &part2)
 

Protected Attributes

DalitzEventPattern _pat
 
std::vector< TLorentzVector > _p
 
double _rememberPhaseSpace
 
RememberAnythingFast< std::complex< double > > _rememberComplexFast
 
RememberAnythingFast< double > _rememberDoubleFast
 
double _aValue
 
double _weight
 
double _generatorPdfRelativeToPhaseSpace
 
std::vector< double > _vectorOfValues
 
std::vector< double > _vectorOfWeights
 
std::vector< std::vector< double > > _s
 
std::vector< std::vector< double > > _t
 
std::map< std::vector< int >, double > _sijMap
 
int _permutationIndex
 
Permutator _perm
 

Static Protected Attributes

static const char ntpNameChars [] = { '#', '~', '\0' }
 
static const char prtNameChars [] = { '+', '-', '\0' }
 
static long int _eventCounter =0
 
static long int _rememberVectorCounter =0
 

Detailed Description

Definition at line 27 of file DalitzEvent.h.

Constructor & Destructor Documentation

◆ DalitzEvent() [1/12]

DalitzEvent::DalitzEvent ( )

Definition at line 33 of file DalitzEvent.cpp.

34  : _rememberPhaseSpace(-9999.)
35  , _aValue(-9999.)
36  , _weight(1)
38  , _vectorOfValues()
41  , _perm()
42 {
43  _eventCounter++;
44 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
Permutator _perm
Definition: DalitzEvent.h:63
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [2/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat,
const std::vector< TLorentzVector > &  mumAndDgtr_p4 
)

Definition at line 45 of file DalitzEvent.cpp.

48  : _pat(pat)
49  , _rememberPhaseSpace(-9999.)
50  , _aValue(-9999.)
51  , _weight(1)
53  , _vectorOfValues()
55  , _s(pat.size())
56  , _t(pat.size())
58  , _perm(pat)
59 {
60  setMomenta(mumAndDgtrs_p4);
61  if(_pat.size() < 3 || _pat.size() != _p.size()) shoutAndKill();
62  resetST();
63  _eventCounter++;
64 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool setMomenta(const std::vector< TLorentzVector > &mumAndDgtr_p4)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [3/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat,
double  t01_in,
double  s12_in,
double  s23_in,
double  s34_in,
double  t40_in 
)

Constructor for a 4-body decay at a given point in phase space.

Definition at line 66 of file DalitzEvent.cpp.

71  : _pat(pat)
72  , _rememberPhaseSpace(-9999.)
73  , _aValue(-9999.)
74  , _weight(1)
76  , _vectorOfValues()
78  , _s(pat.size())
79  , _t(pat.size())
81  , _perm(pat)
82 {
83  bool dbThis=false;
84 
85  if(_pat.size() != 5){
86  cout << "DalitzEvent::DalitzEvent ERROR "
87  << " event pattern " << pat
88  << " doesn't fit 4-body constructor"
89  << endl;
90  throw "shit happens";
91  }
92  Calculate4BodyProps c4bp(t01_in, s12_in, s23_in, s34_in, t40_in
93  , _pat[0].mass()
94  , _pat[1].mass(), _pat[2].mass()
95  , _pat[3].mass(), _pat[4].mass()
96  );
97  if(dbThis){
98  cout << " c4bp phase space " << c4bp.phaseSpaceFactor() << endl;
99  }
100  std::vector<TLorentzVector> p4List(5);
101  if(c4bp.phaseSpaceFactor() > 0){
102  for(int i=0; i<=4; i++){
103  p4List[i] = c4bp.pVec(i);
104  if(dbThis)cout << "p4List[" << i << "] " << p4List[i] << endl;
105  }
106  }else{
107  for(int i=0; i<=4; i++){
108  TLorentzVector p4(-9999, -9999, -9999, -9999);
109  p4List[i] = p4;
110  }
111  }
112  setMomenta(p4List);
113  if(_pat.size() < 3 || _pat.size() != _p.size()) shoutAndKill();
114  resetST();
115 
116  if(dbThis){
117  cout << "constructed event from this pattern: " << _pat << "\n and "
118  << t01_in << ", " << s12_in << ", " << s23_in << ", " << s34_in << ", " << t40_in
119  << endl;
120 
121  cout << " and got out: "
122  << t_intern(0,1) << ", " << s_intern(1,2) << ", " << s_intern(2,3)
123  << ", " << s_intern(3,4) << ", " << t_intern(4,0)
124  << endl;
125  }
126  _eventCounter++;
127 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool setMomenta(const std::vector< TLorentzVector > &mumAndDgtr_p4)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
double t_intern(unsigned int i, unsigned int j) const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double s_intern(unsigned int i, unsigned int j) const
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [4/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat,
const double  s13,
const double  s23 
)

Constructor for a 3-body decay at a given point in phase space.

Definition at line 129 of file DalitzEvent.cpp.

132  : _pat(pat)
133  , _rememberPhaseSpace(-9999.)
134  , _aValue(-9999.)
135  , _weight(1)
137  , _vectorOfValues()
138  , _vectorOfWeights()
139  , _s(pat.size())
140  , _t(pat.size())
141  , _permutationIndex(0)
142  , _perm(pat)
143 {
144  bool dbThis=false;
145 
146  if(_pat.size() != 4){
147  cout << "DalitzEvent::DalitzEvent ERROR "
148  << " event pattern " << pat
149  << " doesn't fit 4-body constructor"
150  << endl;
151  throw "shit happens";
152  }
153 
154  double m0 = pat[0].mass() ;
155  double m1 = pat[1].mass() ;
156  double m2 = pat[2].mass() ;
157  double m3 = pat[3].mass() ;
158 
159  // Calculate energies of pi+, pi- from known parameters
160  double E1 = ( m0 * m0 + m1 * m1 - s23 ) / (2 * m0) ;
161  double E2 = ( m0 * m0 + m2 * m2 - s13 ) / (2 * m0) ;
162 
163  // Get E3 from conservation of energy
164  double E3 = m0 - E1 - E2 ;
165 
166  // Reject events where kinematics are impossible
167  double pz1Sq = E1 * E1 - m1 * m1 ;
168  if(pz1Sq < 0){
169  cout << "Invalid phase space point for the given pattern (s13, s23) = ("
170  << s13 << ", " << s23 << ")" << endl ;
171  throw out_of_range("") ;
172  }
173 
174  // Get pz1 from mass/energy relation (co-ordinate system defined so that px1 = py1 = 0)
175  double pz1 = sqrt(pz1Sq) ;
176 
177  // Calculate pz3 from conservation of momentum
178  double pz3 = ( E2 * E2 - m2 * m2 - E3 * E3 + m3 * m3 - pz1 * pz1 ) / ( 2*pz1 ) ;
179 
180  // Use pz3 to get py3 (co-ordinate system designed so that px3 = 0)
181  double py3Sq = E3 * E3 - m3 * m3 - pz3 * pz3 ;
182 
183  // Reject events where kinematics are impossible
184  if(py3Sq < 0){
185  cout << "Invalid phase space point for the given pattern (s13, s23) = ("
186  << s13 << ", " << s23 << ")" << endl ;
187  throw out_of_range("") ;
188  }
189 
190  double py3 = sqrt(py3Sq) ;
191  // Conservation of momentum => py2 = -py3 and pz2 = -(pz1 + pz3)
192  double py2 = -1 * py3 ;
193  double pz2 = -1 * (pz1 + pz3) ;
194  std::vector<TLorentzVector> p4List(4, TLorentzVector(0, 0, 0, m0));
195  p4List[1] = TLorentzVector(0., 0., pz1, E1) ;
196  p4List[2] = TLorentzVector(0., py2, pz2, E2) ;
197  p4List[3] = TLorentzVector(0., py3, pz3, E3) ;
198 
199  setMomenta(p4List);
200  if(_pat.size() < 3 || _pat.size() != _p.size()) shoutAndKill();
201  resetST();
202 
203  if(!kinematicallyAllowed()){
204  cout << "Event not kinematically allowed! Was given (s13, s23) = ("
205  << s13 << ", " << s23 << "), got (" << s_intern(1, 3) << ", "
206  << s_intern(2, 3) << ")" << endl ;
207  throw out_of_range("") ;
208  }
209  _eventCounter++;
210 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool setMomenta(const std::vector< TLorentzVector > &mumAndDgtr_p4)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
double m2(unsigned int i) const
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double s_intern(unsigned int i, unsigned int j) const
static const double m3
bool kinematicallyAllowed(double epsilon=1.e-9) const
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [5/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat,
const std::vector< TVector3 > &  mumAndDgtr_p3 
)

Definition at line 212 of file DalitzEvent.cpp.

215  : _pat(pat)
216  , _rememberPhaseSpace(-9999.)
217  , _aValue(-9999.)
218  , _weight(1)
220  , _vectorOfValues()
221  , _vectorOfWeights()
222  , _s(pat.size())
223  , _t(pat.size())
224  , _permutationIndex(0)
225  , _perm(pat)
226 {
227  setMomenta(mumAndDgtrs_p3);
228  if(_pat.size() < 3 || _pat.size() != mumAndDgtrs_p3.size()) shoutAndKill();
229  resetST();
230  _eventCounter++;
231 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool setMomenta(const std::vector< TLorentzVector > &mumAndDgtr_p4)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [6/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat)

Definition at line 232 of file DalitzEvent.cpp.

233  : _pat(pat)
234  , _rememberPhaseSpace(-9999.)
235  , _aValue(-9999.)
236  , _weight(1)
238  , _vectorOfValues()
239  , _vectorOfWeights()
240  , _s(pat.size())
241  , _t(pat.size())
242  , _permutationIndex(0)
243  , _perm(pat)
244 {
245  if (_pat.size() < 3) shoutAndKill();
246  _p.resize(_pat.size());
247  resetST();
248  _eventCounter++;
249 } // generate or add momemta later.
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [7/12]

DalitzEvent::DalitzEvent ( const DalitzEventPattern pat,
TRandom *  rnd 
)

Definition at line 251 of file DalitzEvent.cpp.

252  : _pat(pat)
253  , _rememberPhaseSpace(-9999.)
254  , _aValue(-9999.)
255  , _weight(1)
257  , _vectorOfValues()
258  , _vectorOfWeights()
259  , _s(pat.size())
260  , _t(pat.size())
261  , _permutationIndex(0)
262  , _perm(pat)
263 {
264  if (_pat.size() < 3) shoutAndKill();
265  resetST();
267  _eventCounter++;
268 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
bool shoutAndKill()
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
void generateThisToPhaseSpace(TRandom *rnd=0)
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [8/12]

DalitzEvent::DalitzEvent ( TNtupleD *  ntp)

Definition at line 365 of file DalitzEvent.cpp.

366  : _rememberPhaseSpace(-9999.)
367  , _aValue(-9999.)
368  , _weight(1)
370  , _vectorOfValues()
371  , _vectorOfWeights()
372  , _permutationIndex(0)
373 {
374  if(! fromTree(ntp)){
375  cout << "ERROR in DalitzEvent constructor from ntuple"
376  << " something went wrong!"
377  << endl;
378  }
379  _eventCounter++;
380 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
bool fromTree(TTree *tree)
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [9/12]

DalitzEvent::DalitzEvent ( const IDalitzEvent ievt)

Definition at line 293 of file DalitzEvent.cpp.

294  : _pat(other.eventPattern())
295  , _rememberPhaseSpace(other.phaseSpace())
298  , _aValue(other.getAValue())
299  , _weight(other.getWeight())
300  , _generatorPdfRelativeToPhaseSpace(other.getGeneratorPdfRelativeToPhaseSpace())
301  , _vectorOfValues(other.getVectorOfValues())
302  , _vectorOfWeights(other.getVectorOfWeights())
303  , _permutationIndex(0)
304  , _perm(other.eventPattern())
305 {
306  int np = other.eventPattern().size();
307  _p.resize(np);
308  for(int i=0; i < np; i++) _p[i] = other.p(i);
309  resetST();
310  _eventCounter++;
311 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [10/12]

DalitzEvent::DalitzEvent ( const IDalitzEvent ievt)

Definition at line 270 of file DalitzEvent.cpp.

271  : _pat(other->eventPattern())
272  // , _p(other.)
273  , _rememberPhaseSpace(other->phaseSpace())
276  , _aValue(other->getAValue())
277  , _weight(other->getWeight())
278  , _generatorPdfRelativeToPhaseSpace(other->getGeneratorPdfRelativeToPhaseSpace())
279  , _vectorOfValues(other->getVectorOfValues())
280  , _vectorOfWeights(other->getVectorOfWeights())
281  // , _s()
282  //, _t(other._t)
283  , _permutationIndex(0)
284  , _perm(other->eventPattern())
285 {
286  int np = other->eventPattern().size();
287  _p.resize(np);
288  for(int i=0; i < np; i++) _p[i] = other->p(i);
289  resetST();
290  _eventCounter++;
291 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [11/12]

DalitzEvent::DalitzEvent ( const DalitzEvent other)

Definition at line 315 of file DalitzEvent.cpp.

316  : IDalitzEvent() // the compiler wants it.
317  , _pat(other->_pat)
318  , _p(other->_p)
322  , _aValue(other->_aValue)
323  , _weight(other->_weight)
327  , _s(other->_s)
328  , _t(other->_t)
329  , _sijMap(other->_sijMap)
331  , _perm(other->_perm)
332 {
333  bool dbThis=false;
334  if(dbThis) cout << " copy constructor of DalitzEvent called."
335  << " weight before " << other->getWeight()
336  << " and after " << this->getWeight()
337  << endl;
338  _eventCounter++;
339 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48
Permutator _perm
Definition: DalitzEvent.h:63
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
virtual double getWeight() const
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
std::map< std::vector< int >, double > _sijMap
Definition: DalitzEvent.h:59
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ DalitzEvent() [12/12]

DalitzEvent::DalitzEvent ( const DalitzEvent other)

Definition at line 340 of file DalitzEvent.cpp.

341  : IDalitzEvent() // the compiler wants it.
342  , _pat(other._pat)
343  , _p(other._p)
347  , _aValue(other._aValue)
348  , _weight(other._weight)
352  , _s(other._s)
353  , _t(other._t)
354  , _sijMap(other._sijMap)
356  , _perm(other._perm)
357 {
358  bool dbThis=false;
359  if(dbThis) cout << " copy constructor of DalitzEvent called."
360  << " weight before " << other.getWeight()
361  << " and after " << this->getWeight()
362  << endl;
363  _eventCounter++;
364 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48
Permutator _perm
Definition: DalitzEvent.h:63
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
virtual double getWeight() const
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
int _permutationIndex
Definition: DalitzEvent.h:62
static long int _eventCounter
Definition: DalitzEvent.h:37
std::map< std::vector< int >, double > _sijMap
Definition: DalitzEvent.h:59
double _aValue
Definition: DalitzEvent.h:50
double _weight
Definition: DalitzEvent.h:51

◆ ~DalitzEvent()

DalitzEvent::~DalitzEvent ( )
virtual

Definition at line 382 of file DalitzEvent.cpp.

382  {
383  _eventCounter--;
384  // cout << "DalitzEvent destroyed\n" << endl;
385 }
static long int _eventCounter
Definition: DalitzEvent.h:37

Member Function Documentation

◆ assignUniqueRememberNumber()

long int DalitzEvent::assignUniqueRememberNumber ( )
static

Definition at line 29 of file DalitzEvent.cpp.

29  {
30  return _rememberVectorCounter++;
31 }
static long int _rememberVectorCounter
Definition: DalitzEvent.h:39

◆ BDet()

double DalitzEvent::BDet ( ) const
protected

Definition at line 883 of file DalitzEvent.cpp.

883  {
884  /*
885  cout << " for sij " << t(0,1) << ", " << s(1,2) << ", " << s(2,3)
886  << ", " << s(3,4) << ", " << t(4,0) << endl;
887  cout << "\n " << m2(0) << ", " << m2(1) << ", " << m2(2)
888  << ", " << m2(3) << m2(4) << endl;
889  */
890 
891  if(nDgtr() < 4){
892  cout << "ERROR in DalitzEvent::BDet(): nDgtr() = "
893  << nDgtr() << " < 4. How did this happen?"
894  << endl;
895  return +9999;
896  }
897 
898  TMatrixDSym B(6);
899  B(0,0)=0;
900  B(1,0)=B(0,1)=m2(2); B(1,1) = 0;
901  B(2,0)=B(0,2)=s(2,3) ; B(2,1)=B(1,2) = m2(3); B(2,2)=0;
902  B(3,0)=B(0,3)=t(0,1) ; B(3,1)=B(1,3) = s(3,4) ; B(3,2)=B(2,3)=m2(4); B(3,3)=0;
903  B(4,0)=B(0,4)=m2(1); B(4,1)=B(1,4) = s(1,2) ; B(4,2)=B(2,4)=t(4,0) ; B(4,3)=B(3,4)=m2(0); B(4,4)=0;
904  B(5,0)=B(0,5)=B(5,1)=B(1,5)=B(5,2)=B(2,5)=B(5,3)=B(3,5)=B(5,4)=B(4,5)=1;
905  B(5,5)=0;
906 
907  /*
908  double Bdet = B.Determinant();
909  cout << "Bdet = " << Bdet << endl;
910  */
911 
912  return B.Determinant();
913 }
double m2(unsigned int i) const
int nDgtr() const
Definition: DalitzEvent.h:222
virtual double t(unsigned int i, unsigned int j) const
virtual double s(unsigned int i, unsigned int j) const

◆ C_conjugateYourself()

void DalitzEvent::C_conjugateYourself ( )
virtual

Definition at line 808 of file DalitzEvent.cpp.

808  {
810 }
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
DalitzEventPattern makeCPConjugate() const

◆ clone()

IDalitzEvent * DalitzEvent::clone ( ) const
virtual

Implements IDalitzEvent.

Definition at line 387 of file DalitzEvent.cpp.

387  {
388  return new DalitzEvent(*this);
389 }

◆ CP_conjugateYourself()

void DalitzEvent::CP_conjugateYourself ( )
virtual

Definition at line 811 of file DalitzEvent.cpp.

811  {
814 }
virtual void P_conjugateYourself()
virtual void C_conjugateYourself()

◆ currentPermutation()

const Permutation& DalitzEvent::currentPermutation ( ) const
inline

Definition at line 148 of file DalitzEvent.h.

148  {
149  return _perm[_permutationIndex];
150  }
Permutator _perm
Definition: DalitzEvent.h:63
int _permutationIndex
Definition: DalitzEvent.h:62

◆ Delta4()

double DalitzEvent::Delta4 ( ) const
protected

Definition at line 914 of file DalitzEvent.cpp.

914  {
915  return -BDet()/16;
916 }
double BDet() const

◆ evalsij_intern()

double DalitzEvent::evalsij_intern ( const std::vector< int > &  intern_indices) const
protected

Definition at line 699 of file DalitzEvent.cpp.

700  {
701  TLorentzVector psum(0.0, 0.0, 0.0, 0.0);
702  for(unsigned int i=0; i< intern_indices.size(); i++){
703  //cout << " p(" << i << ") " << p(indices[i]);
704  psum += p_intern(intern_indices[i]);
705  }
706  return psum.Mag2();
707 }
const TLorentzVector & p_intern(unsigned int i) const

◆ eventCounter()

static long int DalitzEvent::eventCounter ( )
inlinestatic

Definition at line 125 of file DalitzEvent.h.

125 { return _eventCounter;}
static long int _eventCounter
Definition: DalitzEvent.h:37

◆ eventPattern()

virtual const DalitzEventPattern& DalitzEvent::eventPattern ( ) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 173 of file DalitzEvent.h.

173 {return _pat;}
DalitzEventPattern _pat
Definition: DalitzEvent.h:42

◆ fillNtupleVarArray()

bool DalitzEvent::fillNtupleVarArray ( std::vector< Double_t > &  array) const

Definition at line 1014 of file DalitzEvent.cpp.

1014  {
1015 
1016  if(eventPattern().empty()){
1017  return false;
1018  }
1019  int counter=0;
1020  int maxCounter = (array.size() - 2) - (singleParticleNtpArraySize()-1);
1021  for(unsigned int i=0; i< _p.size() && counter < maxCounter; i++){
1022  array.at(counter++) = _pat[i];
1023  array.at(counter++) = _p[i].E();
1024  array.at(counter++) = _p[i].X();
1025  array.at(counter++) = _p[i].Y();
1026  array.at(counter++) = _p[i].Z();
1027  }
1028  array.at(counter++) = getWeight();
1029  array.at(counter++) = getGeneratorPdfRelativeToPhaseSpace();
1030  return true;
1031 }
virtual const DalitzEventPattern & eventPattern() const
Definition: DalitzEvent.h:173
virtual double getWeight() const
T & at(unsigned int i)
virtual double getGeneratorPdfRelativeToPhaseSpace() const
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
static int singleParticleNtpArraySize()

◆ fromNtuple()

bool DalitzEvent::fromNtuple ( TNtupleD *  ntp)

Definition at line 1061 of file DalitzEvent.cpp.

1061  {
1062  // assumes ntuple is set to the correct
1063  // entry (using ntp->GetEntry(int i))
1064  bool dbThis=false;
1065  if(dbThis) cout << "DalitzEvent::fromNtuple(" << ntp->GetName() << ") called" << endl;
1066  int arraySize = ntp->GetNvar();
1067  // int arraySize = ntp->GetNbranches();
1068  if(dbThis) cout << " arraySize " << arraySize << endl;
1069  bool newFormat = false;
1070  if(0 == (arraySize-2)%singleParticleNtpArraySize() ){
1071  if(dbThis) cout << "New Format" << endl;
1072  newFormat = true; // includes weight & generatorPdf
1073  }else if( 0 == (arraySize)%singleParticleNtpArraySize()){
1074  newFormat = false; // excludes weight & generatorPdf
1075  }else{
1076  cout << "ERROR in DalitzEvent::fromNtuple length of array"
1077  << " in ntuple = " << arraySize
1078  << "\n > not a multiple of # entries per particle = "
1080  << endl;
1081  return false;
1082  }
1083  const Double_t* array = ntp->GetArgs();
1084 
1085  if(dbThis) cout << " array ptr " << array << endl;
1086  int counter = 0;
1087  int maxCounter = arraySize - (singleParticleNtpArraySize()-1);
1088  if(newFormat) maxCounter -= 2;
1089  int nParticles;
1090  if(newFormat){
1091  nParticles = (arraySize-2)/singleParticleNtpArraySize();
1092  }else{
1093  nParticles = arraySize/singleParticleNtpArraySize();
1094  }
1095  if(dbThis){
1096  cout << "nParticles = " << nParticles
1097  << " arraySize = " << arraySize
1098  << endl;
1099  }
1100 
1101  _pat.resize(nParticles);
1102  _p.resize(nParticles);
1103  for(int i=0; i< nParticles && counter < maxCounter; i++){
1104  if(dbThis) cout << "filling " << i << " th particle: " << endl;
1105  _pat[i] = nearestInt(array[counter++]);
1106  _p[i].SetE(array[counter++]);
1107  _p[i].SetX(array[counter++]);
1108  _p[i].SetY(array[counter++]);
1109  _p[i].SetZ(array[counter++]);
1110 
1111  // check units (should be MeV, but will recognise if its GeV and
1112  // translate it to MeV)
1113 
1114  double units(1);
1115  double dMeV = abs(_pat[i].mass() - _p[i].M()*MeV);
1116  double dGeV = abs(_pat[i].mass() - _p[i].M()*GeV);
1117  if(dMeV < dGeV) units=MeV;
1118  else units=GeV;
1119  _p[i] *= units;
1120 
1121  if(dbThis){
1122  cout << "\t" << _pat[i] << " " << _p[i]
1123  << ", m= " << _p[i].M() << endl;
1124  }
1125  }
1126  if(newFormat){
1127  setWeight(array[counter++]);
1128  setGeneratorPdfRelativeToPhaseSpace(array[counter++]);
1129  }
1130  if(dbThis){
1131  cout << " weight, pdf "
1132  << getWeight() << ", "
1134  }
1135  resetST();
1137  if(dbThis){
1138  cout << "this is me now " << endl;
1139  this->print();
1140  cout << "done it all, returning" << endl;
1141  }
1142  return true;
1143 }
virtual void setWeight(double w)
void setPattern(const DalitzEventPattern &pat)
Definition: Permutator.cpp:23
void resize(unsigned int N)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
static const double MeV
int nearestInt(double f)
Definition: Utils.cpp:26
virtual double getWeight() const
virtual double getGeneratorPdfRelativeToPhaseSpace() const
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
static const double GeV
void print(std::ostream &os=std::cout) const
virtual void setGeneratorPdfRelativeToPhaseSpace(double gpdf)
static int singleParticleNtpArraySize()

◆ fromParasTree()

bool DalitzEvent::fromParasTree ( TTree *  ntp)

Definition at line 1271 of file DalitzEvent.cpp.

1271  {
1272  bool dbThis=false;
1273  bool forcePDGMassForFinalState=true;
1274 
1275  if(dbThis) cout << "DalitzEvent::fromTREE(" << ntp << ") called" << endl;
1276 
1277  int arraySize = ntp->GetNbranches();
1278  if(dbThis) cout << " arraySize " << arraySize << endl;
1279  const TObjArray* branchArray = ntp->GetListOfBranches();
1280 
1281  map<string, map<std::string, double> > prtMap;
1282 
1283  double weight=1;
1284  double genPdf=1;
1285  for(int i=0; i< arraySize; i++){
1286  TObjArray* leafArray=
1287  ((TBranch*) (*branchArray)[i])->GetListOfLeaves();
1288  double leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1289  std::string leafName = ((TLeaf*) (*leafArray)[0])->GetName();
1290 
1291  std::string parName, element;
1292  parseNtpEntryName(leafName, parName, element);
1293  if(dbThis){
1294  cout << " leaf[ " << i << "] "
1295  << leafName << "; " << leafVal << endl;
1296  }
1297 
1298  if("weight" == leafName) weight=leafVal;
1299  else if("genPdf" == leafName) genPdf = leafVal;
1300  else{
1301  std::string parName, element;
1302  if(!parseNtpEntryName(leafName, parName, element)){
1303  // leafName doesn't contain an underscore.
1304  continue ;
1305  }
1306  (prtMap[parName])[element] = leafVal;
1307  if(dbThis){
1308  cout << "\t leaf[ " << i << "] name parsed: "
1309  << parName << ", " << element
1310  << endl;
1311  }
1312  }
1313  }
1314  int nParticles = prtMap.size();
1315  if(dbThis) cout << "nParticles = " << nParticles << endl;
1316  _pat.resize(nParticles);
1317  _p.resize(nParticles);
1318 
1319  int counter=0;
1320  for(map<string, map<string, double> >::iterator it =
1321  prtMap.begin();
1322  it != prtMap.end();
1323  it++, counter++){
1324  _pat[counter] = (it->second)["pdg"];
1325  if(dbThis){
1326  cout << "just filled _pat[" << counter << "] with: "
1327  << (it->second)["pdg"] << "; see: " << _pat[counter]
1328  << endl;
1329  }
1330  _p[counter].SetE((it->second)["E"]);
1331  _p[counter].SetPx((it->second)["Px"]);
1332  _p[counter].SetPy((it->second)["Py"]);
1333  _p[counter].SetPz((it->second)["Pz"]);
1334 
1335  // check units (Paras usually GeV, else usually MeV)
1336  double units;
1337  double dMeV = abs(_pat[counter].mass() - _p[counter].M()*MeV);
1338  double dGeV = abs(_pat[counter].mass() - _p[counter].M()*GeV);
1339  if(dMeV < dGeV) units=MeV;
1340  else units=GeV;
1341  _p[counter] *= units;
1342 
1343  if(forcePDGMassForFinalState && counter > 0){
1344  double m = _pat[counter].mass();
1345  double p3 = _p[counter].Rho();
1346  _p[counter].SetE(sqrt(m*m + p3*p3));
1347  }
1348  }
1349  setWeight(weight);
1351  if(dbThis){
1352  cout << " weight, pdf "
1353  << getWeight() << ", "
1355  for(int i=0; i < nParticles; i++){
1356  cout << _pat[i] << " ) " << _p[i] << " m = " << _p[i].M() << endl;
1357  }
1358  }
1359  resetST();
1361  if(dbThis) cout << "done it all, returning" << endl;
1362  return true;
1363  }
virtual void setWeight(double w)
void setPattern(const DalitzEventPattern &pat)
Definition: Permutator.cpp:23
void resize(unsigned int N)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
static const double MeV
static bool parseNtpEntryName(const std::string &entry, std::string &part1, std::string &part2)
virtual double getWeight() const
virtual double getGeneratorPdfRelativeToPhaseSpace() const
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
static const double GeV
double m(unsigned int i) const
virtual void setGeneratorPdfRelativeToPhaseSpace(double gpdf)

◆ fromParasTreeOld()

bool DalitzEvent::fromParasTreeOld ( TTree *  ntp)

Definition at line 1144 of file DalitzEvent.cpp.

1144  {
1145  bool dbThis=false;
1146  bool forcePDGMassForFinalState=true;
1147 
1148  if(dbThis) cout << "DalitzEvent::fromTREE(" << ntp << ") OLD called" << endl;
1149 
1150  int arraySize = ntp->GetNbranches();
1151  if(dbThis) cout << " arraySize " << arraySize << endl;
1152  bool newFormat = false;
1153  if(0 == (arraySize-2)%singleParticleNtpArraySize() ){
1154  newFormat = true; // includes weight & generatorPdf
1155  if(dbThis) cout << " new Format" << endl;
1156  }else if( 0 == (arraySize)%singleParticleNtpArraySize()){
1157  newFormat = false; // excludes weight & generatorPdf
1158  if(dbThis) cout << " old Format" << endl;
1159  }else{
1160  cout << "ERROR in DalitzEvent::fromTree length of array"
1161  << " in ntuple = " << arraySize
1162  << "\n > not a multiple of # entries per particle = "
1164  << endl;
1165  return false;
1166  }
1167  const TObjArray* branchArray = ntp->GetListOfBranches();
1168 
1169  if(dbThis) cout << " array ptr " << branchArray << endl;
1170  int counter = 0;
1171  int maxCounter = arraySize - (singleParticleNtpArraySize()-1);
1172  if(newFormat) maxCounter -= 2;
1173  int nParticles;
1174  if(newFormat){
1175  nParticles = (arraySize-2)/singleParticleNtpArraySize();
1176  }else{
1177  nParticles = arraySize/singleParticleNtpArraySize();
1178  }
1179  if(dbThis){
1180  cout << "nParticles = " << nParticles
1181  << " arraySize = " << arraySize
1182  << endl;
1183  }
1184 
1185  _pat.resize(nParticles);
1186  _p.resize(nParticles);
1187  for(int i=0; i< nParticles && counter < maxCounter; i++){
1188  if(dbThis) cout << "filling " << i << " th particle: " << endl;
1189  TObjArray* leafArray;
1190  double leafVal = -9999;
1191  std::string leafName;
1192  leafArray=
1193  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1194  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1195  leafName = ((TLeaf*) (*leafArray)[0])->GetName();
1196  if(dbThis) cout << " leaf[ " << counter-1 << "] "
1197  << leafVal << endl;
1198  _p[i].SetE(leafVal*GeV);
1199  leafArray =
1200  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1201  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1202  leafName = ((TLeaf*) (*leafArray)[0])->GetName();
1203  if(dbThis){
1204  cout << " leaf[ " << counter-1 << "] " << leafName
1205  << ": " << leafVal << endl;
1206  }
1207  _p[i].SetX(leafVal*GeV);
1208  leafArray =
1209  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1210  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1211  if(dbThis) cout << " leaf[ " << counter-1 << "] " << leafVal << endl;
1212  _p[i].SetY(leafVal*GeV);
1213  leafArray =
1214  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1215  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1216  if(dbThis) cout << " leaf[ " << counter-1 << "] " << leafVal << endl;
1217  _p[i].SetZ(leafVal*GeV);
1218  leafArray =
1219  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1220  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1221  if(dbThis) cout << " leaf[ " << counter-1 << "] " << leafVal << endl;
1222  _pat[i] = nearestInt(leafVal);
1223 
1224  if(forcePDGMassForFinalState && i > 0){
1225  double m=_pat[i].mass();
1226  _p[i].SetE(sqrt(m*m + _p[i].Rho()*_p[i].Rho()));
1227  }
1228  if(dbThis) cout << " got pattern " << (int) _pat[i] << endl;
1229 
1230  if(dbThis){
1231  cout << "\t" << _pat[i] << " " << _p[i]
1232  << ", m= " << _p[i].M() << endl;
1233  }
1234  }
1235  if(newFormat){
1236  TObjArray* leafArray;
1237  double leafVal = -9999;
1238  leafArray =
1239  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1240  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1241  if(dbThis) cout << " leaf[ " << counter-1 << "] " << leafVal << endl;
1242  setWeight(leafVal);
1243 
1244  leafArray =
1245  ((TBranch*) (*branchArray)[counter++])->GetListOfLeaves();
1246  leafVal = ((TLeaf*) (*leafArray)[0])->GetValue();
1248  }
1249  if(dbThis){
1250  cout << " weight, pdf "
1251  << getWeight() << ", "
1253  }
1254  resetST();
1256  if(dbThis) cout << "done it all, returning" << endl;
1257  return true;
1258 }
virtual void setWeight(double w)
void setPattern(const DalitzEventPattern &pat)
Definition: Permutator.cpp:23
void resize(unsigned int N)
bool resetST()
Permutator _perm
Definition: DalitzEvent.h:63
int nearestInt(double f)
Definition: Utils.cpp:26
virtual double getWeight() const
virtual double getGeneratorPdfRelativeToPhaseSpace() const
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
static const double GeV
double m(unsigned int i) const
virtual void setGeneratorPdfRelativeToPhaseSpace(double gpdf)
static int singleParticleNtpArraySize()

◆ fromTree()

bool DalitzEvent::fromTree ( TTree *  tree)

Definition at line 1047 of file DalitzEvent.cpp.

1047  {
1048  bool dbThis=false;
1049  // figure out format: TNtuple or what?
1050  if(0 == tree) return 0;
1051  std::string treeName(tree->ClassName());
1052  if(dbThis) cout << "ClassName \"" << treeName << "\"" << endl;
1053  if(treeName == (std::string) "TNtupleD"){
1054  if(dbThis)cout << "DalitzEvent::fromTree: tree is TNtupleD" << endl;
1055  return fromNtuple((TNtupleD*)tree);
1056  }else{
1057  if(dbThis) cout << "DalitzEvent::fromTree: tree is ParasTree" << endl;
1058  return fromParasTree(tree);
1059  }
1060 }
bool fromParasTree(TTree *ntp)
bool fromNtuple(TNtupleD *ntp)

◆ generateThisToPhaseSpace() [1/3]

void DalitzEvent::generateThisToPhaseSpace ( TRandom *  rnd = 0)

Definition at line 546 of file DalitzEvent.cpp.

546  {
547  generateThisToPhaseSpace(0.0, rnd);
548 }
void generateThisToPhaseSpace(TRandom *rnd=0)

◆ generateThisToPhaseSpace() [2/3]

void DalitzEvent::generateThisToPhaseSpace ( double  p_of_D,
TRandom *  rnd = 0 
)

Definition at line 549 of file DalitzEvent.cpp.

549  {
550  TVector3 p3(0,0,p);
551  generateThisToPhaseSpace(p3, rnd);
552 }
void generateThisToPhaseSpace(TRandom *rnd=0)
virtual const TLorentzVector & p(unsigned int i) const

◆ generateThisToPhaseSpace() [3/3]

void DalitzEvent::generateThisToPhaseSpace ( const TVector3 &  pVec_of_D,
TRandom *  rnd = 0 
)

Definition at line 553 of file DalitzEvent.cpp.

554  {
555  bool dbthis=false;
556  if(dbthis) cout << "got called" << endl;
557  TRandom* remember_gRandom = gRandom;
558  bool needToReset_gRandom=false;
559  if(0 != rnd){
560  needToReset_gRandom=true;
561  gRandom = rnd; // horrible fix because
562  // TGenPhaseSpace can only use gRandom
563  }
564  if(dbthis) cout << "got so far" << endl;
565  TGenPhaseSpace gps; // uses gRandom, so no worries
566  // about re-starting with the same seed (won't happen),
567  // but worries about what kind of
568  // random number generator gRandom is.
569 
570  TLorentzVector mumsP4;
571  mumsP4.SetXYZM(mumsp3.X()
572  , mumsp3.Y()
573  , mumsp3.Z()
574  , _pat[0].mass());
575 
576  Double_t* mdgt = new Double_t[nDgtr()];
577  if(dbthis) cout << "booked " << nDgtr() << " daughters " << endl;
578  for(int i=1; i<= nDgtr(); i++){
579  mdgt[i-1] = this->m(i);
580  }
581 
582  gps.SetDecay(mumsP4, nDgtr() , mdgt);
583  double maxWeight = 1.0;//gps.GetWtMax();
584  double weight=0;
585  if(dbthis) cout << " maxWeight " << maxWeight << endl;
586  //normal:
587  do{
588  weight=gps.Generate();
589  }while(weight < gRandom->Rndm() * maxWeight);
590  // for debug:
591  // weight=gps.Generate();
592 
593  if(needToReset_gRandom){
594  gRandom = remember_gRandom;
595  }
596  if(weight > maxWeight){
597  cout << "ERROR in DalitzEvent::generateThisToPhaseSpace"
598  << " weight = " << weight << " > " << maxWeight << " = maxWeight"
599  << "\n > this is a programming error that needs to be fixed."
600  << endl;
601  throw "...so I got maxWeight wrong.";
602  }
603  if(dbthis) cout << " this weight " << weight << endl;
604  if(dbthis) cout << " max weight " << gps.GetWtMax() << endl;
605 
606  _p.resize(nDgtr()+1);
607  if(dbthis) cout << "generated event" << endl;
608  setP(0, mumsP4);
609  if(dbthis)cout << "set p0" << endl;
610  for(unsigned int i=1; i< _pat.size(); i++){
611  if(dbthis)cout << "setting P (" << i << ") "
612  << (*(gps.GetDecay(i-1)))
613  << "..." << endl;
614  setP(i, (*(gps.GetDecay(i-1))));
615  if(dbthis)cout << " and look what I've set: " << p(i) << endl;
616  }
617  delete[] mdgt;
618  _weight = 1.0;
620  return;
621 }
int nDgtr() const
Definition: DalitzEvent.h:222
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
virtual const TLorentzVector & p(unsigned int i) const
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
double m(unsigned int i) const
double _weight
Definition: DalitzEvent.h:51
void setP(unsigned int i, const TLorentzVector &p4)

◆ getAValue()

virtual double DalitzEvent::getAValue ( ) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 171 of file DalitzEvent.h.

171 {return _aValue;}
double _aValue
Definition: DalitzEvent.h:50

◆ getGeneratorPdfRelativeToPhaseSpace()

double DalitzEvent::getGeneratorPdfRelativeToPhaseSpace ( ) const
virtual

Implements IDalitzEvent.

Definition at line 396 of file DalitzEvent.cpp.

396  {
398 }
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52

◆ getValueFromVector()

double DalitzEvent::getValueFromVector ( unsigned int  i) const
virtual

Implements IDalitzEvent.

Definition at line 425 of file DalitzEvent.cpp.

425  {
426  return _vectorOfValues[i];
427 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54

◆ getVectorOfValues() [1/2]

const std::vector< double > & DalitzEvent::getVectorOfValues ( ) const
virtual

Implements IDalitzEvent.

Definition at line 403 of file DalitzEvent.cpp.

403  {
404  return _vectorOfValues;}
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54

◆ getVectorOfValues() [2/2]

std::vector< double > & DalitzEvent::getVectorOfValues ( )
virtual

Implements IDalitzEvent.

Definition at line 406 of file DalitzEvent.cpp.

406  {
407  return _vectorOfValues;}
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54

◆ getVectorOfWeights() [1/2]

const std::vector< double > & DalitzEvent::getVectorOfWeights ( ) const
virtual

Implements IDalitzEvent.

Definition at line 409 of file DalitzEvent.cpp.

409  {
410  return _vectorOfWeights;}
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55

◆ getVectorOfWeights() [2/2]

std::vector< double > & DalitzEvent::getVectorOfWeights ( )
virtual

Implements IDalitzEvent.

Definition at line 412 of file DalitzEvent.cpp.

412  {
413  return _vectorOfWeights;}
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55

◆ getWeight()

double DalitzEvent::getWeight ( ) const
virtual

Implements IDalitzEvent.

Definition at line 390 of file DalitzEvent.cpp.

390  {
391  return _weight;
392 }
double _weight
Definition: DalitzEvent.h:51

◆ getWeightFromVector()

double DalitzEvent::getWeightFromVector ( unsigned int  i) const
virtual

Implements IDalitzEvent.

Definition at line 429 of file DalitzEvent.cpp.

429  {
430  return _vectorOfWeights[i];
431 }
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55

◆ kinematicallyAllowed()

bool DalitzEvent::kinematicallyAllowed ( double  epsilon = 1.e-9) const

Definition at line 466 of file DalitzEvent.cpp.

466  {
467  bool dbThis=false;
468 
469  TLorentzVector pdiff = (p(0) - p_dgtr_sum());
470  TLorentzVector psum = (p(0) + p_dgtr_sum());
471 
472  double diffsum = 0;
473  double sumsum = 0;
474  for(unsigned int i=0; i<4; i++){
475  diffsum += pdiff[i]*pdiff[i];
476  sumsum += psum[i] *psum[i];
477  }
478  double avsum = sumsum/2.0;
479 
480  if(diffsum > epsilon*epsilon * avsum){
481  if(dbThis){
482  cout << "DalitzEvent::kinematicallyAllowed? "
483  << " failing diffsum: " << diffsum
484  << " max allowed: " << epsilon*epsilon*avsum << endl;
485  }
486  return false;
487  }
488 
489  for(unsigned int i = 0; i<_pat.size(); i++){
490  double mreco = p(i).M();
491  if(mreco < 0){
492  if(dbThis) cout << "failing mreco " << mreco << " < 0" << endl;
493  return false;
494  }
495  double avmass = (mreco + m(i))/2.0;
496  if( fabs(mreco - m(i)) > epsilon*avmass ){
497  if(dbThis){
498  cout << " failing mreco - m(i) fabs("
499  << mreco << " - " << m(i) << ") = "
500  << fabs(mreco - m(i)) << " > " << epsilon*avmass
501  << endl;
502  }
503  return false;
504  }
505  }
506  return true;
507 }
TLorentzVector p_dgtr_sum() const
virtual const TLorentzVector & p(unsigned int i) const
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
double m(unsigned int i) const

◆ m()

double DalitzEvent::m ( unsigned int  i) const

Definition at line 433 of file DalitzEvent.cpp.

433  {
434  // easy access to the nominal masses
435  if(i >= _pat.size()){
436  cout << "ERROR in DalitzEvent::m(unsigned int i)"
437  << " index i= " << i << " out of range!"
438  << " allowed range: [0," << _pat.size()-1 << "]"
439  << " returning -9999"
440  << endl;
441  return -9999;
442  }
443  return _pat[i].mass();
444 }
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const

◆ m2()

double DalitzEvent::m2 ( unsigned int  i) const

Definition at line 445 of file DalitzEvent.cpp.

445  {
446  if(i >= _pat.size()){
447  cout << "ERROR in DalitzEvent::m2(unsigned int i)"
448  << " index i= " << i << " out of range!"
449  << " allowed range: [0," << _pat.size()-1 << "]"
450  << " returning -9999"
451  << endl;
452  return -9999;
453  }
454  return m(i)*m(i);
455 }
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
double m(unsigned int i) const

◆ makeCalculate4BodyProps()

Calculate4BodyProps DalitzEvent::makeCalculate4BodyProps ( ) const

Definition at line 751 of file DalitzEvent.cpp.

751  {
752  if(nDgtr() != 4){
753  cout << "ERROR in DalitzEvent::makeCalculate4BodyProps()"
754  << "\n > You can't expect a sensible Calculate4BodyProps"
755  << " Object for a " << nDgtr() << " body decay!"
756  << endl;
757  return Calculate4BodyProps(0,0,0,0,0, 0,0,0,0,0);
758  }
759  return Calculate4BodyProps(t(0,1), s(1,2), s(2,3), s(3,4), t(4,0)
760  , m(0), m(1), m(2), m(3), m(4)
761  );
762 }
int nDgtr() const
Definition: DalitzEvent.h:222
virtual double t(unsigned int i, unsigned int j) const
double m(unsigned int i) const
virtual double s(unsigned int i, unsigned int j) const

◆ makeNtupleVarnames()

std::string DalitzEvent::makeNtupleVarnames ( ) const

Definition at line 987 of file DalitzEvent.cpp.

987  {
988  std::string s="";
989  if(eventPattern().empty()){
990  return s;
991  }
992 
993  for(unsigned int i=0; i< _pat.size(); i++){
994  std::string name = _pat[i].name();
995  name = "_" + anythingToString((int) i) + "_" + name;
996  name = prtToNtpName(name);
997 
998  s += name + "_pdg:";
999  s += name + "_E:";
1000  s += name + "_Px:";
1001  s += name + "_Py:";
1002  s += name + "_Pz";
1003  // if(i != _pat.size() -1) s+= ":";
1004  s+= ":";
1005  }
1006  s += "weight:";
1007  s += "genPdf";
1008  return s;
1009 }
static std::string prtToNtpName(const std::string &s_in)
std::string name() const
virtual const DalitzEventPattern & eventPattern() const
Definition: DalitzEvent.h:173
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
std::string anythingToString(const T &anything)
Definition: Utils.h:62
virtual double s(unsigned int i, unsigned int j) const

◆ nDgtr()

int DalitzEvent::nDgtr ( ) const
inline

Definition at line 222 of file DalitzEvent.h.

222  {
223  return _pat.size()-1;
224  }
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const

◆ ntpToPrtName()

std::string DalitzEvent::ntpToPrtName ( const std::string &  s_in)
staticprotected

Definition at line 979 of file DalitzEvent.cpp.

979  {
980  std::string s(s_in);
981  for(int i=0; ntpNameChars[i] != '\0'; i++){
982  replace(s.begin(), s.end(), ntpNameChars[i], prtNameChars[i]);
983  }
984  return s;
985 }
static const char prtNameChars[]
Definition: DalitzEvent.h:31
virtual double s(unsigned int i, unsigned int j) const
static const char ntpNameChars[]
Definition: DalitzEvent.h:30

◆ ntupleVarArraySize()

unsigned int DalitzEvent::ntupleVarArraySize ( ) const

Definition at line 1032 of file DalitzEvent.cpp.

1032  {
1033  return _p.size() * singleParticleNtpArraySize() + 2;
1034 }
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
static int singleParticleNtpArraySize()

◆ numPermutations()

int DalitzEvent::numPermutations ( ) const
virtual

Implements IDalitzEvent.

Definition at line 843 of file DalitzEvent.cpp.

843  {
844  return _perm.size();
845 }
Permutator _perm
Definition: DalitzEvent.h:63
unsigned int size() const

◆ p() [1/2]

const TLorentzVector & DalitzEvent::p ( unsigned int  i) const
virtual

Implements IDalitzEvent.

Definition at line 673 of file DalitzEvent.cpp.

673  {
674  bool dbThis=false;
675  if(dbThis){
676  cout << "p(" << i << ") called. Current permutation: "
677  << currentPermutation()[i]
678  << endl;
679  }
680  return p_intern(currentPermutation()[i]);
681 }
const TLorentzVector & p_intern(unsigned int i) const
const Permutation & currentPermutation() const
Definition: DalitzEvent.h:148

◆ p() [2/2]

TLorentzVector & DalitzEvent::p ( unsigned int  i)
virtual

Definition at line 682 of file DalitzEvent.cpp.

682  {
683  bool dbThis=false;
684  if(dbThis){
685  cout << "p(" << i << ") called. Current permutation: "
686  << currentPermutation()[i]
687  << endl;
688  }
689  return p_intern(currentPermutation()[i]);
690 }
const TLorentzVector & p_intern(unsigned int i) const
const Permutation & currentPermutation() const
Definition: DalitzEvent.h:148

◆ P_conjugateYourself()

void DalitzEvent::P_conjugateYourself ( )
virtual

Definition at line 795 of file DalitzEvent.cpp.

795  {
796  // p-conjugates in mums restframe, but keeps old D momentum.
797 
798  resetST(); // shouldn't be necessary, but to be save.
799  if(_p.empty()) return;
800  TVector3 mums3Momentum(_p[0].X(), _p[0].Y(), _p[0].Z());
801  for(unsigned int i=0; i < _p.size(); i++){
802  _p[i].SetX( - _p[i].X() );
803  _p[i].SetY( - _p[i].Y() );
804  _p[i].SetZ( - _p[i].Z() );
805  }
806  setMothers3Momentum(mums3Momentum);
807 }
virtual void setMothers3Momentum(const TVector3 &mp3)
bool resetST()
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ p_dgtr_sum()

TLorentzVector DalitzEvent::p_dgtr_sum ( ) const

Definition at line 457 of file DalitzEvent.cpp.

457  {
458  TLorentzVector pd(0.0, 0.0, 0.0, 0.0);
459 
460  for(unsigned int i=1; i< _p.size(); i++){
461  pd += _p[i];
462  }
463  return pd;
464 }
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ p_intern() [1/2]

const TLorentzVector & DalitzEvent::p_intern ( unsigned int  i) const
protected

Definition at line 623 of file DalitzEvent.cpp.

623  {
624  if(i >= _p.size()){
625  cout << "in TLorentzVector& DalitzEvent::p(unsigned int i):"
626  << " index out of range. DalitzEventPattern:"
627  << "\n > \"" << eventPattern() << "\""
628  << "\n > num mother + daughters = " << _p.size()
629  << endl;
630  throw "index out or range in TLorentzVector& DalitzEvent::p(unsigned int i)";
631  }
632  return _p[i];
633 }
virtual const DalitzEventPattern & eventPattern() const
Definition: DalitzEvent.h:173
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ p_intern() [2/2]

TLorentzVector & DalitzEvent::p_intern ( unsigned int  i)
protected

Definition at line 635 of file DalitzEvent.cpp.

635  {
636  if(i >= _p.size()){
637  cout << "in TLorentzVector& DalitzEvent::p(unsigned int i):"
638  << " index out of range. DalitzEventPattern:"
639  << "\n > \"" << eventPattern() << "\""
640  << "\n > num mother + daughters = " << _p.size()
641  << endl;
642  throw "index out or range in TLorentzVector& DalitzEvent::p(unsigned int i)";
643  }
644  return _p[i];
645 }
virtual const DalitzEventPattern & eventPattern() const
Definition: DalitzEvent.h:173
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ parseNtpEntryName()

bool DalitzEvent::parseNtpEntryName ( const std::string &  entry,
std::string &  part1,
std::string &  part2 
)
staticprotected

Definition at line 1260 of file DalitzEvent.cpp.

1262  {
1263  int posOf_ = entry.find_last_of('_');
1264  if(std::string::npos == posOf_)
1265  return false ;
1266  part1 = entry.substr(0, posOf_);
1267  part2 = entry.substr(posOf_ +1);
1268  return true;
1269 }

◆ permutationIndex()

int DalitzEvent::permutationIndex ( ) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 154 of file DalitzEvent.h.

154  {
155  return _permutationIndex;
156  }
int _permutationIndex
Definition: DalitzEvent.h:62

◆ phaseSpace()

double DalitzEvent::phaseSpace ( ) const
virtual

Implements IDalitzEvent.

Definition at line 510 of file DalitzEvent.cpp.

510  {
511  if(nDgtr() <= 2){
513  return 1;
514  }
515  if(_p.size() != _pat.size()){
516  cout << "ERROR in DalitzEvent::phaseSpace: inconsistency!"
517  << " pattern array size != momentum array size"
518  << " " << _pat.size() << " != " << _p.size()
519  << ". returning -9999"
520  << endl;
521  return -9999;
522  }
523 
524  if(_rememberPhaseSpace > -9998){
525  // cout << " remembering... " << _rememberPhaseSpace << endl;
526  return _rememberPhaseSpace;
527  }
528 
529  if(nDgtr() == 3){
530  _rememberPhaseSpace = phaseSpace3(); // one or zero
531  return _rememberPhaseSpace;
532  }else if(nDgtr() == 4){
534  return _rememberPhaseSpace;
535  }else{
536  cout << "ERROR in DalitzEvent::phaseSpace inconsistency!"
537  << "\n > Sorry, don't know how to calculate phase-space"
538  << " factor for nDgtr = " << nDgtr()
539  << "\n > - can only do up to nDgtr <= 4"
540  << ". returning -9999"
541  << endl;
542  return -9999;
543  }
544 }
double phaseSpace4() const
int nDgtr() const
Definition: DalitzEvent.h:222
double phaseSpace3(double epsilon=1.e-9) const
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const

◆ phaseSpace3()

double DalitzEvent::phaseSpace3 ( double  epsilon = 1.e-9) const
protected

Definition at line 952 of file DalitzEvent.cpp.

952  {
953  bool dbThis=false;
954  // return either one or zero.
955  // one if it's kinematically allowed, zero otherwise.
956  if(kinematicallyAllowed(epsilon)){
957  if(dbThis) {
958  cout << "DalitzEvent::phaseSpace3(" << epsilon << ")"
959  << " returning 1" << endl;
960  }
961  //return pi*pi/(4.0*m(0)*m(0));
962  return 1;
963  }else{
964  if(dbThis) {
965  cout << "DalitzEvent::phaseSpace3(" << epsilon << ")"
966  << " returning 0" << endl;
967  }
968  return 0;
969  }
970 }
bool kinematicallyAllowed(double epsilon=1.e-9) const

◆ phaseSpace4()

double DalitzEvent::phaseSpace4 ( ) const
protected

Definition at line 917 of file DalitzEvent.cpp.

917  {
918  bool dbThis=false;
919  if(dbThis)cout << "\n\n called phaseSpace4() " << endl;
920  if(nDgtr() < 4){
921  cout << "ERROR in DalitzEvent::phaseSpace4: nDgtr() = "
922  << nDgtr() << " < 4. How did this happen?"
923  << endl;
924  return -9999;
925  }
926  if(dbThis){
927  Calculate4BodyProps c4bp(t_intern(0,1), s_intern(1,2), s_intern(2,3)
928  , s_intern(3,4), t_intern(4,0)
929  , m(0), m(1), m(2), m(3), m(4)
930  );
931 
932  cout << "4-body probs gets this phase space: "
933  << c4bp.phaseSpaceFactor() << endl;
934  cout << " ... for t01 <<... = "
935  << t_intern(0,1) << ", " << s_intern(1,2) << ", " << s_intern(2,3)
936  << ", " << s_intern(3,4) << ", " << t_intern(4,0)
937  << endl;
938  }
939 
940  // Note that this calculation assumes that this is a "physical" event.
941  // It doesn't check if it is outside the physically allowed
942  // phase space area, for which there are more requirements than
943  // only Delta4 < 0.
944  double invPhsp2 = -Delta4();
945  if(invPhsp2 <= 0) return 0;
946  double phsp = pi*pi/(32.0*m(0)*m(0)) * 1./sqrt(invPhsp2);
947  if(dbThis) cout << " ... I get: phsp = " << phsp << endl;
948  if(phsp < 0) return 0;
949  else return phsp;
950 }
static const double pi
double Delta4() const
int nDgtr() const
Definition: DalitzEvent.h:222
double t_intern(unsigned int i, unsigned int j) const
double m(unsigned int i) const
double s_intern(unsigned int i, unsigned int j) const

◆ print()

void DalitzEvent::print ( std::ostream &  os = std::cout) const
virtual

Implements IDalitzEvent.

Definition at line 764 of file DalitzEvent.cpp.

764  {
765  os << "Dalitz event for pattern "
766  << _pat;
767  for(unsigned int i=0; i<_p.size(); i++){
768  cout << "\n\t p(" << i << ") " << p(i);
769  }
770  os << "\n\t m(i) ";
771  for(unsigned int i=0; i<_p.size(); i++){
772  cout << m(i) << " ; ";
773  }
774  if(nDgtr() == 4){
775  os << "\n\t sij :"
776  << t(0,1) << ", " << s(1,2) << ", " << s(2,3)
777  << ", " << s(3,4) << ", " << t(4, 0);
778  }
779  os << "\n\t ps = " << this->phaseSpace();
780  os << "\t weight = " << this->getWeight();
781  os << endl;
782 }
virtual double phaseSpace() const
int nDgtr() const
Definition: DalitzEvent.h:222
virtual double getWeight() const
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
virtual const TLorentzVector & p(unsigned int i) const
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
virtual double t(unsigned int i, unsigned int j) const
double m(unsigned int i) const
virtual double s(unsigned int i, unsigned int j) const

◆ prtToNtpName()

std::string DalitzEvent::prtToNtpName ( const std::string &  s_in)
staticprotected

Definition at line 972 of file DalitzEvent.cpp.

972  {
973  std::string s(s_in);
974  for(int i=0; prtNameChars[i] != '\0'; i++){
975  replace(s.begin(), s.end(), prtNameChars[i], ntpNameChars[i]);
976  }
977  return s;
978 }
static const char prtNameChars[]
Definition: DalitzEvent.h:31
virtual double s(unsigned int i, unsigned int j) const
static const char ntpNameChars[]
Definition: DalitzEvent.h:30

◆ resetST()

bool DalitzEvent::resetST ( )

Definition at line 816 of file DalitzEvent.cpp.

816  {
817  // resets (possibly previously calculated) values for sij, tij
818 
819  _s.resize(_pat.size());
820  for(unsigned int i=0; i< _s.size(); i++){
821  _s[i].clear();
822  _s[i].resize(_pat.size());
823  for(unsigned int j=0; j< _s[i].size(); j++){
824  _s[i][j] = -1;
825  }
826  }
827 
828  _t.resize(_pat.size());
829  for(unsigned int i=0; i< _t.size(); i++){
830  _t[i].clear();
831  _t[i].resize(_pat.size());
832  for(unsigned int j=0; j< _t[i].size(); j++){
833  _t[i][j] = -1;
834  }
835  }
836  return true;
837 }
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const

◆ retrieveValue() [1/2]

virtual bool DalitzEvent::retrieveValue ( int  i,
std::complex< double > &  value,
long int  configNumber = 1 
)
inlinevirtual

Implements IDalitzEvent.

Definition at line 187 of file DalitzEvent.h.

188  {
189  /*std::cout << "I remember these numbers: " << std::endl;
190  for(unsigned int j=0; j < _rememberComplexFast.size(); j++){
191  std::complex<double> tempVal;
192  _rememberComplexFast.get(j, tempVal, configNumber);
193  std::cout << j << ") " << tempVal << std::endl;
194  }
195  std::cout << "and I am returning this: " <<
196  _rememberComplexFast.get(i, value, configNumber);
197  std::cout << " " << i << ") " << value << std::endl;
198  */
199  return _rememberComplexFast.get(i, value, configNumber);
200  }
bool get(int i, T &value, long int configNumber=1)
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47

◆ retrieveValue() [2/2]

virtual bool DalitzEvent::retrieveValue ( int  i,
double  value,
long int  configNumber = 1 
)
inlinevirtual

Implements IDalitzEvent.

Definition at line 206 of file DalitzEvent.h.

206  {
207  return _rememberDoubleFast.get(i, value, configNumber);
208  }
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48
bool get(int i, T &value, long int configNumber=1)

◆ s()

double DalitzEvent::s ( unsigned int  i,
unsigned int  j 
) const
virtual

Implements IDalitzEvent.

Definition at line 691 of file DalitzEvent.cpp.

691  {
693 }
const Permutation & currentPermutation() const
Definition: DalitzEvent.h:148
double s_intern(unsigned int i, unsigned int j) const

◆ s_intern()

double DalitzEvent::s_intern ( unsigned int  i,
unsigned int  j 
) const
protected

Definition at line 647 of file DalitzEvent.cpp.

647  {
648  if(_s[i][j] < 0) _s[i][j] = (p_intern(i) + p_intern(j)).Mag2();
649  return _s[i][j];
650 }
const TLorentzVector & p_intern(unsigned int i) const
std::vector< std::vector< double > > _s
Definition: DalitzEvent.h:57

◆ setAValue()

virtual void DalitzEvent::setAValue ( double  a)
inlinevirtual

Implements IDalitzEvent.

Definition at line 170 of file DalitzEvent.h.

170 {_aValue=a;}
double _aValue
Definition: DalitzEvent.h:50

◆ setGeneratorPdfRelativeToPhaseSpace()

void DalitzEvent::setGeneratorPdfRelativeToPhaseSpace ( double  gpdf)
virtual

Implements IDalitzEvent.

Definition at line 399 of file DalitzEvent.cpp.

399  {
401 }
double _generatorPdfRelativeToPhaseSpace
Definition: DalitzEvent.h:52

◆ setMomenta() [1/2]

bool DalitzEvent::setMomenta ( const std::vector< TLorentzVector > &  mumAndDgtr_p4)
protected

Definition at line 847 of file DalitzEvent.cpp.

847  {
848  _rememberPhaseSpace = -9999;
849  if(_pat.size() != p4List.size()){
850  return false;
851  }
852  _p.resize(p4List.size());
853  for(unsigned int i=0; i<p4List.size(); i++){
854  _p[i] = p4List[i];
855  }
856  return true;
857 }
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const

◆ setMomenta() [2/2]

bool DalitzEvent::setMomenta ( const std::vector< TVector3 > &  mumAndDgtr_p3)
protected

Definition at line 859 of file DalitzEvent.cpp.

859  {
860  _rememberPhaseSpace = -9999;
861 
862  _p.resize(p3List.size());
863  if(_pat.size() != p3List.size()){
864  return false;
865  }
866  for(unsigned int i=0; i< _pat.size(); i++){
867  const TVector3& p3 = p3List[i];
868  double m = _pat[i].mass();
869  _p[i].SetXYZM(p3.X(), p3.Y(), p3.Z(), m);
870  }
871  return true;
872 }
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42
unsigned int size() const
double m(unsigned int i) const

◆ setMothers3Momentum()

void DalitzEvent::setMothers3Momentum ( const TVector3 &  mp3)
virtual

Implements IDalitzEvent.

Definition at line 737 of file DalitzEvent.cpp.

737  {
738  TLorentzVector newMum4;
739  newMum4.SetXYZM(mp3.X(), mp3.Y(), mp3.Z(), _p[0].M());
740 
741  TVector3 boostToOldMumsRest(- _p[0].BoostVector());
742  TVector3 boostToNewMumsLab(newMum4.BoostVector());
743  for(unsigned int i=0; i < _p.size(); i++){
744  _p[i].Boost(boostToOldMumsRest); // to D0 rest frame
745  _p[i].Boost(boostToNewMumsLab); // to lab frame with new D0 momentum
746  }
747 
748 }
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ setP()

void DalitzEvent::setP ( unsigned int  i,
const TLorentzVector &  p4 
)
protected

Definition at line 874 of file DalitzEvent.cpp.

874  {
875  _rememberPhaseSpace = -9999;
876  if(i >= _p.size()){
877  _p.resize(i+1);
878  // cout << "had to resize" << endl;
879  }
880  _p[i] = p4;
881 }
double _rememberPhaseSpace
Definition: DalitzEvent.h:45
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43

◆ setPermutationIndex()

void DalitzEvent::setPermutationIndex ( int  i)
virtual

Implements IDalitzEvent.

Definition at line 839 of file DalitzEvent.cpp.

839  {
840  if(i >= numPermutations()) return;
842 }
int numPermutations() const
int _permutationIndex
Definition: DalitzEvent.h:62

◆ setValue() [1/2]

virtual void DalitzEvent::setValue ( int  i,
const std::complex< double > &  value,
long int  configNumber = 1 
)
inlinevirtual

Implements IDalitzEvent.

Definition at line 201 of file DalitzEvent.h.

202  {
203  //std::cout << "setting: " << i << ") " << value << std::endl;
204  _rememberComplexFast.set(i,value, configNumber);
205  }
void set(unsigned int i, const T &value, long int configNumber=1)
RememberAnythingFast< std::complex< double > > _rememberComplexFast
Definition: DalitzEvent.h:47

◆ setValue() [2/2]

virtual void DalitzEvent::setValue ( int  i,
double  value,
long int  configNumber = 1 
)
inlinevirtual

Implements IDalitzEvent.

Definition at line 209 of file DalitzEvent.h.

209  {
210  _rememberDoubleFast.set(i,value, configNumber);
211  }
void set(unsigned int i, const T &value, long int configNumber=1)
RememberAnythingFast< double > _rememberDoubleFast
Definition: DalitzEvent.h:48

◆ setValueInVector()

void DalitzEvent::setValueInVector ( unsigned int  i,
double  value 
)
virtual

Implements IDalitzEvent.

Definition at line 415 of file DalitzEvent.cpp.

415  {
416  if(_vectorOfValues.size() <= i) _vectorOfValues.resize(i+1);
417  _vectorOfValues[i]=value;
418 }
std::vector< double > _vectorOfValues
Definition: DalitzEvent.h:54

◆ setWeight()

void DalitzEvent::setWeight ( double  w)
virtual

Implements IDalitzEvent.

Definition at line 393 of file DalitzEvent.cpp.

393  {
394  _weight = w;
395 }
double _weight
Definition: DalitzEvent.h:51

◆ setWeightInVector()

void DalitzEvent::setWeightInVector ( unsigned int  i,
double  weight 
)
virtual

Implements IDalitzEvent.

Definition at line 420 of file DalitzEvent.cpp.

420  {
421  if(_vectorOfWeights.size() <= i) _vectorOfWeights.resize(i+1);
422  _vectorOfWeights[i]=weight;
423 }
std::vector< double > _vectorOfWeights
Definition: DalitzEvent.h:55

◆ shoutAndKill()

bool DalitzEvent::shoutAndKill ( )
protected

Definition at line 787 of file DalitzEvent.cpp.

787  {
788  cout << "FATAL ERROR in DalitzEvent:"
789  << "\n > pattern: \"" << _pat << "\""
790  << "\n > 4-vector list size: " << _p.size()
791  << endl;
792  throw "probably insonsistent pattern in DalitzEvent";
793 }
std::vector< TLorentzVector > _p
Definition: DalitzEvent.h:43
DalitzEventPattern _pat
Definition: DalitzEvent.h:42

◆ sij() [1/2]

double DalitzEvent::sij ( const std::vector< int > &  indices) const
virtual

Definition at line 709 of file DalitzEvent.cpp.

709  {
710  bool dbThis=false;
711 
712  if(1 == indices.size()){
713  return p(indices[0]).Mag2();
714  }else if(2 == indices.size()){
715  return s(indices[0], indices[1]);
716  }else{
717  if(dbThis){
718  TLorentzVector psum(0.0, 0.0, 0.0, 0.0);
719  for(unsigned int i=0; i< indices.size(); i++){
720  //cout << " p(" << i << ") " << p(indices[i]);
721  psum += p(indices[i]);
722  }
723  cout << "DalitzEvent sij: compare " << psum.Mag2()
724  << " = " << sijMap(indices) << endl;
725  }
726  return sijMap(indices);
727  }
728 }
virtual const TLorentzVector & p(unsigned int i) const
virtual double s(unsigned int i, unsigned int j) const
double sijMap(const std::vector< int > &indices) const

◆ sij() [2/2]

virtual double DalitzEvent::sij ( const MINT::PolymorphVector< int > &  indices) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 179 of file DalitzEvent.h.

179  {
180  return sij(indices.theVector());
181  }
virtual double sij(const std::vector< int > &indices) const
std::vector< T > & theVector()

◆ sijMap()

double DalitzEvent::sijMap ( const std::vector< int > &  indices) const
protected

Definition at line 665 of file DalitzEvent.cpp.

665  {
666  vector<int> intern_indices(indices);
667  for(unsigned int i=0; i < intern_indices.size(); i++){
668  intern_indices[i] = currentPermutation()[ indices[i] ];
669  }
670  return sijMap_intern(intern_indices);
671 }
const Permutation & currentPermutation() const
Definition: DalitzEvent.h:148
double sijMap_intern(const std::vector< int > &intern_indices) const

◆ sijMap_intern()

double DalitzEvent::sijMap_intern ( const std::vector< int > &  intern_indices) const
protected

Definition at line 655 of file DalitzEvent.cpp.

655  {
656  bool successFlag=false;
657  double returnVal = MINT::keyFinder(intern_indices, _sijMap
658  , (double) 0, successFlag);
659  if(! successFlag){
660  returnVal = evalsij_intern(intern_indices);
661  _sijMap[intern_indices] = returnVal;
662  }
663  return returnVal;
664 }
const Val & keyFinder(const Key &k, const std::map< Key, Val > &m, const Val &dummy, bool &successFlag)
Definition: Utils.h:20
double evalsij_intern(const std::vector< int > &intern_indices) const
std::map< std::vector< int >, double > _sijMap
Definition: DalitzEvent.h:59

◆ sijMax() [1/2]

double DalitzEvent::sijMax ( const std::vector< int > &  indices) const

Definition at line 733 of file DalitzEvent.cpp.

733  {
734  return _pat.sijMax(indices);
735 }
double sijMax(const MINT::PolymorphVector< int > &indices) const
DalitzEventPattern _pat
Definition: DalitzEvent.h:42

◆ sijMax() [2/2]

double DalitzEvent::sijMax ( const MINT::PolymorphVector< int > &  indices) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 219 of file DalitzEvent.h.

219  {
220  return sijMax(indices.theVector());}
std::vector< T > & theVector()
double sijMax(const std::vector< int > &indices) const

◆ sijMin() [1/2]

double DalitzEvent::sijMin ( const std::vector< int > &  indices) const

Definition at line 730 of file DalitzEvent.cpp.

730  {
731  return _pat.sijMin(indices);
732 }
double sijMin(const MINT::PolymorphVector< int > &indices) const
DalitzEventPattern _pat
Definition: DalitzEvent.h:42

◆ sijMin() [2/2]

double DalitzEvent::sijMin ( const MINT::PolymorphVector< int > &  indices) const
inlinevirtual

Implements IDalitzEvent.

Definition at line 215 of file DalitzEvent.h.

215  {
216  return sijMin(indices.theVector());}
double sijMin(const std::vector< int > &indices) const
std::vector< T > & theVector()

◆ singleParticleNtpArraySize()

int DalitzEvent::singleParticleNtpArraySize ( )
staticprotected

Definition at line 1011 of file DalitzEvent.cpp.

1011  {
1012  return 5;
1013 }

◆ t()

double DalitzEvent::t ( unsigned int  i,
unsigned int  j 
) const
virtual

Implements IDalitzEvent.

Definition at line 694 of file DalitzEvent.cpp.

694  {
696 }
const Permutation & currentPermutation() const
Definition: DalitzEvent.h:148
double t_intern(unsigned int i, unsigned int j) const

◆ t_intern()

double DalitzEvent::t_intern ( unsigned int  i,
unsigned int  j 
) const
protected

Definition at line 651 of file DalitzEvent.cpp.

651  {
652  if(_t[i][j] < 0) _t[i][j] = (p_intern(i) - p_intern(j)).Mag2();
653  return _t[i][j];
654 }
std::vector< std::vector< double > > _t
Definition: DalitzEvent.h:58
const TLorentzVector & p_intern(unsigned int i) const

Member Data Documentation

◆ _aValue

double DalitzEvent::_aValue
protected

Definition at line 50 of file DalitzEvent.h.

◆ _eventCounter

long int DalitzEvent::_eventCounter =0
staticprotected

Definition at line 37 of file DalitzEvent.h.

◆ _generatorPdfRelativeToPhaseSpace

double DalitzEvent::_generatorPdfRelativeToPhaseSpace
protected

Definition at line 52 of file DalitzEvent.h.

◆ _p

std::vector<TLorentzVector> DalitzEvent::_p
protected

Definition at line 43 of file DalitzEvent.h.

◆ _pat

DalitzEventPattern DalitzEvent::_pat
protected

Definition at line 42 of file DalitzEvent.h.

◆ _perm

Permutator DalitzEvent::_perm
protected

Definition at line 63 of file DalitzEvent.h.

◆ _permutationIndex

int DalitzEvent::_permutationIndex
protected

Definition at line 62 of file DalitzEvent.h.

◆ _rememberComplexFast

RememberAnythingFast<std::complex<double> > DalitzEvent::_rememberComplexFast
protected

Definition at line 47 of file DalitzEvent.h.

◆ _rememberDoubleFast

RememberAnythingFast<double > DalitzEvent::_rememberDoubleFast
protected

Definition at line 48 of file DalitzEvent.h.

◆ _rememberPhaseSpace

double DalitzEvent::_rememberPhaseSpace
mutableprotected

Definition at line 45 of file DalitzEvent.h.

◆ _rememberVectorCounter

long int DalitzEvent::_rememberVectorCounter =0
staticprotected

Definition at line 39 of file DalitzEvent.h.

◆ _s

std::vector< std::vector<double> > DalitzEvent::_s
mutableprotected

Definition at line 57 of file DalitzEvent.h.

◆ _sijMap

std::map<std::vector<int>, double> DalitzEvent::_sijMap
mutableprotected

Definition at line 59 of file DalitzEvent.h.

◆ _t

std::vector< std::vector<double> > DalitzEvent::_t
mutableprotected

Definition at line 58 of file DalitzEvent.h.

◆ _vectorOfValues

std::vector<double> DalitzEvent::_vectorOfValues
protected

Definition at line 54 of file DalitzEvent.h.

◆ _vectorOfWeights

std::vector<double> DalitzEvent::_vectorOfWeights
protected

Definition at line 55 of file DalitzEvent.h.

◆ _weight

double DalitzEvent::_weight
protected

Definition at line 51 of file DalitzEvent.h.

◆ ntpNameChars

const char DalitzEvent::ntpNameChars = { '#', '~', '\0' }
staticprotected

Definition at line 30 of file DalitzEvent.h.

◆ prtNameChars

const char DalitzEvent::prtNameChars = { '+', '-', '\0' }
staticprotected

Definition at line 31 of file DalitzEvent.h.


The documentation for this class was generated from the following files: