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

Abstract base class for DIS Monte Carlo events. More...

#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 erhic::EventBeagle

Public Member Functions

 EventMC ()
 Constructor. More...
 
virtual ~EventMC ()
 Destructor. More...
 
virtual bool RequiresEaParticleFields ()
 
virtual ULong64_t GetN () const
 Returns a unique identifier for this event. More...
 
virtual Int_t GetProcess () const
 Returns a code describing the production process of this event. More...
 
virtual UInt_t GetNTracks () const
 Returns the number of tracks in the event. More...
 
virtual const ParticleMCGetTrack (UInt_t) const
 Returns the nth track. More...
 
virtual ParticleMCGetTrack (UInt_t)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. More...
 
virtual const ParticleMCBeamLepton () const
 Returns a pointer to the incident lepton beam particle. More...
 
virtual const ParticleMCBeamHadron () const
 Returns a pointer to the incident hadron beam particle. More...
 
virtual const ParticleMCExchangeBoson () const
 Returns a pointer to the exchanged boson. More...
 
virtual const ParticleMCScatteredLepton () const
 Returns a pointer to the lepton beam particle after scattering. More...
 
virtual bool Parse (const std::string &)=0
 Populates the event-wise variables from a string. More...
 
virtual void AddLast (ParticleMC *track)
 Add a copy of a track argument to the end of the track list. More...
 
virtual void Reset ()
 Resets event properties to defaults. More...
 
virtual void Clear (Option_t *="")
 Clears event contents. More...
 
virtual void SetProcess (int code)
 Sets the code describing the production process of this event. More...
 
virtual void SetN (int n)
 Sets the unique identifier for this event. More...
 
virtual void SetNTracks (int n)
 Sets the track count for this event. More...
 
virtual void SetELeptonInNuclearFrame (double energy)
 Set incident lepton energy in the nuclear rest frame. More...
 
virtual void SetEScatteredInNuclearFrame (double energy)
 Set scattered lepton energy in the nuclear rest frame. More...
 
void FinalState (ParticlePtrList &particles) const
 Stores pointers to all final state particles in the list. More...
 
void HadronicFinalState (ParticlePtrList &) const
 Yields all particles that belong to the hadronic final state. More...
 
TLorentzVector FinalStateMomentum () const
 Returns the total momentum of the final state in GeV/c. More...
 
TLorentzVector HadronicFinalStateMomentum () const
 Returns the total momentum of the hadronic final state in GeV/c. More...
 
Double_t FinalStateCharge () const
 Returns the total charge of the final state in units of e. More...
 
std::vector< const VirtualParticle * > GetTracks () const
 Returns pointers to all tracks in the event. More...
 
- Public Member Functions inherited from erhic::EventDis
virtual ~EventDis ()
 Destructor. More...
 
 EventDis ()
 Default constructor. More...
 
 EventDis (const EventDis &)
 Constructor copying another event's kinematics. More...
 
EventDisoperator= (const EventDis &)
 Assign another event's kinematics to this EventDis. More...
 
virtual Double_t GetX () const
 Returns Bjorken-x of the event. More...
 
virtual Double_t GetQ2 () const
 Returns the four-momentum transfer (exchange boson mass) Q2. More...
 
virtual Double_t GetY () const
 Returns the event inelasticity. More...
 
virtual Double_t GetYPlus () const
 Returns Y+ = y2 / (1 + (1-y)2) More...
 
virtual Double_t GetW2 () const
 Returns the invariant mass of the hadronic final state. More...
 
virtual Double_t GetNu () const
 Returns the exchange boson energy in the beam hadron rest frame. More...
 
virtual double GetXDoubleAngle () const
 Returns Bjorken x computed via the double-angle method. More...
 
virtual double GetQ2DoubleAngle () const
 Returns Q-squared computed via the double-angle method. More...
 
virtual double GetYDoubleAngle () const
 Returns inelasticity computed via the double-angle method. More...
 
virtual double GetW2DoubleAngle () const
 Returns W-squared computed via the double-angle method. More...
 
virtual double GetXJacquetBlondel () const
 Returns Bjorken x computed via the Jacquet-Blondel method. More...
 
virtual double GetQ2JacquetBlondel () const
 Returns Q-squared computed via the Jacquet-Blondel method. More...
 
virtual double GetYJacquetBlondel () const
 Returns inelasticity computed via the Jacquet-Blondel method. More...
 
virtual double GetW2JacquetBlondel () const
 Returns W-squared computed via the Jacquet-Blondel method. More...
 
virtual void SetLeptonKinematics (const DisKinematics &)
 Set the kinematics computed from the scattered lepton. More...
 
virtual void SetJacquetBlondelKinematics (const DisKinematics &)
 Set the kinematics computed from the Jacquet-Blondel method. More...
 
virtual void SetDoubleAngleKinematics (const DisKinematics &)
 Set the kinematics computed from the double-angle method. More...
 
virtual void CopyKinematics (const EventDis &)
 Copy the kinematics from another event to this event. More...
 
- Public Member Functions inherited from erhic::VirtualEvent
virtual ~VirtualEvent ()
 Destructor. More...
 

Protected Attributes

Int_t number
 Event number. More...
 
Int_t process
 PYTHIA code for the physics process producing the event. More...
 
Int_t nTracks
 Number of Particles in the event (intermediate + final) More...
 
Double32_t ELeptonInNucl
 Incident lepton energy in the nuclear rest frame. More...
 
Double32_t ELeptonOutNucl
 Scattered lepton energy in the nuclear rest frame. More...
 
TClonesArray particles
 Particle list. More...
 

Additional Inherited Members

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

Detailed Description

Abstract base class for DIS Monte Carlo events.

Implements common event properties and methods.

Definition at line 30 of file erhic/EventMC.h.

Constructor & Destructor Documentation

◆ EventMC()

erhic::EventMC::EventMC ( )

Constructor.

Definition at line 34 of file erhic/EventMC.cxx.

◆ ~EventMC()

erhic::EventMC::~EventMC ( )
virtual

Destructor.

Definition at line 43 of file erhic/EventMC.cxx.

Member Function Documentation

◆ AddLast()

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 erhic/EventMC.cxx.

◆ BeamHadron()

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 erhic/EventMC.cxx.

◆ BeamLepton()

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 erhic/EventMC.cxx.

◆ Clear()

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 erhic/EventMC.cxx.

◆ ExchangeBoson()

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 erhic/EventMC.cxx.

◆ FinalState()

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 erhic/EventMC.cxx.

◆ FinalStateCharge()

Double_t erhic::EventMC::FinalStateCharge ( ) const

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

Definition at line 106 of file erhic/EventMC.cxx.

◆ FinalStateMomentum()

TLorentzVector erhic::EventMC::FinalStateMomentum ( ) const

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

Definition at line 84 of file erhic/EventMC.cxx.

◆ GetN()

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

Returns a unique identifier for this event.

Definition at line 215 of file erhic/EventMC.h.

◆ GetNTracks()

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

Returns the number of tracks in the event.

Implements erhic::VirtualEvent.

Definition at line 223 of file erhic/EventMC.h.

◆ GetProcess()

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

Returns a code describing the production process of this event.

Definition at line 219 of file erhic/EventMC.h.

◆ GetTrack() [1/2]

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 235 of file erhic/EventMC.h.

◆ GetTrack() [2/2]

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 227 of file erhic/EventMC.h.

◆ GetTracks()

TrackVector erhic::EventMC::GetTracks ( ) const

Returns pointers to all tracks in the event.

Do not delete the pointers.

Definition at line 48 of file erhic/EventMC.cxx.

◆ HadronicFinalState()

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 erhic/EventMC.cxx.

◆ HadronicFinalStateMomentum()

TLorentzVector erhic::EventMC::HadronicFinalStateMomentum ( ) const

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

Definition at line 95 of file erhic/EventMC.cxx.

◆ Parse()

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::EventBeagle, erhic::EventPythia, erhic::EventGmcTrans, erhic::EventDjangoh, erhic::EventMilou, erhic::EventPepsi, erhic::EventDpmjet, and erhic::EventRapgap.

◆ RequiresEaParticleFields()

virtual bool erhic::EventMC::RequiresEaParticleFields ( )
inlinevirtual

Reimplemented in erhic::EventBeagle.

Definition at line 42 of file erhic/EventMC.h.

◆ Reset()

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 erhic/EventMC.cxx.

◆ ScatteredLepton()

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 erhic/EventMC.cxx.

◆ SetELeptonInNuclearFrame()

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

Set incident lepton energy in the nuclear rest frame.

Definition at line 255 of file erhic/EventMC.h.

◆ SetEScatteredInNuclearFrame()

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

Set scattered lepton energy in the nuclear rest frame.

Definition at line 259 of file erhic/EventMC.h.

◆ SetN()

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

Sets the unique identifier for this event.

Parameters
[in]nThe identifying number, an integer

Definition at line 247 of file erhic/EventMC.h.

◆ SetNTracks()

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

Sets the track count for this event.

Parameters
[in]nThe track count, an integer

Definition at line 251 of file erhic/EventMC.h.

◆ SetProcess()

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 243 of file erhic/EventMC.h.

Member Data Documentation

◆ ELeptonInNucl

Double32_t erhic::EventMC::ELeptonInNucl
protected

Incident lepton energy in the nuclear rest frame.

Definition at line 206 of file erhic/EventMC.h.

◆ ELeptonOutNucl

Double32_t erhic::EventMC::ELeptonOutNucl
protected

Scattered lepton energy in the nuclear rest frame.

Definition at line 208 of file erhic/EventMC.h.

◆ nTracks

Int_t erhic::EventMC::nTracks
protected

Number of Particles in the event (intermediate + final)

Definition at line 205 of file erhic/EventMC.h.

◆ number

Int_t erhic::EventMC::number
protected

Event number.

Definition at line 203 of file erhic/EventMC.h.

◆ particles

TClonesArray erhic::EventMC::particles
protected

Particle list.

Definition at line 210 of file erhic/EventMC.h.

◆ process

Int_t erhic::EventMC::process
protected

PYTHIA code for the physics process producing the event.

Definition at line 204 of file erhic/EventMC.h.


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