EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ParticleMC.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file ParticleMC.cxx
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(0)
57 , KS(0)
58 , id(0)
59 , orig(0)
60 , orig1(0)
61 , daughter(0)
62 , ldaughter(0)
63 , px(0.)
64 , py(0.)
65 , pz(0.)
66 , E(0.)
67 , m(0.)
68 , pt(0.)
69 , xv(0.)
70 , yv(0.)
71 , zv(0.)
72 , parentId(std::numeric_limits<Int_t>::min())
73 , p(0.)
74 , theta(0.)
75 , phi(0.)
76 , rapidity(0.)
77 , eta(0.)
78 , z(0.)
79 , xFeynman(0.)
80 , thetaGamma(0.)
81 , ptVsGamma(0.)
82 , thetaGammaHCM(0.)
83 , ptVsGammaHCM(0.)
84 , phiPrf(0.)
85  {
86  }
87 
88  ParticleMC::ParticleMC(const std::string &line, bool eAflag): ParticleMCbase(), eA(0)
89 {
90  // Initialise to nonsense values to make input errors easy to spot
91  if (!line.empty()) {
92  static std::stringstream ss;
93  ss.str("");
94  ss.clear();
95  ss << line;
96  if (eAflag) {
97  eA = new ParticleMCeA();
98 
99  ss >>
100  //changed by liang to add particle mother1
101  I >> KS >> id >> orig1 >> orig >> daughter >> ldaughter >>
102  px >> py >> pz >> E >> m >> xv >> yv >> zv
103  //added by liang to include add particle data structure
104  >> eA->massNum >> eA->charge >> eA->NoBam;
105 
106  //eA->pz = pz; orig1 = eA->orig1;
107  } else {
108  ss >>
109  I >> KS >> id >> orig >> daughter >> ldaughter >>
110  px >> py >> pz >> E >> m >> xv >> yv >> zv;
111  orig1 = 0;
112  } //if
113  // We should have no stream errors and should have exhausted
114  // the whole of the stream filling the particle.
115  if (ss.fail() || !ss.eof()) {
116  throw std::runtime_error("Bad particle input: " + line);
117  } // if
119  } // if
120 }
121 
123 {
124  if (eA) delete eA;
125 }
126 
127  // FIXME: may also want to print out ParticleMCeA variables?;
128 void ParticleMCbase::Print(Option_t* /* option */) const {
129  std::cout << I << '\t' << KS << '\t' << id << '\t' << orig << '\t' <<
130  daughter << '\t' << ldaughter << '\t' << px << '\t' << py << '\t' << pz
131  << '\t' << E << '\t' << m << '\t' << xv << '\t' << yv << '\t' << zv <<
132  std::endl;
133 }
134 
136  // Calculate quantities that depend only on the properties already read.
137  pt = sqrt(pow(px, 2.) + pow(py, 2.));
138  p = sqrt(pow(pt, 2.) + pow(pz, 2.));
139  // Rapidity and pseudorapidity
140  Double_t Epluspz = E + pz;
141  Double_t Eminuspz = E - pz;
142  Double_t Ppluspz = p + pz;
143  Double_t Pminuspz = p - pz;
144  if (Eminuspz <= 0.0 || Pminuspz == 0.0 ||
145  Ppluspz == 0.0 || Epluspz <= 0.0) {
146  // Dummy values to avoid zero or infinite arguments in calculations
147  rapidity = -19.;
148  eta = -19.;
149  } else {
150  rapidity = 0.5 * log(Epluspz / Eminuspz);
151  eta = 0.5 * log(Ppluspz / Pminuspz);
152  } // if
153  theta = atan2(pt, pz);
154  phi = TVector2::Phi_0_2pi(atan2(py, px));
155 }
156 
158  try {
159  // Get the beam hadon, beam lepton and exchange boson.
160  const TLorentzVector& hadron = event.BeamHadron()->Get4Vector();
161  const TLorentzVector& lepton = event.ScatteredLepton()->Get4Vector();
162  // Determine the exchange boson 4-vector from the scattered lepton,
163  // since we're not always guaranteed to have one (weak neutral current e.g.)
164  // const TLorentzVector& boson = event.ExchangeBoson()->Get4Vector();
165  const TLorentzVector boson = event.BeamLepton()->Get4Vector() - lepton;
166 
167  // Calculate z using the 4-vector definition,
168  // so we don't care about frame of reference.
169  z = hadron.Dot(Get4Vector()) / hadron.Dot(boson);
170  // Calculate properties in the proton rest frame.
171  // We want pT and angle with respect to the virtual photon,
172  // so use that to define the z axis.
173  TLorentzRotation toHadronRest = computeBoost(hadron, &boson);
174  // Boost this particle to the proton rest frame and calculate its
175  // pT and angle with respect to the virtual photon:
176  TLorentzVector p = (Get4Vector() *= toHadronRest);
177  thetaGamma = p.Theta();
178  ptVsGamma = p.Pt();
179  // Calculate phi angle around virtual photon according
180  // to the HERMES convention.
181  TLorentzVector bosonPrf = (TLorentzVector(boson) *= toHadronRest);
182  TLorentzVector leptonPrf = (TLorentzVector(lepton) *= toHadronRest);
183  phiPrf = computeHermesPhiH(p, leptonPrf, bosonPrf);
184  // Feynman x with xF = 2 * pz / W in the boson-hadron CM frame.
185  // First boost to boson-hadron centre-of-mass frame.
186  // Use the photon to define the z direction.
187  TLorentzRotation boost = computeBoost(boson + hadron, &boson);
188  xFeynman = 2. * (Get4Vector() *= boost).Pz() / sqrt(event.GetW2());
189 
190  thetaGammaHCM = (Get4Vector() *= boost).Theta();
191  ptVsGammaHCM = (Get4Vector() *= boost).Pt();
192 
193  // Determine the PDG code of the parent particle, if the particle
194  // has a parent and the parent is present in the particle array.
195  // The index of the particles from the Monte Carlo runs from [1,N]
196  // while the index in the array runs from [0,N-1], so subtract 1
197  // from the parent index to find its position.
198  if (event.GetNTracks() > unsigned(orig - 1)) {
199  parentId = event.GetTrack(orig - 1)->Id();
200  } // if
201  } // try
202  catch(std::exception& e) {
203  std::cerr <<
204  "Exception in Particle::ComputeEventDependentQuantities: " <<
205  e.what() << std::endl;
206  } // catch
207 }
208 
209 TLorentzVector ParticleMCbase::Get4Vector() const {
210  return TLorentzVector(px, py, pz, E);
211 }
212 
214  return static_cast<EventMC*>(event.GetObject());
215 }
216 
217 const ParticleMC* ParticleMC::GetChild(UShort_t u) const {
218  // Look up this particle's child via the event
219  // containing it and the index of the child in that event.
220  if (!GetEvent()) {
221  return NULL;
222  } // if
223  // index is in the range [1,N]
224  unsigned idx = daughter + u;
225  if (daughter < 1 || // If first daughter index = 0, it has no children
226  u >= GetNChildren()) { // Insufficient children
227  return NULL;
228  } // if
229  --idx; // Convert [1,N] --> [0,N)
230  const ParticleMC* p(NULL);
231  // Check this index is within the # of particles in the event
232  if (idx < GetEvent()->GetNTracks()) {
233  p = GetEvent()->GetTrack(idx);
234  } // if
235  return p;
236 }
237 
239  // Look up this particle's parent via the event
240  // containing it and the index of the parent in that event.
241  const ParticleMC* p(NULL);
242  if (GetEvent()) {
243  if (GetEvent()->GetNTracks() >= GetParentIndex()) {
244  p = GetEvent()->GetTrack(GetParentIndex() - 1);
245  } // if
246  } // if
247  return p;
248 }
249 
250 Bool_t ParticleMC::HasChild(Int_t pdg) const {
251  for (UInt_t i(0); i < GetNChildren(); ++i) {
252  if (!GetChild(i)) {
253  continue;
254  } // if
255  if (pdg == GetChild(i)->Id()) {
256  return true;
257  } // if
258  } // for
259  return false;
260 }
261 
263  double p_(0.), e_(0.), px_(ptVsGammaHCM), py_(0.), pz_(0.);
264  // Calculate mangitude of momentum from pT and polar angle in
265  // hadron-boson frame. If theta is ~parallel to the beam just set
266  // p to whatever pT is (not that it makes all that much sense).
267  if (thetaGammaHCM > 1.e-6) {
268  p_ = ptVsGammaHCM / sin(thetaGammaHCM);
269  } // if
270  // Deal with virtual particles later, so check if particle is off mass-shell
271  if (!(m < 0.)) {
272  e_ = sqrt(pow(p_, 2.) + pow(m, 2.));
273  } else {
274  e_ = sqrt(pow(p_, 2.) - pow(m, 2.));
275  if (TMath::IsNaN(e_)) {
276  e_ = 0.;
277  } // if
278  } // if
279  // Calculate pZ from pT and theta, unless it's very close to the beam,
280  // in which case set it to p.
281  if (thetaGammaHCM > 1.e-6) {
282  pz_ = ptVsGammaHCM / tan(thetaGammaHCM);
283  } // if
284  // Calculate px and py, unless a dummy phi value is present.
285  if (phiPrf > -100.) {
286  px_ = ptVsGammaHCM * cos(phiPrf);
287  py_ = ptVsGammaHCM * sin(phiPrf);
288  } // if
289  // If we ended up with no energy, it's likely ths is the exchange boson,
290  // as nothing will have happened above. If so, try to reference the event
291  // record to get the necessary information, as we don't have enough
292  // in the particle itself.
293  // Note that this appears not to work in a TTree as ROOT doesn't read
294  // the event when it reads the particle, though I'm not certain of the
295  // exact cause.
296  if (m < 0. && GetEvent()) {
297  if (GetEvent()->ExchangeBoson()) {
298  // Don't check if GetEvent()->ExchangeBoson() == this, in case this
299  // particle is a copy of the event track that isn't part of a copied
300  // event.
301  if (GetEvent()->ExchangeBoson()->GetIndex() == GetIndex()) {
302  // If it's the exchange boson just set E == nu.
303  e_ = GetEvent()->GetNu();
304  // Calculate p, careful about negative mass.
305  p_ = sqrt(pow(e_, 2.) + pow(m, 2.));
306  pz_ = p_ * cos(thetaGammaHCM);
307  } // if
308  } // if
309  } // if
310  return TLorentzVector(px_, py_, pz_, e_);
311 }
312 
314  event = e;
315 }
316 
317 void ParticleMCbase::Set4Vector(const TLorentzVector& v) {
318  E = v.Energy();
319  px = v.Px();
320  py = v.Py();
321  pz = v.Pz();
322  m = v.M();
323  ComputeDerivedQuantities(); // Rapidity etc
324  // If an event reference is set, recalculate event-dependent quantities
325  // like z, xF
326  EventMC* ev = static_cast<EventMC*>(event.GetObject());
327  if (ev) {
329  } // if
330 }
331 
332 void ParticleMCbase::SetVertex(const TVector3& v) {
333  xv = v.X();
334  yv = v.Y();
335  zv = v.Z();
336 }
337 } // namespace erhic