StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main85.cc
1 // main85.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2011 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 // This program is written by Stefan Prestel.
7 // It illustrates how to do CKKW-L merging,
8 // see the Matrix Element Merging page in the online manual.
9 
10 #include "Pythia8/Pythia.h"
11 #include "Pythia8/Pythia8ToHepMC.h"
12 #include <unistd.h>
13 
14 #include "HepMC/GenEvent.h"
15 #include "HepMC/IO_GenEvent.h"
16 // Following line to be used with HepMC 2.04 onwards.
17 #include "HepMC/Units.h"
18 
19 using namespace Pythia8;
20 
21 //==========================================================================
22 
23 // Example main programm to illustrate merging
24 
25 int main( int argc, char* argv[] ){
26 
27  // Check that correct number of command-line arguments
28  if (argc != 4) {
29  cerr << " Unexpected number of command-line arguments ("<<argc<<"). \n"
30  << " You are expected to provide the arguments" << endl
31  << " 1. Input file for settings" << endl
32  << " 2. Name of the input LHE file (with path), up to the '_tree'"
33  << " identifier" << endl
34  << " 3. Output hepmc file name" << endl
35  << " Program stopped. " << endl;
36  return 1;
37  }
38 
39  Pythia pythia;
40 
41  // Input parameters:
42  pythia.readFile(argv[1]);
43  // Interface for conversion from Pythia8::Event to HepMC one.
44  HepMC::Pythia8ToHepMC ToHepMC;
45  // Specify file where HepMC events will be stored.
46  HepMC::IO_GenEvent ascii_io(argv[3], std::ios::out);
47  // Switch off warnings for parton-level events.
48  ToHepMC.set_print_inconsistency(false);
49  ToHepMC.set_free_parton_warnings(false);
50  // Do not store cross section information, as this will be done manually.
51  ToHepMC.set_store_pdf(false);
52  ToHepMC.set_store_proc(false);
53  ToHepMC.set_store_xsec(false);
54 
55  // Path to input events, with name up to the "_tree" identifier included.
56  string iPath = string(argv[2]);
57 
58  // Number of events.
59  int nEvent = pythia.mode("Main:numberOfEvents");
60  // Maximal number of additional LO jets.
61  int nMaxLO = pythia.mode("Merging:nJetMax");
62 
63 //-----------------------------------------------------------------------------
64 //-----------------------------------------------------------------------------
65 
66  // Switch off all showering and MPI when extimating the cross section after
67  // the merging scale cut.
68  bool fsr = pythia.flag("PartonLevel:FSR");
69  bool isr = pythia.flag("PartonLevel:ISR");
70  bool mpi = pythia.flag("PartonLevel:MPI");
71  bool had = pythia.flag("HadronLevel:all");
72  pythia.settings.flag("PartonLevel:FSR",false);
73  pythia.settings.flag("PartonLevel:ISR",false);
74  pythia.settings.flag("HadronLevel:all",false);
75  pythia.settings.flag("PartonLevel:MPI",false);
76 
77  // Switch on cross section estimation procedure.
78  pythia.settings.flag("Merging:doXSectionEstimate", true);
79 
80  int njetcounterLO = nMaxLO;
81  string iPathTree = iPath + "_tree";
82 
83  // Save estimates in vectors.
84  vector<double> xsecLO;
85  vector<double> nAcceptLO;
86 
87  cout << endl << endl << endl;
88  cout << "Start estimating ckkwl tree level cross section" << endl;
89 
90  while(njetcounterLO >= 0) {
91 
92  // From njetcounter, choose LHE file
93  stringstream in;
94  in << "_" << njetcounterLO << ".lhe";
95 #ifdef GZIPSUPPORT
96  if(access( (iPathTree+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
97 #endif
98  string LHEfile = iPathTree + in.str();
99  LHAupLHEF lhareader((char*)(LHEfile).c_str());
100  pythia.settings.mode("Merging:nRequested", njetcounterLO);
101  pythia.settings.word("Beams:LHEF", LHEfile);
102  pythia.init(&lhareader);
103 
104  // Start generation loop
105  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
106  // Generate next event
107  if( !pythia.next() ) {
108  if( pythia.info.atEndOfFile() ){
109  break;
110  }
111  else continue;
112  }
113  } // end loop over events to generate
114 
115  // print cross section, errors
116  pythia.stat();
117 
118  xsecLO.push_back(pythia.info.sigmaGen());
119  nAcceptLO.push_back(pythia.info.nAccepted());
120 
121  // Restart with ME of a reduced the number of jets
122  if( njetcounterLO > 0 )
123  njetcounterLO--;
124  else
125  break;
126 
127  } // end loop over different jet multiplicities
128 
129 //-----------------------------------------------------------------------------
130 //-----------------------------------------------------------------------------
131 
132  // Switch off cross section estimation.
133  pythia.settings.flag("Merging:doXSectionEstimate", false);
134 
135  // Switch showering and multiple interaction back on.
136  pythia.settings.flag("PartonLevel:FSR",fsr);
137  pythia.settings.flag("PartonLevel:ISR",isr);
138  pythia.settings.flag("HadronLevel:all",had);
139  pythia.settings.flag("PartonLevel:MPI",mpi);
140 
141  // Declare sample cross section for output.
142  double sigmaTemp = 0.;
143  vector<double> sampleXStree;
144  // Cross section an error.
145  double sigmaTotal = 0.;
146  double errorTotal = 0.;
147 
148  int sizeLO = int(xsecLO.size());
149  njetcounterLO = nMaxLO;
150  iPathTree = iPath + "_tree";
151 
152  while(njetcounterLO >= 0){
153 
154  // From njetcounter, choose LHE file
155  stringstream in;
156  in << "_" << njetcounterLO << ".lhe";
157 #ifdef GZIPSUPPORT
158  if(access( (iPathTree+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
159 #endif
160  string LHEfile = iPathTree + in.str();
161  LHAupLHEF lhareader((char*)(LHEfile).c_str());
162 
163  cout << endl << endl << endl
164  << "Start tree level treatment for " << njetcounterLO << " jets"
165  << endl;
166 
167  pythia.settings.mode("Merging:nRequested", njetcounterLO);
168  pythia.settings.word("Beams:LHEF", LHEfile);
169  pythia.init(&lhareader);
170 
171  // Remember position in vector of cross section estimates.
172  int iNow = sizeLO-1-njetcounterLO;
173 
174  // Start generation loop
175  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
176 
177  // Generate next event
178  if( !pythia.next() ) {
179  if( pythia.info.atEndOfFile() ) break;
180  else continue;
181  }
182 
183  // Get event weight(s).
184  double weight = pythia.info.mergingWeight();
185  double evtweight = pythia.info.weight();
186  weight *= evtweight;
187 
188  // Do not print zero-weight events.
189  if ( weight == 0. ) continue;
190 
191  // Construct new empty HepMC event.
192  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
193  // Get correct cross section from previous estimate.
194  double normhepmc = xsecLO[iNow] / nAcceptLO[iNow];
195  // Set event weight
196  hepmcevt->weights().push_back(weight*normhepmc);
197  // Fill HepMC event
198  ToHepMC.fill_next_event( pythia, hepmcevt );
199  // Add the weight of the current event to the cross section.
200  sigmaTotal += weight*normhepmc;
201  sigmaTemp += weight*normhepmc;
202  errorTotal += pow2(weight*normhepmc);
203  // Report cross section to hepmc
205  xsec.set_cross_section( sigmaTotal*1e9, pythia.info.sigmaErr()*1e9 );
206  hepmcevt->set_cross_section( xsec );
207  // Write the HepMC event to file. Done with it.
208  ascii_io << hepmcevt;
209  delete hepmcevt;
210 
211  } // end loop over events to generate
212 
213  // print cross section, errors
214  pythia.stat();
215  // Save sample cross section for output.
216  sampleXStree.push_back(sigmaTemp);
217  sigmaTemp = 0.;
218 
219  // Restart with ME of a reduced the number of jets
220  if( njetcounterLO > 0 )
221  njetcounterLO--;
222  else
223  break;
224 
225  }
226 
227  // Print cross section information.
228  cout << endl << endl;
229  cout << " *---------------------------------------------------*" << endl;
230  cout << " | |" << endl;
231  cout << " | Sample cross sections after CKKW-L merging |" << endl;
232  cout << " | |" << endl;
233  cout << " | Leading order cross sections (mb): |" << endl;
234  for (int i = 0; i < int(sampleXStree.size()); ++i)
235  cout << " | " << sampleXStree.size()-1-i << "-jet: "
236  << setw(17) << scientific << setprecision(6)
237  << sampleXStree[i] << " |" << endl;
238  cout << " | |" << endl;
239  cout << " |---------------------------------------------------|" << endl;
240  cout << " |---------------------------------------------------|" << endl;
241  cout << " | Inclusive cross sections: |" << endl;
242  cout << " | |" << endl;
243  cout << " | CKKW-L merged inclusive cross section: |" << endl;
244  cout << " | " << setw(17) << scientific << setprecision(6)
245  << sigmaTotal << " +- " << setw(17) << sqrt(errorTotal) << " mb "
246  << " |" << endl;
247  cout << " | |" << endl;
248  cout << " | LO inclusive cross section: |" << endl;
249  cout << " | " << setw(17) << scientific << setprecision(6)
250  << xsecLO.back() << " mb |" << endl;
251  cout << " | |" << endl;
252  cout << " *---------------------------------------------------*" << endl;
253  cout << endl << endl;
254 
255  // Done
256  return 0;
257 
258 }
void set_cross_section(double xs, double xs_err)
Set cross section and error in pb.
The GenCrossSection class stores the generated cross section.
void set_cross_section(const GenCrossSection &)
provide a pointer to the GenCrossSection container
Definition: GenEvent.h:752
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
WeightContainer & weights()
direct access to WeightContainer
Definition: GenEvent.h:699