EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EventMC.h
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file EventMC.h
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 
114  virtual const ParticleMC* ScatteredLepton() const;
115 
121  virtual bool Parse(const std::string&) = 0;
122 
127  virtual void AddLast(ParticleMC* track);
128 
133  virtual void Reset();
134 
138  void Print( const Option_t *option="" ) const;
139 
145  virtual void Clear(Option_t* = "");
146 
151  virtual void SetProcess(int code);
152 
157  virtual void SetN(int n);
158 
163  virtual void SetNTracks(int n);
164 
168  virtual void SetELeptonInNuclearFrame(double energy);
169 
173  virtual void SetEScatteredInNuclearFrame(double energy);
174 
181  void FinalState(ParticlePtrList& particles) const;
182 
188  void HadronicFinalState(ParticlePtrList&) const;
189 
193  TLorentzVector FinalStateMomentum() const;
194 
198  TLorentzVector HadronicFinalStateMomentum() const;
199 
203  Double_t FinalStateCharge() const;
204 
209  std::vector<const VirtualParticle*> GetTracks() const;
210 
211  protected:
212  Int_t number;
213  Int_t process;
214  Int_t nTracks;
215  Double32_t ELeptonInNucl;
216 
217  Double32_t ELeptonOutNucl;
218 
219  TClonesArray particles;
220 
221  ClassDef(erhic::EventMC, 2)
222 };
223 
224 inline ULong64_t EventMC::GetN() const {
225  return number;
226 }
227 
228 inline Int_t EventMC::GetProcess() const {
229  return process;
230 }
231 
232 inline UInt_t EventMC::GetNTracks() const {
233  return particles.GetEntries();
234 }
235 
236 inline const ParticleMC* EventMC::GetTrack(UInt_t u) const {
237  if (u < (UInt_t)particles.GetEntries()) {
238  return static_cast<ParticleMC*>(particles.At(u));
239  } else {
240  return NULL;
241  } // if
242 }
243 
244 inline ParticleMC* EventMC::GetTrack(UInt_t u) {
245  if (u < (UInt_t)particles.GetEntries()) {
246  return static_cast<ParticleMC*>(particles.At(u));
247  } else {
248  return NULL;
249  } // if
250 }
251 
252 inline void EventMC::SetProcess(int code) {
253  process = code;
254 }
255 
256 inline void EventMC::SetN(int n) {
257  number = n;
258 }
259 
260 inline void EventMC::SetNTracks(int n) {
261  nTracks = n;
262 }
263 
265  ELeptonInNucl = e;
266 }
267 
269  ELeptonOutNucl = e;
270 }
271 
275 class Reader {
276  public:
280  explicit Reader(const std::string& treeName = "EICTree");
281 
285  virtual ~Reader() { }
286 
292  EventMC* Read(Long64_t number);
293 
299  EventMC* operator()(Long64_t number);
300 
304  TTree* GetTree();
305 
307  TTree* mTree;
308 
309  ClassDef(erhic::Reader, 1)
310 };
311 
312 inline EventMC* Reader::operator()(Long64_t i) {
313  return Read(i);
314 }
315 
316 inline TTree* Reader::GetTree() {
317  return mTree;
318 }
319 
320 } // namespace erhic
321 
322 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_