Belle II Software  release-06-00-14
EvtGenInterface.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 <framework/gearbox/Unit.h>
10 #include <framework/gearbox/Const.h>
11 #include <framework/logging/Logger.h>
12 #include <framework/particledb/EvtGenDatabasePDG.h>
13 #include <framework/utilities/FileSystem.h>
14 #include <generators/evtgen/EvtGenInterface.h>
15 #include <mdst/dataobjects/MCParticleGraph.h>
16 #include <generators/evtgen/EvtGenModelRegister.h>
17 
18 #include <memory>
19 #include <string>
20 #include <queue>
21 
22 #include <EvtGenExternal/EvtExternalGenList.hh>
23 #include <EvtGenBase/EvtAbsRadCorr.hh>
24 #include <EvtGenBase/EvtCPUtil.hh>
25 #include <EvtGenBase/EvtDecayTable.hh>
26 #include <EvtGenBase/EvtDecayBase.hh>
27 #include <EvtGenBase/EvtParticleFactory.hh>
28 #include <EvtGenBase/EvtPDL.hh>
29 #include <EvtGenBase/EvtRandom.hh>
30 
31 #include <TLorentzVector.h>
32 
33 using namespace std;
34 using namespace Belle2;
35 
36 EvtGenFwRandEngine EvtGenInterface::m_eng;
37 
38 EvtGenInterface::~EvtGenInterface()
39 {
40  EvtDecayTable* evtDecayTable = EvtDecayTable::getInstance();
41  for (unsigned int i = 0; i < EvtPDL::entries(); ++i) {
42  for (int j = 0; j < evtDecayTable->getNMode(i); ++j) {
43  delete evtDecayTable->getDecay(i, j);
44  }
45  }
46  if (m_Generator) delete m_Generator;
47 }
48 
49 EvtGen* EvtGenInterface::createEvtGen(const std::string& DECFileName, const bool coherentMixing)
50 {
51  IOIntercept::OutputToLogMessages initLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Debug, 100, 100);
52  initLogCapture.start();
53  EvtRandom::setRandomEngine(&EvtGenInterface::m_eng);
54 
55  // Official BelleII models
56  std::list<EvtDecayBase*> extraModels = EvtGenModelRegister::getModels();
57 
58  // Fill model list with pythia, photos etc.
59  // The parameter 'false' means that Pythia codes must not be converted,
60  // since the conversion is applied to all the decfiles since release-06.
61  EvtExternalGenList genList{false};
62  EvtAbsRadCorr* radCorrEngine = genList.getPhotosModel();
63  list<EvtDecayBase*> modelList = genList.getListOfModels();
64  extraModels.insert(extraModels.end(), modelList.begin(), modelList.end());
65 
67  EvtGenDatabasePDG::Instance()->WriteEvtGenTable(tmp);
68 
69  auto mixingMode = EvtCPUtil::Incoherent;
70  if (coherentMixing)
71  mixingMode = EvtCPUtil::Coherent;
72  else {
73  B2WARNING("Evtgen has been set to decay the B mesons incoherently. This is useful as a workaround only to generate Y(5S, 6S) -> BBar for QCD studies.");
74  B2WARNING("If you are generating Y(4S) events, you _really_ must use the coherent decay mode.");
75  }
76 
77  return new EvtGen(DECFileName.c_str(), tmp.getName().c_str(), &EvtGenInterface::m_eng,
78  radCorrEngine, &extraModels, mixingMode);
79 
80  initLogCapture.finish();
81 }
82 
83 int EvtGenInterface::setup(const std::string& DECFileName, const std::string& parentParticle,
84  const std::string& userFileName, const bool coherentMixing)
85 {
86  B2DEBUG(150, "Begin initialisation of EvtGen Interface.");
87 
88  //tauola prints normal things to stderr.. oh well.
89  IOIntercept::OutputToLogMessages initLogCapture("EvtGen", LogConfig::c_Debug, LogConfig::c_Debug, 100, 100);
90  initLogCapture.start();
91  if (!m_Generator) {
92  m_Generator = createEvtGen(DECFileName, coherentMixing);
93  }
94  if (!userFileName.empty()) {
95  m_Generator->readUDecay(userFileName.c_str());
96  }
97 
98  // Setup Parent Particle in rest frame
99  if (parentParticle != "") {
100  m_ParentInitialized = true;
101  m_ParentParticle = EvtPDL::getId(parentParticle);
102  }
103  initLogCapture.finish();
104 
105  B2DEBUG(150, "End initialisation of EvtGen Interface.");
106 
107  return 0;
108 }
109 
110 
111 int EvtGenInterface::simulateEvent(MCParticleGraph& graph, TLorentzVector pParentParticle, TVector3 pPrimaryVertex,
112  int inclusiveType, const std::string& inclusiveParticle)
113 {
114  EvtId inclusiveParticleID, inclusiveAntiParticleID;
115 
116  if (!m_ParentInitialized)
117  B2FATAL("Parent particle is not initialized.");
118  //Init evtgen
119  m_pinit.set(pParentParticle.E(), pParentParticle.X(), pParentParticle.Y(), pParentParticle.Z());
120 
121  if (inclusiveType != 0) {
122  inclusiveParticleID = EvtPDL::getId(inclusiveParticle);
123  if (inclusiveParticleID.getId() < 0)
124  B2FATAL("Incorrect inclusive particle " << inclusiveParticle);
125  inclusiveAntiParticleID = EvtPDL::chargeConj(inclusiveParticleID);
126  if (inclusiveAntiParticleID.getId() < 0) {
127  B2FATAL("Cannot find the charge-conjugate particle for "
128  << inclusiveParticle);
129  }
130  }
131 
132  bool we_got_inclusive_particle = false;
133  do {
134  m_logCapture.start();
135  m_parent = EvtParticleFactory::particleFactory(m_ParentParticle, m_pinit);
136  m_parent->setVectorSpinDensity();
137  m_Generator->generateDecay(m_parent);
138  m_logCapture.finish();
139 
140  if (inclusiveType != 0) {
141  EvtParticle* p = m_parent;
142  // following loop will go through generated event and check it for
143  // presense of inclusive particle
144  do {
145  //for (int ii = 0; ii < iPart ; ii++) {
146  //std::cout << p->getPDGId() << std::endl;
147  if (p->getId() == inclusiveParticleID ||
148  (inclusiveType == 2 && p->getId() == inclusiveAntiParticleID)) {
149  we_got_inclusive_particle = true;
150  break;
151  }
152  p = p->nextIter(m_parent);
153  } while (p != 0);
154 
155  if (!we_got_inclusive_particle) {
156  m_parent->deleteTree();
157  }
158  } else {
159  // we don't do inclusive skimming so any generated event will do
160  we_got_inclusive_particle = true;
161  }
162  } while (!we_got_inclusive_particle);
163 
164  // B2INFO("after generate Decay.");
165 
166  int iPart = addParticles2Graph(m_parent, graph, pPrimaryVertex, NULL);
167  graph.generateList("", MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
168 
169  m_parent->deleteTree();
170 
171  return iPart; //returns the number of generated particles from evtgen
172 }
173 
174 int EvtGenInterface::simulateDecay(MCParticleGraph& graph,
176 {
177  int pdg;
178  EvtId id;
179  TLorentzVector momentum = parent.get4Vector();
180  TVector3 vertex = parent.getVertex();
181  m_pinit.set(momentum.E(), momentum.X(), momentum.Y(), momentum.Z());
182  m_logCapture.start();
183  // we want to decay the particle so the decay time in the tree needs to be lower
184  // than whatever the daughters will get
185  parent.setDecayTime(-std::numeric_limits<float>::infinity());
186  pdg = parent.getPDG();
187  id = EvtPDL::evtIdFromStdHep(pdg);
188  m_parent = EvtParticleFactory::particleFactory(id, m_pinit);
189  if (pdg == 10022) // Virtual photon
190  m_parent->setVectorSpinDensity();
191  else
192  m_parent->setDiagonalSpinDensity();
193  m_Generator->generateDecay(m_parent);
194  m_logCapture.finish();
195  int iPart = addParticles2Graph(m_parent, graph, vertex, &parent, parent.getProductionTime());
196  m_parent->deleteTree();
197  return iPart;
198 }
199 
200 int EvtGenInterface::addParticles2Graph(EvtParticle* top, MCParticleGraph& graph, TVector3 pPrimaryVertex,
201  MCParticleGraph::GraphParticle* parent, double timeOffset)
202 {
203  //Fill top particle in the tree & starting the queue:
204  const int existingParticles = graph.size();
206  if (parent == NULL)
207  p = &graph.addParticle();
208  else
209  p = parent;
210  updateGraphParticle(top, p, pPrimaryVertex);
211 
212  typedef pair<MCParticleGraph::GraphParticle*, EvtParticle*> halfFamily;
213  halfFamily currFamily;
214  queue < halfFamily > heritancesQueue;
215 
216  for (uint idaughter = 0; idaughter < top->getNDaug(); idaughter++) {
217  currFamily.first = p;
218  currFamily.second = top->getDaug(idaughter);
219  heritancesQueue.push(currFamily);
220  }
221 
222  //now we can go through the queue:
223  while (!heritancesQueue.empty()) {
224  currFamily = heritancesQueue.front(); //get the first entry from the queue
225  heritancesQueue.pop(); //remove the entry.
226 
227  MCParticleGraph::GraphParticle* currMother = currFamily.first;
228  EvtParticle* currDaughter = currFamily.second;
229 
230  //putting the daughter in the graph:
231  MCParticleGraph::GraphParticle* graphDaughter = &graph.addParticle();
232  updateGraphParticle(currDaughter, graphDaughter, pPrimaryVertex, timeOffset);
233 
234  //add relation between mother and daughter to graph:
235  currMother->decaysInto((*graphDaughter));
236 
237  int nGrandChildren = currDaughter->getNDaug();
238 
239  if (nGrandChildren == 0)
240  graphDaughter->addStatus(MCParticle::c_StableInGenerator);
241  else {
242  for (int igrandchild = 0; igrandchild < nGrandChildren; igrandchild++) {
243  currFamily.first = graphDaughter;
244  currFamily.second = currDaughter->getDaug(igrandchild);
245  heritancesQueue.push(currFamily);
246  }
247  }
248  }
249 
250  return graph.size() - existingParticles;
251 }
252 
253 
254 void EvtGenInterface::updateGraphParticle(EvtParticle* eParticle, MCParticleGraph::GraphParticle* gParticle,
255  TVector3 pPrimaryVertex, double timeOffset)
256 {
257  //updating the GraphParticle information from the EvtParticle information
258 
259  gParticle->setStatus(MCParticle::c_PrimaryParticle);
260  gParticle->setMass(eParticle->mass());
261  gParticle->setPDG(EvtPDL::getStdHep(eParticle->getId()));
262 
263  EvtVector4R EvtP4 = eParticle->getP4Lab();
264  TLorentzVector p4(EvtP4.get(1), EvtP4.get(2), EvtP4.get(3), EvtP4.get(0));
265  gParticle->set4Vector(p4);
266 
267  EvtVector4R Evtpos = eParticle->get4Pos();
268 
269  TVector3 pVertex(Evtpos.get(1)*Unit::mm, Evtpos.get(2)*Unit::mm, Evtpos.get(3)*Unit::mm);
270  pVertex = pVertex + pPrimaryVertex;
271 
272  gParticle->setProductionVertex(pVertex(0), pVertex(1), pVertex(2));
273  gParticle->setProductionTime((Evtpos.get(0)*Unit::mm / Const::speedOfLight) + timeOffset);
274  gParticle->setValidVertex(true);
275 
276  //add PHOTOS flag
277  if (eParticle->getAttribute("FSR")) {
278  gParticle->addStatus(MCParticle::c_IsPHOTOSPhoton);
279  }
280 
281 }
Helper file to create a temporary file and ensure deletion if object goes out of scope.
Definition: FileSystem.h:128
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.
void decaysInto(GraphParticle &daughter)
Tells the graph that this particle decays into daughter.
Class to build, validate and sort a particle decay chain.
size_t size() const
Return the number of particles in 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.
void setDecayTime(float time)
Set decay time.
Definition: MCParticle.h:387
void set4Vector(const TLorentzVector &p4)
Sets the 4Vector of particle.
Definition: MCParticle.h:435
void setMass(float mass)
Set particle mass.
Definition: MCParticle.h:363
void setProductionVertex(const TVector3 &vertex)
Set production vertex position.
Definition: MCParticle.h:393
void addStatus(unsigned short int bitmask)
Add bitmask to current status.
Definition: MCParticle.h:350
void setValidVertex(bool valid)
Set indication wether vertex and time information is valid or just default.
Definition: MCParticle.h:375
TVector3 getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition: MCParticle.h:183
void setPDG(int pdg)
Set PDG code of the particle.
Definition: MCParticle.h:332
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
TLorentzVector get4Vector() const
Return 4Vector of particle.
Definition: MCParticle.h:207
void setStatus(unsigned short int status)
Set Status code for the particle.
Definition: MCParticle.h:343
void setProductionTime(float time)
Set production time.
Definition: MCParticle.h:381
Evtgen random generator.
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.