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-
The simplest interface is to call:
cs = jet_reconstruct(particles::AbstractArray{T, 1}; algorithm = JetAlgorithm.AntiKt, R = 1.0, [p = -1,] [recombine = +,] [strategy = RecoStrategy.Best])
-
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
orEEjet
types from this package, a 4-vector fromLorentzVectorHEP
, or aReconstructedParticle
from EDM4hep are suitable (and have the appropriate definitions)
- Any type that supplies the methods
-
algorithm
is the name of the jet algorithm to be used (from theJetAlgorithm
enum)-
JetAlgorithm.AntiKt
anti-${k}_\text{T}$ clustering (default) -
JetAlgorithm.CA
Cambridge/Aachen clustering -
JetAlgorithm.Kt
inclusive$k_\text{T}$ -
JetAlgorithm.GenKt
generalised$k_\text{T}$ (which also requires specification ofp
) -
JetAlgorithm.Durham
the$e^+e-$ $k_\text{T}$ algorithm, also known as the Durham algorithm -
JetAlgorithm.EEKt
the$e^+e-$ generalised$k_\text{T}$ algorithm
-
-
R
- the cone size parameter; no particles more geometrically distance thanR
will be merged (default 1.0; note this parameter is ignored for the Durham algorithm) -
recombine
- the function used to merge two pseudojets (default is a simple 4-vector addition of$(E, \mathbf{p})$ ) -
strategy
- the algorithm strategy to adopt, as described below (defaultRecoStrategy.Best
)
The object returned is a ClusterSequence
, which internally tracks all merge steps.
Alternatively, for pp reconstruction, one can swap the algorithm=...
parameter for the value of p
, the transverse momentum power used in the
-
-1
gives anti-${k}_\text{T}$ clustering (default) -
0
gives Cambridge/Aachen -
1
gives inclusive$k_\text{T}$
Note, for the GenKt
and EEKt
algorithms the p
value must also be given to specify the algorithm fully.
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 Vector{LorentzVectorHEP}
, but different return types can be specified (e.g., T = EEjet
).
As sorting vectors is trivial in Julia, no special sorting methods are provided. As an example, to sort exclusive jets of
sorted_jets = sort!(inclusive_jets(cs::ClusterSequence; ptmin=5.0), by=JetReconstruction.energy, rev=true)
Three strategies are available for the different algorithms:
Strategy Name | Notes | Interface |
---|---|---|
RecoStrategy.Best |
Dynamically switch strategy based on input particle density | jet_reconstruct |
RecoStrategy.N2Plain |
Global matching of particles at each interation (works well for low |
plain_jet_reconstruct |
RecoStrategy.N2Tiled |
Use tiles of radius |
tiled_jet_reconstruct |
Generally one can use the jet_reconstruct
interface, shown above, as the Best strategy safely as the overhead is extremely low. 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::Vector{T}; algorithm = JetAlgorithm.AntiKt, R = 1.0, recombine = +)
Note that there is no strategy
option in these interfaces.
In the examples directory there are a number of example scripts.
See the jetreco.jl
script for an example of how to call jet reconstruction.
julia --project=examples examples/jetreco.jl --algorithm=AntiKt test/data/events.pp13TeV.hepmc3.gz
...
julia --project=examples examples/jetreco.jl --algorithm=Durham test/data/events.eeH.hepmc3.gz
...
julia --project=examples examples/jetreco.jl --maxevents=10 --strategy=N2Plain --algorithm=Kt --exclusive-njets=3 test/data/events.pp13TeV.hepmc3.gz
...
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.
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.
The package also provides methods such as loadjets
, loadjets!
, and
savejets
that one can use to save and load objects on/from disk easily in a
very flexible format. See documentation for more.
Although it has been developed further since the CHEP2023 conference, the CHEP conference proceedings, arXiv:2309.17309, should be cited if you use this package:
@misc{stewart2023polyglot,
title={Polyglot Jet Finding},
author={Graeme Andrew Stewart and Philippe Gras and Benedikt Hegner and Atell Krasnopolski},
year={2023},
eprint={2309.17309},
archivePrefix={arXiv},
primaryClass={hep-ex}
}
Code in this package is authored by:
- Atell Krasnopolski delta_atell@protonmail.com
- Graeme A Stewart graeme.andrew.stewart@cern.ch
- Philippe Gras philippe.gras@cern.ch
and is Copyright 2022-2024 The Authors, CERN.
The code is under the MIT License.