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