PYTHIA  8.303
Public Member Functions | Public Attributes | List of all members
Pythia Class Reference

The Pythia class contains the top-level routines to generate an event. More...

#include <Pythia.h>

Public Member Functions

 Pythia (string xmlDir="../share/Pythia8/xmldoc", bool printBanner=true)
 Constructor. (See Pythia.cc file.) More...
 
 Pythia (Settings &settingsIn, ParticleData &particleDataIn, bool printBanner=true)
 Constructor from pre-initialised ParticleData and Settings objects. More...
 
 Pythia (istream &settingsStrings, istream &particleDataStrings, bool printBanner=true)
 Constructor taking input from streams instead of files. More...
 
 ~Pythia ()
 Destructor.
 
 Pythia (const Pythia &)=delete
 Copy and = constructors cannot be used.
 
Pythiaoperator= (const Pythia &)=delete
 
bool checkVersion ()
 Check consistency of version numbers (called by constructors). More...
 
bool readString (string, bool warn=true)
 Read in one update for a setting or particle data from a single line. More...
 
bool readFile (string fileName, bool warn=true, int subrun=SUBRUNDEFAULT)
 Read in updates for settings or particle data from user-defined file. More...
 
bool readFile (string fileName, int subrun)
 
bool readFile (istream &is=cin, bool warn=true, int subrun=SUBRUNDEFAULT)
 
bool readFile (istream &is, int subrun)
 
bool setPDFPtr (PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn, PDFPtr pdfHardAPtrIn=nullptr, PDFPtr pdfHardBPtrIn=nullptr, PDFPtr pdfPomAPtrIn=nullptr, PDFPtr pdfPomBPtrIn=nullptr, PDFPtr pdfGamAPtrIn=nullptr, PDFPtr pdfGamBPtrIn=nullptr, PDFPtr pdfHardGamAPtrIn=nullptr, PDFPtr pdfHardGamBPtrIn=nullptr, PDFPtr pdfUnresAPtrIn=nullptr, PDFPtr pdfUnresBPtrIn=nullptr, PDFPtr pdfUnresGamAPtrIn=nullptr, PDFPtr pdfUnresGamBPtrIn=nullptr, PDFPtr pdfVMDAPtrIn=nullptr, PDFPtr pdfVMDBPtrIn=nullptr)
 Possibility to pass in pointers to PDF's. More...
 
bool setPDFAPtr (PDFPtr pdfAPtrIn)
 Routine to pass in pointers to PDF's. Usage optional. More...
 
bool setPDFBPtr (PDFPtr pdfBPtrIn)
 Routine to pass in pointers to PDF's. Usage optional. More...
 
bool setPhotonFluxPtr (PDFPtr photonFluxAIn, PDFPtr photonFluxBIn)
 Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
 
bool setLHAupPtr (LHAupPtr lhaUpPtrIn)
 Possibility to pass in pointer to external LHA-interfaced generator.
 
bool setDecayPtr (DecayHandlerPtr decayHandlePtrIn, vector< int > handledParticlesIn)
 Possibility to pass in pointer for external handling of some decays.
 
bool setRndmEnginePtr (RndmEngine *rndmEnginePtrIn)
 Possibility to pass in pointer for external random number generation.
 
bool setUserHooksPtr (UserHooksPtr userHooksPtrIn)
 Possibility to pass in pointer for user hooks.
 
bool addUserHooksPtr (UserHooksPtr userHooksPtrIn)
 Possibility to add further pointers to allow multiple user hooks.
 
bool setMergingPtr (MergingPtr mergingPtrIn)
 Possibility to pass in pointer for full merging class.
 
bool setMergingHooksPtr (MergingHooksPtr mergingHooksPtrIn)
 Possibility to pass in pointer for merging hooks.
 
bool setBeamShapePtr (BeamShapePtr beamShapePtrIn)
 Possibility to pass in pointer for beam shape.
 
bool setSigmaPtr (SigmaProcess *sigmaPtrIn, PhaseSpace *phaseSpacePtrIn=0)
 
bool setResonancePtr (ResonanceWidths *resonancePtrIn)
 Possibility to pass in pointer(s) for external resonance.
 
bool setShowerModelPtr (ShowerModelPtr showerModelPtrIn)
 Possibility to pass in pointer for external showers.
 
bool setHeavyIonsPtr (HeavyIonsPtr heavyIonsPtrIn)
 Possibility to pass in pointer for modelling of heavy ion collisions.
 
bool setHIHooks (HIUserHooksPtr hiHooksPtrIn)
 
HeavyIonsPtr getHeavyIonsPtr ()
 
ShowerModelPtr getShowerModelPtr ()
 Possibility to get the pointer to the parton-shower model.
 
bool setPartonVertexPtr (PartonVertexPtr partonVertexPtrIn)
 Possibility to pass in pointer for setting of parton space-time vertices.
 
bool init ()
 Initialize. More...
 
bool next ()
 Generate the next event. More...
 
bool next (double eCMin)
 Generate the next event, either with new energies or new beam momenta. More...
 
bool next (double eAin, double eBin)
 Variant of the main event-generation routine, for variable beam energies. More...
 
bool next (double pxAin, double pyAin, double pzAin, double pxBin, double pyBin, double pzBin)
 Variant of the main event-generation routine, for variable beam momenta. More...
 
int forceTimeShower (int iBeg, int iEnd, double pTmax, int nBranchMax=0)
 Generate only a single timelike shower as in a decay.
 
bool forceHadronLevel (bool findJunctions=true)
 Generate only the hadronization/decay stage. More...
 
bool moreDecays ()
 Special routine to allow more decays if on/off switches changed.
 
bool forceRHadronDecays ()
 Special routine to force R-hadron decay when not done before.
 
bool doLowEnergyProcess (int i1, int i2, int type)
 Do a low-energy collision between two hadrons in the event record.
 
double getLowEnergySigma (int i1, int i2, int type=0)
 Get low-energy cross section for two hadrons in the event record.
 
double getLowEnergySigma (int id1, int id2, double eCM12, double m1, double m2, int type=0)
 Get low-energy cross section for two hadrons standalone.
 
double getLowEnergySigma (int id1, int id2, double eCM12, int type=0)
 
double getLowEnergySlope (int id1, int id2, double eCM12, double m1, double m2, int type=2)
 Get b slope in elastic and diffractive interactions standalone.
 
void LHAeventList ()
 List the current Les Houches event.
 
bool LHAeventSkip (int nSkip)
 Skip a number of Les Houches events at input.
 
void stat ()
 Main routine to provide final statistics on generation. More...
 
bool flag (string key)
 Read in settings values: shorthand, not new functionality.
 
int mode (string key)
 
double parm (string key)
 
string word (string key)
 
PDFPtr getPDFPtr (int idIn, int sequence=1, string beam="A", bool resolved=true)
 Auxiliary to set parton densities among list of possibilities. More...
 

Public Attributes

Event process = {}
 The event record for the parton-level central process.
 
Event event = {}
 The event record for the complete event history.
 
const Infoinfo = infoPrivate
 Public information and statistic on the generation.
 
Settings settings = {}
 Settings: databases of flags/modes/parms/words to control run.
 
ParticleData particleData = {}
 ParticleData: the particle data table/database.
 
Rndm rndm = {}
 Random number generator.
 
CoupSM coupSM = {}
 Standard Model couplings, including alphaS and alphaEM.
 
CoupSUSY coupSUSY = {}
 SUSY couplings.
 
SLHAinterface slhaInterface = {}
 SLHA Interface.
 
PartonSystems partonSystems = {}
 The partonic content of each subcollision system (auxiliary to event).
 
MergingPtr mergingPtr = {}
 Merging object as wrapper for matrix element merging routines.
 
MergingHooksPtr mergingHooksPtr
 
HeavyIonsPtr heavyIonsPtr = {}
 Pointer to a HeavyIons object for generating heavy ion collisions.
 
HIUserHooksPtr hiHooksPtr = {}
 Pointer to a HIUserHooks object to modify heavy ion modelling.
 
HadronWidths hadronWidths = {}
 HadronWidths: the hadron widths data table/database.
 
BeamParticle beamA = {}
 The two incoming beams.
 
BeamParticle beamB = {}
 

Detailed Description

The Pythia class contains the top-level routines to generate an event.

Constructor & Destructor Documentation

Pythia ( string  xmlDir = "../share/Pythia8/xmldoc",
bool  printBanner = true 
)

Constructor. (See Pythia.cc file.)

Constructor.

Initialise / reset pointers and global variables.

Find path to data files, i.e. xmldoc directory location. Environment variable takes precedence, then constructor input, and finally the pre-processor constant XMLDIR.

Read in files with all flags, modes, parms and words.

Save XML path in settings.

Check that XML and header version numbers match code version number.

Read in files with all particle data.

Write the Pythia banner to output.

Not initialized until at the end of the init() call.

Pythia ( Settings settingsIn,
ParticleData particleDataIn,
bool  printBanner = true 
)

Constructor from pre-initialised ParticleData and Settings objects.

Constructor to copy settings and particle database from another Pythia object instead of XML files (to speed up multiple initialisations).

Initialise / reset pointers and global variables.

Copy XML path from existing Settings database.

Copy settings database and redirect pointers.

Check XML and header version numbers match code version number.

Copy particleData database and redirect pointers.

Write the Pythia banner to output.

Not initialized until at the end of the init() call.

Pythia ( istream &  settingsStrings,
istream &  particleDataStrings,
bool  printBanner = true 
)

Constructor taking input from streams instead of files.

Constructor from string streams.

Initialise / reset pointers and global variables.

Copy settings database

Reset pointers to pertain to this PYTHIA object.

Check XML and header version numbers match code version number.

Read in files with all particle data.

Write the Pythia banner to output.

Not initialized until at the end of the init() call.

Member Function Documentation

bool checkVersion ( )

Check consistency of version numbers (called by constructors).

Check for consistency of version numbers (called by constructors).

Check that XML version number matches code version number.

Check that header version number matches code version number.

All is well that ends well.

bool forceHadronLevel ( bool  findJunctions = true)

Generate only the hadronization/decay stage.

Generate only the hadronization/decay stage, using internal machinery. The "event" instance should already contain a parton-level configuration.

Can only generate event if initialization worked.

Check whether any junctions in system. (Normally done in ProcessLevel.) Avoid it if there are no final-state coloured partons.

Allow for CR before the hadronization.

Setup parton system for SK-I and SK-II colour reconnection. Require all final state particles to have the Ws as mothers.

save spare copy of event in case of failure.

Allow up to ten tries for CR.

Save spare copy of event in case of failure.

Allow up to ten tries for hadron-level processing.

Check whether any resonances need to be handled at process level.

Allow for showers if decays happened at process level.

Hadron-level: hadronization, decays.

If failure then warn, restore original configuration and try again.

Done for simpler option.

Optionally check final event for problems.

Done.

HeavyIonsPtr getHeavyIonsPtr ( )
inline

Possibility to get the pointer to a object modelling heavy ion collisions.

PDFPtr getPDFPtr ( int  idIn,
int  sequence = 1,
string  beam = "A",
bool  resolved = true 
)

Auxiliary to set parton densities among list of possibilities.

Routine to set up a PDF pointer.

Temporary pointer to be returned.

Data file directory.

One option is to treat a Pomeron like a pi0.

Check if photon beam inside proton.

Proton beam, normal or hard choice. Also used for neutron.

Use internal LHAgrid1 implementation for LHAPDF6 files.

Use sets from LHAPDF.

Use internal sets.

Quasi-real photons inside a (anti-)proton beam.

Find the resolved photon PDF to combine with the flux.

Set up the combination of flux and PDF for resolved photons.

Find the pre-set photon PDF, hard or normal.

Set the photon flux pointer and construct approximation. Use the existing machinery for external fluxes.

Externally provided flux.

Find the correct flux for given beam set with setPhotonFluxPtr().

Check that external flux exist and complain if not.

Classic EPA proton by Budnev, Ginzburg, Meledin and Serbo.

Check if Q^2 sampling on and turn off if necessary.

EPA approximation by Drees and Zeppenfeld.

Construct flux object when pointer succesfully created.

Pion beam (or, in one option, Pomeron beam).

If VMD process then scale PDF accordingly: f_a^VMD = alphaEM * (1/f_rho^2 + 1/f_omega^2 + 1/f_phi^2 + 1/f_J/psi)

  • f_a^pi0. COR: New value here includes J/psi

Use internal LHAgrid1 implementation for LHAPDF6 files.

Use sets from LHAPDF.

Use internal set.

Pomeron beam, if not treated like a pi0 beam.

Use internal LHAgrid1 implementation for LHAPDF6 files.

Use sets from LHAPDF.

A generic Q2-independent parametrization.

The H1 Q2-dependent parametrizations. Initialization requires files.

The parametrizations of Alvero, Collins, Terron and Whitmore.

Set up nuclear PDFs.

Which nPDF set to use.

Temporary pointer for storing proton PDF pointer.

Photon beam, either point-like (unresolved) or resolved.

For unresolved beam use the point-like PDF.

Point-like beam if unresolved photons.

Use different PDFs for hard process.

Find the name or number of the hard PDF set.

Use sets from LHAPDF. Only available for hard processes.

Or set up an internal set (though currently only one).

Set up the PDF.

Lepton beam: neutrino, resolved charged lepton or unresolved ditto. Also photon inside lepton PDFs.

Set up resolved photon inside lepton for beam A.

Find the pre-set photon PDF, hard or normal.

Get the mass of lepton and maximum virtuality of the photon.

Initialize the gamma-inside-lepton PDFs with internal photon flux.

Initialize the gamma-inside-lepton PDFs with external photon flux. Requires that the pointer to the flux set.

Set up resolved photon inside lepton for beam B.

Find the pre-set photon PDF, hard or normal.

Get the mass of lepton and maximum virtuality of the photon.

Initialize the gamma-inside-lepton PDFs with internal photon flux.

Initialize the gamma-inside-lepton PDFs with external photon flux.

Usual lepton PDFs.

External photon flux for direct-photon processes.

Dark matter beam set up as pointlike lepton.

Optionally allow extrapolation beyond x and Q2 limits.

Done.

bool init ( )

Initialize.

Routine to initialize with the variable values of the Beams kind.

Check that constructor worked.

Early catching of heavy ion mode.

Master choice of shower model.

Early readout, if return false or changed when no beams.

Check that changes in Settings and ParticleData have worked.

Initialize the random number generator.

Find which frame type to use.

Already get the Les Houches File name here, to have it for later new merging initialization (e.g. by new showers).

Set up values related to CKKW-L merging.

Set up values related to unitarised CKKW merging

Set up values related to NL3 NLO merging

Set up values related to unitarised NLO merging

Set/Reset the weights

Set up MergingHooks object for simple shower model.

Set up Merging object for simple shower model.

Initialization with internal processes: read in and set values.

Initialization with a Les Houches Event File or an LHAup object.

For file input: renew file stream or (re)new Les Houches object.

Header is optional, so use NULL pointer to indicate no value.

Check that file was properly opened.

For object input: at least check that not null pointer.

LHAup object generic abort using fileFound() routine.

Send in pointer to info. Store or replace LHA pointer in other classes.

If second time around, only with new file, then simplify. Optionally skip ahead a number of events at beginning of file.

Set up values related to merging hooks.

Store the name of the input LHEF for merging.

Set LHAinit information (in some external program).

Extract beams from values set in an LHAinit object.

Optionally skip ahead a number of events at beginning of file.

Set up values related to user hooks.

Setup objects for string interactions (swing, colour reconnection, shoving and rope hadronisation).

Back to common initialization. Reset error counters.

Initialize data members extracted from database.

Warn/abort for unallowed process and beam combinations.

Find out if beams are or have a resolved photon beam. The PDF:lepton2gamma is kept for backwards compatibility, now beamA2gamma and beamB2gamma are the master switches.

Check if resolved photons are needed.

Check if unresolved photons are needed.

Check if VMD sampling is required for beam A and/or B.

Turn off central diffraction for VMD processes.

Initialise merging hooks.

Check that combinations of settings are allowed; change if not.

Initialize the SM couplings (needed to initialize resonances).

Initialize SLHA interface (including SLHA/BSM couplings).

Reset couplingsPtr to the correct memory address.

Set headers to distinguish the two event listing kinds.

Final setup stage of particle data, notably resonance widths.

Read in files with particle widths.

Set up R-hadrons particle data, where relevant.

Set up and initialize setting of parton production vertices.

Prepare for low-energy QCD processes.

Set up and initialize the ShowerModel (if not provided by user).

Register shower model as physicsBase object (also sets pointers)

Initialise shower model

Get relevant pointers from shower models.

Initialize pointers in showers.

At this point, the mergingPtr should be set. If that is not the case, then the initialization should be aborted.

Set up values related to beam shape.

Check that beams and beam combination can be handled.

Simplified beam setup when no process level.

Full beam setup: first beam kinematics.

Set up pointers to PDFs.

Set up the two beams and the common remnant system.

Pass information whether the beam will contain a photon beam.

Init also unresolved PDF pointers for photon beams when needed.

Optionally set up new alternative beams for these Pomerons.

Initialise VMD beams from gammas (in leptons). Use pion PDF for VMDs.

Optionally set up photon beams from lepton beams if resolved photons.

Send info/pointers to process level for initialization.

Initialize shower handlers using intialized beams.

Initialize timelike showers already here, since needed in decays. The pointers to the beams are needed by some external plugin showers.

Alternatively only initialize resonance decays.

Send info/pointers to parton level for initialization.

Make pointer to shower available for merging machinery.

Alternatively only initialize final-state showers in resonance decays.

Set up shower variation groups if enabled

Send info/pointers to parton level for trial shower initialization.

Initialise the merging wrapper class.

Send info/pointers to hadron level for initialization. Note: forceHadronLevel() can come, so we must always initialize.

Optionally check particle data table for inconsistencies.

Optionally show settings and particle data, changed or all.

Listing options for the next() routine.

Init junction splitting.

Flags for colour reconnection.

Succeeded.

bool next ( )

Generate the next event.

Main routine to generate the next event, using internal machinery.

Check that constructor worked.

Flexible-use call at the beginning of each new event.

Check if the generation is taken over by the HeavyIons object. Allows HeavyIons::next to call next for this Pythia object without going into a loop.

Regularly print how many events have been generated.

Set/reset info counters specific to each event.

Simpler option when no hard process, i.e. mainly hadron level.

Optionally fetch in resonance decays from LHA interface.

Reset info and partonSystems arrays (while event record contains data).

Set correct energy for system.

Generate parton showers where appropriate.

Generate hadronization and decays.

Reset arrays.

Pick current beam valence flavours (for pi0, K0S, K0L, Pomeron).

Can only generate event if initialization worked.

Pick beam momentum spread and beam vertex.

Recalculate kinematics when beam momentum spread.

Simplified special treatment for low-energy nonperturbative collisions.

Optionally check final event for problems.

Outer loop over hard processes; only relevant for user-set vetoes.

Provide the hard process that starts it off. Only one try.

Reset the event information. Necessary if the previous event was read from LHEF, while the current event is not read from LHEF.

Update tried and selected events immediately after next event was generated. Note: This does not accumulate cross section.

Possibility for a user veto of the process-level event.

Possibility to perform matrix element merging for this event.

Apply possible merging scale cut.

Exit because of vanishing no-emission probability.

Redo resonance decays after the merging, in case the resonance structure has been changed because of reclusterings.

Possibility to stop the generation at this stage.

Save spare copy of process record in case of problems.

Allow up to ten tries for parton- and hadron-level processing.

Restore original process record if problems.

Reset event record and (extracted partons from) beam remnants.

Parton-level evolution: ISR, FSR, MPI.

Abort event generation if parton level is set to abort.

Skip to next hard process for failure owing to deliberate veto, or alternatively retry for the same hard process.

If hard diffractive event has been discarded retry partonLevel.

Else make a new try for other failures.

Possibility for a user veto of the parton-level event.

Boost to lab frame (before decays, for vertices).

Possibility to stop the generation at this stage.

Optionally check final event for problems.

Hadron-level: hadronization, decays.

If R-hadrons have been formed, then (optionally) let them decay.

Optionally check final event for problems.

Stop parton- and hadron-level looping if you got this far.

If event vetoed then to make a new try.

If event failed any other way (after ten tries) then give up.

Process- and parton-level statistics. Event scale.

End of outer loop over hard processes. Done with normal option.

List events.

Done.

bool next ( double  eCMin)

Generate the next event, either with new energies or new beam momenta.

Variant of the main event-generation routine, for variable CM energies.

Check that constructor worked.

Check that generation has been initialized for variable energy.

Check that the frameType matches the input provided.

Save input value.

Call regular next method for event generation.

bool next ( double  eAin,
double  eBin 
)

Variant of the main event-generation routine, for variable beam energies.

Check that constructor worked.

Check that generation has been initialized for variable energy.

Check that the frameType matches the input provided.

Save input values.

Call regular next method for event generation.

bool next ( double  pxAin,
double  pyAin,
double  pzAin,
double  pxBin,
double  pyBin,
double  pzBin 
)

Variant of the main event-generation routine, for variable beam momenta.

Check that constructor worked.

Check that generation has been initialized for variable energy.

Check that the frameType matches the input provided.

Save input value.

Call regular next method for event generation.

bool readFile ( string  fileName,
bool  warn = true,
int  subrun = SUBRUNDEFAULT 
)

Read in updates for settings or particle data from user-defined file.

Check that constructor worked.

Open file for reading.

Hand over real work to next method.

bool readFile ( istream &  is = cin,
bool  warn = true,
int  subrun = SUBRUNDEFAULT 
)

Read in updates for settings or particle data from user-defined stream (or file).

Check that constructor worked.

Read in one line at a time.

Check whether entering, leaving or inside commented-commands section.

Check whether entered new subrun.

Process the line if in correct subrun.

Reached end of input file.

bool readString ( string  line,
bool  warn = true 
)

Read in one update for a setting or particle data from a single line.

Check that constructor worked.

If empty line then done.

If Settings input stretches over several lines then continue with it.

If first character is not a letter/digit, then taken to be a comment.

Send on particle data to the ParticleData database.

Everything else sent on to Settings.

bool setHIHooks ( HIUserHooksPtr  hiHooksPtrIn)
inline

Possibility to pass a HIUserHooks pointer for modifying the behavior of the heavy ion modelling.

bool setPDFAPtr ( PDFPtr  pdfAPtrIn)

Routine to pass in pointers to PDF's. Usage optional.

Reset pointers to be empty.

Switch off external PDF's by zero as input.

Save pointers.

By default same pointers for hard-process PDF's.

Done.

bool setPDFBPtr ( PDFPtr  pdfBPtrIn)

Routine to pass in pointers to PDF's. Usage optional.

Reset pointers to be empty.

Switch off external PDF's by zero as input.

Save pointers.

By default same pointers for hard-process PDF's.

Done.

bool setPDFPtr ( PDFPtr  pdfAPtrIn,
PDFPtr  pdfBPtrIn,
PDFPtr  pdfHardAPtrIn = nullptr,
PDFPtr  pdfHardBPtrIn = nullptr,
PDFPtr  pdfPomAPtrIn = nullptr,
PDFPtr  pdfPomBPtrIn = nullptr,
PDFPtr  pdfGamAPtrIn = nullptr,
PDFPtr  pdfGamBPtrIn = nullptr,
PDFPtr  pdfHardGamAPtrIn = nullptr,
PDFPtr  pdfHardGamBPtrIn = nullptr,
PDFPtr  pdfUnresAPtrIn = nullptr,
PDFPtr  pdfUnresBPtrIn = nullptr,
PDFPtr  pdfUnresGamAPtrIn = nullptr,
PDFPtr  pdfUnresGamBPtrIn = nullptr,
PDFPtr  pdfVMDAPtrIn = nullptr,
PDFPtr  pdfVMDBPtrIn = nullptr 
)

Possibility to pass in pointers to PDF's.

Routine to pass in pointers to PDF's. Usage optional.

Switch off external PDF's by zero as input.

The two PDF objects cannot be one and the same.

Save pointers.

By default same pointers for hard-process PDF's.

Optionally allow separate pointers for hard process.

Optionally allow pointers for Pomerons in the proton.

Optionally allow pointers for Gammas in the leptons.

Optionally allow pointers for Hard PDFs for photons in the leptons.

Optionally allow pointers for unresolved PDFs.

Optionally allow pointers for unresolved PDFs for photons from leptons.

Optionally allow pointers for VMD in the gamma.

Done.

bool setSigmaPtr ( SigmaProcess sigmaPtrIn,
PhaseSpace phaseSpacePtrIn = 0 
)
inline

Possibility to pass in pointer(s) for external cross section, with option to include external phase-space generator(s).

void stat ( )

Main routine to provide final statistics on generation.

Print statistics on event generation.

Read out settings for what to include.

Statistics on cross section and number of events.

Statistics from other classes, currently multiparton interactions.

Merging statistics.

Summary of which and how many warnings/errors encountered.

Loop through all PhysicsBase-derived objects.

Member Data Documentation

MergingHooksPtr mergingHooksPtr

Pointer to MergingHooks object for user interaction with the merging. MergingHooks also more generally steers the matrix element merging.


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