Loading [MathJax]/extensions/tex2jax.js
 |
ATLAS Offline Software
|
Go to the documentation of this file.
5 #ifndef XAODPFLOW_VERSIONS_FLOWELEMENT_V1_H
6 #define XAODPFLOW_VERSIONS_FLOWELEMENT_V1_H
9 #include "AthLinks/ElementLink.h"
69 virtual double pt()
const override;
70 virtual double eta()
const override;
71 virtual double phi()
const override;
72 virtual double m()
const override;
73 virtual double e()
const override;
74 virtual double rapidity()
const override;
133 std::vector<const xAOD::IParticle*>
otherObjects()
const ;
void setOtherObjectLinks(const std::vector< ElementLink< IParticleContainer >> &elV, const std::vector< float > &wV)
virtual double rapidity() const override
The true rapidity (y) of the particle.
virtual double m() const override
The invariant mass of the particle.
vertex_t vertexType() const
std::size_t nOtherObjects() const
void setSignalType(signal_t t)
DATAVECTOR_BASE(xAOD::FlowElement_v1, xAOD::IParticle)
std::vector< std::pair< const xAOD::IParticle *, float > > otherObjectsAndWeights() const
const std::vector< ElementLink< IParticleContainer > > & chargedObjectLinks() const
Access to the EL.
void setVertexType(vertex_t t)
const std::vector< ElementLink< IParticleContainer > > & otherObjectLinks() const
SignalType
Enum to encode the nature of the object this FlowElement represents.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
virtual double phi() const override
The azimuthal angle ( ) of the particle.
virtual double pt() const override
void setOtherObjectLinks(const std::vector< ElementLink< IParticleContainer >> &elV)
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
Class providing the definition of the 4-vector interface.
std::size_t nChargedObjects() const
bool isMatchedToPV(MatchedPVType vxtype=HardScatter) const
TLorentzVector FourMom_t
Definition of the 4-momentum type.
std::pair< const xAOD::IParticle *, float > otherObjectAndWeight(std::size_t i) const
std::vector< std::pair< const xAOD::IParticle *, float > > chargedObjectsAndWeights() const
Description of a calorimeter cluster.
MatchedPVType
Enum to encode high-level information on the vertex associated to this FlowElement.
signal_t signalType() const
const xAOD::IParticle * chargedObject(std::size_t i) const
void setChargedObjectLinks(const std::vector< ElementLink< IParticleContainer >> &elV)
void setChargedObjectLinks(const std::vector< ElementLink< IParticleContainer >> &elV, const std::vector< float > &wV)
virtual double e() const override
The total energy of the particle.
std::vector< const xAOD::IParticle * > chargedObjects() const
An STL vector of pointers that by default owns its pointed-to elements.
virtual double eta() const override
The pseudorapidity ( ) of the particle.
std::pair< const xAOD::IParticle *, float > chargedObjectAndWeight(std::size_t i) const
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
const std::vector< float > & chargedObjectWeights() const
const xAOD::IParticle * otherObject(std::size_t i) const
void setP4(float pt, float eta, float phi, float m)
const std::vector< float > & otherObjectWeights() const
ObjectType
Type of objects that have a representation in the xAOD EDM.
virtual Type::ObjectType type() const override
The type of the object as a simple enumeration.
std::vector< const xAOD::IParticle * > otherObjects() const
A detector object made of other lower level object(s)