CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFDisplacedVertex.h
Go to the documentation of this file.
1 #ifndef DataFormat_ParticleFlowReco_PFDisplacedVertex_h
2 #define DataFormat_ParticleFlowReco_PFDisplacedVertex_h
3 
7 #include "TVector3.h"
8 
9 #include <vector>
10 #include <string>
11 #include <iostream>
12 
13 namespace reco {
14 
15 
17 
25  class PFDisplacedVertex : public Vertex{
26 
27  public:
28 
30  typedef std::pair <unsigned int, unsigned int> PFTrackHitInfo;
31  typedef std::pair <PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
32 
34  enum M_Hypo {
40  };
41 
42 
52  };
53 
56  enum VertexType {
57  ANY = 0,
58  FAKE = 1,
59  LOOPER = 2,
60  NUCL = 10,
61  NUCL_LOOSE = 11,
62  NUCL_KINK = 12,
63  CONVERSION = 20,
66  K0_DECAY = 30,
73  BSM_VERTEX = 100
74  };
75 
76 
79 
82 
84  void addElement( const TrackBaseRef & r, const Track & refTrack,
85  const PFTrackHitFullInfo& hitInfo ,
86  VertexTrackType trackType = T_NOT_FROM_VERTEX, float w=1.0 );
87 
89  void cleanTracks();
90 
93 
96  void setPrimaryDirection(const math::XYZPoint& pvtx);
97 
100 
101  const std::vector < PFTrackHitFullInfo > trackHitFullInfos() const
102  {return trackHitFullInfos_;}
103 
104  const std::vector <VertexTrackType> trackTypes() const
105  {return trackTypes_;}
106 
107 
109 
111  const bool isTherePrimaryTracks() const
112  {return isThereKindTracks(T_TO_VERTEX);}
113 
115  const bool isThereMergedTracks() const
116  {return isThereKindTracks(T_MERGED);}
117 
119  const bool isThereSecondaryTracks() const
121 
123  const bool isThereNotFromVertexTracks() const
125 
126 
127 
128 
131  {
132  size_t itrk = trackPosition(originalTrack);
133  return isTrack(itrk, T_TO_VERTEX);
134  }
135 
138  {
139  size_t itrk = trackPosition(originalTrack);
140  return isTrack(itrk, T_FROM_VERTEX);
141  }
142 
145  {
146  size_t itrk = trackPosition(originalTrack);
147  return isTrack(itrk, T_MERGED);
148  }
149 
151  size_t itrk = trackPosition(originalTrack);
152  return trackHitFullInfos_[itrk];
153  }
154 
155 
156 
159  {
160  size_t itrk = trackPosition(originalTrack);
161  return isTrack(itrk, T_MERGED) || isTrack(itrk, T_TO_VERTEX);
162  }
163 
166  {
167  size_t itrk = trackPosition(originalTrack);
168  return isTrack(itrk, T_FROM_VERTEX);
169  }
170 
171 
173  const int nPrimaryTracks() const
174  {return nKindTracks(T_TO_VERTEX);}
175 
177  const int nMergedTracks() const
178  {return nKindTracks(T_MERGED);}
179 
181  const int nSecondaryTracks() const
182  {return nKindTracks(T_FROM_VERTEX);}
183 
185  const int nNotFromVertexTracks() const
186  {return nKindTracks(T_NOT_FROM_VERTEX);}
187 
189  const int nTracks() const {return trackTypes_.size();}
190 
191  // const reco::VertexTrackType vertexTrackType(reco::TrackBaseRef tkRef) const;
192 
197  secondaryMomentum(std::string massHypo = "PI",
198  bool useRefitted = true, double mass = 0.0) const
199  {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
200 
203  primaryMomentum(std::string massHypo = "PI",
204  bool useRefitted = true, double mass = 0.0) const
205  {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
206 
207 
208 
214  bool useRefitted = true, double mass = 0.0) const
215  {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
216 
220  bool useRefitted = true, double mass = 0.0) const
221  {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
222 
223 
225  defaultPrimaryMomentum_ = momentum("PI", T_TO_VERTEX, false, 0.0);
227  }
228 
229 
230  const double secondaryPt() const
231  {return defaultPrimaryMomentum_.Pt();}
232 
234  const double primaryPt() const
235  {return defaultSecondaryMomentum_.Pt();}
236 
237 
238 
240  const int totalCharge() const;
241 
244  const double angle_io() const;
245 
246 
248  const math::XYZVector primaryDirection() const;
249 
250 
251 
252  bool isFake() const { return vertexType_ == FAKE;}
253  bool isLooper() const { return vertexType_ == LOOPER;}
254  bool isNucl() const { return vertexType_ == NUCL;}
255  bool isNucl_Loose() const { return vertexType_ == NUCL_LOOSE;}
256  bool isNucl_Kink() const { return vertexType_ == NUCL_KINK;}
257  bool isConv() const { return vertexType_ == CONVERSION;}
258  bool isConv_Loose() const { return vertexType_ == CONVERSION_LOOSE;}
259  bool isConvertedBremm() const { return vertexType_ == CONVERTED_BREMM;}
260  bool isK0() const { return vertexType_ == K0_DECAY;}
261  bool isLambda() const { return vertexType_ == LAMBDA_DECAY;}
262  bool isLambdaBar() const { return vertexType_ == LAMBDABAR_DECAY;}
263  bool isKplus() const { return vertexType_ == KPLUS_DECAY;}
264  bool isKminus() const { return vertexType_ == KMINUS_DECAY;}
265  bool isKplus_Loose() const { return vertexType_ == KPLUS_DECAY_LOOSE;}
266  bool isKminus_Loose() const { return vertexType_ == KMINUS_DECAY_LOOSE;}
267  bool isBSM() const { return vertexType_ == BSM_VERTEX;}
268 
269 
270  std::string nameVertexType() const;
271 
273  void Dump(std::ostream& out = std::cout) const;
274 
275  private:
276 
278 
280  const bool isThereKindTracks(VertexTrackType) const;
281 
283  const int nKindTracks(VertexTrackType) const;
284 
286  const math::XYZTLorentzVector momentum(std::string,
288  bool, double mass) const;
289 
291  const math::XYZTLorentzVector momentum(M_Hypo massHypo,
293  bool, double mass) const;
294 
296  const double getMass2(M_Hypo, double) const;
297 
298  const size_t trackPosition(const reco::TrackBaseRef& originalTrack) const;
299 
300  const bool isTrack(size_t itrk, VertexTrackType T) const {
301  return trackTypes_[itrk] == T;
302  }
303 
304 
306 
309 
311  std::vector < VertexTrackType > trackTypes_;
312 
314  std::vector < PFTrackHitFullInfo > trackHitFullInfos_;
315 
317 
320 
321 
322  };
323 }
324 
325 #endif
326 
327 
328 
const double getMass2(M_Hypo, double) const
Get the mass with a given hypothesis.
const bool isSecondaryTrack(const reco::TrackBaseRef &originalTrack) const
Is a secondary track was identified.
const int nKindTracks(VertexTrackType) const
Common tool used to get the number of tracks of a given Kind.
const double angle_io() const
const bool isTherePrimaryTracks() const
--—— Provide useful information --—— ///
void setPrimaryDirection(const math::XYZPoint &pvtx)
const std::vector< PFTrackHitFullInfo > trackHitFullInfos() const
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:86
const math::XYZTLorentzVector secondaryMomentum(std::string massHypo="PI", bool useRefitted=true, double mass=0.0) const
const int totalCharge() const
Total Charge.
const bool isIncomingTrack(const reco::TrackBaseRef &originalTrack) const
Is primary or merged track.
void Dump(std::ostream &out=std::cout) const
cout function
const bool isMergedTrack(const reco::TrackBaseRef &originalTrack) const
Is a secondary track was identified.
const PFTrackHitFullInfo trackHitFullInfo(const reco::TrackBaseRef &originalTrack) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
const bool isThereKindTracks(VertexTrackType) const
---—— TOOLS --------—— ///
const int nNotFromVertexTracks() const
Number of tracks 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.
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 int nPrimaryTracks() const
Number of primary tracks was identified.
const std::vector< VertexTrackType > trackTypes() const
const math::XYZVector primaryDirection() const
Primary Direction.
const bool isPrimaryTrack(const reco::TrackBaseRef &originalTrack) const
Is a primary track was identified.
const double primaryPt() const
Momentum of primary or merged track calculated with a mass hypothesis.
tuple out
Definition: dbtoconf.py:99
VertexType vertexType_
--—— MEMBERS --—— ///
const size_t trackPosition(const reco::TrackBaseRef &originalTrack) const
Block of elements.
PFDisplacedVertex()
Default constructor.
std::vector< PFTrackHitFullInfo > trackHitFullInfos_
Information on the distance between track&#39;s hits and the Vertex.
const bool isThereSecondaryTracks() const
If a secondary track was identified.
std::pair< PFTrackHitInfo, PFTrackHitInfo > PFTrackHitFullInfo
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
const double secondaryPt() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
const bool isOutgoingTrack(const reco::TrackBaseRef &originalTrack) const
Is secondary track.
const bool isTrack(size_t itrk, VertexTrackType T) const
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.
VertexType vertexType()
Get the type of this vertex.
const bool isThereNotFromVertexTracks() const
If there is a track which was not identified.
const math::XYZTLorentzVector momentum(std::string, VertexTrackType, bool, double mass) const
Common tool to calculate the momentum vector of tracks with a given Kind.
math::XYZVector primaryDirection_
void cleanTracks()
Clean the tracks collection and all the associated collections.
std::string nameVertexType() const
tuple mass
Definition: scaleCards.py:27
const int nMergedTracks() const
Number of merged tracks was identified.
math::XYZTLorentzVector defaultPrimaryMomentum_
math::XYZTLorentzVector defaultSecondaryMomentum_
const math::XYZTLorentzVector secondaryMomentum(M_Hypo massHypo, bool useRefitted=true, double mass=0.0) const
tuple cout
Definition: gather_cfg.py:121
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 bool isThereMergedTracks() const
If a merged track was identified.
long double T
const int nSecondaryTracks() const
Number of secondary tracks was identified.
M_Hypo
Mass hypothesis enum.
std::vector< VertexTrackType > trackTypes_
Types of the tracks associated to the vertex.
T w() const
const int nTracks() const
Number of tracks.