00001 #ifndef DataFormat_ParticleFlowReco_PFDisplacedVertex_h
00002 #define DataFormat_ParticleFlowReco_PFDisplacedVertex_h
00003
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 #include "DataFormats/VertexReco/interface/Vertex.h"
00006 #include "DataFormats/Math/interface/LorentzVector.h"
00007 #include "TVector3.h"
00008
00009 #include <vector>
00010 #include <string>
00011 #include <iostream>
00012
00013 namespace reco {
00014
00015
00017
00025 class PFDisplacedVertex : public Vertex{
00026
00027 public:
00028
00030 typedef std::pair <unsigned int, unsigned int> PFTrackHitInfo;
00031 typedef std::pair <PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
00032
00034 enum M_Hypo {
00035 M_CUSTOM,
00036 M_MASSLESS,
00037 M_PION,
00038 M_KAON,
00039 M_LAMBDA
00040 };
00041
00042
00047 enum VertexTrackType {
00048 T_NOT_FROM_VERTEX,
00049 T_TO_VERTEX,
00050 T_FROM_VERTEX,
00051 T_MERGED
00052 };
00053
00056 enum VertexType {
00057 ANY = 0,
00058 FAKE = 1,
00059 LOOPER = 2,
00060 NUCL = 10,
00061 NUCL_LOOSE = 11,
00062 NUCL_KINK = 12,
00063 CONVERSION = 20,
00064 CONVERSION_LOOSE = 21,
00065 CONVERTED_BREMM = 22,
00066 K0_DECAY = 30,
00067 LAMBDA_DECAY = 31,
00068 LAMBDABAR_DECAY = 32,
00069 KPLUS_DECAY = 40,
00070 KMINUS_DECAY = 41,
00071 KPLUS_DECAY_LOOSE = 42,
00072 KMINUS_DECAY_LOOSE = 43,
00073 BSM_VERTEX = 100
00074 };
00075
00076
00078 PFDisplacedVertex();
00079
00081 PFDisplacedVertex(reco::Vertex&);
00082
00084 void addElement( const TrackBaseRef & r, const Track & refTrack,
00085 const PFTrackHitFullInfo& hitInfo ,
00086 VertexTrackType trackType = T_NOT_FROM_VERTEX, float w=1.0 );
00087
00089 void cleanTracks();
00090
00092 void setVertexType(VertexType vertexType) {vertexType_ = vertexType;}
00093
00096 void setPrimaryDirection(const math::XYZPoint& pvtx);
00097
00099 VertexType vertexType(){return vertexType_;}
00100
00101 const std::vector < PFTrackHitFullInfo > trackHitFullInfos() const
00102 {return trackHitFullInfos_;}
00103
00104 const std::vector <VertexTrackType> trackTypes() const
00105 {return trackTypes_;}
00106
00107
00109
00111 const bool isTherePrimaryTracks() const
00112 {return isThereKindTracks(T_TO_VERTEX);}
00113
00115 const bool isThereMergedTracks() const
00116 {return isThereKindTracks(T_MERGED);}
00117
00119 const bool isThereSecondaryTracks() const
00120 {return isThereKindTracks(T_FROM_VERTEX);}
00121
00123 const bool isThereNotFromVertexTracks() const
00124 {return isThereKindTracks(T_NOT_FROM_VERTEX);}
00125
00126
00127
00128
00130 const bool isPrimaryTrack(const reco::TrackBaseRef& originalTrack) const
00131 {
00132 size_t itrk = trackPosition(originalTrack);
00133 return isTrack(itrk, T_TO_VERTEX);
00134 }
00135
00137 const bool isSecondaryTrack(const reco::TrackBaseRef& originalTrack) const
00138 {
00139 size_t itrk = trackPosition(originalTrack);
00140 return isTrack(itrk, T_FROM_VERTEX);
00141 }
00142
00144 const bool isMergedTrack(const reco::TrackBaseRef& originalTrack) const
00145 {
00146 size_t itrk = trackPosition(originalTrack);
00147 return isTrack(itrk, T_MERGED);
00148 }
00149
00150 const PFTrackHitFullInfo trackHitFullInfo(const reco::TrackBaseRef& originalTrack) const{
00151 size_t itrk = trackPosition(originalTrack);
00152 return trackHitFullInfos_[itrk];
00153 }
00154
00155
00156
00158 const bool isIncomingTrack(const reco::TrackBaseRef& originalTrack) const
00159 {
00160 size_t itrk = trackPosition(originalTrack);
00161 return isTrack(itrk, T_MERGED) || isTrack(itrk, T_TO_VERTEX);
00162 }
00163
00165 const bool isOutgoingTrack(const reco::TrackBaseRef& originalTrack) const
00166 {
00167 size_t itrk = trackPosition(originalTrack);
00168 return isTrack(itrk, T_FROM_VERTEX);
00169 }
00170
00171
00173 const int nPrimaryTracks() const
00174 {return nKindTracks(T_TO_VERTEX);}
00175
00177 const int nMergedTracks() const
00178 {return nKindTracks(T_MERGED);}
00179
00181 const int nSecondaryTracks() const
00182 {return nKindTracks(T_FROM_VERTEX);}
00183
00185 const int nNotFromVertexTracks() const
00186 {return nKindTracks(T_NOT_FROM_VERTEX);}
00187
00189 const int nTracks() const {return trackTypes_.size();}
00190
00191
00192
00196 const math::XYZTLorentzVector
00197 secondaryMomentum(std::string massHypo = "PI",
00198 bool useRefitted = true, double mass = 0.0) const
00199 {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
00200
00202 const math::XYZTLorentzVector
00203 primaryMomentum(std::string massHypo = "PI",
00204 bool useRefitted = true, double mass = 0.0) const
00205 {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
00206
00207
00208
00212 const math::XYZTLorentzVector
00213 secondaryMomentum(M_Hypo massHypo,
00214 bool useRefitted = true, double mass = 0.0) const
00215 {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
00216
00218 const math::XYZTLorentzVector
00219 primaryMomentum(M_Hypo massHypo,
00220 bool useRefitted = true, double mass = 0.0) const
00221 {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
00222
00223
00224 void calcKinematics(){
00225 defaultPrimaryMomentum_ = momentum("PI", T_TO_VERTEX, false, 0.0);
00226 defaultSecondaryMomentum_ = momentum("PI", T_FROM_VERTEX, true, 0.0);
00227 }
00228
00229
00230 const double secondaryPt() const
00231 {return defaultPrimaryMomentum_.Pt();}
00232
00234 const double primaryPt() const
00235 {return defaultSecondaryMomentum_.Pt();}
00236
00237
00238
00240 const int totalCharge() const;
00241
00244 const double angle_io() const;
00245
00246
00248 const math::XYZVector primaryDirection() const;
00249
00250
00251
00252 bool isFake() const { return vertexType_ == FAKE;}
00253 bool isLooper() const { return vertexType_ == LOOPER;}
00254 bool isNucl() const { return vertexType_ == NUCL;}
00255 bool isNucl_Loose() const { return vertexType_ == NUCL_LOOSE;}
00256 bool isNucl_Kink() const { return vertexType_ == NUCL_KINK;}
00257 bool isConv() const { return vertexType_ == CONVERSION;}
00258 bool isConv_Loose() const { return vertexType_ == CONVERSION_LOOSE;}
00259 bool isConvertedBremm() const { return vertexType_ == CONVERTED_BREMM;}
00260 bool isK0() const { return vertexType_ == K0_DECAY;}
00261 bool isLambda() const { return vertexType_ == LAMBDA_DECAY;}
00262 bool isLambdaBar() const { return vertexType_ == LAMBDABAR_DECAY;}
00263 bool isKplus() const { return vertexType_ == KPLUS_DECAY;}
00264 bool isKminus() const { return vertexType_ == KMINUS_DECAY;}
00265 bool isKplus_Loose() const { return vertexType_ == KPLUS_DECAY_LOOSE;}
00266 bool isKminus_Loose() const { return vertexType_ == KMINUS_DECAY_LOOSE;}
00267 bool isBSM() const { return vertexType_ == BSM_VERTEX;}
00268
00269
00270 std::string nameVertexType() const;
00271
00273 void Dump(std::ostream& out = std::cout) const;
00274
00275 private:
00276
00278
00280 const bool isThereKindTracks(VertexTrackType) const;
00281
00283 const int nKindTracks(VertexTrackType) const;
00284
00286 const math::XYZTLorentzVector momentum(std::string,
00287 VertexTrackType,
00288 bool, double mass) const;
00289
00291 const math::XYZTLorentzVector momentum(M_Hypo massHypo,
00292 VertexTrackType,
00293 bool, double mass) const;
00294
00296 const double getMass2(M_Hypo, double) const;
00297
00298 const size_t trackPosition(const reco::TrackBaseRef& originalTrack) const;
00299
00300 const bool isTrack(size_t itrk, VertexTrackType T) const {
00301 return trackTypes_[itrk] == T;
00302 }
00303
00304
00306
00308 VertexType vertexType_;
00309
00311 std::vector < VertexTrackType > trackTypes_;
00312
00314 std::vector < PFTrackHitFullInfo > trackHitFullInfos_;
00315
00316 math::XYZVector primaryDirection_;
00317
00318 math::XYZTLorentzVector defaultPrimaryMomentum_;
00319 math::XYZTLorentzVector defaultSecondaryMomentum_;
00320
00321
00322 };
00323 }
00324
00325 #endif
00326
00327
00328