StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ProcessContainer.cc
1 // ProcessContainer.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 // ProcessContainer and SetupContainers classes.
8 
9 #include "Pythia8/ProcessContainer.h"
10 
11 // Internal headers for special processes.
12 #include "Pythia8/SigmaCompositeness.h"
13 #include "Pythia8/SigmaEW.h"
14 #include "Pythia8/SigmaExtraDim.h"
15 #include "Pythia8/SigmaGeneric.h"
16 #include "Pythia8/SigmaHiggs.h"
17 #include "Pythia8/SigmaLeftRightSym.h"
18 #include "Pythia8/SigmaLeptoquark.h"
19 #include "Pythia8/SigmaNewGaugeBosons.h"
20 #include "Pythia8/SigmaQCD.h"
21 #include "Pythia8/SigmaSUSY.h"
22 
23 namespace Pythia8 {
24 
25 //==========================================================================
26 
27 // ProcessContainer class.
28 // Information allowing the generation of a specific process.
29 
30 //--------------------------------------------------------------------------
31 
32 // Constants: could be changed here if desired, but normally should not.
33 // These are of technical nature, as described for each.
34 
35 // Number of event tries to check maximization finding reliability.
36 const int ProcessContainer::N12SAMPLE = 100;
37 
38 // Ditto, but increased for 2 -> 3 processes.
39 const int ProcessContainer::N3SAMPLE = 1000;
40 
41 //--------------------------------------------------------------------------
42 
43 // Initialize phase space and counters.
44 // Argument isFirst distinguishes two hard processes in same event.
45 
46 bool ProcessContainer::init(bool isFirst, Info* infoPtrIn,
47  Settings& settings, ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
48  BeamParticle* beamAPtr, BeamParticle* beamBPtr, Couplings* couplingsPtr,
49  SigmaTotal* sigmaTotPtr, ResonanceDecays* resDecaysPtrIn,
50  SLHAinterface* slhaInterfacePtr, UserHooks* userHooksPtrIn) {
51 
52  // Extract info about current process from SigmaProcess object.
53  isLHA = sigmaProcessPtr->isLHA();
54  isNonDiff = sigmaProcessPtr->isNonDiff();
55  isResolved = sigmaProcessPtr->isResolved();
56  isDiffA = sigmaProcessPtr->isDiffA();
57  isDiffB = sigmaProcessPtr->isDiffB();
58  isDiffC = sigmaProcessPtr->isDiffC();
59  isQCD3body = sigmaProcessPtr->isQCD3body();
60  int nFin = sigmaProcessPtr->nFinal();
61  lhaStrat = (isLHA) ? lhaUpPtr->strategy() : 0;
62  lhaStratAbs = abs(lhaStrat);
63  allowNegSig = sigmaProcessPtr->allowNegativeSigma();
64 
65  // Flag for maximum violation handling.
66  increaseMaximum = settings.flag("PhaseSpace:increaseMaximum");
67 
68  // Pick and create phase space generator. Send pointers where required.
69  if (isLHA) phaseSpacePtr = new PhaseSpaceLHA();
70  else if (isNonDiff) phaseSpacePtr = new PhaseSpace2to2nondiffractive();
71  else if (!isResolved && !isDiffA && !isDiffB && !isDiffC )
72  phaseSpacePtr = new PhaseSpace2to2elastic();
73  else if (!isResolved && !isDiffA && !isDiffB && isDiffC)
74  phaseSpacePtr = new PhaseSpace2to3diffractive();
75  else if (!isResolved) phaseSpacePtr = new PhaseSpace2to2diffractive(
76  isDiffA, isDiffB);
77  else if (nFin == 1) phaseSpacePtr = new PhaseSpace2to1tauy();
78  else if (nFin == 2) phaseSpacePtr = new PhaseSpace2to2tauyz();
79  else if (isQCD3body) phaseSpacePtr = new PhaseSpace2to3yyycyl();
80  else phaseSpacePtr = new PhaseSpace2to3tauycyl();
81 
82  // Store pointers and perform simple initialization.
83  infoPtr = infoPtrIn;
84  particleDataPtr = particleDataPtrIn;
85  rndmPtr = rndmPtrIn;
86  resDecaysPtr = resDecaysPtrIn;
87  userHooksPtr = userHooksPtrIn;
88  canVetoResDecay = (userHooksPtr != 0)
89  ? userHooksPtr->canVetoResonanceDecays() : false;
90  if (isLHA) {
91  sigmaProcessPtr->setLHAPtr(lhaUpPtr);
92  phaseSpacePtr->setLHAPtr(lhaUpPtr);
93  }
94  sigmaProcessPtr->init(infoPtr, &settings, particleDataPtr, rndmPtr,
95  beamAPtr, beamBPtr, couplingsPtr, sigmaTotPtr, slhaInterfacePtr);
96  phaseSpacePtr->init( isFirst, sigmaProcessPtr, infoPtr, &settings,
97  particleDataPtr, rndmPtr, beamAPtr, beamBPtr, couplingsPtr, sigmaTotPtr,
98  userHooksPtr);
99 
100  // Reset cross section statistics.
101  nTry = 0;
102  nSel = 0;
103  nAcc = 0;
104  nTryStat = 0;
105  sigmaMx = 0.;
106  sigmaSum = 0.;
107  sigma2Sum = 0.;
108  sigmaNeg = 0.;
109  sigmaAvg = 0.;
110  sigmaFin = 0.;
111  deltaFin = 0.;
112  wtAccSum = 0.;
113 
114  // Initialize process and allowed incoming partons.
115  sigmaProcessPtr->initProc();
116  if (!sigmaProcessPtr->initFlux()) return false;
117 
118  // Find maximum of differential cross section * phasespace.
119  bool physical = phaseSpacePtr->setupSampling();
120  sigmaMx = phaseSpacePtr->sigmaMax();
121  double sigmaHalfWay = sigmaMx;
122 
123  // Separate signed maximum needed for LHA with negative weight.
124  sigmaSgn = phaseSpacePtr->sigmaSumSigned();
125 
126  // Check maximum by a few events, and extrapolate a further increase.
127  if (physical & !isLHA) {
128  int nSample = (nFin < 3) ? N12SAMPLE : N3SAMPLE;
129  for (int iSample = 0; iSample < nSample; ++iSample) {
130  bool test = false;
131  while (!test) test = phaseSpacePtr->trialKin(false);
132  if (iSample == nSample/2) sigmaHalfWay = phaseSpacePtr->sigmaMax();
133  }
134  double sigmaFullWay = phaseSpacePtr->sigmaMax();
135  sigmaMx = (sigmaHalfWay > 0.) ? pow2(sigmaFullWay) / sigmaHalfWay
136  : sigmaFullWay;
137  phaseSpacePtr->setSigmaMax(sigmaMx);
138  }
139 
140  // Allow Pythia to overwrite incoming beams or Les Houches lifetime input.
141  idRenameBeams = settings.mode("LesHouches:idRenameBeams");
142  setLifetime = settings.mode("LesHouches:setLifetime");
143 
144  // Done.
145  return physical;
146 }
147 
148 //--------------------------------------------------------------------------
149 
150 // Generate a trial event; selected or not.
151 
152 bool ProcessContainer::trialProcess() {
153 
154  // Loop over tries only occurs for Les Houches strategy = +-2.
155  for (int iTry = 0; ; ++iTry) {
156 
157  // Generate a trial phase space point, if meaningful.
158  if (sigmaMx == 0.) return false;
159  infoPtr->setEndOfFile(false);
160  bool repeatSame = (iTry > 0);
161  bool physical = phaseSpacePtr->trialKin(true, repeatSame);
162 
163  // Note if at end of Les Houches file, else do statistics.
164  if (isLHA && !physical) infoPtr->setEndOfFile(true);
165  else {
166  ++nTry;
167  // Statistics for separate Les Houches process codes. Insert new codes.
168  if (isLHA) {
169  int codeLHANow = lhaUpPtr->idProcess();
170  int iFill = -1;
171  for (int i = 0; i < int(codeLHA.size()); ++i)
172  if (codeLHANow == codeLHA[i]) iFill = i;
173  if (iFill >= 0) {
174  ++nTryLHA[iFill];
175  } else {
176  codeLHA.push_back(codeLHANow);
177  nTryLHA.push_back(1);
178  nSelLHA.push_back(0);
179  nAccLHA.push_back(0);
180  for (int i = int(codeLHA.size()) - 1; i > 0; --i) {
181  if (codeLHA[i] < codeLHA[i - 1]) {
182  swap(codeLHA[i], codeLHA[i - 1]);
183  swap(nTryLHA[i], nTryLHA[i - 1]);
184  swap(nSelLHA[i], nSelLHA[i - 1]);
185  swap(nAccLHA[i], nAccLHA[i - 1]);
186  }
187  else break;
188  }
189  }
190  }
191  }
192 
193  // Possibly fail, else cross section.
194  if (!physical) return false;
195  double sigmaNow = phaseSpacePtr->sigmaNow();
196 
197  // Tell if this event comes with weight from cross section.
198  double sigmaWeight = 1.;
199  if (!isLHA && !increaseMaximum && sigmaNow > sigmaMx)
200  sigmaWeight = sigmaNow / sigmaMx;
201  if ( lhaStrat < 0 && sigmaNow < 0.) sigmaWeight = -1.;
202  if ( lhaStratAbs == 4) sigmaWeight = sigmaNow;
203 
204  // Also compensating weight from biased phase-space selection.
205  double biasWeight = phaseSpacePtr->biasSelectionWeight();
206  weightNow = sigmaWeight * biasWeight;
207  infoPtr->setWeight( weightNow, lhaStrat);
208 
209  // Check that not negative cross section when not allowed.
210  if (!allowNegSig) {
211  if (sigmaNow < sigmaNeg) {
212  infoPtr->errorMsg("Warning in ProcessContainer::trialProcess: neg"
213  "ative cross section set 0", "for " + sigmaProcessPtr->name() );
214  sigmaNeg = sigmaNow;
215  }
216  if (sigmaNow < 0.) sigmaNow = 0.;
217  }
218 
219  // Cross section of process may come with a weight. Update sum.
220  double sigmaAdd = sigmaNow * biasWeight;
221  if (lhaStratAbs == 2 || lhaStratAbs == 3) sigmaAdd = sigmaSgn;
222  sigmaSum += sigmaAdd;
223  sigma2Sum += pow2(sigmaAdd);
224 
225  // Check if maximum violated.
226  newSigmaMx = phaseSpacePtr->newSigmaMax();
227  if (newSigmaMx) sigmaMx = phaseSpacePtr->sigmaMax();
228 
229  // Select or reject trial point. Statistics.
230  bool select = true;
231  if (lhaStratAbs < 3) select
232  = (newSigmaMx || rndmPtr->flat() * abs(sigmaMx) < abs(sigmaNow));
233  if (select) {
234  ++nSel;
235  if (isLHA) {
236  int codeLHANow = lhaUpPtr->idProcess();
237  int iFill = -1;
238  for (int i = 0; i < int(codeLHA.size()); ++i)
239  if (codeLHANow == codeLHA[i]) iFill = i;
240  if (iFill >= 0) ++nSelLHA[iFill];
241  }
242  }
243  if (select || lhaStratAbs != 2) return select;
244 
245  }
246 
247 }
248 
249 //--------------------------------------------------------------------------
250 
251 // Accumulate statistics after user veto, including LHA code.
252 
253 void ProcessContainer::accumulate() {
254 
255  ++nAcc;
256  wtAccSum += weightNow;
257  if (isLHA) {
258  int codeLHANow = lhaUpPtr->idProcess();
259  int iFill = -1;
260  for (int i = 0; i < int(codeLHA.size()); ++i)
261  if (codeLHANow == codeLHA[i]) iFill = i;
262  if (iFill >= 0) ++nAccLHA[iFill];
263  }
264 
265 }
266 
267 //--------------------------------------------------------------------------
268 
269 // Pick flavours and colour flow of process.
270 
271 bool ProcessContainer::constructState() {
272 
273  // Construct flavour and colours for selected event.
274  if (isResolved && !isNonDiff) sigmaProcessPtr->pickInState();
275  sigmaProcessPtr->setIdColAcol();
276 
277  // Done.
278  return true;
279 
280 }
281 
282 //--------------------------------------------------------------------------
283 
284 // Give the hard subprocess with kinematics.
285 
286 bool ProcessContainer::constructProcess( Event& process, bool isHardest) {
287 
288  // Construct kinematics from selected phase space point.
289  if (!phaseSpacePtr->finalKin()) return false;
290  int nFin = sigmaProcessPtr->nFinal();
291 
292  // Basic info on process.
293  if (isHardest) infoPtr->setType( name(), code(), nFin, isNonDiff,
294  isResolved, isDiffA, isDiffB, isDiffC, isLHA);
295 
296  // Let hard process record begin with the event as a whole and
297  // the two incoming beam particles.
298  int idA = infoPtr->idA();
299  if (abs(idA) == idRenameBeams) idA = 16;
300  int idB = infoPtr->idB();
301  if (abs(idB) == idRenameBeams) idB = -16;
302  process.append( 90, -11, 0, 0, 0, 0, 0, 0,
303  Vec4(0., 0., 0., infoPtr->eCM()), infoPtr->eCM(), 0. );
304  process.append( idA, -12, 0, 0, 0, 0, 0, 0,
305  Vec4(0., 0., infoPtr->pzA(), infoPtr->eA()), infoPtr->mA(), 0. );
306  process.append( idB, -12, 0, 0, 0, 0, 0, 0,
307  Vec4(0., 0., infoPtr->pzB(), infoPtr->eB()), infoPtr->mB(), 0. );
308 
309  // For nondiffractive process no interaction selected so far, so done.
310  if (isNonDiff) return true;
311 
312  // Entries 3 and 4, now to be added, come from 1 and 2.
313  process[1].daughter1(3);
314  process[2].daughter1(4);
315  double scale = 0.;
316  double scalup = 0.;
317 
318  // For DiffC entries 3 - 5 come jointly from 1 and 2 (to keep HepMC happy).
319  if (isDiffC) {
320  process[1].daughters(3, 5);
321  process[2].daughters(3, 5);
322  }
323 
324  // Insert the subprocess partons - resolved processes.
325  int idRes = sigmaProcessPtr->idSChannel();
326  if (isResolved && !isLHA) {
327 
328  // NOAM: Mothers and daughters without/with intermediate state.
329  int m_M1 = 3;
330  int m_M2 = 4;
331  int m_D1 = 5;
332  int m_D2 = (nFin == 1) ? 0 : 4 + nFin;
333  if (idRes != 0) {
334  m_M1 = 5;
335  m_M2 = 0;
336  m_D1 = 5;
337  m_D2 = 0;
338  }
339 
340  // Find scale from which to begin MPI/ISR/FSR evolution.
341  scale = sqrt(Q2Fac());
342  process.scale( scale );
343 
344  // Loop over incoming and outgoing partons.
345  int colOffset = process.lastColTag();
346  for (int i = 1; i <= 2 + nFin; ++i) {
347 
348  // Read out particle info from SigmaProcess object.
349  int id = sigmaProcessPtr->id(i);
350  int status = (i <= 2) ? -21 : 23;
351  int mother1 = (i <= 2) ? i : m_M1 ;
352  int mother2 = (i <= 2) ? 0 : m_M2 ;
353  int daughter1 = (i <= 2) ? m_D1 : 0;
354  int daughter2 = (i <= 2) ? m_D2 : 0;
355  int col = sigmaProcessPtr->col(i);
356  if (col > 0) col += colOffset;
357  else if (col < 0) col -= colOffset;
358  int acol = sigmaProcessPtr->acol(i);
359  if (acol > 0) acol += colOffset;
360  else if (acol < 0) acol -= colOffset;
361 
362  // Append to process record.
363  int iNow = process.append( id, status, mother1, mother2,
364  daughter1, daughter2, col, acol, phaseSpacePtr->p(i),
365  phaseSpacePtr->m(i), scale);
366 
367  // NOAM: If there is an intermediate state, insert the it in
368  // the process record after the two incoming particles.
369  if (i == 2 && idRes != 0) {
370 
371  // Sign of intermediate state: go by charge.
372  if (particleDataPtr->hasAnti(idRes)
373  && process[3].chargeType() + process[4].chargeType() < 0)
374  idRes *= -1;
375 
376  // The colour configuration of the intermediate state has to be
377  // resolved separately.
378  col = 0;
379  acol = 0;
380  int m_col1 = sigmaProcessPtr->col(1);
381  int m_acol1 = sigmaProcessPtr->acol(1);
382  int m_col2 = sigmaProcessPtr->col(2);
383  int m_acol2 = sigmaProcessPtr->acol(2);
384  if (m_col1 == m_acol2 && m_col2 != m_acol1) {
385  col = m_col2;
386  acol = m_acol1;
387  } else if (m_col2 == m_acol1 && m_col1 != m_acol2) {
388  col = m_col1;
389  acol = m_acol2;
390  }
391  if ( col > 0) col += colOffset;
392  else if ( col < 0) col -= colOffset;
393  if (acol > 0) acol += colOffset;
394  else if (acol < 0) acol -= colOffset;
395 
396  // Insert the intermediate state into the event record.
397  Vec4 pIntMed = phaseSpacePtr->p(1) + phaseSpacePtr->p(2);
398  process.append( idRes, -22, 3, 4, 6, 5 + nFin, col, acol,
399  pIntMed, pIntMed.mCalc(), scale);
400  }
401 
402  // Pick lifetime where relevant, else not.
403  if (process[iNow].tau0() > 0.) process[iNow].tau(
404  process[iNow].tau0() * rndmPtr->exp() );
405  }
406  }
407 
408  // Insert the outgoing particles - unresolved processes.
409  else if (!isLHA) {
410  int id3 = sigmaProcessPtr->id(3);
411  int status3 = (id3 == process[1].id()) ? 14 : 15;
412  process.append( id3, status3, 1, 0, 0, 0, 0, 0,
413  phaseSpacePtr->p(3), phaseSpacePtr->m(3));
414  int id4 = sigmaProcessPtr->id(4);
415  int status4 = (id4 == process[2].id()) ? 14 : 15;
416  process.append( id4, status4, 2, 0, 0, 0, 0, 0,
417  phaseSpacePtr->p(4), phaseSpacePtr->m(4));
418 
419  // For central diffraction: two scattered protons inserted so far,
420  // but modify mothers, add also centrally-produced hadronic system.
421  if (isDiffC) {
422  process[3].mothers( 1, 2);
423  process[4].mothers( 1, 2);
424  int id5 = sigmaProcessPtr->id(5);
425  int status5 = 15;
426  process.append( id5, status5, 1, 2, 0, 0, 0, 0,
427  phaseSpacePtr->p(5), phaseSpacePtr->m(5));
428  }
429  }
430 
431  // Insert the outgoing particles - Les Houches Accord processes.
432  else {
433 
434  // Since LHA partons may be out of order, determine correct one.
435  // (Recall that zeroth particle is empty.)
436  vector<int> newPos;
437  newPos.reserve(lhaUpPtr->sizePart());
438  newPos.push_back(0);
439  for (int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
440  // For iNew == 0 look for the two incoming partons, then for
441  // partons having them as mothers, and so on layer by layer.
442  for (int i = 1; i < lhaUpPtr->sizePart(); ++i)
443  if (lhaUpPtr->mother1(i) == newPos[iNew]) newPos.push_back(i);
444  if (int(newPos.size()) <= iNew) break;
445  }
446 
447  // Find scale from which to begin MPI/ISR/FSR evolution.
448  scalup = lhaUpPtr->scale();
449  scale = scalup;
450  double scalePr = (scale < 0.) ? sqrt(Q2Fac()) : scale;
451  process.scale( scalePr);
452 
453  // Copy over info from LHA event to process, in proper order.
454  for (int i = 1; i < lhaUpPtr->sizePart(); ++i) {
455  int iOld = newPos[i];
456  int id = lhaUpPtr->id(iOld);
457  if (i == 1 && abs(id) == idRenameBeams) id = 16;
458  if (i == 2 && abs(id) == idRenameBeams) id = -16;
459 
460  // Translate from LHA status codes.
461  int lhaStatus = lhaUpPtr->status(iOld);
462  int status = -21;
463  if (lhaStatus == 2 || lhaStatus == 3) status = -22;
464  if (lhaStatus == 1) status = 23;
465 
466  // Find where mothers have been moved by reordering.
467  int mother1Old = lhaUpPtr->mother1(iOld);
468  int mother2Old = lhaUpPtr->mother2(iOld);
469  int mother1 = 0;
470  int mother2 = 0;
471  for (int im = 1; im < i; ++im) {
472  if (mother1Old == newPos[im]) mother1 = im + 2;
473  if (mother2Old == newPos[im]) mother2 = im + 2;
474  }
475  if (i <= 2) mother1 = i;
476 
477  // Ensure that second mother = 0 except for bona fide carbon copies.
478  if (mother1 > 0 && mother2 == mother1) {
479  int sister1 = process[mother1].daughter1();
480  int sister2 = process[mother1].daughter2();
481  if (sister2 != sister1 && sister2 != 0) mother2 = 0;
482  }
483 
484  // Find daughters and where they have been moved by reordering.
485  // (Values shifted two steps to account for inserted beams.)
486  int daughter1 = 0;
487  int daughter2 = 0;
488  for (int im = i + 1; im < lhaUpPtr->sizePart(); ++im) {
489  if (lhaUpPtr->mother1(newPos[im]) == iOld
490  || lhaUpPtr->mother2(newPos[im]) == iOld) {
491  if (daughter1 == 0 || im + 2 < daughter1) daughter1 = im + 2;
492  if (daughter2 == 0 || im + 2 > daughter2) daughter2 = im + 2;
493  }
494  }
495  // For 2 -> 1 hard scatterings reset second daughter to 0.
496  if (daughter2 == daughter1) daughter2 = 0;
497 
498  // Colour trivial, except reset irrelevant colour indices.
499  int colType = particleDataPtr->colType(id);
500  int col1 = (colType == 1 || colType == 2 || abs(colType) == 3)
501  ? lhaUpPtr->col1(iOld) : 0;
502  int col2 = (colType == -1 || colType == 2 || abs(colType) == 3)
503  ? lhaUpPtr->col2(iOld) : 0;
504 
505  // Momentum trivial.
506  double px = lhaUpPtr->px(iOld);
507  double py = lhaUpPtr->py(iOld);
508  double pz = lhaUpPtr->pz(iOld);
509  double e = lhaUpPtr->e(iOld);
510  double m = lhaUpPtr->m(iOld);
511 
512  // Polarization.
513  double pol = lhaUpPtr->spin(iOld);
514 
515  // Allow scale setting for generic partons.
516  double scaleShow = lhaUpPtr->scale(iOld);
517 
518  // For resonance decay products use resonance mass as scale.
519  double scaleNow = scalePr;
520  if (mother1 > 4) scaleNow = process[mother1].m();
521  if (scaleShow >= 0.0) scaleNow = scaleShow;
522 
523  // Store Les Houches Accord partons.
524  int iNow = process.append( id, status, mother1, mother2, daughter1,
525  daughter2, col1, col2, Vec4(px, py, pz, e), m, scaleNow, pol);
526 
527  // Check if need to store lifetime.
528  double tau = lhaUpPtr->tau(iOld);
529  if ( (setLifetime == 1 && abs(id) == 15) || setLifetime == 2)
530  tau = process[iNow].tau0() * rndmPtr->exp();
531  if (tau > 0.) process[iNow].tau(tau);
532  }
533  }
534 
535  // Loop through decay chains and set secondary vertices when needed.
536  for (int i = 3; i < process.size(); ++i) {
537  int iMother = process[i].mother1();
538 
539  // If sister to already assigned vertex then assign same.
540  if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex() )
541  process[i].vProd( process[i - 1].vProd() );
542 
543  // Else if mother already has vertex and/or lifetime then assign.
544  else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
545  process[i].vProd( process[iMother].vDec() );
546  }
547 
548  // Further info on process. Reset quantities that may or may not be known.
549  int id1Now = process[3].id();
550  int id2Now = process[4].id();
551  int id1pdf = 0;
552  int id2pdf = 0;
553  double x1pdf = 0.;
554  double x2pdf = 0.;
555  double pdf1 = 0.;
556  double pdf2 = 0.;
557  double tHat = 0.;
558  double uHat = 0.;
559  double pTHatL = 0.;
560  double m3 = 0.;
561  double m4 = 0.;
562  double theta = 0.;
563  double phi = 0.;
564  double x1Now, x2Now, Q2FacNow, alphaEM, alphaS, Q2Ren, sHat;
565 
566  // Internally generated and stored information.
567  if (!isLHA) {
568  id1pdf = id1Now;
569  id2pdf = id2Now;
570  x1Now = phaseSpacePtr->x1();
571  x2Now = phaseSpacePtr->x2();
572  x1pdf = x1Now;
573  x2pdf = x2Now;
574  pdf1 = sigmaProcessPtr->pdf1();
575  pdf2 = sigmaProcessPtr->pdf2();
576  Q2FacNow = sigmaProcessPtr->Q2Fac();
577  alphaEM = sigmaProcessPtr->alphaEMRen();
578  alphaS = sigmaProcessPtr->alphaSRen();
579  Q2Ren = sigmaProcessPtr->Q2Ren();
580  sHat = phaseSpacePtr->sHat();
581  tHat = phaseSpacePtr->tHat();
582  uHat = phaseSpacePtr->uHat();
583  pTHatL = phaseSpacePtr->pTHat();
584  m3 = phaseSpacePtr->m(3);
585  m4 = phaseSpacePtr->m(4);
586  theta = phaseSpacePtr->thetaHat();
587  phi = phaseSpacePtr->phiHat();
588  }
589 
590  // Les Houches Accord process partly available, partly to be constructed.
591  else {
592  x1Now = 2. * process[3].e() / infoPtr->eCM();
593  x2Now = 2. * process[4].e() / infoPtr->eCM();
594  Q2FacNow = (scale < 0.) ? sigmaProcessPtr->Q2Fac() : pow2(scale);
595  alphaEM = lhaUpPtr->alphaQED();
596  if (alphaEM < 0.001) alphaEM = sigmaProcessPtr->alphaEMRen();
597  alphaS = lhaUpPtr->alphaQCD();
598  if (alphaS < 0.001) alphaS = sigmaProcessPtr->alphaSRen();
599  Q2Ren = (scale < 0.) ? sigmaProcessPtr->Q2Ren() : pow2(scale);
600  Vec4 pSum = process[3].p() + process[4].p();
601  sHat = pSum * pSum;
602  int nFinLH = 0;
603  double pTLH = 0.;
604  for (int i = 5; i < process.size(); ++i)
605  if (process[i].mother1() == 3 && process[i].mother2() == 4) {
606  ++nFinLH;
607  pTLH += process[i].pT();
608  }
609  if (nFinLH > 0) pTHatL = pTLH / nFinLH;
610 
611  // Read info on parton densities if provided.
612  id1pdf = lhaUpPtr->id1pdf();
613  id2pdf = lhaUpPtr->id2pdf();
614  x1pdf = lhaUpPtr->x1pdf();
615  x2pdf = lhaUpPtr->x2pdf();
616  if (lhaUpPtr->pdfIsSet()) {
617  pdf1 = lhaUpPtr->pdf1();
618  pdf2 = lhaUpPtr->pdf2();
619  Q2FacNow = pow2(lhaUpPtr->scalePDF());
620  }
621 
622  // Reconstruct kinematics of 2 -> 2 processes from momenta.
623  if (nFin == 2) {
624  Vec4 pDifT = process[3].p() - process[5].p();
625  tHat = pDifT * pDifT;
626  Vec4 pDifU = process[3].p() - process[6].p();
627  uHat = pDifU * pDifU;
628  pTHatL = process[5].pT();
629  m3 = process[5].m();
630  m4 = process[6].m();
631  Vec4 p5 = process[5].p();
632  p5.bstback(pSum);
633  theta = p5.theta();
634  phi = process[5].phi();
635  }
636  }
637 
638  // Store information.
639  if (isHardest) {
640  infoPtr->setPDFalpha( 0, id1pdf, id2pdf, x1pdf, x2pdf, pdf1, pdf2,
641  Q2FacNow, alphaEM, alphaS, Q2Ren, scalup);
642  infoPtr->setKin( 0, id1Now, id2Now, x1Now, x2Now, sHat, tHat, uHat,
643  pTHatL, m3, m4, theta, phi);
644  }
645  infoPtr->setTypeMPI( code(), pTHatL);
646 
647  // For Les Houches event store subprocess classification.
648  if (isLHA) {
649  int codeSub = lhaUpPtr->idProcess();
650  ostringstream nameSub;
651  nameSub << "user process " << codeSub;
652  infoPtr->setSubType( 0, nameSub.str(), codeSub, nFin);
653  }
654 
655  // Done.
656  return true;
657 
658 }
659 
660 //--------------------------------------------------------------------------
661 
662 // Give the hard resonance decay chain from Les Houches input.
663 // Note: mildly modified copy of some constructProcess code; to streamline.
664 
665 bool ProcessContainer::constructDecays( Event& process) {
666 
667  // Let hard process record begin with the event as a whole.
668  process.clear();
669  process.append( 90, -11, 0, 0, 0, 0, 0, 0, Vec4(0., 0., 0., 0.), 0., 0. );
670 
671  // Since LHA partons may be out of order, determine correct one.
672  // (Recall that zeroth particle is empty.)
673  vector<int> newPos;
674  newPos.reserve(lhaUpPtr->sizePart());
675  newPos.push_back(0);
676  for (int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
677  // For iNew == 0 look for the two incoming partons, then for
678  // partons having them as mothers, and so on layer by layer.
679  for (int i = 1; i < lhaUpPtr->sizePart(); ++i)
680  if (lhaUpPtr->mother1(i) == newPos[iNew]) newPos.push_back(i);
681  if (int(newPos.size()) <= iNew) break;
682  }
683 
684  // Find scale from which to begin MPI/ISR/FSR evolution.
685  double scale = lhaUpPtr->scale();
686  process.scale( scale);
687  Vec4 pSum;
688 
689  // Copy over info from LHA event to process, in proper order.
690  for (int i = 1; i < lhaUpPtr->sizePart(); ++i) {
691  int iOld = newPos[i];
692  int id = lhaUpPtr->id(iOld);
693 
694  // Translate from LHA status codes.
695  int lhaStatus = lhaUpPtr->status(iOld);
696  int status = -21;
697  if (lhaStatus == 2 || lhaStatus == 3) status = -22;
698  if (lhaStatus == 1) status = 23;
699 
700  // Find where mothers have been moved by reordering.
701  int mother1Old = lhaUpPtr->mother1(iOld);
702  int mother2Old = lhaUpPtr->mother2(iOld);
703  int mother1 = 0;
704  int mother2 = 0;
705  for (int im = 1; im < i; ++im) {
706  if (mother1Old == newPos[im]) mother1 = im;
707  if (mother2Old == newPos[im]) mother2 = im;
708  }
709 
710  // Ensure that second mother = 0 except for bona fide carbon copies.
711  if (mother1 > 0 && mother2 == mother1) {
712  int sister1 = process[mother1].daughter1();
713  int sister2 = process[mother1].daughter2();
714  if (sister2 != sister1 && sister2 != 0) mother2 = 0;
715  }
716 
717  // Find daughters and where they have been moved by reordering.
718  int daughter1 = 0;
719  int daughter2 = 0;
720  for (int im = i + 1; im < lhaUpPtr->sizePart(); ++im) {
721  if (lhaUpPtr->mother1(newPos[im]) == iOld
722  || lhaUpPtr->mother2(newPos[im]) == iOld) {
723  if (daughter1 == 0 || im < daughter1) daughter1 = im;
724  if (daughter2 == 0 || im > daughter2) daughter2 = im;
725  }
726  }
727  // For 2 -> 1 hard scatterings reset second daughter to 0.
728  if (daughter2 == daughter1) daughter2 = 0;
729 
730  // Colour trivial, except reset irrelevant colour indices.
731  int colType = particleDataPtr->colType(id);
732  int col1 = (colType == 1 || colType == 2 || abs(colType) == 3)
733  ? lhaUpPtr->col1(iOld) : 0;
734  int col2 = (colType == -1 || colType == 2 || abs(colType) == 3)
735  ? lhaUpPtr->col2(iOld) : 0;
736 
737  // Momentum trivial.
738  double px = lhaUpPtr->px(iOld);
739  double py = lhaUpPtr->py(iOld);
740  double pz = lhaUpPtr->pz(iOld);
741  double e = lhaUpPtr->e(iOld);
742  double m = lhaUpPtr->m(iOld);
743  if (status > 0) pSum += Vec4( px, py, pz, e);
744 
745  // Polarization.
746  double pol = lhaUpPtr->spin(iOld);
747 
748  // For resonance decay products use resonance mass as scale.
749  double scaleNow = scale;
750  if (mother1 > 0) scaleNow = process[mother1].m();
751 
752  // Store Les Houches Accord partons.
753  int iNow = process.append( id, status, mother1, mother2, daughter1,
754  daughter2, col1, col2, Vec4(px, py, pz, e), m, scaleNow, pol);
755 
756  // Check if need to store lifetime.
757  double tau = lhaUpPtr->tau(iOld);
758  if ( (setLifetime == 1 && abs(id) == 15) || setLifetime == 2)
759  tau = process[iNow].tau0() * rndmPtr->exp();
760  if (tau > 0.) process[iNow].tau(tau);
761  }
762 
763  // Update four-momentum of system as a whole.
764  process[0].p( pSum);
765  process[0].m( pSum.mCalc());
766 
767  // Loop through decay chains and set secondary vertices when needed.
768  for (int i = 1; i < process.size(); ++i) {
769  int iMother = process[i].mother1();
770 
771  // If sister to already assigned vertex then assign same.
772  if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex()
773  && i > 1) process[i].vProd( process[i - 1].vProd() );
774 
775  // Else if mother already has vertex and/or lifetime then assign.
776  else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
777  process[i].vProd( process[iMother].vDec() );
778  }
779 
780  // Done.
781  return true;
782 
783 }
784 
785 //--------------------------------------------------------------------------
786 
787 // Handle resonance decays.
788 
789 bool ProcessContainer::decayResonances( Event& process) {
790 
791  // Save current event-record size and status codes.
792  process.saveSize();
793  vector<int> statusSave( process.size());
794  for (int i = 0; i < process.size(); ++i)
795  statusSave[i] = process[i].status();
796  bool physical = true;
797  bool newChain = false;
798  bool newFlavours = false;
799 
800  // Do loop over user veto.
801  do {
802 
803  // Do sequential chain of uncorrelated isotropic decays.
804  do {
805  physical = resDecaysPtr->next( process);
806  if (!physical) return false;
807 
808  // Check whether flavours should be correlated.
809  // (Currently only relevant for f fbar -> gamma*/Z0 gamma*/Z0.)
810  newFlavours = ( sigmaProcessPtr->weightDecayFlav( process)
811  < rndmPtr->flat() );
812 
813  // Reset the decay chains if have to redo.
814  if (newFlavours) {
815  process.restoreSize();
816  for (int i = 0; i < process.size(); ++i)
817  process[i].status( statusSave[i]);
818  }
819 
820  // Loop back where required to generate new decays with new flavours.
821  } while (newFlavours);
822 
823  // Correct to nonisotropic decays.
824  phaseSpacePtr->decayKinematics( process);
825 
826  // Optionally user hooks check/veto on decay chain.
827  if (canVetoResDecay)
828  newChain = userHooksPtr->doVetoResonanceDecays( process);
829 
830  // Reset the decay chains if have to redo.
831  if (newChain) {
832  process.restoreSize();
833  for (int i = 0; i < process.size(); ++i)
834  process[i].status( statusSave[i]);
835  }
836 
837  // Loop back where required to generate new decay chain.
838  } while(newChain);
839 
840  // Done.
841  return true;
842 
843 }
844 
845 //--------------------------------------------------------------------------
846 
847 // Reset event generation statistics; but NOT maximum of cross section.
848 
849 void ProcessContainer::reset() {
850 
851  nTry = 0;
852  nSel = 0;
853  nAcc = 0;
854  nTryStat = 0;
855  sigmaSum = 0.;
856  sigma2Sum = 0.;
857  sigmaNeg = 0.;
858  sigmaAvg = 0.;
859  sigmaFin = 0.;
860  deltaFin = 0.;
861  wtAccSum = 0.;
862 
863 }
864 
865 //--------------------------------------------------------------------------
866 
867 // Estimate integrated cross section and its uncertainty.
868 
869 void ProcessContainer::sigmaDelta() {
870 
871  // Initial values. No analysis meaningful unless accepted events.
872  nTryStat = nTry;
873  sigmaAvg = 0.;
874  sigmaFin = 0.;
875  deltaFin = 0.;
876  if (nAcc == 0) return;
877 
878  // Average value. No error analysis unless at least two events.
879  double nTryInv = 1. / nTry;
880  double nSelInv = 1. / nSel;
881  double nAccInv = 1. / nAcc;
882  sigmaAvg = sigmaSum * nTryInv;
883  double fracAcc = nAcc * nSelInv;
884  sigmaFin = sigmaAvg * fracAcc;
885  deltaFin = sigmaFin;
886  if (nAcc == 1) return;
887 
888  // Estimated error. Quadratic sum of cross section term and
889  // binomial from accept/reject step.
890  double delta2Sig = (lhaStratAbs != 3)
891  ? (sigma2Sum * nTryInv - pow2(sigmaAvg)) * nTryInv / pow2(sigmaAvg)
892  : pow2( lhaUpPtr->xErrSum() / lhaUpPtr->xSecSum());
893  double delta2Veto = (nSel - nAcc) * nAccInv * nSelInv;
894  double delta2Sum = delta2Sig + delta2Veto;
895  deltaFin = sqrtpos(delta2Sum) * sigmaFin;
896 
897 }
898 
899 //==========================================================================
900 
901 // SetupContainer class.
902 // Turns list of user-desired processes into a vector of containers.
903 
904 //--------------------------------------------------------------------------
905 
906 // Main routine to initialize list of processes.
907 
908 bool SetupContainers::init(vector<ProcessContainer*>& containerPtrs,
909  Info *infoPtr, Settings& settings, ParticleData* particleDataPtr,
910  Couplings* couplings) {
911 
912  // Reset process list, if filled in previous subrun.
913  if (containerPtrs.size() > 0) {
914  for (int i = 0; i < int(containerPtrs.size()); ++i)
915  delete containerPtrs[i];
916  containerPtrs.clear();
917  }
918  SigmaProcess* sigmaPtr;
919 
920  // Set up requested objects for soft QCD processes.
921  bool softQCD = settings.flag("SoftQCD:all");
922  bool inelastic = settings.flag("SoftQCD:inelastic");
923  if (softQCD || inelastic || settings.flag("SoftQCD:nonDiffractive")) {
924  sigmaPtr = new Sigma0nonDiffractive;
925  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
926  }
927  if (softQCD || settings.flag("SoftQCD:elastic")) {
928  sigmaPtr = new Sigma0AB2AB;
929  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
930  }
931  if (softQCD || inelastic || settings.flag("SoftQCD:singleDiffractive")) {
932  sigmaPtr = new Sigma0AB2XB;
933  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
934  sigmaPtr = new Sigma0AB2AX;
935  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
936  }
937  if (softQCD || inelastic || settings.flag("SoftQCD:doubleDiffractive")) {
938  sigmaPtr = new Sigma0AB2XX;
939  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
940  }
941  if (softQCD || inelastic || settings.flag("SoftQCD:centralDiffractive")) {
942  sigmaPtr = new Sigma0AB2AXB;
943  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
944  }
945 
946  // Set up requested objects for hard QCD processes.
947  bool hardQCD = settings.flag("HardQCD:all");
948  if (hardQCD || settings.flag("HardQCD:gg2gg")) {
949  sigmaPtr = new Sigma2gg2gg;
950  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
951  }
952  if (hardQCD || settings.flag("HardQCD:gg2qqbar")) {
953  sigmaPtr = new Sigma2gg2qqbar;
954  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
955  }
956  if (hardQCD || settings.flag("HardQCD:qg2qg")) {
957  sigmaPtr = new Sigma2qg2qg;
958  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
959  }
960  if (hardQCD || settings.flag("HardQCD:qq2qq")) {
961  sigmaPtr = new Sigma2qq2qq;
962  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
963  }
964  if (hardQCD || settings.flag("HardQCD:qqbar2gg")) {
965  sigmaPtr = new Sigma2qqbar2gg;
966  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
967  }
968  if (hardQCD || settings.flag("HardQCD:qqbar2qqbarNew")) {
969  sigmaPtr = new Sigma2qqbar2qqbarNew;
970  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
971  }
972 
973  // Set up requested objects for c cbar and b bbar, also hard QCD.
974  bool hardccbar = settings.flag("HardQCD:hardccbar");
975  if (hardQCD || hardccbar || settings.flag("HardQCD:gg2ccbar")) {
976  sigmaPtr = new Sigma2gg2QQbar(4, 121);
977  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
978  }
979  if (hardQCD || hardccbar || settings.flag("HardQCD:qqbar2ccbar")) {
980  sigmaPtr = new Sigma2qqbar2QQbar(4, 122);
981  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
982  }
983  bool hardbbbar = settings.flag("HardQCD:hardbbbar");
984  if (hardQCD || hardbbbar || settings.flag("HardQCD:gg2bbbar")) {
985  sigmaPtr = new Sigma2gg2QQbar(5, 123);
986  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
987  }
988  if (hardQCD || hardbbbar || settings.flag("HardQCD:qqbar2bbbar")) {
989  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
990  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
991  }
992 
993  // Set up requested objects for hard QCD 2 -> 3 processes.
994  bool hardQCD3parton = settings.flag("HardQCD:3parton");
995  if (hardQCD3parton || settings.flag("HardQCD:gg2ggg")) {
996  sigmaPtr = new Sigma3gg2ggg;
997  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
998  }
999  if (hardQCD3parton || settings.flag("HardQCD:qqbar2ggg")) {
1000  sigmaPtr = new Sigma3qqbar2ggg;
1001  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1002  }
1003  if (hardQCD3parton || settings.flag("HardQCD:qg2qgg")) {
1004  sigmaPtr = new Sigma3qg2qgg;
1005  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1006  }
1007  if (hardQCD3parton || settings.flag("HardQCD:qq2qqgDiff")) {
1008  sigmaPtr = new Sigma3qq2qqgDiff;
1009  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1010  }
1011  if (hardQCD3parton || settings.flag("HardQCD:qq2qqgSame")) {
1012  sigmaPtr = new Sigma3qq2qqgSame;
1013  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1014  }
1015  if (hardQCD3parton || settings.flag("HardQCD:qqbar2qqbargDiff")) {
1016  sigmaPtr = new Sigma3qqbar2qqbargDiff;
1017  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1018  }
1019  if (hardQCD3parton || settings.flag("HardQCD:qqbar2qqbargSame")) {
1020  sigmaPtr = new Sigma3qqbar2qqbargSame;
1021  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1022  }
1023  if (hardQCD3parton || settings.flag("HardQCD:gg2qqbarg")) {
1024  sigmaPtr = new Sigma3gg2qqbarg;
1025  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1026  }
1027  if (hardQCD3parton || settings.flag("HardQCD:qg2qqqbarDiff")) {
1028  sigmaPtr = new Sigma3qg2qqqbarDiff;
1029  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1030  }
1031  if (hardQCD3parton || settings.flag("HardQCD:qg2qqqbarSame")) {
1032  sigmaPtr = new Sigma3qg2qqqbarSame;
1033  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1034  }
1035 
1036  // Set up requested objects for prompt photon processes.
1037  bool promptPhotons = settings.flag("PromptPhoton:all");
1038  if (promptPhotons
1039  || settings.flag("PromptPhoton:qg2qgamma")) {
1040  sigmaPtr = new Sigma2qg2qgamma;
1041  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1042  }
1043  if (promptPhotons
1044  || settings.flag("PromptPhoton:qqbar2ggamma")) {
1045  sigmaPtr = new Sigma2qqbar2ggamma;
1046  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1047  }
1048  if (promptPhotons
1049  || settings.flag("PromptPhoton:gg2ggamma")) {
1050  sigmaPtr = new Sigma2gg2ggamma;
1051  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1052  }
1053  if (promptPhotons
1054  || settings.flag("PromptPhoton:ffbar2gammagamma")) {
1055  sigmaPtr = new Sigma2ffbar2gammagamma;
1056  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1057  }
1058  if (promptPhotons
1059  || settings.flag("PromptPhoton:gg2gammagamma")) {
1060  sigmaPtr = new Sigma2gg2gammagamma;
1061  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1062  }
1063 
1064  // Set up requested objects for weak gauge boson t-channel exchange.
1065  bool weakBosonExchanges = settings.flag("WeakBosonExchange:all");
1066  if (weakBosonExchanges
1067  || settings.flag("WeakBosonExchange:ff2ff(t:gmZ)")) {
1068  sigmaPtr = new Sigma2ff2fftgmZ;
1069  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1070  }
1071  if (weakBosonExchanges
1072  || settings.flag("WeakBosonExchange:ff2ff(t:W)")) {
1073  sigmaPtr = new Sigma2ff2fftW;
1074  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1075  }
1076 
1077  // Set up requested objects for weak gauge boson processes.
1078  bool weakSingleBosons = settings.flag("WeakSingleBoson:all");
1079  if (weakSingleBosons
1080  || settings.flag("WeakSingleBoson:ffbar2gmZ")) {
1081  sigmaPtr = new Sigma1ffbar2gmZ;
1082  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1083  }
1084  if (weakSingleBosons
1085  || settings.flag("WeakSingleBoson:ffbar2W")) {
1086  sigmaPtr = new Sigma1ffbar2W;
1087  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1088  }
1089 
1090  // Set up requested object for s-channel gamma exchange.
1091  // Subset of gamma*/Z0 above, intended for multiparton interactions.
1092  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:gm)")) {
1093  sigmaPtr = new Sigma2ffbar2ffbarsgm;
1094  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1095  }
1096 
1097  // Set up requested object for s-channel gamma*/Z0 or W+- exchange.
1098  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:gmZ)")) {
1099  sigmaPtr = new Sigma2ffbar2ffbarsgmZ;
1100  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1101  }
1102  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:W)")) {
1103  sigmaPtr = new Sigma2ffbar2ffbarsW;
1104  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1105  }
1106 
1107  // Set up requested objects for weak gauge boson pair processes.
1108  bool weakDoubleBosons = settings.flag("WeakDoubleBoson:all");
1109  if (weakDoubleBosons
1110  || settings.flag("WeakDoubleBoson:ffbar2gmZgmZ")) {
1111  sigmaPtr = new Sigma2ffbar2gmZgmZ;
1112  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1113  }
1114  if (weakDoubleBosons
1115  || settings.flag("WeakDoubleBoson:ffbar2ZW")) {
1116  sigmaPtr = new Sigma2ffbar2ZW;
1117  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1118  }
1119  if (weakDoubleBosons
1120  || settings.flag("WeakDoubleBoson:ffbar2WW")) {
1121  sigmaPtr = new Sigma2ffbar2WW;
1122  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1123  }
1124 
1125  // Set up requested objects for weak gauge boson + parton processes.
1126  bool weakBosonAndPartons = settings.flag("WeakBosonAndParton:all");
1127  if (weakBosonAndPartons
1128  || settings.flag("WeakBosonAndParton:qqbar2gmZg")) {
1129  sigmaPtr = new Sigma2qqbar2gmZg;
1130  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1131  }
1132  if (weakBosonAndPartons
1133  || settings.flag("WeakBosonAndParton:qg2gmZq")) {
1134  sigmaPtr = new Sigma2qg2gmZq;
1135  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1136  }
1137  if (weakBosonAndPartons
1138  || settings.flag("WeakBosonAndParton:ffbar2gmZgm")) {
1139  sigmaPtr = new Sigma2ffbar2gmZgm;
1140  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1141  }
1142  if (weakBosonAndPartons
1143  || settings.flag("WeakBosonAndParton:fgm2gmZf")) {
1144  sigmaPtr = new Sigma2fgm2gmZf;
1145  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1146  }
1147  if (weakBosonAndPartons
1148  || settings.flag("WeakBosonAndParton:qqbar2Wg")) {
1149  sigmaPtr = new Sigma2qqbar2Wg;
1150  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1151  }
1152  if (weakBosonAndPartons
1153  || settings.flag("WeakBosonAndParton:qg2Wq")) {
1154  sigmaPtr = new Sigma2qg2Wq;
1155  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1156  }
1157  if (weakBosonAndPartons
1158  || settings.flag("WeakBosonAndParton:ffbar2Wgm")) {
1159  sigmaPtr = new Sigma2ffbar2Wgm;
1160  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1161  }
1162  if (weakBosonAndPartons
1163  || settings.flag("WeakBosonAndParton:fgm2Wf")) {
1164  sigmaPtr = new Sigma2fgm2Wf;
1165  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1166  }
1167 
1168  // Set up requested objects for photon collision processes.
1169  bool photonCollisions = settings.flag("PhotonCollision:all");
1170  if (photonCollisions || settings.flag("PhotonCollision:gmgm2qqbar")) {
1171  sigmaPtr = new Sigma2gmgm2ffbar(1, 261);
1172  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1173  }
1174  if (photonCollisions || settings.flag("PhotonCollision:gmgm2ccbar")) {
1175  sigmaPtr = new Sigma2gmgm2ffbar(4, 262);
1176  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1177  }
1178  if (photonCollisions || settings.flag("PhotonCollision:gmgm2bbbar")) {
1179  sigmaPtr = new Sigma2gmgm2ffbar(5, 263);
1180  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1181  }
1182  if (photonCollisions || settings.flag("PhotonCollision:gmgm2ee")) {
1183  sigmaPtr = new Sigma2gmgm2ffbar(11, 264);
1184  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1185  }
1186  if (photonCollisions || settings.flag("PhotonCollision:gmgm2mumu")) {
1187  sigmaPtr = new Sigma2gmgm2ffbar(13, 265);
1188  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1189  }
1190  if (photonCollisions || settings.flag("PhotonCollision:gmgm2tautau")) {
1191  sigmaPtr = new Sigma2gmgm2ffbar(15, 266);
1192  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1193  }
1194 
1195  // Set up requested objects for onia production.
1196  charmonium = SigmaOniaSetup(infoPtr, &settings, particleDataPtr, 4);
1197  bottomonium = SigmaOniaSetup(infoPtr, &settings, particleDataPtr, 5);
1198  vector<SigmaProcess*> charmoniumSigmaPtrs, bottomoniumSigmaPtrs;
1199  charmonium.setupSigma2gg(charmoniumSigmaPtrs);
1200  charmonium.setupSigma2qg(charmoniumSigmaPtrs);
1201  charmonium.setupSigma2qq(charmoniumSigmaPtrs);
1202  bottomonium.setupSigma2gg(bottomoniumSigmaPtrs);
1203  bottomonium.setupSigma2qg(bottomoniumSigmaPtrs);
1204  bottomonium.setupSigma2qq(bottomoniumSigmaPtrs);
1205  for (unsigned int i = 0; i < charmoniumSigmaPtrs.size(); ++i)
1206  containerPtrs.push_back( new ProcessContainer(charmoniumSigmaPtrs[i]) );
1207  for (unsigned int i = 0; i < bottomoniumSigmaPtrs.size(); ++i)
1208  containerPtrs.push_back( new ProcessContainer(bottomoniumSigmaPtrs[i]) );
1209 
1210  // Set up requested objects for top production.
1211  bool tops = settings.flag("Top:all");
1212  if (tops || settings.flag("Top:gg2ttbar")) {
1213  sigmaPtr = new Sigma2gg2QQbar(6, 601);
1214  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1215  }
1216  if (tops || settings.flag("Top:qqbar2ttbar")) {
1217  sigmaPtr = new Sigma2qqbar2QQbar(6, 602);
1218  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1219  }
1220  if (tops || settings.flag("Top:qq2tq(t:W)")) {
1221  sigmaPtr = new Sigma2qq2QqtW(6, 603);
1222  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1223  }
1224  if (tops || settings.flag("Top:ffbar2ttbar(s:gmZ)")) {
1225  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(6, 604);
1226  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1227  }
1228  if (tops || settings.flag("Top:ffbar2tqbar(s:W)")) {
1229  sigmaPtr = new Sigma2ffbar2FfbarsW(6, 0, 605);
1230  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1231  }
1232  if (tops || settings.flag("Top:gmgm2ttbar")) {
1233  sigmaPtr = new Sigma2gmgm2ffbar(6, 606);
1234  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1235  }
1236 
1237  // Set up requested objects for fourth-generation b' production
1238  bool bPrimes = settings.flag("FourthBottom:all");
1239  if (bPrimes || settings.flag("FourthBottom:gg2bPrimebPrimebar")) {
1240  sigmaPtr = new Sigma2gg2QQbar(7, 801);
1241  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1242  }
1243  if (bPrimes || settings.flag("FourthBottom:qqbar2bPrimebPrimebar")) {
1244  sigmaPtr = new Sigma2qqbar2QQbar(7, 802);
1245  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1246  }
1247  if (bPrimes || settings.flag("FourthBottom:qq2bPrimeq(t:W)")) {
1248  sigmaPtr = new Sigma2qq2QqtW(7, 803);
1249  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1250  }
1251  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimebPrimebar(s:gmZ)")) {
1252  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(7, 804);
1253  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1254  }
1255  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimeqbar(s:W)")) {
1256  sigmaPtr = new Sigma2ffbar2FfbarsW(7, 0, 805);
1257  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1258  }
1259  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimetbar(s:W)")) {
1260  sigmaPtr = new Sigma2ffbar2FfbarsW(7, 6, 806);
1261  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1262  }
1263 
1264  // Set up requested objects for fourth-generation t' production
1265  bool tPrimes = settings.flag("FourthTop:all");
1266  if (tPrimes || settings.flag("FourthTop:gg2tPrimetPrimebar")) {
1267  sigmaPtr = new Sigma2gg2QQbar(8, 821);
1268  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1269  }
1270  if (tPrimes || settings.flag("FourthTop:qqbar2tPrimetPrimebar")) {
1271  sigmaPtr = new Sigma2qqbar2QQbar(8, 822);
1272  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1273  }
1274  if (tPrimes || settings.flag("FourthTop:qq2tPrimeq(t:W)")) {
1275  sigmaPtr = new Sigma2qq2QqtW(8, 823);
1276  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1277  }
1278  if (tPrimes || settings.flag("FourthTop:ffbar2tPrimetPrimebar(s:gmZ)")) {
1279  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(8, 824);
1280  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1281  }
1282  if (tPrimes || settings.flag("FourthTop:ffbar2tPrimeqbar(s:W)")) {
1283  sigmaPtr = new Sigma2ffbar2FfbarsW(8, 0, 825);
1284  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1285  }
1286 
1287  // Set up requested objects for two different fourth-generation fermions.
1288  if (bPrimes || tPrimes
1289  || settings.flag("FourthPair:ffbar2tPrimebPrimebar(s:W)")) {
1290  sigmaPtr = new Sigma2ffbar2FfbarsW(8, 7, 841);
1291  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1292  }
1293  if (settings.flag("FourthPair:ffbar2tauPrimenuPrimebar(s:W)")) {
1294  sigmaPtr = new Sigma2ffbar2FfbarsW(17, 18, 842);
1295  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1296  }
1297 
1298  // Flag for global choice between SM and BSM Higgses.
1299  bool useBSMHiggses = settings.flag("Higgs:useBSM");
1300 
1301  // Set up requested objects for Standard-Model Higgs production.
1302  if (!useBSMHiggses) {
1303  bool HiggsesSM = settings.flag("HiggsSM:all");
1304  if (HiggsesSM || settings.flag("HiggsSM:ffbar2H")) {
1305  sigmaPtr = new Sigma1ffbar2H(0);
1306  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1307  }
1308  if (HiggsesSM || settings.flag("HiggsSM:gg2H")) {
1309  sigmaPtr = new Sigma1gg2H(0);
1310  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1311  }
1312  if (HiggsesSM || settings.flag("HiggsSM:gmgm2H")) {
1313  sigmaPtr = new Sigma1gmgm2H(0);
1314  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1315  }
1316  if (HiggsesSM || settings.flag("HiggsSM:ffbar2HZ")) {
1317  sigmaPtr = new Sigma2ffbar2HZ(0);
1318  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1319  }
1320  if (HiggsesSM || settings.flag("HiggsSM:ffbar2HW")) {
1321  sigmaPtr = new Sigma2ffbar2HW(0);
1322  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1323  }
1324  if (HiggsesSM || settings.flag("HiggsSM:ff2Hff(t:ZZ)")) {
1325  sigmaPtr = new Sigma3ff2HfftZZ(0);
1326  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1327  }
1328  if (HiggsesSM || settings.flag("HiggsSM:ff2Hff(t:WW)")) {
1329  sigmaPtr = new Sigma3ff2HfftWW(0);
1330  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1331  }
1332  if (HiggsesSM || settings.flag("HiggsSM:gg2Httbar")) {
1333  sigmaPtr = new Sigma3gg2HQQbar(6,0);
1334  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1335  }
1336  if (HiggsesSM || settings.flag("HiggsSM:qqbar2Httbar")) {
1337  sigmaPtr = new Sigma3qqbar2HQQbar(6,0);
1338  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1339  }
1340 
1341  // Further Standard-Model Higgs processes, not included in "all".
1342  if (settings.flag("HiggsSM:qg2Hq")) {
1343  sigmaPtr = new Sigma2qg2Hq(4,0);
1344  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1345  sigmaPtr = new Sigma2qg2Hq(5,0);
1346  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1347  }
1348  if (settings.flag("HiggsSM:gg2Hbbbar")) {
1349  sigmaPtr = new Sigma3gg2HQQbar(5,0);
1350  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1351  }
1352  if (settings.flag("HiggsSM:qqbar2Hbbbar")) {
1353  sigmaPtr = new Sigma3qqbar2HQQbar(5,0);
1354  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1355  }
1356  if (settings.flag("HiggsSM:gg2Hg(l:t)")) {
1357  sigmaPtr = new Sigma2gg2Hglt(0);
1358  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1359  }
1360  if (settings.flag("HiggsSM:qg2Hq(l:t)")) {
1361  sigmaPtr = new Sigma2qg2Hqlt(0);
1362  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1363  }
1364  if (settings.flag("HiggsSM:qqbar2Hg(l:t)")) {
1365  sigmaPtr = new Sigma2qqbar2Hglt(0);
1366  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1367  }
1368  }
1369 
1370  // Common switch for the group of Higgs production BSM.
1371  if (useBSMHiggses) {
1372  bool HiggsesBSM = settings.flag("HiggsBSM:all");
1373 
1374  // Set up requested objects for BSM H1 production.
1375  bool HiggsesH1 = settings.flag("HiggsBSM:allH1");
1376  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1")) {
1377  sigmaPtr = new Sigma1ffbar2H(1);
1378  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1379  }
1380  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gg2H1")) {
1381  sigmaPtr = new Sigma1gg2H(1);
1382  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1383  }
1384  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gmgm2H1")) {
1385  sigmaPtr = new Sigma1gmgm2H(1);
1386  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1387  }
1388  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1Z")) {
1389  sigmaPtr = new Sigma2ffbar2HZ(1);
1390  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1391  }
1392  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1W")) {
1393  sigmaPtr = new Sigma2ffbar2HW(1);
1394  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1395  }
1396  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ff2H1ff(t:ZZ)")) {
1397  sigmaPtr = new Sigma3ff2HfftZZ(1);
1398  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1399  }
1400  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ff2H1ff(t:WW)")) {
1401  sigmaPtr = new Sigma3ff2HfftWW(1);
1402  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1403  }
1404  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gg2H1ttbar")) {
1405  sigmaPtr = new Sigma3gg2HQQbar(6,1);
1406  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1407  }
1408  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:qqbar2H1ttbar")) {
1409  sigmaPtr = new Sigma3qqbar2HQQbar(6,1);
1410  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1411  }
1412 
1413  // Further BSM H1 processes, not included in "all".
1414  if (settings.flag("HiggsBSM:qg2H1q")) {
1415  sigmaPtr = new Sigma2qg2Hq(4,1);
1416  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1417  sigmaPtr = new Sigma2qg2Hq(5,1);
1418  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1419  }
1420  if (settings.flag("HiggsBSM:gg2H1bbbar")) {
1421  sigmaPtr = new Sigma3gg2HQQbar(5,1);
1422  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1423  }
1424  if (settings.flag("HiggsBSM:qqbar2H1bbbar")) {
1425  sigmaPtr = new Sigma3qqbar2HQQbar(5,1);
1426  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1427  }
1428  if (settings.flag("HiggsBSM:gg2H1g(l:t)")) {
1429  sigmaPtr = new Sigma2gg2Hglt(1);
1430  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1431  }
1432  if (settings.flag("HiggsBSM:qg2H1q(l:t)")) {
1433  sigmaPtr = new Sigma2qg2Hqlt(1);
1434  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1435  }
1436  if (settings.flag("HiggsBSM:qqbar2H1g(l:t)")) {
1437  sigmaPtr = new Sigma2qqbar2Hglt(1);
1438  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1439  }
1440 
1441  // Set up requested objects for BSM H2 production.
1442  bool HiggsesH2 = settings.flag("HiggsBSM:allH2");
1443  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2")) {
1444  sigmaPtr = new Sigma1ffbar2H(2);
1445  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1446  }
1447  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gg2H2")) {
1448  sigmaPtr = new Sigma1gg2H(2);
1449  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1450  }
1451  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gmgm2H2")) {
1452  sigmaPtr = new Sigma1gmgm2H(2);
1453  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1454  }
1455  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2Z")) {
1456  sigmaPtr = new Sigma2ffbar2HZ(2);
1457  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1458  }
1459  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2W")) {
1460  sigmaPtr = new Sigma2ffbar2HW(2);
1461  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1462  }
1463  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ff2H2ff(t:ZZ)")) {
1464  sigmaPtr = new Sigma3ff2HfftZZ(2);
1465  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1466  }
1467  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ff2H2ff(t:WW)")) {
1468  sigmaPtr = new Sigma3ff2HfftWW(2);
1469  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1470  }
1471  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gg2H2ttbar")) {
1472  sigmaPtr = new Sigma3gg2HQQbar(6,2);
1473  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1474  }
1475  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:qqbar2H2ttbar")) {
1476  sigmaPtr = new Sigma3qqbar2HQQbar(6,2);
1477  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1478  }
1479 
1480  // Further BSM H2 processes, not included in "all".
1481  if (settings.flag("HiggsBSM:qg2H2q")) {
1482  sigmaPtr = new Sigma2qg2Hq(4,2);
1483  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1484  sigmaPtr = new Sigma2qg2Hq(5,2);
1485  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1486  }
1487  if (settings.flag("HiggsBSM:gg2H2bbbar")) {
1488  sigmaPtr = new Sigma3gg2HQQbar(5,2);
1489  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1490  }
1491  if (settings.flag("HiggsBSM:qqbar2H2bbbar")) {
1492  sigmaPtr = new Sigma3qqbar2HQQbar(5,2);
1493  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1494  }
1495  if (settings.flag("HiggsBSM:gg2H2g(l:t)")) {
1496  sigmaPtr = new Sigma2gg2Hglt(2);
1497  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1498  }
1499  if (settings.flag("HiggsBSM:qg2H2q(l:t)")) {
1500  sigmaPtr = new Sigma2qg2Hqlt(2);
1501  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1502  }
1503  if (settings.flag("HiggsBSM:qqbar2H2g(l:t)")) {
1504  sigmaPtr = new Sigma2qqbar2Hglt(2);
1505  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1506  }
1507 
1508  // Set up requested objects for BSM A3 production.
1509  bool HiggsesA3 = settings.flag("HiggsBSM:allA3");
1510  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3")) {
1511  sigmaPtr = new Sigma1ffbar2H(3);
1512  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1513  }
1514  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gg2A3")) {
1515  sigmaPtr = new Sigma1gg2H(3);
1516  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1517  }
1518  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gmgm2A3")) {
1519  sigmaPtr = new Sigma1gmgm2H(3);
1520  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1521  }
1522  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3Z")) {
1523  sigmaPtr = new Sigma2ffbar2HZ(3);
1524  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1525  }
1526  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3W")) {
1527  sigmaPtr = new Sigma2ffbar2HW(3);
1528  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1529  }
1530  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ff2A3ff(t:ZZ)")) {
1531  sigmaPtr = new Sigma3ff2HfftZZ(3);
1532  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1533  }
1534  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ff2A3ff(t:WW)")) {
1535  sigmaPtr = new Sigma3ff2HfftWW(3);
1536  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1537  }
1538  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gg2A3ttbar")) {
1539  sigmaPtr = new Sigma3gg2HQQbar(6,3);
1540  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1541  }
1542  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:qqbar2A3ttbar")) {
1543  sigmaPtr = new Sigma3qqbar2HQQbar(6,3);
1544  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1545  }
1546 
1547  // Further BSM A3 processes, not included in "all".
1548  if (settings.flag("HiggsBSM:qg2A3q")) {
1549  sigmaPtr = new Sigma2qg2Hq(4,3);
1550  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1551  sigmaPtr = new Sigma2qg2Hq(5,3);
1552  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1553  }
1554  if (settings.flag("HiggsBSM:gg2A3bbbar")) {
1555  sigmaPtr = new Sigma3gg2HQQbar(5,3);
1556  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1557  }
1558  if (settings.flag("HiggsBSM:qqbar2A3bbbar")) {
1559  sigmaPtr = new Sigma3qqbar2HQQbar(5,3);
1560  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1561  }
1562  if (settings.flag("HiggsBSM:gg2A3g(l:t)")) {
1563  sigmaPtr = new Sigma2gg2Hglt(3);
1564  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1565  }
1566  if (settings.flag("HiggsBSM:qg2A3q(l:t)")) {
1567  sigmaPtr = new Sigma2qg2Hqlt(3);
1568  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1569  }
1570  if (settings.flag("HiggsBSM:qqbar2A3g(l:t)")) {
1571  sigmaPtr = new Sigma2qqbar2Hglt(3);
1572  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1573  }
1574 
1575  // Set up requested objects for Charged Higgs production
1576  bool HiggsesChg = settings.flag("HiggsBSM:allH+-");
1577  if (HiggsesBSM || HiggsesChg || settings.flag("HiggsBSM:ffbar2H+-")) {
1578  sigmaPtr = new Sigma1ffbar2Hchg;
1579  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1580  }
1581  if (HiggsesBSM || HiggsesChg || settings.flag("HiggsBSM:bg2H+-t")) {
1582  sigmaPtr = new Sigma2qg2Hchgq(6, 1062, "b g -> H+- t");
1583  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1584  }
1585 
1586  // Set up requested objects for Higgs pair-production
1587  bool HiggsesPairs = settings.flag("HiggsBSM:allHpair");
1588  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2A3H1")) {
1589  sigmaPtr = new Sigma2ffbar2A3H12(1);
1590  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1591  }
1592  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2A3H2")) {
1593  sigmaPtr = new Sigma2ffbar2A3H12(2);
1594  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1595  }
1596  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+-H1")) {
1597  sigmaPtr = new Sigma2ffbar2HchgH12(1);
1598  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1599  }
1600  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+-H2")) {
1601  sigmaPtr = new Sigma2ffbar2HchgH12(2);
1602  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1603  }
1604  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+H-")) {
1605  sigmaPtr = new Sigma2ffbar2HposHneg();
1606  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1607  }
1608  }
1609 
1610  // Set up requested objects for SUSY pair processes.
1611  if (couplings->isSUSY) {
1612  CoupSUSY* coupSUSY = (CoupSUSY *) couplings;
1613 
1614  bool SUSYs = settings.flag("SUSY:all");
1615  bool nmssm = settings.flag("SLHA:NMSSM");
1616 
1617  // Preselected SUSY codes.
1618  setupIdVecs( settings);
1619 
1620  // MSSM: 4 neutralinos
1621  int nNeut = 4;
1622  if (nmssm) nNeut = 5;
1623 
1624  // Gluino-gluino
1625  if (SUSYs || settings.flag("SUSY:gg2gluinogluino")) {
1626  // Skip if outgoing codes not asked for
1627  if (allowIdVals( 1000021, 1000021)) {
1628  sigmaPtr = new Sigma2gg2gluinogluino();
1629  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1630  }
1631  }
1632  if (SUSYs || settings.flag("SUSY:qqbar2gluinogluino")) {
1633  // Skip if outgoing codes not asked for
1634  if (allowIdVals( 1000021, 1000021)) {
1635  sigmaPtr = new Sigma2qqbar2gluinogluino();
1636  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1637  }
1638  }
1639 
1640  // Gluino-squark
1641  if (SUSYs || settings.flag("SUSY:qg2squarkgluino")) {
1642  int iproc = 1202;
1643  for (int idx = 1; idx <= 6; ++idx) {
1644  for (int iso = 1; iso <= 2; ++iso) {
1645  iproc++;
1646  int id3 = iso + ((idx <= 3)
1647  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1648  int id4 = 1000021;
1649  // Skip if outgoing codes not asked for
1650  if (!allowIdVals( id3, id4)) continue;
1651  sigmaPtr = new Sigma2qg2squarkgluino(id3,iproc);
1652  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1653  }
1654  }
1655  }
1656 
1657  // Squark-antisquark (gg initiated)
1658  if (SUSYs || settings.flag("SUSY:gg2squarkantisquark")) {
1659  int iproc = 1214;
1660  for (int idx = 1; idx <= 6; ++idx) {
1661  for (int iso = 1; iso <= 2; ++iso) {
1662  iproc++;
1663  int id = iso + ((idx <= 3)
1664  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1665  // Skip if outgoing codes not asked for
1666  if (!allowIdVals( id, id)) continue;
1667  sigmaPtr = new Sigma2gg2squarkantisquark(id,iproc);
1668  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1669  }
1670  }
1671  }
1672 
1673  // Squark-antisquark (qqbar initiated)
1674  if (SUSYs || settings.flag("SUSY:qqbar2squarkantisquark")) {
1675  int iproc = 1230;
1676  for (int idx = 1; idx <= 6; ++idx) {
1677  for (int iso = 1; iso <= 2; ++iso) {
1678  for (int jso = iso; jso >= 1; --jso) {
1679  for (int jdx = 1; jdx <= 6; ++jdx) {
1680  if (iso == jso && jdx < idx) continue;
1681  int id1 = iso + ((idx <= 3) ? 1000000+2*(idx-1)
1682  : 2000000+2*(idx-4));
1683  int id2 = jso + ((jdx <= 3) ? 1000000+2*(jdx-1)
1684  : 2000000+2*(jdx-4));
1685  // Update process number counter (for ~q~q, +2 if not self-conj)
1686  //if (iproc == 1302) iproc=1310;
1687  iproc++;
1688  if (iso == jso && id1 != id2) iproc++;
1689  // Skip if outgoing codes not asked for
1690  if (!allowIdVals( id1, id2)) continue;
1691  if (iso == jso && id1 != id2) {
1692  sigmaPtr = new Sigma2qqbar2squarkantisquark(id1,-id2,iproc-1);
1693  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1694  sigmaPtr = new Sigma2qqbar2squarkantisquark(id2,-id1,iproc);
1695  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1696  } else {
1697  sigmaPtr = new Sigma2qqbar2squarkantisquark(id1,-id2,iproc);
1698  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1699  }
1700  }
1701  }
1702  }
1703  }
1704  }
1705 
1706  // Squark-squark
1707  if (SUSYs || settings.flag("SUSY:qq2squarksquark")) {
1708  int iproc = 1350;
1709  for (int idx = 1; idx <= 6; ++idx) {
1710  for (int iso = 1; iso <= 2; ++iso) {
1711  for (int jso = iso; jso >= 1; jso--) {
1712  for (int jdx = 1; jdx <= 6; ++jdx) {
1713  if (iso == jso && jdx < idx) continue;
1714  iproc++;
1715  int id1 = iso + ((idx <= 3)
1716  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1717  int id2 = jso + ((jdx <= 3)
1718  ? 1000000+2*(jdx-1) : 2000000+2*(jdx-4));
1719  // Skip if outgoing codes not asked for
1720  if (!allowIdVals( id1, id2)) continue;
1721  sigmaPtr = new Sigma2qq2squarksquark(id1,id2,iproc);
1722  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1723  }
1724  }
1725  }
1726  }
1727  }
1728 
1729  // Neutralino + squark
1730  if (SUSYs || settings.flag("SUSY:qg2chi0squark")) {
1731  int iproc = 1430;
1732  for (int iNeut= 1; iNeut <= nNeut; iNeut++) {
1733  for (int idx = 1; idx <= 6; idx++) {
1734  bool isUp = false;
1735  for (int iso = 1; iso <= 2; iso++) {
1736  if (iso == 2) isUp = true;
1737  iproc++;
1738  int id3 = coupSUSY->idNeut(iNeut);
1739  int id4 = iso + ((idx <= 3)
1740  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1741  // Skip if outgoing codes not asked for
1742  if (!allowIdVals( id3, id4)) continue;
1743  sigmaPtr = new Sigma2qg2chi0squark(iNeut,idx,isUp,iproc);
1744  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1745  }
1746  }
1747  }
1748  }
1749 
1750  // Chargino + squark
1751  if (SUSYs || settings.flag("SUSY:qg2chi+-squark")) {
1752  int iproc = 1490;
1753  for (int iChar = 1; iChar <= 2; iChar++) {
1754  for (int idx = 1; idx <= 6; idx++) {
1755  bool isUp = false;
1756  for (int iso = 1; iso <= 2; iso++) {
1757  if (iso == 2) isUp = true;
1758  iproc++;
1759  int id3 = coupSUSY->idChar(iChar);
1760  int id4 = iso + ((idx <= 3)
1761  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1762  // Skip if outgoing codes not asked for
1763  if (!allowIdVals( id3, id4)) continue;
1764  sigmaPtr = new Sigma2qg2charsquark(iChar,idx,isUp,iproc);
1765  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1766  }
1767  }
1768  }
1769  }
1770 
1771  // Neutralino pairs
1772  if (SUSYs || settings.flag("SUSY:qqbar2chi0chi0")) {
1773  int iproc = 1550;
1774  for (int iNeut2 = 1; iNeut2 <= nNeut; iNeut2++) {
1775  for (int iNeut1 = 1; iNeut1 <= iNeut2; iNeut1++) {
1776  iproc++;
1777  // Skip if outgoing codes not asked for
1778  if (!allowIdVals( coupSUSY->idNeut(iNeut1),
1779  coupSUSY->idNeut(iNeut2) ) ) continue;
1780  sigmaPtr = new Sigma2qqbar2chi0chi0(iNeut1, iNeut2,iproc);
1781  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1782  }
1783  }
1784  }
1785 
1786  // Neutralino-Chargino
1787  if (SUSYs || settings.flag("SUSY:qqbar2chi+-chi0")) {
1788  int iproc = 1570;
1789  for (int iNeut = 1; iNeut <= nNeut; iNeut++) {
1790  for (int iChar = 1; iChar <= 2; ++iChar) {
1791  iproc += 2;
1792  // Skip if outgoing codes not asked for
1793  if (!allowIdVals( coupSUSY->idNeut(iNeut),
1794  coupSUSY->idChar(iChar) ) ) continue;
1795  sigmaPtr = new Sigma2qqbar2charchi0( iChar, iNeut, iproc-1);
1796  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1797  sigmaPtr = new Sigma2qqbar2charchi0(-iChar, iNeut, iproc);
1798  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1799  }
1800  }
1801  }
1802 
1803  // Chargino-Chargino
1804  if (SUSYs || settings.flag("SUSY:qqbar2chi+chi-")) {
1805  int iproc = 1590;
1806  for (int i = 1; i <= 2; ++i) {
1807  for (int j = 1; j <= 2; ++j) {
1808  iproc++;
1809  // Skip if outgoing codes not asked for
1810  if (!allowIdVals( coupSUSY->idChar(i),
1811  coupSUSY->idChar(j) ) ) continue;
1812  sigmaPtr = new Sigma2qqbar2charchar( i,-j, iproc);
1813  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1814  }
1815  }
1816  }
1817 
1818  // RPV squark production
1819  if(SUSYs || settings.flag("SUSY:qq2antisquark")) {
1820  for (int idx = 1; idx <= 6; ++idx) {
1821  for (int iso = 1; iso <= 2; ++iso) {
1822  int id1 = iso + ((idx <= 3) ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
1823  // Skip if outgoing code not asked for
1824  if (!allowIdVals( id1, 0)) continue;
1825  sigmaPtr = new Sigma1qq2antisquark(id1);
1826  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1827  }
1828  }
1829  }
1830 
1831  // Neutralino-gluino
1832  if (SUSYs || settings.flag("SUSY:qqbar2chi0gluino")) {
1833  int iproc = 1600;
1834  for (int iNeut = 1; iNeut <= nNeut; iNeut++) {
1835  iproc++;
1836  // Skip if outgoing codes not asked for
1837  if (!allowIdVals( coupSUSY->idNeut(iNeut), 1000021)) continue;
1838  sigmaPtr = new Sigma2qqbar2chi0gluino(iNeut, iproc);
1839  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1840  }
1841  }
1842 
1843  // Chargino-Gluino
1844  if (SUSYs || settings.flag("SUSY:qqbar2chi+-gluino")) {
1845  int iproc = 1620;
1846  for (int iChar = 1; iChar <= 2; ++iChar) {
1847  iproc ++;
1848  // Skip if outgoing codes not asked for
1849  if (!allowIdVals( coupSUSY->idChar(iChar), 1000021)) continue;
1850  sigmaPtr = new Sigma2qqbar2chargluino( iChar, iproc);
1851  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1852  }
1853  }
1854 
1855  // Slepton-antislepton (qqbar initiated); Currently no RH sneutrinos
1856  if (SUSYs || settings.flag("SUSY:qqbar2sleptonantislepton")) {
1857  int iproc = 1650;
1858  for (int idx = 1; idx <= 6; ++idx) {
1859  for (int iso = 1; iso <= 2; ++iso) {
1860  for (int jso = iso; jso >= 1; --jso) {
1861  for (int jdx = 1; jdx <= 6; ++jdx) {
1862  if (iso == jso && jdx < idx) continue;
1863  int id1 = iso + ((idx <= 3) ? 1000010+2*(idx-1)
1864  : 2000010+2*(idx-4));
1865  int id2 = jso + ((jdx <= 3) ? 1000010+2*(jdx-1)
1866  : 2000010+2*(jdx-4));
1867  // Update process number counter
1868  iproc++;
1869  if (iso == jso && id1 != id2) iproc++;
1870  // Skip if outgoing codes not asked for
1871  if (!allowIdVals( id1, id2)) continue;
1872  if (iso == jso && id1 != id2) {
1873  sigmaPtr
1874  = new Sigma2qqbar2sleptonantislepton(id1,-id2,iproc-1);
1875  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1876  sigmaPtr = new Sigma2qqbar2sleptonantislepton(id2,-id1,iproc);
1877  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1878  } else {
1879  sigmaPtr = new Sigma2qqbar2sleptonantislepton(id1,-id2,iproc);
1880  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1881  }
1882  }
1883  }
1884  }
1885  }
1886  }
1887 
1888  } // End of SUSY processes.
1889 
1890  // Set up requested objects for New-Gauge-Boson processes.
1891  if (settings.flag("NewGaugeBoson:ffbar2gmZZprime")) {
1892  sigmaPtr = new Sigma1ffbar2gmZZprime();
1893  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1894  }
1895  if (settings.flag("NewGaugeBoson:ffbar2Wprime")) {
1896  sigmaPtr = new Sigma1ffbar2Wprime();
1897  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1898  }
1899  if (settings.flag("NewGaugeBoson:ffbar2R0")) {
1900  sigmaPtr = new Sigma1ffbar2Rhorizontal();
1901  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1902  }
1903 
1904  // Set up requested objects for Left-Right-Symmetry processes.
1905  bool leftrights = settings.flag("LeftRightSymmmetry:all");
1906  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2ZR")) {
1907  sigmaPtr = new Sigma1ffbar2ZRight();
1908  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1909  }
1910  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2WR")) {
1911  sigmaPtr = new Sigma1ffbar2WRight();
1912  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1913  }
1914  if (leftrights || settings.flag("LeftRightSymmmetry:ll2HL")) {
1915  sigmaPtr = new Sigma1ll2Hchgchg(1);
1916  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1917  }
1918  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLe")) {
1919  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 11);
1920  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1921  }
1922  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLmu")) {
1923  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 13);
1924  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1925  }
1926  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLtau")) {
1927  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 15);
1928  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1929  }
1930  if (leftrights || settings.flag("LeftRightSymmmetry:ff2HLff")) {
1931  sigmaPtr = new Sigma3ff2HchgchgfftWW(1);
1932  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1933  }
1934  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2HLHL")) {
1935  sigmaPtr = new Sigma2ffbar2HchgchgHchgchg(1);
1936  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1937  }
1938  if (leftrights || settings.flag("LeftRightSymmmetry:ll2HR")) {
1939  sigmaPtr = new Sigma1ll2Hchgchg(2);
1940  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1941  }
1942  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRe")) {
1943  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 11);
1944  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1945  }
1946  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRmu")) {
1947  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 13);
1948  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1949  }
1950  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRtau")) {
1951  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 15);
1952  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1953  }
1954  if (leftrights || settings.flag("LeftRightSymmmetry:ff2HRff")) {
1955  sigmaPtr = new Sigma3ff2HchgchgfftWW(2);
1956  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1957  }
1958  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2HRHR")) {
1959  sigmaPtr = new Sigma2ffbar2HchgchgHchgchg(2);
1960  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1961  }
1962 
1963  // Set up requested objects for leptoquark LQ processes.
1964  bool leptoquarks = settings.flag("LeptoQuark:all");
1965  if (leptoquarks || settings.flag("LeptoQuark:ql2LQ")) {
1966  sigmaPtr = new Sigma1ql2LeptoQuark;
1967  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1968  }
1969  if (leptoquarks || settings.flag("LeptoQuark:qg2LQl")) {
1970  sigmaPtr = new Sigma2qg2LeptoQuarkl;
1971  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1972  }
1973  if (leptoquarks || settings.flag("LeptoQuark:gg2LQLQbar")) {
1974  sigmaPtr = new Sigma2gg2LQLQbar;
1975  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1976  }
1977  if (leptoquarks || settings.flag("LeptoQuark:qqbar2LQLQbar")) {
1978  sigmaPtr = new Sigma2qqbar2LQLQbar;
1979  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1980  }
1981 
1982  // Set up requested objects for excited-fermion processes.
1983  bool excitedfermions = settings.flag("ExcitedFermion:all");
1984  if (excitedfermions || settings.flag("ExcitedFermion:dg2dStar")) {
1985  sigmaPtr = new Sigma1qg2qStar(1);
1986  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1987  }
1988  if (excitedfermions || settings.flag("ExcitedFermion:ug2uStar")) {
1989  sigmaPtr = new Sigma1qg2qStar(2);
1990  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1991  }
1992  if (excitedfermions || settings.flag("ExcitedFermion:sg2sStar")) {
1993  sigmaPtr = new Sigma1qg2qStar(3);
1994  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1995  }
1996  if (excitedfermions || settings.flag("ExcitedFermion:cg2cStar")) {
1997  sigmaPtr = new Sigma1qg2qStar(4);
1998  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1999  }
2000  if (excitedfermions || settings.flag("ExcitedFermion:bg2bStar")) {
2001  sigmaPtr = new Sigma1qg2qStar(5);
2002  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2003  }
2004  if (excitedfermions || settings.flag("ExcitedFermion:egm2eStar")) {
2005  sigmaPtr = new Sigma1lgm2lStar(11);
2006  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2007  }
2008  if (excitedfermions || settings.flag("ExcitedFermion:mugm2muStar")) {
2009  sigmaPtr = new Sigma1lgm2lStar(13);
2010  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2011  }
2012  if (excitedfermions || settings.flag("ExcitedFermion:taugm2tauStar")) {
2013  sigmaPtr = new Sigma1lgm2lStar(15);
2014  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2015  }
2016  if (excitedfermions || settings.flag("ExcitedFermion:qq2dStarq")) {
2017  sigmaPtr = new Sigma2qq2qStarq(1);
2018  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2019  }
2020  if (excitedfermions || settings.flag("ExcitedFermion:qq2uStarq")) {
2021  sigmaPtr = new Sigma2qq2qStarq(2);
2022  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2023  }
2024  if (excitedfermions || settings.flag("ExcitedFermion:qq2sStarq")) {
2025  sigmaPtr = new Sigma2qq2qStarq(3);
2026  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2027  }
2028  if (excitedfermions || settings.flag("ExcitedFermion:qq2cStarq")) {
2029  sigmaPtr = new Sigma2qq2qStarq(4);
2030  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2031  }
2032  if (excitedfermions || settings.flag("ExcitedFermion:qq2bStarq")) {
2033  sigmaPtr = new Sigma2qq2qStarq(5);
2034  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2035  }
2036  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2eStare")) {
2037  sigmaPtr = new Sigma2qqbar2lStarlbar(11);
2038  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2039  }
2040  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2nueStarnue")) {
2041  sigmaPtr = new Sigma2qqbar2lStarlbar(12);
2042  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2043  }
2044  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2muStarmu")) {
2045  sigmaPtr = new Sigma2qqbar2lStarlbar(13);
2046  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2047  }
2048  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2numuStarnumu")) {
2049  sigmaPtr = new Sigma2qqbar2lStarlbar(14);
2050  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2051  }
2052  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2tauStartau")) {
2053  sigmaPtr = new Sigma2qqbar2lStarlbar(15);
2054  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2055  }
2056  if (excitedfermions
2057  || settings.flag("ExcitedFermion:qqbar2nutauStarnutau")) {
2058  sigmaPtr = new Sigma2qqbar2lStarlbar(16);
2059  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2060  }
2061 
2062  // Set up requested objects for contact interaction processes.
2063  if (settings.flag("ContactInteractions:QCqq2qq")) {
2064  sigmaPtr = new Sigma2QCqq2qq();
2065  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2066  }
2067  if (settings.flag("ContactInteractions:QCqqbar2qqbar")) {
2068  sigmaPtr = new Sigma2QCqqbar2qqbar();
2069  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2070  }
2071  if (settings.flag("ContactInteractions:QCffbar2eebar")) {
2072  sigmaPtr = new Sigma2QCffbar2llbar(11, 4003);
2073  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2074  }
2075  if (settings.flag("ContactInteractions:QCffbar2mumubar")) {
2076  sigmaPtr = new Sigma2QCffbar2llbar(13, 4004);
2077  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2078  }
2079  if (settings.flag("ContactInteractions:QCffbar2tautaubar")) {
2080  sigmaPtr = new Sigma2QCffbar2llbar(15, 4005);
2081  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2082  }
2083 
2084  // Set spin of particles in the Hidden Valley scenario.
2085  int spinFv = settings.mode("HiddenValley:spinFv");
2086  for (int i = 1; i < 7; ++i) {
2087  if (particleDataPtr->spinType( 4900000 + i) != spinFv + 1)
2088  particleDataPtr->spinType( 4900000 + i, spinFv + 1);
2089  if (particleDataPtr->spinType( 4900010 + i) != spinFv + 1)
2090  particleDataPtr->spinType( 4900010 + i, spinFv + 1);
2091  }
2092  if (spinFv != 1) {
2093  if (particleDataPtr->spinType( 4900101) != 2)
2094  particleDataPtr->spinType( 4900101, 2);
2095  } else {
2096  int spinqv = settings.mode("HiddenValley:spinqv");
2097  if (particleDataPtr->spinType( 4900101) != 2 * spinqv + 1)
2098  particleDataPtr->spinType( 4900101, 2 * spinqv + 1);
2099  }
2100 
2101  // Set up requested objects for HiddenValley processes.
2102  bool hiddenvalleys = settings.flag("HiddenValley:all");
2103  if (hiddenvalleys || settings.flag("HiddenValley:gg2DvDvbar")) {
2104  sigmaPtr = new Sigma2gg2qGqGbar( 4900001, 4901, spinFv,
2105  "g g -> Dv Dvbar");
2106  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2107  }
2108  if (hiddenvalleys || settings.flag("HiddenValley:gg2UvUvbar")) {
2109  sigmaPtr = new Sigma2gg2qGqGbar( 4900002, 4902, spinFv,
2110  "g g -> Uv Uvbar");
2111  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2112  }
2113  if (hiddenvalleys || settings.flag("HiddenValley:gg2SvSvbar")) {
2114  sigmaPtr = new Sigma2gg2qGqGbar( 4900003, 4903, spinFv,
2115  "g g -> Sv Svbar");
2116  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2117  }
2118  if (hiddenvalleys || settings.flag("HiddenValley:gg2CvCvbar")) {
2119  sigmaPtr = new Sigma2gg2qGqGbar( 4900004, 4904, spinFv,
2120  "g g -> Cv Cvbar");
2121  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2122  }
2123  if (hiddenvalleys || settings.flag("HiddenValley:gg2BvBvbar")) {
2124  sigmaPtr = new Sigma2gg2qGqGbar( 4900005, 4905, spinFv,
2125  "g g -> Bv Bvbar");
2126  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2127  }
2128  if (hiddenvalleys || settings.flag("HiddenValley:gg2TvTvbar")) {
2129  sigmaPtr = new Sigma2gg2qGqGbar( 4900006, 4906, spinFv,
2130  "g g -> Tv Tvbar");
2131  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2132  }
2133  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2DvDvbar")) {
2134  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900001, 4911, spinFv,
2135  "q qbar -> Dv Dvbar");
2136  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2137  }
2138  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2UvUvbar")) {
2139  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900002, 4912, spinFv,
2140  "q qbar -> Uv Uvbar");
2141  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2142  }
2143  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2SvSvbar")) {
2144  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900003, 4913, spinFv,
2145  "q qbar -> Sv Svbar");
2146  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2147  }
2148  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2CvCvbar")) {
2149  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900004, 4914, spinFv,
2150  "q qbar -> Cv Cvbar");
2151  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2152  }
2153  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2BvBvbar")) {
2154  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900005, 4915, spinFv,
2155  "q qbar -> Bv Bvbar");
2156  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2157  }
2158  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2TvTvbar")) {
2159  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900006, 4916, spinFv,
2160  "q qbar -> Tv Tvbar");
2161  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2162  }
2163  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2DvDvbar")) {
2164  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900001, 4921, spinFv,
2165  "f fbar -> Dv Dvbar");
2166  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2167  }
2168  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2UvUvbar")) {
2169  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900002, 4922, spinFv,
2170  "f fbar -> Uv Uvbar");
2171  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2172  }
2173  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2SvSvbar")) {
2174  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900003, 4923, spinFv,
2175  "f fbar -> Sv Svbar");
2176  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2177  }
2178  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2CvCvbar")) {
2179  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900004, 4924, spinFv,
2180  "f fbar -> Cv Cvbar");
2181  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2182  }
2183  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2BvBvbar")) {
2184  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900005, 4925, spinFv,
2185  "f fbar -> Bv Bvbar");
2186  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2187  }
2188  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2TvTvbar")) {
2189  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900006, 4926, spinFv,
2190  "f fbar -> Tv Tvbar");
2191  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2192  }
2193  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2EvEvbar")) {
2194  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900011, 4931, spinFv,
2195  "f fbar -> Ev Evbar");
2196  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2197  }
2198  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuEvnuEvbar")) {
2199  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900012, 4932, spinFv,
2200  "f fbar -> nuEv nuEvbar");
2201  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2202  }
2203  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2MUvMUvbar")) {
2204  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900013, 4933, spinFv,
2205  "f fbar -> MUv MUvbar");
2206  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2207  }
2208  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuMUvnuMUvbar")) {
2209  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900014, 4934, spinFv,
2210  "f fbar -> nuMUv nuMUvbar");
2211  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2212  }
2213  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2TAUvTAUvbar")) {
2214  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900015, 4935, spinFv,
2215  "f fbar -> TAUv TAUvbar");
2216  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2217  }
2218  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuTAUvnuTAUvbar")) {
2219  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900016, 4936, spinFv,
2220  "f fbar -> nuTAUv nuTAUvbar");
2221  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2222  }
2223  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2Zv")) {
2224  sigmaPtr = new Sigma1ffbar2Zv();
2225  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2226  }
2227 
2228  // Set up requested objects for RS extra-dimensional G* processes.
2229  bool extraDimGstars = settings.flag("ExtraDimensionsG*:all");
2230  if (extraDimGstars || settings.flag("ExtraDimensionsG*:gg2G*")) {
2231  sigmaPtr = new Sigma1gg2GravitonStar;
2232  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2233  }
2234  if (extraDimGstars || settings.flag("ExtraDimensionsG*:ffbar2G*")) {
2235  sigmaPtr = new Sigma1ffbar2GravitonStar;
2236  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2237  }
2238  if (settings.flag("ExtraDimensionsG*:gg2G*g")) {
2239  sigmaPtr = new Sigma2gg2GravitonStarg;
2240  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2241  }
2242  if (settings.flag("ExtraDimensionsG*:qg2G*q")) {
2243  sigmaPtr = new Sigma2qg2GravitonStarq;
2244  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2245  }
2246  if (settings.flag("ExtraDimensionsG*:qqbar2G*g")) {
2247  sigmaPtr = new Sigma2qqbar2GravitonStarg;
2248  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2249  }
2250 
2251  // Set up requested objects for RS extra-dimensional KKgluon processes.
2252  if (settings.flag("ExtraDimensionsG*:qqbar2KKgluon*")) {
2253  sigmaPtr = new Sigma1qqbar2KKgluonStar;
2254  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2255  }
2256 
2257  // NOAM: Set up requested objects for TEV extra-dimensional processes.
2258  if (settings.flag("ExtraDimensionsTEV:ffbar2ddbar")) {
2259  sigmaPtr = new Sigma2ffbar2TEVffbar(1, 5061);
2260  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2261  }
2262  if (settings.flag("ExtraDimensionsTEV:ffbar2uubar")) {
2263  sigmaPtr = new Sigma2ffbar2TEVffbar(2, 5062);
2264  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2265  }
2266  if (settings.flag("ExtraDimensionsTEV:ffbar2ssbar")) {
2267  sigmaPtr = new Sigma2ffbar2TEVffbar(3, 5063);
2268  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2269  }
2270  if (settings.flag("ExtraDimensionsTEV:ffbar2ccbar")) {
2271  sigmaPtr = new Sigma2ffbar2TEVffbar(4, 5064);
2272  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2273  }
2274  if (settings.flag("ExtraDimensionsTEV:ffbar2bbbar")) {
2275  sigmaPtr = new Sigma2ffbar2TEVffbar(5, 5065);
2276  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2277  }
2278  if (settings.flag("ExtraDimensionsTEV:ffbar2ttbar")) {
2279  sigmaPtr = new Sigma2ffbar2TEVffbar(6, 5066);
2280  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2281  }
2282  if (settings.flag("ExtraDimensionsTEV:ffbar2e+e-")) {
2283  sigmaPtr = new Sigma2ffbar2TEVffbar(11, 5071);
2284  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2285  }
2286  if (settings.flag("ExtraDimensionsTEV:ffbar2nuenuebar")) {
2287  sigmaPtr = new Sigma2ffbar2TEVffbar(12, 5072);
2288  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2289  }
2290  if (settings.flag("ExtraDimensionsTEV:ffbar2mu+mu-")) {
2291  sigmaPtr = new Sigma2ffbar2TEVffbar(13, 5073);
2292  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2293  }
2294  if (settings.flag("ExtraDimensionsTEV:ffbar2numunumubar")) {
2295  sigmaPtr = new Sigma2ffbar2TEVffbar(14, 5074);
2296  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2297  }
2298  if (settings.flag("ExtraDimensionsTEV:ffbar2tau+tau-")) {
2299  sigmaPtr = new Sigma2ffbar2TEVffbar(15, 5075);
2300  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2301  }
2302  if (settings.flag("ExtraDimensionsTEV:ffbar2nutaunutaubar")) {
2303  sigmaPtr = new Sigma2ffbar2TEVffbar(16, 5076);
2304  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2305  }
2306 
2307  // Set up requested objects for large extra-dimensional G processes.
2308  bool extraDimLEDmono = settings.flag("ExtraDimensionsLED:monojet");
2309  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:gg2Gg")) {
2310  sigmaPtr = new Sigma2gg2LEDUnparticleg( true );
2311  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2312  }
2313  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:qg2Gq")) {
2314  sigmaPtr = new Sigma2qg2LEDUnparticleq( true );
2315  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2316  }
2317  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:qqbar2Gg")) {
2318  sigmaPtr = new Sigma2qqbar2LEDUnparticleg( true );
2319  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2320  }
2321  if (settings.flag("ExtraDimensionsLED:ffbar2GZ")) {
2322  sigmaPtr = new Sigma2ffbar2LEDUnparticleZ( true );
2323  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2324  }
2325  if (settings.flag("ExtraDimensionsLED:ffbar2Ggamma")) {
2326  sigmaPtr = new Sigma2ffbar2LEDUnparticlegamma( true );
2327  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2328  }
2329  if (settings.flag("ExtraDimensionsLED:ffbar2gammagamma")) {
2330  sigmaPtr = new Sigma2ffbar2LEDgammagamma( true );
2331  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2332  }
2333  if (settings.flag("ExtraDimensionsLED:gg2gammagamma")) {
2334  sigmaPtr = new Sigma2gg2LEDgammagamma( true );
2335  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2336  }
2337  if (settings.flag("ExtraDimensionsLED:ffbar2llbar")) {
2338  sigmaPtr = new Sigma2ffbar2LEDllbar( true );
2339  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2340  }
2341  if (settings.flag("ExtraDimensionsLED:gg2llbar")) {
2342  sigmaPtr = new Sigma2gg2LEDllbar( true );
2343  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2344  }
2345  bool extraDimLEDdij = settings.flag("ExtraDimensionsLED:dijets");
2346  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:gg2DJgg")) {
2347  sigmaPtr = new Sigma2gg2LEDgg;
2348  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2349  }
2350  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:gg2DJqqbar")) {
2351  sigmaPtr = new Sigma2gg2LEDqqbar;
2352  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2353  }
2354  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qg2DJqg")) {
2355  sigmaPtr = new Sigma2qg2LEDqg;
2356  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2357  }
2358  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qq2DJqq")) {
2359  sigmaPtr = new Sigma2qq2LEDqq;
2360  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2361  }
2362  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qqbar2DJgg")) {
2363  sigmaPtr = new Sigma2qqbar2LEDgg;
2364  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2365  }
2366  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qqbar2DJqqbarNew")) {
2367  sigmaPtr = new Sigma2qqbar2LEDqqbarNew;
2368  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2369  }
2370 
2371  // Set up requested objects for unparticle processes.
2372  bool extraDimUnpartmono = settings.flag("ExtraDimensionsUnpart:monojet");
2373  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:gg2Ug")) {
2374  sigmaPtr = new Sigma2gg2LEDUnparticleg( false );
2375  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2376  }
2377  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:qg2Uq")) {
2378  sigmaPtr = new Sigma2qg2LEDUnparticleq( false );
2379  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2380  }
2381  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:qqbar2Ug")) {
2382  sigmaPtr = new Sigma2qqbar2LEDUnparticleg( false );
2383  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2384  }
2385  if (settings.flag("ExtraDimensionsUnpart:ffbar2UZ")) {
2386  sigmaPtr = new Sigma2ffbar2LEDUnparticleZ( false );
2387  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2388  }
2389  if (settings.flag("ExtraDimensionsUnpart:ffbar2Ugamma")) {
2390  sigmaPtr = new Sigma2ffbar2LEDUnparticlegamma( false );
2391  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2392  }
2393  if (settings.flag("ExtraDimensionsUnpart:ffbar2gammagamma")) {
2394  sigmaPtr = new Sigma2ffbar2LEDgammagamma( false );
2395  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2396  }
2397  if (settings.flag("ExtraDimensionsUnpart:gg2gammagamma")) {
2398  sigmaPtr = new Sigma2gg2LEDgammagamma( false );
2399  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2400  }
2401  if (settings.flag("ExtraDimensionsUnpart:ffbar2llbar")) {
2402  sigmaPtr = new Sigma2ffbar2LEDllbar( false );
2403  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2404  }
2405  if (settings.flag("ExtraDimensionsUnpart:gg2llbar")) {
2406  sigmaPtr = new Sigma2gg2LEDllbar( false );
2407  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2408  }
2409 
2410  // Done.
2411  return true;
2412 
2413 }
2414 
2415 //--------------------------------------------------------------------------
2416 
2417 // Routine to initialize list of second hard processes.
2418 
2419 bool SetupContainers::init2(vector<ProcessContainer*>& container2Ptrs,
2420  Settings& settings) {
2421 
2422  // Reset process list, if filled in previous subrun.
2423  if (container2Ptrs.size() > 0) {
2424  for (int i = 0; i < int(container2Ptrs.size()); ++i)
2425  delete container2Ptrs[i];
2426  container2Ptrs.clear();
2427  }
2428  SigmaProcess* sigmaPtr;
2429 
2430  // Two hard QCD jets.
2431  if (settings.flag("SecondHard:TwoJets")) {
2432  sigmaPtr = new Sigma2gg2gg;
2433  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2434  sigmaPtr = new Sigma2gg2qqbar;
2435  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2436  sigmaPtr = new Sigma2qg2qg;
2437  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2438  sigmaPtr = new Sigma2qq2qq;
2439  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2440  sigmaPtr = new Sigma2qqbar2gg;
2441  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2442  sigmaPtr = new Sigma2qqbar2qqbarNew;
2443  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2444  sigmaPtr = new Sigma2gg2QQbar(4, 121);
2445  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2446  sigmaPtr = new Sigma2qqbar2QQbar(4, 122);
2447  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2448  sigmaPtr = new Sigma2gg2QQbar(5, 123);
2449  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2450  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
2451  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2452  }
2453 
2454  // A prompt photon and a hard jet.
2455  if (settings.flag("SecondHard:PhotonAndJet")) {
2456  sigmaPtr = new Sigma2qg2qgamma;
2457  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2458  sigmaPtr = new Sigma2qqbar2ggamma;
2459  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2460  sigmaPtr = new Sigma2gg2ggamma;
2461  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2462  }
2463 
2464  // Two prompt photons.
2465  if (settings.flag("SecondHard:TwoPhotons")) {
2466  sigmaPtr = new Sigma2ffbar2gammagamma;
2467  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2468  sigmaPtr = new Sigma2gg2gammagamma;
2469  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2470  }
2471 
2472  // Charmonium.
2473  if (settings.flag("SecondHard:Charmonium")) {
2474  vector<SigmaProcess*> charmoniumSigmaPtrs;
2475  charmonium.setupSigma2gg(charmoniumSigmaPtrs, true);
2476  charmonium.setupSigma2qg(charmoniumSigmaPtrs, true);
2477  charmonium.setupSigma2qq(charmoniumSigmaPtrs, true);
2478  for (unsigned int i = 0; i < charmoniumSigmaPtrs.size(); ++i)
2479  container2Ptrs.push_back( new ProcessContainer(charmoniumSigmaPtrs[i]) );
2480  }
2481 
2482  // Bottomonium.
2483  if (settings.flag("SecondHard:Bottomonium")) {
2484  vector<SigmaProcess*> bottomoniumSigmaPtrs;
2485  bottomonium.setupSigma2gg(bottomoniumSigmaPtrs, true);
2486  bottomonium.setupSigma2qg(bottomoniumSigmaPtrs, true);
2487  bottomonium.setupSigma2qq(bottomoniumSigmaPtrs, true);
2488  for (unsigned int i = 0; i < bottomoniumSigmaPtrs.size(); ++i)
2489  container2Ptrs.push_back( new ProcessContainer(bottomoniumSigmaPtrs[i]) );
2490  }
2491 
2492  // A single gamma*/Z0.
2493  if (settings.flag("SecondHard:SingleGmZ")) {
2494  sigmaPtr = new Sigma1ffbar2gmZ;
2495  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2496  }
2497 
2498  // A single W+-.
2499  if (settings.flag("SecondHard:SingleW")) {
2500  sigmaPtr = new Sigma1ffbar2W;
2501  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2502  }
2503 
2504  // A gamma*/Z0 and a hard jet.
2505  if (settings.flag("SecondHard:GmZAndJet")) {
2506  sigmaPtr = new Sigma2qqbar2gmZg;
2507  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2508  sigmaPtr = new Sigma2qg2gmZq;
2509  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2510  }
2511 
2512  // A W+- and a hard jet.
2513  if (settings.flag("SecondHard:WAndJet")) {
2514  sigmaPtr = new Sigma2qqbar2Wg;
2515  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2516  sigmaPtr = new Sigma2qg2Wq;
2517  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2518  }
2519 
2520  // Top pair production.
2521  if (settings.flag("SecondHard:TopPair")) {
2522  sigmaPtr = new Sigma2gg2QQbar(6, 601);
2523  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2524  sigmaPtr = new Sigma2qqbar2QQbar(6, 602);
2525  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2526  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(6, 604);
2527  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2528  }
2529 
2530  // Single top production.
2531  if (settings.flag("SecondHard:SingleTop")) {
2532  sigmaPtr = new Sigma2qq2QqtW(6, 603);
2533  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2534  sigmaPtr = new Sigma2ffbar2FfbarsW(6, 0, 605);
2535  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2536  }
2537 
2538  // Two b jets - already part of TwoJets sample above.
2539  if (settings.flag("SecondHard:TwoBJets")) {
2540  sigmaPtr = new Sigma2gg2QQbar(5, 123);
2541  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2542  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
2543  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
2544  }
2545 
2546  // Done.
2547  return true;
2548 
2549 }
2550 
2551 //--------------------------------------------------------------------------
2552 
2553 // Set up arrays of allowed outgoing SUSY particles.
2554 
2555 void SetupContainers::setupIdVecs( Settings& settings) {
2556 
2557  // First array either none, one or many particles.
2558  idVecA.resize(0);
2559  if (settings.mode("SUSY:idA") != 0) {
2560  idVecA.push_back( abs(settings.mode("SUSY:idA")) );
2561  } else {
2562  vector<int> idTmpA = settings.mvec("SUSY:idVecA");
2563  for (int i = 0; i < int(idTmpA.size()); ++i)
2564  if (idTmpA[i] != 0) idVecA.push_back( abs(idTmpA[i]) );
2565  }
2566  nVecA = idVecA.size();
2567 
2568  // Second array either none, one or many particles.
2569  idVecB.resize(0);
2570  if (settings.mode("SUSY:idB") != 0) {
2571  idVecB.push_back( abs(settings.mode("SUSY:idB")) );
2572  } else {
2573  vector<int> idTmpB = settings.mvec("SUSY:idVecB");
2574  for (int i = 0; i < int(idTmpB.size()); ++i)
2575  if (idTmpB[i] != 0) idVecB.push_back( abs(idTmpB[i]) );
2576  }
2577  nVecB = idVecB.size();
2578 
2579 }
2580 
2581 //--------------------------------------------------------------------------
2582 
2583 // Check final state for allowed outgoing SUSY particles.
2584 // Normally check two codes, but allow for only one.
2585 
2586 bool SetupContainers::allowIdVals( int idCheck1, int idCheck2) {
2587 
2588  // If empty arrays or id's no need for checks. Else need absolute values.
2589  if (nVecA == 0 && nVecB == 0) return true;
2590  if (idCheck1 == 0 && idCheck2 == 0) return true;
2591  int idChk1 = abs(idCheck1);
2592  int idChk2 = abs(idCheck2);
2593 
2594  // If only one outgoing particle then check idVecA and idVecB.
2595  if (idChk1 == 0) swap(idChk1, idChk2);
2596  if (idChk2 == 0) {
2597  for (int i = 0; i < nVecA; ++i) if (idChk1 == idVecA[i]) return true;
2598  for (int i = 0; i < nVecB; ++i) if (idChk1 == idVecB[i]) return true;
2599  return false;
2600  }
2601 
2602  // If empty array idVecB then compare with idVecA.
2603  if (nVecB == 0) {
2604  for (int i = 0; i < nVecA; ++i)
2605  if (idChk1 == idVecA[i] || idChk2 == idVecA[i]) return true;
2606  return false;
2607  }
2608 
2609  // If empty array idVecA then compare with idVecB.
2610  if (nVecA == 0) {
2611  for (int i = 0; i < nVecB; ++i)
2612  if (idChk1 == idVecB[i] || idChk2 == idVecB[i]) return true;
2613  return false;
2614  }
2615 
2616  // Else check that pair matches allowed combinations.
2617  for (int i = 0; i < nVecA; ++i)
2618  for (int j = 0; j < nVecB; ++j)
2619  if ( (idChk1 == idVecA[i] && idChk2 == idVecB[j])
2620  || (idChk2 == idVecA[i] && idChk1 == idVecB[j]) ) return true;
2621  return false;
2622 
2623 }
2624 
2625 //==========================================================================
2626 
2627 } // end namespace Pythia8
Definition: AgUStep.h:26