rivet is hosted by Hepforge, IPPP Durham
Rivet  2.7.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 
14 namespace Rivet {
15 
16 
18  class Jet : public ParticleBase {
19  public:
20 
22 
23 
25  Jet(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles()) {
26  setState(pj, particles, tags);
27  }
28 
30  Jet(const FourMomentum& pjet, const Particles& particles=Particles(), const Particles& tags=Particles()) {
31  setState(pjet, particles, tags);
32  }
33 
35  Jet() { clear(); }
36 
38 
39 
41 
42 
44  size_t size() const { return _particles.size(); }
45 
47  Particles& particles() { return _particles; }
49  const Particles& particles() const { return _particles; }
51  const Particles particles(const Cut& c) const { return filter_select(_particles, c); }
53  const Particles particles(const ParticleSelector& s) const { return filter_select(_particles, s); }
54 
56  Particles& constituents() { return particles(); }
58  const Particles& constituents() const { return particles(); }
60  const Particles constituents(const Cut& c) const { return particles(c); }
62  const Particles constituents(const ParticleSelector& s) const { return particles(s); }
63 
65  bool containsParticle(const Particle& particle) const;
67  bool containsPID(const Particle& particle) const { return containsParticle(particle); }
68 
70  bool containsParticleId(PdgId pid) const;
72  bool containsPID(PdgId pid) const { return containsParticleId(pid); }
73 
75  bool containsParticleId(const vector<PdgId>& pids) const;
77  bool containsPID(const vector<PdgId>& pids) const { return containsParticleId(pids); }
78 
80 
81 
86 
87 
89  Particles& tags() { return _tags; }
91  const Particles& tags() const { return _tags; }
95  Particles tags(const ParticleSelector& f) const { return filter_select(tags(), f); }
99  Particles tags(const Cut& c) const;
100 
101 
105  Particles bTags(const Cut& c=Cuts::open()) const;
107  Particles bTags(const ParticleSelector& f) const { return filter_select(bTags(), f); }
108 
110  bool bTagged(const Cut& c=Cuts::open()) const { return !bTags(c).empty(); }
112  bool bTagged(const ParticleSelector& f) const { return !bTags(f).empty(); }
113 
114 
118  Particles cTags(const Cut& c=Cuts::open()) const;
120  Particles cTags(const ParticleSelector& f) const { return filter_select(cTags(), f); }
121 
123  bool cTagged(const Cut& c=Cuts::open()) const { return !cTags(c).empty(); }
125  bool cTagged(const ParticleSelector& f) const { return !cTags(f).empty(); }
126 
127 
131  Particles tauTags(const Cut& c=Cuts::open()) const;
133  Particles tauTags(const ParticleSelector& f) const { return filter_select(tauTags(), f); }
134 
136  bool tauTagged(const Cut& c=Cuts::open()) const { return !tauTags(c).empty(); }
138  bool tauTagged(const ParticleSelector& f) const { return !tauTags(f).empty(); }
139 
140 
153  DEPRECATED("Prefer the bTags() or bTagged() function")
154  bool containsBottom(bool include_decay_products=true) const;
155 
168  DEPRECATED("Prefer the cTags() or cTagged() function")
169  bool containsCharm(bool include_decay_products=true) const;
170 
172 
173 
175 
176 
178  const FourMomentum& momentum() const { return _momentum; }
179 
183  Jet& transformBy(const LorentzTransform& lt);
184 
186  double totalEnergy() const { return momentum().E(); }
187 
189  double neutralEnergy() const;
190 
192  double hadronicEnergy() const;
193 
195 
196 
198 
199 
201  const fastjet::PseudoJet& pseudojet() const { return _pseudojet; }
202 
204  operator const fastjet::PseudoJet& () const { return pseudojet(); }
205 
207 
208 
210 
211 
216  Jet& setState(const fastjet::PseudoJet& pj, const Particles& particles=Particles(), const Particles& tags=Particles());
217 
219  Jet& setState(const FourMomentum& mom, const Particles& particles, const Particles& tags=Particles());
220 
224  Jet& setParticles(const Particles& particles);
225  Jet& setConstituents(const Particles& particles) { return setParticles(particles); }
226 
228  Jet& clear();
229 
231 
232 
233  private:
234 
236  fastjet::PseudoJet _pseudojet;
237 
240  Particles _particles;
241 
243  Particles _tags;
244 
246  mutable FourMomentum _momentum;
247 
248  };
249 
250 
252 
253 
255  std::ostream& operator << (std::ostream& os, const Jet& j);
256 
258 
259 
260 }
261 
262 
263 #include "Rivet/Tools/JetUtils.hh"
264 
265 #endif
Definition: ALICE_2010_I880049.cc:13
Jet & setParticles(const Particles &particles)
Set the particles collection with full particle information.
Definition: Jet.cc:55
size_t size() const
Number of particles in this jet.
Definition: Jet.hh:44
bool containsPID(PdgId pid) const
Nicer alias for containsParticleId.
Definition: Jet.hh:72
Jet & transformBy(const LorentzTransform &lt)
Definition: Jet.cc:91
Particles tags(const ParticleSelector &f) const
Particles which have been tag-matched to this jet and pass a selector function.
Definition: Jet.hh:95
bool bTagged(const ParticleSelector &f) const
Does this jet have at least one b-tag (that passes the supplied selector function)?
Definition: Jet.hh:112
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:110
const FourMomentum & mom() const
Get equivalent single momentum four-vector (const) (alias).
Definition: ParticleBase.hh:39
PdgIdPair pids(const ParticlePair &pp)
Definition: ParticleUtils.hh:735
DEPRECATED("Prefer the bTags() or bTagged() function") bool containsBottom(bool include_decay_products
Check whether this jet contains a bottom-flavoured hadron.
bool containsParticleId(PdgId pid) const
Check whether this jet contains a certain particle type.
Definition: Jet.cc:70
double neutralEnergy() const
Get the energy carried in this jet by neutral particles.
Definition: Jet.cc:100
Particles & particles()
Get the particles in this jet.
Definition: Jet.hh:47
Particle representation, either from a HepMC::GenEvent or reconstructed.
Definition: Particle.hh:18
int pid(const Particle &p)
Unbound function access to PID code.
Definition: ParticleUtils.hh:20
Jet(const FourMomentum &pjet, const Particles &particles=Particles(), const Particles &tags=Particles())
Set the jet data, with optional full particle information.
Definition: Jet.hh:30
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:136
Jet & clear()
Reset this jet as empty.
Definition: Jet.cc:10
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:143
const Particles constituents(const Cut &c) const
Get the particles in this jet which pass a cut (FastJet-like alias, const)
Definition: Jet.hh:60
bool containsPID(const vector< PdgId > &pids) const
Nicer alias for containsParticleId.
Definition: Jet.hh:77
const fastjet::PseudoJet & pseudojet() const
Access the internal FastJet3 PseudoJet (as a const reference)
Definition: Jet.hh:201
bool cTagged(const ParticleSelector &f) const
Does this jet have at least one c-tag (that passes the supplied selector function)?
Definition: Jet.hh:125
Particles bTags(const Cut &c=Cuts::open()) const
b particles which have been tag-matched to this jet (and pass an optional Cut)
Definition: Jet.cc:166
bool tauTagged(const ParticleSelector &f) const
Does this jet have at least one tau-tag (that passes the supplied selector function)?
Definition: Jet.hh:138
Particles tauTags(const Cut &c=Cuts::open()) const
Tau particles which have been tag-matched to this jet (and pass an optional Cut)
Definition: Jet.cc:183
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) ...
Definition: Jet.cc:174
const FourMomentum & momentum() const
Get equivalent single momentum four-vector.
Definition: Jet.hh:178
const Cut & open()
Fully open cut singleton, accepts everything.
Definition: Cuts.cc:81
Particles tauTags(const ParticleSelector &f) const
Tau particles which have been tag-matched to this jet and pass a selector function.
Definition: Jet.hh:133
Object implementing Lorentz transform calculations and boosts.
Definition: LorentzTrans.hh:21
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:120
bool containsParticle(const Particle &particle) const
Check whether this jet contains a particular particle.
Definition: Jet.cc:61
Base class for particle-like things like Particle and Jet.
Definition: ParticleBase.hh:13
const Particles & constituents() const
Get the particles in this jet (FastJet-like alias, const version)
Definition: Jet.hh:58
double E() const
Get energy (time component of momentum).
Definition: Vector4.hh:538
Particles & tags()
Particles which have been tag-matched to this jet.
Definition: Jet.hh:89
const Particles particles(const Cut &c) const
Get the particles in this jet which pass a cut (const)
Definition: Jet.hh:51
const Particles particles(const ParticleSelector &s) const
Get the particles in this jet which pass a filtering functor (const)
Definition: Jet.hh:53
Representation of a clustered jet of particles.
Definition: Jet.hh:18
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:62
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...
Definition: Jet.cc:28
const Particles & tags() const
Particles which have been tag-matched to this jet (const version)
Definition: Jet.hh:91
const Particles & particles() const
Get the particles in this jet (const version)
Definition: Jet.hh:49
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:301
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:123
double hadronicEnergy() const
Get the energy carried in this jet by hadrons.
Definition: Jet.cc:112
double totalEnergy() const
Get the total energy of this jet.
Definition: Jet.hh:186
Jet()
Default constructor – only for STL storability.
Definition: Jet.hh:35
bool containsPID(const Particle &particle) const
Nicer alias for containsParticleId.
Definition: Jet.hh:67
Particles & constituents()
Get the particles in this jet (FastJet-like alias)
Definition: Jet.hh:56
Particles bTags(const ParticleSelector &f) const
b particles which have been tag-matched to this jet and pass a selector function
Definition: Jet.hh:107
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:25