3 <title>The Event Record</title>
4 <link rel="stylesheet" type="text/css" href="pythia.css"/>
5 <link rel="shortcut icon" href="pythia32.gif"/>
9 <h2>The Event Record</h2>
11 A <code>Pythia</code> instance contains two members of the
12 <code>Event</code> class. The one called <code>process</code> provides
13 a brief summary of the main steps of the hard process, while the
14 one called <code>event</code> contains the full history. The
15 user would normally interact mainly with the second one, so
16 we will examplify primarily with that one.
19 The <code>Event</code> class to first approximation is a vector of
20 <code>Particle</code>s, so that it can expand to fit the current
21 event size. The index operator is overloaded, so that e.g.
22 <code>event[i]</code> corresponds to the <i>i</i>'th particle
23 of the object <code>event</code>. Thus <code>event[i].id()</code>
24 returns the identity of the <i>i</i>'th particle, and so on.
25 Therefore the methods of the
26 <code><a href="ParticleProperties.html" target="page">Particle</a></code> class
27 are at least as essential as those of the <code>Event</code> class
31 As used inside PYTHIA, some conventions are imposed on the structure
32 of the event record. Entry 0 of the <code>vector<Particle></code>
33 is used to represent the event as a whole, with its total four-momentum
34 and invariant mass, but does not form part of the event history.
35 Lines 1 and 2 contains the two incoming beams, and only from here on
36 history tracing works as could be expected. That way unassigned mother
37 and daughter indices can be put 0 without ambiguity. Depending on the
38 task at hand, a loop may therefore start at index 1 rather than 0
39 without any loss. Specifically, for translation to other event record
40 formats such as HepMC [<a href="Bibliography.html" target="page">Dob01</a>], where the first index is 1, the
41 Pythia entry 0 definitely ought to be skipped in order to minimize the
42 danger of indexing errors.
45 In the following we will list the methods available.
46 Only a few of them have a function to fill in normal user code.
48 <h3>Basic output methods</h3>
50 Some methods are available to read out information on the
53 <a name="method1"></a>
54 <p/><strong>Particle& Event::operator[](int i) </strong> <br/>
56 <strong>const Particle& Event::operator[](int i) </strong> <br/>
57 returns a (<code>const</code>) reference to the <i>i</i>'th particle
58 in the event record, which can be used to get (or set) all the
59 <a href="ParticleProperties.html" target="page">properties</a> of this particle.
62 <a name="method2"></a>
63 <p/><strong>int Event::size() </strong> <br/>
64 The event size, i.e. the sie of the <code>vector<Particle></code>.
65 Thus valid particles, to be accessed by the above indexing operator,
66 are stored in the range <i>0 <= i < size()</i>. See comment
67 above about the (ir)relevance of entry 0.
70 <a name="method3"></a>
71 <p/><strong>void Event::list() </strong> <br/>
73 <strong>void Event::list(ostream& os) </strong> <br/>
75 <strong>void Event::list(bool showScaleAndVertex, bool showMothersAndDaughters = false) </strong> <br/>
77 <strong>void Event::list(bool showScaleAndVertex, bool showMothersAndDaughters, ostream& os) </strong> <br/>
78 Provide a listing of the whole event, i.e. of the
79 <code>vector<Particle></code>. The methods with fewer arguments
80 call the final one with the respective default values, and are
81 non-inlined so they can be used in a debugger. The basic identity
82 code, status, mother, daughter, colour, four-momentum and mass data
83 are always given, but the methods can also be called with a few
84 optional arguments for further information:
85 <br/><code>argument</code><strong> showScaleAndVertex </strong> (<code>default = <strong>false</strong></code>) : optionally give a
86 second line for each particle, with the production scale (in GeV), the
87 production vertex (in mm or mm/c) and the invariant lifetime
90 <br/><code>argument</code><strong> showMothersAndDaughters </strong> (<code>default = <strong>false</strong></code>) :
91 gives a list of all daughters and mothers of a particle, as defined by
92 the <code>motherList(i)</code> and <code>daughterList(i)</code> methods
93 described below. It is mainly intended for debug purposes.
95 <br/><code>argument</code><strong> os </strong> (<code>default = <strong>cout</strong></code>) : a reference to the <code>ostream</code>
96 object to which the event listing will be directed.
102 Each <code>Particle</code> has two mother and two daughter indices.
103 These may be used to encode zero, one, two or more mothers/daughters,
104 depending on the combination of values and status code, according to
105 well-defined <a href="ParticleProperties.html" target="page">rules</a>. The
106 two methods below can do this job easier for you.
108 <a name="method4"></a>
109 <p/><strong>vector<int> Event::motherList(int i) </strong> <br/>
110 returns a vector of all the mother indices of the particle at index
111 <i>i</i>. This list is empty for entries 0, 1 and 2,
112 i.e. the "system" in line 0 is not counted as part of the history.
113 Normally the list contains one or two mothers, but it can also be more,
114 e.g. in string fragmentation the whole fragmenting system is counted
115 as mothers to the primary hadrons. Many particles may have the same
116 <code>motherList</code>. Mothers are listed in ascending order.
119 <a name="method5"></a>
120 <p/><strong>vector<int> Event::daughterList(int i) </strong> <br/>
121 returns a vector of all the daughter indices of the particle at index
122 <i>i</i>. This list is empty for a particle that did
123 not decay (or, if the evolution is stopped early enough, a parton
124 that did not branch), while otherwise it can contain a list of
125 varying length, from one to many. For the two incoming beam particles,
126 all shower initiators and beam remnants are counted as daughters,
127 with the one in slot 0 being the one leading up to the hardest
128 interaction. The "system" in line 0 does not have any daughters,
129 i.e. is not counted as part of the history. Many partons may have the
130 same <code>daughterList</code>. Daughters are listed in ascending order.
133 <a name="method6"></a>
134 <p/><strong>int Event::statusHepMC(int i) </strong> <br/>
135 returns the status code according to the HepMC conventions agreed in
136 February 2009. This convention does not preserve the full information
137 provided by the internal PYTHIA status code, as obtained by
138 <code>Particle::status()</code>, but comes reasonably close.
139 The allowed output values are:
141 <li>0 : an empty entry, with no meaningful information and therefore
142 to be skipped unconditionally (should not occur in PYTHIA);</li>
143 <li>1 : a final-state particle, i.e. a particle that is not decayed
144 further by the generator (may also include unstable particles that
145 are to be decayed later, as part of the detector simulation);</li>
146 <li>2 : a decayed Standard Model hadron or tau or mu lepton, excepting
147 virtual intermediate states thereof (i.e. the particle must undergo
148 a normal decay, not e.g. a shower branching);</li>
149 <li>3 : a documentation entry (not used in PYTHIA);</li>
150 <li>4 : an incoming beam particle;</li>
151 <li>11 - 200 : an intermediate (decayed/branched/...) particle that does
152 not fulfill the criteria of status code 2, with a generator-dependent
153 classification of its nature; in PYTHIA the absolute value of the normal
154 status code is used.</li>
159 <h3>Further output methods</h3>
161 The above methods are the main ones that a normal user would make
162 frequent use of. There are some further methods that also could come
163 in handy, in the exploration of the history of an event, but where
164 the outcome is not always obvious if one is not familiar with the
165 detailed structure of an event record.
167 <a name="method7"></a>
168 <p/><strong>int Event::iTopCopy(int i) </strong> <br/>
170 <strong>int Event::iBotCopy(int i) </strong> <br/>
171 are used to trace carbon copies of the particle at index <i>i</i> up
172 to its top mother or down to its bottom daughter. If there are no such
173 carbon copies, <i>i</i> itself will be returned. A carbon copy is
174 when the "same" particle appears several times in the event record, but
175 with changed momentum owing to recoil effects.
178 <a name="method8"></a>
179 <p/><strong>int Event::iTopCopyId(int i) </strong> <br/>
181 <strong>int Event::iBotCopyId(int i) </strong> <br/>
182 also trace top mother and bottom daughter, but do not require carbon
183 copies, only that one can find an unbroken chain, of mothers or daughters,
184 with the same flavour <code>id</code> code. When it encounters ambiguities,
185 say a <i>g -> g g</i> branching or a <i>u u -> u u</i> hard scattering,
186 it will stop the tracing and return the current position. It can be confused
187 by nontrivial flavour changes, e.g. a hard process <i>u d -> d u</i>
188 by <i>W^+-</i> exchange will give the wrong answer. These methods
189 therefore are of limited use for common particles, in particular for the
190 gluon, but should work well for "rare" particles.
193 <a name="method9"></a>
194 <p/><strong>vector<int> Event::sisterList(int i) </strong> <br/>
195 returns a vector of all the sister indices of the particle at index
196 <i>i</i>, i.e. all the daughters of the first mother, except the
200 <a name="method10"></a>
201 <p/><strong>vector<int> Event::sisterListTopBot(int i,bool widenSearch = true) </strong> <br/>
202 returns a vector of all the sister indices of the particle at index
203 <i>i</i>, tracking up and back down through carbon copies
204 if required. That is, the particle is first traced up with
205 <code>iTopCopy()</code> before its mother is found, and then all
206 the particles in the <code>daughterList()</code> of this mother are
207 traced down with <code>iBotCopy()</code>, omitting the original
208 particle itself. Any non-final particles are removed from the list.
209 Should this make the list empty the search criterion is widened so that
210 all final daughters are allowed, not only carbon-copy ones. A second
211 argument <code>false</code> inhibits the second step, and increases
212 the risk that an empty list is returned. A typical example of this
213 is for ISR cascades, e.g. <i>e -> e gamma</i> where the photon
214 may not have any obvious sister in the final state if the bottom copy
215 of the photon is an electron that annihilates and thus is not part of
219 <a name="method11"></a>
220 <p/><strong>bool Event::isAncestor(int i, int iAncestor) </strong> <br/>
221 traces the particle <i>i</i> upwards through mother, grandmother,
222 and so on, until either <i>iAncestor</i> is found or the top of
223 the record is reached. Normally one unique mother is required,
224 as is the case e.g. in decay chains or in parton showers, so that
225 e.g. the tracing through a hard scattering would not work. For
226 hadronization, first-rank hadrons are identified with the respective
227 string endpoint quark, which may be useful e.g. for <i>b</i> physics,
228 while higher-rank hadrons give <code>false</code>. Currently also
229 ministrings that collapsed to one single hadron and junction topologies
230 give <code>false</code>.
234 One data member in an <code>Event</code> object is used to keep track
235 of the largest <code>col()</code> or <code>acol()</code> colour tag set
236 so far, so that new ones do not clash.
238 <p/><code>mode </code><strong> Event:startColTag </strong>
239 (<code>default = <strong>100</strong></code>; <code>minimum = 0</code>; <code>maximum = 1000</code>)<br/>
240 This sets the initial colour tag value used, so that the first one
241 assigned is <code>startColTag + 1</code>, etc. The Les Houches accord
242 [<a href="Bibliography.html" target="page">Boo01</a>] suggests this number to be 500, but 100 works equally
246 <a name="method12"></a>
247 <p/><strong>void Event::initColTag(int colTag = 0) </strong> <br/>
248 forces the current colour tag value to be the larger of the input
249 <code>colTag</code> and the above <code>Event:startColTag</code>
253 <a name="method13"></a>
254 <p/><strong>int Event::lastColTag() </strong> <br/>
255 returns the current maximum colour tag.
258 <a name="method14"></a>
259 <p/><strong>int Event::nextColTag() </strong> <br/>
260 increases the current maximum colour tag by one and returns this
261 new value. This method is used whenever a new colour tag is needed.
265 Many event properties are accessible via the <code>Info</code> class,
266 <a href="EventInformation.html" target="page">see here</a>. Since they are used
267 directly in the event generation, a few are stored directly in the
268 <code>Event</code> class, however.
270 <a name="method15"></a>
271 <p/><strong>void Event::scale( double scaleIn) </strong> <br/>
273 <strong>double Event::scale() </strong> <br/>
274 set or get the scale (in GeV) of the hardest process in the event.
275 Matches the function of the <code>scale</code> variable in the
276 <a href="LesHouchesAccord.html" target="page">Les Houches Accord</a>.
279 <a name="method16"></a>
280 <p/><strong>void Event::scaleSecond( double scaleSecondIn) </strong> <br/>
282 <strong>double Event::scaleSecond() </strong> <br/>
283 set or get the scale (in GeV) of a second hard process in the event,
284 in those cases where such a one
285 <a href="SecondHardProcess.html" target="page">has been requested</a>.
288 <h3>Constructors and modifications of the event record</h3>
290 Although you would not normally need to create your own
291 <code>Event</code> instance, there may be times where that
292 could be convenient. The typical exampel would be if you want to
293 create a new event record that is the sum of a few different ones,
294 e.g. if you want to simulate pileup events. There may also be cases
295 where you want to add one or a few particles to an existing event
298 <a name="method17"></a>
299 <p/><strong>Event::Event(int capacity = 100) </strong> <br/>
300 creates an empty event record, but with a reserved size
301 <i>capacity</i> for the <code>Particle</code> vector.
304 <a name="method18"></a>
305 <p/><strong>Event& Event::operator=(const Event& oldEvent) </strong> <br/>
306 copies the input event record.
309 <a name="method19"></a>
310 <p/><strong>Event& Event::operator+=(const Event& addEvent) </strong> <br/>
311 appends an event to an existing one. For the appended particles
312 mother, daughter and colour tags are shifted to make a consistent
313 record. The zeroth particle of the appended event is not copied,
314 but the zeroth particle of the combined event is updated to the
315 full energy-momentum content.
318 <a name="method20"></a>
319 <p/><strong>void Event::init(string headerIn = "", ParticleData* particleDataPtrIn = 0, int startColTagIn = 100) </strong> <br/>
320 initializes colour, the pointer to the particle database, and the
321 header specification used for the event listing. We remind that a
322 <code>Pythia</code> object contains two event records
323 <code>process</code> and <code>event</code>. Thus one may e.g.
324 call either <code>pythia.process.list()</code> or
325 <code>pythia.event.list()</code>. To distinguish those two rapidly
326 at visual inspection, the <code>"Pythia Event Listing"</code> header
327 is printed out differently, in one case adding
328 <code>"(hard process)"</code> and in the other
329 <code>"(complete event)"</code>. When <code>+=</code> is used to
330 append an event, the modified event is printed with
331 <code>"(combination of several events)"</code> as a reminder.
334 <a name="method21"></a>
335 <p/><strong>void Event::clear() </strong> <br/>
336 empties event record. Specifically the <code>Particle</code> vector
337 size is reset to zero.
340 <a name="method22"></a>
341 <p/><strong>void Event::reset() </strong> <br/>
342 empties the event record, as <code>clear()</code> above, but then
343 fills the zero entry of the <code>Particle</code> vector with the
344 pseudoparticle used to represent the event as a whole. At this point
345 the pseudoparticle is not assigned any momentum or mass.
348 <a name="method23"></a>
349 <p/><strong>void Event::popBack(int n = 1) </strong> <br/>
350 removes the last <i>n</i> particle entries; must be a positive
354 <a name="method24"></a>
355 <p/><strong>int Event::append(Particle entryIn) </strong> <br/>
356 appends a particle to the bottom of the event record and
357 returns the index of this position.
360 <a name="method25"></a>
361 <p/><strong>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.) </strong> <br/>
362 appends a particle to the bottom of the event record and
363 returns the index of this position;
364 <a href="ParticleProperties.html" target="page">see here</a> for the meaning
365 of the various particle properties.
368 <a name="method26"></a>
369 <p/><strong>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.) </strong> <br/>
370 appends a particle to the bottom of the event record and
371 returns the index of this position, as above but with four-momentum
372 as a <code>Vec4</code>.
375 <a name="method27"></a>
376 <p/><strong>int Event::append(int id, int status, int col, int acol, double px, double py, double pz, double e, double m = 0.) </strong> <br/>
378 <strong>int Event::append(int id, int status, int col, int acol, Vec4 p, double m = 0.) </strong> <br/>
379 appends a particle to the bottom of the event record and
380 returns the index of this position, as above but with vanishing
381 (i.e. zero) mother and daughter indices.
384 <a name="method28"></a>
385 <p/><strong>int Event::setPDTPtr(int iSet = -1) </strong> <br/>
386 send in a pointer to the <code>ParticleData</code> database for
387 particle <code>iSet</iset>, by default the most recently appended
388 particle. Also generates a pointer to the
389 <code>ParticleDataEntry</code> object of the identity code
393 <a name="method29"></a>
394 <p/><strong>int Event::copy(int iCopy, int newStatus = 0) </strong> <br/>
395 copies the existing particle in entry <code>iCopy</code> to the
396 bottom of the event record and returns the index of this position.
397 By default, i.e. with <code>newStatus = 0</code>, everything is
398 copied precisely as it is, which means that history information
399 has to be modified further by hand to make sense. With a positive
400 <code>newStatus</code>, the new copy is set up to be the daughter of
401 the old, with status code <code>newStatus</code>, while the status
402 code of <code>iCopy</code> is negated. With a negative
403 <code>newStatus</code>, the new copy is instead set up to be the
404 mother of <code>iCopy</code>.
407 <a name="method30"></a>
408 <p/><strong>Particle& Event::back() </strong> <br/>
409 returns a reference to the last particle in the event record.
412 <a name="method31"></a>
413 <p/><strong>void Event::restorePtrs() </strong> <br/>
414 each particle in the event record has a pointer to the whole database
415 and another to the particle species itself, used to find some particle
416 properties. The latter pointer is automatically set/changed whenever
417 the particle identity is set/changed by one of the normal methods.
418 (It is the "changed" part that prompts the inclusion of a pointer
419 to the whole database.) Of course the pointer values are specific to
420 the memory locations of the current run, and so it has no sense to
421 save them if events are written to file. Should you use some
422 persistency scheme that bypasses the normal methods when the event is
423 read back in, you can use <code>restorePtrs()</code> afterwards to set
424 these pointers appropriately.
428 A few methods exist to rotate and boost events. These derive from the
429 <a href="FourVectors.html" target="page">Vec4</a> methods, and affect both the
430 momentum and the vertex (position) components of all particles.
432 <a name="method32"></a>
433 <p/><strong>void Event::rot(double theta, double phi) </strong> <br/>
434 rotate all particles in the event by this polar and azimuthal angle
435 (expressed in radians).
438 <a name="method33"></a>
439 <p/><strong>void Event::bst(double betaX, double betaY, double betaZ) </strong> <br/>
441 <strong>void Event::bst(double betaX, double betaY, double betaZ, double gamma) </strong> <br/>
443 <strong>void Event::bst(const Vec4& vec) </strong> <br/>
444 boost all particles in the event by this three-vector.
445 Optionally you may provide the <i>gamma</i> value as a fourth argument,
446 which may help avoid roundoff errors for big boosts. You may alternatively
447 supply a <code>Vec4</code> four-vector, in which case the boost vector
448 becomes <i>beta = p/E</i>.
451 <a name="method34"></a>
452 <p/><strong>void Event::rotbst(const RotBstMatrix& M) </strong> <br/>
453 rotate and boost by the combined action encoded in the
454 <code><a href="FourVectors.html" target="page">RotBstMatrix</a> M</code>.
457 <h3>The Junction Class</h3>
459 The event record also contains a vector of junctions, which often
460 is empty or else contains only a very few per event. Methods are
461 available to add further junctions or query the current junction list.
462 This is only for the expert user, however, and is not discussed
463 further here, but only the main points.
466 A junction stores the properites associated with a baryon number that
467 is fully resolved, i.e. where three different colour indices are
468 involved. There are two main applications,
470 <li>baryon beams, where at least two valence quarks are kicked out,
471 and so the motion of the baryon number is notrivial;</li>
472 <li>baryon-number violating processes, e.g. in SUSY with broken
473 <i>R</i>-parity.</li>
475 Information on junctions is set, partly in the process generation,
476 partly in the beam remnants machinery, and used by the fragmentation
477 routines, but the normal user does not have to know the details.
480 For each junction, information is stored on the kind of junction, and
481 on the three (anti)colour indices that are involved in the junction.
482 The possibilities foreseen are:
484 <li><code>kind = 1</code> : incoming colourless particle to three
485 outgoing colours (e.g. baryon beam remnant or
486 <i>neutralino -> q q q</i>);</li>
487 <li><code>kind = 2</code> : incoming colourless particle to three
488 outgoing anticolours;</li>
489 <li><code>kind = 3</code> : one incoming anticolor (stored first)
490 and two outgoing colours (e.g. antisquark decaying to quark);</li>
491 <li><code>kind = 4</code> : one incoming color (stored first) and two
492 outgoing anticolours;</li>
493 <li><code>kind = 5</code> : incoming colour octet to three colours,
494 where the incoming colour passes through unchanged and so need not
495 be bokkept here, while the incoming anticolor (stored first) and the
496 two outgoing colours are (e.g. gluino decay to three quarks);</li>
497 <li><code>kind = 6</code> : incoming colour octet to three anticolours,
498 where the incoming anticolour passes through unchanged and so need not
499 be bookkept here, while the incoming color (stored first) and the two
500 outgoing colours are.</li>
502 The odd (even) <code>kind</code> codes corresponds to a +1 (-1) change in
503 baryon number across the junction.
504 <br/><b>Warning:</b> Currently only <code>kind = 1, 2</code> are
508 The kind and colour information in the list of junctions can be set
509 or read with methods of the <code>Event</code> class, but are not of
510 common interest and so not described here.
513 A listing of current junctions can be obtained with the
514 <code>listJunctions()</code> method.
518 Separate from the event record as such, but closely tied to it is the
519 <code><a href="AdvancedUsage.html" target="page">PartonSystems</a></code> class,
520 which mainly stores the parton indices of incoming and outgoing partons,
521 classified by collision subsystem. Such information is needed to
522 interleave multiple interactions, initial-state showers and final-state
523 showers, and append beam remnants. It could also be used in other places.
524 It is intended to be accessed only by experts, such as implementors of
525 <a href="ImplementNewShowers.html" target="page">new showering models</a>.
530 <!-- Copyright (C) 2010 Torbjorn Sjostrand -->