StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main88.cc
1 // main88.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 UNLOPS merging,
8 // see the NLO 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 UNLOPS 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  << " or '_powheg' identifiers" << endl
34  << " 3. Path for output histogram files" << endl
35  << " Program stopped. " << endl;
36  return 1;
37  }
38 
39  Pythia pythia;
40 
41  // Input parameters:
42  // 1. Input file for settings
43  // 2. Path to input LHE file
44  // 3. OUtput histogram path
45  pythia.readFile(argv[1]);
46 
47  // Interface for conversion from Pythia8::Event to HepMC one.
48  HepMC::Pythia8ToHepMC ToHepMC;
49  // Specify file where HepMC events will be stored.
50  HepMC::IO_GenEvent ascii_io(argv[3], std::ios::out);
51  // Switch off warnings for parton-level events.
52  ToHepMC.set_print_inconsistency(false);
53  ToHepMC.set_free_parton_warnings(false);
54  // Do not store cross section information, as this will be done manually.
55  ToHepMC.set_store_pdf(false);
56  ToHepMC.set_store_proc(false);
57  ToHepMC.set_store_xsec(false);
58 
59  // Path to input events, with name up to the "_tree", "_powheg" identifier
60  // included.
61  string iPath = string(argv[2]);
62 
63  // Number of events
64  int nEvent = pythia.mode("Main:numberOfEvents");
65  // Maximal number of additional LO jets.
66  int nMaxLO = pythia.mode("Merging:nJetMax");
67  // maximal number of additional NLO jets.
68  int nMaxNLO = pythia.mode("Merging:nJetMaxNLO");
69 
70 //----------------------------------------------------------------------------
71 //----------------------------------------------------------------------------
72 
73  // Switch off all showering and MPI when extimating the cross section after
74  // the merging scale cut.
75  bool fsr = pythia.flag("PartonLevel:FSR");
76  bool isr = pythia.flag("PartonLevel:ISR");
77  bool mpi = pythia.flag("PartonLevel:MPI");
78  bool had = pythia.flag("HadronLevel:all");
79  pythia.settings.flag("PartonLevel:FSR",false);
80  pythia.settings.flag("PartonLevel:ISR",false);
81  pythia.settings.flag("HadronLevel:all",false);
82  pythia.settings.flag("PartonLevel:MPI",false);
83 
84  // Switch on cross section estimation procedure.
85  pythia.settings.flag("Merging:doXSectionEstimate", true);
86  pythia.settings.flag("Merging:doUNLOPSTree", true);
87 
88  int njetcounterLO = nMaxLO;
89  string iPathTree = iPath + "_tree";
90 
91  // Save estimates in vectors.
92  vector<double> xsecLO;
93  vector<double> nSelectedLO;
94  vector<double> nAcceptLO;
95  vector<int> strategyLO;
96 
97  cout << endl << endl << endl;
98  cout << "Start estimating unlops tree level cross section" << endl;
99 
100  while(njetcounterLO >= 0){
101 
102  // From njetcounter, choose LHE file
103  stringstream in;
104  in << "_" << njetcounterLO << ".lhe";
105 #ifdef GZIPSUPPORT
106  if(access( (iPathTree+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
107 #endif
108  string LHEfile = iPathTree + in.str();
109  LHAupLHEF lhareader((char*)(LHEfile).c_str());
110  pythia.settings.mode("Merging:nRequested", njetcounterLO);
111  pythia.settings.word("Beams:LHEF", LHEfile);
112  pythia.init(&lhareader);
113 
114  // Start generation loop
115  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
116  // Generate next event
117  if( !pythia.next() ) {
118  if( pythia.info.atEndOfFile() ){
119  break;
120  }
121  else continue;
122  }
123  } // end loop over events to generate
124 
125  // print cross section, errors
126  pythia.stat();
127 
128  xsecLO.push_back(pythia.info.sigmaGen());
129  nSelectedLO.push_back(pythia.info.nSelected());
130  nAcceptLO.push_back(pythia.info.nAccepted());
131  strategyLO.push_back(pythia.info.lhaStrategy());
132 
133  // Restart with ME of a reduced the number of jets
134  if( njetcounterLO > 0 )
135  njetcounterLO--;
136  else
137  break;
138 
139  }
140 
141 //-----------------------------------------------------------------------------
142 //-----------------------------------------------------------------------------
143 
144  cout << endl << endl << endl;
145  cout << "Start estimating unlops virtual corrections cross section" << endl;
146 
147  pythia.settings.flag("Merging:doUNLOPSTree",false);
148  pythia.settings.flag("Merging:doUNLOPSLoop", true);
149 
150  int njetcounterNLO = nMaxNLO;
151  string iPathLoop = iPath + "_powheg";
152 
153  // Save estimates in vectors.
154  vector<double> xsecNLO;
155  vector<double> nSelectedNLO;
156  vector<double> nAcceptNLO;
157  vector<int> strategyNLO;
158 
159  while(njetcounterNLO >= 0){
160 
161  // From njetcounter, choose LHE file
162  stringstream in;
163  in << "_" << njetcounterNLO << ".lhe";
164 #ifdef GZIPSUPPORT
165  if(access( (iPathLoop+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
166 #endif
167  string LHEfile = iPathLoop + in.str();
168  LHAupLHEF lhareader((char*)(LHEfile).c_str());
169  pythia.settings.mode("Merging:nRequested", njetcounterNLO);
170  pythia.settings.word("Beams:LHEF", LHEfile);
171  pythia.init(&lhareader);
172 
173  // Start generation loop
174  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
175  // Generate next event
176  if( !pythia.next() ) {
177  if( pythia.info.atEndOfFile() ){
178  break;
179  }
180  else continue;
181  }
182  } // end loop over events to generate
183 
184  // print cross section, errors
185  pythia.stat();
186 
187  xsecNLO.push_back(pythia.info.sigmaGen());
188  nSelectedNLO.push_back(pythia.info.nSelected());
189  nAcceptNLO.push_back(pythia.info.nAccepted());
190  strategyNLO.push_back(pythia.info.lhaStrategy());
191 
192  // Restart with ME of a reduced the number of jets
193  if( njetcounterNLO > 0 )
194  njetcounterNLO--;
195  else
196  break;
197 
198  }
199 
200  int sizeLO = int(xsecLO.size());
201  int sizeNLO = int(xsecNLO.size());
202 
203  // Switch off cross section estimation.
204  pythia.settings.flag("Merging:doXSectionEstimate", false);
205 
206  // Switch showering and multiple interaction back on.
207  pythia.settings.flag("PartonLevel:FSR",fsr);
208  pythia.settings.flag("PartonLevel:ISR",isr);
209  pythia.settings.flag("HadronLevel:all",had);
210  pythia.settings.flag("PartonLevel:MPI",mpi);
211 
212 //-----------------------------------------------------------------------------
213 //-----------------------------------------------------------------------------
214 
215  // Declare sample cross section for output.
216  double sigmaTemp = 0.;
217  vector<double> sampleXStree;
218  vector<double> sampleXSvirt;
219  vector<double> sampleXSsubtTree;
220  vector<double> sampleXSsubtVirt;
221  // Cross section an error.
222  double sigmaTotal = 0.;
223  double errorTotal = 0.;
224 
225  // Switch on tree-level processing.
226  pythia.settings.flag("Merging:doUNLOPSTree",true);
227  pythia.settings.flag("Merging:doUNLOPSLoop",false);
228  pythia.settings.flag("Merging:doUNLOPSSubt",false);
229  pythia.settings.flag("Merging:doUNLOPSSubtNLO",false);
230  pythia.settings.mode("Merging:nRecluster",0);
231 
232  // Start looping through input event files.
233  njetcounterLO = nMaxLO;
234  iPathTree = iPath + "_tree";
235 
236  while(njetcounterLO >= 0){
237 
238  // From njetcounter, choose LHE file
239  stringstream in;
240  in << "_" << njetcounterLO << ".lhe";
241 #ifdef GZIPSUPPORT
242  if(access( (iPathTree+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
243 #endif
244  string LHEfile = iPathTree + in.str();
245  LHAupLHEF lhareader((char*)(LHEfile).c_str());
246 
247  cout << endl << endl << endl
248  << "Start tree level treatment for " << njetcounterLO << " jets"
249  << endl;
250 
251  // UNLOPS does not contain a zero-jet tree-level sample.
252  if ( njetcounterLO == 0 ) break;
253  pythia.settings.mode("Merging:nRequested", njetcounterLO);
254  pythia.settings.word("Beams:LHEF", LHEfile);
255  pythia.init(&lhareader);
256 
257  // Remember position in vector of cross section estimates.
258  int iNow = sizeLO-1-njetcounterLO;
259 
260  // Start generation loop
261  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
262 
263  // Generate next event
264  if( !pythia.next() ) {
265  if( pythia.info.atEndOfFile() ) break;
266  else continue;
267  }
268 
269  // Get event weight(s).
270  double weightNLO = pythia.info.mergingWeightNLO();
271  double evtweight = pythia.info.weight();
272  weightNLO *= evtweight;
273  // Do not print zero-weight events.
274  if ( weightNLO == 0. ) continue;
275 
276  // Construct new empty HepMC event.
277  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
278  // Get correct cross section from previous estimate.
279  double normhepmc = xsecLO[iNow] / nAcceptLO[iNow];
280  // Set hepmc event weight.
281  hepmcevt->weights().push_back(weightNLO*normhepmc);
282  // Fill HepMC event.
283  ToHepMC.fill_next_event( pythia, hepmcevt );
284  // Add the weight of the current event to the cross section.
285  sigmaTotal += weightNLO*normhepmc;
286  sigmaTemp += weightNLO*normhepmc;
287  errorTotal += pow2(weightNLO*normhepmc);
288  // Report cross section to hepmc.
290  xsec.set_cross_section( sigmaTotal*1e9, pythia.info.sigmaErr()*1e9 );
291  hepmcevt->set_cross_section( xsec );
292  // Write the HepMC event to file. Done with it.
293  ascii_io << hepmcevt;
294  delete hepmcevt;
295  } // end loop over events to generate
296 
297  // print cross section, errors
298  pythia.stat();
299 
300  // Restart with ME of a reduced the number of jets
301  if( njetcounterLO > 0 )
302  njetcounterLO--;
303  else
304  break;
305 
306  // Save sample cross section for output.
307  sampleXStree.push_back(sigmaTemp);
308  sigmaTemp = 0.;
309 
310  }
311 
312 //-----------------------------------------------------------------------------
313 //-----------------------------------------------------------------------------
314 
315  cout << endl << endl << endl;
316  cout << "Start unlops virtual corrections part" << endl;
317 
318  // Switch on loop-level processing.
319  pythia.settings.flag("Merging:doUNLOPSTree",false);
320  pythia.settings.flag("Merging:doUNLOPSLoop",true);
321  pythia.settings.flag("Merging:doUNLOPSSubt",false);
322  pythia.settings.flag("Merging:doUNLOPSSubtNLO",false);
323  pythia.settings.mode("Merging:nRecluster",0);
324 
325  njetcounterNLO = nMaxNLO;
326  iPathLoop= iPath + "_powheg";
327 
328  while(njetcounterNLO >= 0){
329 
330  // From njetcounter, choose LHE file
331  stringstream in;
332  in << "_" << njetcounterNLO << ".lhe";
333 #ifdef GZIPSUPPORT
334  if(access( (iPathLoop+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
335 #endif
336  string LHEfile = iPathLoop + in.str();
337  LHAupLHEF lhareader((char*)(LHEfile).c_str());
338 
339  cout << endl << endl << endl
340  << "Start loop level treatment for " << njetcounterNLO << " jets"
341  << endl;
342 
343  pythia.settings.mode("Merging:nRequested", njetcounterNLO);
344  pythia.settings.word("Beams:LHEF", LHEfile);
345  pythia.init(&lhareader);
346 
347  // Remember position in vector of cross section estimates.
348  int iNow = sizeNLO-1-njetcounterNLO;
349 
350  // Start generation loop
351  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
352 
353  // Generate next event
354  if( !pythia.next() ) {
355  if( pythia.info.atEndOfFile() ) break;
356  else continue;
357  }
358 
359  // Get event weight(s).
360  double weightNLO = pythia.info.mergingWeightNLO();
361  double evtweight = pythia.info.weight();
362  weightNLO *= evtweight;
363  // Do not print zero-weight events.
364  if ( weightNLO == 0. ) continue;
365 
366  // Construct new empty HepMC event.
367  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
368  // Get correct cross section from previous estimate.
369  double normhepmc = xsecNLO[iNow] / nAcceptNLO[iNow];
370  // powheg weighted events
371  if( abs(strategyNLO[iNow]) == 4)
372  normhepmc = 1. / (1e9*nSelectedNLO[iNow]);
373  // Set hepmc event weight.
374  hepmcevt->weights().push_back(weightNLO*normhepmc);
375  // Fill HepMC event.
376  ToHepMC.fill_next_event( pythia, hepmcevt );
377  // Add the weight of the current event to the cross section.
378  sigmaTotal += weightNLO*normhepmc;
379  sigmaTemp += weightNLO*normhepmc;
380  errorTotal += pow2(weightNLO*normhepmc);
381  // Report cross section to hepmc
383  xsec.set_cross_section( sigmaTotal*1e9, pythia.info.sigmaErr()*1e9 );
384  hepmcevt->set_cross_section( xsec );
385  // Write the HepMC event to file. Done with it.
386  ascii_io << hepmcevt;
387  delete hepmcevt;
388 
389  } // end loop over events to generate
390 
391  // print cross section, errors
392  pythia.stat();
393  // Save sample cross section for output.
394  sampleXSvirt.push_back(sigmaTemp);
395  sigmaTemp = 0.;
396 
397  // Restart with ME of a reduced the number of jets
398  if( njetcounterNLO > 0)
399  njetcounterNLO--;
400  else
401  break;
402 
403  }
404 
405 //-----------------------------------------------------------------------------
406 //-----------------------------------------------------------------------------
407 
408  cout << endl << endl << endl;
409  cout << "Shower subtractive events" << endl;
410 
411  // Switch on processing of counter-events.
412  pythia.settings.flag("Merging:doUNLOPSTree",false);
413  pythia.settings.flag("Merging:doUNLOPSLoop",false);
414  pythia.settings.flag("Merging:doUNLOPSSubt",true);
415  pythia.settings.flag("Merging:doUNLOPSSubtNLO",false);
416  pythia.settings.mode("Merging:nRecluster",1);
417 
418  int nMaxCT = nMaxLO;
419  int njetcounterCT = nMaxCT;
420  string iPathSubt= iPath + "_tree";
421 
422  while(njetcounterCT >= 1){
423 
424  // From njetcounter, choose LHE file
425  stringstream in;
426  in << "_" << njetcounterCT << ".lhe";
427 #ifdef GZIPSUPPORT
428  if(access( (iPathSubt+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
429 #endif
430  string LHEfile = iPathSubt + in.str();
431  LHAupLHEF lhareader((char*)(LHEfile).c_str());
432 
433  cout << endl << endl << endl
434  << "Start subtractive treatment for " << njetcounterCT << " jets"
435  << endl;
436 
437  pythia.settings.mode("Merging:nRequested", njetcounterCT);
438  pythia.settings.word("Beams:LHEF", LHEfile);
439  pythia.init(&lhareader);
440 
441  // Remember position in vector of cross section estimates.
442  int iNow = sizeLO-1-njetcounterCT;
443 
444  // Start generation loop
445  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
446 
447  // Generate next event
448  if( !pythia.next() ) {
449  if( pythia.info.atEndOfFile() ) break;
450  else continue;
451  }
452 
453  // Get event weight(s).
454  double weightNLO = pythia.info.mergingWeightNLO();
455  double evtweight = pythia.info.weight();
456  weightNLO *= evtweight;
457  // Do not print zero-weight events.
458  if ( weightNLO == 0. ) continue;
459 
460  // Construct new empty HepMC event.
461  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
462  // Get correct cross section from previous estimate.
463  double normhepmc = -1*xsecLO[iNow] / nAcceptLO[iNow];
464  // Set hepmc event weight.
465  hepmcevt->weights().push_back(weightNLO*normhepmc);
466  // Fill HepMC event.
467  ToHepMC.fill_next_event( pythia, hepmcevt );
468  // Add the weight of the current event to the cross section.
469  sigmaTotal += weightNLO*normhepmc;
470  sigmaTemp += weightNLO*normhepmc;
471  errorTotal += pow2(weightNLO*normhepmc);
472  // Report cross section to hepmc.
474  xsec.set_cross_section( sigmaTotal*1e9, pythia.info.sigmaErr()*1e9 );
475  hepmcevt->set_cross_section( xsec );
476  // Write the HepMC event to file. Done with it.
477  ascii_io << hepmcevt;
478  delete hepmcevt;
479 
480  } // end loop over events to generate
481 
482  // print cross section, errors
483  pythia.stat();
484  // Save sample cross section for output.
485  sampleXSsubtTree.push_back(sigmaTemp);
486  sigmaTemp = 0.;
487 
488  // Restart with ME of a reduced the number of jets
489  if( njetcounterCT > 1 )
490  njetcounterCT--;
491  else
492  break;
493 
494  }
495 
496 //-----------------------------------------------------------------------------
497 //-----------------------------------------------------------------------------
498 
499  cout << endl << endl << endl;
500  cout << "Shower subtractive events" << endl;
501 
502  pythia.settings.flag("Merging:doUNLOPSTree",false);
503  pythia.settings.flag("Merging:doUNLOPSLoop",false);
504  pythia.settings.flag("Merging:doUNLOPSSubt",false);
505  pythia.settings.flag("Merging:doUNLOPSSubtNLO",true);
506  pythia.settings.mode("Merging:nRecluster",1);
507 
508  nMaxCT = nMaxNLO;
509  njetcounterCT = nMaxCT;
510  iPathSubt= iPath + "_powheg";
511 
512  while(njetcounterCT >= 1){
513 
514  // From njetcounter, choose LHE file
515  stringstream in;
516  in << "_" << njetcounterCT << ".lhe";
517 #ifdef GZIPSUPPORT
518  if(access( (iPathSubt+in.str()+".gz").c_str(), F_OK) != -1) in << ".gz";
519 #endif
520  string LHEfile = iPathSubt + in.str();
521  LHAupLHEF lhareader((char*)(LHEfile).c_str());
522 
523  cout << endl << endl << endl
524  << "Start subtractive treatment for " << njetcounterCT << " nlo jets"
525  << endl;
526 
527  pythia.settings.mode("Merging:nRequested", njetcounterCT);
528  pythia.settings.word("Beams:LHEF", LHEfile);
529  pythia.init(&lhareader);
530 
531  // Remember position in vector of cross section estimates.
532  int iNow = sizeNLO-1-njetcounterCT;
533 
534  // Start generation loop
535  for( int iEvent=0; iEvent<nEvent; ++iEvent ){
536 
537  // Generate next event
538  if( !pythia.next() ) {
539  if( pythia.info.atEndOfFile() ) break;
540  else continue;
541  }
542 
543  // Get event weight(s).
544  double weightNLO = pythia.info.mergingWeightNLO();
545  double evtweight = pythia.info.weight();
546  weightNLO *= evtweight;
547  // Do not print zero-weight events.
548  if ( weightNLO == 0. ) continue;
549 
550  // Construct new empty HepMC event.
551  HepMC::GenEvent* hepmcevt = new HepMC::GenEvent();
552  // Get correct cross section from previous estimate.
553  double normhepmc = -1*xsecNLO[iNow] / nAcceptNLO[iNow];
554  // powheg weighted events
555  if( abs(strategyNLO[iNow]) == 4)
556  normhepmc = -1. / (1e9*nSelectedNLO[iNow]);
557  // Set hepmc event weight.
558  hepmcevt->weights().push_back(weightNLO*normhepmc);
559  // Fill HepMC event.
560  ToHepMC.fill_next_event( pythia, hepmcevt );
561  // Add the weight of the current event to the cross section.
562  sigmaTotal += weightNLO*normhepmc;
563  sigmaTemp += weightNLO*normhepmc;
564  errorTotal += pow2(weightNLO*normhepmc);
565  // Report cross section to hepmc.
567  xsec.set_cross_section( sigmaTotal*1e9, pythia.info.sigmaErr()*1e9 );
568  hepmcevt->set_cross_section( xsec );
569  // Write the HepMC event to file. Done with it.
570  ascii_io << hepmcevt;
571  delete hepmcevt;
572 
573  } // end loop over events to generate
574 
575  // print cross section, errors
576  pythia.stat();
577  // Save sample cross section for output.
578  sampleXSsubtVirt.push_back(sigmaTemp);
579  sigmaTemp = 0.;
580 
581  // Restart with ME of a reduced the number of jets
582  if( njetcounterCT > 1 )
583  njetcounterCT--;
584  else
585  break;
586 
587  }
588 
589  // Print cross section information.
590  cout << endl << endl;
591  cout << " *---------------------------------------------------*" << endl;
592  cout << " | |" << endl;
593  cout << " | Sample cross sections after UNLOPS merging |" << endl;
594  cout << " | |" << endl;
595  cout << " | Leading order cross sections (mb): |" << endl;
596  for (int i = 0; i < int(sampleXStree.size()); ++i)
597  cout << " | " << sampleXStree.size()-1-i+1 << "-jet: "
598  << setw(17) << scientific << setprecision(6)
599  << sampleXStree[i] << " |" << endl;
600  cout << " | (No 0-jet tree-level sample in UNLOPS) |" << endl;
601  cout << " | |" << endl;
602  cout << " | NLO order cross sections (mb): |" << endl;
603  for (int i = 0; i < int(sampleXSvirt.size()); ++i)
604  cout << " | " << sampleXSvirt.size()-1-i << "-jet: "
605  << setw(17) << scientific << setprecision(6)
606  << sampleXSvirt[i] << " |" << endl;
607  cout << " | |" << endl;
608  cout << " | Leading-order subtractive cross sections (mb): |" << endl;
609  for (int i = 0; i < int(sampleXSsubtTree.size()); ++i)
610  cout << " | " << sampleXSsubtTree.size()-1-i+1 << "-jet: "
611  << setw(17) << scientific << setprecision(6)
612  << sampleXSsubtTree[i] << " |" << endl;
613  cout << " | |" << endl;
614  if ( sampleXSsubtVirt.size() > 0) {
615  cout << " | NLO subtractive cross sections (mb): |" << endl;
616  for (int i = 0; i < int(sampleXSsubtVirt.size()); ++i)
617  cout << " | " << sampleXSsubtVirt.size()-1-i+1 << "-jet: "
618  << setw(17) << scientific << setprecision(6)
619  << sampleXSsubtVirt[i] << " |" << endl;
620  cout << " | |" << endl;
621  }
622  cout << " |---------------------------------------------------|" << endl;
623  cout << " |---------------------------------------------------|" << endl;
624  cout << " | Inclusive cross sections: |" << endl;
625  cout << " | |" << endl;
626  cout << " | UNLOPS merged inclusive cross section: |" << endl;
627  cout << " | " << setw(17) << scientific << setprecision(6)
628  << sigmaTotal << " +- " << setw(17) << sqrt(errorTotal) << " mb "
629  << " |" << endl;
630  cout << " | |" << endl;
631  cout << " | NLO inclusive cross section: |" << endl;
632  cout << " | " << setw(17) << scientific << setprecision(6)
633  << xsecNLO.back() << " mb |" << endl;
634  cout << " | |" << endl;
635  cout << " | LO inclusive cross section: |" << endl;
636  cout << " | " << setw(17) << scientific << setprecision(6)
637  << xsecLO.back() << " mb |" << endl;
638  cout << " | |" << endl;
639  cout << " *---------------------------------------------------*" << endl;
640  cout << endl << endl;
641 
642  // Done
643  return 0;
644 
645 }
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