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.FourMomentumJetReconstruction.HistoryElementJetReconstruction.HistoryElementJetReconstruction.JetWithAncestorsJetReconstruction.SurroundingJetReconstruction.TiledJetJetReconstruction.TiledJetJetReconstruction.TilingJetReconstruction.TilingJetReconstruction.TilingDefJetReconstruction.neighbour_tilesJetReconstruction.rightmost_tilesBase.:+Base.copyBase.isapproxBase.isvalidBase.isvalidBase.iterateBase.iterateBase.iterateBase.showBase.showJetReconstruction.CosThetaJetReconstruction.EJetReconstruction._addjets_with_scaleJetReconstruction._ee_genkt_algorithm!JetReconstruction._plain_jet_reconstruct!JetReconstruction._tiled_jet_reconstruct!JetReconstruction._tj_diJJetReconstruction._tj_distJetReconstruction.add_step_to_history!JetReconstruction.add_untagged_neighbours_to_tile_unionJetReconstruction.addjets_eschemeJetReconstruction.angular_distanceJetReconstruction.cluster_hist_indexJetReconstruction.copy_to_slot!JetReconstruction.declusterJetReconstruction.deltaRJetReconstruction.delta_phiJetReconstruction.deltarJetReconstruction.detach!JetReconstruction.determine_rapidity_extentJetReconstruction.dijJetReconstruction.dij_distJetReconstruction.distJetReconstruction.do_iB_recombination_step!JetReconstruction.energyJetReconstruction.etaJetReconstruction.fast_findminJetReconstruction.find_tile_neighbours!JetReconstruction.geometric_distanceJetReconstruction.get_algorithm_powerJetReconstruction.get_all_ancestorsJetReconstruction.get_dij_distJetReconstruction.get_tileJetReconstruction.get_tile_cartesian_indicesJetReconstruction.initial_historyJetReconstruction.insert!JetReconstruction.is_eeJetReconstruction.is_ppJetReconstruction.jet_ranksJetReconstruction.kt2JetReconstruction.mJetReconstruction.m2JetReconstruction.magJetReconstruction.massJetReconstruction.mass2JetReconstruction.merge_stepsJetReconstruction.nxJetReconstruction.nyJetReconstruction.nzJetReconstruction.open_with_streamJetReconstruction.pJetReconstruction.p2JetReconstruction.p2JetReconstruction.phiJetReconstruction.phiJetReconstruction.ptJetReconstruction.ptJetReconstruction.pt2JetReconstruction.pt2JetReconstruction.pxJetReconstruction.pyJetReconstruction.pzJetReconstruction.rapidityJetReconstruction.rapidityJetReconstruction.reclusterJetReconstruction.reco_stateJetReconstruction.rightneighboursJetReconstruction.set_nearest_neighbours!JetReconstruction.setup_tilingJetReconstruction.surroundingJetReconstruction.tile_indexJetReconstruction.tile_indexJetReconstruction.tiledjet_remove_from_tiles!JetReconstruction.tiledjet_set_jetinfo!JetReconstruction.upd_nn_crosscheck!JetReconstruction.upd_nn_nocross!JetReconstruction.upd_nn_step!JetReconstruction.η
Internal Methods and Types
Base.:+ — Method+(jet1::T, jet2::T) where {T <: FourMomentum}Adds two four-momentum vectors together, returning a new jet.
Details
This addition operation will return a jet with the cluster history index set to
- *This means that this jet cannot be used, or be part of, any clustering
history.*
Base.copy — Methodcopy(j::TiledJet)Create a copy of a TiledJet object.
Arguments
j::TiledJet: TheTiledJetobject to be copied.
Returns
A new TiledJet object with the same attributes as the input object.
Base.isapprox — Methodisapprox(fj1::FinalJet, fj2::FinalJet; kwargs...)Compare two FinalJet objects for approximate equality based on their rapidity, azimuthal angle, and transverse momentum. kwargs are passed to the isapprox function, e.g., tol.
Base.isvalid — Methodisvalid(t::TiledJet)Check if a TiledJet is valid, by seeing if it is not the noTiledJet object.
Arguments
t::TiledJet: TheTiledJetobject to check.
Returns
Bool:trueif theTiledJetobject is valid,falseotherwise.
Base.isvalid — Methodisvalid(j::PseudoJet)Function to check whether a given PseudoJet object is non-zero or not. Primarily to use for checking the validity of outputs of substructure tagging.
Returns
Bool:trueif thePseudoJetobject is non-zero (valid),falseotherwise.
Base.iterate — FunctionBase.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 — FunctionBase.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 — MethodBase.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: TheTiledJetobject to start to iterate over.
Base.show — Methodshow(io::IO, jet::FourMomentum)Print core information about the four-momentum vector of jet to the given IO stream.
Base.show — Methodshow(io::IO, sk::SoftKiller)Pretty-print the SoftKiller grid configuration.
JetReconstruction.CosTheta — MethodCosTheta(jet::T) where {T <: FourMomentum}Compute the cosine of the angle between the momentum vector of jet and the z-axis.
Returns
- The cosine of the angle between the jet and the z-axis.
JetReconstruction.E — FunctionAlias for energy function
JetReconstruction._addjets_with_scale — Method_addjets_with_scale(scale1::Real, scale2::Real, jet1::T, jet2::T, cluster_hist_index::Int) where {T <: FourMomentum}Combine two jets as massless objects using the given scale factors for each jet, and return the new jet with the cluster history index set.
JetReconstruction._ee_genkt_algorithm! — Method_ee_genkt_algorithm!(particles::AbstractVector{EEJet};
algorithm::JetAlgorithm.Algorithm, p::Real, R = 4.0,
recombine = addjets)This function is the internal implementation of the e+e- jet clustering algorithm. It takes a vector of EEJet particles representing the input particles and reconstructs jets based on the specified parameters.
Users of the package should use the ee_genkt_algorithm function as their entry point to this jet reconstruction.
Arguments
particles::AbstractVector{EEJet}: A vector ofEEJetparticles used as input for jet reconstruction. This vector must supply the correctcluster_hist_indexvalues and will be mutated as part of the returnedClusterSequence.algorithm::JetAlgorithm.Algorithm: The jet reconstruction algorithm to use.p::Real: The power to which the transverse momentum (pt) of each particle is raised.R = 4.0: The jet radius parameter.recombine = addjets: The recombination function used to merge two jets.
Returns
clusterseq: The resultingClusterSequenceobject representing the reconstructed jets.
JetReconstruction._plain_jet_reconstruct! — Method_plain_jet_reconstruct!(particles::AbstractVector{PseudoJet};
algorithm::JetAlgorithm.Algorithm, p::Real, R = 1.0,
recombine = addjets)This is the internal implementation of jet reconstruction using the plain algorithm. It takes a vector of PseudoJet particles representing the input particles and reconstructs jets based on the specified parameters.
Users of the package should use the plain_jet_reconstruct function as their entry point to this jet reconstruction.
Arguments
particles::AbstractVector{PseudoJet}: A vector ofPseudoJetparticles used as input for jet reconstruction. This vector must supply the correctcluster_hist_indexvalues and will be mutated as part of the returnedClusterSequence.algorithm::JetAlgorithm.Algorithm: The jet reconstruction algorithm to use.p::Real: The power to which the transverse momentum (pt) of each particle is raised.R = 1.0: The jet radius parameter.recombine = addjets: The recombination scheme to use.
Returns
clusterseq: The resultingClusterSequenceobject representing the reconstructed jets.
JetReconstruction._tiled_jet_reconstruct! — Method_tiled_jet_reconstruct!(particles::AbstractVector{PseudoJet};
algorithm::JetAlgorithm.Algorithm,
p::Real, R = 1.0, recombine = addjets)Main jet internal reconstruction algorithm entry point for reconstructing jets once preprocessing of data types are done. The algorithm parameter must be consistent with the power parameter.
Arguments
particles::AbstractVector{PseudoJet}: A vector ofPseudoJetparticles used as input for jet reconstruction. This vector must supply the correctcluster_hist_indexvalues and will be mutated as part of the returnedClusterSequence.algorithm::JetAlgorithm.Algorithm: The jet reconstruction algorithm to use.p::Real: The power parameter for the jet reconstruction algorithm, thus switching between different algorithms.R = 1.0: The jet radius parameter for the jet reconstruction algorithm.recombine::Function = addjets: The recombination function used for combining pseudojets.
Returns
clusterseq: The resultingClusterSequenceobject representing the reconstructed jets.
Example
_tiled_jet_reconstruct!(particles::Vector{PseudoJet}; algorithm = JetAlgorithm.Kt, p = 1, R = 0.4)JetReconstruction._tj_diJ — Method_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
JetReconstruction._tj_dist — Method_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
JetReconstruction.add_step_to_history! — Methodadd_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.
JetReconstruction.add_untagged_neighbours_to_tile_union — Methodadd_untagged_neighbours_to_tile_union(center_index, tile_union, n_near_tiles, tiling)Adds to the vector tileunion the tiles that are in the neighbourhood of the specified centerindex, including itself and whose tagged status are false - start adding from position nneartiles-1, and increase nneartiles. 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.
JetReconstruction.addjets_escheme — Functionconst addjets_escheme = addjetsAlias for the addjets function, which is the default jet recombination scheme.
JetReconstruction.angular_distance — Methodangular_distance(eereco, i, j) -> Float64Calculate the angular distance between two jets i and j using the formula $1 - cos(θ_{ij})$.
Arguments
eereco: The array ofEERecoJetobjects.i: The first jet.j: The second jet.
Returns
Float64: The angular distance betweeniandj, which is $1 - cos heta$.
JetReconstruction.cluster_hist_index — Methodcluster_hist_index(j::FourMomentum)Return the cluster history index of the jet j.
JetReconstruction.copy_to_slot! — Methodcopy_to_slot!(eereco, i, j)Copy the contents of slot i in the eereco array to slot j.
JetReconstruction.decluster — Methoddecluster(jet::T, clusterseq::ClusterSequence{T}) where {T <: FourMomentum}Given a jet and its clustering sequence, this function attempts to decluster it into its two parent subjets. If both parents exist, it returns them ordered by descending pt².
Returns:
(j1, j2)wherej1is the harder subjet.
JetReconstruction.deltaR — MethoddeltaR(jet1::T, jet2::T) where T <: FourMomentumFunction to calculate the distance in the y-ϕ plane between two jets jet1 and jet2 (that is using rapidity and azimuthal angle).
Returns
- The Euclidean distance in the y-ϕ plane for the two jets.
JetReconstruction.delta_phi — Methoddelta_phi(jet1::T, jet2::T) where {T <: FourMomentum}Computes the difference in azimuthal angle φ between two jets, wrapped into the range [-π, π].
Returns
δφas a Float64 in radians.
JetReconstruction.deltar — Methoddeltar(jet1::T, jet2::T) where T <: FourMomentumFunction to calculate the distance in the η-ϕ plane between two jets jet1 and jet2 (that is, using the pseudorapidity and azimuthal angle).
Returns
- The Euclidean distance in the η-ϕ plane for the two jets.
JetReconstruction.detach! — Methoddetach!(jet::TiledJet)Detach a TiledJet from its linked list by updating the previous and next pointers.
Arguments
jet::TiledJet: TheTiledJetobject to detach.
JetReconstruction.determine_rapidity_extent — Methoddetermine_rapidity_extent(eta::Vector{T}) where T <: AbstractFloatCalculate 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.
JetReconstruction.dij — Methoddij(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.
JetReconstruction.dij_dist — Methoddij_dist(eereco, i, j, dij_factor)Calculate the dij distance between two $e^+e^-$jets.
Arguments
eereco: The array ofEERecoJetobjects.i: The first jet.j: The second jet.dij_factor: The scaling factor to multiply the dij distance by.
Returns
- The dij distance between
iandj.
JetReconstruction.dist — Methoddist(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_arrayandphi_array).j::Int: Index of the second point to consider (indexes intorapidity_arrayandphi_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.
JetReconstruction.do_iB_recombination_step! — Methoddo_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.
JetReconstruction.energy — Methodenergy(j::FourMomentum)Return the energy component of the four-momentum vector of j.
JetReconstruction.eta — Methodeta(jet::T) where {T <: FourMomentum}Compute the pseudorapidity (η) of a jet.
Returns
- The pseudorapidity (η) of the jet.
JetReconstruction.fast_findmin — Methodfast_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 significant performance boost.
Arguments
dij: An array of values.n: The number of elements to consider in thedijarray.
Returns
dij_min: The minimum value in the firstnelements of thedijarray.best: The index of the minimum value in thedijarray.
JetReconstruction.find_tile_neighbours! — Methodfind_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.
JetReconstruction.geometric_distance — Methodgeometric_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.
JetReconstruction.get_algorithm_power — Methodget_algorithm_power(; algorithm::JetAlgorithm.Algorithm, p::Union{Real, Nothing}) -> RealGet the algorithm power
This function returns appropriate power value for the specified jet algorithm if the algorithm is a fixed power algorithm (like AntiKt, CA, Kt, or Durham). If the algorithm is generalized (like GenKt or EEKt), it requires a power value to be specified and the function returns the same value.
Arguments
p::Union{Real, Nothing}: The power value.algorithm::JetAlgorithm.Algorithm: The algorithm.
Returns
Resolved algorithm power value.
Throws
ArgumentError: If no power is specified for a generalized algorithm
JetReconstruction.get_all_ancestors — Methodget_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: TheClusterSequenceobject containing the jet history.
Returns
An array of indices representing the ancestors of the given jet.
JetReconstruction.get_dij_dist — Methodget_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 adjacent to the beam.
JetReconstruction.get_tile — Methodget_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: ATilingDefobject 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.
JetReconstruction.get_tile_cartesian_indices — Methodget_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.
JetReconstruction.initial_history — Methodinitial_history(particles)Create an initial history for the given particles.
Arguments
particles: The initial vector of stable particles.
Returns
history: An array ofHistoryElementobjects.Qtot: The total energy in the event.
JetReconstruction.insert! — Methodinsert!(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: TheTiledJetobject after whichjettomoveshould be inserted.jettomove::TiledJet: TheTiledJetobject to be inserted.
Example
JetReconstruction.is_ee — Methodis_ee(algorithm::JetAlgorithm.Algorithm)Check if the algorithm is a e+e- reconstruction algorithm.
Returns
true if the algorithm is a e+e- reconstruction algorithm, false otherwise.
JetReconstruction.is_pp — Methodis_pp(algorithm::JetAlgorithm.Algorithm)Check if the algorithm is a pp reconstruction algorithm.
Returns
true if the algorithm is a pp reconstruction algorithm, false otherwise.
JetReconstruction.jet_ranks — Methodjet_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: TheClusterSequenceobject 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.
JetReconstruction.kt2 — FunctionAlias for pt2 function
JetReconstruction.m — FunctionAlias for mass function
JetReconstruction.m2 — FunctionAlias for mass2 function
JetReconstruction.mag — Methodmag(jet::T) where {T <: FourMomentum}Return the magnitude of the momentum of a jet, |p|.
Returns
The magnitude of the jet.
JetReconstruction.mass — Methodmass(j::FourMomentum)Return the invariant mass of a four momentum j. By convention if $m^2 < 0$, then $-sqrt{(-m^2)}$ is returned.
JetReconstruction.mass2 — Methodmass2(j::FourMomentum)Return the invariant mass squared of the four-momentum vector j.
JetReconstruction.merge_steps — Methodmerge_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.
JetReconstruction.nx — Methodnx(eej::EEJet)Return the x-component of the unit vector aligned with the momentum of eej
JetReconstruction.ny — Methodny(eej::EEJet)Return the y-component of the unit vector aligned with the momentum of eej
JetReconstruction.nz — Methodnz(eej::EEJet)Return the z-component of the unit vector aligned with the momentum of eej
JetReconstruction.open_with_stream — Methodopen_with_stream(fname::AbstractString)Open a file with a stream decompressor if it is compressed with gzip or zstd, otherwise as a normal file.
JetReconstruction.p — Methodp(j::FourMomentum)Return the momentum of the four-momentum vector of j.
JetReconstruction.p2 — Methodp2(eej::EEJet)Return the squared momentum of the EEJet object eej.
JetReconstruction.p2 — Methodp2(j::FourMomentum)Return the squared momentum of the four-momentum vector of j.
JetReconstruction.phi — Methodphi(j::FourMomentum)Return the azimuthal angle, ϕ, of the four momentum j in the range [0, 2π).
JetReconstruction.phi — Methodphi(p::PseudoJet)Return the azimuthal angle, ϕ, of a PseudoJet object p in the range [0, 2π).
JetReconstruction.pt — Methodpt(j::FourMomentum)Return the momentum of the four-momentum vector of j.
JetReconstruction.pt — Methodpt(p::PseudoJet)Compute the scalar transverse momentum (pt) of a PseudoJet.
Returns
- The transverse momentum (pt) of the PseudoJet.
JetReconstruction.pt2 — Methodpt2(j::FourMomentum)Return the squared transverse momentum of the four-momentum vector of j.
JetReconstruction.pt2 — Methodpt2(p::PseudoJet)Get the squared transverse momentum of a PseudoJet.
Returns
- The squared transverse momentum of the PseudoJet.
JetReconstruction.px — Methodpx(j::FourMomentum)Return the x-component of the four-momentum vector of j.
JetReconstruction.py — Methodpy(j::FourMomentum)Return the y-component of the four-momentum vector of j.
JetReconstruction.pz — Methodpz(j::FourMomentum)Return the z-component of the four-momentum vector of j.
JetReconstruction.rapidity — Methodrapidity(j::FourMomentum)Return the rapidity of the four momentum j.
JetReconstruction.rapidity — Methodrapidity(p::PseudoJet)Compute the rapidity of a PseudoJet object.
Returns
The rapidity of the PseudoJet object.
JetReconstruction.recluster — Methodrecluster(jet::PseudoJet, clusterseq::ClusterSequence{PseudoJet}; R = 1.0,
algorithm::JetAlgorithm.Algorithm = JetAlgorithm.CA) -> PseudoJetReclusters the constituents of a given jet jet with a different clustering algorithm algorithm and different jet radius R.
Arguments
jet::PseudoJet: The jet whose constituents are to be reclustered.clusterseq::ClusterSequence{PseudoJet}: The cluster sequence from which the original jet is obtained.R = 1.0: The new jet radius.algorithm::JetAlgorithm.Algorithm = JetAlgorithm.CA: The new clustering method.
Returns
ClusterSequence: The new cluster sequence.
JetReconstruction.reco_state — Methodreco_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: TheClusterSequenceobject to update.ranks: The ranks of the original clusters, that are inherited by pseudojets
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.
JetReconstruction.rightneighbours — Methodrightneighbours(center::Int, tiling::Tiling)Compute the indices of the right neighbors of a given center index in a tiling. This is used in the initial 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.
JetReconstruction.set_nearest_neighbours! — Methodset_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.
JetReconstruction.setup_tiling — Methodsetup_tiling(eta::Vector{T}, Rparam::AbstractFloat) where T <: AbstractFloatThis 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: ATilingDefobject 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.
JetReconstruction.surrounding — Methodsurrounding(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.
JetReconstruction.tile_index — Methodtile_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.
JetReconstruction.tile_index — Methodtile_index(sk::SoftKiller, p::PseudoJet)Return the tile index for a given PseudoJet p in the SoftKiller grid sk. Returns -1 if the jet is outside the grid bounds.
JetReconstruction.tiledjet_remove_from_tiles! — Methodtiledjet_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.
JetReconstruction.tiledjet_set_jetinfo! — Methodtiledjet_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, jetsindex, NNdist, NN, tile_index, previous, and next fields of the TiledJet object.
Returns:
nothing
JetReconstruction.upd_nn_crosscheck! — Methodupd_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.
JetReconstruction.upd_nn_nocross! — Methodupd_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.
JetReconstruction.upd_nn_step! — Methodupd_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 valid 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.
JetReconstruction.η — Functionconst η = etaAlias for the pseudorapidity function, eta.
JetReconstruction.FourMomentum — Typeabstract type FourMomentum endInterface for composite types that includes fields px, py, py, and E that represents the components of a four-momentum vector. All concrete jets that are used in the package are subtypes of this type.
JetReconstruction.HistoryElement — Typestruct HistoryElementA 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.
JetReconstruction.HistoryElement — MethodHistoryElement(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.
JetReconstruction.JetWithAncestors — Typestruct JetWithAncestors{T <: FourMomentum}A struct representing a jet with its origin ancestors.
Fields
self::T: The jet object itself.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.
JetReconstruction.Surrounding — Typestruct Surrounding{N}Structure used for iterating over neighbour tiles.
Fields
indices::NTuple{N, Int}: A tuple ofNintegers representing the indices.
JetReconstruction.TiledJet — Typemutable 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.
JetReconstruction.TiledJet — MethodTiledJet(id)Constructs a TiledJet object with the given id and initializes its properties to zero.
Arguments
id: The ID of theTiledJetobject.
Returns
A TiledJet object with the specified id and values set to zero or noTiledJet.
JetReconstruction.Tiling — Typestruct TilingThe 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 tofalseinitially, thentruewhen the tile has been setup properly).
JetReconstruction.Tiling — MethodTiling(setup::TilingDef)Constructs a initial Tiling object based on the provided setup parameters.
Arguments
setup::TilingDef: The setup parameters for the tiling.
Returns
A Tiling object.
JetReconstruction.TilingDef — Typestruct TilingDefA struct representing the definition of a specific 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.
JetReconstruction.neighbour_tiles — Typestruct neighbour_tilesA 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
XXXNote, 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
JetReconstruction.rightmost_tiles — Typestruct rightmost_tilesA 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
OOONote, 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