StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main16.cc
1 // main16.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 // This is a simple test program.
7 // It illustrates (a) how to collect the analysis code in a separate class
8 // and (b) how to provide the .cmnd filename on the command line
9 
10 // Once you have linked the main program you can run it with a command line
11 // ./main16.exe main16.cmnd > out16
12 
13 #include "Pythia8/Pythia.h"
14 
15 using namespace Pythia8;
16 
17 //==========================================================================
18 
19 // Put all your own analysis code in the myAnalysis class.
20 
21 class MyAnalysis {
22 
23 public:
24 
25  // Constructor can be empty.
26  MyAnalysis() {}
27 
28  // Initialization actions.
29  void init();
30 
31  // Analysis of each new event.
32  void analyze(Event& event);
33 
34  // Show final results.
35  void finish();
36 
37 private:
38 
39  // Declare variables and objects that span init - analyze - finish.
40  int nEvt;
41  Hist brH, yH, etaChg, mult;
42 
43 };
44 
45 //--------------------------------------------------------------------------
46 
47 // The initialization code.
48 
49 void MyAnalysis::init() {
50 
51  // Initialize counter for number of events.
52  nEvt = 0;
53 
54  // Book histograms.
55  brH.book("Higgs branching ratios by flavour", 30, -0.5, 29.5);
56  yH.book("Higgs rapidity", 100, -10., 10.);
57  etaChg.book("charged pseudorapidity", 100, -10., 10.);
58  mult.book( "charged multiplicity", 100, -0.5, 799.5);
59 
60 }
61 
62 //--------------------------------------------------------------------------
63 
64 // The event analysis code.
65 
66 void MyAnalysis::analyze(Event& event) {
67 
68  // Increase counter.
69  ++nEvt;
70 
71  // Find latest copy of Higgs and plot its rapidity.
72  int iH = 0;
73  for (int i = 0; i < event.size(); ++i)
74  if (event[i].id() == 25) iH = i;
75  yH.fill( event[iH].y() );
76 
77  // Plot flavour of decay channel.
78  int idDau1 = event[ event[iH].daughter1() ].idAbs();
79  int idDau2 = event[ event[iH].daughter2() ].idAbs();
80  int iChan = 29;
81  if (idDau2 == idDau1 && idDau1 < 25) iChan = idDau1;
82  if (min( idDau1, idDau2) == 22 && max( idDau1, idDau2) == 23) iChan = 26;
83  brH.fill( iChan);
84 
85  // Plot pseudorapidity distribution. Sum up charged multiplicity.
86  int nChg = 0;
87  for (int i = 0; i < event.size(); ++i)
88  if (event[i].isFinal() && event[i].isCharged()) {
89  etaChg.fill( event[i].eta() );
90  ++nChg;
91  }
92  mult.fill( nChg );
93 
94 }
95 
96 //--------------------------------------------------------------------------
97 
98 // The finishing code.
99 
100 void MyAnalysis::finish() {
101 
102  // Normalize histograms.
103  double binFactor = 5. / nEvt;
104  yH *= binFactor;
105  etaChg *= binFactor;
106 
107  // Print histograms.
108  cout << brH << yH << etaChg << mult;
109 
110 }
111 
112 //==========================================================================
113 
114 // You should not need to touch the main program: its actions are
115 // determined by the .cmnd file and the rest belongs in MyAnalysis.
116 
117 int main(int argc, char* argv[]) {
118 
119  // Check that correct number of command-line arguments
120  if (argc != 2) {
121  cerr << " Unexpected number of command-line arguments. \n"
122  << " You are expected to provide a file name and nothing else. \n"
123  << " Program stopped! " << endl;
124  return 1;
125  }
126 
127  // Check that the provided file name corresponds to an existing file.
128  ifstream is(argv[1]);
129  if (!is) {
130  cerr << " Command-line file " << argv[1] << " was not found. \n"
131  << " Program stopped! " << endl;
132  return 1;
133  }
134 
135  // Confirm that external file will be used for settings..
136  cout << " PYTHIA settings will be read from file " << argv[1] << endl;
137 
138  // Declare generator. Read in commands from external file.
139  Pythia pythia;
140  pythia.readFile(argv[1]);
141 
142  // Initialization.
143  pythia.init();
144 
145  // Declare user analysis class. Do initialization part of it.
146  MyAnalysis myAnalysis;
147  myAnalysis.init();
148 
149  // Read in number of event and maximal number of aborts.
150  int nEvent = pythia.mode("Main:numberOfEvents");
151  int nAbort = pythia.mode("Main:timesAllowErrors");
152  bool hasPL = pythia.flag("PartonLevel:all");
153 
154  // Begin event loop.
155  int iAbort = 0;
156  for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
157 
158  // Generate events. Quit if too many failures.
159  if (!pythia.next()) {
160  if (++iAbort < nAbort) continue;
161  cout << " Event generation aborted prematurely, owing to error!\n";
162  break;
163  }
164 
165  // User Analysis of current event.
166  myAnalysis.analyze( (hasPL ? pythia.event : pythia.process) );
167 
168  // End of event loop.
169  }
170 
171  // Final statistics.
172  pythia.stat();
173 
174  // User finishing.
175  myAnalysis.finish();
176 
177  // Done.
178  return 0;
179 }