Event Analysis
Introduction
The routines in this section are intended to be used to analyze
event properties. As such they are not part of the main event
generation chain, but can be used in comparisons between Monte
Carlo events and real data. They are rather free-standing, but
assume that input is provided in the PYTHIA 8
Event
format, and use a few basic facilities such
as four-vectors. Their ordering is mainly by history; for current
LHC applications the final one, SlowJet
, is of
special interest.
In addition to the methods presented here, there is also the
possibility to make use of external
jet finders .
Sphericity
The standard sphericity tensor is
S^{ab} = (sum_i p_i^a p_i^b) / (sum_i p_i^2)
where the sum i runs over the particles in the event,
a, b = x, y, z, and p without such an index is
the absolute size of the three-momentum . This tensor can be
diagonalized to find eigenvalues and eigenvectors.
The above tensor can be generalized by introducing a power
r, such that
S^{ab} = (sum_i p_i^a p_i^b p_i^{r-2}) / (sum_i p_i^r)
In particular, r = 1 gives a linear dependence on momenta
and thus a collinear safe definition, unlike sphericity.
To do sphericity analyses you have to set up a Sphericity
instance, and then feed in events to it, one at a time. The results
for the latest event are available as output from a few methods.
Sphericity::Sphericity(double power = 2., int select = 2)
create a sphericity analysis object, where
argument
power (default = 2.
) :
is the power r defined above, i.e.
argumentoption
2. : gives Sphericity, and
argumentoption
1. : gives the linear form.
argument
select (default = 2
) :
tells which particles are analyzed,
argumentoption
1 : all final-state particles,
argumentoption
2 : all observable final-state particles,
i.e. excluding neutrinos and other particles without strong or
electromagnetic interactions (the isVisible()
particle method), and
argumentoption
3 : only charged final-state particles.
bool Sphericity::analyze( const Event& event)
perform a sphericity analysis, where
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
If the routine returns false
the
analysis failed, e.g. if too few particles are present to analyze.
After the analysis has been performed, a few methods are available
to return the result of the analysis of the latest event:
double Sphericity::sphericity()
gives the sphericity (or equivalent if r is not 2),
double Sphericity::aplanarity()
gives the aplanarity (with the same comment),
double Sphericity::eigenValue(int i)
gives one of the three eigenvalues for i = 1, 2 or 3, in
descending order,
Vec4 Sphericity::eventAxis(i)
gives the matching normalized eigenvector, as a Vec4
with vanishing time/energy component.
void Sphericity::list()
provides a listing of the above information.
There is also one method that returns information accumulated for all
the events analyzed so far.
int Sphericity::nError()
tells the number of times analyze(...)
failed to analyze
events, i.e. returned false
.
Thrust
Thrust is obtained by varying the thrust axis so that the longitudinal
momentum component projected onto it is maximized, and thrust itself is
then defined as the sum of absolute longitudinal momenta divided by
the sum of absolute momenta. The major axis is found correspondingly
in the plane transverse to thrust, and the minor one is then defined
to be transverse to both. Oblateness is the difference between the major
and the minor values.
The calculation of thrust is more computer-time-intensive than e.g.
linear sphericity, introduced above, and has no specific advantages except
historical precedent. In the PYTHIA 6 implementation the search was
sped up at the price of then not being guaranteed to hit the absolute
maximum. The current implementation studies all possibilities, but at
the price of being slower, with time consumption for an event with
n particles growing like n^3.
To do thrust analyses you have to set up a Thrust
instance, and then feed in events to it, one at a time. The results
for the latest event are available as output from a few methods.
Thrust::Thrust(int select = 2)
create a thrust analysis object, where
argument
select (default = 2
) :
tells which particles are analyzed,
argumentoption
1 : all final-state particles,
argumentoption
2 : all observable final-state particles,
i.e. excluding neutrinos and other particles without strong or
electromagnetic interactions (the isVisible()
particle method), and
argumentoption
3 : only charged final-state particles.
bool Thrust::analyze( const Event& event)
perform a thrust analysis, where
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
If the routine returns false
the
analysis failed, e.g. if too few particles are present to analyze.
After the analysis has been performed, a few methods are available
to return the result of the analysis of the latest event:
double Thrust::thrust()
double Thrust::tMajor()
double Thrust::tMinor()
double Thrust::oblateness()
gives the thrust, major, minor and oblateness values, respectively,
Vec4 Thrust::eventAxis(int i)
gives the matching normalized event-axis vectors, for i = 1, 2 or 3
corresponding to thrust, major or minor, as a Vec4
with
vanishing time/energy component.
void Thrust::list()
provides a listing of the above information.
There is also one method that returns information accumulated for all
the events analyzed so far.
int Thrust::nError()
tells the number of times analyze(...)
failed to analyze
events, i.e. returned false
.
ClusterJet
ClusterJet
(a.k.a. LUCLUS
and
PYCLUS
) is a clustering algorithm of the type used for
analyses of e^+e^- events, see the PYTHIA 6 manual. All
visible particles in the events are clustered into jets. A few options
are available for some well-known distance measures. Cutoff
distances can either be given in terms of a scaled quadratic quantity
like y = pT^2/E^2 or an unscaled linear one like pT.
Note that we have deliberately chosen not to include the e^+e^-
equivalents of the Cambridge/Aachen and anti-kRT algorithms.
These tend to be good at clustering the densely populated (in angle)
cores of jets, but less successful for the sparsely populated transverse
regions, where many jets may come to consist of a single low-momentum
particle. In hadron collisions such jets could easily be disregarded,
while in e^+e^- annihilation all particles derive back from the
hard process.
To do jet finding analyses you have to set up a ClusterJet
instance, and then feed in events to it, one at a time. The results
for the latest event are available as output from a few methods.
ClusterJet::ClusterJet(string measure = "Lund", int select = 2, int massSet = 2, bool precluster = false, bool reassign = false)
create a ClusterJet
instance, where
argument
measure (default = "Lund"
) : distance measure,
to be provided as a character string (actually, only the first character
is necessary)
argumentoption
"Lund" : the Lund pT distance,
argumentoption
"JADE" : the JADE mass distance, and
argumentoption
"Durham" : the Durham kT measure.
argument
select (default = 2
) :
tells which particles are analyzed,
argumentoption
1 : all final-state particles,
argumentoption
2 : all observable final-state particles,
i.e. excluding neutrinos and other particles without strong or
electromagnetic interactions (the isVisible()
particle
method), and
argumentoption
3 : only charged final-state particles.
argument
massSet (default = 2
) : masses assumed for the particles
used in the analysis
argumentoption
0 : all massless,
argumentoption
1 : photons are massless while all others are
assigned the pi+- mass, and
argumentoption
2 : all given their correct masses.
argument
precluster (default = off
) :
perform or not a early preclustering step, where nearby particles
are lumped together so as to speed up the subsequent normal clustering.
argument
reassign (default = off
) :
reassign all particles to the nearest jet each time after two jets
have been joined.
ClusterJet::analyze( const Event& event, double yScale, double pTscale, int nJetMin = 1, int nJetMax = 0)
performs a jet finding analysis, where
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
argument
yScale :
is the cutoff joining scale, below which jets are joined. Is given
in quadratic dimensionless quantities. Either yScale
or pTscale
must be set nonvanishing, and the larger of
the two dictates the actual value.
argument
pTscale :
is the cutoff joining scale, below which jets are joined. Is given
in linear quantities, such as pT or m depending on
the measure used, but always in units of GeV. Either yScale
or pTscale
must be set nonvanishing, and the larger of
the two dictates the actual value.
argument
nJetMin (default = 1
) :
the minimum number of jets to be reconstructed. If used, it can override
the yScale
and pTscale
values.
argument
nJetMax (default = 0
) :
the maximum number of jets to be reconstructed. Is not used if below
nJetMin
. If used, it can override the yScale
and pTscale
values. Thus e.g.
nJetMin = nJetMax = 3
can be used to reconstruct exactly
3 jets.
If the routine returns false
the analysis failed,
e.g. because the number of particles was smaller than the minimum number
of jets requested.
After the analysis has been performed, a few ClusterJet
class methods are available to return the result of the analysis:
int ClusterJet::size()
gives the number of jets found, with jets numbered 0 through
size() - 1
.
Vec4 ClusterJet::p(int i)
gives a Vec4
corresponding to the four-momentum defined by
the sum of all the contributing particles to the i'th jet.
int ClusterJet::mult(int i)
the number of particles that have been clustered into the i'th jet.
int ClusterJet::jetAssignment(int i)
gives the index of the jet that the particle i of the event
record belongs to,
void ClusterJet::list()
provides a listing of the reconstructed jets.
int ClusterJet::distanceSize()
the number of most recent clustering scales that have been stored
for readout with the next method. Normally this would be five,
but less if fewer clustering steps occurred.
double ClusterJet::distance(int i)
clustering scales, with distance(0)
being the most
recent one, i.e. normally the highest, up to distance(4)
being the fifth most recent. That is, with n being the final
number of jets, ClusterJet::size()
, the scales at which
n+1 jets become n, n+2 become n+1,
and so on till n+5 become n+4. Nonexisting clustering
scales are returned as zero. The physical interpretation of a scale is
as provided by the respective distance measure (Lund, JADE, Durham).
There is also one method that returns information accumulated for all
the events analyzed so far.
int ClusterJet::nError()
tells the number of times analyze(...)
failed to analyze
events, i.e. returned false
.
CellJet
CellJet
(a.k.a. PYCELL
) is a simple cone jet
finder in the UA1 spirit, see the PYTHIA 6 manual. It works in an
(eta, phi, eT) space, where eta is pseudorapidity,
phi azimuthal angle and eT transverse energy.
It will draw cones in R = sqrt(Delta-eta^2 + Delta-phi^2)
around seed cells. If the total eT inside the cone exceeds
the threshold, a jet is formed, and the cells are removed from further
analysis. There are no split or merge procedures, so later-found jets
may be missing some of the edge regions already used up by previous
ones. Not all particles in the event are assigned to jets; leftovers
may be viewed as belonging to beam remnants or the underlying event.
It is not used by any experimental collaboration, but is closely
related to the more recent and better theoretically motivated
anti-kT algorithm [Cac08].
To do jet finding analyses you have to set up a CellJet
instance, and then feed in events to it, one at a time. Especially note
that, if you want to use the options where energies are smeared in
order so emulate detector imperfections, you must hand in an external
random number generator, preferably the one residing in the
Pythia
class. The results for the latest event are
available as output from a few methods.
CellJet::CellJet(double etaMax = 5., int nEta = 50, int nPhi = 32, int select = 2, int smear = 0, double resolution = 0.5, double upperCut = 2., double threshold = 0., Rndm* rndmPtr = 0)
create a CellJet
instance, where
argument
etaMax (default = 5.
) :
the maximum +-pseudorapidity that the detector is assumed to cover.
argument
nEta (default = 50
) :
the number of equal-sized bins that the +-etaMax range
is assumed to be divided into.
argument
nPhi (default = 32
) :
the number of equal-sized bins that the phi range
+-pi is assumed to be divided into.
argument
select (default = 2
) :
tells which particles are analyzed,
argumentoption
1 : all final-state particles,
argumentoption
2 : all observable final-state particles,
i.e. excluding neutrinos and other particles without strong or
electromagnetic interactions (the isVisible()
particle
method),
and
argumentoption
3 : only charged final-state particles.
argument
smear (default = 0
) :
strategy to smear the actual eT bin by bin,
argumentoption
0 : no smearing,
argumentoption
1 : smear the eT according to a Gaussian
with width resolution * sqrt(eT), with the Gaussian truncated
at 0 and upperCut * eT,
argumentoption
2 : smear the e = eT * cosh(eta) according
to a Gaussian with width resolution * sqrt(e), with the
Gaussian truncated at 0 and upperCut * e.
argument
resolution (default = 0.5
) :
see above.
argument
upperCut (default = 2.
) :
see above.
argument
threshold (default = 0 GeV
) :
completely neglect all bins with an eT < threshold.
argument
rndmPtr (default = 0
) :
the random-number generator used to select the smearing described
above. Must be handed in for smearing to be possible. If your
Pythia
class instance is named pythia
,
then &pythia.rndm
would be the logical choice.
bool CellJet::analyze( const Event& event, double eTjetMin = 20., double coneRadius = 0.7, double eTseed = 1.5)
performs a jet finding analysis, where
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
argument
eTjetMin (default = 20. GeV
) :
is the minimum transverse energy inside a cone for this to be
accepted as a jet.
argument
coneRadius (default = 0.7
) :
is the size of the cone in (eta, phi) space drawn around
the geometric center of the jet.
argument
eTseed (default = 1.5 GeV
) :
the minimum eT in a cell for this to be acceptable as
the trial center of a jet.
If the routine returns false
the analysis failed,
but currently this is not foreseen ever to happen.
After the analysis has been performed, a few CellJet
class methods are available to return the result of the analysis:
int CellJet::size()
gives the number of jets found, with jets numbered 0 through
size() - 1
,
double CellJet::eT(i)
gives the eT of the i'th jet, where jets have been
ordered with decreasing eT values,
double CellJet::etaCenter(int i)
double CellJet::phiCenter(int i)
gives the eta and phi coordinates of the geometrical
center of the i'th jet,
double CellJet::etaWeighted(int i)
double CellJet::phiWeighted(int i)
gives the eta and phi coordinates of the
eT-weighted center of the i'th jet,
int CellJet::multiplicity(int i)
gives the number of particles clustered into the i'th jet,
Vec4 CellJet::pMassless(int i)
gives a Vec4
corresponding to the four-momentum defined
by the eT and the weighted center of the i'th jet,
Vec4 CellJet::pMassive(int i)
gives a Vec4
corresponding to the four-momentum defined by
the sum of all the contributing cells to the i'th jet, where
each cell contributes a four-momentum as if all the eT is
deposited in the center of the cell,
double CellJet::m(int i)
gives the invariant mass of the i'th jet, defined by the
pMassive
above,
void CellJet::list()
provides a listing of the above information (except pMassless
,
for reasons of space).
There is also one method that returns information accumulated for all
the events analyzed so far.
int CellJet::nError()
tells the number of times analyze(...)
failed to analyze
events, i.e. returned false
.
SlowJet
SlowJet
is a simple program for doing jet finding according
to either of the kT, anti-kT, and Cambridge/Aachen
algorithms, in a cylindrical coordinate frame. The name is obviously
an homage to the FastJet
program [Cac06, Cac12].
That package contains many more algorithms, with many more options,
and, above all, is much faster. Therefore SlowJet
is not so much intended for massive processing of data or Monte Carlo
files as for simple first tests. Nevertheless, within the advertised
capabilities of SlowJet
, it has been checked to find
identically the same jets as FastJet
. The time consumption
typically is around or below that to generate an LHC pp event
in the first place, so is not prohibitive. But the time rises rapidly
for large multiplicities, so obviously SlowJet
can not
be used for tricks like distributing a dense grid of pseudoparticles
to be able to define jet areas, like FastJet
can, and also
not for events with much pileup or other noise.
The recent introduction of fjcore
, containing the core
functionality of FastJet
in a very much smaller package,
has changed the conditions. It now is possible (even encouraged by the
authors) to distribute the two fjcore
files as part of the
PYTHIA package. Therefore the SlowJet
class doubles as a
convenient front end to fjcore
, managing the conversion
back and forth between PYTHIA and FastJet
variables. Some
applications may still benefit from using the native codes, but the default
now is to let SlowJet
call on fjcore
for the
jet finding.
The first step is to decide which particles should be included in the
analysis, and with what four-momenta. The SlowJet
constructor
allows to pick a maximum pseudorapidity defined by the extent of the
assumed detector, to pick among some standard options of which particles
to analyze, and to allow for some standard mass assumptions, like that
all charged particles have the pion mass. Obviously this is only a
restricted set of possibilities.
Full flexibility can be obtained by deriving from the base class
SlowJetHook
to write your own include
method.
This will be presented with one final-state particle at a time, and
should return true
for those particles that should be
analyzed. It is also possible to return modified four-momenta and masses,
to take into account detector smearing effects or particle identity
misassignments, but you must respect E^2 - p^2 = m^2.
Alternatively you can modify the event record itself, or a copy of it
(if you want to keep the original intact). For instance, only final
particles are considered in the analysis, i.e. particles with positive
status code, so negative status code should then be assigned to those
particles that you do not want to see analyzed. Again four-momenta and
masses can be modified, subject to E^2 - p^2 = m^2.
The jet reconstructions is then based on sequential recombination with
progressive removal, using the E recombination scheme. To be
more specific, the algorithm works as follows.
- Each particle to be analyzed defines an original cluster. It has a
well-defined four-momentum and mass at input. From this information the
triplet (pT, y, phi) is calculated, i.e. the transverse momentum,
rapidity and azimuthal angle of the cluster.
- Define distance measures of all clusters i to the beam
d_iB = pT_i^2p
and of all pairs (i,j) relative to each other
d_ij = min( pT_i^2p, pT_j^2p) DeltaR_ij^2 / R^2
where
DeltaR_ij^2 = (y_i - y_j)^2 + (phi_i - phi_j)^2.
The jet algorithm is determined by the user-selected p value,
where p = -1 corresponds to the anti-kT one,
p = 0 to the Cambridge/Aachen one and p = 1 to the
kT one. Also R is chosen by the user, to give an
approximate measure of the size of jets. However, note that jets need
not have a circular shape in (y, phi) space, so R
can not straight off be interpreted as a jet radius.
- Find the smallest of all d_iB and d_ij.
- If this is a d_iB then cluster i is removed from
the clusters list and instead added to the jets list.
Optionally, a pTjetMin requirement is imposed, where only
clusters with pT > pTjetMin are added to the jets list.
If so, some of the analyzed particles will not be part of any final
jet.
- If instead the smallest measure is a d_ij then the
four-momenta of the i and j clusters are added
to define a single new cluster. Convert this four-momentum to a new
(pT, y, phi) triplet and update the list of d_iB
and d_ij.
- Return to step 3 until no clusters remain.
To do jet finding analyses you first have to set up a SlowJet
instance, where the arguments of the constructor specifies the details
of the subsequent analyses. Thereafter you can feed in events to it,
one at a time, and have them analyzed by the analyze
method.
Information on the resulting jets can be extracted by a few different methods.
The minimal procedure only requires one call per event to do the analysis.
We will begin by presenting it, and only afterwards some extensions.
SlowJet::SlowJet(int power, double R, double pTjetMin = 0., double etaMax = 25., int select = 2, int massSet = 2, SlowJetHook* sjHookPtr = 0, bool useFJcore = true, bool useStandardR = true)
create a SlowJet
instance, where
argument
power :
tells (half) the power of the transverse-momentum dependence of the
distance measure,
argumentoption
-1 : the anti-kT algorithm,
argumentoption
0 : the Cambridge/Aachen algorithm, and
argumentoption
1 : the kT algorithm.
argument
R :
the R size parameter, which is crudely related to the radius of
the jet cone in (y, phi) space around the center of the jet.
argument
pTjetMin (default = 0.0 GeV
) :
the minimum transverse momentum required for a cluster
to become a jet. By default all clusters become jets, and therefore
all analyzed particles are assigned to a jet.
For comparisons with perturbative QCD, however, it is only meaningful
to consider jets with a significant pT.
argument
etaMax (default = 25.
) :
the maximum +-pseudorapidity that the detector is assumed to cover.
If you pick a value above 20 there is assumed to be full coverage
(obviously only meaningful for theoretical studies).
argument
select (default = 2
) :
tells which particles are analyzed,
argumentoption
1 : all final-state particles,
argumentoption
2 : all observable final-state particles,
i.e. excluding neutrinos and other particles without strong or
electromagnetic interactions (the isVisible()
particle
method),
and
argumentoption
3 : only charged final-state particles.
argument
massSet (default = 2
) : masses assumed for the particles
used in the analysis
argumentoption
0 : all massless,
argumentoption
1 : photons are massless while all others are
assigned the pi+- mass, and
argumentoption
2 : all given their correct masses.
argument
sjHookPtr (default = 0
) :
gives the possibility to send in your own selection routine for which
particles should be part of the analysis; see further below on the
SlowJetHook
class. If this pointer is sent in nonzero,
etaMax
and massSet
are disregarded,
and select
only gives the basic selection, to which
the user can add further requirements.
argument
useFJcore (default = on
) : choice of code used for finding
the jets. Does not affect the outcome of the analysis, but only the speed,
and some more specialized options.
argumentoption
on : use the fjcore
package of
FastJet 3.0.5
.
argumentoption
off : use the native SlowJet
implementation,
which gives a slower jet finding, but allows some extra options of
step-by-step jet joining.
argument
useStandardR (default = on
) : definition of R
distance between two jets. This switch is only meaningful for
useFJcore = false
; within the fjcore
package
the standard option below is always used.
argumentoption
on : standard, as described above,
DeltaR_ij^2 = (y_i - y_j)^2 + (phi_i - phi_j)^2.
argumentoption
off : alternative,
DeltaR_ij^2 = 2 (cosh(y_i - y_j) - cos(phi_i - phi_j)),
which corresponds to the rim of the "deformed cone" giving a constant
invariant mass between the two partons considered (for fixed
masses and transverse momenta).
bool SlowJet::analyze( const Event& event)
performs a jet finding analysis, where
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
If the routine returns false
the analysis failed,
but currently this is not foreseen ever to happen.
After the analysis has been performed, a few SlowJet
class methods are available to return the result of the analysis:
int SlowJet::sizeOrig()
gives the original number of particles (and thus clusters) that the
analysis starts with.
int SlowJet::sizeJet()
gives the number of jets found, with jets numbered 0 through
sizeJet() - 1
, and ordered in terms of decreasing
transverse momentum values w.r.t. the beam axis,
double SlowJet::pT(i)
gives the transverse momentum pT of the i'th jet,
double SlowJet::y(int i)
double SlowJet::phi(int i)
gives the rapidity y and azimuthal angle phi
of the center of the i'th jet (defined by the vector sum
of constituent four-momenta),
Vec4 SlowJet::p(int i)
double SlowJet::m(int i)
gives a Vec4
corresponding to the four-momentum
sum of the particles assigned to the i'th jet, and
the invariant mass of this four-vector,
int SlowJet::multiplicity(int i)
gives the number of particles clustered into the i'th jet,
vector<int> SlowJet::constituents(int i)
gives a list of the indices of the particles that have been
clustered into the i'th jet,
vector<int> SlowJet::clusConstituents(int i)
gives a list of the indices of the particles that have been
clustered into the i'th cluster, at the current stage of
the clustering process,
int SlowJet::jetAssignment(int i)
gives the index of the jet that the particle i of the event
record belongs to, or -1 if there is no jet containing particle
i,
void SlowJet::removeJet(int i)
removes the i'th jet,
void SlowJet::list()
provides a listing of the basic jet information from above.
These are the basic methods. For more sophisticated usage
it is possible to trace the clustering, one step at a time.
It requires the native jet finding code, useFJcore = false
in the constructor. If so, the setup
method should be used
to read in the event and find the initial smallest distance. Each subsequent
doStep
will then do one cluster joining and find the new
smallest distance. You can therefore interrogate which clusters will be
joined next before the joining actually is performed. Alternatively you
can take several steps in one go, or take steps down to a predetermined
number of jets plus remaining clusters.
bool SlowJet::setup( const Event& event)
selects the particles to be analyzed, calculates initial distances,
and finds the initial smallest distance.
argument
event : is an object of the Event
class,
most likely the pythia.event
one.
If the routine returns false
the setup failed,
but currently this is not foreseen ever to happen.
bool SlowJet::doStep()
do the next step of the clustering. This can either be that two
clusters are joined to one, or that a cluster is promoted to a jet
(which is discarded if its pT value is below
pTjetMin
).
The routine will only return false
if there are no
clusters left, or if useFJcore = true
.
bool SlowJet::doNSteps(int nStep)
calls the doStep()
method nStep
times,
if possible. Will return false
if the list of clusters
is emptied before then. The stored jet information is still perfectly
fine; it is only the number of steps that is wrong. Will also return
false
if useFJcore = true
.
bool SlowJet::stopAtN(int nStop)
calls the doStep()
method until a total of nStop
jet and cluster objects remain. Will return false
if this
is not possible, specifically if the number of objects already is smaller
than nStop
when the method is called. The stored jet and
cluster information is still perfectly fine; it only does not have the
expected multiplicity. Will also return false
if
useFJcore = true
.
int SlowJet::sizeAll()
gives the total current number of jets and clusters. The jets are
numbered 0 through sizeJet() - 1
, while the clusters
are numbered sizeJet()
through sizeAll() - 1
.
(Internally jets and clusters are represented by two separate arrays,
but are here presented in one flat range.) Note that the jets are ordered
in decreasing pT, while the clusters are not ordered.
When useFJcore = true
there are no intermediate steps, and
thus the number of clusters is zero (after jet finding).
With this extension, the methods double pT(int i)
,
double y(int i)
, double phi(int i)
,
Vec4 p(int i)
, double m(int i)
and
int multiplicity(int i)
can be used as before.
Furthermore, list()
generalizes
void SlowJet::list(bool listAll = false)
provides a listing of the above information.
argument
listAll : lists both jets and clusters if true
,
else only jets.
Three further methods can be used to check what will happen next.
int SlowJet::iNext()
int SlowJet::jNext()
double SlowJet::dNext()
if the next step is to join two clusters, then the methods give
the (i,j, d_ij) values, if instead to promote
a cluster to a jet then (i, -1, d_iB).
If no clusters remain then (-1, -1, 0.). Note that
the cluster numbers are offset as described above, i.e. they begin at
sizeJet()
, which of course easily could be subtracted off.
Also note that the jet and cluster lists are (moderately) reshuffled
in each new step. When useFJcore = true
there are no
intermediate steps, and thus these methods do not return meaningul
information.
Finally, and separately, the SlowJetHook
class can be used
for a more smart selection of which particles to include in the analysis.
For instance, isolated and/or high-pT muons, electrons and
photons should presumably be identified separately at an early stage,
and then not clustered to jets.
Technically, it works like with User Hooks.
That is, PYTHIA contains the base class. You write a derived class.
In the main program you create an instance of this class, and hand it
in to SlowJet
; in this case already as part of the
constructor.
The following methods should be defined in your derived class.
SlowJetHook::SlowJetHook()
virtual SlowJetHook::~SlowJetHook()
the constructor and destructor need not do anything, and if so you
need not write your own destructor.
virtual bool SlowJetHook::include(int iSel, const Event& event, Vec4& pSel, double& mSel)
is the main method that you will need to write. It will be called
once for each final-state particle in an event, subject to the
value of the select
switch in the SlowJet
constructor. The value select = 2
may be convenient
since then you do not need to remove e.g. neutrinos yourself, but
use select = 1
for full control. The method should then
return true
if you want to see particle included in the
jet clustering, and false
if not.
argument
iSel : is the index in the event record of the
currently studied particle.
argument
event : is an object of the Event
class,
most likely the pythia.event
one, where the currently
studied particle is found.
argument
pSel : is at input the four-momentum of the
currently studied particle. You can change the values, e.g. to take
into account energy smearing in the detector, to define the initial
cluster value, without corrupting the event record itself.
argument
mSel : is at input the mass of the currently studied
particle. You can change the value, e.g. to take into account
particle misidentification, to define the initial cluster value,
without corrupting the event record itself. Note that the changes of
pSel
and mSel
must be coordinated such that
E^2 - p^2 = m^2 holds.
It is also possible to define further methods of your own.
One such could e.g. be called directly in the main program before the
analyze
method is called, to identify and bookkeep
some event properties you may not want to reanalyze for each
individual particle.