1 <chapter name="Hadron-Level Standalone">
3 <h2>Hadron-Level Standalone</h2>
5 The Les Houches Accord allows external process-level configurations
6 to be fed in, for subsequent parton-level and hadron-level generation
7 to be handled internally by PYTHIA. There is no correspondingly
8 standardized interface if you have external events that have also
9 been generated through the parton-level stage, so that only the
10 hadron-level remains to be handled. A non-standard way to achieve this
11 exists, however, and can be useful both for real applications and
12 for various tests of the hadronization model on its own.
15 The key trick is to set the flag <code>ProcessLevel:all = off</code>.
16 When <code>pythia.next()</code> is called it then does not try to
17 generate a hard process, and therefore also cannot do anything on the
18 parton level. Instead only the <code>HadronLevel</code> methods are
19 called, to take the current content of the event record stored in
20 <code>pythia.event</code> as a starting point for any hadronization
21 and decays that are allowed by the normal parameters of this step.
22 Often the input would consist solely of partons grouped into colour
23 singlets, but also (colour-singlet) particles are allowed.
26 To set up all the parameters, a <code>pythia.init()</code> call has
27 to be used, without any arguments. In brief, the structure of the
28 main program therefore should be something like
30 Pythia pythia; // Declare generator.
31 Event& event = pythia.event // Convenient shorthand.
32 pythia.readString("ProcessLevel:all = off"); // The trick!
33 pythia.init(); // Initialization.
34 for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
35 // Insert filling of event here!
36 pythia.next(); // Do the hadron level.
39 Of course this should be supplemented by analysis of events, error checks,
40 and so on, as for a normal PYTHIA run. The unique aspect is how to fill
41 the <code>event</code> inside the loop, before <code>pythia.next()</code>
44 <h3>Input configuration</h3>
46 To set up a new configuration the first step is to throw away the current
47 one, with <code>event.reset()</code>. This routine will also reserve
48 the zeroth entry in the even record to represent the event as a whole.
51 With the <code>event.append(...)</code> methods a new entry is added at the
52 bottom of the current record, i.e. the first time it is called entry
53 number 1 is filled, and so on. The <code>append</code> method basically
54 exists in four variants, either without or with history information,
55 and with four-momentum provided either as a <code>Vec4</code> four-vector
56 or as four individual components:
58 append( id, status, col, acol, p, m)
59 append( id, status, col, acol, px, py, pz, e, m)
60 append( id, status, mother1, mother2, daughter1, daughter2, col, acol, p, m)
61 append( id, status, mother1, mother2, daughter1, daughter2, col, acol, px, py, pz, e, m)
63 The methods return the index at which the entry has been stored,
64 but normally you would not use this feature.
67 You can find descriptions of the input variables
68 <aloc href="ParticleProperties">here</aloc>.
69 The PDG particle code <code>id</code> and the Les Houches Accord colour
70 <code>col</code> and anticolour <code>acol</code> tags must be set
71 correctly. The four-momentum and mass have to be provided in units of GeV;
72 if you omit the mass it defaults to 0.
75 Outgoing particles that should hadronize should be given status code 23.
76 Often this is the only status code you need. You could e.g. also fill in
77 incoming partons with -21 and intermediate ones with -22, if you so wish.
78 Usually the choice of status codes is not crucial, so long as you recall
79 that positive numbers correspond to particles that are still around, while
80 negative numbers denote ones that already hadronized or decayed. However,
81 so as not to run into contradictions with the internal PYTHIA checks
82 (when <code>Check:event = on</code>), or with external formats such as
83 HepMC, we do recommend the above codes. When <code>pythia.next()</code>
84 is called the positive-status particles that hadronize/decay get the sign
85 of the status code flipped to negative but the absolute value is retained.
86 The new particles are added with normal PYTHIA status codes.
89 For normal hadronization/decays in <code>pythia.next()</code> the
90 history encoded in the mother and daughter indices is not used.
91 Therefore the first two <code>append</code> methods, which set all these
92 indices vanishing, should suffice. The subsequent hadronization/decays
93 will still be properly documented.
96 The exception is when you want to include junctions in your string
97 topology, i.e. have three string pieces meet. Then you must insert in
98 your event record the (decayed) particle that is the reason for the
99 presence of a junction, e.g. a baryon beam remnant from which several
100 valence quarks have been kicked out, or a neutralino that underwent a
101 baryon-number-violating decay. This particle must have as daughters
102 the three partons that together carry the baryon number.
105 The sample program in <code>main21.cc</code> illustrates how you can work
106 with this facility, both for simple parton configurations and for more
107 complicated ones with junctions.
109 <h3>Repeated hadronization or decay</h3>
111 An alternative approach is possible with the
112 <code>pythia.forceHadronLevel()</code> routine. This method does
113 a call to the <code>HadronLevel</code> methods, irrespective of the
114 value of the <code>HadronLevel:all</code> flag. If you hadronize
115 externally generated events it is equivalent to a
116 <code>pythia.next()</code> call with
117 <code>ProcessLevel:all = off</code>.
120 The real application instead is for repeated hadronization of the same
121 PYTHIA process- and parton-level event. This may for some studies
122 help to save time, given that these two first step are more
123 time-consuming than the hadronization one.
126 For repeated hadronization you should first generate an event as usual,
127 but with <code>HadronLevel:all = off</code>. This event you can save
128 in a temporary copy, e.g. <code>Event savedEvent = pythia.event</code>.
129 Inside a loop you copy back with <code>pythia.event = savedEvent</code>,
130 and call <code>pythia.forceHadronLevel()</code> to obtain a new
131 hadronization history.
134 A more limited form of repetition is if you want to decay a given kind
135 of particle repeatedly, without having to generate the rest of the event
136 anew. This could be the case e.g. in <ei>B</ei> physics applications.
137 Then you can use the <code>pythia.moreDecays()</code> method, which
138 decays all particles in the event record that have not been decayed
139 but should have been done so. The
140 <code>pythia.particleData.mayDecay( id, false/true)</code> method
141 may be used to switch off/on the decays of a particle species
142 <code>id</code>, so that it is not decayed in the
143 <code>pythia.next()</code> call but only inside a loop over a number of
147 Between each loop the newly produced decay products must be
148 removed and the decayed particle status restored to undecayed.
149 The former is simple, since the new products are appended to the
150 end of the event record: <code>event.saveSize()</code> saves the
151 initial size of the event record, and <code>event.restoreSize()</code>
152 can later be used repeatedly to restore this original size, which means
153 that the new particles at the end are thrown away. The latter is more
154 complicated, and requires the user to identify the positions of all
155 particles of the species and restore a positive status code with
156 <code>event[i].statusPos()</code>.
159 The <code>main15.cc</code> program illustrates both these methods,
160 i.e. either repeated hadronization or repeated decay of PYTHIA
165 <!-- Copyright (C) 2010 Torbjorn Sjostrand -->