1 |
grimes |
1.1 |
#ifndef SimDataFormats_TrackingParticle_h
|
2 |
|
|
#define SimDataFormats_TrackingParticle_h
|
3 |
|
|
|
4 |
|
|
/** Concrete TrackingParticle.
|
5 |
|
|
* All track parameters are passed in the constructor and stored internally.
|
6 |
|
|
*/
|
7 |
|
|
|
8 |
|
|
#include <map>
|
9 |
|
|
|
10 |
grimes |
1.1.2.1 |
#include "DataFormats/Math/interface/Point3D.h"
|
11 |
|
|
#include "DataFormats/Math/interface/Vector3D.h"
|
12 |
|
|
#include "DataFormats/Math/interface/LorentzVector.h"
|
13 |
|
|
#include "Rtypes.h"
|
14 |
|
|
|
15 |
grimes |
1.1 |
#include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
|
16 |
grimes |
1.1.2.1 |
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
|
17 |
|
|
#include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
|
18 |
|
|
#include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
|
19 |
grimes |
1.1 |
#include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
|
20 |
|
|
#include "SimDataFormats/Track/interface/SimTrackContainer.h"
|
21 |
|
|
#include "DataFormats/DetId/interface/DetId.h"
|
22 |
|
|
|
23 |
|
|
class TrackingVertex;
|
24 |
|
|
|
25 |
grimes |
1.1.2.1 |
class TrackingParticle
|
26 |
grimes |
1.1 |
{
|
27 |
grimes |
1.1.2.1 |
friend std::ostream& operator<< (std::ostream& s, TrackingParticle const& tp);
|
28 |
grimes |
1.1 |
|
29 |
|
|
public:
|
30 |
|
|
|
31 |
grimes |
1.1.2.1 |
/// electric charge type
|
32 |
|
|
typedef int Charge;
|
33 |
|
|
/// Lorentz vector
|
34 |
|
|
typedef math::XYZTLorentzVectorD LorentzVector;
|
35 |
|
|
/// Lorentz vector
|
36 |
|
|
typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
|
37 |
|
|
/// point in the space
|
38 |
|
|
typedef math::XYZPointD Point;
|
39 |
|
|
/// point in the space
|
40 |
|
|
typedef math::XYZVectorD Vector;
|
41 |
|
|
|
42 |
|
|
/// reference to reco::GenParticle
|
43 |
|
|
typedef reco::GenParticleRefVector::iterator genp_iterator;
|
44 |
|
|
typedef std::vector<SimTrack>::const_iterator g4t_iterator;
|
45 |
grimes |
1.1 |
|
46 |
|
|
/// default constructor
|
47 |
|
|
TrackingParticle() {}
|
48 |
grimes |
1.1.2.1 |
|
49 |
|
|
/// constructor from SimTrack and SimVertex
|
50 |
|
|
TrackingParticle(const SimTrack& simtrk, const TrackingVertexRef& simvtx) {
|
51 |
|
|
addG4Track(simtrk);
|
52 |
|
|
setParentVertex(simvtx);
|
53 |
|
|
}
|
54 |
grimes |
1.1 |
|
55 |
|
|
// destructor
|
56 |
grimes |
1.1.2.1 |
~TrackingParticle() {}
|
57 |
grimes |
1.1 |
|
58 |
|
|
/// PDG id, signal source, crossing number
|
59 |
|
|
int pdgId() const
|
60 |
|
|
{
|
61 |
sarkar |
1.1.2.2 |
return g4Tracks_.at(0).type();
|
62 |
grimes |
1.1 |
}
|
63 |
|
|
EncodedEventId eventId() const
|
64 |
|
|
{
|
65 |
grimes |
1.1.2.1 |
return g4Tracks_.at(0).eventId();
|
66 |
grimes |
1.1 |
}
|
67 |
|
|
|
68 |
grimes |
1.1.2.1 |
// Setters for G4 and reco::GenParticle
|
69 |
|
|
void addGenParticle( const reco::GenParticleRef& ref) {
|
70 |
|
|
genParticles_.push_back(ref);
|
71 |
|
|
}
|
72 |
|
|
void addG4Track( const SimTrack& t) {
|
73 |
|
|
g4Tracks_.push_back(t);
|
74 |
|
|
}
|
75 |
|
|
/// iterators
|
76 |
|
|
genp_iterator genParticle_begin() const {
|
77 |
|
|
return genParticles_.begin();
|
78 |
|
|
}
|
79 |
|
|
genp_iterator genParticle_end() const {
|
80 |
|
|
return genParticles_.end();
|
81 |
|
|
}
|
82 |
|
|
g4t_iterator g4Track_begin() const {
|
83 |
|
|
return g4Tracks_.begin();
|
84 |
|
|
}
|
85 |
|
|
g4t_iterator g4Track_end() const {
|
86 |
|
|
return g4Tracks_.end();
|
87 |
|
|
}
|
88 |
|
|
void setParentVertex(const TrackingVertexRef& ref) {
|
89 |
|
|
parentVertex_ = ref;
|
90 |
|
|
}
|
91 |
|
|
void addDecayVertex(const TrackingVertexRef& ref) {
|
92 |
|
|
decayVertices_.push_back(ref);
|
93 |
|
|
}
|
94 |
|
|
void clearParentVertex() {
|
95 |
|
|
parentVertex_ = TrackingVertexRef();
|
96 |
|
|
}
|
97 |
|
|
void clearDecayVertices() {
|
98 |
|
|
decayVertices_.clear();
|
99 |
|
|
}
|
100 |
|
|
void setMatchedHit(const int& hitnumb) {
|
101 |
|
|
matchedHit_ = hitnumb;
|
102 |
|
|
}
|
103 |
grimes |
1.1 |
// Getters for Embd and Sim Tracks
|
104 |
grimes |
1.1.2.1 |
const reco::GenParticleRefVector& genParticles() const {
|
105 |
grimes |
1.1 |
return genParticles_;
|
106 |
|
|
}
|
107 |
grimes |
1.1.2.1 |
const std::vector<SimTrack>& g4Tracks() const {
|
108 |
|
|
return g4Tracks_;
|
109 |
grimes |
1.1 |
}
|
110 |
grimes |
1.1.2.1 |
const TrackingVertexRef& parentVertex() const {
|
111 |
grimes |
1.1 |
return parentVertex_;
|
112 |
|
|
}
|
113 |
|
|
|
114 |
|
|
// Accessors for vector of decay vertices
|
115 |
grimes |
1.1.2.1 |
const TrackingVertexRefVector& decayVertices() const {
|
116 |
grimes |
1.1 |
return decayVertices_;
|
117 |
|
|
}
|
118 |
grimes |
1.1.2.1 |
tv_iterator decayVertices_begin() const {
|
119 |
grimes |
1.1 |
return decayVertices_.begin();
|
120 |
|
|
}
|
121 |
grimes |
1.1.2.1 |
tv_iterator decayVertices_end() const {
|
122 |
grimes |
1.1 |
return decayVertices_.end();
|
123 |
|
|
}
|
124 |
grimes |
1.1.2.1 |
int matchedHit() const {
|
125 |
grimes |
1.1 |
return matchedHit_;
|
126 |
|
|
}
|
127 |
grimes |
1.1.2.1 |
/// electric charge
|
128 |
|
|
int charge() const {
|
129 |
|
|
return g4Tracks_.at(0).charge() / 3;
|
130 |
|
|
}
|
131 |
|
|
/// electric charge
|
132 |
|
|
int threeCharge() const {
|
133 |
|
|
return g4Tracks_.at(0).charge();
|
134 |
|
|
}
|
135 |
|
|
/// four-momentum Lorentz vector
|
136 |
|
|
const LorentzVector& p4() const {
|
137 |
|
|
return g4Tracks_.at(0).momentum();
|
138 |
|
|
}
|
139 |
|
|
#if 0
|
140 |
|
|
/// four-momentum Lorentz vector
|
141 |
|
|
const PolarLorentzVector& polarP4() const
|
142 |
|
|
{
|
143 |
|
|
}
|
144 |
|
|
#endif
|
145 |
|
|
/// spatial momentum vector
|
146 |
|
|
Vector momentum() const {
|
147 |
|
|
return p4().Vect();
|
148 |
|
|
}
|
149 |
|
|
/// boost vector to boost a Lorentz vector
|
150 |
|
|
/// to the particle center of mass system
|
151 |
|
|
Vector boostToCM() const {
|
152 |
|
|
return p4().BoostToCM();
|
153 |
|
|
}
|
154 |
|
|
/// magnitude of momentum vector
|
155 |
|
|
double p() const {
|
156 |
|
|
return p4().P();
|
157 |
|
|
}
|
158 |
|
|
/// energy
|
159 |
|
|
double energy() const {
|
160 |
|
|
return p4().E();
|
161 |
|
|
}
|
162 |
|
|
/// transverse energy
|
163 |
|
|
double et() const {
|
164 |
|
|
return p4().Et();
|
165 |
|
|
}
|
166 |
|
|
/// mass
|
167 |
|
|
double mass() const {
|
168 |
|
|
return p4().M();
|
169 |
|
|
}
|
170 |
|
|
/// mass squared
|
171 |
|
|
double massSqr() const {
|
172 |
|
|
return pow(mass(), 2);
|
173 |
|
|
}
|
174 |
|
|
/// transverse mass
|
175 |
|
|
double mt() const {
|
176 |
|
|
return p4().Mt();
|
177 |
|
|
}
|
178 |
|
|
/// transverse mass squared
|
179 |
|
|
double mtSqr() const {
|
180 |
|
|
return p4().Mt2();
|
181 |
|
|
}
|
182 |
|
|
/// x coordinate of momentum vector
|
183 |
|
|
double px() const {
|
184 |
|
|
return p4().Px();
|
185 |
|
|
}
|
186 |
|
|
/// y coordinate of momentum vector
|
187 |
|
|
double py() const {
|
188 |
|
|
return p4().Py();
|
189 |
|
|
}
|
190 |
|
|
/// z coordinate of momentum vector
|
191 |
|
|
double pz() const {
|
192 |
|
|
return p4().Pz();
|
193 |
|
|
}
|
194 |
|
|
/// transverse momentum
|
195 |
|
|
double pt() const {
|
196 |
|
|
return p4().Pt();
|
197 |
|
|
}
|
198 |
|
|
/// momentum azimuthal angle
|
199 |
|
|
double phi() const {
|
200 |
|
|
return p4().Phi();
|
201 |
|
|
}
|
202 |
|
|
/// momentum polar angle
|
203 |
|
|
double theta() const {
|
204 |
|
|
return p4().Theta();
|
205 |
|
|
}
|
206 |
|
|
/// momentum pseudorapidity
|
207 |
|
|
double eta() const {
|
208 |
|
|
return p4().Eta();
|
209 |
|
|
}
|
210 |
|
|
/// rapidity
|
211 |
|
|
double rapidity() const {
|
212 |
|
|
return p4().Rapidity();
|
213 |
|
|
}
|
214 |
|
|
/// rapidity
|
215 |
|
|
double y() const {
|
216 |
|
|
return rapidity();
|
217 |
|
|
}
|
218 |
|
|
/// vertex position
|
219 |
|
|
Point vertex() const {
|
220 |
|
|
return Point(vx(), vy(), vz());
|
221 |
|
|
}
|
222 |
|
|
/// x coordinate of vertex position
|
223 |
|
|
double vx() const {
|
224 |
|
|
const TrackingVertex& r = (*parentVertex_);
|
225 |
|
|
return r.position().X();
|
226 |
|
|
}
|
227 |
|
|
/// y coordinate of vertex position
|
228 |
|
|
double vy() const {
|
229 |
|
|
const TrackingVertex& r = (*parentVertex_);
|
230 |
|
|
return r.position().Y();
|
231 |
|
|
}
|
232 |
|
|
/// z coordinate of vertex position
|
233 |
|
|
double vz() const {
|
234 |
|
|
const TrackingVertex& r = (*parentVertex_);
|
235 |
|
|
return r.position().Z();
|
236 |
|
|
}
|
237 |
|
|
/// status word
|
238 |
|
|
int status() const {
|
239 |
|
|
reco::GenParticleRefVector::const_iterator it = genParticle_begin();
|
240 |
|
|
return (*it)->status();
|
241 |
|
|
}
|
242 |
|
|
/// long lived flag
|
243 |
|
|
static const unsigned int longLivedTag;
|
244 |
|
|
/// is long lived?
|
245 |
|
|
bool longLived() const {
|
246 |
|
|
return status() & longLivedTag;
|
247 |
|
|
}
|
248 |
grimes |
1.1 |
|
249 |
|
|
private:
|
250 |
|
|
|
251 |
|
|
/// production time
|
252 |
|
|
float t_;
|
253 |
|
|
|
254 |
|
|
/// Total Number of Hits belonging to the TrackingParticle
|
255 |
|
|
int matchedHit_;
|
256 |
|
|
|
257 |
grimes |
1.1.2.1 |
/// references to G4 and reco::GenParticle tracks
|
258 |
grimes |
1.1 |
std::vector<SimTrack> g4Tracks_;
|
259 |
grimes |
1.1.2.1 |
reco::GenParticleRefVector genParticles_;
|
260 |
grimes |
1.1 |
|
261 |
|
|
// Source and decay vertices
|
262 |
|
|
TrackingVertexRef parentVertex_;
|
263 |
|
|
TrackingVertexRefVector decayVertices_;
|
264 |
|
|
};
|
265 |
|
|
|
266 |
|
|
#endif // SimDataFormats_TrackingParticle_H
|