External Decays
DecayHandler is a base class for the external handling of
decays. It is intended for normal particle decays, primarily
B mesons and tau, and cannot be used to redirect
decays of heavy resonances like t or Z^0.
The user-written derived class is called if a pointer to it has
been given with the
pythia.decayPtr()
method, where it also is specified which particles it will be called for.
This particle information is accessible with the
doExternalDecay()
method.
There is only one pure virtual method in DecayHandler,
to do the decay:
virtual bool DecayHandler::decay(vector<int>& idProd, vector<double>& mProd, vector<Vec4>& pProd, int iDec, const Event& event)
where
argument idProd : is a list of particle PDG identity codes,
argument mProd : is a list of their respective masses (in GeV), and
argument pProd : is a list of their respective four-momenta.
At input, these vectors each have size one, so that idProd[0],
mProd[0] and pProd[0] contain information on the
particle that is to be decayed. At output, the vectors should have
increased by the addition of all the decay products. Even if initially
defined in the rest frame of the mother, the products should have been
boosted so that their four-momenta add up to the pProd[0] of
the decaying particle.
Should it be of interest to know the prehistory of the decaying
particle, e.g. to set some helicity information affecting the
decay angular distribution, the full event record is available
read-only, with info in which slot iDec the decaying particle
is stored.
The routine should return true if it managed the decay and
false otherwise, in which case Pythia will try
to do the decay itself. This e.g. means you can choose to do some decay
channels yourself, and leave others to Pythia. To avoid
double-counting, the channels you want to handle should be switched off
in the Pythia particle database. In the beginning of the
external decay method you should then return
false with a probability given by the sum of the branching
ratios for those channels you do not want to handle yourself.
Note that the decay vertex is always set by Pythia, and that
B-Bbar oscillations have already been taken into account,
if they were switched on. Thus idProd[0] may be the opposite
of event[iDec].id(), where the latter provides the code at
production.
A sample test program is available in main17.cc, providing
a simple example of how to use this facility.