Program Flow

Recall that, to first order, the event generation process can be subdivided into three stages:
  1. Initialization.
  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.

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. (Hypothetically you could use several Pythia objects sequentially, but if done in parallel the static character of some program elements is likely not to give the desired behaviour.)
    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 (by an optional last argument, a reference to an ostream, usually not explicitly written out here).
  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 static classes, and 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 dedicated methods of each class to change the database default values to fit the needs of your current run. However, the

        pythia.readString(string);
    
    method provides a covenient uniform interface to all of them. 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), unless an optional second argument false is used to switch off warnings.
    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 of input files that contain commands to several different 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 ParticleDataTable. Again, an optional second argument false allows you to switch off warning messages for unknown variables.
    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.
  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 [Alw06] 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 [Boo01] 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, but also caused by a failure of Pythia to generate an event, or that an event was generated but something strange was detected in it. 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 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. With optional argument true also further statistics is printed. Currently this means the number of different subprocesses generated in the multiple-interactions framework.

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 to the .xml files, originally stored in 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 (further instructions).

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

  1. PYTHIA has not been written for threadsafe execution. As a rule, you should not have more than one Pythia object in an executable at any time. For 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.

    In some cases it is possible to use more than one Pythia object. The key example would be the simultaneous generation of signal and pileup events, see main19.cc. Here all signal processes must be switched on before the first initialization, and then switched off and replaced by the background ones before the second initialization. Also most other settings can be changed consistently in between the two initializations, but in a few cases the last value will be used. Particle data is always based on the latest information. As a rule, however, it is safer to use two separate runs to store events on disk, in two separate files, and mix afterwards.

  2. 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 a pythia object once 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), which as second argument can take a non-negative subrun number. (Or, alternatively use the longer form pythia.readFile(fileName, warn, subrun).) 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.
  3. 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.