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 {
62  if( genParticles_.empty() ) return g4Tracks_[0].type();
63  else return (*genParticles_.begin())->pdgId();
64  }
65 
71  return g4Tracks_[0].eventId();
72  }
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 
97 
99  float charge() const { return g4Tracks_[0].charge(); }
101  int threeCharge() const { return lrintf(3.f*charge()); }
102 
104  const LorentzVector& p4() const {
105  return g4Tracks_[0].momentum();
106  }
107 
109  Vector momentum() const {
110  return p4().Vect();
111  }
112 
114  Vector boostToCM() const {
115  return p4().BoostToCM();
116  }
117 
119  double p() const {
120  return p4().P();
121  }
122 
124  double energy() const {
125  return p4().E();
126  }
127 
129  double et() const {
130  return p4().Et();
131  }
132 
134  double mass() const {
135  return p4().M();
136  }
137 
139  double massSqr() const {
140  return pow( mass(), 2 );
141  }
142 
144  double mt() const {
145  return p4().Mt();
146  }
147 
149  double mtSqr() const {
150  return p4().Mt2();
151  }
152 
154  double px() const {
155  return p4().Px();
156  }
157 
159  double py() const {
160  return p4().Py();
161  }
162 
164  double pz() const {
165  return p4().Pz();
166  }
167 
169  double pt() const {
170  return p4().Pt();
171  }
172 
174  double phi() const {
175  return p4().Phi();
176  }
177 
179  double theta() const {
180  return p4().Theta();
181  }
182 
184  double eta() const {
185  return p4().Eta();
186  }
187 
189  double rapidity() const {
190  return p4().Rapidity();
191  }
192 
194  double y() const {
195  return rapidity();
196  }
197 
199  Point vertex() const {
200  const TrackingVertex::LorentzVector & p = (*parentVertex_).position();
201  return Point(p.x(),p.y(),p.z());
202  }
203 
205  double vx() const {
206  const TrackingVertex& r=( *parentVertex_);
207  return r.position().X();
208  }
209 
211  double vy() const {
212  const TrackingVertex& r=( *parentVertex_);
213  return r.position().Y();
214  }
215  // @brief z coordinate of parent vertex position
216  double vz() const {
217  const TrackingVertex& r=( *parentVertex_);
218  return r.position().Z();
219  }
220 
224  int status() const {
225  return genParticles_.empty() ? -99 : (*genParticles_[0]).status();
226  }
227 
228  static const unsigned int longLivedTag;
229 
231  bool longLived() const { return status()&longLivedTag;}
232 
236  int numberOfHits() const {return numberOfHits_;}
237 
242 
245  int matchedHit() const;
246 
251 
252  void setNumberOfHits( int numberOfHits );
255 private:
259 
261  std::vector<SimTrack> g4Tracks_;
263 
264  // Source and decay vertices
267 };
268 
269 #endif // SimDataFormats_TrackingParticle_H
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
type
Definition: HCALResponse.h:21
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
genp_iterator genParticle_begin() const
iterators
tv_iterator decayVertices_end() const
tv_iterator decayVertices_begin() const
const std::vector< SimTrack > & g4Tracks() const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
TrackingVertexRefVector decayVertices_
Vector boostToCM() const
Vector to boost to the particle centre of mass frame.
Vector momentum() const
spatial momentum vector
const reco::GenParticleRefVector & genParticles() const
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.
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.
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:96
int status() const
Status word.
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
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.
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 TrackingVertexRef & parentVertex() 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)
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.
void setNumberOfTrackerHits(int numberOfTrackerHits)
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
const TrackingVertexRefVector & decayVertices() 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.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
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
const LorentzVector & position() const
TrackingParticle::g4t_iterator g4t_iterator