11 #ifndef Pythia8_HepMC2_H
12 #define Pythia8_HepMC2_H
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"
32 virtual const char* what()
const throw() {
return "Pythia8ToHepMCException";}
49 ss <<
"Bad end vertex at " << i <<
" for flavour " << pdg_idIn;
55 virtual const char* what()
const throw() {
return msg.c_str();}
58 int index()
const {
return iSave;}
59 int pdg_id()
const {
return idSave;}
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() {;}
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);}
91 GenParticle* rootParticle = 0,
int iBarcode = -1);
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;}
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;}
112 virtual bool fill_next_event( GenEvent* ) {
return 0; }
113 virtual void write_event(
const GenEvent* ) {;}
116 Pythia8ToHepMC(
const Pythia8ToHepMC& ) : IO_BaseClass() {;}
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;
131 inline bool Pythia8ToHepMC::fill_next_event(
Pythia8::Event& pyev,
133 bool append, GenParticle* rootParticle,
int iBarcode) {
137 std::cout <<
" Pythia8ToHepMC::fill_next_event error: passed null event."
145 evt->set_event_number(ievnum);
146 m_internal_event_number = ievnum;
148 evt->set_event_number(m_internal_event_number);
149 ++m_internal_event_number;
154 double momFac = HepMC::Units::conversion_factor(HepMC::Units::GEV,
155 evt->momentum_unit());
156 double lenFac = HepMC::Units::conversion_factor(HepMC::Units::MM,
164 std::cout <<
" Pythia8ToHepMC::fill_next_event error: passed null "
165 <<
"root particle in append mode." << std::endl;
169 newBarcode = (iBarcode > -1) ? iBarcode : evt->particles_size();
171 GenVertex* prod_vtx0 =
new GenVertex();
172 prod_vtx0->add_particle_in( rootParticle );
173 evt->add_vertex( prod_vtx0 );
178 std::vector<GenParticle*> hepevt_particles( pyev.size() );
179 for (
int i = iStart; i < pyev.size(); ++i) {
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() );
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());
201 evt->set_beam_particles( hepevt_particles[1], hepevt_particles[2] );
206 for (
int i = iStart; i < pyev.size(); ++i) {
207 GenParticle* p = hepevt_particles[i];
210 std::vector<int> mothers = pyev[i].motherList();
211 unsigned int imother = 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;
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 );
233 if ( prod_vtx && prod_vtx->position() == FourVector() )
234 prod_vtx->set_position( prod_pos );
239 if ( !mothers.empty() ) mother = mothers[imother];
240 while ( prod_vtx && mother > 0 ) {
243 GenParticle* ppp = (append && mother == 1) ? rootParticle
244 : hepevt_particles[mother];
245 if ( !ppp->end_vertex() ) {
246 prod_vtx->add_particle_in( ppp );
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;
263 mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
268 bool doHadr = (pyset == 0) ? m_free_parton_exception
269 : pyset->flag(
"HadronLevel:all") && pyset->flag(
"HadronLevel:Hadronize");
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 );
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);
294 if (append)
return true;
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;
307 evt->set_pdf_info( PdfInfo( id1pdf, id2pdf, pyinfo->x1pdf(),
308 pyinfo->x2pdf(), pyinfo->QFac(), pyinfo->pdf1(), pyinfo->pdf2() ) );
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() );
321 if (m_store_xsec && pyinfo != 0) {
324 pyinfo->sigmaErr() * 1e9);
325 evt->set_cross_section(xsec);
326 evt->weights().push_back( pyinfo->weight() );
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