rivet is hosted by Hepforge, IPPP Durham
Rivet  2.7.0
Rivet::Particle Class Reference

Particle representation, either from a HepMC::GenEvent or reconstructed. More...

#include <Particle.hh>

Inheritance diagram for Rivet::Particle:
Rivet::ParticleBase Rivet::DressedLepton

Public Member Functions

Constructors
 Particle ()
 
 Particle (PdgId pid, const FourMomentum &mom, const FourVector &pos=FourVector())
 Constructor without GenParticle.
 
 Particle (const GenParticle *gp)
 Constructor from a HepMC GenParticle pointer.
 
 Particle (const GenParticle &gp)
 Constructor from a HepMC GenParticle reference.
 
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.
 
Other representations and implicit casts to momentum-like objects
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
 
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
 
Charge
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 isCharged () const
 Is this Particle charged?
 
Particle species
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 isChargedLepton () const
 Is this a charged 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?
 
bool isParton () const
 Is this a parton? (Hopefully not very often... fiducial FTW)
 
Constituents (for composite particles)
virtual void setConstituents (const Particles &cs, bool setmom=false)
 Set direct constituents of this particle.
 
virtual void addConstituent (const Particle &c, bool addmom=false)
 Add a single direct constituent to this particle.
 
virtual void addConstituents (const Particles &cs, bool addmom=false)
 Add direct constituents to this particle.
 
bool isComposite () const
 Determine if this Particle is a composite of other Rivet Particles.
 
const Particles & constituents () const
 Direct constituents of this particle, returned by reference. More...
 
const Particles constituents (const ParticleSorter &sorter) const
 Direct constituents of this particle, sorted by a functor. More...
 
const Particles constituents (const Cut &c) const
 Direct constituents of this particle, filtered by a Cut. More...
 
const Particles constituents (const Cut &c, const ParticleSorter &sorter) const
 Direct constituents of this particle, sorted by a functor. More...
 
const Particles constituents (const ParticleSelector &selector) const
 Direct constituents of this particle, filtered by a selection functor. More...
 
const Particles constituents (const ParticleSelector &selector, const ParticleSorter &sorter) const
 Direct constituents of this particle, filtered and sorted by functors. More...
 
Particles rawConstituents () const
 Fundamental constituents of this particle. More...
 
const Particles rawConstituents (const ParticleSorter &sorter) const
 Fundamental constituents of this particle, sorted by a functor. More...
 
const Particles rawConstituents (const Cut &c) const
 Fundamental constituents of this particle, filtered by a Cut. More...
 
const Particles rawConstituents (const Cut &c, const ParticleSorter &sorter) const
 Fundamental constituents of this particle, sorted by a functor. More...
 
const Particles rawConstituents (const ParticleSelector &selector) const
 Fundamental constituents of this particle, filtered by a selection functor. More...
 
const Particles rawConstituents (const ParticleSelector &selector, const ParticleSorter &sorter) const
 Fundamental constituents of this particle, filtered and sorted by functors. More...
 
Ancestry (for fundamental particles with a HepMC link)
Particles parents (const Cut &c=Cuts::OPEN) const
 
Particles parents (const ParticleSelector &f) const
 
bool hasParentWith (const ParticleSelector &f) const
 
bool hasParentWith (const Cut &c) const
 
bool hasParentWithout (const ParticleSelector &f) const
 
bool hasParentWithout (const Cut &c) const
 
bool hasParent (PdgId pid) const
 
Particles ancestors (const Cut &c=Cuts::OPEN, bool only_physical=true) const
 
Particles ancestors (const ParticleSelector &f, bool only_physical=true) const
 
bool hasAncestorWith (const ParticleSelector &f, bool only_physical=true) const
 
bool hasAncestorWith (const Cut &c, bool only_physical=true) const
 
bool hasAncestorWithout (const ParticleSelector &f, bool only_physical=true) const
 
bool hasAncestorWithout (const Cut &c, bool only_physical=true) const
 
bool hasAncestor (PdgId pid, bool only_physical=true) const
 
bool fromBottom () const
 Determine whether the particle is from a b-hadron decay. More...
 
bool fromCharm () const
 Determine whether the particle is from a c-hadron decay. More...
 
bool fromHadron () const
 Determine whether the particle is from a hadron decay. More...
 
bool fromTau (bool prompt_taus_only=false) const
 Determine whether the particle is from a tau decay. More...
 
bool fromPromptTau () const
 Determine whether the particle is from a prompt tau decay. More...
 
bool fromHadronicTau (bool prompt_taus_only=false) const
 Determine whether the particle is from a tau which decayed hadronically. More...
 
bool fromDecay () const
 Determine whether the particle is from a hadron or tau decay. More...
 
bool isDirect (bool allow_from_direct_tau=false, bool allow_from_direct_mu=false) const
 Shorthand definition of 'promptness' based on set definition flags. More...
 
bool isPrompt (bool allow_from_prompt_tau=false, bool allow_from_prompt_mu=false) const
 Alias for isDirect.
 
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) More...
 
Particles children (const ParticleSelector &f) const
 Get a list of the direct descendants from the current particle (with selector function)
 
bool hasChildWith (const ParticleSelector &f) const
 
bool hasChildWith (const Cut &c) const
 
bool hasChildWithout (const ParticleSelector &f) const
 
bool hasChildWithout (const Cut &c) const
 
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) More...
 
Particles allDescendants (const ParticleSelector &f, bool remove_duplicates=true) const
 Get a list of all the descendants from the current particle (with selector function)
 
bool hasDescendantWith (const ParticleSelector &f, bool remove_duplicates=true) const
 
bool hasDescendantWith (const Cut &c, bool remove_duplicates=true) const
 
bool hasDescendantWithout (const ParticleSelector &f, bool remove_duplicates=true) const
 
bool hasDescendantWithout (const Cut &c, bool remove_duplicates=true) const
 
Particles stableDescendants (const Cut &c=Cuts::OPEN) const
 
Particles stableDescendants (const ParticleSelector &f) const
 Get a list of all the stable descendants from the current particle (with selector function)
 
bool hasStableDescendantWith (const ParticleSelector &f) const
 
bool hasStableDescendantWith (const Cut &c) const
 
bool hasStableDescendantWithout (const ParticleSelector &f) const
 
bool hasStableDescendantWithout (const Cut &c) const
 
double flightLength () const
 Flight length (divide by mm or cm to get the appropriate units)
 
Duplicate testing
bool isFirstWith (const ParticleSelector &f) const
 Determine whether a particle is the first in a decay chain to meet the function requirement.
 
bool isFirstWithout (const ParticleSelector &f) const
 Determine whether a particle is the first in a decay chain not to meet the function requirement.
 
bool isLastWith (const ParticleSelector &f) const
 Determine whether a particle is the last in a decay chain to meet the function requirement.
 
bool isLastWithout (const ParticleSelector &f) const
 Determine whether a particle is the last in a decay chain not to meet the function requirement.
 
- Public Member Functions inherited from Rivet::ParticleBase
 ParticleBase ()
 Default constructor.
 
virtual ~ParticleBase ()
 Virtual destructor.
 
const FourMomentummom () const
 Get equivalent single momentum four-vector (const) (alias).
 
 operator const FourMomentum & () const
 Cast operator for conversion to FourMomentum.
 
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.
 
Vector3 ptvec () const
 Get the transverse 3-momentum directly.
 
Vector3 pTvec () const
 Get the transverse 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 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)
 

Detailed Description

Particle representation, either from a HepMC::GenEvent or reconstructed.

Constructor & Destructor Documentation

◆ Particle()

Rivet::Particle::Particle ( )
inline

Default constructor.

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

Member Function Documentation

◆ allDescendants()

Particles Rivet::Particle::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)

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

References genParticle(), isStable(), and Rivet::ParticleBase::p().

Referenced by allDescendants(), hasChildWithout(), and hasDescendantWith().

◆ ancestors() [1/2]

Particles Rivet::Particle::ancestors ( const Cut &  c = Cuts::OPEN,
bool  only_physical = true 
) const

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

Note
By default only physical ancestors, with status=2, are returned.
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

References genParticle(), and Rivet::ParticleBase::p().

Referenced by ancestors(), hasAncestorWith(), and hasParentWithout().

◆ ancestors() [2/2]

Particles Rivet::Particle::ancestors ( const ParticleSelector &  f,
bool  only_physical = true 
) const
inline

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

Note
By default only physical ancestors, with status=2, are returned.
This is valid in MC, but may not be answerable experimentally – use this function with care when replicating experimental analyses!

References ancestors(), and Rivet::filter_select().

◆ children()

Particles Rivet::Particle::children ( const Cut &  c = Cuts::OPEN) const

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

Todo:
isDecayed? How to restrict to physical particles?
Todo:
Something going wrong with taus -> GenParticle nullptr?
Todo:
Remove this const mess crap when HepMC doesn't suck
Todo:
Would like to do this, but the range objects are broken

References genParticle(), isStable(), and Rivet::ParticleBase::p().

Referenced by Rivet::FinalPartons::accept(), children(), hasChildWith(), isLastWith(), isPrompt(), Rivet::PartonicTops::project(), Rivet::TAU_EFF_ATLAS_RUN1(), and Rivet::TAU_EFF_ATLAS_RUN2().

◆ constituents() [1/6]

const Particles& Rivet::Particle::constituents ( ) const
inline

Direct constituents of this particle, returned by reference.

The returned vector will be empty if this particle is non-composite, and its entries may themselves be composites.

Referenced by Rivet::DressedLepton::bareLepton(), Rivet::ZFinder::constituentLeptons(), constituents(), isComposite(), Rivet::DressedLepton::photons(), Rivet::ZFinder::project(), Rivet::WFinder::project(), and rawConstituents().

◆ constituents() [2/6]

const Particles Rivet::Particle::constituents ( const ParticleSorter &  sorter) const
inline

Direct constituents of this particle, sorted by a functor.

Note
Returns a copy, thanks to the sorting

References constituents(), and Rivet::sortBy().

◆ constituents() [3/6]

const Particles Rivet::Particle::constituents ( const Cut &  c) const
inline

Direct constituents of this particle, filtered by a Cut.

Note
Returns a copy, thanks to the filtering

References constituents(), and Rivet::filter_select().

◆ constituents() [4/6]

const Particles Rivet::Particle::constituents ( const Cut &  c,
const ParticleSorter &  sorter 
) const
inline

Direct constituents of this particle, sorted by a functor.

Note
Returns a copy, thanks to the filtering and sorting

References constituents(), and Rivet::sortBy().

◆ constituents() [5/6]

const Particles Rivet::Particle::constituents ( const ParticleSelector &  selector) const
inline

Direct constituents of this particle, filtered by a selection functor.

Note
Returns a copy, thanks to the filtering

References constituents(), and Rivet::filter_select().

◆ constituents() [6/6]

const Particles Rivet::Particle::constituents ( const ParticleSelector &  selector,
const ParticleSorter &  sorter 
) const
inline

Direct constituents of this particle, filtered and sorted by functors.

Note
Returns a copy, thanks to the filtering and sorting

References constituents(), rawConstituents(), and Rivet::sortBy().

◆ fromBottom()

bool Rivet::Particle::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!

References genParticle(), hasAncestorWith(), hasBottom(), isHadron(), and Rivet::ParticleBase::p().

Referenced by Rivet::fromBottom(), and hasAncestorWithout().

◆ fromCharm()

bool Rivet::Particle::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!

References genParticle(), hasAncestorWith(), hasCharm(), isHadron(), and Rivet::ParticleBase::p().

Referenced by Rivet::fromCharm(), and hasAncestorWithout().

◆ fromDecay()

bool Rivet::Particle::fromDecay ( ) const
inline

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!
Deprecated:
Too vague: use fromHadron or fromHadronicTau

References fromHadron(), fromPromptTau(), and isDirect().

Referenced by Rivet::FinalPartons::accept(), Rivet::fromDecay(), Rivet::MC_ParticleAnalysis::init(), and Rivet::FastJets::project().

◆ fromHadron()

bool Rivet::Particle::fromHadron ( ) const

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!

References genParticle(), hasAncestorWith(), isHadron(), and Rivet::ParticleBase::p().

Referenced by fromDecay(), Rivet::fromHadron(), fromTau(), and hasAncestorWithout().

◆ fromHadronicTau()

bool Rivet::Particle::fromHadronicTau ( bool  prompt_taus_only = false) const

Determine whether the particle is from a tau which decayed hadronically.

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

References genParticle(), hasAncestorWith(), Rivet::hasHadronicDecay(), isPrompt(), and Rivet::ParticleBase::p().

Referenced by fromPromptTau().

◆ fromPromptTau()

bool Rivet::Particle::fromPromptTau ( ) const
inline

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!

References fromHadronicTau(), and fromTau().

Referenced by fromDecay(), and Rivet::fromPromptTau().

◆ fromTau()

bool Rivet::Particle::fromTau ( bool  prompt_taus_only = false) 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!

References fromHadron(), genParticle(), hasAncestorWith(), and Rivet::ParticleBase::p().

Referenced by fromPromptTau(), Rivet::fromTau(), and hasAncestorWithout().

◆ hasAncestor()

bool Rivet::Particle::hasAncestor ( PdgId  pid,
bool  only_physical = true 
) const

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!
Deprecated:
Prefer hasAncestorWith(Cuts::pid == pid) etc.

References hasAncestorWith().

Referenced by Rivet::hasAncestor(), hasAncestorWithout(), and Rivet::PartonicTops::project().

◆ hasAncestorWith() [1/2]

bool Rivet::Particle::hasAncestorWith ( const ParticleSelector &  f,
bool  only_physical = true 
) const
inline

Check whether any 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!

References ancestors().

Referenced by fromBottom(), fromCharm(), fromHadron(), fromHadronicTau(), fromTau(), hasAncestor(), Rivet::hasAncestorWith(), hasAncestorWith(), and hasAncestorWithout().

◆ hasAncestorWith() [2/2]

bool Rivet::Particle::hasAncestorWith ( const Cut &  c,
bool  only_physical = true 
) const

Check whether any 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!

References hasAncestorWith(), and Rivet::ParticleBase::p().

◆ hasAncestorWithout() [1/2]

bool Rivet::Particle::hasAncestorWithout ( const ParticleSelector &  f,
bool  only_physical = true 
) const
inline

Check whether any particle in the particle's ancestor list does not have 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!

References fromBottom(), fromCharm(), fromHadron(), fromTau(), hasAncestor(), hasAncestorWith(), and Rivet::ParticleBase::p().

Referenced by Rivet::hasAncestorWithout().

◆ hasAncestorWithout() [2/2]

bool Rivet::Particle::hasAncestorWithout ( const Cut &  c,
bool  only_physical = true 
) const

Check whether any particle in the particle's ancestor list does not have 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!

◆ hasChildWith() [1/2]

bool Rivet::Particle::hasChildWith ( const ParticleSelector &  f) const
inline

Check whether any direct child of this particle 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!

References children().

Referenced by Rivet::hasChildWith(), hasChildWith(), hasChildWithout(), Rivet::hasHadronicDecay(), and Rivet::hasLeptonicDecay().

◆ hasChildWith() [2/2]

bool Rivet::Particle::hasChildWith ( const Cut &  c) const

Check whether any direct child of this particle 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!

References hasChildWith(), and Rivet::ParticleBase::p().

◆ hasChildWithout() [1/2]

bool Rivet::Particle::hasChildWithout ( const ParticleSelector &  f) const
inline

Check whether any direct child of this particle does not have 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!

References allDescendants(), hasChildWith(), and Rivet::ParticleBase::p().

Referenced by Rivet::hasChildWithout().

◆ hasChildWithout() [2/2]

bool Rivet::Particle::hasChildWithout ( const Cut &  c) const

Check whether any direct child of this particle does not have 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!

◆ hasDescendantWith() [1/2]

bool Rivet::Particle::hasDescendantWith ( const ParticleSelector &  f,
bool  remove_duplicates = true 
) const
inline

Check whether any descendant of this particle 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!

References allDescendants().

Referenced by Rivet::hasDescendantWith(), hasDescendantWith(), and hasDescendantWithout().

◆ hasDescendantWith() [2/2]

bool Rivet::Particle::hasDescendantWith ( const Cut &  c,
bool  remove_duplicates = true 
) const

Check whether any descendant of this particle 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!

References hasDescendantWith(), and Rivet::ParticleBase::p().

◆ hasDescendantWithout() [1/2]

bool Rivet::Particle::hasDescendantWithout ( const ParticleSelector &  f,
bool  remove_duplicates = true 
) const
inline

Check whether any descendant of this particle does not have 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!

References hasDescendantWith(), Rivet::ParticleBase::p(), and stableDescendants().

Referenced by Rivet::hasDescendantWithout().

◆ hasDescendantWithout() [2/2]

bool Rivet::Particle::hasDescendantWithout ( const Cut &  c,
bool  remove_duplicates = true 
) const

Check whether any descendant of this particle does not have 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!

◆ hasParent()

bool Rivet::Particle::hasParent ( PdgId  pid) const

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. hasParentWith(Cut::pid == 123)

References hasParentWith().

Referenced by hasParentWithout().

◆ hasParentWith() [1/2]

bool Rivet::Particle::hasParentWith ( const ParticleSelector &  f) const
inline

Check whether any 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!

References parents().

Referenced by hasParent(), hasParentWith(), Rivet::hasParentWith(), and hasParentWithout().

◆ hasParentWith() [2/2]

bool Rivet::Particle::hasParentWith ( const Cut &  c) const

Check whether any 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!

References hasParentWith(), and Rivet::ParticleBase::p().

◆ hasParentWithout() [1/2]

bool Rivet::Particle::hasParentWithout ( const ParticleSelector &  f) const
inline

Check whether any particle in the particle's parent list does not have 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!

References ancestors(), hasParent(), hasParentWith(), and Rivet::ParticleBase::p().

Referenced by Rivet::hasParentWithout().

◆ hasParentWithout() [2/2]

bool Rivet::Particle::hasParentWithout ( const Cut &  c) const

Check whether any particle in the particle's parent list does not have 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!

◆ hasStableDescendantWith() [1/2]

bool Rivet::Particle::hasStableDescendantWith ( const ParticleSelector &  f) const
inline

Check whether any stable descendant of this particle 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!

References stableDescendants().

Referenced by Rivet::hasStableDescendantWith(), hasStableDescendantWith(), and hasStableDescendantWithout().

◆ hasStableDescendantWith() [2/2]

bool Rivet::Particle::hasStableDescendantWith ( const Cut &  c) const

Check whether any stable descendant of this particle 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!

References hasStableDescendantWith(), and Rivet::ParticleBase::p().

◆ hasStableDescendantWithout() [1/2]

bool Rivet::Particle::hasStableDescendantWithout ( const ParticleSelector &  f) const
inline

Check whether any stable descendant of this particle does not have 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!

References flightLength(), hasStableDescendantWith(), and Rivet::ParticleBase::p().

Referenced by Rivet::hasStableDescendantWithout().

◆ hasStableDescendantWithout() [2/2]

bool Rivet::Particle::hasStableDescendantWithout ( const Cut &  c) const

Check whether any stable descendant of this particle does not have 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!

◆ isDirect()

bool Rivet::Particle::isDirect ( bool  allow_from_direct_tau = false,
bool  allow_from_direct_mu = false 
) const

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

A "direct" particle is one directly connected to the hard process. It is a preferred alias for "prompt", since it has no confusing implications about distinguishability by timing information.

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

Note
This one doesn't make any judgements about final-stateness

<

Todo:
Replace awkward caching with C++17 std::optional
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

References abspid(), Rivet::beams(), genParticle(), isHadron(), isParton(), Rivet::particles(), and pid().

Referenced by fromDecay(), Rivet::isDirect(), and isPrompt().

◆ operator const GenParticle *()

Rivet::Particle::operator const GenParticle * ( ) const
inline

Cast operator for conversion to GenParticle*

Todo:
This one's a bad idea since it enables accidental Particle comparisons

References genParticle().

◆ parents() [1/2]

Particles Rivet::Particle::parents ( const Cut &  c = Cuts::OPEN) const

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

References genParticle(), and Rivet::ParticleBase::p().

Referenced by hasParentWith(), isFirstWith(), parents(), and rawConstituents().

◆ parents() [2/2]

Particles Rivet::Particle::parents ( const ParticleSelector &  f) const
inline

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!

References Rivet::filter_select(), and parents().

◆ pdgId()

PdgId Rivet::Particle::pdgId ( ) const
inline

This Particle's PDG ID code (alias).

Deprecated:
Prefer the pid/abspid form

Referenced by Rivet::ChargedFinalState::compare().

◆ rawConstituents() [1/6]

Particles Rivet::Particle::rawConstituents ( ) const

Fundamental constituents of this particle.

Note
Returns {{*this}} if this particle is non-composite.

References constituents(), isComposite(), and Rivet::ParticleBase::p().

Referenced by constituents(), and rawConstituents().

◆ rawConstituents() [2/6]

const Particles Rivet::Particle::rawConstituents ( const ParticleSorter &  sorter) const
inline

Fundamental constituents of this particle, sorted by a functor.

Note
Returns a copy, thanks to the sorting

References rawConstituents(), and Rivet::sortBy().

◆ rawConstituents() [3/6]

const Particles Rivet::Particle::rawConstituents ( const Cut &  c) const
inline

Fundamental constituents of this particle, filtered by a Cut.

Note
Returns a copy, thanks to the filtering

References Rivet::filter_select(), and rawConstituents().

◆ rawConstituents() [4/6]

const Particles Rivet::Particle::rawConstituents ( const Cut &  c,
const ParticleSorter &  sorter 
) const
inline

Fundamental constituents of this particle, sorted by a functor.

Note
Returns a copy, thanks to the filtering and sorting

References rawConstituents(), and Rivet::sortBy().

◆ rawConstituents() [5/6]

const Particles Rivet::Particle::rawConstituents ( const ParticleSelector &  selector) const
inline

Fundamental constituents of this particle, filtered by a selection functor.

Note
Returns a copy, thanks to the filtering

References Rivet::filter_select(), and rawConstituents().

◆ rawConstituents() [6/6]

const Particles Rivet::Particle::rawConstituents ( const ParticleSelector &  selector,
const ParticleSorter &  sorter 
) const
inline

Fundamental constituents of this particle, filtered and sorted by functors.

Note
Returns a copy, thanks to the filtering and sorting

References parents(), rawConstituents(), and Rivet::sortBy().

◆ stableDescendants()

Particles Rivet::Particle::stableDescendants ( const Cut &  c = Cuts::OPEN) const

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

References genParticle(), isStable(), and Rivet::ParticleBase::p().

Referenced by hasDescendantWithout(), hasStableDescendantWith(), and stableDescendants().

◆ threeCharge()

int Rivet::Particle::threeCharge ( ) const
inline

Alias for charge3

Deprecated:
Use charge3

References pid().


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