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 ()
 Jet (const vector< Particle > &particles, const FourMomentum &pjet)
 Set all the jet data, with full particle information.
Access jet constituents
Todo:
Add a constructor from PseudoJet
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)
bool containsParticle (const Particle &particle) const
 Check whether this jet contains a particular particle.
bool containsParticleId (PdgId pid) const
 Check whether this jet contains a certain particle type.
bool containsParticleId (const vector< PdgId > &pids) const
 Check whether this jet contains at least one of certain particle types.
bool containsCharm () const
 Check whether this jet contains a charm-flavoured hadron (or decay products from one).
bool containsBottom () const
 Check whether this jet contains a bottom-flavoured hadron (or decay products from one).
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.
Set the jet constituents and properties
JetsetState (const vector< Particle > &particles, const FourMomentum &pjet)
 Set all the jet data, with full particle information.
JetsetMomentum (const FourMomentum &momentum)
 Set the effective 4-momentum of the jet.
JetsetParticles (const vector< Particle > &particles)
 Set the particles collection with full particle information.
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
Todo:
Add a converter and cast operator to FJ3 PseudoJet
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

Particles _particles
 Full particle information including tracks, ID etc. (caching PseudoJet properties)
FourMomentum _momentum
 Effective jet 4-vector (caching PseudoJet properties)

Detailed Description

Representation of a clustered jet of particles.

Definition at line 15 of file Jet.hh.


Constructor & Destructor Documentation

Jet ( ) [inline]

Definition at line 21 of file Jet.hh.

References Jet::clear().

: ParticleBase() { clear(); }
Jet ( const vector< Particle > &  particles,
const FourMomentum pjet 
) [inline]

Set all the jet data, with full particle information.

Definition at line 24 of file Jet.hh.

References Jet::setState().

      : ParticleBase() {
      setState(particles, pjet);
    }

Member Function Documentation

double abspseudorapidity ( ) const [inline, inherited]

Get the $ |\eta| $ directly.

Definition at line 75 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 84 of file ParticleBase.hh.

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

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

Angle between this vector and another.

Definition at line 112 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 114 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 116 of file ParticleBase.hh.

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

{ return momentum().angle(v3); }
Jet & clear ( )

Reset this jet as empty.

Definition at line 117 of file Jet.cc.

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

Referenced by Jet::Jet().

                  {
    //_momenta.clear();
    _particles.clear();
    _momentum = FourMomentum();
    return *this;
  }
bool containsBottom ( ) const

Check whether this jet contains a bottom-flavoured hadron (or decay products from one).

Definition at line 100 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().

Referenced by CDF_2008_S7782535::analyze(), EXAMPLE::analyze(), ATLAS_2012_I1126136::analyze(), MC_TTBAR::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_I1095236::analyze(), and MC_VH2BB::analyze().

                                 {
    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;
      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 ( ) const

Check whether this jet contains a charm-flavoured hadron (or decay products from one).

Definition at line 83 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;
      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 29 of file Jet.cc.

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

                                                           {
    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 38 of file Jet.cc.

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

                                              {
    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 46 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;
  }
double energy ( ) const [inline, inherited]

Get the energy directly (alias).

Definition at line 46 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 73 of file ParticleBase.hh.

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

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

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

Get the energy carried in this jet by hadrons.

Definition at line 71 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 mass ( ) const [inline, inherited]
double mass2 ( ) const [inline, inherited]

Get the mass**2 directly.

Definition at line 68 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().

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

Get equivalent single momentum four-vector.

Implements ParticleBase.

Definition at line 69 of file Jet.hh.

References Jet::_momentum.

Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), STAR_2006_S6870392::analyze(), D0_1996_S3324664::analyze(), D0_2008_S6879055::analyze(), CDF_1993_S2742446::analyze(), ATLAS_2012_I1188891::analyze(), ATLAS_2011_I930220::analyze(), CDF_2008_S7782535::analyze(), CMS_2013_I1272853::analyze(), CDF_1996_S3108457::analyze(), CDF_1997_S3541940::analyze(), ATLAS_2010_S8817804::analyze(), ATLAS_2012_I1082009::analyze(), CMS_2011_S9120041::analyze(), ATLAS_2013_I1243871::analyze(), ATLAS_2012_I1082936::analyze(), CDF_2008_S7541902::analyze(), ATLAS_2010_S8919674::analyze(), MC_PHOTONJETUE::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9128077::analyze(), ATLAS_2011_S9212183::analyze(), MC_TTBAR::analyze(), ATLAS_2013_I1190187::analyze(), ATLAS_2012_CONF_2012_105::analyze(), CDF_1996_S3349578::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2013_I1230812::analyze(), ATLAS_2012_CONF_2012_153::analyze(), MC_SUSY::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2012_I943401::analyze(), MC_VH2BB::analyze(), ATLAS_2011_S9225137::analyze(), ATLAS_2011_S9041966::analyze(), ATLAS_2012_I1094568::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2011_S9126244::analyze(), ATLAS_2012_I1203852::analyze(), JetShape::calc(), JetAlg::jets(), CMS_2013_I1258128::makePhotonCut(), CMS_2013_I1258128::makeZCut(), ATLAS_2011_I945498::selectJets(), Jet::setMomentum(), 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 59 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 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 45 of file Jet.hh.

References Jet::_particles.

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

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

Definition at line 60 of file ParticleBase.hh.

References ParticleBase::pt2().

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

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

Definition at line 107 of file ParticleBase.hh.

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

{ return momentum().polarAngle(); }
double pseudorapidity ( ) const [inline, inherited]

Get the $ \eta $ directly.

Definition at line 71 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 49 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 51 of file ParticleBase.hh.

References ParticleBase::pt().

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

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

Get the $ p_T^2 $ directly.

Definition at line 56 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 58 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]

z component of momentum.

Definition at line 104 of file ParticleBase.hh.

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

Referenced by SFM_1984_S1178091::analyze(), ATLAS_2011_S9108483::analyze(), and ALEPH_1996_S3196992::particleInJet().

{ return momentum().z(); }
double rap ( ) const [inline, inherited]

Get the $ y $ directly (alias).

Definition at line 82 of file ParticleBase.hh.

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

{ return momentum().rapidity(); }
Jet & setMomentum ( const FourMomentum momentum)

Set the effective 4-momentum of the jet.

Definition at line 17 of file Jet.cc.

References Jet::_momentum, and Jet::momentum().

Referenced by Jet::setState().

                                                    {
    _momentum = momentum;
    return *this;
  }
Jet & setParticles ( const vector< Particle > &  particles)

Set the particles collection with full particle information.

Definition at line 23 of file Jet.cc.

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

Referenced by Jet::setState().

                                                          {
    _particles = particles;
    return *this;
  }
Jet & setState ( const vector< Particle > &  particles,
const FourMomentum pjet 
)

Set all the jet data, with full particle information.

Definition at line 10 of file Jet.cc.

References Jet::setMomentum(), and Jet::setParticles().

Referenced by Jet::Jet().

                                                                                {
    setParticles(particles);
    setMomentum(pjet);
    return *this;
  }
size_t size ( ) const [inline]

Number of particles in this jet.

Definition at line 39 of file Jet.hh.

References Jet::_particles.

Referenced by MC_TTBAR::analyze().

{ return _particles.size(); }
double theta ( ) const [inline, inherited]

Synonym for polarAngle.

Definition at line 109 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 72 of file Jet.hh.

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

Referenced by STAR_2009_UE_HELEN::analyze().

{ return momentum().E(); }

Member Data Documentation

Effective jet 4-vector (caching PseudoJet properties)

Definition at line 122 of file Jet.hh.

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

Particles _particles [private]

Full particle information including tracks, ID etc. (caching PseudoJet properties)

Todo:
Add a FJ3 PseudoJet member to unify PseudoJet and Jet

Definition at line 119 of file Jet.hh.

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


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