rivet is hosted by Hepforge, IPPP Durham

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

#include <DressedLeptons.hh>

List of all members.

Public Member Functions

 DressedLepton (const Particle &lepton)
void addPhoton (const Particle &p, bool cluster)
const ParticleconstituentLepton () const
const ParticlesconstituentPhotons () const
Other representations and implicit casts
virtual fastjet::PseudoJet pseudojet () const
 Converter to FastJet3 PseudoJet.
 operator PseudoJet () const
 Cast operator to FastJet3 PseudoJet.
const GenParticle * genParticle () const
 Get a const pointer to the original GenParticle.
 operator const GenParticle * () const
 Cast operator for conversion to GenParticle*.
Kinematic properties
const FourMomentummomentum () const
 The momentum.
ParticlesetMomentum (const FourMomentum &momentum)
 Set the momentum.
ParticlesetMomentum (double E, double px, double py, double pz)
 Set the momentum via components.
ParticletransformBy (const LorentzTransform &lt)
 Apply an active Lorentz transform to this particle.
Positional properties
const FourVectororigin () const
 The origin position.
ParticlesetOrigin (const FourVector &position)
 Set the origin position.
ParticlesetOrigin (double t, double x, double y, double z)
 Set the origin position via components.
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.
double abscharge () const
 The absolute charge of this Particle.
int charge3 () const
 Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
int threeCharge () const
int abscharge3 () const
 Three times the absolute 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?
bool isVisible () const
 Is this particle potentially visible in a detector?
Ancestry properties
Particles parents (const Cut &c=Cuts::OPEN) const
template<typename FN >
Particles parents (const FN &f) const
bool hasParent (PdgId pid) const
template<typename FN >
bool hasParentWith (const FN &f) const
bool hasParentWith (const Cut &c) const
bool hasAncestor (PdgId pid) const
template<typename FN >
bool hasAncestorWith (const FN &f) const
bool hasAncestorWith (const Cut &c) const
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 fromHadron () const
 Determine whether the particle is from a hadron decay.
bool fromTau (bool prompt_taus_only=false) const
 Determine whether the particle is from a tau decay.
bool fromPromptTau () const
 Determine whether the particle is from a prompt tau decay.
bool fromDecay () const
 Determine whether the particle is from a hadron or tau decay.
bool isPrompt (bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const
 Shorthand definition of 'promptness' based on set definition flags.
Decay info
bool isStable () const
 Whether this particle is stable according to the generator.
Particles children (const Cut &c=Cuts::OPEN) const
 Get a list of the direct descendants from the current particle (with optional selection Cut)
template<typename FN >
Particles children (const FN &f) const
 Get a list of the direct descendants from the current particle (with selector function)
Particles allDescendants (const Cut &c=Cuts::OPEN, bool remove_duplicates=true) const
 Get a list of all the descendants from the current particle (with optional selection Cut)
template<typename FN >
Particles allDescendants (const FN &f, bool remove_duplicates=true) const
 Get a list of all the descendants from the current particle (with selector function)
Particles stableDescendants (const Cut &c=Cuts::OPEN) const
template<typename FN >
Particles stableDescendants (const FN &f) const
 Get a list of all the stable descendants from the current particle (with selector function)
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 E2 () const
 Get the energy-squared.
double energy2 () const
 Get the energy-squared (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 Et2 () const
 Get the $ E_T^2 $ 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 p () const
 Get the 3-momentum magnitude directly.
double p2 () const
 Get the 3-momentum magnitude-squared directly.
double px () const
 x component of momentum.
double py () const
 y component of momentum.
double pz () const
 z component of momentum.
double px2 () const
 x component of momentum, squared.
double py2 () const
 y component of momentum, squared.
double pz2 () const
 z component of momentum, squared.
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)

Protected Member Functions

template<typename FN >
bool _hasRelativeWith (HepMC::IteratorRange relation, const FN &f) const

Protected 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 particle.
FourVector _origin
 The creation position of this particle.

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 14 of file DressedLeptons.hh.


Constructor & Destructor Documentation

DressedLepton ( const Particle lepton) [inline]

Definition at line 17 of file DressedLeptons.hh.

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

Member Function Documentation

bool _hasRelativeWith ( HepMC::IteratorRange  relation,
const FN &  f 
) const [inline, protected, inherited]

Definition at line 371 of file Particle.hh.

                                                                          {
      for (const GenParticle* ancestor : particles(genParticle(), relation)) {
        if (f(Particle(ancestor))) return true;
      }
      return false;
    }
double abscharge ( ) const [inline, inherited]

The absolute charge of this Particle.

Definition at line 153 of file Particle.hh.

{ return PID::abscharge(pid()); }
int abscharge3 ( ) const [inline, inherited]

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

Definition at line 163 of file Particle.hh.

{ return PID::abscharge3(pid()); }
double abseta ( ) const [inline, inherited]

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

Definition at line 81 of file ParticleBase.hh.

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

Absolute value of the PDG ID code.

Definition at line 138 of file Particle.hh.

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

Get the $ |\eta| $ directly.

Definition at line 79 of file ParticleBase.hh.

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

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

Definition at line 90 of file ParticleBase.hh.

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

Get the $ |y| $ directly.

Definition at line 88 of file ParticleBase.hh.

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

Definition at line 21 of file DressedLeptons.hh.

                                                    {
      _constituentPhotons.push_back(p);
      if (cluster) setMomentum(momentum() + p.momentum());
    }
vector< Particle > allDescendants ( const Cut c = Cuts::OPEN,
bool  remove_duplicates = true 
) const [inherited]

Get a list of all the descendants from the current particle (with optional selection Cut)

Todo:

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

Use recursion through replica-avoiding functions to avoid bookkeeping duplicates

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

Definition at line 72 of file Particle.cc.

                                                                                      {
    vector<Particle> rtn;
    if (isStable()) return rtn;
    /// @todo Remove this const mess crap when HepMC doesn't suck
    GenVertexPtr gv = const_cast<GenVertexPtr>( genParticle()->end_vertex() );
    if (gv == NULL) return rtn;
    /// @todo Would like to do this, but the range objects are broken
    // foreach (const GenParticlePtr gp, gv->particles(HepMC::descendants))
    for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it) {
      const Particle p(*it);
      if (c != Cuts::OPEN && !c->accept(p)) continue;
      if (remove_duplicates && (*it)->end_vertex() != NULL) {
        // size_t n = 0; ///< @todo Only remove 1-to-1 duplicates?
        bool dup = false;
        /// @todo Yuck, HepMC
        for (GenVertex::particle_iterator it2 = (*it)->end_vertex()->particles_begin(HepMC::children); it2 != (*it)->end_vertex()->particles_end(HepMC::children); ++it2) {
          // n += 1; if (n > 1) break;
          if ((*it)->pdg_id() == (*it2)->pdg_id()) { dup = true; break; }
        }
        if (dup) continue;
      }
      rtn += p;
    }
    return rtn;
  }
Particles allDescendants ( const FN &  f,
bool  remove_duplicates = true 
) const [inline, inherited]

Get a list of all the descendants from the current particle (with selector function)

Definition at line 346 of file Particle.hh.

                                                                             {
      return filter_select(allDescendants(Cuts::OPEN, remove_duplicates), f);
    }
double angle ( const ParticleBase v) const [inline, inherited]

Angle between this vector and another.

Definition at line 125 of file ParticleBase.hh.

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

Angle between this vector and another.

Definition at line 127 of file ParticleBase.hh.

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

Angle between this vector and another (3-vector)

Definition at line 129 of file ParticleBase.hh.

{ return momentum().angle(v3); }
double azimuthalAngle ( const PhiMapping  mapping = ZERO_2PI) const [inline, inherited]

Azimuthal angle $ \phi $.

Definition at line 93 of file ParticleBase.hh.

{ return momentum().azimuthalAngle(mapping); }
double charge ( ) const [inline, inherited]

The charge of this Particle.

Definition at line 150 of file Particle.hh.

{ return PID::charge(pid()); }
int charge3 ( ) const [inline, inherited]

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

Definition at line 156 of file Particle.hh.

{ return PID::charge3(pid()); }
vector< Particle > children ( const Cut c = Cuts::OPEN) const [inherited]

Get a list of the direct descendants from the current particle (with optional selection Cut)

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 52 of file Particle.cc.

                                                        {
    vector<Particle> rtn;
    if (isStable()) return rtn;
    /// @todo Remove this const mess crap when HepMC doesn't suck
    GenVertexPtr gv = const_cast<GenVertexPtr>( genParticle()->end_vertex() );
    if (gv == NULL) return rtn;
    /// @todo Would like to do this, but the range objects are broken
    // foreach (const GenParticlePtr 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) {
      const Particle p(*it);
      if (c != Cuts::OPEN && !c->accept(p)) continue;
      rtn += p;
    }
    return rtn;
  }
Particles children ( const FN &  f) const [inline, inherited]

Get a list of the direct descendants from the current particle (with selector function)

Definition at line 337 of file Particle.hh.

                                          {
      return filter_select(children(), f);
    }
const Particle& constituentLepton ( ) const [inline]

Definition at line 26 of file DressedLeptons.hh.

{ return _constituentLepton; }
const Particles& constituentPhotons ( ) const [inline]

Definition at line 27 of file DressedLeptons.hh.

{ return _constituentPhotons; }
double E ( ) const [inline, inherited]

Get the energy directly.

Definition at line 41 of file ParticleBase.hh.

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

Get the energy-squared.

Definition at line 46 of file ParticleBase.hh.

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

Get the energy directly (alias).

Definition at line 43 of file ParticleBase.hh.

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

Get the energy-squared (alias).

Definition at line 48 of file ParticleBase.hh.

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

Get the $ E_T $ directly.

Definition at line 65 of file ParticleBase.hh.

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

Get the $ E_T^2 $ directly.

Definition at line 67 of file ParticleBase.hh.

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

Get the $ \eta $ directly (alias).

Definition at line 77 of file ParticleBase.hh.

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

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

Definition at line 119 of file Particle.cc.

                                      {
    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!

Definition at line 147 of file Particle.cc.

                                  {
    return _hasRelativeWith(HepMC::ancestors, [](const Particle& p){
        return p.genParticle()->status() == 2 && p.isHadron() && p.hasBottom();
      });
    // const GenVertexPtr prodVtx = genParticle()->production_vertex();
    // if (prodVtx == NULL) return false;
    // foreach (const GenParticlePtr 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!

Definition at line 161 of file Particle.cc.

                                 {
    return _hasRelativeWith(HepMC::ancestors, [](const Particle& p){
        return p.genParticle()->status() == 2 && p.isHadron() && p.hasCharm();
      });
    // const GenVertexPtr prodVtx = genParticle()->production_vertex();
    // if (prodVtx == NULL) return false;
    // foreach (const GenParticlePtr 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 [inline, 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!

Definition at line 313 of file Particle.hh.

{ return fromHadron() || fromPromptTau(); }
bool fromHadron ( ) const [inherited]

Determine whether the particle is from a 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!

Definition at line 175 of file Particle.cc.

                                  {
    return _hasRelativeWith(HepMC::ancestors, [](const Particle& p){
        return p.genParticle()->status() == 2 && p.isHadron();
      });
    // const GenVertexPtr prodVtx = genParticle()->production_vertex();
    // if (prodVtx == NULL) return false;
    // foreach (const GenParticlePtr ancestor, particles(prodVtx, HepMC::ancestors)) {
    //   const PdgId pid = ancestor->pdg_id();
    //   if (ancestor->status() == 2 && PID::isHadron(pid)) return true;
    // }
    // return false;
  }
bool fromPromptTau ( ) const [inline, inherited]

Determine whether the particle is from a prompt 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!

Definition at line 303 of file Particle.hh.

{ return fromTau(true); }
bool fromTau ( bool  prompt_taus_only = false) 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!

Definition at line 189 of file Particle.cc.

                                                    {
    if (prompt_taus_only && fromHadron()) return false;
    return _hasRelativeWith(HepMC::ancestors, [](const Particle& p){
        return p.genParticle()->status() == 2 && isTau(p);
      });
    // const GenVertexPtr prodVtx = genParticle()->production_vertex();
    // if (prodVtx == NULL) return false;
    // foreach (const GenParticlePtr 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 pointer to the original GenParticle.

Definition at line 75 of file Particle.hh.

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

Check whether a given PID is found in the particle'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!

Definition at line 138 of file Particle.cc.

                                            {
    return _hasRelativeWith(HepMC::ancestors, hasPID(pid));
  }
bool hasAncestorWith ( const FN &  f) const [inline, inherited]

Check whether a particle in the particle's ancestor list has the requested property

Note:
This question is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!

Definition at line 252 of file Particle.hh.

                                            {
      return _hasRelativeWith(HepMC::ancestors, f);
    }
bool hasAncestorWith ( const Cut c) const [inherited]

Definition at line 142 of file Particle.cc.

                                                   {
    return hasAncestorWith([&](const Particle& p){return c->accept(p);});
  }
bool hasBottom ( ) const [inline, inherited]

Does this (hadron) contain a b quark?

Definition at line 181 of file Particle.hh.

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

Does this (hadron) contain a c quark?

Definition at line 184 of file Particle.hh.

{ return PID::hasCharm(pid()); }
bool hasParent ( PdgId  pid) const [inherited]

Check whether a given PID is found in the particle's parent list

Note:
This question is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!
Deprecated:
Prefer e.g. parents(Cut::pid == 123).size()

Definition at line 129 of file Particle.cc.

                                          {
    return _hasRelativeWith(HepMC::parents, hasPID(pid));
  }
bool hasParentWith ( const FN &  f) const [inline, inherited]

Check whether a particle in the particle's parent list has the requested property

Note:
This question is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!
Deprecated:
Prefer parents(Cut) or parents(FN) methods and .empty()

Definition at line 234 of file Particle.hh.

                                          {
      return _hasRelativeWith(HepMC::parents, f);
    }
bool hasParentWith ( const Cut c) const [inherited]

Definition at line 133 of file Particle.cc.

                                                 {
    return hasParentWith([&](const Particle& p){return c->accept(p);});
  }
bool isBaryon ( ) const [inline, inherited]

Is this a baryon?

Definition at line 172 of file Particle.hh.

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

Is this a hadron?

Definition at line 166 of file Particle.hh.

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

Is this a lepton?

Definition at line 175 of file Particle.hh.

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

Is this a meson?

Definition at line 169 of file Particle.hh.

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

Is this a neutrino?

Definition at line 178 of file Particle.hh.

{ return PID::isNeutrino(pid()); }
bool isPrompt ( bool  allow_from_prompt_tau = false,
bool  allow_from_prompt_mu = false 
) const [inherited]

Shorthand definition of 'promptness' based on set definition flags.

The boolean arguments allow a decay lepton to be considered prompt if its parent was a "real" prompt lepton.

Note:
This one doesn't make any judgements about final-stateness
Todo:
Would be nicer to be able to write this recursively up the chain, exiting as soon as a parton or string/cluster is seen

Definition at line 215 of file Particle.cc.

                                                                                     {
    if (genParticle() == NULL) return false; // no HepMC connection, give up! Throw UserError exception?
    const GenVertexPtr prodVtx = genParticle()->production_vertex();
    if (prodVtx == NULL) return false; // orphaned particle, has to be assume false
    const pair<GenParticlePtr, GenParticlePtr> beams = prodVtx->parent_event()->beam_particles();

    /// @todo Would be nicer to be able to write this recursively up the chain, exiting as soon as a parton or string/cluster is seen
    for (const GenParticlePtr ancestor : Rivet::particles(prodVtx, HepMC::ancestors)) {
      const PdgId pid = ancestor->pdg_id();
      if (ancestor->status() != 2) continue; // no non-standard statuses or beams to be used in decision making
      if (ancestor == beams.first || ancestor == beams.second) continue; // PYTHIA6 uses status 2 for beams, I think... (sigh)
      if (PID::isParton(pid)) continue; // PYTHIA6 also uses status 2 for some partons, I think... (sigh)
      if (PID::isHadron(pid)) return false; // prompt particles can't be from hadron decays
      if (abs(pid) == PID::TAU && abspid() != PID::TAU && !allow_from_prompt_tau) return false; // allow or ban particles from tau decays (permitting tau copies)
      if (abs(pid) == PID::MUON && abspid() != PID::MUON && !allow_from_prompt_mu) return false; // allow or ban particles from muon decays (permitting muon copies)
    }
    return true;
  }
bool isStable ( ) const [inherited]

Whether this particle is stable according to the generator.

Definition at line 28 of file Particle.cc.

                                {
    return genParticle() != NULL &&
      genParticle()->status() == 1 &&
      genParticle()->end_vertex() == NULL;
  }
bool isVisible ( ) const [inherited]

Is this particle potentially visible in a detector?

Definition at line 14 of file Particle.cc.

                                 {
    // Charged particles are visible
    if ( PID::threeCharge(pid()) != 0 ) return true;
    // Neutral hadrons are visible
    if ( PID::isHadron(pid()) ) return true;
    // Photons are visible
    if ( pid() == PID::PHOTON ) return true;
    // Gluons are visible (for parton level analyses)
    if ( pid() == PID::GLUON ) return true;
    // Everything else is invisible
    return false;
  }
double mass ( ) const [inline, inherited]

Get the mass directly.

Definition at line 70 of file ParticleBase.hh.

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

Get the mass**2 directly.

Definition at line 72 of file ParticleBase.hh.

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

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

Definition at line 29 of file ParticleBase.hh.

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

The momentum.

Implements ParticleBase.

Definition at line 89 of file Particle.hh.

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

Cast operator for conversion to FourMomentum.

Definition at line 32 of file ParticleBase.hh.

{ return momentum(); }
operator const GenParticle * ( ) const [inline, inherited]

Cast operator for conversion to GenParticle*.

Definition at line 80 of file Particle.hh.

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

Cast operator to FastJet3 PseudoJet.

Definition at line 71 of file Particle.hh.

{ return pseudojet(); }
const FourVector& origin ( ) const [inline, inherited]

The origin position.

Definition at line 115 of file Particle.hh.

                                     {
      return _origin;
    }
double p ( ) const [inline, inherited]

Get the 3-momentum magnitude directly.

Definition at line 101 of file ParticleBase.hh.

{ return momentum().p(); }
double p2 ( ) const [inline, inherited]

Get the 3-momentum magnitude-squared directly.

Definition at line 103 of file ParticleBase.hh.

{ return momentum().p2(); }
Vector3 p3 ( ) const [inline, inherited]

Get the 3-momentum directly.

Definition at line 99 of file ParticleBase.hh.

{ return momentum().vector3(); }
vector< Particle > parents ( const Cut c = Cuts::OPEN) const [inherited]
Todo:
Add physicalAncestors, allAncestors?

Get a list of the direct parents of the current particle (with optional selection Cut)

Note:
This is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!
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 35 of file Particle.cc.

                                                       {
    vector<Particle> rtn;
    /// @todo Remove this const mess crap when HepMC doesn't suck
    GenVertexPtr gv = const_cast<GenVertexPtr>( genParticle()->production_vertex() );
    if (gv == NULL) return rtn;
    /// @todo Would like to do this, but the range objects are broken
    // foreach (const GenParticlePtr gp, gv->particles(HepMC::children))
    //   rtn += Particle(gp);
    for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::parents); it != gv->particles_end(HepMC::parents); ++it) {
      const Particle p(*it);
      if (c != Cuts::OPEN && !c->accept(p)) continue;
      rtn += p;
    }
    return rtn;
  }
Particles parents ( const FN &  f) const [inline, inherited]

Get a list of the direct parents of the current particle (with selector function)

Note:
This is valid in MC, but may not be answerable experimentally -- use this function with care when replicating experimental analyses!

Definition at line 213 of file Particle.hh.

                                         {
      return filter_select(parents(), f);
    }
PdgId pdgId ( ) const [inline, inherited]

This Particle's PDG ID code (alias).

Deprecated:
Prefer the pid/abspid form

Definition at line 141 of file Particle.hh.

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

Get the $ p_T $ directly (alias).

Definition at line 55 of file ParticleBase.hh.

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

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

Definition at line 62 of file ParticleBase.hh.

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

Get the $ \phi $ directly.

Definition at line 95 of file ParticleBase.hh.

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

This Particle's PDG ID code.

Definition at line 136 of file Particle.hh.

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

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

Definition at line 120 of file ParticleBase.hh.

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

Converter to FastJet3 PseudoJet.

Definition at line 66 of file Particle.hh.

                                               {
      return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
    }
double pseudorapidity ( ) const [inline, inherited]

Get the $ \eta $ directly.

Definition at line 75 of file ParticleBase.hh.

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

Get the $ p_T $ directly.

Definition at line 51 of file ParticleBase.hh.

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

Get the $ p_T $ directly (alias).

Definition at line 53 of file ParticleBase.hh.

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

Get the $ p_T^2 $ directly.

Definition at line 58 of file ParticleBase.hh.

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

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

Definition at line 60 of file ParticleBase.hh.

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

x component of momentum.

Definition at line 106 of file ParticleBase.hh.

{ return momentum().x(); }
double px2 ( ) const [inline, inherited]

x component of momentum, squared.

Definition at line 113 of file ParticleBase.hh.

{ return momentum().x2(); }
double py ( ) const [inline, inherited]

y component of momentum.

Definition at line 108 of file ParticleBase.hh.

{ return momentum().y(); }
double py2 ( ) const [inline, inherited]

y component of momentum, squared.

Definition at line 115 of file ParticleBase.hh.

{ return momentum().y2(); }
double pz ( ) const [inline, inherited]

z component of momentum.

Definition at line 110 of file ParticleBase.hh.

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

z component of momentum, squared.

Definition at line 117 of file ParticleBase.hh.

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

Get the $ y $ directly (alias).

Definition at line 86 of file ParticleBase.hh.

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

Get the $ y $ directly.

Definition at line 84 of file ParticleBase.hh.

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

Set the momentum.

Definition at line 94 of file Particle.hh.

                                                        {
      _momentum = momentum;
      return *this;
    }
Particle& setMomentum ( double  E,
double  px,
double  py,
double  pz 
) [inline, inherited]

Set the momentum via components.

Definition at line 100 of file Particle.hh.

                                                                     {
      _momentum = FourMomentum(E, px, py, pz);
      return *this;
    }
Particle& setOrigin ( const FourVector position) [inline, inherited]

Set the origin position.

Definition at line 119 of file Particle.hh.

                                                    {
      _origin = position;
      return *this;
    }
Particle& setOrigin ( double  t,
double  x,
double  y,
double  z 
) [inline, inherited]

Set the origin position via components.

Definition at line 124 of file Particle.hh.

                                                                {
      _origin = FourMomentum(t, x, y, z);
      return *this;
    }
vector< Particle > stableDescendants ( const Cut c = Cuts::OPEN) const [inherited]

Get a list of all the stable descendants from the current particle (with optional selection Cut)

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:
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 100 of file Particle.cc.

                                                                 {
    vector<Particle> rtn;
    if (isStable()) return rtn;
    /// @todo Remove this const mess crap when HepMC doesn't suck
    GenVertexPtr gv = const_cast<GenVertexPtr>( genParticle()->end_vertex() );
    if (gv == NULL) return rtn;
    /// @todo Would like to do this, but the range objects are broken
    // foreach (const GenParticlePtr gp, gv->particles(HepMC::descendants))
    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) continue;
      const Particle p(*it);
      if (!p.isStable()) continue;
      if (c != Cuts::OPEN && !c->accept(p)) continue;
      rtn += p;
    }
    return rtn;
  }
Particles stableDescendants ( const FN &  f) const [inline, inherited]

Get a list of all the stable descendants from the current particle (with selector function)

Definition at line 358 of file Particle.hh.

                                                   {
      return filter_select(stableDescendants(), f);
    }
double theta ( ) const [inline, inherited]

Synonym for polarAngle.

Definition at line 122 of file ParticleBase.hh.

{ return momentum().theta(); }
int threeCharge ( ) const [inline, inherited]

Alias for charge3

Deprecated:
Use charge3

Definition at line 160 of file Particle.hh.

{ return PID::threeCharge(pid()); }
Particle & transformBy ( const LorentzTransform lt) [inherited]

Apply an active Lorentz transform to this particle.

Definition at line 8 of file Particle.cc.

                                                            {
    _momentum = lt.transform(_momentum);
    return *this;
  }

Member Data Documentation

Definition at line 32 of file DressedLeptons.hh.

Definition at line 31 of file DressedLeptons.hh.

PdgId _id [protected, inherited]

The PDG ID code for this Particle.

Definition at line 382 of file Particle.hh.

FourMomentum _momentum [protected, inherited]

The momentum of this particle.

Definition at line 385 of file Particle.hh.

FourVector _origin [protected, inherited]

The creation position of this particle.

Definition at line 388 of file Particle.hh.

const GenParticle* _original [protected, inherited]

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

Definition at line 379 of file Particle.hh.


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