-
Notifications
You must be signed in to change notification settings - Fork 496
Expand file tree
/
Copy pathPrimaryGenerator.cxx
More file actions
391 lines (329 loc) · 12.4 KB
/
PrimaryGenerator.cxx
File metadata and controls
391 lines (329 loc) · 12.4 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
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// 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 - June 2017
#include "Generators/PrimaryGenerator.h"
#include <Generators/PrimaryGeneratorParam.h>
#include "Generators/Generator.h"
#include "SimConfig/InteractionDiamondParam.h"
#include "SimulationDataFormat/MCEventHeader.h"
#include "SimulationDataFormat/MCGenProperties.h"
#include "DataFormatsCalibration/MeanVertexObject.h"
#include "DetectorsBase/Stack.h"
#include <fairlogger/Logger.h>
#include "FairGenericStack.h"
#include "TFile.h"
#include "TTree.h"
#include "TDatabasePDG.h"
#include "TVirtualMC.h"
using o2::dataformats::MCEventHeader;
namespace o2
{
namespace eventgen
{
/*****************************************************************/
PrimaryGenerator::~PrimaryGenerator()
{
/** destructor **/
LOG(info) << "Destructing PrimaryGenerator";
if (mEmbedFile && mEmbedFile->IsOpen()) {
mEmbedFile->Close();
delete mEmbedFile;
}
if (mEmbedEvent) {
delete mEmbedEvent;
}
}
/*****************************************************************/
Bool_t PrimaryGenerator::Init()
{
/** init **/
LOG(info) << "Initialising primary generator";
// set generator ID and description
auto& params = PrimaryGeneratorParam::Instance();
setGeneratorId(params.id);
setGeneratorDescription(params.description);
/** embedding **/
if (mEmbedTree) {
LOG(info) << "Embedding into: " << mEmbedFile->GetName()
<< " (" << mEmbedEntries << " events)";
return FairPrimaryGenerator::Init();
}
/** base class init **/
return FairPrimaryGenerator::Init();
}
/*****************************************************************/
Bool_t PrimaryGenerator::GenerateEvent(FairGenericStack* pStack)
{
/** generate event **/
/** normal generation if no embedding **/
if (!mEmbedTree || mEmbedIndex < 0) {
fixInteractionVertex(); // <-- always fixes vertex outside of FairROOT
auto ret = FairPrimaryGenerator::GenerateEvent(pStack);
if (ret) {
setGeneratorInformation();
}
return ret;
}
/** this is for embedding **/
if (mEmbedIndex >= 0) {
/** setup interaction vertex **/
mEmbedTree->GetEntry(mEmbedIndex);
setInteractionVertex(mEmbedEvent);
/** notify event generators **/
auto genList = GetListOfGenerators();
for (int igen = 0; igen < genList->GetEntries(); ++igen) {
auto o2gen = dynamic_cast<Generator*>(genList->At(igen));
if (o2gen) {
o2gen->notifyEmbedding(mEmbedEvent);
}
}
}
/** generate event **/
if (!FairPrimaryGenerator::GenerateEvent(pStack)) {
return kFALSE;
}
/** add embedding info to event header **/
auto o2event = dynamic_cast<MCEventHeader*>(fEvent);
if (o2event) {
o2event->setEmbeddingFileName(mEmbedFile->GetName());
o2event->setEmbeddingEventIndex(mEmbedIndex);
}
setGeneratorInformation();
/** increment embedding counter **/
mEmbedIndex++;
mEmbedIndex %= mEmbedEntries;
/** success **/
return kTRUE;
}
/*****************************************************************/
void PrimaryGenerator::AddTrack(Int_t pdgid, Double_t px, Double_t py, Double_t pz,
Double_t vx, Double_t vy, Double_t vz,
Int_t mother1, Int_t mother2,
Int_t daughter1, Int_t daughter2,
Bool_t wanttracking,
Double_t e, Double_t tof,
Double_t weight, TMCProcess proc, Int_t generatorStatus)
{
/** add track **/
// check the status encoding
if (!mcgenstatus::isEncoded(generatorStatus) && proc == TMCProcess::kPPrimary) {
LOG(fatal) << "Generator status " << generatorStatus << " of particle is not encoded properly.";
}
/** add event vertex to track vertex **/
vx += fVertex.X();
vy += fVertex.Y();
vz += fVertex.Z();
/** check if particle to be tracked exists in PDG database **/
auto particlePDG = TDatabasePDG::Instance()->GetParticle(pdgid);
if (wanttracking && !particlePDG) {
LOG(warn) << "Particle to be tracked is not defined in PDG: pdg = " << pdgid << " (disabling tracking)";
wanttracking = false;
}
/** set all other parameters required by PushTrack **/
Int_t doTracking = 0; // Go to tracking
if (fdoTracking && wanttracking) {
doTracking = 1;
}
Int_t dummyparent = -1; // Primary particle (now the value is -1 by default)
Double_t polx = 0.; // Polarisation
Double_t poly = 0.;
Double_t polz = 0.;
Int_t ntr = 0; // Track number; to be filled by the stack
Int_t status = generatorStatus; // Generation status
// correct for tracks which are in list before generator is called
if (mother1 != -1) {
mother1 += fMCIndexOffset;
}
if (mother2 != -1) {
mother2 += fMCIndexOffset;
}
if (daughter1 != -1) {
daughter1 += fMCIndexOffset;
}
if (daughter2 != -1) {
daughter2 += fMCIndexOffset;
}
/** if it is a K0/antiK0 to be tracked, convert it into K0s/K0L.
NOTE: we could think of pushing the K0/antiK0 without tracking first
and then push she K0s/K0L for tracking.
In this way we would properly keep track of this conversion,
but there is the risk of messing up with the indices, so this
is not done for the time being.
**/
if (abs(pdgid) == 311 && doTracking) {
LOG(warn) << "K0/antiK0 requested for tracking: converting into K0s/K0L";
pdgid = gRandom->Uniform() < 0.5 ? 310 : 130;
}
/** compute particle energy if negative **/
if (e < 0) {
double mass = particlePDG ? particlePDG->Mass() : 0.;
e = std::sqrt(mass * mass + px * px + py * py + pz * pz);
}
/** add track to stack **/
auto stack = dynamic_cast<o2::data::Stack*>(fStack);
if (!stack) {
LOG(fatal) << "Stack must be an o2::data:Stack";
return; // must be the o2 stack
}
stack->PushTrack(doTracking, mother1, pdgid, px, py, pz,
e, vx, vy, vz, tof, polx, poly, polz, TMCProcess::kPPrimary, ntr,
weight, status, mother2, daughter1, daughter2, proc);
fNTracks++;
}
/*****************************************************************/
void PrimaryGenerator::AddTrack(Int_t pdgid, Double_t px, Double_t py,
Double_t pz, Double_t vx, Double_t vy,
Double_t vz, Int_t parent, Bool_t wanttracking,
Double_t e, Double_t tof, Double_t weight, TMCProcess proc)
{
// Do this to encode status code correctly. In FairRoot's PrimaryGenerator, this is simply one number that is set to 0.
// So we basically do the same.
// Assuming that this is treated as the HepMC status code down the line (as it used to be).
AddTrack(pdgid, px, py, pz, vx, vy, vz, parent, -1, -1, -1, wanttracking,
e, tof, weight, proc, mcgenstatus::MCGenStatusEncoding(0, 0).fullEncoding);
}
/*****************************************************************/
void PrimaryGenerator::setInteractionVertex(const MCEventHeader* event)
{
if (!mApplyVertex) {
return;
}
/** set interaction vertex **/
Double_t xyz[3] = {event->GetX(), event->GetY(), event->GetZ()};
SetBeam(xyz[0], xyz[1], 0., 0.);
SetTarget(xyz[2], 0.);
SmearVertexXY(false);
SmearVertexZ(false);
SmearGausVertexXY(false);
SmearGausVertexZ(false);
}
/*****************************************************************/
void PrimaryGenerator::setExternalVertexForNextEvent(double x, double y, double z)
{
if (!mApplyVertex) {
return;
}
mExternalVertexX = x;
mExternalVertexY = y;
mExternalVertexZ = z;
mHaveExternalVertex = true;
}
/*****************************************************************/
void PrimaryGenerator::setVertexMode(o2::conf::VertexMode const& mode, o2::dataformats::MeanVertexObject const* v)
{
mVertexMode = mode;
if (mode == o2::conf::VertexMode::kCCDB) {
if (!v) {
LOG(fatal) << "A valid MeanVertexObject needs to be passed with option o2::conf::VertexMode::kCCDB";
}
mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(*v)));
LOG(info) << "The mean vertex is set to :";
mMeanVertex->print();
}
if (mVertexMode == o2::conf::VertexMode::kNoVertex) {
setApplyVertex(false);
LOG(info) << "Disabling vertexing";
mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(0, 0, 0, 0, 0, 0, 0, 0)));
LOG(info) << "The mean vertex is set to :";
mMeanVertex->print();
}
}
/*****************************************************************/
void PrimaryGenerator::fixInteractionVertex()
{
if (!mApplyVertex) {
SetBeam(0., 0., 0., 0.);
SetTarget(0., 0.);
return;
}
// if someone gave vertex from outside; we will take it
if (mHaveExternalVertex) {
SetBeam(mExternalVertexX, mExternalVertexY, 0., 0.);
SetTarget(mExternalVertexZ, 0.);
mHaveExternalVertex = false; // the vertex is now consumed
return;
}
// sampling a vertex and fixing for next event; no smearing will be done
// inside FairPrimaryGenerator;
SmearVertexXY(false);
SmearVertexZ(false);
SmearGausVertexXY(false);
SmearGausVertexZ(false);
// we use the mMeanVertexObject if initialized (initialize first)
if (mMeanVertex.get() == nullptr) {
if (mVertexMode == o2::conf::VertexMode::kDiamondParam) {
auto const& param = InteractionDiamondParam::Instance();
const auto& xyz = param.position;
const auto& sigma = param.width;
mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(xyz[0], xyz[1], xyz[2], sigma[0], sigma[1], sigma[2], param.slopeX, param.slopeY)));
}
if (mVertexMode == o2::conf::VertexMode::kNoVertex) {
mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(0, 0, 0, 0, 0, 0, 0, 0)));
}
if (mVertexMode == o2::conf::VertexMode::kCCDB) {
// fatal.. then the object should have been passed with setting
LOG(fatal) << "MeanVertexObject is null ... but mode is kCCDB. Please inject the valid CCDB object via setVertexMode";
}
}
auto sampledvertex = mMeanVertex->sample();
if (PrimaryGeneratorParam::Instance().verbose) {
LOG(info) << "Sampled interacting vertex " << sampledvertex;
}
SetBeam(sampledvertex.X(), sampledvertex.Y(), 0., 0.);
SetTarget(sampledvertex.Z(), 0.);
}
/*****************************************************************/
Bool_t PrimaryGenerator::embedInto(TString fname)
{
/** embed into **/
/** check if a file is already open **/
if (mEmbedFile && mEmbedFile->IsOpen()) {
LOG(error) << "Another embedding file is currently open";
return kFALSE;
}
/** open file **/
mEmbedFile = TFile::Open(fname);
if (!mEmbedFile || !mEmbedFile->IsOpen()) {
LOG(error) << "Cannot open file for embedding: " << fname;
return kFALSE;
}
/** get tree **/
mEmbedTree = (TTree*)mEmbedFile->Get("o2sim");
if (!mEmbedTree) {
LOG(error) << R"(Cannot find "o2sim" tree for embedding in )" << fname;
return kFALSE;
}
/** get entries **/
mEmbedEntries = mEmbedTree->GetEntries();
if (mEmbedEntries <= 0) {
LOG(error) << "Invalid number of entries found in tree for embedding: " << mEmbedEntries;
return kFALSE;
}
/** connect MC event header **/
mEmbedEvent = new MCEventHeader;
mEmbedTree->SetBranchAddress("MCEventHeader.", &mEmbedEvent);
/** success **/
return kTRUE;
}
/*****************************************************************/
void PrimaryGenerator::setGeneratorInformation()
{
auto o2event = dynamic_cast<MCEventHeader*>(fEvent);
if (o2event) {
o2event->putInfo<int>(o2::mcgenid::GeneratorProperty::GENERATORID, mGeneratorId);
o2event->putInfo<std::string>(o2::mcgenid::GeneratorProperty::GENERATORDESCRIPTION, mGeneratorDescription);
}
}
/*****************************************************************/
/*****************************************************************/
} /* namespace eventgen */
} /* namespace o2 */
ClassImp(o2::eventgen::PrimaryGenerator);