|  | ATLAS Offline Software
    | 
 
 
 
Go to the documentation of this file.
   12 #ifndef PARTICLEEVENT_PARTICLEIMPL_H 
   13 #define PARTICLEEVENT_PARTICLEIMPL_H 1 
   24 #include "AthLinks/ElementLink.h" 
   35 template< 
class INavigable_t,
 
  107           const std::any& 
weight ) 
const;
 
  212   CLHEP::HepLorentzVector 
hlv() 
const;
 
  221   std::ostream& 
dump( std::ostream& 
out ) 
const;
 
  321 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  336 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  348   m_part      ( rhs.m_part )
 
  351 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  358   if ( 
this != &rhs ) {
 
  371 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  375          IParticle_t>::~ParticleImpl()
 
  378 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  383          IParticle_t>::navigableBase()
 const 
  388 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  393          IParticle_t>::momentumBase()
 const 
  398 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  403          IParticle_t>::particleBase()
 const 
  408 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  418 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  424                       const std::any& 
weight )
 const 
  429 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  439 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  449 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  459 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  469 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  479 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  489 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  499 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  504          IParticle_t>::rapidity()
 const 
  509 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  519 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  529 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  539 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  549 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  559 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  564          IParticle_t>::iPt()
 const 
  569 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  574          IParticle_t>::cosPhi()
 const 
  579 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  584          IParticle_t>::sinPhi()
 const 
  589 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  594          IParticle_t>::tanTh()
 const 
  596   return m_mom.
tanTh();
 
  599 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  604          IParticle_t>::cosTh()
 const 
  606   return m_mom.
cosTh();
 
  609 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  614          IParticle_t>::sinTh()
 const 
  616   return m_mom.
sinTh();
 
  619 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  624          IParticle_t>::cotTh()
 const 
  626   return m_mom.
cotTh();
 
  630 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  632 CLHEP::HepLorentzVector
 
  635          IParticle_t>::hlv()
 const 
  640 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  647   return m_mom.errors();
 
  650 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  655          IParticle_t>::kind()
 const 
  660 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  673 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  680   return m_part.dataType();
 
  683 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  688          IParticle_t>::origin()
 const 
  693 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  698              IParticle_t>::originLink()
 const 
  703 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  708          IParticle_t>::hasCharge()
  const 
  710   return m_part.hasCharge();
 
  713 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  720   return m_part.charge();
 
  723 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  728          IParticle_t>::hasPdgId()
 const 
  730   return m_part.hasPdgId();
 
  733 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  740   return m_part.pdgId();
 
  748 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t> 
 
  753          IParticle_t>::getAthenaBarCode()
 const  
  755   return particleBase().getAthenaBarCodeImpl().getAthenaBarCode();
 
  758 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t> 
 
  765   particleBase().getAthenaBarCodeImpl().setAthenaBarCode(
id);
 
  769 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  776   return particleBase().getAthenaBarCodeImpl().hasSameAthenaBarCode(bc);
 
  779 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  786   return particleBase().getAthenaBarCodeImpl().hasSameAthenaBarCodeExceptVersion(bc);
 
  789 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  794          IParticle_t>::getVersion()
 const  
  796   return particleBase().getAthenaBarCodeImpl().getVersion();
 
  799 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t> 
 
  804          IParticle_t>::newVersion() 
 
  806   particleBase().getAthenaBarCodeImpl().newVersion();
 
  809 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  816   particleBase().getAthenaBarCodeImpl().setVersion(newversion);
 
  837 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  842          IParticle_t>::navigableBase()
 
  847 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  852          IParticle_t>::momentumBase()
 
  857 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  862          IParticle_t>::particleBase()
 
  870 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  877   return m_mom.set4Mom( p4 );
 
  880 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  887   return m_mom.set4Mom( p4 );
 
  891 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  896          IParticle_t>::set4Mom( 
const CLHEP::HepLorentzVector& hlv )
 
  898   return m_mom.set4Mom( hlv );
 
  901 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  908   m_part.set_dataType(
x);
 
  911 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  918   m_part.set_charge(
x);
 
  921 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  926          IParticle_t>::set_pdgId( 
int x )
 
  931 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  939   m_part.set_origin(theContainer, 
index);
 
  942 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  950   m_part.set_origin(theContainer, 
vertex);
 
  953 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  960   m_part.set_origin(origin);
 
  963 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  970 template< 
class INavigable_t, 
class I4Momentum_t, 
class IParticle_t>
 
  978   return p4.dump( 
out );
 
  981 #endif //> PARTICLEEVENT_PARTICLEIMPL_H 
  
virtual double sinTh() const
sinus theta
INavigable_t navigable_type
publish the type of the base class (ie: 'traits-itize' it)
virtual double pt() const
transverse momentum
Scalar phi() const
phi method
virtual int pdgId() const
Return enum indicating particle id the enum file is available in Event/EventKernel/PdtPdg....
virtual const I4MomentumError * errors() const
Access to errors, if available; returns 0 if no errors.
AthenaBarCodeVersion_t getVersion() const
Scalar eta() const
pseudorapidity method
virtual double tanTh() const
tan theta
virtual void fillToken(INavigationToken &thisToken, const std::any &weight) const
virtual const Trk::RecVertex * origin() const
Return a RecVertex corresponding to particle Origin
particle_type & particleBase()
access to underlying base type (IParticle-like)
Dummy type needed fro specialized implementation.
void setAthenaBarCode(AthenaBarCode_t)
virtual ~ParticleImpl()
Default constructor:
virtual double cosPhi() const
cosinus phi
virtual bool hasCharge() const
method to check if charge information is available
AccessorTemplate & operator=(AccessorTemplate &&that)
virtual void set_origin(const VxContainer *theContainer, const Trk::VxCandidate *vertex)
momentum_type & momentumBase()
access to underlying base type (I4Momentum-like)
virtual double cotTh() const
cottan theta
virtual double rapidity() const
rapidity
I4Momentum_t momentum_type
publish the type of the base class (ie: 'traits-itize' it)
ParticleImpl & operator=(const ParticleImpl &rhs)
Assignment operator.
Trk::RecVertex inherits from Trk::Vertex.
void setVersion(AthenaBarCodeVersion_t newversion)
virtual double e() const
energy
virtual double px() const
x component of momentum
virtual double m() const
mass
INavigable_t m_nav
The navigable part.
void set_origin(const ElementLink< VxContainer > &origin)
virtual double p() const
mass momentum magnitude
uint64_t AthenaBarCode_t
barcode for all INav4Mom classes
virtual const ElementLink< VxContainer > & originLink() const
Return an ElementLink corresponding to particle's Origin.
virtual ChargeType charge() const
returns charge as a typedef ChargeType currently Charge Type is a double for jets this may be changed...
virtual void set_pdgId(int x)
IParticle_t m_part
The particle-id part.
virtual std::ostream & dump(std::ostream &out) const
Print I4Momentum content.
virtual void set4Mom(const I4Momentum *const p4)
set 4Momentum (will throw exception if cannot be implemented)
AthenaBarCode_t getAthenaBarCode() const
ParticleImpl< INavigable_t, I4Momentum_t, IParticle_t > self_type
publish the type of the base class (ie: 'traits-itize' it)
virtual void set4Mom(const CLHEP::HepLorentzVector &hlv)
set 4Momentum (will throw exception if cannot be implemented)
virtual double phi() const
phi in [-pi,pi[
virtual ParticleDataType::DataType dataType() const
Return enum indicating real data, fast, or full simulation Return Type has a DataType enum with the f...
virtual I4Momentum::Kind kind() const
add Kind which tells what kind (P4XYZ) the underlying implementation has
virtual double py() const
y component of momentum
bool hasSameAthenaBarCodeExceptVersion(const IAthenaBarCode &) const
virtual double eta() const
pseudo rapidity
const momentum_type & momentumBase() const
access to underlying base type (I4Momentum-like)
virtual CLHEP::HepLorentzVector hlv() const
CLHEP HepLorentzVector.
virtual double pz() const
z component of momentum
virtual void set_charge(ChargeType x)
ParticleImpl()
Default constructor.
IParticle_t particle_type
publish the type of the base class (ie: 'traits-itize' it)
const particle_type & particleBase() const
access to underlying base type (IParticle-like)
const navigable_type & navigableBase() const
access to underlying base type (INavigable-like)
double charge(const T &p)
virtual double p2() const
square of momentum magnitude
std::ostream & operator<<(std::ostream &out, const ParticleImpl< INavigable_t, I4Momentum_t, IParticle_t > &p4)
double ChargeType
typedef ChargeType used to anticipate changes here
virtual double et() const
transverse energy defined to be e*sin(theta)
virtual double m2() const
mass squared
virtual double sinPhi() const
sinus phi
virtual void set4Mom(const I4Momentum &p4)
set 4Momentum (will throw exception if cannot be implemented)
virtual void set_origin(const VxContainer *theContainer, int index)
virtual double iPt() const
inverse of transverse momentum
virtual void fillToken(INavigationToken &thisToken) const
AthenaBarCode_t AthenaBarCodeVersion_t
bool hasSameAthenaBarCode(const IAthenaBarCode &) const
virtual void set_dataType(ParticleDataType::DataType x)
I4Momentum_t m_mom
The 4-momentum part.
navigable_type & navigableBase()
access to underlying base type (INavigable-like)
virtual double cosTh() const
cosinus theta
ParticleImpl(const ParticleImpl &rhs)
Copy constructor.
virtual bool hasPdgId() const
method to check if particle id information is available