eic-smear  1.0.3
A collection of ROOT classes for Monte Carlo events and a fast-smearing code simulating detector effects for the Electron-Ion Collider task force
ParticleMC.cxx
Go to the documentation of this file.
1 
11 
12 #include <cmath>
13 #include <limits>
14 
15 #include <TMCParticle.h>
16 #include <TVector2.h>
17 
18 namespace erhic {
19 namespace hadronic {
20 
22 : KS(-1)
23 , orig(-1)
24 , id(std::numeric_limits<Int_t>::min())
25 , px(NAN)
26 , py(NAN)
27 , pz(NAN)
28 , E(NAN)
29 , p(NAN)
30 , m(NAN)
31 , pt(NAN)
32 , theta(NAN)
33 , phi(NAN)
34 , rapidity(NAN)
35 , eta(NAN)
36 , xFeynman(NAN)
37 , xv(NAN)
38 , yv(NAN)
39 , zv(NAN) {
40 }
41 
42 ParticleMC::ParticleMC(const TMCParticle& mc)
43 : KS(mc.GetKS())
44 , orig(-1)
45 , id(mc.GetKF())
46 , px(mc.GetPx())
47 , py(mc.GetPy())
48 , pz(mc.GetPz())
49 , E(mc.GetEnergy())
50 , m(mc.GetMass())
51 , pt(NAN)
52 , theta(NAN)
53 , rapidity(NAN)
54 , eta(NAN)
55 , xFeynman(NAN)
56 , xv(mc.GetVx())
57 , yv(mc.GetVy())
58 , zv(mc.GetVz()) {
59  TLorentzVector v(mc.GetPx(), mc.GetPy(), mc.GetPz(), mc.GetEnergy());
60  p = v.P();
61  pt = v.Pt();
62  theta = v.Theta();
63  phi = TVector2::Phi_0_2pi(v.Phi());
64  rapidity = v.Rapidity();
65  // Protect against attempting pseudorapidity calculation for
66  // particles with angle close to the beam.
67  if (v.Pt() > 1.e-5) {
68  eta = v.PseudoRapidity();
69  } else {
70  eta = std::numeric_limits<Double32_t>::infinity();
71  } // if
72 }
73 
74 ParticleMC::ParticleMC(const TLorentzVector& ep, const TVector3& v, int pdg,
75  int status, int parent)
76 : KS(status)
77 , orig(parent)
78 , id(pdg)
79 , px(ep.Px())
80 , py(ep.Py())
81 , pz(ep.Pz())
82 , E(ep.E())
83 , p(ep.P())
84 , m(ep.M())
85 , pt(ep.Pt())
86 , theta(ep.Theta())
87 , phi(TVector2::Phi_0_2pi(ep.Phi()))
88 , rapidity(ep.Rapidity())
89 , eta(ep.PseudoRapidity())
90 , xFeynman(NAN)
91 , xv(v.X())
92 , yv(v.Y())
93 , zv(v.Z()) {
94 }
95 
96 
97 void ParticleMC::SetParentIndex(UShort_t i) {
98  orig = i;
99 }
100 
101 void ParticleMC::Set4Vector(const TLorentzVector& v) {
102  E = v.Energy();
103  px = v.Px();
104  py = v.Py();
105  pz = v.Pz();
106  p = v.P();
107  m = v.M();
108  pt = v.Pt();
109  theta = v.Theta();
110  phi = TVector2::Phi_0_2pi(v.Phi());
111  rapidity = v.Rapidity();
112  eta = v.PseudoRapidity();
113 }
114 
115 } // namespace hadronic
116 } // namespace erhic
Double32_t p
Magnitude of momentum (GeV/c)
Definition: ParticleMC.h:177
UShort_t orig
I of parent particle.
Definition: ParticleMC.h:171
Double32_t px
x component of momentum (GeV/c)
Definition: ParticleMC.h:173
Double32_t E
Total energy (GeV)
Definition: ParticleMC.h:176
Double32_t py
y component of momentum (GeV/c)
Definition: ParticleMC.h:174
Double32_t phi
Angle of azimuth (radians [0, 2pi])
Definition: ParticleMC.h:181
Double32_t theta
Polar angle (radians [0, pi])
Definition: ParticleMC.h:180
virtual void SetParentIndex(UShort_t)
Definition: ParticleMC.cxx:97
Double32_t pz
z component of momentum (GeV/c)
Definition: ParticleMC.h:175
virtual void Set4Vector(const TLorentzVector &)
Definition: ParticleMC.cxx:101
Double32_t pt
Momentum transverse to the beam direction (GeV/c)
Definition: ParticleMC.h:179
Double32_t rapidity
Rapidity.
Definition: ParticleMC.h:182
Double32_t m
Invariant mass (GeV/c2)
Definition: ParticleMC.h:178
Double32_t eta
Pseudorapidity.
Definition: ParticleMC.h:183