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
ParticleMCS.cxx
Go to the documentation of this file.
1 
11 
12 #include <iomanip>
13 #include <iostream>
14 
15 #include <TMath.h>
16 
17 namespace Smear {
18 
20 : status(0)
21 , id(0)
22 , px(0.)
23 , py(0.)
24 , pz(0.)
25 , E(0.)
26 , pt(0.)
27 , p(0.)
28 , theta(0.)
29 , phi(0.) {
30 }
31 
32 ParticleMCS::ParticleMCS(const TLorentzVector& ep, int pdg, int stat)
33 : status(stat)
34 , id(pdg)
35 , px(ep.Px())
36 , py(ep.Py())
37 , pz(ep.Pz())
38 , E(ep.E())
39 , pt(ep.Pt())
40 , p(ep.P())
41 , theta(ep.Theta())
42 , phi(ep.Phi()) {
43 }
44 
46 }
47 
48 TLorentzVector ParticleMCS::Get4Vector() const {
49  return TLorentzVector(px, py, pz, E);
50 }
51 
52 void ParticleMCS::Print(Option_t* /* unused */) const {
53  // Store initial flags for cout so we can restore them later.
54  std::ios_base::fmtflags ff = std::cout.flags();
55  // Configure flags for output
56  std::cout.setf(std::ios_base::fixed, std::ios_base::basefield);
57  std::cout.precision(4);
58  // Output values.
59  std::cout <<
60  std::setw(3) << status <<
61  std::setw(12) << id <<
62  std::setw(10) << px <<
63  std::setw(10) << py <<
64  std::setw(10) << pz <<
65  std::setw(10) << E << std::endl;
66  // Restore initial cout flags.
67  std::cout.flags(ff);
68 }
69 
70 Double_t ParticleMCS::GetEta() const {
71  // The default value of -19 is used in the main eRHIC code,
72  // so use that for consistency.
73  double eta(-19.);
74  const double theta = GetTheta();
75  if (theta > 0. && theta < TMath::Pi() && !TMath::IsNaN(theta)) {
76  eta = -log(tan(theta / 2.));
77  } // if
78  return eta;
79 }
80 
81 Double_t ParticleMCS::GetRapidity() const {
82  // The default value of -19 is used in the main eRHIC code,
83  // so use that for consistency.
84  double y(-19.);
85  // Be careful here, as the energy and momentum can become
86  // smeared such that the stored E < pz, giving NaN when
87  // taking a log of a negative.
88  // In that case, return the default value.
89  // E > pz already takes care of the case if E is NaN, so
90  // no need to check that again.
91  if (E > pz && !TMath::IsNaN(pz)) {
92  y = 0.5 * log((E + pz) / (E - pz));
93  } // if
94  return y;
95 }
96 
97 } // namespace Smear
virtual ~ParticleMCS()
Definition: ParticleMCS.cxx:45
virtual TLorentzVector Get4Vector() const
Definition: ParticleMCS.cxx:48
Double32_t E
Energy of particle.
Definition: ParticleMCS.h:176
virtual void Print(Option_t *="") const
Definition: ParticleMCS.cxx:52
Double32_t py
y component of particle momentum
Definition: ParticleMCS.h:174
Double32_t px
x component of particle momentum
Definition: ParticleMCS.h:173
virtual Double_t GetTheta() const
Definition: ParticleMCS.h:217
virtual Double_t GetRapidity() const
Definition: ParticleMCS.cxx:81
UShort_t status
Status code.
Definition: ParticleMCS.h:171
Double32_t pz
z component of particle momentum
Definition: ParticleMCS.h:175
Double32_t theta
Polar angle.
Definition: ParticleMCS.h:179
virtual Double_t GetEta() const
Definition: ParticleMCS.cxx:70