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 multiple 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 occured.

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 diffractively excited?

bool Info::isMinBias()  
is the process a minimum-bias one?

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 minbias 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 occured 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 parton densities and scales

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.

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.

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

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.

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.

Hard process kinematics

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.

Event weight and activity

double Info::weight()  
weight assigned to the current event. Is normally 1 and thus uninteresting. However, for Les Houches events some strategies allow negative weights, which then after unweighting lead to events with weight -1. There are also strategies where no unweighting is done, and therefore a nontrivial event weight must be used e.g. when filling histograms.

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 decys, and in the final-state showering inside resonance decays, respectively.

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

Multiple interactions

double Info::bMI()  
the impact parameter b assumed for the current collision when multiple 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::enhanceMI()  
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 be unity for minimum-bias events (meaning more than that for events with hard processes).

int Info::nMI()  
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 multiple interactions.

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

int Info::iAMI(i)  
int Info::iBMI(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. iAMI and iBMI then denote partons coming from the first or second beam, respectively.

Cross sections

Here are the currently available methods related to the event sample as a whole. 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.

long Info::nTried()  
long Info::nSelected()  
long Info::nAccepted()  
the total number of tried phase-space points, selected hard processes and finally accepted events, summed over all allowed subprocesses. 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()  
double Info::sigmaErr()  
the estimated cross section and its estimated error, summed over all allowed subprocesses, 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 MI/ISR/FSR evolution downwards in pT for the current system.
argumentoption 23 : the number of time MI has been selected for the downwards step above.
argumentoption 24 : the number of time ISR has been selected for the downwards step above.
argumentoption 25 : the number of time FSR has been selected for the downwards step above.
argumentoption 26 : the number of time MI has been accepted as the downwards step above, after the vetoes.
argumentoption 27 : the number of time ISR has been accepted as the downwards step above, after the vetoes.
argumentoption 28 : the number of time 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 time FSR has been selected for the downwards step above.
argumentoption 31 : the number of time 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 30 - 39.
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 30 - 39.
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.