next up previous contents
Next: Event Study and Analysis Up: The Fragmentation and Decay Previous: Parameter values   Contents

Examples

A 10 GeV $\u $ quark jet going out along the $+z$ axis is generated with

      CALL PY1ENT(0,2,10D0,0D0,0D0)
Note that such a single jet is not required to conserve energy, momentum or flavour. In the generation scheme, particles with negative $p_z$ are produced as well, but these are automatically rejected unless MSTJ(3) = -1. While frequently used in former days, the one-jet generation option is not of much current interest.

In e.g. a leptoproduction event a typical situation could be a $\u $ quark going out in the $+z$ direction and a $\u\d _0$ target remnant essentially at rest. (Such a process can be simulated by PYTHIA, but here we illustrate how to do part of it yourself.) The simplest procedure is probably to treat the process in the c.m. frame and boost it to the lab frame afterwards. Hence, if the c.m. energy is 20 GeV and the boost $\beta_z = 0.996$ (corresponding to $x_B = 0.045$), then

      CALL PY2ENT(0,2,2101,20D0)
      CALL PYROBO(0,0,0D0,0D0,0D0,0D0,0.996D0)
The jets could of course also be defined and allowed to fragment in the lab frame with
      CALL PY1ENT(-1,2,223.15D0,0D0,0D0)
      CALL PY1ENT(2,12,0.6837D0,3.1416D0,0D0)
      CALL PYEXEC
Note here that the target diquark is required to move in the backwards direction with $E-p_z = m_{\mathrm{p}}(1-x_B)$ to obtain the correct invariant mass for the system. This is, however, only an artifact of using a fixed diquark mass to represent a varying target remnant mass, and is of no importance for the fragmentation. If one wants a nicer-looking event record, it is possible to use the following
      CALL PY1ENT(-1,2,223.15D0,0D0,0D0)
      MSTU(10)=1
      P(2,5)=0.938D0*(1D0-0.045D0)
      CALL PY1ENT(2,2101,0D0,0D0,0D0)
      MSTU(10)=2
      CALL PYEXEC

A 30 GeV $\u\overline{\mathrm{u}}\mathrm{g}$ event with $E_{\u } = 8$ GeV and $E_{\overline{\mathrm{u}}} = 14$ GeV is simulated with

      CALL PY3ENT(0,2,21,-2,30D0,2D0*8D0/30D0,2D0*14D0/30D0)
The event will be given in a standard orientation with the $\u $ quark along the $+z$ axis and the $\overline{\mathrm{u}}$ in the $\pm z, +x$ half plane. Note that the flavours of the three partons have to be given in the order they are found along a string, if string fragmentation options are to work. Also note that, for 3-jet events, and particularly 4-jet ones, not all setups of kinematical variables $x$ lie within the kinematically allowed regions of phase space.

All common-block variables can obviously be changed by including the corresponding common block in the user-written main program. Alternatively, the routine PYGIVE can be used to feed in values, with some additional checks on array bounds then performed. A call

      CALL PYGIVE('MSTJ(21)=3;PMAS(C663,1)=210.;CHAF(401,1)=funnyino;'//
     &'PMAS(21,4)=')
will thus change the value of MSTJ(21) to 3, the value of PMAS(PYCOMP(663),1) = PMAS(136,1) to 210., the value of CHAF(401,1) to 'funnyino', and print the current value of PMAS(21,4). Since old and new values of parameters changed are written to output, this may offer a convenient way of documenting non-default values used in a given run. On the other hand, if a variable is changed back and forth frequently, the resulting voluminous output may be undesirable, and a direct usage of the common blocks is then to be recommended (the output can also be switched off, see MSTU(13)).

A general rule of thumb is that none of the physics routines (PYSTRF, PYINDF, PYDECY, etc.) should ever be called directly, but only via PYEXEC. This routine may be called repeatedly for one single event. At each call only those entries that are allowed to fragment or decay, and have not yet done so, are treated. Thus

      CALL PY2ENT(1,1,-1,20D0)       ! fill 2 partons without fragmenting
      MSTJ(1)=0                      ! inhibit jet fragmentation
      MSTJ(21)=0                     ! inhibit particle decay
      MDCY(PYCOMP(111),1)=0          ! inhibit pi0 decay
      CALL PYEXEC                    ! will not do anything
      MSTJ(1)=1                      !
      CALL PYEXEC                    ! partons will fragment, no decays
      MSTJ(21)=2                     !
      CALL PYEXEC                    ! particles decay, except pi0
      CALL PYEXEC                    ! nothing new can happen
      MDCY(PYCOMP(111),1)=1          !
      CALL PYEXEC                    ! pi0's decay

A partial exception to the rule above is PYSHOW. Its main application is for internal use by PYEEVT, PYDECY, and PYEVNT, but it can also be directly called by you. Note that a special format for storing colour-flow information in K(I,4) and K(I,5) must then be used. For simple cases, the PY2ENT can be made to take care of that automatically, by calling with the first argument negative.

      CALL PY2ENT(-1,1,-2,80D0)      ! store d ubar with colour flow
      CALL PYSHOW(1,2,80D0)          ! shower partons
      CALL PYEXEC                    ! subsequent fragmentation/decay
For more complicated configurations, PYJOIN should be used.

It is always good practice to list one or a few events during a run to check that the program is working as intended. With

      CALL PYLIST(1)
all particles will be listed and in addition total charge, momentum and energy of stable entries will be given. For string fragmentation these quantities should be conserved exactly (up to machine precision errors), and the same goes when running independent fragmentation with one of the momentum conservation options. PYLIST(1) gives a format that comfortably fits in an 80 column window, at the price of not giving the complete story. With PYLIST(2) a more extensive listing is obtained, and PYLIST(3) also gives vertex information. Further options are available, like PYLIST(12), which gives a list of particle data.

An event, as stored in the PYJETS common block, will contain the original partons and the whole decay chain, i.e. also particles which subsequently decayed. If parton showers are used, the amount of parton information is also considerable: first the on-shell partons before showers have been considered, then a K(I,1) = 22 line with total energy of the showering subsystem, after that the complete shower history tree-like structure, starting off with the same initial partons (now off-shell), and finally the end products of the shower rearranged along the string directions. This detailed record is useful in many connections, but if one only wants to retain the final particles, superfluous information may be removed with PYEDIT. Thus e.g.

      CALL PYEDIT(2)
will leave you with the final charged and neutral particles, except for neutrinos (and some other weakly interacting, neutral particles).

The information in PYJETS may be used directly to study an event. Some useful additional quantities derived from these, such as charge and rapidity, may easily be found via the PYK and PYP functions. Thus electric charge = PYP(I,6) (as integer, three times charge = PYK(I,6)) and true rapidity $y$ with respect to the $z$ axis = PYP(I,17).

A number of utility (MSTU, PARU) and physics (MSTJ, PARJ) switches and parameters are available in common block PYDAT1. All of these have sensible default values. Particle data is stored in common blocks PYDAT2, PYDAT3 and PYDAT4. Note that the data in the arrays KCHG, PMAS, MDCY CHAF and MWID is not stored by KF code, but by the compressed code KC. This code is not to be learnt by heart, but instead accessed via the conversion function PYCOMP, KC = PYCOMP(KF).

In the particle tables, the following particles are considered stable: the photon, $\mathrm{e}^{\pm}$, $\mu^{\pm}$, $\pi^{\pm}$, $\mathrm{K}^{\pm}$, $\mathrm{K}_{\mathrm{L}}^0$, $\mathrm{p}$, $\overline{\mathrm{p}}$, $\mathrm{n}$, $\overline{\mathrm{n}}$ and all the neutrinos. It is, however, always possible to inhibit the decay of any given particle by putting the corresponding MDCY value zero or negative, e.g. MDCY(PYCOMP(310),1) = 0 makes $\mathrm{K}_{\mathrm{S}}^0$ and MDCY(PYCOMP(3122),1) = 0 $\Lambda$ stable. It is also possible to select stability based on the average lifetime (see MSTJ(22)), or based on whether the decay takes place within a given spherical or cylindrical volume around the origin.

The Field-Feynman jet model [Fie78] is available in the program by changing the following values: MSTJ(1) = 2 (independent fragmentation), MSTJ(3) = -1 (retain particles with $p_z < 0$; is not mandatory), MSTJ(11) = 2 (choice of longitudinal fragmentation function, with the $a$ parameter stored in PARJ(51) - PARJ(53)), MSTJ(12) = 0 (no baryon production), MSTJ(13) = 1 (give endpoint quarks $p_{\perp}$ as quarks created in the field), MSTJ(24) = 0 (no mass broadening of resonances), PARJ(2) = 0.5 ($\mathrm{s}/ \u $ ratio for the production of new $\mathrm{q}\overline{\mathrm{q}}$ pairs), PARJ(11) = PARJ(12) = 0.5 (probability for mesons to have spin 1) and PARJ(21) = 0.35 (width of Gaussian transverse momentum distribution). In addition only $\d $, $\u $ and $\mathrm{s}$ single quark jets may be generated following the FF recipe. Today the FF `standard jet' concept is probably dead and buried, so the numbers above should more be taken as an example of the flexibility of the program, than as something to apply in practice.

A wide range of independent fragmentation options are implemented, to be accessed with the master switch MSTJ(1) = 2. In particular, with MSTJ(2) = 1 a gluon jet is assumed to fragment like a random $\d $, $\overline{\mathrm{d}}$, $\u $, $\overline{\mathrm{u}}$, $\mathrm{s}$ or $\overline{\mathrm{s}}$ jet, while with MSTJ(2) = 3 the gluon is split into a $\d\overline{\mathrm{d}}$, $\u\overline{\mathrm{u}}$ or $\mathrm{s}\overline{\mathrm{s}}$ pair of jets sharing the energy according to the Altarelli-Parisi splitting function. Whereas energy, momentum and flavour is not explicitly conserved in independent fragmentation, a number of options are available in MSTJ(3) to ensure this `post facto', e.g. MSTJ(3) = 1 will boost the event to ensure momentum conservation and then (in the c.m. frame) rescale momenta by a common factor to obtain energy conservation, whereas MSTJ(3) = 4 rather uses a method of stretching the jets in longitudinal momentum along the respective jet axis to keep angles between jets fixed.


next up previous contents
Next: Event Study and Analysis Up: The Fragmentation and Decay Previous: Parameter values   Contents
Stephen_Mrenna 2012-10-24