StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MultipartonInteractions.cc
1 // MultipartonInteractions.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 // SigmaMultiparton and MultipartonInteractions classes.
8 
9 #include "Pythia8/MultipartonInteractions.h"
10 
11 // Internal headers for special processes.
12 #include "Pythia8/SigmaQCD.h"
13 #include "Pythia8/SigmaEW.h"
14 #include "Pythia8/SigmaOnia.h"
15 
16 namespace Pythia8 {
17 
18 //==========================================================================
19 
20 // The SigmaMultiparton class.
21 
22 //--------------------------------------------------------------------------
23 
24 // Constants: could be changed here if desired, but normally should not.
25 // These are of technical nature, as described for each.
26 
27 // The sum of outgoing masses must not be too close to the cm energy.
28 const double SigmaMultiparton::MASSMARGIN = 0.1;
29 
30 // Fraction of time not the dominant "gluon t-channel" process is picked.
31 const double SigmaMultiparton::OTHERFRAC = 0.2;
32 
33 //--------------------------------------------------------------------------
34 
35 // Initialize the generation process for given beams.
36 
37 bool SigmaMultiparton::init(int inState, int processLevel, Info* infoPtr,
38  Settings* settingsPtr, ParticleData* particleDataPtr, Rndm* rndmPtrIn,
39  BeamParticle* beamAPtr, BeamParticle* beamBPtr, Couplings* couplingsPtr) {
40 
41  // Store input pointer for future use.
42  rndmPtr = rndmPtrIn;
43 
44  // Reset vector sizes (necessary in case of re-initialization).
45  if (sigmaT.size() > 0) {
46  for (int i = 0; i < int(sigmaT.size()); ++i) delete sigmaT[i];
47  sigmaT.resize(0);
48  }
49  if (sigmaU.size() > 0) {
50  for (int i = 0; i < int(sigmaU.size()); ++i) delete sigmaU[i];
51  sigmaU.resize(0);
52  }
53 
54  // Always store mimimal set of processes:QCD 2 -> 2 t-channel.
55 
56  // Gluon-gluon instate.
57  if (inState == 0) {
58  sigmaT.push_back( new Sigma2gg2gg() );
59  sigmaU.push_back( new Sigma2gg2gg() );
60 
61  // Quark-gluon instate.
62  } else if (inState == 1) {
63  sigmaT.push_back( new Sigma2qg2qg() );
64  sigmaU.push_back( new Sigma2qg2qg() );
65 
66  // Quark-(anti)quark instate.
67  } else {
68  sigmaT.push_back( new Sigma2qq2qq() );
69  sigmaU.push_back( new Sigma2qq2qq() );
70  }
71 
72  // Normally store QCD processes to new flavour.
73  if (processLevel > 0) {
74  if (inState == 0) {
75  sigmaT.push_back( new Sigma2gg2qqbar() );
76  sigmaU.push_back( new Sigma2gg2qqbar() );
77  sigmaT.push_back( new Sigma2gg2QQbar(4, 121) );
78  sigmaU.push_back( new Sigma2gg2QQbar(4, 121) );
79  sigmaT.push_back( new Sigma2gg2QQbar(5, 123) );
80  sigmaU.push_back( new Sigma2gg2QQbar(5, 123) );
81  } else if (inState == 2) {
82  sigmaT.push_back( new Sigma2qqbar2gg() );
83  sigmaU.push_back( new Sigma2qqbar2gg() );
84  sigmaT.push_back( new Sigma2qqbar2qqbarNew() );
85  sigmaU.push_back( new Sigma2qqbar2qqbarNew() );
86  sigmaT.push_back( new Sigma2qqbar2QQbar(4, 122) );
87  sigmaU.push_back( new Sigma2qqbar2QQbar(4, 122) );
88  sigmaT.push_back( new Sigma2qqbar2QQbar(5, 124) );
89  sigmaU.push_back( new Sigma2qqbar2QQbar(5, 124) );
90  }
91  }
92 
93  // Optionally store electroweak processes, mainly photon production.
94  if (processLevel > 1) {
95  if (inState == 0) {
96  sigmaT.push_back( new Sigma2gg2ggamma() );
97  sigmaU.push_back( new Sigma2gg2ggamma() );
98  sigmaT.push_back( new Sigma2gg2gammagamma() );
99  sigmaU.push_back( new Sigma2gg2gammagamma() );
100  } else if (inState == 1) {
101  sigmaT.push_back( new Sigma2qg2qgamma() );
102  sigmaU.push_back( new Sigma2qg2qgamma() );
103  } else if (inState == 2) {
104  sigmaT.push_back( new Sigma2qqbar2ggamma() );
105  sigmaU.push_back( new Sigma2qqbar2ggamma() );
106  sigmaT.push_back( new Sigma2ffbar2gammagamma() );
107  sigmaU.push_back( new Sigma2ffbar2gammagamma() );
108  sigmaT.push_back( new Sigma2ffbar2ffbarsgm() );
109  sigmaU.push_back( new Sigma2ffbar2ffbarsgm() );
110  }
111  if (inState >= 2) {
112  sigmaT.push_back( new Sigma2ff2fftgmZ() );
113  sigmaU.push_back( new Sigma2ff2fftgmZ() );
114  sigmaT.push_back( new Sigma2ff2fftW() );
115  sigmaU.push_back( new Sigma2ff2fftW() );
116  }
117  }
118 
119  // Optionally store charmonium and bottomonium production.
120  if (processLevel > 2) {
121  SigmaOniaSetup charmonium(infoPtr, settingsPtr, particleDataPtr, 4);
122  SigmaOniaSetup bottomonium(infoPtr, settingsPtr, particleDataPtr, 5);
123  if (inState == 0) {
124  charmonium.setupSigma2gg(sigmaT, true);
125  charmonium.setupSigma2gg(sigmaU, true);
126  bottomonium.setupSigma2gg(sigmaT, true);
127  bottomonium.setupSigma2gg(sigmaU, true);
128  } else if (inState == 1) {
129  charmonium.setupSigma2qg(sigmaT, true);
130  charmonium.setupSigma2qg(sigmaU, true);
131  bottomonium.setupSigma2qg(sigmaT, true);
132  bottomonium.setupSigma2qg(sigmaU, true);
133  } else if (inState == 2) {
134  charmonium.setupSigma2qq(sigmaT, true);
135  charmonium.setupSigma2qq(sigmaU, true);
136  bottomonium.setupSigma2qq(sigmaT, true);
137  bottomonium.setupSigma2qq(sigmaU, true);
138  }
139  }
140 
141  // Resize arrays to match sizes above.
142  nChan = sigmaT.size();
143  needMasses.resize(nChan);
144  m3Fix.resize(nChan);
145  m4Fix.resize(nChan);
146  sHatMin.resize(nChan);
147  sigmaTval.resize(nChan);
148  sigmaUval.resize(nChan);
149 
150  // Initialize the processes.
151  for (int i = 0; i < nChan; ++i) {
152  sigmaT[i]->init( infoPtr, settingsPtr, particleDataPtr, rndmPtr,
153  beamAPtr, beamBPtr, couplingsPtr);
154  sigmaT[i]->initProc();
155  sigmaU[i]->init( infoPtr, settingsPtr, particleDataPtr, rndmPtr,
156  beamAPtr, beamBPtr, couplingsPtr);
157  sigmaU[i]->initProc();
158 
159  // Prepare for massive kinematics (but fixed masses!) where required.
160  needMasses[i] = false;
161  int id3Mass = sigmaT[i]->id3Mass();
162  int id4Mass = sigmaT[i]->id4Mass();
163  m3Fix[i] = 0.;
164  m4Fix[i] = 0.;
165  if (id3Mass > 0 || id4Mass > 0) {
166  needMasses[i] = true;
167  m3Fix[i] = particleDataPtr->m0(id3Mass);
168  m4Fix[i] = particleDataPtr->m0(id4Mass);
169  }
170  sHatMin[i] = pow2( m3Fix[i] + m4Fix[i] + MASSMARGIN);
171  }
172 
173  // Done.
174  return true;
175 
176 }
177 
178 //--------------------------------------------------------------------------
179 
180 // Calculate cross section summed over possibilities.
181 
182 double SigmaMultiparton::sigma( int id1, int id2, double x1, double x2,
183  double sHat, double tHat, double uHat, double alpS, double alpEM,
184  bool restore, bool pickOtherIn) {
185 
186  // Choose either the dominant process (in slot 0) or the rest of them.
187  if (restore) pickOther = pickOtherIn;
188  else pickOther = (rndmPtr->flat() < OTHERFRAC);
189 
190  // Iterate over all subprocesses.
191  sigmaTsum = 0.;
192  sigmaUsum = 0.;
193  for (int i = 0; i < nChan; ++i) {
194  sigmaTval[i] = 0.;
195  sigmaUval[i] = 0.;
196 
197  // Skip the not chosen processes.
198  if (i == 0 && pickOther) continue;
199  if (i > 0 && !pickOther) continue;
200 
201  // t-channel-sampling contribution.
202  if (sHat > sHatMin[i]) {
203  sigmaT[i]->set2KinMPI( x1, x2, sHat, tHat, uHat,
204  alpS, alpEM, needMasses[i], m3Fix[i], m4Fix[i]);
205  sigmaTval[i] = sigmaT[i]->sigmaHatWrap(id1, id2);
206  sigmaT[i]->pickInState(id1, id2);
207  // Correction factor for tHat rescaling in massive kinematics.
208  if (needMasses[i]) sigmaTval[i] *= sigmaT[i]->sHBetaMPI() / sHat;
209  sigmaTsum += sigmaTval[i];
210  }
211 
212  // u-channel-sampling contribution.
213  if (sHat > sHatMin[i]) {
214  sigmaU[i]->set2KinMPI( x1, x2, sHat, uHat, tHat,
215  alpS, alpEM, needMasses[i], m3Fix[i], m4Fix[i]);
216  sigmaUval[i] = sigmaU[i]->sigmaHatWrap( id1, id2);
217  sigmaU[i]->pickInState(id1, id2);
218  // Correction factor for tHat rescaling in massive kinematics.
219  if (needMasses[i]) sigmaUval[i] *= sigmaU[i]->sHBetaMPI() / sHat;
220  sigmaUsum += sigmaUval[i];
221  }
222 
223  // Average of t- and u-channel sampling; corrected for not selected channels.
224  }
225  double sigmaAvg = 0.5 * (sigmaTsum + sigmaUsum);
226  if (pickOther) sigmaAvg /= OTHERFRAC;
227  if (!pickOther) sigmaAvg /= (1. - OTHERFRAC);
228  return sigmaAvg;
229 
230 }
231 
232 //--------------------------------------------------------------------------
233 
234 // Return one subprocess, picked according to relative cross sections.
235 
236 SigmaProcess* SigmaMultiparton::sigmaSel() {
237 
238  // Decide between t- and u-channel-sampled kinematics.
239  pickedU = (rndmPtr->flat() * (sigmaTsum + sigmaUsum) < sigmaUsum);
240 
241  // Pick one of t-channel-sampled processes.
242  if (!pickedU) {
243  double sigmaRndm = sigmaTsum * rndmPtr->flat();
244  int iPick = -1;
245  do sigmaRndm -= sigmaTval[++iPick];
246  while (sigmaRndm > 0.);
247  return sigmaT[iPick];
248 
249  // Pick one of u-channel-sampled processes.
250  } else {
251  double sigmaRndm = sigmaUsum * rndmPtr->flat();
252  int iPick = -1;
253  do sigmaRndm -= sigmaUval[++iPick];
254  while (sigmaRndm > 0.);
255  return sigmaU[iPick];
256  }
257 
258 }
259 
260 //==========================================================================
261 
262 // The MultipartonInteractions class.
263 
264 //--------------------------------------------------------------------------
265 
266 // Constants: could be changed here if desired, but normally should not.
267 // These are of technical nature, as described for each.
268 
269 // Factorization scale pT2 by default, but could be shifted to pT2 + pT02.
270 // (A priori possible, but problems for flavour threshold interpretation.)
271 const bool MultipartonInteractions::SHIFTFACSCALE = false;
272 
273 // Pick one parton to represent rescattering cross section to speed up.
274 const bool MultipartonInteractions::PREPICKRESCATTER = true;
275 
276 // Naive upper estimate of cross section too pessimistic, so reduce by this.
277 const double MultipartonInteractions::SIGMAFUDGE = 0.8;
278 
279 // The r value above, picked to allow a flatter correct/trial cross section.
280 const double MultipartonInteractions::RPT20 = 0.25;
281 
282 // Reduce pT0 by factor pT0STEP if sigmaInt < SIGMASTEP * sigmaND.
283 const double MultipartonInteractions::PT0STEP = 0.9;
284 const double MultipartonInteractions::SIGMASTEP = 1.1;
285 
286 // Stop if pT0 or pTmin fall below this, or alpha_s blows up.
287 const double MultipartonInteractions::PT0MIN = 0.2;
288 
289 // Refuse too low expPow in impact parameter profile.
290 const double MultipartonInteractions::EXPPOWMIN = 0.4;
291 
292 // Define low-b region by interaction rate above given number.
293 const double MultipartonInteractions::PROBATLOWB = 0.6;
294 
295 // Basic step size for b integration; sometimes modified.
296 const double MultipartonInteractions::BSTEP = 0.01;
297 
298 // Stop b integration when integrand dropped enough.
299 const double MultipartonInteractions::BMAX = 1e-8;
300 
301 // Do not allow too large argument to exp function.
302 const double MultipartonInteractions::EXPMAX = 50.;
303 
304 // Convergence criterion for k iteration.
305 const double MultipartonInteractions::KCONVERGE = 1e-7;
306 
307 // Conversion of GeV^{-2} to mb for cross section.
308 const double MultipartonInteractions::CONVERT2MB = 0.389380;
309 
310 // Stay away from division by zero in Jacobian for tHat -> pT2.
311 const double MultipartonInteractions::ROOTMIN = 0.01;
312 
313 // No need to reinitialize parameters if energy close to previous.
314 const double MultipartonInteractions::ECMDEV = 0.01;
315 
316 // Settings for x-dependent matter profile:
317 // Number of bins in b (with these settings, no bStep increase and
318 // reintegration needed with a1 ~ 0.20 up to ECM ~ 40TeV).
319 const int MultipartonInteractions::XDEP_BBIN = 500;
320 // Initial value of a0.
321 const double MultipartonInteractions::XDEP_A0 = 1.0;
322 // Width of form ( XDEP_A1 + a1 * log(1 / x) ).
323 const double MultipartonInteractions::XDEP_A1 = 1.0;
324 // Initial step size in b and increment.
325 const double MultipartonInteractions::XDEP_BSTEP = 0.02;
326 const double MultipartonInteractions::XDEP_BSTEPINC = 0.01;
327 // Overlap-weighted cross section in last bin of b must be beneath
328 // XDEP_CUTOFF * sigmaInt.
329 const double MultipartonInteractions::XDEP_CUTOFF = 1e-4;
330 // a0 is calculated in units of sqrt(mb), so convert to fermi.
331 const double MultipartonInteractions::XDEP_SMB2FM = sqrt(0.1);
332 
333 // Only write warning when weight clearly above unity.
334 const double MultipartonInteractions::WTACCWARN = 1.1;
335 
336 //--------------------------------------------------------------------------
337 
338 // Initialize the generation process for given beams.
339 
340 bool MultipartonInteractions::init( bool doMPIinit, int iDiffSysIn,
341  Info* infoPtrIn, Settings& settings, ParticleData* particleDataPtr,
342  Rndm* rndmPtrIn, BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
343  Couplings* couplingsPtrIn, PartonSystems* partonSystemsPtrIn,
344  SigmaTotal* sigmaTotPtrIn, UserHooks* userHooksPtrIn, ostream& os) {
345 
346  // Store input pointers for future use. Done if no initialization.
347  iDiffSys = iDiffSysIn;
348  infoPtr = infoPtrIn;
349  rndmPtr = rndmPtrIn;
350  beamAPtr = beamAPtrIn;
351  beamBPtr = beamBPtrIn;
352  couplingsPtr = couplingsPtrIn;
353  partonSystemsPtr = partonSystemsPtrIn;
354  sigmaTotPtr = sigmaTotPtrIn;
355  userHooksPtr = userHooksPtrIn;
356  if (!doMPIinit) return false;
357 
358  // If both beams are baryons then softer PDF's than for mesons/Pomerons.
359  hasBaryonBeams = ( beamAPtr->isBaryon() && beamBPtr->isBaryon() );
360 
361  // Matching in pT of hard interaction to further interactions.
362  pTmaxMatch = settings.mode("MultipartonInteractions:pTmaxMatch");
363 
364  // Parameters of alphaStrong generation.
365  alphaSvalue = settings.parm("MultipartonInteractions:alphaSvalue");
366  alphaSorder = settings.mode("MultipartonInteractions:alphaSorder");
367  alphaSnfmax = settings.mode("StandardModel:alphaSnfmax");
368 
369  // Parameters of alphaEM generation.
370  alphaEMorder = settings.mode("MultipartonInteractions:alphaEMorder");
371 
372  // Parameters of cross section generation.
373  Kfactor = settings.parm("MultipartonInteractions:Kfactor");
374 
375  // Regularization of QCD evolution for pT -> 0.
376  pT0Ref = settings.parm("MultipartonInteractions:pT0Ref");
377  ecmRef = settings.parm("MultipartonInteractions:ecmRef");
378  ecmPow = settings.parm("MultipartonInteractions:ecmPow");
379  pTmin = settings.parm("MultipartonInteractions:pTmin");
380 
381  // Impact parameter profile: nondiffractive topologies.
382  if (iDiffSys == 0) {
383  bProfile = settings.mode("MultipartonInteractions:bProfile");
384  coreRadius = settings.parm("MultipartonInteractions:coreRadius");
385  coreFraction = settings.parm("MultipartonInteractions:coreFraction");
386  expPow = settings.parm("MultipartonInteractions:expPow");
387  expPow = max(EXPPOWMIN, expPow);
388  // x-dependent impact parameter profile.
389  a1 = settings.parm("MultipartonInteractions:a1");
390 
391  // Impact parameter profile: diffractive topologies.
392  } else {
393  bProfile = settings.mode("Diffraction:bProfile");
394  coreRadius = settings.parm("Diffraction:coreRadius");
395  coreFraction = settings.parm("Diffraction:coreFraction");
396  expPow = settings.parm("Diffraction:expPow");
397  expPow = max(EXPPOWMIN, expPow);
398  }
399 
400  // Common choice of "pT" scale for determining impact parameter.
401  bSelScale = settings.mode("MultipartonInteractions:bSelScale");
402 
403  // Process sets to include in machinery.
404  processLevel = settings.mode("MultipartonInteractions:processLevel");
405 
406  // Parameters of rescattering description.
407  allowRescatter = settings.flag("MultipartonInteractions:allowRescatter");
408  allowDoubleRes
409  = settings.flag("MultipartonInteractions:allowDoubleRescatter");
410  rescatterMode = settings.mode("MultipartonInteractions:rescatterMode");
411  ySepResc = settings.parm("MultipartonInteractions:ySepRescatter");
412  deltaYResc = settings.parm("MultipartonInteractions:deltaYRescatter");
413 
414  // Rescattering not yet implemented for x-dependent impact profile.
415  if (bProfile == 4) allowRescatter = false;
416 
417  // A global recoil FSR stategy restricts rescattering.
418  globalRecoilFSR = settings.flag("TimeShower:globalRecoil");
419  nMaxGlobalRecoilFSR = settings.mode("TimeShower:nMaxGlobalRecoil");
420 
421  // Various other parameters.
422  nQuarkIn = settings.mode("MultipartonInteractions:nQuarkIn");
423  nSample = settings.mode("MultipartonInteractions:nSample");
424 
425  // Optional dampening at small pT's when large multiplicities.
426  enhanceScreening = settings.mode("MultipartonInteractions:enhanceScreening");
427 
428  // Parameters for diffractive systems.
429  sigmaPomP = settings.parm("Diffraction:sigmaRefPomP");
430  mPomP = settings.parm("Diffraction:mRefPomP");
431  pPomP = settings.parm("Diffraction:mPowPomP");
432  mMinPertDiff = settings.parm("Diffraction:mMinPert");
433 
434  // Possibility to allow user veto of MPI
435  canVetoMPI = (userHooksPtr != 0) ? userHooksPtr->canVetoMPIEmission()
436  : false;
437 
438  // Some common combinations for double Gaussian, as shorthand.
439  if (bProfile == 2) {
440  fracA = pow2(1. - coreFraction);
441  fracB = 2. * coreFraction * (1. - coreFraction);
442  fracC = pow2(coreFraction);
443  radius2B = 0.5 * (1. + pow2(coreRadius));
444  radius2C = pow2(coreRadius);
445 
446  // Some common combinations for exp(b^pow), as shorthand.
447  } else if (bProfile == 3) {
448  hasLowPow = (expPow < 2.);
449  expRev = 2. / expPow - 1.;
450  }
451 
452  // Initialize alpha_strong generation.
453  alphaS.init( alphaSvalue, alphaSorder, alphaSnfmax, false);
454  double Lambda3 = alphaS.Lambda3();
455 
456  // Initialize alphaEM generation.
457  alphaEM.init( alphaEMorder, &settings);
458 
459  // Attach matrix-element calculation objects.
460  sigma2gg.init( 0, processLevel, infoPtr, &settings, particleDataPtr,
461  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
462  sigma2qg.init( 1, processLevel, infoPtr, &settings, particleDataPtr,
463  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
464  sigma2qqbarSame.init( 2, processLevel, infoPtr, &settings, particleDataPtr,
465  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
466  sigma2qq.init( 3, processLevel, infoPtr, &settings, particleDataPtr,
467  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
468 
469  // Calculate invariant mass of system.
470  eCM = infoPtr->eCM();
471  sCM = eCM * eCM;
472  mMaxPertDiff = eCM;
473  eCMsave = eCM;
474 
475  // Get the total inelastic and nondiffractive cross section.
476  if (!sigmaTotPtr->hasSigmaTot()) return false;
477  bool isNonDiff = (iDiffSys == 0);
478  sigmaND = sigmaTotPtr->sigmaND();
479  double sigmaMaxViol = 0.;
480 
481  // Output initialization info - first part.
482  bool showMPI = settings.flag("Init:showMultipartonInteractions");
483  if (showMPI) {
484  os << "\n *------- PYTHIA Multiparton Interactions Initialization "
485  << "---------* \n"
486  << " | "
487  << " | \n";
488  if (isNonDiff)
489  os << " | sigmaNonDiffractive = " << fixed
490  << setprecision(2) << setw(7) << sigmaND << " mb | \n";
491  else if (iDiffSys == 1)
492  os << " | diffraction XB "
493  << " | \n";
494  else if (iDiffSys == 2)
495  os << " | diffraction AX "
496  << " | \n";
497  else if (iDiffSys == 3)
498  os << " | diffraction AXB "
499  << " | \n";
500  os << " | "
501  << " | \n";
502  }
503 
504  // For diffraction need to cover range of diffractive masses.
505  nStep = (iDiffSys == 0) ? 1 : 5;
506  eStepSize = (nStep < 2) ? 1.
507  : log(mMaxPertDiff / mMinPertDiff) / (nStep - 1.);
508  for (int iStep = 0; iStep < nStep; ++iStep) {
509 
510  // Update and output current diffractive mass and
511  // fictitious Pomeron-proton cross section for normalization.
512  if (nStep > 1) {
513  eCM = mMinPertDiff * pow( mMaxPertDiff / mMinPertDiff,
514  iStep / (nStep - 1.) );
515  sCM = eCM * eCM;
516  sigmaND = sigmaPomP * pow( eCM / mPomP, pPomP);
517  if (showMPI) os << " | diffractive mass = " << scientific
518  << setprecision(3) << setw(9) << eCM << " GeV and sigmaNorm = "
519  << fixed << setw(6) << sigmaND << " mb | \n";
520  }
521 
522  // Set current pT0 scale.
523  pT0 = pT0Ref * pow(eCM / ecmRef, ecmPow);
524 
525  // The pT0 value may need to be decreased, if sigmaInt < sigmaND.
526  double pT4dSigmaMaxBeg = 0.;
527  for ( ; ; ) {
528 
529  // Derived pT kinematics combinations.
530  pT20 = pT0*pT0;
531  pT2min = pTmin*pTmin;
532  pTmax = 0.5*eCM;
533  pT2max = pTmax*pTmax;
534  pT20R = RPT20 * pT20;
535  pT20minR = pT2min + pT20R;
536  pT20maxR = pT2max + pT20R;
537  pT20min0maxR = pT20minR * pT20maxR;
538  pT2maxmin = pT2max - pT2min;
539 
540  // Provide upper estimate of interaction rate d(Prob)/d(pT2).
541  upperEnvelope();
542 
543  // Setup binning in b for x-dependent matter profile.
544  if (bProfile == 4) {
545  sigmaIntWgt.resize(XDEP_BBIN);
546  sigmaSumWgt.resize(XDEP_BBIN);
547  bstepNow = XDEP_BSTEP;
548  }
549 
550  // Integrate the parton-parton interaction cross section.
551  pT4dSigmaMaxBeg = pT4dSigmaMax;
552  jetCrossSection();
553 
554  // If the overlap-weighted cross section has not fallen below
555  // cutoff, then increase bin size in b and reintegrate.
556  while (bProfile == 4
557  && sigmaIntWgt[XDEP_BBIN - 1] > XDEP_CUTOFF * sigmaInt) {
558  bstepNow += XDEP_BSTEPINC;
559  jetCrossSection();
560  }
561 
562  // Sufficiently big SigmaInt or reduce pT0; maybe also pTmin.
563  if (sigmaInt > SIGMASTEP * sigmaND) break;
564  if (showMPI) os << fixed << setprecision(2) << " | pT0 = "
565  << setw(5) << pT0 << " gives sigmaInteraction = " << setw(7)
566  << sigmaInt << " mb: rejected | \n";
567  if (pTmin > pT0) pTmin *= PT0STEP;
568  pT0 *= PT0STEP;
569 
570  // Give up if pT0 and pTmin fall too low.
571  if ( max(pT0, pTmin) < max(PT0MIN, Lambda3) ) {
572  infoPtr->errorMsg("Error in MultipartonInteractions::init:"
573  " failed to find acceptable pT0 and pTmin");
574  infoPtr->setTooLowPTmin(true);
575  return false;
576  }
577  }
578 
579  // Output for accepted pT0.
580  if (showMPI) os << fixed << setprecision(2) << " | pT0 = "
581  << setw(5) << pT0 << " gives sigmaInteraction = "<< setw(7)
582  << sigmaInt << " mb: accepted | \n";
583 
584  // Calculate factor relating matter overlap and interaction rate.
585  overlapInit();
586 
587  // Maximum violation relative to first estimate.
588  sigmaMaxViol = max( sigmaMaxViol, pT4dSigmaMax / pT4dSigmaMaxBeg);
589 
590  // Save values calculated.
591  if (nStep > 1) {
592  pT0Save[iStep] = pT0;
593  pT4dSigmaMaxSave[iStep] = pT4dSigmaMax;
594  pT4dProbMaxSave[iStep] = pT4dProbMax;
595  sigmaIntSave[iStep] = sigmaInt;
596  for (int j = 0; j <= 100; ++j) sudExpPTSave[iStep][j] = sudExpPT[j];
597  zeroIntCorrSave[iStep] = zeroIntCorr;
598  normOverlapSave[iStep] = normOverlap;
599  kNowSave[iStep] = kNow;
600  bAvgSave[iStep] = bAvg;
601  bDivSave[iStep] = bDiv;
602  probLowBSave[iStep] = probLowB;
603  fracAhighSave[iStep] = fracAhigh;
604  fracBhighSave[iStep] = fracBhigh;
605  fracChighSave[iStep] = fracBhigh;
606  fracABChighSave[iStep] = fracABChigh;
607  cDivSave[iStep] = cDiv;
608  cMaxSave[iStep] = cMax;
609  }
610 
611  // End of loop over diffractive masses.
612  }
613 
614  // Output details for x-dependent matter profile.
615  if (bProfile == 4 && showMPI)
616  os << " | "
617  << " | \n"
618  << fixed << setprecision(2)
619  << " | x-dependent matter profile: a1 = " << a1 << ", "
620  << "a0 = " << a0now * XDEP_SMB2FM << ", bStep = "
621  << bstepNow << " | \n";
622 
623  // End initialization printout.
624  if (showMPI) os << " | "
625  << " | \n"
626  << " *------- End PYTHIA Multiparton Interactions Initialization"
627  << " -----* " << endl;
628 
629  // Amount of violation from upperEnvelope to jetCrossSection.
630  if (sigmaMaxViol > 1.) {
631  ostringstream osWarn;
632  osWarn << "by factor " << fixed << setprecision(3) << sigmaMaxViol;
633  infoPtr->errorMsg("Warning in MultipartonInteractions::init:"
634  " maximum increased", osWarn.str());
635  }
636 
637  // Reset statistics.
638  SigmaMultiparton* dSigma;
639  for (int i = 0; i < 4; ++i) {
640  if (i == 0) dSigma = &sigma2gg;
641  else if (i == 1) dSigma = &sigma2qg;
642  else if (i == 2) dSigma = &sigma2qqbarSame;
643  else dSigma = &sigma2qq;
644  int nProc = dSigma->nProc();
645  for (int iProc = 0; iProc < nProc; ++iProc)
646  nGen[ dSigma->codeProc(iProc) ] = 0;
647  }
648 
649  // Additional setup for x-dependent matter profile.
650  if (bProfile == 4) {
651  sigmaIntWgt.clear();
652  sigmaSumWgt.clear();
653  }
654  // No preselection of sea/valence content and initialise a0.
655  vsc1 = 0;
656  vsc2 = 0;
657 
658  // Done.
659  return true;
660 }
661 
662 //--------------------------------------------------------------------------
663 
664 // Reset impact parameter choice and update the CM energy.
665 // For diffraction also interpolate parameters to current CM energy.
666 
667 void MultipartonInteractions::reset( ) {
668 
669  // Reset impact parameter choice.
670  bIsSet = false;
671  bSetInFirst = false;
672 
673  // Update CM energy. Done if not diffraction and not new energy.
674  eCM = infoPtr->eCM();
675  sCM = eCM * eCM;
676  if (nStep == 1 || abs( eCM / eCMsave - 1.) < ECMDEV) return;
677 
678  // Set fictitious Pomeron-proton cross section for diffractive system.
679  sigmaND = sigmaPomP * pow( eCM / mPomP, pPomP);
680 
681  // Current interpolation point.
682  eCMsave = eCM;
683  eStepSave = log(eCM / mMinPertDiff) / eStepSize;
684  iStepFrom = max( 0, min( nStep - 2, int( eStepSave) ) );
685  iStepTo = iStepFrom + 1;
686  eStepTo = max( 0., min( 1., eStepSave - iStepFrom) );
687  eStepFrom = 1. - eStepTo;
688 
689  // Update pT0 and combinations derived from it.
690  pT0 = eStepFrom * pT0Save[iStepFrom]
691  + eStepTo * pT0Save[iStepTo];
692  pT20 = pT0*pT0;
693  pT2min = pTmin*pTmin;
694  pTmax = 0.5*eCM;
695  pT2max = pTmax*pTmax;
696  pT20R = RPT20 * pT20;
697  pT20minR = pT2min + pT20R;
698  pT20maxR = pT2max + pT20R;
699  pT20min0maxR = pT20minR * pT20maxR;
700  pT2maxmin = pT2max - pT2min;
701 
702  // Update other parameters used in pT choice.
703  pT4dSigmaMax = eStepFrom * pT4dSigmaMaxSave[iStepFrom]
704  + eStepTo * pT4dSigmaMaxSave[iStepTo];
705  pT4dProbMax = eStepFrom * pT4dProbMaxSave[iStepFrom]
706  + eStepTo * pT4dProbMaxSave[iStepTo];
707  sigmaInt = eStepFrom * sigmaIntSave[iStepFrom]
708  + eStepTo * sigmaIntSave[iStepTo];
709  for (int j = 0; j <= 100; ++j)
710  sudExpPT[j] = eStepFrom * sudExpPTSave[iStepFrom][j]
711  + eStepTo * sudExpPTSave[iStepTo][j];
712 
713  // Update parameters related to the impact-parameter picture.
714  zeroIntCorr = eStepFrom * zeroIntCorrSave[iStepFrom]
715  + eStepTo * zeroIntCorrSave[iStepTo];
716  normOverlap = eStepFrom * normOverlapSave[iStepFrom]
717  + eStepTo * normOverlapSave[iStepTo];
718  kNow = eStepFrom * kNowSave[iStepFrom]
719  + eStepTo * kNowSave[iStepTo];
720  bAvg = eStepFrom * bAvgSave[iStepFrom]
721  + eStepTo * bAvgSave[iStepTo];
722  bDiv = eStepFrom * bDivSave[iStepFrom]
723  + eStepTo * bDivSave[iStepTo];
724  probLowB = eStepFrom * probLowBSave[iStepFrom]
725  + eStepTo * probLowBSave[iStepTo];
726  fracAhigh = eStepFrom * fracAhighSave[iStepFrom]
727  + eStepTo * fracAhighSave[iStepTo];
728  fracBhigh = eStepFrom * fracBhighSave[iStepFrom]
729  + eStepTo * fracBhighSave[iStepTo];
730  fracChigh = eStepFrom * fracChighSave[iStepFrom]
731  + eStepTo * fracChighSave[iStepTo];
732  fracABChigh = eStepFrom * fracABChighSave[iStepFrom]
733  + eStepTo * fracABChighSave[iStepTo];
734  cDiv = eStepFrom * cDivSave[iStepFrom]
735  + eStepTo * cDivSave[iStepTo];
736  cMax = eStepFrom * cMaxSave[iStepFrom]
737  + eStepTo * cMaxSave[iStepTo];
738 
739 }
740 
741 //--------------------------------------------------------------------------
742 
743 // Select first = hardest pT in nondiffractive process.
744 // Requires separate treatment at low and high b values.
745 
746 void MultipartonInteractions::pTfirst() {
747 
748  // Pick impact parameter and thereby interaction rate enhancement.
749  // This is not used for the x-dependent matter profile, which
750  // instead uses trial interactions.
751  if (bProfile == 4) isAtLowB = false;
752  else overlapFirst();
753  bSetInFirst = true;
754  double WTacc;
755 
756  // At low b values evolve downwards with Sudakov.
757  if (isAtLowB) {
758  pT2 = pT2max;
759  do {
760 
761  // Pick a pT using a quick-and-dirty cross section estimate.
762  pT2 = fastPT2(pT2);
763 
764  // If fallen below lower cutoff then need to restart at top.
765  if (pT2 < pT2min) {
766  pT2 = pT2max;
767  WTacc = 0.;
768 
769  // Else pick complete kinematics and evaluate cross-section correction.
770  } else {
771  WTacc = sigmaPT2scatter(true) / dSigmaApprox;
772  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
773  "MultipartonInteractions::pTfirst: weight above unity");
774  }
775 
776  // Loop until acceptable pT and acceptable kinematics.
777  } while (WTacc < rndmPtr->flat() || !dSigmaDtSel->final2KinMPI());
778 
779  // At high b values make preliminary pT choice without Sudakov factor.
780  } else {
781 
782  while (true) {
783  do {
784  pT2 = pT20min0maxR / (pT20minR + rndmPtr->flat() * pT2maxmin) - pT20R;
785 
786  // Evaluate upper estimate of cross section for this pT2 choice.
787  dSigmaApprox = pT4dSigmaMax / pow2(pT2 + pT20R);
788 
789  // Pick complete kinematics and evaluate cross-section correction.
790  WTacc = sigmaPT2scatter(true) / dSigmaApprox;
791 
792  // Evaluate and include Sudakov factor.
793  if (bProfile != 4) WTacc *= sudakov( pT2, enhanceB);
794 
795  // Warn for weight above unity
796  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
797  "MultipartonInteractions::pTfirst: weight above unity");
798 
799  // Loop until acceptable pT and acceptable kinematics.
800  } while (WTacc < rndmPtr->flat() || !dSigmaDtSel->final2KinMPI());
801 
802  // For x-dependent matter profile, use trial interactions to
803  // generate Sudakov, otherwise done.
804  if (bProfile != 4) break;
805  else {
806  // Save details of the original hard interaction.
807  pT2Save = pT2;
808  id1Save = id1;
809  id2Save = id2;
810  x1Save = x1;
811  x2Save = x2;
812  sHatSave = sHat;
813  tHatSave = tHat;
814  uHatSave = uHat;
815  alpSsave = alpS;
816  alpEMsave = alpEM;
817  pT2FacSave = pT2Fac;
818  pT2RenSave = pT2Ren;
819  xPDF1nowSave = xPDF1now;
820  xPDF2nowSave = xPDF2now;
821  // Save accepted kinematics and pointer to SigmaProcess.
822  dSigmaDtSel->saveKin();
823  dSigmaDtSelSave = dSigmaDtSel;
824 
825  // Put x1, x2 information into beam pointers to get correct
826  // PDF rescaling in trial interaction (for hard process, can
827  // be sea or valence, not companion).
828  beamAPtr->append( 0, id1, x1);
829  beamAPtr->xfISR( 0, id1, x1, pT2Fac * pT2Fac);
830  vsc1 = beamAPtr->pickValSeaComp();
831  beamBPtr->append( 0, id2, x2);
832  beamBPtr->xfISR( 0, id2, x2, pT2Fac * pT2Fac);
833  vsc2 = beamBPtr->pickValSeaComp();
834 
835  // Pick b according to O(b, x1, x2).
836  double w1 = XDEP_A1 + a1 * log(1. / x1);
837  double w2 = XDEP_A1 + a1 * log(1. / x2);
838  double fac = a02now * (w1 * w1 + w2 * w2);
839  double expb2 = rndmPtr->flat();
840  b2now = - fac * log(expb2);
841  bNow = sqrt(b2now);
842 
843  // Enhancement factor for the hard process and overestimate
844  // for fastPT2. Note that existing framework has a (1. / sigmaND)
845  // present.
846  enhanceB = sigmaND / M_PI / fac * expb2;
847  enhanceBmax = sigmaND / 2. / M_PI / a02now *
848  exp( -b2now / 2. / a2max );
849 
850  // Trial interaction with dummy event.
851  Event evDummy;
852  double pTtrial = pTnext(pTmax, pTmin, evDummy);
853 
854  // Restore beams.
855  beamAPtr->clear();
856  beamBPtr->clear();
857 
858  // Accept if fallen beneath factorisation scale.
859  if (pTtrial < sqrt(pT2FacSave)) {
860  // Restore previous values and original sigma.
861  swap(pT2, pT2Save);
862  swap(pT2Fac, pT2FacSave);
863  swap(pT2Ren, pT2RenSave);
864  swap(id1, id1Save);
865  swap(id2, id2Save);
866  swap(x1, x1Save);
867  swap(x2, x2Save);
868  swap(sHat, sHatSave);
869  swap(tHat, tHatSave);
870  swap(uHat, uHatSave);
871  swap(alpS, alpSsave);
872  swap(alpEM, alpEMsave);
873  swap(xPDF1now, xPDF1nowSave);
874  swap(xPDF2now, xPDF2nowSave);
875  if (dSigmaDtSel == dSigmaDtSelSave) dSigmaDtSel->swapKin();
876  else swap(dSigmaDtSel, dSigmaDtSelSave);
877 
878  // Accept.
879  bIsSet = true;
880  break;
881  }
882  } // if (bProfile == 4)
883  } // while (true)
884 
885  // End handling for high b.
886  }
887 
888 }
889 
890 //--------------------------------------------------------------------------
891 
892 // Set up kinematics for first = hardest pT in nondiffractive process.
893 
894 void MultipartonInteractions::setupFirstSys( Event& process) {
895 
896  // Last beam-status particles. Offset relative to normal beam locations.
897  int sizeProc = process.size();
898  int nBeams = 3;
899  for (int i = 3; i < sizeProc; ++i)
900  if (process[i].statusAbs() < 20) nBeams = i + 1;
901  int nOffset = nBeams - 3;
902 
903  // Remove any partons of previous failed interactions.
904  if (sizeProc > nBeams) {
905  process.popBack( sizeProc - nBeams);
906  process.initColTag();
907  }
908 
909  // Entries 3 and 4, now to be added, come from 1 and 2.
910  process[1 + nOffset].daughter1(3 + nOffset);
911  process[2 + nOffset].daughter1(4 + nOffset);
912 
913  // Negate beam status, if not already done. (Case with offset beams.)
914  process[1 + nOffset].statusNeg();
915  process[2 + nOffset].statusNeg();
916 
917  // Loop over four partons and offset info relative to subprocess itself.
918  int colOffset = process.lastColTag();
919  for (int i = 1; i <= 4; ++i) {
920  Particle parton = dSigmaDtSel->getParton(i);
921  if (i <= 2) parton.status(-21);
922  else parton.status(23);
923  if (i <= 2) parton.mothers( i + nOffset, 0);
924  else parton.mothers( 3 + nOffset, 4 + nOffset);
925  if (i <= 2) parton.daughters( 5 + nOffset, 6 + nOffset);
926  else parton.daughters( 0, 0);
927  int col = parton.col();
928  if (col > 0) parton.col( col + colOffset);
929  int acol = parton.acol();
930  if (acol > 0) parton.acol( acol + colOffset);
931 
932  // Put the partons into the event record.
933  process.append(parton);
934  }
935 
936  // Set scale from which to begin evolution.
937  process.scale( sqrt(pT2Fac) );
938 
939  // Info on subprocess - specific to mimimum-bias events.
940  string nameSub = dSigmaDtSel->name();
941  int codeSub = dSigmaDtSel->code();
942  int nFinalSub = dSigmaDtSel->nFinal();
943  double pTMPI = dSigmaDtSel->pTMPIFin();
944  infoPtr->setSubType( iDiffSys, nameSub, codeSub, nFinalSub);
945  if (iDiffSys == 0) infoPtr->setTypeMPI( codeSub, pTMPI, 0, 0,
946  enhanceB / zeroIntCorr);
947 
948  // Further standard info on process.
949  infoPtr->setPDFalpha( iDiffSys, id1, id2, x1, x2, xPDF1now, xPDF2now,
950  pT2Fac, alpEM, alpS, pT2Ren, 0.);
951  double m3 = dSigmaDtSel->m(3);
952  double m4 = dSigmaDtSel->m(4);
953  double theta = dSigmaDtSel->thetaMPI();
954  double phi = dSigmaDtSel->phiMPI();
955  infoPtr->setKin( iDiffSys, id1, id2, x1, x2, sHat, tHat, uHat, sqrt(pT2),
956  m3, m4, theta, phi);
957 
958 }
959 
960 //--------------------------------------------------------------------------
961 
962 // Find whether to limit maximum scale of emissions.
963 
964 bool MultipartonInteractions::limitPTmax( Event& event) {
965 
966  // User-set cases.
967  if (pTmaxMatch == 1) return true;
968  if (pTmaxMatch == 2) return false;
969 
970  // Always restrict SoftQCD processes.
971  if (infoPtr->isNonDiffractive() || infoPtr->isDiffractiveA()
972  || infoPtr->isDiffractiveB() || infoPtr->isDiffractiveC() ) return true;
973 
974  // Look if only quarks (u, d, s, c, b), gluons and photons in final state.
975  bool onlyQGP1 = true;
976  bool onlyQGP2 = true;
977  int n21 = 0;
978  for (int i = 5; i < event.size(); ++i) {
979  if (event[i].status() == -21) ++n21;
980  else if (n21 == 0) {
981  int idAbs = event[i].idAbs();
982  if (idAbs > 5 && idAbs != 21 && idAbs != 22) onlyQGP1 = false;
983  } else if (n21 == 2) {
984  int idAbs = event[i].idAbs();
985  if (idAbs > 5 && idAbs != 21 && idAbs != 22) onlyQGP2 = false;
986  }
987  }
988 
989  // If two hard interactions then limit if one only contains q/g/gamma.
990  bool onlyQGP = (n21 == 2) ? (onlyQGP1 || onlyQGP2) : onlyQGP1;
991  return (onlyQGP);
992 
993 }
994 
995 //--------------------------------------------------------------------------
996 
997 // Select next pT in downwards evolution.
998 
999 double MultipartonInteractions::pTnext( double pTbegAll, double pTendAll,
1000  Event& event) {
1001 
1002  // Initial values.
1003  bool pickRescatter, acceptKin;
1004  double dSigmaScatter, dSigmaRescatter, WTacc;
1005  double pT2end = pow2( max(pTmin, pTendAll) );
1006 
1007  // With the x-dependent matter profile, and minimum bias or diffraction,
1008  // it is possible to reuse values that have been stored during the trial
1009  // interactions for a slight speedup.
1010  // bIsSet is false during trial interactions, counter 21 in case partonLevel
1011  // is retried and counter 22 for the first pass through partonLevel.
1012  if (bProfile == 4 && bIsSet && bSetInFirst && infoPtr->getCounter(21) == 1
1013  && infoPtr->getCounter(22) == 1) {
1014  if (pT2Save < pT2end) return 0.;
1015  pT2 = pT2Save;
1016  pT2Fac = pT2FacSave;
1017  pT2Ren = pT2RenSave;
1018  id1 = id1Save;
1019  id2 = id2Save;
1020  x1 = x1Save;
1021  x2 = x2Save;
1022  sHat = sHatSave;
1023  tHat = tHatSave;
1024  uHat = uHatSave;
1025  alpS = alpSsave;
1026  alpEM = alpEMsave;
1027  xPDF1now = xPDF1nowSave;
1028  xPDF2now = xPDF2nowSave;
1029  if (dSigmaDtSel == dSigmaDtSelSave) dSigmaDtSel->swapKin();
1030  else dSigmaDtSel = dSigmaDtSelSave;
1031  return sqrt(pT2);
1032  }
1033 
1034  // Do not allow rescattering while still FSR with global recoil.
1035  bool allowRescatterNow = allowRescatter;
1036  if (globalRecoilFSR && partonSystemsPtr->sizeOut(0) <= nMaxGlobalRecoilFSR)
1037  allowRescatterNow = false;
1038 
1039  // Initial pT2 value.
1040  pT2 = pow2(pTbegAll);
1041 
1042  // Find the set of already scattered partons on the two sides.
1043  if (allowRescatterNow) findScatteredPartons( event);
1044 
1045  // Pick a pT2 using a quick-and-dirty cross section estimate.
1046  do {
1047  do {
1048  pT2 = fastPT2(pT2);
1049  if (pT2 < pT2end) return 0.;
1050 
1051  // Initial values: no rescattering.
1052  i1Sel = 0;
1053  i2Sel = 0;
1054  dSigmaSum = 0.;
1055  pickRescatter = false;
1056 
1057  // Pick complete kinematics and evaluate interaction cross-section.
1058  dSigmaScatter = sigmaPT2scatter(false);
1059 
1060  // Also cross section from rescattering if allowed.
1061  dSigmaRescatter = (allowRescatterNow) ? sigmaPT2rescatter( event) : 0.;
1062 
1063  // Normalize to dSigmaApprox, which was set in fastPT2 above.
1064  WTacc = (dSigmaScatter + dSigmaRescatter) / dSigmaApprox;
1065  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
1066  "MultipartonInteractions::pTnext: weight above unity");
1067 
1068  // Idea suggested by Gosta Gustafson: increased screening in events
1069  // with large activity can be simulated by pT0_eff = sqrt(n) * pT0.
1070  if (enhanceScreening > 0) {
1071  int nSysNow = infoPtr->nMPI() + 1;
1072  if (enhanceScreening == 2) nSysNow += infoPtr->nISR();
1073  double WTscreen = pow2( (pT2 + pT20) / (pT2 + nSysNow * pT20) );
1074  WTacc *= WTscreen;
1075  }
1076 
1077  // x-dependent matter profile overlap weighting.
1078  if (bProfile == 4) {
1079  double w1 = XDEP_A1 + a1 * log(1. / x1);
1080  double w2 = XDEP_A1 + a1 * log(1. / x2);
1081  double fac = a02now * (w1 * w1 + w2 * w2);
1082  // Correct enhancement factor and weight
1083  enhanceBnow = sigmaND / M_PI / fac * exp( - b2now / fac);
1084  double oWgt = enhanceBnow / enhanceBmax;
1085  if (oWgt > 1.0000000001) infoPtr->errorMsg("Warning in Multiparton"
1086  "Interactions::pTnext: overlap weight above unity");
1087  WTacc *= oWgt;
1088  }
1089 
1090  // Decide whether to keep the event based on weight.
1091  } while (WTacc < rndmPtr->flat());
1092 
1093  // When rescattering possible: new interaction or rescattering?
1094  if (allowRescatterNow) {
1095  pickRescatter = (i1Sel > 0 || i2Sel > 0);
1096 
1097  // Restore kinematics for selected scattering/rescattering.
1098  id1 = id1Sel;
1099  id2 = id2Sel;
1100  x1 = x1Sel;
1101  x2 = x2Sel;
1102  sHat = sHatSel;
1103  tHat = tHatSel;
1104  uHat = uHatSel;
1105  sigma2Sel->sigma( id1, id2, x1, x2, sHat, tHat, uHat, alpS, alpEM,
1106  true, pickOtherSel);
1107  }
1108 
1109  // Pick one of the possible channels summed above.
1110  dSigmaDtSel = sigma2Sel->sigmaSel();
1111  if (sigma2Sel->swapTU()) swap( tHat, uHat);
1112 
1113  // Decide to keep event based on kinematics (normally always OK).
1114  // Rescattering: need to provide incoming four-vectors and masses.
1115  if (pickRescatter) {
1116  Vec4 p1Res = (i1Sel == 0) ? 0.5 * eCM * x1Sel * Vec4( 0., 0., 1., 1.)
1117  : event[i1Sel].p();
1118  Vec4 p2Res = (i2Sel == 0) ? 0.5 * eCM * x2Sel * Vec4( 0., 0., -1., 1.)
1119  : event[i2Sel].p();
1120  double m1Res = (i1Sel == 0) ? 0. : event[i1Sel].m();
1121  double m2Res = (i2Sel == 0) ? 0. : event[i2Sel].m();
1122  acceptKin = dSigmaDtSel->final2KinMPI( i1Sel, i2Sel, p1Res, p2Res,
1123  m1Res, m2Res);
1124  // New interaction: already stored values enough.
1125  } else acceptKin = dSigmaDtSel->final2KinMPI();
1126  } while (!acceptKin);
1127 
1128  // Done.
1129  return sqrt(pT2);
1130 
1131 }
1132 
1133 //--------------------------------------------------------------------------
1134 
1135 // Set up the kinematics of the 2 -> 2 scattering process,
1136 // and store the scattering in the event record.
1137 
1138 bool MultipartonInteractions::scatter( Event& event) {
1139 
1140  // Last beam-status particles. Offset relative to normal beam locations.
1141  int sizeProc = event.size();
1142  int nBeams = 3;
1143  for (int i = 3; i < sizeProc; ++i)
1144  if (event[i].statusAbs() < 20) nBeams = i + 1;
1145  int nOffset = nBeams - 3;
1146 
1147  // Loop over four partons and offset info relative to subprocess itself.
1148  int motherOffset = event.size();
1149  int colOffset = event.lastColTag();
1150  for (int i = 1; i <= 4; ++i) {
1151  Particle parton = dSigmaDtSel->getParton(i);
1152  if (i <= 2 ) parton.mothers( i + nOffset, 0);
1153  else parton.mothers( motherOffset, motherOffset + 1);
1154  if (i <= 2 ) parton.daughters( motherOffset + 2, motherOffset + 3);
1155  else parton.daughters( 0, 0);
1156  int col = parton.col();
1157  if (col > 0) parton.col( col + colOffset);
1158  int acol = parton.acol();
1159  if (acol > 0) parton.acol( acol + colOffset);
1160 
1161  // Put the partons into the event record.
1162  event.append(parton);
1163  }
1164 
1165  // Allow veto of MPI. If so restore event record to before scatter.
1166  if (canVetoMPI && userHooksPtr->doVetoMPIEmission(sizeProc, event)) {
1167  event.popBack(event.size() - sizeProc);
1168  return false;
1169  }
1170 
1171  // Store participating partons as a new set in list of all systems.
1172  int iSys = partonSystemsPtr->addSys();
1173  partonSystemsPtr->setInA(iSys, motherOffset);
1174  partonSystemsPtr->setInB(iSys, motherOffset + 1);
1175  partonSystemsPtr->addOut(iSys, motherOffset + 2);
1176  partonSystemsPtr->addOut(iSys, motherOffset + 3);
1177  partonSystemsPtr->setSHat(iSys, sHat);
1178 
1179  // Tag double rescattering graphs that annihilate one initial colour.
1180  bool annihil1 = false;
1181  bool annihil2 = false;
1182  if (i1Sel > 0 && i2Sel > 0) {
1183  if (event[motherOffset].col() == event[motherOffset + 1].acol()
1184  && event[motherOffset].col() > 0) annihil1 = true;
1185  if (event[motherOffset].acol() == event[motherOffset + 1].col()
1186  && event[motherOffset].acol() > 0) annihil2 = true;
1187  }
1188 
1189  // Beam remnant A: add scattered partons to list.
1190  BeamParticle& beamA = *beamAPtr;
1191  int iA = beamA.append( motherOffset, id1, x1);
1192 
1193  // Find whether incoming partons are valence or sea, so prepared for ISR.
1194  if (i1Sel == 0) {
1195  beamA.xfISR( iA, id1, x1, pT2);
1196  beamA.pickValSeaComp();
1197 
1198  // Remove rescattered parton from final state and change history.
1199  // Propagate existing colour labels throught graph.
1200  } else {
1201  beamA[iA].companion(-10);
1202  event[i1Sel].statusNeg();
1203  event[i1Sel].daughters( motherOffset, motherOffset);
1204  event[motherOffset].mothers( i1Sel, i1Sel);
1205  int colOld = event[i1Sel].col();
1206  if (colOld > 0) {
1207  int colNew = event[motherOffset].col();
1208  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1209  if (event[i].col() == colNew) event[i].col( colOld);
1210  if (event[i].acol() == colNew) event[i].acol( colOld);
1211  }
1212  }
1213  int acolOld = event[i1Sel].acol();
1214  if (acolOld > 0) {
1215  int acolNew = event[motherOffset].acol();
1216  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1217  if (event[i].col() == acolNew) event[i].col( acolOld);
1218  if (event[i].acol() == acolNew) event[i].acol( acolOld);
1219  }
1220  }
1221  }
1222 
1223  // Beam remnant B: add scattered partons to list.
1224  BeamParticle& beamB = *beamBPtr;
1225  int iB = beamB.append( motherOffset + 1, id2, x2);
1226 
1227  // Find whether incoming partons are valence or sea, so prepared for ISR.
1228  if (i2Sel == 0) {
1229  beamB.xfISR( iB, id2, x2, pT2);
1230  beamB.pickValSeaComp();
1231 
1232  // Remove rescattered parton from final state and change history.
1233  // Propagate existing colour labels throught graph.
1234  } else {
1235  beamB[iB].companion(-10);
1236  event[i2Sel].statusNeg();
1237  event[i2Sel].daughters( motherOffset + 1, motherOffset + 1);
1238  event[motherOffset + 1].mothers( i2Sel, i2Sel);
1239  int colOld = event[i2Sel].col();
1240  if (colOld > 0) {
1241  int colNew = event[motherOffset + 1].col();
1242  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1243  if (event[i].col() == colNew) event[i].col( colOld);
1244  if (event[i].acol() == colNew) event[i].acol( colOld);
1245  }
1246  }
1247  int acolOld = event[i2Sel].acol();
1248  if (acolOld > 0) {
1249  int acolNew = event[motherOffset + 1].acol();
1250  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1251  if (event[i].col() == acolNew) event[i].col( acolOld);
1252  if (event[i].acol() == acolNew) event[i].acol( acolOld);
1253  }
1254  }
1255  }
1256 
1257  // Annihilating colour in double rescattering requires relabelling
1258  // of one colour into the other in the whole preceding event.
1259  if (annihil1 || annihil2) {
1260  int colLeft = (annihil1) ? event[motherOffset].col()
1261  : event[motherOffset].acol();
1262  int mother1 = event[motherOffset].mother1();
1263  int mother2 = event[motherOffset + 1].mother1();
1264  int colLost = (annihil1)
1265  ? event[mother1].col() + event[mother2].acol() - colLeft
1266  : event[mother1].acol() + event[mother2].col() - colLeft;
1267  for (int i = 0; i < motherOffset; ++i) {
1268  if (event[i].col() == colLost) event[i].col( colLeft );
1269  if (event[i].acol() == colLost) event[i].acol( colLeft );
1270  }
1271  }
1272 
1273  // Store info on subprocess code and rescattered partons.
1274  int codeMPI = dSigmaDtSel->code();
1275  double pTMPI = dSigmaDtSel->pTMPIFin();
1276  if (iDiffSys == 0) infoPtr->setTypeMPI( codeMPI, pTMPI, i1Sel, i2Sel,
1277  enhanceBnow / zeroIntCorr);
1278  partonSystemsPtr->setPTHat(iSys, pTMPI);
1279 
1280  // Done.
1281  return true;
1282 }
1283 
1284 //--------------------------------------------------------------------------
1285 
1286 // Determine constant in d(Prob)/d(pT2) < const / (pT2 + r * pT20)^2.
1287 
1288 void MultipartonInteractions::upperEnvelope() {
1289 
1290  // Initially determine constant in jet cross section upper estimate
1291  // d(sigma_approx)/d(pT2) < const / (pT2 + r * pT20)^2.
1292  pT4dSigmaMax = 0.;
1293 
1294  // Loop thorough allowed pT range logarithmically evenly.
1295  for (int iPT = 0; iPT < 100; ++iPT) {
1296  double pT = pTmin * pow( pTmax/pTmin, 0.01 * (iPT + 0.5) );
1297  pT2 = pT*pT;
1298  pT2shift = pT2 + pT20;
1299  pT2Ren = pT2shift;
1300  pT2Fac = (SHIFTFACSCALE) ? pT2shift : pT2;
1301  xT = 2. * pT / eCM;
1302 
1303  // Evaluate parton density sums at x1 = x2 = xT.
1304  double xPDF1sumMax = (9./4.) * beamAPtr->xf(21, xT, pT2Fac);
1305  for (int id = 1; id <= nQuarkIn; ++id)
1306  xPDF1sumMax += beamAPtr->xf( id, xT, pT2Fac)
1307  + beamAPtr->xf(-id, xT, pT2Fac);
1308  double xPDF2sumMax = (9./4.) * beamBPtr->xf(21, xT, pT2Fac);
1309  for (int id = 1; id <= nQuarkIn; ++id)
1310  xPDF2sumMax += beamBPtr->xf( id, xT, pT2Fac)
1311  + beamBPtr->xf(-id, xT, pT2Fac);
1312 
1313  // Evaluate alpha_strong and _EM, matrix element and phase space volume.
1314  alpS = alphaS.alphaS(pT2Ren);
1315  alpEM = alphaEM.alphaEM(pT2Ren);
1316  double dSigmaPartonApprox = CONVERT2MB * Kfactor * 0.5 * M_PI
1317  * pow2(alpS / pT2shift);
1318  double yMax = log(1./xT + sqrt(1./(xT*xT) - 1.));
1319  double volumePhSp = pow2(2. * yMax);
1320 
1321  // Final comparison to determine upper estimate.
1322  double dSigmaApproxNow = SIGMAFUDGE * xPDF1sumMax * xPDF2sumMax
1323  * dSigmaPartonApprox * volumePhSp;
1324  double pT4dSigmaNow = pow2(pT2 + pT20R) * dSigmaApproxNow;
1325  if ( pT4dSigmaNow > pT4dSigmaMax) pT4dSigmaMax = pT4dSigmaNow;
1326  }
1327 
1328  // Get wanted constant by dividing by the nondiffractive cross section.
1329  pT4dProbMax = pT4dSigmaMax / sigmaND;
1330 
1331 }
1332 
1333 //--------------------------------------------------------------------------
1334 
1335 // Integrate the parton-parton interaction cross section,
1336 // using stratified Monte Carlo sampling.
1337 // Store result in pT bins for use as Sudakov form factors.
1338 
1339 void MultipartonInteractions::jetCrossSection() {
1340 
1341  // Common factor from bin size in dpT2 / (pT2 + r * pT20)^2 and statistics.
1342  double sigmaFactor = (1. / pT20minR - 1. / pT20maxR) / (100. * nSample);
1343 
1344  // Reset overlap-weighted cross section for x-dependent matter profile.
1345  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++)
1346  sigmaIntWgt[bBin] = 0.;
1347 
1348  // Loop through allowed pT range evenly in dpT2 / (pT2 + r * pT20)^2.
1349  sigmaInt = 0.;
1350  double dSigmaMax = 0.;
1351  sudExpPT[100] = 0.;
1352 
1353  for (int iPT = 99; iPT >= 0; --iPT) {
1354  double sigmaSum = 0.;
1355 
1356  // Reset pT-binned overlap-weigted integration.
1357  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++)
1358  sigmaSumWgt[bBin] = 0.;
1359 
1360  // In each pT bin sample a number of random pT values.
1361  for (int iSample = 0; iSample < nSample; ++iSample) {
1362  double mappedPT2 = 1. - 0.01 * (iPT + rndmPtr->flat());
1363  pT2 = pT20min0maxR / (pT20minR + mappedPT2 * pT2maxmin) - pT20R;
1364 
1365  // Evaluate cross section dSigma/dpT2 in phase space point.
1366  double dSigma = sigmaPT2scatter(true);
1367 
1368  // Multiply by (pT2 + r * pT20)^2 to compensate for pT sampling. Sum.
1369  dSigma *= pow2(pT2 + pT20R);
1370  sigmaSum += dSigma;
1371  if (dSigma > dSigmaMax) dSigmaMax = dSigma;
1372 
1373  // Overlap-weighted cross section for x-dependent matter profile.
1374  // Note that dSigma can be 0. when points are rejected.
1375  if (bProfile == 4 && dSigma > 0.) {
1376  double w1 = XDEP_A1 + a1 * log(1. / x1);
1377  double w2 = XDEP_A1 + a1 * log(1. / x2);
1378  double fac = XDEP_A0 * XDEP_A0 * (w1 * w1 + w2 * w2);
1379  double b = 0.5 * bstepNow;
1380  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
1381  double wgt = exp( - b * b / fac ) / fac / M_PI;
1382  sigmaSumWgt[bBin] += dSigma * wgt;
1383  b += bstepNow;
1384  }
1385  }
1386  }
1387 
1388  // Store total cross section and exponent of Sudakov.
1389  sigmaSum *= sigmaFactor;
1390  sigmaInt += sigmaSum;
1391  sudExpPT[iPT] = sudExpPT[iPT + 1] + sigmaSum / sigmaND;
1392 
1393  // Sum overlap-weighted cross section
1394  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
1395  sigmaSumWgt[bBin] *= sigmaFactor;
1396  sigmaIntWgt[bBin] += sigmaSumWgt[bBin];
1397  }
1398 
1399  // End of loop over pT values.
1400  }
1401 
1402  // Update upper estimate of differential cross section. Done.
1403  if (dSigmaMax > pT4dSigmaMax) {
1404  pT4dSigmaMax = dSigmaMax;
1405  pT4dProbMax = dSigmaMax / sigmaND;
1406  }
1407 
1408 }
1409 
1410 //--------------------------------------------------------------------------
1411 
1412 // Evaluate "Sudakov form factor" for not having a harder interaction
1413 // at the selected b value, given the pT scale of the event.
1414 
1415 double MultipartonInteractions::sudakov(double pT2sud, double enhance) {
1416 
1417  // Find bin the pT2 scale falls in.
1418  double xBin = (pT2sud - pT2min) * pT20maxR
1419  / (pT2maxmin * (pT2sud + pT20R));
1420  xBin = max(1e-6, min(100. - 1e-6, 100. * xBin) );
1421  int iBin = int(xBin);
1422 
1423  // Interpolate inside bin. Optionally include enhancement factor.
1424  double sudExp = sudExpPT[iBin] + (xBin - iBin)
1425  * (sudExpPT[iBin + 1] - sudExpPT[iBin]);
1426  return exp( -enhance * sudExp);
1427 
1428 }
1429 
1430 //--------------------------------------------------------------------------
1431 
1432 // Pick a trial next pT, based on a simple upper estimate of the
1433 // d(sigma)/d(pT2) spectrum.
1434 
1435 double MultipartonInteractions::fastPT2( double pT2beg) {
1436 
1437  // Use d(Prob)/d(pT2) < pT4dProbMax / (pT2 + r * pT20)^2.
1438  double pT20begR = pT2beg + pT20R;
1439  double pT4dProbMaxNow = pT4dProbMax * enhanceBmax;
1440  double pT2try = pT4dProbMaxNow * pT20begR
1441  / (pT4dProbMaxNow - pT20begR * log(rndmPtr->flat())) - pT20R;
1442 
1443  // Save cross section associated with ansatz above. Done.
1444  dSigmaApprox = pT4dSigmaMax / pow2(pT2try + pT20R);
1445  return pT2try;
1446 
1447 }
1448 
1449 //--------------------------------------------------------------------------
1450 
1451 // Calculate the actual cross section to decide whether fast choice is OK.
1452 // Select flavours and kinematics for interaction at given pT.
1453 // Slightly different treatment for first interaction and subsequent ones.
1454 
1455 double MultipartonInteractions::sigmaPT2scatter(bool isFirst) {
1456 
1457  // Derive recormalization and factorization scales, amd alpha_strong/em.
1458  pT2shift = pT2 + pT20;
1459  pT2Ren = pT2shift;
1460  pT2Fac = (SHIFTFACSCALE) ? pT2shift : pT2;
1461  alpS = alphaS.alphaS(pT2Ren);
1462  alpEM = alphaEM.alphaEM(pT2Ren);
1463 
1464  // Derive rapidity limits from chosen pT2.
1465  xT = 2. * sqrt(pT2) / eCM;
1466  if (xT >= 1.) return 0.;
1467  xT2 = xT*xT;
1468  double yMax = log(1./xT + sqrt(1./xT2 - 1.));
1469 
1470  // Select rapidities y3 and y4 of the two produced partons.
1471  double y3 = yMax * (2. * rndmPtr->flat() - 1.);
1472  double y4 = yMax * (2. * rndmPtr->flat() - 1.);
1473  y = 0.5 * (y3 + y4);
1474 
1475  // Reject some events at large rapidities to improve efficiency.
1476  // (Works for baryons, not pions or Pomerons if they have hard PDF's.)
1477  double WTy = (hasBaryonBeams)
1478  ? (1. - pow2(y3/yMax)) * (1. - pow2(y4/yMax)) : 1.;
1479  if (WTy < rndmPtr->flat()) return 0.;
1480 
1481  // Failure if x1 or x2 exceed what is left in respective beam.
1482  x1 = 0.5 * xT * (exp(y3) + exp(y4));
1483  x2 = 0.5 * xT * (exp(-y3) + exp(-y4));
1484  if (isFirst) {
1485  if (x1 > 1. || x2 > 1.) return 0.;
1486  } else {
1487  if (x1 > beamAPtr->xMax() || x2 > beamBPtr->xMax()) return 0.;
1488  }
1489  tau = x1 * x2;
1490 
1491  // Begin evaluate parton densities at actual x1 and x2.
1492  double xPDF1[21];
1493  double xPDF1sum = 0.;
1494  double xPDF2[21];
1495  double xPDF2sum = 0.;
1496 
1497  // For first interaction use normal densities.
1498  if (isFirst) {
1499  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1500  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xf(21, x1, pT2Fac);
1501  else xPDF1[id+10] = beamAPtr->xf(id, x1, pT2Fac);
1502  xPDF1sum += xPDF1[id+10];
1503  }
1504  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1505  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xf(21, x2, pT2Fac);
1506  else xPDF2[id+10] = beamBPtr->xf(id, x2, pT2Fac);
1507  xPDF2sum += xPDF2[id+10];
1508  }
1509 
1510  // For subsequent interactions use rescaled densities.
1511  } else {
1512  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1513  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xfMPI(21, x1, pT2Fac);
1514  else xPDF1[id+10] = beamAPtr->xfMPI(id, x1, pT2Fac);
1515  xPDF1sum += xPDF1[id+10];
1516  }
1517  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1518  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xfMPI(21, x2, pT2Fac);
1519  else xPDF2[id+10] = beamBPtr->xfMPI(id, x2, pT2Fac);
1520  xPDF2sum += xPDF2[id+10];
1521  }
1522  }
1523 
1524  // Select incoming flavours according to actual PDF's.
1525  id1 = -nQuarkIn - 1;
1526  double temp = xPDF1sum * rndmPtr->flat();
1527  do { xPDF1now = xPDF1[(++id1) + 10]; temp -= xPDF1now; }
1528  while (temp > 0. && id1 < nQuarkIn);
1529  if (id1 == 0) id1 = 21;
1530  id2 = -nQuarkIn-1;
1531  temp = xPDF2sum * rndmPtr->flat();
1532  do { xPDF2now = xPDF2[(++id2) + 10]; temp -= xPDF2now;}
1533  while (temp > 0. && id2 < nQuarkIn);
1534  if (id2 == 0) id2 = 21;
1535 
1536  // Assign pointers to processes relevant for incoming flavour choice:
1537  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1538  // Factor 4./9. per incoming gluon to compensate for preweighting.
1539  SigmaMultiparton* sigma2Tmp;
1540  double gluFac = 1.;
1541  if (id1 == 21 && id2 == 21) {
1542  sigma2Tmp = &sigma2gg;
1543  gluFac = 16. / 81.;
1544  } else if (id1 == 21 || id2 == 21) {
1545  sigma2Tmp = &sigma2qg;
1546  gluFac = 4. / 9.;
1547  } else if (id1 == -id2) sigma2Tmp = &sigma2qqbarSame;
1548  else sigma2Tmp = &sigma2qq;
1549 
1550  // Prepare to generate differential cross sections.
1551  sHat = tau * sCM;
1552  double root = sqrtpos(1. - xT2 / tau);
1553  tHat = -0.5 * sHat * (1. - root);
1554  uHat = -0.5 * sHat * (1. + root);
1555 
1556  // Evaluate cross sections, include possibility of K factor.
1557  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1558  // (Not necessary, but makes for better MC efficiency.)
1559  double dSigmaPartonCorr = Kfactor * gluFac
1560  * sigma2Tmp->sigma( id1, id2, x1, x2, sHat, tHat, uHat, alpS, alpEM);
1561 
1562  // Combine cross section, pdf's and phase space integral.
1563  double volumePhSp = pow2(2. * yMax) / WTy;
1564  double dSigmaScat = dSigmaPartonCorr * xPDF1sum * xPDF2sum * volumePhSp;
1565 
1566  // Dampen cross section at small pT values; part of formalism.
1567  // Use pT2 corrected for massive kinematics at this step.??
1568  // double pT2massive = dSigmaDtSel->pT2MPI();
1569  double pT2massive = pT2;
1570  dSigmaScat *= pow2( pT2massive / (pT2massive + pT20) );
1571 
1572  // Sum up total contribution for all scatterings and rescatterings.
1573  dSigmaSum += dSigmaScat;
1574 
1575  // Save values for comparison with rescattering processes.
1576  i1Sel = 0;
1577  i2Sel = 0;
1578  id1Sel = id1;
1579  id2Sel = id2;
1580  x1Sel = x1;
1581  x2Sel = x2;
1582  sHatSel = sHat;
1583  tHatSel = tHat;
1584  uHatSel = uHat;
1585  sigma2Sel = sigma2Tmp;
1586  pickOtherSel = sigma2Tmp->pickedOther();
1587 
1588  // For first interaction: pick one of the possible channels summed above.
1589  if (isFirst) {
1590  dSigmaDtSel = sigma2Tmp->sigmaSel();
1591  if (sigma2Tmp->swapTU()) swap( tHat, uHat);
1592  }
1593 
1594  // Done.
1595  return dSigmaScat;
1596 }
1597 
1598 //--------------------------------------------------------------------------
1599 
1600 // Find the partons that are allowed to rescatter.
1601 
1602 void MultipartonInteractions::findScatteredPartons( Event& event) {
1603 
1604  // Reset arrays.
1605  scatteredA.resize(0);
1606  scatteredB.resize(0);
1607  double yTmp, probA, probB;
1608 
1609  // Loop though the event record and catch "final" partons.
1610  for (int i = 0; i < event.size(); ++i)
1611  if (event[i].isFinal() && (event[i].idAbs() <= nQuarkIn
1612  || event[i].id() == 21)) {
1613  yTmp = event[i].y();
1614 
1615  // Different strategies to determine which partons may rescatter.
1616  switch(rescatterMode) {
1617 
1618  // Case 0: step function at origin
1619  case 0:
1620  if ( yTmp > 0.) scatteredA.push_back( i);
1621  if (-yTmp > 0.) scatteredB.push_back( i);
1622  break;
1623 
1624  // Case 1: step function as ySepResc.
1625  case 1:
1626  if ( yTmp > ySepResc) scatteredA.push_back( i);
1627  if (-yTmp > ySepResc) scatteredB.push_back( i);
1628  break;
1629 
1630  // Case 2: linear rise from ySep - deltaY to ySep + deltaY.
1631  case 2:
1632  probA = 0.5 * (1. + ( yTmp - ySepResc) / deltaYResc);
1633  if (probA > rndmPtr->flat()) scatteredA.push_back( i);
1634  probB = 0.5 * (1. + (-yTmp - ySepResc) / deltaYResc);
1635  if (probB > rndmPtr->flat()) scatteredB.push_back( i);
1636  break;
1637 
1638  // Case 3: rise like (1/2) * ( 1 + tanh((y - ySep) / deltaY) ).
1639  // Use that (1/2) (1 + tanh(x)) = 1 / (1 + exp(-2x)).
1640  case 3:
1641  probA = 1. / (1. + exp(-2. * ( yTmp - ySepResc) / deltaYResc));
1642  if (probA > rndmPtr->flat()) scatteredA.push_back( i);
1643  probB = 1. / (1. + exp(-2. * (-yTmp - ySepResc) / deltaYResc));
1644  if (probB > rndmPtr->flat()) scatteredB.push_back( i);
1645  break;
1646 
1647  // Case 4 and undefined values: all partons can rescatter.
1648  default:
1649  scatteredA.push_back( i);
1650  scatteredB.push_back( i);
1651  break;
1652 
1653  // End of loop over partons. Done.
1654  }
1655  }
1656 
1657 }
1658 
1659 //--------------------------------------------------------------------------
1660 
1661 // Rescattering contribution summed over all already scattered partons.
1662 // Calculate the actual cross section to decide whether fast choice is OK.
1663 // Select flavours and kinematics for interaction at given pT.
1664 
1665 double MultipartonInteractions::sigmaPT2rescatter( Event& event) {
1666 
1667  // Derive xT scale from chosen pT2.
1668  xT = 2. * sqrt(pT2) / eCM;
1669  if (xT >= 1.) return 0.;
1670  xT2 = xT*xT;
1671 
1672  // Pointer to cross section and total rescatter contribution.
1673  SigmaMultiparton* sigma2Tmp;
1674  double dSigmaResc = 0.;
1675 
1676  // Normally save time by picking one random scattered parton from side A.
1677  int nScatA = scatteredA.size();
1678  int iScatA = -1;
1679  if (PREPICKRESCATTER && nScatA > 0) iScatA = min( nScatA - 1,
1680  int( rndmPtr->flat() * double(nScatA) ) );
1681 
1682  // Loop over all already scattered partons from side A.
1683  for (int iScat = 0; iScat < nScatA; ++iScat) {
1684  if (PREPICKRESCATTER && iScat != iScatA) continue;
1685  int iA = scatteredA[iScat];
1686  int id1Tmp = event[iA].id();
1687 
1688  // Calculate maximum possible sHat and check whether big enough.
1689  double x1Tmp = event[iA].pPos() / eCM;
1690  double sHatMax = x1Tmp * beamBPtr->xMax() * sCM;
1691  if (sHatMax < 4. * pT2) continue;
1692 
1693  // Select rapidity separation between two produced partons.
1694  double dyMax = log( sqrt(0.25 * sHatMax / pT2)
1695  * (1. + sqrtpos(1. - 4. * pT2 / sHatMax)) );
1696  double dy = dyMax * (2. * rndmPtr->flat() - 1.);
1697 
1698  // Reconstruct kinematical variables, especially x2.
1699  // Incoming c/b masses handled better if tau != x1 * x2.
1700  double m2Tmp = event[iA].m2();
1701  double sHatTmp = m2Tmp + 4. * pT2 * pow2(cosh(dy));
1702  double x2Tmp = (sHatTmp - m2Tmp) /(x1Tmp * sCM);
1703  double tauTmp = sHatTmp / sCM;
1704  double root = sqrtpos(1. - xT2 / tauTmp);
1705  double tHatTmp = -0.5 * sHatTmp * (1. - root);
1706  double uHatTmp = -0.5 * sHatTmp * (1. + root);
1707 
1708  // Begin evaluate parton densities at actual x2.
1709  double xPDF2[21];
1710  double xPDF2sum = 0.;
1711 
1712  // Use rescaled densities, with preweighting 9/4 for gluons.
1713  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1714  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xfMPI(21, x2Tmp, pT2Fac);
1715  else xPDF2[id+10] = beamBPtr->xfMPI(id, x2Tmp, pT2Fac);
1716  xPDF2sum += xPDF2[id+10];
1717  }
1718 
1719  // Select incoming flavour according to actual PDF's.
1720  int id2Tmp = -nQuarkIn - 1;
1721  double temp = xPDF2sum * rndmPtr->flat();
1722  do { xPDF2now = xPDF2[(++id2Tmp) + 10]; temp -= xPDF2now;}
1723  while (temp > 0. && id2Tmp < nQuarkIn);
1724  if (id2Tmp == 0) id2Tmp = 21;
1725 
1726  // Assign pointers to processes relevant for incoming flavour choice:
1727  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1728  // Factor 4./9. for incoming gluon 2 to compensate for preweighting.
1729  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
1730  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
1731  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
1732  else sigma2Tmp = &sigma2qq;
1733  double gluFac = (id2Tmp == 21) ? 4. / 9. : 1.;
1734 
1735  // Evaluate cross sections, include possibility of K factor.
1736  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1737  // (Not necessary, but makes for better MC efficiency.)
1738  double dSigmaPartonCorr = Kfactor * gluFac
1739  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
1740  uHatTmp, alpS, alpEM);
1741 
1742  // Combine cross section, pdf's and phase space integral.
1743  double volumePhSp = 4. *dyMax;
1744  double dSigmaCorr = dSigmaPartonCorr * xPDF2sum * volumePhSp;
1745 
1746  // Dampen cross section at small pT values; part of formalism.
1747  // Use pT2 corrected for massive kinematics at this step.
1748  //?? double pT2massive = dSigmaDtSel->pT2MPI();
1749  double pT2massive = pT2;
1750  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
1751 
1752  // Compensate for prepicked rescattering if appropriate.
1753  if (PREPICKRESCATTER) dSigmaCorr *= nScatA;
1754 
1755  // Sum up total contribution for all scatterings or rescattering only.
1756  dSigmaSum += dSigmaCorr;
1757  dSigmaResc += dSigmaCorr;
1758 
1759  // Determine whether current rescattering should be selected.
1760  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
1761  i1Sel = iA;
1762  i2Sel = 0;
1763  id1Sel = id1Tmp;
1764  id2Sel = id2Tmp;
1765  x1Sel = x1Tmp;
1766  x2Sel = x2Tmp;
1767  sHatSel = sHatTmp;
1768  tHatSel = tHatTmp;
1769  uHatSel = uHatTmp;
1770  sigma2Sel = sigma2Tmp;
1771  pickOtherSel = sigma2Tmp->pickedOther();
1772  }
1773  }
1774 
1775  // Normally save time by picking one random scattered parton from side B.
1776  int nScatB = scatteredB.size();
1777  int iScatB = -1;
1778  if (PREPICKRESCATTER && nScatB > 0) iScatB = min( nScatB - 1,
1779  int( rndmPtr->flat() * double(nScatB) ) );
1780 
1781  // Loop over all already scattered partons from side B.
1782  for (int iScat = 0; iScat < nScatB; ++iScat) {
1783  if (PREPICKRESCATTER && iScat != iScatB) continue;
1784  int iB = scatteredB[iScat];
1785  int id2Tmp = event[iB].id();
1786 
1787  // Calculate maximum possible sHat and check whether big enough.
1788  double x2Tmp = event[iB].pNeg() / eCM;
1789  double sHatMax = beamAPtr->xMax() * x2Tmp * sCM;
1790  if (sHatMax < 4. * pT2) continue;
1791 
1792  // Select rapidity separation between two produced partons.
1793  double dyMax = log( sqrt(0.25 * sHatMax / pT2)
1794  * (1. + sqrtpos(1. - 4. * pT2 / sHatMax)) );
1795  double dy = dyMax * (2. * rndmPtr->flat() - 1.);
1796 
1797  // Reconstruct kinematical variables, especially x1.
1798  // Incoming c/b masses handled better if tau != x1 * x2.
1799  double m2Tmp = event[iB].m2();
1800  double sHatTmp = m2Tmp + 4. * pT2 * pow2(cosh(dy));
1801  double x1Tmp = (sHatTmp - m2Tmp) /(x2Tmp * sCM);
1802  double tauTmp = sHatTmp / sCM;
1803  double root = sqrtpos(1. - xT2 / tauTmp);
1804  double tHatTmp = -0.5 * sHatTmp * (1. - root);
1805  double uHatTmp = -0.5 * sHatTmp * (1. + root);
1806 
1807  // Begin evaluate parton densities at actual x1.
1808  double xPDF1[21];
1809  double xPDF1sum = 0.;
1810 
1811  // Use rescaled densities, with preweighting 9/4 for gluons.
1812  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1813  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xfMPI(21, x1Tmp, pT2Fac);
1814  else xPDF1[id+10] = beamAPtr->xfMPI(id, x1Tmp, pT2Fac);
1815  xPDF1sum += xPDF1[id+10];
1816  }
1817 
1818  // Select incoming flavour according to actual PDF's.
1819  int id1Tmp = -nQuarkIn - 1;
1820  double temp = xPDF1sum * rndmPtr->flat();
1821  do { xPDF1now = xPDF1[(++id1Tmp) + 10]; temp -= xPDF1now;}
1822  while (temp > 0. && id1Tmp < nQuarkIn);
1823  if (id1Tmp == 0) id1Tmp = 21;
1824 
1825  // Assign pointers to processes relevant for incoming flavour choice:
1826  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1827  // Factor 4./9. for incoming gluon 2 to compensate for preweighting.
1828  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
1829  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
1830  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
1831  else sigma2Tmp = &sigma2qq;
1832  double gluFac = (id1Tmp == 21) ? 4. / 9. : 1.;
1833 
1834  // Evaluate cross sections, include possibility of K factor.
1835  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1836  // (Not necessary, but makes for better MC efficiency.)
1837  double dSigmaPartonCorr = Kfactor * gluFac
1838  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
1839  uHatTmp, alpS, alpEM);
1840 
1841  // Combine cross section, pdf's and phase space integral.
1842  double volumePhSp = 4. *dyMax;
1843  double dSigmaCorr = dSigmaPartonCorr * xPDF1sum * volumePhSp;
1844 
1845  // Dampen cross section at small pT values; part of formalism.
1846  // Use pT2 corrected for massive kinematics at this step.
1847  //?? double pT2massive = dSigmaDtSel->pT2MPI();
1848  double pT2massive = pT2;
1849  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
1850 
1851  // Compensate for prepicked rescattering if appropriate.
1852  if (PREPICKRESCATTER) dSigmaCorr *= nScatB;
1853 
1854  // Sum up total contribution for all scatterings or rescattering only.
1855  dSigmaSum += dSigmaCorr;
1856  dSigmaResc += dSigmaCorr;
1857 
1858  // Determine whether current rescattering should be selected.
1859  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
1860  i1Sel = 0;
1861  i2Sel = iB;
1862  id1Sel = id1Tmp;
1863  id2Sel = id2Tmp;
1864  x1Sel = x1Tmp;
1865  x2Sel = x2Tmp;
1866  sHatSel = sHatTmp;
1867  tHatSel = tHatTmp;
1868  uHatSel = uHatTmp;
1869  sigma2Sel = sigma2Tmp;
1870  pickOtherSel = sigma2Tmp->pickedOther();
1871  }
1872  }
1873 
1874  // Double rescatter: loop over already scattered partons from both sides.
1875  // As before, allow to use only one parton per side to speed up.
1876  if (allowDoubleRes) {
1877  for (int iScat1 = 0; iScat1 < nScatA; ++iScat1) {
1878  if (PREPICKRESCATTER && iScat1 != iScatA) continue;
1879  int iA = scatteredA[iScat1];
1880  int id1Tmp = event[iA].id();
1881  for (int iScat2 = 0; iScat2 < nScatB; ++iScat2) {
1882  if (PREPICKRESCATTER && iScat2 != iScatB) continue;
1883  int iB = scatteredB[iScat2];
1884  int id2Tmp = event[iB].id();
1885 
1886  // Calculate current sHat and check whether big enough.
1887  double sHatTmp = (event[iA].p() + event[iB].p()).m2Calc();
1888  if (sHatTmp < 4. * pT2) continue;
1889 
1890  // Reconstruct other kinematical variables. (x values not needed.)
1891  double x1Tmp = event[iA].pPos() / eCM;
1892  double x2Tmp = event[iB].pNeg() / eCM;
1893  double tauTmp = sHatTmp / sCM;
1894  double root = sqrtpos(1. - xT2 / tauTmp);
1895  double tHatTmp = -0.5 * sHatTmp * (1. - root);
1896  double uHatTmp = -0.5 * sHatTmp * (1. + root);
1897 
1898  // Assign pointers to processes relevant for incoming flavour choice:
1899  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1900  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
1901  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
1902  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
1903  else sigma2Tmp = &sigma2qq;
1904 
1905  // Evaluate cross sections, include possibility of K factor.
1906  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1907  // (Not necessary, but makes for better MC efficiency.)
1908  double dSigmaPartonCorr = Kfactor
1909  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
1910  uHatTmp, alpS, alpEM);
1911 
1912  // Combine cross section and Jacobian tHat -> pT2
1913  // (with safety minvalue).
1914  double dSigmaCorr = dSigmaPartonCorr / max(ROOTMIN, root);
1915 
1916  // Dampen cross section at small pT values; part of formalism.
1917  // Use pT2 corrected for massive kinematics at this step.
1918  //?? double pT2massive = dSigmaDtSel->pT2MPI();
1919  double pT2massive = pT2;
1920  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
1921 
1922  // Compensate for prepicked rescattering if appropriate.
1923  if (PREPICKRESCATTER) dSigmaCorr *= nScatA * nScatB;
1924 
1925  // Sum up total contribution for all scatterings or rescattering only.
1926  dSigmaSum += dSigmaCorr;
1927  dSigmaResc += dSigmaCorr;
1928 
1929  // Determine whether current rescattering should be selected.
1930  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
1931  i1Sel = iA;
1932  i2Sel = iB;
1933  id1Sel = id1Tmp;
1934  id2Sel = id2Tmp;
1935  x1Sel = x1Tmp;
1936  x2Sel = x2Tmp;
1937  sHatSel = sHatTmp;
1938  tHatSel = tHatTmp;
1939  uHatSel = uHatTmp;
1940  sigma2Sel = sigma2Tmp;
1941  pickOtherSel = sigma2Tmp->pickedOther();
1942  }
1943  }
1944  }
1945  }
1946 
1947  // Done.
1948  return dSigmaResc;
1949 }
1950 
1951 
1952 //--------------------------------------------------------------------------
1953 
1954 // Calculate factor relating matter overlap and interaction rate,
1955 // i.e. k in <n_interaction(b)> = k * overlap(b) (neglecting
1956 // n_int = 0 corrections and energy-momentum conservation effects).
1957 
1958 void MultipartonInteractions::overlapInit() {
1959 
1960  // Initial values for iteration. Step size of b integration.
1961  nAvg = sigmaInt / sigmaND;
1962  kNow = 0.5;
1963  int stepDir = 1;
1964  double deltaB = BSTEP;
1965  if (bProfile == 2) deltaB *= min( 0.5, 2.5 * coreRadius);
1966  if (bProfile == 3) deltaB *= max(1., pow(2. / expPow, 1. / expPow));
1967 
1968  // Further variables, with dummy initial values.
1969  double nNow = 0.;
1970  double kLow = 0.;
1971  double nLow = 0.;
1972  double kHigh = 0.;
1973  double nHigh = 0.;
1974  double overlapNow = 0.;
1975  double probNow = 0.;
1976  double overlapInt = 0.5;
1977  double probInt = 0.;
1978  double probOverlapInt = 0.;
1979  double bProbInt = 0.;
1980  normPi = 1. / (2. * M_PI);
1981 
1982  // Subdivision into low-b and high-b region by interaction rate.
1983  bool pastBDiv = false;
1984  double overlapHighB = 0.;
1985 
1986  // For x-dependent matter profile, try to find a0 rather than k.
1987  // Existing framework is used, but with substitutions:
1988  // a0 tuned according to Int( Pint(b), d^2b ) = sigmaND,
1989  // nAvg = sigmaND, kNow = a0now, kLow = a0low, kHigh = a0high,
1990  // nNow = probInt, nLow = probIntLow, nHigh = probIntHigh.
1991  double rescale2 = 1.;
1992  if (bProfile == 4) {
1993  nAvg = sigmaND;
1994  kNow = XDEP_A0 / 2.0;
1995  }
1996 
1997  // First close k into an interval by binary steps,
1998  // then find k by successive interpolation.
1999  do {
2000  if (stepDir == 1) kNow *= 2.;
2001  else if (stepDir == -1) kNow *= 0.5;
2002  else kNow = kLow + (nAvg - nLow) * (kHigh - kLow) / (nHigh - nLow);
2003 
2004  // Overlap trivial if no impact parameter dependence.
2005  if (bProfile <= 0 || bProfile > 4) {
2006  probInt = 0.5 * M_PI * (1. - exp(-kNow));
2007  probOverlapInt = probInt / M_PI;
2008  bProbInt = probInt;
2009 
2010  // Ratio of b-integrated k * overlap / (1 - exp( - k * overlap)).
2011  nNow = M_PI * kNow * overlapInt / probInt;
2012 
2013  // Else integrate overlap over impact parameter.
2014  } else if (bProfile < 4) {
2015 
2016  // Reset integrals.
2017  overlapInt = (bProfile == 3) ? 0. : 0.5;
2018  probInt = 0.;
2019  probOverlapInt = 0.;
2020  bProbInt = 0.;
2021  pastBDiv = false;
2022  overlapHighB = 0.;
2023 
2024  // Step in b space.
2025  double b = -0.5 * deltaB;
2026  double bArea = 0.;
2027  do {
2028  b += deltaB;
2029  bArea = 2. * M_PI * b * deltaB;
2030 
2031  // Evaluate overlap at current b value.
2032  if (bProfile == 1) {
2033  overlapNow = normPi * exp( -b*b);
2034  } else if (bProfile == 2) {
2035  overlapNow = normPi * ( fracA * exp( -min(EXPMAX, b*b))
2036  + fracB * exp( -min(EXPMAX, b*b / radius2B)) / radius2B
2037  + fracC * exp( -min(EXPMAX, b*b / radius2C)) / radius2C );
2038  } else {
2039  overlapNow = normPi * exp( -pow( b, expPow));
2040  overlapInt += bArea * overlapNow;
2041  }
2042  if (pastBDiv) overlapHighB += bArea * overlapNow;
2043 
2044  // Calculate interaction probability and integrate.
2045  probNow = 1. - exp( -min(EXPMAX, M_PI * kNow * overlapNow));
2046  probInt += bArea * probNow;
2047  probOverlapInt += bArea * overlapNow * probNow;
2048  bProbInt += b * bArea * probNow;
2049 
2050  // Check when interaction probability has dropped sufficiently.
2051  if (!pastBDiv && probNow < PROBATLOWB) {
2052  bDiv = b + 0.5 * deltaB;
2053  pastBDiv = true;
2054  }
2055 
2056  // Continue out in b until overlap too small.
2057  } while (b < 1. || b * probNow > BMAX);
2058 
2059  // Ratio of b-integrated k * overlap / (1 - exp( - k * overlap)).
2060  nNow = M_PI * kNow * overlapInt / probInt;
2061 
2062  // x-dependent matter profile.
2063  } else if (bProfile == 4) {
2064  rescale2 = pow2(kNow / XDEP_A0);
2065  probInt = 0.;
2066  double b = 0.5 * bstepNow;
2067  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
2068  double bArea = 2. * M_PI * b * bstepNow;
2069  double pIntNow = 1 - exp( -min(EXPMAX, sigmaIntWgt[bBin] / rescale2) );
2070  probInt += bArea * rescale2 * pIntNow;
2071  b += bstepNow;
2072  }
2073  nNow = probInt;
2074  }
2075 
2076  // Replace lower or upper limit of k.
2077  if (nNow < nAvg) {
2078  kLow = kNow;
2079  nLow = nNow;
2080  if (stepDir == -1) stepDir = 0;
2081  } else {
2082  kHigh = kNow;
2083  nHigh = nNow;
2084  if (stepDir == 1) stepDir = -1;
2085  }
2086 
2087  // Continue iteration until convergence.
2088  } while (abs(nNow - nAvg) > KCONVERGE * nAvg);
2089 
2090  // Save relevant final numbers for overlap values.
2091  if (bProfile >= 0 && bProfile < 4) {
2092  double avgOverlap = probOverlapInt / probInt;
2093  zeroIntCorr = probOverlapInt / overlapInt;
2094  normOverlap = normPi * zeroIntCorr / avgOverlap;
2095  bAvg = bProbInt / probInt;
2096 
2097  // Values for x-dependent matter profile.
2098  } else if (bProfile == 4) {
2099  // bAvg = Int(b * Pint(b), d2b) / sigmaND.
2100  // zeroIntCorr = Int(<n(b)> * Pint(b), d2b) / sigmaInt.
2101  bAvg = 0.;
2102  zeroIntCorr = 0.;
2103  double b = 0.5 * bstepNow;
2104  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
2105  double bArea = 2. * M_PI * b * bstepNow;
2106  double pIntNow = 1. - exp( -min(EXPMAX, sigmaIntWgt[bBin] / rescale2) );
2107  bAvg += sqrt(rescale2) * b * bArea * rescale2 * pIntNow;
2108  zeroIntCorr += bArea * sigmaIntWgt[bBin] * pIntNow;
2109  b += bstepNow;
2110  }
2111  bAvg /= nNow;
2112  zeroIntCorr /= sigmaInt;
2113 
2114  // Other required values.
2115  a0now = kNow;
2116  infoPtr->seta0MPI(a0now * XDEP_SMB2FM);
2117  a02now = a0now * a0now;
2118  double xMin = 2. * pTmin / infoPtr->eCM();
2119  a2max = a0now * (XDEP_A1 + a1 * log(1. / xMin));
2120  a2max *= a2max;
2121  }
2122 
2123  // Relative rates for preselection of low-b and high-b region.
2124  // Other useful combinations for subsequent selection.
2125  if (bProfile > 0 && bProfile <= 3) {
2126  probLowB = M_PI * bDiv*bDiv;
2127  double probHighB = M_PI * kNow * overlapHighB;
2128  if (bProfile == 1) probHighB = M_PI * kNow * 0.5 * exp( -bDiv*bDiv);
2129  else if (bProfile == 2) {
2130  fracAhigh = fracA * exp( -bDiv*bDiv);
2131  fracBhigh = fracB * exp( -bDiv*bDiv / radius2B);
2132  fracChigh = fracC * exp( -bDiv*bDiv / radius2C);
2133  fracABChigh = fracAhigh + fracBhigh + fracChigh;
2134  probHighB = M_PI * kNow * 0.5 * fracABChigh;
2135  } else {
2136  cDiv = pow( bDiv, expPow);
2137  cMax = max(2. * expRev, cDiv);
2138  }
2139  probLowB /= (probLowB + probHighB);
2140  }
2141 
2142 }
2143 
2144 //--------------------------------------------------------------------------
2145 
2146 // Pick impact parameter and interaction rate enhancement beforehand,
2147 // i.e. before even the hardest interaction for minimum-bias events.
2148 
2149 void MultipartonInteractions::overlapFirst() {
2150 
2151  // Trivial values if no impact parameter dependence.
2152  if (bProfile <= 0 || bProfile > 4) {
2153  bNow = bAvg;
2154  enhanceB = zeroIntCorr;
2155  bIsSet = true;
2156  isAtLowB = true;
2157  return;
2158  }
2159 
2160  // Preliminary choice between and inside low-b and high-b regions.
2161  double overlapNow = 0.;
2162  double probAccept = 0.;
2163  do {
2164 
2165  // Treatment in low-b region: pick b flat in area.
2166  if (rndmPtr->flat() < probLowB) {
2167  isAtLowB = true;
2168  bNow = bDiv * sqrt(rndmPtr->flat());
2169 
2170  // Evaluate overlap and from that acceptance probability.
2171  if (bProfile == 1) overlapNow = normPi * exp( -bNow*bNow);
2172  else if (bProfile == 2) overlapNow = normPi *
2173  ( fracA * exp( -bNow*bNow)
2174  + fracB * exp( -bNow*bNow / radius2B) / radius2B
2175  + fracC * exp( -bNow*bNow / radius2C) / radius2C );
2176  else overlapNow = normPi * exp( -pow( bNow, expPow));
2177  probAccept = 1. - exp( -min(EXPMAX, M_PI * kNow * overlapNow));
2178 
2179  // Treatment in high-b region: pick b according to overlap.
2180  } else {
2181  isAtLowB = false;
2182 
2183  // For simple and double Gaussian pick b according to exp(-b^2 / r^2).
2184  if (bProfile == 1) {
2185  bNow = sqrt(bDiv*bDiv - log(rndmPtr->flat()));
2186  overlapNow = normPi * exp( -min(EXPMAX, bNow*bNow));
2187  } else if (bProfile == 2) {
2188  double pickFrac = rndmPtr->flat() * fracABChigh;
2189  if (pickFrac < fracAhigh)
2190  bNow = sqrt(bDiv*bDiv - log(rndmPtr->flat()));
2191  else if (pickFrac < fracAhigh + fracBhigh)
2192  bNow = sqrt(bDiv*bDiv - radius2B * log(rndmPtr->flat()));
2193  else bNow = sqrt(bDiv*bDiv - radius2C * log(rndmPtr->flat()));
2194  overlapNow = normPi * ( fracA * exp( -min(EXPMAX, bNow*bNow))
2195  + fracB * exp( -min(EXPMAX, bNow*bNow / radius2B)) / radius2B
2196  + fracC * exp( -min(EXPMAX, bNow*bNow / radius2C)) / radius2C );
2197 
2198  // For exp( - b^expPow) transform to variable c = b^expPow so that
2199  // f(b) = b * exp( - b^expPow) -> f(c) = c^r * exp(-c) with r = expRev.
2200  // case hasLowPow: expPow < 2 <=> r > 0: preselect according to
2201  // f(c) < N exp(-c/2) and then accept with N' * c^r * exp(-c/2).
2202  } else if (hasLowPow) {
2203  double cNow, acceptC;
2204  do {
2205  cNow = cDiv - 2. * log(rndmPtr->flat());
2206  acceptC = pow(cNow / cMax, expRev) * exp( -0.5 * (cNow - cMax));
2207  } while (acceptC < rndmPtr->flat());
2208  bNow = pow( cNow, 1. / expPow);
2209  overlapNow = normPi * exp( -cNow);
2210 
2211  // case !hasLowPow: expPow >= 2 <=> - 1 < r < 0: preselect according to
2212  // f(c) < N exp(-c) and then accept with N' * c^r.
2213  } else {
2214  double cNow, acceptC;
2215  do {
2216  cNow = cDiv - log(rndmPtr->flat());
2217  acceptC = pow(cNow / cDiv, expRev);
2218  } while (acceptC < rndmPtr->flat());
2219  bNow = pow( cNow, 1. / expPow);
2220  overlapNow = normPi * exp( -cNow);
2221  }
2222  double temp = M_PI * kNow *overlapNow;
2223  probAccept = (1. - exp( -min(EXPMAX, temp))) / temp;
2224  }
2225 
2226  // Confirm choice of b value. Derive enhancement factor.
2227  } while (probAccept < rndmPtr->flat());
2228 
2229  // Same enhancement for hardest process and all subsequent MPI
2230  enhanceB = enhanceBmax = enhanceBnow = (normOverlap / normPi) * overlapNow;
2231 
2232  // Done.
2233  bIsSet = true;
2234 
2235 }
2236 
2237 //--------------------------------------------------------------------------
2238 
2239 // Pick impact parameter and interaction rate enhancement afterwards,
2240 // i.e. after a hard interaction is known but before rest of MPI treatment.
2241 
2242 void MultipartonInteractions::overlapNext(Event& event, double pTscale) {
2243 
2244  // Default, valid for bProfile = 0. Also initial Sudakov.
2245  enhanceB = zeroIntCorr;
2246  if (bProfile <= 0 || bProfile > 4) return;
2247 
2248  // Alternative choices of event scale for Sudakov in (pT, b) space.
2249  if (bSelScale == 1) {
2250  vector<double> mmT;
2251  for (int i = 5; i < event.size(); ++i) if (event[i].isFinal()) {
2252  mmT.push_back( event[i].m() + event[i].mT() );
2253  for (int j = int(mmT.size()) - 1; j > 0; --j)
2254  if (mmT[j] > mmT[j - 1]) swap( mmT[j], mmT[j - 1] );
2255  }
2256  pTscale = 0.5 * mmT[0];
2257  for (int j = 1; j < int(mmT.size()); ++j) pTscale += mmT[j] / (j + 1.);
2258  } else if (bSelScale == 2) pTscale = event.scale();
2259  double pT2scale = pTscale*pTscale;
2260 
2261  // Use trial interaction for x-dependent matter profile.
2262  if (bProfile == 4) {
2263  double pTtrial = 0.;
2264  do {
2265  // Pick b according to wanted O(b, x1, x2).
2266  double expb2 = rndmPtr->flat();
2267  double w1 = XDEP_A1 + a1 * log(1. / infoPtr->x1());
2268  double w2 = XDEP_A1 + a1 * log(1. / infoPtr->x2());
2269  double fac = a02now * (w1 * w1 + w2 * w2);
2270  b2now = - fac * log(expb2);
2271  bNow = sqrt(b2now);
2272 
2273  // Enhancement factor for the hard process and overestimate
2274  // for fastPT2. Note that existing framework has a (1. / sigmaND)
2275  // present.
2276  enhanceB = sigmaND / M_PI / fac * expb2;
2277  enhanceBmax = sigmaND / 2. / M_PI / a02now
2278  * exp( -b2now / 2. / a2max );
2279 
2280  // Trial interaction. Keep going until pTtrial < pTscale.
2281  pTtrial = pTnext(pTmax, pTmin, event);
2282  } while (pTtrial > pTscale);
2283  bIsSet = true;
2284  return;
2285  }
2286 
2287  // Begin loop over pT-dependent rejection of b value.
2288  do {
2289 
2290  // Flat enhancement distribution for simple Gaussian.
2291  if (bProfile == 1) {
2292  double expb2 = rndmPtr->flat();
2293  // Same enhancement for hardest process and all subsequent MPI.
2294  enhanceB = enhanceBmax = enhanceBnow = normOverlap * expb2;
2295  bNow = sqrt( -log(expb2));
2296 
2297  // For double Gaussian go the way via b, according to each Gaussian.
2298  } else if (bProfile == 2) {
2299  double bType = rndmPtr->flat();
2300  double b2 = -log( rndmPtr->flat() );
2301  if (bType < fracA) ;
2302  else if (bType < fracA + fracB) b2 *= radius2B;
2303  else b2 *= radius2C;
2304  // Same enhancement for hardest process and all subsequent MPI.
2305  enhanceB = enhanceBmax = enhanceBnow = normOverlap *
2306  ( fracA * exp( -min(EXPMAX, b2))
2307  + fracB * exp( -min(EXPMAX, b2 / radius2B)) / radius2B
2308  + fracC * exp( -min(EXPMAX, b2 / radius2C)) / radius2C );
2309  bNow = sqrt(b2);
2310 
2311  // For exp( - b^expPow) transform to variable c = b^expPow so that
2312  // f(b) = b * exp( - b^expPow) -> f(c) = c^r * exp(-c) with r = expRev.
2313  // case hasLowPow: expPow < 2 <=> r > 0:
2314  // f(c) < r^r exp(-r) for c < 2r, < (2r)^r exp(-r-c/2) for c > 2r.
2315  } else if (bProfile == 3 && hasLowPow) {
2316  double cNow, acceptC;
2317  double probLowC = expRev / (expRev + pow(2., expRev) * exp( - expRev));
2318  do {
2319  if (rndmPtr->flat() < probLowC) {
2320  cNow = 2. * expRev * rndmPtr->flat();
2321  acceptC = pow( cNow / expRev, expRev) * exp(expRev - cNow);
2322  } else {
2323  cNow = 2. * (expRev - log( rndmPtr->flat() ));
2324  acceptC = pow(0.5 * cNow / expRev, expRev)
2325  * exp(expRev - 0.5 * cNow);
2326  }
2327  } while (acceptC < rndmPtr->flat());
2328  // Same enhancement for hardest process and all subsequent MPI.
2329  enhanceB = enhanceBmax = enhanceBnow = normOverlap *exp(-cNow);
2330  bNow = pow( cNow, 1. / expPow);
2331 
2332  // case !hasLowPow: expPow >= 2 <=> - 1 < r < 0:
2333  // f(c) < c^r for c < 1, f(c) < exp(-c) for c > 1.
2334  } else if (bProfile == 3 && !hasLowPow) {
2335  double cNow, acceptC;
2336  double probLowC = expPow / (2. * exp(-1.) + expPow);
2337  do {
2338  if (rndmPtr->flat() < probLowC) {
2339  cNow = pow( rndmPtr->flat(), 0.5 * expPow);
2340  acceptC = exp(-cNow);
2341  } else {
2342  cNow = 1. - log( rndmPtr->flat() );
2343  acceptC = pow( cNow, expRev);
2344  }
2345  } while (acceptC < rndmPtr->flat());
2346  // Same enhancement for hardest process and all subsequent MPI.
2347  enhanceB = enhanceBmax = enhanceBnow = normOverlap * exp(-cNow);
2348  bNow = pow( cNow, 1. / expPow);
2349  }
2350 
2351  // Evaluate "Sudakov form factor" for not having a harder interaction.
2352  } while (sudakov(pT2scale, enhanceB) < rndmPtr->flat());
2353 
2354  // Done.
2355  bIsSet = true;
2356 }
2357 
2358 //--------------------------------------------------------------------------
2359 
2360 // Printe statistics on number of multiparton-interactions processes.
2361 
2362 void MultipartonInteractions::statistics(bool resetStat, ostream& os) {
2363 
2364  // Header.
2365  os << "\n *------- PYTHIA Multiparton Interactions Statistics -----"
2366  << "---*\n"
2367  << " | "
2368  << " |\n"
2369  << " | Note: excludes hardest subprocess if already listed above "
2370  << " |\n"
2371  << " | "
2372  << " |\n"
2373  << " | Subprocess Code | Times"
2374  << " |\n"
2375  << " | | "
2376  << " |\n"
2377  << " |------------------------------------------------------------"
2378  << "-|\n"
2379  << " | | "
2380  << " |\n";
2381 
2382  // Loop over existing processes. Sum of all subprocesses.
2383  int numberSum = 0;
2384  for ( map<int, int>::iterator iter = nGen.begin(); iter != nGen.end();
2385  ++iter) {
2386  int code = iter -> first;
2387  int number = iter->second;
2388  numberSum += number;
2389 
2390  // Find process name that matches code.
2391  string name = " ";
2392  bool foundName = false;
2393  SigmaMultiparton* dSigma;
2394  for (int i = 0; i < 4; ++i) {
2395  if (i == 0) dSigma = &sigma2gg;
2396  else if (i == 1) dSigma = &sigma2qg;
2397  else if (i == 2) dSigma = &sigma2qqbarSame;
2398  else dSigma = &sigma2qq;
2399  int nProc = dSigma->nProc();
2400  for (int iProc = 0; iProc < nProc; ++iProc)
2401  if (dSigma->codeProc(iProc) == code) {
2402  name = dSigma->nameProc(iProc);
2403  foundName = true;
2404  }
2405  if (foundName) break;
2406  }
2407 
2408  // Print individual process info.
2409  os << " | " << left << setw(40) << name << right << setw(5) << code
2410  << " | " << setw(11) << number << " |\n";
2411  }
2412 
2413  // Print summed process info.
2414  os << " | "
2415  << " |\n"
2416  << " | " << left << setw(45) << "sum" << right << " | " << setw(11)
2417  << numberSum << " |\n";
2418 
2419  // Listing finished.
2420  os << " | | "
2421  << " |\n"
2422  << " *------- End PYTHIA Multiparton Interactions Statistics ----"
2423  << "-*" << endl;
2424 
2425  // Optionally reset statistics contents.
2426  if (resetStat) resetStatistics();
2427 
2428 }
2429 
2430 //==========================================================================
2431 
2432 } // end namespace Pythia8
Definition: AgUStep.h:26