18 #ifndef Pythia8_JetMatching_H
19 #define Pythia8_JetMatching_H
22 #include "Pythia8/Pythia.h"
23 #include "GeneratorInput.h"
24 using namespace Pythia8;
37 JetMatching() : cellJet(NULL), slowJet(NULL), slowJetHard(NULL) {}
39 if (cellJet)
delete cellJet;
40 if (slowJet)
delete slowJet;
41 if (slowJetHard)
delete slowJetHard;
45 virtual bool initAfterBeams() = 0;
48 bool canVetoProcessLevel() {
return doMerge; }
49 bool doVetoProcessLevel(
Event& process) {
50 eventProcessOrig = process;
55 bool canVetoPartonLevelEarly() {
return doMerge; }
56 bool doVetoPartonLevelEarly(
const Event& event);
59 int numberVetoStep() {
return 1;}
60 bool canVetoStep() {
return false; }
61 bool doVetoStep(
int,
int,
int,
const Event& ) {
return false; }
66 static const bool MATCHINGDEBUG, MATCHINGCHECK;
69 virtual void sortIncomingProcess(
const Event &)=0;
70 virtual void jetAlgorithmInput(
const Event &,
int)=0;
71 virtual void runJetAlgorithm()=0;
72 virtual bool matchPartonsToJets(
int)=0;
73 virtual int matchPartonsToJetsLight()=0;
74 virtual int matchPartonsToJetsHeavy()=0;
76 enum vetoStatus { NONE, LESS_JETS, MORE_JETS, HARD_JET, UNMATCHED_PARTON };
77 enum partonTypes { ID_CHARM=4, ID_BOT=5, ID_TOP=6, ID_LEPMIN=11,
78 ID_LEPMAX=16, ID_GLUON=21, ID_PHOTON=22 };
91 double eTjetMin, coneRadius, etaJetMax, etaJetMaxAlgo;
105 Event eventProcessOrig, eventProcess, workEventJet;
108 vector<int> typeIdx[3];
113 vector<Vec4> jetMomenta;
117 double eTseed, eTthreshold;
120 int jetAllow, jetMatch, exclusiveMode;
121 double coneMatchLight, coneRadiusHeavy, coneMatchHeavy;
142 bool initAfterBeams();
147 void sortIncomingProcess(
const Event &);
148 void jetAlgorithmInput(
const Event &,
int);
149 void runJetAlgorithm();
150 bool matchPartonsToJets(
int);
151 int matchPartonsToJetsLight();
152 int matchPartonsToJetsHeavy();
155 void sortTypeIdx(vector < int > &vecIn);
158 static const double GHOSTENERGY, ZEROTHRESHOLD;
175 bool initAfterBeams();
178 int numberVetoStep() {
return 1;}
179 bool canVetoStep() {
return doShowerKt; }
180 bool doVetoStep(
int,
int,
int,
const Event& );
185 void sortIncomingProcess(
const Event &);
186 void jetAlgorithmInput(
const Event &,
int);
187 void runJetAlgorithm();
188 bool matchPartonsToJets(
int);
189 int matchPartonsToJetsLight();
190 int matchPartonsToJetsHeavy();
191 bool doShowerKtVeto(
double pTfirst);
194 vector<int> origTypeIdx[3];
196 double qCut, qCutSq, clFact;
199 double qCutME, qCutMESq;
212 const bool JetMatching::MATCHINGDEBUG =
false;
213 const bool JetMatching::MATCHINGCHECK =
false;
219 bool JetMatching::doVetoPartonLevelEarly(
const Event& event) {
231 sortIncomingProcess(event);
235 if ( doShowerKt )
return false;
240 cout << endl <<
"-------- Begin Madgraph Debug --------" << endl;
242 cout << endl <<
"Original incoming process:";
243 eventProcessOrig.list();
245 cout << endl <<
"Final-state incoming process:";
248 for (
size_t i = 0; i < typeIdx[0].size(); i++)
249 cout << ((i == 0) ?
"Light jets: " :
", ") << setw(3) << typeIdx[0][i];
250 if( typeIdx[0].size()== 0 )
251 cout <<
"Light jets: None";
253 for (
size_t i = 0; i < typeIdx[1].size(); i++)
254 cout << ((i == 0) ?
"\nHeavy jets: " :
", ") << setw(3) << typeIdx[1][i];
255 for (
size_t i = 0; i < typeIdx[2].size(); i++)
256 cout << ((i == 0) ?
"\nOther: " :
", ") << setw(3) << typeIdx[2][i];
258 cout << endl << endl <<
"Event:";
261 cout << endl <<
"Work event:";
266 int iTypeEnd = (typeIdx[1].empty()) ? 1 : 2;
267 for (
int iType = 0; iType < iTypeEnd; iType++) {
271 jetAlgorithmInput(event, iType);
276 cout << endl <<
"Jet algorithm event (iType = " << iType <<
"):";
285 if (matchPartonsToJets(iType) ==
true) {
288 cout << endl <<
"Event vetoed" << endl
289 <<
"---------- End MLM Debug ----------" << endl;
297 cout << endl <<
"Event accepted" << endl
298 <<
"---------- End MLM Debug ----------" << endl;
319 const double JetMatchingAlpgen::GHOSTENERGY = 1e-15;
322 const double JetMatchingAlpgen::ZEROTHRESHOLD = 1e-10;
330 void JetMatchingAlpgen::sortTypeIdx(vector < int > &vecIn) {
331 for (
size_t i = 0; i < vecIn.size(); i++) {
333 double vMax = (jetAlgorithm == 1) ?
334 eventProcess[vecIn[i]].eT() :
335 eventProcess[vecIn[i]].pT();
336 for (
size_t j = i + 1; j < vecIn.size(); j++) {
337 double vNow = (jetAlgorithm == 1)
338 ? eventProcess[vecIn[j]].eT() : eventProcess[vecIn[j]].pT();
344 if (jMax != i) swap(vecIn[i], vecIn[jMax]);
353 bool JetMatchingAlpgen::initAfterBeams() {
356 doMerge = settingsPtr->flag(
"JetMatching:merge");
357 jetAlgorithm = settingsPtr->mode(
"JetMatching:jetAlgorithm");
358 nJet = settingsPtr->mode(
"JetMatching:nJet");
359 nJetMax = settingsPtr->mode(
"JetMatching:nJetMax");
360 eTjetMin = settingsPtr->parm(
"JetMatching:eTjetMin");
361 coneRadius = settingsPtr->parm(
"JetMatching:coneRadius");
362 etaJetMax = settingsPtr->parm(
"JetMatching:etaJetMax");
363 doShowerKt = settingsPtr->flag(
"JetMatching:doShowerKt");
366 etaJetMaxAlgo = etaJetMax + coneRadius;
369 nEta = settingsPtr->mode(
"JetMatching:nEta");
370 nPhi = settingsPtr->mode(
"JetMatching:nPhi");
371 eTseed = settingsPtr->parm(
"JetMatching:eTseed");
372 eTthreshold = settingsPtr->parm(
"JetMatching:eTthreshold");
375 slowJetPower = settingsPtr->mode(
"JetMatching:slowJetPower");
376 coneMatchLight = settingsPtr->parm(
"JetMatching:coneMatchLight");
377 coneRadiusHeavy = settingsPtr->parm(
"JetMatching:coneRadiusHeavy");
378 if (coneRadiusHeavy < 0.) coneRadiusHeavy = coneRadius;
379 coneMatchHeavy = settingsPtr->parm(
"JetMatching:coneMatchHeavy");
382 jetAllow = settingsPtr->mode(
"JetMatching:jetAllow");
383 jetMatch = settingsPtr->mode(
"JetMatching:jetMatch");
384 exclusiveMode = settingsPtr->mode(
"JetMatching:exclusive");
387 if (!doMerge)
return true;
390 if (exclusiveMode == 2) {
393 if (nJet < 0 || nJetMax < 0) {
394 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:init: "
395 "missing jet multiplicity information; running in exclusive mode");
400 exclusive = (nJet == nJetMax) ?
false :
true;
405 exclusive = (exclusiveMode == 0) ?
false :
true;
409 if (jetAlgorithm == 1) {
414 int nSel = 2, smear = 0;
415 double resolution = 0.5, upperCut = 2.;
416 cellJet =
new CellJet(etaJetMaxAlgo, nEta, nPhi, nSel,
417 smear, resolution, upperCut, eTthreshold);
420 }
else if (jetAlgorithm == 2) {
421 slowJet =
new SlowJet(slowJetPower, coneRadius, eTjetMin, etaJetMaxAlgo);
425 if (jetAlgorithm == 1 && jetMatch == 2) {
426 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:init: "
427 "jetMatch = 2 only valid with SlowJet algorithm. "
428 "Reverting to jetMatch = 1");
433 eventProcessOrig.init(
"(eventProcessOrig)", particleDataPtr);
434 eventProcess.init(
"(eventProcess)", particleDataPtr);
435 workEventJet.init(
"(workEventJet)", particleDataPtr);
438 string jetStr = (jetAlgorithm == 1) ?
"CellJet" :
439 (slowJetPower == -1) ?
"anti-kT" :
440 (slowJetPower == 0) ?
"C/A" :
441 (slowJetPower == 1) ?
"kT" :
"unknown";
442 string modeStr = (exclusive) ?
"exclusive" :
"inclusive";
443 stringstream nJetStr, nJetMaxStr;
444 if (nJet >= 0) nJetStr << nJet;
else nJetStr <<
"unknown";
445 if (nJetMax >= 0) nJetMaxStr << nJetMax;
else nJetMaxStr <<
"unknown";
447 <<
" *------- MLM matching parameters -------*" << endl
448 <<
" | nJet | " << setw(14)
449 << nJetStr.str() <<
" |" << endl
450 <<
" | nJetMax | " << setw(14)
451 << nJetMaxStr.str() <<
" |" << endl
452 <<
" | Jet algorithm | " << setw(14)
453 << jetStr <<
" |" << endl
454 <<
" | eTjetMin | " << setw(14)
455 << eTjetMin <<
" |" << endl
456 <<
" | coneRadius | " << setw(14)
457 << coneRadius <<
" |" << endl
458 <<
" | etaJetMax | " << setw(14)
459 << etaJetMax <<
" |" << endl
460 <<
" | jetAllow | " << setw(14)
461 << jetAllow <<
" |" << endl
462 <<
" | jetMatch | " << setw(14)
463 << jetMatch <<
" |" << endl
464 <<
" | coneMatchLight | " << setw(14)
465 << coneMatchLight <<
" |" << endl
466 <<
" | coneRadiusHeavy | " << setw(14)
467 << coneRadiusHeavy <<
" |" << endl
468 <<
" | coneMatchHeavy | " << setw(14)
469 << coneMatchHeavy <<
" |" << endl
470 <<
" | Mode | " << setw(14)
471 << modeStr <<
" |" << endl
472 <<
" *-----------------------------------------*" << endl;
481 void JetMatchingAlpgen::sortIncomingProcess(
const Event &event) {
485 omitResonanceDecays(eventProcessOrig,
true);
486 eventProcess = workEvent;
497 for (
int i = 0; i < 3; i++) {
501 for (
int i = 0; i < eventProcess.size(); i++) {
503 if (!eventProcess[i].isFinal())
continue;
507 if (eventProcess[i].
id() == ID_GLUON
508 || (eventProcess[i].idAbs() <= ID_BOT
509 && abs(eventProcess[i].m()) < ZEROTHRESHOLD)) idx = 0;
512 else if (eventProcess[i].idAbs() >= ID_CHARM
513 && eventProcess[i].idAbs() <= ID_TOP) idx = 1;
516 typeIdx[idx].push_back(i);
517 typeSet[idx].insert(eventProcess[i].daughter1());
529 void JetMatchingAlpgen::jetAlgorithmInput(
const Event &event,
int iType) {
532 workEventJet = workEvent;
535 for (
int i = 0; i < workEventJet.size(); ++i) {
536 if (!workEventJet[i].isFinal())
continue;
543 int id = workEventJet[i].idAbs();
544 if ( (
id >= ID_LEPMIN &&
id <= ID_LEPMAX) ||
id == ID_TOP
545 ||
id == ID_PHOTON) {
546 workEventJet[i].statusNeg();
552 int idx = workEventJet[i].daughter1();
561 if (typeSet[1].find(idx) != typeSet[1].end() ||
562 typeSet[2].find(idx) != typeSet[2].end()) {
563 workEventJet[i].statusNeg();
570 idx =
event[idx].mother1();
573 }
else if (iType == 1) {
576 if (typeSet[1].find(idx) != typeSet[1].end())
break;
581 workEventJet[i].statusNeg();
586 idx =
event[idx].mother1();
595 for (
int i = 0; i < int(typeIdx[iType].size()); i++) {
597 Vec4 pIn = eventProcess[typeIdx[iType][i]].p();
598 double y = pIn.rap();
599 double phi = pIn.phi();
602 double e = GHOSTENERGY;
603 double e2y = exp(2. * y);
604 double pz = e * (e2y - 1.) / (e2y + 1.);
605 double pt = sqrt(e*e - pz*pz);
606 double px = pt * cos(phi);
607 double py = pt * sin(phi);
608 workEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0, px, py, pz, e);
613 int lastIdx = workEventJet.size() - 1;
614 if (abs(y - workEventJet[lastIdx].y()) > ZEROTHRESHOLD ||
615 abs(phi - workEventJet[lastIdx].phi()) > ZEROTHRESHOLD)
616 infoPtr->errorMsg(
"Warning in JetMatchingAlpgen:jetAlgorithmInput: "
617 "ghost particle y/phi mismatch");
628 void JetMatchingAlpgen::runJetAlgorithm() {
631 if (jetAlgorithm == 1)
632 cellJet->analyze(workEventJet, eTjetMin, coneRadius, eTseed);
634 slowJet->analyze(workEventJet);
640 int iJet = (jetAlgorithm == 1) ? cellJet->size() - 1:
641 slowJet->sizeJet() - 1;
642 for (
int i = iJet; i > -1; i--) {
643 Vec4 jetMom = (jetAlgorithm == 1) ? cellJet->pMassive(i) :
645 double eta = jetMom.eta();
647 if (abs(eta) > etaJetMax) {
648 if (jetAlgorithm == 2) slowJet->removeJet(i);
651 jetMomenta.push_back(jetMom);
655 reverse(jetMomenta.begin(), jetMomenta.end());
662 bool JetMatchingAlpgen::matchPartonsToJets(
int iType) {
666 if (iType == 0)
return (matchPartonsToJetsLight() > 0);
667 else return (matchPartonsToJetsHeavy() > 0);
683 int JetMatchingAlpgen::matchPartonsToJetsLight() {
686 if (jetMomenta.size() < typeIdx[0].size())
return LESS_JETS;
688 if (exclusive && jetMomenta.size() > typeIdx[0].size())
return MORE_JETS;
691 sortTypeIdx(typeIdx[0]);
694 int nParton = typeIdx[0].size();
697 vector < bool > jetAssigned;
698 jetAssigned.assign(jetMomenta.size(),
false);
704 for (
int i = 0; i < nParton; i++) {
705 Vec4 p1 = eventProcess[typeIdx[0][i]].p();
712 for (
int j = 0; j < int(jetMomenta.size()); j++) {
713 if (jetAssigned[j])
continue;
716 double dR = (jetAlgorithm == 1)
717 ? REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
718 if (jMin < 0 || dR < dRmin) {
725 if (jMin >= 0 && dRmin < coneRadius * coneMatchLight) {
730 if (jMin >= nParton)
return HARD_JET;
733 jetAssigned[jMin] =
true;
736 }
else return UNMATCHED_PARTON;
744 for (
int i = workEventJet.size() - nParton;
745 i < workEventJet.size(); i++) {
746 int jMin = slowJet->jetAssignment(i);
752 if (jMin >= nParton)
return HARD_JET;
753 if (jMin < 0 || jetAssigned[jMin])
return UNMATCHED_PARTON;
756 jetAssigned[jMin] =
true;
764 eTpTlightMin = (jetAlgorithm == 1) ? jetMomenta[nParton - 1].eT()
765 : jetMomenta[nParton - 1].pT();
783 int JetMatchingAlpgen::matchPartonsToJetsHeavy() {
786 if (jetMomenta.empty())
return NONE;
789 int nParton = typeIdx[1].size();
792 set < int > removeJets;
798 for (
int i = 0; i < nParton; i++) {
799 Vec4 p1 = eventProcess[typeIdx[1][i]].p();
802 for (
int j = 0; j < int(jetMomenta.size()); j++) {
803 double dR = (jetAlgorithm == 1) ?
804 REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
805 if (dR < coneRadiusHeavy * coneMatchHeavy)
806 removeJets.insert(j);
816 for (
int i = workEventJet.size() - nParton;
817 i < workEventJet.size(); i++) {
818 int jMin = slowJet->jetAssignment(i);
819 if (jMin >= 0) removeJets.insert(jMin);
825 for (set < int >::reverse_iterator it = removeJets.rbegin();
826 it != removeJets.rend(); it++)
827 jetMomenta.erase(jetMomenta.begin() + *it);
830 if (!jetMomenta.empty()) {
833 if (exclusive)
return MORE_JETS;
836 else if (eTpTlightMin >= 0.)
837 for (
size_t j = 0; j < jetMomenta.size(); j++) {
839 if ( (jetAlgorithm == 1 && jetMomenta[j].eT() > eTpTlightMin) ||
840 (jetAlgorithm == 2 && jetMomenta[j].pT() > eTpTlightMin) )
861 bool JetMatchingMadgraph::initAfterBeams() {
864 bool setMad = settingsPtr->flag(
"JetMatching:setMad");
868 string parStr = infoPtr->header(
"MGRunCard");
869 if (!parStr.empty()) {
876 if ( par.haveParam(
"xqcut") && par.haveParam(
"maxjetflavor")
877 && par.haveParam(
"alpsfact") && par.haveParam(
"ickkw") ) {
878 settingsPtr->flag(
"JetMatching:merge", par.getParam(
"ickkw"));
879 settingsPtr->parm(
"JetMatching:qCut", par.getParam(
"xqcut"));
880 settingsPtr->mode(
"JetMatching:nQmatch",
881 par.getParamAsInt(
"maxjetflavor"));
882 settingsPtr->parm(
"JetMatching:clFact",
883 clFact = par.getParam(
"alpsfact"));
884 if (par.getParamAsInt(
"ickkw") == 0)
885 infoPtr->errorMsg(
"Error in JetMatchingMadgraph:init: "
886 "Madgraph file parameters are not set for merging");
890 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:init: "
891 "Madgraph merging parameters not found");
892 if (!par.haveParam(
"xqcut")) infoPtr->errorMsg(
"Warning in "
893 "JetMatchingMadgraph:init: No xqcut");
894 if (!par.haveParam(
"ickkw")) infoPtr->errorMsg(
"Warning in "
895 "JetMatchingMadgraph:init: No ickkw");
896 if (!par.haveParam(
"maxjetflavor")) infoPtr->errorMsg(
"Warning in "
897 "JetMatchingMadgraph:init: No maxjetflavor");
898 if (!par.haveParam(
"alpsfact")) infoPtr->errorMsg(
"Warning in "
899 "JetMatchingMadgraph:init: No alpsfact");
904 doFxFx = settingsPtr->flag(
"JetMatching:doFxFx");
905 nPartonsNow = settingsPtr->mode(
"JetMatching:nPartonsNow");
906 qCutME = settingsPtr->parm(
"JetMatching:qCutME");
907 qCutMESq = pow(qCutME,2);
910 doMerge = settingsPtr->flag(
"JetMatching:merge");
911 doShowerKt = settingsPtr->flag(
"JetMatching:doShowerKt");
912 qCut = settingsPtr->parm(
"JetMatching:qCut");
913 nQmatch = settingsPtr->mode(
"JetMatching:nQmatch");
914 clFact = settingsPtr->parm(
"JetMatching:clFact");
917 jetAlgorithm = settingsPtr->mode(
"JetMatching:jetAlgorithm");
918 nJetMax = settingsPtr->mode(
"JetMatching:nJetMax");
919 eTjetMin = settingsPtr->parm(
"JetMatching:eTjetMin");
920 coneRadius = settingsPtr->parm(
"JetMatching:coneRadius");
921 etaJetMax = settingsPtr->parm(
"JetMatching:etaJetMax");
922 slowJetPower = settingsPtr->mode(
"JetMatching:slowJetPower");
925 jetAllow = settingsPtr->mode(
"JetMatching:jetAllow");
926 exclusiveMode = settingsPtr->mode(
"JetMatching:exclusive");
927 qCutSq = pow(qCut,2);
928 etaJetMaxAlgo = etaJetMax;
931 if (!doMerge)
return true;
934 if (exclusiveMode == 2) {
938 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:init: "
939 "missing jet multiplicity information; running in exclusive mode");
949 slowJet =
new SlowJet(slowJetPower, coneRadius, eTjetMin,
950 etaJetMaxAlgo, 2, 2, NULL,
false);
955 slowJetHard =
new SlowJet(slowJetPower, coneRadius, qCutME,
956 etaJetMaxAlgo, 2, 2, NULL,
false);
959 eventProcessOrig.init(
"(eventProcessOrig)", particleDataPtr);
960 eventProcess.init(
"(eventProcess)", particleDataPtr);
961 workEventJet.init(
"(workEventJet)", particleDataPtr);
964 string jetStr = (jetAlgorithm == 1) ?
"CellJet" :
965 (slowJetPower == -1) ?
"anti-kT" :
966 (slowJetPower == 0) ?
"C/A" :
967 (slowJetPower == 1) ?
"kT" :
"unknown";
968 string modeStr = (exclusiveMode) ?
"exclusive" :
"inclusive";
970 <<
" *----- Madgraph matching parameters -----*" << endl
971 <<
" | qCut | " << setw(14)
972 << qCut <<
" |" << endl
973 <<
" | nQmatch | " << setw(14)
974 << nQmatch <<
" |" << endl
975 <<
" | clFact | " << setw(14)
976 << clFact <<
" |" << endl
977 <<
" | Jet algorithm | " << setw(14)
978 << jetStr <<
" |" << endl
979 <<
" | eTjetMin | " << setw(14)
980 << eTjetMin <<
" |" << endl
981 <<
" | etaJetMax | " << setw(14)
982 << etaJetMax <<
" |" << endl
983 <<
" | jetAllow | " << setw(14)
984 << jetAllow <<
" |" << endl
985 <<
" | Mode | " << setw(14)
986 << modeStr <<
" |" << endl
987 <<
" *-----------------------------------------*" << endl;
994 bool JetMatchingMadgraph::doVetoStep(
int iPos,
int nISR,
int nFSR,
995 const Event& event) {
998 if ( !doShowerKt )
return false;
1002 if (
int(typeIdx[0].size()) > nJetMax )
return false;
1005 if ( nISR + nFSR > 1 )
return false;
1008 if (iPos == 5)
return false;
1012 sortIncomingProcess(event);
1015 double pTfirst = 0.;
1016 for (
int i = 0; i < workEvent.size(); i++){
1019 if ( workEvent[i].isFinal()
1020 && (workEvent[i].statusAbs()==43 || workEvent[i].statusAbs()==51)) {
1023 bool QCDemission =
true;
1024 while ( workEvent[jPos].statusAbs() > 23 ) {
1025 if ( workEvent[jPos].
id() == 22 || workEvent[jPos].id() == 23
1026 || workEvent[jPos].idAbs() == 24){
1027 QCDemission =
false;
1030 jPos = workEvent[jPos].mother1();
1034 pTfirst = workEvent[i].pT();
1041 if ( doShowerKtVeto(pTfirst) )
return true;
1050 bool JetMatchingMadgraph::doShowerKtVeto(
double pTfirst) {
1053 if ( !doShowerKt )
return false;
1056 bool doVeto =
false;
1060 int nParton = typeIdx[0].size();
1061 double pTminME=1e10;
1062 for (
int i = 0; i < nParton; ++i)
1063 pTminME = min(pTminME,eventProcess[typeIdx[0][i]].pT());
1066 if ( nParton > 0 && pow(pTminME,2) < qCutSq ) doVeto =
true;
1070 if ( exclusive && pow(pTfirst,2) > qCutSq ) {
1074 }
else if ( !exclusive && nParton > 0 && pTfirst > pTminME ) {
1087 void JetMatchingMadgraph::sortIncomingProcess(
const Event &event) {
1091 omitResonanceDecays(eventProcessOrig,
true);
1097 eventProcess.clear();
1098 workEventJet.clear();
1099 for(
int i=0; i < workEvent.size(); ++i) {
1102 int id = workEvent[i].idAbs();
1103 if ((
id >= ID_LEPMIN &&
id <= ID_LEPMAX) ||
id == ID_TOP
1104 ||
id == ID_PHOTON ||
id == 23 ||
id == 24 ||
id == 25) {
1105 eventProcess.append(workEvent[i]);
1107 workEventJet.append(workEvent[i]);
1112 if (!slowJetHard->setup(workEventJet) ) {
1113 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:sortIncomingProcess"
1114 ": the SlowJet algorithm failed on setup");
1119 double localQcutSq = qCutMESq;
1121 while ( slowJetHard->sizeAll() - slowJetHard->sizeJet() > 0 ) {
1123 if( slowJetHard->dNext() > localQcutSq )
break;
1125 if( slowJetHard->sizeAll()-slowJetHard->sizeJet() <= nPartonsNow)
break;
1126 slowJetHard->doStep();
1132 int nJets = slowJetHard->sizeJet();
1133 int nClus = slowJetHard->sizeAll();
1135 for (
int i = nJets; i < nClus; ++i) {
1137 if (i < nClus-nJets) parts = slowJetHard->clusConstituents(i);
1138 else parts = slowJetHard->constituents(nClus-nJets-i);
1139 int flavour = ID_GLUON;
1140 for(
int j=0; j < int(parts.size()); ++j)
1141 if (workEventJet[parts[j]].
id() == ID_BOT)
1143 eventProcess.append( flavour, 98,
1144 workEventJet[parts.back()].mother1(),
1145 workEventJet[parts.back()].mother2(),
1146 workEventJet[parts.back()].daughter1(),
1147 workEventJet[parts.back()].daughter2(),
1148 0, 0, slowJetHard->p(i).px(), slowJetHard->p(i).py(),
1149 slowJetHard->p(i).pz(), slowJetHard->p(i).e() );
1154 workEventJet.clear();
1159 eventProcess = workEvent;
1171 for (
int i = 0; i < 3; i++) {
1174 origTypeIdx[i].clear();
1176 for (
int i = 0; i < eventProcess.size(); i++) {
1178 if (!eventProcess[i].isFinal())
continue;
1183 if (eventProcess[i].
id() == ID_GLUON
1184 || (eventProcess[i].idAbs() <= nQmatch) ) {
1188 if ( eventProcess[i].scale() < 1.999*sqrt(infoPtr->eA()*infoPtr->eB()) )
1193 else if (eventProcess[i].idAbs() > nQmatch
1194 && eventProcess[i].idAbs() <= ID_TOP) {
1204 typeIdx[idx].push_back(i);
1205 typeSet[idx].insert(eventProcess[i].daughter1());
1206 origTypeIdx[orig_idx].push_back(i);
1210 if (exclusiveMode == 2) {
1213 int nParton = origTypeIdx[0].size();
1214 exclusive = (nParton == nJetMax) ?
false :
true;
1218 exclusive = (exclusiveMode == 0) ?
false :
true;
1230 void JetMatchingMadgraph::jetAlgorithmInput(
const Event &event,
int iType) {
1233 workEventJet = workEvent;
1236 for (
int i = 0; i < workEventJet.size(); ++i) {
1237 if (!workEventJet[i].isFinal())
continue;
1240 if (jetAllow == 1) {
1244 int id = workEventJet[i].idAbs();
1245 if ((
id >= ID_LEPMIN &&
id <= ID_LEPMAX) ||
id == ID_TOP
1246 ||
id == ID_PHOTON || (
id > nQmatch &&
id!=21)) {
1247 workEventJet[i].statusNeg();
1253 int idx = workEventJet[i].daughter1();
1262 if (typeSet[1].find(idx) != typeSet[1].end() ||
1263 typeSet[2].find(idx) != typeSet[2].end()) {
1264 workEventJet[i].statusNeg();
1269 if (idx == 0)
break;
1271 idx =
event[idx].mother1();
1274 }
else if (iType == 1) {
1277 if (typeSet[1].find(idx) != typeSet[1].end())
break;
1282 workEventJet[i].statusNeg();
1287 idx =
event[idx].mother1();
1300 void JetMatchingMadgraph::runJetAlgorithm() {; }
1306 bool JetMatchingMadgraph::matchPartonsToJets(
int iType) {
1310 if (iType == 0)
return (matchPartonsToJetsLight() > 0);
1311 else return (matchPartonsToJetsHeavy() > 0);
1327 int JetMatchingMadgraph::matchPartonsToJetsLight() {
1331 int nParton = typeIdx[0].size();
1334 if (!slowJet->setup(workEventJet) ) {
1335 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1336 "Light: the SlowJet algorithm failed on setup");
1339 double localQcutSq = qCutSq;
1342 while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1343 if( slowJet->dNext() > localQcutSq )
break;
1344 dOld = slowJet->dNext();
1347 int nJets = slowJet->sizeJet();
1348 int nClus = slowJet->sizeAll();
1351 if (MATCHINGDEBUG) slowJet->list(
true);
1354 int nCLjets = nClus - nJets;
1356 int nRequested = (doFxFx) ? nPartonsNow : nParton;
1359 if ( nCLjets < nRequested )
return LESS_JETS;
1362 if ( exclusive && !doFxFx ) {
1363 if ( nCLjets > nRequested )
return MORE_JETS;
1369 if ( doFxFx && nRequested < nJetMax && nCLjets > nRequested )
1376 if (!slowJet->setup(workEventJet) ) {
1377 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1378 "Light: the SlowJet algorithm failed on setup");
1385 while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1386 if( slowJet->dNext() > localQcutSq )
break;
1392 while ( slowJet->sizeAll() - slowJet->sizeJet() > nParton )
1399 if ( clFact >= 0. && nParton > 0 ) {
1400 vector<double> partonPt;
1401 for (
int i = 0; i < nParton; ++i)
1402 partonPt.push_back( eventProcess[typeIdx[0][i]].pT2() );
1403 sort( partonPt.begin(), partonPt.end());
1404 localQcutSq = max( qCutSq, partonPt[0]);
1406 nJets = slowJet->sizeJet();
1407 nClus = slowJet->sizeAll();
1410 if ( clFact != 0. ) localQcutSq *= pow2(clFact);
1413 tempEvent.init(
"(tempEvent)", particleDataPtr);
1415 double pTminEstimate = -1.;
1419 for (
int i = nJets; i < nClus; ++i) {
1420 tempEvent.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0, slowJet->p(i).px(),
1421 slowJet->p(i).py(), slowJet->p(i).pz(), slowJet->p(i).e() );
1423 pTminEstimate = max( pTminEstimate, slowJet->pT(i));
1424 if(nPass == nRequested)
break;
1427 int tempSize = tempEvent.size();
1429 vector<bool> jetAssigned;
1430 jetAssigned.assign( tempSize,
false);
1434 vector< vector<bool> > partonMatchesJet;
1435 for (
int i=0; i < nParton; ++i )
1436 partonMatchesJet.push_back( vector<bool>(tempEvent.size(),
false) );
1451 while ( doFxFx && iNow < tempSize ) {
1455 tempEventJet.init(
"(tempEventJet)", particleDataPtr);
1456 for (
int i=0; i < nParton; ++i ) {
1463 tempEventJet.clear();
1464 tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1465 tempEvent[iNow].px(), tempEvent[iNow].py(),
1466 tempEvent[iNow].pz(), tempEvent[iNow].e() );
1468 Vec4 pIn = eventProcess[typeIdx[0][i]].p();
1469 tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1470 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1473 if ( !slowJet->setup(tempEventJet) ) {
1474 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1475 "Light: the SlowJet algorithm failed on setup");
1481 if ( slowJet->iNext() == tempEventJet.size() - 1
1482 && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1483 jetAssigned[iNow] =
true;
1484 partonMatchesJet[i][iNow] =
true;
1490 if ( !jetAssigned[iNow] )
return UNMATCHED_PARTON;
1504 while (!doFxFx && iNow < nParton ) {
1506 tempEventJet.init(
"(tempEventJet)", particleDataPtr);
1507 for (
int i = 0; i < tempSize; ++i) {
1508 if (jetAssigned[i])
continue;
1509 Vec4 pIn = tempEvent[i].p();
1511 tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1512 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1515 Vec4 pIn = eventProcess[typeIdx[0][iNow]].p();
1517 tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1518 pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1519 if ( !slowJet->setup(tempEventJet) ) {
1520 infoPtr->errorMsg(
"Warning in JetMatchingMadgraph:matchPartonsToJets"
1521 "Light: the SlowJet algorithm failed on setup");
1526 if ( slowJet->iNext() == tempEventJet.size() - 1
1527 && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1529 for (
int i = 0; i != tempSize; ++i) {
1530 if (jetAssigned[i])
continue;
1533 if (iKnt == slowJet->jNext() ) jetAssigned[i] =
true;
1536 return UNMATCHED_PARTON;
1544 if (nParton > 0 && pTminEstimate > 0) eTpTlightMin = pTminEstimate;
1545 else eTpTlightMin = -1.;
1561 int JetMatchingMadgraph::matchPartonsToJetsHeavy() {
1565 if (jetMomenta.empty())
return NONE;
1573 #endif // end Pythia8_JetMatching_H