11 #include "Pythia8/JunctionSplitting.h"
25 const double JunctionSplitting::JJSTRINGM2MAX = 25.;
26 const double JunctionSplitting::JJSTRINGM2FRAC = 0.1;
29 const double JunctionSplitting::CONVJNREST = 1e-5;
30 const int JunctionSplitting::NTRYJNREST = 20;
33 const double JunctionSplitting::MTHAD = 0.9;
36 const double JunctionSplitting::MINANGLE = 1e-7;
42 void JunctionSplitting::init( Info* infoPtrIn, Settings& settings,
43 Rndm* rndmPtrIn, ParticleData* particleDataPtrIn) {
49 colTrace.init(infoPtrIn);
50 stringLength.init(infoPtrIn, settings);
53 flavSel.init(settings, particleDataPtrIn, rndmPtr, infoPtr);
54 pTSel.init( settings, particleDataPtrIn, rndmPtr, infoPtr);
55 zSel.init( settings, *particleDataPtrIn, rndmPtr, infoPtr);
58 stringFrag.init(infoPtr, settings, particleDataPtrIn, rndmPtr,
59 &flavSel, &pTSel, &zSel);
62 eNormJunction = settings.parm(
"StringFragmentation:eNormJunction");
63 allowDoubleJunRem = settings.flag(
"ColourReconnection:allowDoubleJunRem");
71 bool JunctionSplitting::checkColours(
Event& event) {
74 for (
int i = 0; i <
event.size(); ++i)
75 if (abs(event[i].px()) >= 0. && abs(event[i].py()) >= 0.
76 && abs(event[i].pz()) >= 0. && abs(event[i].e()) >= 0.
77 && abs(event[i].m()) >= 0.);
79 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
80 "not-a-number energy/momentum/mass");
85 for (
int i = 0; i <
event.size(); ++i) {
86 if (event[i].isFinal() &&
event[i].col() != 0 &&
87 event[i].col() ==
event[i].acol()) {
88 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
89 "Made a gluon colour singlet; redoing colours");
95 colTrace.setupColList(event);
97 vector<vector <int > > iPartonJun, iPartonAntiJun;
98 getPartonLists(event, iPartonJun, iPartonAntiJun);
101 if (!splitJunGluons(event, iPartonJun, iPartonAntiJun) ) {
102 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
103 "Not possible to split junctions; making new colours");
108 if (!splitJunChains(event) ) {
109 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
110 "Not possible to split junctions; making new colours");
115 getPartonLists(event, iPartonJun, iPartonAntiJun);
116 if (!splitJunPairs(event, iPartonJun, iPartonAntiJun) ) {
117 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
118 "Not possible to split junctions; making new colours");
133 bool JunctionSplitting::splitJunGluons(
Event& event,
134 vector<vector< int > >& iPartonJun, vector<vector< int > >& iPartonAntiJun) {
137 for (
int iJun = 0; iJun < int(iPartonJun.size()); ++iJun) {
140 vector<vector <int> > iJunLegs;
143 for (
int i = 0; i < int(iPartonJun[iJun].size()); ++i) {
144 if ( iPartonJun[iJun][i]/10 == iPartonJun[iJun][0]/10) ++leg;
145 iJunLegs[leg].push_back(iPartonJun[iJun][i]);
149 for (
int i = 0;i < int(iJunLegs.size()); ++i) {
152 if (iJunLegs[i].back() > 0)
154 int identJun = iJunLegs[i][0];
156 if (iJunLegs[i].size() == 2)
159 int identAntiJun = 0, iAntiLeg = -1;
162 int colQ = 0, acolQ = 0;
163 int idQ = flavSel.pickLightQ();
167 if ( iJunLegs[i].size() == 3) {
170 if (event[ iJunLegs[i][1] ].idAbs() != 21)
continue;
173 colQ =
event[ iJunLegs[i][1] ].col();
174 acolQ =
event[ iJunLegs[i][1] ].acol();
175 Vec4 pQ = 0.5 *
event[ iJunLegs[i][1] ].p();
176 double mQ = 0.5 *
event[ iJunLegs[i][1] ].m();
177 int iQ =
event.append( idQ, 75, iJunLegs[i][1], 0, 0, 0, colQ, 0,
179 int iQbar =
event.append( -idQ, 75, iJunLegs[i][1], 0, 0, 0, 0, acolQ,
183 event[ iJunLegs[i][1] ].statusNeg();
184 event[ iJunLegs[i][1] ].daughters( iQ, iQbar);
187 identAntiJun = iJunLegs[i].back();
188 int iOld = iJunLegs[i][1];
189 bool erasing =
false;
190 for (
int j = 0; j < int(iPartonJun[iJun].size()); ++j) {
191 if (iPartonJun[iJun][j] == iOld)
193 if (iPartonJun[iJun][j] == identAntiJun) {
194 iPartonJun[iJun][j] = iQ;
198 iPartonJun[iJun].erase(iPartonJun[iJun].begin() + j);
205 for (
int j = 0; j < int(iPartonAntiJun.size()); j++)
206 if ( iPartonAntiJun[j][0]/10 == identAntiJun/10) {
211 if (iAntiJun == -1) {
212 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
213 "Something went wrong in finding antijunction");
219 for (
int j = 0; j < int(iPartonAntiJun[iAntiJun].size()); ++j) {
220 if ( iPartonAntiJun[iAntiJun][j] / 10 == identAntiJun / 10)
222 if ( iPartonAntiJun[iAntiJun][j] == identJun) {
223 iPartonAntiJun[iAntiJun][j] = iQbar;
229 else if (iJunLegs[i].size() > 3) {
232 vector<double> m2Pair;
234 for (
int j = 1; j < int(iJunLegs[i].size()) - 2; ++j) {
235 double m2Now = 0.5 *
event[ iJunLegs[i][j] ].p()
236 *
event[ iJunLegs[i][j + 1] ].p();
237 m2Pair.push_back(m2Now);
242 double m2Reg = m2Sum * rndmPtr->flat();
244 do m2Reg -= m2Pair[++iReg];
245 while (m2Reg > 0. && iReg <
int(m2Pair.size()) - 1);
246 m2Reg = m2Pair[iReg];
252 double m2Temp = min( JJSTRINGM2MAX, JJSTRINGM2FRAC * m2Reg);
256 double zTemp = zSel.zFrag( idQ, 0, m2Temp);
258 xNeg = m2Temp / (zTemp * m2Reg);
260 if (rndmPtr->flat() > 0.5) swap(xPos, xNeg);
263 if ( event[ iJunLegs[i][iReg] ].idAbs() != 21
264 || event[ iJunLegs[i][iReg + 1] ].idAbs() != 21 )
continue;
267 Particle& gJun =
event[ iJunLegs[i][iReg] ];
268 Particle& gAnti =
event[ iJunLegs[i][iReg + 1] ];
271 int dau1 =
event.size();
272 gJun.daughters(dau1, dau1 + 3);
273 gAnti.daughters(dau1, dau1 + 3);
274 int mother1 = min( iJunLegs[i][iReg], iJunLegs[i][iReg + 1]);
275 int mother2 = max( iJunLegs[i][iReg], iJunLegs[i][iReg + 1]);
279 int gJunCol = gJun.col();
280 int gJunAcol = gJun.acol();
281 int gAntiAcol = gAnti.acol();
282 Vec4 gJunP = gJun.p();
283 Vec4 gAntiP = gAnti.p();
284 double gJunM = gJun.m();
285 double gAntiM = gAnti.m();
289 acolQ =
event.nextColTag();
292 int iGjun =
event.append( 21, 75, mother1, mother2, 0, 0,
293 gJunCol, gJunAcol, (1. - 0.5 * xPos) * gJunP,
294 (1. - 0.5 * xPos) * gJunM);
295 event.append( 21, 75, mother1, mother2, 0, 0,
296 acolQ, gAntiAcol, (1. - 0.5 * xNeg) * gAntiP,
297 (1. - 0.5 * xNeg) * gAntiM);
300 int iQ =
event.append( idQ, 75, mother1, mother2, 0, 0,
301 colQ, 0, 0.5 * xNeg * gAntiP, 0.5 * xNeg * gAntiM );
302 int iQbar =
event.append( -idQ, 75, mother1, mother2, 0, 0,
303 0, acolQ, 0.5 * xPos * gJunP, 0.5 * xPos * gJunM );
306 identAntiJun = iJunLegs[i].back();
307 bool erasing =
false;
308 for (
int j = 0; j < int(iPartonJun[iJun].size()); ++j) {
309 if (iPartonJun[iJun][j] == mother1 ||
310 iPartonJun[iJun][j] == mother2)
313 if ( iPartonJun[iJun][j] == identAntiJun) {
314 iPartonJun[iJun][j] = iQ;
315 iPartonJun[iJun].insert(iPartonJun[iJun].begin() + j, iGjun);
319 iPartonJun[iJun].erase(iPartonJun[iJun].begin() + j);
326 for (
int j = 0; j < int(iPartonAntiJun.size());j++)
327 if ( iPartonAntiJun[j][0]/10 == identAntiJun/10) {
332 if (iAntiJun == -1) {
333 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
334 "Something went wrong in finding antijunction");
339 for (
int j = 0; j < int(iPartonAntiJun[iAntiJun].size()); ++j) {
340 if ( iPartonAntiJun[iAntiJun][j] / 10 == identAntiJun / 10)
342 if (iPartonAntiJun[iAntiJun][j] == identJun) {
343 iPartonAntiJun[iAntiJun][j] = iQbar;
350 event.endColJunction((-identJun)/10 - 1, i, colQ);
351 event.endColJunction((-identAntiJun)/10 - 1, iAntiLeg, acolQ);
365 bool JunctionSplitting::splitJunChains(
Event& event) {
368 event.saveJunctionSize();
369 vector<vector<int> > junChains = colTrace.getJunChains(event);
373 for (
int i = 0; i < int(junChains.size()); ++i) {
374 if (junChains[i].size() < 3)
377 vector<int> cols, acols;
378 for (
int j = 0; j < int(junChains[i].size()); ++j) {
380 junRem.push_back(junChains[i][j]);
381 if (event.kindJunction(junChains[i][j]) % 2 == 0)
382 for (
int jLeg = 0; jLeg < 3; ++jLeg)
383 acols.push_back(event.colJunction(junChains[i][j],jLeg));
385 for (
int jLeg = 0; jLeg < 3; ++jLeg)
386 cols.push_back(event.colJunction(junChains[i][j],jLeg));
389 for (
int j = 0; j < int(cols.size()); ++j)
390 for (
int k = 0; k < int(acols.size()); ++k)
391 if (cols[j] == acols[k]) {
392 cols.erase(cols.begin() + j);
393 acols.erase(acols.begin() + k);
399 while (cols.size() > acols.size()) {
400 int i1 = int(rndmPtr->flat() *cols.size());
402 cols.erase(cols.begin() + i1);
403 int i2 = int(rndmPtr->flat() *cols.size());
405 cols.erase(cols.begin() + i2);
406 int i3 = int(rndmPtr->flat() *cols.size());
408 cols.erase(cols.begin() + i3);
409 event.appendJunction(1, col1, col2, col3);
413 while (acols.size() > cols.size()) {
414 int i1 = int(rndmPtr->flat() *acols.size());
415 int acol1 = acols[i1];
416 acols.erase(acols.begin() + i1);
417 int i2 = int(rndmPtr->flat() *acols.size());
418 int acol2 = acols[i2];
419 acols.erase(acols.begin() + i2);
420 int i3 = int(rndmPtr->flat() *acols.size());
421 int acol3 = acols[i3];
422 acols.erase(acols.begin() + i3);
423 event.appendJunction(2,acol1,acol2,acol3);
428 while (
int(acols.size()) > 1) {
429 int i1 = int(rndmPtr->flat() *cols.size());
431 cols.erase(cols.begin() + i1);
432 int i2 = int(rndmPtr->flat() *cols.size());
434 cols.erase(cols.begin() + i2);
435 int i3 = int(rndmPtr->flat() *acols.size());
436 int acol1 = acols[i3];
437 acols.erase(acols.begin() + i3);
438 int i4 = int(rndmPtr->flat() *acols.size());
439 int acol2 = acols[i4];
440 acols.erase(acols.begin() + i4);
441 int newCol =
event.nextColTag();
442 event.appendJunction(1, col1, col2, newCol);
443 event.appendJunction(2, acol1, acol2, newCol);
447 if (
int(acols.size()) == 1) {
449 for (
int iPar = 0; iPar <
event.size(); ++iPar)
450 if (event[iPar].isFinal() &&
event[iPar].col() == cols[0])
453 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
454 "Splitting multiple directly connected junctions failed");
457 int iNew =
event.copy(iCol, 76);
458 event[iNew].col(acols[0]);
463 sort(junRem.begin(),junRem.end());
464 reverse(junRem.begin(),junRem.end());
465 for (
int i = 0; i < int(junRem.size()); ++i)
466 event.eraseJunction(junRem[i]);
467 event.saveJunctionSize();
479 bool JunctionSplitting::splitJunPairs(
Event& event,
480 vector<vector< int > >& iPartonJun, vector<vector< int > >& iPartonAntiJun) {
483 event.saveJunctionSize();
487 vector<vector<int> > junChains = colTrace.getJunChains(event);
489 for (
int i = 0; i < int(junChains.size()); ++i) {
490 if (junChains[i].size() == 2) {
491 vector<pair<int,int> > matchedLegs;
492 for (
int j = 0; j < 3; ++j)
493 for (
int k = 0; k < 3; ++k)
494 if (event.colJunction(junChains[i][0],j) ==
495 event.colJunction(junChains[i][1],k))
496 matchedLegs.push_back(make_pair(j,k));
500 if (matchedLegs.size() == 3) {
501 junRem.push_back(junChains[i][0]);
502 junRem.push_back(junChains[i][1]);
507 else if (matchedLegs.size() == 2) {
511 if (matchedLegs[0].first != 1 && matchedLegs[1].first != 1) i1 = 1;
512 if (matchedLegs[0].first != 2 && matchedLegs[1].first != 2) i1 = 2;
516 if (matchedLegs[0].second != 1 && matchedLegs[1].second != 1) j1 = 1;
517 if (matchedLegs[0].second != 2 && matchedLegs[1].second != 2) j1 = 2;
520 int col =
event.colJunction(junChains[i][0],i1);
521 int acol =
event.colJunction(junChains[i][1],j1);
522 if (event.kindJunction(junChains[i][1]) % 2 == 1)
527 for (
int j = 0;j <
event.size();++j)
528 if (event[j].isFinal() &&
event[j].acol() == acol) {
533 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
534 "Anti colour not found when combining two junctions to a string");
539 int iNew =
event.copy(iAcol,66);
540 event[iNew].acol(col);
543 junRem.push_back(junChains[i][0]);
544 junRem.push_back(junChains[i][1]);
548 else if (matchedLegs.size() == 1) {
551 int iJun = junChains[i][0];
552 int iAnti = junChains[i][1];
553 int iLeg = matchedLegs[0].first;
554 int iAntiLeg = matchedLegs[0].second;
555 if (event.kindJunction(iAnti) % 2 == 1) {
557 swap(iLeg, iAntiLeg);
561 int iJunList = -1, iAntiList = -1;
562 for (
int l = 0;l < int(iPartonJun.size()); ++l)
563 if (- iPartonJun[l][0]/10 - 1 == iJun) {
567 for (
int l = 0;l < int(iPartonAntiJun.size()); ++l)
568 if (- iPartonAntiJun[l][0]/10 - 1 == iAnti) {
572 if (iJunList == -1 || iAntiList == -1) {
573 infoPtr->errorMsg(
"Error in JunctionSplitting::SplitJunChain:"
574 " failed to find junctions in the parton list");
579 vector<vector <int> > iJunLegs;
583 for (
int l = 0; l < int(iPartonJun[iJunList].size()); ++l) {
584 if ( iPartonJun[iJunList][l]/10 == iPartonJun[iJunList][0]/10) ++leg;
585 iJunLegs[leg].push_back(iPartonJun[iJunList][l]);
589 vector<vector <int> > iAntiLegs;
592 for (
int l = 0; l < int(iPartonAntiJun[iAntiList].size()); ++l) {
593 if ( iPartonAntiJun[iAntiList][l]/10
594 == iPartonAntiJun[iAntiList][0]/10) ++leg;
595 iAntiLegs[leg].push_back(iPartonAntiJun[iAntiList][l]);
599 vector<int>& iJunLeg0 = (iLeg == 0) ? iJunLegs[1] : iJunLegs[0];
600 vector<int>& iJunLeg1 = (iLeg == 2) ? iJunLegs[1] : iJunLegs[2];
601 vector<int>& iAntiLeg0 = (iAntiLeg == 0) ? iAntiLegs[1] : iAntiLegs[0];
602 vector<int>& iAntiLeg1 = (iAntiLeg == 2) ? iAntiLegs[1] : iAntiLegs[2];
606 if (iAntiLeg0[1] < 0 || iAntiLeg1[1] < 0)
continue;
609 Vec4 pJunLeg0 =
event[ iJunLeg0[1] ].p();
610 Vec4 pJunLeg1 =
event[ iJunLeg1[1] ].p();
611 Vec4 pAntiLeg0 =
event[ iAntiLeg0[1] ].p();
612 Vec4 pAntiLeg1 =
event[ iAntiLeg1[1] ].p();
615 if ( theta(pJunLeg0,pJunLeg1) < MINANGLE
616 || theta(pAntiLeg0,pAntiLeg1) < MINANGLE
617 || theta(pJunLeg0,pAntiLeg0) < MINANGLE
618 || theta(pJunLeg0,pAntiLeg1) < MINANGLE
619 || theta(pJunLeg1,pAntiLeg0) < MINANGLE
620 || theta(pJunLeg1,pAntiLeg1) < MINANGLE) {
621 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
622 "parallel junction state not allowed.");
627 Vec4 pStart = pJunLeg0 / pJunLeg0.e() + pJunLeg1 / pJunLeg1.e()
628 + pAntiLeg0 / pAntiLeg0.e() + pAntiLeg1 / pAntiLeg1.e();
631 RotBstMatrix MtoJRF, MtoARF;
632 Vec4 pInJRF[3], pInARF[3];
633 for (
int iJunLocal = 0; iJunLocal < 2; ++iJunLocal) {
634 int offset = (iJunLocal == 0) ? 0 : 2;
637 RotBstMatrix MtoRF, Mstep;
638 MtoRF.bstback(pStart);
646 pInRF[0 + offset] = pJunLeg0;
647 pInRF[1 + offset] = pJunLeg1;
648 pInRF[2 - offset] = pAntiLeg0;
649 pInRF[3 - offset] = pAntiLeg1;
650 for (
int l = 0; l < 4; ++l) pInRF[l].rotbst(MtoRF);
653 double wt2 = 1. - exp( -pInRF[2].e() / eNormJunction);
654 double wt3 = 1. - exp( -pInRF[3].e() / eNormJunction);
655 pInRF[2] = wt2 * pInRF[2] + wt3 * pInRF[3];
658 Mstep = stringFrag.junctionRestFrame( pInRF[0], pInRF[1], pInRF[2]);
659 MtoRF.rotbst( Mstep );
660 }
while (iter < 3 || (Mstep.deviation() > CONVJNREST
661 && iter < NTRYJNREST) );
664 if (iJunLocal == 0) {
666 for (
int l = 0; l < 3; ++l) pInJRF[l] = pInRF[l];
669 for (
int l = 0; l < 3; ++l) pInARF[l] = pInRF[l];
674 RotBstMatrix MfromJRF = MtoJRF;
676 RotBstMatrix MfromARF = MtoARF;
680 Vec4 vJun(0., 0., 0., 1.);
681 vJun.rotbst(MfromJRF);
682 Vec4 vAnti(0., 0., 0., 1.);
683 vAnti.rotbst(MfromARF);
684 Vec4 pLabJ[3], pLabA[3];
685 for (
int l = 0; l < 3; ++l) {
686 pLabJ[l] = pInJRF[l];
687 pLabJ[l].rotbst(MfromJRF);
688 pLabA[l] = pInARF[l];
689 pLabA[l].rotbst(MfromARF);
693 double vJvA = vJun * vAnti;
694 double vJvAe2y = vJvA + sqrt(vJvA*vJvA - 1.);
695 double lambdaJA = stringLength.getJuncLength(pInJRF[0], pInJRF[1],
696 pInARF[0], pInARF[1]);
698 double lambda00 = stringLength.getStringLength(pLabJ[0], pLabA[0]) +
699 stringLength.getStringLength(pLabJ[1], pLabA[1]);
701 double lambda01 = stringLength.getStringLength(pLabJ[0], pLabA[1]) +
702 stringLength.getStringLength(pLabJ[1], pLabA[0]);
705 if (lambdaJA > min( lambda00, lambda01) && allowDoubleJunRem) {
708 int iCol1 = iJunLeg0[1];
709 int iCol2 = iJunLeg1[1];
710 int iAcol1 = iAntiLeg0[1];
711 int iAcol2 = iAntiLeg1[1];
712 if (lambda00 > lambda01)
713 swap(iAcol1, iAcol2);
716 int iNew1 =
event.copy(iAcol1, 66);
717 event[iNew1].acol(event[iCol1].col());
718 int iNew2 =
event.copy(iAcol2, 66);
719 event[iNew2].acol(event[iCol2].col());
720 junRem.push_back(junChains[i][0]);
721 junRem.push_back(junChains[i][1]);
732 int idJ0Abs =
event[ iJunLeg0[1] ].idAbs();
733 bool acceptJ0 = idJ0Abs < 6 || idJ0Abs == 21 ||
734 (idJ0Abs > 1000 && idJ0Abs < 6000 && (idJ0Abs / 10) % 10 == 0);
735 int idJ1Abs =
event[ iJunLeg1[1] ].idAbs();
736 bool acceptJ1 = idJ1Abs < 6 || idJ1Abs == 21 ||
737 (idJ1Abs > 1000 && idJ1Abs < 6000 && (idJ1Abs / 10) % 10 == 0);
738 int idA0Abs =
event[ iAntiLeg0[1] ].idAbs();
739 bool acceptA0 = idA0Abs < 6 || idA0Abs == 21 ||
740 (idA0Abs > 1000 && idA0Abs < 6000 && (idA0Abs / 10) % 10 == 0);
741 int idA1Abs =
event[ iAntiLeg1[1] ].idAbs();
742 bool acceptA1 = idA1Abs < 6 || idA1Abs == 21 ||
743 (idA1Abs > 1000 && idA1Abs < 6000 && (idA1Abs / 10) % 10 == 0);
745 if ( !(acceptJ0 || acceptJ1)) {
746 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
747 "No light quarks available in junction split.");
751 if ( !(acceptA0 || acceptA1)) {
752 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
753 "No light quarks available in junction split.");
757 double eShift = MTHAD / (3. * sqrt(vJvAe2y));
758 double fracJ0 = 0, fracJ1 = 0, fracA0 = 0, fracA1 = 0;
760 fracJ0 = min(0.5, eShift / pInJRF[0].e());
762 fracJ1 = min(0.5, eShift / pInJRF[1].e());
763 Vec4 pFromJun = fracJ0 * pJunLeg0 + fracJ1 * pJunLeg1;
765 fracA0 = min(0.5, eShift / pInARF[0].e());
767 fracA1 = min(0.5, eShift / pInARF[1].e());
768 Vec4 pFromAnti = fracA0 * pAntiLeg0 + fracA1 * pAntiLeg1;
771 int idQ = flavSel.pickLightQ();
773 int junMother1 = min(iJunLeg0[1], iJunLeg1[1]);
774 int junMother2 = max(iJunLeg0[1], iJunLeg1[1]);
775 int antiMother1 = min(iAntiLeg0[1], iAntiLeg1[1]);
776 int antiMother2 = max(iAntiLeg0[1], iAntiLeg1[1]);
779 int iJunNew1 =
event.copy(iJunLeg0[1], 76);
780 event[iJunNew1].rescale5(1. - fracJ0);
781 iJunLeg0[1] = iJunNew1;
782 event[iJunNew1].mothers(junMother2, junMother1);
784 int iJunNew2 =
event.copy(iJunLeg1[1], 76);
785 event[iJunNew2].rescale5(1. - fracJ1);
786 iJunLeg1[1] = iJunNew2;
787 event[iJunNew2].mothers(junMother2, junMother1);
791 int acolQ =
event.nextColTag();
792 event.endColJunction(iAnti, iAntiLeg, acolQ);
793 event.colJunction(iAnti, iAntiLeg, acolQ);
794 int iNewA =
event.append( -idQ, 76, junMother2, junMother1, 0, 0,
795 0, acolQ, pFromJun, pFromJun.mCalc() );
798 int iAntiNew1 =
event.copy(iAntiLeg0[1], 76);
799 event[iAntiNew1].rescale5(1. - fracA0);
800 iAntiLeg0[1] = iAntiNew1;
801 event[iAntiNew1].mothers(antiMother2, antiMother1);
803 int iAntiNew2 =
event.copy(iAntiLeg1[1], 76);
804 event[iAntiNew2].rescale5(1. - fracA1);
805 iAntiLeg1[1] = iAntiNew2;
806 event[iAntiNew2].mothers(antiMother2, antiMother1);
809 int colQ =
event.nextColTag();
810 event.endColJunction(iJun, iLeg, colQ);
811 event.colJunction(iJun, iLeg, colQ);
812 int iNewJ =
event.append( idQ, 76, antiMother2, antiMother1, 0, 0,
813 colQ, 0, pFromAnti, pFromAnti.mCalc() );
816 event[
event[iJunNew1].mother1()].daughters( iJunNew1, iNewA);
817 event[
event[iJunNew1].mother2()].daughters( iJunNew1, iNewA);
818 event[
event[iAntiNew1].mother1()].daughters( iAntiNew1, iNewJ);
819 event[
event[iAntiNew1].mother2()].daughters( iAntiNew1, iNewJ);
827 sort(junRem.begin(),junRem.end());
828 reverse(junRem.begin(),junRem.end());
829 for (
int i = 0; i < int(junRem.size()); ++i)
830 event.eraseJunction(junRem[i]);
831 event.saveJunctionSize();
843 bool JunctionSplitting::getPartonLists(
Event& event,
844 vector<vector< int > > & iPartonJun, vector<vector<int > >& iPartonAntiJun) {
847 colTrace.setupColList(event);
850 iPartonAntiJun.clear();
856 for (
int iLoop = 0; iLoop < 2*
event.sizeJunction(); ++iLoop) {
858 int iJun = iLoop %
event.sizeJunction();
859 if ( !event.remainsJunction(iJun))
continue;
862 int kindJun =
event.kindJunction(iJun);
863 if ( iLoop < event.sizeJunction() && kindJun%2 == 0)
continue;
864 else if ( iLoop >= event.sizeJunction() && kindJun%2 == 1)
continue;
872 for (
int iCol = 0; iCol < 3; ++iCol) {
873 int indxCol =
event.colJunction(iJun, iCol);
874 iParton.push_back( -(10 + 10 * iJun + iCol) );
876 if ( kindJun%2 == 1 && !colTrace.traceFromAcol(indxCol,event, iJun,
877 iCol, iParton) )
return false;
879 else if ( kindJun%2 == 0 && !colTrace.traceFromCol(indxCol,event, iJun,
880 iCol, iParton) )
return false;
885 for (
int i = 0; i < int(iParton.size()); ++i)
if (iParton[i] < 0) ++nNeg;
888 if ( kindJun%2 == 1 ) iPartonJun.push_back(iParton);
889 else iPartonAntiJun.push_back(iParton);
903 bool JunctionSplitting::setAcol(
Event& event,
int col,
int acol) {
906 for (
int j = 0;j <
event.size(); ++j)
907 if (event[j].isFinal() &&
event[j].acol() == acol) {
908 int iNew =
event.copy(j,66);
909 event[iNew].acol(col);
913 for (
int j = 0;j <
event.sizeJunction(); ++j)
914 for (
int jLeg = 0;jLeg < 3; ++jLeg)
915 if (event.colJunction(j, jLeg) == acol) {
916 event.colJunction(j, jLeg, col);
921 infoPtr->errorMsg(
"Warning in JunctionSplitting::setAcol:"
922 "Anti colour not found when combing two junctions to a string");