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) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, 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 BeamRemnants class.
16 
17 //--------------------------------------------------------------------------
18 
19 // Constants: could be changed here if desired, but normally should not.
20 // These are of technical nature, as described for each.
21 
22 // If same (anti)colour appears twice in final state, repair or reject.
23 const bool BeamRemnants::ALLOWCOLOURTWICE = true;
24 
25 // Maximum number of tries to match colours and kinematics in the event.
26 const int BeamRemnants::NTRYCOLMATCH = 10;
27 const int BeamRemnants::NTRYKINMATCH = 10;
28 
29 // Overall correction step for energy-momentum conservation; only
30 // becomes relevant in rescattering scenarios when FSR dipole emissions
31 // and primordial kT is added. Should hopefully not be needed.
32 const bool BeamRemnants::CORRECTMISMATCH = false;
33 
34 //--------------------------------------------------------------------------
35 
36 // Initialization.
37 
38 bool BeamRemnants::init( PartonVertexPtr partonVertexPtrIn,
39  ColRecPtr colourReconnectionPtrIn) {
40 
41  // Save pointers.
42  partonVertexPtr = partonVertexPtrIn;
43  colourReconnectionPtr = colourReconnectionPtrIn;
44 
45  // Width of primordial kT distribution.
46  doPrimordialKT = flag("BeamRemnants:primordialKT");
47  primordialKTsoft = parm("BeamRemnants:primordialKTsoft");
48  primordialKThard = parm("BeamRemnants:primordialKThard");
49  primordialKTremnant = parm("BeamRemnants:primordialKTremnant");
50  halfScaleForKT = parm("BeamRemnants:halfScaleForKT");
51  halfMassForKT = parm("BeamRemnants:halfMassForKT");
52  reducedKTatHighY = parm("BeamRemnants:reducedKTatHighY");
53 
54  // Handling of rescattering kinematics uncertainties from primodial kT.
55  allowRescatter = flag("MultipartonInteractions:allowRescatter");
56  doRescatterRestoreY = flag("BeamRemnants:rescatterRestoreY");
57 
58  // Choice of beam remnant and colour reconnection scenarios.
59  remnantMode = mode("BeamRemnants:remnantMode");
60  doReconnect = flag("ColourReconnection:reconnect");
61  reconnectMode = mode("ColourReconnection:mode");
62 
63  // Do multiparton interactions.
64  doMPI = flag("PartonLevel:MPI");
65 
66  // Flags for photoproduction from either (or both) side.
67  beamA2gamma = flag("PDF:beamA2gamma");
68  beamB2gamma = flag("PDF:beamB2gamma");
69 
70  // Check that remnant model and colour reconnection model work together.
71  if (remnantMode == 1 && reconnectMode == 0) {
72  infoPtr->errorMsg("Abort from BeamRemnants::init: The remnant model"
73  " and colour reconnection model does not work together");
74  return false;
75  }
76 
77  // Total and squared CM energy at nominal energy.
78  eCM = infoPtr->eCM();
79  sCM = eCM * eCM;
80 
81  // Initialize junction splitting class.
82  junctionSplitting.init();
83 
84  // Possibility to set parton vertex information.
85  doPartonVertex = flag("PartonVertex:setVertex")
86  && (partonVertexPtr != 0);
87 
88  // Done.
89  return true;
90 
91 }
92 
93 //--------------------------------------------------------------------------
94 
95 // Select the flavours/kinematics/colours of the two beam remnants.
96 // Notation: iPar = all partons, iSys = matched systems of two beams,
97 // iRem = additional partons in remnants.
98 
99 bool BeamRemnants::add( Event& event, int iFirst, bool doDiffCR) {
100 
101  // Update to current CM energy.
102  eCM = infoPtr->eCM();
103  sCM = eCM * eCM;
104 
105  // Check that flavour bookkept in event and in beam particles agree.
106  for (int i = 0; i < beamAPtr->size(); ++i) {
107  int j = (*beamAPtr)[i].iPos();
108  if ((*beamAPtr)[i].id() != event[j].id()) {
109  infoPtr->errorMsg("Error in BeamRemnants::add: "
110  "event and beam flavours do not match");
111  return false;
112  }
113  }
114  for (int i = 0; i < beamBPtr->size(); ++i) {
115  int j = (*beamBPtr)[i].iPos();
116  if ((*beamBPtr)[i].id() != event[j].id()) {
117  infoPtr->errorMsg("Error in BeamRemnants::add: "
118  "event and beam flavours do not match");
119  return false;
120  }
121  }
122 
123  // Deeply inelastic scattering needs special remnant handling.
124  // Here getGammaMode separates from photoproduction.
125  isDIS = (beamAPtr->isLepton() && !beamBPtr->isLepton()
126  && beamAPtr->getGammaMode() == 0)
127  || (beamBPtr->isLepton() && !beamAPtr->isLepton()
128  && beamBPtr->getGammaMode() == 0);
129 
130  // Number of scattering subsystems. Size of event record before treatment.
131  nSys = partonSystemsPtr->sizeSys();
132  oldSize = event.size();
133 
134  // Store event as it was before adding anything.
135  Event eventSave = event;
136  BeamParticle beamAsave = (*beamAPtr);
137  BeamParticle beamBsave = (*beamBPtr);
138  PartonSystems partonSystemsSave = (*partonSystemsPtr);
139 
140  // Two different methods to add the beam remnants.
141  if (remnantMode == 0) {
142  if (!addOld(event)) return false;
143  } else
144  if (!addNew(event)) return false;
145  if (isDIS) return true;
146 
147  // Store event before doing colour reconnections.
148  Event eventTmpSave = event;
149  bool colCorrect = false;
150  for (int i = 0; i < 10; ++i) {
151  if (doReconnect && doDiffCR
152  && (reconnectMode == 1 || reconnectMode == 2)) {
153  colourReconnectionPtr->next(event, iFirst);
154 
155  // Check that the new colour structure is physical.
156  if (!junctionSplitting.checkColours(event))
157  event = eventTmpSave;
158  else {
159  colCorrect = true;
160  break;
161  }
162  // If no colour reconnection, just check the colour configuration once.
163  } else {
164  if (junctionSplitting.checkColours(event))
165  colCorrect = true;
166  break;
167  }
168  }
169 
170  // Possibility to add vertex information to beam particles and remnants.
171  if (doPartonVertex) for (int iBeam = 0; iBeam < 2; ++iBeam) {
172  BeamParticle& beamNow = (iBeam == 0) ? *beamAPtr : *beamBPtr;
173  vector<int> iRemn, iInit;
174  for (int i = beamNow.sizeInit(); i < beamNow.size(); ++i)
175  iRemn.push_back( beamNow[i].iPos() );
176  for (int i = 0; i < beamNow.sizeInit(); ++i)
177  iInit.push_back( beamNow[i].iPos() );
178  partonVertexPtr->vertexBeam(iBeam, iRemn, iInit, event);
179  }
180 
181  // Restore event and return false if colour reconnection failed.
182  if (!colCorrect) {
183  event = eventSave;
184  (*beamAPtr) = beamAsave;
185  (*beamBPtr) = beamBsave;
186  (*partonSystemsPtr) = partonSystemsSave;
187  infoPtr->errorMsg("Error in BeamRemnants::add: "
188  "failed to find physical colour state after colour reconnection");
189  return false;
190  }
191 
192  // Done.
193  return true;
194 }
195 
196 //--------------------------------------------------------------------------
197 
198 // Old function for adding beam remnant.
199 
200 bool BeamRemnants::addOld( Event& event) {
201 
202  // Add required extra remnant flavour content.
203  // Start all over if fails (in option where junctions not allowed).
204  if ( !beamAPtr->remnantFlavours(event, isDIS)
205  || !beamBPtr->remnantFlavours(event, isDIS) ) {
206  infoPtr->errorMsg("Error in BeamRemnants::add:"
207  " remnant flavour setup failed");
208  return false;
209  }
210 
211  // Do the kinematics of the collision subsystems and two beam remnants.
212  if (!setKinematics(event)) return false;
213 
214  // Allow colour reconnections.
215  if (doReconnect && reconnectMode == 0 && !isDIS)
216  colourReconnectionPtr->next(event, oldSize);
217 
218  // Save current modifiable colour configuration for fast restoration.
219  vector<int> colSave;
220  vector<int> acolSave;
221  for (int i = oldSize; i < event.size(); ++i) {
222  colSave.push_back( event[i].col() );
223  acolSave.push_back( event[i].acol() );
224  }
225  event.saveJunctionSize();
226 
227  // Allow several tries to match colours of initiators and remnants.
228  // Frequent "failures" since shortcutting colours separately on
229  // the two event sides may give "colour singlet gluons" etc.
230  bool physical = true;
231  for (int iTry = 0; iTry < NTRYCOLMATCH; ++iTry) {
232  physical = true;
233 
234  // Reset list of colour "collapses" (transformations).
235  colFrom.resize(0);
236  colTo.resize(0);
237 
238  // First process each set of beam colours on its own.
239  if (!beamAPtr->remnantColours(event, colFrom, colTo))
240  physical = false;
241  if (!beamBPtr->remnantColours(event, colFrom, colTo))
242  physical = false;
243 
244  // Then check that colours and anticolours are matched in whole event.
245  if ( physical && !checkColours(event) ) physical = false;
246 
247  // If no problems then done, else restore and loop.
248  if (physical) break;
249  for (int i = oldSize; i < event.size(); ++i)
250  event[i].cols( colSave[i - oldSize], acolSave[i - oldSize] );
251  event.restoreJunctionSize();
252  infoPtr->errorMsg("Warning in BeamRemnants::addOld:"
253  " colour tracing failed; will try again");
254  }
255 
256  // If no solution after several tries then failed.
257  if (!physical) {
258  infoPtr->errorMsg("Error in BeamRemnants::addOld:"
259  " colour tracing failed after several attempts");
260  return false;
261  }
262 
263  // Done.
264  return true;
265 }
266 
267 //--------------------------------------------------------------------------
268 
269 // New function for adding beam remnant.
270 
271 bool BeamRemnants::addNew( Event& event) {
272 
273  // Start by saving a copy of the event, if the beam remnant fails.
274  Event eventSave = event;
275  BeamParticle beamAsave = (*beamAPtr);
276  BeamParticle beamBsave = (*beamBPtr);
277  PartonSystems partonSystemsSave = (*partonSystemsPtr);
278 
279  // Do several tries in case an unphysical colour contruction is made.
280  bool beamRemnantFound = false;
281  int nMaxTries = 10;
282 
283  for (int iTries = 0;iTries < nMaxTries; ++iTries) {
284 
285  // Set the initial colours.
286  beamAPtr->setInitialCol(event);
287  beamBPtr->setInitialCol(event);
288 
289  // Find colour state of outgoing partons and reconnect colours to match it.
290  beamAPtr->findColSetup(event);
291  beamBPtr->updateCol(beamAPtr->getColUpdates());
292 
293  beamBPtr->findColSetup(event);
294  beamAPtr->updateCol(beamBPtr->getColUpdates());
295 
296  // Add beam remnants.
297  beamAPtr->remnantFlavoursNew(event);
298  beamBPtr->remnantFlavoursNew(event);
299 
300  // It is possible junctions were added, so update list.
301  event.saveJunctionSize();
302 
303  // Do the kinematics of the collision subsystems and two beam remnants.
304  if (!setKinematics(event)) {
305  // If it does not work, try parton level again.
306  event = eventSave;
307  (*beamAPtr) = beamAsave;
308  (*beamBPtr) = beamBsave;
309  (*partonSystemsPtr) = partonSystemsSave;
310  return false;
311  }
312 
313  // Update the colour changes for all final state particles.
314  updateColEvent(event, beamAPtr->getColUpdates());
315  updateColEvent(event, beamBPtr->getColUpdates());
316 
317  // Check whether the new colour structure can be accepted.
318  if (junctionSplitting.checkColours(event)) {
319  beamRemnantFound = true;
320  break;
321  }
322 
323  // If failed, restore earlier configuration and try to find new
324  // colour structure.
325  else {
326  event = eventSave;
327  (*beamAPtr) = beamAsave;
328  (*beamBPtr) = beamBsave;
329  (*partonSystemsPtr) = partonSystemsSave;
330  }
331  }
332 
333  // Return if it was not possible to find physical colour structure.
334  if (!beamRemnantFound) {
335  infoPtr->errorMsg("Error in BeamRemnants::addNew: "
336  "failed to find physical colour structure");
337  // Restore event to previous state.
338  event = eventSave;
339  (*beamAPtr) = beamAsave;
340  (*beamBPtr) = beamBsave;
341  (*partonSystemsPtr) = partonSystemsSave;
342  return false;
343  }
344 
345  // Done.
346  return true;
347 }
348 
349 //--------------------------------------------------------------------------
350 
351 // Set up trial transverse and longitudinal kinematics for each beam
352 // separately. Final decisions involve comparing the two beams.
353 
354 bool BeamRemnants::setKinematics( Event& event) {
355 
356  // References to beams to simplify indexing.
357  BeamParticle& beamA = *beamAPtr;
358  BeamParticle& beamB = *beamBPtr;
359 
360  // Simple handling of lepton-lepton scattering.
361  if (beamA.isLepton() && beamB.isLepton()) {
362  if (!beamA.isUnresolvedLepton()) {
363  double eGamma = max(0., event[1].e() - event[event[1].daughter1()].e() );
364  event.append( 22, 63, 1, 0, 0, 0, 0, 0, 0., 0., eGamma, eGamma, 0.);
365  }
366  if (!beamB.isUnresolvedLepton()) {
367  double eGamma = max(0., event[2].e() - event[event[2].daughter1()].e() );
368  event.append( 22, 63, 2, 0, 0, 0, 0, 0, 0., 0., -eGamma, eGamma, 0.);
369  }
370  return true;
371  }
372 
373  // Photon from photon beam is unresolved.
374  if ( beamA.isGamma() && beamA[0].id() == 22 ) {
375  beamA.resolvedGamma(false);
376  }
377  if ( beamB.isGamma() && beamB[0].id() == 22 ) {
378  beamB.resolvedGamma(false);
379  }
380 
381  // Check if photon beams need remnants.
382  bool beamAisGamma = beamA.isGamma();
383  bool beamBisGamma = beamB.isGamma();
384  bool gammaAResolved = beamA.resolvedGamma();
385  bool gammaBResolved = beamB.resolvedGamma();
386  bool gammaOneResolved = false;
387 
388  // Do nothing if ISR have found the original beam photons and only one MPI.
389  if ( ( beamAisGamma && !gammaAResolved && beamBisGamma && !gammaBResolved )
390  && infoPtr->nMPI() == 1 )
391  return true;
392 
393  // Check that has not already used up beams. Unless direct photon in photon.
394  if ( beamAisGamma && beamA[0].id() == 22 ) { }
395  else if ( beamBisGamma && beamB[0].id() == 22 ) { }
396  else if ( ( !(beamA.isLepton() || (beamAisGamma && !gammaAResolved) )
397  && beamA.xMax(-1) <= 0.) ||
398  ( !(beamB.isLepton() || (beamBisGamma && !gammaBResolved) )
399  && beamB.xMax(-1) <= 0.) ) {
400  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
401  " no momentum left for beam remnants");
402  return false;
403  }
404 
405  // Check if exactly one remnant for photon-photon collisions.
406  if ( beamAisGamma && beamBisGamma && gammaAResolved != gammaBResolved)
407  gammaOneResolved = true;
408 
409  // Unresolved photon + hadron.
410  if ( ( (beamAisGamma && !gammaAResolved)
411  || (beamBisGamma && !gammaBResolved) )
412  && !(beamAisGamma && beamBisGamma) )
413  gammaOneResolved = true;
414 
415  // Unresolved photon on one side.
416  if ( ( beamA.getGammaMode() == 2 && !beamB2gamma )
417  || ( beamB.getGammaMode() == 2 && !beamA2gamma ) )
418  gammaOneResolved = true;
419 
420  // Special kinematics setup for one-remnant systems (DIS).
421  if ( (gammaOneResolved && infoPtr->nMPI() == 1) || isDIS )
422  return setOneRemnKinematics(event);
423 
424  // Last beam-status particles. Offset relative to normal beam locations.
425  int nBeams = 3;
426  for (int i = 3; i < event.size(); ++i)
427  if (event[i].statusAbs() < 20) nBeams = i + 1;
428  int nOffset = nBeams - 3;
429 
430  // If extra photons in event fix the offset.
431  if ( beamA2gamma || beamB2gamma ) {
432  if (beamA.hasResGamma()) --nOffset;
433  if (beamB.hasResGamma()) --nOffset;
434  }
435 
436  // Reserve space for extra information on the systems and beams.
437  int nMaxBeam = max( beamA.size(), beamB.size() );
438  vector<double> sHatSys(nMaxBeam);
439  vector<double> kTwidth(nMaxBeam);
440  vector<double> kTcomp(nMaxBeam);
441  vector<RotBstMatrix> Msys(nSys);
442 
443  // Loop over all subsystems. Default values. Find invariant mass.
444  double kTcompSumA = 0.;
445  double kTcompSumB = 0.;
446  for (int iSys = 0; iSys < nSys; ++iSys) {
447  // Safety check. All parton systems now remaining should be 2->N.
448  if (!partonSystemsPtr->hasInAB(iSys)) {
449  stringstream ss;
450  ss << "iSys = "<<iSys;
451  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
452  " Encountered parton system without beam partons.",ss.str());
453  return false;
454  }
455  double kTwidthNow = 0.;
456  double mHatDamp = 1.;
457  int iInA = partonSystemsPtr->getInA(iSys);
458  int iInB = partonSystemsPtr->getInB(iSys);
459  double sHatNow = (event[iInA].p() + event[iInB].p()).m2Calc();
460 
461  // Set width of primordial kT distribution.
462  if (doPrimordialKT) {
463  // Les Houches events use primordialKThard.
464  if (iSys == 0 && infoPtr->isLHA()) {
465  kTwidthNow = primordialKThard;
466  }
467  // Internal processes and MPI use pT-dependent interpolation between
468  // primordialKThard and primordialKTsoft.
469  else {
470  // For hardest interaction pT -> renormalization scale so also 2 -> 1.
471  double scale = (iSys == 0) ? infoPtr->QRen(iDS)
472  : partonSystemsPtr->getPTHat(iSys);
473  kTwidthNow = (halfScaleForKT * primordialKTsoft
474  + scale * primordialKThard) / (halfScaleForKT + scale);
475  }
476  // Dampen primordial kT width for very low masses / extreme rapidities.
477  double mHat = sqrt(sHatNow);
478  double yDamp =
479  pow( (event[iInA].e() + event[iInB].e()) / mHat, reducedKTatHighY );
480  mHatDamp = mHat / (mHat + halfMassForKT * yDamp);
481  kTwidthNow *= mHatDamp;
482  }
483 
484  // Store properties of compensation systems and total compensation power.
485  // Rescattered partons do not compensate, but may be massive.
486  sHatSys[iSys] = sHatNow;
487  kTwidth[iSys] = kTwidthNow ;
488  kTcomp[iSys] = mHatDamp;
489  if (beamA[iSys].isFromBeam()) kTcompSumA += mHatDamp;
490  else beamA[iSys].m( event[iInA].m() );
491  if (beamB[iSys].isFromBeam()) kTcompSumB += mHatDamp;
492  else beamB[iSys].m( event[iInB].m() );
493  }
494 
495  // Primordial kT and compensation power among remnants.
496  double kTwidthNow = (doPrimordialKT) ? primordialKTremnant : 0.;
497  for (int iRem = nSys; iRem < nMaxBeam; ++iRem) {
498  sHatSys[iRem] = 0.;
499  kTwidth[iRem] = kTwidthNow ;
500  kTcomp[iRem] = 1.;
501  }
502  kTcompSumA += beamA.size() - nSys;
503  kTcompSumB += beamB.size() - nSys;
504 
505  // Allow ten tries to construct kinematics (but normally works first).
506  bool physical;
507  double xSum[2], xInvM[2], w2Beam[2], wPosRem, wNegRem, w2Rem;
508  for (int iTry = 0; iTry < NTRYKINMATCH; ++iTry) {
509  physical = true;
510 
511  // Loop over the two beams.
512  for (int iBeam = 0; iBeam < 2; ++iBeam) {
513  BeamParticle& beam = (iBeam == 0) ? beamA : beamB;
514  int nPar = beam.size();
515 
516  // Generate Gaussian pT for initiator partons inside hadrons.
517  // Store/restore rescattered parton momenta before primordial kT.
518  if ( (beam.isHadron() || beam.isGamma()) && doPrimordialKT ) {
519  double pxSum = 0.;
520  double pySum = 0.;
521  for (int iPar = 0; iPar < nPar; ++iPar) {
522  // Do nothing if initiator parton from gamma->qqbar
523  // where gamma is the original beam particle or for
524  // the other valence arising from the same splitting.
525  if ( iPar == beam.gamVal() || ( beam.gamVal() >= 0 &&
526  ( iPar != beam.gamVal() && beam[iPar].isValence()) ) ) ;
527  else if ( beam[iPar].isFromBeam() ) {
528  pair<double, double> gauss2 = rndmPtr->gauss2();
529  double px = kTwidth[iPar] * gauss2.first;
530  double py = kTwidth[iPar] * gauss2.second;
531  beam[iPar].px(px);
532  beam[iPar].py(py);
533  pxSum += px;
534  pySum += py;
535  } else {
536  int iInAB = (iBeam == 0) ? partonSystemsPtr->getInA(iPar)
537  : partonSystemsPtr->getInB(iPar);
538  beam[iPar].p( event[iInAB].p() );
539  }
540  }
541 
542  // Share recoil between all initiator partons, rescatterers excluded.
543  double kTcompSum = (iBeam == 0) ? kTcompSumA : kTcompSumB;
544 
545  // If parton from gamma->qqbar in ISR, set pT accordingly.
546  if ( beam.gamVal() >= 0 ) {
547 
548  // Get the pT2 value for the given gamma->qqbar splitting.
549  double pT2corr = beam.pT2gamma2qqbar();
550 
551  // Sample the direction and set px and py.
552  double phi = 2*M_PI*rndmPtr->flat();
553  double px = cos(phi) * sqrt(pT2corr);
554  double py = sin(phi) * sqrt(pT2corr);
555  beam[beam.gamVal()].px(px);
556  beam[beam.gamVal()].py(py);
557  pxSum += px;
558  pySum += py;
559 
560  // Do not change the momentum of these partons when sharing recoil.
561  kTcompSum -= kTcomp[beam.gamVal()];
562 
563  // Other valence parton takes (most of) the recoil but this will
564  // still be part of pT balancing.
565  for (int iPar = 0; iPar < nPar; ++iPar)
566  if ( iPar != beam.gamVal() && beam[iPar].isValence() ) {
567  beam[iPar].px(-px);
568  beam[iPar].py(-py);
569  pxSum += -px;
570  pySum += -py;
571  break;
572  }
573  }
574 
575  for (int iPar = 0; iPar < nPar; ++iPar) {
576  if (beam[iPar].isFromBeam() && iPar != beam.gamVal() ) {
577  beam[iPar].px( beam[iPar].px() - pxSum * kTcomp[iPar] / kTcompSum);
578  beam[iPar].py( beam[iPar].py() - pySum * kTcomp[iPar] / kTcompSum);
579  }
580  }
581 
582  // Without primordial kT: still need to store rescattered partons.
583  } else if ( beam.isHadron() || beam.isGamma() ) {
584  for (int iPar = 0; iPar < nPar; ++iPar) {
585  if ( !beam[iPar].isFromBeam() ) {
586  int iInAB = (iBeam == 0) ? partonSystemsPtr->getInA(iPar)
587  : partonSystemsPtr->getInB(iPar);
588  beam[iPar].p( event[iInAB].p() );
589  }
590  }
591 
592  // And set the pT in case gamma->qqbar splitting with MPIs.
593  if ( beam.gamVal() >= 0 ) {
594  double phi = 2*M_PI*rndmPtr->flat();
595  double px = cos(phi) * sqrt(beam.pT2gamma2qqbar());
596  double py = sin(phi) * sqrt(beam.pT2gamma2qqbar());
597  beam[beam.gamVal()].px(px);
598  beam[beam.gamVal()].py(py);
599 
600  // Other valence parton takes the recoil but this will still be
601  // part of pT balancing.
602  for (int iPar = 0; iPar < nPar; ++iPar) {
603  if ( iPar != beam.gamVal() && beam[iPar].isValence() ) {
604  beam[iPar].px(-px);
605  beam[iPar].py(-py);
606  break;
607  }
608  }
609  }
610 
611  }
612 
613  // Pick unrescaled x values for remnants. Sum up (unscaled) p+ and p-.
614  xSum[iBeam] = 0.;
615  xInvM[iBeam] = 0.;
616  for (int iRem = nSys; iRem < nPar; ++iRem) {
617  double xPrel = beam.xRemnant( iRem);
618  beam[iRem].x(xPrel);
619  xSum[iBeam] += xPrel;
620  xInvM[iBeam] += beam[iRem].mT2()/xPrel;
621  }
622 
623  // Squared transverse mass for each beam, using lightcone x.
624  w2Beam[iBeam] = xSum[iBeam] * xInvM[iBeam];
625 
626  // End separate treatment of the two beams.
627  }
628 
629  // Recalculate kinematics of initiator systems with primordial kT.
630  wPosRem = eCM;
631  wNegRem = eCM;
632  for (int iSys = 0; iSys < nSys; ++iSys) {
633  int iA = beamA[iSys].iPos();
634  int iB = beamB[iSys].iPos();
635  double sHat = sHatSys[iSys];
636 
637  // Classify system: rescattering on either or both sides?
638  bool normalA = beamA[iSys].isFromBeam();
639  bool normalB = beamB[iSys].isFromBeam();
640  bool normalSys = normalA && normalB;
641  bool doubleRes = !normalA && !normalB;
642 
643  // Check whether final-state system momentum matches initial-state one.
644  if (allowRescatter && CORRECTMISMATCH) {
645  Vec4 pInitial = event[iA].p() + event[iB].p();
646  Vec4 pFinal;
647  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
648  int iAB = partonSystemsPtr->getOut(iSys, iMem);
649  if (event[iAB].isFinal()) pFinal += event[iAB].p();
650  }
651 
652  // Scale down primordial kT from side A if p+ increased.
653  if (normalA && pFinal.pPos() > pInitial.pPos())
654  beamA[iSys].scalePT( pInitial.pPos() / pFinal.pPos() );
655 
656  // Scale down primordial kT from side B if p- increased.
657  if (normalB && pFinal.pNeg() > pInitial.pNeg())
658  beamB[iSys].scalePT( pInitial.pNeg() / pFinal.pNeg() );
659  }
660 
661  // Rescatter: possible change in sign of lightcone momentum of a
662  // rescattered parton. If this happens, try to pick
663  // new primordial kT values
664  if (allowRescatter
665  && (event[iA].pPos() / beamA[iSys].pPos() < 0
666  || event[iB].pNeg() / beamB[iSys].pNeg() < 0) ) {
667  infoPtr->errorMsg("Warning in BeamRemnants::setKinematics:"
668  " change in lightcone momentum sign; retrying kinematics");
669  physical = false;
670  break;
671  }
672 
673  // Begin kinematics of partons after primordial kT has been added.
674  double sHatTAft = sHat + pow2( beamA[iSys].px() + beamB[iSys].px())
675  + pow2( beamA[iSys].py() + beamB[iSys].py());
676  double w2A = beamA[iSys].mT2();
677  double w2B = beamB[iSys].mT2();
678  double w2Diff = sHatTAft - w2A - w2B;
679  double lambda = pow2(w2Diff) - 4. * w2A * w2B;
680 
681  // Too large transverse momenta means that kinematics will not work.
682  if (lambda <= 0.) {
683  physical = false;
684  break;
685  }
686  double lamRoot = sqrtpos( lambda );
687 
688  // Mirror solution if the two incoming have reverse rapidity ordering.
689  if (allowRescatter && doubleRes && (event[iA].pPos() * event[iB].pNeg()
690  < event[iA].pNeg() * event[iB].pPos()) ) lamRoot = -lamRoot;
691 
692  // Two procedures, which agree for normal scattering, separate here.
693  // First option keeps rapidity (and mass) of system unchanged by
694  // primordial kT, by modification of rescattered parton.
695  if (normalSys || doRescatterRestoreY || doubleRes) {
696 
697  // Find kinematics of initial system: transverse mass, and
698  // longitudinal momentum carried by all or rescattered partons.
699  double sHatTBef = sHat;
700  double wPosBef, wNegBef, wPosBefRes, wNegBefRes;
701  // Normal scattering.
702  if (normalSys) {
703  wPosBef = beamA[iSys].x() * eCM;
704  wNegBef = beamB[iSys].x() * eCM;
705  wPosBefRes = 0.;
706  wNegBefRes = 0.;
707  // Rescattering on side A.
708  } else if (normalB) {
709  sHatTBef += event[iA].pT2();
710  wPosBef = event[iA].pPos();
711  wNegBef = event[iA].pNeg() + beamB[iSys].x() * eCM;
712  wPosBefRes = beamA[iSys].pPos();
713  wNegBefRes = beamA[iSys].pNeg();
714  // Rescattering on side B.
715  } else if (normalA) {
716  sHatTBef += event[iB].pT2();
717  wPosBef = beamA[iSys].x() * eCM + event[iB].pPos();
718  wNegBef = event[iB].pNeg();
719  wPosBefRes = beamB[iSys].pPos();
720  wNegBefRes = beamB[iSys].pNeg();
721  // Rescattering on both sides.
722  } else {
723  sHatTBef += pow2( event[iA].px() + event[iB].px())
724  + pow2( event[iA].py() + event[iB].py());
725  wPosBef = event[iA].pPos() + event[iB].pPos();
726  wNegBef = event[iA].pNeg() + event[iB].pNeg();
727  wPosBefRes = beamA[iSys].pPos() + beamB[iSys].pPos();
728  wNegBefRes = beamA[iSys].pNeg() + beamB[iSys].pNeg();
729  }
730 
731  // Rescale outgoing momenta to keep same mass and rapidity of system
732  // as originally, and solve for kinematics.
733  double rescale = sqrt(sHatTAft / sHatTBef);
734  double wPosAft = rescale * wPosBef;
735  double wNegAft = rescale * wNegBef;
736  wPosRem -= wPosAft - wPosBefRes;
737  wNegRem -= wNegAft - wNegBefRes;
738  double wPosA = 0.5 * (sHatTAft + w2A - w2B + lamRoot) / wNegAft;
739  double wNegB = 0.5 * (sHatTAft + w2B - w2A + lamRoot) / wPosAft;
740 
741  // Store modified beam parton momenta.
742  beamA[iSys].e( 0.5 * (wPosA + w2A / wPosA) );
743  beamA[iSys].pz( 0.5 * (wPosA - w2A / wPosA) );
744  beamB[iSys].e( 0.5 * (w2B / wNegB + wNegB) );
745  beamB[iSys].pz( 0.5 * (w2B / wNegB - wNegB) );
746 
747  // Second option keeps rescattered parton (and system mass) unchanged
748  // by primordial kT, by modification of system rapidity.
749  } else {
750 
751  // Rescattering on side A: preserve already scattered parton.
752  if (normalB) {
753  double wPosA = beamA[iSys].pPos();
754  double wNegB = 0.5 * (w2Diff + lamRoot) / wPosA;
755  beamB[iSys].e( 0.5 * (w2B / wNegB + wNegB) );
756  beamB[iSys].pz( 0.5 * (w2B / wNegB - wNegB) );
757  wPosRem -= w2B / wNegB;
758  wNegRem -= wNegB;
759 
760 
761  // Rescattering on side B: preserve already scattered parton.
762  } else if (normalA) {
763  double wNegB = beamB[iSys].pNeg();
764  double wPosA = 0.5 * (w2Diff + lamRoot) / wNegB;
765  beamA[iSys].e( 0.5 * (wPosA + w2A / wPosA) );
766  beamA[iSys].pz( 0.5 * (wPosA - w2A / wPosA) );
767  wPosRem -= wPosA;
768  wNegRem -= w2A / wPosA;
769 
770  // Primordial kT in double rescattering does change the mass of
771  // the system without any possible fix in this framework, which
772  // leads to incorrect boosts. Current choice is therefore always
773  // to handle it with the first procedure, where mass is retained.
774  } else {
775  }
776  }
777 
778  // Construct system rotation and boost caused by primordial kT.
779  Msys[iSys].reset();
780  Msys[iSys].toCMframe( event[iA].p(), event[iB].p() );
781  Msys[iSys].fromCMframe( beamA[iSys].p(), beamB[iSys].p() );
782 
783  // Boost rescattered partons in subsequent beam A list.
784  for (int iSys2 = iSys + 1; iSys2 < nSys; ++iSys2) {
785  if (!beamA[iSys2].isFromBeam()) {
786  int iBefResc = event[ beamA[iSys2].iPos() ].mother1();
787  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem)
788  if (partonSystemsPtr->getOut(iSys, iMem) == iBefResc) {
789  Vec4 pTemp = event[iBefResc].p();
790  pTemp.rotbst( Msys[iSys] );
791  beamA[iSys2].p( pTemp );
792  }
793  }
794 
795  // Boost rescattered partons in subsequent beam B list.
796  if (!beamB[iSys2].isFromBeam()) {
797  int iBefResc = event[ beamB[iSys2].iPos() ].mother1();
798  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem)
799  if (partonSystemsPtr->getOut(iSys, iMem) == iBefResc) {
800  Vec4 pTemp = event[iBefResc].p();
801  pTemp.rotbst( Msys[iSys] );
802  beamB[iSys2].p( pTemp );
803  }
804  }
805  }
806  }
807 
808  // Check that remaining momentum is enough for remnants.
809  if (wPosRem < 0. || wNegRem < 0.) physical = false;
810  w2Rem = wPosRem * wNegRem;
811  if (sqrtpos(w2Rem) < sqrt(w2Beam[0]) + sqrt(w2Beam[1]))
812  physical = false;
813 
814  // End of loop over ten tries. Do not loop when solution found.
815  if (physical) break;
816  }
817 
818  // If no solution after ten tries then failed.
819  if (!physical) {
820  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
821  " kinematics construction failed");
822  return false;
823  }
824 
825  // For successful initiator kinematics process whole systems.
826  Vec4 pSumOut;
827  for (int iSys = 0; iSys < nSys; ++iSys) {
828 
829  // Copy initiators and their systems and boost them accordingly.
830  // Update subsystem and beams info on new positions of partons.
831  // Update daughter info of mothers, i.e. of beams, for hardest interaction.
832  if (beamA[iSys].isFromBeam()) {
833  int iA = beamA[iSys].iPos();
834  int iAcopy = event.copy(iA, -61);
835  event[iAcopy].rotbst(Msys[iSys], false);
836  partonSystemsPtr->setInA(iSys, iAcopy);
837  beamA[iSys].iPos( iAcopy);
838  if (iSys == 0) {
839  int mother = event[iAcopy].mother1();
840  event[mother].daughter1(iAcopy);
841  }
842  }
843  if (beamB[iSys].isFromBeam()) {
844  int iB = beamB[iSys].iPos();
845  int iBcopy = event.copy(iB, -61);
846  event[iBcopy].rotbst(Msys[iSys], false);
847  partonSystemsPtr->setInB(iSys, iBcopy);
848  beamB[iSys].iPos( iBcopy);
849  if (iSys == 0) {
850  int mother = event[iBcopy].mother1();
851  event[mother].daughter1(iBcopy);
852  }
853  }
854 
855  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
856  int iAB = partonSystemsPtr->getOut(iSys, iMem);
857  if (event[iAB].isFinal()) {
858  int iABcopy = event.copy(iAB, 62);
859  event[iABcopy].rotbst(Msys[iSys], false);
860  partonSystemsPtr->setOut(iSys, iMem, iABcopy);
861  pSumOut += event[iABcopy].p();
862  }
863  }
864 
865  }
866 
867  // Colour dipoles spanning systems gives mismatch between FSR recoils
868  // and primordial kT boosts.
869  if (allowRescatter && CORRECTMISMATCH) {
870 
871  // Find summed pT of beam remnants = - wanted pT of systems.
872  double pxBeams = 0.;
873  double pyBeams = 0.;
874  for (int iRem = nSys; iRem < beamA.size(); ++iRem) {
875  pxBeams += beamA[iRem].px();
876  pyBeams += beamA[iRem].py();
877  }
878  for (int iRem = nSys; iRem < beamB.size(); ++iRem) {
879  pxBeams += beamB[iRem].px();
880  pyBeams += beamB[iRem].py();
881  }
882 
883  // Boost all final partons in systems transversely, and also their sum.
884  Vec4 pSumTo( -pxBeams, -pyBeams, pSumOut.pz(), sqrt( pow2(pxBeams)
885  + pow2(pyBeams) + pow2(pSumOut.pz()) + pSumOut.m2Calc() ) );
886  RotBstMatrix Mmismatch;
887  Mmismatch.bst( pSumOut, pSumTo);
888  for (int iSys = 0; iSys < nSys; ++iSys)
889  for (int iMem = 0; iMem < partonSystemsPtr->sizeOut(iSys); ++iMem) {
890  int iAB = partonSystemsPtr->getOut(iSys, iMem);
891  if (event[iAB].isFinal()) event[iAB].rotbst(Mmismatch, false);
892  }
893  pSumOut.rotbst(Mmismatch);
894 
895  // Reset energy and momentum sum, to be compensated by beam remnants.
896  wPosRem = eCM - (pSumOut.e() + pSumOut.pz());
897  wNegRem = eCM - (pSumOut.e() - pSumOut.pz());
898  w2Rem = wPosRem * wNegRem;
899  if ( wPosRem < 0. || wNegRem < 0.
900  || sqrtpos(w2Rem) < sqrt(w2Beam[0]) + sqrt(w2Beam[1])) {
901  infoPtr->errorMsg("Error in BeamRemnants::setKinematics:"
902  " kinematics construction failed owing to recoil mismatch");
903  return false;
904  }
905  }
906 
907  // Construct x rescaling factors for the two remants.
908  double lambdaRoot = sqrtpos( pow2(w2Rem - w2Beam[0] - w2Beam[1])
909  - 4. * w2Beam[0] * w2Beam[1] );
910  double rescaleA = (w2Rem + w2Beam[0] - w2Beam[1] + lambdaRoot)
911  / (2. * w2Rem * xSum[0]) ;
912  double rescaleB = (w2Rem + w2Beam[1] - w2Beam[0] + lambdaRoot)
913  / (2. * w2Rem * xSum[1]) ;
914 
915  // Construct energy and pz for remnants in first beam.
916  for (int iRem = nSys; iRem < beamA.size(); ++iRem) {
917  double pPos = rescaleA * beamA[iRem].x() * wPosRem;
918  double pNeg = beamA[iRem].mT2() / pPos;
919  beamA[iRem].e( 0.5 * (pPos + pNeg) );
920  beamA[iRem].pz( 0.5 * (pPos - pNeg) );
921 
922  // Add these partons to the normal event record.
923  int iNew = event.append( beamA[iRem].id(), 63, 1 + nOffset, 0, 0, 0,
924  beamA[iRem].col(), beamA[iRem].acol(), beamA[iRem].p(),
925  beamA[iRem].m() );
926  beamA[iRem].iPos( iNew);
927  }
928 
929  // Construct energy and pz for remnants in second beam.
930  for (int iRem = nSys; iRem < beamB.size(); ++iRem) {
931  double pNeg = rescaleB * beamB[iRem].x() * wNegRem;
932  double pPos = beamB[iRem].mT2() / pNeg;
933  beamB[iRem].e( 0.5 * (pPos + pNeg) );
934  beamB[iRem].pz( 0.5 * (pPos - pNeg) );
935 
936  // Add these partons to the normal event record.
937  int iNew = event.append( beamB[iRem].id(), 63, 2 + nOffset, 0, 0, 0,
938  beamB[iRem].col(), beamB[iRem].acol(), beamB[iRem].p(),
939  beamB[iRem].m() );
940  beamB[iRem].iPos( iNew);
941  }
942 
943  // Done.
944  return true;
945 
946 }
947 
948 //--------------------------------------------------------------------------
949 
950 // Special beam remnant kinematics for gamma-gamma collision with one
951 // remnant system, other created by ISR, and for Deeply Inelastic Scattering.
952 // Currently assumes unresolved lepton.
953 
954 bool BeamRemnants::setOneRemnKinematics( Event& event) {
955 
956  // Identify beams with and without remnant.
957  int iBeamHad;
958  if ( beamAPtr->isGamma() && beamBPtr->isGamma() )
959  iBeamHad = beamAPtr->resolvedGamma() ? 1 : 2;
960  else iBeamHad = beamAPtr->isHadron() ? 1 : 2;
961  BeamParticle& beamHad = (iBeamHad == 1) ? *beamAPtr : *beamBPtr;
962  BeamParticle& beamOther = (iBeamHad == 2) ? *beamAPtr : *beamBPtr;
963 
964  // Normal beam particle locations to start with e.g. for DIS.
965  int beamOffset = 0;
966  int iBeamA = 1;
967  int iBeamB = 2;
968 
969  // Check first that event contains a reasonable number of particles
970  // as a safety measure to avoid calling container out of its range.
971  if (event.size() < 4) {
972  infoPtr->errorMsg("Error in BeamRemnants::setOneRemnKinematics:"
973  " unexpected number of particles in the event record");
974  return false;
975  }
976 
977  // Find the current beam locations from event record.
978  // Offset possible in photonproduction and (hard) diffraction.
979  for (int i = 3; i < event.size(); ++i) {
980  if ( abs(event[i].status()) > 20 ) {
981  beamOffset = i - 3;
982  iBeamA = i - 2;
983  iBeamB = i - 1;
984  break;
985  }
986  }
987 
988  // Identify remnant-side hadronic four-momentum and scattered lepton if DIS.
989  int iLepScat = isDIS ? (beamOther[0].iPos() + 2) : -1;
990 
991  // Definition of scattered lepton more intricate in the presence of QED
992  // radiation. For now, define highest energy lepton as scattered.
993  int iLepOut = particleDataPtr->isLepton(beamAPtr->id())
994  ? iBeamA + 2 + 2 : iBeamB + 2 + 2;
995  if ( !beamOther.isGamma()
996  && (iLepScat > event.size()-1
997  || !event[iLepOut].isLepton()
998  || !event[iLepScat].isLepton()
999  || !event[iLepScat].isFinal())) {
1000  double eMax = -1.0;
1001  for (int i = event.size()-1; i > 0 ; --i) {
1002  if ( event[i].isFinal()
1003  && event[iLepOut].id() == event[i].id()
1004  && event[i].e() > eMax) {
1005  iLepScat = i;
1006  eMax = event[i].e();
1007  break;
1008  }
1009  }
1010  }
1011 
1012  // In the presence of QED radiation from the beam particles, we should now
1013  // find beam particles and transform to CM frame of the final beam
1014  // configuration.
1015  if (beamAPtr->isLepton() && beamBPtr->isHadron())
1016  { iBeamA = (*beamAPtr)[0].iPos(); }
1017  if (beamBPtr->isLepton() && beamAPtr->isHadron())
1018  { iBeamB = (*beamBPtr)[0].iPos(); }
1019  double etot = (event[iBeamA].p() + event[iBeamB].p()).mCalc();
1020  RotBstMatrix MtoRest;
1021  MtoRest.bst( event[iBeamA].p() + event[iBeamB].p(), Vec4(0.,0.,0.,etot));
1022  if ( abs(event[iBeamA].pz() + event[iBeamB].pz()) > 1e-5)
1023  event.rotbst(MtoRest);
1024 
1025  Vec4 pHadScat;
1026  for (int i = 5 + beamOffset; i < event.size(); ++i)
1027  if (event[i].isFinal() && i != iLepScat) pHadScat += event[i].p();
1028 
1029  // Set scattered lepton momentum if DIS and find the hadronic system.
1030  Vec4 pLepScat = isDIS ? event[iLepScat].p() : Vec4();
1031  Vec4 pHadTot = event[iBeamA].p() + event[iBeamB].p();
1032  if ( isDIS ) pHadTot -= pLepScat;
1033  Vec4 pRemnant = pHadTot - pHadScat;
1034  double w2Tot = pHadTot.m2Calc();
1035  double w2Scat = pHadScat.m2Calc();
1036 
1037  // Find a boost to the hadronic rest frame.
1038  RotBstMatrix MtoHadRest;
1039  // Do not flip the direction with boosts.
1040  if (iBeamHad == 1) MtoHadRest.toCMframe( pRemnant, pHadScat);
1041  else MtoHadRest.toCMframe( pHadScat, pRemnant);
1042 
1043  // Boost the event to hadronic rest frame.
1044  event.rotbst( MtoHadRest);
1045  pHadScat.rotbst( MtoHadRest);
1046 
1047  // Set the variables for kT-generation.
1048  int nBeam = beamHad.size();
1049  vector<double> kTwidth(nBeam);
1050  vector<double> kTcomp(nBeam);
1051  // Currently no MPI's for one-remnant systems.
1052  int iSys = 0;
1053  double kTcompSum = 0;
1054  if (doPrimordialKT) {
1055 
1056  // Calculate the kT width from hard process.
1057  int iInA = partonSystemsPtr->getInA(iSys);
1058  int iInB = partonSystemsPtr->getInB(iSys);
1059  double sHatNow = (event[iInA].p() + event[iInB].p()).m2Calc();
1060  double mHat = sqrt(sHatNow);
1061  double mHatDamp = mHat / (mHat + halfMassForKT);
1062  double scale = infoPtr->QRen(iDS);
1063  kTwidth[0] = ( (halfScaleForKT * primordialKTsoft
1064  + scale * primordialKThard) / (halfScaleForKT + scale) ) * mHatDamp;
1065  kTcomp[0] = mHatDamp;
1066  kTcompSum = mHatDamp + nBeam - 1.;
1067 
1068  // Calculate the kT width for remnant partons.
1069  for ( int iRem = 1; iRem < nBeam; ++iRem) {
1070  kTwidth[iRem] = primordialKTremnant;
1071  kTcomp[iRem] = 1.;
1072  }
1073 
1074  }
1075 
1076  // Allow ten tries to construct kinematics (but normally works first).
1077  bool isPhysical = true;
1078  double xSum, xInvM, w2Remn, wDiff;
1079  for (int iTry = 0; iTry < NTRYKINMATCH; ++iTry) {
1080  isPhysical = true;
1081 
1082  // Add primordial kT for one-remnant system.
1083  if (doPrimordialKT) {
1084 
1085  // Generate Gaussian pT for initiator partons inside hadron.
1086  double pxSum = 0.;
1087  double pySum = 0.;
1088  for (int iPar = 0; iPar < nBeam; ++iPar) {
1089  pair<double, double> gauss2 = rndmPtr->gauss2();
1090  double px = kTwidth[iPar] * gauss2.first;
1091  double py = kTwidth[iPar] * gauss2.second;
1092  beamHad[iPar].px(px);
1093  beamHad[iPar].py(py);
1094  pxSum += px;
1095  pySum += py;
1096  }
1097 
1098  // Share recoil between all initiator partons.
1099  for (int iPar = 0; iPar < nBeam; ++iPar) {
1100  beamHad[iPar].px( beamHad[iPar].px() - pxSum*kTcomp[iPar]/kTcompSum );
1101  beamHad[iPar].py( beamHad[iPar].py() - pySum*kTcomp[iPar]/kTcompSum );
1102  }
1103  }
1104 
1105  // Pick unrescaled x values for remnants. Sum up (unscaled) p+ and p-.
1106  xSum = 0.;
1107  xInvM = 0.;
1108  for (int iRem = 1; iRem < beamHad.size(); ++iRem) {
1109  double xPrel = beamHad.xRemnant( iRem);
1110  beamHad[iRem].x(xPrel);
1111  xSum += xPrel;
1112  xInvM += beamHad[iRem].mT2() / xPrel;
1113  }
1114 
1115  // Include the primordial kT for scattered system if any.
1116  double mT2Scat = w2Scat;
1117  if (doPrimordialKT) mT2Scat += pow2( beamHad[0].pT());
1118 
1119  // Squared transverse mass for remnant, may give failure.
1120  w2Remn = xSum * xInvM;
1121  wDiff = sqrt(w2Tot) - sqrtpos(mT2Scat) - sqrtpos(w2Remn);
1122  if (wDiff < 0.) isPhysical = false;
1123  if (isPhysical) break;
1124  }
1125 
1126  if (!isPhysical) {
1127  // Boost back event.
1128  MtoHadRest.invert();
1129  event.rotbst( MtoHadRest);
1130  infoPtr->errorMsg("Error in BeamRemnants::setOneRemnKinematics:"
1131  " too big beam remnant invariant mass");
1132  return false;
1133  }
1134 
1135  // Calculate kinematics for the scattered and remnant system with kT.
1136  if (doPrimordialKT) {
1137 
1138  int iInA = partonSystemsPtr->getInA(iSys);
1139  int iInB = partonSystemsPtr->getInB(iSys);
1140  int iInit = (iBeamHad == 1) ? iInA : iInB;
1141 
1142  // Beam momentum for boosts.
1143  Vec4 pBeam = beamOther.p();
1144 
1145  // Calculate the energy and pz of initiator with added kT.
1146  double w2 = (beamOther.p() + event[iInit].p()).m2Calc();
1147  double m2 = pow2(beamHad.m());
1148  double mT2 = pow2(beamHad[iSys].pT()) + m2;
1149  double eInitNew = beamOther.e()*mT2/(w2-m2) + 0.25*(w2-m2)/beamOther.e();
1150  double pzInitNew = beamOther.e()*mT2/(w2-m2) - 0.25*(w2-m2)/beamOther.e();
1151  if(iBeamHad == 1) pzInitNew *= -1.0;
1152  beamHad[iSys].pz( pzInitNew);
1153  beamHad[iSys].e( eInitNew);
1154 
1155  // Define the boost from event rest frame to frame where initiator
1156  // has some kT.
1157  RotBstMatrix MtoInitWithkT;
1158 
1159  // Find the rest frame of beam particle and initiator.
1160  if (iBeamHad == 1) {
1161  MtoInitWithkT.toCMframe( beamHad[iSys].p(), beamOther.p());
1162  } else {
1163  MtoInitWithkT.toCMframe( beamOther.p(), beamHad[iSys].p());
1164  }
1165 
1166  // Boost to frame where beam particle has its original momentum.
1167  pBeam.rotbst(MtoInitWithkT);
1168  MtoInitWithkT.bst( pBeam, beamOther.p());
1169 
1170  // Invert the boost.
1171  MtoInitWithkT.invert();
1172 
1173  // Add the initiator with kT to the event record.
1174  int iBeam = beamHad[iSys].iPos();
1175  int iCopy = event.copy(iBeam, -61);
1176  event[iCopy].rotbst(MtoInitWithkT);
1177  if (iBeamHad == 1) partonSystemsPtr->setInA(iSys, iCopy);
1178  else partonSystemsPtr->setInB(iSys, iCopy);
1179  beamHad[iSys].iPos( iCopy);
1180  int mother = event[iCopy].mother1();
1181  event[mother].daughter1(iCopy);
1182 
1183  // Calculate transverse mass for scattered system with kT.
1184  w2Scat += pow2( beamHad[iSys].pT());
1185  }
1186 
1187  // Boost of scattered system to compensate for remnant mass.
1188  double lambda = pow2( w2Tot - w2Scat - w2Remn) - 4. * w2Scat * w2Remn;
1189  double pzNew = 0.5 * sqrt( lambda / w2Tot);
1190  // Keep the direction of momentum to ensure stable boost.
1191  if(iBeamHad == 1) pzNew *= -1.0;
1192  double eNewScat = 0.5 * (w2Tot + w2Scat - w2Remn) / sqrt(w2Tot);
1193  Vec4 pNewScat( beamHad[iSys].px(), beamHad[iSys].py(), pzNew, eNewScat);
1194  RotBstMatrix MforScat;
1195  MforScat.bst( pHadScat, pNewScat);
1196  int sizeSave = event.size();
1197  for (int i = 5 + beamOffset; i < sizeSave; ++i)
1198  if ( i != iLepScat && event[i].isFinal() ) {
1199  int iNew = event.copy( i, 62);
1200  event[iNew].rotbst( MforScat, false);
1201  }
1202 
1203  // Calculate kinematics of remnants and insert into event record.
1204  double eNewRemn = 0.5 * (w2Tot + w2Remn - w2Scat) / sqrt(w2Tot);
1205  double wNewRemn = eNewRemn + pzNew;
1206  for (int iRem = 1; iRem < beamHad.size(); ++iRem) {
1207  double wNegNow = wNewRemn * beamHad[iRem].x() / xSum;
1208  double wPosNow = beamHad[iRem].mT2() / wNegNow;
1209  beamHad[iRem].pz(-0.5 * (wNegNow - wPosNow));
1210  beamHad[iRem].e ( 0.5 * (wPosNow + wNegNow));
1211  int iNew = event.append( beamHad[iRem].id(), 63, iBeamHad + beamOffset,
1212  0, 0, 0, beamHad[iRem].col(), beamHad[iRem].acol(), beamHad[iRem].p(),
1213  beamHad[iRem].m() );
1214  beamHad[iRem].iPos( iNew);
1215  }
1216 
1217  // Boost back event.
1218  MtoHadRest.invert();
1219  event.rotbst( MtoHadRest);
1220 
1221  // Done.
1222  return true;
1223 }
1224 
1225 //--------------------------------------------------------------------------
1226 
1227 // Collapse colours and check that they are consistent.
1228 
1229 bool BeamRemnants::checkColours( Event& event) {
1230 
1231  // No colours in lepton beams so no need to do anything.
1232  if (beamAPtr->isLepton() && beamBPtr->isLepton()) return true;
1233 
1234  // Remove ambiguities when one colour collapses two ways.
1235  // Resolve chains where one colour is mapped to another.
1236  for (int iCol = 1; iCol < int(colFrom.size()); ++iCol)
1237  for (int iColRef = 0; iColRef < iCol; ++iColRef) {
1238  if (colFrom[iCol] == colFrom[iColRef]) {
1239  colFrom[iCol] = colTo[iCol];
1240  colTo[iCol] = colTo[iColRef];
1241  }
1242  if (colTo[iCol] == colFrom[iColRef]) colTo[iCol] = colTo[iColRef];
1243  }
1244 
1245  // Transform event record colours from beam remnant colour collapses.
1246  for (int i = oldSize; i < event.size(); ++i) {
1247  int col = event[i].col();
1248  int acol = event[i].acol();
1249  for (int iCol = 0; iCol < int(colFrom.size()); ++iCol) {
1250  if (col == colFrom[iCol]) {col = colTo[iCol]; event[i].col(col);}
1251  if (acol == colFrom[iCol]) {acol = colTo[iCol]; event[i].acol(acol);}
1252  // Sextets have extra, negative, tags.
1253  if (col == -colFrom[iCol]) {col = -colTo[iCol]; event[i].col(col);}
1254  if (acol == -colFrom[iCol]) {acol = -colTo[iCol]; event[i].acol(acol);}
1255  }
1256  }
1257 
1258  // Transform junction colours from beam remnant colour collapses.
1259  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun)
1260  for (int leg = 0; leg < 3; ++leg) {
1261  int col = event.colJunction(iJun, leg);
1262  for (int iCol = 0; iCol < int(colFrom.size()); ++iCol) {
1263  if (col == colFrom[iCol]) {
1264  col = colTo[iCol];
1265  event.colJunction(iJun, leg, col);
1266  }
1267  }
1268  }
1269 
1270  // Arrays for current colours and anticolours, and for singlet gluons.
1271  vector<int> colList;
1272  vector<int> acolList;
1273  vector<int> iSingletGluon;
1274 
1275  // Find current colours and anticolours in the event record.
1276  for (int i = oldSize; i < event.size(); ++i)
1277  if (event[i].isFinal()) {
1278  int id = event[i].id();
1279  int col = event[i].col();
1280  int acol = event[i].acol();
1281  int colType = event[i].colType();
1282 
1283  // Quarks must have colour set, antiquarks anticolour, gluons both.
1284  if ( (id > 0 && id < 9 && (col <= 0 || acol != 0) )
1285  || (id < 0 && id > -9 && (col != 0 || acol <= 0) )
1286  || (id == 21 && (col <= 0 || acol <= 0) ) ) {
1287  infoPtr->errorMsg("Error in BeamRemnants::checkColours: "
1288  "q/qbar/g has wrong colour slots set");
1289  return false;
1290  }
1291 
1292  // Sextets must have one positive and one negative tag
1293  if ( (colType == 3 && (col <= 0 || acol >= 0))
1294  || (colType == -3 && (col >= 0 || acol <= 0)) ) {
1295  infoPtr->errorMsg("Error in BeamRemnants::checkColours: "
1296  "sextet has wrong colours");
1297  }
1298 
1299  // Save colours/anticolours, and position of colour singlet gluons.
1300  if ( col > 0) colList.push_back( col );
1301  if (acol > 0) acolList.push_back( acol );
1302  if (col > 0 && acol == col) iSingletGluon.push_back(i);
1303  // Colour sextets
1304  if ( col < 0) acolList.push_back( -col );
1305  if (acol < 0) colList.push_back( -acol );
1306  }
1307 
1308  // Run through list of singlet gluons and put them on final-state dipole
1309  // (i,j) that offers smallest (p_g p_i) * (p_g p_j) / (p_i p_j).
1310  for (int iS = 0; iS < int(iSingletGluon.size()); ++iS) {
1311  int iGlu = iSingletGluon[iS];
1312  int iAcolDip = -1;
1313  double pT2DipMin = sCM;
1314  for (int iC = oldSize; iC < event.size(); ++iC)
1315  if (iC != iGlu && event[iC].isFinal()) {
1316  int colDip = event[iC].col();
1317  if (colDip > 0 && event[iC].acol() !=colDip)
1318  for (int iA = oldSize; iA < event.size(); ++iA)
1319  if (iA != iGlu && iA != iC && event[iA].isFinal()
1320  && event[iA].acol() == colDip && event[iA].col() !=colDip) {
1321  double pT2Dip = (event[iGlu].p() * event[iC].p())
1322  * (event[iGlu].p() * event[iA].p())
1323  / (event[iC].p() * event[iA].p());
1324  if (pT2Dip < pT2DipMin) {
1325  iAcolDip = iA;
1326  pT2DipMin = pT2Dip;
1327  }
1328  }
1329  }
1330 
1331  // Fail if no dipole. Else insert singlet gluon onto relevant dipole.
1332  if (iAcolDip == -1) return false;
1333  event[iGlu].acol( event[iAcolDip].acol() );
1334  event[iAcolDip].acol( event[iGlu].col() );
1335 
1336  // Update any junction legs that match reconnected dipole.
1337  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
1338 
1339  // Only junctions need to be updated, not antijunctions.
1340  if (event.kindJunction(iJun) % 2 == 0) continue;
1341  for (int leg = 0; leg < 3; ++leg) {
1342  int col = event.colJunction(iJun, leg);
1343  if (col == event[iGlu].acol())
1344  event.colJunction(iJun, leg, event[iGlu].col());
1345  }
1346  }
1347 
1348  }
1349 
1350  // Check that not the same colour or anticolour appears twice.
1351  for (int iCol = 0; iCol < int(colList.size()) - 1; ++iCol) {
1352  int col = colList[iCol];
1353  for (int iCol2 = iCol + 1; iCol2 < int(colList.size()); ++iCol2)
1354  if (colList[iCol2] == col) {
1355  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1356  " colour appears twice");
1357  if (!ALLOWCOLOURTWICE) return false;
1358  }
1359  }
1360  for (int iAcol = 0; iAcol < int(acolList.size()) - 1; ++iAcol) {
1361  int acol = acolList[iAcol];
1362  for (int iAcol2 = iAcol + 1; iAcol2 < int(acolList.size()); ++iAcol2)
1363  if (acolList[iAcol2] == acol) {
1364  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1365  " anticolour appears twice");
1366  if (!ALLOWCOLOURTWICE) return false;
1367  }
1368  }
1369 
1370  // Remove all matching colour-anticolour pairs.
1371  bool foundPair = true;
1372  while (foundPair && colList.size() > 0 && acolList.size() > 0) {
1373  foundPair = false;
1374  for (int iCol = 0; iCol < int(colList.size()); ++iCol) {
1375  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol) {
1376  if (acolList[iAcol] == colList[iCol]) {
1377  colList[iCol] = colList.back(); colList.pop_back();
1378  acolList[iAcol] = acolList.back(); acolList.pop_back();
1379  foundPair = true;
1380  break;
1381  }
1382  }
1383  if (foundPair) break;
1384  }
1385  }
1386 
1387  // Check that remaining (anti)colours are accounted for by junctions.
1388  for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) {
1389  int kindJun = event.kindJunction(iJun);
1390  for (int leg = 0; leg < 3; ++leg) {
1391  int colEnd = event.colJunction(iJun, leg);
1392 
1393  // Junction connected to three colours.
1394  if (kindJun == 1) {
1395  for (int iCol = 0; iCol < int(colList.size()); ++iCol)
1396  if (colList[iCol] == colEnd) {
1397  // Found colour match: remove and exit.
1398  colList[iCol] = colList.back();
1399  colList.pop_back();
1400  break;
1401  }
1402  }
1403 
1404  // Junction connected to three anticolours.
1405  else if (kindJun == 2) {
1406  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol)
1407  if (acolList[iAcol] == colEnd) {
1408  // Found colour match: remove and exit.
1409  acolList[iAcol] = acolList.back();
1410  acolList.pop_back();
1411  break;
1412  }
1413  }
1414 
1415  // Other junction types
1416  else if ( kindJun == 3 || kindJun == 5) {
1417  for (int iCol = 0; iCol < int(colList.size()); ++iCol)
1418  if (colList[iCol] == colEnd) {
1419  // Found colour match: remove and exit.
1420  colList[iCol] = colList.back();
1421  colList.pop_back();
1422  break;
1423  }
1424  }
1425 
1426  // Other antijunction types
1427  else if ( kindJun == 4 || kindJun == 6) {
1428  for (int iAcol = 0; iAcol < int(acolList.size()); ++iAcol)
1429  if (acolList[iAcol] == colEnd) {
1430  // Found colour match: remove and exit.
1431  acolList[iAcol] = acolList.back();
1432  acolList.pop_back();
1433  break;
1434  }
1435  }
1436 
1437  // End junction check.
1438  }
1439  }
1440 
1441 
1442  // Repair step - sometimes needed when rescattering allowed.
1443  if (colList.size() > 0 || acolList.size() > 0) {
1444  infoPtr->errorMsg("Warning in BeamRemnants::checkColours:"
1445  " need to repair unmatched colours");
1446  }
1447  while (colList.size() > 0 && acolList.size() > 0) {
1448 
1449  // Replace one colour and one anticolour index by a new common one.
1450  int colMatch = colList.back();
1451  int acolMatch = acolList.back();
1452  int colNew = event.nextColTag();
1453  colList.pop_back();
1454  acolList.pop_back();
1455  for (int i = oldSize; i < event.size(); ++i) {
1456  if (event[i].isFinal() && event[i].col() == colMatch) {
1457  event[i].col( colNew);
1458  break;
1459  }
1460  else if (event[i].isFinal() && event[i].acol() == -colMatch) {
1461  event[i].acol( -colNew);
1462  break;
1463  }
1464  }
1465  for (int i = oldSize; i < event.size(); ++i) {
1466  if (event[i].isFinal() && event[i].acol() == acolMatch) {
1467  event[i].acol( colNew);
1468  break;
1469  }
1470  if (event[i].isFinal() && event[i].col() == -acolMatch) {
1471  event[i].col( -colNew);
1472  break;
1473  }
1474  }
1475  }
1476 
1477  // Done.
1478  return (colList.size() == 0 && acolList.size() == 0);
1479 
1480 }
1481 
1482 //--------------------------------------------------------------------------
1483 
1484 // Update colours of outgoing particles in the event record.
1485 
1486 void BeamRemnants::updateColEvent( Event& event,
1487  vector<pair <int,int> > colChanges) {
1488 
1489  for (int iCol = 0; iCol < int(colChanges.size()); ++iCol) {
1490 
1491  int oldCol = colChanges[iCol].first;
1492  int newCol = colChanges[iCol].second;
1493  if (oldCol == newCol)
1494  continue;
1495 
1496  // Add a copy of final particles with old colour and change the colour.
1497  for (int j = 0; j < event.size(); ++j) {
1498  if (event[j].isFinal() && event[j].col() == oldCol)
1499  event[event.copy(j, 64)].col(newCol);
1500  if (event[j].isFinal() && event[j].acol() == -oldCol)
1501  event[event.copy(j, 64)].acol(-newCol);
1502 
1503  if (event[j].isFinal() && event[j].acol() == oldCol)
1504  event[event.copy(j,64)].acol(newCol);
1505  if (event[j].isFinal() && event[j].col() == -oldCol)
1506  event[event.copy(j,64)].col(-newCol);
1507  }
1508 
1509  // Update junction.
1510  for (int j = 0;j < event.sizeJunction(); ++j)
1511  for (int jCol = 0; jCol < 3; ++jCol)
1512  if (event.colJunction(j,jCol) == oldCol)
1513  event.colJunction(j,jCol,newCol);
1514  }
1515 
1516 }
1517 
1518 //==========================================================================
1519 
1520 } // end namespace Pythia8
Definition: beam.h:43
Definition: AgUStep.h:26