rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0
Rivet::Jet Class Reference

Representation of a clustered jet of particles. More...

#include <Jet.hh>

Inheritance diagram for Rivet::Jet:
Rivet::ParticleBase

Public Member Functions

Constructors
 Jet (const fastjet::PseudoJet &pj, const Particles &particles=Particles(), const Particles &tags=Particles())
 Constructor from a FastJet PseudoJet, with optional full particle constituents information.
 
 Jet (const FourMomentum &pjet, const Particles &particles=Particles(), const Particles &tags=Particles())
 Set the jet data, with optional full particle information.
 
 Jet ()
 Default constructor – only for STL storability.
 
Access jet constituents
size_t size () const
 Number of particles in this jet.
 
Particlesparticles ()
 Get the particles in this jet.
 
const Particlesparticles () const
 Get the particles in this jet (const version)
 
const Particles particles (const Cut &c) const
 Get the particles in this jet which pass a cut (const)
 
const Particles particles (const ParticleSelector &s) const
 Get the particles in this jet which pass a filtering functor (const)
 
Particlesconstituents ()
 Get the particles in this jet (FastJet-like alias)
 
const Particlesconstituents () const
 Get the particles in this jet (FastJet-like alias, const version)
 
const Particles constituents (const Cut &c) const
 Get the particles in this jet which pass a cut (FastJet-like alias, const)
 
const Particles constituents (const ParticleSelector &s) const
 Get the particles in this jet which pass a filtering functor (FastJet-like alias, const)
 
bool containsParticle (const Particle &particle) const
 Check whether this jet contains a particular particle.
 
bool containsPID (const Particle &particle) const
 Nicer alias for containsParticleId.
 
bool containsParticleId (PdgId pid) const
 Check whether this jet contains a certain particle type.
 
bool containsPID (PdgId pid) const
 Nicer alias for containsParticleId.
 
bool containsParticleId (const vector< PdgId > &pids) const
 Check whether this jet contains at least one of certain particle types.
 
bool containsPID (const vector< PdgId > &pids) const
 Nicer alias for containsParticleId.
 
Tagging
Note
General sources of tag particles are planned. The default jet finding adds b-hadron, c-hadron, and tau tags by ghost association.
Particlestags ()
 Particles which have been tag-matched to this jet.
 
const Particlestags () const
 Particles which have been tag-matched to this jet (const version)
 
Particles tags (const ParticleSelector &f) const
 Particles which have been tag-matched to this jet and pass a selector function.
 
Particles tags (const Cut &c) const
 Particles which have been tag-matched to this jet and pass a Cut.
 
Particles bTags (const Cut &c=Cuts::open()) const
 b particles which have been tag-matched to this jet (and pass an optional Cut)
 
Particles bTags (const ParticleSelector &f) const
 b particles which have been tag-matched to this jet and pass a selector function
 
bool bTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one b-tag (that passes an optional Cut)?
 
bool bTagged (const ParticleSelector &f) const
 Does this jet have at least one b-tag (that passes the supplied selector function)?
 
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)
 
Particles cTags (const ParticleSelector &f) const
 c (and not b) particles which have been tag-matched to this jet and pass a selector function
 
bool cTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one c-tag (that passes an optional Cut)?
 
bool cTagged (const ParticleSelector &f) const
 Does this jet have at least one c-tag (that passes the supplied selector function)?
 
Particles tauTags (const Cut &c=Cuts::open()) const
 Tau particles which have been tag-matched to this jet (and pass an optional Cut)
 
Particles tauTags (const ParticleSelector &f) const
 Tau particles which have been tag-matched to this jet and pass a selector function.
 
bool tauTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one tau-tag (that passes an optional Cut)?
 
bool tauTagged (const ParticleSelector &f) const
 Does this jet have at least one tau-tag (that passes the supplied selector function)?
 
Effective jet 4-vector properties
const FourMomentummomentum () const
 Get equivalent single momentum four-vector.
 
JettransformBy (const LorentzTransform &lt)
 
double totalEnergy () const
 Get the total energy of this jet.
 
double neutralEnergy () const
 Get the energy carried in this jet by neutral particles.
 
double hadronicEnergy () const
 Get the energy carried in this jet by hadrons.
 
Interaction with FastJet
const fastjet::PseudoJet & pseudojet () const
 Access the internal FastJet3 PseudoJet (as a const reference)
 
 operator const fastjet::PseudoJet & () const
 Cast operator to FastJet3 PseudoJet (as a const reference)
 
Set the jet constituents and properties
JetsetState (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.
 
JetsetState (const FourMomentum &mom, const Particles &particles, const Particles &tags=Particles())
 Set all the jet data, with optional full particle constituent and tag information.
 
JetsetParticles (const Particles &particles)
 Set the particles collection with full particle information.
 
JetsetConstituents (const Particles &particles)
 
Jetclear ()
 Reset this jet as empty.
 
Effective momentum accessors
const FourMomentummom () const
 Get equivalent single momentum four-vector (const) (alias).
 
 operator const FourMomentum & () const
 Cast operator for conversion to FourMomentum.
 
Convenience access to the effective 4-vector properties
double E () const
 Get the energy directly.
 
double energy () const
 Get the energy directly (alias).
 
double E2 () const
 Get the energy-squared.
 
double energy2 () const
 Get the energy-squared (alias).
 
double pt () const
 Get the \( p_T \) directly.
 
double pT () const
 Get the \( p_T \) directly (alias).
 
double perp () const
 Get the \( p_T \) directly (alias).
 
double pt2 () const
 Get the \( p_T^2 \) directly.
 
double pT2 () const
 Get the \( p_T^2 \) directly (alias).
 
double perp2 () const
 Get the \( p_T^2 \) directly (alias).
 
double Et () const
 Get the \( E_T \) directly.
 
double Et2 () const
 Get the \( E_T^2 \) directly.
 
double mass () const
 Get the mass directly.
 
double mass2 () const
 Get the mass**2 directly.
 
double pseudorapidity () const
 Get the \( \eta \) directly.
 
double eta () const
 Get the \( \eta \) directly (alias).
 
double abspseudorapidity () const
 Get the \( |\eta| \) directly.
 
double abseta () const
 Get the \( |\eta| \) directly (alias).
 
double rapidity () const
 Get the \( y \) directly.
 
double rap () const
 Get the \( y \) directly (alias).
 
double absrapidity () const
 Get the \( |y| \) directly.
 
double absrap () const
 Get the \( |y| \) directly (alias).
 
double azimuthalAngle (const PhiMapping mapping=ZERO_2PI) const
 Azimuthal angle \( \phi \).
 
double phi (const PhiMapping mapping=ZERO_2PI) const
 Get the \( \phi \) directly.
 
ThreeMomentum p3 () const
 Get the 3-momentum directly.
 
double p () const
 Get the 3-momentum magnitude directly.
 
double p2 () const
 Get the 3-momentum magnitude-squared directly.
 
Vector3 ptvec () const
 Get the transverse 3-momentum directly.
 
Vector3 pTvec () const
 Get the transverse 3-momentum directly.
 
double px () const
 x component of momentum.
 
double py () const
 y component of momentum.
 
double pz () const
 z component of momentum.
 
double px2 () const
 x component of momentum, squared.
 
double py2 () const
 y component of momentum, squared.
 
double pz2 () const
 z component of momentum, squared.
 
double polarAngle () const
 Angle subtended by the 3-vector and the z-axis.
 
double theta () const
 Synonym for polarAngle.
 
double angle (const ParticleBase &v) const
 Angle between this vector and another.
 
double angle (const FourVector &v) const
 Angle between this vector and another.
 
double angle (const Vector3 &v3) const
 Angle between this vector and another (3-vector)
 
double dot (const ParticleBase &v) const
 Lorentz dot product between this 4-vector and another.
 
double dot (const FourVector &v) const
 Angle between this 4-vector and another.
 

Detailed Description

Representation of a clustered jet of particles.

Member Function Documentation

◆ bTags()

Particles Rivet::Jet::bTags ( const Cut c = Cuts::open()) const

b particles which have been tag-matched to this jet (and pass an optional Cut)

The default jet finding adds b-hadron tags by ghost association.

Referenced by bTagged(), and bTagged().

◆ cTags()

Particles Rivet::Jet::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)

The default jet finding adds c-hadron tags by ghost association.

Referenced by cTagged(), and cTagged().

◆ momentum()

const FourMomentum & Rivet::Jet::momentum ( ) const
inlinevirtual

Get equivalent single momentum four-vector.

Implements Rivet::ParticleBase.

Referenced by totalEnergy().

◆ setParticles()

Jet & Rivet::Jet::setParticles ( const Particles particles)

Set the particles collection with full particle information.

If set, this overrides particle info extracted from the PseudoJet

◆ setState()

Jet & Rivet::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.

Note
The particles() list will be extracted from PseudoJet constituents by default, making use of an attached user info if one is found.

Referenced by Jet(), and Jet().

◆ tags() [1/2]

Particles Rivet::Jet::tags ( const Cut c) const

Particles which have been tag-matched to this jet and pass a Cut.

Note
Note the less efficient return by value, due to the cut-pass filtering.

◆ tags() [2/2]

Particles Rivet::Jet::tags ( const ParticleSelector f) const
inline

Particles which have been tag-matched to this jet and pass a selector function.

Note
Note the less efficient return by value, due to the filtering.

References Rivet::select(), and tags().

Referenced by tags().

◆ tauTags()

Particles Rivet::Jet::tauTags ( const Cut c = Cuts::open()) const

Tau particles which have been tag-matched to this jet (and pass an optional Cut)

The default jet finding adds tau tags by ghost association.

Referenced by Rivet::TAUJET_EFF_ATLAS_RUN1(), Rivet::TAUJET_EFF_ATLAS_RUN2(), tauTagged(), and tauTagged().

◆ transformBy()

Jet & Rivet::Jet::transformBy ( const LorentzTransform lt)

Apply an active Lorentz transform to this jet

Note
The Rivet jet momentum, constituent particles, and tag particles will be modified.
Warning
The FastJet cluster sequence and pseudojets will not be modified: don't use them after transformation!

The documentation for this class was generated from the following file:
  • /Users/chrisg/software/rivet/include/Rivet/Jet.hh