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 ()
 
const TrackingVertexRefVectordecayVertices () const
 
tv_iterator decayVertices_begin () const
 
tv_iterator decayVertices_end () const
 
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
 
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 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 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
 
double y () const
 Same as rapidity(). 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

electric charge type

Definition at line 33 of file TrackingParticle.h.

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

Definition at line 41 of file TrackingParticle.h.

reference to reco::GenParticle

Definition at line 40 of file TrackingParticle.h.

Lorentz vector.

Definition at line 34 of file TrackingParticle.h.

point in the space

Definition at line 36 of file TrackingParticle.h.

Lorentz vector.

Definition at line 35 of file TrackingParticle.h.

point in the space

Definition at line 37 of file TrackingParticle.h.

Constructor & Destructor Documentation

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.

11 {
12  // No operation
13 }
TrackingParticle::TrackingParticle ( const SimTrack simtrk,
const TrackingVertexRef parentVertex 
)

Definition at line 15 of file TrackingParticle.cc.

References addG4Track(), and setParentVertex().

16 {
17  addG4Track( simtrk );
18  setParentVertex( parentVertex );
19 }
void setParentVertex(const TrackingVertexRef &ref)
void addG4Track(const SimTrack &t)
TrackingParticle::~TrackingParticle ( )

Definition at line 21 of file TrackingParticle.cc.

22 {
23 }

Member Function Documentation

void TrackingParticle::addDecayVertex ( const TrackingVertexRef ref)

Definition at line 60 of file TrackingParticle.cc.

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

Referenced by eventId(), and TrackingTruthAccumulator::fillSimHits().

61 {
63 }
TrackingVertexRefVector decayVertices_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
void TrackingParticle::addG4Track ( const SimTrack t)

Definition at line 30 of file TrackingParticle.cc.

References g4Tracks_.

Referenced by eventId(), TrackingTruthAccumulator::fillSimHits(), and TrackingParticle().

31 {
32  g4Tracks_.push_back( t );
33 }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
void TrackingParticle::addGenParticle ( const reco::GenParticleRef ref)

Definition at line 25 of file TrackingParticle.cc.

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

Referenced by eventId(), and TrackingTruthAccumulator::fillSimHits().

26 {
27  genParticles_.push_back( ref );
28 }
reco::GenParticleRefVector genParticles_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
Vector TrackingParticle::boostToCM ( ) const
inline

Vector to boost to the particle centre of mass frame.

Definition at line 114 of file TrackingParticle.h.

References p4().

114  {
115  return p4().BoostToCM();
116  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
float TrackingParticle::charge ( void  ) const
inline
void TrackingParticle::clearDecayVertices ( )

Definition at line 70 of file TrackingParticle.cc.

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

Referenced by eventId().

71 {
73 }
TrackingVertexRefVector decayVertices_
void clear()
Clear the vector.
Definition: RefVector.h:147
void TrackingParticle::clearParentVertex ( )

Definition at line 65 of file TrackingParticle.cc.

References parentVertex_.

Referenced by eventId().

66 {
68 }
edm::Ref< TrackingVertexCollection > TrackingVertexRef
TrackingVertexRef parentVertex_
const TrackingVertexRefVector& TrackingParticle::decayVertices ( ) const
inline

Definition at line 93 of file TrackingParticle.h.

References decayVertices_.

Referenced by operator<<().

93 { return decayVertices_; }
TrackingVertexRefVector decayVertices_
tv_iterator TrackingParticle::decayVertices_begin ( ) const
inline

Definition at line 94 of file TrackingParticle.h.

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

Referenced by operator<<().

94 { return decayVertices_.begin(); }
TrackingVertexRefVector decayVertices_
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
tv_iterator TrackingParticle::decayVertices_end ( ) const
inline

Definition at line 95 of file TrackingParticle.h.

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

Referenced by operator<<().

95 { return decayVertices_.end(); }
TrackingVertexRefVector decayVertices_
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
double TrackingParticle::energy ( void  ) const
inline

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

Definition at line 124 of file TrackingParticle.h.

References p4().

Referenced by Jet.Jet::rawEnergy().

124  {
125  return p4().E();
126  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::et ( ) const
inline

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

Definition at line 129 of file TrackingParticle.h.

References p4().

129  {
130  return p4().Et();
131  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::eta ( void  ) const
inline

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

Definition at line 184 of file TrackingParticle.h.

References p4().

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

184  {
185  return p4().Eta();
186  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
EncodedEventId TrackingParticle::eventId ( ) const
inline
TrackingParticle::g4t_iterator TrackingParticle::g4Track_begin ( ) const

Definition at line 45 of file TrackingParticle.cc.

References g4Tracks_.

Referenced by eventId(), operator<<(), and QuickTrackAssociatorByHitsImpl::trackingParticleContainsIdentifier().

46 {
47  return g4Tracks_.begin();
48 }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
TrackingParticle::g4t_iterator TrackingParticle::g4Track_end ( ) const

Definition at line 50 of file TrackingParticle.cc.

References g4Tracks_.

Referenced by eventId(), operator<<(), and QuickTrackAssociatorByHitsImpl::trackingParticleContainsIdentifier().

51 {
52  return g4Tracks_.end();
53 }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
const std::vector<SimTrack>& TrackingParticle::g4Tracks ( ) const
inline

Definition at line 89 of file TrackingParticle.h.

References g4Tracks_.

Referenced by TrackingParticleNumberOfLayers::calculate(), TrackingTruthAccumulator::fillSimHits(), and TTClusterAssociator< T >::produce().

89 { return g4Tracks_; }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
TrackingParticle::genp_iterator TrackingParticle::genParticle_begin ( ) const

iterators

Definition at line 35 of file TrackingParticle.cc.

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

Referenced by eventId(), TrackingParticleSelector::operator()(), and operator<<().

36 {
37  return genParticles_.begin();
38 }
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
reco::GenParticleRefVector genParticles_
TrackingParticle::genp_iterator TrackingParticle::genParticle_end ( ) const

Definition at line 40 of file TrackingParticle.cc.

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

Referenced by eventId(), TrackingParticleSelector::operator()(), and operator<<().

41 {
42  return genParticles_.end();
43 }
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
reco::GenParticleRefVector genParticles_
const reco::GenParticleRefVector& TrackingParticle::genParticles ( ) const
inline

Definition at line 88 of file TrackingParticle.h.

References genParticles_.

Referenced by TrackingTruthAccumulator::fillSimHits().

88 { return genParticles_; }
reco::GenParticleRefVector genParticles_
bool TrackingParticle::longLived ( ) const
inline

is long lived?

Definition at line 231 of file TrackingParticle.h.

References longLivedTag, and status().

231 { return status()&longLivedTag;}
int status() const
Status word.
static const unsigned int longLivedTag
long lived flag
double TrackingParticle::mass ( ) const
inline

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

Definition at line 134 of file TrackingParticle.h.

References p4().

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

134  {
135  return p4().M();
136  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::massSqr ( ) const
inline

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

Definition at line 139 of file TrackingParticle.h.

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

139  {
140  return pow( mass(), 2 );
141  }
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:40
int TrackingParticle::matchedHit ( ) const

Definition at line 75 of file TrackingParticle.cc.

References numberOfTrackerLayers_.

Referenced by numberOfTrackerHits().

76 {
77  edm::LogWarning("TrackingParticle") << "The method matchedHit() has been deprecated. Use numberOfTrackerLayers() instead.";
79 }
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
Vector TrackingParticle::momentum ( ) const
inline

spatial momentum vector

Definition at line 109 of file TrackingParticle.h.

References p4().

Referenced by MuonTrackValidator::analyze(), MultiTrackValidator::dqmAnalyze(), RecoMuonValidator::MuonME::fill(), MultiTrackValidator::tpParametersAndSelection(), and track_associator::trackAssociationChi2().

109  {
110  return p4().Vect();
111  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::mt ( ) const
inline

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

Definition at line 144 of file TrackingParticle.h.

References p4().

144  {
145  return p4().Mt();
146  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::mtSqr ( ) const
inline

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

Definition at line 149 of file TrackingParticle.h.

References p4().

149  {
150  return p4().Mt2();
151  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
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 236 of file TrackingParticle.h.

References numberOfHits_.

Referenced by TrackingTruthAccumulator::fillSimHits(), TrackAssociatorByHitsImpl::getShared(), numberOfTrackerLayers(), and setNumberOfHits().

236 {return numberOfHits_;}
int numberOfHits_
The total number of hits.
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 241 of file TrackingParticle.h.

References matchedHit(), and numberOfTrackerHits_.

Referenced by MultiTrackValidator::dqmAnalyze(), TrackingTruthAccumulator::fillSimHits(), numberOfTrackerLayers(), and setNumberOfTrackerHits().

241 {return numberOfTrackerHits_;}
int numberOfTrackerHits_
The number of tracker only hits.
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 250 of file TrackingParticle.h.

References numberOfHits(), numberOfTrackerHits(), numberOfTrackerLayers_, setNumberOfHits(), setNumberOfTrackerHits(), and setNumberOfTrackerLayers().

Referenced by TrackingTruthAccumulator::fillSimHits(), TrackingParticleSelector::operator()(), and setNumberOfTrackerLayers().

250 {return numberOfTrackerLayers_;}
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
double TrackingParticle::p ( ) const
inline

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

Definition at line 119 of file TrackingParticle.h.

References p4().

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

119  {
120  return p4().P();
121  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
const LorentzVector& TrackingParticle::p4 ( ) const
inline

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

Definition at line 104 of file TrackingParticle.h.

References g4Tracks_.

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

104  {
105  return g4Tracks_[0].momentum();
106  }
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
const TrackingVertexRef& TrackingParticle::parentVertex ( ) const
inline
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__(), TrackingParticleNumberOfLayers::calculate(), MuonSimClassifier::convertAndPush(), MuonMCClassifier::convertAndPush(), TrackTimeValueMapProducer::extractTrackVertexTime(), TrackingTruthAccumulator::fillSimHits(), TrackingNtuple::fillSimHits(), TrackingParticleSelector::operator()(), operator<<(), and ShallowSimTracksProducer::produce().

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

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

Definition at line 174 of file TrackingParticle.h.

References p4().

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

174  {
175  return p4().Phi();
176  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::pt ( ) const
inline
double TrackingParticle::px ( ) const
inline

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

Definition at line 154 of file TrackingParticle.h.

References p4().

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

154  {
155  return p4().Px();
156  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::py ( ) const
inline

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

Definition at line 159 of file TrackingParticle.h.

References p4().

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

159  {
160  return p4().Py();
161  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::pz ( ) const
inline

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

Definition at line 164 of file TrackingParticle.h.

References p4().

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

164  {
165  return p4().Pz();
166  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
double TrackingParticle::rapidity ( ) const
inline

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

Definition at line 189 of file TrackingParticle.h.

References p4().

Referenced by y().

189  {
190  return p4().Rapidity();
191  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
void TrackingParticle::setNumberOfHits ( int  numberOfHits)

Definition at line 82 of file TrackingParticle.cc.

References numberOfHits(), and numberOfHits_.

Referenced by TrackingTruthAccumulator::fillSimHits(), and numberOfTrackerLayers().

83 {
85 }
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...
void TrackingParticle::setNumberOfTrackerHits ( int  numberOfTrackerHits)

Definition at line 87 of file TrackingParticle.cc.

References numberOfTrackerHits(), and numberOfTrackerHits_.

Referenced by TrackingTruthAccumulator::fillSimHits(), and numberOfTrackerLayers().

88 {
90 }
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.
void TrackingParticle::setNumberOfTrackerLayers ( const int  numberOfTrackerLayers)

Definition at line 92 of file TrackingParticle.cc.

References numberOfTrackerLayers(), and numberOfTrackerLayers_.

Referenced by TrackingTruthAccumulator::fillSimHits(), and numberOfTrackerLayers().

93 {
95 }
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.
void TrackingParticle::setParentVertex ( const TrackingVertexRef ref)

Definition at line 55 of file TrackingParticle.cc.

References parentVertex_.

Referenced by eventId(), TrackingTruthAccumulator::fillSimHits(), and TrackingParticle().

56 {
57  parentVertex_=ref;
58 }
TrackingVertexRef parentVertex_
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 224 of file TrackingParticle.h.

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

Referenced by MuonSimClassifier::convertAndPush(), MuonMCClassifier::convertAndPush(), longLived(), and TrackingParticleSelector::operator()().

224  {
225  return genParticles_.empty() ? -99 : (*genParticles_[0]).status();
226  }
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
int status() const
Status word.
reco::GenParticleRefVector genParticles_
double TrackingParticle::theta ( void  ) const
inline

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

Definition at line 179 of file TrackingParticle.h.

References p4().

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

179  {
180  return p4().Theta();
181  }
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
int TrackingParticle::threeCharge ( ) const
inline

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

Definition at line 101 of file TrackingParticle.h.

References charge(), and f.

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

Parent vertex position.

Definition at line 199 of file TrackingParticle.h.

References p().

Referenced by MuonTrackValidator::analyze(), MuonSimClassifier::convertAndPush(), MuonMCClassifier::convertAndPush(), MultiTrackValidator::dqmAnalyze(), Tau.Tau::dxy(), RecoMuonValidator::MuonME::fill(), TrackingParticleSelector::operator()(), operator<<(), MultiTrackValidator::tpParametersAndSelection(), and track_associator::trackAssociationChi2().

199  {
200  const TrackingVertex::LorentzVector & p = (*parentVertex_).position();
201  return Point(p.x(),p.y(),p.z());
202  }
math::XYZPointD Point
point in the space
math::XYZTLorentzVectorD LorentzVector
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
double TrackingParticle::vx ( ) const
inline

x coordinate of parent vertex position

Definition at line 205 of file TrackingParticle.h.

References TrackingVertex::position(), and alignCSCRings::r.

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

205  {
206  const TrackingVertex& r=( *parentVertex_);
207  return r.position().X();
208  }
const LorentzVector & position() const
double TrackingParticle::vy ( ) const
inline

y coordinate of parent vertex position

Definition at line 211 of file TrackingParticle.h.

References TrackingVertex::position(), and alignCSCRings::r.

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

211  {
212  const TrackingVertex& r=( *parentVertex_);
213  return r.position().Y();
214  }
const LorentzVector & position() const
double TrackingParticle::vz ( ) const
inline

Definition at line 216 of file TrackingParticle.h.

References TrackingVertex::position(), and alignCSCRings::r.

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

216  {
217  const TrackingVertex& r=( *parentVertex_);
218  return r.position().Z();
219  }
const LorentzVector & position() const
double TrackingParticle::y ( ) const
inline

Friends And Related Function Documentation

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

Definition at line 97 of file TrackingParticle.cc.

98 {
99  s << "TP momentum, q, ID, & Event #: "
100  << tp.p4() << " " << tp.charge() << " " << tp.pdgId() << " "
101  << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
102 
103  for (TrackingParticle::genp_iterator hepT = tp.genParticle_begin(); hepT != tp.genParticle_end(); ++hepT)
104  {
105  s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
106  }
107 
108  for (TrackingParticle::g4t_iterator g4T = tp.g4Track_begin(); g4T != tp.g4Track_end(); ++g4T)
109  {
110  s << " Geant Track Momentum " << g4T->momentum() << std::endl;
111  s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
112  if (g4T->type() != tp.pdgId())
113  {
114  s << " Mismatch b/t TrackingParticle and Geant types" << std::endl;
115  }
116  }
117  // Loop over decay vertices
118  s << " TP Vertex " << tp.vertex() << std::endl;
119  s << " Source vertex: " << tp.parentVertex()->position() << std::endl;
120  s << " " << tp.decayVertices().size() << " Decay vertices" << std::endl;
121  for (tv_iterator iTV = tp.decayVertices_begin(); iTV != tp.decayVertices_end(); ++iTV)
122  {
123  s << " Decay vertices: " << (**iTV).position() << std::endl;
124  }
125 
126  return s;
127 }
std::vector< SimTrack >::const_iterator g4t_iterator

Member Data Documentation

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

references to G4 and reco::GenParticle tracks

Definition at line 261 of file TrackingParticle.h.

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

reco::GenParticleRefVector TrackingParticle::genParticles_
private
const unsigned int TrackingParticle::longLivedTag = 65536
static

long lived flag

Definition at line 228 of file TrackingParticle.h.

Referenced by longLived().

int TrackingParticle::numberOfHits_
private

The total number of hits.

Definition at line 256 of file TrackingParticle.h.

Referenced by numberOfHits(), and setNumberOfHits().

int TrackingParticle::numberOfTrackerHits_
private

The number of tracker only hits.

Definition at line 257 of file TrackingParticle.h.

Referenced by numberOfTrackerHits(), and setNumberOfTrackerHits().

int TrackingParticle::numberOfTrackerLayers_
private

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

Definition at line 258 of file TrackingParticle.h.

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

TrackingVertexRef TrackingParticle::parentVertex_
private

Definition at line 265 of file TrackingParticle.h.

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