9 #include "Pythia8/BeamParticle.h"
25 vector< pair<int,int> > lastSteps;
40 const double BeamParticle::XMINUNRESOLVED = 1. - 1e-10;
43 const double BeamParticle::POMERONMASS = 1.;
46 const double BeamParticle::XMAXCOMPANION = 0.99;
49 const double BeamParticle::TINYZREL = 1e-8;
52 const int BeamParticle::NMAX = 1000;
56 const int BeamParticle::NRANDOMTRIES = 1000;
62 void BeamParticle::init(
int idIn,
double pzIn,
double eIn,
double mIn,
64 Rndm* rndmPtrIn,
PDF* pdfInPtr,
PDF* pdfHardInPtr,
bool isUnresolvedIn,
69 particleDataPtr = particleDataPtrIn;
71 pdfBeamPtr = pdfInPtr;
72 pdfHardBeamPtr = pdfHardInPtr;
73 isUnresolvedBeam = isUnresolvedIn;
74 flavSelPtr = flavSelPtrIn;
78 pdfBeamPtrSave = pdfInPtr;
79 pdfHardBeamPtrSave = pdfHardInPtr;
82 bool beamHasGamma = settings.flag(
"PDF:lepton2gamma");
85 hasResGammaInBeam = beamHasGamma;
88 maxValQuark = settings.mode(
"BeamRemnants:maxValQuark");
91 valencePowerMeson = settings.parm(
"BeamRemnants:valencePowerMeson");
92 valencePowerUinP = settings.parm(
"BeamRemnants:valencePowerUinP");
93 valencePowerDinP = settings.parm(
"BeamRemnants:valencePowerDinP");
96 valenceDiqEnhance = settings.parm(
"BeamRemnants:valenceDiqEnhance");
99 companionPower = settings.mode(
"BeamRemnants:companionPower");
102 gluonPower = settings.parm(
"BeamRemnants:gluonPower");
103 xGluonCutoff = settings.parm(
"BeamRemnants:xGluonCutoff");
106 allowJunction = settings.flag(
"BeamRemnants:allowJunction");
110 beamJunction = settings.flag(
"beamRemnants:beamJunction");
113 allowBeamJunctions = settings.flag(
"beamRemnants:allowBeamJunction");
116 pickQuarkNorm = settings.parm(
"Diffraction:pickQuarkNorm");
117 pickQuarkPower = settings.parm(
"Diffraction:pickQuarkPower");
120 beamSat = settings.parm(
"BeamRemnants:saturation");
123 diffPrimKTwidth = settings.parm(
"Diffraction:primKTwidth");
126 diffLargeMassSuppress = settings.parm(
"Diffraction:largeMassSuppress");
129 doND = settings.flag(
"SoftQCD:nonDiffractive");
130 doISR = settings.flag(
"PartonLevel:ISR");
131 doMPI = settings.flag(
"PartonLevel:MPI");
132 pTminISR = settings.parm(
"SpaceShower:pTmin");
137 pBeam =
Vec4( 0., 0., pzIn, eIn);
142 resetGammaInLepton();
154 void BeamParticle::initBeamKind() {
157 idBeamAbs = abs(idBeam);
158 isLeptonBeam =
false;
159 isHadronBeam =
false;
161 isBaryonBeam =
false;
170 if ( (idBeamAbs > 10 && idBeamAbs < 17)
171 || (idBeamAbs > 50 && idBeamAbs < 60) ) {
179 if (idBeamAbs == 22) {
189 if (idBeamAbs < 101 || idBeamAbs > 9999)
return;
192 if (idBeamAbs == 990) {
200 }
else if (idBeamAbs < 1000) {
201 int id1 = idBeamAbs/100;
202 int id2 = (idBeamAbs/10)%10;
203 if ( id1 < 1 || id1 > maxValQuark
204 || id2 < 1 || id2 > maxValQuark )
return;
222 int id1 = idBeamAbs/1000;
223 int id2 = (idBeamAbs/100)%10;
224 int id3 = (idBeamAbs/10)%10;
225 if ( id1 < 1 || id1 > maxValQuark || id2 < 1 || id2 > maxValQuark
226 || id3 < 1 || id3 > maxValQuark)
return;
227 if (id2 > id1 || id3 > id1)
return;
231 nValKinds = 1; idVal[0] = id1; nVal[0] = 1;
232 if (id2 == id1) ++nVal[0];
238 if (id3 == id1) ++nVal[0];
239 else if (id3 == id2) ++nVal[1];
241 idVal[nValKinds] = id3;
248 if (idBeam < 0)
for (
int i = 0; i < nValKinds; ++i) idVal[i] = -idVal[i];
258 void BeamParticle::initUnres(PDF* pdfUnresInPtr) {
261 pdfUnresBeamPtr = pdfUnresInPtr;
262 isResUnres = (pdfUnresBeamPtr != 0 ) ?
true :
false;
270 void BeamParticle::newValenceContent() {
273 if (idBeam == 111 || idBeam == 113 || idBeam == 223) {
274 idVal[0] = (rndmPtr->flat() < 0.5) ? 1 : 2;
275 idVal[1] = -idVal[0];
278 }
else if (idBeam == 130 || idBeam == 310) {
279 idVal[0] = (rndmPtr->flat() < 0.5) ? 1 : 3;
280 idVal[1] = (idVal[0] == 1) ? -3 : -1;
283 }
else if (idBeam == 990) {
284 idVal[0] = (rndmPtr->flat() < 0.5) ? 1 : 2;
285 idVal[1] = -idVal[0];
288 }
else if (idBeam == 22) {
292 if (idTmp == 113 || idTmp == 223) {
293 idVal[0] = (rndmPtr->flat() < 0.5) ? 1 : 2;
294 idVal[1] = -idVal[0];
296 }
else if (idTmp == 333) {
304 idVal[1] = -idVal[0];
308 }
else if (idBeam == 333) {
310 idVal[1] = -idVal[0];
316 pdfBeamPtr->newValenceContent( idVal[0], idVal[1]);
317 if (pdfHardBeamPtr != pdfBeamPtr && pdfHardBeamPtr != 0)
318 pdfHardBeamPtr->newValenceContent( idVal[0], idVal[1]);
324 double BeamParticle::xMax(
int iSkip) {
328 if (idBeam == 990) xLeft -= POMERONMASS / e();
329 else if (isHadron()) xLeft -= m() / e();
330 if (size() == 0)
return xLeft;
333 for (
int i = 0; i < size(); ++i)
334 if (i != iSkip && resolved[i].isFromBeam()) xLeft -= resolved[i].x();
345 double BeamParticle::xfModified(
int iSkip,
int idIn,
double x,
double Q2) {
356 if (x >= 1.)
return 0.;
357 bool canBeVal =
false;
358 for (
int i = 0; i < nValKinds; ++i)
359 if (idIn == idVal[i]) canBeVal =
true;
361 xqVal = xfVal( idIn, x, Q2);
362 xqgSea = xfSea( idIn, x, Q2);
364 else xqgSea = xf( idIn, x, Q2);
371 for (
int i = 0; i < size(); ++i)
372 if (i != iSkip) xUsed += resolved[i].x();
373 double xLeft = 1. - xUsed;
374 if (x >= xLeft)
return 0.;
375 double xRescaled = x / xLeft;
379 double xValLeft = 0.;
380 for (
int i = 0; i < nValKinds; ++i) {
381 nValLeft[i] = nVal[i];
382 for (
int j = 0; j < size(); ++j)
383 if (j != iSkip && resolved[j].isValence()
384 && resolved[j].id() == idVal[i]) --nValLeft[i];
385 double xValNow = xValFrac(i, Q2);
386 xValTot += nVal[i] * xValNow;
387 xValLeft += nValLeft[i] * xValNow;
391 double xCompAdded = 0.;
392 for (
int i = 0; i < size(); ++i)
393 if (i != iSkip && resolved[i].isUnmatched()) xCompAdded
394 += xCompFrac( resolved[i].x() / (xLeft + resolved[i].x()) )
398 * (1. + resolved[i].x() / xLeft);
401 double rescaleGS = max( 0., (1. - xValLeft - xCompAdded)
403 xqgSea = rescaleGS * xfSea( idIn, xRescaled, Q2);
406 for (
int i = 0; i < nValKinds; ++i)
407 if (idIn == idVal[i] && nValLeft[i] > 0)
408 xqVal = xfVal( idIn, xRescaled, Q2)
409 * double(nValLeft[i]) / double(nVal[i]);
412 for (
int i = 0; i < size(); ++i)
413 if (i != iSkip && resolved[i].
id() == -idIn
414 && resolved[i].isUnmatched()) {
415 double xsRescaled = resolved[i].x() / (xLeft + resolved[i].x());
416 double xcRescaled = x / (xLeft + resolved[i].x());
417 double xqCompNow = xCompDist( xcRescaled, xsRescaled);
420 if ( isGamma() ) xqCompNow *= xIntegratedPDFs(Q2);
421 resolved[i].xqCompanion( xqCompNow);
422 xqCompSum += xqCompNow;
428 xqgTot = xqVal + xqgSea + xqCompSum;
431 if (isGammaBeam && doISR)
return xqgTot;
434 if (resolved[iSkip].isValence())
return xqVal;
435 if (resolved[iSkip].isUnmatched())
return xqgSea + xqCompSum;
447 int BeamParticle::pickValSeaComp() {
450 int oldCompanion = resolved[iSkipSave].companion();
451 if (oldCompanion >= 0) resolved[oldCompanion].companion(-2);
457 if (idSave == 21 || idSave == 22) vsc = -1;
460 else if (isLeptonBeam && idSave == idBeam) vsc = -3;
465 double xqRndm = xqgTot * rndmPtr->flat();
466 if (xqRndm < xqVal && !isGammaBeam ) vsc = -3;
467 else if (xqRndm < xqVal + xqgSea) vsc = -2;
471 xqRndm -= xqVal + xqgSea;
472 for (
int i = 0; i < size(); ++i)
473 if (i != iSkipSave && resolved[i].
id() == -idSave
474 && resolved[i].isUnmatched()) {
475 xqRndm -= resolved[i].xqCompanion();
476 if (xqRndm < 0.) vsc = i;
483 resolved[iSkipSave].companion(vsc);
484 if (vsc >= 0) resolved[vsc].companion(iSkipSave);
496 double BeamParticle::xValFrac(
int j,
double Q2) {
499 if (Q2 != Q2ValFracSav) {
503 double llQ2 = log( log( max( 1., Q2) / 0.04 ));
506 uValInt = 0.48 / (1. + 1.56 * llQ2);
507 dValInt = 0.385 / (1. + 1.60 * llQ2);
511 if (isBaryonBeam && nValKinds == 3)
return (2. * uValInt + dValInt) / 3.;
514 if (isBaryonBeam && nVal[j] == 1)
return dValInt;
515 if (isBaryonBeam && nVal[j] == 2)
return uValInt;
518 return 0.5 * (2. * uValInt + dValInt);
528 double BeamParticle::xCompFrac(
double xs) {
531 if (xs > XMAXCOMPANION)
return 0.;
534 switch (companionPower) {
537 return xs * ( 5. + xs * (-9. - 2. * xs * (-3. + xs)) + 3. * log(xs) )
538 / ( (-1. + xs) * (2. + xs * (-1. + 2. * xs)) );
541 return -1. -3. * xs + ( 2. * pow2(-1. + xs) * (1. + xs + xs*xs))
542 / ( 2. + xs*xs * (xs - 3.) + 3. * xs * log(xs) );
545 return xs * ( (1. - xs) * (19. + xs * (43. + 4. * xs))
546 + 6. * log(xs) * (1. + 6. * xs + 4.*xs*xs) ) /
547 ( 4. * ( (xs - 1.) * (1. + xs * (4. + xs) )
548 - 3. * xs * log(xs) * (1 + xs) ) );
551 return 3. * xs * ( (xs - 1.) * (7. + xs * (28. + 13. * xs))
552 - 2. * log(xs) * (1. + xs * (9. + 2. * xs * (6. + xs))) )
553 / ( 4. + 27. * xs - 31. * pow3(xs)
554 + 6. * xs * log(xs) * (3. + 2. * xs * (3.+xs)) );
557 return ( -9. * xs * (xs*xs - 1.) * (5. + xs * (24. + xs)) + 12. * xs
558 * log(xs) * (1. + 2. * xs) * (1. + 2. * xs * (5. + 2. * xs)) )
559 / ( 8. * (1. + 2. * xs) * ((xs - 1.) * (1. + xs * (10. + xs))
560 - 6. * xs * log(xs) * (1. + xs)) );
571 double BeamParticle::xCompDist(
double xc,
double xs) {
574 if (xs > XMAXCOMPANION)
return 0.;
578 if (xg > 1.)
return 0.;
582 double fac = 3. * xc * xs * (xc*xc + xs*xs) / pow4(xg);
585 switch (companionPower) {
588 return fac / ( 2. - xs * (3. - xs * (3. - 2. * xs)) );
591 return fac * (1. - xg) / ( 2. + xs*xs * (-3. + xs) + 3. * xs * log(xs) );
594 return fac * pow2(1. - xg) / ( 2. * ((1. - xs) * (1. + xs * (4. + xs))
595 + 3. * xs * (1. + xs) * log(xs)) );
598 return fac * pow3(1. - xg) * 2. / ( 4. + 27. * xs - 31. * pow3(xs)
599 + 6. * xs * log(xs) * (3. + 2. * xs * (3. + xs)) );
602 return fac * pow4(1. - xg) / ( 2. * (1. + 2. * xs) * ((1. - xs)
603 * (1. + xs * (10. + xs)) + 6. * xs * log(xs) * (1. + xs)) );
612 void BeamParticle::setGammaMode(
int gammaModeIn) {
617 pdfBeamPtr = pdfBeamPtrSave;
618 pdfHardBeamPtr = pdfHardBeamPtrSave;
619 hasResGammaInBeam =
false;
620 isResolvedGamma =
false;
625 gammaMode = gammaModeIn;
628 if (gammaMode == 2 && isResUnres) {
629 pdfBeamPtr = pdfUnresBeamPtr;
630 pdfHardBeamPtr = pdfUnresBeamPtr;
631 isResolvedGamma =
false;
632 hasResGammaInBeam =
false;
635 if ( isGamma()) isUnresolvedBeam =
true;
639 pdfBeamPtr = pdfBeamPtrSave;
640 pdfHardBeamPtr = pdfHardBeamPtrSave;
641 isUnresolvedBeam =
false;
642 if ( isGamma()) isResolvedGamma =
true;
643 else isResolvedGamma =
false;
644 if ( isLepton() && gammaMode == 1 ) hasResGammaInBeam =
true;
645 else hasResGammaInBeam =
false;
654 bool BeamParticle::gammaInitiatorIsVal(
int iResolved,
double Q2) {
655 return gammaInitiatorIsVal( iResolved, resolved[iResolved].
id(),
656 resolved[iResolved].x(), Q2 );
664 bool BeamParticle::gammaInitiatorIsVal(
int iResolved,
int idInit,
665 double x,
double Q2) {
671 if ( idInit == 0 || abs(idInit) == 21 ) {
672 idVal[0] = pdfBeamPtr->sampleGammaValFlavor(Q2);
673 idVal[1] = -idVal[0];
681 pdfBeamPtr->newValenceContent( idVal[0], idVal[1]);
684 if ( iResolved == iGamVal ) {
689 }
else if ( Q2 < pdfBeamPtr->gammaPDFRefScale(idInit) ) {
695 double xVal = xfVal( idInit, x, Q2);
696 double xSea = xfSea( idInit, x, Q2);
697 if ( rndmPtr->flat() < xVal/( xVal + xSea ) ) {
703 idVal[0] = pdfBeamPtr->sampleGammaValFlavor(Q2);
704 idVal[1] = -idVal[0];
715 int BeamParticle::gammaValSeaComp(
int iResolved) {
721 if ( resolved[iResolved].
id() == 21 ||
722 resolved[iResolved].
id() == 22 ) vsc = -1;
725 else if (iResolved == iPosVal) vsc = -3;
726 resolved[iResolved].companion(vsc);
735 bool BeamParticle::remnantFlavours(
Event& event,
bool isDIS) {
738 if (isHadronBeam && isUnresolvedBeam) {
739 append( 0, idBeam, 0., -1);
740 resolved[1].m( particleDataPtr->m0( idBeamAbs ) );
745 hasJunctionBeam = (isBaryon());
749 if (isDIS && nInit != 1)
return false;
755 if (resolved[0].
id() == 22)
return true;
758 if ( doISR && !doMPI ) {
762 if ( isResolvedGamma ) {
763 gammaInitiatorIsVal(0, pow2(pTminISR));
769 }
else if ( doMPI || doND ) {
773 double pTmin = doISR ? pTminISR : pTminMPI;
777 if ( iGamVal >= 0 ) {
778 gammaInitiatorIsVal(iGamVal, pow2(pTmin));
779 int valComp = resolved[iGamVal].companion();
780 if ( valComp >= 0 ) resolved[valComp].companion(-3);
781 gammaValSeaComp(iGamVal);
787 for (
int i = 0; i < size(); ++i) {
788 bool isValence = gammaInitiatorIsVal(i, pow2(pTminISR));
789 int origComp = resolved[i].companion();
794 if ( origComp >= 0 ) resolved[origComp].companion(-3);
803 else if ( isLeptonBeam && gammaMode == 2)
return true;
806 for (
int i = 0; i < nValKinds; ++i) {
807 nValLeft[i] = nVal[i];
808 for (
int j = 0; j < nInit; ++j)
if (resolved[j].isValence()
809 && resolved[j].
id() == idVal[i]) --nValLeft[i];
811 if ( isGammaBeam && doISR && !isResolvedGamma && !doMPI ) nValLeft[i] = 0;
813 for (
int k = 0; k < nValLeft[i]; ++k) append(0, idVal[i], 0., -3);
817 int nInitPlusVal = size();
818 if (isBaryon() && nInitPlusVal - nInit >= 2) {
823 if (nInitPlusVal - nInit == 3) {
824 double pickDq = 3. * rndmPtr->flat();
825 if (pickDq > 1.) iQ2 = nInit + 2;
826 if (pickDq > 2.) iQ1 = nInit + 1;
830 int idDq = flavSelPtr->makeDiquark( resolved[iQ1].
id(),
831 resolved[iQ2].
id(), idBeam);
834 resolved[iQ1].id(idDq);
835 if (nInitPlusVal - nInit == 3 && iQ2 == nInit + 1)
836 resolved[nInit + 1].id( resolved[nInit + 2].
id() );
838 hasJunctionBeam =
false;
842 for (
int i = 0; i < nInit; ++i)
843 if (resolved[i].isUnmatched()) {
846 append(0, -resolved[i].
id(), 0., i);
847 resolved[i].companion(size() - 1);
852 if ( size() == nInit && !isUnresolvedBeam &&
853 (!isGammaBeam || isResolvedGamma || doMPI) ) {
854 int idRemnant = (isHadronBeam || isGammaBeam) ? 21 : 22;
855 append(0, idRemnant, 0., -1);
859 if (isHadronBeam && isDIS && size() > 2 && resolved[0].
id() != 21) {
861 infoPtr->errorMsg(
"Error in BeamParticle::remnantFlavours: "
862 "unexpected number of beam remnants for DIS");
867 int colTypeComp = particleDataPtr->colType( resolved[3].
id() );
868 int colType1 = particleDataPtr->colType( resolved[1].
id() );
869 int i12 = (colType1 == -colTypeComp) ? 1 : 2;
872 int idHad = flavSelPtr->combineId( resolved[i12].
id(), resolved[3].
id() );
874 infoPtr->errorMsg(
"Error in BeamParticle::remnantFlavours: "
875 "failed to combine hadron for DIS");
880 resolved[i12].id(idHad);
882 resolved[0].companion(-3);
886 for (
int i = 0; i < size(); ++i) {
887 if (i < nInit) resolved[i].m(0.);
888 else resolved[i].m( particleDataPtr->m0( resolved[i].id() ) );
892 if (hasJunctionBeam && !allowJunction)
return false;
895 for (
int i = nInit; i < size(); ++i) {
896 int colType = particleDataPtr->colType( resolved[i].
id() );
897 int col = (colType == 1 || colType == 2) ? event.nextColTag() : 0;
898 int acol = (colType == -1 || colType == 2) ? event.nextColTag() : 0;
899 resolved[i].cols( col, acol);
911 bool BeamParticle::remnantColours(
Event& event, vector<int>& colFrom,
912 vector<int>& colTo) {
915 if (isLeptonBeam)
return true;
918 for (
int i = 0; i < size(); ++i) {
919 int j = resolved[i].iPos();
920 resolved[i].cols( event[j].col(), event[j].acol());
928 for (
int i = 0; i < size(); ++i)
929 if (resolved[i].isFromBeam()) {
930 if ( resolved[i].isValence() ) iVal.push_back(i);
931 else if ( resolved[i].isCompanion() && resolved[i].companion() > i )
933 else if ( resolved[i].
id() == 21
934 && resolved[i].col() != resolved[i].acol() ) iGlu.push_back(i);
940 if (iVal.size() != 0) {
942 if (iVal.size() == 2) {
943 if ( abs(resolved[iValSel].
id()) > 10 ) iValSel = iVal[1];
945 double rndmValSel = 3. * rndmPtr->flat();
946 if (rndmValSel > 1.) iValSel= iVal[1];
947 if (rndmValSel > 2.) iValSel= iVal[2];
953 bool hasCol = (resolved[iBeg].col() > 0);
954 int begCol = (hasCol) ? resolved[iBeg].col() : resolved[iBeg].acol();
957 vector<int> iGluRndm;
958 for (
int i = 0; i < int(iGlu.size()); ++i)
959 iGluRndm.push_back( iGlu[i] );
960 for (
int iOrder = 0; iOrder < int(iGlu.size()); ++iOrder) {
961 int iRndm = int(
double(iGluRndm.size()) * rndmPtr->flat());
962 int iGluSel = iGluRndm[iRndm];
963 iGluRndm[iRndm] = iGluRndm[iGluRndm.size() - 1];
968 int endCol = (hasCol) ? resolved[iEnd].acol() : resolved[iEnd].col();
971 iEnd = resolved[iEnd].companion();
972 endCol = (hasCol) ? resolved[iEnd].acol() : resolved[iEnd].col();
976 if (begCol < endCol) {
977 if (hasCol) resolved[iEnd].acol(begCol);
978 else resolved[iEnd].col(begCol);
979 colFrom.push_back(endCol);
980 colTo.push_back(begCol);
982 if (hasCol) resolved[iBeg].col(endCol);
983 else resolved[iBeg].acol(endCol);
984 colFrom.push_back(begCol);
985 colTo.push_back(endCol);
990 begCol = (hasCol) ? resolved[iBeg].col() : resolved[iBeg].acol();
993 iBeg = resolved[iBeg].companion();
994 begCol = (hasCol) ? resolved[iBeg].col() : resolved[iBeg].acol();
1002 vector<int> colList;
1003 vector<int> acolList;
1004 for (
int i = 0; i < size(); ++i)
1005 if ( resolved[i].isFromBeam() )
1006 if ( resolved[i].col() != resolved[i].acol() ) {
1007 if (resolved[i].col() > 0) colList.push_back( resolved[i].col() );
1008 if (resolved[i].acol() > 0) acolList.push_back( resolved[i].acol() );
1012 bool foundPair =
true;
1013 while (foundPair && colList.size() > 0 && acolList.size() > 0) {
1015 for (
int iCol = 0; iCol < int(colList.size()); ++iCol) {
1016 for (
int iAcol = 0; iAcol < int(acolList.size()); ++iAcol) {
1017 if (acolList[iAcol] == colList[iCol]) {
1018 colList[iCol] = colList.back();
1020 acolList[iAcol] = acolList.back();
1021 acolList.pop_back();
1025 }
if (foundPair)
break;
1030 if (colList.size() == 1 && acolList.size() == 1) {
1031 int finalFrom = max( colList[0], acolList[0]);
1032 int finalTo = min( colList[0], acolList[0]);
1033 for (
int i = 0; i < size(); ++i)
1034 if ( resolved[i].isFromBeam() ) {
1035 if (resolved[i].col() == finalFrom) resolved[i].col(finalTo);
1036 if (resolved[i].acol() == finalFrom) resolved[i].acol(finalTo);
1038 colFrom.push_back(finalFrom);
1039 colTo.push_back(finalTo);
1042 }
else if (hasJunctionBeam && colList.size() == 3
1043 && acolList.size() == 0) {
1044 event.appendJunction( 1, colList[0], colList[1], colList[2]);
1045 junCol[0] = colList[0];
1046 junCol[1] = colList[1];
1047 junCol[2] = colList[2];
1048 }
else if (hasJunctionBeam && acolList.size() == 3
1049 && colList.size() == 0) {
1050 event.appendJunction( 2, acolList[0], acolList[1], acolList[2]);
1051 junCol[0] = acolList[0];
1052 junCol[1] = acolList[1];
1053 junCol[2] = acolList[2];
1056 }
else if (colList.size() > 0 || acolList.size() > 0) {
1057 infoPtr->errorMsg(
"Error in BeamParticle::remnantColours: "
1058 "leftover unmatched colours");
1063 for (
int i = nInit; i < size(); ++i)
1064 event[resolved[i].iPos()].cols( resolved[i].col(), resolved[i].acol() );
1075 double BeamParticle::xRemnant(
int i) {
1080 int idAbs = abs(resolved[i].
id());
1081 if (idAbs > 100 && (idAbs/10)%10 != 0) {
1085 }
else if (resolved[i].isValence()) {
1088 int id1 = resolved[i].id();
1090 if (abs(id1) > 1000) {
1091 id2 = (id1 > 0) ? (id1/100)%10 : -(((-id1)/100)%10);
1092 id1 = (id1 > 0) ? id1/1000 : -((-id1)/1000);
1096 for (
int iId = 0; iId < 2; ++iId) {
1097 int idNow = (iId == 0) ? id1 : id2;
1098 if (idNow == 0)
break;
1102 double xPow = valencePowerMeson;
1104 if (nValKinds == 3 || nValKinds == 1)
1105 xPow = (3. * rndmPtr->flat() < 2.)
1106 ? valencePowerUinP : valencePowerDinP ;
1107 else if (nValence(idNow) == 2) xPow = valencePowerUinP;
1108 else xPow = valencePowerDinP;
1110 do xPart = pow2( rndmPtr->flat() );
1111 while ( pow(1. - xPart, xPow) < rndmPtr->flat() );
1116 if (id2 != 0) x *= valenceDiqEnhance;
1119 }
else if (resolved[i].isCompanion()) {
1123 for (
int iInit = 0; iInit < nInit; ++iInit)
1124 if (resolved[iInit].isFromBeam()) xLeft -= resolved[iInit].x();
1125 double xCompanion = resolved[ resolved[i].companion() ].x();
1126 xCompanion /= (xLeft + xCompanion);
1129 do x = pow( xCompanion, rndmPtr->flat()) - xCompanion;
1130 while ( pow( (1. - x - xCompanion) / (1. - xCompanion), companionPower)
1131 * (pow2(x) + pow2(xCompanion)) / pow2(x + xCompanion)
1132 < rndmPtr->flat() );
1137 do x = pow(xGluonCutoff, 1 - rndmPtr->flat());
1138 while ( pow(1. - x, gluonPower) < rndmPtr->flat() );
1149 double BeamParticle::remnantMass(
int idIn) {
1155 double mRem = particleDataPtr->m0(
id());
1156 int valSign1 = (nValence(idIn) > 0) ? -1 : 1;
1157 mRem += valSign1 * particleDataPtr->m0(idIn);
1162 }
else if ( isGamma() ) {
1163 if ( isUnresolved() )
return 0.;
1164 double mRem = (idIn == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1165 particleDataPtr->m0( idIn );
1176 bool BeamParticle::roomFor1Remnant(
double eCM) {
1179 if (!isResolvedGamma)
return true;
1182 double x1 = resolved[0].x();
1183 int id1 = resolved[0].id();
1184 return roomFor1Remnant(id1, x1, eCM);
1191 bool BeamParticle::roomFor1Remnant(
int id1,
double x1,
double eCM) {
1198 double mRemnant = (id1 == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1199 particleDataPtr->m0( id1 );
1200 return ( mRemnant < eCM*( 1 - sqrt(x1) ) );
1207 bool BeamParticle::roomFor2Remnants(
int id1,
double x1,
double eCM) {
1211 double id2 = resolved[0].id();
1212 double x2 = resolved[0].x();
1216 double m1 = (id1 == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1217 particleDataPtr->m0( id1 );
1218 double m2 = (id2 == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1219 particleDataPtr->m0( id2 );
1220 return ( (m1 + m2) < eCM*sqrt( (1.0 - x1)*(1.0 - x2) ) );
1228 bool BeamParticle::roomForRemnants(BeamParticle beamOther) {
1231 double xLeftA = this->xMax(-1);
1232 double xLeftB = beamOther.xMax(-1);
1233 double eCM = infoPtr->eCM();
1234 double Wleft = eCM * sqrt(xLeftA * xLeftB);
1237 bool allGluonsA =
true;
1238 bool allGluonsB =
true;
1241 for (
int i = 0; i < this->size(); ++i)
1242 if ( resolved[i].
id() != 21 ) {
1245 if ( resolved[i].companion() < 0 && resolved[i].companion() != -3 )
1246 mRemA += particleDataPtr->m0( resolved[i].id() );
1248 for (
int i = 0; i < beamOther.size(); ++i)
1249 if ( beamOther[i].
id() != 21 ) {
1251 if ( beamOther[i].companion() < 0 && beamOther[i].companion() != -3 )
1252 mRemB += particleDataPtr->m0( beamOther[i].id() );
1257 if ( allGluonsA) mRemA = this->isGamma() ? 2*particleDataPtr->m0(2) : 0.;
1258 if ( allGluonsB) mRemB = beamOther.isGamma() ? 2*particleDataPtr->m0(2) : 0.;
1261 if ( Wleft < mRemA + mRemB )
return false;
1269 void BeamParticle::list()
const {
1272 cout <<
"\n -------- PYTHIA Partons resolved in beam -----------------"
1273 <<
"-------------------------------------------------------------\n"
1274 <<
"\n i iPos id x comp xqcomp pTfact "
1275 <<
"colours p_x p_y p_z e m \n";
1280 for (
int i = 0; i < size(); ++i) {
1281 ResolvedParton res = resolved[i];
1282 cout << fixed << setprecision(6) << setw(5) << i << setw(6) << res.iPos()
1283 << setw(8) << res.id() << setw(10) << res.x() << setw(6)
1284 << res.companion() << setw(10) << res.xqCompanion() << setw(10)
1285 << res.pTfactor() << setprecision(3) << setw(6) << res.col()
1286 << setw(6) << res.acol() << setw(11) << res.px() << setw(11)
1287 << res.py() << setw(11) << res.pz() << setw(11) << res.e()
1288 << setw(11) << res.m() <<
"\n";
1291 if (res.companion() != -10) {
1298 cout << setprecision(6) <<
" x sum:" << setw(10) << xSum
1299 << setprecision(3) <<
" p sum:"
1300 << setw(11) << pSum.px() << setw(11) << pSum.py() << setw(11)
1301 << pSum.pz() << setw(11) << pSum.e()
1302 <<
"\n\n -------- End PYTHIA Partons resolved in beam -----------"
1303 <<
"---------------------------------------------------------------"
1311 bool BeamParticle::isUnresolvedLepton() {
1314 if (!isLeptonBeam || resolved.size() > 2 || resolved[1].id() != 22
1315 || resolved[0].x() < XMINUNRESOLVED)
return false;
1324 bool BeamParticle::pickGluon(
double mDiff) {
1327 double probPickQuark = pickQuarkNorm / pow( mDiff, pickQuarkPower);
1328 return ( (1. + probPickQuark) * rndmPtr->flat() < 1. );
1336 int BeamParticle::pickValence() {
1339 int nTotVal = (isBaryonBeam) ? 3 : 2;
1340 double rnVal = rndmPtr->flat() * nTotVal;
1341 int iVal = (rnVal < 1.) ? 1 : ( (rnVal < 2.) ? 2 : 3 );
1348 for (
int i = 0; i < nValKinds; ++i)
1349 for (
int j = 0; j < nVal[i]; ++j) {
1351 if (iNow == iVal) idVal1 = idVal[i];
1352 else if ( idVal2 == 0) idVal2 = idVal[i];
1353 else idVal3 = idVal[i];
1357 if (idVal3 != 0) idVal2 = flavSelPtr->makeDiquark( idVal2, idVal3, idBeam);
1368 double BeamParticle::zShare(
double mDiff,
double m1,
double m2) {
1371 append(0, idVal1, 0., -3);
1372 append(0, idVal2, 0., -3);
1373 double m2Diff = mDiff*mDiff;
1378 double x1 = xRemnant(0);
1379 double x2 = xRemnant(0);
1380 zRel = max( TINYZREL, min( 1. - TINYZREL, x1 / (x1 + x2) ) );
1381 pair<double, double> gauss2 = rndmPtr->gauss2();
1382 pxRel = diffPrimKTwidth * gauss2.first;
1383 pyRel = diffPrimKTwidth * gauss2.second;
1386 double mTS1 = m1*m1 + pxRel*pxRel + pyRel*pyRel;
1387 double mTS2 = m2*m2 + pxRel*pxRel + pyRel*pyRel;
1388 double m2Sys = mTS1 / zRel + mTS2 / (1. - zRel);
1389 wtAcc = (m2Sys < m2Diff)
1390 ? pow( 1. - m2Sys / m2Diff, diffLargeMassSuppress) : 0.;
1391 }
while (wtAcc < rndmPtr->flat());
1403 void BeamParticle::findColSetup(
Event& event) {
1406 colSetup = make_pair(0,0);
1414 vector<vector <vector <ColState> > > colStates;
1415 colStates.resize(size() + 1);
1416 for (
int i = 0; i < size() + 1; ++i) {
1417 colStates[i].resize(2*(i + 1));
1418 for (
int j = 0; j < 2*(i+1); ++j)
1419 colStates[i][j].resize(2*(i+1));
1421 colStates[0][0][0].nTotal = 1.;
1423 bool noColouredParticles =
true;
1425 for (
int i = 0; i < size(); ++i) {
1426 for (
int j = 0; j < int(colStates[i].size()); ++j) {
1427 for (
int k = 0; k < int(colStates[i][j].size()); ++k) {
1428 if (colStates[i][j][k].nTotal < 0.5)
continue;
1429 int idParton = resolved[i].id();
1432 if (idParton > 0 && idParton < 9) {
1433 colStates[i+1][j+1][k].lastSteps.push_back(make_pair(j,k));
1434 colStates[i+1][j+1][k].nTotal += colStates[i][j][k].nTotal;
1436 colStates[i+1][j][k -1].lastSteps.push_back(make_pair(j,k));
1437 colStates[i+1][j][k -1].nTotal += colStates[i][j][k].nTotal;
1441 if (j > 0 && allowBeamJunctions) {
1442 colStates[i+1][j - 1][k + 1].lastSteps.push_back(make_pair(j,k));
1443 colStates[i+1][j - 1][k + 1].nTotal += colStates[i][j][k].nTotal;
1448 if (idParton < 0 && idParton > -9) {
1449 colStates[i+1][j][k + 1].lastSteps.push_back(make_pair(j,k));
1450 colStates[i+1][j][k + 1].nTotal += colStates[i][j][k].nTotal;
1452 colStates[i+1][j - 1][k].lastSteps.push_back(make_pair(j,k));
1453 colStates[i+1][j - 1][k].nTotal += colStates[i][j][k].nTotal;
1457 if (k > 0 && allowBeamJunctions) {
1458 colStates[i+1][j + 1][k - 1].lastSteps.push_back(make_pair(j,k));
1459 colStates[i+1][j + 1][k - 1].nTotal += colStates[i][j][k].nTotal;
1464 if (idParton == 21) {
1465 colStates[i+1][j + 1][k + 1].lastSteps.push_back(make_pair(j,k));
1466 colStates[i+1][j + 1][k + 1].nTotal += colStates[i][j][k].nTotal;
1468 colStates[i+1][j][k].lastSteps.push_back(make_pair(j,k));
1469 colStates[i+1][j][k].nTotal += colStates[i][j][k].nTotal;
1472 colStates[i+1][j][k].lastSteps.push_back(make_pair(j,k));
1473 colStates[i+1][j][k].nTotal += colStates[i][j][k].nTotal;
1475 if (j > 0 && k > 0) {
1476 colStates[i+1][j - 1][k - 1].lastSteps.push_back(make_pair(j,k));
1477 colStates[i+1][j - 1][k - 1].nTotal += colStates[i][j][k].nTotal;
1481 if (k > 0 && allowBeamJunctions) {
1482 colStates[i+1][j + 2][k - 1].lastSteps.push_back(make_pair(j,k));
1483 colStates[i+1][j + 2][k - 1].nTotal += colStates[i][j][k].nTotal;
1485 if (j > 0 && allowBeamJunctions) {
1486 colStates[i+1][j - 1][k + 2].lastSteps.push_back(make_pair(j,k));
1487 colStates[i+1][j - 1][k + 2].nTotal += colStates[i][j][k].nTotal;
1489 if (j > 1 && allowBeamJunctions) {
1490 colStates[i+1][j - 2][k + 1].lastSteps.push_back(make_pair(j,k));
1491 colStates[i+1][j - 2][k + 1].nTotal += colStates[i][j][k].nTotal;
1493 if (k > 1 && allowBeamJunctions) {
1494 colStates[i+1][j + 1][k - 2].lastSteps.push_back(make_pair(j,k));
1495 colStates[i+1][j + 1][k - 2].nTotal += colStates[i][j][k].nTotal;
1499 if (idParton != 21 && abs(idParton) > 9) {
1500 colStates[i+1][j][k].lastSteps.push_back(make_pair(j,k));
1501 colStates[i+1][j][k].nTotal += colStates[i][j][k].nTotal;
1503 noColouredParticles =
false;
1510 double totalSize = 0;
1511 for (
int i = 0;i < int(colStates[size()].size()); ++i) {
1512 for (
int j = 0;j < int(colStates[size()][i].size()); ++j) {
1515 if (i == 0 && j == 0 && !noColouredParticles)
continue;
1517 double multipletSize = (i + 1) * (j + 1) * (i + j + 2) / 2.;
1518 totalSize += colStates[size()][i][j].nTotal *
1519 multipletSize * exp(-multipletSize / beamSat);
1525 double chosenSize = rndmPtr->flat() * totalSize;
1526 for (
int i = 0;i < int(colStates[size()].size()); ++i) {
1527 for (
int j = 0;j < int(colStates[size()][i].size()); ++j) {
1530 if (i == 0 && j == 0 && !noColouredParticles)
continue;
1533 double multipletSize = (i + 1) * (j + 1) * (i + j + 2) / 2.;
1534 curSize += colStates[size()][i][j].nTotal *
1535 multipletSize * exp(-multipletSize / beamSat);
1536 if (curSize > chosenSize) {
1538 colSetup.second = j;
1542 if (curSize > chosenSize)
break;
1546 vector<pair<int, int> > colSetupChain;
1547 colSetupChain.resize(size() + 1);
1548 pair<int,int> curColSetup = make_pair(colSetup.first, colSetup.second);
1549 for (
int i = size(); i > 0; --i) {
1550 colSetupChain[i] = curColSetup;
1551 int curColSize = colStates[i][curColSetup.first][curColSetup.second].
1553 int iCurCol = int(rndmPtr->flat() * curColSize);
1554 curColSetup = colStates[i][curColSetup.first][curColSetup.second].
1557 colSetupChain[0] = curColSetup;
1561 for (
int i = 0; i < size() ; ++i) {
1564 if (resolved[i].
id() > 0 && resolved[i].id() < 9) {
1566 if (colSetupChain[i].first + 1 == colSetupChain[i + 1].first)
1567 cols.push_back(resolved[i].col());
1570 else if (colSetupChain[i].second - 1 == colSetupChain[i + 1].second) {
1571 int iAcol = int(acols.size() * rndmPtr->flat());
1572 int acol = acols[iAcol];
1573 acols.erase(acols.begin() + iAcol);
1574 updateSingleCol(acol, resolved[i].col());
1579 int iCol = int(cols.size() * rndmPtr->flat());
1580 int juncCol =
event.nextColTag();
1581 event.appendJunction(1, cols[iCol], resolved[i].col(), juncCol);
1582 event.saveJunctionSize();
1583 acols.push_back(juncCol);
1584 cols.erase(cols.begin() + iCol);
1590 if (resolved[i].
id() < 0 && resolved[i].id() > -9) {
1592 if (colSetupChain[i].second + 1 == colSetupChain[i + 1].second)
1593 acols.push_back(resolved[i].acol());
1596 else if (colSetupChain[i].first - 1 == colSetupChain[i + 1].first) {
1597 int iCol = int(cols.size() * rndmPtr->flat());
1598 int col = cols[iCol];
1599 cols.erase(cols.begin() + iCol);
1600 updateSingleCol(col, resolved[i].acol());
1605 int iAcol = int(acols.size() * rndmPtr->flat());
1606 int juncCol =
event.nextColTag();
1607 event.appendJunction(2, acols[iAcol], resolved[i].acol(), juncCol);
1608 event.saveJunctionSize();
1609 cols.push_back(juncCol);
1610 acols.erase(acols.begin() + iAcol);
1616 if (resolved[i].
id() == 21) {
1618 if (colSetupChain[i].first + 1 == colSetupChain[i + 1].first &&
1619 colSetupChain[i].second + 1 == colSetupChain[i + 1].second ) {
1620 acols.push_back(resolved[i].acol());
1621 cols.push_back(resolved[i].col());
1625 if (colSetupChain[i].first - 1 == colSetupChain[i + 1].first &&
1626 colSetupChain[i].second - 1 == colSetupChain[i + 1].second ) {
1628 int iCol = int(cols.size() * rndmPtr->flat());
1629 int col = cols[iCol];
1630 cols.erase(cols.begin() + iCol);
1631 updateSingleCol(col, resolved[i].acol());
1633 int iAcol = int(acols.size() * rndmPtr->flat());
1634 int acol = acols[iAcol];
1635 acols.erase(acols.begin() + iAcol);
1636 updateSingleCol(acol, resolved[i].col());
1641 if (colSetupChain[i].first == colSetupChain[i + 1].first &&
1642 colSetupChain[i].second == colSetupChain[i + 1].second ) {
1643 bool removeColour =
true;
1644 if (cols.size() > 0 && acols.size() > 0)
1645 removeColour = (rndmPtr->flat() > 0.5);
1646 else if (acols.size() > 0)
1647 removeColour =
false;
1650 int iCol = int(cols.size() * rndmPtr->flat());
1651 int col = cols[iCol];
1652 cols.erase(cols.begin() + iCol);
1653 cols.push_back(resolved[i].col());
1654 updateSingleCol(col, resolved[i].acol());
1658 int iAcol = int(acols.size() * rndmPtr->flat());
1659 int acol = acols[iAcol];
1660 acols.erase(acols.begin() + iAcol);
1661 acols.push_back(resolved[i].acol());
1662 updateSingleCol(acol, resolved[i].col());
1669 if (colSetupChain[i].first + 2 == colSetupChain[i + 1].first &&
1670 colSetupChain[i].second - 1 == colSetupChain[i + 1].second ) {
1671 int iAcol = int(acols.size() * rndmPtr->flat());
1672 int acol = acols[iAcol];
1673 acols.erase(acols.begin() + iAcol);
1674 int acol3 =
event.nextColTag();
1675 event.appendJunction(2,acol,resolved[i].acol(),acol3);
1676 cols.push_back(resolved[i].col());
1677 cols.push_back(acol3);
1682 if (colSetupChain[i].first - 1 == colSetupChain[i + 1].first &&
1683 colSetupChain[i].second + 2 == colSetupChain[i + 1].second ) {
1684 int iCol = int(cols.size() * rndmPtr->flat());
1685 int col = cols[iCol];
1686 cols.erase(cols.begin() + iCol);
1687 int col3 =
event.nextColTag();
1688 event.appendJunction(1,col,resolved[i].col(),col3);
1689 acols.push_back(resolved[i].acol());
1690 acols.push_back(col3);
1695 if (colSetupChain[i].first + 1 == colSetupChain[i + 1].first &&
1696 colSetupChain[i].second - 2 == colSetupChain[i + 1].second ) {
1697 int iAcol = int(acols.size() * rndmPtr->flat());
1698 int acol = acols[iAcol];
1699 acols.erase(acols.begin() + iAcol);
1700 int iAcol2 = int(acols.size() * rndmPtr->flat());
1701 int acol2 = acols[iAcol2];
1702 acols.erase(acols.begin() + iAcol2);
1703 event.appendJunction(2,acol,resolved[i].acol(),acol2);
1704 cols.push_back(resolved[i].col());
1709 if (colSetupChain[i].first - 2 == colSetupChain[i + 1].first &&
1710 colSetupChain[i].second + 1 == colSetupChain[i + 1].second ) {
1711 int iCol = int(cols.size() * rndmPtr->flat());
1712 int col = cols[iCol];
1713 cols.erase(cols.begin() + iCol);
1714 int iCol2 = int(cols.size() * rndmPtr->flat());
1715 int col2 = cols[iCol2];
1716 cols.erase(cols.begin() + iCol2);
1717 event.appendJunction(1,col,resolved[i].col(),col2);
1718 acols.push_back(resolved[i].acol());
1731 bool BeamParticle::remnantFlavoursNew(
Event& event) {
1734 hasJunctionBeam = (isBaryon());
1740 for (
int i = 0; i < nValKinds; ++i) {
1741 nValLeft[i] = nVal[i];
1742 for (
int j = 0; j < nInit; ++j)
if (resolved[j].isValence()
1743 && resolved[j].
id() == idVal[i]) --nValLeft[i];
1745 for (
int k = 0; k < nValLeft[i]; ++k) append(0, idVal[i], 0., -3);
1747 int nInitPlusVal = size();
1750 for (
int i = 0; i < nInit; ++i)
1751 if (resolved[i].isUnmatched()) {
1754 append(0, -resolved[i].
id(), 0., i);
1755 resolved[i].companion(size() - 1);
1758 if (isBaryon()) beamJunc = 1;
1759 if (
id() < 0) beamJunc = -beamJunc;
1762 int nGluons = (colSetup.first + colSetup.second - (size() - nInit)
1763 +abs( (nJuncs - nAjuncs) - beamJunc)) / 2;
1764 for (
int i = 0; i < nGluons; i++) append(0,21,0.,-1);
1768 if (size() == nInit) {
1770 append(0, 22, 0., -1);
1772 int idRemnant = int(3*rndmPtr->flat())+1;
1773 append(0, -idRemnant, 0., -1);
1774 append(0, idRemnant, 0., -1);
1775 resolved[size()-2].companion(size() - 1);
1776 resolved[size()-1].companion(size() - 2);
1780 usedCol = vector<bool>(size(),
false);
1781 usedAcol = vector<bool>(size(),
false);
1785 nDiffJuncs = nJuncs - nAjuncs - beamJunc;
1786 if (isBaryon() && nInitPlusVal - nInit >= 2 && (
1787 (nDiffJuncs > 0 && beamJunc < 0) ||
1788 (nDiffJuncs < 0 && beamJunc > 0)) ) {
1792 int iQ2 = nInit + 1;
1793 if (nInitPlusVal - nInit == 3) {
1794 double pickDq = 3. * rndmPtr->flat();
1795 if (pickDq > 1.) iQ2 = nInit + 2;
1796 if (pickDq > 2.) iQ1 = nInit + 1;
1802 if (resolved[iQ1].
id() < 0) {
1805 usedAcol[iQ1] =
true;
1806 usedAcol[iQ2] =
true;
1809 int acol = findSingleCol(event,
true,
true);
1810 if ( acol == 0)
return false;
1813 int newCol1 =
event.nextColTag();
1814 int newCol2 =
event.nextColTag();
1815 resolved[iQ1].acol(newCol1);
1816 resolved[iQ2].acol(newCol2);
1817 event.appendJunction(2, resolved[iQ1].acol(), resolved[iQ2].acol(),
1825 usedCol[iQ1] =
true;
1826 usedCol[iQ2] =
true;
1829 int col = findSingleCol(event,
false,
true);
1830 if (col == 0)
return false;
1833 int newCol1 =
event.nextColTag();
1834 int newCol2 =
event.nextColTag();
1835 resolved[iQ1].col(newCol1);
1836 resolved[iQ2].col(newCol2);
1837 event.appendJunction(1,resolved[iQ1].col(),resolved[iQ2].col(),col);
1845 int idDq = flavSelPtr->makeDiquark( resolved[iQ1].
id(),
1846 resolved[iQ2].
id(), idBeam);
1849 if (nInitPlusVal - nInit == 3)
1850 resolved[nInit + 2].id( resolved[3 * nInit + 3 - iQ2 - iQ1].
id() );
1851 resolved[nInit].id(idDq);
1852 resolved.erase(resolved.begin() + nInit + 1);
1853 hasJunctionBeam =
false;
1856 if (idDq > 0) nDiffJuncs++;
1862 while (nDiffJuncs > 0) {
1863 int acol1 = findSingleCol(event,
true,
false);
1864 int acol2 = findSingleCol(event,
true,
false);
1865 int acol3 = findSingleCol(event,
true,
true);
1866 event.appendJunction(2,acol1,acol2,acol3);
1870 while (nDiffJuncs < 0) {
1871 int col1 = findSingleCol(event,
false,
false);
1872 int col2 = findSingleCol(event,
false,
false);
1873 int col3 = findSingleCol(event,
false,
true);
1874 event.appendJunction(1,col1,col2,col3);
1879 for (
int j = 0;j < NRANDOMTRIES; ++j) {
1880 int i = int(rndmPtr->flat() * (size() - nInit) + nInit );
1882 if ( resolved[i].hasCol() && !usedCol[i]) {
1884 int acol = findSingleCol(event,
true,
true);
1885 if ( acol == 0)
return false;
1886 resolved[i].col(acol);
1889 if ( resolved[i].hasAcol() && !usedAcol[i]) {
1891 int col = findSingleCol(event,
false,
true);
1892 if (col == 0)
return false;
1893 resolved[i].acol(col);
1898 for (
int i = nInit;i < size();i++) {
1900 if ( resolved[i].hasCol() && !usedCol[i]) {
1902 int acol = findSingleCol(event,
true,
true);
1903 if ( acol == 0)
return false;
1904 resolved[i].col(acol);
1907 if ( resolved[i].hasAcol() && !usedAcol[i]) {
1909 int col = findSingleCol(event,
false,
true);
1910 if (col == 0)
return false;
1911 resolved[i].acol(col);
1916 if (cols.size() != 0 || acols.size() != 0) {
1917 infoPtr->errorMsg(
"Error in BeamParticle::RemnantFlavours: "
1918 "Colour not conserved in beamRemnants");
1923 for (
int i = 0; i < size(); ++i) {
1924 if (i < nInit) resolved[i].m(0.);
1925 else resolved[i].m( particleDataPtr->m0( resolved[i].id() ) );
1929 if (hasJunctionBeam && !allowJunction)
return false;
1940 void BeamParticle::setInitialCol(
Event& event) {
1943 for (
int i = 0;i < size(); ++i) {
1944 if (event[resolved[i].iPos()].col() != 0)
1945 resolved[i].col(event[resolved[i].iPos()].col());
1946 if (event[resolved[i].iPos()].acol() != 0)
1947 resolved[i].acol(event[resolved[i].iPos()].acol());
1956 int BeamParticle::findSingleCol(
Event& event,
bool isAcol,
1957 bool useHardScatters) {
1960 if (useHardScatters) {
1962 if (acols.size() > 0) {
1963 int iAcol = int(acols.size() * rndmPtr->flat());
1964 int acol = acols[iAcol];
1965 acols.erase(acols.begin() + iAcol);
1969 if (
int(cols.size()) > 0) {
1970 int iCol = int(cols.size() * rndmPtr->flat());
1971 int col = cols[iCol];
1972 cols.erase(cols.begin() + iCol);
1980 for (
int i = 0;i < NMAX; ++i) {
1981 int iBeam = int((size() - nInit) * rndmPtr->flat()) + nInit;
1982 if (resolved[iBeam].hasAcol() && !usedAcol[iBeam]) {
1983 int acol =
event.nextColTag();
1984 resolved[iBeam].acol(acol);
1985 usedAcol[iBeam] =
true;
1990 for (
int i = 0; i < NMAX; ++i) {
1991 int iBeam = int((size() - nInit) * rndmPtr->flat()) + nInit;
1992 if (resolved[iBeam].hasCol() && !usedCol[iBeam]) {
1993 int col =
event.nextColTag();
1994 resolved[iBeam].col(col);
1995 usedCol[iBeam] =
true;
2002 infoPtr->errorMsg(
"Error in BeamParticle::findSingleCol: "
2003 "could not find matching anti colour");
2011 void BeamParticle::updateCol(vector<pair<int,int> > colourChanges) {
2013 for (
int iCol = 0;iCol < int(colourChanges.size()); ++iCol) {
2014 int oldCol = colourChanges[iCol].first;
2015 int newCol = colourChanges[iCol].second;
2018 for (
int i = 0;i < int(cols.size()); ++i)
2019 if (cols[i] == oldCol) cols[i] = newCol;
2020 for (
int i = 0;i < int(acols.size()); ++i)
2021 if (acols[i] == oldCol) acols[i] = newCol;
2024 for (
int i = 0;i < int(resolved.size()); ++i) {
2025 if (resolved[i].acol() == oldCol) resolved[i].acol(newCol);
2026 if (resolved[i].col() == oldCol) resolved[i].col(newCol);
2034 void BeamParticle::updateSingleCol(
int oldCol,
int newCol) {
2037 for (
int i = 0; i < int(cols.size()); ++i)
2038 if (cols[i] == oldCol) cols[i] = newCol;
2040 for (
int i = 0; i < int(acols.size()); ++i)
2041 if (acols[i] == oldCol) acols[i] = newCol;
2044 for (
int i = 0; i < int(resolved.size()); ++i) {
2045 if (resolved[i].acol() == oldCol) resolved[i].acol(newCol);
2046 if (resolved[i].col() == oldCol) resolved[i].col(newCol);
2049 colUpdates.push_back(make_pair(oldCol,newCol));