CMS 3D CMS Logo

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  friend std::ostream& operator<<(std::ostream& s, TrackingParticle const& tp);
31 
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 {
62  if (genParticles_.empty())
63  return g4Tracks_[0].type();
64  else
65  return (*genParticles_.begin())->pdgId();
66  }
67 
72  EncodedEventId eventId() const { return g4Tracks_[0].eventId(); }
73 
74  // Setters for G4 and reco::GenParticle
75  void addGenParticle(const reco::GenParticleRef& ref);
76  void addG4Track(const SimTrack& t);
81  g4t_iterator g4Track_end() const;
82  void setParentVertex(const TrackingVertexRef& ref);
83  void addDecayVertex(const TrackingVertexRef& ref);
84  void clearParentVertex();
85  void clearDecayVertices();
86 
87  // Getters for Embd and Sim Tracks
89  const std::vector<SimTrack>& g4Tracks() const { return g4Tracks_; }
90  const TrackingVertexRef& parentVertex() const { return parentVertex_; }
91 
92  // Accessors for vector of decay vertices
96 
98  float charge() const { return g4Tracks_[0].charge(); }
100  int threeCharge() const { return lrintf(3.f * charge()); }
101 
103  const LorentzVector& p4() const { return g4Tracks_[0].momentum(); }
104 
106  Vector momentum() const { return p4().Vect(); }
107 
109  Vector boostToCM() const { return p4().BoostToCM(); }
110 
112  double p() const { return p4().P(); }
113 
115  double qoverp() const { return charge() / p(); }
116 
118  double energy() const { return p4().E(); }
119 
121  double et() const { return p4().Et(); }
122 
124  double mass() const { return p4().M(); }
125 
127  double massSqr() const { return pow(mass(), 2); }
128 
130  double mt() const { return p4().Mt(); }
131 
133  double mtSqr() const { return p4().Mt2(); }
134 
136  double px() const { return p4().Px(); }
137 
139  double py() const { return p4().Py(); }
140 
142  double pz() const { return p4().Pz(); }
143 
145  double pt() const { return p4().Pt(); }
146 
148  double phi() const { return p4().Phi(); }
149 
151  double theta() const { return p4().Theta(); }
152 
154  double eta() const { return p4().Eta(); }
155 
157  double lambda() const { return M_PI_2 - theta(); }
158 
160  double tanl() const { return tan(lambda()); }
161 
163  double rapidity() const { return p4().Rapidity(); }
164 
166  double y() const { return rapidity(); }
167 
169  Point vertex() const {
170  const TrackingVertex::LorentzVector& p = (*parentVertex_).position();
171  return Point(p.x(), p.y(), p.z());
172  }
173 
175  double vx() const {
176  const TrackingVertex& r = (*parentVertex_);
177  return r.position().X();
178  }
179 
181  double vy() const {
182  const TrackingVertex& r = (*parentVertex_);
183  return r.position().Y();
184  }
185 
187  double vz() const {
188  const TrackingVertex& r = (*parentVertex_);
189  return r.position().Z();
190  }
191 
193  double dxy() const { return (-vx() * py() + vy() * px()) / pt(); }
194 
196  double d0() const { return -dxy(); }
197 
199  double dz() const { return vz() - (vx() * px() + vy() * py()) * pz() / p4().Perp2(); }
200 
202  double z0() const { return dz(); }
203 
207  int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
208 
209  static const unsigned int longLivedTag;
210 
212  bool longLived() const { return status() & longLivedTag; }
213 
217  int numberOfHits() const { return numberOfHits_; }
218 
223 
226  int matchedHit() const;
227 
232 
233  void setNumberOfHits(int numberOfHits);
236 
237 private:
241 
243  std::vector<SimTrack> g4Tracks_;
245 
246  // Source and decay vertices
249 };
250 
251 #endif // SimDataFormats_TrackingParticle_H
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
double dxy() const
dxy parameter.
const reco::GenParticleRefVector & genParticles() const
int pdgId() const
PDG ID.
EncodedEventId eventId() const
Signal source, crossing number.
Vector momentum() const
spatial momentum vector
int status() const
Status word.
double mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
double vx() const
x coordinate of parent vertex position
TrackingVertexRefVector decayVertices_
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
tv_iterator decayVertices_end() const
TrackingParticle()
Default constructor. Note that the object will be useless until it is provided with a SimTrack and pa...
double et() const
Transverse energy. Note this is taken from the first SimTrack only.
#define M_PI_2
double rapidity() const
Rapidity. Note this is taken from the first SimTrack only.
int numberOfTrackerHits_
The number of tracker only hits.
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
double y() const
Same as rapidity().
bool longLived() const
is long lived?
const std::vector< SimTrack > & g4Tracks() const
g4t_iterator g4Track_end() const
void setParentVertex(const TrackingVertexRef &ref)
math::XYZPointD Point
point in the space
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
double energy() const
Energy. Note this is taken from the first SimTrack only.
void addGenParticle(const reco::GenParticleRef &ref)
math::XYZTLorentzVectorD LorentzVector
double lambda() const
Lambda angle. Note this is taken from the first SimTrack only.
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
reco::GenParticleRefVector genParticles_
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
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...
genp_iterator genParticle_begin() const
iterators
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
double mass() const
Mass. Note this is taken from the first SimTrack only.
void addG4Track(const SimTrack &t)
double mt() const
Transverse mass. Note this is taken from the first SimTrack only.
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
int numberOfHits_
The total number of hits.
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
double f[11][100]
int matchedHit() const
double massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
std::vector< SimTrack >::const_iterator g4t_iterator
double eta() const
Momentum pseudorapidity. Note this is taken from the first SimTrack only.
g4t_iterator g4Track_begin() const
const TrackingVertexRef & parentVertex() 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.
double tanl() const
tangent of the lambda angle. Note this is taken from the first SimTrack only.
void addDecayVertex(const TrackingVertexRef &ref)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
double vy() const
y coordinate of parent vertex position
const TrackingVertexRefVector & decayVertices() const
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...
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. ...
void setNumberOfTrackerHits(int numberOfTrackerHits)
double phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
Structure Point Contains parameters of Gaussian fits to DMRs.
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
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
genp_iterator genParticle_end() const
Monte Carlo truth information used for tracking validation.
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
math::XYZVectorD Vector
point in the space
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
double qoverp() const
Quotient of the electric charge over the magnitude of the momentum vector. Note this is taken from th...
Vector boostToCM() const
Vector to boost to the particle centre of mass frame.
TrackingVertexRef parentVertex_
tv_iterator decayVertices_begin() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
friend std::ostream & operator<<(std::ostream &s, TrackingParticle const &tp)
Point vertex() const
Parent vertex position.
void setNumberOfHits(int numberOfHits)
double z0() const
z0 parameter
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
TrackingParticle::g4t_iterator g4t_iterator
double vz() const
z coordinate of parent vertex position