main134
Back to index.
// main134.cc is a part of the PYTHIA event generator.
// Copyright (C) 2024 Torbjorn Sjostrand.
// PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
// Please respect the MCnet Guidelines, see GUIDELINES for details.
// Author: Stefan Prestel
// Contact: Christian Preuss <preuss@uni-wuppertal.de>
// Keywords:
// LHE file
// Hepmc
// This program illustrates how a file with HepMC2 or HepMC3 events
// can be generated by Pythia8. Input and output files are specified
// on the command line, e.g. like
// ./main134 main134.cmnd main134.hepmc > main134.log
#include "Pythia8/Pythia.h"
#ifndef HEPMC2
#include "Pythia8Plugins/HepMC3.h"
#else
#include "Pythia8Plugins/HepMC2.h"
#endif
#include <unistd.h>
using namespace Pythia8;
//==========================================================================
// Example main programm to illustrate merging.
int main( int argc, char* argv[] ) {
// Check that correct number of command-line arguments
if (argc != 3) {
cerr << " Unexpected number of command-line arguments ("<<argc<<"). \n"
<< " You are expected to provide the arguments" << endl
<< " 1. Input file for settings" << endl
<< " 2. Output file for HepMC events" << endl
<< " Program stopped. " << endl;
return 1;
}
// Input parameters.
Pythia pythia;
pythia.readFile(argv[1],0);
// Interface for conversion from Pythia8::Event to HepMC one.
// Specify file where HepMC events will be stored.
Pythia8ToHepMC toHepMC(argv[2]);
// Allow abort of run if many errors.
int nAbort = pythia.mode("Main:timesAllowErrors");
int iAbort = 0;
bool doAbort = false;
// Read in loop parameters.
cout << endl << endl << endl;
cout << "Start generating events" << endl;
long nEvent = pythia.settings.mode("Main:numberOfEvents");
int nRuns = pythia.mode("Main:numberOfSubruns");
double sigmaTotal(0.), errorTotal(0.);
// Loop over subruns with varying number of jets.
for (int iRuns = 0; iRuns < nRuns; ++iRuns) {
double sigmaSample = 0., errorSample = 0.;
// Read in name of LHE file for current subrun and initialize.
pythia.readFile(argv[1], iRuns);
// If Pythia fails to initialize, exit with error.
if (!pythia.init()) return 1;
// Get the inclusive x-section by summing over all process x-sections.
double xs = 0.;
for (int i=0; i < pythia.info.nProcessesLHEF(); ++i)
xs += pythia.info.sigmaLHEF(i);
// Start generation loop.
while( pythia.info.nSelected() < nEvent ) {
// Generate next event.
if( !pythia.next() ) {
if ( pythia.info.atEndOfFile() ) break;
else if (++iAbort > nAbort) {doAbort = true; break;}
else continue;
}
// Get event weight(s).
double evtweight = pythia.info.weight();
// Do not print zero-weight events.
if ( evtweight == 0. ) continue;
// Inform HepMC3 about the naming of the weights.
#ifndef HEPMC2
toHepMC.setWeightNames(pythia.info.weightNameVector());
#endif
// Work with weighted (LHA strategy=-4) events.
double normhepmc = 1.;
if (abs(pythia.info.lhaStrategy()) == 4)
normhepmc = 1. / double(1e9*nEvent);
// Work with unweighted events.
else
normhepmc = xs / double(1e9*nEvent);
// Set event weights (optional).
#ifndef HEPMC2
// hepmcevt.weights().push_back(evtweight*normhepmc);
#endif
// Fill a new HepMC event.
toHepMC.fillNextEvent( pythia );
// Add the weight of the current event to the cross section.
sigmaTotal += evtweight*normhepmc;
sigmaSample += evtweight*normhepmc;
errorTotal += pow2(evtweight*normhepmc);
errorSample += pow2(evtweight*normhepmc);
// Report cross section to HepMC.
toHepMC.setXSec( sigmaTotal*1e9, errorTotal*1e9 );
// Write the HepMC event to file. Done with it.
toHepMC.writeEvent();
} // End loop over events to generate.
if (doAbort) break;
// Print cross section and errors.
pythia.stat();
cout << endl << " Contribution of sample " << iRuns
<< " to the inclusive cross section : "
<< scientific << setprecision(8)
<< sigmaSample << " +- " << sqrt(errorSample) << endl;
}
// Abort warning.
cout << endl << endl << endl;
if (doAbort)
cout << " Run was not completed owing to too many aborted events" << endl;
cout << endl << endl << endl;
// Done
return 0;
}