00001 #include "GeneratorInterface/Pythia8Interface/plugins/EmissionVetoHook1.h"
00002
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 void EmissionVetoHook1::fatalEmissionVeto(string message) {
00005 throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
00006 << "EmissionVeto: " << message << endl;
00007 }
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 double EmissionVetoHook1::pTpythia(const Pythia8::Event &e, int RadAfterBranch,
00018 int EmtAfterBranch, int RecAfterBranch, bool FSR) {
00019
00020
00021 Pythia8::Vec4 radVec = e[RadAfterBranch].p();
00022 Pythia8::Vec4 emtVec = e[EmtAfterBranch].p();
00023 Pythia8::Vec4 recVec = e[RecAfterBranch].p();
00024 int radID = e[RadAfterBranch].id();
00025
00026
00027 double sign = (FSR) ? 1. : -1.;
00028 Pythia8::Vec4 Q(radVec + sign * emtVec);
00029 double Qsq = sign * Q.m2Calc();
00030
00031
00032 double m2Rad = (abs(radID) >= 4 && abs(radID) < 7) ?
00033 Pythia8::pow2(particleDataPtr->m0(radID)) : 0.;
00034
00035
00036 double z, pTnow;
00037 if (FSR) {
00038
00039 Pythia8::Vec4 sum = radVec + recVec + emtVec;
00040 double m2Dip = sum.m2Calc();
00041 double x1 = 2. * (sum * radVec) / m2Dip;
00042 double x3 = 2. * (sum * emtVec) / m2Dip;
00043 z = x1 / (x1 + x3);
00044 pTnow = z * (1. - z);
00045
00046 } else {
00047
00048 Pythia8::Vec4 qBR(radVec - emtVec + recVec);
00049 Pythia8::Vec4 qAR(radVec + recVec);
00050 z = qBR.m2Calc() / qAR.m2Calc();
00051 pTnow = (1. - z);
00052 }
00053
00054
00055 pTnow *= (Qsq - sign * m2Rad);
00056
00057
00058 if (pTnow < 0.) {
00059 cout << "Warning: pTpythia was negative" << endl;
00060 return -1.;
00061 }
00062
00063 #ifdef DBGOUTPUT
00064 cout << "pTpythia: rad = " << RadAfterBranch << ", emt = "
00065 << EmtAfterBranch << ", rec = " << RecAfterBranch
00066 << ", pTnow = " << sqrt(pTnow) << endl;
00067 #endif
00068
00069
00070 return sqrt(pTnow);
00071 }
00072
00073
00074 double EmissionVetoHook1::pTpowheg(const Pythia8::Event &e, int i, int j, bool FSR) {
00075
00076
00077 double pTnow = 0.;
00078 if (FSR) {
00079
00080
00081
00082 int iInA = partonSystemsPtr->getInA(0);
00083 int iInB = partonSystemsPtr->getInB(0);
00084 double betaZ = - ( e[iInA].pz() + e[iInB].pz() ) /
00085 ( e[iInA].e() + e[iInB].e() );
00086 Pythia8::Vec4 iVecBst(e[i].p()), jVecBst(e[j].p());
00087 iVecBst.bst(0., 0., betaZ);
00088 jVecBst.bst(0., 0., betaZ);
00089 pTnow = sqrt( (iVecBst + jVecBst).m2Calc() *
00090 iVecBst.e() * jVecBst.e() /
00091 Pythia8::pow2(iVecBst.e() + jVecBst.e()) );
00092
00093 } else {
00094
00095 pTnow = e[j].pT();
00096 }
00097
00098
00099 if (pTnow < 0.) {
00100 cout << "Warning: pTpowheg was negative" << endl;
00101 return -1.;
00102 }
00103
00104 #ifdef DBGOUTPUT
00105 cout << "pTpowheg: i = " << i << ", j = " << j
00106 << ", pTnow = " << pTnow << endl;
00107 #endif
00108
00109 return pTnow;
00110 }
00111
00112
00113
00114
00115
00116
00117 double EmissionVetoHook1::pTcalc(const Pythia8::Event &e, int i, int j, int k, int r, int xSRin) {
00118
00119
00120 double pTemt = -1., pTnow;
00121 int xSR1 = (xSRin == -1) ? 0 : xSRin;
00122 int xSR2 = (xSRin == -1) ? 2 : xSRin + 1;
00123 for (int xSR = xSR1; xSR < xSR2; xSR++) {
00124
00125 bool FSR = (xSR == 0) ? false : true;
00126
00127
00128
00129 if ((pTdefMode == 0 || pTdefMode == 1) && i > 0 && j > 0) {
00130 pTemt = pTpowheg(e, i, j, (pTdefMode == 0) ? false : FSR);
00131
00132
00133 } else if (!FSR && pTdefMode == 2 && i > 0 && j > 0 && r > 0) {
00134 pTemt = pTpythia(e, i, j, r, FSR);
00135
00136
00137 } else if (FSR && pTdefMode == 2 && j > 0 && k > 0 && r > 0) {
00138 pTemt = pTpythia(e, k, j, r, FSR);
00139
00140
00141 } else {
00142
00143
00144
00145
00146 int iInA = partonSystemsPtr->getInA(0);
00147 int iInB = partonSystemsPtr->getInB(0);
00148 while (e[iInA].mother1() != 1) { iInA = e[iInA].mother1(); }
00149 while (e[iInB].mother1() != 2) { iInB = e[iInB].mother1(); }
00150
00151
00152 int jNow = (j > 0) ? j : 0;
00153 int jMax = (j > 0) ? j + 1 : e.size();
00154 for (; jNow < jMax; jNow++) {
00155
00156
00157 if (!e[jNow].isFinal()) continue;
00158 if (e[jNow].colType() == 0 && e[jNow].id() != 22) continue;
00159
00160
00161 if (pTdefMode == 0 || pTdefMode == 1) {
00162
00163
00164 if (!FSR) {
00165 pTnow = pTpowheg(e, iInA, jNow, (pTdefMode == 0) ? false : FSR);
00166 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
00167
00168
00169
00170
00171 } else {
00172
00173 int outSize = partonSystemsPtr->sizeOut(0);
00174 for (int iMem = 0; iMem < outSize; iMem++) {
00175 int iNow = partonSystemsPtr->getOut(0, iMem);
00176
00177
00178 if (iNow == jNow || e[iNow].colType() == 0) continue;
00179 if (jNow == e[iNow].daughter1()
00180 && jNow == e[iNow].daughter2()) continue;
00181
00182 pTnow = pTpowheg(e, iNow, jNow, (pTdefMode == 0)
00183 ? false : FSR);
00184 if (pTnow > 0.) pTemt = (pTemt < 0)
00185 ? pTnow : min(pTemt, pTnow);
00186 }
00187
00188 }
00189
00190
00191 } else if (pTdefMode == 2) {
00192
00193
00194 if (!FSR) {
00195 pTnow = pTpythia(e, iInA, jNow, iInB, FSR);
00196 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
00197 pTnow = pTpythia(e, iInB, jNow, iInA, FSR);
00198 if (pTnow > 0.) pTemt = (pTemt < 0) ? pTnow : min(pTemt, pTnow);
00199
00200
00201
00202 } else {
00203 for (int kNow = 0; kNow < e.size(); kNow++) {
00204 if (kNow == jNow || !e[kNow].isFinal() ||
00205 e[kNow].colType() == 0) continue;
00206
00207
00208
00209 pTnow = pTpythia(e, kNow, jNow, iInA, FSR);
00210 if (pTnow > 0.) pTemt = (pTemt < 0)
00211 ? pTnow : min(pTemt, pTnow);
00212 pTnow = pTpythia(e, kNow, jNow, iInB, FSR);
00213 if (pTnow > 0.) pTemt = (pTemt < 0)
00214 ? pTnow : min(pTemt, pTnow);
00215
00216
00217 for (int rNow = 0; rNow < e.size(); rNow++) {
00218 if (rNow == kNow || rNow == jNow ||
00219 !e[rNow].isFinal() || e[rNow].colType() == 0) continue;
00220 pTnow = pTpythia(e, kNow, jNow, rNow, FSR);
00221 if (pTnow > 0.) pTemt = (pTemt < 0)
00222 ? pTnow : min(pTemt, pTnow);
00223 }
00224
00225 }
00226 }
00227 }
00228 }
00229 }
00230 }
00231
00232 #ifdef DBGOUTPUT
00233 cout << "pTcalc: i = " << i << ", j = " << j << ", k = " << k
00234 << ", r = " << r << ", xSR = " << xSRin
00235 << ", pTemt = " << pTemt << endl;
00236 #endif
00237
00238 return pTemt;
00239 }
00240
00241
00242
00243
00244
00245
00246
00247
00248 bool EmissionVetoHook1::doVetoMPIStep(int nMPI, const Pythia8::Event &e) {
00249
00250 if (nMPI > 1) return false;
00251
00252
00253
00254
00255 int count = 0, inonfinal = 0;
00256 double pT1 = 0., pTsum = 0.;
00257 for (int i = e.size() - 1; i > 0; i--) {
00258 inonfinal = i;
00259 if (e[i].isFinal()) {
00260 count++;
00261 pT1 = e[i].pT();
00262 pTsum += e[i].pT();
00263 } else break;
00264 }
00265
00266 nFinal = nFinalExt;
00267
00268 if (nFinal < 0) {
00269 int first = -1, myid;
00270 int last = -1;
00271 for(int ip = 2; ip < e.size(); ip++) {
00272 myid = e[ip].id();
00273 if(abs(myid) < 6 || abs(myid) == 21) continue;
00274 first = ip;
00275 break;
00276 }
00277 if(first < 0) fatalEmissionVeto(string("signal particles not found"));
00278 for(int ip = first; ip < e.size(); ip++) {
00279 myid = e[ip].id();
00280 if(abs(myid) < 6 || abs(myid) == 21) continue;
00281 last = ip;
00282 }
00283 nFinal = last - inonfinal;
00284 }
00285
00286
00287 if (count != nFinal && count != nFinal + 1)
00288 fatalEmissionVeto(string("Wrong number of final state particles in event"));
00289
00290
00291 bool isEmt = (count == nFinal) ? false : true;
00292 int iEmt = (isEmt) ? e.size() - 1 : -1;
00293
00294
00295 if (!isEmt || pThardMode == 0) {
00296 pThard = infoPtr->QRen();
00297
00298
00299
00300 } else if (pThardMode == 1) {
00301 pThard = pTcalc(e, -1, iEmt, -1, -1, -1);
00302
00303
00304
00305
00306 } else if (pThardMode == 2) {
00307 pThard = pTcalc(e, -1, -1, -1, -1, -1);
00308
00309 }
00310
00311
00312 if (MPIvetoOn) {
00313 pTMPI = (isEmt) ? pTsum / 2. : pT1;
00314 }
00315
00316 if(Verbosity)
00317 cout << "doVetoMPIStep: QFac = " << infoPtr->QFac()
00318 << ", QRen = " << infoPtr->QRen()
00319 << ", pThard = " << pThard << endl << endl;
00320
00321
00322 accepted = false;
00323 nAcceptSeq = 0;
00324
00325
00326 return false;
00327 }
00328
00329
00330
00331
00332
00333 bool EmissionVetoHook1::doVetoISREmission(int, const Pythia8::Event &e, int iSys) {
00334
00335 if (iSys != 0) return false;
00336
00337
00338 if (vetoOn && nAcceptSeq >= vetoCount) return false;
00339
00340
00341 int iRadAft = -1, iEmt = -1, iRecAft = -1;
00342 for (int i = e.size() - 1; i > 0; i--) {
00343 if (iRadAft == -1 && e[i].status() == -41) iRadAft = i;
00344 else if (iEmt == -1 && e[i].status() == 43) iEmt = i;
00345 else if (iRecAft == -1 && e[i].status() == -42) iRecAft = i;
00346 if (iRadAft != -1 && iEmt != -1 && iRecAft != -1) break;
00347 }
00348 if (iRadAft == -1 || iEmt == -1 || iRecAft == -1) {
00349 e.list();
00350 fatalEmissionVeto(string("Couldn't find Pythia ISR emission"));
00351 }
00352
00353
00354
00355
00356 int xSR = (pTemtMode == 0) ? 0 : -1;
00357 int i = (pTemtMode == 0) ? iRadAft : -1;
00358 int j = (pTemtMode != 2) ? iEmt : -1;
00359 int k = -1;
00360 int r = (pTemtMode == 0) ? iRecAft : -1;
00361 double pTemt = pTcalc(e, i, j, k, r, xSR);
00362
00363 #ifdef DBGOUTPUT
00364 cout << "doVetoISREmission: pTemt = " << pTemt << endl << endl;
00365 #endif
00366
00367
00368 if (pTemt > pThard) {
00369 nAcceptSeq = 0;
00370 nISRveto++;
00371 return true;
00372 }
00373
00374
00375 nAcceptSeq++;
00376 accepted = true;
00377 return false;
00378 }
00379
00380
00381
00382
00383
00384 bool EmissionVetoHook1::doVetoFSREmission(int, const Pythia8::Event &e, int iSys, bool) {
00385
00386 if (iSys != 0) return false;
00387
00388
00389 if (vetoOn && nAcceptSeq >= vetoCount) return false;
00390
00391
00392 int iRecAft = e.size() - 1;
00393 int iEmt = e.size() - 2;
00394 int iRadAft = e.size() - 3;
00395 int iRadBef = e[iEmt].mother1();
00396 if ( (e[iRecAft].status() != 52 && e[iRecAft].status() != -53) ||
00397 e[iEmt].status() != 51 || e[iRadAft].status() != 51) {
00398 e.list();
00399 fatalEmissionVeto(string("Couldn't find Pythia FSR emission"));
00400 }
00401
00402
00403
00404
00405
00406 int xSR = (pTemtMode == 0) ? 1 : -1;
00407 int i = (pTemtMode == 0) ? iRadBef : -1;
00408 int k = (pTemtMode == 0) ? iRadAft : -1;
00409 int r = (pTemtMode == 0) ? iRecAft : -1;
00410
00411
00412 double pTemt = -1.;
00413 if (pTemtMode == 0 || pTemtMode == 1) {
00414
00415
00416
00417
00418
00419 int j = iRadAft;
00420 if (emittedMode == 0 || (emittedMode == 2 && rndmPtr->flat() < 0.5)) j++;
00421
00422 for (int jLoop = 0; jLoop < 2; jLoop++) {
00423 if (jLoop == 0) pTemt = pTcalc(e, i, j, k, r, xSR);
00424 else if (jLoop == 1) pTemt = min(pTemt, pTcalc(e, i, j, k, r, xSR));
00425
00426
00427 if (emittedMode != 3) break;
00428 if (k != -1) swap(j, k); else j = iEmt;
00429 }
00430
00431
00432 } else if (pTemtMode == 2) {
00433 pTemt = pTcalc(e, i, -1, k, r, xSR);
00434
00435 }
00436
00437 #ifdef DBGOUTPUT
00438 cout << "doVetoFSREmission: pTemt = " << pTemt << endl << endl;
00439 #endif
00440
00441
00442 if (pTemt > pThard) {
00443 nAcceptSeq = 0;
00444 nFSRveto++;
00445 return true;
00446 }
00447
00448
00449 nAcceptSeq++;
00450 accepted = true;
00451 return false;
00452 }
00453
00454
00455
00456
00457
00458 bool EmissionVetoHook1::doVetoMPIEmission(int, const Pythia8::Event &e) {
00459 if (MPIvetoOn) {
00460 if (e[e.size() - 1].pT() > pTMPI) {
00461 #ifdef DBGOUTPUT
00462 cout << "doVetoMPIEmission: pTnow = " << e[e.size() - 1].pT()
00463 << ", pTMPI = " << pTMPI << endl << endl;
00464 #endif
00465 return true;
00466 }
00467 }
00468 return false;
00469 }
00470