StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main08.cc
1 // main08.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 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 is a simple test program.
7 // It illustrates how to combine subruns in pT bins.
8 
9 #include "Pythia.h"
10 
11 using namespace Pythia8;
12 
13 int main() {
14 
15  // Two different modes are illustrated for setting the pT ranges.
16  // 1 : Hardcoded in the main program.
17  // 2 : Using the Main:subrun keyword in a separate command file.
18  int mode = 2;
19 
20  // Number of events to generate per bin, and to list.
21  int nEvent = 10000;
22 
23  // Book histograms.
24  Hist pTraw("pTHat distribution, unweighted", 100, 0., 1000.);
25  Hist pTnorm("pTHat distribution, weighted", 100, 0., 1000.);
26  Hist pTpow3("pTHat distribution, pT3*weighted", 100, 0., 1000.);
27  Hist pTpow5("pTHat distribution, pT5*weighted", 100, 0., 1000.);
28  Hist pTnormPart("pTHat distribution, weighted", 100, 0., 1000.);
29  Hist pTpow3Part("pTHat distribution, pT3*weighted", 100, 0., 1000.);
30  Hist pTpow5Part("pTHat distribution, pT5*weighted", 100, 0., 1000.);
31 
32  // Generator.
33  Pythia pythia;
34 
35  // Shorthand for some public members of pythia (also static ones).
36  Settings& settings = pythia.settings;
37  Info& info = pythia.info;
38 
39  // Set up to generate QCD jets, but only the hard process itself.
40  pythia.readString("HardQCD:all = on");
41  pythia.readString("PartonLevel:all = off");
42 
43  // Number of bins to use. In mode 2 read from main08.cmnd file.
44  int nBin = 5;
45  if (mode == 2) {
46  pythia.readFile("main08.cmnd");
47  nBin = pythia.mode("Main:numberOfSubruns");
48  }
49 
50  // Mode 1: set up five pT bins - last one open-ended.
51  double pTlimit[6] = {100., 150., 250., 400., 600., 0.};
52 
53  // Loop over number of bins, i.e. number of subruns.
54  for (int iBin = 0; iBin < nBin; ++iBin) {
55 
56  // Mode 1: hardcoded here. Using settings.parm allows non-string input.
57  if (mode <= 1) {
58  settings.parm("PhaseSpace:pTHatMin", pTlimit[iBin]);
59  settings.parm("PhaseSpace:pTHatMax", pTlimit[iBin + 1]);
60  }
61 
62  // Mode 2: subruns stored in the main08.cmnd file.
63  else pythia.readFile("main08.cmnd", iBin);
64 
65  // Initialize for LHC at 14 TeV.
66  pythia.readString("Beams:eCM = 14000.");
67  pythia.init();
68 
69  // List changed data.
70  settings.listChanged();
71 
72  // Reset local histograms (that need to be rescaled before added).
73  pTnormPart.null();
74  pTpow3Part.null();
75  pTpow5Part.null();
76 
77  // Begin event loop.
78  for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
79 
80  // Generate events. Quit if failure.
81  if (!pythia.next()) break;
82 
83  // Fill hard scale of event.
84  double pTHat = info. pTHat();
85  pTraw.fill( pTHat );
86  pTnormPart.fill( pTHat );
87  pTpow3Part.fill( pTHat, pow3(pTHat) );
88  pTpow5Part.fill( pTHat, pow5(pTHat) );
89 
90  // End of event loop. Statistics.
91  }
92  pythia.stat();
93 
94  // Normalize each case to cross section/(bin * event), and add to sum.
95  double sigmaNorm = info.sigmaGen() / (10. * nEvent);
96  pTnorm += sigmaNorm * pTnormPart;
97  pTpow3 += sigmaNorm * pTpow3Part;
98  pTpow5 += sigmaNorm * pTpow5Part;
99 
100  // End of pT-bin loop.
101  }
102 
103  // Output histograms.
104  cout << pTraw << pTnorm << pTpow3 << pTpow5;
105 
106  // Done.
107  return 0;
108 }