CMS 3D CMS Logo

SimTrack.h
Go to the documentation of this file.
1 #ifndef SimTrack_H
2 #define SimTrack_H
3 
8 
9 class SimTrack : public CoreSimTrack {
10 public:
11  typedef CoreSimTrack Core;
12 
14  SimTrack();
15  SimTrack(int ipart, const math::XYZTLorentzVectorD& p);
16 
20  SimTrack(int ipart, const math::XYZTLorentzVectorD& p, int iv, int ig);
21 
22  SimTrack(int ipart,
24  int iv,
25  int ig,
26  const math::XYZVectorD& tkp,
27  const math::XYZTLorentzVectorD& tkm);
28 
30  SimTrack(const CoreSimTrack& t, int iv, int ig);
31 
33  int vertIndex() const { return ivert; }
34  bool noVertex() const { return ivert == -1; }
35 
37  int genpartIndex() const { return igenpart; }
38  bool noGenpart() const { return igenpart == -1; }
39 
41 
43 
44  inline void setTkPosition(const math::XYZVectorD& pos) { tkposition = pos; }
45 
46  inline void setTkMomentum(const math::XYZTLorentzVectorD& mom) { tkmomentum = mom; }
47 
48  inline void setVertexIndex(const int v) { ivert = v; }
49 
51  int idAtBoundary,
52  math::XYZTLorentzVectorF positionAtBoundary,
53  math::XYZTLorentzVectorF momentumAtBoundary) {
55  idAtBoundary_ = idAtBoundary;
56  positionAtBoundary_ = positionAtBoundary;
57  momentumAtBoundary_ = momentumAtBoundary;
58  }
59  bool crossedBoundary() const { return crossedBoundary_; }
62  int getIDAtBoundary() const { return idAtBoundary_; }
63 
64 private:
65  int ivert;
66  int igenpart;
67 
70 
75 };
76 
77 #include <iosfwd>
78 std::ostream& operator<<(std::ostream& o, const SimTrack& t);
79 
80 #endif
std::ostream & operator<<(std::ostream &o, const SimTrack &t)
Definition: SimTrack.cc:21
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:40
CoreSimTrack Core
Definition: SimTrack.h:11
int getIDAtBoundary() const
Definition: SimTrack.h:62
SimTrack()
constructor
Definition: SimTrack.cc:3
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:33
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
int igenpart
Definition: SimTrack.h:66
const math::XYZTLorentzVectorF & getPositionAtBoundary() const
Definition: SimTrack.h:60
math::XYZTLorentzVectorF momentumAtBoundary_
Definition: SimTrack.h:74
bool crossedBoundary() const
Definition: SimTrack.h:59
const math::XYZTLorentzVectorF & getMomentumAtBoundary() const
Definition: SimTrack.h:61
void setVertexIndex(const int v)
Definition: SimTrack.h:48
void setCrossedBoundaryVars(bool crossedBoundary, int idAtBoundary, math::XYZTLorentzVectorF positionAtBoundary, math::XYZTLorentzVectorF momentumAtBoundary)
Definition: SimTrack.h:50
int idAtBoundary_
Definition: SimTrack.h:72
bool noVertex() const
Definition: SimTrack.h:34
math::XYZTLorentzVectorD tkmomentum
Definition: SimTrack.h:69
void setTkPosition(const math::XYZVectorD &pos)
Definition: SimTrack.h:44
bool noGenpart() const
Definition: SimTrack.h:38
int ivert
Definition: SimTrack.h:65
math::XYZVectorD tkposition
Definition: SimTrack.h:68
void setTkMomentum(const math::XYZTLorentzVectorD &mom)
Definition: SimTrack.h:46
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:42
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
bool crossedBoundary_
Definition: SimTrack.h:71
math::XYZTLorentzVectorF positionAtBoundary_
Definition: SimTrack.h:73
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:37