9 #include "Pythia8/ProcessContainer.h"
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"
38 const int ProcessContainer::N12SAMPLE = 100;
41 const int ProcessContainer::N3SAMPLE = 1000;
48 bool ProcessContainer::init(
bool isFirst, ResonanceDecays* resDecaysPtrIn,
49 SLHAinterface* slhaInterfacePtr, GammaKinematics* gammaKinPtrIn) {
51 registerSubObject(*sigmaProcessPtr);
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();
69 useStrictLHEFscales = flag(
"Beams:strictLHEFscale");
73 nTryRequested = mode(
"Main:numberOfTriedEvents");
74 nSelRequested = mode(
"Main:numberOfSelectedEvents");
75 nAccRequested = mode(
"Main:numberOfAcceptedEvents");
78 increaseMaximum = flag(
"PhaseSpace:increaseMaximum");
81 gammaKinPtr = gammaKinPtrIn;
84 approximatedGammaFlux = beamAPtr->hasApproxGammaFlux() ||
85 beamBPtr->hasApproxGammaFlux();
88 bool beamAhasGamma = flag(
"PDF:beamA2gamma");
89 bool beamBhasGamma = flag(
"PDF:beamB2gamma");
90 beamHasGamma = beamAhasGamma || beamBhasGamma;
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(
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();
108 resDecaysPtr = resDecaysPtrIn;
109 canVetoResDecay = (userHooksPtr != 0)
110 ? userHooksPtr->canVetoResonanceDecays() :
false;
112 sigmaProcessPtr->setLHAPtr(lhaUpPtr);
113 phaseSpacePtr->setLHAPtr(lhaUpPtr);
116 double relEnergyDiff = ( lhaUpPtr->eBeamA() - lhaUpPtr->eBeamB() )
117 / ( lhaUpPtr->eBeamA() + lhaUpPtr->eBeamB() );
118 isAsymLHA = ( abs(relEnergyDiff) > 1.e-10 )
119 || ( lhaUpPtr->idBeamA() != lhaUpPtr->idBeamB() );
121 betazLHA = (sqrtpos(pow2(lhaUpPtr->eBeamA()) - pow2(infoPtr->mA()))
122 - sqrtpos(pow2(lhaUpPtr->eBeamB()) - pow2(infoPtr->mB())))
123 / (lhaUpPtr->eBeamA() + lhaUpPtr->eBeamB());
126 sigmaProcessPtr->init(beamAPtr, beamBPtr, slhaInterfacePtr);
130 string inState = sigmaProcessPtr->inFlux();
134 if ( beamAPtr->isGamma() || beamAhasGamma ) {
135 if ( inState ==
"gmg" || inState ==
"gmq" || inState ==
"gmgm" )
137 else beamAgammaMode = 1;
139 if ( beamBPtr->isGamma() || beamBhasGamma ) {
140 if ( inState ==
"ggm" || inState ==
"qgm" || inState ==
"gmgm" )
142 else beamBgammaMode = 1;
147 if ( ( beamAPtr->isGamma() || beamBPtr->isGamma() ) || beamHasGamma )
148 setBeamModes(isSoftQCD(),
false);
151 beamAhasResGamma = beamAPtr->hasResGamma();
152 beamBhasResGamma = beamBPtr->hasResGamma();
153 beamHasResGamma = beamAhasResGamma || beamBhasResGamma;
156 if (phaseSpacePtr) registerSubObject(*phaseSpacePtr);
157 phaseSpacePtr->init( isFirst, sigmaProcessPtr);
160 if ( beamAhasResGamma || beamBhasResGamma )
161 phaseSpacePtr->setGammaKinPtr( gammaKinPtr);
182 normVar3 = (lhaStratAbs == 3)
183 ? pow2( lhaUpPtr->xErrSum() / lhaUpPtr->xSecSum()) : 0.;
186 sigmaProcessPtr->initProc();
187 if (!sigmaProcessPtr->initFlux())
return false;
190 bool physical = phaseSpacePtr->setupSampling();
191 sigmaMx = phaseSpacePtr->sigmaMax();
192 double sigmaHalfWay = sigmaMx;
195 sigmaSgn = phaseSpacePtr->sigmaSumSigned();
198 if (physical & !isLHA) {
199 int nSample = (nFin < 3) ? N12SAMPLE : N3SAMPLE;
200 for (
int iSample = 0; iSample < nSample; ++iSample) {
202 while (!test) test = phaseSpacePtr->trialKin(
false);
203 if (iSample == nSample/2) sigmaHalfWay = phaseSpacePtr->sigmaMax();
205 double sigmaFullWay = phaseSpacePtr->sigmaMax();
206 sigmaMx = (sigmaHalfWay > 0.) ? pow2(sigmaFullWay) / sigmaHalfWay
208 phaseSpacePtr->setSigmaMax(sigmaMx);
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]);
230 bool ProcessContainer::trialProcess() {
238 if ( beamAPtr->isGamma() || beamBPtr->isGamma() || beamHasGamma )
242 for (
int iTry = 0; ; ++iTry) {
245 if (sigmaMx == 0.)
return false;
246 infoPtr->setEndOfFile(
false);
247 bool repeatSame = (iTry > 0);
248 bool physical = phaseSpacePtr->trialKin(
true, repeatSame);
251 if ( physical && !isSoftQCD() && beamHasGamma ) {
254 if ( !beamAhasResGamma ) beamAPtr->xGamma( x1());
255 if ( !beamBhasResGamma ) beamBPtr->xGamma( x2());
258 if ( !gammaKinPtr->sampleKTgamma() ) physical =
false;
261 if ( physical && !(beamAhasResGamma && beamBhasResGamma) ) {
262 double sHatNew = gammaKinPtr->calcNewSHat( phaseSpacePtr->sHat() );
263 phaseSpacePtr->rescaleSigma( sHatNew);
267 if ( physical && beamHasGamma && approximatedGammaFlux ) {
268 wtPDF = phaseSpacePtr->weightGammaPDFApprox();
269 wtFlux = gammaKinPtr->fluxWeight();
278 bool doTryNext =
true;
281 if (isLHA && !physical) infoPtr->setEndOfFile(
true);
285 if (nSelRequested > 0) doTryNext = (nTry < nSelRequested);
288 if (doTryNext) ++nTry;
292 int codeLHANow = lhaUpPtr->idProcess();
294 for (
int i = 0; i < int(codeLHA.size()); ++i)
295 if (codeLHANow == codeLHA[i]) iFill = i;
298 if (doTryNext) ++nTryLHA[iFill];
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]);
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");
325 if (!physical)
return false;
326 double sigmaNow = phaseSpacePtr->sigmaNow();
329 if (beamHasGamma && approximatedGammaFlux) sigmaNow *= wtFlux * wtPDF;
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;
339 double biasWeight = phaseSpacePtr->biasSelectionWeight();
340 weightNow = sigmaWeight * biasWeight;
343 if (!doTryNext) weightNow = 0.0;
345 infoPtr->setWeight( weightNow, lhaStrat);
349 if (sigmaNow < sigmaNeg) {
350 infoPtr->errorMsg(
"Warning in ProcessContainer::trialProcess: neg"
351 "ative cross section set 0",
"for " + sigmaProcessPtr->name() );
354 if (sigmaNow < 0.) sigmaNow = 0.;
358 double sigmaAdd = sigmaNow * biasWeight;
359 if (lhaStratAbs == 2 || lhaStratAbs == 3) sigmaAdd = sigmaSgn;
371 if (lhaStratAbs >= 3 ) {
372 sigmaTemp = sigmaAdd;
373 sigma2Temp = pow2(sigmaAdd);
375 sigmaTemp += sigmaAdd;
376 sigma2Temp += pow2(sigmaAdd);
380 newSigmaMx = phaseSpacePtr->newSigmaMax();
381 if (newSigmaMx) sigmaMx = phaseSpacePtr->sigmaMax();
385 if (lhaStratAbs < 3) select
386 = (newSigmaMx || rndmPtr->flat() * abs(sigmaMx) < abs(sigmaNow));
389 if (doTryNext) ++nSel;
392 int codeLHANow = lhaUpPtr->idProcess();
394 for (
int i = 0; i < int(codeLHA.size()); ++i)
395 if (codeLHANow == codeLHA[i]) iFill = i;
397 if (iFill >= 0 && doTryNext) ++nSelLHA[iFill];
400 if (select || lhaStratAbs != 2)
return select;
410 void ProcessContainer::accumulate() {
414 double wgtNow = infoPtr->weight();
419 if (lhaStratAbs == 4) wgtNow /= 1e9;
423 int codeLHANow = lhaUpPtr->idProcess();
425 for (
int i = 0; i < int(codeLHA.size()); ++i)
426 if (codeLHANow == codeLHA[i]) iFill = i;
427 if (iFill >= 0) ++nAccLHA[iFill];
437 bool ProcessContainer::constructState() {
440 if (isResolved && !isNonDiff) sigmaProcessPtr->pickInState();
441 sigmaProcessPtr->setIdColAcol();
445 if ( beamAPtr->isResolvedUnresolved() || beamBPtr->isResolvedUnresolved() )
457 void ProcessContainer::setBeamModes(
bool setVMD,
bool isSampled) {
460 beamAPtr->setGammaMode(beamAgammaMode);
461 beamBPtr->setGammaMode(beamBgammaMode);
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;
471 infoPtr->setGammaMode(gammaModeEvent);
474 if (setVMD && !isSampled) {
475 if (beamAgammaMode > 0) infoPtr->setVMDstateA(
true, 22, 0.,0.);
476 if (beamBgammaMode > 0) infoPtr->setVMDstateB(
true, 22, 0.,0.);
481 if (infoPtr->isVMDstateA()) beamAPtr->setVMDstate(
true,
482 infoPtr->idVMDA(), infoPtr->mVMDA(), infoPtr->scaleVMDA(),
484 if (infoPtr->isVMDstateB()) beamBPtr->setVMDstate(
true,
485 infoPtr->idVMDB(), infoPtr->mVMDB(), infoPtr->scaleVMDB(),
494 bool ProcessContainer::constructProcess(
Event& process,
bool isHardest) {
497 if (!phaseSpacePtr->finalKin())
return false;
498 int nFin = sigmaProcessPtr->nFinal();
501 if (isHardest) infoPtr->setType( name(), code(), nFin, isNonDiff,
502 isResolved, isDiffA, isDiffB, isDiffC, isLHA);
505 if ( beamHasGamma && !isSoftQCD() ) gammaKinPtr->finalize();
508 if ( beamHasGamma && !(beamAhasResGamma && beamBhasResGamma) ) {
509 double sHatNew = infoPtr->sHatNew();
510 phaseSpacePtr->rescaleMomenta( sHatNew);
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. );
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. );
538 process.append( 22, -13, 1, 0, 0, 0, 0, 0,
539 Vec4(0., 0., xGm1*infoPtr->pzA(), xGm1*infoPtr->eA()), 0, 0. );
541 process[1].daughter1(3);
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. );
548 process.append( 22, -13, 2, 0, 0, 0, 0, 0,
549 Vec4(0., 0., xGm2*infoPtr->pzB(), xGm2*infoPtr->eB()), 0, 0. );
551 process[1 + nOffsetGamma].daughter1(3 + nOffsetGamma);
556 if (isNonDiff)
return true;
560 process[1 + nOffsetGamma].daughter1(3 + nOffsetGamma);
561 process[2 + nOffsetGamma].daughter1(4 + nOffsetGamma);
567 process[1].daughters(3, 5);
568 process[2].daughters(3, 5);
572 int idRes = sigmaProcessPtr->idSChannel();
573 if (isResolved && !isLHA) {
579 int m_D2 = (nFin == 1) ? 0 : 4 + nFin;
588 scale = sqrt(Q2Fac());
589 process.scale( scale );
592 int colOffset = process.lastColTag();
593 for (
int i = 1; i <= 2 + nFin; ++i) {
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;
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;
620 int iNow = process.append(
id, status, mother1, mother2,
621 daughter1, daughter2, col, acol, phaseSpacePtr->p(i),
622 phaseSpacePtr->m(i), scale);
626 if (i == 2 && idRes != 0) {
629 if (particleDataPtr->hasAnti(idRes)
630 && process[3].chargeType() + process[4].chargeType() < 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) {
644 }
else if (m_col2 == m_acol1 && m_col1 != m_acol2) {
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;
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);
660 if (process[iNow].tau0() > 0.) process[iNow].tau(
661 process[iNow].tau0() * rndmPtr->exp() );
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));
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));
699 process[3].mothers( 1, 2);
700 process[4].mothers( 1, 2);
701 int id5 = sigmaProcessPtr->id(5);
703 process.append( id5, status5, 1, 2, 0, 0, 0, 0,
704 phaseSpacePtr->p(5), phaseSpacePtr->m(5));
714 int nOffsetSecond = 0;
715 for (
int i = 1; i < lhaUpPtr->sizePart(); ++i)
716 if (lhaUpPtr->status(i) == -1) {
718 if (n21 == 3) nOffsetSecond = i - 1;
720 bool twoHard = (n21 == 4);
725 newPos.reserve(lhaUpPtr->sizePart());
727 for (
int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
729 if (twoHard) newPos.push_back(iNew + 1);
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;
740 scalup = lhaUpPtr->scale();
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) );
749 double scalePr2 = process.scaleSecond();
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;
762 int lhaStatus = lhaUpPtr->status(iOld);
764 if (lhaStatus == 2 || lhaStatus == 3) status = -22;
765 if (lhaStatus == 1) status = 23;
768 int mother1Old = lhaUpPtr->mother1(iOld);
769 int mother2Old = lhaUpPtr->mother2(iOld);
772 for (
int im = 1; im < i; ++im) {
773 if (mother1Old == newPos[im]) mother1 = im + 2;
774 if (mother2Old == newPos[im]) mother2 = im + 2;
778 mother1 = 1 + (iIn - 1)%2;
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;
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;
800 if (daughter2 == daughter1) daughter2 = 0;
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;
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);
820 double pol = lhaUpPtr->spin(iOld);
823 double scaleShow = lhaUpPtr->scale(iOld);
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;
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);
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);
844 if (status == 23) iFinal.push_back( iNow);
848 int iFinalSz = iFinal.size();
849 for (
int iF = 0; iF < iFinalSz; ++iF) {
850 int iMod = iFinal[iF];
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]) )
862 if (iQLmod == 0)
continue;
863 double mNew = mNewM[iQLmod];
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;
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;
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;
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;
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;
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;}
916 bool doneShift =
false;
917 Vec4 pMod = process[iMod].p();
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);
927 if (!process[iRec].isFinal()) {
928 vector<int> iDauRec = process[iRec].daughterListRecursive();
930 Mdau.bst( pRecOld, pRecNew);
931 for (
int iD = 0; iD < int(iDauRec.size()); ++iD)
932 process[iDauRec[iD]].rotbst( Mdau);
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");
949 if (matchInOut && !twoHard) {
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");
960 e1 = min( e1, 0.5 * process[0].e());
961 e2 = min( e2, 0.5 * process[0].e());
972 for (
int i = 3; i < process.size(); ++i) {
973 int iMother = process[i].mother1();
976 if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex() )
977 process[i].vProd( process[i - 1].vProd() );
980 else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
981 process[i].vProd( process[iMother].vDec() );
985 int id1Now = process[3 + nOffsetGamma].id();
986 int id2Now = process[4 + nOffsetGamma].id();
1000 double x1Now, x2Now, Q2FacNow, alphaEM, alphaS, Q2Ren, sHat;
1006 x1Now = phaseSpacePtr->x1();
1007 x2Now = phaseSpacePtr->x2();
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();
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();
1040 for (
int i = 5; i < process.size(); ++i)
1041 if (process[i].mother1() == 3 && process[i].mother2() == 4) {
1043 pTLH += process[i].pT();
1045 if (nFinLH > 0) pTHatL = pTLH / nFinLH;
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());
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();
1070 phi = process[5].phi();
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);
1081 infoPtr->setTypeMPI( code(), pTHatL);
1085 int codeSub = lhaUpPtr->idProcess();
1086 ostringstream nameSub;
1087 nameSub <<
"user process " << codeSub;
1088 infoPtr->setSubType( 0, nameSub.str(), codeSub, nFin);
1101 bool ProcessContainer::constructDecays(
Event& process) {
1105 process.append( 90, -11, 0, 0, 0, 0, 0, 0, Vec4(0., 0., 0., 0.), 0., 0. );
1110 newPos.reserve(lhaUpPtr->sizePart());
1111 newPos.push_back(0);
1112 for (
int iNew = 0; iNew < lhaUpPtr->sizePart(); ++iNew) {
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;
1121 double scale = lhaUpPtr->scale();
1122 process.scale( scale);
1126 for (
int i = 1; i < lhaUpPtr->sizePart(); ++i) {
1127 int iOld = newPos[i];
1128 int id = lhaUpPtr->id(iOld);
1131 int lhaStatus = lhaUpPtr->status(iOld);
1133 if (lhaStatus == 2 || lhaStatus == 3) status = -22;
1134 if (lhaStatus == 1) status = 23;
1137 int mother1Old = lhaUpPtr->mother1(iOld);
1138 int mother2Old = lhaUpPtr->mother2(iOld);
1141 for (
int im = 1; im < i; ++im) {
1142 if (mother1Old == newPos[im]) mother1 = im;
1143 if (mother2Old == newPos[im]) mother2 = im;
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;
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;
1164 if (daughter2 == daughter1) daughter2 = 0;
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;
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);
1182 double pol = lhaUpPtr->spin(iOld);
1185 double scaleNow = scale;
1186 if (mother1 > 0) scaleNow = process[mother1].m();
1189 int iNow = process.append(
id, status, mother1, mother2, daughter1,
1190 daughter2, col1, col2, Vec4(px, py, pz, e), m, scaleNow, pol);
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);
1200 process[0].p( pSum);
1201 process[0].m( pSum.mCalc());
1204 for (
int i = 1; i < process.size(); ++i) {
1205 int iMother = process[i].mother1();
1208 if ( process[i - 1].mother1() == iMother && process[i - 1].hasVertex()
1209 && i > 1) process[i].vProd( process[i - 1].vProd() );
1212 else if ( process[iMother].hasVertex() || process[iMother].tau() > 0.)
1213 process[i].vProd( process[iMother].vDec() );
1225 bool ProcessContainer::decayResonances(
Event& process) {
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;
1241 physical = resDecaysPtr->next( process);
1242 if (!physical)
return false;
1246 newFlavours = ( sigmaProcessPtr->weightDecayFlav( process)
1247 < rndmPtr->flat() );
1251 process.restoreSize();
1252 for (
int i = 0; i < process.size(); ++i)
1253 process[i].status( statusSave[i]);
1257 }
while (newFlavours);
1260 phaseSpacePtr->decayKinematics( process);
1263 if (canVetoResDecay)
1264 newChain = userHooksPtr->doVetoResonanceDecays( process);
1268 process.restoreSize();
1269 for (
int i = 0; i < process.size(); ++i)
1270 process[i].status( statusSave[i]);
1285 void ProcessContainer::reset() {
1307 void ProcessContainer::sigmaDelta() {
1314 if (nAcc == 0)
return;
1318 double wgtNow = infoPtr->weight();
1320 if (lhaStrat == 0) wgtNow = sigmaTemp;
1321 if (lhaStratAbs == 3) wgtNow *= sigmaTemp;
1322 if (lhaStratAbs == 4) wgtNow /= 1e9;
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;
1335 double nTryInv = 1. / nTry;
1336 double nSelInv = 1. / nSel;
1337 double nAccInv = 1. / nAcc;
1338 sigmaAvg = sigmaSum * nTryInv;
1339 double fracAcc = nAcc * nSelInv;
1343 if (lhaStratAbs >= 3 ) fracAcc = 1.;
1344 sigmaFin = sigmaAvg * fracAcc;
1345 deltaFin = sigmaFin;
1346 if (nAcc == 1)
return;
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;
1367 bool SetupContainers::init(vector<ProcessContainer*>& containerPtrs,
1371 if (containerPtrs.size() > 0) {
1372 for (
int i = 0; i < int(containerPtrs.size()); ++i)
1373 delete containerPtrs[i];
1374 containerPtrs.clear();
1376 SigmaProcess* sigmaPtr;
1379 Settings& settings = *infoPtr->settingsPtr;
1380 ParticleData* particleDataPtr = infoPtr->particleDataPtr;
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) );
1389 if (softQCD || settings.flag(
"SoftQCD:elastic")) {
1390 sigmaPtr =
new Sigma0AB2AB;
1391 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1399 if (softQCD || inelastic || settings.flag(
"SoftQCD:doubleDiffractive")) {
1400 sigmaPtr =
new Sigma0AB2XX;
1401 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1403 if (softQCD || inelastic || settings.flag(
"SoftQCD:centralDiffractive")) {
1404 sigmaPtr =
new Sigma0AB2AXB;
1405 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1414 if (hardQCD || settings.flag(
"HardQCD:gg2qqbar")) {
1415 sigmaPtr =
new Sigma2gg2qqbar;
1416 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1418 if (hardQCD || settings.flag(
"HardQCD:qg2qg")) {
1419 sigmaPtr =
new Sigma2qg2qg;
1420 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1422 if (hardQCD || settings.flag(
"HardQCD:qq2qq")) {
1423 sigmaPtr =
new Sigma2qq2qq;
1424 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1426 if (hardQCD || settings.flag(
"HardQCD:qqbar2gg")) {
1427 sigmaPtr =
new Sigma2qqbar2gg;
1428 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1430 if (hardQCD || settings.flag(
"HardQCD:qqbar2qqbarNew")) {
1431 sigmaPtr =
new Sigma2qqbar2qqbarNew;
1432 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1441 if (hardQCD || hardccbar || settings.flag(
"HardQCD:qqbar2ccbar")) {
1442 sigmaPtr =
new Sigma2qqbar2QQbar(4, 122);
1443 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1450 if (hardQCD || hardbbbar || settings.flag(
"HardQCD:qqbar2bbbar")) {
1451 sigmaPtr =
new Sigma2qqbar2QQbar(5, 124);
1452 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1461 if (hardQCD3parton || settings.flag(
"HardQCD:qqbar2ggg")) {
1462 sigmaPtr =
new Sigma3qqbar2ggg;
1463 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1465 if (hardQCD3parton || settings.flag(
"HardQCD:qg2qgg")) {
1466 sigmaPtr =
new Sigma3qg2qgg;
1467 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1469 if (hardQCD3parton || settings.flag(
"HardQCD:qq2qqgDiff")) {
1470 sigmaPtr =
new Sigma3qq2qqgDiff;
1471 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1473 if (hardQCD3parton || settings.flag(
"HardQCD:qq2qqgSame")) {
1474 sigmaPtr =
new Sigma3qq2qqgSame;
1475 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1477 if (hardQCD3parton || settings.flag(
"HardQCD:qqbar2qqbargDiff")) {
1478 sigmaPtr =
new Sigma3qqbar2qqbargDiff;
1479 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1481 if (hardQCD3parton || settings.flag(
"HardQCD:qqbar2qqbargSame")) {
1482 sigmaPtr =
new Sigma3qqbar2qqbargSame;
1483 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1485 if (hardQCD3parton || settings.flag(
"HardQCD:gg2qqbarg")) {
1486 sigmaPtr =
new Sigma3gg2qqbarg;
1487 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1489 if (hardQCD3parton || settings.flag(
"HardQCD:qg2qqqbarDiff")) {
1490 sigmaPtr =
new Sigma3qg2qqqbarDiff;
1491 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1493 if (hardQCD3parton || settings.flag(
"HardQCD:qg2qqqbarSame")) {
1494 sigmaPtr =
new Sigma3qg2qqqbarSame;
1495 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1499 bool promptPhotons = settings.flag(
"PromptPhoton:all");
1501 || settings.flag(
"PromptPhoton:qg2qgamma")) {
1502 sigmaPtr =
new Sigma2qg2qgamma;
1503 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1506 || settings.flag(
"PromptPhoton:qqbar2ggamma")) {
1507 sigmaPtr =
new Sigma2qqbar2ggamma;
1508 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1511 || settings.flag(
"PromptPhoton:gg2ggamma")) {
1512 sigmaPtr =
new Sigma2gg2ggamma;
1513 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1516 || settings.flag(
"PromptPhoton:ffbar2gammagamma")) {
1517 sigmaPtr =
new Sigma2ffbar2gammagamma;
1518 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1521 || settings.flag(
"PromptPhoton:gg2gammagamma")) {
1522 sigmaPtr =
new Sigma2gg2gammagamma;
1523 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1533 if (weakBosonExchanges
1534 || settings.flag(
"WeakBosonExchange:ff2ff(t:W)")) {
1535 sigmaPtr =
new Sigma2ff2fftW;
1536 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1546 if (weakSingleBosons
1547 || settings.flag(
"WeakSingleBoson:ffbar2W")) {
1548 sigmaPtr =
new Sigma1ffbar2W;
1549 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1554 if (settings.flag(
"WeakSingleBoson:ffbar2ffbar(s:gm)")) {
1555 sigmaPtr =
new Sigma2ffbar2ffbarsgm;
1556 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1560 if (settings.flag(
"WeakSingleBoson:ffbar2ffbar(s:gmZ)")) {
1561 sigmaPtr =
new Sigma2ffbar2ffbarsgmZ;
1562 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1564 if (settings.flag(
"WeakSingleBoson:ffbar2ffbar(s:W)")) {
1565 sigmaPtr =
new Sigma2ffbar2ffbarsW;
1566 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1576 if (weakDoubleBosons
1577 || settings.flag(
"WeakDoubleBoson:ffbar2ZW")) {
1578 sigmaPtr =
new Sigma2ffbar2ZW;
1579 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1581 if (weakDoubleBosons
1582 || settings.flag(
"WeakDoubleBoson:ffbar2WW")) {
1583 sigmaPtr =
new Sigma2ffbar2WW;
1584 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1594 if (weakBosonAndPartons
1595 || settings.flag(
"WeakBosonAndParton:qg2gmZq")) {
1596 sigmaPtr =
new Sigma2qg2gmZq;
1597 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1599 if (weakBosonAndPartons
1600 || settings.flag(
"WeakBosonAndParton:ffbar2gmZgm")) {
1601 sigmaPtr =
new Sigma2ffbar2gmZgm;
1602 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1604 if (weakBosonAndPartons
1605 || settings.flag(
"WeakBosonAndParton:fgm2gmZf")) {
1606 sigmaPtr =
new Sigma2fgm2gmZf;
1607 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1609 if (weakBosonAndPartons
1610 || settings.flag(
"WeakBosonAndParton:qqbar2Wg")) {
1611 sigmaPtr =
new Sigma2qqbar2Wg;
1612 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1614 if (weakBosonAndPartons
1615 || settings.flag(
"WeakBosonAndParton:qg2Wq")) {
1616 sigmaPtr =
new Sigma2qg2Wq;
1617 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1619 if (weakBosonAndPartons
1620 || settings.flag(
"WeakBosonAndParton:ffbar2Wgm")) {
1621 sigmaPtr =
new Sigma2ffbar2Wgm;
1622 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1624 if (weakBosonAndPartons
1625 || settings.flag(
"WeakBosonAndParton:fgm2Wf")) {
1626 sigmaPtr =
new Sigma2fgm2Wf;
1627 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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);
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);
1652 bool initGammaGamma = ( (initGammaA || !gammaA) && (initGammaB || !gammaB) )
1653 || (photonMode == 4 && gammaA && gammaB );
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) );
1662 if ( ( photonCollisions || settings.flag(
"PhotonCollision:gmgm2ccbar") )
1663 && initGammaGamma ) {
1664 sigmaPtr =
new Sigma2gmgm2ffbar(4, 262);
1665 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1667 if ( ( photonCollisions || settings.flag(
"PhotonCollision:gmgm2bbbar") )
1668 && initGammaGamma ) {
1669 sigmaPtr =
new Sigma2gmgm2ffbar(5, 263);
1670 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1672 if ( ( photonCollisions || settings.flag(
"PhotonCollision:gmgm2ee") )
1673 && initGammaGamma ) {
1674 sigmaPtr =
new Sigma2gmgm2ffbar(11, 264);
1675 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1677 if ( ( photonCollisions || settings.flag(
"PhotonCollision:gmgm2mumu") )
1678 && initGammaGamma ) {
1679 sigmaPtr =
new Sigma2gmgm2ffbar(13, 265);
1680 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1682 if ( ( photonCollisions || settings.flag(
"PhotonCollision:gmgm2tautau") )
1683 && initGammaGamma ) {
1684 sigmaPtr =
new Sigma2gmgm2ffbar(15, 266);
1685 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1690 if (photonParton || settings.flag(
"PhotonParton:ggm2qqbar")) {
1692 sigmaPtr =
new Sigma2ggm2qqbar(1, 271,
"ggm");
1693 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1696 sigmaPtr =
new Sigma2ggm2qqbar(1, 281,
"gmg");
1697 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1700 if (photonParton || settings.flag(
"PhotonParton:ggm2ccbar")) {
1702 sigmaPtr =
new Sigma2ggm2qqbar(4, 272,
"ggm");
1703 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1706 sigmaPtr =
new Sigma2ggm2qqbar(4, 282,
"gmg");
1707 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1710 if (photonParton || settings.flag(
"PhotonParton:ggm2bbbar")) {
1712 sigmaPtr =
new Sigma2ggm2qqbar(5, 273,
"ggm");
1713 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1716 sigmaPtr =
new Sigma2ggm2qqbar(5, 283,
"gmg");
1717 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1720 if (photonParton || settings.flag(
"PhotonParton:qgm2qg")) {
1722 sigmaPtr =
new Sigma2qgm2qg(274,
"qgm");
1723 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1726 sigmaPtr =
new Sigma2qgm2qg(284,
"gmq");
1727 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1730 if (settings.flag(
"PhotonParton:qgm2qgm")) {
1732 sigmaPtr =
new Sigma2qgm2qgm(275,
"qgm");
1733 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1736 sigmaPtr =
new Sigma2qgm2qgm(285,
"gmq");
1737 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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]) );
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) );
1764 if (tops || settings.flag(
"Top:qqbar2ttbar")) {
1765 sigmaPtr =
new Sigma2qqbar2QQbar(6, 602);
1766 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1768 if (tops || settings.flag(
"Top:qq2tq(t:W)")) {
1769 sigmaPtr =
new Sigma2qq2QqtW(6, 603);
1770 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1772 if (tops || settings.flag(
"Top:ffbar2ttbar(s:gmZ)")) {
1773 sigmaPtr =
new Sigma2ffbar2FFbarsgmZ(6, 604);
1774 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1776 if (tops || settings.flag(
"Top:ffbar2tqbar(s:W)")) {
1777 sigmaPtr =
new Sigma2ffbar2FfbarsW(6, 0, 605);
1778 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1780 if (tops || settings.flag(
"Top:gmgm2ttbar")) {
1781 sigmaPtr =
new Sigma2gmgm2ffbar(6, 606);
1782 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1784 if (tops || settings.flag(
"Top:ggm2ttbar")) {
1786 sigmaPtr =
new Sigma2ggm2qqbar(6, 607,
"ggm");
1787 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1790 sigmaPtr =
new Sigma2ggm2qqbar(6, 617,
"gmg");
1791 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1801 if (bPrimes || settings.flag(
"FourthBottom:qqbar2bPrimebPrimebar")) {
1802 sigmaPtr =
new Sigma2qqbar2QQbar(7, 802);
1803 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1805 if (bPrimes || settings.flag(
"FourthBottom:qq2bPrimeq(t:W)")) {
1806 sigmaPtr =
new Sigma2qq2QqtW(7, 803);
1807 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1809 if (bPrimes || settings.flag(
"FourthBottom:ffbar2bPrimebPrimebar(s:gmZ)")) {
1810 sigmaPtr =
new Sigma2ffbar2FFbarsgmZ(7, 804);
1811 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1813 if (bPrimes || settings.flag(
"FourthBottom:ffbar2bPrimeqbar(s:W)")) {
1814 sigmaPtr =
new Sigma2ffbar2FfbarsW(7, 0, 805);
1815 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1817 if (bPrimes || settings.flag(
"FourthBottom:ffbar2bPrimetbar(s:W)")) {
1818 sigmaPtr =
new Sigma2ffbar2FfbarsW(7, 6, 806);
1819 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1828 if (tPrimes || settings.flag(
"FourthTop:qqbar2tPrimetPrimebar")) {
1829 sigmaPtr =
new Sigma2qqbar2QQbar(8, 822);
1830 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1832 if (tPrimes || settings.flag(
"FourthTop:qq2tPrimeq(t:W)")) {
1833 sigmaPtr =
new Sigma2qq2QqtW(8, 823);
1834 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1836 if (tPrimes || settings.flag(
"FourthTop:ffbar2tPrimetPrimebar(s:gmZ)")) {
1837 sigmaPtr =
new Sigma2ffbar2FFbarsgmZ(8, 824);
1838 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1840 if (tPrimes || settings.flag(
"FourthTop:ffbar2tPrimeqbar(s:W)")) {
1841 sigmaPtr =
new Sigma2ffbar2FfbarsW(8, 0, 825);
1842 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1851 if (settings.flag(
"FourthPair:ffbar2tauPrimenuPrimebar(s:W)")) {
1852 sigmaPtr =
new Sigma2ffbar2FfbarsW(17, 18, 842);
1853 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1857 bool useBSMHiggses = settings.flag(
"Higgs:useBSM");
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) );
1866 if (HiggsesSM || settings.flag(
"HiggsSM:gg2H")) {
1867 sigmaPtr =
new Sigma1gg2H(0);
1868 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1870 if (HiggsesSM || settings.flag(
"HiggsSM:gmgm2H")) {
1871 sigmaPtr =
new Sigma1gmgm2H(0);
1872 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1874 if (HiggsesSM || settings.flag(
"HiggsSM:ffbar2HZ")) {
1875 sigmaPtr =
new Sigma2ffbar2HZ(0);
1876 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1878 if (HiggsesSM || settings.flag(
"HiggsSM:ffbar2HW")) {
1879 sigmaPtr =
new Sigma2ffbar2HW(0);
1880 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1882 if (HiggsesSM || settings.flag(
"HiggsSM:ff2Hff(t:ZZ)")) {
1883 sigmaPtr =
new Sigma3ff2HfftZZ(0);
1884 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1886 if (HiggsesSM || settings.flag(
"HiggsSM:ff2Hff(t:WW)")) {
1887 sigmaPtr =
new Sigma3ff2HfftWW(0);
1888 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1890 if (HiggsesSM || settings.flag(
"HiggsSM:gg2Httbar")) {
1891 sigmaPtr =
new Sigma3gg2HQQbar(6,0);
1892 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1894 if (HiggsesSM || settings.flag(
"HiggsSM:qqbar2Httbar")) {
1895 sigmaPtr =
new Sigma3qqbar2HQQbar(6,0);
1896 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1906 if (settings.flag(
"HiggsSM:gg2Hbbbar")) {
1907 sigmaPtr =
new Sigma3gg2HQQbar(5,0);
1908 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1910 if (settings.flag(
"HiggsSM:qqbar2Hbbbar")) {
1911 sigmaPtr =
new Sigma3qqbar2HQQbar(5,0);
1912 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1914 if (settings.flag(
"HiggsSM:gg2Hg(l:t)")) {
1915 sigmaPtr =
new Sigma2gg2Hglt(0);
1916 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1918 if (settings.flag(
"HiggsSM:qg2Hq(l:t)")) {
1919 sigmaPtr =
new Sigma2qg2Hqlt(0);
1920 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1922 if (settings.flag(
"HiggsSM:qqbar2Hg(l:t)")) {
1923 sigmaPtr =
new Sigma2qqbar2Hglt(0);
1924 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1929 if (useBSMHiggses) {
1930 bool HiggsesBSM = settings.flag(
"HiggsBSM:all");
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) );
1938 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:gg2H1")) {
1939 sigmaPtr =
new Sigma1gg2H(1);
1940 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1942 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:gmgm2H1")) {
1943 sigmaPtr =
new Sigma1gmgm2H(1);
1944 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1946 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:ffbar2H1Z")) {
1947 sigmaPtr =
new Sigma2ffbar2HZ(1);
1948 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1950 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:ffbar2H1W")) {
1951 sigmaPtr =
new Sigma2ffbar2HW(1);
1952 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1954 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:ff2H1ff(t:ZZ)")) {
1955 sigmaPtr =
new Sigma3ff2HfftZZ(1);
1956 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1958 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:ff2H1ff(t:WW)")) {
1959 sigmaPtr =
new Sigma3ff2HfftWW(1);
1960 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1962 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:gg2H1ttbar")) {
1963 sigmaPtr =
new Sigma3gg2HQQbar(6,1);
1964 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1966 if (HiggsesBSM || HiggsesH1 || settings.flag(
"HiggsBSM:qqbar2H1ttbar")) {
1967 sigmaPtr =
new Sigma3qqbar2HQQbar(6,1);
1968 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
1978 if (settings.flag(
"HiggsBSM:gg2H1bbbar")) {
1979 sigmaPtr =
new Sigma3gg2HQQbar(5,1);
1980 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1982 if (settings.flag(
"HiggsBSM:qqbar2H1bbbar")) {
1983 sigmaPtr =
new Sigma3qqbar2HQQbar(5,1);
1984 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1986 if (settings.flag(
"HiggsBSM:gg2H1g(l:t)")) {
1987 sigmaPtr =
new Sigma2gg2Hglt(1);
1988 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1990 if (settings.flag(
"HiggsBSM:qg2H1q(l:t)")) {
1991 sigmaPtr =
new Sigma2qg2Hqlt(1);
1992 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
1994 if (settings.flag(
"HiggsBSM:qqbar2H1g(l:t)")) {
1995 sigmaPtr =
new Sigma2qqbar2Hglt(1);
1996 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2005 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:gg2H2")) {
2006 sigmaPtr =
new Sigma1gg2H(2);
2007 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2009 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:gmgm2H2")) {
2010 sigmaPtr =
new Sigma1gmgm2H(2);
2011 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2013 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:ffbar2H2Z")) {
2014 sigmaPtr =
new Sigma2ffbar2HZ(2);
2015 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2017 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:ffbar2H2W")) {
2018 sigmaPtr =
new Sigma2ffbar2HW(2);
2019 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2021 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:ff2H2ff(t:ZZ)")) {
2022 sigmaPtr =
new Sigma3ff2HfftZZ(2);
2023 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2025 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:ff2H2ff(t:WW)")) {
2026 sigmaPtr =
new Sigma3ff2HfftWW(2);
2027 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2029 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:gg2H2ttbar")) {
2030 sigmaPtr =
new Sigma3gg2HQQbar(6,2);
2031 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2033 if (HiggsesBSM || HiggsesH2 || settings.flag(
"HiggsBSM:qqbar2H2ttbar")) {
2034 sigmaPtr =
new Sigma3qqbar2HQQbar(6,2);
2035 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2045 if (settings.flag(
"HiggsBSM:gg2H2bbbar")) {
2046 sigmaPtr =
new Sigma3gg2HQQbar(5,2);
2047 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2049 if (settings.flag(
"HiggsBSM:qqbar2H2bbbar")) {
2050 sigmaPtr =
new Sigma3qqbar2HQQbar(5,2);
2051 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2053 if (settings.flag(
"HiggsBSM:gg2H2g(l:t)")) {
2054 sigmaPtr =
new Sigma2gg2Hglt(2);
2055 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2057 if (settings.flag(
"HiggsBSM:qg2H2q(l:t)")) {
2058 sigmaPtr =
new Sigma2qg2Hqlt(2);
2059 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2061 if (settings.flag(
"HiggsBSM:qqbar2H2g(l:t)")) {
2062 sigmaPtr =
new Sigma2qqbar2Hglt(2);
2063 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2072 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:gg2A3")) {
2073 sigmaPtr =
new Sigma1gg2H(3);
2074 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2076 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:gmgm2A3")) {
2077 sigmaPtr =
new Sigma1gmgm2H(3);
2078 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2080 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:ffbar2A3Z")) {
2081 sigmaPtr =
new Sigma2ffbar2HZ(3);
2082 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2084 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:ffbar2A3W")) {
2085 sigmaPtr =
new Sigma2ffbar2HW(3);
2086 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2088 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:ff2A3ff(t:ZZ)")) {
2089 sigmaPtr =
new Sigma3ff2HfftZZ(3);
2090 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2092 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:ff2A3ff(t:WW)")) {
2093 sigmaPtr =
new Sigma3ff2HfftWW(3);
2094 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2096 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:gg2A3ttbar")) {
2097 sigmaPtr =
new Sigma3gg2HQQbar(6,3);
2098 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2100 if (HiggsesBSM || HiggsesA3 || settings.flag(
"HiggsBSM:qqbar2A3ttbar")) {
2101 sigmaPtr =
new Sigma3qqbar2HQQbar(6,3);
2102 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2112 if (settings.flag(
"HiggsBSM:gg2A3bbbar")) {
2113 sigmaPtr =
new Sigma3gg2HQQbar(5,3);
2114 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2116 if (settings.flag(
"HiggsBSM:qqbar2A3bbbar")) {
2117 sigmaPtr =
new Sigma3qqbar2HQQbar(5,3);
2118 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2120 if (settings.flag(
"HiggsBSM:gg2A3g(l:t)")) {
2121 sigmaPtr =
new Sigma2gg2Hglt(3);
2122 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2124 if (settings.flag(
"HiggsBSM:qg2A3q(l:t)")) {
2125 sigmaPtr =
new Sigma2qg2Hqlt(3);
2126 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2128 if (settings.flag(
"HiggsBSM:qqbar2A3g(l:t)")) {
2129 sigmaPtr =
new Sigma2qqbar2Hglt(3);
2130 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
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) );
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) );
2150 if (HiggsesBSM || HiggsesPairs || settings.flag(
"HiggsBSM:ffbar2A3H2")) {
2151 sigmaPtr =
new Sigma2ffbar2A3H12(2);
2152 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2154 if (HiggsesBSM || HiggsesPairs || settings.flag(
"HiggsBSM:ffbar2H+-H1")) {
2155 sigmaPtr =
new Sigma2ffbar2HchgH12(1);
2156 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2158 if (HiggsesBSM || HiggsesPairs || settings.flag(
"HiggsBSM:ffbar2H+-H2")) {
2159 sigmaPtr =
new Sigma2ffbar2HchgH12(2);
2160 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2162 if (HiggsesBSM || HiggsesPairs || settings.flag(
"HiggsBSM:ffbar2H+H-")) {
2163 sigmaPtr =
new Sigma2ffbar2HposHneg();
2164 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2169 CoupSUSY* coupSUSYPtr = infoPtr->coupSUSYPtr;
2170 if (coupSUSYPtr->isSUSY) {
2172 bool SUSYs = settings.flag(
"SUSY:all");
2173 bool nmssm = settings.flag(
"SLHA:NMSSM");
2176 setupIdVecs( settings);
2180 if (nmssm) nNeut = 5;
2183 if (SUSYs || settings.flag(
"SUSY:gg2gluinogluino")) {
2185 if (allowIdVals( 1000021, 1000021)) {
2186 sigmaPtr =
new Sigma2gg2gluinogluino();
2187 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2190 if (SUSYs || settings.flag(
"SUSY:qqbar2gluinogluino")) {
2192 if (allowIdVals( 1000021, 1000021)) {
2193 sigmaPtr =
new Sigma2qqbar2gluinogluino();
2194 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2199 if (SUSYs || settings.flag(
"SUSY:qg2squarkgluino")) {
2201 for (
int idx = 1; idx <= 6; ++idx) {
2202 for (
int iso = 1; iso <= 2; ++iso) {
2204 int id3 = iso + ((idx <= 3)
2205 ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
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) );
2218 if (SUSYs || settings.flag(
"SUSY:gg2squarkantisquark")) {
2220 for (
int idx = 1; idx <= 6; ++idx) {
2221 for (
int iso = 1; iso <= 2; ++iso) {
2223 int id = iso + ((idx <= 3)
2224 ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2226 if (!allowIdVals(
id,
id))
continue;
2227 sigmaPtr =
new Sigma2gg2squarkantisquark(
id,iproc);
2228 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2234 if (SUSYs || settings.flag(
"SUSY:qqbar2squarkantisquark")) {
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));
2248 if (iso == jso && id1 != id2) iproc++;
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) );
2257 sigmaPtr =
new Sigma2qqbar2squarkantisquark(id1,-id2,iproc);
2258 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2267 if (SUSYs || settings.flag(
"SUSY:qq2squarksquark")) {
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;
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));
2280 if (!allowIdVals( id1, id2))
continue;
2281 sigmaPtr =
new Sigma2qq2squarksquark(id1,id2,iproc);
2282 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2290 if (SUSYs || settings.flag(
"SUSY:qg2chi0squark")) {
2292 for (
int iNeut= 1; iNeut <= nNeut; iNeut++) {
2293 for (
int idx = 1; idx <= 6; idx++) {
2295 for (
int iso = 1; iso <= 2; iso++) {
2296 if (iso == 2) isUp =
true;
2298 int id3 = coupSUSYPtr->idNeut(iNeut);
2299 int id4 = iso + ((idx <= 3)
2300 ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2302 if (!allowIdVals( id3, id4))
continue;
2303 sigmaPtr =
new Sigma2qg2chi0squark(iNeut,idx,isUp,iproc);
2304 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2311 if (SUSYs || settings.flag(
"SUSY:qg2chi+-squark")) {
2313 for (
int iChar = 1; iChar <= 2; iChar++) {
2314 for (
int idx = 1; idx <= 6; idx++) {
2316 for (
int iso = 1; iso <= 2; iso++) {
2317 if (iso == 2) isUp =
true;
2319 int id3 = coupSUSYPtr->idChar(iChar);
2320 int id4 = iso + ((idx <= 3)
2321 ? 1000000+2*(idx-1) : 2000000+2*(idx-4));
2323 if (!allowIdVals( id3, id4))
continue;
2324 sigmaPtr =
new Sigma2qg2charsquark(iChar,idx,isUp,iproc);
2325 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2332 if (SUSYs || settings.flag(
"SUSY:qqbar2chi0chi0")) {
2334 for (
int iNeut2 = 1; iNeut2 <= nNeut; iNeut2++) {
2335 for (
int iNeut1 = 1; iNeut1 <= iNeut2; iNeut1++) {
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) );
2347 if (SUSYs || settings.flag(
"SUSY:qqbar2chi+-chi0")) {
2349 for (
int iNeut = 1; iNeut <= nNeut; iNeut++) {
2350 for (
int iChar = 1; iChar <= 2; ++iChar) {
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) );
2364 if (SUSYs || settings.flag(
"SUSY:qqbar2chi+chi-")) {
2366 for (
int i = 1; i <= 2; ++i) {
2367 for (
int j = 1; j <= 2; ++j) {
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) );
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));
2384 if (!allowIdVals( id1, 0))
continue;
2385 sigmaPtr =
new Sigma1qq2antisquark(id1);
2386 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2392 if (SUSYs || settings.flag(
"SUSY:qqbar2chi0gluino")) {
2394 for (
int iNeut = 1; iNeut <= nNeut; iNeut++) {
2397 if (!allowIdVals( coupSUSYPtr->idNeut(iNeut), 1000021))
continue;
2398 sigmaPtr =
new Sigma2qqbar2chi0gluino(iNeut, iproc);
2399 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2404 if (SUSYs || settings.flag(
"SUSY:qqbar2chi+-gluino")) {
2406 for (
int iChar = 1; iChar <= 2; ++iChar) {
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) );
2419 if (SUSYs || settings.flag(
"SUSY:qqbar2sleptonantislepton")) {
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));
2432 if (iso == jso && id1 != id2) iproc++;
2434 if (abs(id1) >= 2000012 && id1 % 2 == 0)
continue;
2435 if (abs(id2) >= 2000012 && id2 % 2 == 0)
continue;
2437 if (!allowIdVals( id1, id2))
continue;
2438 if (iso == jso && id1 != id2) {
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) );
2445 sigmaPtr =
new Sigma2qqbar2sleptonantislepton(id1,-id2,iproc);
2446 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2457 if (settings.flag(
"NewGaugeBoson:ffbar2gmZZprime")) {
2458 sigmaPtr =
new Sigma1ffbar2gmZZprime();
2459 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2461 if (settings.flag(
"NewGaugeBoson:ffbar2Wprime")) {
2462 sigmaPtr =
new Sigma1ffbar2Wprime();
2463 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2465 if (settings.flag(
"NewGaugeBoson:ffbar2R0")) {
2466 sigmaPtr =
new Sigma1ffbar2Rhorizontal();
2467 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2476 if (leftrights || settings.flag(
"LeftRightSymmmetry:ffbar2WR")) {
2477 sigmaPtr =
new Sigma1ffbar2WRight();
2478 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2480 if (leftrights || settings.flag(
"LeftRightSymmmetry:ll2HL")) {
2481 sigmaPtr =
new Sigma1ll2Hchgchg(1);
2482 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2484 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HLe")) {
2485 sigmaPtr =
new Sigma2lgm2Hchgchgl(1, 11);
2486 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2488 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HLmu")) {
2489 sigmaPtr =
new Sigma2lgm2Hchgchgl(1, 13);
2490 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2492 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HLtau")) {
2493 sigmaPtr =
new Sigma2lgm2Hchgchgl(1, 15);
2494 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2496 if (leftrights || settings.flag(
"LeftRightSymmmetry:ff2HLff")) {
2497 sigmaPtr =
new Sigma3ff2HchgchgfftWW(1);
2498 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2500 if (leftrights || settings.flag(
"LeftRightSymmmetry:ffbar2HLHL")) {
2501 sigmaPtr =
new Sigma2ffbar2HchgchgHchgchg(1);
2502 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2504 if (leftrights || settings.flag(
"LeftRightSymmmetry:ll2HR")) {
2505 sigmaPtr =
new Sigma1ll2Hchgchg(2);
2506 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2508 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HRe")) {
2509 sigmaPtr =
new Sigma2lgm2Hchgchgl(2, 11);
2510 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2512 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HRmu")) {
2513 sigmaPtr =
new Sigma2lgm2Hchgchgl(2, 13);
2514 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2516 if (leftrights || settings.flag(
"LeftRightSymmmetry:lgm2HRtau")) {
2517 sigmaPtr =
new Sigma2lgm2Hchgchgl(2, 15);
2518 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2520 if (leftrights || settings.flag(
"LeftRightSymmmetry:ff2HRff")) {
2521 sigmaPtr =
new Sigma3ff2HchgchgfftWW(2);
2522 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2524 if (leftrights || settings.flag(
"LeftRightSymmmetry:ffbar2HRHR")) {
2525 sigmaPtr =
new Sigma2ffbar2HchgchgHchgchg(2);
2526 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2535 if (leptoquarks || settings.flag(
"LeptoQuark:qg2LQl")) {
2536 sigmaPtr =
new Sigma2qg2LeptoQuarkl;
2537 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2539 if (leptoquarks || settings.flag(
"LeptoQuark:gg2LQLQbar")) {
2540 sigmaPtr =
new Sigma2gg2LQLQbar;
2541 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2543 if (leptoquarks || settings.flag(
"LeptoQuark:qqbar2LQLQbar")) {
2544 sigmaPtr =
new Sigma2qqbar2LQLQbar;
2545 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2554 if (excitedfermions || settings.flag(
"ExcitedFermion:ug2uStar")) {
2555 sigmaPtr =
new Sigma1qg2qStar(2);
2556 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2558 if (excitedfermions || settings.flag(
"ExcitedFermion:sg2sStar")) {
2559 sigmaPtr =
new Sigma1qg2qStar(3);
2560 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2562 if (excitedfermions || settings.flag(
"ExcitedFermion:cg2cStar")) {
2563 sigmaPtr =
new Sigma1qg2qStar(4);
2564 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2566 if (excitedfermions || settings.flag(
"ExcitedFermion:bg2bStar")) {
2567 sigmaPtr =
new Sigma1qg2qStar(5);
2568 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2570 if (excitedfermions || settings.flag(
"ExcitedFermion:egm2eStar")) {
2571 sigmaPtr =
new Sigma1lgm2lStar(11);
2572 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2574 if (excitedfermions || settings.flag(
"ExcitedFermion:mugm2muStar")) {
2575 sigmaPtr =
new Sigma1lgm2lStar(13);
2576 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2578 if (excitedfermions || settings.flag(
"ExcitedFermion:taugm2tauStar")) {
2579 sigmaPtr =
new Sigma1lgm2lStar(15);
2580 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2582 if (excitedfermions || settings.flag(
"ExcitedFermion:qq2dStarq")) {
2583 sigmaPtr =
new Sigma2qq2qStarq(1);
2584 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2586 if (excitedfermions || settings.flag(
"ExcitedFermion:qq2uStarq")) {
2587 sigmaPtr =
new Sigma2qq2qStarq(2);
2588 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2590 if (excitedfermions || settings.flag(
"ExcitedFermion:qq2sStarq")) {
2591 sigmaPtr =
new Sigma2qq2qStarq(3);
2592 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2594 if (excitedfermions || settings.flag(
"ExcitedFermion:qq2cStarq")) {
2595 sigmaPtr =
new Sigma2qq2qStarq(4);
2596 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2598 if (excitedfermions || settings.flag(
"ExcitedFermion:qq2bStarq")) {
2599 sigmaPtr =
new Sigma2qq2qStarq(5);
2600 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2602 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2eStare")) {
2603 sigmaPtr =
new Sigma2qqbar2lStarlbar(11);
2604 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2606 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2nueStarnue")) {
2607 sigmaPtr =
new Sigma2qqbar2lStarlbar(12);
2608 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2610 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2muStarmu")) {
2611 sigmaPtr =
new Sigma2qqbar2lStarlbar(13);
2612 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2614 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2numuStarnumu")) {
2615 sigmaPtr =
new Sigma2qqbar2lStarlbar(14);
2616 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2618 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2tauStartau")) {
2619 sigmaPtr =
new Sigma2qqbar2lStarlbar(15);
2620 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2623 || settings.flag(
"ExcitedFermion:qqbar2nutauStarnutau")) {
2624 sigmaPtr =
new Sigma2qqbar2lStarlbar(16);
2625 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2627 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2eStareStar")) {
2628 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(11);
2629 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2632 || settings.flag(
"ExcitedFermion:qqbar2nueStarnueStar")) {
2633 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(12);
2634 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2636 if (excitedfermions || settings.flag(
"ExcitedFermion:qqbar2muStarmuStar")) {
2637 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(13);
2638 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2641 || settings.flag(
"ExcitedFermion:qqbar2numuStarnumuStar")) {
2642 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(14);
2643 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2646 || settings.flag(
"ExcitedFermion:qqbar2tauStartauStar")) {
2647 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(15);
2648 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2651 || settings.flag(
"ExcitedFermion:qqbar2nutauStarnutauStar")) {
2652 sigmaPtr =
new Sigma2qqbar2lStarlStarBar(16);
2653 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2657 if (settings.flag(
"ContactInteractions:QCqq2qq")) {
2658 sigmaPtr =
new Sigma2QCqq2qq();
2659 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2661 if (settings.flag(
"ContactInteractions:QCqqbar2qqbar")) {
2662 sigmaPtr =
new Sigma2QCqqbar2qqbar();
2663 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2665 if (settings.flag(
"ContactInteractions:QCffbar2eebar")) {
2666 sigmaPtr =
new Sigma2QCffbar2llbar(11, 4203);
2667 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2669 if (settings.flag(
"ContactInteractions:QCffbar2mumubar")) {
2670 sigmaPtr =
new Sigma2QCffbar2llbar(13, 4204);
2671 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2673 if (settings.flag(
"ContactInteractions:QCffbar2tautaubar")) {
2674 sigmaPtr =
new Sigma2QCffbar2llbar(15, 4205);
2675 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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);
2687 if (particleDataPtr->spinType( 4900101) != 2)
2688 particleDataPtr->spinType( 4900101, 2);
2690 int spinqv = settings.mode(
"HiddenValley:spinqv");
2691 if (particleDataPtr->spinType( 4900101) != 2 * spinqv + 1)
2692 particleDataPtr->spinType( 4900101, 2 * spinqv + 1);
2696 bool hiddenvalleys = settings.flag(
"HiddenValley:all");
2697 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2DvDvbar")) {
2698 sigmaPtr =
new Sigma2gg2qGqGbar( 4900001, 4901, spinFv,
2700 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2702 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2UvUvbar")) {
2703 sigmaPtr =
new Sigma2gg2qGqGbar( 4900002, 4902, spinFv,
2705 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2707 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2SvSvbar")) {
2708 sigmaPtr =
new Sigma2gg2qGqGbar( 4900003, 4903, spinFv,
2710 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2712 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2CvCvbar")) {
2713 sigmaPtr =
new Sigma2gg2qGqGbar( 4900004, 4904, spinFv,
2715 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2717 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2BvBvbar")) {
2718 sigmaPtr =
new Sigma2gg2qGqGbar( 4900005, 4905, spinFv,
2720 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2722 if (hiddenvalleys || settings.flag(
"HiddenValley:gg2TvTvbar")) {
2723 sigmaPtr =
new Sigma2gg2qGqGbar( 4900006, 4906, spinFv,
2725 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
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) );
2817 if (hiddenvalleys || settings.flag(
"HiddenValley:ffbar2Zv")) {
2818 sigmaPtr =
new Sigma1ffbar2Zv();
2819 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2828 if (extraDimGstars || settings.flag(
"ExtraDimensionsG*:ffbar2G*")) {
2829 sigmaPtr =
new Sigma1ffbar2GravitonStar;
2830 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2832 if (settings.flag(
"ExtraDimensionsG*:gg2G*g")) {
2833 sigmaPtr =
new Sigma2gg2GravitonStarg;
2834 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2836 if (settings.flag(
"ExtraDimensionsG*:qg2G*q")) {
2837 sigmaPtr =
new Sigma2qg2GravitonStarq;
2838 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2840 if (settings.flag(
"ExtraDimensionsG*:qqbar2G*g")) {
2841 sigmaPtr =
new Sigma2qqbar2GravitonStarg;
2842 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2846 if (settings.flag(
"ExtraDimensionsG*:qqbar2KKgluon*")) {
2847 sigmaPtr =
new Sigma1qqbar2KKgluonStar;
2848 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2852 if (settings.flag(
"ExtraDimensionsTEV:ffbar2ddbar")) {
2853 sigmaPtr =
new Sigma2ffbar2TEVffbar(1, 5061);
2854 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2856 if (settings.flag(
"ExtraDimensionsTEV:ffbar2uubar")) {
2857 sigmaPtr =
new Sigma2ffbar2TEVffbar(2, 5062);
2858 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2860 if (settings.flag(
"ExtraDimensionsTEV:ffbar2ssbar")) {
2861 sigmaPtr =
new Sigma2ffbar2TEVffbar(3, 5063);
2862 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2864 if (settings.flag(
"ExtraDimensionsTEV:ffbar2ccbar")) {
2865 sigmaPtr =
new Sigma2ffbar2TEVffbar(4, 5064);
2866 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2868 if (settings.flag(
"ExtraDimensionsTEV:ffbar2bbbar")) {
2869 sigmaPtr =
new Sigma2ffbar2TEVffbar(5, 5065);
2870 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2872 if (settings.flag(
"ExtraDimensionsTEV:ffbar2ttbar")) {
2873 sigmaPtr =
new Sigma2ffbar2TEVffbar(6, 5066);
2874 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2876 if (settings.flag(
"ExtraDimensionsTEV:ffbar2e+e-")) {
2877 sigmaPtr =
new Sigma2ffbar2TEVffbar(11, 5071);
2878 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2880 if (settings.flag(
"ExtraDimensionsTEV:ffbar2nuenuebar")) {
2881 sigmaPtr =
new Sigma2ffbar2TEVffbar(12, 5072);
2882 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2884 if (settings.flag(
"ExtraDimensionsTEV:ffbar2mu+mu-")) {
2885 sigmaPtr =
new Sigma2ffbar2TEVffbar(13, 5073);
2886 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2888 if (settings.flag(
"ExtraDimensionsTEV:ffbar2numunumubar")) {
2889 sigmaPtr =
new Sigma2ffbar2TEVffbar(14, 5074);
2890 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2892 if (settings.flag(
"ExtraDimensionsTEV:ffbar2tau+tau-")) {
2893 sigmaPtr =
new Sigma2ffbar2TEVffbar(15, 5075);
2894 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2896 if (settings.flag(
"ExtraDimensionsTEV:ffbar2nutaunutaubar")) {
2897 sigmaPtr =
new Sigma2ffbar2TEVffbar(16, 5076);
2898 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2907 if (extraDimLEDmono || settings.flag(
"ExtraDimensionsLED:qg2Gq")) {
2908 sigmaPtr =
new Sigma2qg2LEDUnparticleq(
true );
2909 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2911 if (extraDimLEDmono || settings.flag(
"ExtraDimensionsLED:qqbar2Gg")) {
2912 sigmaPtr =
new Sigma2qqbar2LEDUnparticleg(
true );
2913 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2915 if (settings.flag(
"ExtraDimensionsLED:ffbar2GZ")) {
2916 sigmaPtr =
new Sigma2ffbar2LEDUnparticleZ(
true );
2917 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2919 if (settings.flag(
"ExtraDimensionsLED:ffbar2Ggamma")) {
2920 sigmaPtr =
new Sigma2ffbar2LEDUnparticlegamma(
true );
2921 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2923 if (settings.flag(
"ExtraDimensionsLED:ffbar2gammagamma")) {
2924 sigmaPtr =
new Sigma2ffbar2LEDgammagamma(
true );
2925 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2927 if (settings.flag(
"ExtraDimensionsLED:gg2gammagamma")) {
2928 sigmaPtr =
new Sigma2gg2LEDgammagamma(
true );
2929 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2931 if (settings.flag(
"ExtraDimensionsLED:ffbar2llbar")) {
2932 sigmaPtr =
new Sigma2ffbar2LEDllbar(
true );
2933 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2935 if (settings.flag(
"ExtraDimensionsLED:gg2llbar")) {
2936 sigmaPtr =
new Sigma2gg2LEDllbar(
true );
2937 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2944 if (extraDimLEDdij || settings.flag(
"ExtraDimensionsLED:gg2DJqqbar")) {
2945 sigmaPtr =
new Sigma2gg2LEDqqbar;
2946 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2948 if (extraDimLEDdij || settings.flag(
"ExtraDimensionsLED:qg2DJqg")) {
2949 sigmaPtr =
new Sigma2qg2LEDqg;
2950 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2952 if (extraDimLEDdij || settings.flag(
"ExtraDimensionsLED:qq2DJqq")) {
2953 sigmaPtr =
new Sigma2qq2LEDqq;
2954 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2956 if (extraDimLEDdij || settings.flag(
"ExtraDimensionsLED:qqbar2DJgg")) {
2957 sigmaPtr =
new Sigma2qqbar2LEDgg;
2958 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2960 if (extraDimLEDdij || settings.flag(
"ExtraDimensionsLED:qqbar2DJqqbarNew")) {
2961 sigmaPtr =
new Sigma2qqbar2LEDqqbarNew;
2962 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
2971 if (extraDimUnpartmono || settings.flag(
"ExtraDimensionsUnpart:qg2Uq")) {
2972 sigmaPtr =
new Sigma2qg2LEDUnparticleq(
false );
2973 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2975 if (extraDimUnpartmono || settings.flag(
"ExtraDimensionsUnpart:qqbar2Ug")) {
2976 sigmaPtr =
new Sigma2qqbar2LEDUnparticleg(
false );
2977 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2979 if (settings.flag(
"ExtraDimensionsUnpart:ffbar2UZ")) {
2980 sigmaPtr =
new Sigma2ffbar2LEDUnparticleZ(
false );
2981 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2983 if (settings.flag(
"ExtraDimensionsUnpart:ffbar2Ugamma")) {
2984 sigmaPtr =
new Sigma2ffbar2LEDUnparticlegamma(
false );
2985 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2987 if (settings.flag(
"ExtraDimensionsUnpart:ffbar2gammagamma")) {
2988 sigmaPtr =
new Sigma2ffbar2LEDgammagamma(
false );
2989 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2991 if (settings.flag(
"ExtraDimensionsUnpart:gg2gammagamma")) {
2992 sigmaPtr =
new Sigma2gg2LEDgammagamma(
false );
2993 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2995 if (settings.flag(
"ExtraDimensionsUnpart:ffbar2llbar")) {
2996 sigmaPtr =
new Sigma2ffbar2LEDllbar(
false );
2997 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
2999 if (settings.flag(
"ExtraDimensionsUnpart:gg2llbar")) {
3000 sigmaPtr =
new Sigma2gg2LEDllbar(
false );
3001 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3005 if (settings.flag(
"DM:ffbar2Zp2XX")) {
3006 sigmaPtr =
new Sigma1ffbar2Zp2XX();
3007 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3009 if (settings.flag(
"DM:ffbar2Zp2XXj")) {
3010 sigmaPtr =
new Sigma2qqbar2Zpg2XXj();
3011 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3013 if (settings.flag(
"DM:ffbar2ZpH")) {
3014 sigmaPtr =
new Sigma2ffbar2ZpH();
3015 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3017 if (settings.flag(
"DM:gg2S2XX")) {
3018 sigmaPtr =
new Sigma1gg2S2XX();
3019 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3021 if (settings.flag(
"DM:gg2S2XXj")) {
3022 sigmaPtr =
new Sigma2gg2Sg2XXj();
3023 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3025 if (settings.flag(
"DM:qqbar2DY")) {
3026 sigmaPtr =
new Sigma2qqbar2DY();
3027 containerPtrs.push_back(
new ProcessContainer(sigmaPtr) );
3030 for ( ProcessContainer * cont : containerPtrs )
3031 cont->initInfoPtr(*infoPtr);
3042 bool SetupContainers::init2(vector<ProcessContainer*>& container2Ptrs,
3045 Settings& settings = *infoPtr->settingsPtr;
3048 if (container2Ptrs.size() > 0) {
3049 for (
int i = 0; i < int(container2Ptrs.size()); ++i)
3050 delete container2Ptrs[i];
3051 container2Ptrs.clear();
3053 SigmaProcess* sigmaPtr;
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) );
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) );
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) );
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]));
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]));
3118 if (settings.flag(
"SecondHard:SingleGmZ")) {
3119 sigmaPtr =
new Sigma1ffbar2gmZ;
3120 container2Ptrs.push_back(
new ProcessContainer(sigmaPtr) );
3124 if (settings.flag(
"SecondHard:SingleW")) {
3125 sigmaPtr =
new Sigma1ffbar2W;
3126 container2Ptrs.push_back(
new ProcessContainer(sigmaPtr) );
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) );
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) );
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) );
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) );
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) );
3171 for ( ProcessContainer * cont : container2Ptrs )
3172 cont->initInfoPtr(*infoPtr);
3183 void SetupContainers::setupIdVecs( Settings& settings) {
3187 if (settings.mode(
"SUSY:idA") != 0) {
3188 idVecA.push_back( abs(settings.mode(
"SUSY:idA")) );
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]) );
3194 nVecA = idVecA.size();
3198 if (settings.mode(
"SUSY:idB") != 0) {
3199 idVecB.push_back( abs(settings.mode(
"SUSY:idB")) );
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]) );
3205 nVecB = idVecB.size();
3214 bool SetupContainers::allowIdVals(
int idCheck1,
int idCheck2) {
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);
3223 if (idChk1 == 0) swap(idChk1, idChk2);
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;
3232 for (
int i = 0; i < nVecA; ++i)
3233 if (idChk1 == idVecA[i] || idChk2 == idVecA[i])
return true;
3239 for (
int i = 0; i < nVecB; ++i)
3240 if (idChk1 == idVecB[i] || idChk2 == idVecB[i])
return true;
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;