rivet is hosted by Hepforge, IPPP Durham
Rivet 3.1.6
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
23 class Jets : public std::vector<Jet> {
24 public:
25 using base = std::vector<Jet>; //< using-declarations don't like template syntax
26 using base::base; //< import base-class constructors
27 Jets();
28 Jets(const std::vector<Jet>& vjs);
29 FourMomenta moms() const;
30 PseudoJets pseudojets() const;
31 operator FourMomenta () const { return moms(); }
32 operator PseudoJets () const { return pseudojets(); }
33 Jets& operator += (const Jet& j);
34 Jets& operator += (const Jets& js);
35 };
36
37 Jets operator + (const Jets& a, const Jets& b);
38
40
41
42
44
45
46
48 class Jet : public ParticleBase {
49 public:
50
52
53
55 Jet(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles()) {
57 }
58
61 setState(pjet, particles, tags);
62 }
63
65 Jet() { clear(); }
66
68
69
71
72
74 size_t size() const { return _particles.size(); }
75
77 Particles& particles() { return _particles; }
79 const Particles& particles() const { return _particles; }
81 const Particles particles(const Cut& c) const { return filter_select(_particles, c); }
83 const Particles particles(const ParticleSelector& s) const { return filter_select(_particles, s); }
84
88 const Particles& constituents() const { return particles(); }
90 const Particles constituents(const Cut& c) const { return particles(c); }
92 const Particles constituents(const ParticleSelector& s) const { return particles(s); }
93
95 bool containsParticle(const Particle& particle) const;
97 bool containsPID(const Particle& particle) const { return containsParticle(particle); }
98
100 bool containsParticleId(PdgId pid) const;
102 bool containsPID(PdgId pid) const { return containsParticleId(pid); }
103
105 bool containsParticleId(const vector<PdgId>& pids) const;
107 bool containsPID(const vector<PdgId>& pids) const { return containsParticleId(pids); }
108
110
111
116
117
119 Particles& tags() { return _tags; }
121 const Particles& tags() const { return _tags; }
125 Particles tags(const ParticleSelector& f) const { return filter_select(tags(), f); }
129 Particles tags(const Cut& c) const;
130
131
135 Particles bTags(const Cut& c=Cuts::open()) const;
137 Particles bTags(const ParticleSelector& f) const { return filter_select(bTags(), f); }
138
140 bool bTagged(const Cut& c=Cuts::open()) const { return !bTags(c).empty(); }
142 bool bTagged(const ParticleSelector& f) const { return !bTags(f).empty(); }
143
144
148 Particles cTags(const Cut& c=Cuts::open()) const;
150 Particles cTags(const ParticleSelector& f) const { return filter_select(cTags(), f); }
151
153 bool cTagged(const Cut& c=Cuts::open()) const { return !cTags(c).empty(); }
155 bool cTagged(const ParticleSelector& f) const { return !cTags(f).empty(); }
156
157
161 Particles tauTags(const Cut& c=Cuts::open()) const;
163 Particles tauTags(const ParticleSelector& f) const { return filter_select(tauTags(), f); }
164
166 bool tauTagged(const Cut& c=Cuts::open()) const { return !tauTags(c).empty(); }
168 bool tauTagged(const ParticleSelector& f) const { return !tauTags(f).empty(); }
169
171
172
174
175
177 const FourMomentum& momentum() const { return _momentum; }
178
183
185 double totalEnergy() const { return momentum().E(); }
186
188 double neutralEnergy() const;
189
191 double hadronicEnergy() const;
192
194
195
197
198
200 const fastjet::PseudoJet& pseudojet() const { return _pseudojet; }
201
203 operator const fastjet::PseudoJet& () const { return pseudojet(); }
204
206
207
209
210
215 Jet& setState(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles());
216
219
224 Jet& setConstituents(const Particles& particles) { return setParticles(particles); }
225
228
230
231
232 private:
233
235 fastjet::PseudoJet _pseudojet;
236
239 Particles _particles;
240
242 Particles _tags;
243
245 mutable FourMomentum _momentum;
246
247 // /// Provide but hide the equality operators, to avoid implicit comparison via fastjet::PseudoJet
248 // bool operator == (const Jet&) const;
249 // bool operator != (const Jet&) const;
250
251 };
252
253
255
256
258 std::ostream& operator << (std::ostream& os, const Jet& j);
259
261
262
263}
264
265
266#include "Rivet/Tools/JetUtils.hh"
267
268#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:306
Representation of a clustered jet of particles.
Definition: Jet.hh:48
Particles cTags(const ParticleSelector &f) const
c (and not b) particles which have been tag-matched to this jet and pass a selector function
Definition: Jet.hh:150
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.
Particles bTags(const Cut &c=Cuts::open()) const
b particles which have been tag-matched to this jet (and pass an optional Cut)
const fastjet::PseudoJet & pseudojet() const
Access the internal FastJet3 PseudoJet (as a const reference)
Definition: Jet.hh:200
const Particles particles(const ParticleSelector &s) const
Get the particles in this jet which pass a filtering functor (const)
Definition: Jet.hh:83
Particles tauTags(const Cut &c=Cuts::open()) const
Tau particles which have been tag-matched to this jet (and pass an optional Cut)
Jet & transformBy(const LorentzTransform &lt)
size_t size() const
Number of particles in this jet.
Definition: Jet.hh:74
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:55
const Particles constituents(const Cut &c) const
Get the particles in this jet which pass a cut (FastJet-like alias, const)
Definition: Jet.hh:90
bool tauTagged(const ParticleSelector &f) const
Does this jet have at least one tau-tag (that passes the supplied selector function)?
Definition: Jet.hh:168
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 bTags(const ParticleSelector &f) const
b particles which have been tag-matched to this jet and pass a selector function
Definition: Jet.hh:137
const Particles & constituents() const
Get the particles in this jet (FastJet-like alias, const version)
Definition: Jet.hh:88
Jet(const FourMomentum &pjet, const Particles &particles=Particles(), const Particles &tags=Particles())
Set the jet data, with optional full particle information.
Definition: Jet.hh:60
Jet & clear()
Reset this jet as empty.
bool cTagged(const Cut &c=Cuts::open()) const
Does this jet have at least one c-tag (that passes an optional Cut)?
Definition: Jet.hh:153
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:92
const Particles & tags() const
Particles which have been tag-matched to this jet (const version)
Definition: Jet.hh:121
bool bTagged(const Cut &c=Cuts::open()) const
Does this jet have at least one b-tag (that passes an optional Cut)?
Definition: Jet.hh:140
Particles & particles()
Get the particles in this jet.
Definition: Jet.hh:77
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:102
const Particles particles(const Cut &c) const
Get the particles in this jet which pass a cut (const)
Definition: Jet.hh:81
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.
bool bTagged(const ParticleSelector &f) const
Does this jet have at least one b-tag (that passes the supplied selector function)?
Definition: Jet.hh:142
Particles cTags(const Cut &c=Cuts::open()) const
c (and not b) particles which have been tag-matched to this jet (and pass an optional Cut)
const FourMomentum & momentum() const
Get equivalent single momentum four-vector.
Definition: Jet.hh:177
Jet()
Default constructor – only for STL storability.
Definition: Jet.hh:65
Particles tags(const ParticleSelector &f) const
Particles which have been tag-matched to this jet and pass a selector function.
Definition: Jet.hh:125
Particles & constituents()
Get the particles in this jet (FastJet-like alias)
Definition: Jet.hh:86
bool containsPID(const Particle &particle) const
Nicer alias for containsParticleId.
Definition: Jet.hh:97
Particles & tags()
Particles which have been tag-matched to this jet.
Definition: Jet.hh:119
bool tauTagged(const Cut &c=Cuts::open()) const
Does this jet have at least one tau-tag (that passes an optional Cut)?
Definition: Jet.hh:166
Particles tauTags(const ParticleSelector &f) const
Tau particles which have been tag-matched to this jet and pass a selector function.
Definition: Jet.hh:163
const Particles & particles() const
Get the particles in this jet (const version)
Definition: Jet.hh:79
double neutralEnergy() const
Get the energy carried in this jet by neutral particles.
bool cTagged(const ParticleSelector &f) const
Does this jet have at least one c-tag (that passes the supplied selector function)?
Definition: Jet.hh:155
double totalEnergy() const
Get the total energy of this jet.
Definition: Jet.hh:185
bool containsParticle(const Particle &particle) const
Check whether this jet contains a particular particle.
bool containsPID(const vector< PdgId > &pids) const
Nicer alias for containsParticleId.
Definition: Jet.hh:107
Particles tags(const Cut &c) const
Particles which have been tag-matched to this jet and pass a Cut.
Specialised vector of Jet objects.
Definition: Jet.hh:23
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
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition: Particle.hh:53
Specialised vector of Particle objects.
Definition: Particle.hh:25
Jets filter_select(const Jets &jets, const Cut &c)
Filter a jet collection in-place to the subset that passes the supplied Cut.
Definition: JetUtils.hh:157
double E() const
Get energy (time component of momentum).
Definition: Vector4.hh:543
PdgIdPair pids(const ParticlePair &pp)
Definition: ParticleUtils.hh:749
int pid(const Particle &p)
Unbound function access to PID code.
Definition: ParticleUtils.hh:23
const Cut & open()
Fully open cut singleton, accepts everything.
Definition: MC_Cent_pPb.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