StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main15.cc
1 // main15.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 how either
8 // (a) B decays (sections marked "Repeated decays:), or
9 // (b) all hadronization (sections marked "Repeated hadronization:")
10 // could be repeated a number of times for each event,
11 // to improve statistics when this could be a problem.
12 // Option (a) is faster than (b), but less generic.
13 
14 // Note 1: the compartmentalization of hadronization in forceHadronLevel
15 // from the rest of the event processing somewhat limits the ways the
16 // program can retry in case of problems, and so an occasional abort
17 // may occur more easily than normally.
18 
19 // Note 2: for simple cases, where it is only one particle that is to be
20 // decayed repeatedly, the event[i].undoDecay() method is handy.
21 // When used for several particles, remember that the position of
22 // some particles may be moved by the undoDecay step.
23 
24 #include "Pythia8/Pythia.h"
25 using namespace Pythia8;
26 
27 int main() {
28 
29  // Main switches: redo B decays only or redo all hadronization, but not both.
30  bool redoBDecays = false;
31  bool redoHadrons = true;
32  if (redoHadrons) redoBDecays = false;
33 
34  // Number of events. Number to list redone events.
35  int nEvent = 100;
36  int nListRedo = 1;
37 
38  // Number of times decays/hadronization should be redone for each event.
39  int nRepeat = 10;
40  if (!redoBDecays && !redoHadrons) nRepeat = 1;
41 
42  // Generator. Shorthand for event.
43  Pythia pythia;
44  Event& event = pythia.event;
45 
46  // Simulate b production above given pTmin scale.
47  // Warning: these processes do not catch all possible production modes.
48  // You would need to use HardQCD:all or even SoftQCD:nonDiffractive for that.
49  pythia.readString("HardQCD:gg2bbbar = on");
50  pythia.readString("HardQCD:qqbar2bbbar = on");
51  pythia.readString("PhaseSpace:pTHatMin = 50.");
52 
53  // Repeated decays: list of weakly decaying B hadrons.
54  // Note: this list is overkill; some will never be produced.
55  int bCodes[28] = {511, 521, 531, 541, 5122, 5132, 5142, 5232, 5242,
56  5332, 5342, 5412, 5414, 5422, 5424, 5432, 5434, 5442, 5444, 5512,
57  5514, 5522, 5524, 5532, 5534, 5542, 5544, 5544 };
58  int nCodes = 28;
59 
60  // Repeated decays: location of B handrons.
61  vector<int> iBHad;
62  int nBHad = 0;
63 
64  // Repeated hadronization: spare copy of event.
65  Event savedEvent;
66 
67  // Repeated hadronization: switch off normal HadronLevel call.
68  if (redoHadrons) pythia.readString("HadronLevel:all = off");
69 
70  // Initialize for LHC energies; default 14 TeV
71  pythia.init();
72 
73  // Histogram invariant mass of muon pairs.
74  Hist nBperEvent("number of b quarks in an event", 10, -0.5, 9.5);
75  Hist nSameEvent("number of times same event is used", 10, -0.5, 9.5);
76  Hist oppSignMass("mass of opposite-sign muon pair", 100, 0.0, 100.0);
77  Hist sameSignMass("mass of same-sign muon pair", 100, 0.0, 100.0);
78 
79  // Begin event loop.
80  for (int iEvent = 0; iEvent < nEvent; ++iEvent) {
81 
82  // Repeated decays: switch off decays of weakly decaying B hadrons.
83  // (More compact solution than repeated readString(..).)
84  if (redoBDecays) for (int iC = 0; iC < nCodes; ++iC)
85  pythia.particleData.mayDecay( bCodes[iC], false);
86 
87  // Generate event. Skip it if error.
88  if (!pythia.next()) continue;
89 
90  // Find and histogram number of b quarks.
91  int nBquark = 0;
92  int stat;
93  for (int i = 0; i < event.size(); ++i) {
94  stat = event[i].statusAbs();
95  if (event[i].idAbs() == 5 && (stat == 62 || stat == 63)) ++nBquark;
96  }
97  nBperEvent.fill( nBquark );
98 
99  // Repeated decays: find all locations where B hadrons are stored.
100  if (redoBDecays) {
101  iBHad.resize(0);
102  for (int i = 0; i < event.size(); ++i) {
103  int idAbs = event[i].idAbs();
104  for (int iC = 0; iC < 28; ++iC)
105  if (idAbs == bCodes[iC]) {
106  iBHad.push_back(i);
107  break;
108  }
109  }
110 
111  // Repeated decays: check that #b = #B.
112  nBHad = iBHad.size();
113  if (nBquark != nBHad) cout << " Warning: " << nBquark
114  << " b quarks but " << nBHad << " B hadrons" << endl;
115 
116  // Repeated decays: store size of current event.
117  event.saveSize();
118 
119  // Repeated decays: switch back on weakly decaying B hadrons.
120  for (int iC = 0; iC < nCodes; ++iC)
121  pythia.particleData.mayDecay( bCodes[iC], true);
122 
123  // Repeated hadronization: copy event into spare position.
124  } else if (redoHadrons) {
125  savedEvent = event;
126  }
127 
128  // Begin loop over rounds of decays / hadronization for same event.
129  int nWithPair = 0;
130  for (int iRepeat = 0; iRepeat < nRepeat; ++iRepeat) {
131 
132  // Repeated decays: remove B decay products from previous round.
133  if (redoBDecays) {
134  if (iRepeat > 0) {
135  event.restoreSize();
136 
137  // Repeated decays: mark decayed B hadrons as undecayed.
138  for (int iB = 0; iB < nBHad; ++iB) event[ iBHad[iB] ].statusPos();
139  }
140 
141  // Repeated decays: do decays of B hadrons, sequentially for products.
142  // Note: modeDecays does not work for bottomonium (or heavier) states,
143  // since there decays like Upsilon -> g g g also need hadronization.
144  // Also, there is no provision for Bose-Einstein effects.
145  if (!pythia.moreDecays()) continue;
146 
147 
148  // Repeated hadronization: restore saved event record.
149  } else if (redoHadrons) {
150  if (iRepeat > 0) event = savedEvent;
151 
152  // Repeated hadronization: do HadronLevel (repeatedly).
153  // Note: argument false needed owing to bug in junction search??
154  if (!pythia.forceHadronLevel(false)) continue;
155  }
156 
157  // List last repetition of first few events.
158  if ( (redoBDecays || redoHadrons) && iEvent < nListRedo
159  && iRepeat == nRepeat - 1) event.list();
160 
161  // Look for muons among decay products (also from charm/tau/...).
162  vector<int> iMuNeg, iMuPos;
163  for (int i = 0; i < event.size(); ++i) {
164  int id = event[i].id();
165  if (id == 13) iMuNeg.push_back(i);
166  if (id == -13) iMuPos.push_back(i);
167  }
168 
169  // Check whether pair(s) present.
170  int nMuNeg = iMuNeg.size();
171  int nMuPos = iMuPos.size();
172  if (nMuNeg + nMuPos > 1) {
173  ++nWithPair;
174 
175  // Fill masses of opposite-sign pairs.
176  for (int iN = 0; iN < nMuNeg; ++iN)
177  for (int iP = 0; iP < nMuPos; ++iP)
178  oppSignMass.fill(
179  (event[iMuNeg[iN]].p() + event[iMuPos[iP]].p()).mCalc() );
180 
181  // Fill masses of same-sign pairs.
182  for (int i1 = 0; i1 < nMuNeg - 1; ++i1)
183  for (int i2 = i1 + 1; i2 < nMuNeg; ++i2)
184  sameSignMass.fill(
185  (event[iMuNeg[i1]].p() + event[iMuNeg[i2]].p()).mCalc() );
186  for (int i1 = 0; i1 < nMuPos - 1; ++i1)
187  for (int i2 = i1 + 1; i2 < nMuPos; ++i2)
188  sameSignMass.fill(
189  (event[iMuPos[i1]].p() + event[iMuPos[i2]].p()).mCalc() );
190 
191  // Finished analysis of current round.
192  }
193 
194  // End of loop over many rounds. fill number of rounds with pairs.
195  }
196  nSameEvent.fill( nWithPair );
197 
198  // End of event loop.
199  }
200 
201  // Statistics. Histograms.
202  pythia.stat();
203  cout << nBperEvent << nSameEvent << oppSignMass << sameSignMass << endl;
204 
205  // Done.
206  return 0;
207 }