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/EventMC.cxx
Go to the documentation of this file.
1 
10 #include "eicsmear/erhic/EventMC.h"
11 
12 #include <iostream>
13 #include <list>
14 #include <string>
15 #include <vector>
16 
17 #include <TCollection.h>
18 #include <TDatabasePDG.h>
19 #include <TDirectory.h>
20 #include <TParticlePDG.h>
21 #include <TTree.h>
22 
25 
26 namespace erhic {
27 
28 // We use these vectors/iterators a lot, so define
29 // some typedefs for convenience.
30 typedef std::vector<const VirtualParticle*> TrackVector;
31 typedef TrackVector::iterator TrackVectorIter;
32 typedef TrackVector::const_iterator TrackVectorCIter;
33 
35 : number(-1)
36 , process(-1)
37 , nTracks(-1)
38 , ELeptonInNucl(NAN)
39 , ELeptonOutNucl(NAN)
40 , particles("erhic::ParticleMC", 100) {
41 }
42 
44  // No memory to clear. The TClonesArray takes care of the
45  // particles when it is destroyed.
46 }
47 
49  TrackVector tracks;
50  TObject* object(NULL);
51  TIter next(&particles);
52  while ((object = next())) {
53  tracks.push_back(static_cast<ParticleMC*>(object));
54  } // while
55  return tracks;
56 }
57 
59  // This is simple but not very efficient.
60  // Copy vector to a list and use list::remove to get rid
61  // of the scattered lepton. Then copy this list back to
62  // the vector.
63  FinalState(final_);
64  std::list<const VirtualParticle*> plist(final_.begin(),
65  final_.end());
66  plist.remove(ScatteredLepton());
67  final_ = TrackVector(plist.begin(), plist.end());
68 }
69 
70 // Get the particles that belong to the hadronic final state.
71 // The stored Particle* are pointers to the original particles in the event
72 // so don't delete them!
73 void EventMC::FinalState(TrackVector& final_) const {
74  final_.clear();
75  TIter next(&particles);
76  ParticleMC* p(NULL);
77  while ((p = static_cast<ParticleMC*>(next()))) {
78  if (1 == p->GetStatus()) {
79  final_.push_back(p);
80  } // if
81  } // while
82 }
83 
84 TLorentzVector EventMC::FinalStateMomentum() const {
85  TrackVector final_;
86  FinalState(final_);
87  TLorentzVector mom;
88  for (TrackVectorCIter i = final_.begin(); i != final_.end(); ++i) {
89  mom += (*i)->Get4Vector();
90  } // for
91 
92  return mom;
93 }
94 
95 TLorentzVector EventMC::HadronicFinalStateMomentum() const {
96  TrackVector final_;
97  HadronicFinalState(final_);
98  TLorentzVector mom;
99  for (TrackVectorCIter i = final_.begin(); i != final_.end(); ++i) {
100  mom += (*i)->Get4Vector();
101  } // for
102 
103  return mom;
104 }
105 
106 Double_t EventMC::FinalStateCharge() const {
107  TrackVector final_;
108  FinalState(final_);
109  Double_t charge(0);
110  TDatabasePDG* pdg = TDatabasePDG::Instance();
111  for (TrackVectorCIter i = final_.begin(); i != final_.end(); ++i) {
112  TParticlePDG* part = pdg->GetParticle((*i)->Id());
113  if (part) {
114  charge += part->Charge() / 3.;
115  } else {
116  std::cout << "Unknown particle: " << (*i)->Id() << std::endl;
117  } // if
118  } // for
119  return charge;
120 }
121 
123  return GetTrack(0);
124 }
125 
127  return GetTrack(1);
128 }
129 
131  return GetTrack(3);
132 }
133 
135  return GetTrack(2);
136 }
137 
139  number = -1;
140  process = -1;
141  nTracks = -1;
142  x = QSquared = y = WSquared = nu = ELeptonInNucl = ELeptonOutNucl = NAN;
143 }
144 
145 void EventMC::Clear(Option_t* /*option*/) {
146  Reset();
147  particles.Clear();
148 }
149 
151  new(particles[particles.GetEntries()]) ParticleMC(*track);
152  nTracks = particles.GetEntries();
153 }
154 
155 //
156 // class Reader
157 //
158 
159 Reader::Reader(const std::string& treeName)
160 : mEvent(NULL)
161 , mTree(NULL) {
162  if (gDirectory) gDirectory->GetObject(treeName.c_str(), mTree);
163  if (mTree) mTree->SetBranchAddress("event", &mEvent);
164 }
165 
166 EventMC* Reader::Read(Long64_t i) {
167  EventMC* event(NULL);
168  if (mTree) {
169  mTree->GetEntry(i);
170  event = mEvent;
171  } // if
172  return event;
173 }
174 
175 } // namespace erhic
erhic::EventDis::WSquared
Double32_t WSquared
Invariant mass of the hadronic system.
Definition: EventDis.h:184
erhic::EventMC::ExchangeBoson
virtual const ParticleMC * ExchangeBoson() const
Returns a pointer to the exchanged boson.
Definition: erhic/EventMC.cxx:130
erhic::EventMC::EventMC
EventMC()
Constructor.
Definition: erhic/EventMC.cxx:34
erhic
Definition: EventDis.cxx:14
erhic::EventDis::QSquared
Double32_t QSquared
Q2 calculated from scattered electron.
Definition: EventDis.h:182
erhic::Reader::mEvent
EventMC * mEvent
The last event read.
Definition: erhic/EventMC.h:297
erhic::EventMC::process
Int_t process
PYTHIA code for the physics process producing the event.
Definition: erhic/EventMC.h:204
erhic::EventMC::Clear
virtual void Clear(Option_t *="")
Clears event contents.
Definition: erhic/EventMC.cxx:145
erhic::EventMC::FinalStateMomentum
TLorentzVector FinalStateMomentum() const
Returns the total momentum of the final state in GeV/c.
Definition: erhic/EventMC.cxx:84
erhic::Reader::Read
EventMC * Read(Long64_t number)
Read and return the numbered event from the tree.
Definition: erhic/EventMC.cxx:166
EventMC.h
erhic::ParticleMCbase::GetStatus
virtual UShort_t GetStatus() const
Returns the status of the particle.
Definition: erhic/ParticleMC.h:461
erhic::Reader::mTree
TTree * mTree
The tree being read.
Definition: erhic/EventMC.h:298
erhic::EventMC::ScatteredLepton
virtual const ParticleMC * ScatteredLepton() const
Returns a pointer to the lepton beam particle after scattering.
Definition: erhic/EventMC.cxx:134
erhic::EventMC::FinalState
void FinalState(ParticlePtrList &particles) const
Stores pointers to all final state particles in the list.
Definition: erhic/EventMC.cxx:73
erhic::EventMC::BeamHadron
virtual const ParticleMC * BeamHadron() const
Returns a pointer to the incident hadron beam particle.
Definition: erhic/EventMC.cxx:126
erhic::EventMC::ELeptonOutNucl
Double32_t ELeptonOutNucl
Scattered lepton energy in the nuclear rest frame.
Definition: erhic/EventMC.h:208
erhic::EventMC::ELeptonInNucl
Double32_t ELeptonInNucl
Incident lepton energy in the nuclear rest frame.
Definition: erhic/EventMC.h:206
erhic::EventMC::BeamLepton
virtual const ParticleMC * BeamLepton() const
Returns a pointer to the incident lepton beam particle.
Definition: erhic/EventMC.cxx:122
erhic::EventDis::nu
Double32_t nu
Energy transfer from the electron.
Definition: EventDis.h:185
erhic::TrackVectorCIter
TrackVector::const_iterator TrackVectorCIter
Definition: erhic/EventMC.cxx:32
erhic::EventDis::x
Double32_t x
Bjorken scaling variable.
Definition: EventDis.h:181
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::TrackVector
std::vector< const VirtualParticle * > TrackVector
Definition: erhic/EventMC.cxx:30
erhic::EventMC::number
Int_t number
Event number.
Definition: erhic/EventMC.h:203
erhic::TrackVectorIter
TrackVector::iterator TrackVectorIter
Definition: erhic/EventMC.cxx:31
erhic::Reader::Reader
Reader(const std::string &treeName="EICTree")
Construct a Reader for the named TTree.
Definition: erhic/EventMC.cxx:159
erhic::EventMC::AddLast
virtual void AddLast(ParticleMC *track)
Add a copy of a track argument to the end of the track list.
Definition: erhic/EventMC.cxx:150
erhic::EventMC::nTracks
Int_t nTracks
Number of Particles in the event (intermediate + final)
Definition: erhic/EventMC.h:205
erhic::EventDis::y
Double32_t y
Inelasticity.
Definition: EventDis.h:183
erhic::EventMC::GetTrack
virtual const ParticleMC * GetTrack(UInt_t) const
Returns the nth track.
Definition: erhic/EventMC.h:227
ParticleMC.h
erhic::EventMC::HadronicFinalState
void HadronicFinalState(ParticlePtrList &) const
Yields all particles that belong to the hadronic final state.
Definition: erhic/EventMC.cxx:58
erhic::EventMC::particles
TClonesArray particles
Particle list.
Definition: erhic/EventMC.h:210
erhic::EventMC::HadronicFinalStateMomentum
TLorentzVector HadronicFinalStateMomentum() const
Returns the total momentum of the hadronic final state in GeV/c.
Definition: erhic/EventMC.cxx:95
erhic::EventMC::FinalStateCharge
Double_t FinalStateCharge() const
Returns the total charge of the final state in units of e.
Definition: erhic/EventMC.cxx:106
ParticleIdentifier.h
erhic::EventMC::GetTracks
std::vector< const VirtualParticle * > GetTracks() const
Returns pointers to all tracks in the event.
Definition: erhic/EventMC.cxx:48
erhic::EventMC
Abstract base class for DIS Monte Carlo events.
Definition: erhic/EventMC.h:30
erhic::EventMC::~EventMC
virtual ~EventMC()
Destructor.
Definition: erhic/EventMC.cxx:43
erhic::EventMC::Reset
virtual void Reset()
Resets event properties to defaults.
Definition: erhic/EventMC.cxx:138