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