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
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 the effective jet 4-vector properties
const FourMomentummomentum () const
 Get equivalent single momentum four-vector.
double eta () const
 Get the unweighted average $ \eta $ for this jet. (caches)
double phi () const
 Get the unweighted average $ \phi $ for this jet. (caches)
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.
double ptSum () const
 Get the sum of the $ p_T $ values of the constituent tracks. (caches)
double EtSum () const
 Get the sum of the $ E_T $ values of the constituent tracks. (caches)
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.

Private Attributes

ParticleVector _particles
 Full particle information including tracks, ID etc.
FourMomentum _momentum
 Effective jet 4-vector.

Detailed Description

Representation of a clustered jet of particles.

Definition at line 12 of file Jet.hh.


Constructor & Destructor Documentation

Jet ( ) [inline]

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

References Jet::setState().

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

Member Function Documentation

Jet & clear ( )

Reset this jet as empty.

Definition at line 146 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 129 of file Jet.cc.

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

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.pdgId();
      if (abs(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 112 of file Jet.cc.

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

                                {
    foreach (const Particle& p, particles()) {
      const PdgId pid = p.pdgId();
      if (abs(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 58 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 67 of file Jet.cc.

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

                                              {
    foreach (const Particle& p, particles()) {
      if (p.pdgId() == 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 75 of file Jet.cc.

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

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

Get the unweighted average $ \eta $ for this jet. (caches)

Definition at line 108 of file Jet.hh.

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

Referenced by STAR_2006_S6870392::analyze(), ATLAS_2012_I1082009::analyze(), and ATLAS_2010_S8919674::analyze().

{ return momentum().eta(); }
double EtSum ( ) const [inline]

Get the sum of the $ E_T $ values of the constituent tracks. (caches)

Definition at line 126 of file Jet.hh.

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

Referenced by CDF_2004_S5839831::analyze(), Rivet::cmpJetsByEt(), and Rivet::cmpJetsByEtDesc().

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

Get the energy carried in this jet by hadrons.

Definition at line 100 of file Jet.cc.

References FourMomentum::E(), Rivet::PID::isHadron(), Particle::momentum(), Jet::particles(), and Particle::pdgId().

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

Get equivalent single momentum four-vector.

Implements ParticleBase.

Definition at line 105 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(), MC_JetAnalysis::analyze(), CMS_2012_I1087342::analyze(), CMS_2011_S9088458::analyze(), CDF_2006_S6450792::analyze(), CMS_2011_S9086218::analyze(), CMS_2011_S8957746::analyze(), STAR_2006_S6870392::analyze(), CDF_2007_S7057202::analyze(), D0_1996_S3324664::analyze(), CDF_1998_S3618439::analyze(), CDF_2001_S4563131::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CDF_2008_S7828950::analyze(), D0_2008_S6879055::analyze(), ZEUS_2001_S4815815::analyze(), CDF_1993_S2742446::analyze(), CDF_2008_S7540469::analyze(), MC_DIJET::analyze(), CDF_2008_S7782535::analyze(), ATLAS_2011_S8971293::analyze(), CDF_1996_S3108457::analyze(), CDF_1997_S3541940::analyze(), D0_2008_S7863608::analyze(), ATLAS_2010_S8817804::analyze(), D0_2008_S7662670::analyze(), ATLAS_2012_I1082009::analyze(), D0_2009_S8349509::analyze(), STAR_2009_UE_HELEN::analyze(), D0_2009_S8202443::analyze(), ATLAS_2012_I1082936::analyze(), MC_PHOTONJETUE::analyze(), CDF_2008_S7541902::analyze(), D0_1996_S3214044::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1186556::analyze(), MC_TTBAR::analyze(), ATLAS_2011_S9212183::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_2011_S9128077::analyze(), ATLAS_2012_I1125961::analyze(), MC_ZZJETS::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1117704::analyze(), MC_WWJETS::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2012_CONF_2012_153::analyze(), MC_SUSY::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2011_S9041966::analyze(), ATLAS_2012_I1094568::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2011_S9126244::analyze(), ATLAS_2011_I919017::analyze(), JetShape::calc(), Rivet::cmpJetsByAscAbsPseudorapidity(), Rivet::cmpJetsByAscAbsRapidity(), Rivet::cmpJetsByAscE(), Rivet::cmpJetsByAscP(), Rivet::cmpJetsByAscPseudorapidity(), Rivet::cmpJetsByAscRapidity(), Rivet::cmpJetsByDescAbsPseudorapidity(), Rivet::cmpJetsByDescAbsRapidity(), Rivet::cmpJetsByDescPseudorapidity(), Rivet::cmpJetsByDescRapidity(), Rivet::cmpJetsByE(), Rivet::cmpJetsByP(), Rivet::deltaEta(), Rivet::deltaPhi(), Rivet::deltaR(), Jet::eta(), Jet::EtSum(), JetAlg::jets(), Jet::phi(), Jet::ptSum(), 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 88 of file Jet.cc.

References FourMomentum::E(), Particle::momentum(), Jet::particles(), Particle::pdgId(), 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.pdgId();
      if (PID::threeCharge(pid) == 0) {
        e_neutral += p.momentum().E();
      }
    }
    return e_neutral;
  }
const vector<Particle>& particles ( ) const [inline]

Get the particles in this jet (const version)

Definition at line 81 of file Jet.hh.

References Jet::_particles.

{ return _particles; }
double phi ( ) const [inline]

Get the unweighted average $ \phi $ for this jet. (caches)

Definition at line 111 of file Jet.hh.

References Jet::momentum(), and FourVector::phi().

Referenced by CDF_2001_S4751469::analyze().

{ return momentum().phi(); }
double ptSum ( ) const [inline]

Get the sum of the $ p_T $ values of the constituent tracks. (caches)

Definition at line 123 of file Jet.hh.

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

Referenced by CDF_2001_S4751469::analyze(), Rivet::cmpJetsByAscPt(), and Rivet::cmpJetsByPt().

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

Set the effective 4-momentum of the jet.

Definition at line 24 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 30 of file Jet.cc.

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

Referenced by Jet::setState().

                                                         {
   _particles = particles;
   // foreach (const Particle& p, particles) {
   //   _momenta.push_back(p.momentum());
   // }
   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 totalEnergy ( ) const [inline]

Get the total energy of this jet.

Definition at line 114 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.

Definition at line 171 of file Jet.hh.

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

Full particle information including tracks, ID etc.

Definition at line 164 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: