ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/grimes/TrackingTruthCode/SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h
Revision: 1.1.2.1
Committed: Wed Mar 20 17:40:03 2013 UTC (12 years, 1 month ago) by grimes
Content type: text/plain
Branch: NewTrackingParticles
Changes since 1.1: +202 -84 lines
Log Message:
First commit of the new TrackingParticle implementation

File Contents

# User Rev Content
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 grimes 1.1.2.1 if( genParticles_.empty() ) return 0;
62    
63     reco::GenParticleRefVector::const_iterator it = genParticle_begin();
64     return (*it)->pdgId();
65 grimes 1.1 }
66     EncodedEventId eventId() const
67     {
68 grimes 1.1.2.1 return g4Tracks_.at(0).eventId();
69 grimes 1.1 }
70    
71 grimes 1.1.2.1 // Setters for G4 and reco::GenParticle
72     void addGenParticle( const reco::GenParticleRef& ref) {
73     genParticles_.push_back(ref);
74     }
75     void addG4Track( const SimTrack& t) {
76     g4Tracks_.push_back(t);
77     }
78     /// iterators
79     genp_iterator genParticle_begin() const {
80     return genParticles_.begin();
81     }
82     genp_iterator genParticle_end() const {
83     return genParticles_.end();
84     }
85     g4t_iterator g4Track_begin() const {
86     return g4Tracks_.begin();
87     }
88     g4t_iterator g4Track_end() const {
89     return g4Tracks_.end();
90     }
91     void setParentVertex(const TrackingVertexRef& ref) {
92     parentVertex_ = ref;
93     }
94     void addDecayVertex(const TrackingVertexRef& ref) {
95     decayVertices_.push_back(ref);
96     }
97     void clearParentVertex() {
98     parentVertex_ = TrackingVertexRef();
99     }
100     void clearDecayVertices() {
101     decayVertices_.clear();
102     }
103     void setMatchedHit(const int& hitnumb) {
104     matchedHit_ = hitnumb;
105     }
106 grimes 1.1 // Getters for Embd and Sim Tracks
107 grimes 1.1.2.1 const reco::GenParticleRefVector& genParticles() const {
108 grimes 1.1 return genParticles_;
109     }
110 grimes 1.1.2.1 const std::vector<SimTrack>& g4Tracks() const {
111     return g4Tracks_;
112 grimes 1.1 }
113 grimes 1.1.2.1 const TrackingVertexRef& parentVertex() const {
114 grimes 1.1 return parentVertex_;
115     }
116    
117     // Accessors for vector of decay vertices
118 grimes 1.1.2.1 const TrackingVertexRefVector& decayVertices() const {
119 grimes 1.1 return decayVertices_;
120     }
121 grimes 1.1.2.1 tv_iterator decayVertices_begin() const {
122 grimes 1.1 return decayVertices_.begin();
123     }
124 grimes 1.1.2.1 tv_iterator decayVertices_end() const {
125 grimes 1.1 return decayVertices_.end();
126     }
127 grimes 1.1.2.1 int matchedHit() const {
128 grimes 1.1 return matchedHit_;
129     }
130 grimes 1.1.2.1 /// electric charge
131     int charge() const {
132     return g4Tracks_.at(0).charge() / 3;
133     }
134     /// electric charge
135     int threeCharge() const {
136     return g4Tracks_.at(0).charge();
137     }
138     /// four-momentum Lorentz vector
139     const LorentzVector& p4() const {
140     return g4Tracks_.at(0).momentum();
141     }
142     #if 0
143     /// four-momentum Lorentz vector
144     const PolarLorentzVector& polarP4() const
145     {
146     }
147     #endif
148     /// spatial momentum vector
149     Vector momentum() const {
150     return p4().Vect();
151     }
152     /// boost vector to boost a Lorentz vector
153     /// to the particle center of mass system
154     Vector boostToCM() const {
155     return p4().BoostToCM();
156     }
157     /// magnitude of momentum vector
158     double p() const {
159     return p4().P();
160     }
161     /// energy
162     double energy() const {
163     return p4().E();
164     }
165     /// transverse energy
166     double et() const {
167     return p4().Et();
168     }
169     /// mass
170     double mass() const {
171     return p4().M();
172     }
173     /// mass squared
174     double massSqr() const {
175     return pow(mass(), 2);
176     }
177     /// transverse mass
178     double mt() const {
179     return p4().Mt();
180     }
181     /// transverse mass squared
182     double mtSqr() const {
183     return p4().Mt2();
184     }
185     /// x coordinate of momentum vector
186     double px() const {
187     return p4().Px();
188     }
189     /// y coordinate of momentum vector
190     double py() const {
191     return p4().Py();
192     }
193     /// z coordinate of momentum vector
194     double pz() const {
195     return p4().Pz();
196     }
197     /// transverse momentum
198     double pt() const {
199     return p4().Pt();
200     }
201     /// momentum azimuthal angle
202     double phi() const {
203     return p4().Phi();
204     }
205     /// momentum polar angle
206     double theta() const {
207     return p4().Theta();
208     }
209     /// momentum pseudorapidity
210     double eta() const {
211     return p4().Eta();
212     }
213     /// rapidity
214     double rapidity() const {
215     return p4().Rapidity();
216     }
217     /// rapidity
218     double y() const {
219     return rapidity();
220     }
221     /// vertex position
222     Point vertex() const {
223     return Point(vx(), vy(), vz());
224     }
225     /// x coordinate of vertex position
226     double vx() const {
227     const TrackingVertex& r = (*parentVertex_);
228     return r.position().X();
229     }
230     /// y coordinate of vertex position
231     double vy() const {
232     const TrackingVertex& r = (*parentVertex_);
233     return r.position().Y();
234     }
235     /// z coordinate of vertex position
236     double vz() const {
237     const TrackingVertex& r = (*parentVertex_);
238     return r.position().Z();
239     }
240     /// status word
241     int status() const {
242     reco::GenParticleRefVector::const_iterator it = genParticle_begin();
243     return (*it)->status();
244     }
245     /// long lived flag
246     static const unsigned int longLivedTag;
247     /// is long lived?
248     bool longLived() const {
249     return status() & longLivedTag;
250     }
251 grimes 1.1
252     private:
253    
254     /// production time
255     float t_;
256    
257     /// Total Number of Hits belonging to the TrackingParticle
258     int matchedHit_;
259    
260 grimes 1.1.2.1 /// references to G4 and reco::GenParticle tracks
261 grimes 1.1 std::vector<SimTrack> g4Tracks_;
262 grimes 1.1.2.1 reco::GenParticleRefVector genParticles_;
263 grimes 1.1
264     // Source and decay vertices
265     TrackingVertexRef parentVertex_;
266     TrackingVertexRefVector decayVertices_;
267     };
268    
269     #endif // SimDataFormats_TrackingParticle_H