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.h
Go to the documentation of this file.
1 
10 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
11 #define INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
12 
13 #include <string>
14 #include <vector>
15 
16 #include <TClonesArray.h>
17 #include <TLorentzVector.h>
18 
21 
22 class TTree;
23 
24 namespace erhic {
25 
30 class EventMC : public EventDis {
31  public:
35  EventMC();
36 
40  virtual ~EventMC();
41 
42  virtual bool RequiresEaParticleFields() { return false; };
43 
47  virtual ULong64_t GetN() const;
48 
52  virtual Int_t GetProcess() const;
53 
57  virtual UInt_t GetNTracks() const;
58 
64  virtual const ParticleMC* GetTrack(UInt_t) const;
65 
69  virtual ParticleMC* GetTrack(UInt_t);
70 
82  virtual const ParticleMC* BeamLepton() const;
83 
91  virtual const ParticleMC* BeamHadron() const;
92 
100  virtual const ParticleMC* ExchangeBoson() const;
101 
110  virtual const ParticleMC* ScatteredLepton() const;
111 
117  virtual bool Parse(const std::string&) = 0;
118 
123  virtual void AddLast(ParticleMC* track);
124 
129  virtual void Reset();
130 
136  virtual void Clear(Option_t* = "");
137 
142  virtual void SetProcess(int code);
143 
148  virtual void SetN(int n);
149 
154  virtual void SetNTracks(int n);
155 
159  virtual void SetELeptonInNuclearFrame(double energy);
160 
164  virtual void SetEScatteredInNuclearFrame(double energy);
165 
172  void FinalState(ParticlePtrList& particles) const;
173 
179  void HadronicFinalState(ParticlePtrList&) const;
180 
184  TLorentzVector FinalStateMomentum() const;
185 
189  TLorentzVector HadronicFinalStateMomentum() const;
190 
194  Double_t FinalStateCharge() const;
195 
200  std::vector<const VirtualParticle*> GetTracks() const;
201 
202  protected:
203  Int_t number;
204  Int_t process;
205  Int_t nTracks;
206  Double32_t ELeptonInNucl;
207  Double32_t ELeptonOutNucl;
209  TClonesArray particles;
211 
212  ClassDef(erhic::EventMC, 2)
213 };
214 
215 inline ULong64_t EventMC::GetN() const {
216  return number;
217 }
218 
219 inline Int_t EventMC::GetProcess() const {
220  return process;
221 }
222 
223 inline UInt_t EventMC::GetNTracks() const {
224  return particles.GetEntries();
225 }
226 
227 inline const ParticleMC* EventMC::GetTrack(UInt_t u) const {
228  if (u < (UInt_t)particles.GetEntries()) {
229  return static_cast<ParticleMC*>(particles.At(u));
230  } else {
231  return NULL;
232  } // if
233 }
234 
235 inline ParticleMC* EventMC::GetTrack(UInt_t u) {
236  if (u < (UInt_t)particles.GetEntries()) {
237  return static_cast<ParticleMC*>(particles.At(u));
238  } else {
239  return NULL;
240  } // if
241 }
242 
243 inline void EventMC::SetProcess(int code) {
244  process = code;
245 }
246 
247 inline void EventMC::SetN(int n) {
248  number = n;
249 }
250 
251 inline void EventMC::SetNTracks(int n) {
252  nTracks = n;
253 }
254 
255 inline void EventMC::SetELeptonInNuclearFrame(double e) {
256  ELeptonInNucl = e;
257 }
258 
260  ELeptonOutNucl = e;
261 }
262 
266 class Reader {
267  public:
271  explicit Reader(const std::string& treeName = "EICTree");
272 
276  virtual ~Reader() { }
277 
283  EventMC* Read(Long64_t number);
284 
290  EventMC* operator()(Long64_t number);
291 
295  TTree* GetTree();
296 
298  TTree* mTree;
299 
300  ClassDef(erhic::Reader, 1)
301 };
302 
303 inline EventMC* Reader::operator()(Long64_t i) {
304  return Read(i);
305 }
306 
307 inline TTree* Reader::GetTree() {
308  return mTree;
309 }
310 
311 } // namespace erhic
312 
313 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
erhic::Reader::~Reader
virtual ~Reader()
Destructor.
Definition: erhic/EventMC.h:276
erhic::EventMC::Parse
virtual bool Parse(const std::string &)=0
Populates the event-wise variables from a string.
erhic::EventMC::SetProcess
virtual void SetProcess(int code)
Sets the code describing the production process of this event.
Definition: erhic/EventMC.h:243
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::VirtualEvent::ParticlePtrList
std::vector< const erhic::VirtualParticle * > ParticlePtrList
typedef for a track pointer collection.
Definition: VirtualEvent.h:52
erhic::EventMC::GetN
virtual ULong64_t GetN() const
Returns a unique identifier for this event.
Definition: erhic/EventMC.h:215
erhic::Reader::mEvent
EventMC * mEvent
The last event read.
Definition: erhic/EventMC.h:297
erhic::EventMC::SetN
virtual void SetN(int n)
Sets the unique identifier for this event.
Definition: erhic/EventMC.h:247
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
EventDis.h
erhic::Reader::Read
EventMC * Read(Long64_t number)
Read and return the numbered event from the tree.
Definition: erhic/EventMC.cxx:166
erhic::EventMC::GetProcess
virtual Int_t GetProcess() const
Returns a code describing the production process of this event.
Definition: erhic/EventMC.h:219
erhic::Reader::mTree
TTree * mTree
The tree being read.
Definition: erhic/EventMC.h:298
erhic::Reader::operator()
EventMC * operator()(Long64_t number)
Read and return the numbered event from the tree.
Definition: erhic/EventMC.h:303
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::RequiresEaParticleFields
virtual bool RequiresEaParticleFields()
Definition: erhic/EventMC.h:42
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::SetNTracks
virtual void SetNTracks(int n)
Sets the track count for this event.
Definition: erhic/EventMC.h:251
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
A deeply inelastic scattering event.
Definition: EventDis.h:37
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::EventMC::number
Int_t number
Event number.
Definition: erhic/EventMC.h:203
erhic::EventMC::SetELeptonInNuclearFrame
virtual void SetELeptonInNuclearFrame(double energy)
Set incident lepton energy in the nuclear rest frame.
Definition: erhic/EventMC.h:255
erhic::Reader::GetTree
TTree * GetTree()
Return the tree read by this reader.
Definition: erhic/EventMC.h:307
erhic::Reader::Reader
Reader(const std::string &treeName="EICTree")
Construct a Reader for the named TTree.
Definition: erhic/EventMC.cxx:159
erhic::Reader
Wrapper for getting tree from file and event from tree.
Definition: erhic/EventMC.h:266
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::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::SetEScatteredInNuclearFrame
virtual void SetEScatteredInNuclearFrame(double energy)
Set scattered lepton energy in the nuclear rest frame.
Definition: erhic/EventMC.h:259
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
erhic::EventMC::GetTracks
std::vector< const VirtualParticle * > GetTracks() const
Returns pointers to all tracks in the event.
Definition: erhic/EventMC.cxx:48
erhic::EventMC::GetNTracks
virtual UInt_t GetNTracks() const
Returns the number of tracks in the event.
Definition: erhic/EventMC.h:223
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