PYTHIA  8.313
HepMC3.h
1 // HepMC3.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 //
6 // Author: HepMC 3 Collaboration, hepmc-dev@.cern.ch
7 // Based on the HepMC2 interface by Mikhail Kirsanov, Mikhail.Kirsanov@cern.ch.
8 // Header file and function definitions for the Pythia8ToHepMC class,
9 // which converts a PYTHIA event record to the standard HepMC format.
10 
11 #ifndef Pythia8_HepMC3_H
12 #define Pythia8_HepMC3_H
13 
14 #ifdef Pythia8_HepMC2_H
15 #error Cannot include HepMC3.h if HepMC2.h has already been included.
16 #endif
17 
18 #include <vector>
19 #include "Pythia8/Pythia.h"
20 #include "Pythia8/HIInfo.h"
21 #include "HepMC3/GenVertex.h"
22 #include "HepMC3/GenParticle.h"
23 #include "HepMC3/GenEvent.h"
24 #include "HepMC3/WriterAscii.h"
25 #include "HepMC3/WriterAsciiHepMC2.h"
26 #include "HepMC3/GenHeavyIon.h"
27 #include "HepMC3/GenPdfInfo.h"
28 
29 namespace HepMC3 {
30 
31 //==========================================================================
32 
34 
35 public:
36 
37  // Constructor and destructor
38  Pythia8ToHepMC3(): m_internal_event_number(0), m_print_inconsistency(true),
39  m_free_parton_warnings(true), m_crash_on_problem(false),
40  m_convert_gluon_to_0(false), m_store_pdf(true), m_store_proc(true),
41  m_store_xsec(true), m_store_weights(true) {}
42  virtual ~Pythia8ToHepMC3() {}
43 
44  // The recommended method to convert Pythia events into HepMC3 ones.
45  bool fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt,
46  int ievnum = -1 ) { return fill_next_event( pythia.event, evt,
47  ievnum, &pythia.info, &pythia.settings); }
48  bool fill_next_event( Pythia8::Pythia& pythia, GenEvent& evt) {
49  return fill_next_event( pythia, &evt); }
50 
51  // Alternative method to convert Pythia events into HepMC3 ones.
52  bool fill_next_event( Pythia8::Event& pyev, GenEvent&evt, int ievnum = -1,
53  const Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0) {
54  return fill_next_event(pyev, &evt, ievnum, pyinfo, pyset); }
55  bool fill_next_event( Pythia8::Event& pyev, GenEvent* evt, int ievnum = -1,
56  const Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0) {
57 
58  // 1. Error if no event passed.
59  if (evt == nullptr) return warning(pyinfo,
60  "Pythia8ToHepMC::fill_next_event", "passed null event");
61 
62  // Event number counter.
63  if ( ievnum >= 0 ) {
64  evt->set_event_number(ievnum);
65  m_internal_event_number = ievnum;
66  }
67  else {
68  evt->set_event_number(m_internal_event_number);
69  ++m_internal_event_number;
70  }
71 
72  // Set units to be GeV and mm, to agree with Pythia ones.
73  evt->set_units(Units::GEV,Units::MM);
74 
75  // 1a. If there is a HIInfo object fill info from that.
76  if ( pyinfo && pyinfo->hiInfo ) {
77  auto ion = make_shared<HepMC3::GenHeavyIon>();
78  ion->Ncoll_hard = pyinfo->hiInfo->nCollND();
79  ion->Ncoll = pyinfo->hiInfo->nCollTot();
80  ion->Npart_proj = pyinfo->hiInfo->nAbsProj() +
81  pyinfo->hiInfo->nDiffProj();
82  ion->Npart_targ = pyinfo->hiInfo->nAbsTarg() +
83  pyinfo->hiInfo->nDiffTarg();
84  ion->impact_parameter = pyinfo->hiInfo->b();
85  evt->set_heavy_ion(ion);
86  }
87 
88  // 2. Fill particle information.
89  std::vector<GenParticlePtr> hepevt_particles;
90  hepevt_particles.reserve( pyev.size() );
91  for(int i = 0; i < pyev.size(); ++i) {
92  hepevt_particles.push_back( std::make_shared<GenParticle>(
93  FourVector( pyev[i].px(), pyev[i].py(), pyev[i].pz(), pyev[i].e() ),
94  pyev[i].id(), pyev[i].statusHepMC() ) );
95  hepevt_particles[i]->set_generated_mass( pyev[i].m() );
96  }
97 
98  // 3. Fill vertex information.
99  std::vector<GenVertexPtr> vertex_cache;
100  vector<GenParticlePtr> beam_particles;
101  for (int i = 1; i < pyev.size(); ++i) {
102  vector<int> mothers = pyev[i].motherList();
103  sort(mothers.begin(),mothers.end());
104  for (;;) {
105  if (!mothers.empty() && mothers.front() == 0)
106  mothers.erase(mothers.begin());
107  else break;
108  }
109  if (mothers.size()) {
110  GenVertexPtr prod_vtx = hepevt_particles[mothers[0]]->end_vertex();
111  if (!prod_vtx) {
112  prod_vtx = make_shared<GenVertex>();
113  vertex_cache.push_back(prod_vtx);
114  for (unsigned int j = 0; j < mothers.size(); ++j)
115  prod_vtx->add_particle_in( hepevt_particles[mothers[j]] );
116  }
117  FourVector prod_pos( pyev[i].xProd(), pyev[i].yProd(),pyev[i].zProd(),
118  pyev[i].tProd() );
119 
120  // Update vertex position if necessary.
121  if (!prod_pos.is_zero() && prod_vtx->position().is_zero())
122  prod_vtx->set_position( prod_pos );
123  prod_vtx->add_particle_out( hepevt_particles[i] );
124  } else beam_particles.push_back(hepevt_particles[i]);
125  }
126 
127  // Reserve memory for the event.
128  evt->reserve( hepevt_particles.size(), vertex_cache.size() );
129 
130  // Add particles and vertices in topological order.
131  evt->add_tree( beam_particles );
132 
133  // Attributes should be set after adding the particles to event.
134  for (int i = 0; i < pyev.size(); ++i) {
135  /* TODO: Set polarization */
136  // Colour flow uses index 1 and 2.
137  int colType = pyev[i].colType();
138  if (colType == -1 ||colType == 1 || colType == 2) {
139  int flow1 = 0, flow2 = 0;
140  if (colType == 1 || colType == 2) flow1 = pyev[i].col();
141  if (colType == -1 || colType == 2) flow2 = pyev[i].acol();
142  hepevt_particles[i]->add_attribute("flow1",
143  make_shared<IntAttribute>(flow1));
144  hepevt_particles[i]->add_attribute("flow2",
145  make_shared<IntAttribute>(flow2));
146  }
147  }
148 
149  // If hadronization switched on then no final coloured particles.
150  bool doHadr = (pyset == 0) ? m_free_parton_warnings
151  : pyset->flag("HadronLevel:all") && pyset->flag("HadronLevel:Hadronize");
152 
153  // 4. Check for particles which come from nowhere, i.e. are without
154  // mothers or daughters. These need to be attached to a vertex, or else
155  // they will never become part of the event.
156  for (int i = 1; i < pyev.size(); ++i) {
157 
158  // Check for particles not added to the event.
159  // NOTE: We have to check if this step makes any sense in
160  // the HepMC event standard.
161  if ( hepevt_particles[i] == nullptr ||
162  !hepevt_particles[i]->in_event()) {
163  warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
164  "found orphan particle", "i = " + Pythia8::toString(i));
165  GenVertexPtr prod_vtx = make_shared<GenVertex>();
166  prod_vtx->add_particle_out( hepevt_particles[i] );
167  evt->add_vertex(prod_vtx);
168  }
169 
170  // Also check for free partons (= gluons and quarks; not diquarks?).
171  if ( doHadr && m_free_parton_warnings ) {
172  if ( hepevt_particles[i]->pid() == 21
173  && hepevt_particles[i]->end_vertex() == nullptr ) {
174  warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
175  "found gluon without end vertex", "i = " + Pythia8::toString(i));
176  if ( m_crash_on_problem ) exit(1);
177  }
178  if ( abs(hepevt_particles[i]->pid()) <= 6
179  && hepevt_particles[i]->end_vertex() == nullptr ) {
180  warning(pyinfo, "Pythia8ToHepMC::fill_next_event",
181  "found quark without end vertex", "i = " + Pythia8::toString(i));
182  if ( m_crash_on_problem ) exit(1);
183  }
184  }
185  }
186 
187  // 5. Store PDF, weight, cross section and other event information.
188  // Flavours of incoming partons.
189  if (m_store_pdf && pyinfo != 0) {
190  int id1pdf = pyinfo->id1pdf();
191  int id2pdf = pyinfo->id2pdf();
192  if ( m_convert_gluon_to_0 ) {
193  if (id1pdf == 21) id1pdf = 0;
194  if (id2pdf == 21) id2pdf = 0;
195  }
196 
197  // Store PDF information.
198  GenPdfInfoPtr pdfinfo = make_shared<GenPdfInfo>();
199  pdfinfo->set(id1pdf, id2pdf, pyinfo->x1pdf(), pyinfo->x2pdf(),
200  pyinfo->QFac(), pyinfo->pdf1(), pyinfo->pdf2() );
201  evt->set_pdf_info( pdfinfo );
202  }
203 
204  // Store process code, scale, alpha_em, alpha_s.
205  if (m_store_proc && pyinfo != 0) {
206  evt->add_attribute("signal_process_id",
207  std::make_shared<IntAttribute>( pyinfo->code()));
208  evt->add_attribute("mpi",
209  std::make_shared<IntAttribute>( pyinfo->nMPI()));
210  evt->add_attribute("event_scale",
211  std::make_shared<DoubleAttribute>(pyinfo->QRen()));
212  evt->add_attribute("alphaQCD",
213  std::make_shared<DoubleAttribute>(pyinfo->alphaS()));
214  evt->add_attribute("alphaQED",
215  std::make_shared<DoubleAttribute>(pyinfo->alphaEM()));
216  }
217 
218  // Store event weights.
219  if (m_store_weights && pyinfo != 0) {
220  evt->weights().clear();
221  for (int iWeight = 0; iWeight < pyinfo->numberOfWeights(); ++iWeight)
222  evt->weights().push_back(pyinfo->weightValueByIndex(iWeight));
223  }
224 
225  // Store cross-section information in pb.
226  if (m_store_xsec && pyinfo != 0) {
227  // First set atribute to event, such that
228  // GenCrossSection::set_cross_section knows how many weights the
229  // event has and sets the number of cross sections accordingly.
230  GenCrossSectionPtr xsec = make_shared<GenCrossSection>();
231  evt->set_cross_section(xsec);
232  xsec->set_cross_section( pyinfo->sigmaGen() * 1e9,
233  pyinfo->sigmaErr() * 1e9);
234  // If multiweights with possibly different xsec, overwrite central value
235  vector<double> xsecVec = pyinfo->weightContainerPtr->getTotalXsec();
236  if (xsecVec.size() > 0) {
237  for (unsigned int iXsec = 0; iXsec < xsecVec.size(); ++iXsec) {
238  xsec->set_xsec(iXsec, xsecVec[iXsec]*1e9);
239  }
240  }
241  }
242 
243  // Done.
244  return true;
245  }
246 
247  // Read out values for some switches.
248  bool print_inconsistency() const { return m_print_inconsistency; }
249  bool free_parton_warnings() const { return m_free_parton_warnings; }
250  bool crash_on_problem() const { return m_crash_on_problem; }
251  bool convert_gluon_to_0() const { return m_convert_gluon_to_0; }
252  bool store_pdf() const { return m_store_pdf; }
253  bool store_proc() const { return m_store_proc; }
254  bool store_xsec() const { return m_store_xsec; }
255  bool store_weights() const { return m_store_weights; }
256 
257  // Set values for some switches.
258  void set_print_inconsistency(bool b = true) { m_print_inconsistency = b; }
259  void set_free_parton_warnings(bool b = true) { m_free_parton_warnings = b; }
260  void set_crash_on_problem(bool b = false) { m_crash_on_problem = b; }
261  void set_convert_gluon_to_0(bool b = false) { m_convert_gluon_to_0 = b; }
262  void set_store_pdf(bool b = true) { m_store_pdf = b; }
263  void set_store_proc(bool b = true) { m_store_proc = b; }
264  void set_store_xsec(bool b = true) { m_store_xsec = b; }
265  void set_store_weights(bool b = true) { m_store_weights = b; }
266 
267 private:
268 
269  // Try to send warning message to the logger if present, otherwise
270  // send it to cout if print_inconsistency().
271  bool warning(const Pythia8::Info * pyinfo, string loc,
272  string message, string extraInfo = "") {
273  if ( pyinfo )
274  pyinfo->loggerPtr->warningMsg(loc, message, extraInfo);
275  else if ( print_inconsistency() )
276  cout << "Warning in " << loc << ": " << message << extraInfo << endl;
277  return false;
278  }
279 
280  // Following methods are not implemented for this class.
281  virtual bool fill_next_event( GenEvent* ) { return 0; }
282  virtual void write_event( const GenEvent* ) {}
283 
284  // Use of copy constructor is not allowed.
285  Pythia8ToHepMC3( const Pythia8ToHepMC3& ) {}
286 
287  // Data members.
288  int m_internal_event_number;
289  bool m_print_inconsistency, m_free_parton_warnings, m_crash_on_problem,
290  m_convert_gluon_to_0, m_store_pdf, m_store_proc, m_store_xsec,
291  m_store_weights;
292 
293 };
294 
295 //==========================================================================
296 
297 } // end namespace HepMC3
298 
299 namespace Pythia8 {
300 
301 //==========================================================================
302 // This a wrapper around HepMC::Pythia8ToHepMC in the Pythia8
303 // namespace that simplify the most common use cases. It stores the
304 // current GenEvent and output stream internally to avoid cluttering
305 // of user code. This class is also defined in HepMC2.h with the same
306 // signatures, and the user can therefore switch between HepMC version
307 // 2 and 3, by simply changing the include file.
308 class Pythia8ToHepMC : public HepMC3::Pythia8ToHepMC3 {
309 
310 public:
311 
312  // We can either have standard ascii output version 2 or three or
313  // none at all.
314  enum OutputType { none, ascii2, ascii3 };
315 
316  // Typedef for the version 3 specific classes used.
317  typedef HepMC3::GenEvent GenEvent;
318  typedef shared_ptr<GenEvent> EventPtr;
319  typedef HepMC3::Writer Writer;
320  typedef shared_ptr<Writer> WriterPtr;
321 
322  // The empty constructor does not creat an aoutput stream.
323  Pythia8ToHepMC() : runinfo(make_shared<HepMC3::GenRunInfo>()) {}
324 
325  // Construct an object with an internal output stream.
326  Pythia8ToHepMC(string filename, OutputType ft = ascii3)
327  : runinfo(make_shared<HepMC3::GenRunInfo>()) {
328  setNewFile(filename, ft);
329  }
330 
331  // Open a new external output stream.
332  bool setNewFile(string filename, OutputType ft = ascii3) {
333  switch ( ft ) {
334  case ascii3:
335  writerPtr = make_shared<HepMC3::WriterAscii>(filename);
336  break;
337  case ascii2:
338  writerPtr = make_shared<HepMC3::WriterAsciiHepMC2>(filename);
339  break;
340  case none:
341  break;
342  }
343  return writerPtr != nullptr;
344  }
345 
346  // Create a new GenEvent object and fill it with information from
347  // the given Pythia object.
348  bool fillNextEvent(Pythia & pythia) {
349  geneve = make_shared<HepMC3::GenEvent>(runinfo);
350  if (runinfo->weight_names().size() == 0)
351  setWeightNames(pythia.info.weightNameVector());
352  return fill_next_event(pythia, *geneve);
353  }
354 
355  // Write out the current GenEvent to the internal stream.
356  void writeEvent() {
357  writerPtr->write_event(*geneve);
358  }
359 
360  // Create a new GenEvent object and fill it with information from
361  // the given Pythia object and write it out directly to the
362  // internal stream.
363  bool writeNextEvent(Pythia & pythia) {
364  if ( !fillNextEvent(pythia) ) return false;
365  writeEvent();
366  return !writerPtr->failed();
367  }
368 
369  // Get a reference to the current GenEvent.
370  GenEvent & event() {
371  return *geneve;
372  }
373 
374  // Get a pointer to the current GenEvent.
375  EventPtr getEventPtr() {
376  return geneve;
377  }
378 
379  // Get a reference to the internal stream.
380  Writer & output() {
381  return *writerPtr;
382  }
383 
384  // Get a pointer to the internal stream.
385  WriterPtr outputPtr() {
386  return writerPtr;
387  }
388 
389  // Set cross section information in the current GenEvent.
390  void setXSec(double xsec, double xsecerr) {
391  auto xsecptr = geneve->cross_section();
392  if ( !xsecptr ) {
393  xsecptr = make_shared<HepMC3::GenCrossSection>();
394  geneve->set_cross_section(xsecptr);
395  }
396  xsecptr->set_cross_section(xsec, xsecerr);
397  }
398 
399  // Update all weights in the current GenEvent.
400  void setWeights(const vector<double> & wv) {
401  geneve->weights() = wv;
402  }
403 
404  // Set all weight names in the current run.
405  void setWeightNames(const vector<string> &wnv) {
406  runinfo->set_weight_names(wnv);
407  }
408 
409  // Update the PDF information in the current GenEvent
410  void setPdfInfo(int id1, int id2, double x1, double x2,
411  double scale, double xf1, double xf2,
412  int pdf1 = 0, int pdf2 = 0) {
413  auto pdf = make_shared<HepMC3::GenPdfInfo>();
414  pdf->set(id1, id2, x1, x2, scale, xf1, xf2, pdf1, pdf2);
415  geneve->set_pdf_info(pdf);
416  }
417 
418  // Add an additional attribute derived from HepMC3::Attribute
419  // to the current event.
420  template<class T>
421  void addAtribute(const string& name, T& attribute) {
422  shared_ptr<HepMC3::Attribute> att = make_shared<T>(attribute);
423  geneve->add_attribute(name, att);
424  }
425 
426  // Add an attribute of double type.
427  template<class T=double>
428  void addAttribute(const string& name, double& attribute) {
429  auto dAtt = HepMC3::DoubleAttribute(attribute);
430  shared_ptr<HepMC3::Attribute> att =
431  make_shared<HepMC3::DoubleAttribute>(dAtt);
432  geneve->add_attribute(name, att);
433  }
434 
435  // Remove an attribute from the current event.
436  void removeAttribute(const string& name) {
437  geneve->remove_attribute(name);
438  }
439 
440 private:
441 
442  // The current GenEvent
443  EventPtr geneve = nullptr;
444 
445  // The output stream.
446  WriterPtr writerPtr = nullptr;
447 
448  // The current run info.
449  shared_ptr<HepMC3::GenRunInfo> runinfo;
450 
451 };
452 
453 }
454 
455 #endif // end Pythia8_HepMC3_H
void addAttribute(const string &name, double &attribute)
Add an attribute of double type.
Definition: HepMC3.h:428
Definition: HepMC3.h:29
void setWeightNames(const vector< string > &wnv)
Set all weight names in the current run.
Definition: HepMC3.h:405
GenEvent & event()
Get a reference to the current GenEvent.
Definition: HepMC3.h:370
Definition: Info.h:45
The Event class holds all info on the generated event.
Definition: Event.h:408
void setPdfInfo(int id1, int id2, double x1, double x2, double scale, double xf1, double xf2, int pdf1=0, int pdf2=0)
Update the PDF information in the current GenEvent.
Definition: HepMC3.h:410
bool fill_next_event(Pythia8::Event &pyev, GenEvent *evt, int ievnum=-1, const Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0)
Definition: HepMC3.h:55
bool setNewFile(string filename, OutputType ft=ascii3)
Open a new external output stream.
Definition: HepMC3.h:332
Logger * loggerPtr
Pointer to the logger.
Definition: Info.h:86
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:374
void removeAttribute(const string &name)
Remove an attribute from the current event.
Definition: HepMC3.h:436
bool writeNextEvent(Pythia &pythia)
Definition: HepMC3.h:363
Definition: HepMC3.h:33
Event event
The event record for the complete event history.
Definition: Pythia.h:365
void setXSec(double xsec, double xsecerr)
Set cross section information in the current GenEvent.
Definition: HepMC3.h:390
Pythia8ToHepMC3()
Constructor and destructor.
Definition: HepMC3.h:38
void warningMsg(string loc, string message, string extraInfo="", bool showAlways=false)
Warnings indicate that there might be an issue, but the run will continue.
Definition: Logger.cc:32
void writeEvent()
Write out the current GenEvent to the internal stream.
Definition: HepMC3.h:356
Pythia8ToHepMC(string filename, OutputType ft=ascii3)
Construct an object with an internal output stream.
Definition: HepMC3.h:326
EventPtr getEventPtr()
Get a pointer to the current GenEvent.
Definition: HepMC3.h:375
void set_print_inconsistency(bool b=true)
Set values for some switches.
Definition: HepMC3.h:258
Pythia8ToHepMC()
The empty constructor does not creat an aoutput stream.
Definition: HepMC3.h:323
void addAtribute(const string &name, T &attribute)
Definition: HepMC3.h:421
string toString(bool val)
Convert a boolean to a string.
Definition: PythiaStdlib.h:217
int size() const
Event record size.
Definition: Event.h:459
HepMC3::GenEvent GenEvent
Typedef for the version 3 specific classes used.
Definition: HepMC3.h:317
WriterPtr outputPtr()
Get a pointer to the internal stream.
Definition: HepMC3.h:385
void setWeights(const vector< double > &wv)
Update all weights in the current GenEvent.
Definition: HepMC3.h:400
bool fill_next_event(Pythia8::Pythia &pythia, GenEvent *evt, int ievnum=-1)
The recommended method to convert Pythia events into HepMC3 ones.
Definition: HepMC3.h:45
bool fill_next_event(Pythia8::Event &pyev, GenEvent &evt, int ievnum=-1, const Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0)
Alternative method to convert Pythia events into HepMC3 ones.
Definition: HepMC3.h:52
OutputType
We can either have standard ascii output or none at all.
Definition: HepMC2.h:414
bool print_inconsistency() const
Read out values for some switches.
Definition: HepMC3.h:248
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:71
bool fillNextEvent(Pythia &pythia)
Definition: HepMC3.h:348
const Info & info
Public information and statistic on the generation.
Definition: Pythia.h:368
Writer & output()
Get a reference to the internal stream.
Definition: HepMC3.h:380
Definition: Settings.h:196