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
|
Go to the documentation of this file.
18 #include <TLorentzRotation.h>
19 #include <TParticlePDG.h>
20 #include <TRotation.h>
35 TLorentzRotation computeBoost(
const TLorentzVector& rest,
36 const TLorentzVector* z) {
37 TLorentzRotation toRest(-(rest.BoostVector()));
40 TLorentzVector boostedZ(*z);
42 rotate.SetZAxis(boostedZ.Vect());
45 rotate = rotate.Inverse();
46 toRest.Transform(rotate);
58 , id(std::numeric_limits<Int_t>::min())
71 , parentId(std::numeric_limits<Int_t>::min())
89 static std::stringstream ss;
113 if (ss.fail() || !ss.eof()) {
114 throw std::runtime_error(
"Bad particle input: " + line);
127 std::cout <<
I <<
'\t' <<
KS <<
'\t' <<
id <<
'\t' <<
orig <<
'\t' <<
129 <<
'\t' <<
E <<
'\t' <<
m <<
'\t' <<
xv <<
'\t' <<
yv <<
'\t' <<
zv <<
135 pt = sqrt(pow(
px, 2.) + pow(
py, 2.));
136 p = sqrt(pow(
pt, 2.) + pow(
pz, 2.));
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) {
148 rapidity = 0.5 * log(Epluspz / Eminuspz);
149 eta = 0.5 * log(Ppluspz / Pminuspz);
152 phi = TVector2::Phi_0_2pi(atan2(
py,
px));
158 const TLorentzVector& hadron =
event.BeamHadron()->Get4Vector();
159 const TLorentzVector& lepton =
event.ScatteredLepton()->Get4Vector();
160 const TLorentzVector& boson =
event.ExchangeBoson()->Get4Vector();
167 TLorentzRotation toHadronRest = computeBoost(hadron, &boson);
175 TLorentzVector bosonPrf = (TLorentzVector(boson) *= toHadronRest);
176 TLorentzVector leptonPrf = (TLorentzVector(lepton) *= toHadronRest);
181 TLorentzRotation boost = computeBoost(boson + hadron, &boson);
188 if (
event.GetNTracks() >
unsigned(
orig - 1)) {
192 catch(std::exception& e) {
194 "Exception in Particle::ComputeEventDependentQuantities: " <<
195 e.what() << std::endl;
200 return TLorentzVector(
px,
py,
pz,
E);
204 return static_cast<EventMC*
>(
event.GetObject());
222 if (idx <
GetEvent()->GetNTracks()) {
253 double p_(0.), e_(0.), px_(
ptVsGamma), py_(0.), pz_(0.);
262 e_ = sqrt(pow(p_, 2.) + pow(
m, 2.));
264 e_ = sqrt(pow(p_, 2.) - pow(
m, 2.));
265 if (TMath::IsNaN(e_)) {
295 p_ = sqrt(pow(e_, 2.) + pow(
m, 2.));
300 return TLorentzVector(px_, py_, pz_, e_);
363 static std::stringstream ss;
373 >> massNum >> charge >> NoBam;
376 if (ss.fail() || !ss.eof()) {
377 throw std::runtime_error(
"Bad particle input: " + line);
Int_t KS
Particle status code: see PYTHIA manual.
virtual Pid Id() const
Returns the ID of the particle.
virtual const ParticleMC * GetParent() const
Returns a pointer to the parent of this particle.
Int_t NoBam
0, 2 create in hard collisions not affected by INC: 0 inside nucleus, 2 outside nucleus;
Double32_t ptVsGamma
pt w.r.t.
Double_t yv
y coordinate of particle production vertex
Double32_t phiPrf
Azimuthal angle around virtual photon in hadron beam rest frame.
virtual TLorentzVector Get4Vector() const
Returns the (E,p) 4-vector in the lab frame.
virtual Double_t GetNu() const
Returns the exchange boson energy in the beam hadron rest frame.
ParticleMC()
Default constructor.
Double_t py
y component of particle momentum
virtual UInt_t GetNChildren() const
Returns the number of children of this particle.
Double_t pz
z component of particle momentum
A particle produced by a Monte Carlo generator.
Double32_t thetaGamma
Angle between particle and the exchange boson in the hadron beam rest frame.
Double32_t m
Invariant mass of particle.
Double32_t rapidity
Rapidity of particle.
UShort_t orig
I of parent particle.
UShort_t daughter
I of first child particle.
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.
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...
Int_t parentId
PDG code of this particle's parent.
UShort_t orig1
I of parent particle1.
virtual UShort_t GetParentIndex() const
Returns the index of this particle's parent in an event.
Double32_t p
Total momentum of particle.
virtual TLorentzVector Get4VectorInHadronBosonFrame() const
Returns the (E,p) 4-vector in the hadron-boson frame.
virtual UInt_t GetIndex() const
Returns the particle index in an event, in the range [1, N].
virtual void ComputeEventDependentQuantities(EventMC &)
Sets quantities that depend on the properties of the event or associations of one particle with anoth...
Double32_t E
Energy of particle.
virtual void Print(Option_t *="") const
Print the contents of Particle to standard output.
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,...
Double32_t xFeynman
Feynman x = pz/(2sqrt(s))
Double_t px
x component of particle momentum
virtual void SetVertex(const TVector3 &)
Sets the origin coordinates.
const EventMC * GetEvent() const
Returns a pointer to the event containing this particle.
Double_t zv
z coordinate of particle production vertex
UShort_t ldaughter
I of last child particle.
Double32_t pt
Transverse momentum of particle.
virtual void ComputeDerivedQuantities()
Sets quantities derived from the four-momentum (E, px, py, pz), namely.
Double32_t theta
Polar angle.
UShort_t I
Particle index in event.
virtual const ParticleMC * GetTrack(UInt_t) const
Returns the nth track.
Double_t xv
x coordinate of particle production vertex
Double32_t eta
Pseudorapidity of particle.
void SetEvent(EventMC *event)
Set the event with which to associate this particle.
virtual void Set4Vector(const TLorentzVector &)
Sets the four-momentum of the particle.
Double32_t z
Fraction of virtual photon energy carried by particle.
TRef event
Persistent reference to the event containing this particle.
Double32_t phi
Azimuthal angle.
Abstract base class for DIS Monte Carlo events.