CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HepMCA2.h
Go to the documentation of this file.
1 // HepMC2.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Author: Mikhail Kirsanov, Mikhail.Kirsanov@cern.ch
7 // Header file and function definitions for the Pythia8ToHepMC class,
8 // which converts a PYTHIA event record to the standard HepMC format.
9 
10 #ifndef Pythia8_HepMCA2_H
11 #define Pythia8_HepMCA2_H
12 
13 #include <vector>
14 #include "HepMC/IO_BaseClass.h"
15 #include "HepMC/IO_GenEvent.h"
16 #include "HepMC/GenEvent.h"
17 #include "HepMC/Units.h"
18 #include "Pythia8/Pythia.h"
19 
20 namespace HepMC {
21 
22 //==========================================================================
23 
24 // The Pythia8ToHepMCA class.
25 
26 class Pythia8ToHepMCA : public IO_BaseClass {
27 
28 public:
29 
30  // Constructor and destructor.
34  virtual ~Pythia8ToHepMCA() {;}
35 
36  // Alternative method to convert Pythia events into HepMC ones.
37  bool append_event( Pythia8::Event& pyev, GenEvent* evt, GenParticle* rootpart,
38  int ibarcode = -1, Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0);
39 
40  // Read out values for some switches.
43  bool crash_on_problem() const {return m_crash_on_problem;}
45 
46  // Set values for some switches.
49  void set_crash_on_problem(bool b = false) {m_crash_on_problem = b;}
51 
52 private:
53 
54  // Following methods are not implemented for this class.
55  virtual bool fill_next_event( GenEvent* ) { return 0; }
56  virtual void write_event( const GenEvent* ) {;}
57 
58  // Use of copy constructor is not allowed.
59  Pythia8ToHepMCA( const Pythia8ToHepMCA& ) : IO_BaseClass() {;}
60 
61  // Data members.
65 
66 };
67 
68 //==========================================================================
69 
70 // Main method to append PYTHIA event to HepMC event.
71 // Read one event from Pythia8, append GenEvent
72 // and return T/F = success/failure.
73 
74 inline bool Pythia8ToHepMCA::append_event( Pythia8::Event& pyev, GenEvent* evt, GenParticle* rootpart,
75  int ibarcode, Pythia8::Info* pyinfo, Pythia8::Settings* pyset) {
76 
77  // 1. Error if no event passed.
78  if (!evt) {
79  std::cerr << "Pythia8ToHepMCA::fill_next_event error - passed null event."
80  << std::endl;
81  return 0;
82  }
83 
84  // Conversion factors from Pythia units GeV and mm to HepMC ones.
85  double momFac = HepMC::Units::conversion_factor(HepMC::Units::GEV,
86  evt->momentum_unit());
87  double lenFac = HepMC::Units::conversion_factor(HepMC::Units::MM,
88  evt->length_unit());
89 
90  int NewBarcode = evt->particles_size();
91  if (ibarcode > -1) NewBarcode = ibarcode;
92 
93  GenVertex* prod_vtx0 = new GenVertex();
94  prod_vtx0->add_particle_in( rootpart );
95  evt->add_vertex( prod_vtx0 );
96 
97  // 2. Create a particle instance for each entry and fill a map, and
98  // a vector which maps from the particle index to the GenParticle address.
99  std::vector<GenParticle*> hepevt_particles( pyev.size() );
100  for (int i = 2; i < pyev.size(); ++i) {
101 
102  // Fill the particle.
103 
104  hepevt_particles[i] = new GenParticle(
105  FourVector( momFac * pyev[i].px(), momFac * pyev[i].py(),
106  momFac * pyev[i].pz(), momFac * pyev[i].e() ),
107  pyev[i].id(), pyev[i].statusHepMC() );
108  if (ibarcode !=0) NewBarcode++;
109  hepevt_particles[i]->suggest_barcode(NewBarcode);
110  hepevt_particles[i]->set_generated_mass( momFac * pyev[i].m() );
111 
112  // Colour flow uses index 1 and 2.
113  int colType = pyev[i].colType();
114  if (colType == 1 || colType == 2)
115  hepevt_particles[i]->set_flow(1, pyev[i].col());
116  if (colType == -1 || colType == 2)
117  hepevt_particles[i]->set_flow(2, pyev[i].acol());
118  }
119 
120  // 3. Loop over particles AGAIN, this time creating vertices.
121  // We build the production vertex for each entry in hepevt.
122  // The HEPEVT pointers are bi-directional, so gives decay vertices as well.
123  for (int i = 2; i < pyev.size(); ++i) {
124  GenParticle* p = hepevt_particles[i];
125 
126  // 3a. Search to see if a production vertex already exists.
127  std::vector<int> mothers = pyev[i].motherList();
128 
129  unsigned int imother = 0;
130  int mother = -1; // note that in Pythia8 there is a particle number 0!
131  if ( !mothers.empty() ) mother = mothers[imother];
132  GenVertex* prod_vtx = p->production_vertex();
133  while ( !prod_vtx && mother > 0 ) {
134  if(mother == 1) {
135  prod_vtx = rootpart->end_vertex();
136  } else {
137  prod_vtx = hepevt_particles[mother]->end_vertex();
138  }
139  if ( prod_vtx ) prod_vtx->add_particle_out( p );
140  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
141  }
142 
143  // 3b. If no suitable production vertex exists - and the particle has
144  // at least one mother or position information to store - make one.
145  FourVector prod_pos( lenFac * pyev[i].xProd(), lenFac * pyev[i].yProd(),
146  lenFac * pyev[i].zProd(), lenFac * pyev[i].tProd() );
147  if ( !prod_vtx && ( mothers.size() > 0 || prod_pos != FourVector() ) ) {
148  prod_vtx = new GenVertex();
149  prod_vtx->add_particle_out( p );
150  evt->add_vertex( prod_vtx );
151  }
152 
153  // 3c. If prod_vtx doesn't already have position specified, fill it.
154  if ( prod_vtx && prod_vtx->position() == FourVector() )
155  prod_vtx->set_position( prod_pos );
156 
157  // 3d. loop over mothers to make sure their end_vertices are consistent.
158  imother = 0;
159  mother = -1;
160  if ( !mothers.empty() ) mother = mothers[imother];
161  while ( prod_vtx && mother > 0 ) {
162 
163  // If end vertex of the mother isn't specified, do it now.
164  GenParticle* ppp;
165  if(mother == 1) ppp = rootpart;
166  if(mother > 1) ppp = hepevt_particles[mother];
167 
168  if ( !ppp->end_vertex() ) {
169  prod_vtx->add_particle_in( ppp );
170 
171  // Problem scenario: the mother already has a decay vertex which
172  // differs from the daughter's production vertex. This means there is
173  // internal inconsistency in the HEPEVT event record. Print an error.
174  // Note: we could provide a fix by joining the two vertices with a
175  // dummy particle if the problem arises often.
176  } else if (ppp->end_vertex() != prod_vtx ) {
178  << "HepMC::Pythia8ToHepMC: inconsistent mother/daugher "
179  << "information in Pythia8 event " << std::endl
180  << "i = " << i << " mother = " << mother
181  << "\n This warning can be turned off with the "
182  << "Pythia8ToHepMC::print_inconsistency switch." << std::endl;
183  }
184 
185  // End of vertex-setting loops.
186  mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
187  }
188  }
189 
190  // If hadronization switched on then no final coloured particles.
191  bool doHadr = (pyset == 0) ? m_free_parton_warnings
192  : pyset->flag("HadronLevel:all") && pyset->flag("HadronLevel:Hadronize");
193 
194  // 4. Check for particles which come from nowhere, i.e. are without
195  // mothers or daughters. These need to be attached to a vertex, or else
196  // they will never become part of the event.
197  for (int i = 2; i < pyev.size(); ++i) {
198  if ( !hepevt_particles[i]->end_vertex() &&
199  !hepevt_particles[i]->production_vertex() ) {
200  std::cerr << "hanging particle " << i << std::endl;
201  GenVertex* prod_vtx = new GenVertex();
202  prod_vtx->add_particle_out( hepevt_particles[i] );
203  evt->add_vertex( prod_vtx );
204  }
205 
206  // Also check for free partons (= gluons and quarks; not diquarks?).
207  if ( doHadr && m_free_parton_warnings ) {
208  if ( hepevt_particles[i]->pdg_id() == 21 &&
209  !hepevt_particles[i]->end_vertex() ) {
210  std::cerr << "gluon without end vertex " << i << std::endl;
211  if ( m_crash_on_problem ) exit(1);
212  }
213  if ( abs(hepevt_particles[i]->pdg_id()) <= 6 &&
214  !hepevt_particles[i]->end_vertex() ) {
215  std::cerr << "quark without end vertex " << i << std::endl;
216  if ( m_crash_on_problem ) exit(1);
217  }
218  }
219  }
220 
221  // Done.
222  return true;
223 
224 }
225 
226 //==========================================================================
227 
228 } // end namespace HepMC
229 
230 #endif // end Pythia8_HepMCA2_H
int i
Definition: DBlmapReader.cc:9
bool m_free_parton_warnings
Definition: HepMCA2.h:63
virtual ~Pythia8ToHepMCA()
Definition: HepMCA2.h:34
void set_crash_on_problem(bool b=false)
Definition: HepMCA2.h:49
void set_print_inconsistency(bool b=true)
Definition: HepMCA2.h:47
bool crash_on_problem() const
Definition: HepMCA2.h:43
bool print_inconsistency() const
Definition: HepMCA2.h:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool m_print_inconsistency
Definition: HepMCA2.h:63
virtual void write_event(const GenEvent *)
Definition: HepMCA2.h:56
double b
Definition: hdecay.h:120
bool free_parton_warnings() const
Definition: HepMCA2.h:42
void set_convert_gluon_to_0(bool b=false)
Definition: HepMCA2.h:50
bool convert_gluon_to_0() const
Definition: HepMCA2.h:44
void set_free_parton_warnings(bool b=true)
Definition: HepMCA2.h:48
volatile std::atomic< bool > shutdown_flag false
Pythia8ToHepMCA(const Pythia8ToHepMCA &)
Definition: HepMCA2.h:59
int col
Definition: cuy.py:1008
bool append_event(Pythia8::Event &pyev, GenEvent *evt, GenParticle *rootpart, int ibarcode=-1, Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0)
Definition: HepMCA2.h:74
virtual bool fill_next_event(GenEvent *)
Definition: HepMCA2.h:55