Jet Reconstruction Public Documentation
Documentation for JetReconstruction.jl's public interfaces.
Index
JetReconstruction.ClusterSequenceJetReconstruction.ClusterSequenceJetReconstruction.EEJetJetReconstruction.EEJetJetReconstruction.EEJetJetReconstruction.EEJetJetReconstruction.EEJetJetReconstruction.FinalJetJetReconstruction.FinalJetJetReconstruction.FinalJetsJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.PseudoJetJetReconstruction.SoftKillerJetReconstruction.addjetsJetReconstruction.addjets_pt2schemeJetReconstruction.addjets_ptschemeJetReconstruction.constituent_indexesJetReconstruction.constituentsJetReconstruction.ee_genkt_algorithmJetReconstruction.exclusive_jetsJetReconstruction.final_jetsJetReconstruction.generate_lund_emissionsJetReconstruction.inclusive_jetsJetReconstruction.jet_filteringJetReconstruction.jet_reconstructJetReconstruction.jet_trimmingJetReconstruction.kt_scaleJetReconstruction.lorentzvectorJetReconstruction.lorentzvector_cylJetReconstruction.mass_dropJetReconstruction.n_exclusive_jetsJetReconstruction.parent_jetsJetReconstruction.plain_jet_reconstructJetReconstruction.preprocess_pt2schemeJetReconstruction.preprocess_ptschemeJetReconstruction.preprocess_ptschemeJetReconstruction.preprocess_ptschemeJetReconstruction.pt_fractionJetReconstruction.read_final_state_particlesJetReconstruction.select_ABS_RAP_maxJetReconstruction.soft_dropJetReconstruction.softkillerJetReconstruction.tiled_jet_reconstruct
Public Methods and Types
JetReconstruction.addjets — Methodaddjets(jet1::T, jet2::T; cluster_hist_index::Int) where {T <: FourMomentum}Add jets' four momenta together, returning a new jet of type T with the specified cluster history index.
Details
This method is also known as the E_scheme in Fastjet.
JetReconstruction.addjets_pt2scheme — Methodaddjets_pt2scheme(jet1::T, jet2::T, cluster_hist_index::Int) where {T <: FourMomentum}Use the massless $p_T^2$ scheme for combining two jets, setting the appropriate cluster history index for the new jet.
JetReconstruction.addjets_ptscheme — Methodaddjets_ptscheme(jet1::T, jet2::T, cluster_hist_index::Int) where {T <: FourMomentum}Use the massless $p_T$ scheme for combining two jets, setting the appropriate cluster history index for the new jet.
JetReconstruction.constituent_indexes — Methodconstituent_indexes(jet::T, cs::ClusterSequence{T}) where T <: FourMomentumReturn 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.
JetReconstruction.constituents — Methodconstituents(jet::T, cs::ClusterSequence{T}) where T <: FourMomentumGet a copy of 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) copied from the constituents of the given jet, with reset cluster history indexes.
JetReconstruction.ee_genkt_algorithm — Methodee_genkt_algorithm(particles::AbstractVector{T}; algorithm::JetAlgorithm.Algorithm,
p::Union{Real, Nothing} = nothing, R = 4.0, recombine = addjets,
preprocess = nothing) where {T}Run an e+e- reconstruction algorithm on a set of initial particles.
Arguments
particles::AbstractVector{T}: A vector of particles to be clustered.algorithm::JetAlgorithm.Algorithm: The jet algorithm to use.p::Union{Real, Nothing} = nothing: The power parameter for the algorithm. This is not required for theDurhamalgorithm, but must be specified for theEEKt` algorithm.R = 4.0: The jet radius parameter. Not required and ignored for theDurhamalgorithm.recombine = addjets: The recombination scheme to use.preprocess = nothing: Preprocessing function for input particles.
Returns
- The result of the jet clustering as a
ClusterSequenceobject.
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, R is nominally set to 4. If the algorithm is EEkt, power p must also be specified.
JetReconstruction.exclusive_jets — Methodexclusive_jets(clusterseq::ClusterSequence{U}, ::Type{T} = LorentzVector{Float64}; dcut = nothing, njets = nothing) where {T, 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: TheClusterSequenceobject containing the clustering history and jets.::Type{T} = LorentzVectorCyl{Float64}: The return type used for the selected jets.dcut::Union{Nothing, Real}: The distance parameter used to define the exclusive jets. Ifdcutis 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. Ifnjetsis provided, the distance parameterdcutwill be calculated based on this number.
Note: Either dcut or njets must be provided (but not both).
Returns
- An array of
Tobjects representing the exclusive jets.
Valid return types are LorentzVector LorentzVectorCyl or the jet type of the input clusterseq (U - either PseudoJet or EEJet depending which algorithm was used).
Exceptions
ArgumentError: If neitherdcutnornjetsis provided.ArgumentError: If the algorithm used in theClusterSequenceobject 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, PseudoJet, njets = 3)JetReconstruction.final_jets — Methodfinal_jets(jets::Vector{T}, ptmin::AbstractFloat=0.0)This function takes a vector of T objects, which should be jets, and a minimum transverse momentum ptmin as input. It returns a vector of FinalJet objects that satisfy the transverse momentum condition.
Arguments
jets::Vector{T}: A vector ofTobjects 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.
JetReconstruction.generate_lund_emissions — Methodgenerate_lund_emissions(jet::PseudoJet, cs::ClusterSequence{PseudoJet})Generates the Lund plane emissions for a given jet. The jet is reclustered using the CA algorithm with a very large R to fully capture the jet structure.
Returns:
lundPoints: A vector of named tuples, each representing one step in the declustering with the following fields:h_pt: harder branch pts_pt: softer branch ptz: momentum fraction of the softer branchdelta: angular distance between brancheskt: transverse momentum of the softer branch relative to the harderpsi: azimuthal angle between brancheskappa: z * ΔR
JetReconstruction.inclusive_jets — Methodinclusive_jets(clusterseq::ClusterSequence{U}, ::Type{T} = LorentzVectorCyl{Float64}; ptmin = 0.0) where {T, U}Return all inclusive jets of a ClusterSequence with pt > ptmin.
Arguments
clusterseq::ClusterSequence: TheClusterSequenceobject containing the clustering history and jets.::Type{T} = LorentzVectorCyl{Float64}: The return type used for the selected jets.ptmin::Float64 = 0.0: The minimum transverse momentum (pt) threshold for the inclusive 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 LorentzVector LorentzVectorCyl or the jet type of the input clusterseq (U - either PseudoJet or EEJet depending which algorithm was used).
Example
inclusive_jets(clusterseq; ptmin = 10.0)JetReconstruction.jet_filtering — Methodjet_filtering(jet::PseudoJet, clusterseq::ClusterSequence{PseudoJet}; radius::Real,
hardest_jets::Integer) -> PseudoJetFilters 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{PseudoJet}: ClusterSequence containing jet history.radius::Real: Radius parameter to recluster subjets.hardest_jets::Integer: Number of hardest jets to retain in the filtered result.
Returns:
PseudoJet: Filtered jet composed of the hardest subjets.
JetReconstruction.jet_reconstruct — Methodjet_reconstruct(particles::AbstractVector; algorithm::JetAlgorithm.Algorithm,
p::Union{Real, Nothing} = nothing, R = 1.0,
recombine = addjets, preprocess = nothing,
strategy::RecoStrategy.Strategy = RecoStrategy.Best)Reconstructs jets from a collection of particles using a specified algorithm and strategy.
Arguments
particles::AbstractVector: A collection of particles used for jet reconstruction.algorithm::JetAlgorithm.Algorithm: The algorithm to use for jet reconstruction.p::Union{Real, Nothing} = nothing: The power value used for the distance measure for generalised k_T algorithms (GenKt, EEKt). Other algorithms will ignore this value.R = 1.0: The jet radius parameter.recombine = addjets: The recombination scheme used for combining particles.preprocess = nothing: The function to preprocess the particles before reconstruction (e.g., for massless schemes).nothingmeans the particles are not preprocessed.strategy::RecoStrategy.Strategy = RecoStrategy.Best: The jet reconstruction strategy to use.RecoStrategy.Bestmakes a dynamic decision based on the number of starting particles.
Note that p must be specified for GenKt and EEKt algorithms, other algorithms will ignore its value. When an algorithm has no R dependence the R parameter is ignored.
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, PseudoJet and EEJet, for which these methods are already predefined in the JetReconstruction namespace.
Note when using PseudoJet or EEJet, the history indices (_cluster_hist_index) must be set correctly. For initial jets, this means assigning to each jet its index in the vector if you construct them manually. When using other jet types correct indices are set automatically internally.
recombine argument
The recombine argument is the function used to merge pairs of particles. The default is addjets, which uses 4-momenta addition (a.k.a. the E-scheme).
preprocess argument
The preprocess argument is a function that will be called for all original input particles and which returns a new particle, usually matching a non-standard recombination scheme, e.g., massless particles for $p_T$ or $p_T^2$ recombination. nothing means no preprocessing is done.
Example
jet_reconstruct(particles; algorithm = JetAlgorithm.AntiKt, 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)
jet_reconstruct(particles; algorithm = JetAlgorithm.AntiKt, R = 1.0, preprocess = preprocess_ptscheme,
recombine = addjets_ptscheme)JetReconstruction.jet_trimming — Methodjet_trimming(jet::PseudoJet, clusterseq::ClusterSequence{PseudoJet}; radius::Real,
fraction::Real, recluster_method::JetAlgorithm.Algorithm) -> PseudoJetTrims a jet by removing subjets with transverse momentum below a specified fraction.
Arguments:
jet::PseudoJet: PseudoJet instance representing the jet to trim.clusterseq::ClusterSequence{PseudoJet}: ClusterSequence containing jet history.radius::Real: Radius used for reclustering in trimming.fraction::Real: Minimum momentum fraction for retained subjets.recluster_method::JetAlgorithm.Algorithm: Method identifier for reclustering.
Returns:
PseudoJet: Trimmed jet composed of retained subjets.
JetReconstruction.kt_scale — Methodkt_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.
JetReconstruction.lorentzvector — Methodlorentzvector(jet::T) where {T <: FourMomentum} -> -> LorentzVectorReturn a cartesian LorentzVector from a jet.
JetReconstruction.lorentzvector_cyl — Methodlorentzvector_cyl(jet::T) where T <: FourMomentum -> LorentzVectorHEP.LorentzVectorCylReturn a cylindrical LorentzVectorCyl from a jet.
JetReconstruction.mass_drop — Methodmass_drop(jet::PseudoJet, clusterseq::ClusterSequence{PseudoJet}; mu::Real,
y::Real) -> PseudoJetIdentifies 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{PseudoJet}: ClusterSequence with jet clustering history.mu::Real: Maximum allowed mass ratio for a jet to pass tagging.y::Real: Minimum kT distance threshold for parent separation.
Returns:
PseudoJet: The jet (or subjet) satisfying the mass drop conditions, if tagging is successful, otherwiseinvalid_pseudojetobject
JetReconstruction.n_exclusive_jets — Methodn_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: TheClusterSequenceobject 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)JetReconstruction.parent_jets — Methodparent_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).
JetReconstruction.plain_jet_reconstruct — Methodplain_jet_reconstruct(particles::AbstractVector{T};
algorithm::JetAlgorithm.Algorithm,
p::Union{Real, Nothing} = nothing, R = 1.0,
recombine = addjets, preprocess = nothing) where {T}Perform pp jet reconstruction using the plain algorithm.
The power value maps to specific pp jet reconstruction algorithms, but can be omitted when the algorithm implies the power value to use. It must be specified for the GenKt algorithm.
Arguments
particles::AbstractVector{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.algorithm::JetAlgorithm.Algorithm: The jet algorithm to use.p::Union{Real, Nothing} = nothing: The power value used for jet reconstruction. Must be specified for GenKt algorithm. Other algorithms will ignore this value.R = 1.0: The radius parameter used for jet reconstruction.recombine::Function = addjets: The recombination function used to combine particles into a new jet.preprocess::Function = nothing: A function to preprocess the input particles.
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.
Returns
clusterseq: The resultingClusterSequenceobject representing the reconstructed jets.
Example
jets = plain_jet_reconstruct(particles; algorithm = JetAlgorithm.Kt, R = 1.0)
jets = plain_jet_reconstruct(particles; algorithm = JetAlgorithm.GenKt, p = -0.5, R = 0.4)JetReconstruction.preprocess_pt2scheme — Functionconst preprocess_pt2scheme = preprocess_ptschemePreprocessing for p_T and p_T^2 schemes are identical.
JetReconstruction.preprocess_ptscheme — Methodpreprocess_ptscheme(jet::T; cluster_hist_index::Int = 0) -> T where {T <: FourMomentum}Jet preprocessor for the massless $p_T$ schemes, resetting the energy of the jet to be equal to the 3-momentum of the input jet.
JetReconstruction.preprocess_ptscheme — Methodpreprocess_ptscheme(jet::T, ::Type{OutputT};
cluster_hist_index::Int = 0) -> OutputT where {T <: FourMomentum,
OutputT <: FourMomentum}Jet preprocessor for the massless $p_T$ schemes, resetting the energy of the jet to be equal to the 3-momentum of the input jet.
JetReconstruction.preprocess_ptscheme — Methodpreprocess_ptscheme(particle::Union{LorentzVector, LorentzVectorCyl},
::Type{OutputT}= PseudoJet,;
cluster_hist_index::Int = 0) -> OutputT where {OutputT <: FourMomentum}Jet preprocessor for the massless $p_T$ schemes, resetting the energy of the jet to be equal to the 3-momentum of the input jet (generic particle type).
Details
This function is used to convert a particle of type LorentzVector or LorentzVectorCyl into a jet_type object, which is a subtype of FourMomentum. (This is a work around until LorentzVectorBase can be used, which will make the accessors uniform.)
JetReconstruction.pt_fraction — Methodpt_fraction(jet1::T, jet2::T) where T <: FourMomentumComputes the transverse momentum fraction of the softer of two jets.
Returns
- The transverse momentum fraction of the softer of the two jets.
JetReconstruction.read_final_state_particles — Methodread_final_state_particles(fname, ::Type{T} = PseudoJet; maxevents = -1, skipevents = 0) where {T}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.::Type{T}=PseudoJet: The type of object to construct and return.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.
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).
JetReconstruction.select_ABS_RAP_max — Methodselect_ABS_RAP_max(event, absrapmax)Filter a collection of PseudoJets, returning only those with absolute rapidity less than or equal to absrapmax.
JetReconstruction.soft_drop — Methodsoft_drop(jet::PseudoJet, clusterseq::ClusterSequence{PseudoJet}; zcut::Real,
beta::Real, radius::Real = 1.0) -> PseudoJetApplies 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{PseudoJet}: ClusterSequence containing jet history.zcut::Real: Minimum allowed energy fraction for subjets.beta::Real: Angular exponent controlling soft radiation suppression.radius::Real: The new radius that will be used to recluster the components of the jet, by default set to 1.0.
Returns:
PseudoJet: Groomed jet orinvalid_pseudojetobject if grooming fails.
JetReconstruction.softkiller — Methodsoftkiller(sk::SoftKiller, event::Vector{PseudoJet})Apply the SoftKiller algorithm to an event (a vector of PseudoJets). Returns a tuple (reduced_event, pt_threshold), where reduced_event is the filtered event and pt_threshold is the computed pt threshold.
JetReconstruction.tiled_jet_reconstruct — Methodtiled_jet_reconstruct(particles::AbstractVector{T};
algorithm::JetAlgorithm.Algorithm,
p::Union{Real, Nothing} = nothing, R = 1.0,
recombine = addjets, preprocess = nothing) where {T}Main jet reconstruction algorithm entry point for reconstructing jets using the tiled strategy for generic jet type T.
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::AbstractVector{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).algorithm::JetAlgorithm.Algorithm: The jet algorithm to use.p::Union{Real, Nothing} = nothing: The power value used for jet reconstruction. Must be specified for GenKt algorithm. Other algorithms will ignore this value.R = 1.0: The jet radius parameter for the jet reconstruction algorithm.recombine::Function = addjets: The recombination function used to combine particles into a new jet.preprocess::Function = nothing: A function to preprocess the input particles.
Returns
ClusterSequence: The resultingClusterSequenceobject representing the reconstructed jets.
Example
tiled_jet_reconstruct(particles::Vector{LorentzVectorHEP}; algorithm = JetAlgorithm.GenKt, p = 0.5, R = 0.4)
tiled_jet_reconstruct(particles::Vector{LorentzVectorHEP}; algorithm = JetAlgorithm.AntiKt, R = 0.4)JetReconstruction.ClusterSequence — Typestruct ClusterSequenceA struct holding the full history of a jet clustering sequence, including the final jets.
Fields
algorithm::JetAlgorithm.Algorithm: The algorithm 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.strategy::RecoStrategy.Strategy: The strategy used for clustering.jets::Vector{T}: The actual jets in the cluster sequence, which are of typeT <: 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::Float64: The total energy of the event.
JetReconstruction.ClusterSequence — MethodClusterSequence(algorithm::JetAlgorithm.Algorithm, p, R,
strategy::RecoStrategy.Strategy, jets::Vector{T}, history,
Qtot) where {T <: FourMomentum}Construct a ClusterSequence object.
Arguments
algorithm::JetAlgorithm.Algorithm: The algorithm used for clustering.p: The power value used for the clustering algorithm.R: 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 <: FourMomentumhistory::Vector{HistoryElement}: The branching history of the cluster sequence.Qtot: The total energy of the event.
JetReconstruction.EEJet — Typestruct EEJet <: FourMomentumThe 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.
JetReconstruction.EEJet — MethodEEJet(px::Real, py::Real, pz::Real, E::Real, cluster_hist_index::Int)Constructs an EEJet object from the given momentum components, energy, and cluster history index.
Details
If the default value of cluster_hist_index=0 is used, the EEJet cannot be used in a reconstruction sequence.
JetReconstruction.EEJet — MethodEEJet(jet::Any; cluster_hist_index::Int = 0)Construct a EEJet from a generic object jet with the given cluster index. This functions also for LorentzVectorCyl objects.
Details
This function is used to convert a generic object jet into an EEJet, for this to work the object must have the methods px, py, pz, and energy defined, which are used to extract the four-momentum components of the object.
The cluster_hist_index is optional, but needed if the jet is part of a reconstruction sequence. If not provided, it defaults to 0 as an "invalid" value.
JetReconstruction.EEJet — MethodEEJet(jet::LorentzVector; cluster_hist_index::Int = 0)Construct a EEJet from a LorentzVector object with optional cluster index.
JetReconstruction.EEJet — MethodEEJet(jet::PseudoJet; cluster_hist_index::Int=0)Constructs an EEJet from a PseudoJet.
Details
The cluster_hist_index is set to the value of the cluster_hist_index of the PseudoJet if 0 is passed. Otherwise it is set to the value, >0, passed in.
JetReconstruction.FinalJet — Typestruct FinalJetA 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.
JetReconstruction.FinalJet — MethodFinalJet(jet::T)Convert a jet of type T to a FinalJet object. T must be a type that supports the methods rapidity, phi, and pt.
JetReconstruction.FinalJets — Typestruct FinalJetsA 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 ofFinalJetobjects representing the jets.
JetReconstruction.PseudoJet — Typestruct PseudoJet <: FourMomentumThe 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.
JetReconstruction.PseudoJet — MethodPseudoJet(px::Real, py::Real, pz::Real, E::Real; cluster_hist_index::Int = 0)Construct a PseudoJet from a four momentum (px, py, pz, E)with cluster indexclusterhistindex`.
Details
If the (default) value of cluster_hist_index=0 is used, the PseudoJet cannot be used in a reconstruction sequence.
JetReconstruction.PseudoJet — MethodPseudoJet(jet::Any; cluster_hist_index::Int = 0)Construct a PseudoJet from a generic object jet with the given cluster index.
Details
This function is used to convert a generic object jet into a PseudoJet, for this to work the object must have the methods px, py, pz, and energy defined, which are used to extract the four-momentum components of the object.
The cluster_hist_index is optional, but needed if the jet is part of a reconstruction sequence. If not provided, it defaults to 0 as an "invalid" value.
JetReconstruction.PseudoJet — MethodPseudoJet(jet::EEJet)Constructs a PseudoJet object from an EEJet object, with the same four momentum and cluster history index.
JetReconstruction.PseudoJet — MethodPseudoJet(jet::LorentzVectorCyl; cluster_hist_index::Int = 0)Construct a PseudoJet from a LorentzVectorCyl object with the given cluster index.
JetReconstruction.PseudoJet — MethodPseudoJet(jet::LorentzVector; cluster_hist_index::Int = 0)Construct a PseudoJet from a LorentzVector object with the cluster index.
JetReconstruction.PseudoJet — MethodPseudoJet(;pt::Real, rap::Real, phi::Real, m::Real = 0, cluster_hist_index::Int = 0)Construct a PseudoJet from (pt, y, ϕ, m) with the cluster index cluster_hist_index.
Details
If the (default) value of cluster_hist_index=0 is used, the PseudoJet cannot be used in a reconstruction sequence.
JetReconstruction.SoftKiller — Typestruct SoftKillerImplements the SoftKiller pileup mitigation algorithm, inspired by the FastJet contrib package.
The original algorithm is described in: https://arxiv.org/abs/1407.0408 by Matteo Cacciari, Gavin P. Salam, Gregory Soyez
This version inspired by the SoftKiller implementation in the FastJet contrib package Original C++ code: https://fastjet.hepforge.org/contrib/
The SoftKiller algorithm divides the rapidity-phi plane into a grid of tiles and determines a dynamic pt threshold by finding the median of the maximum pt in each tile. Particles with pt below this threshold are removed from the event.
Fields
_ymax::Float64: Maximum rapidity of the grid._ymin::Float64: Minimum rapidity of the grid._requested_drap::Float64: Requested grid spacing in rapidity._requested_dphi::Float64: Requested grid spacing in phi._ntotal::Int: Total number of tiles._dy::Float64: Actual grid spacing in rapidity._dphi::Float64: Actual grid spacing in phi._cell_area::Float64: Area of a single tile._inverse_dy::Float64: Inverse of rapidity grid spacing._inverse_dphi::Float64: Inverse of phi grid spacing._ny::Int: Number of tiles in rapidity._nphi::Int: Number of tiles in phi.
Constructors
SoftKiller(rapmin::Float64, rapmax::Float64, drap::Float64, dphi::Float64): Construct a grid fromrapmintorapmaxin rapidity, with tile sizesdrapanddphi.SoftKiller(rapmax::Float64, grid_size::Float64): Construct a square grid from-rapmaxtorapmaxin rapidity, with tile sizegrid_size.