FairShip
Loading...
Searching...
No Matches
ShipParticle.h
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-3.0-or-later
2// SPDX-FileCopyrightText: Copyright CERN for the benefit of the SHiP
3// Collaboration
4
5#ifndef SHIPDATA_SHIPPARTICLE_H_
6#define SHIPDATA_SHIPPARTICLE_H_
7
8#include <array> // for std::array
9
10#include "Rtypes.h" // for Double_t, Int_t, Double32_t, etc
11#include "ShipDetectorList.h" // for DetectorId
12#include "TDatabasePDG.h" // for TDatabasePDG
13#include "TLorentzVector.h" // for TLorentzVector
14#include "TMath.h" // for Sqrt
15#include "TObject.h" // for TObject
16#include "TVector3.h" // for TVector3
17
18class ShipParticle : public TObject {
19 public:
22 ShipParticle(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2,
23 Int_t daughter1, Int_t daughter2, Double_t px, Double_t py,
24 Double_t pz, Double_t etot, Double_t vx, Double_t vy,
25 Double_t vz, Double_t time);
26 ShipParticle(Int_t pdg, Int_t status, Int_t mother1, Int_t mother2,
27 Int_t daughter1, Int_t daughter2, const TLorentzVector& p,
28 const TLorentzVector& v);
30 virtual ~ShipParticle();
31
33 ShipParticle(const ShipParticle& particle) = default;
34 ShipParticle& operator=(const ShipParticle& particle) = default;
35
37 void Print(Int_t iTrack = 0) const;
38
40 Int_t GetPdgCode() const { return fPdgCode; }
41 Int_t GetStatusCode() const { return fStatusCode; }
42 Int_t GetMother(Int_t i) const { return fMother[i]; }
43 Int_t GetFirstMother() const { return fMother[0]; }
44 Int_t GetSecondMother() const { return fMother[1]; }
45 Int_t GetDaughter(Int_t i) const { return fDaughter[i]; }
46 Int_t GetFirstDaughter() const { return fDaughter[0]; }
47 Int_t GetLastDaughter() const { return fDaughter[1]; }
48 const char* GetName() const override;
49
51 Double_t Px() const { return fPx; }
52 Double_t Py() const { return fPy; }
53 Double_t Pz() const { return fPz; }
54 Double_t P() const { return TMath::Sqrt(fPx * fPx + fPy * fPy + fPz * fPz); }
55 Double_t Pt() const { return TMath::Sqrt(fPx * fPx + fPy * fPy); }
56 Double_t Energy() const { return fE; }
57 void GetMomentum(TLorentzVector& momentum) const;
58 TLorentzVector GetMomentum() const {
59 return TLorentzVector(fPx, fPy, fPz, fE);
60 }
61 void Momentum(TLorentzVector& momentum) const { GetMomentum(momentum); }
62 Double_t GetMass() const;
63
65 Double_t Vx() const { return fVx; }
66 Double_t Vy() const { return fVy; }
67 Double_t Vz() const { return fVz; }
68 Double_t T() const { return fVt; }
69 void GetVertex(TVector3& vertex) const;
70 TVector3 GetVertex() const { return TVector3(fVx, fVy, fVz); }
71 void ProductionVertex(TLorentzVector& vertex) const {
72 vertex.SetXYZT(fVx, fVy, fVz, fVt);
73 }
74
76 // Get momentum covariance matrix elements (4x4 symmetric, 10 unique values)
77 const std::array<Double_t, 10>& GetCovP() const { return fCovP; }
78 void SetCovP(const Double_t* covElements);
79 void SetCovP(const std::array<Double_t, 10>& covElements) {
80 fCovP = covElements;
81 }
82
83 // Get vertex covariance matrix elements (3x3 symmetric, 6 unique values)
84 const std::array<Double_t, 6>& GetCovV() const { return fCovV; }
85 void SetCovV(const Double_t* covElements);
86 void SetCovV(const std::array<Double_t, 6>& covElements) {
87 fCovV = covElements;
88 }
89
91 Double_t GetDoca() const { return fDoca; }
92 void SetDoca(Double_t x) { fDoca = x; }
93
94 private:
95 // Particle identification and status
96 Int_t fPdgCode; // PDG code of the particle
97 Int_t fStatusCode; // generation status code
98 Int_t fMother[2]; // Indices of the mother particles
99 Int_t fDaughter[2]; // Indices of the daughter particles
100
101 // 4-momentum (px, py, pz, E)
102 Double_t fPx; // x component of momentum [GeV]
103 Double_t fPy; // y component of momentum [GeV]
104 Double_t fPz; // z component of momentum [GeV]
105 Double_t fE; // Energy [GeV]
106
107 // 4-vertex (vx, vy, vz, t)
108 Double_t fVx; // x of production vertex [cm]
109 Double_t fVy; // y of production vertex [cm]
110 Double_t fVz; // z of production vertex [cm]
111 Double_t fVt; // t of production vertex [ns]
112
113 // Covariance matrices stored as upper triangular elements
114 // fCovP: 4x4 symmetric momentum covariance (10 elements: 0,0; 0,1; 0,2; 0,3;
115 // 1,1; 1,2; 1,3; 2,2; 2,3; 3,3)
116 std::array<Double_t, 10> fCovP;
117 // fCovV: 3x3 symmetric vertex covariance (6 elements: 0,0; 0,1; 0,2; 1,1;
118 // 1,2; 2,2)
119 std::array<Double_t, 6> fCovV;
120
121 // Distance of closest approach
122 Double_t fDoca;
123
125};
126
127#endif // SHIPDATA_SHIPPARTICLE_H_
void SetDoca(Double_t x)
Definition: ShipParticle.h:92
Double_t Vx() const
Definition: ShipParticle.h:65
const std::array< Double_t, 6 > & GetCovV() const
Definition: ShipParticle.h:84
TLorentzVector GetMomentum() const
Definition: ShipParticle.h:58
void Momentum(TLorentzVector &momentum) const
Definition: ShipParticle.h:61
Double_t fE
Definition: ShipParticle.h:105
Int_t GetFirstDaughter() const
Definition: ShipParticle.h:46
Int_t GetDaughter(Int_t i) const
Definition: ShipParticle.h:45
void SetCovP(const Double_t *covElements)
virtual ~ShipParticle()
ShipParticle & operator=(const ShipParticle &particle)=default
Double_t Py() const
Definition: ShipParticle.h:52
Double_t Vz() const
Definition: ShipParticle.h:67
Double_t fPx
Definition: ShipParticle.h:102
Int_t fDaughter[2]
Definition: ShipParticle.h:99
std::array< Double_t, 10 > fCovP
Definition: ShipParticle.h:116
Double_t fVz
Definition: ShipParticle.h:110
Double_t GetDoca() const
Definition: ShipParticle.h:91
Double_t P() const
Definition: ShipParticle.h:54
Double_t fVy
Definition: ShipParticle.h:109
Int_t fStatusCode
Definition: ShipParticle.h:97
Int_t GetPdgCode() const
Definition: ShipParticle.h:40
Double_t fPz
Definition: ShipParticle.h:104
Double_t Pt() const
Definition: ShipParticle.h:55
Int_t GetLastDaughter() const
Definition: ShipParticle.h:47
Int_t GetSecondMother() const
Definition: ShipParticle.h:44
Int_t GetMother(Int_t i) const
Definition: ShipParticle.h:42
Double_t Px() const
Definition: ShipParticle.h:51
Double_t GetMass() const
TVector3 GetVertex() const
Definition: ShipParticle.h:70
const std::array< Double_t, 10 > & GetCovP() const
Definition: ShipParticle.h:77
void ProductionVertex(TLorentzVector &vertex) const
Definition: ShipParticle.h:71
Int_t fMother[2]
Definition: ShipParticle.h:98
void SetCovP(const std::array< Double_t, 10 > &covElements)
Definition: ShipParticle.h:79
Double_t fDoca
Definition: ShipParticle.h:122
Double_t T() const
Definition: ShipParticle.h:68
Int_t GetStatusCode() const
Definition: ShipParticle.h:41
std::array< Double_t, 6 > fCovV
Definition: ShipParticle.h:119
Double_t fVx
Definition: ShipParticle.h:108
Double_t Energy() const
Definition: ShipParticle.h:56
void SetCovV(const Double_t *covElements)
Int_t fPdgCode
Definition: ShipParticle.h:96
Int_t GetFirstMother() const
Definition: ShipParticle.h:43
Double_t fPy
Definition: ShipParticle.h:103
ShipParticle(const ShipParticle &particle)=default
Double_t Pz() const
Definition: ShipParticle.h:53
Double_t fVt
Definition: ShipParticle.h:111
ClassDefOverride(ShipParticle, 3)
void SetCovV(const std::array< Double_t, 6 > &covElements)
Definition: ShipParticle.h:86
void Print(Int_t iTrack=0) const
Double_t Vy() const
Definition: ShipParticle.h:66
const char * GetName() const override