Belle II Software  release-06-00-14
FragmentationModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <generators/modules/fragmentation/FragmentationModule.h>
10 
11 #include <generators/evtgen/EvtGenInterface.h>
12 #include <generators/utilities/GeneratorConst.h>
13 
14 #include <framework/gearbox/Unit.h>
15 #include <framework/gearbox/Const.h>
16 #include <framework/particledb/EvtGenDatabasePDG.h>
17 #include <framework/utilities/FileSystem.h>
18 
19 #include <TRandom3.h>
20 
21 #include <mdst/dataobjects/MCParticleGraph.h>
22 
23 #include <framework/logging/Logger.h>
24 #include <framework/utilities/IOIntercept.h>
25 
26 #include <string>
27 
28 using namespace std;
29 using namespace Belle2;
30 using namespace Pythia8;
31 
32 //-----------------------------------------------------------------
33 // Register the Module
34 //-----------------------------------------------------------------
35 REG_MODULE(Fragmentation)
36 
37 //-----------------------------------------------------------------
38 // Implementation
39 //-----------------------------------------------------------------
41 {
42  //Set module properties
43  setDescription("Fragmention of (u/d/s/c) quarks using PYTHIA8");
44 
45  //Parameter definition
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);
54 
55  //initialize member variables
56  evtgen = 0;
57  nAdded = 0;
58  nQuarks = 0;
59  nVpho = 0;
60  nAll = 0;
61  nGood = 0;
62 }
63 
64 
65 FragmentationModule::~FragmentationModule()
66 {
67 
68 }
69 
70 void FragmentationModule::terminate()
71 {
72 
73  // print internal pythia error statistics
74  IOIntercept::OutputToLogMessages statLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Error, 50, 100);
75  statLogCapture.start();
76  m_Pythia->stat();
77  statLogCapture.finish();
78 
79  if (nAll != nGood) {
80  double ratio = 0.; //ratio of good over all events
81  if (nAll) ratio = 100.0 * nGood / nAll;
82 
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
85  << "%");
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.");
88  }
89 }
90 
91 //-----------------------------------------------------------------
92 // Initialize
93 //-----------------------------------------------------------------
94 void FragmentationModule::initialize()
95 {
96  m_mcparticles.isRequired(m_particleList);
97 
98  B2DEBUG(150, "Initialize PYTHIA8");
99 
100  // Generator and the shorthand m_PythiaEvent = pythia->event are declared in .h file
101  // A simple way to collect all the changes is to store the parameter values in a separate file,
102  // with one line per change. This should be done between the creation of the Pythia object
103  // and the init call for it.
104 
105  IOIntercept::OutputToLogMessages initLogCapture("PYTHIA", LogConfig::c_Debug, LogConfig::c_Warning, 100, 100);
106  initLogCapture.start();
107 
108  m_Pythia = new Pythia;
109  m_PythiaEvent = &m_Pythia->event;
110  (*m_PythiaEvent) = 0;
111 
112  // Load EvtGen particle data.
113  EvtGen* evtGen = EvtGenInterface::createEvtGen(m_DecFile, m_coherentMixing);
114  loadEvtGenParticleData(m_Pythia);
115 
116  // Switch off ProcessLevel
117  m_Pythia->readString("ProcessLevel:all = off");
118 
119  // Read the PYTHIA input file, overrides parameters
120  if (!m_Pythia->readFile(m_parameterfile))
121  B2FATAL("Cannot read Pythia parameter file.");
122 
123  // Set framework generator
124  FragmentationRndm* fragRndm = new FragmentationRndm();
125  m_Pythia->setRndmEnginePtr(fragRndm);
126 
127  // Initialize PYTHIA
128  m_Pythia->init();
129 
130  // Set EvtGen (after m_Pythia->init())
131  evtgen = 0;
132 
133  if (m_useEvtGen) {
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");
138  return;
139  }
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 << "\"");
144  }
145  evtgen = new EvtGenDecays(m_Pythia, evtGen);
146  evtgen->readDecayFile(m_UserDecFile);
147  // Workaround for Pythia bug. It is the only way to call
148  // EvtGenDecays::updateData(true) to disable particle decays
149  // for all particles from EvtGen decay table. Thus, EvtGen generation
150  // has to be called before generation of the first Pythia event.
151  // Since Pythia event is currently empty, it actually only updates
152  // the particle properties (exactly what is necessary).
153  evtgen->decay();
154  }
155 
156  // List variable(s) that differ from their defaults
157  m_Pythia->settings.listChanged();
158 
159  initLogCapture.finish();
160 }
161 
162 //-----------------------------------------------------------------
163 // Event
164 //-----------------------------------------------------------------
165 void FragmentationModule::event()
166 {
167  // Reset the indices of the graph
168  mcParticleGraph.clear();
169  mcParticleGraph.loadList(m_particleList);
170 
171  // Reset PYTHIA event record to allow for new event
172  IOIntercept::OutputToLogMessages resetLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Error, 100, 100);
173  resetLogCapture.start();
174  m_PythiaEvent->reset();
175  resetLogCapture.finish();
176 
177  // Reset counter for added quarks and vphos
178  nAdded = 0;
179  nQuarks = 0;
180  nVpho = 0;
181 
182  // Store which MCParticle index belongs to which Pythia index
183  std::map<int, int> indexPYTHIA;
184 
185  // Store which Pythia index belongs to which MCParticle index
186  std::map<int, int> indexMCGraph;
187 
188  // Store position of the quark (mother of hadronized final state)
189  int quarkPosition = 0;
190 
191  // Loop over all particles to find the quark pair
192  int nPart = m_mcparticles.getEntries();
193  for (int iPart = 0; iPart < nPart; iPart++) {
194  MCParticle* currParticle = m_mcparticles[iPart];
195 
196  //returns quark id if it finds a quark, zero otherwise
197  //increments a counter for the number of found quarks
198  int pythiaIndex = addParticleToPYTHIA(*currParticle);
199 
200  if (pythiaIndex != 0) {
201  indexPYTHIA[nAdded] = iPart;
202  if (pythiaIndex > 0) quarkPosition = iPart;
203  }
204  }
205 
206  // Check needed if virtual exchange boson and two quarks are present
207  if (nQuarks != 2) {
208  B2FATAL("Invalid number of quarks: " << nQuarks << " (should be 2)!");
209  }
210 
211  if (nVpho != 1) {
212  B2WARNING("No virtual exchange particle given, no PYTHIA FSR in Decays");
213  } else {
214  // Adding QCD and QED FSR
215  m_Pythia->forceTimeShower(2, 3, 20.00);
216  }
217 
218  // Check needed if event is energetically possible
219  // ...
220 
221  // Do the fragmentation using PYTHIA
222  setReturnValue(1); //return value is 1...
223  nAll = nAll + 1;
224 
225  IOIntercept::OutputToLogMessages eventLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Error, 50, 100);
226  eventLogCapture.start();
227  int success = m_Pythia->next();
228  eventLogCapture.finish();
229 
230  if (!success) {
231  IOIntercept::OutputToLogMessages listLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Error, 50, 100);
232  listLogCapture.start();
233  m_PythiaEvent->list();
234  listLogCapture.finish();
235 
236  setReturnValue(-1); //return value becomes -1 if trials were not successfull
237  } else {
238  nGood = nGood + 1;
239  }
240 
241  // use evtgen to perform the decay
242  if (m_useEvtGen) {
243  IOIntercept::OutputToLogMessages decayLogCapture("PYTHIA", LogConfig::c_Debug, LogConfig::c_Warning, 100, 100);
244  decayLogCapture.start();
245  evtgen->decay();
246  decayLogCapture.finish();
247  }
248 
249  // Loop over the PYTHIA list and assign the mother-daughter relation
250  // Might not work if the mother appear below the daughter in the event record
251  for (int iPythiaPart = 0; iPythiaPart < m_Pythia->event.size(); ++iPythiaPart) {
252  auto oldindex = indexPYTHIA.find(iPythiaPart);
253 
254  //skip "system" particle generated by PYTHIA
255  if (m_Pythia->event[iPythiaPart].id() == 90) continue;
256 
257  if (oldindex == end(indexPYTHIA)) {
258  // --> new particle
259 
260  // Add to particle grapg
261  int position = mcParticleGraph.size();
262  mcParticleGraph.addParticle();
263  indexMCGraph[iPythiaPart] = position;
264 
265  MCParticleGraph::GraphParticle* p = &mcParticleGraph[position];
266 
267  // from PYTHIA manual: "<1: an empty entry, with no meaningful information and
268  // therefore to be skipped unconditionally (should not occur in PYTHIA)"
269  if (m_Pythia->event[iPythiaPart].statusHepMC() < 1) continue;
270 
271  // Set PDG code
272  p->setPDG(m_Pythia->event[iPythiaPart].id());
273 
274  // Set four vector
275  TLorentzVector p4(m_Pythia->event[iPythiaPart].px(), m_Pythia->event[iPythiaPart].py(), m_Pythia->event[iPythiaPart].pz(),
276  m_Pythia->event[iPythiaPart].e());
277  p->set4Vector(p4);
278  p->setMass(m_Pythia->event[iPythiaPart].m());
279 
280  // Set vertex
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);
285 
286  // Set all(!) particles from the generator to primary
287  p->addStatus(MCParticleGraph::GraphParticle::c_PrimaryParticle);
288 
289  // Set FSR flag from PYTHIA TimeShower:QEDshowerByQ
290  if (m_Pythia->event[iPythiaPart].status() == 51 && m_Pythia->event[iPythiaPart].id() == 22) {
291  p->addStatus(MCParticleGraph::GraphParticle::c_IsFSRPhoton);
292  }
293 
294  // Set PHOTOS flag from PYTHIA-EvtGen
295  if (m_Pythia->event[iPythiaPart].status() == GeneratorConst::FSR_STATUS_CODE && m_Pythia->event[iPythiaPart].id() == 22) {
296  p->addStatus(MCParticleGraph::GraphParticle::c_IsPHOTOSPhoton);
297  }
298 
299  // Set stable at generator level
300  if (m_Pythia->event[iPythiaPart].statusHepMC() == 1) {
301  p->addStatus(MCParticleGraph::GraphParticle::c_StableInGenerator);
302  }
303 
304  //set mother
305  const int motherid = m_Pythia->event[iPythiaPart].mother1();
306 
307  //check if mother exists in indexMCGraph
308  auto motherindex = indexMCGraph.find(motherid);
309 
310  if (motherindex != end(indexMCGraph)) {
311  int motheridingraph = indexMCGraph[motherid];
312  MCParticleGraph::GraphParticle* q = &mcParticleGraph[motheridingraph];
313  p->comesFrom(*q);
314  } else {
315  // Particle has no mother from PYTHIA, add quark as mother
316  MCParticleGraph::GraphParticle* q = &mcParticleGraph[quarkPosition];
317  p->comesFrom(*q);
318  }
319 
320  } else {
321  // particle is already in the graph
322  // modify here if needed
323  // ...
324  }
325  }
326 
327  // Print original PYTHIA list
328  if (m_listEvent) m_PythiaEvent->list();
329 
330  // Create new ParticleGraph
331  mcParticleGraph.generateList(m_particleList,
332  MCParticleGraph::c_clearParticles | MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
333 
334 }
335 
336 //-----------------------------------------------------------------
337 // addParticleToPYTHIA
338 //-----------------------------------------------------------------
339 int FragmentationModule::addParticleToPYTHIA(const MCParticle& mcParticle)
340 {
341  //get PDG code
342  const int id = mcParticle.getPDG();
343 
344  //check that this particle is a quark or a virtual photon(-like)
345  bool isVPho = false;
346  bool isQuark = false;
347  if (abs(id) >= 1 && abs(id) <= 5) isQuark = true;
348  if (id == 23) isVPho = true;
349 
350  if (!(isVPho || isQuark)) return 0;
351 
352  //check that there is no daughter for the quarks
353  if (isQuark && mcParticle.getDaughters().size()) {
354  B2WARNING("Quark already has a daughter!");
355  return 0;
356  }
357 
358  //get some basic kinematics
359  const double mass = mcParticle.getMass();
360  const TVector3& p = mcParticle.getMomentum();
361  const double energy = sqrt(mass * mass + p.Mag2());
362 
363  //add this (anti)quark to the m_PythiaEvent
364  if (id == 23) {
365  m_PythiaEvent->append(23, -22, 0, 0, 2, 3, 0, 0, p.X(), p.Y(), p.Z(), energy, mass);
366  nVpho++;
367  } else if (id > 0) {
368  m_PythiaEvent->append(id, 23, 1, 0, 0, 0, 101, 0, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
369  nQuarks++;
370  } else if (id < 0) {
371  m_PythiaEvent->append(id, 23, 1, 0, 0, 0, 0, 101, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
372  nQuarks++;
373  }
374 
375  nAdded++;
376 
377  return id;
378 }
379 
380 void FragmentationModule::loadEvtGenParticleData(Pythia8::Pythia* pythia)
381 {
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);
386  /*
387  * Pythia uses absolute value of the PDG code, thus, only positive codes
388  * are necessary. Only leptons (11 <= pdg <= 20) and hadrons (pdg > 100)
389  * are updated.
390  */
391  int pdg = EvtPDL::getStdHep(evtgenParticle);
392  if (pdg <= 10 || (pdg > 20 && pdg <= 100))
393  continue;
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),
401  // colType == 0 for uncolored particles.
402  0,
403  EvtPDL::getMeanMass(evtgenParticle),
404  EvtPDL::getWidth(evtgenParticle),
405  EvtPDL::getMinMass(evtgenParticle),
406  EvtPDL::getMaxMass(evtgenParticle),
407  EvtPDL::getctau(evtgenParticle));
408  }
409  }
410 }
411 
412 //-----------------------------------------------------------------
413 // random generator for PYTHIA
414 //-----------------------------------------------------------------
415 FragmentationRndm::FragmentationRndm() : Pythia8::RndmEngine()
416 {
417 
418 }
419 
421 {
422  double value = gRandom->Rndm();
423  return value;
424 }
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.
Definition: IOIntercept.h:155
Capture stdout and stderr and convert into log messages.
Definition: IOIntercept.h:226
bool finish()
Finish the capture and emit the message if output has appeard on stdout or stderr.
Definition: IOIntercept.cc:235
Class to represent Particle data in graph.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:50
float getMass() const
Return the particle mass in GeV.
Definition: MCParticle.h:135
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:198
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
A class to perform decays via the external EvtGen decay program, see http://evtgen....
Definition: EvtGenDecays.h:88
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.