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 LOG(info) << "** Generator to use: '" << genconfig << "'";
97 if (genconfig.compare("boxgen") == 0) {
98 // a simple "box" generator configurable via BoxGunparam
99 auto& boxparam = BoxGunParam::Instance();
100 LOG(info) << "Init generic box generator with following parameters";
101 LOG(info) << boxparam;
102 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);
103 primGen->AddGenerator(boxGen);
104 } else if (genconfig.compare("fwmugen") == 0) {
105 // a simple "box" generator for forward muons
106 LOG(info) << "Init box forward muons generator";
107 auto boxGen = makeBoxGen(13, 1, -4, -2.5, 50., 50., 0., 360);
108 primGen->AddGenerator(boxGen);
109 } else if (genconfig.compare("hmpidgun") == 0) {
110 // a simple "box" generator for forward muons
111 LOG(info) << "Init hmpid gun generator";
112 auto boxGen = makeBoxGen(-211, 100, -0.5, -0.5, 2, 5, -5, 60);
113 primGen->AddGenerator(boxGen);
114 } else if (genconfig.compare("fwpigen") == 0) {
115 // a simple "box" generator for forward pions
116 LOG(info) << "Init box forward pions generator";
117 auto boxGen = makeBoxGen(-211, 10, -4, -2.5, 7, 7, 0, 360);
118 primGen->AddGenerator(boxGen);
119 } else if (genconfig.compare("fwrootino") == 0) {
120 // a simple "box" generator for forward rootinos
121 LOG(info) << "Init box forward rootinos generator";
122 auto boxGen = makeBoxGen(0, 1, -4, -2.5, 1, 5, 0, 360);
123 primGen->AddGenerator(boxGen);
124 } else if (genconfig.compare("zdcgen") == 0) {
125 // a simple "box" generator for forward neutrons
126 LOG(info) << "Init box forward/backward zdc generator";
127 auto boxGenC = makeBoxGen(2112 /*neutrons*/, 1, -8, -9999, 500, 1000, 0., 360.);
128 auto boxGenA = makeBoxGen(2112 /*neutrons*/, 1, 8, 9999, 500, 1000, 0., 360.);
129 primGen->AddGenerator(boxGenC);
130 primGen->AddGenerator(boxGenA);
131 } else if (genconfig.compare("emcgenele") == 0) {
132 // box generator with one electron per event
133 LOG(info) << "Init box generator for electrons in EMCAL";
134 // using phi range of emcal
135 auto elecgen = makeBoxGen(11, 1, -0.67, 0.67, 15, 15, 80, 187);
136 primGen->AddGenerator(elecgen);
137 } else if (genconfig.compare("emcgenphoton") == 0) {
138 LOG(info) << "Init box generator for photons in EMCAL";
139 auto photongen = makeBoxGen(22, 1, -0.67, 0.67, 15, 15, 80, 187);
140 primGen->AddGenerator(photongen);
141 } else if (genconfig.compare("fddgen") == 0) {
142 LOG(info) << "Init box FDD generator";
143 auto boxGenFDC = makeBoxGen(13, 1000, -7, -4.8, 10, 500, 0, 360.);
144 auto boxGenFDA = makeBoxGen(13, 1000, 4.9, 6.3, 10, 500, 0., 360);
145 primGen->AddGenerator(boxGenFDA);
146 primGen->AddGenerator(boxGenFDC);
147 } else if (genconfig.compare("extkin") == 0) {
148 // external kinematics
149 // needs precense of a kinematics file "Kinematics.root"
150 // TODO: make this configurable and check for presence
151 auto extGen = new o2::eventgen::GeneratorFromFile(conf.getExtKinematicsFileName().c_str());
152 extGen->SetStartEvent(conf.getStartEvent());
153 primGen->AddGenerator(extGen);
154 LOG(info) << "using external kinematics";
155 } else if (genconfig.compare("extkinO2") == 0) {
156 // external kinematics from previous O2 output
157 auto& singleton = GeneratorFromO2KineParam::Instance();
158 auto name1 = singleton.fileName;
159 auto name2 = conf.getExtKinematicsFileName();
160 auto pars = O2KineGenConfig{
161 .skipNonTrackable = singleton.skipNonTrackable,
162 .continueMode = singleton.continueMode,
163 .roundRobin = singleton.roundRobin,
164 .randomize = singleton.randomize,
165 .rngseed = singleton.rngseed,
166 .randomphi = singleton.randomphi,
167 .fileName = name1.size() > 0 ? name1.c_str() : name2.c_str()};
168 auto extGen = new o2::eventgen::GeneratorFromO2Kine(pars);
169 extGen->SetStartEvent(conf.getStartEvent());
170 primGen->AddGenerator(extGen);
171 if (pars.continueMode) {
172 auto o2PrimGen = dynamic_cast<o2::eventgen::PrimaryGenerator*>(primGen);
173 if (o2PrimGen) {
174 o2PrimGen->setApplyVertex(false);
175 }
176 }
177 LOG(info) << "using external O2 kinematics";
178 } else if (genconfig.compare("evtpool") == 0) {
179 // case of an "event-pool" which is a specialization of extkinO2
180 // with some additional logic in file management and less configurability
181 // and not features such as "continue transport"
183 primGen->AddGenerator(extGen);
184 LOG(info) << "using the eventpool generator";
185 } else if (genconfig.compare("tparticle") == 0) {
186 // External ROOT file(s) with tree of TParticle in clones array,
187 // or external program generating such a file
188 auto& param0 = GeneratorFileOrCmdParam::Instance();
190 LOG(info) << "Init 'GeneratorTParticle' with the following parameters";
191 LOG(info) << param0;
192 LOG(info) << param;
193 auto tgen = new o2::eventgen::GeneratorTParticle();
194 tgen->setup(param0, param, conf);
195 primGen->AddGenerator(tgen);
196#ifdef GENERATORS_WITH_HEPMC3
197 } else if (genconfig.compare("hepmc") == 0) {
198 // external HepMC file, or external program writing HepMC event
199 // records to standard output.
200 auto& param0 = GeneratorFileOrCmdParam::Instance();
202 LOG(info) << "Init \'GeneratorHepMC\' with following parameters";
203 LOG(info) << param0;
204 LOG(info) << param;
205 auto hepmcGen = new o2::eventgen::GeneratorHepMC();
206 hepmcGen->setup(param0, param, conf);
207 primGen->AddGenerator(hepmcGen);
208#endif
209#ifdef GENERATORS_WITH_PYTHIA8
210 } else if (genconfig.compare("alldets") == 0) {
211 // a simple generator for test purposes - making sure to generate hits
212 // in all detectors
213 // I compose it of:
214 // 1) pythia8
215 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_inel.cfg";
216 auto py8 = makePythia8Gen(py8config);
217 primGen->AddGenerator(py8);
218 // 2) forward muons
219 auto muon = makeBoxGen(13, 100, -2.5, -4.0, 100, 100, 0., 360);
220 primGen->AddGenerator(muon);
221 } else if (genconfig.compare("pythia8") == 0) {
222 auto py8config = std::string();
223 auto py8 = makePythia8Gen(py8config);
224 primGen->AddGenerator(py8);
225 } else if (genconfig.compare("pythia8pp") == 0) {
226 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_inel.cfg";
227 auto py8 = makePythia8Gen(py8config);
228 primGen->AddGenerator(py8);
229 } else if (genconfig.compare("pythia8hf") == 0) {
230 // pythia8 pp (HF production)
231 // configures pythia for HF production in pp collisions at 14 TeV
232 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_hf.cfg";
233 auto py8 = makePythia8Gen(py8config);
234 primGen->AddGenerator(py8);
235 } else if (genconfig.compare("pythia8hi") == 0) {
236 // pythia8 heavy-ion
237 // exploits pythia8 heavy-ion machinery (available from v8.230)
238 // configures pythia for min.bias Pb-Pb collisions at 5.52 TeV
239 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_hi.cfg";
240 auto py8 = makePythia8Gen(py8config);
241 primGen->AddGenerator(py8);
242 } else if (genconfig.compare("pythia8powheg") == 0) {
243 // pythia8 with powheg
244 auto py8config = std::string(std::getenv("O2_ROOT")) + "/share/Generators/egconfig/pythia8_powheg.cfg";
245 auto py8 = makePythia8Gen(py8config);
246 primGen->AddGenerator(py8);
247#endif
248 } else if (genconfig.compare("external") == 0 || genconfig.compare("extgen") == 0) {
249 // external generator via configuration macro
251 LOG(info) << "Setting up external generator with following parameters";
252 LOG(info) << params;
253 auto extgen_filename = params.fileName;
254 auto extgen_func = params.funcName;
255 auto extgen = o2::conf::GetFromMacro<FairGenerator*>(extgen_filename, extgen_func, "FairGenerator*", "extgen");
256 if (!extgen) {
257 LOG(fatal) << "Failed to retrieve \'extgen\': problem with configuration ";
258 }
259 primGen->AddGenerator(extgen);
260 } else if (genconfig.compare("toftest") == 0) { // 1 muon per sector and per module
261 LOG(info) << "Init tof test generator -> 1 muon per sector and per module";
262 for (int i = 0; i < 18; i++) {
263 for (int j = 0; j < 5; j++) {
264 auto boxGen = new FairBoxGenerator(13, 1); /*protons*/
265 boxGen->SetEtaRange(-0.8 + 0.32 * j + 0.15, -0.8 + 0.32 * j + 0.17);
266 boxGen->SetPRange(9, 10);
267 boxGen->SetPhiRange(10 + 20. * i - 1, 10 + 20. * i + 1);
268 boxGen->SetDebug(kTRUE);
269 primGen->AddGenerator(boxGen);
270 }
271 }
272#if defined(GENERATORS_WITH_PYTHIA8) && defined(GENERATORS_WITH_HEPMC3)
273 } else if (genconfig.compare("hybrid") == 0) { // hybrid using multiple generators
274 LOG(info) << "Init hybrid generator";
275 auto& hybridparam = GeneratorHybridParam::Instance();
276 std::string config = hybridparam.configFile;
277 // check if config string points to an existing and not empty file
278 if (config.empty()) {
279 LOG(fatal) << "No configuration file provided for hybrid generator";
280 return;
281 }
282 // check if file named config exists and it's not empty
283 else if (gSystem->AccessPathName(config.c_str())) {
284 LOG(fatal) << "Configuration file for hybrid generator does not exist";
285 return;
286 }
287 auto hybrid = new o2::eventgen::GeneratorHybrid(config);
288 primGen->AddGenerator(hybrid);
289#endif
290 } else {
291 LOG(fatal) << "Invalid generator";
292 }
293
295 // to be set via GeneratorFactory only if generator is not hybrid
296 // external settings via JSON are supported in the latter
297
298 Trigger trigger = nullptr;
299 DeepTrigger deeptrigger = nullptr;
300 if (!(genconfig.compare("hybrid") == 0)) {
301 auto trgconfig = conf.getTrigger();
302 if (trgconfig.empty()) {
303 return;
304 } else if (trgconfig.compare("particle") == 0) {
306 } else if (trgconfig.compare("external") == 0) {
307 // external trigger via configuration macro
309 LOG(info) << "Setting up external trigger with following parameters";
310 LOG(info) << params;
311 auto external_trigger_filename = params.fileName;
312 auto external_trigger_func = params.funcName;
313 trigger = o2::conf::GetFromMacro<o2::eventgen::Trigger>(external_trigger_filename, external_trigger_func, "o2::eventgen::Trigger", "trigger");
314 if (!trigger) {
315 LOG(info) << "Trying to retrieve a \'o2::eventgen::DeepTrigger\' type" << std::endl;
316 deeptrigger = o2::conf::GetFromMacro<o2::eventgen::DeepTrigger>(external_trigger_filename, external_trigger_func, "o2::eventgen::DeepTrigger", "deeptrigger");
317 }
318 if (!trigger && !deeptrigger) {
319 LOG(fatal) << "Failed to retrieve \'external trigger\': problem with configuration ";
320 }
321 } else {
322 LOG(fatal) << "Invalid trigger";
323 }
324
326 auto generators = primGen->GetListOfGenerators();
327 for (int igen = 0; igen < generators->GetEntries(); ++igen) {
328 auto generator = dynamic_cast<o2::eventgen::Generator*>(generators->At(igen));
329 if (!generator) {
330 LOG(fatal) << "request to add a trigger to an unsupported generator";
331 return;
332 }
334 if (trigger) {
335 generator->addTrigger(trigger);
336 }
337 if (deeptrigger) {
338 generator->addDeepTrigger(deeptrigger);
339 }
340 }
341 }
342}
343
344} // end namespace eventgen
345} // end namespace o2
default_random_engine gen(dev())
int32_t i
uint32_t j
Definition RawData.h:0
std::ostringstream debug
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
void setTriggerMode(ETriggerMode_t val)
Definition Generator.h:85
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"