Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2024-11-16 09:02:05

0001 /**
0002  \file
0003  Declaration of class erhic::EventMC.
0004  
0005  \author    Thomas Burton
0006  \date      2011-10-10
0007  \copyright 2011 Brookhaven National Lab
0008  */
0009 
0010 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
0011 #define INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
0012 
0013 #include <string>
0014 #include <vector>
0015 
0016 #include <TClonesArray.h>
0017 #include <TLorentzVector.h>
0018 
0019 #include "eicsmear/erhic/EventDis.h"
0020 #include "eicsmear/erhic/ParticleMC.h"
0021 
0022 class TTree;
0023 
0024 namespace erhic {
0025 
0026 /**
0027  Abstract base class for DIS Monte Carlo events.
0028  Implements common event properties and methods.
0029  */
0030 class EventMC : public EventDis {
0031  public:
0032   /**
0033    Constructor.
0034    */
0035   EventMC();
0036 
0037   /**
0038    Destructor.
0039    */
0040   virtual ~EventMC();
0041 
0042   virtual bool RequiresEaParticleFields() { return false; };
0043 
0044   /**
0045    Returns a unique identifier for this event.
0046    */
0047   virtual ULong64_t GetN() const;
0048 
0049   /**
0050    Returns a code describing the production process of this event.
0051    */
0052   virtual Int_t GetProcess() const;
0053 
0054   /**
0055    Returns the number of tracks in the event.
0056    */
0057   virtual UInt_t GetNTracks() const;
0058 
0059   /**
0060    Returns the nth track.
0061    Returns NULL if the track number is out of the range [0, GetNTracks()).
0062    @param [in] The track index, in the range [0, GetNTracks()).
0063    */
0064   virtual const ParticleMC* GetTrack(UInt_t) const;
0065 
0066   /**
0067    \overload const ParticleMC* GetTrack(UInt_t) const
0068    */
0069   virtual ParticleMC* GetTrack(UInt_t);
0070 
0071   /**
0072    Returns a pointer to the incident lepton beam particle.
0073    Returns a NULL pointer if the particle cannot be located in the event.
0074    IMPORTANT - DO NOT DELETE THE POINTER OR BAD THINGS WILL HAPPEN!
0075    
0076    In the standard eRHIC Monte Carlo format, the incident lepton beam
0077    is assumed to be the first particle in the particle list.
0078    This is the behaviour implemented here.
0079    Derived classes can implement other selection mechanisms depending on
0080    their data format.
0081    */
0082   virtual const ParticleMC* BeamLepton() const;
0083 
0084   /**
0085    Returns a pointer to the incident hadron beam particle.
0086    See also notes in BeamLepton().
0087    
0088    In the standard eRHIC Monte Carlo format, the incident hadron beam
0089    is assumed to be the second particle in the particle list.
0090    */
0091   virtual const ParticleMC* BeamHadron() const;
0092 
0093   /**
0094    Returns a pointer to the exchanged boson.
0095    See also notes in BeamLepton().
0096    
0097    In the standard eRHIC Monte Carlo format, the exchanged boson
0098    is assumed to be the fourth particle in the particle list.
0099    */
0100   virtual const ParticleMC* ExchangeBoson() const;
0101 
0102   /**
0103    Returns a pointer to the lepton beam particle after scattering.
0104    See also notes in BeamLepton().
0105    
0106    In the standard eRHIC Monte Carlo format, the scattered lepton beam
0107    is assumed to be the first final-state particle in the particle list
0108    with the same PDG code as the incident lepton beam.
0109 
0110    Please overwrite this method accordingly!
0111    By default, it will simply use the fourth particle in the particle list.
0112    See e.g. EventPythia, EventSimple.
0113    */
0114   virtual const ParticleMC* ScatteredLepton() const;
0115 
0116   /**
0117    Populates the event-wise variables from a string.
0118    Does not populate the particle list or compute derived quantities.
0119    See also Compute().
0120    */
0121   virtual bool Parse(const std::string&) = 0;
0122 
0123   /**
0124    Add a copy of a track argument to the end of the track list.
0125    @param [in] Pointer to the track to add.
0126    */
0127   virtual void AddLast(ParticleMC* track);
0128 
0129   /**
0130    Resets event properties to defaults.
0131    Does not clear particle list - use Clear() for that.
0132    */
0133   virtual void Reset();
0134 
0135   /**
0136      Quick Event list
0137    */
0138   void Print( const Option_t *option="" ) const;
0139 
0140   /**
0141    Clears event contents.
0142    Event properties are reset to defaults and track list
0143    is deleted.
0144    */
0145   virtual void Clear(Option_t* = "");
0146 
0147   /**
0148    Sets the code describing the production process of this event.
0149    @param [in] code The identifying code, an integer
0150    */
0151   virtual void SetProcess(int code);
0152 
0153   /**
0154    Sets the unique identifier for this event.
0155    @param [in] n The identifying number, an integer
0156    */
0157   virtual void SetN(int n);
0158 
0159   /**
0160    Sets the track count for this event.
0161    @param [in] n The track count, an integer
0162    */
0163   virtual void SetNTracks(int n);
0164 
0165   /**
0166    Set incident lepton energy in the nuclear rest frame.
0167    */
0168   virtual void SetELeptonInNuclearFrame(double energy);
0169 
0170   /**
0171    Set scattered lepton energy in the nuclear rest frame.
0172    */
0173   virtual void SetEScatteredInNuclearFrame(double energy);
0174 
0175   /**
0176    Stores pointers to all final state particles in the list.
0177    These pointers should not be deleted by the user.
0178    Any existing entries in the list are not changed.
0179    @param [out] particles The list in which to store particles.
0180    */
0181   void FinalState(ParticlePtrList& particles) const;
0182 
0183   /**
0184    Yields all particles that belong to the hadronic final state.
0185    This is the same as the result of FinalState(), minus the scattered
0186    beam lepton (i.e. including leptons and bosons).
0187    */
0188   void HadronicFinalState(ParticlePtrList&) const;
0189 
0190   /**
0191    Returns the total momentum of the final state in GeV/c.
0192    */
0193   TLorentzVector FinalStateMomentum() const;
0194 
0195   /**
0196    Returns the total momentum of the hadronic final state in GeV/c.
0197    */
0198   TLorentzVector HadronicFinalStateMomentum() const;
0199 
0200   /**
0201    Returns the total charge of the final state in units of e.
0202    */
0203   Double_t FinalStateCharge() const;
0204 
0205   /**
0206    Returns pointers to all tracks in the event.
0207    Do not delete the pointers.
0208    */
0209   std::vector<const VirtualParticle*> GetTracks() const;
0210 
0211  protected:
0212   Int_t number;  ///< Event number
0213   Int_t process;  ///< PYTHIA code for the physics process producing the event
0214   Int_t nTracks;  ///< Number of Particles in the event (intermediate + final)
0215   Double32_t ELeptonInNucl;  ///< Incident lepton energy in the
0216                              ///< nuclear rest frame
0217   Double32_t ELeptonOutNucl;  ///< Scattered lepton energy in the
0218                               ///< nuclear rest frame
0219   TClonesArray particles;  ///< Particle list
0220 
0221   ClassDef(erhic::EventMC, 2)
0222 };
0223 
0224 inline ULong64_t EventMC::GetN() const {
0225   return number;
0226 }
0227 
0228 inline Int_t EventMC::GetProcess() const {
0229   return process;
0230 }
0231 
0232 inline UInt_t EventMC::GetNTracks() const {
0233   return particles.GetEntries();
0234 }
0235 
0236 inline const ParticleMC* EventMC::GetTrack(UInt_t u) const {
0237   if (u < (UInt_t)particles.GetEntries()) {
0238     return static_cast<ParticleMC*>(particles.At(u));
0239   } else {
0240     return NULL;
0241   }  // if
0242 }
0243 
0244 inline ParticleMC* EventMC::GetTrack(UInt_t u) {
0245   if (u < (UInt_t)particles.GetEntries()) {
0246     return static_cast<ParticleMC*>(particles.At(u));
0247   } else {
0248     return NULL;
0249   }  // if
0250 }
0251 
0252 inline void EventMC::SetProcess(int code) {
0253   process = code;
0254 }
0255 
0256 inline void EventMC::SetN(int n) {
0257   number = n;
0258 }
0259 
0260 inline void EventMC::SetNTracks(int n) {
0261   nTracks = n;
0262 }
0263 
0264 inline void EventMC::SetELeptonInNuclearFrame(double e) {
0265   ELeptonInNucl = e;
0266 }
0267 
0268 inline void EventMC::SetEScatteredInNuclearFrame(double e) {
0269   ELeptonOutNucl = e;
0270 }
0271 
0272 /**
0273  Wrapper for getting tree from file and event from tree.
0274  */
0275 class Reader {
0276  public:
0277   /**
0278    Construct a Reader for the named TTree.
0279    */
0280   explicit Reader(const std::string& treeName = "EICTree");
0281 
0282   /**
0283    Destructor.
0284    */
0285   virtual ~Reader() { }
0286 
0287   /**
0288    Read and return the numbered event from the tree.
0289    
0290    The event number should lie in the range [0, numEvents).
0291    */
0292   EventMC* Read(Long64_t number);
0293 
0294   /**
0295    Read and return the numbered event from the tree.
0296    
0297    The event number should lie in the range [0, numEvents).
0298    */
0299   EventMC* operator()(Long64_t number);
0300 
0301   /**
0302    Return the tree read by this reader.
0303    */
0304   TTree* GetTree();
0305 
0306   EventMC* mEvent;  ///< The last event read
0307   TTree* mTree;   ///< The tree being read
0308 
0309   ClassDef(erhic::Reader, 1)
0310 };
0311 
0312 inline EventMC* Reader::operator()(Long64_t i) {
0313   return Read(i);
0314 }
0315 
0316 inline TTree* Reader::GetTree() {
0317   return mTree;
0318 }
0319 
0320 }  // namespace erhic
0321 
0322 #endif  // INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_