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 {
44  kE, kP, kTheta, kPhi, kPz, kPt, kInvalidKinType
45 };
46 
48 enum EGenre {
49  kAll = 0, kElectromagnetic = 1, kHadronic = 2
50 };
51 
53 enum ECharge {
54  kNeutral, kCharged, kAllCharges
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_
virtual Double_t GetPhi() const =0
virtual Double_t GetP() const =0
virtual Double_t GetE() const =0
virtual UShort_t GetStatus() const =0
virtual Double_t GetPz() const =0
virtual Double_t GetPt() const =0
virtual Pid Id() const =0
virtual Double_t GetTheta() const =0
Abstract base class for a general particle.