StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HepMC2.h
1 // HepMC2.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 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: Mikhail Kirsanov, Mikhail.Kirsanov@cern.ch
7 // Eexception classes provided by James Monk, with minor changes.
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_HepMC2_H
12 #define Pythia8_HepMC2_H
13 
14 #include <exception>
15 #include <sstream>
16 #include <vector>
17 #include "HepMC/IO_BaseClass.h"
18 #include "HepMC/IO_GenEvent.h"
19 #include "HepMC/GenEvent.h"
20 #include "HepMC/Units.h"
21 #include "Pythia8/Pythia.h"
22 
23 namespace HepMC {
24 
25 //==========================================================================
26 
27 // Base exception for all exceptions that Pythia8ToHepMC might throw.
28 
29 class Pythia8ToHepMCException : public std::exception {
30 
31 public:
32  virtual const char* what() const throw() { return "Pythia8ToHepMCException";}
33 
34 };
35 
36 //--------------------------------------------------------------------------
37 
38 // Exception thrown when an undecayed parton is written into the record.
39 
41 
42 public:
43 
44  // Constructor and destructor.
45  PartonEndVertexException(int i, int pdg_idIn) : Pythia8ToHepMCException() {
46  iSave = i;
47  idSave = pdg_idIn;
48  std::stringstream ss;
49  ss << "Bad end vertex at " << i << " for flavour " << pdg_idIn;
50  msg = ss.str();
51  }
52  virtual ~PartonEndVertexException() throw() {}
53 
54  // Throw exception.
55  virtual const char* what() const throw() {return msg.c_str();}
56 
57  // Return info on location and flavour of bad parton.
58  int index() const {return iSave;}
59  int pdg_id() const {return idSave;}
60 
61 protected:
62 
63  std::string msg;
64  int iSave, idSave;
65 };
66 
67 //==========================================================================
68 
69 // The Pythia8ToHepMC class.
70 
71 class Pythia8ToHepMC : public IO_BaseClass {
72 
73 public:
74 
75  // Constructor and destructor.
76  Pythia8ToHepMC() : m_internal_event_number(0), m_print_inconsistency(true),
77  m_free_parton_exception(true), m_convert_gluon_to_0(false),
78  m_store_pdf(true), m_store_proc(true), m_store_xsec(true) {;}
79  virtual ~Pythia8ToHepMC() {;}
80 
81  // The recommended method to convert Pythia events into HepMC ones.
82  bool fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt,
83  int ievnum = -1, bool append = false, GenParticle* rootParticle = 0,
84  int iBarcode = -1 ) {return fill_next_event( pythia.event, evt, ievnum,
85  &pythia.info, &pythia.settings, append, rootParticle, iBarcode);}
86 
87  // Alternative method to convert Pythia events into HepMC ones.
88  bool fill_next_event( Pythia8::Event& pyev, GenEvent* evt,
89  int ievnum = -1, Pythia8::Info* pyinfo = 0,
90  Pythia8::Settings* pyset = 0, bool append = false,
91  GenParticle* rootParticle = 0, int iBarcode = -1);
92 
93  // Read out values for some switches.
94  bool print_inconsistency() const {return m_print_inconsistency;}
95  bool free_parton_exception() const {return m_free_parton_exception;}
96  bool convert_gluon_to_0() const {return m_convert_gluon_to_0;}
97  bool store_pdf() const {return m_store_pdf;}
98  bool store_proc() const {return m_store_proc;}
99  bool store_xsec() const {return m_store_xsec;}
100 
101  // Set values for some switches.
102  void set_print_inconsistency(bool b = true) {m_print_inconsistency = b;}
103  void set_free_parton_exception(bool b = true) {m_free_parton_exception = b;}
104  void set_convert_gluon_to_0(bool b = false) {m_convert_gluon_to_0 = b;}
105  void set_store_pdf(bool b = true) {m_store_pdf = b;}
106  void set_store_proc(bool b = true) {m_store_proc = b;}
107  void set_store_xsec(bool b = true) {m_store_xsec = b;}
108 
109 private:
110 
111  // Following methods are not implemented for this class.
112  virtual bool fill_next_event( GenEvent* ) { return 0; }
113  virtual void write_event( const GenEvent* ) {;}
114 
115  // Use of copy constructor is not allowed.
116  Pythia8ToHepMC( const Pythia8ToHepMC& ) : IO_BaseClass() {;}
117 
118  // Data members.
119  int m_internal_event_number;
120  bool m_print_inconsistency, m_free_parton_exception, m_convert_gluon_to_0,
121  m_store_pdf, m_store_proc, m_store_xsec;
122 
123 };
124 
125 //==========================================================================
126 
127 // Main method for conversion from PYTHIA event to HepMC event.
128 // Read one event from Pythia8 and fill a new GenEvent, alternatively
129 // append to an existing GenEvent, and return T/F = success/failure.
130 
131 inline bool Pythia8ToHepMC::fill_next_event( Pythia8::Event& pyev,
132  GenEvent* evt, int ievnum, Pythia8::Info* pyinfo, Pythia8::Settings* pyset,
133  bool append, GenParticle* rootParticle, int iBarcode) {
134 
135  // 1. Error if no event passed.
136  if (!evt) {
137  std::cout << " Pythia8ToHepMC::fill_next_event error: passed null event."
138  << std::endl;
139  return 0;
140  }
141 
142  // Update event number counter.
143  if (!append) {
144  if (ievnum >= 0) {
145  evt->set_event_number(ievnum);
146  m_internal_event_number = ievnum;
147  } else {
148  evt->set_event_number(m_internal_event_number);
149  ++m_internal_event_number;
150  }
151  }
152 
153  // Conversion factors from Pythia units GeV and mm to HepMC ones.
154  double momFac = HepMC::Units::conversion_factor(HepMC::Units::GEV,
155  evt->momentum_unit());
156  double lenFac = HepMC::Units::conversion_factor(HepMC::Units::MM,
157  evt->length_unit());
158 
159  // Set up for alternative to append to an existing event.
160  int iStart = 1;
161  int newBarcode = 0;
162  if (append) {
163  if (!rootParticle) {
164  std::cout << " Pythia8ToHepMC::fill_next_event error: passed null "
165  << "root particle in append mode." << std::endl;
166  return 0;
167  }
168  iStart = 2;
169  newBarcode = (iBarcode > -1) ? iBarcode : evt->particles_size();
170  // New vertex associated with appended particles.
171  GenVertex* prod_vtx0 = new GenVertex();
172  prod_vtx0->add_particle_in( rootParticle );
173  evt->add_vertex( prod_vtx0 );
174  }
175 
176  // 2. Create a particle instance for each entry and fill a map, and
177  // a vector which maps from the particle index to the GenParticle address.
178  std::vector<GenParticle*> hepevt_particles( pyev.size() );
179  for (int i = iStart; i < pyev.size(); ++i) {
180 
181  // Fill the particle.
182  hepevt_particles[i] = new GenParticle(
183  FourVector( momFac * pyev[i].px(), momFac * pyev[i].py(),
184  momFac * pyev[i].pz(), momFac * pyev[i].e() ),
185  pyev[i].id(), pyev[i].statusHepMC() );
186  if (iBarcode != 0) ++newBarcode;
187  hepevt_particles[i]->suggest_barcode( (append) ? newBarcode : i );
188  hepevt_particles[i]->set_generated_mass( momFac * pyev[i].m() );
189 
190  // Colour flow uses index 1 and 2.
191  int colType = pyev[i].colType();
192  if (colType == 1 || colType == 2)
193  hepevt_particles[i]->set_flow(1, pyev[i].col());
194  if (colType == -1 || colType == 2)
195  hepevt_particles[i]->set_flow(2, pyev[i].acol());
196  }
197 
198  // Here we assume that the first two particles in the list
199  // are the incoming beam particles.
200  if (!append)
201  evt->set_beam_particles( hepevt_particles[1], hepevt_particles[2] );
202 
203  // 3. Loop over particles AGAIN, this time creating vertices.
204  // We build the production vertex for each entry in hepevt.
205  // The HEPEVT pointers are bi-directional, so gives decay vertices as well.
206  for (int i = iStart; i < pyev.size(); ++i) {
207  GenParticle* p = hepevt_particles[i];
208 
209  // 3a. Search to see if a production vertex already exists.
210  std::vector<int> mothers = pyev[i].motherList();
211  unsigned int imother = 0;
212  int mother = -1; // note that in Pythia8 there is a particle number 0!
213  if ( !mothers.empty() ) mother = mothers[imother];
214  GenVertex* prod_vtx = p->production_vertex();
215  while ( !prod_vtx && mother > 0 ) {
216  prod_vtx = (append && mother == 1) ? rootParticle->end_vertex()
217  : hepevt_particles[mother]->end_vertex();
218  if (prod_vtx) prod_vtx->add_particle_out( p );
219  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
220  }
221 
222  // 3b. If no suitable production vertex exists - and the particle has
223  // at least one mother or position information to store - make one.
224  FourVector prod_pos( lenFac * pyev[i].xProd(), lenFac * pyev[i].yProd(),
225  lenFac * pyev[i].zProd(), lenFac * pyev[i].tProd() );
226  if ( !prod_vtx && ( mothers.size() > 0 || prod_pos != FourVector() ) ) {
227  prod_vtx = new GenVertex();
228  prod_vtx->add_particle_out( p );
229  evt->add_vertex( prod_vtx );
230  }
231 
232  // 3c. If prod_vtx doesn't already have position specified, fill it.
233  if ( prod_vtx && prod_vtx->position() == FourVector() )
234  prod_vtx->set_position( prod_pos );
235 
236  // 3d. loop over mothers to make sure their end_vertices are consistent.
237  imother = 0;
238  mother = -1;
239  if ( !mothers.empty() ) mother = mothers[imother];
240  while ( prod_vtx && mother > 0 ) {
241 
242  // If end vertex of the mother isn't specified, do it now.
243  GenParticle* ppp = (append && mother == 1) ? rootParticle
244  : hepevt_particles[mother];
245  if ( !ppp->end_vertex() ) {
246  prod_vtx->add_particle_in( ppp );
247 
248  // Problem scenario: the mother already has a decay vertex which
249  // differs from the daughter's production vertex. This means there is
250  // internal inconsistency in the HEPEVT event record. Print an error.
251  // Note: we could provide a fix by joining the two vertices with a
252  // dummy particle if the problem arises often.
253  } else if (ppp->end_vertex() != prod_vtx ) {
254  if ( m_print_inconsistency ) std::cout
255  << " Pythia8ToHepMC::fill_next_event: inconsistent mother/daugher "
256  << "information in Pythia8 event " << std::endl
257  << "i = " << i << " mother = " << mother
258  << "\n This warning can be turned off with the "
259  << "Pythia8ToHepMC::print_inconsistency switch." << std::endl;
260  }
261 
262  // End of vertex-setting loops.
263  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
264  }
265  }
266 
267  // If hadronization switched on then no final coloured particles.
268  bool doHadr = (pyset == 0) ? m_free_parton_exception
269  : pyset->flag("HadronLevel:all") && pyset->flag("HadronLevel:Hadronize");
270 
271  // 4. Check for particles which come from nowhere, i.e. are without
272  // mothers or daughters. These need to be attached to a vertex, or else
273  // they will never become part of the event.
274  for (int i = iStart; i < pyev.size(); ++i) {
275  if ( !hepevt_particles[i]->end_vertex() &&
276  !hepevt_particles[i]->production_vertex() ) {
277  std::cout << " Pythia8ToHepMC::fill_next_event error: "
278  << "hanging particle " << i << std::endl;
279  GenVertex* prod_vtx = new GenVertex();
280  prod_vtx->add_particle_out( hepevt_particles[i] );
281  evt->add_vertex( prod_vtx );
282  }
283 
284  // Also check for free partons (= gluons and quarks; not diquarks?).
285  if ( doHadr && m_free_parton_exception ) {
286  int pdg_tmp = hepevt_particles[i]->pdg_id();
287  if ( (abs(pdg_tmp) <= 6 || pdg_tmp == 21)
288  && !hepevt_particles[i]->end_vertex() )
289  throw PartonEndVertexException(i, pdg_tmp);
290  }
291  }
292 
293  // Done if only appending to already existing event.
294  if (append) return true;
295 
296  // 5. Store PDF, weight, cross section and other event information.
297  // Flavours of incoming partons.
298  if (m_store_pdf && pyinfo != 0) {
299  int id1pdf = pyinfo->id1pdf();
300  int id2pdf = pyinfo->id2pdf();
301  if ( m_convert_gluon_to_0 ) {
302  if (id1pdf == 21) id1pdf = 0;
303  if (id2pdf == 21) id2pdf = 0;
304  }
305 
306  // Store PDF information.
307  evt->set_pdf_info( PdfInfo( id1pdf, id2pdf, pyinfo->x1pdf(),
308  pyinfo->x2pdf(), pyinfo->QFac(), pyinfo->pdf1(), pyinfo->pdf2() ) );
309  }
310 
311  // Store process code, scale, alpha_em, alpha_s.
312  if (m_store_proc && pyinfo != 0) {
313  evt->set_signal_process_id( pyinfo->code() );
314  evt->set_event_scale( pyinfo->QRen() );
315  if (evt->alphaQED() <= 0) evt->set_alphaQED( pyinfo->alphaEM() );
316  if (evt->alphaQCD() <= 0) evt->set_alphaQCD( pyinfo->alphaS() );
317  }
318 
319  // Store cross-section information in pb and event weight. The latter is
320  // usually dimensionless, but in units of pb for Les Houches strategies +-4.
321  if (m_store_xsec && pyinfo != 0) {
323  xsec.set_cross_section( pyinfo->sigmaGen() * 1e9,
324  pyinfo->sigmaErr() * 1e9);
325  evt->set_cross_section(xsec);
326  evt->weights().push_back( pyinfo->weight() );
327  }
328 
329  // Done for new event.
330  return true;
331 
332 }
333 
334 //==========================================================================
335 
336 } // end namespace HepMC
337 
338 #endif // end Pythia8_HepMC2_H
void set_cross_section(double xs, double xs_err)
Set cross section and error in pb.
The GenCrossSection class stores the generated cross section.
all input/output classes inherit from IO_BaseClass
Definition: IO_BaseClass.h:34