Skip to content

Commit ebf8670

Browse files
jackal1-66sawenzel
authored andcommitted
Generator example for quick HepMC extraction from Pythia8
1 parent 1e69f5e commit ebf8670

File tree

4 files changed

+80
-0
lines changed

4 files changed

+80
-0
lines changed
Lines changed: 44 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,44 @@
1+
/// \author Marco Giacalone - March 2025
2+
3+
// A simple wrapper and demonstrator around Pythia8 for extracting HepMC3 files.
4+
5+
#include "Pythia8/Pythia.h"
6+
#include "Pythia8Plugins/HepMC3.h"
7+
8+
using namespace o2::eventgen;
9+
10+
class HepMC3_Pythia8Wrapper : public GeneratorPythia8
11+
{
12+
public:
13+
HepMC3_Pythia8Wrapper(std::string filename = "pythia8.hepmc") : GeneratorPythia8(), mFileName(filename)
14+
{
15+
// HepMC conversion object.
16+
mToHepMC = std::make_unique<Pythia8::Pythia8ToHepMC>();
17+
mToHepMC->setNewFile((filename == "" ? "pythia.hepmc" : filename));
18+
};
19+
~HepMC3_Pythia8Wrapper() = default;
20+
21+
bool importParticles() override
22+
{
23+
// events are written after the importParticles step
24+
// since some filtering is happening there
25+
auto ret = GeneratorPythia8::importParticles();
26+
if (ret) {
27+
LOG(info) << "Writing event to HepMC3 format";
28+
mToHepMC->writeNextEvent(mPythia);
29+
}
30+
return ret;
31+
};
32+
33+
private:
34+
std::string mFileName = "pythia8.hepmc";
35+
std::unique_ptr<Pythia8::Pythia8ToHepMC> mToHepMC;
36+
};
37+
38+
FairGenerator*
39+
hepmc_pythia8(std::string filename = "pythia8.hepmc")
40+
{
41+
std::cout << "HepMC3_Pythia8Wrapper initialising with filename: " << filename << std::endl;
42+
auto py8 = new HepMC3_Pythia8Wrapper(filename);
43+
return py8;
44+
}
Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,16 @@
1+
<!-- doxy
2+
\page refrunSimExamplesPythiaHepMCWrapper Example showing easy HepMC extraction using GeneratorPythia8
3+
/doxy -->
4+
5+
This example demonstrates how we can extend GeneratorPythia8 in a user-defined macro (or external generator),
6+
to achieve additional HepMC3 export of generated Pythia8 events.
7+
8+
The example provides a small utility for poeple in need to obtain HepMC files from Pythia8.
9+
Note that many other methods to achieve this are possible (See original Pythia8 example).
10+
11+
The example provides:
12+
13+
- The external generator implementation `Pythia8HepMC3.C`
14+
- a `run.sh` script demonstrating it's usage and a check feeding back the generated hepmc into the simulation
15+
16+
Lines changed: 19 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,19 @@
1+
#!/bin/bash
2+
3+
#
4+
# Script doing Pythia8 event generation and writing these events into HepMC3 files
5+
# (next to generating the usual MCTrack kinematics output).
6+
#
7+
# The script also performs a second event generation based on the generated HepMC3 files.
8+
# In principle it should yield identical kinematics files.
9+
#
10+
11+
NEVENTS=1000
12+
SEED=11
13+
14+
o2-sim -j 1 -g external --configKeyValues 'GeneratorExternal.fileName=Pythia8HepMC3.macro;GeneratorExternal.funcName=hepmc_pythia8("skimmed.hepmc");GeneratorPythia8.config=${O2_ROOT}/share/Generators/egconfig/pythia8_inel.cfg' --seed ${SEED} --noGeant -o pythia8_skimmed -n ${NEVENTS}
15+
o2-sim -j 1 -g external --configKeyValues 'GeneratorExternal.fileName=Pythia8HepMC3.macro;GeneratorExternal.funcName=hepmc_pythia8("unskimmed.hepmc");GeneratorPythia8.config=${O2_ROOT}/share/Generators/egconfig/pythia8_inel.cfg;GeneratorPythia8.includePartonEvent=true' --seed ${SEED} --noGeant -o pythia8_unskimmed -n ${NEVENTS}
16+
17+
# propagate generated hepmc file; it should produce the same kinematics as the original Pythia8
18+
o2-sim -j 1 -g hepmc --configKeyValues="GeneratorFileOrCmd.fileNames=skimmed.hepmc" --vertexMode kNoVertex --noGeant -o fromhepmc_skimmed -n ${NEVENTS} --seed ${SEED}
19+
o2-sim -j 1 -g hepmc --configKeyValues="GeneratorFileOrCmd.fileNames=unskimmed.hepmc" --vertexMode kNoVertex --noGeant -o fromhepmc_unskimmed -n ${NEVENTS} --seed ${SEED}

run/SimExamples/README.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,7 @@
66

77
<!-- doxy
88
* \subpage refrunSimExamplesPythia8
9+
* \subpage refrunSimExamplesPythiaHepMCWrapper
910
* \subpage refrunSimExamplesHF_Embedding_Pythia8
1011
* \subpage refrunSimExamplesSignal_ImpactB
1112
* \subpage refrunSimExamplesTrigger_ImpactB_Pythia8

0 commit comments

Comments
 (0)