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