Jet Reconstruction Public Documentation

Documentation for JetReconstruction.jl's public interfaces.

Index

Public Methods and Types

JetReconstruction.constituent_indexesMethod
constituent_indexes(jet::T, cs::ClusterSequence{T}) where T <: FourMomentum

Return the indexes of the original particles which are the constituents of the given jet.

Arguments

  • jet::T: The jet for which to retrieve the constituents.
  • cs::ClusterSequence{T}: The cluster sequence object.

Returns

An vector of indices representing the original constituents of the given jet.

source
JetReconstruction.constituentsMethod
constituents(jet::T, cs::ClusterSequence{T}) where T <: FourMomentum

Get the constituents of a given jet in a cluster sequence.

Arguments

  • cs::ClusterSequence{T}: The cluster sequence object.
  • jet::T: The jet for which to retrieve the constituents.

Returns

An array of jet objects (which are of the same type as the input jet) representing the constituents of the given jet,

source
JetReconstruction.ee_genkt_algorithmMethod
ee_genkt_algorithm(particles::Vector{T}; p = -1, R = 4.0,
                   algorithm::JetAlgorithm.Algorithm = JetAlgorithm.Durham,
                   recombine = +) where {T}

Run an e+e- reconstruction algorithm on a set of initial particles.

Arguments

  • particles::Vector{T}: A vector of particles to be clustered.
  • p = 1: The power parameter for the algorithm. Not required / ignored for the Durham algorithm when it is set to 1.
  • R = 4.0: The jet radius parameter. Not required / ignored for the Durham algorithm.
  • algorithm::JetAlgorithm.Algorithm = JetAlgorithm.Durham: The specific jet algorithm to use.
  • recombine: The recombination scheme to use. Defaults to +.

Returns

  • The result of the jet clustering as a ClusterSequence object.

Notes

This is the public interface to the e+e- jet clustering algorithm. The function will check for consistency between the algorithm and the power parameter as needed. It will then prepare the internal EDM particles for the clustering itself, and call the actual reconstruction method _ee_genkt_algorithm.

If the algorithm is Durham, p is set to 1 and R is nominally set to 4.

Note that unlike pp reconstruction the algorithm has to be specified explicitly.

source
JetReconstruction.exclusive_jetsMethod
exclusive_jets(clusterseq::ClusterSequence{U}; dcut = nothing, njets = nothing, T = LorentzVectorCyl) where {U}

Return all exclusive jets of a ClusterSequence, with either a specific number of jets or a cut on the maximum distance parameter.

Arguments

  • clusterseq::ClusterSequence: The ClusterSequence object containing the clustering history and jets.
  • dcut::Union{Nothing, Real}: The distance parameter used to define the exclusive jets. If dcut is provided, the number of exclusive jets will be calculated based on this parameter.
  • njets::Union{Nothing, Integer}: The number of exclusive jets to be calculated. If njets is provided, the distance parameter dcut will be calculated based on this number.
  • T = LorentzVectorCyl: The return type used for the selected jets.

Note: Either dcut or njets must be provided (but not both).

Returns

  • An array of T objects representing the exclusive jets.

Valid return types are LorentzVectorCyl and the jet type of the input clusterseq (U - either PseudoJet or EEjet depending which algorithm was used) (N.B. this will evolve in the future to be any subtype of FourMomentumBase; currently unrecognised types will return LorentzVectorCyl).

Exceptions

  • ArgumentError: If neither dcut nor njets is provided.
  • ArgumentError: If the algorithm used in the ClusterSequence object is not suitable for exclusive jets.
  • ErrorException: If the cluster sequence is incomplete and exclusive jets are unavailable.

Examples

exclusive_jets(clusterseq, dcut = 20.0)
exclusive_jets(clusterseq, njets = 3, T = PseudoJet)
source
JetReconstruction.final_jetsFunction
final_jets(jets::Vector{PseudoJet}, ptmin::AbstractFloat=0.0)

This function takes a vector of PseudoJet objects and a minimum transverse momentum ptmin as input. It returns a vector of FinalJet objects that satisfy the transverse momentum condition.

Arguments

  • jets::Vector{PseudoJet}: A vector of PseudoJet objects representing the input jets.
  • ptmin::AbstractFloat=0.0: The minimum transverse momentum required for a jet to be included in the final jets vector.

Returns

A vector of FinalJet objects that satisfy the transverse momentum condition.

source
JetReconstruction.inclusive_jetsMethod
inclusive_jets(clusterseq::ClusterSequence{U}; ptmin = 0.0, T = LorentzVectorCyl) where {U}

Return all inclusive jets of a ClusterSequence with pt > ptmin.

Arguments

  • clusterseq::ClusterSequence: The ClusterSequence object containing the clustering history and jets.
  • ptmin::Float64 = 0.0: The minimum transverse momentum (pt) threshold for the inclusive jets.
  • T = LorentzVectorCyl: The return type used for the selected jets.

Returns

An array of T objects representing the inclusive jets.

Description

This function computes the inclusive jets from a given ClusterSequence object. It iterates over the clustering history and checks the transverse momentum of each parent jet. If the transverse momentum is greater than or equal to ptmin, the jet is added to the array of inclusive jets.

Valid return types are LorentzVectorCyl and the jet type of the input clusterseq (U - either PseudoJet or EEjet depending which algorithm was used) (N.B. this will evolve in the future to be any subtype of FourMomentumBase; currently unrecognised types will return LorentzVectorCyl).

Example

inclusive_jets(clusterseq; ptmin = 10.0)
source
JetReconstruction.jet_filteringMethod
jet_filtering(jet, clusterseq, filter) -> PseudoJet

Filters a jet to retain only the hardest subjets based on a specified radius and number.

Arguments:

  • jet::PseudoJet: PseudoJet instance representing the jet to filter.
  • clusterseq::ClusterSequence: ClusterSequence containing jet history.
  • filter::JetFilter: Filter instance specifying radius and number of subjets.

Returns:

  • PseudoJet: Filtered jet composed of the hardest subjets.
source
JetReconstruction.jet_reconstructMethod
jet_reconstruct(particles; p::Union{Real, Nothing} = nothing,
                     algorithm::Union{JetAlgorithm.Algorithm, Nothing} = nothing,
                     R = 1.0, recombine = +,
                     strategy::RecoStrategy.Strategy = RecoStrategy.Best)

Reconstructs jets from a collection of particles using a specified algorithm and strategy.

Arguments

  • particles: A collection of particles used for jet reconstruction.
  • p::Union{Real, Nothing} = nothing: The power value used for the distance measure for generalised k_T, which maps to a particular reconstruction algorithm (-1 = AntiKt, 0 = Cambridge/Aachen, 1 = Kt).
  • algorithm::Union{JetAlgorithm.Algorithm, Nothing} = nothing: The algorithm to use for jet reconstruction.
  • R = 1.0: The jet radius parameter.
  • recombine = +: The recombination scheme used for combining particles.
  • strategy::RecoStrategy.Strategy = RecoStrategy.Best: The jet reconstruction strategy to use. RecoStrategy.Best makes a dynamic decision based on the number of starting particles.

Note that one of p or algorithm must be specified, with algorithm preferred.

Returns

A cluster sequence object containing the reconstructed jets and the merging history.

Details

particles argument

Any type that supplies the methods pt2(), phi(), rapidity(), px(), py(), pz(), energy() (in the JetReconstruction namespace) can be used. This includes LorentzVector, LorentzVectorCyl, and PseudoJet, for which these methods are already predefined in the JetReconstruction namespace.

recombine argument

The recombine argument is the function used to merge pairs of particles. The default is simply +(jet1,jet2), i.e. 4-momenta addition or the E-scheme.

Consistency of p, algorithm and R arguments

If an algorithm is explicitly specified the p value should be consistent with it or nothing. If the algorithm is one where p can vary, then it has to be given, along with the algorithm.

If the p parameter is passed and algorithm=nothing, then pp-type reconstruction is implied (i.e., AntiKt, CA, Kt or GenKt will be used, depending on the value of p).

When an algorithm has no R dependence the R parameter is ignored.

Example

jet_reconstruct(particles; p = -1, R = 0.4)
jet_reconstruct(particles; algorithm = JetAlgorithm.Kt, R = 1.0)
jet_reconstruct(particles; algorithm = JetAlgorithm.Durham)
jet_reconstruct(particles; algorithm = JetAlgorithm.GenKt, p = 0.5, R = 1.0)
source
JetReconstruction.jet_trimmingMethod
jet_trimming(jet, clusterseq, trim) -> PseudoJet

Trims a jet by removing subjets with transverse momentum below a specified fraction.

Arguments:

  • jet::PseudoJet: PseudoJet instance representing the jet to trim.
  • clusterseq::ClusterSequence: ClusterSequence containing jet history.
  • trim::JetTrim: Trim instance specifying trimming parameters.

Returns:

  • PseudoJet: Trimmed jet composed of retained subjets.
source
JetReconstruction.kt_scaleMethod
kt_scale(jet1::T, jet2::T) where {T <: FourMomentum}

Computes the transverse momentum scale as the product of the minimum pt and the angular separation in the η-ϕ plane (using pseudorapidity).

Returns

  • The transverse momentum scale of the two jets.
source
JetReconstruction.loadjets!Method
loadjets!(filename, jets; splitby=isspace, constructor=(px,py,pz,E)->LorentzVectorHEP(E,px,py,pz), dtype=Float64)

Loads the jets from a file. Ignores lines that start with '#'. Each line gets processed in the following way: the line is split using split(line, splitby) or simply split(line) by default. Every value in this line is then converted to the dtype (which is Float64 by default). These values are then used as arguments for the constructor function which should produce individual jets. By default, the constructor constructs Lorentz vectors.

Everything that was already in jets is not affected as we only use push! on it.

Example

# Load jets from two files into one array
jets = []
loadjets!("myjets1.dat", jets)
loadjets!("myjets2.dat", jets)
source
JetReconstruction.loadjetsMethod
loadjets(filename; splitby=isspace, constructor=(px,py,pz,E)->LorentzVectorHEP(E,px,py,pz), VT=LorentzVector)

Load jets from a file.

Arguments

  • filename: The name of the file to load jets from.
  • splitby: The delimiter used to split the data in the file. Default is isspace.
  • constructor: A function that constructs a VT object from the jet data. Default is (px,py,pz,E)->LorentzVector(E,px,py,pz).
  • VT: The type of the vector used to store the jet data. Default is LorentzVector.

Returns

  • A vector of VT objects representing the loaded jets.
source
JetReconstruction.mass_dropMethod
mass_drop(jet, clusterseq, tag) -> PseudoJet

Identifies subjets in a jet that pass the mass drop tagging condition. The method stops at the first jet satisfying the mass and distance thresholds.

Arguments:

  • jet::PseudoJet: PseudoJet instance representing the jet to tag.
  • clusterseq::ClusterSequence: ClusterSequence with jet clustering history.
  • tag::MassDropTagger: MassDropTagger instance providing mass drop parameters.

Returns:

  • PseudoJet: The jet (or subjet) satisfying the mass drop conditions, if tagging is successful, otherwise a zero-momentum PseudoJet
source
JetReconstruction.n_exclusive_jetsMethod
n_exclusive_jets(clusterseq::ClusterSequence; dcut::AbstractFloat)

Return the number of exclusive jets of a ClusterSequence that are above a certain dcut value.

Arguments

  • clusterseq::ClusterSequence: The ClusterSequence object containing the clustering history.
  • dcut::AbstractFloat: The maximum value for the distance parameter in the reconstruction.

Returns

The number of exclusive jets in the ClusterSequence object.

Example

n_exclusive_jets(clusterseq, dcut = 20.0)
source
JetReconstruction.parent_jetsMethod
parent_jets(jet::T, cs::ClusterSequence{T})::Tuple{Union{Nothing, T}, Union{Nothing, T}} where {T <: FourMomentum}

Find the parent jets of a given jet in a cluster sequence.

Arguments

  • jet::T: The jet for which to find the parent jets.
  • cs::ClusterSequence: The cluster sequence object.

Returns

A tuple of two elements, each of which is either the parent jet object or nothing (if the jet has no parent).

source
JetReconstruction.plain_jet_reconstructMethod
plain_jet_reconstruct(particles::Vector{T}; p::Union{Real, Nothing} = -1,
                           algorithm::Union{JetAlgorithm.Algorithm, Nothing} = nothing,
                           R = 1.0, recombine = +) where {T}

Perform pp jet reconstruction using the plain algorithm.

Arguments

  • particles::Vector{T}: A vector of particles used for jet reconstruction, any array of particles, which supports suitable 4-vector methods, viz. pt2(), phi(), rapidity(), px(), py(), pz(), energy(), can be used. for each element.
  • p::Union{Real, Nothing} = -1: The power value used for jet reconstruction.
  • algorithm::Union{JetAlgorithm, Nothing} = nothing: The explicit jet algorithm to use.
  • R::Float64 = 1.0: The radius parameter used for jet reconstruction.
  • recombine::Function = +: The recombination function used for jet reconstruction.

Note for the particles argument, the 4-vector methods need to exist in the JetReconstruction package namespace.

This code will use the k_t algorithm types, operating in (rapidity, φ) space.

It is not necessary to specify both the algorithm and the p (power) value. If both are given they must be consistent or an exception is thrown.

Returns

  • Vector{PseudoJet}: A vector of reconstructed jets.

Example

jets = plain_jet_reconstruct(particles; p = -1, R = 0.4)
jets = plain_jet_reconstruct(particles; algorithm = JetAlgorithm.Kt, R = 1.0)
source
JetReconstruction.pt_fractionMethod
pt_fraction(jet1::T, jet2::T) where T <: FourMomentum

Computes the transverse momentum fraction of the softer of two jets.

Returns

  • The transverse momentum fraction of the softer of the two jets.
source
JetReconstruction.read_final_state_particlesMethod
read_final_state_particles(fname; maxevents = -1, skipevents = 0, T=PseudoJet)

Reads final state particles from a file and returns them as a vector of type T.

Arguments

  • fname: The name of the HepMC3 ASCII file to read particles from. If the file is gzipped, the function will automatically decompress it.
  • maxevents=-1: The maximum number of events to read. -1 means all events will be read.
  • skipevents=0: The number of events to skip before an event is included.
  • T=PseudoJet: The type of object to construct and return.

Returns

A vector of vectors of T objects, where each inner vector represents all the particles of a particular event. In particular T can be PseudoJet or a LorentzVector type. Note, if T is not PseudoJet, the order of the arguments in the constructor must be (t, x, y, z).

source
JetReconstruction.savejetsMethod
savejets(filename, jets; format="px py pz E")

Save jet data to a file.

Arguments

  • filename: The name of the file to save the jet data to.
  • jets: An array of jet objects to save.
  • format="px py pz E": (optional) A string specifying the format of the jet data to save. The default format is "px py pz E".

Details

This function saves jet data to a file in a specific format. Each line in the file represents a jet and contains the information about the jet in the specified format. The format string can include the following placeholders:

  • "E" or "energy": Jet energy
  • "px": Momentum along the x-axis
  • "py": Momentum along the y-axis
  • "pz": Momentum along the z-axis
  • "pt2": Square of the transverse momentum
  • "phi": Azimuth angle
  • "rapidity": Rapidity

Lines starting with '#' are treated as comments and are ignored.

It is strongly NOT recommended to put something other than values and (possibly custom) separators in the format string.

source
JetReconstruction.soft_dropMethod
soft_drop(jet, clusterseq, tag) -> PseudoJet

Applies soft-drop grooming to remove soft, wide-angle radiation from jets. This function reclusters the jet and iteratively checks the soft-drop condition on subjets.

Arguments:

  • jet::PseudoJet: PseudoJet instance to groom.
  • clusterseq::ClusterSequence: ClusterSequence containing jet history.
  • tag::SoftDropTagger: SoftDropTagger instance with soft-drop parameters.

Returns:

  • PseudoJet: Groomed jet or zero-momentum PseudoJet if grooming fails.
source
JetReconstruction.tiled_jet_reconstructMethod
tiled_jet_reconstruct(particles::Vector{T}; p::Union{Real, Nothing} = -1,
                           algorithm::Union{JetAlgorithm.Algorithm, Nothing} = nothing,
                           R = 1.0, recombine = +) where {T}

Main jet reconstruction algorithm entry point for reconstructing jets using the tiled strategy for generic jet type T.

Note - if a non-standard recombination is used, it must be defined for JetReconstruction.PseudoJet, as this struct is used internally.

This code will use the k_t algorithm types, operating in (rapidity, φ) space.

It is not necessary to specify both the algorithm and the p (power) value. If both are given they must be consistent or an exception is thrown.

Arguments

  • particles::Vector{T}: A vector of particles used as input for jet reconstruction. T must support methods px, py, pz and energy (defined in the JetReconstruction namespace)
  • p::Union{Real, Nothing} = -1: The power parameter for the jet reconstruction algorithm, thus switching between different algorithms.
  • algorithm::Union{JetAlgorithm.Algorithm, Nothing} = nothing: The explicit jet algorithm to use.
  • R::Float64 = 1.0: The jet radius parameter for the jet reconstruction algorithm.
  • recombine::Function = +: The recombination function used for combining pseudojets.

Returns

  • Vector{PseudoJet}: A vector of reconstructed jets.

Example

tiled_jet_reconstruct(particles::Vector{LorentzVectorHEP}; p = -1, R = 0.4, recombine = +)
source
JetReconstruction.ClusterSequenceType
struct ClusterSequence

A struct holding the full history of a jet clustering sequence, including the final jets.

Fields

  • algorithm::JetAlgorithm.Algorithm: The algorithm used for clustering.
  • strategy::RecoStrategy.Strategy: The strategy used for clustering.
  • power::Float64: The power value used for the clustering algorithm (not that this value is always stored as a Float64 to be type stable)
  • R::Float64: The R parameter used for the clustering algorithm.
  • jets::Vector{T}: The actual jets in the cluster sequence, which are of type T <: FourMomentum.
  • n_initial_jets::Int: The initial number of particles used for exclusive jets.
  • history::Vector{HistoryElement}: The branching history of the cluster sequence. Each stage in the history indicates where to look in the jets vector to get the physical PseudoJet.
  • Qtot::Any: The total energy of the event.
source
JetReconstruction.ClusterSequenceMethod
ClusterSequence(algorithm::JetAlgorithm.Algorithm, p::Real, R::Float64, strategy::RecoStrategy.Strategy, jets::T, history, Qtot) where T <: FourMomentum

Construct a ClusterSequence object.

Arguments

  • algorithm::JetAlgorithm.Algorithm: The algorithm used for clustering.
  • p::Real: The power value used for the clustering algorithm.
  • R::Float64: The R parameter used for the clustering algorithm.
  • strategy::RecoStrategy.Strategy: The strategy used for clustering.
  • jets::Vector{T}: The jets in the cluster sequence, which are of T <: FourMomentum
  • history::Vector{HistoryElement}: The branching history of the cluster sequence.
  • Qtot::Any: The total energy of the event.
source
JetReconstruction.EEjetType
mutable struct EEjet

The EEjet struct is a 4-momentum object used for the e+e jet reconstruction routines.

Fields

  • px::Float64: The x-component of the jet momentum.
  • py::Float64: The y-component of the jet momentum.
  • pz::Float64: The z-component of the jet momentum.
  • E::Float64: The energy of the jet.
  • _cluster_hist_index::Int: The index of the cluster histogram.
  • _p2::Float64: The squared momentum of the jet.
  • _inv_p::Float64: The inverse momentum of the jet.
source
JetReconstruction.FinalJetType
struct FinalJet

A struct representing the final properties of a jet, used for JSON serialisation.

Fields

  • rap::Float64: The rapidity of the jet.
  • phi::Float64: The azimuthal angle of the jet.
  • pt::Float64: The transverse momentum of the jet.
source
JetReconstruction.FinalJetsType
struct FinalJets

A struct with the vector of all jets for a certain jet identifier, used for JSON serialisation.

Fields

  • jetid::Int64: The ID of the jet.
  • jets::Vector{FinalJet}: A vector of FinalJet objects representing the jets.
source
JetReconstruction.JetFilterType
struct JetFilter

Filters jets based on radius and number of hardest subjets, reducing contamination.

Fields:

  • filter_radius::Float64: Radius parameter to recluster subjets.
  • num_hardest_jets::Int64: Number of hardest jets to retain in the filtered result.
source
JetReconstruction.JetTrimType
struct JetTrim

Trims soft, large-angle components from jets based on fraction and radius.

Fields:

  • trim_radius::Float64: Radius used for reclustering in trimming.
  • trim_fraction::Float64: Minimum momentum fraction for retained subjets.
  • recluster_method::JetAlgorithm.Algorithm: Method identifier for reclustering.
source
JetReconstruction.MassDropTaggerType
struct MassDropTagger

Used for tagging jets that undergo mass drop, a common technique in jet substructure.

Fields:

  • mu::Float64: Maximum allowed mass ratio for a jet to pass tagging.
  • y::Float64: Minimum kT distance threshold for parent separation.
source
JetReconstruction.PseudoJetType
mutable struct PseudoJet <: FourMomentum

The PseudoJet struct represents a pseudojet, a four-momentum object used in jet reconstruction algorithms. Additional information for the link back into the history of the clustering is stored in the _cluster_hist_index field. There is caching of the more expensive calculations for rapidity and azimuthal angle.

Fields

  • px::Float64: The x-component of the momentum.
  • py::Float64: The y-component of the momentum.
  • pz::Float64: The z-component of the momentum.
  • E::Float64: The energy component of the momentum.
  • _cluster_hist_index::Int: The index of the cluster history.
  • _pt2::Float64: The squared transverse momentum.
  • _inv_pt2::Float64: The inverse squared transverse momentum.
  • _rap::Float64: The rapidity.
  • _phi::Float64: The azimuthal angle.
source
JetReconstruction.PseudoJetMethod
PseudoJet(px::Real, py::Real, pz::Real, E::Real)

Constructs a PseudoJet object with the given momentum components and energy.

Arguments

  • px::Real: The x-component of the momentum.
  • py::Real: The y-component of the momentum.
  • pz::Real: The z-component of the momentum.
  • E::Real: The energy.

Returns

A PseudoJet object.

source
JetReconstruction.PseudoJetMethod
PseudoJet(px::Real, py::Real, pz::Real, E::Real,
    _cluster_hist_index::Int,
    pt2::Real)

Constructs a PseudoJet object with the given momentum components and energy and history index.

Arguments

  • px::Real: The x-component of the momentum.
  • py::Real: The y-component of the momentum.
  • pz::Real: The z-component of the momentum.
  • E::Real: The energy.
  • _cluster_hist_index::Int: The cluster history index.
  • pt2::Real: The transverse momentum squared.

Returns

A PseudoJet object.

source
JetReconstruction.SoftDropTaggerType
struct SoftDropTagger

Applies a soft-drop condition on jets, trimming away soft, wide-angle radiation.

Fields:

  • zcut::Float64: Minimum allowed energy fraction for subjets.
  • b::Float64: Angular exponent controlling soft radiation suppression.
  • cluster_rad::Float64: The new radius that will be used to recluster the components of the jet, by default set to 1.0.
source