Program Flow

Recall that, to first order, the event generation process can be subdivided into three stages:
  1. Initializaion.
  2. The event loop.
  3. Finishing.
This is reflected in how the top-level Pythia class should be used in the user-supplied main program, further outlined in the following. Since the nature of the run is defined at the initialization stage, this is where most of the PYTHIA user code has to be written. So as not to confuse the reader unduly, the description of initialization options has been subdivided into what would normally be used and what is intended for more special applications.

At the bottom of this webpge is a complete survey of all public Pythia methods and data members, in a more formal style than the task-oriented descriptions found in the preceding sections. This offers complementary information.

Initialization - normal usage

  1. Already at the top of the main program file, you need to include the proper header file
        #include "Pythia.h"
    
    To simplify typing, it also makes sense to declare
        using namespace Pythia8; 
    
  2. The first step is to create a generator object, e.g. with
         Pythia pythia;
    
    It is this object that we will use from now on. Normally a run will only contain one Pythia object. (But you can use several Pythia objects, which then will be independent of each other.)
    By default all output from Pythia will be on the cout stream, but the list methods below do allow output to alternative streams or files.
  3. You next want to set up the character of the run. The pages under the "Setup Run Tasks" heading in the index describe all the options available (with some very few exceptions, found on the other pages). The default values and your modifications are stored in two databases, one for generic settings and one for particle data. Both of these are initialized with their default values by the Pythia constructor. The default values can then be changed, primarily by one of the two ways below, or by a combination of them.

    a) You can use the

        pythia.readString(string);
    
    method repeatedly to do a change of a property at a time. The information in the string is case-insensitive, but upper- and lowercase can be combined for clarity. The rules are that
    (i) if the first nonblank character of the string is a letter it is assumed to contain a setting, and is sent on to pythia.settings.readString(string);
    (ii) if instead the string begins with a digit it is assumed to contain particle data updates, and so sent on to pythia.particleData.readString(string);
    (iii) if none of the above, the string is assumed to be a comment, i.e. nothing will be done.
    In the former two cases, a warning is issued whenever a string cannot be recognized (maybe because of a spelling mistake).
    Some examples would be
        pythia.readString("TimeShower:pTmin = 1.0");
        pythia.readString("111:mayDecay = false");
    
    The readString(string) method is intended primarily for a few changes. It can also be useful if you want to construct a parser for input files that contain commands both to PYTHIA and to other libraries.

    b) You can read in a file containing a list of those variables you want to see changed, with a

        pythia.readFile(fileName);
    
    Each line in this file with be processes by the readString(string) method introduced above. You can thus freely mix comment lines and lines handed on to Settings or to ParticleData.
    This approach is better suited for more extensive changes than a direct usage of readString(string), and can also avoid having to recompile and relink your main program between runs.
    It is also possible to read input from an istream, by default cin, rather than from a file. This may be convenient if information is generated on-the-fly, within the same run.
  4. Next comes the initialization stage, where all remaining details of the generation are to be specified. The init(...) method allows a few different input formats, so you can pick the one convenient for you:

    a) pythia.init( idA, idB, eCM);
    lets you specify the identities and the CM energy of the two incoming beam particles, with A (B) assumed moving in the +z (-z) direction.

    b) pythia.init( idA, idB, eA, eB);
    is similar, but the two beam energies can be different, so the collisions do not occur in the CM frame. If one of the beam energies is below the particle mass you obtain a fixed-target topology.

    c) pythia.init( idA, idB, pxA, pyA, pzA, pxB, pyB, pzB);
    is similar, but here you provide the three-momenta (p_x, p_y, p_z) of the two incoming particles, to allow for arbitrary beam directions.

    d) pythia.init(fileName);
    assumes a file in the Les Houches Event File format is provided.

    e) pythia.init();
    with no arguments will read the beam parameters from the Main group of variables, which provides you with the same possibilities as the above options a, b, c and d. If you don't change any of those you will default to proton-proton collisions at 14 TeV, i.e. the nominal LHC values.

    f) pythia.init( LHAup*);
    assumes Les Houches Accord initialization and event information is available in an LHAup class object, and that a pointer to this object is handed in.

  5. If you want to have a list of the generator and particle data used, either only what has been changed or everything, you can use
        pythia.settings.listChanged();
        pythia.settings.listAll();
        pythia.particleData.listChanged(); 
        pythia.particleData.listAll(); 
    

The event loop

  1. Inside the event generation loop you generate the next event using the next() method,
        pythia.next();
    
    This method takes no arguments; everything has already been specified. It does return a bool value, however, false when the generation failed. This can be a "programmed death" when the supply of input parton-level configurations on file is exhausted. It can alternatively signal a failure of Pythia to generate an event, or unphysical features in the event record at the end of the generation step. It makes sense to allow a few false values before a run is aborted, so long as the related faulty events are skipped.
  2. The generated event is now stored in the event object, of type Event, which is a public member of pythia. You therefore have access to all the tools described on the pages under the "Study Output" header in the index. For instance, an event can be listed with pythia.event.list(), the identity of the i'th particle is given by pythia.event[i].id(), and so on.
    The hard process - roughly the information normally stored in the Les Houches Accord event record - is available as a second object, process, also of type Event.
    A third useful public object is info, which offers a set of one-of-a kind pieces of information about the most recent event.

Finishing

  1. At the end of the generation process, you can call
        pythia.statistics(); 
    
    to get some run statistics, on cross sections and the number of errors and warnings encountered.

Advanced usage, mainly for initialization

A) Necessary data are automatically loaded when you use the default PYTHIA installation directory structure and run the main programs in the examples subdirectory. However, in the general case, you must provide the path of the xmldoc directory, where default settings and particle data are found. This can be done in two ways.
  1. You can set the environment variable PYTHIA8DATA to contain the location of the xmldoc directory. In the csh and tcsh shells this could e.g. be
         setenv PYTHIA8DATA /home/myname/pythia81xx/xmldoc
    
    while in other shells it could be
         export PYTHIA8DATA=/home/myname/pythia81xx/xmldoc
    
    where xx is the subversion number.
    Recall that environment variables set locally are only defined in the current instance of the shell. The above lines should go into your .cshrc and .bashrc files, respectively, if you want a more permanant assignment.
  2. You can provide the path as argument to the Pythia constructor, e.g.
         Pythia pythia("/home/myname/pythia81xx/xmldoc");
    
where again xx is the subversion number.
When PYTHIA8DATA is set it takes precedence, else the path in the constructor is used, else one defaults to the ../xmldoc directory.

B) You can override the default behaviour of PYTHIA not only by the settings and particle data, but also by replacing some of the PYTHIA standard routines by ones of your own. Of course, this is only possible if your routines fit into the general PYTHIA framework. Therefore they must be coded according to the the rules relevant in each case, as a derived class of a PYTHIA base class, and a pointer to such an object must be handed in by one of the methods below. These calls must be made before the pythia.init(...) call.

  1. If you are not satisfied with the list of parton density functions that are implemented internally or available via the LHAPDF interface (see the PDF Selection page), you can suppy your own by a call to the setPDFPtr(...) method
          pythia.setPDFptr( pdfAPtr, pdfBPtr); 
    
    where pdfAPtr and pdfBPtr are pointers to two Pythia PDF objects. Note that pdfAPtr and pdfBPtr cannot point to the same object; even if the PDF set is the same, two copies are needed to keep track of two separate sets of x and density values.
    If you further wish to use separate PDF's for the hard process of an event than the ones being used for everything else, the extended form
          pythia.setPDFptr( pdfAPtr, pdfBPtr, pdfHardAPtr, pdfHardBPtr); 
    
    allows you to specify those separately, and then the first two sets would only be used for the showers and for multiple interactions.
  2. If you want to perform some particle decays with an external generator, you can call the setDecayPtr(...) method
          pythia.setDecayPtr( decayHandlePtr, particles);
    
    where the decayHandlePtr derives from the DecayHandler base class and particles is a vector of particle codes to be handled.
  3. If you want to use an external random number generator, you can call the setRndmEnginePtr(...) method
          pythia.setRndmEnginePtr( rndmEnginePtr); 
    
    where rndmEnginePtr derives from the RndmEngine base class. The Pythia default random number generator is perfectly good, so this is only intended for consistency in bigger frameworks.
  4. If you want to interrupt the evolution at various stages, to interrogate the event and possibly veto it, or you want to reweight the cross section, you can use
          pythia.setUserHooksPtr( userHooksPtr); 
    
    where userHooksPtr derives from the UserHooks base class.
  5. If you want to use your own parametrization of beam momentum spread and interaction vertex, rather than the provided simple Gaussian parametrization (off by default), you can call
          pythia.setBeamShapePtr( beamShapePtr); 
    
    where beamShapePtr derives from the BeamShape base class.
  6. If you want to implement a cross section of your own, but still make use of the built-in phase space selection machinery, you can use
          pythia.setSigmaPtr( sigmaPtr);
    
    where sigmaPtr of type SigmaProcess* is an instance of a class derived from one of the Sigma1Process, Sigma2Process and Sigma3Process base classes in their turn derived from SigmaProcess. This call can be used repeatedly to hand in several different processes.
  7. If your cross section contains the production of a new resonance with known analytical expression for all the relevant partial widths, you can make this resonance available to the program with
          pythia.setResonancePtr( resonancePtr);
    
    where resonancePtr of type ResonanceWidths* is an instance of a class derived from the ResonanceWidths base class. In addition you need to add the particle to the normal particle and decay database. This procedure can be used repeatedly to hand in several different resonances.
  8. If you are a real expert and want to replace the PYTHIA initial- and final-state showers, you can use
          pythia.setShowerPtr( timesDecPtr, timesPtr, spacePtr);
    
    where timesDecPtr and timesPtr derive from the TimeShower base class, and spacePtr from SpaceShower.
  9. The initTunes(...) method allows you to initialize the settings of a specific tune at an early stage, so that you can change some of the tune values between the initTunes(...) and the init(...) calls. This method should not be used if you want to use an existing tune as is, or make your own changes altogether. See the initTunes(...) description below for more details.

C) Some comments on collecting several tasks in the same run.

  1. PYTHIA has not been written for threadsafe execution on multicore processors. If you want to use all cores, the most efficient way presumably is to start correspondingly many jobs, with different random number seeds, and add the statistics at the end.
  2. In some cases it is convenient to use more than one Pythia object. The key example would be the simultaneous generation of signal and pileup events, see main19.cc. The two objects are then set up and initialized separately, and generate events completely independently of each other. It is only afterwards that the event records are combined into one single super-event per beam crossing.
  3. When time is not an issue, it may be that you want to perform several separate subruns sequentially inside a run, e.g. to combine results for several kinematical regions or to compare results for some different tunes of the underlying event. One way to go is to create (and destroy) one pythia object for each subrun, in which case they are completely separate. You can also use the same pythia object, only doing a new init(...) call for each subrun. In that case, the settings and particle databases remain as they were in the previous subrun, only affected by the specific changes you introduced in the meantime. You can put those changes in the main program, with pythia.readString(string), using your own logic to decide which ones to execute in which subrun. A corresponding possibility exists with pythia.readFile(fileName, subrun) (or an istream instead of a fileName), which as second argument can take a non-negative subrun number. Then only those sections of the file before any Main:subrun = ... line or with matching subrun number will be read. That is, the file could have a structure like
        ( lines always read, i.e. "default values" always (re)set )
        Main:subrun = 1
        ( lines only read with readFile(fileName, 1) )
        Main:subrun = 2
        ( lines only read with readFile(fileName, 2) )
    
    Both of these possibilities are illustrated in main08.cc.
  4. When working with Les Houches Event Files, it may well be that your intended input event sample is spread over several files, that you all want to turn into complete events in one and the same run. There is no problem with looping over several subruns, where each new subrun is initialized with a new file. However, in that case you will do a complete re-initialization each time around. If you want to avoid this, note that there is an optional second argument for LHEF initialization: pythia.init(fileName, skipInit). Alternatively, the tag Main:LHEFskipInit can be put in a file of commands to obtain the same effect. Here skipInit defaults to false, but if set true then the new file will be simulated with the same initialization data as already set in a previous pythia.init(...) call. The burden rests on you to ensure that this is indeed correct, e.g. that the two event samples have not been generated for different beam energies.

The Pythia Class

Here follows the complete survey of all public Pythia methods and data members.

Constructor and destructor

Pythia::Pythia(string xmlDir = "../xmldoc")  
creates an instance of the Pythia event generators, and sets initial default values, notably for all settings and particle data. You may use several Pythia instances in the same run; only when you want to access external static libraries could this cause problems. (This includes in particular Fortran libraries such as LHAPDF.)
argument xmlDir (default = ../xmldoc) : allows you to choose from which directory the default settings and particle data values are read in. If the PYTHIA8DATA environment variable has been set it takes precedence. Else this optional argument allows you to choose another directory location than the default one. Note that it is only the directory location you can change, its contents must be the ones of the xmldoc directory in the standard distribution.

Pythia::~Pythia  
the destructor deletes the objects created by the constructor.

Set up run

bool Pythia::readString(string line, bool warn = true)  
reads in a single string, that is interpreted as an instruction to modify the value of a setting or particle data, as already described above.
argument line : the string to be interpreted as an instruction.
argument warn (default = true) : write a warning message or not whenever the instruction does not make sense, e.g. if the variable does not exist in the databases.
Note: the method returns false if it fails to make sense out of the string.

bool Pythia::readFile(string fileName, bool warn = true, int subrun = SUBRUNDEFAULT)  
bool Pythia::readFile(string fileName, int subrun = SUBRUNDEFAULT)  
bool Pythia::readFile(istream& inStream = cin, bool warn = true, int subrun = SUBRUNDEFAULT)  
bool Pythia::readFile(istream& inStream = cin, int subrun = SUBRUNDEFAULT)  
reads in a whole file, where each line is interpreted as an instruction to modify the value of a setting or particle data, cf. the above readString method. All four forms of the readFile command share code for actually reading a file.
argument fileName : the file from which instructions are read.
argument inStream : an istream from which instructions are read.
argument warn (default = true) : write a warning message or not whenever the instruction does not make sense, e.g. if the variable does not exist in the databases. In the command forms where warn is omitted it is true.
argument subrun : allows you have several optional sets of commands within the same file. Only those sections of the file before any Main:subrun = ... line or following such a line with matching subrun number will be read. The subrun number should not be negative; negative codes like SUBRUNDEFAULT corresponds to no specific subrun.
Note: the method returns false if it fails to make sense out of any one line.

bool Pythia::setPDFPtr( PDF* pdfAPtr, PDF* pdfBPtr, PDF* pdfHardAPtr = 0, PDF* pdfHardBPtr = 0)  
offers the possibility to link in external PDF sets for usage inside the program. The rules for constructing your own class from the PDF base class are described here.
argument pdfAPtr, pdfBPtr : pointers to two PDF-derived objects, one for each of the incoming beams. The two objects have to be instantiated by you in your program. Even if the two beam particles are the same (protons, say) two separate instances are required, since current information is cached in the objects. If both arguments are zero then any previous linkage to external PDF's is disconnected, see further Note 2 below.
argument pdfHardAPtr, pdfHardBPtr (default = 0) : pointers to two further PDF-derived objects, one for each of the incoming beams. Normally only the first two arguments above would be used, and then the same PDF sets would be invoked everywhere. If you provide these two further pointers then two different sets of PDF's are used. This second set is then exclusively for the generation of the hard process from the process matrix elements library. The first set above is for everything else, notably parton showers and multiple interactions.
Note 1: The method returns false if the input is obviously incorrect, e.g. if two (nonzero) pointers agree.
Note 2: If you want to combine several subruns you can call setPDFPtr with new arguments before each Pythia::init(...) call. To revert from external PDF's to the normal internal PDF selection you must call setPDFPtr(0, 0) before Pythia::init(...).

bool Pythia::setDecayPtr( DecayHandler* decayHandlePtr, vector<int> handledParticles)  
offers the possibility to link to an external program that can do some of the particle decays, instad of using the internal decay machinery. With particles we here mean the normal hadrons and leptons, not top quarks, electroweak bosons or new particles in BSM scenarios. The rules for constructing your own class from the DecayHandler base class are described here. Note that you can only provide one external object, but this object in its turn could very well hand on different particles to separate decay libraries.
argument decayHandlePtr : pointer to a DecayHandler-derived object. This object must be instantiated by you in your program.
argument handledParticles : vector with the PDG identity codes of the particles that should be handled by the external decay package. You should only give the particle (positive) codes; the respective antiparticle is always included as well.
Note: The method currently always returns true.

bool Pythia::setRndmEnginePtr( RndmEngine* rndmEnginePtr)  
offers the possibility to link to an external random number generator. The rules for constructing your own class from the RndmEngine base class are described here.
argument rndmEnginePtr : pointer to a RndmEngine-derived object. This object must be instantiated by you in your program.
Note: The method returns true if the pointer is different from 0.

bool Pythia::setUserHooksPtr( UserHooks* userHooksPtr)  
offers the possibility to interact with the generation process at a few different specified points, e.g. to reject undesirable events at an early stage to save computer time. The rules for constructing your own class from the UserHooks base class are described here. You can only hand in one such pointer, but this may be to a class that implements several of the different allowed possibilities.
argument userHooksPtr : pointer to a userHooks-derived object. This object must be instantiated by you in your program.
Note: The method currently always returns true.

bool Pythia::setBeamShapePtr( BeamShape* beamShapePtr)  
offers the possibility to provide your own shape of the momentum and space-time spread of the incoming beams. The rules for constructing your own class from the BeamShape base class are described here.
argument BeamShapePtr : pointer to a BeamShape-derived object. This object must be instantiated by you in your program.
Note: The method currently always returns true.

bool Pythia::setSigmaPtr( SigmaProcess* sigmaPtr)  
offers the possibility to link your own implementation of a process and its cross section, to make it a part of the normal process generation machinery, without having to recompile the Pythia library itself. The rules for constructing your own class from the SigmaProcess base class are described here. You may call this routine repeatedly, to add as many new processes as you wish.
argument sigmaPtr : pointer to a SigmaProcess-derived object. This object must be instantiated by you in your program.
Note: The method currently always returns true.

bool Pythia::setResonancePtr( ResonanceWidths* resonancePtr)  
offers the possibility to link your own implementation of the calculation of partial resonance widths, to make it a part of the normal process generation machinery, without having to recompile the Pythia library itself. This allows the decay of new resonances to be handled internally, when combined with new particle data. Note that the decay of normal hadrons cannot be modelled here; this is for New Physics resonances. The rules for constructing your own class from the ResonanceWidths base class are described here. You may call this routine repeatedly, to add as many new resonances as you wish.
argument resonancePtr : pointer to a ResonanceWidths-derived object. This object must be instantiated by you in your program.
Note: The method currently always returns true.

bool Pythia::setShowerPtr( TimeShower* timesDecPtr, TimeShower* timesPtr = 0, SpaceShower* spacePtr = 0)  
offers the possibility to link your own parton shower routines as replacements for the default ones. This is much more complicated since the showers are so central and are so interlinked with other parts of the program. Therefore it is also possible to do the replacement in stages, from the more independent to the more intertwined. The rules for constructing your own classes from the TimeShower and SpaceShowerbase classes are described here. These objects must be instantiated by you in your program.
argument timesDecPtr : pointer to a TimeShower-derived object for doing timelike shower evolution in resonance decays, e.g. of a Z^0. This is decoupled from beam remnants and parton distributions, and is therefore the simplest kind of shower to write. If you provide a value 0 then the internal shower routine will be used.
argument timesPtr (default = 0) : pointer to a TimeShower-derived object for doing all other timelike shower evolution, which is normally interleaved with multiple interactions and spacelike showers, introducing both further physics and further technical issues. If you retain the default value 0 then the internal shower routine will be used. You are allowed to use the same pointer as above for the timesDecPtr if the same shower can fulfill both tasks.
argument spacePtr (default = 0) : pointer to a SpaceShower-derived object for doing all spacelike shower evolution, which is normally interleaved with multiple interactions and timelike showers. If you retain the default value 0 then the internal shower routine will be used.
Note: The method currently always returns true.

void initTunes(int eeTune = 0, int ppTune = 0)  
prepackages the parameter settings of some PYTHIA tunes. While the readString(...) and readFile(...) methods can be used to set the parameters of some already existing tune, one by one, this is time-consuming and error-prone. Therefore the initTunes(...) method provides a set of tunes, hopefully growing with time, that can be selected by the two switches Tune:ee and Tune:pp. Default is 0, meaning no changes. Nornally initTunes(...) is automatically called from inside the init(...) method, so you would not call it yourself. This can be constraining if you want to use almost, but not quite, an existing tune. Specifically, any attempts to set one or a few of the parameters differently by readString(...) or readFile(...) calls before initialization would not work, since they would be overwritten during init(...). This is where the possibility to call initTunes(...) directly comes to the rescue. You could then by hand pick the desired tune, afterwards change whatever settings you wish, and retain these changes during initialization by having Tune:ee = Tune:pp = 0.
argument eeTune (default = 0) : the settings obtained from e^+e^- data, see Tune:ee for a list of valid values. If you call initTunes(...) with this option nonzero you would want to keep Tune:ee equal to zero, or else the whole point is lost.
argument eeTune (default = 0) : the additional settings that can only be extracted from pp/ppbar data see Tune:pp for a list of valid values. Comment as above.

Initialize

At the initialization stage all the information provided above is processed, and the stage is set up for the subsequent generation of events. Several alterative forms of the init method are available for this stage; pick the one most convenient.

bool Pythia::init( int idA, int idB, double eCM)  
initialize for collisions in the center-of-mass frame, with the beams moving in the +-z directions.
argument idA, idB : particle identity code for the two incoming beams.
argument eCM : the CM energy of the collisions.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

bool Pythia::init( int idA, int idB, double eA, double eB)  
initialize for collisions with back-to-back beams, moving in the +-z directions, but with different energies.
argument idA, idB : particle identity code for the two incoming beams.
argument eA, eB : the energies of the two beams. If an energy is set to be below the mass of the respective beam particle that particle is taken to be at rest. This offers a simple possibility to simulate fixed-target collisions.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

bool Pythia::init( int idA, int idB, double pxA, double pyA, double pzA, double pxB, double pyB, double pzB)  
initialize for collisions with arbitrary beam directions.
argument idA, idB : particle identity code for the two incoming beams.
argument pxA, pyA, pzA : the three-momntum vector (p_x, p_y, p_z) of the first incoming beam.
argument pxB, pyB, pzB : the three-momntum vector (p_x, p_y, p_z) of the second incoming beam.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

bool Pythia::init( string LesHouchesEventFile, bool skipInit = false)  
initialize for hard-process collisions fed in from an external file with events, written according to the Les Houches Event File standard.
argument LesHouchesEventFile : the file name (including path, where required) where the events are stored, including relevant information on beam identities and energies.
argument skipInit (default = false) : By default this method does a complete reinitialization of the generation process. If you set this argument to true then no reinitialization will occur, only the pointer to the event file is updated. This may come in handy if the full event sample is split across several files generated under the same conditions (except random numbers, of course). You then do the first initialization with the default, and all subsequent ones with true. Note that things may go wrong if the files are not created under the same conditions.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

bool Pythia::init()  
initialize for collisions, in any of the four above possibilities. In this option the beams are not specified by input arguments, but instead by the settings in the Beam Parameters section. This allows the beams to be specified in the same file as other run instructions. The default settings give pp collisions at 14 TeV.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

bool Pythia::init( LHAup* lhaUpPtr)  
initialize for hard-process collisions fed in from an external source of events, consistent with the Les Houches Accord standard. The rules for constructing your own class from the LHAup base class are described here. This class is also required to provide the beam parameters.
argument lhaUpPtr : pointer to a LHAup-derived object. This object must be instantiated by you in your program.
Note: The method returns false if the initialization fails. It is then not possible to generate any events.

Generate events

The next() method is the main one to generate events. In this section we also put a few other specialized methods that may be useful in some circumstances.

bool Pythia::next()  
generate the next event. No input parameters are required; all instructions have already been set up in the initialization stage.
Note: The method returns false if the event generation fails. The event record is then not consistent and should not be studied. When reading in hard collisions from a Les Houches Event File the problem may be that the end of the file has been reached. This can be checked with the Info::atEndOfFile() method.

bool Pythia::forceHadronLevel()  
hadronize the existing event record, i.e. perform string fragmentation and particle decays. There are two main applications. Firstly, you can use the same parton-level content as a basis for repeated hadronization attempts, in schemes intended to save computer time. Secondly, you may have an external program that can simulate the full partonic level of the event - hard process, parton showers, multiple interactions, beam remnants, colour flow, and so on - but not hadronization. Further details are found here.
Note: The method returns false if the hadronization fails. The event record is then not consistent and should not be studied.

bool Pythia::moreDecays()  
perform decays of all particles in the event record that have not been decayed but should have been done so. This can be used e.g. for repeated decay attempts, in schemes intended to save computer time. Further details are found here.
Note: The method returns false if the decays fail. The event record is then not consistent and should not be studied.

void Pythia::LHAeventList(ostream& os = cout)  
list the Les Houches Accord information on the current event, see LHAup::listEvent(...). (Other listings are available via the class members below, so this listing is a special case that would not fit elsewhere.)
argument os (default = cout) : output stream where the listing occurs.

bool Pythia::LHAeventSkip(int nSkip)  
skip ahead a number of events in the Les Houches generation sequence, without doing anything further with them, see LHAup::skipEvent(nSkip). Mainly intended for debug purposes, e.g. when an event at a known location in a Les Houches Event File is causing problems.
argument nSkip : number of events to skip.
Note: The method returns false if the operation fails, specifically if the end of a LHEF has been reached, cf. next() above.

Finalize

There is no required finalization step; you can stop generating events when and how you want. It is still recommended that you make it a routine to call the following method at the end.

void Pythia::statistics(bool all = false, bool reset = false)  
list statistics on the event generation, specifically total and partial cross sections and the number of different errors. For more details see here.
argument all (default = false) : if true also statistics on multiple interactions is shown, by default not.
argument reset (default = false) : if true then all counters, e.g on events generated and errors experienced, are reset to zero whenever the routine is called. The default instead is that all stored statistics information is unaffected by the call. Counters are automatically reset in each new Pythia::init(...) call, however, so the only time the reset option makes a difference is if statistics(...) is called several times in a (sub)run.

Interrogate settings

Normally settings are used in the setup and initialization stages to determine the character of a run, e.g. read from a file with the above-described Pythia::readFile(...) method. There is no strict need for a user to interact with the Settings database in any other way. However, as an option, some settings variables have been left free for the user to set in such a file, and then use in the main program to directly affect the performance of that program, see here. A typical example would be the number of events to generate. For such applications the following shortcuts to some Settings methods may be convenient.

bool Pythia::flag(string key)  
read in a boolean variable from the Settings database.
argument key : the name of the variable to be read.

int Pythia::mode(string key)  
read in an integer variable from the Settings database.
argument key : the name of the variable to be read.

double Pythia::parm(string key)  
read in a double-precision variable from the Settings database.
argument key : the name of the variable to be read.

string Pythia::word(string key)  
read in a string variable from the Settings database.
argument key : the name of the variable to be read.

Data members

The Pythia class contains a few public data members, several of which play a central role. We list them here, with links to the places where they are further described.

Event Pythia::process  
the hard-process event record, see here for further details.

Event Pythia::event  
the complete event record, see here for further details.

Info Pythia::info  
further information on the event-generation process, see here for further details.

Settings Pythia::settings  
the settings database, see here for further details.

ParticleData Pythia::particleData  
the particle properties and decay tables database, see here for further details.

Rndm Pythia::rndm  
the random number generator, see here and here for further details.

CoupSM Pythia::coupSM  
Standard Model couplings and mixing matrices, see here for further details.

SusyLesHouches Pythia::slha  
parameters and particle data in the context of supersymmetric models, see here for further details.

PartonSystems Pythia::partonSystems  
a grouping of the partons in the event record by subsystem, see here for further details.