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_ERHIC_PARTICLEMC_H_
11
#define INCLUDE_EICSMEAR_ERHIC_PARTICLEMC_H_
12
13
#include <string>
14
15
#include <TLorentzVector.h>
16
#include <TRef.h>
17
#include <TVector3.h>
18
19
#include "
eicsmear/erhic/Pid.h
"
20
#include "
eicsmear/erhic/VirtualParticle.h
"
21
22
namespace
erhic {
23
24
class
ParticleMCeA
:
public
TObject
{
25
public
:
26
explicit
ParticleMCeA
() {};
27
virtual
~ParticleMCeA
() {};
28
29
// Let them all be public;
30
//added by liang to include add on particle data structure
31
Int_t
charge
;
32
Int_t
massNum
;
33
Int_t
NoBam
;
34
35
36
ClassDef(
ParticleMCeA
, 1)
37
};
38
39
class
EventMC
;
43
class
ParticleMCbase
:
public
VirtualParticle
{
44
public
:
45
explicit
ParticleMCbase
();
49
virtual
~ParticleMCbase
() {};
50
57
virtual
void
Print
(Option_t* =
""
)
const
;
58
62
virtual
UInt_t
GetIndex
()
const
;
63
69
virtual
UShort_t
GetStatus
()
const
;
70
74
virtual
UShort_t
GetParentIndex
()
const
;
75
virtual
UShort_t
GetParentIndex1
()
const
;
76
81
virtual
UShort_t
GetChild1Index
()
const
;
82
87
virtual
UShort_t
GetChildNIndex
()
const
;
88
93
virtual
UInt_t
GetNChildren
()
const
;
94
98
virtual
Double_t
GetPx
()
const
;
99
103
virtual
Double_t
GetPy
()
const
;
104
108
virtual
Double_t
GetPz
()
const
;
109
113
virtual
Double_t
GetM
()
const
;
114
118
virtual
Double_t
GetPt
()
const
;
119
124
virtual
TVector3
GetVertex
()
const
;
125
130
virtual
Pid
GetParentId
()
const
;
131
135
virtual
Double_t
GetP
()
const
;
136
140
virtual
Double_t
GetTheta
()
const
;
141
145
virtual
Double_t
GetPhi
()
const
;
146
150
virtual
Double_t
GetRapidity
()
const
;
151
155
virtual
Double_t
GetEta
()
const
;
156
161
virtual
Double_t
GetZ
()
const
;
162
167
virtual
Double_t
GetXFeynman
()
const
;
168
174
virtual
Double_t
GetThetaVsGamma
()
const
;
175
180
virtual
Double_t
GetPtVsGamma
()
const
;
181
187
const
EventMC
*
GetEvent
()
const
;
188
192
void
SetEvent
(
EventMC
*
event
);
193
197
virtual
TLorentzVector
Get4Vector
()
const
;
198
202
virtual
TLorentzVector
PxPyPzE
()
const
{
return
Get4Vector
(); }
203
221
virtual
TLorentzVector
Get4VectorInHadronBosonFrame
()
const
;
222
226
virtual
Double_t
GetE
()
const
;
227
228
virtual
void
SetE
(Double_t);
229
230
virtual
void
SetM
(Double_t);
231
232
virtual
void
SetP
(Double_t);
233
234
virtual
void
SetPt
(Double_t);
235
236
virtual
void
SetPz
(Double_t);
237
238
virtual
void
SetPhi
(Double_t);
239
240
virtual
void
SetTheta
(Double_t);
241
242
virtual
void
SetStatus
(UShort_t);
243
247
virtual
Pid
Id
()
const
;
248
252
virtual
Pid
GetPdgCode
()
const
{
return
Id
(); }
253
267
virtual
void
ComputeDerivedQuantities
();
268
283
virtual
void
ComputeEventDependentQuantities
(
EventMC
&);
284
290
virtual
void
SetIndex
(
int
i) {
I
= i; }
291
294
virtual
void
SetStatus
(
int
i) {
KS
= i; }
295
299
virtual
void
SetId
(
int
i) {
id
= i; }
300
303
virtual
void
SetParentIndex
(
int
i) {
orig
= i; }
304
307
virtual
void
SetParentIndex1
(
int
i) {
orig1
= i; }
308
311
virtual
void
SetChild1Index
(
int
i) {
daughter
= i; }
312
315
virtual
void
SetChildNIndex
(
int
i) {
ldaughter
= i; }
316
321
virtual
void
Set4Vector
(
const
TLorentzVector&);
322
326
virtual
void
SetVertex
(
const
TVector3&);
327
331
virtual
void
SetParentId
(
int
i) {
parentId
= i; }
332
333
// protected:
334
335
UShort_t
I
;
336
//UShort_t KS; ///< Particle status code: see PYTHIA manual
337
//modified by liang to include negative KS values
338
Int_t
KS
;
339
Int_t
id
;
340
// Looks strange (see also access methods); keep like this for sort of
341
// backward compatibility;
342
UShort_t
orig
;
343
UShort_t
orig1
;
344
UShort_t
daughter
;
345
UShort_t
ldaughter
;
346
347
// AYK, 2017/04/02: changed {px,py,pz} & {xv,yv,zv} variable types
348
// from Double32_t (float) to Double_t (double);
349
Double_t
px
;
350
Double_t
py
;
351
Double_t
pz
;
352
353
// Reducing file size.
354
// The following variables could be removed, at the cost of increased
355
// running time to calculate them on-the-fly:
356
// p - compute from px/y/z
357
// m - accessible via PID with TDatabasePDG
358
// E - compute from p and m
359
// parentId - accessible via event record and parent index IF a
360
// persistent pointer to the containing event can be stored
361
// daughter - IF daughter particles always start immediately after this one
362
// pt - from px, py
363
// theta - from pt, pz
364
// phi - from px, py
365
// rapidity - from E, pz
366
// eta - from theta
367
// xFeynman - if event is accessible
368
// thetaGamma - duplicated in pHadronBoson
369
// ptVsGamma - duplicated in pHadronBoson
370
// phiPrf - duplicated in pHadronBoson
371
372
Double32_t
E
;
373
Double32_t
m
;
374
Double32_t
pt
;
375
Double_t
xv
;
376
Double_t
yv
;
377
Double_t
zv
;
378
379
// Can be deprecated if parent particle itself is available
380
Int_t
parentId
;
381
382
Double32_t
p
;
383
Double32_t
theta
;
384
Double32_t
phi
;
385
Double32_t
rapidity
;
386
Double32_t
eta
;
387
Double32_t
z
;
388
389
Double32_t
xFeynman
;
390
Double32_t
thetaGamma
;
391
392
Double32_t
ptVsGamma
;
393
394
Double32_t
thetaGammaHCM
;
395
396
Double32_t
ptVsGammaHCM
;
397
398
Double32_t
phiPrf
;
399
400
401
TRef
event
;
402
403
404
ClassDef(
ParticleMCbase
, 5)
405
};
406
407
// The unfortunate consequence of adding 'eA' pointer to the ParticleMC
408
// class is that default copy ctor does not work; the easiest workaround
409
// (unless somebody would like to copy all the other variables by hand)
410
// is to offload these variables to an intermediate base class;
411
class
ParticleMC
:
public
ParticleMCbase
{
412
public
:
419
explicit
ParticleMC
():
eA
(0) {};
420
explicit
ParticleMC
(
const
std::string&,
bool
eAflag);
421
// So this is the evil copy ctor, which is damn simple now;
422
ParticleMC
(
const
ParticleMC
&src):
ParticleMCbase
(src) {
423
eA
= src.
eA
?
new
ParticleMCeA
(*src.
eA
) : 0;
//orig1 = src.orig1;
424
};
425
426
~ParticleMC
();
427
435
virtual
const
ParticleMC
*
GetParent
()
const
;
436
447
virtual
const
ParticleMC
*
GetChild
(UShort_t)
const
;
448
455
virtual
Bool_t
HasChild
(Int_t)
const
;
456
457
// FIXME: let it be public?; NB: TRef (like for ParticleMCbase.event) is
458
// not needed here, since it is just a POD instance with a single reference;
459
ParticleMCeA
*
eA
;
460
//UShort_t orig1; ///< I of parent particle1
461
462
ClassDef(
ParticleMC
, 2)
463
};
464
465
inline
UInt_t
ParticleMCbase::GetIndex
()
const
{
466
return
I
;
467
}
468
469
inline
UShort_t
ParticleMCbase::GetStatus
()
const
{
470
return
KS
;
471
}
472
473
inline
UShort_t
ParticleMCbase::GetParentIndex
()
const
{
474
return
orig
;
475
}
476
inline
UShort_t
ParticleMCbase::GetParentIndex1
()
const
{
477
return
orig1
;
478
}
479
480
inline
UShort_t
ParticleMCbase::GetChild1Index
()
const
{
481
return
daughter
;
482
}
483
484
inline
UShort_t
ParticleMCbase::GetChildNIndex
()
const
{
485
return
ldaughter
;
486
}
487
488
inline
Double_t
ParticleMCbase::GetPx
()
const
{
489
return
px
;
490
}
491
492
inline
Double_t
ParticleMCbase::GetPy
()
const
{
493
return
py
;
494
}
495
496
inline
Double_t
ParticleMCbase::GetPz
()
const
{
497
return
pz
;
498
}
499
500
inline
Double_t
ParticleMCbase::GetM
()
const
{
501
return
m
;
502
}
503
504
inline
Double_t
ParticleMCbase::GetPt
()
const
{
505
return
pt
;
506
}
507
508
inline
TVector3
ParticleMCbase::GetVertex
()
const
{
509
return
TVector3(
xv
,
yv
,
zv
);
510
}
511
512
inline
erhic::Pid
ParticleMCbase::GetParentId
()
const
{
513
return
erhic::Pid
(
parentId
);
514
}
515
516
inline
Double_t
ParticleMCbase::GetP
()
const
{
517
return
p
;
518
}
519
520
inline
Double_t
ParticleMCbase::GetTheta
()
const
{
521
return
theta
;
522
}
523
524
inline
Double_t
ParticleMCbase::GetPhi
()
const
{
525
return
phi
;
526
}
527
528
inline
Double_t
ParticleMCbase::GetRapidity
()
const
{
529
return
rapidity
;
530
}
531
532
inline
Double_t
ParticleMCbase::GetEta
()
const
{
533
return
eta
;
534
}
535
536
inline
Double_t
ParticleMCbase::GetZ
()
const
{
537
return
z
;
538
}
539
540
inline
Double_t
ParticleMCbase::GetXFeynman
()
const
{
541
return
xFeynman
;
542
}
543
544
inline
Double_t
ParticleMCbase::GetThetaVsGamma
()
const
{
545
return
thetaGamma
;
546
}
547
548
inline
Double_t
ParticleMCbase::GetPtVsGamma
()
const
{
549
return
ptVsGamma
;
550
}
551
552
inline
Pid
ParticleMCbase::Id
()
const
{
553
return
Pid
(
id
);
554
}
555
556
inline
UInt_t
ParticleMCbase::GetNChildren
()
const
{
557
if
(0 ==
daughter
)
return
0;
558
if
(0 ==
ldaughter
)
return
1;
559
return
ldaughter
-
daughter
+ 1;
560
}
561
562
inline
Double_t
ParticleMCbase::GetE
()
const
{
563
return
E
;
564
}
565
566
inline
void
ParticleMCbase::SetE
(Double_t
e
) {
567
E
=
e
;
568
}
569
570
inline
void
ParticleMCbase::SetM
(Double_t
mass
) {
571
m
=
mass
;
572
}
573
574
inline
void
ParticleMCbase::SetP
(Double_t
momentum
) {
575
p
=
momentum
;
576
}
577
578
inline
void
ParticleMCbase::SetPt
(Double_t
momentum
) {
579
pt
=
momentum
;
580
}
581
582
inline
void
ParticleMCbase::SetPz
(Double_t
momentum
) {
583
pz
=
momentum
;
584
}
585
586
inline
void
ParticleMCbase::SetPhi
(Double_t
value
) {
587
phi
=
value
;
588
}
589
590
inline
void
ParticleMCbase::SetTheta
(Double_t
value
) {
591
theta
=
value
;
592
}
593
594
inline
void
ParticleMCbase::SetStatus
(UShort_t status) {
595
KS
= status;
596
}
597
598
}
// namespace erhic
599
600
#endif // INCLUDE_EICSMEAR_ERHIC_PARTICLEMC_H_
eic-smear
blob
master
include
eicsmear
erhic
ParticleMC.h
Built by
Jin Huang
. updated:
Mon Jan 22 2024 12:43:32
using
1.8.2 with
EIC GitHub integration