rivet is hosted by Hepforge, IPPP Durham
Particle Class Reference

Representation of particles from a HepMC::GenEvent. More...

#include <Particle.hh>

Inheritance diagram for Particle:
Collaboration diagram for Particle:

List of all members.

Public Member Functions

 Particle ()
 Particle (PdgId pid, const FourMomentum &mom)
 Constructor without GenParticle.
 Particle (const GenParticle &gp)
 Constructor from a HepMC GenParticle.
 Particle (const GenParticle *gp)
 Constructor from a HepMC GenParticle pointer.
Basic particle specific properties
const GenParticle * genParticle () const
 Get a const reference to the original GenParticle.
 operator const GenParticle * () const
 Cast operator for conversion to GenParticle*.
const FourMomentummomentum () const
 The momentum.
ParticlesetMomentum (const FourMomentum &momentum)
 Set the momentum.
Particle ID code accessors
PdgId pid () const
 This Particle's PDG ID code.
PdgId abspid () const
 Absolute value of the PDG ID code.
PdgId pdgId () const
Particle species properties
double charge () const
 The charge of this Particle.
int threeCharge () const
 Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
bool isHadron () const
 Is this a hadron?
bool isMeson () const
 Is this a meson?
bool isBaryon () const
 Is this a baryon?
bool isLepton () const
 Is this a lepton?
bool isNeutrino () const
 Is this a neutrino?
bool hasBottom () const
 Does this (hadron) contain a b quark?
bool hasCharm () const
 Does this (hadron) contain a c quark?
Ancestry properties
bool hasAncestor (PdgId pdg_id) const
bool fromDecay () const
 Determine whether the particle is from a hadron or tau decay.
bool fromBottom () const
 Determine whether the particle is from a b-hadron decay.
bool fromCharm () const
 Determine whether the particle is from a c-hadron decay.
bool fromTau () const
 Determine whether the particle is from a tau decay.
Effective momentum accessors
 operator const FourMomentum & () const
 Cast operator for conversion to FourMomentum.
Convenience access to the effective 4-vector properties
Todo:
Add a 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 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 phi () const
 Get the $ \phi $ directly.

Private Attributes

const GenParticle * _original
 A pointer to the original GenParticle from which this Particle is projected.
PdgId _id
 The PDG ID code for this Particle.
FourMomentum _momentum
 The momentum of this projection of the Particle.

Detailed Description

Representation of particles from a HepMC::GenEvent.

Definition at line 17 of file Particle.hh.


Constructor & Destructor Documentation

Particle ( ) [inline]

Default constructor.

Deprecated:
A particle without info is useless. This only exists to keep STL containers happy.

Definition at line 22 of file Particle.hh.

      : ParticleBase(),
        _original(0), _id(0), _momentum()
    { }
Particle ( PdgId  pid,
const FourMomentum mom 
) [inline]

Constructor without GenParticle.

Definition at line 28 of file Particle.hh.

      : ParticleBase(),
        _original(0), _id(pid), _momentum(mom)
    { }
Particle ( const GenParticle &  gp) [inline]

Constructor from a HepMC GenParticle.

Definition at line 34 of file Particle.hh.

      : ParticleBase(),
        _original(&gp), _id(gp.pdg_id()),
        _momentum(gp.momentum())
    { }
Particle ( const GenParticle *  gp) [inline]

Constructor from a HepMC GenParticle pointer.

Definition at line 41 of file Particle.hh.

      : ParticleBase(),
        _original(gp), _id(gp->pdg_id()),
        _momentum(gp->momentum())
    { }

Member Function Documentation

double abseta ( ) const [inline, inherited]

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

Definition at line 68 of file ParticleBase.hh.

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

Referenced by ATLAS_2011_S9225137::analyze().

{ return momentum().abseta(); }
PdgId abspid ( ) const [inline]

Absolute value of the PDG ID code.

Definition at line 80 of file Particle.hh.

References Particle::_id.

Referenced by LHCB_2013_I1218996::analyze().

{ return _id; }
double abspseudorapidity ( ) const [inline, inherited]

Get the $ |\eta| $ directly.

Definition at line 66 of file ParticleBase.hh.

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

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

Get the $ |y| $ directly (alias).

Definition at line 77 of file ParticleBase.hh.

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

Referenced by ATLAS_2011_I930220::analyze(), and LHCB_2013_I1218996::analyze().

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

Get the $ |y| $ directly.

Definition at line 75 of file ParticleBase.hh.

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

{ return momentum().absrapidity(); }
double charge ( ) const [inline]

The charge of this Particle.

Definition at line 93 of file Particle.hh.

References Particle::pdgId().

Referenced by NeutralFinalState::project().

                          {
      return PID::charge(pdgId());
    }
double E ( ) const [inline, inherited]

Get the energy directly.

Definition at line 42 of file ParticleBase.hh.

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

Referenced by Jet::hadronicEnergy(), and Jet::neutralEnergy().

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

Get the energy directly (alias).

Definition at line 44 of file ParticleBase.hh.

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

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

Get the $ E_T $ directly.

Definition at line 54 of file ParticleBase.hh.

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

Referenced by CDF_2004_S5839831::analyze(), CentralEtHCM::project(), and NeutralFinalState::project().

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

Get the $ \eta $ directly (alias).

Definition at line 64 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_I1087342::analyze(), CMS_2010_S8656010::analyze(), CMS_2012_I1184941::analyze(), TOTEM_2012_I1115294::analyze(), ATLAS_2010_S8591806::analyze(), CMS_2011_S9088458::analyze(), STAR_2008_S7993412::analyze(), ALICE_2012_I1181770::analyze(), CMS_2011_S8957746::analyze(), STAR_2006_S6870392::analyze(), MC_DIPHOTON::analyze(), MC_HINC::analyze(), CDF_2005_S6080774::analyze(), CMS_2010_S8547297::analyze(), ATLAS_2011_S8994773::analyze(), D0_1996_S3324664::analyze(), CDF_2001_S4563131::analyze(), MC_ZINC::analyze(), ATLAS_2010_CONF_2010_049::analyze(), D0_2008_S6879055::analyze(), CDF_1993_S2742446::analyze(), CDF_2009_S8436959::analyze(), CMS_2011_S8973270::analyze(), CMS_2011_S9215166::analyze(), D0_2010_S8570965::analyze(), MC_WINC::analyze(), MC_IDENTIFIED::analyze(), ATLAS_2012_I1119557::analyze(), CDF_1997_S3541940::analyze(), CMS_2011_S8884919::analyze(), ATLAS_2012_I1082009::analyze(), D0_2009_S8202443::analyze(), CMS_2011_S9120041::analyze(), ATLAS_2012_I946427::analyze(), STAR_2009_UE_HELEN::analyze(), UA1_1990_S2044935::analyze(), ATLAS_2010_S8919674::analyze(), CMS_2013_I1218372::analyze(), ATLAS_2012_I1183818::analyze(), D0_1996_S3214044::analyze(), ATLAS_2012_CONF_2012_104::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2011_S9120807::analyze(), ATLAS_2011_S9128077::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2010_S8914702::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1125961::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_I943401::analyze(), ATLAS_2012_I1093738::analyze(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2011_S9041966::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2012_I1084540::fillMap(), ATLAS_2010_S8918562::fillPtEtaNch(), ATLAS_2012_I1091481::getSeta(), TriggerCDFRun0Run1::project(), TriggerCDFRun2::project(), NeutralFinalState::project(), TriggerUA5::project(), and FinalState::project().

{ return momentum().eta(); }
bool fromBottom ( ) const

Determine whether the particle is from a b-hadron decay.

Note:
This question is valid in MC, but may not be perfectly answerable experimentally -- use this function with care when replicating experimental analyses!
Todo:
Shouldn't a const vertex be being returned? Ah, HepMC...

Definition at line 34 of file Particle.cc.

References Particle::genParticle(), Particle::hasBottom(), Particle::isHadron(), Rivet::particles(), and Particle::pid().

Referenced by LHCB_2013_I1218996::analyze().

                                  {
    /// @todo Shouldn't a const vertex be being returned? Ah, HepMC...
    GenVertex* prodVtx = genParticle()->production_vertex();
    if (prodVtx == NULL) return false;
    foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
      const PdgId pid = ancestor->pdg_id();
      if (ancestor->status() == 2 && (PID::isHadron(pid) && PID::hasBottom(pid))) return true;
    }
    return false;
  }
bool fromCharm ( ) const

Determine whether the particle is from a c-hadron decay.

Note:
If a hadron contains b and c quarks it is considered a bottom hadron and NOT a charm hadron.
This question is valid in MC, but may not be perfectly answerable experimentally -- use this function with care when replicating experimental analyses!
Todo:
Shouldn't a const vertex be being returned? Ah, HepMC...

Definition at line 46 of file Particle.cc.

References Particle::genParticle(), Particle::hasBottom(), Particle::hasCharm(), Particle::isHadron(), Rivet::particles(), and Particle::pid().

                                 {
    /// @todo Shouldn't a const vertex be being returned? Ah, HepMC...
    GenVertex* prodVtx = genParticle()->production_vertex();
    if (prodVtx == NULL) return false;
    foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
      const PdgId pid = ancestor->pdg_id();
      if (ancestor->status() == 2 && (PID::isHadron(pid) && PID::hasCharm(pid) && !PID::hasBottom(pid))) return true;
    }
    return false;
  }
bool fromDecay ( ) const

Determine whether the particle is from a hadron or tau decay.

Specifically, walk up the ancestor chain until a status 2 hadron or tau is found, if at all.

Note:
This question is valid in MC, but may not be perfectly answerable experimentally -- use this function with care when replicating experimental analyses!
Todo:
Shouldn't a const vertex be being returned? Ah, HepMC...

Definition at line 22 of file Particle.cc.

References Particle::genParticle(), Particle::isHadron(), Rivet::particles(), Particle::pid(), and Rivet::PID::TAU.

Referenced by DressedLeptons::project().

                                 {
    /// @todo Shouldn't a const vertex be being returned? Ah, HepMC...
    GenVertex* prodVtx = genParticle()->production_vertex();
    if (prodVtx == NULL) return false;
    foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
      const PdgId pid = ancestor->pdg_id();
      if (ancestor->status() == 2 && (PID::isHadron(pid) || abs(pid) == PID::TAU)) return true;
    }
    return false;
  }
bool fromTau ( ) const

Determine whether the particle is from a tau decay.

Note:
This question is valid in MC, but may not be perfectly answerable experimentally -- use this function with care when replicating experimental analyses!
Todo:
Shouldn't a const vertex be being returned? Ah, HepMC...

Definition at line 58 of file Particle.cc.

References Particle::genParticle(), Rivet::particles(), Particle::pid(), and Rivet::PID::TAU.

                               {
    /// @todo Shouldn't a const vertex be being returned? Ah, HepMC...
    GenVertex* prodVtx = genParticle()->production_vertex();
    if (prodVtx == NULL) return false;
    foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
      const PdgId pid = ancestor->pdg_id();
      if (ancestor->status() == 2 && abs(pid) == PID::TAU) return true;
    }
    return false;
  }
bool hasAncestor ( PdgId  pdg_id) const

Check whether a given PID is found in the GenParticle's ancestor list

Note:
This question is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!
Todo:
Neaten this up with C++11, via one waliker function and several uses with lamba tests
Todo:
Shouldn't a const vertex be being returned? Ah, HepMC...

Definition at line 11 of file Particle.cc.

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

Referenced by ALICE_2011_S8909580::analyze(), ALICE_2011_S8945144::analyze(), and CMS_2011_S8978280::analyze().

                                               {
    /// @todo Shouldn't a const vertex be being returned? Ah, HepMC...
    GenVertex* prodVtx = genParticle()->production_vertex();
    if (prodVtx == 0) return false;
    foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
      if (ancestor->pdg_id() == pdg_id) return true;
    }
    return false;
  }
bool hasBottom ( ) const [inline]

Does this (hadron) contain a b quark?

Definition at line 117 of file Particle.hh.

References Particle::pdgId().

Referenced by Particle::fromBottom(), and Particle::fromCharm().

{ return PID::hasBottom(pdgId()); }
bool hasCharm ( ) const [inline]

Does this (hadron) contain a c quark?

Definition at line 120 of file Particle.hh.

References Particle::pdgId().

Referenced by LHCB_2013_I1218996::analyze(), and Particle::fromCharm().

{ return PID::hasCharm(pdgId()); }
bool isBaryon ( ) const [inline]

Is this a baryon?

Definition at line 108 of file Particle.hh.

References Particle::pdgId().

{ return PID::isBaryon(pdgId()); }
bool isHadron ( ) const [inline]

Is this a hadron?

Definition at line 102 of file Particle.hh.

References Particle::pdgId().

Referenced by LHCB_2013_I1218996::analyze(), Particle::fromBottom(), Particle::fromCharm(), and Particle::fromDecay().

{ return PID::isHadron(pdgId()); }
bool isLepton ( ) const [inline]

Is this a lepton?

Definition at line 111 of file Particle.hh.

References Particle::pdgId().

{ return PID::isLepton(pdgId()); }
bool isMeson ( ) const [inline]

Is this a meson?

Definition at line 105 of file Particle.hh.

References Particle::pdgId().

{ return PID::isMeson(pdgId()); }
bool isNeutrino ( ) const [inline]

Is this a neutrino?

Definition at line 114 of file Particle.hh.

References Particle::pdgId().

{ return PID::isNeutrino(pdgId()); }
double mass ( ) const [inline, inherited]

Get the mass directly.

Definition at line 57 of file ParticleBase.hh.

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

Referenced by ALICE_2011_S8945144::analyze().

{ return momentum().mass(); }
double mass2 ( ) const [inline, inherited]

Get the mass**2 directly.

Definition at line 59 of file ParticleBase.hh.

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

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

The momentum.

Implements ParticleBase.

Definition at line 62 of file Particle.hh.

References Particle::_momentum.

Referenced by CDF_2004_S5839831::_calcTransCones(), ClusteredLepton::addPhoton(), CLEO_2004_S5809304::analyze(), BABAR_2003_I593379::analyze(), BABAR_2007_S6895344::analyze(), BABAR_2005_S6181155::analyze(), BELLE_2006_S6265367::analyze(), BELLE_2001_S4598261::analyze(), ARGUS_1993_S2653028::analyze(), ATLAS_2011_I894867::analyze(), OPAL_1994_S2927284::analyze(), H1_2000_S4129130::analyze(), OPAL_1993_S2692198::analyze(), CMS_2012_I1193338::analyze(), SLD_2004_S5693039::analyze(), OPAL_1998_S3780481::analyze(), CMS_2012_I1184941::analyze(), STAR_2008_S7993412::analyze(), SLD_1999_S3743934::analyze(), OPAL_1995_S3198391::analyze(), OPAL_1997_S3608263::analyze(), OPAL_2000_S4418603::analyze(), ARGUS_1993_S2669951::analyze(), ALICE_2012_I1181770::analyze(), ALEPH_2002_S4823664::analyze(), H1_1994_S2919893::analyze(), OPAL_1996_S3257789::analyze(), OPAL_1998_S3702294::analyze(), DELPHI_1999_S3960137::analyze(), ATLAS_2011_S9002537::analyze(), CMS_QCD_10_024::analyze(), BABAR_2007_S7266081::analyze(), CDF_2005_S6080774::analyze(), UA5_1982_S875503::analyze(), DELPHI_1995_S3137023::analyze(), EXAMPLE_CUTS::analyze(), MC_DIPHOTON::analyze(), ATLAS_2011_S8994773::analyze(), CDF_1990_S2089246::analyze(), ATLAS_2011_S9035664::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CDF_1993_S2742446::analyze(), CDF_2009_S8436959::analyze(), CDF_2008_S7540469::analyze(), OPAL_2002_S5361494::analyze(), SLD_1996_S3398250::analyze(), CDF_2012_NOTE10874::analyze(), DELPHI_2000_S4328825::analyze(), CMS_2011_S8973270::analyze(), CMS_2011_S9215166::analyze(), MC_PHOTONINC::analyze(), OPAL_1998_S3749908::analyze(), SFM_1984_S1178091::analyze(), D0_2010_S8570965::analyze(), LHCB_2011_I919315::analyze(), MC_PHOTONKTSPLITTINGS::analyze(), MC_WINC::analyze(), OPAL_1997_S3396100::analyze(), D0_2006_S6438750::analyze(), ALICE_2010_S8625980::analyze(), MC_LEADJETUE::analyze(), ALEPH_1996_S3196992::analyze(), CDF_2009_S8233977::analyze(), UA5_1986_S1583476::analyze(), MC_PHOTONJETS::analyze(), MC_PHOTONS::analyze(), ATLAS_2012_I1082009::analyze(), CDF_2008_S8095620::analyze(), CMS_2011_S9120041::analyze(), LHCB_2010_S8758301::analyze(), ATLAS_2012_I946427::analyze(), CDF_2006_S6653332::analyze(), STAR_2009_UE_HELEN::analyze(), CMS_2012_I1107658::analyze(), ARGUS_1993_S2789213::analyze(), ATLAS_2013_I1243871::analyze(), LHCB_2013_I1208105::analyze(), ALEPH_1999_S4193598::analyze(), CDF_2008_S7541902::analyze(), CMS_2013_I1218372::analyze(), MC_PHOTONJETUE::analyze(), ATLAS_2012_I1183818::analyze(), CDF_2010_S8591881_QCD::analyze(), LHCB_2011_I917009::analyze(), ATLAS_2012_CONF_2012_104::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2011_S9108483::analyze(), TASSO_1990_S2148048::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2011_S9120807::analyze(), ATLAS_2012_CONF_2012_105::analyze(), MC_TTBAR::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1112263::analyze(), CDF_2001_S4751469::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2010_S8914702::analyze(), ATLAS_2012_I1095236::analyze(), D0_2008_S7719523::analyze(), ATLAS_2012_I1125961::analyze(), CMS_2013_I1224539_WJET::analyze(), CMS_2013_I1224539_ZJET::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_CONF_2012_109::analyze(), LHCB_2012_I1119400::analyze(), ATLAS_2012_I1180197::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(), MC_SUSY::analyze(), ATLAS_2011_S9212353::analyze(), ATLAS_2012_I943401::analyze(), ATLAS_2012_I1093738::analyze(), ATLAS_2011_S9225137::analyze(), ATLAS_2011_S9041966::analyze(), MC_VH2BB::analyze(), ATLAS_2012_I1094568::analyze(), ALEPH_1996_S3486095::analyze(), ALEPH_2004_S5765862::analyze(), DELPHI_1996_S3430090::analyze(), CDF_2004_S5839831::analyze(), ATLAS_2012_I1093734::analyze(), BeamThrust::calc(), InvMassFinalState::calc(), FParameter::calc(), JetShape::calc(), Spherocity::calc(), Thrust::calc(), Sphericity::calc(), FastJets::calc(), Rivet::deltaEta(), Rivet::deltaPhi(), Rivet::deltaR(), ATLAS_2012_I1091481::getPionEnergy(), ATLAS_2012_I1091481::getSeta(), CMS_2013_I1258128::makePhotonCut(), CMS_2013_I1258128::makeZCut(), ALEPH_1996_S3196992::particleInJet(), DISLepton::project(), DISKinematics::project(), DISFinalState::project(), DressedLeptons::project(), MissingMomentum::project(), Hemispheres::project(), ZFinder::project(), WFinder::project(), and Particle::setMomentum().

                                         {
      return _momentum;
    }
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(); }
operator const GenParticle * ( ) const [inline]

Cast operator for conversion to GenParticle*.

Definition at line 59 of file Particle.hh.

References Particle::genParticle().

{ return genParticle(); }
PdgId pdgId ( ) const [inline]

This Particle's PDG ID code (alias). pid/abspid form is nicer (don't need to remember lower/upper case, doesn't visually stick out in the code, etc. ...)

Definition at line 84 of file Particle.hh.

References Particle::_id.

Referenced by Rivet::PID::abspid(), BABAR_2003_I593379::analyze(), BABAR_2005_S6181155::analyze(), BABAR_2007_S6895344::analyze(), BELLE_2001_S4598261::analyze(), BELLE_2006_S6265367::analyze(), CLEO_2004_S5809304::analyze(), ARGUS_1993_S2653028::analyze(), OPAL_1994_S2927284::analyze(), PDG_HADRON_MULTIPLICITIES::analyze(), PDG_HADRON_MULTIPLICITIES_RATIOS::analyze(), OPAL_1993_S2692198::analyze(), TOTEM_2012_002::analyze(), SLD_2004_S5693039::analyze(), CMS_2010_S8656010::analyze(), OPAL_1998_S3780481::analyze(), CMS_2012_PAS_QCD_11_010::analyze(), LHCF_2012_I1115479::analyze(), ALICE_2011_S8945144::analyze(), ALICE_2011_S8909580::analyze(), SLD_1999_S3743934::analyze(), ARGUS_1993_S2669951::analyze(), OPAL_1995_S3198391::analyze(), OPAL_1997_S3608263::analyze(), OPAL_2000_S4418603::analyze(), ALEPH_2002_S4823664::analyze(), H1_1994_S2919893::analyze(), ALICE_2012_I1181770::analyze(), OPAL_1996_S3257789::analyze(), OPAL_1998_S3702294::analyze(), DELPHI_1999_S3960137::analyze(), BABAR_2007_S7266081::analyze(), DELPHI_1995_S3137023::analyze(), CMS_2010_S8547297::analyze(), ATLAS_2011_S9035664::analyze(), CDF_1993_S2742446::analyze(), CMS_2011_S8973270::analyze(), OPAL_2002_S5361494::analyze(), DELPHI_2000_S4328825::analyze(), SLD_1996_S3398250::analyze(), CDF_2008_S7540469::analyze(), OPAL_1998_S3749908::analyze(), LHCB_2011_I919315::analyze(), MC_WINC::analyze(), OPAL_1997_S3396100::analyze(), MC_IDENTIFIED::analyze(), STAR_2006_S6500200::analyze(), CMS_2011_S8978280::analyze(), STAR_2006_S6860818::analyze(), CMS_2011_S8884919::analyze(), ATLAS_2012_I1082009::analyze(), LHCB_2010_S8758301::analyze(), ARGUS_1993_S2789213::analyze(), ALEPH_1999_S4193598::analyze(), CDF_2008_S7541902::analyze(), ATLAS_2012_I1183818::analyze(), LHCB_2011_I917009::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1126136::analyze(), STAR_2008_S7869363::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_CONF_2012_109::analyze(), LHCB_2012_I1119400::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_153::analyze(), MC_SUSY::analyze(), ATLAS_2011_I944826::analyze(), ALEPH_1996_S3486095::analyze(), DELPHI_1996_S3430090::analyze(), InvMassFinalState::calc(), Particle::charge(), Jet::containsBottom(), Jet::containsCharm(), Jet::containsParticleId(), Rivet::hadronFilter(), Jet::hadronicEnergy(), Rivet::PID::hasBottom(), Particle::hasBottom(), Rivet::PID::hasCharm(), Particle::hasCharm(), Rivet::PID::hasDown(), Rivet::PID::hasStrange(), Rivet::PID::hasTop(), Rivet::PID::hasUp(), LeadingParticlesFinalState::inList(), Rivet::PID::isBaryon(), Particle::isBaryon(), Rivet::PID::isDiQuark(), Rivet::PID::isHadron(), Particle::isHadron(), Rivet::isInvisibleFilter(), Rivet::PID::isLepton(), Particle::isLepton(), Rivet::PID::isMeson(), Particle::isMeson(), Particle::isNeutrino(), Rivet::PID::isNucleus(), Rivet::PID::isPentaquark(), Rivet::PID::isRhadron(), Rivet::PID::isSUSY(), Rivet::PID::jSpin(), Rivet::PID::lSpin(), Jet::neutralEnergy(), Rivet::nonHadronFilter(), ChargedLeptons::project(), DISLepton::project(), PrimaryHadrons::project(), NeutralFinalState::project(), HeavyHadrons::project(), ZFinder::project(), IdentifiedFinalState::project(), WFinder::project(), FinalState::project(), VetoedFinalState::project(), Rivet::PID::sSpin(), Rivet::PID::threeCharge(), and Particle::threeCharge().

{ return _id; }
double perp ( ) const [inline, inherited]

Get the $ p_T $ directly (alias).

Definition at line 51 of file ParticleBase.hh.

References ParticleBase::pt().

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

Get the $ \phi $ directly.

Definition at line 80 of file ParticleBase.hh.

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

Referenced by CDF_2001_S4751469::analyze().

{ return momentum().phi(); }
PdgId pid ( ) const [inline]

This Particle's PDG ID code.

Definition at line 78 of file Particle.hh.

References Particle::_id.

Referenced by Particle::fromBottom(), Particle::fromCharm(), Particle::fromDecay(), and Particle::fromTau().

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

Get the $ \eta $ directly.

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

References ParticleBase::pt().

Referenced by 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(), ATLAS_2011_S8994773::analyze(), MC_ZINC::analyze(), CDF_1988_S1865951::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CDF_2008_S7828950::analyze(), ZEUS_2001_S4815815::analyze(), ALICE_2010_S8706239::analyze(), ATLAS_2012_I1188891::analyze(), CDF_2008_S7540469::analyze(), ATLAS_2011_I930220::analyze(), CDF_2012_NOTE10874::analyze(), CMS_2011_S8973270::analyze(), MC_WINC::analyze(), MC_LEADJETUE::analyze(), MC_DIJET::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(), ATLAS_2010_S8817804::analyze(), D0_2008_S7662670::analyze(), MC_WWJETS::analyze(), MC_ZZJETS::analyze(), CMS_2011_S9120041::analyze(), STAR_2009_UE_HELEN::analyze(), CMS_2012_I1107658::analyze(), UA1_1990_S2044935::analyze(), LHCB_2013_I1218996::analyze(), ATLAS_2012_I1082936::analyze(), MC_PHOTONJETUE::analyze(), CDF_2010_S8591881_QCD::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1190891::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), STAR_2008_S7869363::analyze(), MC_TTBAR::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_CONF_2012_103::analyze(), CDF_2001_S4751469::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2011_I926145::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2012_I1083318::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(), ATLAS_2011_S9225137::analyze(), MC_VH2BB::analyze(), ATLAS_2011_S9041966::analyze(), ATLAS_2012_I1094568::analyze(), ATLAS_2011_I944826::analyze(), ATLAS_2011_I919017::analyze(), JetShape::calc(), ATLAS_2012_I1084540::fillMap(), ATLAS_2010_S8918562::fillPtEtaNch(), STARRandomFilter::operator()(), LeadingParticlesFinalState::project(), FinalState::project(), and VetoedFinalState::project().

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

Get the $ y $ directly (alias).

Definition at line 73 of file ParticleBase.hh.

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

{ return momentum().rapidity(); }
Particle& setMomentum ( const FourMomentum momentum) [inline]

Set the momentum.

Definition at line 66 of file Particle.hh.

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

Referenced by ClusteredLepton::addPhoton(), and DISFinalState::project().

                                                        {
      _momentum = momentum;
      return *this;
    }
int threeCharge ( ) const [inline]

Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).

Definition at line 97 of file Particle.hh.

References Particle::pdgId().

                            {
      return PID::threeCharge(pdgId());
    }

Member Data Documentation

PdgId _id [private]

The PDG ID code for this Particle.

Definition at line 191 of file Particle.hh.

Referenced by Particle::abspid(), Particle::pdgId(), and Particle::pid().

The momentum of this projection of the Particle.

Definition at line 194 of file Particle.hh.

Referenced by Particle::momentum(), and Particle::setMomentum().

const GenParticle* _original [private]

A pointer to the original GenParticle from which this Particle is projected.

Definition at line 188 of file Particle.hh.

Referenced by Particle::genParticle().


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