StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Merging.cc
1 // MergingHooks.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This file is written by Stefan Prestel.
7 // Function definitions (not found in the header) for the Merging class.
8 
9 #include "Pythia8/Merging.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The Merging class.
16 
17 //--------------------------------------------------------------------------
18 
19 // Factor by which the maximal value of the merging scale can deviate before
20 // a warning is printed.
21 const double Merging::TMSMISMATCH = 1.5;
22 
23 //--------------------------------------------------------------------------
24 
25 // Initialise Merging class
26 
27 void Merging::init(){
28 
29  // Reset minimal tms value.
30  tmsNowMin = infoPtr->eCM();
31 
32 }
33 
34 //--------------------------------------------------------------------------
35 
36 // Function to print information.
37 void Merging::statistics() {
38 
39  // Recall switch to enfore merging scale cut.
40  bool enforceCutOnLHE = settingsPtr->flag("Merging:enforceCutOnLHE");
41  // Recall merging scale value.
42  double tmsval = mergingHooksPtr->tms();
43  bool printBanner = enforceCutOnLHE && tmsNowMin > TMSMISMATCH*tmsval;
44  // Reset minimal tms value.
45  tmsNowMin = infoPtr->eCM();
46 
47  if (!printBanner) return;
48 
49  // Header.
50  cout << "\n *------- PYTHIA Matrix Element Merging Information ------"
51  << "-------------------------------------------------------*\n"
52  << " | "
53  << " |\n";
54  // Print warning if the minimal tms value of any event was significantly
55  // above the desired merging scale value.
56  cout << " | Warning in Merging::statistics: All Les Houches events"
57  << " significantly above Merging:TMS cut. Please check. |\n";
58 
59  // Listing finished.
60  cout << " | "
61  << " |\n"
62  << " *------- End PYTHIA Matrix Element Merging Information -----"
63  << "-----------------------------------------------------*" << endl;
64 }
65 
66 //--------------------------------------------------------------------------
67 
68 // Function to steer different merging prescriptions.
69 
70 int Merging::mergeProcess(Event& process){
71 
72  int vetoCode = 1;
73 
74  // Reinitialise hard process.
75  mergingHooksPtr->hardProcess->clear();
76  mergingHooksPtr->processSave = settingsPtr->word("Merging:Process");
77  mergingHooksPtr->hardProcess->initOnProcess(
78  settingsPtr->word("Merging:Process"), particleDataPtr);
79 
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");
110 
111  mergingHooksPtr->hasJetMaxLocal = false;
112  mergingHooksPtr->nJetMaxLocal
113  = mergingHooksPtr->nJetMaxSave;
114  mergingHooksPtr->nJetMaxNLOLocal
115  = mergingHooksPtr->nJetMaxNLOSave;
116  mergingHooksPtr->nRequestedSave
117  = settingsPtr->mode("Merging:nRequested");
118 
119  // Ensure that merging weight is not counted twice.
120  bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
121 
122  // Possibility to apply merging scale to an input event.
123  bool applyTMSCut = settingsPtr->flag("Merging:doXSectionEstimate");
124  if ( applyTMSCut && cutOnProcess(process) ) {
125  if (includeWGT) infoPtr->updateWeight(0.);
126  return -1;
127  }
128  // Done if only a cut should be applied.
129  if ( applyTMSCut ) return 1;
130 
131  // Possibility to perform CKKW-L merging on this event.
132  if ( mergingHooksPtr->doCKKWLMerging() )
133  vetoCode = mergeProcessCKKWL(process);
134 
135  // Possibility to perform UMEPS merging on this event.
136  if ( mergingHooksPtr->doUMEPSMerging() )
137  vetoCode = mergeProcessUMEPS(process);
138 
139  // Possibility to perform NL3 NLO merging on this event.
140  if ( mergingHooksPtr->doNL3Merging() )
141  vetoCode = mergeProcessNL3(process);
142 
143  // Possibility to perform UNLOPS merging on this event.
144  if ( mergingHooksPtr->doUNLOPSMerging() )
145  vetoCode = mergeProcessUNLOPS(process);
146 
147  return vetoCode;
148 
149 }
150 
151 //--------------------------------------------------------------------------
152 
153 // Function to perform CKKW-L merging on this event.
154 
155 int Merging::mergeProcessCKKWL( Event& process) {
156 
157  // Ensure that merging hooks to not veto events in the trial showers.
158  mergingHooksPtr->doIgnoreStep(true);
159  // For pp > h, allow cut on state, so that underlying processes
160  // can be clustered to gg > h
161  if ( mergingHooksPtr->getProcessString().compare("pp>h") == 0 )
162  mergingHooksPtr->allowCutOnRecState(true);
163  // For now, prefer construction of ordered histories.
164  mergingHooksPtr->orderHistories(true);
165 
166  // Ensure that merging weight is not counted twice.
167  bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
168 
169  // Reset weight of the event.
170  double wgt = 1.0;
171  mergingHooksPtr->setWeightCKKWL(1.);
172  mergingHooksPtr->muMI(-1.);
173 
174  // Prepare process record for merging. If Pythia has already decayed
175  // resonances used to define the hard process, remove resonance decay
176  // products.
177  Event newProcess( mergingHooksPtr->bareEvent( process, true) );
178  // Reset any incoming spins for W+-.
179  if (mergingHooksPtr->doWeakClustering())
180  for (int i = 0;i < newProcess.size();++i)
181  newProcess[i].pol(9);
182  // Store candidates for the splitting V -> qqbar'.
183  mergingHooksPtr->storeHardProcessCandidates( newProcess);
184 
185  // Check if event passes the merging scale cut.
186  double tmsval = mergingHooksPtr->tms();
187  // Get merging scale in current event.
188  double tmsnow = mergingHooksPtr->tmsNow( newProcess );
189  // Calculate number of clustering steps.
190  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess, true);
191 
192  // Check if hard event cut should be applied later.
193  bool allowReject = settingsPtr->flag("Merging:applyVeto");
194 
195  // Too few steps can be possible if a chain of resonance decays has been
196  // removed. In this case, reject this event, since it will be handled in
197  // lower-multiplicity samples.
198  int nRequested = mergingHooksPtr->nRequested();
199 
200  // Store hard event cut information, reset veto information.
201  mergingHooksPtr->setHardProcessInfo(nSteps, tmsnow);
202  mergingHooksPtr->setEventVetoInfo(-1, -1.);
203 
204  if (nSteps < nRequested && allowReject) {
205  if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
206  if ( includeWGT) infoPtr->updateWeight(0.);
207  return -1;
208  }
209 
210  // Reset the minimal tms value, if necessary.
211  tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
212 
213  // Get random number to choose a path.
214  double RN = rndmPtr->flat();
215 
216  // Set dummy process scale.
217  newProcess.scale(0.0);
218  // Generate all histories.
219  History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
220  (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
221  trialPartonLevelPtr, coupSMPtr, true, true, true, true, 1.0, 0);
222 
223  // Project histories onto desired branches, e.g. only ordered paths.
224  FullHistory.projectOntoDesiredHistories();
225 
226  // Setup the selected path. Needed for
227  FullHistory.select(RN)->setSelectedChild();
228 
229  // Do not apply cut if the configuration could not be projected onto an
230  // underlying born configuration.
231  bool applyCut = allowReject
232  && nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
233 
234  // Enfore merging scale cut if the event did not pass the merging scale
235  // criterion.
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.);
243  return -1;
244  }
245 
246  // Check if more steps should be taken.
247  int nFinalP = 0;
248  int nFinalW = 0;
249  Event coreProcess = Event();
250  coreProcess.clear();
251  coreProcess.init( "(hard process-modified)", particleDataPtr );
252  coreProcess.clear();
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 )
257  nFinalP++;
258  if ( coreProcess[i].idAbs() == 24 )
259  nFinalW++;
260  }
261 
262  bool complete = (FullHistory.select(RN)->nClusterings() == nSteps) ||
263  ( mergingHooksPtr->doWeakClustering() && nFinalP == 2 && nFinalW == 0 );
264 
265  if ( !complete ) {
266  string message="Warning in Merging::mergeProcessCKKWL: No clusterings";
267  message+=" found. History incomplete.";
268  infoPtr->errorMsg(message);
269  }
270 
271  // Calculate CKKWL weight:
272  // Perform reweighting with Sudakov factors, save alpha_s ratios and
273  // PDF ratio weights.
274  wgt = FullHistory.weightTREE( trialPartonLevelPtr,
275  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
276  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
277 
278  // Event with production scales set for further (trial) showering
279  // and starting conditions for the shower.
280  FullHistory.getStartingConditions( RN, process );
281  // If necessary, reattach resonance decay products.
282  mergingHooksPtr->reattachResonanceDecays(process);
283 
284  // Allow to dampen histories in which the lowest multiplicity reclustered
285  // state does not pass the lowest multiplicity cut of the matrix element.
286  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
287  FullHistory.lowestMultProc(RN) );
288  // Save the weight of the event for histogramming. Only change the
289  // event weight after trial shower on the matrix element
290  // multiplicity event (= in doVetoStep).
291  wgt *= dampWeight;
292 
293  // Save the weight of the event for histogramming.
294  if (!includeWGT) mergingHooksPtr->setWeightCKKWL(wgt);
295 
296  // Update the event weight.
297  double norm = (abs(infoPtr->lhaStrategy()) == 4) ? 1/1e9 : 1.;
298  if ( includeWGT) infoPtr->updateWeight(infoPtr->weight()*wgt*norm);
299 
300  // Allow merging hooks to veto events from now on.
301  mergingHooksPtr->doIgnoreStep(false);
302 
303  // If no-emission probability is zero.
304  if ( allowReject && wgt == 0. ) return 0;
305 
306  // Done
307  return 1;
308 
309 }
310 
311 //--------------------------------------------------------------------------
312 
313 // Function to perform UMEPS merging on this event.
314 
315 int Merging::mergeProcessUMEPS( Event& process) {
316 
317  // Initialise which part of UMEPS merging is applied.
318  bool doUMEPSTree = settingsPtr->flag("Merging:doUMEPSTree");
319  bool doUMEPSSubt = settingsPtr->flag("Merging:doUMEPSSubt");
320  // Save number of looping steps
321  mergingHooksPtr->nReclusterSave = settingsPtr->mode("Merging:nRecluster");
322  int nRecluster = settingsPtr->mode("Merging:nRecluster");
323 
324  // Ensure that merging hooks does not remove emissions.
325  mergingHooksPtr->doIgnoreEmissions(true);
326  // For pp > h, allow cut on state, so that underlying processes
327  // can be clustered to gg > h
328  if ( mergingHooksPtr->getProcessString().compare("pp>h") == 0 )
329  mergingHooksPtr->allowCutOnRecState(true);
330  // For now, prefer construction of ordered histories.
331  mergingHooksPtr->orderHistories(true);
332 
333  // Ensure that merging weight is not counted twice.
334  bool includeWGT = mergingHooksPtr->includeWGTinXSEC();
335 
336  // Reset any incoming spins for W+-.
337  if (mergingHooksPtr->doWeakClustering())
338  for (int i = 0;i < process.size();++i)
339  process[i].pol(9);
340 
341  // Reset weights of the event.
342  double wgt = 1.;
343  mergingHooksPtr->setWeightCKKWL(1.);
344  mergingHooksPtr->muMI(-1.);
345 
346  // Prepare process record for merging. If Pythia has already decayed
347  // resonances used to define the hard process, remove resonance decay
348  // products.
349  Event newProcess( mergingHooksPtr->bareEvent( process, true) );
350  // Store candidates for the splitting V -> qqbar'.
351  mergingHooksPtr->storeHardProcessCandidates( newProcess );
352 
353  // Check if event passes the merging scale cut.
354  double tmsval = mergingHooksPtr->tms();
355  // Get merging scale in current event.
356  double tmsnow = mergingHooksPtr->tmsNow( newProcess );
357  // Calculate number of clustering steps.
358  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess, true);
359  int nRequested = mergingHooksPtr->nRequested();
360 
361  // Too few steps can be possible if a chain of resonance decays has been
362  // removed. In this case, reject this event, since it will be handled in
363  // lower-multiplicity samples.
364  if (nSteps < nRequested) {
365  if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
366  if ( includeWGT) infoPtr->updateWeight(0.);
367  return -1;
368  }
369 
370  // Reset the minimal tms value, if necessary.
371  tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
372 
373  // Get random number to choose a path.
374  double RN = rndmPtr->flat();
375  // Set dummy process scale.
376  newProcess.scale(0.0);
377  // Generate all histories.
378  History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
379  (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
380  trialPartonLevelPtr, coupSMPtr, true, true, true, true, 1.0, 0);
381  // Project histories onto desired branches, e.g. only ordered paths.
382  FullHistory.projectOntoDesiredHistories();
383 
384  // Do not apply cut if the configuration could not be projected onto an
385  // underlying born configuration.
386  bool applyCut = nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
387 
388  // Enfore merging scale cut if the event did not pass the merging scale
389  // criterion.
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.);
397  return -1;
398  }
399 
400  // Check reclustering steps to correctly apply MPI.
401  int nPerformed = 0;
402  if ( nSteps > 0 && doUMEPSSubt
403  && !FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster,
404  newProcess, nPerformed, false ) ) {
405  // Discard if the state could not be reclustered to a state above TMS.
406  if (!includeWGT) mergingHooksPtr->setWeightCKKWL(0.);
407  if ( includeWGT) infoPtr->updateWeight(0.);
408  return -1;
409  }
410 
411  mergingHooksPtr->nMinMPI(nSteps - nPerformed);
412 
413  // Calculate CKKWL weight:
414  // Perform reweighting with Sudakov factors, save alpha_s ratios and
415  // PDF ratio weights.
416  if ( doUMEPSTree ) {
417  wgt = FullHistory.weight_UMEPS_TREE( trialPartonLevelPtr,
418  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
419  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
420  } else {
421  wgt = FullHistory.weight_UMEPS_SUBT( trialPartonLevelPtr,
422  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
423  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN);
424  }
425 
426  // Event with production scales set for further (trial) showering
427  // and starting conditions for the shower.
428  if ( doUMEPSTree ) FullHistory.getStartingConditions( RN, process );
429  // Do reclustering (looping) steps.
430  else FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster, process,
431  nPerformed, true );
432 
433  // Allow to dampen histories in which the lowest multiplicity reclustered
434  // state does not pass the lowest multiplicity cut of the matrix element
435  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
436  FullHistory.lowestMultProc(RN) );
437  // Save the weight of the event for histogramming. Only change the
438  // event weight after trial shower on the matrix element
439  // multiplicity event (= in doVetoStep)
440  wgt *= dampWeight;
441 
442  // Save the weight of the event for histogramming.
443  if (!includeWGT) mergingHooksPtr->setWeightCKKWL(wgt);
444 
445  // Update the event weight.
446  double norm = (abs(infoPtr->lhaStrategy()) == 4) ? 1/1e9 : 1.;
447  if ( includeWGT) infoPtr->updateWeight(infoPtr->weight()*wgt*norm);
448 
449  // Set QCD 2->2 starting scale different from arbitrary scale in LHEF!
450  // --> Set to minimal mT of partons.
451  int nFinal = 0;
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 ) ) {
456  nFinal++;
457  muf = min( muf, abs(process[i].mT()) );
458  }
459 
460  // For pure QCD dijet events (only!), set the process scale to the
461  // transverse momentum of the outgoing partons.
462  // Calculate number of clustering steps.
463  int nStepsNew = mergingHooksPtr->getNumberOfClusteringSteps( process );
464  if ( nStepsNew == 0
465  && ( mergingHooksPtr->getProcessString().compare("pp>jj") == 0
466  || mergingHooksPtr->getProcessString().compare("pp>aj") == 0) )
467  process.scale(muf);
468 
469  // Reset hard process candidates (changed after clustering a parton).
470  mergingHooksPtr->storeHardProcessCandidates( process );
471  // If necessary, reattach resonance decay products.
472  mergingHooksPtr->reattachResonanceDecays(process);
473 
474  // Allow merging hooks to remove emissions from now on.
475  mergingHooksPtr->doIgnoreEmissions(false);
476 
477  // If no-emission probability is zero.
478  if ( wgt == 0. ) return 0;
479 
480  // Done
481  return 1;
482 
483 }
484 
485 //--------------------------------------------------------------------------
486 
487 // Function to perform NL3 NLO merging on this event.
488 
489 int Merging::mergeProcessNL3( Event& process) {
490 
491  // Initialise which part of NL3 merging is applied.
492  bool doNL3Tree = settingsPtr->flag("Merging:doNL3Tree");
493  bool doNL3Loop = settingsPtr->flag("Merging:doNL3Loop");
494  bool doNL3Subt = settingsPtr->flag("Merging:doNL3Subt");
495 
496  // Ensure that hooks (NL3 part) to not remove emissions.
497  mergingHooksPtr->doIgnoreEmissions(true);
498  // Ensure that hooks (CKKWL part) to not veto events in trial showers.
499  mergingHooksPtr->doIgnoreStep(true);
500  // For pp > h, allow cut on state, so that underlying processes
501  // can be clustered to gg > h
502  if ( mergingHooksPtr->getProcessString().compare("pp>h") == 0)
503  mergingHooksPtr->allowCutOnRecState(true);
504  // For now, prefer construction of ordered histories.
505  mergingHooksPtr->orderHistories(true);
506 
507  // Reset weight of the event
508  double wgt = 1.;
509  mergingHooksPtr->setWeightCKKWL(1.);
510  // Reset the O(alphaS)-term of the CKKW-L weight.
511  double wgtFIRST = 0.;
512  mergingHooksPtr->setWeightFIRST(0.);
513  mergingHooksPtr->muMI(-1.);
514 
515  // Prepare process record for merging. If Pythia has already decayed
516  // resonances used to define the hard process, remove resonance decay
517  // products.
518  Event newProcess( mergingHooksPtr->bareEvent( process, true) );
519  // Store candidates for the splitting V -> qqbar'
520  mergingHooksPtr->storeHardProcessCandidates( newProcess);
521 
522  // Check if event passes the merging scale cut.
523  double tmsval = mergingHooksPtr->tms();
524  // Get merging scale in current event.
525  double tmsnow = mergingHooksPtr->tmsNow( newProcess );
526  // Calculate number of clustering steps
527  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess, true);
528  int nRequested = mergingHooksPtr->nRequested();
529 
530  // Too few steps can be possible if a chain of resonance decays has been
531  // removed. In this case, reject this event, since it will be handled in
532  // lower-multiplicity samples.
533  if (nSteps < nRequested) {
534  mergingHooksPtr->setWeightCKKWL(0.);
535  mergingHooksPtr->setWeightFIRST(0.);
536  return -1;
537  }
538 
539  // Reset the minimal tms value, if necessary.
540  tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
541 
542  // Enfore merging scale cut if the event did not pass the merging scale
543  // criterion.
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.);
552  return -1;
553  }
554 
555  // Get random number to choose a path.
556  double RN = rndmPtr->flat();
557  // Set dummy process scale.
558  newProcess.scale(0.0);
559  // Generate all histories
560  History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
561  (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
562  trialPartonLevelPtr, coupSMPtr, true, true, true, true, 1.0, 0);
563  // Project histories onto desired branches, e.g. only ordered paths.
564  FullHistory.projectOntoDesiredHistories();
565 
566  // Discard states that cannot be projected unto a state with one less jet.
567  if ( nSteps > 0 && doNL3Subt
568  && FullHistory.select(RN)->nClusterings() == 0 ){
569  mergingHooksPtr->setWeightCKKWL(0.);
570  mergingHooksPtr->setWeightFIRST(0.);
571  return -1;
572  }
573 
574  // Potentially recluster real emission jets for powheg input containing
575  // "too many" jets, i.e. real-emission kinematics.
576  bool containsRealKin = nSteps > nRequested && nSteps > 0;
577 
578  // Perform one reclustering for real emission kinematics, then apply merging
579  // scale cut on underlying Born kinematics.
580  if ( containsRealKin ) {
581  Event dummy = Event();
582  // Initialise temporary output of reclustering.
583  dummy.clear();
584  dummy.init( "(hard process-modified)", particleDataPtr );
585  dummy.clear();
586  // Recluster once.
587  if ( !FullHistory.getClusteredEvent( RN, nSteps, dummy )) {
588  mergingHooksPtr->setWeightCKKWL(0.);
589  mergingHooksPtr->setWeightFIRST(0.);
590  return -1;
591  }
592  double tnowNew = mergingHooksPtr->tmsNow( dummy );
593  // Veto if underlying Born kinematics do not pass merging scale cut.
594  if ( enforceCutOnLHE && nSteps > 0 && nRequested > 0
595  && tnowNew < tmsval ) {
596  mergingHooksPtr->setWeightCKKWL(0.);
597  mergingHooksPtr->setWeightFIRST(0.);
598  return -1;
599  }
600  }
601 
602  // Remember number of jets, to include correct MPI no-emission probabilities.
603  if ( doNL3Subt || containsRealKin ) mergingHooksPtr->nMinMPI(nSteps - 1);
604  else mergingHooksPtr->nMinMPI(nSteps);
605 
606  // Calculate weight
607  // Do LO or first part of NLO tree-level reweighting
608  if( doNL3Tree ) {
609  // Perform reweighting with Sudakov factors, save as ratios and
610  // PDF ratio weights
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 ) {
615  // No reweighting, just set event scales properly and incorporate MPI
616  // no-emission probabilities.
617  wgt = FullHistory.weightLOOP( trialPartonLevelPtr, RN);
618  }
619 
620  // Event with production scales set for further (trial) showering
621  // and starting conditions for the shower
622  if ( !doNL3Subt && !containsRealKin )
623  FullHistory.getStartingConditions(RN, process);
624  // For sutraction of nSteps-additional resolved partons from
625  // the nSteps-1 parton phase space, recluster the last parton
626  // in nSteps-parton events, and sutract later
627  else {
628  // Function to return the reclustered event
629  if ( !FullHistory.getClusteredEvent( RN, nSteps, process )) {
630  mergingHooksPtr->setWeightCKKWL(0.);
631  mergingHooksPtr->setWeightFIRST(0.);
632  return -1;
633  }
634  }
635 
636  // Allow to dampen histories in which the lowest multiplicity reclustered
637  // state does not pass the lowest multiplicity cut of the matrix element
638  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
639  FullHistory.lowestMultProc(RN) );
640  // Save the weight of the event for histogramming. Only change the
641  // event weight after trial shower on the matrix element
642  // multiplicity event (= in doVetoStep)
643  wgt *= dampWeight;
644 
645  // For tree level samples in NL3, rescale with k-Factor
646  if (doNL3Tree ){
647  // Find k-factor
648  double kFactor = 1.;
649  if( nSteps > mergingHooksPtr->nMaxJetsNLO() )
650  kFactor = mergingHooksPtr->kFactor( mergingHooksPtr->nMaxJetsNLO() );
651  else kFactor = mergingHooksPtr->kFactor(nSteps);
652  // For NLO merging, rescale CKKW-L weight with k-factor
653  wgt *= kFactor;
654  }
655 
656  // Save the weight of the event for histogramming
657  mergingHooksPtr->setWeightCKKWL(wgt);
658 
659  // Check if we need to subtract the O(\alpha_s)-term. If the number
660  // of additional partons is larger than the number of jets for
661  // which loop matrix elements are available, do standard CKKW-L
662  bool doOASTree = doNL3Tree && nSteps <= mergingHooksPtr->nMaxJetsNLO();
663 
664  // Now begin NLO part for tree-level events
665  if ( doOASTree ) {
666  // Calculate the O(\alpha_s)-term of the CKKWL weight
667  wgtFIRST = FullHistory.weightFIRST( trialPartonLevelPtr,
668  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
669  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(), RN,
670  rndmPtr );
671  // If necessary, also dampen the O(\alpha_s)-term
672  wgtFIRST *= dampWeight;
673  // Set the subtractive weight to the value calculated so far
674  mergingHooksPtr->setWeightFIRST(wgtFIRST);
675  // Subtract the O(\alpha_s)-term from the CKKW-L weight
676  // If PDF contributions have not been included, subtract these later
677  wgt = wgt - wgtFIRST;
678  }
679 
680  // Set qcd 2->2 starting scale different from arbirtrary scale in LHEF!
681  // --> Set to pT of partons
682  double pT = 0.;
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));
686  break;
687  }
688  // For pure QCD dijet events (only!), set the process scale to the
689  // transverse momentum of the outgoing partons.
690  if ( nSteps == 0
691  && mergingHooksPtr->getProcessString().compare("pp>jj") == 0)
692  process.scale(pT);
693 
694  // Reset hard process candidates (changed after clustering a parton).
695  mergingHooksPtr->storeHardProcessCandidates( process );
696  // If necessary, reattach resonance decay products.
697  mergingHooksPtr->reattachResonanceDecays(process);
698 
699  // Allow merging hooks (NL3 part) to remove emissions from now on.
700  mergingHooksPtr->doIgnoreEmissions(false);
701  // Allow merging hooks (CKKWL part) to veto events from now on.
702  mergingHooksPtr->doIgnoreStep(false);
703 
704  // Done
705  return 1;
706 
707 }
708 
709 //--------------------------------------------------------------------------
710 
711 // Function to perform UNLOPS merging on this event.
712 
713 int Merging::mergeProcessUNLOPS( Event& process) {
714 
715  // Initialise which part of UNLOPS merging is applied.
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");
721  // Save number of looping steps
722  mergingHooksPtr->nReclusterSave = settingsPtr->mode("Merging:nRecluster");
723  int nRecluster = settingsPtr->mode("Merging:nRecluster");
724 
725  // Ensure that merging hooks to not remove emissions
726  mergingHooksPtr->doIgnoreEmissions(true);
727  // For now, prefer construction of ordered histories.
728  mergingHooksPtr->orderHistories(true);
729  // For pp > h, allow cut on state, so that underlying processes
730  // can be clustered to gg > h
731  if ( mergingHooksPtr->getProcessString().compare("pp>h") == 0)
732  mergingHooksPtr->allowCutOnRecState(true);
733 
734  // Reset weight of the event.
735  double wgt = 1.;
736  mergingHooksPtr->setWeightCKKWL(1.);
737  // Reset the O(alphaS)-term of the UMEPS weight.
738  double wgtFIRST = 0.;
739  mergingHooksPtr->setWeightFIRST(0.);
740  mergingHooksPtr->muMI(-1.);
741 
742  // Prepare process record for merging. If Pythia has already decayed
743  // resonances used to define the hard process, remove resonance decay
744  // products.
745  Event newProcess( mergingHooksPtr->bareEvent( process, true) );
746  // Store candidates for the splitting V -> qqbar'
747  mergingHooksPtr->storeHardProcessCandidates( newProcess );
748 
749  // Check if event passes the merging scale cut.
750  double tmsval = mergingHooksPtr->tms();
751  // Get merging scale in current event.
752  double tmsnow = mergingHooksPtr->tmsNow( newProcess );
753  // Calculate number of clustering steps
754  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess, true);
755  int nRequested = mergingHooksPtr->nRequested();
756 
757  // Too few steps can be possible if a chain of resonance decays has been
758  // removed. In this case, reject this event, since it will be handled in
759  // lower-multiplicity samples.
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.);
766  return -1;
767  }
768 
769  // Reset the minimal tms value, if necessary.
770  tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
771 
772  // Get random number to choose a path.
773  double RN = rndmPtr->flat();
774  // Set dummy process scale.
775  newProcess.scale(0.0);
776  // Generate all histories
777  History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
778  (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
779  trialPartonLevelPtr, coupSMPtr, true, true, true, true, 1.0, 0);
780  // Project histories onto desired branches, e.g. only ordered paths.
781  FullHistory.projectOntoDesiredHistories();
782 
783  // Do not apply cut if the configuration could not be projected onto an
784  // underlying born configuration.
785  bool applyCut = nSteps > 0 && FullHistory.select(RN)->nClusterings() > 0;
786 
787  // Enfore merging scale cut if the event did not pass the merging scale
788  // criterion.
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.);
797  return -1;
798  }
799 
800  // Potentially recluster real emission jets for powheg input containing
801  // "too many" jets, i.e. real-emission kinematics.
802  bool containsRealKin = nSteps > nRequested && nSteps > 0;
803  if ( containsRealKin ) nRecluster += nSteps - nRequested;
804 
805  // Remove real emission events without underlying Born configuration from
806  // the loop sample, since such states will be taken care of by tree-level
807  // samples.
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.);
814  return -1;
815  }
816 
817  // Discard if the state could not be reclustered to any state above TMS.
818  int nPerformed = 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.);
825  return -1;
826  }
827 
828  // Check reclustering steps to correctly apply MPI.
829  mergingHooksPtr->nMinMPI(nSteps - nPerformed);
830 
831  // Perform one reclustering for real emission kinematics, then apply
832  // merging scale cut on underlying Born kinematics.
833  if ( containsRealKin ) {
834  Event dummy = Event();
835  // Initialise temporary output of reclustering.
836  dummy.clear();
837  dummy.init( "(hard process-modified)", particleDataPtr );
838  dummy.clear();
839  // Recluster once.
840  FullHistory.getClusteredEvent( RN, nSteps, dummy );
841  double tnowNew = mergingHooksPtr->tmsNow( dummy );
842  // Veto if underlying Born kinematics do not pass merging scale cut.
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.);
850  return -1;
851  }
852  }
853 
854  // New UNLOPS strategy based on UN2LOPS.
855  bool doUNLOPS2 = false;
856  int depth = (!doUNLOPS2) ? -1 : ( (containsRealKin) ? nSteps-1 : nSteps);
857 
858  // Calculate weights.
859  // Do LO or first part of NLO tree-level reweighting
860  if( doUNLOPSTree ) {
861  // Perform reweighting with Sudakov factors, save as ratios and
862  // PDF ratio weights
863  wgt = FullHistory.weight_UNLOPS_TREE( trialPartonLevelPtr,
864  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
865  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
866  RN, depth);
867  } else if( doUNLOPSLoop ) {
868  // Set event scales properly, reweight for new UNLOPS
869  wgt = FullHistory.weight_UNLOPS_LOOP( trialPartonLevelPtr,
870  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
871  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
872  RN, depth);
873  } else if( doUNLOPSSubtNLO ) {
874  // Set event scales properly, reweight for new UNLOPS
875  wgt = FullHistory.weight_UNLOPS_SUBTNLO( trialPartonLevelPtr,
876  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
877  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
878  RN, depth);
879  } else if( doUNLOPSSubt ) {
880  // Perform reweighting with Sudakov factors, save as ratios and
881  // PDF ratio weights
882  wgt = FullHistory.weight_UNLOPS_SUBT( trialPartonLevelPtr,
883  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(),
884  mergingHooksPtr->AlphaEM_FSR(), mergingHooksPtr->AlphaEM_ISR(),
885  RN, depth);
886  }
887 
888  // Event with production scales set for further (trial) showering
889  // and starting conditions for the shower.
890  if (!doUNLOPSSubt && !doUNLOPSSubtNLO && !containsRealKin )
891  FullHistory.getStartingConditions(RN, process);
892  // Do reclustering (looping) steps.
893  else FullHistory.getFirstClusteredEventAboveTMS( RN, nRecluster, process,
894  nPerformed, true );
895 
896  // Allow to dampen histories in which the lowest multiplicity reclustered
897  // state does not pass the lowest multiplicity cut of the matrix element
898  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
899  FullHistory.lowestMultProc(RN) );
900  // Save the weight of the event for histogramming. Only change the
901  // event weight after trial shower on the matrix element
902  // multiplicity event (= in doVetoStep)
903  wgt *= dampWeight;
904 
905  // For tree-level or subtractive sammples, rescale with k-Factor
906  if ( doUNLOPSTree || doUNLOPSSubt ){
907  // Find k-factor
908  double kFactor = 1.;
909  if ( nSteps > mergingHooksPtr->nMaxJetsNLO() )
910  kFactor = mergingHooksPtr->kFactor( mergingHooksPtr->nMaxJetsNLO() );
911  else kFactor = mergingHooksPtr->kFactor(nSteps);
912  // For NLO merging, rescale CKKW-L weight with k-factor
913  wgt *= (nRecluster == 2 && nloTilde) ? 1. : kFactor;
914  }
915 
916  // Save the weight of the event for histogramming
917  mergingHooksPtr->setWeightCKKWL(wgt);
918 
919  // Check if we need to subtract the O(\alpha_s)-term. If the number
920  // of additional partons is larger than the number of jets for
921  // which loop matrix elements are available, do standard UMEPS.
922  int nMaxNLO = mergingHooksPtr->nMaxJetsNLO();
923  bool doOASTree = doUNLOPSTree && nSteps <= nMaxNLO;
924  bool doOASSubt = doUNLOPSSubt && nSteps <= nMaxNLO+1 && nSteps > 0;
925 
926  // Now begin NLO part for tree-level events
927  if ( doOASTree || doOASSubt ) {
928 
929  // Decide on which order to expand to.
930  int order = ( nSteps > 0 && nSteps <= nMaxNLO) ? 1 : -1;
931 
932  // Exclusive inputs:
933  // Subtract only the O(\alpha_s^{n+0})-term from the tree-level
934  // subtraction, if we're at the highest NLO multiplicity (nMaxNLO).
935  if ( nloTilde && doUNLOPSSubt && nRecluster == 1
936  && nSteps == nMaxNLO+1 ) order = 0;
937 
938  // Exclusive inputs:
939  // Do not remove the O(as)-term if the number of reclusterings
940  // exceeds the number of NLO jets, or if more clusterings have
941  // been performed.
942  if (nloTilde && doUNLOPSSubt && ( nSteps > nMaxNLO+1
943  || (nSteps == nMaxNLO+1 && nPerformed != nRecluster) ))
944  order = -1;
945 
946  // Calculate terms in expansion of the CKKW-L weight.
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 );
951 
952  // Exclusive inputs:
953  // Subtract the O(\alpha_s^{n+1})-term from the tree-level
954  // subtraction, not the O(\alpha_s^{n+0})-terms.
955  if ( nloTilde && doUNLOPSSubt && nRecluster == 1
956  && nPerformed == nRecluster && nSteps <= nMaxNLO )
957  wgtFIRST += 1.;
958 
959  // If necessary, also dampen the O(\alpha_s)-term
960  wgtFIRST *= dampWeight;
961 
962  // Set the subtractive weight to the value calculated so far
963  mergingHooksPtr->setWeightFIRST(wgtFIRST);
964  // Subtract the O(\alpha_s)-term from the CKKW-L weight
965  // If PDF contributions have not been included, subtract these later
966  // New UNLOPS based on UN2LOPS.
967  if (doUNLOPS2 && order > -1) wgt = -wgt*(wgtFIRST-1.);
968  else if (order > -1) wgt = wgt - wgtFIRST;
969 
970  }
971 
972  // Set QCD 2->2 starting scale different from arbitrary scale in LHEF!
973  // --> Set to minimal mT of partons.
974  int nFinal = 0;
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 ) ) {
979  nFinal++;
980  muf = min( muf, abs(process[i].mT()) );
981  }
982  // For pure QCD dijet events (only!), set the process scale to the
983  // transverse momentum of the outgoing partons.
984  if ( nSteps == 0 && nFinal == 2
985  && ( mergingHooksPtr->getProcessString().compare("pp>jj") == 0
986  || mergingHooksPtr->getProcessString().compare("pp>aj") == 0) )
987  process.scale(muf);
988 
989  // Reset hard process candidates (changed after clustering a parton).
990  mergingHooksPtr->storeHardProcessCandidates( process );
991 
992  // Check if resonance structure has been changed
993  // (e.g. because of clustering W/Z/gluino)
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());
1002  // Compare old and new resonances
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;
1007  break;
1008  }
1009  bool hasNewResonances = (newResonance.size() != oldResonance.size());
1010  for ( int i=0; i < int(oldResonance.size()); ++i )
1011  hasNewResonances = (oldResonance[i] != 99);
1012 
1013  // If necessary, reattach resonance decay products.
1014  if (!hasNewResonances) mergingHooksPtr->reattachResonanceDecays(process);
1015 
1016  // Allow merging hooks to remove emissions from now on.
1017  mergingHooksPtr->doIgnoreEmissions(false);
1018 
1019  // If no-emission probability is zero.
1020  if ( wgt == 0. ) return 0;
1021 
1022  // If the resonance structure of the process has changed due to reclustering,
1023  // redo the resonance decays in Pythia::next()
1024  if (hasNewResonances) return 2;
1025 
1026  // Done
1027  return 1;
1028 
1029 }
1030 
1031 //--------------------------------------------------------------------------
1032 
1033 // Function to apply the merging scale cut on an input event.
1034 
1035 bool Merging::cutOnProcess( Event& process) {
1036 
1037  // Save number of looping steps
1038  mergingHooksPtr->nReclusterSave = settingsPtr->mode("Merging:nRecluster");
1039 
1040  // For now, prefer construction of ordered histories.
1041  mergingHooksPtr->orderHistories(true);
1042  // For pp > h, allow cut on state, so that underlying processes
1043  // can be clustered to gg > h
1044  if ( mergingHooksPtr->getProcessString().compare("pp>h") == 0)
1045  mergingHooksPtr->allowCutOnRecState(true);
1046 
1047  // Reset any incoming spins for W+-.
1048  if (mergingHooksPtr->doWeakClustering())
1049  for (int i = 0;i < process.size();++i)
1050  process[i].pol(9);
1051 
1052  // Prepare process record for merging. If Pythia has already decayed
1053  // resonances used to define the hard process, remove resonance decay
1054  // products.
1055  Event newProcess( mergingHooksPtr->bareEvent( process, true) );
1056  // Store candidates for the splitting V -> qqbar'
1057  mergingHooksPtr->storeHardProcessCandidates( newProcess );
1058 
1059  // Check if event passes the merging scale cut.
1060  double tmsval = mergingHooksPtr->tms();
1061  // Get merging scale in current event.
1062  double tmsnow = mergingHooksPtr->tmsNow( newProcess );
1063  // Calculate number of clustering steps
1064  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( newProcess, true);
1065 
1066  // Too few steps can be possible if a chain of resonance decays has been
1067  // removed. In this case, reject this event, since it will be handled in
1068  // lower-multiplicity samples.
1069  int nRequested = mergingHooksPtr->nRequested();
1070  if (nSteps < nRequested) return true;
1071 
1072  // Reset the minimal tms value, if necessary.
1073  tmsNowMin = (nSteps == 0) ? 0. : min(tmsNowMin, tmsnow);
1074 
1075  // Potentially recluster real emission jets for powheg input containing
1076  // "too many" jets, i.e. real-emission kinematics.
1077  bool containsRealKin = nSteps > nRequested && nSteps > 0;
1078 
1079  // Get random number to choose a path.
1080  double RN = rndmPtr->flat();
1081  // Set dummy process scale.
1082  newProcess.scale(0.0);
1083  // Generate all histories
1084  History FullHistory( nSteps, 0.0, newProcess, Clustering(), mergingHooksPtr,
1085  (*beamAPtr), (*beamBPtr), particleDataPtr, infoPtr,
1086  trialPartonLevelPtr, coupSMPtr, true, true, true, true, 1.0, 0);
1087  // Project histories onto desired branches, e.g. only ordered paths.
1088  FullHistory.projectOntoDesiredHistories();
1089 
1090  // Remove real emission events without underlying Born configuration from
1091  // the loop sample, since such states will be taken care of by tree-level
1092  // samples.
1093  bool allowIncompleteReal =
1094  settingsPtr->flag("Merging:allowIncompleteHistoriesInReal");
1095  if ( containsRealKin && !allowIncompleteReal
1096  && FullHistory.select(RN)->nClusterings() == 0 )
1097  return true;
1098 
1099  // Cut if no history passes the cut on the lowest-multiplicity state.
1100  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
1101  FullHistory.lowestMultProc(RN) );
1102  if ( dampWeight == 0. ) return true;
1103 
1104  // Do not apply cut if the configuration could not be projected onto an
1105  // underlying born configuration.
1106  if ( nSteps > 0 && FullHistory.select(RN)->nClusterings() == 0 )
1107  return false;
1108 
1109  // Now enfore merging scale cut if the event did not pass the merging scale
1110  // criterion.
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);
1115  return true;
1116  }
1117 
1118  // Check if more steps should be taken.
1119  int nFinalP = 0;
1120  int nFinalW = 0;
1121  Event coreProcess = Event();
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 )
1129  nFinalP++;
1130  if ( coreProcess[i].idAbs() == 24 )
1131  nFinalW++;
1132  }
1133 
1134  bool complete = (FullHistory.select(RN)->nClusterings() == nSteps) ||
1135  ( mergingHooksPtr->doWeakClustering() && nFinalP == 2 && nFinalW == 0 );
1136 
1137  if ( !complete ) {
1138  string message="Warning in Merging::cutOnProcess: No clusterings";
1139  message+=" found. History incomplete.";
1140  infoPtr->errorMsg(message);
1141  }
1142 
1143  // Done if no real-emission jets are present.
1144  if ( !containsRealKin ) return false;
1145 
1146  // Now cut on events that contain an additional real-emission jet.
1147  // Perform one reclustering for real emission kinematics, then apply merging
1148  // scale cut on underlying Born kinematics.
1149  if ( containsRealKin ) {
1150  Event dummy = Event();
1151  // Initialise temporary output of reclustering.
1152  dummy.clear();
1153  dummy.init( "(hard process-modified)", particleDataPtr );
1154  dummy.clear();
1155  // Recluster once.
1156  FullHistory.getClusteredEvent( RN, nSteps, dummy );
1157  double tnowNew = mergingHooksPtr->tmsNow( dummy );
1158  // Veto if underlying Born kinematics do not pass merging scale cut.
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);
1163  return true;
1164  }
1165  }
1166 
1167  // Done if only interested in cross section estimate after cuts.
1168  return false;
1169 
1170 }
1171 
1172 //==========================================================================
1173 
1174 } // end namespace Pythia8
Definition: AgUStep.h:26