1 <chapter name="Accessing Pythia6 Processes">
3 <h2>Access PYTHIA 6 Processes</h2>
5 Gradually all relevant processes from PYTHIA 6 are being
6 re-implemented in PYTHIA 8, but this transition is not finished.
7 For a while longer it may therefore at times be convenient to
8 access the Fortran PYTHIA 6 process library. In order to give
9 this access at runtime, and not only by writing/reading event files,
10 an interface is provided to C++. This interface is residing in
11 <code>Pythia6Interface.h</code>, and in addition the PYTHIA 6 library
12 must be linked. The latter should normally be the most recent
13 Fortran PYTHIA version, but must be at least 6.314, since this is
14 the first version that allows processes to be output in the Les
15 Houches format (and not only input).
18 The routines interfaced are
20 <li><code>pygive(command)</code> to give in a command to change a
21 setting, e.g. to decide which processes should be generated,
22 <li><code>pyinit( frame, beam, target, wIn)</code> to initialize
23 the event generation chain,</li>
24 <li><code>pyupin()</code> to fill this initialization information
25 in the <code>HEPRUP</code> commonblock,</li>
26 <li><code>pyupev()</code> to generate the next hard process and
27 fill the event information in the <code>HEPEUP</code> commonblock,</li>
28 <li><code>pylist( mode)</code> to list the event at the process
30 <li><code>pystat( mode)</code> to print statistics on the event
31 generation process.</li>
33 Details on allowed arguments are given in the PYTHIA 6.4 manual.
36 These methods can be used in context of the
37 <aloc href="LesHouchesAccord"><code>LHAupFortran</code></aloc> class.
38 The existing code there takes care of converting <code>HEPRUP</code>
39 and <code>HEPEUP</code> commonblock information from Fortran to C++,
40 and of making it available to the PYTHIA 8 methods. What needs to be
41 supplied are the two <code>LHAupFortran::fillHepRup()</code> and
42 <code>LHAupFortran::fillHepEup()</code> methods. The first can
43 contain an arbitrary number of <code>pygive(...)</code>, followed
44 by <code>pyinit(...)</code> and <code>pyupin()</code> in that order.
45 The second only needs to contain <code>pyupev()</code>. Finally,
46 the use of <code>pylist(...)</code> and <code>pystat(...)</code>
47 is entirely optional, and calls are best put directly in the
51 This means that all other Fortran routines have not been interfaced
52 and cannot be accessed directly from the C++ code; there is no need
53 for them in the current setup.
56 PYTHIA 6.4 does its own rejection of events internally, according to
57 the strategy option 3. However, the current runtime interface does not
58 take cross-section information from PYTHIA 6.4. This means that both
59 the initial maxima and the final cross sections printed by the PYTHIA 8
60 routines are irrelevant in this case. Instead you have to study the
61 standard PYTHIA 6.4 initialization printout and call on
62 <code>pystat(...)</code> at the end of the run to obtain this information.
63 It also means that you cannot mix with internally generated events,
64 unlike what could have been allowed for strategy 3. Should a strong need
65 arise, PYTHIA 6.4 could be modified to work with strategy option 1 and
66 thus allow a mixing with internal processes, but we do not expect this
70 An example of a <code>fillHepRup()</code> method to set up
71 <ei>Z^0</ei> production at LHC, with masses above 50 GeV, would be
73 bool LHAinitFortran::fillHepRup() {
74 Pythia6Interface::pygive("msel = 11");
75 Pythia6Interface::pygive("ckin(1) = 50.");
76 Pythia6Interface::pyinit("cms","p","p",14000.);
77 Pythia6Interface::pyupin();
81 and the process-generic <code>fillHepEup()</code> method would be
83 bool LHAupFortran::fillHepEup() {
84 Pythia6Interface::pyupev();
88 Note that, of all parameters that could be set with the
89 <code>PYGIVE</code>, only those that influence the generation of the
90 hard processes have any impact, since this is the only part of the
91 Fortran code that is used. Also, if you want to set many parameters,
92 you can easily collect them in one file (separate from PYTHIA 8
93 input) and parse that file.
96 All hard PYTHIA 6 processes should be available for full generation
97 in PYTHIA 8, at least to the extent that they are defined for beams
98 of protons and antiprotons, which is the key application for PYTHIA 8
99 so far. Soft processes, i.e. elastic and diffractive scattering, as well
100 as minimum-bias events, require a different kinematics machinery, and
101 can only be generated with the internal PYTHIA 8 processes.
104 A simple example is found in <code>main51.cc</code>, another with parsed
105 input in <code>main52.cc</code> and a third with HepMC output in
106 <code>main54.cc</code>.
110 <!-- Copyright (C) 2008 Torbjorn Sjostrand -->