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
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 
45  virtual ULong64_t GetN() const;
46 
50  virtual Int_t GetProcess() const;
51 
55  virtual UInt_t GetNTracks() const;
56 
62  virtual const ParticleMC* GetTrack(UInt_t) const;
63 
67  virtual ParticleMC* GetTrack(UInt_t);
68 
80  virtual const ParticleMC* BeamLepton() const;
81 
89  virtual const ParticleMC* BeamHadron() const;
90 
98  virtual const ParticleMC* ExchangeBoson() const;
99 
108  virtual const ParticleMC* ScatteredLepton() const;
109 
115  virtual bool Parse(const std::string&) = 0;
116 
121  virtual void AddLast(ParticleMC* track);
122 
127  virtual void Reset();
128 
134  virtual void Clear(Option_t* = "");
135 
140  virtual void SetProcess(int code);
141 
146  virtual void SetN(int n);
147 
152  virtual void SetNTracks(int n);
153 
157  virtual void SetELeptonInNuclearFrame(double energy);
158 
162  virtual void SetEScatteredInNuclearFrame(double energy);
163 
170  void FinalState(ParticlePtrList& particles) const;
171 
177  void HadronicFinalState(ParticlePtrList&) const;
178 
182  TLorentzVector FinalStateMomentum() const;
183 
187  TLorentzVector HadronicFinalStateMomentum() const;
188 
192  Double_t FinalStateCharge() const;
193 
198  std::vector<const VirtualParticle*> GetTracks() const;
199 
200  protected:
201  Int_t number;
202  Int_t process;
203  Int_t nTracks;
204  Double32_t ELeptonInNucl;
205  Double32_t ELeptonOutNucl;
207  TClonesArray particles;
209 
210  ClassDef(erhic::EventMC, 2)
211 };
212 
213 inline ULong64_t EventMC::GetN() const {
214  return number;
215 }
216 
217 inline Int_t EventMC::GetProcess() const {
218  return process;
219 }
220 
221 inline UInt_t EventMC::GetNTracks() const {
222  return particles.GetEntries();
223 }
224 
225 inline const ParticleMC* EventMC::GetTrack(UInt_t u) const {
226  if (u < (UInt_t)particles.GetEntries()) {
227  return static_cast<ParticleMC*>(particles.At(u));
228  } else {
229  return NULL;
230  } // if
231 }
232 
233 inline ParticleMC* EventMC::GetTrack(UInt_t u) {
234  if (u < (UInt_t)particles.GetEntries()) {
235  return static_cast<ParticleMC*>(particles.At(u));
236  } else {
237  return NULL;
238  } // if
239 }
240 
241 inline void EventMC::SetProcess(int code) {
242  process = code;
243 }
244 
245 inline void EventMC::SetN(int n) {
246  number = n;
247 }
248 
249 inline void EventMC::SetNTracks(int n) {
250  nTracks = n;
251 }
252 
253 inline void EventMC::SetELeptonInNuclearFrame(double e) {
254  ELeptonInNucl = e;
255 }
256 
258  ELeptonOutNucl = e;
259 }
260 
264 class Reader {
265  public:
269  explicit Reader(const std::string& treeName = "EICTree");
270 
274  virtual ~Reader() { }
275 
281  EventMC* Read(Long64_t number);
282 
288  EventMC* operator()(Long64_t number);
289 
293  TTree* GetTree();
294 
296  TTree* mTree;
297 
298  ClassDef(erhic::Reader, 1)
299 };
300 
301 inline EventMC* Reader::operator()(Long64_t i) {
302  return Read(i);
303 }
304 
305 inline TTree* Reader::GetTree() {
306  return mTree;
307 }
308 
309 } // namespace erhic
310 
311 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
virtual const ParticleMC * ScatteredLepton() const
Definition: EventMC.cxx:134
virtual void Clear(Option_t *="")
Definition: EventMC.cxx:145
virtual void SetN(int n)
Definition: EventMC.h:245
virtual void SetNTracks(int n)
Definition: EventMC.h:249
virtual bool Parse(const std::string &)=0
std::vector< const erhic::VirtualParticle * > ParticlePtrList
Definition: VirtualEvent.h:52
Double_t FinalStateCharge() const
Definition: EventMC.cxx:106
virtual UInt_t GetNTracks() const
Definition: EventMC.h:221
virtual void SetEScatteredInNuclearFrame(double energy)
Definition: EventMC.h:257
Int_t nTracks
Number of Particles in the event (intermediate + final)
Definition: EventMC.h:203
virtual ~Reader()
Definition: EventMC.h:274
virtual void Reset()
Definition: EventMC.cxx:138
TTree * GetTree()
Definition: EventMC.h:305
Int_t number
Event number.
Definition: EventMC.h:201
TClonesArray particles
Particle list.
Definition: EventMC.h:208
EventMC * mEvent
The last event read.
Definition: EventMC.h:295
Int_t process
PYTHIA code for the physics process producing the event.
Definition: EventMC.h:202
void HadronicFinalState(ParticlePtrList &) const
Definition: EventMC.cxx:58
Reader(const std::string &treeName="EICTree")
Definition: EventMC.cxx:159
virtual ULong64_t GetN() const
Definition: EventMC.h:213
virtual void AddLast(ParticleMC *track)
Definition: EventMC.cxx:150
virtual const ParticleMC * ExchangeBoson() const
Definition: EventMC.cxx:130
virtual const ParticleMC * GetTrack(UInt_t) const
Definition: EventMC.h:225
EventMC * Read(Long64_t number)
Definition: EventMC.cxx:166
TTree * mTree
The tree being read.
Definition: EventMC.h:296
virtual const ParticleMC * BeamLepton() const
Definition: EventMC.cxx:122
TLorentzVector HadronicFinalStateMomentum() const
Definition: EventMC.cxx:95
EventMC * operator()(Long64_t number)
Definition: EventMC.h:301
virtual const ParticleMC * BeamHadron() const
Definition: EventMC.cxx:126
virtual void SetProcess(int code)
Definition: EventMC.h:241
TLorentzVector FinalStateMomentum() const
Definition: EventMC.cxx:84
void FinalState(ParticlePtrList &particles) const
Definition: EventMC.cxx:73
virtual Int_t GetProcess() const
Definition: EventMC.h:217
virtual void SetELeptonInNuclearFrame(double energy)
Definition: EventMC.h:253
Double32_t ELeptonOutNucl
Definition: EventMC.h:206
virtual ~EventMC()
Definition: EventMC.cxx:43
std::vector< const VirtualParticle * > GetTracks() const
Definition: EventMC.cxx:48
Double32_t ELeptonInNucl
Definition: EventMC.h:204