11 #include <generators/modules/fragmentation/FragmentationModule.h>
13 #include <generators/evtgen/EvtGenInterface.h>
14 #include <generators/utilities/GeneratorConst.h>
16 #include <framework/gearbox/Unit.h>
17 #include <framework/gearbox/Const.h>
18 #include <framework/particledb/EvtGenDatabasePDG.h>
19 #include <framework/utilities/FileSystem.h>
23 #include <mdst/dataobjects/MCParticleGraph.h>
25 #include <framework/logging/Logger.h>
26 #include <framework/utilities/IOIntercept.h>
32 using namespace Pythia8;
45 setDescription(
"Fragmention of (u/d/s/c) quarks using PYTHIA8");
48 addParam(
"ParameterFile", m_parameterfile,
"Input parameter file for PYTHIA",
49 FileSystem::findFile(
"generators/modules/fragmentation/data/pythia_belle2.dat"));
50 addParam(
"ListPYTHIAEvent", m_listEvent,
"List event record of PYTHIA after hadronization", 0);
51 addParam(
"UseEvtGen", m_useEvtGen,
"Use EvtGen for specific decays", 1);
52 addParam(
"DecFile", m_DecFile,
"EvtGen decay file (DECAY.DEC)",
53 FileSystem::findFile(
"decfiles/dec/DECAY_BELLE2.DEC",
true));
54 addParam(
"UserDecFile", m_UserDecFile,
"User EvtGen decay file", std::string(
""));
55 addParam(
"CoherentMixing", m_coherentMixing,
"Decay the B0-B0bar coherently (should always be true)",
true);
67 FragmentationModule::~FragmentationModule()
72 void FragmentationModule::terminate()
77 statLogCapture.
start();
83 if (nAll) ratio = 100.0 * nGood / nAll;
85 B2WARNING(
"Not all events could be fragmented: " << nAll - nGood <<
" events failed.");
86 B2WARNING(
"Total number of events: " << nAll <<
", of these fragmented: " << nGood <<
", success-ratio (should be >97%): " << ratio
88 B2WARNING(
"Please contact the generator librarian if the success ratio is below 97%.");
89 B2WARNING(
"Please treat the success-ratio as correction of the effective cross section due to unphysical events.");
96 void FragmentationModule::initialize()
98 m_mcparticles.isRequired(m_particleList);
100 B2DEBUG(150,
"Initialize PYTHIA8");
108 initLogCapture.
start();
110 m_Pythia =
new Pythia;
111 m_PythiaEvent = &m_Pythia->event;
112 (*m_PythiaEvent) = 0;
115 EvtGen* evtGen = EvtGenInterface::createEvtGen(m_DecFile, m_coherentMixing);
116 loadEvtGenParticleData(m_Pythia);
119 m_Pythia->readString(
"ProcessLevel:all = off");
122 if (!m_Pythia->readFile(m_parameterfile))
123 B2FATAL(
"Cannot read Pythia parameter file.");
127 m_Pythia->setRndmEnginePtr(fragRndm);
136 B2INFO(
"Using PYTHIA EvtGen Interface");
137 const std::string defaultDecFile = FileSystem::findFile(
"decfiles/dec/DECAY_BELLE2.DEC",
true);
138 if (m_DecFile.empty()) {
139 B2ERROR(
"No global decay file defined, please make sure the parameter 'DecFile' is set correctly");
142 if (defaultDecFile.empty()) {
143 B2WARNING(
"Cannot find default decay file");
144 }
else if (defaultDecFile != m_DecFile) {
145 B2INFO(
"Using non-standard DECAY file \"" << m_DecFile <<
"\"");
148 evtgen->readDecayFile(m_UserDecFile);
159 m_Pythia->settings.listChanged();
167 void FragmentationModule::event()
170 mcParticleGraph.clear();
171 mcParticleGraph.loadList(m_particleList);
175 resetLogCapture.
start();
176 m_PythiaEvent->reset();
185 std::map<int, int> indexPYTHIA;
188 std::map<int, int> indexMCGraph;
191 int quarkPosition = 0;
194 int nPart = m_mcparticles.getEntries();
195 for (
int iPart = 0; iPart < nPart; iPart++) {
196 MCParticle* currParticle = m_mcparticles[iPart];
200 int pythiaIndex = addParticleToPYTHIA(*currParticle);
202 if (pythiaIndex != 0) {
203 indexPYTHIA[nAdded] = iPart;
204 if (pythiaIndex > 0) quarkPosition = iPart;
210 B2FATAL(
"Invalid number of quarks: " << nQuarks <<
" (should be 2)!");
214 B2WARNING(
"No virtual exchange particle given, no PYTHIA FSR in Decays");
217 m_Pythia->forceTimeShower(2, 3, 20.00);
228 eventLogCapture.
start();
229 int success = m_Pythia->next();
234 listLogCapture.
start();
235 m_PythiaEvent->list();
246 decayLogCapture.
start();
253 for (
int iPythiaPart = 0; iPythiaPart < m_Pythia->event.size(); ++iPythiaPart) {
254 auto oldindex = indexPYTHIA.find(iPythiaPart);
257 if (m_Pythia->event[iPythiaPart].id() == 90)
continue;
259 if (oldindex == end(indexPYTHIA)) {
263 int position = mcParticleGraph.size();
264 mcParticleGraph.addParticle();
265 indexMCGraph[iPythiaPart] = position;
271 if (m_Pythia->event[iPythiaPart].statusHepMC() < 1)
continue;
274 p->setPDG(m_Pythia->event[iPythiaPart].id());
277 TLorentzVector p4(m_Pythia->event[iPythiaPart].px(), m_Pythia->event[iPythiaPart].py(), m_Pythia->event[iPythiaPart].pz(),
278 m_Pythia->event[iPythiaPart].e());
280 p->setMass(m_Pythia->event[iPythiaPart].m());
283 p->setProductionVertex(m_Pythia->event[iPythiaPart].xProd() * Unit::mm, m_Pythia->event[iPythiaPart].yProd() * Unit::mm,
284 m_Pythia->event[iPythiaPart].zProd() * Unit::mm);
285 p->setProductionTime(m_Pythia->event[iPythiaPart].tProd() * Unit::mm / Const::speedOfLight);
286 p->setValidVertex(
true);
289 p->addStatus(MCParticleGraph::GraphParticle::c_PrimaryParticle);
292 if (m_Pythia->event[iPythiaPart].status() == 51 && m_Pythia->event[iPythiaPart].id() == 22) {
293 p->addStatus(MCParticleGraph::GraphParticle::c_IsFSRPhoton);
297 if (m_Pythia->event[iPythiaPart].status() == GeneratorConst::FSR_STATUS_CODE && m_Pythia->event[iPythiaPart].id() == 22) {
298 p->addStatus(MCParticleGraph::GraphParticle::c_IsPHOTOSPhoton);
302 if (m_Pythia->event[iPythiaPart].statusHepMC() == 1) {
303 p->addStatus(MCParticleGraph::GraphParticle::c_StableInGenerator);
307 const int motherid = m_Pythia->event[iPythiaPart].mother1();
310 auto motherindex = indexMCGraph.find(motherid);
312 if (motherindex != end(indexMCGraph)) {
313 int motheridingraph = indexMCGraph[motherid];
330 if (m_listEvent) m_PythiaEvent->list();
333 mcParticleGraph.generateList(m_particleList,
334 MCParticleGraph::c_clearParticles | MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
341 int FragmentationModule::addParticleToPYTHIA(
const MCParticle& mcParticle)
344 const int id = mcParticle.
getPDG();
348 bool isQuark =
false;
349 if (abs(
id) >= 1 && abs(
id) <= 5) isQuark =
true;
350 if (
id == 23) isVPho =
true;
352 if (!(isVPho || isQuark))
return 0;
356 B2WARNING(
"Quark already has a daughter!");
361 const double mass = mcParticle.
getMass();
363 const double energy = sqrt(mass * mass + p.Mag2());
367 m_PythiaEvent->append(23, -22, 0, 0, 2, 3, 0, 0, p.X(), p.Y(), p.Z(), energy, mass);
370 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 101, 0, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
373 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 0, 101, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
382 void FragmentationModule::loadEvtGenParticleData(Pythia8::Pythia* pythia)
384 Pythia8::ParticleData* particleData = &pythia->particleData;
385 int nParticles = EvtPDL::entries();
386 for (
int i = 0; i < nParticles; ++i) {
387 EvtId evtgenParticle = EvtPDL::getEntry(i);
393 int pdg = EvtPDL::getStdHep(evtgenParticle);
394 if (pdg <= 10 || (pdg > 20 && pdg <= 100))
396 EvtId evtgenAntiParticle = EvtPDL::chargeConj(evtgenParticle);
397 if (particleData->isParticle(pdg)) {
398 particleData->setAll(pdg,
399 EvtPDL::name(evtgenParticle),
400 EvtPDL::name(evtgenAntiParticle),
401 EvtPDL::getSpinType(evtgenParticle),
402 EvtPDL::chg3(evtgenParticle),
405 EvtPDL::getMeanMass(evtgenParticle),
406 EvtPDL::getWidth(evtgenParticle),
407 EvtPDL::getMinMass(evtgenParticle),
408 EvtPDL::getMaxMass(evtgenParticle),
409 EvtPDL::getctau(evtgenParticle));
417 FragmentationRndm::FragmentationRndm() : Pythia8::RndmEngine()
424 double value = gRandom->Rndm();