Project
Loading...
Searching...
No Matches
GeneratorPythia8.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
17#include <fairlogger/Logger.h>
18#include "TParticle.h"
19#include "TF1.h"
20#include "TRandom.h"
24#if PYTHIA_VERSION_INTEGER >= 8310
25#include "Pythia8/HIInfo.h"
26#else
27#include "Pythia8/HIUserHooks.h"
28#endif
29#include "Pythia8Plugins/PowhegHooks.h"
30#include "TSystem.h"
33#include <filesystem>
35
36#include <iostream>
37#include <unordered_map>
38#include <numeric>
39
40namespace o2
41{
42namespace eventgen
43{
44
46
47/*****************************************************************/
48/*****************************************************************/
49
50// the default construct uses the GeneratorPythia8Param singleton to extract a config and delegates
51// to the proper constructor
53{
54 LOG(info) << "GeneratorPythia8 constructed from GeneratorPythia8Param ConfigurableParam";
55}
56
57/*****************************************************************/
58
59GeneratorPythia8::GeneratorPythia8(Pythia8GenConfig const& config) : Generator("ALICEo2", "ALICEo2 Pythia8 Generator")
60{
64
65 mInterface = reinterpret_cast<void*>(&mPythia);
66 mInterfaceName = "pythia8";
67
68 LOG(info) << "Instance \'Pythia8\' generator with following parameters";
69 LOG(info) << "config: " << config.config;
70 LOG(info) << "hooksFileName: " << config.hooksFileName;
71 LOG(info) << "hooksFuncName: " << config.hooksFuncName;
72
73 mGenConfig = config;
74
78}
79
80/*****************************************************************/
81
82GeneratorPythia8::GeneratorPythia8(const Char_t* name, const Char_t* title) : Generator(name, title)
83{
86 mInterface = reinterpret_cast<void*>(&mPythia);
87 mInterfaceName = "pythia8";
88}
89
91{
92 // check first of all if Init not yet called and seed not <0
93 if (mIsInitialized) {
94 return false;
95 }
96 if (seed < 0) {
97 return false;
98 }
99 // sets the initial seed and applies the correct Pythia
100 // range
101 mInitialRNGSeed = seed % (MAX_SEED + 1);
102 LOG(info) << "GeneratorPythia8: Setting initial seed to " << mInitialRNGSeed;
103 return true;
104}
105
106/*****************************************************************/
108{
113
116
117 auto seed = mInitialRNGSeed;
118 if (seed == -1) {
119 // Will use the mInitialRNGSeed if it was set.
120 // Otherwise will seed the generator with the state of
121 // TRandom::GetSeed. This is the seed that is influenced from
122 // SimConfig --seed command line options options.
123 seed = gRandom->TRandom::GetSeed(); // this uses the "original" seed
124 // we advance the seed by one so that the next Pythia8 generator gets a different value
125 if (mThisPythia8InstanceID > 0) {
126 gRandom->Rndm();
127 LOG(info) << "Multiple Pythia8 generator instances detected .. automatically adjusting seed further to avoid overlap ";
128 seed = seed ^ gRandom->GetSeed(); // this uses the "current" seed
129 }
130 // apply max seed cuttof
131 seed = seed % (MAX_SEED + 1);
132 LOG(info) << "GeneratorPythia8: Using random seed from gRandom % 900000001: " << seed;
133 }
134 mPythia.readString("Random:setSeed on");
135 mPythia.readString("Random:seed " + std::to_string(seed));
136}
137
138/*****************************************************************/
139
141{
146
152
154 if (!mConfig.empty()) {
155 std::stringstream ss(mConfig);
156 std::string config;
157 while (getline(ss, config, ' ')) {
158 config = gSystem->ExpandPathName(config.c_str());
159 LOG(info) << "Reading configuration from file: " << config;
160 if (!mPythia.readFile(config, true)) {
161 LOG(fatal) << "Failed to init \'Pythia8\': problems with configuration file "
162 << config;
163 return false;
164 }
165 }
166 }
167
169 if (!mHooksFileName.empty()) {
170 LOG(info) << "Applying \'Pythia8\' user hooks: " << mHooksFileName << " -> " << mHooksFuncName;
171 auto hooks = o2::conf::GetFromMacro<Pythia8::UserHooks*>(mHooksFileName, mHooksFuncName, "Pythia8::UserHooks*", "pythia8_user_hooks");
172 if (!hooks) {
173 LOG(fatal) << "Failed to init \'Pythia8\': problem with user hooks configuration ";
174 return false;
175 }
176 setUserHooks(hooks);
177 }
178
179#if PYTHIA_VERSION_INTEGER < 8300
188 mPythia.readString("HadronLevel:Decay off");
189#endif
190 if (mPythia.settings.mode("Beams:frameType") == 4) {
191 // Hook for POWHEG
192 // Read in key POWHEG merging settings
193 int vetoMode = mPythia.settings.mode("POWHEG:veto");
194 int MPIvetoMode = mPythia.settings.mode("POWHEG:MPIveto");
195 bool loadHooks = (vetoMode > 0 || MPIvetoMode > 0);
196 // Add in user hooks for shower vetoing
197 std::shared_ptr<Pythia8::PowhegHooks> powhegHooks;
198 if (loadHooks) {
199 // Set ISR and FSR to start at the kinematical limit
200 if (vetoMode > 0) {
201 mPythia.readString("SpaceShower:pTmaxMatch = 2");
202 mPythia.readString("TimeShower:pTmaxMatch = 2");
203 }
204 // Set MPI to start at the kinematical limit
205 if (MPIvetoMode > 0) {
206 mPythia.readString("MultipartonInteractions:pTmaxMatch = 2");
207 }
208 powhegHooks = std::make_shared<Pythia8::PowhegHooks>();
209 mPythia.setUserHooksPtr((Pythia8::UserHooksPtr)powhegHooks);
210 }
211 }
213 if (!mPythia.init()) {
214 LOG(fatal) << "Failed to init \'Pythia8\': init returned with error";
215 return false;
216 }
217
219
220 mIsInitialized = true;
221
223 return true;
224}
225
226/*****************************************************************/
227void GeneratorPythia8::setUserHooks(Pythia8::UserHooks* hooks)
228{
229#if PYTHIA_VERSION_INTEGER < 8300
230 mPythia.setUserHooksPtr(hooks);
231#else
232 mPythia.setUserHooksPtr(std::shared_ptr<Pythia8::UserHooks>(hooks));
233#endif
234}
235
236/*****************************************************************/
237
238Bool_t
240{
242 if (!mPythia.next()) {
243 return false;
244 }
245
246#if PYTHIA_VERSION_INTEGER < 8300
261 auto nParticles = mPythia.event.size();
262 for (int iparticle = 0; iparticle < nParticles; iparticle++) {
263 auto& aParticle = mPythia.event[iparticle];
264 aParticle.xProd(0.);
265 aParticle.yProd(0.);
266 aParticle.zProd(0.);
267 aParticle.tProd(0.);
268 }
269
271 if (!mPythia.moreDecays())
272 return false;
273#endif
274
276 return true;
277}
278
279/*****************************************************************/
281 const std::vector<int>& old2New,
282 size_t index,
283 std::vector<bool>& done,
284 GetRelatives getter,
285 SetRelatives setter,
286 FirstLastRelative firstLast,
287 const std::string& what,
288 const std::string& ind)
289{
290 // Utility to find new index, or -1 if not found
291 auto findNew = [old2New](size_t old) -> int {
292 return old2New[old];
293 };
294 int newIdx = findNew(index);
295 int hepmc = event[index].statusHepMC();
296
297 LOG(debug) << ind
298 << index << " -> "
299 << newIdx << " (" << hepmc << ") ";
300 if (done[index]) {
301 LOG(debug) << ind << " already done";
302 return;
303 }
304
305 // Our list of new relatives
306 using IdList = std::pair<int, int>;
307 constexpr int invalid = 0xFFFFFFF;
308 IdList newRelatives = std::make_pair(invalid, -invalid);
309
310 // Utility to add id
311 auto addId = [](IdList& l, size_t id) {
312 l.first = std::min(int(id), l.first);
313 l.second = std::max(int(id), l.second);
314 };
315
316 // Get particle and relatives
317 auto& particle = event[index];
318 auto relatives = getter(particle);
319
320 LOG(debug) << ind << " Check " << what << "s ["
321 << std::setw(3) << firstLast(particle).first << ","
322 << std::setw(3) << firstLast(particle).second << "] "
323 << relatives.size();
324
325 for (auto relativeIdx : relatives) {
326 int newRelative = findNew(relativeIdx);
327 if (newRelative >= 0) {
328 // If this relative is to be kept, then append to list of new
329 // relatives.
330 LOG(debug) << ind << " "
331 << what << " "
332 << relativeIdx << " -> "
333 << newRelative << " to be kept" << std::endl;
334 addId(newRelatives, newRelative);
335 continue;
336 }
337 LOG(debug) << ind << " "
338 << what << " "
339 << relativeIdx << " not to be kept "
340 << (done[relativeIdx] ? "already done" : "to be done")
341 << std::endl;
342
343 // Below is code for when the relative is not to be kept
344 auto& relative = event[relativeIdx];
345 if (not done[relativeIdx]) {
346 // IF the relative hasn't been processed yet, do so now
348 old2New, // Map from old to new
349 relativeIdx, // current particle index
350 done, // cache flag
351 getter, // get mother relatives
352 setter, // set mother relatives
353 firstLast, // get first and last
354 what, // what we're looking at
355 ind + " "); // Logging indent
356 }
357
358 // If this relative was already done, then get its relatives and
359 // add them to the list of new relatives.
360 auto grandRelatives = firstLast(relative);
361 int grandRelative1 = grandRelatives.first;
362 int grandRelative2 = grandRelatives.second;
363 assert(grandRelative1 != invalid);
364 assert(grandRelative2 != -invalid);
365 if (grandRelative1 > 0) {
366 addId(newRelatives, grandRelative1);
367 }
368 if (grandRelative2 > 0) {
369 addId(newRelatives, grandRelative2);
370 }
371 LOG(debug) << ind << " "
372 << what << " "
373 << relativeIdx << " gave new relatives "
374 << grandRelative1 << " -> " << grandRelative2;
375 }
376 LOG(debug) << ind << " Got "
377 << (newRelatives.second - newRelatives.first + 1) << " new "
378 << what << "s ";
379
380 if (newRelatives.first != invalid) {
381 // If the first relative is not invalid, then the second isn't
382 // either (possibly the same though).
383 int newRelative1 = newRelatives.first;
384 int newRelative2 = newRelatives.second;
385 setter(particle, newRelative1, newRelative2);
386 LOG(debug) << ind << " " << what << "s: "
387 << firstLast(particle).first << " ("
388 << newRelative1 << "),"
389 << firstLast(particle).second << " ("
390 << newRelative2 << ")";
391
392 } else {
393 setter(particle, 0, 0);
394 }
395 done[index] = true;
396}
397
398/*****************************************************************/
400{
401 // Mapping from old to new index.
402 std::vector<int> old2new(event.size(), -1);
403
404 // Particle 0 is a system particle, and we will skip that in the
405 // following.
406 size_t newId = 0;
407
408 // Loop over particles and store those we need
409 for (size_t i = 1; i < event.size(); i++) {
410 auto& particle = event[i];
411 if (select(particle)) {
412 ++newId;
413 old2new[i] = newId;
414 }
415 }
416 // Utility to find new index, or -1 if not found
417 auto findNew = [old2new](size_t old) -> int {
418 return old2new[old];
419 };
420
421 // First loop, investigate mothers - from the bottom
422 auto getMothers = [](const Pythia8::Particle& particle) { return particle.motherList(); };
423 auto setMothers = [](Pythia8::Particle& particle, int m1, int m2) { particle.mothers(m1, m2); };
424 auto firstLastMothers = [](const Pythia8::Particle& particle) { return std::make_pair(particle.mother1(), particle.mother2()); };
425
426 std::vector<bool> motherDone(event.size(), false);
427 for (size_t i = 1; i < event.size(); ++i) {
429 old2new, // Map from old to new
430 i, // current particle index
431 motherDone, // cache flag
432 getMothers, // get mother relatives
433 setMothers, // set mother relatives
434 firstLastMothers, // get first and last
435 "mother"); // what we're looking at
436 }
437
438 // Second loop, investigate daughters - from the top
439 auto getDaughters = [](const Pythia8::Particle& particle) {
440 // In case of |status|==13 (diffractive), we cannot use
441 // Pythia8::Particle::daughterList as it will give more than
442 // just the immediate daughters. In that cae, we do it
443 // ourselves.
444 if (std::abs(particle.status()) == 13) {
445 int d1 = particle.daughter1();
446 int d2 = particle.daughter2();
447 if (d1 == 0 and d2 == 0) {
448 return std::vector<int>();
449 }
450 if (d2 == 0) {
451 return std::vector<int>{d1};
452 }
453 if (d2 > d1) {
454 std::vector<int> ret(d2-d1+1);
455 std::iota(ret.begin(), ret.end(), d1);
456 return ret;
457 }
458 return std::vector<int>{d2,d1};
459 }
460 return particle.daughterList(); };
461 auto setDaughters = [](Pythia8::Particle& particle, int d1, int d2) { particle.daughters(d1, d2); };
462 auto firstLastDaughters = [](const Pythia8::Particle& particle) { return std::make_pair(particle.daughter1(), particle.daughter2()); };
463
464 std::vector<bool> daughterDone(event.size(), false);
465 for (size_t i = event.size() - 1; i > 0; --i) {
467 old2new, // Map from old to new
468 i, // current particle index
469 daughterDone, // cache flag
470 getDaughters, // get mother relatives
471 setDaughters, // set mother relatives
472 firstLastDaughters, // get first and last
473 "daughter"); // what we're looking at
474 }
475
476 // Make a pruned event
477 Pythia8::Event pruned;
478 pruned.init("Pruned event", &mPythia.particleData);
479 pruned.reset();
480
481 for (size_t i = 1; i < event.size(); i++) {
482 int newIdx = findNew(i);
483 if (newIdx < 0) {
484 continue;
485 }
486
487 auto particle = event[i];
488 int realIdx = pruned.append(particle);
489 assert(realIdx == newIdx);
490 }
491
492 // We may have that two or more mothers share some daughters, but
493 // that one or more mothers have more daughters than the other
494 // mothers, and hence not all daughters point back to all mothers.
495 // This can happen, for example, if a beam particle radiates
496 // on-shell particles before an interaction with any daughters
497 // from the other mothers. Thus, we need to take care of that or
498 // the event record will be corrupted.
499 //
500 // What we do is that for all particles, we look up the daughters.
501 // Then for each daughter, we check the mothers of those
502 // daughters. If this list of mothers include other mothers than
503 // the currently investigated mother, we must change the mothers
504 // of the currently investigated daughters.
505 using IdList = std::pair<int, int>;
506 // Utility to add id
507 auto addId = [](IdList& l, size_t id) {
508 l.first = std::min(int(id), l.first);
509 l.second = std::max(int(id), l.second);
510 };
511 constexpr int invalid = 0xFFFFFFF;
512
513 std::vector<bool> shareDone(pruned.size(), false);
514 for (size_t i = 1; i < pruned.size(); i++) {
515 if (shareDone[i]) {
516 continue;
517 }
518
519 auto& particle = pruned[i];
520 auto daughters = particle.daughterList();
521 IdList allDaughters = std::make_pair(invalid, -invalid);
522 IdList allMothers = std::make_pair(invalid, -invalid);
523 addId(allMothers, i);
524 for (auto daughterIdx : daughters) {
525 // Add this daughter to set of all daughters
526 addId(allDaughters, daughterIdx);
527 auto& daughter = pruned[daughterIdx];
528 auto otherMothers = daughter.motherList();
529 for (auto otherMotherIdx : otherMothers) {
530 // Add this mother to set of all mothers. That is, take all
531 // mothers of the current daughter of the current particle
532 // and store that. In this way, we register mothers that
533 // share a daughter with the current particle.
534 addId(allMothers, otherMotherIdx);
535 // We also need to take all the daughters of this shared
536 // mother and reister those.
537 auto& otherMother = pruned[otherMotherIdx];
538 int otherDaughter1 = otherMother.daughter1();
539 int otherDaughter2 = otherMother.daughter2();
540 if (otherDaughter1 > 0) {
541 addId(allDaughters, otherDaughter1);
542 }
543 if (otherDaughter2 > 0) {
544 addId(allDaughters, otherDaughter2);
545 }
546 }
547 // At this point, we have added all mothers of current
548 // daughter, and all daughters of those mothers.
549 }
550 // At this point, we have all mothers that share daughters with
551 // the current particle, and we have all of the daughters
552 // too.
553 //
554 // We can now update the daughter information on all mothers
555 int minDaughter = allDaughters.first;
556 int maxDaughter = allDaughters.second;
557 int minMother = allMothers.first;
558 int maxMother = allMothers.second;
559 if (minMother != invalid) {
560 // If first mother isn't invalid, then second isn't either
561 for (size_t motherIdx = minMother; motherIdx <= maxMother; //
562 motherIdx++) {
563 shareDone[motherIdx] = true;
564 if (minDaughter == invalid) {
565 pruned[motherIdx].daughters(0, 0);
566 } else {
567 pruned[motherIdx].daughters(minDaughter, maxDaughter);
568 }
569 }
570 }
571 if (minDaughter != invalid) {
572 // If least mother isn't invalid, then largest mother will not
573 // be invalid either.
574 for (size_t daughterIdx = minDaughter; daughterIdx <= maxDaughter; //
575 daughterIdx++) {
576 if (minMother == invalid) {
577 pruned[daughterIdx].mothers(0, 0);
578 } else {
579 pruned[daughterIdx].mothers(minMother, maxMother);
580 }
581 }
582 }
583 }
584 if (mGenConfig.verbose) {
585 LOG(info) << "Pythia event was pruned from " << event.size()
586 << " to " << pruned.size() << " particles";
587 }
588 // Assign our pruned event to the event passed in
589 event = pruned;
590}
591
592/*****************************************************************/
594{
595 mUserFilterFcn = [](Pythia8::Particle const&) -> bool { return true; };
596
597 std::string filter = mGenConfig.particleFilter;
598 if (filter.size() > 0) {
599 LOG(info) << "Initializing the callback for user-based particle pruning " << filter;
600 auto expandedFileName = o2::utils::expandShellVarsInFileName(filter);
601 if (std::filesystem::exists(expandedFileName)) {
602 // if the filter is in a file we will compile the hook on the fly
603 mUserFilterFcn = o2::conf::GetFromMacro<UserFilterFcn>(expandedFileName, "filterPythia()", "o2::eventgen::GeneratorPythia8::UserFilterFcn", "o2mc_pythia8_userfilter_hook");
604 LOG(info) << "Hook initialized from file " << expandedFileName;
605 } else {
606 // if it's not a file we interpret it as a C++ lambda string and JIT it directly;
607 LOG(error) << "Did not find a file " << expandedFileName << " ; Will not execute hook";
608 }
609 mApplyPruning = true;
610 }
611}
612
613/*****************************************************************/
614
615Bool_t
617{
620 // The right moment to filter out unwanted stuff (like parton-level
621 // event information) Here, we aim to filter out everything before
622 // hadronization with the motivation to reduce the size of the MC
623 // event record in the AOD.
624
625 std::function<bool(const Pythia8::Particle&)> partonSelect = [](const Pythia8::Particle&) { return true; };
626 bool includeParton = mGenConfig.includePartonEvent;
627 if (not includeParton) {
628
629 // Select pythia particles
630 partonSelect = [](const Pythia8::Particle& particle) {
631 switch (particle.statusHepMC()) {
632 case 1: // Final st
633 case 2: // Decayed
634 case 4: // Beam
635 return true;
636 }
637 // For example to keep diffractive particles
638 // if (particle.id() == 9902210) return true;
639 return false;
640 };
641 mApplyPruning = true;
642 }
643
644 if (mApplyPruning) {
645 auto finalSelect = [partonSelect, this](const Pythia8::Particle& p) { return partonSelect(p) && mUserFilterFcn(p); };
646 pruneEvent(event, finalSelect);
647 }
648
649 /* loop over particles */
650 auto nParticles = event.size();
651 for (Int_t iparticle = 1; iparticle < nParticles; iparticle++) {
652 // first particle is system
653 auto particle = event[iparticle];
654 auto pdg = particle.id();
655 auto st = o2::mcgenstatus::MCGenStatusEncoding(particle.statusHepMC(), //
656 particle.status()) //
658 mParticles.push_back(TParticle(particle.id(), // Particle type
659 st, // status
660 particle.mother1() - 1, // first mother
661 particle.mother2() - 1, // second mother
662 particle.daughter1() - 1, // first daughter
663 particle.daughter2() - 1, // second daughter
664 particle.px(), // X-momentum
665 particle.py(), // Y-momentum
666 particle.pz(), // Z-momentum
667 particle.e(), // Energy
668 particle.xProd(), // Production X
669 particle.yProd(), // Production Y
670 particle.zProd(), // Production Z
671 particle.tProd())); // Production t
672 mParticles.back().SetBit(ParticleStatus::kToBeDone, //
673 particle.statusHepMC() == 1);
674 }
675
677 return kTRUE;
678}
679
680/*****************************************************************/
681
683{
686
687 eventHeader->putInfo<std::string>(Key::generator, "pythia8");
688 eventHeader->putInfo<int>(Key::generatorVersion, PYTHIA_VERSION_INTEGER);
689 eventHeader->putInfo<std::string>(Key::processName, mPythia.info.name());
690 eventHeader->putInfo<int>(Key::processCode, mPythia.info.code());
691 eventHeader->putInfo<float>(Key::weight, mPythia.info.weight());
692
693 auto& info = mPythia.info;
694
695 eventHeader->putInfo<int>(Key::acceptedEvents, info.nAccepted());
696 eventHeader->putInfo<int>(Key::attemptedEvents, info.nTried());
697
698 // Set PDF information
699 eventHeader->putInfo<int>(Key::pdfParton1Id, info.id1pdf());
700 eventHeader->putInfo<int>(Key::pdfParton2Id, info.id2pdf());
701 eventHeader->putInfo<float>(Key::pdfX1, info.x1pdf());
702 eventHeader->putInfo<float>(Key::pdfX2, info.x2pdf());
703 eventHeader->putInfo<float>(Key::pdfScale, info.QFac());
704 eventHeader->putInfo<float>(Key::pdfXF1, info.pdf1());
705 eventHeader->putInfo<float>(Key::pdfXF2, info.pdf2());
706
707 // Set cross section
708 eventHeader->putInfo<float>(Key::xSection, info.sigmaGen() * 1e9);
709 eventHeader->putInfo<float>(Key::xSectionError, info.sigmaErr() * 1e9);
710
711 // Set event scale and nMPI
712 eventHeader->putInfo<float>(Key::eventScale, info.QRen());
713 eventHeader->putInfo<int>(Key::mpi, info.nMPI());
714
715 // Set weights (overrides cross-section for each weight)
716 size_t iw = 0;
717 auto xsecErr = info.weightContainerPtr->getTotalXsecErr();
718 for (auto w : info.weightContainerPtr->getTotalXsec()) {
719 std::string post = (iw == 0 ? "" : "_" + std::to_string(iw));
720 eventHeader->putInfo<float>(Key::weight + post, info.weightValueByIndex(iw));
721 eventHeader->putInfo<float>(Key::xSection + post, w * 1e9);
722 eventHeader->putInfo<float>(Key::xSectionError + post, xsecErr[iw] * 1e9);
723 iw++;
724 }
725
726#if PYTHIA_VERSION_INTEGER < 8300
727 auto hiinfo = mPythia.info.hiinfo;
728#else
729 auto hiinfo = mPythia.info.hiInfo;
730#endif
731
732 if (hiinfo) {
734 eventHeader->SetB(hiinfo->b());
735 eventHeader->putInfo<double>(Key::impactParameter, hiinfo->b());
737#if PYTHIA_VERSION_INTEGER >= 8310
738 eventHeader->putInfo<double>(Key::planeAngle, hiinfo->phi());
739#endif
740 auto bImp = hiinfo->b();
742 int nColl, nPart;
743 int nPartProtonProj, nPartNeutronProj, nPartProtonTarg, nPartNeutronTarg;
744 int nRemnProtonProj, nRemnNeutronProj, nRemnProtonTarg, nRemnNeutronTarg;
745 int nFreeNeutronProj, nFreeProtonProj, nFreeNeutronTarg, nFreeProtonTarg;
746 getNcoll(nColl);
747 getNpart(nPart);
748 getNpart(nPartProtonProj, nPartNeutronProj, nPartProtonTarg, nPartNeutronTarg);
749 getNremn(nRemnProtonProj, nRemnNeutronProj, nRemnProtonTarg, nRemnNeutronTarg);
750 getNfreeSpec(nFreeNeutronProj, nFreeProtonProj, nFreeNeutronTarg, nFreeProtonTarg);
751 eventHeader->putInfo<int>(Key::nColl, nColl);
752 // These are all non-HepMC3 fields - of limited use
753 eventHeader->putInfo<int>("Npart", nPart);
754 eventHeader->putInfo<int>("Npart_proj_p", nPartProtonProj);
755 eventHeader->putInfo<int>("Npart_proj_n", nPartNeutronProj);
756 eventHeader->putInfo<int>("Npart_targ_p", nPartProtonTarg);
757 eventHeader->putInfo<int>("Npart_targ_n", nPartNeutronTarg);
758 eventHeader->putInfo<int>("Nremn_proj_p", nRemnProtonProj);
759 eventHeader->putInfo<int>("Nremn_proj_n", nRemnNeutronProj);
760 eventHeader->putInfo<int>("Nremn_targ_p", nRemnProtonTarg);
761 eventHeader->putInfo<int>("Nremn_targ_n", nRemnNeutronTarg);
762 eventHeader->putInfo<int>("Nfree_proj_n", nFreeNeutronProj);
763 eventHeader->putInfo<int>("Nfree_proj_p", nFreeProtonProj);
764 eventHeader->putInfo<int>("Nfree_targ_n", nFreeNeutronTarg);
765 eventHeader->putInfo<int>("Nfree_targ_p", nFreeProtonTarg);
766
767 // --- HepMC3 conforming information ---
768 // This is how the Pythia authors define Ncoll
769 // eventHeader->putInfo<int>(Key::nColl,
770 // hiinfo->nAbsProj() + hiinfo->nDiffProj() +
771 // hiinfo->nAbsTarg() + hiinfo->nDiffTarg() -
772 // hiiinfo->nCollND() - hiinfo->nCollDD());
773 eventHeader->putInfo<int>(Key::nPartProjectile,
774 hiinfo->nAbsProj() + hiinfo->nDiffProj());
775 eventHeader->putInfo<int>(Key::nPartTarget,
776 hiinfo->nAbsTarg() + hiinfo->nDiffTarg());
777 eventHeader->putInfo<int>(Key::nCollHard, hiinfo->nCollNDTot());
778 }
779}
780
781/*****************************************************************/
782
783void GeneratorPythia8::selectFromAncestor(int ancestor, Pythia8::Event& inputEvent, Pythia8::Event& outputEvent)
784{
785
790 // recursive selection via lambda function
791 std::set<int> selected;
792 std::function<void(int)> select;
793 select = [&](int i) {
794 selected.insert(i);
795 auto dl = inputEvent[i].daughterList();
796 for (auto j : dl) {
797 select(j);
798 }
799 };
800 select(ancestor);
801
802 // map selected particle index to output index
803 std::map<int, int> indexMap;
804 int index = outputEvent.size();
805 for (auto i : selected) {
806 indexMap[i] = index++;
807 }
808
809 // adjust mother/daughter indices and append to output event
810 for (auto i : selected) {
811 auto p = mPythia.event[i];
812 auto m1 = indexMap[p.mother1()];
813 auto m2 = indexMap[p.mother2()];
814 auto d1 = indexMap[p.daughter1()];
815 auto d2 = indexMap[p.daughter2()];
816 p.mothers(m1, m2);
817 p.daughters(d1, d2);
818
819 outputEvent.append(p);
820 }
821}
822
823/*****************************************************************/
824
825void GeneratorPythia8::getNcoll(const Pythia8::Info& info, int& nColl)
826{
827
830#if PYTHIA_VERSION_INTEGER < 8300
831 auto hiinfo = info.hiinfo;
832#else
833 auto hiinfo = info.hiInfo;
834#endif
835 nColl = 0;
836 if (!hiinfo) {
837 LOG(warn) << "No heavy-ion information from Pythia";
838 return;
839 }
840
841 // This is how the Pythia authors define Ncoll
842 nColl = (hiinfo->nAbsProj() + hiinfo->nDiffProj() +
843 hiinfo->nAbsTarg() + hiinfo->nDiffTarg() -
844 hiinfo->nCollND() - hiinfo->nCollDD());
845
846 if (not hiinfo->subCollisionsPtr()) {
847#if PYTHIA_VERSION_INTEGER < 8310
848 LOG(fatal) << "No sub-collision pointer from Pythia";
849#endif
850 return;
851 }
852
853 // loop over sub-collisions
854 auto scptr = hiinfo->subCollisionsPtr();
855 for (auto sc : *scptr) {
856
857 // wounded nucleon flag in projectile/target
858 auto pW = sc.proj->status() == Pythia8::Nucleon::ABS; // according to C.Bierlich this should be == Nucleon::ABS
859 auto tW = sc.targ->status() == Pythia8::Nucleon::ABS;
860
861 // increase number of collisions if both are wounded
862 if (pW && tW) {
863 nColl++;
864 }
865 }
866}
867
868/*****************************************************************/
869
870void GeneratorPythia8::getNpart(const Pythia8::Info& info, int& nPart)
871{
872
875#if PYTHIA_VERSION_INTEGER < 8300
876 auto hiinfo = info.hiinfo;
877#else
878 auto hiinfo = info.hiInfo;
879#endif
880 nPart = 0;
881 if (not hiinfo) {
882 return;
883 }
884
885 // This is how the Pythia authors calculate Npart
886 nPart = (hiinfo->nAbsProj() + hiinfo->nDiffProj() +
887 hiinfo->nAbsTarg() + hiinfo->nDiffTarg());
888 if (not hiinfo->subCollisionsPtr()) {
889 return;
890 }
891
892 int nProtonProj, nNeutronProj, nProtonTarg, nNeutronTarg;
893 getNpart(info, nProtonProj, nNeutronProj, nProtonTarg, nNeutronTarg);
894 nPart = nProtonProj + nNeutronProj + nProtonTarg + nNeutronTarg;
895}
896
897/*****************************************************************/
898
899void GeneratorPythia8::getNpart(const Pythia8::Info& info, int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg)
900{
901
904#if PYTHIA_VERSION_INTEGER < 8300
905 auto hiinfo = info.hiinfo;
906#else
907 auto hiinfo = info.hiInfo;
908#endif
909
910 nProtonProj = nNeutronProj = nProtonTarg = nNeutronTarg = 0;
911 if (!hiinfo) {
912 return;
913 }
914
915 nProtonProj = hiinfo->nAbsProj() + hiinfo->nDiffProj();
916 nProtonTarg = hiinfo->nAbsTarg() + hiinfo->nDiffTarg();
917 if (not hiinfo->subCollisionsPtr()) {
918#if PYTHIA_VERSION_INTEGER < 8310
919 LOG(fatal) << "No sub-collision pointer from Pythia";
920#endif
921 return;
922 }
923
924 // keep track of wounded nucleons
925 std::vector<Pythia8::Nucleon*> projW;
926 std::vector<Pythia8::Nucleon*> targW;
927
928 // loop over sub-collisions
929 auto scptr = hiinfo->subCollisionsPtr();
930 for (auto sc : *scptr) {
931
932 // wounded nucleon flag in projectile/target
933 auto pW = sc.proj->status() == Pythia8::Nucleon::ABS || sc.proj->status() == Pythia8::Nucleon::DIFF; // according to C.Bierlich this should be == Nucleon::ABS || Nucleon::DIFF
934 auto tW = sc.targ->status() == Pythia8::Nucleon::ABS || sc.targ->status() == Pythia8::Nucleon::DIFF;
935
936 // increase number of wounded projectile nucleons if not yet in the wounded vector
937 if (pW && std::find(projW.begin(), projW.end(), sc.proj) == projW.end()) {
938 projW.push_back(sc.proj);
939 if (sc.proj->id() == 2212) {
940 nProtonProj++;
941 } else if (sc.proj->id() == 2112) {
942 nNeutronProj++;
943 }
944 }
945
946 // increase number of wounded target nucleons if not yet in the wounded vector
947 if (tW && std::find(targW.begin(), targW.end(), sc.targ) == targW.end()) {
948 targW.push_back(sc.targ);
949 if (sc.targ->id() == 2212) {
950 nProtonTarg++;
951 } else if (sc.targ->id() == 2112) {
952 nNeutronTarg++;
953 }
954 }
955 }
956}
957
958/*****************************************************************/
959
960void GeneratorPythia8::getNremn(const Pythia8::Event& event, int& nProtonProj, int& nNeutronProj, int& nProtonTarg, int& nNeutronTarg)
961{
962
965 // reset
966 nProtonProj = nNeutronProj = nProtonTarg = nNeutronTarg = 0;
967 auto nNucRem = 0;
968
969 // particle loop
970 auto nparticles = event.size();
971 for (int ipa = 0; ipa < nparticles; ++ipa) {
972 const auto particle = event[ipa];
973 auto pdg = particle.id();
974
975 // nuclear remnants have pdg code = ±10LZZZAAA9
976 if (pdg < 1000000000) {
977 continue; // must be nucleus
978 }
979 if (pdg % 10 != 9) {
980 continue; // first digit must be 9
981 }
982 nNucRem++;
983
984 // extract A, Z and L from pdg code
985 pdg /= 10;
986 auto A = pdg % 1000;
987 pdg /= 1000;
988 auto Z = pdg % 1000;
989 pdg /= 1000;
990 auto L = pdg % 10;
991
992 if (particle.pz() > 0.) {
993 nProtonProj = Z;
994 nNeutronProj = A - Z;
995 }
996 if (particle.pz() < 0.) {
997 nProtonTarg = Z;
998 nNeutronTarg = A - Z;
999 }
1000
1001 } // end of particle loop
1002
1003 if (nNucRem > 2) {
1004 LOG(warning) << " GeneratorPythia8: found more than two nuclear remnants (weird)";
1005 }
1006}
1007/*****************************************************************/
1008
1009/*****************************************************************/
1010
1011void GeneratorPythia8::getNfreeSpec(const Pythia8::Info& info, int& nFreenProj, int& nFreepProj, int& nFreenTarg, int& nFreepTarg)
1012{
1015#if PYTHIA_VERSION_INTEGER < 8300
1016 auto hiinfo = info.hiinfo;
1017#else
1018 auto hiinfo = info.hiInfo;
1019#endif
1020
1021 if (!hiinfo) {
1022 return;
1023 }
1024
1025 double b = hiinfo->b();
1026
1027 static o2::zdc::FragmentParam frag; // data-driven model to get free spectators given impact parameter
1028
1029 TF1 const& fneutrons = frag.getfNeutrons();
1030 TF1 const& fsigman = frag.getsigmaNeutrons();
1031 TF1 const& fprotons = frag.getfProtons();
1032 TF1 const& fsigmap = frag.getsigmaProtons();
1033
1034 // Calculating no. of free spectators from parametrization
1035 int nneu[2] = {0, 0};
1036 for (int i = 0; i < 2; i++) {
1037 float nave = fneutrons.Eval(b);
1038 float sigman = fsigman.Eval(b);
1039 float nfree = gRandom->Gaus(nave, 0.68 * sigman * nave);
1040 nneu[i] = (int)nfree;
1041 if (nave < 0 || nneu[i] < 0) {
1042 nneu[i] = 0;
1043 }
1044 if (nneu[i] > 126) {
1045 nneu[i] = 126;
1046 }
1047 }
1048 //
1049 int npro[2] = {0, 0};
1050 for (int i = 0; i < 2; i++) {
1051 float pave = fprotons.Eval(b);
1052 float sigmap = fsigman.Eval(b);
1053 float pfree = gRandom->Gaus(pave, 0.68 * sigmap * pave) / 0.7;
1054 npro[i] = (int)pfree;
1055 if (pave < 0 || npro[i] < 0) {
1056 npro[i] = 0;
1057 }
1058 if (npro[i] > 82) {
1059 npro[i] = 82;
1060 }
1061 }
1062
1063 nFreenProj = nneu[0];
1064 nFreenTarg = nneu[1];
1065 nFreepProj = npro[0];
1066 nFreepTarg = npro[1];
1067 /*****************************************************************/
1068}
1069
1070} /* namespace eventgen */
1071} /* namespace o2 */
o2::monitoring::tags::Key Key
int32_t i
bool done
@ kToBeDone
uint32_t j
Definition RawData.h:0
std::ostringstream debug
benchmark::State & st
Definition A.h:16
void putInfo(std::string const &key, T const &value)
void setHooksFuncName(std::string val)
void selectFromAncestor(int ancestor, Pythia8::Event &inputEvent, Pythia8::Event &outputEvent)
void pruneEvent(Pythia8::Event &event, Select select)
static std::atomic< int > Pythia8InstanceCounter
void setHooksFileName(std::string val)
void investigateRelatives(Pythia8::Event &event, const std::vector< int > &old2New, size_t index, std::vector< bool > &done, GetRelatives getter, SetRelatives setter, FirstLastRelative firstLast, const std::string &what, const std::string &ind="")
void getNremn(int &nProtonProj, int &nNeutronProj, int &nProtonTarg, int &nNeutronTarg)
void seedGenerator()
performs seeding of the random state of Pythia (called from Init)
void setUserHooks(Pythia8::UserHooks *hooks)
void updateHeader(o2::dataformats::MCEventHeader *eventHeader) override
void getNfreeSpec(int &nFreenProj, int &nFreepProj, int &nFreenTarg, int &nFreepTarg)
void setConfig(std::string val)
std::string mInterfaceName
Definition Generator.h:120
std::vector< TParticle > mParticles
Definition Generator.h:140
Bool_t Init() override
Definition Generator.cxx:57
TF1 const & getfProtons() const
TF1 const & getfNeutrons() const
TF1 const & getsigmaNeutrons() const
TF1 const & getsigmaProtons() const
struct _cl_event * event
Definition glcorearb.h:2982
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLint GLint GLint GLint GLint GLint GLint GLbitfield GLenum filter
Definition glcorearb.h:1308
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLuint id
Definition glcorearb.h:650
std::vector< InputSpec > select(char const *matcher="")
std::string expandShellVarsInFileName(std::string const &input)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"