StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main62.cc
1 // main62.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Author: Mikhail Kirsanov, Mikhail.Kirsanov@cern.ch.
7 // This program illustrates how a file with HepMC events
8 // can be generated by Pythia8.
9 // It is similar to main42 and main61, except that it is compiled with LHAPDF
10 // and allows for several subruns, e.g. from related LHEF.
11 // Input and output files are specified on the command line, e.g. like
12 // ./main62.exe main62.cmnd hepmcout62.dat > out
13 // The main program contains no analysis; this is intended to happen later.
14 // It therefore "never" has to be recompiled to handle different tasks.
15 
16 // WARNING: typically one needs 25 MB/100 events at the LHC.
17 // Therefore large event samples may be impractical.
18 
19 #include "Pythia8/Pythia.h"
20 #include "Pythia8/Pythia8ToHepMC.h"
21 #include "HepMC/GenEvent.h"
22 #include "HepMC/IO_GenEvent.h"
23 
24 using namespace Pythia8;
25 
26 int main(int argc, char* argv[]) {
27 
28  // Check that correct number of command-line arguments
29  if (argc != 3) {
30  cerr << " Unexpected number of command-line arguments. \n You are"
31  << " expected to provide one input and one output file name. \n"
32  << " Program stopped! " << endl;
33  return 1;
34  }
35 
36  // Check that the provided input name corresponds to an existing file.
37  ifstream is(argv[1]);
38  if (!is) {
39  cerr << " Command-line file " << argv[1] << " was not found. \n"
40  << " Program stopped! " << endl;
41  return 1;
42  }
43 
44  // Confirm that external files will be used for input and output.
45  cout << "\n >>> PYTHIA settings will be read from file " << argv[1]
46  << " <<< \n >>> HepMC events will be written to file "
47  << argv[2] << " <<< \n" << endl;
48 
49  // Interface for conversion from Pythia8::Event to HepMC event.
50  HepMC::Pythia8ToHepMC ToHepMC;
51 
52  // Specify file where HepMC events will be stored.
53  HepMC::IO_GenEvent ascii_io(argv[2], std::ios::out);
54 
55  // Generator.
56  Pythia pythia;
57 
58  // Read in subrun-independent commands from external file.
59  pythia.readFile( argv[1]);
60 
61  // Extract data to be used in main program. Set counters.
62  int nSubrun = pythia.mode("Main:numberOfSubruns");
63  int nAbort = pythia.mode("Main:timesAllowErrors");
64  int iAbort = 0;
65 
66  // Begin loop over subruns.
67  for (int iSubrun = 1; iSubrun <= nSubrun; ++iSubrun) {
68 
69  // Read in subrun-specific data from external file.
70  pythia.readFile( argv[1], iSubrun);
71 
72  // Initialization.
73  pythia.init();
74 
75  // Print name of Les Houches Event File.
76  string lheFile = pythia.settings.word("Beams:LHEF");
77  cout << "\n >>> Now begin subrun " << iSubrun
78  << " with events from file " << lheFile << " <<< \n"
79  << endl;
80 
81  // Begin infinite event loop - to be exited at end of file.
82  for (int iEvent = 0; ; ++iEvent) {
83 
84  // Generate event.
85  if (!pythia.next()) {
86 
87  // Leave event loop if at end of file.
88  if (pythia.info.atEndOfFile()) break;
89 
90  // First few failures write off as "acceptable" errors, then quit.
91  if (++iAbort < nAbort) continue;
92  cout << " Event generation aborted prematurely, owing to error!\n";
93  break;
94  }
95 
96  // Construct new empty HepMC event and fill it.
97  // Units will be as chosen for HepMC build, but can be changed
98  // by arguments, e.g. GenEvt( HepMC::Units::GEV, HepMC::Units::MM)
99  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
100  ToHepMC.fill_next_event( pythia, hepmcevt );
101 
102  // Write the HepMC event to file. Done with it.
103  ascii_io << hepmcevt;
104  delete hepmcevt;
105 
106  // End of event loop.
107  }
108 
109  // End of subrun loop.
110  }
111 
112  // Statistics.
113  pythia.stat();
114 
115  // Done.
116  return 0;
117 }
The GenEvent class is the core of HepMC.
Definition: GenEvent.h:155
IO_GenEvent also deals with HeavyIon and PdfInfo.
Definition: IO_GenEvent.h:63