Particle Data Scheme

  1. Databases
  2. Stored properties for particles
  3. Stored properties for decays
  4. Operation
  5. The public methods
The particle data scheme may take somewhat longer to understand than the settings one. In particular the set of methods to access information is rather more varied, to allow better functionality for advanced usage. However, PYTHIA does come with a sensible default set of particle properties and decay tables. Thus there is no need to learn any of the methods on this page to get going. Only when you perceive a specific need does it make sense to learn the basics.

The central section on this page is the Operation one. The preceding sections are there mainly to introduce the basic structure and the set of properties that can be accessed. The subsequent sections provide a complete listing of the existing public methods, which most users probably will have little interaction with.

Databases

The management of particle data is based on three classes: The objects of these classes together form a database that is continuously being used as the program has to assign particle masses, select decay modes, etc.

Each Pythia object has a public member particleData of the ParticleData class. Therefore you access the particle data methods as pythia.particleData.command(argument), assuming that pythia is an instance of the Pythia class. Further, for some of the most frequent user tasks, Pythia methods have been defined, so that pythia.command(argument) would work, see further below.

A fundamental difference between the particle data classes and the settings ones is that the former are accessed regularly during the event generation process, as a new particle is produced and its mass need to be set, e.g., while the latter are mainly/only used at the initialization stage. Nevertheless, it is not a good idea to change data in either of them in mid-run, since this may lead to inconsistencies.

Stored properties for particles

The main properties stored for each particle are as follows. Different ways to set and get these properties will be described further down.

Stored properties for decays

An unstable particle has a decay table consisting of one or more decay channel. The following properties are stored for each such channel. Again different ways to set and get these properties will be described further down.

Operation

The normal flow of the particle data operations is:
  1. When a Pythia object pythia is created, the pythia.particleData member is asked to scan the ParticleData.xml file.

    All lines beginning with <particle are scanned for information on a particle species, and all lines beginning with <channel are assumed to contain a decay channel of the enclosing particle. In both cases XML syntax is used, with attributes used to identify the stored properties, and with omitted properties defaulting back to 0 where meaningful. The particle and channel information may be split over several lines, up to the > endtoken. The format of a <particle tag is:

     
        <particle id="..." name="..." antiName="..." spinType="..." chargeType="..." colType="..." 
           m0="..." mWidth="..." mMin="..." mMax="..." tau0="..."> 
        </particle> 
    
    where the fields are the properties already introduced above. Note that tauCalc, isResonance, mayDecay, doExternalDecay, isVisible and doForceWidth are not set here, but are provided with default values by the rules described above. Once initialized, also these latter properties can be changed, see below.
    The format of a <channel> tag is:
     
        <channel onMode="..." bRatio="..." meMode="..." products="..." /> 
    
    again see properties above. The products are given as a blank-separated list of id codes.
    Important: the values in the .xml file should not be changed, except by the PYTHIA authors. Any changes should be done with the help of the methods described below.
  2. Between the creation of the Pythia object and the init call for it, you may use the methods of the ParticleData class to modify some of the default values. Several different approaches can be chosen for this.

    a) Inside your main program you can directly set values with

     
        pythia.readString(string); 
    
    where both the variable name and the value are contained inside the character string, separated by blanks and/or a =, e.g.
     
        pythia.readString("111:mayDecay = off"); 
    
    switches off the decays of the pi^0.
    The particle id (> 0) and the property to be changed must be given, separated by a colon.
    The allowed properties are: name, antiName, spinType, chargeType, colType, m0, mWidth, mMin, mMax, tau0, tauCalc, isResonance, mayDecay, doExternalDecay, isVisible and doForceWidth. All of these names are case-insensitive. Names that do not match an existing variable are ignored.
    Strings beginning with a non-alphanumeric character, like # or !, are assumed to be comments and are not processed at all. For bool values, the following notation may be used interchangeably: true = on = yes = ok = 1, while everything else gives false (including but not limited to false, off, no and 0).

    Particle data often comes in sets of closely related information. Therefore some properties expect the value to consist of several numbers. These can then be separated by blanks (or by commas). A simple example is names, which expects both the name and antiname to be given. A more interesting one is the all property,

     
        id:all = name antiName spinType chargeType colType m0 mWidth mMin mMax tau0 
    
    where all the current information on the particle itself is replaced, but any decay channels are kept unchanged. Using new instead of all also removes any previous decay channels. If the string contains fewer fields than expected the trailing properties are set to vanish ("void", 0 or 0.). Note that such a truncated string should not be followed by a comment, since this comment would then be read in as if it contained the missing properties. The truncation can be done anywhere, specifically a string with only id:new defines a new "empty" particle. As before, tauCalc, isResonance, mayDecay, doExternalDecay, isVisible and doForceWidth are (re)set to their default values, and would have to be changed separately if required.

    A further command is rescaleBR, which rescales each of the existing branching ratios with a common factor, such that their new sum is the provided value. This may be a first step towards adding new decay channels, see further below.

    Alternatively the id code may be followed by another integer, which then gives the decay channel number. This then has to be followed by the property specific to this channel, either onMode, bRatio, meMode or products. In the latter case all the products of the channel should be given:

     
        id:channel:products =  product1 product2 .... 
    
    The line will be scanned until the end of the line, or until a non-number word is encountered, or until the maximum allowed number of eight products is encountered, whichever happens first. (Thus the multiplicity of a decay channel need not be input; it is automatically calculated from the products list.) It is also possible to replace all the properties of a channel in a similar way:
     
        id:channel:all = onMode bRatio meMode product1 product2 .... 
    
    To add a new channel at the end, use
     
        id:addChannel = onMode bRatio meMode product1 product2 .... 
    

    It is currently not possible to remove a channel selectively, but setting its branching ratio vanishing is as effective. If you want to remove all existing channels and force decays into one new channel you can use

     
        id:oneChannel = onMode bRatio meMode product1 product2 .... 
    
    A first oneChannel command could be followed by several subsequent addChannel ones, to build up a completely new decay table for an existing particle.

    When adding new channels or changing branching ratios in general, note that, once a particle is to be decayed, the sum of branching ratios is always rescaled to unity. Beforehand, rescaleBR may be used to rescale an existing branching ratio by the given factor.

    There are a few commands that will study all the decay channels of the given particle, to switch them on or off as desired. The

     
        id:onMode = onMode 
    
    will set the onMode property of all channels to the desired value. The
     
        id:offIfAny   = product1 product2 .... 
        id:onIfAny    = product1 product2 .... 
        id:onPosIfAny = product1 product2 .... 
        id:onNegIfAny = product1 product2 .... 
    
    will set the onMode 0, 1, 2 or 3, respectively, for all channels which contain any of the enumerated products, where the matching to these products is done without distinction of particles and antiparticles. Note that "Pos" and "Neg" are slightly misleading since it refers to the particle and antiparticle of the id species rather than charge, but should still be simpler to remember and understand than alternative notations. Correspondingly
     
        id:offIfAll   = product1 product2 .... 
        id:onIfAll    = product1 product2 .... 
        id:onPosIfAll = product1 product2 .... 
        id:onNegIfAll = product1 product2 .... 
    
    will set the onMode 0, 1, 2 or 3, respectively, for all channels which contain all of the enumerated products, again without distinction of particles and antiparticles. If the same product appears twice in the list it must also appear twice in the decay channel, and so on. The decay channel is allowed to contain further particles, beyond the product list. By contrast,
     
        id:offIfMatch   = product1 product2 .... 
        id:onIfMatch    = product1 product2 .... 
        id:onPosIfMatch = product1 product2 .... 
        id:onNegIfMatch = product1 product2 .... 
    
    requires the decay-channel multiplicity to agree with that of the product list, but otherwise works as the onIfAll/offIfAll methods.

    Note that the action of several of the commands depends on the order in which they are executed, as one would logically expect. For instance, id:oneChannel removes all decay channels of id and thus all previous changes in this decay table, while subsequent additions or changes would still take effect. Another example would be that 23:onMode = off followed by 23:onIfAny = 1 2 3 4 5 would let the Z^0 decay to quarks, while no decays would be allowed if the order were to be reversed.

    b) The Pythia readString(string) method actually does not do changes itself, but sends on the string either to the ParticleData class or to the Settings one, depending on whether the string begins with a digit or a letter. If desired, it is possible to communicate directly with the corresponding ParticleData method:

     
        pythia.particleData.readString("111:mayDecay = off"); 
        pythia.particleData.readString("15:2:products = 16 -211"); 
    
    In this case, changes intended for Settings would not be understood.

    c) Underlying this are commands for all the individual properties in the ParticleData class, one for each. These are further described below. Thus, an example now reads

     
        pythia.particleData.mayDecay(111, false); 
    
    Boolean values should here be given as true or false.

    d) A simpler and more useful way is to collect all your changes in a separate file, with one line per change, e.g.

     
        111:mayDecay = off 
    
    The file can be read by the
     
        pythia.readFile(fileName); 
    
    method, where fileName is a string, e.g. pythia.readFile("main.cmnd") (or an istream instead of a fileName). Each line is processed as described for the string in 2a). This file can freely mix commands to the Settings and ParticleData classes.
  3. A routine reInit(fileName) is provided, and can be used to zero the particle data table and reinitialize it from scratch. Such a call might be useful if several subruns are to be made with widely different particle data - normally the maps are only built from scratch once, namely when the Pythia() object is created. Also, there is no other possibility to restore the default values, unlike for the settings.

  4. You may at any time obtain a listing of all the particle data by calling

     
        pythia.particleData.listAll(); 
    
    The listing is by increasing id number. It shows the basic quantities introduced above. Some are abbreviated in the header to fit on the lines: spn = spinType, chg = chargeType, col = colType, res = isResonance, dec = mayDecay && canDecay (the latter checks that decay channels have been defined), ext = doExternalDecay, vis = isVisible and wid = doForceWidth.
    To list only those particles that were changed (one way or another, the listing will not tell what property or decay channel was changed), instead use
     
        pythia.particleData.listChanged(); 
    
    (This info is based on a further hasChanged flag of a particle or a channel, set true whenever any of the changing methods are used. It is possible to manipulate this value, but this is not recommended.) By default the internal initialization of the widths of resonances such as gamma^*/Z^0, W^+-, t/tbar, H^0 do not count as changes; if you want to list also those changes instead call listChanged(true).
    To list only one particle, give its id code as argument to the list(...) function.. To list a restricted set of particles, give in their id codes to list(...) as a vector<int>.
  5. For wholesale changes of particle properties all available data can be written out, edited, and then read back in again. These methods are mainly intended for expert users. You can choose between two alternative syntaxes.

    a) XML syntax, using the <particle and <channel lines already described. You use the method particleData.listXML(fileName) to produce such an XML file and particleData.readXML(fileName) to read it back in after editing.

    b) Fixed/free format, using exactly the same information as illustrated for the <particle and <channel lines above, but now without any tags. This means that all information fields must be provided (if there is no antiparticle then write void), in the correct order (while the order is irrelevant with XML syntax), and all on one line. Information is written out in properly lined-up columns, but the reading is done using free format, so fields need only be separated by at least one blank. Each new particle is supposed to be separated by (at least) one blank line, whereas no blank lines are allowed between the particle line and the subsequent decay channel lines, if any. You use the method particleData.listFF(fileName) to produce such a fixed/free file and particleData.readFF(fileName) to read it back in after editing.

    As an alternative to the readXML and readFF methods you can also use the particleData.reInit(fileName, xmlFormat) method, where xmlFormat = true (default) corresponds to reading an XML file and xmlFormat = false to a fixed/free format one.

    To check that the new particle and decay tables makes sense, you can use the particleData.checkTable() method, either directly or by switching it on among the standard error checks.



The public methods

In the following we present briefly the public methods in the three classes used to build up the particle database. The order is top-down, i.e from the full table of all particles to a single particle to a single channel. Note that these methods usually are less elegant and safe than the input methods outlined above. If you use any of these methods, it is likely to be the ones in the full database, i.e. the first ones to be covered in the following.

For convenience, we have grouped related input and output methods together. It should be obvious from the context which is which: the input is of type void and has an extra last argument, namely is the input value, while the output method returns a quantity of the expected type.

The ParticleData methods

ParticleData::ParticleData()  
the constructor has no arguments and does not do anything. Internal.

ParticleData& operator=( const ParticleData& particleDataIn)  
copy the database from an existing ParticleData object. Can be useful when running with multiple Pythia instances. Does not include the links to the resonance decay handling set up by initWidths(...).

void ParticleData::initPtr(Info* infoPtr, Settings* settingsPtrIn, Rndm* rndmPtrIn, CoupSM* coupSMPtrIn)  
initialize pointers to a few other classes. Internal.

bool ParticleData::init(string startFile = "../share/Pythia8/xmldoc/ParticleData.xml")  
read in an XML-style file with particle data and initialize the particle data tables accordingly. This command is executed in the Pythia constructor, i.e. is mainly for internal use.
argument startFile : the name of the data file to be read. When called from the Pythia constructor the directory is provided by the PYTHIA8DATA environment variable, if set, else by the argument of this constructor, which has the default value "../share/Pythia8/xmldoc".

bool ParticleData::init(const ParticleData& particleDataIn  
copy particle data from information stored in an existing ParticleData instance.

bool ParticleData::init(istream& is  
copy particle data from a stream (rather than from a file).

bool ParticleData::reInit(string startFile, bool xmlFormat = true)  
overwrite the existing database by reading from the specified file. Unlike init above this method is not called by the Pythia constructor, but is entirely intended for users who want to replace the existing particle data with their own.
argument startFile : the path and name of file to be read.
argument xmlFormat : if true read the same kind of XML-style file as used by init, if not use an alternative "free format" file (i.e. without any XML tags, but with well-defined rules specifying in which order properties are stored).

void ParticleData::initWidths( vector<ResonanceWidthsPtr> resonancePtrs)  
initialize Breit-Wigner shape parameters for all particles, and the detailed handling of resonances, i.e. particles with perturbatively calculable partial widths, which can be used to obtain a mass-dependent Breit-Wigner and a dynamic choice of decay channels. Called from Pythia::init().

bool ParticleData::readXML(string inFile, bool reset = true)  
void ParticleData::listXML(string outFile)  
read in XML-style data from a file or write it out to a file. For the former one can also decide whether to reset all particles to scratch, or only overwrite those particles in the file. The former method is used by init and reInit above.

bool ParticleData::readFF(string inFile, bool reset = true)  
void ParticleData::listFF(string outFile)  
read in free-format-style data from a file or write it out to a file. For the former one can also decide whether to reset all particles to scratch, or only overwrite those particles in the file. The former method is used by reInit above.

bool ParticleData::readString(string line, bool warn = true)  
read in a string and interpret is as a new or changed particle data. The possibilities are extensively described above. It is normally used indirectly, via Pythia::readString(...) and Pythia::readFile(...).
argument line : the string to be interpreted as an instruction.
argument warn (default = on) : write a warning message or not whenever the instruction does not make sense, e.g. if the particle does not exist in the database.
Note: the method returns false if it fails to make sense out of the input string.

void ParticleData::listAll()  
void ParticleData::listChanged(bool changedRes = false)  
void ParticleData::list(bool changedOnly = false, bool changedRes = true)  
void ParticleData::listAll(ostream& stream)  
void ParticleData::list(ostream& stream, bool changedOnly = false, bool changedRes = true)  
methods intended to present a listing of particle data in a readable format. The first two are special cases of the third. The first lists all particle data, the second only data for those particles that were changed after the original creation of the particle data table. Resonances are a special case since they can get their data changed by being linked to an object that does the calculation of branching ratios. By default the second method does not count such resonances as changed, whereas the third does and thus lists all resonances. The last two methods are simple variants where the output can be redirected.

void ParticleData::list(int idList)  
void ParticleData::list(vector<int> idList)  
list particle data for one single particle, with the identity code as input, or for a set of particles, with an input vector of identity codes.

vector<string> ParticleData::getReadHistory()  
Method to retrieve the history of readString commands that have been processed by the ParticleData instance, e.g. for inspection. Note that ParticleData commands read by Pythia::readString and Pythia::readFile are interpreted by readString and thus also are listed.

vector<string> ParticleData::getReadHistory(int subrun)  
Method to retrieve the history of readString commands that have been processed by the ParticleData instance, for a specific subrun (see the section on Main-Program Settings). For subrun = -1, returns the readString history common to all subruns. For subrun >= 0, returns the history of readString commands for that specific subrun (omitting the common part).

void ParticleData::checkTable(int verbosity)  
check that the particle decay table makes sense, especially for decays.
argument verbosity : level of checks. 0 is only minimal, e.g. if a particle has no open decay channels. 1, which is the level of the first method, provides warning if any individual channel is closed, except for resonances. 2 also prints the branching-ratio-averaged threshold mass. 11 and 12 are like 1 and 2, but also include resonances in the detailed checks.

void ParticleData::addParticle(int id, string name = " ", int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0., double tau0 = 0.)  
void ParticleData::addParticle(int id, string name, string antiName, int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0., double tau0 = 0.)  
add a particle to the decay table; in the first form a particle which is its own antiparticle, in the second where a separate antiparticle exists.

void ParticleData::setAll(int id, string name, string antiName, int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0.,double tau0 = 0.)  
change all the properties of the particle associated with a given identity code.

bool ParticleData::isParticle(int id)  
query whether the particle data table contains the particle of the identity code.

ParticleDataEntry* ParticleData::findParticle(int id)  
const ParticleDataEntry* ParticleData::findParticle(int id)  
query whether the particle data table contains the particle of the identity code, and if so return a (const) iterator to it.

int ParticleData::nextId(int id)  
return the identity code of the sequentially next particle stored in table.

map<int, ParticleDataEntry>::iterator begin()  
map<int, ParticleDataEntry>::iterator end()  
iterators over all the entries in the particle data table.

bool ParticleData::hasAnti(int id)  
bool whether a distinct antiparticle exists or not. Is true if an antiparticle name has been set (and is different from void).

int ParticleData::antiId(int id)  
gets the id of the antiparticle of this particle, which is -id if a distinct antiparticle exists, and id otherwise.

void ParticleData::name(int id, string name)  
void ParticleData::antiName(int id, string antiName)  
void ParticleData::names(int id, string name, string antiName)  
string ParticleData::name(int id)  
particle and antiparticle names are stored separately, the sign of id determines which of the two is returned, with void used to indicate the absence of an antiparticle.

void ParticleData::spinType(int id, int spinType)  
int ParticleData::spinType(int id)  
the spin type, of the form 2 s + 1, with special code 0 for entries of unknown or indeterminate spin.

void ParticleData::chargeType(int id, int chargeType)  
int ParticleData::chargeType(int id)  
three times the charge (to make it an integer), taking into account the sign of id.

double ParticleData::charge(int id)  
the electrical charge of a particle, equal to chargeType(id)/3.

void ParticleData::colType(int id, int colType)  
int ParticleData::colType(int id)  
the colour type, with 0 uncoloured, 1 triplet, -1 antitriplet and 2 octet, taking into account the sign of id.

void ParticleData::m0(int id, double m0)  
double ParticleData::m0(int id)  
the nominal mass m_0 (in GeV).

void ParticleData::mWidth(int id, double mWidth)  
double ParticleData::mWidth(int id)  
the width Gamma of the Breit-Wigner distribution (in GeV).

void ParticleData::mMin(int id, double mMin)  
double ParticleData::mMin(int id)  
the lower limit of the allowed mass range generated by the Breit-Wigner (in GeV). Has no meaning for particles without width, and would typically be 0 there.

void ParticleData::mMax(int id, double mMax)  
double ParticleData::mMax(int id)  
the upper limit of the allowed mass range generated by the Breit-Wigner (in GeV). If mMax < mMin then no upper limit is imposed. Has no meaning for particles without width, and would typically be 0 there.

double ParticleData::m0Min(int id)  
similar to mMin() above, except that for particles with no width the m0(id) value is returned.

double ParticleData::m0Max(int id)  
similar to mMax() above, except that for particles with no width the m0(id) value is returned.

void ParticleData::tau0(int id, double tau0)  
double ParticleData::tau0(int id)  
the nominal proper lifetime tau_0 (in mm/c).

void ParticleData::tauCalc(int id, bool tauCalc)  
bool ParticleData::tauCalc(int id)  
a flag telling whether the nominal proper lifetime tau_0 should be calculated from the particle width, mWidth, at initialization (the default behavior). If set false, then the proper lifetime can be varied independently from the particle width. This capability is useful for simulating models of long-lived particles where the dependence of a signature on the particle width is negligible.

void ParticleData::isResonance(int id, bool isResonance)  
bool ParticleData::isResonance(int id)  
a flag telling whether a particle species are considered as a resonance or not. Here "resonance" is used as shorthand for any massive particle where the decay process should be counted as part of the hard process itself, and thus be performed before showers and other event aspects are added. Restrictions on allowed decay channels is also directly reflected in the cross section of simulated processes, while those of normal hadrons and other light particles are not. In practice, it is reserved for states above the b bbar bound systems in mass, i.e. for W, Z, t, Higgs states, supersymmetric states and (most?) other states in any new theory. All particles with m0 above 20 GeV are by default initialized to be considered as resonances.

void ParticleData::mayDecay(int id, bool mayDecay)  
bool ParticleData::mayDecay(int id)  
a flag telling whether a particle species may decay or not, offering the main user switch. Whether a given particle of this kind then actually will decay also depends on it having allowed decay channels, and on other flags for particle decays (or resonance decays). All particles with tau0 below 1000 mm are by default initialized to allow decays.

void ParticleData::doExternalDecays(int id, bool doExternalDecays)  
bool ParticleData::doExternalDecay(int id)  
a flag telling whether a particle should be handled by an external decay package or not, with the latter default. Can be manipulated as described on this page, but should normally not be. Instead the pythia.decayPtr method should be provided with the list of relevant particles.

void ParticleData::isVisible(int id, bool isVisible)  
bool ParticleData::isVisible(int id)  
a flag telling whether a particle species is to be considered as visible in a detector or not, as used e.g. in analysis routines. By default the invisibles include neutrinos, Dark Matter particles (codes 51 - 60) and a few BSM particles (gravitino, sneutrinos, neutralinos) that have neither strong nor electromagnetic charge, and are not made up of constituents that have it. The value of this flag is only relevant if a particle is long-lived enough actually to make it to a detector.

void ParticleData::doForceWidth(int id, bool doForceWidth)  
bool ParticleData::doForceWidth(int id)  
a flag applicable only for resonances (see isResonance above), whereby it is possible to force resonances to retain their assigned widths, whatever that is, see Resonance Decays for details. The normal behaviour is false, i.e. the width is based on hardcoded calculations whenever available.

void ParticleData::hasChanged(int id, bool hasChanged)  
bool ParticleData::hasChanged(int id)  
keep track of whether the data for a particle has been changed in any respect between initialization and the current status. Is used e.g. by the listChanged method to determine which particles to list.

bool ParticleData::useBreitWigner(int id)  
tells whether a particle will have a Breit-Wigner mass distribution or not. Is determined by an internal logic based on the particle width and on the value of the ParticleData:modeBreitWigner switch.

bool ParticleData::varWidth(int id)  
returns whether this particle has a defined mass-dependent width. These widths are parametrised in the HadronWidths class.

double ParticleData::constituentMass(int id)  
is the constituent mass for a quark, hardcoded as m_u = m_d = 0.325, m_s = 0.50, m_c = 1.60 and m_b = 5.0 GeV, for a diquark the sum of quark constituent masses, and for everything else the same as the ordinary mass.

double ParticleData::mSel(int id)  
returns a mass distributed according to a truncated Breit-Wigner, with parameters as described here. Is equal to m0(id) for particles without width.

double ParticleData::mRun(int id, double mH)  
calculate the running mass of species id when probed at a hard mass scale of mH. Only applied to obtain the running quark masses; for all other particle the normal fixed mass is used.

bool ParticleData::canDecay(int id)  
true for a particle with at least one decay channel defined.

bool ParticleData::isLepton(int id)  
true for a lepton or an antilepton (including neutrinos).

bool ParticleData::isQuark(int id)  
true for a quark or an antiquark.

bool ParticleData::isGluon(int id)  
true for a gluon.

bool ParticleData::isDiquark(int id)  
true for a diquark or antidiquark.

bool ParticleData::isParton()  
true for a gluon, a quark or antiquark up to the b (but excluding top), and a diquark or antidiquark consisting of quarks up to the b.

bool ParticleData::isHadron(int id)  
true for a hadron (made up out of normal quarks and gluons, i.e. not for R-hadrons and other exotic states).

bool ParticleData::isMeson(int id)  
true for a meson.

bool ParticleData::isBaryon(int id)  
true for a baryon or antibaryon.

bool ParticleData::isOnium(int id)  
true for a charmonium, bottomonium or (hypothetical) toponium state, i.e any state composed of c cbar, b bbar or t tbar irrespective of colour (singlet or octet) and spin, angular momentum or radial excitation.

bool ParticleData::isExotic(int id)  
true for an exotic hadron (the 9th PDG ID digit is 9 and the remaining digits are consistent with a hadron).

bool ParticleData::isOctetHadron(int id)  
true for an intermediate hadron-like state with a colour octet charge as used in the colour octet model for onia production.

int ParticleData::heaviestQuark(int id)  
extracts the heaviest quark or antiquark, i.e. one with largest id number, for a hadron.

int ParticleData::baryonNumberType(int id)  
is 1 for a quark, 2 for a diquark, 3 for a baryon, the same with a minus sign for antiparticles, and else zero.

void ParticleData::rescaleBR(int id, double newSumBR = 1.)  
rescales all partial branching ratios by a common factor, such that the sum afterward becomes newSumBR.

void setResonancePtr(int id, ResonanceWidthsPtr resonancePtr)  
set a shared pointer for a particle kind to a ResonanceWidths object. This is done, from inside ParticleData::initWidths, only for resonances, i.e. for particles such as Z^0, W^+-, top, Higgs, and new unstable states beyond the Standard Model. The presence of such an object will allow a more dynamic calculation of partial and total widths, as illustrated by the following methods.

void ParticleData::resInit(int id)  
initialize the treatment of a resonance.

double ParticleData::resWidth(int id, double mHat, int idInFlav = 0, bool openOnly = false, bool setBR = false)  
calculate the total with for a resonance of a given current mass, optionally including coupling to incoming flavour state (consider the gamma*/Z^0 combination), optionally excluding decay channels that have been closed by the user, and optionally storing the results in the normal decay table.

double ParticleData::resWidthOpen(int id, double mHat, int idInFlav = 0)  
special case of resWidth, where only open channels are included, but results are not stored in the normal decay table.

double ParticleData::resWidthStore(int id, double mHat, int idInFlav = 0)  
special case of resWidth, where only open channels are included, and results are stored in the normal decay table.

double ParticleData::resOpenFrac(int id1, int id2 = 0, int id3 = 0)  
calculate the fraction of the full branching ratio that is left open by the user choice of allowed decay channels. Can be applied to a final state with up to three resonances. Since the procedure is multiplicative, it would be easy to generalize also to more.

double ParticleData::resWidthRescaleFactor(int id)  
the factor used to rescale all partial widths in case the total width is being forced to a specific value by the user.

double ParticleData::resWidthChan(int id, double mHat, int idAbs1 = 0, int idAbs2 = 0)  
special case to calculate one final-state width; currently only used for Higgs decay to q qbar, g g or gamma gamma.

ParticleDataEntry* ParticleData::particleDataEntryPtr(int id)  
returns a pointer to the ParticleDataEntry object. The methods in the next section can then be used to manipulate this object.

bool ParticleData::getIsInit()  
return true if the database has been initialized, else false.

The ParticleDataEntry methods

Most of the methods that can be applied to a single ParticleDataEntry object are almost identical with those used above for the ParticleData, except that the id argument is no longer needed to find the right entry in the table. By and large, this makes direct access to the ParticleDataEntry methods superfluous. There are a few methods that are unique to each class, however. Furthermore, to avoid some naming ambiguities, many methods that set values begin with set.

ParticleDataEntry::ParticleDataEntry(int id = 0, string name = " ", int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0., double tau0 = 0.)  
ParticleDataEntry::ParticleDataEntry(int id, string name, string antiName, int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0., double tau0 = 0.)  
there are two alternative constructors, that both expect the properties of a particle as input. The first assumes that there is only one particle, the latter that there is a particle-antiparticle pair (but if the antiparticle name is void one reverts back to the particle-only case).

ParticleDataEntry& operator=( const ParticleDataEntry& particleDataEntryIn)  
copy the values stored in an existing ParticleDataEntry object.

void ParticleDataEntry::setDefaults()  
initialize some particle flags with default values, e.g. whether a particle is a resonance, may decay, or is visible. Is called from the constructors and from setAll.

void ParticleDataEntry::initPtr(ParticleData* particleDataPtrIn)  
initialize pointer back to the whole database (so that masses of decay products can be accessed, e.g.).

void ParticleDataEntry::setAll( string name, string antiName, int spinType = 0, int chargeType = 0, int colType = 0, double m0 = 0., double mWidth = 0., double mMin = 0., double mMax = 0.,double tau0 = 0.)  
change all the properties of the particle associated with a given identity code.

int ParticleDataEntry::id()  
the PDG identity code.

bool ParticleDataEntry::hasAnti()  
tell whether a separate antiparticle exists.

void ParticleDataEntry::setName(string name)  
void ParticleDataEntry::setAntiName(string antiName)  
void ParticleDataEntry::setNames(string name, string antiName)  
string ParticleDataEntry::name(int id = 1)  
set or get the particle or antiparticle name. Only the sign of id is needed to distinguish particle/antiparticle.

void ParticleDataEntry::setSpinType(int spinType)  
int ParticleDataEntry::spinType()  
set or get the particle spin type, i.e. 2 s + 1, or 0 in some special cases.

void ParticleDataEntry::setChargeType(int chargeType)  
int ParticleDataEntry::chargeType(int id = 1)  
double ParticleDataEntry::charge(int id = 1)  
set or get the particle charge type, i.e. three times the charge, or the charge itself. Only the sign of id is needed to distinguish particle/antiparticle.

void ParticleDataEntry::setColType(int colType)  
int ParticleDataEntry::colType(int id = 1)  
set or get the particle colour type, 0 for singlet, 1 for triplet, -1 for antitriplet, 2 for octet. Only the sign of id is needed to distinguish particle/antiparticle.

void ParticleDataEntry::setM0(double m0)  
double ParticleDataEntry::m0()  
the nominal mass m_0 (in GeV).

void ParticleDataEntry::setMWidth(double mWidth)  
double ParticleDataEntry::mWidth()  
the width Gamma of the Breit-Wigner distribution (in GeV).

void ParticleDataEntry::setMMin(double mMin)  
double ParticleDataEntry::mMin()  
the lower limit of the allowed mass range generated by the Breit-Wigner (in GeV). Has no meaning for particles without width, and would typically be 0 there.

void ParticleDataEntry::setMMax(double mMax)  
double ParticleDataEntry::mMax()  
the upper limit of the allowed mass range generated by the Breit-Wigner (in GeV). If mMax < mMin then no upper limit is imposed. Has no meaning for particles without width, and would typically be 0 there.

double ParticleDataEntry::m0Min()  
similar to mMin() above, except that for particles with no width the m0(id) value is returned.

double ParticleDataEntry::m0Max()  
similar to mMax() above, except that for particles with no width the m0(id) value is returned.

void ParticleDataEntry::setTau0(double tau0, bool countAsChanged = true)  
double ParticleDataEntry::tau0()  
the nominal proper lifetime tau_0 (in mm/c). Second argument of input method for internal use only.

void ParticleDataEntry::setTauCalc(bool tauCalc)  
bool ParticleDataEntry::tauCalc()  
a flag telling whether the nominal proper lifetime tau_0 should be calculated from the particle width, mWidth, at initialization (the default behavior). If set false, then the proper lifetime can be varied independently from the particle width. This capability is useful for simulating models of long-lived particles where the dependence of a signature on the particle width is negligible.

void ParticleDataEntry::setIsResonance(bool isResonance)  
bool ParticleDataEntry::isResonance()  
a flag telling whether a particle species are considered as a resonance or not. Here "resonance" is used as shorthand for any massive particle where the decay process should be counted as part of the hard process itself, and thus be performed before showers and other event aspects are added. Restrictions on allowed decay channels is also directly reflected in the cross section of simulated processes, while those of normal hadrons and other light particles are not. In practice, it is reserved for states above the b bbar bound systems in mass, i.e. for W, Z, t, Higgs states, supersymmetric states and (most?) other states in any new theory. All particles with m0 above 20 GeV are by default initialized to be considered as resonances.

void ParticleDataEntry::setMayDecay(bool mayDecay)  
bool ParticleDataEntry::mayDecay()  
a flag telling whether a particle species may decay or not, offering the main user switch. Whether a given particle of this kind then actually will decay also depends on it having allowed decay channels, and on other flags for particle decays (or resonance decays). All particles with tau0 below 1000 mm are by default initialized to allow decays.

void ParticleDataEntry::setDoExternalDecays(bool doExternalDecays)  
bool ParticleDataEntry::doExternalDecay()  
a flag telling whether a particle should be handled by an external decay package or not, with the latter default. Can be manipulated as described on this page, but should normally not be. Instead the pythia.decayPtr method should be provided with the list of relevant particles.

void ParticleDataEntry::setIsVisible(bool isVisible)  
bool ParticleDataEntry::isVisible()  
a flag telling whether a particle species is to be considered as visible in a detector or not, as used e.g. in analysis routines. By default the invisibles include neutrinos, Dark Matter particles (codes 51 - 60) and a few BSM particles (gravitino, sneutrinos, neutralinos) that have neither strong nor electromagnetic charge, and are not made up of constituents that have it. The value of this flag is only relevant if a particle is long-lived enough actually to make it to a detector.

void ParticleDataEntry::setDoForceWidth(bool doForceWidth)  
bool ParticleDataEntry::doForceWidth()  
a flag applicable only for resonances (see isResonance above), whereby it is possible to force resonances to retain their assigned widths, whatever that is, see Resonance Decays for details. The normal behaviour is false, i.e. the width is based on hardcoded calculations whenever available.

void ParticleDataEntry::setVarWidth(bool varWidth)  
bool ParticleDataEntry::varWidth()  
sets or returns whether this particle has a defined mass-dependent width. These widths are parametrised in the HadronWidths class.

void ParticleDataEntry::setHasChanged(bool hasChanged)  

void ParticleDataEntry::hasChanged(bool hasChanged)  
keep track of whether the data for a particle has been changed in any respect between initialization and the current status. Is used e.g. by the ParticleData::listChanged method to determine which particles to list.

void ParticleDataEntry::initBWmass()  
Prepare the Breit-Wigner mass selection by precalculating frequently-used expressions.

double ParticleDataEntry::constituentMass()  
is the constituent mass for a quark, hardcoded as m_u = m_d = 0.325, m_s = 0.50, m_c = 1.60 and m_b = 5.0 GeV, for a diquark the sum of quark constituent masses, and for everything else the same as the ordinary mass.

double ParticleDataEntry::mSel()  
give the mass of a particle, either at the nominal value or picked according to a (linear or quadratic) Breit-Wigner.

double ParticleDataEntry::mRun(double mH)  
calculate the running quark mass at a hard scale mH. For other particles the on-shell mass is given.

bool ParticleDataEntry::useBreitWigner()  
tells whether a particle will have a Breit-Wigner mass distribution or not. Is determined by an internal logic based on the particle width and on the value of the ParticleData:modeBreitWigner switch.

bool ParticleDataEntry::canDecay(int id)  
true for a particle with at least one decay channel defined.

bool ParticleDataEntry::isLepton()  
true for a lepton or an antilepton (including neutrinos).

bool ParticleDataEntry::isQuark()  
true for a quark or an antiquark.

bool ParticleDataEntry::isGluon()  
true for a gluon.

bool ParticleDataEntry::isDiquark()  
true for a diquark or antidiquark.

bool ParticleDataEntry::isParton()  
true for a gluon, a quark or antiquark up to the b (but excluding top), and a diquark or antidiquark consisting of quarks up to the b.

bool ParticleDataEntry::isHadron()  
true for a hadron (made up out of normal quarks and gluons, i.e. not for R-hadrons and other exotic states).

bool ParticleDataEntry::isMeson()  
true for a meson.

bool ParticleDataEntry::isBaryon()  
true for a baryon or antibaryon.

bool ParticleDataEntry::isOnium()  
true for a charmonium, bottomonium or (hypothetical) toponium state, i.e any state composed of c cbar, b bbar or t tbar irrespective of colour (singlet or octet) and spin, angular momentum or radial excitation.

bool ParticleDataEntry::isExotic()  
true for an exotic hadron (the 9th PDG ID digit is 9 and the remaining digits are consistent with a hadron).

bool ParticleDataEntry::isOctetHadron()  
true for an intermediate hadron-like state with a colour octet charge as used in the colour octet model for onia production.

int ParticleDataEntry::heaviestQuark(int id)  
extracts the heaviest quark or antiquark, i.e. one with largest id number, for a hadron. Only the sign of the input argument is relevant.

int ParticleDataEntry::baryonNumberType(int id)  
is 1 for a quark, 2 for a diquark, 3 for a baryon, the same with a minus sign for antiparticles, and else zero. Only the sign of the input argument is relevant.

void ParticleDataEntry::clearChannels()  
resets to an empty decay table.

void ParticleDataEntry::addChannel(int onMode = 0, double bRatio = 0., int meMode = 0, int prod0 = 0, int prod1 = 0, int prod2 = 0, int prod3 = 0, int prod4 = 0, int prod5 = 0, int prod6 = 0, int prod7 = 0,)  
adds a decay channel with up to 8 products.

int ParticleDataEntry::sizeChannels()  
returns the number of decay channels for a particle.

DecayChannel& ParticleDataEntry::channel(int i)  
const DecayChannel& ParticleDataEntry::channel(int i)  
gain access to a specified channel in the decay table.

void ParticleDataEntry::rescaleBR(double newSumBR = 1.)  
rescales all partial branching ratios by a common factor, such that the sum afterward becomes newSumBR.

bool ParticleDataEntry::preparePick(int idSgn, double mHat = 0., int idInFlav = 0)  
prepare to pick a decay channel.

DecayChannel& ParticleDataEntry::pickChannel()  
pick a decay channel according to branching ratios from preparePick.

void ParticleDataEntry::setResonancePtr(ResonanceWidthsPtr resonancePtr)  
ResonanceWidthsPtr ParticleDataEntry::getResonancePtr()  
set or get a shared pointer to an object that can be used for dynamic calculation of partial and total resonance widths. Here a resonance is a particle such as top, Z^0, W^+-, Higgs, and new unstable states beyond the Standard Model.

void ParticleDataEntry::resInit(Info* infoPtrIn, Settings* settingsPtrIn, ParticleData* particleDataPtrIn, CoupSM* coupSMPtrIn)  
initialize the treatment of a resonance.

double ParticleDataEntry::resWidth(int idSgn, double mHat, int idInFlav = 0, bool openOnly = false, bool setBR = false)  
calculate the total with for a resonance of a given current mass, optionally including coupling to incoming flavour state (consider the gamma*/Z^0 combination), optionally excluding decay channels that have been closed by the user, and optionally storing the results in the normal decay table. For the first argument only the sign is relevant.

double ParticleDataEntry::resWidthOpen(int idSgn, double mHat, int idInFlav = 0)  
special case of resWidth, where only open channels are included, but results are not stored in the normal decay table.

double ParticleDataEntry::resWidthStore(int idSgn, double mHat, int idInFlav = 0)  
special case of resWidth, where only open channels are included, and results are stored in the normal decay table.

double ParticleDataEntry::resOpenFrac(int idSgn)  
calculate the fraction of the full branching ratio that is left open by the user choice of allowed decay channels.

double ParticleDataEntry::resWidthRescaleFactor()  
the factor used to rescale all partial widths in case the total width is being forced to a specific value by the user.

double ParticleDataEntry::resWidthChan(double mHat, int idAbs1 = 0, int idAbs2 = 0)  
special case to calculate one final-state width; currently only used for Higgs decay to q qbar, g g or gamma gamma.

The DecayChannel methods

The properties stored in an individual decay channel can be set or get by the methods in this section.

DecayChannel::DecayChannel(int onMode = 0, double bRatio = 0., int meMode = 0, int prod0 = 0, int prod1 = 0, int prod2 = 0, int prod3 = 0, int prod4 = 0, int prod5 = 0, int prod6 = 0, int prod7 = 0)  
the constructor for a decay channel. Internal.

DecayChannel& operator=( const DecayChannel& decayChannelIn)  
copy the values stored in an existing DecayChannel object.

void DecayChannel::onMode(int onMode)  
int DecayChannel::onMode()  
set or get the onMode of a decay channel,
0 if a channel is off,
1 if on,
2 if on for a particle but off for an antiparticle,
3 if on for an antiparticle but off for a particle.
If a particle is its own antiparticle then 2 is on and 3 off but, of course, for such particles it is much simpler and safer to use only 1 and 0.
The 2 and 3 options can be used e.g. to encode CP violation in B decays, or to let the W's in a q qbar → W^+ W^- process decay in different channels.

void DecayChannel::bRatio(double bRatio, bool countAsChanged = true)  
double DecayChannel::bRatio()  
set or get the branching ratio of the channel. Second argument only for internal use.

void DecayChannel::rescaleBR(double fac)  
multiply the current branching ratio by fac.

void DecayChannel::meMode(int meMode)  
int DecayChannel::meMode()  
set or get the mode of processing this channel, possibly with matrix elements (see the particle decays and resonance decays descriptions).

void DecayChannel::multiplicity(int multiplicity)  
int DecayChannel::multiplicity()  
set or get the number of decay products in a channel, at most 8. (Is normally not to be set by hand, since it is automatically updated whenever the products list is changed.)

void DecayChannel::product(int i, int product)  
int DecayChannel::product(int i)  
set or get a list of the decay products, 8 products 0 <= i < 8, with trailing unused ones set to 0.

void DecayChannel::setHasChanged(bool hasChanged)  
bool DecayChannel::hasChanged()  
used for internal purposes, to know which decay modes have been changed.

bool DecayChannel::contains(int id1)  
bool DecayChannel::contains(int id1, int id2)  
bool DecayChannel::contains(int id1, int id2, int id3)  
find if the decay product list contains the one, two or three particle identities provided. If the same code is repeated then so must it be in the products list. Matching also requires correct sign.

void DecayChannel::currentBR(double currentBR)  
double DecayChannel::currentBR()  
set or get the current branching ratio, taking into account on/off switches and dynamic width for resonances. For internal use.

void DecayChannel::onShellWidth(double onShellWidth)  
double DecayChannel::onShellWidth()  
set or get the current partial width of the channel; intended for resonances where the widths are recalculated based on the current resonance mass. For internal use.

void DecayChannel::onShellWidthFactor(double factor)  
multiply the current partial width by factor.

void DecayChannel::openSec(int idSgn, double openSecIn)  
double DecayChannel::openSec(nt idSgn)  
set or get the fraction of secondary open widths, separately for positive and negative particles. For internal use.