Hi Andre,
since you are on 6.11, you can take advantage of TDataFrame. Basically, TDataFrame is a tool that allows you to formulate your analysis or MC generation in a declarative manner. TDataFrame is rather powerful: it aims to take away cumbersome tasks from users and get out of their way. It allows also to seamlessly write parallelised code. In the following example, I provide a recipe that will allow you:
- Download Pythia8
- Build a dictionary for the Event class
- Build a program to create a generated dataset in root format containing pythia events. Two parameters are available: number of cpus and total number of events to be generated.
The program tdfGeneration.cpp: As you can see the overhead of TDataFrame in terms of lines of code is minimal:
#include <ROOT/TDataFrame.hxx>
#include <Pythia8/Pythia.h>
int main(int argc, char** argv){
int nworkers = 1;
unsigned nevents = 1000;
if (argc == 2) {
nworkers = std::atoi(argv[1]);
}
if (argc == 3) {
nworkers = std::atoi(argv[1]);
nevents = std::atoi(argv[2]);
}
// Generate the Pythias. The random seed is different for each instance
Pythia8::Pythia pythias[nworkers];
auto seed = 1;
for (auto&& pythia : pythias) {
pythia.readString("Beams:idA = 2212"); // Specifies Proton Beam
pythia.readString("Beams:idB = 2212"); // Specifies Proton Beam
pythia.readString("HardQCD:all = on");
pythia.readString("PhaseSpace:pTHatMin = 80.");
pythia.readString("Beams:eCM = 7000.");
pythia.readString("Random:setSeed = on");
pythia.readString("Random:seed = " + std::to_string(seed++));
pythia.init();
}
// go parallel if needed
if (nworkers!=1) ROOT::EnableImplicitMT(nworkers);
// The "generator function"
auto genFunc = [&](unsigned int slot) {
while (!pythias[slot].next()) continue;
return &pythias[slot].event;
};
// Use the TDF to write out
ROOT::Experimental::TDataFrame tdf(nevents);
tdf.DefineSlot("event", genFunc)
.Snapshot<Pythia8::Event*>("tree", "hardQCD.root", {"event"});
}
Now, the instruction to get Pythia8:
wget http://home.thep.lu.se/~torbjorn/pythia8/pythia8226.tgz
tar -zxf pythia8226.tgz
mkdir pythia
cd pythia8226
./configure --with-root=$ROOTSYS --enable-shared --prefix=`pwd`/../pythia
make -j 6
make install
cd ..
And now, the instruction to setup the environment and build the dictionary:
export LD_LIBRARY_PATH=pythia/lib/:$LD_LIBRARY_PATH
export ROOT_INCLUDE_PATH=pythia/include
rootcling -f pythiaDict.cc -rml libpythiaDict.so -rmf libpythiaDict.rootmap pythia8226/examples/main92.h ./pythia8226/examples/main92LinkDef.h
g++ -fPIC -shared -o libpythiaDict.so pythiaDict.cc -I pythia/include/ -L pythia/lib/ -l pythia8 `root-config --cflags --libs`
g++ tdfGeneration.cpp -o tdfGeneration -I pythia/include/ -L pythia/lib/ -l pythia8 `root-config --cflags --libs` -lTreePlayer
Now, to generate 100 pp collisions resulting in hard QCD processes on 8 cpus in parallel, just type:
./tdfGeneration 8 100
This should allow you to profitably use all of the cores at your disposal.
I hope that helps.
Cheers,
D