CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingVertex.h
Go to the documentation of this file.
1 #ifndef SimDataFormats_TrackingVertex_h
2 #define SimDataFormats_TrackingVertex_h
3 
15 
21 
23 {
24 
25  friend std::ostream& operator<< (std::ostream& s, const TrackingVertex & tv);
26 
27 public:
28 
33  typedef std::vector<SimVertex>::const_iterator g4v_iterator;
35 
36 // Default constructor and constructor from values
38  TrackingVertex(const LorentzVector &position, const bool inVolume,
39  const EncodedEventId e = EncodedEventId(0));
40 
41 // Setters
43  {
44  eId_=e;
45  };
46 
47 // Track and vertex iterators
48  genv_iterator genVertices_begin() const; // Ref's to HepMC and Geant4
49  genv_iterator genVertices_end() const; // vertices associated with
50  g4v_iterator g4Vertices_begin() const; // this vertex, respectively
51  g4v_iterator g4Vertices_end() const; // ....
52 
53  tp_iterator daughterTracks_begin() const; // Ref's to daughter and source
54  tp_iterator daughterTracks_end() const; // tracks associated with
55  tp_iterator sourceTracks_begin() const; // this vertex, respectively
56  tp_iterator sourceTracks_end() const; // ....
57 
58  unsigned int nG4Vertices() const
59  {
60  return g4Vertices_.size();
61  };
62  unsigned int nGenVertices() const
63  {
64  return genVertices_.size();
65  };
66  unsigned int nDaughterTracks() const
67  {
68  return daughterTracks_.size();
69  };
70  unsigned int nSourceTracks() const
71  {
72  return sourceTracks_.size();
73  };
74 
75 // Add references to TrackingParticles, Geant4, and HepMC vertices to containers
76  void addG4Vertex( const SimVertex& );
77  void addGenVertex( const GenVertexRef& );
79  void addParentTrack( const TrackingParticleRef&);
80  void clearDaughterTracks();
81  void clearParentTracks();
82 
83 // Getters for RefVectors
84  const std::vector<SimVertex>& g4Vertices() const;
85  const GenVertexRefVector& genVertices() const;
88 
89 // Getters for other info
90  const LorentzVector& position() const
91  {
92  return position_;
93  };
94  const EncodedEventId& eventId() const
95  {
96  return eId_;
97  };
98  const bool inVolume() const
99  {
100  return inVolume_;
101  };
102 
103 private:
104 
105  LorentzVector position_; // Vertex position and time
106  bool inVolume_; // Is it inside tracker volume?
108 
109 // References to G4 and generator vertices and TrackingParticles
110 
111  std::vector<SimVertex> g4Vertices_;
115 };
116 
117 #endif
edm::RefVector< edm::HepMCProduct, HepMC::GenVertex > GenVertexRefVector
tp_iterator daughterTracks_begin() const
const TrackingParticleRefVector & sourceTracks() const
void setEventId(EncodedEventId e)
GenVertexRefVector::iterator genv_iterator
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
TrackingVertex::g4v_iterator g4v_iterator
std::vector< SimVertex >::const_iterator g4v_iterator
LorentzVector position_
const std::vector< SimVertex > & g4Vertices() const
g4v_iterator g4Vertices_end() const
const GenVertexRefVector & genVertices() const
void addDaughterTrack(const TrackingParticleRef &)
EncodedEventId eId_
void clearDaughterTracks()
void addParentTrack(const TrackingParticleRef &)
const bool inVolume() const
g4v_iterator g4Vertices_begin() const
void addG4Vertex(const SimVertex &)
math::XYZTLorentzVectorD LorentzVector
genv_iterator genVertices_end() const
TrackingParticleRefVector daughterTracks_
friend std::ostream & operator<<(std::ostream &s, const TrackingVertex &tv)
GenVertexRefVector genVertices_
tp_iterator daughterTracks_end() const
unsigned int nG4Vertices() const
tp_iterator sourceTracks_begin() const
unsigned int nGenVertices() const
unsigned int nDaughterTracks() const
tp_iterator sourceTracks_end() const
void clearParentTracks()
std::vector< SimVertex > g4Vertices_
void addGenVertex(const GenVertexRef &)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > LorentzVector
Definition: analysisEnums.h:9
const EncodedEventId & eventId() const
const TrackingParticleRefVector & daughterTracks() const
TrackingParticleRefVector::iterator tp_iterator
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
TrackingParticleRefVector sourceTracks_
edm::Ref< edm::HepMCProduct, HepMC::GenVertex > GenVertexRef
genv_iterator genVertices_begin() const
unsigned int nSourceTracks() const
const LorentzVector & position() const