File indexing completed on 2024-04-06 12:04:51
0001 #ifndef DataFormat_ParticleFlowReco_PFDisplacedVertex_h
0002 #define DataFormat_ParticleFlowReco_PFDisplacedVertex_h
0003
0004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0005 #include "DataFormats/VertexReco/interface/Vertex.h"
0006 #include "DataFormats/Math/interface/LorentzVector.h"
0007
0008 #include <vector>
0009 #include <string>
0010 #include <iostream>
0011
0012 namespace reco {
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023 class PFDisplacedVertex : public Vertex {
0024 public:
0025
0026 typedef std::pair<unsigned int, unsigned int> PFTrackHitInfo;
0027 typedef std::pair<PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
0028
0029
0030 enum M_Hypo { M_CUSTOM, M_MASSLESS, M_PION, M_KAON, M_LAMBDA };
0031
0032
0033
0034
0035
0036 enum VertexTrackType { T_NOT_FROM_VERTEX, T_TO_VERTEX, T_FROM_VERTEX, T_MERGED };
0037
0038
0039
0040 enum VertexType {
0041 ANY = 0,
0042 FAKE = 1,
0043 LOOPER = 2,
0044 NUCL = 10,
0045 NUCL_LOOSE = 11,
0046 NUCL_KINK = 12,
0047 CONVERSION = 20,
0048 CONVERSION_LOOSE = 21,
0049 CONVERTED_BREMM = 22,
0050 K0_DECAY = 30,
0051 LAMBDA_DECAY = 31,
0052 LAMBDABAR_DECAY = 32,
0053 KPLUS_DECAY = 40,
0054 KMINUS_DECAY = 41,
0055 KPLUS_DECAY_LOOSE = 42,
0056 KMINUS_DECAY_LOOSE = 43,
0057 BSM_VERTEX = 100
0058 };
0059
0060
0061 PFDisplacedVertex();
0062
0063
0064 PFDisplacedVertex(reco::Vertex&);
0065
0066
0067 void addElement(const TrackBaseRef& r,
0068 const Track& refTrack,
0069 const PFTrackHitFullInfo& hitInfo,
0070 VertexTrackType trackType = T_NOT_FROM_VERTEX,
0071 float w = 1.0);
0072
0073
0074 void cleanTracks();
0075
0076
0077 void setVertexType(VertexType vertexType) { vertexType_ = vertexType; }
0078
0079
0080
0081 void setPrimaryDirection(const math::XYZPoint& pvtx);
0082
0083
0084 VertexType vertexType() { return vertexType_; }
0085
0086 const std::vector<PFTrackHitFullInfo> trackHitFullInfos() const { return trackHitFullInfos_; }
0087
0088 const std::vector<VertexTrackType> trackTypes() const { return trackTypes_; }
0089
0090
0091
0092
0093 const bool isTherePrimaryTracks() const { return isThereKindTracks(T_TO_VERTEX); }
0094
0095
0096 const bool isThereMergedTracks() const { return isThereKindTracks(T_MERGED); }
0097
0098
0099 const bool isThereSecondaryTracks() const { return isThereKindTracks(T_FROM_VERTEX); }
0100
0101
0102 const bool isThereNotFromVertexTracks() const { return isThereKindTracks(T_NOT_FROM_VERTEX); }
0103
0104
0105 const bool isPrimaryTrack(const reco::TrackBaseRef& originalTrack) const {
0106 size_t itrk = trackPosition(originalTrack);
0107 return isTrack(itrk, T_TO_VERTEX);
0108 }
0109
0110
0111 const bool isSecondaryTrack(const reco::TrackBaseRef& originalTrack) const {
0112 size_t itrk = trackPosition(originalTrack);
0113 return isTrack(itrk, T_FROM_VERTEX);
0114 }
0115
0116
0117 const bool isMergedTrack(const reco::TrackBaseRef& originalTrack) const {
0118 size_t itrk = trackPosition(originalTrack);
0119 return isTrack(itrk, T_MERGED);
0120 }
0121
0122 const PFTrackHitFullInfo trackHitFullInfo(const reco::TrackBaseRef& originalTrack) const {
0123 size_t itrk = trackPosition(originalTrack);
0124 return trackHitFullInfos_[itrk];
0125 }
0126
0127
0128 const bool isIncomingTrack(const reco::TrackBaseRef& originalTrack) const {
0129 size_t itrk = trackPosition(originalTrack);
0130 return isTrack(itrk, T_MERGED) || isTrack(itrk, T_TO_VERTEX);
0131 }
0132
0133
0134 const bool isOutgoingTrack(const reco::TrackBaseRef& originalTrack) const {
0135 size_t itrk = trackPosition(originalTrack);
0136 return isTrack(itrk, T_FROM_VERTEX);
0137 }
0138
0139
0140 const int nPrimaryTracks() const { return nKindTracks(T_TO_VERTEX); }
0141
0142
0143 const int nMergedTracks() const { return nKindTracks(T_MERGED); }
0144
0145
0146 const int nSecondaryTracks() const { return nKindTracks(T_FROM_VERTEX); }
0147
0148
0149 const int nNotFromVertexTracks() const { return nKindTracks(T_NOT_FROM_VERTEX); }
0150
0151
0152 const int nTracks() const { return trackTypes_.size(); }
0153
0154
0155
0156
0157
0158
0159 const math::XYZTLorentzVector secondaryMomentum(std::string massHypo = "PI",
0160 bool useRefitted = true,
0161 double mass = 0.0) const {
0162 return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);
0163 }
0164
0165
0166 const math::XYZTLorentzVector primaryMomentum(std::string massHypo = "PI",
0167 bool useRefitted = true,
0168 double mass = 0.0) const {
0169 return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);
0170 }
0171
0172
0173
0174
0175 const math::XYZTLorentzVector secondaryMomentum(M_Hypo massHypo, bool useRefitted = true, double mass = 0.0) const {
0176 return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);
0177 }
0178
0179
0180 const math::XYZTLorentzVector primaryMomentum(M_Hypo massHypo, bool useRefitted = true, double mass = 0.0) const {
0181 return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);
0182 }
0183
0184 void calcKinematics() {
0185 defaultPrimaryMomentum_ = momentum("PI", T_TO_VERTEX, false, 0.0);
0186 defaultSecondaryMomentum_ = momentum("PI", T_FROM_VERTEX, true, 0.0);
0187 }
0188
0189 const double secondaryPt() const { return defaultPrimaryMomentum_.Pt(); }
0190
0191
0192 const double primaryPt() const { return defaultSecondaryMomentum_.Pt(); }
0193
0194
0195 const int totalCharge() const;
0196
0197
0198
0199 const double angle_io() const;
0200
0201
0202 const math::XYZVector primaryDirection() const;
0203
0204 bool isFake() const { return vertexType_ == FAKE; }
0205 bool isLooper() const { return vertexType_ == LOOPER; }
0206 bool isNucl() const { return vertexType_ == NUCL; }
0207 bool isNucl_Loose() const { return vertexType_ == NUCL_LOOSE; }
0208 bool isNucl_Kink() const { return vertexType_ == NUCL_KINK; }
0209 bool isConv() const { return vertexType_ == CONVERSION; }
0210 bool isConv_Loose() const { return vertexType_ == CONVERSION_LOOSE; }
0211 bool isConvertedBremm() const { return vertexType_ == CONVERTED_BREMM; }
0212 bool isK0() const { return vertexType_ == K0_DECAY; }
0213 bool isLambda() const { return vertexType_ == LAMBDA_DECAY; }
0214 bool isLambdaBar() const { return vertexType_ == LAMBDABAR_DECAY; }
0215 bool isKplus() const { return vertexType_ == KPLUS_DECAY; }
0216 bool isKminus() const { return vertexType_ == KMINUS_DECAY; }
0217 bool isKplus_Loose() const { return vertexType_ == KPLUS_DECAY_LOOSE; }
0218 bool isKminus_Loose() const { return vertexType_ == KMINUS_DECAY_LOOSE; }
0219 bool isBSM() const { return vertexType_ == BSM_VERTEX; }
0220
0221 std::string nameVertexType() const;
0222
0223
0224 void Dump(std::ostream& out = std::cout) const;
0225
0226 private:
0227
0228
0229
0230 const bool isThereKindTracks(VertexTrackType) const;
0231
0232
0233 const int nKindTracks(VertexTrackType) const;
0234
0235
0236 const math::XYZTLorentzVector momentum(std::string, VertexTrackType, bool, double mass) const;
0237
0238
0239 const math::XYZTLorentzVector momentum(M_Hypo massHypo, VertexTrackType, bool, double mass) const;
0240
0241
0242 const double getMass2(M_Hypo, double) const;
0243
0244 const size_t trackPosition(const reco::TrackBaseRef& originalTrack) const;
0245
0246 const bool isTrack(size_t itrk, VertexTrackType T) const { return trackTypes_[itrk] == T; }
0247
0248
0249
0250
0251 VertexType vertexType_;
0252
0253
0254 std::vector<VertexTrackType> trackTypes_;
0255
0256
0257 std::vector<PFTrackHitFullInfo> trackHitFullInfos_;
0258
0259 math::XYZVector primaryDirection_;
0260
0261 math::XYZTLorentzVector defaultPrimaryMomentum_;
0262 math::XYZTLorentzVector defaultSecondaryMomentum_;
0263 };
0264 }
0265
0266 #endif