1 | <html> |
---|
2 | <head> |
---|
3 | <title>Hadron-Level Standalone</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 | <script language=javascript type=text/javascript> |
---|
10 | function stopRKey(evt) { |
---|
11 | var evt = (evt) ? evt : ((event) ? event : null); |
---|
12 | var node = (evt.target) ? evt.target :((evt.srcElement) ? evt.srcElement : null); |
---|
13 | if ((evt.keyCode == 13) && (node.type=="text")) |
---|
14 | {return false;} |
---|
15 | } |
---|
16 | |
---|
17 | document.onkeypress = stopRKey; |
---|
18 | </script> |
---|
19 | <?php |
---|
20 | if($_POST['saved'] == 1) { |
---|
21 | if($_POST['filepath'] != "files/") { |
---|
22 | echo "<font color='red'>SETTINGS SAVED TO FILE</font><br/><br/>"; } |
---|
23 | else { |
---|
24 | echo "<font color='red'>NO FILE SELECTED YET.. PLEASE DO SO </font><a href='SaveSettings.php'>HERE</a><br/><br/>"; } |
---|
25 | } |
---|
26 | ?> |
---|
27 | |
---|
28 | <form method='post' action='HadronLevelStandalone.php'> |
---|
29 | |
---|
30 | <h2>Hadron-Level Standalone</h2> |
---|
31 | |
---|
32 | The Les Houches Accord allows external process-level configurations |
---|
33 | to be fed in, for subsequent parton-level and hadron-level generation |
---|
34 | to be handled internally by PYTHIA. There is no correspondingly |
---|
35 | standardized interface if you have external events that have also |
---|
36 | been generated through the parton-level stage, so that only the |
---|
37 | hadron-level remains to be handled. A non-standard way to achieve this |
---|
38 | exists, however, and can be useful both for real applications and |
---|
39 | for various tests of the hadronization model on its own. |
---|
40 | |
---|
41 | <p/> |
---|
42 | The key trick is to set the flag <code>ProcessLevel:all = off</code>. |
---|
43 | When <code>pythia.next()</code> is called it then does not try to |
---|
44 | generate a hard process. Since there are no beams, it is also not |
---|
45 | possible to perform the normal <code>PartonLevel</code> step. |
---|
46 | (It is still possible to generate final-state radiation, but this |
---|
47 | is not automatic. It would have to be done by hand, using the |
---|
48 | <code>pythia.forceTimeShower(...)</code> method, before |
---|
49 | <code>pythia.next()</code> is called.) |
---|
50 | Thus only the <code>HadronLevel</code> methods are |
---|
51 | called, to take the current content of the event record stored in |
---|
52 | <code>pythia.event</code> as a starting point for any hadronization |
---|
53 | and decays that are allowed by the normal parameters of this step. |
---|
54 | Often the input would consist solely of partons grouped into colour |
---|
55 | singlets, but also (colour-singlet) particles are allowed. |
---|
56 | |
---|
57 | <p/> |
---|
58 | To set up all the parameters, a <code>pythia.init()</code> call has |
---|
59 | to be used, without any arguments. In brief, the structure of the |
---|
60 | main program therefore should be something like |
---|
61 | <pre> |
---|
62 | Pythia pythia; // Declare generator. |
---|
63 | Event& event = pythia.event // Convenient shorthand. |
---|
64 | pythia.readString("ProcessLevel:all = off"); // The trick! |
---|
65 | pythia.init(); // Initialization. |
---|
66 | for (int iEvent = 0; iEvent < nEvent; ++iEvent) { |
---|
67 | // Insert filling of event here! |
---|
68 | pythia.next(); // Do the hadron level. |
---|
69 | } |
---|
70 | </pre> |
---|
71 | Of course this should be supplemented by analysis of events, error checks, |
---|
72 | and so on, as for a normal PYTHIA run. The unique aspect is how to fill |
---|
73 | the <code>event</code> inside the loop, before <code>pythia.next()</code> |
---|
74 | is called. |
---|
75 | |
---|
76 | <h3>Input configuration</h3> |
---|
77 | |
---|
78 | To set up a new configuration the first step is to throw away the current |
---|
79 | one, with <code>event.reset()</code>. This routine will also reserve |
---|
80 | the zeroth entry in the even record to represent the event as a whole. |
---|
81 | |
---|
82 | <p/> |
---|
83 | With the <code>event.append(...)</code> methods a new entry is added at the |
---|
84 | bottom of the current record, i.e. the first time it is called entry |
---|
85 | number 1 is filled, and so on. The <code>append</code> method basically |
---|
86 | exists in four variants, either without or with history information, |
---|
87 | and with four-momentum provided either as a <code>Vec4</code> four-vector |
---|
88 | or as four individual components: |
---|
89 | <pre> |
---|
90 | append( id, status, col, acol, p, m) |
---|
91 | append( id, status, col, acol, px, py, pz, e, m) |
---|
92 | append( id, status, mother1, mother2, daughter1, daughter2, col, acol, p, m) |
---|
93 | append( id, status, mother1, mother2, daughter1, daughter2, col, acol, px, py, pz, e, m) |
---|
94 | </pre> |
---|
95 | The methods return the index at which the entry has been stored, |
---|
96 | but normally you would not use this feature. |
---|
97 | |
---|
98 | <p/> |
---|
99 | You can find descriptions of the input variables |
---|
100 | <?php $filepath = $_GET["filepath"]; |
---|
101 | echo "<a href='ParticleProperties.php?filepath=".$filepath."' target='page'>";?>here</a>. |
---|
102 | The PDG particle code <code>id</code> and the Les Houches Accord colour |
---|
103 | <code>col</code> and anticolour <code>acol</code> tags must be set |
---|
104 | correctly. The four-momentum and mass have to be provided in units of GeV; |
---|
105 | if you omit the mass it defaults to 0. |
---|
106 | |
---|
107 | <p/> |
---|
108 | Outgoing particles that should hadronize should be given status code 23. |
---|
109 | Often this is the only status code you need. You could e.g. also fill in |
---|
110 | incoming partons with -21 and intermediate ones with -22, if you so wish. |
---|
111 | Usually the choice of status codes is not crucial, so long as you recall |
---|
112 | that positive numbers correspond to particles that are still around, while |
---|
113 | negative numbers denote ones that already hadronized or decayed. However, |
---|
114 | so as not to run into contradictions with the internal PYTHIA checks |
---|
115 | (when <code>Check:event = on</code>), or with external formats such as |
---|
116 | HepMC, we do recommend the above codes. When <code>pythia.next()</code> |
---|
117 | is called the positive-status particles that hadronize/decay get the sign |
---|
118 | of the status code flipped to negative but the absolute value is retained. |
---|
119 | The new particles are added with normal PYTHIA status codes. |
---|
120 | |
---|
121 | <p/> |
---|
122 | For normal hadronization/decays in <code>pythia.next()</code> the |
---|
123 | history encoded in the mother and daughter indices is not used. |
---|
124 | Therefore the first two <code>append</code> methods, which set all these |
---|
125 | indices vanishing, should suffice. The subsequent hadronization/decays |
---|
126 | will still be properly documented. |
---|
127 | |
---|
128 | <p/> |
---|
129 | The exception is when you want to include junctions in your string |
---|
130 | topology, i.e. have three string pieces meet. Then you must insert in |
---|
131 | your event record the (decayed) particle that is the reason for the |
---|
132 | presence of a junction, e.g. a baryon beam remnant from which several |
---|
133 | valence quarks have been kicked out, or a neutralino that underwent a |
---|
134 | baryon-number-violating decay. This particle must have as daughters |
---|
135 | the three partons that together carry the baryon number. |
---|
136 | |
---|
137 | <p/> |
---|
138 | The sample program in <code>main21.cc</code> illustrates how you can work |
---|
139 | with this facility, both for simple parton configurations and for more |
---|
140 | complicated ones with junctions. |
---|
141 | |
---|
142 | <p/> |
---|
143 | As an alternative to setting up a topology with the methods above, |
---|
144 | a <?php $filepath = $_GET["filepath"]; |
---|
145 | echo "<a href='LesHouchesAccord.php?filepath=".$filepath."' target='page'>";?>Les Houches Event File</a> (LHEF) |
---|
146 | can also provide the configurations. Since no beams or processes are |
---|
147 | defined, only the <code><event>....</event></code> blocks |
---|
148 | need to be present, one for each event, so strictly it is not a true |
---|
149 | LHEF. You need to select <code>Beams:frameType = 4</code>, provide |
---|
150 | the file name in <code>Beams:LHEF</code> and, as above, set |
---|
151 | <code>ProcessLevel:all = off</code>. Needless to say, an externally |
---|
152 | linked <code>LHAup</code> class works as well as an LHEF, |
---|
153 | with <code>Beams:frameType = 5</code>. |
---|
154 | |
---|
155 | <p/> |
---|
156 | The event information to store in the LHEF, or provide by the |
---|
157 | <code>LHAup</code>, is essentially the same as above. The only |
---|
158 | difference is in status codes: outgoing particles should have 1 |
---|
159 | instead of 23, and intermediate resonances 2 instead of -22. |
---|
160 | Incoming partons, if any, are -1 instead of -21. |
---|
161 | |
---|
162 | <h3>Extensions to resonance decays</h3> |
---|
163 | |
---|
164 | With the above scheme, <code>pythia.next()</code> will generate |
---|
165 | hadronization, i.e. string fragmentation and subsequent decays of |
---|
166 | normal unstable particles. It will not decay |
---|
167 | <?php $filepath = $_GET["filepath"]; |
---|
168 | echo "<a href='ResonanceDecays.php?filepath=".$filepath."' target='page'>";?>resonances</a>, i.e. |
---|
169 | <i>W, Z</i>, top, Higgs, SUSY and other massive particles. |
---|
170 | |
---|
171 | <p/> |
---|
172 | If the decay products are already provided, of course those |
---|
173 | products will be hadronized, but without any of the showers |
---|
174 | that one would expect in <i>Z^0 -> q qbar</i>, say. That is, the |
---|
175 | presence of a decayed resonance in the event record can be nice |
---|
176 | for documentation purposes, but otherwise it plays no role. |
---|
177 | It is possible to change this behaviour with the following flag. |
---|
178 | |
---|
179 | <br/><br/><strong>Standalone:allowResDec</strong> <input type="radio" name="1" value="on"><strong>On</strong> |
---|
180 | <input type="radio" name="1" value="off" checked="checked"><strong>Off</strong> |
---|
181 | (<code>default = <strong>off</strong></code>)<br/> |
---|
182 | If off, then resonances are stable if final-state particles, |
---|
183 | and irrelevant if intermediate particles. If on, then |
---|
184 | resonances are decayed, sequentially if necessary. After each |
---|
185 | decay step a parton shower is applied, and subsequent decay |
---|
186 | kinematics is shifted by the recoil of such branchings. |
---|
187 | If the decay chain and kinematics of the resonance is provided |
---|
188 | at input, this information is used, otherwise it is generated |
---|
189 | internally according to built-in branching ratios etc. |
---|
190 | |
---|
191 | |
---|
192 | <p/> |
---|
193 | The input configuration has to follow the rules described above, |
---|
194 | with <code>ProcessLevel:all = off</code>. (Terminology not quite |
---|
195 | consistent, since resonance decays normally is part of the |
---|
196 | process-level step.) It is possible to combine several resonances, |
---|
197 | and other coloured or uncoloured particles into the same event. |
---|
198 | |
---|
199 | <p/> |
---|
200 | Final-state radiation is applied to each step of the decay |
---|
201 | sequence by default, but can be switched off with |
---|
202 | <code>PartonLevel:FSR = off</code> or |
---|
203 | <code>PartonLevel:FSRinResonances</code>. |
---|
204 | |
---|
205 | <h3>Repeated hadronization or decay</h3> |
---|
206 | |
---|
207 | An alternative approach is possible with the |
---|
208 | <code>pythia.forceHadronLevel()</code> routine. This method does |
---|
209 | a call to the <code>HadronLevel</code> methods, irrespective of the |
---|
210 | value of the <code>HadronLevel:all</code> flag. If you hadronize |
---|
211 | externally generated events it is equivalent to a |
---|
212 | <code>pythia.next()</code> call with |
---|
213 | <code>ProcessLevel:all = off</code>. |
---|
214 | |
---|
215 | <p/> |
---|
216 | This method truly sticks to the hadron level, and thus cannot handle |
---|
217 | resonance decays. You therefore <b>must not</b> mix it with the |
---|
218 | <code>Standalone:allowResDec = on</code> framework. |
---|
219 | |
---|
220 | <p/> |
---|
221 | The similarity of names indicates that |
---|
222 | <code>pythia.forceTimeShower( int iBeg, int iEnd, double pTmax, |
---|
223 | int nBranchMax = 0)</code> is intended to belong to the same set of |
---|
224 | work-by-hand methods. Here <code>iBeg</code> and <code>iEnd</code> |
---|
225 | give the range of partons that should be allowed to shower, |
---|
226 | <code>pTmax</code> the maximum <i>pT</i> scale of emissions, |
---|
227 | and a nonzero <code>nBranchMax</code> a maximum number of allowed |
---|
228 | branchings. Additionally, a <code>scale</code> has to be set for each |
---|
229 | parton that should shower, which requires an additional final argument |
---|
230 | to the <code>append</code> methods above. This scale limits the maximum |
---|
231 | <i>pT</i> allowed for each parton, in addition to the global |
---|
232 | <code>pTmax</code>. When not set the scale defaults to 0, meaning no |
---|
233 | radiation for that parton. |
---|
234 | |
---|
235 | <p/> |
---|
236 | The real application instead is for repeated hadronization of the same |
---|
237 | PYTHIA process- and parton-level event. This may for some studies |
---|
238 | help to save time, given that these two first step are more |
---|
239 | time-consuming than the hadronization one. |
---|
240 | |
---|
241 | <p/> |
---|
242 | For repeated hadronization you should first generate an event as usual, |
---|
243 | but with <code>HadronLevel:all = off</code>. This event you can save |
---|
244 | in a temporary copy, e.g. <code>Event savedEvent = pythia.event</code>. |
---|
245 | Inside a loop you copy back with <code>pythia.event = savedEvent</code>, |
---|
246 | and call <code>pythia.forceHadronLevel()</code> to obtain a new |
---|
247 | hadronization history. |
---|
248 | |
---|
249 | <p/> |
---|
250 | A more limited form of repetition is if you want to decay a given kind |
---|
251 | of particle repeatedly, without having to generate the rest of the event |
---|
252 | anew. This could be the case e.g. in <i>B</i> physics applications. |
---|
253 | Then you can use the <code>pythia.moreDecays()</code> method, which |
---|
254 | decays all particles in the event record that have not been decayed |
---|
255 | but should have been done so. The |
---|
256 | <code>pythia.particleData.mayDecay( id, false/true)</code> method |
---|
257 | may be used to switch off/on the decays of a particle species |
---|
258 | <code>id</code>, so that it is not decayed in the |
---|
259 | <code>pythia.next()</code> call but only inside a loop over a number of |
---|
260 | tries. |
---|
261 | |
---|
262 | <p/> |
---|
263 | Between each loop the newly produced decay products must be |
---|
264 | removed and the decayed particle status restored to undecayed. |
---|
265 | The former is simple, since the new products are appended to the |
---|
266 | end of the event record: <code>event.saveSize()</code> saves the |
---|
267 | initial size of the event record, and <code>event.restoreSize()</code> |
---|
268 | can later be used repeatedly to restore this original size, which means |
---|
269 | that the new particles at the end are thrown away. The latter is more |
---|
270 | complicated, and requires the user to identify the positions of all |
---|
271 | particles of the species and restore a positive status code with |
---|
272 | <code>event[i].statusPos()</code>. |
---|
273 | |
---|
274 | <p/> |
---|
275 | The <code>main15.cc</code> program illustrates both these methods, |
---|
276 | i.e. either repeated hadronization or repeated decay of PYTHIA |
---|
277 | events. |
---|
278 | |
---|
279 | <input type="hidden" name="saved" value="1"/> |
---|
280 | |
---|
281 | <?php |
---|
282 | echo "<input type='hidden' name='filepath' value='".$_GET["filepath"]."'/>"?> |
---|
283 | |
---|
284 | <table width="100%"><tr><td align="right"><input type="submit" value="Save Settings" /></td></tr></table> |
---|
285 | </form> |
---|
286 | |
---|
287 | <?php |
---|
288 | |
---|
289 | if($_POST["saved"] == 1) |
---|
290 | { |
---|
291 | $filepath = $_POST["filepath"]; |
---|
292 | $handle = fopen($filepath, 'a'); |
---|
293 | |
---|
294 | if($_POST["1"] != "off") |
---|
295 | { |
---|
296 | $data = "Standalone:allowResDec = ".$_POST["1"]."\n"; |
---|
297 | fwrite($handle,$data); |
---|
298 | } |
---|
299 | fclose($handle); |
---|
300 | } |
---|
301 | |
---|
302 | ?> |
---|
303 | </body> |
---|
304 | </html> |
---|
305 | |
---|
306 | <!-- Copyright (C) 2012 Torbjorn Sjostrand --> |
---|