File indexing completed on 2024-04-06 12:29:43
0001 #ifndef SimDataFormats_TrackingParticle_h
0002 #define SimDataFormats_TrackingParticle_h
0003
0004 #include <vector>
0005 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
0006 #include "DataFormats/Math/interface/Point3D.h"
0007 #include "DataFormats/Math/interface/Vector3D.h"
0008 #include "DataFormats/Math/interface/LorentzVector.h"
0009 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
0010 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0011
0012
0013
0014
0015 class TrackingVertex;
0016 class SimTrack;
0017 class EncodedEventId;
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 class TrackingParticle {
0030 friend std::ostream& operator<<(std::ostream& s, TrackingParticle const& tp);
0031
0032 public:
0033 typedef int Charge;
0034 typedef math::XYZTLorentzVectorD LorentzVector;
0035 typedef math::PtEtaPhiMLorentzVector PolarLorentzVector;
0036 typedef math::XYZPointD Point;
0037 typedef math::XYZVectorD Vector;
0038
0039
0040 typedef reco::GenParticleRefVector::iterator genp_iterator;
0041 typedef std::vector<SimTrack>::const_iterator g4t_iterator;
0042
0043
0044
0045
0046
0047
0048
0049
0050 TrackingParticle();
0051
0052 TrackingParticle(const SimTrack& simtrk, const TrackingVertexRef& parentVertex);
0053
0054
0055 ~TrackingParticle();
0056
0057
0058
0059
0060
0061 int pdgId() const {
0062 if (genParticles_.empty())
0063 return g4Tracks_[0].type();
0064 else
0065 return (*genParticles_.begin())->pdgId();
0066 }
0067
0068
0069
0070
0071
0072 EncodedEventId eventId() const { return g4Tracks_[0].eventId(); }
0073
0074
0075 void addGenParticle(const reco::GenParticleRef& ref);
0076 void addG4Track(const SimTrack& t);
0077
0078 genp_iterator genParticle_begin() const;
0079 genp_iterator genParticle_end() const;
0080 g4t_iterator g4Track_begin() const;
0081 g4t_iterator g4Track_end() const;
0082 void setParentVertex(const TrackingVertexRef& ref);
0083 void addDecayVertex(const TrackingVertexRef& ref);
0084 void clearParentVertex();
0085 void clearDecayVertices();
0086
0087
0088 const reco::GenParticleRefVector& genParticles() const { return genParticles_; }
0089 const std::vector<SimTrack>& g4Tracks() const { return g4Tracks_; }
0090 const TrackingVertexRef& parentVertex() const { return parentVertex_; }
0091
0092
0093 const TrackingVertexRefVector& decayVertices() const { return decayVertices_; }
0094 tv_iterator decayVertices_begin() const { return decayVertices_.begin(); }
0095 tv_iterator decayVertices_end() const { return decayVertices_.end(); }
0096
0097
0098 float charge() const { return g4Tracks_[0].charge(); }
0099
0100 int threeCharge() const { return lrintf(3.f * charge()); }
0101
0102
0103 const LorentzVector& p4() const { return g4Tracks_[0].momentum(); }
0104
0105
0106 Vector momentum() const { return p4().Vect(); }
0107
0108
0109 Vector boostToCM() const { return p4().BoostToCM(); }
0110
0111
0112 double p() const { return p4().P(); }
0113
0114
0115 double qoverp() const { return charge() / p(); }
0116
0117
0118 double energy() const { return p4().E(); }
0119
0120
0121 double et() const { return p4().Et(); }
0122
0123
0124 double mass() const { return p4().M(); }
0125
0126
0127 double massSqr() const { return pow(mass(), 2); }
0128
0129
0130 double mt() const { return p4().Mt(); }
0131
0132
0133 double mtSqr() const { return p4().Mt2(); }
0134
0135
0136 double px() const { return p4().Px(); }
0137
0138
0139 double py() const { return p4().Py(); }
0140
0141
0142 double pz() const { return p4().Pz(); }
0143
0144
0145 double pt() const { return p4().Pt(); }
0146
0147
0148 double phi() const { return p4().Phi(); }
0149
0150
0151 double theta() const { return p4().Theta(); }
0152
0153
0154 double eta() const { return p4().Eta(); }
0155
0156
0157 double lambda() const { return M_PI_2 - theta(); }
0158
0159
0160 double tanl() const { return tan(lambda()); }
0161
0162
0163 double rapidity() const { return p4().Rapidity(); }
0164
0165
0166 double y() const { return rapidity(); }
0167
0168
0169 Point vertex() const {
0170 const TrackingVertex::LorentzVector& p = (*parentVertex_).position();
0171 return Point(p.x(), p.y(), p.z());
0172 }
0173
0174
0175 double vx() const {
0176 const TrackingVertex& r = (*parentVertex_);
0177 return r.position().X();
0178 }
0179
0180
0181 double vy() const {
0182 const TrackingVertex& r = (*parentVertex_);
0183 return r.position().Y();
0184 }
0185
0186
0187 double vz() const {
0188 const TrackingVertex& r = (*parentVertex_);
0189 return r.position().Z();
0190 }
0191
0192
0193 double dxy() const { return (-vx() * py() + vy() * px()) / pt(); }
0194
0195
0196 double d0() const { return -dxy(); }
0197
0198
0199 double dz() const { return vz() - (vx() * px() + vy() * py()) * pz() / p4().Perp2(); }
0200
0201
0202 double z0() const { return dz(); }
0203
0204
0205
0206
0207 int status() const { return genParticles_.empty() ? -99 : (*genParticles_[0]).status(); }
0208
0209 static const unsigned int longLivedTag;
0210
0211
0212 bool longLived() const { return status() & longLivedTag; }
0213
0214
0215
0216
0217 int numberOfHits() const { return numberOfHits_; }
0218
0219
0220
0221
0222 int numberOfTrackerHits() const { return numberOfTrackerHits_; }
0223
0224
0225
0226 int matchedHit() const;
0227
0228
0229
0230
0231 int numberOfTrackerLayers() const { return numberOfTrackerLayers_; }
0232
0233 void setNumberOfHits(int numberOfHits);
0234 void setNumberOfTrackerHits(int numberOfTrackerHits);
0235 void setNumberOfTrackerLayers(const int numberOfTrackerLayers);
0236
0237 private:
0238 int numberOfHits_;
0239 int numberOfTrackerHits_;
0240 int numberOfTrackerLayers_;
0241
0242
0243 std::vector<SimTrack> g4Tracks_;
0244 reco::GenParticleRefVector genParticles_;
0245
0246
0247 TrackingVertexRef parentVertex_;
0248 TrackingVertexRefVector decayVertices_;
0249 };
0250
0251 #endif