Project
Loading...
Searching...
No Matches
GeneratorFactory.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
16#include "FairGenerator.h"
17#include "FairBoxGenerator.h"
18#include <fairlogger/Logger.h>
19#include <SimConfig/SimConfig.h>
23#ifdef GENERATORS_WITH_PYTHIA8
26#endif
30#ifdef GENERATORS_WITH_HEPMC3
33#endif
34#if defined(GENERATORS_WITH_PYTHIA8) && defined(GENERATORS_WITH_HEPMC3)
37#endif
44
45#include "TRandom.h"
46
47namespace o2
48{
49namespace eventgen
50{
51
52// reusable helper class
53// main purpose is to init a FairPrimGen given some (Sim)Config
55{
56 if (!primGen) {
57 LOG(warning) << "No primary generator instance; Cannot setup";
58 return;
59 }
60
61 auto primGenO2 = dynamic_cast<PrimaryGenerator*>(primGen);
62
63 auto makeBoxGen = [](int pdgid, int mult, double etamin, double etamax, double pmin, double pmax, double phimin, double phimax, bool debug = false) {
64 auto gen = new FairBoxGenerator(pdgid, mult);
65 gen->SetEtaRange(etamin, etamax);
66 gen->SetPRange(pmin, pmax);
67 gen->SetPhiRange(phimin, phimax);
68 gen->SetDebug(debug);
69 return gen;
70 };
71
72#ifdef GENERATORS_WITH_PYTHIA8
73 auto makePythia8Gen = [](std::string& config) {
74 auto& singleton = GeneratorPythia8Param::Instance();
76 .config = config.size() > 0 ? config : singleton.config,
77 .hooksFileName = singleton.hooksFileName,
78 .hooksFuncName = singleton.hooksFuncName,
79 .includePartonEvent = singleton.includePartonEvent,
80 .particleFilter = singleton.particleFilter,
81 .verbose = singleton.verbose,
82 };
83 auto gen = new o2::eventgen::GeneratorPythia8(pars);
84 if (!config.empty()) {
85 LOG(info) << "Setting \'Pythia8\' base configuration: " << config << std::endl;
86 gen->setConfig(config); // assign config; will be executed in Init function
87 }
88 return gen;
89 };
90#endif
91
94 o2::O2DatabasePDG::addALICEParticles(TDatabasePDG::Instance());
95 auto genconfig = conf.getGenerator();
96#if defined(GENERATORS_WITH_PYTHIA8) && defined(GENERATORS_WITH_HEPMC3)
97 if (GeneratorHybridParam::Instance().switchExtToHybrid && (genconfig.compare("external") == 0 || genconfig.compare("extgen") == 0)) {
98 LOG(info) << "Switching external generator to hybrid mode";
99 genconfig = "hybrid";
100 }
101#endif
102 LOG(info) << "** Generator to use: '" << genconfig << "'";
103 if (genconfig.compare("boxgen") == 0) {
104 // a simple "box" generator configurable via BoxGunparam
105 auto& boxparam = BoxGunParam::Instance();
106 LOG(info) << "Init generic box generator with following parameters";
107 LOG(info) << boxparam;
108 auto boxGen = makeBoxGen(boxparam.pdg, boxparam.number, boxparam.eta[0], boxparam.eta[1], boxparam.prange[0], boxparam.prange[1], boxparam.phirange[0], boxparam.phirange[1], boxparam.debug);
109 primGen->AddGenerator(boxGen);
110 } else if (genconfig.compare("fwmugen") == 0) {
111 // a simple "box" generator for forward muons
112 LOG(info) << "Init box forward muons generator";
113 auto boxGen = makeBoxGen(13, 1, -4, -2.5, 50., 50., 0., 360);
114 primGen->AddGenerator(boxGen);
115 } else if (genconfig.compare("hmpidgun") == 0) {
116 // a simple "box" generator for forward muons
117 LOG(info) << "Init hmpid gun generator";
118 auto boxGen = makeBoxGen(-211, 100, -0.5, -0.5, 2, 5, -5, 60);
119 primGen->AddGenerator(boxGen);
120 } else if (genconfig.compare("fwpigen") == 0) {
121 // a simple "box" generator for forward pions
122 LOG(info) << "Init box forward pions generator";
123 auto boxGen = makeBoxGen(-211, 10, -4, -2.5, 7, 7, 0, 360);
124 primGen->AddGenerator(boxGen);
125 } else if (genconfig.compare("fwrootino") == 0) {
126 // a simple "box" generator for forward rootinos
127 LOG(info) << "Init box forward rootinos generator";
128 auto boxGen = makeBoxGen(0, 1, -4, -2.5, 1, 5, 0, 360);
129 primGen->AddGenerator(boxGen);
130 } else if (genconfig.compare("zdcgen") == 0) {
131 // a simple "box" generator for forward neutrons
132 LOG(info) << "Init box forward/backward zdc generator";
133 auto boxGenC = makeBoxGen(2112 /*neutrons*/, 1, -8, -9999, 500, 1000, 0., 360.);
134 auto boxGenA = makeBoxGen(2112 /*neutrons*/, 1, 8, 9999, 500, 1000, 0., 360.);
135 primGen->AddGenerator(boxGenC);
136 primGen->AddGenerator(boxGenA);
137 } else if (genconfig.compare("emcgenele") == 0) {
138 // box generator with one electron per event
139 LOG(info) << "Init box generator for electrons in EMCAL";
140 // using phi range of emcal
141 auto elecgen = makeBoxGen(11, 1, -0.67, 0.67, 15, 15, 80, 187);
142 primGen->AddGenerator(elecgen);
143 } else if (genconfig.compare("emcgenphoton") == 0) {
144 LOG(info) << "Init box generator for photons in EMCAL";
145 auto photongen = makeBoxGen(22, 1, -0.67, 0.67, 15, 15, 80, 187);
146 primGen->AddGenerator(photongen);
147 } else if (genconfig.compare("fddgen") == 0) {
148 LOG(info) << "Init box FDD generator";
149 auto boxGenFDC = makeBoxGen(13, 1000, -7, -4.8, 10, 500, 0, 360.);
150 auto boxGenFDA = makeBoxGen(13, 1000, 4.9, 6.3, 10, 500, 0., 360);
151 primGen->AddGenerator(boxGenFDA);
152 primGen->AddGenerator(boxGenFDC);
153 } else if (genconfig.compare("extkin") == 0) {
154 // external kinematics
155 // needs precense of a kinematics file "Kinematics.root"
156 // TODO: make this configurable and check for presence
157 auto extGen = new o2::eventgen::GeneratorFromFile(conf.getExtKinematicsFileName().c_str());
158 extGen->SetStartEvent(conf.getStartEvent());
159 primGen->AddGenerator(extGen);
160 LOG(info) << "using external kinematics";
161 } else if (genconfig.compare("extkinO2") == 0) {
162 // external kinematics from previous O2 output
163 auto& singleton = GeneratorFromO2KineParam::Instance();
164 auto name1 = singleton.fileName;
165 auto name2 = conf.getExtKinematicsFileName();
166 auto pars = O2KineGenConfig{
167 .skipNonTrackable = singleton.skipNonTrackable,
168 .continueMode = singleton.continueMode,
169 .roundRobin = singleton.roundRobin,
170 .randomize = singleton.randomize,
171 .rngseed = singleton.rngseed,
172 .randomphi = singleton.randomphi,
173 .fileName = name1.size() > 0 ? name1.c_str() : name2.c_str()};
174 auto extGen = new o2::eventgen::GeneratorFromO2Kine(pars);
175 extGen->SetStartEvent(conf.getStartEvent());
176 primGen->AddGenerator(extGen);
177 if (pars.continueMode) {
178 auto o2PrimGen = dynamic_cast<o2::eventgen::PrimaryGenerator*>(primGen);
179 if (o2PrimGen) {
180 o2PrimGen->setApplyVertex(false);
181 }
182 }
183 LOG(info) << "using external O2 kinematics";
184 } else if (genconfig.compare("evtpool") == 0) {
185 // case of an "event-pool" which is a specialization of extkinO2
186 // with some additional logic in file management and less configurability
187 // and not features such as "continue transport"
189 primGen->AddGenerator(extGen);
190 LOG(info) << "using the eventpool generator";
191 } else if (genconfig.compare("tparticle") == 0) {
192 // External ROOT file(s) with tree of TParticle in clones array,
193 // or external program generating such a file
194 auto& param0 = GeneratorFileOrCmdParam::Instance();
196 LOG(info) << "Init 'GeneratorTParticle' with the following parameters";
197 LOG(info) << param0;
198 LOG(info) << param;
199 auto tgen = new o2::eventgen::GeneratorTParticle();
200 tgen->setup(param0, param, conf);
201 primGen->AddGenerator(tgen);
202#ifdef GENERATORS_WITH_HEPMC3
203 } else if (genconfig.compare("hepmc") == 0) {
204 // external HepMC file, or external program writing HepMC event
205 // records to standard output.
206 auto& param0 = GeneratorFileOrCmdParam::Instance();
208 LOG(info) << "Init \'GeneratorHepMC\' with following parameters";
209 LOG(info) << param0;
210 LOG(info) << param;
211 auto hepmcGen = new o2::eventgen::GeneratorHepMC();
212 hepmcGen->setup(param0, param, conf);
213 primGen->AddGenerator(hepmcGen);
214#endif
215#ifdef GENERATORS_WITH_PYTHIA8
216 } else if (genconfig.compare("alldets") == 0) {
217 // a simple generator for test purposes - making sure to generate hits
218 // in all detectors
219 // I compose it of:
220 // 1) pythia8
221 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_inel.cfg";
222 auto py8 = makePythia8Gen(py8config);
223 primGen->AddGenerator(py8);
224 // 2) forward muons
225 auto muon = makeBoxGen(13, 100, -2.5, -4.0, 100, 100, 0., 360);
226 primGen->AddGenerator(muon);
227 } else if (genconfig.compare("pythia8") == 0) {
228 auto py8config = std::string();
229 auto py8 = makePythia8Gen(py8config);
230 primGen->AddGenerator(py8);
231 } else if (genconfig.compare("pythia8pp") == 0) {
232 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_inel.cfg";
233 auto py8 = makePythia8Gen(py8config);
234 primGen->AddGenerator(py8);
235 } else if (genconfig.compare("pythia8hf") == 0) {
236 // pythia8 pp (HF production)
237 // configures pythia for HF production in pp collisions at 14 TeV
238 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_hf.cfg";
239 auto py8 = makePythia8Gen(py8config);
240 primGen->AddGenerator(py8);
241 } else if (genconfig.compare("pythia8hi") == 0) {
242 // pythia8 heavy-ion
243 // exploits pythia8 heavy-ion machinery (available from v8.230)
244 // configures pythia for min.bias Pb-Pb collisions at 5.52 TeV
245 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_hi.cfg";
246 auto py8 = makePythia8Gen(py8config);
247 primGen->AddGenerator(py8);
248 } else if (genconfig.compare("pythia8powheg") == 0) {
249 // pythia8 with powheg
250 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_powheg.cfg";
251 auto py8 = makePythia8Gen(py8config);
252 primGen->AddGenerator(py8);
253#endif
254 } else if (genconfig.compare("external") == 0 || genconfig.compare("extgen") == 0) {
255 // external generator via configuration macro
257 LOG(info) << "Setting up external generator with following parameters";
258 LOG(info) << params;
259 auto extgen_filename = params.fileName;
260 auto extgen_func = params.funcName;
261 auto extgen = o2::conf::GetFromMacro<FairGenerator*>(extgen_filename, extgen_func, "FairGenerator*", "extgen");
262 if (!extgen) {
263 LOG(fatal) << "Failed to retrieve \'extgen\': problem with configuration ";
264 }
265 primGen->AddGenerator(extgen);
266 } else if (genconfig.compare("toftest") == 0) { // 1 muon per sector and per module
267 LOG(info) << "Init tof test generator -> 1 muon per sector and per module";
268 for (int i = 0; i < 18; i++) {
269 for (int j = 0; j < 5; j++) {
270 auto boxGen = new FairBoxGenerator(13, 1); /*protons*/
271 boxGen->SetEtaRange(-0.8 + 0.32 * j + 0.15, -0.8 + 0.32 * j + 0.17);
272 boxGen->SetPRange(9, 10);
273 boxGen->SetPhiRange(10 + 20. * i - 1, 10 + 20. * i + 1);
274 boxGen->SetDebug(kTRUE);
275 primGen->AddGenerator(boxGen);
276 }
277 }
278#if defined(GENERATORS_WITH_PYTHIA8) && defined(GENERATORS_WITH_HEPMC3)
279 } else if (genconfig.compare("hybrid") == 0) { // hybrid using multiple generators
280 LOG(info) << "Init hybrid generator";
281 auto& hybridparam = GeneratorHybridParam::Instance();
282 std::string config = hybridparam.configFile;
283 // check if config string points to an existing and not empty file
284 if (config.empty()) {
285 LOG(fatal) << "No configuration file provided for hybrid generator";
286 return;
287 }
288 auto& hybrid = o2::eventgen::GeneratorHybrid::Instance(config);
289 primGen->AddGenerator(&hybrid);
290#endif
291 } else {
292 LOG(fatal) << "Invalid generator";
293 }
294
296 // to be set via GeneratorFactory only if generator is not hybrid
297 // external settings via JSON are supported in the latter
298
299 Trigger trigger = nullptr;
300 DeepTrigger deeptrigger = nullptr;
301 if (!(genconfig.compare("hybrid") == 0)) {
302 auto trgconfig = conf.getTrigger();
303 if (trgconfig.empty()) {
304 return;
305 } else if (trgconfig.compare("particle") == 0) {
307 } else if (trgconfig.compare("external") == 0) {
308 // external trigger via configuration macro
310 LOG(info) << "Setting up external trigger with following parameters";
311 LOG(info) << params;
312 auto external_trigger_filename = params.fileName;
313 auto external_trigger_func = params.funcName;
314 trigger = o2::conf::GetFromMacro<o2::eventgen::Trigger>(external_trigger_filename, external_trigger_func, "o2::eventgen::Trigger", "trigger");
315 if (!trigger) {
316 LOG(info) << "Trying to retrieve a \'o2::eventgen::DeepTrigger\' type" << std::endl;
317 deeptrigger = o2::conf::GetFromMacro<o2::eventgen::DeepTrigger>(external_trigger_filename, external_trigger_func, "o2::eventgen::DeepTrigger", "deeptrigger");
318 }
319 if (!trigger && !deeptrigger) {
320 LOG(fatal) << "Failed to retrieve \'external trigger\': problem with configuration ";
321 }
322 } else {
323 LOG(fatal) << "Invalid trigger";
324 }
325
327 auto generators = primGen->GetListOfGenerators();
328 for (int igen = 0; igen < generators->GetEntries(); ++igen) {
329 auto generator = dynamic_cast<o2::eventgen::Generator*>(generators->At(igen));
330 if (!generator) {
331 LOG(fatal) << "request to add a trigger to an unsupported generator";
332 return;
333 }
335 if (trigger) {
336 generator->addTrigger(trigger);
337 }
338 if (deeptrigger) {
339 generator->addDeepTrigger(deeptrigger);
340 }
341 }
342 }
343}
344
345} // end namespace eventgen
346} // end namespace o2
default_random_engine gen(dev())
std::ostringstream debug
int32_t i
uint32_t j
Definition RawData.h:0
static void addALICEParticles(TDatabasePDG *db=TDatabasePDG::Instance())
std::string getExtKinematicsFileName() const
Definition SimConfig.h:158
unsigned int getStartEvent() const
Definition SimConfig.h:160
std::string getTrigger() const
Definition SimConfig.h:155
std::string getGenerator() const
Definition SimConfig.h:154
static GeneratorHybrid & Instance(const std::string &inputgens="")
void setTriggerMode(ETriggerMode_t val)
Definition Generator.h:94
GLenum const GLfloat * params
Definition glcorearb.h:272
GLenum GLfloat param
Definition glcorearb.h:271
std::function< bool(void *, std::string)> DeepTrigger
Definition Trigger.h:26
std::function< bool(const std::vector< TParticle > &)> Trigger
Definition Trigger.h:25
Trigger TriggerParticle(const TriggerParticleParam &param)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
static void setPrimaryGenerator(o2::conf::SimConfig const &, FairPrimaryGenerator *)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"