9 #include <generators/modules/fragmentation/FragmentationModule.h>
11 #include <generators/evtgen/EvtGenInterface.h>
12 #include <generators/utilities/GeneratorConst.h>
14 #include <framework/gearbox/Unit.h>
15 #include <framework/gearbox/Const.h>
16 #include <framework/particledb/EvtGenDatabasePDG.h>
17 #include <framework/utilities/FileSystem.h>
21 #include <mdst/dataobjects/MCParticleGraph.h>
23 #include <framework/logging/Logger.h>
24 #include <framework/utilities/IOIntercept.h>
30 using namespace Pythia8;
43 setDescription(
"Fragmention of (u/d/s/c) quarks using PYTHIA8");
46 addParam(
"ParameterFile", m_parameterfile,
"Input parameter file for PYTHIA",
47 FileSystem::findFile(
"generators/modules/fragmentation/data/pythia_belle2.dat"));
48 addParam(
"ListPYTHIAEvent", m_listEvent,
"List event record of PYTHIA after hadronization", 0);
49 addParam(
"UseEvtGen", m_useEvtGen,
"Use EvtGen for specific decays", 1);
50 addParam(
"DecFile", m_DecFile,
"EvtGen decay file (DECAY.DEC)",
51 FileSystem::findFile(
"decfiles/dec/DECAY_BELLE2.DEC",
true));
52 addParam(
"UserDecFile", m_UserDecFile,
"User EvtGen decay file", std::string(
""));
53 addParam(
"CoherentMixing", m_coherentMixing,
"Decay the B0-B0bar coherently (should always be true)",
true);
65 FragmentationModule::~FragmentationModule()
70 void FragmentationModule::terminate()
75 statLogCapture.
start();
81 if (nAll) ratio = 100.0 * nGood / nAll;
83 B2WARNING(
"Not all events could be fragmented: " << nAll - nGood <<
" events failed.");
84 B2WARNING(
"Total number of events: " << nAll <<
", of these fragmented: " << nGood <<
", success-ratio (should be >97%): " << ratio
86 B2WARNING(
"Please contact the generator librarian if the success ratio is below 97%.");
87 B2WARNING(
"Please treat the success-ratio as correction of the effective cross section due to unphysical events.");
94 void FragmentationModule::initialize()
96 m_mcparticles.isRequired(m_particleList);
98 B2DEBUG(150,
"Initialize PYTHIA8");
106 initLogCapture.
start();
108 m_Pythia =
new Pythia;
109 m_PythiaEvent = &m_Pythia->event;
110 (*m_PythiaEvent) = 0;
113 EvtGen* evtGen = EvtGenInterface::createEvtGen(m_DecFile, m_coherentMixing);
114 loadEvtGenParticleData(m_Pythia);
117 m_Pythia->readString(
"ProcessLevel:all = off");
120 if (!m_Pythia->readFile(m_parameterfile))
121 B2FATAL(
"Cannot read Pythia parameter file.");
125 m_Pythia->setRndmEnginePtr(fragRndm);
134 B2INFO(
"Using PYTHIA EvtGen Interface");
135 const std::string defaultDecFile = FileSystem::findFile(
"decfiles/dec/DECAY_BELLE2.DEC",
true);
136 if (m_DecFile.empty()) {
137 B2ERROR(
"No global decay file defined, please make sure the parameter 'DecFile' is set correctly");
140 if (defaultDecFile.empty()) {
141 B2WARNING(
"Cannot find default decay file");
142 }
else if (defaultDecFile != m_DecFile) {
143 B2INFO(
"Using non-standard DECAY file \"" << m_DecFile <<
"\"");
146 evtgen->readDecayFile(m_UserDecFile);
157 m_Pythia->settings.listChanged();
165 void FragmentationModule::event()
168 mcParticleGraph.clear();
169 mcParticleGraph.loadList(m_particleList);
173 resetLogCapture.
start();
174 m_PythiaEvent->reset();
183 std::map<int, int> indexPYTHIA;
186 std::map<int, int> indexMCGraph;
189 int quarkPosition = 0;
192 int nPart = m_mcparticles.getEntries();
193 for (
int iPart = 0; iPart < nPart; iPart++) {
194 MCParticle* currParticle = m_mcparticles[iPart];
198 int pythiaIndex = addParticleToPYTHIA(*currParticle);
200 if (pythiaIndex != 0) {
201 indexPYTHIA[nAdded] = iPart;
202 if (pythiaIndex > 0) quarkPosition = iPart;
208 B2FATAL(
"Invalid number of quarks: " << nQuarks <<
" (should be 2)!");
212 B2WARNING(
"No virtual exchange particle given, no PYTHIA FSR in Decays");
215 m_Pythia->forceTimeShower(2, 3, 20.00);
226 eventLogCapture.
start();
227 int success = m_Pythia->next();
232 listLogCapture.
start();
233 m_PythiaEvent->list();
244 decayLogCapture.
start();
251 for (
int iPythiaPart = 0; iPythiaPart < m_Pythia->event.size(); ++iPythiaPart) {
252 auto oldindex = indexPYTHIA.find(iPythiaPart);
255 if (m_Pythia->event[iPythiaPart].id() == 90)
continue;
257 if (oldindex == end(indexPYTHIA)) {
261 int position = mcParticleGraph.size();
262 mcParticleGraph.addParticle();
263 indexMCGraph[iPythiaPart] = position;
269 if (m_Pythia->event[iPythiaPart].statusHepMC() < 1)
continue;
272 p->setPDG(m_Pythia->event[iPythiaPart].id());
275 TLorentzVector p4(m_Pythia->event[iPythiaPart].px(), m_Pythia->event[iPythiaPart].py(), m_Pythia->event[iPythiaPart].pz(),
276 m_Pythia->event[iPythiaPart].e());
278 p->setMass(m_Pythia->event[iPythiaPart].m());
281 p->setProductionVertex(m_Pythia->event[iPythiaPart].xProd() * Unit::mm, m_Pythia->event[iPythiaPart].yProd() * Unit::mm,
282 m_Pythia->event[iPythiaPart].zProd() * Unit::mm);
283 p->setProductionTime(m_Pythia->event[iPythiaPart].tProd() * Unit::mm / Const::speedOfLight);
284 p->setValidVertex(
true);
287 p->addStatus(MCParticleGraph::GraphParticle::c_PrimaryParticle);
290 if (m_Pythia->event[iPythiaPart].status() == 51 && m_Pythia->event[iPythiaPart].id() == 22) {
291 p->addStatus(MCParticleGraph::GraphParticle::c_IsFSRPhoton);
295 if (m_Pythia->event[iPythiaPart].status() == GeneratorConst::FSR_STATUS_CODE && m_Pythia->event[iPythiaPart].id() == 22) {
296 p->addStatus(MCParticleGraph::GraphParticle::c_IsPHOTOSPhoton);
300 if (m_Pythia->event[iPythiaPart].statusHepMC() == 1) {
301 p->addStatus(MCParticleGraph::GraphParticle::c_StableInGenerator);
305 const int motherid = m_Pythia->event[iPythiaPart].mother1();
308 auto motherindex = indexMCGraph.find(motherid);
310 if (motherindex != end(indexMCGraph)) {
311 int motheridingraph = indexMCGraph[motherid];
328 if (m_listEvent) m_PythiaEvent->list();
331 mcParticleGraph.generateList(m_particleList,
332 MCParticleGraph::c_clearParticles | MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
339 int FragmentationModule::addParticleToPYTHIA(
const MCParticle& mcParticle)
342 const int id = mcParticle.
getPDG();
346 bool isQuark =
false;
347 if (abs(
id) >= 1 && abs(
id) <= 5) isQuark =
true;
348 if (
id == 23) isVPho =
true;
350 if (!(isVPho || isQuark))
return 0;
354 B2WARNING(
"Quark already has a daughter!");
359 const double mass = mcParticle.
getMass();
361 const double energy = sqrt(mass * mass + p.Mag2());
365 m_PythiaEvent->append(23, -22, 0, 0, 2, 3, 0, 0, p.X(), p.Y(), p.Z(), energy, mass);
368 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 101, 0, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
371 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 0, 101, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
380 void FragmentationModule::loadEvtGenParticleData(Pythia8::Pythia* pythia)
382 Pythia8::ParticleData* particleData = &pythia->particleData;
383 int nParticles = EvtPDL::entries();
384 for (
int i = 0; i < nParticles; ++i) {
385 EvtId evtgenParticle = EvtPDL::getEntry(i);
391 int pdg = EvtPDL::getStdHep(evtgenParticle);
392 if (pdg <= 10 || (pdg > 20 && pdg <= 100))
394 EvtId evtgenAntiParticle = EvtPDL::chargeConj(evtgenParticle);
395 if (particleData->isParticle(pdg)) {
396 particleData->setAll(pdg,
397 EvtPDL::name(evtgenParticle),
398 EvtPDL::name(evtgenAntiParticle),
399 EvtPDL::getSpinType(evtgenParticle),
400 EvtPDL::chg3(evtgenParticle),
403 EvtPDL::getMeanMass(evtgenParticle),
404 EvtPDL::getWidth(evtgenParticle),
405 EvtPDL::getMinMass(evtgenParticle),
406 EvtPDL::getMaxMass(evtgenParticle),
407 EvtPDL::getctau(evtgenParticle));
415 FragmentationRndm::FragmentationRndm() : Pythia8::RndmEngine()
422 double value = gRandom->Rndm();
The Fragmentation module.
Minimal class for external random generator to be used in PYTHIA.
double flat()
flat random generator.
bool start()
Start intercepting the output.
Capture stdout and stderr and convert into log messages.
bool finish()
Finish the capture and emit the message if output has appeard on stdout or stderr.
Class to represent Particle data in graph.
A Class to store the Monte Carlo particle information.
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
float getMass() const
Return the particle mass in GeV.
TVector3 getMomentum() const
Return momentum.
int getPDG() const
Return PDG code of particle.
A class to perform decays via the external EvtGen decay program, see http://evtgen....
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.