xAOD::PFO_v1 Class Reference

#include <PFO_v1.h>

Inheritance diagram for xAOD::PFO_v1:
xAOD::IParticle SG::AuxElement SG::IAuxElement

List of all members.

Public Member Functions

 PFO_v1 ()
 Default constructor.
 PFO_v1 (const PFO_v1 &other)
void setP4 (const FourMom_t &vec)
 set the 4-vec
void setP4 (float pt, float eta, float phi, float m=0.0)
 set the 4-vec
const FourMom_tp4EM () const
void setP4EM (const FourMom_t &p4EM)
void setP4EM (float pt, float eta, float phi, float m)
virtual double ptEM () const
virtual double etaEM () const
virtual double phiEM () const
virtual double mEM () const
virtual double eEM () const
float bdtPi0Score () const
void setBDTPi0Score (float BDTPi0Score)
float centerMag () const
void setCenterMag (float CenterMag)
float charge () const
void setCharge (float charge)
template<class T >
void setAttribute (PFODetails::PFOAttributes AttributeType, const T &anAttribute)
template<class T >
bool attribute (PFODetails::PFOAttributes AttributeType, T &anAttribute) const
template<class T >
void setAttribute (const std::string &AttributeType, const T &anAttribute)
template<class T >
bool attribute (const std::string &AttributeType, T &anAttribute) const
bool getClusterMoment (float &theMoment, xAOD::CaloCluster::MomentType momentType) const
const CaloClustercluster (unsigned int index) const
const TrackParticletrack (unsigned int index) const
const xAOD::Vertexvertex () const
bool setVertexLink (const ElementLink< xAOD::VertexContainer > &theVertexLink)
bool setTrackLink (const ElementLink< xAOD::TrackParticleContainer > &theTrack)
bool setClusterLink (const ElementLink< xAOD::CaloClusterContainer > &theCluster)
bool setAssociatedParticleLink (PFODetails::PFOParticleType ParticleType, const ElementLink< IParticleContainer > &theParticle)
void setAssociatedParticleLink (const std::string &ParticleType, const ElementLink< IParticleContainer > &theParticle)
bool setAssociatedParticleLinks (PFODetails::PFOParticleType ParticleType, const std::vector< ElementLink< IParticleContainer > > &theParticles)
bool associatedParticles (PFODetails::PFOParticleType ParticleType, std::vector< const IParticle * > &theParticles) const
void setAssociatedParticleLinks (const std::string &ParticleType, const std::vector< ElementLink< IParticleContainer > > &theParticles)
bool associatedParticles (const std::string &ParticleType, std::vector< const IParticle * > &theParticles) const
TLorentzVector GetVertexCorrectedFourVec (const xAOD::Vertex &vertexToCorrectTo) const
TLorentzVector GetVertexCorrectedFourVec (const TVector3 &vertexToCorrectTo) const
TLorentzVector GetVertexCorrectedEMFourVec (const xAOD::Vertex &vertexToCorrectTo) const
TLorentzVector GetVertexCorrectedEMFourVec (const TVector3 &vertexToCorrectTo) const
void toPersistent ()
template<>
void setAttribute (const std::string &AttributeType, const xAOD::PFODetails::PFOLeptonType &anAttribute)
template<>
bool attribute (const std::string &AttributeType, xAOD::PFODetails::PFOLeptonType &anAttribute) const
template<>
void setAttribute (PFODetails::PFOAttributes AttributeType, const float &anAttribute)
template<>
bool attribute (PFODetails::PFOAttributes AttributeType, float &anAttribute) const
template<>
void setAttribute (PFODetails::PFOAttributes AttributeType, const double &anAttribute)
template<>
bool attribute (PFODetails::PFOAttributes AttributeType, double &anAttribute) const
template<>
void setAttribute (const std::string &AttributeType, const double &anAttribute)
template<>
bool attribute (const std::string &AttributeType, double &anAttribute) const
xAOD::IParticle functions



virtual double pt () const
 The transverse momentum ($p_T$) of the particle.
virtual double eta () const
 The pseudorapidity ($\eta$) of the particle.
virtual double phi () const
 The azimuthal angle ($\phi$) of the particle.
virtual double m () const
 The invariant mass of the particle.
virtual double e () const
 The total energy of the particle.
virtual double rapidity () const
 The true rapidity (y) of the particle.
virtual const FourMom_tp4 () const
 The full 4-momentum of the particle.
virtual Type::ObjectType type () const
 The type of the object as a simple enumeration.

Detailed Description

Class describing a particle flow object


Constructor & Destructor Documentation

xAOD::PFO_v1::PFO_v1 ( const PFO_v1 other  ) 

Copy Constructor


Member Function Documentation

bool xAOD::PFO_v1::associatedParticles ( const std::string &  ParticleType,
std::vector< const IParticle * > &  theParticles 
) const

get a vector of PFO constituent particle types via string

bool xAOD::PFO_v1::associatedParticles ( PFODetails::PFOParticleType  ParticleType,
std::vector< const IParticle * > &  theParticles 
) const

get a vector of PFO constituent particle types via enum

template<class T >
bool xAOD::PFO_v1::attribute ( const std::string &  AttributeType,
T anAttribute 
) const [inline]

Get a PFO Variable via string

template<class T >
bool xAOD::PFO_v1::attribute ( PFODetails::PFOAttributes  AttributeType,
T anAttribute 
) const [inline]

get a PFO Variable via enum

float xAOD::PFO_v1::bdtPi0Score (  )  const

get BDT Score used to classify clusters as Pi0 like or not

float xAOD::PFO_v1::centerMag (  )  const

get CenterMag moment needed for vertex correction

float xAOD::PFO_v1::charge (  )  const

get charge of PFO

const CaloCluster * xAOD::PFO_v1::cluster ( unsigned int  index  )  const

Retrieve a const pointer to a CaloCluster

double xAOD::PFO_v1::eEM (  )  const [virtual]

get EM scale energy

double xAOD::PFO_v1::etaEM (  )  const [virtual]

get EM scale eta

bool xAOD::PFO_v1::getClusterMoment ( float &  theMoment,
xAOD::CaloCluster::MomentType  momentType 
) const

Accessor for cluster moments

TLorentzVector xAOD::PFO_v1::GetVertexCorrectedEMFourVec ( const TVector3 &  vertexToCorrectTo  )  const

Correct EM scale 4-vector to point at a vertex

TLorentzVector xAOD::PFO_v1::GetVertexCorrectedEMFourVec ( const xAOD::Vertex vertexToCorrectTo  )  const

Correct EM scale 4-vector to point at a vertex

TLorentzVector xAOD::PFO_v1::GetVertexCorrectedFourVec ( const TVector3 &  vertexToCorrectTo  )  const

Correct 4-vector to point at a vertex

TLorentzVector xAOD::PFO_v1::GetVertexCorrectedFourVec ( const xAOD::Vertex vertexToCorrectTo  )  const

Correct 4-vector to point at a vertex

double xAOD::PFO_v1::mEM (  )  const [virtual]

get EM scale mass

const PFO_v1::FourMom_t & xAOD::PFO_v1::p4EM (  )  const

get EM scale 4-vector

double xAOD::PFO_v1::phiEM (  )  const [virtual]

get EM scale phi

double xAOD::PFO_v1::ptEM (  )  const [virtual]

get EM scale pt

void xAOD::PFO_v1::setAssociatedParticleLink ( const std::string &  ParticleType,
const ElementLink< IParticleContainer > &  theParticle 
)

Set an IParticle constituent via string

bool xAOD::PFO_v1::setAssociatedParticleLink ( PFODetails::PFOParticleType  ParticleType,
const ElementLink< IParticleContainer > &  theParticle 
)

Set an IParticle constituent via enum

void xAOD::PFO_v1::setAssociatedParticleLinks ( const std::string &  ParticleType,
const std::vector< ElementLink< IParticleContainer > > &  theParticles 
)

Set a vector of PFO constituent particle types via string - overwrite is allowed

bool xAOD::PFO_v1::setAssociatedParticleLinks ( PFODetails::PFOParticleType  ParticleType,
const std::vector< ElementLink< IParticleContainer > > &  theParticles 
)

Set a vector of PFO constituent particle types via enum - overwrite is allowed

template<>
void xAOD::PFO_v1::setAttribute ( PFODetails::PFOAttributes  AttributeType,
const double &  anAttribute 
) [inline]

special implementations for doubles to prevent user from putting doubles in the aux store - convert to float in this case

template<>
void xAOD::PFO_v1::setAttribute ( const std::string &  AttributeType,
const xAOD::PFODetails::PFOLeptonType &  anAttribute 
) [inline]

specaial implementations for floats, for eflowRec JetETMiss variables, to reduce disk space usage

template<class T >
void xAOD::PFO_v1::setAttribute ( const std::string &  AttributeType,
const T anAttribute 
) [inline]

Set a PFO Variable via string - overwrite is allowed

template<class T >
void xAOD::PFO_v1::setAttribute ( PFODetails::PFOAttributes  AttributeType,
const T anAttribute 
) [inline]

Set a PFO Variable via enum - overwrite is allowed

void xAOD::PFO_v1::setBDTPi0Score ( float  BDTPi0Score  ) 

set BDT Score used to classify clusters as Pi0 like or not

void xAOD::PFO_v1::setCenterMag ( float  CenterMag  ) 

set CenterMag moment needed for vertex correction

void xAOD::PFO_v1::setCharge ( float  charge  ) 

set charge of PFO

bool xAOD::PFO_v1::setClusterLink ( const ElementLink< xAOD::CaloClusterContainer > &  theCluster  ) 

Set a cluster constituent

void xAOD::PFO_v1::setP4EM ( float  pt,
float  eta,
float  phi,
float  m 
)

set EM scale 4-vector

void xAOD::PFO_v1::setP4EM ( const FourMom_t p4EM  ) 

set EM scale 4-vector

bool xAOD::PFO_v1::setTrackLink ( const ElementLink< xAOD::TrackParticleContainer > &  theTrack  ) 

Set a track constituent

bool xAOD::PFO_v1::setVertexLink ( const ElementLink< xAOD::VertexContainer > &  theVertexLink  ) 

Set a vertex link

void xAOD::PFO_v1::toPersistent (  ) 

prepare all links for persistification

const TrackParticle * xAOD::PFO_v1::track ( unsigned int  index  )  const

Retrieve a const pointer to a Rec::TrackParticle

const xAOD::Vertex * xAOD::PFO_v1::vertex (  )  const

Retrieve a const pointer to the xAOD::Vertex a charged PFO is associated to


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Generated on 1 Dec 2017 for RootCore Packages by  doxygen 1.6.1