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
erhic/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 
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  {
83  }
84 
85  ParticleMC::ParticleMC(const std::string &line, bool eAflag): ParticleMCbase(), eA(0)
86 {
87  // Initialise to nonsense values to make input errors easy to spot
88  if (!line.empty()) {
89  static std::stringstream ss;
90  ss.str("");
91  ss.clear();
92  ss << line;
93  if (eAflag) {
94  eA = new ParticleMCeA();
95 
96  ss >>
97  //changed by liang to add particle mother1
98  // I >> KS >> id >> orig >> daughter >> ldaughter >>
99  I >> KS >> id >> orig1 >> orig >> daughter >> ldaughter >>
100  px >> py >> pz >> E >> m >> xv >> yv >> zv
101  //added by liang to include add particle data structure
102  >> eA->massNum >> eA->charge >> eA->NoBam;
103 
104  //eA->pz = pz; orig1 = eA->orig1;
105  } else {
106  ss >>
107  I >> KS >> id >> orig >> daughter >> ldaughter >>
108  px >> py >> pz >> E >> m >> xv >> yv >> zv;
109  orig1 = 0;
110  } //if
111  // We should have no stream errors and should have exhausted
112  // the whole of the stream filling the particle.
113  if (ss.fail() || !ss.eof()) {
114  throw std::runtime_error("Bad particle input: " + line);
115  } // if
117  } // if
118 }
119 
121 {
122  if (eA) delete eA;
123 }
124 
125  // FIXME: may also want to print out ParticleMCeA variables?;
126 void ParticleMCbase::Print(Option_t* /* option */) const {
127  std::cout << I << '\t' << KS << '\t' << id << '\t' << orig << '\t' <<
128  daughter << '\t' << ldaughter << '\t' << px << '\t' << py << '\t' << pz
129  << '\t' << E << '\t' << m << '\t' << xv << '\t' << yv << '\t' << zv <<
130  std::endl;
131 }
132 
134  // Calculate quantities that depend only on the properties already read.
135  pt = sqrt(pow(px, 2.) + pow(py, 2.));
136  p = sqrt(pow(pt, 2.) + pow(pz, 2.));
137  // Rapidity and pseudorapidity
138  Double_t Epluspz = E + pz;
139  Double_t Eminuspz = E - pz;
140  Double_t Ppluspz = p + pz;
141  Double_t Pminuspz = p - pz;
142  if (Eminuspz <= 0.0 || Pminuspz == 0.0 ||
143  Ppluspz == 0.0 || Epluspz <= 0.0) {
144  // Dummy values to avoid zero or infinite arguments in calculations
145  rapidity = -19.;
146  eta = -19.;
147  } else {
148  rapidity = 0.5 * log(Epluspz / Eminuspz);
149  eta = 0.5 * log(Ppluspz / Pminuspz);
150  } // if
151  theta = atan2(pt, pz);
152  phi = TVector2::Phi_0_2pi(atan2(py, px));
153 }
154 
156  try {
157  // Get the beam hadon, beam lepton and exchange boson.
158  const TLorentzVector& hadron = event.BeamHadron()->Get4Vector();
159  const TLorentzVector& lepton = event.ScatteredLepton()->Get4Vector();
160  const TLorentzVector& boson = event.ExchangeBoson()->Get4Vector();
161  // Calculate z using the 4-vector definition,
162  // so we don't care about frame of reference.
163  z = hadron.Dot(Get4Vector()) / hadron.Dot(boson);
164  // Calculate properties in the proton rest frame.
165  // We want pT and angle with respect to the virtual photon,
166  // so use that to define the z axis.
167  TLorentzRotation toHadronRest = computeBoost(hadron, &boson);
168  // Boost this particle to the proton rest frame and calculate its
169  // pT and angle with respect to the virtual photon:
170  TLorentzVector p = (Get4Vector() *= toHadronRest);
171  thetaGamma = p.Theta();
172  ptVsGamma = p.Pt();
173  // Calculate phi angle around virtual photon according
174  // to the HERMES convention.
175  TLorentzVector bosonPrf = (TLorentzVector(boson) *= toHadronRest);
176  TLorentzVector leptonPrf = (TLorentzVector(lepton) *= toHadronRest);
177  phiPrf = computeHermesPhiH(p, leptonPrf, bosonPrf);
178  // Feynman x with xF = 2 * pz / W in the boson-hadron CM frame.
179  // First boost to boson-hadron centre-of-mass frame.
180  // Use the photon to define the z direction.
181  TLorentzRotation boost = computeBoost(boson + hadron, &boson);
182  xFeynman = 2. * (Get4Vector() *= boost).Pz() / sqrt(event.GetW2());
183  // Determine the PDG code of the parent particle, if the particle
184  // has a parent and the parent is present in the particle array.
185  // The index of the particles from the Monte Carlo runs from [1,N]
186  // while the index in the array runs from [0,N-1], so subtract 1
187  // from the parent index to find its position.
188  if (event.GetNTracks() > unsigned(orig - 1)) {
189  parentId = event.GetTrack(orig - 1)->Id();
190  } // if
191  } // try
192  catch(std::exception& e) {
193  std::cerr <<
194  "Exception in Particle::ComputeEventDependentQuantities: " <<
195  e.what() << std::endl;
196  } // catch
197 }
198 
199 TLorentzVector ParticleMCbase::Get4Vector() const {
200  return TLorentzVector(px, py, pz, E);
201 }
202 
204  return static_cast<EventMC*>(event.GetObject());
205 }
206 
207 const ParticleMC* ParticleMC::GetChild(UShort_t u) const {
208  // Look up this particle's child via the event
209  // containing it and the index of the child in that event.
210  if (!GetEvent()) {
211  return NULL;
212  } // if
213  // index is in the range [1,N]
214  unsigned idx = daughter + u;
215  if (daughter < 1 || // If first daughter index = 0, it has no children
216  u >= GetNChildren()) { // Insufficient children
217  return NULL;
218  } // if
219  --idx; // Convert [1,N] --> [0,N)
220  const ParticleMC* p(NULL);
221  // Check this index is within the # of particles in the event
222  if (idx < GetEvent()->GetNTracks()) {
223  p = GetEvent()->GetTrack(idx);
224  } // if
225  return p;
226 }
227 
229  // Look up this particle's parent via the event
230  // containing it and the index of the parent in that event.
231  const ParticleMC* p(NULL);
232  if (GetEvent()) {
233  if (GetEvent()->GetNTracks() >= GetParentIndex()) {
234  p = GetEvent()->GetTrack(GetParentIndex() - 1);
235  } // if
236  } // if
237  return p;
238 }
239 
240 Bool_t ParticleMC::HasChild(Int_t pdg) const {
241  for (UInt_t i(0); i < GetNChildren(); ++i) {
242  if (!GetChild(i)) {
243  continue;
244  } // if
245  if (pdg == GetChild(i)->Id()) {
246  return true;
247  } // if
248  } // for
249  return false;
250 }
251 
253  double p_(0.), e_(0.), px_(ptVsGamma), py_(0.), pz_(0.);
254  // Calculate mangitude of momentum from pT and polar angle in
255  // hadron-boson frame. If theta is ~parallel to the beam just set
256  // p to whatever pT is (not that it makes all that much sense).
257  if (thetaGamma > 1.e-6) {
258  p_ = ptVsGamma / sin(thetaGamma);
259  } // if
260  // Deal with virtual particles later, so check if particle is off mass-shell
261  if (!(m < 0.)) {
262  e_ = sqrt(pow(p_, 2.) + pow(m, 2.));
263  } else {
264  e_ = sqrt(pow(p_, 2.) - pow(m, 2.));
265  if (TMath::IsNaN(e_)) {
266  e_ = 0.;
267  } // if
268  } // if
269  // Calculate pZ from pT and theta, unless it's very close to the beam,
270  // in which case set it to p.
271  if (thetaGamma > 1.e-6) {
272  pz_ = ptVsGamma / tan(thetaGamma);
273  } // if
274  // Calculate px and py, unless a dummy phi value is present.
275  if (phiPrf > -100.) {
276  px_ = ptVsGamma * cos(phiPrf);
277  py_ = ptVsGamma * sin(phiPrf);
278  } // if
279  // If we ended up with no energy, it's likely ths is the exchange boson,
280  // as nothing will have happened above. If so, try to reference the event
281  // record to get the necessary information, as we don't have enough
282  // in the particle itself.
283  // Note that this appears not to work in a TTree as ROOT doesn't read
284  // the event when it reads the particle, though I'm not certain of the
285  // exact cause.
286  if (m < 0. && GetEvent()) {
287  if (GetEvent()->ExchangeBoson()) {
288  // Don't check if GetEvent()->ExchangeBoson() == this, in case this
289  // particle is a copy of the event track that isn't part of a copied
290  // event.
291  if (GetEvent()->ExchangeBoson()->GetIndex() == GetIndex()) {
292  // If it's the exchange boson just set E == nu.
293  e_ = GetEvent()->GetNu();
294  // Calculate p, careful about negative mass.
295  p_ = sqrt(pow(e_, 2.) + pow(m, 2.));
296  pz_ = p_ * cos(thetaGamma);
297  } // if
298  } // if
299  } // if
300  return TLorentzVector(px_, py_, pz_, e_);
301 }
302 
304  event = e;
305 }
306 
307 void ParticleMCbase::Set4Vector(const TLorentzVector& v) {
308  E = v.Energy();
309  px = v.Px();
310  py = v.Py();
311  pz = v.Pz();
312  m = v.M();
313  ComputeDerivedQuantities(); // Rapidity etc
314  // If an event reference is set, recalculate event-dependent quantities
315  // like z, xF
316  EventMC* ev = static_cast<EventMC*>(event.GetObject());
317  if (ev) {
319  } // if
320 }
321 
322 void ParticleMCbase::SetVertex(const TVector3& v) {
323  xv = v.X();
324  yv = v.Y();
325  zv = v.Z();
326 }
327 
328 #if _OLD_
329 ParticleMCeA::ParticleMCeA(const std::string& line)
330 : /*I(-1)
331 , KS(-1)
332 , id(std::numeric_limits<Int_t>::min())
333 , orig(-1)
334 , daughter(-1)
335 , ldaughter(-1)
336 , px(0.)
337 , py(0.)
338 , pz(0.)
339 , E(0.)
340 , m(0.)
341 , pt(0.)
342 , xv(0.)
343 , yv(0.)
344 , zv(0.)
345 , parentId(std::numeric_limits<Int_t>::min())
346 , p(0.)
347 , theta(0.)
348 , phi(0.)
349 , rapidity(0.)
350 , eta(0.)
351 , z(0.)
352 , xFeynman(0.)
353 , thetaGamma(0.)
354 , ptVsGamma(0.)
355 , phiPrf(0.)
356 //added by liang to include add particle data structure
357 ,*/ orig1(-1)
358 , charge(-999)
359 , massNum(-999)
360 , NoBam(-999) {
361  // Initialise to nonsense values to make input errors easy to spot
362  if (!line.empty()) {
363  static std::stringstream ss;
364  ss.str("");
365  ss.clear();
366  ss << line;
367  ss >>
368  //changed by liang to add particle mother1
369 // I >> KS >> id >> orig >> daughter >> ldaughter >>
370  I >> KS >> id >> orig1 >> orig >> daughter >> ldaughter >>
371  px >> py >> pz >> E >> m >> xv >> yv >> zv
372  //added by liang to include add particle data structure
373  >> massNum >> charge >> NoBam;
374  // We should have no stream errors and should have exhausted
375  // the whole of the stream filling the particle.
376  if (ss.fail() || !ss.eof()) {
377  throw std::runtime_error("Bad particle input: " + line);
378  } // if
380  } // if
381 }
382 #endif
383 
384 } // namespace erhic
erhic::ParticleMC::~ParticleMC
~ParticleMC()
Definition: erhic/ParticleMC.cxx:120
erhic::ParticleMCbase::KS
Int_t KS
Particle status code: see PYTHIA manual.
Definition: erhic/ParticleMC.h:334
erhic::ParticleMCeA
Definition: erhic/ParticleMC.h:24
erhic::ParticleMCbase::Id
virtual Pid Id() const
Returns the ID of the particle.
Definition: erhic/ParticleMC.h:544
erhic::ParticleMC::GetParent
virtual const ParticleMC * GetParent() const
Returns a pointer to the parent of this particle.
Definition: erhic/ParticleMC.cxx:228
erhic::ParticleMCeA::NoBam
Int_t NoBam
0, 2 create in hard collisions not affected by INC: 0 inside nucleus, 2 outside nucleus;
Definition: erhic/ParticleMC.h:33
erhic::ParticleMCbase::ptVsGamma
Double32_t ptVsGamma
pt w.r.t.
Definition: erhic/ParticleMC.h:388
erhic
Definition: EventDis.cxx:14
erhic::ParticleMCbase::yv
Double_t yv
y coordinate of particle production vertex
Definition: erhic/ParticleMC.h:372
erhic::ParticleMCbase::phiPrf
Double32_t phiPrf
Azimuthal angle around virtual photon in hadron beam rest frame.
Definition: erhic/ParticleMC.h:390
erhic::ParticleMCbase::Get4Vector
virtual TLorentzVector Get4Vector() const
Returns the (E,p) 4-vector in the lab frame.
Definition: erhic/ParticleMC.cxx:199
erhic::EventDis::GetNu
virtual Double_t GetNu() const
Returns the exchange boson energy in the beam hadron rest frame.
Definition: EventDis.h:202
erhic::ParticleMC::ParticleMC
ParticleMC()
Default constructor.
Definition: erhic/ParticleMC.h:411
erhic::ParticleMCbase::py
Double_t py
y component of particle momentum
Definition: erhic/ParticleMC.h:346
erhic::ParticleMCbase::GetNChildren
virtual UInt_t GetNChildren() const
Returns the number of children of this particle.
Definition: erhic/ParticleMC.h:548
erhic::ParticleMCbase::pz
Double_t pz
z component of particle momentum
Definition: erhic/ParticleMC.h:347
erhic::ParticleMCbase::ParticleMCbase
ParticleMCbase()
Definition: erhic/ParticleMC.cxx:55
EventBase.h
erhic::ParticleMCbase
A particle produced by a Monte Carlo generator.
Definition: erhic/ParticleMC.h:43
erhic::ParticleMCbase::thetaGamma
Double32_t thetaGamma
Angle between particle and the exchange boson in the hadron beam rest frame.
Definition: erhic/ParticleMC.h:386
erhic::ParticleMCeA::ParticleMCeA
ParticleMCeA()
Definition: erhic/ParticleMC.h:26
functions.h
erhic::ParticleMCbase::m
Double32_t m
Invariant mass of particle.
Definition: erhic/ParticleMC.h:369
erhic::ParticleMCbase::rapidity
Double32_t rapidity
Rapidity of particle.
Definition: erhic/ParticleMC.h:381
erhic::ParticleMCbase::orig
UShort_t orig
I of parent particle.
Definition: erhic/ParticleMC.h:338
erhic::ParticleMCbase::daughter
UShort_t daughter
I of first child particle.
Definition: erhic/ParticleMC.h:340
erhic::ParticleMC::HasChild
virtual Bool_t HasChild(Int_t) const
Returns true if n in the range [0, N), where N is the number of children of this particle.
Definition: erhic/ParticleMC.cxx:240
computeHermesPhiH
double computeHermesPhiH(const TLorentzVector &hadronInPrf, const TLorentzVector &leptonInPrf, const TLorentzVector &photonInPrf)
Calculate the hadron azimuthal angle around the virtual photon direction with respect to the lepton s...
Definition: functions.cxx:52
erhic::ParticleMCbase::parentId
Int_t parentId
PDG code of this particle's parent.
Definition: erhic/ParticleMC.h:376
erhic::ParticleMCbase::orig1
UShort_t orig1
I of parent particle1.
Definition: erhic/ParticleMC.h:339
erhic::ParticleMCbase::GetParentIndex
virtual UShort_t GetParentIndex() const
Returns the index of this particle's parent in an event.
Definition: erhic/ParticleMC.h:465
erhic::ParticleMCbase::p
Double32_t p
Total momentum of particle.
Definition: erhic/ParticleMC.h:378
erhic::ParticleMCbase::Get4VectorInHadronBosonFrame
virtual TLorentzVector Get4VectorInHadronBosonFrame() const
Returns the (E,p) 4-vector in the hadron-boson frame.
Definition: erhic/ParticleMC.cxx:252
erhic::ParticleMC::eA
ParticleMCeA * eA
Definition: erhic/ParticleMC.h:451
erhic::ParticleMCbase::GetIndex
virtual UInt_t GetIndex() const
Returns the particle index in an event, in the range [1, N].
Definition: erhic/ParticleMC.h:457
erhic::ParticleMCbase::ComputeEventDependentQuantities
virtual void ComputeEventDependentQuantities(EventMC &)
Sets quantities that depend on the properties of the event or associations of one particle with anoth...
Definition: erhic/ParticleMC.cxx:155
erhic::ParticleMCbase::E
Double32_t E
Energy of particle.
Definition: erhic/ParticleMC.h:368
erhic::ParticleMCbase::Print
virtual void Print(Option_t *="") const
Print the contents of Particle to standard output.
Definition: erhic/ParticleMC.cxx:126
erhic::ParticleMC::GetChild
virtual const ParticleMC * GetChild(UShort_t) const
Returns a pointer to the nth child particle of this particle, where n is in the range [0,...
Definition: erhic/ParticleMC.cxx:207
erhic::ParticleMCbase::xFeynman
Double32_t xFeynman
Feynman x = pz/(2sqrt(s))
Definition: erhic/ParticleMC.h:385
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::ParticleMCeA::massNum
Int_t massNum
Definition: erhic/ParticleMC.h:32
erhic::ParticleMCbase::px
Double_t px
x component of particle momentum
Definition: erhic/ParticleMC.h:345
erhic::ParticleMCbase::SetVertex
virtual void SetVertex(const TVector3 &)
Sets the origin coordinates.
Definition: erhic/ParticleMC.cxx:322
erhic::ParticleMCbase::GetEvent
const EventMC * GetEvent() const
Returns a pointer to the event containing this particle.
Definition: erhic/ParticleMC.cxx:203
erhic::ParticleMCbase::zv
Double_t zv
z coordinate of particle production vertex
Definition: erhic/ParticleMC.h:373
erhic::ParticleMCbase::ldaughter
UShort_t ldaughter
I of last child particle.
Definition: erhic/ParticleMC.h:341
erhic::ParticleMCbase::pt
Double32_t pt
Transverse momentum of particle.
Definition: erhic/ParticleMC.h:370
erhic::ParticleMCbase::ComputeDerivedQuantities
virtual void ComputeDerivedQuantities()
Sets quantities derived from the four-momentum (E, px, py, pz), namely.
Definition: erhic/ParticleMC.cxx:133
erhic::ParticleMCbase::theta
Double32_t theta
Polar angle.
Definition: erhic/ParticleMC.h:379
erhic::ParticleMCbase::I
UShort_t I
Particle index in event.
Definition: erhic/ParticleMC.h:331
erhic::ParticleMCeA::charge
Int_t charge
Definition: erhic/ParticleMC.h:27
erhic::EventMC::GetTrack
virtual const ParticleMC * GetTrack(UInt_t) const
Returns the nth track.
Definition: erhic/EventMC.h:227
erhic::ParticleMCbase::xv
Double_t xv
x coordinate of particle production vertex
Definition: erhic/ParticleMC.h:371
ParticleMC.h
erhic::ParticleMCbase::eta
Double32_t eta
Pseudorapidity of particle.
Definition: erhic/ParticleMC.h:382
erhic::ParticleMCbase::SetEvent
void SetEvent(EventMC *event)
Set the event with which to associate this particle.
Definition: erhic/ParticleMC.cxx:303
erhic::ParticleMCbase::Set4Vector
virtual void Set4Vector(const TLorentzVector &)
Sets the four-momentum of the particle.
Definition: erhic/ParticleMC.cxx:307
erhic::ParticleMCbase::z
Double32_t z
Fraction of virtual photon energy carried by particle.
Definition: erhic/ParticleMC.h:383
erhic::ParticleMCbase::event
TRef event
Persistent reference to the event containing this particle.
Definition: erhic/ParticleMC.h:393
erhic::ParticleMCbase::phi
Double32_t phi
Azimuthal angle.
Definition: erhic/ParticleMC.h:380
erhic::EventMC
Abstract base class for DIS Monte Carlo events.
Definition: erhic/EventMC.h:30