Event Information

The Info class collects various one-of-a-kind information, some relevant for all events and others for the current event. An object info is a public member of the Pythia class, so if you e.g. have declared Pythia pythia, the Info methods can be accessed by pythia.info.method(). Most of this is information that could also be obtained e.g. from the event record, but is here more directly available. It is primarily intended for processes generated internally in PYTHIA, but many of the methods would work also for events fed in via the Les Houches Accord.

List information

void Info::list()  
a listing of most of the information set for the current event.

The beams

int Info::idA()  
int Info::idB()  
the identities of the two beam particles.

double Info::pzA()  
double Info::pzB()  
the longitudinal momenta of the two beam particles.

double Info::eA()  
double Info::eB()  
the energies of the two beam particles.

double Info::mA()  
double Info::mB()  
the masses of the two beam particles.

double Info::eCM()  
double Info::s()  
the CM energy and its square for the two beams.

Initialization

bool Info::tooLowPTmin()  
normally false, but true if the proposed pTmin scale was too low in timelike or spacelike showers, or in multiparton interactions. In the former case the pTmin is raised to some minimal value, in the latter the initialization fails (it is impossible to obtain a minijet cross section bigger than the nondiffractive one by reducing pTmin).

The event type

string Info::name()  
int Info::code()  
the name and code of the process that occurred.

int Info::nFinal()  
the number of final-state partons in the hard process.

bool Info::isResolved()  
are beam particles resolved, i.e. were PDF's used for the process?

bool Info::isDiffractiveA()  
bool Info::isDiffractiveB()  
is either beam soft diffractively excited?

bool Info::isDiffractiveC()  
is there soft central diffraction (a.k.a. double Pomeron exchange)?

bool Info::isHardDiffractiveA()  
bool Info::isHardDiffractiveB()  
is either beam hard diffractively excited?

bool Info::isNonDiffractive()  
is the process the SoftQCD:nonDiffractive one, i.e. corresponding to the full inelastic nondiffractive part of the total cross section. (Note that a hard process, say Z^0 production, normally is nondiffractive, but this is not what we aim at here, and so the method would return false, unless the Z^0 had been generated as part of the MPI machinery for the SoftQCD:nonDiffractive component.)

bool Info::isMinBias()  
the same as above, retained for backwards compatibility, but to be removed in PYTHIA 8.2.

bool Info::isLHA()  
has the process been generated from external Les Houches Accord information?

bool Info::atEndOfFile()  
true if a linked Les Houches class refuses to return any further events, presumably because it has reached the end of the file from which events have been read in.

bool Info::hasSub()  
does the process have a subprocess classification? Currently only true for nondiffractive and Les Houches events, where it allows the hardest collision to be identified.

string Info::nameSub()  
int Info::codeSub()  
int Info::nFinalSub()  
the name, code and number of final-state partons in the subprocess that occurred when hasSub() is true. For a minimum-bias event the code would always be 101, while codeSub() would vary depending on the actual hardest interaction, e.g. 111 for g g → g g. For a Les Houches event the code would always be 9999, while codeSub() would be the external user-defined classification code. The methods below would also provide information for such particular subcollisions.

Hard process initiators

The methods in this sections refer to the two initial partons of the hard 2 → n process (diffraction excluded; see below).

int Info::id1()  
int Info::id2()  
the identities of the two partons coming in to the hard process.

double Info::x1()  
double Info::x2()  
x fractions of the two partons coming in to the hard process.

double Info::y()  
double Info::tau()  
rapidity and scaled mass-squared of the hard-process subsystem, as defined by the above x values.

bool Info::isValence1()  
bool Info::isValence2()  
true if the two hard incoming partons have been picked to belong to the valence piece of the parton-density distribution, else false. Should be interpreted with caution. Information is not set if you switch off parton-level processing.

Hard process parton densities and scales

The methods in this section refer to the partons for which parton densities have been defined, in order to calculate the cross section of the hard process (diffraction excluded; see below).

These partons would normally agree with the ones above, the initiators of the 2 → n process, but it does not have to be so. Currently the one counterexample is POWHEG events [Ali10]. Here the original hard process could be 2 → (n-1). The NLO machinery at times would add an initial-state branching to give a 2 → n process with a changed initial state. In that case the values in this section refer to the original 2 → (n-1) state and the initiator ones above to the complete2 → n process. The Info::list() printout will contain a warning in such cases.

For external events in the Les Houches format, the pdf information is obtained from the optional #pdf line. When this information is absent, the parton identities and x values agree with the initiator ones above, while the pdf values are unknown and therefore set to vanish. The alpha_s and alpha_em values are part of the compulsory information. The factorization and renormalization scales are both equated with the one compulsory scale value in the Les Houches standard, except when a #pdf line provides the factorization scale separately. If alpha_s, alpha_em or the compulsory scale value are negative at input then new values are defined as for internal processes.

int Info::id1pdf()  
int Info::id2pdf()  
the identities of the two partons for which parton density values are defined.

double Info::x1pdf()  
double Info::x2pdf()  
x fractions of the two partons for which parton density values are defined.

double Info::pdf1()  
double Info::pdf2()  
parton densities x*f(x,Q^2) evaluated for the two incoming partons; could be used e.g. for reweighting purposes in conjunction with the idpdf, xpdf and QFac methods. Events obtained from external programs or files may not contain this information and, if so, 0 is returned.

double Info::QFac()  
double Info::Q2Fac()  
the Q or Q^2 factorization scale at which the densities were evaluated.

double Info::alphaS()  
double Info::alphaEM()  
the alpha_strong and alpha_electromagnetic values used for the hard process.

double Info::QRen()  
double Info::Q2Ren()  
the Q or Q^2 renormalization scale at which alpha_strong and alpha_electromagnetic were evaluated.

double Info::scalup()  
returns the stored SCALUP value for Les Houches events, and else zero. It may agree with both the QFac() and QRen() values, as explained above. However, to repeat, should the input SCALUP scale be negative, separate positive factorization and renormalization scales are calculated and set as for internally generated events. Furthermore, when PDF info is supplied for the Les Houches event, the factorization scale is set by this PDF info (scalePDF), which can disagree with SCALUP.

Hard process kinematics

The methods in this section provide info on the kinematics of the hard processes, with special emphasis on 2 → 2 (diffraction excluded; see below).

double Info::mHat()  
double Info::sHat()  
the invariant mass and its square for the hard process.

double Info::tHat()  
double Info::uHat()  
the remaining two Mandelstam variables; only defined for 2 → 2 processes.

double Info::pTHat()  
double Info::pT2Hat()  
transverse momentum and its square in the rest frame of a 2 → 2 processes.

double Info::m3Hat()  
double Info::m4Hat()  
the masses of the two outgoing particles in a 2 → 2 processes.

double Info::thetaHat()  
double Info::phiHat()  
the polar and azimuthal scattering angles in the rest frame of a 2 → 2 process.

Soft Diffraction

Information on the primary elastic or diffractive process (A B → A B, X1 B, A X2, X1 X2, A X B) can be obtained with the methods in the "Hard process kinematics" section above. The variables here obviously are s, t, u, ... rather than sHat, tHat, uHat, ..., but the method names remain to avoid unnecessary duplication. Most other methods are irrelevant for a primary elastic/diffractive process.

Central diffraction A B → A X B is a 2 → 3 process, and therefore most of the 2 → 2 variables are no longer relevant. The tHat() and uHat() methods instead return the two t values at the A → A and B → B vertices, and pTHat() the average transverse momentum of the three outgoing "particles", while thetaHat() and phiHat() are undefined.

While the primary interaction does not contain a hard process, the diffractive subsystems can contain them, but need not. Specifically, double diffraction can contain two separate hard subprocesses, which breaks the methods above. Most of them have been expanded with an optional argument to address properties of diffractive subsystems. This argument can take four values:

The argument is defined for all of the methods in the three sections above, "Hard process initiators", "Hard process parton densities and scales" and "Hard process kinematics", with the exception of the isValence methods. Also the four final methods of "The event type" section, the ...Sub() methods, take this argument. But recall that they will only provide meaningful answers, firstly if there is a system of the requested type, and secondly if there is a hard subprocess in this system. A simple check for this is that id1() has to be nonvanishing. The methods below this section do not currently provide information specific to diffractive subsystems, e.g. the MPI information is not bookkept in such cases.

Hard Diffraction

Information on the momentum fraction taken from the beam and the momentum transfer in the hard diffractive process. Note that when side A is diffractively exited, then the Pomeron has been taken from side B and vice versa.

double Info::xPomeronA()  
double Info::xPomeronB()  
x fractions of momenta carried by the Pomeron in the hard diffractive process.

double Info::tPomeronA()  
double Info::tPomeronB()  
The momentum transfer t in the hard diffractive process.

Photons from lepton beams

Information about the kinematics of photon-photon collisions from lepton beams.

double Info::eCMsub()  
Collision energy of the gamma-gamma sub-system.

double Info::xGammaA()  
double Info::xGammaB()  
x fractions of lepton momenta carried by the photons.

double Info::Q2GammaA()  
double Info::Q2GammaB()  
Virtualities of the photons emitted by the leptons.

double Info::thetaScatLepA()  
double Info::thetaScatLepB()  
Scattering angles of the leptons wrt. the beam direction.

int Info::photonMode()  
Type of photon process, see Photoproduction for details.

Event weight and activity

double Info::weight()  
weight assigned to the current event. Is normally 1 and thus uninteresting. However, there are several cases where one may have nontrivial event weights. These weights must the be used e.g. when filling histograms.
(i) In the PhaseSpace:increaseMaximum = off default strategy, an event with a differential cross-section above the assumed one (in a given phase-space point) is assigned a weight correspondingly above unity. This should happen only very rarely, if at all, and so could normally be disregarded.
(ii) The User Hooks class offers the possibility to bias the selection of phase space points, which means that events come with a compensating weight, stored here.
(iii) For Les Houches events some strategies allow negative weights, which then after unweighting lead to events with weight -1. There are also Les Houches strategies where no unweighting is done, so events come with a weight. Specifically, for strategies +4 and -4, the event weight is in units of pb. (Internally in mb, but converted at output.)

double Info::weightSum()  
Sum of weights accumulated during the run. For unweighted events this agrees with the number of generated events. In order to obtain histograms normalized "per event", at the end of a run, histogram contents should be divided by this weight. (And additionally divided by the bin width.) Normalization to cross section also required multiplication by sigmaGen() below.

int Info::lhaStrategy()  
normally 0, but if Les Houches events are input then it gives the event weighting strategy, see Les Houches Accord.

int Info::nISR()  
int Info::nFSRinProc()  
int Info::nFSRinRes()  
the number of emissions in the initial-state showering, in the final-state showering excluding resonance decays, and in the final-state showering inside resonance decays, respectively.

double Info::pTmaxMPI()  
double Info::pTmaxISR()  
double Info::pTmaxFSR()  
Maximum pT scales set for MPI, ISR and FSR, given the process type and scale choice for the hard interactions. The actual evolution will run down from these scales.

double Info::pTnow()  
The current pT scale in the combined MPI, ISR and FSR evolution. Useful for classification in user hooks, but not once the event has been evolved.

double Info::mergingWeight()  
combined leading-order merging weight assigned to the current event, if tree-level multi-jet merging (i.e. CKKW-L or UMEPS merging) is attempted. If tree-level multi-jet merging is performed, all histograms should be filled with this weight, as discussed in CKKW-L Merging and UMEPS Merging.

double Info::mergingWeightNLO()  
combined NLO merging weight assigned to the current event, if NLO multi-jet merging (i.e. NL3 or UNLOPS merging) is attempted. If NLO multi-jet merging is performed, all histograms should be filled with this weight, as discussed in NLO Merging.

Multiparton interactions

As already noted, these methods do not make sense for diffractive topologies, and should not be used there. Partly this is physics, but mainly it is for technical reasons, e.g. that double diffraction involves two separate systems that would have to be bookkept as such.

double Info::a0MPI()  
The value of a0 when an x-dependent matter profile is used, MultipartonInteractions:bProfile = 4.

double Info::bMPI()  
The impact parameter b assumed for the current collision when multiparton interactions are simulated. Is not expressed in any physical size (like fm), but only rescaled so that the average should be unity for minimum-bias events (meaning less than that for events with hard processes).

double Info::enhanceMPI()  
The choice of impact parameter implies an enhancement or depletion of the rate of subsequent interactions, as given by this number. Again the average is normalized to be unity for minimum-bias events (meaning more than that for events with hard processes).

double Info::enhanceMPIavg()  
The average enhancement factor expected for hard processes in those cases where it can be calculated already at initialization, i.e. excluding the x-dependent b profile.

int Info::nMPI()  
The number of hard interactions in the current event. Is 0 for elastic and diffractive events, and else at least 1, with more possible from multiparton interactions.

int Info::codeMPI(int i)  
double Info::pTMPI(int i)  
the process code and transverse momentum of the i'th subprocess, with i in the range from 0 to nMPI() - 1. The values for subprocess 0 is redundant with information already provided above.

int Info::iAMPI(int i)  
int Info::iBMPI(int i)  
are normally zero. However, if the i'th subprocess is a rescattering, i.e. either or both incoming partons come from the outgoing state of previous scatterings, they give the position in the event record of the outgoing-state parton that rescatters. iAMPI and iBMPI then denote partons coming from the first or second beam, respectively.

double Info::eMPI(int i)  
The enhancement or depletion of the rate of the i'th subprocess. Is primarily of interest for the MultipartonInteractions:bProfile = 4 option, where the size of the proton depends on the x values of the colliding partons. Note that eMPI(0) = enhanceMPI().

double Info::bMPIold()  
double Info::enhanceMPIold()  
double Info::enhanceMPIoldavg()  
These methods are only relevant for hard diffraction with the requirement of no MPI in the hadron-hadron collision. Then an impact parameter and associated enhancement factor is picked for this collision, but afterwards overwritten when the Pomeron-hadron subcollision is considered. In such cases the old hadron-hadron values can be found here, while bMPI, enhanceMPI and enhanceMPIavg provide the new Pomeron-hadron ones.

Cross sections

Here are the currently available methods related to the event sample as a whole, for the default value i = 0, and otherwise for the specific process code provided as argument. This is the number obtained with Info::code(), while the further subdivision given by Info::codeSub() is not bookkept. While continuously updated during the run, it is recommended only to study these properties at the end of the event generation, when the full statistics is available. The individual process results are not available if a second hard process has been chosen, but can be gleaned from the pythia.stat() output.

vector<int> Info::codesHard()  
returns a vector with all the process codes set up for the current run, i.e. the valid nonzero arguments for the five methods below.

string Info::nameProc(int i = 0)  
returns the process name for process code i.

long Info::nTried(int i = 0)  
long Info::nSelected(int i = 0)  
long Info::nAccepted(int i = 0)  
the total number of tried phase-space points, selected hard processes and finally accepted events, summed over all allowed processes (i = 0) or for the given process. The first number is only intended for a study of the phase-space selection efficiency. The last two numbers usually only disagree if the user introduces some veto during the event-generation process; then the former is the number of acceptable events found by PYTHIA and the latter the number that also were approved by the user. If you set a second hard process there may also be a mismatch.

double Info::sigmaGen(int i = 0)  
double Info::sigmaErr(int i = 0)  
the estimated cross section and its estimated error, summed over all allowed processes (i = 0) or for the given process, in units of mb. The numbers refer to the accepted event sample above, i.e. after any user veto.

Loop counters

Mainly for internal/debug purposes, a number of loop counters from various parts of the program are stored in the Info class, so that one can keep track of how the event generation is progressing. This may be especially useful in the context of the User Hooks facility.

int Info::getCounter(int i)  
the method that gives you access to the value of the various loop counters.
argument i : the counter number you want to access:
argumentoption 0 - 9 : counters that refer to the run as a whole, i.e. are set 0 at the beginning of the run and then only can increase.
argumentoption 0 : the number of successful constructor calls for the Pythia class (can only be 0 or 1).
argumentoption 1 : the number of times a Pythia::init() call has been begun.
argumentoption 2 : the number of times a Pythia::init() call has been completed successfully.
argumentoption 3 : the number of times a Pythia::next() call has been begun.
argumentoption 4 : the number of times a Pythia::next() call has been completed successfully.
argumentoption 10 - 19 : counters that refer to each individual event, and are reset and updated in the top-level Pythia::next() method.
argumentoption 10 : the number of times the selection of a new hard process has been begun. Normally this should only happen once, unless a user veto is set to abort the current process and try a new one.
argumentoption 11 : the number of times the selection of a new hard process has been completed successfully.
argumentoption 12 : as 11, but additionally the process should survive any user veto and go on to the parton- and hadron-level stages.
argumentoption 13 : as 11, but additionally the process should survive the parton- and hadron-level stage and any user cuts.
argumentoption 14 : the number of times the loop over parton- and hadron-level processing has begun for a hard process. Is reset each time counter 12 above is reached.
argumentoption 15 : the number of times the above loop has successfully completed the parton-level step.
argumentoption 16 : the number of times the above loop has successfully completed the checks and user vetoes after the parton-level step.
argumentoption 17 : the number of times the above loop has successfully completed the hadron-level step.
argumentoption 18 : the number of times the above loop has successfully completed the checks and user vetoes after the hadron-level step.
argumentoption 20 - 39 : counters that refer to a local part of the individual event, and are reset at the beginning of this part.
argumentoption 20 : the current system being processed in PartonLevel::next(). Is almost always 1, but for double diffraction the two diffractive systems are 1 and 2, respectively.
argumentoption 21 : the number of times the processing of the current system (see above) has begun.
argumentoption 22 : the number of times a step has begun in the combined MPI/ISR/FSR evolution downwards in pT for the current system.
argumentoption 23 : the number of times MPI has been selected for the downwards step above.
argumentoption 24 : the number of times ISR has been selected for the downwards step above.
argumentoption 25 : the number of times FSR has been selected for the downwards step above.
argumentoption 26 : the number of times MPI has been accepted as the downwards step above, after the vetoes.
argumentoption 27 : the number of times ISR has been accepted as the downwards step above, after the vetoes.
argumentoption 28 : the number of times FSR has been accepted as the downwards step above, after the vetoes.
argumentoption 29 : the number of times a step has begun in the separate (optional) FSR evolution downwards in pT for the current system.
argumentoption 30 : the number of times FSR has been selected for the downwards step above.
argumentoption 31 : the number of times FSR has been accepted as the downwards step above, after the vetoes.
argumentoption 40 - 49 : counters that are unused (currently), and that therefore are free to use, with the help of the two methods below.

void Info::setCounter(int i, int value = 0)  
set the above counters to a given value. Only to be used by you for the unassigned counters 40 - 49.
argument i : the counter number, see above.
argument value (default = 0) : set the counter to this number; normally the default value is what you want.

void Info::addCounter(int i, int value = 0)  
increase the above counters by a given amount. Only to be used by you for the unassigned counters 40 - 49.
argument i : the counter number, see above.
argument value (default = 1) : increase the counter by this amount; normally the default value is what you want.

Parton shower history

The following methods are mainly intended for internal use, e.g. for matrix-element matching.

void Info::hasHistory(bool hasHistoryIn)  
bool Info::hasHistory()  
set/get knowledge whether the likely shower history of an event has been traced.

void Info::zNowISR(bool zNowIn)  
double Info::zNowISR()  
set/get value of z in latest ISR branching.

void Info::pT2NowISR(bool pT2NowIn)  
double Info::pT2NowISR()  
set/get value of pT^2 in latest ISR branching.

Les Houches Event File 3.0 information

Les Houches Event files can conform to version 1.0 and version 3.0 of the standard (version 2.0 having been extended to 3.0). The LHEF version of an input file can can be accessed by
int Info::LHEFversion()  

The Info class also provides a suitable interface to the information stored after reading Les Houches Event files in the updated format [But14]. An example main program using LHEF 3.0 information is main38.cc.

LHEF 3.0 offers new features both in the initialisation and the event sections of the input files. Possible information include extended use of XML tags in the <header> and <init> blocks. The LHEF 3.0 information is stored in a series struct's:

--     The <initrwgt> tag is a container tag for weight and weightgroup tags. This information is stored internally in LHAinitrwgt. Currently, there is no dedicated output for this tag. However, all the information stored in the tag can be retrieved by using the Info class member pointer LHAinitrwgt Info::initrwgt.

--     Multiple <weightgroup> tags: Container tag for weight tags. Currently, there is no dedicated output for this tag. However, all the information stored in the tag can be retrieved by using the Info class member pointer vector<LHAweightgroups> * Info::weightgroups.

--     Multiple <weight> tags: Tag defining auxiliary information on an event weight, e.g. the identifier and information on what the weight represents. All the information stored in the tag can be retrieved by using the Info class member pointer vector<LHAweightgroups> * Info::init_weights. This vector contains all <weight> tags in the <initrwgt> container and its subcontainer <weightgroup> tags. The size of the vector can be accessed through the method
int Info::getInitrwgtSize()  

--     Multiple <generator> tags: Store information on the generators used in the event generation. All the information stored in the tag can be retrieved by using the Info class member pointer vector<LHAgenerators> * Info::generators. More easy-to-use output functions are available. The size of this vector can be obtained from
int Info::getGeneratorSize()  

The complete header can be obtained with the Info class member string getHeaderBlock().

The contents of a <generator> tag can be accessed through the method
string Info::getGeneratorValue(unsigned int n = 0)  
Return the contents of the n'th <generator> tag in the vector of tags.

Attributes of the <generator> tag (e.g. the generator name and version) can be accessed via
string Info::getGeneratorAttribute(unsigned int n, string key, bool doRemoveWhitespace = false)  
Return the value of the generator attribute named key for the n'th generator in the vector. Setting doRemoveWhitespace to true will return the value, stripped of any whitespace. An empty string is returned if the attribute named key does not exist.

To obtain information on cross sections, the following two methods can be used
int Info::nProcessesLHEF()  
return the number of processes for which the cross section is stored.
double Info::sigmaLHEF(int iProcess)  
return the cross section of the iProcess'th process.

Possible information also includes extended use of XML tags in the <event> blocks:

--     The <rwgt> tag is a container tag for wgt tags. Currently, there is no dedicated output for this tag. It can however be retrieved by using the Info class member pointer LHArwgt Info::rwgt.

--     Multiple <wgt> tags: Tag defining the event weight in the detailed version of LHEF 3.0. All the information stored in the tag can be retrieved by using the Info class member pointer vector<LHAwgt> * Info::weights_detailed. More easy-to-use output functions are available. The size of this vector can be obtained from
unsigned int Info::getWeightsDetailedSize()  

A convenient access point to the information stored in the <wgt> tags is the Info class member vector<double> Info::weights_detailed_vector. The entries of this vector are ordered according to how <wgt> tags appear in the event block.

The contents of a <wgt> tag can be accessed through the method
double Info::getWeightsDetailedValue(string n)  
Return the value of the n'th <wgt> tag in the event.

Attributes of the <wgt> tag (e.g. the weight id) can be accessed via
string Info::getWeightsDetailedAttribute(string n, string key, bool doRemoveWhitespace = false)  
Return the value of the wgt attribute named key for the n'th wgt in the vector. Setting doRemoveWhitespace to true will return the value, stripped of any whitespace. An empty string is returned if the attribute named key does not exist.

--     The <weights> tag: Tag containing a vector of double entries for weights in the compressed version of LHEF 3.0. All the information stored in the tag can be retrieved by using the Info class member pointer LHAweights * Info::weights and the vector vector<double> Info::weights_compressed. More easy-to-use output functions are available. The size of this vector can be obtained from
unsigned int Info::getWeightsCompressedSize()  

The n'th weight can be accessed through the method
double Info::getWeightsCompressedValue(unsigned int n)  

Attributes of the <weights> tag (not normally used) can be accessed via
string Info::getWeightsCompressedAttribute(string key, bool doRemoveWhitespace = false)  
Return the value of the <weights> tag's attribute named key. Setting doRemoveWhitespace to true will return the value, stripped of any whitespace. An empty string is returned if the attribute named key does not exist.

--     The <scales> tag: Contains information on different scales used by the matrix element generator. All the information stored in the tag can be retrieved by using the Info class member pointer LHAweights * Info::scales. More easy-to-use output functions are available. The contents of the scales tag can be obtained from
string Info::getScalesValue()  

However, note that the actual scale values are stored as attributes (called e.g. muf or mur). Attributes of the <scales> tag can be accessed via
double Info::getScalesAttribute(string key)  
Return the value of the <scales> tag's attribute named key. Not-a-number will be returned if the attribute named key does not exist.

Finally, arbitrary attributes of the <event> tag are supported. Attributes of the <event> tag can be accessed by
string Info::getEventAttribute(string key, bool doRemoveWhitespace = false)  
return the value of the event attribute named key. Setting doRemoveWhitespace to true will return the value, stripped of any whitespace. An empty string is returned if the attribute named key does not exist.

Additional comments appearing in the <event> tag can be obtained with the Info class member string getEventComments().

Header information

A simple string key/value store, mainly intended for accessing information that is stored in the header block of Les Houches Event (LHE) files. In principle, any LHAup derived class can set this header information, which can then be read out later. Although the naming convention is arbitrary, in practice, it is dictated by the XML-like format of LHE files, see Les Houches Accord for more details.

string Info::header(string key)  
return the header named key

vector <string> Info::headerKeys()  
return a vector of all header key names

void Info::setHeader(string key, string val)  
set the header named key with the contents of val