Jet Reconstruction Public Documentation
Documentation for JetReconstruction.jl
's public interfaces.
Index
JetReconstruction.ClusterSequence
JetReconstruction.ClusterSequence
JetReconstruction.ClusterSequence
JetReconstruction.FinalJet
JetReconstruction.FinalJets
JetReconstruction.PseudoJet
JetReconstruction.PseudoJet
JetReconstruction.PseudoJet
JetReconstruction.exclusive_jets
JetReconstruction.final_jets
JetReconstruction.final_jets
JetReconstruction.final_jets
JetReconstruction.inclusive_jets
JetReconstruction.jet_reconstruct
JetReconstruction.loadjets
JetReconstruction.loadjets!
JetReconstruction.n_exclusive_jets
JetReconstruction.plain_jet_reconstruct
JetReconstruction.read_final_state_particles
JetReconstruction.read_final_state_particles_lv
JetReconstruction.savejets
JetReconstruction.tiled_jet_reconstruct
JetReconstruction.tiled_jet_reconstruct
Public Methods and Types
exclusive_jets(clusterseq::ClusterSequence; dcut = nothing, njets = nothing)
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
: TheClusterSequence
object containing the clustering history and jets.dcut::Union{Nothing, Real}
: The distance parameter used to define the exclusive jets. Ifdcut
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. Ifnjets
is provided, the distance parameterdcut
will be calculated based on this number.
Note: Either dcut
or njets
must be provided (but not both).
Returns
- An array of
LorentzVectorCyl
objects representing the exclusive jets.
Exceptions
ArgumentError
: If neitherdcut
nornjets
is provided.ArgumentError
: If the algorithm used in theClusterSequence
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)
Specialisation for final jets from LorentzVectors (TODO: merge into more general function)
Specialisation for final jets from LorentzVectorCyl (TODO: merge into more general function)
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 ofPseudoJet
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.
inclusive_jets(clusterseq::ClusterSequence, ptmin = 0.0)
Return all inclusive jets of a ClusterSequence with pt > ptmin.
Arguments
clusterseq::ClusterSequence
: TheClusterSequence
object containing the clustering history and jets.ptmin::Float64 = 0.0
: The minimum transverse momentum (pt) threshold for the inclusive jets.
Returns
An array of LorentzVectorCyl
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.
Example
inclusive_jets(clusterseq, ptmin = 10.0)
jet_reconstruct(particles; p = -1, R = 1.0, recombine = +, 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=-1
: 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).R=1.0
: The jet radius parameter.recombine=+
: The recombination scheme used for combining particles.strategy=RecoStrategy.Best
: The jet reconstruction strategy to use.RecoStrategy.Best
makes a dynamic decision based on the number of starting particles.
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.
Example
jet_reconstruct(particles; p = -1, R = 0.4)
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)
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 isisspace
.constructor
: A function that constructs aVT
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 isLorentzVector
.
Returns
- A vector of
VT
objects representing the loaded jets.
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
: TheClusterSequence
object containing the clustering history.dcut::AbstractFloat
: The maximum calue 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)
plain_jet_reconstruct(particles::Vector{T}; p = -1, R = 1.0, recombine = +) where T
Perform 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::Int=-1
: The integer value used for jet reconstruction.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.
Returns
Vector{PseudoJet}
: A vector of reconstructed jets.
Example
jets = plain_jet_reconstruct(particles; p = -1, R = 1.0)
read_final_state_particles(fname; maxevents = -1, skipevents = 0)
Reads final state particles from a file and returns them as a vector of PseudoJet objects.
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
: The number of events to skip before an event is included. Default is 0.
Returns
A vector of vectors of PseudoJet objects, where each inner vector represents all the particles of a particular event.
read_final_state_particles_lv(fname; maxevents = -1, skipevents = 0)
Reads final state particles from a file and returns them as a vector of LorentzVector objects.
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
: The number of events to skip before an event is included. Default is 0.
Returns
A vector of vectors of LorentzVector objects, where each inner vector represents all the particles of a particular event.
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.
tiled_jet_reconstruct(particles::Vector{PseudoJet}; p = -1, R = 1.0, recombine = +) where {T}
Main jet reconstruction algorithm entry point for reconstructing jets using the tiled stragegy for the prefered internal jet type, PseudoJet
.
Arguments
particles::Vector{PseudoJet}
: A vector ofPseudoJet
particles used as input for jet reconstruction.p::Int = -1
: The power parameter for the jet reconstruction algorithm, thus swiching between different algorithms.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{PseudoJet}; p = 1, R = 1.0, recombine = +)
tiled_jet_reconstruct(particles::Vector{T}; p = -1, R = 1.0, recombine = +) where {T}
Main jet reconstruction algorithm entry point for reconstructing jets using the tiled stragegy 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.
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::Int = -1
: The power parameter for the jet reconstruction algorithm, thus swiching between different algorithms.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 = +)
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.jets::Vector{PseudoJet}
: The physical PseudoJets in the cluster sequence. Each PseudoJet corresponds to a position in the history.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.
ClusterSequence(p::Int, strategy::RecoStrategy.Strategy, jets, history, Qtot)
Constructs a ClusterSequence
object with a specific power value mapped to a pp reconstruction algorithm.
Arguments
p::Int
: The power value for the algorithm.strategy::RecoStrategy.Strategy
: The reconstruction strategy.jets
: The jets to be clustered.history
: The length of the jets.Qtot
: The total energy of the jets.
Returns
A ClusterSequence
object.
ClusterSequence(alg::JetAlgorithm.Algorithm, strategy::RecoStrategy.Strategy, jets, history, Qtot)
Constructs a ClusterSequence
object with a specific algorithm spcified.
Arguments
alg::JetAlgorithm.Algorithm
: The algorithm.strategy::RecoStrategy.Strategy
: The reconstruction strategy.jets
: The jets to be clustered.history
: The length of the jets.Qtot
: The total energy of the jets.
Returns
A ClusterSequence
object.
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.
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 ofFinalJet
objects representing the jets.
mutable struct PseudoJet <: FourMomentum
The PseudoJet
struct represents a pseudojet, a four-momentum object used in jet reconstruction algorithms. Additonal 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.
PseudoJet(px::Float64, py::Float64, pz::Float64, E::Float64)
Constructs a PseudoJet object with the given momentum components and energy.
Arguments
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.
Returns
A PseudoJet object.
PseudoJet(px::Float64, py::Float64, pz::Float64, E::Float64,
_cluster_hist_index::Int,
pt2::Float64)
Constructs a PseudoJet object with the given momentum components and energy and history index.
Arguments
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._cluster_hist_index::Int
: The cluster history index.pt2::Float64
: The transverse momentum squared.
Returns
A PseudoJet
object.