3 <title>Program Flow</title>
4 <link rel="stylesheet" type="text/css" href="pythia.css"/>
5 <link rel="shortcut icon" href="pythia32.gif"/>
11 Recall that, to first order, the event generation process can be
12 subdivided into three stages:
14 <li>Initializaion.</li>
15 <li>The event loop.</li>
18 This is reflected in how the top-level <code>Pythia</code> class should
19 be used in the user-supplied main program, further outlined in the
20 following. Since the nature of the run is defined at the initialization
21 stage, this is where most of the PYTHIA user code has to be written.
22 So as not to confuse the reader unduly, the description of initialization
23 options has been subdivided into what would normally be used and what is
24 intended for more special applications.
27 At the bottom of this webpage is a complete survey of all public
28 <code>Pythia</code> methods and data members, in a more formal style
29 than the task-oriented descriptions found in the preceding sections.
30 This offers complementary information.
32 <h3>Initialization - normal usage</h3>
37 Already at the top of the main program file, you need to include the proper
42 To simplify typing, it also makes sense to declare
44 using namespace Pythia8;
50 The first step is to create a generator object,
55 It is this object that we will use from now on. Normally a run
56 will only contain one <code>Pythia</code> object. (But you can
57 use several <code>Pythia</code> objects, which then will be
58 independent of each other.)<br/>
59 By default all output from <code>Pythia</code> will be on the
60 <code>cout</code> stream, but the <code>list</code> methods below do
61 allow output to alternative streams or files.
66 You next want to set up the character of the run.
67 The pages under the "Setup Run Tasks" heading in the index
68 describe all the options available (with some very few exceptions,
69 found on the other pages).
70 The default values and your modifications are stored in two databases,
71 one for <a href="SettingsScheme.html" target="page">generic settings</a>
72 and one for <a href="ParticleDataScheme.html" target="page">particle data</a>.
73 Both of these are initialized with their default values by the
74 <code>Pythia</code> constructor. The default values can then be
75 changed, primarily by one of the two ways below, or by a combination
81 pythia.readString(string);
83 method repeatedly to do a change of a property at a time.
84 The information in the string is case-insensitive, but upper- and
85 lowercase can be combined for clarity. The rules are that<br/>
86 (i) if the first nonblank character of the string is a letter
87 it is assumed to contain a setting, and is sent on to
88 <code>pythia.settings.readString(string)</code>;<br/>
89 (ii) if instead the string begins with a digit it is assumed to
90 contain particle data updates, and so sent on to
91 <code>pythia.particleData.readString(string)</code>;<br/>
92 (iii) if none of the above, the string is assumed to be a comment,
93 i.e. nothing will be done.<br/>
94 In the former two cases, a warning is issued whenever a string
95 cannot be recognized (maybe because of a spelling mistake).<br/>
96 Some examples would be
98 pythia.readString("TimeShower:pTmin = 1.0");
99 pythia.readString("111:mayDecay = false");
101 The <code>readString(string)</code> method is intended primarily for
102 a few changes. It can also be useful if you want to construct a
103 parser for input files that contain commands both to PYTHIA and to
104 other libraries.<br/>
107 b) You can read in a file containing a list of those variables
108 you want to see changed, with a
110 pythia.readFile(fileName);
112 Each line in this file with be processes by the
113 <code>readString(string)</code> method introduced above. You can thus
114 freely mix comment lines and lines handed on to <code>Settings</code>
115 or to <code>ParticleData</code>.<br/>
116 This approach is better suited for more extensive changes than a direct
117 usage of <code>readString(string)</code>, and can also avoid having to
118 recompile and relink your main program between runs.<br/>
119 It is also possible to read input from an <code>istream</code>, by
120 default <code>cin</code>, rather than from a file. This may be convenient
121 if information is generated on-the-fly, within the same run.
124 Changes are made sequentially in the order the commands are encountered
125 during execution, meaning that if a parameter is changed several times
126 it is the last one that counts. The two special
127 <code><a href="Tunes.html" target="page">Tune:ee</a></code> and
128 <code><a href="Tunes.html" target="page">Tune:pp</a></code>
129 modes are expanded to change several settings in one go, but these obey
130 the same ordering rules.
136 Next comes the initialization stage, where all
137 remaining details of the generation are to be specified.
138 There is one standard method to use for this
141 <code>pythia.init();</code><br/>
142 with no arguments will read all relevant information from the
143 <code><a href="SettingsScheme.html" target="page">Settings</a></code>
144 and <code><a href="ParticleDataScheme.html" target="page">ParticleData</a></code>
145 databases. Specifically the setup of incoming beams and energies
146 is governed by the the beam parameters from the
147 <code><a href="BeamParameters.html" target="page">Beams</a></code>
148 group of variables. If you don't change any of those you will
149 default to proton-proton collisions at 14 TeV, i.e. the nominal LHC
153 A few alternative forms are available, where the arguments of the
154 <code>init(...)</code> call can be used to set the beam parameters.
155 These alternatives are now deprecated, and will be removed for
159 a) <code>pythia.init( idA, idB, eCM);</code><br/>
160 lets you specify the identities and the CM energy of the two incoming
161 beam particles, with A (B) assumed moving in the <i>+z (-z)</i>
165 b) <code>pythia.init( idA, idB, eA, eB);</code><br/>
166 is similar, but the two beam energies can be different, so the
167 collisions do not occur in the CM frame. If one of the beam energies
168 is below the particle mass you obtain a fixed-target topology.
171 c) <code>pythia.init( idA, idB, pxA, pyA, pzA, pxB, pyB, pzB);</code><br/>
172 is similar, but here you provide the three-momenta
173 <i>(p_x, p_y, p_z)</i> of the two incoming particles,
174 to allow for arbitrary beam directions.
177 d) <code>pythia.init(fileName);</code> <br/>
178 assumes a file in the <a href="LesHouchesAccord.html" target="page">Les Houches
179 Event File</a> format is provided.
182 e) <code>pythia.init( LHAup*);</code> <br/>
183 assumes <a href="LesHouchesAccord.html" target="page">Les Houches Accord</a>
184 initialization and event information is available in an <code>LHAup</code>
185 class object, and that a pointer to this object is handed in.
189 If you want to have a list of the generator and particle data used,
190 either only what has been changed or everything, you can use
192 pythia.settings.listChanged();
193 pythia.settings.listAll();
194 pythia.particleData.listChanged();
195 pythia.particleData.listAll();
201 <h3>The event loop</h3>
206 Inside the event generation loop you generate the
207 next event using the <code>next()</code> method,
211 This method takes no arguments; everything has already been specified.
212 It does return a bool value, however, <code>false</code> when the
213 generation failed. This can be a "programmed death" when the
214 supply of input parton-level configurations on file is exhausted.
215 It can alternatively signal a failure of <code>Pythia</code> to
216 generate an event, or unphysical features in the event record at the
217 end of the generation step. It makes sense to allow a few <code>false</code>
218 values before a run is aborted, so long as the related faulty
224 The generated event is now stored in the <code>event</code>
225 object, of type <code><a href="EventRecord.html" target="page">Event</a></code>,
226 which is a public member of <code>pythia</code>. You therefore have
227 access to all the tools described on the pages under the "Study Output"
228 header in the index. For instance, an event can be listed with
229 <code>pythia.event.list()</code>, the identity of the <i>i</i>'th
230 <a href="ParticleProperties.html" target="page">particle</a> is given by
231 <code>pythia.event[i].id()</code>, and so on.<br/>
232 The hard process - roughly the information normally stored in the
233 Les Houches Accord event record - is available as a second object,
234 <code>process</code>, also of type <code>Event</code>.<br/>
235 A third useful public object is
236 <code><a href="EventInformation.html" target="page">info</a></code>, which offers
237 a set of one-of-a kind pieces of information about the most recent
247 <li>At the end of the generation process, you can call
251 to get some run statistics, on cross sections and the number of errors
252 and warnings encountered. The alternative
253 <code>pythia.statistics(...);</code> is equivalent but deprecated.
258 <h3>Advanced usage, mainly for initialization</h3>
260 A) Necessary data are automatically loaded when you use the
261 default PYTHIA installation directory structure and run the main
262 programs in the <code>examples</code> subdirectory. However, in the
263 general case, you must provide the path of the <code>xmldoc</code>
264 directory, where default settings and particle data are found.
265 This can be done in two ways.
270 You can set the environment variable <code>PYTHIA8DATA</code> to
271 contain the location of the <code>xmldoc</code> directory. In the
272 <code>csh</code> and <code>tcsh</code> shells this could e.g. be
274 setenv PYTHIA8DATA /home/myname/pythia81xx/xmldoc
276 while in other shells it could be
278 export PYTHIA8DATA=/home/myname/pythia81xx/xmldoc
280 where xx is the subversion number.<br/>
281 Recall that environment variables set locally are only defined in the
282 current instance of the shell. The above lines should go into your
283 <code>.cshrc</code> and <code>.bashrc</code> files, respectively,
284 if you want a more permanent assignment.
289 You can provide the path as argument to the <code>Pythia</code>
292 Pythia pythia("/home/myname/pythia81xx/xmldoc");
296 where again xx is the subversion number.<br/>
297 When <code>PYTHIA8DATA</code> is set it takes precedence, else
298 the path in the constructor is used, else one defaults to the
299 <code>../xmldoc</code> directory.
302 B) You can override the default behaviour of PYTHIA not only by the
303 settings and particle data, but also by replacing some of the
304 PYTHIA standard routines by ones of your own. Of course, this is only
305 possible if your routines fit into the general PYTHIA framework.
306 Therefore they must be coded according to the the rules relevant
307 in each case, as a derived class of a PYTHIA base class, and a pointer
308 to such an object must be handed in by one of the methods below.
309 These calls must be made before the <code>pythia.init(...)</code> call.
314 If you are not satisfied with the list of parton density functions that
315 are implemented internally or available via the LHAPDF interface
316 (see the <a href="PDFSelection.html" target="page">PDF Selection</a> page), you
317 can supply your own by a call to the <code>setPDFPtr(...)</code> method
319 pythia.setPDFptr( pdfAPtr, pdfBPtr);
321 where <code>pdfAPtr</code> and <code>pdfBPtr</code> are pointers to
322 two <code>Pythia</code> <a href="PartonDistributions.html" target="page">PDF
323 objects</a>. Note that <code>pdfAPtr</code> and <code>pdfBPtr</code>
324 cannot point to the same object; even if the PDF set is the same,
325 two copies are needed to keep track of two separate sets of <i>x</i>
326 and density values.<br/>
327 If you further wish to use separate PDF's for the hard process of an
328 event than the ones being used for everything else, the extended form
330 pythia.setPDFptr( pdfAPtr, pdfBPtr, pdfHardAPtr, pdfHardBPtr);
332 allows you to specify those separately, and then the first two sets
333 would only be used for the showers and for multiparton interactions.
338 If you want to link to an external generator that feeds in events
339 in the LHA format, you can call the <code>setLHAupPtr(...)</code>
342 pythia.setLHAupPtr( lhaUpPtr);
344 where the <code>lhaUpPtr</code> derives from the
345 <a href="LesHouchesAccord.html" target="page">LHAup</a> base class.
350 If you want to perform some particle decays with an
351 external generator, you can call the <code>setDecayPtr(...)</code>
354 pythia.setDecayPtr( decayHandlePtr, particles);
356 where the <code>decayHandlePtr</code> derives from the
357 <code><a href="ExternalDecays.html" target="page">DecayHandler</a></code> base
358 class and <code>particles</code> is a vector of particle codes to be
364 If you want to use an external random number generator,
365 you can call the <code>setRndmEnginePtr(...)</code> method
367 pythia.setRndmEnginePtr( rndmEnginePtr);
369 where <code>rndmEnginePtr</code> derives from the
370 <code><a href="RandomNumbers.html" target="page">RndmEngine</a></code> base class.
371 The <code>Pythia</code> default random number generator is perfectly
372 good, so this is only intended for consistency in bigger frameworks.
377 If you want to interrupt the evolution at various stages,
378 to interrogate the event and possibly veto it, or you want to
379 reweight the cross section, you can use
381 pythia.setUserHooksPtr( userHooksPtr);
383 where <code>userHooksPtr</code> derives from the
384 <code><a href="UserHooks.html" target="page">UserHooks</a></code> base class.
389 If you want to use your own merging scale definition for
390 matrix element + parton shower merging, you can call
392 pythia.setMergingHooksPtr( mergingHooksPtr);
394 where <code>mergingHooksPtr</code> derives from the
395 <code><a href="MatrixElementMerging.html" target="page">MergingHooks</a></code> base class.
400 If you want to use your own parametrization of beam momentum spread and
401 interaction vertex, rather than the provided simple Gaussian
402 parametrization (off by default), you can call
404 pythia.setBeamShapePtr( beamShapePtr);
406 where <code>beamShapePtr</code> derives from the
407 <code><a href="BeamShape.html" target="page">BeamShape</a></code> base class.
412 If you want to implement a cross section of your own, but still make use
413 of the built-in phase space selection machinery, you can use
415 pythia.setSigmaPtr( sigmaPtr);
417 where <code>sigmaPtr</code> of type <code>SigmaProcess*</code> is an
418 instance of a class derived from one of the <code>Sigma1Process</code>,
419 <code>Sigma2Process</code> and <code>Sigma3Process</code> base classes
420 in their turn derived from
421 <code><a href="SemiInternalProcesses.html" target="page">SigmaProcess</a></code>.
422 This call can be used repeatedly to hand in several different processes.
427 If your cross section contains the production of a new resonance
428 with known analytical expression for all the relevant partial widths,
429 you can make this resonance available to the program with
431 pythia.setResonancePtr( resonancePtr);
433 where <code>resonancePtr</code> of type <code>ResonanceWidths*</code>
434 is an instance of a class derived from the
435 <code><a href="SemiInternalResonances.html" target="page">ResonanceWidths</a></code>
436 base class. In addition you need to add the particle to the normal
437 <a href="ParticleDataScheme.html" target="page">particle and decay database</a>.
438 This procedure can be used repeatedly to hand in several different
444 If you are a real expert and want to <a href="ImplementNewShowers.html" target="page">replace
445 the PYTHIA initial- and final-state showers</a>, you can use
447 pythia.setShowerPtr( timesDecPtr, timesPtr, spacePtr);
449 where <code>timesDecPtr</code> and <code>timesPtr</code>
450 derive from the <code>TimeShower</code> base class, and
451 <code>spacePtr</code> from <code>SpaceShower</code>.
457 C) Some comments on collecting several tasks in the same run.
461 PYTHIA has not been written for threadsafe execution on multicore
462 processors. If you want to use all cores,
463 the most efficient way presumably is to start correspondingly many jobs,
464 with different random number seeds, and add the statistics at the end.
465 However, note that several instances can be set up in the same main
466 program, since instances are completely independent of each other,
467 so each instance could be run inside a separate thread.
472 In some cases it is convenient to use more than one <code>Pythia</code>
473 object. The key example would be the simultaneous generation of signal
474 and pileup events, see <code>main19.cc</code>. The two objects are then
475 set up and initialized separately, and generate events completely
476 independently of each other. It is only afterwards that the event records
477 are combined into one single super-event per beam crossing.
482 When time is not an issue, it may be that you want to perform several
483 separate subruns sequentially inside a run, e.g. to combine results for
484 several kinematical regions or to compare results for some different
485 tunes of the underlying event. One way to go is to create (and destroy)
486 one <code>pythia</code> object for each subrun, in which case they are
487 completely separate. You can also use the same <code>pythia</code> object,
488 only doing a new <code>init(...)</code> call for each subrun. In that
489 case, the settings and particle databases remain as they were in the
490 previous subrun, only affected by the specific changes you introduced in
491 the meantime. You can put those changes in the main program, with
492 <code>pythia.readString(string)</code>, using your own logic to decide
493 which ones to execute in which subrun. A corresponding possibility
494 exists with <code>pythia.readFile(fileName, subrun)</code> (or an
495 <code>istream</code> instead of a <code>fileName</code>), which as second
496 argument can take a non-negative subrun number. Then only those
497 sections of the file before any <code>Main:subrun = ...</code> line
498 or with matching <code>subrun</code> number will be read. That is, the
499 file could have a structure like
501 ( lines always read, i.e. "default values" always (re)set )
503 ( lines only read with readFile(fileName, 1) )
505 ( lines only read with readFile(fileName, 2) )
507 Both of these possibilities are illustrated in <code>main08.cc</code>.
512 When working with Les Houches Event Files, it may well be that your
513 intended input event sample is spread over several files, that you all
514 want to turn into complete events in one and the same run. There is no
515 problem with looping over several subruns, where each new subrun
516 is initialized with a new file, with name set in <code>Beams:LHEF</code>.
517 However, in that case you will do a complete re-initialization each time
518 around. If you want to avoid this, note that the flag
519 <code>Beams:newLHEFsameInit = true</code> can be set for the second and
520 subsequent subruns. Then the new file will be simulated with the same
521 initialization data as already set in a previous
522 <code>pythia.init()</code> call. The burden rests on you to ensure
523 that this is indeed correct, e.g. that the two event samples have not
524 been generated for different beam energies. Also note that cross
525 sections for processes will be based on the information in the
526 first-read file, when the full initialization is performed.
531 <h2>The Pythia Class</h2>
533 Here follows the complete survey of all public <code>Pythia</code>
534 methods and data members.
536 <h3>Constructor and destructor</h3>
538 <a name="method1"></a>
539 <p/><strong>Pythia::Pythia(string xmlDir = "../xmldoc",bool printBanner = true) </strong> <br/>
540 creates an instance of the <code>Pythia</code> event generators,
541 and sets initial default values, notably for all settings and
542 particle data. You may use several <code>Pythia</code> instances
543 in the same run; only when you want to access external static
544 libraries could this cause problems. (This includes in particular
545 Fortran libraries such as <a href="PDFSelection.html" target="page">LHAPDF</a>.)
546 <br/><code>argument</code><strong> xmlDir </strong> (<code>default = <strong>../xmldoc</strong></code>) : allows you to choose
547 from which directory the default settings and particle data values
548 are read in. If the <code>PYTHIA8DATA</code> environment variable
549 has been set it takes precedence. Else this optional argument allows
550 you to choose another directory location than the default one. Note
551 that it is only the directory location you can change, its contents
552 must be the ones of the <code>xmldoc</code> directory in the
553 standard distribution.
555 <br/><code>argument</code><strong> printBanner </strong> (<code>default = <strong>true</strong></code>) : can be set
556 <code>false</code> to stop the program from printing a banner.
557 The banner contains useful information, so this option is only
558 intended for runs with multiple <code>Pythia</code> instances,
559 where output needs to be restricted.
563 <a name="method2"></a>
564 <p/><strong>Pythia::~Pythia </strong> <br/>
565 the destructor deletes the objects created by the constructor.
569 <a name="method3"></a>
570 <p/><strong>bool Pythia::readString(string line, bool warn = true) </strong> <br/>
571 reads in a single string, that is interpreted as an instruction to
572 modify the value of a <a href="SettingsScheme.html" target="page">setting</a> or
573 <a href="ParticleDataScheme.html" target="page">particle data</a>, as already described
575 <br/><code>argument</code><strong> line </strong> :
576 the string to be interpreted as an instruction.
578 <br/><code>argument</code><strong> warn </strong> (<code>default = <strong>true</strong></code>) :
579 write a warning message or not whenever the instruction does not make
580 sense, e.g. if the variable does not exist in the databases.
582 <br/><b>Note:</b> the method returns false if it fails to
583 make sense out of the string.
586 <a name="method4"></a>
587 <p/><strong>bool Pythia::readFile(string fileName, bool warn = true, int subrun = SUBRUNDEFAULT) </strong> <br/>
589 <strong>bool Pythia::readFile(string fileName, int subrun = SUBRUNDEFAULT) </strong> <br/>
591 <strong>bool Pythia::readFile(istream& inStream = cin, bool warn = true, int subrun = SUBRUNDEFAULT) </strong> <br/>
593 <strong>bool Pythia::readFile(istream& inStream = cin, int subrun = SUBRUNDEFAULT) </strong> <br/>
594 reads in a whole file, where each line is interpreted as an instruction
595 to modify the value of a <a href="SettingsScheme.html" target="page">setting</a> or
596 <a href="ParticleDataScheme.html" target="page">particle data</a>, cf. the above
597 <code>readString</code> method. All four forms of the
598 <code>readFile</code> command share code for actually reading a file.
599 <br/><code>argument</code><strong> fileName </strong> :
600 the file from which instructions are read.
602 <br/><code>argument</code><strong> inStream </strong> :
603 an istream from which instructions are read.
605 <br/><code>argument</code><strong> warn </strong> (<code>default = <strong>true</strong></code>) :
606 write a warning message or not whenever the instruction does not make
607 sense, e.g. if the variable does not exist in the databases. In the
608 command forms where <code>warn</code> is omitted it is true.
610 <br/><code>argument</code><strong> subrun </strong> :
611 allows you have several optional sets of commands within the same file.
612 Only those sections of the file before any <code>Main:subrun = ...</code>
613 line or following such a line with matching subrun number will be read.
614 The subrun number should not be negative; negative codes like
615 <code>SUBRUNDEFAULT</code> corresponds to no specific subrun.
617 <br/><b>Note:</b> the method returns false if it fails to
618 make sense out of any one line.
621 <a name="method5"></a>
622 <p/><strong>bool Pythia::setPDFPtr( PDF* pdfAPtr, PDF* pdfBPtr, PDF* pdfHardAPtr = 0, PDF* pdfHardBPtr = 0) </strong> <br/>
623 offers the possibility to link in external PDF sets for usage inside
624 the program. The rules for constructing your own class from
625 the <code>PDF</code> base class are described
626 <a href="PartonDistributions.html" target="page">here</a>.
627 <br/><code>argument</code><strong> pdfAPtr, pdfBPtr </strong> :
628 pointers to two <code>PDF</code>-derived objects, one for each of
629 the incoming beams. The two objects have to be instantiated by you
630 in your program. Even if the two beam particles are the same
631 (protons, say) two separate instances are required, since current
632 information is cached in the objects. If both arguments are zero
633 then any previous linkage to external PDF's is disconnected,
634 see further Note 2 below.
636 <br/><code>argument</code><strong> pdfHardAPtr, pdfHardBPtr </strong> (<code>default = <strong>0</strong></code>) :
637 pointers to two further <code>PDF</code>-derived objects, one for each
638 of the incoming beams. Normally only the first two arguments above would
639 be used, and then the same PDF sets would be invoked everywhere. If you
640 provide these two further pointers then two different sets of PDF's are
641 used. This second set is then exclusively for the generation of the hard
642 process from the process matrix elements library. The first set above
643 is for everything else, notably parton showers and multiparton interactions.
645 <br/><b>Note 1:</b> The method returns false if the input is obviously
646 incorrect, e.g. if two (nonzero) pointers agree.
647 <br/><b>Note 2:</b> If you want to combine several subruns you can
648 call <code>setPDFPtr</code> with new arguments before each
649 <code>Pythia::init(...)</code> call. To revert from external PDF's
650 to the normal internal PDF selection you must call
651 <code>setPDFPtr(0, 0)</code> before <code>Pythia::init(...)</code>.
654 <a name="method6"></a>
655 <p/><strong>bool Pythia::setLHAupPtr( LHAup* lhaUpPtrIn) </strong> <br/>
656 offers linkage to an external generator that feeds in events
657 in the LHA format, see
658 <a href="LesHouchesAccord.html" target="page">Les Houches Accord</a>,
660 <code><a href="BeamParameters.html" target="page">Beams:frameType = 5</a></code>
662 <br/><code>argument</code><strong> lhaUpPtrIn </strong> :
663 pointer to a <code>LHAup</code>-derived object.
665 <br/><b>Note:</b> The method currently always returns true.
668 <a name="method7"></a>
669 <p/><strong>bool Pythia::setDecayPtr( DecayHandler* decayHandlePtr, vector<int> handledParticles) </strong> <br/>
670 offers the possibility to link to an external program that can do some
671 of the particle decays, instead of using the internal decay machinery.
672 With particles we here mean the normal hadrons and leptons, not
673 top quarks, electroweak bosons or new particles in BSM scenarios.
674 The rules for constructing your own class from the
675 <code>DecayHandler</code> base class are described
676 <a href="ExternalDecays.html" target="page">here</a>. Note that you can only
677 provide one external object, but this object in its turn could
678 very well hand on different particles to separate decay libraries.
679 <br/><code>argument</code><strong> decayHandlePtr </strong> :
680 pointer to a <code>DecayHandler</code>-derived object. This object
681 must be instantiated by you in your program.
683 <br/><code>argument</code><strong> handledParticles </strong> : vector with the PDG identity codes
684 of the particles that should be handled by the external decay package.
685 You should only give the particle (positive) codes; the respective
686 antiparticle is always included as well.
688 <br/><b>Note:</b> The method currently always returns true.
691 <a name="method8"></a>
692 <p/><strong>bool Pythia::setRndmEnginePtr( RndmEngine* rndmEnginePtr) </strong> <br/>
693 offers the possibility to link to an external random number generator.
694 The rules for constructing your own class from the
695 <code>RndmEngine</code> base class are described
696 <a href="RandomNumbers.html" target="page">here</a>.
697 <br/><code>argument</code><strong> rndmEnginePtr </strong> :
698 pointer to a <code>RndmEngine</code>-derived object. This object
699 must be instantiated by you in your program.
701 <br/><b>Note:</b> The method returns true if the pointer is different
705 <a name="method9"></a>
706 <p/><strong>bool Pythia::setUserHooksPtr( UserHooks* userHooksPtr) </strong> <br/>
707 offers the possibility to interact with the generation process at
708 a few different specified points, e.g. to reject undesirable events
709 at an early stage to save computer time. The rules for constructing
710 your own class from the <code>UserHooks</code> base class are described
711 <a href="UserHooks.html" target="page">here</a>. You can only hand in one such
712 pointer, but this may be to a class that implements several of the
713 different allowed possibilities.
714 <br/><code>argument</code><strong> userHooksPtr </strong> :
715 pointer to a <code>userHooks</code>-derived object. This object
716 must be instantiated by you in your program.
718 <br/><b>Note:</b> The method currently always returns true.
721 <a name="method10"></a>
722 <p/><strong>bool Pythia::setBeamShapePtr( BeamShape* beamShapePtr) </strong> <br/>
723 offers the possibility to provide your own shape of the momentum and
724 space-time spread of the incoming beams. The rules for constructing
725 your own class from the <code>BeamShape</code> base class are described
726 <a href="BeamShape.html" target="page">here</a>.
727 <br/><code>argument</code><strong> BeamShapePtr </strong> :
728 pointer to a <code>BeamShape</code>-derived object. This object
729 must be instantiated by you in your program.
731 <br/><b>Note:</b> The method currently always returns true.
734 <a name="method11"></a>
735 <p/><strong>bool Pythia::setSigmaPtr( SigmaProcess* sigmaPtr) </strong> <br/>
736 offers the possibility to link your own implementation of a process
737 and its cross section, to make it a part of the normal process
738 generation machinery, without having to recompile the
739 <code>Pythia</code> library itself. The rules for constructing your
740 own class from the <code>SigmaProcess</code> base class are described
741 <a href="SemiInternalProcesses.html" target="page">here</a>. You may call this
742 routine repeatedly, to add as many new processes as you wish.
743 <br/><code>argument</code><strong> sigmaPtr </strong> :
744 pointer to a <code>SigmaProcess</code>-derived object. This object
745 must be instantiated by you in your program.
747 <br/><b>Note:</b> The method currently always returns true.
750 <a name="method12"></a>
751 <p/><strong>bool Pythia::setResonancePtr( ResonanceWidths* resonancePtr) </strong> <br/>
752 offers the possibility to link your own implementation of the
753 calculation of partial resonance widths, to make it a part of the
754 normal process generation machinery, without having to recompile the
755 <code>Pythia</code> library itself. This allows the decay of new
756 resonances to be handled internally, when combined with new particle
757 data. Note that the decay of normal hadrons cannot be modeled here;
758 this is for New Physics resonances. The rules for constructing your
759 own class from the <code>ResonanceWidths</code> base class are described
760 <a href="SemiInternalResonances.html" target="page">here</a>. You may call this
761 routine repeatedly, to add as many new resonances as you wish.
762 <br/><code>argument</code><strong> resonancePtr </strong> :
763 pointer to a <code>ResonanceWidths</code>-derived object. This object
764 must be instantiated by you in your program.
766 <br/><b>Note:</b> The method currently always returns true.
769 <a name="method13"></a>
770 <p/><strong>bool Pythia::setShowerPtr( TimeShower* timesDecPtr, TimeShower* timesPtr = 0, SpaceShower* spacePtr = 0) </strong> <br/>
771 offers the possibility to link your own parton shower routines as
772 replacements for the default ones. This is much more complicated
773 since the showers are so central and are so interlinked with other
774 parts of the program. Therefore it is also possible to do the
775 replacement in stages, from the more independent to the more
776 intertwined. The rules for constructing your own classes from the
777 <code>TimeShower</code> and <code>SpaceShower</code>base classes
778 are described <a href="ImplementNewShowers.html" target="page">here</a>. These
779 objects must be instantiated by you in your program.
780 <br/><code>argument</code><strong> timesDecPtr </strong> :
781 pointer to a <code>TimeShower</code>-derived object for doing
782 timelike shower evolution in resonance decays, e.g. of a
783 <i>Z^0</i>. This is decoupled from beam remnants and parton
784 distributions, and is therefore the simplest kind of shower
785 to write. If you provide a value 0 then the internal shower
786 routine will be used.
788 <br/><code>argument</code><strong> timesPtr </strong> (<code>default = <strong>0</strong></code>) :
789 pointer to a <code>TimeShower</code>-derived object for doing
790 all other timelike shower evolution, which is normally interleaved
791 with multiparton interactions and spacelike showers, introducing
792 both further physics and further technical issues. If you retain
793 the default value 0 then the internal shower routine will be used.
794 You are allowed to use the same pointer as above for the
795 <code>timesDecPtr</code> if the same shower can fulfill both tasks.
797 <br/><code>argument</code><strong> spacePtr </strong> (<code>default = <strong>0</strong></code>) :
798 pointer to a <code>SpaceShower</code>-derived object for doing
799 all spacelike shower evolution, which is normally interleaved
800 with multiparton interactions and timelike showers. If you retain
801 the default value 0 then the internal shower routine will be used.
803 <br/><b>Note:</b> The method currently always returns true.
808 At the initialization stage all the information provided above is
809 processed, and the stage is set up for the subsequent generation
810 of events. Currently several alternative forms of the <code>init</code>
811 method are available for this stage, but only the first one is
814 <a name="method14"></a>
815 <p/><strong>bool Pythia::init() </strong> <br/>
816 initialize for collisions, in any of the five separate possibilities
817 below. In this option the beams are not specified by input arguments,
818 but instead by the settings in the
819 <a href="BeamParameters.html" target="page">Beam Parameters</a> section.
820 This allows the beams to be specified in the same file as other
821 run instructions. The default settings give pp collisions at 14 TeV.
822 <br/><b>Note:</b> The method returns false if the
823 initialization fails. It is then not possible to generate any
827 <a name="method15"></a>
828 <p/><strong>bool Pythia::init( int idA, int idB, double eCM) </strong> <br/>
829 initialize for collisions in the center-of-mass frame, with the
830 beams moving in the <i>+-z</i> directions.
831 <br/><code>argument</code><strong> idA, idB </strong> :
832 particle identity code for the two incoming beams.
834 <br/><code>argument</code><strong> eCM </strong> :
835 the CM energy of the collisions.
837 <br/><b>Notes:</b> Deprecated. The method returns false if the
838 initialization fails. It is then not possible to generate any
842 <a name="method16"></a>
843 <p/><strong>bool Pythia::init( int idA, int idB, double eA, double eB) </strong> <br/>
844 initialize for collisions with back-to-back beams,
845 moving in the <i>+-z</i> directions, but with different energies.
846 <br/><code>argument</code><strong> idA, idB </strong> :
847 particle identity code for the two incoming beams.
849 <br/><code>argument</code><strong> eA, eB </strong> :
850 the energies of the two beams. If an energy is set to be below
851 the mass of the respective beam particle that particle is taken to
852 be at rest. This offers a simple possibility to simulate
853 fixed-target collisions.
855 <br/><b>Notes:</b> Deprecated. The method returns false if the
856 initialization fails. It is then not possible to generate any
860 <a name="method17"></a>
861 <p/><strong>bool Pythia::init( int idA, int idB, double pxA, double pyA, double pzA, double pxB, double pyB, double pzB) </strong> <br/>
862 initialize for collisions with arbitrary beam directions.
863 <br/><code>argument</code><strong> idA, idB </strong> :
864 particle identity code for the two incoming beams.
866 <br/><code>argument</code><strong> pxA, pyA, pzA </strong> :
867 the three-momentum vector <i>(p_x, p_y, p_z)</i> of the first
870 <br/><code>argument</code><strong> pxB, pyB, pzB </strong> :
871 the three-momentum vector <i>(p_x, p_y, p_z)</i> of the second
874 <br/><b>Notes:</b> Deprecated. The method returns false if the
875 initialization fails. It is then not possible to generate any
879 <a name="method18"></a>
880 <p/><strong>bool Pythia::init( string LesHouchesEventFile, bool skipInit = false) </strong> <br/>
881 initialize for hard-process collisions fed in from an external file
882 with events, written according to the
883 <a href="LesHouchesAccord.html" target="page">Les Houches Event File</a>
885 <br/><code>argument</code><strong> LesHouchesEventFile </strong> :
886 the file name (including path, where required) where the
887 events are stored, including relevant information on beam
888 identities and energies.
890 <br/><code>argument</code><strong> skipInit </strong> (<code>default = <strong>false</strong></code>) :
891 By default this method does a complete reinitialization of the
892 generation process. If you set this argument to true then
893 no reinitialization will occur, only the pointer to the event
894 file is updated. This may come in handy if the full event sample
895 is split across several files generated under the same conditions
896 (except random numbers, of course). You then do the first
897 initialization with the default, and all subsequent ones with
898 true. Note that things may go wrong if the files are not created
899 under the same conditions.
901 <br/><b>Notes:</b> Deprecated. The method returns false if the
902 initialization fails. It is then not possible to generate any
906 <a name="method19"></a>
907 <p/><strong>bool Pythia::init( LHAup* lhaUpPtr) </strong> <br/>
908 initialize for hard-process collisions fed in from an external
909 source of events, consistent with the Les Houches Accord standard.
910 The rules for constructing your own class from the <code>LHAup</code>
911 base class are described <a href="LesHouchesAccord.html" target="page">here</a>.
912 This class is also required to provide the beam parameters.
913 <br/><code>argument</code><strong> lhaUpPtr </strong> :
914 pointer to a <code>LHAup</code>-derived object. This object
915 must be instantiated by you in your program.
917 <br/><b>Notes:</b> Deprecated. The method returns false if the
918 initialization fails. It is then not possible to generate any
922 <h3>Generate events</h3>
924 The <code>next()</code> method is the main one to generate events.
925 In this section we also put a few other specialized methods that
926 may be useful in some circumstances.
928 <a name="method20"></a>
929 <p/><strong>bool Pythia::next() </strong> <br/>
930 generate the next event. No input parameters are required; all
931 instructions have already been set up in the initialization stage.
932 <br/><b>Note:</b> The method returns false if the event generation
933 fails. The event record is then not consistent and should not be
934 studied. When reading in hard collisions from a Les Houches Event File
935 the problem may be that the end of the file has been reached. This
936 can be checked with the
937 <code><a href="EventInformation.html" target="page">Info::atEndOfFile()</a></code>
941 <a name="method21"></a>
942 <p/><strong>int Pythia::forceTimeShower( int iBeg, int iEnd, double pTmax, int nBranchMax = 0) </strong> <br/>
943 perform a final-state shower evolution on partons in the
944 <code>event</code> event record. This could be used for externally
945 provided simple events, or even parts of events, for which
946 a complete generation is not foreseen. Since the mother source of
947 the parton system is not known, one cannot expect as good accuracy
948 as in a normal generation. When two different timelike shower
949 instances are set up, it is the one used for showering in resonance
950 decays that is used here. The <code>forceTimeShower</code> method
951 can be used in conjunction with the <code>forceHadronLevel</code>
952 one below. Further comments are found
953 <a href="HadronLevelStandalone.html" target="page">here</a>.
954 <br/><code>argument</code><strong> iBeg, iEnd </strong> : the first and last entry of the event
955 record to be affected by the call.
957 <br/><code>argument</code><strong> pTmax </strong> : the maximum <i>pT</i> scale of emissions.
958 Additionally, as always, the <code>scale</code> variable of each parton
959 sets the maximum <i>pT</i> scale of branchings of this parton.
960 Recall that this scale defaults to 0 if not set, so that no radiation
963 <br/><code>argument</code><strong> nBranchMax </strong> (<code>default = <strong>0</strong></code>) : when positive, it sets the
964 maximum number of branchings that are allowed to occur in the shower,
965 i.e. the shower may stop evolving before reaching the lower cutoff.
966 The argument has no effect when zero or negative, i.e. then the shower
967 will continue to the lower cutoff.
969 <br/><b>Note:</b> The method returns the number of branchings that
973 <a name="method22"></a>
974 <p/><strong>bool Pythia::forceHadronLevel(bool findJunctions = true) </strong> <br/>
975 hadronize the existing event record, i.e. perform string fragmentation
976 and particle decays. There are two main applications. Firstly,
977 you can use the same parton-level content as a basis for repeated
978 hadronization attempts, in schemes intended to save computer time.
979 Secondly, you may have an external program that can simulate the full
980 partonic level of the event - hard process, parton showers, multiparton
981 interactions, beam remnants, colour flow, and so on - but not
982 hadronization. Further details are found
983 <a href="HadronLevelStandalone.html" target="page">here</a>.
984 <br/><code>argument</code><strong> findJunctions </strong> (<code>default = <strong>true</strong></code>) :
985 normally this routine will search through the event record and try to
986 figure out if any colour junctions are present. If so, the colour
987 topology of such junctions must be sorted out. In tricky cases this
988 might fail, and then hadronization will not work. A user who is
989 aware of this and knows the intended colour flow can set up the
990 junction information in the event record, and then call
991 <code>forceHadronLevel(false)</code> so as not to have this information
993 <br/><b>Note:</b> The method returns false if the hadronization
994 fails. The event record is then not consistent and should not be
998 <a name="method23"></a>
999 <p/><strong>bool Pythia::moreDecays() </strong> <br/>
1000 perform decays of all particles in the event record that have not been
1001 decayed but should have been done so. This can be used e.g. for
1002 repeated decay attempts, in schemes intended to save computer time.
1003 Further details are found <a href="HadronLevelStandalone.html" target="page">here</a>.
1004 <br/><b>Note:</b> The method returns false if the decays fail. The
1005 event record is then not consistent and should not be studied.
1008 <a name="method24"></a>
1009 <p/><strong>bool Pythia::forceRHadronDecays() </strong> <br/>
1010 perform decays of R-hadrons that were previously considered stable.
1011 This could be if an R-hadron is sufficiently long-lived that
1012 it may interact in the detector between production and decay, so that
1013 its four-momentum is changed. Further details are found
1014 <a href="RHadrons.html" target="page">here</a>.
1015 <br/><b>Note:</b> The method returns false if the decays fail. The
1016 event record is then not consistent and should not be studied.
1019 <a name="method25"></a>
1020 <p/><strong>void Pythia::LHAeventList(ostream& os = cout) </strong> <br/>
1021 list the Les Houches Accord information on the current event, see
1022 <code><a href="LesHouchesAccord.html" target="page">LHAup::listEvent(...)</a></code>.
1023 (Other listings are available via the class members below, so this
1024 listing is a special case that would not fit elsewhere.)
1025 <br/><code>argument</code><strong> os </strong> (<code>default = <strong>cout</strong></code>) :
1026 output stream where the listing occurs.
1030 <a name="method26"></a>
1031 <p/><strong>bool Pythia::LHAeventSkip(int nSkip) </strong> <br/>
1032 skip ahead a number of events in the Les Houches generation
1033 sequence, without doing anything further with them, see
1034 <code><a href="LesHouchesAccord.html" target="page">LHAup::skipEvent(nSkip)</a></code>.
1035 Mainly intended for debug purposes, e.g. when an event at a known
1036 location in a Les Houches Event File is causing problems.
1037 <br/><code>argument</code><strong> nSkip </strong> :
1038 number of events to skip.
1040 <br/><b>Note:</b> The method returns false if the operation fails,
1041 specifically if the end of a LHEF has been reached, cf.
1042 <code>next()</code> above.
1047 There is no required finalization step; you can stop generating events
1048 when and how you want. It is still recommended that you make it a
1049 routine to call the following method at the end. A second method provides
1050 a deprecated alternative.
1052 <a name="method27"></a>
1053 <p/><strong>void Pythia::stat() </strong> <br/>
1054 list statistics on the event generation, specifically total and partial
1055 cross sections and the number of different errors. For more details see
1056 <a href="EventStatistics.html" target="page">here</a> and for available options
1057 <a href="MainProgramSettings.html" target="page">here</a>.
1060 <a name="method28"></a>
1061 <p/><strong>void Pythia::statistics(bool all = false, bool reset = false) </strong> <br/>
1062 list statistics on the event generation, specifically total and partial
1063 cross sections and the number of different errors. For more details see
1064 <a href="EventStatistics.html" target="page">here</a>.
1065 <br/><code>argument</code><strong> all </strong> (<code>default = <strong>false</strong></code>) :
1066 if true also statistics on multiparton interactions is shown, by default not.
1068 <br/><code>argument</code><strong> reset </strong> (<code>default = <strong>false</strong></code>) : if true then all counters,
1069 e.g on events generated and errors experienced, are reset to zero
1070 whenever the routine is called. The default instead is that
1071 all stored statistics information is unaffected by the call. Counters
1072 are automatically reset in each new <code>Pythia::init(...)</code>
1073 call, however, so the only time the <code>reset</code> option makes a
1074 difference is if <code>statistics(...)</code> is called several times
1077 <br/><b>Note:</b> Deprecated.
1080 <h3>Interrogate settings</h3>
1082 Normally settings are used in the setup and initialization stages
1083 to determine the character of a run, e.g. read from a file with the
1084 above-described <code>Pythia::readFile(...)</code> method.
1085 There is no strict need for a user to interact with the
1086 <code>Settings</code> database in any other way. However, as an option,
1087 some settings variables have been left free for the user to set in
1088 such a file, and then use in the main program to directly affect the
1089 performance of that program, see
1090 <a href="MainProgramSettings.html" target="page">here</a>. A typical example would
1091 be the number of events to generate. For such applications the
1092 following shortcuts to some <code>Settings</code> methods may be
1095 <a name="method29"></a>
1096 <p/><strong>bool Pythia::flag(string key) </strong> <br/>
1097 read in a boolean variable from the <code>Settings</code> database.
1098 <br/><code>argument</code><strong> key </strong> :
1099 the name of the variable to be read.
1103 <a name="method30"></a>
1104 <p/><strong>int Pythia::mode(string key) </strong> <br/>
1105 read in an integer variable from the <code>Settings</code> database.
1106 <br/><code>argument</code><strong> key </strong> :
1107 the name of the variable to be read.
1111 <a name="method31"></a>
1112 <p/><strong>double Pythia::parm(string key) </strong> <br/>
1113 read in a double-precision variable from the <code>Settings</code>
1115 <br/><code>argument</code><strong> key </strong> :
1116 the name of the variable to be read.
1120 <a name="method32"></a>
1121 <p/><strong>string Pythia::word(string key) </strong> <br/>
1122 read in a string variable from the <code>Settings</code> database.
1123 <br/><code>argument</code><strong> key </strong> :
1124 the name of the variable to be read.
1128 <h3>Get a PDF set</h3>
1130 <code>Pythia</code> contains an number of parton density sets
1131 internally, plus an interface to LHAPDF. With the method below,
1132 this machinery is also made available for external usage.
1134 <a name="method33"></a>
1135 <p/><strong>PDF* getPDFPtr(int id, int sequence = 1) </strong> <br/>
1136 get a pointer to a PDF object. Which PDF is returned depends on the
1137 <a href="PDFSelection.html" target="page">PDF Selection</a> settings.
1138 <br/><code>argument</code><strong> id </strong> :
1139 the identity code of the incoming particle.
1141 <br/><code>argument</code><strong> sequence </strong> :
1142 should normally be 1, but 2 can be used for protons to let the PDF
1143 selection be determined by the special settings for hard processes
1144 (<code>PDF:useHard</code> etc.).
1148 <h3>Data members</h3>
1150 The <code>Pythia</code> class contains a few public data members,
1151 several of which play a central role. We list them here, with
1152 links to the places where they are further described.
1154 <a name="method34"></a>
1155 <p/><strong>Event Pythia::process </strong> <br/>
1156 the hard-process event record, see <a href="EventRecord.html" target="page">here</a>
1157 for further details.
1160 <a name="method35"></a>
1161 <p/><strong>Event Pythia::event </strong> <br/>
1162 the complete event record, see <a href="EventRecord.html" target="page">here</a>
1163 for further details.
1166 <a name="method36"></a>
1167 <p/><strong>Info Pythia::info </strong> <br/>
1168 further information on the event-generation process, see
1169 <a href="EventInformation.html" target="page">here</a> for further details.
1172 <a name="method37"></a>
1173 <p/><strong>Settings Pythia::settings </strong> <br/>
1174 the settings database, see <a href="SettingsScheme.html" target="page">here</a>
1175 for further details.
1178 <a name="method38"></a>
1179 <p/><strong>ParticleData Pythia::particleData </strong> <br/>
1180 the particle properties and decay tables database, see
1181 <a href="ParticleDataScheme.html" target="page">here</a> for further details.
1184 <a name="method39"></a>
1185 <p/><strong>Rndm Pythia::rndm </strong> <br/>
1186 the random number generator, see <a href="RandomNumberSeed.html" target="page">here</a>
1187 and <a href="RandomNumbers.html" target="page">here</a> for further details.
1190 <a name="method40"></a>
1191 <p/><strong>CoupSM Pythia::coupSM </strong> <br/>
1192 Standard Model couplings and mixing matrices, see
1193 <a href="StandardModelParameters.html" target="page">here</a> for further details.
1196 <a name="method41"></a>
1197 <p/><strong>SusyLesHouches Pythia::slha </strong> <br/>
1198 parameters and particle data in the context of supersymmetric models,
1199 see <a href="SUSYLesHouchesAccord.html" target="page">here</a> for further details.
1202 <a name="method42"></a>
1203 <p/><strong>PartonSystems Pythia::partonSystems </strong> <br/>
1204 a grouping of the partons in the event record by subsystem,
1205 see <a href="AdvancedUsage.html" target="page">here</a> for further details.
1211 <!-- Copyright (C) 2013 Torbjorn Sjostrand -->