Jet Reconstruction Internal Documentation
Documentation for JetReconstruction.jl
's internal methods and types.
N.B. no guarantee is made of stability of these interfaces or types.
Index
JetReconstruction.FourMomentum
JetReconstruction.HistoryElement
JetReconstruction.HistoryElement
JetReconstruction.JetWithAncestors
JetReconstruction.Surrounding
JetReconstruction.TiledJet
JetReconstruction.TiledJet
JetReconstruction.Tiling
JetReconstruction.Tiling
JetReconstruction.TilingDef
JetReconstruction.neighbour_tiles
JetReconstruction.rightmost_tiles
Base.:+
Base.copy
Base.iterate
Base.iterate
Base.iterate
Base.show
Base.tryparse
JetReconstruction.CosTheta
JetReconstruction._ensure_valid_rap_phi
JetReconstruction._plain_jet_reconstruct
JetReconstruction._set_rap_phi!
JetReconstruction._tj_diJ
JetReconstruction._tj_dist
JetReconstruction.add_step_to_history!
JetReconstruction.add_untagged_neighbours_to_tile_union
JetReconstruction.detach!
JetReconstruction.determine_rapidity_extent
JetReconstruction.dij
JetReconstruction.dist
JetReconstruction.do_iB_recombination_step!
JetReconstruction.do_ij_recombination_step!
JetReconstruction.energy
JetReconstruction.eta
JetReconstruction.fast_findmin
JetReconstruction.find_tile_neighbours!
JetReconstruction.geometric_distance
JetReconstruction.get_all_ancestors
JetReconstruction.get_dij_dist
JetReconstruction.get_tile
JetReconstruction.get_tile_cartesian_indices
JetReconstruction.initial_history
JetReconstruction.insert!
JetReconstruction.isvalid
JetReconstruction.jet_ranks
JetReconstruction.m
JetReconstruction.m2
JetReconstruction.mag
JetReconstruction.mass
JetReconstruction.mass2
JetReconstruction.merge_steps
JetReconstruction.phi
JetReconstruction.phi_02pi
JetReconstruction.pt
JetReconstruction.pt2
JetReconstruction.px
JetReconstruction.py
JetReconstruction.pz
JetReconstruction.rapidity
JetReconstruction.reco_state
JetReconstruction.rightneighbours
JetReconstruction.set_momentum!
JetReconstruction.set_nearest_neighbours!
JetReconstruction.setup_tiling
JetReconstruction.surrounding
JetReconstruction.tile_index
JetReconstruction.tiledjet_remove_from_tiles!
JetReconstruction.tiledjet_set_jetinfo!
JetReconstruction.upd_nn_crosscheck!
JetReconstruction.upd_nn_nocross!
JetReconstruction.upd_nn_step!
JetReconstruction.η
Public Interface
+(j1::PseudoJet, j2::PseudoJet)
Addition operator for PseudoJet
objects.
Arguments
j1::PseudoJet
: The firstPseudoJet
object.j2::PseudoJet
: The secondPseudoJet
object.
Returns
A new PseudoJet
object with the sum of the momenta and energy of j1
and j2
.
copy(j::TiledJet)
Create a copy of a TiledJet
object.
Arguments
j::TiledJet
: TheTiledJet
object to be copied.
Returns
A new TiledJet
object with the same attributes as the input object.
Base.iterate(t::neighbour_tiles, state=1)
Iterate over the neighbour_tiles
object, returning all the neighbour tiles for a given Cartesian tile index.
Base.iterate(t::rightmost_tiles, state=1)
Iterate over the rightmost_tiles
object, returning all the rightmost tiles for a given Cartesian tile index.
Base.iterate(tj::TiledJet)
Iterate over a TiledJet
object's linked list, walking over all jets until the end (then the next jet is invalid).
Arguments
tj::TiledJet
: TheTiledJet
object to start to iterate over.
show(io::IO, jet::PseudoJet)
Print a PseudoJet
object to the specified IO stream.
Arguments
io::IO
: The IO stream to which the information will be printed.jet::PseudoJet
: ThePseudoJet
object whose information will be printed.
Base.tryparse(E::Type{<:Enum}, str::String)
Parser that converts a string to an enum value if it exists, otherwise returns nothing.
CosTheta(p::PseudoJet)
Compute the cosine of the angle between the momentum vector p
and the z-axis.
Arguments
p::PseudoJet
: The PseudoJet object representing the momentum vector.
Returns
- The cosine of the angle between
p
and the z-axis.
_ensure_valid_rap_phi(p::PseudoJet)
Ensure that the rapidity and azimuthal angle of the PseudoJet p
are valid. If the azimuthal angle is invalid (used as a proxy for both variables), they are set to a valid value using _set_rap_phi!
.
Arguments
p::PseudoJet
: The PseudoJet object to ensure valid rapidity and azimuthal angle for.
_plain_jet_reconstruct(; particles::Vector{PseudoJet}, p = -1, R = 1.0, recombine = +)
This is the internal implementation of jet reconstruction using the plain algorithm. It takes a vector of particles
representing the input particles and reconstructs jets based on the specified parameters. Here the particles must be of type PseudoJet
.
Users of the package should use the plain_jet_reconstruct
function as their entry point to this jet reconstruction.
The power value maps to specific pp jet reconstruction algorithms: -1 = AntiKt, 0 = Cambridge/Aachen, 1 = Inclusive Kt.
Arguments
particles
: A vector ofPseudoJet
objects representing the input particles.p=-1
: The power to which the transverse momentum (pt
) of each particle is raised.R=1.0
: The jet radius parameter.recombine
: The recombination function used to merge two jets. Default is+
(additive recombination).
Returns
clusterseq
: The resultingClusterSequence
object representing the reconstructed jets.
_set_rap_phi!(p::PseudoJet)
Set the rapidity and azimuthal angle of the PseudoJet p
.
Arguments
p::PseudoJet
: The PseudoJet object for which to set the rapidity and azimuthal angle.
Description
This function calculates and sets the rapidity and azimuthal angle of the PseudoJet p
based on its momentum components. The rapidity is calculated in a way that is insensitive to roundoff errors when the momentum components are large. If the PseudoJet represents a point with infinite rapidity, a large number is assigned to the rapidity in order to lift the degeneracy between different zero-pt momenta.
Note - the ϕ angle is calculated in the range [0, 2π).
_tj_diJ(jet)
Compute the dij metric value for a given jet.
Arguments
jet
: The input jet.
Returns
- The dij value for the jet.
Example
_tj_dist(jetA, jetB)
Compute the geometric distance in the (y, ϕ)-plane between two jets in the TiledAlgoLL module.
Arguments
jetA
: The first jet.jetB
: The second jet.
Returns
The squared distance between jetA
and jetB
.
Examples
add_step_to_history!(clusterseq::ClusterSequence, parent1, parent2, jetp_index, dij)
Add a new jet's history into the recombination sequence.
Arguments:
clusterseq::ClusterSequence
: The cluster sequence object.parent1
: The index of the first parent.parent2
: The index of the second parent.jetp_index
: The index of the jet.dij
: The dij value.
This function adds a new HistoryElement
to the history
vector of the clusterseq
object. The HistoryElement
contains information about the parents, child, jet index, dij value, and the maximum dij value so far. It also updates the child index of the parent elements.
If the parent1
or parent2
have already been recombined, an InternalError
is thrown. The jetp_index
is used to update the _cluster_hist_index
of the corresponding PseudoJet
object.
add_untagged_neighbours_to_tile_union(center_index, tile_union, n_near_tiles, tiling)
Adds to the vector tile_union the tiles that are in the neighbourhood of the specified center_index, including itself and whose tagged status are false - start adding from position n_near_tiles-1, and increase n_near_tiles. When a neighbour is added its tagged status is set to true.
Arguments
center_index
: The index of the center tile.tile_union
: An array to store the indices of neighbouring tiles.n_near_tiles
: The number of neighbouring tiles.tiling
: The tiling object containing the tile tags.
Returns
The updated number of near tiles.
detach!(jet::TiledJet)
Detach a TiledJet
from its linked list by updating the previous
and next
pointers.
Arguments
jet::TiledJet
: TheTiledJet
object to detach.
determine_rapidity_extent(eta::Vector{T}) where T <: AbstractFloat
Calculate the minimum and maximum rapidities based on the input vector eta
. The function determines the rapidity extent by binning the multiplicities as a function of rapidity and finding the minimum and maximum rapidities such that the edge bins contain a certain fraction (~1/4) of the busiest bin and a minimum number of particles.
This is the heuristic which is used by FastJet (inline comments are from FastJet).
Arguments
eta::Vector{T}
: A vector of rapidity values.
Returns
minrap::T
: The minimum rapidity value.maxrap::T
: The maximum rapidity value.
dij(i, kt2_array, nn, nndist)
Compute the dij value for a given index i
to its nearest neighbor. The nearest neighbor is determined from nn[i]
, and the metric distance to the nearest neighbor is given by the distance nndist[i]
applying the lower of the kt2_array
values for the two particles.ßß
Arguments
i
: The index of the element.kt2_array
: An array of kt2 values.nn
: An array of nearest neighbors.nndist
: An array of nearest neighbor distances.
Returns
- The computed dij value.
dist(i, j, rapidity_array, phi_array)
Compute the distance between points in a 2D space defined by rapidity and phi coordinates.
Arguments
i::Int
: Index of the first point to consider (indexes intorapidity_array
andphi_array
).j::Int
: Index of the second point to consider (indexes intorapidity_array
andphi_array
).rapidity_array::Vector{Float64}
: Array of rapidity coordinates.phi_array::Vector{Float64}
: Array of phi coordinates.
Returns
distance::Float64
: The distance between the two points.
do_iB_recombination_step!(clusterseq::ClusterSequence, jet_i, diB)
Bookkeeping for recombining a jet with the beam (i.e., finalising the jet) by adding a step to the history of the cluster sequence.
Arguments
clusterseq::ClusterSequence
: The cluster sequence object.jet_i
: The index of the jet.diB
: The diB value.
do_ij_recombination_step!(clusterseq::ClusterSequence, jet_i, jet_j, dij, recombine=+)
Perform the bookkeeping associated with the step of recombining jet_i and jet_j (assuming a distance dij).
Arguments
clusterseq::ClusterSequence
: The cluster sequence object.jet_i
: The index of the first jet to be recombined.jet_j
: The index of the second jet to be recombined.dij
: The distance between the two jets.recombine=+
: The recombination function to be used. Default is addition.
Returns
newjet_k
: The index of the newly created jet.
Description
This function performs the i-j recombination step in the cluster sequence. It creates a new jet by recombining the first two jets using the specified recombination function. The new jet is then added to the cluster sequence. The function also updates the indices and history information of the new jet and sorts out the history.
energy(p::PseudoJet)
Return the energy of a PseudoJet
.
Arguments
p::PseudoJet
: ThePseudoJet
object.
Returns
- The energy of the
PseudoJet
.
eta(p::PseudoJet)
Compute the pseudorapidity (η) of a PseudoJet.
Arguments
p::PseudoJet
: The PseudoJet object for which to compute the pseudorapidity.
Returns
- The pseudorapidity (η) of the PseudoJet.
fast_findmin(dij, n)
Find the minimum value and its index in the first n
elements of the dij
array. The use of @turbo
macro gives a significiant performance boost.
Arguments
dij
: An array of values.n
: The number of elements to consider in thedij
array.
Returns
dij_min
: The minimum value in the firstn
elements of thedij
array.best
: The index of the minimum value in thedij
array.
find_tile_neighbours!(tile_union, jetA, jetB, oldB, tiling)
Find the union of neighbouring tiles of jetA
, jetB
, and oldB
and add them to the tile_union
. This established the set of tiles over which searches for updated and new nearest-neighbours must be run
Arguments
tile_union
: The tile union to which the neighbouring tiles will be added.jetA
: The first jet.jetB
: The second jet.oldB
: The old second jet.tiling
: The tiling information.
Returns
The number of neighbouring tiles added to the tile_union
.
geometric_distance(eta1::AbstractFloat, phi1::AbstractFloat, eta2::AbstractFloat, phi2::AbstractFloat)
Compute the geometric distance between two points in the rap-phi plane.
Arguments
eta1::AbstractFloat
: The eta coordinate of the first point.phi1::AbstractFloat
: The phi coordinate of the first point.eta2::AbstractFloat
: The eta coordinate of the second point.phi2::AbstractFloat
: The phi coordinate of the second point.
Returns
distance::Float64
: The geometric distance between the two points.
get_all_ancestors(idx, cs::ClusterSequence)
Recursively finds all ancestors of a given index in a ClusterSequence
object.
Arguments
idx
: The index of the jet for which to find ancestors.cs
: TheClusterSequence
object containing the jet history.
Returns
An array of indices representing the ancestors of the given jet.
get_dij_dist(nn_dist, kt2_1, kt2_2, R2)
Compute the dij metric distance between two jets.
Arguments
nn_dist
: The nearest-neighbor distance between two jets.kt2_1
: The squared momentum metric value of the first jet.kt2_2
: The squared momentum metric value of the second jet.R2
: The jet radius parameter squared.
Returns
The distance between the two jets.
If kt2_2
is equal to 0.0, then the first jet doesn't actually have a valid neighbour, so it's treated as a single jet adjecent to the beam.
get_tile(tiling_setup::TilingDef, eta::AbstractFloat, phi::AbstractFloat)
Given a tiling_setup
object, eta
and phi
values, this function calculates the tile indices for the given eta
and phi
values.
Arguments
tiling_setup
: ATilingDef
object that contains the tiling setup parameters.eta
: The eta value for which to calculate the tile index.phi
: The phi value for which to calculate the tile index.
Returns
ieta
: The tile index along the eta direction.iphi
: The tile index along the phi direction.
get_tile_linear_index(tiling_setup::TilingDef, i_η::Int, i_ϕ::Int)
Compute the linear index of a tile in a tiled setup. This is much faster in this function than using the LinearIndices construct (like x100, which is bonkers, but there you go...)
Arguments
tiling_setup::TilingDef
: The tiling setup defining the number of tiles in each dimension.i_η::Int
: The index of the tile in the η dimension.i_ϕ::Int
: The index of the tile in the ϕ dimension.
Returns
- The linear index of the tile.
initial_history(particles)
Create an initial history for the given particles.
Arguments
particles
: The initial vector of stable particles.
Returns
history
: An array ofHistoryElement
objects.Qtot
: The total energy in the event.
insert!(nextjet::TiledJet, jettomove::TiledJet)
Inserts a TiledJet
object into the linked list of TiledJet
objects, before the nextjet
object. The jet to move can be an isolated jet, a jet from another list or a jet from the same list
Arguments
nextjet::TiledJet
: TheTiledJet
object after whichjettomove
should be inserted.jettomove::TiledJet
: TheTiledJet
object to be inserted.
Example
isvalid(t::TiledJet)
Check if a TiledJet
is valid, by seeing if it is not the noTiledJet
object.
Arguments
t::TiledJet
: TheTiledJet
object to check.
Returns
Bool
:true
if theTiledJet
object is valid,false
otherwise.
jet_ranks(clusterseq::ClusterSequence; compare_fn = JetReconstruction.pt)
Compute the ranks of jets in a given ClusterSequence
object based on a specified comparison function.
Arguments
clusterseq::ClusterSequence
: TheClusterSequence
object containing the jets to rank.compare_fn = JetReconstruction.pt
: The comparison function used to determine the order of the jets. Defaults toJetReconstruction.pt
, which compares jets based on their transverse momentum.
Returns
A dictionary mapping each jet index to its rank.
Note
This is a utility function that can be used to rank initial clusters based on a specified jet property. It can be used to assign a consistent "rank" to each reconstructed jet in the cluster sequence, which is useful for stable plotting of jet outputs.
m(p::PseudoJet)
Compute the invariant mass of a PseudoJet
object. By convention if m^2 < 0, then -sqrt{(-m^2)} is returned.
Arguments
p::PseudoJet
: ThePseudoJet
object for which to compute the invariant mass.
Returns
The invariant mass of the PseudoJet
object.
m2(p::PseudoJet)
Calculate the invariant mass squared (m^2) of a PseudoJet.
Arguments
p::PseudoJet
: The PseudoJet object for which to calculate the invariant mass squared.
Returns
- The invariant mass squared (m^2) of the PseudoJet.
mag(p::PseudoJet)
Return the magnitude of the momentum of a PseudoJet
, |p|
.
Arguments
p::PseudoJet
: ThePseudoJet
object for which to compute the magnitude.
Returns
The magnitude of the PseudoJet
object.
mass(p::PseudoJet)
Compute the invariant mass (alias for m(p)
).
Arguments
p::PseudoJet
: The PseudoJet object for which to compute the mass.
Returns
- The mass of the PseudoJet.
merge_steps(clusterseq::ClusterSequence)
Compute the number of jet-jet merge steps in a cluster sequence. This is useful to give the number of meaningful recombination steps in a jet reconstruction sequence (beam merge steps are not counted).
Arguments
clusterseq::ClusterSequence
: The cluster sequence object.
Returns
merge_steps::Int
: The number of merge steps.
phi(p::PseudoJet)
Compute the ϕ angle of a PseudoJet
object p
.
Note this function is a wrapper for phi_02pi(p)
.
Arguments
p::PseudoJet
: ThePseudoJet
object for which to compute the azimuthal angle.
Returns
- The azimuthal angle of
p
in the range [0, 2π).
phi_02pi(p::PseudoJet)
Compute the azimuthal angle of a PseudoJet
object p
in the range [0, 2π).
Arguments
p::PseudoJet
: ThePseudoJet
object for which to compute the azimuthal angle.
Returns
- The azimuthal angle of
p
in the range [0, 2π).
pt(p::PseudoJet)
Compute the scalar transverse momentum (pt) of a PseudoJet.
Arguments
p::PseudoJet
: The PseudoJet object for which to compute the transverse momentum.
Returns
- The transverse momentum (pt) of the PseudoJet.
pt2(p::PseudoJet)
Get the squared transverse momentum of a PseudoJet.
Arguments
p::PseudoJet
: The PseudoJet object.
Returns
- The squared transverse momentum of the PseudoJet.
px(p::PseudoJet)
Return the x-component of the momentum of a PseudoJet
.
Arguments
p::PseudoJet
: ThePseudoJet
object.
Returns
- The x-component of the momentum of the
PseudoJet
.
py(p::PseudoJet)
Return the y-component of the momentum of a PseudoJet
.
Arguments
p::PseudoJet
: ThePseudoJet
object.
Returns
- The y-component of the momentum of the
PseudoJet
.
pz(p::PseudoJet)
Return the z-component of the momentum of a PseudoJet
.
Arguments
p::PseudoJet
: ThePseudoJet
object.
Returns
- The z-component of the momentum of the
PseudoJet
.
rapidity(p::PseudoJet)
Compute the rapidity of a PseudoJet
object.
Arguments
p::PseudoJet
: ThePseudoJet
object for which to compute the rapidity.
Returns
The rapidity of the PseudoJet
object.
reco_state(cs::ClusterSequence, pt_ranks; iteration=0)
This function returns the reconstruction state of a ClusterSequence
object based on a given iteration number in the reconstruction.
Arguments
cs::ClusterSequence
: TheClusterSequence
object to update.ranks
: The ranks of the original clusters, that are inherited by peudojets
during the reconstruction process.
iteration=0
: The iteration number to consider for updating the reconstruction state (0 represents the initial state).ignore_beam_merge=true
: Ignore beam merging steps in the reconstruction (which produce no change in status).
Returns
A dictionary representing a snapshot of the reconstruction state.
Details
The function starts by initializing the reconstruction state with the initial particles. Then, it walks over the iteration sequence and updates the reconstruction state based on the history of recombination and finalization/beam merger steps.
rightneighbours(center::Int, tiling::Tiling)
Compute the indices of the right neighbors of a given center index in a tiling. This is used in the inital sweep to calculate the nearest neighbors, where the search between jets for the nearest neighbour is bi-directional, thus when a tile is considered only the right neighbours are needed to compare jet distances as the left-hand tiles have been done from that tile already.
Arguments
center::Int
: The center index.tiling::Tiling
: The tiling object.
Returns
Surrounding
: An object containing the indices of the right neighbors.
set_momentum!(j::PseudoJet, px, py, pz, E)
Set the momentum components and energy of a PseudoJet
object.
Arguments
j::PseudoJet
: ThePseudoJet
object to set the momentum for.px
: The x-component of the momentum.py
: The y-component of the momentum.pz
: The z-component of the momentum.E
: The energy of the particle.
set_nearest_neighbours!(clusterseq::ClusterSequence, tiling::Tiling, tiledjets::Vector{TiledJet})
This function sets the nearest neighbor information for all jets in the tiledjets
vector.
Arguments
clusterseq::ClusterSequence
: The cluster sequence object.tiling::Tiling
: The tiling object.tiledjets::Vector{TiledJet}
: The vector of tiled jets.
Returns
NNs::Vector{TiledJet}
: The vector of nearest neighbor jets.diJ::Vector{Float64}
: The vector of diJ values.
The function iterates over each tile in the tiling
and sets the nearest neighbor information for each jet in the tile. It then looks for neighbor jets in the neighboring tiles and updates the nearest neighbor information accordingly. Finally, it creates the diJ table and returns the vectors of nearest neighbor jets and diJ values.
Note: The diJ values are calculated as the kt distance multiplied by R^2.
setup_tiling(eta::Vector{T}, Rparam::AbstractFloat) where T <: AbstractFloat
This function sets up the tiling parameters for a reconstruction given a vector of rapidities eta
and a radius parameter Rparam
.
Arguments
eta::Vector{T}
: A vector of rapidities.Rparam::AbstractFloat
: The jet radius parameter.
Returns
tiling_setup
: ATilingDef
object containing the tiling setup parameters.
Description
The function first decides the tile sizes based on the Rparam
value. It then determines the number of tiles in the phi direction (n_tiles_phi
) based on the tile size. Next, it determines the rapidity extent of the input eta
vector and adjusts the values accordingly. Finally, it creates a TilingDef
object with the calculated tiling parameters and returns it.
surrounding(center::Int, tiling::Tiling)
Compute the surrounding indices of a given center index in a tiling.
Arguments
center::Int
: The center index.tiling::Tiling
: The tiling object.
Returns
Surrounding
: An object containing the surrounding indices.
tile_index(tiling_setup, eta::Float64, phi::Float64)
Compute the tile index for a given (eta, phi) coordinate.
Arguments
tiling_setup
: The tiling setup object containing the tile size and number of tiles.eta::Float64
: The eta coordinate.phi::Float64
: The phi coordinate.
Returns
The tile index corresponding to the (eta, phi) coordinate.
tiledjet_remove_from_tiles!(tiling, jet)
Remove a jet from the given tiling structure.
Arguments
tiling
: The tiling structure from which the jet will be removed.jet
: The jet to be removed from the tiling structure.
Description
This function removes a jet from the tiling structure. It adjusts the linked list to be consistent with the removal of the jet.
tiledjet_set_jetinfo!(jet::TiledJet, clusterseq::ClusterSequence, tiling::Tiling, jets_index, R2, p)
Initialise a tiled jet from a PseudoJet (using an index into our ClusterSequence)
Arguments:
jet::TiledJet
: The TiledJet object to set the information for.clusterseq::ClusterSequence
: The ClusterSequence object containing the jets.tiling::Tiling
: The Tiling object containing the tile information.jets_index
: The index of the jet in the ClusterSequence.R2
: The jet radius parameter squared.p
: The power to raise the pt2 value to.
This function sets the eta, phi, kt2, jets_index, NN_dist, NN, tile_index, previous, and next fields of the TiledJet object.
Returns:
nothing
upd_nn_crosscheck!(i, from, to, rapidity_array, phi_array, R2, nndist, nn)
Update the nearest neighbor information for a given particle index i
against all particles in the range indexes from
to to
. The function updates the nndist
and nn
arrays with the nearest neighbor distance and index respectively, both for particle i
and the checked particles [from:to]
(hence crosscheck).
Arguments
i::Int
: The index of the particle to update and check against.from::Int
: The starting index of the range of particles to check against.to::Int
: The ending index of the range of particles to check against.rapidity_array
: An array containing the rapidity values of all particles.phi_array
: An array containing the phi values of the all particles.R2
: The squared jet distance threshold for considering a particle as a neighbour.nndist
: The array that stores the nearest neighbor distances.nn
: The array that stores the nearest neighbor indices.
upd_nn_nocross!(i, from, to, rapidity_array, phi_array, R2, nndist, nn)
Update the nearest neighbor information for a given particle index i
against all particles in the range indexes from
to to
. The function updates the nndist
and nn
arrays with the nearest neighbor distance and index respectively, only for particle i
(hence nocross).
Arguments
i::Int
: The index of the particle to update and check against.from::Int
: The starting index of the range of particles to check against.to::Int
: The ending index of the range of particles to check against.rapidity_array
: An array containing the rapidity values of all particles.phi_array
: An array containing the phi values of the all particles.R2
: The squared jet distance threshold for considering a particle as a neighbour.nndist
: The array that stores the nearest neighbor distances.nn
: The array that stores the nearest neighbor indices.
upd_nn_step!(i, j, k, N, Nn, kt2_array, rapidity_array, phi_array, R2, nndist, nn, nndij)
Update the nearest neighbor information after a jet merge step.
Arguments:
i
: Index of the first particle in the last merge step.j
: Index of the second particle in the last merge step.k
: Index of the current particle for which the nearest neighbour will be updated.N
: Total number of particles (currently vaild array indexes are[1:N]
).Nn
: Number of nearest neighbors to consider.kt2_array
: Array of transverse momentum squared values.rapidity_array
: Array of rapidity values.phi_array
: Array of azimuthal angle values.R2
: Distance threshold squared for nearest neighbors.nndist
: Array of nearest neighbor geometric distances.nn
: Array of nearest neighbor indices.nndij
: Array of metric distances between particles.
This function updates the nearest neighbor information for the current particle k
by considering the distances to particles i
and j
. It checks if the distance between k
and i
is smaller than the current nearest neighbor distance for k
, and updates the nearest neighbor information accordingly. It also updates the nearest neighbor information for i
if the distance between k
and i
is smaller than the current nearest neighbor distance for i
. Finally, it checks if the nearest neighbor of k
is the total number of particles Nn
and updates it to j
if necessary.
Interface for composite types that includes fields px, py, py, and E that represents the components of a four-momentum vector.
struct HistoryElement
A struct holding a record of jet mergers and finalisations
Fields:
parent1
: Index in history where first parent of this jet was created (NonexistentParent if this jet is an original particle)parent2
: Index in history where second parent of this jet was created (NonexistentParent if this jet is an original particle); BeamJet if this history entry just labels the fact that the jet has recombined with the beam)child
: Index in history where the current jet is recombined with another jet to form its child. It is Invalid if this jet does not further recombine.jetp_index
: Index in the jets vector where we will find the PseudoJet object corresponding to this jet (i.e. the jet created at this entry of the history). NB: if this element of the history corresponds to a beam recombination, thenjetp_index=Invalid
.dij
: The distance corresponding to the recombination at this stage of the clustering.max_dij_so_far
: The largest recombination distance seen so far in the clustering history.
HistoryElement(jetp_index)
Constructs a HistoryElement
object with the given jetp_index
, used for initialising the history with original particles.
Arguments
jetp_index
: The index of the jetp.
Returns
A HistoryElement
object.
struct JetWithAncestors
A struct representing a jet with its origin ancestors.
Fields
self::PseudoJet
: The PseudoJet object for this jet.jetp_index::Int
: The index of the jet in the corresponding cluster sequence.ancestors::Set{Int}
: A set of indices representing the jetp_indexes of ancestors of the jet (in the cluster sequence).jet_rank::Int
: The rank of the jet based on a comparison of all of the jet's ancestors
Note
This structure needs its associated cluster sequence origin to be useful.
struct Surrounding{N}
Structure used for iterating over neighbour tiles.
Fields
indices::NTuple{N, Int}
: A tuple ofN
integers representing the indices.
struct TiledJet
TiledJet represents a jet in a tiled algorithm for jet reconstruction, with additional information to track the jet's position in the tiled structures.
Fields
id::Int
: The ID of the jet.eta::Float64
: The rapidity of the jet.phi::Float64
: The azimuthal angle of the jet.kt2::Float64
: The transverse momentum squared of the jet.NN_dist::Float64
: The distance to the nearest neighbor.jets_index::Int
: The index of the jet in the jet array.tile_index::Int
: The index of the tile in the tile array.dij_posn::Int
: The position of this jet in the dij compact array.NN::TiledJet
: The nearest neighbor.previous::TiledJet
: The previous jet.next::TiledJet
: The next jet.
TiledJet(id)
Constructs a TiledJet
object with the given id
and initializes its properties to zero.
Arguments
id
: The ID of theTiledJet
object.
Returns
A TiledJet
object with the specified id
and values set to zero or noTiledJet.
struct Tiling
The Tiling
struct represents a tiling configuration for jet reconstruction.
Fields
setup::TilingDef
: The tiling definition used for the configuration.tiles::Matrix{TiledJet}
: A matrix of tiled jets, containing the first jet in each tile (then the linked list of the first jet is followed to get access to all jets in this tile).positions::Matrix{Int}
: Used to track tiles that are on the edge of ϕ array, where neighbours need to be wrapped around.tags::Matrix{Bool}
: The matrix of tags indicating whether a tile is valid or not (set tofalse
initially, thentrue
when the tile has been setup properly).
Tiling(setup::TilingDef)
Constructs a intial Tiling
object based on the provided setup
parameters.
Arguments
setup::TilingDef
: The setup parameters for the tiling.
Returns
A Tiling
object.
struct TilingDef
A struct representing the definition of a spcific tiling scheme.
Fields
_tiles_eta_min::Float64
: The minimum rapidity of the tiles._tiles_eta_max::Float64
: The maximum rapidity of the tiles._tile_size_eta::Float64
: The size of a tile in rapidity (usually R^2)._tile_size_phi::Float64
: The size of a tile in phi (usually a bit more than R^2)._n_tiles_eta::Int
: The number of tiles across rapidity._n_tiles_phi::Int
: The number of tiles across phi._n_tiles::Int
: The total number of tiles._tiles_ieta_min::Int
: The minimum rapidity tile index._tiles_ieta_max::Int
: The maximum rapidity tile index.
Constructor
TilingDef(_tiles_eta_min, _tiles_eta_max, _tile_size_eta, _tile_size_phi,
_n_tiles_eta, _n_tiles_phi, _tiles_ieta_min, _tiles_ieta_max)
Constructs a TilingDef
object with the given parameters.
struct neighbour_tiles
A struct representing the neighbouring tiles.
A struct for iterating over all neighbour tiles for a given Cartesian tile index. These are the tiles above and to the right of the given tile (X=included, O=not included):
XXX
X.X
XXX
Note, rapidity coordinate must be in range, ϕ coordinate wraps
Fields
n_η::Int
: Number of η tilesn_ϕ::Int
: Number of ϕ tilesstart_η::Int
: Centre η tile coordinatestart_ϕ::Int
: Centre ϕ tile coordinate
struct rightmost_tiles
A struct for iterating over rightmost tiles for a given Cartesian tile index. These are the tiles above and to the right of the given tile (X=included, O=not included):
XXX
O.X
OOO
Note, rapidity coordinate must be in range, ϕ coordinate wraps
Fields
n_η::Int
: Number of η tilesn_ϕ::Int
: Number of ϕ tilesstart_η::Int
: Centre η tile coordinatestart_ϕ::Int
: Centre ϕ tile coordinate