CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Attributes | Private Attributes | Friends
TrackingParticle Class Reference

Monte Carlo truth information used for tracking validation. More...

#include <TrackingParticle.h>

Public Types

typedef int Charge
 electric charge type More...
 
typedef std::vector< SimTrack >::const_iterator g4t_iterator
 
typedef reco::GenParticleRefVector::iterator genp_iterator
 reference to reco::GenParticle More...
 
typedef math::XYZTLorentzVectorD LorentzVector
 Lorentz vector. More...
 
typedef math::XYZPointD Point
 point in the space More...
 
typedef math::PtEtaPhiMLorentzVector PolarLorentzVector
 Lorentz vector. More...
 
typedef math::XYZVectorD Vector
 point in the space More...
 

Public Member Functions

void addDecayVertex (const TrackingVertexRef &ref)
 
void addG4Track (const SimTrack &t)
 
void addGenParticle (const reco::GenParticleRef &ref)
 
Vector boostToCM () const
 Vector to boost to the particle centre of mass frame. More...
 
float charge () const
 Electric charge. Note this is taken from the first SimTrack only. More...
 
void clearDecayVertices ()
 
void clearParentVertex ()
 
double d0 () const
 dxy parameter in perigee convention (d0 = -dxy) More...
 
const TrackingVertexRefVectordecayVertices () const
 
tv_iterator decayVertices_begin () const
 
tv_iterator decayVertices_end () const
 
double dxy () const
 dxy parameter. More...
 
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 (0,0,0). More...
 
double energy () const
 Energy. Note this is taken from the first SimTrack only. More...
 
double et () const
 Transverse energy. Note this is taken from the first SimTrack only. More...
 
double eta () const
 Momentum pseudorapidity. Note this is taken from the first SimTrack only. More...
 
EncodedEventId eventId () const
 Signal source, crossing number. More...
 
g4t_iterator g4Track_begin () const
 
g4t_iterator g4Track_end () const
 
const std::vector< SimTrack > & g4Tracks () const
 
genp_iterator genParticle_begin () const
 iterators More...
 
genp_iterator genParticle_end () const
 
const reco::GenParticleRefVectorgenParticles () const
 
double lambda () const
 Lambda angle. Note this is taken from the first SimTrack only. More...
 
bool longLived () const
 is long lived? More...
 
double mass () const
 Mass. Note this is taken from the first SimTrack only. More...
 
double massSqr () const
 Mass squared. Note this is taken from the first SimTrack only. More...
 
int matchedHit () const
 
Vector momentum () const
 spatial momentum vector More...
 
double mt () const
 Transverse mass. Note this is taken from the first SimTrack only. More...
 
double mtSqr () const
 Transverse mass squared. Note this is taken from the first SimTrack only. More...
 
int numberOfHits () const
 Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separately. More...
 
int numberOfTrackerHits () const
 The number of hits in the tracker. Hits on overlaps in the same layer count separately. More...
 
int numberOfTrackerLayers () const
 The number of tracker layers with a hit. More...
 
double p () const
 Magnitude of momentum vector. Note this is taken from the first SimTrack only. More...
 
const LorentzVectorp4 () const
 Four-momentum Lorentz vector. Note this is taken from the first SimTrack only. More...
 
const TrackingVertexRefparentVertex () const
 
int pdgId () const
 PDG ID. More...
 
double phi () const
 Momentum azimuthal angle. Note this is taken from the first SimTrack only. More...
 
double pt () const
 Transverse momentum. Note this is taken from the first SimTrack only. More...
 
double px () const
 x coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
double py () const
 y coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
double pz () const
 z coordinate of momentum vector. Note this is taken from the first SimTrack only. More...
 
double qoverp () const
 Quotient of the electric charge over the magnitude of the momentum vector. Note this is taken from the first SimTrack only. More...
 
double rapidity () const
 Rapidity. Note this is taken from the first SimTrack only. More...
 
void setNumberOfHits (int numberOfHits)
 
void setNumberOfTrackerHits (int numberOfTrackerHits)
 
void setNumberOfTrackerLayers (const int numberOfTrackerLayers)
 
void setParentVertex (const TrackingVertexRef &ref)
 
int status () const
 Status word. More...
 
double tanl () const
 tangent of the lambda angle. Note this is taken from the first SimTrack only. More...
 
double theta () const
 Momentum polar angle. Note this is taken from the first SimTrack only. More...
 
int threeCharge () const
 Gives charge in unit of quark charge (should be 3 times "charge()") More...
 
 TrackingParticle ()
 Default constructor. Note that the object will be useless until it is provided with a SimTrack and parent TrackingVertex. More...
 
 TrackingParticle (const SimTrack &simtrk, const TrackingVertexRef &parentVertex)
 
Point vertex () const
 Parent vertex position. More...
 
double vx () const
 x coordinate of parent vertex position More...
 
double vy () const
 y coordinate of parent vertex position More...
 
double vz () const
 z coordinate of parent vertex position More...
 
double y () const
 Same as rapidity(). More...
 
double z0 () const
 z0 parameter More...
 
 ~TrackingParticle ()
 

Static Public Attributes

static const unsigned int longLivedTag = 65536
 long lived flag More...
 

Private Attributes

TrackingVertexRefVector decayVertices_
 
std::vector< SimTrackg4Tracks_
 references to G4 and reco::GenParticle tracks More...
 
reco::GenParticleRefVector genParticles_
 
int numberOfHits_
 The total number of hits. More...
 
int numberOfTrackerHits_
 The number of tracker only hits. More...
 
int numberOfTrackerLayers_
 The number of tracker layers with hits. Equivalent to the old matchedHit. More...
 
TrackingVertexRef parentVertex_
 

Friends

std::ostream & operator<< (std::ostream &s, TrackingParticle const &tp)
 

Detailed Description

Monte Carlo truth information used for tracking validation.

Object with references to the original SimTrack and parent and daughter TrackingVertices. Simulation with high (~100) pileup was taking too much memory so the class was slimmed down and copies of the SimHits were removed.

Author
original author unknown, re-engineering and slimming by Subir Sarkar (subir.nosp@m..sar.nosp@m.kar@c.nosp@m.ern..nosp@m.ch), some tweaking and documentation by Mark Grimes (mark..nosp@m.grim.nosp@m.es@br.nosp@m.isto.nosp@m.l.ac..nosp@m.uk).
Date
original date unknown, re-engineering Jan-May 2013

Definition at line 29 of file TrackingParticle.h.

Member Typedef Documentation

◆ Charge

electric charge type

Definition at line 33 of file TrackingParticle.h.

◆ g4t_iterator

typedef std::vector<SimTrack>::const_iterator TrackingParticle::g4t_iterator

Definition at line 41 of file TrackingParticle.h.

◆ genp_iterator

reference to reco::GenParticle

Definition at line 40 of file TrackingParticle.h.

◆ LorentzVector

Lorentz vector.

Definition at line 34 of file TrackingParticle.h.

◆ Point

point in the space

Definition at line 36 of file TrackingParticle.h.

◆ PolarLorentzVector

Lorentz vector.

Definition at line 35 of file TrackingParticle.h.

◆ Vector

point in the space

Definition at line 37 of file TrackingParticle.h.

Constructor & Destructor Documentation

◆ TrackingParticle() [1/2]

TrackingParticle::TrackingParticle ( )

Default constructor. Note that the object will be useless until it is provided with a SimTrack and parent TrackingVertex.

Most of the methods assume there is a SimTrack and parent TrackingVertex set, so will either crash or give undefined results if this isn't true. This constructor should only be used to create a placeholder until setParentVertex() and addG4Track() can be called.

Definition at line 10 of file TrackingParticle.cc.

10  {
11  // No operation
12 }

◆ TrackingParticle() [2/2]

TrackingParticle::TrackingParticle ( const SimTrack simtrk,
const TrackingVertexRef parentVertex 
)

Definition at line 14 of file TrackingParticle.cc.

References addG4Track(), parentVertex(), and setParentVertex().

14  {
15  addG4Track(simtrk);
17 }
void setParentVertex(const TrackingVertexRef &ref)
void addG4Track(const SimTrack &t)
const TrackingVertexRef & parentVertex() const

◆ ~TrackingParticle()

TrackingParticle::~TrackingParticle ( )

Definition at line 19 of file TrackingParticle.cc.

19 {}

Member Function Documentation

◆ addDecayVertex()

void TrackingParticle::addDecayVertex ( const TrackingVertexRef ref)

Definition at line 35 of file TrackingParticle.cc.

References decayVertices_, and edm::RefVector< C, T, F >::push_back().

35 { decayVertices_.push_back(ref); }
TrackingVertexRefVector decayVertices_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67

◆ addG4Track()

void TrackingParticle::addG4Track ( const SimTrack t)

Definition at line 23 of file TrackingParticle.cc.

References g4Tracks_, and submitPVValidationJobs::t.

Referenced by TrackingParticle().

23 { g4Tracks_.push_back(t); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ addGenParticle()

void TrackingParticle::addGenParticle ( const reco::GenParticleRef ref)

Definition at line 21 of file TrackingParticle.cc.

References genParticles_, and edm::RefVector< C, T, F >::push_back().

21 { genParticles_.push_back(ref); }
reco::GenParticleRefVector genParticles_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67

◆ boostToCM()

Vector TrackingParticle::boostToCM ( ) const
inline

Vector to boost to the particle centre of mass frame.

Definition at line 109 of file TrackingParticle.h.

References p4().

109 { return p4().BoostToCM(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ charge()

float TrackingParticle::charge ( void  ) const
inline

◆ clearDecayVertices()

void TrackingParticle::clearDecayVertices ( )

Definition at line 39 of file TrackingParticle.cc.

References edm::RefVector< C, T, F >::clear(), and decayVertices_.

39 { decayVertices_.clear(); }
TrackingVertexRefVector decayVertices_
void clear()
Clear the vector.
Definition: RefVector.h:142

◆ clearParentVertex()

void TrackingParticle::clearParentVertex ( )

Definition at line 37 of file TrackingParticle.cc.

References parentVertex_.

edm::Ref< TrackingVertexCollection > TrackingVertexRef
TrackingVertexRef parentVertex_

◆ d0()

double TrackingParticle::d0 ( ) const
inline

dxy parameter in perigee convention (d0 = -dxy)

Definition at line 196 of file TrackingParticle.h.

References dxy().

196 { return -dxy(); }
double dxy() const
dxy parameter.

◆ decayVertices()

const TrackingVertexRefVector& TrackingParticle::decayVertices ( ) const
inline

Definition at line 93 of file TrackingParticle.h.

References decayVertices_.

93 { return decayVertices_; }
TrackingVertexRefVector decayVertices_

◆ decayVertices_begin()

tv_iterator TrackingParticle::decayVertices_begin ( ) const
inline

Definition at line 94 of file TrackingParticle.h.

References edm::RefVector< C, T, F >::begin(), and decayVertices_.

94 { return decayVertices_.begin(); }
TrackingVertexRefVector decayVertices_
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ decayVertices_end()

tv_iterator TrackingParticle::decayVertices_end ( ) const
inline

Definition at line 95 of file TrackingParticle.h.

References decayVertices_, and edm::RefVector< C, T, F >::end().

95 { return decayVertices_.end(); }
TrackingVertexRefVector decayVertices_
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228

◆ dxy()

double TrackingParticle::dxy ( ) const
inline

dxy parameter.

Definition at line 193 of file TrackingParticle.h.

References pt(), px(), py(), vx(), and vy().

Referenced by d0(), and ntupleDataFormat.Track::dxyPull().

193 { return (-vx() * py() + vy() * px()) / pt(); }
double vx() const
x coordinate of parent vertex position
double vy() const
y coordinate of parent vertex position
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.

◆ dz()

double TrackingParticle::dz ( ) const
inline

dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to (0,0,0).

Definition at line 199 of file TrackingParticle.h.

References p4(), px(), py(), pz(), vx(), vy(), and vz().

Referenced by ntupleDataFormat.Track::dzPull(), and z0().

199 { return vz() - (vx() * px() + vy() * py()) * pz() / p4().Perp2(); }
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double vx() const
x coordinate of parent vertex position
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.
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double vz() const
z coordinate of parent vertex position

◆ energy()

double TrackingParticle::energy ( void  ) const
inline

Energy. Note this is taken from the first SimTrack only.

Definition at line 118 of file TrackingParticle.h.

References p4().

Referenced by Jet.Jet::rawEnergy().

118 { return p4().E(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ et()

double TrackingParticle::et ( ) const
inline

Transverse energy. Note this is taken from the first SimTrack only.

Definition at line 121 of file TrackingParticle.h.

References p4().

121 { return p4().Et(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ eta()

double TrackingParticle::eta ( void  ) const
inline

Momentum pseudorapidity. Note this is taken from the first SimTrack only.

Definition at line 154 of file TrackingParticle.h.

References p4().

Referenced by Particle.Particle::__str__(), datamodel.Object::DeltaR(), RecoMuonValidator::MuonME::fill(), Jet.Jet::jetID(), datamodel.Object::p4(), ShallowSimTracksProducer::produce(), and Jet.Jet::puJetId().

154 { return p4().Eta(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ eventId()

EncodedEventId TrackingParticle::eventId ( ) const
inline

Signal source, crossing number.

Note this is taken from the first SimTrack only, but there shouldn't be any SimTracks from different crossings in the TrackingParticle.

Definition at line 72 of file TrackingParticle.h.

References g4Tracks_.

Referenced by MultiTrackValidator::dqmAnalyze(), ntupleDataFormat.Event::eventIdStr(), StubsSimHitsMatcher::match(), TTClusterAssociator< T >::produce(), trackingParticleIsMuonInOmtfBx0(), and trackingParticleIsMuonInOmtfEvent0().

72 { return g4Tracks_[0].eventId(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ g4Track_begin()

TrackingParticle::g4t_iterator TrackingParticle::g4Track_begin ( ) const

Definition at line 29 of file TrackingParticle.cc.

References g4Tracks_.

Referenced by QuickTrackAssociatorByHitsImpl::trackingParticleContainsIdentifier().

29 { return g4Tracks_.begin(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ g4Track_end()

TrackingParticle::g4t_iterator TrackingParticle::g4Track_end ( ) const

Definition at line 31 of file TrackingParticle.cc.

References g4Tracks_.

Referenced by QuickTrackAssociatorByHitsImpl::trackingParticleContainsIdentifier().

31 { return g4Tracks_.end(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ g4Tracks()

const std::vector<SimTrack>& TrackingParticle::g4Tracks ( ) const
inline

Definition at line 89 of file TrackingParticle.h.

References g4Tracks_.

Referenced by StubsSimHitsMatcher::match(), TTClusterAssociator< T >::produce(), l1tVertexFinder::TP::TP(), tmtt::TP::TP(), and trackingParticleIsMuonInOmtfBx0().

89 { return g4Tracks_; }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ genParticle_begin()

TrackingParticle::genp_iterator TrackingParticle::genParticle_begin ( ) const

iterators

Definition at line 25 of file TrackingParticle.cc.

References edm::RefVector< C, T, F >::begin(), and genParticles_.

25 { return genParticles_.begin(); }
reco::GenParticleRefVector genParticles_
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ genParticle_end()

TrackingParticle::genp_iterator TrackingParticle::genParticle_end ( ) const

Definition at line 27 of file TrackingParticle.cc.

References edm::RefVector< C, T, F >::end(), and genParticles_.

27 { return genParticles_.end(); }
reco::GenParticleRefVector genParticles_
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228

◆ genParticles()

const reco::GenParticleRefVector& TrackingParticle::genParticles ( ) const
inline

Definition at line 88 of file TrackingParticle.h.

References genParticles_.

88 { return genParticles_; }
reco::GenParticleRefVector genParticles_

◆ lambda()

double TrackingParticle::lambda ( ) const
inline

Lambda angle. Note this is taken from the first SimTrack only.

Definition at line 157 of file TrackingParticle.h.

References M_PI_2, and theta().

Referenced by tanl().

157 { return M_PI_2 - theta(); }
#define M_PI_2
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.

◆ longLived()

bool TrackingParticle::longLived ( ) const
inline

is long lived?

Definition at line 212 of file TrackingParticle.h.

References longLivedTag, and status().

212 { return status() & longLivedTag; }
int status() const
Status word.
static const unsigned int longLivedTag
long lived flag

◆ mass()

double TrackingParticle::mass ( ) const
inline

Mass. Note this is taken from the first SimTrack only.

Definition at line 124 of file TrackingParticle.h.

References p4().

Referenced by Particle.Particle::__str__(), DiObject.DiMuon::__str__(), massSqr(), and datamodel.Object::p4().

124 { return p4().M(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ massSqr()

double TrackingParticle::massSqr ( ) const
inline

Mass squared. Note this is taken from the first SimTrack only.

Definition at line 127 of file TrackingParticle.h.

References mass(), and funct::pow().

127 { return pow(mass(), 2); }
double mass() const
Mass. Note this is taken from the first SimTrack only.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ matchedHit()

int TrackingParticle::matchedHit ( ) const

Definition at line 41 of file TrackingParticle.cc.

References numberOfTrackerLayers_.

41  {
42  edm::LogWarning("TrackingParticle")
43  << "The method matchedHit() has been deprecated. Use numberOfTrackerLayers() instead.";
45 }
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
Log< level::Warning, false > LogWarning

◆ momentum()

Vector TrackingParticle::momentum ( ) const
inline

spatial momentum vector

Definition at line 106 of file TrackingParticle.h.

References p4().

Referenced by MuonTrackValidator::analyze(), RecoMuonValidator::MuonME::fill(), CandidateSimMuonMatcher::match(), MatchingResult::MatchingResult(), CandidateSimMuonMatcher::propagate(), CandidateSimMuonMatcher::simTrackToFts(), track_associator::trackAssociationChi2(), and trackingParticleIsMuonInOmtfBx0().

106 { return p4().Vect(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ mt()

double TrackingParticle::mt ( ) const
inline

Transverse mass. Note this is taken from the first SimTrack only.

Definition at line 130 of file TrackingParticle.h.

References p4().

130 { return p4().Mt(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ mtSqr()

double TrackingParticle::mtSqr ( ) const
inline

Transverse mass squared. Note this is taken from the first SimTrack only.

Definition at line 133 of file TrackingParticle.h.

References p4().

133 { return p4().Mt2(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ numberOfHits()

int TrackingParticle::numberOfHits ( ) const
inline

Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separately.

Equivalent to trackPSimHit().size() in the old TrackingParticle implementation.

Definition at line 217 of file TrackingParticle.h.

References numberOfHits_.

Referenced by setNumberOfHits().

217 { return numberOfHits_; }
int numberOfHits_
The total number of hits.

◆ numberOfTrackerHits()

int TrackingParticle::numberOfTrackerHits ( ) const
inline

The number of hits in the tracker. Hits on overlaps in the same layer count separately.

Equivalent to trackPSimHit(DetId::Tracker).size() in the old TrackingParticle implementation.

Definition at line 222 of file TrackingParticle.h.

References numberOfTrackerHits_.

Referenced by setNumberOfTrackerHits().

222 { return numberOfTrackerHits_; }
int numberOfTrackerHits_
The number of tracker only hits.

◆ numberOfTrackerLayers()

int TrackingParticle::numberOfTrackerLayers ( ) const
inline

The number of tracker layers with a hit.

Different from numberOfTrackerHits because this method counts multiple hits on overlaps in the layer as one hit.

Definition at line 231 of file TrackingParticle.h.

References numberOfTrackerLayers_.

Referenced by setNumberOfTrackerLayers().

231 { return numberOfTrackerLayers_; }
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.

◆ p()

double TrackingParticle::p ( ) const
inline

Magnitude of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 112 of file TrackingParticle.h.

References p4().

Referenced by RecoMuonValidator::MuonME::fill(), ShallowSimTracksProducer::produce(), Electron.Electron::ptErr(), qoverp(), and vertex().

112 { return p4().P(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ p4()

const LorentzVector& TrackingParticle::p4 ( ) const
inline

Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

Definition at line 103 of file TrackingParticle.h.

References g4Tracks_.

Referenced by OuterTrackerMonitorTrackingParticles::analyze(), boostToCM(), Tau.Tau::dxy_approx(), Tau.Tau::dz(), dz(), energy(), et(), eta(), TTClusterAssociationMap< T >::isGenuine(), mass(), momentum(), mt(), mtSqr(), p(), Lepton.Lepton::p4WithFSR(), phi(), pt(), px(), py(), pz(), rapidity(), and theta().

103 { return g4Tracks_[0].momentum(); }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks

◆ parentVertex()

const TrackingVertexRef& TrackingParticle::parentVertex ( ) const
inline

◆ pdgId()

int TrackingParticle::pdgId ( ) const
inline

PDG ID.

Returns the PDG ID of the first associated gen particle. If there are no gen particles associated then it returns type() from the first SimTrack.

Definition at line 61 of file TrackingParticle.h.

References edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::empty(), g4Tracks_, and genParticles_.

Referenced by Particle.Particle::__str__(), OuterTrackerMonitorTrackingParticles::analyze(), StubsSimHitsMatcher::match(), CandidateSimMuonMatcher::match(), MatchingResult::MatchingResult(), ShallowSimTracksProducer::produce(), CandidateSimMuonMatcher::simTrackToFts(), and trackingParticleIsMuonInOmtfBx0().

61  {
62  if (genParticles_.empty())
63  return g4Tracks_[0].type();
64  else
65  return (*genParticles_.begin())->pdgId();
66  }
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
int pdgId() const
PDG ID.
reco::GenParticleRefVector genParticles_
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223

◆ phi()

double TrackingParticle::phi ( void  ) const
inline

Momentum azimuthal angle. Note this is taken from the first SimTrack only.

Definition at line 148 of file TrackingParticle.h.

References p4().

Referenced by Particle.Particle::__str__(), datamodel.Object::DeltaR(), RecoMuonValidator::MuonME::fill(), datamodel.Object::p4(), ntupleDataFormat.Track::phiPull(), and ShallowSimTracksProducer::produce().

148 { return p4().Phi(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ pt()

double TrackingParticle::pt ( ) const
inline

◆ px()

double TrackingParticle::px ( ) const
inline

x coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 136 of file TrackingParticle.h.

References p4().

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), dxy(), and dz().

136 { return p4().Px(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ py()

double TrackingParticle::py ( ) const
inline

y coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 139 of file TrackingParticle.h.

References p4().

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), dxy(), and dz().

139 { return p4().Py(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ pz()

double TrackingParticle::pz ( ) const
inline

z coordinate of momentum vector. Note this is taken from the first SimTrack only.

Definition at line 142 of file TrackingParticle.h.

References p4().

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), and dz().

142 { return p4().Pz(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ qoverp()

double TrackingParticle::qoverp ( ) const
inline

Quotient of the electric charge over the magnitude of the momentum vector. Note this is taken from the first SimTrack only.

Definition at line 115 of file TrackingParticle.h.

References charge(), and p().

115 { return charge() / p(); }
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
float charge() const
Electric charge. Note this is taken from the first SimTrack only.

◆ rapidity()

double TrackingParticle::rapidity ( ) const
inline

Rapidity. Note this is taken from the first SimTrack only.

Definition at line 163 of file TrackingParticle.h.

References p4().

Referenced by y().

163 { return p4().Rapidity(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ setNumberOfHits()

void TrackingParticle::setNumberOfHits ( int  numberOfHits)

Definition at line 47 of file TrackingParticle.cc.

References numberOfHits(), and numberOfHits_.

int numberOfHits_
The total number of hits.
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...

◆ setNumberOfTrackerHits()

void TrackingParticle::setNumberOfTrackerHits ( int  numberOfTrackerHits)

Definition at line 49 of file TrackingParticle.cc.

References numberOfTrackerHits(), and numberOfTrackerHits_.

int numberOfTrackerHits_
The number of tracker only hits.
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.

◆ setNumberOfTrackerLayers()

void TrackingParticle::setNumberOfTrackerLayers ( const int  numberOfTrackerLayers)

Definition at line 51 of file TrackingParticle.cc.

References numberOfTrackerLayers(), and numberOfTrackerLayers_.

51  {
53 }
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.

◆ setParentVertex()

void TrackingParticle::setParentVertex ( const TrackingVertexRef ref)

Definition at line 33 of file TrackingParticle.cc.

References parentVertex_.

Referenced by TrackingParticle().

33 { parentVertex_ = ref; }
TrackingVertexRef parentVertex_

◆ status()

int TrackingParticle::status ( void  ) const
inline

Status word.

Returns status() from the first gen particle, or -99 if there are no gen particles attached.

Definition at line 207 of file TrackingParticle.h.

References edm::RefVector< C, T, F >::empty(), genParticles_, and status().

Referenced by longLived(), and status().

207 { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
int status() const
Status word.
reco::GenParticleRefVector genParticles_

◆ tanl()

double TrackingParticle::tanl ( ) const
inline

tangent of the lambda angle. Note this is taken from the first SimTrack only.

Definition at line 160 of file TrackingParticle.h.

References lambda(), and funct::tan().

160 { return tan(lambda()); }
double lambda() const
Lambda angle. Note this is taken from the first SimTrack only.
Tan< T >::type tan(const T &t)
Definition: Tan.h:22

◆ theta()

double TrackingParticle::theta ( void  ) const
inline

Momentum polar angle. Note this is taken from the first SimTrack only.

Definition at line 151 of file TrackingParticle.h.

References p4().

Referenced by lambda(), ShallowSimTracksProducer::produce(), and Tau.Tau::zImpact().

151 { return p4().Theta(); }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.

◆ threeCharge()

int TrackingParticle::threeCharge ( ) const
inline

Gives charge in unit of quark charge (should be 3 times "charge()")

Definition at line 100 of file TrackingParticle.h.

References charge(), and f.

100 { return lrintf(3.f * charge()); }
double f[11][100]
float charge() const
Electric charge. Note this is taken from the first SimTrack only.

◆ vertex()

Point TrackingParticle::vertex ( ) const
inline

Parent vertex position.

Definition at line 169 of file TrackingParticle.h.

References p().

Referenced by Tau.Tau::dxy(), RecoMuonValidator::MuonME::fill(), and track_associator::trackAssociationChi2().

169  {
170  const TrackingVertex::LorentzVector& p = (*parentVertex_).position();
171  return Point(p.x(), p.y(), p.z());
172  }
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
math::XYZPointD Point
point in the space
math::XYZTLorentzVectorD LorentzVector

◆ vx()

double TrackingParticle::vx ( ) const
inline

x coordinate of parent vertex position

Definition at line 175 of file TrackingParticle.h.

References alignCSCRings::r.

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), dxy(), dz(), and CandidateSimMuonMatcher::simTrackToFts().

175  {
176  const TrackingVertex& r = (*parentVertex_);
177  return r.position().X();
178  }

◆ vy()

double TrackingParticle::vy ( ) const
inline

y coordinate of parent vertex position

Definition at line 181 of file TrackingParticle.h.

References alignCSCRings::r.

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), dxy(), dz(), and CandidateSimMuonMatcher::simTrackToFts().

181  {
182  const TrackingVertex& r = (*parentVertex_);
183  return r.position().Y();
184  }

◆ vz()

double TrackingParticle::vz ( ) const
inline

z coordinate of parent vertex position

Definition at line 187 of file TrackingParticle.h.

References alignCSCRings::r.

Referenced by FWTrackingParticleProxyBuilder::build(), FWTrackingParticleProxyBuilderFullFramework::build(), dz(), and CandidateSimMuonMatcher::simTrackToFts().

187  {
188  const TrackingVertex& r = (*parentVertex_);
189  return r.position().Z();
190  }

◆ y()

double TrackingParticle::y ( ) const
inline

Same as rapidity().

Definition at line 166 of file TrackingParticle.h.

References rapidity().

Referenced by svgfig.Ellipse::__repr__(), geometryXMLparser.Alignable::pos(), ntupleDataFormat._HitObject::r(), and ntupleDataFormat._HitObject::r3D().

166 { return rapidity(); }
double rapidity() const
Rapidity. Note this is taken from the first SimTrack only.

◆ z0()

double TrackingParticle::z0 ( ) const
inline

z0 parameter

Definition at line 202 of file TrackingParticle.h.

References dz().

202 { return dz(); }
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...

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  s,
TrackingParticle const &  tp 
)
friend

Definition at line 55 of file TrackingParticle.cc.

55  {
56  s << "TP momentum, q, ID, & Event #: " << tp.p4() << " " << tp.charge() << " " << tp.pdgId() << " "
57  << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
58 
59  for (TrackingParticle::genp_iterator hepT = tp.genParticle_begin(); hepT != tp.genParticle_end(); ++hepT) {
60  s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
61  }
62 
63  for (TrackingParticle::g4t_iterator g4T = tp.g4Track_begin(); g4T != tp.g4Track_end(); ++g4T) {
64  s << " Geant Track Momentum " << g4T->momentum() << std::endl;
65  s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
66  if (g4T->type() != tp.pdgId()) {
67  s << " Mismatch b/t TrackingParticle and Geant types" << std::endl;
68  }
69  }
70  // Loop over decay vertices
71  s << " TP Vertex " << tp.vertex() << std::endl;
72  s << " Source vertex: " << tp.parentVertex()->position() << std::endl;
73  s << " " << tp.decayVertices().size() << " Decay vertices" << std::endl;
74  for (tv_iterator iTV = tp.decayVertices_begin(); iTV != tp.decayVertices_end(); ++iTV) {
75  s << " Decay vertices: " << (**iTV).position() << std::endl;
76  }
77 
78  return s;
79 }
std::vector< SimTrack >::const_iterator g4t_iterator

Member Data Documentation

◆ decayVertices_

TrackingVertexRefVector TrackingParticle::decayVertices_
private

◆ g4Tracks_

std::vector<SimTrack> TrackingParticle::g4Tracks_
private

references to G4 and reco::GenParticle tracks

Definition at line 243 of file TrackingParticle.h.

Referenced by addG4Track(), charge(), eventId(), g4Track_begin(), g4Track_end(), g4Tracks(), p4(), and pdgId().

◆ genParticles_

reco::GenParticleRefVector TrackingParticle::genParticles_
private

◆ longLivedTag

const unsigned int TrackingParticle::longLivedTag = 65536
static

long lived flag

Definition at line 209 of file TrackingParticle.h.

Referenced by longLived().

◆ numberOfHits_

int TrackingParticle::numberOfHits_
private

The total number of hits.

Definition at line 238 of file TrackingParticle.h.

Referenced by numberOfHits(), and setNumberOfHits().

◆ numberOfTrackerHits_

int TrackingParticle::numberOfTrackerHits_
private

The number of tracker only hits.

Definition at line 239 of file TrackingParticle.h.

Referenced by numberOfTrackerHits(), and setNumberOfTrackerHits().

◆ numberOfTrackerLayers_

int TrackingParticle::numberOfTrackerLayers_
private

The number of tracker layers with hits. Equivalent to the old matchedHit.

Definition at line 240 of file TrackingParticle.h.

Referenced by matchedHit(), numberOfTrackerLayers(), and setNumberOfTrackerLayers().

◆ parentVertex_

TrackingVertexRef TrackingParticle::parentVertex_
private

Definition at line 247 of file TrackingParticle.h.

Referenced by clearParentVertex(), parentVertex(), and setParentVertex().