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