18 #ifndef Pythia8_JetMatching_H
19 #define Pythia8_JetMatching_H
22 #include "Pythia8/Pythia.h"
23 #include "Pythia8Plugins/GeneratorInput.h"
32 HJSlowJet(
int powerIn,
double Rin,
double pTjetMinIn = 0.,
33 double etaMaxIn = 25.,
int selectIn = 1,
int massSetIn = 2,
34 SlowJetHook* sjHookPtrIn = 0,
bool useFJcoreIn =
false,
35 bool useStandardRin =
true) :
36 SlowJet(powerIn, Rin, pTjetMinIn, etaMaxIn, selectIn, massSetIn,
37 sjHookPtrIn, useFJcoreIn, useStandardRin) {}
49 void HJSlowJet::findNext() {
57 for (
int i = 1; i < clSize; ++i) {
58 for (
int j = 0; j < i; ++j) {
59 if (dij[i*(i-1)/2 + j] < dMin) {
62 dMin = dij[i*(i-1)/2 + j];
87 JetMatching() : cellJet(NULL), slowJet(NULL), slowJetHard(NULL),
90 if (cellJet)
delete cellJet;
91 if (slowJet)
delete slowJet;
92 if (slowJetHard)
delete slowJetHard;
93 if (hjSlowJet)
delete hjSlowJet;
97 virtual bool initAfterBeams() = 0;
100 bool canVetoProcessLevel() {
return doMerge; }
101 bool doVetoProcessLevel(
Event& process) {
102 eventProcessOrig = process;
107 bool canVetoPartonLevelEarly() {
return doMerge; }
108 bool doVetoPartonLevelEarly(
const Event& event);
111 int numberVetoStep() {
return 1;}
112 bool canVetoStep() {
return false; }
113 bool doVetoStep(
int,
int,
int,
const Event& ) {
return false; }
118 static const bool MATCHINGDEBUG, MATCHINGCHECK;
121 virtual void sortIncomingProcess(
const Event &)=0;
122 virtual void jetAlgorithmInput(
const Event &,
int)=0;
123 virtual void runJetAlgorithm()=0;
124 virtual bool matchPartonsToJets(
int)=0;
125 virtual int matchPartonsToJetsLight()=0;
126 virtual int matchPartonsToJetsHeavy()=0;
128 enum vetoStatus { NONE, LESS_JETS, MORE_JETS, HARD_JET,
129 UNMATCHED_PARTON, INCLUSIVE_VETO };
130 enum partonTypes { ID_CHARM=4, ID_BOT=5, ID_TOP=6, ID_LEPMIN=11,
131 ID_LEPMAX=16, ID_GLUON=21, ID_PHOTON=22 };
144 double eTjetMin, coneRadius, etaJetMax, etaJetMaxAlgo;
159 Event eventProcessOrig, eventProcess, workEventJet;
162 vector<int> typeIdx[3];
167 vector<Vec4> jetMomenta;
171 double eTseed, eTthreshold;
174 int jetAllow, jetMatch, exclusiveMode;
175 double coneMatchLight, coneRadiusHeavy, coneMatchHeavy;
196 bool initAfterBeams();
201 void sortIncomingProcess(
const Event &);
202 void jetAlgorithmInput(
const Event &,
int);
203 void runJetAlgorithm();
204 bool matchPartonsToJets(
int);
205 int matchPartonsToJetsLight();
206 int matchPartonsToJetsHeavy();
209 void sortTypeIdx(vector < int > &vecIn);
212 static const double GHOSTENERGY, ZEROTHRESHOLD;
229 bool initAfterBeams();
232 bool canVetoProcessLevel() {
return doMerge; }
233 bool doVetoProcessLevel(
Event& process);
236 int numberVetoStep() {
return 1;}
237 bool canVetoStep() {
return doShowerKt; }
238 bool doVetoStep(
int,
int,
int,
const Event& );
245 vector<double> getDJR() {
return DJR;}
246 pair<int,int> nMEpartons() {
return nMEpartonsSave;}
252 Event getWorkEventJet() {
return workEventJetSave; }
253 Event getProcessSubset() {
return processSubsetSave; }
254 bool getExclusive() {
return exclusive; }
255 double getPTfirst() {
return pTfirstSave; }
261 Event processSubsetSave;
262 Event workEventJetSave;
270 void sortIncomingProcess(
const Event &);
271 void jetAlgorithmInput(
const Event &,
int);
272 void runJetAlgorithm();
273 bool matchPartonsToJets(
int);
274 int matchPartonsToJetsLight();
275 int matchPartonsToJetsHeavy();
276 int matchPartonsToJetsOther();
277 bool doShowerKtVeto(
double pTfirst);
280 void clearDJR() { DJR.resize(0);}
281 void setDJR(
const Event& event);
283 void clear_nMEpartons() { nMEpartonsSave.first = nMEpartonsSave.second =-1;}
284 void set_nMEpartons(
const int nOrig,
const int nMatch) {
286 nMEpartonsSave.first = nOrig;
287 nMEpartonsSave.second = nMatch;
291 vector<int> origTypeIdx[3];
293 double qCut, qCutSq, clFact;
296 double qCutME, qCutMESq;
302 pair<int,int> nMEpartonsSave;
319 const bool JetMatching::MATCHINGDEBUG =
false;
320 const bool JetMatching::MATCHINGCHECK =
false;
326 inline bool JetMatching::doVetoPartonLevelEarly(
const Event& event) {
338 sortIncomingProcess(event);
342 if ( doShowerKt )
return false;
347 cout << endl <<
"-------- Begin Madgraph Debug --------" << endl;
349 cout << endl <<
"Original incoming process:";
350 eventProcessOrig.list();
352 cout << endl <<
"Final-state incoming process:";
355 for (
size_t i = 0; i < typeIdx[0].size(); i++)
356 cout << ((i == 0) ?
"Light jets: " :
", ") << setw(3) << typeIdx[0][i];
357 if( typeIdx[0].size()== 0 )
358 cout <<
"Light jets: None";
360 for (
size_t i = 0; i < typeIdx[1].size(); i++)
361 cout << ((i == 0) ?
"\nHeavy jets: " :
", ") << setw(3) << typeIdx[1][i];
362 for (
size_t i = 0; i < typeIdx[2].size(); i++)
363 cout << ((i == 0) ?
"\nOther: " :
", ") << setw(3) << typeIdx[2][i];
365 cout << endl << endl <<
"Event:";
368 cout << endl <<
"Work event:";
373 int iTypeEnd = (typeIdx[2].empty()) ? 2 : 3;
374 for (
int iType = 0; iType < iTypeEnd; iType++) {
378 jetAlgorithmInput(event, iType);
383 cout << endl <<
"Jet algorithm event (iType = " << iType <<
"):";
392 if (matchPartonsToJets(iType) ==
true) {
395 cout << endl <<
"Event vetoed" << endl
396 <<
"---------- End MLM Debug ----------" << endl;
404 cout << endl <<
"Event accepted" << endl
405 <<
"---------- End MLM Debug ----------" << endl;
426 const double JetMatchingAlpgen::GHOSTENERGY = 1e-15;
429 const double JetMatchingAlpgen::ZEROTHRESHOLD = 1e-10;
437 inline void JetMatchingAlpgen::sortTypeIdx(vector < int > &vecIn) {
438 for (
size_t i = 0; i < vecIn.size(); i++) {
440 double vMax = (jetAlgorithm == 1) ?
441 eventProcess[vecIn[i]].eT() :
442 eventProcess[vecIn[i]].pT();
443 for (
size_t j = i + 1; j < vecIn.size(); j++) {
444 double vNow = (jetAlgorithm == 1)
445 ? eventProcess[vecIn[j]].eT() : eventProcess[vecIn[j]].pT();
451 if (jMax != i) swap(vecIn[i], vecIn[jMax]);
460 inline bool JetMatchingAlpgen::initAfterBeams() {
463 doMerge = settingsPtr->flag(
"JetMatching:merge");
464 jetAlgorithm = settingsPtr->mode(
"JetMatching:jetAlgorithm");
465 nJet = settingsPtr->mode(
"JetMatching:nJet");
466 nJetMax = settingsPtr->mode(
"JetMatching:nJetMax");
467 eTjetMin = settingsPtr->parm(
"JetMatching:eTjetMin");
468 coneRadius = settingsPtr->parm(
"JetMatching:coneRadius");
469 etaJetMax = settingsPtr->parm(
"JetMatching:etaJetMax");
470 doShowerKt = settingsPtr->flag(
"JetMatching:doShowerKt");
473 etaJetMaxAlgo = etaJetMax + coneRadius;
476 nEta = settingsPtr->mode(
"JetMatching:nEta");
477 nPhi = settingsPtr->mode(
"JetMatching:nPhi");
478 eTseed = settingsPtr->parm(
"JetMatching:eTseed");
479 eTthreshold = settingsPtr->parm(
"JetMatching:eTthreshold");
482 slowJetPower = settingsPtr->mode(
"JetMatching:slowJetPower");
483 coneMatchLight = settingsPtr->parm(
"JetMatching:coneMatchLight");
484 coneRadiusHeavy = settingsPtr->parm(
"JetMatching:coneRadiusHeavy");
485 if (coneRadiusHeavy < 0.) coneRadiusHeavy = coneRadius;
486 coneMatchHeavy = settingsPtr->parm(
"JetMatching:coneMatchHeavy");
489 jetAllow = settingsPtr->mode(
"JetMatching:jetAllow");
490 jetMatch = settingsPtr->mode(
"JetMatching:jetMatch");
491 exclusiveMode = settingsPtr->mode(
"JetMatching:exclusive");
494 if (!doMerge)
return true;
497 if (exclusiveMode == 2) {
500 if (nJet < 0 || nJetMax < 0) {
501 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:init: "
502 "missing jet multiplicity information; running in exclusive mode");
507 exclusive = (nJet == nJetMax) ?
false :
true;
512 exclusive = (exclusiveMode == 0) ?
false :
true;
516 if (jetAlgorithm == 1) {
521 int nSel = 2, smear = 0;
522 double resolution = 0.5, upperCut = 2.;
523 cellJet =
new CellJet(etaJetMaxAlgo, nEta, nPhi, nSel,
524 smear, resolution, upperCut, eTthreshold);
527 }
else if (jetAlgorithm == 2) {
528 slowJet =
new SlowJet(slowJetPower, coneRadius, eTjetMin, etaJetMaxAlgo);
532 if (jetAlgorithm == 1 && jetMatch == 2) {
533 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:init: "
534 "jetMatch = 2 only valid with SlowJet algorithm. "
535 "Reverting to jetMatch = 1");
540 eventProcessOrig.init(
"(eventProcessOrig)", particleDataPtr);
541 eventProcess.init(
"(eventProcess)", particleDataPtr);
542 workEventJet.init(
"(workEventJet)", particleDataPtr);
545 string jetStr = (jetAlgorithm == 1) ?
"CellJet" :
546 (slowJetPower == -1) ?
"anti-kT" :
547 (slowJetPower == 0) ?
"C/A" :
548 (slowJetPower == 1) ?
"kT" :
"unknown";
549 string modeStr = (exclusive) ?
"exclusive" :
"inclusive";
550 stringstream nJetStr, nJetMaxStr;
551 if (nJet >= 0) nJetStr << nJet;
else nJetStr <<
"unknown";
552 if (nJetMax >= 0) nJetMaxStr << nJetMax;
else nJetMaxStr <<
"unknown";
554 <<
" *------- MLM matching parameters -------*" << endl
555 <<
" | nJet | " << setw(14)
556 << nJetStr.str() <<
" |" << endl
557 <<
" | nJetMax | " << setw(14)
558 << nJetMaxStr.str() <<
" |" << endl
559 <<
" | Jet algorithm | " << setw(14)
560 << jetStr <<
" |" << endl
561 <<
" | eTjetMin | " << setw(14)
562 << eTjetMin <<
" |" << endl
563 <<
" | coneRadius | " << setw(14)
564 << coneRadius <<
" |" << endl
565 <<
" | etaJetMax | " << setw(14)
566 << etaJetMax <<
" |" << endl
567 <<
" | jetAllow | " << setw(14)
568 << jetAllow <<
" |" << endl
569 <<
" | jetMatch | " << setw(14)
570 << jetMatch <<
" |" << endl
571 <<
" | coneMatchLight | " << setw(14)
572 << coneMatchLight <<
" |" << endl
573 <<
" | coneRadiusHeavy | " << setw(14)
574 << coneRadiusHeavy <<
" |" << endl
575 <<
" | coneMatchHeavy | " << setw(14)
576 << coneMatchHeavy <<
" |" << endl
577 <<
" | Mode | " << setw(14)
578 << modeStr <<
" |" << endl
579 <<
" *-----------------------------------------*" << endl;
588 inline void JetMatchingAlpgen::sortIncomingProcess(
const Event &event) {
592 omitResonanceDecays(eventProcessOrig,
true);
593 eventProcess = workEvent;
604 for (
int i = 0; i < 3; i++) {
608 for (
int i = 0; i < eventProcess.size(); i++) {
610 if (!eventProcess[i].isFinal())
continue;
614 if (eventProcess[i].
id() == ID_GLUON
615 || (eventProcess[i].idAbs() <= ID_BOT
616 && abs(eventProcess[i].m()) < ZEROTHRESHOLD)) idx = 0;
619 else if (eventProcess[i].idAbs() >= ID_CHARM
620 && eventProcess[i].idAbs() <= ID_TOP) idx = 1;
623 typeIdx[idx].push_back(i);
624 typeSet[idx].insert(eventProcess[i].daughter1());
636 inline void JetMatchingAlpgen::jetAlgorithmInput(
const Event &event,
640 workEventJet = workEvent;
643 for (
int i = 0; i < workEventJet.size(); ++i) {
644 if (!workEventJet[i].isFinal())
continue;
651 int id = workEventJet[i].idAbs();
652 if ( (
id >= ID_LEPMIN &&
id <= ID_LEPMAX) ||
id == ID_TOP
653 ||
id == ID_PHOTON) {
654 workEventJet[i].statusNeg();
660 int idx = workEventJet[i].daughter1();
669 if (typeSet[1].find(idx) != typeSet[1].end() ||
670 typeSet[2].find(idx) != typeSet[2].end()) {
671 workEventJet[i].statusNeg();
678 idx =
event[idx].mother1();
681 }
else if (iType == 1) {
684 if (typeSet[1].find(idx) != typeSet[1].end())
break;
689 workEventJet[i].statusNeg();
694 idx =
event[idx].mother1();
697 }
else if (iType == 2) {
700 if (typeSet[2].find(idx) != typeSet[2].end())
break;
705 workEventJet[i].statusNeg();
710 idx =
event[idx].mother1();
719 for (
int i = 0; i < int(typeIdx[iType].size()); i++) {
721 Vec4 pIn = eventProcess[typeIdx[iType][i]].p();
722 double y = pIn.rap();
723 double phi = pIn.phi();
726 double e = GHOSTENERGY;
727 double e2y = exp(2. * y);
728 double pz = e * (e2y - 1.) / (e2y + 1.);
729 double pt = sqrt(e*e - pz*pz);
730 double px = pt * cos(phi);
731 double py = pt * sin(phi);
732 workEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0, px, py, pz, e);
737 int lastIdx = workEventJet.size() - 1;
738 if (abs(y - workEventJet[lastIdx].y()) > ZEROTHRESHOLD ||
739 abs(phi - workEventJet[lastIdx].phi()) > ZEROTHRESHOLD)
740 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:jetAlgorithmInput: "
741 "ghost particle y/phi mismatch");
752 inline void JetMatchingAlpgen::runJetAlgorithm() {
755 if (jetAlgorithm == 1)
756 cellJet->analyze(workEventJet, eTjetMin, coneRadius, eTseed);
758 slowJet->analyze(workEventJet);
764 int iJet = (jetAlgorithm == 1) ? cellJet->size() - 1:
765 slowJet->sizeJet() - 1;
766 for (
int i = iJet; i > -1; i--) {
767 Vec4 jetMom = (jetAlgorithm == 1) ? cellJet->pMassive(i) :
769 double eta = jetMom.eta();
771 if (abs(eta) > etaJetMax) {
772 if (jetAlgorithm == 2) slowJet->removeJet(i);
775 jetMomenta.push_back(jetMom);
779 reverse(jetMomenta.begin(), jetMomenta.end());
786 inline bool JetMatchingAlpgen::matchPartonsToJets(
int iType) {
790 if (iType == 0)
return (matchPartonsToJetsLight() > 0);
791 else if (iType == 1)
return (matchPartonsToJetsHeavy() > 0);
792 else if (iType == 2)
return false;
809 inline int JetMatchingAlpgen::matchPartonsToJetsLight() {
812 if (jetMomenta.size() < typeIdx[0].size())
return LESS_JETS;
814 if (exclusive && jetMomenta.size() > typeIdx[0].size())
return MORE_JETS;
817 sortTypeIdx(typeIdx[0]);
820 int nParton = typeIdx[0].size();
823 vector < bool > jetAssigned;
824 jetAssigned.assign(jetMomenta.size(),
false);
830 for (
int i = 0; i < nParton; i++) {
831 Vec4 p1 = eventProcess[typeIdx[0][i]].p();
838 for (
int j = 0; j < int(jetMomenta.size()); j++) {
839 if (jetAssigned[j])
continue;
842 double dR = (jetAlgorithm == 1)
843 ? REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
844 if (jMin < 0 || dR < dRmin) {
851 if (jMin >= 0 && dRmin < coneRadius * coneMatchLight) {
856 if (jMin >= nParton)
return HARD_JET;
859 jetAssigned[jMin] =
true;
862 }
else return UNMATCHED_PARTON;
870 for (
int i = workEventJet.size() - nParton;
871 i < workEventJet.size(); i++) {
872 int jMin = slowJet->jetAssignment(i);
878 if (jMin >= nParton)
return HARD_JET;
879 if (jMin < 0 || jetAssigned[jMin])
return UNMATCHED_PARTON;
882 jetAssigned[jMin] =
true;
890 eTpTlightMin = (jetAlgorithm == 1) ? jetMomenta[nParton - 1].eT()
891 : jetMomenta[nParton - 1].pT();
909 inline int JetMatchingAlpgen::matchPartonsToJetsHeavy() {
912 if (jetMomenta.empty())
return NONE;
915 int nParton = typeIdx[1].size();
918 set < int > removeJets;
924 for (
int i = 0; i < nParton; i++) {
925 Vec4 p1 = eventProcess[typeIdx[1][i]].p();
928 for (
int j = 0; j < int(jetMomenta.size()); j++) {
929 double dR = (jetAlgorithm == 1) ?
930 REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
931 if (dR < coneRadiusHeavy * coneMatchHeavy)
932 removeJets.insert(j);
942 for (
int i = workEventJet.size() - nParton;
943 i < workEventJet.size(); i++) {
944 int jMin = slowJet->jetAssignment(i);
945 if (jMin >= 0) removeJets.insert(jMin);
951 for (set < int >::reverse_iterator it = removeJets.rbegin();
952 it != removeJets.rend(); it++)
953 jetMomenta.erase(jetMomenta.begin() + *it);
956 if (!jetMomenta.empty()) {
959 if (exclusive)
return MORE_JETS;
962 else if (eTpTlightMin >= 0.)
963 for (
size_t j = 0; j < jetMomenta.size(); j++) {
965 if ( (jetAlgorithm == 1 && jetMomenta[j].eT() > eTpTlightMin) ||
966 (jetAlgorithm == 2 && jetMomenta[j].pT() > eTpTlightMin) )
987 inline bool JetMatchingMadgraph::initAfterBeams() {
991 processSubsetSave.init(
"(eventProcess)", particleDataPtr);
992 workEventJetSave.init(
"(workEventJet)", particleDataPtr);
995 bool setMad = settingsPtr->flag(
"JetMatching:setMad");
999 string parStr = infoPtr->header(
"MGRunCard");
1000 if (!parStr.empty()) {
1007 if ( par.haveParam(
"xqcut") && par.haveParam(
"maxjetflavor")
1008 && par.haveParam(
"alpsfact") && par.haveParam(
"ickkw") ) {
1009 settingsPtr->flag(
"JetMatching:merge", par.getParam(
"ickkw"));
1010 settingsPtr->parm(
"JetMatching:qCut", par.getParam(
"xqcut"));
1011 settingsPtr->mode(
"JetMatching:nQmatch",
1012 par.getParamAsInt(
"maxjetflavor"));
1013 settingsPtr->parm(
"JetMatching:clFact",
1014 clFact = par.getParam(
"alpsfact"));
1015 if (par.getParamAsInt(
"ickkw") == 0)
1016 infoPtr->errorMsg(
"Error in JetMatchingMadgraph:init: "
1017 "Madgraph file parameters are not set for merging");
1021 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:init: "
1022 "Madgraph merging parameters not found");
1023 if (!par.haveParam(
"xqcut")) infoPtr->errorMsg(
"Warning in "
1024 "JetMatchingMadgraph:init: No xqcut");
1025 if (!par.haveParam(
"ickkw")) infoPtr->errorMsg(
"Warning in "
1026 "JetMatchingMadgraph:init: No ickkw");
1027 if (!par.haveParam(
"maxjetflavor")) infoPtr->errorMsg(
"Warning in "
1028 "JetMatchingMadgraph:init: No maxjetflavor");
1029 if (!par.haveParam(
"alpsfact")) infoPtr->errorMsg(
"Warning in "
1030 "JetMatchingMadgraph:init: No alpsfact");
1035 doFxFx = settingsPtr->flag(
"JetMatching:doFxFx");
1036 nPartonsNow = settingsPtr->mode(
"JetMatching:nPartonsNow");
1037 qCutME = settingsPtr->parm(
"JetMatching:qCutME");
1038 qCutMESq = pow(qCutME,2);
1041 doMerge = settingsPtr->flag(
"JetMatching:merge");
1042 doShowerKt = settingsPtr->flag(
"JetMatching:doShowerKt");
1043 qCut = settingsPtr->parm(
"JetMatching:qCut");
1044 nQmatch = settingsPtr->mode(
"JetMatching:nQmatch");
1045 clFact = settingsPtr->parm(
"JetMatching:clFact");
1048 jetAlgorithm = settingsPtr->mode(
"JetMatching:jetAlgorithm");
1049 nJetMax = settingsPtr->mode(
"JetMatching:nJetMax");
1050 eTjetMin = settingsPtr->parm(
"JetMatching:eTjetMin");
1051 coneRadius = settingsPtr->parm(
"JetMatching:coneRadius");
1052 etaJetMax = settingsPtr->parm(
"JetMatching:etaJetMax");
1053 slowJetPower = settingsPtr->mode(
"JetMatching:slowJetPower");
1056 jetAllow = settingsPtr->mode(
"JetMatching:jetAllow");
1057 exclusiveMode = settingsPtr->mode(
"JetMatching:exclusive");
1058 qCutSq = pow(qCut,2);
1059 etaJetMaxAlgo = etaJetMax;
1062 performVeto = settingsPtr->flag(
"JetMatching:doVeto");
1065 if (!doMerge)
return true;
1068 if (exclusiveMode == 2) {
1072 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:init: "
1073 "missing jet multiplicity information; running in exclusive mode");
1083 slowJet =
new SlowJet(slowJetPower, coneRadius, eTjetMin,
1084 etaJetMaxAlgo, 2, 2, NULL,
false);
1089 slowJetHard =
new SlowJet(slowJetPower, coneRadius, qCutME,
1090 etaJetMaxAlgo, 2, 2, NULL,
false);
1093 slowJetDJR =
new SlowJet(slowJetPower, coneRadius, qCutME,
1094 etaJetMaxAlgo, 2, 2, NULL,
false);
1097 hjSlowJet =
new HJSlowJet(slowJetPower, coneRadius, 0.0,
1098 100.0, 1, 2, NULL,
false,
true);
1101 eventProcessOrig.init(
"(eventProcessOrig)", particleDataPtr);
1102 eventProcess.init(
"(eventProcess)", particleDataPtr);
1103 workEventJet.init(
"(workEventJet)", particleDataPtr);
1106 string jetStr = (jetAlgorithm == 1) ?
"CellJet" :
1107 (slowJetPower == -1) ?
"anti-kT" :
1108 (slowJetPower == 0) ?
"C/A" :
1109 (slowJetPower == 1) ?
"kT" :
"unknown";
1110 string modeStr = (exclusiveMode) ?
"exclusive" :
"inclusive";
1112 <<
" *----- Madgraph matching parameters -----*" << endl
1113 <<
" | qCut | " << setw(14)
1114 << qCut <<
" |" << endl
1115 <<
" | nQmatch | " << setw(14)
1116 << nQmatch <<
" |" << endl
1117 <<
" | clFact | " << setw(14)
1118 << clFact <<
" |" << endl
1119 <<
" | Jet algorithm | " << setw(14)
1120 << jetStr <<
" |" << endl
1121 <<
" | eTjetMin | " << setw(14)
1122 << eTjetMin <<
" |" << endl
1123 <<
" | etaJetMax | " << setw(14)
1124 << etaJetMax <<
" |" << endl
1125 <<
" | jetAllow | " << setw(14)
1126 << jetAllow <<
" |" << endl
1127 <<
" | Mode | " << setw(14)
1128 << modeStr <<
" |" << endl
1129 <<
" *-----------------------------------------*" << endl;
1138 inline bool JetMatchingMadgraph::doVetoProcessLevel(
Event& process) {
1140 eventProcessOrig = process;
1147 sortIncomingProcess(process);
1150 if ( !doFxFx &&
int(typeIdx[0].size()) > nJetMax )
1152 if ( doFxFx && npNLO() < nJetMax &&
int(typeIdx[0].size()) > nJetMax )
1162 inline bool JetMatchingMadgraph::doVetoStep(
int iPos,
int nISR,
int nFSR,
1163 const Event& event) {
1166 if ( !doShowerKt )
return false;
1169 if ( nISR + nFSR > 1 )
return false;
1172 if (iPos == 5)
return false;
1176 sortIncomingProcess(event);
1179 double pTfirst = 0.;
1182 vector<int> weakBosons;
1183 for (
int i = 0; i <
event.size(); i++) {
1184 if ( event[i].
id() == 22
1185 &&
event[i].id() == 23
1186 &&
event[i].idAbs() == 24)
1187 weakBosons.push_back(i);
1190 for (
int i = workEvent.size()-1; i > 0; --i) {
1191 if ( workEvent[i].isFinal() && workEvent[i].colType() != 0
1192 && (workEvent[i].statusAbs() == 43 || workEvent[i].statusAbs() == 51)) {
1196 bool QCDemission =
true;
1199 int iPosOld = workEvent[i].daughter1();
1200 for (
int j = 0; i < int(weakBosons.size()); ++i)
1201 if ( event[iPosOld].isAncestor(j)) {
1202 QCDemission =
false;
1207 pTfirst = workEvent[i].pT();
1214 pTfirstSave = pTfirst;
1216 if (!performVeto)
return false;
1219 if ( doShowerKtVeto(pTfirst) )
return true;
1228 inline bool JetMatchingMadgraph::doShowerKtVeto(
double pTfirst) {
1231 if ( !doShowerKt )
return false;
1234 bool doVeto =
false;
1238 int nParton = typeIdx[0].size();
1239 double pTminME=1e10;
1240 for (
int i = 0; i < nParton; ++i)
1241 pTminME = min(pTminME,eventProcess[typeIdx[0][i]].pT());
1244 if ( nParton > 0 && pow(pTminME,2) < qCutSq ) doVeto =
true;
1248 if ( exclusive && pow(pTfirst,2) > qCutSq ) {
1252 }
else if ( !exclusive && nParton > 0 && pTfirst > pTminME ) {
1265 inline void JetMatchingMadgraph::setDJR(
const Event& event) {
1269 vector<double> result;
1272 if (!slowJetDJR->setup(event) ) {
1273 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:setDJR"
1274 ": the SlowJet algorithm failed on setup");
1279 while ( slowJetDJR->sizeAll() - slowJetDJR->sizeJet() > 0 ) {
1281 result.push_back(sqrt(slowJetDJR->dNext()));
1283 slowJetDJR->doStep();
1287 for (
int i=
int(result.size())-1; i >= 0; --i)
1288 DJR.push_back(result[i]);
1297 inline int JetMatchingMadgraph::npNLO(){
1298 string npIn = infoPtr->getEventAttribute(
"npNLO",
true);
1299 int np = (npIn !=
"") ? atoi((
char*)npIn.c_str()) : -1;
1309 inline void JetMatchingMadgraph::sortIncomingProcess(
const Event &event) {
1313 omitResonanceDecays(eventProcessOrig,
true);
1321 eventProcess.clear();
1322 workEventJet.clear();
1323 for(
int i=0; i < workEvent.size(); ++i) {
1326 int id = workEvent[i].idAbs();
1327 if ((
id >= ID_LEPMIN &&
id <= ID_LEPMAX) ||
id == ID_TOP
1328 ||
id == ID_PHOTON ||
id == 23 ||
id == 24 ||
id == 25) {
1329 eventProcess.append(workEvent[i]);
1331 workEventJet.append(workEvent[i]);
1336 if (!slowJetHard->setup(workEventJet) ) {
1337 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:sortIncomingProcess"
1338 ": the SlowJet algorithm failed on setup");
1343 double localQcutSq = qCutMESq;
1345 while ( slowJetHard->sizeAll() - slowJetHard->sizeJet() > 0 ) {
1347 if( slowJetHard->dNext() > localQcutSq )
break;
1349 if( slowJetHard->sizeAll()-slowJetHard->sizeJet() <= npNLO())
break;
1350 slowJetHard->doStep();
1356 int nJets = slowJetHard->sizeJet();
1357 int nClus = slowJetHard->sizeAll();
1359 for (
int i = nJets; i < nClus; ++i) {
1361 if (i < nClus-nJets) parts = slowJetHard->clusConstituents(i);
1362 else parts = slowJetHard->constituents(nClus-nJets-i);
1363 int flavour = ID_GLUON;
1364 for(
int j=0; j < int(parts.size()); ++j)
1365 if (workEventJet[parts[j]].
id() == ID_BOT)
1367 eventProcess.append( flavour, 98,
1368 workEventJet[parts.back()].mother1(),
1369 workEventJet[parts.back()].mother2(),
1370 workEventJet[parts.back()].daughter1(),
1371 workEventJet[parts.back()].daughter2(),
1372 0, 0, slowJetHard->p(i).px(), slowJetHard->p(i).py(),
1373 slowJetHard->p(i).pz(), slowJetHard->p(i).e() );
1378 workEventJet.clear();
1383 eventProcess = workEvent;
1395 for (
int i = 0; i < 3; i++) {
1398 origTypeIdx[i].clear();
1400 for (
int i = 0; i < eventProcess.size(); i++) {
1402 if (!eventProcess[i].isFinal())
continue;
1407 if (eventProcess[i].isGluon()
1408 || (eventProcess[i].idAbs() <= nQmatch) ) {
1412 idx = ( eventProcess[i].scale() < 1.999 * sqrt(infoPtr->eA()
1413 * infoPtr->eB()) ) ? 0 : 2;
1417 else if (eventProcess[i].idAbs() > nQmatch
1418 && eventProcess[i].idAbs() <= ID_TOP) {
1422 }
else if (eventProcess[i].colType() != 0
1423 && eventProcess[i].idAbs() > ID_TOP) {
1427 if( idx < 0 )
continue;
1429 typeIdx[idx].push_back(i);
1430 typeSet[idx].insert(eventProcess[i].daughter1());
1431 origTypeIdx[orig_idx].push_back(i);
1435 if (exclusiveMode == 2) {
1438 int nParton = origTypeIdx[0].size();
1439 exclusive = (nParton == nJetMax) ?
false :
true;
1443 exclusive = (exclusiveMode == 0) ?
false :
true;
1451 int nParton = typeIdx[0].size();
1452 processSubsetSave.clear();
1453 for (
int i = 0; i < nParton; ++i)
1454 processSubsetSave.append( eventProcess[typeIdx[0][i]] );
1462 inline void JetMatchingMadgraph::jetAlgorithmInput(
const Event &event,
1466 workEventJet = workEvent;
1469 for (
int i = 0; i < workEventJet.size(); ++i) {
1470 if (!workEventJet[i].isFinal())
continue;
1473 if (jetAllow == 1) {
1475 if( workEventJet[i].colType() == 0 ) {
1476 workEventJet[i].statusNeg();
1482 int idx = workEventJet[i].daughter1();
1491 if (typeSet[1].find(idx) != typeSet[1].end() ||
1492 typeSet[2].find(idx) != typeSet[2].end()) {
1493 workEventJet[i].statusNeg();
1498 if (idx == 0)
break;
1500 idx =
event[idx].mother1();
1503 }
else if (iType == 1) {
1506 if (typeSet[1].find(idx) != typeSet[1].end())
break;
1511 workEventJet[i].statusNeg();
1516 idx =
event[idx].mother1();
1519 }
else if (iType == 2) {
1522 if (typeSet[2].find(idx) != typeSet[2].end())
break;
1527 workEventJet[i].statusNeg();
1532 idx =
event[idx].mother1();
1545 inline void JetMatchingMadgraph::runJetAlgorithm() {; }
1551 inline bool JetMatchingMadgraph::matchPartonsToJets(
int iType) {
1558 setDJR(workEventJet);
1559 set_nMEpartons(origTypeIdx[0].size(), typeIdx[0].size());
1561 return (matchPartonsToJetsLight() > 0);
1562 }
else if (iType == 1) {
1563 return (matchPartonsToJetsHeavy() > 0);
1565 return (matchPartonsToJetsOther() > 0);
1583 inline int JetMatchingMadgraph::matchPartonsToJetsLight() {
1586 workEventJetSave = workEventJet;
1588 if (!performVeto)
return false;
1591 int nParton = typeIdx[0].size();
1594 if (!slowJet->setup(workEventJet) ) {
1595 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1596 "Light: the SlowJet algorithm failed on setup");
1599 double localQcutSq = qCutSq;
1602 while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1603 if( slowJet->dNext() > localQcutSq )
break;
1604 dOld = slowJet->dNext();
1607 int nJets = slowJet->sizeJet();
1608 int nClus = slowJet->sizeAll();
1611 if (MATCHINGDEBUG) slowJet->list(
true);
1614 int nCLjets = nClus - nJets;
1616 int nRequested = (doFxFx) ? npNLO() : nParton;
1619 if ( nCLjets < nRequested )
return LESS_JETS;
1622 if ( exclusive && !doFxFx ) {
1623 if ( nCLjets > nRequested )
return MORE_JETS;
1629 if ( doFxFx && nRequested < nJetMax && nCLjets > nRequested )
1636 if (!slowJet->setup(workEventJet) ) {
1637 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1638 "Light: the SlowJet algorithm failed on setup");
1645 while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1646 if( slowJet->dNext() > localQcutSq )
break;
1652 while ( slowJet->sizeAll() - slowJet->sizeJet() > nParton )
1659 if ( clFact >= 0. && nParton > 0 ) {
1660 vector<double> partonPt;
1661 for (
int i = 0; i < nParton; ++i)
1662 partonPt.push_back( eventProcess[typeIdx[0][i]].pT2() );
1663 sort( partonPt.begin(), partonPt.end());
1664 localQcutSq = max( qCutSq, partonPt[0]);
1666 nJets = slowJet->sizeJet();
1667 nClus = slowJet->sizeAll();
1670 if ( clFact != 0. ) localQcutSq *= pow2(clFact);
1673 tempEvent.init(
"(tempEvent)", particleDataPtr);
1675 double pTminEstimate = -1.;
1679 for (
int i = nJets; i < nClus; ++i) {
1680 tempEvent.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0, slowJet->p(i).px(),
1681 slowJet->p(i).py(), slowJet->p(i).pz(), slowJet->p(i).e() );
1683 pTminEstimate = max( pTminEstimate, slowJet->pT(i));
1684 if(nPass == nRequested)
break;
1687 int tempSize = tempEvent.size();
1689 vector<bool> jetAssigned;
1690 jetAssigned.assign( tempSize,
false);
1694 vector< vector<bool> > partonMatchesJet;
1695 for (
int i=0; i < nParton; ++i )
1696 partonMatchesJet.push_back( vector<bool>(tempEvent.size(),
false) );
1712 while ( doFxFx && iNow < tempSize ) {
1716 tempEventJet.init(
"(tempEventJet)", particleDataPtr);
1717 for (
int i=0; i < nParton; ++i ) {
1724 tempEventJet.clear();
1725 tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1726 tempEvent[iNow].px(), tempEvent[iNow].py(),
1727 tempEvent[iNow].pz(), tempEvent[iNow].e() );
1729 Vec4 pIn = eventProcess[typeIdx[0][i]].p();
1730 tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1731 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1734 if ( !slowJet->setup(tempEventJet) ) {
1735 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1736 "Light: the SlowJet algorithm failed on setup");
1742 if ( slowJet->iNext() == tempEventJet.size() - 1
1743 && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1744 jetAssigned[iNow] =
true;
1745 partonMatchesJet[i][iNow] =
true;
1751 if ( jetAssigned[iNow] ) nMatched++;
1759 if ( nRequested < nJetMax && nMatched != nRequested )
1760 return UNMATCHED_PARTON;
1761 if ( nRequested == nJetMax && nMatched < nRequested )
1762 return UNMATCHED_PARTON;
1773 while (!doFxFx && iNow < nParton ) {
1775 tempEventJet.init(
"(tempEventJet)", particleDataPtr);
1776 for (
int i = 0; i < tempSize; ++i) {
1777 if (jetAssigned[i])
continue;
1778 Vec4 pIn = tempEvent[i].p();
1780 tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1781 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1784 Vec4 pIn = eventProcess[typeIdx[0][iNow]].p();
1786 tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1787 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1788 if ( !slowJet->setup(tempEventJet) ) {
1789 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1790 "Light: the SlowJet algorithm failed on setup");
1795 if ( slowJet->iNext() == tempEventJet.size() - 1
1796 && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1798 for (
int i = 0; i != tempSize; ++i) {
1799 if (jetAssigned[i])
continue;
1802 if (iKnt == slowJet->jNext() ) jetAssigned[i] =
true;
1805 return UNMATCHED_PARTON;
1813 if (nParton > 0 && pTminEstimate > 0) eTpTlightMin = pTminEstimate;
1814 else eTpTlightMin = -1.;
1817 setDJR(workEventJet);
1833 inline int JetMatchingMadgraph::matchPartonsToJetsHeavy() {
1845 int nParton = typeIdx[1].size();
1847 Event tempEventJet(workEventJet);
1853 for(
int i=0; i<nParton; ++i) {
1854 scaleF = eventProcessOrig[0].e()/workEventJet[typeIdx[1][i]].pT();
1855 tempEventJet[typeIdx[1][i]].rescale5(scaleF);
1858 if (!hjSlowJet->setup(tempEventJet) ) {
1859 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1860 "Heavy: the SlowJet algorithm failed on setup");
1865 while ( hjSlowJet->sizeAll() - hjSlowJet->sizeJet() > 0 ) {
1866 if( hjSlowJet->dNext() > qCutSq )
break;
1867 hjSlowJet->doStep();
1873 for(
int idx=0 ; idx< hjSlowJet->sizeAll(); ++idx) {
1874 if( hjSlowJet->pT(idx) > sqrt(qCutSq) ) nCLjets++;
1878 if (MATCHINGDEBUG) hjSlowJet->list(
true);
1883 int nRequested = nParton;
1886 if ( nCLjets < nRequested ) {
1887 if (MATCHINGDEBUG) cout <<
"veto : hvy LESS_JETS " << endl;
1888 if (MATCHINGDEBUG) cout <<
"nCLjets = " << nCLjets <<
"; nRequest = "
1889 << nRequested << endl;
1895 if ( nCLjets > nRequested ) {
1896 if (MATCHINGDEBUG) cout <<
"veto : excl hvy MORE_JETS " << endl;
1915 inline int JetMatchingMadgraph::matchPartonsToJetsOther() {
1927 int nParton = typeIdx[2].size();
1929 Event tempEventJet(workEventJet);
1935 for(
int i=0; i<nParton; ++i) {
1936 scaleF = eventProcessOrig[0].e()/workEventJet[typeIdx[2][i]].pT();
1937 tempEventJet[typeIdx[2][i]].rescale5(scaleF);
1940 if (!hjSlowJet->setup(tempEventJet) ) {
1941 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1942 "Heavy: the SlowJet algorithm failed on setup");
1947 while ( hjSlowJet->sizeAll() - hjSlowJet->sizeJet() > 0 ) {
1948 if( hjSlowJet->dNext() > qCutSq )
break;
1949 hjSlowJet->doStep();
1955 for(
int idx=0 ; idx< hjSlowJet->sizeAll(); ++idx) {
1956 if( hjSlowJet->pT(idx) > sqrt(qCutSq) ) nCLjets++;
1960 if (MATCHINGDEBUG) hjSlowJet->list(
true);
1965 int nRequested = nParton;
1968 if ( nCLjets < nRequested ) {
1969 if (MATCHINGDEBUG) cout <<
"veto : other LESS_JETS " << endl;
1970 if (MATCHINGDEBUG) cout <<
"nCLjets = " << nCLjets <<
"; nRequest = "
1971 << nRequested << endl;
1977 if ( nCLjets > nRequested ) {
1978 if (MATCHINGDEBUG) cout <<
"veto : excl other MORE_JETS" << endl;
1991 #endif // end Pythia8_JetMatching_H