20 using namespace Pythia8;
38 Event&
event = pythia.event;
41 pythia.readString(
"HardQCD:all = on");
42 pythia.readString(
"PhaseSpace:pTHatMin = 200.");
45 pythia.readString(
"Next:numberShowInfo = 0");
46 pythia.readString(
"Next:numberShowProcess = 0");
47 pythia.readString(
"Next:numberShowEvent = 0");
50 pythia.readString(
"Beams:eCM = 14000.");
54 SlowJet slowJet( power, R, pTMin, etaMax, select, massSet);
65 fastjet::JetDefinition jetDef(fastjet::genkt_algorithm, R, power);
66 std::vector <fastjet::PseudoJet> fjInputs;
69 Hist nJetsS(
"number of jets, SlowJet", 100, -0.5, 99.5);
70 Hist nJetsF(
"number of jets, FastJet", 100, -0.5, 99.5);
71 Hist nJetsD(
"number of jets, SlowJet-FastJet", 99, -49.5, 49.5);
72 Hist pTjetsS(
"pT for jets, SlowJet", 100, 0., 500.);
73 Hist pTjetsF(
"pT for jets, FastJet", 100, 0., 500.);
74 Hist pTjetsD(
"pT for jets, SlowJet - FastJet", 100, 0., 500.);
75 Hist RdistD(
"R distance SlowJet to FastJet", 100, 0., 1.);
76 Hist etaJets(
"eta for jets", 100, -5., 5.);
77 Hist phiJets(
"phi for jets", 100, -M_PI, M_PI);
78 Hist distJets(
"R distance between jets", 100, 0., 10.);
79 Hist pTdiff(
"pT difference between consecutive jets", 100, -100., 400.);
80 Hist nAna(
"multiplicity of analyzed event", 100, -0.5, 999.5);
81 Hist tGen(
"generation time as fn of multiplicity", 100, -0.5, 999.5);
82 Hist tSlow(
"SlowJet time as fn of multiplicity", 100, -0.5, 999.5);
83 Hist tFast(
"FastJet time as fn of multiplicity", 100, -0.5, 999.5);
84 Hist tSlowGen(
"SlowJet/generation time as fn of multiplicity",
86 Hist tFastGen(
"FastJet/generation time as fn of multiplicity",
90 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
91 clock_t befGen = clock();
92 if (!pythia.next())
continue;
93 clock_t aftGen = clock();
96 clock_t befSlow = clock();
97 slowJet.analyze( pythia.event );
98 clock_t aftSlow = clock();
99 if (iEvent < nListJets) slowJet.list();
102 int nSlow = slowJet.sizeJet();
103 nJetsS.fill( nSlow );
104 for (
int i = 0; i < nSlow; ++i) {
105 pTjetsS.fill( slowJet.pT(i) );
106 etaJets.fill( slowJet.y(i) );
107 phiJets.fill( slowJet.phi(i) );
111 for (
int i = 0; i < nSlow - 1; ++i)
112 for (
int j = i + 1; j < nSlow; ++j) {
113 double dY = slowJet.y(i) - slowJet.y(j);
114 double dPhi = abs( slowJet.phi(i) - slowJet.phi(j) );
115 if (dPhi > M_PI) dPhi = 2. * M_PI - dPhi;
116 double dR = sqrt( pow2(dY) + pow2(dPhi) );
121 for (
int i = 1; i < nSlow; ++i)
122 pTdiff.fill( slowJet.pT(i-1)- slowJet.pT(i) );
125 clock_t befFast = clock();
130 for (
int i = 0; i <
event.size(); ++i)
if (event[i].isFinal()) {
133 if (select > 2 && event[i].isNeutral() )
continue;
134 else if (select == 2 && !event[i].isVisible() )
continue;
135 if (etaMax < 20. && abs(event[i].eta()) > etaMax)
continue;
138 fastjet::PseudoJet particleTemp =
event[i];
141 pTemp =
event[i].p();
142 mTemp =
event[i].m();
144 mTemp = (massSet == 0 ||
event[i].id() == 22) ? 0. : 0.13957;
145 pTemp.e( sqrt(pTemp.pAbs2() + mTemp*mTemp) );
146 particleTemp.reset_momentum( pTemp.px(), pTemp.py(),
147 pTemp.pz(), pTemp.e() );
153 fjInputs.push_back( particleTemp);
158 vector <fastjet::PseudoJet> inclusiveJets, sortedJets;
159 fastjet::ClusterSequence clustSeq(fjInputs, jetDef);
160 inclusiveJets = clustSeq.inclusive_jets(pTMin);
161 sortedJets = sorted_by_pt(inclusiveJets);
162 clock_t aftFast = clock();
167 if (iEvent < nListJets) {
168 cout <<
"\n -------- FastJet jets, p = " << setw(2) << power
169 <<
" --------------------------------------------------\n\n "
170 <<
" i pT y phi mult chgmult photons"
171 <<
" hardest pT in neutral " << endl
173 <<
" constituent hadrons " << endl;
174 for (
int i = 0; i < int(sortedJets.size()); ++i) {
175 vector<fastjet::PseudoJet> constituents
176 = sortedJets[i].constituents();
177 fastjet::PseudoJet hardest
178 = fastjet::SelectorNHardest(1)(constituents)[0];
179 vector<fastjet::PseudoJet> neutral_hadrons
180 = ( fastjet::SelectorIsHadron()
181 && fastjet::SelectorIsNeutral())(constituents);
182 double neutral_hadrons_pt = join(neutral_hadrons).perp();
183 cout << setw(4) << i << fixed << setprecision(3) << setw(11)
184 << sortedJets[i].perp() << setw(9) << sortedJets[i].rap()
185 << setw(9) << sortedJets[i].phi_std()
186 << setw(6) << constituents.size()
187 << setw(8) << fastjet::SelectorIsCharged().count(constituents)
188 << setw(8) << fastjet::SelectorId(22).count(constituents)
189 << setw(13) << hardest.user_info<
Particle>().name()
190 <<
" " << setw(10) << neutral_hadrons_pt << endl;
192 cout <<
"\n -------- End FastJet Listing ------------------"
193 <<
"---------------------------------" << endl;
197 int nFast = sortedJets.size();
198 nJetsF.fill( nFast );
199 for (
int i = 0; i < int(sortedJets.size()); ++i)
200 pTjetsF.fill( sortedJets[i].perp() );
203 nJetsD.fill( nSlow - nFast);
204 if (nFast == nSlow) {
205 for (
int i = 0; i < nSlow; ++i) {
206 double dist2 = pow2( slowJet.y(i) - sortedJets[i].rap())
207 + pow2( slowJet.phi(i) - sortedJets[i].phi_std());
208 RdistD.fill( sqrt(dist2) );
213 nAna.fill( nAnalyze);
214 tGen.fill( nAnalyze, aftGen - befGen);
215 tSlow.fill( nAnalyze, aftSlow - befSlow);
216 tFast.fill( nAnalyze, aftFast - befFast);
223 pTjetsD = pTjetsS - pTjetsF;
224 tSlowGen = tSlow / tGen;
225 tFastGen = tFast / tGen;
230 cout << nJetsS << nJetsF << nJetsD << pTjetsS << pTjetsF << pTjetsD
231 << RdistD << etaJets << phiJets << distJets << pTdiff
232 << nAna << tGen << tSlow << tFast << tSlowGen << tFastGen;