CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingParticle.h
Go to the documentation of this file.
1 #ifndef SimDataFormats_TrackingParticle_h
2 #define SimDataFormats_TrackingParticle_h
3 
4 #include <vector>
11 
12 //
13 // Forward declarations
14 //
15 class TrackingVertex;
16 class SimTrack;
17 class EncodedEventId;
18 
30 {
31  friend std::ostream& operator<< (std::ostream& s, TrackingParticle const& tp);
32 public:
33  typedef int Charge;
38 
41  typedef std::vector<SimTrack>::const_iterator g4t_iterator;
42 
51 
53 
54  // destructor
56 
61  int pdgId() const;
66  EncodedEventId eventId() const;
67 
68  // Setters for G4 and reco::GenParticle
69  void addGenParticle( const reco::GenParticleRef& ref);
70  void addG4Track( const SimTrack& t);
75  g4t_iterator g4Track_end() const;
76  void setParentVertex(const TrackingVertexRef& ref);
77  void addDecayVertex(const TrackingVertexRef& ref);
78  void clearParentVertex();
79  void clearDecayVertices();
80  // Getters for Embd and Sim Tracks
82  const std::vector<SimTrack>& g4Tracks() const;
83  const TrackingVertexRef& parentVertex() const;
84 
85  // Accessors for vector of decay vertices
89 
90 
92  float charge() const { return g4Tracks_[0].charge(); }
94  int threeCharge() const { return lrintf(3.f*charge()); }
95 
96  const LorentzVector& p4() const;
97 
98  Vector momentum() const;
99 
100  Vector boostToCM() const;
101 
102  double p() const;
103  double energy() const;
104  double et() const;
105  double mass() const;
106  double massSqr() const;
107  double mt() const;
108  double mtSqr() const;
109  double px() const;
110  double py() const;
111  double pz() const;
112  double pt() const;
113  double phi() const;
114  double theta() const;
115  double eta() const;
116  double rapidity() const;
117  double y() const;
118 
120  Point vertex() const {
121  const TrackingVertex::LorentzVector & p = (*parentVertex_).position();
122  return Point(p.x(),p.y(),p.z());
123  }
124 
125  double vx() const;
126  double vy() const;
127  double vz() const;
128 
131  int status() const {
132  return genParticles_.empty() ? -99 : (*genParticles_[0]).status();
133  }
134 
135  static const unsigned int longLivedTag;
136 
138  bool longLived() const { return status()&longLivedTag;}
139 
143  int numberOfHits() const {return numberOfHits_;}
144 
149 
152  int matchedHit() const;
157 
158  void setNumberOfHits( int numberOfHits );
161 private:
165 
167  std::vector<SimTrack> g4Tracks_;
169 
170  // Source and decay vertices
173 };
174 
175 #endif // SimDataFormats_TrackingParticle_H
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
genp_iterator genParticle_begin() const
iterators
tv_iterator decayVertices_end() const
tv_iterator decayVertices_begin() const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
TrackingVertexRefVector decayVertices_
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
g4t_iterator g4Track_begin() const
int pdgId() const
PDG ID.
TrackingParticle()
Default constructor. Note that the object will be useless until it is provided with a SimTrack and pa...
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
Vector boostToCM() const
Vector to boost to the particle centre of mass frame.
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
int matchedHit() const
int numberOfTrackerHits_
The number of tracker only hits.
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:86
int status() const
Status word.
void setParentVertex(const TrackingVertexRef &ref)
math::XYZPointD Point
point in the space
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void addGenParticle(const reco::GenParticleRef &ref)
math::XYZTLorentzVectorD LorentzVector
reco::GenParticleRefVector genParticles_
double pz() const
z coordinate 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.
const reco::GenParticleRefVector & genParticles() const
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
double mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
bool longLived() const
is long lived?
double vy() const
y coordinate of parent vertex position
double y() const
Same as rapidity().
double massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
void addG4Track(const SimTrack &t)
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
int numberOfHits_
The total number of hits.
double mt() const
Transverse mass. Note this is taken from the first SimTrack only.
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
double f[11][100]
const TrackingVertexRefVector & decayVertices() const
std::vector< SimTrack >::const_iterator g4t_iterator
double mass() const
Mass. Note this is taken from the first SimTrack only.
genp_iterator genParticle_end() const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
void addDecayVertex(const TrackingVertexRef &ref)
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
const std::vector< SimTrack > & g4Tracks() const
double vx() const
x coordinate of parent vertex position
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...
double eta() const
Momentum pseudorapidity. Note this is taken from the first SimTrack only.
double rapidity() const
Rapidity. Note this is taken from the first SimTrack only.
Point vertex() const
Parent vertex position.
EncodedEventId eventId() const
Signal source, crossing number.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > LorentzVector
Definition: analysisEnums.h:9
void setNumberOfTrackerHits(int numberOfTrackerHits)
const TrackingVertexRef & parentVertex() const
Vector momentum() const
spatial momentum vector
int Charge
electric charge type
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
static const unsigned int longLivedTag
long lived flag
Monte Carlo truth information used for tracking validation.
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
math::XYZVectorD Vector
point in the space
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times &quot;charge()&quot;)
g4t_iterator g4Track_end() const
double et() const
Transverse energy. Note this is taken from the first SimTrack only.
TrackingVertexRef parentVertex_
double energy() const
Energy. Note this is taken from the first SimTrack only.
double phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
friend std::ostream & operator<<(std::ostream &s, TrackingParticle const &tp)
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
void setNumberOfHits(int numberOfHits)
double vz() const
z coordinate of parent vertex position