The Event Record

  1. Basic output methods
  2. Further output methods
  3. Constructors and modifications of the event record
  4. The Junction Class
  5. The HVcols Class
  6. Subsystems
A Pythia instance contains two members of the Event class. The one called process provides a brief summary of the main steps of the hard process, while the one called event contains the full history. The user would normally interact mainly with the second one, so we will exemplify primarily with that one.

The Event class to first approximation is a vector of Particles, so that it can expand to fit the current event size. The index operator is overloaded, so that e.g. event[i] corresponds to the i'th particle of the object event. Thus event[i].id() returns the identity of the i'th particle, and so on. Therefore the methods of the Particle class are at least as essential as those of the Event class itself.

As used inside PYTHIA, some conventions are imposed on the structure of the event record. Entry 0 of the vector<Particle> is used to represent the event as a whole, with its total four-momentum and invariant mass, but does not form part of the event history. Lines 1 and 2 contains the two incoming beams, and only from here on history tracing works as could be expected. That way unassigned mother and daughter indices can be put 0 without ambiguity. Depending on the task at hand, a loop may therefore start at index 1 rather than 0 without any loss. Specifically, for translation to other event record formats such as HepMC [Dob01], where the first index is 1, the Pythia entry 0 definitely ought to be skipped in order to minimize the danger of indexing errors.

In the following we will list the methods available. Only a few of them have a function to fill in normal user code.

Basic output methods

Some methods are available to read out information on the current event record:

Particle& Event::operator[](int i)  
const Particle& Event::operator[](int i)  
returns a (const) reference to the i'th particle in the event record, which can be used to get (or set) all the properties of this particle. const vector<Particle>* particles()  
alternatively, this method directly returns the (const) vector of particles stored in the event.

Particle& Event::front()  
Particle& Event::at(int i)  
Particle& Event::back()  
returns a reference to the zeroth, i'th or last particle in the event record, as an alternative to the methods above.

int Event::size()  
The event size, i.e. the size of the vector<Particle>. Thus valid particles, to be accessed by the above indexing operator, are stored in the range 0 <= i < size(). See comment above about the (ir)relevance of entry 0.

void Event::list(bool showScaleAndVertex = false, bool showMothersAndDaughters = false, int precision = 3)  
Provide a listing of the whole event, i.e. of the vector<Particle>. The basic identity code, status, mother, daughter, colour, four-momentum and mass data are always given, but the methods can also be called with a few optional arguments for further information:
argument showScaleAndVertex (default = off) : optionally give a second line for each particle, with the production scale (in GeV), the particle polarization (dimensionless), the production vertex (in mm or mm/c) and the invariant lifetime (also in mm/c).
argument showMothersAndDaughters (default = off) : gives a list of all daughters and mothers of a particle, as defined by the motherList(i) and daughterList(i) methods described below. It is mainly intended for debug purposes.
argument precision (default = 3) : the number of digits to the right of the decimal point shown for momenta, energies andf masses. Can be set above 3, but reducing it below 3 will have no effect. This option is intended for expert users, e.g. for debugging purposes, and so no effort has been made to stretch header and footer to match.

Further output methods

Many event properties are accessible via the Info class, see here. Since they are used directly in the event generation, a few are stored directly in the Event class, however.

void Event::scale( double scaleIn)  
double Event::scale()  
set or get the scale (in GeV) of the hardest process in the event. Matches the function of the scale variable in the Les Houches Accord.

void Event::scaleSecond( double scaleSecondIn)  
double Event::scaleSecond()  
set or get the scale (in GeV) of a second hard process in the event, in those cases where such a one has been requested.

One data member in an Event object is used to keep track of the largest col() or acol() colour tag set so far, so that new ones do not clash.

mode  Event:startColTag   (default = 100; minimum = 0; maximum = 1000)
This sets the initial colour tag value used, so that the first one assigned is startColTag + 1, etc. The Les Houches accord [Boo01] suggests this number to be 500, but 100 works equally well.

void Event::initColTag(int colTag = 0)  
forces the current colour tag value to be the larger of the input colTag and the above Event:startColTag values.

int Event::lastColTag()  
returns the current maximum colour tag.

int Event::nextColTag()  
increases the current maximum colour tag by one and returns this new value. This method is used whenever a new colour tag is needed.

Constructors and modifications of the event record

Although you would not normally need to create your own Event instance, there may be times where that could be convenient. The typical example would be if you want to create a new event record that is the sum of a few different ones, e.g. if you want to simulate pileup events. There may also be cases where you want to add one or a few particles to an existing event record.

Event::Event(int capacity = 100)  
creates an empty event record, but with a reserved size capacity for the Particle vector.

Event& Event::operator=(const Event& oldEvent)  
copies the input event record.

Event& Event::operator+=(const Event& addEvent)  
appends an event to an existing one. For the appended particles mother, daughter and colour tags are shifted to make a consistent record. The zeroth particle of the appended event is not copied, but the zeroth particle of the combined event is updated to the full energy-momentum content.

void Event::init(string headerIn = "", ParticleData* particleDataPtrIn = 0, int startColTagIn = 100)  
initializes colour, the pointer to the particle database, and the header specification used for the event listing. We remind that a Pythia object contains two event records process and event. Thus one may e.g. call either pythia.process.list() or pythia.event.list(). To distinguish those two rapidly at visual inspection, the "Pythia Event Listing" header is printed out differently, in one case adding "(hard process)" and in the other "(complete event)". When += is used to append an event, the modified event is printed with "(combination of several events)" as a reminder.

void Event::clear()  
empties event record. Specifically the Particle vector size is reset to zero.

void Event::free()  
empties event record, like clear() above, but also frees up the memory of the Particle vector.

void Event::reset()  
empties the event record, as clear() above, but then fills the zero entry of the Particle vector with the pseudoparticle used to represent the event as a whole. At this point the pseudoparticle is not assigned any momentum or mass.

void Event::popBack(int n = 1)  
removes the last n particle entries; must be a positive number. History (and other) information of remaning entries is untouched, and so may be internally inconsistent.

void Event::remove(int iFirst, int iLast, bool shiftHistory = true)  
removes particles in the range between indices iFirst and iLast, including the endpoints. By default all mother and daughter indices above the removed range are shifted down by the number of removed entries, while indices in the removed range are put zero. Optionally these shifts can be omitted. Other information remains unchanged, which may lead to inconsistencies. If the decay products of a particle are removed, e.g., the mother particle status should be set positive, cf. Particle::undoDecay().

int Event::append(Particle entryIn)  
appends a particle to the bottom of the event record and returns the index of this position.

int Event::append(int id, int status, int mother1, int mother2, int daughter1, int daughter2, int col, int acol, double px, double py, double pz, double e, double m = 0., double scale = 0., double pol = 9.)  
appends a particle to the bottom of the event record and returns the index of this position; see here for the meaning of the various particle properties.

int Event::append(int id, int status, int mother1, int mother2, int daughter1, int daughter2, int col, int acol, Vec4 p, double m = 0., double scale = 0., double pol = 9.)  
appends a particle to the bottom of the event record and returns the index of this position, as above but with four-momentum as a Vec4.

int Event::append(int id, int status, int col, int acol, double px, double py, double pz, double e, double m = 0., double scale = 0., double pol = 9.)  
int Event::append(int id, int status, int col, int acol, Vec4 p, double m = 0., double scale = 0., double pol = 9.)  
appends a particle to the bottom of the event record and returns the index of this position, as above but with vanishing (i.e. zero) mother and daughter indices.

int Event::setEvtPtr(int iSet = -1)  
send in the this pointer of the current Event itself to the particle iSet, by default the most recently appended particle. Also generates a pointer to the ParticleDataEntry object of the identity code of the particle.

int Event::copy(int iCopy, int newStatus = 0)  
copies the existing particle in entry iCopy to the bottom of the event record and returns the index of this position. By default, i.e. with newStatus = 0, everything is copied precisely as it is, which means that history information has to be modified further by hand to make sense. With a positive newStatus, the new copy is set up to be the daughter of the old, with status code newStatus, while the status code of iCopy is negated. With a negative newStatus, the new copy is instead set up to be the mother of iCopy. An attempt to copy an out-of-range entry will return -1.

void Event::restorePtrs()  
each particle in the event record has a pointer to the event itself and another to the particle species it belongs to. The latter pointer is automatically set/changed whenever the particle identity is set/changed by one of the normal methods. Of course the pointer values are specific to the memory locations of the current run, and so it has no sense to save them if events are written to file. Should you use some persistency scheme that bypasses the normal methods when the event is read back in, you can use restorePtrs() afterwards to set these pointers appropriately.

int Event::nFinal(bool chargedOnly = false)  
return the number of final-state particles in the current event record, optionally counting charged particles only.

double Event::dyAbs(int i1, int i2)  
double Event::detaAbs(int i1, int i2)  
double Event::dphiAbs(int i1, int i2)  
double Event::RRapPhi(int i1, int i2)  
double Event::REtaPhi(int i1, int i2)  
return the separation between two particles in the event record, in true rapidity, in pseudorapidity, in phi angle, and in the R combination of true or pseudorapidity with phi, when required with absolute sign so as to avoid negative numbers.

A few methods exist to rotate and boost events. These derive from the Vec4 methods, and affect both the momentum and the vertex (position) components of all particles.

void Event::rot(double theta, double phi)  
rotate all particles in the event by this polar and azimuthal angle (expressed in radians).

void Event::bst(double betaX, double betaY, double betaZ)  
void Event::bst(double betaX, double betaY, double betaZ, double gamma)  
void Event::bst(const Vec4& vec)  
boost all particles in the event by this three-vector. Optionally you may provide the gamma value as a fourth argument, which may help avoid roundoff errors for big boosts. You may alternatively supply a Vec4 four-vector, in which case the boost vector becomes beta = p/E.

void Event::rotbst(const RotBstMatrix& M, bool boostVertices = true)  
rotate and boost by the combined action encoded in the RotBstMatrix M. If the optional second argument is false only the four-momenta are boosted, and not the production vertices.

The Junction Class

The event record also contains a vector of junctions, which often is empty or else contains only a very few per event. Methods are available to add further junctions or query the current junction list. This is only for the expert user, however, and is not discussed further here, but only the main points.

A junction stores the properties associated with a baryon number that is fully resolved, i.e. where three different colour indices are involved. There are two main applications,

  1. baryon beams, where at least two valence quarks are kicked out, and so the motion of the baryon number is nontrivial;
  2. baryon-number violating processes, e.g. in SUSY with broken R-parity.
Information on junctions is set, partly in the process generation, partly in the beam remnants machinery, and used by the fragmentation routines, but the normal user does not have to know the details.

For each junction, information is stored on the kind of junction, and on the three (anti)colour indices that are involved in the junction. The possibilities foreseen are:

The odd (even) kind codes corresponds to a +1 (-1) change in baryon number across the junction.

The kind and colour information in the list of junctions can be set or read with methods of the Event class, but are not of common interest and so not described here.

A listing of current junctions can be obtained with the listJunctions() method.

The HVcols Class

In Hidden Valley SU(N) scenarios, HV colour and anticolour tags need to be bookkept. So as not to have to expand the information on all particles by two more integers, a separate bookkeeping is set up only for HV-coloured particles. It is based on the HVcols class, which only stores three integers: the index of a HV-coloured particle in the full event record, and the related HV colour and anticolour tags. A vector of such objects contains the list of all HV-coloured particles in the full event, and would thus normally be empty. This information can be accessed for each particle by the pointer back to the full event, as described on the Particle Properties page. In addition some public methods exist for the event, notably

bool Event::hasHVcols()  
tell whether the event has any HV-coloured particles or not.

void Event::listHVcols()  
list the indices of particles that have HV colour, along with their respective HV colour and HV anticolour tags.

There are also some further methods mainly for internal use

int Event::maxHVcols()  
the maximum HV-colour tag used in the current event.

void Event::saveVcolsSize()  
void Event::restoreHVcolsSize()  
save the current number of HV-coloured particles, such that the vector of HV-coloured particles can be restored to this size, in case recent additions need to be undone.

Subsystems

Separate from the event record as such, but closely tied to it is the PartonSystems class, which mainly stores the parton indices of incoming and outgoing partons, classified by collision subsystem. Such information is needed to interleave multiparton interactions, initial-state showers and final-state showers, and append beam remnants. It could also be used in other places. It is intended to be accessed only by experts, such as implementors of new showering models.