CMS 3D CMS Logo

PFDisplacedVertex.h
Go to the documentation of this file.
1 #ifndef DataFormat_ParticleFlowReco_PFDisplacedVertex_h
2 #define DataFormat_ParticleFlowReco_PFDisplacedVertex_h
3 
7 
8 #include <vector>
9 #include <string>
10 #include <iostream>
11 
12 namespace reco {
13 
15 
23  class PFDisplacedVertex : public Vertex {
24  public:
26  typedef std::pair<unsigned int, unsigned int> PFTrackHitInfo;
27  typedef std::pair<PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
28 
31 
37 
40  enum VertexType {
41  ANY = 0,
42  FAKE = 1,
43  LOOPER = 2,
44  NUCL = 10,
45  NUCL_LOOSE = 11,
46  NUCL_KINK = 12,
47  CONVERSION = 20,
50  K0_DECAY = 30,
57  BSM_VERTEX = 100
58  };
59 
62 
65 
67  void addElement(const TrackBaseRef& r,
68  const Track& refTrack,
69  const PFTrackHitFullInfo& hitInfo,
71  float w = 1.0);
72 
74  void cleanTracks();
75 
78 
81  void setPrimaryDirection(const math::XYZPoint& pvtx);
82 
85 
86  const std::vector<PFTrackHitFullInfo> trackHitFullInfos() const { return trackHitFullInfos_; }
87 
88  const std::vector<VertexTrackType> trackTypes() const { return trackTypes_; }
89 
91 
93  const bool isTherePrimaryTracks() const { return isThereKindTracks(T_TO_VERTEX); }
94 
96  const bool isThereMergedTracks() const { return isThereKindTracks(T_MERGED); }
97 
99  const bool isThereSecondaryTracks() const { return isThereKindTracks(T_FROM_VERTEX); }
100 
103 
106  size_t itrk = trackPosition(originalTrack);
107  return isTrack(itrk, T_TO_VERTEX);
108  }
109 
112  size_t itrk = trackPosition(originalTrack);
113  return isTrack(itrk, T_FROM_VERTEX);
114  }
115 
117  const bool isMergedTrack(const reco::TrackBaseRef& originalTrack) const {
118  size_t itrk = trackPosition(originalTrack);
119  return isTrack(itrk, T_MERGED);
120  }
121 
123  size_t itrk = trackPosition(originalTrack);
124  return trackHitFullInfos_[itrk];
125  }
126 
129  size_t itrk = trackPosition(originalTrack);
130  return isTrack(itrk, T_MERGED) || isTrack(itrk, T_TO_VERTEX);
131  }
132 
135  size_t itrk = trackPosition(originalTrack);
136  return isTrack(itrk, T_FROM_VERTEX);
137  }
138 
140  const int nPrimaryTracks() const { return nKindTracks(T_TO_VERTEX); }
141 
143  const int nMergedTracks() const { return nKindTracks(T_MERGED); }
144 
146  const int nSecondaryTracks() const { return nKindTracks(T_FROM_VERTEX); }
147 
149  const int nNotFromVertexTracks() const { return nKindTracks(T_NOT_FROM_VERTEX); }
150 
152  const int nTracks() const { return trackTypes_.size(); }
153 
154  // const reco::VertexTrackType vertexTrackType(reco::TrackBaseRef tkRef) const;
155 
160  bool useRefitted = true,
161  double mass = 0.0) const {
162  return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);
163  }
164 
167  bool useRefitted = true,
168  double mass = 0.0) const {
169  return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);
170  }
171 
175  const math::XYZTLorentzVector secondaryMomentum(M_Hypo massHypo, bool useRefitted = true, double mass = 0.0) const {
176  return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);
177  }
178 
180  const math::XYZTLorentzVector primaryMomentum(M_Hypo massHypo, bool useRefitted = true, double mass = 0.0) const {
181  return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);
182  }
183 
184  void calcKinematics() {
185  defaultPrimaryMomentum_ = momentum("PI", T_TO_VERTEX, false, 0.0);
187  }
188 
189  const double secondaryPt() const { return defaultPrimaryMomentum_.Pt(); }
190 
192  const double primaryPt() const { return defaultSecondaryMomentum_.Pt(); }
193 
195  const int totalCharge() const;
196 
199  const double angle_io() const;
200 
202  const math::XYZVector primaryDirection() const;
203 
204  bool isFake() const { return vertexType_ == FAKE; }
205  bool isLooper() const { return vertexType_ == LOOPER; }
206  bool isNucl() const { return vertexType_ == NUCL; }
207  bool isNucl_Loose() const { return vertexType_ == NUCL_LOOSE; }
208  bool isNucl_Kink() const { return vertexType_ == NUCL_KINK; }
209  bool isConv() const { return vertexType_ == CONVERSION; }
210  bool isConv_Loose() const { return vertexType_ == CONVERSION_LOOSE; }
211  bool isConvertedBremm() const { return vertexType_ == CONVERTED_BREMM; }
212  bool isK0() const { return vertexType_ == K0_DECAY; }
213  bool isLambda() const { return vertexType_ == LAMBDA_DECAY; }
214  bool isLambdaBar() const { return vertexType_ == LAMBDABAR_DECAY; }
215  bool isKplus() const { return vertexType_ == KPLUS_DECAY; }
216  bool isKminus() const { return vertexType_ == KMINUS_DECAY; }
217  bool isKplus_Loose() const { return vertexType_ == KPLUS_DECAY_LOOSE; }
218  bool isKminus_Loose() const { return vertexType_ == KMINUS_DECAY_LOOSE; }
219  bool isBSM() const { return vertexType_ == BSM_VERTEX; }
220 
221  std::string nameVertexType() const;
222 
224  void Dump(std::ostream& out = std::cout) const;
225 
226  private:
228 
230  const bool isThereKindTracks(VertexTrackType) const;
231 
233  const int nKindTracks(VertexTrackType) const;
234 
237 
239  const math::XYZTLorentzVector momentum(M_Hypo massHypo, VertexTrackType, bool, double mass) const;
240 
242  const double getMass2(M_Hypo, double) const;
243 
244  const size_t trackPosition(const reco::TrackBaseRef& originalTrack) const;
245 
246  const bool isTrack(size_t itrk, VertexTrackType T) const { return trackTypes_[itrk] == T; }
247 
249 
252 
254  std::vector<VertexTrackType> trackTypes_;
255 
257  std::vector<PFTrackHitFullInfo> trackHitFullInfos_;
258 
260 
263  };
264 } // namespace reco
265 
266 #endif
const math::XYZTLorentzVector primaryMomentum(std::string massHypo="PI", bool useRefitted=true, double mass=0.0) const
Momentum of primary or merged track calculated with a mass hypothesis.
const bool isPrimaryTrack(const reco::TrackBaseRef &originalTrack) const
Is a primary track was identified.
std::string nameVertexType() const
void setPrimaryDirection(const math::XYZPoint &pvtx)
const bool isTherePrimaryTracks() const
--—— Provide useful information --—— ///
T w() const
const int nMergedTracks() const
Number of merged tracks was identified.
void Dump(std::ostream &out=std::cout) const
cout function
const int nPrimaryTracks() const
Number of primary tracks was identified.
const bool isIncomingTrack(const reco::TrackBaseRef &originalTrack) const
Is primary or merged track.
const int nSecondaryTracks() const
Number of secondary tracks was identified.
const double secondaryPt() const
const int nNotFromVertexTracks() const
Number of tracks which was not identified.
const std::vector< VertexTrackType > trackTypes() const
const math::XYZVector primaryDirection() const
Primary Direction.
const int nTracks() const
Number of tracks.
const math::XYZTLorentzVector secondaryMomentum(std::string massHypo="PI", bool useRefitted=true, double mass=0.0) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const bool isThereNotFromVertexTracks() const
If there is a track which was not identified.
void addElement(const TrackBaseRef &r, const Track &refTrack, const PFTrackHitFullInfo &hitInfo, VertexTrackType trackType=T_NOT_FROM_VERTEX, float w=1.0)
Add a new track to the vertex.
std::vector< PFTrackHitFullInfo > trackHitFullInfos_
Information on the distance between track&#39;s hits and the Vertex.
void setVertexType(VertexType vertexType)
Set the type of this vertex.
std::pair< unsigned int, unsigned int > PFTrackHitInfo
Information on the distance between track&#39;s hits and the Vertex.
const bool isThereSecondaryTracks() const
If a secondary track was identified.
const int totalCharge() const
Total Charge.
const std::vector< PFTrackHitFullInfo > trackHitFullInfos() const
const size_t trackPosition(const reco::TrackBaseRef &originalTrack) const
const bool isThereMergedTracks() const
If a merged track was identified.
VertexType vertexType_
--—— MEMBERS --—— ///
Block of elements.
PFDisplacedVertex()
Default constructor.
std::pair< PFTrackHitInfo, PFTrackHitInfo > PFTrackHitFullInfo
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
const PFTrackHitFullInfo trackHitFullInfo(const reco::TrackBaseRef &originalTrack) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< VertexTrackType > trackTypes_
Types of the tracks associated to the vertex.
VertexType vertexType()
Get the type of this vertex.
const int nKindTracks(VertexTrackType) const
Common tool used to get the number of tracks of a given Kind.
const double getMass2(M_Hypo, double) const
Get the mass with a given hypothesis.
math::XYZVector primaryDirection_
const math::XYZTLorentzVector primaryMomentum(M_Hypo massHypo, bool useRefitted=true, double mass=0.0) const
Momentum of primary or merged track calculated with a mass hypothesis.
const math::XYZTLorentzVector momentum(std::string, VertexTrackType, bool, double mass) const
Common tool to calculate the momentum vector of tracks with a given Kind.
const bool isThereKindTracks(VertexTrackType) const
---—— TOOLS --------—— ///
fixed size matrix
void cleanTracks()
Clean the tracks collection and all the associated collections.
const bool isOutgoingTrack(const reco::TrackBaseRef &originalTrack) const
Is secondary track.
math::XYZTLorentzVector defaultPrimaryMomentum_
math::XYZTLorentzVector defaultSecondaryMomentum_
const double angle_io() const
const bool isMergedTrack(const reco::TrackBaseRef &originalTrack) const
Is a secondary track was identified.
const double primaryPt() const
Momentum of primary or merged track calculated with a mass hypothesis.
const bool isSecondaryTrack(const reco::TrackBaseRef &originalTrack) const
Is a secondary track was identified.
long double T
const math::XYZTLorentzVector secondaryMomentum(M_Hypo massHypo, bool useRefitted=true, double mass=0.0) const
M_Hypo
Mass hypothesis enum.
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:81
const bool isTrack(size_t itrk, VertexTrackType T) const