rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.2
FastJets.hh
1// -*- C++ -*-
2#ifndef RIVET_FastJets_HH
3#define RIVET_FastJets_HH
4
5#include "Rivet/Jet.hh"
6#include "Rivet/Particle.hh"
7#include "Rivet/Projection.hh"
8#include "Rivet/Projections/JetFinder.hh"
9#include "Rivet/Projections/FinalState.hh"
10#include "Rivet/Tools/RivetFastJet.hh"
11
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"
20
21#include "Rivet/Projections/PxConePlugin.hh"
22#include "Rivet/Tools/TypeTraits.hh"
23
24namespace Rivet {
25
26
28 class FastJets : public JetFinder {
29 public:
30
31 using JetFinder::operator =;
32
33
36
40 FastJets(const FinalState& fsp,
41 const fastjet::JetDefinition& jdef,
42 JetMuons usemuons=JetMuons::ALL,
43 JetInvisibles useinvis=JetInvisibles::NONE,
44 fastjet::AreaDefinition* adef=nullptr)
45 : JetFinder(fsp, usemuons, useinvis), _jdef(jdef),
46 _adef(adef), _cuts(Cuts::OPEN)
47 {
48 _initBase();
49 }
50
54 FastJets(const FinalState& fsp,
55 const fastjet::JetDefinition& jdef,
56 const Cut& c,
57 JetMuons usemuons=JetMuons::ALL,
58 JetInvisibles useinvis=JetInvisibles::NONE,
59 fastjet::AreaDefinition* adef=nullptr)
60 : JetFinder(fsp, usemuons, useinvis), _jdef(jdef),
61 _adef(adef), _cuts(c)
62 {
63 _initBase();
64 }
65
69 FastJets(const FinalState& fsp,
70 const fastjet::JetDefinition& jdef,
71 fastjet::AreaDefinition* adef,
72 JetMuons usemuons=JetMuons::ALL,
73 JetInvisibles useinvis=JetInvisibles::NONE)
74 : FastJets(fsp, jdef, usemuons, useinvis, adef)
75 { }
76
81 FastJets(const FinalState& fsp,
82 const fastjet::JetDefinition& jdef,
83 fastjet::AreaDefinition* adef,
84 const Cut& c,
85 JetMuons usemuons=JetMuons::ALL,
86 JetInvisibles useinvis=JetInvisibles::NONE)
87 : FastJets(fsp, jdef, c, usemuons, useinvis, adef)
88 { }
89
90
94 FastJets(const FinalState& fsp,
95 fastjet::JetAlgorithm type,
96 fastjet::RecombinationScheme recom, double rparameter,
97 JetMuons usemuons=JetMuons::ALL,
98 JetInvisibles useinvis=JetInvisibles::NONE,
99 fastjet::AreaDefinition* adef=nullptr)
100 : FastJets(fsp, fastjet::JetDefinition(type, rparameter, recom), usemuons, useinvis, adef)
101 { }
102
107 fastjet::JetAlgorithm type,
108 fastjet::RecombinationScheme recom, double rparameter,
109 const Cut& c,
110 JetMuons usemuons=JetMuons::ALL,
111 JetInvisibles useinvis=JetInvisibles::NONE,
112 fastjet::AreaDefinition* adef=nullptr)
113 : FastJets(fsp, fastjet::JetDefinition(type, rparameter, recom), c, usemuons, useinvis, adef)
114 { }
115
120 fastjet::JetAlgorithm type,
121 fastjet::RecombinationScheme recom, double rparameter,
122 fastjet::AreaDefinition* adef,
123 JetMuons usemuons=JetMuons::ALL,
124 JetInvisibles useinvis=JetInvisibles::NONE)
125 : FastJets(fsp, type, recom, rparameter, usemuons, useinvis, adef)
126 { }
127
132 fastjet::JetAlgorithm type,
133 fastjet::RecombinationScheme recom, double rparameter,
134 fastjet::AreaDefinition* adef,
135 const Cut& c,
136 JetMuons usemuons=JetMuons::ALL,
137 JetInvisibles useinvis=JetInvisibles::NONE)
138 : FastJets(fsp, type, recom, rparameter, c, usemuons, useinvis, adef)
139 { }
140
145 fastjet::JetDefinition::Plugin* plugin,
146 JetMuons usemuons=JetMuons::ALL,
147 JetInvisibles useinvis=JetInvisibles::NONE,
148 fastjet::AreaDefinition* adef=nullptr)
149 : FastJets(fsp, fastjet::JetDefinition(plugin), usemuons, useinvis, adef)
150 {
151 _plugin.reset(plugin);
152 }
153
158 fastjet::JetDefinition::Plugin* plugin,
159 const Cut& c,
160 JetMuons usemuons=JetMuons::ALL,
161 JetInvisibles useinvis=JetInvisibles::NONE,
162 fastjet::AreaDefinition* adef=nullptr)
163 : FastJets(fsp, fastjet::JetDefinition(plugin), c, usemuons, useinvis, adef)
164 {
165 _plugin.reset(plugin);
166 }
167
172 fastjet::JetDefinition::Plugin* plugin,
173 fastjet::AreaDefinition* adef,
174 JetMuons usemuons=JetMuons::ALL,
175 JetInvisibles useinvis=JetInvisibles::NONE)
176 : FastJets(fsp, plugin, usemuons, useinvis, adef)
177 { }
178
183 fastjet::JetDefinition::Plugin* plugin,
184 fastjet::AreaDefinition* adef,
185 const Cut& c,
186 JetMuons usemuons=JetMuons::ALL,
187 JetInvisibles useinvis=JetInvisibles::NONE)
188 : FastJets(fsp, plugin, c, usemuons, useinvis, adef)
189 { }
190
199 JetAlg alg, double rparameter,
200 JetMuons usemuons=JetMuons::ALL,
201 JetInvisibles useinvis=JetInvisibles::NONE,
202 fastjet::AreaDefinition* adef=nullptr,
203 double seed_threshold=1.0)
204 : JetFinder(fsp, usemuons, useinvis),
205 _adef(adef), _cuts(Cuts::OPEN)
206 {
207 _initBase();
208 _initJdef(alg, rparameter, seed_threshold);
209 }
210
219 JetAlg alg, double rparameter,
220 const Cut& c,
221 JetMuons usemuons=JetMuons::ALL,
222 JetInvisibles useinvis=JetInvisibles::NONE,
223 fastjet::AreaDefinition* adef=nullptr,
224 double seed_threshold=1.0)
225 : JetFinder(fsp, usemuons, useinvis), _adef(adef), _cuts(c)
226 {
227 _initBase();
228 _initJdef(alg, rparameter, seed_threshold);
229 }
230
231
232
235
237
238
240 using Projection::operator =;
241
242
245
247 static PseudoJets mkClusterInputs(const Particles& fsparticles, const Particles& tagparticles=Particles());
249 static Jet mkJet(const PseudoJet& pj, const Particles& fsparticles, const Particles& tagparticles=Particles());
251 static Jets mkJets(const PseudoJets& pjs, const Particles& fsparticles, const Particles& tagparticles=Particles());
252
254
255
257 void reset();
258
259
262
267 void useJetArea(fastjet::AreaDefinition* adef) {
268 _adef.reset(adef);
269 }
270
273 _adef.reset();
274 }
275
277
278
281
286 void addTrf(fastjet::Transformer* trf) {
287 _trfs.push_back(shared_ptr<fastjet::Transformer>(trf));
288 }
289
294 template<typename TRFS, typename TRF=typename TRFS::value_type>
295 typename std::enable_if<Derefable<TRF>::value, void>::type
296 addTrfs(const TRFS& trfs) {
297 for (auto& trf : trfs) addTrf(trf);
298 }
299
301 void clearTrfs() {
302 _trfs.clear();
303 }
304
306
307
310
312 Jets _jets() const;
313
316 PseudoJets pseudojets(double ptmin=0.0) const;
317
319 PseudoJets pseudojetsByPt(double ptmin=0.0) const {
320 return sorted_by_pt(pseudojets(ptmin));
321 }
322
324 PseudoJets pseudojetsByE(double ptmin=0.0) const {
325 return sorted_by_E(pseudojets(ptmin));
326 }
327
329 PseudoJets pseudojetsByRapidity(double ptmin=0.0) const {
330 return sorted_by_rapidity(pseudojets(ptmin));
331 }
332
334
335
338
341 const shared_ptr<fastjet::ClusterSequence> clusterSeq() const {
342 return _cseq;
343 }
344
347 const shared_ptr<fastjet::ClusterSequenceArea> clusterSeqArea() const {
348 return areaDef() ? dynamic_pointer_cast<fastjet::ClusterSequenceArea>(_cseq) : nullptr;
349 }
350
352 const fastjet::JetDefinition& jetDef() const {
353 return _jdef;
354 }
355
360 const shared_ptr<fastjet::AreaDefinition> areaDef() const {
361 return _adef;
362 }
363
365
366
367 protected:
368
370 void _initBase();
371 void _initJdef(JetAlg alg, double rparameter, double seed_threshold);
372
373 protected:
374
376 void project(const Event& e);
377
379 CmpState compare(const Projection& p) const;
380
381 public:
382
384 void calc(const Particles& fsparticles, const Particles& tagparticles=Particles());
385
386
387 protected:
388
390 fastjet::JetDefinition _jdef;
391
393 std::shared_ptr<fastjet::AreaDefinition> _adef;
394
396 std::shared_ptr<fastjet::ClusterSequence> _cseq;
397
399 Cut _cuts;
400
402 std::shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
403
405 std::vector< std::shared_ptr<fastjet::Transformer> > _trfs;
406
408 mutable std::map<int, vector<double> > _yscales;
409
411 Particles _fsparticles, _tagparticles;
412
413 };
414
415}
416
417#endif
Representation of a HepMC event, and enabler of Projection caching.
Definition Event.hh:22
Find jets using jet algorithms via the FastJet package.
Definition FastJets.hh:28
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...
CmpState compare(const Projection &p) const
Compare projections.
PseudoJets pseudojetsByE(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition FastJets.hh:324
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of ...
Definition FastJets.hh:171
void clearTrfs()
Don't apply any jet transformers.
Definition FastJets.hh:301
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:69
PseudoJets pseudojetsByRapidity(double ptmin=0.0) const
Get the pseudo jets, ordered by rapidity.
Definition FastJets.hh:329
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:40
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:54
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Explicitly pass in an externally-constructed plugin.
Definition FastJets.hh:144
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:94
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Definition FastJets.hh:106
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr)
Explicitly pass in an externally-constructed plugin.
Definition FastJets.hh:157
std::enable_if< Derefable< TRF >::value, void >::type addTrfs(const TRFS &trfs)
Add a list of grooming transformers.
Definition FastJets.hh:296
FastJets(const FinalState &fsp, JetAlg alg, double rparameter, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE, fastjet::AreaDefinition *adef=nullptr, double seed_threshold=1.0)
Convenience constructor using Cut argument and Rivet enums for most common jet algs (including some p...
Definition FastJets.hh:218
FastJets(const FinalState &fsp, const fastjet::JetDefinition &jdef, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:81
void clearJetArea()
Don't calculate a jet area.
Definition FastJets.hh:272
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.
const shared_ptr< fastjet::ClusterSequence > clusterSeq() const
Definition FastJets.hh:341
PseudoJets pseudojets(double ptmin=0.0) const
FastJets(const FinalState &fsp, JetAlg alg, double rparameter, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::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:198
void calc(const Particles &fsparticles, const Particles &tagparticles=Particles())
Do the calculation locally (no caching).
const shared_ptr< fastjet::ClusterSequenceArea > clusterSeqArea() const
Definition FastJets.hh:347
void project(const Event &e)
Perform the projection on the Event.
RIVET_DEFAULT_PROJ_CLONE(FastJets)
Clone on the heap.
void useJetArea(fastjet::AreaDefinition *adef)
Use provided jet area definition.
Definition FastJets.hh:267
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...
const fastjet::JetDefinition & jetDef() const
Return the jet definition.
Definition FastJets.hh:352
void addTrf(fastjet::Transformer *trf)
Add a grooming transformer (base class of fastjet::Filter, etc.)
Definition FastJets.hh:286
FastJets(const FinalState &fsp, fastjet::JetDefinition::Plugin *plugin, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Explicitly pass in an externally-constructed plugin, with reordered args for easier specification of ...
Definition FastJets.hh:182
PseudoJets pseudojetsByPt(double ptmin=0.0) const
Get the pseudo jets, ordered by .
Definition FastJets.hh:319
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:119
void reset()
Reset the projection. Jet def, etc. are unchanged.
const shared_ptr< fastjet::AreaDefinition > areaDef() const
Return the area definition.
Definition FastJets.hh:360
FastJets(const FinalState &fsp, fastjet::JetAlgorithm type, fastjet::RecombinationScheme recom, double rparameter, fastjet::AreaDefinition *adef, const Cut &c, JetMuons usemuons=JetMuons::ALL, JetInvisibles useinvis=JetInvisibles::NONE)
Definition FastJets.hh:131
Project out all final-state particles in an event. Probably the most important projection in Rivet!
Definition FinalState.hh:12
Abstract base class for projections which can return a set of Jets.
Definition JetFinder.hh:23
Representation of a clustered jet of particles.
Definition Jet.hh:42
Specialised vector of Jet objects.
Definition Jet.hh:21
Specialised vector of Particle objects.
Definition Particle.hh:21
Base class for all Rivet projections.
Definition Projection.hh:29
Definition MC_CENT_PPB_Projections.hh:10
JetInvisibles
Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding.
Definition JetFinder.hh:18
std::vector< PseudoJet > PseudoJets
Definition RivetFastJet.hh:30
JetMuons
Enum for the treatment of muons: whether to include all, some, or none in jet-finding.
Definition JetFinder.hh:15