Belle II Software  release-08-01-10
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 using namespace std;
32 using namespace Belle2;
33 
34 EvtGenFwRandEngine EvtGenInterface::m_eng;
35 
36 EvtGenInterface::~EvtGenInterface()
37 {
38  EvtDecayTable* evtDecayTable = EvtDecayTable::getInstance();
39  for (unsigned int i = 0; i < EvtPDL::entries(); ++i) {
40  for (int j = 0; j < evtDecayTable->getNMode(i); ++j) {
41  delete evtDecayTable->getDecay(i, j);
42  }
43  }
44  if (m_Generator) delete m_Generator;
45 }
46 
47 EvtGen* EvtGenInterface::createEvtGen(const std::string& DECFileName, const bool coherentMixing)
48 {
49  // Tauola prints normal things to stderr.. oh well.
50  IOIntercept::OutputToLogMessages initLogCapture("EvtGen", LogConfig::c_Info, LogConfig::c_Info);
51  initLogCapture.setIndent(" ");
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(20, "Begin initialisation of EvtGen Interface.");
87 
88  // Tauola prints normal things to stderr.. oh well.
89  IOIntercept::OutputToLogMessages initLogCapture("EvtGen", LogConfig::c_Info, LogConfig::c_Info);
90  initLogCapture.setIndent(" ");
91  initLogCapture.start();
92  if (!m_Generator) {
93  m_Generator = createEvtGen(DECFileName, coherentMixing);
94  }
95  if (!userFileName.empty()) {
96  m_Generator->readUDecay(userFileName.c_str());
97  }
98 
99  // Setup Parent Particle in rest frame
100  if (parentParticle != "") {
101  m_ParentInitialized = true;
102  m_ParentParticle = EvtPDL::getId(parentParticle);
103  }
104  initLogCapture.finish();
105 
106  B2DEBUG(20, "End initialisation of EvtGen Interface.");
107 
108  return 0;
109 }
110 
111 
112 int EvtGenInterface::simulateEvent(MCParticleGraph& graph, ROOT::Math::PxPyPzEVector pParentParticle,
113  ROOT::Math::XYZVector pPrimaryVertex,
114  int inclusiveType, const std::string& inclusiveParticle)
115 {
116  EvtId inclusiveParticleID, inclusiveAntiParticleID;
117 
118  if (!m_ParentInitialized)
119  B2FATAL("Parent particle is not initialized.");
120  //Init evtgen
121  m_pinit.set(pParentParticle.E(), pParentParticle.X(), pParentParticle.Y(), pParentParticle.Z());
122 
123  if (inclusiveType != 0) {
124  inclusiveParticleID = EvtPDL::getId(inclusiveParticle);
125  if (inclusiveParticleID.getId() < 0)
126  B2FATAL("Incorrect inclusive particle " << inclusiveParticle);
127  inclusiveAntiParticleID = EvtPDL::chargeConj(inclusiveParticleID);
128  if (inclusiveAntiParticleID.getId() < 0) {
129  B2FATAL("Cannot find the charge-conjugate particle for "
130  << inclusiveParticle);
131  }
132  }
133 
134  bool we_got_inclusive_particle = false;
135  do {
136  m_logCaptureDebug.start();
137  m_parent = EvtParticleFactory::particleFactory(m_ParentParticle, m_pinit);
138  m_parent->setVectorSpinDensity();
139  m_Generator->generateDecay(m_parent);
140  m_logCaptureDebug.finish();
141 
142  if (inclusiveType != 0) {
143  EvtParticle* p = m_parent;
144  // following loop will go through generated event and check it for
145  // presense of inclusive particle
146  do {
147  //for (int ii = 0; ii < iPart ; ii++) {
148  //std::cout << p->getPDGId() << std::endl;
149  if (p->getId() == inclusiveParticleID ||
150  (inclusiveType == 2 && p->getId() == inclusiveAntiParticleID)) {
151  we_got_inclusive_particle = true;
152  break;
153  }
154  p = p->nextIter(m_parent);
155  } while (p != 0);
156 
157  if (!we_got_inclusive_particle) {
158  m_parent->deleteTree();
159  }
160  } else {
161  // we don't do inclusive skimming so any generated event will do
162  we_got_inclusive_particle = true;
163  }
164  } while (!we_got_inclusive_particle);
165 
166  // B2INFO("after generate Decay.");
167 
168  int iPart = addParticles2Graph(m_parent, graph, pPrimaryVertex, NULL);
169  graph.generateList("", MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
170 
171  m_parent->deleteTree();
172 
173  return iPart; //returns the number of generated particles from evtgen
174 }
175 
176 int EvtGenInterface::simulateDecay(MCParticleGraph& graph,
178 {
179  int pdg;
180  EvtId id;
181  ROOT::Math::PxPyPzEVector momentum = parent.get4Vector();
182  ROOT::Math::XYZVector vertex = parent.getVertex();
183  m_pinit.set(momentum.E(), momentum.X(), momentum.Y(), momentum.Z());
184  m_logCaptureDebug.start();
185  // we want to decay the particle so the decay time in the tree needs to be lower
186  // than whatever the daughters will get
187  parent.setDecayTime(-std::numeric_limits<float>::infinity());
188  pdg = parent.getPDG();
189  id = EvtPDL::evtIdFromStdHep(pdg);
190  m_parent = EvtParticleFactory::particleFactory(id, m_pinit);
191  if (pdg == 10022) // Virtual photon
192  m_parent->setVectorSpinDensity();
193  else
194  m_parent->setDiagonalSpinDensity();
195  m_Generator->generateDecay(m_parent);
196  m_logCaptureDebug.finish();
197  int iPart = addParticles2Graph(m_parent, graph, vertex, &parent, parent.getProductionTime());
198  m_parent->deleteTree();
199  return iPart;
200 }
201 
202 int EvtGenInterface::addParticles2Graph(EvtParticle* top, MCParticleGraph& graph, ROOT::Math::XYZVector pPrimaryVertex,
203  MCParticleGraph::GraphParticle* parent, double timeOffset)
204 {
205  //Fill top particle in the tree & starting the queue:
206  const int existingParticles = graph.size();
208  if (parent == NULL)
209  p = &graph.addParticle();
210  else
211  p = parent;
212  updateGraphParticle(top, p, pPrimaryVertex);
213 
214  typedef pair<MCParticleGraph::GraphParticle*, EvtParticle*> halfFamily;
215  halfFamily currFamily;
216  queue < halfFamily > heritancesQueue;
217 
218  for (uint idaughter = 0; idaughter < top->getNDaug(); idaughter++) {
219  currFamily.first = p;
220  currFamily.second = top->getDaug(idaughter);
221  heritancesQueue.push(currFamily);
222  }
223 
224  //now we can go through the queue:
225  while (!heritancesQueue.empty()) {
226  currFamily = heritancesQueue.front(); //get the first entry from the queue
227  heritancesQueue.pop(); //remove the entry.
228 
229  MCParticleGraph::GraphParticle* currMother = currFamily.first;
230  EvtParticle* currDaughter = currFamily.second;
231 
232  //putting the daughter in the graph:
233  MCParticleGraph::GraphParticle* graphDaughter = &graph.addParticle();
234  updateGraphParticle(currDaughter, graphDaughter, pPrimaryVertex, timeOffset);
235 
236  //add relation between mother and daughter to graph:
237  currMother->decaysInto((*graphDaughter));
238 
239  int nGrandChildren = currDaughter->getNDaug();
240 
241  if (nGrandChildren == 0)
242  graphDaughter->addStatus(MCParticle::c_StableInGenerator);
243  else {
244  for (int igrandchild = 0; igrandchild < nGrandChildren; igrandchild++) {
245  currFamily.first = graphDaughter;
246  currFamily.second = currDaughter->getDaug(igrandchild);
247  heritancesQueue.push(currFamily);
248  }
249  }
250  }
251 
252  return graph.size() - existingParticles;
253 }
254 
255 
256 void EvtGenInterface::updateGraphParticle(EvtParticle* eParticle, MCParticleGraph::GraphParticle* gParticle,
257  ROOT::Math::XYZVector pPrimaryVertex, double timeOffset)
258 {
259  //updating the GraphParticle information from the EvtParticle information
260 
261  gParticle->setStatus(MCParticle::c_PrimaryParticle);
262  gParticle->setMass(eParticle->mass());
263  gParticle->setPDG(EvtPDL::getStdHep(eParticle->getId()));
264 
265  EvtVector4R EvtP4 = eParticle->getP4Lab();
266  ROOT::Math::PxPyPzEVector p4(EvtP4.get(1), EvtP4.get(2), EvtP4.get(3), EvtP4.get(0));
267  gParticle->set4Vector(p4);
268 
269  EvtVector4R Evtpos = eParticle->get4Pos();
270 
271  ROOT::Math::XYZVector pVertex(Evtpos.get(1)*Unit::mm, Evtpos.get(2)*Unit::mm, Evtpos.get(3)*Unit::mm);
272  pVertex = pVertex + pPrimaryVertex;
273 
274  gParticle->setProductionVertex(pVertex.x(), pVertex.y(), pVertex.z());
275  gParticle->setProductionTime((Evtpos.get(0)*Unit::mm / Const::speedOfLight) + timeOffset);
276  gParticle->setValidVertex(true);
277 
278  //add PHOTOS flag
279  if (eParticle->getAttribute("FSR")) {
280  gParticle->addStatus(MCParticle::c_IsPHOTOSPhoton);
281  }
282 
283 }
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
void setIndent(const std::string &indent)
Set the indent for each line of the output, default is the supplied name + ": "
Definition: IOIntercept.h:257
bool finish()
Finish the capture and emit the message if output has appeard on stdout or stderr.
Definition: IOIntercept.cc:236
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:390
void setMass(float mass)
Set particle mass.
Definition: MCParticle.h:366
void addStatus(unsigned short int bitmask)
Add bitmask to current status.
Definition: MCParticle.h:353
ROOT::Math::XYZVector getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition: MCParticle.h:183
void setValidVertex(bool valid)
Set indication wether vertex and time information is valid or just default.
Definition: MCParticle.h:378
void setProductionVertex(const ROOT::Math::XYZVector &vertex)
Set production vertex position.
Definition: MCParticle.h:396
ROOT::Math::PxPyPzEVector get4Vector() const
Return 4Vector of particle.
Definition: MCParticle.h:207
void setPDG(int pdg)
Set PDG code of the particle.
Definition: MCParticle.h:335
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets the 4Vector of particle.
Definition: MCParticle.h:438
void setStatus(unsigned short int status)
Set Status code for the particle.
Definition: MCParticle.h:346
void setProductionTime(float time)
Set production time.
Definition: MCParticle.h:384
Evtgen random generator.
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.