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
Smear.h
Go to the documentation of this file.
1 
10 #ifndef INCLUDE_EICSMEAR_SMEAR_SMEAR_H_
11 #define INCLUDE_EICSMEAR_SMEAR_SMEAR_H_
12 
13 #include <cmath>
14 #include <fstream>
15 #include <iostream>
16 #include <memory>
17 #include <sstream>
18 
19 #include <TF1.h>
20 #include <TF2.h>
21 #include <TLorentzVector.h>
22 #include <TMath.h>
23 #include <TRandom3.h>
24 #include <TROOT.h>
25 #include <TString.h>
26 #include <TSystem.h>
27 #include <TVector2.h>
28 
34 
35 namespace Smear {
36 
43 enum KinType {
45 };
46 
48 enum EGenre {
50 };
51 
53 enum ECharge {
55 };
56 
68 inline int PGenre(const erhic::VirtualParticle& prt) {
69  int genre(kAll);
70  const int id = abs(prt.Id()); // Sign doesn't matter
71  if (1 == prt.GetStatus()) { // Only check stable particles.
72  if (id == 11 || id == 22) {
73  genre = kElectromagnetic;
74  } else if (id >110) {
75  genre = kHadronic;
76  } // if
77  } // if
78  return genre;
79 }
80 
85 inline double FixTheta(double theta) {
86  while (theta < 0. || theta > TMath::Pi()) {
87  if (theta < 0.) {
88  theta = -theta;
89  } // if
90  if (theta > TMath::Pi()) {
91  theta = TMath::TwoPi() - theta;
92  } // if
93  }
94  return theta;
95 }
96 
101 inline double FixPhi(double phi) {
102  return TVector2::Phi_0_2pi(phi);
103 }
104 
108 inline double GetVariable(const erhic::VirtualParticle& prt, KinType kin) {
109  double z(0.);
110  switch (kin) {
111  case kE:
112  z = prt.GetE(); break;
113  case kP:
114  z = prt.GetP(); break;
115  case kTheta:
116  z = prt.GetTheta(); break;
117  case kPhi:
118  z = prt.GetPhi(); break;
119  case kPz:
120  z = prt.GetPz(); break;
121  case kPt:
122  z = prt.GetPt(); break;
123  default:
124  break;
125  } // switch
126  return z;
127 }
128 
132 inline void SetVariable(ParticleMCS &prt, double z, KinType kin) {
133  switch (kin) {
134  case kE:
135  prt.SetE(z); break;
136  case kP:
137  prt.SetP(z); break;
138  case kTheta:
139  prt.SetTheta(z); break;
140  case kPhi:
141  prt.SetPhi(z); break;
142  case kPz:
143  prt.SetPz(z); break;
144  case kPt:
145  prt.SetPt(z); break;
146  default:
147  break;
148  } // switch
149 }
150 
155 inline void HandleBogusValues(ParticleMCS &prt, KinType kin) {
156  double fault(0.);
157  if (kE == kin && prt.GetE() < 0.) {
158  prt.SetE(fault);
159  } else if (kP == kin && prt.GetP() < 0.) {
160  prt.SetP(fault);
161  } else if (kPt == kin && prt.GetPt() < 0.) {
162  prt.SetPt(fault);
163  } else if (kPz == kin && prt.GetPz() < 0.) {
164  prt.SetPz(fault);
165  } // if
166 }
167 
168 inline void HandleBogusValues(ParticleMCS& prt) {
169  double fault = NAN; // -999.;
170  if (prt.GetE() < 0.) {
171  prt.SetE(fault);
172  } // if
173  if (prt.GetP() < 0.) {
174  prt.SetP(fault);
175  } // if
176  if (prt.GetPt() < 0.) {
177  prt.SetPt(fault);
178  } // if
179  if (prt.GetPz() < 0.) {
180  prt.SetPz(fault);
181  } // if
182 }
183 
184 inline bool IsCoreType(KinType kin) {
185  if (kin == kE || kin == kP || kin == kTheta || kin == kPhi) return true;
186  return false;
187 }
188 
189 int ParseInputFunction(TString &s, KinType &kin1, KinType &kin2);
190 
191 } // namespace Smear
192 
193 #endif // INCLUDE_EICSMEAR_SMEAR_SMEAR_H_
Smear::ParseInputFunction
int ParseInputFunction(TString &s, KinType &kin1, KinType &kin2)
Definition: Smear.cxx:18
Smear::kPz
@ kPz
Definition: Smear.h:44
Smear
Definition: Acceptance.cxx:16
Smear::ParticleMCS::SetPhi
virtual void SetPhi(Double_t)
Definition: ParticleMCS.h:245
Smear::FixPhi
double FixPhi(double phi)
Fix an azimuthal angle so that it lies within [0,2*pi).
Definition: Smear.h:101
erhic::VirtualParticle::GetPz
virtual Double_t GetPz() const =0
Returns the z component of 3-momentum.
Smear::ParticleMCS::SetPt
virtual void SetPt(Double_t)
Definition: ParticleMCS.h:237
Smear::ParticleMCS::SetPz
virtual void SetPz(Double_t)
Definition: ParticleMCS.h:241
Smear::kAllCharges
@ kAllCharges
Definition: Smear.h:54
Smear::HandleBogusValues
void HandleBogusValues(ParticleMCS &prt, KinType kin)
This dictates how the namespace deals with positive definite variables which have been smeared to neg...
Definition: Smear.h:155
Smear::ParticleMCS::GetP
virtual Double_t GetP() const
Returns the total momentum (GeV).
Definition: ParticleMCS.h:213
Smear::kAll
@ kAll
Definition: Smear.h:49
ParticleMCS.h
Smear::kE
@ kE
Definition: Smear.h:44
erhic::VirtualParticle::GetP
virtual Double_t GetP() const =0
Returns the magnitude of 3-momentum (GeV).
Smear::kPt
@ kPt
Definition: Smear.h:44
Smear::kTheta
@ kTheta
Definition: Smear.h:44
Smear::ParticleMCS::GetE
virtual Double_t GetE() const
Returns the energy of the particle in the lab frame.
Definition: ParticleMCS.h:197
Smear::ParticleMCS::SetP
virtual void SetP(Double_t)
Definition: ParticleMCS.h:233
Smear::ParticleMCS::SetE
virtual void SetE(Double_t)
Definition: ParticleMCS.h:229
erhic::VirtualParticle::GetStatus
virtual UShort_t GetStatus() const =0
A general "status" code for the particle (definition depends on implementation).
erhic::VirtualParticle::GetE
virtual Double_t GetE() const =0
Returns total energy.
Smear::SetVariable
void SetVariable(ParticleMCS &prt, double z, KinType kin)
Stores z in the ParticleS.K where K is the kinematic variable associated with kin.
Definition: Smear.h:132
Smear::ParticleMCS::GetPz
virtual Double_t GetPz() const
Returns the z component of 3-momentum.
Definition: ParticleMCS.h:193
Smear::FixTheta
double FixTheta(double theta)
Fix a polar angle so that it lies within [0,pi].
Definition: Smear.h:85
VirtualParticle.h
erhic::VirtualParticle
Abstract base class for a general particle.
Definition: VirtualParticle.h:23
Smear::IsCoreType
bool IsCoreType(KinType kin)
Definition: Smear.h:184
Smear::kCharged
@ kCharged
Definition: Smear.h:54
Particle.h
Smear::kPhi
@ kPhi
Definition: Smear.h:44
Smear::ParticleMCS
A smeared Monte Carlo particle.
Definition: ParticleMCS.h:27
Smear::ParticleMCS::SetTheta
virtual void SetTheta(Double_t)
Definition: ParticleMCS.h:249
Smear::kP
@ kP
Definition: Smear.h:44
erhic::VirtualParticle::GetPt
virtual Double_t GetPt() const =0
Returns momentum perpendicular to the beam direction.
erhic::VirtualParticle::GetTheta
virtual Double_t GetTheta() const =0
Returns the polar angle in the range [0, pi] radians.
Smear::kInvalidKinType
@ kInvalidKinType
Definition: Smear.h:44
Smear::ECharge
ECharge
Particle charged.
Definition: Smear.h:53
Smear::KinType
KinType
Enumerator listing particle wise kinematic variables.
Definition: Smear.h:43
Smear::PGenre
int PGenre(const erhic::VirtualParticle &prt)
Determine particle "genre".
Definition: Smear.h:68
erhic::VirtualParticle::GetPhi
virtual Double_t GetPhi() const =0
Returns the polar angle in the range [0, 2pi] radians.
Smear::GetVariable
double GetVariable(const erhic::VirtualParticle &prt, KinType kin)
Returns the kinematic variable associated with kin from the input particle.
Definition: Smear.h:108
ParticleIdentifier.h
Smear::kHadronic
@ kHadronic
Definition: Smear.h:49
Smear::EGenre
EGenre
Classes of particles.
Definition: Smear.h:48
Smear::kElectromagnetic
@ kElectromagnetic
Definition: Smear.h:49
Kinematics.h
Smear::ParticleMCS::GetPt
virtual Double_t GetPt() const
Returns momentum transverse to the beam direction.
Definition: ParticleMCS.h:205
Smear::kNeutral
@ kNeutral
Definition: Smear.h:54
erhic::VirtualParticle::Id
virtual Pid Id() const =0
Returns identity information for the Particle species.