Hadron-Level Standalone
The Les Houches Accord allows external process-level configurations
to be fed in, for subsequent parton-level and hadron-level generation
to be handled internally by PYTHIA. There is no correspondingly
standardized interface if you have external events that have also
been generated through the parton-level stage, so that only the
hadron-level remains to be handled. A non-standard way to achieve this
exists, however, and can be useful both for real applications and
for various tests of the hadronization model on its own.
The key trick is to set the flag ProcessLevel:all = off
.
When pythia.next()
is called it then does not try to
generate a hard process. Since there are no beams, it is also not
possible to perform the normal PartonLevel
step.
(It is still possible to generate final-state radiation, but this
is not automatic. It would have to be done by hand, using the
pythia.forceTimeShower(...)
method, before
pythia.next()
is called.)
Thus only the HadronLevel
methods are
called, to take the current content of the event record stored in
pythia.event
as a starting point for any hadronization
and decays that are allowed by the normal parameters of this step.
Often the input would consist solely of partons grouped into colour
singlets, but also (colour-singlet) particles are allowed.
To set up all the parameters, a pythia.init()
call has
to be used, without any arguments. In brief, the structure of the
main program therefore should be something like
Pythia pythia; // Declare generator.
Event& event = pythia.event // Convenient shorthand.
pythia.readString("ProcessLevel:all = off"); // The trick!
pythia.init(); // Initialization.
for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
// Insert filling of event here!
pythia.next(); // Do the hadron level.
}
Of course this should be supplemented by analysis of events, error checks,
and so on, as for a normal PYTHIA run. The unique aspect is how to fill
the event
inside the loop, before pythia.next()
is called.
Input configuration
To set up a new configuration the first step is to throw away the current
one, with event.reset()
. This routine will also reserve
the zeroth entry in the even record to represent the event as a whole.
With the event.append(...)
methods a new entry is added at the
bottom of the current record, i.e. the first time it is called entry
number 1 is filled, and so on. The append
method basically
exists in four variants, either without or with history information,
and with four-momentum provided either as a Vec4
four-vector
or as four individual components:
append( id, status, col, acol, p, m)
append( id, status, col, acol, px, py, pz, e, m)
append( id, status, mother1, mother2, daughter1, daughter2, col, acol, p, m)
append( id, status, mother1, mother2, daughter1, daughter2, col, acol, px, py, pz, e, m)
The methods return the index at which the entry has been stored,
but normally you would not use this feature.
You can find descriptions of the input variables
here.
The PDG particle code id
and the Les Houches Accord colour
col
and anticolour acol
tags must be set
correctly. The four-momentum and mass have to be provided in units of GeV;
if you omit the mass it defaults to 0.
Outgoing particles that should hadronize should be given status code 23.
Often this is the only status code you need. You could e.g. also fill in
incoming partons with -21 and intermediate ones with -22, if you so wish.
Usually the choice of status codes is not crucial, so long as you recall
that positive numbers correspond to particles that are still around, while
negative numbers denote ones that already hadronized or decayed. However,
so as not to run into contradictions with the internal PYTHIA checks
(when Check:event = on
), or with external formats such as
HepMC, we do recommend the above codes. When pythia.next()
is called the positive-status particles that hadronize/decay get the sign
of the status code flipped to negative but the absolute value is retained.
The new particles are added with normal PYTHIA status codes.
For normal hadronization/decays in pythia.next()
the
history encoded in the mother and daughter indices is not used.
Therefore the first two append
methods, which set all these
indices vanishing, should suffice. The subsequent hadronization/decays
will still be properly documented.
The exception is when you want to include junctions in your string
topology, i.e. have three string pieces meet. Then you must insert in
your event record the (decayed) particle that is the reason for the
presence of a junction, e.g. a baryon beam remnant from which several
valence quarks have been kicked out, or a neutralino that underwent a
baryon-number-violating decay. This particle must have as daughters
the three partons that together carry the baryon number.
The sample program in main21.cc
illustrates how you can work
with this facility, both for simple parton configurations and for more
complicated ones with junctions.
As an alternative to setting up a topology with the methods above,
a Les Houches Event File (LHEF)
can also provide the configurations. Since no beams or processes are
defined, only the <event>....</event>
blocks
need to be present, one for each event, so strictly it is not a true
LHEF. You need to select Beams:frameType = 4
, provide
the file name in Beams:LHEF
and, as above, set
ProcessLevel:all = off
. Needless to say, an externally
linked LHAup
class works as well as an LHEF,
with Beams:frameType = 5
.
The event information to store in the LHEF, or provide by the
LHAup
, is essentially the same as above. The only
difference is in status codes: outgoing particles should have 1
instead of 23, and intermediate resonances 2 instead of -22.
Incoming partons, if any, are -1 instead of -21.
Extensions to resonance decays
With the above scheme, pythia.next()
will generate
hadronization, i.e. string fragmentation and subsequent decays of
normal unstable particles. It will not decay
resonances, i.e.
W, Z, top, Higgs, SUSY and other massive particles.
If the decay products are already provided, of course those
products will be hadronized, but without any of the showers
that one would expect in Z^0 -> q qbar, say. That is, the
presence of a decayed resonance in the event record can be nice
for documentation purposes, but otherwise it plays no role.
It is possible to change this behaviour with the following flag.
flag
Standalone:allowResDec
(default = off
)
If off, then resonances are stable if final-state particles,
and irrelevant if intermediate particles. If on, then
resonances are decayed, sequentially if necessary. After each
decay step a parton shower is applied, and subsequent decay
kinematics is shifted by the recoil of such branchings.
If the decay chain and kinematics of the resonance is provided
at input, this information is used, otherwise it is generated
internally according to built-in branching ratios etc.
The input configuration has to follow the rules described above,
with ProcessLevel:all = off
. (Terminology not quite
consistent, since resonance decays normally is part of the
process-level step.) It is possible to combine several resonances,
and other coloured or uncoloured particles into the same event.
Final-state radiation is applied to each step of the decay
sequence by default, but can be switched off with
PartonLevel:FSR = off
or
PartonLevel:FSRinResonances
.
Repeated hadronization or decay
An alternative approach is possible with the
pythia.forceHadronLevel()
routine. This method does
a call to the HadronLevel
methods, irrespective of the
value of the HadronLevel:all
flag. If you hadronize
externally generated events it is equivalent to a
pythia.next()
call with
ProcessLevel:all = off
.
This method truly sticks to the hadron level, and thus cannot handle
resonance decays. You therefore must not mix it with the
Standalone:allowResDec = on
framework.
The similarity of names indicates that
pythia.forceTimeShower( int iBeg, int iEnd, double pTmax,
int nBranchMax = 0)
is intended to belong to the same set of
work-by-hand methods. Here iBeg
and iEnd
give the range of partons that should be allowed to shower,
pTmax
the maximum pT scale of emissions,
and a nonzero nBranchMax
a maximum number of allowed
branchings. Additionally, a scale
has to be set for each
parton that should shower, which requires an additional final argument
to the append
methods above. This scale limits the maximum
pT allowed for each parton, in addition to the global
pTmax
. When not set the scale defaults to 0, meaning no
radiation for that parton.
The real application instead is for repeated hadronization of the same
PYTHIA process- and parton-level event. This may for some studies
help to save time, given that these two first step are more
time-consuming than the hadronization one.
For repeated hadronization you should first generate an event as usual,
but with HadronLevel:all = off
. This event you can save
in a temporary copy, e.g. Event savedEvent = pythia.event
.
Inside a loop you copy back with pythia.event = savedEvent
,
and call pythia.forceHadronLevel()
to obtain a new
hadronization history.
A more limited form of repetition is if you want to decay a given kind
of particle repeatedly, without having to generate the rest of the event
anew. This could be the case e.g. in B physics applications.
Then you can use the pythia.moreDecays()
method, which
decays all particles in the event record that have not been decayed
but should have been done so. The
pythia.particleData.mayDecay( id, false/true)
method
may be used to switch off/on the decays of a particle species
id
, so that it is not decayed in the
pythia.next()
call but only inside a loop over a number of
tries.
Between each loop the newly produced decay products must be
removed and the decayed particle status restored to undecayed.
The former is simple, since the new products are appended to the
end of the event record: event.saveSize()
saves the
initial size of the event record, and event.restoreSize()
can later be used repeatedly to restore this original size, which means
that the new particles at the end are thrown away. The latter is more
complicated, and requires the user to identify the positions of all
particles of the species and restore a positive status code with
event[i].statusPos()
.
The main15.cc
program illustrates both these methods,
i.e. either repeated hadronization or repeated decay of PYTHIA
events.