13 #ifndef Pythia8_HepMC2_H 14 #define Pythia8_HepMC2_H 16 #ifdef Pythia8_HepMC3_H 17 #error Cannot include HepMC2.h if HepMC3.h has already been included. 23 #include "HepMC/IO_BaseClass.h" 24 #include "HepMC/IO_GenEvent.h" 25 #include "HepMC/GenEvent.h" 26 #include "HepMC/Units.h" 27 #include "Pythia8/Pythia.h" 28 #include "Pythia8/HIInfo.h" 39 virtual const char* what()
const throw() {
return "Pythia8ToHepMCException";}
56 ss <<
"Bad end vertex at " << i <<
" for flavour " << pdg_idIn;
62 virtual const char*
what()
const throw() {
return msg.c_str();}
65 int index()
const {
return iSave;}
66 int pdg_id()
const {
return idSave;}
84 m_free_parton_exception(true), m_convert_gluon_to_0(false),
85 m_store_pdf(true), m_store_proc(true), m_store_xsec(true),
86 m_store_weights(true) {;}
91 int ievnum = -1,
bool append =
false, GenParticle* rootParticle = 0,
92 int iBarcode = -1 ) {
return fill_next_event( pythia.
event, &evt, ievnum,
93 &pythia.
info, &pythia.
settings, append, rootParticle, iBarcode);}
95 int ievnum = -1,
bool append =
false, GenParticle* rootParticle = 0,
96 int iBarcode = -1 ) {
return fill_next_event( pythia.
event, evt, ievnum,
97 &pythia.
info, &pythia.
settings, append, rootParticle, iBarcode);}
103 GenParticle* rootParticle = 0,
int iBarcode = -1) {
104 return fill_next_event(pyev, &evt, ievnum, pyinfo, pyset,
105 append, rootParticle, iBarcode); }
110 GenParticle* rootParticle = 0,
int iBarcode = -1);
114 bool free_parton_exception()
const {
return m_free_parton_exception;}
115 bool free_parton_warnings()
const {
return m_free_parton_exception;}
116 bool convert_gluon_to_0()
const {
return m_convert_gluon_to_0;}
117 bool store_pdf()
const {
return m_store_pdf;}
118 bool store_proc()
const {
return m_store_proc;}
119 bool store_xsec()
const {
return m_store_xsec;}
120 bool store_weights()
const {
return m_store_weights;}
124 void set_free_parton_exception(
bool b =
true) {m_free_parton_exception = b;}
125 void set_free_parton_warnings(
bool b =
true) {m_free_parton_exception = b;}
126 void set_convert_gluon_to_0(
bool b =
false) {m_convert_gluon_to_0 = b;}
127 void set_store_pdf(
bool b =
true) {m_store_pdf = b;}
128 void set_store_proc(
bool b =
true) {m_store_proc = b;}
129 void set_store_xsec(
bool b =
true) {m_store_xsec = b;}
130 void set_store_weights(
bool b =
true) {m_store_weights = b;}
135 virtual bool fill_next_event( GenEvent* ) {
return 0; }
136 virtual void write_event(
const GenEvent* ) {;}
142 int m_internal_event_number;
143 bool m_print_inconsistency, m_free_parton_exception, m_convert_gluon_to_0,
144 m_store_pdf, m_store_proc, m_store_xsec, m_store_weights;
161 std::cout <<
" Pythia8ToHepMC::fill_next_event error: passed null event." 169 evt->set_event_number(ievnum);
170 m_internal_event_number = ievnum;
172 evt->set_event_number(m_internal_event_number);
173 ++m_internal_event_number;
178 double momFac = HepMC::Units::conversion_factor(HepMC::Units::GEV,
179 evt->momentum_unit());
180 double lenFac = HepMC::Units::conversion_factor(HepMC::Units::MM,
188 std::cout <<
" Pythia8ToHepMC::fill_next_event error: passed null " 189 <<
"root particle in append mode." << std::endl;
193 newBarcode = (iBarcode > -1) ? iBarcode : evt->particles_size();
195 GenVertex* prod_vtx0 =
new GenVertex();
196 prod_vtx0->add_particle_in( rootParticle );
197 evt->add_vertex( prod_vtx0 );
201 if ( pyinfo && pyinfo->
hiInfo ) {
214 ion.set_impact_parameter(pyinfo->
hiInfo->
b());
215 evt->set_heavy_ion(ion);
220 std::vector<GenParticle*> hepevt_particles( pyev.
size() );
221 for (
int i = iStart; i < pyev.
size(); ++i) {
224 hepevt_particles[i] =
new GenParticle(
225 FourVector( momFac * pyev[i].px(), momFac * pyev[i].py(),
226 momFac * pyev[i].pz(), momFac * pyev[i].e() ),
227 pyev[i].
id(), pyev[i].statusHepMC() );
228 if (iBarcode != 0) ++newBarcode;
229 hepevt_particles[i]->suggest_barcode( (append) ? newBarcode : i );
230 hepevt_particles[i]->set_generated_mass( momFac * pyev[i].m() );
233 int colType = pyev[i].colType();
234 if (colType == 1 || colType == 2)
235 hepevt_particles[i]->set_flow(1, pyev[i].col());
236 if (colType == -1 || colType == 2)
237 hepevt_particles[i]->set_flow(2, pyev[i].acol());
243 evt->set_beam_particles( hepevt_particles[1], hepevt_particles[2] );
248 for (
int i = iStart; i < pyev.
size(); ++i) {
249 GenParticle* p = hepevt_particles[i];
252 std::vector<int> mothers = pyev[i].motherList();
253 unsigned int imother = 0;
255 if ( !mothers.empty() ) mother = mothers[imother];
256 GenVertex* prod_vtx = p->production_vertex();
257 while ( !prod_vtx && mother > 0 ) {
258 prod_vtx = (append && mother == 1) ? rootParticle->end_vertex()
259 : hepevt_particles[mother]->end_vertex();
260 if (prod_vtx) prod_vtx->add_particle_out( p );
261 mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
266 FourVector prod_pos( lenFac * pyev[i].xProd(), lenFac * pyev[i].yProd(),
267 lenFac * pyev[i].zProd(), lenFac * pyev[i].tProd() );
268 if ( !prod_vtx && ( mothers.size() > 0 || prod_pos != FourVector() ) ) {
269 prod_vtx =
new GenVertex();
270 prod_vtx->add_particle_out( p );
271 evt->add_vertex( prod_vtx );
275 if ( prod_vtx && prod_vtx->position() == FourVector() )
276 prod_vtx->set_position( prod_pos );
281 if ( !mothers.empty() ) mother = mothers[imother];
282 while ( prod_vtx && mother > 0 ) {
285 GenParticle* ppp = (append && mother == 1) ? rootParticle
286 : hepevt_particles[mother];
287 if ( !ppp->end_vertex() ) {
288 prod_vtx->add_particle_in( ppp );
289 }
else if (ppp->end_vertex() != prod_vtx ) {
295 if ( (pyev[i].statusAbs() == 43 || pyev[i].statusAbs() == 51
296 || pyev[i].statusAbs() == 53) && mother >= 1)
break;
301 if ( m_print_inconsistency ) std::cout
302 <<
" Pythia8ToHepMC::fill_next_event: inconsistent mother/daugher " 303 <<
"information in Pythia8 event " << std::endl
304 <<
"i = " << i <<
" mother = " << mother
305 <<
"\n This warning can be turned off with the " 306 <<
"Pythia8ToHepMC::print_inconsistency switch." << std::endl;
310 mother = ( ++imother < mothers.size() ) ? mothers[imother] : -1;
315 bool doHadr = (pyset == 0) ? m_free_parton_exception
316 : pyset->
flag(
"HadronLevel:all") && pyset->
flag(
"HadronLevel:Hadronize");
321 for (
int i = iStart; i < pyev.
size(); ++i) {
322 if ( !hepevt_particles[i]->end_vertex() &&
323 !hepevt_particles[i]->production_vertex() ) {
324 std::cout <<
" Pythia8ToHepMC::fill_next_event error: " 325 <<
"hanging particle " << i << std::endl;
326 GenVertex* prod_vtx =
new GenVertex();
327 prod_vtx->add_particle_out( hepevt_particles[i] );
328 evt->add_vertex( prod_vtx );
332 if ( doHadr && m_free_parton_exception ) {
333 int pdg_tmp = hepevt_particles[i]->pdg_id();
334 if ( (abs(pdg_tmp) <= 6 || pdg_tmp == 21)
335 && !hepevt_particles[i]->end_vertex() )
341 if (append)
return true;
345 if (m_store_pdf && pyinfo != 0) {
346 int id1pdf = pyinfo->
id1pdf();
347 int id2pdf = pyinfo->id2pdf();
348 if ( m_convert_gluon_to_0 ) {
349 if (id1pdf == 21) id1pdf = 0;
350 if (id2pdf == 21) id2pdf = 0;
354 evt->set_pdf_info( PdfInfo( id1pdf, id2pdf, pyinfo->x1pdf(),
355 pyinfo->x2pdf(), pyinfo->QFac(), pyinfo->pdf1(), pyinfo->pdf2() ) );
359 if (m_store_proc && pyinfo != 0) {
360 evt->set_signal_process_id( pyinfo->code() );
361 evt->set_event_scale( pyinfo->QRen() );
362 if (evt->alphaQED() <= 0) evt->set_alphaQED( pyinfo->alphaEM() );
363 if (evt->alphaQCD() <= 0) evt->set_alphaQCD( pyinfo->alphaS() );
368 if (m_store_xsec && pyinfo != 0) {
369 HepMC::GenCrossSection xsec;
370 xsec.set_cross_section( pyinfo->
sigmaGen() * 1e9,
372 evt->set_cross_section(xsec);
374 std::vector<double> xsecVec = pyinfo->weightContainerPtr->
getTotalXsec();
375 if (xsecVec.size() > 0) {
376 xsec.set_cross_section(xsecVec[0]*1e9);
377 evt->set_cross_section(xsec);
381 if (m_store_weights && pyinfo != 0) {
382 for (
int iweight = 0; iweight < pyinfo->numberOfWeights();
384 std::string name = pyinfo->weightNameByIndex(iweight);
385 double value = pyinfo->weightValueByIndex(iweight);
386 evt->weights()[name] = value;
417 typedef shared_ptr<GenEvent> EventPtr;
418 typedef HepMC::IO_GenEvent Writer;
419 typedef shared_ptr<Writer> WriterPtr;
426 setNewFile(filename, ft);
433 writerPtr = make_shared<HepMC::IO_GenEvent>(filename);
439 return writerPtr !=
nullptr;
445 geneve = make_shared<GenEvent>();
446 return fill_next_event(pythia, *geneve);
451 writerPtr->write_event(&*geneve);
458 if ( !fillNextEvent(pythia) )
return false;
485 HepMC::GenCrossSection xs;
486 xs.set_cross_section(xsec, xsecerr);
487 geneve->set_cross_section(xs);
492 if ( wv.size() != weightNames.size() )
493 geneve->weights() = wv;
495 for (
int i= 0, N = wv.size(); i < N; ++i )
496 geneve->weights()[weightNames[i]] = wv[i];
506 double scale,
double xf1,
double xf2,
507 int pdf1 = 0,
int pdf2 = 0) {
508 HepMC::PdfInfo pdf(id1, id2, x1, x2, scale, xf1, xf2, pdf1, pdf2);
509 geneve->set_pdf_info(pdf);
515 EventPtr geneve =
nullptr;
518 WriterPtr writerPtr =
nullptr;
521 vector<string> weightNames;
Pythia8ToHepMC()
Constructor and destructor.
Definition: HepMC2.h:83
void setWeightNames(const vector< string > &wnv)
Set all weight names in the current run.
Definition: HepMC2.h:500
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1478
GenEvent & event()
Get a reference to the current GenEvent.
Definition: HepMC2.h:464
int nCollND() const
Definition: HIInfo.h:94
The Event class holds all info on the generated event.
Definition: Event.h:453
Exception thrown when an undecayed parton is written into the record.
Definition: HepMC2.h:47
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: HepMC2.h:505
Pythia8ToHepMC(string filename, OutputType ft=ascii2)
Construct an object with an internal output stream.
Definition: HepMC2.h:425
vector< double > getTotalXsec()
Return cross section estimate for total run.
Definition: Weights.cc:1025
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:332
double b() const
The impact-parameter distance in the current event.
Definition: HIInfo.h:43
bool writeNextEvent(Pythia &pythia)
Definition: HepMC2.h:457
The Pythia8ToHepMC class.
Definition: HepMC2.h:78
double sigmaErr(int i=0) const
The uncertainty on the estimated cross-section in units of mb.
Definition: Info.h:314
Event event
The event record for the complete event history.
Definition: Pythia.h:323
Base exception for all exceptions that Pythia8ToHepMC might throw.
Definition: HepMC2.h:36
PartonEndVertexException(int i, int pdg_idIn)
Constructor and destructor.
Definition: HepMC2.h:52
int id1pdf(int i=0) const
Hard process flavours, x values, parton densities, couplings, Q2 scales.
Definition: Info.h:169
void setXSec(double xsec, double xsecerr)
Set cross section information in the current GenEvent.
Definition: HepMC2.h:484
void writeEvent()
Write out the current GenEvent to the internal stream.
Definition: HepMC2.h:450
virtual const char * what() const
Throw exception.
Definition: HepMC2.h:62
EventPtr getEventPtr()
Get a pointer to the current GenEvent.
Definition: HepMC2.h:469
Pythia8ToHepMC()
The empty constructor does not creat an aoutput stream.
Definition: HepMC2.h:422
HIInfo * hiInfo
Definition: Info.h:113
int nCollDD() const
Definition: HIInfo.h:109
bool fill_next_event(Pythia8::Event &pyev, GenEvent &evt, int ievnum=-1, const Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0, bool append=false, GenParticle *rootParticle=0, int iBarcode=-1)
Alternative method to convert Pythia events into HepMC ones.
Definition: HepMC2.h:100
double sigmaGen(int i=0) const
The estimated cross-section in units of mb.
Definition: Info.h:306
bool setNewFile(string filename, OutputType ft=ascii2)
Open a new external output stream.
Definition: HepMC2.h:430
int size() const
Event record size.
Definition: Event.h:504
int nAbsTarg() const
Definition: HIInfo.h:139
WriterPtr outputPtr()
Get a pointer to the internal stream.
Definition: HepMC2.h:479
bool print_inconsistency() const
Read out values for some switches.
Definition: HepMC2.h:113
int nAbsProj() const
Definition: HIInfo.h:123
void setWeights(const vector< double > &wv)
Update all weights in the current GenEvent.
Definition: HepMC2.h:491
HepMC::GenEvent GenEvent
Typedef for the version 2 specific classes used.
Definition: HepMC2.h:416
void set_print_inconsistency(bool b=true)
Set values for some switches.
Definition: HepMC2.h:123
int nCollNDTot() const
The total number of non-diffractive sub collisions in the current event.
Definition: HIInfo.h:97
OutputType
We can either have standard ascii output or none at all.
Definition: HepMC2.h:413
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: HepMC2.h:444
const Info & info
Public information and statistic on the generation.
Definition: Pythia.h:326
Writer & output()
Get a reference to the internal stream.
Definition: HepMC2.h:474
int nDiffTarg() const
Definition: HIInfo.h:143
int index() const
Return info on location and flavour of bad parton.
Definition: HepMC2.h:65
int nDiffProj() const
Definition: HIInfo.h:127
bool fill_next_event(Pythia8::Pythia &pythia, GenEvent &evt, int ievnum=-1, bool append=false, GenParticle *rootParticle=0, int iBarcode=-1)
The recommended method to convert Pythia events into HepMC ones.
Definition: HepMC2.h:90
Definition: Settings.h:195