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 <iostream>
13 #include <limits>
14 #include <sstream>
15 #include <stdexcept>
16 #include <string>
17 
18 #include <TLorentzRotation.h>
19 #include <TParticlePDG.h>
20 #include <TRotation.h>
21 
23 #include "eicsmear/functions.h"
24 
25 namespace {
26 
27 /*
28  Returns the boost to transform to the rest frame of the vector "rest".
29  If z is non-NULL, rotate the frame so that z AFTER BOOSTING
30  defines the positive z direction of that frame.
31  e.g. To boost a gamma-p interaction from the lab frame to the
32  proton rest frame with the virtual photon defining z:
33  computeBoost(protonLab, photonLab);
34  */
35 TLorentzRotation computeBoost(const TLorentzVector& rest,
36  const TLorentzVector* z) {
37  TLorentzRotation toRest(-(rest.BoostVector()));
38  if (z) {
39  TRotation rotate;
40  TLorentzVector boostedZ(*z);
41  boostedZ *= toRest;
42  rotate.SetZAxis(boostedZ.Vect());
43  // We need the rotation of the frame, so take the inverse.
44  // See the TRotation documentation for more explanation.
45  rotate = rotate.Inverse();
46  toRest.Transform(rotate);
47  } // if
48  return toRest;
49 }
50 
51 } // anonymous namespace
52 
53 namespace erhic {
54 
55 ParticleMC::ParticleMC(const std::string& line)
56 : I(-1)
57 , KS(-1)
58 , id(std::numeric_limits<Int_t>::min())
59 , orig(-1)
60 , daughter(-1)
61 , ldaughter(-1)
62 , px(0.)
63 , py(0.)
64 , pz(0.)
65 , E(0.)
66 , m(0.)
67 , pt(0.)
68 , xv(0.)
69 , yv(0.)
70 , zv(0.)
71 , parentId(std::numeric_limits<Int_t>::min())
72 , p(0.)
73 , theta(0.)
74 , phi(0.)
75 , rapidity(0.)
76 , eta(0.)
77 , z(0.)
78 , xFeynman(0.)
79 , thetaGamma(0.)
80 , ptVsGamma(0.)
81 , phiPrf(0.) {
82  // Initialise to nonsense values to make input errors easy to spot
83  if (!line.empty()) {
84  static std::stringstream ss;
85  ss.str("");
86  ss.clear();
87  ss << line;
88  ss >>
89  I >> KS >> id >> orig >> daughter >> ldaughter >>
90  px >> py >> pz >> E >> m >> xv >> yv >> zv;
91  // We should have no stream errors and should have exhausted
92  // the whole of the stream filling the particle.
93  if (ss.fail() || !ss.eof()) {
94  throw std::runtime_error("Bad particle input: " + line);
95  } // if
97  } // if
98 }
99 
101 }
102 
103 void ParticleMC::Print(Option_t* /* option */) const {
104  std::cout << I << '\t' << KS << '\t' << id << '\t' << orig << '\t' <<
105  daughter << '\t' << ldaughter << '\t' << px << '\t' << py << '\t' << pz
106  << '\t' << E << '\t' << m << '\t' << xv << '\t' << yv << '\t' << zv <<
107  std::endl;
108 }
109 
111  // Calculate quantities that depend only on the properties already read.
112  pt = sqrt(pow(px, 2.) + pow(py, 2.));
113  p = sqrt(pow(pt, 2.) + pow(pz, 2.));
114  // Rapidity and pseudorapidity
115  Double_t Epluspz = E + pz;
116  Double_t Eminuspz = E - pz;
117  Double_t Ppluspz = p + pz;
118  Double_t Pminuspz = p - pz;
119  if (Eminuspz <= 0.0 || Pminuspz == 0.0 ||
120  Ppluspz == 0.0 || Epluspz <= 0.0) {
121  // Dummy values to avoid zero or infinite arguments in calculations
122  rapidity = -19.;
123  eta = -19.;
124  } else {
125  rapidity = 0.5 * log(Epluspz / Eminuspz);
126  eta = 0.5 * log(Ppluspz / Pminuspz);
127  } // if
128  theta = atan2(pt, pz);
129  phi = TVector2::Phi_0_2pi(atan2(py, px));
130 }
131 
133  try {
134  // Get the beam hadon, beam lepton and exchange boson.
135  const TLorentzVector& hadron = event.BeamHadron()->Get4Vector();
136  const TLorentzVector& lepton = event.ScatteredLepton()->Get4Vector();
137  const TLorentzVector& boson = event.ExchangeBoson()->Get4Vector();
138  // Calculate z using the 4-vector definition,
139  // so we don't care about frame of reference.
140  z = hadron.Dot(Get4Vector()) / hadron.Dot(boson);
141  // Calculate properties in the proton rest frame.
142  // We want pT and angle with respect to the virtual photon,
143  // so use that to define the z axis.
144  TLorentzRotation toHadronRest = computeBoost(hadron, &boson);
145  // Boost this particle to the proton rest frame and calculate its
146  // pT and angle with respect to the virtual photon:
147  TLorentzVector p = (Get4Vector() *= toHadronRest);
148  thetaGamma = p.Theta();
149  ptVsGamma = p.Pt();
150  // Calculate phi angle around virtual photon according
151  // to the HERMES convention.
152  TLorentzVector bosonPrf = (TLorentzVector(boson) *= toHadronRest);
153  TLorentzVector leptonPrf = (TLorentzVector(lepton) *= toHadronRest);
154  phiPrf = computeHermesPhiH(p, leptonPrf, bosonPrf);
155  // Feynman x with xF = 2 * pz / W in the boson-hadron CM frame.
156  // First boost to boson-hadron centre-of-mass frame.
157  // Use the photon to define the z direction.
158  TLorentzRotation boost = computeBoost(boson + hadron, &boson);
159  xFeynman = 2. * (Get4Vector() *= boost).Pz() / sqrt(event.GetW2());
160  // Determine the PDG code of the parent particle, if the particle
161  // has a parent and the parent is present in the particle array.
162  // The index of the particles from the Monte Carlo runs from [1,N]
163  // while the index in the array runs from [0,N-1], so subtract 1
164  // from the parent index to find its position.
165  if (event.GetNTracks() > unsigned(orig - 1)) {
166  parentId = event.GetTrack(orig - 1)->Id();
167  } // if
168  } // try
169  catch(std::exception& e) {
170  std::cerr <<
171  "Exception in Particle::ComputeEventDependentQuantities: " <<
172  e.what() << std::endl;
173  } // catch
174 }
175 
176 TLorentzVector ParticleMC::Get4Vector() const {
177  return TLorentzVector(px, py, pz, E);
178 }
179 
181  return static_cast<EventMC*>(event.GetObject());
182 }
183 
184 const ParticleMC* ParticleMC::GetChild(UShort_t u) const {
185  // Look up this particle's child via the event
186  // containing it and the index of the child in that event.
187  if (!GetEvent()) {
188  return NULL;
189  } // if
190  // index is in the range [1,N]
191  unsigned idx = daughter + u;
192  if (daughter < 1 || // If first daughter index = 0, it has no children
193  u >= GetNChildren()) { // Insufficient children
194  return NULL;
195  } // if
196  --idx; // Convert [1,N] --> [0,N)
197  const ParticleMC* p(NULL);
198  // Check this index is within the # of particles in the event
199  if (idx < GetEvent()->GetNTracks()) {
200  p = GetEvent()->GetTrack(idx);
201  } // if
202  return p;
203 }
204 
206  // Look up this particle's parent via the event
207  // containing it and the index of the parent in that event.
208  const ParticleMC* p(NULL);
209  if (GetEvent()) {
210  if (GetEvent()->GetNTracks() >= GetParentIndex()) {
211  p = GetEvent()->GetTrack(GetParentIndex() - 1);
212  } // if
213  } // if
214  return p;
215 }
216 
217 Bool_t ParticleMC::HasChild(Int_t pdg) const {
218  for (UInt_t i(0); i < GetNChildren(); ++i) {
219  if (!GetChild(i)) {
220  continue;
221  } // if
222  if (pdg == GetChild(i)->Id()) {
223  return true;
224  } // if
225  } // for
226  return false;
227 }
228 
230  double p_(0.), e_(0.), px_(ptVsGamma), py_(0.), pz_(0.);
231  // Calculate mangitude of momentum from pT and polar angle in
232  // hadron-boson frame. If theta is ~parallel to the beam just set
233  // p to whatever pT is (not that it makes all that much sense).
234  if (thetaGamma > 1.e-6) {
235  p_ = ptVsGamma / sin(thetaGamma);
236  } // if
237  // Deal with virtual particles later, so check if particle is off mass-shell
238  if (!(m < 0.)) {
239  e_ = sqrt(pow(p_, 2.) + pow(m, 2.));
240  } else {
241  e_ = sqrt(pow(p_, 2.) - pow(m, 2.));
242  if (TMath::IsNaN(e_)) {
243  e_ = 0.;
244  } // if
245  } // if
246  // Calculate pZ from pT and theta, unless it's very close to the beam,
247  // in which case set it to p.
248  if (thetaGamma > 1.e-6) {
249  pz_ = ptVsGamma / tan(thetaGamma);
250  } // if
251  // Calculate px and py, unless a dummy phi value is present.
252  if (phiPrf > -100.) {
253  px_ = ptVsGamma * cos(phiPrf);
254  py_ = ptVsGamma * sin(phiPrf);
255  } // if
256  // If we ended up with no energy, it's likely ths is the exchange boson,
257  // as nothing will have happened above. If so, try to reference the event
258  // record to get the necessary information, as we don't have enough
259  // in the particle itself.
260  // Note that this appears not to work in a TTree as ROOT doesn't read
261  // the event when it reads the particle, though I'm not certain of the
262  // exact cause.
263  if (m < 0. && GetEvent()) {
264  if (GetEvent()->ExchangeBoson()) {
265  // Don't check if GetEvent()->ExchangeBoson() == this, in case this
266  // particle is a copy of the event track that isn't part of a copied
267  // event.
268  if (GetEvent()->ExchangeBoson()->GetIndex() == GetIndex()) {
269  // If it's the exchange boson just set E == nu.
270  e_ = GetEvent()->GetNu();
271  // Calculate p, careful about negative mass.
272  p_ = sqrt(pow(e_, 2.) + pow(m, 2.));
273  pz_ = p_ * cos(thetaGamma);
274  } // if
275  } // if
276  } // if
277  return TLorentzVector(px_, py_, pz_, e_);
278 }
279 
281  event = e;
282 }
283 
284 void ParticleMC::Set4Vector(const TLorentzVector& v) {
285  E = v.Energy();
286  px = v.Px();
287  py = v.Py();
288  pz = v.Pz();
289  m = v.M();
290  ComputeDerivedQuantities(); // Rapidity etc
291  // If an event reference is set, recalculate event-dependent quantities
292  // like z, xF
293  EventMC* ev = static_cast<EventMC*>(event.GetObject());
294  if (ev) {
296  } // if
297 }
298 
299 void ParticleMC::SetVertex(const TVector3& v) {
300  xv = v.X();
301  yv = v.Y();
302  zv = v.Z();
303 }
304 
305 } // namespace erhic
Double32_t thetaGamma
Definition: ParticleMC.h:399
Double32_t m
Invariant mass of particle.
Definition: ParticleMC.h:382
UShort_t ldaughter
I of last child particle.
Definition: ParticleMC.h:356
Double32_t px
x component of particle momentum
Definition: ParticleMC.h:358
Double32_t zv
z coordinate of particle production vertex
Definition: ParticleMC.h:386
Double32_t E
Energy of particle.
Definition: ParticleMC.h:381
void SetEvent(EventMC *event)
Definition: ParticleMC.cxx:280
Double32_t eta
Pseudorapidity of particle.
Definition: ParticleMC.h:395
Double32_t ptVsGamma
Definition: ParticleMC.h:401
Double32_t py
y component of particle momentum
Definition: ParticleMC.h:359
UShort_t I
Particle index in event.
Definition: ParticleMC.h:351
virtual void ComputeDerivedQuantities()
Definition: ParticleMC.cxx:110
double computeHermesPhiH(const TLorentzVector &hadronInPrf, const TLorentzVector &leptonInPrf, const TLorentzVector &photonInPrf)
Definition: functions.cxx:52
virtual UInt_t GetNTracks() const
Definition: EventMC.h:221
Double32_t xv
x coordinate of particle production vertex
Definition: ParticleMC.h:384
virtual TLorentzVector Get4Vector() const
Definition: ParticleMC.cxx:176
ParticleMC(const std::string &="")
Definition: ParticleMC.cxx:55
const EventMC * GetEvent() const
Definition: ParticleMC.cxx:180
virtual void Print(Option_t *="") const
Definition: ParticleMC.cxx:103
UShort_t orig
I of parent particle.
Definition: ParticleMC.h:354
Double32_t pz
z component of particle momentum
Definition: ParticleMC.h:360
virtual void SetVertex(const TVector3 &)
Definition: ParticleMC.cxx:299
Double32_t yv
y coordinate of particle production vertex
Definition: ParticleMC.h:385
Double32_t phi
Azimuthal angle.
Definition: ParticleMC.h:393
Double32_t pt
Transverse momentum of particle.
Definition: ParticleMC.h:383
Double32_t rapidity
Rapidity of particle.
Definition: ParticleMC.h:394
virtual const ParticleMC * GetParent() const
Definition: ParticleMC.cxx:205
virtual ~ParticleMC()
Definition: ParticleMC.cxx:100
virtual void Set4Vector(const TLorentzVector &)
Definition: ParticleMC.cxx:284
Double32_t theta
Polar angle.
Definition: ParticleMC.h:392
virtual const ParticleMC * GetChild(UShort_t) const
Definition: ParticleMC.cxx:184
virtual Bool_t HasChild(Int_t) const
Definition: ParticleMC.cxx:217
virtual const ParticleMC * GetTrack(UInt_t) const
Definition: EventMC.h:225
virtual void ComputeEventDependentQuantities(EventMC &)
Definition: ParticleMC.cxx:132
virtual TLorentzVector Get4VectorInHadronBosonFrame() const
Definition: ParticleMC.cxx:229
Double32_t phiPrf
Definition: ParticleMC.h:403
UShort_t daughter
I of first child particle.
Definition: ParticleMC.h:355
virtual UInt_t GetNChildren() const
Definition: ParticleMC.h:500
Double32_t p
Total momentum of particle.
Definition: ParticleMC.h:391
virtual Pid Id() const
Definition: ParticleMC.h:496
Int_t parentId
PDG code of this particle's parent.
Definition: ParticleMC.h:389
virtual UInt_t GetIndex() const
Definition: ParticleMC.h:412
virtual Double_t GetNu() const
Definition: EventDis.h:202
UShort_t KS
Particle status code: see PYTHIA manual.
Definition: ParticleMC.h:352
Double32_t xFeynman
Feynman x = pz/(2sqrt(s))
Definition: ParticleMC.h:398
virtual Double_t GetW2() const
Definition: EventDis.h:210
virtual UShort_t GetParentIndex() const
Definition: ParticleMC.h:420