-
Notifications
You must be signed in to change notification settings - Fork 496
Expand file tree
/
Copy pathGeneratorHepMC.cxx
More file actions
179 lines (143 loc) · 4.82 KB
/
GeneratorHepMC.cxx
File metadata and controls
179 lines (143 loc) · 4.82 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
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
// 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.
/// \author R+Preghenella - August 2017
#include "Generators/GeneratorHepMC.h"
#include "Generators/GeneratorHepMCParam.h"
#include "HepMC3/ReaderAscii.h"
#include "HepMC3/ReaderAsciiHepMC2.h"
#include "HepMC3/GenEvent.h"
#include "HepMC3/GenParticle.h"
#include "HepMC3/GenVertex.h"
#include "HepMC3/FourVector.h"
#include "TParticle.h"
#include "TSystem.h"
#include "FairLogger.h"
#include "FairPrimaryGenerator.h"
#include <cmath>
namespace o2
{
namespace eventgen
{
/*****************************************************************/
/*****************************************************************/
GeneratorHepMC::GeneratorHepMC()
: Generator("ALICEo2", "ALICEo2 HepMC Generator"), mStream(), mFileName(), mVersion(3), mReader(nullptr), mEvent(nullptr)
{
/** default constructor **/
mEvent = new HepMC3::GenEvent();
mInterface = reinterpret_cast<void*>(mEvent);
mInterfaceName = "hepmc";
}
/*****************************************************************/
GeneratorHepMC::GeneratorHepMC(const Char_t* name, const Char_t* title)
: Generator(name, title), mStream(), mFileName(), mVersion(3), mReader(nullptr), mEvent(nullptr)
{
/** constructor **/
mEvent = new HepMC3::GenEvent();
mInterface = reinterpret_cast<void*>(mEvent);
mInterfaceName = "hepmc";
}
/*****************************************************************/
GeneratorHepMC::~GeneratorHepMC()
{
/** default destructor **/
if (mStream.is_open())
mStream.close();
if (mReader) {
mReader->close();
delete mReader;
}
if (mEvent)
delete mEvent;
}
/*****************************************************************/
Bool_t GeneratorHepMC::generateEvent()
{
/** generate event **/
/** clear and read event **/
mEvent->clear();
mReader->read_event(*mEvent);
if (mReader->failed())
return kFALSE;
/** set units to desired output **/
mEvent->set_units(HepMC3::Units::GEV, HepMC3::Units::MM);
/** success **/
return kTRUE;
}
/*****************************************************************/
Bool_t GeneratorHepMC::importParticles()
{
/** import particles **/
/** loop over particles **/
auto particles = mEvent->particles();
for (int i = 0; i < particles.size(); ++i) {
/** get particle information **/
auto particle = particles.at(i);
auto pdg = particle->pid();
auto st = particle->status();
auto momentum = particle->momentum();
auto vertex = particle->production_vertex()->position();
auto parents = particle->parents(); // less efficient than via vertex
auto children = particle->children(); // less efficient than via vertex
/** get momentum information **/
auto px = momentum.x();
auto py = momentum.y();
auto pz = momentum.z();
auto et = momentum.t();
/** get vertex information **/
auto vx = vertex.x();
auto vy = vertex.y();
auto vz = vertex.z();
auto vt = vertex.t();
/** get mother information **/
auto m1 = parents.empty() ? -1 : parents.front()->id() - 1;
auto m2 = parents.empty() ? -1 : parents.back()->id() - 1;
/** get daughter information **/
auto d1 = children.empty() ? -1 : children.front()->id() - 1;
auto d2 = children.empty() ? -1 : children.back()->id() - 1;
/** add to particle vector **/
mParticles.push_back(TParticle(pdg, st, m1, m2, d1, d2, px, py, pz, et, vx, vy, vz, vt));
} /** end of loop over particles **/
/** success **/
return kTRUE;
}
/*****************************************************************/
Bool_t GeneratorHepMC::Init()
{
/** init **/
/** init base class **/
Generator::Init();
/** open file **/
std::string filename = gSystem->ExpandPathName(mFileName.c_str());
mStream.open(filename);
if (!mStream.is_open()) {
LOG(FATAL) << "Cannot open input file: " << filename << std::endl;
return kFALSE;
}
/** create reader according to HepMC version **/
switch (mVersion) {
case 2:
mStream.close();
mReader = new HepMC3::ReaderAsciiHepMC2(filename);
break;
case 3:
mReader = new HepMC3::ReaderAscii(mStream);
break;
default:
LOG(FATAL) << "Unsupported HepMC version: " << mVersion << std::endl;
return kFALSE;
}
/** success **/
return !mReader->failed();
}
/*****************************************************************/
/*****************************************************************/
} /* namespace eventgen */
} /* namespace o2 */