rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.3
Jet.hh
1// -*- C++ -*-
2#ifndef RIVET_Jet_HH
3#define RIVET_Jet_HH
4
5#include "Rivet/Config/RivetCommon.hh"
6#include "Rivet/Jet.fhh"
7#include "Rivet/Particle.hh"
8#include "Rivet/Tools/Cuts.hh"
9#include "Rivet/Tools/Utils.hh"
10#include "Rivet/Tools/RivetFastJet.hh"
11#include "Rivet/Math/LorentzTrans.hh"
12#include <numeric>
13
14namespace Rivet {
15
16
21 class Jets : public std::vector<Jet> {
22 public:
23 using base = std::vector<Jet>; //< using-declarations don't like template syntax
24 using base::base; //< import base-class constructors
25 Jets();
26 Jets(const std::vector<Jet>& vjs);
27 FourMomenta moms() const;
28 PseudoJets pseudojets() const;
29 operator FourMomenta () const { return moms(); }
30 operator PseudoJets () const { return pseudojets(); }
31 Jets& operator += (const Jet& j);
32 Jets& operator += (const Jets& js);
33 };
34
35 Jets operator + (const Jets& a, const Jets& b);
36
37
39
40
42 class Jet : public ParticleBase {
43 public:
44
47
49 Jet(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles()) {
51 }
52
57
59 Jet() { clear(); }
60
62
63
66
68 size_t size() const { return _particles.size(); }
69
71 Particles& particles() { return _particles; }
73 const Particles& particles() const { return _particles; }
75 const Particles particles(const Cut& c) const { return select(_particles, c); }
77 const Particles particles(const ParticleSelector& s) const { return select(_particles, s); }
78
82 const Particles& constituents() const { return particles(); }
84 const Particles constituents(const Cut& c) const { return particles(c); }
86 const Particles constituents(const ParticleSelector& s) const { return particles(s); }
87
91 bool containsPID(const Particle& particle) const { return containsParticle(particle); }
92
96 bool containsPID(PdgId pid) const { return containsParticleId(pid); }
97
99 bool containsParticleId(const vector<PdgId>& pids) const;
101 bool containsPID(const vector<PdgId>& pids) const { return containsParticleId(pids); }
102
104
105
111
113 Particles& tags() { return _tags; }
115 const Particles& tags() const { return _tags; }
119 Particles tags(const ParticleSelector& f, double dRmax=-1) const {
120 const Particles filttags = select(tags(), f);
121 // Return, via dR filter if requested
122 return dRmax < 0 ? filttags : select(filttags, deltaRLess(this->mom(), dRmax));
123 }
127 Particles tags(const Cut& c, double dRmax=-1) const;
128
129
136 Particles bTags(const Cut& c=Cuts::open(), double dRmax=-1) const;
137
144 Particles bTags(const ParticleSelector& f, double dRmax=-1) const { return select(bTags(), f); }
145
147 bool bTagged(const Cut& c=Cuts::open(), double dRmax=-1) const { return !bTags(c).empty(); }
148
150 bool bTagged(const ParticleSelector& f, double dRmax=-1) const { return !bTags(f).empty(); }
151
152
159 Particles cTags(const Cut& c=Cuts::open(), double dRmax=-1) const;
160
167 Particles cTags(const ParticleSelector& f, double dRmax=-1) const { return select(cTags(), f); }
168
170 bool cTagged(const Cut& c=Cuts::open(), double dRmax=-1) const { return !cTags(c).empty(); }
171
173 bool cTagged(const ParticleSelector& f, double dRmax=-1) const { return !cTags(f).empty(); }
174
175
182 Particles tauTags(const Cut& c=Cuts::open(), double dRmax=-1) const;
183
190 Particles tauTags(const ParticleSelector& f, double dRmax=-1) const { return select(tauTags(), f); }
191
193 bool tauTagged(const Cut& c=Cuts::open(), double dRmax=-1) const { return !tauTags(c).empty(); }
194
196 bool tauTagged(const ParticleSelector& f, double dRmax=-1) const { return !tauTags(f).empty(); }
197
198
200
202
203
206
208 const FourMomentum& momentum() const { return _momentum; }
209
214
216 double totalEnergy() const { return momentum().E(); }
217
219 double neutralEnergy() const;
220
222 double hadronicEnergy() const;
223
225
226
229
231 const fastjet::PseudoJet& pseudojet() const { return _pseudojet; }
232
234 operator const fastjet::PseudoJet& () const { return pseudojet(); }
235
237
238
241
246 Jet& setState(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles());
247
250
255 Jet& setConstituents(const Particles& particles) { return setParticles(particles); }
256
259
261
262
263 private:
264
266 fastjet::PseudoJet _pseudojet;
267
270 Particles _particles;
271
273 Particles _tags;
274
276 mutable FourMomentum _momentum;
277
278 // /// Provide but hide the equality operators, to avoid implicit comparison via fastjet::PseudoJet
279 // bool operator == (const Jet&) const;
280 // bool operator != (const Jet&) const;
281
282 };
283
284
287
289 std::ostream& operator << (std::ostream& os, const Jet& j);
290
292
293
294}
295
296
297#include "Rivet/Tools/JetUtils.hh"
298
299#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:316
double E() const
Get energy (time component of momentum).
Definition Vector4.hh:553
Representation of a clustered jet of particles.
Definition Jet.hh:42
double hadronicEnergy() const
Get the energy carried in this jet by hadrons.
bool containsParticleId(const vector< PdgId > &pids) const
Check whether this jet contains at least one of certain particle types.
const fastjet::PseudoJet & pseudojet() const
Access the internal FastJet3 PseudoJet (as a const reference)
Definition Jet.hh:231
const Particles particles(const ParticleSelector &s) const
Get the particles in this jet which pass a filtering functor (const)
Definition Jet.hh:77
Jet & transformBy(const LorentzTransform &lt)
size_t size() const
Number of particles in this jet.
Definition Jet.hh:68
Jet(const fastjet::PseudoJet &pj, const Particles &particles=Particles(), const Particles &tags=Particles())
Constructor from a FastJet PseudoJet, with optional full particle constituents information.
Definition Jet.hh:49
const Particles constituents(const Cut &c) const
Get the particles in this jet which pass a cut (FastJet-like alias, const)
Definition Jet.hh:84
bool cTagged(const Cut &c=Cuts::open(), double dRmax=-1) const
Does this jet have at least one c-tag? (with optional Cut and dR restriction)
Definition Jet.hh:170
Jet & setParticles(const Particles &particles)
Set the particles collection with full particle information.
bool containsParticleId(PdgId pid) const
Check whether this jet contains a certain particle type.
Particles tags(const Cut &c, double dRmax=-1) const
Particles which have been tag-matched to this jet and pass a Cut or dR requirement.
const Particles & constituents() const
Get the particles in this jet (FastJet-like alias, const version)
Definition Jet.hh:82
Jet(const FourMomentum &pjet, const Particles &particles=Particles(), const Particles &tags=Particles())
Set the jet data, with optional full particle information.
Definition Jet.hh:54
Jet & clear()
Reset this jet as empty.
bool tauTagged(const ParticleSelector &f, double dRmax=-1) const
Does this jet have at least one tau-tag (with optional selector function and dR restriction)
Definition Jet.hh:196
Particles cTags(const ParticleSelector &f, double dRmax=-1) const
Get the c (and not b) particles which have been tag-matched to this jet.
Definition Jet.hh:167
const Particles constituents(const ParticleSelector &s) const
Get the particles in this jet which pass a filtering functor (FastJet-like alias, const)
Definition Jet.hh:86
const Particles & tags() const
Particles which have been tag-matched to this jet (const version)
Definition Jet.hh:115
Particles & particles()
Get the particles in this jet.
Definition Jet.hh:71
Jet & setState(const FourMomentum &mom, const Particles &particles, const Particles &tags=Particles())
Set all the jet data, with optional full particle constituent and tag information.
bool containsPID(PdgId pid) const
Nicer alias for containsParticleId.
Definition Jet.hh:96
bool tauTagged(const Cut &c=Cuts::open(), double dRmax=-1) const
Does this jet have at least one tau-tag (with optional Cut and dR restriction)
Definition Jet.hh:193
const Particles particles(const Cut &c) const
Get the particles in this jet which pass a cut (const)
Definition Jet.hh:75
Jet & setState(const fastjet::PseudoJet &pj, const Particles &particles=Particles(), const Particles &tags=Particles())
Set the jet data from a FastJet PseudoJet, with optional particle constituents and tags lists.
const FourMomentum & momentum() const
Get equivalent single momentum four-vector.
Definition Jet.hh:208
Jet()
Default constructor – only for STL storability.
Definition Jet.hh:59
Particles bTags(const Cut &c=Cuts::open(), double dRmax=-1) const
Get the b particles tag-matched to this jet.
Particles & constituents()
Get the particles in this jet (FastJet-like alias)
Definition Jet.hh:80
bool containsPID(const Particle &particle) const
Nicer alias for containsParticleId.
Definition Jet.hh:91
Particles & tags()
Particles which have been tag-matched to this jet.
Definition Jet.hh:113
Particles bTags(const ParticleSelector &f, double dRmax=-1) const
Get the b particles tag-matched to this jet (with optional selector function and dR restriction)
Definition Jet.hh:144
bool bTagged(const Cut &c=Cuts::open(), double dRmax=-1) const
Does this jet have at least one b-tag? (with optional Cut and dR restriction)
Definition Jet.hh:147
Particles tauTags(const Cut &c=Cuts::open(), double dRmax=-1) const
Get the tau particles tag-matched to this jet.
Particles tags(const ParticleSelector &f, double dRmax=-1) const
Particles which have been tag-matched to this jet and pass a selector function or dR requirement.
Definition Jet.hh:119
bool cTagged(const ParticleSelector &f, double dRmax=-1) const
Does this jet have at least one c-tag? (with optional selector function and dR restriction)
Definition Jet.hh:173
Particles cTags(const Cut &c=Cuts::open(), double dRmax=-1) const
Get the c (and not b) particles tag-matched to this jet.
const Particles & particles() const
Get the particles in this jet (const version)
Definition Jet.hh:73
double neutralEnergy() const
Get the energy carried in this jet by neutral particles.
bool bTagged(const ParticleSelector &f, double dRmax=-1) const
Does this jet have at least one b-tag? (with optional selector function and dR restriction)
Definition Jet.hh:150
double totalEnergy() const
Get the total energy of this jet.
Definition Jet.hh:216
bool containsParticle(const Particle &particle) const
Check whether this jet contains a particular particle.
Particles tauTags(const ParticleSelector &f, double dRmax=-1) const
Get the tau particles tag-matched to this jet.
Definition Jet.hh:190
bool containsPID(const vector< PdgId > &pids) const
Nicer alias for containsParticleId.
Definition Jet.hh:101
Specialised vector of Jet objects.
Definition Jet.hh:21
Object implementing Lorentz transform calculations and boosts.
Definition LorentzTrans.hh:21
Base class for particle-like things like Particle and Jet.
Definition ParticleBase.hh:13
const FourMomentum & mom() const
Get equivalent single momentum four-vector (const) (alias).
Definition ParticleBase.hh:39
double eta() const
Get the directly (alias).
Definition ParticleBase.hh:87
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition Particle.hh:45
Specialised vector of Particle objects.
Definition Particle.hh:21
Jets select(const Jets &jets, const Cut &c)
Filter a jet collection in-place to the subset that passes the supplied Cut.
Definition JetUtils.hh:152
int pid(const Particle &p)
Unbound function access to PID code.
Definition ParticleUtils.hh:23
PdgIdPair pids(const ParticlePair &pp)
Get the PDG ID codes of a ParticlePair.
Definition ParticleUtils.hh:717
const Cut & open()
Fully open cut singleton, accepts everything.
Definition MC_CENT_PPB_Projections.hh:10
std::ostream & operator<<(std::ostream &os, const AnalysisInfo &ai)
Stream an AnalysisInfo as a text description.
Definition AnalysisInfo.hh:362
std::vector< PseudoJet > PseudoJets
Definition RivetFastJet.hh:30
(with respect to another 4-momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:205