CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Decay3Body.cc
Go to the documentation of this file.
1 #include <cmath>
3 
4 #include "CLHEP/Units/GlobalSystemOfUnits.h"
5 #include "CLHEP/Units/GlobalPhysicalConstants.h"
8 #include "Math/GenVector/Rotation3D.h"
9 #include "Math/GenVector/EulerAngles.h"
10 #include "Math/GenVector/Boost.h"
11 #include <vector>
12 #include <cmath>
13 
14 #include "Randomize.hh"
15 
17 }
18 
20 }
21 
22 
23 void Decay3Body::doDecay(const G4LorentzVector & mother,
24  G4LorentzVector & daughter1,
25  G4LorentzVector & daughter2,
26  G4LorentzVector & daughter3) {
27 
28  double m0 = mother.m();
29  double m1 = daughter1.m();
30  double m2 = daughter2.m();
31  double m3 = daughter3.m();
32  double sumM2 = m0*m0 + m1*m1 + m2*m2 + m3*m3;
33  double tolerance = 1.0e-9;
34 
35  math::XYZTLorentzVectorD mmm(mother.px(),mother.py(),mother.pz(),mother.e());
36 
37  if (m0 < m1+m2+m3) {
38  std::cout << "Error: Daughters too heavy!" << std::endl;
39  std::cout << "M: " << m0/GeV <<
40  " < m1+m2+m3: " << m1/GeV + m2/GeV + m3/GeV << std::endl;
41  return;
42  } else {
43  double m2_12max = sqr(m0-m3);
44  double m2_12min = sqr(m1+m2);
45  double m2_23max = sqr(m0-m1);
46  double m2_23min = sqr(m2+m3);
47 
48  double x1,x2;
49  double m2_12 = 0.0;
50  double m2_23 = 0.0;
51  double E2_12,E3_12;
52  double m2_23max_12,m2_23min_12;
53 
54  do {
55 // Pick values for m2_12 and m2_23 uniformly:
56  x1 = G4UniformRand();
57  m2_12 = m2_12min + x1*(m2_12max-m2_12min);
58  x2 = G4UniformRand();
59  m2_23 = m2_23min + x2*(m2_23max-m2_23min);
60 
61 // From the allowed range of m2_23 (given m2_12), determine if the point is valid:
62 // (formulae taken from PDG booklet 2004 kinematics, page 305, Eqs. 38.22a+b)
63  E2_12 = (m2_12 - m1*m1 + m2*m2)/(2.0*sqrt(m2_12));
64  E3_12 = (m0*m0 - m2_12 - m3*m3)/(2.0*sqrt(m2_12));
65  m2_23max_12 = sqr(E2_12+E3_12)-sqr(sqrt(sqr(E2_12)-m2*m2)-sqrt(sqr(E3_12)-m3*m3));
66  m2_23min_12 = sqr(E2_12+E3_12)-sqr(sqrt(sqr(E2_12)-m2*m2)+sqrt(sqr(E3_12)-m3*m3));
67  } while ((m2_23 > m2_23max_12) || (m2_23 < m2_23min_12));
68 
69 // Determine the value of the third invariant mass squared:
70  double m2_13 = sumM2 - m2_12 - m2_23;
71 
72 // Calculate the energy and size of the momentum of the three daughters:
73  double e1 = (m0*m0 + m1*m1 - m2_23)/(2.0*m0);
74  double e2 = (m0*m0 + m2*m2 - m2_13)/(2.0*m0);
75  double e3 = (m0*m0 + m3*m3 - m2_12)/(2.0*m0);
76  double p1 = sqrt(e1*e1 - m1*m1);
77  double p2 = sqrt(e2*e2 - m2*m2);
78  double p3 = sqrt(e3*e3 - m3*m3);
79 
80 // Calculate cosine of the relative angles between the three daughters:
81  double cos12 = (m1*m1 + m2*m2 + 2.0*e1*e2 - m2_12)/(2.0*p1*p2);
82  double cos13 = (m1*m1 + m3*m3 + 2.0*e1*e3 - m2_13)/(2.0*p1*p3);
83  double cos23 = (m2*m2 + m3*m3 + 2.0*e2*e3 - m2_23)/(2.0*p2*p3);
84  if (fabs(cos12) > 1.0) std::cout << "Error: Undefined angle12!" << std::endl;
85  if (fabs(cos13) > 1.0) std::cout << "Error: Undefined angle13!" << std::endl;
86  if (fabs(cos23) > 1.0) std::cout << "Error: Undefined angle23!" << std::endl;
87 
88 // Find the four vectors of the particles in a chosen (i.e. simple) frame:
89  double xi = 2.0 * pi * G4UniformRand();
90  math::XYZVectorD q1(0.0,0.0,p1);
91  math::XYZVectorD q2( sin(acos(cos12))*cos(xi)*p2, sin(acos(cos12))*sin(xi)*p2,cos12*p2);
92  math::XYZVectorD q3(-sin(acos(cos13))*cos(xi)*p3,-sin(acos(cos13))*sin(xi)*p3,cos13*p3);
93 
94 // Rotate all three daughters momentum with the angles theta and phi:
95  double theta = acos(2.0 * G4UniformRand() - 1.0);
96  double phi = 2.0 * pi * G4UniformRand();
97  double psi = 2.0 * pi * G4UniformRand();
98 
99  ROOT::Math::EulerAngles ang(phi,theta,psi);
100  ROOT::Math::Rotation3D rot(ang);
101 
102  math::XYZVectorD q1rot = rot*q1;
103  math::XYZVectorD q2rot = rot*q2;
104  math::XYZVectorD q3rot = rot*q3;
105 
106  math::XYZTLorentzVectorD daughter1_orig(q1rot.X(),q1rot.Y(),q1rot.Z(),e1);
107  math::XYZTLorentzVectorD daughter2_orig(q2rot.X(),q2rot.Y(),q2rot.Z(),e2);
108  math::XYZTLorentzVectorD daughter3_orig(q3rot.X(),q3rot.Y(),q3rot.Z(),e3);
109 
110  ROOT::Math::Boost cmboost(mmm.BoostToCM());
111 
112  // Check of total angle and momentum:
113  if (acos(cos12)+acos(cos13)+acos(cos23)-2.0*pi > tolerance)
114  std::cout << "Error: Total angle not 2pi! " <<
115  acos(cos12)+acos(cos13)+acos(cos23)-2.0*pi << std::endl;
116  if (fabs(daughter1_orig.px()+daughter2_orig.px()+daughter3_orig.px())/GeV > tolerance)
117  std::cout << "Error: Total 3B Px not conserved! " <<
118  (daughter1_orig.px()+daughter2_orig.px()+daughter3_orig.px())/GeV << std::endl;
119  if (fabs(daughter1_orig.py()+daughter2_orig.py()+daughter3_orig.py())/GeV > tolerance)
120  std::cout << "Error: Total 3B Py not conserved! " <<
121  (daughter1_orig.py()+daughter2_orig.py()+daughter3_orig.py())/GeV << std::endl;
122  if (fabs(daughter1_orig.pz()+daughter2_orig.pz()+daughter3_orig.pz())/GeV > tolerance)
123  std::cout << "Error: Total 3B Pz not conserved! " <<
124  (daughter1.pz()+daughter2.pz()+daughter3.pz())/GeV << std::endl;
125 
126  // Boost the daughters back to the frame of the mother:
127 
128  math::XYZTLorentzVectorD temp1(cmboost(daughter1_orig));
129  math::XYZTLorentzVectorD temp2(cmboost(daughter2_orig));
130  math::XYZTLorentzVectorD temp3(cmboost(daughter3_orig));
131 
132  daughter1.setPx(temp1.Px());
133  daughter1.setPy(temp1.Py());
134  daughter1.setPz(temp1.Pz());
135  daughter1.setE(temp1.e());
136 
137  daughter2.setPx(temp2.Px());
138  daughter2.setPy(temp2.Py());
139  daughter2.setPz(temp2.Pz());
140  daughter2.setE(temp2.e());
141 
142  daughter3.setPx(temp3.Px());
143  daughter3.setPy(temp3.Py());
144  daughter3.setPz(temp3.Pz());
145  daughter3.setE(temp3.e());
146 
147  return;
148  }
149 }
150 
151 
152 double Decay3Body::sqr(double a) {
153  return a*a;
154 }
155 
void doDecay(const G4LorentzVector &mother, G4LorentzVector &daughter1, G4LorentzVector &daughter2, G4LorentzVector &daughter3)
Definition: Decay3Body.cc:23
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double sqr(double a)
Definition: Decay3Body.cc:152
Geom::Theta< T > theta() const
double q2[4]
Definition: TauolaWrapper.h:88
std::map< std::string, int, std::less< std::string > > psi
const Double_t pi
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:9
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double p2[4]
Definition: TauolaWrapper.h:90
AlgebraicVector EulerAngles
Definition: Definitions.h:36
double q1[4]
Definition: TauolaWrapper.h:87
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
double p3[4]
Definition: TauolaWrapper.h:91
Definition: DDAxes.h:10