StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
BeamRemnants.cc
1 // BeamRemnants.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // BeamRemnants class.
8 
9 #include "Pythia8/BeamRemnants.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The BeamDipole class is purely internal to reconnectColours.
16 
17 class BeamDipole {
18 
19 public:
20 
21  // Constructor.
22  BeamDipole( int colIn = 0, int iColIn = 0, int iAcolIn = 0)
23  : col(colIn), iCol(iColIn), iAcol(iAcolIn) {}
24 
25  // Members.
26  int col, iCol, iAcol;
27  double p1p2;
28 
29 };
30 
31 //==========================================================================
32 
33 // The BeamRemnants class.
34 
35 //--------------------------------------------------------------------------
36 
37 // Constants: could be changed here if desired, but normally should not.
38 // These are of technical nature, as described for each.
39 
40 // If same (anti)colour appears twice in final state, repair or reject.
41 const bool BeamRemnants::ALLOWCOLOURTWICE = true;
42 
43 // Maximum number of tries to match colours and kinematics in the event.
44 const int BeamRemnants::NTRYCOLMATCH = 10;
45 const int BeamRemnants::NTRYKINMATCH = 10;
46 
47 // Overall correction step for energy-momentum conservation; only
48 // becomes relevant in rescattering scenarios when FSR dipole emissions
49 // and primordial kT is added. Should hopefully not be needed.
50 const bool BeamRemnants::CORRECTMISMATCH = false;
51 
52 //--------------------------------------------------------------------------
53 
54 // Initialization.
55 
56 bool BeamRemnants::init( Info* infoPtrIn, Settings& settings, Rndm* rndmPtrIn,
57  BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
58  PartonSystems* partonSystemsPtrIn) {
59 
60  // Save pointers.
61  infoPtr = infoPtrIn;
62  rndmPtr = rndmPtrIn;
63  beamAPtr = beamAPtrIn;
64  beamBPtr = beamBPtrIn;
65  partonSystemsPtr = partonSystemsPtrIn;
66 
67  // Width of primordial kT distribution.
68  doPrimordialKT = settings.flag("BeamRemnants:primordialKT");
69  primordialKTsoft = settings.parm("BeamRemnants:primordialKTsoft");
70  primordialKThard = settings.parm("BeamRemnants:primordialKThard");
71  primordialKTremnant = settings.parm("BeamRemnants:primordialKTremnant");
72  halfScaleForKT = settings.parm("BeamRemnants:halfScaleForKT");
73  halfMassForKT = settings.parm("BeamRemnants:halfMassForKT");
74 
75  // Handling of rescattering kinematics uncertainties from primodial kT.
76  allowRescatter = settings.flag("MultipartonInteractions:allowRescatter");
77  doRescatterRestoreY = settings.flag("BeamRemnants:rescatterRestoreY");
78 
79  // Parameters for colour reconnection scenario, partly borrowed from
80  // multiparton interactions not to introduce too many new ones.
81  doReconnect = settings.flag("BeamRemnants:reconnectColours");
82  reconnectRange = settings.parm("BeamRemnants:reconnectRange");
83  pT0Ref = settings.parm("MultipartonInteractions:pT0Ref");
84  ecmRef = settings.parm("MultipartonInteractions:ecmRef");
85  ecmPow = settings.parm("MultipartonInteractions:ecmPow");
86 
87  // Total and squared CM energy at nominal energy.
88  eCM = infoPtr->eCM();
89  sCM = eCM * eCM;
90 
91  // The MPI pT0 smoothening scale and its reconnection-strength combination.
92  pT0 = pT0Ref * pow(eCM / ecmRef, ecmPow);
93  pT20Rec = pow2(reconnectRange * pT0);
94 
95  // Done.
96  return true;
97 
98 }
99 
100 //--------------------------------------------------------------------------
101 
102 // Select the flavours/kinematics/colours of the two beam remnants.
103 // Notation: iPar = all partons, iSys = matched systems of two beams,
104 // iRem = additional partons in remnants.
105 
106 bool BeamRemnants::add( Event& event) {
107 
108  // Update to current CM energy.
109  eCM = infoPtr->eCM();
110  sCM = eCM * eCM;
111 
112  // Check that flavour bookkept in event and in beam particles agree.
113  for (int i = 0; i < beamAPtr->size(); ++i) {
114  int j = (*beamAPtr)[i].iPos();
115  if ((*beamAPtr)[i].id() != event[j].id()) {
116  infoPtr->errorMsg("Error in BeamRemnants::add: "
117  "event and beam flavours do not match");
118  return false;
119  }
120  }
121  for (int i = 0; i < beamBPtr->size(); ++i) {
122  int j = (*beamBPtr)[i].iPos();
123  if ((*beamBPtr)[i].id() != event[j].id()) {
124  infoPtr->errorMsg("Error in BeamRemnants::add: "
125  "event and beam flavours do not match");
126  return false;
127  }
128  }
129 
130  // Number of scattering subsystems. Size of event record before treatment.
131  nSys = partonSystemsPtr->sizeSys();
132  oldSize = event.size();
133 
134  // Add required extra remnant flavour content.
135  // Start all over if fails (in option where junctions not allowed).
136  if ( !beamAPtr->remnantFlavours(event)
137  || !beamBPtr->remnantFlavours(event) ) {
138  infoPtr->errorMsg("Error in BeamRemnants::add:"
139  " remnant flavour setup failed");
140  return false;
141  }
142 
143  // Do the kinematics of the collision subsystems and two beam remnants.
144  if (!setKinematics(event)) return false;
145 
146  // Allow colour reconnections.
147  if (doReconnect) reconnectColours(event);
148 
149  // Save current modifiable colour configuration for fast restoration.
150  vector<int> colSave;
151  vector<int> acolSave;
152  for (int i = oldSize; i < event.size(); ++i) {
153  colSave.push_back( event[i].col() );
154  acolSave.push_back( event[i].acol() );
155  }
156  event.saveJunctionSize();
157 
158  // Allow several tries to match colours of initiators and remnants.
159  // Frequent "failures" since shortcutting colours separately on
160  // the two event sides may give "colour singlet gluons" etc.
161  bool physical = true;
162  for (int iTry = 0; iTry < NTRYCOLMATCH; ++iTry) {
163  physical = true;
164 
165  // Reset list of colour "collapses" (transformations).
166  colFrom.resize(0);
167  colTo.resize(0);
168 
169  // First process each set of beam colours on its own.
170  if (!beamAPtr->remnantColours(event, colFrom, colTo))
171  physical = false;
172  if (!beamBPtr->remnantColours(event, colFrom, colTo))
173  physical = false;
174 
175  // Then check that colours and anticolours are matched in whole event.
176  if ( physical && !checkColours(event) ) physical = false;
177 
178  // If no problems then done, else restore and loop.
179  if (physical) break;
180  for (int i = oldSize; i < event.size(); ++i)
181  event[i].cols( colSave[i - oldSize], acolSave[i - oldSize] );
182  event.restoreJunctionSize();
183  infoPtr->errorMsg("Warning in BeamRemnants::add:"
184  " colour tracing failed; will try again");
185  }
186 
187  // If no solution after several tries then failed.
188  if (!physical) {
189  infoPtr->errorMsg("Error in BeamRemnants::add:"
190  " colour tracing failed after several attempts");
191  return false;
192  }
193 
194  // Done.
195  return true;
196 
197 }
198 
199 //--------------------------------------------------------------------------
200 
201 // Set up trial transverse and longitudinal kinematics for each beam
202 // separately. Final decisions involve comparing the two beams.
203 
204 bool BeamRemnants::setKinematics( Event& event) {
205 
206  // References to beams to simplify indexing.
207  BeamParticle& beamA = *beamAPtr;
208  BeamParticle& beamB = *beamBPtr;
209 
210  // Nothing to do for lepton-lepton scattering with all energy already used.
211  if ( beamA.isUnresolvedLepton() && beamB.isUnresolvedLepton() )
212  return true;
213 
214  // Check that has not already used up beams.
215  if ( (!beamA.isLepton() && beamA.xMax(-1) <= 0.)
216  || (!beamB.isLepton() && beamB.xMax(-1) <= 0.) ) {
217  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
218  " no momentum left for beam remnants");
219  return false;
220  }
221 
222  // Last beam-status particles. Offset relative to normal beam locations.
223  int nBeams = 3;
224  for (int i = 3; i < event.size(); ++i)
225  if (event[i].statusAbs() < 20) nBeams = i + 1;
226  int nOffset = nBeams - 3;
227 
228  // Reserve space for extra information on the systems and beams.
229  int nMaxBeam = max( beamA.size(), beamB.size() );
230  vector<double> sHatSys(nMaxBeam);
231  vector<double> kTwidth(nMaxBeam);
232  vector<double> kTcomp(nMaxBeam);
233  vector<RotBstMatrix> Msys(nSys);
234 
235  // Loop over all subsystems. Default values. Find invariant mass.
236  double kTcompSumA = 0.;
237  double kTcompSumB = 0.;
238  for (int iSys = 0; iSys < nSys; ++iSys) {
239  double kTwidthNow = 0.;
240  double mHatDamp = 1.;
241  int iInA = partonSystemsPtr->getInA(iSys);
242  int iInB = partonSystemsPtr->getInB(iSys);
243  double sHatNow = (event[iInA].p() + event[iInB].p()).m2Calc();
244 
245  // Allow primordial kT reduction for small-mass and small-pT systems
246  // (for hardest interaction pT -> renormalization scale so also 2 -> 1).
247  if (doPrimordialKT) {
248  double mHat = sqrt(sHatNow);
249  mHatDamp = mHat / (mHat + halfMassForKT);
250  double scale = (iSys == 0) ? infoPtr->QRen(iDS)
251  : partonSystemsPtr->getPTHat(iSys);
252  kTwidthNow = ( (halfScaleForKT * primordialKTsoft
253  + scale * primordialKThard) / (halfScaleForKT + scale) ) * mHatDamp;
254  }
255 
256  // Store properties of compensation systems and total compensation power.
257  // Rescattered partons do not compensate, but may be massive.
258  sHatSys[iSys] = sHatNow;
259  kTwidth[iSys] = kTwidthNow ;
260  kTcomp[iSys] = mHatDamp;
261  if (beamA[iSys].isFromBeam()) kTcompSumA += mHatDamp;
262  else beamA[iSys].m( event[iInA].m() );
263  if (beamB[iSys].isFromBeam()) kTcompSumB += mHatDamp;
264  else beamB[iSys].m( event[iInB].m() );
265  }
266 
267  // Primordial kT and compensation power among remnants.
268  double kTwidthNow = (doPrimordialKT) ? primordialKTremnant : 0.;
269  for (int iRem = nSys; iRem < nMaxBeam; ++iRem) {
270  sHatSys[iRem] = 0.;
271  kTwidth[iRem] = kTwidthNow ;
272  kTcomp[iRem] = 1.;
273  }
274  kTcompSumA += beamA.size() - nSys;
275  kTcompSumB += beamB.size() - nSys;
276 
277  // Allow ten tries to construct kinematics (but normally works first).
278  bool physical;
279  double xSum[2], xInvM[2], w2Beam[2], wPosRem, wNegRem, w2Rem;
280  for (int iTry = 0; iTry < NTRYKINMATCH; ++iTry) {
281  physical = true;
282 
283  // Loop over the two beams.
284  for (int iBeam = 0; iBeam < 2; ++iBeam) {
285  BeamParticle& beam = (iBeam == 0) ? beamA : beamB;
286  int nPar = beam.size();
287 
288  // Generate Gaussian pT for initiator partons inside hadrons.
289  // Store/restore rescattered parton momenta before primordial kT.
290  if (beam.isHadron() && doPrimordialKT) {
291  double pxSum = 0.;
292  double pySum = 0.;
293  for (int iPar = 0; iPar < nPar; ++iPar) {
294  if ( beam[iPar].isFromBeam() ) {
295  pair<double, double> gauss2 = rndmPtr->gauss2();
296  double px = kTwidth[iPar] * gauss2.first;
297  double py = kTwidth[iPar] * gauss2.second;
298  beam[iPar].px(px);
299  beam[iPar].py(py);
300  pxSum += px;
301  pySum += py;
302  } else {
303  int iInAB = (iBeam == 0) ? partonSystemsPtr->getInA(iPar)
304  : partonSystemsPtr->getInB(iPar);
305  beam[iPar].p( event[iInAB].p() );
306  }
307  }
308 
309  // Share recoil between all initiator partons, rescatterers excluded.
310  double kTcompSum = (iBeam == 0) ? kTcompSumA : kTcompSumB;
311  for (int iPar = 0; iPar < nPar; ++iPar)
312  if (beam[iPar].isFromBeam() ) {
313  beam[iPar].px( beam[iPar].px() - pxSum * kTcomp[iPar] / kTcompSum );
314  beam[iPar].py( beam[iPar].py() - pySum * kTcomp[iPar] / kTcompSum );
315  }
316 
317  // Without primordial kT: still need to store rescattered partons.
318  } else if (beam.isHadron()) {
319  for (int iPar = 0; iPar < nPar; ++iPar)
320  if ( !beam[iPar].isFromBeam() ) {
321  int iInAB = (iBeam == 0) ? partonSystemsPtr->getInA(iPar)
322  : partonSystemsPtr->getInB(iPar);
323  beam[iPar].p( event[iInAB].p() );
324  }
325  }
326 
327  // Pick unrescaled x values for remnants. Sum up (unscaled) p+ and p-.
328  xSum[iBeam] = 0.;
329  xInvM[iBeam] = 0.;
330  for (int iRem = nSys; iRem < nPar; ++iRem) {
331  double xPrel = beam.xRemnant( iRem);
332  beam[iRem].x(xPrel);
333  xSum[iBeam] += xPrel;
334  xInvM[iBeam] += beam[iRem].mT2()/xPrel;
335  }
336 
337  // Squared transverse mass for each beam, using lightcone x.
338  w2Beam[iBeam] = xSum[iBeam] * xInvM[iBeam];
339 
340  // End separate treatment of the two beams.
341  }
342 
343  // Recalculate kinematics of initiator systems with primordial kT.
344  wPosRem = eCM;
345  wNegRem = eCM;
346  for (int iSys = 0; iSys < nSys; ++iSys) {
347  int iA = beamA[iSys].iPos();
348  int iB = beamB[iSys].iPos();
349  double sHat = sHatSys[iSys];
350 
351  // Classify system: rescattering on either or both sides?
352  bool normalA = beamA[iSys].isFromBeam();
353  bool normalB = beamB[iSys].isFromBeam();
354  bool normalSys = normalA && normalB;
355  bool doubleRes = !normalA && !normalB;
356 
357  // Check whether final-state system momentum matches initial-state one.
358  if (allowRescatter && CORRECTMISMATCH) {
359  Vec4 pInitial = event[iA].p() + event[iB].p();
360  Vec4 pFinal;
361  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
362  int iAB = partonSystemsPtr->getOut(iSys, iMem);
363  if (event[iAB].isFinal()) pFinal += event[iAB].p();
364  }
365 
366  // Scale down primordial kT from side A if p+ increased.
367  if (normalA && pFinal.pPos() > pInitial.pPos())
368  beamA[iSys].scalePT( pInitial.pPos() / pFinal.pPos() );
369 
370  // Scale down primordial kT from side B if p- increased.
371  if (normalB && pFinal.pNeg() > pInitial.pNeg())
372  beamB[iSys].scalePT( pInitial.pNeg() / pFinal.pNeg() );
373  }
374 
375  // Rescatter: possible change in sign of lightcone momentum of a
376  // rescattered parton. If this happens, try to pick
377  // new primordial kT values
378  if (allowRescatter
379  && (event[iA].pPos() / beamA[iSys].pPos() < 0
380  || event[iB].pNeg() / beamB[iSys].pNeg() < 0) ) {
381  infoPtr->errorMsg("Warning in BeamRemnants::setKinematics:"
382  " change in lightcone momentum sign; retrying kinematics");
383  physical = false;
384  break;
385  }
386 
387  // Begin kinematics of partons after primordial kT has been added.
388  double sHatTAft = sHat + pow2( beamA[iSys].px() + beamB[iSys].px())
389  + pow2( beamA[iSys].py() + beamB[iSys].py());
390  double w2A = beamA[iSys].mT2();
391  double w2B = beamB[iSys].mT2();
392  double w2Diff = sHatTAft - w2A - w2B;
393  double lambda = pow2(w2Diff) - 4. * w2A * w2B;
394 
395  // Too large transverse momenta means that kinematics will not work.
396  if (lambda <= 0.) {
397  physical = false;
398  break;
399  }
400  double lamRoot = sqrtpos( lambda );
401 
402  // Mirror solution if the two incoming have reverse rapidity ordering.
403  if (allowRescatter && doubleRes && (event[iA].pPos() * event[iB].pNeg()
404  < event[iA].pNeg() * event[iB].pPos()) ) lamRoot = -lamRoot;
405 
406  // Two procedures, which agree for normal scattering, separate here.
407  // First option keeps rapidity (and mass) of system unchanged by
408  // primordial kT, by modification of rescattered parton.
409  if (normalSys || doRescatterRestoreY || doubleRes) {
410 
411  // Find kinematics of initial system: transverse mass, and
412  // longitudinal momentum carried by all or rescattered partons.
413  double sHatTBef = sHat;
414  double wPosBef, wNegBef, wPosBefRes, wNegBefRes;
415  // Normal scattering.
416  if (normalSys) {
417  wPosBef = beamA[iSys].x() * eCM;
418  wNegBef = beamB[iSys].x() * eCM;
419  wPosBefRes = 0.;
420  wNegBefRes = 0.;
421  // Rescattering on side A.
422  } else if (normalB) {
423  sHatTBef += event[iA].pT2();
424  wPosBef = event[iA].pPos();
425  wNegBef = event[iA].pNeg() + beamB[iSys].x() * eCM;
426  wPosBefRes = beamA[iSys].pPos();
427  wNegBefRes = beamA[iSys].pNeg();
428  // Rescattering on side B.
429  } else if (normalA) {
430  sHatTBef += event[iB].pT2();
431  wPosBef = beamA[iSys].x() * eCM + event[iB].pPos();
432  wNegBef = event[iB].pNeg();
433  wPosBefRes = beamB[iSys].pPos();
434  wNegBefRes = beamB[iSys].pNeg();
435  // Rescattering on both sides.
436  } else {
437  sHatTBef += pow2( event[iA].px() + event[iB].px())
438  + pow2( event[iA].py() + event[iB].py());
439  wPosBef = event[iA].pPos() + event[iB].pPos();
440  wNegBef = event[iA].pNeg() + event[iB].pNeg();
441  wPosBefRes = beamA[iSys].pPos() + beamB[iSys].pPos();
442  wNegBefRes = beamA[iSys].pNeg() + beamB[iSys].pNeg();
443  }
444 
445  // Rescale outgoing momenta to keep same mass and rapidity of system
446  // as originally, and solve for kinematics.
447  double rescale = sqrt(sHatTAft / sHatTBef);
448  double wPosAft = rescale * wPosBef;
449  double wNegAft = rescale * wNegBef;
450  wPosRem -= wPosAft - wPosBefRes;
451  wNegRem -= wNegAft - wNegBefRes;
452  double wPosA = 0.5 * (sHatTAft + w2A - w2B + lamRoot) / wNegAft;
453  double wNegB = 0.5 * (sHatTAft + w2B - w2A + lamRoot) / wPosAft;
454 
455  // Store modified beam parton momenta.
456  beamA[iSys].e( 0.5 * (wPosA + w2A / wPosA) );
457  beamA[iSys].pz( 0.5 * (wPosA - w2A / wPosA) );
458  beamB[iSys].e( 0.5 * (w2B / wNegB + wNegB) );
459  beamB[iSys].pz( 0.5 * (w2B / wNegB - wNegB) );
460 
461  // Second option keeps rescattered parton (and system mass) unchanged
462  // by primordial kT, by modification of system rapidity.
463  } else {
464 
465  // Rescattering on side A: preserve already scattered parton.
466  if (normalB) {
467  double wPosA = beamA[iSys].pPos();
468  double wNegB = 0.5 * (w2Diff + lamRoot) / wPosA;
469  beamB[iSys].e( 0.5 * (w2B / wNegB + wNegB) );
470  beamB[iSys].pz( 0.5 * (w2B / wNegB - wNegB) );
471  wPosRem -= w2B / wNegB;
472  wNegRem -= wNegB;
473 
474 
475  // Rescattering on side B: preserve already scattered parton.
476  } else if (normalA) {
477  double wNegB = beamB[iSys].pNeg();
478  double wPosA = 0.5 * (w2Diff + lamRoot) / wNegB;
479  beamA[iSys].e( 0.5 * (wPosA + w2A / wPosA) );
480  beamA[iSys].pz( 0.5 * (wPosA - w2A / wPosA) );
481  wPosRem -= wPosA;
482  wNegRem -= w2A / wPosA;
483 
484  // Primordial kT in double rescattering does change the mass of
485  // the system without any possible fix in this framework, which
486  // leads to incorrect boosts. Current choice is therefore always
487  // to handle it with the first procedure, where mass is retained.
488  } else {
489  }
490  }
491 
492  // Construct system rotation and boost caused by primordial kT.
493  Msys[iSys].reset();
494  Msys[iSys].toCMframe( event[iA].p(), event[iB].p() );
495  Msys[iSys].fromCMframe( beamA[iSys].p(), beamB[iSys].p() );
496 
497  // Boost rescattered partons in subsequent beam A list.
498  for (int iSys2 = iSys + 1; iSys2 < nSys; ++iSys2) {
499  if (!beamA[iSys2].isFromBeam()) {
500  int iBefResc = event[ beamA[iSys2].iPos() ].mother1();
501  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem)
502  if (partonSystemsPtr->getOut(iSys, iMem) == iBefResc) {
503  Vec4 pTemp = event[iBefResc].p();
504  pTemp.rotbst( Msys[iSys] );
505  beamA[iSys2].p( pTemp );
506  }
507  }
508 
509  // Boost rescattered partons in subsequent beam B list.
510  if (!beamB[iSys2].isFromBeam()) {
511  int iBefResc = event[ beamB[iSys2].iPos() ].mother1();
512  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem)
513  if (partonSystemsPtr->getOut(iSys, iMem) == iBefResc) {
514  Vec4 pTemp = event[iBefResc].p();
515  pTemp.rotbst( Msys[iSys] );
516  beamB[iSys2].p( pTemp );
517  }
518  }
519  }
520  }
521 
522  // Check that remaining momentum is enough for remnants.
523  if (wPosRem < 0. || wNegRem < 0.) physical = false;
524  w2Rem = wPosRem * wNegRem;
525  if (sqrtpos(w2Rem) < sqrt(w2Beam[0]) + sqrt(w2Beam[1]))
526  physical = false;
527 
528  // End of loop over ten tries. Do not loop when solution found.
529  if (physical) break;
530  }
531 
532  // If no solution after ten tries then failed.
533  if (!physical) {
534  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
535  " kinematics construction failed");
536  return false;
537  }
538 
539  // For successful initiator kinematics process whole systems.
540  Vec4 pSumOut;
541  for (int iSys = 0; iSys < nSys; ++iSys) {
542 
543  // Copy initiators and their systems and boost them accordingly.
544  // Update subsystem and beams info on new positions of partons.
545  // Update daughter info of mothers, i.e. of beams, for hardest interaction.
546  if (beamA[iSys].isFromBeam()) {
547  int iA = beamA[iSys].iPos();
548  int iAcopy = event.copy(iA, -61);
549  event[iAcopy].rotbst(Msys[iSys]);
550  partonSystemsPtr->setInA(iSys, iAcopy);
551  beamA[iSys].iPos( iAcopy);
552  if (iSys == 0) {
553  int mother = event[iAcopy].mother1();
554  event[mother].daughter1(iAcopy);
555  }
556  }
557  if (beamB[iSys].isFromBeam()) {
558  int iB = beamB[iSys].iPos();
559  int iBcopy = event.copy(iB, -61);
560  event[iBcopy].rotbst(Msys[iSys]);
561  partonSystemsPtr->setInB(iSys, iBcopy);
562  beamB[iSys].iPos( iBcopy);
563  if (iSys == 0) {
564  int mother = event[iBcopy].mother1();
565  event[mother].daughter1(iBcopy);
566  }
567  }
568 
569  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
570  int iAB = partonSystemsPtr->getOut(iSys, iMem);
571  if (event[iAB].isFinal()) {
572  int iABcopy = event.copy(iAB, 62);
573  event[iABcopy].rotbst(Msys[iSys]);
574  partonSystemsPtr->setOut(iSys, iMem, iABcopy);
575  pSumOut += event[iABcopy].p();
576  }
577  }
578 
579  }
580 
581  // Colour dipoles spanning systems gives mismatch between FSR recoils
582  // and primordial kT boosts.
583  if (allowRescatter && CORRECTMISMATCH) {
584 
585  // Find summed pT of beam remnants = - wanted pT of systems.
586  double pxBeams = 0.;
587  double pyBeams = 0.;
588  for (int iRem = nSys; iRem < beamA.size(); ++iRem) {
589  pxBeams += beamA[iRem].px();
590  pyBeams += beamA[iRem].py();
591  }
592  for (int iRem = nSys; iRem < beamB.size(); ++iRem) {
593  pxBeams += beamB[iRem].px();
594  pyBeams += beamB[iRem].py();
595  }
596 
597  // Boost all final partons in systems transversely, and also their sum.
598  Vec4 pSumTo( -pxBeams, -pyBeams, pSumOut.pz(), sqrt( pow2(pxBeams)
599  + pow2(pyBeams) + pow2(pSumOut.pz()) + pSumOut.m2Calc() ) );
600  RotBstMatrix Mmismatch;
601  Mmismatch.bst( pSumOut, pSumTo);
602  for (int iSys = 0; iSys < nSys; ++iSys)
603  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
604  int iAB = partonSystemsPtr->getOut(iSys, iMem);
605  if (event[iAB].isFinal()) event[iAB].rotbst(Mmismatch);
606  }
607  pSumOut.rotbst(Mmismatch);
608 
609  // Reset energy and momentum sum, to be compensated by beam remnants.
610  wPosRem = eCM - (pSumOut.e() + pSumOut.pz());
611  wNegRem = eCM - (pSumOut.e() - pSumOut.pz());
612  w2Rem = wPosRem * wNegRem;
613  if ( wPosRem < 0. || wNegRem < 0.
614  || sqrtpos(w2Rem) < sqrt(w2Beam[0]) + sqrt(w2Beam[1])) {
615  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
616  " kinematics construction failed owing to recoil mismatch");
617  return false;
618  }
619  }
620 
621  // Construct x rescaling factors for the two remants.
622  double lambdaRoot = sqrtpos( pow2(w2Rem - w2Beam[0] - w2Beam[1])
623  - 4. * w2Beam[0] * w2Beam[1] );
624  double rescaleA = (w2Rem + w2Beam[0] - w2Beam[1] + lambdaRoot)
625  / (2. * w2Rem * xSum[0]) ;
626  double rescaleB = (w2Rem + w2Beam[1] - w2Beam[0] + lambdaRoot)
627  / (2. * w2Rem * xSum[1]) ;
628 
629  // Construct energy and pz for remnants in first beam.
630  for (int iRem = nSys; iRem < beamA.size(); ++iRem) {
631  double pPos = rescaleA * beamA[iRem].x() * wPosRem;
632  double pNeg = beamA[iRem].mT2() / pPos;
633  beamA[iRem].e( 0.5 * (pPos + pNeg) );
634  beamA[iRem].pz( 0.5 * (pPos - pNeg) );
635 
636  // Add these partons to the normal event record.
637  int iNew = event.append( beamA[iRem].id(), 63, 1 + nOffset, 0, 0, 0,
638  beamA[iRem].col(), beamA[iRem].acol(), beamA[iRem].p(),
639  beamA[iRem].m() );
640  beamA[iRem].iPos( iNew);
641  }
642 
643  // Construct energy and pz for remnants in second beam.
644  for (int iRem = nSys; iRem < beamB.size(); ++iRem) {
645  double pNeg = rescaleB * beamB[iRem].x() * wNegRem;
646  double pPos = beamB[iRem].mT2() / pNeg;
647  beamB[iRem].e( 0.5 * (pPos + pNeg) );
648  beamB[iRem].pz( 0.5 * (pPos - pNeg) );
649 
650  // Add these partons to the normal event record.
651  int iNew = event.append( beamB[iRem].id(), 63, 2 + nOffset, 0, 0, 0,
652  beamB[iRem].col(), beamB[iRem].acol(), beamB[iRem].p(),
653  beamB[iRem].m() );
654  beamB[iRem].iPos( iNew);
655  }
656 
657  // Done.
658  return true;
659 
660 }
661 
662 //--------------------------------------------------------------------------
663 
664 // Allow colour reconnections by mergings of collision subsystems.
665 // iRec is system that may be reconnected, by moving its gluons to iSys,
666 // where minimal pT (or equivalently Lambda) is used to pick location.
667 // Therefore all dipoles in iSys have to be found, and all gluons in iRec.
668 // Matching q-qbar pairs are treated by analogy with gluons.
669 // Note: owing to rescatterings some outgoing partons must be skipped.
670 
671 bool BeamRemnants::reconnectColours( Event& event) {
672 
673  // References to beams to simplify indexing.
674  BeamParticle& beamA = *beamAPtr;
675  BeamParticle& beamB = *beamBPtr;
676 
677  // Prepare record of which systems should be merged onto another.
678  // The iSys system must have colour in final state to attach to it.
679  vector<int> iMerge(nSys);
680  vector<bool> hasColour(nSys);
681  for (int iSys = 0; iSys < nSys; ++iSys) {
682  iMerge[iSys] = iSys;
683  bool hasCol = false;
684  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
685  int iNow = partonSystemsPtr->getOut( iSys, iMem);
686  if (event[iNow].isFinal() && (event[iNow].col() > 0
687  || event[iNow].acol() > 0) ) {
688  hasCol = true;
689  break;
690  }
691  }
692  hasColour[iSys] = hasCol;
693  }
694 
695  // Loop over systems to decide which should be reconnected.
696  for (int iRec = nSys - 1; iRec > 0; --iRec) {
697 
698  // Determine reconnection strength from pT scale of system.
699  double pT2Rec = pow2( partonSystemsPtr->getPTHat(iRec) );
700  double probRec = pT20Rec / (pT20Rec + pT2Rec);
701 
702  // Loop over other systems iSys at higher pT scale and
703  // decide whether to reconnect the iRec gluons onto one of them.
704  for (int iSys = iRec - 1; iSys >= 0; --iSys)
705  if (hasColour[iSys] && probRec > rndmPtr->flat()) {
706 
707  // The iRec system and all merged with it to be merged with iSys.
708  iMerge[iRec] = iSys;
709  for (int iRec2 = iRec + 1; iRec2 < nSys; ++iRec2)
710  if (iMerge[iRec2] == iRec) iMerge[iRec2] = iSys;
711 
712  // Once a system has been merged do not test it anymore.
713  break;
714  }
715  }
716 
717  // Loop over systems. Check whether other systems to be merged with it.
718  for (int iSys = 0; iSys < nSys; ++iSys) {
719  int nMerge = 0;
720  for (int iRec = iSys + 1; iRec < nSys; ++iRec)
721  if (iMerge[iRec] == iSys) ++nMerge;
722  if (nMerge == 0) continue;
723 
724  // Incoming partons not counted if rescattered.
725  int iInASys = partonSystemsPtr->getInA(iSys);
726  bool hasInA = (beamA[iSys].isFromBeam());
727  int iInBSys = partonSystemsPtr->getInB(iSys);
728  bool hasInB = (beamB[iSys].isFromBeam());
729 
730  // Begin find dipoles in iSys system.
731  vector<BeamDipole> dipoles;
732  int sizeOut = partonSystemsPtr->sizeOut(iSys);
733  for (int iMem = 0; iMem < sizeOut; ++iMem) {
734 
735  // Find colour dipoles to beam remnant.
736  int iNow = partonSystemsPtr->getOut( iSys, iMem);
737  if (!event[iNow].isFinal()) continue;
738  int col = event[iNow].col();
739  if (col > 0) {
740  if (hasInA && event[iInASys].col() == col)
741  dipoles.push_back( BeamDipole( col, iNow, iInASys ) );
742  else if (hasInB && event[iInBSys].col() == col)
743  dipoles.push_back( BeamDipole( col, iNow, iInBSys ) );
744 
745  // Find colour dipole between final-state partons.
746  else for (int iMem2 = 0; iMem2 < sizeOut; ++iMem2)
747  if (iMem2 != iMem) {
748  int iNow2 = partonSystemsPtr->getOut( iSys, iMem2);
749  if (!event[iNow2].isFinal()) continue;
750  if (event[iNow2].acol() == col) {
751  dipoles.push_back( BeamDipole( col, iNow, iNow2) );
752  break;
753  }
754  }
755  }
756 
757  // Find anticolour dipoles to beam remnant.
758  int acol = event[iNow].acol();
759  if (acol > 0) {
760  if (hasInA && event[iInASys].acol() == acol)
761  dipoles.push_back( BeamDipole( acol, iInASys, iNow ) );
762  else if (hasInB && event[iInBSys].acol() == acol)
763  dipoles.push_back( BeamDipole( acol, iInBSys, iNow ) );
764  }
765  }
766 
767  // Skip mergings if no dipoles found.
768  if (dipoles.size() == 0) continue;
769 
770  // Find dipole sizes.
771  for (int iDip = 0; iDip < int(dipoles.size()); ++iDip)
772  dipoles[iDip].p1p2 = event[dipoles[iDip].iCol].p()
773  * event[dipoles[iDip].iAcol].p();
774 
775  // Loop over systems iRec to be merged with iSys.
776  for (int iRec = iSys + 1; iRec < nSys; ++iRec) {
777  if (iMerge[iRec] != iSys) continue;
778 
779  // Information on iRec. Vectors for gluons and anything else.
780  int sizeRec = partonSystemsPtr->sizeOut(iRec);
781  int iInARec = partonSystemsPtr->getInA(iRec);
782  int iInBRec = partonSystemsPtr->getInB(iRec);
783  int nGluRec = 0;
784  vector<int> iGluRec;
785  vector<double> pT2GluRec;
786  int nAnyRec = 0;
787  vector<int> iAnyRec;
788  vector<bool> freeAnyRec;
789 
790  // Copy of gluon positions in descending order.
791  for (int iMem = 0; iMem < sizeRec; ++iMem) {
792  int iNow = partonSystemsPtr->getOut( iRec, iMem);
793  if (!event[iNow].isFinal()) continue;
794  if (event[iNow].isGluon()) {
795  ++nGluRec;
796  iGluRec.push_back( iNow );
797  pT2GluRec.push_back( event[iNow].pT2() );
798  for (int i = nGluRec - 1; i > 1; --i) {
799  if (pT2GluRec[i - 1] > pT2GluRec[i]) break;
800  swap( iGluRec[i - 1], iGluRec[i] );
801  swap( pT2GluRec[i - 1], pT2GluRec[i] );
802  }
803  // Copy of anything else, mainly quarks, in no particular order.
804  } else {
805  ++nAnyRec;
806  iAnyRec.push_back( iNow );
807  freeAnyRec.push_back( true );
808  }
809  }
810 
811  // For each gluon in iRec now find the dipole that gives the smallest
812  // (pGlu * pI) (pGlu * pJ) / (pI * pJ), i.e. minimal pT (and Lambda).
813  for (int iGRec = 0; iGRec < nGluRec; ++iGRec) {
814  int iGlu = iGluRec[iGRec];
815  Vec4 pGlu = event[iGlu].p();
816  int iDipMin = 0;
817  double pT2DipMin = sCM;
818  for (int iDip = 0; iDip < int(dipoles.size()); ++iDip) {
819  double pT2Dip = (pGlu * event[dipoles[iDip].iCol].p())
820  * (pGlu * event[dipoles[iDip].iAcol].p()) / dipoles[iDip].p1p2;
821  if (pT2Dip < pT2DipMin) {
822  iDipMin = iDip;
823  pT2DipMin = pT2Dip;
824  }
825  }
826 
827  // Attach the gluon to the dipole, i.e. split the dipole in two.
828  int colGlu = event[iGlu].col();
829  int acolGlu = event[iGlu].acol();
830  int colDip = dipoles[iDipMin].col;
831  int iColDip = dipoles[iDipMin].iCol;
832  int iAcolDip = dipoles[iDipMin].iAcol;
833  event[iGlu].acol( colDip );
834  if (event[iAcolDip].acol() == colDip)
835  event[iAcolDip].acol( colGlu );
836  else event[iAcolDip].col( colGlu );
837  dipoles[iDipMin].iAcol = iGlu;
838  dipoles[iDipMin].p1p2 = event[iColDip].p() * pGlu;
839  dipoles.push_back( BeamDipole( colGlu, iGlu, iAcolDip ) );
840  dipoles.back().p1p2 = pGlu * event[iAcolDip].p();
841 
842  // Remove gluon from old system: reconnect colours.
843  for (int i = oldSize; i < event.size(); ++i)
844  if (i != iGlu && i != iAcolDip) {
845  if (event[i].isFinal()) {
846  if (event[i].acol() == colGlu) event[i].acol( acolGlu );
847  } else {
848  if (event[i].col() == colGlu) event[i].col( acolGlu );
849  }
850  }
851 
852  // Update any junction legs that match reconnected dipole.
853  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
854 
855  // Only junctions need to be updated, not antijunctions.
856  if (event.kindJunction(iJun) % 2 == 0) continue;
857  for (int leg = 0; leg < 3; ++leg) {
858  int col = event.colJunction(iJun, leg);
859  if (col == colDip)
860  event.colJunction(iJun, leg, colGlu);
861  }
862  }
863 
864  }
865 
866  // See if any matching quark-antiquark pairs among the rest.
867  for (int iQRec = 0; iQRec < nAnyRec; ++iQRec) {
868  int iQ = iAnyRec[iQRec];
869  int idQ = event[iQ].id();
870  if (freeAnyRec[iQRec] && idQ > 0 && idQ < 6)
871  for (int iQbarRec = 0; iQbarRec < nAnyRec; ++iQbarRec) {
872  int iQbar = iAnyRec[iQbarRec];
873  if (freeAnyRec[iQbarRec] && event[iQbar].id() == -idQ) {
874 
875  // Check that these can be traced back to same gluon splitting.
876  // For now also avoid qqbar pairs produced in rescatterings.??
877  int iTopQ = event.iTopCopyId(iQ);
878  int iTopQbar = event.iTopCopyId(iQbar);
879  int iMother = event[iTopQ].mother1();
880  if (event[iTopQbar].mother1() == iMother
881  && event[iMother].isGluon() && event[iMother].status() != -34
882  && event[iMother + 1].status() != -34 ) {
883 
884  // Now find the dipole that gives the smallest
885  // ((pQ + pQbar) * pI) ((pQ + pQbar) * pJ) / (pI * pJ).
886  Vec4 pGlu = event[iQ].p() + event[iQbar].p();
887  int iDipMin = 0;
888  double pT2DipMin = sCM;
889  for (int iDip = 0; iDip < int(dipoles.size()); ++iDip) {
890  double pT2Dip = (pGlu * event[dipoles[iDip].iCol].p())
891  * (pGlu * event[dipoles[iDip].iAcol].p())
892  / dipoles[iDip].p1p2;
893  if (pT2Dip < pT2DipMin) {
894  iDipMin = iDip;
895  pT2DipMin = pT2Dip;
896  }
897  }
898 
899  // Attach the q-qbar pair to the dipole, i.e. split the dipole.
900  int colGlu = event[iQ].col();
901  int acolGlu = event[iQbar].acol();
902  int colDip = dipoles[iDipMin].col;
903  int iColDip = dipoles[iDipMin].iCol;
904  int iAcolDip = dipoles[iDipMin].iAcol;
905  event[iQbar].acol( colDip );
906  if (event[iAcolDip].acol() == colDip)
907  event[iAcolDip].acol( colGlu );
908  else event[iAcolDip].col( colGlu );
909  dipoles[iDipMin].iAcol = iQbar;
910  dipoles[iDipMin].p1p2 = event[iColDip].p() * event[iQbar].p();
911  dipoles.push_back( BeamDipole( colGlu, iQ, iAcolDip ) );
912  dipoles.back().p1p2 = event[iQ].p() * event[iAcolDip].p();
913 
914  // Remove q-qbar pair from old system: reconnect colours.
915  freeAnyRec[iQRec] = false;
916  freeAnyRec[iQbarRec] = false;
917  for (int i = oldSize; i < event.size(); ++i)
918  if (i != iQRec && i != iQbarRec && i != iColDip
919  && i != iAcolDip) {
920  if (event[i].isFinal()) {
921  if (event[i].acol() == colGlu) event[i].acol( acolGlu );
922  } else {
923  if (event[i].col() == colGlu) event[i].col( acolGlu );
924  }
925  }
926 
927  // Update any junction legs that match reconnected dipole.
928  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
929 
930  // Only junctions need to be updated, not antijunctions.
931  if (event.kindJunction(iJun) % 2 == 0) continue;
932  for (int leg = 0; leg < 3; ++leg) {
933  int col = event.colJunction(iJun, leg);
934  if (col == colDip)
935  event.colJunction(iJun, leg, colGlu);
936  }
937  }
938 
939  // Done with processing of q-qbar pairs.
940  }
941  }
942  }
943  }
944 
945  // If only two beam gluons left of system, set their colour = anticolour.
946  // Used by BeamParticle::remnantColours to skip irrelevant gluons.
947  if ( event[iInARec].isGluon() && !event[iInARec].isRescatteredIncoming()
948  && event[iInBRec].isGluon() && !event[iInBRec].isRescatteredIncoming()
949  && event[iInARec].col() == event[iInBRec].acol()
950  && event[iInARec].acol() == event[iInBRec].col() ) {
951  event[iInARec].acol( event[iInARec].col() );
952  event[iInBRec].acol( event[iInBRec].col() );
953  }
954 
955  // End of loops over iRec and iSys systems.
956  }
957  }
958 
959  // Done.
960  return true;
961 
962 }
963 
964 //--------------------------------------------------------------------------
965 
966 // Collapse colours and check that they are consistent.
967 
968 bool BeamRemnants::checkColours( Event& event) {
969 
970  // No colours in lepton beams so no need to do anything.
971  if (beamAPtr->isLepton() && beamBPtr->isLepton()) return true;
972 
973  // Remove ambiguities when one colour collapses two ways.
974  // Resolve chains where one colour is mapped to another.
975  for (int iCol = 1; iCol < int(colFrom.size()); ++iCol)
976  for (int iColRef = 0; iColRef < iCol; ++iColRef) {
977  if (colFrom[iCol] == colFrom[iColRef]) {
978  colFrom[iCol] = colTo[iCol];
979  colTo[iCol] = colTo[iColRef];
980  }
981  if (colTo[iCol] == colFrom[iColRef]) colTo[iCol] = colTo[iColRef];
982  }
983 
984  // Transform event record colours from beam remnant colour collapses.
985  for (int i = oldSize; i < event.size(); ++i) {
986  int col = event[i].col();
987  int acol = event[i].acol();
988  for (int iCol = 0; iCol < int(colFrom.size()); ++iCol) {
989  if (col == colFrom[iCol]) {col = colTo[iCol]; event[i].col(col);}
990  if (acol == colFrom[iCol]) {acol = colTo[iCol]; event[i].acol(acol);}
991  }
992  }
993 
994  // Transform junction colours from beam remnant colour collapses.
995  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun)
996  for (int leg = 0; leg < 3; ++leg) {
997  int col = event.colJunction(iJun, leg);
998  for (int iCol = 0; iCol < int(colFrom.size()); ++iCol) {
999  if (col == colFrom[iCol]) {
1000  col = colTo[iCol];
1001  event.colJunction(iJun, leg, col);
1002  }
1003  }
1004  }
1005 
1006  // Arrays for current colours and anticolours, and for singlet gluons.
1007  vector<int> colList;
1008  vector<int> acolList;
1009  vector<int> iSingletGluon;
1010 
1011  // Find current colours and anticolours in the event record.
1012  for (int i = oldSize; i < event.size(); ++i)
1013  if (event[i].isFinal()) {
1014  int id = event[i].id();
1015  int col = event[i].col();
1016  int acol = event[i].acol();
1017  int colType = event[i].colType();
1018 
1019  // Quarks must have colour set, antiquarks anticolour, gluons both.
1020  if ( (id > 0 && id < 9 && (col <= 0 || acol != 0) )
1021  || (id < 0 && id > -9 && (col != 0 || acol <= 0) )
1022  || (id == 21 && (col <= 0 || acol <= 0) ) ) {
1023  infoPtr->errorMsg("Error in BeamRemnants::checkColours: "
1024  "q/qbar/g has wrong colour slots set");
1025  return false;
1026  }
1027 
1028  // Sextets must have one positive and one negative tag
1029  if ( (colType == 3 && (col <= 0 || acol >= 0))
1030  || (colType == -3 && (col >= 0 || acol <= 0)) ) {
1031  infoPtr->errorMsg("Error in BeamRemnants::checkColours: "
1032  "sextet has wrong colours");
1033  }
1034 
1035  // Save colours/anticolours, and position of colour singlet gluons.
1036  if ( col > 0) colList.push_back( col );
1037  if (acol > 0) acolList.push_back( acol );
1038  if (col > 0 && acol == col) iSingletGluon.push_back(i);
1039  // Colour sextets
1040  if ( col < 0) acolList.push_back( -col );
1041  if (acol < 0) colList.push_back( -acol );
1042  }
1043 
1044  // Run though list of singlet gluons and put them on final-state dipole
1045  // (i,j) that offers smallest (p_g p_i) * (p_g p_j) / (p_i p_j).
1046  for (int iS = 0; iS < int(iSingletGluon.size()); ++iS) {
1047  int iGlu = iSingletGluon[iS];
1048  int iAcolDip = -1;
1049  double pT2DipMin = sCM;
1050  for (int iC = oldSize; iC < event.size(); ++iC)
1051  if (iC != iGlu && event[iC].isFinal()) {
1052  int colDip = event[iC].col();
1053  if (colDip > 0 && event[iC].acol() !=colDip)
1054  for (int iA = oldSize; iA < event.size(); ++iA)
1055  if (iA != iGlu && iA != iC && event[iA].isFinal()
1056  && event[iA].acol() == colDip && event[iA].col() !=colDip) {
1057  double pT2Dip = (event[iGlu].p() * event[iC].p())
1058  * (event[iGlu].p() * event[iA].p())
1059  / (event[iC].p() * event[iA].p());
1060  if (pT2Dip < pT2DipMin) {
1061  iAcolDip = iA;
1062  pT2DipMin = pT2Dip;
1063  }
1064  }
1065  }
1066 
1067  // Fail if no dipole. Else insert singlet gluon onto relevant dipole.
1068  if (iAcolDip == -1) return false;
1069  event[iGlu].acol( event[iAcolDip].acol() );
1070  event[iAcolDip].acol( event[iGlu].col() );
1071 
1072  // Update any junction legs that match reconnected dipole.
1073  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
1074 
1075  // Only junctions need to be updated, not antijunctions.
1076  if (event.kindJunction(iJun) % 2 == 0) continue;
1077  for (int leg = 0; leg < 3; ++leg) {
1078  int col = event.colJunction(iJun, leg);
1079  if (col == event[iGlu].acol())
1080  event.colJunction(iJun, leg, event[iGlu].col());
1081  }
1082  }
1083 
1084  }
1085 
1086  // Check that not the same colour or anticolour appears twice.
1087  for (int iCol = 0; iCol < int(colList.size()) - 1; ++iCol) {
1088  int col = colList[iCol];
1089  for (int iCol2 = iCol + 1; iCol2 < int(colList.size()); ++iCol2)
1090  if (colList[iCol2] == col) {
1091  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1092  " colour appears twice");
1093  if (!ALLOWCOLOURTWICE) return false;
1094  }
1095  }
1096  for (int iAcol = 0; iAcol < int(acolList.size()) - 1; ++iAcol) {
1097  int acol = acolList[iAcol];
1098  for (int iAcol2 = iAcol + 1; iAcol2 < int(acolList.size()); ++iAcol2)
1099  if (acolList[iAcol2] == acol) {
1100  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1101  " anticolour appears twice");
1102  if (!ALLOWCOLOURTWICE) return false;
1103  }
1104  }
1105 
1106  // Remove all matching colour-anticolour pairs.
1107  bool foundPair = true;
1108  while (foundPair && colList.size() > 0 && acolList.size() > 0) {
1109  foundPair = false;
1110  for (int iCol = 0; iCol < int(colList.size()); ++iCol) {
1111  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol) {
1112  if (acolList[iAcol] == colList[iCol]) {
1113  colList[iCol] = colList.back(); colList.pop_back();
1114  acolList[iAcol] = acolList.back(); acolList.pop_back();
1115  foundPair = true;
1116  break;
1117  }
1118  }
1119  if (foundPair) break;
1120  }
1121  }
1122 
1123  // Check that remaining (anti)colours are accounted for by junctions.
1124  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
1125  int kindJun = event.kindJunction(iJun);
1126  for (int leg = 0; leg < 3; ++leg) {
1127  int colEnd = event.colJunction(iJun, leg);
1128 
1129  // Junction connected to three colours.
1130  if (kindJun == 1) {
1131  for (int iCol = 0; iCol < int(colList.size()); ++iCol)
1132  if (colList[iCol] == colEnd) {
1133  // Found colour match: remove and exit.
1134  colList[iCol] = colList.back();
1135  colList.pop_back();
1136  break;
1137  }
1138  }
1139 
1140  // Junction connected to three anticolours.
1141  else if (kindJun == 2) {
1142  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol)
1143  if (acolList[iAcol] == colEnd) {
1144  // Found colour match: remove and exit.
1145  acolList[iAcol] = acolList.back();
1146  acolList.pop_back();
1147  break;
1148  }
1149  }
1150 
1151  // Other junction types
1152  else if ( kindJun == 3 || kindJun == 5) {
1153  for (int iCol = 0; iCol < int(colList.size()); ++iCol)
1154  if (colList[iCol] == colEnd) {
1155  // Found colour match: remove and exit.
1156  colList[iCol] = colList.back();
1157  colList.pop_back();
1158  break;
1159  }
1160  }
1161 
1162  // Other antijunction types
1163  else if ( kindJun == 4 || kindJun == 6) {
1164  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol)
1165  if (acolList[iAcol] == colEnd) {
1166  // Found colour match: remove and exit.
1167  acolList[iAcol] = acolList.back();
1168  acolList.pop_back();
1169  break;
1170  }
1171  }
1172 
1173  // End junction check.
1174  }
1175  }
1176 
1177 
1178  // Repair step - sometimes needed when rescattering allowed.
1179  if (colList.size() > 0 || acolList.size() > 0) {
1180  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1181  " need to repair unmatched colours");
1182  }
1183  while (colList.size() > 0 && acolList.size() > 0) {
1184 
1185  // Replace one colour and one anticolour index by a new common one.
1186  int colMatch = colList.back();
1187  int acolMatch = acolList.back();
1188  int colNew = event.nextColTag();
1189  colList.pop_back();
1190  acolList.pop_back();
1191  for (int i = oldSize; i < event.size(); ++i)
1192  if (event[i].isFinal() && event[i].col() == colMatch) {
1193  event[i].col( colNew);
1194  break;
1195  }
1196  for (int i = oldSize; i < event.size(); ++i)
1197  if (event[i].isFinal() && event[i].acol() == acolMatch) {
1198  event[i].acol( colNew);
1199  break;
1200  }
1201  }
1202 
1203  // Done.
1204  return (colList.size() == 0 && acolList.size() == 0);
1205 
1206 }
1207 
1208 //==========================================================================
1209 
1210 } // end namespace Pythia8
Definition: beam.h:43
Definition: AgUStep.h:26