rivet is hosted by Hepforge, IPPP Durham
DressedLepton Class Reference

A charged lepton meta-particle created by clustering photons close to the bare lepton. More...

#include <DressedLeptons.hh>

Inheritance diagram for DressedLepton:
Collaboration diagram for DressedLepton:

List of all members.

Public Member Functions

 DressedLepton (const Particle &lepton)
void addPhoton (const Particle &p, bool cluster)
const ParticleconstituentLepton () const
const ParticlesconstituentPhotons () const
virtual fastjet::PseudoJet pseudojet () const
 Converter to FastJet3 PseudoJet.
 operator PseudoJet () const
 Cast operator to FastJet3 PseudoJet.
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.
Decay info
bool isStable () const
 Whether this particle is stable according to the generator.
vector< Particlechildren () const
 Get a list of the direct descendants from the current particle.
vector< ParticleallDescendants () const
vector< ParticlestableDescendants () const
double flightLength () const
 Flight length (divide by mm or cm to get the appropriate units)
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

Particles _constituentPhotons
Particle _constituentLepton

Detailed Description

A charged lepton meta-particle created by clustering photons close to the bare lepton.

Definition at line 18 of file DressedLeptons.hh.


Constructor & Destructor Documentation

DressedLepton ( const Particle lepton) [inline]

Definition at line 21 of file DressedLeptons.hh.

                                          :
      Particle(lepton.pid(), lepton.momentum()),
      _constituentLepton(lepton) {}

Member Function Documentation

PdgId abspid ( ) const [inline, inherited]

Absolute value of the PDG ID code.

Definition at line 90 of file Particle.hh.

References Particle::_id.

Referenced by BABAR_2005_S6181155::analyze(), BABAR_2007_S6895344::analyze(), CLEO_2004_S5809304::analyze(), BELLE_2006_S6265367::analyze(), BABAR_2013_I1238276::analyze(), BELLE_2013_I1216515::analyze(), OPAL_1994_S2927284::analyze(), PDG_HADRON_MULTIPLICITIES::analyze(), PDG_HADRON_MULTIPLICITIES_RATIOS::analyze(), BELLE_2008_I786560::analyze(), SLD_2004_S5693039::analyze(), CMS_2012_PAS_QCD_11_010::analyze(), SLD_1999_S3743934::analyze(), OPAL_1997_S3608263::analyze(), OPAL_2000_S4418603::analyze(), ARGUS_1993_S2669951::analyze(), OPAL_1995_S3198391::analyze(), ALEPH_2002_S4823664::analyze(), OPAL_1996_S3257789::analyze(), OPAL_1998_S3702294::analyze(), DELPHI_1999_S3960137::analyze(), BABAR_2007_S7266081::analyze(), DELPHI_1995_S3137023::analyze(), ATLAS_2011_S9035664::analyze(), CMS_2011_S8973270::analyze(), OPAL_1998_S3749908::analyze(), LHCB_2011_I919315::analyze(), OPAL_1997_S3396100::analyze(), MC_IDENTIFIED::analyze(), CMS_2011_S8978280::analyze(), ATLAS_2012_I1082009::analyze(), ARGUS_1993_S2789213::analyze(), ALEPH_1999_S4193598::analyze(), CDF_2008_S7541902::analyze(), LHCB_2013_I1218996::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2013_I1190187::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2012_CONF_2012_153::analyze(), ATLAS_2011_I944826::analyze(), DELPHI_1996_S3430090::analyze(), ATLAS_2012_I1204447::get_tau_neutrino_mom(), Rivet::identifyZstates(), PromptFinalState::isPrompt(), HasPID::operator()(), and TauFinder::project().

{ return std::abs(_id); }
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().

{ return momentum().absrapidity(); }
void addPhoton ( const Particle p,
bool  cluster 
) [inline]

Definition at line 25 of file DressedLeptons.hh.

References DressedLepton::_constituentPhotons, Particle::momentum(), and Particle::setMomentum().

                                                    {
      _constituentPhotons.push_back(p);
      if (cluster) setMomentum(momentum() + p.momentum());
    }
vector<Particle> allDescendants ( ) const [inline, inherited]

Get a list of all the descendants (including duplication of parents and children) from the current particle

Todo:

Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates

Insist that the current particle is post-hadronization, otherwise throw an exception?

Todo:
Remove this const mess crap when HepMC doesn't suck
Todo:
Would like to do this, but the range objects are broken

Definition at line 220 of file Particle.hh.

References Particle::genParticle(), Particle::isStable(), and Particle::Particle().

                                            {
      vector<Particle> rtn;
      if (isStable()) return rtn;
      /// @todo Remove this const mess crap when HepMC doesn't suck
      HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
      /// @todo Would like to do this, but the range objects are broken
      // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
      //   rtn += Particle(gp);
      for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
        rtn += Particle(*it);
      return rtn;
    }
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); }
double charge ( ) const [inline, inherited]

The charge of this Particle.

Definition at line 103 of file Particle.hh.

References Particle::pid().

Referenced by ATLAS_2014_I1282441::analyze(), NeutralFinalState::project(), and PromptFinalState::project().

                          {
      return PID::charge(pid());
    }
vector<Particle> children ( ) const [inline, inherited]

Get a list of the direct descendants from the current particle.

Todo:
Remove this const mess crap when HepMC doesn't suck
Todo:
Would like to do this, but the range objects are broken

Definition at line 204 of file Particle.hh.

References Particle::genParticle(), Particle::isStable(), and Particle::Particle().

                                      {
      vector<Particle> rtn;
      if (isStable()) return rtn;
      /// @todo Remove this const mess crap when HepMC doesn't suck
      HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
      /// @todo Would like to do this, but the range objects are broken
      // foreach (const GenParticle* gp, gv->particles(HepMC::children))
      //   rtn += Particle(gp);
      for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::children); it != gv->particles_end(HepMC::children); ++it)
        rtn += Particle(*it);
      return rtn;
    }
const Particles& constituentPhotons ( ) const [inline]
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_I1087342::analyze(), CMS_2012_I1184941::analyze(), CMS_2010_S8656010::analyze(), ATLAS_2010_S8591806::analyze(), STAR_2008_S7993412::analyze(), CMSTOTEM_2014_I1294140::analyze(), ALICE_2012_I1181770::analyze(), STAR_2006_S6870392::analyze(), H1_1994_S2919893::analyze(), CMS_QCD_10_024::analyze(), CDF_2005_S6080774::analyze(), CMS_2010_S8547297::analyze(), D0_1996_S3324664::analyze(), ATLAS_2011_S8994773::analyze(), MC_DIPHOTON::analyze(), MC_HINC::analyze(), CDF_1990_S2089246::analyze(), D0_2008_S6879055::analyze(), CDF_1993_S2742446::analyze(), MC_ZINC::analyze(), CDF_2009_S8436959::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_WINC::analyze(), ATLAS_2012_I1119557::analyze(), MC_DIJET::analyze(), ATLAS_2014_I1298811::analyze(), CMS_2011_S8884919::analyze(), ATLAS_2012_I946427::analyze(), LHCB_2013_I1208105::analyze(), ATLAS_2014_I1307756::analyze(), UA1_1990_S2044935::analyze(), CMS_2013_I1218372::analyze(), ATLAS_2012_I1183818::analyze(), ATLAS_2012_I1199269::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1190891::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9120807::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2011_CONF_2011_098::analyze(), D0_2008_S7719523::analyze(), ATLAS_2010_S8914702::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1095236::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(), ATLAS_2011_S9041966::analyze(), MC_VH2BB::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(), FinalState::project(), TriggerCDFRun2::project(), TriggerCDFRun0Run1::project(), NeutralFinalState::project(), and TriggerUA5::project().

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

Flight length (divide by mm or cm to get the appropriate units)

Definition at line 252 of file Particle.hh.

References Particle::genParticle(), Particle::isStable(), and Rivet::sqr().

                                {
      if (isStable()) return -1;
      if (genParticle() == NULL) return 0;
      if (genParticle()->production_vertex() == NULL) return 0;
      const HepMC::FourVector v1 = genParticle()->production_vertex()->position();
      const HepMC::FourVector v2 = genParticle()->end_vertex()->position();
      return sqrt(sqr(v2.x()-v1.x()) + sqr(v2.y()-v1.y()) + sqr(v2.z()-v1.z()));
    }
bool fromBottom ( ) const [inherited]

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 [inherited]

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 [inherited]

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 MC_ParticleAnalysis::_analyze(), and 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 [inherited]

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;
  }
const GenParticle* genParticle ( ) const [inline, inherited]

Get a const reference to the original GenParticle.

Definition at line 64 of file Particle.hh.

References Particle::_original.

Referenced by FinalState::accept(), Particle::allDescendants(), BELLE_2001_S4598261::analyze(), BABAR_2003_I593379::analyze(), BELLE_2006_S6265367::analyze(), BABAR_2005_S6181155::analyze(), BABAR_2013_I1238276::analyze(), ARGUS_1993_S2653028::analyze(), H1_2000_S4129130::analyze(), BELLE_2008_I786560::analyze(), ARGUS_1993_S2669951::analyze(), H1_1994_S2919893::analyze(), BABAR_2007_S7266081::analyze(), ATLAS_2011_S8994773::analyze(), ATLAS_2011_S9035664::analyze(), CDF_2008_S7540469::analyze(), D0_2010_S8570965::analyze(), CMS_2013_I1256943::analyze(), ARGUS_1993_S2789213::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2010_S8894728::analyze(), InvMassFinalState::calc(), Particle::children(), Jet::containsBottom(), Jet::containsCharm(), Jet::containsParticle(), ATLAS_2011_I944826::daughtersSurviveCuts(), Particle::flightLength(), Particle::fromBottom(), Particle::fromCharm(), Particle::fromDecay(), Particle::fromTau(), ATLAS_2012_I1204447::get_tau_neutrino_mom(), LHCB_2010_S8758301::getLongestLivedAncestor(), LHCB_2011_I917009::getMotherLifeTimeSum(), LHCB_2012_I1119400::getMotherLifeTimeSum(), ATLAS_2011_I944826::getPerpFlightDistance(), Particle::hasAncestor(), hasDecayedTo(), PromptFinalState::isPrompt(), Particle::isStable(), Particle::operator const GenParticle *(), DISLepton::project(), MergedFinalState::project(), DISFinalState::project(), PrimaryHadrons::project(), HeavyHadrons::project(), VetoedFinalState::project(), and Particle::stableDescendants().

                                           {
      return _original;
    }
bool hasAncestor ( PdgId  pdg_id) const [inherited]

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(), ATLAS_2014_I1282441::analyze(), CMS_2011_S8978280::analyze(), and ATLAS_2012_I1204447::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, inherited]

Does this (hadron) contain a b quark?

Definition at line 127 of file Particle.hh.

References Particle::pid().

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

{ return PID::hasBottom(pid()); }
bool hasCharm ( ) const [inline, inherited]

Does this (hadron) contain a c quark?

Definition at line 130 of file Particle.hh.

References Particle::pid().

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

{ return PID::hasCharm(pid()); }
bool isBaryon ( ) const [inline, inherited]

Is this a baryon?

Definition at line 118 of file Particle.hh.

References Particle::pid().

{ return PID::isBaryon(pid()); }
bool isHadron ( ) const [inline, inherited]

Is this a hadron?

Definition at line 112 of file Particle.hh.

References Particle::pid().

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

{ return PID::isHadron(pid()); }
bool isLepton ( ) const [inline, inherited]

Is this a lepton?

Definition at line 121 of file Particle.hh.

References Particle::pid().

{ return PID::isLepton(pid()); }
bool isMeson ( ) const [inline, inherited]

Is this a meson?

Definition at line 115 of file Particle.hh.

References Particle::pid().

{ return PID::isMeson(pid()); }
bool isNeutrino ( ) const [inline, inherited]

Is this a neutrino?

Definition at line 124 of file Particle.hh.

References Particle::pid().

Referenced by ATLAS_2014_I1307756::analyze().

{ return PID::isNeutrino(pid()); }
bool isStable ( ) const [inline, inherited]

Whether this particle is stable according to the generator.

Definition at line 199 of file Particle.hh.

References Particle::genParticle().

Referenced by Particle::allDescendants(), Particle::children(), Particle::flightLength(), and Particle::stableDescendants().

                          {
      return genParticle() != NULL && genParticle()->status() == 1 && genParticle()->end_vertex() == NULL;
    }
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 Hemispheres::calc(), Particle::pseudojet(), Jet::setState(), and Variables::Variables().

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

The momentum.

Implements ParticleBase.

Definition at line 72 of file Particle.hh.

References Particle::_momentum.

Referenced by CDF_2004_S5839831::_calcTransCones(), DressedLepton::addPhoton(), BABAR_2003_I593379::analyze(), BELLE_2006_S6265367::analyze(), BABAR_2007_S6895344::analyze(), CLEO_2004_S5809304::analyze(), BELLE_2001_S4598261::analyze(), BABAR_2005_S6181155::analyze(), ARGUS_1993_S2653028::analyze(), BABAR_2013_I1238276::analyze(), BELLE_2013_I1216515::analyze(), ATLAS_2011_I894867::analyze(), CMS_2012_I1193338::analyze(), OPAL_1993_S2692198::analyze(), H1_2000_S4129130::analyze(), BELLE_2008_I786560::analyze(), ARGUS_1993_S2669951::analyze(), ALICE_2012_I1181770::analyze(), H1_1994_S2919893::analyze(), ATLAS_2011_S9002537::analyze(), BABAR_2007_S7266081::analyze(), CDF_2005_S6080774::analyze(), EXAMPLE_CUTS::analyze(), MC_DIPHOTON::analyze(), CDF_2009_S8436959::analyze(), CDF_2008_S7540469::analyze(), MC_PHOTONINC::analyze(), MC_PHOTONKTSPLITTINGS::analyze(), MC_WINC::analyze(), MC_HFJETS::analyze(), MC_PHOTONJETS::analyze(), LHCB_2010_S8758301::analyze(), MC_PHOTONS::analyze(), CMS_2013_I1256943::analyze(), ATLAS_2012_I946427::analyze(), LHCB_2011_I917009::analyze(), ARGUS_1993_S2789213::analyze(), ATLAS_2013_I1243871::analyze(), ATLAS_2014_I1307756::analyze(), CDF_2008_S7541902::analyze(), ATLAS_2012_I1183818::analyze(), ATLAS_2012_I1199269::analyze(), ATLAS_2011_I921594::analyze(), ATLAS_2013_I1263495::analyze(), ATLAS_2011_S9108483::analyze(), ATLAS_2012_CONF_2012_104::analyze(), ATLAS_2012_I1126136::analyze(), ATLAS_2012_I1190891::analyze(), LHCB_2012_I1119400::analyze(), ATLAS_2012_I1186556::analyze(), ATLAS_2011_S9120807::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), ATLAS_2013_I1190187::analyze(), MC_TTBAR::analyze(), ATLAS_2011_CONF_2011_098::analyze(), ATLAS_2012_I1112263::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2010_S8914702::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2012_CONF_2012_001::analyze(), ATLAS_2012_I1180197::analyze(), ATLAS_2012_CONF_2012_109::analyze(), ATLAS_2011_S8983313::analyze(), ATLAS_2011_S9019561::analyze(), ATLAS_2011_CONF_2011_090::analyze(), ATLAS_2014_I1304688::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(), CDF_2004_S5839831::analyze(), ATLAS_2012_I1203852::analyze(), ATLAS_2012_I1093734::analyze(), BeamThrust::calc(), InvMassFinalState::calc(), FParameter::calc(), JetShape::calc(), Spherocity::calc(), Thrust::calc(), Sphericity::calc(), FastJets::calc(), ATLAS_2012_I1204447::isonZ(), CMS_2013_I1258128::makePhotonCut(), CMS_2013_I1258128::makeZCut(), DISLepton::project(), DISKinematics::project(), DISFinalState::project(), MissingMomentum::project(), DressedLeptons::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, inherited]

Cast operator for conversion to GenParticle*.

Definition at line 69 of file Particle.hh.

References Particle::genParticle().

{ return genParticle(); }
operator PseudoJet ( ) const [inline, inherited]

Cast operator to FastJet3 PseudoJet.

Definition at line 55 of file Particle.hh.

References Particle::pseudojet().

{ return pseudojet(); }
PdgId pdgId ( ) const [inline, inherited]

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 94 of file Particle.hh.

References Particle::_id.

Referenced by ATLAS_2014_I1298811::analyze(), and FinalState::project().

{ return _id; }
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]
PdgId pid ( ) const [inline, inherited]

This Particle's PDG ID code.

Definition at line 88 of file Particle.hh.

References Particle::_id.

Referenced by BABAR_2003_I593379::analyze(), BELLE_2001_S4598261::analyze(), ARGUS_1993_S2653028::analyze(), OPAL_1993_S2692198::analyze(), TOTEM_2012_002::analyze(), SLD_2004_S5693039::analyze(), BELLE_2008_I786560::analyze(), CMS_2010_S8656010::analyze(), OPAL_1998_S3780481::analyze(), LHCF_2012_I1115479::analyze(), SLD_1999_S3743934::analyze(), ALICE_2011_S8909580::analyze(), ALICE_2011_S8945144::analyze(), ARGUS_1993_S2669951::analyze(), ALICE_2012_I1181770::analyze(), H1_1994_S2919893::analyze(), BABAR_2007_S7266081::analyze(), CMS_2010_S8547297::analyze(), CDF_1993_S2742446::analyze(), SLD_1996_S3398250::analyze(), CDF_2008_S7540469::analyze(), DELPHI_2000_S4328825::analyze(), OPAL_2002_S5361494::analyze(), MC_WINC::analyze(), STAR_2006_S6500200::analyze(), MC_HFJETS::analyze(), MC_IDENTIFIED::analyze(), STAR_2006_S6860818::analyze(), CMS_2011_S8884919::analyze(), LHCB_2010_S8758301::analyze(), CMS_2013_I1256943::analyze(), LHCB_2011_I917009::analyze(), ARGUS_1993_S2789213::analyze(), ATLAS_2012_I1183818::analyze(), ATLAS_2012_I1126136::analyze(), LHCB_2012_I1119400::analyze(), STAR_2008_S7869363::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_CONF_2012_153::analyze(), MC_SUSY::analyze(), ALEPH_1996_S3486095::analyze(), ATLAS_2011_I944826::analyze(), ATLAS_2012_I1203852::analyze(), InvMassFinalState::calc(), Particle::charge(), Jet::containsBottom(), Jet::containsCharm(), Jet::containsParticleId(), Particle::fromBottom(), Particle::fromCharm(), Particle::fromDecay(), Particle::fromTau(), Rivet::hadronFilter(), Jet::hadronicEnergy(), Particle::hasBottom(), Particle::hasCharm(), Rivet::identifyZstates(), LeadingParticlesFinalState::inList(), Particle::isBaryon(), Particle::isHadron(), Rivet::isInvisibleFilter(), Particle::isLepton(), Particle::isMeson(), Particle::isNeutrino(), ATLAS_2012_I1204447::isonZ(), Jet::neutralEnergy(), Rivet::nonHadronFilter(), HasPID::operator()(), ChargedLeptons::project(), DISLepton::project(), NeutralFinalState::project(), PrimaryHadrons::project(), PromptFinalState::project(), HeavyHadrons::project(), ZFinder::project(), WFinder::project(), IdentifiedFinalState::project(), VetoedFinalState::project(), and Particle::threeCharge().

{ return _id; }
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(); }
virtual fastjet::PseudoJet pseudojet ( ) const [inline, virtual, inherited]

Converter to FastJet3 PseudoJet.

Definition at line 50 of file Particle.hh.

References ParticleBase::E(), ParticleBase::mom(), ParticleBase::px(), ParticleBase::py(), and ParticleBase::pz().

Referenced by Particle::operator PseudoJet().

                                               {
      return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
    }
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(), TOTEM_2012_002::analyze(), CMS_2012_I1087342::analyze(), CMS_2010_S8656010::analyze(), CMS_2012_PAS_QCD_11_010::analyze(), CDF_2006_S6450792::analyze(), CMS_2011_S9088458::analyze(), ATLAS_2010_S8591806::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(), ATLAS_2014_I1282441::analyze(), CMS_2010_S8547297::analyze(), ATLAS_2011_S8994773::analyze(), CMS_2013_I1273574::analyze(), MC_HINC::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(), MC_ZINC::analyze(), CMS_2011_S8973270::analyze(), ATLAS_2011_I930220::analyze(), CDF_2008_S7540469::analyze(), CDF_2012_NOTE10874::analyze(), SFM_1984_S1178091::analyze(), MC_WINC::analyze(), CDF_1996_S3108457::analyze(), MC_DIJET::analyze(), MC_HFJETS::analyze(), MC_LEADJETUE::analyze(), STAR_2006_S6500200::analyze(), ATLAS_2014_I1268975::analyze(), ATLAS_2014_I1298811::analyze(), CDF_2009_S8233977::analyze(), CMS_2011_S8978280::analyze(), STAR_2006_S6860818::analyze(), CMS_2011_S8884919::analyze(), ATLAS_2010_S8817804::analyze(), D0_2008_S7662670::analyze(), CMS_2011_S9120041::analyze(), MC_PHOTONS::analyze(), MC_WWJETS::analyze(), MC_ZZJETS::analyze(), STAR_2009_UE_HELEN::analyze(), LHCB_2013_I1208105::analyze(), ATLAS_2013_I1243871::analyze(), CMS_2012_I1107658::analyze(), UA1_1990_S2044935::analyze(), LHCB_2013_I1218996::analyze(), ATLAS_2012_I1082936::analyze(), MC_PHOTONJETUE::analyze(), ATLAS_2012_I1124167::analyze(), ATLAS_2011_I921594::analyze(), CDF_2010_S8591881_QCD::analyze(), ATLAS_2013_I1263495::analyze(), ATLAS_2012_CONF_2012_104::analyze(), MC_GENERIC::analyze(), ATLAS_2012_I1190891::analyze(), ATLAS_2012_I1186556::analyze(), STAR_2008_S7869363::analyze(), ATLAS_2011_S9212183::analyze(), ATLAS_2012_CONF_2012_105::analyze(), MC_TTBAR::analyze(), ATLAS_2013_I1190187::analyze(), ATLAS_2012_I1112263::analyze(), CDF_2001_S4751469::analyze(), ATLAS_2012_CONF_2012_103::analyze(), ATLAS_2012_I1117704::analyze(), ATLAS_2012_I1095236::analyze(), ATLAS_2012_I1125961::analyze(), ATLAS_2012_I1204447::analyze(), ATLAS_2011_I926145::analyze(), CMS_2013_I1224539_WJET::analyze(), CMS_2013_I1224539_ZJET::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_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(), ATLAS_2011_S9041966::analyze(), MC_VH2BB::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()(), 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().

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

Set the momentum.

Definition at line 76 of file Particle.hh.

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

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

                                                        {
      _momentum = momentum;
      return *this;
    }
vector<Particle> stableDescendants ( ) const [inline, inherited]

Get a list of all the stable descendants from the current particle

Todo:

Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates

Insist that the current particle is post-hadronization, otherwise throw an exception?

Todo:
Remove this const mess crap when HepMC doesn't suck
Todo:
Would like to do this, but the range objects are broken

Definition at line 236 of file Particle.hh.

References Particle::genParticle(), Particle::isStable(), and Particle::Particle().

                                               {
      vector<Particle> rtn;
      if (isStable()) return rtn;
      /// @todo Remove this const mess crap when HepMC doesn't suck
      HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
      /// @todo Would like to do this, but the range objects are broken
      // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
      //   if (gp->status() == 1 && gp->end_vertex() == NULL)
      //     rtn += Particle(gp);
      for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
        if ((*it)->status() == 1 && (*it)->end_vertex() == NULL)
          rtn += Particle(*it);
      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(); }
int threeCharge ( ) const [inline, inherited]

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

Definition at line 107 of file Particle.hh.

References Particle::pid().

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

Member Data Documentation

Definition at line 36 of file DressedLeptons.hh.

Referenced by DressedLepton::constituentLepton().


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