PYTHIA  8.312
Public Types | Public Member Functions | List of all members
Pythia8ToHepMC Class Reference

#include <HepMC2.h>

Inheritance diagram for Pythia8ToHepMC:
Pythia8ToHepMC3 Pythia8ToHepMC

Public Types

enum  OutputType {
  none, ascii2, none, ascii2,
  ascii3
}
 We can either have standard ascii output or none at all.
 
enum  OutputType {
  none, ascii2, none, ascii2,
  ascii3
}
 
typedef HepMC::GenEvent GenEvent
 Typedef for the version 2 specific classes used.
 
typedef shared_ptr< GenEventEventPtr
 
typedef HepMC::IO_GenEvent Writer
 
typedef shared_ptr< Writer > WriterPtr
 
typedef HepMC3::GenEvent GenEvent
 Typedef for the version 3 specific classes used.
 
typedef shared_ptr< GenEventEventPtr
 
typedef HepMC3::Writer Writer
 
typedef shared_ptr< Writer > WriterPtr
 

Public Member Functions

 Pythia8ToHepMC ()
 The empty constructor does not creat an aoutput stream.
 
 Pythia8ToHepMC (string filename, OutputType ft=ascii2)
 Construct an object with an internal output stream.
 
bool setNewFile (string filename, OutputType ft=ascii2)
 Open a new external output stream.
 
bool fillNextEvent (Pythia &pythia)
 
void writeEvent ()
 Write out the current GenEvent to the internal stream.
 
bool writeNextEvent (Pythia &pythia)
 
GenEventevent ()
 Get a reference to the current GenEvent.
 
EventPtr getEventPtr ()
 Get a pointer to the current GenEvent.
 
Writer & output ()
 Get a reference to the internal stream.
 
WriterPtr outputPtr ()
 Get a pointer to the internal stream.
 
void setXSec (double xsec, double xsecerr)
 Set cross section information in the current GenEvent.
 
void setWeights (const vector< double > &wv)
 Update all weights in the current GenEvent.
 
void setWeightNames (const vector< string > &wnv)
 Set all weight names in the current run.
 
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.
 
 Pythia8ToHepMC ()
 The empty constructor does not creat an aoutput stream.
 
 Pythia8ToHepMC (string filename, OutputType ft=ascii3)
 Construct an object with an internal output stream.
 
bool setNewFile (string filename, OutputType ft=ascii3)
 Open a new external output stream.
 
bool fillNextEvent (Pythia &pythia)
 
void writeEvent ()
 Write out the current GenEvent to the internal stream.
 
bool writeNextEvent (Pythia &pythia)
 
GenEventevent ()
 Get a reference to the current GenEvent.
 
EventPtr getEventPtr ()
 Get a pointer to the current GenEvent.
 
Writer & output ()
 Get a reference to the internal stream.
 
WriterPtr outputPtr ()
 Get a pointer to the internal stream.
 
void setXSec (double xsec, double xsecerr)
 Set cross section information in the current GenEvent.
 
void setWeights (const vector< double > &wv)
 Update all weights in the current GenEvent.
 
void setWeightNames (const vector< string > &wnv)
 Set all weight names in the current run.
 
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.
 
template<class T >
void addAtribute (const string &name, T &attribute)
 
template<class T = double>
void addAttribute (const string &name, double &attribute)
 Add an attribute of double type.
 
void removeAttribute (const string &name)
 Remove an attribute from the current event.
 
- Public Member Functions inherited from Pythia8ToHepMC3
 Pythia8ToHepMC3 ()
 Constructor and destructor.
 
bool fill_next_event (Pythia8::Pythia &pythia, GenEvent *evt, int ievnum=-1)
 The recommended method to convert Pythia events into HepMC3 ones.
 
bool fill_next_event (Pythia8::Pythia &pythia, GenEvent &evt)
 
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.
 
bool fill_next_event (Pythia8::Event &pyev, GenEvent *evt, int ievnum=-1, const Pythia8::Info *pyinfo=0, Pythia8::Settings *pyset=0)
 
bool print_inconsistency () const
 Read out values for some switches.
 
bool free_parton_warnings () const
 
bool crash_on_problem () const
 
bool convert_gluon_to_0 () const
 
bool store_pdf () const
 
bool store_proc () const
 
bool store_xsec () const
 
bool store_weights () const
 
void set_print_inconsistency (bool b=true)
 Set values for some switches.
 
void set_free_parton_warnings (bool b=true)
 
void set_crash_on_problem (bool b=false)
 
void set_convert_gluon_to_0 (bool b=false)
 
void set_store_pdf (bool b=true)
 
void set_store_proc (bool b=true)
 
void set_store_xsec (bool b=true)
 
void set_store_weights (bool b=true)
 
- Public Member Functions inherited from Pythia8ToHepMC
 Pythia8ToHepMC ()
 Constructor and destructor.
 
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.
 
bool fill_next_event (Pythia8::Pythia &pythia, GenEvent *evt, int ievnum=-1, bool append=false, GenParticle *rootParticle=0, int iBarcode=-1)
 
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.
 
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)
 
bool print_inconsistency () const
 Read out values for some switches.
 
bool free_parton_exception () const
 
bool free_parton_warnings () const
 
bool convert_gluon_to_0 () const
 
bool store_pdf () const
 
bool store_proc () const
 
bool store_xsec () const
 
bool store_weights () const
 
void set_print_inconsistency (bool b=true)
 Set values for some switches.
 
void set_free_parton_exception (bool b=true)
 
void set_free_parton_warnings (bool b=true)
 
void set_convert_gluon_to_0 (bool b=false)
 
void set_store_pdf (bool b=true)
 
void set_store_proc (bool b=true)
 
void set_store_xsec (bool b=true)
 
void set_store_weights (bool b=true)
 

Detailed Description

This a wrapper around HepMC::Pythia8ToHepMC in the Pythia8 namespace that simplify the most common use cases. It stores the current GenEvent and output stream internally to avoid cluttering of user code. This class is also defined in HepMC3.h with the same signatures, and the user can therefore switch between HepMC version 2 and 3, by simply changing the include file.

This a wrapper around HepMC::Pythia8ToHepMC in the Pythia8 namespace that simplify the most common use cases. It stores the current GenEvent and output stream internally to avoid cluttering of user code. This class is also defined in HepMC2.h with the same signatures, and the user can therefore switch between HepMC version 2 and 3, by simply changing the include file.

Member Enumeration Documentation

enum OutputType

We can either have standard ascii output version 2 or three or none at all.

Member Function Documentation

void addAtribute ( const string &  name,
T &  attribute 
)
inline

Add an additional attribute derived from HepMC3::Attribute to the current event.

bool fillNextEvent ( Pythia pythia)
inline

Create a new GenEvent object and fill it with information from the given Pythia object.

bool fillNextEvent ( Pythia pythia)
inline

Create a new GenEvent object and fill it with information from the given Pythia object.

bool writeNextEvent ( Pythia pythia)
inline

Create a new GenEvent object and fill it with information from the given Pythia object and write it out directly to the internal stream.

bool writeNextEvent ( Pythia pythia)
inline

Create a new GenEvent object and fill it with information from the given Pythia object and write it out directly to the internal stream.


The documentation for this class was generated from the following files: