![]() |
|
|||
File indexing completed on 2025-08-03 08:20:09
0001 /******************************************************************************* 0002 * Copyright (c) The JETSCAPE Collaboration, 2018 0003 * 0004 * Modular, task-based framework for simulating all aspects of heavy-ion collisions 0005 * 0006 * For the list of contributors see AUTHORS. 0007 * 0008 * Report issues at https://github.com/JETSCAPE/JETSCAPE/issues 0009 * 0010 * or via email to bugs.jetscape@gmail.com 0011 * 0012 * Distributed under the GNU General Public License 3.0 (GPLv3 or later). 0013 * See COPYING for details. 0014 ******************************************************************************/ 0015 0016 /******************************************************************************* 0017 * kk Sep 22, 2020: Significant updates to status codes 0018 * to conform with Appendix A in https://arxiv.org/pdf/1912.08005.pdf 0019 * We will follow Pythia's example of preserving as much of the internal status code as possible 0020 * Specifically: 0021 - beam particles. We don't have those, yet. If they ever become part of initial 0022 state moduls, those module writers MUST set the status to 4, all we do here 0023 is respect that code (for the future). 0024 However, this particle won't be a parton, so we only check that for hadrons. 0025 In practice, this will probably require a revamp of the graph structure, both 0026 in the framework and in here. Then, probably also use GenEvent::add_beam_particle() 0027 - decayed particle. We don't have those either yet, but it's a feature that may 0028 come pretty soon. 0029 Here, we use that exclusively to mean decayed unstable hadrons, 0030 e.g. K0S -> pi pi 0031 In this case, the module that produced the hadron list MUST ensure to set the 0032 status of K0s to 2. The pions are final particles, see below. 0033 - Final hadrons will all be forced to status 1 0034 - Partons: By default, we will accept and use any code provided by the framework 0035 within 11<=status<=200 and use the absolute value. 0036 Parton exceptions: 0037 1) If the code is not HepMC-legal (|status|<11 or >200, often 0), 0038 we assign 12 to most and 11 to the final partons before hadronization 0039 (to mimic the 1, 2 scheme) 0040 2) If there are NO hadrons in the event, we assume the user would like to treat 0041 the final partons (like 11 from above) as final particles and assign 1 0042 ******************************************************************************/ 0043 0044 0045 #include "JetScapeWriterRootHepMC.h" 0046 #include "JetScapeLogger.h" 0047 #include "HardProcess.h" 0048 #include "JetScapeSignalManager.h" 0049 #include "GTL/node.h" 0050 #include <GTL/topsort.h> 0051 0052 using HepMC3::Units; 0053 0054 namespace Jetscape { 0055 0056 JetScapeWriterRootHepMC::~JetScapeWriterRootHepMC() { 0057 if (GetActive()) 0058 Close(); 0059 } 0060 0061 void JetScapeWriterRootHepMC::WriteHeaderToFile() { 0062 // Create event here - not actually writing 0063 // TODO: GeV seems right, but I don't think we actually measure in mm 0064 // Should multiply all lengths by 1e-12 probably 0065 evt = GenEvent(Units::GEV, Units::MM); 0066 0067 // Expects pb, pythia delivers mb 0068 auto xsec = make_shared<HepMC3::GenCrossSection>(); 0069 xsec->set_cross_section(GetHeader().GetSigmaGen() * 1e9, 0); 0070 xsec->set_cross_section(GetHeader().GetSigmaGen() * 1e9, 0071 GetHeader().GetSigmaErr() * 1e9); 0072 evt.set_cross_section(xsec); 0073 evt.weights().push_back(GetHeader().GetEventWeight()); 0074 0075 auto heavyion = make_shared<HepMC3::GenHeavyIon>(); 0076 // see https://gitlab.cern.ch/hepmc/HepMC3/blob/master/include/HepMC/GenHeavyIon.h 0077 if (GetHeader().GetNpart() > -1) { 0078 // Not clear what the difference is... 0079 heavyion->Ncoll_hard = GetHeader().GetNcoll(); 0080 heavyion->Ncoll = GetHeader().GetNcoll(); 0081 } 0082 if (GetHeader().GetNcoll() > -1) { 0083 // Hepmc separates into target and projectile. 0084 // Set one? Which? Both? half to each? setting projectile for now. 0085 // setting both might lead to weird problems when they get added up 0086 heavyion->Npart_proj = GetHeader().GetNpart(); 0087 } 0088 if (GetHeader().GetTotalEntropy() > -1) { 0089 // nothing good in the HepMC standard. Something related to mulitplicity would work 0090 } 0091 0092 if (GetHeader().GetEventPlaneAngle() > -999) { 0093 heavyion->event_plane_angle = GetHeader().GetEventPlaneAngle(); 0094 } 0095 0096 evt.set_heavy_ion(heavyion); 0097 0098 // also a good moment to initialize the hadron boolean 0099 } 0100 0101 void JetScapeWriterRootHepMC::WriteEvent() { 0102 VERBOSE(1) << "Run JetScapeWriterHepMC: Write event # " << GetCurrentEvent(); 0103 0104 // Have collected all vertices now. 0105 // Add all vertices to the event 0106 for (auto v : vertices){ 0107 evt.add_vertex(v); 0108 } 0109 0110 VERBOSE(1) << " found " << vertices.size() << " vertices in the list"; 0111 0112 // If there are no hadrons, promote final partons 0113 // The graph support of hepmc is a bit rudimentary. 0114 // easiest is to just check all childless particles 0115 // Note, one could just check for status==11, 0116 // but modules are allowed to assign that number to non-final partons 0117 if ( !hashadrons ) { 0118 VERBOSE(1) << " found no hadrons, promoting final partons to status 1"; 0119 for ( auto p : evt.particles() ){ 0120 if ( p->children().size() == 0 ){ 0121 if ( p->status() !=11 ){ 0122 JSWARN << "Found a final parton with status!=11 : status=" << p->status() << ". This should not happen"; 0123 } 0124 p->set_status(1); 0125 } 0126 } 0127 } 0128 evt.set_event_number(GetCurrentEvent()); 0129 write_event(evt); 0130 vertices.clear(); 0131 hadronizationvertex = 0; 0132 } 0133 0134 //This function dumps the particles in a specific parton shower to the event 0135 void JetScapeWriterRootHepMC::Write(weak_ptr<PartonShower> ps) { 0136 shared_ptr<PartonShower> pShower = ps.lock(); 0137 if (!pShower) 0138 return; 0139 0140 // Need topological order, see 0141 // https://hepmc.web.cern.ch/hepmc/differences.html 0142 // That means if parton p1 comes into vertex v, and p2 goes out of v, 0143 // then p1 has to be created (bestowed an id) before p2 0144 // Take inspiration from hepmc3.0.0/interfaces/pythia8/src/Pythia8ToHepMC3.cc 0145 // But pythia showers are different from our existing graph structure, 0146 // So instead try to modify the first attempt to respect top. order 0147 // and don't create vertices and particles more than once 0148 0149 // Using GTL's topsort 0150 // 1. Check that our graph is sane 0151 if (!pShower->is_acyclic()) 0152 throw std::runtime_error( 0153 "PROBLEM in JetScapeWriterHepMC: Graph is not acyclic."); 0154 0155 // 2. 0156 topsort topsortsearch; 0157 topsortsearch.scan_whole_graph(true); 0158 topsortsearch.start_node(); // defaults to first node 0159 topsortsearch.run(*pShower); 0160 auto nEnd = 0161 topsortsearch.top_order_end(); // this is a topsort::topsort_iterator 0162 0163 // Need to keep track of already created ones 0164 map<int, GenParticlePtr> CreatedPartons; 0165 0166 // probably unnecessary, used for consistency checks 0167 map<int, bool> vused; 0168 bool foundRoot = false; 0169 for (auto nIt = topsortsearch.top_order_begin(); nIt != nEnd; ++nIt) { 0170 // cout << *nIt << " " << nIt->indeg() << " " << nIt->outdeg() << endl; 0171 0172 // Should be the only time we see this node. 0173 if (vused.find(nIt->id()) != vused.end()){ 0174 throw std::runtime_error( "PROBLEM in JetScapeWriterHepMC: Reusing a vertex."); 0175 } 0176 vused[nIt->id()] = true; 0177 0178 // 0. No incoming edges? 0179 // --------------------- 0180 // This is typically a shower initiator. 0181 // HepMC needs an incoming and outgoing particle for every vertex. 0182 // That could be a place to attach partons or ions. 0183 // We could also do both, but as of now, JETSCAPE actually attaches a dummy vertex 0184 // to the start of the initiator, that can safely go away. 0185 // Previously, we attached a dummy or clone of the outgoing one here. 0186 // Instead, we can just skip the vertex. Its outgoing edges will be picked up 0187 // as incomers in a later vertex. 0188 // Note that the [0]=>[1] connection in JETSCAPE 0189 // already uses a dummy node[0], and [1] is at time t=0; removing that seems correct. 0190 if (nIt->indeg() == 0) continue; 0191 0192 // 1. Create a new vertex. 0193 // -------------------------------------------- 0194 auto v = castVtxToHepMC(pShower->GetVertex(*nIt)); 0195 0196 // 2. Incoming edges 0197 // ----------------- 0198 // In the current framework, it should only be one. 0199 // So we will catch anything more but provide a mechanism that should work anyway. 0200 if (nIt->indeg() > 1) { 0201 JSWARN << "Found more than one mother parton! Should only happen if we " 0202 "added medium particles. " 0203 << "The code should work, but proceed with caution"; 0204 } 0205 0206 auto inIt = nIt->in_edges_begin(); 0207 auto inEnd = nIt->in_edges_end(); 0208 for (/* nop */; inIt != inEnd; ++inIt) { 0209 auto phepin = CreatedPartons.find(inIt->id()); 0210 if (phepin != CreatedPartons.end()) { 0211 // We should already have one! 0212 v->add_particle_in(phepin->second); 0213 } else { 0214 // This indicates we skipped an earlier vertex without incomers. 0215 // JSWARN << "Incoming particle out of nowhere. This could maybe happen " 0216 // "if we pick up medium particles " 0217 // << " but is probably a topsort problem. Try using the code " 0218 // "after this throw() but be very careful."; 0219 // throw std::runtime_error("PROBLEM in JetScapeWriterHepMC: Incoming " 0220 // "particle out of nowhere."); 0221 0222 auto in = pShower->GetParton(*inIt); 0223 auto hepin = castPartonToHepMC(in); 0224 auto status = std::abs(hepin->status()); 0225 if ( status < 11 || status > 200) { 0226 // incoming edge can't be final 0227 status = 12; 0228 } 0229 hepin->set_status(status); 0230 CreatedPartons[inIt->id()] = hepin; 0231 v->add_particle_in(hepin); 0232 0233 if ( nIt->outdeg() == 0 ) { 0234 // However, motherless AND childless particles do exist 0235 // I.e., a shower initiator that never actually showers 0236 // For this, we need an out going clone, much like 3) below 0237 auto hepout = castPartonToHepMC(in); 0238 // Since the status information is preserved in the incomer, we'll force 11 0239 // Note: if we later see in WriteEvent() that there are no hadrons, this will be overwritten to 1 0240 hepout->set_status( 11 ); 0241 // Note that we do not register this particle. Since it's pointing nowhere it can never be reused. 0242 v->add_particle_out(hepout); 0243 } 0244 } 0245 } 0246 0247 // 3. Outgoing edges? 0248 // -------------------------------------------- 0249 // 3.1: No. Need to create one. 0250 // We'll use this opportunity to copy the incomer but give it a final code 0251 if (nIt->outdeg() == 0) { 0252 if (nIt->indeg() != 1) { 0253 // This won't work with multiple incomers (but that's pretty unphysical) 0254 throw std::runtime_error("PROBLEM in JetScapeWriterHepMC: Need exactly " 0255 "one parent to clone final state partons."); 0256 } 0257 auto in = pShower->GetParton(*(nIt->in_edges_begin())); 0258 auto hepout = castPartonToHepMC(in); 0259 // an outgoing edge without terminator is "final" 0260 // Since the status information is preserved in the incomer, we'll force 11 0261 // Note: if we later see in WriteEvent() that there are no hadrons, this will be overwritten to 1 0262 hepout->set_status( 11 ); 0263 v->add_particle_out(hepout); 0264 // Note that we do not register this particle. Since it's pointing nowhere it can never be reused. 0265 } 0266 0267 // 3.2: Otherwise use and register the outgoing edge 0268 if (nIt->outdeg() > 0) { 0269 auto outIt = nIt->out_edges_begin(); 0270 auto outEnd = nIt->out_edges_end(); 0271 for (/* nop */; outIt != outEnd; ++outIt) { 0272 if (CreatedPartons.find(outIt->id()) != CreatedPartons.end()) { 0273 throw std::runtime_error("PROBLEM in JetScapeWriterHepMC: Trying to " 0274 "recreate a preexisting GenParticle."); 0275 } 0276 auto out = pShower->GetParton(*outIt); 0277 auto hepout = castPartonToHepMC(out); 0278 if ( !hepout->status()) { 0279 // incoming and outgoing -> status 12 0280 hepout->set_status(12); 0281 } 0282 0283 CreatedPartons[outIt->id()] = hepout; 0284 v->add_particle_out(hepout); 0285 } 0286 } 0287 0288 vertices.push_back(v); 0289 } 0290 } 0291 0292 void JetScapeWriterRootHepMC::Write(weak_ptr<Hadron> h) { 0293 auto hadron = h.lock(); 0294 if (!hadron) 0295 return; 0296 0297 // No clear source for most hadrons 0298 // Also, a graph with e.g. recombination hadrons would have loops, 0299 // (though the direction should still make it acyclic?) 0300 // Not sure how this is supposed to be done in HepMC3 0301 // Our solution: Attach all hadrons to one dedicated hadronization vertex. 0302 // Future option: Have separate shower and bulk vertices? 0303 0304 // Create if it doesn't exist yet 0305 if (!hadronizationvertex) { 0306 // dummy position 0307 HepMC3::FourVector vtxPosition(0, 0, 0, 100); // set it to a late time... 0308 hadronizationvertex = make_shared<GenVertex>(vtxPosition); 0309 0310 // dummy mother -- could also maybe use the first/hardest shower initiator 0311 HepMC3::FourVector pmom(0, 0, 0, 0); 0312 make_shared<GenParticle>(pmom, 0, 0); 0313 hadronizationvertex->add_particle_in(make_shared<GenParticle>(pmom, 0, 0)); 0314 0315 vertices.push_back(hadronizationvertex); 0316 hashadrons=true; 0317 } 0318 0319 // now attach 0320 auto hepmc = castHadronToHepMC(hadron); 0321 if ( !hepmc->status() ) { 0322 // unless otherwise specified, all hadrons get status 1 0323 // TODO: Need to better account for short-lived hadrons 0324 hepmc->set_status(1); 0325 } 0326 hadronizationvertex->add_particle_out(hepmc); 0327 } 0328 0329 void JetScapeWriterRootHepMC::Init() { 0330 if (GetActive()) { 0331 JSINFO << "JetScape HepMC Writer initialized with output file = " 0332 << GetOutputFileName(); 0333 } 0334 } 0335 0336 void JetScapeWriterRootHepMC::Exec() { 0337 // Nothing to do 0338 } 0339 } // end namespace Jetscape
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |