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>
30using namespace Pythia8;
48 addParam(
"ListPYTHIAEvent",
m_listEvent,
"List event record of PYTHIA after hadronization", 0);
75 statLogCapture.
start();
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.");
98 B2DEBUG(150,
"Initialize PYTHIA8");
106 initLogCapture.
start();
110 (*m_PythiaEvent) = 0;
117 m_Pythia->readString(
"ProcessLevel:all = off");
121 B2FATAL(
"Cannot read Pythia parameter file.");
125 m_Pythia->setRndmEnginePtr(fragRndm);
134 B2INFO(
"Using PYTHIA EvtGen Interface");
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 <<
"\"");
177 resetLogCapture.
start();
187 std::map<int, int> indexPYTHIA;
190 std::map<int, int> indexMCGraph;
193 int quarkPosition = 0;
197 for (
int iPart = 0; iPart < nPart; iPart++) {
204 if (pythiaIndex != 0) {
205 indexPYTHIA[
nAdded] = iPart;
206 if (pythiaIndex > 0) quarkPosition = iPart;
212 B2FATAL(
"Invalid number of quarks: " <<
nQuarks <<
" (should be 2)!");
216 B2WARNING(
"No virtual exchange particle given, no PYTHIA FSR in Decays");
219 m_Pythia->forceTimeShower(2, 3, 20.00);
230 eventLogCapture.
start();
236 listLogCapture.
start();
248 decayLogCapture.
start();
255 for (
int iPythiaPart = 0; iPythiaPart <
m_Pythia->event.size(); ++iPythiaPart) {
256 auto oldindex = indexPYTHIA.find(iPythiaPart);
259 if (
m_Pythia->event[iPythiaPart].id() == 90)
continue;
261 if (oldindex == end(indexPYTHIA)) {
267 indexMCGraph[iPythiaPart] = position;
273 if (
m_Pythia->event[iPythiaPart].statusHepMC() < 1)
continue;
276 p->setPDG(
m_Pythia->event[iPythiaPart].id());
279 ROOT::Math::PxPyPzEVector p4(
m_Pythia->event[iPythiaPart].px(),
m_Pythia->event[iPythiaPart].py(),
283 p->setMass(
m_Pythia->event[iPythiaPart].m());
289 p->setValidVertex(
true);
295 if (
m_Pythia->event[iPythiaPart].status() == 51 &&
m_Pythia->event[iPythiaPart].id() == 22) {
300 if (
m_Pythia->event[iPythiaPart].status() == GeneratorConst::FSR_STATUS_CODE &&
m_Pythia->event[iPythiaPart].id() == 22) {
305 if (
m_Pythia->event[iPythiaPart].statusHepMC() == 1) {
310 const int motherid =
m_Pythia->event[iPythiaPart].mother1();
313 auto motherindex = indexMCGraph.find(motherid);
315 if (motherindex != end(indexMCGraph)) {
316 int motheridingraph = indexMCGraph[motherid];
347 const int id = mcParticle.
getPDG();
351 bool isQuark =
false;
352 if (abs(
id) >= 1 && abs(
id) <= 5) isQuark =
true;
353 if (
id == 23) isVPho =
true;
355 if (!(isVPho || isQuark))
return 0;
359 B2WARNING(
"Quark already has a daughter!");
364 const double mass = mcParticle.
getMass();
365 const ROOT::Math::XYZVector& p = mcParticle.
getMomentum();
366 const double energy =
sqrt(mass * mass + p.Mag2());
370 m_PythiaEvent->append(23, -22, 0, 0, 2, 3, 0, 0, p.X(), p.Y(), p.Z(), energy, mass);
373 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 101, 0, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
376 m_PythiaEvent->append(
id, 23, 1, 0, 0, 0, 0, 101, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
387 Pythia8::ParticleData* particleData = &pythia->particleData;
388 int nParticles = EvtPDL::entries();
389 for (
int i = 0; i < nParticles; ++i) {
390 EvtId evtgenParticle = EvtPDL::getEntry(i);
396 int pdg = EvtPDL::getStdHep(evtgenParticle);
397 if (pdg <= 10 || (pdg > 20 && pdg <= 100))
399 EvtId evtgenAntiParticle = EvtPDL::chargeConj(evtgenParticle);
400 if (particleData->isParticle(pdg)) {
401 particleData->setAll(pdg,
402 EvtPDL::name(evtgenParticle),
403 EvtPDL::name(evtgenAntiParticle),
404 EvtPDL::getSpinType(evtgenParticle),
405 EvtPDL::chg3(evtgenParticle),
408 EvtPDL::getMeanMass(evtgenParticle),
409 EvtPDL::getWidth(evtgenParticle),
410 EvtPDL::getMinMass(evtgenParticle),
411 EvtPDL::getMaxMass(evtgenParticle),
412 EvtPDL::getctau(evtgenParticle));
420FragmentationRndm::FragmentationRndm() : Pythia8::RndmEngine()
427 double value = gRandom->Rndm();
static const double speedOfLight
[cm/ns]
static EvtGen * createEvtGen(const std::string &decayFileName, bool coherentMixing)
Create and initialize an EvtGen instance:
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
EvtGenDecays * evtgen
EvtGen decay engine inside PYTHIA8.
int nGood
number of events with successful fragmentation.
FragmentationModule()
Constructor.
void loadEvtGenParticleData(Pythia8::Pythia *pythia)
Load EvtGen particle data.
bool m_coherentMixing
decay the B0-B0bar coherently.
int nVpho
number of virtual exchange particles.
virtual void initialize() override
Initialize the module.
virtual void event() override
Event method (process events)
Pythia8::Pythia * m_Pythia
Pythia generator.
StoreArray< MCParticle > m_mcparticles
store array for the MCParticles
virtual ~FragmentationModule()
Destructor.
virtual void terminate() override
terminate the module
Pythia8::Event * m_PythiaEvent
Pythia event.
int addParticleToPYTHIA(const MCParticle &mcParticle)
Add particle to Pythia event.
std::string m_DecFile
EvtGen decay file.
int m_useEvtGen
use EvtGen for some decays.
int nAll
number of events created.
std::string m_UserDecFile
User EvtGen decay file.
std::string m_particleList
The name of the MCParticle collection.
int nQuarks
number of quarks.
MCParticleGraph mcParticleGraph
An instance of the MCParticle graph.
std::string m_parameterfile
Module parameters.
int nAdded
number of added particles.
int m_listEvent
list event generated by PYTHIA.
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.
@ c_Error
Error: for things that went wrong and have to be fixed.
@ c_Debug
Debug: for code development.
@ c_Warning
Warning: for potential problems that the user should pay attention to.
Class to represent Particle data in graph.
@ c_checkCyclic
Check for cyclic dependencies.
@ c_clearParticles
Clear the particle list before adding the graph.
@ c_setDecayInfo
Set decay time and vertex.
size_t size() const
Return the number of particles in the graph.
void loadList(const std::string &name="")
Load the MCParticle list given by name into the Graph.
void generateList(const std::string &name="", int options=c_setNothing)
Generates the MCParticle list and stores it in the StoreArray with the given name.
A Class to store the Monte Carlo particle information.
@ c_IsFSRPhoton
bit 7: Particle is from finial state radiation
@ c_IsPHOTOSPhoton
bit 8: Particle is an radiative photon from PHOTOS
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
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.
int getPDG() const
Return PDG code of particle.
ROOT::Math::XYZVector getMomentum() const
Return momentum.
void setDescription(const std::string &description)
Sets the description of the module.
void setReturnValue(int value)
Sets the return value for this module as integer.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
int getEntries() const
Get the number of objects in the array.
static const double mm
[millimeters]
A class to perform decays via the external EvtGen decay program, see http://evtgen....
void readDecayFile(std::string decayFile, bool xml=false)
Read an EvtGen user decay file.
std::string signalSuffix
The suffix indicating an EvtGen particle or alias is signal.
double decay()
Perform all decays and return the event weight.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
double sqrt(double a)
sqrt for double
void clear()
Reset particles and decay information to make the class reusable.
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.