2 #ifndef RIVET_FastJets_HH 3 #define RIVET_FastJets_HH 5 #include "Rivet/Jet.hh" 6 #include "Rivet/Particle.hh" 7 #include "Rivet/Projection.hh" 8 #include "Rivet/Projections/JetAlg.hh" 9 #include "Rivet/Projections/FinalState.hh" 10 #include "Rivet/Tools/RivetFastJet.hh" 12 #include "fastjet/SISConePlugin.hh" 13 #include "fastjet/ATLASConePlugin.hh" 14 #include "fastjet/CMSIterativeConePlugin.hh" 15 #include "fastjet/CDFJetCluPlugin.hh" 16 #include "fastjet/CDFMidPointPlugin.hh" 17 #include "fastjet/D0RunIIConePlugin.hh" 18 #include "fastjet/TrackJetPlugin.hh" 19 #include "fastjet/JadePlugin.hh" 21 #include "Rivet/Projections/PxConePlugin.hh" 22 #include "Rivet/Tools/TypeTraits.hh" 35 CA=2, CAM=2, CAMBRIDGE=2,
38 CDFJETCLU, CDFMIDPOINT, D0ILCONE,
39 JADE, DURHAM, TRACKJET, GENKTEE };
49 const fastjet::JetDefinition& jdef,
52 fastjet::AreaDefinition* adef=
nullptr)
53 :
JetAlg(fsp, usemuons, useinvis), _jdef(jdef), _adef(adef)
62 const fastjet::JetDefinition& jdef,
63 fastjet::AreaDefinition* adef,
66 :
FastJets(fsp, jdef, usemuons, useinvis, adef)
73 fastjet::JetAlgorithm type,
74 fastjet::RecombinationScheme recom,
double rparameter,
77 fastjet::AreaDefinition* adef=
nullptr)
78 :
FastJets(fsp,
fastjet::JetDefinition(type, rparameter, recom), usemuons, useinvis, adef)
85 fastjet::JetAlgorithm type,
86 fastjet::RecombinationScheme recom,
double rparameter,
87 fastjet::AreaDefinition* adef,
90 :
FastJets(fsp, type, recom, rparameter, usemuons, useinvis, adef)
97 fastjet::JetDefinition::Plugin* plugin,
100 fastjet::AreaDefinition* adef=
nullptr)
101 :
FastJets(fsp,
fastjet::JetDefinition(plugin), usemuons, useinvis, adef)
103 _plugin.
reset(plugin);
110 fastjet::JetDefinition::Plugin* plugin,
111 fastjet::AreaDefinition* adef,
114 :
FastJets(fsp, plugin, usemuons, useinvis, adef)
125 Algo alg,
double rparameter,
128 fastjet::AreaDefinition* adef=
nullptr,
129 double seed_threshold=1.0)
130 :
JetAlg(fsp, usemuons, useinvis)
133 _initJdef(alg, rparameter, seed_threshold);
163 static Jet mkJet(
const PseudoJet& pj,
const Particles& fsparticles,
const Particles& tagparticles=Particles());
165 static Jets
mkJets(
const PseudoJets& pjs,
const Particles& fsparticles,
const Particles& tagparticles=Particles());
201 _trfs.push_back(shared_ptr<fastjet::Transformer>(trf));
208 template<
typename TRFS,
typename TRF=
typename TRFS::value_type>
209 typename std::enable_if<Derefable<TRF>::value,
void>::type
211 for (
auto& trf : trfs)
addTrf(trf);
222 Jet trimJet(
const Jet& input,
const fastjet::Filter& trimmer)
const;
267 const shared_ptr<fastjet::ClusterSequence>
clusterSeq()
const {
274 return areaDef() ? dynamic_pointer_cast<fastjet::ClusterSequenceArea>(_cseq) :
nullptr;
278 const fastjet::JetDefinition&
jetDef()
const {
286 const shared_ptr<fastjet::AreaDefinition>
areaDef()
const {
297 void _initJdef(
Algo alg,
double rparameter,
double seed_threshold);
310 void calc(
const Particles& fsparticles,
const Particles& tagparticles=Particles());
316 fastjet::JetDefinition _jdef;
319 std::shared_ptr<fastjet::AreaDefinition> _adef;
322 std::shared_ptr<fastjet::ClusterSequence> _cseq;
325 std::shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
328 std::vector< std::shared_ptr<fastjet::Transformer> > _trfs;
331 mutable std::map<int, vector<double> > _yscales;
334 Particles _fsparticles, _tagparticles;
Definition: MC_Cent_pPb.hh:10
void clearJetArea()
Don't calculate a jet area.
Definition: FastJets.hh:186
const shared_ptr< fastjet::ClusterSequenceArea > clusterSeqArea() const
Definition: FastJets.hh:273
PseudoJets pseudojetsByE(double ptmin=0.0) const
Alias.
Definition: FastJets.hh:250
DEFAULT_RIVET_PROJ_CLONE(FastJets)
Clone on the heap.
PseudoJets pseudojetsByRapidity(double ptmin=0.0) const
Alias.
Definition: FastJets.hh:257
const shared_ptr< fastjet::ClusterSequence > clusterSeq() const
Definition: FastJets.hh:267
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, fastjet::AreaDefinition *adef, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE)
Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of ...
Definition: FastJets.hh:109
void project(const Event &e)
Perform the projection on the Event.
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition: FastJets.hh:48
void useJetArea(fastjet::AreaDefinition *adef)
Use provided jet area definition.
Definition: FastJets.hh:181
PseudoJets pseudoJets(double ptmin=0.0) const
Get the pseudo jets (unordered).
const shared_ptr< fastjet::AreaDefinition > areaDef() const
Return the area definition.
Definition: FastJets.hh:286
Project out jets found using the FastJet package jet algorithms.
Definition: FastJets.hh:28
static PseudoJets mkClusterInputs(const Particles &fsparticles, const Particles &tagparticles=Particles())
Make PseudoJets for input to a ClusterSequence, with user_index codes for constituent- and tag-partic...
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition: FastJets.hh:72
Representation of a HepMC event, and enabler of Projection caching.
Definition: Event.hh:22
PseudoJets pseudoJetsByE(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition: FastJets.hh:246
Algo
Definition: FastJets.hh:33
Abstract base class for projections which can return a set of Jets.
Definition: JetFinder.hh:15
Definition: RivetFastJet.hh:14
PseudoJets pseudojets(double ptmin=0.0) const
Alias.
Definition: FastJets.hh:236
CmpState compare(const Projection &p) const
Compare projections.
Invisibles
Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding...
Definition: JetFinder.hh:22
static Jets mkJets(const PseudoJets &pjs, const Particles &fsparticles, const Particles &tagparticles=Particles())
Convert a whole list of PseudoJets to a list of Jets, with mkJet-style unpacking. ...
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE)
Definition: FastJets.hh:61
Project out all final-state particles in an event. Probably the most important projection in Rivet! ...
Definition: FinalState.hh:12
static Jet mkJet(const PseudoJet &pj, const Particles &fsparticles, const Particles &tagparticles=Particles())
Make a Rivet Jet from a PseudoJet holding a user_index code for lookup of Rivet fsparticle or tagpart...
PseudoJets pseudojetsByPt(double ptmin=0.0) const
Alias.
Definition: FastJets.hh:243
std::vector< PseudoJet > PseudoJets
Definition: RivetFastJet.hh:29
const fastjet::JetDefinition & jetDef() const
Return the jet definition.
Definition: FastJets.hh:278
void reset()
Reset the projection. Jet def, etc. are unchanged.
void calc(const Particles &fsparticles, const Particles &tagparticles=Particles())
Do the calculation locally (no caching).
Representation of a clustered jet of particles.
Definition: Jet.hh:18
void clearTrfs()
Don't apply any jet transformers.
Definition: FastJets.hh:215
PseudoJets pseudoJetsByRapidity(double ptmin=0.0) const
Get the pseudo jets, ordered by rapidity.
Definition: FastJets.hh:253
std::enable_if< Derefable< TRF >::value, void >::type addTrfs(const TRFS &trfs)
Add a list of grooming transformers.
Definition: FastJets.hh:210
Base class for all Rivet projections.
Definition: Projection.hh:29
PseudoJets pseudoJetsByPt(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition: FastJets.hh:239
Jet trimJet(const Jet &input, const fastjet::Filter &trimmer) const
Trim (filter) a jet, keeping tag and constituent info in the resulting jet.
FastJets(const FinalState &fsp, Algo alg, double rparameter, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE, fastjet::AreaDefinition *adef=nullptr, double seed_threshold=1.0)
Convenience constructor using Rivet enums for most common jet algs (including some plugins)...
Definition: FastJets.hh:124
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE)
Definition: FastJets.hh:84
void addTrf(fastjet::Transformer *trf)
Add a grooming transformer (base class of fastjet::Filter, etc.)
Definition: FastJets.hh:200
Muons
Enum for the treatment of muons: whether to include all, some, or none in jet-finding.
Definition: JetFinder.hh:19
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, JetAlg::Muons usemuons=JetAlg::Muons::ALL, JetAlg::Invisibles useinvis=JetAlg::Invisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Explicitly pass in an externally-constructed plugin.
Definition: FastJets.hh:96