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 Class Referenceabstract

#include <EventMC.h>

Inheritance diagram for erhic::EventMC:
erhic::EventDis erhic::VirtualEvent erhic::EventDjangoh erhic::EventDpmjet erhic::EventGmcTrans erhic::EventMilou erhic::EventPepsi erhic::EventPythia erhic::EventRapgap

Public Member Functions

 EventMC ()
 
virtual ~EventMC ()
 
virtual ULong64_t GetN () const
 
virtual Int_t GetProcess () const
 
virtual UInt_t GetNTracks () const
 
virtual const ParticleMCGetTrack (UInt_t) const
 
virtual ParticleMCGetTrack (UInt_t)
 
virtual const ParticleMCBeamLepton () const
 
virtual const ParticleMCBeamHadron () const
 
virtual const ParticleMCExchangeBoson () const
 
virtual const ParticleMCScatteredLepton () const
 
virtual bool Parse (const std::string &)=0
 
virtual void AddLast (ParticleMC *track)
 
virtual void Reset ()
 
virtual void Clear (Option_t *="")
 
virtual void SetProcess (int code)
 
virtual void SetN (int n)
 
virtual void SetNTracks (int n)
 
virtual void SetELeptonInNuclearFrame (double energy)
 
virtual void SetEScatteredInNuclearFrame (double energy)
 
void FinalState (ParticlePtrList &particles) const
 
void HadronicFinalState (ParticlePtrList &) const
 
TLorentzVector FinalStateMomentum () const
 
TLorentzVector HadronicFinalStateMomentum () const
 
Double_t FinalStateCharge () const
 
std::vector< const
VirtualParticle * > 
GetTracks () const
 
- Public Member Functions inherited from erhic::EventDis
virtual ~EventDis ()
 
 EventDis ()
 
 EventDis (const EventDis &)
 
EventDisoperator= (const EventDis &)
 
virtual Double_t GetX () const
 
virtual Double_t GetQ2 () const
 
virtual Double_t GetY () const
 
virtual Double_t GetYPlus () const
 
virtual Double_t GetW2 () const
 
virtual Double_t GetNu () const
 
virtual double GetXDoubleAngle () const
 
virtual double GetQ2DoubleAngle () const
 
virtual double GetYDoubleAngle () const
 
virtual double GetW2DoubleAngle () const
 
virtual double GetXJacquetBlondel () const
 
virtual double GetQ2JacquetBlondel () const
 
virtual double GetYJacquetBlondel () const
 
virtual double GetW2JacquetBlondel () const
 
virtual void SetLeptonKinematics (const DisKinematics &)
 
virtual void SetJacquetBlondelKinematics (const DisKinematics &)
 
virtual void SetDoubleAngleKinematics (const DisKinematics &)
 
virtual void CopyKinematics (const EventDis &)
 
- Public Member Functions inherited from erhic::VirtualEvent
virtual ~VirtualEvent ()
 

Protected Attributes

Int_t number
 Event number.
 
Int_t process
 PYTHIA code for the physics process producing the event.
 
Int_t nTracks
 Number of Particles in the event (intermediate + final)
 
Double32_t ELeptonInNucl
 
Double32_t ELeptonOutNucl
 
TClonesArray particles
 Particle list.
 

Additional Inherited Members

- Public Types inherited from erhic::VirtualEvent
typedef std::vector< const
erhic::VirtualParticle * > 
ParticlePtrList
 
- Public Attributes inherited from erhic::EventDis
Double32_t x
 Bjorken scaling variable.
 
Double32_t QSquared
 Q2 calculated from scattered electron.
 
Double32_t y
 Inelasticity.
 
Double32_t WSquared
 Invariant mass of the hadronic system.
 
Double32_t nu
 Energy transfer from the electron.
 
Double32_t yJB
 y calculated via the Jacquet-Blondel method
 
Double32_t QSquaredJB
 Q2 calculated via the Jacquet-Blondel method.
 
Double32_t xJB
 x calculated via the Jacquet-Blondel method
 
Double32_t WSquaredJB
 W2 calculated via the Jacquet-Blondel method.
 
Double32_t yDA
 y calculated via the double-angle method
 
Double32_t QSquaredDA
 Q2 calculated via the double-angle method.
 
Double32_t xDA
 x calculated via the double-angle method
 
Double32_t WSquaredDA
 W2 calculated via the double-angle method.
 

Detailed Description

Abstract base class for DIS Monte Carlo events. Implements common event properties and methods.

Definition at line 30 of file EventMC.h.

Constructor & Destructor Documentation

erhic::EventMC::EventMC ( )

Constructor.

Definition at line 34 of file EventMC.cxx.

erhic::EventMC::~EventMC ( )
virtual

Destructor.

Definition at line 43 of file EventMC.cxx.

Member Function Documentation

void erhic::EventMC::AddLast ( ParticleMC track)
virtual

Add a copy of a track argument to the end of the track list.

Parameters
[in]Pointerto the track to add.

Definition at line 150 of file EventMC.cxx.

const ParticleMC * erhic::EventMC::BeamHadron ( ) const
virtual

Returns a pointer to the incident hadron beam particle. See also notes in BeamLepton().

In the standard eRHIC Monte Carlo format, the incident hadron beam is assumed to be the second particle in the particle list.

Implements erhic::EventDis.

Definition at line 126 of file EventMC.cxx.

const ParticleMC * erhic::EventMC::BeamLepton ( ) const
virtual

Returns a pointer to the incident lepton beam particle. Returns a NULL pointer if the particle cannot be located in the event. IMPORTANT - DO NOT DELETE THE POINTER OR BAD THINGS WILL HAPPEN!

In the standard eRHIC Monte Carlo format, the incident lepton beam is assumed to be the first particle in the particle list. This is the behaviour implemented here. Derived classes can implement other selection mechanisms depending on their data format.

Implements erhic::EventDis.

Definition at line 122 of file EventMC.cxx.

void erhic::EventMC::Clear ( Option_t *  = "")
virtual

Clears event contents. Event properties are reset to defaults and track list is deleted.

Definition at line 145 of file EventMC.cxx.

const ParticleMC * erhic::EventMC::ExchangeBoson ( ) const
virtual

Returns a pointer to the exchanged boson. See also notes in BeamLepton().

In the standard eRHIC Monte Carlo format, the exchanged boson is assumed to be the third particle in the particle list.

Implements erhic::EventDis.

Reimplemented in erhic::EventDjangoh, and erhic::EventPepsi.

Definition at line 130 of file EventMC.cxx.

void erhic::EventMC::FinalState ( ParticlePtrList particles) const

Stores pointers to all final state particles in the list. These pointers should not be deleted by the user. Any existing entries in the list are not changed.

Parameters
[out]particlesThe list in which to store particles.

Definition at line 73 of file EventMC.cxx.

Double_t erhic::EventMC::FinalStateCharge ( ) const

Returns the total charge of the final state in units of e.

Definition at line 106 of file EventMC.cxx.

TLorentzVector erhic::EventMC::FinalStateMomentum ( ) const

Returns the total momentum of the final state in GeV/c.

Definition at line 84 of file EventMC.cxx.

ULong64_t erhic::EventMC::GetN ( ) const
inlinevirtual

Returns a unique identifier for this event.

Definition at line 213 of file EventMC.h.

UInt_t erhic::EventMC::GetNTracks ( ) const
inlinevirtual

Returns the number of tracks in the event.

Implements erhic::VirtualEvent.

Definition at line 221 of file EventMC.h.

Int_t erhic::EventMC::GetProcess ( ) const
inlinevirtual

Returns a code describing the production process of this event.

Definition at line 217 of file EventMC.h.

const ParticleMC * erhic::EventMC::GetTrack ( UInt_t  u) const
inlinevirtual

Returns the nth track. Returns NULL if the track number is out of the range [0, GetNTracks()).

Parameters
[in]Thetrack index, in the range [0, GetNTracks()).

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Implements erhic::VirtualEvent.

Definition at line 225 of file EventMC.h.

ParticleMC * erhic::EventMC::GetTrack ( UInt_t  )
inlinevirtual

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Implements erhic::VirtualEvent.

Definition at line 233 of file EventMC.h.

TrackVector erhic::EventMC::GetTracks ( ) const

Returns pointers to all tracks in the event. Do not delete the pointers.

Definition at line 48 of file EventMC.cxx.

void erhic::EventMC::HadronicFinalState ( ParticlePtrList ) const
virtual

Yields all particles that belong to the hadronic final state. This is the same as the result of FinalState(), minus the scattered beam lepton.

Reimplemented from erhic::VirtualEvent.

Definition at line 58 of file EventMC.cxx.

TLorentzVector erhic::EventMC::HadronicFinalStateMomentum ( ) const

Returns the total momentum of the hadronic final state in GeV/c.

Definition at line 95 of file EventMC.cxx.

virtual bool erhic::EventMC::Parse ( const std::string &  )
pure virtual

Populates the event-wise variables from a string. Does not populate the particle list or compute derived quantities. See also Compute().

Implemented in erhic::EventPythia, erhic::EventGmcTrans, erhic::EventDjangoh, erhic::EventMilou, erhic::EventPepsi, erhic::EventDpmjet, and erhic::EventRapgap.

void erhic::EventMC::Reset ( )
virtual

Resets event properties to defaults. Does not clear particle list - use Clear() for that.

Definition at line 138 of file EventMC.cxx.

const ParticleMC * erhic::EventMC::ScatteredLepton ( ) const
virtual

Returns a pointer to the lepton beam particle after scattering. See also notes in BeamLepton().

In the standard eRHIC Monte Carlo format, the scattered lepton beam is assumed to be the first final-state particle in the particle list with the same PDG code as the incident lepton beam.

Implements erhic::EventDis.

Reimplemented in erhic::EventPythia, erhic::EventDjangoh, and erhic::EventPepsi.

Definition at line 134 of file EventMC.cxx.

void erhic::EventMC::SetELeptonInNuclearFrame ( double  energy)
inlinevirtual

Set incident lepton energy in the nuclear rest frame.

Definition at line 253 of file EventMC.h.

void erhic::EventMC::SetEScatteredInNuclearFrame ( double  energy)
inlinevirtual

Set scattered lepton energy in the nuclear rest frame.

Definition at line 257 of file EventMC.h.

void erhic::EventMC::SetN ( int  n)
inlinevirtual

Sets the unique identifier for this event.

Parameters
[in]nThe identifying number, an integer

Definition at line 245 of file EventMC.h.

void erhic::EventMC::SetNTracks ( int  n)
inlinevirtual

Sets the track count for this event.

Parameters
[in]nThe track count, an integer

Definition at line 249 of file EventMC.h.

void erhic::EventMC::SetProcess ( int  code)
inlinevirtual

Sets the code describing the production process of this event.

Parameters
[in]codeThe identifying code, an integer

Definition at line 241 of file EventMC.h.

Member Data Documentation

Double32_t erhic::EventMC::ELeptonInNucl
protected

Incident lepton energy in the nuclear rest frame

Definition at line 204 of file EventMC.h.

Double32_t erhic::EventMC::ELeptonOutNucl
protected

Scattered lepton energy in the nuclear rest frame

Definition at line 206 of file EventMC.h.


The documentation for this class was generated from the following files: