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() {
45 colTrace.init(infoPtr);
46 stringLength.init(infoPtr, *settingsPtr);
54 stringFrag.init(&flavSel, &pTSel, &zSel);
57 eNormJunction = parm(
"StringFragmentation:eNormJunction");
58 allowDoubleJunRem = flag(
"ColourReconnection:allowDoubleJunRem");
66 bool JunctionSplitting::checkColours(
Event& event) {
69 for (
int i = 0; i <
event.size(); ++i)
70 if (abs(event[i].px()) >= 0. && abs(event[i].py()) >= 0.
71 && abs(event[i].pz()) >= 0. && abs(event[i].e()) >= 0.
72 && abs(event[i].m()) >= 0.);
74 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
75 "not-a-number energy/momentum/mass");
80 for (
int i = 0; i <
event.size(); ++i) {
81 if (event[i].isFinal() &&
event[i].col() != 0 &&
82 event[i].col() ==
event[i].acol()) {
83 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
84 "Made a gluon colour singlet; redoing colours");
90 colTrace.setupColList(event);
92 vector<vector <int > > iPartonJun, iPartonAntiJun;
93 getPartonLists(event, iPartonJun, iPartonAntiJun);
96 if (!splitJunGluons(event, iPartonJun, iPartonAntiJun) ) {
97 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
98 "Not possible to split junctions; making new colours");
103 if (!splitJunChains(event) ) {
104 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
105 "Not possible to split junctions; making new colours");
110 getPartonLists(event, iPartonJun, iPartonAntiJun);
111 if (!splitJunPairs(event, iPartonJun, iPartonAntiJun) ) {
112 infoPtr->errorMsg(
"Warning in JunctionSplitting::CheckColours: "
113 "Not possible to split junctions; making new colours");
128 bool JunctionSplitting::splitJunGluons(
Event& event,
129 vector<vector< int > >& iPartonJun, vector<vector< int > >& iPartonAntiJun) {
132 for (
int iJun = 0; iJun < int(iPartonJun.size()); ++iJun) {
135 vector<vector <int> > iJunLegs;
138 for (
int i = 0; i < int(iPartonJun[iJun].size()); ++i) {
139 if ( iPartonJun[iJun][i]/10 == iPartonJun[iJun][0]/10) ++leg;
140 iJunLegs[leg].push_back(iPartonJun[iJun][i]);
144 for (
int i = 0;i < int(iJunLegs.size()); ++i) {
147 if (iJunLegs[i].back() > 0)
149 int identJun = iJunLegs[i][0];
151 if (iJunLegs[i].size() == 2)
154 int identAntiJun = 0, iAntiLeg = -1;
157 int colQ = 0, acolQ = 0;
158 int idQ = flavSel.pickLightQ();
162 if ( iJunLegs[i].size() == 3) {
165 if (event[ iJunLegs[i][1] ].idAbs() != 21)
continue;
168 colQ =
event[ iJunLegs[i][1] ].col();
169 acolQ =
event[ iJunLegs[i][1] ].acol();
170 Vec4 pQ = 0.5 *
event[ iJunLegs[i][1] ].p();
171 double mQ = 0.5 *
event[ iJunLegs[i][1] ].m();
172 int iQ =
event.append( idQ, 75, iJunLegs[i][1], 0, 0, 0, colQ, 0,
174 int iQbar =
event.append( -idQ, 75, iJunLegs[i][1], 0, 0, 0, 0, acolQ,
178 event[ iJunLegs[i][1] ].statusNeg();
179 event[ iJunLegs[i][1] ].daughters( iQ, iQbar);
182 identAntiJun = iJunLegs[i].back();
183 int iOld = iJunLegs[i][1];
184 bool erasing =
false;
185 for (
int j = 0; j < int(iPartonJun[iJun].size()); ++j) {
186 if (iPartonJun[iJun][j] == iOld)
188 if (iPartonJun[iJun][j] == identAntiJun) {
189 iPartonJun[iJun][j] = iQ;
193 iPartonJun[iJun].erase(iPartonJun[iJun].begin() + j);
200 for (
int j = 0; j < int(iPartonAntiJun.size()); j++)
201 if ( iPartonAntiJun[j][0]/10 == identAntiJun/10) {
206 if (iAntiJun == -1) {
207 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
208 "Something went wrong in finding antijunction");
214 for (
int j = 0; j < int(iPartonAntiJun[iAntiJun].size()); ++j) {
215 if ( iPartonAntiJun[iAntiJun][j] / 10 == identAntiJun / 10)
217 if ( iPartonAntiJun[iAntiJun][j] == identJun) {
218 iPartonAntiJun[iAntiJun][j] = iQbar;
224 else if (iJunLegs[i].size() > 3) {
227 vector<double> m2Pair;
229 for (
int j = 1; j < int(iJunLegs[i].size()) - 2; ++j) {
230 double m2Now = 0.5 *
event[ iJunLegs[i][j] ].p()
231 *
event[ iJunLegs[i][j + 1] ].p();
232 m2Pair.push_back(m2Now);
237 double m2Reg = m2Sum * rndmPtr->flat();
239 do m2Reg -= m2Pair[++iReg];
240 while (m2Reg > 0. && iReg <
int(m2Pair.size()) - 1);
241 m2Reg = m2Pair[iReg];
247 double m2Temp = min( JJSTRINGM2MAX, JJSTRINGM2FRAC * m2Reg);
251 double zTemp = zSel.zFrag( idQ, 0, m2Temp);
253 xNeg = m2Temp / (zTemp * m2Reg);
255 if (rndmPtr->flat() > 0.5) swap(xPos, xNeg);
258 if ( event[ iJunLegs[i][iReg] ].idAbs() != 21
259 || event[ iJunLegs[i][iReg + 1] ].idAbs() != 21 )
continue;
262 Particle& gJun =
event[ iJunLegs[i][iReg] ];
263 Particle& gAnti =
event[ iJunLegs[i][iReg + 1] ];
266 int dau1 =
event.size();
267 gJun.daughters(dau1, dau1 + 3);
268 gAnti.daughters(dau1, dau1 + 3);
269 int mother1 = min( iJunLegs[i][iReg], iJunLegs[i][iReg + 1]);
270 int mother2 = max( iJunLegs[i][iReg], iJunLegs[i][iReg + 1]);
274 int gJunCol = gJun.col();
275 int gJunAcol = gJun.acol();
276 int gAntiAcol = gAnti.acol();
277 Vec4 gJunP = gJun.p();
278 Vec4 gAntiP = gAnti.p();
279 double gJunM = gJun.m();
280 double gAntiM = gAnti.m();
284 acolQ =
event.nextColTag();
287 int iGjun =
event.append( 21, 75, mother1, mother2, 0, 0,
288 gJunCol, gJunAcol, (1. - 0.5 * xPos) * gJunP,
289 (1. - 0.5 * xPos) * gJunM);
290 event.append( 21, 75, mother1, mother2, 0, 0,
291 acolQ, gAntiAcol, (1. - 0.5 * xNeg) * gAntiP,
292 (1. - 0.5 * xNeg) * gAntiM);
295 int iQ =
event.append( idQ, 75, mother1, mother2, 0, 0,
296 colQ, 0, 0.5 * xNeg * gAntiP, 0.5 * xNeg * gAntiM );
297 int iQbar =
event.append( -idQ, 75, mother1, mother2, 0, 0,
298 0, acolQ, 0.5 * xPos * gJunP, 0.5 * xPos * gJunM );
301 identAntiJun = iJunLegs[i].back();
302 bool erasing =
false;
303 for (
int j = 0; j < int(iPartonJun[iJun].size()); ++j) {
304 if (iPartonJun[iJun][j] == mother1 ||
305 iPartonJun[iJun][j] == mother2)
308 if ( iPartonJun[iJun][j] == identAntiJun) {
309 iPartonJun[iJun][j] = iQ;
310 iPartonJun[iJun].insert(iPartonJun[iJun].begin() + j, iGjun);
314 iPartonJun[iJun].erase(iPartonJun[iJun].begin() + j);
321 for (
int j = 0; j < int(iPartonAntiJun.size());j++)
322 if ( iPartonAntiJun[j][0]/10 == identAntiJun/10) {
327 if (iAntiJun == -1) {
328 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
329 "Something went wrong in finding antijunction");
334 for (
int j = 0; j < int(iPartonAntiJun[iAntiJun].size()); ++j) {
335 if ( iPartonAntiJun[iAntiJun][j] / 10 == identAntiJun / 10)
337 if (iPartonAntiJun[iAntiJun][j] == identJun) {
338 iPartonAntiJun[iAntiJun][j] = iQbar;
345 event.endColJunction((-identJun)/10 - 1, i, colQ);
346 event.endColJunction((-identAntiJun)/10 - 1, iAntiLeg, acolQ);
360 bool JunctionSplitting::splitJunChains(
Event& event) {
363 event.saveJunctionSize();
364 vector<vector<int> > junChains = colTrace.getJunChains(event);
368 for (
int i = 0; i < int(junChains.size()); ++i) {
369 if (junChains[i].size() < 3)
372 vector<int> cols, acols;
373 for (
int j = 0; j < int(junChains[i].size()); ++j) {
375 junRem.push_back(junChains[i][j]);
376 if (event.kindJunction(junChains[i][j]) % 2 == 0)
377 for (
int jLeg = 0; jLeg < 3; ++jLeg)
378 acols.push_back(event.colJunction(junChains[i][j],jLeg));
380 for (
int jLeg = 0; jLeg < 3; ++jLeg)
381 cols.push_back(event.colJunction(junChains[i][j],jLeg));
384 for (
int j = 0; j < int(cols.size()); ++j)
385 for (
int k = 0; k < int(acols.size()); ++k)
386 if (cols[j] == acols[k]) {
387 cols.erase(cols.begin() + j);
388 acols.erase(acols.begin() + k);
394 while (cols.size() > acols.size()) {
395 int i1 = int(rndmPtr->flat() *cols.size());
397 cols.erase(cols.begin() + i1);
398 int i2 = int(rndmPtr->flat() *cols.size());
400 cols.erase(cols.begin() + i2);
401 int i3 = int(rndmPtr->flat() *cols.size());
403 cols.erase(cols.begin() + i3);
404 event.appendJunction(1, col1, col2, col3);
408 while (acols.size() > cols.size()) {
409 int i1 = int(rndmPtr->flat() *acols.size());
410 int acol1 = acols[i1];
411 acols.erase(acols.begin() + i1);
412 int i2 = int(rndmPtr->flat() *acols.size());
413 int acol2 = acols[i2];
414 acols.erase(acols.begin() + i2);
415 int i3 = int(rndmPtr->flat() *acols.size());
416 int acol3 = acols[i3];
417 acols.erase(acols.begin() + i3);
418 event.appendJunction(2,acol1,acol2,acol3);
423 while (
int(acols.size()) > 1) {
424 int i1 = int(rndmPtr->flat() *cols.size());
426 cols.erase(cols.begin() + i1);
427 int i2 = int(rndmPtr->flat() *cols.size());
429 cols.erase(cols.begin() + i2);
430 int i3 = int(rndmPtr->flat() *acols.size());
431 int acol1 = acols[i3];
432 acols.erase(acols.begin() + i3);
433 int i4 = int(rndmPtr->flat() *acols.size());
434 int acol2 = acols[i4];
435 acols.erase(acols.begin() + i4);
436 int newCol =
event.nextColTag();
437 event.appendJunction(1, col1, col2, newCol);
438 event.appendJunction(2, acol1, acol2, newCol);
442 if (
int(acols.size()) == 1) {
444 for (
int iPar = 0; iPar <
event.size(); ++iPar)
445 if (event[iPar].isFinal() &&
event[iPar].col() == cols[0])
448 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
449 "Splitting multiple directly connected junctions failed");
452 int iNew =
event.copy(iCol, 76);
453 event[iNew].col(acols[0]);
458 sort(junRem.begin(),junRem.end());
459 reverse(junRem.begin(),junRem.end());
460 for (
int i = 0; i < int(junRem.size()); ++i)
461 event.eraseJunction(junRem[i]);
462 event.saveJunctionSize();
474 bool JunctionSplitting::splitJunPairs(
Event& event,
475 vector<vector< int > >& iPartonJun, vector<vector< int > >& iPartonAntiJun) {
478 event.saveJunctionSize();
482 vector<vector<int> > junChains = colTrace.getJunChains(event);
484 for (
int i = 0; i < int(junChains.size()); ++i) {
485 if (junChains[i].size() == 2) {
486 vector<pair<int,int> > matchedLegs;
487 for (
int j = 0; j < 3; ++j)
488 for (
int k = 0; k < 3; ++k)
489 if (event.colJunction(junChains[i][0],j) ==
490 event.colJunction(junChains[i][1],k))
491 matchedLegs.push_back(make_pair(j,k));
495 if (matchedLegs.size() == 3) {
496 junRem.push_back(junChains[i][0]);
497 junRem.push_back(junChains[i][1]);
502 else if (matchedLegs.size() == 2) {
506 if (matchedLegs[0].first != 1 && matchedLegs[1].first != 1) i1 = 1;
507 if (matchedLegs[0].first != 2 && matchedLegs[1].first != 2) i1 = 2;
511 if (matchedLegs[0].second != 1 && matchedLegs[1].second != 1) j1 = 1;
512 if (matchedLegs[0].second != 2 && matchedLegs[1].second != 2) j1 = 2;
515 int col =
event.colJunction(junChains[i][0],i1);
516 int acol =
event.colJunction(junChains[i][1],j1);
517 if (event.kindJunction(junChains[i][1]) % 2 == 1)
522 for (
int j = 0;j <
event.size();++j)
523 if (event[j].isFinal() &&
event[j].acol() == acol) {
528 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunChain:"
529 "Anti colour not found when combining two junctions to a string");
534 int iNew =
event.copy(iAcol,66);
535 event[iNew].acol(col);
538 junRem.push_back(junChains[i][0]);
539 junRem.push_back(junChains[i][1]);
543 else if (matchedLegs.size() == 1) {
546 int iJun = junChains[i][0];
547 int iAnti = junChains[i][1];
548 int iLeg = matchedLegs[0].first;
549 int iAntiLeg = matchedLegs[0].second;
550 if (event.kindJunction(iAnti) % 2 == 1) {
552 swap(iLeg, iAntiLeg);
556 int iJunList = -1, iAntiList = -1;
557 for (
int l = 0;l < int(iPartonJun.size()); ++l)
558 if (- iPartonJun[l][0]/10 - 1 == iJun) {
562 for (
int l = 0;l < int(iPartonAntiJun.size()); ++l)
563 if (- iPartonAntiJun[l][0]/10 - 1 == iAnti) {
567 if (iJunList == -1 || iAntiList == -1) {
568 infoPtr->errorMsg(
"Error in JunctionSplitting::SplitJunChain:"
569 " failed to find junctions in the parton list");
574 vector<vector <int> > iJunLegs;
578 for (
int l = 0; l < int(iPartonJun[iJunList].size()); ++l) {
579 if ( iPartonJun[iJunList][l]/10 == iPartonJun[iJunList][0]/10) ++leg;
580 iJunLegs[leg].push_back(iPartonJun[iJunList][l]);
584 vector<vector <int> > iAntiLegs;
587 for (
int l = 0; l < int(iPartonAntiJun[iAntiList].size()); ++l) {
588 if ( iPartonAntiJun[iAntiList][l]/10
589 == iPartonAntiJun[iAntiList][0]/10) ++leg;
590 iAntiLegs[leg].push_back(iPartonAntiJun[iAntiList][l]);
594 vector<int>& iJunLeg0 = (iLeg == 0) ? iJunLegs[1] : iJunLegs[0];
595 vector<int>& iJunLeg1 = (iLeg == 2) ? iJunLegs[1] : iJunLegs[2];
596 vector<int>& iAntiLeg0 = (iAntiLeg == 0) ? iAntiLegs[1] : iAntiLegs[0];
597 vector<int>& iAntiLeg1 = (iAntiLeg == 2) ? iAntiLegs[1] : iAntiLegs[2];
601 if (iAntiLeg0[1] < 0 || iAntiLeg1[1] < 0)
continue;
604 Vec4 pJunLeg0 =
event[ iJunLeg0[1] ].p();
605 Vec4 pJunLeg1 =
event[ iJunLeg1[1] ].p();
606 Vec4 pAntiLeg0 =
event[ iAntiLeg0[1] ].p();
607 Vec4 pAntiLeg1 =
event[ iAntiLeg1[1] ].p();
610 if ( theta(pJunLeg0,pJunLeg1) < MINANGLE
611 || theta(pAntiLeg0,pAntiLeg1) < MINANGLE
612 || theta(pJunLeg0,pAntiLeg0) < MINANGLE
613 || theta(pJunLeg0,pAntiLeg1) < MINANGLE
614 || theta(pJunLeg1,pAntiLeg0) < MINANGLE
615 || theta(pJunLeg1,pAntiLeg1) < MINANGLE) {
616 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
617 "parallel junction state not allowed.");
622 Vec4 pStart = pJunLeg0 / pJunLeg0.e() + pJunLeg1 / pJunLeg1.e()
623 + pAntiLeg0 / pAntiLeg0.e() + pAntiLeg1 / pAntiLeg1.e();
626 RotBstMatrix MtoJRF, MtoARF;
627 Vec4 pInJRF[3], pInARF[3];
628 for (
int iJunLocal = 0; iJunLocal < 2; ++iJunLocal) {
629 int offset = (iJunLocal == 0) ? 0 : 2;
632 RotBstMatrix MtoRF, Mstep;
633 MtoRF.bstback(pStart);
641 pInRF[0 + offset] = pJunLeg0;
642 pInRF[1 + offset] = pJunLeg1;
643 pInRF[2 - offset] = pAntiLeg0;
644 pInRF[3 - offset] = pAntiLeg1;
645 for (
int l = 0; l < 4; ++l) pInRF[l].rotbst(MtoRF);
648 double wt2 = 1. - exp( -pInRF[2].e() / eNormJunction);
649 double wt3 = 1. - exp( -pInRF[3].e() / eNormJunction);
650 pInRF[2] = wt2 * pInRF[2] + wt3 * pInRF[3];
653 Mstep = stringFrag.junctionRestFrame( pInRF[0], pInRF[1], pInRF[2]);
654 MtoRF.rotbst( Mstep );
655 }
while (iter < 3 || (Mstep.deviation() > CONVJNREST
656 && iter < NTRYJNREST) );
659 if (iJunLocal == 0) {
661 for (
int l = 0; l < 3; ++l) pInJRF[l] = pInRF[l];
664 for (
int l = 0; l < 3; ++l) pInARF[l] = pInRF[l];
669 RotBstMatrix MfromJRF = MtoJRF;
671 RotBstMatrix MfromARF = MtoARF;
675 Vec4 vJun(0., 0., 0., 1.);
676 vJun.rotbst(MfromJRF);
677 Vec4 vAnti(0., 0., 0., 1.);
678 vAnti.rotbst(MfromARF);
679 Vec4 pLabJ[3], pLabA[3];
680 for (
int l = 0; l < 3; ++l) {
681 pLabJ[l] = pInJRF[l];
682 pLabJ[l].rotbst(MfromJRF);
683 pLabA[l] = pInARF[l];
684 pLabA[l].rotbst(MfromARF);
688 double vJvA = vJun * vAnti;
689 double vJvAe2y = vJvA + sqrt(vJvA*vJvA - 1.);
690 double lambdaJA = stringLength.getJuncLength(pInJRF[0], pInJRF[1],
691 pInARF[0], pInARF[1]);
693 double lambda00 = stringLength.getStringLength(pLabJ[0], pLabA[0]) +
694 stringLength.getStringLength(pLabJ[1], pLabA[1]);
696 double lambda01 = stringLength.getStringLength(pLabJ[0], pLabA[1]) +
697 stringLength.getStringLength(pLabJ[1], pLabA[0]);
700 if (lambdaJA > min( lambda00, lambda01) && allowDoubleJunRem) {
703 int iCol1 = iJunLeg0[1];
704 int iCol2 = iJunLeg1[1];
705 int iAcol1 = iAntiLeg0[1];
706 int iAcol2 = iAntiLeg1[1];
707 if (lambda00 > lambda01)
708 swap(iAcol1, iAcol2);
711 int iNew1 =
event.copy(iAcol1, 66);
712 event[iNew1].acol(event[iCol1].col());
713 int iNew2 =
event.copy(iAcol2, 66);
714 event[iNew2].acol(event[iCol2].col());
715 junRem.push_back(junChains[i][0]);
716 junRem.push_back(junChains[i][1]);
727 int idJ0Abs =
event[ iJunLeg0[1] ].idAbs();
728 bool acceptJ0 = idJ0Abs < 6 || idJ0Abs == 21 ||
729 (idJ0Abs > 1000 && idJ0Abs < 6000 && (idJ0Abs / 10) % 10 == 0);
730 int idJ1Abs =
event[ iJunLeg1[1] ].idAbs();
731 bool acceptJ1 = idJ1Abs < 6 || idJ1Abs == 21 ||
732 (idJ1Abs > 1000 && idJ1Abs < 6000 && (idJ1Abs / 10) % 10 == 0);
733 int idA0Abs =
event[ iAntiLeg0[1] ].idAbs();
734 bool acceptA0 = idA0Abs < 6 || idA0Abs == 21 ||
735 (idA0Abs > 1000 && idA0Abs < 6000 && (idA0Abs / 10) % 10 == 0);
736 int idA1Abs =
event[ iAntiLeg1[1] ].idAbs();
737 bool acceptA1 = idA1Abs < 6 || idA1Abs == 21 ||
738 (idA1Abs > 1000 && idA1Abs < 6000 && (idA1Abs / 10) % 10 == 0);
740 if ( !(acceptJ0 || acceptJ1)) {
741 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
742 "No light quarks available in junction split.");
746 if ( !(acceptA0 || acceptA1)) {
747 infoPtr->errorMsg(
"Warning in JunctionSplitting::SplitJunPairs: "
748 "No light quarks available in junction split.");
752 double eShift = MTHAD / (3. * sqrt(vJvAe2y));
753 double fracJ0 = 0, fracJ1 = 0, fracA0 = 0, fracA1 = 0;
755 fracJ0 = min(0.5, eShift / pInJRF[0].e());
757 fracJ1 = min(0.5, eShift / pInJRF[1].e());
758 Vec4 pFromJun = fracJ0 * pJunLeg0 + fracJ1 * pJunLeg1;
760 fracA0 = min(0.5, eShift / pInARF[0].e());
762 fracA1 = min(0.5, eShift / pInARF[1].e());
763 Vec4 pFromAnti = fracA0 * pAntiLeg0 + fracA1 * pAntiLeg1;
766 int idQ = flavSel.pickLightQ();
768 int junMother1 = min(iJunLeg0[1], iJunLeg1[1]);
769 int junMother2 = max(iJunLeg0[1], iJunLeg1[1]);
770 int antiMother1 = min(iAntiLeg0[1], iAntiLeg1[1]);
771 int antiMother2 = max(iAntiLeg0[1], iAntiLeg1[1]);
774 int iJunNew1 =
event.copy(iJunLeg0[1], 76);
775 event[iJunNew1].rescale5(1. - fracJ0);
776 iJunLeg0[1] = iJunNew1;
777 event[iJunNew1].mothers(junMother2, junMother1);
779 int iJunNew2 =
event.copy(iJunLeg1[1], 76);
780 event[iJunNew2].rescale5(1. - fracJ1);
781 iJunLeg1[1] = iJunNew2;
782 event[iJunNew2].mothers(junMother2, junMother1);
786 int acolQ =
event.nextColTag();
787 event.endColJunction(iAnti, iAntiLeg, acolQ);
788 event.colJunction(iAnti, iAntiLeg, acolQ);
789 int iNewA =
event.append( -idQ, 76, junMother2, junMother1, 0, 0,
790 0, acolQ, pFromJun, pFromJun.mCalc() );
793 int iAntiNew1 =
event.copy(iAntiLeg0[1], 76);
794 event[iAntiNew1].rescale5(1. - fracA0);
795 iAntiLeg0[1] = iAntiNew1;
796 event[iAntiNew1].mothers(antiMother2, antiMother1);
798 int iAntiNew2 =
event.copy(iAntiLeg1[1], 76);
799 event[iAntiNew2].rescale5(1. - fracA1);
800 iAntiLeg1[1] = iAntiNew2;
801 event[iAntiNew2].mothers(antiMother2, antiMother1);
804 int colQ =
event.nextColTag();
805 event.endColJunction(iJun, iLeg, colQ);
806 event.colJunction(iJun, iLeg, colQ);
807 int iNewJ =
event.append( idQ, 76, antiMother2, antiMother1, 0, 0,
808 colQ, 0, pFromAnti, pFromAnti.mCalc() );
811 event[
event[iJunNew1].mother1()].daughters( iJunNew1, iNewA);
812 event[
event[iJunNew1].mother2()].daughters( iJunNew1, iNewA);
813 event[
event[iAntiNew1].mother1()].daughters( iAntiNew1, iNewJ);
814 event[
event[iAntiNew1].mother2()].daughters( iAntiNew1, iNewJ);
822 sort(junRem.begin(),junRem.end());
823 reverse(junRem.begin(),junRem.end());
824 for (
int i = 0; i < int(junRem.size()); ++i)
825 event.eraseJunction(junRem[i]);
826 event.saveJunctionSize();
838 bool JunctionSplitting::getPartonLists(
Event& event,
839 vector<vector< int > > & iPartonJun, vector<vector<int > >& iPartonAntiJun) {
842 colTrace.setupColList(event);
845 iPartonAntiJun.clear();
851 for (
int iLoop = 0; iLoop < 2*
event.sizeJunction(); ++iLoop) {
853 int iJun = iLoop %
event.sizeJunction();
854 if ( !event.remainsJunction(iJun))
continue;
857 int kindJun =
event.kindJunction(iJun);
858 if ( iLoop < event.sizeJunction() && kindJun%2 == 0)
continue;
859 else if ( iLoop >= event.sizeJunction() && kindJun%2 == 1)
continue;
867 for (
int iCol = 0; iCol < 3; ++iCol) {
868 int indxCol =
event.colJunction(iJun, iCol);
869 iParton.push_back( -(10 + 10 * iJun + iCol) );
871 if ( kindJun%2 == 1 && !colTrace.traceFromAcol(indxCol,event, iJun,
872 iCol, iParton) )
return false;
874 else if ( kindJun%2 == 0 && !colTrace.traceFromCol(indxCol,event, iJun,
875 iCol, iParton) )
return false;
880 for (
int i = 0; i < int(iParton.size()); ++i)
if (iParton[i] < 0) ++nNeg;
883 if ( kindJun%2 == 1 ) iPartonJun.push_back(iParton);
884 else iPartonAntiJun.push_back(iParton);
898 bool JunctionSplitting::setAcol(
Event& event,
int col,
int acol) {
901 for (
int j = 0;j <
event.size(); ++j)
902 if (event[j].isFinal() &&
event[j].acol() == acol) {
903 int iNew =
event.copy(j,66);
904 event[iNew].acol(col);
908 for (
int j = 0;j <
event.sizeJunction(); ++j)
909 for (
int jLeg = 0;jLeg < 3; ++jLeg)
910 if (event.colJunction(j, jLeg) == acol) {
911 event.colJunction(j, jLeg, col);
916 infoPtr->errorMsg(
"Warning in JunctionSplitting::setAcol:"
917 "Anti colour not found when combing two junctions to a string");