3 <title>The Particle Data Scheme</title>
4 <link rel="stylesheet" type="text/css" href="pythia.css"/>
5 <link rel="shortcut icon" href="pythia32.gif"/>
9 <h2>The Particle Data Scheme</h2>
11 The particle data scheme may take somewhat longer to understand than
12 the settings one. In particular the set of methods to access information
13 is rather more varied, to allow better functionality for advanced usage.
14 However, PYTHIA does come with a sensible default set of particle
15 properties and decay tables. Thus there is no need to learn any of the
16 methods on this page to get going. Only when you perceive a specific need
17 does it make sense to learn the basics.
20 The central section on this page is the Operation one. The preceding
21 sections are there mainly to introduce the basic structure and the set
22 of properties that can be accessed. The subsequent sections provide a
23 complete listing of the existing public methods, which most users
24 probably will have little interaction with.
28 The management of particle data is based on three classes:
30 <li><code>ParticleData</code>, which is the top-level class, with
31 methods that can be used to interrogate all particle data. It contains
32 a map of PDG particle identity numbers [<a href="Bibliography.html" target="page">Yao06</a>] onto the relevant
33 <code>ParticleDataEntry</code> objects,</li>
34 <li><code>ParticleDataEntry</code>, which stores the relevant information
35 on an individual particle species, and</li>
36 <li><code>DecayChannel</code>, which stores info on one particular decay
37 mode of a particle.</li>
40 The objects of these classes together form a database that is
41 continuously being used as the program has to assign particle masses,
42 select decay modes, etc.
45 Each <code>Pythia</code> object has a public member
46 <code>particleData</code> of the <code>ParticleData</code> class.
47 Therefore you access the particle data methods as
48 <code>pythia.particleData.command(argument)</code>,
49 assuming that <code>pythia</code> is an instance of the
50 <code>Pythia</code> class. Further, for some of the most frequent user
51 tasks, <code>Pythia</code> methods have been defined, so that
52 <code>pythia.command(argument)</code>
53 would work, see further below.
56 A fundamental difference between the particle data classes and the
57 settings ones is that the former are accessed regularly during the
58 event generation process, as a new particle is produced and its mass
59 need to be set, e.g., while the latter are mainly/only used
60 at the initialization stage. Nevertheless, it is not a good idea to
61 change data in either of them in mid-run, since this may lead to
64 <h3>Stored properties for particles</h3>
66 The main properties stored for each particle are as follows.
67 Different ways to set and get these properties will be described
72 <li><code>name</code>: a character string with the name of the
73 particle. Particle and antiparticle names are stored separately,
74 with <code>void</code> returned when no antiparticle exists.</li>
76 <li><code>spinType</code>: the spin type, of the form <i>2 s + 1</i>,
77 with special code 0 for entries of unknown or indeterminate spin.</li>
79 <li><code>chargeType</code>: three times the charge (to make it an
82 <li><code>colType</code>: the colour type, with 0 uncoloured, 1 triplet,
83 -1 antitriplet and 2 octet.</li>
85 <li><code>m0</code>: the nominal mass <i>m_0</i> (in GeV).</li>
87 <li><code>mWidth</code>: the width <i>Gamma</i> of the Breit-Wigner
88 distribution (in GeV).</li>
90 <li><code>mMin</code>: the lower limit of the allowed mass range
91 generated by the Breit-Wigner (in GeV). Has no meaning for particles
92 without width, and would typically be 0 there.</li>
94 <li><code>mMax</code>: the upper limit of the allowed mass range
95 generated by the Breit-Wigner (in GeV). If <i>mMax < mMin</i> then
96 no upper limit is imposed. Has no meaning for particles without width,
97 and would typically be 0 there.</li>
99 <li><code>tau0</code>: the nominal proper lifetime <i>tau_0</i>
102 <li><code>isResonance</code>: a flag telling whether a particle species
103 is considered as a resonance or not. Here
104 <a href="ResonanceDecays.html" target="page">"resonance"</a> is used as shorthand
105 for any massive particle where the decay process should be counted as part
106 of the hard process itself, and thus be performed before showers and other
107 event aspects are added. Restrictions on allowed decay channels is also
108 directly reflected in the cross section of simulated processes, while
109 those of normal hadrons and other light particles are not.
110 In practice, it is reserved for states above the <i>b bbar</i>
111 bound systems in mass, i.e. for <i>W, Z, t</i>, Higgs states,
112 supersymmetric states and (most?) other states in any new theory.
113 All particles with <code>m0</code> above 20 GeV are by default
114 initialized to be considered as resonances.</li>
116 <li><code>mayDecay</code>: a flag telling whether a particle species
117 may decay or not, offering the main user switch. Whether a given particle
118 of this kind then actually will decay also depends on it having allowed
119 decay channels, and on other flags for
120 <a href="ParticleDecays.html" target="page">particle decays</a>.
121 All particles with <code>tau0</code> below 1000 mm are
122 by default initialized to allow decays.</li>
124 <li><code>doExternalDecays</code>: a flag telling whether a particle
125 should be handled by an external decay package or not, with the latter
126 default. Can be manipulated as described on this page, but should
127 normally not be. Instead the
128 <code><a href="ExternalDecays.html" target="page">Pythia::decayPtr(...)</a></code>
129 method should be provided with the list of relevant particles.</li>
131 <li><code>isVisible</code>: a flag telling whether a particle species
132 is to be considered as visible in a detector or not, as used e.g. in
133 analysis routines. By default this includes neutrinos and a few BSM
134 particles (gravitino, sneutrinos, neutralinos) that have neither strong
135 nor electromagnetic charge, and are not made up of constituents that
136 have it. The value of this flag is only relevant if a particle is
137 long-lived enough actually to make it to a detector.</li>
139 <li><code>doForceWidth</code>: a flag valid only for resonances where
140 PYTHIA contains code to calculate the width of the resonance from
141 encoded matrix-element expressions, i.e. the <i>Z^0</i>, <i>W^+-</i>,
142 <i>t</i>, <i>h^0</i>, and a few more. The normal behaviour
143 (<code>false</code>) is then that the width is calculated from the mass,
144 but it is possible to <a href="ResonanceDecays.html" target="page">force</a> the
145 resonance to retain the nominal width. Branching ratios and the running
146 of the total width are unaffected.</li>
150 <h3>Stored properties for decays</h3>
152 An unstable particle has a decay table consisting of one or more
153 decay channel. The following properties are stored for each such channel.
154 Again different ways to set and get these properties will be described
158 <li><code>onMode</code>: integer code for use or not of channel,<br/>
159 0 if a channel is off,<br/>
161 2 if on for a particle but off for an antiparticle,<br/>
162 3 if on for an antiparticle but off for a particle.<br/>
163 If a particle is its own antiparticle then 2 is on and 3 off
164 but, of course, for such particles it is much simpler and safer
165 to use only 1 and 0.<br/>
166 The 2 and 3 options can be used e.g. to encode CP violation in
167 B decays, or to let the <i>W</i>'s in a <i>q qbar -> W^+ W^-</i>
168 process decay in different channels. </li>
170 <li><code>bRatio</code>: the branching ratio of the channel.</li>
172 <li><code>meMode</code>: the mode of processing this channel, possibly
173 with matrix elements; see the
174 <a href="ParticleDecays.html" target="page">particle decays</a> description</li>
175 for the list of possibilities.
177 <li><code>multiplicity</code>: the number of decay products of the
178 channel. Can be at most 8.</li>
180 <li><code>product(i)</code>: the identity code of the decay products,
181 where <code>i</code> runs between <code>0</code> and
182 <code>multiplicity - 1</code>. Trailing positions are filled with 0.
189 The normal flow of the particle data operations is:
194 When a <code>Pythia</code> object <code>pythia</code> is created, the
195 <code>pythia.particleData</code> member is asked to scan the
196 <code>ParticleData.xml</code> file.
199 All lines beginning with <code><particle</code> are scanned for
200 information on a particle species, and all lines beginning with
201 <code><channel</code> are assumed to contain a decay channel of the
202 enclosing particle. In both cases XML syntax is used, with attributes
203 used to identify the stored properties, and with omitted properties
204 defaulting back to 0 where meaningful. The particle and channel
205 information may be split over several lines, up to the > endtoken.
206 The format of a <code><particle</code> tag is:
208 <particle id="..." name="..." antiName="..." spinType="..." chargeType="..." colType="..."
209 m0="..." mWidth="..." mMin="..." mMax="..." tau0="...">
212 where the fields are the properties already introduced above.
213 Note that <code>isResonance</code>, <code>mayDecay</code>,
214 <code>doExternalDecay</code>, <code>isVisible</code> and
215 <code>doForceWidth</code> are not set here, but are provided with
216 default values by the rules described above. Once initialized, also
217 these latter properties can be changed, see below.<br/>
219 The format of a <code><channel></code> tag is:
221 <channel onMode="..." bRatio="..." meMode="..." products="..." />
223 again see properties above. The products are given as a blank-separated
224 list of <code>id</code> codes.
225 <br/><b>Important</b>: the values in the <code>.xml</code> file should not
226 be changed, except by the PYTHIA authors. Any changes should be done
227 with the help of the methods described below.
231 Between the creation of the <code>Pythia</code> object and the
232 <code>init</code> call for it, you may use the methods of the
233 <code>ParticleData</code> class to modify some of the default values.
234 Several different approaches can be chosen for this.
237 a) Inside your main program you can directly set values with
239 pythia.readString(string);
241 where both the variable name and the value are contained inside
242 the character string, separated by blanks and/or a =, e.g.
244 pythia.readString("111:mayDecay = off");
246 switches off the decays of the <i>pi^0</i>.<br/>
248 The particle id (> 0) and the property to be changed must be given,
249 separated by a colon.<br/>
251 The allowed properties are: <code>name</code>, <code>antiName</code>,
252 <code>spinType</code>, <code>chargeType</code>, <code>colType</code>,
253 <code>m0</code>, <code>mWidth</code>, <code>mMin</code>,
254 <code>mMax</code>, <code>tau0</code>, <code>isResonance</code>,
255 <code>mayDecay</code>, <code>doExternalDecay</code>,
256 <code>isVisible</code> and <code>doForceWidth</code>. All of these
257 names are case-insensitive. Names that do not match an existing
258 variable are ignored.<br/>
260 Strings beginning with a non-alphanumeric character, like # or !,
261 are assumed to be comments and are not processed at all. For
262 <code>bool</code> values, the following notation may be used
263 interchangeably: <code>true = on = yes = ok = 1</code>, while everything
264 else gives <code>false</code> (including but not limited to
265 <code>false</code>, <code>off</code>, <code>no</code> and
269 Particle data often comes in sets of closely related information.
270 Therefore some properties expect the value to consist of several
271 numbers. These can then be separated by blanks (or by commas).
272 A simple example is <code>names</code>, which expects both the
273 name and antiname to be given. A more interesting one is the
274 <code>all</code> property,
276 id:all = name antiName spinType chargeType colType m0 mWidth mMin mMax tau0
278 where all the current information on the particle itself is replaced,
279 but any decay channels are kept unchanged. Using <code>new</code> instead
280 of <code>all</code> also removes any previous decay channels.
281 If the string contains fewer fields than expected the trailing
282 properties are set to vanish ("void", 0 or 0.). Note that such a
283 truncated string should not be followed by a comment, since this
284 comment would then be read in as if it contained the missing properties.
285 The truncation can be done anywhere, specifically a string with only
286 <code>id:new</code> defines a new "empty" particle.
287 As before, <code>isResonance</code>, <code>mayDecay</code>,
288 <code>doExternalDecay</code>, <code>isVisible</code> and
289 <code>doForceWidth</code> are (re)set to their default values, and
290 would have to be changed separately if required.
293 A further command is <code>rescaleBR</code>, which rescales each of the
294 existing branching ratios with a common factor, such that their new
295 sum is the provided value. This may be a first step towards adding
296 new decay channels, see further below.
299 Alternatively the <code>id</code> code may be followed by another integer,
300 which then gives the decay channel number. This then has to be
301 followed by the property specific to this channel, either
302 <code>onMode</code>, <code>bRatio</code>, <code>meMode</code> or
303 <code>products</code>. In the latter case all the products of
304 the channel should be given:
306 id:channel:products = product1 product2 ....
308 The line will be scanned until the end of the line, or until a
309 non-number word is encountered, or until the maximum allowed number
310 of eight products is encountered, whichever happens first. (Thus the
311 multiplicity of a decay channel need not be input; it is automatically
312 calculated from the products list.) It is also possible to replace all
313 the properties of a channel in a similar way:
315 id:channel:all = onMode bRatio meMode product1 product2 ....
317 To add a new channel at the end, use
319 id:addChannel = onMode bRatio meMode product1 product2 ....
323 It is currently not possible to remove a channel selectively, but
324 setting its branching ratio vanishing is as effective. If you want to
325 remove all existing channels and force decays into one new channel
328 id:oneChannel = onMode bRatio meMode product1 product2 ....
330 A first <code>oneChannel</code> command could be followed by
331 several subsequent <code>addChannel</code> ones, to build
332 up a completely new decay table for an existing particle.
335 When adding new channels or changing branching ratios in general,
336 note that, once a particle is to be decayed, the sum of branching
337 ratios is always rescaled to unity. Beforehand, <code>rescaleBR</code>
338 may be used to rescale an existing branching ratio by the given factor.
341 There are a few commands that will study all the decay channels of the
342 given particle, to switch them on or off as desired. The
346 will set the <code>onMode</code> property of all channels to the
349 id:offIfAny = product1 product2 ....
350 id:onIfAny = product1 product2 ....
351 id:onPosIfAny = product1 product2 ....
352 id:onNegIfAny = product1 product2 ....
354 will set the <code>onMode</code> 0, 1, 2 or 3, respectively, for all
355 channels which contain any of the enumerated products, where the matching
356 to these products is done without distinction of particles and
357 antiparticles. Note that "<code>Pos</code>" and "<code>Neg</code>"
358 are slightly misleading since it refers to the particle and antiparticle
359 of the <code>id</code> species rather than charge, but should still be
360 simpler to remember and understand than alternative notations.
363 id:offIfAll = product1 product2 ....
364 id:onIfAll = product1 product2 ....
365 id:onPosIfAll = product1 product2 ....
366 id:onNegIfAll = product1 product2 ....
368 will set the <code>onMode</code> 0, 1, 2 or 3, respectively, for all
369 channels which contain all of the enumerated products, again without
370 distinction of particles and antiparticles. If the same product appears
371 twice in the list it must also appear twice in the decay channel, and
372 so on. The decay channel is allowed to contain further particles,
373 beyond the product list. By contrast,
375 id:offIfMatch = product1 product2 ....
376 id:onIfMatch = product1 product2 ....
377 id:onPosIfMatch = product1 product2 ....
378 id:onPosIfMatch = product1 product2 ....
380 requires the decay-channel multiplicity to agree with that of the product
381 list, but otherwise works as the <code>onIfAll/offIfAll</code> methods.
384 Note that the action of several of the commands depends on the order
385 in which they are executed, as one would logically expect. For instance,
386 <code>id:oneChannel</code> removes all decay channels of <code>id</code>
387 and thus all previous changes in this decay table, while subsequent
388 additions or changes would still take effect. Another example would be that
389 <code>23:onMode = off</code> followed by <code>23:onIfAny = 1 2 3 4 5</code>
390 would let the <i>Z^0</i> decay to quarks, while no decays would be
391 allowed if the order were to be reversed.
394 b) The <code>Pythia</code> <code>readString(string)</code> method actually
395 does not do changes itself, but sends on the string either to the
396 <code>ParticleData</code> class or to the <code>Settings</code> one,
397 depending on whether the string begins with a digit or a letter.
398 If desired, it is possible to communicate directly with the corresponding
399 <code>ParticleData</code> method:
401 pythia.particleData.readString("111:mayDecay = off");
402 pythia.particleData.readString("15:2:products = 16 -211");
404 In this case, changes intended for <code>Settings</code> would not be
408 c) Underlying this are commands for all the individual properties in
409 the <code>ParticleData</code> class, one for each. These are
410 further described below. Thus, an example now reads
412 pythia.particleData.mayDecay(111, false);
414 Boolean values should here be given as <code>true</code> or
418 d) A simpler and more useful way is to collect all your changes
419 in a separate file, with one line per change, e.g.
423 The file can be read by the
425 pythia.readFile(fileName);
427 method, where <code>fileName</code> is a string, e.g.
428 <code>pythia.readFile("main.cmnd")</code> (or an <code>istream</code>
429 instead of a <code>fileName</code>). Each line is processed as
430 described for the string in 2a). This file can freely mix commands to
431 the <code>Settings</code> and <code>ParticleData</code> classes.
435 A routine <code>reInit(fileName)</code> is provided, and can be used to
436 zero the particle data table and reinitialize it from scratch.
437 Such a call might be useful if several subruns are to be made with
438 widely different particle data - normally the maps are only built
439 from scratch once, namely when the <code>Pythia()</code> object is
440 created. Also, there is no other possibility to restore the default
441 values, unlike for the settings.
445 You may at any time obtain a listing of all the particle data by calling
447 pythia.particleData.listAll();
449 The listing is by increasing <code>id</code> number. It shows the basic
450 quantities introduced above. Some are abbreviated in the header to fit on
451 the lines: <code>spn = spinType</code>, <code>chg = chargeType</code>,
452 <code>col = colType</code>, <code>res = isResonance</code>,
453 <code>dec = mayDecay && canDecay</code> (the latter checks that decay
454 channels have been defined), <code>ext = doExternalDecay</code>,
455 <code>vis = isVisible</code> and <code>wid = doForceWidth</code>.<br/>
457 To list only those particles that were changed (one way or another, the
458 listing will not tell what property or decay channel was changed), instead use
460 pythia.particleData.listChanged();
462 (This info is based on a further <code>hasChanged</code> flag of a particle
463 or a channel, set <code>true</code> whenever any of the changing methods are
464 used. It is possible to manipulate this value, but this is not recommended.)
465 By default the internal initialization of the widths of resonances such as
466 <i>gamma^*/Z^0, W^+-, t/tbar, H^0</i> do not count as changes; if you want
467 to list also those changes instead call <code>listChanged(true)</code>.
470 To list only one particle, give its <code>id</code> code as argument to
471 the <code>list(...)</code> function.. To list a restricted set of particles,
472 give in their <code>id</code> codes to <code>list(...)</code> as a
473 <code>vector<int></code>.
477 For wholesale changes of particle properties all available data can be
478 written out, edited, and then read back in again. These methods are
479 mainly intended for expert users. You can choose between two alternative
483 a) XML syntax, using the <code><particle</code> and
484 <code><channel</code> lines already described. You use the method
485 <code>particleData.listXML(fileName)</code> to produce such an XML
486 file and <code>particleData.readXML(fileName)</code> to read it back
490 b) Fixed/free format, using exactly the same information as illustrated
491 for the <code><particle</code> and <code><channel</code> lines
492 above, but now without any tags. This means that all information fields
493 must be provided (if there is no antiparticle then write
494 <code>void</code>), in the correct order (while the order is irrelevant
495 with XML syntax), and all on one line. Information is written out in
496 properly lined-up columns, but the reading is done using free format,
497 so fields need only be separated by at least one blank. Each new particle
498 is supposed to be separated by (at least) one blank line, whereas no
499 blank lines are allowed between the particle line and the subsequent
500 decay channel lines, if any. You use the method
501 <code>particleData.listFF(fileName)</code> to produce such a fixed/free
502 file and <code>particleData.readFF(fileName)</code> to read it back
506 As an alternative to the <code>readXML</code> and <code>readFF</code>
507 methods you can also use the
508 <code>particleData.reInit(fileName, xmlFormat)</code> method, where
509 <code>xmlFormat = true</code> (default) corresponds to reading an XML
510 file and <code>xmlFormat = false</code> to a fixed/free format one.
513 To check that the new particle and decay tables makes sense, you can use
514 the <code>particleData.checkTable()</code> method, either directly or by
515 switching it on among the standard
516 <a href="ErrorChecks.html" target="page">error checks</a>.
521 <h2>The public methods</h2>
523 In the following we present briefly the public methods in the three
524 classes used to build up the particle database. The order
525 is top-down, i.e from the full table of all particles to a single
526 particle to a single channel.
527 Note that these methods usually are less elegant and safe than the
528 input methods outlined above. If you use any of these methods, it is
529 likely to be the ones in the full database, i.e. the first ones to be
530 covered in the following.
533 For convenience, we have grouped related input and output methods
534 together. It should be obvious from the context which is which:
535 the input is of type <code>void</code> and has an extra last argument,
536 namely is the input value, while the output method returns a
537 quantity of the expected type.
539 <h3>The ParticleData methods</h3>
541 <a name="method1"></a>
542 <p/><strong>ParticleData::ParticleData() </strong> <br/>
543 the constructor has no arguments and does not do anything. Internal.
546 <a name="method2"></a>
547 <p/><strong>void ParticleData::initPtr(Info* infoPtr,Settings* settingsPtrIn, Rndm* rndmPtrIn, CoupSM* coupSMPtrIn) </strong> <br/>
548 initialize pointers to a few other classes. Internal.
551 <a name="method3"></a>
552 <p/><strong>bool ParticleData::init(string startFile = "../xmldoc/ParticleData.xml") </strong> <br/>
553 read in an XML-style file with particle data and initialize the
554 particle data tables accordingly. This command is executed
555 in the <code>Pythia</code> constructor, i.e. is mainly for
557 <br/><code>argument</code><strong> startFile </strong> (<code>default = <strong>../xmldoc/ParticleData.xml</strong></code>) :
558 the name of the data file to be read. When called from the
559 <code>Pythia</code> constructor the directory is provided by the
560 <code><a href="ProgramFlow.html" target="page">PYTHIA8DATA</a></code>
561 environment variable, if set, else by the argument of this constructor,
562 which has the default value "../xmldoc".
566 <a name="method4"></a>
567 <p/><strong>bool ParticleData::reInit(string startFile,bool xmlFormat = true) </strong> <br/>
568 overwrite the existing database by reading from the specified file.
569 Unlike <code>init</code> above this method is not called by the
570 <code>Pythia</code> constructor, but is entirely intended for users
571 who want to replace the existing particle data with their own.
572 <br/><code>argument</code><strong> startFile </strong> : the path and name of file to be read.
574 <br/><code>argument</code><strong> xmlFormat </strong> : if true read the same kind of XML-style file
575 as used by <code>init</code>, if not use an alternative "free format"
576 file (i.e. without any XML tags, but with well-defined rules
577 specifying in which order properties are stored).
581 <a name="method5"></a>
582 <p/><strong>void ParticleData::initWidths(vector<ResonanceWidths*> resonancePtrs) </strong> <br/>
583 initialize Breit-Wigner shape parameters for all particles,
584 and the detailed handling of resonances, i.e. particles with
585 perturbatively calculable partial widths, which can be used to
586 obtain a mass-dependent Breit-Wigner and a dynamic choice of
587 decay channels. Called from <code>Pythia::init()</code>.
590 <a name="method6"></a>
591 <p/><strong>bool ParticleData::readXML(string inFile, bool reset = true) </strong> <br/>
593 <strong>void ParticleData::listXML(string outFile) </strong> <br/>
594 read in XML-style data from a file or write it out to a file. For the
595 former one can also decide whether to reset all particles to scratch,
596 or only overwrite those particles in the file. The former method is
597 used by <code>init</code> and <code>reInit</code> above.
600 <a name="method7"></a>
601 <p/><strong>bool ParticleData::readFF(string inFile, bool reset = true) </strong> <br/>
603 <strong>void ParticleData::listFF(string outFile) </strong> <br/>
604 read in free-format-style data from a file or write it out to a file.
605 For the former one can also decide whether to reset all particles to
606 scratch, or only overwrite those particles in the file. The former
607 method is used by <code>reInit</code> above.
610 <a name="method8"></a>
611 <p/><strong>bool ParticleData::readString(string line, bool warn = true, ostream& os = cout) </strong> <br/>
612 read in a string and interpret is as a new or changed particle data.
613 The possibilities are extensively described above. It is normally
614 used indirectly, via <code>Pythia::readString(...)</code> and
615 <code>Pythia::readFile(...)</code>.
616 <br/><code>argument</code><strong> line </strong> :
617 the string to be interpreted as an instruction.
619 <br/><code>argument</code><strong> warn </strong> (<code>default = <strong>true</strong></code>) :
620 write a warning message or not whenever the instruction does not make
621 sense, e.g. if the particle does not exist in the database.
623 <br/><code>argument</code><strong> os </strong> (<code>default = <strong>cout</strong></code>) :
624 stream for error printout.
626 <br/><b>Note:</b> the method returns false if it fails to
627 make sense out of the input string.
630 <a name="method9"></a>
631 <p/><strong>void ParticleData::listAll(ostream& os = cout) </strong> <br/>
633 <strong>void ParticleData::listChanged(ostream& os = cout) </strong> <br/>
635 <strong>void ParticleData::listChangedAndRes(ostream& os = cout) </strong> <br/>
637 <strong>void ParticleData::list(bool changedOnly = false, bool changedRes = true, ostream& os = cout) </strong> <br/>
638 methods intended to present a listing of particle data in a readable
639 format. The first three are special cases of the fourth. The first
640 lists all particle data, the second only data for those particles that
641 were changed after the original creation of the particle data table.
642 Resonances are a special case since they can get their data changed
643 by being linked to an object that does the calculation of branching
644 ratios. The second method does not count such resonances as changed,
645 whereas the third does and thus lists all resonances.
648 <a name="method10"></a>
649 <p/><strong>void ParticleData::list(int idList, ostream& os = cout) </strong> <br/>
651 <strong>void ParticleData::list(vector<int> idList, ostream& os = cout) </strong> <br/>
652 list particle data for one single particle, with the identity code as
653 input, or for a set of particles, with an input vector of identity codes.
656 <a name="method11"></a>
657 <p/><strong>void ParticleData::checkTable(ostream& os = cout) </strong> <br/>
659 <strong>void ParticleData::checkTable(int verbosity,ostream& os = cout) </strong> <br/>
660 check that the particle decay table makes sense, especially for decays.
661 <br/><code>argument</code><strong> verbosity </strong> : level of checks. 0 is only mininal,
662 e.g. if a particle has no open decay channels. 1, which is the level
663 of the first method, provides warning if any individual channel is
664 closed, except for resonances. 2 also prints the
665 branching-ratio-averaged threshold mass. 11 and 12 are like 1 and 2,
666 but also include resonances in the detailed checks.
670 <a name="method12"></a>
671 <p/><strong>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.) </strong> <br/>
673 <strong>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.) </strong> <br/>
674 add a particle to the decay table; in the first form a partcle which is
675 its own antiparticle, in the second where a separate antiparticle exists.
678 <a name="method13"></a>
679 <p/><strong>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.) </strong> <br/>
680 change all the properties of the particle associated with a given
684 <a name="method14"></a>
685 <p/><strong>bool ParticleData::isParticle(int id) </strong> <br/>
686 query whether the particle data table contains the particle of the
690 <a name="method15"></a>
691 <p/><strong>int ParticleData::nextId(int id) </strong> <br/>
692 return the identity code of the sequentially next particle stored in table.
695 <a name="method16"></a>
696 <p/><strong>bool ParticleData::hasAnti(int id) </strong> <br/>
697 bool whether a distinct antiparticle exists or not. Is true if an
698 antiparticle name has been set (and is different from
702 <a name="method17"></a>
703 <p/><strong>void ParticleData::name(int id, string name) </strong> <br/>
705 <strong>void ParticleData::antiName(int id, string antiName) </strong> <br/>
707 <strong>void ParticleData::names(int id, string name, string antiName) </strong> <br/>
709 <<strong>string ParticleData::name(int id) </strong> <br/>
710 particle and antiparticle names are stored separately, the sign of
711 <code>id</code> determines which of the two is returned, with
712 <code>void</code> used to indicate the absence of an antiparticle.
715 <a name="method18"></a>
716 <p/><strong>void ParticleData::spinType(int id, int spinType) </strong> <br/>
718 <strong>int ParticleData::spinType(int id) </strong> <br/>
719 the spin type, of the form <i>2 s + 1</i>, with special code 0
720 for entries of unknown or indeterminate spin.
723 <a name="method19"></a>
724 <p/><strong>void ParticleData::chargeType(int id, int chargeType) </strong> <br/>
726 <strong>int ParticleData::chargeType(int id) </strong> <br/>
727 three times the charge (to make it an integer), taking into account
728 the sign of <code>id</code>.
731 <a name="method20"></a>
732 <p/><strong>double ParticleData::charge(int id) </strong> <br/>
733 the electrical charge of a particle, equal to
734 <code>chargeType(id)/3</code>.
737 <a name="method21"></a>
738 <p/><strong>void ParticleData::colType(int id, int colType) </strong> <br/>
740 <strong>int ParticleData::colType(int id) </strong> <br/>
741 the colour type, with 0 uncoloured, 1 triplet, -1 antitriplet and 2
742 octet, taking into account the sign of <code>id</code>.
745 <a name="method22"></a>
746 <p/><strong>void ParticleData::m0(int id, double m0) </strong> <br/>
748 <strong>double ParticleData::m0(int id) </strong> <br/>
749 the nominal mass <i>m_0</i> (in GeV).
752 <a name="method23"></a>
753 <p/><strong>void ParticleData::mWidth(int id, double mWidth) </strong> <br/>
755 <strong>double ParticleData::mWidth(int id) </strong> <br/>
756 the width <i>Gamma</i> of the Breit-Wigner distribution (in GeV).
759 <a name="method24"></a>
760 <p/><strong>void ParticleData::mMin(int id, double mMin) </strong> <br/>
762 <strong>double ParticleData::mMin(int id) </strong> <br/>
763 the lower limit of the allowed mass range generated by the Breit-Wigner
764 (in GeV). Has no meaning for particles without width, and would
765 typically be 0 there.
768 <a name="method25"></a>
769 <p/><strong>void ParticleData::mMax(int id, double mMax) </strong> <br/>
771 <strong>double ParticleData::mMax(int id) </strong> <br/>
772 the upper limit of the allowed mass range generated by the Breit-Wigner
773 (in GeV). If <i>mMax < mMin</i> then no upper limit is imposed.
774 Has no meaning for particles without width, and would typically
778 <a name="method26"></a>
779 <p/><strong>double ParticleData::m0Min(int id) </strong> <br/>
780 similar to <code>mMin()</code> above, except that for particles with
781 no width the <code>m0(id)</code> value is returned.
784 <a name="method27"></a>
785 <p/><strong>double ParticleData::m0Max(int id) </strong> <br/>
786 similar to <code>mMax()</code> above, except that for particles with
787 no width the <code>m0(id)</code> value is returned.
790 <a name="method28"></a>
791 <p/><strong>void ParticleData::tau0(int id, double tau0) </strong> <br/>
793 <strong>double ParticleData::tau0(int id) </strong> <br/>
794 the nominal proper lifetime <i>tau_0</i> (in mm/c).
797 <a name="method29"></a>
798 <p/><strong>void ParticleData::isResonance(int id, bool isResonance) </strong> <br/>
800 <strong>bool ParticleData::isResonance(int id) </strong> <br/>
801 a flag telling whether a particle species are considered as a resonance
802 or not. Here <a href="ResonanceDecays.html" target="page">"resonance"</a>
803 is used as shorthand for any massive particle
804 where the decay process should be counted as part of the hard process
805 itself, and thus be performed before showers and other event aspects
806 are added. Restrictions on allowed decay channels is also directly
807 reflected in the cross section of simulated processes, while those of
808 normal hadrons and other light particles are not.
809 In practice, it is reserved for states above the <i>b bbar</i>
810 bound systems in mass, i.e. for <i>W, Z, t</i>, Higgs states,
811 supersymmetric states and (most?) other states in any new theory.
812 All particles with <code>m0</code> above 20 GeV are by default
813 initialized to be considered as resonances.
816 <a name="method30"></a>
817 <p/><strong>void ParticleData::mayDecay(int id, bool mayDecay) </strong> <br/>
819 <strong>bool ParticleData::mayDecay(int id) </strong> <br/>
820 a flag telling whether a particle species may decay or not, offering
821 the main user switch. Whether a given particle of this kind then actually
822 will decay also depends on it having allowed decay channels, and on
823 other flags for <a href="ParticleDecays.html" target="page">particle decays</a>.
824 All particles with <code>tau0</code> below 1000 mm are
825 by default initialized to allow decays.
828 <a name="method31"></a>
829 <p/><strong>void ParticleData::doExternalDecays(int id, bool doExternalDecays) </strong> <br/>
831 <strong>bool ParticleData::doExternalDecay(int id) </strong> <br/>
832 a flag telling whether a particle should be handled by an external
833 decay package or not, with the latter default. Can be manipulated as
834 described on this page, but should normally not be. Instead the
835 <code><a href="ExternalDecays.html" target="page">pythia.decayPtr</a></code>
836 method should be provided with the list of relevant particles.
839 <a name="method32"></a>
840 <p/><strong>void ParticleData::isVisible(int id, bool isVisible) </strong> <br/>
842 <strong>bool ParticleData::isVisible(int id) </strong> <br/>
843 a flag telling whether a particle species is to be considered as
844 visible in a detector or not, as used e.g. in analysis routines.
845 By default this includes neutrinos and a few BSM particles
846 (gravitino, sneutrinos, neutralinos) that have neither strong nor
847 electromagnetic charge, and are not made up of constituents that
848 have it. The value of this flag is only relevant if a particle is
849 long-lived enough actually to make it to a detector.
852 <a name="method33"></a>
853 <p/><strong>void ParticleData::doForceWidth(int id, bool doForceWidth) </strong> <br/>
855 <strong>bool ParticleData::doForceWidth(int id) </strong> <br/>
856 a flag valid only for resonances where PYTHIA contains code to
857 calculate the width of the resonance from encoded matrix-element
858 expressions, i.e. the <i>Z^0</i>, <i>W^+-</i>, <i>t</i>,
859 <i>h^0</i>, and a few more. The normal behaviour (<code>false</code>)
860 is then that the width is calculated from the mass, but it is
861 possible to <a href="ResonanceDecays.html" target="page">force</a> the resonance
862 to retain the nominal width. Branching ratios and the running of the
863 total width are unaffected.
866 <a name="method34"></a>
867 <p/><strong>void ParticleData::hasChanged(int id, bool hasChanged) </strong> <br/>
869 <strong>bool ParticleData::hasChanged(int id) </strong> <br/>
870 keep track of whether the data for a particle has been changed
871 in any respect between initialization and the current status.
872 Is used e.g. by the <code>listChanged</code> method to determine
873 which particles to list.
876 <a name="method35"></a>
877 <p/><strong>bool ParticleData::useBreitWigner(int id) </strong> <br/>
878 tells whether a particle will have a Breit-Wigner mass distribution or
879 not. Is determined by an internal logic based on the particle width and
881 <code><a href="ParticleData.html" target="page">ParticleData:modeBreitWigner</a></code>
885 <a name="method36"></a>
886 <p/><strong>double ParticleData::constituentMass(int id) </strong> <br/>
887 is the constituent mass for a quark, hardcoded as
888 <i>m_u = m_d = 0.325</i>, <i>m_s = 0.50</i>, <i>m_c = 1.60</i>
889 and <i>m_b = 5.0</i> GeV, for a diquark the sum of quark constituent
890 masses, and for everything else the same as the ordinary mass.
893 <a name="method37"></a>
894 <p/><strong>double ParticleData::mass(int id) </strong> <br/>
895 returns a mass distributed according to a truncated Breit-Wigner,
896 with parameters as described here. Is equal to <code>m0(id)</code> for
897 particles without width.
900 <a name="method38"></a>
901 <p/><strong>double ParticleData::mRun(int id, double mH) </strong> <br/>
902 calculate the running mass of species <code>id</code> when probed at a
903 hard mass scale of <code>mH</code>. Only applied to obtain the
904 running quark masses; for all other particle the normal fixed mass
908 <a name="method39"></a>
909 <p/><strong>bool ParticleData::canDecay(int id) </strong> <br/>
910 true for a particle with at least one decay channel defined.
913 <a name="method40"></a>
914 <p/><strong>bool ParticleData::isLepton(int id) </strong> <br/>
915 true for a lepton or an antilepton (including neutrinos).
918 <a name="method41"></a>
919 <p/><strong>bool ParticleData::isQuark(int id) </strong> <br/>
920 true for a quark or an antiquark.
923 <a name="method42"></a>
924 <p/><strong>bool ParticleData::isGluon(int id) </strong> <br/>
928 <a name="method43"></a>
929 <p/><strong>bool ParticleData::isDiquark(int id) </strong> <br/>
930 true for a diquark or antidiquark.
933 <a name="method44"></a>
934 <p/><strong>bool ParticleData::isHadron(int id) </strong> <br/>
935 true for a hadron (made up out of normal quarks and gluons,
936 i.e. not for R-hadrons and other exotic states).
939 <a name="method45"></a>
940 <p/><strong>bool ParticleData::isMeson(int id) </strong> <br/>
944 <a name="method46"></a>
945 <p/><strong>bool ParticleData::isBaryon(int id) </strong> <br/>
946 true for a baryon or antibaryon.
949 <a name="method47"></a>
950 <p/><strong>bool ParticleData::isOctetHadron(int id) </strong> <br/>
951 true for an intermediate hadron-like state with a colour octet charge
952 as used in the colour octet model for
953 <a href="OniaProcesses.html" target="page">onia</a> production.
956 <a name="method48"></a>
957 <p/><strong>int ParticleData::heaviestQuark(int id) </strong> <br/>
958 extracts the heaviest quark or antiquark, i.e. one with largest
959 <code>id</code> number, for a hadron.
962 <a name="method49"></a>
963 <p/><strong>int ParticleData::baryonNumberType(int id) </strong> <br/>
964 is 1 for a quark, 2 for a diquark, 3 for a baryon, the same with a
965 minus sign for antiparticles, and else zero.
968 <a name="method50"></a>
969 <p/><strong>void ParticleData::rescaleBR(int id, double newSumBR = 1.) </strong> <br/>
970 rescales all partial branching ratios by a common factor, such that
971 the sum afterward becomes <code>newSumBR</code>.
974 <a name="method51"></a>
975 <p/><strong>void setResonancePtr(int id, ResonanceWidths* resonancePtr) </strong> <br/>
976 set a pointer for a particle kind to a <code>ResonanceWidths</code> object.
977 This is done, from inside <code>ParticleData::initWidths</code>, only for
978 resonances, i.e. for particles such as <i>Z^0</i>, <i>W^+-</i>, top,
979 Higgs, and new unstable states beyond the Standard Model. The presence
980 of such an object will allow a more dynamic calculation of partial and
981 total widths, as illustrated by the following methods.
984 <a name="method52"></a>
985 <p/><strong>void ParticleData::resInit(int id) </strong> <br/>
986 initialize the treatment of a resonance.
989 <a name="method53"></a>
990 <p/><strong>double ParticleData::resWidth(int id, double mHat, int idInFlav = 0, bool openOnly = false, bool setBR = false) </strong> <br/>
991 calculate the total with for a resonance of a given current mass,
992 optionally including coupling to incoming flavour state (consider
993 the <i>gamma*/Z^0</i> combination), optionally excluding decay
994 channels that have been closed by the user, and optionally storing
995 the results in the normal decay table.
998 <a name="method54"></a>
999 <p/><strong>double ParticleData::resWidthOpen(int id, double mHat, int idInFlav = 0) </strong> <br/>
1000 special case of <code>resWidth</code>, where only open channels are
1001 included, but results are not stored in the normal decay table.
1004 <a name="method55"></a>
1005 <p/><strong>double ParticleData::resWidthStore(int id, double mHat, int idInFlav = 0) </strong> <br/>
1006 special case of <code>resWidth</code>, where only open channels are
1007 included, and results are stored in the normal decay table.
1010 <a name="method56"></a>
1011 <p/><strong>double ParticleData::resOpenFrac(int id1, int id2 = 0, int id3 = 0) </strong> <br/>
1012 calculate the fraction of the full branching ratio that is left
1013 open by the user choice of allowed decay channels. Can be applied
1014 to a final state with up to three resonances. Since the procedure
1015 is multiplicative, it would be easy to generalize also to more.
1018 <a name="method57"></a>
1019 <p/><strong>double ParticleData::resWidthRescaleFactor(int id) </strong> <br/>
1020 the factor used to rescale all partial widths in case the total
1021 width is being forced to a specific value by the user.
1024 <a name="method58"></a>
1025 <p/><strong>double ParticleData::resWidthChan(int id,double mHat, int idAbs1 = 0, int idAbs2 = 0) </strong> <br/>
1026 special case to calculate one final-state width; currently only used
1027 for Higgs decay to <i>q qbar</i>, <i>g g</i> or
1031 <a name="method59"></a>
1032 <p/><strong>ParticleDataEntry* ParticleData::particleDataEntryPtr(int id) </strong> <br/>
1033 returns a pointer to the <code>ParticleDataEntry</code> object.
1034 The methods in the next section can then be used to manipulate
1038 <h3>The ParticleDataEntry methods</h3>
1040 Most of the methods that can be applied to a single
1041 <code>ParticleDataEntry</code> object are almost identical with
1042 those used above for the <code>ParticleData</code>, except
1043 that the <code>id</code> argument is no longer needed to find
1044 the right entry in the table. By and large, this makes direct
1045 access to the <code>ParticleDataEntry</code> methods superfluous.
1046 There are a few methods that are unique to each class, however.
1047 Furthermore, to avoid some naming ambiguities, many methods that
1048 set values begin with <code>set</code>.
1050 <a name="method60"></a>
1051 <p/><strong>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.) </strong> <br/>
1053 <strong>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.) </strong> <br/>
1054 there are two alternative constructors, that both expect the
1055 properties of a particle as input. The first assumes that there
1056 is only one particle, thet latter that there is a
1057 particle-antiparticle pair (but if the antiparticle name is
1058 <code>void</code> one reverts back to the particle-only case).
1061 <a name="method61"></a>
1062 <p/><strong>ParticleDataEntry::~ParticleDataEntry </strong> <br/>
1063 the destructor is needed to delete any <code>ResonanceWidths</code>
1064 objects that have been created and linked to the respective particle.
1067 <a name="method62"></a>
1068 <p/><strong>void ParticleDataEntry::setDefaults() </strong> <br/>
1069 initialize some particle flags with default values, e.g. whether
1070 a particle is a resonance, may decay, or is visible. Is called from the
1071 constructors and from <code>setAll</code>.
1074 <a name="method63"></a>
1075 <p/><strong>void ParticleDataEntry::initPtr(ParticleData* particleDataPtrIn) </strong> <br/>
1076 initialize pointer back to the whole database (so that masses of
1077 decay products can be accessed, e.g.).
1080 <a name="method64"></a>
1081 <p/><strong>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.) </strong> <br/>
1082 change all the properties of the particle associated with a given
1086 <a name="method65"></a>
1087 <p/><strong>int ParticleDataEntry::id() </strong> <br/>
1088 the PDG identity code.
1091 <a name="method66"></a>
1092 <p/><strong>bool ParticleDataEntry::hasAnti() </strong> <br/>
1093 tell whether a separate antiparticle exists.
1096 <a name="method67"></a>
1097 <p/><strong>void ParticleDataEntry::setName(string name) </strong> <br/>
1099 <strong>void ParticleDataEntry::setAntiName(string antiName) </strong> <br/>
1101 <strong>void ParticleDataEntry::setNames(string name,string antiName) </strong> <br/>
1103 <strong>string ParticleDataEntry::name(int id = 1) </strong> <br/>
1104 set or get the particle or antiparticle name. Only the sign of
1105 <code>id</code> is needed to distinguish particle/antiparticle.
1108 <a name="method68"></a>
1109 <p/><strong>void ParticleDataEntry::setSpinType(int spinType) </strong> <br/>
1111 <strong>int ParticleDataEntry::spinType() </strong> <br/>
1112 set or get the particle spin type, i.e. <i>2 s + 1</i>, or 0 in some
1116 <a name="method69"></a>
1117 <p/><strong>void ParticleDataEntry::setChargeType(int chargeType) </strong> <br/>
1119 <strong>int ParticleDataEntry::chargeType(int id = 1) </strong> <br/>
1121 <strong>double ParticleDataEntry::charge(int id = 1) </strong> <br/>
1122 set or get the particle charge type, i.e. three times the charge,
1123 or the charge itself. Only the sign of <code>id</code> is needed
1124 to distinguish particle/antiparticle.
1127 <a name="method70"></a>
1128 <p/><strong>void ParticleDataEntry::setColType(int colType) </strong> <br/>
1130 <strong>int ParticleDataEntry::colType(int id = 1) </strong> <br/>
1131 set or get the particle colour type, 0 for singlet, 1 for triplet,
1132 -1 for antitriplet, 2 for octet. Only the sign of <code>id</code>
1133 is needed to distinguish particle/antiparticle.
1136 <a name="method71"></a>
1137 <p/><strong>void ParticleDataEntry::setM0(double m0) </strong> <br/>
1139 <strong>double ParticleDataEntry::m0() </strong> <br/>
1140 the nominal mass <i>m_0</i> (in GeV).
1143 <a name="method72"></a>
1144 <p/><strong>void ParticleDataEntry::setMWidth(double mWidth) </strong> <br/>
1146 <strong>double ParticleDataEntry::mWidth() </strong> <br/>
1147 the width <i>Gamma</i> of the Breit-Wigner distribution (in GeV).
1150 <a name="method73"></a>
1151 <p/><strong>void ParticleDataEntry::setMMin(double mMin) </strong> <br/>
1153 <strong>double ParticleDataEntry::mMin() </strong> <br/>
1154 the lower limit of the allowed mass range generated by the Breit-Wigner
1155 (in GeV). Has no meaning for particles without width, and would
1156 typically be 0 there.
1159 <a name="method74"></a>
1160 <p/><strong>void ParticleDataEntry::setMMax(double mMax) </strong> <br/>
1162 <strong>double ParticleDataEntry::mMax() </strong> <br/>
1163 the upper limit of the allowed mass range generated by the Breit-Wigner
1164 (in GeV). If <i>mMax < mMin</i> then no upper limit is imposed.
1165 Has no meaning for particles without width, and would typically
1169 <a name="method75"></a>
1170 <p/><strong>double ParticleDataEntry::m0Min() </strong> <br/>
1171 similar to <code>mMin()</code> above, except that for particles with
1172 no width the <code>m0(id)</code> value is returned.
1175 <a name="method76"></a>
1176 <p/><strong>double ParticleDataEntry::m0Max() </strong> <br/>
1177 similar to <code>mMax()</code> above, except that for particles with
1178 no width the <code>m0(id)</code> value is returned.
1181 <a name="method77"></a>
1182 <p/><strong>void ParticleDataEntry::setTau0(double tau0) </strong> <br/>
1184 <strong>double ParticleDataEntry::tau0() </strong> <br/>
1185 the nominal proper lifetime <i>tau_0</i> (in mm/c).
1188 <a name="method78"></a>
1189 <p/><strong>void ParticleDataEntry::setIsResonance(bool isResonance) </strong> <br/>
1191 <strong>bool ParticleDataEntry::isResonance() </strong> <br/>
1192 a flag telling whether a particle species are considered as a resonance
1193 or not. Here <a href="ResonanceDecays.html" target="page">"resonance"</a>
1194 is used as shorthand for any massive particle
1195 where the decay process should be counted as part of the hard process
1196 itself, and thus be performed before showers and other event aspects
1197 are added. Restrictions on allowed decay channels is also directly
1198 reflected in the cross section of simulated processes, while those of
1199 normal hadrons and other light particles are not.
1200 In practice, it is reserved for states above the <i>b bbar</i>
1201 bound systems in mass, i.e. for <i>W, Z, t</i>, Higgs states,
1202 supersymmetric states and (most?) other states in any new theory.
1203 All particles with <code>m0</code> above 20 GeV are by default
1204 initialized to be considered as resonances.
1207 <a name="method79"></a>
1208 <p/><strong>void ParticleDataEntry::setMayDecay(bool mayDecay) </strong> <br/>
1210 <strong>bool ParticleDataEntry::mayDecay() </strong> <br/>
1211 a flag telling whether a particle species may decay or not, offering
1212 the main user switch. Whether a given particle of this kind then actually
1213 will decay also depends on it having allowed decay channels, and on
1214 other flags for <a href="ParticleDecays.html" target="page">particle decays</a>.
1215 All particles with <code>tau0</code> below 1000 mm are
1216 by default initialized to allow decays.
1219 <a name="method80"></a>
1220 <p/><strong>void ParticleDataEntry::setDoExternalDecays(bool doExternalDecays) </strong> <br/>
1222 <strong>bool ParticleDataEntry::doExternalDecay() </strong> <br/>
1223 a flag telling whether a particle should be handled by an external
1224 decay package or not, with the latter default. Can be manipulated as
1225 described on this page, but should normally not be. Instead the
1226 <code><a href="ExternalDecays.html" target="page">pythia.decayPtr</a></code>
1227 method should be provided with the list of relevant particles.
1230 <a name="method81"></a>
1231 <p/><strong>void ParticleDataEntry::setIsVisible(bool isVisible) </strong> <br/>
1233 <strong>bool ParticleDataEntry::isVisible() </strong> <br/>
1234 a flag telling whether a particle species is to be considered as
1235 visible in a detector or not, as used e.g. in analysis routines.
1236 By default this includes neutrinos and a few BSM particles
1237 (gravitino, sneutrinos, neutralinos) that have neither strong nor
1238 electromagnetic charge, and are not made up of constituents that
1239 have it. The value of this flag is only relevant if a particle is
1240 long-lived enough actually to make it to a detector.
1243 <a name="method82"></a>
1244 <p/><strong>void ParticleDataEntry::setDoForceWidth(bool doForceWidth) </strong> <br/>
1246 <strong>bool ParticleDataEntry::doForceWidth() </strong> <br/>
1247 a flag valid only for resonances where PYTHIA contains code to
1248 calculate the width of the resonance from encoded matrix-element
1249 expressions, i.e. the <i>Z^0</i>, <i>W^+-</i>, <i>t</i>,
1250 <i>h^0</i>, and a few more. The normal behaviour (<code>false</code>)
1251 is then that the width is calculated from the mass, but it is
1252 possible to <a href="ResonanceDecays.html" target="page">force</a> the resonance
1253 to retain the nominal width. Branching ratios and the running of the
1254 total width are unaffected.
1257 <a name="method83"></a>
1258 <p/><strong>void ParticleDataEntry::setHasChanged(bool hasChanged) </strong> <br/>
1260 <a name="method84"></a>
1261 <p/><strong>void ParticleDataEntry::hasChanged(bool hasChanged) </strong> <br/>
1262 keep track of whether the data for a particle has been changed
1263 in any respect between initialization and the current status.
1264 Is used e.g. by the <code>ParticleData::listChanged</code> method
1265 to determine which particles to list.
1268 <a name="method85"></a>
1269 <p/><strong>void ParticleDataEntry::initBWmass() </strong> <br/>
1270 Prepare the Breit-Wigner mass selection by precalculating
1271 frequently-used expressions.
1274 <a name="method86"></a>
1275 <p/><strong>double ParticleDataEntry::constituentMass() </strong> <br/>
1276 is the constituent mass for a quark, hardcoded as
1277 <i>m_u = m_d = 0.325</i>, <i>m_s = 0.50</i>, <i>m_c = 1.60</i>
1278 and <i>m_b = 5.0</i> GeV, for a diquark the sum of quark constituent
1279 masses, and for everything else the same as the ordinary mass.
1282 <a name="method87"></a>
1283 <p/><strong>double ParticleDataEntry::mass() </strong> <br/>
1284 give the mass of a particle, either at the nominal value
1285 or picked according to a (linear or quadratic) Breit-Wigner.
1288 <a name="method88"></a>
1289 <p/><strong>double ParticleDataEntry::mRun(double mH) </strong> <br/>
1290 calculate the running quark mass at a hard scale <code>mH</code>.
1291 For other particles the on-shell mass is given.
1294 <a name="method89"></a>
1295 <p/><strong>bool ParticleDataEntry::useBreitWigner() </strong> <br/>
1296 tells whether a particle will have a Breit-Wigner mass distribution or
1297 not. Is determined by an internal logic based on the particle width and
1298 on the value of the <code><a href="ParticleData.html" target="page">
1299 ParticleData:modeBreitWigner</a></code> switch.
1302 <a name="method90"></a>
1303 <p/><strong>bool ParticleDataEntry::canDecay(int id) </strong> <br/>
1304 true for a particle with at least one decay channel defined.
1307 <a name="method91"></a>
1308 <p/><strong>bool ParticleDataEntry::isLepton() </strong> <br/>
1309 true for a lepton or an antilepton (including neutrinos).
1312 <a name="method92"></a>
1313 <p/><strong>bool ParticleDataEntry::isQuark() </strong> <br/>
1314 true for a quark or an antiquark.
1317 <a name="method93"></a>
1318 <p/><strong>bool ParticleDataEntry::isGluon() </strong> <br/>
1322 <a name="method94"></a>
1323 <p/><strong>bool ParticleDataEntry::isDiquark() </strong> <br/>
1324 true for a diquark or antidiquark.
1327 <a name="method95"></a>
1328 <p/><strong>bool ParticleDataEntry::isHadron() </strong> <br/>
1329 true for a hadron (made up out of normal quarks and gluons,
1330 i.e. not for R-hadrons and other exotic states).
1333 <a name="method96"></a>
1334 <p/><strong>bool ParticleDataEntry::isMeson() </strong> <br/>
1338 <a name="method97"></a>
1339 <p/><strong>bool ParticleDataEntry::isBaryon() </strong> <br/>
1340 true for a baryon or antibaryon.
1343 <a name="method98"></a>
1344 <p/><strong>bool ParticleDataEntry::isOctetHadron() </strong> <br/>
1345 true for an intermediate hadron-like state with a colour octet charge
1346 as used in the colour octet model for
1347 <a href="OniaProcesses.html" target="page">onia</a> production.
1350 <a name="method99"></a>
1351 <p/><strong>int ParticleDataEntry::heaviestQuark(int id) </strong> <br/>
1352 extracts the heaviest quark or antiquark, i.e. one with largest
1353 <code>id</code> number, for a hadron. Only the sign of the input
1354 argument is relevant.
1357 <a name="method100"></a>
1358 <p/><strong>int ParticleDataEntry::baryonNumberType(int id) </strong> <br/>
1359 is 1 for a quark, 2 for a diquark, 3 for a baryon, the same with a
1360 minus sign for antiparticles, and else zero. Only the sign of the
1361 input argument is relevant.
1364 <a name="method101"></a>
1365 <p/><strong>void ParticleDataEntry::clearChannels() </strong> <br/>
1366 resets to an empty decay table.
1369 <a name="method102"></a>
1370 <p/><strong>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,) </strong> <br/>
1371 adds a decay channel with up to 8 products.
1374 <a name="method103"></a>
1375 <p/><strong>int ParticleDataEntry::sizeChannels() </strong> <br/>
1376 returns the number of decay channels for a particle.
1379 <a name="method104"></a>
1380 <p/><strong>DecayChannel& ParticleDataEntry::channel(int i) </strong> <br/>
1382 <strong>const DecayChannel& ParticleDataEntry::channel(int i) </strong> <br/>
1383 gain access to a specified channel in the decay table.
1386 <a name="method105"></a>
1387 <p/><strong>void ParticleDataEntry::rescaleBR(double newSumBR = 1.) </strong> <br/>
1388 rescales all partial branching ratios by a common factor, such that
1389 the sum afterward becomes <code>newSumBR</code>.
1392 <a name="method106"></a>
1393 <p/><strong>bool ParticleDataEntry::preparePick(int idSgn, double mHat = 0., int idInFlav = 0) </strong> <br/>
1394 prepare to pick a decay channel.
1397 <a name="method107"></a>
1398 <p/><strong>DecayChannel& ParticleDataEntry::pickChannel() </strong> <br/>
1399 pick a decay channel according to branching ratios from
1400 <code>preparePick</code>.
1403 <a name="method108"></a>
1404 <p/><strong>void ParticleDataEntry::setResonancePtr(ResonanceWidths* resonancePtr) </strong> <br/>
1406 <strong>ResonanceWidths* ParticleDataEntry::getResonancePtr() </strong> <br/>
1407 set or get a pointer to an object that can be used for dynamic calculation
1408 of partial and total resonance widths. Here a resonance is a particle
1409 such as top, <i>Z^0</i>, <i>W^+-</i>, Higgs, and new unstable states
1410 beyond the Standard Model.
1413 <a name="method109"></a>
1414 <p/><strong>void ParticleDataEntry::resInit(Info* infoPtrIn, Settings* settingsPtrIn, ParticleData* particleDataPtrIn, CoupSM* coupSMPtrIn) </strong> <br/>
1415 initialize the treatment of a resonance.
1418 <a name="method110"></a>
1419 <p/><strong>double ParticleDataEntry::resWidth(int idSgn,double mHat, int idInFlav = 0, bool openOnly = false, bool setBR = false) </strong> <br/>
1420 calculate the total with for a resonance of a given current mass,
1421 optionally including coupling to incoming flavour state (consider
1422 the <i>gamma*/Z^0</i> combination), optionally excluding decay
1423 channels that have been closed by the user, and optionally storing
1424 the results in the normal decay table. For the first argument only
1425 the sign is relevant.
1428 <a name="method111"></a>
1429 <p/><strong>double ParticleDataEntry::resWidthOpen(int idSgn,double mHat, int idInFlav = 0) </strong> <br/>
1430 special case of <code>resWidth</code>, where only open channels are
1431 included, but results are not stored in the normal decay table.
1434 <a name="method112"></a>
1435 <p/><strong>double ParticleDataEntry::resWidthStore(int idSgn,double mHat, int idInFlav = 0) </strong> <br/>
1436 special case of <code>resWidth</code>, where only open channels are
1437 included, and results are stored in the normal decay table.
1440 <a name="method113"></a>
1441 <p/><strong>double ParticleDataEntry::resOpenFrac(int idSgn) </strong> <br/>
1442 calculate the fraction of the full branching ratio that is left
1443 open by the user choice of allowed decay channels.
1446 <a name="method114"></a>
1447 <p/><strong>double ParticleDataEntry::resWidthRescaleFactor() </strong> <br/>
1448 the factor used to rescale all partial widths in case the total
1449 width is being forced to a specific value by the user.
1452 <a name="method115"></a>
1453 <p/><strong>double ParticleDataEntry::resWidthChan(double mHat, int idAbs1 = 0, int idAbs2 = 0) </strong> <br/>
1454 special case to calculate one final-state width; currently only used
1455 for Higgs decay to <i>q qbar</i>, <i>g g</i> or
1459 <h3>The DecayChannel methods</h3>
1461 The properties stored in an individual decay channel can be set or get
1462 by the methods in this section.
1464 <a name="method116"></a>
1465 <p/><strong>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) </strong> <br/>
1466 the constructor for a decay channel. Internal.
1469 <a name="method117"></a>
1470 <p/><strong>void DecayChannel::onMode(int onMode) </strong> <br/>
1472 <strong>int DecayChannel::onMode() </strong> <br/>
1473 set or get the <code>onMode</code> of a decay channel,<br/>
1474 0 if a channel is off,<br/>
1476 2 if on for a particle but off for an antiparticle,<br/>
1477 3 if on for an antiparticle but off for a particle.<br/>
1478 If a particle is its own antiparticle then 2 is on and 3 off
1479 but, of course, for such particles it is much simpler and safer
1480 to use only 1 and 0.<br/>
1481 The 2 and 3 options can be used e.g. to encode CP violation in
1482 B decays, or to let the <i>W</i>'s in a <i>q qbar -> W^+ W^-</i>
1483 process decay in different channels.
1486 <a name="method118"></a>
1487 <p/><strong>void DecayChannel::bRatio(double bRatio, bool countAsChanged = true) </strong> <br/>
1489 <strong>double DecayChannel::bRatio() </strong> <br/>
1490 set or get the branching ratio of the channel. Second argument only
1494 <a name="method119"></a>
1495 <p/><strong>void DecayChannel::rescaleBR(double fac) </strong> <br/>
1496 multiply the current branching ratio by <code>fac</code>.
1499 <a name="method120"></a>
1500 <p/><strong>void DecayChannel::meMode(int meMode) </strong> <br/>
1502 <strong>int DecayChannel::meMode() </strong> <br/>
1503 set or get the mode of processing this channel, possibly with matrix
1504 elements (see the <a href="ParticleDecays.html" target="page">particle decays</a>
1508 <a name="method121"></a>
1509 <p/><strong>void DecayChannel::multiplicity(int multiplicity) </strong> <br/>
1511 <strong>int DecayChannel::multiplicity() </strong> <br/>
1512 set or get the number of decay products in a channel, at most 8.
1513 (Is normally not to be set by hand, since it is automatically
1514 updated whenever the products list is changed.)
1517 <a name="method122"></a>
1518 <p/><strong>void DecayChannel::product(int i, int product) </strong> <br/>
1520 <strong>int DecayChannel::product(int i) </strong> <br/>
1521 set or get a list of the decay products, 8 products 0 <= i < 8,
1522 with trailing unused ones set to 0.
1525 <a name="method123"></a>
1526 <p/><strong>void DecayChannel::setHasChanged(bool hasChanged) </strong> <br/>
1528 <strong>bool DecayChannel::hasChanged() </strong> <br/>
1529 used for internal purposes, to know which decay modes have been changed.
1532 <a name="method124"></a>
1533 <p/><strong>bool DecayChannel::contains(int id1) </strong> <br/>
1535 <strong>bool DecayChannel::contains(int id1, int id2) </strong> <br/>
1537 <strong>bool DecayChannel::contains(int id1, int id2, int id3) </strong> <br/>
1538 find if the decay product list contains the one, two or three particle
1539 identities provided. If the same code is repeated then so must it be in
1540 the products list. Matching also requires correct sign.
1543 <a name="method125"></a>
1544 <p/><strong>void DecayChannel::currentBR(double currentBR) </strong> <br/>
1546 <strong>double DecayChannel::currentBR() </strong> <br/>
1547 set or get the current branching ratio, taking into account on/off
1548 switches and dynamic width for resonances. For internal use.
1551 <a name="method126"></a>
1552 <p/><strong>void DecayChannel::onShellWidth(double onShellWidth) </strong> <br/>
1554 <strong>double DecayChannel::onShellWidth() </strong> <br/>
1555 set or get the current partial width of the channel; intended for
1556 resonances where the widhts are recalculated based on the current
1557 resonance mass. For internal use.
1560 <a name="method127"></a>
1561 <p/><strong>void DecayChannel::onShellWidthFactor(double factor) </strong> <br/>
1562 multiply the current partial width by <code>factor</code>.
1565 <a name="method128"></a>
1566 <p/><strong>void DecayChannel::openSec(int idSgn, double openSecIn) </strong> <br/>
1568 <strong>double DecayChannel::openSec(nt idSgn) </strong> <br/>
1569 set or get the fraction of secondary open widths, separately for
1570 positive and negative particles. For internal use.
1576 <!-- Copyright (C) 2010 Torbjorn Sjostrand -->