Go to the documentation of this file. 1 #ifndef SimDataFormats_TrackingParticle_h
2 #define SimDataFormats_TrackingParticle_h
112 double p()
const {
return p4().P(); }
121 double et()
const {
return p4().Et(); }
130 double mt()
const {
return p4().Mt(); }
136 double px()
const {
return p4().Px(); }
139 double py()
const {
return p4().Py(); }
142 double pz()
const {
return p4().Pz(); }
145 double pt()
const {
return p4().Pt(); }
148 double phi()
const {
return p4().Phi(); }
154 double eta()
const {
return p4().Eta(); }
177 return r.position().X();
183 return r.position().Y();
189 return r.position().Z();
196 double d0()
const {
return -
dxy(); }
202 double z0()
const {
return dz(); }
251 #endif // SimDataFormats_TrackingParticle_H
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
void addGenParticle(const reco::GenParticleRef &ref)
double eta() const
Momentum pseudorapidity. Note this is taken from the first SimTrack only.
g4t_iterator g4Track_begin() const
double mass() const
Mass. Note this is taken from the first SimTrack only.
double mt() const
Transverse mass. Note this is taken from the first SimTrack only.
double tanl() const
tangent of the lambda angle. Note this is taken from the first SimTrack only.
void clearDecayVertices()
TrackingParticle::g4t_iterator g4t_iterator
int numberOfHits_
The total number of hits.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
void addDecayVertex(const TrackingVertexRef &ref)
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
const_iterator begin() const
Initialize an iterator over the RefVector.
std::vector< SimTrack >::const_iterator g4t_iterator
double vy() const
y coordinate of parent vertex position
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
const TrackingVertexRefVector & decayVertices() const
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only.
const TrackingVertexRef & parentVertex() const
math::XYZVectorD Vector
point in the space
double massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
Vector boostToCM() const
Vector to boost to the particle centre of mass frame.
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
static const unsigned int longLivedTag
long lived flag
TrackingVertexRef parentVertex_
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only.
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
const_iterator end() const
Termination of iteration.
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
double qoverp() const
Quotient of the electric charge over the magnitude of the momentum vector. Note this is taken from th...
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...
bool empty() const
Is the RefVector empty.
genp_iterator genParticle_end() const
double dxy() const
dxy parameter.
void setNumberOfTrackerHits(int numberOfTrackerHits)
Monte Carlo truth information used for tracking validation.
double phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
friend std::ostream & operator<<(std::ostream &s, TrackingParticle const &tp)
Vector momentum() const
spatial momentum vector
Point vertex() const
Parent vertex position.
double mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
int Charge
electric charge type
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
const reco::GenParticleRefVector & genParticles() const
Structure Point Contains parameters of Gaussian fits to DMRs.
double z0() const
z0 parameter
int status() const
Status word.
TrackingVertexRefVector decayVertices_
double vz() const
z coordinate of parent vertex position
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
tv_iterator decayVertices_begin() const
Tan< T >::type tan(const T &t)
tv_iterator decayVertices_end() const
void setNumberOfHits(int numberOfHits)
TrackingParticle()
Default constructor. Note that the object will be useless until it is provided with a SimTrack and pa...
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only.
double rapidity() const
Rapidity. Note this is taken from the first SimTrack only.
double et() const
Transverse energy. Note this is taken from the first SimTrack only.
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
EncodedEventId eventId() const
Signal source, crossing number.
double lambda() const
Lambda angle. Note this is taken from the first SimTrack only.
void setParentVertex(const TrackingVertexRef &ref)
double vx() const
x coordinate of parent vertex position
reco::GenParticleRefVector genParticles_
bool longLived() const
is long lived?
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
g4t_iterator g4Track_end() const
Power< A, B >::type pow(const A &a, const B &b)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
math::XYZPointD Point
point in the space
int numberOfTrackerHits_
The number of tracker only hits.
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
double y() const
Same as rapidity().
double energy() const
Energy. Note this is taken from the first SimTrack only.
void addG4Track(const SimTrack &t)
math::XYZTLorentzVectorD LorentzVector
const std::vector< SimTrack > & g4Tracks() const
genp_iterator genParticle_begin() const
iterators