FairShip
Loading...
Searching...
No Matches
ShipParticle.cxx
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// -------------------------------------------------------------------------
6// ----- ShipParticle source file -----
7// -------------------------------------------------------------------------
8#include "ShipParticle.h"
9
10#include "FairLogger.h" // for FairLogger, etc
11#include "TLorentzVector.h" // for TLorentzVector
12#include "TVector3.h" // for TVector3
13
14// ----- Default constructor -------------------------------------------
16 : TObject(),
17 fPdgCode(0),
18 fStatusCode(0),
19 fMother{0, 0},
20 fDaughter{0, 0},
21 fPx(0.),
22 fPy(0.),
23 fPz(0.),
24 fE(0.),
25 fVx(0.),
26 fVy(0.),
27 fVz(0.),
28 fVt(0.),
29 fCovP{0., 0., 0., 0., 0., 0., 0., 0., 0., 0.},
30 fCovV{0., 0., 0., 0., 0., 0.},
31 fDoca(0.) {}
32// -------------------------------------------------------------------------
33
34// ----- Standard constructor ------------------------------------------
35ShipParticle::ShipParticle(Int_t pdg, Int_t status, Int_t mother1,
36 Int_t mother2, Int_t daughter1, Int_t daughter2,
37 Double_t px, Double_t py, Double_t pz, Double_t etot,
38 Double_t vx, Double_t vy, Double_t vz, Double_t time)
39 : TObject(),
40 fPdgCode(pdg),
41 fStatusCode(status),
42 fMother{mother1, mother2},
43 fDaughter{daughter1, daughter2},
44 fPx(px),
45 fPy(py),
46 fPz(pz),
47 fE(etot),
48 fVx(vx),
49 fVy(vy),
50 fVz(vz),
51 fVt(time),
52 fCovP{0., 0., 0., 0., 0., 0., 0., 0., 0., 0.},
53 fCovV{0., 0., 0., 0., 0., 0.},
54 fDoca(0.) {}
55// -------------------------------------------------------------------------
56
57// ----- Constructor with TLorentzVector -------------------------------
58ShipParticle::ShipParticle(Int_t pdg, Int_t status, Int_t mother1,
59 Int_t mother2, Int_t daughter1, Int_t daughter2,
60 const TLorentzVector& p, const TLorentzVector& v)
61 : TObject(),
62 fPdgCode(pdg),
63 fStatusCode(status),
64 fMother{mother1, mother2},
65 fDaughter{daughter1, daughter2},
66 fPx(p.Px()),
67 fPy(p.Py()),
68 fPz(p.Pz()),
69 fE(p.E()),
70 fVx(v.X()),
71 fVy(v.Y()),
72 fVz(v.Z()),
73 fVt(v.T()),
74 fCovP{0., 0., 0., 0., 0., 0., 0., 0., 0., 0.},
75 fCovV{0., 0., 0., 0., 0., 0.},
76 fDoca(0.) {}
77// -------------------------------------------------------------------------
78
79// ----- Destructor ----------------------------------------------------
81// -------------------------------------------------------------------------
82
83// ----- Public method SetCovP -----------------------------------------
84void ShipParticle::SetCovP(const Double_t* covElements) {
85 // Set 4x4 symmetric momentum covariance matrix from 10 unique elements
86 // Order: (0,0), (0,1), (0,2), (0,3), (1,1), (1,2), (1,3), (2,2), (2,3), (3,3)
87 fCovP[0] = covElements[0]; // (0,0)
88 fCovP[1] = covElements[1]; // (0,1)
89 fCovP[2] = covElements[2]; // (0,2)
90 fCovP[3] = covElements[3]; // (0,3)
91 fCovP[4] = covElements[4]; // (1,1)
92 fCovP[5] = covElements[5]; // (1,2)
93 fCovP[6] = covElements[6]; // (1,3)
94 fCovP[7] = covElements[7]; // (2,2)
95 fCovP[8] = covElements[8]; // (2,3)
96 fCovP[9] = covElements[9]; // (3,3)
97}
98// -------------------------------------------------------------------------
99
100// ----- Public method SetCovV -----------------------------------------
101void ShipParticle::SetCovV(const Double_t* covElements) {
102 // Set 3x3 symmetric vertex covariance matrix from 6 unique elements
103 // Order: (0,0), (0,1), (0,2), (1,1), (1,2), (2,2)
104 fCovV[0] = covElements[0]; // (0,0)
105 fCovV[1] = covElements[1]; // (0,1)
106 fCovV[2] = covElements[2]; // (0,2)
107 fCovV[3] = covElements[3]; // (1,1)
108 fCovV[4] = covElements[4]; // (1,2)
109 fCovV[5] = covElements[5]; // (2,2)
110}
111// -------------------------------------------------------------------------
112
113// ----- Public method GetMomentum -------------------------------------
114void ShipParticle::GetMomentum(TLorentzVector& momentum) const {
115 momentum.SetPxPyPzE(fPx, fPy, fPz, fE);
116}
117// -------------------------------------------------------------------------
118
119// ----- Public method GetVertex ---------------------------------------
120void ShipParticle::GetVertex(TVector3& vertex) const {
121 vertex.SetXYZ(fVx, fVy, fVz);
122}
123// -------------------------------------------------------------------------
124
125// ----- Public method GetMass -----------------------------------------
126Double_t ShipParticle::GetMass() const {
127 return TMath::Sqrt(fE * fE - fPx * fPx - fPy * fPy - fPz * fPz);
128}
129// -------------------------------------------------------------------------
130
131// ----- Public method Print -------------------------------------------
132void ShipParticle::Print(Int_t iTrack) const {
133 LOG(info) << "ShipParticle " << iTrack << ":";
134 LOG(info) << " PDG code: " << fPdgCode << ", Status: " << fStatusCode;
135 LOG(info) << " Mothers: " << fMother[0] << ", " << fMother[1];
136 LOG(info) << " Daughters: " << fDaughter[0] << ", " << fDaughter[1];
137 LOG(info) << " Momentum (" << fPx << ", " << fPy << ", " << fPz << ", " << fE
138 << ") GeV";
139 LOG(info) << " Vertex (" << fVx << ", " << fVy << ", " << fVz
140 << ") cm, t = " << fVt << " ns";
141 LOG(info) << " DOCA: " << fDoca << " cm";
142}
143// -------------------------------------------------------------------------
144
145const char* ShipParticle::GetName() const {
146 // Return particle name based on PDG code
147 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(fPdgCode);
148 if (particle) {
149 return particle->GetName();
150 }
151 return "Unknown";
152}
153// -------------------------------------------------------------------------
Definition: diagrams_e.h:4
TLorentzVector GetMomentum() const
Definition: ShipParticle.h:58
Double_t fE
Definition: ShipParticle.h:105
void SetCovP(const Double_t *covElements)
virtual ~ShipParticle()
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 fVy
Definition: ShipParticle.h:109
Int_t fStatusCode
Definition: ShipParticle.h:97
Double_t fPz
Definition: ShipParticle.h:104
Double_t GetMass() const
TVector3 GetVertex() const
Definition: ShipParticle.h:70
Int_t fMother[2]
Definition: ShipParticle.h:98
Double_t fDoca
Definition: ShipParticle.h:122
std::array< Double_t, 6 > fCovV
Definition: ShipParticle.h:119
Double_t fVx
Definition: ShipParticle.h:108
void SetCovV(const Double_t *covElements)
Int_t fPdgCode
Definition: ShipParticle.h:96
Double_t fPy
Definition: ShipParticle.h:103
Double_t fVt
Definition: ShipParticle.h:111
void Print(Int_t iTrack=0) const
const char * GetName() const override