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) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // 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 #include "Pythia8/SigmaDM.h"
23 #include "Pythia8/SusyCouplings.h"
24 
25 namespace Pythia8 {
26 
27 //==========================================================================
28 
29 // ProcessContainer class.
30 // Information allowing the generation of a specific process.
31 
32 //--------------------------------------------------------------------------
33 
34 // Constants: could be changed here if desired, but normally should not.
35 // These are of technical nature, as described for each.
36 
37 // Number of event tries to check maximization finding reliability.
38 const int ProcessContainer::N12SAMPLE = 100;
39 
40 // Ditto, but increased for 2 -> 3 processes.
41 const int ProcessContainer::N3SAMPLE = 1000;
42 
43 //--------------------------------------------------------------------------
44 
45 // Initialize phase space and counters.
46 // Argument isFirst distinguishes two hard processes in same event.
47 
48 bool ProcessContainer::init(bool isFirst, ResonanceDecays* resDecaysPtrIn,
49  SLHAinterface* slhaInterfacePtr, GammaKinematics* gammaKinPtrIn) {
50 
51  registerSubObject(*sigmaProcessPtr);
52 
53  // Extract info about current process from SigmaProcess object.
54  isLHA = sigmaProcessPtr->isLHA();
55  isNonDiff = sigmaProcessPtr->isNonDiff();
56  isResolved = sigmaProcessPtr->isResolved();
57  isDiffA = sigmaProcessPtr->isDiffA();
58  isDiffB = sigmaProcessPtr->isDiffB();
59  isDiffC = sigmaProcessPtr->isDiffC();
60  isQCD3body = sigmaProcessPtr->isQCD3body();
61  int nFin = sigmaProcessPtr->nFinal();
62  lhaStrat = (isLHA) ? lhaUpPtr->strategy() : 0;
63  lhaStratAbs = abs(lhaStrat);
64  allowNegSig = sigmaProcessPtr->allowNegativeSigma();
65  processCode = sigmaProcessPtr->code();
66 
67  // Switch to ensure that the scales set in a LH input event are
68  // respected, even in resonance showers.
69  useStrictLHEFscales = flag("Beams:strictLHEFscale");
70 
71  // Maximal number of events for more in-depth control of
72  // how many events are selected/tried/accepted.
73  nTryRequested = mode("Main:numberOfTriedEvents");
74  nSelRequested = mode("Main:numberOfSelectedEvents");
75  nAccRequested = mode("Main:numberOfAcceptedEvents");
76 
77  // Flag for maximum violation handling.
78  increaseMaximum = flag("PhaseSpace:increaseMaximum");
79 
80  // Store whether beam particle has a photon and save the mode.
81  gammaKinPtr = gammaKinPtrIn;
82 
83  // Use external photon flux.
84  approximatedGammaFlux = beamAPtr->hasApproxGammaFlux() ||
85  beamBPtr->hasApproxGammaFlux();
86 
87  // Check whether photon sub-beams present.
88  bool beamAhasGamma = flag("PDF:beamA2gamma");
89  bool beamBhasGamma = flag("PDF:beamB2gamma");
90  beamHasGamma = beamAhasGamma || beamBhasGamma;
91 
92  // Pick and create phase space generator. Send pointers where required.
93  if (phaseSpacePtr != 0) ;
94  else if (isLHA) phaseSpacePtr = new PhaseSpaceLHA();
95  else if (isNonDiff) phaseSpacePtr = new PhaseSpace2to2nondiffractive();
96  else if (!isResolved && !isDiffA && !isDiffB && !isDiffC )
97  phaseSpacePtr = new PhaseSpace2to2elastic();
98  else if (!isResolved && !isDiffA && !isDiffB && isDiffC)
99  phaseSpacePtr = new PhaseSpace2to3diffractive();
100  else if (!isResolved) phaseSpacePtr = new PhaseSpace2to2diffractive(
101  isDiffA, isDiffB);
102  else if (nFin == 1) phaseSpacePtr = new PhaseSpace2to1tauy();
103  else if (nFin == 2) phaseSpacePtr = new PhaseSpace2to2tauyz();
104  else if (isQCD3body) phaseSpacePtr = new PhaseSpace2to3yyycyl();
105  else phaseSpacePtr = new PhaseSpace2to3tauycyl();
106 
107  // Store pointers and perform simple initialization.
108  resDecaysPtr = resDecaysPtrIn;
109  canVetoResDecay = (userHooksPtr != 0)
110  ? userHooksPtr->canVetoResonanceDecays() : false;
111  if (isLHA) {
112  sigmaProcessPtr->setLHAPtr(lhaUpPtr);
113  phaseSpacePtr->setLHAPtr(lhaUpPtr);
114 
115  // Check if beams asymmetric and define a boost to CM frame if yes.
116  double relEnergyDiff = ( lhaUpPtr->eBeamA() - lhaUpPtr->eBeamB() )
117  / ( lhaUpPtr->eBeamA() + lhaUpPtr->eBeamB() );
118  isAsymLHA = ( abs(relEnergyDiff) > 1.e-10 )
119  || ( lhaUpPtr->idBeamA() != lhaUpPtr->idBeamB() );
120  if (isAsymLHA) {
121  betazLHA = (sqrtpos(pow2(lhaUpPtr->eBeamA()) - pow2(infoPtr->mA()))
122  - sqrtpos(pow2(lhaUpPtr->eBeamB()) - pow2(infoPtr->mB())))
123  / (lhaUpPtr->eBeamA() + lhaUpPtr->eBeamB());
124  }
125  }
126  sigmaProcessPtr->init(beamAPtr, beamBPtr, slhaInterfacePtr);
127 
128  // Store the state of photon beams using inFlux: 0 = not a photon beam;
129  // 1 = resolved photon; 2 = unresolved photon.
130  string inState = sigmaProcessPtr->inFlux();
131  beamAgammaMode = 0;
132  beamBgammaMode = 0;
133  gammaModeEvent = 0;
134  if ( beamAPtr->isGamma() || beamAhasGamma ) {
135  if ( inState == "gmg" || inState == "gmq" || inState == "gmgm" )
136  beamAgammaMode = 2;
137  else beamAgammaMode = 1;
138  }
139  if ( beamBPtr->isGamma() || beamBhasGamma ) {
140  if ( inState == "ggm" || inState == "qgm" || inState == "gmgm" )
141  beamBgammaMode = 2;
142  else beamBgammaMode = 1;
143  }
144 
145  // Save the photon modes and propagate further.
146  // For soft processes set possible VM states.
147  if ( ( beamAPtr->isGamma() || beamBPtr->isGamma() ) || beamHasGamma )
148  setBeamModes(isSoftQCD(), false);
149 
150  // Save the state of photon beams.
151  beamAhasResGamma = beamAPtr->hasResGamma();
152  beamBhasResGamma = beamBPtr->hasResGamma();
153  beamHasResGamma = beamAhasResGamma || beamBhasResGamma;
154 
155  // Initialize also phaseSpace pointer.
156  if (phaseSpacePtr) registerSubObject(*phaseSpacePtr);
157  phaseSpacePtr->init( isFirst, sigmaProcessPtr);
158 
159  // Send the pointer to gammaKinematics for sampling of soft QCD processes.
160  if ( beamAhasResGamma || beamBhasResGamma )
161  phaseSpacePtr->setGammaKinPtr( gammaKinPtr);
162 
163  // Reset cross section statistics.
164  nTry = 0;
165  nSel = 0;
166  nAcc = 0;
167  nTryStat = 0;
168  sigmaMx = 0.;
169  sigmaSum = 0.;
170  sigma2Sum = 0.;
171  sigmaNeg = 0.;
172  sigmaAvg = 0.;
173  sigmaFin = 0.;
174  deltaFin = 0.;
175  wtAccSum = 0.;
176 
177  // Reset temporary cross section summand.
178  sigmaTemp = 0.;
179  sigma2Temp = 0.;
180 
181  // Set up normalized variance for lhaStratAbs = 3.
182  normVar3 = (lhaStratAbs == 3)
183  ? pow2( lhaUpPtr->xErrSum() / lhaUpPtr->xSecSum()) : 0.;
184 
185  // Initialize process and allowed incoming partons.
186  sigmaProcessPtr->initProc();
187  if (!sigmaProcessPtr->initFlux()) return false;
188 
189  // Find maximum of differential cross section * phasespace.
190  bool physical = phaseSpacePtr->setupSampling();
191  sigmaMx = phaseSpacePtr->sigmaMax();
192  double sigmaHalfWay = sigmaMx;
193 
194  // Separate signed maximum needed for LHA with negative weight.
195  sigmaSgn = phaseSpacePtr->sigmaSumSigned();
196 
197  // Check maximum by a few events, and extrapolate a further increase.
198  if (physical & !isLHA) {
199  int nSample = (nFin < 3) ? N12SAMPLE : N3SAMPLE;
200  for (int iSample = 0; iSample < nSample; ++iSample) {
201  bool test = false;
202  while (!test) test = phaseSpacePtr->trialKin(false);
203  if (iSample == nSample/2) sigmaHalfWay = phaseSpacePtr->sigmaMax();
204  }
205  double sigmaFullWay = phaseSpacePtr->sigmaMax();
206  sigmaMx = (sigmaHalfWay > 0.) ? pow2(sigmaFullWay) / sigmaHalfWay
207  : sigmaFullWay;
208  phaseSpacePtr->setSigmaMax(sigmaMx);
209  }
210 
211  // Allow Pythia to overwrite incoming beams or parts of Les Houches input.
212  idRenameBeams = mode("LesHouches:idRenameBeams");
213  setLifetime = mode("LesHouches:setLifetime");
214  setQuarkMass = mode("LesHouches:setQuarkMass");
215  setLeptonMass = mode("LesHouches:setLeptonMass");
216  mRecalculate = parm("LesHouches:mRecalculate");
217  matchInOut = flag("LesHouches:matchInOut");
218  for (int i = 0; i < 6; ++i) idNewM[i] = i;
219  for (int i = 6; i < 9; ++i) idNewM[i] = 2 * i - 1;
220  for (int i = 1; i < 9; ++i) mNewM[i] = particleDataPtr->m0(idNewM[i]);
221 
222  // Done.
223  return physical;
224 }
225 
226 //--------------------------------------------------------------------------
227 
228 // Generate a trial event; selected or not.
229 
230 bool ProcessContainer::trialProcess() {
231 
232  // Weights for photon flux oversampling with external flux.
233  double wtPDF = 1.;
234  double wtFlux = 1.;
235 
236  // For photon beams set PDF pointer to resolved or unresolved.
237  // Do not reset VM states since they will be sampled in phaseSpace->trialKin.
238  if ( beamAPtr->isGamma() || beamBPtr->isGamma() || beamHasGamma )
239  setBeamModes(false);
240 
241  // Loop over tries only occurs for Les Houches strategy = +-2.
242  for (int iTry = 0; ; ++iTry) {
243 
244  // Generate a trial phase space point, if meaningful.
245  if (sigmaMx == 0.) return false;
246  infoPtr->setEndOfFile(false);
247  bool repeatSame = (iTry > 0);
248  bool physical = phaseSpacePtr->trialKin(true, repeatSame);
249 
250  // For acceptable kinematics sample also kT for photons from leptons.
251  if ( physical && !isSoftQCD() && beamHasGamma ) {
252 
253  // Save the x_gamma values for unresolved photons.
254  if ( !beamAhasResGamma ) beamAPtr->xGamma( x1());
255  if ( !beamBhasResGamma ) beamBPtr->xGamma( x2());
256 
257  // Sample the kinematics of virtual photons.
258  if ( !gammaKinPtr->sampleKTgamma() ) physical = false;
259 
260  // Processes with direct photons rescale momenta and cross section.
261  if ( physical && !(beamAhasResGamma && beamBhasResGamma) ) {
262  double sHatNew = gammaKinPtr->calcNewSHat( phaseSpacePtr->sHat() );
263  phaseSpacePtr->rescaleSigma( sHatNew);
264  }
265 
266  // With external photon flux calculate weight wrt. approximated fluxes.
267  if ( physical && beamHasGamma && approximatedGammaFlux ) {
268  wtPDF = phaseSpacePtr->weightGammaPDFApprox();
269  wtFlux = gammaKinPtr->fluxWeight();
270  } else {
271  wtPDF = 1.;
272  wtFlux = 1.;
273  }
274 
275  }
276 
277  // Flag to check if more events should be generated.
278  bool doTryNext = true;
279 
280  // Note if at end of Les Houches file, else do statistics.
281  if (isLHA && !physical) infoPtr->setEndOfFile(true);
282  else {
283 
284  // Check if the number of selected events should be limited.
285  if (nSelRequested > 0) doTryNext = (nTry < nSelRequested);
286 
287  // Only add to counter if event should count towards cross section.
288  if (doTryNext) ++nTry;
289 
290  // Statistics for separate Les Houches process codes. Insert new codes.
291  if (isLHA) {
292  int codeLHANow = lhaUpPtr->idProcess();
293  int iFill = -1;
294  for (int i = 0; i < int(codeLHA.size()); ++i)
295  if (codeLHANow == codeLHA[i]) iFill = i;
296  if (iFill >= 0) {
297  // Only add to counter if event should count towards cross section.
298  if (doTryNext) ++nTryLHA[iFill];
299  } else {
300  codeLHA.push_back(codeLHANow);
301  nTryLHA.push_back(1);
302  nSelLHA.push_back(0);
303  nAccLHA.push_back(0);
304  for (int i = int(codeLHA.size()) - 1; i > 0; --i) {
305  if (codeLHA[i] < codeLHA[i - 1]) {
306  swap(codeLHA[i], codeLHA[i - 1]);
307  swap(nTryLHA[i], nTryLHA[i - 1]);
308  swap(nSelLHA[i], nSelLHA[i - 1]);
309  swap(nAccLHA[i], nAccLHA[i - 1]);
310  }
311  else break;
312  }
313  }
314  }
315  }
316 
317  // Incoming top beams will not be handled, although allowed by LHAPDF.
318  if (isLHA && (abs(lhaUpPtr->id(1)) == 6 || abs(lhaUpPtr->id(2)) == 6)) {
319  infoPtr->errorMsg("Error in ProcessContainer::trialProcess(): top"
320  " not allowed incoming beam parton; event skipped");
321  return false;
322  }
323 
324  // Possibly fail, else cross section.
325  if (!physical) return false;
326  double sigmaNow = phaseSpacePtr->sigmaNow();
327 
328  // For photons with external flux correct the cross section.
329  if (beamHasGamma && approximatedGammaFlux) sigmaNow *= wtFlux * wtPDF;
330 
331  // Tell if this event comes with weight from cross section.
332  double sigmaWeight = 1.;
333  if (!isLHA && !increaseMaximum && sigmaNow > sigmaMx)
334  sigmaWeight = sigmaNow / sigmaMx;
335  if ( lhaStrat < 0 && sigmaNow < 0.) sigmaWeight = -1.;
336  if ( lhaStratAbs == 4) sigmaWeight = sigmaNow;
337 
338  // Also compensating weight from biased phase-space selection.
339  double biasWeight = phaseSpacePtr->biasSelectionWeight();
340  weightNow = sigmaWeight * biasWeight;
341 
342  // Set event weight to zero if this event should not count.
343  if (!doTryNext) weightNow = 0.0;
344 
345  infoPtr->setWeight( weightNow, lhaStrat);
346 
347  // Check that not negative cross section when not allowed.
348  if (!allowNegSig) {
349  if (sigmaNow < sigmaNeg) {
350  infoPtr->errorMsg("Warning in ProcessContainer::trialProcess: neg"
351  "ative cross section set 0", "for " + sigmaProcessPtr->name() );
352  sigmaNeg = sigmaNow;
353  }
354  if (sigmaNow < 0.) sigmaNow = 0.;
355  }
356 
357  // Cross section of process may come with a weight. Update sum.
358  double sigmaAdd = sigmaNow * biasWeight;
359  if (lhaStratAbs == 2 || lhaStratAbs == 3) sigmaAdd = sigmaSgn;
360 
361  // Do not add event weight to the cross section if it should not count.
362  if (!doTryNext) {
363  sigmaAdd=0.;
364  // Reset stored weight.
365  sigmaTemp = 0.;
366  sigma2Temp = 0.;
367  }
368 
369  // Only update once an event is accepted, as is needed if the event weight
370  // can still change by a finite amount due to reweighting.
371  if (lhaStratAbs >= 3 ) {
372  sigmaTemp = sigmaAdd;
373  sigma2Temp = pow2(sigmaAdd);
374  } else {
375  sigmaTemp += sigmaAdd;
376  sigma2Temp += pow2(sigmaAdd);
377  }
378 
379  // Check if maximum violated.
380  newSigmaMx = phaseSpacePtr->newSigmaMax();
381  if (newSigmaMx) sigmaMx = phaseSpacePtr->sigmaMax();
382 
383  // Select or reject trial point. Statistics.
384  bool select = true;
385  if (lhaStratAbs < 3) select
386  = (newSigmaMx || rndmPtr->flat() * abs(sigmaMx) < abs(sigmaNow));
387  if (select) {
388  // Only add to counter if event should count towards cross section.
389  if (doTryNext) ++nSel;
390 
391  if (isLHA) {
392  int codeLHANow = lhaUpPtr->idProcess();
393  int iFill = -1;
394  for (int i = 0; i < int(codeLHA.size()); ++i)
395  if (codeLHANow == codeLHA[i]) iFill = i;
396  // Only add to counter if event should count towards cross section.
397  if (iFill >= 0 && doTryNext) ++nSelLHA[iFill];
398  }
399  }
400  if (select || lhaStratAbs != 2) return select;
401 
402  }
403 
404 }
405 
406 //--------------------------------------------------------------------------
407 
408 // Accumulate statistics after user veto, including LHA code.
409 
410 void ProcessContainer::accumulate() {
411 
412  // Only update for non-vanishing weight - vanishing weight means
413  // event should not count as accepted.
414  double wgtNow = infoPtr->weight();
415  if (wgtNow != 0.0) {
416  ++nAcc;
417 
418  // infoPtr->weight() performs coversion to pb if lhaStratAbs = 4
419  if (lhaStratAbs == 4) wgtNow /= 1e9;
420 
421  wtAccSum += wgtNow;
422  if (isLHA) {
423  int codeLHANow = lhaUpPtr->idProcess();
424  int iFill = -1;
425  for (int i = 0; i < int(codeLHA.size()); ++i)
426  if (codeLHANow == codeLHA[i]) iFill = i;
427  if (iFill >= 0) ++nAccLHA[iFill];
428  }
429  }
430 
431 }
432 
433 //--------------------------------------------------------------------------
434 
435 // Pick flavours and colour flow of process.
436 
437 bool ProcessContainer::constructState() {
438 
439  // Construct flavour and colours for selected event.
440  if (isResolved && !isNonDiff) sigmaProcessPtr->pickInState();
441  sigmaProcessPtr->setIdColAcol();
442 
443  // Set beam modes correctly for given process.
444  // Propagate the sampled VM states to beams as they are now sampled.
445  if ( beamAPtr->isResolvedUnresolved() || beamBPtr->isResolvedUnresolved() )
446  setBeamModes();
447 
448  // Done.
449  return true;
450 
451 }
452 
453 //--------------------------------------------------------------------------
454 
455 // Set the photon modes according to the process.
456 
457 void ProcessContainer::setBeamModes(bool setVMD, bool isSampled) {
458 
459  // Set the modes for the current beams.
460  beamAPtr->setGammaMode(beamAgammaMode);
461  beamBPtr->setGammaMode(beamBgammaMode);
462 
463  // Set the gammaMode for given process.
464  if (beamAgammaMode <= 1 && beamBgammaMode <= 1) gammaModeEvent = 1;
465  else if (beamAgammaMode <= 1 && beamBgammaMode == 2) gammaModeEvent = 2;
466  else if (beamAgammaMode == 2 && beamBgammaMode <= 1) gammaModeEvent = 3;
467  else if (beamAgammaMode == 2 && beamBgammaMode == 2) gammaModeEvent = 4;
468  else gammaModeEvent = 0;
469 
470  // Propagate gammaMode and VMD state to info pointer.
471  infoPtr->setGammaMode(gammaModeEvent);
472 
473  // Check whether explicit VMD sampling is required later on.
474  if (setVMD && !isSampled) {
475  if (beamAgammaMode > 0) infoPtr->setVMDstateA(true, 22, 0.,0.);
476  if (beamBgammaMode > 0) infoPtr->setVMDstateB(true, 22, 0.,0.);
477  }
478 
479  // Propagate the sampled VMD states to beams.
480  if (isSampled) {
481  if (infoPtr->isVMDstateA()) beamAPtr->setVMDstate(true,
482  infoPtr->idVMDA(), infoPtr->mVMDA(), infoPtr->scaleVMDA(),
483  false);
484  if (infoPtr->isVMDstateB()) beamBPtr->setVMDstate(true,
485  infoPtr->idVMDB(), infoPtr->mVMDB(), infoPtr->scaleVMDB(),
486  false);
487  }
488 }
489 
490 //--------------------------------------------------------------------------
491 
492 // Give the hard subprocess with kinematics.
493 
494 bool ProcessContainer::constructProcess( Event& process, bool isHardest) {
495 
496  // Construct kinematics from selected phase space point.
497  if (!phaseSpacePtr->finalKin()) return false;
498  int nFin = sigmaProcessPtr->nFinal();
499 
500  // Basic info on process.
501  if (isHardest) infoPtr->setType( name(), code(), nFin, isNonDiff,
502  isResolved, isDiffA, isDiffB, isDiffC, isLHA);
503 
504  // Save sampled values for further use, requires info stored by setType.
505  if ( beamHasGamma && !isSoftQCD() ) gammaKinPtr->finalize();
506 
507  // Rescale the momenta again when unresolved photons after finalKin.
508  if ( beamHasGamma && !(beamAhasResGamma && beamBhasResGamma) ) {
509  double sHatNew = infoPtr->sHatNew();
510  phaseSpacePtr->rescaleMomenta( sHatNew);
511  }
512 
513  // Let hard process record begin with the event as a whole and
514  // the two incoming beam particles.
515  int idA = infoPtr->idA();
516  if (abs(idA) == idRenameBeams) idA = 16;
517  int idB = infoPtr->idB();
518  if (abs(idB) == idRenameBeams) idB = -16;
519  process.append( 90, -11, 0, 0, 0, 0, 0, 0,
520  Vec4(0., 0., 0., infoPtr->eCM()), infoPtr->eCM(), 0. );
521  process.append( idA, -12, 0, 0, 0, 0, 0, 0,
522  Vec4(0., 0., infoPtr->pzA(), infoPtr->eA()), infoPtr->mA(), 0. );
523  process.append( idB, -12, 0, 0, 0, 0, 0, 0,
524  Vec4(0., 0., infoPtr->pzB(), infoPtr->eB()), infoPtr->mB(), 0. );
525 
526  // Add intermediate gammas for lepton -> gamma -> parton processes
527  // for both non-diffractive and hard processes, including direct-resolved.
528  // Add a copy of hadron beam when e->gamma + p.
529  int nOffsetGamma = 0;
530  bool isGammaHadronDir = (beamAgammaMode == 2 && beamBgammaMode == 0)
531  || (beamAgammaMode == 0 && beamBgammaMode == 2);
532  if ( beamHasResGamma || (isGammaHadronDir && beamHasGamma) ) {
533  double xGm1 = beamAPtr->xGamma();
534  if ( !(beamAPtr->gammaInBeam()) ) {
535  process.append( beamAPtr->id(), -13, 1, 0, 0, 0, 0, 0,
536  Vec4(0., 0., infoPtr->pzA(), infoPtr->eA()), beamAPtr->m(), 0. );
537  } else {
538  process.append( 22, -13, 1, 0, 0, 0, 0, 0,
539  Vec4(0., 0., xGm1*infoPtr->pzA(), xGm1*infoPtr->eA()), 0, 0. );
540  }
541  process[1].daughter1(3);
542  ++nOffsetGamma;
543  double xGm2 = beamBPtr->xGamma();
544  if ( !(beamBPtr->gammaInBeam()) ) {
545  process.append( beamBPtr->id(), -13, 2, 0, 0, 0, 0, 0,
546  Vec4(0., 0., infoPtr->pzB(), infoPtr->eB()), beamBPtr->m(), 0. );
547  } else {
548  process.append( 22, -13, 2, 0, 0, 0, 0, 0,
549  Vec4(0., 0., xGm2*infoPtr->pzB(), xGm2*infoPtr->eB()), 0, 0. );
550  }
551  process[1 + nOffsetGamma].daughter1(3 + nOffsetGamma);
552  ++nOffsetGamma;
553  }
554 
555  // For nondiffractive process no interaction selected so far, so done.
556  if (isNonDiff) return true;
557 
558  // Entries 3 and 4, now to be added, come from 1 and 2.
559  // Offset from normal locations possible due to intermediate photons.
560  process[1 + nOffsetGamma].daughter1(3 + nOffsetGamma);
561  process[2 + nOffsetGamma].daughter1(4 + nOffsetGamma);
562  double scale = 0.;
563  double scalup = 0.;
564 
565  // For DiffC entries 3 - 5 come jointly from 1 and 2 (to keep HepMC happy).
566  if (isDiffC) {
567  process[1].daughters(3, 5);
568  process[2].daughters(3, 5);
569  }
570 
571  // Insert the subprocess partons - resolved processes.
572  int idRes = sigmaProcessPtr->idSChannel();
573  if (isResolved && !isLHA) {
574 
575  // NOAM: Mothers and daughters without/with intermediate state.
576  int m_M1 = 3;
577  int m_M2 = 4;
578  int m_D1 = 5;
579  int m_D2 = (nFin == 1) ? 0 : 4 + nFin;
580  if (idRes != 0) {
581  m_M1 = 5;
582  m_M2 = 0;
583  m_D1 = 5;
584  m_D2 = 0;
585  }
586 
587  // Find scale from which to begin MPI/ISR/FSR evolution.
588  scale = sqrt(Q2Fac());
589  process.scale( scale );
590 
591  // Loop over incoming and outgoing partons.
592  int colOffset = process.lastColTag();
593  for (int i = 1; i <= 2 + nFin; ++i) {
594 
595  // Read out particle info from SigmaProcess object.
596  int id = sigmaProcessPtr->id(i);
597  int status = (i <= 2) ? -21 : 23;
598  int mother1 = (i <= 2) ? i : m_M1 ;
599  int mother2 = (i <= 2) ? 0 : m_M2 ;
600  int daughter1 = (i <= 2) ? m_D1 : 0;
601  int daughter2 = (i <= 2) ? m_D2 : 0;
602  int col = sigmaProcessPtr->col(i);
603  if (col > 0) col += colOffset;
604  else if (col < 0) col -= colOffset;
605  int acol = sigmaProcessPtr->acol(i);
606  if (acol > 0) acol += colOffset;
607  else if (acol < 0) acol -= colOffset;
608 
609  // If extra photons in event record, offset the mother/daughter list.
610  if ( beamAhasResGamma || beamBhasResGamma
611  || (beamAgammaMode == 2 && beamBgammaMode == 0)
612  || (beamAgammaMode == 0 && beamBgammaMode == 2) ) {
613  mother1 += nOffsetGamma;
614  if (mother2 > 0) mother2 += nOffsetGamma;
615  if (daughter1 > 0) daughter1 += nOffsetGamma;
616  if (daughter2 > 0) daughter2 += nOffsetGamma;
617  }
618 
619  // Append to process record.
620  int iNow = process.append( id, status, mother1, mother2,
621  daughter1, daughter2, col, acol, phaseSpacePtr->p(i),
622  phaseSpacePtr->m(i), scale);
623 
624  // NOAM: If there is an intermediate state, insert the it in
625  // the process record after the two incoming particles.
626  if (i == 2 && idRes != 0) {
627 
628  // Sign of intermediate state: go by charge.
629  if (particleDataPtr->hasAnti(idRes)
630  && process[3].chargeType() + process[4].chargeType() < 0)
631  idRes *= -1;
632 
633  // The colour configuration of the intermediate state has to be
634  // resolved separately.
635  col = 0;
636  acol = 0;
637  int m_col1 = sigmaProcessPtr->col(1);
638  int m_acol1 = sigmaProcessPtr->acol(1);
639  int m_col2 = sigmaProcessPtr->col(2);
640  int m_acol2 = sigmaProcessPtr->acol(2);
641  if (m_col1 == m_acol2 && m_col2 != m_acol1) {
642  col = m_col2;
643  acol = m_acol1;
644  } else if (m_col2 == m_acol1 && m_col1 != m_acol2) {
645  col = m_col1;
646  acol = m_acol2;
647  }
648  if ( col > 0) col += colOffset;
649  else if ( col < 0) col -= colOffset;
650  if (acol > 0) acol += colOffset;
651  else if (acol < 0) acol -= colOffset;
652 
653  // Insert the intermediate state into the event record.
654  Vec4 pIntMed = phaseSpacePtr->p(1) + phaseSpacePtr->p(2);
655  process.append( idRes, -22, 3, 4, 6, 5 + nFin, col, acol,
656  pIntMed, pIntMed.mCalc(), scale);
657  }
658 
659  // Pick lifetime where relevant, else not.
660  if (process[iNow].tau0() > 0.) process[iNow].tau(
661  process[iNow].tau0() * rndmPtr->exp() );
662  }
663  }
664 
665  // Insert the outgoing particles - unresolved processes.
666  else if (!isLHA) {
667 
668  // COR: Special handling of soft QCD with photons. If VMD states,
669  // then outgoing photon has to be changed to VMD
670  if (isSoftQCD() && (infoPtr->isVMDstateA()
671  || infoPtr->isVMDstateB())) {
672  int id3orig = sigmaProcessPtr->id(3);
673  int status3 = (id3orig == process[1+nOffsetGamma].id()) ? 14 : 15;
674  int id3 = (status3 == 14 && infoPtr->isVMDstateA())
675  ? infoPtr->idVMDA() : id3orig;
676  process.append( id3, status3, 1 + nOffsetGamma, 0, 0, 0, 0, 0,
677  phaseSpacePtr->p(3), phaseSpacePtr->m(3));
678  int id4orig = sigmaProcessPtr->id(4);
679  int status4 = (id4orig == process[2+nOffsetGamma].id()) ? 14 : 15;
680  int id4 = (status4 == 14 && infoPtr->isVMDstateB())
681  ? infoPtr->idVMDB() : id4orig;
682  process.append( id4, status4, 2 + nOffsetGamma, 0, 0, 0, 0, 0,
683  phaseSpacePtr->p(4), phaseSpacePtr->m(4));
684 
685  } else {
686  int id3 = sigmaProcessPtr->id(3);
687  int status3 = (id3 == process[1].id()) ? 14 : 15;
688  process.append( id3, status3, 1 + nOffsetGamma, 0, 0, 0, 0, 0,
689  phaseSpacePtr->p(3), phaseSpacePtr->m(3));
690  int id4 = sigmaProcessPtr->id(4);
691  int status4 = (id4 == process[2].id()) ? 14 : 15;
692  process.append( id4, status4, 2 + nOffsetGamma, 0, 0, 0, 0, 0,
693  phaseSpacePtr->p(4), phaseSpacePtr->m(4));
694  }
695 
696  // For central diffraction: two scattered protons inserted so far,
697  // but modify mothers, add also centrally-produced hadronic system.
698  if (isDiffC) {
699  process[3].mothers( 1, 2);
700  process[4].mothers( 1, 2);
701  int id5 = sigmaProcessPtr->id(5);
702  int status5 = 15;
703  process.append( id5, status5, 1, 2, 0, 0, 0, 0,
704  phaseSpacePtr->p(5), phaseSpacePtr->m(5));
705  }
706 
707  }
708 
709  // Insert the outgoing particles - Les Houches Accord processes.
710  else {
711 
712  // Check if second process; if so partons must be in order.
713  int n21 = 0;
714  int nOffsetSecond = 0;
715  for (int i = 1; i < lhaUpPtr->sizePart(); ++i)
716  if (lhaUpPtr->status(i) == -1) {
717  ++n21;
718  if (n21 == 3) nOffsetSecond = i - 1;
719  }
720  bool twoHard = (n21 == 4);
721 
722  // Since LHA partons may be out of order, determine correct one.
723  // (Recall that zeroth particle is empty.)
724  vector<int> newPos;
725  newPos.reserve(lhaUpPtr->sizePart());
726  newPos.push_back(0);
727  for (int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
728  // Simple copy in unchanged order for two hard interactions.
729  if (twoHard) newPos.push_back(iNew + 1);
730  // For iNew == 0 look for the two incoming partons, then for
731  // partons having them as mothers, and so on layer by layer.
732  else {
733  for (int i = 1; i < lhaUpPtr->sizePart(); ++i)
734  if (lhaUpPtr->mother1(i) == newPos[iNew]) newPos.push_back(i);
735  if (int(newPos.size()) <= iNew) break;
736  }
737  }
738 
739  // Find scale from which to begin MPI/ISR/FSR evolution.
740  scalup = lhaUpPtr->scale();
741  scale = scalup;
742  double scalePr = (scale < 0.) ? sqrt(Q2Fac()) : scale;
743  process.scale( scalePr);
744  if (twoHard) process.scaleSecond( scalePr);
745  if ( lhaUpPtr->scaleShowersIsSet() ) {
746  process.scale( lhaUpPtr->scaleShowers(0) );
747  process.scaleSecond( lhaUpPtr->scaleShowers(1) );
748  }
749  double scalePr2 = process.scaleSecond();
750 
751  // Copy over info from LHA event to process, in proper order.
752  vector<int> iFinal;
753  int iIn = 0;
754  for (int i = 1; i < lhaUpPtr->sizePart(); ++i) {
755  int iOld = newPos[i];
756  int id = lhaUpPtr->id(iOld);
757  if (i == 1 && abs(id) == idRenameBeams) id = 16;
758  if (i == 2 && abs(id) == idRenameBeams) id = -16;
759  int idAbs = abs(id);
760 
761  // Translate from LHA status codes.
762  int lhaStatus = lhaUpPtr->status(iOld);
763  int status = -21;
764  if (lhaStatus == 2 || lhaStatus == 3) status = -22;
765  if (lhaStatus == 1) status = 23;
766 
767  // Find where mothers have been moved by reordering.
768  int mother1Old = lhaUpPtr->mother1(iOld);
769  int mother2Old = lhaUpPtr->mother2(iOld);
770  int mother1 = 0;
771  int mother2 = 0;
772  for (int im = 1; im < i; ++im) {
773  if (mother1Old == newPos[im]) mother1 = im + 2;
774  if (mother2Old == newPos[im]) mother2 = im + 2;
775  }
776  if (status == -21) {
777  ++iIn;
778  mother1 = 1 + (iIn - 1)%2;
779  }
780 
781  // Ensure that second mother = 0 except for bona fide carbon copies.
782  if (mother1 > 0 && mother2 == mother1) {
783  int sister1 = process[mother1].daughter1();
784  int sister2 = process[mother1].daughter2();
785  if (sister2 != sister1 && sister2 != 0) mother2 = 0;
786  }
787 
788  // Find daughters and where they have been moved by reordering.
789  // (Values shifted two steps to account for inserted beams.)
790  int daughter1 = 0;
791  int daughter2 = 0;
792  for (int im = i + 1; im < lhaUpPtr->sizePart(); ++im) {
793  if (lhaUpPtr->mother1(newPos[im]) == iOld
794  || lhaUpPtr->mother2(newPos[im]) == iOld) {
795  if (daughter1 == 0 || im + 2 < daughter1) daughter1 = im + 2;
796  if (daughter2 == 0 || im + 2 > daughter2) daughter2 = im + 2;
797  }
798  }
799  // For 2 -> 1 hard scatterings reset second daughter to 0.
800  if (daughter2 == daughter1) daughter2 = 0;
801 
802  // Colour trivial, except reset irrelevant colour indices.
803  int colType = particleDataPtr->colType(id);
804  int col1 = (colType == 1 || colType == 2 || abs(colType) == 3)
805  ? lhaUpPtr->col1(iOld) : 0;
806  int col2 = (colType == -1 || colType == 2 || abs(colType) == 3)
807  ? lhaUpPtr->col2(iOld) : 0;
808 
809  // Copy momentum and ensure consistent (E, p, m) set.
810  double px = lhaUpPtr->px(iOld);
811  double py = lhaUpPtr->py(iOld);
812  double pz = lhaUpPtr->pz(iOld);
813  double e = lhaUpPtr->e(iOld);
814  double m = lhaUpPtr->m(iOld);
815  if (mRecalculate > 0. && m > mRecalculate)
816  m = sqrtpos( e*e - px*px - py*py - pz*pz);
817  else e = sqrt( m*m + px*px + py*py + pz*pz);
818 
819  // Polarization.
820  double pol = lhaUpPtr->spin(iOld);
821 
822  // Allow scale setting for generic partons.
823  double scaleShow = lhaUpPtr->scale(iOld);
824 
825  // For resonance decay products use resonance mass as scale.
826  double scaleNow = (iIn < 3) ? scalePr : scalePr2;
827  int motherBeg = (iIn < 3) ? 4 : 4 + nOffsetSecond;
828  if (mother1 > motherBeg && !useStrictLHEFscales)
829  scaleNow = process[mother1].m();
830  if (scaleShow >= 0.0) scaleNow = scaleShow;
831 
832  // Store Les Houches Accord partons. Boost to CM frame if not already.
833  int iNow = process.append( id, status, mother1, mother2, daughter1,
834  daughter2, col1, col2, Vec4(px, py, pz, e), m, scaleNow, pol);
835  if (isAsymLHA) process[iNow].bst( 0., 0., -betazLHA);
836 
837  // Check if need to store lifetime.
838  double tau = lhaUpPtr->tau(iOld);
839  if ( (setLifetime == 1 && idAbs == 15) || setLifetime == 2)
840  tau = process[iNow].tau0() * rndmPtr->exp();
841  if (tau > 0.) process[iNow].tau(tau);
842 
843  // Track outgoing final particles.
844  if (status == 23) iFinal.push_back( iNow);
845  }
846 
847  // Second pass to catch vanishing final lepton and quark masses.
848  int iFinalSz = iFinal.size();
849  for (int iF = 0; iF < iFinalSz; ++iF) {
850  int iMod = iFinal[iF];
851  int iQLmod = 0;
852  double mOld = process[iMod].m();
853  for (int iQL = 1; iQL < 9; ++iQL)
854  if (process[iMod].idAbs() == idNewM[iQL]) {
855  if ( iQL < 6 && setQuarkMass > 0 && (iQL == 4 || iQL == 5
856  || setQuarkMass == 2) && (mOld < 0.5 * mNewM[iQL]
857  || mOld > 1.5 * mNewM[iQL]) ) iQLmod = iQL;
858  if ( iQL > 5 && setLeptonMass > 0 && (setLeptonMass == 2
859  || mOld < 0.9 * mNewM[iQL] || mOld > 1.1 * mNewM[iQL]) )
860  iQLmod = iQL;
861  }
862  if (iQLmod == 0) continue;
863  double mNew = mNewM[iQLmod];
864 
865  // Find partner to exchange energy and momentum with: general.
866  int iRec = 0;
867  if (iFinalSz == 2) iRec = iFinal[1 - iF];
868  else if (process[iMod].mother1() > 2) {
869  int iMother = process[iMod].mother1();
870  int iMDau1 = process[iMother].daughter1();
871  int iMDau2 = process[iMother].daughter2();
872  if (iMDau2 == iMDau1 + 1 && iMod == iMDau1) iRec = iMDau2;
873  if (iMDau2 == iMDau1 + 1 && iMod == iMDau2) iRec = iMDau1;
874  }
875 
876  // Find partner to exchange energy and momentum with: lepton.
877  if (iRec == 0 && iQLmod > 5) {
878  for (int iR = 0; iR < iFinalSz; ++iR) if (iR != iF) {
879  int iRtmp = iFinal[iR];
880  if (process[iRtmp].idAbs() == idNewM[iQLmod] + 1
881  && process[iRtmp].id() * process[iMod].id() < 0) iRec = iRtmp;
882  }
883  }
884  if (iRec == 0 && iQLmod > 5) {
885  for (int iR = 0; iR < iFinalSz; ++iR) if (iR != iF) {
886  int iRtmp = iFinal[iR];
887  if (process[iRtmp].id() == -process[iMod].id()) iRec = iRtmp;
888  }
889  }
890 
891  // Find partner to exchange energy and momentum with: quark.
892  if (iRec == 0 && iQLmod < 6 && process[iMod].col() != 0) {
893  for (int iR = 0; iR < iFinalSz; ++iR) if (iR != iF) {
894  int iRtmp = iFinal[iR];
895  if (process[iRtmp].acol() == process[iMod].col()) iRec = iRtmp;
896  }
897  }
898  if (iRec == 0 && iQLmod < 6 && process[iMod].acol() != 0) {
899  for (int iR = 0; iR < iFinalSz; ++iR) if (iR != iF) {
900  int iRtmp = iFinal[iR];
901  if (process[iRtmp].col() == process[iMod].acol()) iRec = iRtmp;
902  }
903  }
904 
905  // Find partner to exchange energy and momentum with: largest mass.
906  if (iRec == 0) {
907  double mMax = 0.;
908  for (int iR = 0; iR < iFinalSz; ++iR) if (iR != iF) {
909  int iRtmp = iFinal[iR];
910  double mTmp = m( process[iMod], process[iRtmp]) - process[iRtmp].m();
911  if (mTmp > mMax) { iRec = iRtmp; mMax = mTmp;}
912  }
913  }
914 
915  // Carry out exchange of energy and momentum, if possible.
916  bool doneShift = false;
917  Vec4 pMod = process[iMod].p();
918  if (iRec != 0) {
919  Vec4 pRecOld = process[iRec].p();
920  Vec4 pRecNew = pRecOld;
921  double mRec = process[iRec].m();
922  if ( pShift( pMod, pRecNew, mNew, mRec) ) {
923  process[iMod].p( pMod);
924  process[iMod].m( mNew);
925  process[iRec].p( pRecNew);
926  // Recursive boost of recoiler decay products, if necessary.
927  if (!process[iRec].isFinal()) {
928  vector<int> iDauRec = process[iRec].daughterListRecursive();
929  RotBstMatrix Mdau;
930  Mdau.bst( pRecOld, pRecNew);
931  for (int iD = 0; iD < int(iDauRec.size()); ++iD)
932  process[iDauRec[iD]].rotbst( Mdau);
933  }
934  doneShift = true;
935  }
936  }
937 
938  // Emergency solution: just add energy as needed, => new incoming.
939  if (!doneShift) {
940  pMod.e( sqrtpos( pMod.pAbs2() + mNew * mNew) );
941  process[iMod].p( pMod);
942  process[iMod].m( mNew);
943  infoPtr->errorMsg("Warning in ProcessContainer::constructProcess: "
944  "unsuitable recoiler found");
945  }
946  }
947 
948  // Reassign momenta and masses for incoming partons.
949  if (matchInOut && !twoHard) {
950  Vec4 pSumOut;
951  for (int iF = 0; iF < iFinalSz; ++iF)
952  pSumOut += process[iFinal[iF]].p();
953  double e1 = 0.5 * (pSumOut.e() + pSumOut.pz());
954  double e2 = 0.5 * (pSumOut.e() - pSumOut.pz());
955  if (max(e1, e2) > 0.500001 * process[0].e()) {
956  infoPtr->errorMsg("Error in ProcessContainer::constructProcess: "
957  "setting mass failed");
958  return false;
959  }
960  e1 = min( e1, 0.5 * process[0].e());
961  e2 = min( e2, 0.5 * process[0].e());
962  process[3].pz( e1);
963  process[3].e( e1);
964  process[3].m( 0.);
965  process[4].pz(-e2);
966  process[4].e( e2);
967  process[4].m( 0.);
968  }
969  }
970 
971  // Loop through decay chains and set secondary vertices when needed.
972  for (int i = 3; i < process.size(); ++i) {
973  int iMother = process[i].mother1();
974 
975  // If sister to already assigned vertex then assign same.
976  if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex() )
977  process[i].vProd( process[i - 1].vProd() );
978 
979  // Else if mother already has vertex and/or lifetime then assign.
980  else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
981  process[i].vProd( process[iMother].vDec() );
982  }
983 
984  // Further info on process. Reset quantities that may or may not be known.
985  int id1Now = process[3 + nOffsetGamma].id();
986  int id2Now = process[4 + nOffsetGamma].id();
987  int id1pdf = 0;
988  int id2pdf = 0;
989  double x1pdf = 0.;
990  double x2pdf = 0.;
991  double pdf1 = 0.;
992  double pdf2 = 0.;
993  double tHat = 0.;
994  double uHat = 0.;
995  double pTHatL = 0.;
996  double m3 = 0.;
997  double m4 = 0.;
998  double theta = 0.;
999  double phi = 0.;
1000  double x1Now, x2Now, Q2FacNow, alphaEM, alphaS, Q2Ren, sHat;
1001 
1002  // Internally generated and stored information.
1003  if (!isLHA) {
1004  id1pdf = id1Now;
1005  id2pdf = id2Now;
1006  x1Now = phaseSpacePtr->x1();
1007  x2Now = phaseSpacePtr->x2();
1008  x1pdf = x1Now;
1009  x2pdf = x2Now;
1010  pdf1 = sigmaProcessPtr->pdf1();
1011  pdf2 = sigmaProcessPtr->pdf2();
1012  Q2FacNow = sigmaProcessPtr->Q2Fac();
1013  alphaEM = sigmaProcessPtr->alphaEMRen();
1014  alphaS = sigmaProcessPtr->alphaSRen();
1015  Q2Ren = sigmaProcessPtr->Q2Ren();
1016  sHat = phaseSpacePtr->sHat();
1017  tHat = phaseSpacePtr->tHat();
1018  uHat = phaseSpacePtr->uHat();
1019  pTHatL = phaseSpacePtr->pTHat();
1020  m3 = phaseSpacePtr->m(3);
1021  m4 = phaseSpacePtr->m(4);
1022  theta = phaseSpacePtr->thetaHat();
1023  phi = phaseSpacePtr->phiHat();
1024  }
1025 
1026  // Les Houches Accord process partly available, partly to be constructed.
1027  else {
1028  x1Now = 2. * process[3].e() / infoPtr->eCM();
1029  x2Now = 2. * process[4].e() / infoPtr->eCM();
1030  Q2FacNow = (scale < 0.) ? sigmaProcessPtr->Q2Fac() : pow2(scale);
1031  alphaEM = lhaUpPtr->alphaQED();
1032  if (alphaEM < 0.001) alphaEM = sigmaProcessPtr->alphaEMRen();
1033  alphaS = lhaUpPtr->alphaQCD();
1034  if (alphaS < 0.001) alphaS = sigmaProcessPtr->alphaSRen();
1035  Q2Ren = (scale < 0.) ? sigmaProcessPtr->Q2Ren() : pow2(scale);
1036  Vec4 pSum = process[3].p() + process[4].p();
1037  sHat = pSum * pSum;
1038  int nFinLH = 0;
1039  double pTLH = 0.;
1040  for (int i = 5; i < process.size(); ++i)
1041  if (process[i].mother1() == 3 && process[i].mother2() == 4) {
1042  ++nFinLH;
1043  pTLH += process[i].pT();
1044  }
1045  if (nFinLH > 0) pTHatL = pTLH / nFinLH;
1046 
1047  // Read info on parton densities if provided.
1048  id1pdf = lhaUpPtr->id1pdf();
1049  id2pdf = lhaUpPtr->id2pdf();
1050  x1pdf = lhaUpPtr->x1pdf();
1051  x2pdf = lhaUpPtr->x2pdf();
1052  if (lhaUpPtr->pdfIsSet()) {
1053  pdf1 = lhaUpPtr->pdf1();
1054  pdf2 = lhaUpPtr->pdf2();
1055  Q2FacNow = pow2(lhaUpPtr->scalePDF());
1056  }
1057 
1058  // Reconstruct kinematics of 2 -> 2 processes from momenta.
1059  if (nFin == 2) {
1060  Vec4 pDifT = process[3].p() - process[5].p();
1061  tHat = pDifT * pDifT;
1062  Vec4 pDifU = process[3].p() - process[6].p();
1063  uHat = pDifU * pDifU;
1064  pTHatL = process[5].pT();
1065  m3 = process[5].m();
1066  m4 = process[6].m();
1067  Vec4 p5 = process[5].p();
1068  p5.bstback(pSum);
1069  theta = p5.theta();
1070  phi = process[5].phi();
1071  }
1072  }
1073 
1074  // Store information.
1075  if (isHardest) {
1076  infoPtr->setPDFalpha( 0, id1pdf, id2pdf, x1pdf, x2pdf, pdf1, pdf2,
1077  Q2FacNow, alphaEM, alphaS, Q2Ren, scalup);
1078  infoPtr->setKin( 0, id1Now, id2Now, x1Now, x2Now, sHat, tHat, uHat,
1079  pTHatL, m3, m4, theta, phi);
1080  }
1081  infoPtr->setTypeMPI( code(), pTHatL);
1082 
1083  // For Les Houches event store subprocess classification.
1084  if (isLHA) {
1085  int codeSub = lhaUpPtr->idProcess();
1086  ostringstream nameSub;
1087  nameSub << "user process " << codeSub;
1088  infoPtr->setSubType( 0, nameSub.str(), codeSub, nFin);
1089  }
1090 
1091  // Done.
1092  return true;
1093 
1094 }
1095 
1096 //--------------------------------------------------------------------------
1097 
1098 // Give the hard resonance decay chain from Les Houches input.
1099 // Note: mildly modified copy of some constructProcess code; to streamline.
1100 
1101 bool ProcessContainer::constructDecays( Event& process) {
1102 
1103  // Let hard process record begin with the event as a whole.
1104  process.clear();
1105  process.append( 90, -11, 0, 0, 0, 0, 0, 0, Vec4(0., 0., 0., 0.), 0., 0. );
1106 
1107  // Since LHA partons may be out of order, determine correct one.
1108  // (Recall that zeroth particle is empty.)
1109  vector<int> newPos;
1110  newPos.reserve(lhaUpPtr->sizePart());
1111  newPos.push_back(0);
1112  for (int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
1113  // For iNew == 0 look for the two incoming partons, then for
1114  // partons having them as mothers, and so on layer by layer.
1115  for (int i = 1; i < lhaUpPtr->sizePart(); ++i)
1116  if (lhaUpPtr->mother1(i) == newPos[iNew]) newPos.push_back(i);
1117  if (int(newPos.size()) <= iNew) break;
1118  }
1119 
1120  // Find scale from which to begin MPI/ISR/FSR evolution.
1121  double scale = lhaUpPtr->scale();
1122  process.scale( scale);
1123  Vec4 pSum;
1124 
1125  // Copy over info from LHA event to process, in proper order.
1126  for (int i = 1; i < lhaUpPtr->sizePart(); ++i) {
1127  int iOld = newPos[i];
1128  int id = lhaUpPtr->id(iOld);
1129 
1130  // Translate from LHA status codes.
1131  int lhaStatus = lhaUpPtr->status(iOld);
1132  int status = -21;
1133  if (lhaStatus == 2 || lhaStatus == 3) status = -22;
1134  if (lhaStatus == 1) status = 23;
1135 
1136  // Find where mothers have been moved by reordering.
1137  int mother1Old = lhaUpPtr->mother1(iOld);
1138  int mother2Old = lhaUpPtr->mother2(iOld);
1139  int mother1 = 0;
1140  int mother2 = 0;
1141  for (int im = 1; im < i; ++im) {
1142  if (mother1Old == newPos[im]) mother1 = im;
1143  if (mother2Old == newPos[im]) mother2 = im;
1144  }
1145 
1146  // Ensure that second mother = 0 except for bona fide carbon copies.
1147  if (mother1 > 0 && mother2 == mother1) {
1148  int sister1 = process[mother1].daughter1();
1149  int sister2 = process[mother1].daughter2();
1150  if (sister2 != sister1 && sister2 != 0) mother2 = 0;
1151  }
1152 
1153  // Find daughters and where they have been moved by reordering.
1154  int daughter1 = 0;
1155  int daughter2 = 0;
1156  for (int im = i + 1; im < lhaUpPtr->sizePart(); ++im) {
1157  if (lhaUpPtr->mother1(newPos[im]) == iOld
1158  || lhaUpPtr->mother2(newPos[im]) == iOld) {
1159  if (daughter1 == 0 || im < daughter1) daughter1 = im;
1160  if (daughter2 == 0 || im > daughter2) daughter2 = im;
1161  }
1162  }
1163  // For 2 -> 1 hard scatterings reset second daughter to 0.
1164  if (daughter2 == daughter1) daughter2 = 0;
1165 
1166  // Colour trivial, except reset irrelevant colour indices.
1167  int colType = particleDataPtr->colType(id);
1168  int col1 = (colType == 1 || colType == 2 || abs(colType) == 3)
1169  ? lhaUpPtr->col1(iOld) : 0;
1170  int col2 = (colType == -1 || colType == 2 || abs(colType) == 3)
1171  ? lhaUpPtr->col2(iOld) : 0;
1172 
1173  // Momentum trivial.
1174  double px = lhaUpPtr->px(iOld);
1175  double py = lhaUpPtr->py(iOld);
1176  double pz = lhaUpPtr->pz(iOld);
1177  double e = lhaUpPtr->e(iOld);
1178  double m = lhaUpPtr->m(iOld);
1179  if (status > 0) pSum += Vec4( px, py, pz, e);
1180 
1181  // Polarization.
1182  double pol = lhaUpPtr->spin(iOld);
1183 
1184  // For resonance decay products use resonance mass as scale.
1185  double scaleNow = scale;
1186  if (mother1 > 0) scaleNow = process[mother1].m();
1187 
1188  // Store Les Houches Accord partons.
1189  int iNow = process.append( id, status, mother1, mother2, daughter1,
1190  daughter2, col1, col2, Vec4(px, py, pz, e), m, scaleNow, pol);
1191 
1192  // Check if need to store lifetime.
1193  double tau = lhaUpPtr->tau(iOld);
1194  if ( (setLifetime == 1 && abs(id) == 15) || setLifetime == 2)
1195  tau = process[iNow].tau0() * rndmPtr->exp();
1196  if (tau > 0.) process[iNow].tau(tau);
1197  }
1198 
1199  // Update four-momentum of system as a whole.
1200  process[0].p( pSum);
1201  process[0].m( pSum.mCalc());
1202 
1203  // Loop through decay chains and set secondary vertices when needed.
1204  for (int i = 1; i < process.size(); ++i) {
1205  int iMother = process[i].mother1();
1206 
1207  // If sister to already assigned vertex then assign same.
1208  if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex()
1209  && i > 1) process[i].vProd( process[i - 1].vProd() );
1210 
1211  // Else if mother already has vertex and/or lifetime then assign.
1212  else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
1213  process[i].vProd( process[iMother].vDec() );
1214  }
1215 
1216  // Done.
1217  return true;
1218 
1219 }
1220 
1221 //--------------------------------------------------------------------------
1222 
1223 // Handle resonance decays.
1224 
1225 bool ProcessContainer::decayResonances( Event& process) {
1226 
1227  // Save current event-record size and status codes.
1228  process.saveSize();
1229  vector<int> statusSave( process.size());
1230  for (int i = 0; i < process.size(); ++i)
1231  statusSave[i] = process[i].status();
1232  bool physical = true;
1233  bool newChain = false;
1234  bool newFlavours = false;
1235 
1236  // Do loop over user veto.
1237  do {
1238 
1239  // Do sequential chain of uncorrelated isotropic decays.
1240  do {
1241  physical = resDecaysPtr->next( process);
1242  if (!physical) return false;
1243 
1244  // Check whether flavours should be correlated.
1245  // (Currently only relevant for f fbar -> gamma*/Z0 gamma*/Z0.)
1246  newFlavours = ( sigmaProcessPtr->weightDecayFlav( process)
1247  < rndmPtr->flat() );
1248 
1249  // Reset the decay chains if have to redo.
1250  if (newFlavours) {
1251  process.restoreSize();
1252  for (int i = 0; i < process.size(); ++i)
1253  process[i].status( statusSave[i]);
1254  }
1255 
1256  // Loop back where required to generate new decays with new flavours.
1257  } while (newFlavours);
1258 
1259  // Correct to nonisotropic decays.
1260  phaseSpacePtr->decayKinematics( process);
1261 
1262  // Optionally user hooks check/veto on decay chain.
1263  if (canVetoResDecay)
1264  newChain = userHooksPtr->doVetoResonanceDecays( process);
1265 
1266  // Reset the decay chains if have to redo.
1267  if (newChain) {
1268  process.restoreSize();
1269  for (int i = 0; i < process.size(); ++i)
1270  process[i].status( statusSave[i]);
1271  }
1272 
1273  // Loop back where required to generate new decay chain.
1274  } while(newChain);
1275 
1276  // Done.
1277  return true;
1278 
1279 }
1280 
1281 //--------------------------------------------------------------------------
1282 
1283 // Reset event generation statistics; but NOT maximum of cross section.
1284 
1285 void ProcessContainer::reset() {
1286 
1287  nTry = 0;
1288  nSel = 0;
1289  nAcc = 0;
1290  nTryStat = 0;
1291  sigmaSum = 0.;
1292  sigma2Sum = 0.;
1293  sigmaNeg = 0.;
1294  sigmaAvg = 0.;
1295  sigmaFin = 0.;
1296  deltaFin = 0.;
1297  wtAccSum = 0.;
1298  sigmaTemp = 0.;
1299  sigma2Temp= 0.;
1300 
1301 }
1302 
1303 //--------------------------------------------------------------------------
1304 
1305 // Estimate integrated cross section and its uncertainty.
1306 
1307 void ProcessContainer::sigmaDelta() {
1308 
1309  // Initial values. No analysis meaningful unless accepted events.
1310  nTryStat = nTry;
1311  sigmaAvg = 0.;
1312  sigmaFin = 0.;
1313  deltaFin = 0.;
1314  if (nAcc == 0) return;
1315 
1316  // Only update once an event is accepted, as is needed if the event weight
1317  // can still change by a finite amount due to reweighting.
1318  double wgtNow = infoPtr->weight();
1319  // infoPtr->weight() performs conversion to pb if lhaStratAbs = 4
1320  if (lhaStrat == 0) wgtNow = sigmaTemp;
1321  if (lhaStratAbs == 3) wgtNow *= sigmaTemp;
1322  if (lhaStratAbs == 4) wgtNow /= 1e9;
1323  sigmaSum += wgtNow;
1324 
1325  double wgtNow2 = 1.0;
1326  if (lhaStrat == 0) wgtNow2 = sigma2Temp;
1327  if (lhaStratAbs == 3) wgtNow2 = pow2(wgtNow)*sigma2Temp;
1328  if (lhaStratAbs == 4) wgtNow2 = pow2(wgtNow/1e9);
1329  sigma2Sum += wgtNow2;
1330 
1331  sigmaTemp = 0.0;
1332  sigma2Temp = 0.0;
1333 
1334  // Average value. No error analysis unless at least two events.
1335  double nTryInv = 1. / nTry;
1336  double nSelInv = 1. / nSel;
1337  double nAccInv = 1. / nAcc;
1338  sigmaAvg = sigmaSum * nTryInv;
1339  double fracAcc = nAcc * nSelInv;
1340 
1341  // If Pythia should not perform the unweighting, always simply add accepted
1342  // event weights.
1343  if (lhaStratAbs >= 3 ) fracAcc = 1.;
1344  sigmaFin = sigmaAvg * fracAcc;
1345  deltaFin = sigmaFin;
1346  if (nAcc == 1) return;
1347 
1348  // Estimated error. Quadratic sum of cross section term and
1349  // binomial from accept/reject step.
1350  double delta2Sig = (lhaStratAbs == 3) ? normVar3
1351  : (sigma2Sum * nTryInv - pow2(sigmaAvg)) * nTryInv / pow2(sigmaAvg);
1352  double delta2Veto = (nSel - nAcc) * nAccInv * nSelInv;
1353  double delta2Sum = delta2Sig + delta2Veto;
1354  deltaFin = sqrtpos(delta2Sum) * sigmaFin;
1355 
1356 }
1357 
1358 //==========================================================================
1359 
1360 // SetupContainer class.
1361 // Turns list of user-desired processes into a vector of containers.
1362 
1363 //--------------------------------------------------------------------------
1364 
1365 // Main routine to initialize list of processes.
1366 
1367 bool SetupContainers::init(vector<ProcessContainer*>& containerPtrs,
1368  Info* infoPtr) {
1369 
1370  // Reset process list, if filled in previous subrun.
1371  if (containerPtrs.size() > 0) {
1372  for (int i = 0; i < int(containerPtrs.size()); ++i)
1373  delete containerPtrs[i];
1374  containerPtrs.clear();
1375  }
1376  SigmaProcess* sigmaPtr;
1377 
1378  // Reference to settings and pointer to particle data.
1379  Settings& settings = *infoPtr->settingsPtr;
1380  ParticleData* particleDataPtr = infoPtr->particleDataPtr;
1381 
1382  // Set up requested objects for soft QCD processes.
1383  bool softQCD = settings.flag("SoftQCD:all");
1384  bool inelastic = settings.flag("SoftQCD:inelastic");
1385  if (softQCD || inelastic || settings.flag("SoftQCD:nonDiffractive")) {
1386  sigmaPtr = new Sigma0nonDiffractive;
1387  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1388  }
1389  if (softQCD || settings.flag("SoftQCD:elastic")) {
1390  sigmaPtr = new Sigma0AB2AB;
1391  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1392  }
1393  if (softQCD || inelastic || settings.flag("SoftQCD:singleDiffractive")) {
1394  sigmaPtr = new Sigma0AB2XB;
1395  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1396  sigmaPtr = new Sigma0AB2AX;
1397  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1398  }
1399  if (softQCD || inelastic || settings.flag("SoftQCD:doubleDiffractive")) {
1400  sigmaPtr = new Sigma0AB2XX;
1401  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1402  }
1403  if (softQCD || inelastic || settings.flag("SoftQCD:centralDiffractive")) {
1404  sigmaPtr = new Sigma0AB2AXB;
1405  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1406  }
1407 
1408  // Set up requested objects for hard QCD processes.
1409  bool hardQCD = settings.flag("HardQCD:all");
1410  if (hardQCD || settings.flag("HardQCD:gg2gg")) {
1411  sigmaPtr = new Sigma2gg2gg;
1412  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1413  }
1414  if (hardQCD || settings.flag("HardQCD:gg2qqbar")) {
1415  sigmaPtr = new Sigma2gg2qqbar;
1416  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1417  }
1418  if (hardQCD || settings.flag("HardQCD:qg2qg")) {
1419  sigmaPtr = new Sigma2qg2qg;
1420  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1421  }
1422  if (hardQCD || settings.flag("HardQCD:qq2qq")) {
1423  sigmaPtr = new Sigma2qq2qq;
1424  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1425  }
1426  if (hardQCD || settings.flag("HardQCD:qqbar2gg")) {
1427  sigmaPtr = new Sigma2qqbar2gg;
1428  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1429  }
1430  if (hardQCD || settings.flag("HardQCD:qqbar2qqbarNew")) {
1431  sigmaPtr = new Sigma2qqbar2qqbarNew;
1432  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1433  }
1434 
1435  // Set up requested objects for c cbar and b bbar, also hard QCD.
1436  bool hardccbar = settings.flag("HardQCD:hardccbar");
1437  if (hardQCD || hardccbar || settings.flag("HardQCD:gg2ccbar")) {
1438  sigmaPtr = new Sigma2gg2QQbar(4, 121);
1439  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1440  }
1441  if (hardQCD || hardccbar || settings.flag("HardQCD:qqbar2ccbar")) {
1442  sigmaPtr = new Sigma2qqbar2QQbar(4, 122);
1443  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1444  }
1445  bool hardbbbar = settings.flag("HardQCD:hardbbbar");
1446  if (hardQCD || hardbbbar || settings.flag("HardQCD:gg2bbbar")) {
1447  sigmaPtr = new Sigma2gg2QQbar(5, 123);
1448  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1449  }
1450  if (hardQCD || hardbbbar || settings.flag("HardQCD:qqbar2bbbar")) {
1451  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
1452  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1453  }
1454 
1455  // Set up requested objects for hard QCD 2 -> 3 processes.
1456  bool hardQCD3parton = settings.flag("HardQCD:3parton");
1457  if (hardQCD3parton || settings.flag("HardQCD:gg2ggg")) {
1458  sigmaPtr = new Sigma3gg2ggg;
1459  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1460  }
1461  if (hardQCD3parton || settings.flag("HardQCD:qqbar2ggg")) {
1462  sigmaPtr = new Sigma3qqbar2ggg;
1463  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1464  }
1465  if (hardQCD3parton || settings.flag("HardQCD:qg2qgg")) {
1466  sigmaPtr = new Sigma3qg2qgg;
1467  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1468  }
1469  if (hardQCD3parton || settings.flag("HardQCD:qq2qqgDiff")) {
1470  sigmaPtr = new Sigma3qq2qqgDiff;
1471  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1472  }
1473  if (hardQCD3parton || settings.flag("HardQCD:qq2qqgSame")) {
1474  sigmaPtr = new Sigma3qq2qqgSame;
1475  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1476  }
1477  if (hardQCD3parton || settings.flag("HardQCD:qqbar2qqbargDiff")) {
1478  sigmaPtr = new Sigma3qqbar2qqbargDiff;
1479  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1480  }
1481  if (hardQCD3parton || settings.flag("HardQCD:qqbar2qqbargSame")) {
1482  sigmaPtr = new Sigma3qqbar2qqbargSame;
1483  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1484  }
1485  if (hardQCD3parton || settings.flag("HardQCD:gg2qqbarg")) {
1486  sigmaPtr = new Sigma3gg2qqbarg;
1487  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1488  }
1489  if (hardQCD3parton || settings.flag("HardQCD:qg2qqqbarDiff")) {
1490  sigmaPtr = new Sigma3qg2qqqbarDiff;
1491  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1492  }
1493  if (hardQCD3parton || settings.flag("HardQCD:qg2qqqbarSame")) {
1494  sigmaPtr = new Sigma3qg2qqqbarSame;
1495  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1496  }
1497 
1498  // Set up requested objects for prompt photon processes.
1499  bool promptPhotons = settings.flag("PromptPhoton:all");
1500  if (promptPhotons
1501  || settings.flag("PromptPhoton:qg2qgamma")) {
1502  sigmaPtr = new Sigma2qg2qgamma;
1503  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1504  }
1505  if (promptPhotons
1506  || settings.flag("PromptPhoton:qqbar2ggamma")) {
1507  sigmaPtr = new Sigma2qqbar2ggamma;
1508  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1509  }
1510  if (promptPhotons
1511  || settings.flag("PromptPhoton:gg2ggamma")) {
1512  sigmaPtr = new Sigma2gg2ggamma;
1513  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1514  }
1515  if (promptPhotons
1516  || settings.flag("PromptPhoton:ffbar2gammagamma")) {
1517  sigmaPtr = new Sigma2ffbar2gammagamma;
1518  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1519  }
1520  if (promptPhotons
1521  || settings.flag("PromptPhoton:gg2gammagamma")) {
1522  sigmaPtr = new Sigma2gg2gammagamma;
1523  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1524  }
1525 
1526  // Set up requested objects for weak gauge boson t-channel exchange.
1527  bool weakBosonExchanges = settings.flag("WeakBosonExchange:all");
1528  if (weakBosonExchanges
1529  || settings.flag("WeakBosonExchange:ff2ff(t:gmZ)")) {
1530  sigmaPtr = new Sigma2ff2fftgmZ;
1531  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1532  }
1533  if (weakBosonExchanges
1534  || settings.flag("WeakBosonExchange:ff2ff(t:W)")) {
1535  sigmaPtr = new Sigma2ff2fftW;
1536  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1537  }
1538 
1539  // Set up requested objects for weak gauge boson processes.
1540  bool weakSingleBosons = settings.flag("WeakSingleBoson:all");
1541  if (weakSingleBosons
1542  || settings.flag("WeakSingleBoson:ffbar2gmZ")) {
1543  sigmaPtr = new Sigma1ffbar2gmZ;
1544  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1545  }
1546  if (weakSingleBosons
1547  || settings.flag("WeakSingleBoson:ffbar2W")) {
1548  sigmaPtr = new Sigma1ffbar2W;
1549  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1550  }
1551 
1552  // Set up requested object for s-channel gamma exchange.
1553  // Subset of gamma*/Z0 above, intended for multiparton interactions.
1554  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:gm)")) {
1555  sigmaPtr = new Sigma2ffbar2ffbarsgm;
1556  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1557  }
1558 
1559  // Set up requested object for s-channel gamma*/Z0 or W+- exchange.
1560  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:gmZ)")) {
1561  sigmaPtr = new Sigma2ffbar2ffbarsgmZ;
1562  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1563  }
1564  if (settings.flag("WeakSingleBoson:ffbar2ffbar(s:W)")) {
1565  sigmaPtr = new Sigma2ffbar2ffbarsW;
1566  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1567  }
1568 
1569  // Set up requested objects for weak gauge boson pair processes.
1570  bool weakDoubleBosons = settings.flag("WeakDoubleBoson:all");
1571  if (weakDoubleBosons
1572  || settings.flag("WeakDoubleBoson:ffbar2gmZgmZ")) {
1573  sigmaPtr = new Sigma2ffbar2gmZgmZ;
1574  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1575  }
1576  if (weakDoubleBosons
1577  || settings.flag("WeakDoubleBoson:ffbar2ZW")) {
1578  sigmaPtr = new Sigma2ffbar2ZW;
1579  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1580  }
1581  if (weakDoubleBosons
1582  || settings.flag("WeakDoubleBoson:ffbar2WW")) {
1583  sigmaPtr = new Sigma2ffbar2WW;
1584  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1585  }
1586 
1587  // Set up requested objects for weak gauge boson + parton processes.
1588  bool weakBosonAndPartons = settings.flag("WeakBosonAndParton:all");
1589  if (weakBosonAndPartons
1590  || settings.flag("WeakBosonAndParton:qqbar2gmZg")) {
1591  sigmaPtr = new Sigma2qqbar2gmZg;
1592  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1593  }
1594  if (weakBosonAndPartons
1595  || settings.flag("WeakBosonAndParton:qg2gmZq")) {
1596  sigmaPtr = new Sigma2qg2gmZq;
1597  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1598  }
1599  if (weakBosonAndPartons
1600  || settings.flag("WeakBosonAndParton:ffbar2gmZgm")) {
1601  sigmaPtr = new Sigma2ffbar2gmZgm;
1602  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1603  }
1604  if (weakBosonAndPartons
1605  || settings.flag("WeakBosonAndParton:fgm2gmZf")) {
1606  sigmaPtr = new Sigma2fgm2gmZf;
1607  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1608  }
1609  if (weakBosonAndPartons
1610  || settings.flag("WeakBosonAndParton:qqbar2Wg")) {
1611  sigmaPtr = new Sigma2qqbar2Wg;
1612  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1613  }
1614  if (weakBosonAndPartons
1615  || settings.flag("WeakBosonAndParton:qg2Wq")) {
1616  sigmaPtr = new Sigma2qg2Wq;
1617  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1618  }
1619  if (weakBosonAndPartons
1620  || settings.flag("WeakBosonAndParton:ffbar2Wgm")) {
1621  sigmaPtr = new Sigma2ffbar2Wgm;
1622  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1623  }
1624  if (weakBosonAndPartons
1625  || settings.flag("WeakBosonAndParton:fgm2Wf")) {
1626  sigmaPtr = new Sigma2fgm2Wf;
1627  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1628  }
1629 
1630  // Read in settings related to photon beams.
1631  bool photonParton = settings.flag("PhotonParton:all");
1632  int photonMode = settings.mode("Photon:ProcessType");
1633  bool gammaA = settings.flag("PDF:beamA2gamma")
1634  || (abs(infoPtr->idA()) == 22);
1635  bool gammaB = settings.flag("PDF:beamB2gamma")
1636  || (abs(infoPtr->idB()) == 22);
1637 
1638  // Initialize process for photon on side A only when explicitly direct
1639  // photons requested for side A or when a ''normal'' hadron on side A and
1640  // a resolved photon on side B as important to use the proper flux/PDF.
1641  // Otherwise initialize only one process (with the photon on side B).
1642  bool initGammaA = ( (photonMode == 0) && (gammaA || gammaB) )
1643  || ( (photonMode == 3) && gammaA )
1644  || ( (photonMode == 1) && (gammaB && !gammaA) );
1645  bool initGammaB = (photonMode == 0)
1646  || ( (photonMode == 2) && gammaB )
1647  || ( (photonMode == 1) && (gammaA && !gammaB) )
1648  || (!gammaA && !gammaB);
1649 
1650  // Initialize gm+gm processes only when explicitly needed with photon beams
1651  // or when normal hadronic collisions.
1652  bool initGammaGamma = ( (initGammaA || !gammaA) && (initGammaB || !gammaB) )
1653  || (photonMode == 4 && gammaA && gammaB );
1654 
1655  // Set up requested objects for photon collision processes.
1656  bool photonCollisions = settings.flag("PhotonCollision:all");
1657  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2qqbar") )
1658  && initGammaGamma ) {
1659  sigmaPtr = new Sigma2gmgm2ffbar(1, 261);
1660  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1661  }
1662  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2ccbar") )
1663  && initGammaGamma ) {
1664  sigmaPtr = new Sigma2gmgm2ffbar(4, 262);
1665  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1666  }
1667  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2bbbar") )
1668  && initGammaGamma ) {
1669  sigmaPtr = new Sigma2gmgm2ffbar(5, 263);
1670  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1671  }
1672  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2ee") )
1673  && initGammaGamma ) {
1674  sigmaPtr = new Sigma2gmgm2ffbar(11, 264);
1675  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1676  }
1677  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2mumu") )
1678  && initGammaGamma ) {
1679  sigmaPtr = new Sigma2gmgm2ffbar(13, 265);
1680  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1681  }
1682  if ( ( photonCollisions || settings.flag("PhotonCollision:gmgm2tautau") )
1683  && initGammaGamma ) {
1684  sigmaPtr = new Sigma2gmgm2ffbar(15, 266);
1685  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1686  }
1687 
1688  // Set up requested objects for photon-parton processes.
1689  // In case of unresolved photons insert process with correct initiator.
1690  if (photonParton || settings.flag("PhotonParton:ggm2qqbar")) {
1691  if ( initGammaB ) {
1692  sigmaPtr = new Sigma2ggm2qqbar(1, 271, "ggm");
1693  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1694  }
1695  if ( initGammaA ) {
1696  sigmaPtr = new Sigma2ggm2qqbar(1, 281, "gmg");
1697  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1698  }
1699  }
1700  if (photonParton || settings.flag("PhotonParton:ggm2ccbar")) {
1701  if ( initGammaB ) {
1702  sigmaPtr = new Sigma2ggm2qqbar(4, 272, "ggm");
1703  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1704  }
1705  if ( initGammaA) {
1706  sigmaPtr = new Sigma2ggm2qqbar(4, 282, "gmg");
1707  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1708  }
1709  }
1710  if (photonParton || settings.flag("PhotonParton:ggm2bbbar")) {
1711  if ( initGammaB ) {
1712  sigmaPtr = new Sigma2ggm2qqbar(5, 273, "ggm");
1713  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1714  }
1715  if ( initGammaA ) {
1716  sigmaPtr = new Sigma2ggm2qqbar(5, 283, "gmg");
1717  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1718  }
1719  }
1720  if (photonParton || settings.flag("PhotonParton:qgm2qg")) {
1721  if ( initGammaB ) {
1722  sigmaPtr = new Sigma2qgm2qg(274, "qgm");
1723  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1724  }
1725  if ( initGammaA ) {
1726  sigmaPtr = new Sigma2qgm2qg(284, "gmq");
1727  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1728  }
1729  }
1730  if (settings.flag("PhotonParton:qgm2qgm")) {
1731  if ( initGammaB ) {
1732  sigmaPtr = new Sigma2qgm2qgm(275, "qgm");
1733  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1734  }
1735  if ( initGammaA ) {
1736  sigmaPtr = new Sigma2qgm2qgm(285, "gmq");
1737  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1738  }
1739  }
1740 
1741  // Set up requested objects for onia production.
1742  charmonium = SigmaOniaSetup(infoPtr, 4);
1743  bottomonium = SigmaOniaSetup(infoPtr, 5);
1744  vector<SigmaProcess*> charmoniumSigmaPtrs, bottomoniumSigmaPtrs;
1745  charmonium.setupSigma2gg(charmoniumSigmaPtrs);
1746  charmonium.setupSigma2qg(charmoniumSigmaPtrs);
1747  charmonium.setupSigma2qq(charmoniumSigmaPtrs);
1748  charmonium.setupSigma2dbl(charmoniumSigmaPtrs);
1749  bottomonium.setupSigma2gg(bottomoniumSigmaPtrs);
1750  bottomonium.setupSigma2qg(bottomoniumSigmaPtrs);
1751  bottomonium.setupSigma2qq(bottomoniumSigmaPtrs);
1752  bottomonium.setupSigma2dbl(bottomoniumSigmaPtrs);
1753  for (unsigned int i = 0; i < charmoniumSigmaPtrs.size(); ++i)
1754  containerPtrs.push_back( new ProcessContainer(charmoniumSigmaPtrs[i]) );
1755  for (unsigned int i = 0; i < bottomoniumSigmaPtrs.size(); ++i)
1756  containerPtrs.push_back( new ProcessContainer(bottomoniumSigmaPtrs[i]) );
1757 
1758  // Set up requested objects for top production.
1759  bool tops = settings.flag("Top:all");
1760  if (tops || settings.flag("Top:gg2ttbar")) {
1761  sigmaPtr = new Sigma2gg2QQbar(6, 601);
1762  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1763  }
1764  if (tops || settings.flag("Top:qqbar2ttbar")) {
1765  sigmaPtr = new Sigma2qqbar2QQbar(6, 602);
1766  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1767  }
1768  if (tops || settings.flag("Top:qq2tq(t:W)")) {
1769  sigmaPtr = new Sigma2qq2QqtW(6, 603);
1770  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1771  }
1772  if (tops || settings.flag("Top:ffbar2ttbar(s:gmZ)")) {
1773  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(6, 604);
1774  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1775  }
1776  if (tops || settings.flag("Top:ffbar2tqbar(s:W)")) {
1777  sigmaPtr = new Sigma2ffbar2FfbarsW(6, 0, 605);
1778  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1779  }
1780  if (tops || settings.flag("Top:gmgm2ttbar")) {
1781  sigmaPtr = new Sigma2gmgm2ffbar(6, 606);
1782  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1783  }
1784  if (tops || settings.flag("Top:ggm2ttbar")) {
1785  if ( initGammaB ) {
1786  sigmaPtr = new Sigma2ggm2qqbar(6, 607, "ggm");
1787  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1788  }
1789  if ( initGammaA ) {
1790  sigmaPtr = new Sigma2ggm2qqbar(6, 617, "gmg");
1791  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1792  }
1793  }
1794 
1795  // Set up requested objects for fourth-generation b' production
1796  bool bPrimes = settings.flag("FourthBottom:all");
1797  if (bPrimes || settings.flag("FourthBottom:gg2bPrimebPrimebar")) {
1798  sigmaPtr = new Sigma2gg2QQbar(7, 801);
1799  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1800  }
1801  if (bPrimes || settings.flag("FourthBottom:qqbar2bPrimebPrimebar")) {
1802  sigmaPtr = new Sigma2qqbar2QQbar(7, 802);
1803  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1804  }
1805  if (bPrimes || settings.flag("FourthBottom:qq2bPrimeq(t:W)")) {
1806  sigmaPtr = new Sigma2qq2QqtW(7, 803);
1807  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1808  }
1809  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimebPrimebar(s:gmZ)")) {
1810  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(7, 804);
1811  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1812  }
1813  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimeqbar(s:W)")) {
1814  sigmaPtr = new Sigma2ffbar2FfbarsW(7, 0, 805);
1815  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1816  }
1817  if (bPrimes || settings.flag("FourthBottom:ffbar2bPrimetbar(s:W)")) {
1818  sigmaPtr = new Sigma2ffbar2FfbarsW(7, 6, 806);
1819  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1820  }
1821 
1822  // Set up requested objects for fourth-generation t' production
1823  bool tPrimes = settings.flag("FourthTop:all");
1824  if (tPrimes || settings.flag("FourthTop:gg2tPrimetPrimebar")) {
1825  sigmaPtr = new Sigma2gg2QQbar(8, 821);
1826  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1827  }
1828  if (tPrimes || settings.flag("FourthTop:qqbar2tPrimetPrimebar")) {
1829  sigmaPtr = new Sigma2qqbar2QQbar(8, 822);
1830  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1831  }
1832  if (tPrimes || settings.flag("FourthTop:qq2tPrimeq(t:W)")) {
1833  sigmaPtr = new Sigma2qq2QqtW(8, 823);
1834  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1835  }
1836  if (tPrimes || settings.flag("FourthTop:ffbar2tPrimetPrimebar(s:gmZ)")) {
1837  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(8, 824);
1838  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1839  }
1840  if (tPrimes || settings.flag("FourthTop:ffbar2tPrimeqbar(s:W)")) {
1841  sigmaPtr = new Sigma2ffbar2FfbarsW(8, 0, 825);
1842  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1843  }
1844 
1845  // Set up requested objects for two different fourth-generation fermions.
1846  if (bPrimes || tPrimes
1847  || settings.flag("FourthPair:ffbar2tPrimebPrimebar(s:W)")) {
1848  sigmaPtr = new Sigma2ffbar2FfbarsW(8, 7, 841);
1849  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1850  }
1851  if (settings.flag("FourthPair:ffbar2tauPrimenuPrimebar(s:W)")) {
1852  sigmaPtr = new Sigma2ffbar2FfbarsW(17, 18, 842);
1853  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1854  }
1855 
1856  // Flag for global choice between SM and BSM Higgses.
1857  bool useBSMHiggses = settings.flag("Higgs:useBSM");
1858 
1859  // Set up requested objects for Standard-Model Higgs production.
1860  if (!useBSMHiggses) {
1861  bool HiggsesSM = settings.flag("HiggsSM:all");
1862  if (HiggsesSM || settings.flag("HiggsSM:ffbar2H")) {
1863  sigmaPtr = new Sigma1ffbar2H(0);
1864  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1865  }
1866  if (HiggsesSM || settings.flag("HiggsSM:gg2H")) {
1867  sigmaPtr = new Sigma1gg2H(0);
1868  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1869  }
1870  if (HiggsesSM || settings.flag("HiggsSM:gmgm2H")) {
1871  sigmaPtr = new Sigma1gmgm2H(0);
1872  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1873  }
1874  if (HiggsesSM || settings.flag("HiggsSM:ffbar2HZ")) {
1875  sigmaPtr = new Sigma2ffbar2HZ(0);
1876  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1877  }
1878  if (HiggsesSM || settings.flag("HiggsSM:ffbar2HW")) {
1879  sigmaPtr = new Sigma2ffbar2HW(0);
1880  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1881  }
1882  if (HiggsesSM || settings.flag("HiggsSM:ff2Hff(t:ZZ)")) {
1883  sigmaPtr = new Sigma3ff2HfftZZ(0);
1884  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1885  }
1886  if (HiggsesSM || settings.flag("HiggsSM:ff2Hff(t:WW)")) {
1887  sigmaPtr = new Sigma3ff2HfftWW(0);
1888  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1889  }
1890  if (HiggsesSM || settings.flag("HiggsSM:gg2Httbar")) {
1891  sigmaPtr = new Sigma3gg2HQQbar(6,0);
1892  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1893  }
1894  if (HiggsesSM || settings.flag("HiggsSM:qqbar2Httbar")) {
1895  sigmaPtr = new Sigma3qqbar2HQQbar(6,0);
1896  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1897  }
1898 
1899  // Further Standard-Model Higgs processes, not included in "all".
1900  if (settings.flag("HiggsSM:qg2Hq")) {
1901  sigmaPtr = new Sigma2qg2Hq(4,0);
1902  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1903  sigmaPtr = new Sigma2qg2Hq(5,0);
1904  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1905  }
1906  if (settings.flag("HiggsSM:gg2Hbbbar")) {
1907  sigmaPtr = new Sigma3gg2HQQbar(5,0);
1908  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1909  }
1910  if (settings.flag("HiggsSM:qqbar2Hbbbar")) {
1911  sigmaPtr = new Sigma3qqbar2HQQbar(5,0);
1912  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1913  }
1914  if (settings.flag("HiggsSM:gg2Hg(l:t)")) {
1915  sigmaPtr = new Sigma2gg2Hglt(0);
1916  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1917  }
1918  if (settings.flag("HiggsSM:qg2Hq(l:t)")) {
1919  sigmaPtr = new Sigma2qg2Hqlt(0);
1920  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1921  }
1922  if (settings.flag("HiggsSM:qqbar2Hg(l:t)")) {
1923  sigmaPtr = new Sigma2qqbar2Hglt(0);
1924  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1925  }
1926  }
1927 
1928  // Common switch for the group of Higgs production BSM.
1929  if (useBSMHiggses) {
1930  bool HiggsesBSM = settings.flag("HiggsBSM:all");
1931 
1932  // Set up requested objects for BSM H1 production.
1933  bool HiggsesH1 = settings.flag("HiggsBSM:allH1");
1934  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1")) {
1935  sigmaPtr = new Sigma1ffbar2H(1);
1936  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1937  }
1938  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gg2H1")) {
1939  sigmaPtr = new Sigma1gg2H(1);
1940  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1941  }
1942  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gmgm2H1")) {
1943  sigmaPtr = new Sigma1gmgm2H(1);
1944  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1945  }
1946  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1Z")) {
1947  sigmaPtr = new Sigma2ffbar2HZ(1);
1948  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1949  }
1950  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ffbar2H1W")) {
1951  sigmaPtr = new Sigma2ffbar2HW(1);
1952  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1953  }
1954  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ff2H1ff(t:ZZ)")) {
1955  sigmaPtr = new Sigma3ff2HfftZZ(1);
1956  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1957  }
1958  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:ff2H1ff(t:WW)")) {
1959  sigmaPtr = new Sigma3ff2HfftWW(1);
1960  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1961  }
1962  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:gg2H1ttbar")) {
1963  sigmaPtr = new Sigma3gg2HQQbar(6,1);
1964  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1965  }
1966  if (HiggsesBSM || HiggsesH1 || settings.flag("HiggsBSM:qqbar2H1ttbar")) {
1967  sigmaPtr = new Sigma3qqbar2HQQbar(6,1);
1968  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1969  }
1970 
1971  // Further BSM H1 processes, not included in "all".
1972  if (settings.flag("HiggsBSM:qg2H1q")) {
1973  sigmaPtr = new Sigma2qg2Hq(4,1);
1974  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1975  sigmaPtr = new Sigma2qg2Hq(5,1);
1976  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1977  }
1978  if (settings.flag("HiggsBSM:gg2H1bbbar")) {
1979  sigmaPtr = new Sigma3gg2HQQbar(5,1);
1980  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1981  }
1982  if (settings.flag("HiggsBSM:qqbar2H1bbbar")) {
1983  sigmaPtr = new Sigma3qqbar2HQQbar(5,1);
1984  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1985  }
1986  if (settings.flag("HiggsBSM:gg2H1g(l:t)")) {
1987  sigmaPtr = new Sigma2gg2Hglt(1);
1988  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1989  }
1990  if (settings.flag("HiggsBSM:qg2H1q(l:t)")) {
1991  sigmaPtr = new Sigma2qg2Hqlt(1);
1992  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1993  }
1994  if (settings.flag("HiggsBSM:qqbar2H1g(l:t)")) {
1995  sigmaPtr = new Sigma2qqbar2Hglt(1);
1996  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
1997  }
1998 
1999  // Set up requested objects for BSM H2 production.
2000  bool HiggsesH2 = settings.flag("HiggsBSM:allH2");
2001  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2")) {
2002  sigmaPtr = new Sigma1ffbar2H(2);
2003  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2004  }
2005  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gg2H2")) {
2006  sigmaPtr = new Sigma1gg2H(2);
2007  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2008  }
2009  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gmgm2H2")) {
2010  sigmaPtr = new Sigma1gmgm2H(2);
2011  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2012  }
2013  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2Z")) {
2014  sigmaPtr = new Sigma2ffbar2HZ(2);
2015  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2016  }
2017  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ffbar2H2W")) {
2018  sigmaPtr = new Sigma2ffbar2HW(2);
2019  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2020  }
2021  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ff2H2ff(t:ZZ)")) {
2022  sigmaPtr = new Sigma3ff2HfftZZ(2);
2023  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2024  }
2025  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:ff2H2ff(t:WW)")) {
2026  sigmaPtr = new Sigma3ff2HfftWW(2);
2027  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2028  }
2029  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:gg2H2ttbar")) {
2030  sigmaPtr = new Sigma3gg2HQQbar(6,2);
2031  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2032  }
2033  if (HiggsesBSM || HiggsesH2 || settings.flag("HiggsBSM:qqbar2H2ttbar")) {
2034  sigmaPtr = new Sigma3qqbar2HQQbar(6,2);
2035  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2036  }
2037 
2038  // Further BSM H2 processes, not included in "all".
2039  if (settings.flag("HiggsBSM:qg2H2q")) {
2040  sigmaPtr = new Sigma2qg2Hq(4,2);
2041  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2042  sigmaPtr = new Sigma2qg2Hq(5,2);
2043  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2044  }
2045  if (settings.flag("HiggsBSM:gg2H2bbbar")) {
2046  sigmaPtr = new Sigma3gg2HQQbar(5,2);
2047  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2048  }
2049  if (settings.flag("HiggsBSM:qqbar2H2bbbar")) {
2050  sigmaPtr = new Sigma3qqbar2HQQbar(5,2);
2051  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2052  }
2053  if (settings.flag("HiggsBSM:gg2H2g(l:t)")) {
2054  sigmaPtr = new Sigma2gg2Hglt(2);
2055  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2056  }
2057  if (settings.flag("HiggsBSM:qg2H2q(l:t)")) {
2058  sigmaPtr = new Sigma2qg2Hqlt(2);
2059  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2060  }
2061  if (settings.flag("HiggsBSM:qqbar2H2g(l:t)")) {
2062  sigmaPtr = new Sigma2qqbar2Hglt(2);
2063  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2064  }
2065 
2066  // Set up requested objects for BSM A3 production.
2067  bool HiggsesA3 = settings.flag("HiggsBSM:allA3");
2068  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3")) {
2069  sigmaPtr = new Sigma1ffbar2H(3);
2070  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2071  }
2072  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gg2A3")) {
2073  sigmaPtr = new Sigma1gg2H(3);
2074  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2075  }
2076  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gmgm2A3")) {
2077  sigmaPtr = new Sigma1gmgm2H(3);
2078  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2079  }
2080  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3Z")) {
2081  sigmaPtr = new Sigma2ffbar2HZ(3);
2082  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2083  }
2084  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ffbar2A3W")) {
2085  sigmaPtr = new Sigma2ffbar2HW(3);
2086  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2087  }
2088  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ff2A3ff(t:ZZ)")) {
2089  sigmaPtr = new Sigma3ff2HfftZZ(3);
2090  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2091  }
2092  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:ff2A3ff(t:WW)")) {
2093  sigmaPtr = new Sigma3ff2HfftWW(3);
2094  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2095  }
2096  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:gg2A3ttbar")) {
2097  sigmaPtr = new Sigma3gg2HQQbar(6,3);
2098  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2099  }
2100  if (HiggsesBSM || HiggsesA3 || settings.flag("HiggsBSM:qqbar2A3ttbar")) {
2101  sigmaPtr = new Sigma3qqbar2HQQbar(6,3);
2102  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2103  }
2104 
2105  // Further BSM A3 processes, not included in "all".
2106  if (settings.flag("HiggsBSM:qg2A3q")) {
2107  sigmaPtr = new Sigma2qg2Hq(4,3);
2108  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2109  sigmaPtr = new Sigma2qg2Hq(5,3);
2110  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2111  }
2112  if (settings.flag("HiggsBSM:gg2A3bbbar")) {
2113  sigmaPtr = new Sigma3gg2HQQbar(5,3);
2114  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2115  }
2116  if (settings.flag("HiggsBSM:qqbar2A3bbbar")) {
2117  sigmaPtr = new Sigma3qqbar2HQQbar(5,3);
2118  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2119  }
2120  if (settings.flag("HiggsBSM:gg2A3g(l:t)")) {
2121  sigmaPtr = new Sigma2gg2Hglt(3);
2122  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2123  }
2124  if (settings.flag("HiggsBSM:qg2A3q(l:t)")) {
2125  sigmaPtr = new Sigma2qg2Hqlt(3);
2126  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2127  }
2128  if (settings.flag("HiggsBSM:qqbar2A3g(l:t)")) {
2129  sigmaPtr = new Sigma2qqbar2Hglt(3);
2130  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2131  }
2132 
2133  // Set up requested objects for Charged Higgs production
2134  bool HiggsesChg = settings.flag("HiggsBSM:allH+-");
2135  if (HiggsesBSM || HiggsesChg || settings.flag("HiggsBSM:ffbar2H+-")) {
2136  sigmaPtr = new Sigma1ffbar2Hchg;
2137  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2138  }
2139  if (HiggsesBSM || HiggsesChg || settings.flag("HiggsBSM:bg2H+-t")) {
2140  sigmaPtr = new Sigma2qg2Hchgq(6, 1062, "b g -> H+- t");
2141  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2142  }
2143 
2144  // Set up requested objects for Higgs pair-production
2145  bool HiggsesPairs = settings.flag("HiggsBSM:allHpair");
2146  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2A3H1")) {
2147  sigmaPtr = new Sigma2ffbar2A3H12(1);
2148  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2149  }
2150  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2A3H2")) {
2151  sigmaPtr = new Sigma2ffbar2A3H12(2);
2152  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2153  }
2154  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+-H1")) {
2155  sigmaPtr = new Sigma2ffbar2HchgH12(1);
2156  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2157  }
2158  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+-H2")) {
2159  sigmaPtr = new Sigma2ffbar2HchgH12(2);
2160  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2161  }
2162  if (HiggsesBSM || HiggsesPairs || settings.flag("HiggsBSM:ffbar2H+H-")) {
2163  sigmaPtr = new Sigma2ffbar2HposHneg();
2164  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2165  }
2166  }
2167 
2168  // Set up requested objects for SUSY pair processes.
2169  CoupSUSY* coupSUSYPtr = infoPtr->coupSUSYPtr;
2170  if (coupSUSYPtr->isSUSY) {
2171 
2172  bool SUSYs = settings.flag("SUSY:all");
2173  bool nmssm = settings.flag("SLHA:NMSSM");
2174 
2175  // Preselected SUSY codes.
2176  setupIdVecs( settings);
2177 
2178  // MSSM: 4 neutralinos
2179  int nNeut = 4;
2180  if (nmssm) nNeut = 5;
2181 
2182  // Gluino-gluino
2183  if (SUSYs || settings.flag("SUSY:gg2gluinogluino")) {
2184  // Skip if outgoing codes not asked for
2185  if (allowIdVals( 1000021, 1000021)) {
2186  sigmaPtr = new Sigma2gg2gluinogluino();
2187  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2188  }
2189  }
2190  if (SUSYs || settings.flag("SUSY:qqbar2gluinogluino")) {
2191  // Skip if outgoing codes not asked for
2192  if (allowIdVals( 1000021, 1000021)) {
2193  sigmaPtr = new Sigma2qqbar2gluinogluino();
2194  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2195  }
2196  }
2197 
2198  // Gluino-squark
2199  if (SUSYs || settings.flag("SUSY:qg2squarkgluino")) {
2200  int iproc = 1202;
2201  for (int idx = 1; idx <= 6; ++idx) {
2202  for (int iso = 1; iso <= 2; ++iso) {
2203  iproc += 2;
2204  int id3 = iso + ((idx <= 3)
2205  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2206  int id4 = 1000021;
2207  // Skip if outgoing codes not asked for
2208  if (!allowIdVals( id3, id4)) continue;
2209  sigmaPtr = new Sigma2qg2squarkgluino(id3,iproc-1);
2210  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2211  sigmaPtr = new Sigma2qg2squarkgluino(-id3,iproc);
2212  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2213  }
2214  }
2215  }
2216 
2217  // Squark-antisquark (gg initiated)
2218  if (SUSYs || settings.flag("SUSY:gg2squarkantisquark")) {
2219  int iproc = 1214;
2220  for (int idx = 1; idx <= 6; ++idx) {
2221  for (int iso = 1; iso <= 2; ++iso) {
2222  iproc++;
2223  int id = iso + ((idx <= 3)
2224  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2225  // Skip if outgoing codes not asked for
2226  if (!allowIdVals( id, id)) continue;
2227  sigmaPtr = new Sigma2gg2squarkantisquark(id,iproc);
2228  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2229  }
2230  }
2231  }
2232 
2233  // Squark-antisquark (qqbar initiated)
2234  if (SUSYs || settings.flag("SUSY:qqbar2squarkantisquark")) {
2235  int iproc = 1230;
2236  for (int idx = 1; idx <= 6; ++idx) {
2237  for (int iso = 1; iso <= 2; ++iso) {
2238  for (int jso = iso; jso >= 1; --jso) {
2239  for (int jdx = 1; jdx <= 6; ++jdx) {
2240  if (iso == jso && jdx < idx) continue;
2241  int id1 = iso + ((idx <= 3) ? 1000000+2*(idx-1)
2242  : 2000000+2*(idx-4));
2243  int id2 = jso + ((jdx <= 3) ? 1000000+2*(jdx-1)
2244  : 2000000+2*(jdx-4));
2245  // Update process number counter (for ~q~q, +2 if not self-conj)
2246  //if (iproc == 1302) iproc=1310;
2247  iproc++;
2248  if (iso == jso && id1 != id2) iproc++;
2249  // Skip if outgoing codes not asked for
2250  if (!allowIdVals( id1, id2)) continue;
2251  if (iso == jso && id1 != id2) {
2252  sigmaPtr = new Sigma2qqbar2squarkantisquark(id1,-id2,iproc-1);
2253  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2254  sigmaPtr = new Sigma2qqbar2squarkantisquark(id2,-id1,iproc);
2255  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2256  } else {
2257  sigmaPtr = new Sigma2qqbar2squarkantisquark(id1,-id2,iproc);
2258  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2259  }
2260  }
2261  }
2262  }
2263  }
2264  }
2265 
2266  // Squark-squark
2267  if (SUSYs || settings.flag("SUSY:qq2squarksquark")) {
2268  int iproc = 1350;
2269  for (int idx = 1; idx <= 6; ++idx) {
2270  for (int iso = 1; iso <= 2; ++iso) {
2271  for (int jso = iso; jso >= 1; jso--) {
2272  for (int jdx = 1; jdx <= 6; ++jdx) {
2273  if (iso == jso && jdx < idx) continue;
2274  iproc++;
2275  int id1 = iso + ((idx <= 3)
2276  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2277  int id2 = jso + ((jdx <= 3)
2278  ? 1000000+2*(jdx-1) : 2000000+2*(jdx-4));
2279  // Skip if outgoing codes not asked for.
2280  if (!allowIdVals( id1, id2)) continue;
2281  sigmaPtr = new Sigma2qq2squarksquark(id1,id2,iproc);
2282  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2283  }
2284  }
2285  }
2286  }
2287  }
2288 
2289  // Neutralino + squark
2290  if (SUSYs || settings.flag("SUSY:qg2chi0squark")) {
2291  int iproc = 1430;
2292  for (int iNeut= 1; iNeut <= nNeut; iNeut++) {
2293  for (int idx = 1; idx <= 6; idx++) {
2294  bool isUp = false;
2295  for (int iso = 1; iso <= 2; iso++) {
2296  if (iso == 2) isUp = true;
2297  iproc++;
2298  int id3 = coupSUSYPtr->idNeut(iNeut);
2299  int id4 = iso + ((idx <= 3)
2300  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2301  // Skip if outgoing codes not asked for
2302  if (!allowIdVals( id3, id4)) continue;
2303  sigmaPtr = new Sigma2qg2chi0squark(iNeut,idx,isUp,iproc);
2304  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2305  }
2306  }
2307  }
2308  }
2309 
2310  // Chargino + squark
2311  if (SUSYs || settings.flag("SUSY:qg2chi+-squark")) {
2312  int iproc = 1490;
2313  for (int iChar = 1; iChar <= 2; iChar++) {
2314  for (int idx = 1; idx <= 6; idx++) {
2315  bool isUp = false;
2316  for (int iso = 1; iso <= 2; iso++) {
2317  if (iso == 2) isUp = true;
2318  iproc++;
2319  int id3 = coupSUSYPtr->idChar(iChar);
2320  int id4 = iso + ((idx <= 3)
2321  ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2322  // Skip if outgoing codes not asked for
2323  if (!allowIdVals( id3, id4)) continue;
2324  sigmaPtr = new Sigma2qg2charsquark(iChar,idx,isUp,iproc);
2325  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2326  }
2327  }
2328  }
2329  }
2330 
2331  // Neutralino pairs
2332  if (SUSYs || settings.flag("SUSY:qqbar2chi0chi0")) {
2333  int iproc = 1550;
2334  for (int iNeut2 = 1; iNeut2 <= nNeut; iNeut2++) {
2335  for (int iNeut1 = 1; iNeut1 <= iNeut2; iNeut1++) {
2336  iproc++;
2337  // Skip if outgoing codes not asked for
2338  if (!allowIdVals( coupSUSYPtr->idNeut(iNeut1),
2339  coupSUSYPtr->idNeut(iNeut2) ) ) continue;
2340  sigmaPtr = new Sigma2qqbar2chi0chi0(iNeut1, iNeut2,iproc);
2341  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2342  }
2343  }
2344  }
2345 
2346  // Neutralino-Chargino
2347  if (SUSYs || settings.flag("SUSY:qqbar2chi+-chi0")) {
2348  int iproc = 1570;
2349  for (int iNeut = 1; iNeut <= nNeut; iNeut++) {
2350  for (int iChar = 1; iChar <= 2; ++iChar) {
2351  iproc += 2;
2352  // Skip if outgoing codes not asked for
2353  if (!allowIdVals( coupSUSYPtr->idNeut(iNeut),
2354  coupSUSYPtr->idChar(iChar) ) ) continue;
2355  sigmaPtr = new Sigma2qqbar2charchi0( iChar, iNeut, iproc-1);
2356  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2357  sigmaPtr = new Sigma2qqbar2charchi0(-iChar, iNeut, iproc);
2358  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2359  }
2360  }
2361  }
2362 
2363  // Chargino-Chargino
2364  if (SUSYs || settings.flag("SUSY:qqbar2chi+chi-")) {
2365  int iproc = 1590;
2366  for (int i = 1; i <= 2; ++i) {
2367  for (int j = 1; j <= 2; ++j) {
2368  iproc++;
2369  // Skip if outgoing codes not asked for
2370  if (!allowIdVals( coupSUSYPtr->idChar(i),
2371  coupSUSYPtr->idChar(j) ) ) continue;
2372  sigmaPtr = new Sigma2qqbar2charchar( i,-j, iproc);
2373  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2374  }
2375  }
2376  }
2377 
2378  // RPV squark production
2379  if(SUSYs || settings.flag("SUSY:qq2antisquark")) {
2380  for (int idx = 1; idx <= 6; ++idx) {
2381  for (int iso = 1; iso <= 2; ++iso) {
2382  int id1 = iso + ((idx <= 3) ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2383  // Skip if outgoing code not asked for
2384  if (!allowIdVals( id1, 0)) continue;
2385  sigmaPtr = new Sigma1qq2antisquark(id1);
2386  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2387  }
2388  }
2389  }
2390 
2391  // Neutralino-gluino
2392  if (SUSYs || settings.flag("SUSY:qqbar2chi0gluino")) {
2393  int iproc = 1600;
2394  for (int iNeut = 1; iNeut <= nNeut; iNeut++) {
2395  iproc++;
2396  // Skip if outgoing codes not asked for
2397  if (!allowIdVals( coupSUSYPtr->idNeut(iNeut), 1000021)) continue;
2398  sigmaPtr = new Sigma2qqbar2chi0gluino(iNeut, iproc);
2399  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2400  }
2401  }
2402 
2403  // Chargino-Gluino
2404  if (SUSYs || settings.flag("SUSY:qqbar2chi+-gluino")) {
2405  int iproc = 1620;
2406  for (int iChar = 1; iChar <= 2; ++iChar) {
2407  iproc += 2;
2408  // Skip if outgoing codes not asked for
2409  if (!allowIdVals( coupSUSYPtr->idChar(iChar), 1000021)) continue;
2410  sigmaPtr = new Sigma2qqbar2chargluino( iChar, iproc-1);
2411  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2412  sigmaPtr = new Sigma2qqbar2chargluino( -iChar, iproc);
2413  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2414 
2415  }
2416  }
2417 
2418  // Slepton-antislepton (qqbar initiated); Currently no RH sneutrinos
2419  if (SUSYs || settings.flag("SUSY:qqbar2sleptonantislepton")) {
2420  int iproc = 1650;
2421  for (int idx = 1; idx <= 6; ++idx) {
2422  for (int iso = 1; iso <= 2; ++iso) {
2423  for (int jso = iso; jso >= 1; --jso) {
2424  for (int jdx = 1; jdx <= 6; ++jdx) {
2425  if (iso == jso && jdx < idx) continue;
2426  int id1 = iso + ((idx <= 3) ? 1000010+2*(idx-1)
2427  : 2000010+2*(idx-4));
2428  int id2 = jso + ((jdx <= 3) ? 1000010+2*(jdx-1)
2429  : 2000010+2*(jdx-4));
2430  // Update process number counter
2431  iproc++;
2432  if (iso == jso && id1 != id2) iproc++;
2433  // Exclude RH neutrinos from allowed final states
2434  if (abs(id1) >= 2000012 && id1 % 2 == 0) continue;
2435  if (abs(id2) >= 2000012 && id2 % 2 == 0) continue;
2436  // Skip if outgoing codes not asked for
2437  if (!allowIdVals( id1, id2)) continue;
2438  if (iso == jso && id1 != id2) {
2439  sigmaPtr
2440  = new Sigma2qqbar2sleptonantislepton(id1,-id2,iproc-1);
2441  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2442  sigmaPtr = new Sigma2qqbar2sleptonantislepton(id2,-id1,iproc);
2443  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2444  } else {
2445  sigmaPtr = new Sigma2qqbar2sleptonantislepton(id1,-id2,iproc);
2446  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2447  }
2448  }
2449  }
2450  }
2451  }
2452  }
2453 
2454  } // End of SUSY processes.
2455 
2456  // Set up requested objects for New-Gauge-Boson processes.
2457  if (settings.flag("NewGaugeBoson:ffbar2gmZZprime")) {
2458  sigmaPtr = new Sigma1ffbar2gmZZprime();
2459  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2460  }
2461  if (settings.flag("NewGaugeBoson:ffbar2Wprime")) {
2462  sigmaPtr = new Sigma1ffbar2Wprime();
2463  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2464  }
2465  if (settings.flag("NewGaugeBoson:ffbar2R0")) {
2466  sigmaPtr = new Sigma1ffbar2Rhorizontal();
2467  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2468  }
2469 
2470  // Set up requested objects for Left-Right-Symmetry processes.
2471  bool leftrights = settings.flag("LeftRightSymmmetry:all");
2472  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2ZR")) {
2473  sigmaPtr = new Sigma1ffbar2ZRight();
2474  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2475  }
2476  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2WR")) {
2477  sigmaPtr = new Sigma1ffbar2WRight();
2478  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2479  }
2480  if (leftrights || settings.flag("LeftRightSymmmetry:ll2HL")) {
2481  sigmaPtr = new Sigma1ll2Hchgchg(1);
2482  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2483  }
2484  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLe")) {
2485  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 11);
2486  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2487  }
2488  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLmu")) {
2489  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 13);
2490  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2491  }
2492  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HLtau")) {
2493  sigmaPtr = new Sigma2lgm2Hchgchgl(1, 15);
2494  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2495  }
2496  if (leftrights || settings.flag("LeftRightSymmmetry:ff2HLff")) {
2497  sigmaPtr = new Sigma3ff2HchgchgfftWW(1);
2498  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2499  }
2500  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2HLHL")) {
2501  sigmaPtr = new Sigma2ffbar2HchgchgHchgchg(1);
2502  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2503  }
2504  if (leftrights || settings.flag("LeftRightSymmmetry:ll2HR")) {
2505  sigmaPtr = new Sigma1ll2Hchgchg(2);
2506  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2507  }
2508  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRe")) {
2509  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 11);
2510  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2511  }
2512  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRmu")) {
2513  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 13);
2514  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2515  }
2516  if (leftrights || settings.flag("LeftRightSymmmetry:lgm2HRtau")) {
2517  sigmaPtr = new Sigma2lgm2Hchgchgl(2, 15);
2518  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2519  }
2520  if (leftrights || settings.flag("LeftRightSymmmetry:ff2HRff")) {
2521  sigmaPtr = new Sigma3ff2HchgchgfftWW(2);
2522  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2523  }
2524  if (leftrights || settings.flag("LeftRightSymmmetry:ffbar2HRHR")) {
2525  sigmaPtr = new Sigma2ffbar2HchgchgHchgchg(2);
2526  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2527  }
2528 
2529  // Set up requested objects for leptoquark LQ processes.
2530  bool leptoquarks = settings.flag("LeptoQuark:all");
2531  if (leptoquarks || settings.flag("LeptoQuark:ql2LQ")) {
2532  sigmaPtr = new Sigma1ql2LeptoQuark;
2533  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2534  }
2535  if (leptoquarks || settings.flag("LeptoQuark:qg2LQl")) {
2536  sigmaPtr = new Sigma2qg2LeptoQuarkl;
2537  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2538  }
2539  if (leptoquarks || settings.flag("LeptoQuark:gg2LQLQbar")) {
2540  sigmaPtr = new Sigma2gg2LQLQbar;
2541  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2542  }
2543  if (leptoquarks || settings.flag("LeptoQuark:qqbar2LQLQbar")) {
2544  sigmaPtr = new Sigma2qqbar2LQLQbar;
2545  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2546  }
2547 
2548  // Set up requested objects for excited-fermion processes.
2549  bool excitedfermions = settings.flag("ExcitedFermion:all");
2550  if (excitedfermions || settings.flag("ExcitedFermion:dg2dStar")) {
2551  sigmaPtr = new Sigma1qg2qStar(1);
2552  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2553  }
2554  if (excitedfermions || settings.flag("ExcitedFermion:ug2uStar")) {
2555  sigmaPtr = new Sigma1qg2qStar(2);
2556  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2557  }
2558  if (excitedfermions || settings.flag("ExcitedFermion:sg2sStar")) {
2559  sigmaPtr = new Sigma1qg2qStar(3);
2560  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2561  }
2562  if (excitedfermions || settings.flag("ExcitedFermion:cg2cStar")) {
2563  sigmaPtr = new Sigma1qg2qStar(4);
2564  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2565  }
2566  if (excitedfermions || settings.flag("ExcitedFermion:bg2bStar")) {
2567  sigmaPtr = new Sigma1qg2qStar(5);
2568  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2569  }
2570  if (excitedfermions || settings.flag("ExcitedFermion:egm2eStar")) {
2571  sigmaPtr = new Sigma1lgm2lStar(11);
2572  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2573  }
2574  if (excitedfermions || settings.flag("ExcitedFermion:mugm2muStar")) {
2575  sigmaPtr = new Sigma1lgm2lStar(13);
2576  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2577  }
2578  if (excitedfermions || settings.flag("ExcitedFermion:taugm2tauStar")) {
2579  sigmaPtr = new Sigma1lgm2lStar(15);
2580  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2581  }
2582  if (excitedfermions || settings.flag("ExcitedFermion:qq2dStarq")) {
2583  sigmaPtr = new Sigma2qq2qStarq(1);
2584  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2585  }
2586  if (excitedfermions || settings.flag("ExcitedFermion:qq2uStarq")) {
2587  sigmaPtr = new Sigma2qq2qStarq(2);
2588  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2589  }
2590  if (excitedfermions || settings.flag("ExcitedFermion:qq2sStarq")) {
2591  sigmaPtr = new Sigma2qq2qStarq(3);
2592  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2593  }
2594  if (excitedfermions || settings.flag("ExcitedFermion:qq2cStarq")) {
2595  sigmaPtr = new Sigma2qq2qStarq(4);
2596  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2597  }
2598  if (excitedfermions || settings.flag("ExcitedFermion:qq2bStarq")) {
2599  sigmaPtr = new Sigma2qq2qStarq(5);
2600  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2601  }
2602  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2eStare")) {
2603  sigmaPtr = new Sigma2qqbar2lStarlbar(11);
2604  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2605  }
2606  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2nueStarnue")) {
2607  sigmaPtr = new Sigma2qqbar2lStarlbar(12);
2608  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2609  }
2610  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2muStarmu")) {
2611  sigmaPtr = new Sigma2qqbar2lStarlbar(13);
2612  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2613  }
2614  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2numuStarnumu")) {
2615  sigmaPtr = new Sigma2qqbar2lStarlbar(14);
2616  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2617  }
2618  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2tauStartau")) {
2619  sigmaPtr = new Sigma2qqbar2lStarlbar(15);
2620  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2621  }
2622  if (excitedfermions
2623  || settings.flag("ExcitedFermion:qqbar2nutauStarnutau")) {
2624  sigmaPtr = new Sigma2qqbar2lStarlbar(16);
2625  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2626  }
2627  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2eStareStar")) {
2628  sigmaPtr = new Sigma2qqbar2lStarlStarBar(11);
2629  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2630  }
2631  if (excitedfermions
2632  || settings.flag("ExcitedFermion:qqbar2nueStarnueStar")) {
2633  sigmaPtr = new Sigma2qqbar2lStarlStarBar(12);
2634  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2635  }
2636  if (excitedfermions || settings.flag("ExcitedFermion:qqbar2muStarmuStar")) {
2637  sigmaPtr = new Sigma2qqbar2lStarlStarBar(13);
2638  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2639  }
2640  if (excitedfermions
2641  || settings.flag("ExcitedFermion:qqbar2numuStarnumuStar")) {
2642  sigmaPtr = new Sigma2qqbar2lStarlStarBar(14);
2643  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2644  }
2645  if (excitedfermions
2646  || settings.flag("ExcitedFermion:qqbar2tauStartauStar")) {
2647  sigmaPtr = new Sigma2qqbar2lStarlStarBar(15);
2648  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2649  }
2650  if (excitedfermions
2651  || settings.flag("ExcitedFermion:qqbar2nutauStarnutauStar")) {
2652  sigmaPtr = new Sigma2qqbar2lStarlStarBar(16);
2653  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2654  }
2655 
2656  // Set up requested objects for contact interaction processes.
2657  if (settings.flag("ContactInteractions:QCqq2qq")) {
2658  sigmaPtr = new Sigma2QCqq2qq();
2659  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2660  }
2661  if (settings.flag("ContactInteractions:QCqqbar2qqbar")) {
2662  sigmaPtr = new Sigma2QCqqbar2qqbar();
2663  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2664  }
2665  if (settings.flag("ContactInteractions:QCffbar2eebar")) {
2666  sigmaPtr = new Sigma2QCffbar2llbar(11, 4203);
2667  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2668  }
2669  if (settings.flag("ContactInteractions:QCffbar2mumubar")) {
2670  sigmaPtr = new Sigma2QCffbar2llbar(13, 4204);
2671  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2672  }
2673  if (settings.flag("ContactInteractions:QCffbar2tautaubar")) {
2674  sigmaPtr = new Sigma2QCffbar2llbar(15, 4205);
2675  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2676  }
2677 
2678  // Set spin of particles in the Hidden Valley scenario.
2679  int spinFv = settings.mode("HiddenValley:spinFv");
2680  for (int i = 1; i < 7; ++i) {
2681  if (particleDataPtr->spinType( 4900000 + i) != spinFv + 1)
2682  particleDataPtr->spinType( 4900000 + i, spinFv + 1);
2683  if (particleDataPtr->spinType( 4900010 + i) != spinFv + 1)
2684  particleDataPtr->spinType( 4900010 + i, spinFv + 1);
2685  }
2686  if (spinFv != 1) {
2687  if (particleDataPtr->spinType( 4900101) != 2)
2688  particleDataPtr->spinType( 4900101, 2);
2689  } else {
2690  int spinqv = settings.mode("HiddenValley:spinqv");
2691  if (particleDataPtr->spinType( 4900101) != 2 * spinqv + 1)
2692  particleDataPtr->spinType( 4900101, 2 * spinqv + 1);
2693  }
2694 
2695  // Set up requested objects for HiddenValley processes.
2696  bool hiddenvalleys = settings.flag("HiddenValley:all");
2697  if (hiddenvalleys || settings.flag("HiddenValley:gg2DvDvbar")) {
2698  sigmaPtr = new Sigma2gg2qGqGbar( 4900001, 4901, spinFv,
2699  "g g -> Dv Dvbar");
2700  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2701  }
2702  if (hiddenvalleys || settings.flag("HiddenValley:gg2UvUvbar")) {
2703  sigmaPtr = new Sigma2gg2qGqGbar( 4900002, 4902, spinFv,
2704  "g g -> Uv Uvbar");
2705  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2706  }
2707  if (hiddenvalleys || settings.flag("HiddenValley:gg2SvSvbar")) {
2708  sigmaPtr = new Sigma2gg2qGqGbar( 4900003, 4903, spinFv,
2709  "g g -> Sv Svbar");
2710  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2711  }
2712  if (hiddenvalleys || settings.flag("HiddenValley:gg2CvCvbar")) {
2713  sigmaPtr = new Sigma2gg2qGqGbar( 4900004, 4904, spinFv,
2714  "g g -> Cv Cvbar");
2715  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2716  }
2717  if (hiddenvalleys || settings.flag("HiddenValley:gg2BvBvbar")) {
2718  sigmaPtr = new Sigma2gg2qGqGbar( 4900005, 4905, spinFv,
2719  "g g -> Bv Bvbar");
2720  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2721  }
2722  if (hiddenvalleys || settings.flag("HiddenValley:gg2TvTvbar")) {
2723  sigmaPtr = new Sigma2gg2qGqGbar( 4900006, 4906, spinFv,
2724  "g g -> Tv Tvbar");
2725  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2726  }
2727  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2DvDvbar")) {
2728  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900001, 4911, spinFv,
2729  "q qbar -> Dv Dvbar");
2730  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2731  }
2732  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2UvUvbar")) {
2733  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900002, 4912, spinFv,
2734  "q qbar -> Uv Uvbar");
2735  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2736  }
2737  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2SvSvbar")) {
2738  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900003, 4913, spinFv,
2739  "q qbar -> Sv Svbar");
2740  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2741  }
2742  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2CvCvbar")) {
2743  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900004, 4914, spinFv,
2744  "q qbar -> Cv Cvbar");
2745  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2746  }
2747  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2BvBvbar")) {
2748  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900005, 4915, spinFv,
2749  "q qbar -> Bv Bvbar");
2750  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2751  }
2752  if (hiddenvalleys || settings.flag("HiddenValley:qqbar2TvTvbar")) {
2753  sigmaPtr = new Sigma2qqbar2qGqGbar( 4900006, 4916, spinFv,
2754  "q qbar -> Tv Tvbar");
2755  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2756  }
2757  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2DvDvbar")) {
2758  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900001, 4921, spinFv,
2759  "f fbar -> Dv Dvbar");
2760  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2761  }
2762  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2UvUvbar")) {
2763  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900002, 4922, spinFv,
2764  "f fbar -> Uv Uvbar");
2765  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2766  }
2767  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2SvSvbar")) {
2768  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900003, 4923, spinFv,
2769  "f fbar -> Sv Svbar");
2770  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2771  }
2772  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2CvCvbar")) {
2773  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900004, 4924, spinFv,
2774  "f fbar -> Cv Cvbar");
2775  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2776  }
2777  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2BvBvbar")) {
2778  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900005, 4925, spinFv,
2779  "f fbar -> Bv Bvbar");
2780  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2781  }
2782  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2TvTvbar")) {
2783  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900006, 4926, spinFv,
2784  "f fbar -> Tv Tvbar");
2785  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2786  }
2787  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2EvEvbar")) {
2788  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900011, 4931, spinFv,
2789  "f fbar -> Ev Evbar");
2790  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2791  }
2792  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuEvnuEvbar")) {
2793  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900012, 4932, spinFv,
2794  "f fbar -> nuEv nuEvbar");
2795  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2796  }
2797  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2MUvMUvbar")) {
2798  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900013, 4933, spinFv,
2799  "f fbar -> MUv MUvbar");
2800  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2801  }
2802  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuMUvnuMUvbar")) {
2803  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900014, 4934, spinFv,
2804  "f fbar -> nuMUv nuMUvbar");
2805  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2806  }
2807  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2TAUvTAUvbar")) {
2808  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900015, 4935, spinFv,
2809  "f fbar -> TAUv TAUvbar");
2810  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2811  }
2812  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2nuTAUvnuTAUvbar")) {
2813  sigmaPtr = new Sigma2ffbar2fGfGbar( 4900016, 4936, spinFv,
2814  "f fbar -> nuTAUv nuTAUvbar");
2815  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2816  }
2817  if (hiddenvalleys || settings.flag("HiddenValley:ffbar2Zv")) {
2818  sigmaPtr = new Sigma1ffbar2Zv();
2819  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2820  }
2821 
2822  // Set up requested objects for RS extra-dimensional G* processes.
2823  bool extraDimGstars = settings.flag("ExtraDimensionsG*:all");
2824  if (extraDimGstars || settings.flag("ExtraDimensionsG*:gg2G*")) {
2825  sigmaPtr = new Sigma1gg2GravitonStar;
2826  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2827  }
2828  if (extraDimGstars || settings.flag("ExtraDimensionsG*:ffbar2G*")) {
2829  sigmaPtr = new Sigma1ffbar2GravitonStar;
2830  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2831  }
2832  if (settings.flag("ExtraDimensionsG*:gg2G*g")) {
2833  sigmaPtr = new Sigma2gg2GravitonStarg;
2834  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2835  }
2836  if (settings.flag("ExtraDimensionsG*:qg2G*q")) {
2837  sigmaPtr = new Sigma2qg2GravitonStarq;
2838  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2839  }
2840  if (settings.flag("ExtraDimensionsG*:qqbar2G*g")) {
2841  sigmaPtr = new Sigma2qqbar2GravitonStarg;
2842  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2843  }
2844 
2845  // Set up requested objects for RS extra-dimensional KKgluon processes.
2846  if (settings.flag("ExtraDimensionsG*:qqbar2KKgluon*")) {
2847  sigmaPtr = new Sigma1qqbar2KKgluonStar;
2848  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2849  }
2850 
2851  // NOAM: Set up requested objects for TEV extra-dimensional processes.
2852  if (settings.flag("ExtraDimensionsTEV:ffbar2ddbar")) {
2853  sigmaPtr = new Sigma2ffbar2TEVffbar(1, 5061);
2854  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2855  }
2856  if (settings.flag("ExtraDimensionsTEV:ffbar2uubar")) {
2857  sigmaPtr = new Sigma2ffbar2TEVffbar(2, 5062);
2858  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2859  }
2860  if (settings.flag("ExtraDimensionsTEV:ffbar2ssbar")) {
2861  sigmaPtr = new Sigma2ffbar2TEVffbar(3, 5063);
2862  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2863  }
2864  if (settings.flag("ExtraDimensionsTEV:ffbar2ccbar")) {
2865  sigmaPtr = new Sigma2ffbar2TEVffbar(4, 5064);
2866  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2867  }
2868  if (settings.flag("ExtraDimensionsTEV:ffbar2bbbar")) {
2869  sigmaPtr = new Sigma2ffbar2TEVffbar(5, 5065);
2870  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2871  }
2872  if (settings.flag("ExtraDimensionsTEV:ffbar2ttbar")) {
2873  sigmaPtr = new Sigma2ffbar2TEVffbar(6, 5066);
2874  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2875  }
2876  if (settings.flag("ExtraDimensionsTEV:ffbar2e+e-")) {
2877  sigmaPtr = new Sigma2ffbar2TEVffbar(11, 5071);
2878  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2879  }
2880  if (settings.flag("ExtraDimensionsTEV:ffbar2nuenuebar")) {
2881  sigmaPtr = new Sigma2ffbar2TEVffbar(12, 5072);
2882  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2883  }
2884  if (settings.flag("ExtraDimensionsTEV:ffbar2mu+mu-")) {
2885  sigmaPtr = new Sigma2ffbar2TEVffbar(13, 5073);
2886  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2887  }
2888  if (settings.flag("ExtraDimensionsTEV:ffbar2numunumubar")) {
2889  sigmaPtr = new Sigma2ffbar2TEVffbar(14, 5074);
2890  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2891  }
2892  if (settings.flag("ExtraDimensionsTEV:ffbar2tau+tau-")) {
2893  sigmaPtr = new Sigma2ffbar2TEVffbar(15, 5075);
2894  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2895  }
2896  if (settings.flag("ExtraDimensionsTEV:ffbar2nutaunutaubar")) {
2897  sigmaPtr = new Sigma2ffbar2TEVffbar(16, 5076);
2898  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2899  }
2900 
2901  // Set up requested objects for large extra-dimensional G processes.
2902  bool extraDimLEDmono = settings.flag("ExtraDimensionsLED:monojet");
2903  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:gg2Gg")) {
2904  sigmaPtr = new Sigma2gg2LEDUnparticleg( true );
2905  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2906  }
2907  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:qg2Gq")) {
2908  sigmaPtr = new Sigma2qg2LEDUnparticleq( true );
2909  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2910  }
2911  if (extraDimLEDmono || settings.flag("ExtraDimensionsLED:qqbar2Gg")) {
2912  sigmaPtr = new Sigma2qqbar2LEDUnparticleg( true );
2913  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2914  }
2915  if (settings.flag("ExtraDimensionsLED:ffbar2GZ")) {
2916  sigmaPtr = new Sigma2ffbar2LEDUnparticleZ( true );
2917  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2918  }
2919  if (settings.flag("ExtraDimensionsLED:ffbar2Ggamma")) {
2920  sigmaPtr = new Sigma2ffbar2LEDUnparticlegamma( true );
2921  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2922  }
2923  if (settings.flag("ExtraDimensionsLED:ffbar2gammagamma")) {
2924  sigmaPtr = new Sigma2ffbar2LEDgammagamma( true );
2925  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2926  }
2927  if (settings.flag("ExtraDimensionsLED:gg2gammagamma")) {
2928  sigmaPtr = new Sigma2gg2LEDgammagamma( true );
2929  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2930  }
2931  if (settings.flag("ExtraDimensionsLED:ffbar2llbar")) {
2932  sigmaPtr = new Sigma2ffbar2LEDllbar( true );
2933  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2934  }
2935  if (settings.flag("ExtraDimensionsLED:gg2llbar")) {
2936  sigmaPtr = new Sigma2gg2LEDllbar( true );
2937  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2938  }
2939  bool extraDimLEDdij = settings.flag("ExtraDimensionsLED:dijets");
2940  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:gg2DJgg")) {
2941  sigmaPtr = new Sigma2gg2LEDgg;
2942  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2943  }
2944  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:gg2DJqqbar")) {
2945  sigmaPtr = new Sigma2gg2LEDqqbar;
2946  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2947  }
2948  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qg2DJqg")) {
2949  sigmaPtr = new Sigma2qg2LEDqg;
2950  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2951  }
2952  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qq2DJqq")) {
2953  sigmaPtr = new Sigma2qq2LEDqq;
2954  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2955  }
2956  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qqbar2DJgg")) {
2957  sigmaPtr = new Sigma2qqbar2LEDgg;
2958  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2959  }
2960  if (extraDimLEDdij || settings.flag("ExtraDimensionsLED:qqbar2DJqqbarNew")) {
2961  sigmaPtr = new Sigma2qqbar2LEDqqbarNew;
2962  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2963  }
2964 
2965  // Set up requested objects for unparticle processes.
2966  bool extraDimUnpartmono = settings.flag("ExtraDimensionsUnpart:monojet");
2967  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:gg2Ug")) {
2968  sigmaPtr = new Sigma2gg2LEDUnparticleg( false );
2969  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2970  }
2971  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:qg2Uq")) {
2972  sigmaPtr = new Sigma2qg2LEDUnparticleq( false );
2973  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2974  }
2975  if (extraDimUnpartmono || settings.flag("ExtraDimensionsUnpart:qqbar2Ug")) {
2976  sigmaPtr = new Sigma2qqbar2LEDUnparticleg( false );
2977  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2978  }
2979  if (settings.flag("ExtraDimensionsUnpart:ffbar2UZ")) {
2980  sigmaPtr = new Sigma2ffbar2LEDUnparticleZ( false );
2981  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2982  }
2983  if (settings.flag("ExtraDimensionsUnpart:ffbar2Ugamma")) {
2984  sigmaPtr = new Sigma2ffbar2LEDUnparticlegamma( false );
2985  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2986  }
2987  if (settings.flag("ExtraDimensionsUnpart:ffbar2gammagamma")) {
2988  sigmaPtr = new Sigma2ffbar2LEDgammagamma( false );
2989  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2990  }
2991  if (settings.flag("ExtraDimensionsUnpart:gg2gammagamma")) {
2992  sigmaPtr = new Sigma2gg2LEDgammagamma( false );
2993  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2994  }
2995  if (settings.flag("ExtraDimensionsUnpart:ffbar2llbar")) {
2996  sigmaPtr = new Sigma2ffbar2LEDllbar( false );
2997  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
2998  }
2999  if (settings.flag("ExtraDimensionsUnpart:gg2llbar")) {
3000  sigmaPtr = new Sigma2gg2LEDllbar( false );
3001  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3002  }
3003 
3004  // Set up requested objects for Dark Matter processes.
3005  if (settings.flag("DM:ffbar2Zp2XX")) {
3006  sigmaPtr = new Sigma1ffbar2Zp2XX();
3007  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3008  }
3009  if (settings.flag("DM:ffbar2Zp2XXj")) {
3010  sigmaPtr = new Sigma2qqbar2Zpg2XXj();
3011  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3012  }
3013  if (settings.flag("DM:ffbar2ZpH")) {
3014  sigmaPtr = new Sigma2ffbar2ZpH();
3015  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3016  }
3017  if (settings.flag("DM:gg2S2XX")) {
3018  sigmaPtr = new Sigma1gg2S2XX();
3019  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3020  }
3021  if (settings.flag("DM:gg2S2XXj")) {
3022  sigmaPtr = new Sigma2gg2Sg2XXj();
3023  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3024  }
3025  if (settings.flag("DM:qqbar2DY")) {
3026  sigmaPtr = new Sigma2qqbar2DY();
3027  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
3028  }
3029 
3030  for ( ProcessContainer * cont : containerPtrs )
3031  cont->initInfoPtr(*infoPtr);
3032 
3033  // Done.
3034  return true;
3035 
3036 }
3037 
3038 //--------------------------------------------------------------------------
3039 
3040 // Routine to initialize list of second hard processes.
3041 
3042 bool SetupContainers::init2(vector<ProcessContainer*>& container2Ptrs,
3043  Info* infoPtr) {
3044 
3045  Settings& settings = *infoPtr->settingsPtr;
3046 
3047  // Reset process list, if filled in previous subrun.
3048  if (container2Ptrs.size() > 0) {
3049  for (int i = 0; i < int(container2Ptrs.size()); ++i)
3050  delete container2Ptrs[i];
3051  container2Ptrs.clear();
3052  }
3053  SigmaProcess* sigmaPtr;
3054 
3055  // Two hard QCD jets.
3056  if (settings.flag("SecondHard:TwoJets")) {
3057  sigmaPtr = new Sigma2gg2gg;
3058  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3059  sigmaPtr = new Sigma2gg2qqbar;
3060  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3061  sigmaPtr = new Sigma2qg2qg;
3062  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3063  sigmaPtr = new Sigma2qq2qq;
3064  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3065  sigmaPtr = new Sigma2qqbar2gg;
3066  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3067  sigmaPtr = new Sigma2qqbar2qqbarNew;
3068  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3069  sigmaPtr = new Sigma2gg2QQbar(4, 121);
3070  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3071  sigmaPtr = new Sigma2qqbar2QQbar(4, 122);
3072  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3073  sigmaPtr = new Sigma2gg2QQbar(5, 123);
3074  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3075  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
3076  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3077  }
3078 
3079  // A prompt photon and a hard jet.
3080  if (settings.flag("SecondHard:PhotonAndJet")) {
3081  sigmaPtr = new Sigma2qg2qgamma;
3082  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3083  sigmaPtr = new Sigma2qqbar2ggamma;
3084  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3085  sigmaPtr = new Sigma2gg2ggamma;
3086  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3087  }
3088 
3089  // Two prompt photons.
3090  if (settings.flag("SecondHard:TwoPhotons")) {
3091  sigmaPtr = new Sigma2ffbar2gammagamma;
3092  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3093  sigmaPtr = new Sigma2gg2gammagamma;
3094  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3095  }
3096 
3097  // Charmonium.
3098  if (settings.flag("SecondHard:Charmonium")) {
3099  vector<SigmaProcess*> charmoniumSigmaPtrs;
3100  charmonium.setupSigma2gg(charmoniumSigmaPtrs, true);
3101  charmonium.setupSigma2qg(charmoniumSigmaPtrs, true);
3102  charmonium.setupSigma2qq(charmoniumSigmaPtrs, true);
3103  for (unsigned int i = 0; i < charmoniumSigmaPtrs.size(); ++i)
3104  container2Ptrs.push_back( new ProcessContainer(charmoniumSigmaPtrs[i]));
3105  }
3106 
3107  // Bottomonium.
3108  if (settings.flag("SecondHard:Bottomonium")) {
3109  vector<SigmaProcess*> bottomoniumSigmaPtrs;
3110  bottomonium.setupSigma2gg(bottomoniumSigmaPtrs, true);
3111  bottomonium.setupSigma2qg(bottomoniumSigmaPtrs, true);
3112  bottomonium.setupSigma2qq(bottomoniumSigmaPtrs, true);
3113  for (unsigned int i = 0; i < bottomoniumSigmaPtrs.size(); ++i)
3114  container2Ptrs.push_back( new ProcessContainer(bottomoniumSigmaPtrs[i]));
3115  }
3116 
3117  // A single gamma*/Z0.
3118  if (settings.flag("SecondHard:SingleGmZ")) {
3119  sigmaPtr = new Sigma1ffbar2gmZ;
3120  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3121  }
3122 
3123  // A single W+-.
3124  if (settings.flag("SecondHard:SingleW")) {
3125  sigmaPtr = new Sigma1ffbar2W;
3126  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3127  }
3128 
3129  // A gamma*/Z0 and a hard jet.
3130  if (settings.flag("SecondHard:GmZAndJet")) {
3131  sigmaPtr = new Sigma2qqbar2gmZg;
3132  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3133  sigmaPtr = new Sigma2qg2gmZq;
3134  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3135  }
3136 
3137  // A W+- and a hard jet.
3138  if (settings.flag("SecondHard:WAndJet")) {
3139  sigmaPtr = new Sigma2qqbar2Wg;
3140  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3141  sigmaPtr = new Sigma2qg2Wq;
3142  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3143  }
3144 
3145  // Top pair production.
3146  if (settings.flag("SecondHard:TopPair")) {
3147  sigmaPtr = new Sigma2gg2QQbar(6, 601);
3148  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3149  sigmaPtr = new Sigma2qqbar2QQbar(6, 602);
3150  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3151  sigmaPtr = new Sigma2ffbar2FFbarsgmZ(6, 604);
3152  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3153  }
3154 
3155  // Single top production.
3156  if (settings.flag("SecondHard:SingleTop")) {
3157  sigmaPtr = new Sigma2qq2QqtW(6, 603);
3158  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3159  sigmaPtr = new Sigma2ffbar2FfbarsW(6, 0, 605);
3160  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3161  }
3162 
3163  // Two b jets - already part of TwoJets sample above.
3164  if (settings.flag("SecondHard:TwoBJets")) {
3165  sigmaPtr = new Sigma2gg2QQbar(5, 123);
3166  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3167  sigmaPtr = new Sigma2qqbar2QQbar(5, 124);
3168  container2Ptrs.push_back( new ProcessContainer(sigmaPtr) );
3169  }
3170 
3171  for ( ProcessContainer * cont : container2Ptrs )
3172  cont->initInfoPtr(*infoPtr);
3173 
3174  // Done.
3175  return true;
3176 
3177 }
3178 
3179 //--------------------------------------------------------------------------
3180 
3181 // Set up arrays of allowed outgoing SUSY particles.
3182 
3183 void SetupContainers::setupIdVecs( Settings& settings) {
3184 
3185  // First array either none, one or many particles.
3186  idVecA.resize(0);
3187  if (settings.mode("SUSY:idA") != 0) {
3188  idVecA.push_back( abs(settings.mode("SUSY:idA")) );
3189  } else {
3190  vector<int> idTmpA = settings.mvec("SUSY:idVecA");
3191  for (int i = 0; i < int(idTmpA.size()); ++i)
3192  if (idTmpA[i] != 0) idVecA.push_back( abs(idTmpA[i]) );
3193  }
3194  nVecA = idVecA.size();
3195 
3196  // Second array either none, one or many particles.
3197  idVecB.resize(0);
3198  if (settings.mode("SUSY:idB") != 0) {
3199  idVecB.push_back( abs(settings.mode("SUSY:idB")) );
3200  } else {
3201  vector<int> idTmpB = settings.mvec("SUSY:idVecB");
3202  for (int i = 0; i < int(idTmpB.size()); ++i)
3203  if (idTmpB[i] != 0) idVecB.push_back( abs(idTmpB[i]) );
3204  }
3205  nVecB = idVecB.size();
3206 
3207 }
3208 
3209 //--------------------------------------------------------------------------
3210 
3211 // Check final state for allowed outgoing SUSY particles.
3212 // Normally check two codes, but allow for only one.
3213 
3214 bool SetupContainers::allowIdVals( int idCheck1, int idCheck2) {
3215 
3216  // If empty arrays or id's no need for checks. Else need absolute values.
3217  if (nVecA == 0 && nVecB == 0) return true;
3218  if (idCheck1 == 0 && idCheck2 == 0) return true;
3219  int idChk1 = abs(idCheck1);
3220  int idChk2 = abs(idCheck2);
3221 
3222  // If only one outgoing particle then check idVecA and idVecB.
3223  if (idChk1 == 0) swap(idChk1, idChk2);
3224  if (idChk2 == 0) {
3225  for (int i = 0; i < nVecA; ++i) if (idChk1 == idVecA[i]) return true;
3226  for (int i = 0; i < nVecB; ++i) if (idChk1 == idVecB[i]) return true;
3227  return false;
3228  }
3229 
3230  // If empty array idVecB then compare with idVecA.
3231  if (nVecB == 0) {
3232  for (int i = 0; i < nVecA; ++i)
3233  if (idChk1 == idVecA[i] || idChk2 == idVecA[i]) return true;
3234  return false;
3235  }
3236 
3237  // If empty array idVecA then compare with idVecB.
3238  if (nVecA == 0) {
3239  for (int i = 0; i < nVecB; ++i)
3240  if (idChk1 == idVecB[i] || idChk2 == idVecB[i]) return true;
3241  return false;
3242  }
3243 
3244  // Else check that pair matches allowed combinations.
3245  for (int i = 0; i < nVecA; ++i)
3246  for (int j = 0; j < nVecB; ++j)
3247  if ( (idChk1 == idVecA[i] && idChk2 == idVecB[j])
3248  || (idChk2 == idVecA[i] && idChk1 == idVecB[j]) ) return true;
3249  return false;
3250 
3251 }
3252 
3253 //==========================================================================
3254 
3255 } // end namespace Pythia8
Definition: AgUStep.h:26