StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SigmaProcess.cc
1 // SigmaProcess.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // SigmaProcess class, and classes derived from it.
8 
9 #include "SigmaProcess.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The SigmaProcess class.
16 // Base class for cross sections.
17 
18 //--------------------------------------------------------------------------
19 
20 // Constants: could be changed here if desired, but normally should not.
21 // These are of technical nature, as described for each.
22 
23 // Conversion of GeV^{-2} to mb for cross section.
24 const double SigmaProcess::CONVERT2MB = 0.389380;
25 
26 // The sum of outgoing masses must not be too close to the cm energy.
27 const double SigmaProcess::MASSMARGIN = 0.1;
28 
29 // Parameters of momentum rescaling procedure: maximally allowed
30 // relative energy error and number of iterations.
31 const double SigmaProcess::COMPRELERR = 1e-10;
32 const int SigmaProcess::NCOMPSTEP = 10;
33 
34 //--------------------------------------------------------------------------
35 
36 // Perform simple initialization and store pointers.
37 
38 void SigmaProcess::init(Info* infoPtrIn, Settings* settingsPtrIn,
39  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn, BeamParticle* beamAPtrIn,
40  BeamParticle* beamBPtrIn, Couplings* couplingsPtrIn,
41  SigmaTotal* sigmaTotPtrIn, SusyLesHouches* slhaPtrIn) {
42 
43  // Store pointers.
44  infoPtr = infoPtrIn;
45  settingsPtr = settingsPtrIn;
46  particleDataPtr = particleDataPtrIn;
47  rndmPtr = rndmPtrIn;
48  beamAPtr = beamAPtrIn;
49  beamBPtr = beamBPtrIn;
50  couplingsPtr = couplingsPtrIn;
51  sigmaTotPtr = sigmaTotPtrIn;
52  slhaPtr = slhaPtrIn;
53 
54  // Read out some properties of beams to allow shorthand.
55  idA = (beamAPtr > 0) ? beamAPtr->id() : 0;
56  idB = (beamBPtr > 0) ? beamBPtr->id() : 0;
57  mA = (beamAPtr > 0) ? beamAPtr->m() : 0.;
58  mB = (beamBPtr > 0) ? beamBPtr->m() : 0.;
59  isLeptonA = (beamAPtr > 0) ? beamAPtr->isLepton() : false;
60  isLeptonB = (beamBPtr > 0) ? beamBPtr->isLepton() : false;
61  hasLeptonBeams = isLeptonA || isLeptonB;
62 
63  // K factor, multiplying resolved processes. (But not here for MPI.)
64  Kfactor = settingsPtr->parm("SigmaProcess:Kfactor");
65 
66  // Maximum incoming quark flavour.
67  nQuarkIn = settingsPtr->mode("PDFinProcess:nQuarkIn");
68 
69  // Medium heavy fermion masses set massless or not in ME expressions.
70  mcME = (settingsPtr->flag("SigmaProcess:cMassiveME"))
71  ? particleDataPtr->m0(4) : 0.;
72  mbME = (settingsPtr->flag("SigmaProcess:bMassiveME"))
73  ? particleDataPtr->m0(5) : 0.;
74  mmuME = (settingsPtr->flag("SigmaProcess:muMassiveME"))
75  ? particleDataPtr->m0(13) : 0.;
76  mtauME = (settingsPtr->flag("SigmaProcess:tauMassiveME"))
77  ? particleDataPtr->m0(15) : 0.;
78 
79  // Renormalization scale choice.
80  renormScale1 = settingsPtr->mode("SigmaProcess:renormScale1");
81  renormScale2 = settingsPtr->mode("SigmaProcess:renormScale2");
82  renormScale3 = settingsPtr->mode("SigmaProcess:renormScale3");
83  renormScale3VV = settingsPtr->mode("SigmaProcess:renormScale3VV");
84  renormMultFac = settingsPtr->parm("SigmaProcess:renormMultFac");
85  renormFixScale = settingsPtr->parm("SigmaProcess:renormFixScale");
86 
87  // Factorization scale choice.
88  factorScale1 = settingsPtr->mode("SigmaProcess:factorScale1");
89  factorScale2 = settingsPtr->mode("SigmaProcess:factorScale2");
90  factorScale3 = settingsPtr->mode("SigmaProcess:factorScale3");
91  factorScale3VV = settingsPtr->mode("SigmaProcess:factorScale3VV");
92  factorMultFac = settingsPtr->parm("SigmaProcess:factorMultFac");
93  factorFixScale = settingsPtr->parm("SigmaProcess:factorFixScale");
94 
95  // CP violation parameters for the BSM Higgs sector.
96  higgsH1parity = settingsPtr->mode("HiggsH1:parity");
97  higgsH1eta = settingsPtr->parm("HiggsH1:etaParity");
98  higgsH2parity = settingsPtr->mode("HiggsH2:parity");
99  higgsH2eta = settingsPtr->parm("HiggsH2:etaParity");
100  higgsA3parity = settingsPtr->mode("HiggsA3:parity");
101  higgsA3eta = settingsPtr->parm("HiggsA3:etaParity");
102 
103  // If BSM not switched on then H1 should have SM properties.
104  if (!settingsPtr->flag("Higgs:useBSM")){
105  higgsH1parity = 1;
106  higgsH1eta = 0.;
107  }
108 
109 }
110 
111 //--------------------------------------------------------------------------
112 
113 // Set up allowed flux of incoming partons.
114 // addBeam: set up PDF's that need to be evaluated for the two beams.
115 // addPair: set up pairs of incoming partons from the two beams.
116 
117 bool SigmaProcess::initFlux() {
118 
119  // Reset arrays (in case of several init's in same run).
120  inBeamA.clear();
121  inBeamB.clear();
122  inPair.clear();
123 
124  // Read in process-specific channel information.
125  string fluxType = inFlux();
126 
127  // Case with g g incoming state.
128  if (fluxType == "gg") {
129  addBeamA(21);
130  addBeamB(21);
131  addPair(21, 21);
132  }
133 
134  // Case with q g incoming state.
135  else if (fluxType == "qg") {
136  for (int i = -nQuarkIn; i <= nQuarkIn; ++i) {
137  int idNow = (i == 0) ? 21 : i;
138  addBeamA(idNow);
139  addBeamB(idNow);
140  }
141  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
142  if (idNow != 0) {
143  addPair(idNow, 21);
144  addPair(21, idNow);
145  }
146  }
147 
148  // Case with q q', q qbar' or qbar qbar' incoming state.
149  else if (fluxType == "qq") {
150  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
151  if (idNow != 0) {
152  addBeamA(idNow);
153  addBeamB(idNow);
154  }
155  for (int id1Now = -nQuarkIn; id1Now <= nQuarkIn; ++id1Now)
156  if (id1Now != 0)
157  for (int id2Now = -nQuarkIn; id2Now <= nQuarkIn; ++id2Now)
158  if (id2Now != 0)
159  addPair(id1Now, id2Now);
160  }
161 
162  // Case with q qbar incoming state.
163  else if (fluxType == "qqbarSame") {
164  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
165  if (idNow != 0) {
166  addBeamA(idNow);
167  addBeamB(idNow);
168  }
169  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
170  if (idNow != 0)
171  addPair(idNow, -idNow);
172  }
173 
174  // Case with f f', f fbar', fbar fbar' incoming state.
175  else if (fluxType == "ff") {
176  // If beams are leptons then they are also the colliding partons.
177  if ( isLeptonA && isLeptonB ) {
178  addBeamA(idA);
179  addBeamB(idB);
180  addPair(idA, idB);
181  // First beam is lepton and second is hadron.
182  } else if ( isLeptonA ) {
183  addBeamA(idA);
184  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
185  if (idNow != 0) {
186  addBeamB(idNow);
187  addPair(idA, idNow);
188  }
189  // First beam is hadron and second is lepton.
190  } else if ( isLeptonB ) {
191  addBeamB(idB);
192  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
193  if (idNow != 0) {
194  addBeamA(idNow);
195  addPair(idNow, idB);
196  }
197  // Hadron beams gives quarks.
198  } else {
199  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
200  if (idNow != 0) {
201  addBeamA(idNow);
202  addBeamB(idNow);
203  }
204  for (int id1Now = -nQuarkIn; id1Now <= nQuarkIn; ++id1Now)
205  if (id1Now != 0)
206  for (int id2Now = -nQuarkIn; id2Now <= nQuarkIn; ++id2Now)
207  if (id2Now != 0)
208  addPair(id1Now, id2Now);
209  }
210  }
211 
212  // Case with f fbar incoming state.
213  else if (fluxType == "ffbarSame") {
214  // If beams are antiparticle pair and leptons then also colliding partons.
215  if ( idA + idB == 0 && isLeptonA ) {
216  addBeamA(idA);
217  addBeamB(idB);
218  addPair(idA, idB);
219  // Else assume both to be hadrons, for better or worse.
220  } else {
221  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
222  if (idNow != 0) {
223  addBeamA(idNow);
224  addBeamB(idNow);
225  }
226  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
227  if (idNow != 0)
228  addPair(idNow, -idNow);
229  }
230  }
231 
232  // Case with f fbar' charged(+-1) incoming state.
233  else if (fluxType == "ffbarChg") {
234  // If beams are leptons then also colliding partons.
235  if ( isLeptonA && isLeptonB && abs( particleDataPtr->chargeType(idA)
236  + particleDataPtr->chargeType(idB) ) == 3 ) {
237  addBeamA(idA);
238  addBeamB(idB);
239  addPair(idA, idB);
240  // Hadron beams gives quarks.
241  } else {
242  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
243  if (idNow != 0) {
244  addBeamA(idNow);
245  addBeamB(idNow);
246  }
247  for (int id1Now = -nQuarkIn; id1Now <= nQuarkIn; ++id1Now)
248  if (id1Now != 0)
249  for (int id2Now = -nQuarkIn; id2Now <= nQuarkIn; ++id2Now)
250  if (id2Now != 0 && id1Now * id2Now < 0
251  && (abs(id1Now) + abs(id2Now))%2 == 1) addPair(id1Now, id2Now);
252  }
253  }
254 
255  // Case with f fbar' generic incoming state.
256  else if (fluxType == "ffbar") {
257  // If beams are leptons then also colliding partons.
258  if (isLeptonA && isLeptonB && idA * idB < 0) {
259  addBeamA(idA);
260  addBeamB(idB);
261  addPair(idA, idB);
262  // Hadron beams gives quarks.
263  } else {
264  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
265  if (idNow != 0) {
266  addBeamA(idNow);
267  addBeamB(idNow);
268  }
269  for (int id1Now = -nQuarkIn; id1Now <= nQuarkIn; ++id1Now)
270  if (id1Now != 0)
271  for (int id2Now = -nQuarkIn; id2Now <= nQuarkIn; ++id2Now)
272  if (id2Now != 0 && id1Now * id2Now < 0)
273  addPair(id1Now, id2Now);
274  }
275  }
276 
277  // Case with f gamma incoming state.
278  else if (fluxType == "fgm") {
279  // Fermion from incoming side A.
280  if ( isLeptonA ) {
281  addBeamA(idA);
282  addPair(idA, 22);
283  } else {
284  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
285  if (idNow != 0) {
286  addBeamA(idNow);
287  addPair(idNow, 22);
288  }
289  }
290  // Fermion from incoming side B.
291  if ( isLeptonB ) {
292  addBeamB( idB);
293  addPair(22, idB);
294  } else {
295  for (int idNow = -nQuarkIn; idNow <= nQuarkIn; ++idNow)
296  if (idNow != 0) {
297  addBeamB(idNow);
298  addPair(22, idNow);
299  }
300  }
301  // Photons in the beams.
302  addBeamA(22);
303  addBeamB(22);
304  }
305 
306  // Case with gamma gamma incoming state.
307  else if (fluxType == "ggm") {
308  addBeamA(21);
309  addBeamA(22);
310  addBeamB(21);
311  addBeamB(22);
312  addPair(21, 22);
313  addPair(22, 21);
314  }
315 
316  // Case with gamma gamma incoming state.
317  else if (fluxType == "gmgm") {
318  addBeamA(22);
319  addBeamB(22);
320  addPair(22, 22);
321  }
322 
323  // Unrecognized fluxType is bad sign. Else done.
324  else {
325  infoPtr->errorMsg("Error in SigmaProcess::initFlux: "
326  "unrecognized inFlux type", fluxType);
327  return false;
328  }
329  return true;
330 
331 }
332 
333 //--------------------------------------------------------------------------
334 
335 // Convolute matrix-element expression(s) with parton flux and K factor.
336 
337 double SigmaProcess::sigmaPDF() {
338 
339  // Evaluate and store the required parton densities.
340  for (int j = 0; j < sizeBeamA(); ++j)
341  inBeamA[j].pdf = beamAPtr->xfHard( inBeamA[j].id, x1Save, Q2FacSave);
342  for (int j = 0; j < sizeBeamB(); ++j)
343  inBeamB[j].pdf = beamBPtr->xfHard( inBeamB[j].id, x2Save, Q2FacSave);
344 
345  // Loop over allowed incoming channels.
346  sigmaSumSave = 0.;
347  for (int i = 0; i < sizePair(); ++i) {
348 
349  // Evaluate hard-scattering cross section. Include K factor.
350  inPair[i].pdfSigma = Kfactor
351  * sigmaHatWrap(inPair[i].idA, inPair[i].idB);
352 
353  // Multiply by respective parton densities.
354  for (int j = 0; j < sizeBeamA(); ++j)
355  if (inPair[i].idA == inBeamA[j].id) {
356  inPair[i].pdfA = inBeamA[j].pdf;
357  inPair[i].pdfSigma *= inBeamA[j].pdf;
358  break;
359  }
360  for (int j = 0; j < sizeBeamB(); ++j)
361  if (inPair[i].idB == inBeamB[j].id) {
362  inPair[i].pdfB = inBeamB[j].pdf;
363  inPair[i].pdfSigma *= inBeamB[j].pdf;
364  break;
365  }
366 
367  // Sum for all channels.
368  sigmaSumSave += inPair[i].pdfSigma;
369  }
370 
371  // Done.
372  return sigmaSumSave;
373 
374 }
375 
376 //--------------------------------------------------------------------------
377 
378 // Select incoming parton channel and extract parton densities (resolved).
379 
380 void SigmaProcess::pickInState(int id1in, int id2in) {
381 
382  // Multiparton interactions: partons already selected.
383  if (id1in != 0 && id2in != 0) {
384  id1 = id1in;
385  id2 = id2in;
386  }
387 
388  // Pick channel. Extract channel flavours and pdf's.
389  double sigmaRand = sigmaSumSave * rndmPtr->flat();
390  for (int i = 0; i < sizePair(); ++i) {
391  sigmaRand -= inPair[i].pdfSigma;
392  if (sigmaRand <= 0.) {
393  id1 = inPair[i].idA;
394  id2 = inPair[i].idB;
395  pdf1Save = inPair[i].pdfA;
396  pdf2Save = inPair[i].pdfB;
397  break;
398  }
399  }
400 
401 }
402 
403 //--------------------------------------------------------------------------
404 
405 // Calculate incoming modified masses and four-vectors for matrix elements.
406 
407 bool SigmaProcess::setupForMEin() {
408 
409  // Initially assume it will work out to set up modified kinematics.
410  bool allowME = true;
411 
412  // Correct incoming c, b, mu and tau to be massive or not.
413  mME[0] = 0.;
414  int id1Tmp = abs(id1);
415  if (id1Tmp == 4) mME[0] = mcME;
416  if (id1Tmp == 5) mME[0] = mbME;
417  if (id1Tmp == 13) mME[0] = mmuME;
418  if (id1Tmp == 15) mME[0] = mtauME;
419  mME[1] = 0.;
420  int id2Tmp = abs(id2);
421  if (id2Tmp == 4) mME[1] = mcME;
422  if (id2Tmp == 5) mME[1] = mbME;
423  if (id2Tmp == 13) mME[1] = mmuME;
424  if (id2Tmp == 15) mME[1] = mtauME;
425 
426  // If kinematically impossible return to massless case, but set error.
427  if (mME[0] + mME[1] >= mH) {
428  mME[0] = 0.;
429  mME[1] = 0.;
430  allowME = false;
431  }
432 
433  // Do incoming two-body kinematics for massless or massive cases.
434  if (mME[0] == 0. && mME[1] == 0.) {
435  pME[0] = 0.5 * mH * Vec4( 0., 0., 1., 1.);
436  pME[1] = 0.5 * mH * Vec4( 0., 0., -1., 1.);
437  } else {
438  double e0 = 0.5 * (mH * mH + mME[0] * mME[0] - mME[1] * mME[1]) / mH;
439  double pz0 = sqrtpos(e0 * e0 - mME[0] * mME[0]);
440  pME[0] = Vec4( 0., 0., pz0, e0);
441  pME[1] = Vec4( 0., 0., -pz0, mH - e0);
442  }
443 
444  // Done.
445  return allowME;
446 
447 }
448 
449 //--------------------------------------------------------------------------
450 
451 // Evaluate weight for W decay distribution in t -> W b -> f fbar b.
452 
453 double SigmaProcess::weightTopDecay( Event& process, int iResBeg,
454  int iResEnd) {
455 
456  // If not pair W d/s/b and mother t then return unit weight.
457  if (iResEnd - iResBeg != 1) return 1.;
458  int iW1 = iResBeg;
459  int iB2 = iResBeg + 1;
460  int idW1 = process[iW1].idAbs();
461  int idB2 = process[iB2].idAbs();
462  if (idW1 != 24) {
463  swap(iW1, iB2);
464  swap(idW1, idB2);
465  }
466  if (idW1 != 24 || (idB2 != 1 && idB2 != 3 && idB2 != 5)) return 1.;
467  int iT = process[iW1].mother1();
468  if (iT <= 0 || process[iT].idAbs() != 6) return 1.;
469 
470  // Find sign-matched order of W decay products.
471  int iF = process[iW1].daughter1();
472  int iFbar = process[iW1].daughter2();
473  if (iFbar - iF != 1) return 1.;
474  if (process[iT].id() * process[iF].id() < 0) swap(iF, iFbar);
475 
476  // Weight and maximum weight.
477  double wt = (process[iT].p() * process[iFbar].p())
478  * (process[iF].p() * process[iB2].p());
479  double wtMax = ( pow4(process[iT].m()) - pow4(process[iW1].m()) ) / 8.;
480 
481  // Done.
482  return wt / wtMax;
483 
484 }
485 
486 
487 //--------------------------------------------------------------------------
488 
489 // Evaluate weight for Z0/W+- decay distributions in H -> Z0/W+ Z0/W- -> 4f.
490 
491 double SigmaProcess::weightHiggsDecay( Event& process, int iResBeg,
492  int iResEnd) {
493 
494  // If not pair Z0 Z0 or W+ W- then return unit weight.
495  if (iResEnd - iResBeg != 1) return 1.;
496  int iZW1 = iResBeg;
497  int iZW2 = iResBeg + 1;
498  int idZW1 = process[iZW1].id();
499  int idZW2 = process[iZW2].id();
500  if (idZW1 < 0) {
501  swap(iZW1, iZW2);
502  swap(idZW1, idZW2);
503  }
504  if ( (idZW1 != 23 || idZW2 != 23) && (idZW1 != 24 || idZW2 != -24) )
505  return 1.;
506 
507  // If mother is not Higgs then return unit weight.
508  int iH = process[iZW1].mother1();
509  if (iH <= 0) return 1.;
510  int idH = process[iH].id();
511  if (idH != 25 && idH != 35 && idH !=36) return 1.;
512 
513  // Parameters depend on Higgs type: H0(H_1), H^0(H_2) or A^0(H_3).
514  int higgsParity = higgsH1parity;
515  double higgsEta = higgsH1eta;
516  if (idH == 35) {
517  higgsParity = higgsH2parity;
518  higgsEta = higgsH2eta;
519  } else if (idH == 36) {
520  higgsParity = higgsA3parity;
521  higgsEta = higgsA3eta;
522  }
523 
524  // Option with isotropic decays.
525  if (higgsParity == 0) return 1.;
526 
527  // Maximum and initial weight.
528  double wtMax = pow4(process[iH].m());
529  double wt = wtMax;
530 
531  // Find sign-matched order of Z0/W+- decay products.
532  int i3 = process[iZW1].daughter1();
533  int i4 = process[iZW1].daughter2();
534  if (process[i3].id() < 0) swap( i3, i4);
535  int i5 = process[iZW2].daughter1();
536  int i6 = process[iZW2].daughter2();
537  if (process[i5].id() < 0) swap( i5, i6);
538 
539  // Evaluate four-vector products and find masses..
540  double p35 = 2. * process[i3].p() * process[i5].p();
541  double p36 = 2. * process[i3].p() * process[i6].p();
542  double p45 = 2. * process[i4].p() * process[i5].p();
543  double p46 = 2. * process[i4].p() * process[i6].p();
544  double p34 = 2. * process[i3].p() * process[i4].p();
545  double p56 = 2. * process[i5].p() * process[i6].p();
546  double mZW1 = process[iZW1].m();
547  double mZW2 = process[iZW2].m();
548 
549  // For mixed CP states need epsilon product and gauge boson masses.
550  double epsilonProd = 0.;
551  if (higgsParity == 3) {
552  double p[4][4];
553  for (int i = 0; i < 4; ++i) {
554  int ii = i3;
555  if (i == 1) ii = i4;
556  if (i == 2) ii = i5;
557  if (i == 3) ii = i6;
558  p[i][0] = process[ii].e();
559  p[i][1] = process[ii].px();
560  p[i][2] = process[ii].py();
561  p[i][3] = process[ii].pz();
562  }
563  epsilonProd
564  = p[0][0]*p[1][1]*p[2][2]*p[3][3] - p[0][0]*p[1][1]*p[2][3]*p[3][2]
565  - p[0][0]*p[1][2]*p[2][1]*p[3][3] + p[0][0]*p[1][2]*p[2][3]*p[3][1]
566  + p[0][0]*p[1][3]*p[2][1]*p[3][2] - p[0][0]*p[1][3]*p[2][2]*p[3][1]
567  - p[0][1]*p[1][0]*p[2][2]*p[3][3] + p[0][1]*p[1][0]*p[2][3]*p[3][2]
568  + p[0][1]*p[1][2]*p[2][0]*p[3][3] - p[0][1]*p[1][2]*p[2][3]*p[3][0]
569  - p[0][1]*p[1][3]*p[2][0]*p[3][2] + p[0][1]*p[1][3]*p[2][2]*p[3][0]
570  + p[0][2]*p[1][0]*p[2][1]*p[3][3] - p[0][2]*p[1][0]*p[2][3]*p[3][1]
571  - p[0][2]*p[1][1]*p[2][0]*p[3][3] + p[0][2]*p[1][1]*p[2][3]*p[3][0]
572  + p[0][2]*p[1][3]*p[2][0]*p[3][1] - p[0][2]*p[1][3]*p[2][1]*p[3][0]
573  - p[0][3]*p[1][0]*p[2][1]*p[3][2] + p[0][3]*p[1][0]*p[2][2]*p[3][1]
574  + p[0][3]*p[1][1]*p[2][0]*p[3][2] - p[0][3]*p[1][1]*p[2][2]*p[3][0]
575  - p[0][3]*p[1][2]*p[2][0]*p[3][1] + p[0][3]*p[1][2]*p[2][1]*p[3][0];
576  }
577 
578  // Z0 Z0 decay: vector and axial couplings of two fermion pairs.
579  if (idZW1 == 23) {
580  double vf1 = couplingsPtr->vf(process[i3].idAbs());
581  double af1 = couplingsPtr->af(process[i3].idAbs());
582  double vf2 = couplingsPtr->vf(process[i5].idAbs());
583  double af2 = couplingsPtr->af(process[i5].idAbs());
584  double va12asym = 4. * vf1 * af1 * vf2 * af2
585  / ( (vf1*vf1 + af1*af1) * (vf2*vf2 + af2*af2) );
586  double etaMod = higgsEta / pow2( particleDataPtr->m0(23) );
587 
588  // Normal CP-even decay.
589  if (higgsParity == 1) wt = 8. * (1. + va12asym) * p35 * p46
590  + 8. * (1. - va12asym) * p36 * p45;
591 
592  // CP-odd decay (normal for A0(H_3)).
593  else if (higgsParity == 2) wt = ( pow2(p35 + p46)
594  + pow2(p36 + p45) - 2. * p34 * p56
595  - 2. * pow2(p35 * p46 - p36 * p45) / (p34 * p56)
596  + va12asym * (p35 + p36 - p45 - p46) * (p35 + p45 - p36 - p46) )
597  / (1. + va12asym);
598 
599  // Mixed CP states.
600  else wt = 32. * ( 0.25 * ( (1. + va12asym) * p35 * p46
601  + (1. - va12asym) * p36 * p45 ) - 0.5 * etaMod * epsilonProd
602  * ( (1. + va12asym) * (p35 + p46) - (1. - va12asym) * (p36 + p45) )
603  + 0.0625 * etaMod * etaMod * (-2. * pow2(p34 * p56)
604  - 2. * pow2(p35 * p46 - p36 * p45)
605  + p34 * p56 * (pow2(p35 + p46) + pow2(p36 + p45))
606  + va12asym * p34 * p56 * (p35 + p36 - p45 - p46)
607  * (p35 + p45 - p36 - p46) ) ) / ( 1. + 2. * etaMod * mZW1 * mZW2
608  + 2. * pow2(etaMod * mZW1 * mZW2) * (1. + va12asym) );
609 
610  // W+ W- decay.
611  } else if (idZW1 == 24) {
612  double etaMod = higgsEta / pow2( particleDataPtr->m0(24) );
613 
614  // Normal CP-even decay.
615  if (higgsParity == 1) wt = 16. * p35 * p46;
616 
617  // CP-odd decay (normal for A0(H_3)).
618  else if (higgsParity == 2) wt = 0.5 * ( pow2(p35 + p46)
619  + pow2(p36 + p45) - 2. * p34 * p56
620  - 2. * pow2(p35 * p46 - p36 * p45) / (p34 * p56)
621  + (p35 + p36 - p45 - p46) * (p35 + p45 - p36 - p46) );
622 
623  // Mixed CP states.
624  else wt = 32. * ( 0.25 * 2. * p35 * p46
625  - 0.5 * etaMod * epsilonProd * 2. * (p35 + p46)
626  + 0.0625 * etaMod * etaMod * (-2. * pow2(p34 * p56)
627  - 2. * pow2(p35 * p46 - p36 * p45)
628  + p34 * p56 * (pow2(p35 + p46) + pow2(p36 + p45))
629  + p34 * p56 * (p35 + p36 - p45 - p46) * (p35 + p45 - p36 - p46) ) )
630  / ( 1. * 2. * etaMod * mZW1 * mZW2 + 2. * pow2(etaMod * mZW1 * mZW2) );
631  }
632 
633  // Done.
634  return wt / wtMax;
635 
636 }
637 
638 //==========================================================================
639 
640 // The Sigma1Process class.
641 // Base class for resolved 2 -> 1 cross sections; derived from SigmaProcess.
642 
643 //--------------------------------------------------------------------------
644 
645 // Wrapper to sigmaHat, to (a) store current incoming flavours,
646 // (b) convert from GeV^-2 to mb where required, and
647 // (c) convert from |M|^2 to d(sigmaHat)/d(tHat) where required.
648 
649 double Sigma1Process::sigmaHatWrap(int id1in, int id2in) {
650 
651  id1 = id1in;
652  id2 = id2in;
653  double sigmaTmp = sigmaHat();
654  if (convertM2()) {
655  sigmaTmp /= 2. * sH;
656  // Convert 2 * pi * delta(p^2 - m^2) to Breit-Wigner with same area.
657  int idTmp = resonanceA();
658  double mTmp = particleDataPtr->m0(idTmp);
659  double GamTmp = particleDataPtr->mWidth(idTmp);
660  sigmaTmp *= 2. * mTmp * GamTmp / ( pow2(sH - mTmp * mTmp)
661  + pow2(mTmp * GamTmp) );
662  }
663  if (convert2mb()) sigmaTmp *= CONVERT2MB;
664  return sigmaTmp;
665 
666 }
667 
668 //--------------------------------------------------------------------------
669 
670 // Input and complement kinematics for resolved 2 -> 1 process.
671 
672 void Sigma1Process::store1Kin( double x1in, double x2in, double sHin) {
673 
674  // Default value only sensible for these processes.
675  swapTU = false;
676 
677  // Incoming parton momentum fractions and sHat.
678  x1Save = x1in;
679  x2Save = x2in;
680  sH = sHin;
681  mH = sqrt(sH);
682  sH2 = sH * sH;
683 
684  // Different options for renormalization scale, but normally sHat.
685  Q2RenSave = renormMultFac * sH;
686  if (renormScale1 == 2) Q2RenSave = renormFixScale;
687 
688  // Different options for factorization scale, but normally sHat.
689  Q2FacSave = factorMultFac * sH;
690  if (factorScale1 == 2) Q2FacSave = factorFixScale;
691 
692  // Evaluate alpha_strong and alpha_EM.
693  alpS = couplingsPtr->alphaS(Q2RenSave);
694  alpEM = couplingsPtr->alphaEM(Q2RenSave);
695 
696 }
697 
698 //--------------------------------------------------------------------------
699 
700 // Calculate modified masses and four-vectors for matrix elements.
701 
702 bool Sigma1Process::setupForME() {
703 
704  // Common initial-state handling.
705  bool allowME = setupForMEin();
706 
707  // Final state trivial here.
708  mME[2] = mH;
709  pME[2] = Vec4( 0., 0., 0., mH);
710 
711  // Done.
712  return allowME;
713 
714 }
715 
716 //==========================================================================
717 
718 // The Sigma2Process class.
719 // Base class for resolved 2 -> 2 cross sections; derived from SigmaProcess.
720 
721 //--------------------------------------------------------------------------
722 
723 // Input and complement kinematics for resolved 2 -> 2 process.
724 
725 void Sigma2Process::store2Kin( double x1in, double x2in, double sHin,
726  double tHin, double m3in, double m4in, double runBW3in, double runBW4in) {
727 
728  // Default ordering of particles 3 and 4.
729  swapTU = false;
730 
731  // Incoming parton momentum fractions.
732  x1Save = x1in;
733  x2Save = x2in;
734 
735  // Incoming masses and their squares.
736  bool masslessKin = (id3Mass() == 0) && (id4Mass() == 0);
737  if (masslessKin) {
738  m3 = 0.;
739  m4 = 0.;
740  } else {
741  m3 = m3in;
742  m4 = m4in;
743  }
744  mSave[3] = m3;
745  mSave[4] = m4;
746  s3 = m3 * m3;
747  s4 = m4 * m4;
748 
749  // Standard Mandelstam variables and their squares.
750  sH = sHin;
751  tH = tHin;
752  uH = (masslessKin) ? -(sH + tH) : s3 + s4 - (sH + tH);
753  mH = sqrt(sH);
754  sH2 = sH * sH;
755  tH2 = tH * tH;
756  uH2 = uH * uH;
757 
758  // The nominal Breit-Wigner factors with running width.
759  runBW3 = runBW3in;
760  runBW4 = runBW4in;
761 
762  // Calculate squared transverse momentum.
763  pT2 = (masslessKin) ? tH * uH / sH : (tH * uH - s3 * s4) / sH;
764 
765  // Special case: pick scale as if 2 -> 1 process in disguise.
766  if (isSChannel()) {
767 
768  // Different options for renormalization scale, but normally sHat.
769  Q2RenSave = renormMultFac * sH;
770  if (renormScale1 == 2) Q2RenSave = renormFixScale;
771 
772  // Different options for factorization scale, but normally sHat.
773  Q2FacSave = factorMultFac * sH;
774  if (factorScale1 == 2) Q2FacSave = factorFixScale;
775 
776  // Normal case with "true" 2 -> 2.
777  } else {
778 
779  // Different options for renormalization scale.
780  if (masslessKin) Q2RenSave = (renormScale2 < 4) ? pT2 : sH;
781  else if (renormScale2 == 1) Q2RenSave = pT2 + min(s3, s4);
782  else if (renormScale2 == 2) Q2RenSave = sqrt((pT2 + s3) * (pT2 + s4));
783  else if (renormScale2 == 3) Q2RenSave = pT2 + 0.5 * (s3 + s4);
784  else Q2RenSave = sH;
785  Q2RenSave *= renormMultFac;
786  if (renormScale2 == 5) Q2RenSave = renormFixScale;
787 
788  // Different options for factorization scale.
789  if (masslessKin) Q2FacSave = (factorScale2 < 4) ? pT2 : sH;
790  else if (factorScale2 == 1) Q2FacSave = pT2 + min(s3, s4);
791  else if (factorScale2 == 2) Q2FacSave = sqrt((pT2 + s3) * (pT2 + s4));
792  else if (factorScale2 == 3) Q2FacSave = pT2 + 0.5 * (s3 + s4);
793  else Q2FacSave = sH;
794  Q2FacSave *= factorMultFac;
795  if (factorScale2 == 5) Q2FacSave = factorFixScale;
796  }
797 
798  // Evaluate alpha_strong and alpha_EM.
799  alpS = couplingsPtr->alphaS(Q2RenSave);
800  alpEM = couplingsPtr->alphaEM(Q2RenSave);
801 
802 }
803 
804 //--------------------------------------------------------------------------
805 
806 // As above, special kinematics for multiparton interactions.
807 
808 void Sigma2Process::store2KinMPI( double x1in, double x2in,
809  double sHin, double tHin, double uHin, double alpSin, double alpEMin,
810  bool needMasses, double m3in, double m4in) {
811 
812  // Default ordering of particles 3 and 4.
813  swapTU = false;
814 
815  // Incoming x values.
816  x1Save = x1in;
817  x2Save = x2in;
818 
819  // Standard Mandelstam variables and their squares.
820  sH = sHin;
821  tH = tHin;
822  uH = uHin;
823  mH = sqrt(sH);
824  sH2 = sH * sH;
825  tH2 = tH * tH;
826  uH2 = uH * uH;
827 
828  // Strong and electroweak couplings.
829  alpS = alpSin;
830  alpEM = alpEMin;
831 
832  // Assume vanishing masses. (Will be modified in final kinematics.)
833  m3 = 0.;
834  s3 = 0.;
835  m4 = 0.;
836  s4 = 0.;
837  sHBeta = sH;
838 
839  // Scattering angle.
840  cosTheta = (tH - uH) / sH;
841  sinTheta = 2. * sqrtpos( tH * uH ) / sH;
842 
843  // In some cases must use masses and redefine meaning of tHat and uHat.
844  if (needMasses) {
845  m3 = m3in;
846  s3 = m3 * m3;
847  m4 = m4in;
848  s4 = m4 * m4;
849  sHMass = sH - s3 - s4;
850  sHBeta = sqrtpos(sHMass*sHMass - 4. * s3 * s4);
851  tH = -0.5 * (sHMass - sHBeta * cosTheta);
852  uH = -0.5 * (sHMass + sHBeta * cosTheta);
853  tH2 = tH * tH;
854  uH2 = uH * uH;
855  }
856 
857  // pT2 with masses (at this stage) included.
858  pT2Mass = 0.25 * sHBeta * pow2(sinTheta);
859 
860 }
861 
862 //--------------------------------------------------------------------------
863 
864 // Perform kinematics for a multiparton interaction, including a rescattering.
865 
866 bool Sigma2Process::final2KinMPI( int i1Res, int i2Res, Vec4 p1Res, Vec4 p2Res,
867  double m1Res, double m2Res) {
868 
869  // Have to set flavours and colours.
870  setIdColAcol();
871 
872  // Check that masses of outgoing particles not too big.
873  m3 = particleDataPtr->m0(idSave[3]);
874  m4 = particleDataPtr->m0(idSave[4]);
875  mH = sqrt(sH);
876  if (m3 + m4 + MASSMARGIN > mH) return false;
877  s3 = m3 * m3;
878  s4 = m4 * m4;
879 
880  // Do kinematics of the production; without or with masses.
881  double e1In = 0.5 * mH;
882  double e2In = e1In;
883  double pzIn = e1In;
884  if (i1Res > 0 || i2Res > 0) {
885  double s1 = m1Res * m1Res;
886  double s2 = m2Res * m2Res;
887  e1In = 0.5 * (sH + s1 - s2) / mH;
888  e2In = 0.5 * (sH + s2 - s1) / mH;
889  pzIn = sqrtpos( e1In*e1In - s1 );
890  }
891 
892  // Do kinematics of the decay.
893  double e3 = 0.5 * (sH + s3 - s4) / mH;
894  double e4 = 0.5 * (sH + s4 - s3) / mH;
895  double pAbs = sqrtpos( e3*e3 - s3 );
896  phi = 2. * M_PI * rndmPtr->flat();
897  double pZ = pAbs * cosTheta;
898  pTFin = pAbs * sinTheta;
899  double pX = pTFin * sin(phi);
900  double pY = pTFin * cos(phi);
901  double scale = 0.5 * mH * sinTheta;
902 
903  // Fill particle info.
904  int status1 = (i1Res == 0) ? -31 : -34;
905  int status2 = (i2Res == 0) ? -31 : -34;
906  parton[1] = Particle( idSave[1], status1, 0, 0, 3, 4,
907  colSave[1], acolSave[1], 0., 0., pzIn, e1In, m1Res, scale);
908  parton[2] = Particle( idSave[2], status2, 0, 0, 3, 4,
909  colSave[2], acolSave[2], 0., 0., -pzIn, e2In, m2Res, scale);
910  parton[3] = Particle( idSave[3], 33, 1, 2, 0, 0,
911  colSave[3], acolSave[3], pX, pY, pZ, e3, m3, scale);
912  parton[4] = Particle( idSave[4], 33, 1, 2, 0, 0,
913  colSave[4], acolSave[4], -pX, -pY, -pZ, e4, m4, scale);
914 
915  // Boost particles from subprocess rest frame to event rest frame.
916  // Normal multiparton interaction: only longitudinal boost.
917  if (i1Res == 0 && i2Res == 0) {
918  double betaZ = (x1Save - x2Save) / (x1Save + x2Save);
919  for (int i = 1; i <= 4; ++i) parton[i].bst(0., 0., betaZ);
920  // Rescattering: generic rotation and boost required.
921  } else {
922  RotBstMatrix M;
923  M.fromCMframe( p1Res, p2Res);
924  for (int i = 1; i <= 4; ++i) parton[i].rotbst(M);
925  }
926 
927  // Done.
928  return true;
929 
930 }
931 
932 //--------------------------------------------------------------------------
933 
934 // Calculate modified masses and four-vectors for matrix elements.
935 
936 bool Sigma2Process::setupForME() {
937 
938  // Common initial-state handling.
939  bool allowME = setupForMEin();
940 
941  // Correct outgoing c, b, mu and tau to be massive or not.
942  mME[2] = m3;
943  int id3Tmp = abs(id3Mass());
944  if (id3Tmp == 4) mME[2] = mcME;
945  if (id3Tmp == 5) mME[2] = mbME;
946  if (id3Tmp == 13) mME[2] = mmuME;
947  if (id3Tmp == 15) mME[2] = mtauME;
948  mME[3] = m4;
949  int id4Tmp = abs(id4Mass());
950  if (id4Tmp == 4) mME[3] = mcME;
951  if (id4Tmp == 5) mME[3] = mbME;
952  if (id4Tmp == 13) mME[3] = mmuME;
953  if (id4Tmp == 15) mME[3] = mtauME;
954 
955  // If kinematically impossible turn to massless case, but set error.
956  if (mME[2] + mME[3] >= mH) {
957  mME[2] = 0.;
958  mME[3] = 0.;
959  allowME = false;
960  }
961 
962  // Calculate scattering angle in subsystem rest frame.
963  double sH34 = sqrtpos( pow2(sH - s3 - s4) - 4. * s3 * s4);
964  double cThe = (tH - uH) / sH34;
965  double sThe = sqrtpos(1. - cThe * cThe);
966 
967  // Setup massive kinematics with preserved scattering angle.
968  double s3ME = pow2(mME[2]);
969  double s4ME = pow2(mME[3]);
970  double sH34ME = sqrtpos( pow2(sH - s3ME - s4ME) - 4. * s3ME * s4ME);
971  double pAbsME = 0.5 * sH34ME / mH;
972 
973  // Normally allowed with unequal (or vanishing) masses.
974  if (id3Tmp == 0 || id3Tmp != id4Tmp) {
975  pME[2] = Vec4( pAbsME * sThe, 0., pAbsME * cThe,
976  0.5 * (sH + s3ME - s4ME) / mH);
977  pME[3] = Vec4( -pAbsME * sThe, 0., -pAbsME * cThe,
978  0.5 * (sH + s4ME - s3ME) / mH);
979 
980  // For equal (anti)particles (e.g. W+ W-) use averaged mass.
981  } else {
982  mME[2] = sqrtpos(0.5 * (s3ME + s4ME) - 0.25 * pow2(s3ME - s4ME) / sH);
983  mME[3] = mME[2];
984  pME[2] = Vec4( pAbsME * sThe, 0., pAbsME * cThe, 0.5 * mH);
985  pME[3] = Vec4( -pAbsME * sThe, 0., -pAbsME * cThe, 0.5 * mH);
986  }
987 
988  // Done.
989  return allowME;
990 
991 }
992 
993 //==========================================================================
994 
995 // The Sigma3Process class.
996 // Base class for resolved 2 -> 3 cross sections; derived from SigmaProcess.
997 
998 //--------------------------------------------------------------------------
999 
1000 // Input and complement kinematics for resolved 2 -> 3 process.
1001 
1002 void Sigma3Process::store3Kin( double x1in, double x2in, double sHin,
1003  Vec4 p3cmIn, Vec4 p4cmIn, Vec4 p5cmIn, double m3in, double m4in,
1004  double m5in, double runBW3in, double runBW4in, double runBW5in) {
1005 
1006  // Default ordering of particles 3 and 4 - not relevant here.
1007  swapTU = false;
1008 
1009  // Incoming parton momentum fractions.
1010  x1Save = x1in;
1011  x2Save = x2in;
1012 
1013  // Incoming masses and their squares.
1014  if (id3Mass() == 0 && id4Mass() == 0 && id5Mass() == 0) {
1015  m3 = 0.;
1016  m4 = 0.;
1017  m5 = 0.;
1018  } else {
1019  m3 = m3in;
1020  m4 = m4in;
1021  m5 = m5in;
1022  }
1023  mSave[3] = m3;
1024  mSave[4] = m4;
1025  mSave[5] = m5;
1026  s3 = m3 * m3;
1027  s4 = m4 * m4;
1028  s5 = m5 * m5;
1029 
1030  // Standard Mandelstam variables and four-momenta in rest frame.
1031  sH = sHin;
1032  mH = sqrt(sH);
1033  sH2 = sH * sH;
1034  p3cm = p3cmIn;
1035  p4cm = p4cmIn;
1036  p5cm = p5cmIn;
1037 
1038  // The nominal Breit-Wigner factors with running width.
1039  runBW3 = runBW3in;
1040  runBW4 = runBW4in;
1041  runBW5 = runBW5in;
1042 
1043  // Special case: pick scale as if 2 -> 1 process in disguise.
1044  if (isSChannel()) {
1045 
1046  // Different options for renormalization scale, but normally sHat.
1047  Q2RenSave = renormMultFac * sH;
1048  if (renormScale1 == 2) Q2RenSave = renormFixScale;
1049 
1050  // Different options for factorization scale, but normally sHat.
1051  Q2FacSave = factorMultFac * sH;
1052  if (factorScale1 == 2) Q2RenSave = factorFixScale;
1053 
1054  // "Normal" 2 -> 3 processes, i.e. not vector boson fusion.
1055  } else if ( idTchan1() != 23 && idTchan1() != 24 && idTchan2() != 23
1056  && idTchan2() != 24 ) {
1057  double mT3S = s3 + p3cm.pT2();
1058  double mT4S = s4 + p4cm.pT2();
1059  double mT5S = s5 + p5cm.pT2();
1060 
1061  // Different options for renormalization scale.
1062  if (renormScale3 == 1) Q2RenSave = min( mT3S, min(mT4S, mT5S) );
1063  else if (renormScale3 == 2) Q2RenSave = sqrt( mT3S * mT4S * mT5S
1064  / max( mT3S, max(mT4S, mT5S) ) );
1065  else if (renormScale3 == 3) Q2RenSave = pow( mT3S * mT4S * mT5S,
1066  1./3. );
1067  else if (renormScale3 == 4) Q2RenSave = (mT3S + mT4S + mT5S) / 3.;
1068  else Q2RenSave = sH;
1069  Q2RenSave *= renormMultFac;
1070  if (renormScale3 == 6) Q2RenSave = renormFixScale;
1071 
1072  // Different options for factorization scale.
1073  if (factorScale3 == 1) Q2FacSave = min( mT3S, min(mT4S, mT5S) );
1074  else if (factorScale3 == 2) Q2FacSave = sqrt( mT3S * mT4S * mT5S
1075  / max( mT3S, max(mT4S, mT5S) ) );
1076  else if (factorScale3 == 3) Q2FacSave = pow( mT3S * mT4S * mT5S,
1077  1./3. );
1078  else if (factorScale3 == 4) Q2FacSave = (mT3S + mT4S + mT5S) / 3.;
1079  else Q2FacSave = sH;
1080  Q2FacSave *= factorMultFac;
1081  if (factorScale3 == 6) Q2FacSave = factorFixScale;
1082 
1083  // Vector boson fusion 2 -> 3 processes; recoils in positions 4 and 5.
1084  } else {
1085  double sV4 = pow2( particleDataPtr->m0(idTchan1()) );
1086  double sV5 = pow2( particleDataPtr->m0(idTchan2()) );
1087  double mT3S = s3 + p3cm.pT2();
1088  double mTV4S = sV4 + p4cm.pT2();
1089  double mTV5S = sV5 + p5cm.pT2();
1090 
1091  // Different options for renormalization scale.
1092  if (renormScale3VV == 1) Q2RenSave = max( sV4, sV5);
1093  else if (renormScale3VV == 2) Q2RenSave = sqrt( mTV4S * mTV5S );
1094  else if (renormScale3VV == 3) Q2RenSave = pow( mT3S * mTV4S * mTV5S,
1095  1./3. );
1096  else if (renormScale3VV == 4) Q2RenSave = (mT3S * mTV4S * mTV5S) / 3.;
1097  else Q2RenSave = sH;
1098  Q2RenSave *= renormMultFac;
1099  if (renormScale3VV == 6) Q2RenSave = renormFixScale;
1100 
1101  // Different options for factorization scale.
1102  if (factorScale3VV == 1) Q2FacSave = max( sV4, sV5);
1103  else if (factorScale3VV == 2) Q2FacSave = sqrt( mTV4S * mTV5S );
1104  else if (factorScale3VV == 3) Q2FacSave = pow( mT3S * mTV4S * mTV5S,
1105  1./3. );
1106  else if (factorScale3VV == 4) Q2FacSave = (mT3S * mTV4S * mTV5S) / 3.;
1107  else Q2FacSave = sH;
1108  Q2FacSave *= factorMultFac;
1109  if (factorScale3VV == 6) Q2FacSave = factorFixScale;
1110  }
1111 
1112  // Evaluate alpha_strong and alpha_EM.
1113  alpS = couplingsPtr->alphaS(Q2RenSave);
1114  alpEM = couplingsPtr->alphaEM(Q2RenSave);
1115 
1116 }
1117 
1118 //--------------------------------------------------------------------------
1119 
1120 // Calculate modified masses and four-vectors for matrix elements.
1121 
1122 bool Sigma3Process::setupForME() {
1123 
1124  // Common initial-state handling.
1125  bool allowME = setupForMEin();
1126 
1127  // Correct outgoing c, b, mu and tau to be massive or not.
1128  mME[2] = m3;
1129  int id3Tmp = abs(id3Mass());
1130  if (id3Tmp == 4) mME[2] = mcME;
1131  if (id3Tmp == 5) mME[2] = mbME;
1132  if (id3Tmp == 13) mME[2] = mmuME;
1133  if (id3Tmp == 15) mME[2] = mtauME;
1134  mME[3] = m4;
1135  int id4Tmp = abs(id4Mass());
1136  if (id4Tmp == 4) mME[3] = mcME;
1137  if (id4Tmp == 5) mME[3] = mbME;
1138  if (id4Tmp == 13) mME[3] = mmuME;
1139  if (id4Tmp == 15) mME[3] = mtauME;
1140  mME[4] = m5;
1141  int id5Tmp = abs(id5Mass());
1142  if (id5Tmp == 4) mME[4] = mcME;
1143  if (id5Tmp == 5) mME[4] = mbME;
1144  if (id5Tmp == 13) mME[4] = mmuME;
1145  if (id5Tmp == 15) mME[4] = mtauME;
1146 
1147  // If kinematically impossible turn to massless case, but set error.
1148  if (mME[2] + mME[3] + mME[4] >= mH) {
1149  mME[2] = 0.;
1150  mME[3] = 0.;
1151  mME[4] = 0.;
1152  allowME = false;
1153  }
1154 
1155  // Form new average masses if identical particles.
1156  if (id3Tmp != 0 && id4Tmp == id3Tmp && id5Tmp == id3Tmp) {
1157  double mAvg = (mME[2] + mME[3] + mME[4]) / 3.;
1158  mME[2] = mAvg;
1159  mME[3] = mAvg;
1160  mME[4] = mAvg;
1161  } else if (id3Tmp != 0 && id4Tmp == id3Tmp) {
1162  mME[2] = sqrtpos(0.5 * (pow2(mME[2]) + pow2(mME[3]))
1163  - 0.25 * pow2(pow2(mME[2]) - pow2(mME[3])) / sH);
1164  mME[3] = mME[2];
1165  } else if (id3Tmp != 0 && id5Tmp == id3Tmp) {
1166  mME[2] = sqrtpos(0.5 * (pow2(mME[2]) + pow2(mME[4]))
1167  - 0.25 * pow2(pow2(mME[2]) - pow2(mME[4])) / sH);
1168  mME[4] = mME[2];
1169  } else if (id4Tmp != 0 && id5Tmp == id4Tmp) {
1170  mME[3] = sqrtpos(0.5 * (pow2(mME[3]) + pow2(mME[4]))
1171  - 0.25 * pow2(pow2(mME[3]) - pow2(mME[4])) / sH);
1172  mME[4] = mME[2];
1173  }
1174 
1175  // Iterate rescaled three-momenta until convergence.
1176  double m2ME3 = pow2(mME[2]);
1177  double m2ME4 = pow2(mME[3]);
1178  double m2ME5 = pow2(mME[4]);
1179  double p2ME3 = p3cm.pAbs2();
1180  double p2ME4 = p4cm.pAbs2();
1181  double p2ME5 = p5cm.pAbs2();
1182  double p2sum = p2ME3 + p2ME4 + p2ME5;
1183  double eME3 = sqrt(m2ME3 + p2ME3);
1184  double eME4 = sqrt(m2ME4 + p2ME4);
1185  double eME5 = sqrt(m2ME5 + p2ME5);
1186  double esum = eME3 + eME4 + eME5;
1187  double p2rat = p2ME3 / eME3 + p2ME4 / eME4 + p2ME5 / eME5;
1188  int iStep = 0;
1189  while ( abs(esum - mH) > COMPRELERR * mH && iStep < NCOMPSTEP ) {
1190  ++iStep;
1191  double compFac = 1. + 2. * (mH - esum) / p2rat;
1192  p2ME3 *= compFac;
1193  p2ME4 *= compFac;
1194  p2ME5 *= compFac;
1195  eME3 = sqrt(m2ME3 + p2ME3);
1196  eME4 = sqrt(m2ME4 + p2ME4);
1197  eME5 = sqrt(m2ME5 + p2ME5);
1198  esum = eME3 + eME4 + eME5;
1199  p2rat = p2ME3 / eME3 + p2ME4 / eME4 + p2ME5 / eME5;
1200  }
1201 
1202  // If failed convergence set error flag.
1203  if (abs(esum - mH) > COMPRELERR * mH) allowME = false;
1204 
1205  // Set up accepted kinematics.
1206  double totFac = sqrt( (p2ME3 + p2ME4 + p2ME5) / p2sum);
1207  pME[2] = totFac * p3cm;
1208  pME[2].e( eME3);
1209  pME[3] = totFac * p4cm;
1210  pME[3].e( eME4);
1211  pME[4] = totFac * p5cm;
1212  pME[4].e( eME5);
1213 
1214  // Done.
1215  return allowME;
1216 
1217 }
1218 
1219 //==========================================================================
1220 
1221 // The SigmaLHAProcess class.
1222 // Wrapper for Les Houches Accord external input; derived from SigmaProcess.
1223 // Note: arbitrary subdivision into PhaseSpaceLHA and SigmaLHAProcess tasks.
1224 
1225 //--------------------------------------------------------------------------
1226 
1227 // Set scale, alpha_strong and alpha_EM when not set.
1228 
1229 void SigmaLHAProcess::setScale() {
1230 
1231  // If scale has not been set, then to set.
1232  double scaleLHA = lhaUpPtr->scale();
1233  if (scaleLHA < 0.) {
1234 
1235  // Final-state partons and their invariant mass.
1236  vector<int> iFin;
1237  Vec4 pFinSum;
1238  for (int i = 3; i < lhaUpPtr->sizePart(); ++i)
1239  if (lhaUpPtr->mother1(i) == 1) {
1240  iFin.push_back(i);
1241  pFinSum += Vec4( lhaUpPtr->px(i), lhaUpPtr->py(i),
1242  lhaUpPtr->pz(i), lhaUpPtr->e(i) );
1243  }
1244  int nFin = iFin.size();
1245  sH = pFinSum * pFinSum;
1246  mH = sqrt(sH);
1247  sH2 = sH * sH;
1248 
1249  // If 1 final-state particle then use Sigma1Process logic.
1250  if (nFin == 1) {
1251  Q2RenSave = renormMultFac * sH;
1252  if (renormScale1 == 2) Q2RenSave = renormFixScale;
1253  Q2FacSave = factorMultFac * sH;
1254  if (factorScale1 == 2) Q2FacSave = factorFixScale;
1255 
1256  // If 2 final-state particles then use Sigma2Process logic.
1257  } else if (nFin == 2) {
1258  double s3 = pow2(lhaUpPtr->m(iFin[0]));
1259  double s4 = pow2(lhaUpPtr->m(iFin[1]));
1260  double pT2 = pow2(lhaUpPtr->px(iFin[0])) + pow2(lhaUpPtr->py(iFin[0]));
1261  if (renormScale2 == 1) Q2RenSave = pT2 + min(s3, s4);
1262  else if (renormScale2 == 2) Q2RenSave = sqrt((pT2 + s3) * (pT2 + s4));
1263  else if (renormScale2 == 3) Q2RenSave = pT2 + 0.5 * (s3 + s4);
1264  else Q2RenSave = sH;
1265  Q2RenSave *= renormMultFac;
1266  if (renormScale2 == 5) Q2RenSave = renormFixScale;
1267  if (factorScale2 == 1) Q2FacSave = pT2 + min(s3, s4);
1268  else if (factorScale2 == 2) Q2FacSave = sqrt((pT2 + s3) * (pT2 + s4));
1269  else if (factorScale2 == 3) Q2FacSave = pT2 + 0.5 * (s3 + s4);
1270  else Q2FacSave = sH;
1271  Q2FacSave *= factorMultFac;
1272  if (factorScale2 == 5) Q2FacSave = factorFixScale;
1273 
1274  // If 3 or more final-state particles then use Sigma3Process logic.
1275  } else {
1276  double mTSlow = sH;
1277  double mTSmed = sH;
1278  double mTSprod = 1.;
1279  double mTSsum = 0.;
1280  for (int i = 0; i < nFin; ++i) {
1281  double mTSnow = pow2(lhaUpPtr->m(iFin[i]))
1282  + pow2(lhaUpPtr->px(iFin[i])) + pow2(lhaUpPtr->py(iFin[i]));
1283  if (mTSnow < mTSlow) {mTSmed = mTSlow; mTSlow = mTSnow;}
1284  else if (mTSnow < mTSmed) mTSmed = mTSnow;
1285  mTSprod *= mTSnow;
1286  mTSsum += mTSnow;
1287  }
1288  if (renormScale3 == 1) Q2RenSave = mTSlow;
1289  else if (renormScale3 == 2) Q2RenSave = sqrt(mTSlow * mTSmed);
1290  else if (renormScale3 == 3) Q2RenSave = pow(mTSprod, 1. / nFin);
1291  else if (renormScale3 == 4) Q2RenSave = mTSsum / nFin;
1292  else Q2RenSave = sH;
1293  Q2RenSave *= renormMultFac;
1294  if (renormScale3 == 6) Q2RenSave = renormFixScale;
1295  if (factorScale3 == 1) Q2FacSave = mTSlow;
1296  else if (factorScale3 == 2) Q2FacSave = sqrt(mTSlow * mTSmed);
1297  else if (factorScale3 == 3) Q2FacSave = pow(mTSprod, 1. / nFin);
1298  else if (factorScale3 == 4) Q2FacSave = mTSsum / nFin;
1299  else Q2FacSave = sH;
1300  Q2FacSave *= factorMultFac;
1301  if (factorScale3 == 6) Q2FacSave = factorFixScale;
1302  }
1303  }
1304 
1305  // If alpha_strong and alpha_EM have not been set, then set them.
1306  if (lhaUpPtr->alphaQCD() < 0.001) {
1307  double Q2RenNow = (scaleLHA < 0.) ? Q2RenSave : pow2(scaleLHA);
1308  alpS = couplingsPtr->alphaS(Q2RenNow);
1309  }
1310  if (lhaUpPtr->alphaQED() < 0.001) {
1311  double Q2RenNow = (scaleLHA < 0.) ? Q2RenSave : pow2(scaleLHA);
1312  alpEM = couplingsPtr->alphaEM(Q2RenNow);
1313  }
1314 
1315 }
1316 
1317 //--------------------------------------------------------------------------
1318 
1319 // Obtain number of final-state partons from LHA object.
1320 
1321 int SigmaLHAProcess::nFinal() const {
1322 
1323  // At initialization size unknown, so return 0.
1324  if (lhaUpPtr->sizePart() <= 0) return 0;
1325 
1326  // Sum up all particles that has first mother = 1.
1327  int nFin = 0;
1328  for (int i = 3; i < lhaUpPtr->sizePart(); ++i)
1329  if (lhaUpPtr->mother1(i) == 1) ++nFin;
1330  return nFin;
1331 
1332 }
1333 
1334 //==========================================================================
1335 
1336 } // end namespace Pythia8
Definition: AgUStep.h:26