9 #include "Pythia8/Analysis.h"
10 #include "Pythia8/FJcore.h"
25 const int Sphericity::NSTUDYMIN = 2;
28 const int Sphericity::TIMESTOPRINT = 1;
31 const double Sphericity::P2MIN = 1e-20;
34 const double Sphericity::EIGENVALUEMIN = 1e-10;
40 bool Sphericity::analyze(
const Event& event) {
43 eVal1 = eVal2 = eVal3 = 0.;
44 eVec1 = eVec2 = eVec3 = 0.;
46 for (
int j = 1; j < 4; ++j)
47 for (
int k = j; k < 4; ++k) tt[j][k] = 0.;
52 for (
int i = 0; i <
event.size(); ++i)
53 if (event[i].isFinal()) {
54 if (select > 2 && event[i].isNeutral() )
continue;
55 if (select == 2 && !event[i].isVisible() )
continue;
60 pNow[1] =
event[i].px();
61 pNow[2] =
event[i].py();
62 pNow[3] =
event[i].pz();
63 double p2Now = pNow[1]*pNow[1] + pNow[2]*pNow[2] + pNow[3]*pNow[3];
65 if (powerInt == 1) pWeight = 1. / sqrt(max(P2MIN, p2Now));
66 else if (powerInt == 0) pWeight = pow( max(P2MIN, p2Now), powerMod);
67 for (
int j = 1; j < 4; ++j)
68 for (
int k = j; k < 4; ++k) tt[j][k] += pWeight * pNow[j] * pNow[k];
69 denom += pWeight * p2Now;
73 if (nStudy < NSTUDYMIN) {
74 if (nFew < TIMESTOPRINT) cout <<
" PYTHIA Error in "
75 <<
"Sphericity::analyze: too few particles" << endl;
81 for (
int j = 1; j < 4; ++j)
82 for (
int k = j; k < 4; ++k) tt[j][k] /= denom;
85 double qCoef = ( tt[1][1] * tt[2][2] + tt[1][1] * tt[3][3]
86 + tt[2][2] * tt[3][3] - pow2(tt[1][2]) - pow2(tt[1][3])
87 - pow2(tt[2][3]) ) / 3. - 1./9.;
88 double qCoefRt = sqrt( -qCoef);
89 double rCoef = -0.5 * ( qCoef + 1./9. + tt[1][1] * pow2(tt[2][3])
90 + tt[2][2] * pow2(tt[1][3]) + tt[3][3] * pow2(tt[1][2])
91 - tt[1][1] * tt[2][2] * tt[3][3] )
92 + tt[1][2] * tt[1][3] * tt[2][3] + 1./27.;
93 double pTemp = max( min( rCoef / pow3(qCoefRt), 1.), -1.);
94 double pCoef = cos( acos(pTemp) / 3.);
95 double pCoefRt = sqrt( 3. * (1. - pow2(pCoef)) );
96 eVal1 = 1./3. + qCoefRt * max( 2. * pCoef, pCoefRt - pCoef);
97 eVal3 = 1./3. + qCoefRt * min( 2. * pCoef, -pCoefRt - pCoef);
98 eVal2 = 1. - eVal1 - eVal3;
101 for (
int iVal = 0; iVal < 2; ++iVal) {
102 double eVal = (iVal == 0) ? eVal1 : eVal3;
105 if (iVal > 0 && eVal2 < EIGENVALUEMIN) {
106 if ( abs(eVec1.pz()) > 0.5) eVec3 = Vec4( 1., 0., 0., 0.);
107 else eVec3 = Vec4( 0., 0., 1., 0.);
108 eVec3 -= dot3( eVec1, eVec3) * eVec1;
109 eVec3 /= eVec3.pAbs();
110 eVec2 = cross3( eVec1, eVec3);
116 for (
int j = 1; j < 4; ++j) {
117 dd[j][j] = tt[j][j] - eVal;
118 for (
int k = j + 1; k < 4; ++k) {
128 for (
int j = 1; j < 4; ++j)
129 for (
int k = 1; k < 4; ++k)
130 if (abs(dd[j][k]) > ddMax) {
133 ddMax = abs(dd[j][k]);
139 for (
int j = 1; j < 4; ++j)
141 double pivot = dd[j][kMax] / dd[jMax][kMax];
142 for (
int k = 1; k < 4; ++k) {
143 dd[j][k] -= pivot * dd[jMax][k];
144 if (abs(dd[j][k]) > ddMax) {
146 ddMax = abs(dd[j][k]);
152 int k1 = kMax + 1;
if (k1 > 3) k1 -= 3;
153 int k2 = kMax + 2;
if (k2 > 3) k2 -= 3;
155 eVec[k1] = -dd[jMax2][k2];
156 eVec[k2] = dd[jMax2][k1];
157 eVec[kMax] = (dd[jMax][k1] * dd[jMax2][k2]
158 - dd[jMax][k2] * dd[jMax2][k1]) / dd[jMax][kMax];
159 double length = sqrt( pow2(eVec[1]) + pow2(eVec[2])
163 if (iVal == 0) eVec1 = Vec4( eVec[1] / length,
164 eVec[2] / length, eVec[3] / length, 0.);
165 else eVec3 = Vec4( eVec[1] / length,
166 eVec[2] / length, eVec[3] / length, 0.);
170 eVec2 = cross3( eVec1, eVec3);
181 void Sphericity::list()
const {
184 cout <<
"\n -------- PYTHIA Sphericity Listing -------- \n";
185 if (powerInt !=2) cout<<
" Nonstandard momentum power = "
186 << fixed << setprecision(3) << setw(6) << power <<
"\n";
187 cout <<
"\n no lambda e_x e_y e_z \n";
190 cout << setprecision(5);
191 cout <<
" 1" << setw(11) << eVal1 << setw(11) << eVec1.px()
192 << setw(10) << eVec1.py() << setw(10) << eVec1.pz() <<
"\n";
193 cout <<
" 2" << setw(11) << eVal2 << setw(11) << eVec2.px()
194 << setw(10) << eVec2.py() << setw(10) << eVec2.pz() <<
"\n";
195 cout <<
" 3" << setw(11) << eVal3 << setw(11) << eVec3.px()
196 << setw(10) << eVec3.py() << setw(10) << eVec3.pz() <<
"\n";
199 cout <<
"\n -------- End PYTHIA Sphericity Listing ----" << endl;
215 const int Thrust::NSTUDYMIN = 2;
218 const int Thrust::TIMESTOPRINT = 1;
221 const double Thrust::MAJORMIN = 1e-10;
227 bool Thrust::analyze(
const Event& event) {
230 eVal1 = eVal2 = eVal3 = 0.;
231 eVec1 = eVec2 = eVec3 = 0.;
234 Vec4 pSum, nRef, pPart, pFull, pMax;
237 for (
int i = 0; i <
event.size(); ++i)
238 if (event[i].isFinal()) {
239 if (select > 2 && event[i].isNeutral() )
continue;
240 if (select == 2 && !event[i].isVisible() )
continue;
244 Vec4 pNow =
event[i].p();
247 pOrder.push_back(pNow);
251 if (nStudy < NSTUDYMIN) {
252 if (nFew < TIMESTOPRINT) cout <<
" PYTHIA Error in "
253 <<
"Thrust::analyze: too few particles" << endl;
259 for (
int i1 = 0; i1 < nStudy - 1; ++i1)
260 for (
int i2 = i1 + 1; i2 < nStudy; ++i2) {
261 nRef = cross3( pOrder[i1], pOrder[i2]);
266 for (
int i = 0; i < nStudy; ++i)
if (i != i1 && i != i2) {
267 if (dot3(pOrder[i], nRef) > 0.) pPart += pOrder[i];
268 else pPart -= pOrder[i];
270 for (
int j = 0; j < 4; ++j) {
271 if (j == 0) pFull = pPart + pOrder[i1] + pOrder[i2];
272 else if (j == 1) pFull = pPart + pOrder[i1] - pOrder[i2];
273 else if (j == 2) pFull = pPart - pOrder[i1] + pOrder[i2];
274 else pFull = pPart - pOrder[i1] - pOrder[i2];
275 pFull.e(pFull.pAbs());
276 if (pFull.e() > pMax.e()) pMax = pFull;
281 eVal1 = pMax.e() / pSum.e();
282 eVec1 = pMax / pMax.e();
287 for (
int i = 0; i < nStudy; ++i) {
288 pOrder[i] -= dot3( eVec1, pOrder[i]) * eVec1;
289 pOrder[i].e(pOrder[i].pAbs());
290 pAbsSum += pOrder[i].e();
294 if (pAbsSum < MAJORMIN * pSum.e()) {
295 if ( abs(eVec1.pz()) > 0.5) eVec2 = Vec4( 1., 0., 0., 0.);
296 else eVec2 = Vec4( 0., 0., 1., 0.);
297 eVec2 -= dot3( eVec1, eVec2) * eVec1;
298 eVec2 /= eVec2.pAbs();
299 eVec3 = cross3( eVec1, eVec2);
305 for (
int i1 = 0; i1 < nStudy; ++i1) {
306 nRef = cross3( pOrder[i1], eVec1);
311 for (
int i = 0; i < nStudy; ++i)
if (i != i1) {
312 if (dot3(pOrder[i], nRef) > 0.) pPart += pOrder[i];
313 else pPart -= pOrder[i];
315 pFull = pPart + pOrder[i1];
316 pFull.e(pFull.pAbs());
317 if (pFull.e() > pMax.e()) pMax = pFull;
318 pFull = pPart - pOrder[i1];
319 pFull.e(pFull.pAbs());
320 if (pFull.e() > pMax.e()) pMax = pFull;
324 eVal2 = pMax.e() / pSum.e();
325 eVec2 = pMax / pMax.e();
329 eVec3 = cross3( eVec1, eVec2);
331 for (
int i = 0; i < nStudy; ++i)
332 pAbsSum += abs( dot3(eVec3, pOrder[i]) );
333 eVal3 = pAbsSum / pSum.e();
344 void Thrust::list()
const {
347 cout <<
"\n -------- PYTHIA Thrust Listing ------------ \n"
348 <<
"\n value e_x e_y e_z \n";
351 cout << setprecision(5);
352 cout <<
" Thr" << setw(11) << eVal1 << setw(11) << eVec1.px()
353 << setw(10) << eVec1.py() << setw(10) << eVec1.pz() <<
"\n";
354 cout <<
" Maj" << setw(11) << eVal2 << setw(11) << eVec2.px()
355 << setw(10) << eVec2.py() << setw(10) << eVec2.pz() <<
"\n";
356 cout <<
" Min" << setw(11) << eVal3 << setw(11) << eVec3.px()
357 << setw(10) << eVec3.py() << setw(10) << eVec3.pz() <<
"\n";
360 cout <<
"\n -------- End PYTHIA Thrust Listing --------" << endl;
375 const double SingleClusterJet::PABSMIN = 1e-10;
381 double dist2Fun(
int measure,
const SingleClusterJet& j1,
382 const SingleClusterJet& j2) {
385 if (measure == 2)
return 2. * j1.pJet.e() * j2.pJet.e()
386 * (1. - dot3( j1.pJet, j2.pJet) / (j1.pAbs * j2.pAbs) );
389 if (measure == 3)
return 2. * pow2( min( j1.pJet.e(), j2.pJet.e() ) )
390 * (1. - dot3( j1.pJet, j2.pJet) / (j1.pAbs * j2.pAbs) );
393 return (j1.pAbs * j2.pAbs - dot3( j1.pJet, j2.pJet))
394 * 2. * j1.pAbs * j2.pAbs / pow2(j1.pAbs + j2.pAbs);
410 const int ClusterJet::TIMESTOPRINT = 1;
413 const double ClusterJet::PIMASS = 0.13957;
416 const double ClusterJet::PABSMIN = 1e-10;
419 const double ClusterJet::PRECLUSTERFRAC = 0.1;
422 const double ClusterJet::PRECLUSTERSTEP = 0.8;
428 bool ClusterJet::analyze(
const Event& event,
double yScaleIn,
429 double pTscaleIn,
int nJetMinIn,
int nJetMaxIn) {
442 for (
int i = 0; i <
event.size(); ++i)
443 if (event[i].isFinal()) {
444 if (select > 2 && event[i].isNeutral() )
continue;
445 if (select == 2 && !event[i].isVisible() )
continue;
448 Vec4 pTemp =
event[i].p();
449 if (massSet == 0 || massSet == 1) {
450 double mTemp = (massSet == 0 ||
event[i].id() == 22)
452 double eTemp = sqrt(pTemp.pAbs2() + pow2(mTemp));
455 particles.push_back( SingleClusterJet(pTemp, i) );
460 nParticles = particles.size();
461 if (nParticles < nJetMin) {
462 if (nFew < TIMESTOPRINT) cout <<
" PYTHIA Error in "
463 <<
"ClusterJet::analyze: too few particles" << endl;
469 double p2Sum = pSum.m2Calc();
470 dist2Join = max( yScale * p2Sum, pow2(pTscale));
471 dist2BigMin = 2. * max( dist2Join, p2Sum);
474 if (doPrecluster && nParticles > nJetMin + 2) {
476 if (doReassign) reassign();
480 else for (
int i = 0; i < nParticles; ++i) {
481 jets.push_back( SingleClusterJet(particles[i]) );
482 particles[i].daughter = i;
489 double dist2Min = dist2BigMin;
492 for (
int j = 0; j < int(jets.size()) - 1; ++j)
493 for (
int k = j + 1; k < int(jets.size()); ++k) {
494 double dist2 = dist2Fun( measure, jets[j], jets[k]);
495 if (dist2 < dist2Min) {
503 if ( dist2Min > dist2Join
504 && (nJetMax < nJetMin ||
int(jets.size()) <= nJetMax) )
break;
507 if (
int(jets.size()) <= nJetMin)
break;
510 jets[jMin].pJet += jets[kMin].pJet;
511 jets[jMin].pAbs = max( PABSMIN, jets[jMin].pJet.pAbs());
512 jets[jMin].multiplicity += jets[kMin].multiplicity;
513 for (
int i = 0; i < nParticles; ++i)
514 if (particles[i].daughter == kMin) particles[i].daughter = jMin;
517 distances.push_front(dist2Min);
518 if (distances.size() > 5) distances.pop_back();
521 jets[kMin] = jets.back();
523 int iEnd = jets.size();
524 for (
int i = 0; i < nParticles; ++i)
525 if (particles[i].daughter == iEnd) particles[i].daughter = kMin;
528 if (doReassign) reassign();
532 for (
int j = 0; j < int(jets.size()) - 1; ++j)
533 for (
int k =
int(jets.size()) - 1; k > j; --k)
534 if (jets[k].pJet.e() > jets[k-1].pJet.e()) {
535 swap( jets[k], jets[k-1]);
536 for (
int i = 0; i < nParticles; ++i) {
537 if (particles[i].daughter == k) particles[i].daughter = k-1;
538 else if (particles[i].daughter == k-1) particles[i].daughter = k;
550 void ClusterJet::precluster() {
553 distPre = PRECLUSTERFRAC * sqrt(dist2Join) / PRECLUSTERSTEP;
555 distPre *= PRECLUSTERSTEP;
556 dist2Pre = pow2(distPre);
557 for (
int i = 0; i < nParticles; ++i) {
558 particles[i].daughter = -1;
559 particles[i].isAssigned =
false;
565 for (
int i = 0; i < nParticles; ++i)
566 if (particles[i].pAbs < 2. * distPre) {
567 pCentral += particles[i].pJet;
568 multCentral += particles[i].multiplicity;
569 particles[i].isAssigned =
true;
571 if (pCentral.pAbs() > 2. * distPre) {
572 jets.push_back( SingleClusterJet(pCentral) );
573 jets.back().multiplicity = multCentral;
574 for (
int i = 0; i < nParticles; ++i)
575 if (particles[i].isAssigned) particles[i].daughter = 0;
582 for (
int i = 0; i < nParticles; ++i)
583 if ( !particles[i].isAssigned && particles[i].pAbs > pMax) {
585 pMax = particles[i].pAbs;
587 if (iMax == -1)
break;
593 for (
int i = 0; i < nParticles; ++i)
594 if ( !particles[i].isAssigned) {
595 double dist2 = dist2Fun( measure, particles[iMax],
597 if (dist2 < dist2Pre) {
598 pPre += particles[i].pJet;
600 particles[i].isAssigned =
true;
601 particles[i].daughter = jets.size();
604 jets.push_back( SingleClusterJet(pPre) );
605 jets.back().multiplicity = multPre;
608 if (
int(jets.size()) + nRemain < nJetMin)
break;
610 if (
int(jets.size()) >= nJetMin)
break;
619 void ClusterJet::reassign() {
622 for (
int j = 0; j < int(jets.size()); ++j) {
624 jets[j].multiplicity = 0;
628 for (
int i = 0; i < nParticles; ++i) {
629 particles[i].daughter = -1;
630 double dist2Min = dist2BigMin;
632 for (
int j = 0; j < int(jets.size()); ++j) {
633 double dist2 = dist2Fun( measure, particles[i], jets[j]);
634 if (dist2 < dist2Min) {
639 jets[jMin].pTemp += particles[i].pJet;
640 ++jets[jMin].multiplicity;
641 particles[i].daughter = jMin;
645 for (
int j = 0; j < int(jets.size()); ++j) {
646 jets[j].pJet = jets[j].pTemp;
647 jets[j].pAbs = max( PABSMIN, jets[j].pJet.pAbs());
655 for (
int j = 0; j < int(jets.size()); ++j)
656 if (jets[j].multiplicity == 0) jEmpty = j;
657 if (jEmpty == -1)
return;
661 double dist2Max = 0.;
662 for (
int i = 0; i < nParticles; ++i) {
663 int j = particles[i].daughter;
664 double dist2 = dist2Fun( measure, particles[i], jets[j]);
665 if (dist2 > dist2Max) {
672 int jSplit = particles[iSplit].daughter;
673 jets[jEmpty] = SingleClusterJet( particles[iSplit].pJet );
674 jets[jSplit].pJet -= particles[iSplit].pJet;
675 jets[jSplit].pAbs = max( PABSMIN,jets[jSplit].pJet.pAbs());
676 particles[iSplit].daughter = jEmpty;
677 --jets[jSplit].multiplicity;
686 void ClusterJet::list()
const {
689 string method = (measure == 1) ?
"Lund pT"
690 : ( (measure == 2) ?
"JADE m" :
"Durham kT" ) ;
691 cout <<
"\n -------- PYTHIA ClusterJet Listing, " << setw(9) << method
692 <<
" =" << fixed << setprecision(3) << setw(7) << sqrt(dist2Join)
693 <<
" GeV --- \n \n no mult p_x p_y p_z "
697 for (
int i = 0; i < int(jets.size()); ++i) {
698 cout << setw(4) << i << setw(6) << jets[i].multiplicity << setw(11)
699 << jets[i].pJet.px() << setw(11) << jets[i].pJet.py()
700 << setw(11) << jets[i].pJet.pz() << setw(11)
701 << jets[i].pJet.e() << setw(11) << jets[i].pJet.mCalc()
706 cout <<
"\n -------- End PYTHIA ClusterJet Listing ---------------"
707 <<
"--------" << endl;
721 const int CellJet::TIMESTOPRINT = 1;
727 bool CellJet::analyze(
const Event& event,
double eTjetMinIn,
728 double coneRadiusIn,
double eTseedIn) {
731 eTjetMin = eTjetMinIn;
732 coneRadius = coneRadiusIn;
735 vector<SingleCell> cells;
738 for (
int i = 0; i <
event.size(); ++i)
739 if (event[i].isFinal()) {
740 if (select > 2 && event[i].isNeutral() )
continue;
741 if (select == 2 && !event[i].isVisible() )
continue;
744 double etaNow =
event[i].eta();
745 if (abs(etaNow) > etaMax)
continue;
746 double phiNow =
event[i].phi();
747 double pTnow =
event[i].pT();
748 int iEtaNow = max(1, min( nEta, 1 +
int(nEta * 0.5
749 * (1. + etaNow / etaMax) ) ) );
750 int iPhiNow = max(1, min( nPhi, 1 +
int(nPhi * 0.5
751 * (1. + phiNow / M_PI) ) ) );
752 int iCell = nPhi * iEtaNow + iPhiNow;
756 for (
int j = 0; j < int(cells.size()); ++j) {
757 if (iCell == cells[j].iCell) {
759 ++cells[j].multiplicity;
760 cells[j].eTcell += pTnow;
765 double etaCell = (etaMax / nEta) * (2 * iEtaNow - 1 - nEta);
766 double phiCell = (M_PI / nPhi) * (2 * iPhiNow - 1 - nPhi);
767 cells.push_back( SingleCell( iCell, etaCell, phiCell, pTnow, 1) );
772 if (smear > 0 && rndmPtr != 0)
773 for (
int j = 0; j < int(cells.size()); ++j) {
774 double eTeConv = (smear < 2) ? 1. : cosh( cells[j].etaCell );
775 double eBef = cells[j].eTcell * eTeConv;
777 do eAft = eBef + resolution * sqrt(eBef) * rndmPtr->gauss();
778 while (eAft < 0 || eAft > upperCut * eBef);
779 cells[j].eTcell = eAft / eTeConv;
783 for (
int j = 0; j < int(cells.size()); ++j) {
784 if (cells[j].eTcell < eTseed) cells[j].canBeSeed =
false;
785 if (cells[j].eTcell < threshold) cells[j].isUsed =
true;
792 for (
int j = 0; j < int(cells.size()); ++j)
793 if (cells[j].canBeSeed && cells[j].eTcell > eTmax) {
795 eTmax = cells[j].eTcell;
799 if (eTmax < eTseed)
break;
800 double etaCenterNow = cells[jMax].etaCell;
801 double phiCenterNow = cells[jMax].phiCell;
802 double eTjetNow = 0.;
805 for (
int j = 0; j < int(cells.size()); ++j) {
806 if (cells[j].isUsed)
continue;
807 double dEta = abs( cells[j].etaCell - etaCenterNow );
808 if (dEta > coneRadius)
continue;
809 double dPhi = abs( cells[j].phiCell - phiCenterNow );
810 if (dPhi > M_PI) dPhi = 2. * M_PI - dPhi;
811 if (dPhi > coneRadius)
continue;
812 if (pow2(dEta) + pow2(dPhi) > pow2(coneRadius))
continue;
813 cells[j].isAssigned =
true;
814 eTjetNow += cells[j].eTcell;
818 if (eTjetNow < eTjetMin) {
819 cells[jMax].canBeSeed =
false;
820 for (
int j = 0; j < int(cells.size()); ++j)
821 cells[j].isAssigned =
false;
825 double etaWeightedNow = 0.;
826 double phiWeightedNow = 0.;
827 int multiplicityNow = 0;
829 for (
int j = 0; j < int(cells.size()); ++j)
830 if (cells[j].isAssigned) {
831 cells[j].canBeSeed =
false;
832 cells[j].isUsed =
true;
833 cells[j].isAssigned =
false;
834 etaWeightedNow += cells[j].eTcell * cells[j].etaCell;
835 double phiCell = cells[j].phiCell;
836 if (abs(phiCell - phiCenterNow) > M_PI)
837 phiCell += (phiCenterNow > 0.) ? 2. * M_PI : -2. * M_PI;
838 phiWeightedNow += cells[j].eTcell * phiCell;
839 multiplicityNow += cells[j].multiplicity;
840 pMassiveNow += cells[j].eTcell * Vec4(
841 cos(cells[j].phiCell), sin(cells[j].phiCell),
842 sinh(cells[j].etaCell), cosh(cells[j].etaCell) );
844 etaWeightedNow /= eTjetNow;
845 phiWeightedNow /= eTjetNow;
848 jets.push_back( SingleCellJet( eTjetNow, etaCenterNow, phiCenterNow,
849 etaWeightedNow, phiWeightedNow, multiplicityNow, pMassiveNow) );
850 for (
int i =
int(jets.size()) - 1; i > 0; --i) {
851 if (jets[i-1].eTjet > jets[i].eTjet)
break;
852 swap( jets[i-1], jets[i]);
865 void CellJet::list()
const {
868 cout <<
"\n -------- PYTHIA CellJet Listing, eTjetMin = "
869 << fixed << setprecision(3) << setw(8) << eTjetMin
870 <<
", coneRadius = " << setw(5) << coneRadius
871 <<
" ------------------------------ \n \n no "
872 <<
" eTjet etaCtr phiCtr etaWt phiWt mult p_x"
873 <<
" p_y p_z e m \n";
876 for (
int i = 0; i < int(jets.size()); ++i) {
877 cout << setw(4) << i << setw(10) << jets[i].eTjet << setw(8)
878 << jets[i].etaCenter << setw(8) << jets[i].phiCenter << setw(8)
879 << jets[i].etaWeighted << setw(8) << jets[i].phiWeighted
880 << setw(5) << jets[i].multiplicity << setw(11)
881 << jets[i].pMassive.px() << setw(11) << jets[i].pMassive.py()
882 << setw(11) << jets[i].pMassive.pz() << setw(11)
883 << jets[i].pMassive.e() << setw(11)
884 << jets[i].pMassive.mCalc() <<
"\n";
888 cout <<
"\n -------- End PYTHIA CellJet Listing ------------------"
889 <<
"-------------------------------------------------"
904 const int SlowJet::TIMESTOPRINT = 1;
907 const double SlowJet::PIMASS = 0.13957;
910 const double SlowJet::TINY = 1e-20;
916 bool SlowJet::setup(
const Event& event) {
925 double mTemp, pT2Temp, mTTemp, yTemp, phiTemp;
926 for (
int i = 0; i <
event.size(); ++i)
927 if (event[i].isFinal()) {
930 if (chargedOnly && event[i].isNeutral() )
continue;
931 else if (visibleOnly && !event[i].isVisible() )
continue;
937 if (cutInEta && abs(event[i].eta()) > etaMax)
continue;
940 pTemp =
event[i].p();
941 mTemp =
event[i].m();
943 mTemp = (massSet == 0 ||
event[i].id() == 22) ? 0. : PIMASS;
944 pTemp.e( sqrt(pTemp.pAbs2() + mTemp*mTemp) );
950 pTemp =
event[i].p();
951 mTemp =
event[i].m();
952 if ( !sjHookPtr->include( i, event, pTemp, mTemp) )
continue;
956 pT2Temp = max( TINY*TINY, pTemp.pT2());
957 mTTemp = sqrt( mTemp*mTemp + pT2Temp);
958 yTemp = (pTemp.pz() > 0)
959 ? log( max( TINY, pTemp.e() + pTemp.pz() ) / mTTemp )
960 : log( mTTemp / max( TINY, pTemp.e() - pTemp.pz() ) );
961 phiTemp = pTemp.phi();
962 clusters.push_back( SingleSlowJet(pTemp, pT2Temp, yTemp, phiTemp, i) );
964 origSize = clusters.size();
967 if (useFJcore)
return true;
973 dij.resize(clSize * (clSize - 1) / 2);
976 for (
int i = 0; i < clSize; ++i) {
977 if (isAnti) diB[i] = 1. / clusters[i].pT2;
978 else if (isKT) diB[i] = clusters[i].pT2;
982 for (
int j = 0; j < i; ++j) {
983 dPhi = abs( clusters[i].phi - clusters[j].phi );
984 if (dPhi > M_PI) dPhi = 2. * M_PI - dPhi;
985 dijTemp = (useStandardR)
986 ? (pow2( clusters[i].y - clusters[j].y) + dPhi*dPhi) / R2
987 : 2. * (cosh( clusters[i].y - clusters[j].y) - cos(dPhi)) / R2 ;
988 if (isAnti) dijTemp /= max(clusters[i].pT2, clusters[j].pT2);
989 else if (isKT) dijTemp *= min(clusters[i].pT2, clusters[j].pT2);
990 dij[i*(i-1)/2 + j] = dijTemp;
1008 bool SlowJet::doStep() {
1011 if (useFJcore)
return false;
1012 if (clSize == 0)
return false;
1018 if (clusters[iMin].pT2 > pT2jetMin) {
1019 jets.push_back( SingleSlowJet(clusters[iMin]) );
1023 for (
int i = jtSize - 1; i > 0; --i) {
1024 if (jets[i].pT2 < jets[i-1].pT2)
break;
1025 swap( jets[i], jets[i-1]);
1034 clusters[jMin].p += clusters[iMin].p;
1035 clusters[jMin].pT2 = max( TINY*TINY, clusters[jMin].p.pT2());
1036 double mTTemp = sqrt(clusters[jMin].p.m2Calc() + clusters[jMin].pT2);
1037 clusters[jMin].y = (clusters[jMin].p.pz() > 0)
1038 ? log( max( TINY, clusters[jMin].p.e() + clusters[jMin].p.pz() )
1039 / mTTemp ) : log( mTTemp
1040 / max( TINY, clusters[jMin].p.e() - clusters[jMin].p.pz() ) );
1041 clusters[jMin].phi = clusters[jMin].p.phi();
1042 clusters[jMin].mult += clusters[iMin].mult;
1043 clusters[jMin].idx.insert(clusters[iMin].idx.begin(),
1044 clusters[iMin].idx.end());
1047 if (isAnti) diB[jMin] = 1. / clusters[jMin].pT2;
1048 else if (isKT) diB[jMin] = clusters[jMin].pT2;
1049 else diB[jMin] = 1.;
1050 for (
int i = 0; i < clSize; ++i)
if (i != jMin && i != iMin) {
1051 dPhi = abs( clusters[i].phi - clusters[jMin].phi );
1052 if (dPhi > M_PI) dPhi = 2. * M_PI - dPhi;
1053 dijTemp = (useStandardR)
1054 ? (pow2( clusters[i].y - clusters[jMin].y) + dPhi*dPhi) / R2
1055 : 2. * (cosh( clusters[i].y - clusters[jMin].y) - cos(dPhi)) / R2 ;
1056 if (isAnti) dijTemp /= max(clusters[i].pT2, clusters[jMin].pT2);
1057 else if (isKT) dijTemp *= min(clusters[i].pT2, clusters[jMin].pT2);
1058 if (i < jMin) dij[jMin*(jMin-1)/2 + i] = dijTemp;
1059 else dij[i*(i-1)/2 + jMin] = dijTemp;
1064 if (iMin < clLast) {
1065 clusters[iMin] = clusters[clLast];
1066 diB[iMin] = diB[clLast];
1067 for (
int j = 0; j < iMin; ++j)
1068 dij[iMin*(iMin-1)/2 + j] = dij[clLast*(clLast-1)/2 + j];
1069 for (
int j = iMin + 1; j < clLast; ++j)
1070 dij[j*(j-1)/2 + iMin] = dij[clLast*(clLast-1)/2 + j];
1074 clusters.pop_back();
1090 void SlowJet::list(
bool listAll)
const {
1093 if (useFJcore) cout <<
"\n -- PYTHIA SlowJet(fjcore) Listing, p = ";
1094 else cout <<
"\n -- PYTHIA SlowJet(native) Listing, p = ";
1095 cout << setw(2) << power <<
", R = " << fixed << setprecision(3) << setw(5)
1096 << R <<
", pTjetMin =" << setw(8) << pTjetMin <<
", etaMax = "
1097 << setw(6) << etaMax <<
" -- \n \n no pTjet y phi"
1098 <<
" mult p_x p_y p_z e m \n";
1101 for (
int i = 0; i < jtSize; ++i) {
1102 cout << setw(5) << i << setw(11) << sqrt(jets[i].pT2) << setw(9)
1103 << jets[i].y << setw(9) << jets[i].phi << setw(6)
1104 << jets[i].mult << setw(11) << jets[i].p.px() << setw(11)
1105 << jets[i].p.py() << setw(11) << jets[i].p.pz() << setw(11)
1106 << jets[i].p.e() << setw(11) << jets[i].p.mCalc() <<
"\n";
1110 if (listAll && clSize > 0) {
1111 cout <<
" -------- Below this line follows remaining clusters,"
1112 <<
" still pT-unordered -------------------\n";
1113 for (
int i = 0; i < clSize; ++i) {
1114 cout << setw(5) << i + jtSize << setw(11) << sqrt(clusters[i].pT2)
1115 << setw(9) << clusters[i].y << setw(9) << clusters[i].phi
1116 << setw(6) << clusters[i].mult << setw(11) << clusters[i].p.px()
1117 << setw(11) << clusters[i].p.py() << setw(11) << clusters[i].p.pz()
1118 << setw(11) << clusters[i].p.e() << setw(11)
1119 << clusters[i].p.mCalc() <<
"\n";
1124 cout <<
"\n -------- End PYTHIA SlowJet Listing ------------------"
1125 <<
"--------------------------------------" << endl;
1133 void SlowJet::findNext() {
1140 for (
int i = 1; i < clSize; ++i) {
1141 if (diB[i] < dMin) {
1146 for (
int j = 0; j < i; ++j) {
1147 if (dij[i*(i-1)/2 + j] < dMin) {
1150 dMin = dij[i*(i-1)/2 + j];
1168 bool SlowJet::clusterFJ() {
1171 vector<fjcore::PseudoJet> inFJ;
1172 for (
int i = 0; i < int(clusters.size()); ++i) {
1173 inFJ.push_back( fjcore::PseudoJet( clusters[i].p.px(),
1174 clusters[i].p.py(), clusters[i].p.pz(), clusters[i].p.e() ) );
1175 set<int>::iterator idxTmp = clusters[i].idx.begin();
1176 inFJ[i].set_user_index( *idxTmp);
1180 fjcore::JetAlgorithm jetAlg = fjcore::cambridge_algorithm;
1181 if (isAnti) jetAlg = fjcore::antikt_algorithm;
1182 if (isKT) jetAlg = fjcore::kt_algorithm;
1183 fjcore::JetDefinition jetDef( jetAlg, R);
1186 fjcore::ClusterSequence clust_seq( inFJ, jetDef);
1189 vector<fjcore::PseudoJet> outFJ
1190 = sorted_by_pt( clust_seq.inclusive_jets( pTjetMin) );
1194 double pT2Temp, yTemp, phiTemp;
1195 for (
int i = 0; i < int(outFJ.size()); ++i) {
1196 pTemp = Vec4( outFJ[i].px(), outFJ[i].py(), outFJ[i].pz(), outFJ[i].e() );
1197 pT2Temp = outFJ[i].pt2();
1198 yTemp = outFJ[i].rap();
1199 phiTemp = outFJ[i].phi_std();
1200 jets.push_back( SingleSlowJet(pTemp, pT2Temp, yTemp, phiTemp, 0) );
1203 jets[i].idx.clear();
1204 vector<fjcore::PseudoJet> constFJ = outFJ[i].constituents();
1205 for (
int j = 0; j < int(constFJ.size()); ++j)
1206 jets[i].idx.insert( constFJ[j].user_index() );
1207 jets[i].mult = constFJ.size();
1213 jtSize = outFJ.size();