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) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // 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 // Limit below which scientific notation is used for printing.
337 const double MultipartonInteractions::SIGMAMBLIMIT = 1.;
338 
339 //--------------------------------------------------------------------------
340 
341 // Initialize the generation process for given beams.
342 
343 bool MultipartonInteractions::init( bool doMPIinit, int iDiffSysIn,
344  Info* infoPtrIn, Settings& settings, ParticleData* particleDataPtr,
345  Rndm* rndmPtrIn, BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
346  Couplings* couplingsPtrIn, PartonSystems* partonSystemsPtrIn,
347  SigmaTotal* sigmaTotPtrIn, UserHooks* userHooksPtrIn,
348  PartonVertex* partonVertexPtrIn, bool hasGammaIn) {
349 
350  // Store input pointers for future use. Done if no initialization.
351  iDiffSys = iDiffSysIn;
352  infoPtr = infoPtrIn;
353  rndmPtr = rndmPtrIn;
354  beamAPtr = beamAPtrIn;
355  beamBPtr = beamBPtrIn;
356  couplingsPtr = couplingsPtrIn;
357  partonSystemsPtr = partonSystemsPtrIn;
358  sigmaTotPtr = sigmaTotPtrIn;
359  userHooksPtr = userHooksPtrIn;
360  partonVertexPtr = partonVertexPtrIn;
361  hasGamma = hasGammaIn;
362  if (!doMPIinit) return false;
363 
364  // If both beams are baryons then softer PDF's than for mesons/Pomerons.
365  hasBaryonBeams = ( beamAPtr->isBaryon() && beamBPtr->isBaryon() );
366  hasPomeronBeams = ( beamAPtr->id() == 990 || beamBPtr->id() == 990 );
367 
368  // Matching in pT of hard interaction to further interactions.
369  pTmaxMatch = settings.mode("MultipartonInteractions:pTmaxMatch");
370 
371  // Parameters of alphaStrong generation.
372  alphaSvalue = settings.parm("MultipartonInteractions:alphaSvalue");
373  alphaSorder = settings.mode("MultipartonInteractions:alphaSorder");
374  alphaSnfmax = settings.mode("StandardModel:alphaSnfmax");
375 
376  // Parameters of alphaEM generation.
377  alphaEMorder = settings.mode("MultipartonInteractions:alphaEMorder");
378 
379  // Parameters of cross section generation.
380  Kfactor = settings.parm("MultipartonInteractions:Kfactor");
381 
382  // Check if photon-photon or photon-hadron collision.
383  isGammaGamma = beamAPtr->isGamma() && beamBPtr->isGamma();
384  isGammaHadron = beamAPtr->isGamma() && beamBPtr->isHadron();
385  isHadronGamma = beamAPtr->isHadron() && beamBPtr->isGamma();
386 
387  // Regularization of QCD evolution for pT -> 0.
388  // Separate default parameters for photon-photon collisions.
389  if (isGammaGamma) {
390  pT0paramMode = settings.mode("PhotonPhoton:pT0parametrization");
391  pT0Ref = settings.parm("PhotonPhoton:pT0Ref");
392  ecmRef = settings.parm("PhotonPhoton:ecmRef");
393  ecmPow = settings.parm("PhotonPhoton:ecmPow");
394  pTmin = settings.parm("PhotonPhoton:pTmin");
395  } else {
396  pT0paramMode = settings.mode("MultipartonInteractions:pT0parametrization");
397  pT0Ref = settings.parm("MultipartonInteractions:pT0Ref");
398  ecmRef = settings.parm("MultipartonInteractions:ecmRef");
399  ecmPow = settings.parm("MultipartonInteractions:ecmPow");
400  pTmin = settings.parm("MultipartonInteractions:pTmin");
401  }
402 
403  // Impact parameter profile: nondiffractive topologies.
404  if (iDiffSys == 0) {
405  bProfile = settings.mode("MultipartonInteractions:bProfile");
406  coreRadius = settings.parm("MultipartonInteractions:coreRadius");
407  coreFraction = settings.parm("MultipartonInteractions:coreFraction");
408  expPow = settings.parm("MultipartonInteractions:expPow");
409  expPow = max(EXPPOWMIN, expPow);
410  // x-dependent impact parameter profile.
411  a1 = settings.parm("MultipartonInteractions:a1");
412 
413  // Impact parameter profile: diffractive topologies.
414  } else {
415  bProfile = settings.mode("Diffraction:bProfile");
416  coreRadius = settings.parm("Diffraction:coreRadius");
417  coreFraction = settings.parm("Diffraction:coreFraction");
418  expPow = settings.parm("Diffraction:expPow");
419  expPow = max(EXPPOWMIN, expPow);
420  }
421 
422  // No x-dependent impact-parameter profile for diffraction.
423  if ((iDiffSys > 0 || settings.flag("Diffraction:doHard")) && bProfile == 4) {
424  infoPtr->errorMsg("Error in MultipartonInteractions::init:"
425  " chosen b profile not allowed for diffraction");
426  return false;
427  }
428 
429  // Common choice of "pT" scale for determining impact parameter.
430  bSelScale = settings.mode("MultipartonInteractions:bSelScale");
431 
432  // Process sets to include in machinery.
433  processLevel = settings.mode("MultipartonInteractions:processLevel");
434 
435  // Parameters of rescattering description.
436  allowRescatter = settings.flag("MultipartonInteractions:allowRescatter");
437  allowDoubleRes
438  = settings.flag("MultipartonInteractions:allowDoubleRescatter");
439  rescatterMode = settings.mode("MultipartonInteractions:rescatterMode");
440  ySepResc = settings.parm("MultipartonInteractions:ySepRescatter");
441  deltaYResc = settings.parm("MultipartonInteractions:deltaYRescatter");
442 
443  // Rescattering not yet implemented for x-dependent impact profile.
444  if (bProfile == 4) allowRescatter = false;
445 
446  // A global recoil FSR stategy restricts rescattering.
447  globalRecoilFSR = settings.flag("TimeShower:globalRecoil");
448  nMaxGlobalRecoilFSR = settings.mode("TimeShower:nMaxGlobalRecoil");
449 
450  // Various other parameters.
451  nQuarkIn = settings.mode("MultipartonInteractions:nQuarkIn");
452  nSample = settings.mode("MultipartonInteractions:nSample");
453 
454  // Optional dampening at small pT's when large multiplicities.
455  enhanceScreening = settings.mode("MultipartonInteractions:enhanceScreening");
456 
457  // Parameters for diffractive systems.
458  sigmaPomP = settings.parm("Diffraction:sigmaRefPomP");
459  mPomP = settings.parm("Diffraction:mRefPomP");
460  pPomP = settings.parm("Diffraction:mPowPomP");
461  mMinPertDiff = settings.parm("Diffraction:mMinPert");
462  bSelHard = settings.mode("Diffraction:bSelHard");
463 
464  // Beam particles might not be found from the usual positions.
465  beamOffset = 0;
466 
467  // Possibility to allow user veto of MPI.
468  canVetoMPI = (userHooksPtr != 0) ? userHooksPtr->canVetoMPIEmission()
469  : false;
470 
471  // Possibility to set parton vertex information.
472  doPartonVertex = settings.flag("PartonVertex:setVertex")
473  && (partonVertexPtr != 0);
474 
475  // Some common combinations for double Gaussian, as shorthand.
476  if (bProfile == 2) {
477  fracA = pow2(1. - coreFraction);
478  fracB = 2. * coreFraction * (1. - coreFraction);
479  fracC = pow2(coreFraction);
480  radius2B = 0.5 * (1. + pow2(coreRadius));
481  radius2C = pow2(coreRadius);
482 
483  // Some common combinations for exp(b^pow), as shorthand.
484  } else if (bProfile == 3) {
485  hasLowPow = (expPow < 2.);
486  expRev = 2. / expPow - 1.;
487  }
488  enhanceBavg = 1.;
489 
490  // Initialize alpha_strong generation.
491  alphaS.init( alphaSvalue, alphaSorder, alphaSnfmax, false);
492  double Lambda3 = alphaS.Lambda3();
493 
494  // Initialize alphaEM generation.
495  alphaEM.init( alphaEMorder, &settings);
496 
497  // Attach matrix-element calculation objects.
498  sigma2gg.init( 0, processLevel, infoPtr, &settings, particleDataPtr,
499  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
500  sigma2qg.init( 1, processLevel, infoPtr, &settings, particleDataPtr,
501  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
502  sigma2qqbarSame.init( 2, processLevel, infoPtr, &settings, particleDataPtr,
503  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
504  sigma2qq.init( 3, processLevel, infoPtr, &settings, particleDataPtr,
505  rndmPtr, beamAPtr, beamBPtr, couplingsPtr);
506 
507  // Calculate invariant mass of system.
508  eCM = infoPtr->eCM();
509  sCM = eCM * eCM;
510  mMaxPertDiff = eCM;
511  eCMsave = eCM;
512 
513  // Limits on invariant mass of gm+gm system.
514  mGmGmMin = settings.parm("Photon:Wmin");
515  mGmGmMax = settings.parm("Photon:Wmax");
516  if ( mGmGmMax < mGmGmMin ) mGmGmMax = eCM;
517 
518  // Get the total inelastic and nondiffractive cross section.
519  // Ensure correct cross sections for VMD photons.
520  if (infoPtr->isVMDstateA() || infoPtr->isVMDstateB())
521  sigmaTotPtr->calc(infoPtr->idA(), infoPtr->idB(), infoPtr->eCM());
522  // Ensure correct cross sections also for non-VMD photon beams.
523  else if ( (isGammaGamma || isGammaHadron || isHadronGamma) && !hasGamma)
524  sigmaTotPtr->calc(infoPtr->idA(), infoPtr->idB(), infoPtr->eCM());
525  if (!sigmaTotPtr->hasSigmaTot()) return false;
526  bool isNonDiff = (iDiffSys == 0);
527  sigmaND = sigmaTotPtr->sigmaND();
528  double sigmaMaxViol = 0.;
529 
530  // Output initialization info - first part.
531  bool showMPI = settings.flag("Init:showMultipartonInteractions");
532  if (showMPI) {
533  cout << "\n *------- PYTHIA Multiparton Interactions Initialization "
534  << "---------* \n"
535  << " | "
536  << " | \n";
537  if (isNonDiff && !hasGamma)
538  cout << " | sigmaNonDiffractive = " << setprecision(2)
539  << ((sigmaND > 1.) ? fixed : scientific) << setw(8) << sigmaND
540  << " mb | \n";
541  else if (iDiffSys == 1)
542  cout << " | diffraction XB "
543  << " | \n";
544  else if (iDiffSys == 2)
545  cout << " | diffraction AX "
546  << " | \n";
547  else if (iDiffSys == 3)
548  cout << " | diffraction AXB "
549  << " | \n";
550  else if ( hasGamma && isGammaGamma )
551  cout << " | l+l- -> gamma+gamma -> X "
552  << " | \n";
553  else if ( hasGamma && isGammaHadron )
554  cout << " | lepton+hadron -> gamma+hadron -> X "
555  << " | \n";
556  else if ( hasGamma && isHadronGamma )
557  cout << " | hadron+lepton -> hadron+gamma -> X "
558  << " | \n";
559  cout << " | "
560  << " | \n";
561  }
562 
563  // For diffraction need to cover range of diffractive masses.
564  nStep = (iDiffSys == 0) ? 1 : 5;
565  eStepSize = (nStep < 2) ? 1.
566  : log(mMaxPertDiff / mMinPertDiff) / (nStep - 1.);
567 
568  // For photons from lepton cover range of gm+gm invariant masses.
569  if (hasGamma) {
570  nStep = 5;
571  eStepSize = log(mGmGmMax / mGmGmMin) / (nStep - 1.);
572  }
573 
574  for (int iStep = 0; iStep < nStep; ++iStep) {
575 
576  // Update and output current diffractive mass and
577  // fictitious Pomeron-proton cross section for normalization.
578  if (nStep > 1) {
579  if (!hasGamma) eCM = mMinPertDiff * pow( mMaxPertDiff / mMinPertDiff,
580  iStep / (nStep - 1.) );
581  else eCM = mGmGmMin * pow( mGmGmMax / mGmGmMin, iStep / (nStep - 1.) );
582  sCM = eCM * eCM;
583 
584 
585  // MPI for diffractive events. Rescale Pom/p flux to use for Pom/gamma.
586  if (hasPomeronBeams) {
587  double gamPomRatio = 1.;
588  if (hasGamma) {
589  sigmaTotPtr->calc(22, 2212, eCM);
590  double sigGamP = sigmaTotPtr->sigmaTot();
591  sigmaTotPtr->calc(2212, 2212, eCM);
592  double sigPP = sigmaTotPtr->sigmaTot();
593  gamPomRatio = sigGamP / sigPP;
594  }
595  sigmaND = gamPomRatio * sigmaPomP * pow( eCM / mPomP, pPomP);
596  if (showMPI) cout << " | diffractive mass = " << scientific
597  << setprecision(2) << setw(8) << eCM << " GeV and sigmaNorm = "
598  << ((sigmaND > SIGMAMBLIMIT) ? fixed : scientific)
599  << setw(8) << sigmaND << " mb | \n";
600 
601  // Keep track of pomeron momentum fraction.
602  if ( beamAPtr->id() == 990 && beamBPtr->id() == 990 ) {
603  beamAPtr->xPom(eCM/eCMsave);
604  beamBPtr->xPom(eCM/eCMsave);
605  }
606  else if ( beamAPtr->id() == 990 )
607  beamAPtr->xPom(pow2(eCM/eCMsave));
608  else if ( beamBPtr->id() == 990 )
609  beamBPtr->xPom(pow2(eCM/eCMsave));
610  // MPI with photons from leptons.
611  } else {
612 
613  // Hadron-photon case.
614  if ( isHadronGamma ) {
615  sigmaTotPtr->calc( beamAPtr->id(), 22, eCM );
616  sigmaND = sigmaTotPtr->sigmaND();
617  if (showMPI) cout << " | hadron+gamma eCM = " << scientific
618  << setprecision(2) << setw(8) << eCM << " GeV and sigmaNorm = "
619  << ((sigmaND > SIGMAMBLIMIT) ? fixed : scientific)
620  << setw(8) << sigmaND << " mb | \n";
621 
622  // Photon-hadron case.
623  } else if ( isGammaHadron ) {
624  sigmaTotPtr->calc( 22, beamBPtr->id(), eCM );
625  sigmaND = sigmaTotPtr->sigmaND();
626  if (showMPI) cout << " | gamma+hadron eCM = " << scientific
627  << setprecision(2) << setw(8) << eCM << " GeV and sigmaNorm = "
628  << ((sigmaND > SIGMAMBLIMIT) ? fixed : scientific)
629  << setw(8) << sigmaND << " mb | \n";
630 
631  // Photon-photon case.
632  } else {
633  sigmaTotPtr->calc( 22, 22, eCM );
634  sigmaND = sigmaTotPtr->sigmaND();
635  if (showMPI) cout << " | gamma+gamma eCM = " << scientific
636  << setprecision(2) << setw(8) << eCM << " GeV and sigmaNorm = "
637  << ((sigmaND > SIGMAMBLIMIT) ? fixed : scientific)
638  << setw(8) << sigmaND << " mb | \n";
639  }
640  }
641 
642  }
643 
644  // Set current pT0 scale according to chosed parametrization.
645  if (pT0paramMode == 0) pT0 = pT0Ref * pow(eCM / ecmRef, ecmPow);
646  else pT0 = pT0Ref + ecmPow * log (eCM / ecmRef);
647 
648  // The pT0 value may need to be decreased, if sigmaInt < sigmaND.
649  double pT4dSigmaMaxBeg = 0.;
650  for ( ; ; ) {
651 
652  // Derived pT kinematics combinations.
653  pT20 = pT0*pT0;
654  pT2min = pTmin*pTmin;
655  pTmax = 0.5*eCM;
656  pT2max = pTmax*pTmax;
657  pT20R = RPT20 * pT20;
658  pT20minR = pT2min + pT20R;
659  pT20maxR = pT2max + pT20R;
660  pT20min0maxR = pT20minR * pT20maxR;
661  pT2maxmin = pT2max - pT2min;
662 
663  // Provide upper estimate of interaction rate d(Prob)/d(pT2).
664  upperEnvelope();
665 
666  // Setup binning in b for x-dependent matter profile.
667  if (bProfile == 4) {
668  sigmaIntWgt.resize(XDEP_BBIN);
669  sigmaSumWgt.resize(XDEP_BBIN);
670  bstepNow = XDEP_BSTEP;
671  }
672 
673  // Integrate the parton-parton interaction cross section.
674  pT4dSigmaMaxBeg = pT4dSigmaMax;
675  jetCrossSection();
676 
677  // If the overlap-weighted cross section has not fallen below
678  // cutoff, then increase bin size in b and reintegrate.
679  while (bProfile == 4
680  && sigmaIntWgt[XDEP_BBIN - 1] > XDEP_CUTOFF * sigmaInt) {
681  bstepNow += XDEP_BSTEPINC;
682  jetCrossSection();
683  }
684 
685  // Sufficiently big SigmaInt or reduce pT0; maybe also pTmin.
686  if (sigmaInt > SIGMASTEP * sigmaND) break;
687  if (showMPI) cout << fixed << setprecision(2) << " | pT0 = "
688  << setw(5) << pT0 << " gives sigmaInteraction = " << setw(8)
689  << ((sigmaInt > SIGMAMBLIMIT) ? fixed : scientific) << sigmaInt
690  << " mb: rejected | \n";
691  if (pTmin > pT0) pTmin *= PT0STEP;
692  pT0 *= PT0STEP;
693 
694  // Give up if pT0 and pTmin fall too low.
695  if ( max(pT0, pTmin) < max(PT0MIN, Lambda3) ) {
696  infoPtr->errorMsg("Error in MultipartonInteractions::init:"
697  " failed to find acceptable pT0 and pTmin");
698  infoPtr->setTooLowPTmin(true);
699  return false;
700  }
701  }
702 
703  // Output for accepted pT0.
704  if (showMPI) cout << fixed << setprecision(2) << " | pT0 = "
705  << setw(5) << pT0 << " gives sigmaInteraction = "<< setw(8)
706  << ((sigmaInt > SIGMAMBLIMIT) ? fixed : scientific) << sigmaInt
707  << " mb: accepted | \n";
708 
709  // Calculate factor relating matter overlap and interaction rate.
710  overlapInit();
711 
712  // Maximum violation relative to first estimate.
713  sigmaMaxViol = max( sigmaMaxViol, pT4dSigmaMax / pT4dSigmaMaxBeg);
714 
715  // Save values calculated.
716  if (nStep > 1) {
717  pT0Save[iStep] = pT0;
718  pT4dSigmaMaxSave[iStep] = pT4dSigmaMax;
719  pT4dProbMaxSave[iStep] = pT4dProbMax;
720  sigmaIntSave[iStep] = sigmaInt;
721  for (int j = 0; j <= 100; ++j) sudExpPTSave[iStep][j] = sudExpPT[j];
722  zeroIntCorrSave[iStep] = zeroIntCorr;
723  normOverlapSave[iStep] = normOverlap;
724  kNowSave[iStep] = kNow;
725  bAvgSave[iStep] = bAvg;
726  bDivSave[iStep] = bDiv;
727  probLowBSave[iStep] = probLowB;
728  fracAhighSave[iStep] = fracAhigh;
729  fracBhighSave[iStep] = fracBhigh;
730  fracChighSave[iStep] = fracBhigh;
731  fracABChighSave[iStep] = fracABChigh;
732  cDivSave[iStep] = cDiv;
733  cMaxSave[iStep] = cMax;
734  }
735 
736  // End of loop over diffractive/invariant gamma+gamma masses.
737  }
738 
739  // Reset pomeron momentum fraction.
740  beamAPtr->xPom();
741  beamBPtr->xPom();
742 
743  // Output details for x-dependent matter profile.
744  if (bProfile == 4 && showMPI)
745  cout << " | "
746  << " | \n"
747  << fixed << setprecision(2)
748  << " | x-dependent matter profile: a1 = " << a1 << ", "
749  << "a0 = " << a0now * XDEP_SMB2FM << ", bStep = "
750  << bstepNow << " | \n";
751 
752  // End initialization printout.
753  if (showMPI) cout << " | "
754  << " | \n"
755  << " *------- End PYTHIA Multiparton Interactions Initialization"
756  << " -----* " << endl;
757 
758  // Amount of violation from upperEnvelope to jetCrossSection.
759  if (sigmaMaxViol > 1.) {
760  ostringstream osWarn;
761  osWarn << "by factor " << fixed << setprecision(3) << sigmaMaxViol;
762  infoPtr->errorMsg("Warning in MultipartonInteractions::init:"
763  " maximum increased", osWarn.str());
764  }
765 
766  // Reset statistics.
767  SigmaMultiparton* dSigma;
768  for (int i = 0; i < 4; ++i) {
769  if (i == 0) dSigma = &sigma2gg;
770  else if (i == 1) dSigma = &sigma2qg;
771  else if (i == 2) dSigma = &sigma2qqbarSame;
772  else dSigma = &sigma2qq;
773  int nProc = dSigma->nProc();
774  for (int iProc = 0; iProc < nProc; ++iProc)
775  nGen[ dSigma->codeProc(iProc) ] = 0;
776  }
777 
778  // Additional setup for x-dependent matter profile.
779  if (bProfile == 4) {
780  sigmaIntWgt.clear();
781  sigmaSumWgt.clear();
782  }
783  // No preselection of sea/valence content and initialise a0.
784  vsc1 = 0;
785  vsc2 = 0;
786 
787  // Done.
788  return true;
789 }
790 
791 //--------------------------------------------------------------------------
792 
793 // Reset impact parameter choice and update the CM energy.
794 // For diffraction also interpolate parameters to current CM energy.
795 
796 void MultipartonInteractions::reset( ) {
797 
798  // Reset impact parameter choice.
799  bIsSet = false;
800  bSetInFirst = false;
801 
802  // Update CM energy. Done if not diffraction and not new energy.
803  eCM = infoPtr->eCM();
804  sCM = eCM * eCM;
805  if (nStep == 1 || abs( eCM / eCMsave - 1.) < ECMDEV) return;
806 
807  // Set fictitious Pomeron-proton cross section for diffractive system.
808  if (!hasGamma) sigmaND = sigmaPomP * pow( eCM / mPomP, pPomP);
809  // For photons from leptons calculate sigmaND at updated CM energy.
810  else {
811  sigmaTotPtr->calc( beamAPtr->id(), beamBPtr->id(), eCM );
812  sigmaND = sigmaTotPtr->sigmaND();
813  }
814 
815  // Current interpolation point.
816  eCMsave = eCM;
817  if (!hasGamma) eStepSave = log(eCM / mMinPertDiff) / eStepSize;
818  else eStepSave = log(eCM / mGmGmMin) / eStepSize;
819  iStepFrom = max( 0, min( nStep - 2, int( eStepSave) ) );
820  iStepTo = iStepFrom + 1;
821  eStepTo = max( 0., min( 1., eStepSave - iStepFrom) );
822  eStepFrom = 1. - eStepTo;
823 
824  // Update pT0 and combinations derived from it.
825  pT0 = eStepFrom * pT0Save[iStepFrom]
826  + eStepTo * pT0Save[iStepTo];
827  pT20 = pT0*pT0;
828  pT2min = pTmin*pTmin;
829  pTmax = 0.5*eCM;
830  pT2max = pTmax*pTmax;
831  pT20R = RPT20 * pT20;
832  pT20minR = pT2min + pT20R;
833  pT20maxR = pT2max + pT20R;
834  pT20min0maxR = pT20minR * pT20maxR;
835  pT2maxmin = pT2max - pT2min;
836 
837  // Update other parameters used in pT choice.
838  pT4dSigmaMax = eStepFrom * pT4dSigmaMaxSave[iStepFrom]
839  + eStepTo * pT4dSigmaMaxSave[iStepTo];
840  pT4dProbMax = eStepFrom * pT4dProbMaxSave[iStepFrom]
841  + eStepTo * pT4dProbMaxSave[iStepTo];
842  sigmaInt = eStepFrom * sigmaIntSave[iStepFrom]
843  + eStepTo * sigmaIntSave[iStepTo];
844  for (int j = 0; j <= 100; ++j)
845  sudExpPT[j] = eStepFrom * sudExpPTSave[iStepFrom][j]
846  + eStepTo * sudExpPTSave[iStepTo][j];
847 
848  // Update parameters related to the impact-parameter picture.
849  zeroIntCorr = eStepFrom * zeroIntCorrSave[iStepFrom]
850  + eStepTo * zeroIntCorrSave[iStepTo];
851  normOverlap = eStepFrom * normOverlapSave[iStepFrom]
852  + eStepTo * normOverlapSave[iStepTo];
853  kNow = eStepFrom * kNowSave[iStepFrom]
854  + eStepTo * kNowSave[iStepTo];
855  bAvg = eStepFrom * bAvgSave[iStepFrom]
856  + eStepTo * bAvgSave[iStepTo];
857  bDiv = eStepFrom * bDivSave[iStepFrom]
858  + eStepTo * bDivSave[iStepTo];
859  probLowB = eStepFrom * probLowBSave[iStepFrom]
860  + eStepTo * probLowBSave[iStepTo];
861  fracAhigh = eStepFrom * fracAhighSave[iStepFrom]
862  + eStepTo * fracAhighSave[iStepTo];
863  fracBhigh = eStepFrom * fracBhighSave[iStepFrom]
864  + eStepTo * fracBhighSave[iStepTo];
865  fracChigh = eStepFrom * fracChighSave[iStepFrom]
866  + eStepTo * fracChighSave[iStepTo];
867  fracABChigh = eStepFrom * fracABChighSave[iStepFrom]
868  + eStepTo * fracABChighSave[iStepTo];
869  cDiv = eStepFrom * cDivSave[iStepFrom]
870  + eStepTo * cDivSave[iStepTo];
871  cMax = eStepFrom * cMaxSave[iStepFrom]
872  + eStepTo * cMaxSave[iStepTo];
873 
874 }
875 
876 //--------------------------------------------------------------------------
877 
878 // Select first = hardest pT in nondiffractive process.
879 // Requires separate treatment at low and high b values.
880 
881 void MultipartonInteractions::pTfirst() {
882 
883  // Pick impact parameter and thereby interaction rate enhancement.
884  // This is not used for the x-dependent matter profile, which
885  // instead uses trial interactions.
886  if (bProfile == 4) isAtLowB = false;
887  else overlapFirst();
888  bSetInFirst = true;
889  double WTacc;
890 
891  // At low b values evolve downwards with Sudakov.
892  if (isAtLowB) {
893  pT2 = pT2max;
894  do {
895 
896  // Pick a pT using a quick-and-dirty cross section estimate.
897  pT2 = fastPT2(pT2);
898 
899  // If fallen below lower cutoff then need to restart at top.
900  if (pT2 < pT2min) {
901  pT2 = pT2max;
902  WTacc = 0.;
903 
904  // Else pick complete kinematics and evaluate cross-section correction.
905  } else {
906  WTacc = sigmaPT2scatter(true) / dSigmaApprox;
907  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
908  "MultipartonInteractions::pTfirst: weight above unity");
909  }
910 
911  // Loop until acceptable pT and acceptable kinematics.
912  } while (WTacc < rndmPtr->flat() || !dSigmaDtSel->final2KinMPI());
913 
914  // At high b values make preliminary pT choice without Sudakov factor.
915  } else {
916 
917  while (true) {
918  do {
919  pT2 = pT20min0maxR / (pT20minR + rndmPtr->flat() * pT2maxmin) - pT20R;
920 
921  // Evaluate upper estimate of cross section for this pT2 choice.
922  dSigmaApprox = pT4dSigmaMax / pow2(pT2 + pT20R);
923 
924  // Pick complete kinematics and evaluate cross-section correction.
925  WTacc = sigmaPT2scatter(true) / dSigmaApprox;
926 
927  // Evaluate and include Sudakov factor.
928  if (bProfile != 4) WTacc *= sudakov( pT2, enhanceB);
929 
930  // Warn for weight above unity
931  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
932  "MultipartonInteractions::pTfirst: weight above unity");
933 
934  // Loop until acceptable pT and acceptable kinematics.
935  } while (WTacc < rndmPtr->flat() || !dSigmaDtSel->final2KinMPI());
936 
937  // For x-dependent matter profile, use trial interactions to
938  // generate Sudakov, otherwise done.
939  if (bProfile != 4) break;
940  else {
941  // Save details of the original hard interaction.
942  pT2Save = pT2;
943  id1Save = id1;
944  id2Save = id2;
945  x1Save = x1;
946  x2Save = x2;
947  sHatSave = sHat;
948  tHatSave = tHat;
949  uHatSave = uHat;
950  alpSsave = alpS;
951  alpEMsave = alpEM;
952  pT2FacSave = pT2Fac;
953  pT2RenSave = pT2Ren;
954  xPDF1nowSave = xPDF1now;
955  xPDF2nowSave = xPDF2now;
956  // Save accepted kinematics and pointer to SigmaProcess.
957  dSigmaDtSel->saveKin();
958  dSigmaDtSelSave = dSigmaDtSel;
959 
960  // Put x1, x2 information into beam pointers to get correct
961  // PDF rescaling in trial interaction (for hard process, can
962  // be sea or valence, not companion).
963  beamAPtr->append( 0, id1, x1);
964  beamAPtr->xfISR( 0, id1, x1, pT2Fac * pT2Fac);
965  vsc1 = beamAPtr->pickValSeaComp();
966  beamBPtr->append( 0, id2, x2);
967  beamBPtr->xfISR( 0, id2, x2, pT2Fac * pT2Fac);
968  vsc2 = beamBPtr->pickValSeaComp();
969 
970  // Pick b according to O(b, x1, x2).
971  double w1 = XDEP_A1 + a1 * log(1. / x1);
972  double w2 = XDEP_A1 + a1 * log(1. / x2);
973  double fac = a02now * (w1 * w1 + w2 * w2);
974  double expb2 = 1.;
975  if ( userHooksPtr && userHooksPtr->canSetImpactParameter() ) {
976  bNow = userHooksPtr->doSetImpactParameter() * bAvg;
977  b2now = pow2(bNow);
978  expb2 = exp(-b2now / fac);
979  } else {
980  expb2 = rndmPtr->flat();
981  b2now = - fac * log(expb2);
982  bNow = sqrt(b2now);
983  }
984 
985  // Enhancement factor for the hard process and overestimate
986  // for fastPT2. Note that existing framework has a (1. / sigmaND)
987  // present.
988  enhanceB = sigmaND / M_PI / fac * expb2;
989  enhanceBmax = sigmaND / 2. / M_PI / a02now *
990  exp( -b2now / 2. / a2max );
991 
992  // Trial interaction with dummy event.
993  Event evDummy;
994  double pTtrial = pTnext(pTmax, pTmin, evDummy);
995 
996  // Restore beams.
997  beamAPtr->clear();
998  beamBPtr->clear();
999 
1000  // Accept if fallen beneath factorisation scale.
1001  if (pTtrial < sqrt(pT2FacSave)) {
1002  // Restore previous values and original sigma.
1003  swap(pT2, pT2Save);
1004  swap(pT2Fac, pT2FacSave);
1005  swap(pT2Ren, pT2RenSave);
1006  swap(id1, id1Save);
1007  swap(id2, id2Save);
1008  swap(x1, x1Save);
1009  swap(x2, x2Save);
1010  swap(sHat, sHatSave);
1011  swap(tHat, tHatSave);
1012  swap(uHat, uHatSave);
1013  swap(alpS, alpSsave);
1014  swap(alpEM, alpEMsave);
1015  swap(xPDF1now, xPDF1nowSave);
1016  swap(xPDF2now, xPDF2nowSave);
1017  if (dSigmaDtSel == dSigmaDtSelSave) dSigmaDtSel->swapKin();
1018  else swap(dSigmaDtSel, dSigmaDtSelSave);
1019 
1020  // Accept.
1021  bNow /= bAvg;
1022  bIsSet = true;
1023  break;
1024  }
1025  } // if (bProfile == 4)
1026  } // while (true)
1027 
1028  // End handling for high b.
1029  }
1030 
1031 }
1032 
1033 //--------------------------------------------------------------------------
1034 
1035 // Set up kinematics for first = hardest pT in nondiffractive process.
1036 
1037 void MultipartonInteractions::setupFirstSys( Event& process) {
1038 
1039  // Last beam-status particles. Offset relative to normal beam locations.
1040  int sizeProc = process.size();
1041  int nBeams = 3;
1042  for (int i = 3; i < sizeProc; ++i)
1043  if (process[i].statusAbs() < 20) nBeams = i + 1;
1044  int nOffset = nBeams - 3;
1045 
1046  // Remove any partons of previous failed interactions.
1047  if (sizeProc > nBeams) {
1048  process.popBack( sizeProc - nBeams);
1049  process.initColTag();
1050  }
1051 
1052  // Entries 3 and 4, now to be added, come from 1 and 2.
1053  process[1 + nOffset].daughter1(3 + nOffset);
1054  process[2 + nOffset].daughter1(4 + nOffset);
1055 
1056  // Negate beam status, if not already done. (Case with offset beams.)
1057  process[1 + nOffset].statusNeg();
1058  process[2 + nOffset].statusNeg();
1059 
1060  // Loop over four partons and offset info relative to subprocess itself.
1061  int colOffset = process.lastColTag();
1062  for (int i = 1; i <= 4; ++i) {
1063  Particle parton = dSigmaDtSel->getParton(i);
1064  if (i <= 2) parton.status(-21);
1065  else parton.status(23);
1066  if (i <= 2) parton.mothers( i + nOffset, 0);
1067  else parton.mothers( 3 + nOffset, 4 + nOffset);
1068  if (i <= 2) parton.daughters( 5 + nOffset, 6 + nOffset);
1069  else parton.daughters( 0, 0);
1070  int col = parton.col();
1071  if (col > 0) parton.col( col + colOffset);
1072  int acol = parton.acol();
1073  if (acol > 0) parton.acol( acol + colOffset);
1074 
1075  // Put the partons into the event record.
1076  process.append(parton);
1077  }
1078 
1079  // Set scale from which to begin evolution.
1080  process.scale( sqrt(pT2Fac) );
1081 
1082  // Info on subprocess - specific to mimimum-bias events.
1083  string nameSub = dSigmaDtSel->name();
1084  int codeSub = dSigmaDtSel->code();
1085  int nFinalSub = dSigmaDtSel->nFinal();
1086  double pTMPI = dSigmaDtSel->pTMPIFin();
1087  infoPtr->setSubType( iDiffSys, nameSub, codeSub, nFinalSub);
1088  if (iDiffSys == 0) infoPtr->setTypeMPI( codeSub, pTMPI, 0, 0,
1089  enhanceB / zeroIntCorr);
1090 
1091  // Further standard info on process.
1092  infoPtr->setPDFalpha( iDiffSys, id1, id2, x1, x2,
1093  (id1 == 21 ? 4./9. : 1.) * xPDF1now, (id2 == 21 ? 4./9. : 1.) * xPDF2now,
1094  pT2Fac, alpEM, alpS, pT2Ren, 0.);
1095  double m3 = dSigmaDtSel->m(3);
1096  double m4 = dSigmaDtSel->m(4);
1097  double theta = dSigmaDtSel->thetaMPI();
1098  double phi = dSigmaDtSel->phiMPI();
1099  infoPtr->setKin( iDiffSys, id1, id2, x1, x2, sHat, tHat, uHat, sqrt(pT2),
1100  m3, m4, theta, phi);
1101 
1102 }
1103 
1104 //--------------------------------------------------------------------------
1105 
1106 // Find whether to limit maximum scale of emissions.
1107 
1108 bool MultipartonInteractions::limitPTmax( Event& event) {
1109 
1110  // User-set cases.
1111  if (pTmaxMatch == 1) return true;
1112  if (pTmaxMatch == 2) return false;
1113 
1114  // Always restrict SoftQCD processes.
1115  if (infoPtr->isNonDiffractive() || infoPtr->isDiffractiveA()
1116  || infoPtr->isDiffractiveB() || infoPtr->isDiffractiveC() ) return true;
1117 
1118  // Look if only quarks (u, d, s, c, b), gluons and photons in final state.
1119  bool onlyQGP1 = true;
1120  bool onlyQGP2 = true;
1121  double scaleLimit1 = 0.;
1122  double scaleLimit2 = 0.;
1123  int n21 = 0;
1124  int iBegin = 5 + beamOffset;
1125  for (int i = iBegin; i < event.size(); ++i) {
1126  if (event[i].status() == -21) ++n21;
1127  else if (n21 == 0) {
1128  scaleLimit1 += 0.5 * event[i].pT();
1129  int idAbs = event[i].idAbs();
1130  if (idAbs > 5 && idAbs != 21 && idAbs != 22) onlyQGP1 = false;
1131  } else if (n21 == 2) {
1132  scaleLimit2 += 0.5 * event[i].pT();
1133  int idAbs = event[i].idAbs();
1134  if (idAbs > 5 && idAbs != 21 && idAbs != 22) onlyQGP2 = false;
1135  }
1136  }
1137 
1138  // If two hard interactions then limit if one only contains q/g/gamma.
1139  scaleLimitPTsave = (n21 == 2) ? min( scaleLimit1, scaleLimit2) : scaleLimit1;
1140  bool onlyQGP = (n21 == 2) ? (onlyQGP1 || onlyQGP2) : onlyQGP1;
1141  return (onlyQGP);
1142 
1143 }
1144 
1145 //--------------------------------------------------------------------------
1146 
1147 // Select next pT in downwards evolution.
1148 
1149 double MultipartonInteractions::pTnext( double pTbegAll, double pTendAll,
1150  Event& event) {
1151 
1152  // Initial values.
1153  bool pickRescatter, acceptKin;
1154  double dSigmaScatter, dSigmaRescatter, WTacc;
1155  double pT2end = pow2( max(pTmin, pTendAll) );
1156 
1157  // With the x-dependent matter profile, and minimum bias or diffraction,
1158  // it is possible to reuse values that have been stored during the trial
1159  // interactions for a slight speedup.
1160  // bIsSet is false during trial interactions, counter 21 in case partonLevel
1161  // is retried and counter 22 for the first pass through partonLevel.
1162  if (bProfile == 4 && bIsSet && bSetInFirst && infoPtr->getCounter(21) == 1
1163  && infoPtr->getCounter(22) == 1) {
1164  if (pT2Save < pT2end) return 0.;
1165  pT2 = pT2Save;
1166  pT2Fac = pT2FacSave;
1167  pT2Ren = pT2RenSave;
1168  id1 = id1Save;
1169  id2 = id2Save;
1170  x1 = x1Save;
1171  x2 = x2Save;
1172  sHat = sHatSave;
1173  tHat = tHatSave;
1174  uHat = uHatSave;
1175  alpS = alpSsave;
1176  alpEM = alpEMsave;
1177  xPDF1now = xPDF1nowSave;
1178  xPDF2now = xPDF2nowSave;
1179  if (dSigmaDtSel == dSigmaDtSelSave) dSigmaDtSel->swapKin();
1180  else dSigmaDtSel = dSigmaDtSelSave;
1181  return sqrt(pT2);
1182  }
1183 
1184  // Do not allow rescattering while still FSR with global recoil.
1185  bool allowRescatterNow = allowRescatter;
1186  if (globalRecoilFSR && partonSystemsPtr->sizeOut(0) <= nMaxGlobalRecoilFSR)
1187  allowRescatterNow = false;
1188 
1189  // Initial pT2 value.
1190  pT2 = pow2(pTbegAll);
1191 
1192  // Find the set of already scattered partons on the two sides.
1193  if (allowRescatterNow) findScatteredPartons( event);
1194 
1195  // Pick a pT2 using a quick-and-dirty cross section estimate.
1196  do {
1197  do {
1198  pT2 = fastPT2(pT2);
1199  if (pT2 < pT2end) return 0.;
1200 
1201  // Initial values: no rescattering.
1202  i1Sel = 0;
1203  i2Sel = 0;
1204  dSigmaSum = 0.;
1205  pickRescatter = false;
1206 
1207  // Pick complete kinematics and evaluate interaction cross-section.
1208  dSigmaScatter = sigmaPT2scatter(false);
1209 
1210  // Also cross section from rescattering if allowed.
1211  dSigmaRescatter = (allowRescatterNow) ? sigmaPT2rescatter( event) : 0.;
1212 
1213  // Normalize to dSigmaApprox, which was set in fastPT2 above.
1214  WTacc = (dSigmaScatter + dSigmaRescatter) / dSigmaApprox;
1215  if (WTacc > WTACCWARN) infoPtr->errorMsg("Warning in "
1216  "MultipartonInteractions::pTnext: weight above unity");
1217 
1218  // Idea suggested by Gosta Gustafson: increased screening in events
1219  // with large activity can be simulated by pT0_eff = sqrt(n) * pT0.
1220  if (enhanceScreening > 0) {
1221  int nSysNow = infoPtr->nMPI() + 1;
1222  if (enhanceScreening == 2) nSysNow += infoPtr->nISR();
1223  double WTscreen = pow2( (pT2 + pT20) / (pT2 + nSysNow * pT20) );
1224  WTacc *= WTscreen;
1225  }
1226 
1227  // x-dependent matter profile overlap weighting.
1228  if (bProfile == 4) {
1229  double w1 = XDEP_A1 + a1 * log(1. / x1);
1230  double w2 = XDEP_A1 + a1 * log(1. / x2);
1231  double fac = a02now * (w1 * w1 + w2 * w2);
1232  // Correct enhancement factor and weight
1233  enhanceBnow = sigmaND / M_PI / fac * exp( - b2now / fac);
1234  double oWgt = enhanceBnow / enhanceBmax;
1235  if (oWgt > 1.0000000001) infoPtr->errorMsg("Warning in Multiparton"
1236  "Interactions::pTnext: overlap weight above unity");
1237  WTacc *= oWgt;
1238  }
1239 
1240  // Decide whether to keep the event based on weight.
1241  } while (WTacc < rndmPtr->flat());
1242 
1243  // When rescattering possible: new interaction or rescattering?
1244  if (allowRescatterNow) {
1245  pickRescatter = (i1Sel > 0 || i2Sel > 0);
1246 
1247  // Restore kinematics for selected scattering/rescattering.
1248  id1 = id1Sel;
1249  id2 = id2Sel;
1250  x1 = x1Sel;
1251  x2 = x2Sel;
1252  sHat = sHatSel;
1253  tHat = tHatSel;
1254  uHat = uHatSel;
1255  sigma2Sel->sigma( id1, id2, x1, x2, sHat, tHat, uHat, alpS, alpEM,
1256  true, pickOtherSel);
1257  }
1258 
1259  // Pick one of the possible channels summed above.
1260  dSigmaDtSel = sigma2Sel->sigmaSel();
1261  if (sigma2Sel->swapTU()) swap( tHat, uHat);
1262 
1263  // Decide to keep event based on kinematics (normally always OK).
1264  // Rescattering: need to provide incoming four-vectors and masses.
1265  if (pickRescatter) {
1266  Vec4 p1Res = (i1Sel == 0) ? 0.5 * eCM * x1Sel * Vec4( 0., 0., 1., 1.)
1267  : event[i1Sel].p();
1268  Vec4 p2Res = (i2Sel == 0) ? 0.5 * eCM * x2Sel * Vec4( 0., 0., -1., 1.)
1269  : event[i2Sel].p();
1270  double m1Res = (i1Sel == 0) ? 0. : event[i1Sel].m();
1271  double m2Res = (i2Sel == 0) ? 0. : event[i2Sel].m();
1272  acceptKin = dSigmaDtSel->final2KinMPI( i1Sel, i2Sel, p1Res, p2Res,
1273  m1Res, m2Res);
1274  // New interaction: already stored values enough.
1275  } else acceptKin = dSigmaDtSel->final2KinMPI();
1276  } while (!acceptKin);
1277 
1278  // Done.
1279  return sqrt(pT2);
1280 
1281 }
1282 
1283 //--------------------------------------------------------------------------
1284 
1285 // Set up the kinematics of the 2 -> 2 scattering process,
1286 // and store the scattering in the event record.
1287 
1288 bool MultipartonInteractions::scatter( Event& event) {
1289 
1290  // Last beam-status particles. Offset relative to normal beam locations.
1291  int sizeProc = event.size();
1292  int nBeams = 3;
1293  for (int i = 3; i < sizeProc; ++i)
1294  if (event[i].statusAbs() < 20) nBeams = i + 1;
1295  int nOffset = nBeams - 3;
1296 
1297  // Loop over four partons and offset info relative to subprocess itself.
1298  int motherOffset = event.size();
1299  int colOffset = event.lastColTag();
1300  for (int i = 1; i <= 4; ++i) {
1301  Particle parton = dSigmaDtSel->getParton(i);
1302  if (i <= 2 ) parton.mothers( i + nOffset, 0);
1303  else parton.mothers( motherOffset, motherOffset + 1);
1304  if (i <= 2 ) parton.daughters( motherOffset + 2, motherOffset + 3);
1305  else parton.daughters( 0, 0);
1306  int col = parton.col();
1307  if (col > 0) parton.col( col + colOffset);
1308  int acol = parton.acol();
1309  if (acol > 0) parton.acol( acol + colOffset);
1310 
1311  // Put the partons into the event record.
1312  event.append(parton);
1313  }
1314 
1315  // Allow setting of new parton production vertices.
1316  if (doPartonVertex)
1317  partonVertexPtr->vertexMPI( sizeProc, 4, bNow, event);
1318 
1319  // Allow veto of MPI. If so restore event record to before scatter.
1320  if (canVetoMPI && userHooksPtr->doVetoMPIEmission(sizeProc, event)) {
1321  event.popBack(event.size() - sizeProc);
1322  return false;
1323  }
1324 
1325  // Store participating partons as a new set in list of all systems.
1326  int iSys = partonSystemsPtr->addSys();
1327  partonSystemsPtr->setInA(iSys, motherOffset);
1328  partonSystemsPtr->setInB(iSys, motherOffset + 1);
1329  partonSystemsPtr->addOut(iSys, motherOffset + 2);
1330  partonSystemsPtr->addOut(iSys, motherOffset + 3);
1331  partonSystemsPtr->setSHat(iSys, sHat);
1332 
1333  // Tag double rescattering graphs that annihilate one initial colour.
1334  bool annihil1 = false;
1335  bool annihil2 = false;
1336  if (i1Sel > 0 && i2Sel > 0) {
1337  if (event[motherOffset].col() == event[motherOffset + 1].acol()
1338  && event[motherOffset].col() > 0) annihil1 = true;
1339  if (event[motherOffset].acol() == event[motherOffset + 1].col()
1340  && event[motherOffset].acol() > 0) annihil2 = true;
1341  }
1342 
1343  // Beam remnant A: add scattered partons to list.
1344  BeamParticle& beamA = *beamAPtr;
1345  int iA = beamA.append( motherOffset, id1, x1);
1346 
1347  // Find whether incoming partons are valence or sea, so prepared for ISR.
1348  if (i1Sel == 0) {
1349  beamA.xfISR( iA, id1, x1, pT2);
1350  beamA.pickValSeaComp();
1351 
1352  // Remove rescattered parton from final state and change history.
1353  // Propagate existing colour labels throught graph.
1354  } else {
1355  beamA[iA].companion(-10);
1356  event[i1Sel].statusNeg();
1357  event[i1Sel].daughters( motherOffset, motherOffset);
1358  event[motherOffset].mothers( i1Sel, i1Sel);
1359  int colOld = event[i1Sel].col();
1360  if (colOld > 0) {
1361  int colNew = event[motherOffset].col();
1362  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1363  if (event[i].col() == colNew) event[i].col( colOld);
1364  if (event[i].acol() == colNew) event[i].acol( colOld);
1365  }
1366  }
1367  int acolOld = event[i1Sel].acol();
1368  if (acolOld > 0) {
1369  int acolNew = event[motherOffset].acol();
1370  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1371  if (event[i].col() == acolNew) event[i].col( acolOld);
1372  if (event[i].acol() == acolNew) event[i].acol( acolOld);
1373  }
1374  }
1375  }
1376 
1377  // Beam remnant B: add scattered partons to list.
1378  BeamParticle& beamB = *beamBPtr;
1379  int iB = beamB.append( motherOffset + 1, id2, x2);
1380 
1381  // Find whether incoming partons are valence or sea, so prepared for ISR.
1382  if (i2Sel == 0) {
1383  beamB.xfISR( iB, id2, x2, pT2);
1384  beamB.pickValSeaComp();
1385 
1386  // Remove rescattered parton from final state and change history.
1387  // Propagate existing colour labels throught graph.
1388  } else {
1389  beamB[iB].companion(-10);
1390  event[i2Sel].statusNeg();
1391  event[i2Sel].daughters( motherOffset + 1, motherOffset + 1);
1392  event[motherOffset + 1].mothers( i2Sel, i2Sel);
1393  int colOld = event[i2Sel].col();
1394  if (colOld > 0) {
1395  int colNew = event[motherOffset + 1].col();
1396  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1397  if (event[i].col() == colNew) event[i].col( colOld);
1398  if (event[i].acol() == colNew) event[i].acol( colOld);
1399  }
1400  }
1401  int acolOld = event[i2Sel].acol();
1402  if (acolOld > 0) {
1403  int acolNew = event[motherOffset + 1].acol();
1404  for (int i = motherOffset; i < motherOffset + 4; ++i) {
1405  if (event[i].col() == acolNew) event[i].col( acolOld);
1406  if (event[i].acol() == acolNew) event[i].acol( acolOld);
1407  }
1408  }
1409  }
1410 
1411  // Annihilating colour in double rescattering requires relabelling
1412  // of one colour into the other in the whole preceding event.
1413  if (annihil1 || annihil2) {
1414  int colLeft = (annihil1) ? event[motherOffset].col()
1415  : event[motherOffset].acol();
1416  int mother1 = event[motherOffset].mother1();
1417  int mother2 = event[motherOffset + 1].mother1();
1418  int colLost = (annihil1)
1419  ? event[mother1].col() + event[mother2].acol() - colLeft
1420  : event[mother1].acol() + event[mother2].col() - colLeft;
1421  for (int i = 0; i < motherOffset; ++i) {
1422  if (event[i].col() == colLost) event[i].col( colLeft );
1423  if (event[i].acol() == colLost) event[i].acol( colLeft );
1424  }
1425  }
1426 
1427  // With gamma+gamma check that room for beam remnants for current scattering.
1428  // Otherwise take the partons out from event record.
1429  // roomForRemnants treats both beam equally so need to do only once.
1430  if ( beamAPtr->isGamma() || beamBPtr->isGamma() ) {
1431  if ( !beamAPtr->roomForRemnants(*beamBPtr) ) {
1432  // Remove the partons associated to the latest scattering from the
1433  // event record.
1434  event.popBack(4);
1435  beamAPtr->popBack();
1436  beamBPtr->popBack();
1437  partonSystemsPtr->popBack();
1438 
1439  infoPtr->errorMsg("Warning in MultipartonInteractions::scatter:"
1440  " No room for remnants for given scattering");
1441  return false;
1442  }
1443  }
1444 
1445  // Store the pT value for valence decision of resolved photons.
1446  beamA.pTMPI(sqrtpos(pT2));
1447  beamB.pTMPI(sqrtpos(pT2));
1448 
1449  // Store info on subprocess code and rescattered partons.
1450  int codeMPI = dSigmaDtSel->code();
1451  double pTMPI = dSigmaDtSel->pTMPIFin();
1452  if (iDiffSys == 0) infoPtr->setTypeMPI( codeMPI, pTMPI, i1Sel, i2Sel,
1453  enhanceBnow / zeroIntCorr);
1454  partonSystemsPtr->setPTHat(iSys, pTMPI);
1455 
1456  // Done.
1457  return true;
1458 }
1459 
1460 //--------------------------------------------------------------------------
1461 
1462 // Determine constant in d(Prob)/d(pT2) < const / (pT2 + r * pT20)^2.
1463 
1464 void MultipartonInteractions::upperEnvelope() {
1465 
1466  // Initially determine constant in jet cross section upper estimate
1467  // d(sigma_approx)/d(pT2) < const / (pT2 + r * pT20)^2.
1468  pT4dSigmaMax = 0.;
1469 
1470  // Loop thorough allowed pT range logarithmically evenly.
1471  for (int iPT = 0; iPT < 100; ++iPT) {
1472  double pT = pTmin * pow( pTmax/pTmin, 0.01 * (iPT + 0.5) );
1473  pT2 = pT*pT;
1474  pT2shift = pT2 + pT20;
1475  pT2Ren = pT2shift;
1476  pT2Fac = (SHIFTFACSCALE) ? pT2shift : pT2;
1477  xT = 2. * pT / eCM;
1478 
1479  // Evaluate parton density sums at x1 = x2 = xT.
1480  double xPDF1sumMax = (9./4.) * beamAPtr->xf(21, xT, pT2Fac);
1481  for (int id = 1; id <= nQuarkIn; ++id)
1482  xPDF1sumMax += beamAPtr->xf( id, xT, pT2Fac)
1483  + beamAPtr->xf(-id, xT, pT2Fac);
1484  double xPDF2sumMax = (9./4.) * beamBPtr->xf(21, xT, pT2Fac);
1485  for (int id = 1; id <= nQuarkIn; ++id)
1486  xPDF2sumMax += beamBPtr->xf( id, xT, pT2Fac)
1487  + beamBPtr->xf(-id, xT, pT2Fac);
1488 
1489  // Evaluate alpha_strong and _EM, matrix element and phase space volume.
1490  alpS = alphaS.alphaS(pT2Ren);
1491  alpEM = alphaEM.alphaEM(pT2Ren);
1492  double dSigmaPartonApprox = CONVERT2MB * Kfactor * 0.5 * M_PI
1493  * pow2(alpS / pT2shift);
1494  double yMax = log(1./xT + sqrt(1./(xT*xT) - 1.));
1495  double volumePhSp = pow2(2. * yMax);
1496 
1497  // Final comparison to determine upper estimate.
1498  double dSigmaApproxNow = SIGMAFUDGE * xPDF1sumMax * xPDF2sumMax
1499  * dSigmaPartonApprox * volumePhSp;
1500  double pT4dSigmaNow = pow2(pT2 + pT20R) * dSigmaApproxNow;
1501  if ( pT4dSigmaNow > pT4dSigmaMax) pT4dSigmaMax = pT4dSigmaNow;
1502  }
1503 
1504  // Get wanted constant by dividing by the nondiffractive cross section.
1505  pT4dProbMax = pT4dSigmaMax / sigmaND;
1506 
1507 }
1508 
1509 //--------------------------------------------------------------------------
1510 
1511 // Integrate the parton-parton interaction cross section,
1512 // using stratified Monte Carlo sampling.
1513 // Store result in pT bins for use as Sudakov form factors.
1514 
1515 void MultipartonInteractions::jetCrossSection() {
1516 
1517  // Common factor from bin size in dpT2 / (pT2 + r * pT20)^2 and statistics.
1518  double sigmaFactor = (1. / pT20minR - 1. / pT20maxR) / (100. * nSample);
1519 
1520  // Reset overlap-weighted cross section for x-dependent matter profile.
1521  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++)
1522  sigmaIntWgt[bBin] = 0.;
1523 
1524  // Loop through allowed pT range evenly in dpT2 / (pT2 + r * pT20)^2.
1525  sigmaInt = 0.;
1526  double dSigmaMax = 0.;
1527  sudExpPT[100] = 0.;
1528 
1529  for (int iPT = 99; iPT >= 0; --iPT) {
1530  double sigmaSum = 0.;
1531 
1532  // Reset pT-binned overlap-weigted integration.
1533  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++)
1534  sigmaSumWgt[bBin] = 0.;
1535 
1536  // In each pT bin sample a number of random pT values.
1537  for (int iSample = 0; iSample < nSample; ++iSample) {
1538  double mappedPT2 = 1. - 0.01 * (iPT + rndmPtr->flat());
1539  pT2 = pT20min0maxR / (pT20minR + mappedPT2 * pT2maxmin) - pT20R;
1540 
1541  // Evaluate cross section dSigma/dpT2 in phase space point.
1542  double dSigma = sigmaPT2scatter(true);
1543 
1544  // Multiply by (pT2 + r * pT20)^2 to compensate for pT sampling. Sum.
1545  dSigma *= pow2(pT2 + pT20R);
1546  sigmaSum += dSigma;
1547  if (dSigma > dSigmaMax) dSigmaMax = dSigma;
1548 
1549  // Overlap-weighted cross section for x-dependent matter profile.
1550  // Note that dSigma can be 0. when points are rejected.
1551  if (bProfile == 4 && dSigma > 0.) {
1552  double w1 = XDEP_A1 + a1 * log(1. / x1);
1553  double w2 = XDEP_A1 + a1 * log(1. / x2);
1554  double fac = XDEP_A0 * XDEP_A0 * (w1 * w1 + w2 * w2);
1555  double b = 0.5 * bstepNow;
1556  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
1557  double wgt = exp( - b * b / fac ) / fac / M_PI;
1558  sigmaSumWgt[bBin] += dSigma * wgt;
1559  b += bstepNow;
1560  }
1561  }
1562  }
1563 
1564  // Store total cross section and exponent of Sudakov.
1565  sigmaSum *= sigmaFactor;
1566  sigmaInt += sigmaSum;
1567  sudExpPT[iPT] = sudExpPT[iPT + 1] + sigmaSum / sigmaND;
1568 
1569  // Sum overlap-weighted cross section
1570  if (bProfile == 4) for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
1571  sigmaSumWgt[bBin] *= sigmaFactor;
1572  sigmaIntWgt[bBin] += sigmaSumWgt[bBin];
1573  }
1574 
1575  // End of loop over pT values.
1576  }
1577 
1578  // Update upper estimate of differential cross section. Done.
1579  if (dSigmaMax > pT4dSigmaMax) {
1580  pT4dSigmaMax = dSigmaMax;
1581  pT4dProbMax = dSigmaMax / sigmaND;
1582  }
1583 
1584 }
1585 
1586 //--------------------------------------------------------------------------
1587 
1588 // Evaluate "Sudakov form factor" for not having a harder interaction
1589 // at the selected b value, given the pT scale of the event.
1590 
1591 double MultipartonInteractions::sudakov(double pT2sud, double enhance) {
1592 
1593  // Find bin the pT2 scale falls in.
1594  double xBin = (pT2sud - pT2min) * pT20maxR
1595  / (pT2maxmin * (pT2sud + pT20R));
1596  xBin = max(1e-6, min(100. - 1e-6, 100. * xBin) );
1597  int iBin = int(xBin);
1598 
1599  // Interpolate inside bin. Optionally include enhancement factor.
1600  double sudExp = sudExpPT[iBin] + (xBin - iBin)
1601  * (sudExpPT[iBin + 1] - sudExpPT[iBin]);
1602  return exp( -enhance * sudExp);
1603 
1604 }
1605 
1606 //--------------------------------------------------------------------------
1607 
1608 // Pick a trial next pT, based on a simple upper estimate of the
1609 // d(sigma)/d(pT2) spectrum.
1610 
1611 double MultipartonInteractions::fastPT2( double pT2beg) {
1612 
1613  // Use d(Prob)/d(pT2) < pT4dProbMax / (pT2 + r * pT20)^2.
1614  double pT20begR = pT2beg + pT20R;
1615  double pT4dProbMaxNow = pT4dProbMax * enhanceBmax;
1616  double pT2try = pT4dProbMaxNow * pT20begR
1617  / (pT4dProbMaxNow - pT20begR * log(rndmPtr->flat())) - pT20R;
1618 
1619  if ( pT2try + pT20R <= 0.0 ) return 0.0;
1620 
1621  // Save cross section associated with ansatz above. Done.
1622  dSigmaApprox = pT4dSigmaMax / pow2(pT2try + pT20R);
1623  return pT2try;
1624 
1625 }
1626 
1627 //--------------------------------------------------------------------------
1628 
1629 // Calculate the actual cross section to decide whether fast choice is OK.
1630 // Select flavours and kinematics for interaction at given pT.
1631 // Slightly different treatment for first interaction and subsequent ones.
1632 
1633 double MultipartonInteractions::sigmaPT2scatter(bool isFirst) {
1634 
1635  // Derive recormalization and factorization scales, amd alpha_strong/em.
1636  pT2shift = pT2 + pT20;
1637  pT2Ren = pT2shift;
1638  pT2Fac = (SHIFTFACSCALE) ? pT2shift : pT2;
1639  alpS = alphaS.alphaS(pT2Ren);
1640  alpEM = alphaEM.alphaEM(pT2Ren);
1641 
1642  // Derive rapidity limits from chosen pT2.
1643  xT = 2. * sqrt(pT2) / eCM;
1644  if (xT >= 1.) return 0.;
1645  xT2 = xT*xT;
1646  double yMax = log(1./xT + sqrt(1./xT2 - 1.));
1647 
1648  // Select rapidities y3 and y4 of the two produced partons.
1649  double y3 = yMax * (2. * rndmPtr->flat() - 1.);
1650  double y4 = yMax * (2. * rndmPtr->flat() - 1.);
1651  y = 0.5 * (y3 + y4);
1652 
1653  // Failure if x1 or x2 exceed what is left in respective beam.
1654  x1 = 0.5 * xT * (exp(y3) + exp(y4));
1655  x2 = 0.5 * xT * (exp(-y3) + exp(-y4));
1656  if (isFirst && iDiffSys == 0) {
1657  if (x1 > 1. || x2 > 1.) return 0.;
1658  } else {
1659  if (x1 > beamAPtr->xMax() || x2 > beamBPtr->xMax()) return 0.;
1660  }
1661  tau = x1 * x2;
1662 
1663  // Begin evaluate parton densities at actual x1 and x2.
1664  double xPDF1[21];
1665  double xPDF1sum = 0.;
1666  double xPDF2[21];
1667  double xPDF2sum = 0.;
1668 
1669  // For first interaction use normal densities.
1670  if (isFirst) {
1671  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1672  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xf(21, x1, pT2Fac);
1673  else xPDF1[id+10] = beamAPtr->xf(id, x1, pT2Fac);
1674  xPDF1sum += xPDF1[id+10];
1675  }
1676  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1677  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xf(21, x2, pT2Fac);
1678  else xPDF2[id+10] = beamBPtr->xf(id, x2, pT2Fac);
1679  xPDF2sum += xPDF2[id+10];
1680  }
1681 
1682  // For subsequent interactions use rescaled densities.
1683  } else {
1684  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1685  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xfMPI(21, x1, pT2Fac);
1686  else xPDF1[id+10] = beamAPtr->xfMPI(id, x1, pT2Fac);
1687  xPDF1sum += xPDF1[id+10];
1688  }
1689  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1690  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xfMPI(21, x2, pT2Fac);
1691  else xPDF2[id+10] = beamBPtr->xfMPI(id, x2, pT2Fac);
1692  xPDF2sum += xPDF2[id+10];
1693  }
1694  }
1695 
1696  // Select incoming flavours according to actual PDF's.
1697  id1 = -nQuarkIn - 1;
1698  double temp = xPDF1sum * rndmPtr->flat();
1699  do { xPDF1now = xPDF1[(++id1) + 10]; temp -= xPDF1now; }
1700  while (temp > 0. && id1 < nQuarkIn);
1701  if (id1 == 0) id1 = 21;
1702  id2 = -nQuarkIn-1;
1703  temp = xPDF2sum * rndmPtr->flat();
1704  do { xPDF2now = xPDF2[(++id2) + 10]; temp -= xPDF2now;}
1705  while (temp > 0. && id2 < nQuarkIn);
1706  if (id2 == 0) id2 = 21;
1707 
1708  // Check whether room for remnants left after scattering with photon beams.
1709  if ( isFirst && ( beamAPtr->isGamma() || beamBPtr->isGamma() ) ) {
1710  double mTRem = eCM * sqrt( (1 - x1) * (1 - x2) );
1711  double m1 = beamAPtr->remnantMass(id1);
1712  double m2 = beamBPtr->remnantMass(id2);
1713  if (mTRem < m1 + m2) return 0.;
1714  }
1715 
1716  // Assign pointers to processes relevant for incoming flavour choice:
1717  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1718  // Factor 4./9. per incoming gluon to compensate for preweighting.
1719  SigmaMultiparton* sigma2Tmp;
1720  double gluFac = 1.;
1721  if (id1 == 21 && id2 == 21) {
1722  sigma2Tmp = &sigma2gg;
1723  gluFac = 16. / 81.;
1724  } else if (id1 == 21 || id2 == 21) {
1725  sigma2Tmp = &sigma2qg;
1726  gluFac = 4. / 9.;
1727  } else if (id1 == -id2) sigma2Tmp = &sigma2qqbarSame;
1728  else sigma2Tmp = &sigma2qq;
1729 
1730  // Prepare to generate differential cross sections.
1731  sHat = tau * sCM;
1732  double root = sqrtpos(1. - xT2 / tau);
1733  tHat = -0.5 * sHat * (1. - root);
1734  uHat = -0.5 * sHat * (1. + root);
1735 
1736  // Evaluate cross sections, include possibility of K factor.
1737  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1738  // (Not necessary, but makes for better MC efficiency.)
1739  double dSigmaPartonCorr = Kfactor * gluFac
1740  * sigma2Tmp->sigma( id1, id2, x1, x2, sHat, tHat, uHat, alpS, alpEM);
1741 
1742  // Combine cross section, pdf's and phase space integral.
1743  double volumePhSp = pow2(2. * yMax);
1744  double dSigmaScat = dSigmaPartonCorr * xPDF1sum * 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  dSigmaScat *= pow2( pT2massive / (pT2massive + pT20) );
1751 
1752  // Sum up total contribution for all scatterings and rescatterings.
1753  dSigmaSum += dSigmaScat;
1754 
1755  // Save values for comparison with rescattering processes.
1756  i1Sel = 0;
1757  i2Sel = 0;
1758  id1Sel = id1;
1759  id2Sel = id2;
1760  x1Sel = x1;
1761  x2Sel = x2;
1762  sHatSel = sHat;
1763  tHatSel = tHat;
1764  uHatSel = uHat;
1765  sigma2Sel = sigma2Tmp;
1766  pickOtherSel = sigma2Tmp->pickedOther();
1767 
1768  // For first interaction: pick one of the possible channels summed above.
1769  if (isFirst) {
1770  dSigmaDtSel = sigma2Tmp->sigmaSel();
1771  if (sigma2Tmp->swapTU()) swap( tHat, uHat);
1772  }
1773 
1774  // Done.
1775  return dSigmaScat;
1776 }
1777 
1778 //--------------------------------------------------------------------------
1779 
1780 // Find the partons that are allowed to rescatter.
1781 
1782 void MultipartonInteractions::findScatteredPartons( Event& event) {
1783 
1784  // Reset arrays.
1785  scatteredA.resize(0);
1786  scatteredB.resize(0);
1787  double yTmp, probA, probB;
1788 
1789  // Loop though the event record and catch "final" partons.
1790  for (int i = 0; i < event.size(); ++i)
1791  if (event[i].isFinal() && (event[i].idAbs() <= nQuarkIn
1792  || event[i].id() == 21)) {
1793  yTmp = event[i].y();
1794 
1795  // Different strategies to determine which partons may rescatter.
1796  switch(rescatterMode) {
1797 
1798  // Case 0: step function at origin
1799  case 0:
1800  if ( yTmp > 0.) scatteredA.push_back( i);
1801  if (-yTmp > 0.) scatteredB.push_back( i);
1802  break;
1803 
1804  // Case 1: step function as ySepResc.
1805  case 1:
1806  if ( yTmp > ySepResc) scatteredA.push_back( i);
1807  if (-yTmp > ySepResc) scatteredB.push_back( i);
1808  break;
1809 
1810  // Case 2: linear rise from ySep - deltaY to ySep + deltaY.
1811  case 2:
1812  probA = 0.5 * (1. + ( yTmp - ySepResc) / deltaYResc);
1813  if (probA > rndmPtr->flat()) scatteredA.push_back( i);
1814  probB = 0.5 * (1. + (-yTmp - ySepResc) / deltaYResc);
1815  if (probB > rndmPtr->flat()) scatteredB.push_back( i);
1816  break;
1817 
1818  // Case 3: rise like (1/2) * ( 1 + tanh((y - ySep) / deltaY) ).
1819  // Use that (1/2) (1 + tanh(x)) = 1 / (1 + exp(-2x)).
1820  case 3:
1821  probA = 1. / (1. + exp(-2. * ( yTmp - ySepResc) / deltaYResc));
1822  if (probA > rndmPtr->flat()) scatteredA.push_back( i);
1823  probB = 1. / (1. + exp(-2. * (-yTmp - ySepResc) / deltaYResc));
1824  if (probB > rndmPtr->flat()) scatteredB.push_back( i);
1825  break;
1826 
1827  // Case 4 and undefined values: all partons can rescatter.
1828  default:
1829  scatteredA.push_back( i);
1830  scatteredB.push_back( i);
1831  break;
1832 
1833  // End of loop over partons. Done.
1834  }
1835  }
1836 
1837 }
1838 
1839 //--------------------------------------------------------------------------
1840 
1841 // Rescattering contribution summed over all already scattered partons.
1842 // Calculate the actual cross section to decide whether fast choice is OK.
1843 // Select flavours and kinematics for interaction at given pT.
1844 
1845 double MultipartonInteractions::sigmaPT2rescatter( Event& event) {
1846 
1847  // Derive xT scale from chosen pT2.
1848  xT = 2. * sqrt(pT2) / eCM;
1849  if (xT >= 1.) return 0.;
1850  xT2 = xT*xT;
1851 
1852  // Pointer to cross section and total rescatter contribution.
1853  SigmaMultiparton* sigma2Tmp;
1854  double dSigmaResc = 0.;
1855 
1856  // Normally save time by picking one random scattered parton from side A.
1857  int nScatA = scatteredA.size();
1858  int iScatA = -1;
1859  if (PREPICKRESCATTER && nScatA > 0) iScatA = min( nScatA - 1,
1860  int( rndmPtr->flat() * double(nScatA) ) );
1861 
1862  // Loop over all already scattered partons from side A.
1863  for (int iScat = 0; iScat < nScatA; ++iScat) {
1864  if (PREPICKRESCATTER && iScat != iScatA) continue;
1865  int iA = scatteredA[iScat];
1866  int id1Tmp = event[iA].id();
1867 
1868  // Calculate maximum possible sHat and check whether big enough.
1869  double x1Tmp = event[iA].pPos() / eCM;
1870  double sHatMax = x1Tmp * beamBPtr->xMax() * sCM;
1871  if (sHatMax < 4. * pT2) continue;
1872 
1873  // Select rapidity separation between two produced partons.
1874  double dyMax = log( sqrt(0.25 * sHatMax / pT2)
1875  * (1. + sqrtpos(1. - 4. * pT2 / sHatMax)) );
1876  double dy = dyMax * (2. * rndmPtr->flat() - 1.);
1877 
1878  // Reconstruct kinematical variables, especially x2.
1879  // Incoming c/b masses handled better if tau != x1 * x2.
1880  double m2Tmp = event[iA].m2();
1881  double sHatTmp = m2Tmp + 4. * pT2 * pow2(cosh(dy));
1882  double x2Tmp = (sHatTmp - m2Tmp) /(x1Tmp * sCM);
1883  double tauTmp = sHatTmp / sCM;
1884  double root = sqrtpos(1. - xT2 / tauTmp);
1885  double tHatTmp = -0.5 * sHatTmp * (1. - root);
1886  double uHatTmp = -0.5 * sHatTmp * (1. + root);
1887 
1888  // Begin evaluate parton densities at actual x2.
1889  double xPDF2[21];
1890  double xPDF2sum = 0.;
1891 
1892  // Use rescaled densities, with preweighting 9/4 for gluons.
1893  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1894  if (id == 0) xPDF2[10] = (9./4.) * beamBPtr->xfMPI(21, x2Tmp, pT2Fac);
1895  else xPDF2[id+10] = beamBPtr->xfMPI(id, x2Tmp, pT2Fac);
1896  xPDF2sum += xPDF2[id+10];
1897  }
1898 
1899  // Select incoming flavour according to actual PDF's.
1900  int id2Tmp = -nQuarkIn - 1;
1901  double temp = xPDF2sum * rndmPtr->flat();
1902  do { xPDF2now = xPDF2[(++id2Tmp) + 10]; temp -= xPDF2now;}
1903  while (temp > 0. && id2Tmp < nQuarkIn);
1904  if (id2Tmp == 0) id2Tmp = 21;
1905 
1906  // Assign pointers to processes relevant for incoming flavour choice:
1907  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
1908  // Factor 4./9. for incoming gluon 2 to compensate for preweighting.
1909  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
1910  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
1911  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
1912  else sigma2Tmp = &sigma2qq;
1913  double gluFac = (id2Tmp == 21) ? 4. / 9. : 1.;
1914 
1915  // Evaluate cross sections, include possibility of K factor.
1916  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
1917  // (Not necessary, but makes for better MC efficiency.)
1918  double dSigmaPartonCorr = Kfactor * gluFac
1919  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
1920  uHatTmp, alpS, alpEM);
1921 
1922  // Combine cross section, pdf's and phase space integral.
1923  double volumePhSp = 4. *dyMax;
1924  double dSigmaCorr = dSigmaPartonCorr * xPDF2sum * volumePhSp;
1925 
1926  // Dampen cross section at small pT values; part of formalism.
1927  // Use pT2 corrected for massive kinematics at this step.
1928  //?? double pT2massive = dSigmaDtSel->pT2MPI();
1929  double pT2massive = pT2;
1930  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
1931 
1932  // Compensate for prepicked rescattering if appropriate.
1933  if (PREPICKRESCATTER) dSigmaCorr *= nScatA;
1934 
1935  // Sum up total contribution for all scatterings or rescattering only.
1936  dSigmaSum += dSigmaCorr;
1937  dSigmaResc += dSigmaCorr;
1938 
1939  // Determine whether current rescattering should be selected.
1940  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
1941  i1Sel = iA;
1942  i2Sel = 0;
1943  id1Sel = id1Tmp;
1944  id2Sel = id2Tmp;
1945  x1Sel = x1Tmp;
1946  x2Sel = x2Tmp;
1947  sHatSel = sHatTmp;
1948  tHatSel = tHatTmp;
1949  uHatSel = uHatTmp;
1950  sigma2Sel = sigma2Tmp;
1951  pickOtherSel = sigma2Tmp->pickedOther();
1952  }
1953  }
1954 
1955  // Normally save time by picking one random scattered parton from side B.
1956  int nScatB = scatteredB.size();
1957  int iScatB = -1;
1958  if (PREPICKRESCATTER && nScatB > 0) iScatB = min( nScatB - 1,
1959  int( rndmPtr->flat() * double(nScatB) ) );
1960 
1961  // Loop over all already scattered partons from side B.
1962  for (int iScat = 0; iScat < nScatB; ++iScat) {
1963  if (PREPICKRESCATTER && iScat != iScatB) continue;
1964  int iB = scatteredB[iScat];
1965  int id2Tmp = event[iB].id();
1966 
1967  // Calculate maximum possible sHat and check whether big enough.
1968  double x2Tmp = event[iB].pNeg() / eCM;
1969  double sHatMax = beamAPtr->xMax() * x2Tmp * sCM;
1970  if (sHatMax < 4. * pT2) continue;
1971 
1972  // Select rapidity separation between two produced partons.
1973  double dyMax = log( sqrt(0.25 * sHatMax / pT2)
1974  * (1. + sqrtpos(1. - 4. * pT2 / sHatMax)) );
1975  double dy = dyMax * (2. * rndmPtr->flat() - 1.);
1976 
1977  // Reconstruct kinematical variables, especially x1.
1978  // Incoming c/b masses handled better if tau != x1 * x2.
1979  double m2Tmp = event[iB].m2();
1980  double sHatTmp = m2Tmp + 4. * pT2 * pow2(cosh(dy));
1981  double x1Tmp = (sHatTmp - m2Tmp) /(x2Tmp * sCM);
1982  double tauTmp = sHatTmp / sCM;
1983  double root = sqrtpos(1. - xT2 / tauTmp);
1984  double tHatTmp = -0.5 * sHatTmp * (1. - root);
1985  double uHatTmp = -0.5 * sHatTmp * (1. + root);
1986 
1987  // Begin evaluate parton densities at actual x1.
1988  double xPDF1[21];
1989  double xPDF1sum = 0.;
1990 
1991  // Use rescaled densities, with preweighting 9/4 for gluons.
1992  for (int id = -nQuarkIn; id <= nQuarkIn; ++id) {
1993  if (id == 0) xPDF1[10] = (9./4.) * beamAPtr->xfMPI(21, x1Tmp, pT2Fac);
1994  else xPDF1[id+10] = beamAPtr->xfMPI(id, x1Tmp, pT2Fac);
1995  xPDF1sum += xPDF1[id+10];
1996  }
1997 
1998  // Select incoming flavour according to actual PDF's.
1999  int id1Tmp = -nQuarkIn - 1;
2000  double temp = xPDF1sum * rndmPtr->flat();
2001  do { xPDF1now = xPDF1[(++id1Tmp) + 10]; temp -= xPDF1now;}
2002  while (temp > 0. && id1Tmp < nQuarkIn);
2003  if (id1Tmp == 0) id1Tmp = 21;
2004 
2005  // Assign pointers to processes relevant for incoming flavour choice:
2006  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
2007  // Factor 4./9. for incoming gluon 2 to compensate for preweighting.
2008  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
2009  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
2010  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
2011  else sigma2Tmp = &sigma2qq;
2012  double gluFac = (id1Tmp == 21) ? 4. / 9. : 1.;
2013 
2014  // Evaluate cross sections, include possibility of K factor.
2015  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
2016  // (Not necessary, but makes for better MC efficiency.)
2017  double dSigmaPartonCorr = Kfactor * gluFac
2018  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
2019  uHatTmp, alpS, alpEM);
2020 
2021  // Combine cross section, pdf's and phase space integral.
2022  double volumePhSp = 4. *dyMax;
2023  double dSigmaCorr = dSigmaPartonCorr * xPDF1sum * volumePhSp;
2024 
2025  // Dampen cross section at small pT values; part of formalism.
2026  // Use pT2 corrected for massive kinematics at this step.
2027  //?? double pT2massive = dSigmaDtSel->pT2MPI();
2028  double pT2massive = pT2;
2029  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
2030 
2031  // Compensate for prepicked rescattering if appropriate.
2032  if (PREPICKRESCATTER) dSigmaCorr *= nScatB;
2033 
2034  // Sum up total contribution for all scatterings or rescattering only.
2035  dSigmaSum += dSigmaCorr;
2036  dSigmaResc += dSigmaCorr;
2037 
2038  // Determine whether current rescattering should be selected.
2039  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
2040  i1Sel = 0;
2041  i2Sel = iB;
2042  id1Sel = id1Tmp;
2043  id2Sel = id2Tmp;
2044  x1Sel = x1Tmp;
2045  x2Sel = x2Tmp;
2046  sHatSel = sHatTmp;
2047  tHatSel = tHatTmp;
2048  uHatSel = uHatTmp;
2049  sigma2Sel = sigma2Tmp;
2050  pickOtherSel = sigma2Tmp->pickedOther();
2051  }
2052  }
2053 
2054  // Double rescatter: loop over already scattered partons from both sides.
2055  // As before, allow to use only one parton per side to speed up.
2056  if (allowDoubleRes) {
2057  for (int iScat1 = 0; iScat1 < nScatA; ++iScat1) {
2058  if (PREPICKRESCATTER && iScat1 != iScatA) continue;
2059  int iA = scatteredA[iScat1];
2060  int id1Tmp = event[iA].id();
2061  for (int iScat2 = 0; iScat2 < nScatB; ++iScat2) {
2062  if (PREPICKRESCATTER && iScat2 != iScatB) continue;
2063  int iB = scatteredB[iScat2];
2064  int id2Tmp = event[iB].id();
2065 
2066  // Calculate current sHat and check whether big enough.
2067  double sHatTmp = (event[iA].p() + event[iB].p()).m2Calc();
2068  if (sHatTmp < 4. * pT2) continue;
2069 
2070  // Reconstruct other kinematical variables. (x values not needed.)
2071  double x1Tmp = event[iA].pPos() / eCM;
2072  double x2Tmp = event[iB].pNeg() / eCM;
2073  double tauTmp = sHatTmp / sCM;
2074  double root = sqrtpos(1. - xT2 / tauTmp);
2075  double tHatTmp = -0.5 * sHatTmp * (1. - root);
2076  double uHatTmp = -0.5 * sHatTmp * (1. + root);
2077 
2078  // Assign pointers to processes relevant for incoming flavour choice:
2079  // g + g, q + g, q + qbar (same flavour), q + q(bar) (the rest).
2080  if (id1Tmp == 21 && id2Tmp == 21) sigma2Tmp = &sigma2gg;
2081  else if (id1Tmp == 21 || id2Tmp == 21) sigma2Tmp = &sigma2qg;
2082  else if (id1Tmp == -id2Tmp) sigma2Tmp = &sigma2qqbarSame;
2083  else sigma2Tmp = &sigma2qq;
2084 
2085  // Evaluate cross sections, include possibility of K factor.
2086  // Use kinematics for two symmetrical configurations (tHat <-> uHat).
2087  // (Not necessary, but makes for better MC efficiency.)
2088  double dSigmaPartonCorr = Kfactor
2089  * sigma2Tmp->sigma( id1Tmp, id2Tmp, x1Tmp, x2Tmp, sHatTmp, tHatTmp,
2090  uHatTmp, alpS, alpEM);
2091 
2092  // Combine cross section and Jacobian tHat -> pT2
2093  // (with safety minvalue).
2094  double dSigmaCorr = dSigmaPartonCorr / max(ROOTMIN, root);
2095 
2096  // Dampen cross section at small pT values; part of formalism.
2097  // Use pT2 corrected for massive kinematics at this step.
2098  //?? double pT2massive = dSigmaDtSel->pT2MPI();
2099  double pT2massive = pT2;
2100  dSigmaCorr *= pow2( pT2massive / (pT2massive + pT20) );
2101 
2102  // Compensate for prepicked rescattering if appropriate.
2103  if (PREPICKRESCATTER) dSigmaCorr *= nScatA * nScatB;
2104 
2105  // Sum up total contribution for all scatterings or rescattering only.
2106  dSigmaSum += dSigmaCorr;
2107  dSigmaResc += dSigmaCorr;
2108 
2109  // Determine whether current rescattering should be selected.
2110  if (dSigmaCorr > rndmPtr->flat() * dSigmaSum) {
2111  i1Sel = iA;
2112  i2Sel = iB;
2113  id1Sel = id1Tmp;
2114  id2Sel = id2Tmp;
2115  x1Sel = x1Tmp;
2116  x2Sel = x2Tmp;
2117  sHatSel = sHatTmp;
2118  tHatSel = tHatTmp;
2119  uHatSel = uHatTmp;
2120  sigma2Sel = sigma2Tmp;
2121  pickOtherSel = sigma2Tmp->pickedOther();
2122  }
2123  }
2124  }
2125  }
2126 
2127  // Done.
2128  return dSigmaResc;
2129 }
2130 
2131 
2132 //--------------------------------------------------------------------------
2133 
2134 // Calculate factor relating matter overlap and interaction rate,
2135 // i.e. k in <n_interaction(b)> = k * overlap(b) (neglecting
2136 // n_int = 0 corrections and energy-momentum conservation effects).
2137 
2138 void MultipartonInteractions::overlapInit() {
2139 
2140  // Initial values for iteration. Step size of b integration.
2141  nAvg = sigmaInt / sigmaND;
2142  kNow = 0.5;
2143  int stepDir = 1;
2144  double deltaB = BSTEP;
2145  if (bProfile == 2) deltaB *= min( 0.5, 2.5 * coreRadius);
2146  if (bProfile == 3) deltaB *= max(1., pow(2. / expPow, 1. / expPow));
2147 
2148  // Further variables, with dummy initial values.
2149  double nNow = 0.;
2150  double kLow = 0.;
2151  double nLow = 0.;
2152  double kHigh = 0.;
2153  double nHigh = 0.;
2154  double overlapNow = 0.;
2155  double probNow = 0.;
2156  double overlapInt = 0.5;
2157  double overlap2Int = 0.;
2158  double probInt = 0.;
2159  double probOverlapInt = 0.;
2160  double bProbInt = 0.;
2161  normPi = 1. / (2. * M_PI);
2162 
2163  // Subdivision into low-b and high-b region by interaction rate.
2164  bool pastBDiv = false;
2165  double overlapHighB = 0.;
2166 
2167  // For x-dependent matter profile, try to find a0 rather than k.
2168  // Existing framework is used, but with substitutions:
2169  // a0 tuned according to Int( Pint(b), d^2b ) = sigmaND,
2170  // nAvg = sigmaND, kNow = a0now, kLow = a0low, kHigh = a0high,
2171  // nNow = probInt, nLow = probIntLow, nHigh = probIntHigh.
2172  double rescale2 = 1.;
2173  if (bProfile == 4) {
2174  nAvg = sigmaND;
2175  kNow = XDEP_A0 / 2.0;
2176  }
2177 
2178  // First close k into an interval by binary steps,
2179  // then find k by successive interpolation.
2180  do {
2181  if (stepDir == 1) kNow *= 2.;
2182  else if (stepDir == -1) kNow *= 0.5;
2183  else kNow = kLow + (nAvg - nLow) * (kHigh - kLow) / (nHigh - nLow);
2184 
2185  // Overlap trivial if no impact parameter dependence.
2186  if (bProfile <= 0 || bProfile > 4) {
2187  probInt = 0.5 * M_PI * (1. - exp(-kNow));
2188  probOverlapInt = probInt / M_PI;
2189  bProbInt = probInt;
2190 
2191  // Ratio of b-integrated k * overlap / (1 - exp( - k * overlap)).
2192  nNow = M_PI * kNow * overlapInt / probInt;
2193 
2194  // Else integrate overlap over impact parameter.
2195  } else if (bProfile < 4) {
2196 
2197  // Reset integrals.
2198  overlapInt = (bProfile == 3) ? 0. : 0.5;
2199  overlap2Int = 0.;
2200  probInt = 0.;
2201  probOverlapInt = 0.;
2202  bProbInt = 0.;
2203  pastBDiv = false;
2204  overlapHighB = 0.;
2205 
2206  // Step in b space.
2207  double b = -0.5 * deltaB;
2208  double bArea = 0.;
2209  do {
2210  b += deltaB;
2211  bArea = 2. * M_PI * b * deltaB;
2212 
2213  // Evaluate overlap at current b value.
2214  if (bProfile == 1) {
2215  overlapNow = normPi * exp( -b*b);
2216  } else if (bProfile == 2) {
2217  overlapNow = normPi * ( fracA * exp( -min(EXPMAX, b*b))
2218  + fracB * exp( -min(EXPMAX, b*b / radius2B)) / radius2B
2219  + fracC * exp( -min(EXPMAX, b*b / radius2C)) / radius2C );
2220  } else {
2221  overlapNow = normPi * exp( -pow( b, expPow));
2222  overlapInt += bArea * overlapNow;
2223  }
2224  if (pastBDiv) overlapHighB += bArea * overlapNow;
2225 
2226  // Calculate interaction probability and integrate.
2227  probNow = 1. - exp( -min(EXPMAX, M_PI * kNow * overlapNow));
2228  overlap2Int += bArea * pow2(overlapNow);
2229  probInt += bArea * probNow;
2230  probOverlapInt += bArea * overlapNow * probNow;
2231  bProbInt += b * bArea * probNow;
2232 
2233  // Check when interaction probability has dropped sufficiently.
2234  if (!pastBDiv && probNow < PROBATLOWB) {
2235  bDiv = b + 0.5 * deltaB;
2236  pastBDiv = true;
2237  }
2238 
2239  // Continue out in b until overlap too small.
2240  } while (b < 1. || b * probNow > BMAX);
2241 
2242  // Ratio of b-integrated k * overlap / (1 - exp( - k * overlap)).
2243  nNow = M_PI * kNow * overlapInt / probInt;
2244 
2245  // x-dependent matter profile.
2246  } else if (bProfile == 4) {
2247  rescale2 = pow2(kNow / XDEP_A0);
2248  probInt = 0.;
2249  double b = 0.5 * bstepNow;
2250  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
2251  double bArea = 2. * M_PI * b * bstepNow;
2252  double pIntNow = 1 - exp( -min(EXPMAX, sigmaIntWgt[bBin] / rescale2) );
2253  probInt += bArea * rescale2 * pIntNow;
2254  b += bstepNow;
2255  }
2256  nNow = probInt;
2257  }
2258 
2259  // Replace lower or upper limit of k.
2260  if (nNow < nAvg) {
2261  kLow = kNow;
2262  nLow = nNow;
2263  if (stepDir == -1) stepDir = 0;
2264  } else {
2265  kHigh = kNow;
2266  nHigh = nNow;
2267  if (stepDir == 1) stepDir = -1;
2268  }
2269 
2270  // Continue iteration until convergence.
2271  } while (abs(nNow - nAvg) > KCONVERGE * nAvg);
2272 
2273  // Save relevant final numbers for overlap values.
2274  if (bProfile >= 0 && bProfile < 4) {
2275  double avgOverlap = probOverlapInt / probInt;
2276  zeroIntCorr = probOverlapInt / overlapInt;
2277  normOverlap = normPi * zeroIntCorr / avgOverlap;
2278  bAvg = bProbInt / probInt;
2279  enhanceBavg = (overlap2Int * probInt) / pow2(overlapInt);
2280 
2281  // Values for x-dependent matter profile.
2282  } else if (bProfile == 4) {
2283  // bAvg = Int(b * Pint(b), d2b) / sigmaND.
2284  // zeroIntCorr = Int(<n(b)> * Pint(b), d2b) / sigmaInt.
2285  bAvg = 0.;
2286  zeroIntCorr = 0.;
2287  double b = 0.5 * bstepNow;
2288  for (int bBin = 0; bBin < XDEP_BBIN; bBin++) {
2289  double bArea = 2. * M_PI * b * bstepNow;
2290  double pIntNow = 1. - exp( -min(EXPMAX, sigmaIntWgt[bBin] / rescale2) );
2291  bAvg += sqrt(rescale2) * b * bArea * rescale2 * pIntNow;
2292  zeroIntCorr += bArea * sigmaIntWgt[bBin] * pIntNow;
2293  b += bstepNow;
2294  }
2295  bAvg /= nNow;
2296  zeroIntCorr /= sigmaInt;
2297 
2298  // Other required values.
2299  a0now = kNow;
2300  infoPtr->seta0MPI(a0now * XDEP_SMB2FM);
2301  a02now = a0now * a0now;
2302  double xMin = 2. * pTmin / infoPtr->eCM();
2303  a2max = a0now * (XDEP_A1 + a1 * log(1. / xMin));
2304  a2max *= a2max;
2305  }
2306 
2307  // Relative rates for preselection of low-b and high-b region.
2308  // Other useful combinations for subsequent selection.
2309  if (bProfile > 0 && bProfile <= 3) {
2310  probLowB = M_PI * bDiv*bDiv;
2311  double probHighB = M_PI * kNow * overlapHighB;
2312  if (bProfile == 1) probHighB = M_PI * kNow * 0.5 * exp( -bDiv*bDiv);
2313  else if (bProfile == 2) {
2314  fracAhigh = fracA * exp( -bDiv*bDiv);
2315  fracBhigh = fracB * exp( -bDiv*bDiv / radius2B);
2316  fracChigh = fracC * exp( -bDiv*bDiv / radius2C);
2317  fracABChigh = fracAhigh + fracBhigh + fracChigh;
2318  probHighB = M_PI * kNow * 0.5 * fracABChigh;
2319  } else {
2320  cDiv = pow( bDiv, expPow);
2321  cMax = max(2. * expRev, cDiv);
2322  }
2323  probLowB /= (probLowB + probHighB);
2324  }
2325 
2326 }
2327 
2328 //--------------------------------------------------------------------------
2329 
2330 // Pick impact parameter and interaction rate enhancement beforehand,
2331 // i.e. before even the hardest interaction for minimum-bias events.
2332 
2333 void MultipartonInteractions::overlapFirst() {
2334 
2335  // Trivial values if no impact parameter dependence.
2336  if (bProfile <= 0 || bProfile > 4) {
2337  bNow = 1.;
2338  enhanceB = enhanceBmax = enhanceBnow = zeroIntCorr;
2339  bIsSet = true;
2340  isAtLowB = true;
2341  return;
2342  }
2343 
2344  // Possibility for user to set impact parameter. Evaluate overlap.
2345  double overlapNow = 0.;
2346  if ( userHooksPtr && userHooksPtr->canSetImpactParameter() ) {
2347  bNow = userHooksPtr->doSetImpactParameter() * bAvg;
2348  isAtLowB = ( bNow < bDiv );
2349 
2350  if (bProfile == 1) overlapNow = normPi * exp( -min(EXPMAX, bNow*bNow));
2351  else if (bProfile == 2) overlapNow = normPi *
2352  ( fracA * exp( -min(EXPMAX, bNow*bNow))
2353  + fracB * exp( -min(EXPMAX, bNow*bNow / radius2B)) / radius2B
2354  + fracC * exp( -min(EXPMAX, bNow*bNow / radius2C)) / radius2C );
2355  else overlapNow = normPi * exp( -pow( bNow, expPow));
2356  // Same enhancement for hardest process and all subsequent MPI.
2357  enhanceB = enhanceBmax = enhanceBnow = (normOverlap / normPi) * overlapNow;
2358 
2359  // Done.
2360  bNow /= bAvg;
2361  bIsSet = true;
2362  return;
2363  }
2364 
2365  // Preliminary choice between and inside low-b and high-b regions.
2366  double probAccept = 0.;
2367  do {
2368 
2369  // Treatment in low-b region: pick b flat in area.
2370  if (rndmPtr->flat() < probLowB) {
2371  isAtLowB = true;
2372  bNow = bDiv * sqrt(rndmPtr->flat());
2373 
2374  // Evaluate overlap and from that acceptance probability.
2375  if (bProfile == 1) overlapNow = normPi * exp( -bNow*bNow);
2376  else if (bProfile == 2) overlapNow = normPi *
2377  ( fracA * exp( -bNow*bNow)
2378  + fracB * exp( -bNow*bNow / radius2B) / radius2B
2379  + fracC * exp( -bNow*bNow / radius2C) / radius2C );
2380  else overlapNow = normPi * exp( -pow( bNow, expPow));
2381  probAccept = 1. - exp( -min(EXPMAX, M_PI * kNow * overlapNow));
2382 
2383  // Treatment in high-b region: pick b according to overlap.
2384  } else {
2385  isAtLowB = false;
2386 
2387  // For simple and double Gaussian pick b according to exp(-b^2 / r^2).
2388  if (bProfile == 1) {
2389  bNow = sqrt(bDiv*bDiv - log(rndmPtr->flat()));
2390  overlapNow = normPi * exp( -min(EXPMAX, bNow*bNow));
2391  } else if (bProfile == 2) {
2392  double pickFrac = rndmPtr->flat() * fracABChigh;
2393  if (pickFrac < fracAhigh)
2394  bNow = sqrt(bDiv*bDiv - log(rndmPtr->flat()));
2395  else if (pickFrac < fracAhigh + fracBhigh)
2396  bNow = sqrt(bDiv*bDiv - radius2B * log(rndmPtr->flat()));
2397  else bNow = sqrt(bDiv*bDiv - radius2C * log(rndmPtr->flat()));
2398  overlapNow = normPi * ( fracA * exp( -min(EXPMAX, bNow*bNow))
2399  + fracB * exp( -min(EXPMAX, bNow*bNow / radius2B)) / radius2B
2400  + fracC * exp( -min(EXPMAX, bNow*bNow / radius2C)) / radius2C );
2401 
2402  // For exp( - b^expPow) transform to variable c = b^expPow so that
2403  // f(b) = b * exp( - b^expPow) -> f(c) = c^r * exp(-c) with r = expRev.
2404  // case hasLowPow: expPow < 2 <=> r > 0: preselect according to
2405  // f(c) < N exp(-c/2) and then accept with N' * c^r * exp(-c/2).
2406  } else if (hasLowPow) {
2407  double cNow, acceptC;
2408  do {
2409  cNow = cDiv - 2. * log(rndmPtr->flat());
2410  acceptC = pow(cNow / cMax, expRev) * exp( -0.5 * (cNow - cMax));
2411  } while (acceptC < rndmPtr->flat());
2412  bNow = pow( cNow, 1. / expPow);
2413  overlapNow = normPi * exp( -cNow);
2414 
2415  // case !hasLowPow: expPow >= 2 <=> - 1 < r < 0: preselect according to
2416  // f(c) < N exp(-c) and then accept with N' * c^r.
2417  } else {
2418  double cNow, acceptC;
2419  do {
2420  cNow = cDiv - log(rndmPtr->flat());
2421  acceptC = pow(cNow / cDiv, expRev);
2422  } while (acceptC < rndmPtr->flat());
2423  bNow = pow( cNow, 1. / expPow);
2424  overlapNow = normPi * exp( -cNow);
2425  }
2426  double temp = M_PI * kNow *overlapNow;
2427  probAccept = (1. - exp( -min(EXPMAX, temp))) / temp;
2428  }
2429 
2430  // Confirm choice of b value. Derive enhancement factor.
2431  } while (probAccept < rndmPtr->flat());
2432 
2433  // Same enhancement for hardest process and all subsequent MPI
2434  enhanceB = enhanceBmax = enhanceBnow = (normOverlap / normPi) * overlapNow;
2435 
2436  // Done.
2437  bNow /= bAvg;
2438  bIsSet = true;
2439 
2440 }
2441 
2442 //--------------------------------------------------------------------------
2443 
2444 // Pick impact parameter and interaction rate enhancement afterwards,
2445 // i.e. after a hard interaction is known but before rest of MPI treatment.
2446 
2447 void MultipartonInteractions::overlapNext(Event& event, double pTscale,
2448  bool rehashB) {
2449 
2450  // Special case for hard diffraction if unchanged/related b.
2451  if (rehashB && bSelHard < 3) {
2452 
2453  // One option: bring b closer to its average value.
2454  bNow = infoPtr->bMPI();
2455  if (bSelHard == 2) bNow = sqrt(bNow);
2456  bNow *= bAvg;
2457  double b2 = bNow * bNow;
2458 
2459  // Caclulate new overlap enhancement factor.
2460  if (bProfile == 1) {
2461  double expb2 = exp( -min(EXPMAX, b2));
2462  enhanceB = enhanceBmax = enhanceBnow = normOverlap * expb2;
2463  } else if (bProfile == 2) {
2464  enhanceB = enhanceBmax = enhanceBnow = normOverlap *
2465  ( fracA * exp( -min(EXPMAX, b2))
2466  + fracB * exp( -min(EXPMAX, b2 / radius2B)) / radius2B
2467  + fracC * exp( -min(EXPMAX, b2 / radius2C)) / radius2C );
2468  } else {
2469  double cNow = pow( bNow, expPow);
2470  enhanceB = enhanceBmax = enhanceBnow = normOverlap * exp(-cNow);
2471  }
2472 
2473  // Done for simple cases.
2474  bNow /= bAvg;
2475  bIsSet = true;
2476  return;
2477  }
2478 
2479  // Default, valid for bProfile = 0. Also initial Sudakov.
2480  enhanceB = enhanceBmax = enhanceBnow = zeroIntCorr;
2481  if (bProfile <= 0 || bProfile > 4) return;
2482 
2483  // Alternative choices of event scale for Sudakov in (pT, b) space.
2484  if (bSelScale == 1) {
2485  vector<double> mmT;
2486  for (int i = 5; i < event.size(); ++i) if (event[i].isFinal()) {
2487  mmT.push_back( event[i].m() + event[i].mT() );
2488  for (int j = int(mmT.size()) - 1; j > 0; --j)
2489  if (mmT[j] > mmT[j - 1]) swap( mmT[j], mmT[j - 1] );
2490  }
2491  pTscale = 0.5 * mmT[0];
2492  for (int j = 1; j < int(mmT.size()); ++j) pTscale += mmT[j] / (j + 1.);
2493  } else if (bSelScale == 2) pTscale = event.scale();
2494  double pT2scale = pTscale*pTscale;
2495 
2496  // Use trial interaction for x-dependent matter profile.
2497  if (bProfile == 4) {
2498  double pTtrial = 0.;
2499  do {
2500  // Pick b according to wanted O(b, x1, x2).
2501  double expb2 = rndmPtr->flat();
2502  double w1 = XDEP_A1 + a1 * log(1. / infoPtr->x1());
2503  double w2 = XDEP_A1 + a1 * log(1. / infoPtr->x2());
2504  double fac = a02now * (w1 * w1 + w2 * w2);
2505  b2now = - fac * log(expb2);
2506  bNow = sqrt(b2now);
2507 
2508  // Enhancement factor for the hard process and overestimate
2509  // for fastPT2. Note that existing framework has a (1. / sigmaND)
2510  // present.
2511  enhanceB = sigmaND / M_PI / fac * expb2;
2512  enhanceBmax = sigmaND / 2. / M_PI / a02now
2513  * exp( -b2now / 2. / a2max );
2514 
2515  // Trial interaction. Keep going until pTtrial < pTscale.
2516  pTtrial = pTnext(pTmax, pTmin, event);
2517  } while (pTtrial > pTscale);
2518  bNow /= bAvg;
2519  bIsSet = true;
2520  return;
2521  }
2522 
2523  // Begin loop over pT-dependent rejection of b value.
2524  do {
2525 
2526  // Flat enhancement distribution for simple Gaussian.
2527  if (bProfile == 1) {
2528  double expb2 = rndmPtr->flat();
2529  // Same enhancement for hardest process and all subsequent MPI.
2530  enhanceB = enhanceBmax = enhanceBnow = normOverlap * expb2;
2531  bNow = sqrt( -log(expb2));
2532 
2533  // For double Gaussian go the way via b, according to each Gaussian.
2534  } else if (bProfile == 2) {
2535  double bType = rndmPtr->flat();
2536  double b2 = -log( rndmPtr->flat() );
2537  if (bType < fracA) ;
2538  else if (bType < fracA + fracB) b2 *= radius2B;
2539  else b2 *= radius2C;
2540  // Same enhancement for hardest process and all subsequent MPI.
2541  enhanceB = enhanceBmax = enhanceBnow = normOverlap *
2542  ( fracA * exp( -min(EXPMAX, b2))
2543  + fracB * exp( -min(EXPMAX, b2 / radius2B)) / radius2B
2544  + fracC * exp( -min(EXPMAX, b2 / radius2C)) / radius2C );
2545  bNow = sqrt(b2);
2546 
2547  // For exp( - b^expPow) transform to variable c = b^expPow so that
2548  // f(b) = b * exp( - b^expPow) -> f(c) = c^r * exp(-c) with r = expRev.
2549  // case hasLowPow: expPow < 2 <=> r > 0:
2550  // f(c) < r^r exp(-r) for c < 2r, < (2r)^r exp(-r-c/2) for c > 2r.
2551  } else if (bProfile == 3 && hasLowPow) {
2552  double cNow, acceptC;
2553  double probLowC = expRev / (expRev + pow(2., expRev) * exp( - expRev));
2554  do {
2555  if (rndmPtr->flat() < probLowC) {
2556  cNow = 2. * expRev * rndmPtr->flat();
2557  acceptC = pow( cNow / expRev, expRev) * exp(expRev - cNow);
2558  } else {
2559  cNow = 2. * (expRev - log( rndmPtr->flat() ));
2560  acceptC = pow(0.5 * cNow / expRev, expRev)
2561  * exp(expRev - 0.5 * cNow);
2562  }
2563  } while (acceptC < rndmPtr->flat());
2564  // Same enhancement for hardest process and all subsequent MPI.
2565  enhanceB = enhanceBmax = enhanceBnow = normOverlap *exp(-cNow);
2566  bNow = pow( cNow, 1. / expPow);
2567 
2568  // case !hasLowPow: expPow >= 2 <=> - 1 < r < 0:
2569  // f(c) < c^r for c < 1, f(c) < exp(-c) for c > 1.
2570  } else if (bProfile == 3 && !hasLowPow) {
2571  double cNow, acceptC;
2572  double probLowC = expPow / (2. * exp(-1.) + expPow);
2573  do {
2574  if (rndmPtr->flat() < probLowC) {
2575  cNow = pow( rndmPtr->flat(), 0.5 * expPow);
2576  acceptC = exp(-cNow);
2577  } else {
2578  cNow = 1. - log( rndmPtr->flat() );
2579  acceptC = pow( cNow, expRev);
2580  }
2581  } while (acceptC < rndmPtr->flat());
2582  // Same enhancement for hardest process and all subsequent MPI.
2583  enhanceB = enhanceBmax = enhanceBnow = normOverlap * exp(-cNow);
2584  bNow = pow( cNow, 1. / expPow);
2585  }
2586 
2587  // Evaluate "Sudakov form factor" for not having a harder interaction.
2588  } while (sudakov(pT2scale, enhanceB) < rndmPtr->flat());
2589 
2590  // Done.
2591  bNow /= bAvg;
2592  bIsSet = true;
2593 }
2594 
2595 //--------------------------------------------------------------------------
2596 
2597 // Print statistics on number of multiparton-interactions processes.
2598 
2599 void MultipartonInteractions::statistics(bool resetStat) {
2600 
2601  // Header.
2602  cout << "\n *------- PYTHIA Multiparton Interactions Statistics -----"
2603  << "---*\n"
2604  << " | "
2605  << " |\n"
2606  << " | Note: excludes hardest subprocess if already listed above "
2607  << " |\n"
2608  << " | "
2609  << " |\n"
2610  << " | Subprocess Code | Times"
2611  << " |\n"
2612  << " | | "
2613  << " |\n"
2614  << " |------------------------------------------------------------"
2615  << "-|\n"
2616  << " | | "
2617  << " |\n";
2618 
2619  // Loop over existing processes. Sum of all subprocesses.
2620  int numberSum = 0;
2621  for ( map<int, int>::iterator iter = nGen.begin(); iter != nGen.end();
2622  ++iter) {
2623  int code = iter -> first;
2624  int number = iter->second;
2625  numberSum += number;
2626 
2627  // Find process name that matches code.
2628  string name = " ";
2629  bool foundName = false;
2630  SigmaMultiparton* dSigma;
2631  for (int i = 0; i < 4; ++i) {
2632  if (i == 0) dSigma = &sigma2gg;
2633  else if (i == 1) dSigma = &sigma2qg;
2634  else if (i == 2) dSigma = &sigma2qqbarSame;
2635  else dSigma = &sigma2qq;
2636  int nProc = dSigma->nProc();
2637  for (int iProc = 0; iProc < nProc; ++iProc)
2638  if (dSigma->codeProc(iProc) == code) {
2639  name = dSigma->nameProc(iProc);
2640  foundName = true;
2641  }
2642  if (foundName) break;
2643  }
2644 
2645  // Print individual process info.
2646  cout << " | " << left << setw(40) << name << right << setw(5) << code
2647  << " | " << setw(11) << number << " |\n";
2648  }
2649 
2650  // Print summed process info.
2651  cout << " | "
2652  << " |\n"
2653  << " | " << left << setw(45) << "sum" << right << " | " << setw(11)
2654  << numberSum << " |\n";
2655 
2656  // Listing finished.
2657  cout << " | | "
2658  << " |\n"
2659  << " *------- End PYTHIA Multiparton Interactions Statistics ----"
2660  << "-*" << endl;
2661 
2662  // Optionally reset statistics contents.
2663  if (resetStat) resetStatistics();
2664 
2665 }
2666 
2667 //==========================================================================
2668 
2669 } // end namespace Pythia8
Definition: AgUStep.h:26