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
|
Go to the documentation of this file.
10 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
11 #define INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
16 #include <TClonesArray.h>
17 #include <TLorentzVector.h>
47 virtual ULong64_t
GetN()
const;
117 virtual bool Parse(
const std::string&) = 0;
129 virtual void Reset();
136 virtual void Clear(Option_t* =
"");
148 virtual void SetN(
int n);
200 std::vector<const VirtualParticle*>
GetTracks()
const;
228 if (u < (UInt_t)
particles.GetEntries()) {
236 if (u < (UInt_t)
particles.GetEntries()) {
271 explicit Reader(
const std::string& treeName =
"EICTree");
313 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTMC_H_
virtual ~Reader()
Destructor.
virtual bool Parse(const std::string &)=0
Populates the event-wise variables from a string.
virtual void SetProcess(int code)
Sets the code describing the production process of this event.
virtual const ParticleMC * ExchangeBoson() const
Returns a pointer to the exchanged boson.
std::vector< const erhic::VirtualParticle * > ParticlePtrList
typedef for a track pointer collection.
virtual ULong64_t GetN() const
Returns a unique identifier for this event.
EventMC * mEvent
The last event read.
virtual void SetN(int n)
Sets the unique identifier for this event.
Int_t process
PYTHIA code for the physics process producing the event.
virtual void Clear(Option_t *="")
Clears event contents.
TLorentzVector FinalStateMomentum() const
Returns the total momentum of the final state in GeV/c.
EventMC * Read(Long64_t number)
Read and return the numbered event from the tree.
virtual Int_t GetProcess() const
Returns a code describing the production process of this event.
TTree * mTree
The tree being read.
EventMC * operator()(Long64_t number)
Read and return the numbered event from the tree.
virtual const ParticleMC * ScatteredLepton() const
Returns a pointer to the lepton beam particle after scattering.
virtual bool RequiresEaParticleFields()
void FinalState(ParticlePtrList &particles) const
Stores pointers to all final state particles in the list.
virtual void SetNTracks(int n)
Sets the track count for this event.
virtual const ParticleMC * BeamHadron() const
Returns a pointer to the incident hadron beam particle.
Double32_t ELeptonOutNucl
Scattered lepton energy in the nuclear rest frame.
Double32_t ELeptonInNucl
Incident lepton energy in the nuclear rest frame.
virtual const ParticleMC * BeamLepton() const
Returns a pointer to the incident lepton beam particle.
A deeply inelastic scattering event.
Int_t number
Event number.
virtual void SetELeptonInNuclearFrame(double energy)
Set incident lepton energy in the nuclear rest frame.
TTree * GetTree()
Return the tree read by this reader.
Reader(const std::string &treeName="EICTree")
Construct a Reader for the named TTree.
Wrapper for getting tree from file and event from tree.
virtual void AddLast(ParticleMC *track)
Add a copy of a track argument to the end of the track list.
Int_t nTracks
Number of Particles in the event (intermediate + final)
virtual const ParticleMC * GetTrack(UInt_t) const
Returns the nth track.
void HadronicFinalState(ParticlePtrList &) const
Yields all particles that belong to the hadronic final state.
TClonesArray particles
Particle list.
virtual void SetEScatteredInNuclearFrame(double energy)
Set scattered lepton energy in the nuclear rest frame.
TLorentzVector HadronicFinalStateMomentum() const
Returns the total momentum of the hadronic final state in GeV/c.
Double_t FinalStateCharge() const
Returns the total charge of the final state in units of e.
std::vector< const VirtualParticle * > GetTracks() const
Returns pointers to all tracks in the event.
virtual UInt_t GetNTracks() const
Returns the number of tracks in the event.
Abstract base class for DIS Monte Carlo events.
virtual ~EventMC()
Destructor.
virtual void Reset()
Resets event properties to defaults.