Belle II Software development
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
28using namespace std;
29using namespace Belle2;
30using namespace Pythia8;
31
32//-----------------------------------------------------------------
33// Register the Module
34//-----------------------------------------------------------------
35REG_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
66{
67
68}
69
71{
72
73 // print internal pythia error statistics
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//-----------------------------------------------------------------
95{
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
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.
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 // Set signal particle suffix.
147 // Signal particle names must end in "sig" in user decay file.
148 evtgen->signalSuffix = "sig";
149 // Read the user decay file.
151 // Workaround for Pythia bug. It is the only way to call
152 // EvtGenDecays::updateData(true) to disable particle decays
153 // for all particles from EvtGen decay table. Thus, EvtGen generation
154 // has to be called before generation of the first Pythia event.
155 // Since Pythia event is currently empty, it actually only updates
156 // the particle properties (exactly what is necessary).
157 evtgen->decay();
158 }
159
160 // List variable(s) that differ from their defaults
161 m_Pythia->settings.listChanged();
162
163 initLogCapture.finish();
164}
165
166//-----------------------------------------------------------------
167// Event
168//-----------------------------------------------------------------
170{
171 // Reset the indices of the graph
174
175 // Reset PYTHIA event record to allow for new event
176 IOIntercept::OutputToLogMessages resetLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Error, 100, 100);
177 resetLogCapture.start();
178 m_PythiaEvent->reset();
179 resetLogCapture.finish();
180
181 // Reset counter for added quarks and vphos
182 nAdded = 0;
183 nQuarks = 0;
184 nVpho = 0;
185
186 // Store which MCParticle index belongs to which Pythia index
187 std::map<int, int> indexPYTHIA;
188
189 // Store which Pythia index belongs to which MCParticle index
190 std::map<int, int> indexMCGraph;
191
192 // Store position of the quark (mother of hadronized final state)
193 int quarkPosition = 0;
194
195 // Loop over all particles to find the quark pair
196 int nPart = m_mcparticles.getEntries();
197 for (int iPart = 0; iPart < nPart; iPart++) {
198 MCParticle* currParticle = m_mcparticles[iPart];
199
200 //returns quark id if it finds a quark, zero otherwise
201 //increments a counter for the number of found quarks
202 int pythiaIndex = addParticleToPYTHIA(*currParticle);
203
204 if (pythiaIndex != 0) {
205 indexPYTHIA[nAdded] = iPart;
206 if (pythiaIndex > 0) quarkPosition = iPart;
207 }
208 }
209
210 // Check needed if virtual exchange boson and two quarks are present
211 if (nQuarks != 2) {
212 B2FATAL("Invalid number of quarks: " << nQuarks << " (should be 2)!");
213 }
214
215 if (nVpho != 1) {
216 B2WARNING("No virtual exchange particle given, no PYTHIA FSR in Decays");
217 } else {
218 // Adding QCD and QED FSR
219 m_Pythia->forceTimeShower(2, 3, 20.00);
220 }
221
222 // Check needed if event is energetically possible
223 // ...
224
225 // Do the fragmentation using PYTHIA
226 setReturnValue(1); //return value is 1...
227 nAll = nAll + 1;
228
230 eventLogCapture.start();
231 int success = m_Pythia->next();
232 eventLogCapture.finish();
233
234 if (!success) {
236 listLogCapture.start();
237 m_PythiaEvent->list();
238 listLogCapture.finish();
239
240 setReturnValue(-1); //return value becomes -1 if trials were not successfull
241 } else {
242 nGood = nGood + 1;
243 }
244
245 // use evtgen to perform the decay
246 if (m_useEvtGen) {
248 decayLogCapture.start();
249 evtgen->decay();
250 decayLogCapture.finish();
251 }
252
253 // Loop over the PYTHIA list and assign the mother-daughter relation
254 // Might not work if the mother appear below the daughter in the event record
255 for (int iPythiaPart = 0; iPythiaPart < m_Pythia->event.size(); ++iPythiaPart) {
256 auto oldindex = indexPYTHIA.find(iPythiaPart);
257
258 //skip "system" particle generated by PYTHIA
259 if (m_Pythia->event[iPythiaPart].id() == 90) continue;
260
261 if (oldindex == end(indexPYTHIA)) {
262 // --> new particle
263
264 // Add to particle grapg
265 int position = mcParticleGraph.size();
267 indexMCGraph[iPythiaPart] = position;
268
270
271 // from PYTHIA manual: "<1: an empty entry, with no meaningful information and
272 // therefore to be skipped unconditionally (should not occur in PYTHIA)"
273 if (m_Pythia->event[iPythiaPart].statusHepMC() < 1) continue;
274
275 // Set PDG code
276 p->setPDG(m_Pythia->event[iPythiaPart].id());
277
278 // Set four vector
279 ROOT::Math::PxPyPzEVector p4(m_Pythia->event[iPythiaPart].px(), m_Pythia->event[iPythiaPart].py(),
280 m_Pythia->event[iPythiaPart].pz(),
281 m_Pythia->event[iPythiaPart].e());
282 p->set4Vector(p4);
283 p->setMass(m_Pythia->event[iPythiaPart].m());
284
285 // Set vertex
286 p->setProductionVertex(m_Pythia->event[iPythiaPart].xProd() * Unit::mm, m_Pythia->event[iPythiaPart].yProd() * Unit::mm,
287 m_Pythia->event[iPythiaPart].zProd() * Unit::mm);
288 p->setProductionTime(m_Pythia->event[iPythiaPart].tProd() * Unit::mm / Const::speedOfLight);
289 p->setValidVertex(true);
290
291 // Set all(!) particles from the generator to primary
293
294 // Set FSR flag from PYTHIA TimeShower:QEDshowerByQ
295 if (m_Pythia->event[iPythiaPart].status() == 51 && m_Pythia->event[iPythiaPart].id() == 22) {
297 }
298
299 // Set PHOTOS flag from PYTHIA-EvtGen
300 if (m_Pythia->event[iPythiaPart].status() == GeneratorConst::FSR_STATUS_CODE && m_Pythia->event[iPythiaPart].id() == 22) {
302 }
303
304 // Set stable at generator level
305 if (m_Pythia->event[iPythiaPart].statusHepMC() == 1) {
307 }
308
309 //set mother
310 const int motherid = m_Pythia->event[iPythiaPart].mother1();
311
312 //check if mother exists in indexMCGraph
313 auto motherindex = indexMCGraph.find(motherid);
314
315 if (motherindex != end(indexMCGraph)) {
316 int motheridingraph = indexMCGraph[motherid];
317 MCParticleGraph::GraphParticle* q = &mcParticleGraph[motheridingraph];
318 p->comesFrom(*q);
319 } else {
320 // Particle has no mother from PYTHIA, add quark as mother
322 p->comesFrom(*q);
323 }
324
325 } else {
326 // particle is already in the graph
327 // modify here if needed
328 // ...
329 }
330 }
331
332 // Print original PYTHIA list
333 if (m_listEvent) m_PythiaEvent->list();
334
335 // Create new ParticleGraph
338
339}
340
341//-----------------------------------------------------------------
342// addParticleToPYTHIA
343//-----------------------------------------------------------------
345{
346 //get PDG code
347 const int id = mcParticle.getPDG();
348
349 //check that this particle is a quark or a virtual photon(-like)
350 bool isVPho = false;
351 bool isQuark = false;
352 if (abs(id) >= 1 && abs(id) <= 5) isQuark = true;
353 if (id == 23) isVPho = true;
354
355 if (!(isVPho || isQuark)) return 0;
356
357 //check that there is no daughter for the quarks
358 if (isQuark && mcParticle.getDaughters().size()) {
359 B2WARNING("Quark already has a daughter!");
360 return 0;
361 }
362
363 //get some basic kinematics
364 const double mass = mcParticle.getMass();
365 const ROOT::Math::XYZVector& p = mcParticle.getMomentum();
366 const double energy = sqrt(mass * mass + p.Mag2());
367
368 //add this (anti)quark to the m_PythiaEvent
369 if (id == 23) {
370 m_PythiaEvent->append(23, -22, 0, 0, 2, 3, 0, 0, p.X(), p.Y(), p.Z(), energy, mass);
371 nVpho++;
372 } else if (id > 0) {
373 m_PythiaEvent->append(id, 23, 1, 0, 0, 0, 101, 0, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
374 nQuarks++;
375 } else if (id < 0) {
376 m_PythiaEvent->append(id, 23, 1, 0, 0, 0, 0, 101, p.X(), p.Y(), p.Z(), energy, mass, 20.0);
377 nQuarks++;
378 }
379
380 nAdded++;
381
382 return id;
383}
384
386{
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);
391 /*
392 * Pythia uses absolute value of the PDG code, thus, only positive codes
393 * are necessary. Only leptons (11 <= pdg <= 20) and hadrons (pdg > 100)
394 * are updated.
395 */
396 int pdg = EvtPDL::getStdHep(evtgenParticle);
397 if (pdg <= 10 || (pdg > 20 && pdg <= 100))
398 continue;
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),
406 // colType == 0 for uncolored particles.
407 0,
408 EvtPDL::getMeanMass(evtgenParticle),
409 EvtPDL::getWidth(evtgenParticle),
410 EvtPDL::getMinMass(evtgenParticle),
411 EvtPDL::getMaxMass(evtgenParticle),
412 EvtPDL::getctau(evtgenParticle));
413 }
414 }
415}
416
417//-----------------------------------------------------------------
418// random generator for PYTHIA
419//-----------------------------------------------------------------
420FragmentationRndm::FragmentationRndm() : Pythia8::RndmEngine()
421{
422
423}
424
426{
427 double value = gRandom->Rndm();
428 return value;
429}
static const double speedOfLight
[cm/ns]
Definition: Const.h:695
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...
Definition: FileSystem.cc:151
EvtGenDecays * evtgen
EvtGen decay engine inside PYTHIA8.
int nGood
number of events with successful fragmentation.
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.
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:236
@ c_Error
Error: for things that went wrong and have to be fixed.
Definition: LogConfig.h:30
@ c_Debug
Debug: for code development.
Definition: LogConfig.h:26
@ c_Warning
Warning: for potential problems that the user should pay attention to.
Definition: LogConfig.h:29
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.
Definition: MCParticle.h:32
@ c_IsFSRPhoton
bit 7: Particle is from finial state radiation
Definition: MCParticle.h:61
@ c_IsPHOTOSPhoton
bit 8: Particle is an radiative photon from PHOTOS
Definition: MCParticle.h:63
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
Definition: MCParticle.h:49
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:52
float getMass() const
Return the particle mass in GeV.
Definition: MCParticle.h:135
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
ROOT::Math::XYZVector getMomentum() const
Return momentum.
Definition: MCParticle.h:198
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setReturnValue(int value)
Sets the return value for this module as integer.
Definition: Module.cc:220
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.
Definition: StoreArray.h:216
static const double mm
[millimeters]
Definition: Unit.h:70
A class to perform decays via the external EvtGen decay program, see http://evtgen....
Definition: EvtGenDecays.h:88
void readDecayFile(std::string decayFile, bool xml=false)
Read an EvtGen user decay file.
Definition: EvtGenDecays.h:131
std::string signalSuffix
The suffix indicating an EvtGen particle or alias is signal.
Definition: EvtGenDecays.h:154
double decay()
Perform all decays and return the event weight.
Definition: EvtGenDecays.h:327
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
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.
STL namespace.