forked from mcoquet642/AliceO2
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathGeneratorTGenerator.cxx
More file actions
133 lines (106 loc) · 3.77 KB
/
GeneratorTGenerator.cxx
File metadata and controls
133 lines (106 loc) · 3.77 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
// Copyright CERN and copyright holders of ALICE O2. This software is
// distributed under the terms of the GNU General Public License v3 (GPL
// Version 3), copied verbatim in the file "COPYING".
//
// See http://alice-o2.web.cern.ch/license for full licensing information.
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
#include "Generators/GeneratorTGenerator.h"
#include "FairLogger.h"
#include "FairPrimaryGenerator.h"
#include "TGenerator.h"
#include "TClonesArray.h"
#include "TParticle.h"
namespace o2
{
namespace eventgen
{
/*****************************************************************/
/*****************************************************************/
GeneratorTGenerator::GeneratorTGenerator() : Generator("ALICEo2", "ALICEo2 TGenerator Generator"),
mTGenerator(nullptr),
mParticles(nullptr)
{
/** default constructor **/
}
/*****************************************************************/
GeneratorTGenerator::GeneratorTGenerator(const Char_t* name, const Char_t* title) : Generator(name, title),
mTGenerator(nullptr),
mParticles(nullptr)
{
/** constructor **/
}
/*****************************************************************/
GeneratorTGenerator::~GeneratorTGenerator()
{
/** default destructor **/
if (mParticles)
delete mParticles;
}
/*****************************************************************/
Bool_t
GeneratorTGenerator::generateEvent()
{
/** generate event **/
mTGenerator->GenerateEvent();
mTGenerator->ImportParticles(mParticles, "Final");
/** success **/
return kTRUE;
}
/*****************************************************************/
Bool_t
GeneratorTGenerator::boostEvent(Double_t boost)
{
/** boost event **/
LOG(WARNING) << "Boost not implemented yet" << std::endl;
return kTRUE;
}
/*****************************************************************/
Bool_t
GeneratorTGenerator::addTracks(FairPrimaryGenerator* primGen) const
{
/** add tracks **/
/* loop over particles */
Int_t nParticles = mParticles->GetEntries();
TParticle* particle = nullptr;
for (Int_t iparticle = 0; iparticle < nParticles; iparticle++) {
particle = (TParticle*)mParticles->At(iparticle);
if (!particle)
continue;
primGen->AddTrack(particle->GetPdgCode(),
particle->Px() * mMomentumUnit,
particle->Py() * mMomentumUnit,
particle->Pz() * mMomentumUnit,
particle->Vx() * mPositionUnit,
particle->Vy() * mPositionUnit,
particle->Vz() * mPositionUnit,
particle->GetMother(0),
particle->GetStatusCode() == 1,
particle->Energy() * mEnergyUnit,
particle->T() * mTimeUnit,
particle->GetWeight());
}
/** success **/
return kTRUE;
}
/*****************************************************************/
Bool_t
GeneratorTGenerator::Init()
{
/** init **/
if (!mTGenerator) {
LOG(FATAL) << "No TGenerator inteface assigned" << std::endl;
return kFALSE;
}
/** array of generated particles **/
mParticles = new TClonesArray("TParticle");
mParticles->SetOwner(kTRUE);
/** success **/
return kTRUE;
}
/*****************************************************************/
/*****************************************************************/
} /* namespace eventgen */
} /* namespace o2 */