1 | |
---|
2 | <chapter name="ALPGEN and MLM merging"> |
---|
3 | |
---|
4 | <h2>ALPGEN and MLM merging</h2> |
---|
5 | |
---|
6 | This manual page describes the ALPGEN <ref>Man03</ref> and MLM merging |
---|
7 | <ref>Man02</ref> interfaces for PYTHIA 8. While future versions of |
---|
8 | ALPGEN will be able to write out events in LHEF format, previous |
---|
9 | versions always output events in an ALPGEN native format (a combination |
---|
10 | of a ".unw" and a "_unw.par" file). The ALPGEN component of this code |
---|
11 | contains a reader for this native format (for unweighted events), as |
---|
12 | well as parameter reading for both ALPGEN native and LHE file formats. |
---|
13 | Although designed to work together with the MLM component, they are |
---|
14 | implemented entirely separately and it is possible to use one without |
---|
15 | the other. |
---|
16 | |
---|
17 | <p/> |
---|
18 | It should be noted that all the functionality described here is provided |
---|
19 | through external routines, and therefore the presence of these features is |
---|
20 | dependent on the main program being used. This structure allows for the |
---|
21 | easy extensibility of the merging scheme. The files of interest are located |
---|
22 | in the <code>examples/</code> subdirectory: |
---|
23 | <ul> |
---|
24 | <li> |
---|
25 | <code>MLMhooks.h</code> : provides a |
---|
26 | <aloc href="UserHooks"><code>UserHooks</code></aloc> derived |
---|
27 | class for performing the MLM merging procedure. All <code>MLM:*</code> |
---|
28 | options, described below, are implemented in this file. Further |
---|
29 | technical details of this class are given at the end of this manual |
---|
30 | page. |
---|
31 | </li> |
---|
32 | <li> |
---|
33 | <code>LHAupAlpgen.h</code> : provides three classes for the reading of |
---|
34 | ALPGEN event and parameter files. <code>LHAupAlpgen</code> is an |
---|
35 | <aloc href="LesHouchesAccord"><code>LHAup</code></aloc> derived |
---|
36 | class for reading in ALPGEN native format event files. |
---|
37 | <code>AlpgenPar</code> is a class for the parsing of ALPGEN parameter |
---|
38 | files, making the information available through a simple interface. |
---|
39 | <code>AlpgenHooks</code> is a |
---|
40 | <aloc href="UserHooks"><code>UserHooks</code></aloc> derived class that |
---|
41 | provides the <code>Alpgen:*</code> options, described below. Further |
---|
42 | technical details of these classes are given at the end of this manual |
---|
43 | page. |
---|
44 | </li> |
---|
45 | <li> |
---|
46 | <code>main32.cc, main32.cmnd</code> : a sample main program and card |
---|
47 | file showing the usage of the previous two files. It combines the Alpgen |
---|
48 | and MLM UserHooks classes together, such that the functionality of both |
---|
49 | is available, and reads in a sample ALPGEN event file while performing |
---|
50 | the MLM merging procedure. Some commented out sets of options are |
---|
51 | provided in the card file, which can be activated to try different |
---|
52 | merging setups. |
---|
53 | </li> |
---|
54 | <li> |
---|
55 | <code>main32.unw, main32_unw.par</code> : an ALPGEN format event and |
---|
56 | parameter file containing 100 W + 3 jet events. It is not feasible |
---|
57 | to package large event files with the PYTHIA distribution, but this |
---|
58 | sample is enough to show the different components in action. |
---|
59 | </li> |
---|
60 | </ul> |
---|
61 | |
---|
62 | <h2>ALPGEN</h2> |
---|
63 | |
---|
64 | <b>NB: these options are provided by the AlpgenHooks class, |
---|
65 | which must be loaded for this functionality to be present</b> |
---|
66 | |
---|
67 | <p/> |
---|
68 | ALPGEN event files that have been written out in LHEF format should be |
---|
69 | read in through the normal LHEF machinery |
---|
70 | (see <aloc href="BeamParameters">beam parameters</aloc>). Files in |
---|
71 | ALPGEN's native format, instead, may be processed using the |
---|
72 | <code>Alpgen:file</code> option below. When using this option, the |
---|
73 | ALPGEN parameter file is stored in the PYTHIA Info object under the key |
---|
74 | <code>AlpgenPar</code>, see the "Header information" section of the |
---|
75 | <aloc href="EventInformation">Event Information</aloc> manual page for |
---|
76 | more details. Processes not implemented by the PYTHIA 6 interface |
---|
77 | supplied with ALPGEN are also not implemented here. |
---|
78 | |
---|
79 | <p/> |
---|
80 | When reading in ALPGEN native event files, some momenta are shifted by |
---|
81 | the file reader to ensure energy-momentum conservation. The magnitude of |
---|
82 | these shifts should be small (around the MeV level in the worst case) |
---|
83 | and warnings will be produced if they are above a set threshold. A large |
---|
84 | number of warnings may signify unexpected behaviour and should |
---|
85 | potentially be investigated. It is also known that certain event |
---|
86 | classes, for example an event with both light and heavy <ei>b</ei> |
---|
87 | quarks may give rise to these warnings. |
---|
88 | |
---|
89 | <p/> |
---|
90 | The ALPGEN file reader supports the reading of the event and parameter |
---|
91 | files in gzip format with file extensions ".unw.gz" and "_unw.par.gz" |
---|
92 | respectively. This requires the use of external libraries, however, and |
---|
93 | the <code>README</code> file in the main directory contains instructions |
---|
94 | on how to enable this. |
---|
95 | |
---|
96 | <p/> |
---|
97 | All other <code>Alpgen:*</code> options apply to both LHE and native |
---|
98 | file formats, and include options to guide the MLM merging procedure |
---|
99 | based on the parameters that are read in with the events file. |
---|
100 | |
---|
101 | <word name="Alpgen:file" default="void"> |
---|
102 | This option is used to read in ALPGEN format event files. Using this option |
---|
103 | overrides any previously set beam options inside PYTHIA. The path to the |
---|
104 | files, not including any file extension, should be provided e.g. for input |
---|
105 | files <ei>input_unw.par</ei> and <ei>input.unw</ei>, the value |
---|
106 | <ei>input</ei> should be used. |
---|
107 | </word> |
---|
108 | |
---|
109 | <flag name="Alpgen:setMasses" default="on"> |
---|
110 | When switched on, any particle masses provided by ALPGEN are set in |
---|
111 | the PYTHIA <aloc href="ParticleDataScheme">particle database</aloc>. |
---|
112 | </flag> |
---|
113 | |
---|
114 | <flag name="Alpgen:setMLM" default="on"> |
---|
115 | When switched on, the merging parameters (see below) are set according to |
---|
116 | the ALPGEN hard process cuts: |
---|
117 | <ul> |
---|
118 | <li> <ei>MLM:eTjetMin = min(ptjmin + 5., 1.2 * ptjmin), </ei> </li> |
---|
119 | <li> <ei>MLM:coneRadius = drjmin, </ei> |
---|
120 | <li> <ei>MLM:etaJetMax = etajmax, </ei> |
---|
121 | </ul> |
---|
122 | where the <code>ptjmin</code>, <code>drjmin</code> and |
---|
123 | <code>etajmax</code> are the incoming ALPGEN parameters. Note that any |
---|
124 | existing values of these parameters are overwritten. |
---|
125 | </flag> |
---|
126 | |
---|
127 | <flag name="Alpgen:setNjet" default="on"> |
---|
128 | When switched on, the <code>MLM:nJet</code> parameter (see below) is set |
---|
129 | to the incoming <code>njet</code> ALPGEN parameter. Note that any |
---|
130 | existing value of this parameter is overwritten. |
---|
131 | </flag> |
---|
132 | |
---|
133 | |
---|
134 | <h2>MLM Merging</h2> |
---|
135 | |
---|
136 | <b>NB: these options are provided by the MLMhooks class, |
---|
137 | which must be loaded for this functionality to be present</b> |
---|
138 | |
---|
139 | <p/> |
---|
140 | This section describes the MLM merging interface for PYTHIA 8. The most |
---|
141 | common reference to the algorithm is <ref>Man02</ref>. In many respects, |
---|
142 | however, the implementation provided in the ALPGEN package should be |
---|
143 | considered the official description of the MLM merging procedure. |
---|
144 | Although designed primarily to work with events generated with ALPGEN, it |
---|
145 | can in principle also be used with events from a different source. This |
---|
146 | should not be done without thought, however, and it is up to the user to |
---|
147 | understand the details of the algorithm and the implications of using a |
---|
148 | different hard process generator. |
---|
149 | |
---|
150 | <p/> |
---|
151 | A brief description of the MLM procedure, as implemented, is given here. |
---|
152 | First, either the CellJet or SlowJet jet algorithm is chosen, see the |
---|
153 | <aloc href="EventAnalysis">Event Analysis</aloc> page for more details. |
---|
154 | Both of these algorithms have an <ei>R</ei> and <ei>etaMax</ei> |
---|
155 | parameter. In addition, CellJet has an <ei>eTmin</ei> and SlowJet has a |
---|
156 | <ei>pTmin</ei> parameter. These are the primary three parameters of the |
---|
157 | merging procedure, and in practice are set dependent on the cuts applied |
---|
158 | to the matrix element (ME) generation. We stress that the merging |
---|
159 | procedure is not tied to the geometry of a specific physical detector, |
---|
160 | but only to the match between the original partons and the resulting |
---|
161 | jets, using standard jet algorithms in the phase space region where |
---|
162 | partons have been generated. |
---|
163 | |
---|
164 | <p/> |
---|
165 | ME samples with different jet multiplicities are run through the event |
---|
166 | generator, and the generation interrupted after parton showers have been |
---|
167 | applied, but before resonance decays and beam remnants have been |
---|
168 | processed. Note in particular that top quarks will not yet |
---|
169 | be decayed, which may lead to slight differences with the PYTHIA 6 |
---|
170 | interface included with the ALPGEN package. In what follows, the |
---|
171 | hardness measure of jets/partons is taken to be <ei>eT</ei> when CellJet |
---|
172 | is used and <ei>pT</ei> when SlowJet is used. The hard system (ignoring |
---|
173 | all MPI systems) is then analysed: |
---|
174 | <ul> |
---|
175 | <li> |
---|
176 | The particles in the original matrix element process are sorted into |
---|
177 | light partons, heavy partons and other particles. For backwards |
---|
178 | compatibility, a light parton is defined as the set <ei>(d, u, s, c, |
---|
179 | b, g)</ei> with zero mass. A heavy parton is defined as the set |
---|
180 | <ei>(c, b, t)</ei> with non-zero mass. |
---|
181 | </li> |
---|
182 | <li> |
---|
183 | All particles not originating from the heavy partons or other |
---|
184 | particles are passed to the jet algorithm and clustered. |
---|
185 | </li> |
---|
186 | <li> |
---|
187 | Clustered jets are matched to the light partons in the original ME |
---|
188 | process. There are two different methods which can be used: |
---|
189 | <ul> |
---|
190 | <li> |
---|
191 | Method 1: The following is done for each parton, in order |
---|
192 | of decreasing hardness. The <ei>delta R</ei> between the parton |
---|
193 | and all jets is calculated and the smallest value taken. If |
---|
194 | this is less than the jet <ei>R</ei> parameter, possibly |
---|
195 | multiplied by a constant, the jet and parton are considered to |
---|
196 | match, and the jet is removed from further consideration. |
---|
197 | Note that for CellJet the <ei>delta R</ei> measure is in |
---|
198 | <ei>(eta, phi)</ei>, while for SlowJet, it is in |
---|
199 | <ei>(y, phi)</ei>. |
---|
200 | </li> |
---|
201 | <li> |
---|
202 | Method 2: This method is only possible when using the SlowJet |
---|
203 | algorithm. Before the clustering is performed, extremely soft |
---|
204 | "ghost" particles are added to the event at the |
---|
205 | <ei>(y, phi)</ei> coordinates of the original matrix element |
---|
206 | partons. If such a particle is clustered into a jet, the parton |
---|
207 | and jet are considered to match. The idea of "ghost" particles |
---|
208 | was originally introduced by FastJet as a way to measure jet |
---|
209 | areas <ref>Cac06</ref> and should not affect clustering with an |
---|
210 | infrared-safe jet algorithm. |
---|
211 | </li> |
---|
212 | </ul> |
---|
213 | <li> |
---|
214 | If there is a light ME parton remaining which has not been matched |
---|
215 | to a jet, then the event is vetoed. If all ME partons have been |
---|
216 | matched to a jet, but there are still some extra jets remaining, |
---|
217 | then two options are possible: |
---|
218 | <ul> |
---|
219 | <li> |
---|
220 | Exclusive mode: the event is vetoed. This is typically used when |
---|
221 | there are ME samples with higher jet multiplicities, which would |
---|
222 | fill in the extra jets. |
---|
223 | </li> |
---|
224 | <li> |
---|
225 | Inclusive mode: the event is retained if the extra jets are softer |
---|
226 | than the softest matched jet. This is typically used when |
---|
227 | there is no ME sample with higher jet multiplicity, so the parton |
---|
228 | shower should be allowed to give extra jets. |
---|
229 | </ul> |
---|
230 | </li> |
---|
231 | <li> |
---|
232 | All particles originating from the heavy partons are passed to the |
---|
233 | jet algorithm and clustered. |
---|
234 | </li> |
---|
235 | <li> |
---|
236 | The clustered jets are again matched to the original partons, but |
---|
237 | there is no requirement for a match to be present; all matched jets |
---|
238 | are immediately discarded. The matching procedure is much the same |
---|
239 | as for light partons, but with two differences when <ei>delta R</ei> |
---|
240 | matching is used. First, a different <ei>R</ei> parameter than that |
---|
241 | used by the jet algorithm may optionally be given. Second, all jets |
---|
242 | that are within the given radius of the parton are matched, not |
---|
243 | just the one with the smallest <ei>delta R</ei> measure. If there |
---|
244 | are still extra jets remaining then in exclusive mode the event is |
---|
245 | immediately vetoed, while in inclusive mode the event is retained if |
---|
246 | the extra jets are softer than the softest <em>light</em> matched jet. |
---|
247 | </li> |
---|
248 | </ul> |
---|
249 | |
---|
250 | <p/> |
---|
251 | Some different options are provided, specified further below. These |
---|
252 | are set so that, by default, the algorithm closely follows the official |
---|
253 | MLM interface provided in the ALPGEN package. All vetoing of events |
---|
254 | is done through the usual |
---|
255 | <aloc href="UserHooks"><code>UserHooks</code></aloc> machinery, and is |
---|
256 | therefore already taken into account in the cross section. |
---|
257 | In the output from |
---|
258 | <code><aloc href="EventStatistics">Pythia::statistics()</aloc></code>, |
---|
259 | the difference between the "Selected" and "Accepted" columns gives the |
---|
260 | number of events that have not survived the vetoing procedure. It is |
---|
261 | still the responsibility of the user to add together the results from |
---|
262 | runs with different jet multiplicities. In the simplest case, when |
---|
263 | ALPGEN input is used and the hard process parameters are used to guide |
---|
264 | the merging procedure, it is enough to set the <code>MLM:nJetMax</code> |
---|
265 | parameter (see below). |
---|
266 | |
---|
267 | |
---|
268 | <h3>Merging</h3> |
---|
269 | |
---|
270 | <flag name="MLM:merge" default="off"> |
---|
271 | Master switch to activate MLM merging. |
---|
272 | </flag> |
---|
273 | |
---|
274 | |
---|
275 | <h3>Exclusive mode</h3> |
---|
276 | |
---|
277 | The following settings determine whether clustered jets which do not |
---|
278 | match an original hard parton are allowed. They are typically permitted |
---|
279 | in the highest jet multiplicity sample, where the parton shower may |
---|
280 | produce extra hard jets, without risk of double counting. Any |
---|
281 | extra jet produced by the shower must be softer than any matched light |
---|
282 | jet, or else the event is vetoed. |
---|
283 | |
---|
284 | <modepick name="MLM:exclusive" default="2" min="0" max="2"> |
---|
285 | Exclusive or inclusive merging. |
---|
286 | <option value="0"> |
---|
287 | The merging is run in inclusive mode. |
---|
288 | </option> |
---|
289 | <option value="1"> |
---|
290 | The merging is run in exclusive mode. |
---|
291 | </option> |
---|
292 | <option value="2"> |
---|
293 | If <ei>MLM:nJet < MLM:nJetMax</ei>, then the merging is run in |
---|
294 | exclusive mode, otherwise it is run in inclusive mode. If |
---|
295 | <ei>MLM:nJetMax < 0</ei> or <ei>MLM:nJet < 0</ei>, then the |
---|
296 | algorithm defaults to exclusive mode. |
---|
297 | </option> |
---|
298 | |
---|
299 | <modepick name="MLM:nJet" default="-1" min="-1"> |
---|
300 | The number of additional light jets in the incoming process. |
---|
301 | This is used when <ei>MLM:exclusive = 2</ei>. |
---|
302 | Note that for ALPGEN input, this value may be automatically set. |
---|
303 | </modepick> |
---|
304 | |
---|
305 | <modepick name="MLM:nJetMax" default="-1" min="-1"> |
---|
306 | This parameter is used to indicate the maximum number of jets that will be |
---|
307 | matched. This is used when <ei>MLM:exclusive = 2</ei>. |
---|
308 | </modepick> |
---|
309 | |
---|
310 | |
---|
311 | <h3>Jet algorithm</h3> |
---|
312 | |
---|
313 | The choice of jet algorithm and associated parameters can be adjusted with |
---|
314 | the settings below. The PYTHIA 8 internal <code>CellJet</code> and |
---|
315 | <code>SlowJet</code> routines are used, see the |
---|
316 | <aloc href="EventAnalysis">Event Analysis</aloc> page for more details. |
---|
317 | |
---|
318 | <modepick name="MLM:jetAlgorithm" default="1" min="1" max="2"> |
---|
319 | The choice of jet algorithm to use when merging against hard partons. |
---|
320 | <option value="1">The CellJet algorithm.</option> |
---|
321 | <option value="2">The SlowJet algorithm.</option> |
---|
322 | </modepick> |
---|
323 | |
---|
324 | <modepick name="MLM:nEta" default="100" min="50"> |
---|
325 | Specific to the CellJet algorithm, the number of bins in pseudorapidity. |
---|
326 | </modepick> |
---|
327 | |
---|
328 | <modepick name="MLM:nPhi" default="60" min="30"> |
---|
329 | Specific to the CellJet algorithm, the number of bins in phi. |
---|
330 | </modepick> |
---|
331 | |
---|
332 | <parm name="MLM:eTseed" default="1.5" min="0.0"> |
---|
333 | Specific to the CellJet algorithm, the minimum <ei>eT</ei> for a cell |
---|
334 | to be acceptable as the trial center of a jet. |
---|
335 | </parm> |
---|
336 | |
---|
337 | <parm name="MLM:eTthreshold" default="0.1"> |
---|
338 | Specific to the CellJet algorithm, cells with <ei>eT < eTthreshold</ei> |
---|
339 | are completely neglected by the jet algorithm. |
---|
340 | </parm> |
---|
341 | |
---|
342 | <modepick name="MLM:slowJetPower" default="-1" min="-1" max="1"> |
---|
343 | The power to use in the SlowJet algorithm. |
---|
344 | <option value="-1">The anti-kT algorithm.</option> |
---|
345 | <option value="0">The Cambridge/Aachen algorithm.</option> |
---|
346 | <option value="1">The kT algorithm.</option> |
---|
347 | </modepick> |
---|
348 | |
---|
349 | |
---|
350 | <h3>Merging parameters</h3> |
---|
351 | |
---|
352 | The following options are the three main parameters for the merging |
---|
353 | procedure. Although here they are in principle free parameters, they should |
---|
354 | be heavily influenced by the hard process generation cuts. For ALPGEN |
---|
355 | input, these values may be set automatically. |
---|
356 | |
---|
357 | <parm name="MLM:eTjetMin" default="20.0" min="5.0"> |
---|
358 | For the CellJet algorithm, this gives the minimum transverse energy |
---|
359 | inside a cone for a jet to be accepted. For the SlowJet algorithm, this |
---|
360 | is instead the minimum transverse momentum required for a cluster to be |
---|
361 | accepted as a jet. |
---|
362 | </parm> |
---|
363 | |
---|
364 | <parm name="MLM:coneRadius" default="0.7" min="0.1"> |
---|
365 | For the CellJet algorithm, this gives the size of the cone in (eta, phi) |
---|
366 | space drawn around the geometric center of the jet. For the SlowJet |
---|
367 | algorithm, this gives the <ei>R</ei> parameter. |
---|
368 | </parm> |
---|
369 | |
---|
370 | <parm name="MLM:etaJetMax" default="2.5" min="0.1"> |
---|
371 | For both jet algorithms, this defines the maximum pseudorapidity that |
---|
372 | the detector is assumed to cover. In this context, however, it is tied |
---|
373 | to the phase space region in which partons have been generated. |
---|
374 | Specifically, particles within <ei>etaJetMax + coneRadius</ei> are |
---|
375 | passed to the jet algorithm, with only jets within <ei>etaJetMax</ei> |
---|
376 | retained in the merging. |
---|
377 | </parm> |
---|
378 | |
---|
379 | |
---|
380 | <h3>Jet matching</h3> |
---|
381 | |
---|
382 | The following parameters control the criteria for matching a clustered jet |
---|
383 | to a hard parton. |
---|
384 | |
---|
385 | <modepick name="MLM:jetAllow" default="1" min="1" max="2"> |
---|
386 | Controls which particles are clustered by the jet algorithm. |
---|
387 | <option value="1"> |
---|
388 | This option explicitly disallows top quarks, leptons and photons. All other |
---|
389 | particle types are passed to the jet algorithm. |
---|
390 | </option> |
---|
391 | <option value="2"> |
---|
392 | No extra particles are disallowed. |
---|
393 | </option> |
---|
394 | </modepick> |
---|
395 | |
---|
396 | <modepick name="MLM:jetMatch" default="1" min="1" max="2"> |
---|
397 | Criteria for matching a clustered jet to a parton. |
---|
398 | <option value="1"> |
---|
399 | This option can be used with both the CellJet and SlowJet |
---|
400 | algorithms. The <ei>delta R</ei> between each parton and |
---|
401 | jet is taken, and the minimal value compared against |
---|
402 | <ei>MLM:coneMatchLight * MLM:coneRadius</ei> for light jets or |
---|
403 | <ei>MLM:coneMatchHeavy * MLM:coneRadiusHeavy</ei> for heavy jets. |
---|
404 | Note that by default <ei>MLM:coneRadiusHeavy = MLM:coneRadius</ei>, |
---|
405 | see below. If below this value, the parton and jet are considered |
---|
406 | to match. With CellJet, the <ei>delta R</ei> measure is in |
---|
407 | <ei>(eta, phi)</ei>, while with SlowJet it is in <ei>(y, phi)</ei>. |
---|
408 | </option> |
---|
409 | <option value="2"> |
---|
410 | This option can only be used with the SlowJet algorithm. The hard |
---|
411 | partons are inserted into the parton level event as "ghost" particles, |
---|
412 | but at the correct <ei>(y, phi)</ei> position. If this particle is then |
---|
413 | clustered into a jet, it is considered a match. |
---|
414 | </option> |
---|
415 | </modepick> |
---|
416 | |
---|
417 | <parm name="MLM:coneMatchLight" default="1.5" min="0.1"> |
---|
418 | The <code>coneMatchLight</code> parameter used when |
---|
419 | <ei>MLM:jetMatch = 1</ei>. |
---|
420 | </parm> |
---|
421 | |
---|
422 | <parm name="MLM:coneRadiusHeavy" default="-1.0"> |
---|
423 | The <code>coneRadiusHeavy</code> parameter used when |
---|
424 | <ei>MLM:jetMatch = 1</ei>. When assigned a negative value, |
---|
425 | the value of <code>MLM:coneRadius</code> is used. |
---|
426 | </parm> |
---|
427 | |
---|
428 | <parm name="MLM:coneMatchHeavy" default="1.0" min="0.1"> |
---|
429 | The <code>coneMatchHeavy</code> parameter used when |
---|
430 | <ei>MLM:jetMatch = 1</ei>. |
---|
431 | </parm> |
---|
432 | |
---|
433 | |
---|
434 | <h2>Class information</h2> |
---|
435 | |
---|
436 | Some more technical information about the different classes is given |
---|
437 | below. For clarity, some limited information on certain private methods |
---|
438 | is provided. |
---|
439 | |
---|
440 | <h3>LHAupAlpgen</h3> |
---|
441 | |
---|
442 | This class is derived from the |
---|
443 | <aloc href="LesHouchesAccord"><code>LHAup</code></aloc> base class, and |
---|
444 | uses the standard machinery to pass initialisation and event data to |
---|
445 | PYTHIA. These standard functions are not documented here. The complete |
---|
446 | parameter file is stored in the PYTHIA Info object, if given, under the |
---|
447 | key <code>AlpgenPar</code>. |
---|
448 | |
---|
449 | <method name="LHAupAlpgen::LHAupAlpgen(const char *baseFNin, |
---|
450 | Info *infoPtrIn = NULL)"> |
---|
451 | The constructor for the class takes the base filename for the ALPGEN |
---|
452 | format files (without file extensions) and optionally a pointer to a |
---|
453 | PYTHIA Info class, used for warning/error message printing and for |
---|
454 | storing the ALPGEN parameter file. The event and |
---|
455 | parameter files are opened immediately, with the <code>AlpgenPar</code> |
---|
456 | class, described below, used to parse the parameter file. |
---|
457 | </method> |
---|
458 | |
---|
459 | <method name="bool LHAupAlpgen::addResonances()"> |
---|
460 | This is a private method used when an event is read in. The information |
---|
461 | read from the event file does not always contain a complete listing of |
---|
462 | all particles and four-momenta, and so various details must be |
---|
463 | reconstructed. Exactly which details are filled in can vary based on the |
---|
464 | ALPGEN process in question. |
---|
465 | </method> |
---|
466 | |
---|
467 | <method name="bool LHAupAlpgen::rescaleMomenta()"> |
---|
468 | This is another private method used when an event is read in. |
---|
469 | It shuffles and rescales momenta in an event to ensure energy-momentum |
---|
470 | conservation. First, <ei>pT</ei> is made to balance by splitting any |
---|
471 | imbalance between all outgoing particles with their energies also |
---|
472 | scaled. Second, the <ei>e/pZ</ei> of the two incoming particles are |
---|
473 | scaled to balance the outgoing particles. Finally, any intermediate |
---|
474 | resonances are recalculated from their decay products. |
---|
475 | </method> |
---|
476 | |
---|
477 | <h3>AlpgenPar</h3> |
---|
478 | |
---|
479 | This class parses an ALPGEN parameter file and makes the information |
---|
480 | available through a simple interface. The information is stored |
---|
481 | internally in key/value (string/double) format. All lines prior to: |
---|
482 | <pre> ************** run parameters </pre> |
---|
483 | are ignored, and in the general case, a line e.g. |
---|
484 | <pre> 10 3.00000000000000 ! njets</pre> |
---|
485 | would be stored with key "njets" and value "3.0". The following lines |
---|
486 | are special cases where the line may be split or the key translated: |
---|
487 | <pre> |
---|
488 | 3 ! hard process code |
---|
489 | 0.000 4.700 174.300 80.419 91.188 120.000 ! mc,mb,mt,mw,mz,mh |
---|
490 | 912.905 0.0914176 ! Crosssection +- error (pb) |
---|
491 | 100 29787.4 ! unwtd events, lum (pb-1) Njob= 2 |
---|
492 | </pre> |
---|
493 | In the first line, the key "hard process code" is translated to |
---|
494 | "hpc". In the second, the mass values are split and each given an entry |
---|
495 | in the internal store. In the third, the cross section and cross section |
---|
496 | error are stored under the keys "xsecup" and "xerrup" respectively. |
---|
497 | Finally, the number of events and luminosity are stored under the keys |
---|
498 | "nevent" and "lum" respectively. In the event that a duplicate key is |
---|
499 | present, with differing values, the stored value is overwritten and a |
---|
500 | warning given. |
---|
501 | |
---|
502 | <method name="AlpgenPar::AlpgenPar(Info *infoPtrIn = NULL)"> |
---|
503 | The constructor does nothing except for store the PYTHIA Info |
---|
504 | pointer, if given. This is used for warning/error message printing. |
---|
505 | </method> |
---|
506 | |
---|
507 | <method name="bool AlpgenPar::parse(const string paramStr)"> |
---|
508 | This method parses an ALPGEN parameter file. The parameter file is |
---|
509 | passed as a single string, mainly intended to be read out from the |
---|
510 | PYTHIA Info object using the header information methods. |
---|
511 | </method> |
---|
512 | |
---|
513 | <method name="bool AlpgenPar::haveParam(const string &paramIn)"> |
---|
514 | Method to check if a parameter with key <code>paramIn</code> is present. |
---|
515 | Returns true if present, else false. |
---|
516 | </method> |
---|
517 | |
---|
518 | <method name="double AlpgenPar::getParam(const string &paramIn)"> |
---|
519 | </method> |
---|
520 | <methodmore name="int AlpgenPar::getParamAsInt(const string &paramIn)"> |
---|
521 | Return the parameter with key <code>paramIn</code> as a double or |
---|
522 | integer. The presence of a parameter should have already been checked |
---|
523 | using the <code>haveParam()</code> function above. If the parameter is |
---|
524 | not present, 0 is returned. |
---|
525 | </methodmore> |
---|
526 | |
---|
527 | <method name="void AlpgenPar::void printParams()"> |
---|
528 | Method to print a list of stored parameters. |
---|
529 | </method> |
---|
530 | |
---|
531 | <h3>AlpgenHooks</h3> |
---|
532 | |
---|
533 | This <aloc href="UserHooks"><code>UserHooks</code></aloc> derived class |
---|
534 | provides all the <code>Alpgen:*</code> options. It is provided as a |
---|
535 | UserHooks class such that the code works regardless of whether ALPGEN |
---|
536 | native or LHE file formats are used. It is declared with virtual |
---|
537 | inheritance so that it may be combine with other UserHooks classes, see |
---|
538 | the "Combining UserHooks" section below. |
---|
539 | |
---|
540 | <method name="AlpgenHooks(Pythia &pythia)"> |
---|
541 | The constructor takes a PYTHIA object as input, so that the beam |
---|
542 | parameter settings can be overridden if the <code>Alpgen:file</code> |
---|
543 | option is given. If this is the case, an <code>LHAupAlpgen</code> |
---|
544 | instance is automatically created and passed to PYTHIA. |
---|
545 | </method> |
---|
546 | |
---|
547 | <method name="bool initAfterBeams()"> |
---|
548 | This is the only UserHooks method that is overridden. It is called |
---|
549 | directly after PYTHIA has initialised the beams, and therefore the |
---|
550 | header information should be present in the PYTHIA Info object. The |
---|
551 | <code>AlpgenPar</code> class is used to parse ALPGEN parameters, if |
---|
552 | present, which are then used to set further PYTHIA settings. |
---|
553 | </method> |
---|
554 | |
---|
555 | <h3>MLMhooks</h3> |
---|
556 | |
---|
557 | This <aloc href="UserHooks"><code>UserHooks</code></aloc> derived class |
---|
558 | provides all the <code>MLM:*</code> options. It is also declared with |
---|
559 | virtual inheritance for combining with other UserHooks classes. It |
---|
560 | uses standard UserHooks methods to perform the merging procedure |
---|
561 | outlined previously in this manual page, and therefore detailed method |
---|
562 | information is not given. |
---|
563 | |
---|
564 | <h3>Combining UserHooks</h3> |
---|
565 | |
---|
566 | It is possible to combine multiple UserHooks classes together, such that |
---|
567 | the functionality of both is present. A prerequisite is that the |
---|
568 | different UserHooks classes should be declared with virtual inheritance, |
---|
569 | e.g. |
---|
570 | <pre> |
---|
571 | class MLMhooks : virtual public UserHooks |
---|
572 | </pre> |
---|
573 | Without this option, when combining two UserHooks derived classes |
---|
574 | together, two copies of the base UserHooks class would be created |
---|
575 | leading to ambiguity. |
---|
576 | |
---|
577 | <p/> |
---|
578 | An example combining the AlpgenHooks and MLMhooks classes together is |
---|
579 | given in <code>main32.cc</code>: |
---|
580 | <pre> |
---|
581 | class AlpgenAndMLMhooks : public AlpgenHooks, public MLMhooks { |
---|
582 | public: |
---|
583 | AlpgenAndMLMhooks(Pythia &pythia) |
---|
584 | : AlpgenHooks(pythia), MLMhooks() {} |
---|
585 | |
---|
586 | virtual bool initAfterBeams() { |
---|
587 | // Call init of both AlpgenHooks and MLMhooks (order important) |
---|
588 | if (!AlpgenHooks::initAfterBeams()) return false; |
---|
589 | if (!MLMhooks::initAfterBeams()) return false; |
---|
590 | return true; |
---|
591 | } |
---|
592 | }; |
---|
593 | </pre> |
---|
594 | This class inherits from both AlpgenHooks and MLMhooks. Any functions |
---|
595 | which are present in both classes should be overridden with a function |
---|
596 | that calls the different parent methods in the desired order. In the |
---|
597 | above example, the only shared methods are the constructor and |
---|
598 | <code>initAfterBeams()</code>. |
---|
599 | |
---|
600 | </chapter> |
---|
601 | |
---|
602 | <!-- Copyright (C) 2012 Torbjorn Sjostrand --> |
---|