The Event Record
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 examplify 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. 
  
int Event::size()   
The event size, i.e. the sie 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()   
  
void Event::list(ostream& os)   
  
void Event::list(bool showScaleAndVertex, bool showMothersAndDaughters = false)   
  
void Event::list(bool showScaleAndVertex, bool showMothersAndDaughters, ostream& os)   
Provide a listing of the whole event, i.e. of the 
vector<Particle>. The methods with fewer arguments 
call the final one with the respective default values, and are 
non-inlined so they can be used in a debugger. 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 = false) :  optionally give a 
second line for each particle, with the production scale (in GeV), the
production vertex (in mm or mm/c) and the invariant lifetime
(also in mm/c).
  
argument showMothersAndDaughters  (default = false) : 
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 os  (default = cout) :  a reference to the ostream
object to which the event listing will be directed.
  
  
Each Particle has two mother and two daughter indices.
These may be used to encode zero, one, two or more mothers/daughters,
depending on the combination of values and status code, according to 
well-defined rules. The 
two methods below can do this job easier for you.
vector<int> Event::motherList(int i)   
returns a vector of all the mother indices of the particle at index 
i. This list is empty for entries 0, 1 and 2,
i.e. the "system" in line 0 is not counted as part of the history. 
Normally the list contains one or two mothers, but it can also be more, 
e.g. in string fragmentation the whole fragmenting system is counted 
as mothers to the primary hadrons. Many particles may have the same
motherList. Mothers are listed in ascending order.
  
vector<int> Event::daughterList(int i)   
returns a vector of all the daughter indices of the particle at index 
i. This list is empty for a particle that did 
not decay (or, if the evolution is stopped early enough, a parton
that did not branch), while otherwise it can contain a list of 
varying length, from one to many. For the two incoming beam particles, 
all shower initiators and beam remnants are counted as daughters, 
with the one in slot 0 being the one leading up to the hardest 
interaction. The "system" in line 0 does not have any daughters, 
i.e. is not counted as part of the history. Many partons may have the 
same daughterList. Daughters are listed in ascending order.
  
int Event::statusHepMC(int i)   
returns the status code according to the HepMC conventions agreed in
February 2009. This convention does not preserve the full information
provided by the internal PYTHIA status code, as obtained by
Particle::status(), but comes reasonably close. 
The allowed output values are:
- 0 : an empty entry, with no meaningful information and therefore
to be skipped unconditionally (should not occur in PYTHIA);
  
- 1 : a final-state particle, i.e. a particle that is not decayed
further by the generator (may also include unstable particles that 
are to be decayed later, as part of the detector simulation);
  
- 2 : a decayed Standard Model hadron or tau or mu lepton, excepting 
virtual intermediate states thereof (i.e. the particle must undergo
a normal decay, not e.g. a shower branching);
 
- 3 : a documentation entry (not used in PYTHIA);
  
- 4 : an incoming beam particle;
  
- 11 - 200 : an intermediate (decayed/branched/...) particle that does 
not fulfill the criteria of status code 2, with a generator-dependent 
classification of its nature; in PYTHIA the absolute value of the normal 
status code is used.
  
  
Further output methods
The above methods are the main ones that a normal user would make
frequent use of. There are some further methods that also could come
in handy, in the exploration of the history of an event, but where 
the outcome is not always obvious if one is not familiar with the
detailed structure of an event record.
int Event::iTopCopy(int i)   
  
int Event::iBotCopy(int i)   
are used to trace carbon copies of the particle at index i up 
to its top mother or down to its bottom daughter. If there are no such 
carbon copies, i itself will be returned. A carbon copy is 
when the "same" particle appears several times in the event record, but 
with changed momentum owing to recoil effects. 
  
int Event::iTopCopyId(int i)   
  
int Event::iBotCopyId(int i)   
also trace top mother and bottom daughter, but do not require carbon 
copies, only that one can find an unbroken chain, of mothers or daughters, 
with the same flavour id code. When it encounters ambiguities,
say a g -> g g branching or a u u -> u u hard scattering,
it will stop the tracing and return the current position. It can be confused
by nontrivial flavour changes, e.g. a hard process u d -> d u  
by W^+- exchange will give the wrong answer. These methods
therefore are of limited use for common particles, in particular for the
gluon, but should work well for "rare" particles. 
  
vector<int> Event::sisterList(int i)   
returns a vector of all the sister indices of the particle at index 
i, i.e. all the daughters of the first mother, except the 
particle itself. 
  
vector<int> Event::sisterListTopBot(int i,bool widenSearch = true)   
returns a vector of all the sister indices of the particle at index 
i, tracking up and back down through carbon copies 
if required. That is, the particle is first traced up with 
iTopCopy() before its mother is found, and then all 
the particles in the daughterList() of this mother are 
traced down with iBotCopy(), omitting the original 
particle itself. Any non-final particles are removed from the list.
Should this make the list empty the search criterion is widened so that
all final daughters are allowed, not only carbon-copy ones. A second
argument false inhibits the second step, and increases 
the risk that an empty list is returned. A typical example of this
is for ISR cascades, e.g. e -> e gamma where the photon 
may not have any obvious sister in the final state if the bottom copy 
of the photon is an electron that annihilates and thus is not part of 
the final state.  
  
bool Event::isAncestor(int i, int iAncestor)   
traces the particle i upwards through mother, grandmother, 
and so on, until either iAncestor is found or the top of 
the record is reached. Normally one unique mother is required,
as is the case e.g. in decay chains or in parton showers, so that
e.g. the tracing through a hard scattering would not work. For
hadronization, first-rank hadrons are identified with the respective 
string endpoint quark, which may be useful e.g. for b physics, 
while higher-rank hadrons give false. Currently also 
ministrings that collapsed to one single hadron and junction topologies 
give false.  
  
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.
   
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.
  
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 exampel 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::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.
  
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.)   
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.)   
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.)   
  
int Event::append(int id, int status, int col, int acol, Vec4 p, double m = 0.)   
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::setPDTPtr(int iSet = -1)   
send in a pointer to the ParticleData database for 
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.
  
Particle& Event::back()   
returns a reference to the last particle in the event record.
  
void Event::restorePtrs()   
each particle in the event record has a pointer to the whole database
and another to the particle species itself, used to find some particle
properties. The latter pointer is automatically set/changed whenever 
the particle identity is set/changed by one of the normal methods. 
(It is the "changed" part that prompts the inclusion of a pointer 
to the whole database.) 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.
  
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)   
rotate and boost by the combined action encoded in the 
RotBstMatrix M.
  
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 properites associated with a baryon number that
is fully resolved, i.e. where three different colour indices are 
involved. There are two main applications,
 
- baryon beams, where at least two valence quarks are kicked out,
and so the motion of the baryon number is notrivial;
 
- 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:
kind = 1 : incoming colourless particle to three 
outgoing colours (e.g. baryon beam remnant or 
neutralino -> q q q); 
kind = 2 : incoming colourless particle to three 
outgoing anticolours; 
kind = 3 : one incoming anticolor (stored first) 
and two outgoing  colours (e.g. antisquark decaying to quark); 
kind = 4 : one incoming color (stored first) and two 
outgoing anticolours; 
kind = 5 : incoming colour octet to three colours, 
where the incoming colour passes through unchanged and so need not 
be bokkept here, while the incoming anticolor (stored first) and the 
two outgoing colours are (e.g. gluino decay to three quarks); 
kind = 6 : incoming colour octet to three anticolours, 
where the incoming anticolour passes through unchanged and so need not 
be bookkept here, while the incoming color (stored first) and the two 
outgoing colours are. 
The odd (even) kind codes corresponds to a +1 (-1) change in 
baryon number across the junction.
Warning: Currently only kind = 1, 2 are 
implemented.
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.
 
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 multiple 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.