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/EventPythia.cxx
Go to the documentation of this file.
1 
11 
12 #include <cmath>
13 #include <limits>
14 #include <sstream>
15 #include <string>
16 #include <vector>
17 
18 namespace erhic {
19 
20 EventPythia::EventPythia(const std::string& /* Unused */)
21 : nucleon(std::numeric_limits<Int_t>::max())
22 , tgtparton(std::numeric_limits<Int_t>::max())
23 , beamparton(std::numeric_limits<Int_t>::max())
24 , genevent(-1)
25 , xtgtparton(NAN)
26 , xbeamparton(NAN)
27 , thetabeamparton(NAN)
28 , leptonphi(NAN)
29 , F1(NAN)
30 , sigma_rad(NAN)
31 , t_hat(NAN)
32 , u_hat(NAN)
33 , Q2_hat(NAN)
34 , SigRadCor(NAN)
35 , EBrems(NAN)
36 , photonflux(NAN)
37 , trueY(NAN)
38 , trueQ2(NAN)
39 , trueX(NAN)
40 , trueW2(NAN)
41 , trueNu(NAN)
42 , F2(NAN)
43 , R(NAN)
44 , pt2_hat(NAN)
45 , sHat(NAN) {
46 }
47 
49 
50 bool EventPythia::Parse(const std::string& line) {
51  static std::stringstream ss;
52  ss.str("");
53  ss.clear();
54  ss << line;
55  ss >>
56  number >> number >> // Skip first int in the line
59  trueX >> trueW2 >> trueNu >> leptonphi >> sHat >> t_hat >> u_hat >>
60  pt2_hat >> Q2_hat >> F2 >> F1 >> R >> sigma_rad >> SigRadCor >> EBrems >>
62  // Protect against errors in the input file or the stream
63  return !ss.fail();
64 }
65 
66 // Look for the scattered lepton in the event record.
67 // This is the first (only?) particle that matches the following:
68 // 1) pdg code equals that of incident lepton beam.
69 // 2) status code is 1 i.e. it's a stable/final-state particle.
70 // 3) the parent is track three (counting from 1).
72  // Look for the lepton beam to get the species.
73  // If we don't get it we can't find the scattered
74  // lepton so return NULL.
75  const VirtualParticle* beam = BeamLepton();
76  if (!beam) {
77  return NULL;
78  } // if
79  const int species = beam->Id().Code();
80  // Get the final state particles and search them for
81  // the scattered lepton.
82  std::vector<const VirtualParticle*> final;
83  FinalState(final);
84  std::vector<const VirtualParticle*>::const_iterator iter;
85  for (iter = final.begin(); iter != final.end(); ++iter) {
86  // We already know the particle is final state, so
87  // check its species and parent index.
88  if ((*iter)->GetParentIndex() == 3 &&
89  (*iter)->Id().Code() == species) {
90  // Found it, cast to required particle type and return.
91  return static_cast<const ParticleMC*>(*iter);
92  } // if
93  } // for
94  // No luck, couldn't find the scattered lepton.
95  return NULL;
96 }
97 
98  EventBeagle::EventBeagle(const std::string& str/* Unused */):
99  EventPythia(str)
100  , lepton(std::numeric_limits<Int_t>::max())
101  , Atarg(std::numeric_limits<Int_t>::max())
102  , Ztarg(std::numeric_limits<Int_t>::max())
103  , pzlep(NAN)
104  , pztarg(NAN)
105  , pznucl(NAN)
106  , crang(NAN)
107  , crori(NAN)
108  , b(NAN)
109  , Phib(NAN)
110  , Thickness(NAN)
111  , ThickScl(NAN)
112  , Ncollt(std::numeric_limits<Int_t>::max())
113  , Ncolli(std::numeric_limits<Int_t>::max())
114  , Nwound(std::numeric_limits<Int_t>::max())
115  , Nwdch(std::numeric_limits<Int_t>::max())
116  , Nnevap(std::numeric_limits<Int_t>::max())
117  , Npevap(std::numeric_limits<Int_t>::max())
118  , Aremn(std::numeric_limits<Int_t>::max())
119  , NINC(std::numeric_limits<Int_t>::max())
120  , NINCch(std::numeric_limits<Int_t>::max())
121  , d1st(NAN)
122  , davg(NAN)
123  , pxf(NAN)
124  , pyf(NAN)
125  , pzf(NAN)
126  , Eexc(NAN)
127  , RAevt(NAN)
128  , User1(NAN)
129  , User2(NAN)
130  , User3(NAN)
131  {
132 
133  } // EventBeagle::EventBeagle()
134 
136 
137 bool EventBeagle::Parse(const std::string& line) {
138  static std::stringstream ss;
139  ss.str("");
140  ss.clear();
141  ss << line;
142  ss >>
143 // number >> number >> // Skip first int in the line
144 // genevent >> process >> nucleon >> tgtparton >> xtgtparton >>
145 // beamparton >> xbeamparton >> thetabeamparton >> trueY >> trueQ2 >>
146 // trueX >> trueW2 >> trueNu >> leptonphi >> sHat >> t_hat >> u_hat >>
147 // pt2_hat >> Q2_hat >> F2 >> F1 >> R >> sigma_rad >> SigRadCor >> EBrems >>
148 // photonflux >> nTracks;
149  number >> number >> // Skip first int in the line
150  genevent >>
151  lepton >> Atarg >> Ztarg >> pzlep >> pztarg >> pznucl >> crang >> crori >>//added variables
152  process >> nucleon >> tgtparton >> xtgtparton >>
154  trueX >> trueW2 >> trueNu >> leptonphi >> sHat >> t_hat >> u_hat >>
155  pt2_hat >> Q2_hat >> F2 >> F1 >> R >> sigma_rad >> SigRadCor >> EBrems >>
156  photonflux >>
157  b >> Phib >> Thickness >> ThickScl >> Ncollt >> Ncolli >> Nwound >> Nwdch >> Nnevap >> Npevap >> Aremn >> //added variables
158  NINC >> NINCch >> d1st >> davg >> pxf >> pyf >> pzf >> Eexc >> RAevt >> User1 >> User2 >> User3 >> //added variables
159  nTracks;
160  // Protect against errors in the input file or the stream
161  return !ss.fail();
162 }
163 } // namespace erhic
erhic::EventPythia::t_hat
Double32_t t_hat
Mandelstam t of the hard subprocess, see PARI(15)
Definition: erhic/EventPythia.h:331
erhic::EventPythia::xtgtparton
Double32_t xtgtparton
Momentum fraction taken by the target parton, see PARI(34)
Definition: erhic/EventPythia.h:320
erhic::EventBeagle::Ncollt
Int_t Ncollt
Number of collisions in target.
Definition: erhic/EventPythia.h:570
erhic::EventBeagle::RAevt
Double32_t RAevt
Nuclear PDF ratio for the up sea for the given event kinematics.
Definition: erhic/EventPythia.h:585
erhic::EventPythia::sigma_rad
Double32_t sigma_rad
Value used for radiative corrections.
Definition: erhic/EventPythia.h:330
erhic::EventBeagle::lepton
Int_t lepton
=================additional variables for BeAGLE==================
Definition: erhic/EventPythia.h:558
erhic::EventPythia::EventPythia
EventPythia(const std::string &str="")
Constructor.
Definition: erhic/EventPythia.cxx:20
erhic
Definition: EventDis.cxx:14
erhic::EventBeagle::Npevap
Int_t Npevap
Number of protons from evaporation.
Definition: erhic/EventPythia.h:575
erhic::EventPythia::sHat
Double32_t sHat
Mandelstam s of the hard subprocess, see PARI(14)
Definition: erhic/EventPythia.h:352
EventPythia.h
erhic::EventBeagle::pyf
Double32_t pyf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:582
erhic::EventPythia::nucleon
Int_t nucleon
PDG code of the hadron beam, see MSTI(12)
Definition: erhic/EventPythia.h:312
erhic::EventPythia::u_hat
Double32_t u_hat
Mandelstam u of the hard subprocess, see PARI(16)
Definition: erhic/EventPythia.h:333
erhic::EventBeagle::Ztarg
Int_t Ztarg
charge number of target beam
Definition: erhic/EventPythia.h:560
erhic::EventPythia::Q2_hat
Double32_t Q2_hat
Q2 of the hard subprocess, see PARI(22)
Definition: erhic/EventPythia.h:335
erhic::EventMC::process
Int_t process
PYTHIA code for the physics process producing the event.
Definition: erhic/EventMC.h:204
erhic::EventBeagle::crang
Double32_t crang
crossing angle (mr). crang=1000*atan(px/pz), one beam px=py=0, the other py=0
Definition: erhic/EventPythia.h:564
erhic::EventBeagle::NINC
Int_t NINC
Number of stable hadrons from intranuclear cascade.
Definition: erhic/EventPythia.h:577
erhic::EventPythia::~EventPythia
virtual ~EventPythia()
Destructor.
Definition: erhic/EventPythia.cxx:48
erhic::EventBeagle::Atarg
Int_t Atarg
mass number of target beam
Definition: erhic/EventPythia.h:559
erhic::EventBeagle::d1st
Double32_t d1st
density-weighted distance from first collision to the edge of the nucleus (amount of material travers...
Definition: erhic/EventPythia.h:579
erhic::EventPythia::trueQ2
Double32_t trueQ2
Generated Q2 of the event, see VINT(307)
Definition: erhic/EventPythia.h:343
erhic::EventBeagle::User1
Double32_t User1
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:586
erhic::EventBeagle::~EventBeagle
~EventBeagle()
Definition: erhic/EventPythia.cxx:135
erhic::EventBeagle::NINCch
Int_t NINCch
Number of charged stable hadrons from intranuclear cascade.
Definition: erhic/EventPythia.h:578
erhic::EventPythia::R
Double32_t R
Value used for radiative corrections.
Definition: erhic/EventPythia.h:349
erhic::EventPythia::thetabeamparton
Double32_t thetabeamparton
Polar angle of the beam parton in the cm frame, between 0 and pi radians, see PARI(53)
Definition: erhic/EventPythia.h:324
erhic::EventPythia::ScatteredLepton
virtual const ParticleMC * ScatteredLepton() const
Returns the scattered lepton.
Definition: erhic/EventPythia.cxx:71
erhic::EventPythia::trueX
Double32_t trueX
Generated x of the event.
Definition: erhic/EventPythia.h:345
erhic::EventPythia
Pythia 6 DIS event.
Definition: erhic/EventPythia.h:28
erhic::Pid::Code
Int_t Code() const
Returns the PDG code, an integer.
Definition: Pid.h:105
erhic::EventPythia::SigRadCor
Double32_t SigRadCor
Value used for radiative corrections.
Definition: erhic/EventPythia.h:337
erhic::EventBeagle::Nwdch
Int_t Nwdch
Number of wounded proton including those in INC.
Definition: erhic/EventPythia.h:573
erhic::EventMC::FinalState
void FinalState(ParticlePtrList &particles) const
Stores pointers to all final state particles in the list.
Definition: erhic/EventMC.cxx:73
erhic::EventBeagle::User2
Double32_t User2
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:587
erhic::EventBeagle::b
Double32_t b
impact parameter value
Definition: erhic/EventPythia.h:566
erhic::EventBeagle::pxf
Double32_t pxf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:581
erhic::EventPythia::genevent
Int_t genevent
Trials required for this event.
Definition: erhic/EventPythia.h:319
erhic::EventBeagle::Nnevap
Int_t Nnevap
Number of neutrons from evaporation.
Definition: erhic/EventPythia.h:574
erhic::VirtualParticle
Abstract base class for a general particle.
Definition: VirtualParticle.h:23
erhic::EventPythia::F2
Double32_t F2
Value used for radiative corrections.
Definition: erhic/EventPythia.h:348
erhic::EventBeagle::crori
Double32_t crori
crossing angle orientation, +-1 lepton beam along +-z, +-2 hadron beam along +-z, 0 meas no crossing ...
Definition: erhic/EventPythia.h:565
erhic::EventBeagle::davg
Double32_t davg
Average density-weighted distance from all inelastic collisions to the edge of the nucleus.
Definition: erhic/EventPythia.h:580
erhic::EventPythia::trueNu
Double32_t trueNu
Generated nu of the event.
Definition: erhic/EventPythia.h:347
erhic::EventBeagle::User3
Double32_t User3
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:588
erhic::EventPythia::xbeamparton
Double32_t xbeamparton
Momentum fraction taken by the beam parton, see PARI(33)
Definition: erhic/EventPythia.h:322
erhic::EventBeagle::pznucl
Double32_t pznucl
target nucleon momentum
Definition: erhic/EventPythia.h:563
erhic::EventMC::BeamLepton
virtual const ParticleMC * BeamLepton() const
Returns a pointer to the incident lepton beam particle.
Definition: erhic/EventMC.cxx:122
erhic::EventBeagle::pztarg
Double32_t pztarg
target beam momentum
Definition: erhic/EventPythia.h:562
erhic::EventPythia::leptonphi
Double32_t leptonphi
Azimuthal angle of the scattered lepton in the cm frame.
Definition: erhic/EventPythia.h:327
erhic::EventPythia::EBrems
Double32_t EBrems
Energy per radiative photon in the nuclear rest frame.
Definition: erhic/EventPythia.h:338
erhic::EventPythia::tgtparton
Int_t tgtparton
PDG code of the struck parton in the hadron beam, see MSTI(16)
Definition: erhic/EventPythia.h:314
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::EventBeagle::Phib
Double32_t Phib
phi of impact parameter vector
Definition: erhic/EventPythia.h:567
erhic::EventBeagle::Ncolli
Int_t Ncolli
Number of inelastic collisions in target.
Definition: erhic/EventPythia.h:571
erhic::EventBeagle::ThickScl
Double32_t ThickScl
T(b)/rho0 in fm.
Definition: erhic/EventPythia.h:569
erhic::EventBeagle::Nwound
Int_t Nwound
Number of wounded nucleon including those in INC.
Definition: erhic/EventPythia.h:572
erhic::EventPythia::photonflux
Double32_t photonflux
Flux factor, see VINT(319)
Definition: erhic/EventPythia.h:340
erhic::EventPythia::Parse
virtual bool Parse(const std::string &)
Parses the event information from a text string.
Definition: erhic/EventPythia.cxx:50
erhic::EventMC::number
Int_t number
Event number.
Definition: erhic/EventMC.h:203
erhic::EventPythia::trueY
Double32_t trueY
Generated y of the event, see VINT(309)
Definition: erhic/EventPythia.h:341
erhic::EventBeagle::Aremn
Int_t Aremn
A of the nuclear remnant after evaporation and breakup.
Definition: erhic/EventPythia.h:576
erhic::EventBeagle::Eexc
Double32_t Eexc
Excitation energy in the nuclear remnant before evaporation and breakup.
Definition: erhic/EventPythia.h:584
erhic::EventBeagle::EventBeagle
EventBeagle(const std::string &str="")
Definition: erhic/EventPythia.cxx:98
erhic::EventPythia::F1
Double32_t F1
Value used for radiative corrections.
Definition: erhic/EventPythia.h:329
erhic::EventMC::nTracks
Int_t nTracks
Number of Particles in the event (intermediate + final)
Definition: erhic/EventMC.h:205
erhic::EventBeagle::pzf
Double32_t pzf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:583
erhic::EventBeagle::Thickness
Double32_t Thickness
T(b) in nucleons/fm^2.
Definition: erhic/EventPythia.h:568
erhic::EventBeagle::Parse
bool Parse(const std::string &)
Parses the event information from a text string.
Definition: erhic/EventPythia.cxx:137
erhic::EventPythia::trueW2
Double32_t trueW2
Generated W2 of the event.
Definition: erhic/EventPythia.h:346
erhic::EventPythia::pt2_hat
Double32_t pt2_hat
Squared pT of the hard subprocess, see PARI(18)
Definition: erhic/EventPythia.h:350
erhic::EventPythia::beamparton
Int_t beamparton
Parton interacting with the hadron beam in the case of resolved photon processes and soft VMD,...
Definition: erhic/EventPythia.h:316
erhic::EventBeagle::pzlep
Double32_t pzlep
lepton beam momentum
Definition: erhic/EventPythia.h:561
erhic::VirtualParticle::Id
virtual Pid Id() const =0
Returns identity information for the Particle species.