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