29 #include "fastjet/PseudoJet.hh"
30 #include "fastjet/ClusterSequence.hh"
32 using namespace Pythia8;
36 const double expCrossSec[] = { 798.0, 53.5, 6.8, 0.84, 0.074 };
47 pythia.readString(
"WeakSingleBoson:ffbar2W = on");
49 pythia.readString(
"24:onMode = off");
50 pythia.readString(
"24:onIfAny = 11 12");
52 if (doMPI ==
false) pythia.readString(
"PartonLevel:MPI = off");
55 pythia.readString(
"Beams:idB = -2212");
56 pythia.readString(
"Beams:eCM = 1960.");
60 Hist dSigma1(
"1-jet cross-section (E_jet1 > 20 GeV)", 70, 0.0, 350.0);
61 Hist dSigma2(
"2-jet cross-section (E_jet2 > 20 GeV)", 38, 0.0, 190.0);
62 Hist dSigma3(
"3-jet cross-section (E_jet3 > 20 GeV)", 16, 0.0, 80.0);
63 Hist dSigma4(
"4-jet cross-section (E_jet4 > 20 GeV)", 7, 0.0, 35.0);
64 Hist *dSigmaHist[5] = { NULL, &dSigma1, &dSigma2, &dSigma3, &dSigma4 };
65 double dSigmaBin[5] = { 0.0, 350.0 / 70.0, 190.0 / 38.0,
66 80.0 / 16.0, 35.0 / 7.0 };
70 fastjet::Strategy strategy = fastjet::Best;
71 fastjet::RecombinationScheme recombScheme = fastjet::E_scheme;
72 fastjet::JetDefinition *jetDef = NULL;
73 jetDef =
new fastjet::JetDefinition(fastjet::kt_algorithm, Rparam,
74 recombScheme, strategy);
77 std::vector <fastjet::PseudoJet> fjInputs;
80 int nEventAccept25[5] = { 0, 0, 0, 0, 0 };
81 int vetoCount[4] = { 0, 0, 0, 0 };
82 const char *vetoStr[] = {
"ET(elec)",
"|eta(elec)|",
83 "ET(missing)",
"deltaR(elec, jet)" };
84 bool firstEvent =
true;
87 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
88 if (!pythia.next())
continue;
93 for (
int i = pythia.event.size() - 1; i > 0; i--) {
94 if (pythia.event[i].idAbs() == 24) {
100 cout <<
"Error: Could not find W" << endl;
107 int daughter = pythia.event[idxElec].daughter1();
108 if (daughter == 0)
break;
109 else idxElec = daughter;
111 if (pythia.event[idxElec].idAbs() != 11 ||
112 !pythia.event[idxElec].isFinal()) {
113 cout <<
"Error: Found incorrect decay product of the W" << endl;
118 if (pythia.event[idxElec].pT() < 20.0) {
122 if (abs(pythia.event[idxElec].eta()) > 1.1) {
134 for (
int i = 0; i < pythia.event.size(); ++i) {
136 if (!pythia.event[i].isFinal())
continue;
139 if (pythia.event[i].idAbs() == 12 || pythia.event[i].idAbs() == 14 ||
140 pythia.event[i].idAbs() == 16)
continue;
143 if (fabs(pythia.event[i].eta()) > 3.6)
continue;
146 missingETvec += pythia.event[i].p();
149 if (i == idxElec)
continue;
152 fjInputs.push_back( fastjet::PseudoJet( pythia.event[i].px(),
153 pythia.event[i].py(), pythia.event[i].pz(), pythia.event[i].e() ) );
156 if (fjInputs.size() == 0) {
157 cout <<
"Error: event with no final state particles" << endl;
162 vector <fastjet::PseudoJet> inclusiveJets, sortedJets;
163 fastjet::ClusterSequence clustSeq(fjInputs, *jetDef);
167 cout <<
"Ran " << jetDef->description() << endl;
168 cout <<
"Strategy adopted by FastJet was "
169 << clustSeq.strategy_string() << endl << endl;
174 inclusiveJets = clustSeq.inclusive_jets(20.0);
175 sortedJets = sorted_by_pt(inclusiveJets);
178 double missingET = missingETvec.pT();
179 if (missingET < 30.0) {
185 int jetCount20 = 0, jetCount25 = 0;
187 bool vetoEvent =
false;
188 fastjet::PseudoJet fjElec(pythia.event[idxElec].px(),
189 pythia.event[idxElec].py(),
190 pythia.event[idxElec].pz(),
191 pythia.event[idxElec].e());
193 for (
unsigned int i = 0; i < sortedJets.size(); i++) {
195 if (fabs(sortedJets[i].rap()) > 2.0)
continue;
197 if (fjElec.squared_distance(sortedJets[i]) < 0.52 * 0.52)
198 { vetoEvent =
true;
break; }
201 if (sortedJets[i].perp() > 25.0)
205 dSigmaHist[++jetCount20]->fill(sortedJets[i].perp());
207 if (vetoEvent) { vetoCount[3]++;
continue; }
209 if (jetCount25 > 4) jetCount25 = 4;
210 for (
int i = jetCount25; i >= 0; i--)
220 double sigmapb = pythia.info.sigmaGen() * 1.0E9;
222 for (
int i = 1; i <= 4; i++)
223 (*dSigmaHist[i]) = ((*dSigmaHist[i]) * sigmapb) / nEvent / dSigmaBin[i];
224 cout << dSigma1 << dSigma2 << dSigma3 << dSigma4 << endl;
227 cout <<
"Jet algorithm is kT" << endl;
228 cout <<
"Multiparton interactions are switched "
229 << ( (doMPI) ?
"on" :
"off" ) << endl;
230 cout << endl << nEvent <<
" events generated. " << nEventAccept25[0]
231 <<
" events passed cuts." << endl;
232 cout <<
"Vetos:" << endl;
233 for (
int i = 0; i < 4; i++)
234 cout <<
" " << vetoStr[i] <<
" = " << vetoCount[i] << endl;
236 cout << endl <<
"Inclusive cross-sections (pb):" << endl;
237 for (
int i = 0; i < 5; i++) {
238 cout << scientific << setprecision(3)
239 <<
" " << i <<
"-jet - Pythia = "
240 << ((double) nEventAccept25[i] / (
double) nEvent) * sigmapb;
241 cout << ", Experimental = " << expCrossSec[i];
243 cout << scientific << setprecision(3)
244 <<
", Pythia ratio to " << i - 1 <<
"-jet = "
245 << ((double) nEventAccept25[i] / (
double) nEventAccept25[i - 1]);
246 cout << scientific << setprecision(3)
247 << ", Experimental ratio to " << i - 1 << "-jet = "
248 << expCrossSec[i] / expCrossSec[i - 1];