15 using namespace Pythia8;
20 bool redoBDecays =
false;
21 bool redoHadrons =
true;
22 if (redoHadrons) redoBDecays =
false;
30 if (!redoBDecays && !redoHadrons) nRepeat = 1;
34 Event&
event = pythia.event;
39 pythia.readString(
"HardQCD:gg2bbbar = on");
40 pythia.readString(
"HardQCD:qqbar2bbbar = on");
41 pythia.readString(
"PhaseSpace:pTHatMin = 50.");
45 int bCodes[28] = {511, 521, 531, 541, 5122, 5132, 5142, 5232, 5242,
46 5332, 5342, 5412, 5414, 5422, 5424, 5432, 5434, 5442, 5444, 5512,
47 5514, 5522, 5524, 5532, 5534, 5542, 5544, 5544 };
58 if (redoHadrons) pythia.readString(
"HadronLevel:all = off");
64 Hist nBperEvent(
"number of b quarks in an event", 10, -0.5, 9.5);
65 Hist nSameEvent(
"number of times same event is used", 10, -0.5, 9.5);
66 Hist oppSignMass(
"mass of opposite-sign muon pair", 100, 0.0, 100.0);
67 Hist sameSignMass(
"mass of same-sign muon pair", 100, 0.0, 100.0);
70 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
74 if (redoBDecays)
for (
int iC = 0; iC < nCodes; ++iC)
75 pythia.particleData.mayDecay( bCodes[iC],
false);
78 if (!pythia.next())
continue;
83 for (
int i = 0; i <
event.size(); ++i) {
84 stat =
event[i].statusAbs();
85 if (event[i].idAbs() == 5 && (stat == 62 || stat == 63)) ++nBquark;
87 nBperEvent.fill( nBquark );
92 for (
int i = 0; i <
event.size(); ++i) {
93 int idAbs =
event[i].idAbs();
94 for (
int iC = 0; iC < 28; ++iC)
95 if (idAbs == bCodes[iC]) {
102 nBHad = iBHad.size();
103 if (nBquark != nBHad) cout <<
" Warning: " << nBquark
104 <<
" b quarks but " << nBHad <<
" B hadrons" << endl;
110 for (
int iC = 0; iC < nCodes; ++iC)
111 pythia.particleData.mayDecay( bCodes[iC],
true);
114 }
else if (redoHadrons) {
120 for (
int iRepeat = 0; iRepeat < nRepeat; ++iRepeat) {
128 for (
int iB = 0; iB < nBHad; ++iB) event[ iBHad[iB] ].statusPos();
135 if (!pythia.moreDecays())
continue;
139 }
else if (redoHadrons) {
140 if (iRepeat > 0)
event = savedEvent;
144 if (!pythia.forceHadronLevel(
false))
continue;
148 if ( (redoBDecays || redoHadrons) && iEvent < nListRedo
149 && iRepeat == nRepeat - 1) event.list();
152 vector<int> iMuNeg, iMuPos;
153 for (
int i = 0; i <
event.size(); ++i) {
154 int id =
event[i].id();
155 if (
id == 13) iMuNeg.push_back(i);
156 if (
id == -13) iMuPos.push_back(i);
160 int nMuNeg = iMuNeg.size();
161 int nMuPos = iMuPos.size();
162 if (nMuNeg + nMuPos > 1) {
166 for (
int iN = 0; iN < nMuNeg; ++iN)
167 for (
int iP = 0; iP < nMuPos; ++iP)
169 (event[iMuNeg[iN]].p() +
event[iMuPos[iP]].p()).mCalc() );
172 for (
int i1 = 0; i1 < nMuNeg - 1; ++i1)
173 for (
int i2 = i1 + 1; i2 < nMuNeg; ++i2)
175 (event[iMuNeg[i1]].p() +
event[iMuNeg[i2]].p()).mCalc() );
176 for (
int i1 = 0; i1 < nMuPos - 1; ++i1)
177 for (
int i2 = i1 + 1; i2 < nMuPos; ++i2)
179 (event[iMuPos[i1]].p() +
event[iMuPos[i2]].p()).mCalc() );
186 nSameEvent.fill( nWithPair );
193 cout << nBperEvent << nSameEvent << oppSignMass << sameSignMass << endl;