14 #ifndef Pythia8_PowhegHooks_H
15 #define Pythia8_PowhegHooks_H
18 #include "Pythia8/Pythia.h"
37 bool initAfterBeams() {
38 nFinal = settingsPtr->mode(
"POWHEG:nFinal");
39 vetoMode = settingsPtr->mode(
"POWHEG:veto");
40 vetoCount = settingsPtr->mode(
"POWHEG:vetoCount");
41 pThardMode = settingsPtr->mode(
"POWHEG:pThard");
42 pTemtMode = settingsPtr->mode(
"POWHEG:pTemt");
43 emittedMode = settingsPtr->mode(
"POWHEG:emitted");
44 pTdefMode = settingsPtr->mode(
"POWHEG:pTdef");
45 MPIvetoMode = settingsPtr->mode(
"POWHEG:MPIveto");
46 QEDvetoMode = settingsPtr->mode(
"POWHEG:QEDveto");
58 inline double pTpythia(
const Event &e,
int RadAfterBranch,
59 int EmtAfterBranch,
int RecAfterBranch,
bool FSR) {
62 Vec4 radVec = e[RadAfterBranch].p();
63 Vec4 emtVec = e[EmtAfterBranch].p();
64 Vec4 recVec = e[RecAfterBranch].p();
65 int radID = e[RadAfterBranch].id();
68 double sign = (FSR) ? 1. : -1.;
69 Vec4 Q(radVec + sign * emtVec);
70 double Qsq = sign * Q.m2Calc();
73 double m2Rad = (abs(radID) >= 4 && abs(radID) < 7) ?
74 pow2(particleDataPtr->m0(radID)) : 0.;
80 Vec4 sum = radVec + recVec + emtVec;
81 double m2Dip = sum.m2Calc();
82 double x1 = 2. * (sum * radVec) / m2Dip;
83 double x3 = 2. * (sum * emtVec) / m2Dip;
89 Vec4 qBR(radVec - emtVec + recVec);
90 Vec4 qAR(radVec + recVec);
91 z = qBR.m2Calc() / qAR.m2Calc();
96 pTnow *= (Qsq - sign * m2Rad);
100 cout <<
"Warning: pTpythia was negative" << endl;
109 inline double pTpowheg(
const Event &e,
int i,
int j,
bool FSR) {
117 int iInA = partonSystemsPtr->getInA(0);
118 int iInB = partonSystemsPtr->getInB(0);
119 double betaZ = - ( e[iInA].pz() + e[iInB].pz() ) /
120 ( e[iInA].e() + e[iInB].e() );
121 Vec4 iVecBst(e[i].p()), jVecBst(e[j].p());
122 iVecBst.bst(0., 0., betaZ);
123 jVecBst.bst(0., 0., betaZ);
124 pTnow = sqrt( (iVecBst + jVecBst).m2Calc() *
125 iVecBst.e() * jVecBst.e() /
126 pow2(iVecBst.e() + jVecBst.e()) );
135 cout <<
"Warning: pTpowheg was negative" << endl;
147 inline double pTcalc(
const Event &e,
int i,
int j,
int k,
int r,
int xSRin) {
150 double pTemt = -1., pTnow;
151 int xSR1 = (xSRin == -1) ? 0 : xSRin;
152 int xSR2 = (xSRin == -1) ? 2 : xSRin + 1;
153 for (
int xSR = xSR1; xSR < xSR2; xSR++) {
155 bool FSR = (xSR == 0) ?
false :
true;
159 if ((pTdefMode == 0 || pTdefMode == 1) && i > 0 && j > 0) {
160 pTemt = pTpowheg(e, i, j, (pTdefMode == 0) ?
false : FSR);
163 }
else if (!FSR && pTdefMode == 2 && i > 0 && j > 0 && r > 0) {
164 pTemt = pTpythia(e, i, j, r, FSR);
167 }
else if (FSR && pTdefMode == 2 && j > 0 && k > 0 && r > 0) {
168 pTemt = pTpythia(e, k, j, r, FSR);
176 int iInA = partonSystemsPtr->getInA(0);
177 int iInB = partonSystemsPtr->getInB(0);
178 while (e[iInA].mother1() != 1) { iInA = e[iInA].mother1(); }
179 while (e[iInB].mother1() != 2) { iInB = e[iInB].mother1(); }
182 int jNow = (j > 0) ? j : 0;
183 int jMax = (j > 0) ? j + 1 : e.size();
184 for (; jNow < jMax; jNow++) {
187 if (!e[jNow].isFinal())
continue;
189 if (QEDvetoMode==0 && e[jNow].colType() == 0)
continue;
192 if (pTdefMode == 0 || pTdefMode == 1) {
196 pTnow = pTpowheg(e, iInA, jNow, (pTdefMode == 0) ?
false : FSR);
197 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
204 int outSize = partonSystemsPtr->sizeOut(0);
205 for (
int iMem = 0; iMem < outSize; iMem++) {
206 int iNow = partonSystemsPtr->getOut(0, iMem);
209 if (iNow == jNow )
continue;
211 if (QEDvetoMode==0 && e[iNow].colType() == 0)
continue;
212 if (jNow == e[iNow].daughter1()
213 && jNow == e[iNow].daughter2())
continue;
215 pTnow = pTpowheg(e, iNow, jNow, (pTdefMode == 0)
217 if (pTnow > 0.) pTemt = (pTemt < 0)
218 ? pTnow : min(pTemt, pTnow);
224 }
else if (pTdefMode == 2) {
228 pTnow = pTpythia(e, iInA, jNow, iInB, FSR);
229 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
230 pTnow = pTpythia(e, iInB, jNow, iInA, FSR);
231 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
236 for (
int kNow = 0; kNow < e.size(); kNow++) {
237 if (kNow == jNow || !e[kNow].isFinal())
continue;
238 if (QEDvetoMode==0 && e[kNow].colType() == 0)
continue;
242 pTnow = pTpythia(e, kNow, jNow, iInA, FSR);
243 if (pTnow > 0.) pTemt = (pTemt < 0)
244 ? pTnow : min(pTemt, pTnow);
245 pTnow = pTpythia(e, kNow, jNow, iInB, FSR);
246 if (pTnow > 0.) pTemt = (pTemt < 0)
247 ? pTnow : min(pTemt, pTnow);
250 for (
int rNow = 0; rNow < e.size(); rNow++) {
251 if (rNow == kNow || rNow == jNow ||
252 !e[rNow].isFinal())
continue;
253 if(QEDvetoMode==0 && e[rNow].colType() == 0)
continue;
254 pTnow = pTpythia(e, kNow, jNow, rNow, FSR);
255 if (pTnow > 0.) pTemt = (pTemt < 0)
256 ? pTnow : min(pTemt, pTnow);
281 inline bool canVetoMPIStep() {
return true; }
282 inline int numberVetoMPIStep() {
return 1; }
283 inline bool doVetoMPIStep(
int nMPI,
const Event &e) {
285 if (nMPI > 1)
return false;
291 double pT1 = 0., pTsum = 0.;
292 for (
int i = e.size() - 1; i > 0; i--) {
293 if (e[i].isFinal()) {
300 if (count != nFinal && count != nFinal + 1) {
301 cout <<
"Error: wrong number of final state particles in event" << endl;
305 isEmt = (count == nFinal) ?
false :
true;
306 int iEmt = (isEmt) ? e.size() - 1 : -1;
309 if (!isEmt || pThardMode == 0) {
310 pThard = infoPtr->scalup();
314 }
else if (pThardMode == 1) {
315 pThard = pTcalc(e, -1, iEmt, -1, -1, -1);
320 }
else if (pThardMode == 2) {
321 pThard = pTcalc(e, -1, -1, -1, -1, -1);
325 if (MPIvetoMode == 1) {
326 pTMPI = (isEmt) ? pTsum / 2. : pT1;
331 nAcceptSeq = nISRveto = nFSRveto = 0;
341 inline bool canVetoISREmission() {
return (vetoMode == 0) ?
false :
true; }
342 inline bool doVetoISREmission(
int,
const Event &e,
int iSys) {
344 if (iSys != 0)
return false;
347 if (vetoMode == 1 && nAcceptSeq >= vetoCount)
return false;
350 int iRadAft = -1, iEmt = -1, iRecAft = -1;
351 for (
int i = e.size() - 1; i > 0; i--) {
352 if (iRadAft == -1 && e[i].status() == -41) iRadAft = i;
353 else if (iEmt == -1 && e[i].status() == 43) iEmt = i;
354 else if (iRecAft == -1 && e[i].status() == -42) iRecAft = i;
355 if (iRadAft != -1 && iEmt != -1 && iRecAft != -1)
break;
357 if (iRadAft == -1 || iEmt == -1 || iRecAft == -1) {
359 cout <<
"Error: couldn't find Pythia ISR emission" << endl;
366 int xSR = (pTemtMode == 0) ? 0 : -1;
367 int i = (pTemtMode == 0) ? iRadAft : -1;
368 int j = (pTemtMode != 2) ? iEmt : -1;
370 int r = (pTemtMode == 0) ? iRecAft : -1;
371 double pTemt = pTcalc(e, i, j, k, r, xSR);
375 bool vetoParton = (!isEmt && e[iEmt].colType()==0 && QEDvetoMode==2)
379 if (pTemt > pThard) {
382 nAcceptSeq = vetoCount-1;
400 inline bool canVetoFSREmission() {
return (vetoMode == 0) ?
false :
true; }
401 inline bool doVetoFSREmission(
int,
const Event &e,
int iSys,
bool) {
403 if (iSys != 0)
return false;
406 if (vetoMode == 1 && nAcceptSeq >= vetoCount)
return false;
409 int iRecAft = e.size() - 1;
410 int iEmt = e.size() - 2;
411 int iRadAft = e.size() - 3;
412 int iRadBef = e[iEmt].mother1();
413 if ( (e[iRecAft].status() != 52 && e[iRecAft].status() != -53) ||
414 e[iEmt].status() != 51 || e[iRadAft].status() != 51) {
416 cout <<
"Error: couldn't find Pythia FSR emission" << endl;
424 int xSR = (pTemtMode == 0) ? 1 : -1;
425 int i = (pTemtMode == 0) ? iRadBef : -1;
426 int k = (pTemtMode == 0) ? iRadAft : -1;
427 int r = (pTemtMode == 0) ? iRecAft : -1;
431 if (pTemtMode == 0 || pTemtMode == 1) {
438 if (emittedMode == 0 || (emittedMode == 2 && rndmPtr->flat() < 0.5)) j++;
440 for (
int jLoop = 0; jLoop < 2; jLoop++) {
441 if (jLoop == 0) pTemt = pTcalc(e, i, j, k, r, xSR);
442 else if (jLoop == 1) pTemt = min(pTemt, pTcalc(e, i, j, k, r, xSR));
445 if (emittedMode != 3)
break;
446 if (k != -1) swap(j, k);
else j = iEmt;
450 }
else if (pTemtMode == 2) {
451 pTemt = pTcalc(e, i, -1, k, r, xSR);
457 bool vetoParton = (!isEmt && e[iEmt].colType()==0 && QEDvetoMode==2)
461 if (pTemt > pThard) {
464 nAcceptSeq = vetoCount-1;
482 inline bool canVetoMPIEmission() {
return (MPIvetoMode == 0) ?
false :
true;}
483 inline bool doVetoMPIEmission(
int,
const Event &e) {
484 if (MPIvetoMode == 1) {
485 if (e[e.size() - 1].pT() > pTMPI)
return true;
494 inline int getNISRveto() {
return nISRveto; }
495 inline int getNFSRveto() {
return nFSRveto; }
500 int nFinal, vetoMode, vetoCount, pThardMode, pTemtMode,
501 emittedMode, pTdefMode, MPIvetoMode, QEDvetoMode;
502 double pThard, pTMPI;
503 bool accepted, isEmt;
508 unsigned long int nISRveto, nFSRveto;
516 #endif // end Pythia8_PowhegHooks_H