9 #include "Pythia8/Merging.h"
21 const double Merging::TMSMISMATCH = 1.5;
30 tmsNowMin = infoPtr->eCM();
37 void Merging::statistics() {
40 bool enforceCutOnLHE = settingsPtr->flag(
"Merging:enforceCutOnLHE");
42 double tmsval = mergingHooksPtr->tms();
43 bool printBanner = enforceCutOnLHE && tmsNowMin > TMSMISMATCH*tmsval;
45 tmsNowMin = infoPtr->eCM();
47 if (!printBanner)
return;
50 cout <<
"\n *------- PYTHIA Matrix Element Merging Information ------"
51 <<
"-------------------------------------------------------*\n"
56 cout <<
" | Warning in Merging::statistics: All Les Houches events"
57 <<
" significantly above Merging:TMS cut. Please check. |\n";
62 <<
" *------- End PYTHIA Matrix Element Merging Information -----"
63 <<
"-----------------------------------------------------*" << endl;
70 int Merging::mergeProcess(
Event& process){
75 mergingHooksPtr->hardProcess->clear();
76 mergingHooksPtr->processSave = settingsPtr->word(
"Merging:Process");
77 mergingHooksPtr->hardProcess->initOnProcess(
78 settingsPtr->word(
"Merging:Process"), particleDataPtr);
80 mergingHooksPtr->doUserMergingSave
81 = settingsPtr->flag(
"Merging:doUserMerging");
82 mergingHooksPtr->doMGMergingSave
83 = settingsPtr->flag(
"Merging:doMGMerging");
84 mergingHooksPtr->doKTMergingSave
85 = settingsPtr->flag(
"Merging:doKTMerging");
86 mergingHooksPtr->doPTLundMergingSave
87 = settingsPtr->flag(
"Merging:doPTLundMerging");
88 mergingHooksPtr->doCutBasedMergingSave
89 = settingsPtr->flag(
"Merging:doCutBasedMerging");
90 mergingHooksPtr->doNL3TreeSave
91 = settingsPtr->flag(
"Merging:doNL3Tree");
92 mergingHooksPtr->doNL3LoopSave
93 = settingsPtr->flag(
"Merging:doNL3Loop");
94 mergingHooksPtr->doNL3SubtSave
95 = settingsPtr->flag(
"Merging:doNL3Subt");
96 mergingHooksPtr->doUNLOPSTreeSave
97 = settingsPtr->flag(
"Merging:doUNLOPSTree");
98 mergingHooksPtr->doUNLOPSLoopSave
99 = settingsPtr->flag(
"Merging:doUNLOPSLoop");
100 mergingHooksPtr->doUNLOPSSubtSave
101 = settingsPtr->flag(
"Merging:doUNLOPSSubt");
102 mergingHooksPtr->doUNLOPSSubtNLOSave
103 = settingsPtr->flag(
"Merging:doUNLOPSSubtNLO");
104 mergingHooksPtr->doUMEPSTreeSave
105 = settingsPtr->flag(
"Merging:doUMEPSTree");
106 mergingHooksPtr->doUMEPSSubtSave
107 = settingsPtr->flag(
"Merging:doUMEPSSubt");
108 mergingHooksPtr->nReclusterSave
109 = settingsPtr->mode(
"Merging:nRecluster");
111 mergingHooksPtr->hasJetMaxLocal =
false;
112 mergingHooksPtr->nJetMaxLocal
113 = mergingHooksPtr->nJetMaxSave;
114 mergingHooksPtr->nJetMaxNLOLocal
115 = mergingHooksPtr->nJetMaxNLOSave;
116 mergingHooksPtr->nRequestedSave
117 = settingsPtr->mode(
"Merging:nRequested");
120 bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
123 bool applyTMSCut = settingsPtr->flag(
"Merging:doXSectionEstimate");
124 if ( applyTMSCut && cutOnProcess(process) ) {
125 if (includeWGT) infoPtr->updateWeight(0.);
129 if ( applyTMSCut )
return 1;
132 if ( mergingHooksPtr->doCKKWLMerging() )
133 vetoCode = mergeProcessCKKWL(process);
136 if ( mergingHooksPtr->doUMEPSMerging() )
137 vetoCode = mergeProcessUMEPS(process);
140 if ( mergingHooksPtr->doNL3Merging() )
141 vetoCode = mergeProcessNL3(process);
144 if ( mergingHooksPtr->doUNLOPSMerging() )
145 vetoCode = mergeProcessUNLOPS(process);
155 int Merging::mergeProcessCKKWL(
Event& process) {
158 mergingHooksPtr->doIgnoreStep(
true);
161 if ( mergingHooksPtr->getProcessString().compare(
"pp>h") == 0 )
162 mergingHooksPtr->allowCutOnRecState(
true);
164 mergingHooksPtr->orderHistories(
true);
167 bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
171 mergingHooksPtr->setWeightCKKWL(1.);
172 mergingHooksPtr->muMI(-1.);
177 Event newProcess( mergingHooksPtr->bareEvent( process,
true) );
179 if (mergingHooksPtr->doWeakClustering())
180 for (
int i = 0;i < newProcess.size();++i)
181 newProcess[i].pol(9);
183 mergingHooksPtr->storeHardProcessCandidates( newProcess);
186 double tmsval = mergingHooksPtr->tms();
188 double tmsnow = mergingHooksPtr->tmsNow( newProcess );
190 int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess,
true);
193 bool allowReject = settingsPtr->flag(
"Merging:applyVeto");
198 int nRequested = mergingHooksPtr->nRequested();
201 mergingHooksPtr->setHardProcessInfo(nSteps, tmsnow);
202 mergingHooksPtr->setEventVetoInfo(-1, -1.);
204 if (nSteps < nRequested && allowReject) {
205 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
206 if ( includeWGT) infoPtr->updateWeight(0.);
211 tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
214 double RN = rndmPtr->flat();
217 newProcess.scale(0.0);
219 History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
220 (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
221 trialPartonLevelPtr, coupSMPtr,
true,
true,
true,
true, 1.0, 0);
224 FullHistory.projectOntoDesiredHistories();
227 FullHistory.select(RN)->setSelectedChild();
231 bool applyCut = allowReject
232 && nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
236 bool enforceCutOnLHE = settingsPtr->flag(
"Merging:enforceCutOnLHE");
237 if ( enforceCutOnLHE && applyCut && tmsnow < tmsval ) {
238 string message=
"Warning in Merging::mergeProcessCKKWL: Les Houches Event";
239 message+=
" fails merging scale cut. Reject event.";
240 infoPtr->errorMsg(message);
241 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
242 if ( includeWGT) infoPtr->updateWeight(0.);
251 coreProcess.init(
"(hard process-modified)", particleDataPtr );
253 coreProcess = FullHistory.lowestMultProc(RN);
254 for (
int i = 0; i < coreProcess.size(); ++i )
255 if ( coreProcess[i].isFinal() ) {
256 if ( coreProcess[i].colType() != 0 )
258 if ( coreProcess[i].idAbs() == 24 )
262 bool complete = (FullHistory.select(RN)->nClusterings() == nSteps) ||
263 ( mergingHooksPtr->doWeakClustering() && nFinalP == 2 && nFinalW == 0 );
266 string message=
"Warning in Merging::mergeProcessCKKWL: No clusterings";
267 message+=
" found. History incomplete.";
268 infoPtr->errorMsg(message);
274 wgt = FullHistory.weightTREE( trialPartonLevelPtr,
275 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
276 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
280 FullHistory.getStartingConditions( RN, process );
282 mergingHooksPtr->reattachResonanceDecays(process);
286 double dampWeight = mergingHooksPtr->dampenIfFailCuts(
287 FullHistory.lowestMultProc(RN) );
294 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(wgt);
297 double norm = (abs(infoPtr->lhaStrategy()) == 4) ? 1/1e9 : 1.;
298 if ( includeWGT) infoPtr->updateWeight(infoPtr->weight()*wgt*norm);
301 mergingHooksPtr->doIgnoreStep(
false);
304 if ( allowReject && wgt == 0. )
return 0;
315 int Merging::mergeProcessUMEPS(
Event& process) {
318 bool doUMEPSTree = settingsPtr->flag(
"Merging:doUMEPSTree");
319 bool doUMEPSSubt = settingsPtr->flag(
"Merging:doUMEPSSubt");
321 mergingHooksPtr->nReclusterSave = settingsPtr->mode(
"Merging:nRecluster");
322 int nRecluster = settingsPtr->mode(
"Merging:nRecluster");
325 mergingHooksPtr->doIgnoreEmissions(
true);
328 if ( mergingHooksPtr->getProcessString().compare(
"pp>h") == 0 )
329 mergingHooksPtr->allowCutOnRecState(
true);
331 mergingHooksPtr->orderHistories(
true);
334 bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
337 if (mergingHooksPtr->doWeakClustering())
338 for (
int i = 0;i < process.size();++i)
343 mergingHooksPtr->setWeightCKKWL(1.);
344 mergingHooksPtr->muMI(-1.);
349 Event newProcess( mergingHooksPtr->bareEvent( process,
true) );
351 mergingHooksPtr->storeHardProcessCandidates( newProcess );
354 double tmsval = mergingHooksPtr->tms();
356 double tmsnow = mergingHooksPtr->tmsNow( newProcess );
358 int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess,
true);
359 int nRequested = mergingHooksPtr->nRequested();
364 if (nSteps < nRequested) {
365 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
366 if ( includeWGT) infoPtr->updateWeight(0.);
371 tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
374 double RN = rndmPtr->flat();
376 newProcess.scale(0.0);
378 History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
379 (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
380 trialPartonLevelPtr, coupSMPtr,
true,
true,
true,
true, 1.0, 0);
382 FullHistory.projectOntoDesiredHistories();
386 bool applyCut = nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
390 bool enforceCutOnLHE = settingsPtr->flag(
"Merging:enforceCutOnLHE");
391 if ( enforceCutOnLHE && applyCut && tmsnow < tmsval ) {
392 string message=
"Warning in Merging::mergeProcessUMEPS: Les Houches Event";
393 message+=
" fails merging scale cut. Reject event.";
394 infoPtr->errorMsg(message);
395 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
396 if ( includeWGT) infoPtr->updateWeight(0.);
402 if ( nSteps > 0 && doUMEPSSubt
403 && !FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster,
404 newProcess, nPerformed,
false ) ) {
406 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
407 if ( includeWGT) infoPtr->updateWeight(0.);
411 mergingHooksPtr->nMinMPI(nSteps - nPerformed);
417 wgt = FullHistory.weight_UMEPS_TREE( trialPartonLevelPtr,
418 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
419 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
421 wgt = FullHistory.weight_UMEPS_SUBT( trialPartonLevelPtr,
422 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
423 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
428 if ( doUMEPSTree ) FullHistory.getStartingConditions( RN, process );
430 else FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster, process,
435 double dampWeight = mergingHooksPtr->dampenIfFailCuts(
436 FullHistory.lowestMultProc(RN) );
443 if (!includeWGT) mergingHooksPtr->setWeightCKKWL(wgt);
446 double norm = (abs(infoPtr->lhaStrategy()) == 4) ? 1/1e9 : 1.;
447 if ( includeWGT) infoPtr->updateWeight(infoPtr->weight()*wgt*norm);
452 double muf = process[0].e();
453 for (
int i=0; i < process.size(); ++i )
454 if ( process[i].isFinal()
455 && (process[i].colType() != 0 || process[i].id() == 22 ) ) {
457 muf = min( muf, abs(process[i].mT()) );
463 int nStepsNew = mergingHooksPtr->getNumberOfClusteringSteps( process );
465 && ( mergingHooksPtr->getProcessString().compare(
"pp>jj") == 0
466 || mergingHooksPtr->getProcessString().compare(
"pp>aj") == 0) )
470 mergingHooksPtr->storeHardProcessCandidates( process );
472 mergingHooksPtr->reattachResonanceDecays(process);
475 mergingHooksPtr->doIgnoreEmissions(
false);
478 if ( wgt == 0. )
return 0;
489 int Merging::mergeProcessNL3(
Event& process) {
492 bool doNL3Tree = settingsPtr->flag(
"Merging:doNL3Tree");
493 bool doNL3Loop = settingsPtr->flag(
"Merging:doNL3Loop");
494 bool doNL3Subt = settingsPtr->flag(
"Merging:doNL3Subt");
497 mergingHooksPtr->doIgnoreEmissions(
true);
499 mergingHooksPtr->doIgnoreStep(
true);
502 if ( mergingHooksPtr->getProcessString().compare(
"pp>h") == 0)
503 mergingHooksPtr->allowCutOnRecState(
true);
505 mergingHooksPtr->orderHistories(
true);
509 mergingHooksPtr->setWeightCKKWL(1.);
511 double wgtFIRST = 0.;
512 mergingHooksPtr->setWeightFIRST(0.);
513 mergingHooksPtr->muMI(-1.);
518 Event newProcess( mergingHooksPtr->bareEvent( process,
true) );
520 mergingHooksPtr->storeHardProcessCandidates( newProcess);
523 double tmsval = mergingHooksPtr->tms();
525 double tmsnow = mergingHooksPtr->tmsNow( newProcess );
527 int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess,
true);
528 int nRequested = mergingHooksPtr->nRequested();
533 if (nSteps < nRequested) {
534 mergingHooksPtr->setWeightCKKWL(0.);
535 mergingHooksPtr->setWeightFIRST(0.);
540 tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
544 bool enforceCutOnLHE = settingsPtr->flag(
"Merging:enforceCutOnLHE");
545 if ( enforceCutOnLHE && nSteps > 0 && nSteps == nRequested
546 && tmsnow < tmsval ) {
547 string message=
"Warning in Merging::mergeProcessNL3: Les Houches Event";
548 message+=
" fails merging scale cut. Reject event.";
549 infoPtr->errorMsg(message);
550 mergingHooksPtr->setWeightCKKWL(0.);
551 mergingHooksPtr->setWeightFIRST(0.);
556 double RN = rndmPtr->flat();
558 newProcess.scale(0.0);
560 History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
561 (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
562 trialPartonLevelPtr, coupSMPtr,
true,
true,
true,
true, 1.0, 0);
564 FullHistory.projectOntoDesiredHistories();
567 if ( nSteps > 0 && doNL3Subt
568 && FullHistory.select(RN)->nClusterings() == 0 ){
569 mergingHooksPtr->setWeightCKKWL(0.);
570 mergingHooksPtr->setWeightFIRST(0.);
576 bool containsRealKin = nSteps > nRequested && nSteps > 0;
580 if ( containsRealKin ) {
584 dummy.init(
"(hard process-modified)", particleDataPtr );
587 if ( !FullHistory.getClusteredEvent( RN, nSteps, dummy )) {
588 mergingHooksPtr->setWeightCKKWL(0.);
589 mergingHooksPtr->setWeightFIRST(0.);
592 double tnowNew = mergingHooksPtr->tmsNow( dummy );
594 if ( enforceCutOnLHE && nSteps > 0 && nRequested > 0
595 && tnowNew < tmsval ) {
596 mergingHooksPtr->setWeightCKKWL(0.);
597 mergingHooksPtr->setWeightFIRST(0.);
603 if ( doNL3Subt || containsRealKin ) mergingHooksPtr->nMinMPI(nSteps - 1);
604 else mergingHooksPtr->nMinMPI(nSteps);
611 wgt = FullHistory.weightTREE( trialPartonLevelPtr,
612 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
613 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
614 }
else if( doNL3Loop || doNL3Subt ) {
617 wgt = FullHistory.weightLOOP( trialPartonLevelPtr, RN);
622 if ( !doNL3Subt && !containsRealKin )
623 FullHistory.getStartingConditions(RN, process);
629 if ( !FullHistory.getClusteredEvent( RN, nSteps, process )) {
630 mergingHooksPtr->setWeightCKKWL(0.);
631 mergingHooksPtr->setWeightFIRST(0.);
638 double dampWeight = mergingHooksPtr->dampenIfFailCuts(
639 FullHistory.lowestMultProc(RN) );
649 if( nSteps > mergingHooksPtr->nMaxJetsNLO() )
650 kFactor = mergingHooksPtr->kFactor( mergingHooksPtr->nMaxJetsNLO() );
651 else kFactor = mergingHooksPtr->kFactor(nSteps);
657 mergingHooksPtr->setWeightCKKWL(wgt);
662 bool doOASTree = doNL3Tree && nSteps <= mergingHooksPtr->nMaxJetsNLO();
667 wgtFIRST = FullHistory.weightFIRST( trialPartonLevelPtr,
668 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
669 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN,
672 wgtFIRST *= dampWeight;
674 mergingHooksPtr->setWeightFIRST(wgtFIRST);
677 wgt = wgt - wgtFIRST;
683 for(
int i=0; i < process.size(); ++i)
684 if(process[i].isFinal() && process[i].colType() != 0) {
685 pT = sqrt(pow(process[i].px(),2) + pow(process[i].py(),2));
691 && mergingHooksPtr->getProcessString().compare(
"pp>jj") == 0)
695 mergingHooksPtr->storeHardProcessCandidates( process );
697 mergingHooksPtr->reattachResonanceDecays(process);
700 mergingHooksPtr->doIgnoreEmissions(
false);
702 mergingHooksPtr->doIgnoreStep(
false);
713 int Merging::mergeProcessUNLOPS(
Event& process) {
716 bool nloTilde = settingsPtr->flag(
"Merging:doUNLOPSTilde");
717 bool doUNLOPSTree = settingsPtr->flag(
"Merging:doUNLOPSTree");
718 bool doUNLOPSLoop = settingsPtr->flag(
"Merging:doUNLOPSLoop");
719 bool doUNLOPSSubt = settingsPtr->flag(
"Merging:doUNLOPSSubt");
720 bool doUNLOPSSubtNLO = settingsPtr->flag(
"Merging:doUNLOPSSubtNLO");
722 mergingHooksPtr->nReclusterSave = settingsPtr->mode(
"Merging:nRecluster");
723 int nRecluster = settingsPtr->mode(
"Merging:nRecluster");
726 mergingHooksPtr->doIgnoreEmissions(
true);
728 mergingHooksPtr->orderHistories(
true);
731 if ( mergingHooksPtr->getProcessString().compare(
"pp>h") == 0)
732 mergingHooksPtr->allowCutOnRecState(
true);
736 mergingHooksPtr->setWeightCKKWL(1.);
738 double wgtFIRST = 0.;
739 mergingHooksPtr->setWeightFIRST(0.);
740 mergingHooksPtr->muMI(-1.);
745 Event newProcess( mergingHooksPtr->bareEvent( process,
true) );
747 mergingHooksPtr->storeHardProcessCandidates( newProcess );
750 double tmsval = mergingHooksPtr->tms();
752 double tmsnow = mergingHooksPtr->tmsNow( newProcess );
754 int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess,
true);
755 int nRequested = mergingHooksPtr->nRequested();
760 if (nSteps < nRequested) {
761 string message=
"Warning in Merging::mergeProcessUNLOPS: Les Houches Event";
762 message+=
" after removing decay products does not contain enough partons.";
763 infoPtr->errorMsg(message);
764 mergingHooksPtr->setWeightCKKWL(0.);
765 mergingHooksPtr->setWeightFIRST(0.);
770 tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
773 double RN = rndmPtr->flat();
775 newProcess.scale(0.0);
777 History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
778 (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
779 trialPartonLevelPtr, coupSMPtr,
true,
true,
true,
true, 1.0, 0);
781 FullHistory.projectOntoDesiredHistories();
785 bool applyCut = nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
789 bool enforceCutOnLHE = settingsPtr->flag(
"Merging:enforceCutOnLHE");
790 if ( enforceCutOnLHE && applyCut && nSteps == nRequested
791 && tmsnow < tmsval ) {
792 string message=
"Warning in Merging::mergeProcessUNLOPS: Les Houches";
793 message+=
" Event fails merging scale cut. Reject event.";
794 infoPtr->errorMsg(message);
795 mergingHooksPtr->setWeightCKKWL(0.);
796 mergingHooksPtr->setWeightFIRST(0.);
802 bool containsRealKin = nSteps > nRequested && nSteps > 0;
803 if ( containsRealKin ) nRecluster += nSteps - nRequested;
808 bool allowIncompleteReal =
809 settingsPtr->flag(
"Merging:allowIncompleteHistoriesInReal");
810 if ( doUNLOPSLoop && containsRealKin && !allowIncompleteReal
811 && FullHistory.select(RN)->nClusterings() == 0 ) {
812 mergingHooksPtr->setWeightCKKWL(0.);
813 mergingHooksPtr->setWeightFIRST(0.);
819 if ( nSteps > 0 && !allowIncompleteReal
820 && ( doUNLOPSSubt || doUNLOPSSubtNLO || containsRealKin )
821 && !FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster,
822 newProcess, nPerformed,
false ) ) {
823 mergingHooksPtr->setWeightCKKWL(0.);
824 mergingHooksPtr->setWeightFIRST(0.);
829 mergingHooksPtr->nMinMPI(nSteps - nPerformed);
833 if ( containsRealKin ) {
837 dummy.init(
"(hard process-modified)", particleDataPtr );
840 FullHistory.getClusteredEvent( RN, nSteps, dummy );
841 double tnowNew = mergingHooksPtr->tmsNow( dummy );
843 if ( enforceCutOnLHE && nSteps > 0 && nRequested > 0
844 && tnowNew < tmsval ) {
845 string message=
"Warning in Merging::mergeProcessUNLOPS: Les Houches";
846 message+=
" Event fails merging scale cut. Reject event.";
847 infoPtr->errorMsg(message);
848 mergingHooksPtr->setWeightCKKWL(0.);
849 mergingHooksPtr->setWeightFIRST(0.);
855 bool doUNLOPS2 =
false;
856 int depth = (!doUNLOPS2) ? -1 : ( (containsRealKin) ? nSteps-1 : nSteps);
863 wgt = FullHistory.weight_UNLOPS_TREE( trialPartonLevelPtr,
864 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
865 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
867 }
else if( doUNLOPSLoop ) {
869 wgt = FullHistory.weight_UNLOPS_LOOP( trialPartonLevelPtr,
870 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
871 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
873 }
else if( doUNLOPSSubtNLO ) {
875 wgt = FullHistory.weight_UNLOPS_SUBTNLO( trialPartonLevelPtr,
876 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
877 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
879 }
else if( doUNLOPSSubt ) {
882 wgt = FullHistory.weight_UNLOPS_SUBT( trialPartonLevelPtr,
883 mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
884 mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
890 if (!doUNLOPSSubt && !doUNLOPSSubtNLO && !containsRealKin )
891 FullHistory.getStartingConditions(RN, process);
893 else FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster, process,
898 double dampWeight = mergingHooksPtr->dampenIfFailCuts(
899 FullHistory.lowestMultProc(RN) );
906 if ( doUNLOPSTree || doUNLOPSSubt ){
909 if ( nSteps > mergingHooksPtr->nMaxJetsNLO() )
910 kFactor = mergingHooksPtr->kFactor( mergingHooksPtr->nMaxJetsNLO() );
911 else kFactor = mergingHooksPtr->kFactor(nSteps);
913 wgt *= (nRecluster == 2 && nloTilde) ? 1. : kFactor;
917 mergingHooksPtr->setWeightCKKWL(wgt);
922 int nMaxNLO = mergingHooksPtr->nMaxJetsNLO();
923 bool doOASTree = doUNLOPSTree && nSteps <= nMaxNLO;
924 bool doOASSubt = doUNLOPSSubt && nSteps <= nMaxNLO+1 && nSteps > 0;
927 if ( doOASTree || doOASSubt ) {
930 int order = ( nSteps > 0 && nSteps <= nMaxNLO) ? 1 : -1;
935 if ( nloTilde && doUNLOPSSubt && nRecluster == 1
936 && nSteps == nMaxNLO+1 ) order = 0;
942 if (nloTilde && doUNLOPSSubt && ( nSteps > nMaxNLO+1
943 || (nSteps == nMaxNLO+1 && nPerformed != nRecluster) ))
947 wgtFIRST = FullHistory.weight_UNLOPS_CORRECTION( order,
948 trialPartonLevelPtr, mergingHooksPtr->AlphaS_FSR(),
949 mergingHooksPtr->AlphaS_ISR(), mergingHooksPtr->AlphaEM_FSR(),
950 mergingHooksPtr->AlphaEM_ISR(), RN, rndmPtr );
955 if ( nloTilde && doUNLOPSSubt && nRecluster == 1
956 && nPerformed == nRecluster && nSteps <= nMaxNLO )
960 wgtFIRST *= dampWeight;
963 mergingHooksPtr->setWeightFIRST(wgtFIRST);
967 if (doUNLOPS2 && order > -1) wgt = -wgt*(wgtFIRST-1.);
968 else if (order > -1) wgt = wgt - wgtFIRST;
975 double muf = process[0].e();
976 for (
int i=0; i < process.size(); ++i )
977 if ( process[i].isFinal()
978 && (process[i].colType() != 0 || process[i].id() == 22 ) ) {
980 muf = min( muf, abs(process[i].mT()) );
984 if ( nSteps == 0 && nFinal == 2
985 && ( mergingHooksPtr->getProcessString().compare(
"pp>jj") == 0
986 || mergingHooksPtr->getProcessString().compare(
"pp>aj") == 0) )
990 mergingHooksPtr->storeHardProcessCandidates( process );
994 vector <int> oldResonance;
995 for (
int i=0; i < newProcess.size(); ++i )
996 if ( newProcess[i].status() == 22 )
997 oldResonance.push_back(newProcess[i].id());
998 vector <int> newResonance;
999 for (
int i=0; i < process.size(); ++i )
1000 if ( process[i].status() == 22 )
1001 newResonance.push_back(process[i].id());
1003 for (
int i=0; i < int(oldResonance.size()); ++i )
1004 for (
int j=0; j < int(newResonance.size()); ++j )
1005 if ( newResonance[j] == oldResonance[i] ) {
1006 oldResonance[i] = 99;
1009 bool hasNewResonances = (newResonance.size() != oldResonance.size());
1010 for (
int i=0; i < int(oldResonance.size()); ++i )
1011 hasNewResonances = (oldResonance[i] != 99);
1014 if (!hasNewResonances) mergingHooksPtr->reattachResonanceDecays(process);
1017 mergingHooksPtr->doIgnoreEmissions(
false);
1020 if ( wgt == 0. )
return 0;
1024 if (hasNewResonances)
return 2;
1035 bool Merging::cutOnProcess(
Event& process) {
1038 mergingHooksPtr->nReclusterSave = settingsPtr->mode(
"Merging:nRecluster");
1041 mergingHooksPtr->orderHistories(
true);
1044 if ( mergingHooksPtr->getProcessString().compare(
"pp>h") == 0)
1045 mergingHooksPtr->allowCutOnRecState(
true);
1048 if (mergingHooksPtr->doWeakClustering())
1049 for (
int i = 0;i < process.size();++i)
1055 Event newProcess( mergingHooksPtr->bareEvent( process,
true) );
1057 mergingHooksPtr->storeHardProcessCandidates( newProcess );
1060 double tmsval = mergingHooksPtr->tms();
1062 double tmsnow = mergingHooksPtr->tmsNow( newProcess );
1064 int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess,
true);
1069 int nRequested = mergingHooksPtr->nRequested();
1070 if (nSteps < nRequested)
return true;
1073 tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
1077 bool containsRealKin = nSteps > nRequested && nSteps > 0;
1080 double RN = rndmPtr->flat();
1082 newProcess.scale(0.0);
1084 History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
1085 (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
1086 trialPartonLevelPtr, coupSMPtr,
true,
true,
true,
true, 1.0, 0);
1088 FullHistory.projectOntoDesiredHistories();
1093 bool allowIncompleteReal =
1094 settingsPtr->flag(
"Merging:allowIncompleteHistoriesInReal");
1095 if ( containsRealKin && !allowIncompleteReal
1096 && FullHistory.select(RN)->nClusterings() == 0 )
1100 double dampWeight = mergingHooksPtr->dampenIfFailCuts(
1101 FullHistory.lowestMultProc(RN) );
1102 if ( dampWeight == 0. )
return true;
1106 if ( nSteps > 0 && FullHistory.select(RN)->nClusterings() == 0 )
1111 if ( nSteps > 0 && nSteps == nRequested && tmsnow < tmsval ) {
1112 string message=
"Warning in Merging::cutOnProcess: Les Houches Event";
1113 message+=
" fails merging scale cut. Reject event.";
1114 infoPtr->errorMsg(message);
1122 coreProcess.clear();
1123 coreProcess.init(
"(hard process-modified)", particleDataPtr );
1124 coreProcess.clear();
1125 coreProcess = FullHistory.lowestMultProc(RN);
1126 for (
int i = 0; i < coreProcess.size(); ++i )
1127 if ( coreProcess[i].isFinal() ) {
1128 if ( coreProcess[i].colType() != 0 )
1130 if ( coreProcess[i].idAbs() == 24 )
1134 bool complete = (FullHistory.select(RN)->nClusterings() == nSteps) ||
1135 ( mergingHooksPtr->doWeakClustering() && nFinalP == 2 && nFinalW == 0 );
1138 string message=
"Warning in Merging::cutOnProcess: No clusterings";
1139 message+=
" found. History incomplete.";
1140 infoPtr->errorMsg(message);
1144 if ( !containsRealKin )
return false;
1149 if ( containsRealKin ) {
1153 dummy.init(
"(hard process-modified)", particleDataPtr );
1156 FullHistory.getClusteredEvent( RN, nSteps, dummy );
1157 double tnowNew = mergingHooksPtr->tmsNow( dummy );
1159 if ( nSteps > 0 && nRequested > 0 && tnowNew < tmsval ) {
1160 string message=
"Warning in Merging::cutOnProcess: Les Houches Event";
1161 message+=
" fails merging scale cut. Reject event.";
1162 infoPtr->errorMsg(message);