rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0
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) const { return select(tags(), f); }
123 Particles tags(const Cut& c) const;
124
125
129 Particles bTags(const Cut& c=Cuts::open()) const;
131 Particles bTags(const ParticleSelector& f) const { return select(bTags(), f); }
132
134 bool bTagged(const Cut& c=Cuts::open()) const { return !bTags(c).empty(); }
136 bool bTagged(const ParticleSelector& f) const { return !bTags(f).empty(); }
137
138
142 Particles cTags(const Cut& c=Cuts::open()) const;
144 Particles cTags(const ParticleSelector& f) const { return select(cTags(), f); }
145
147 bool cTagged(const Cut& c=Cuts::open()) const { return !cTags(c).empty(); }
149 bool cTagged(const ParticleSelector& f) const { return !cTags(f).empty(); }
150
151
157 Particles tauTags(const ParticleSelector& f) const { return select(tauTags(), f); }
158
160 bool tauTagged(const Cut& c=Cuts::open()) const { return !tauTags(c).empty(); }
162 bool tauTagged(const ParticleSelector& f) const { return !tauTags(f).empty(); }
163
165
166
169
171 const FourMomentum& momentum() const { return _momentum; }
172
177
179 double totalEnergy() const { return momentum().E(); }
180
182 double neutralEnergy() const;
183
185 double hadronicEnergy() const;
186
188
189
192
194 const fastjet::PseudoJet& pseudojet() const { return _pseudojet; }
195
197 operator const fastjet::PseudoJet& () const { return pseudojet(); }
198
200
201
204
209 Jet& setState(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles());
210
213
218 Jet& setConstituents(const Particles& particles) { return setParticles(particles); }
219
222
224
225
226 private:
227
229 fastjet::PseudoJet _pseudojet;
230
233 Particles _particles;
234
236 Particles _tags;
237
239 mutable FourMomentum _momentum;
240
241 // /// Provide but hide the equality operators, to avoid implicit comparison via fastjet::PseudoJet
242 // bool operator == (const Jet&) const;
243 // bool operator != (const Jet&) const;
244
245 };
246
247
250
252 std::ostream& operator << (std::ostream& os, const Jet& j);
253
255
256
257}
258
259
260#include "Rivet/Tools/JetUtils.hh"
261
262#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
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:144
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:194
const Particles particles(const ParticleSelector &s) const
Get the particles in this jet which pass a filtering functor (const)
Definition Jet.hh:77
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: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 tauTagged(const ParticleSelector &f) const
Does this jet have at least one tau-tag (that passes the supplied selector function)?
Definition Jet.hh:162
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:131
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 cTagged(const Cut &c=Cuts::open()) const
Does this jet have at least one c-tag (that passes an optional Cut)?
Definition Jet.hh:147
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
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:134
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
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.
bool bTagged(const ParticleSelector &f) const
Does this jet have at least one b-tag (that passes the supplied selector function)?
Definition Jet.hh:136
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:171
Jet()
Default constructor – only for STL storability.
Definition Jet.hh:59
Particles tags(const ParticleSelector &f) const
Particles which have been tag-matched to this jet and pass a selector function.
Definition Jet.hh:119
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
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:160
Particles tauTags(const ParticleSelector &f) const
Tau particles which have been tag-matched to this jet and pass a selector function.
Definition Jet.hh:157
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 cTagged(const ParticleSelector &f) const
Does this jet have at least one c-tag (that passes the supplied selector function)?
Definition Jet.hh:149
double totalEnergy() const
Get the total energy of this jet.
Definition Jet.hh:179
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:101
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: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