User Hooks
Sometimes it may be convenient to step in during the generation
process: to modify the built-in cross sections, to veto undesirable
events or simply to collect statistics at various stages of the
evolution. There is a base class UserHooks
that gives
you this access at a few selected places. This class in itself does
nothing; the idea is that you should write your own derived class
for your task. One simple derived class (SuppressSmallPT
)
comes with the program, mainly as illustration, and the
main10.cc
program provides a complete (toy) example how
a derived class could be set up and used.
There are six sets of routines, that give you different kinds of
freedom. They are, in no particular order:
(i) Ones that give you access to the event record in between
the process-level and parton-level steps, or in between the
parton-level and hadron-level ones. You can study the event record
and decide whether to veto this event.
(ii) Ones that allow you to set a scale at which the combined
parton-level MI+ISR+FSR downwards evolution in pT is
temporarily interrupted, so the event can be studied and either
vetoed or allowed to continue the evolution.
(iii) Ones that allow you to to study the event after the first
few ISR/FSR emissions, or first few MI, so the event can be vetoed
or allowed to continue the evolution.
(iv) Ones that allow you to study the latest initial- or
final-state emission and veto that emission, without vetoing the
event as a whole.
(v) Ones that give you access to the properties of the trial
hard process, so that you can modify the internal Pythia cross section
by your own correction factors.
(vi) Ones that let you set the scale of shower evolution,
specifically for matching in resonance decays.
They are described further in the following numbered subsections.
All the possibilities above can be combined freely and also be combined
with the standard flags. An event would then survive only if it survived
each of the possible veto methods. There are no hidden interdependencies
in this game, but of course some combinations may not be particularly
meaningful. For instance, if you set PartonLevel:all = off
then the doVetoPT(...)
and doVetoPartonLevel(...)
locations in the code are not even reached, so they would never be called.
The effect of the vetoes of types (i), (ii) and (iii) can be studied
in the output of the
Pythia::statistics()
method. The "Selected" column represents the number of events that were
found acceptable by the internal Pythia machinery, whereas the "Accepted"
one are the events that also survived the user cuts. The cross section
is based on the latter number, and so is reduced by the amount associated
by the vetoed events. Also type (v) modifies the cross section, while
types (iv) and (vi) do not.
The basic components
For a derived UserHooks
class to be called during the
execution, a pointer to an object of this class should be handed in
with the
Pythia::setUserHooksPtr( UserHooks*)
method. The first step therefore is to construct your own derived
class, of course. This must contain a constructor and a destructor. The
initPtr
method comes "for free", and is set up without
any intervention from you.
UserHooks::UserHooks()
virtual UserHooks::~UserHooks()
The constructor and destructor do not need to do anything.
void UserHooks::initPtr( Info* infoPtr, Settings* settingsPtr, ParticleData* particleDataPtr,PartonSystems* partonSystemsPtr)
this (non-virtual) method is automatically called during the initialization
stage to set four useful pointers, and to set up the workEvent
below. The corresponding objects can later be used to extract some useful
information.
Info:
general event and run information, including some loop counters.
Settings:
the settings used to determine the character of the run.
Particle Data:
the particle data used in the event record
(including workEvent
below).
PartonSystems:
the list of partons that belong to each individual subcollision system.
Next you overload the desired methods listed in the sections below.
These often come in pairs or triplets, where the first must return
true for the last method to be called. This latter method typically
hands you a reference to the event record, which you then can use to
decide whether or not to veto. Often the event record can be quite
lengthy and difficult to overview. The following methods and data member
can then come in handy.
void UserHooks::omitResonanceDecays(const Event& process)
is a protected method that you can make use of in your own methods to
extract a simplified list of the hard process, where all resonance decay
chains are omitted. Intended for the can/doVetoProcessLevel
routines. Note that the normal process-level generation does include
resonance decays. That is, if a top quark is produced in the hard process,
then also decays such as t -> b W+, W+ -> u dbar will be generated
and stored in process
. The omitResonanceDecays
routine will take the input process
and copy it to
workEvent
(see below), minus the resonance decay chains.
All particles produced in the hard process, such as the top, will be
considered final-state ones, with positive status and no daughters,
just as it is before resonances are allowed to decay.
(In the PartonLevel
routines, these decay chains will
initially not be copied from process
to event
.
Instead the combined MI, ISR and FSR evolution is done with the top
above as final particle. Only afterwards will the resonance decay chains
be copied over, with kinematics changes reflecting those of the top, and
showers in the decays carried out.)
void UserHooks::subEvent(const Event& event, bool isHardest = true)
is a protected method that you can make use of in your own methods to
extract a brief list of the current
partons of interest, with all irrelevant ones omitted. For the default
isHardest = true
only the outgoing partons from the hardest
interaction (including the partons added to it by ISR and FSR) are
extracted, as relevant e.g. for doVetoPT( iPos, event)
with
iPos = 0 - 4
. With isHardest = false
instead
the outgoing partons of the latest "subprocess" are extracted, as relevant
when iPos = 5
, where it corresponds to the outgoing partons
in the currently considered decay. The resut is stored in
workEvent
below. The partons are stored in slots 0 through
workEvent.size() - 1
. The mother1()
and
mother2()
indices of a particle in this list both return 0.
The daughter1()
and daughter2()
indices both
return the index of the same parton in the original event record
(event
; possibly process
),
so that you can trace the full history, if of interest.
Event UserHooks::workEvent
This protected class member contains the outcome of the above
omitResonanceDecays(...)
and
subEvent(...)
methods. Alternatively you can use it for
whatever temporary purposes you wish. You are free to use standard
operations, e.g. to boost the event to its rest frame before analysis,
or remove particles that should not be analyzed.
The workEvent
can also be sent on to a
jet clustering algorithm.
(i) Interrupt between the main generation levels
virtual bool UserHooks::canVetoProcessLevel()
In the base class this method returns false. If you redefine it
to return true then the method doVetoProcessLevel(...)
will be called immediately after a hard process (and associated
resonance decays) has been selected and stored in the
process
event record.
At this stage, the process
record typically contains
the two beams in slots 1 and 2, the two incoming partons to the hard
process in slots 3 and 4, the N (usually 1, 2 or 3) primary produced
particles in slots 5 through 4 + N, and thereafter recursively the
resonance decay chains, if any. Use the method
omitResonanceDecays(...)
if you want to skip these
decay chains. There are exceptions to this structure,
for soft QCD processes (where
the partonic process may not yet have been selected at this stage),
and when a second hard process has
been requested (where two hard processes are bookkept). In general
it is useful to begin the development work by listing a few
process
records, to clarify what the structure is for
the cases of interest.
virtual bool UserHooks::doVetoProcessLevel(const Event& process)
can optionally be called, as described above. You can study, but not
modify, the process
event record of the hard process.
Based on that you can decide whether to veto the event, true, or let
it continue to evolve, false. If you veto, then this event is not
counted among the accepted ones, and does not contribute to the estimated
cross section. The Pytha::next()
method will begin a
completely new event, so the vetoed event will not appear in the
output of Pythia::next()
.
Note: the above veto is different from setting the flag
PartonLevel:all = off
.
Also in the latter case the event generation will stop after the process
level, but an event generated up to this point is considered perfectly
acceptable. It can be studied and it contributes to the cross section.
That is, PartonLevel:all = off
is intended for simple studies
of hard processes, where one can save a lot of time by not generating
the rest of the story. By contrast, the doVetoProcessLevel()
method allows you to throw away uninteresting events at an early stage
to save time, but those events that do survive the veto are allowed to
develop into complete final states (unless flags have been set otherwise).
virtual bool UserHooks::canVetoPartonLevel()
In the base class this method returns false. If you redefine it
to return true then the method doVetoPartonLevel(...)
will be called immediately after the parton level has been generated
and stored in the event
event record. Thus showers, multiple interactions and beam remnants
have been set up, but hadronization and decays have not yet been
performed. This is already a fairly complete event, possibly with quite
a complex parton-level history. Therefore it is usually only meaningful
to study the hardest interaction, e.g. using subEvent(...)
introduced above, or fairly generic properties, such as the parton-level
jet structure.
virtual bool UserHooks::doVetoPartonLevel(const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the event, true, or let
it continue to evolve, false. If you veto, then this event is not
counted among the accepted ones, and does not contribute to the estimated
cross section. The Pytha::next()
method will begin a
completely new event, so the vetoed event will not appear in the
output of Pythia::next()
.
Note: the above veto is different from setting the flag
HadronLevel:all = off
.
Also in the latter case the event generation will stop after the parton
level, but an event generated up to this point is considered perfectly
acceptable. It can be studied and it contributes to the cross section.
That is, HadronLevel:all = off
is intended for simple
studies of complete partonic states, where one can save time by not
generating the complete hadronic final state. By contrast, the
doVetoPartonLevel()
method allows you to throw away
uninteresting events to save time that way, but those events that
do survive the veto are allowed to develop into complete final states
(unless flags have been set otherwise).
(ii) Interrupt during the parton-level evolution, at a
pT scale
During the parton-level evolution, multiple interactions (MI),
initial-state radiation (ISR) and final-state radiation (FSR)
are normally evolved downwards in
one interleaved evolution sequence of decreasing pT values.
For some applications, e.g matrix-element-matching approaches, it
may be convenient to stop the evolution temporarily when the "hard"
emissions have been considered, but before continuing with the more
time-consuming soft activity. Based on these hard partons one can make
a decision whether the event at all falls in the intended event class,
e.g. has the "right" number of parton-level jets. If yes then, as for
the methods above, the evolution will continue all the way up to a
complete event. Also as above, if no, then the event will not be
considered in the final cross section.
Recall that the new or modified partons resulting from a MI, ISR or FSR
step are always appended to the end of the then-current event record.
Previously existing partons are not touched, except for the
status, mother and daughter
values, which are updated to reflect the modified history. It is
therefore straightforward to find the partons associated with the most
recent occurence.
An MI results in four new partons being appended, two incoming
and two outgoing ones.
An ISR results in the whole affected system being copied down,
with one of the two incoming partons being replaced by a new one, and
one more outgoing parton.
An FSR results in three new partons, two that come from the
branching and one that takes the recoil.
The story becomes more messy when rescattering is allowed as part
of the MI machinery. Then there will not only be a new system, as
outlined above, but additionally some existing systems will undergo
cascade effects, and be copied down with changed kinematics.
In this subsection we outline the possibility to interrupt at a given
pT scale, in the next to interrupt after a given number of
emissions.
virtual bool UserHooks::canVetoPT()
In the base class this method returns false. If you redefine it
to return true then the method doVetoPT(...)
will
interrupt the downward evolution at scaleVetoPT()
.
virtual double UserHooks::scaleVetoPT()
In the base class this method returns 0. You should redefine it
to return the pT scale at which you want to study the event.
virtual bool UserHooks::doVetoPT(int iPos, const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the event, true, or let
it continue to evolve, false. If you veto, then this event is not
counted among the accepted ones, and does not contribute to the estimated
cross section. The Pytha::next()
method will begin a
completely new event, so the vetoed event will not appear in the
output of Pythia::next()
.
argument
iPos : is the position/status when the routine is
called, information that can help you decide your course of action:
argumentoption
0 : when no MI, ISR or FSR occured above the veto scale;
argumentoption
1 : when inside the interleaved MI + ISR + FSR evolution,
after an MI process;
argumentoption
2 : when inside the interleaved MI + ISR + FSR evolution,
after an ISR emission;
argumentoption
3 : when inside the interleaved MI + ISR + FSR evolution,
after an FSR emission;
argumentoption
4 : for the optional case where FSR is deferred from the
interleaved evolution and only considered separately afterward (then
alternative 3 would never occur);
argumentoption
5 : is for subsequent resonance decays, and is called once
for each decaying resonance in a chain such as t -> b W, W -> u dbar.
argument
event : the event record contains a list of all partons
generated so far, also including intermediate ones not part of the
"current final state", and also those from further multiple interactions.
This may not be desirable for comparisons with matrix-element calculations.
You may want to make use of the subEvent(...)
method below to
obtain a simplified event record workEvent
.
(iii) Interrupt during the parton-level evolution, after a step
These options are closely related to the ones above in section (ii), so
we do not repeat the introduction, nor the possibilities to study the
event record, also by using subEvent(...)
and
workEvent
.
What is different is that the methods in this section give access to the
event as it looks like after each of the first few steps in the downwards
evolution, irrespectively of the pT scales of these branchings.
Furthermore, it is here assumed that the focus normally is on the hardest
subprocess, so that ISR/FSR emissions associated with additional MI's
are not considered. For MI studies, however, a separate simpler
alternative is offered to consider the event after a given number
of interactions.
virtual bool UserHooks::canVetoStep()
In the base class this method returns false. If you redefine it
to return true then the method doVetoStep(...)
will
interrupt the downward ISR and FSR evolution the first
numberVetoStep()
times.
virtual int UserHooks::numberVetoStep()
Returns the number of steps n each of ISR and FSR, for the
hardest interaction, that you want to be able to study. That is,
the method will be called after the first n ISR emissions,
irrespective of the number of FSR ones at the time, and after the
first n FSR emissions, irespective of the number of ISR ones.
The number of steps defaults to the first one only, but you are free
to pick another value. Note that double diffraction is handled as two
separate Pomeron-proton collisions, and thus has two sequences of
emissions.
virtual bool UserHooks::doVetoStep(int iPos, int nISR, int nFSR, const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the event, true, or let
it continue to evolve, false. If you veto, then this event is not
counted among the accepted ones, and does not contribute to the estimated
cross section. The Pytha::next()
method will begin a
completely new event, so the vetoed event will not appear in the
output of Pythia::next()
.
argument
iPos : is the position/status when the routine is
called, information that can help you decide your course of action.
Agrees with options 2 - 5 of the doVetoPT(...)
routine
above, while options 0 and 1 are not relevant here.
argument
nISR : is the number of ISR emissions in the hardest
process so far. For resonance decays, iPos = 5
, it is 0.
argument
nFSR : is the number of FSR emissions in the hardest
process so far. For resonance decays, iPos = 5
, it is the
number of emissions in the currently studied system.
argument
event : the event record contains a list of all partons
generated so far, also including intermediate ones not part of the
"current final state", and also those from further multiple interactions.
This may not be desirable for comparisons with matrix-element calculations.
You may want to make use of the subEvent(...)
method above to
obtain a simplified event record.
virtual bool UserHooks::canVetoMIStep()
In the base class this method returns false. If you redefine it
to return true then the method doVetoMIStep(...)
will
interrupt the downward MI evolution the first
numberVetoMIStep()
times.
virtual int UserHooks::numberVetoMIStep()
Returns the number of steps in the MI evolution that you want to be
able to study, right after each new step has been taken and the
subcollision has been added to the event record. The number of steps
defaults to the first one only, but you are free to pick another value.
Note that the hardest interaction of an events counts as the first
multiple interaction. For most hard processes it thus at the first
step offers nothing not available with the VetoProcessLevel
functionality above. For the minimum-bias and diffractive systems the
hardest interaction is not selected at the process level, however, so
there a check after the first multiple interaction offers new
functionality. Note that double diffraction is handled as two separate
Pomeron-proton collisions, and thus has two sequences of interactions.
Also, if you have set up a second hard process then a check is made
after these first two, and the first interaction coming from the MI
machinery would have sequence number 3.
virtual bool UserHooks::doVetoMIStep(int nMI,const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the event, true, or let
it continue to evolve, false. If you veto, then this event is not
counted among the accepted ones, and does not contribute to the estimated
cross section. The Pytha::next()
method will begin a
completely new event, so the vetoed event will not appear in the
output of Pythia::next()
.
argument
nMI : is the number of MI subprocesses has occured
so far.
argument
event : the event record contains a list of all partons
generated so far, also including intermediate ones not part of the
"current final state", e.g. leftovers from the ISR and FSR evolution
of previously generated systems. The most recently added one has not
had time to radiate, of course.
(iv) Veto emissions
The methods in this group are intended to allow the veto of an emission
in ISR or FSR, without affecting the evolution in any other way.
If an emission is vetoed, the event record is "rolled back" to the
way it was before the emission occured, and the evolution in pT
is continued downwards from the rejected value. The decision can be
based on full knowledge of the kinematics of the branching.
virtual bool UserHooks::canVetoISREmission()
In the base class this method returns false. If you redefine it
to return true then the method doVetoISREmission(...)
will interrupt the initial-state shower immediately after each
emission and allow that emission to be vetoed.
virtual bool UserHooks::doVetoISREmission( const int sizeOld, const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the emission, true, or
not, false. If you veto, then the latest emission is removed from
the event record. In either case the evolution of the shower will
continue from the point where it was left off.
argument
sizeOld : is the size of the event record before the
latest emission was added to it. It will also become the new size if
the emission is vetoed.
argument
event : the event record contains a list of all partons
generated so far. Of special interest are the ones associated with the
most recent emission, which are stored in entries from sizeOld
through event.size() - 1
inclusive. If you veto the emission
these entries will be removed, and the history info in the remaining
partons will be restored to a state as if the emission had never occured.
virtual bool UserHooks::canVetoFSREmission()
In the base class this method returns false. If you redefine it
to return true then the method doVetoFSREmission(...)
will interrupt the final-state shower immediately after each
emission and allow that emission to be vetoed.
virtual bool UserHooks::doVetoFSREmission( const int sizeOld, const Event& event)
can optionally be called, as described above. You can study, but not
modify, the event
event record of the partonic process.
Based on that you can decide whether to veto the emission, true, or
not, false. If you veto, then the latest emission is removed from
the event record. In either case the evolution of the shower will
continue from the point where it was left off.
argument
sizeOld : is the size of the event record before the
latest emission was added to it. It will also become the new size if
the emission is vetoed.
argument
event : the event record contains a list of all partons
generated so far. Of special interest are the ones associated with the
most recent emission, which are stored in entries from sizeOld
through event.size() - 1
inclusive. If you veto the emission
these entries will be removed, and the history info in the remaining
partons will be restored to a state as if the emission had never occured.
(v) Modify cross-sections
virtual bool UserHooks::canModifySigma()
In the base class this method returns false. If you redefine it
to return true then the method multiplySigmaBy(...)
will
allow you to modify the cross section weight assigned to the current
event.
virtual double UserHooks::multiplySigmaBy(const SigmaProcess* sigmaProcessPtr, const PhaseSpace* phaseSpacePtr, bool inEvent)
when called this method should provide the factor by which you want to
see the cross section weight of the current event modified by. If you
return unity then the normal cross section is obtained. Note that, unlike
the methods above, these modifications do not lead to a difference between
the number of "selected" events and the number of "accepted" ones,
since the modifications occur already before the "selected" level.
The integrated cross section of a process is modified, of course.
Note that the cross section is only modifiable for normal hard processes.
It does not affect the cross section in further multiple interactions,
nor in elastic/diffractive/minimum-bias events.
argument
sigmaProcessPtr, phaseSpacePtr : :
what makes this routine somewhat tricky to write is that the
hard-process event has not yet been constructed, so one is restricted
to use the information available in the phase-space and cross-section
objects currently being accessed. Which of their methods are applicable
depends on the process, in particular the number of final-state particles.
The multiplySigmaBy
code in UserHooks.cc
contains explicit instructions about which methods provide meaningful
information, and so offers a convenient starting point.
argument
inEvent : : this flag is true when the method is
called from within the event-generation machinery and false
when it is called at the initialization stage of the run, when the
cross section is explored to find a maximum for later Monte Carlo usage.
Cross-section modifications should be independent of this flag,
for consistency, but if multiplySigmaBy(...)
is used to
collect statistics on the original kinematics distributions before cuts,
then it is important to be able to exclude the initialization stage
from comparisons.
One derived class is supplied as an example how this facility can be used
to reweight cross sections in the same spirit as is done with QCD cross
sections for the minimum-bias/underlying-event description:
class
SuppressSmallPT : public UserHooks
suppress small-pT production for 2 -> 2 processes
only, while leaving other processes unaffected. The basic suppression
factor is pT^4 / ((k*pT0)^2 + pT^2)^2, where pT
refers to the current hard subprocess and pT0 is the same
energy-dependent dampening scale as used for
multiple interactions.
This class contains canModifySigma()
and
multiplySigmaBy()
methods that overload the base class ones.
SuppressSmallPT::SuppressSmallPT( double pT0timesMI = 1., int numberAlphaS = 0, bool useSameAlphaSasMI = true)
The optional arguments of the constructor provides further variability.
argument
pT0timesMI :
corresponds to the additional factor k in the above formula.
It is by default equal to 1 but can be used to explore deviations from
the expected value.
argument
numberAlphaS :
if this number n is bigger than the default 0, the
corresponding number of alpha_strong factors is also
reweighted from the normal renormalization scale to a modified one,
i.e. a further suppression factor
( alpha_s((k*pT0)^2 + Q^2_ren) / alpha_s(Q^2_ren) )^n
is introduced.
argument
useSameAlphaSasMI :
regulates which kind of new alpha_strong value is evaluated
for the numerator in the above expression. It is by default the same
as set for multiple interactions (i.e. same starting value at
M_Z and same order of running), but if false
instead the one for hard subprocesses. The denominator
alpha_s(Q^2_ren) is always the value used for the "original",
unweighted cross section.
(vi) Modify scale in shower evolution
The choice of maximum shower scale in resonance decays is normally not a
big issue, since the shower here is expected to cover the full phase
space. In some special cases a matching scheme is intended, where hard
radiation is covered by matrix elements, and only softer by showers. The
below two methods support such an approach. Note that the two methods
are not used in the TimeShower
class itself, but when
showers are called from the PartonLevel
generation. Thus
user calls directly to TimeShower
are not affected.
virtual bool UserHooks::canSetResonanceScale()
In the base class this method returns false. If you redefine it
to return true then the method scaleResonance(...)
will set the initial scale of downwards shower evolution.
virtual double UserHooks::scaleResonance( const int iRes, const Event& event)
can optionally be called, as described above. You should return the maximum
scale, in GeV, from which the shower evolution will begin. The base class
method returns 0, i.e. gives no shower evolution at all.
You can study, but not modify, the event
event record
of the partonic process to check which resonance is decaying, and into what.
argument
iRes : is the location in the event record of the
resonance that decayed to the particles that now will shower.
argument
event : the event record contains a list of all partons
generated so far, specifically the decaying resonance and its immediate
decay products.