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>
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>
36 EvtGenInterface::~EvtGenInterface()
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);
44 if (m_Generator)
delete m_Generator;
47 EvtGen* EvtGenInterface::createEvtGen(
const std::string& DECFileName,
const bool coherentMixing)
52 initLogCapture.
start();
53 EvtRandom::setRandomEngine(&EvtGenInterface::m_eng);
56 std::list<EvtDecayBase*> extraModels = EvtGenModelRegister::getModels();
61 EvtExternalGenList genList{
false};
62 EvtAbsRadCorr* radCorrEngine = genList.getPhotosModel();
63 list<EvtDecayBase*> modelList = genList.getListOfModels();
64 extraModels.insert(extraModels.end(), modelList.begin(), modelList.end());
67 EvtGenDatabasePDG::Instance()->WriteEvtGenTable(tmp);
69 auto mixingMode = EvtCPUtil::Incoherent;
71 mixingMode = EvtCPUtil::Coherent;
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.");
77 return new EvtGen(DECFileName.c_str(), tmp.getName().c_str(), &EvtGenInterface::m_eng,
78 radCorrEngine, &extraModels, mixingMode);
83 int EvtGenInterface::setup(
const std::string& DECFileName,
const std::string& parentParticle,
84 const std::string& userFileName,
const bool coherentMixing)
86 B2DEBUG(20,
"Begin initialisation of EvtGen Interface.");
91 initLogCapture.
start();
93 m_Generator = createEvtGen(DECFileName, coherentMixing);
95 if (!userFileName.empty()) {
96 m_Generator->readUDecay(userFileName.c_str());
100 if (parentParticle !=
"") {
101 m_ParentInitialized =
true;
102 m_ParentParticle = EvtPDL::getId(parentParticle);
106 B2DEBUG(20,
"End initialisation of EvtGen Interface.");
112 int EvtGenInterface::simulateEvent(
MCParticleGraph& graph, ROOT::Math::PxPyPzEVector pParentParticle,
113 ROOT::Math::XYZVector pPrimaryVertex,
114 int inclusiveType,
const std::string& inclusiveParticle)
116 EvtId inclusiveParticleID, inclusiveAntiParticleID;
118 if (!m_ParentInitialized)
119 B2FATAL(
"Parent particle is not initialized.");
121 m_pinit.set(pParentParticle.E(), pParentParticle.X(), pParentParticle.Y(), pParentParticle.Z());
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);
134 bool we_got_inclusive_particle =
false;
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();
142 if (inclusiveType != 0) {
143 EvtParticle* p = m_parent;
149 if (p->getId() == inclusiveParticleID ||
150 (inclusiveType == 2 && p->getId() == inclusiveAntiParticleID)) {
151 we_got_inclusive_particle =
true;
154 p = p->nextIter(m_parent);
157 if (!we_got_inclusive_particle) {
158 m_parent->deleteTree();
162 we_got_inclusive_particle =
true;
164 }
while (!we_got_inclusive_particle);
168 int iPart = addParticles2Graph(m_parent, graph, pPrimaryVertex, NULL);
169 graph.
generateList(
"", MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
171 m_parent->deleteTree();
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();
187 parent.
setDecayTime(-std::numeric_limits<float>::infinity());
189 id = EvtPDL::evtIdFromStdHep(pdg);
190 m_parent = EvtParticleFactory::particleFactory(
id, m_pinit);
192 m_parent->setVectorSpinDensity();
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();
202 int EvtGenInterface::addParticles2Graph(EvtParticle* top,
MCParticleGraph& graph, ROOT::Math::XYZVector pPrimaryVertex,
206 const int existingParticles = graph.
size();
212 updateGraphParticle(top, p, pPrimaryVertex);
214 typedef pair<MCParticleGraph::GraphParticle*, EvtParticle*> halfFamily;
215 halfFamily currFamily;
216 queue < halfFamily > heritancesQueue;
218 for (uint idaughter = 0; idaughter < top->getNDaug(); idaughter++) {
219 currFamily.first = p;
220 currFamily.second = top->getDaug(idaughter);
221 heritancesQueue.push(currFamily);
225 while (!heritancesQueue.empty()) {
226 currFamily = heritancesQueue.front();
227 heritancesQueue.pop();
230 EvtParticle* currDaughter = currFamily.second;
234 updateGraphParticle(currDaughter, graphDaughter, pPrimaryVertex, timeOffset);
239 int nGrandChildren = currDaughter->getNDaug();
241 if (nGrandChildren == 0)
242 graphDaughter->
addStatus(MCParticle::c_StableInGenerator);
244 for (
int igrandchild = 0; igrandchild < nGrandChildren; igrandchild++) {
245 currFamily.first = graphDaughter;
246 currFamily.second = currDaughter->getDaug(igrandchild);
247 heritancesQueue.push(currFamily);
252 return graph.
size() - existingParticles;
257 ROOT::Math::XYZVector pPrimaryVertex,
double timeOffset)
261 gParticle->
setStatus(MCParticle::c_PrimaryParticle);
262 gParticle->
setMass(eParticle->mass());
263 gParticle->
setPDG(EvtPDL::getStdHep(eParticle->getId()));
265 EvtVector4R EvtP4 = eParticle->getP4Lab();
266 ROOT::Math::PxPyPzEVector p4(EvtP4.get(1), EvtP4.get(2), EvtP4.get(3), EvtP4.get(0));
269 EvtVector4R Evtpos = eParticle->get4Pos();
271 ROOT::Math::XYZVector pVertex(Evtpos.get(1)*Unit::mm, Evtpos.get(2)*Unit::mm, Evtpos.get(3)*Unit::mm);
272 pVertex = pVertex + pPrimaryVertex;
275 gParticle->
setProductionTime((Evtpos.get(0)*Unit::mm / Const::speedOfLight) + timeOffset);
279 if (eParticle->getAttribute(
"FSR")) {
280 gParticle->
addStatus(MCParticle::c_IsPHOTOSPhoton);
Helper file to create a temporary file and ensure deletion if object goes out of scope.
bool start()
Start intercepting the output.
Capture stdout and stderr and convert into log messages.
void setIndent(const std::string &indent)
Set the indent for each line of the output, default is the supplied name + ": "
bool finish()
Finish the capture and emit the message if output has appeard on stdout or stderr.
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.
void setMass(float mass)
Set particle mass.
void addStatus(unsigned short int bitmask)
Add bitmask to current status.
ROOT::Math::XYZVector getVertex() const
Return production vertex position, shorthand for getProductionVertex().
void setValidVertex(bool valid)
Set indication wether vertex and time information is valid or just default.
void setProductionVertex(const ROOT::Math::XYZVector &vertex)
Set production vertex position.
ROOT::Math::PxPyPzEVector get4Vector() const
Return 4Vector of particle.
void setPDG(int pdg)
Set PDG code of the particle.
int getPDG() const
Return PDG code of particle.
float getProductionTime() const
Return production time in ns.
void set4Vector(const ROOT::Math::PxPyPzEVector &p4)
Sets the 4Vector of particle.
void setStatus(unsigned short int status)
Set Status code for the particle.
void setProductionTime(float time)
Set production time.
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.