StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SpaceShower.cc
1 // SpaceShower.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 // SpaceShower class.
8 
9 #include "Pythia8/SpaceShower.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The SpaceShower 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 // Leftover companion can give PDF > 0 at small Q2 where other PDF's = 0,
23 // and then one can end in infinite loop of impossible kinematics.
24 const int SpaceShower::MAXLOOPTINYPDF = 10;
25 
26 // Minimal allowed c and b quark masses, for flavour thresholds.
27 const double SpaceShower::MCMIN = 1.2;
28 const double SpaceShower::MBMIN = 4.0;
29 
30 // Switch to alternative (but equivalent) backwards evolution for
31 // g -> Q Qbar (Q = c or b) when below QTHRESHOLD * mQ2.
32 const double SpaceShower::CTHRESHOLD = 2.0;
33 const double SpaceShower::BTHRESHOLD = 2.0;
34 
35 // Renew evaluation of PDF's when the pT2 step is bigger than this
36 // (in addition to initial scale and c and b thresholds.)
37 const double SpaceShower::EVALPDFSTEP = 0.1;
38 
39 // Lower limit on PDF value in order to avoid division by zero.
40 const double SpaceShower::TINYPDF = 1e-10;
41 
42 // Lower limit on estimated evolution rate, below which stop.
43 const double SpaceShower::TINYKERNELPDF = 1e-6;
44 
45 // Lower limit on pT2, below which branching is rejected.
46 const double SpaceShower::TINYPT2 = 0.25e-6;
47 
48 // No attempt to do backwards evolution of a heavy (c or b) quark
49 // if evolution starts at a scale pT2 < HEAVYPT2EVOL * mQ2.
50 const double SpaceShower::HEAVYPT2EVOL = 1.1;
51 
52 // No attempt to do backwards evolution of a heavy (c or b) quark
53 // if evolution starts at a x > HEAVYXEVOL * x_max, where
54 // x_max is the largest possible x value for a g -> Q Qbar branching.
55 const double SpaceShower::HEAVYXEVOL = 0.9;
56 
57 // When backwards evolution Q -> g + Q creates a heavy quark Q,
58 // an earlier branching g -> Q + Qbar will restrict kinematics
59 // to M_{Q Qbar}^2 > EXTRASPACEQ * 4 m_Q^2. (Smarter to be found??)
60 const double SpaceShower::EXTRASPACEQ = 2.0;
61 
62 // Never pick pT so low that alphaS is evaluated too close to Lambda_3.
63 const double SpaceShower::LAMBDA3MARGIN = 1.1;
64 
65 // Do not warn for large PDF ratios at small pT2 scales.
66 const double SpaceShower::PT2MINWARN = 1.;
67 
68 // Cutoff for f_e^e at x < 1 - 10^{-10} to be used in z selection.
69 // Note: the x_min quantity come from 1 - x_max.
70 const double SpaceShower::LEPTONXMIN = 1e-10;
71 const double SpaceShower::LEPTONXMAX = 1. - 1e-10;
72 
73 // Stop l -> l gamma evolution slightly above m2l.
74 const double SpaceShower::LEPTONPT2MIN = 1.2;
75 
76 // Enhancement of l -> l gamma trial rate to compensate imperfect modelling.
77 const double SpaceShower::LEPTONFUDGE = 10.;
78 
79 // Overestimation extra factor for t-channel weak ME corrections.
80 const double SpaceShower::WEAKPSWEIGHT = 5.;
81 
82 // Overestimation extra factors by branching type
83 const double SpaceShower::HEADROOMQ2G = 1.35;
84 const double SpaceShower::HEADROOMQ2Q = 1.15;
85 const double SpaceShower::HEADROOMG2Q = 1.35;
86 const double SpaceShower::HEADROOMG2G = 1.35;
87 
88 //--------------------------------------------------------------------------
89 
90 // Initialize alphaStrong, alphaEM and related pTmin parameters.
91 
92 void SpaceShower::init( BeamParticle* beamAPtrIn,
93  BeamParticle* beamBPtrIn) {
94 
95  // Store input pointers for future use.
96  beamAPtr = beamAPtrIn;
97  beamBPtr = beamBPtrIn;
98 
99  // Main flags to switch on and off branchings.
100  doQCDshower = settingsPtr->flag("SpaceShower:QCDshower");
101  doQEDshowerByQ = settingsPtr->flag("SpaceShower:QEDshowerByQ");
102  doQEDshowerByL = settingsPtr->flag("SpaceShower:QEDshowerByL");
103  doWeakShower = settingsPtr->flag("SpaceShower:WeakShower");
104 
105  // Matching in pT of hard interaction to shower evolution.
106  pTmaxMatch = settingsPtr->mode("SpaceShower:pTmaxMatch");
107  pTdampMatch = settingsPtr->mode("SpaceShower:pTdampMatch");
108  pTmaxFudge = settingsPtr->parm("SpaceShower:pTmaxFudge");
109  pTmaxFudgeMPI = settingsPtr->parm("SpaceShower:pTmaxFudgeMPI");
110  pTdampFudge = settingsPtr->parm("SpaceShower:pTdampFudge");
111 
112  // Optionally force emissions to be ordered in rapidity/angle.
113  doRapidityOrder = settingsPtr->flag("SpaceShower:rapidityOrder");
114 
115  // Charm, bottom and lepton mass thresholds.
116  mc = max( MCMIN, particleDataPtr->m0(4));
117  mb = max( MBMIN, particleDataPtr->m0(5));
118  m2c = pow2(mc);
119  m2b = pow2(mb);
120 
121  // Parameters of scale choices.
122  renormMultFac = settingsPtr->parm("SpaceShower:renormMultFac");
123  factorMultFac = settingsPtr->parm("SpaceShower:factorMultFac");
124  useFixedFacScale = settingsPtr->flag("SpaceShower:useFixedFacScale");
125  fixedFacScale2 = pow2(settingsPtr->parm("SpaceShower:fixedFacScale"));
126 
127  // Parameters of alphaStrong generation.
128  alphaSvalue = settingsPtr->parm("SpaceShower:alphaSvalue");
129  alphaSorder = settingsPtr->mode("SpaceShower:alphaSorder");
130  alphaSnfmax = settingsPtr->mode("StandardModel:alphaSnfmax");
131  alphaSuseCMW = settingsPtr->flag("SpaceShower:alphaSuseCMW");
132  alphaS2pi = 0.5 * alphaSvalue / M_PI;
133 
134  // Initialize alpha_strong generation.
135  alphaS.init( alphaSvalue, alphaSorder, alphaSnfmax, alphaSuseCMW);
136 
137  // Lambda for 5, 4 and 3 flavours.
138  Lambda5flav = alphaS.Lambda5();
139  Lambda4flav = alphaS.Lambda4();
140  Lambda3flav = alphaS.Lambda3();
141  Lambda5flav2 = pow2(Lambda5flav);
142  Lambda4flav2 = pow2(Lambda4flav);
143  Lambda3flav2 = pow2(Lambda3flav);
144 
145  // Regularization of QCD evolution for pT -> 0. Can be taken
146  // same as for multiparton interactions, or be set separately.
147  useSamePTasMPI = settingsPtr->flag("SpaceShower:samePTasMPI");
148  if (useSamePTasMPI) {
149  pT0Ref = settingsPtr->parm("MultipartonInteractions:pT0Ref");
150  ecmRef = settingsPtr->parm("MultipartonInteractions:ecmRef");
151  ecmPow = settingsPtr->parm("MultipartonInteractions:ecmPow");
152  pTmin = settingsPtr->parm("MultipartonInteractions:pTmin");
153  } else {
154  pT0Ref = settingsPtr->parm("SpaceShower:pT0Ref");
155  ecmRef = settingsPtr->parm("SpaceShower:ecmRef");
156  ecmPow = settingsPtr->parm("SpaceShower:ecmPow");
157  pTmin = settingsPtr->parm("SpaceShower:pTmin");
158  }
159 
160  // Calculate nominal invariant mass of events. Set current pT0 scale.
161  sCM = m2( beamAPtr->p(), beamBPtr->p());
162  eCM = sqrt(sCM);
163  pT0 = pT0Ref * pow(eCM / ecmRef, ecmPow);
164 
165  // Restrict pTmin to ensure that alpha_s(pTmin^2 + pT_0^2) does not blow up.
166  double pTminAbs = sqrtpos(pow2(LAMBDA3MARGIN) * Lambda3flav2 / renormMultFac
167  - pT0*pT0);
168  if (pTmin < pTminAbs) {
169  pTmin = pTminAbs;
170  ostringstream newPTmin;
171  newPTmin << fixed << setprecision(3) << pTmin;
172  infoPtr->errorMsg("Warning in SpaceShower::init: pTmin too low",
173  ", raised to " + newPTmin.str() );
174  infoPtr->setTooLowPTmin(true);
175  }
176 
177  // Parameters of alphaEM generation.
178  alphaEMorder = settingsPtr->mode("SpaceShower:alphaEMorder");
179 
180  // Initialize alphaEM generation.
181  alphaEM.init( alphaEMorder, settingsPtr);
182 
183  // Parameters of QED evolution.
184  pTminChgQ = settingsPtr->parm("SpaceShower:pTminchgQ");
185  pTminChgL = settingsPtr->parm("SpaceShower:pTminchgL");
186 
187  // Derived parameters of QCD evolution.
188  pT20 = pow2(pT0);
189  pT2min = pow2(pTmin);
190  pT2minChgQ = pow2(pTminChgQ);
191  pT2minChgL = pow2(pTminChgL);
192 
193  // Parameters of weak evolution.
194  weakMode = settingsPtr->mode("SpaceShower:weakShowerMode");
195  pTweakCut = settingsPtr->parm("SpaceShower:pTminWeak");
196  pT2weakCut = pow2(pTweakCut);
197  weakEnhancement = settingsPtr->parm("WeakShower:enhancement");
198  singleWeakEmission = settingsPtr->flag("WeakShower:singleEmission");
199  vetoWeakJets = settingsPtr->flag("WeakShower:vetoWeakJets");
200  vetoWeakDeltaR2 = pow2(settingsPtr->parm("weakShower:vetoWeakDeltaR"));
201 
202  // Various other parameters.
203  doMEcorrections = settingsPtr->flag("SpaceShower:MEcorrections");
204  doMEafterFirst = settingsPtr->flag("SpaceShower:MEafterFirst");
205  doPhiPolAsym = settingsPtr->flag("SpaceShower:phiPolAsym");
206  doPhiIntAsym = settingsPtr->flag("SpaceShower:phiIntAsym");
207  strengthIntAsym = settingsPtr->parm("SpaceShower:strengthIntAsym");
208  nQuarkIn = settingsPtr->mode("SpaceShower:nQuarkIn");
209 
210  // Z0 and W+- properties needed for weak showers.
211  mZ = particleDataPtr->m0(23);
212  gammaZ = particleDataPtr->mWidth(23);
213  thetaWRat = 1. / (16. * coupSMPtr->sin2thetaW()
214  * coupSMPtr->cos2thetaW());
215  mW = particleDataPtr->m0(24);
216  gammaW = particleDataPtr->mWidth(24);
217 
218  // Possibility of two predetermined hard emissions in event.
219  doSecondHard = settingsPtr->flag("SecondHard:generate");
220 
221  // Optional dampening at small pT's when large multiplicities.
222  enhanceScreening
223  = settingsPtr->mode("MultipartonInteractions:enhanceScreening");
224  if (!useSamePTasMPI) enhanceScreening = 0;
225 
226  // Possibility to allow user veto of emission step.
227  canVetoEmission = (userHooksPtr != 0)
228  ? userHooksPtr->canVetoISREmission() : false;
229 
230  // Default values for the weak shower.
231  hasWeaklyRadiated = false;
232  weakMaxWt = 1.;
233 
234 }
235 
236 //--------------------------------------------------------------------------
237 
238 // Find whether to limit maximum scale of emissions.
239 // Also allow for dampening at factorization or renormalization scale.
240 
241 bool SpaceShower::limitPTmax( Event& event, double Q2Fac, double Q2Ren) {
242 
243  // Find whether to limit pT. Begin by user-set cases.
244  bool dopTlimit = false;
245  dopTlimit1 = dopTlimit2 = false;
246  if (pTmaxMatch == 1) dopTlimit = dopTlimit1 = dopTlimit2 = true;
247  else if (pTmaxMatch == 2) dopTlimit = dopTlimit1 = dopTlimit2 = false;
248 
249  // Always restrict SoftQCD processes.
250  else if (infoPtr->isNonDiffractive() || infoPtr->isDiffractiveA()
251  || infoPtr->isDiffractiveB() || infoPtr->isDiffractiveC() )
252  dopTlimit = dopTlimit1 = dopTlimit2 = true;
253 
254  // Look if any quark (u, d, s, c, b), gluon or photon in final state.
255  else {
256  int n21 = 0;
257  for (int i = 5; i < event.size(); ++i) {
258  if (event[i].status() == -21) ++n21;
259  else if (n21 == 0) {
260  int idAbs = event[i].idAbs();
261  if (idAbs <= 5 || idAbs == 21 || idAbs == 22) dopTlimit1 = true;
262  } else if (n21 == 2) {
263  int idAbs = event[i].idAbs();
264  if (idAbs <= 5 || idAbs == 21 || idAbs == 22) dopTlimit2 = true;
265  }
266  }
267  dopTlimit = (doSecondHard) ? (dopTlimit1 && dopTlimit2) : dopTlimit1;
268  }
269 
270  // Dampening at factorization or renormalization scale; only for hardest.
271  dopTdamp = false;
272  pT2damp = 0.;
273  if ( !dopTlimit1 && (pTdampMatch == 1 || pTdampMatch == 2) ) {
274  dopTdamp = true;
275  pT2damp = pow2(pTdampFudge) * ((pTdampMatch == 1) ? Q2Fac : Q2Ren);
276  }
277 
278  // Done.
279  return dopTlimit;
280 
281 }
282 
283 //--------------------------------------------------------------------------
284 
285 // Prepare system for evolution; identify ME.
286 // Routine may be called after multiparton interactions, for a new subystem.
287 
288 void SpaceShower::prepare( int iSys, Event& event, bool limitPTmaxIn) {
289 
290  // Reset W/Z radiation flag and counters at first call for new event.
291  if (iSys == 0) {
292  nRadA.clear();
293  nRadB.clear();
294  hasWeaklyRadiated = false;
295  }
296 
297  // Find positions of incoming colliding partons.
298  int in1 = partonSystemsPtr->getInA(iSys);
299  int in2 = partonSystemsPtr->getInB(iSys);
300 
301  // Rescattered partons cannot radiate.
302  bool canRadiate1 = !(event[in1].isRescatteredIncoming());
303  bool canRadiate2 = !(event[in2].isRescatteredIncoming());
304 
305  // Reset dipole-ends list for first interaction. Also resonances.
306  if (iSys == 0) dipEnd.resize(0);
307  if (iSys == 0) idResFirst = 0;
308  if (iSys == 1) idResSecond = 0;
309 
310  // Find matrix element corrections for system.
311  int MEtype = findMEtype( iSys, event, false);
312 
313  // In case of DPS overwrite limitPTmaxIn by saved value.
314  if (doSecondHard && iSys == 0) limitPTmaxIn = dopTlimit1;
315  if (doSecondHard && iSys == 1) limitPTmaxIn = dopTlimit2;
316 
317  // Maximum pT scale for dipole ends.
318  double pTmax1 = (limitPTmaxIn) ? event[in1].scale() : eCM;
319  double pTmax2 = (limitPTmaxIn) ? event[in2].scale() : eCM;
320  if ( limitPTmaxIn && (iSys == 0 || (iSys == 1 && doSecondHard)) ) {
321  pTmax1 *= pTmaxFudge;
322  pTmax2 *= pTmaxFudge;
323  } else if (limitPTmaxIn && iSys > 0) {
324  pTmax1 *= pTmaxFudgeMPI;
325  pTmax2 *= pTmaxFudgeMPI;
326  }
327 
328  // Find dipole ends for QCD radiation.
329  // Note: colour type can change during evolution, so book also if zero.
330  if (doQCDshower) {
331  int colType1 = event[in1].colType();
332  if (canRadiate1) dipEnd.push_back( SpaceDipoleEnd( iSys, 1,
333  in1, in2, pTmax1, colType1, 0, 0, MEtype, canRadiate2) );
334  int colType2 = event[in2].colType();
335  if (canRadiate2) dipEnd.push_back( SpaceDipoleEnd( iSys, 2,
336  in2, in1, pTmax2, colType2, 0, 0, MEtype, canRadiate1) );
337  }
338 
339  // Find dipole ends for QED radiation.
340  // Note: charge type can change during evolution, so book also if zero.
341  if (doQEDshowerByQ || doQEDshowerByL) {
342  int chgType1 = ( (event[in1].isQuark() && doQEDshowerByQ)
343  || (event[in1].isLepton() && doQEDshowerByL) )
344  ? event[in1].chargeType() : 0;
345  // Special: photons have charge zero, but can evolve (only off Q for now)
346  if (event[in1].id() == 22 && doQEDshowerByQ) chgType1 = 22 ;
347  if (canRadiate1) dipEnd.push_back( SpaceDipoleEnd( iSys, -1,
348  in1, in2, pTmax1, 0, chgType1, 0, MEtype, canRadiate2) );
349  int chgType2 = ( (event[in2].isQuark() && doQEDshowerByQ)
350  || (event[in2].isLepton() && doQEDshowerByL) )
351  ? event[in2].chargeType() : 0;
352  // Special: photons have charge zero, but can evolve (only off Q for now)
353  if (event[in2].id() == 22 && doQEDshowerByQ) chgType2 = 22 ;
354  if (canRadiate2) dipEnd.push_back( SpaceDipoleEnd( iSys, -2,
355  in2, in1, pTmax2, 0, chgType2, 0, MEtype, canRadiate1) );
356  }
357 
358  // Find dipole ends for weak radiation. No right-handed W emission.
359  // Currently leptons are not allow to emit W bosons and only
360  // emissions from the hard process are included.
361  if (doWeakShower && iSys == 0) {
362 
363  // Determine what type of 2 -> 2 process it is.
364  int MEtypeWeak = findMEtype( iSys, event, true);
365  if (MEtypeWeak == 201 || MEtypeWeak == 202 || MEtypeWeak == 203 ||
366  MEtypeWeak == 206 || MEtypeWeak == 207 || MEtypeWeak == 208) {
367 
368  // Nonidentical incoming flavours.
369  if (event[in1].id() != event[in2].id()) {
370  if (event[in1].id() == event[in1 + 2].id()) tChannel = true;
371  else if (event[in2].id() == event[in1 + 2].id()) tChannel = false;
372  // No quark matches the outgoing, choose randomly.
373  else tChannel = (rndmPtr->flat() < 0.5);
374  // In case of same quark flavours, choose randomly.
375  } else tChannel = (rndmPtr->flat() < 0.5);
376  }
377 
378  // Set up weak dipole ends for first incoming parton.
379  int weakPol = (rndmPtr->flat() > 0.5) ? -1 : 1;
380  if (event[in1].idAbs() < 20) event[in1].pol(weakPol);
381  if (canRadiate1) {
382  if ( (weakMode == 0 || weakMode == 1) && weakPol == -1
383  && event[in1].isQuark() )
384  dipEnd.push_back( SpaceDipoleEnd( iSys, 1, in1, in2, pTmax1,
385  0, 0, 1, MEtypeWeak, canRadiate2, weakPol) );
386  if ( (weakMode == 0 || weakMode == 2)
387  && (event[in1].isQuark() || event[in1].isLepton()) )
388  dipEnd.push_back( SpaceDipoleEnd( iSys, 1, in1, in2, pTmax1,
389  0, 0, 2, MEtypeWeak + 5, canRadiate2, weakPol) );
390  }
391 
392  // Set up weak dipole ends for second incoming parton.
393  if (event[in1].id() != - event[in2].id())
394  weakPol = (rndmPtr->flat() > 0.5) ? -1 : 1;
395  if (event[in2].idAbs() < 20) event[in2].pol(weakPol);
396  if (canRadiate2) {
397  if ( (weakMode == 0 || weakMode == 1) && weakPol == -1
398  && event[in2].isQuark())
399  dipEnd.push_back( SpaceDipoleEnd( iSys, 2, in2, in1, pTmax2,
400  0, 0, 1, MEtypeWeak, canRadiate1, weakPol) );
401  if ( (weakMode == 0 || weakMode == 2) &&
402  (event[in2].isQuark() || event[in2].isLepton()) )
403  dipEnd.push_back( SpaceDipoleEnd( iSys, 2, in2, in1, pTmax2,
404  0, 0, 2, MEtypeWeak + 5, canRadiate1, weakPol) );
405  }
406  }
407 
408  // Store the z and pT2 values of the last previous splitting
409  // when an event history has already been constructed.
410  if (iSys == 0 && infoPtr->hasHistory()) {
411  double zNow = infoPtr->zNowISR();
412  double pT2Now = infoPtr->pT2NowISR();
413  for (int iDipEnd = 0; iDipEnd < int(dipEnd.size()); ++iDipEnd) {
414  dipEnd[iDipEnd].zOld = zNow;
415  dipEnd[iDipEnd].pT2Old = pT2Now;
416  ++dipEnd[iDipEnd].nBranch;
417  }
418  }
419 
420 }
421 
422 //--------------------------------------------------------------------------
423 
424 // Select next pT in downwards evolution of the existing dipoles.
425 
426 double SpaceShower::pTnext( Event& event, double pTbegAll, double pTendAll,
427  int nRadIn) {
428 
429  // Current cm energy, in case it varies between events.
430  sCM = m2( beamAPtr->p(), beamBPtr->p());
431  eCM = sqrt(sCM);
432  pTbegRef = pTbegAll;
433 
434  // Starting values: no radiating dipole found.
435  nRad = nRadIn;
436  double pT2sel = pow2(pTendAll);
437  iDipSel = 0;
438  iSysSel = 0;
439  dipEndSel = 0;
440 
441  // Loop over all possible dipole ends.
442  for (int iDipEnd = 0; iDipEnd < int(dipEnd.size()); ++iDipEnd) {
443  iDipNow = iDipEnd;
444  dipEndNow = &dipEnd[iDipEnd];
445  iSysNow = dipEndNow->system;
446  dipEndNow->pT2 = 0.;
447  double pTbegDip = min( pTbegAll, dipEndNow->pTmax );
448 
449  // Check whether dipole end should be allowed to shower.
450  double pT2begDip = pow2(pTbegDip);
451  if (pT2begDip > pT2sel && ( dipEndNow->colType != 0
452  || dipEndNow->chgType != 0 || dipEndNow->weakType != 0) ) {
453  double pT2endDip = 0.;
454 
455  // Determine lower cut for evolution, for QCD or QED (q or l).
456  if (dipEndNow->colType != 0)
457  pT2endDip = max( pT2sel, pT2min );
458  else if (abs(dipEndNow->weakType) != 0)
459  pT2endDip = max( pT2sel, pT2weakCut);
460  else if (abs(dipEndNow->chgType) != 3 && dipEndNow->chgType != 0)
461  pT2endDip = max( pT2sel, pT2minChgQ );
462  else
463  pT2endDip = max( pT2sel, pT2minChgL );
464 
465  // Find properties of dipole and radiating dipole end.
466  sideA = ( abs(dipEndNow->side) == 1 );
467  BeamParticle& beamNow = (sideA) ? *beamAPtr : *beamBPtr;
468  BeamParticle& beamRec = (sideA) ? *beamBPtr : *beamAPtr;
469  iNow = beamNow[iSysNow].iPos();
470  iRec = beamRec[iSysNow].iPos();
471  idDaughter = beamNow[iSysNow].id();
472  xDaughter = beamNow[iSysNow].x();
473  x1Now = (sideA) ? xDaughter : beamRec[iSysNow].x();
474  x2Now = (sideA) ? beamRec[iSysNow].x() : xDaughter;
475 
476  // Note dipole mass correction when recoiler is a rescatter.
477  m2Rec = (dipEndNow->normalRecoil) ? 0. : event[iRec].m2();
478  m2Dip = x1Now * x2Now * sCM + m2Rec;
479 
480  // Now do evolution in pT2, for QCD, QED or weak.
481  if (pT2begDip > pT2endDip) {
482  if (dipEndNow->colType != 0) pT2nextQCD( pT2begDip, pT2endDip);
483  else if (dipEndNow->chgType != 0 || idDaughter == 22)
484  pT2nextQED( pT2begDip, pT2endDip);
485  else if (dipEndNow->weakType != 0) pT2nextWeak( pT2begDip, pT2endDip);
486 
487  // Update if found larger pT than current maximum.
488  if (dipEndNow->pT2 > pT2sel) {
489  pT2sel = dipEndNow->pT2;
490  iDipSel = iDipNow;
491  iSysSel = iSysNow;
492  dipEndSel = dipEndNow;
493  }
494  }
495  }
496  // End loop over dipole ends.
497  }
498 
499  // Return nonvanishing value if found pT is bigger than already found.
500  return (dipEndSel == 0) ? 0. : sqrt(pT2sel);
501 }
502 
503 //--------------------------------------------------------------------------
504 
505 // Evolve a QCD dipole end.
506 
507 void SpaceShower::pT2nextQCD( double pT2begDip, double pT2endDip) {
508 
509  // Some properties and kinematical starting values.
510  BeamParticle& beam = (sideA) ? *beamAPtr : *beamBPtr;
511  bool isGluon = (idDaughter == 21);
512  bool isValence = beam[iSysNow].isValence();
513  int MEtype = dipEndNow->MEtype;
514  double pT2 = pT2begDip;
515  double xMaxAbs = beam.xMax(iSysNow);
516  double zMinAbs = xDaughter / xMaxAbs;
517  if (xMaxAbs < 0.) {
518  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQCD: "
519  "xMaxAbs negative");
520  return;
521  }
522 
523  // Starting values for handling of massive quarks (c/b), if any.
524  double idMassive = 0;
525  if ( abs(idDaughter) == 4 ) idMassive = 4;
526  if ( abs(idDaughter) == 5 ) idMassive = 5;
527  bool isMassive = (idMassive > 0);
528  double m2Massive = 0.;
529  double mRatio = 0.;
530  double zMaxMassive = 1.;
531  double m2Threshold = pT2;
532 
533  // Evolution below scale of massive quark or at large x is impossible.
534  if (isMassive) {
535  m2Massive = (idMassive == 4) ? m2c : m2b;
536  if (pT2 < HEAVYPT2EVOL * m2Massive) return;
537  mRatio = sqrt( m2Massive / m2Dip );
538  zMaxMassive = (1. - mRatio) / ( 1. + mRatio * (1. - mRatio) );
539  if (xDaughter > HEAVYXEVOL * zMaxMassive * xMaxAbs) return;
540 
541  // Find threshold scale below which only g -> Q + Qbar will be allowed.
542  m2Threshold = (idMassive == 4) ? min( pT2, CTHRESHOLD * m2c)
543  : min( pT2, BTHRESHOLD * m2b);
544  }
545 
546  // Variables used inside evolution loop. (Mainly dummy starting values.)
547  int nFlavour = 3;
548  double b0 = 4.5;
549  double Lambda2 = Lambda3flav2;
550  double pT2minNow = pT2endDip;
551  int idMother = 0;
552  int idSister = 0;
553  double z = 0.;
554  double zMaxAbs = 0.;
555  double zRootMax = 0.;
556  double zRootMin = 0.;
557  double g2gInt = 0.;
558  double q2gInt = 0.;
559  double q2qInt = 0.;
560  double g2qInt = 0.;
561  double g2Qenhance = 0.;
562  double xPDFdaughter = 0.;
563  double xPDFmother[21] = {0.};
564  double xPDFgMother = 0.;
565  double xPDFmotherSum = 0.;
566  double kernelPDF = 0.;
567  double xMother = 0.;
568  double wt = 0.;
569  double Q2 = 0.;
570  double mSister = 0.;
571  double m2Sister = 0.;
572  double pT2corr = 0.;
573  double pT2PDF = pT2;
574  bool needNewPDF = true;
575 
576  // Begin evolution loop towards smaller pT values.
577  int loopTinyPDFdau = 0;
578  bool hasTinyPDFdau = false;
579  do {
580  wt = 0.;
581 
582  // Bad sign if repeated looping with small daughter PDF, so fail.
583  // (Example: if all PDF's = 0 below Q_0, except for c/b companion.)
584  if (hasTinyPDFdau) ++loopTinyPDFdau;
585  if (loopTinyPDFdau > MAXLOOPTINYPDF) {
586  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQCD: "
587  "small daughter PDF");
588  return;
589  }
590 
591  // Initialize integrals of splitting kernels and evaluate parton
592  // densities at the beginning. Reinitialize after long evolution
593  // in pT2 or when crossing c and b flavour thresholds.
594  if (needNewPDF || pT2 < EVALPDFSTEP * pT2PDF) {
595  pT2PDF = pT2;
596  hasTinyPDFdau = false;
597 
598  // Determine overestimated z range; switch at c and b masses.
599  if (pT2 > m2b) {
600  nFlavour = 5;
601  pT2minNow = max( m2b, pT2endDip);
602  b0 = 23./6.;
603  Lambda2 = Lambda5flav2;
604  } else if (pT2 > m2c) {
605  nFlavour = 4;
606  pT2minNow = max( m2c, pT2endDip);
607  b0 = 25./6.;
608  Lambda2 = Lambda4flav2;
609  } else {
610  nFlavour = 3;
611  pT2minNow = pT2endDip;
612  b0 = 27./6.;
613  Lambda2 = Lambda3flav2;
614  }
615  // A change of renormalization scale expressed by a change of Lambda.
616  Lambda2 /= renormMultFac;
617  zMaxAbs = 1. - 0.5 * (pT2minNow / m2Dip) *
618  ( sqrt( 1. + 4. * m2Dip / pT2minNow ) - 1. );
619  if (isMassive) zMaxAbs = min( zMaxAbs, zMaxMassive);
620 
621  // Go to another z range with lower mass scale if current is closed.
622  if (zMinAbs > zMaxAbs) {
623  if (nFlavour == 3 || (idMassive == 4 && nFlavour == 4)
624  || idMassive == 5) return;
625  pT2 = (nFlavour == 4) ? m2c : m2b;
626  continue;
627  }
628 
629  // Parton density of daughter at current scale.
630  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
631  xPDFdaughter = beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2);
632  if (xPDFdaughter < TINYPDF) {
633  xPDFdaughter = TINYPDF;
634  hasTinyPDFdau = true;
635  }
636 
637  // Integrals of splitting kernels for gluons: g -> g, q -> g.
638  if (isGluon) {
639  g2gInt = HEADROOMG2G * 6.
640  * log(zMaxAbs * (1.-zMinAbs) / (zMinAbs * (1.-zMaxAbs)));
641  if (doMEcorrections) g2gInt *= calcMEmax(MEtype, 21, 21);
642  q2gInt = HEADROOMQ2G * (16./3.)
643  * (1./sqrt(zMinAbs) - 1./sqrt(zMaxAbs));
644  if (doMEcorrections) q2gInt *= calcMEmax(MEtype, 1, 21);
645 
646  // Parton density of potential quark mothers to a g.
647  xPDFmotherSum = 0.;
648  for (int i = -nQuarkIn; i <= nQuarkIn; ++i) {
649  if (i == 0) {
650  xPDFmother[10] = 0.;
651  } else {
652  xPDFmother[i+10] = beam.xfISR(iSysNow, i, xDaughter, pdfScale2);
653  xPDFmotherSum += xPDFmother[i+10];
654  }
655  }
656 
657  // Total QCD evolution coefficient for a gluon.
658  kernelPDF = g2gInt + q2gInt * xPDFmotherSum / xPDFdaughter;
659 
660  // For valence quark only need consider q -> q g branchings.
661  // Introduce an extra factor sqrt(z) to smooth bumps.
662  } else if (isValence) {
663  zRootMin = (1. + sqrt(zMinAbs)) / (1. - sqrt(zMinAbs));
664  zRootMax = (1. + sqrt(zMaxAbs)) / (1. - sqrt(zMaxAbs));
665  q2qInt = (8./3.) * log( zRootMax / zRootMin );
666  if (doMEcorrections) q2qInt *= calcMEmax(MEtype, 1, 1);
667  kernelPDF = q2qInt;
668 
669  // Integrals of splitting kernels for quarks: q -> q, g -> q.
670  } else {
671  q2qInt = HEADROOMQ2Q * (8./3.)
672  * log( (1. - zMinAbs) / (1. - zMaxAbs) );
673  if (doMEcorrections) q2qInt *= calcMEmax(MEtype, 1, 1);
674  g2qInt = HEADROOMG2Q * 0.5 * (zMaxAbs - zMinAbs);
675  if (doMEcorrections) g2qInt *= calcMEmax(MEtype, 21, 1);
676 
677  // Increase estimated upper weight for g -> Q + Qbar.
678  if (isMassive) {
679  if (alphaSorder == 0) g2Qenhance = log(pT2/m2Massive)
680  / log(m2Threshold/m2Massive);
681  else {
682  double m2log = log( m2Massive / Lambda2);
683  g2Qenhance = log( log(pT2/Lambda2) / m2log )
684  / log( log(m2Threshold/Lambda2) / m2log );
685  }
686  g2qInt *= g2Qenhance;
687  }
688 
689  // Parton density of a potential gluon mother to a q.
690  xPDFgMother = beam.xfISR(iSysNow, 21, xDaughter, pdfScale2);
691 
692  // Total QCD evolution coefficient for a quark.
693  kernelPDF = q2qInt + g2qInt * xPDFgMother / xPDFdaughter;
694  }
695 
696  // End evaluation of splitting kernels and parton densities.
697  needNewPDF = false;
698  }
699  if (kernelPDF < TINYKERNELPDF) return;
700 
701  // Pick pT2 (in overestimated z range), for one of three different cases.
702  // Assume form alphas(pT0^2 + pT^2) * dpT^2/(pT0^2 + pT^2).
703  double Q2alphaS;
704 
705  // Fixed alpha_strong.
706  if (alphaSorder == 0) {
707  pT2 = (pT2 + pT20) * pow( rndmPtr->flat(),
708  1. / (alphaS2pi * kernelPDF)) - pT20;
709 
710  // First-order alpha_strong.
711  } else if (alphaSorder == 1) {
712  pT2 = Lambda2 * pow( (pT2 + pT20) / Lambda2,
713  pow(rndmPtr->flat(), b0 / kernelPDF) ) - pT20;
714 
715  // For second order reject by second term in alpha_strong expression.
716  } else {
717  do {
718  pT2 = Lambda2 * pow( (pT2 + pT20) / Lambda2,
719  pow(rndmPtr->flat(), b0 / kernelPDF) ) - pT20;
720  Q2alphaS = renormMultFac * max( pT2 + pT20,
721  pow2(LAMBDA3MARGIN) * Lambda3flav2);
722  } while (alphaS.alphaS2OrdCorr(Q2alphaS) < rndmPtr->flat()
723  && pT2 > pT2minNow);
724  }
725 
726  // Check for pT2 values that prompt special action.
727 
728  // If fallen into b threshold region, force g -> b + bbar.
729  if (idMassive == 5 && pT2 < m2Threshold) {
730  pT2nearQCDthreshold( beam, m2Massive, m2Threshold, xMaxAbs,
731  zMinAbs, zMaxMassive );
732  return;
733 
734  // If crossed b threshold, continue evolution from this threshold.
735  } else if (nFlavour == 5 && pT2 < m2b) {
736  needNewPDF = true;
737  pT2 = m2b;
738  continue;
739 
740  // If fallen into c threshold region, force g -> c + cbar.
741  } else if (idMassive == 4 && pT2 < m2Threshold) {
742  pT2nearQCDthreshold( beam, m2Massive, m2Threshold, xMaxAbs,
743  zMinAbs, zMaxMassive );
744  return;
745 
746  // If crossed c threshold, continue evolution from this threshold.
747  } else if (nFlavour == 4 && pT2 < m2c) {
748  needNewPDF = true;
749  pT2 = m2c;
750  continue;
751 
752  // Abort evolution if below cutoff scale, or below another branching.
753  } else if (pT2 < pT2endDip) return;
754 
755  // Select z value of branching to g, and corrective weight.
756  if (isGluon) {
757  // g -> g (+ g).
758  if (rndmPtr->flat() * kernelPDF < g2gInt) {
759  idMother = 21;
760  idSister = 21;
761  z = 1. / ( 1. + ((1. - zMinAbs) / zMinAbs) * pow( (zMinAbs *
762  (1. - zMaxAbs)) / (zMaxAbs * (1. - zMinAbs)), rndmPtr->flat() ) );
763  wt = pow2( 1. - z * (1. - z));
764  // Account for headroom factor used to enhance trial probability
765  wt /= HEADROOMG2G;
766  } else {
767  // q -> g (+ q): also select flavour.
768  double temp = xPDFmotherSum * rndmPtr->flat();
769  idMother = -nQuarkIn - 1;
770  do { temp -= xPDFmother[(++idMother) + 10]; }
771  while (temp > 0. && idMother < nQuarkIn);
772  idSister = idMother;
773  z = (zMinAbs * zMaxAbs) / pow2( sqrt(zMinAbs) + rndmPtr->flat()
774  * ( sqrt(zMaxAbs)- sqrt(zMinAbs) ));
775  wt = 0.5 * (1. + pow2(1. - z)) * sqrt(z)
776  * xPDFdaughter / xPDFmother[idMother + 10];
777  // Account for headroom factor used to enhance trial probability
778  wt /= HEADROOMQ2G;
779  }
780 
781  // Select z value of branching to q, and corrective weight.
782  // Include massive kernel corrections for c and b quarks.
783  } else {
784  // q -> q (+ g).
785  if (isValence || rndmPtr->flat() * kernelPDF < q2qInt) {
786  idMother = idDaughter;
787  idSister = 21;
788  // Valence more peaked at large z.
789  if (isValence) {
790  double zTmp = zRootMin * pow(zRootMax / zRootMin, rndmPtr->flat() );
791  z = pow2( (1. - zTmp) / (1. + zTmp) );
792  } else {
793  z = 1. - (1. - zMinAbs) * pow( (1. - zMaxAbs) / (1. - zMinAbs),
794  rndmPtr->flat() );
795  }
796  if (!isMassive) {
797  wt = 0.5 * (1. + pow2(z));
798  } else {
799  //?? Bug?? should be 2 more for massive part??
800  // wt = 0.5 * (1. + pow2(z) - z * pow2(1.-z) * m2Massive / pT2);
801  wt = 0.5 * (1. + pow2(z)) - z * pow2(1.-z) * m2Massive / pT2;
802  }
803  if (isValence) wt *= sqrt(z);
804  // Account for headroom factor for sea quarks
805  else wt /= HEADROOMQ2Q;
806  // g -> q (+ qbar).
807  } else {
808  idMother = 21;
809  idSister = - idDaughter;
810  z = zMinAbs + rndmPtr->flat() * (zMaxAbs - zMinAbs);
811  if (!isMassive) {
812  wt = (pow2(z) + pow2(1.-z)) * xPDFdaughter / xPDFgMother;
813  } else {
814  wt = (pow2(z) + pow2(1.-z) + 2. * z * (1.-z) * m2Massive / pT2)
815  * xPDFdaughter / (xPDFgMother * g2Qenhance) ;
816  }
817  // Account for headroom factor for gluons
818  wt /= HEADROOMG2Q;
819  }
820  }
821 
822  // Derive Q2 and x of mother from pT2 and z.
823  Q2 = pT2 / (1. - z);
824  xMother = xDaughter / z;
825  // Correction to x for massive recoiler from rescattering.
826  if (!dipEndNow->normalRecoil) {
827  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
828  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
829  }
830  if (xMother > xMaxAbs) { wt = 0.; continue; }
831 
832  // Forbidden emission if outside allowed z range for given pT2.
833  mSister = particleDataPtr->m0(idSister);
834  m2Sister = pow2(mSister);
835  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
836  if (pT2corr < TINYPT2) { wt = 0.; continue; }
837 
838  // Optionally veto emissions not ordered in rapidity (= angle).
839  if ( doRapidityOrder && dipEndNow->nBranch > 0
840  && pT2 > pow2( (1. - z) / (z * (1. - dipEndNow->zOld)) )
841  * dipEndNow->pT2Old ) { wt = 0.; continue; }
842 
843  // If creating heavy quark by Q -> g + Q then next need g -> Q + Qbar.
844  // So minimum total mass2 is 4 * m2Sister, but use more to be safe.
845  if ( isGluon && ( abs(idMother) == 4 || abs(idMother) == 5 )) {
846  double m2QQsister = EXTRASPACEQ * 4. * m2Sister;
847  double pT2QQcorr = Q2 - z * (m2Dip + Q2) * (Q2 + m2QQsister) / m2Dip;
848  if (pT2QQcorr < TINYPT2) { wt = 0.; continue; }
849  }
850 
851  // Evaluation of ME correction.
852  if (doMEcorrections) wt *= calcMEcorr(MEtype, idMother, idDaughter,
853  m2Dip, z, Q2, m2Sister) / calcMEmax(MEtype, idMother, idDaughter);
854 
855  // Optional dampening of large pT values in first radiation.
856  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
857  wt *= pT2damp / (pT2 + pT2damp);
858 
859  // Idea suggested by Gosta Gustafson: increased screening in events
860  // with large activity can be simulated by pT0_eff = sqrt(n) * pT0.
861  if (enhanceScreening == 2) {
862  int nSysNow = infoPtr->nMPI() + infoPtr->nISR() + 1;
863  double WTscreen = pow2( (pT2 + pT20) / (pT2 + nSysNow * pT20) );
864  wt *= WTscreen;
865  }
866 
867  // Evaluation of new daughter and mother PDF's.
868  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
869  double xPDFdaughterNew = max ( TINYPDF,
870  beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2) );
871  double xPDFmotherNew =
872  beam.xfISR(iSysNow, idMother, xMother, pdfScale2);
873  wt *= xPDFmotherNew / xPDFdaughterNew;
874 
875  // Check that valence step does not cause problem.
876  if (wt > 1. && pT2 > PT2MINWARN) infoPtr->errorMsg("Warning in "
877  "SpaceShower::pT2nextQCD: weight above unity");
878 
879  // Iterate until acceptable pT (or have fallen below pTmin).
880  } while (wt < rndmPtr->flat()) ;
881 
882  // Save values for (so far) acceptable branching.
883  dipEndNow->store( idDaughter,idMother, idSister, x1Now, x2Now, m2Dip,
884  pT2, z, xMother, Q2, mSister, m2Sister, pT2corr);
885 
886 }
887 
888 //--------------------------------------------------------------------------
889 
890 // Evolve a QCD dipole end near threshold, with g -> Q + Qbar enforced.
891 // Note: No explicit Sudakov factor formalism here. Instead use that
892 // df_Q(x, pT2) = (alpha_s/2pi) * (dT2/pT2) * ((gluon) * (splitting)).
893 // This implies that effects of Q -> Q + g are neglected in this range.
894 
895 void SpaceShower::pT2nearQCDthreshold( BeamParticle& beam,
896  double m2Massive, double m2Threshold, double xMaxAbs,
897  double zMinAbs, double zMaxMassive) {
898 
899  // Initial values, to be used in kinematics and weighting.
900  double Lambda2 = (abs(idDaughter) == 4) ? Lambda4flav2 : Lambda5flav2;
901  Lambda2 /= renormMultFac;
902  double logM2Lambda2 = (alphaSorder > 0) ? log( m2Massive / Lambda2 ) : 1.;
903  pdfScale2 = (useFixedFacScale) ? fixedFacScale2
904  : factorMultFac * m2Threshold;
905  double xPDFmotherOld = beam.xfISR(iSysNow, 21, xDaughter, pdfScale2);
906 
907  // Variables used inside evolution loop. (Mainly dummy start values.)
908  int loop = 0;
909  double wt = 0.;
910  double pT2 = 0.;
911  double z = 0.;
912  double Q2 = 0.;
913  double pT2corr = 0.;
914  double xMother = 0.;
915 
916  // Begin loop over tries to find acceptable g -> Q + Qbar branching.
917  do {
918  wt = 0.;
919 
920  // Check that not caught in infinite loop with impossible kinematics.
921  if (++loop > 100) {
922  infoPtr->errorMsg("Error in SpaceShower::pT2nearQCDthreshold: "
923  "stuck in loop");
924  return;
925  }
926 
927  // Pick dpT2/pT2 in range [m2Massive,thresholdRatio * m2Massive].
928  pT2 = m2Massive * pow( m2Threshold / m2Massive, rndmPtr->flat() );
929 
930  // Pick z flat in allowed range.
931  z = zMinAbs + rndmPtr->flat() * (zMaxMassive - zMinAbs);
932 
933  // Check that kinematically possible choice.
934  Q2 = pT2 / (1.-z) - m2Massive;
935  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Massive) / m2Dip;
936  if (pT2corr < TINYPT2) continue;
937 
938  // Correction factor for running alpha_s. (Only first order for now.)
939  wt = (alphaSorder > 0) ? logM2Lambda2 / log( pT2 / Lambda2 ) : 1.;
940 
941  // Correction factor for splitting kernel.
942  wt *= pow2(z) + pow2(1.-z) + 2. * z * (1.-z) * m2Massive / pT2;
943 
944  // x, including correction for massive recoiler from rescattering.
945  xMother = xDaughter / z;
946  if (!dipEndNow->normalRecoil) {
947  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
948  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
949  }
950  if (xMother > xMaxAbs) { wt = 0.; continue; }
951 
952  // Correction factor for gluon density.
953  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
954  double xPDFmotherNew = beam.xfISR(iSysNow, 21, xMother, pdfScale2);
955  wt *= xPDFmotherNew / xPDFmotherOld;
956 
957  // Iterate until acceptable pT and z.
958  } while (wt < rndmPtr->flat()) ;
959 
960  // Save values for (so far) acceptable branching.
961  double mSister = (abs(idDaughter) == 4) ? mc : mb;
962  dipEndNow->store( idDaughter, 21, -idDaughter, x1Now, x2Now, m2Dip,
963  pT2, z, xMother, Q2, mSister, pow2(mSister), pT2corr);
964 
965 }
966 
967 //--------------------------------------------------------------------------
968 
969 // Evolve a QED dipole end.
970 
971 void SpaceShower::pT2nextQED( double pT2begDip, double pT2endDip) {
972 
973  // Type of dipole and starting values.
974  BeamParticle& beam = (sideA) ? *beamAPtr : *beamBPtr;
975  bool isLeptonBeam = beam.isLepton();
976  int MEtype = dipEndNow->MEtype;
977  bool isPhoton = (idDaughter == 22);
978  double pT2 = pT2begDip;
979  double m2Lepton = (isLeptonBeam) ? pow2(beam.m()) : 0.;
980  if (isLeptonBeam && pT2begDip < m2Lepton) return;
981 
982  // Currently no f -> gamma branching implemented for lepton beams
983  if (isPhoton && isLeptonBeam) return;
984 
985  // alpha_em at maximum scale provides upper estimate.
986  double alphaEMmax = alphaEM.alphaEM(renormMultFac * pT2begDip);
987  double alphaEM2pi = alphaEMmax / (2. * M_PI);
988 
989  // Maximum x of mother implies minimum z = xDaughter / xMother.
990  double xMaxAbs = (isLeptonBeam) ? LEPTONXMAX : beam.xMax(iSysNow);
991  double zMinAbs = xDaughter / xMaxAbs;
992  if (xMaxAbs < 0.) {
993  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQED: "
994  "xMaxAbs negative");
995  return;
996  }
997 
998  // Maximum z from minimum pT and, for lepton, from minimum x_gamma.
999  double zMaxAbs = 1. - 0.5 * (pT2endDip / m2Dip) *
1000  ( sqrt( 1. + 4. * m2Dip / pT2endDip ) - 1. );
1001  if (isLeptonBeam) {
1002  double zMaxLepton = xDaughter / (xDaughter + LEPTONXMIN);
1003  if (zMaxLepton < zMaxAbs) zMaxAbs = zMaxLepton;
1004  }
1005  if (zMaxAbs < zMinAbs) return;
1006 
1007  // Variables used inside evolution loop. (Mainly dummy start values.)
1008  int idMother = 0;
1009  int idSister = 22;
1010  double z = 0.;
1011  double xMother = 0.;
1012  double wt = 0.;
1013  double Q2 = 0.;
1014  double mSister = 0.;
1015  double m2Sister = 0.;
1016  double pT2corr = 0.;
1017 
1018  // QED evolution of fermions
1019  if (!isPhoton) {
1020 
1021  // Integrals of splitting kernels for fermions: f -> f. Use 1 + z^2 < 2.
1022  // Ansatz f(z) = 2 / (1 - z), with + 2 / (z - xDaughter) for lepton.
1023  double f2fInt = 0.;
1024  double f2fIntA = 2. * log( (1. - zMinAbs) / (1. - zMaxAbs) );
1025  double f2fIntB = 0.;
1026  if (isLeptonBeam) {
1027  f2fIntB = 2. * log( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter) );
1028  f2fInt = f2fIntA + f2fIntB;
1029  } else f2fInt = pow2(dipEndNow->chgType / 3.) * f2fIntA;
1030 
1031  // Upper estimate for evolution equation, including fudge factor.
1032  if (doMEcorrections) f2fInt *= calcMEmax(MEtype, 1, 1);
1033  double kernelPDF = alphaEM2pi * f2fInt;
1034  double fudge = (isLeptonBeam) ? LEPTONFUDGE * log(m2Dip/m2Lepton) : 1.;
1035  kernelPDF *= fudge;
1036  if (kernelPDF < TINYKERNELPDF) return;
1037 
1038  // Begin evolution loop towards smaller pT values.
1039  do {
1040 
1041  // Pick pT2 (in overestimated z range).
1042  // For l -> l gamma include extrafactor 1 / ln(pT2 / m2l) in evolution.
1043  double shift = pow(rndmPtr->flat(), 1. / kernelPDF);
1044  if (isLeptonBeam) pT2 = m2Lepton * pow( pT2 / m2Lepton, shift);
1045  else pT2 = pT2 * shift;
1046 
1047  // Abort evolution if below cutoff scale, or below another branching.
1048  if (pT2 < pT2endDip) return;
1049  if (isLeptonBeam && pT2 < LEPTONPT2MIN * m2Lepton) return;
1050 
1051  // Select z value of branching f -> f + gamma, and corrective weight.
1052  idMother = idDaughter;
1053  wt = 1.;
1054  if (isLeptonBeam) {
1055  if (f2fIntA > rndmPtr->flat() * (f2fIntA + f2fIntB)) {
1056  z = 1. - (1. - zMinAbs)
1057  * pow( (1. - zMaxAbs) / (1. - zMinAbs), rndmPtr->flat() );
1058  } else {
1059  z = xDaughter + (zMinAbs - xDaughter)
1060  * pow( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter),
1061  rndmPtr->flat() );
1062  }
1063  wt *= (z - xDaughter) / (1. - xDaughter);
1064  } else {
1065  z = 1. - (1. - zMinAbs)
1066  * pow( (1. - zMaxAbs) / (1. - zMinAbs), rndmPtr->flat() );
1067  }
1068  wt *= 0.5 * (1. + pow2(z));
1069 
1070  // Derive Q2 and x of mother from pT2 and z.
1071  Q2 = pT2 / (1. - z);
1072  xMother = xDaughter / z;
1073  // Correction to x for massive recoiler from rescattering.
1074  if (!dipEndNow->normalRecoil) {
1075  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
1076  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
1077  }
1078  if (xMother > xMaxAbs) { wt = 0.; continue; }
1079 
1080  // Forbidden emission if outside allowed z range for given pT2.
1081  mSister = 0.;
1082  m2Sister = 0.;
1083  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
1084  if (pT2corr < TINYPT2) { wt = 0.; continue; }
1085 
1086  // Correct by ln(pT2 / m2l) and fudge factor.
1087  if (isLeptonBeam) wt *= log(pT2 / m2Lepton) / fudge;
1088 
1089  // Evaluation of ME correction.
1090  if (doMEcorrections) wt *= calcMEcorr(MEtype, idMother, idDaughter,
1091  m2Dip, z, Q2, m2Sister) / calcMEmax(MEtype, idMother, idDaughter);
1092 
1093  // Extra QED correction for f fbar -> W+- gamma. Debug??
1094  if (doMEcorrections && MEtype == 1 && idDaughter == idMother
1095  && ( (iSysNow == 0 && idResFirst == 24)
1096  || (iSysNow == 1 && idResSecond == 24) ) ) {
1097  double tHe = -Q2;
1098  double uHe = Q2 - m2Dip * (1. - z) / z;
1099  double chg1 = abs(dipEndNow->chgType / 3.);
1100  double chg2 = 1. - chg1;
1101  wt *= pow2(chg1 * uHe - chg2 * tHe)
1102  / ( (tHe + uHe) * (pow2(chg1) * uHe + pow2(chg2) * tHe) );
1103  }
1104 
1105  // Optional dampening of large pT values in first radiation.
1106  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
1107  wt *= pT2damp / (pT2 + pT2damp);
1108 
1109  // Correct to current value of alpha_EM.
1110  double alphaEMnow = alphaEM.alphaEM(renormMultFac * pT2);
1111  wt *= (alphaEMnow / alphaEMmax);
1112 
1113  // Evaluation of new daughter and mother PDF's.
1114  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
1115  double xPDFdaughterNew = max ( TINYPDF,
1116  beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2) );
1117  double xPDFmotherNew =
1118  beam.xfISR(iSysNow, idMother, xMother, pdfScale2);
1119  wt *= xPDFmotherNew / xPDFdaughterNew;
1120 
1121  // Iterate until acceptable pT (or have fallen below pTmin).
1122  } while (wt < rndmPtr->flat()) ;
1123  }
1124 
1125  // QED evolution of photons (so far only for hadron beams).
1126  else {
1127 
1128  // Initial values
1129  int nFlavour = 3;
1130  double kernelPDF = 0.0;
1131  double xPDFdaughter = 0.;
1132  double xPDFmother[21] = {0.};
1133  double xPDFmotherSum = 0.0;
1134  double pT2PDF = pT2;
1135  double pT2minNow = pT2endDip;
1136  bool needNewPDF = true;
1137 
1138  // Begin evolution loop towards smaller pT values.
1139  int loopTinyPDFdau = 0;
1140  bool hasTinyPDFdau = false;
1141  do {
1142  wt = 0.;
1143 
1144  // Bad sign if repeated looping with small daughter PDF, so fail.
1145  if (hasTinyPDFdau) ++loopTinyPDFdau;
1146  if (loopTinyPDFdau > MAXLOOPTINYPDF) {
1147  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQED: "
1148  "small daughter PDF");
1149  return;
1150  }
1151 
1152  // Initialize integrals of splitting kernels and evaluate parton
1153  // densities at the beginning. Reinitialize after long evolution
1154  // in pT2 or when crossing c and b flavour thresholds.
1155  if (needNewPDF || pT2 < EVALPDFSTEP * pT2PDF) {
1156 
1157  pT2PDF = pT2;
1158  hasTinyPDFdau = false;
1159 
1160  // Determine overestimated z range; switch at c and b masses.
1161  if (pT2 > m2b && nQuarkIn >= 5) {
1162  nFlavour = 5;
1163  pT2minNow = max( m2b, pT2endDip);
1164  } else if (pT2 > m2c && nQuarkIn >= 4) {
1165  nFlavour = 4;
1166  pT2minNow = max( m2c, pT2endDip);
1167  } else {
1168  nFlavour = 3;
1169  pT2minNow = pT2endDip;
1170  }
1171 
1172  // Compute upper z limit
1173  zMaxAbs = 1. - 0.5 * (pT2minNow / m2Dip) *
1174  ( sqrt( 1. + 4. * m2Dip / pT2minNow ) - 1. );
1175 
1176  // Parton density of daughter at current scale.
1177  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
1178  xPDFdaughter = beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2);
1179  if (xPDFdaughter < TINYPDF) {
1180  xPDFdaughter = TINYPDF;
1181  hasTinyPDFdau = true;
1182  }
1183 
1184  // Integral over f -> gamma f splitting kernel.
1185  // Normalized so: 4/3 aS/2pi P(z) -> eq^2 * aEM/2pi P(z).
1186  // (Charge-weighting happens below.)
1187  double q2gInt = 4. * (1./sqrt(zMinAbs) - 1./sqrt(zMaxAbs));
1188 
1189  // Charge-weighted Parton density of potential quark mothers.
1190  xPDFmotherSum = 0.;
1191  for (int i = -nFlavour; i <= nFlavour; ++i) {
1192  if (i == 0) {
1193  xPDFmother[10] = 0.;
1194  } else {
1195  xPDFmother[i+10] = pow2((abs(i+1) % 2 + 1)/3.0)
1196  * beam.xfISR(iSysNow, i, xDaughter, pdfScale2);
1197  xPDFmotherSum += xPDFmother[i+10];
1198  }
1199  }
1200 
1201  // Total QED evolution coefficient for a photon.
1202  kernelPDF = q2gInt * xPDFmotherSum / xPDFdaughter;
1203 
1204  // End evaluation of splitting kernels and parton densities.
1205  needNewPDF = false;
1206  }
1207  if (kernelPDF < TINYKERNELPDF) return;
1208 
1209  // Select pT2 for next trial branching
1210  pT2 *= pow( rndmPtr->flat(), 1. / (alphaEM2pi * kernelPDF));
1211 
1212  // If crossed b threshold, continue evolution from this threshold.
1213  if (nFlavour == 5 && pT2 < m2b) {
1214  needNewPDF = true;
1215  pT2 = m2b;
1216  continue;
1217  }
1218 
1219  // If crossed c threshold, continue evolution from this threshold.
1220  else if (nFlavour == 4 && pT2 < m2c) {
1221  needNewPDF = true;
1222  pT2 = m2c;
1223  continue;
1224  }
1225 
1226  // Abort evolution if below cutoff scale, or below another branching.
1227  else if (pT2 < pT2endDip) return;
1228 
1229  // Select flavour for trial branching
1230  double temp = xPDFmotherSum * rndmPtr->flat();
1231  idMother = -nQuarkIn - 1;
1232  do {
1233  temp -= xPDFmother[(++idMother) + 10];
1234  } while (temp > 0. && idMother < nQuarkIn);
1235 
1236  // Sister is same as mother, but can have m2 > 0
1237  idSister = idMother;
1238  mSister = particleDataPtr->m0(idSister);
1239  m2Sister = pow2(mSister);
1240 
1241  // Select z for trial branching
1242  z = (zMinAbs * zMaxAbs) / pow2( sqrt(zMinAbs) + rndmPtr->flat()
1243  * ( sqrt(zMaxAbs)- sqrt(zMinAbs) ));
1244 
1245  // Trial weight: splitting kernel
1246  wt = 0.5 * (1. + pow2(1. - z)) * sqrt(z);
1247 
1248  // Trial weight: running alpha_EM
1249  double alphaEMnow = alphaEM.alphaEM(renormMultFac * pT2);
1250  wt *= (alphaEMnow / alphaEMmax);
1251 
1252  // Derive Q2 and x of mother from pT2 and z
1253  Q2 = pT2 / (1. - z);
1254  xMother = xDaughter / z;
1255  // Correction to x for massive recoiler from rescattering.
1256  if (!dipEndNow->normalRecoil) {
1257  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
1258  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
1259  }
1260 
1261  // Compute pT2corr
1262  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
1263  if (pT2corr < TINYPT2) { wt = 0.; continue; }
1264 
1265  // If creating heavy quark by Q -> gamma + Q then next g -> Q + Qbar.
1266  // So minimum total mass2 is 4 * m2Sister, but use more to be safe.
1267  if ( abs(idMother) == 4 || abs(idMother) == 5 ) {
1268  double m2QQsister = EXTRASPACEQ * 4. * m2Sister;
1269  double pT2QQcorr = Q2 - z * (m2Dip + Q2) * (Q2 + m2QQsister) / m2Dip;
1270  if (pT2QQcorr < TINYPT2) { wt = 0.; continue; }
1271  }
1272 
1273  // Optional dampening of large pT values in first radiation.
1274  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
1275  wt *= pT2damp / (pT2 + pT2damp);
1276 
1277  // Evaluation of new daughter PDF
1278  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
1279  double xPDFdaughterNew = beam.xfISR(iSysNow, idDaughter, xDaughter,
1280  pdfScale2);
1281  if (xPDFdaughterNew < TINYPDF) {
1282  xPDFdaughterNew = TINYPDF;
1283  }
1284 
1285  // Evaluation of new charge-weighted mother PDF
1286  double xPDFmotherNew = pow2( (abs(idMother+1) % 2 + 1)/3.0 )
1287  * beam.xfISR(iSysNow, idMother, xMother, pdfScale2);
1288 
1289  // Trial weight: divide out old pdf ratio
1290  wt *= xPDFdaughter / xPDFmother[idMother + 10];
1291 
1292  // Trial weight: new pdf ratio
1293  wt *= xPDFmotherNew / xPDFdaughterNew;
1294 
1295  // Iterate until acceptable pT (or have fallen below pTmin).
1296  } while (wt < rndmPtr->flat()) ;
1297  }
1298 
1299  // Save values for (so far) acceptable branching.
1300  dipEndNow->store( idDaughter, idMother, idSister, x1Now, x2Now, m2Dip,
1301  pT2, z, xMother, Q2, mSister, m2Sister, pT2corr);
1302 
1303 }
1304 
1305 //--------------------------------------------------------------------------
1306 
1307 void SpaceShower::pT2nextWeak( double pT2begDip, double pT2endDip) {
1308 
1309  // Type of dipole and starting values.
1310  BeamParticle& beam = (sideA) ? *beamAPtr : *beamBPtr;
1311  bool isLeptonBeam = beam.isLepton();
1312  bool isValence = beam[iSysNow].isValence();
1313  int MEtype = dipEndNow->MEtype;
1314  double pT2 = pT2begDip;
1315  double m2Lepton = (isLeptonBeam) ? pow2(beam.m()) : 0.;
1316  if (isLeptonBeam && pT2begDip < m2Lepton) return;
1317 
1318  // alpha_em at maximum scale provides upper estimate.
1319  double alphaEMmax = alphaEM.alphaEM(pT2begDip);
1320  double alphaEM2pi = alphaEMmax / (2. * M_PI);
1321 
1322  // Maximum x of mother implies minimum z = xDaughter / xMother.
1323  double xMaxAbs = (isLeptonBeam) ? LEPTONXMAX : beam.xMax(iSysNow);
1324  double zMinAbs = xDaughter / xMaxAbs;
1325  if (xMaxAbs < 0.) {
1326  infoPtr->errorMsg("Warning in SpaceShower::pT2nextWeak: "
1327  "xMaxAbs negative");
1328  return;
1329  }
1330 
1331  // Maximum z from minimum pT and, for lepton, from minimum x_gamma.
1332  double zMaxAbs = 1. - 0.5 * (pT2endDip / m2Dip) *
1333  ( sqrt( 1. + 4. * m2Dip / pT2endDip ) - 1. );
1334  if (isLeptonBeam) {
1335  double zMaxLepton = xDaughter / (xDaughter + LEPTONXMIN);
1336  if (zMaxLepton < zMaxAbs) zMaxAbs = zMaxLepton;
1337  }
1338  if (zMaxAbs < zMinAbs) return;
1339 
1340  // Variables used inside evolution loop. (Mainly dummy start values.)
1341  int idMother = 0;
1342  int idSister = 0;
1343  // Check whether emission of W+, W- or Z0.
1344  if (dipEndNow->weakType == 1) {
1345  idSister = (idDaughter > 0) ? -24 : 24;
1346  if (abs(idDaughter) % 2 == 1) idSister = -idSister;
1347  } else if (dipEndNow->weakType == 2) idSister = 23;
1348  double z = 0.;
1349  double xMother = 0.;
1350  double wt = 0.;
1351  double Q2 = 0.;
1352  double mSister = particleDataPtr->mSel(idSister);
1353  double m2Sister = pow2(mSister);
1354  double pT2corr = 0.;
1355 
1356  // Find maximum z due to massive sister.
1357  // Still need to prove that this actually is an overestimate.
1358  double mRatio = mSister/ sqrt(m2Dip);
1359  double m2R1 = 1. + pow2(mRatio);
1360  double zMaxMassive = 1. / (m2R1 + pT2endDip/m2Dip);
1361  if (zMaxMassive < zMaxAbs) zMaxAbs = zMaxMassive;
1362  if (zMaxAbs < zMinAbs) return;
1363 
1364  // Weak evolution of fermions.
1365  // Integrals of splitting kernels for fermions: f -> f.
1366  // Old ansatz f(z) = 2 / (1 - z), with + 2 / (z - xDaughter) for lepton.
1367  // New Ansatz f(z) = (1 + (1+r^2)^2) / (1 - z * (1 + r^2)).
1368  // This should always be an overestimate for massive emissions.
1369  // Not yet implemented correctly for lepton beam.
1370  double f2fInt = 0.;
1371  double f2fIntA = (1. + pow2(zMaxAbs * m2R1)) / m2R1
1372  * log( (1. - m2R1 * zMinAbs) / (1. - m2R1 * zMaxAbs) );
1373  double f2fIntB = 0.;
1374  double zRootMin = (1. + sqrt(m2R1 * zMinAbs)) / (1. - sqrt(m2R1 * zMinAbs));
1375  double zRootMax = (1. + sqrt(m2R1 * zMaxAbs)) / (1. - sqrt(m2R1 * zMaxAbs));
1376  if (isLeptonBeam) {
1377  f2fIntB = 2. * log( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter) );
1378  f2fInt = f2fIntA + f2fIntB;
1379  } else if (isValence)
1380  f2fInt = (1. + pow2(zMaxAbs) * pow2(m2R1))/ sqrt(m2R1)
1381  * log( zRootMax / zRootMin );
1382  else
1383  f2fInt = f2fIntA;
1384 
1385  // Calculate the weak coupling: separate for left- and right-handed fermions.
1386  double weakCoupling = 0;
1387  if (dipEndNow->weakType == 1)
1388  weakCoupling = 2. * alphaEM2pi / (4. * coupSMPtr->sin2thetaW());
1389  else if (dipEndNow->weakType == 2 && dipEndNow->weakPol == -1)
1390  weakCoupling = alphaEM2pi * thetaWRat
1391  * pow2(2. * coupSMPtr->lf( abs(idDaughter) ));
1392  else
1393  weakCoupling = alphaEM2pi * thetaWRat
1394  * pow2(2. * coupSMPtr->rf(abs(idDaughter) ));
1395 
1396  // Find the possible mothers for a W emission. So far quarks only.
1397  vector<int> possibleMothers;
1398  if (abs(idDaughter) % 2 == 0) {
1399  possibleMothers.push_back(1);
1400  possibleMothers.push_back(3);
1401  possibleMothers.push_back(5);
1402  } else {
1403  possibleMothers.push_back(2);
1404  possibleMothers.push_back(4);
1405  }
1406  if (idDaughter < 0)
1407  for (unsigned int i = 0;i < possibleMothers.size();i++)
1408  possibleMothers[i] = - possibleMothers[i];
1409 
1410  // Check if daughter estimate is 0, return in that case.
1411  // Only write warning if u, d or g is the daughter.
1412  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2begDip;
1413  double xPDFdaughter = beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2);
1414  if (xPDFdaughter < TINYPDF) {
1415  if (abs(idDaughter) == 1 || abs(idDaughter) == 2 || abs(idDaughter) == 21)
1416  infoPtr->errorMsg("Warning in SpaceShower::pT2nextWeak: "
1417  "very small PDF");
1418  return;
1419  }
1420 
1421  // PDF and CKM upper estimate needed for W emission.
1422  double overEstimatePDFCKM = 0.;
1423  if (dipEndNow->weakType == 1) {
1424  for (unsigned int i = 0; i < possibleMothers.size(); i++)
1425  overEstimatePDFCKM += coupSMPtr->V2CKMid(idDaughter, possibleMothers[i])
1426  * beam.xfISR(iSysNow, possibleMothers[i], xDaughter, pdfScale2)
1427  / xPDFdaughter;
1428  }
1429  if (dipEndNow->weakType == 2) overEstimatePDFCKM = 1.;
1430 
1431  // Upper estimate for evolution equation, including fudge factor.
1432  if (doMEcorrections) f2fInt *= calcMEmax(MEtype, 1, 1);
1433  double kernelPDF = weakCoupling * f2fInt * weakEnhancement;
1434 
1435  // PDF and CKM overestimate.
1436  kernelPDF *= overEstimatePDFCKM;
1437  double fudge = (isLeptonBeam) ? LEPTONFUDGE * log(m2Dip/m2Lepton) : 1.;
1438  kernelPDF *= fudge;
1439  if (kernelPDF < TINYKERNELPDF) return;
1440 
1441  // Begin evolution loop towards smaller pT values.
1442  do {
1443 
1444  // Pick pT2 (in overestimated z range).
1445  // For l -> l gamma include extrafactor 1 / ln(pT2 / m2l) in evolution.
1446  double shift = pow(rndmPtr->flat(), 1. / kernelPDF);
1447  if (isLeptonBeam) pT2 = m2Lepton * pow( pT2 / m2Lepton, shift);
1448  else pT2 = pT2 * shift;
1449 
1450  // Abort evolution if below cutoff scale, or below another branching.
1451  if (pT2 < pT2endDip) return;
1452  if (isLeptonBeam && pT2 < LEPTONPT2MIN * m2Lepton) return;
1453 
1454  // Abort evolution if below mass treshold.
1455  if (pT2 < HEAVYPT2EVOL * pow2(particleDataPtr->m0(idDaughter))) return;
1456 
1457  // Set the id for the mother particle to be equal to the daughter
1458  // particle. This is correct for Z, and it will later be changed for W.
1459  idMother = idDaughter;
1460 
1461  // Select z value of branching f -> f + Z/W, and corrective weight.
1462  wt = 1.0;
1463  if (isLeptonBeam) {
1464  if (f2fIntA > rndmPtr->flat() * (f2fIntA + f2fIntB)) {
1465  z = 1. - (1. - zMinAbs)
1466  * pow( (1. - zMaxAbs) / (1. - zMinAbs), rndmPtr->flat() );
1467  } else {
1468  z = xDaughter + (zMinAbs - xDaughter)
1469  * pow( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter),
1470  rndmPtr->flat() );
1471  }
1472  wt *= (z - xDaughter) / (1. - xDaughter);
1473  } else if (isValence) {
1474  // Valence more peaked at large z.
1475  double zTmp = zRootMin * pow(zRootMax / zRootMin, rndmPtr->flat() );
1476  z = pow2( (1. - zTmp) / (1. + zTmp) ) / m2R1;
1477  wt *= sqrt(z);
1478  } else {
1479  z = (1. - (1. - zMinAbs * m2R1) * pow( (1. - zMaxAbs * m2R1)
1480  / (1. - zMinAbs * m2R1), rndmPtr->flat() ) ) / m2R1;
1481  }
1482  wt *= (1. + pow2(z * m2R1)) / (1. + pow2(zMaxAbs * m2R1));
1483 
1484  // Derive Q2 and x of mother from pT2 and z.
1485  Q2 = pT2 / (1. - z);
1486  xMother = xDaughter / z;
1487 
1488  // Correction to x for massive recoiler from rescattering.
1489  if (!dipEndNow->normalRecoil) {
1490  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
1491  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
1492  }
1493  if (xMother > xMaxAbs) { wt = 0.; continue; }
1494 
1495  // Forbidden emission if outside allowed z range for given pT2.
1496  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
1497  if (pT2corr < TINYPT2) { wt = 0.; continue; }
1498 
1499  // Correct by ln(pT2 / m2l) and fudge factor.
1500  if (isLeptonBeam) wt *= log(pT2 / m2Lepton) / fudge;
1501 
1502  // Evaluation of ME correction.
1503  if (doMEcorrections) wt *= calcMEcorr(MEtype, idMother, idDaughter,
1504  m2Dip, z, Q2, m2Sister) / calcMEmax(MEtype, idMother, idDaughter);
1505 
1506  // Optional dampening of large pT values in first radiation.
1507  // Allow damping also for corrected matrix elements
1508  if (dopTdamp && iSysNow == 0 && nRad == 0)
1509  wt *= pT2damp / (pT2 + pT2damp);
1510 
1511  // Correct to current value of alpha_EM.
1512  double alphaEMnow = alphaEM.alphaEM(pT2);
1513  wt *= (alphaEMnow / alphaEMmax);
1514 
1515  // Evaluation of new daughter and mother PDF's for Z.
1516  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
1517  double xPDFdaughterNew = max ( TINYPDF,
1518  beam.xfISR(iSysNow, idDaughter, xDaughter, pdfScale2) );
1519  if (dipEndNow->weakType == 2) {
1520  double xPDFmotherNew
1521  = beam.xfISR(iSysNow, idMother, xMother, pdfScale2);
1522  wt *= xPDFmotherNew / xPDFdaughterNew;
1523  }
1524 
1525  // Evaluation of daughter and mother PDF's for W.
1526  if (dipEndNow->weakType == 1) {
1527  double valPDFCKM[3] = {0.};
1528  double sumPDFCKM = 0.;
1529  for (unsigned int i = 0; i < possibleMothers.size(); i++) {
1530  valPDFCKM[i] = coupSMPtr->V2CKMid(idDaughter,possibleMothers[i])
1531  * beam.xfISR(iSysNow, possibleMothers[i], xMother, pdfScale2)
1532  / xPDFdaughterNew;
1533  sumPDFCKM += valPDFCKM[i];
1534  }
1535  wt *= sumPDFCKM / overEstimatePDFCKM;
1536 
1537  // Choose id for mother in case of W.
1538  double pickId = sumPDFCKM * rndmPtr->flat();
1539  int iId = -1;
1540  int pmSize = possibleMothers.size();
1541  do pickId -= valPDFCKM[++iId];
1542  while (pickId > 0. && iId < pmSize);
1543  idMother = possibleMothers[iId];
1544  }
1545 
1546  // Warn if too big weight.
1547  if (wt > 1.) infoPtr->errorMsg("Warning in SpaceShower::pT2nextWeak: "
1548  "weight is above unity.");
1549 
1550  // Iterate until acceptable pT (or have fallen below pTmin).
1551  } while (wt < rndmPtr->flat()) ;
1552 
1553  // Save values for (so far) acceptable branching.
1554  dipEndNow->store( idDaughter, idMother, idSister, x1Now, x2Now, m2Dip,
1555  pT2, z, xMother, Q2, mSister, m2Sister, pT2corr);
1556 
1557 }
1558 
1559 //--------------------------------------------------------------------------
1560 
1561 // Kinematics of branching.
1562 // Construct mother -> daughter + sister, with recoiler on other side.
1563 
1564 bool SpaceShower::branch( Event& event) {
1565 
1566  // Side on which branching occured.
1567  int side = abs(dipEndSel->side);
1568  double sideSign = (side == 1) ? 1. : -1.;
1569 
1570  // Read in flavour and colour variables.
1571  int iDaughter = partonSystemsPtr->getInA(iSysSel);
1572  int iRecoiler = partonSystemsPtr->getInB(iSysSel);
1573  if (side == 2) swap(iDaughter, iRecoiler);
1574  int idDaughterNow = dipEndSel->idDaughter;
1575  int idMother = dipEndSel->idMother;
1576  int idSister = dipEndSel->idSister;
1577  int idRecoiler = event[iRecoiler].id();
1578  int colDaughter = event[iDaughter].col();
1579  int acolDaughter = event[iDaughter].acol();
1580 
1581  // Recoil parton may be rescatterer, requiring special processing.
1582  bool normalRecoil = dipEndSel->normalRecoil;
1583  int iRecoilMother = event[iRecoiler].mother1();
1584 
1585  // Read in kinematical variables.
1586  double x1 = dipEndSel->x1;
1587  double x2 = dipEndSel->x2;
1588  double xMo = dipEndSel->xMo;
1589  double m2 = dipEndSel->m2Dip;
1590  double m = sqrt(m2);
1591  double pT2 = dipEndSel->pT2;
1592  double z = dipEndSel->z;
1593  double Q2 = dipEndSel->Q2;
1594  double mSister = dipEndSel->mSister;
1595  double m2Sister = dipEndSel->m2Sister;
1596  double pT2corr = dipEndSel->pT2corr;
1597  double x1New = (side == 1) ? xMo : x1;
1598  double x2New = (side == 2) ? xMo : x2;
1599 
1600  // Read in MEtype:
1601  int MEtype = dipEndSel->MEtype;
1602 
1603  // Rescatter: kinematics may fail; use the rescatterFail flag to tell
1604  // parton level to try again.
1605  rescatterFail = false;
1606 
1607  // Construct kinematics of mother, sister and recoiler in old rest frame.
1608  // Normally both mother and recoiler are taken massless.
1609  double eNewRec, pzNewRec, pTbranch, pzMother;
1610  if (normalRecoil) {
1611  eNewRec = 0.5 * (m2 + Q2) / m;
1612  pzNewRec = -sideSign * eNewRec;
1613  pTbranch = sqrt(pT2corr) * m2 / ( z * (m2 + Q2) );
1614  pzMother = sideSign * 0.5 * m * ( (m2 - Q2) / ( z * (m2 + Q2) )
1615  + (Q2 + m2Sister) / m2 );
1616  // More complicated kinematics when recoiler not massless. May fail.
1617  } else {
1618  m2Rec = event[iRecoiler].m2();
1619  double s1Tmp = m2 + Q2 - m2Rec;
1620  double s3Tmp = m2 / z - m2Rec;
1621  double r1Tmp = sqrt(s1Tmp * s1Tmp + 4. * Q2 * m2Rec);
1622  eNewRec = 0.5 * (m2 + m2Rec + Q2) / m;
1623  pzNewRec = -sideSign * 0.5 * r1Tmp / m;
1624  double pT2br = Q2 * s3Tmp * (m2 / z - m2 - Q2)
1625  - m2Sister * s1Tmp * s3Tmp - m2Rec * pow2(Q2 + m2Sister);
1626  if (pT2br <= 0.) return false;
1627  pTbranch = sqrt(pT2br) / r1Tmp;
1628  pzMother = sideSign * (m * s3Tmp
1629  - eNewRec * (m2 / z - Q2 - m2Rec - m2Sister)) / r1Tmp;
1630  }
1631  // Common final kinematics steps for both normal and rescattering.
1632  double eMother = sqrt( pow2(pTbranch) + pow2(pzMother) );
1633  double pzSister = pzMother + pzNewRec;
1634  double eSister = sqrt( pow2(pTbranch) + pow2(pzSister) + m2Sister );
1635  Vec4 pMother( pTbranch, 0., pzMother, eMother );
1636  Vec4 pSister( pTbranch, 0., pzSister, eSister );
1637  Vec4 pNewRec( 0., 0., pzNewRec, eNewRec );
1638 
1639  // Current event and subsystem size.
1640  int eventSizeOld = event.size();
1641  int systemSizeOld = partonSystemsPtr->sizeAll(iSysSel);
1642 
1643  // Save properties to be restored in case of user-hook veto of emission.
1644  int beamOff1 = 1 + beamOffset;
1645  int beamOff2 = 2 + beamOffset;
1646  int ev1Dau1V = event[beamOff1].daughter1();
1647  int ev2Dau1V = event[beamOff2].daughter1();
1648  vector<int> statusV, mother1V, mother2V, daughter1V, daughter2V;
1649 
1650  // Check if the first emission shoild be checked for removal
1651  bool canMergeFirst = (mergingHooksPtr != 0)
1652  ? mergingHooksPtr->canVetoEmission() : false;
1653  if (canVetoEmission || canMergeFirst || doWeakShower) {
1654  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1655  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1656  statusV.push_back( event[iOldCopy].status());
1657  mother1V.push_back( event[iOldCopy].mother1());
1658  mother2V.push_back( event[iOldCopy].mother2());
1659  daughter1V.push_back( event[iOldCopy].daughter1());
1660  daughter2V.push_back( event[iOldCopy].daughter2());
1661  }
1662  }
1663 
1664  // Take copy of existing system, to be given modified kinematics.
1665  // Incoming negative status. Rescattered also negative, but after copy.
1666  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1667  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1668  int statusOld = event[iOldCopy].status();
1669  int statusNew = (iOldCopy == iDaughter
1670  || iOldCopy == iRecoiler) ? statusOld : 44;
1671  int iNewCopy = event.copy(iOldCopy, statusNew);
1672  if (statusOld < 0) event[iNewCopy].statusNeg();
1673  }
1674 
1675  // Define colour flow in branching.
1676  // Default corresponds to f -> f + gamma.
1677  int colMother = colDaughter;
1678  int acolMother = acolDaughter;
1679  int colSister = 0;
1680  int acolSister = 0;
1681  if (idSister == 22 || idSister == 23 || abs(idSister) == 24) ;
1682  // q -> q + g and 50% of g -> g + g; need new colour.
1683  else if (idSister == 21 && ( (idMother > 0 && idMother < 9)
1684  || (idMother == 21 && rndmPtr->flat() < 0.5) ) ) {
1685  colMother = event.nextColTag();
1686  colSister = colMother;
1687  acolSister = colDaughter;
1688  // qbar -> qbar + g and other 50% of g -> g + g; need new colour.
1689  } else if (idSister == 21) {
1690  acolMother = event.nextColTag();
1691  acolSister = acolMother;
1692  colSister = acolDaughter;
1693  // q -> g + q.
1694  } else if (idDaughterNow == 21 && idMother > 0) {
1695  colMother = colDaughter;
1696  acolMother = 0;
1697  colSister = acolDaughter;
1698  // qbar -> g + qbar
1699  } else if (idDaughterNow == 21) {
1700  acolMother = acolDaughter;
1701  colMother = 0;
1702  acolSister = colDaughter;
1703  // g -> q + qbar.
1704  } else if (idDaughterNow > 0 && idDaughterNow < 9) {
1705  acolMother = event.nextColTag();
1706  acolSister = acolMother;
1707  // g -> qbar + q.
1708  } else if (idDaughterNow < 0 && idDaughterNow > -9) {
1709  colMother = event.nextColTag();
1710  colSister = colMother;
1711  // q -> gamma + q.
1712  } else if (idDaughterNow == 22 && idMother > 0) {
1713  colMother = event.nextColTag();
1714  colSister = colMother;
1715  // qbar -> gamma + qbar.
1716  } else if (idDaughterNow == 22) {
1717  acolMother = event.nextColTag();
1718  acolSister = acolMother;
1719  }
1720 
1721  // Indices of partons involved. Add new sister.
1722  int iMother = eventSizeOld + side - 1;
1723  int iNewRecoiler = eventSizeOld + 2 - side;
1724  int iSister = event.append( idSister, 43, iMother, 0, 0, 0,
1725  colSister, acolSister, pSister, mSister, sqrt(pT2) );
1726 
1727  // References to the partons involved.
1728  Particle& daughter = event[iDaughter];
1729  Particle& mother = event[iMother];
1730  Particle& newRecoiler = event[iNewRecoiler];
1731  Particle& sister = event.back();
1732 
1733  // Replace old by new mother; update new recoiler.
1734  mother.id( idMother );
1735  mother.status( -41);
1736  mother.cols( colMother, acolMother);
1737  mother.p( pMother);
1738  if (mother.idAbs() == 21 || mother.idAbs() == 22) mother.pol(9);
1739  newRecoiler.status( (normalRecoil) ? -42 : -46 );
1740  newRecoiler.p( pNewRec);
1741  if (!normalRecoil) newRecoiler.m( event[iRecoiler].m() );
1742 
1743  // Update mother and daughter pointers; also for beams.
1744  daughter.mothers( iMother, 0);
1745  mother.daughters( iSister, iDaughter);
1746  if (iSysSel == 0) {
1747  event[beamOff1].daughter1( (side == 1) ? iMother : iNewRecoiler );
1748  event[beamOff2].daughter1( (side == 2) ? iMother : iNewRecoiler );
1749  }
1750 
1751  // Special checks to set weak particles status equal to 47.
1752  if (sister.idAbs() == 23 || sister.idAbs() == 24) sister.status(47);
1753 
1754  // Find boost to old rest frame.
1755  RotBstMatrix Mtot;
1756  if (normalRecoil) Mtot.bst(0., 0., (x2 - x1) / (x1 + x2) );
1757  else if (side == 1)
1758  Mtot.toCMframe( event[iDaughter].p(), event[iRecoiler].p() );
1759  else Mtot.toCMframe( event[iRecoiler].p(), event[iDaughter].p() );
1760 
1761  // Initially select phi angle of branching at random.
1762  double phi = 2. * M_PI * rndmPtr->flat();
1763 
1764  // Evaluate coefficient of azimuthal asymmetry from gluon polarization.
1765  findAsymPol( event, dipEndSel);
1766  int iFinPol = dipEndSel->iFinPol;
1767  double cPol = dipEndSel->asymPol;
1768  double phiPol = (iFinPol == 0) ? 0. : event[iFinPol].phi();
1769 
1770  // If interference: try to match sister (anti)colour to final state.
1771  int iFinInt = 0;
1772  double cInt = 0.;
1773  double phiInt = 0.;
1774  if (doPhiIntAsym) {
1775  for (int i = 0; i < partonSystemsPtr->sizeOut(iSysSel); ++ i) {
1776  int iOut = partonSystemsPtr->getOut(iSysSel, i);
1777  if ( (acolSister != 0 && event[iOut].col() == acolSister)
1778  || (colSister != 0 && event[iOut].acol() == colSister) )
1779  iFinInt = iOut;
1780  }
1781  if (iFinInt != 0) {
1782  // Boost final-state parton to current frame of new kinematics.
1783  Vec4 pFinTmp = event[iFinInt].p();
1784  pFinTmp.rotbst(Mtot);
1785  double theFin = pFinTmp.theta();
1786  if (side == 2) theFin = M_PI - theFin;
1787  double theSis = pSister.theta();
1788  if (side == 2) theSis = M_PI - theSis;
1789  cInt = strengthIntAsym * 2. * theSis * theFin
1790  / (pow2(theSis) + pow2(theFin));
1791  phiInt = event[iFinInt].phi();
1792  }
1793  }
1794 
1795  // Bias phi distribution for polarization and interference.
1796  if (iFinPol != 0 || iFinInt != 0) {
1797  double cPhiPol, cPhiInt, weight;
1798  do {
1799  phi = 2. * M_PI * rndmPtr->flat();
1800  weight = 1.;
1801  if (iFinPol !=0 ) {
1802  cPhiPol = cos(phi - phiPol);
1803  weight *= ( 1. + cPol * (2. * pow2(cPhiPol) - 1.) )
1804  / ( 1. + abs(cPol) );
1805  }
1806  if (iFinInt !=0 ) {
1807  cPhiInt = cos(phi - phiInt);
1808  weight *= (1. - cInt) * (1. - cInt * cPhiInt)
1809  / (1. + pow2(cInt) - 2. * cInt * cPhiInt);
1810  }
1811  } while (weight < rndmPtr->flat());
1812  }
1813 
1814  // Include rotation -phi on boost to old rest frame.
1815  Mtot.rot(0., -phi);
1816 
1817  // Find boost from old rest frame to event cm frame.
1818  RotBstMatrix MfromRest;
1819  // The boost to the new rest frame.
1820  Vec4 sumNew = pMother + pNewRec;
1821  double betaX = sumNew.px() / sumNew.e();
1822  double betaZ = sumNew.pz() / sumNew.e();
1823  MfromRest.bst( -betaX, 0., -betaZ);
1824  // Alignment of radiator + recoiler to +- z axis, and rotation +phi.
1825  // Note: with spacelike (E < 0) recoiler p'_x_mother < 0 can happen!
1826  pMother.rotbst(MfromRest);
1827  double theta = pMother.theta();
1828  if (pMother.px() < 0.) theta = -theta;
1829  if (side == 2) theta += M_PI;
1830  MfromRest.rot(-theta, phi);
1831  // Boost to radiator + recoiler in event cm frame.
1832  if (normalRecoil) {
1833  MfromRest.bst( 0., 0., (x1New - x2New) / (x1New + x2New) );
1834  } else if (side == 1) {
1835  Vec4 pMotherWanted( 0., 0., 0.5 * eCM * x1New, 0.5 * eCM * x1New);
1836  MfromRest.fromCMframe( pMotherWanted, event[iRecoiler].p() );
1837 
1838  } else {
1839  Vec4 pMotherWanted( 0., 0., -0.5 * eCM * x2New, 0.5 * eCM * x2New);
1840  MfromRest.fromCMframe( event[iRecoiler].p(), pMotherWanted );
1841  }
1842  Mtot.rotbst(MfromRest);
1843 
1844  // ME correction for weak emissions in the t-channel.
1845  double wt;
1846  if (MEtype == 201 || MEtype == 202 || MEtype == 203 ||
1847  MEtype == 206 || MEtype == 207 || MEtype == 208) {
1848 
1849  // Start by finding the correct outgoing particles
1850  // to use in the ME correction.
1851  Vec4 pA0 = mother.p();
1852  Vec4 pB = newRecoiler.p();
1853  bool sideRad = (abs(side) == 1);
1854  Vec4 p1 = event[5].p();
1855  Vec4 p2 = event[6].p();
1856  int id1 = event[5].id();
1857  int id2 = event[6].id();
1858  if (!tChannel) {swap(p1,p2); swap(id1,id2);}
1859  if (!sideRad) {swap(p1,p2); swap(id1,id2);}
1860 
1861  // Rotate with -phi to keep correct for the later +phi rotation.
1862  p1.rot(0., -phi);
1863  p2.rot(0., -phi);
1864 
1865  // Calculate the actual weight.
1866  wt = calcMEcorrWeak(MEtype, m2, z, pT2, pA0, pB,
1867  event[iDaughter].p(), event[iRecoiler].p(), p1, p2, sister.p());
1868  if (wt > weakMaxWt) {
1869  weakMaxWt = wt;
1870  infoPtr->errorMsg("Warning in SpaceShower::Branch: "
1871  "weight is above unity for weak emission.");
1872  }
1873 
1874  // If weighting fails then restore event record to state before emission.
1875  if (wt < rndmPtr->flat()) {
1876  event.popBack( event.size() - eventSizeOld);
1877  event[beamOff1].daughter1( ev1Dau1V);
1878  event[beamOff2].daughter1( ev2Dau1V);
1879  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1880  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1881  event[iOldCopy].status( statusV[iCopy]);
1882  event[iOldCopy].mothers( mother1V[iCopy], mother2V[iCopy]);
1883  event[iOldCopy].daughters( daughter1V[iCopy], daughter2V[iCopy]);
1884  }
1885  return false;
1886  }
1887  }
1888 
1889  // Perform cumulative rotation/boost operation.
1890  // Mother, recoiler and sister from old rest frame to event cm frame.
1891  mother.rotbst(MfromRest);
1892  newRecoiler.rotbst(MfromRest);
1893  sister.rotbst(MfromRest);
1894  // The rest from (and to) event cm frame.
1895  for ( int i = eventSizeOld + 2; i < eventSizeOld + systemSizeOld; ++i)
1896  event[i].rotbst(Mtot);
1897 
1898  // Remove double counting. Only implemented for QCD hard processes
1899  // and for the first emission.
1900  if (infoPtr->nISR() + infoPtr->nFSRinProc() == 0
1901  && infoPtr->code() > 110 && infoPtr->code() < 130
1902  && MEtype >= 200 && MEtype < 210 && vetoWeakJets) {
1903 
1904  // Find outgoing particles.
1905  int iP1 = event[5].daughter1();
1906  int iP2 = event[6].daughter1();
1907  Vec4 pP1 = event[iP1].p();
1908  Vec4 pP2 = event[iP2].p();
1909 
1910  // Set start pT2 as pT2 of emitted particle and therefore no cut.
1911  double d = sister.pT2();
1912  bool cut = false;
1913 
1914  if (pP1.pT2() < d) {d = pP1.pT2(); cut = true;}
1915  if (pP2.pT2() < d) {d = pP2.pT2(); cut = true;}
1916 
1917  // Check for angle between weak boson and quarks
1918  // (require final state particle to be a fermion).
1919  if (event[iP1].idAbs() < 20) {
1920  double dij = min(pP1.pT2(),sister.pT2())
1921  * pow2(RRapPhi(pP1,sister.p()))/vetoWeakDeltaR2;
1922  if (dij < d) {
1923  d = dij;
1924  cut = false;
1925  }
1926  }
1927 
1928  if (event[iP2].idAbs() < 20) {
1929  double dij = min(pP2.pT2(),sister.pT2())
1930  * pow2(RRapPhi(pP2,sister.p()))/vetoWeakDeltaR2;
1931  if (dij < d) {
1932  d = dij;
1933  cut = false;
1934  }
1935  }
1936 
1937  // Check for angle between recoiler and radiator, if quark anti-quark pair,
1938  // or if the recoiler is a gluon.
1939  if (event[iP1].idAbs() == 21 || event[iP2].idAbs() == 21 ||
1940  event[iP1].id() == - event[iP2].id()) {
1941  double dij = min(pP1.pT2(),pP2.pT2())
1942  * pow2(RRapPhi(pP1,pP2))/vetoWeakDeltaR2;
1943  if (dij < d) {
1944  d = dij;
1945  cut = true;
1946  }
1947  }
1948 
1949  // Clean up event if the emission should be removed.
1950  if (cut) {
1951  event.popBack( event.size() - eventSizeOld);
1952  event[beamOff1].daughter1( ev1Dau1V);
1953  event[beamOff2].daughter1( ev2Dau1V);
1954  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1955  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1956  event[iOldCopy].status( statusV[iCopy]);
1957  event[iOldCopy].mothers( mother1V[iCopy], mother2V[iCopy]);
1958  event[iOldCopy].daughters( daughter1V[iCopy], daughter2V[iCopy]);
1959  }
1960  return false;
1961  }
1962  }
1963 
1964  // Allow veto of branching. If so restore event record to before emission.
1965  if ( (canVetoEmission
1966  && userHooksPtr->doVetoISREmission(eventSizeOld, event, iSysSel))
1967  || (canMergeFirst
1968  && mergingHooksPtr->doVetoEmission( event )) ) {
1969  event.popBack( event.size() - eventSizeOld);
1970  event[beamOff1].daughter1( ev1Dau1V);
1971  event[beamOff2].daughter1( ev2Dau1V);
1972  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1973  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1974  event[iOldCopy].status( statusV[iCopy]);
1975  event[iOldCopy].mothers( mother1V[iCopy], mother2V[iCopy]);
1976  event[iOldCopy].daughters( daughter1V[iCopy], daughter2V[iCopy]);
1977  }
1978  return false;
1979  }
1980 
1981  // Update list of partons in system; adding newly produced one.
1982  partonSystemsPtr->setInA(iSysSel, eventSizeOld);
1983  partonSystemsPtr->setInB(iSysSel, eventSizeOld + 1);
1984  for (int iCopy = 2; iCopy < systemSizeOld; ++iCopy)
1985  partonSystemsPtr->setOut(iSysSel, iCopy - 2, eventSizeOld + iCopy);
1986  partonSystemsPtr->addOut(iSysSel, eventSizeOld + systemSizeOld);
1987  partonSystemsPtr->setSHat(iSysSel, m2 / z);
1988 
1989  // Add dipoles for q -> g q, where the daughter is the gluon.
1990  if (idDaughter == 21 && idMother != 21) {
1991  if (doQEDshowerByQ) {
1992  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
1993  iNewRecoiler, pT2, 0, mother.chargeType(), 0, 0, normalRecoil) );
1994  }
1995  if (doWeakShower && iSysSel == 0) {
1996  int MEtypeNew = 203;
1997  if (idRecoiler == 21) MEtypeNew = 201;
1998  if (idRecoiler == idMother) MEtypeNew = 202;
1999  // If original was a Drell-Yan, keep as Drell-Yan.
2000  if( event[3].id() == - event[4].id()) MEtypeNew = 200;
2001  int weakPol = (rndmPtr->flat() > 0.5) ? -1 : 1;
2002  event[iMother].pol(weakPol);
2003  if ((weakMode == 0 || weakMode == 1) && weakPol == -1)
2004  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2005  iNewRecoiler, pT2, 0, 0, 1, MEtypeNew, normalRecoil, weakPol) );
2006  if (weakMode == 0 || weakMode == 2)
2007  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2008  iNewRecoiler, pT2, 0, 0, 2, MEtypeNew + 5, normalRecoil, weakPol) );
2009  }
2010  }
2011 
2012  // Add dipoles for q -> q gamma, where the daughter is the gamma.
2013  if (idDaughter == 22 && idMother != 22) {
2014  if (doQCDshower && mother.colType() != 0) {
2015  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2016  iNewRecoiler, pT2, mother.colType(), 0, 0, 0, normalRecoil) );
2017  }
2018  if (doQEDshowerByQ && mother.chargeType() != 3) {
2019  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2020  iNewRecoiler, pT2, 0, mother.chargeType(), 0, 0, normalRecoil) );
2021  }
2022  if (doQEDshowerByL && mother.chargeType() == 3) {
2023  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2024  iNewRecoiler, pT2, 0, mother.chargeType(), 0, 0, normalRecoil) );
2025  }
2026  if (doWeakShower && iSysSel == 0) {
2027  int MEtypeNew = 203;
2028  if (idRecoiler == 21) MEtypeNew = 201;
2029  if (idRecoiler == idMother) MEtypeNew = 202;
2030  // If original was a Drell-Yan, keep as Drell-Yan.
2031  if( event[3].id() == - event[4].id()) MEtypeNew = 200;
2032  int weakPol = (rndmPtr->flat() > 0.5) ? -1 : 1;
2033  event[iMother].pol(weakPol);
2034  if ((weakMode == 0 || weakMode == 1) && weakPol == -1)
2035  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2036  iNewRecoiler, pT2, 0, 0, 1, MEtypeNew, normalRecoil, weakPol) );
2037  if (weakMode == 0 || weakMode == 2)
2038  dipEnd.push_back( SpaceDipoleEnd( iSysSel, side, iMother,
2039  iNewRecoiler, pT2, 0, 0, 2, MEtypeNew + 5, normalRecoil, weakPol) );
2040  }
2041  }
2042 
2043  // dipEnd array may have expanded and been moved, so regenerate dipEndSel.
2044  dipEndSel = &dipEnd[iDipSel];
2045 
2046  // Set flag to tell that a weak emission has happened.
2047  if (dipEndSel->weakType != 0) hasWeaklyRadiated = true;
2048 
2049  // Update list of QCD emissions in side A and B in given iSysSel
2050  // This is used to veto jets in W/z events.
2051  while (iSysSel >= int(nRadA.size()) || iSysSel >= int(nRadB.size())) {
2052  nRadA.push_back(0);
2053  nRadB.push_back(0);
2054  }
2055  if (dipEndSel->colType != 0 && side == 1) ++nRadA[iSysSel];
2056  else if (dipEndSel->colType != 0) ++nRadB[iSysSel];
2057 
2058  // Update info on radiating dipole ends (QCD, QED or weak).
2059  for (int iDip = 0; iDip < int(dipEnd.size()); ++iDip)
2060  if ( dipEnd[iDip].system == iSysSel) {
2061  if (abs(dipEnd[iDip].side) == side) {
2062  dipEnd[iDip].iRadiator = iMother;
2063  dipEnd[iDip].iRecoiler = iNewRecoiler;
2064  if (dipEnd[iDip].colType != 0)
2065  dipEnd[iDip].colType = mother.colType();
2066  else if (dipEnd[iDip].chgType != 0) {
2067  dipEnd[iDip].chgType = 0;
2068  if ( (mother.isQuark() && doQEDshowerByQ)
2069  || (mother.isLepton() && doQEDshowerByL) )
2070  dipEnd[iDip].chgType = mother.chargeType();
2071  }
2072  else if (dipEnd[iDip].weakType != 0) {
2073  // Kill weak dipole if mother becomes gluon / photon.
2074  if (!(mother.isLepton() || mother.isQuark()))
2075  dipEnd[iDip].weakType = 0;
2076  if (singleWeakEmission && hasWeaklyRadiated)
2077  dipEnd[iDip].weakType = 0;
2078  }
2079 
2080  // Kill ME corrections after first emission for everything
2081  // but weak showers.
2082  if (dipEnd[iDip].weakType == 0) dipEnd[iDip].MEtype = 0;
2083 
2084  // Update info on recoiling dipole ends (QCD or QED).
2085  } else {
2086  dipEnd[iDip].iRadiator = iNewRecoiler;
2087  dipEnd[iDip].iRecoiler = iMother;
2088  // Optionally also kill recoiler ME corrections after first emission.
2089  if (!doMEafterFirst && dipEnd[iDip].weakType == 0)
2090  dipEnd[iDip].MEtype = 0;
2091  // Remove weak dipoles if we only want a single emission.
2092  if (dipEnd[iDip].weakType != 0 && singleWeakEmission
2093  && hasWeaklyRadiated) dipEnd[iDip].weakType = 0;
2094  }
2095  }
2096 
2097  // Set polarisation of mother for weak emissions.
2098  if (dipEndSel->weakType != 0) mother.pol(dipEndSel->weakPol);
2099 
2100  // Update info on beam remnants.
2101  BeamParticle& beamNow = (side == 1) ? *beamAPtr : *beamBPtr;
2102  double xNew = (side == 1) ? x1New : x2New;
2103  beamNow[iSysSel].update( iMother, idMother, xNew);
2104  // Redo choice of companion kind whenever new flavour.
2105  if (idMother != idDaughterNow) {
2106  pdfScale2 = (useFixedFacScale) ? fixedFacScale2 : factorMultFac * pT2;
2107  beamNow.xfISR( iSysSel, idMother, xNew, pdfScale2);
2108  beamNow.pickValSeaComp();
2109  }
2110  BeamParticle& beamRec = (side == 1) ? *beamBPtr : *beamAPtr;
2111  beamRec[iSysSel].iPos( iNewRecoiler);
2112 
2113  // Store branching values of current dipole. (For rapidity ordering.)
2114  ++dipEndSel->nBranch;
2115  dipEndSel->pT2Old = pT2;
2116  dipEndSel->zOld = z;
2117 
2118  // Update history if recoiler rescatters.
2119  if (!normalRecoil)
2120  event[iRecoilMother].daughters( iNewRecoiler, iNewRecoiler);
2121 
2122  // Start list of rescatterers that force changed kinematics.
2123  vector<int> iRescatterer;
2124  for ( int i = 0; i < systemSizeOld - 2; ++i) {
2125  int iOutNew = partonSystemsPtr->getOut( iSysSel, i);
2126  if (!event[iOutNew].isFinal()) iRescatterer.push_back(iOutNew);
2127  }
2128 
2129  // Start iterate over list of such rescatterers.
2130  int iRescNow = -1;
2131  while (++iRescNow < int(iRescatterer.size())) {
2132 
2133  // Identify partons that induce or are affected by rescatter shift.
2134  // In following Old is before change of kinematics, New after,
2135  // Out scatterer in outstate and In in instate of another system.
2136  // Daughter sequence is (iOutOld ->) iOutNew -> iInNew -> iInOld.
2137  int iOutNew = iRescatterer[iRescNow];
2138  int iInOld = event[iOutNew].daughter1();
2139  int iSysResc = partonSystemsPtr->getSystemOf(iInOld, true);
2140 
2141  // Copy incoming partons of rescattered system and hook them up.
2142  int iOldA = partonSystemsPtr->getInA(iSysResc);
2143  int iOldB = partonSystemsPtr->getInB(iSysResc);
2144  bool rescSideA = event[iOldA].isRescatteredIncoming();
2145  int statusNewA = (rescSideA) ? -45 : -42;
2146  int statusNewB = (rescSideA) ? -42 : -45;
2147  int iNewA = event.copy(iOldA, statusNewA);
2148  int iNewB = event.copy(iOldB, statusNewB);
2149 
2150  // Copy outgoing partons of rescattered system and hook them up.
2151  int eventSize = event.size();
2152  int sizeOutAB = partonSystemsPtr->sizeOut(iSysResc);
2153  int iOldAB, statusOldAB, iNewAB;
2154  for (int iOutAB = 0; iOutAB < sizeOutAB; ++iOutAB) {
2155  iOldAB = partonSystemsPtr->getOut(iSysResc, iOutAB);
2156  statusOldAB = event[iOldAB].status();
2157  iNewAB = event.copy(iOldAB, 44);
2158  // Status could be negative for parton that rescatters in its turn.
2159  if (statusOldAB < 0) {
2160  event[iNewAB].statusNeg();
2161  iRescatterer.push_back(iNewAB);
2162  }
2163  }
2164 
2165  // Hook up new outgoing with new incoming parton.
2166  int iInNew = (rescSideA) ? iNewA : iNewB;
2167  event[iOutNew].daughters( iInNew, iInNew);
2168  event[iInNew].mothers( iOutNew, iOutNew);
2169 
2170  // Rescale recoiling incoming parton for correct invariant mass.
2171  event[iInNew].p( event[iOutNew].p() );
2172  double momFac = (rescSideA)
2173  ? event[iInOld].pPos() / event[iInNew].pPos()
2174  : event[iInOld].pNeg() / event[iInNew].pNeg();
2175  int iInRec = (rescSideA) ? iNewB : iNewA;
2176 
2177  // Rescatter: A previous boost may cause the light cone momentum of a
2178  // rescattered parton to change sign. If this happens, tell
2179  // parton level to try again.
2180  if (momFac < 0.0) {
2181  infoPtr->errorMsg("Warning in SpaceShower::branch: "
2182  "change in lightcone momentum sign; retrying parton level");
2183  rescatterFail = true;
2184  return false;
2185  }
2186  event[iInRec].rescale4( momFac);
2187 
2188  // Boost outgoing partons to new frame of incoming.
2189  RotBstMatrix MmodResc;
2190  MmodResc.toCMframe( event[iOldA].p(), event[iOldB].p());
2191  MmodResc.fromCMframe(event[iNewA].p(), event[iNewB].p());
2192  for (int iOutAB = 0; iOutAB < sizeOutAB; ++iOutAB)
2193  event[eventSize + iOutAB].rotbst(MmodResc);
2194 
2195  // Update list of partons in system.
2196  partonSystemsPtr->setInA(iSysResc, iNewA);
2197  partonSystemsPtr->setInB(iSysResc, iNewB);
2198  for (int iCopy = 0; iCopy < sizeOutAB; ++iCopy)
2199  partonSystemsPtr->setOut(iSysResc, iCopy, eventSize + iCopy);
2200 
2201  // Update info on radiating dipole ends (QCD or QED).
2202  for (int iDip = 0; iDip < int(dipEnd.size()); ++iDip)
2203  if ( dipEnd[iDip].system == iSysResc) {
2204  bool sideAnow = (abs(dipEnd[iDip].side) == 1);
2205  dipEnd[iDip].iRadiator = (sideAnow) ? iNewA : iNewB;
2206  dipEnd[iDip].iRecoiler = (sideAnow) ? iNewB : iNewA;
2207  }
2208 
2209  // Update info on beam remnants.
2210  BeamParticle& beamResc = (rescSideA) ? *beamAPtr : *beamBPtr;
2211  beamResc[iSysResc].iPos( iInNew);
2212  beamResc[iSysResc].p( event[iInNew].p() );
2213  beamResc[iSysResc].scaleX( 1. / momFac );
2214  BeamParticle& beamReco = (rescSideA) ? *beamBPtr : *beamAPtr;
2215  beamReco[iSysResc].iPos( iInRec);
2216  beamReco[iSysResc].scaleX( momFac);
2217 
2218  // End iterate over list of rescatterers.
2219  }
2220 
2221  // Check that beam momentum not used up by rescattered-system boosts.
2222  if (beamAPtr->xMax(-1) < 0.0 || beamBPtr->xMax(-1) < 0.0) {
2223  infoPtr->errorMsg("Warning in SpaceShower::branch: "
2224  "used up beam momentum; retrying parton level");
2225  rescatterFail = true;
2226  return false;
2227  }
2228 
2229  // Done without any errors.
2230  return true;
2231 
2232 }
2233 
2234 //--------------------------------------------------------------------------
2235 
2236 // Find class of ME correction.
2237 
2238  int SpaceShower::findMEtype( int iSys, Event& event, bool weakRadiation) {
2239 
2240  // Default values and no action.
2241  int MEtype = 0;
2242  if (!doMEcorrections) return MEtype;
2243 
2244  // Identify systems producing a single resonance.
2245  if (partonSystemsPtr->sizeOut( iSys) == 1 && !weakRadiation) {
2246  int idIn1 = event[partonSystemsPtr->getInA(iSys)].id();
2247  int idIn2 = event[partonSystemsPtr->getInA(iSys)].id();
2248  int idRes = event[partonSystemsPtr->getOut(iSys, 0)].id();
2249  if (iSys == 0) idResFirst = abs(idRes);
2250  if (iSys == 1) idResSecond = abs(idRes);
2251 
2252  // f + fbar -> vector boson.
2253  if ( (idRes == 23 || abs(idRes) == 24 || idRes == 32
2254  || idRes == 33 || abs(idRes) == 34 || abs(idRes) == 41)
2255  && abs(idIn1) < 20 && abs(idIn2) < 20 ) MEtype = 1;
2256 
2257  // g + g, gamma + gamma -> Higgs boson.
2258  if ( (idRes == 25 || idRes == 35 || idRes == 36)
2259  && ( ( idIn1 == 21 && idIn2 == 21 )
2260  || ( idIn1 == 22 && idIn2 == 22 ) ) ) MEtype = 2;
2261 
2262  // f + fbar -> Higgs boson.
2263  if ( (idRes == 25 || idRes == 35 || idRes == 36)
2264  && abs(idIn1) < 20 && abs(idIn2) < 20 ) MEtype = 3;
2265  }
2266 
2267  // Weak ME corrections.
2268  if (weakRadiation) {
2269  if (event[3].id() == -event[4].id()
2270  || event[event[3].daughter1()].idAbs() == 24 || infoPtr->nFinal() != 2)
2271  MEtype = 200;
2272  else if (event[3].idAbs() == 21 || event[4].idAbs() == 21) MEtype = 201;
2273  else if (event[3].id() == event[4].id()) MEtype = 202;
2274  else MEtype = 203;
2275  }
2276 
2277  // Done.
2278  return MEtype;
2279 
2280 }
2281 
2282 //--------------------------------------------------------------------------
2283 
2284 // Provide maximum of expected ME weight; for preweighting of evolution.
2285 
2286 double SpaceShower::calcMEmax( int MEtype, int idMother, int idDaughterIn) {
2287 
2288  // Main non-unity case: g(gamma) f -> V f'.
2289  if (MEtype == 1 && idMother > 20 && idDaughterIn < 20) return 3.;
2290 
2291  // Added a case for t-channel W/Z exchange, since the PS is not an
2292  // overestimate. This does not help fully, but it should only be small
2293  // pT quarks / gluons that break the overscattering.
2294  if ( MEtype == 201 || MEtype == 202 || MEtype == 203
2295  || MEtype == 206 || MEtype == 207 || MEtype == 208) return WEAKPSWEIGHT;
2296 
2297  // Default.
2298  return 1.;
2299 
2300 }
2301 
2302 //--------------------------------------------------------------------------
2303 
2304 // Provide actual ME weight for current branching.
2305 // Note: currently ME corrections are only allowed for first branching
2306 // on each side, so idDaughter is essentially known and checks overkill.
2307 
2308 double SpaceShower::calcMEcorr(int MEtype, int idMother, int idDaughterIn,
2309  double M2, double z, double Q2, double m2Sister) {
2310 
2311  // Convert to Mandelstam variables. Sometimes may need to swap later.
2312  double sH = M2 / z;
2313  double tH = -Q2;
2314  double uH = Q2 - M2 * (1. - z) / z;
2315  int idMabs = abs(idMother);
2316  int idDabs = abs(idDaughterIn);
2317 
2318  // Corrections for f + fbar -> s-channel vector boson.
2319  if (MEtype == 1) {
2320  if (idMabs < 20 && idDabs < 20) {
2321  return (tH*tH + uH*uH + 2. * M2 * sH) / (sH*sH + M2*M2);
2322  } else if (idDabs < 20) {
2323  // g(gamma) f -> V f': -Q2 = (p_g - p_f')^2 in PS while
2324  // tHat = (p_f - p_f')^2 in ME so need to swap tHat <-> uHat.
2325  swap( tH, uH);
2326  return (sH*sH + uH*uH + 2. * M2 * tH) / (pow2(sH - M2) + M2*M2);
2327  }
2328 
2329  // Corrections for g + g -> Higgs boson.
2330  } else if (MEtype == 2) {
2331  if (idMabs < 20 && idDabs > 20) {
2332  return (sH*sH + uH*uH) / (sH*sH + pow2(sH - M2));
2333  } else if (idDabs > 20) {
2334  return 0.5 * (pow4(sH) + pow4(tH) + pow4(uH) + pow4(M2))
2335  / pow2(sH*sH - M2 * (sH - M2));
2336  }
2337 
2338  // Corrections for f + fbar -> Higgs boson (f = b mainly).
2339  } else if (MEtype == 3) {
2340  if (idMabs < 20 && idDabs < 20) {
2341  // The PS and ME answers agree for f fbar -> H g/gamma.
2342  return 1.;
2343  } else if (idDabs < 20) {
2344  // Need to swap tHat <-> uHat, cf. vector-boson production above.
2345  swap( tH, uH);
2346  return (sH*sH + uH*uH + 2. * (M2 - uH) * (M2 - sH))
2347  / (pow2(sH - M2) + M2*M2);
2348  }
2349 
2350  // Corrections for f -> f' + W/Z (s-channel).
2351  } else if (MEtype == 200 || MEtype == 205) {
2352  // Need to redo calculations of uH since we now emit a massive particle.
2353  uH += m2Sister;
2354  double wtME = (uH*uH + tH*tH + 2 * sH * (m2Sister + M2)) / (uH*tH)
2355  - M2 * m2Sister * (1/(tH*tH) + 1/(uH*uH));
2356  double wtPS = (sH*sH + pow2(M2 + m2Sister)) / (tH*uH);
2357  return wtME / wtPS;
2358  } else if (MEtype == 201 || MEtype == 202 || MEtype == 203 ||
2359  MEtype == 206 || MEtype == 207 || MEtype == 208)
2360  return calcMEmax(MEtype, 0, 0);
2361 
2362  // Default.
2363  return 1.;
2364 
2365 }
2366 
2367 //--------------------------------------------------------------------------
2368 
2369 // Provide actual ME weight for current branching for weak t-channel emissions.
2370 
2371 double SpaceShower::calcMEcorrWeak(int MEtype, double m2, double z,
2372  double pT2, Vec4 pMother, Vec4 pB, Vec4 pDaughter,
2373  Vec4 pB0, Vec4 p1, Vec4 p2, Vec4 pSister) {
2374 
2375  // Find daughter four-momentum in current frame.
2376  Vec4 pA = pMother - pSister;
2377 
2378  // Scale outgoing vectors to conserve energy / momentum.
2379  double scaleFactor2 = (pA + pB).m2Calc() / (p1 + p2).m2Calc();
2380  double scaleFactor = sqrt(scaleFactor2);
2381  RotBstMatrix rot2to2frame;
2382  rot2to2frame.bstback(p1 + p2);
2383  p1.rotbst(rot2to2frame);
2384  p2.rotbst(rot2to2frame);
2385  p1 *= scaleFactor;
2386  p2 *= scaleFactor;
2387 
2388  // Find 2 to 2 rest frame for incoming particles.
2389  // This is done before one of the two are made virtual (Q^2 mass).
2390  RotBstMatrix rot2to2frameInc;
2391  rot2to2frameInc.bstback(pDaughter + pB0);
2392  pDaughter.rotbst(rot2to2frameInc);
2393  pB0.rotbst(rot2to2frameInc);
2394  double sHat = (p1 + p2).m2Calc();
2395  double tHat = (p1 - pDaughter).m2Calc();
2396  double uHat = (p1 - pB0).m2Calc();
2397 
2398  // Calculate the weak t-channel correction.
2399  double m2R1 = 1. + pSister.m2Calc() / m2;
2400  double wt = 4. * sHat / (pMother + pB).m2Calc() * pT2 * ( 1. - z * m2R1)
2401  / (1. + pow2(z * m2R1)) / (1.-z);
2402  if (MEtype == 201 || MEtype == 206)
2403  wt *= weakShowerMEs.getTchanneluGuGZME(pMother, pB, p2, pSister, p1)
2404  / weakShowerMEs.getTchanneluGuGME(sHat, tHat, uHat);
2405  else if (MEtype == 202 || MEtype == 207)
2406  wt *= weakShowerMEs.getTchannelududZME(pMother, pB, pSister, p2, p1)
2407  / weakShowerMEs.getTchanneluuuuME(sHat, tHat, uHat);
2408  else if (MEtype == 203 || MEtype == 208)
2409  wt *= weakShowerMEs.getTchannelududZME(pMother, pB, pSister, p2, p1)
2410  / weakShowerMEs.getTchannelududME(sHat, tHat, uHat);
2411 
2412  // Split of ME into an ISR part and FSR part.
2413  wt *= (pSister + p1).m2Calc() / ( (pSister + p1).m2Calc()
2414  + abs((-pMother + pSister).m2Calc()) );
2415 
2416  // Remove the addition weight that was used to get an overestimate.
2417  wt /= calcMEmax(MEtype, 0, 0);
2418 
2419  return wt;
2420 }
2421 
2422 //--------------------------------------------------------------------------
2423 
2424 // Find coefficient of azimuthal asymmetry from gluon polarization.
2425 
2426 void SpaceShower::findAsymPol( Event& event, SpaceDipoleEnd* dip) {
2427 
2428  // Default is no asymmetry. Only gluons are studied.
2429  dip->iFinPol = 0;
2430  dip->asymPol = 0.;
2431  int iRad = dip->iRadiator;
2432  if (!doPhiPolAsym || dip->idDaughter != 21) return;
2433 
2434  // At least two particles in final state, whereof at least one coloured.
2435  int systemSizeOut = partonSystemsPtr->sizeOut( iSysSel);
2436  if (systemSizeOut < 2) return;
2437  bool foundColOut = false;
2438  for (int ii = 0; ii < systemSizeOut; ++ii) {
2439  int i = partonSystemsPtr->getOut( iSysSel, ii);
2440  if (event[i].col() != 0 || event[i].acol() != 0) foundColOut = true;
2441  }
2442  if (!foundColOut) return;
2443 
2444  // Check if granddaughter in final state of hard scattering.
2445  // (May need to trace across carbon copies to find granddaughters.)
2446  // If so, only accept 2 -> 2 scatterings with gg or qq in final state.
2447  int iGrandD1 = event[iRad].daughter1();
2448  int iGrandD2 = event[iRad].daughter2();
2449  bool traceCopy = false;
2450  do {
2451  traceCopy = false;
2452  if (iGrandD1 > 0 && iGrandD2 == iGrandD1) {
2453  iGrandD1 = event[iGrandD2].daughter1();
2454  iGrandD2 = event[iGrandD2].daughter2();
2455  traceCopy = true;
2456  }
2457  } while (traceCopy);
2458  int statusGrandD1 = event[ iGrandD1 ].statusAbs();
2459  bool isHardProc = (statusGrandD1 == 23 || statusGrandD1 == 33);
2460  if (isHardProc) {
2461  if (iGrandD2 != iGrandD1 + 1) return;
2462  if (event[iGrandD1].isGluon() && event[iGrandD2].isGluon());
2463  else if (event[iGrandD1].isQuark() && event[iGrandD2].isQuark());
2464  else return;
2465  }
2466  dip->iFinPol = iGrandD1;
2467 
2468  // Coefficient from gluon production.
2469  if (dip->idMother == 21) dip->asymPol = pow2( (1. - dip->z)
2470  / (1. - dip->z * (1. - dip->z) ) );
2471  else dip->asymPol = 2. * (1. - dip->z) / (1. + pow2(1. - dip->z) );
2472 
2473  // Coefficients from gluon decay. Put z = 1/2 for hard process.
2474  double zDau = (isHardProc) ? 0.5 : dip->zOld;
2475  if (event[iGrandD1].isGluon()) dip->asymPol *= pow2( (1. - zDau)
2476  / (1. - zDau * (1. - zDau) ) );
2477  else dip->asymPol *= -2. * zDau *( 1. - zDau )
2478  / (1. - 2. * zDau * (1. - zDau) );
2479 
2480 }
2481 
2482 //--------------------------------------------------------------------------
2483 
2484 // Remove weak dipoles if FSR already emitted a W/Z
2485 // and only a single weak emission is permited.
2486 
2487 void SpaceShower::update(int , Event &, bool hasWeakRad) {
2488 
2489  if (hasWeakRad && singleWeakEmission)
2490  for (int i = 0; i < int(dipEnd.size()); i++)
2491  if (dipEnd[i].weakType != 0) dipEnd[i].weakType = 0;
2492  if (hasWeakRad) hasWeaklyRadiated = true;
2493 
2494 }
2495 
2496 //-------------------------------------------------------------------------
2497 
2498 // Print the list of dipoles.
2499 
2500 void SpaceShower::list(ostream& os) const {
2501 
2502  // Header.
2503  os << "\n -------- PYTHIA SpaceShower Dipole Listing -------------- \n"
2504  << "\n i syst side rad rec pTmax col chg ME rec \n"
2505  << fixed << setprecision(3);
2506 
2507  // Loop over dipole list and print it.
2508  for (int i = 0; i < int(dipEnd.size()); ++i)
2509  os << setw(5) << i << setw(6) << dipEnd[i].system
2510  << setw(6) << dipEnd[i].side << setw(6) << dipEnd[i].iRadiator
2511  << setw(6) << dipEnd[i].iRecoiler << setw(12) << dipEnd[i].pTmax
2512  << setw(5) << dipEnd[i].colType << setw(5) << dipEnd[i].chgType
2513  << setw(5) << dipEnd[i].MEtype << setw(4)
2514  << dipEnd[i].normalRecoil << "\n";
2515 
2516  // Done.
2517  os << "\n -------- End PYTHIA SpaceShower Dipole Listing ----------"
2518  << endl;
2519 
2520 
2521 }
2522 
2523 //==========================================================================
2524 
2525 } // end namespace Pythia8
2526 
Definition: beam.h:43
Definition: AgUStep.h:26