rivet is hosted by Hepforge, IPPP Durham
Jet Class Reference

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

#include <Jet.hh>

Inheritance diagram for Jet:
Collaboration diagram for Jet:

List of all members.

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 (const Particles &particles, const FourMomentum &pjet)
 Jet ()
 Default constructor -- only for STL storability.
Access jet constituents
size_t size () const
 Number of particles in this jet.
vector< Particle > & particles ()
 Get the particles in this jet.
const vector< Particle > & particles () const
 Get the particles in this jet (const version)
vector< Particle > & constituents ()
 Get the particles in this jet (FastJet-like alias)
const vector< Particle > & constituents () const
 Get the particles in this jet (FastJet-like alias, const version)
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 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)
bool bTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one b-tag (that passes an optional Cut)?
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)
bool cTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one c-tag (that passes an optional Cut)?
Particles tauTags (const Cut &c=Cuts::open()) const
 Tau particles which have been tag-matched to this jet (and pass an optional Cut)
bool tauTagged (const Cut &c=Cuts::open()) const
 Does this jet have at least one tau-tag (that passes an optional Cut)?
bool containsBottom (bool include_decay_products=true) const
 Check whether this jet contains a bottom-flavoured hadron.
bool containsCharm (bool include_decay_products=true) const
 Check whether this jet contains a charm-flavoured hadron.
Access additional effective jet 4-vector properties
const FourMomentummomentum () const
 Get equivalent single momentum four-vector.
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.
JetsetState (const Particles &particles, const FourMomentum &mom)
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 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 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.
Vector3 p3 () const
 Get the 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 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)

Private Attributes

fastjet::PseudoJet _pseudojet
 FJ3 PseudoJet member to unify PseudoJet and Jet.
Particles _particles
Particles _tags
 Particles used to tag this jet (can be anything, but c and b hadrons are the most common)
FourMomentum _momentum
 Effective jet 4-vector (just for caching)

Detailed Description

Representation of a clustered jet of particles.

Definition at line 17 of file Jet.hh.


Constructor & Destructor Documentation

Jet ( const fastjet::PseudoJet &  pj,
const Particles particles = Particles(),
const Particles tags = Particles() 
) [inline]

Constructor from a FastJet PseudoJet, with optional full particle constituents information.

Definition at line 24 of file Jet.hh.

References Jet::particles(), Jet::setState(), and Jet::tags().

                                                                                                               {
      setState(pj, particles, tags);
    }
Jet ( const FourMomentum pjet,
const Particles particles = Particles(),
const Particles tags = Particles() 
) [inline]

Set the jet data, with optional full particle information.

Definition at line 29 of file Jet.hh.

References Jet::particles(), Jet::setState(), and Jet::tags().

                                                                                                             {
      setState(pjet, particles, tags);
    }
Jet ( const Particles particles,
const FourMomentum pjet 
) [inline]

Set all the jet data, with full particle information.

Deprecated:
Prefer the form where the 4-vec comes first and the particles list is optional.

Definition at line 36 of file Jet.hh.

References Jet::setState().

                                                              {
      setState(pjet, particles);
    }
Jet ( ) [inline]

Default constructor -- only for STL storability.

Definition at line 41 of file Jet.hh.

References Jet::clear().

{ clear(); }

Member Function Documentation

double abseta ( ) const [inline, inherited]

Get the $ |\eta| $ directly (alias).

Definition at line 73 of file ParticleBase.hh.

References FourVector::abseta(), and ParticleBase::momentum().

Referenced by CMS_2012_I1087342::analyze(), TOTEM_2012_I1115294::analyze(), TOTEM_2014_I1328627::analyze(), CMS_2011_S9088458::analyze(), CMSTOTEM_2014_I1294140::analyze(), CMS_2011_S8957746::analyze(), UA5_1982_S875503::analyze(), CDF_2001_S4563131::analyze(), CMS_2013_I1273574::analyze(), D0_1996_S3324664::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CDF_2008_S7540469::analyze(), CMS_2011_S9215166::analyze(), MC_IDENTIFIED::analyze(), CDF_1997_S3541940::analyze(), D0_2008_S7863608::analyze(), UA5_1986_S1583476::analyze(), D0_2009_S8202443::analyze(), ATLAS_2013_I1219109::analyze(), CMS_2015_I1385107::analyze(), ATLAS_2012_I1082009::analyze(), CMS_2011_S9120041::analyze(), CMS_2015_I1310737::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1126136::analyze(), D0_2009_S8349509::analyze(), STAR_2009_UE_HELEN::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2012_I1183818::analyze(), ATLAS_2012_I1199269::analyze(), ATLAS_2011_I921594::analyze(), ATLAS_2013_I1217863_Z::analyze(), ATLAS_2013_I1244522::analyze(), ATLAS_2014_I1306294::analyze(), ATLAS_2013_I1217863_W::analyze(), ATLAS_2013_I1263495::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2011_S9128077::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2012_CONF_2012_153::analyze(), ATLAS_2014_I1327229::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2011_S9225137::analyze(), and ATLAS_2012_I1203852::analyze().

{ return momentum().abseta(); }
double abspseudorapidity ( ) const [inline, inherited]

Get the $ |\eta| $ directly.

Definition at line 71 of file ParticleBase.hh.

References FourVector::abspseudorapidity(), and ParticleBase::momentum().

{ return momentum().abspseudorapidity(); }
double absrapidity ( ) const [inline, inherited]

Get the $ |y| $ directly.

Definition at line 80 of file ParticleBase.hh.

References FourMomentum::absrapidity(), and ParticleBase::momentum().

Referenced by CMS_2014_I1298810::analyze().

{ return momentum().absrapidity(); }
double angle ( const ParticleBase v) const [inline, inherited]

Angle between this vector and another.

Definition at line 108 of file ParticleBase.hh.

References FourVector::angle(), and ParticleBase::momentum().

Referenced by H1_2000_S4129130::analyze(), and H1_1994_S2919893::analyze().

{ return momentum().angle(v.momentum()); }
double angle ( const FourVector v) const [inline, inherited]

Angle between this vector and another.

Definition at line 110 of file ParticleBase.hh.

References FourVector::angle(), and ParticleBase::momentum().

{ return momentum().angle(v); }
double angle ( const Vector3 v3) const [inline, inherited]

Angle between this vector and another (3-vector)

Definition at line 112 of file ParticleBase.hh.

References FourVector::angle(), and ParticleBase::momentum().

{ return momentum().angle(v3); }
bool bTagged ( const Cut c = Cuts::open()) const [inline]
Particles 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.

Definition at line 162 of file Jet.cc.

References Rivet::PID::hasBottom(), and Jet::tags().

Referenced by MC_JETTAGS::analyze(), ATLAS_2014_I1304688::analyze(), ATLAS_2015_I1345452::analyze(), and Jet::bTagged().

                                         {
    Particles rtn;
    foreach (const Particle& tp, tags()) {
      if (hasBottom(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
Jet & clear ( )

Reset this jet as empty.

Definition at line 10 of file Jet.cc.

References Jet::_momentum, Jet::_particles, and Jet::_pseudojet.

Referenced by Jet::Jet(), and Jet::setState().

                  {
    _momentum = FourMomentum();
    _pseudojet.reset(0,0,0,0);
    _particles.clear();
    return *this;
  }
vector<Particle>& constituents ( ) [inline]

Get the particles in this jet (FastJet-like alias)

Definition at line 58 of file Jet.hh.

References Jet::particles().

{ return particles(); }
const vector<Particle>& constituents ( ) const [inline]

Get the particles in this jet (FastJet-like alias, const version)

Definition at line 60 of file Jet.hh.

References Jet::particles().

{ return particles(); }
bool containsBottom ( bool  include_decay_products = true) const

Check whether this jet contains a bottom-flavoured hadron.

Deprecated:
The bTags() or bTagged() function is probably what you want for tagging. This one ignores the tags() list and draws conclusions based directly on the jet constituents; the other gives a much better match to typical experimental methods.
Note:
The decision is made by first trying to find a bottom-flavoured particle in the particles list. Most likely this will fail unless bottom hadrons are set stable. If include_decay_products is true (the default), a fallback is attempted, using the post-hadronization ancestor history of all constituents.

Definition at line 134 of file Jet.cc.

References Rivet::PID::BQUARK, Particle::genParticle(), Rivet::PID::hasBottom(), Rivet::PID::isHadron(), Rivet::particles(), Jet::particles(), Rivet::pi, and Particle::pid().

                                                            {
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (abs(pid) == PID::BQUARK) return true;
      if (PID::isHadron(pid) && PID::hasBottom(pid)) return true;
      if (include_decay_products) {
        const HepMC::GenVertex* gv = p.genParticle()->production_vertex();
        if (gv) {
          foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) {
            const PdgId pid2 = pi->pdg_id();
            if (PID::isHadron(pid2) && PID::hasBottom(pid2)) return true;
          }
        }
      }
    }
    return false;
  }
bool containsCharm ( bool  include_decay_products = true) const

Check whether this jet contains a charm-flavoured hadron.

Deprecated:
The cTags() or cTagged() function is probably what you want for tagging. This one ignores the tags() list and draws conclusions based directly on the jet constituents; the other gives a much better match to typical experimental methods.
Note:
The decision is made by first trying to find a charm-flavoured particle in the particles list. Most likely this will fail unless charmed hadrons are set stable. If include_decay_products is true (the default), a fallback is attempted, using the post-hadronization ancestor history of all constituents.

Definition at line 115 of file Jet.cc.

References Rivet::PID::CQUARK, Particle::genParticle(), Rivet::PID::hasCharm(), Rivet::PID::isHadron(), Rivet::particles(), Jet::particles(), Rivet::pi, and Particle::pid().

                                                           {
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (abs(pid) == PID::CQUARK) return true;
      if (PID::isHadron(pid) && PID::hasCharm(pid)) return true;
      if (include_decay_products) {
        const HepMC::GenVertex* gv = p.genParticle()->production_vertex();
        if (gv) {
          foreach (const GenParticle* pi, Rivet::particles(gv, HepMC::ancestors)) {
            const PdgId pid2 = pi->pdg_id();
            if (PID::isHadron(pid2) && PID::hasCharm(pid2)) return true;
          }
        }
      }
    }
    return false;
  }
bool containsParticle ( const Particle particle) const

Check whether this jet contains a particular particle.

Definition at line 61 of file Jet.cc.

References Particle::genParticle(), and Jet::particles().

Referenced by Jet::containsPID().

                                                           {
    const int barcode = particle.genParticle()->barcode();
    foreach (const Particle& p, particles()) {
      if (p.genParticle()->barcode() == barcode) return true;
    }
    return false;
  }
bool containsParticleId ( PdgId  pid) const

Check whether this jet contains a certain particle type.

Definition at line 70 of file Jet.cc.

References Jet::particles(), and Particle::pid().

Referenced by Jet::containsPID().

                                              {
    foreach (const Particle& p, particles()) {
      if (p.pid() == pid) return true;
    }
    return false;
  }
bool containsParticleId ( const vector< PdgId > &  pids) const

Check whether this jet contains at least one of certain particle types.

Definition at line 78 of file Jet.cc.

References Jet::particles(), and Particle::pid().

                                                              {
    foreach (const Particle& p, particles()) {
      foreach (PdgId pid, pids) {
        if (p.pid() == pid) return true;
      }
    }
    return false;
  }
bool containsPID ( const Particle particle) const [inline]

Nicer alias for containsParticleId.

Definition at line 65 of file Jet.hh.

References Jet::containsParticle().

{ return containsParticle(particle); }
bool containsPID ( PdgId  pid) const [inline]

Nicer alias for containsParticleId.

Definition at line 70 of file Jet.hh.

References Jet::containsParticleId().

{ return containsParticleId(pid); }
bool containsPID ( const vector< PdgId > &  pids) const [inline]

Nicer alias for containsParticleId.

Definition at line 75 of file Jet.hh.

References Jet::containsParticleId().

{ return containsParticleId(pids); }
bool cTagged ( const Cut c = Cuts::open()) const [inline]

Does this jet have at least one c-tag (that passes an optional Cut)?

Definition at line 109 of file Jet.hh.

References Jet::cTags().

Referenced by ATLAS_2012_I1188891::analyze().

{ return !cTags(c).empty(); }
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)

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

Todo:
Is making b and c tags exclusive the right thing to do?

Definition at line 171 of file Jet.cc.

References Rivet::PID::hasBottom(), Rivet::PID::hasCharm(), and Jet::tags().

Referenced by MC_JETTAGS::analyze(), and Jet::cTagged().

                                         {
    Particles rtn;
    foreach (const Particle& tp, tags()) {
      /// @todo Is making b and c tags exclusive the right thing to do?
      if (hasCharm(tp) && !hasBottom(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
double energy ( ) const [inline, inherited]

Get the energy directly (alias).

Definition at line 42 of file ParticleBase.hh.

References FourMomentum::E(), and ParticleBase::momentum().

{ return momentum().E(); }
double eta ( ) const [inline, inherited]

Get the $ \eta $ directly (alias).

Definition at line 69 of file ParticleBase.hh.

References FourVector::eta(), and ParticleBase::momentum().

Referenced by ATLAS_2012_I1118269::analyze(), ATLAS_2011_I894867::analyze(), LHCB_2010_I867355::analyze(), CMS_2012_I1193338::analyze(), TOTEM_2012_002::analyze(), CMS_2012_I1184941::analyze(), CMS_2010_S8656010::analyze(), ATLAS_2010_S8591806::analyze(), CMS_2015_I1384119::analyze(), STAR_2008_S7993412::analyze(), CMSTOTEM_2014_I1294140::analyze(), STAR_2006_S6870392::analyze(), ALICE_2012_I1181770::analyze(), ATLAS_2015_I1387176::analyze(), H1_1994_S2919893::analyze(), CMS_QCD_10_024::analyze(), MC_HINC::analyze(), CDF_2005_S6080774::analyze(), CMS_2010_S8547297::analyze(), D0_1996_S3324664::analyze(), ATLAS_2011_S8994773::analyze(), MC_DIPHOTON::analyze(), CDF_1990_S2089246::analyze(), D0_2008_S6879055::analyze(), CDF_1993_S2742446::analyze(), CDF_2009_S8436959::analyze(), MC_ZINC::analyze(), CDF_2008_S7540469::analyze(), CMS_2011_S8973270::analyze(), CMS_2011_S9215166::analyze(), D0_2010_S8570965::analyze(), SFM_1984_S1178091::analyze(), ALICE_2010_S8625980::analyze(), D0_2006_S6438750::analyze(), MC_DIJET::analyze(), ATLAS_2012_I1119557::analyze(), ATLAS_2014_I1298811::analyze(), MC_WINC::analyze(), CMS_2011_S8884919::analyze(), ATLAS_2012_I946427::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_CONF_2012_104::analyze(), LHCB_2013_I1208105::analyze(), ATLAS_2014_I1307756::analyze(), ATLAS_2012_I1190891::analyze(), UA1_1990_S2044935::analyze(), ATLAS_2012_I1186556::analyze(), CMS_2013_I1218372::analyze(), ATLAS_2012_I1183818::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2012_I1199269::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2013_I1244522::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2011_S9120807::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2010_S8914702::analyze(), D0_2008_S7719523::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2010_S8894728::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2011_S9225137::analyze(), ATLAS_2012_I1093738::analyze(), MC_VH2BB::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2012_I1093734::analyze(), ATLAS_2014_I1282447::analyze(), ATLAS_2012_I1204447::apply_reco_eff(), ATLAS_2014_I1327229::apply_reco_eff(), CMS_2013_I1261026::eventDecomp(), ATLAS_2012_I1094061::fillHistos(), ATLAS_2012_I1084540::fillMap(), ATLAS_2011_I928289_W::FillPlots1d(), ATLAS_2010_S8918562::fillPtEtaNch(), ATLAS_2012_I1091481::getSeta(), CMS_2013_I1258128::makePhotonCut(), TriggerCDFRun2::project(), TriggerCDFRun0Run1::project(), FinalState::project(), NeutralFinalState::project(), and TriggerUA5::project().

{ return momentum().eta(); }
double hadronicEnergy ( ) const

Get the energy carried in this jet by hadrons.

Definition at line 103 of file Jet.cc.

References ParticleBase::E(), Rivet::PID::isHadron(), Jet::particles(), and Particle::pid().

                                   {
    double e_hadr = 0.0;
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (PID::isHadron(pid)) {
        e_hadr += p.E();
      }
    }
    return e_hadr;
  }
double mass2 ( ) const [inline, inherited]

Get the mass**2 directly.

Definition at line 64 of file ParticleBase.hh.

References FourMomentum::mass2(), and ParticleBase::momentum().

{ return momentum().mass2(); }
const FourMomentum& mom ( ) const [inline, inherited]

Get equivalent single momentum four-vector (const) (alias).

Definition at line 28 of file ParticleBase.hh.

References ParticleBase::momentum().

Referenced by CMS_2015_I1346843::analyze(), Hemispheres::calc(), Particle::pseudojet(), Jet::setState(), and Variables::Variables().

{ return momentum(); };
const FourMomentum& momentum ( ) const [inline, virtual]

Get equivalent single momentum four-vector.

Implements ParticleBase.

Definition at line 157 of file Jet.hh.

References Jet::_momentum.

Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), Variables::_getNumGapJets(), CDF_1996_S3349578::_threeJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), STAR_2006_S6870392::analyze(), ATLAS_2012_I1188891::analyze(), D0_1996_S3324664::analyze(), D0_2008_S6879055::analyze(), CDF_1993_S2742446::analyze(), CMS_2013_I1208923::analyze(), ATLAS_2011_I930220::analyze(), CDF_1996_S3108457::analyze(), CDF_2008_S7782535::analyze(), CMS_2013_I1272853::analyze(), CDF_1997_S3541940::analyze(), MC_HFJETS::analyze(), ATLAS_2014_I1268975::analyze(), CMS_2015_I1385107::analyze(), ATLAS_2010_S8817804::analyze(), ATLAS_2012_I1082009::analyze(), ATLAS_2014_I1326641::analyze(), CMS_2011_S9120041::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2013_I1243871::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2012_I1082936::analyze(), CDF_2008_S7541902::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2013_I1244522::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2013_I1190187::analyze(), ATLAS_2014_I1307243::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2011_S9128077::analyze(), CDF_1996_S3349578::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_001::analyze(), MC_TTBAR::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2013_I1230812::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2014_I1327229::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2011_S9041966::analyze(), MC_SUSY::analyze(), ATLAS_2014_I1306615::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2012_I1094568::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2011_S9126244::analyze(), ATLAS_2012_I1203852::analyze(), JetShape::calc(), ATLAS_2011_I945498::selectJets(), and Jet::totalEnergy().

{ return _momentum; }
double neutralEnergy ( ) const

Get the energy carried in this jet by neutral particles.

Todo:
Jet::containsMatch(Matcher m) { ... if m(pid) return true; ... }

Definition at line 91 of file Jet.cc.

References ParticleBase::E(), Jet::particles(), Particle::pid(), and Rivet::PID::threeCharge().

Referenced by STAR_2009_UE_HELEN::analyze().

                                  {
    double e_neutral = 0.0;
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pid();
      if (PID::threeCharge(pid) == 0) {
        e_neutral += p.E();
      }
    }
    return e_neutral;
  }
operator const fastjet::PseudoJet & ( ) const [inline]

Cast operator to FastJet3 PseudoJet (as a const reference)

Definition at line 178 of file Jet.hh.

References Jet::pseudojet().

{ return pseudojet(); }
operator const FourMomentum & ( ) const [inline, inherited]

Cast operator for conversion to FourMomentum.

Definition at line 31 of file ParticleBase.hh.

References ParticleBase::momentum().

{ return momentum(); }
const vector<Particle>& particles ( ) const [inline]

Get the particles in this jet (const version)

Definition at line 55 of file Jet.hh.

References Jet::_particles.

{ return _particles; }
double perp2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly (alias).

Definition at line 56 of file ParticleBase.hh.

References ParticleBase::pt2().

{ return pt2(); }
double phi ( const PhiMapping  mapping = ZERO_2PI) const [inline, inherited]
double polarAngle ( ) const [inline, inherited]

Angle subtended by the 3-vector and the z-axis.

Definition at line 103 of file ParticleBase.hh.

References ParticleBase::momentum(), and FourVector::polarAngle().

{ return momentum().polarAngle(); }
const fastjet::PseudoJet& pseudojet ( ) const [inline]

Access the internal FastJet3 PseudoJet (as a const reference)

Definition at line 175 of file Jet.hh.

References Jet::_pseudojet.

Referenced by Jet::operator const fastjet::PseudoJet &().

{ return _pseudojet; }
double pseudorapidity ( ) const [inline, inherited]

Get the $ \eta $ directly.

Definition at line 67 of file ParticleBase.hh.

References FourVector::eta(), and ParticleBase::momentum().

{ return momentum().eta(); }
double pt ( ) const [inline, inherited]

Get the $ p_T $ directly.

Definition at line 45 of file ParticleBase.hh.

References ParticleBase::momentum(), and FourMomentum::pt().

Referenced by ParticleBase::perp(), and ParticleBase::pT().

{ return momentum().pt(); }
double pT ( ) const [inline, inherited]

Get the $ p_T $ directly (alias).

Definition at line 47 of file ParticleBase.hh.

References ParticleBase::pt().

Referenced by ATLAS_2012_I1118269::analyze(), MC_JetAnalysis::analyze(), CMS_2012_I1087342::analyze(), TOTEM_2012_002::analyze(), CDF_2006_S6450792::analyze(), CMS_2010_S8656010::analyze(), ATLAS_2010_S8591806::analyze(), CMS_2011_S9088458::analyze(), CMS_2012_PAS_QCD_11_010::analyze(), STAR_2008_S7993412::analyze(), LHCF_2012_I1115479::analyze(), ALICE_2014_I1300380::analyze(), ALICE_2011_S8909580::analyze(), ALICE_2011_S8945144::analyze(), CMS_2011_S9086218::analyze(), CDF_2007_S7057202::analyze(), ATLAS_2012_I1188891::analyze(), ATLAS_2014_I1282441::analyze(), ALICE_2015_I1357424::analyze(), CMS_2010_S8547297::analyze(), MC_HINC::analyze(), ATLAS_2011_S8994773::analyze(), CMS_2013_I1273574::analyze(), CDF_1988_S1865951::analyze(), CDF_2008_S7828950::analyze(), ATLAS_2010_CONF_2010_049::analyze(), ZEUS_2001_S4815815::analyze(), ALICE_2010_S8706239::analyze(), CDF_2008_S7540469::analyze(), CDF_2012_NOTE10874::analyze(), CMS_2011_S8973270::analyze(), ATLAS_2011_I930220::analyze(), MC_ZINC::analyze(), CDF_1996_S3108457::analyze(), SFM_1984_S1178091::analyze(), CMS_2015_I1346843::analyze(), MC_DIJET::analyze(), MC_HFJETS::analyze(), MC_LEADJETUE::analyze(), STAR_2006_S6500200::analyze(), ATLAS_2014_I1325553::analyze(), CDF_2009_S8233977::analyze(), ATLAS_2014_I1268975::analyze(), ATLAS_2014_I1298811::analyze(), MC_WINC::analyze(), CMS_2015_I1385107::analyze(), ATLAS_2013_I1219109::analyze(), STAR_2006_S6860818::analyze(), CMS_2011_S8978280::analyze(), CMS_2011_S8884919::analyze(), CMS_2014_I1298810::analyze(), ATLAS_2010_S8817804::analyze(), D0_2008_S7662670::analyze(), MC_ZZJETS::analyze(), ATLAS_2014_I1326641::analyze(), CMS_2011_S9120041::analyze(), MC_PHOTONS::analyze(), MC_WWJETS::analyze(), CMS_2015_I1310737::analyze(), STAR_2009_UE_HELEN::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_CONF_2012_104::analyze(), LHCB_2013_I1208105::analyze(), ATLAS_2013_I1243871::analyze(), UA1_1990_S2044935::analyze(), CMS_2012_I1107658::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1186556::analyze(), LHCB_2013_I1218996::analyze(), ATLAS_2012_I1082936::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2013_I1244522::analyze(), ATLAS_2012_I1124167::analyze(), CDF_2010_S8591881_QCD::analyze(), ATLAS_2013_I1217863_Z::analyze(), ATLAS_2011_I921594::analyze(), ATLAS_2014_I1306294::analyze(), ATLAS_2013_I1217863_W::analyze(), ATLAS_2013_I1263495::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), MC_GENERIC::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2013_I1190187::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1125961::analyze(), STAR_2008_S7869363::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_CONF_2012_109::analyze(), MC_TTBAR::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2011_I926145::analyze(), CDF_2001_S4751469::analyze(), CMS_2013_I1224539_WJET::analyze(), CMS_2014_I1303894::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2012_I1204447::analyze(), CMS_2013_I1224539_ZJET::analyze(), ATLAS_2012_CONF_2012_153::analyze(), ATLAS_2014_I1327229::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2010_S8894728::analyze(), ATLAS_2012_I1125575::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2014_I1312627::analyze(), MC_SUSY::analyze(), ATLAS_2011_S9041966::analyze(), ATLAS_2015_I1345452::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2011_I944826::analyze(), ATLAS_2012_I1094568::analyze(), ATLAS_2014_I1319490::analyze(), ATLAS_2012_I1093734::analyze(), ATLAS_2014_I1282447::analyze(), ATLAS_2011_I919017::analyze(), ATLAS_2012_I1204447::apply_reco_eff(), ATLAS_2014_I1327229::apply_reco_eff(), HeavyHadrons::bHadrons(), JetShape::calc(), HeavyHadrons::cHadrons(), CMS_2013_I1261026::eventDecomp(), ATLAS_2012_I1084540::fillMap(), ATLAS_2014_I1319490::fillPlots(), ATLAS_2010_S8918562::fillPtEtaNch(), CMS_2013_I1258128::makePhotonCut(), STARRandomFilter::operator()(), FinalState::project(), LeadingParticlesFinalState::project(), VetoedFinalState::project(), and Variables::Variables().

{ return pt(); }
double pt2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly.

Definition at line 52 of file ParticleBase.hh.

References ParticleBase::momentum(), and FourMomentum::pt2().

Referenced by ParticleBase::perp2(), and ParticleBase::pT2().

{ return momentum().pt2(); }
double pT2 ( ) const [inline, inherited]

Get the $ p_T^2 $ directly (alias).

Definition at line 54 of file ParticleBase.hh.

References ParticleBase::pt2().

Referenced by ATLAS_2011_S9108483::analyze().

{ return pt2(); }
double px ( ) const [inline, inherited]
double py ( ) const [inline, inherited]
double pz ( ) const [inline, inherited]
double rap ( ) const [inline, inherited]

Get the $ y $ directly (alias).

Definition at line 78 of file ParticleBase.hh.

References ParticleBase::momentum(), and FourMomentum::rapidity().

Referenced by ATLAS_2013_I1244522::analyze().

{ return momentum().rapidity(); }
Jet& setConstituents ( const Particles particles) [inline]

Definition at line 203 of file Jet.hh.

References Jet::setParticles().

{ return setParticles(particles); }
Jet & setParticles ( const Particles particles)

Set the particles collection with full particle information.

If set, this overrides particle info extracted from the PseudoJet

Definition at line 55 of file Jet.cc.

References Jet::_particles, and Jet::particles().

Referenced by Jet::setConstituents().

                                                          {
    _particles = particles;
    return *this;
  }
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::Jet().

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.

Definition at line 18 of file Jet.cc.

References Jet::_momentum, Jet::_particles, Jet::_pseudojet, Jet::_tags, Jet::clear(), FourMomentum::E(), ParticleBase::mom(), Jet::particles(), FourMomentum::px(), FourMomentum::py(), FourMomentum::pz(), and Jet::tags().

                                                                                               {
    clear();
    _momentum = mom;
    _pseudojet = fastjet::PseudoJet(mom.px(), mom.py(), mom.pz(), mom.E());
    _particles = particles;
    _tags = tags;
    return *this;
  }
Jet& setState ( const Particles particles,
const FourMomentum mom 
) [inline]
Deprecated:
Prefer the 4-mom first-arg versions

Definition at line 197 of file Jet.hh.

References Jet::setState().

Referenced by Jet::setState().

{ return setState(mom, particles); }
size_t size ( ) const [inline]

Number of particles in this jet.

Definition at line 50 of file Jet.hh.

References Jet::_particles.

Referenced by ATLAS_2015_CONF_2015_041::analyze(), and MC_TTBAR::analyze().

{ return _particles.size(); }
Particles& tags ( ) [inline]

Particles which have been tag-matched to this jet.

Definition at line 87 of file Jet.hh.

References Jet::_tags.

Referenced by Jet::bTags(), Jet::cTags(), Jet::Jet(), Jet::setState(), Jet::tags(), and Jet::tauTags().

{ return _tags; }
const Particles& tags ( ) const [inline]

Particles which have been tag-matched to this jet (const version)

Definition at line 89 of file Jet.hh.

References Jet::_tags.

{ return _tags; }
Particles 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.

Definition at line 153 of file Jet.cc.

References Jet::tags().

                                        {
    Particles rtn;
    foreach (const Particle& p, tags()) {
      if (c->accept(p)) rtn.push_back(p);
    }
    return rtn;
  }
bool tauTagged ( const Cut c = Cuts::open()) const [inline]

Does this jet have at least one tau-tag (that passes an optional Cut)?

Definition at line 117 of file Jet.hh.

References Jet::tauTags().

{ return !tauTags(c).empty(); }
Particles 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.

Definition at line 181 of file Jet.cc.

References Rivet::PID::isTau(), and Jet::tags().

Referenced by MC_JETTAGS::analyze(), and Jet::tauTagged().

                                           {
    Particles rtn;
    foreach (const Particle& tp, tags()) {
      if (isTau(tp) && c->accept(tp)) rtn.push_back(tp);
    }
    return rtn;
  }
double theta ( ) const [inline, inherited]

Synonym for polarAngle.

Definition at line 105 of file ParticleBase.hh.

References ParticleBase::momentum(), and FourVector::theta().

Referenced by CMS_2012_I1184941::analyze(), and ALEPH_1996_S3196992::analyze().

{ return momentum().theta(); }
double totalEnergy ( ) const [inline]

Get the total energy of this jet.

Definition at line 160 of file Jet.hh.

References FourMomentum::E(), and Jet::momentum().

Referenced by STAR_2009_UE_HELEN::analyze().

{ return momentum().E(); }

Member Data Documentation

FourMomentum _momentum [mutable, private]

Effective jet 4-vector (just for caching)

Definition at line 224 of file Jet.hh.

Referenced by Jet::clear(), Jet::momentum(), and Jet::setState().

Particles _particles [private]

Full constituent particle information. (Filled from PseudoJet if possible.)

Todo:
Make these mutable or similar? Add a flag to force a cache rebuild?

Definition at line 218 of file Jet.hh.

Referenced by Jet::clear(), Jet::particles(), Jet::setParticles(), Jet::setState(), and Jet::size().

fastjet::PseudoJet _pseudojet [private]

FJ3 PseudoJet member to unify PseudoJet and Jet.

Definition at line 214 of file Jet.hh.

Referenced by Jet::clear(), Jet::pseudojet(), and Jet::setState().

Particles _tags [private]

Particles used to tag this jet (can be anything, but c and b hadrons are the most common)

Definition at line 221 of file Jet.hh.

Referenced by Jet::setState(), and Jet::tags().


The documentation for this class was generated from the following files: