StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
StarGenEvent Class Reference

Base class for event records. More...

#include <StarGenEvent.h>

Inheritance diagram for StarGenEvent:
StarGenAAEvent StarGenEPEvent StarGenPPEvent

Public Types

enum  FilterResult { kUnknown = 0, kAccept = 0x0001, kReject = 0x0002, kFlag = 0x0004 }
 Filter result enumeration.
 

Public Member Functions

 StarGenEvent (const Char_t *name="event", const Char_t *title="")
 Constructor.
 
 ~StarGenEvent ()
 Destructor.
 
void SetGeneratorId (Int_t id)
 Set the generator ID.
 
Int_t GetGeneratorId ()
 Returns the generator ID.
 
void SetProcessId (Int_t id)
 Set the process ID for this event.
 
Int_t GetProcessId ()
 Get the process ID for this event.
 
void SetOffset (Int_t o)
 Set the particle index offset for this event.
 
Int_t GetOffset ()
 Get the particle index offset for this event.
 
Int_t operator++ ()
 Prefix increment operator increases event number.
 
Int_t operator++ (Int_t)
 Postfix increment operator increases event number.
 
Int_t GetEventNumber ()
 Get the event number.
 
void SetRunNumber (Int_t run)
 Set the run number for this event.
 
Int_t GetRunNumber ()
 Get the run number for this event.
 
void SetDaqRunNumber (Int_t run)
 Set the run number for the DAQ file in an embedding job.
 
Int_t GetDaqRunNumber ()
 Get the run number for the DAQ file in this event.
 
void SetDaqFileNumber (Int_t fnum)
 Set the file number in an embedding job.
 
Int_t GetDaqFileNumber ()
 Get the file number in an embedding job.
 
void SetBlue (Int_t id)
 Set the blue beam ID.
 
Int_t GetBlue ()
 Get the blue beam ID.
 
void SetYell (Int_t id)
 Set the yellow beam ID.
 
Int_t GetYell ()
 Get the yellow beam ID.
 
void SetRootS (Double_t rs)
 Set the sqrt(s) of the collision.
 
Double_t GetRootS ()
 Get the sqrt(s) of the collision.
 
StarGenParticleAddParticle ()
 
StarGenParticleAddParticle (Int_t status, Int_t pdg, Int_t m1, Int_t m2, Int_t d1, Int_t d2, Double_t px, Double_t py, Double_t pz, Double_t E, Double_t M, Double_t vx, Double_t vy, Double_t vz, Double_t vt)
 
StarGenParticleAddParticle (StarGenParticle *particle)
 
void Print (const Option_t *opts="head") const
 
TIter IterAll (Bool_t dir=kIterForward)
 
virtual void Clear (const Option_t *opts="part,data")
 Clear the event.
 
StarGenParticleoperator[] (Int_t idx)
 
Int_t GetNumberOfParticles ()
 Obtain the number of particles in the event record.
 
void SetFilterResult (UInt_t result)
 Sets the filter result.
 
void AddFilterResult (UInt_t result)
 
UInt_t GetFilterResult ()
 Returns the filter result.
 
void AddUserWeight (double w)
 Attach a user-weight to the event.
 
const std::vector< double > & GetUserWeights ()
 Retrieve the user weights vector.
 
const std::vector< double > & GetUserWeights () const
 Retrieve the user weights vector.
 

Protected Member Functions

void InitArrays ()
 Initialize clowns arrays.
 
 ClassDef (StarGenEvent, 1)
 

Protected Attributes

TString mName
 
TString mTitle
 
TClonesArray * mParticles
 Array of particles.
 
Int_t mGeneratorId
 Generator Id.
 
Int_t mProcessId
 Event generator process ID.
 
Int_t mOffset
 Event generator offset.
 
Int_t mEventNumber
 Event number.
 
Int_t mRunNumber
 Monte Carlo run number.
 
Int_t mDaqRunNumber
 DAQ run number (for embedding)
 
Int_t mDaqFileNumber
 File number (for embedding)
 
Int_t mBlueId
 PDG for blue beam.
 
Int_t mYellId
 PDG for yellow beam.
 
Double_t mCmsEnergy
 aka sqrt(s)
 
Int_t mNumRejected [3]
 0=total, 1=EG, 2=filter
 
UInt_t mFilterResult
 Result of filter.
 
std::vector< Double_t > mWeights
 User weights.
 
Int_t mNumParticles
 Number of particles in the record.
 

Detailed Description

Base class for event records.

StarGenEvent is the base class for event records in the STAR framework. Concrete event generators will generally use one of the derived classes which represent pp, ep, AA, or eA collisions:

The base class contains common functionality across all typical event generators in STAR:

Derived classes save additional information revelant to the kinematics of the system being studied, such as (in a pp collision)...

StarGenEvent, and its derived classes, are intended to allow the developer to provide a persistent record of the state of an event generator on each event. It provides a record of the primary particles generated, and a record of the final state particles which are pushed to the GEANT stack for simulation.

Description

The particle record is based on the HEPEVT and HepMC standards. Generated particles are created as instances of the StarGenParticle class and added to a list of particles. The mother-daughter relationships of these particles is encoded as integers which refer to the index of the mother (daughter) particle within the list. Developers should take care in implementing their event generators to make sure that the status-code conventions defined in StarGenParticle are understood and respected. The decision whether a particle is stable and should be stacked out for geant simulationn is handled by the framework... but it is based on that convention.

The particle record maintains two lists of particles. The complete event record, listing all particles which have been generated, and the list of particles which have been stacked out to the geant simulation. It is the responsability of the developer to fill in the full list of particles. The framework will take care of the decision of whether to stack the particle out and fill in the reference list of particles whcih have been fully simulated.

Header inforamtion is kept simple, as different event generators will have a different set of variables which a developer will want to be filled. So here we provide only the event number, run number, generator and process ID. Additional variables to store kinematic information are implemented in the derived classes.

Usage

StarPrimaryMaker is setup to save the event record in a TTree. The primary maker aggregates the particles generated by all event generators in the chain, adding them to an instance of StarGenEvent.

Author
Jason C Webb

Definition at line 81 of file StarGenEvent.h.

Member Function Documentation

void StarGenEvent::AddFilterResult ( UInt_t  result)
inline

Performs a bitwise-or between the current filter result and the passed parameter

Definition at line 203 of file StarGenEvent.h.

References mFilterResult.

StarGenParticle * StarGenEvent::AddParticle ( )

Add a particle to the list of particles. StarGenEvent is responsible for cleaning up the memory.

Definition at line 56 of file StarGenEvent.cxx.

References mNumParticles, mParticles, and StarGenParticle::SetIndex().

Referenced by StarKinematics::AddParticle(), and AddParticle().

StarGenParticle * StarGenEvent::AddParticle ( Int_t  status,
Int_t  pdg,
Int_t  m1,
Int_t  m2,
Int_t  d1,
Int_t  d2,
Double_t  px,
Double_t  py,
Double_t  pz,
Double_t  E,
Double_t  M,
Double_t  vx,
Double_t  vy,
Double_t  vz,
Double_t  vt 
)

Add a particle to the list of particles. StarGenEvent is responsible for cleaning up the memory.

Parameters
statusthe status code of the particle using the HepMC standard
pdgthe particle datagroup ID of the particle
m1the first mother in the event record
m2the last mother in the event record. Note: m1 to m2 must be continuous.
d1the first daughter in the event record
d2the last daughter in the event record. Note: d1 to d2 must be continuous.
pxthe x-component of the momentum
pythe y-component of the momentum
pzthe z-component of the momentum
Ethe energy
Mthe mass
vxthe x-component of the production vertex
vythe y-component of the production vertex
vzthe z-component of the production vertex
tofthe time of flight to the production vertex

Definition at line 63 of file StarGenEvent.cxx.

References AddParticle(), StarGenParticle::SetEnergy(), StarGenParticle::SetFirstDaughter(), StarGenParticle::SetFirstMother(), StarGenParticle::SetId(), StarGenParticle::SetLastDaughter(), StarGenParticle::SetLastMother(), StarGenParticle::SetMass(), StarGenParticle::SetPx(), StarGenParticle::SetPy(), StarGenParticle::SetPz(), StarGenParticle::SetStatus(), StarGenParticle::SetTof(), StarGenParticle::SetVx(), StarGenParticle::SetVy(), and StarGenParticle::SetVz().

StarGenParticle * StarGenEvent::AddParticle ( StarGenParticle particle)
TIter StarGenEvent::IterAll ( Bool_t  dir = kIterForward)
inline

Obtain an iterator over all generated particles.

Parameters
dirkIterForward or kIterBackward

Definition at line 173 of file StarGenEvent.h.

References mParticles.

Referenced by FcsDYFilter::Filter(), FcsJetFilter::Filter(), FcsJPsiFilter::Filter(), and StDijetFilter::Filter().

StarGenParticle* StarGenEvent::operator[] ( Int_t  idx)
inline

Obtain a pointer to the particle at index idx

Parameters
idxThe index of the particle. Users are cautioned that no bounds checking is applied. Do not reference idx larger than the number of particles.

Definition at line 182 of file StarGenEvent.h.

void StarGenEvent::Print ( const Option_t *  opts = "head") const

Print the event.

Parameters
optsDefault prints everything. If "simu" option is provided, then only the particles which were stacked for simulation are printed.

Definition at line 98 of file StarGenEvent.cxx.

References mEventNumber, mFilterResult, mGeneratorId, mOffset, mParticles, mRunNumber, StarGenParticle::Print(), and StarGenParticle::Simulate().

Referenced by StarGenEventReader::Generate(), and StarPrimaryMaker::Make().


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