EIC Software
Reference for
EIC
simulation and reconstruction software on GitHub
Home page
Related Pages
Modules
Namespaces
Classes
Files
External Links
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
ParticleMC.h
Go to the documentation of this file.
Or view
the newest version in sPHENIX GitHub for file ParticleMC.h
1
10
#ifndef INCLUDE_EICSMEAR_HADRONIC_PARTICLEMC_H_
11
#define INCLUDE_EICSMEAR_HADRONIC_PARTICLEMC_H_
12
13
#include <iostream>
14
15
#include <TLorentzVector.h>
16
#include <TVector3.h>
17
18
#include "
eicsmear/erhic/Pid.h
"
19
#include "
eicsmear/erhic/VirtualParticle.h
"
20
21
class
TMCParticle;
// ROOT Monte Carlo particle class
22
23
namespace
erhic {
24
namespace
hadronic {
25
30
class
ParticleMC
:
public
erhic::VirtualParticle
{
31
public
:
35
virtual
~ParticleMC
() { }
36
40
ParticleMC
();
41
45
explicit
ParticleMC
(
const
TMCParticle&);
46
51
ParticleMC
(
const
TLorentzVector&,
const
TVector3&,
int
,
int
,
int
);
52
56
virtual
erhic::Pid
Id
()
const
;
57
61
virtual
TLorentzVector
Get4Vector
()
const
;
62
66
virtual
Double_t
GetPx
()
const
;
67
71
virtual
Double_t
GetPy
()
const
;
72
76
virtual
Double_t
GetPz
()
const
;
77
81
virtual
Double_t
GetE
()
const
;
82
86
virtual
Double_t
GetP
()
const
;
87
91
virtual
Double_t
GetM
()
const
;
92
96
virtual
Double_t
GetPt
()
const
;
97
101
virtual
Double_t
GetTheta
()
const
;
102
106
virtual
Double_t
GetPhi
()
const
;
107
111
virtual
Double_t
GetRapidity
()
const
;
112
116
virtual
Double_t
GetEta
()
const
;
117
122
virtual
TVector3
GetVertex
()
const
;
123
128
virtual
UShort_t
GetStatus
()
const
;
129
134
virtual
UShort_t
GetParentIndex
()
const
;
135
140
virtual
Double_t
GetXFeynman
()
const
;
141
145
virtual
void
SetStatus
(UShort_t);
146
151
virtual
void
SetParentIndex
(UShort_t);
152
156
virtual
void
SetXFeynman
(
double
xf);
157
162
virtual
void
Set4Vector
(
const
TLorentzVector&);
163
167
virtual
void
SetVertex
(
const
TVector3&);
168
169
protected
:
170
UShort_t
KS
;
171
UShort_t
orig
;
172
Int_t
id
;
173
Double32_t
px
;
174
Double32_t
py
;
175
Double32_t
pz
;
176
Double32_t
E
;
177
Double32_t
p
;
178
Double32_t
m
;
179
Double32_t
pt
;
180
Double32_t
theta
;
181
Double32_t
phi
;
182
Double32_t
rapidity
;
183
Double32_t
eta
;
184
Double32_t
xFeynman
;
185
Double32_t
xv
;
186
Double32_t
yv
;
187
Double32_t
zv
;
188
189
ClassDef(
erhic::hadronic::ParticleMC
, 1)
190
};
191
192
inline
Pid
ParticleMC::Id
()
const
{
193
return
Pid
(
id
);
194
}
195
196
inline
TLorentzVector
ParticleMC::Get4Vector
()
const
{
197
return
TLorentzVector(
GetPx
(),
GetPy
(),
GetPz
(),
GetE
());
198
}
199
200
inline
Double_t
ParticleMC::GetPx
()
const
{
201
return
px
;
202
}
203
204
inline
Double_t
ParticleMC::GetPy
()
const
{
205
return
py
;
206
}
207
208
inline
Double_t
ParticleMC::GetPz
()
const
{
209
return
pz
;
210
}
211
212
inline
Double_t
ParticleMC::GetE
()
const
{
213
return
E
;
214
}
215
216
inline
Double_t
ParticleMC::GetP
()
const
{
217
return
p
;
218
}
219
220
inline
Double_t
ParticleMC::GetM
()
const
{
221
return
m
;
222
}
223
224
inline
Double_t
ParticleMC::GetPt
()
const
{
225
return
pt
;
226
}
227
228
inline
Double_t
ParticleMC::GetTheta
()
const
{
229
return
theta
;
230
}
231
232
inline
Double_t
ParticleMC::GetPhi
()
const
{
233
return
phi
;
234
}
235
236
inline
Double_t
ParticleMC::GetRapidity
()
const
{
237
return
rapidity
;
238
}
239
240
inline
Double_t
ParticleMC::GetEta
()
const
{
241
return
eta
;
242
}
243
244
inline
TVector3
ParticleMC::GetVertex
()
const
{
245
return
TVector3(
xv
,
yv
,
zv
);
246
}
247
248
inline
UShort_t
ParticleMC::GetStatus
()
const
{
249
return
KS
;
250
}
251
252
inline
UShort_t
ParticleMC::GetParentIndex
()
const
{
253
return
orig
;
254
}
255
256
inline
Double_t
ParticleMC::GetXFeynman
()
const
{
257
return
xFeynman
;
258
}
259
260
inline
void
ParticleMC::SetStatus
(UShort_t i) {
261
KS
= i;
262
}
263
264
inline
void
ParticleMC::SetXFeynman
(
double
xf) {
265
xFeynman
= xf;
266
}
267
268
inline
void
ParticleMC::SetVertex
(
const
TVector3&
v
) {
269
xv
= v.x();
270
yv
= v.y();
271
zv
= v.z();
272
}
273
274
}
// namespace hadronic
275
}
// namespace erhic
276
277
#endif // INCLUDE_EICSMEAR_HADRONIC_PARTICLEMC_H_
eic-smear
blob
master
include
eicsmear
hadronic
ParticleMC.h
Built by
Jin Huang
. updated:
Mon Jan 22 2024 12:43:32
using
1.8.2 with
EIC GitHub integration