Project
Loading...
Searching...
No Matches
GeneratorService.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
14#include "SimConfig/SimConfig.h"
17
18using namespace o2::eventgen;
19
20void GeneratorService::initService(std::string const& genName,
21 std::string const& triggerName,
22 VertexOption const& vxtOption)
23{
24 auto localSimConfig = o2::conf::SimConfig::make();
25 localSimConfig.getConfigData().mGenerator = genName;
26 localSimConfig.getConfigData().mTrigger = triggerName;
27 localSimConfig.getConfigData().mNEvents = o2::eventgen::Generator::getTotalNEvents();
28
30
31 // initialize vertexing based on type
32 if (dynamic_cast<MeanVertexObjectOption const*>(&vxtOption) != nullptr) {
33 auto ccdboption = dynamic_cast<MeanVertexObjectOption const*>(&vxtOption);
34 LOG(info) << "Init prim gen with MeanVertexObject";
35 // assign the mean vertex object
36 if (!ccdboption->meanVertexObject) {
37 LOG(fatal) << "No mean vertex object found - Cannot initialize Event generator";
38 }
39 mPrimGen.setVertexMode(o2::conf::VertexMode::kCCDB, ccdboption->meanVertexObject);
40 } else if (dynamic_cast<NoVertexOption const*>(&vxtOption) != nullptr) {
42 } else if (dynamic_cast<DiamondParamVertexOption const*>(&vxtOption) != nullptr) {
44 } else {
45 LOG(error) << "Unknown VertexOption passed to Generator initialization";
46 }
47
48 mStack.setExternalMode(true);
49 mPrimGen.Init();
50}
51
53{
54 mPrimGen.SetEvent(&header);
55 mStack.Reset();
56 mPrimGen.GenerateEvent(&mStack); // this is the usual FairROOT interface going via stack
57
58 tracks.reserve(mStack.getPrimaries().size());
59 for (auto& tparticle : mStack.getPrimaries()) {
60 tracks.emplace_back(tparticle);
61 }
62}
63
65{
66 std::vector<o2::MCTrack> tracks;
68 generateEvent_MCTracks(tracks, header);
69 return std::pair<std::vector<MCTrack>, o2::dataformats::MCEventHeader>(tracks, header);
70}
71
73{
74 mPrimGen.SetEvent(&header);
75 mStack.Reset();
76 mPrimGen.GenerateEvent(&mStack); // this is the usual FairROOT interface going via stack
77
78 tracks.clear();
79 tracks = mStack.getPrimaries();
80}
static SimConfig make()
Definition SimConfig.h:118
void setExternalMode(bool m)
Definition Stack.h:205
void Reset() override
Resets arrays and stack and deletes particles and tracks.
Definition Stack.cxx:588
const std::vector< TParticle > & getPrimaries() const
Definition Stack.h:187
std::pair< std::vector< MCTrack >, o2::dataformats::MCEventHeader > generateEvent()
void generateEvent_MCTracks(std::vector< MCTrack > &tracks, o2::dataformats::MCEventHeader &header)
void generateEvent_TParticles(std::vector< TParticle > &tparts, o2::dataformats::MCEventHeader &header)
void initService(std::string const &generatorName, std::string const &triggerName, VertexOption const &vtxOption)
static unsigned int getTotalNEvents()
Definition Generator.h:93
Bool_t GenerateEvent(FairGenericStack *pStack) override
void setVertexMode(o2::conf::VertexMode const &mode, o2::dataformats::MeanVertexObject const *obj=nullptr)
static void setPrimaryGenerator(o2::conf::SimConfig const &, FairPrimaryGenerator *)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"