SND@LHC Software
Loading...
Searching...
No Matches
ShipMCTrack.h
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- ShipMCTrack header file -----
3// -------------------------------------------------------------------------
4
5
13#ifndef ShipMCTrack_H
14#define ShipMCTrack_H 1
15
16#include "TObject.h" // for TObject
17#include "ShipDetectorList.h" // for DetectorId
18#include "Rtypes.h" // for Double_t, Int_t, Double32_t, etc
19#include "TLorentzVector.h" // for TLorentzVector
20#include "TMath.h" // for Sqrt
21#include "TVector3.h" // for TVector3
22#include "TMCProcess.h" // enum with process ids
23#include "TString.h"
24class TParticle;
25
26class ShipMCTrack : public TObject
27{
28
29 public:
30
31
34
35
37 ShipMCTrack(Int_t pdgCode, Int_t motherID, Double_t px, Double_t py,
38 Double_t pz, Double_t E, Double_t x, Double_t y, Double_t z,
39 Double_t t, Int_t nPoints, Double_t w);
40
42 ShipMCTrack(const ShipMCTrack& track);
43
44
46 ShipMCTrack(TParticle* particle);
47
48
50 virtual ~ShipMCTrack();
51
52
54 void Print(Int_t iTrack=0) const;
55
56
58 Int_t GetPdgCode() const { return fPdgCode; }
59 Int_t GetMotherId() const { return fMotherId; }
60 Double_t GetPx() const { return fPx; }
61 Double_t GetPy() const { return fPy; }
62 Double_t GetPz() const { return fPz; }
63 Double_t GetStartX() const { return fStartX; }
64 Double_t GetStartY() const { return fStartY; }
65 Double_t GetStartZ() const { return fStartZ; }
66 Double_t GetStartT() const { return fStartT; }
67 void SetProcID(Int_t i) { fProcID = i; }
68 Int_t GetProcID() const { return fProcID; }
69 TString GetProcName() const { return TMCProcessName[fProcID]; }
70 Double_t GetMass() const;
71 Double_t GetEnergy() const;
72 Double_t GetPt() const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
73 Double_t GetP() const { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz); }
74 Double_t GetRapidity() const;
75 void MultiplyWeight(Double_t w) {fW = fW*w;}
76 void SetWeight(Double_t w) {fW = w;}
77 Double_t GetWeight() const;
78 void GetMomentum(TVector3& momentum);
79 void Get4Momentum(TLorentzVector& momentum);
80 void GetStartVertex(TVector3& vertex);
81
82
84 Int_t GetNPoints(DetectorId detId) const;
85
86
88 void SetMotherId(Int_t id) { fMotherId = id; }
89 void SetNPoints(Int_t iDet, Int_t np);
90
91
92
93 private:
94
96 Int_t fPdgCode;
97
99 Int_t fMotherId;
100
102 Double32_t fPx, fPy, fPz, fM;
103
106
108 Double32_t fW;
109
111 Int_t fProcID;
112
128 Int_t fNPoints;
129
130
132
133};
134
135
136
137// ========== Inline functions ========================================
138
139
140inline void ShipMCTrack::GetMomentum(TVector3& momentum)
141{
142 momentum.SetXYZ(fPx,fPy,fPz);
143}
144
145
146inline void ShipMCTrack::Get4Momentum(TLorentzVector& momentum)
147{
148 momentum.SetXYZT(fPx,fPy,fPz,GetEnergy());
149}
150
151
152inline void ShipMCTrack::GetStartVertex(TVector3& vertex)
153{
154 vertex.SetXYZ(fStartX,fStartY,fStartZ);
155}
156
157
158
159
160
161#endif
DetectorId
Double_t GetWeight() const
Double_t GetStartZ() const
Definition ShipMCTrack.h:65
Double_t GetPt() const
Definition ShipMCTrack.h:72
void MultiplyWeight(Double_t w)
Definition ShipMCTrack.h:75
void Print(Int_t iTrack=0) const
Double32_t fPy
ClassDef(ShipMCTrack, 8)
Int_t fPdgCode
Definition ShipMCTrack.h:96
Double_t GetStartX() const
Definition ShipMCTrack.h:63
Double_t GetRapidity() const
Double_t GetMass() const
void SetWeight(Double_t w)
Definition ShipMCTrack.h:76
void GetMomentum(TVector3 &momentum)
Double_t GetStartY() const
Definition ShipMCTrack.h:64
Double32_t fStartX
Int_t GetProcID() const
Definition ShipMCTrack.h:68
Double32_t fW
Double32_t fStartZ
Double32_t fM
Double32_t fStartY
Double_t GetStartT() const
Definition ShipMCTrack.h:66
void GetStartVertex(TVector3 &vertex)
void SetProcID(Int_t i)
Definition ShipMCTrack.h:67
Double32_t fStartT
Double_t GetPy() const
Definition ShipMCTrack.h:61
Double_t GetEnergy() const
Double_t GetPx() const
Definition ShipMCTrack.h:60
Int_t GetNPoints(DetectorId detId) const
TString GetProcName() const
Definition ShipMCTrack.h:69
Int_t GetPdgCode() const
Definition ShipMCTrack.h:58
void Get4Momentum(TLorentzVector &momentum)
void SetNPoints(Int_t iDet, Int_t np)
Int_t GetMotherId() const
Definition ShipMCTrack.h:59
Int_t fMotherId
Definition ShipMCTrack.h:99
Double_t GetP() const
Definition ShipMCTrack.h:73
Double32_t fPz
virtual ~ShipMCTrack()
Double32_t fPx
void SetMotherId(Int_t id)
Definition ShipMCTrack.h:88
Double_t GetPz() const
Definition ShipMCTrack.h:62