File indexing completed on 2025-09-17 09:08:38
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef Pythia8_HepMC3_H
0012 #define Pythia8_HepMC3_H
0013
0014 #ifdef Pythia8_HepMC2_H
0015 #error Cannot include HepMC3.h if HepMC2.h has already been included.
0016 #endif
0017
0018 #include <vector>
0019 #include "Pythia8/Pythia.h"
0020 #include "Pythia8/HIInfo.h"
0021 #include "HepMC3/GenVertex.h"
0022 #include "HepMC3/GenParticle.h"
0023 #include "HepMC3/GenEvent.h"
0024 #include "HepMC3/WriterAscii.h"
0025 #include "HepMC3/WriterAsciiHepMC2.h"
0026 #include "HepMC3/GenHeavyIon.h"
0027 #include "HepMC3/GenPdfInfo.h"
0028
0029 namespace HepMC3 {
0030
0031
0032
0033 class Pythia8ToHepMC3 {
0034
0035 public:
0036
0037
0038 Pythia8ToHepMC3(): m_internal_event_number(0), m_print_inconsistency(true),
0039 m_free_parton_warnings(true), m_crash_on_problem(false),
0040 m_convert_gluon_to_0(false), m_store_pdf(true), m_store_proc(true),
0041 m_store_xsec(true), m_store_weights(true) {}
0042 virtual ~Pythia8ToHepMC3() {}
0043
0044
0045 bool fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt,
0046 int ievnum = -1 ) { return fill_next_event( pythia.event, evt,
0047 ievnum, &pythia.info, &pythia.settings); }
0048 bool fill_next_event( Pythia8::Pythia& pythia, GenEvent& evt) {
0049 return fill_next_event( pythia, &evt); }
0050
0051
0052 bool fill_next_event( Pythia8::Event& pyev, GenEvent&evt, int ievnum = -1,
0053 const Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0) {
0054 return fill_next_event(pyev, &evt, ievnum, pyinfo, pyset); }
0055 bool fill_next_event( Pythia8::Event& pyev, GenEvent* evt, int ievnum = -1,
0056 const Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0) {
0057
0058
0059 if (evt == nullptr) return warning(pyinfo,
0060 "Pythia8ToHepMC::fill_next_event", "passed null event");
0061
0062
0063 if ( ievnum >= 0 ) {
0064 evt->set_event_number(ievnum);
0065 m_internal_event_number = ievnum;
0066 }
0067 else {
0068 evt->set_event_number(m_internal_event_number);
0069 ++m_internal_event_number;
0070 }
0071
0072
0073 evt->set_units(Units::GEV,Units::MM);
0074
0075
0076 if ( pyinfo && pyinfo->hiInfo ) {
0077 auto ion = make_shared<HepMC3::GenHeavyIon>();
0078 ion->Ncoll_hard = pyinfo->hiInfo->nCollND();
0079 ion->Ncoll = pyinfo->hiInfo->nCollTot();
0080 ion->Npart_proj = pyinfo->hiInfo->nAbsProj() +
0081 pyinfo->hiInfo->nDiffProj();
0082 ion->Npart_targ = pyinfo->hiInfo->nAbsTarg() +
0083 pyinfo->hiInfo->nDiffTarg();
0084 ion->impact_parameter = pyinfo->hiInfo->b();
0085 evt->set_heavy_ion(ion);
0086 }
0087
0088
0089 std::vector<GenParticlePtr> hepevt_particles;
0090 hepevt_particles.reserve( pyev.size() );
0091 for(int i = 0; i < pyev.size(); ++i) {
0092 hepevt_particles.push_back( std::make_shared<GenParticle>(
0093 FourVector( pyev[i].px(), pyev[i].py(), pyev[i].pz(), pyev[i].e() ),
0094 pyev[i].id(), pyev[i].statusHepMC() ) );
0095 hepevt_particles[i]->set_generated_mass( pyev[i].m() );
0096 }
0097
0098
0099 std::vector<GenVertexPtr> vertex_cache;
0100 vector<GenParticlePtr> beam_particles;
0101 for (int i = 1; i < pyev.size(); ++i) {
0102 vector<int> mothers = pyev[i].motherList();
0103 sort(mothers.begin(),mothers.end());
0104 for (;;) {
0105 if (!mothers.empty() && mothers.front() == 0)
0106 mothers.erase(mothers.begin());
0107 else break;
0108 }
0109 if (mothers.size()) {
0110 GenVertexPtr prod_vtx = hepevt_particles[mothers[0]]->end_vertex();
0111 if (!prod_vtx) {
0112 prod_vtx = make_shared<GenVertex>();
0113 vertex_cache.push_back(prod_vtx);
0114 for (unsigned int j = 0; j < mothers.size(); ++j)
0115 prod_vtx->add_particle_in( hepevt_particles[mothers[j]] );
0116 }
0117 FourVector prod_pos( pyev[i].xProd(), pyev[i].yProd(),pyev[i].zProd(),
0118 pyev[i].tProd() );
0119
0120
0121 if (!prod_pos.is_zero() && prod_vtx->position().is_zero())
0122 prod_vtx->set_position( prod_pos );
0123 prod_vtx->add_particle_out( hepevt_particles[i] );
0124 } else beam_particles.push_back(hepevt_particles[i]);
0125 }
0126
0127
0128 evt->reserve( hepevt_particles.size(), vertex_cache.size() );
0129
0130
0131 evt->add_tree( beam_particles );
0132
0133
0134 for (int i = 0; i < pyev.size(); ++i) {
0135
0136
0137 int colType = pyev[i].colType();
0138 if (colType == -1 ||colType == 1 || colType == 2) {
0139 int flow1 = 0, flow2 = 0;
0140 if (colType == 1 || colType == 2) flow1 = pyev[i].col();
0141 if (colType == -1 || colType == 2) flow2 = pyev[i].acol();
0142 hepevt_particles[i]->add_attribute("flow1",
0143 make_shared<IntAttribute>(flow1));
0144 hepevt_particles[i]->add_attribute("flow2",
0145 make_shared<IntAttribute>(flow2));
0146 }
0147 }
0148
0149
0150 bool doHadr = (pyset == 0) ? m_free_parton_warnings
0151 : pyset->flag("HadronLevel:all") && pyset->flag("HadronLevel:Hadronize");
0152
0153
0154
0155
0156 for (int i = 1; i < pyev.size(); ++i) {
0157
0158
0159
0160
0161 if ( hepevt_particles[i] == nullptr ||
0162 !hepevt_particles[i]->in_event()) {
0163 warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
0164 "found orphan particle", "i = " + Pythia8::toString(i));
0165 GenVertexPtr prod_vtx = make_shared<GenVertex>();
0166 prod_vtx->add_particle_out( hepevt_particles[i] );
0167 evt->add_vertex(prod_vtx);
0168 }
0169
0170
0171 if ( doHadr && m_free_parton_warnings ) {
0172 if ( hepevt_particles[i]->pid() == 21
0173 && hepevt_particles[i]->end_vertex() == nullptr ) {
0174 warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
0175 "found gluon without end vertex", "i = " + Pythia8::toString(i));
0176 if ( m_crash_on_problem ) exit(1);
0177 }
0178 if ( abs(hepevt_particles[i]->pid()) <= 6
0179 && hepevt_particles[i]->end_vertex() == nullptr ) {
0180 warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
0181 "found quark without end vertex", "i = " + Pythia8::toString(i));
0182 if ( m_crash_on_problem ) exit(1);
0183 }
0184 }
0185 }
0186
0187
0188
0189 if (m_store_pdf && pyinfo != 0) {
0190 int id1pdf = pyinfo->id1pdf();
0191 int id2pdf = pyinfo->id2pdf();
0192 if ( m_convert_gluon_to_0 ) {
0193 if (id1pdf == 21) id1pdf = 0;
0194 if (id2pdf == 21) id2pdf = 0;
0195 }
0196
0197
0198 GenPdfInfoPtr pdfinfo = make_shared<GenPdfInfo>();
0199 pdfinfo->set(id1pdf, id2pdf, pyinfo->x1pdf(), pyinfo->x2pdf(),
0200 pyinfo->QFac(), pyinfo->pdf1(), pyinfo->pdf2() );
0201 evt->set_pdf_info( pdfinfo );
0202 }
0203
0204
0205 if (m_store_proc && pyinfo != 0) {
0206 evt->add_attribute("signal_process_id",
0207 std::make_shared<IntAttribute>( pyinfo->code()));
0208 evt->add_attribute("mpi",
0209 std::make_shared<IntAttribute>( pyinfo->nMPI()));
0210 evt->add_attribute("event_scale",
0211 std::make_shared<DoubleAttribute>(pyinfo->QRen()));
0212 evt->add_attribute("alphaQCD",
0213 std::make_shared<DoubleAttribute>(pyinfo->alphaS()));
0214 evt->add_attribute("alphaQED",
0215 std::make_shared<DoubleAttribute>(pyinfo->alphaEM()));
0216 }
0217
0218
0219 if (m_store_weights && pyinfo != 0) {
0220 evt->weights().clear();
0221 for (int iWeight = 0; iWeight < pyinfo->numberOfWeights(); ++iWeight)
0222 evt->weights().push_back(pyinfo->weightValueByIndex(iWeight));
0223 }
0224
0225
0226 if (m_store_xsec && pyinfo != 0) {
0227
0228
0229
0230 GenCrossSectionPtr xsec = make_shared<GenCrossSection>();
0231 evt->set_cross_section(xsec);
0232 xsec->set_cross_section( pyinfo->sigmaGen() * 1e9,
0233 pyinfo->sigmaErr() * 1e9);
0234
0235 vector<double> xsecVec = pyinfo->weightContainerPtr->getTotalXsec();
0236 if (xsecVec.size() > 0) {
0237 for (unsigned int iXsec = 0; iXsec < xsecVec.size(); ++iXsec) {
0238 xsec->set_xsec(iXsec, xsecVec[iXsec]*1e9);
0239 }
0240 }
0241 }
0242
0243
0244 return true;
0245 }
0246
0247
0248 bool print_inconsistency() const { return m_print_inconsistency; }
0249 bool free_parton_warnings() const { return m_free_parton_warnings; }
0250 bool crash_on_problem() const { return m_crash_on_problem; }
0251 bool convert_gluon_to_0() const { return m_convert_gluon_to_0; }
0252 bool store_pdf() const { return m_store_pdf; }
0253 bool store_proc() const { return m_store_proc; }
0254 bool store_xsec() const { return m_store_xsec; }
0255 bool store_weights() const { return m_store_weights; }
0256
0257
0258 void set_print_inconsistency(bool b = true) { m_print_inconsistency = b; }
0259 void set_free_parton_warnings(bool b = true) { m_free_parton_warnings = b; }
0260 void set_crash_on_problem(bool b = false) { m_crash_on_problem = b; }
0261 void set_convert_gluon_to_0(bool b = false) { m_convert_gluon_to_0 = b; }
0262 void set_store_pdf(bool b = true) { m_store_pdf = b; }
0263 void set_store_proc(bool b = true) { m_store_proc = b; }
0264 void set_store_xsec(bool b = true) { m_store_xsec = b; }
0265 void set_store_weights(bool b = true) { m_store_weights = b; }
0266
0267 private:
0268
0269
0270
0271 bool warning(const Pythia8::Info * pyinfo, string loc,
0272 string message, string extraInfo = "") {
0273 if ( pyinfo )
0274 pyinfo->loggerPtr->warningMsg(loc, message, extraInfo);
0275 else if ( print_inconsistency() )
0276 cout << "Warning in " << loc << ": " << message << extraInfo << endl;
0277 return false;
0278 }
0279
0280
0281 virtual bool fill_next_event( GenEvent* ) { return 0; }
0282 virtual void write_event( const GenEvent* ) {}
0283
0284
0285 Pythia8ToHepMC3( const Pythia8ToHepMC3& ) {}
0286
0287
0288 int m_internal_event_number;
0289 bool m_print_inconsistency, m_free_parton_warnings, m_crash_on_problem,
0290 m_convert_gluon_to_0, m_store_pdf, m_store_proc, m_store_xsec,
0291 m_store_weights;
0292
0293 };
0294
0295
0296
0297 }
0298
0299 namespace Pythia8 {
0300
0301
0302
0303
0304
0305
0306
0307
0308 class Pythia8ToHepMC : public HepMC3::Pythia8ToHepMC3 {
0309
0310 public:
0311
0312
0313
0314 enum OutputType { none, ascii2, ascii3 };
0315
0316
0317 typedef HepMC3::GenEvent GenEvent;
0318 typedef shared_ptr<GenEvent> EventPtr;
0319 typedef HepMC3::Writer Writer;
0320 typedef shared_ptr<Writer> WriterPtr;
0321
0322
0323 Pythia8ToHepMC() : runinfo(make_shared<HepMC3::GenRunInfo>()) {}
0324
0325
0326 Pythia8ToHepMC(string filename, OutputType ft = ascii3)
0327 : runinfo(make_shared<HepMC3::GenRunInfo>()) {
0328 setNewFile(filename, ft);
0329 }
0330
0331
0332 bool setNewFile(string filename, OutputType ft = ascii3) {
0333 switch ( ft ) {
0334 case ascii3:
0335 writerPtr = make_shared<HepMC3::WriterAscii>(filename);
0336 break;
0337 case ascii2:
0338 writerPtr = make_shared<HepMC3::WriterAsciiHepMC2>(filename);
0339 break;
0340 case none:
0341 break;
0342 }
0343 return writerPtr != nullptr;
0344 }
0345
0346
0347
0348 bool fillNextEvent(Pythia & pythia) {
0349 geneve = make_shared<HepMC3::GenEvent>(runinfo);
0350 if (runinfo->weight_names().size() == 0)
0351 setWeightNames(pythia.info.weightNameVector());
0352 return fill_next_event(pythia, *geneve);
0353 }
0354
0355
0356 void writeEvent() {
0357 writerPtr->write_event(*geneve);
0358 }
0359
0360
0361
0362
0363 bool writeNextEvent(Pythia & pythia) {
0364 if ( !fillNextEvent(pythia) ) return false;
0365 writeEvent();
0366 return !writerPtr->failed();
0367 }
0368
0369
0370 GenEvent & event() {
0371 return *geneve;
0372 }
0373
0374
0375 EventPtr getEventPtr() {
0376 return geneve;
0377 }
0378
0379
0380 Writer & output() {
0381 return *writerPtr;
0382 }
0383
0384
0385 WriterPtr outputPtr() {
0386 return writerPtr;
0387 }
0388
0389
0390 void setXSec(double xsec, double xsecerr) {
0391 auto xsecptr = geneve->cross_section();
0392 if ( !xsecptr ) {
0393 xsecptr = make_shared<HepMC3::GenCrossSection>();
0394 geneve->set_cross_section(xsecptr);
0395 }
0396 xsecptr->set_cross_section(xsec, xsecerr);
0397 }
0398
0399
0400 void setWeights(const vector<double> & wv) {
0401 geneve->weights() = wv;
0402 }
0403
0404
0405 void setWeightNames(const vector<string> &wnv) {
0406 runinfo->set_weight_names(wnv);
0407 }
0408
0409
0410 void setPdfInfo(int id1, int id2, double x1, double x2,
0411 double scale, double xf1, double xf2,
0412 int pdf1 = 0, int pdf2 = 0) {
0413 auto pdf = make_shared<HepMC3::GenPdfInfo>();
0414 pdf->set(id1, id2, x1, x2, scale, xf1, xf2, pdf1, pdf2);
0415 geneve->set_pdf_info(pdf);
0416 }
0417
0418
0419
0420 template<class T>
0421 void addAtribute(const string& name, T& attribute) {
0422 shared_ptr<HepMC3::Attribute> att = make_shared<T>(attribute);
0423 geneve->add_attribute(name, att);
0424 }
0425
0426
0427 template<class T=double>
0428 void addAttribute(const string& name, double& attribute) {
0429 auto dAtt = HepMC3::DoubleAttribute(attribute);
0430 shared_ptr<HepMC3::Attribute> att =
0431 make_shared<HepMC3::DoubleAttribute>(dAtt);
0432 geneve->add_attribute(name, att);
0433 }
0434
0435
0436 void removeAttribute(const string& name) {
0437 geneve->remove_attribute(name);
0438 }
0439
0440 private:
0441
0442
0443 EventPtr geneve = nullptr;
0444
0445
0446 WriterPtr writerPtr = nullptr;
0447
0448
0449 shared_ptr<HepMC3::GenRunInfo> runinfo;
0450
0451 };
0452
0453 }
0454
0455 #endif