class Rivet::DressedLepton

Rivet::DressedLepton

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

#include <DressedLeptons.hh>

Inherits from Rivet::Particle, Rivet::ParticleBase

Public Functions

Name
virtual const FourMomentum &momentum() const
The momentum.
Particle &setMomentum(const FourMomentum & momentum)
Set the momentum.
Particle &setMomentum(double E, double px, double py, double pz)
Set the momentum via components.
Particle &transformBy(const LorentzTransform & lt)
Apply an active Lorentz transform to this particle.
const FourVector &origin() const
The origin position (and time).
Particle &setOrigin(const FourVector & position)
Set the origin position.
Particle &setOrigin(double t, double x, double y, double z)
Set the origin position via components.
Vector3closestApproach() const
Find the point of closest approach to the primary vertex.
virtual fastjet::PseudoJetpseudojet() const
Converter to FastJet3 PseudoJet.
operator PseudoJet() const
Cast operator to FastJet3 PseudoJet.
Particle &setGenParticle(ConstGenParticlePtr gp)
Set a const pointer to the original GenParticle.
ConstGenParticlePtrgenParticle() const
Get a const pointer to the original GenParticle.
operator ConstGenParticlePtr() const
PdgIdpid() const
This Particle’s PDG ID code.
PdgIdabspid() const
Absolute value of the PDG ID code.
doublecharge() const
The charge of this Particle.
doubleabscharge() const
The absolute charge of this Particle.
intcharge3() const
Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
intabscharge3() const
Three times the absolute charge of this Particle (i.e. integer multiple of smallest quark charge).
boolisCharged() const
Is this Particle charged?
boolisHadron() const
Is this a hadron?
boolisMeson() const
Is this a meson?
boolisBaryon() const
Is this a baryon?
boolisLepton() const
Is this a lepton?
boolisChargedLepton() const
Is this a charged lepton?
boolisNeutrino() const
Is this a neutrino?
boolhasBottom() const
Does this (hadron) contain a b quark?
boolhasCharm() const
Does this (hadron) contain a c quark?
boolisVisible() const
Is this particle potentially visible in a detector?
boolisParton() const
Is this a parton? (Hopefully not very often… fiducial FTW)
virtual voidsetConstituents(const Particles & cs, bool setmom =false)
Set direct constituents of this particle.
virtual voidaddConstituent(const Particle & c, bool addmom =false)
Add a single direct constituent to this particle.
virtual voidaddConstituents(const Particles & cs, bool addmom =false)
Add direct constituents to this particle.
boolisComposite() const
Determine if this Particle is a composite of other RivetParticles.
const Particles &constituents() const
Direct constituents of this particle, returned by reference.
const Particlesconstituents(const ParticleSorter & sorter) const
Direct constituents of this particle, sorted by a functor.
const Particlesconstituents(const Cut & c) const
Direct constituents of this particle, filtered by a Cut.
const Particlesconstituents(const Cut & c, const ParticleSorter & sorter) const
Direct constituents of this particle, sorted by a functor.
const Particlesconstituents(const ParticleSelector & selector) const
Direct constituents of this particle, filtered by a selection functor.
const Particlesconstituents(const ParticleSelector & selector, const ParticleSorter & sorter) const
Direct constituents of this particle, filtered and sorted by functors.
ParticlesrawConstituents() const
Fundamental constituents of this particle.
const ParticlesrawConstituents(const ParticleSorter & sorter) const
Fundamental constituents of this particle, sorted by a functor.
const ParticlesrawConstituents(const Cut & c) const
Fundamental constituents of this particle, filtered by a Cut.
const ParticlesrawConstituents(const Cut & c, const ParticleSorter & sorter) const
Fundamental constituents of this particle, sorted by a functor.
const ParticlesrawConstituents(const ParticleSelector & selector) const
Fundamental constituents of this particle, filtered by a selection functor.
const ParticlesrawConstituents(const ParticleSelector & selector, const ParticleSorter & sorter) const
Fundamental constituents of this particle, filtered and sorted by functors.
Particlesparents(const Cut & c =Cuts::OPEN) const
Particlesparents(const ParticleSelector & f) const
boolhasParentWith(const ParticleSelector & f) const
boolhasParentWith(const Cut & c) const
boolhasParentWithout(const ParticleSelector & f) const
boolhasParentWithout(const Cut & c) const
boolhasParent(PdgId pid) const
Particlesancestors(const Cut & c =Cuts::OPEN, bool only_physical =true) const
Particlesancestors(const ParticleSelector & f, bool only_physical =true) const
boolhasAncestorWith(const ParticleSelector & f, bool only_physical =true) const
boolhasAncestorWith(const Cut & c, bool only_physical =true) const
boolhasAncestorWithout(const ParticleSelector & f, bool only_physical =true) const
boolhasAncestorWithout(const Cut & c, bool only_physical =true) const
boolhasAncestor(PdgId pid, bool only_physical =true) const
boolfromBottom() const
Determine whether the particle is from a b-hadron decay.
boolfromCharm() const
Determine whether the particle is from a c-hadron decay.
boolfromHadron() const
Determine whether the particle is from a hadron decay.
boolfromTau(bool prompt_taus_only =false) const
Determine whether the particle is from a tau decay.
boolfromPromptTau() const
Determine whether the particle is from a prompt tau decay.
boolfromHadronicTau(bool prompt_taus_only =false) const
Determine whether the particle is from a tau which decayed hadronically.
DEPRECATED(“Too vague: use fromHadron)
boolisDirect(bool allow_from_direct_tau =false, bool allow_from_direct_mu =false) const
Shorthand definition of ‘promptness’ based on set definition flags.
boolisPrompt(bool allow_from_prompt_tau =false, bool allow_from_prompt_mu =false) const
Alias for isDirect.
boolisStable() const
Whether this particle is stable according to the generator.
Particleschildren(const Cut & c =Cuts::OPEN) const
Get a list of the direct descendants from the current particle (with optional selection Cut)
Particleschildren(const ParticleSelector & f) const
Get a list of the direct descendants from the current particle (with selector function)
boolhasChildWith(const ParticleSelector & f) const
boolhasChildWith(const Cut & c) const
boolhasChildWithout(const ParticleSelector & f) const
boolhasChildWithout(const Cut & c) const
ParticlesallDescendants(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)
ParticlesallDescendants(const ParticleSelector & f, bool remove_duplicates =true) const
Get a list of all the descendants from the current particle (with selector function)
boolhasDescendantWith(const ParticleSelector & f, bool remove_duplicates =true) const
boolhasDescendantWith(const Cut & c, bool remove_duplicates =true) const
boolhasDescendantWithout(const ParticleSelector & f, bool remove_duplicates =true) const
boolhasDescendantWithout(const Cut & c, bool remove_duplicates =true) const
ParticlesstableDescendants(const Cut & c =Cuts::OPEN) const
ParticlesstableDescendants(const ParticleSelector & f) const
Get a list of all the stable descendants from the current particle (with selector function)
boolhasStableDescendantWith(const ParticleSelector & f) const
boolhasStableDescendantWith(const Cut & c) const
boolhasStableDescendantWithout(const ParticleSelector & f) const
boolhasStableDescendantWithout(const Cut & c) const
doubleflightLength() const
boolisFirstWith(const ParticleSelector & f) const
Determine whether a particle is the first in a decay chain to meet the function requirement.
boolisFirstWithout(const ParticleSelector & f) const
Determine whether a particle is the first in a decay chain not to meet the function requirement.
boolisLastWith(const ParticleSelector & f) const
Determine whether a particle is the last in a decay chain to meet the function requirement.
boolisLastWithout(const ParticleSelector & f) const
Determine whether a particle is the last in a decay chain not to meet the function requirement.
boolisSame(const Particle & other) const
const FourMomentum &mom() const
Get equivalent single momentum four-vector (const) (alias).
operator const FourMomentum &() const
Cast operator for conversion to FourMomentum.
doubleE() const
Get the energy directly.
doubleenergy() const
Get the energy directly (alias).
doubleE2() const
Get the energy-squared.
doubleenergy2() const
Get the energy-squared (alias).
doublept() const
Get the ( p_T ) directly.
doublepT() const
Get the ( p_T ) directly (alias).
doubleperp() const
Get the ( p_T ) directly (alias).
doublept2() const
Get the ( p_T^2 ) directly.
doublepT2() const
Get the ( p_T^2 ) directly (alias).
doubleperp2() const
Get the ( p_T^2 ) directly (alias).
doubleEt() const
Get the ( E_T ) directly.
doubleEt2() const
Get the ( E_T^2 ) directly.
doublemass() const
Get the mass directly.
doublemass2() const
Get the mass**2 directly.
doublepseudorapidity() const
Get the ( \eta ) directly.
doubleeta() const
Get the ( \eta ) directly (alias).
doubleabspseudorapidity() const
Get the (
doubleabseta() const
Get the (
doublerapidity() const
Get the ( y ) directly.
doublerap() const
Get the ( y ) directly (alias).
doubleabsrapidity() const
Get the (
doubleabsrap() const
Get the (
doubleazimuthalAngle(const PhiMapping mapping =ZERO_2PI) const
Azimuthal angle ( \phi ).
doublephi(const PhiMapping mapping =ZERO_2PI) const
Get the ( \phi ) directly.
Vector3p3() const
Get the 3-momentum directly.
doublep() const
Get the 3-momentum magnitude directly.
doublep2() const
Get the 3-momentum magnitude-squared directly.
Vector3ptvec() const
Get the transverse 3-momentum directly.
Vector3pTvec() const
Get the transverse 3-momentum directly.
doublepx() const
x component of momentum.
doublepy() const
y component of momentum.
doublepz() const
z component of momentum.
doublepx2() const
x component of momentum, squared.
doublepy2() const
y component of momentum, squared.
doublepz2() const
z component of momentum, squared.
doublepolarAngle() const
Angle subtended by the 3-vector and the z-axis.
doubletheta() const
Synonym for polarAngle.
doubleangle(const ParticleBase & v) const
Angle between this vector and another.
doubleangle(const FourVector & v) const
Angle between this vector and another.
doubleangle(const Vector3 & v3) const
Angle between this vector and another (3-vector)
doubledot(const ParticleBase & v) const
Lorentz dot product between this 4-vector and another.
doubledot(const FourVector & v) const
Angle between this 4-vector and another.
DressedLepton(const Particle & dlepton)
Copy constructor (from Particle)
DressedLepton(const Particle & lepton, const Particles & photons, bool momsum =true)
voidaddPhoton(const Particle & p, bool momsum =true)
const Particle &bareLepton() const
Retrieve the bare lepton.
const Particle &constituentLepton() const
const Particlesphotons() const
Retrieve the clustered photons.
const ParticlesconstituentPhotons() const

Additional inherited members

Public Functions inherited from Rivet::Particle

Name
Particle()
Particle(PdgId pid, const FourMomentum & mom, const FourVector & pos =FourVector(), ConstGenParticlePtr gp =nullptr)
Constructor from PID and momentum.
Particle(PdgId pid, const FourMomentum & mom, ConstGenParticlePtr gp, const FourVector & pos =FourVector())
Constructor from PID, momentum, and a GenParticle for relational links.
Particle(ConstGenParticlePtr gp)
Constructor from a HepMC GenParticle pointer.
Particle(const RivetHepMC::GenParticle & gp)
Constructor from a HepMC GenParticle reference.

Public Functions inherited from Rivet::ParticleBase

Name
ParticleBase()
Default constructor.
virtual~ParticleBase()
Virtual destructor.

Detailed Description

class Rivet::DressedLepton;

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

Deprecated:

Just use Particle.constituents() now.

Todo: Remove completely – it’s unnecessary and too confusing (esp. between copying & aggregating)

Public Functions Documentation

function momentum

inline virtual const FourMomentum & momentum() const

The momentum.

Reimplements: Rivet::ParticleBase::momentum

function setMomentum

inline Particle & setMomentum(
    const FourMomentum & momentum
)

Set the momentum.

function setMomentum

inline Particle & setMomentum(
    double E,
    double px,
    double py,
    double pz
)

Set the momentum via components.

function transformBy

Particle & transformBy(
    const LorentzTransform & lt
)

Apply an active Lorentz transform to this particle.

function origin

inline const FourVector & origin() const

The origin position (and time).

function setOrigin

inline Particle & setOrigin(
    const FourVector & position
)

Set the origin position.

function setOrigin

inline Particle & setOrigin(
    double t,
    double x,
    double y,
    double z
)

Set the origin position via components.

function closestApproach

inline Vector3 closestApproach() const

Find the point of closest approach to the primary vertex.

TodoCheck that this works with all angles

function pseudojet

inline virtual fastjet::PseudoJet pseudojet() const

Converter to FastJet3 PseudoJet.

function operator PseudoJet

inline operator PseudoJet() const

Cast operator to FastJet3 PseudoJet.

function setGenParticle

inline Particle & setGenParticle(
    ConstGenParticlePtr gp
)

Set a const pointer to the original GenParticle.

function genParticle

inline ConstGenParticlePtr genParticle() const

Get a const pointer to the original GenParticle.

function operator ConstGenParticlePtr

inline explicit operator ConstGenParticlePtr() const

Note: Not implicit since that would enable accidental Particle::operator== comparisons

Cast operator for conversion to GenParticle*

function pid

inline PdgId pid() const

This Particle’s PDG ID code.

function abspid

inline PdgId abspid() const

Absolute value of the PDG ID code.

function charge

inline double charge() const

The charge of this Particle.

function abscharge

inline double abscharge() const

The absolute charge of this Particle.

function charge3

inline int charge3() const

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

function abscharge3

inline int abscharge3() const

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

function isCharged

inline bool isCharged() const

Is this Particle charged?

function isHadron

inline bool isHadron() const

Is this a hadron?

function isMeson

inline bool isMeson() const

Is this a meson?

function isBaryon

inline bool isBaryon() const

Is this a baryon?

function isLepton

inline bool isLepton() const

Is this a lepton?

function isChargedLepton

inline bool isChargedLepton() const

Is this a charged lepton?

function isNeutrino

inline bool isNeutrino() const

Is this a neutrino?

function hasBottom

inline bool hasBottom() const

Does this (hadron) contain a b quark?

function hasCharm

inline bool hasCharm() const

Does this (hadron) contain a c quark?

function isVisible

bool isVisible() const

Is this particle potentially visible in a detector?

function isParton

inline bool isParton() const

Is this a parton? (Hopefully not very often… fiducial FTW)

function setConstituents

virtual void setConstituents(
    const Particles & cs,
    bool setmom =false
)

Set direct constituents of this particle.

function addConstituent

virtual void addConstituent(
    const Particle & c,
    bool addmom =false
)

Add a single direct constituent to this particle.

function addConstituents

virtual void addConstituents(
    const Particles & cs,
    bool addmom =false
)

Add direct constituents to this particle.

function isComposite

inline bool isComposite() const

Determine if this Particle is a composite of other RivetParticles.

function constituents

inline const Particles & constituents() const

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.

function constituents

inline const Particles constituents(
    const ParticleSorter & sorter
) const

Direct constituents of this particle, sorted by a functor.

Note: Returns a copy, thanks to the sorting

function constituents

inline const Particles constituents(
    const Cut & c
) const

Direct constituents of this particle, filtered by a Cut.

Note: Returns a copy, thanks to the filtering

function constituents

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

Direct constituents of this particle, sorted by a functor.

Note: Returns a copy, thanks to the filtering and sorting

function constituents

inline const Particles constituents(
    const ParticleSelector & selector
) const

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

Note: Returns a copy, thanks to the filtering

function constituents

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

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

Note: Returns a copy, thanks to the filtering and sorting

function rawConstituents

Particles rawConstituents() const

Fundamental constituents of this particle.

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

function rawConstituents

inline const Particles rawConstituents(
    const ParticleSorter & sorter
) const

Fundamental constituents of this particle, sorted by a functor.

Note: Returns a copy, thanks to the sorting

function rawConstituents

inline const Particles rawConstituents(
    const Cut & c
) const

Fundamental constituents of this particle, filtered by a Cut.

Note: Returns a copy, thanks to the filtering

function rawConstituents

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

Fundamental constituents of this particle, sorted by a functor.

Note: Returns a copy, thanks to the filtering and sorting

function rawConstituents

inline const Particles rawConstituents(
    const ParticleSelector & selector
) const

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

Note: Returns a copy, thanks to the filtering

function rawConstituents

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

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

Note: Returns a copy, thanks to the filtering and sorting

function parents

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

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

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

function parents

inline Particles parents(
    const ParticleSelector & f
) const

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

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

function hasParentWith

inline bool hasParentWith(
    const ParticleSelector & f
) const

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

Check whether any particle in the particle’s parent list has the requested property

function hasParentWith

bool hasParentWith(
    const Cut & c
) const

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

Check whether any particle in the particle’s parent list has the requested property

function hasParentWithout

inline bool hasParentWithout(
    const ParticleSelector & f
) const

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

Check whether any particle in the particle’s parent list does not have the requested property

function hasParentWithout

bool hasParentWithout(
    const Cut & c
) const

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

Check whether any particle in the particle’s parent list does not have the requested property

function hasParent

bool hasParent(
    PdgId pid
) const

Deprecated:

Prefer e.g. hasParentWith(Cut::pid == 123)

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

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

function ancestors

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

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!

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

function ancestors

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

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!

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

function hasAncestorWith

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

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

Check whether any particle in the particle’s ancestor list has the requested property

function hasAncestorWith

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

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

Check whether any particle in the particle’s ancestor list has the requested property

function hasAncestorWithout

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

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

Check whether any particle in the particle’s ancestor list does not have the requested property

function hasAncestorWithout

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

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

Check whether any particle in the particle’s ancestor list does not have the requested property

function hasAncestor

bool hasAncestor(
    PdgId pid,
    bool only_physical =true
) const

Deprecated:

Prefer hasAncestorWith(Cuts::pid == pid) etc.

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

Check whether a given PID is found in the particle’s ancestor list

function fromBottom

bool fromBottom() const

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

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

function fromCharm

bool fromCharm() const

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

Note:

  • If a hadron contains b and c quarks it is considered a bottom hadron and NOT a charm hadron.
  • This question is valid in MC, but may not be perfectly answerable experimentally – use this function with care when replicating experimental analyses!

function fromHadron

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

function fromTau

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

function fromPromptTau

inline bool fromPromptTau() const

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!

function fromHadronicTau

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

function DEPRECATED

inline DEPRECATED(
    "Too vague: use  fromHadron) || fromPromptTau(,
    or isDirect()" 
) const

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

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

function isDirect

bool isDirect(
    bool allow_from_direct_tau =false,
    bool allow_from_direct_mu =false
) const

Shorthand definition of ‘promptness’ based on set definition flags.

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

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.

function isPrompt

inline bool isPrompt(
    bool allow_from_prompt_tau =false,
    bool allow_from_prompt_mu =false
) const

Alias for isDirect.

function isStable

bool isStable() const

Whether this particle is stable according to the generator.

function children

Particles 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?

function children

inline Particles children(
    const ParticleSelector & f
) const

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

function hasChildWith

inline bool hasChildWith(
    const ParticleSelector & f
) const

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

Check whether any direct child of this particle has the requested property

function hasChildWith

bool hasChildWith(
    const Cut & c
) const

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

Check whether any direct child of this particle has the requested property

function hasChildWithout

inline bool hasChildWithout(
    const ParticleSelector & f
) const

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

Check whether any direct child of this particle does not have the requested property

function hasChildWithout

bool hasChildWithout(
    const Cut & c
) const

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

Check whether any direct child of this particle does not have the requested property

function allDescendants

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)

function allDescendants

inline Particles allDescendants(
    const ParticleSelector & f,
    bool remove_duplicates =true
) const

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

function hasDescendantWith

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

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

Check whether any descendant of this particle has the requested property

function hasDescendantWith

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

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

Check whether any descendant of this particle has the requested property

function hasDescendantWithout

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

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

Check whether any descendant of this particle does not have the requested property

function hasDescendantWithout

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

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

Check whether any descendant of this particle does not have the requested property

function stableDescendants

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

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

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

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

function stableDescendants

inline Particles stableDescendants(
    const ParticleSelector & f
) const

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

function hasStableDescendantWith

inline bool hasStableDescendantWith(
    const ParticleSelector & f
) const

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

Check whether any stable descendant of this particle has the requested property

function hasStableDescendantWith

bool hasStableDescendantWith(
    const Cut & c
) const

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

Check whether any stable descendant of this particle has the requested property

function hasStableDescendantWithout

inline bool hasStableDescendantWithout(
    const ParticleSelector & f
) const

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

Check whether any stable descendant of this particle does not have the requested property

function hasStableDescendantWithout

bool hasStableDescendantWithout(
    const Cut & c
) const

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

Check whether any stable descendant of this particle does not have the requested property

function flightLength

double flightLength() const

Note: Divide by mm or cm as usual to get the appropriate units.

Flight length of the particle from origin to decay

function isFirstWith

inline bool isFirstWith(
    const ParticleSelector & f
) const

Determine whether a particle is the first in a decay chain to meet the function requirement.

function isFirstWithout

inline bool isFirstWithout(
    const ParticleSelector & f
) const

Determine whether a particle is the first in a decay chain not to meet the function requirement.

function isLastWith

inline bool isLastWith(
    const ParticleSelector & f
) const

Determine whether a particle is the last in a decay chain to meet the function requirement.

function isLastWithout

inline bool isLastWithout(
    const ParticleSelector & f
) const

Determine whether a particle is the last in a decay chain not to meet the function requirement.

function isSame

inline bool isSame(
    const Particle & other
) const

Note: Not a deep comparison: GenParticle ptr and constituents are not used in the comparison

Compare particles, based on “external” characteristics, with a little angular tolerance

function mom

inline const FourMomentum & mom() const

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

function operator const FourMomentum &

inline operator const FourMomentum &() const

Cast operator for conversion to FourMomentum.

function E

inline double E() const

Get the energy directly.

function energy

inline double energy() const

Get the energy directly (alias).

function E2

inline double E2() const

Get the energy-squared.

function energy2

inline double energy2() const

Get the energy-squared (alias).

function pt

inline double pt() const

Get the ( p_T ) directly.

function pT

inline double pT() const

Get the ( p_T ) directly (alias).

function perp

inline double perp() const

Get the ( p_T ) directly (alias).

function pt2

inline double pt2() const

Get the ( p_T^2 ) directly.

function pT2

inline double pT2() const

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

function perp2

inline double perp2() const

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

function Et

inline double Et() const

Get the ( E_T ) directly.

function Et2

inline double Et2() const

Get the ( E_T^2 ) directly.

function mass

inline double mass() const

Get the mass directly.

function mass2

inline double mass2() const

Get the mass**2 directly.

function pseudorapidity

inline double pseudorapidity() const

Get the ( \eta ) directly.

function eta

inline double eta() const

Get the ( \eta ) directly (alias).

function abspseudorapidity

inline double abspseudorapidity() const

Get the ( |\eta| ) directly.

function abseta

inline double abseta() const

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

function rapidity

inline double rapidity() const

Get the ( y ) directly.

function rap

inline double rap() const

Get the ( y ) directly (alias).

function absrapidity

inline double absrapidity() const

Get the ( |y| ) directly.

function absrap

inline double absrap() const

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

function azimuthalAngle

inline double azimuthalAngle(
    const PhiMapping mapping =ZERO_2PI
) const

Azimuthal angle ( \phi ).

function phi

inline double phi(
    const PhiMapping mapping =ZERO_2PI
) const

Get the ( \phi ) directly.

function p3

inline Vector3 p3() const

Get the 3-momentum directly.

function p

inline double p() const

Get the 3-momentum magnitude directly.

function p2

inline double p2() const

Get the 3-momentum magnitude-squared directly.

function ptvec

inline Vector3 ptvec() const

Get the transverse 3-momentum directly.

function pTvec

inline Vector3 pTvec() const

Get the transverse 3-momentum directly.

function px

inline double px() const

x component of momentum.

function py

inline double py() const

y component of momentum.

function pz

inline double pz() const

z component of momentum.

function px2

inline double px2() const

x component of momentum, squared.

function py2

inline double py2() const

y component of momentum, squared.

function pz2

inline double pz2() const

z component of momentum, squared.

function polarAngle

inline double polarAngle() const

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

function theta

inline double theta() const

Synonym for polarAngle.

function angle

inline double angle(
    const ParticleBase & v
) const

Angle between this vector and another.

function angle

inline double angle(
    const FourVector & v
) const

Angle between this vector and another.

function angle

inline double angle(
    const Vector3 & v3
) const

Angle between this vector and another (3-vector)

function dot

inline double dot(
    const ParticleBase & v
) const

Lorentz dot product between this 4-vector and another.

function dot

inline double dot(
    const FourVector & v
) const

Angle between this 4-vector and another.

function DressedLepton

DressedLepton(
    const Particle & dlepton
)

Copy constructor (from Particle)

function DressedLepton

DressedLepton(
    const Particle & lepton,
    const Particles & photons,
    bool momsum =true
)

Note: This is not a copy constructor, hence the explicit second argument even if empty

Components constructor

function addPhoton

void addPhoton(
    const Particle & p,
    bool momsum =true
)

Todo: Deprecate and override add/setConstituents instead?

Add a photon to the dressed lepton

function bareLepton

const Particle & bareLepton() const

Retrieve the bare lepton.

function constituentLepton

inline const Particle & constituentLepton() const

Deprecated:

Prefer the more physicsy bareLepton()

Retrieve the bare lepton (alias)

function photons

inline const Particles photons() const

Retrieve the clustered photons.

function constituentPhotons

inline const Particles constituentPhotons() const

Deprecated:

Prefer the shorter photons()

Retrieve the clustered photons (alias)


Updated on 2022-08-07 at 20:17:16 +0100