5 #include "Rivet/Rivet.hh"
39 size_t size()
const {
return _particles.size(); }
78 vector<Particle>&
particles() {
return _particles; }
81 const vector<Particle>&
particles()
const {
return _particles; }
164 ParticleVector _particles;
271 inline double deltaR(
const Jet& j1,
const Jet& j2,
273 return deltaR(j1.momentum(), j2.momentum(), scheme);
276 inline double deltaR(
const Jet& j,
const Particle& p,
278 return deltaR(j.momentum(), p.momentum(), scheme);
281 inline double deltaR(
const Particle& p,
const Jet& j,
283 return deltaR(p.momentum(), j.momentum(), scheme);
286 inline double deltaR(
const Jet& j,
const FourMomentum& v,
288 return deltaR(j.momentum(), v, scheme);
291 inline double deltaR(
const Jet& j,
const FourVector& v,
293 return deltaR(j.momentum(), v, scheme);
296 inline double deltaR(
const Jet& j,
const Vector3& v) {
297 return deltaR(j.momentum(), v);
300 inline double deltaR(
const Jet& j,
double eta,
double phi) {
301 return deltaR(j.momentum(),
eta,
phi);
304 inline double deltaR(
const FourMomentum& v,
const Jet& j,
306 return deltaR(v, j.momentum(), scheme);
309 inline double deltaR(
const FourVector& v,
const Jet& j,
311 return deltaR(v, j.momentum(), scheme);
314 inline double deltaR(
const Vector3& v,
const Jet& j) {
315 return deltaR(v, j.momentum());
318 inline double deltaR(
double eta,
double phi,
const Jet& j) {
319 return deltaR(eta, phi, j.momentum());
323 inline double deltaPhi(
const Jet& j1,
const Jet& j2) {
324 return deltaPhi(j1.momentum(), j2.momentum());
327 inline double deltaPhi(
const Jet& j,
const Particle& p) {
328 return deltaPhi(j.momentum(), p.momentum());
331 inline double deltaPhi(
const Particle& p,
const Jet& j) {
332 return deltaPhi(p.momentum(), j.momentum());
335 inline double deltaPhi(
const Jet& j,
const FourMomentum& v) {
336 return deltaPhi(j.momentum(), v);
339 inline double deltaPhi(
const Jet& j,
const FourVector& v) {
340 return deltaPhi(j.momentum(), v);
343 inline double deltaPhi(
const Jet& j,
const Vector3& v) {
344 return deltaPhi(j.momentum(), v);
347 inline double deltaPhi(
const Jet& j,
double phi) {
348 return deltaPhi(j.momentum(),
phi);
351 inline double deltaPhi(
const FourMomentum& v,
const Jet& j) {
352 return deltaPhi(v, j.momentum());
355 inline double deltaPhi(
const FourVector& v,
const Jet& j) {
356 return deltaPhi(v, j.momentum());
359 inline double deltaPhi(
const Vector3& v,
const Jet& j) {
360 return deltaPhi(v, j.momentum());
363 inline double deltaPhi(
double phi,
const Jet& j) {
364 return deltaPhi(phi, j.momentum());
368 inline double deltaEta(
const Jet& j1,
const Jet& j2) {
369 return deltaEta(j1.momentum(), j2.momentum());
372 inline double deltaEta(
const Jet& j,
const Particle& p) {
373 return deltaEta(j.momentum(), p.momentum());
376 inline double deltaEta(
const Particle& p,
const Jet& j) {
377 return deltaEta(p.momentum(), j.momentum());
380 inline double deltaEta(
const Jet& j,
const FourMomentum& v) {
381 return deltaEta(j.momentum(), v);
384 inline double deltaEta(
const Jet& j,
const FourVector& v) {
385 return deltaEta(j.momentum(), v);
388 inline double deltaEta(
const Jet& j,
const Vector3& v) {
389 return deltaEta(j.momentum(), v);
392 inline double deltaEta(
const Jet& j,
double eta) {
393 return deltaEta(j.momentum(),
eta);
396 inline double deltaEta(
const FourMomentum& v,
const Jet& j) {
397 return deltaEta(v, j.momentum());
400 inline double deltaEta(
const FourVector& v,
const Jet& j) {
401 return deltaEta(v, j.momentum());
404 inline double deltaEta(
const Vector3& v,
const Jet& j) {
405 return deltaEta(v, j.momentum());
408 inline double deltaEta(
double eta,
const Jet& j) {
409 return deltaEta(eta, j.momentum());
bool cmpJetsByEt(const Jet &a, const Jet &b)
Use this so that highest is at the front of the list.
Definition: Jet.hh:205
Definition: MC_JetAnalysis.hh:9
bool cmpJetsByAscRapidity(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:254
bool cmpJetsByAscE(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:221
bool cmpJetsByE(const Jet &a, const Jet &b)
Compare jets by (descending - usual sorting for HEP) Use this so that highest is at the front of th...
Definition: Jet.hh:216
bool cmpJetsByP(const Jet &a, const Jet &b)
Compare jets by descending momentum, .
Definition: Jet.hh:195
double E() const
Get energy (time component of momentum).
Definition: Vector4.hh:355
double eta() const
Synonym for pseudorapidity.
Definition: Vector4.hh:134
bool cmpJetsByDescAbsRapidity(const Jet &a, const Jet &b)
Compare jets by (descending) Use this so that highest is at the front of the list.
Definition: Jet.hh:260
double mod() const
Calculate the modulus of a vector. .
Definition: VectorN.hh:87
const FourMomentum & momentum() const
Get equivalent single momentum four-vector.
Definition: Jet.hh:105
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition: Vector4.hh:139
double phi(const Vector3 &v, const PhiMapping mapping=ZERO_2PI)
Synonym for azimuthalAngle.
Definition: Vector3.hh:281
std::vector< Jet > Jets
Typedef for a collection of Jet objects.
Definition: Jet.hh:177
vector< Particle > & particles()
Get the particles in this jet.
Definition: Jet.hh:78
Jet(const vector< Particle > &particles, const FourMomentum &pjet)
Set all the jet data, with full particle information.
Definition: Jet.hh:21
double phi(const PhiMapping mapping=ZERO_2PI) const
Synonym for azimuthalAngle.
Definition: Vector4.hh:114
Representation of particles from a HepMC::GenEvent.
Definition: Particle.hh:16
bool containsBottom() const
Check whether this jet contains a bottom-flavoured hadron (or decay products from one)...
Definition: Jet.cc:129
double totalEnergy() const
Get the total energy of this jet.
Definition: Jet.hh:114
Jet & setState(const vector< Particle > &particles, const FourMomentum &pjet)
Set all the jet data, with full particle information.
Definition: Jet.cc:10
bool cmpJetsByEtDesc(const Jet &a, const Jet &b)
Use this so that lowest is at the front of the list.
Definition: Jet.hh:210
bool cmpJetsByAscAbsRapidity(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:265
double pT() const
Calculate the transverse momentum .
Definition: Vector4.hh:410
bool cmpJetsByAscPt(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:190
Jet & clear()
Reset this jet as empty.
Definition: Jet.cc:146
bool cmpJetsByAscP(const Jet &a, const Jet &b)
Compare jets by ascending momentum, .
Definition: Jet.hh:199
double eta() const
Get the unweighted average for this jet. (caches)
Definition: Jet.hh:108
size_t size() const
Number of particles in this jet.
Definition: Jet.hh:39
bool cmpJetsByPt(const Jet &a, const Jet &b)
Compare jets by (descending - usual sorting for HEP) Use this so that highest is at the front of th...
Definition: Jet.hh:185
Jet & setMomentum(const FourMomentum &momentum)
Set the effective 4-momentum of the jet.
Definition: Jet.cc:24
double EtSum() const
Get the sum of the values of the constituent tracks. (caches)
Definition: Jet.hh:126
const vector< Particle > & particles() const
Get the particles in this jet (const version)
Definition: Jet.hh:81
double hadronicEnergy() const
Get the energy carried in this jet by hadrons.
Definition: Jet.cc:100
bool cmpJetsByDescRapidity(const Jet &a, const Jet &b)
Compare jets by (descending) Use this so that highest is at the front of the list.
Definition: Jet.hh:249
double neutralEnergy() const
Get the energy carried in this jet by neutral particles.
Definition: Jet.cc:88
bool cmpJetsByDescPseudorapidity(const Jet &a, const Jet &b)
Compare jets by (descending) Use this so that highest is at the front of the list.
Definition: Jet.hh:227
Base class for particle-like things like Particle and Jet.
Definition: ParticleBase.hh:11
bool containsParticleId(PdgId pid) const
Check whether this jet contains a certain particle type.
Definition: Jet.cc:67
double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components)
Definition: Vector4.hh:129
RapScheme
Enum for rapidity variable to be used in calculating , applying rapidity cuts, etc.
Definition: MathHeader.hh:60
bool containsCharm() const
Check whether this jet contains a charm-flavoured hadron (or decay products from one).
Definition: Jet.cc:112
double Et() const
Calculate the transverse energy .
Definition: Vector4.hh:420
Representation of a clustered jet of particles.
Definition: Jet.hh:12
double eta(const Vector3 &v)
Synonym for pseudorapidity.
Definition: Vector3.hh:299
Jet & setParticles(const vector< Particle > &particles)
Set the particles collection with full particle information.
Definition: Jet.cc:30
double ptSum() const
Get the sum of the values of the constituent tracks. (caches)
Definition: Jet.hh:123
bool cmpJetsByAscAbsPseudorapidity(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:243
double rapidity() const
Calculate the rapidity.
Definition: Vector4.hh:400
bool cmpJetsByAscPseudorapidity(const Jet &a, const Jet &b)
Compare jets by (ascending) Use this so that lowest is at the front of the list.
Definition: Jet.hh:232
Specialized version of the FourVector with momentum/energy functionality.
Definition: Vector4.hh:324
bool containsParticle(const Particle &particle) const
Check whether this jet contains a particular particle.
Definition: Jet.cc:58
bool cmpJetsByDescAbsPseudorapidity(const Jet &a, const Jet &b)
Compare jets by (descending) Use this so that highest is at the front of the list.
Definition: Jet.hh:238
double phi() const
Get the unweighted average for this jet. (caches)
Definition: Jet.hh:111