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) 2012 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 "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 // Turn to true to allow debug printout.
23 const bool SpaceShower::DEBUG = false;
24 
25 // Leftover companion can give PDF > 0 at small Q2 where other PDF's = 0,
26 // and then one can end in infinite loop of impossible kinematics.
27 const int SpaceShower::MAXLOOPTINYPDF = 10;
28 
29 // Switch to alternative (but equivalent) backwards evolution for
30 // g -> Q Qbar (Q = c or b) when below QTHRESHOLD * mQ2.
31 const double SpaceShower::CTHRESHOLD = 2.0;
32 const double SpaceShower::BTHRESHOLD = 2.0;
33 
34 // Renew evaluation of PDF's when the pT2 step is bigger than this
35 // (in addition to initial scale and c and b thresholds.)
36 const double SpaceShower::EVALPDFSTEP = 0.1;
37 
38 // Lower limit on PDF value in order to avoid division by zero.
39 const double SpaceShower::TINYPDF = 1e-10;
40 
41 // Lower limit on estimated evolution rate, below which stop.
42 const double SpaceShower::TINYKERNELPDF = 1e-6;
43 
44 // Lower limit on pT2, below which branching is rejected.
45 const double SpaceShower::TINYPT2 = 0.25e-6;
46 
47 // No attempt to do backwards evolution of a heavy (c or b) quark
48 // if evolution starts at a scale pT2 < HEAVYPT2EVOL * mQ2.
49 const double SpaceShower::HEAVYPT2EVOL = 1.1;
50 
51 // No attempt to do backwards evolution of a heavy (c or b) quark
52 // if evolution starts at a x > HEAVYXEVOL * x_max, where
53 // x_max is the largest possible x value for a g -> Q Qbar branching.
54 const double SpaceShower::HEAVYXEVOL = 0.9;
55 
56 // When backwards evolution Q -> g + Q creates a heavy quark Q,
57 // an earlier branching g -> Q + Qbar will restrict kinematics
58 // to M_{Q Qbar}^2 > EXTRASPACEQ * 4 m_Q^2. (Smarter to be found??)
59 const double SpaceShower::EXTRASPACEQ = 2.0;
60 
61 // Never pick pT so low that alphaS is evaluated too close to Lambda_3.
62 const double SpaceShower::LAMBDA3MARGIN = 1.1;
63 
64 // Cutoff for f_e^e at x < 1 - 10^{-10} to be used in z selection.
65 // Note: the x_min quantity come from 1 - x_max.
66 const double SpaceShower::LEPTONXMIN = 1e-10;
67 const double SpaceShower::LEPTONXMAX = 1. - 1e-10;
68 
69 // Stop l -> l gamma evolution slightly above m2l.
70 const double SpaceShower::LEPTONPT2MIN = 1.2;
71 
72 // Enhancement of l -> l gamma trial rate to compensate imperfect modelling.
73 const double SpaceShower::LEPTONFUDGE = 10.;
74 
75 //--------------------------------------------------------------------------
76 
77 // Initialize alphaStrong, alphaEM and related pTmin parameters.
78 
79 void SpaceShower::init( BeamParticle* beamAPtrIn,
80  BeamParticle* beamBPtrIn) {
81 
82  // Store input pointers for future use.
83  beamAPtr = beamAPtrIn;
84  beamBPtr = beamBPtrIn;
85 
86  // Main flags to switch on and off branchings.
87  doQCDshower = settingsPtr->flag("SpaceShower:QCDshower");
88  doQEDshowerByQ = settingsPtr->flag("SpaceShower:QEDshowerByQ");
89  doQEDshowerByL = settingsPtr->flag("SpaceShower:QEDshowerByL");
90 
91  // Matching in pT of hard interaction to shower evolution.
92  pTmaxMatch = settingsPtr->mode("SpaceShower:pTmaxMatch");
93  pTdampMatch = settingsPtr->mode("SpaceShower:pTdampMatch");
94  pTmaxFudge = settingsPtr->parm("SpaceShower:pTmaxFudge");
95  pTmaxFudgeMPI = settingsPtr->parm("SpaceShower:pTmaxFudgeMPI");
96  pTdampFudge = settingsPtr->parm("SpaceShower:pTdampFudge");
97 
98  // Optionally force emissions to be ordered in rapidity/angle.
99  doRapidityOrder = settingsPtr->flag("SpaceShower:rapidityOrder");
100 
101  // Charm, bottom and lepton mass thresholds.
102  mc = particleDataPtr->m0(4);
103  mb = particleDataPtr->m0(5);
104  m2c = pow2(mc);
105  m2b = pow2(mb);
106 
107  // Parameters of alphaStrong generation.
108  alphaSvalue = settingsPtr->parm("SpaceShower:alphaSvalue");
109  alphaSorder = settingsPtr->mode("SpaceShower:alphaSorder");
110  alphaS2pi = 0.5 * alphaSvalue / M_PI;
111 
112  // Initialize alpha_strong generation.
113  alphaS.init( alphaSvalue, alphaSorder);
114 
115  // Lambda for 5, 4 and 3 flavours.
116  Lambda5flav = alphaS.Lambda5();
117  Lambda4flav = alphaS.Lambda4();
118  Lambda3flav = alphaS.Lambda3();
119  Lambda5flav2 = pow2(Lambda5flav);
120  Lambda4flav2 = pow2(Lambda4flav);
121  Lambda3flav2 = pow2(Lambda3flav);
122 
123  // Regularization of QCD evolution for pT -> 0. Can be taken
124  // same as for multiparton interactions, or be set separately.
125  useSamePTasMPI = settingsPtr->flag("SpaceShower:samePTasMPI");
126  if (useSamePTasMPI) {
127  pT0Ref = settingsPtr->parm("MultipartonInteractions:pT0Ref");
128  ecmRef = settingsPtr->parm("MultipartonInteractions:ecmRef");
129  ecmPow = settingsPtr->parm("MultipartonInteractions:ecmPow");
130  pTmin = settingsPtr->parm("MultipartonInteractions:pTmin");
131  } else {
132  pT0Ref = settingsPtr->parm("SpaceShower:pT0Ref");
133  ecmRef = settingsPtr->parm("SpaceShower:ecmRef");
134  ecmPow = settingsPtr->parm("SpaceShower:ecmPow");
135  pTmin = settingsPtr->parm("SpaceShower:pTmin");
136  }
137 
138  // Calculate nominal invariant mass of events. Set current pT0 scale.
139  sCM = m2( beamAPtr->p(), beamBPtr->p());
140  eCM = sqrt(sCM);
141  pT0 = pT0Ref * pow(eCM / ecmRef, ecmPow);
142 
143  // Restrict pTmin to ensure that alpha_s(pTmin^2 + pT_0^2) does not blow up.
144  double pTminAbs = sqrtpos(pow2(LAMBDA3MARGIN) * Lambda3flav2 - pT0*pT0);
145  if (pTmin < pTminAbs) {
146  pTmin = pTminAbs;
147  ostringstream newPTmin;
148  newPTmin << fixed << setprecision(3) << pTmin;
149  infoPtr->errorMsg("Warning in SpaceShower::init: pTmin too low",
150  ", raised to " + newPTmin.str() );
151  infoPtr->setTooLowPTmin(true);
152  }
153 
154  // Parameters of alphaEM generation.
155  alphaEMorder = settingsPtr->mode("SpaceShower:alphaEMorder");
156 
157  // Initialize alphaEM generation.
158  alphaEM.init( alphaEMorder, settingsPtr);
159 
160  // Parameters of QED evolution.
161  pTminChgQ = settingsPtr->parm("SpaceShower:pTminchgQ");
162  pTminChgL = settingsPtr->parm("SpaceShower:pTminchgL");
163 
164  // Derived parameters of QCD evolution.
165  pT20 = pow2(pT0);
166  pT2min = pow2(pTmin);
167  pT2minChgQ = pow2(pTminChgQ);
168  pT2minChgL = pow2(pTminChgL);
169 
170  // Various other parameters.
171  doMEcorrections = settingsPtr->flag("SpaceShower:MEcorrections");
172  doMEafterFirst = settingsPtr->flag("SpaceShower:MEafterFirst");
173  doPhiPolAsym = settingsPtr->flag("SpaceShower:phiPolAsym");
174  doPhiIntAsym = settingsPtr->flag("SpaceShower:phiIntAsym");
175  strengthIntAsym = settingsPtr->parm("SpaceShower:strengthIntAsym");
176  nQuarkIn = settingsPtr->mode("SpaceShower:nQuarkIn");
177 
178  // Optional dampening at small pT's when large multiplicities.
179  enhanceScreening
180  = settingsPtr->mode("MultipartonInteractions:enhanceScreening");
181  if (!useSamePTasMPI) enhanceScreening = 0;
182 
183  // Possibility to allow user veto of emission step.
184  canVetoEmission = (userHooksPtr > 0) ? userHooksPtr->canVetoISREmission()
185  : false;
186 
187 }
188 
189 //--------------------------------------------------------------------------
190 
191 // Find whether to limit maximum scale of emissions.
192 // Also allow for dampening at factorization or renormalization scale.
193 
194 bool SpaceShower::limitPTmax( Event& event, double Q2Fac, double Q2Ren) {
195 
196  // Find whether to limit pT. Begin by user-set cases.
197  bool dopTlimit = false;
198  if (pTmaxMatch == 1) dopTlimit = true;
199  else if (pTmaxMatch == 2) dopTlimit = false;
200 
201  // Look if any quark (u, d, s, c, b), gluon or photon in final state.
202  else {
203  for (int i = 5; i < event.size(); ++i)
204  if (event[i].status() != -21) {
205  int idAbs = event[i].idAbs();
206  if (idAbs <= 5 || idAbs == 21 || idAbs == 22) dopTlimit = true;
207  }
208  }
209 
210  // Dampening at factorization or renormalization scale.
211  dopTdamp = false;
212  pT2damp = 0.;
213  if ( !dopTlimit && (pTdampMatch == 1 || pTdampMatch == 2) ) {
214  dopTdamp = true;
215  pT2damp = pow2(pTdampFudge) * ((pTdampMatch == 1) ? Q2Fac : Q2Ren);
216  }
217 
218  // Done.
219  return dopTlimit;
220 
221 }
222 
223 //--------------------------------------------------------------------------
224 
225 // Prepare system for evolution; identify ME.
226 // Routine may be called after multiparton interactions, for a new subystem.
227 
228 void SpaceShower::prepare( int iSys, Event& event, bool limitPTmaxIn) {
229 
230  // Find positions of incoming colliding partons.
231  int in1 = partonSystemsPtr->getInA(iSys);
232  int in2 = partonSystemsPtr->getInB(iSys);
233 
234  // Rescattered partons cannot radiate.
235  bool canRadiate1 = !(event[in1].isRescatteredIncoming());
236  bool canRadiate2 = !(event[in2].isRescatteredIncoming());
237 
238  // Reset dipole-ends list for first interaction. Also resonances.
239  if (iSys == 0) dipEnd.resize(0);
240  if (iSys == 0) idResFirst = 0;
241  if (iSys == 1) idResSecond = 0;
242 
243  // Find matrix element corrections for system.
244  int MEtype = findMEtype( iSys, event);
245 
246  // Maximum pT scale for dipole ends.
247  double pTmax1 = (limitPTmaxIn) ? event[in1].scale() : eCM;
248  double pTmax2 = (limitPTmaxIn) ? event[in2].scale() : eCM;
249  if (iSys == 0 && limitPTmaxIn) {
250  pTmax1 *= pTmaxFudge;
251  pTmax2 *= pTmaxFudge;
252  } else if (iSys > 0 && limitPTmaxIn) {
253  pTmax1 *= pTmaxFudgeMPI;
254  pTmax2 *= pTmaxFudgeMPI;
255  }
256 
257  // Find dipole ends for QCD radiation.
258  // Note: colour type can change during evolution, so book also if zero.
259  if (doQCDshower) {
260  int colType1 = event[in1].colType();
261  if (canRadiate1) dipEnd.push_back( SpaceDipoleEnd( iSys, 1,
262  in1, in2, pTmax1, colType1, 0, MEtype, canRadiate2) );
263  int colType2 = event[in2].colType();
264  if (canRadiate2) dipEnd.push_back( SpaceDipoleEnd( iSys, 2,
265  in2, in1, pTmax2, colType2, 0, MEtype, canRadiate1) );
266  }
267 
268  // Find dipole ends for QED radiation.
269  // Note: charge type can change during evolution, so book also if zero.
270  if (doQEDshowerByQ || doQEDshowerByL) {
271  int chgType1 = ( (event[in1].isQuark() && doQEDshowerByQ)
272  || (event[in1].isLepton() && doQEDshowerByL) )
273  ? event[in1].chargeType() : 0;
274  // Special: photons have charge zero, but can evolve (only off Q for now)
275  if (event[in1].id() == 22 && doQEDshowerByQ) chgType1 = 22 ;
276  if (canRadiate1) dipEnd.push_back( SpaceDipoleEnd( iSys, -1,
277  in1, in2, pTmax1, 0, chgType1, MEtype, canRadiate2) );
278  int chgType2 = ( (event[in2].isQuark() && doQEDshowerByQ)
279  || (event[in2].isLepton() && doQEDshowerByL) )
280  ? event[in2].chargeType() : 0;
281  // Special: photons have charge zero, but can evolve (only off Q for now)
282  if (event[in2].id() == 22 && doQEDshowerByQ) chgType2 = 22 ;
283  if (canRadiate2) dipEnd.push_back( SpaceDipoleEnd( iSys, -2,
284  in2, in1, pTmax2, 0, chgType2, MEtype, canRadiate1) );
285  }
286 
287  // Store the z and pT2 values of the last previous splitting
288  // when an event history has already been constructed.
289  if (iSys == 0 && infoPtr->hasHistory()) {
290  double zNow = infoPtr->zNowISR();
291  double pT2Now = infoPtr->pT2NowISR();
292  for (int iDipEnd = 0; iDipEnd < int(dipEnd.size()); ++iDipEnd) {
293  dipEnd[iDipEnd].zOld = zNow;
294  dipEnd[iDipEnd].pT2Old = pT2Now;
295  ++dipEnd[iDipEnd].nBranch;
296  }
297  }
298 
299 }
300 
301 //--------------------------------------------------------------------------
302 
303 // Select next pT in downwards evolution of the existing dipoles.
304 
305 double SpaceShower::pTnext( Event& event, double pTbegAll, double pTendAll,
306  int nRadIn) {
307 
308  // Current cm energy, in case it varies between events.
309  sCM = m2( beamAPtr->p(), beamBPtr->p());
310  eCM = sqrt(sCM);
311  pTbegRef = pTbegAll;
312 
313  // Starting values: no radiating dipole found.
314  nRad = nRadIn;
315  double pT2sel = pow2(pTendAll);
316  iDipSel = 0;
317  iSysSel = 0;
318  dipEndSel = 0;
319 
320  // Loop over all possible dipole ends.
321  for (int iDipEnd = 0; iDipEnd < int(dipEnd.size()); ++iDipEnd) {
322  iDipNow = iDipEnd;
323  dipEndNow = &dipEnd[iDipEnd];
324  iSysNow = dipEndNow->system;
325  dipEndNow->pT2 = 0.;
326 
327  // Check whether dipole end should be allowed to shower.
328  double pT2begDip = pow2( min( pTbegAll, dipEndNow->pTmax ));
329  if (pT2begDip > pT2sel
330  && ( dipEndNow->colType != 0 || dipEndNow->chgType != 0 ) ) {
331  double pT2endDip = 0.;
332 
333  // Determine lower cut for evolution, for QCD or QED (q or l).
334  if (dipEndNow->colType != 0) pT2endDip = max( pT2sel, pT2min );
335  else if (abs(dipEndNow->chgType) != 3) pT2endDip
336  = max( pT2sel, pT2minChgQ );
337  else pT2endDip = max( pT2sel, pT2minChgL );
338 
339  // Find properties of dipole and radiating dipole end.
340  sideA = ( abs(dipEndNow->side) == 1 );
341  BeamParticle& beamNow = (sideA) ? *beamAPtr : *beamBPtr;
342  BeamParticle& beamRec = (sideA) ? *beamBPtr : *beamAPtr;
343  iNow = beamNow[iSysNow].iPos();
344  iRec = beamRec[iSysNow].iPos();
345  idDaughter = beamNow[iSysNow].id();
346  xDaughter = beamNow[iSysNow].x();
347  x1Now = (sideA) ? xDaughter : beamRec[iSysNow].x();
348  x2Now = (sideA) ? beamRec[iSysNow].x() : xDaughter;
349 
350  // Note dipole mass correction when recoiler is a rescatter.
351  m2Rec = (dipEndNow->normalRecoil) ? 0. : event[iRec].m2();
352  m2Dip = x1Now * x2Now * sCM + m2Rec;
353 
354  // Now do evolution in pT2, for QCD or QED
355  if (pT2begDip > pT2endDip) {
356  if (dipEndNow->colType != 0) pT2nextQCD( pT2begDip, pT2endDip);
357  else pT2nextQED( pT2begDip, pT2endDip);
358  }
359 
360  // Update if found larger pT than current maximum.
361  if (dipEndNow->pT2 > pT2sel) {
362  pT2sel = dipEndNow->pT2;
363  iDipSel = iDipNow;
364  iSysSel = iSysNow;
365  dipEndSel = dipEndNow;
366  }
367 
368  // End loop over dipole ends.
369  }
370  }
371 
372  // Return nonvanishing value if found pT is bigger than already found.
373  return (dipEndSel == 0) ? 0. : sqrt(pT2sel);
374 }
375 
376 //--------------------------------------------------------------------------
377 
378 // Evolve a QCD dipole end.
379 
380 void SpaceShower::pT2nextQCD( double pT2begDip, double pT2endDip) {
381 
382  // Some properties and kinematical starting values.
383  BeamParticle& beam = (sideA) ? *beamAPtr : *beamBPtr;
384  bool isGluon = (idDaughter == 21);
385  bool isValence = beam[iSysNow].isValence();
386  int MEtype = dipEndNow->MEtype;
387  double pT2 = pT2begDip;
388  double xMaxAbs = beam.xMax(iSysNow);
389  double zMinAbs = xDaughter / xMaxAbs;
390  if (xMaxAbs < 0.) {
391  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQCD: "
392  "xMaxAbs negative");
393  return;
394  }
395 
396  // Starting values for handling of massive quarks (c/b), if any.
397  double idMassive = 0;
398  if ( abs(idDaughter) == 4 ) idMassive = 4;
399  if ( abs(idDaughter) == 5 ) idMassive = 5;
400  bool isMassive = (idMassive > 0);
401  double m2Massive = 0.;
402  double mRatio = 0.;
403  double zMaxMassive = 1.;
404  double m2Threshold = pT2;
405 
406  // Evolution below scale of massive quark or at large x is impossible.
407  if (isMassive) {
408  m2Massive = (idMassive == 4) ? m2c : m2b;
409  if (pT2 < HEAVYPT2EVOL * m2Massive) return;
410  mRatio = sqrt( m2Massive / m2Dip );
411  zMaxMassive = (1. - mRatio) / ( 1. + mRatio * (1. - mRatio) );
412  if (xDaughter > HEAVYXEVOL * zMaxMassive * xMaxAbs) return;
413 
414  // Find threshold scale below which only g -> Q + Qbar will be allowed.
415  m2Threshold = (idMassive == 4) ? min( pT2, CTHRESHOLD * m2c)
416  : min( pT2, BTHRESHOLD * m2b);
417  }
418 
419  // Variables used inside evolution loop. (Mainly dummy starting values.)
420  int nFlavour = 3;
421  double b0 = 4.5;
422  double Lambda2 = Lambda3flav2;
423  double pT2minNow = pT2endDip;
424  int idMother = 0;
425  int idSister = 0;
426  double z = 0.;
427  double zMaxAbs = 0.;
428  double zRootMax = 0.;
429  double zRootMin = 0.;
430  double g2gInt = 0.;
431  double q2gInt = 0.;
432  double q2qInt = 0.;
433  double g2qInt = 0.;
434  double g2Qenhance = 0.;
435  double xPDFdaughter = 0.;
436  double xPDFmother[21] = {0.};
437  double xPDFgMother = 0.;
438  double xPDFmotherSum = 0.;
439  double kernelPDF = 0.;
440  double xMother = 0.;
441  double wt = 0.;
442  double Q2 = 0.;
443  double mSister = 0.;
444  double m2Sister = 0.;
445  double pT2corr = 0.;
446  double pT2PDF = pT2;
447  bool needNewPDF = true;
448 
449  // Begin evolution loop towards smaller pT values.
450  int loopTinyPDFdau = 0;
451  bool hasTinyPDFdau = false;
452  do {
453  wt = 0.;
454 
455  // Bad sign if repeated looping with small daughter PDF, so fail.
456  // (Example: if all PDF's = 0 below Q_0, except for c/b companion.)
457  if (hasTinyPDFdau) ++loopTinyPDFdau;
458  if (loopTinyPDFdau > MAXLOOPTINYPDF) {
459  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQCD: "
460  "small daughter PDF");
461  return;
462  }
463 
464  // Initialize integrals of splitting kernels and evaluate parton
465  // densities at the beginning. Reinitialize after long evolution
466  // in pT2 or when crossing c and b flavour thresholds.
467  if (needNewPDF || pT2 < EVALPDFSTEP * pT2PDF) {
468  pT2PDF = pT2;
469  hasTinyPDFdau = false;
470 
471  // Determine overestimated z range; switch at c and b masses.
472  if (pT2 > m2b) {
473  nFlavour = 5;
474  pT2minNow = max( m2b, pT2endDip);
475  b0 = 23./6.;
476  Lambda2 = Lambda5flav2;
477  } else if (pT2 > m2c) {
478  nFlavour = 4;
479  pT2minNow = max( m2c, pT2endDip);
480  b0 = 25./6.;
481  Lambda2 = Lambda4flav2;
482  } else {
483  nFlavour = 3;
484  pT2minNow = pT2endDip;
485  b0 = 27./6.;
486  Lambda2 = Lambda3flav2;
487  }
488  zMaxAbs = 1. - 0.5 * (pT2minNow / m2Dip) *
489  ( sqrt( 1. + 4. * m2Dip / pT2minNow ) - 1. );
490  if (isMassive) zMaxAbs = min( zMaxAbs, zMaxMassive);
491 
492  // Go to another z range with lower mass scale if current is closed.
493  if (zMinAbs > zMaxAbs) {
494  if (nFlavour == 3 || (idMassive == 4 && nFlavour == 4)
495  || idMassive == 5) return;
496  pT2 = (nFlavour == 4) ? m2c : m2b;
497  continue;
498  }
499 
500  // Parton density of daughter at current scale.
501  xPDFdaughter = beam.xfISR(iSysNow, idDaughter, xDaughter, pT2);
502  if (xPDFdaughter < TINYPDF) {
503  xPDFdaughter = TINYPDF;
504  hasTinyPDFdau = true;
505  }
506 
507  // Integrals of splitting kernels for gluons: g -> g, q -> g.
508  if (isGluon) {
509  g2gInt = 6. * log(zMaxAbs * (1.-zMinAbs)
510  / (zMinAbs * (1.-zMaxAbs)));
511  if (doMEcorrections) g2gInt *= calcMEmax(MEtype, 21, 21);
512  q2gInt = (16./3.) * (1./sqrt(zMinAbs) - 1./sqrt(zMaxAbs));
513  if (doMEcorrections) q2gInt *= calcMEmax(MEtype, 1, 21);
514 
515  // Parton density of potential quark mothers to a g.
516  xPDFmotherSum = 0.;
517  for (int i = -nQuarkIn; i <= nQuarkIn; ++i) {
518  if (i == 0) {
519  xPDFmother[10] = 0.;
520  } else {
521  xPDFmother[i+10] = beam.xfISR(iSysNow, i, xDaughter, pT2);
522  xPDFmotherSum += xPDFmother[i+10];
523  }
524  }
525 
526  // Total QCD evolution coefficient for a gluon.
527  kernelPDF = g2gInt + q2gInt * xPDFmotherSum / xPDFdaughter;
528 
529  // For valence quark only need consider q -> q g branchings.
530  // Introduce an extra factor sqrt(z) to smooth bumps.
531  } else if (isValence) {
532  zRootMin = (1. + sqrt(zMinAbs)) / (1. - sqrt(zMinAbs));
533  zRootMax = (1. + sqrt(zMaxAbs)) / (1. - sqrt(zMaxAbs));
534  q2qInt = (8./3.) * log( zRootMax / zRootMin );
535  if (doMEcorrections) q2qInt *= calcMEmax(MEtype, 1, 1);
536  kernelPDF = q2qInt;
537 
538  // Integrals of splitting kernels for quarks: q -> q, g -> q.
539  } else {
540  q2qInt = (8./3.) * log( (1. - zMinAbs) / (1. - zMaxAbs) );
541  if (doMEcorrections) q2qInt *= calcMEmax(MEtype, 1, 1);
542  g2qInt = 0.5 * (zMaxAbs - zMinAbs);
543  if (doMEcorrections) g2qInt *= calcMEmax(MEtype, 21, 1);
544 
545  // Increase estimated upper weight for g -> Q + Qbar.
546  if (isMassive) {
547  if (alphaSorder == 0) g2Qenhance = log(pT2/m2Massive)
548  / log(m2Threshold/m2Massive);
549  else {
550  double m2log = log( m2Massive / Lambda2);
551  g2Qenhance = log( log(pT2/Lambda2) / m2log )
552  / log( log(m2Threshold/Lambda2) / m2log );
553  }
554  g2qInt *= g2Qenhance;
555  }
556 
557  // Parton density of a potential gluon mother to a q.
558  xPDFgMother = beam.xfISR(iSysNow, 21, xDaughter, pT2);
559 
560  // Total QCD evolution coefficient for a quark.
561  kernelPDF = q2qInt + g2qInt * xPDFgMother / xPDFdaughter;
562  }
563 
564  // End evaluation of splitting kernels and parton densities.
565  needNewPDF = false;
566  }
567  if (kernelPDF < TINYKERNELPDF) return;
568 
569  // Pick pT2 (in overestimated z range), for one of three different cases.
570  // Assume form alphas(pT0^2 + pT^2) * dpT^2/(pT0^2 + pT^2).
571  double Q2alphaS;
572 
573  // Fixed alpha_strong.
574  if (alphaSorder == 0) {
575  pT2 = (pT2 + pT20) * pow( rndmPtr->flat(),
576  1. / (alphaS2pi * kernelPDF)) - pT20;
577 
578  // First-order alpha_strong.
579  } else if (alphaSorder == 1) {
580  pT2 = Lambda2 * pow( (pT2 + pT20) / Lambda2,
581  pow(rndmPtr->flat(), b0 / kernelPDF) ) - pT20;
582 
583  // For second order reject by second term in alpha_strong expression.
584  } else {
585  do {
586  pT2 = Lambda2 * pow( (pT2 + pT20) / Lambda2,
587  pow(rndmPtr->flat(), b0 / kernelPDF) ) - pT20;
588  Q2alphaS = max(pT2 + pT20, pow2(LAMBDA3MARGIN) * Lambda3flav2);
589  } while (alphaS.alphaS2OrdCorr(Q2alphaS) < rndmPtr->flat()
590  && pT2 > pT2minNow);
591  }
592 
593  // Check for pT2 values that prompt special action.
594 
595  // If fallen into b threshold region, force g -> b + bbar.
596  if (idMassive == 5 && pT2 < m2Threshold) {
597  pT2nearQCDthreshold( beam, m2Massive, m2Threshold, xMaxAbs,
598  zMinAbs, zMaxMassive );
599  return;
600 
601  // If crossed b threshold, continue evolution from this threshold.
602  } else if (nFlavour == 5 && pT2 < m2b) {
603  needNewPDF = true;
604  pT2 = m2b;
605  continue;
606 
607  // If fallen into c threshold region, force g -> c + cbar.
608  } else if (idMassive == 4 && pT2 < m2Threshold) {
609  pT2nearQCDthreshold( beam, m2Massive, m2Threshold, xMaxAbs,
610  zMinAbs, zMaxMassive );
611  return;
612 
613  // If crossed c threshold, continue evolution from this threshold.
614  } else if (nFlavour == 4 && pT2 < m2c) {
615  needNewPDF = true;
616  pT2 = m2c;
617  continue;
618 
619  // Abort evolution if below cutoff scale, or below another branching.
620  } else if (pT2 < pT2endDip) return;
621 
622  // Select z value of branching to g, and corrective weight.
623  if (isGluon) {
624  // g -> g (+ g).
625  if (rndmPtr->flat() * kernelPDF < g2gInt) {
626  idMother = 21;
627  idSister = 21;
628  z = 1. / ( 1. + ((1. - zMinAbs) / zMinAbs) * pow( (zMinAbs *
629  (1. - zMaxAbs)) / (zMaxAbs * (1. - zMinAbs)), rndmPtr->flat() ) );
630  wt = pow2( 1. - z * (1. - z));
631  } else {
632  // q -> g (+ q): also select flavour.
633  double temp = xPDFmotherSum * rndmPtr->flat();
634  idMother = -nQuarkIn - 1;
635  do { temp -= xPDFmother[(++idMother) + 10]; }
636  while (temp > 0. && idMother < nQuarkIn);
637  idSister = idMother;
638  z = (zMinAbs * zMaxAbs) / pow2( sqrt(zMinAbs) + rndmPtr->flat()
639  * ( sqrt(zMaxAbs)- sqrt(zMinAbs) ));
640  wt = 0.5 * (1. + pow2(1. - z)) * sqrt(z)
641  * xPDFdaughter / xPDFmother[idMother + 10];
642  }
643 
644  // Select z value of branching to q, and corrective weight.
645  // Include massive kernel corrections for c and b quarks.
646  } else {
647  // q -> q (+ g).
648  if (isValence || rndmPtr->flat() * kernelPDF < q2qInt) {
649  idMother = idDaughter;
650  idSister = 21;
651  // Valence more peaked at large z.
652  if (isValence) {
653  double zTmp = zRootMin * pow(zRootMax / zRootMin, rndmPtr->flat() );
654  z = pow2( (1. - zTmp) / (1. + zTmp) );
655  } else {
656  z = 1. - (1. - zMinAbs) * pow( (1. - zMaxAbs) / (1. - zMinAbs),
657  rndmPtr->flat() );
658  }
659  if (!isMassive) {
660  wt = 0.5 * (1. + pow2(z));
661  } else {
662  //?? Bug?? should be 2 more for massive part??
663  // wt = 0.5 * (1. + pow2(z) - z * pow2(1.-z) * m2Massive / pT2);
664  wt = 0.5 * (1. + pow2(z)) - z * pow2(1.-z) * m2Massive / pT2;
665  }
666  if (isValence) wt *= sqrt(z);
667  // g -> q (+ qbar).
668  } else {
669  idMother = 21;
670  idSister = - idDaughter;
671  z = zMinAbs + rndmPtr->flat() * (zMaxAbs - zMinAbs);
672  if (!isMassive) {
673  wt = (pow2(z) + pow2(1.-z)) * xPDFdaughter / xPDFgMother ;
674  } else {
675  wt = (pow2(z) + pow2(1.-z) + 2. * z * (1.-z) * m2Massive / pT2)
676  * xPDFdaughter / (xPDFgMother * g2Qenhance) ;
677  }
678  }
679  }
680 
681  // Derive Q2 and x of mother from pT2 and z.
682  Q2 = pT2 / (1. - z);
683  xMother = xDaughter / z;
684  // Correction to x for massive recoiler from rescattering.
685  if (!dipEndNow->normalRecoil) {
686  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
687  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
688  }
689  if(xMother > xMaxAbs) { wt = 0.; continue; }
690 
691  // Forbidden emission if outside allowed z range for given pT2.
692  mSister = particleDataPtr->m0(idSister);
693  m2Sister = pow2(mSister);
694  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
695  if(pT2corr < TINYPT2) { wt = 0.; continue; }
696 
697  // Optionally veto emissions not ordered in rapidity (= angle).
698  if ( doRapidityOrder && dipEndNow->nBranch > 0
699  && pT2 > pow2( (1. - z) / (z * (1. - dipEndNow->zOld)) )
700  * dipEndNow->pT2Old ) { wt = 0.; continue; }
701 
702  // If creating heavy quark by Q -> g + Q then next need g -> Q + Qbar.
703  // So minimum total mass2 is 4 * m2Sister, but use more to be safe.
704  if ( isGluon && ( abs(idMother) == 4 || abs(idMother) == 5 )) {
705  double m2QQsister = EXTRASPACEQ * 4. * m2Sister;
706  double pT2QQcorr = Q2 - z * (m2Dip + Q2) * (Q2 + m2QQsister) / m2Dip;
707  if(pT2QQcorr < TINYPT2) { wt = 0.; continue; }
708  }
709 
710  // Evaluation of ME correction.
711  if (doMEcorrections) wt *= calcMEcorr(MEtype, idMother, idDaughter,
712  m2Dip, z, Q2) / calcMEmax(MEtype, idMother, idDaughter);
713 
714  // Optional dampening of large pT values in first radiation.
715  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
716  wt *= pT2damp / (pT2 + pT2damp);
717 
718  // Idea suggested by Gosta Gustafson: increased screening in events
719  // with large activity can be simulated by pT0_eff = sqrt(n) * pT0.
720  if (enhanceScreening == 2) {
721  int nSysNow = infoPtr->nMPI() + infoPtr->nISR() + 1;
722  double WTscreen = pow2( (pT2 + pT20) / (pT2 + nSysNow * pT20) );
723  wt *= WTscreen;
724  }
725 
726  // Evaluation of new daughter and mother PDF's.
727  double xPDFdaughterNew = max ( TINYPDF,
728  beam.xfISR(iSysNow, idDaughter, xDaughter, pT2) );
729  double xPDFmotherNew = beam.xfISR(iSysNow, idMother, xMother, pT2);
730  wt *= xPDFmotherNew / xPDFdaughterNew;
731 
732  // Check that valence step does not cause problem.
733  if (wt > 1.) infoPtr->errorMsg("Warning in SpaceShower::"
734  "pT2nextQCD: weight above unity");
735 
736  // Iterate until acceptable pT (or have fallen below pTmin).
737  } while (wt < rndmPtr->flat()) ;
738 
739  // Save values for (so far) acceptable branching.
740  dipEndNow->store( idDaughter,idMother, idSister, x1Now, x2Now, m2Dip,
741  pT2, z, xMother, Q2, mSister, m2Sister, pT2corr);
742 
743 }
744 
745 //--------------------------------------------------------------------------
746 
747 // Evolve a QCD dipole end near threshold, with g -> Q + Qbar enforced.
748 // Note: No explicit Sudakov factor formalism here. Instead use that
749 // df_Q(x, pT2) = (alpha_s/2pi) * (dT2/pT2) * ((gluon) * (splitting)).
750 // This implies that effects of Q -> Q + g are neglected in this range.
751 
752 void SpaceShower::pT2nearQCDthreshold( BeamParticle& beam,
753  double m2Massive, double m2Threshold, double xMaxAbs,
754  double zMinAbs, double zMaxMassive) {
755 
756  // Initial values, to be used in kinematics and weighting.
757  double Lambda2 = (abs(idDaughter) == 4) ? Lambda4flav2 : Lambda5flav2;
758  double logM2Lambda2 = log( m2Massive / Lambda2 );
759  double xPDFmotherOld = beam.xfISR(iSysNow, 21, xDaughter, m2Threshold);
760 
761  // Variables used inside evolution loop. (Mainly dummy start values.)
762  int loop = 0;
763  double wt = 0.;
764  double pT2 = 0.;
765  double z = 0.;
766  double Q2 = 0.;
767  double pT2corr = 0.;
768  double xMother = 0.;
769 
770  // Begin loop over tries to find acceptable g -> Q + Qbar branching.
771  do {
772  wt = 0.;
773 
774  // Check that not caught in infinite loop with impossible kinematics.
775  if (++loop > 100) {
776  infoPtr->errorMsg("Error in SpaceShower::pT2nearQCDthreshold: "
777  "stuck in loop");
778  return;
779  }
780 
781  // Pick dpT2/pT2 in range [m2Massive,thresholdRatio * m2Massive].
782  pT2 = m2Massive * pow( m2Threshold / m2Massive, rndmPtr->flat() );
783 
784  // Pick z flat in allowed range.
785  z = zMinAbs + rndmPtr->flat() * (zMaxMassive - zMinAbs);
786 
787  // Check that kinematically possible choice.
788  Q2 = pT2 / (1.-z) - m2Massive;
789  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Massive) / m2Dip;
790  if(pT2corr < TINYPT2) continue;
791 
792  // Correction factor for running alpha_s. ??
793  wt = logM2Lambda2 / log( pT2 / Lambda2 );
794 
795  // Correction factor for splitting kernel.
796  wt *= pow2(z) + pow2(1.-z) + 2. * z * (1.-z) * m2Massive / pT2;
797 
798  // x, including correction for massive recoiler from rescattering.
799  xMother = xDaughter / z;
800  if (!dipEndNow->normalRecoil) {
801  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
802  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
803  }
804  if (xMother > xMaxAbs) { wt = 0.; continue; }
805 
806  // Correction factor for gluon density.
807  double xPDFmotherNew = beam.xfISR(iSysNow, 21, xMother, pT2);
808  wt *= xPDFmotherNew / xPDFmotherOld;
809 
810  // Iterate until acceptable pT and z.
811  } while (wt < rndmPtr->flat()) ;
812 
813  // Save values for (so far) acceptable branching.
814  double mSister = (abs(idDaughter) == 4) ? mc : mb;
815  dipEndNow->store( idDaughter, 21, -idDaughter, x1Now, x2Now, m2Dip,
816  pT2, z, xMother, Q2, mSister, pow2(mSister), pT2corr);
817 
818 }
819 
820 //--------------------------------------------------------------------------
821 
822 // Evolve a QED dipole end.
823 
824 void SpaceShower::pT2nextQED( double pT2begDip, double pT2endDip) {
825 
826  // Type of dipole and starting values.
827  BeamParticle& beam = (sideA) ? *beamAPtr : *beamBPtr;
828  bool isLeptonBeam = beam.isLepton();
829  int MEtype = dipEndNow->MEtype;
830  bool isPhoton = (idDaughter == 22);
831  double pT2 = pT2begDip;
832  double m2Lepton = (isLeptonBeam) ? pow2(beam.m()) : 0.;
833  if (isLeptonBeam && pT2begDip < m2Lepton) return;
834 
835  // Currently no f -> gamma branching implemented for lepton beams
836  if (isPhoton && isLeptonBeam) return;
837 
838  // alpha_em at maximum scale provides upper estimate.
839  double alphaEMmax = alphaEM.alphaEM(pT2begDip);
840  double alphaEM2pi = alphaEMmax / (2. * M_PI);
841 
842  // Maximum x of mother implies minimum z = xDaughter / xMother.
843  double xMaxAbs = (isLeptonBeam) ? LEPTONXMAX : beam.xMax(iSysNow);
844  double zMinAbs = xDaughter / xMaxAbs;
845  if (xMaxAbs < 0.) {
846  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQED: "
847  "xMaxAbs negative");
848  return;
849  }
850 
851  // Maximum z from minimum pT and, for lepton, from minimum x_gamma.
852  double zMaxAbs = 1. - 0.5 * (pT2endDip / m2Dip) *
853  ( sqrt( 1. + 4. * m2Dip / pT2endDip ) - 1. );
854  if (isLeptonBeam) {
855  double zMaxLepton = xDaughter / (xDaughter + LEPTONXMIN);
856  if (zMaxLepton < zMaxAbs) zMaxAbs = zMaxLepton;
857  }
858  if (zMaxAbs < zMinAbs) return;
859 
860  // Variables used inside evolution loop. (Mainly dummy start values.)
861  int idMother = 0;
862  int idSister =22;
863  double z = 0.;
864  double xMother = 0.;
865  double wt = 0.;
866  double Q2 = 0.;
867  double mSister = 0.;
868  double m2Sister = 0.;
869  double pT2corr = 0.;
870 
871  // QED evolution of fermions
872  if (!isPhoton) {
873 
874  // Integrals of splitting kernels for fermions: f -> f. Use 1 + z^2 < 2.
875  // Ansatz f(z) = 2 / (1 - z), with + 2 / (z - xDaughter) for lepton.
876  double f2fInt = 0.;
877  double f2fIntA = 2. * log( (1. - zMinAbs) / (1. - zMaxAbs) );
878  double f2fIntB = 0.;
879  if (isLeptonBeam) {
880  f2fIntB = 2. * log( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter) );
881  f2fInt = f2fIntA + f2fIntB;
882  } else f2fInt = pow2(dipEndNow->chgType / 3.) * f2fIntA;
883 
884  // Upper estimate for evolution equation, including fudge factor.
885  if (doMEcorrections) f2fInt *= calcMEmax(MEtype, 1, 1);
886  double kernelPDF = alphaEM2pi * f2fInt;
887  double fudge = (isLeptonBeam) ? LEPTONFUDGE * log(m2Dip/m2Lepton) : 1.;
888  kernelPDF *= fudge;
889  if (kernelPDF < TINYKERNELPDF) return;
890 
891  // Begin evolution loop towards smaller pT values.
892  do {
893 
894  // Pick pT2 (in overestimated z range).
895  // For l -> l gamma include extrafactor 1 / ln(pT2 / m2l) in evolution.
896  double shift = pow(rndmPtr->flat(), 1. / kernelPDF);
897  if (isLeptonBeam) pT2 = m2Lepton * pow( pT2 / m2Lepton, shift);
898  else pT2 = pT2 * shift;
899 
900  // Abort evolution if below cutoff scale, or below another branching.
901  if (pT2 < pT2endDip) return;
902  if (isLeptonBeam && pT2 < LEPTONPT2MIN * m2Lepton) return;
903 
904  // Select z value of branching f -> f + gamma, and corrective weight.
905  idMother = idDaughter;
906  wt = 0.5 * (1. + pow2(z));
907  if (isLeptonBeam) {
908  if (f2fIntA > rndmPtr->flat() * (f2fIntA + f2fIntB)) {
909  z = 1. - (1. - zMinAbs)
910  * pow( (1. - zMaxAbs) / (1. - zMinAbs), rndmPtr->flat() );
911  } else {
912  z = xDaughter + (zMinAbs - xDaughter)
913  * pow( (zMaxAbs - xDaughter) / (zMinAbs - xDaughter),
914  rndmPtr->flat() );
915  }
916  wt *= (z - xDaughter) / (1. - xDaughter);
917  } else {
918  z = 1. - (1. - zMinAbs)
919  * pow( (1. - zMaxAbs) / (1. - zMinAbs), rndmPtr->flat() );
920  }
921 
922  // Derive Q2 and x of mother from pT2 and z.
923  Q2 = pT2 / (1. - z);
924  xMother = xDaughter / z;
925  // Correction to x for massive recoiler from rescattering.
926  if (!dipEndNow->normalRecoil) {
927  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
928  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
929  }
930  if(xMother > xMaxAbs) { wt = 0.; continue; }
931 
932  // Forbidden emission if outside allowed z range for given pT2.
933  mSister = 0.;
934  m2Sister = 0.;
935  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
936  if(pT2corr < TINYPT2) { wt = 0.; continue; }
937 
938  // Correct by ln(pT2 / m2l) and fudge factor.
939  if (isLeptonBeam) wt *= log(pT2 / m2Lepton) / fudge;
940 
941  // Evaluation of ME correction.
942  if (doMEcorrections) wt *= calcMEcorr(MEtype, idMother, idDaughter,
943  m2Dip, z, Q2) / calcMEmax(MEtype, idMother, idDaughter);
944 
945  // Extra QED correction for f fbar -> W+- gamma. Debug??
946  if (doMEcorrections && MEtype == 1 && idDaughter == idMother
947  && ( (iSysNow == 0 && idResFirst == 24)
948  || (iSysNow == 1 && idResSecond == 24) ) ) {
949  double tHe = -Q2;
950  double uHe = Q2 - m2Dip * (1. - z) / z;
951  double chg1 = abs(dipEndNow->chgType / 3.);
952  double chg2 = 1. - chg1;
953  wt *= pow2(chg1 * uHe - chg2 * tHe)
954  / ( (tHe + uHe) * (pow2(chg1) * uHe + pow2(chg2) * tHe) );
955  }
956 
957  // Optional dampening of large pT values in first radiation.
958  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
959  wt *= pT2damp / (pT2 + pT2damp);
960 
961  // Correct to current value of alpha_EM.
962  double alphaEMnow = alphaEM.alphaEM(pT2);
963  wt *= (alphaEMnow / alphaEMmax);
964 
965  // Evaluation of new daughter and mother PDF's.
966  double xPDFdaughterNew = max ( TINYPDF,
967  beam.xfISR(iSysNow, idDaughter, xDaughter, pT2) );
968  double xPDFmotherNew = beam.xfISR(iSysNow, idMother, xMother, pT2);
969  wt *= xPDFmotherNew / xPDFdaughterNew;
970 
971  // Iterate until acceptable pT (or have fallen below pTmin).
972  } while (wt < rndmPtr->flat()) ;
973  }
974 
975  // QED evolution of photons (so far only for hadron beams).
976  else {
977 
978  // Initial values
979  int nFlavour = 3;
980  double kernelPDF = 0.0;
981  double xPDFdaughter = 0.;
982  double xPDFmother[21] = {0.};
983  double xPDFmotherSum = 0.0;
984  double pT2PDF = pT2;
985  double pT2minNow = pT2endDip;
986  bool needNewPDF = true;
987 
988  // Begin evolution loop towards smaller pT values.
989  int loopTinyPDFdau = 0;
990  bool hasTinyPDFdau = false;
991  do {
992  wt = 0.;
993 
994  // Bad sign if repeated looping with small daughter PDF, so fail.
995  if (hasTinyPDFdau) ++loopTinyPDFdau;
996  if (loopTinyPDFdau > MAXLOOPTINYPDF) {
997  infoPtr->errorMsg("Warning in SpaceShower::pT2nextQED: "
998  "small daughter PDF");
999  return;
1000  }
1001 
1002  // Initialize integrals of splitting kernels and evaluate parton
1003  // densities at the beginning. Reinitialize after long evolution
1004  // in pT2 or when crossing c and b flavour thresholds.
1005  if (needNewPDF || pT2 < EVALPDFSTEP * pT2PDF) {
1006 
1007  pT2PDF = pT2;
1008  hasTinyPDFdau = false;
1009 
1010  // Determine overestimated z range; switch at c and b masses.
1011  if (pT2 > m2b && nQuarkIn >= 5) {
1012  nFlavour = 5;
1013  pT2minNow = max( m2b, pT2endDip);
1014  } else if (pT2 > m2c && nQuarkIn >= 4) {
1015  nFlavour = 4;
1016  pT2minNow = max( m2c, pT2endDip);
1017  } else {
1018  nFlavour = 3;
1019  pT2minNow = pT2endDip;
1020  }
1021 
1022  // Compute upper z limit
1023  zMaxAbs = 1. - 0.5 * (pT2minNow / m2Dip) *
1024  ( sqrt( 1. + 4. * m2Dip / pT2minNow ) - 1. );
1025 
1026  // Parton density of daughter at current scale.
1027  xPDFdaughter = beam.xfISR(iSysNow, idDaughter, xDaughter, pT2);
1028  if (xPDFdaughter < TINYPDF) {
1029  xPDFdaughter = TINYPDF;
1030  hasTinyPDFdau = true;
1031  }
1032 
1033  // Integral over f -> gamma f splitting kernel.
1034  // Normalized so: 4/3 aS/2pi P(z) -> eq^2 * aEM/2pi P(z).
1035  // (Charge-weighting happens below.)
1036  double q2gInt = 4. * (1./sqrt(zMinAbs) - 1./sqrt(zMaxAbs));
1037 
1038  // Charge-weighted Parton density of potential quark mothers.
1039  xPDFmotherSum = 0.;
1040  for (int i = -nFlavour; i <= nFlavour; ++i) {
1041  if (i == 0) {
1042  xPDFmother[10] = 0.;
1043  } else {
1044  xPDFmother[i+10] = pow2((abs(i+1) % 2 + 1)/3.0)
1045  * beam.xfISR(iSysNow, i, xDaughter, pT2);
1046  xPDFmotherSum += xPDFmother[i+10];
1047  }
1048  }
1049 
1050  // Total QED evolution coefficient for a photon.
1051  kernelPDF = q2gInt * xPDFmotherSum / xPDFdaughter;
1052 
1053  // End evaluation of splitting kernels and parton densities.
1054  needNewPDF = false;
1055  }
1056  if (kernelPDF < TINYKERNELPDF) return;
1057 
1058  // Select pT2 for next trial branching
1059  pT2 *= pow( rndmPtr->flat(), 1. / (alphaEM2pi * kernelPDF));
1060 
1061  // If crossed b threshold, continue evolution from this threshold.
1062  if (nFlavour == 5 && pT2 < m2b) {
1063  needNewPDF = true;
1064  pT2 = m2b;
1065  continue;
1066  }
1067 
1068  // If crossed c threshold, continue evolution from this threshold.
1069  else if (nFlavour == 4 && pT2 < m2c) {
1070  needNewPDF = true;
1071  pT2 = m2c;
1072  continue;
1073  }
1074 
1075  // Abort evolution if below cutoff scale, or below another branching.
1076  else if (pT2 < pT2endDip) return;
1077 
1078  // Select flavour for trial branching
1079  double temp = xPDFmotherSum * rndmPtr->flat();
1080  idMother = -nQuarkIn - 1;
1081  do {
1082  temp -= xPDFmother[(++idMother) + 10];
1083  } while (temp > 0. && idMother < nQuarkIn);
1084 
1085  // Sister is same as mother, but can have m2 > 0
1086  idSister = idMother;
1087  mSister = particleDataPtr->m0(idSister);
1088  m2Sister = pow2(mSister);
1089 
1090  // Select z for trial branching
1091  z = (zMinAbs * zMaxAbs) / pow2( sqrt(zMinAbs) + rndmPtr->flat()
1092  * ( sqrt(zMaxAbs)- sqrt(zMinAbs) ));
1093 
1094  // Trial weight: splitting kernel
1095  wt = 0.5 * (1. + pow2(1. - z)) * sqrt(z);
1096 
1097  // Trial weight: running alpha_EM
1098  double alphaEMnow = alphaEM.alphaEM(pT2);
1099  wt *= (alphaEMnow / alphaEMmax);
1100 
1101  // Derive Q2 and x of mother from pT2 and z
1102  Q2 = pT2 / (1. - z);
1103  xMother = xDaughter / z;
1104  // Correction to x for massive recoiler from rescattering.
1105  if (!dipEndNow->normalRecoil) {
1106  if (sideA) xMother += (m2Rec / (x2Now * sCM)) * (1. / z - 1.);
1107  else xMother += (m2Rec / (x1Now * sCM)) * (1. / z - 1.);
1108  }
1109 
1110  // Compute pT2corr
1111  pT2corr = Q2 - z * (m2Dip + Q2) * (Q2 + m2Sister) / m2Dip;
1112  if(pT2corr < TINYPT2) { wt = 0.; continue; }
1113 
1114  // If creating heavy quark by Q -> gamma + Q then next g -> Q + Qbar.
1115  // So minimum total mass2 is 4 * m2Sister, but use more to be safe.
1116  if ( abs(idMother) == 4 || abs(idMother) == 5 ) {
1117  double m2QQsister = EXTRASPACEQ * 4. * m2Sister;
1118  double pT2QQcorr = Q2 - z * (m2Dip + Q2) * (Q2 + m2QQsister) / m2Dip;
1119  if(pT2QQcorr < TINYPT2) { wt = 0.; continue; }
1120  }
1121 
1122  // Optional dampening of large pT values in first radiation.
1123  if (dopTdamp && iSysNow == 0 && MEtype == 0 && nRad == 0)
1124  wt *= pT2damp / (pT2 + pT2damp);
1125 
1126  // Evaluation of new daughter PDF
1127  double xPDFdaughterNew = beam.xfISR(iSysNow, idDaughter, xDaughter, pT2);
1128  if (xPDFdaughterNew < TINYPDF) {
1129  xPDFdaughterNew = TINYPDF;
1130  }
1131 
1132  // Evaluation of new charge-weighted mother PDF
1133  double xPDFmotherNew = pow2( (abs(idMother+1) % 2 + 1)/3.0 )
1134  * beam.xfISR(iSysNow, idMother, xMother, pT2);
1135 
1136  // Trial weight: divide out old pdf ratio
1137  wt *= xPDFdaughter / xPDFmother[idMother + 10];
1138 
1139  // Trial weight: new pdf ratio
1140  wt *= xPDFmotherNew / xPDFdaughterNew;
1141 
1142  // Debug information.
1143  if (DEBUG) cout << " Trial weight is : " << wt << "\n pT2 = "
1144  << pT2 << " z = " << z << " alpha/alphahat = " << alphaEMnow
1145  << "/" << alphaEMmax << " idMother = " << idMother
1146  << " xPDFd/xPDFm = " << xPDFdaughter << "/" << xPDFmother[idMother+10]
1147  << " xPDFmNew/xPDFdNew " << xPDFmotherNew << "/" << xPDFdaughterNew
1148  << " pT2damp = " << pT2damp << " Q2 = " << Q2 << " pT2corr = "
1149  << pT2corr << endl;
1150 
1151  // Iterate until acceptable pT (or have fallen below pTmin).
1152  } while (wt < rndmPtr->flat()) ;
1153  }
1154 
1155  // Save values for (so far) acceptable branching.
1156  dipEndNow->store( idDaughter, idMother, idSister, x1Now, x2Now, m2Dip,
1157  pT2, z, xMother, Q2, mSister, m2Sister, pT2corr);
1158 
1159 }
1160 
1161 //--------------------------------------------------------------------------
1162 
1163 // Kinematics of branching.
1164 // Construct mother -> daughter + sister, with recoiler on other side.
1165 
1166 bool SpaceShower::branch( Event& event) {
1167 
1168  // Side on which branching occured.
1169  int side = abs(dipEndSel->side);
1170  double sideSign = (side == 1) ? 1. : -1.;
1171 
1172  // Read in flavour and colour variables.
1173  int iDaughter = partonSystemsPtr->getInA(iSysSel);
1174  int iRecoiler = partonSystemsPtr->getInB(iSysSel);
1175  if (side == 2) swap(iDaughter, iRecoiler);
1176  int idDaughterNow = dipEndSel->idDaughter;
1177  int idMother = dipEndSel->idMother;
1178  int idSister = dipEndSel->idSister;
1179  int colDaughter = event[iDaughter].col();
1180  int acolDaughter = event[iDaughter].acol();
1181 
1182  // Recoil parton may be rescatterer, requiring special processing.
1183  bool normalRecoil = dipEndSel->normalRecoil;
1184  int iRecoilMother = event[iRecoiler].mother1();
1185 
1186  // Read in kinematical variables.
1187  double x1 = dipEndSel->x1;
1188  double x2 = dipEndSel->x2;
1189  double xMo = dipEndSel->xMo;
1190  double m2 = dipEndSel->m2Dip;
1191  double m = sqrt(m2);
1192  double pT2 = dipEndSel->pT2;
1193  double z = dipEndSel->z;
1194  double Q2 = dipEndSel->Q2;
1195  double mSister = dipEndSel->mSister;
1196  double m2Sister = dipEndSel->m2Sister;
1197  double pT2corr = dipEndSel->pT2corr;
1198  double x1New = (side == 1) ? xMo : x1;
1199  double x2New = (side == 2) ? xMo : x2;
1200 
1201  // Rescatter: kinematics may fail; use the rescatterFail flag to tell
1202  // parton level to try again.
1203  rescatterFail = false;
1204 
1205  // Construct kinematics of mother, sister and recoiler in old rest frame.
1206  // Normally both mother and recoiler are taken massless.
1207  double eNewRec, pzNewRec, pTbranch, pzMother;
1208  if (normalRecoil) {
1209  eNewRec = 0.5 * (m2 + Q2) / m;
1210  pzNewRec = -sideSign * eNewRec;
1211  pTbranch = sqrt(pT2corr) * m2 / ( z * (m2 + Q2) );
1212  pzMother = sideSign * 0.5 * m * ( (m2 - Q2) / ( z * (m2 + Q2) )
1213  + (Q2 + m2Sister) / m2 );
1214  // More complicated kinematics when recoiler not massless. May fail.
1215  } else {
1216  m2Rec = event[iRecoiler].m2();
1217  double s1Tmp = m2 + Q2 - m2Rec;
1218  double s3Tmp = m2 / z - m2Rec;
1219  double r1Tmp = sqrt(s1Tmp * s1Tmp + 4. * Q2 * m2Rec);
1220  eNewRec = 0.5 * (m2 + m2Rec + Q2) / m;
1221  pzNewRec = -sideSign * 0.5 * r1Tmp / m;
1222  double pT2br = Q2 * s3Tmp * (m2 / z - m2 - Q2)
1223  - m2Sister * s1Tmp * s3Tmp - m2Rec * pow2(Q2 + m2Sister);
1224  if (pT2br <= 0.) return false;
1225  pTbranch = sqrt(pT2br) / r1Tmp;
1226  pzMother = sideSign * (m * s3Tmp
1227  - eNewRec * (m2 / z - Q2 - m2Rec - m2Sister)) / r1Tmp;
1228  }
1229  // Common final kinematics steps for both normal and rescattering.
1230  double eMother = sqrt( pow2(pTbranch) + pow2(pzMother) );
1231  double pzSister = pzMother + pzNewRec;
1232  double eSister = sqrt( pow2(pTbranch) + pow2(pzSister) + m2Sister );
1233  Vec4 pMother( pTbranch, 0., pzMother, eMother );
1234  Vec4 pSister( pTbranch, 0., pzSister, eSister );
1235  Vec4 pNewRec( 0., 0., pzNewRec, eNewRec );
1236 
1237  // Current event and subsystem size.
1238  int eventSizeOld = event.size();
1239  int systemSizeOld = partonSystemsPtr->sizeAll(iSysSel);
1240 
1241  // Save properties to be restored in case of user-hook veto of emission.
1242  int ev1dau1V = event[1].daughter1();
1243  int ev2dau1V = event[2].daughter1();
1244  vector<int> statusV, mother1V, mother2V, daughter1V, daughter2V;
1245  if (canVetoEmission) {
1246  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1247  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1248  statusV.push_back( event[iOldCopy].status());
1249  mother1V.push_back( event[iOldCopy].mother1());
1250  mother2V.push_back( event[iOldCopy].mother2());
1251  daughter1V.push_back( event[iOldCopy].daughter1());
1252  daughter2V.push_back( event[iOldCopy].daughter2());
1253  }
1254  }
1255 
1256  // Take copy of existing system, to be given modified kinematics.
1257  // Incoming negative status. Rescattered also negative, but after copy.
1258  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1259  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1260  int statusOld = event[iOldCopy].status();
1261  int statusNew = (iOldCopy == iDaughter
1262  || iOldCopy == iRecoiler) ? statusOld : 44;
1263  int iNewCopy = event.copy(iOldCopy, statusNew);
1264  if (statusOld < 0) event[iNewCopy].statusNeg();
1265  }
1266 
1267  // Define colour flow in branching.
1268  // Default corresponds to f -> f + gamma.
1269  int colMother = colDaughter;
1270  int acolMother = acolDaughter;
1271  int colSister = 0;
1272  int acolSister = 0;
1273  if (idSister == 22) ;
1274  // q -> q + g and 50% of g -> g + g; need new colour.
1275  else if (idSister == 21 && ( (idMother > 0 && idMother < 9)
1276  || (idMother == 21 && rndmPtr->flat() < 0.5) ) ) {
1277  colMother = event.nextColTag();
1278  colSister = colMother;
1279  acolSister = colDaughter;
1280  // qbar -> qbar + g and other 50% of g -> g + g; need new colour.
1281  } else if (idSister == 21) {
1282  acolMother = event.nextColTag();
1283  acolSister = acolMother;
1284  colSister = acolDaughter;
1285  // q -> g + q.
1286  } else if (idDaughterNow == 21 && idMother > 0) {
1287  colMother = colDaughter;
1288  acolMother = 0;
1289  colSister = acolDaughter;
1290  // qbar -> g + qbar
1291  } else if (idDaughterNow == 21) {
1292  acolMother = acolDaughter;
1293  colMother = 0;
1294  acolSister = colDaughter;
1295  // g -> q + qbar.
1296  } else if (idDaughterNow > 0 && idDaughterNow < 9) {
1297  acolMother = event.nextColTag();
1298  acolSister = acolMother;
1299  // g -> qbar + q.
1300  } else if (idDaughterNow < 0 && idDaughterNow > -9) {
1301  colMother = event.nextColTag();
1302  colSister = colMother;
1303  // q -> gamma + q.
1304  } else if (idDaughterNow == 22 && idMother > 0) {
1305  colMother = event.nextColTag();
1306  colSister = colMother;
1307  // qbar -> gamma + qbar.
1308  } else if (idDaughterNow == 22) {
1309  acolMother = event.nextColTag();
1310  acolSister = acolMother;
1311  }
1312 
1313  // Indices of partons involved. Add new sister.
1314  int iMother = eventSizeOld + side - 1;
1315  int iNewRecoiler = eventSizeOld + 2 - side;
1316  int iSister = event.append( idSister, 43, iMother, 0, 0, 0,
1317  colSister, acolSister, pSister, mSister, sqrt(pT2) );
1318 
1319  // References to the partons involved.
1320  Particle& daughter = event[iDaughter];
1321  Particle& mother = event[iMother];
1322  Particle& newRecoiler = event[iNewRecoiler];
1323  Particle& sister = event.back();
1324 
1325  // Replace old by new mother; update new recoiler.
1326  mother.id( idMother );
1327  mother.status( -41);
1328  mother.cols( colMother, acolMother);
1329  mother.p( pMother);
1330  newRecoiler.status( (normalRecoil) ? -42 : -46 );
1331  newRecoiler.p( pNewRec);
1332  if (!normalRecoil) newRecoiler.m( event[iRecoiler].m() );
1333 
1334  // Update mother and daughter pointers; also for beams.
1335  daughter.mothers( iMother, 0);
1336  mother.daughters( iSister, iDaughter);
1337  if (iSysSel == 0) {
1338  event[1].daughter1( (side == 1) ? iMother : iNewRecoiler );
1339  event[2].daughter1( (side == 2) ? iMother : iNewRecoiler );
1340  }
1341 
1342  // Find boost to old rest frame.
1343  RotBstMatrix Mtot;
1344  if (normalRecoil) Mtot.bst(0., 0., (x2 - x1) / (x1 + x2) );
1345  else if (side == 1)
1346  Mtot.toCMframe( event[iDaughter].p(), event[iRecoiler].p() );
1347  else Mtot.toCMframe( event[iRecoiler].p(), event[iDaughter].p() );
1348 
1349  // Initially select phi angle of branching at random.
1350  double phi = 2. * M_PI * rndmPtr->flat();
1351 
1352  // Evaluate coefficient of azimuthal asymmetry from gluon polarization.
1353  findAsymPol( event, dipEndSel);
1354  int iFinPol = dipEndSel->iFinPol;
1355  double cPol = dipEndSel->asymPol;
1356  double phiPol = (iFinPol == 0) ? 0. : event[iFinPol].phi();
1357 
1358  // If interference: try to match sister (anti)colour to final state.
1359  int iFinInt = 0;
1360  double cInt = 0.;
1361  double phiInt = 0.;
1362  if (doPhiIntAsym) {
1363  for (int i = 0; i < partonSystemsPtr->sizeOut(iSysSel); ++ i) {
1364  int iOut = partonSystemsPtr->getOut(iSysSel, i);
1365  if ( (acolSister != 0 && event[iOut].col() == acolSister)
1366  || (colSister != 0 && event[iOut].acol() == colSister) )
1367  iFinInt = iOut;
1368  }
1369  if (iFinInt != 0) {
1370  // Boost final-state parton to current frame of new kinematics.
1371  Vec4 pFinTmp = event[iFinInt].p();
1372  pFinTmp.rotbst(Mtot);
1373  double theFin = pFinTmp.theta();
1374  if (side == 2) theFin = M_PI - theFin;
1375  double theSis = pSister.theta();
1376  if (side == 2) theSis = M_PI - theSis;
1377  cInt = strengthIntAsym * 2. * theSis * theFin
1378  / (pow2(theSis) + pow2(theFin));
1379  phiInt = event[iFinInt].phi();
1380  }
1381  }
1382 
1383  // Bias phi distribution for polarization and interference.
1384  if (iFinPol != 0 || iFinInt != 0) {
1385  double cPhiPol, cPhiInt, weight;
1386  do {
1387  phi = 2. * M_PI * rndmPtr->flat();
1388  weight = 1.;
1389  if (iFinPol !=0 ) {
1390  cPhiPol = cos(phi - phiPol);
1391  weight *= ( 1. + cPol * (2. * pow2(cPhiPol) - 1.) )
1392  / ( 1. + abs(cPol) );
1393  }
1394  if (iFinInt !=0 ) {
1395  cPhiInt = cos(phi - phiInt);
1396  weight *= (1. - cInt) * (1. - cInt * cPhiInt)
1397  / (1. + pow2(cInt) - 2. * cInt * cPhiInt);
1398  }
1399  } while (weight < rndmPtr->flat());
1400  }
1401 
1402  // Include rotation -phi on boost to old rest frame.
1403  Mtot.rot(0., -phi);
1404 
1405  // Find boost from old rest frame to event cm frame.
1406  RotBstMatrix MfromRest;
1407  // The boost to the new rest frame.
1408  Vec4 sumNew = pMother + pNewRec;
1409  double betaX = sumNew.px() / sumNew.e();
1410  double betaZ = sumNew.pz() / sumNew.e();
1411  MfromRest.bst( -betaX, 0., -betaZ);
1412  // Alignment of radiator + recoiler to +- z axis, and rotation +phi.
1413  // Note: with spacelike (E < 0) recoiler p'_x_mother < 0 can happen!
1414  pMother.rotbst(MfromRest);
1415  double theta = pMother.theta();
1416  if (pMother.px() < 0.) theta = -theta;
1417  if (side == 2) theta += M_PI;
1418  MfromRest.rot(-theta, phi);
1419  // Boost to radiator + recoiler in event cm frame.
1420  if (normalRecoil) {
1421  MfromRest.bst( 0., 0., (x1New - x2New) / (x1New + x2New) );
1422  } else if (side == 1) {
1423  Vec4 pMotherWanted( 0., 0., 0.5 * eCM * x1New, 0.5 * eCM * x1New);
1424  MfromRest.fromCMframe( pMotherWanted, event[iRecoiler].p() );
1425 
1426  } else {
1427  Vec4 pMotherWanted( 0., 0., -0.5 * eCM * x2New, 0.5 * eCM * x2New);
1428  MfromRest.fromCMframe( event[iRecoiler].p(), pMotherWanted );
1429  }
1430  Mtot.rotbst(MfromRest);
1431 
1432  // Perform cumulative rotation/boost operation.
1433  // Mother, recoiler and sister from old rest frame to event cm frame.
1434  mother.rotbst(MfromRest);
1435  newRecoiler.rotbst(MfromRest);
1436  sister.rotbst(MfromRest);
1437  // The rest from (and to) event cm frame.
1438  for ( int i = eventSizeOld + 2; i < eventSizeOld + systemSizeOld; ++i)
1439  event[i].rotbst(Mtot);
1440 
1441  // Allow veto of branching. If so restore event record to before emission.
1442  if ( canVetoEmission
1443  && userHooksPtr->doVetoISREmission(eventSizeOld, event, iSysSel) ) {
1444  event.popBack( event.size() - eventSizeOld);
1445  event[1].daughter1( ev1dau1V);
1446  event[2].daughter1( ev2dau1V);
1447  for ( int iCopy = 0; iCopy < systemSizeOld; ++iCopy) {
1448  int iOldCopy = partonSystemsPtr->getAll(iSysSel, iCopy);
1449  event[iOldCopy].status( statusV[iCopy]);
1450  event[iOldCopy].mothers( mother1V[iCopy], mother2V[iCopy]);
1451  event[iOldCopy].daughters( daughter1V[iCopy], daughter2V[iCopy]);
1452  }
1453  return false;
1454  }
1455 
1456  // Update list of partons in system; adding newly produced one.
1457  partonSystemsPtr->setInA(iSysSel, eventSizeOld);
1458  partonSystemsPtr->setInB(iSysSel, eventSizeOld + 1);
1459  for (int iCopy = 2; iCopy < systemSizeOld; ++iCopy)
1460  partonSystemsPtr->setOut(iSysSel, iCopy - 2, eventSizeOld + iCopy);
1461  partonSystemsPtr->addOut(iSysSel, eventSizeOld + systemSizeOld);
1462  partonSystemsPtr->setSHat(iSysSel, m2 / z);
1463 
1464  // Update info on radiating dipole ends (QCD or QED).
1465  for (int iDip = 0; iDip < int(dipEnd.size()); ++iDip)
1466  if ( dipEnd[iDip].system == iSysSel) {
1467  if (abs(dipEnd[iDip].side) == side) {
1468  dipEnd[iDip].iRadiator = iMother;
1469  dipEnd[iDip].iRecoiler = iNewRecoiler;
1470  if (dipEnd[iDip].side > 0)
1471  dipEnd[iDip].colType = mother.colType();
1472  else {
1473  dipEnd[iDip].chgType = 0;
1474  if ( (mother.isQuark() && doQEDshowerByQ)
1475  || (mother.isLepton() && doQEDshowerByL) )
1476  dipEnd[iDip].chgType = mother.chargeType();
1477  }
1478  // Kill ME corrections after first emission.
1479  dipEnd[iDip].MEtype = 0;
1480 
1481  // Update info on recoiling dipole ends (QCD or QED).
1482  } else {
1483  dipEnd[iDip].iRadiator = iNewRecoiler;
1484  dipEnd[iDip].iRecoiler = iMother;
1485  // Optionally also kill recoiler ME corrections after first emission.
1486  if (!doMEafterFirst) dipEnd[iDip].MEtype = 0;
1487  }
1488  }
1489 
1490  // Update info on beam remnants.
1491  BeamParticle& beamNow = (side == 1) ? *beamAPtr : *beamBPtr;
1492  double xNew = (side == 1) ? x1New : x2New;
1493  beamNow[iSysSel].update( iMother, idMother, xNew);
1494  // Redo choice of companion kind whenever new flavour.
1495  if (idMother != idDaughterNow) {
1496  beamNow.xfISR( iSysSel, idMother, xNew, pT2);
1497  beamNow.pickValSeaComp();
1498  }
1499  BeamParticle& beamRec = (side == 1) ? *beamBPtr : *beamAPtr;
1500  beamRec[iSysSel].iPos( iNewRecoiler);
1501 
1502  // Store branching values of current dipole. (For rapidity ordering.)
1503  ++dipEndSel->nBranch;
1504  dipEndSel->pT2Old = pT2;
1505  dipEndSel->zOld = z;
1506 
1507  // Update history if recoiler rescatters.
1508  if (!normalRecoil)
1509  event[iRecoilMother].daughters( iNewRecoiler, iNewRecoiler);
1510 
1511  // Start list of rescatterers that force changed kinematics.
1512  vector<int> iRescatterer;
1513  for ( int i = 0; i < systemSizeOld - 2; ++i) {
1514  int iOutNew = partonSystemsPtr->getOut( iSysSel, i);
1515  if (!event[iOutNew].isFinal()) iRescatterer.push_back(iOutNew);
1516  }
1517 
1518  // Start iterate over list of such rescatterers.
1519  int iRescNow = -1;
1520  while (++iRescNow < int(iRescatterer.size())) {
1521 
1522  // Identify partons that induce or are affected by rescatter shift.
1523  // In following Old is before change of kinematics, New after,
1524  // Out scatterer in outstate and In in instate of another system.
1525  // Daughter sequence is (iOutOld ->) iOutNew -> iInNew -> iInOld.
1526  int iOutNew = iRescatterer[iRescNow];
1527  int iInOld = event[iOutNew].daughter1();
1528  int iSysResc = partonSystemsPtr->getSystemOf(iInOld, true);
1529 
1530  // Copy incoming partons of rescattered system and hook them up.
1531  int iOldA = partonSystemsPtr->getInA(iSysResc);
1532  int iOldB = partonSystemsPtr->getInB(iSysResc);
1533  bool rescSideA = event[iOldA].isRescatteredIncoming();
1534  int statusNewA = (rescSideA) ? -45 : -42;
1535  int statusNewB = (rescSideA) ? -42 : -45;
1536  int iNewA = event.copy(iOldA, statusNewA);
1537  int iNewB = event.copy(iOldB, statusNewB);
1538 
1539  // Copy outgoing partons of rescattered system and hook them up.
1540  int eventSize = event.size();
1541  int sizeOutAB = partonSystemsPtr->sizeOut(iSysResc);
1542  int iOldAB, statusOldAB, iNewAB;
1543  for (int iOutAB = 0; iOutAB < sizeOutAB; ++iOutAB) {
1544  iOldAB = partonSystemsPtr->getOut(iSysResc, iOutAB);
1545  statusOldAB = event[iOldAB].status();
1546  iNewAB = event.copy(iOldAB, 44);
1547  // Status could be negative for parton that rescatters in its turn.
1548  if (statusOldAB < 0) {
1549  event[iNewAB].statusNeg();
1550  iRescatterer.push_back(iNewAB);
1551  }
1552  }
1553 
1554  // Hook up new outgoing with new incoming parton.
1555  int iInNew = (rescSideA) ? iNewA : iNewB;
1556  event[iOutNew].daughters( iInNew, iInNew);
1557  event[iInNew].mothers( iOutNew, iOutNew);
1558 
1559  // Rescale recoiling incoming parton for correct invariant mass.
1560  event[iInNew].p( event[iOutNew].p() );
1561  double momFac = (rescSideA)
1562  ? event[iInOld].pPos() / event[iInNew].pPos()
1563  : event[iInOld].pNeg() / event[iInNew].pNeg();
1564  int iInRec = (rescSideA) ? iNewB : iNewA;
1565 
1566  // Rescatter: A previous boost may cause the light cone momentum of a
1567  // rescattered parton to change sign. If this happens, tell
1568  // parton level to try again.
1569  if (momFac < 0.0) {
1570  infoPtr->errorMsg("Warning in SpaceShower::branch: "
1571  "change in lightcone momentum sign; retrying parton level");
1572  rescatterFail = true;
1573  return false;
1574  }
1575 
1576  event[iInRec].rescale4( momFac);
1577 
1578  // Boost outgoing partons to new frame of incoming.
1579  RotBstMatrix MmodResc;
1580  MmodResc.toCMframe( event[iOldA].p(), event[iOldB].p());
1581  MmodResc.fromCMframe(event[iNewA].p(), event[iNewB].p());
1582  for (int iOutAB = 0; iOutAB < sizeOutAB; ++iOutAB)
1583  event[eventSize + iOutAB].rotbst(MmodResc);
1584 
1585  // Update list of partons in system.
1586  partonSystemsPtr->setInA(iSysResc, iNewA);
1587  partonSystemsPtr->setInB(iSysResc, iNewB);
1588  for (int iCopy = 0; iCopy < sizeOutAB; ++iCopy)
1589  partonSystemsPtr->setOut(iSysResc, iCopy, eventSize + iCopy);
1590 
1591  // Update info on radiating dipole ends (QCD or QED).
1592  for (int iDip = 0; iDip < int(dipEnd.size()); ++iDip)
1593  if ( dipEnd[iDip].system == iSysResc) {
1594  bool sideAnow = (abs(dipEnd[iDip].side) == 1);
1595  dipEnd[iDip].iRadiator = (sideAnow) ? iNewA : iNewB;
1596  dipEnd[iDip].iRecoiler = (sideAnow) ? iNewB : iNewA;
1597  }
1598 
1599  // Update info on beam remnants.
1600  BeamParticle& beamResc = (rescSideA) ? *beamAPtr : *beamBPtr;
1601  beamResc[iSysResc].iPos( iInNew);
1602  beamResc[iSysResc].p( event[iInNew].p() );
1603  beamResc[iSysResc].scaleX( 1. / momFac );
1604  BeamParticle& beamReco = (rescSideA) ? *beamBPtr : *beamAPtr;
1605  beamReco[iSysResc].iPos( iInRec);
1606  beamReco[iSysResc].scaleX( momFac);
1607 
1608  // End iterate over list of rescatterers.
1609  }
1610 
1611  // Check that beam momentum not used up by rescattered-system boosts.
1612  if (beamAPtr->xMax(-1) < 0.0 || beamBPtr->xMax(-1) < 0.0) {
1613  infoPtr->errorMsg("Warning in SpaceShower::branch: "
1614  "used up beam momentum; retrying parton level");
1615  rescatterFail = true;
1616  return false;
1617  }
1618 
1619  // Done without any errors.
1620  return true;
1621 
1622 }
1623 
1624 //--------------------------------------------------------------------------
1625 
1626 // Find class of ME correction.
1627 
1628 int SpaceShower::findMEtype( int iSys, Event& event) {
1629 
1630  // Default values and no action.
1631  int MEtype = 0;
1632  if (!doMEcorrections) ;
1633 
1634  // Identify systems producing a single resonance.
1635  else if (partonSystemsPtr->sizeOut( iSys) == 1) {
1636  int idIn1 = event[partonSystemsPtr->getInA(iSys)].id();
1637  int idIn2 = event[partonSystemsPtr->getInA(iSys)].id();
1638  int idRes = event[partonSystemsPtr->getOut(iSys, 0)].id();
1639  if (iSys == 0) idResFirst = abs(idRes);
1640  if (iSys == 1) idResSecond = abs(idRes);
1641 
1642  // f + fbar -> vector boson.
1643  if ( (idRes == 23 || abs(idRes) == 24 || idRes == 32
1644  || idRes == 33 || abs(idRes) == 34 || abs(idRes) == 41)
1645  && abs(idIn1) < 20 && abs(idIn2) < 20 ) MEtype = 1;
1646 
1647  // g + g, gamma + gamma -> Higgs boson.
1648  if ( (idRes == 25 || idRes == 35 || idRes == 36)
1649  && ( ( idIn1 == 21 && idIn2 == 21 )
1650  || ( idIn1 == 22 && idIn2 == 22 ) ) ) MEtype = 2;
1651 
1652  // f + fbar -> Higgs boson.
1653  if ( (idRes == 25 || idRes == 35 || idRes == 36)
1654  && abs(idIn1) < 20 && abs(idIn2) < 20 ) MEtype = 3;
1655  }
1656 
1657  // Done.
1658  return MEtype;
1659 
1660 }
1661 
1662 //--------------------------------------------------------------------------
1663 
1664 // Provide maximum of expected ME weight; for preweighting of evolution.
1665 
1666 double SpaceShower::calcMEmax( int MEtype, int idMother, int idDaughterIn) {
1667 
1668  // Currently only one non-unity case: g(gamma) f -> V f'.
1669  if (MEtype == 1 && idMother > 20 && idDaughterIn < 20) return 3.;
1670  return 1.;
1671 
1672 }
1673 
1674 //--------------------------------------------------------------------------
1675 
1676 // Provide actual ME weight for current branching.
1677 // Note: currently ME corrections are only allowed for first branching
1678 // on each side, so idDaughter is essentially known and checks overkill.
1679 
1680 double SpaceShower::calcMEcorr(int MEtype, int idMother, int idDaughterIn,
1681  double M2, double z, double Q2) {
1682 
1683  // Convert to Mandelstam variables. Sometimes may need to swap later.
1684  double sH = M2 / z;
1685  double tH = -Q2;
1686  double uH = Q2 - M2 * (1. - z) / z;
1687  int idMabs = abs(idMother);
1688  int idDabs = abs(idDaughterIn);
1689 
1690  // Corrections for f + fbar -> s-channel vector boson.
1691  if (MEtype == 1) {
1692  if (idMabs < 20 && idDabs < 20) {
1693  return (tH*tH + uH*uH + 2. * M2 * sH) / (sH*sH + M2*M2);
1694  } else if (idDabs < 20) {
1695  // g(gamma) f -> V f': -Q2 = (p_g - p_f')^2 in PS while
1696  // tHat = (p_f - p_f')^2 in ME so need to swap tHat <-> uHat.
1697  swap( tH, uH);
1698  return (sH*sH + uH*uH + 2. * M2 * tH) / (pow2(sH - M2) + M2*M2);
1699  }
1700 
1701  // Corrections for g + g -> Higgs boson.
1702  } else if (MEtype == 2) {
1703  if (idMabs < 20 && idDabs > 20) {
1704  return (sH*sH + uH*uH) / (sH*sH + pow2(sH - M2));
1705  } else if (idDabs > 20) {
1706  return 0.5 * (pow4(sH) + pow4(tH) + pow4(uH) + pow4(M2))
1707  / pow2(sH*sH - M2 * (sH - M2));
1708  }
1709 
1710  // Corrections for f + fbar -> Higgs boson (f = b mainly).
1711  } else if (MEtype == 3) {
1712  if (idMabs < 20 && idDabs < 20) {
1713  // The PS and ME answers agree for f fbar -> H g/gamma.
1714  return 1.;
1715  } else if (idDabs < 20) {
1716  // Need to swap tHat <-> uHat, cf. vector-boson production above.
1717  swap( tH, uH);
1718  return (sH*sH + uH*uH + 2. * (M2 - uH) * (M2 - sH))
1719  / (pow2(sH - M2) + M2*M2);
1720  }
1721  }
1722 
1723  return 1.;
1724 
1725 }
1726 
1727 //--------------------------------------------------------------------------
1728 
1729 // Find coefficient of azimuthal asymmetry from gluon polarization.
1730 
1731 void SpaceShower::findAsymPol( Event& event, SpaceDipoleEnd* dip) {
1732 
1733  // Default is no asymmetry. Only gluons are studied.
1734  dip->iFinPol = 0;
1735  dip->asymPol = 0.;
1736  int iRad = dip->iRadiator;
1737  if (!doPhiPolAsym || dip->idDaughter != 21) return;
1738 
1739  // Check if granddaughter in final state of hard scattering.
1740  // (May need to trace across carbon copies to find granddaughters.)
1741  // If so, only accept 2 -> 2 scatterings with gg or qq in final state.
1742  int iGrandD1 = event[iRad].daughter1();
1743  int iGrandD2 = event[iRad].daughter2();
1744  bool traceCopy = false;
1745  do {
1746  traceCopy = false;
1747  if (iGrandD1 > 0 && iGrandD2 == iGrandD1) {
1748  iGrandD1 = event[iGrandD2].daughter1();
1749  iGrandD2 = event[iGrandD2].daughter2();
1750  traceCopy = true;
1751  }
1752  } while (traceCopy);
1753  int statusGrandD1 = event[ iGrandD1 ].statusAbs();
1754  bool isHardProc = (statusGrandD1 == 23 || statusGrandD1 == 33);
1755  if (isHardProc) {
1756  if (iGrandD2 != iGrandD1 + 1) return;
1757  if (event[iGrandD1].isGluon() && event[iGrandD2].isGluon());
1758  else if (event[iGrandD1].isQuark() && event[iGrandD2].isQuark());
1759  else return;
1760  }
1761  dip->iFinPol = iGrandD1;
1762 
1763  // Coefficient from gluon production.
1764  if (dip->idMother == 21) dip->asymPol = pow2( (1. - dip->z)
1765  / (1. - dip->z * (1. - dip->z) ) );
1766  else dip->asymPol = 2. * (1. - dip->z) / (1. + pow2(1. - dip->z) );
1767 
1768  // Coefficients from gluon decay. Put z = 1/2 for hard process.
1769  double zDau = (isHardProc) ? 0.5 : dip->zOld;
1770  if (event[iGrandD1].isGluon()) dip->asymPol *= pow2( (1. - zDau)
1771  / (1. - zDau * (1. - zDau) ) );
1772  else dip->asymPol *= -2. * zDau *( 1. - zDau )
1773  / (1. - 2. * zDau * (1. - zDau) );
1774 
1775 }
1776 
1777 //--------------------------------------------------------------------------
1778 
1779 // Print the list of dipoles.
1780 
1781 void SpaceShower::list(ostream& os) const {
1782 
1783  // Header.
1784  os << "\n -------- PYTHIA SpaceShower Dipole Listing -------------- \n"
1785  << "\n i syst side rad rec pTmax col chg ME rec \n"
1786  << fixed << setprecision(3);
1787 
1788  // Loop over dipole list and print it.
1789  for (int i = 0; i < int(dipEnd.size()); ++i)
1790  os << setw(5) << i << setw(6) << dipEnd[i].system
1791  << setw(6) << dipEnd[i].side << setw(6) << dipEnd[i].iRadiator
1792  << setw(6) << dipEnd[i].iRecoiler << setw(12) << dipEnd[i].pTmax
1793  << setw(5) << dipEnd[i].colType << setw(5) << dipEnd[i].chgType
1794  << setw(5) << dipEnd[i].MEtype << setw(4)
1795  << dipEnd[i].normalRecoil << "\n";
1796 
1797  // Done.
1798  os << "\n -------- End PYTHIA SpaceShower Dipole Listing ----------"
1799  << endl;
1800 
1801 
1802 }
1803 
1804 //==========================================================================
1805 
1806 } // end namespace Pythia8
1807 
Definition: beam.h:43
Definition: AgUStep.h:26