]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PYTHIA8/pythia8145/examples/main08.cc
Use Output directive instead of the old OutputFile and OUtputArchive. Save fileinfo...
[u/mrichter/AliRoot.git] / PYTHIA8 / pythia8145 / examples / main08.cc
1 // main08.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2010 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   int nList = 1;
23
24   // Book histograms.
25   Hist pTraw("pTHat distribution, unweighted", 100, 0., 1000.);
26   Hist pTnorm("pTHat distribution, weighted", 100, 0., 1000.);
27   Hist pTpow3("pTHat distribution, pT3*weighted", 100, 0., 1000.);
28   Hist pTpow5("pTHat distribution, pT5*weighted", 100, 0., 1000.);
29   Hist pTnormPart("pTHat distribution, weighted", 100, 0., 1000.);
30   Hist pTpow3Part("pTHat distribution, pT3*weighted", 100, 0., 1000.);
31   Hist pTpow5Part("pTHat distribution, pT5*weighted", 100, 0., 1000.);
32
33   // Generator.
34   Pythia pythia;
35
36   // Shorthand for some public members of pythia (also static ones).
37   Settings& settings = pythia.settings;
38   Info& info = pythia.info;
39
40   // Set up to generate QCD jets.
41   pythia.readString("HardQCD:all = on");  
42   pythia.readString("PartonLevel:all = off");  
43
44   // Number of bins to use. In mode 2 read from main08.cmnd file.
45   int nBin = 5;
46   if (mode == 2) {
47     pythia.readFile("main08.cmnd"); 
48     nBin = pythia.mode("Main:numberOfSubruns");
49   } 
50
51   // Mode 1: set up five pT bins - last one open-ended.
52   double pTlimit[6] = {100., 150., 250., 400., 600., 0.};
53
54   // Loop over number of bins, i.e. number of subruns.
55   for (int iBin = 0; iBin < nBin; ++iBin) {
56
57      // Mode 1: hardcoded here. Using settings.parm allows non-string input.  
58      if (mode <= 1) { 
59        settings.parm("PhaseSpace:pTHatMin", pTlimit[iBin]);  
60        settings.parm("PhaseSpace:pTHatMax", pTlimit[iBin + 1]);
61      }
62
63      // Mode 2: subruns stored in the main08.cmnd file.
64      else pythia.readFile("main08.cmnd", iBin);  
65
66      // Initialize for LHC.
67      pythia.init( 2212, 2212, 14000.);
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       // List first few events, only hard process.
84       if (iEvent < nList) pythia.process.list();
85
86       // Fill hard scale of event.
87       double pTHat = info. pTHat();
88       pTraw.fill( pTHat );
89       pTnormPart.fill( pTHat );
90       pTpow3Part.fill( pTHat, pow3(pTHat) );
91       pTpow5Part.fill( pTHat, pow5(pTHat) );
92
93
94     // End of event loop. Statistics.
95     }
96     pythia.statistics();
97
98     // Normalize each case to cross section/(bin * event), and add to sum.
99     double sigmaNorm = info.sigmaGen() / (10. * nEvent);
100     pTnorm += sigmaNorm * pTnormPart;
101     pTpow3 += sigmaNorm * pTpow3Part;
102     pTpow5 += sigmaNorm * pTpow5Part;
103
104   // End of pT-bin loop.
105   }
106
107   // Output histograms.
108   cout << pTraw << pTnorm << pTpow3 << pTpow5; 
109
110   // Done.
111   return 0;
112 }