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 energy() const { return p4().E(); }
116 
118  double et() const { return p4().Et(); }
119 
121  double mass() const { return p4().M(); }
122 
124  double massSqr() const { return pow(mass(), 2); }
125 
127  double mt() const { return p4().Mt(); }
128 
130  double mtSqr() const { return p4().Mt2(); }
131 
133  double px() const { return p4().Px(); }
134 
136  double py() const { return p4().Py(); }
137 
139  double pz() const { return p4().Pz(); }
140 
142  double pt() const { return p4().Pt(); }
143 
145  double phi() const { return p4().Phi(); }
146 
148  double theta() const { return p4().Theta(); }
149 
151  double eta() const { return p4().Eta(); }
152 
154  double rapidity() const { return p4().Rapidity(); }
155 
157  double y() const { return rapidity(); }
158 
160  Point vertex() const {
161  const TrackingVertex::LorentzVector& p = (*parentVertex_).position();
162  return Point(p.x(), p.y(), p.z());
163  }
164 
166  double vx() const {
167  const TrackingVertex& r = (*parentVertex_);
168  return r.position().X();
169  }
170 
172  double vy() const {
173  const TrackingVertex& r = (*parentVertex_);
174  return r.position().Y();
175  }
176  // @brief z coordinate of parent vertex position
177  double vz() const {
178  const TrackingVertex& r = (*parentVertex_);
179  return r.position().Z();
180  }
181 
185  int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
186 
187  static const unsigned int longLivedTag;
188 
190  bool longLived() const { return status() & longLivedTag; }
191 
195  int numberOfHits() const { return numberOfHits_; }
196 
201 
204  int matchedHit() const;
205 
210 
211  void setNumberOfHits(int numberOfHits);
214 
215 private:
219 
221  std::vector<SimTrack> g4Tracks_;
223 
224  // Source and decay vertices
227 };
228 
229 #endif // SimDataFormats_TrackingParticle_H
TrackingParticle::numberOfTrackerHits
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
Definition: TrackingParticle.h:200
TrackingParticle::numberOfTrackerLayers
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
Definition: TrackingParticle.h:209
TrackingParticle::addGenParticle
void addGenParticle(const reco::GenParticleRef &ref)
Definition: TrackingParticle.cc:21
TrackingParticle::eta
double eta() const
Momentum pseudorapidity. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:151
TrackingParticle::g4Track_begin
g4t_iterator g4Track_begin() const
Definition: TrackingParticle.cc:29
TrackingParticle::mass
double mass() const
Mass. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:121
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
TrackingParticle::mt
double mt() const
Transverse mass. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:127
TrackingParticle::clearDecayVertices
void clearDecayVertices()
Definition: TrackingParticle.cc:39
g4t_iterator
TrackingParticle::g4t_iterator g4t_iterator
Definition: TrackingTruthValid.cc:24
TrackingParticle::numberOfHits_
int numberOfHits_
The total number of hits.
Definition: TrackingParticle.h:216
math::XYZTLorentzVectorD
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
TrackingParticle::addDecayVertex
void addDecayVertex(const TrackingVertexRef &ref)
Definition: TrackingParticle.cc:35
TrackingParticle::g4Tracks_
std::vector< SimTrack > g4Tracks_
references to G4 and reco::GenParticle tracks
Definition: TrackingParticle.h:221
edm::RefVector::begin
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
TrackingParticle::g4t_iterator
std::vector< SimTrack >::const_iterator g4t_iterator
Definition: TrackingParticle.h:41
TrackingParticle::vy
double vy() const
y coordinate of parent vertex position
Definition: TrackingParticle.h:172
TrackingParticle::p4
const LorentzVector & p4() const
Four-momentum Lorentz vector. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:103
TrackingParticle::decayVertices
const TrackingVertexRefVector & decayVertices() const
Definition: TrackingParticle.h:93
TrackingParticle::py
double py() const
y coordinate of momentum vector. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:136
TrackingParticle::matchedHit
int matchedHit() const
Definition: TrackingParticle.cc:41
TrackingParticle::parentVertex
const TrackingVertexRef & parentVertex() const
Definition: TrackingParticle.h:90
TrackingParticle::Vector
math::XYZVectorD Vector
point in the space
Definition: TrackingParticle.h:37
edm::RefVector< GenParticleCollection >
TrackingParticle::massSqr
double massSqr() const
Mass squared. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:124
TrackingParticle::genp_iterator
reco::GenParticleRefVector::iterator genp_iterator
reference to reco::GenParticle
Definition: TrackingParticle.h:40
TrackingParticle::setNumberOfTrackerLayers
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
Definition: TrackingParticle.cc:51
TrackingParticle::boostToCM
Vector boostToCM() const
Vector to boost to the particle centre of mass frame.
Definition: TrackingParticle.h:109
TrackingParticle::numberOfTrackerLayers_
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
Definition: TrackingParticle.h:218
TrackingParticle::longLivedTag
static const unsigned int longLivedTag
long lived flag
Definition: TrackingParticle.h:187
edm::Ref< TrackingVertexCollection >
TrackingParticle::parentVertex_
TrackingVertexRef parentVertex_
Definition: TrackingParticle.h:225
GenParticle.h
EncodedEventId
Definition: EncodedEventId.h:11
TrackingParticle::charge
float charge() const
Electric charge. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:98
TrackingParticle::px
double px() const
x coordinate of momentum vector. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:133
TrackingParticle::PolarLorentzVector
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: TrackingParticle.h:35
alignCSCRings.s
s
Definition: alignCSCRings.py:92
edm::RefVector::end
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
TrackingParticle::pt
double pt() const
Transverse momentum. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:142
TrackingParticle::numberOfHits
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...
Definition: TrackingParticle.h:195
TrackingVertexContainer.h
edm::RefVector::empty
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
TrackingParticle::genParticle_end
genp_iterator genParticle_end() const
Definition: TrackingParticle.cc:27
TrackingParticle::clearParentVertex
void clearParentVertex()
Definition: TrackingParticle.cc:37
TrackingParticle::setNumberOfTrackerHits
void setNumberOfTrackerHits(int numberOfTrackerHits)
Definition: TrackingParticle.cc:49
TrackingParticle
Monte Carlo truth information used for tracking validation.
Definition: TrackingParticle.h:29
TrackingParticle::phi
double phi() const
Momentum azimuthal angle. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:145
TrackingParticle::theta
double theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:148
TrackingParticle::operator<<
friend std::ostream & operator<<(std::ostream &s, TrackingParticle const &tp)
Definition: TrackingParticle.cc:55
TrackingParticle::momentum
Vector momentum() const
spatial momentum vector
Definition: TrackingParticle.h:106
TrackingParticle::vertex
Point vertex() const
Parent vertex position.
Definition: TrackingParticle.h:160
TrackingParticle::mtSqr
double mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:130
TrackingParticle::Charge
int Charge
electric charge type
Definition: TrackingParticle.h:33
cmsswSequenceInfo.tp
tp
Definition: cmsswSequenceInfo.py:17
math::XYZVectorD
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
TrackingParticle::genParticles
const reco::GenParticleRefVector & genParticles() const
Definition: TrackingParticle.h:88
Point
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
TrackingParticle::status
int status() const
Status word.
Definition: TrackingParticle.h:185
TrackingParticle::decayVertices_
TrackingVertexRefVector decayVertices_
Definition: TrackingParticle.h:226
TrackingParticle::vz
double vz() const
Definition: TrackingParticle.h:177
math::XYZPointD
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:8
TrackingParticle::decayVertices_begin
tv_iterator decayVertices_begin() const
Definition: TrackingParticle.h:94
TrackingParticle::~TrackingParticle
~TrackingParticle()
Definition: TrackingParticle.cc:19
TrackingVertex
Definition: TrackingVertex.h:22
gainCalibHelper::gainCalibPI::type
type
Definition: SiPixelGainCalibHelper.h:39
LorentzVector.h
TrackingParticle::decayVertices_end
tv_iterator decayVertices_end() const
Definition: TrackingParticle.h:95
TrackingParticle::pdgId
int pdgId() const
PDG ID.
Definition: TrackingParticle.h:61
TrackingParticle::setNumberOfHits
void setNumberOfHits(int numberOfHits)
Definition: TrackingParticle.cc:47
TrackingParticle::TrackingParticle
TrackingParticle()
Default constructor. Note that the object will be useless until it is provided with a SimTrack and pa...
Definition: TrackingParticle.cc:10
TrackingParticle::pz
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:139
TrackingParticle::rapidity
double rapidity() const
Rapidity. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:154
TrackingParticle::et
double et() const
Transverse energy. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:118
TrackingParticle::LorentzVector
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
Definition: TrackingParticle.h:34
TrackingParticle::eventId
EncodedEventId eventId() const
Signal source, crossing number.
Definition: TrackingParticle.h:72
alignCSCRings.r
r
Definition: alignCSCRings.py:93
TrackingParticle::setParentVertex
void setParentVertex(const TrackingVertexRef &ref)
Definition: TrackingParticle.cc:33
TrackingParticle::vx
double vx() const
x coordinate of parent vertex position
Definition: TrackingParticle.h:166
TrackingParticle::genParticles_
reco::GenParticleRefVector genParticles_
Definition: TrackingParticle.h:222
TrackingParticle::longLived
bool longLived() const
is long lived?
Definition: TrackingParticle.h:190
TrackingParticleFwd.h
SimTrack
Definition: SimTrack.h:6
TrackingParticle::p
double p() const
Magnitude of momentum vector. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:112
TrackingParticle::threeCharge
int threeCharge() const
Gives charge in unit of quark charge (should be 3 times "charge()")
Definition: TrackingParticle.h:100
TrackingParticle::g4Track_end
g4t_iterator g4Track_end() const
Definition: TrackingParticle.cc:31
edm::RefVectorIterator
Definition: EDProductfwd.h:33
Point3D.h
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
math::PtEtaPhiMLorentzVector
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
TrackingParticle::Point
math::XYZPointD Point
point in the space
Definition: TrackingParticle.h:36
Vector3D.h
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
TrackingParticle::numberOfTrackerHits_
int numberOfTrackerHits_
The number of tracker only hits.
Definition: TrackingParticle.h:217
TrackingParticle::y
double y() const
Same as rapidity().
Definition: TrackingParticle.h:157
TrackingParticle::energy
double energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: TrackingParticle.h:115
TrackingParticle::addG4Track
void addG4Track(const SimTrack &t)
Definition: TrackingParticle.cc:23
TrackingVertex::LorentzVector
math::XYZTLorentzVectorD LorentzVector
Definition: TrackingVertex.h:28
TrackingParticle::g4Tracks
const std::vector< SimTrack > & g4Tracks() const
Definition: TrackingParticle.h:89
TrackingParticle::genParticle_begin
genp_iterator genParticle_begin() const
iterators
Definition: TrackingParticle.cc:25