JetReconstruction.jl

March 13, 2026 · View on GitHub

All Contributors

Build Status DOI Aqua QA CITATION-cff

This package implements sequential Jet Reconstruction (clustering)

Algorithms

Algorithms used are based on the C++ FastJet package (https://fastjet.fr, hep-ph/0512210, arXiv:1111.6097), reimplemented natively in Julia.

The algorithms include anti-kT`{k}_\text{T}`, Cambridge/Aachen, inclusive kT`k_\text{T}`, generalised kT`k_\text{T}` for pp`pp` events; and the Durham algorithm and generalised kT`k_\text{T}` for e+e`e^+e^-`.

Interface

The simplest interface is to call:

cs = jet_reconstruct(particles::AbstractVector{T}; algorithm = JetAlgorithm.AntiKt, R = 1.0)
  • particles - a one dimensional array (vector) of input particles for the clustering
    • Any type that supplies the methods pt2(), phi(), rapidity(), px(), py(), pz(), energy() can be used
    • These methods have to be defined in the namespace of this package, i.e., JetReconstruction.pt2(::T)
    • The PseudoJet or EEJet types from this package, a 4-vector from LorentzVectorHEP, or a ReconstructedParticle from EDM4hep are suitable (and have the appropriate definitions)
  • algorithm is the name of the jet algorithm to be used (from the JetAlgorithm enum)
    • JetAlgorithm.AntiKt anti-kT`{k}_\text{T}` clustering
    • JetAlgorithm.CA Cambridge/Aachen clustering
    • JetAlgorithm.Kt inclusive kTk_\text{T}
    • JetAlgorithm.GenKt generalised kTk_\text{T} (which also requires specification of p)
    • JetAlgorithm.Durham the e+ee^+e- kTk_\text{T} algorithm, also known as the Durham algorithm
    • JetAlgorithm.EEKt the e+ee^+e- generalised kTk_\text{T} algorithm (which also requires specification of p)
    • JetAlgorithm.Valencia the Valencia e+ee^+e- pile-up resistant algorithm (which also requires specification of p (γ) and β)
  • R - the cone size parameter; no particles more geometrically distance than R will be merged (default 1.0; note this parameter is ignored for the Durham algorithm)

The object returned is a ClusterSequence, which internally tracks all merge steps.

For a more complete description of all possible parameters please refer to the documentation.

To obtain the final inclusive jets, use the inclusive_jets method:

final_jets = inclusive_jets(cs::ClusterSequence; ptmin=0.0)

Only jets passing the cut pT>pTminp_T > p_{Tmin} will be returned. The result is returned as a Vector{LorentzVectorHEP}, but different return types can be specified (e.g., inclusive_jets(cs::ClusterSequence, EEJet; ptmin=0.0)).

Sorting

As sorting vectors is trivial in Julia, no special sorting methods are provided. As an example, to sort exclusive jets of >5.0>5.0 (usually GeV, depending on your EDM) from highest energy to lowest:

sorted_jets = sort!(inclusive_jets(cs::ClusterSequence; ptmin=5.0), by=JetReconstruction.energy, rev=true)

Strategy

Three strategies are available for the different algorithms, but only make a difference for pp`pp` events:

Strategy NameNotesInterface
RecoStrategy.BestDynamically switch strategy based on input particle densityjet_reconstruct
RecoStrategy.N2PlainGlobal matching of particles at each interaction (works well for low NN)plain_jet_reconstruct
RecoStrategy.N2TiledUse tiles of radius RR to limit search space (works well for higher NN)tiled_jet_reconstruct

Generally one can use the jet_reconstruct interface, shown above, as the Best strategy safely as the overhead is utterly negligible. That interface supports a strategy option to switch to a different option.

Another option, if one wishes to use a specific strategy, is to call that strategy's interface directly, e.g.,

# For N2Plain strategy called directly
plain_jet_reconstruct(particles::AbstractVector{T}; algorithm = JetAlgorithm.AntiKt, R = 1.0)

Note that there is no strategy option in these interfaces.

Examples

In the examples directory there are a number of example scripts (note there is a specific Project.toml for examples).

See the jetreco.jl script for an example of how to call jet reconstruction.

julia --project jetreco.jl --algorithm=AntiKt ../test/data/events.pp13TeV.hepmc3.zst
...
julia --project jetreco.jl --algorithm=Durham ../test/data/events.eeH.hepmc3.zst
...
julia --project jetreco.jl --maxevents=10 --strategy=N2Plain --algorithm=Kt --exclusive-njets=3 ../test/data/events.pp13TeV.hepmc3.zst
...
julia --project ./jetreco.jl --algorithm=Valencia --gamma=1.2 --beta=1.2 -R 0.8 ../test/data/events.eeH.hepmc3.zst

There are options to explicitly set the algorithm (use --help to see these).

The example also shows how to use JetReconstruction.HepMC3 to read HepMC3 ASCII files (via the read_final_state_particles() wrapper).

Further examples, which show visualisation, timing measurements, profiling, etc. are given - see the README.md file in the examples directory.

Note that due to additional dependencies the Project.toml file for the examples is different from the package itself.

Plotting

illustration

To visualise the clustered jets as a 3d bar plot (see illustration above) we now use Makie.jl. See the jetsplot function in ext/JetVisualisation.jl and its documentation for more. There are two worked examples in the examples directory.

The plotting code is a package extension and will load if the one of the Makie modules is loaded in the environment.

Reference

The current recommended reference for JetReconstruction.jl is:

@article{refId0,
  author = {{Stewart, Graeme Andrew} and {Ganguly, Sanmay} and {Ghosh, Sattwamo} and {Gras, Philippe} and {Krasnopolski, Atell}},
  title = {Fast Jet Finding in Julia},
  DOI= "10.1051/epjconf/202533701067",
  url= "https://doi.org/10.1051/epjconf/202533701067",
  journal = {EPJ Web Conf.},
  year = 2025,
  volume = 337,
  pages = "01067",
}

Also available as arXiv:2503.08146.

Code in this package is Copyright 2022-2026 Graeme A Stewart, Philippe Gras, Atell Krasnopolski, Mateusz Fila, CERN.

The code is under the MIT License.

Contributors ✨

Thanks goes to these contributors to this code (emoji key):

Graeme A Stewart
Graeme A Stewart

🚇 ⚠️ 📖 💻 💡 👀
Philippe Gras
Philippe Gras

⚠️ 📖 💻 💡
Atell Krasnopolski
Atell Krasnopolski

⚠️ 📖 💻 💡
Sattwamo Ghosh
Sattwamo Ghosh

⚠️ 📖 💻 💡
Jerry Ling
Jerry Ling

📖 👀
hegner
hegner

📖
Mateusz Jakub Fila
Mateusz Jakub Fila

📖 🚇 ⚠️ 👀
ExpandingMan
ExpandingMan

🚇
emadmtr
emadmtr

⚠️ 💻 📖 💡

This project follows the all-contributors specification. Contributions of any kind welcome!