CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingParticle.cc
Go to the documentation of this file.
3 
5 
7 
8 const unsigned int TrackingParticle::longLivedTag = 65536;
9 
11 {
12  // No operation
13 }
14 
15 TrackingParticle::TrackingParticle( const SimTrack& simtrk, const TrackingVertexRef& parentVertex )
16 {
17  addG4Track( simtrk );
18  setParentVertex( parentVertex );
19 }
20 
22 {
23 }
24 
26 {
27  if( genParticles_.empty() ) return g4Tracks_[0].type();
28  else return (*genParticles_.begin())->pdgId();
29 }
30 
32 {
33  return g4Tracks_[0].eventId();
34 }
35 
37 {
38  genParticles_.push_back( ref );
39 }
40 
42 {
43  g4Tracks_.push_back( t );
44 }
45 
47 {
48  return genParticles_.begin();
49 }
50 
52 {
53  return genParticles_.end();
54 }
55 
57 {
58  return g4Tracks_.begin();
59 }
60 
62 {
63  return g4Tracks_.end();
64 }
65 
67 {
68  parentVertex_=ref;
69 }
70 
72 {
74 }
75 
77 {
79 }
80 
82 {
84 }
85 
87 {
88  return genParticles_;
89 }
90 
91 const std::vector<SimTrack>& TrackingParticle::g4Tracks() const
92 {
93  return g4Tracks_;
94 }
95 
97 {
98  return parentVertex_;
99 }
100 
102 {
103  return decayVertices_;
104 }
105 
107 {
108  return decayVertices_.begin();
109 }
110 
112 {
113  return decayVertices_.end();
114 }
115 
116 
118 {
119  return g4Tracks_[0].momentum();
120 }
121 
123 {
124  return p4().Vect();
125 }
126 
128 {
129  return p4().BoostToCM();
130 }
131 
132 double TrackingParticle::p() const
133 {
134  return p4().P();
135 }
136 
138 {
139  return p4().E();
140 }
141 
142 double TrackingParticle::et() const
143 {
144  return p4().Et();
145 }
146 
148 {
149  return p4().M();
150 }
151 
153 {
154  return pow( mass(), 2 );
155 }
156 
157 double TrackingParticle::mt() const
158 {
159  return p4().Mt();
160 }
161 
163 {
164  return p4().Mt2();
165 }
166 
167 double TrackingParticle::px() const
168 {
169  return p4().Px();
170 }
171 
172 double TrackingParticle::py() const
173 {
174  return p4().Py();
175 }
176 
177 double TrackingParticle::pz() const
178 {
179  return p4().Pz();
180 }
181 
182 double TrackingParticle::pt() const
183 {
184  return p4().Pt();
185 }
186 
187 double TrackingParticle::phi() const
188 {
189  return p4().Phi();
190 }
191 
193 {
194  return p4().Theta();
195 }
196 
197 double TrackingParticle::eta() const
198 {
199  return p4().Eta();
200 }
201 
203 {
204  return p4().Rapidity();
205 }
206 
207 double TrackingParticle::y() const
208 {
209  return rapidity();
210 }
211 
212 double TrackingParticle::vx() const
213 {
214  const TrackingVertex& r=( *parentVertex_);
215  return r.position().X();
216 }
217 
218 double TrackingParticle::vy() const
219 {
220  const TrackingVertex& r=( *parentVertex_);
221  return r.position().Y();
222 }
223 
224 double TrackingParticle::vz() const
225 {
226  const TrackingVertex& r=( *parentVertex_);
227  return r.position().Z();
228 }
229 
230 
232 {
233  edm::LogWarning("TrackingParticle") << "The method matchedHit() has been deprecated. Use numberOfTrackerLayers() instead.";
234  return numberOfTrackerLayers_;
235 }
236 
237 
238 void TrackingParticle::setNumberOfHits( int numberOfHits )
239 {
241 }
242 
243 void TrackingParticle::setNumberOfTrackerHits( int numberOfTrackerHits )
244 {
246 }
247 
248 void TrackingParticle::setNumberOfTrackerLayers( const int numberOfTrackerLayers )
249 {
251 }
252 
253 std::ostream& operator<< (std::ostream& s, TrackingParticle const & tp)
254 {
255  s << "TP momentum, q, ID, & Event #: "
256  << tp.p4() << " " << tp.charge() << " " << tp.pdgId() << " "
257  << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
258 
259  for (TrackingParticle::genp_iterator hepT = tp.genParticle_begin(); hepT != tp.genParticle_end(); ++hepT)
260  {
261  s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
262  }
263 
264  for (TrackingParticle::g4t_iterator g4T = tp.g4Track_begin(); g4T != tp.g4Track_end(); ++g4T)
265  {
266  s << " Geant Track Momentum " << g4T->momentum() << std::endl;
267  s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
268  if (g4T->type() != tp.pdgId())
269  {
270  s << " Mismatch b/t TrackingParticle and Geant types" << std::endl;
271  }
272  }
273  // Loop over decay vertices
274  s << " TP Vertex " << tp.vertex() << std::endl;
275  s << " Source vertex: " << tp.parentVertex()->position() << std::endl;
276  s << " " << tp.decayVertices().size() << " Decay vertices" << std::endl;
277  for (tv_iterator iTV = tp.decayVertices_begin(); iTV != tp.decayVertices_end(); ++iTV)
278  {
279  s << " Decay vertices: " << (**iTV).position() << std::endl;
280  }
281 
282  return s;
283 }
type
Definition: HCALResponse.h:21
math::XYZTLorentzVectorD LorentzVector
Lorentz vector.
genp_iterator genParticle_begin() const
iterators
tv_iterator decayVertices_end() const
Electric charge. Note this is taken from the first SimTrack only.
int event() const
get the contents of the subdetector field (should be protected?)
tv_iterator decayVertices_begin() const
TrackingVertexRefVector decayVertices_
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.
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:86
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
void setParentVertex(const TrackingVertexRef &ref)
void addGenParticle(const reco::GenParticleRef &ref)
reco::GenParticleRefVector genParticles_
double pz() const
z coordinate of momentum vector. Note this is taken from the first SimTrack only. ...
float charge() const
Gives charge in unit of quark charge (should be 3 time the abaove)
const reco::GenParticleRefVector & genParticles() const
double mtSqr() const
Transverse mass squared. Note this is taken from the first SimTrack only.
double vy() const
y coordinate of parent vertex position
double y() const
Same as rapidity().
int bunchCrossing() const
get the detector field from this detid
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.
const TrackingVertexRefVector & decayVertices() const
std::vector< SimTrack >::const_iterator g4t_iterator
edm::Ref< TrackingVertexCollection > TrackingVertexRef
double mass() const
Mass. Note this is taken from the first SimTrack only.
genp_iterator genParticle_end() const
void clear()
Clear the vector.
Definition: RefVector.h:133
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
EncodedEventId eventId() const
Signal source, crossing number.
void setNumberOfTrackerHits(int numberOfTrackerHits)
const TrackingVertexRef & parentVertex() const
Vector momentum() const
spatial momentum vector
int numberOfTrackerLayers_
The number of tracker layers with hits. Equivalent to the old matchedHit.
static const unsigned int longLivedTag
long lived flag
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
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.
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
math::XYZVectorD Vector
point in the space
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.
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.
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
const LorentzVector & position() const