Project
Loading...
Searching...
No Matches
PrimaryGenerator.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
21#include "DetectorsBase/Stack.h"
22#include <fairlogger/Logger.h>
23
24#include "FairGenericStack.h"
25#include "TFile.h"
26#include "TTree.h"
27
28#include "TDatabasePDG.h"
29#include "TVirtualMC.h"
30
32
33namespace o2
34{
35namespace eventgen
36{
37
38/*****************************************************************/
39
41{
43 LOG(info) << "Destructing PrimaryGenerator";
44 if (mEmbedFile && mEmbedFile->IsOpen()) {
45 mEmbedFile->Close();
46 delete mEmbedFile;
47 }
48 if (mEmbedEvent) {
49 delete mEmbedEvent;
50 }
51}
52
53/*****************************************************************/
54
56{
59 LOG(info) << "Initialising primary generator";
60
61 // set generator ID and description
65
67 if (mEmbedTree) {
68 LOG(info) << "Embedding into: " << mEmbedFile->GetName()
69 << " (" << mEmbedEntries << " events)";
70 return FairPrimaryGenerator::Init();
71 }
72
74 return FairPrimaryGenerator::Init();
75}
76
77/*****************************************************************/
78
80{
84 if (!mEmbedTree) {
85 fixInteractionVertex(); // <-- always fixes vertex outside of FairROOT
86 auto ret = FairPrimaryGenerator::GenerateEvent(pStack);
87 if (ret) {
88 setGeneratorInformation();
89 }
90 return ret;
91 }
92
96 mEmbedTree->GetEntry(mEmbedIndex);
98
100 auto genList = GetListOfGenerators();
101 for (int igen = 0; igen < genList->GetEntries(); ++igen) {
102 auto o2gen = dynamic_cast<Generator*>(genList->At(igen));
103 if (o2gen) {
105 }
106 }
107
109 if (!FairPrimaryGenerator::GenerateEvent(pStack)) {
110 return kFALSE;
111 }
112
114 auto o2event = dynamic_cast<MCEventHeader*>(fEvent);
115 if (o2event) {
116 o2event->setEmbeddingFileName(mEmbedFile->GetName());
117 o2event->setEmbeddingEventIndex(mEmbedIndex);
118 }
119 setGeneratorInformation();
120
122 mEmbedIndex++;
124
126 return kTRUE;
127}
128
129/*****************************************************************/
130
131void PrimaryGenerator::AddTrack(Int_t pdgid, Double_t px, Double_t py, Double_t pz,
132 Double_t vx, Double_t vy, Double_t vz,
133 Int_t mother1, Int_t mother2,
134 Int_t daughter1, Int_t daughter2,
135 Bool_t wanttracking,
136 Double_t e, Double_t tof,
137 Double_t weight, TMCProcess proc, Int_t generatorStatus)
138{
141 // check the status encoding
142 if (!mcgenstatus::isEncoded(generatorStatus) && proc == TMCProcess::kPPrimary) {
143 LOG(fatal) << "Generator status " << generatorStatus << " of particle is not encoded properly.";
144 }
145
147 vx += fVertex.X();
148 vy += fVertex.Y();
149 vz += fVertex.Z();
150
152 auto particlePDG = TDatabasePDG::Instance()->GetParticle(pdgid);
153 if (wanttracking && !particlePDG) {
154 LOG(warn) << "Particle to be tracked is not defined in PDG: pdg = " << pdgid << " (disabling tracking)";
155 wanttracking = false;
156 }
157
159 Int_t doTracking = 0; // Go to tracking
160 if (fdoTracking && wanttracking) {
161 doTracking = 1;
162 }
163 Int_t dummyparent = -1; // Primary particle (now the value is -1 by default)
164 Double_t polx = 0.; // Polarisation
165 Double_t poly = 0.;
166 Double_t polz = 0.;
167 Int_t ntr = 0; // Track number; to be filled by the stack
168 Int_t status = generatorStatus; // Generation status
169
170 // correct for tracks which are in list before generator is called
171 if (mother1 != -1) {
172 mother1 += fMCIndexOffset;
173 }
174 if (mother2 != -1) {
175 mother2 += fMCIndexOffset;
176 }
177 if (daughter1 != -1) {
178 daughter1 += fMCIndexOffset;
179 }
180 if (daughter2 != -1) {
181 daughter2 += fMCIndexOffset;
182 }
183
192 if (abs(pdgid) == 311 && doTracking) {
193 LOG(warn) << "K0/antiK0 requested for tracking: converting into K0s/K0L";
194 pdgid = gRandom->Uniform() < 0.5 ? 310 : 130;
195 }
196
198 if (e < 0) {
199 double mass = particlePDG ? particlePDG->Mass() : 0.;
200 e = std::sqrt(mass * mass + px * px + py * py + pz * pz);
201 }
202
204 auto stack = dynamic_cast<o2::data::Stack*>(fStack);
205 if (!stack) {
206 LOG(fatal) << "Stack must be an o2::data:Stack";
207 return; // must be the o2 stack
208 }
209 stack->PushTrack(doTracking, mother1, pdgid, px, py, pz,
210 e, vx, vy, vz, tof, polx, poly, polz, TMCProcess::kPPrimary, ntr,
211 weight, status, mother2, daughter1, daughter2, proc);
212
213 fNTracks++;
214}
215
216/*****************************************************************/
217
218void PrimaryGenerator::AddTrack(Int_t pdgid, Double_t px, Double_t py,
219 Double_t pz, Double_t vx, Double_t vy,
220 Double_t vz, Int_t parent, Bool_t wanttracking,
221 Double_t e, Double_t tof, Double_t weight, TMCProcess proc)
222{
223 // Do this to encode status code correctly. In FairRoot's PrimaryGenerator, this is simply one number that is set to 0.
224 // So we basically do the same.
225 // Assuming that this is treated as the HepMC status code down the line (as it used to be).
226 AddTrack(pdgid, px, py, pz, vx, vy, vz, parent, -1, -1, -1, wanttracking,
227 e, tof, weight, proc, mcgenstatus::MCGenStatusEncoding(0, 0).fullEncoding);
228}
229
230/*****************************************************************/
231
233{
234 if (!mApplyVertex) {
235 return;
236 }
239 Double_t xyz[3] = {event->GetX(), event->GetY(), event->GetZ()};
240 SetBeam(xyz[0], xyz[1], 0., 0.);
241 SetTarget(xyz[2], 0.);
242 SmearVertexXY(false);
243 SmearVertexZ(false);
244 SmearGausVertexXY(false);
245 SmearGausVertexZ(false);
246}
247
248/*****************************************************************/
250{
251 if (!mApplyVertex) {
252 return;
253 }
257 mHaveExternalVertex = true;
258}
259
260/*****************************************************************/
261
263{
266 if (!v) {
267 LOG(fatal) << "A valid MeanVertexObject needs to be passed with option o2::conf::VertexMode::kCCDB";
268 }
269 mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(*v)));
270 LOG(info) << "The mean vertex is set to :";
271 mMeanVertex->print();
272 }
274 setApplyVertex(false);
275 LOG(info) << "Disabling vertexing";
276 mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(0, 0, 0, 0, 0, 0, 0, 0)));
277 LOG(info) << "The mean vertex is set to :";
278 mMeanVertex->print();
279 }
280}
281
282/*****************************************************************/
283
285{
286 if (!mApplyVertex) {
287 SetBeam(0., 0., 0., 0.);
288 SetTarget(0., 0.);
289 return;
290 }
291
292 // if someone gave vertex from outside; we will take it
294 SetBeam(mExternalVertexX, mExternalVertexY, 0., 0.);
295 SetTarget(mExternalVertexZ, 0.);
296 mHaveExternalVertex = false; // the vertex is now consumed
297 return;
298 }
299
300 // sampling a vertex and fixing for next event; no smearing will be done
301 // inside FairPrimaryGenerator;
302 SmearVertexXY(false);
303 SmearVertexZ(false);
304 SmearGausVertexXY(false);
305 SmearGausVertexZ(false);
306
307 // we use the mMeanVertexObject if initialized (initialize first)
308 if (mMeanVertex.get() == nullptr) {
311 const auto& xyz = param.position;
312 const auto& sigma = param.width;
313 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)));
314 }
316 mMeanVertex = std::move(std::unique_ptr<o2::dataformats::MeanVertexObject>(new o2::dataformats::MeanVertexObject(0, 0, 0, 0, 0, 0, 0, 0)));
317 }
319 // fatal.. then the object should have been passed with setting
320 LOG(fatal) << "MeanVertexObject is null ... but mode is kCCDB. Please inject the valid CCDB object via setVertexMode";
321 }
322 }
323 auto sampledvertex = mMeanVertex->sample();
324
325 if (PrimaryGeneratorParam::Instance().verbose) {
326 LOG(info) << "Sampled interacting vertex " << sampledvertex;
327 }
328 SetBeam(sampledvertex.X(), sampledvertex.Y(), 0., 0.);
329 SetTarget(sampledvertex.Z(), 0.);
330}
331
332/*****************************************************************/
333
334Bool_t PrimaryGenerator::embedInto(TString fname)
335{
339 if (mEmbedFile && mEmbedFile->IsOpen()) {
340 LOG(error) << "Another embedding file is currently open";
341 return kFALSE;
342 }
343
345 mEmbedFile = TFile::Open(fname);
346 if (!mEmbedFile || !mEmbedFile->IsOpen()) {
347 LOG(error) << "Cannot open file for embedding: " << fname;
348 return kFALSE;
349 }
350
352 mEmbedTree = (TTree*)mEmbedFile->Get("o2sim");
353 if (!mEmbedTree) {
354 LOG(error) << R"(Cannot find "o2sim" tree for embedding in )" << fname;
355 return kFALSE;
356 }
357
359 mEmbedEntries = mEmbedTree->GetEntries();
360 if (mEmbedEntries <= 0) {
361 LOG(error) << "Invalid number of entries found in tree for embedding: " << mEmbedEntries;
362 return kFALSE;
363 }
364
367 mEmbedTree->SetBranchAddress("MCEventHeader.", &mEmbedEvent);
368
370 return kTRUE;
371}
372
373/*****************************************************************/
374
375void PrimaryGenerator::setGeneratorInformation()
376{
377 auto o2event = dynamic_cast<MCEventHeader*>(fEvent);
378 if (o2event) {
379 o2event->putInfo<int>(o2::mcgenid::GeneratorProperty::GENERATORID, mGeneratorId);
380 o2event->putInfo<std::string>(o2::mcgenid::GeneratorProperty::GENERATORDESCRIPTION, mGeneratorDescription);
381 }
382}
383
384/*****************************************************************/
385/*****************************************************************/
386
387} /* namespace eventgen */
388} /* namespace o2 */
389
Definition of the Stack class.
ClassImp(o2::eventgen::PrimaryGenerator)
uint32_t stack
Definition RawData.h:1
void PushTrack(Int_t toBeDone, Int_t parentID, Int_t pdgCode, Double_t px, Double_t py, Double_t pz, Double_t e, Double_t vx, Double_t vy, Double_t vz, Double_t time, Double_t polx, Double_t poly, Double_t polz, TMCProcess proc, Int_t &ntr, Double_t weight, Int_t is) override
Definition Stack.cxx:176
void setEmbeddingFileName(std::string const &value)
void putInfo(std::string const &key, T const &value)
virtual void notifyEmbedding(const o2::dataformats::MCEventHeader *eventHeader)
Definition Generator.h:99
void setGeneratorDescription(std::string const &desc)
Bool_t GenerateEvent(FairGenericStack *pStack) override
void setVertexMode(o2::conf::VertexMode const &mode, o2::dataformats::MeanVertexObject const *obj=nullptr)
void setInteractionVertex(const o2::dataformats::MCEventHeader *event)
void 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=-1, Int_t mother2=-1, Int_t daughter1=-1, Int_t daughter2=-1, Bool_t wanttracking=true, Double_t e=-9e9, Double_t tof=0., Double_t weight=0., TMCProcess proc=kPPrimary, Int_t generatorStatus=0)
void setExternalVertexForNextEvent(double x, double y, double z)
o2::dataformats::MCEventHeader * mEmbedEvent
std::unique_ptr< o2::dataformats::MeanVertexObject > mMeanVertex
static constexpr Property GENERATORID
static constexpr Property GENERATORDESCRIPTION
struct _cl_event * event
Definition glcorearb.h:2982
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum mode
Definition glcorearb.h:266
const GLdouble * v
Definition glcorearb.h:832
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLint y
Definition glcorearb.h:270
GLenum const GLfloat * params
Definition glcorearb.h:272
GLenum GLfloat param
Definition glcorearb.h:271
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
bool isEncoded(MCGenStatusEncoding statusCode)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"