StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
JetMatching.h
1 // JetMatching.h 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 // Authors: Richard Corke (implementation of MLM matching as
7 // in Alpgen for Alpgen input)
8 // and Stephen Mrenna (implementation of MLM-style matching as
9 // in Madgraph for Alpgen or Madgraph 5 input.)
10 // and Simon de Visscher, Stefan Prestel (implementation of shower-kT
11 // MLM-style matching and flavour treatment for Madgraph input)
12 // and Stefan Prestel (FxFx NLO jet matching with aMC@NLO.)
13 // This file provides the classes to perform MLM matching of
14 // Alpgen or MadGraph 5 input.
15 // Example usage is shown in main32.cc, and further details
16 // can be found in the 'Jet Matching Style' manual page.
17 
18 #ifndef Pythia8_JetMatching_H
19 #define Pythia8_JetMatching_H
20 
21 // Includes
22 #include "Pythia8/Pythia.h"
23 #include "Pythia8Plugins/GeneratorInput.h"
24 
25 namespace Pythia8 {
26 
27 //==========================================================================
28 
29 class HJSlowJet: public SlowJet {
30 
31  public:
32  HJSlowJet(int powerIn, double Rin, double pTjetMinIn = 0.,
33  double etaMaxIn = 25., int selectIn = 1, int massSetIn = 2,
34  SlowJetHook* sjHookPtrIn = 0, bool useFJcoreIn = false,
35  bool useStandardRin = true) :
36  SlowJet(powerIn, Rin, pTjetMinIn, etaMaxIn, selectIn, massSetIn,
37  sjHookPtrIn, useFJcoreIn, useStandardRin) {}
38 
39  protected:
40 
41  void findNext();
42 
43 };
44 
45 //--------------------------------------------------------------------------
46 
47 // Find next cluster pair to join.
48 
49 void HJSlowJet::findNext() {
50 
51  // Find smallest of diB, dij.
52  if (clSize > 0) {
53  iMin = 0;
54  jMin = -1;
55  dMin = 1.0/TINY;
56  // Remove the possibility of choosing a beam clustering
57  for (int i = 1; i < clSize; ++i) {
58  for (int j = 0; j < i; ++j) {
59  if (dij[i*(i-1)/2 + j] < dMin) {
60  iMin = i;
61  jMin = j;
62  dMin = dij[i*(i-1)/2 + j];
63  }
64  }
65  }
66 
67  // If no clusters left then instead default values.
68  } else {
69  iMin = -1;
70  jMin = -1;
71  dMin = 0.;
72  }
73 
74 }
75 
76 //==========================================================================
77 
78 // Declaration of main JetMatching class to perform MLM matching.
79 // Note that it is defined with virtual inheritance, so that it can
80 // be combined with other UserHooks classes, see e.g. main33.cc.
81 
82 class JetMatching : virtual public UserHooks {
83 
84 public:
85 
86  // Constructor and destructor
87  JetMatching() : cellJet(NULL), slowJet(NULL), slowJetHard(NULL),
88  hjSlowJet(NULL) {}
89  ~JetMatching() {
90  if (cellJet) delete cellJet;
91  if (slowJet) delete slowJet;
92  if (slowJetHard) delete slowJetHard;
93  if (hjSlowJet) delete hjSlowJet;
94  }
95 
96  // Initialisation
97  virtual bool initAfterBeams() = 0;
98 
99  // Process level vetos
100  bool canVetoProcessLevel() { return doMerge; }
101  bool doVetoProcessLevel(Event& process) {
102  eventProcessOrig = process;
103  return false;
104  }
105 
106  // Parton level vetos (before beam remnants and resonance decays)
107  bool canVetoPartonLevelEarly() { return doMerge; }
108  bool doVetoPartonLevelEarly(const Event& event);
109 
110  // Shower step vetoes (after the first emission, for Shower-kT scheme)
111  int numberVetoStep() {return 1;}
112  bool canVetoStep() { return false; }
113  bool doVetoStep(int, int, int, const Event& ) { return false; }
114 
115 protected:
116 
117  // Constants to be changed for debug printout or extra checks.
118  static const bool MATCHINGDEBUG, MATCHINGCHECK;
119 
120  // Different steps of the matching algorithm.
121  virtual void sortIncomingProcess(const Event &)=0;
122  virtual void jetAlgorithmInput(const Event &, int)=0;
123  virtual void runJetAlgorithm()=0;
124  virtual bool matchPartonsToJets(int)=0;
125  virtual int matchPartonsToJetsLight()=0;
126  virtual int matchPartonsToJetsHeavy()=0;
127 
128  enum vetoStatus { NONE, LESS_JETS, MORE_JETS, HARD_JET,
129  UNMATCHED_PARTON, INCLUSIVE_VETO };
130  enum partonTypes { ID_CHARM=4, ID_BOT=5, ID_TOP=6, ID_LEPMIN=11,
131  ID_LEPMAX=16, ID_GLUON=21, ID_PHOTON=22 };
132 
133  // Master switch for merging
134  bool doMerge;
135  // Switch for merging in the shower-kT scheme. Needed here because
136  // the scheme uses different UserHooks functionality.
137  bool doShowerKt;
138 
139  // Maximum and current number of jets
140  int nJetMax, nJet;
141 
142  // Jet algorithm parameters
143  int jetAlgorithm;
144  double eTjetMin, coneRadius, etaJetMax, etaJetMaxAlgo;
145 
146  // Internal jet algorithms
147  CellJet* cellJet;
148  SlowJet* slowJet;
149  SlowJet* slowJetHard;
150  HJSlowJet* hjSlowJet;
151 
152  // SlowJet specific
153  int slowJetPower;
154 
155  // Event records to store original incoming process, final-state of the
156  // incoming process and what will be passed to the jet algorithm.
157  // Not completely necessary to store all steps, but makes tracking the
158  // steps of the algorithm a lot easier.
159  Event eventProcessOrig, eventProcess, workEventJet;
160 
161  // Sort final-state of incoming process into light/heavy jets and 'other'
162  vector<int> typeIdx[3];
163  set<int> typeSet[3];
164 
165  // Momenta output of jet algorithm (to provide same output regardless of
166  // the selected jet algorithm)
167  vector<Vec4> jetMomenta;
168 
169  // CellJet specific
170  int nEta, nPhi;
171  double eTseed, eTthreshold;
172 
173  // Merging procedure parameters
174  int jetAllow, jetMatch, exclusiveMode;
175  double coneMatchLight, coneRadiusHeavy, coneMatchHeavy;
176  bool exclusive;
177 
178  // Store the minimum eT/pT of matched light jets
179  double eTpTlightMin;
180 
181 };
182 
183 //==========================================================================
184 
185 // Declaration of main UserHooks class to perform Alpgen matching.
186 
187 class JetMatchingAlpgen : virtual public JetMatching {
188 
189 public:
190 
191  // Constructor and destructor
192  JetMatchingAlpgen() { }
193  ~JetMatchingAlpgen() { }
194 
195  // Initialisation
196  bool initAfterBeams();
197 
198 private:
199 
200  // Different steps of the matching algorithm.
201  void sortIncomingProcess(const Event &);
202  void jetAlgorithmInput(const Event &, int);
203  void runJetAlgorithm();
204  bool matchPartonsToJets(int);
205  int matchPartonsToJetsLight();
206  int matchPartonsToJetsHeavy();
207 
208  // Sorting utility
209  void sortTypeIdx(vector < int > &vecIn);
210 
211  // Constants
212  static const double GHOSTENERGY, ZEROTHRESHOLD;
213 
214 };
215 
216 //==========================================================================
217 
218 // Declaration of main UserHooks class to perform Madgraph matching.
219 
220 class JetMatchingMadgraph : virtual public JetMatching {
221 
222 public:
223 
224  // Constructor and destructor
225  JetMatchingMadgraph() : slowJetDJR(NULL) { }
226  ~JetMatchingMadgraph() { if (slowJetDJR) delete slowJetDJR; }
227 
228  // Initialisation
229  bool initAfterBeams();
230 
231  // Process level vetos
232  bool canVetoProcessLevel() { return doMerge; }
233  bool doVetoProcessLevel(Event& process);
234 
235  // Shower step vetoes (after the first emission, for Shower-kT scheme)
236  int numberVetoStep() {return 1;}
237  bool canVetoStep() { return doShowerKt; }
238  bool doVetoStep(int, int, int, const Event& );
239 
240  // Jet algorithm to access the jet separations in the cleaned event
241  // after showering.
242  SlowJet* slowJetDJR;
243  // Functions to return the jet clustering scales and number of ME partons.
244  // These are useful to investigate the matching systematics.
245  vector<double> getDJR() { return DJR;}
246  pair<int,int> nMEpartons() { return nMEpartonsSave;}
247 
248  // For systematic variations of the jet matching parameters, it is helpful
249  // to decouple the jet matching veto from the internal book-keeping. The
250  // veto can then be applied in hindsight by an expert user. The functions
251  // below return all the necessary information to do this.
252  Event getWorkEventJet() { return workEventJetSave; }
253  Event getProcessSubset() { return processSubsetSave; }
254  bool getExclusive() { return exclusive; }
255  double getPTfirst() { return pTfirstSave; }
256 
257 protected:
258 
259  // Stored values of all inputs necessary to perform the jet matching, as
260  // needed when the veto is applied externally.
261  Event processSubsetSave;
262  Event workEventJetSave;
263  double pTfirstSave;
264 
265  // Save if code should apply the veto, or simply store the things necessary
266  // to perform the veto externally.
267  bool performVeto;
268 
269  // Different steps of the matching algorithm.
270  void sortIncomingProcess(const Event &);
271  void jetAlgorithmInput(const Event &, int);
272  void runJetAlgorithm();
273  bool matchPartonsToJets(int);
274  int matchPartonsToJetsLight();
275  int matchPartonsToJetsHeavy();
276  int matchPartonsToJetsOther();
277  bool doShowerKtVeto(double pTfirst);
278 
279  // Functions to clear and set the jet clustering scales.
280  void clearDJR() { DJR.resize(0);}
281  void setDJR( const Event& event);
282  // Functions to clear and set the jet clustering scales.
283  void clear_nMEpartons() { nMEpartonsSave.first = nMEpartonsSave.second =-1;}
284  void set_nMEpartons( const int nOrig, const int nMatch) {
285  clear_nMEpartons();
286  nMEpartonsSave.first = nOrig;
287  nMEpartonsSave.second = nMatch;
288  };
289 
290  // Variables.
291  vector<int> origTypeIdx[3];
292  int nQmatch;
293  double qCut, qCutSq, clFact;
294  bool doFxFx;
295  int nPartonsNow;
296  double qCutME, qCutMESq;
297 
298  // Vector to store the jet clustering scales.
299  vector<double> DJR;
300  // Pair of integers giving the number of ME partons read from LHEF and used
301  // in the matching (can be different if some partons should not be matched)
302  pair<int,int> nMEpartonsSave;
303 
304  // Function to get the current number of partons in the Born state, as
305  // read from LHE.
306  int npNLO();
307 
308 };
309 
310 //==========================================================================
311 
312 // Main implementation of JetMatching class.
313 // This may be split out to a separate C++ file if desired,
314 // but currently included here for ease of use.
315 
316 //--------------------------------------------------------------------------
317 
318 // Constants to be changed for debug printout or extra checks.
319 const bool JetMatching::MATCHINGDEBUG = false;
320 const bool JetMatching::MATCHINGCHECK = false;
321 
322 //--------------------------------------------------------------------------
323 
324 // Early parton level veto (before beam remnants and resonance showers)
325 
326 inline bool JetMatching::doVetoPartonLevelEarly(const Event& event) {
327 
328  // 1) Sort the original incoming process. After this step is performed,
329  // the following assignments have been made:
330  // eventProcessOrig - the original incoming process
331  // eventProcess - the final-state of the incoming process with
332  // resonance decays removed (and resonances
333  // themselves now with positive status code)
334  // typeIdx[0/1/2] - Indices into 'eventProcess' of
335  // light jets/heavy jets/other
336  // typeSet[0/1/2] - Indices into 'event' of light jets/heavy jets/other
337  // workEvent - partons from the hardest subsystem + ISR + FSR only
338  sortIncomingProcess(event);
339 
340  // For the shower-kT scheme, do not perform any veto here, as any vetoing
341  // will already have taken place in doVetoStep.
342  if ( doShowerKt ) return false;
343 
344  // Debug printout.
345  if (MATCHINGDEBUG) {
346  // Begin
347  cout << endl << "-------- Begin Madgraph Debug --------" << endl;
348  // Original incoming process
349  cout << endl << "Original incoming process:";
350  eventProcessOrig.list();
351  // Final-state of original incoming process
352  cout << endl << "Final-state incoming process:";
353  eventProcess.list();
354  // List categories of sorted particles
355  for (size_t i = 0; i < typeIdx[0].size(); i++)
356  cout << ((i == 0) ? "Light jets: " : ", ") << setw(3) << typeIdx[0][i];
357  if( typeIdx[0].size()== 0 )
358  cout << "Light jets: None";
359 
360  for (size_t i = 0; i < typeIdx[1].size(); i++)
361  cout << ((i == 0) ? "\nHeavy jets: " : ", ") << setw(3) << typeIdx[1][i];
362  for (size_t i = 0; i < typeIdx[2].size(); i++)
363  cout << ((i == 0) ? "\nOther: " : ", ") << setw(3) << typeIdx[2][i];
364  // Full event at this stage
365  cout << endl << endl << "Event:";
366  event.list();
367  // Work event (partons from hardest subsystem + ISR + FSR)
368  cout << endl << "Work event:";
369  workEvent.list();
370  }
371 
372  // 2) Light/heavy jets: iType = 0 (light jets), 1 (heavy jets)
373  int iTypeEnd = (typeIdx[2].empty()) ? 2 : 3;
374  for (int iType = 0; iType < iTypeEnd; iType++) {
375 
376  // 2a) Find particles which will be passed from the jet algorithm.
377  // Input from 'workEvent' and output in 'workEventJet'.
378  jetAlgorithmInput(event, iType);
379 
380  // Debug printout.
381  if (MATCHINGDEBUG) {
382  // Jet algorithm event
383  cout << endl << "Jet algorithm event (iType = " << iType << "):";
384  workEventJet.list();
385  }
386 
387  // 2b) Run jet algorithm on 'workEventJet'.
388  // Output is stored in jetMomenta.
389  runJetAlgorithm();
390 
391  // 2c) Match partons to jets and decide if veto is necessary
392  if (matchPartonsToJets(iType) == true) {
393  // Debug printout.
394  if (MATCHINGDEBUG) {
395  cout << endl << "Event vetoed" << endl
396  << "---------- End MLM Debug ----------" << endl;
397  }
398  return true;
399  }
400  }
401 
402  // Debug printout.
403  if (MATCHINGDEBUG) {
404  cout << endl << "Event accepted" << endl
405  << "---------- End MLM Debug ----------" << endl;
406  }
407 
408  // If we reached here, then no veto
409  return false;
410 
411 }
412 
413 //==========================================================================
414 
415 // Main implementation of Alpgen UserHooks class.
416 // This may be split out to a separate C++ file if desired,
417 // but currently included here for ease of use.
418 
419 //--------------------------------------------------------------------------
420 
421 // Constants: could be changed here if desired, but normally should not.
422 // These are of technical nature, as described for each.
423 
424 // The energy of ghost particles. For technical reasons, this cannot be
425 // set arbitrarily low, see 'Particle::TINY' in 'Event.cc' for details.
426 const double JetMatchingAlpgen::GHOSTENERGY = 1e-15;
427 
428 // A zero threshold value for double comparisons.
429 const double JetMatchingAlpgen::ZEROTHRESHOLD = 1e-10;
430 
431 //--------------------------------------------------------------------------
432 
433 // Function to sort typeIdx vectors into descending eT/pT order.
434 // Uses a selection sort, as number of partons generally small
435 // and so efficiency not a worry.
436 
437 inline void JetMatchingAlpgen::sortTypeIdx(vector < int > &vecIn) {
438  for (size_t i = 0; i < vecIn.size(); i++) {
439  size_t jMax = i;
440  double vMax = (jetAlgorithm == 1) ?
441  eventProcess[vecIn[i]].eT() :
442  eventProcess[vecIn[i]].pT();
443  for (size_t j = i + 1; j < vecIn.size(); j++) {
444  double vNow = (jetAlgorithm == 1)
445  ? eventProcess[vecIn[j]].eT() : eventProcess[vecIn[j]].pT();
446  if (vNow > vMax) {
447  vMax = vNow;
448  jMax = j;
449  }
450  }
451  if (jMax != i) swap(vecIn[i], vecIn[jMax]);
452  }
453 }
454 
455 //--------------------------------------------------------------------------
456 
457 // Initialisation routine automatically called from Pythia::init().
458 // Setup all parts needed for the merging.
459 
460 inline bool JetMatchingAlpgen::initAfterBeams() {
461 
462  // Read in parameters
463  doMerge = settingsPtr->flag("JetMatching:merge");
464  jetAlgorithm = settingsPtr->mode("JetMatching:jetAlgorithm");
465  nJet = settingsPtr->mode("JetMatching:nJet");
466  nJetMax = settingsPtr->mode("JetMatching:nJetMax");
467  eTjetMin = settingsPtr->parm("JetMatching:eTjetMin");
468  coneRadius = settingsPtr->parm("JetMatching:coneRadius");
469  etaJetMax = settingsPtr->parm("JetMatching:etaJetMax");
470  doShowerKt = settingsPtr->flag("JetMatching:doShowerKt");
471 
472  // Use etaJetMax + coneRadius in input to jet algorithms
473  etaJetMaxAlgo = etaJetMax + coneRadius;
474 
475  // CellJet specific
476  nEta = settingsPtr->mode("JetMatching:nEta");
477  nPhi = settingsPtr->mode("JetMatching:nPhi");
478  eTseed = settingsPtr->parm("JetMatching:eTseed");
479  eTthreshold = settingsPtr->parm("JetMatching:eTthreshold");
480 
481  // SlowJet specific
482  slowJetPower = settingsPtr->mode("JetMatching:slowJetPower");
483  coneMatchLight = settingsPtr->parm("JetMatching:coneMatchLight");
484  coneRadiusHeavy = settingsPtr->parm("JetMatching:coneRadiusHeavy");
485  if (coneRadiusHeavy < 0.) coneRadiusHeavy = coneRadius;
486  coneMatchHeavy = settingsPtr->parm("JetMatching:coneMatchHeavy");
487 
488  // Matching procedure
489  jetAllow = settingsPtr->mode("JetMatching:jetAllow");
490  jetMatch = settingsPtr->mode("JetMatching:jetMatch");
491  exclusiveMode = settingsPtr->mode("JetMatching:exclusive");
492 
493  // If not merging, then done
494  if (!doMerge) return true;
495 
496  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
497  if (exclusiveMode == 2) {
498 
499  // No nJet or nJetMax, so default to exclusive mode
500  if (nJet < 0 || nJetMax < 0) {
501  infoPtr->errorMsg("Warning in JetMatchingAlpgen:init: "
502  "missing jet multiplicity information; running in exclusive mode");
503  exclusive = true;
504 
505  // Inclusive if nJet == nJetMax, exclusive otherwise
506  } else {
507  exclusive = (nJet == nJetMax) ? false : true;
508  }
509 
510  // Otherwise, just set as given
511  } else {
512  exclusive = (exclusiveMode == 0) ? false : true;
513  }
514 
515  // Initialise chosen jet algorithm. CellJet.
516  if (jetAlgorithm == 1) {
517 
518  // Extra options for CellJet. nSel = 1 means that all final-state
519  // particles are taken and we retain control of what to select.
520  // smear/resolution/upperCut are not used and are set to default values.
521  int nSel = 2, smear = 0;
522  double resolution = 0.5, upperCut = 2.;
523  cellJet = new CellJet(etaJetMaxAlgo, nEta, nPhi, nSel,
524  smear, resolution, upperCut, eTthreshold);
525 
526  // SlowJet
527  } else if (jetAlgorithm == 2) {
528  slowJet = new SlowJet(slowJetPower, coneRadius, eTjetMin, etaJetMaxAlgo);
529  }
530 
531  // Check the jetMatch parameter; option 2 only works with SlowJet
532  if (jetAlgorithm == 1 && jetMatch == 2) {
533  infoPtr->errorMsg("Warning in JetMatchingAlpgen:init: "
534  "jetMatch = 2 only valid with SlowJet algorithm. "
535  "Reverting to jetMatch = 1");
536  jetMatch = 1;
537  }
538 
539  // Setup local event records
540  eventProcessOrig.init("(eventProcessOrig)", particleDataPtr);
541  eventProcess.init("(eventProcess)", particleDataPtr);
542  workEventJet.init("(workEventJet)", particleDataPtr);
543 
544  // Print information
545  string jetStr = (jetAlgorithm == 1) ? "CellJet" :
546  (slowJetPower == -1) ? "anti-kT" :
547  (slowJetPower == 0) ? "C/A" :
548  (slowJetPower == 1) ? "kT" : "unknown";
549  string modeStr = (exclusive) ? "exclusive" : "inclusive";
550  stringstream nJetStr, nJetMaxStr;
551  if (nJet >= 0) nJetStr << nJet; else nJetStr << "unknown";
552  if (nJetMax >= 0) nJetMaxStr << nJetMax; else nJetMaxStr << "unknown";
553  cout << endl
554  << " *------- MLM matching parameters -------*" << endl
555  << " | nJet | " << setw(14)
556  << nJetStr.str() << " |" << endl
557  << " | nJetMax | " << setw(14)
558  << nJetMaxStr.str() << " |" << endl
559  << " | Jet algorithm | " << setw(14)
560  << jetStr << " |" << endl
561  << " | eTjetMin | " << setw(14)
562  << eTjetMin << " |" << endl
563  << " | coneRadius | " << setw(14)
564  << coneRadius << " |" << endl
565  << " | etaJetMax | " << setw(14)
566  << etaJetMax << " |" << endl
567  << " | jetAllow | " << setw(14)
568  << jetAllow << " |" << endl
569  << " | jetMatch | " << setw(14)
570  << jetMatch << " |" << endl
571  << " | coneMatchLight | " << setw(14)
572  << coneMatchLight << " |" << endl
573  << " | coneRadiusHeavy | " << setw(14)
574  << coneRadiusHeavy << " |" << endl
575  << " | coneMatchHeavy | " << setw(14)
576  << coneMatchHeavy << " |" << endl
577  << " | Mode | " << setw(14)
578  << modeStr << " |" << endl
579  << " *-----------------------------------------*" << endl;
580 
581  return true;
582 }
583 
584 //--------------------------------------------------------------------------
585 
586 // Step (1): sort the incoming particles
587 
588 inline void JetMatchingAlpgen::sortIncomingProcess(const Event &event) {
589 
590  // Remove resonance decays from original process and keep only final
591  // state. Resonances will have positive status code after this step.
592  omitResonanceDecays(eventProcessOrig, true);
593  eventProcess = workEvent;
594 
595  // Sort original process final state into light/heavy jets and 'other'.
596  // Criteria:
597  // 1 <= ID <= 5 and massless, or ID == 21 --> light jet (typeIdx[0])
598  // 4 <= ID <= 6 and massive --> heavy jet (typeIdx[1])
599  // All else --> other (typeIdx[2])
600  // Note that 'typeIdx' stores indices into 'eventProcess' (after resonance
601  // decays are omitted), while 'typeSet' stores indices into the original
602  // process record, 'eventProcessOrig', but these indices are also valid
603  // in 'event'.
604  for (int i = 0; i < 3; i++) {
605  typeIdx[i].clear();
606  typeSet[i].clear();
607  }
608  for (int i = 0; i < eventProcess.size(); i++) {
609  // Ignore nonfinal and default to 'other'
610  if (!eventProcess[i].isFinal()) continue;
611  int idx = 2;
612 
613  // Light jets
614  if (eventProcess[i].id() == ID_GLUON
615  || (eventProcess[i].idAbs() <= ID_BOT
616  && abs(eventProcess[i].m()) < ZEROTHRESHOLD)) idx = 0;
617 
618  // Heavy jets
619  else if (eventProcess[i].idAbs() >= ID_CHARM
620  && eventProcess[i].idAbs() <= ID_TOP) idx = 1;
621 
622  // Store
623  typeIdx[idx].push_back(i);
624  typeSet[idx].insert(eventProcess[i].daughter1());
625  }
626 
627  // Extract partons from hardest subsystem + ISR + FSR only into
628  // workEvent. Note no resonance showers or MPIs.
629  subEvent(event);
630 }
631 
632 //--------------------------------------------------------------------------
633 
634 // Step (2a): pick which particles to pass to the jet algorithm
635 
636 inline void JetMatchingAlpgen::jetAlgorithmInput(const Event &event,
637  int iType) {
638 
639  // Take input from 'workEvent' and put output in 'workEventJet'
640  workEventJet = workEvent;
641 
642  // Loop over particles and decide what to pass to the jet algorithm
643  for (int i = 0; i < workEventJet.size(); ++i) {
644  if (!workEventJet[i].isFinal()) continue;
645 
646  // jetAllow option to disallow certain particle types
647  if (jetAllow == 1) {
648 
649  // Original AG+Py6 algorithm explicitly excludes tops,
650  // leptons and photons.
651  int id = workEventJet[i].idAbs();
652  if ( (id >= ID_LEPMIN && id <= ID_LEPMAX) || id == ID_TOP
653  || id == ID_PHOTON) {
654  workEventJet[i].statusNeg();
655  continue;
656  }
657  }
658 
659  // Get the index of this particle in original event
660  int idx = workEventJet[i].daughter1();
661 
662  // Start with particle idx, and afterwards track mothers
663  while (true) {
664 
665  // Light jets
666  if (iType == 0) {
667 
668  // Do not include if originates from heavy jet or 'other'
669  if (typeSet[1].find(idx) != typeSet[1].end() ||
670  typeSet[2].find(idx) != typeSet[2].end()) {
671  workEventJet[i].statusNeg();
672  break;
673  }
674 
675  // Made it to start of event record so done
676  if (idx == 0) break;
677  // Otherwise next mother and continue
678  idx = event[idx].mother1();
679 
680  // Heavy jets
681  } else if (iType == 1) {
682 
683  // Only include if originates from heavy jet
684  if (typeSet[1].find(idx) != typeSet[1].end()) break;
685 
686  // Made it to start of event record with no heavy jet mother,
687  // so DO NOT include particle
688  if (idx == 0) {
689  workEventJet[i].statusNeg();
690  break;
691  }
692 
693  // Otherwise next mother and continue
694  idx = event[idx].mother1();
695 
696  // Other jets
697  } else if (iType == 2) {
698 
699  // Only include if originates from other jet
700  if (typeSet[2].find(idx) != typeSet[2].end()) break;
701 
702  // Made it to start of event record with no heavy jet mother,
703  // so DO NOT include particle
704  if (idx == 0) {
705  workEventJet[i].statusNeg();
706  break;
707  }
708 
709  // Otherwise next mother and continue
710  idx = event[idx].mother1();
711 
712  } // if (iType)
713  } // while (true)
714  } // for (i)
715 
716  // For jetMatch = 2, insert ghost particles corresponding to
717  // each hard parton in the original process
718  if (jetMatch == 2) {
719  for (int i = 0; i < int(typeIdx[iType].size()); i++) {
720  // Get y/phi of the parton
721  Vec4 pIn = eventProcess[typeIdx[iType][i]].p();
722  double y = pIn.rap();
723  double phi = pIn.phi();
724 
725  // Create a ghost particle and add to the workEventJet
726  double e = GHOSTENERGY;
727  double e2y = exp(2. * y);
728  double pz = e * (e2y - 1.) / (e2y + 1.);
729  double pt = sqrt(e*e - pz*pz);
730  double px = pt * cos(phi);
731  double py = pt * sin(phi);
732  workEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0, px, py, pz, e);
733 
734  // Extra check on reconstructed y/phi values. If many warnings
735  // of this type, GHOSTENERGY may be set too low.
736  if (MATCHINGCHECK) {
737  int lastIdx = workEventJet.size() - 1;
738  if (abs(y - workEventJet[lastIdx].y()) > ZEROTHRESHOLD ||
739  abs(phi - workEventJet[lastIdx].phi()) > ZEROTHRESHOLD)
740  infoPtr->errorMsg("Warning in JetMatchingAlpgen:jetAlgorithmInput: "
741  "ghost particle y/phi mismatch");
742  }
743 
744  } // for (i)
745  } // if (jetMatch == 2)
746 }
747 
748 //--------------------------------------------------------------------------
749 
750 // Step (2b): run jet algorithm and provide common output
751 
752 inline void JetMatchingAlpgen::runJetAlgorithm() {
753 
754  // Run the jet clustering algorithm
755  if (jetAlgorithm == 1)
756  cellJet->analyze(workEventJet, eTjetMin, coneRadius, eTseed);
757  else
758  slowJet->analyze(workEventJet);
759 
760  // Extract four-momenta of jets with |eta| < etaJetMax and
761  // put into jetMomenta. Note that this is done backwards as
762  // jets are removed with SlowJet.
763  jetMomenta.clear();
764  int iJet = (jetAlgorithm == 1) ? cellJet->size() - 1:
765  slowJet->sizeJet() - 1;
766  for (int i = iJet; i > -1; i--) {
767  Vec4 jetMom = (jetAlgorithm == 1) ? cellJet->pMassive(i) :
768  slowJet->p(i);
769  double eta = jetMom.eta();
770 
771  if (abs(eta) > etaJetMax) {
772  if (jetAlgorithm == 2) slowJet->removeJet(i);
773  continue;
774  }
775  jetMomenta.push_back(jetMom);
776  }
777 
778  // Reverse jetMomenta to restore eT/pT ordering
779  reverse(jetMomenta.begin(), jetMomenta.end());
780 }
781 
782 //--------------------------------------------------------------------------
783 
784 // Step (2c): veto decision (returning true vetoes the event)
785 
786 inline bool JetMatchingAlpgen::matchPartonsToJets(int iType) {
787 
788  // Use two different routines for light/heavy jets as
789  // different veto conditions and for clarity
790  if (iType == 0) return (matchPartonsToJetsLight() > 0);
791  else if (iType == 1) return (matchPartonsToJetsHeavy() > 0);
792  else if (iType == 2) return false;
793  return true;
794 }
795 
796 //--------------------------------------------------------------------------
797 
798 // Step(2c): light jets
799 // Return codes are given indicating the reason for a veto.
800 // Although not currently used, they are a useful debugging tool:
801 // 0 = no veto
802 // 1 = veto as number of jets less than number of partons
803 // 2 = veto as exclusive mode and number of jets greater than
804 // number of partons
805 // 3 = veto as inclusive mode and there would be an extra jet
806 // that is harder than any matched soft jet
807 // 4 = veto as there is a parton which does not match a jet
808 
809 inline int JetMatchingAlpgen::matchPartonsToJetsLight() {
810 
811  // Always veto if number of jets is less than original number of jets
812  if (jetMomenta.size() < typeIdx[0].size()) return LESS_JETS;
813  // Veto if in exclusive mode and number of jets bigger than original
814  if (exclusive && jetMomenta.size() > typeIdx[0].size()) return MORE_JETS;
815 
816  // Sort partons by eT/pT
817  sortTypeIdx(typeIdx[0]);
818 
819  // Number of hard partons
820  int nParton = typeIdx[0].size();
821 
822  // Keep track of which jets have been assigned a hard parton
823  vector < bool > jetAssigned;
824  jetAssigned.assign(jetMomenta.size(), false);
825 
826  // Jet matching procedure: (1) deltaR between partons and jets
827  if (jetMatch == 1) {
828 
829  // Loop over light hard partons and get 4-momentum
830  for (int i = 0; i < nParton; i++) {
831  Vec4 p1 = eventProcess[typeIdx[0][i]].p();
832 
833  // Track which jet has the minimal dR measure with this parton
834  int jMin = -1;
835  double dRmin = 0.;
836 
837  // Loop over all jets (skipping those already assigned).
838  for (int j = 0; j < int(jetMomenta.size()); j++) {
839  if (jetAssigned[j]) continue;
840 
841  // DeltaR between parton/jet and store if minimum
842  double dR = (jetAlgorithm == 1)
843  ? REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
844  if (jMin < 0 || dR < dRmin) {
845  dRmin = dR;
846  jMin = j;
847  }
848  } // for (j)
849 
850  // Check for jet-parton match
851  if (jMin >= 0 && dRmin < coneRadius * coneMatchLight) {
852 
853  // If the matched jet is not one of the nParton hardest jets,
854  // the extra left over jet would be harder than some of the
855  // matched jets. This is disallowed, so veto.
856  if (jMin >= nParton) return HARD_JET;
857 
858  // Mark jet as assigned.
859  jetAssigned[jMin] = true;
860 
861  // If no match, then event will be vetoed in all cases
862  } else return UNMATCHED_PARTON;
863 
864  } // for (i)
865 
866  // Jet matching procedure: (2) ghost particles in SlowJet
867  } else {
868 
869  // Loop over added 'ghost' particles and find if assigned to a jet
870  for (int i = workEventJet.size() - nParton;
871  i < workEventJet.size(); i++) {
872  int jMin = slowJet->jetAssignment(i);
873 
874  // Veto if:
875  // 1) not one of nParton hardest jets
876  // 2) not assigned to a jet
877  // 3) jet has already been assigned
878  if (jMin >= nParton) return HARD_JET;
879  if (jMin < 0 || jetAssigned[jMin]) return UNMATCHED_PARTON;
880 
881  // Mark jet as assigned
882  jetAssigned[jMin] = true;
883 
884  } // for (i)
885  } // if (jetMatch)
886 
887  // Minimal eT/pT (CellJet/SlowJet) of matched light jets. Needed
888  // later for heavy jet vetos in inclusive mode.
889  if (nParton > 0)
890  eTpTlightMin = (jetAlgorithm == 1) ? jetMomenta[nParton - 1].eT()
891  : jetMomenta[nParton - 1].pT();
892  else
893  eTpTlightMin = -1.;
894 
895  // No veto
896  return NONE;
897 }
898 
899 //--------------------------------------------------------------------------
900 
901 // Step(2c): heavy jets
902 // Return codes are given indicating the reason for a veto.
903 // Although not currently used, they are a useful debugging tool:
904 // 0 = no veto as there are no extra jets present
905 // 1 = veto as in exclusive mode and extra jets present
906 // 2 = veto as in inclusive mode and extra jets were harder
907 // than any matched light jet
908 
909 inline int JetMatchingAlpgen::matchPartonsToJetsHeavy() {
910 
911  // If there are no extra jets, then accept
912  if (jetMomenta.empty()) return NONE;
913 
914  // Number of hard partons
915  int nParton = typeIdx[1].size();
916 
917  // Remove jets that are close to heavy quarks
918  set < int > removeJets;
919 
920  // Jet matching procedure: (1) deltaR between partons and jets
921  if (jetMatch == 1) {
922 
923  // Loop over heavy hard partons and get 4-momentum
924  for (int i = 0; i < nParton; i++) {
925  Vec4 p1 = eventProcess[typeIdx[1][i]].p();
926 
927  // Loop over all jets, find dR and mark for removal if match
928  for (int j = 0; j < int(jetMomenta.size()); j++) {
929  double dR = (jetAlgorithm == 1) ?
930  REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
931  if (dR < coneRadiusHeavy * coneMatchHeavy)
932  removeJets.insert(j);
933 
934  } // for (j)
935  } // for (i)
936 
937  // Jet matching procedure: (2) ghost particles in SlowJet
938  } else {
939 
940  // Loop over added 'ghost' particles and if assigned to a jet
941  // then mark this jet for removal
942  for (int i = workEventJet.size() - nParton;
943  i < workEventJet.size(); i++) {
944  int jMin = slowJet->jetAssignment(i);
945  if (jMin >= 0) removeJets.insert(jMin);
946  }
947 
948  }
949 
950  // Remove jets (backwards order to not disturb indices)
951  for (set < int >::reverse_iterator it = removeJets.rbegin();
952  it != removeJets.rend(); it++)
953  jetMomenta.erase(jetMomenta.begin() + *it);
954 
955  // Handle case if there are still extra jets
956  if (!jetMomenta.empty()) {
957 
958  // Exclusive mode, so immediate veto
959  if (exclusive) return MORE_JETS;
960 
961  // Inclusive mode; extra jets must be softer than any matched light jet
962  else if (eTpTlightMin >= 0.)
963  for (size_t j = 0; j < jetMomenta.size(); j++) {
964  // CellJet uses eT, SlowJet uses pT
965  if ( (jetAlgorithm == 1 && jetMomenta[j].eT() > eTpTlightMin) ||
966  (jetAlgorithm == 2 && jetMomenta[j].pT() > eTpTlightMin) )
967  return HARD_JET;
968  }
969 
970  } // if (!jetMomenta.empty())
971 
972  // No extra jets were present so no veto
973  return NONE;
974 }
975 
976 //==========================================================================
977 
978 // Main implementation of Madgraph UserHooks class.
979 // This may be split out to a separate C++ file if desired,
980 // but currently included here for ease of use.
981 
982 //--------------------------------------------------------------------------
983 
984 // Initialisation routine automatically called from Pythia::init().
985 // Setup all parts needed for the merging.
986 
987 inline bool JetMatchingMadgraph::initAfterBeams() {
988 
989  // Initialise values for stored jet matching veto inputs.
990  pTfirstSave = -1.;
991  processSubsetSave.init("(eventProcess)", particleDataPtr);
992  workEventJetSave.init("(workEventJet)", particleDataPtr);
993 
994  // Read in Madgraph specific configuration variables
995  bool setMad = settingsPtr->flag("JetMatching:setMad");
996 
997  // If Madgraph parameters are present, then parse in MadgraphPar object
998  MadgraphPar par(infoPtr);
999  string parStr = infoPtr->header("MGRunCard");
1000  if (!parStr.empty()) {
1001  par.parse(parStr);
1002  par.printParams();
1003  }
1004 
1005  // Set Madgraph merging parameters from the file if requested
1006  if (setMad) {
1007  if ( par.haveParam("xqcut") && par.haveParam("maxjetflavor")
1008  && par.haveParam("alpsfact") && par.haveParam("ickkw") ) {
1009  settingsPtr->flag("JetMatching:merge", par.getParam("ickkw"));
1010  settingsPtr->parm("JetMatching:qCut", par.getParam("xqcut"));
1011  settingsPtr->mode("JetMatching:nQmatch",
1012  par.getParamAsInt("maxjetflavor"));
1013  settingsPtr->parm("JetMatching:clFact",
1014  clFact = par.getParam("alpsfact"));
1015  if (par.getParamAsInt("ickkw") == 0)
1016  infoPtr->errorMsg("Error in JetMatchingMadgraph:init: "
1017  "Madgraph file parameters are not set for merging");
1018 
1019  // Warn if setMad requested, but one or more parameters not present
1020  } else {
1021  infoPtr->errorMsg("Warning in JetMatchingMadgraph:init: "
1022  "Madgraph merging parameters not found");
1023  if (!par.haveParam("xqcut")) infoPtr->errorMsg("Warning in "
1024  "JetMatchingMadgraph:init: No xqcut");
1025  if (!par.haveParam("ickkw")) infoPtr->errorMsg("Warning in "
1026  "JetMatchingMadgraph:init: No ickkw");
1027  if (!par.haveParam("maxjetflavor")) infoPtr->errorMsg("Warning in "
1028  "JetMatchingMadgraph:init: No maxjetflavor");
1029  if (!par.haveParam("alpsfact")) infoPtr->errorMsg("Warning in "
1030  "JetMatchingMadgraph:init: No alpsfact");
1031  }
1032  }
1033 
1034  // Read in FxFx matching parameters
1035  doFxFx = settingsPtr->flag("JetMatching:doFxFx");
1036  nPartonsNow = settingsPtr->mode("JetMatching:nPartonsNow");
1037  qCutME = settingsPtr->parm("JetMatching:qCutME");
1038  qCutMESq = pow(qCutME,2);
1039 
1040  // Read in Madgraph merging parameters
1041  doMerge = settingsPtr->flag("JetMatching:merge");
1042  doShowerKt = settingsPtr->flag("JetMatching:doShowerKt");
1043  qCut = settingsPtr->parm("JetMatching:qCut");
1044  nQmatch = settingsPtr->mode("JetMatching:nQmatch");
1045  clFact = settingsPtr->parm("JetMatching:clFact");
1046 
1047  // Read in jet algorithm parameters
1048  jetAlgorithm = settingsPtr->mode("JetMatching:jetAlgorithm");
1049  nJetMax = settingsPtr->mode("JetMatching:nJetMax");
1050  eTjetMin = settingsPtr->parm("JetMatching:eTjetMin");
1051  coneRadius = settingsPtr->parm("JetMatching:coneRadius");
1052  etaJetMax = settingsPtr->parm("JetMatching:etaJetMax");
1053  slowJetPower = settingsPtr->mode("JetMatching:slowJetPower");
1054 
1055  // Matching procedure
1056  jetAllow = settingsPtr->mode("JetMatching:jetAllow");
1057  exclusiveMode = settingsPtr->mode("JetMatching:exclusive");
1058  qCutSq = pow(qCut,2);
1059  etaJetMaxAlgo = etaJetMax;
1060 
1061  // Read if veto should be performed internally.
1062  performVeto = settingsPtr->flag("JetMatching:doVeto");
1063 
1064  // If not merging, then done
1065  if (!doMerge) return true;
1066 
1067  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
1068  if (exclusiveMode == 2) {
1069 
1070  // No nJet or nJetMax, so default to exclusive mode
1071  if (nJetMax < 0) {
1072  infoPtr->errorMsg("Warning in JetMatchingMadgraph:init: "
1073  "missing jet multiplicity information; running in exclusive mode");
1074  exclusiveMode = 1;
1075  }
1076  }
1077 
1078  // Initialise chosen jet algorithm.
1079  // Currently, this only supports the kT-algorithm in SlowJet.
1080  // Use the QCD distance measure by default.
1081  jetAlgorithm = 2;
1082  slowJetPower = 1;
1083  slowJet = new SlowJet(slowJetPower, coneRadius, eTjetMin,
1084  etaJetMaxAlgo, 2, 2, NULL, false);
1085 
1086  // For FxFx, also initialise jet algorithm to define matrix element jets.
1087  // Currently, this only supports the kT-algorithm in SlowJet.
1088  // Use the QCD distance measure by default.
1089  slowJetHard = new SlowJet(slowJetPower, coneRadius, qCutME,
1090  etaJetMaxAlgo, 2, 2, NULL, false);
1091 
1092  // To access the DJR's
1093  slowJetDJR = new SlowJet(slowJetPower, coneRadius, qCutME,
1094  etaJetMaxAlgo, 2, 2, NULL, false);
1095 
1096  // A special version of SlowJet to handle heavy and other partons
1097  hjSlowJet = new HJSlowJet(slowJetPower, coneRadius, 0.0,
1098  100.0, 1, 2, NULL, false, true);
1099 
1100  // Setup local event records
1101  eventProcessOrig.init("(eventProcessOrig)", particleDataPtr);
1102  eventProcess.init("(eventProcess)", particleDataPtr);
1103  workEventJet.init("(workEventJet)", particleDataPtr);
1104 
1105  // Print information
1106  string jetStr = (jetAlgorithm == 1) ? "CellJet" :
1107  (slowJetPower == -1) ? "anti-kT" :
1108  (slowJetPower == 0) ? "C/A" :
1109  (slowJetPower == 1) ? "kT" : "unknown";
1110  string modeStr = (exclusiveMode) ? "exclusive" : "inclusive";
1111  cout << endl
1112  << " *----- Madgraph matching parameters -----*" << endl
1113  << " | qCut | " << setw(14)
1114  << qCut << " |" << endl
1115  << " | nQmatch | " << setw(14)
1116  << nQmatch << " |" << endl
1117  << " | clFact | " << setw(14)
1118  << clFact << " |" << endl
1119  << " | Jet algorithm | " << setw(14)
1120  << jetStr << " |" << endl
1121  << " | eTjetMin | " << setw(14)
1122  << eTjetMin << " |" << endl
1123  << " | etaJetMax | " << setw(14)
1124  << etaJetMax << " |" << endl
1125  << " | jetAllow | " << setw(14)
1126  << jetAllow << " |" << endl
1127  << " | Mode | " << setw(14)
1128  << modeStr << " |" << endl
1129  << " *-----------------------------------------*" << endl;
1130 
1131  return true;
1132 }
1133 
1134 //--------------------------------------------------------------------------
1135 
1136 // Process level vetos
1137 
1138 inline bool JetMatchingMadgraph::doVetoProcessLevel(Event& process) {
1139 
1140  eventProcessOrig = process;
1141 
1142  // Setup for veto if hard ME has too many partons.
1143  // This is done to achieve consistency with the Pythia6 implementation.
1144 
1145  // Clear the event of MPI systems and resonace decay products. Store trimmed
1146  // event in workEvent.
1147  sortIncomingProcess(process);
1148 
1149  // Veto in case the hard input matrix element already has too many partons.
1150  if ( !doFxFx && int(typeIdx[0].size()) > nJetMax )
1151  return true;
1152  if ( doFxFx && npNLO() < nJetMax && int(typeIdx[0].size()) > nJetMax )
1153  return true;
1154 
1155  // Done
1156  return false;
1157 
1158 }
1159 
1160 //--------------------------------------------------------------------------
1161 
1162 inline bool JetMatchingMadgraph::doVetoStep(int iPos, int nISR, int nFSR,
1163  const Event& event) {
1164 
1165  // Do not perform any veto if not in the Shower-kT scheme.
1166  if ( !doShowerKt ) return false;
1167 
1168  // Do nothing for emissions after the first one.
1169  if ( nISR + nFSR > 1 ) return false;
1170 
1171  // Do nothing in resonance decay showers.
1172  if (iPos == 5) return false;
1173 
1174  // Clear the event of MPI systems and resonace decay products. Store trimmed
1175  // event in workEvent.
1176  sortIncomingProcess(event);
1177 
1178  // Get (kinematical) pT of first emission
1179  double pTfirst = 0.;
1180 
1181  // Get weak bosons, for later checks if the emission is a "QCD emission".
1182  vector<int> weakBosons;
1183  for (int i = 0; i < event.size(); i++) {
1184  if ( event[i].id() == 22
1185  && event[i].id() == 23
1186  && event[i].idAbs() == 24)
1187  weakBosons.push_back(i);
1188  }
1189 
1190  for (int i = workEvent.size()-1; i > 0; --i) {
1191  if ( workEvent[i].isFinal() && workEvent[i].colType() != 0
1192  && (workEvent[i].statusAbs() == 43 || workEvent[i].statusAbs() == 51)) {
1193  // Check if any of the EW bosons are ancestors of this parton. This
1194  // should never happen for the first non-resonance shower emission.
1195  // Check just to be sure.
1196  bool QCDemission = true;
1197  // Get position of this parton in the actual event (workEvent does
1198  // not contain right mother-daughter relations). Stored in daughters.
1199  int iPosOld = workEvent[i].daughter1();
1200  for (int j = 0; i < int(weakBosons.size()); ++i)
1201  if ( event[iPosOld].isAncestor(j)) {
1202  QCDemission = false;
1203  break;
1204  }
1205  // Done for a QCD emission.
1206  if (QCDemission){
1207  pTfirst = workEvent[i].pT();
1208  break;
1209  }
1210  }
1211  }
1212 
1213  // Store things that are necessary to perform the shower-kT veto externally.
1214  pTfirstSave = pTfirst;
1215  // Done if only inputs for an external vetoing procedure should be stored.
1216  if (!performVeto) return false;
1217 
1218  // Check veto.
1219  if ( doShowerKtVeto(pTfirst) ) return true;
1220 
1221  // No veto if come this far.
1222  return false;
1223 
1224 }
1225 
1226 //--------------------------------------------------------------------------
1227 
1228 inline bool JetMatchingMadgraph::doShowerKtVeto(double pTfirst) {
1229 
1230  // Only check veto in the shower-kT scheme.
1231  if ( !doShowerKt ) return false;
1232 
1233  // Reset veto code
1234  bool doVeto = false;
1235 
1236  // Find the (kinematical) pT of the softest (light) parton in the hard
1237  // process.
1238  int nParton = typeIdx[0].size();
1239  double pTminME=1e10;
1240  for ( int i = 0; i < nParton; ++i)
1241  pTminME = min(pTminME,eventProcess[typeIdx[0][i]].pT());
1242 
1243  // Veto if the softest hard process parton is below Qcut.
1244  if ( nParton > 0 && pow(pTminME,2) < qCutSq ) doVeto = true;
1245 
1246  // For non-highest multiplicity, veto if the hardest emission is harder
1247  // than Qcut.
1248  if ( exclusive && pow(pTfirst,2) > qCutSq ) {
1249  doVeto = true;
1250  // For highest multiplicity sample, veto if the hardest emission is harder
1251  // than the hard process parton.
1252  } else if ( !exclusive && nParton > 0 && pTfirst > pTminME ) {
1253  doVeto = true;
1254  }
1255 
1256  // Return veto
1257  return doVeto;
1258 
1259 }
1260 
1261 //--------------------------------------------------------------------------
1262 
1263 // Function to set the jet clustering scales (to be used as output)
1264 
1265 inline void JetMatchingMadgraph::setDJR( const Event& event) {
1266 
1267  // Clear members.
1268  clearDJR();
1269  vector<double> result;
1270 
1271  // Initialize SlowJetDJR jet algorithm with event
1272  if (!slowJetDJR->setup(event) ) {
1273  infoPtr->errorMsg("Warning in JetMatchingMadgraph:setDJR"
1274  ": the SlowJet algorithm failed on setup");
1275  return;
1276  }
1277 
1278  // Cluster in steps to find all hadronic jets
1279  while ( slowJetDJR->sizeAll() - slowJetDJR->sizeJet() > 0 ) {
1280  // Save the next clustering scale.
1281  result.push_back(sqrt(slowJetDJR->dNext()));
1282  // Perform step.
1283  slowJetDJR->doStep();
1284  }
1285 
1286  // Save clustering scales in reserve order.
1287  for (int i=int(result.size())-1; i >= 0; --i)
1288  DJR.push_back(result[i]);
1289 
1290 }
1291 
1292 //--------------------------------------------------------------------------
1293 
1294 // Function to get the current number of partons in the Born state, as
1295 // read from LHE.
1296 
1297 inline int JetMatchingMadgraph::npNLO(){
1298  string npIn = infoPtr->getEventAttribute("npNLO",true);
1299  int np = (npIn != "") ? atoi((char*)npIn.c_str()) : -1;
1300  if ( np < 0 ) { ; }
1301  else return np;
1302  return nPartonsNow;
1303 }
1304 
1305 //--------------------------------------------------------------------------
1306 
1307 // Step (1): sort the incoming particles
1308 
1309 inline void JetMatchingMadgraph::sortIncomingProcess(const Event &event) {
1310 
1311  // Remove resonance decays from original process and keep only final
1312  // state. Resonances will have positive status code after this step.
1313  omitResonanceDecays(eventProcessOrig, true);
1314  clearDJR();
1315  clear_nMEpartons();
1316 
1317  // For FxFx, pre-cluster partons in the event into jets.
1318  if (doFxFx) {
1319 
1320  // Get final state partons
1321  eventProcess.clear();
1322  workEventJet.clear();
1323  for( int i=0; i < workEvent.size(); ++i) {
1324  // Original AG+Py6 algorithm explicitly excludes tops,
1325  // leptons and photons.
1326  int id = workEvent[i].idAbs();
1327  if ((id >= ID_LEPMIN && id <= ID_LEPMAX) || id == ID_TOP
1328  || id == ID_PHOTON || id == 23 || id == 24 || id == 25) {
1329  eventProcess.append(workEvent[i]);
1330  } else {
1331  workEventJet.append(workEvent[i]);
1332  }
1333  }
1334 
1335  // Initialize SlowJetHard jet algorithm with current working event
1336  if (!slowJetHard->setup(workEventJet) ) {
1337  infoPtr->errorMsg("Warning in JetMatchingMadgraph:sortIncomingProcess"
1338  ": the SlowJet algorithm failed on setup");
1339  return;
1340  }
1341 
1342  // Get matrix element cut scale.
1343  double localQcutSq = qCutMESq;
1344  // Cluster in steps to find all hadronic jets at the scale qCutME
1345  while ( slowJetHard->sizeAll() - slowJetHard->sizeJet() > 0 ) {
1346  // Done if next step is above qCut
1347  if( slowJetHard->dNext() > localQcutSq ) break;
1348  // Done if we're at or below the number of partons in the Born state.
1349  if( slowJetHard->sizeAll()-slowJetHard->sizeJet() <= npNLO()) break;
1350  slowJetHard->doStep();
1351  }
1352 
1353  // Construct a master copy of the event containing only the
1354  // hardest nPartonsNow hadronic clusters. While constructing the event,
1355  // the parton type (ID_GLUON) and status (98,99) are arbitrary.
1356  int nJets = slowJetHard->sizeJet();
1357  int nClus = slowJetHard->sizeAll();
1358  int nNow = 0;
1359  for (int i = nJets; i < nClus; ++i) {
1360  vector<int> parts;
1361  if (i < nClus-nJets) parts = slowJetHard->clusConstituents(i);
1362  else parts = slowJetHard->constituents(nClus-nJets-i);
1363  int flavour = ID_GLUON;
1364  for(int j=0; j < int(parts.size()); ++j)
1365  if (workEventJet[parts[j]].id() == ID_BOT)
1366  flavour = ID_BOT;
1367  eventProcess.append( flavour, 98,
1368  workEventJet[parts.back()].mother1(),
1369  workEventJet[parts.back()].mother2(),
1370  workEventJet[parts.back()].daughter1(),
1371  workEventJet[parts.back()].daughter2(),
1372  0, 0, slowJetHard->p(i).px(), slowJetHard->p(i).py(),
1373  slowJetHard->p(i).pz(), slowJetHard->p(i).e() );
1374  nNow++;
1375  }
1376 
1377  // Done. Clean-up
1378  workEventJet.clear();
1379 
1380  // For MLM matching, simply take hard process state from workEvent,
1381  // without any preclustering.
1382  } else {
1383  eventProcess = workEvent;
1384  }
1385 
1386  // Sort original process final state into light/heavy jets and 'other'.
1387  // Criteria:
1388  // 1 <= ID <= nQmatch, or ID == 21 --> light jet (typeIdx[0])
1389  // nQMatch < ID --> heavy jet (typeIdx[1])
1390  // All else that is colored --> other (typeIdx[2])
1391  // Note that 'typeIdx' stores indices into 'eventProcess' (after resonance
1392  // decays are omitted), while 'typeSet' stores indices into the original
1393  // process record, 'eventProcessOrig', but these indices are also valid
1394  // in 'event'.
1395  for (int i = 0; i < 3; i++) {
1396  typeIdx[i].clear();
1397  typeSet[i].clear();
1398  origTypeIdx[i].clear();
1399  }
1400  for (int i = 0; i < eventProcess.size(); i++) {
1401  // Ignore non-final state and default to 'other'
1402  if (!eventProcess[i].isFinal()) continue;
1403  int idx = -1;
1404  int orig_idx = -1;
1405 
1406  // Light jets: all gluons and quarks with id less than or equal to nQmatch
1407  if (eventProcess[i].isGluon()
1408  || (eventProcess[i].idAbs() <= nQmatch) ) {
1409  orig_idx = 0;
1410  // Crucial point: MG puts the scale of a non-QCD particle to eCM. For
1411  // such particles, we should keep the default "2"
1412  idx = ( eventProcess[i].scale() < 1.999 * sqrt(infoPtr->eA()
1413  * infoPtr->eB()) ) ? 0 : 2;
1414  }
1415 
1416  // Heavy jets: all quarks with id greater than nQmatch
1417  else if (eventProcess[i].idAbs() > nQmatch
1418  && eventProcess[i].idAbs() <= ID_TOP) {
1419  idx = 1;
1420  orig_idx = 1;
1421  // Update to include non-SM colored particles
1422  } else if (eventProcess[i].colType() != 0
1423  && eventProcess[i].idAbs() > ID_TOP) {
1424  idx = 1;
1425  orig_idx = 1;
1426  }
1427  if( idx < 0 ) continue;
1428  // Store
1429  typeIdx[idx].push_back(i);
1430  typeSet[idx].insert(eventProcess[i].daughter1());
1431  origTypeIdx[orig_idx].push_back(i);
1432  }
1433 
1434  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
1435  if (exclusiveMode == 2) {
1436 
1437  // Inclusive if nJet == nJetMax, exclusive otherwise
1438  int nParton = origTypeIdx[0].size();
1439  exclusive = (nParton == nJetMax) ? false : true;
1440 
1441  // Otherwise, just set as given
1442  } else {
1443  exclusive = (exclusiveMode == 0) ? false : true;
1444  }
1445 
1446  // Extract partons from hardest subsystem + ISR + FSR only into
1447  // workEvent. Note no resonance showers or MPIs.
1448  subEvent(event);
1449 
1450  // Store things that are necessary to perform the kT-MLM veto externally.
1451  int nParton = typeIdx[0].size();
1452  processSubsetSave.clear();
1453  for ( int i = 0; i < nParton; ++i)
1454  processSubsetSave.append( eventProcess[typeIdx[0][i]] );
1455 
1456 }
1457 
1458 //--------------------------------------------------------------------------
1459 
1460 // Step (2a): pick which particles to pass to the jet algorithm
1461 
1462 inline void JetMatchingMadgraph::jetAlgorithmInput(const Event &event,
1463  int iType) {
1464 
1465  // Take input from 'workEvent' and put output in 'workEventJet'
1466  workEventJet = workEvent;
1467 
1468  // Loop over particles and decide what to pass to the jet algorithm
1469  for (int i = 0; i < workEventJet.size(); ++i) {
1470  if (!workEventJet[i].isFinal()) continue;
1471 
1472  // jetAllow option to disallow certain particle types
1473  if (jetAllow == 1) {
1474  // Remove all non-QCD partons from veto list
1475  if( workEventJet[i].colType() == 0 ) {
1476  workEventJet[i].statusNeg();
1477  continue;
1478  }
1479  }
1480 
1481  // Get the index of this particle in original event
1482  int idx = workEventJet[i].daughter1();
1483 
1484  // Start with particle idx, and afterwards track mothers
1485  while (true) {
1486 
1487  // Light jets
1488  if (iType == 0) {
1489 
1490  // Do not include if originates from heavy jet or 'other'
1491  if (typeSet[1].find(idx) != typeSet[1].end() ||
1492  typeSet[2].find(idx) != typeSet[2].end()) {
1493  workEventJet[i].statusNeg();
1494  break;
1495  }
1496 
1497  // Made it to start of event record so done
1498  if (idx == 0) break;
1499  // Otherwise next mother and continue
1500  idx = event[idx].mother1();
1501 
1502  // Heavy jets
1503  } else if (iType == 1) {
1504 
1505  // Only include if originates from heavy jet
1506  if (typeSet[1].find(idx) != typeSet[1].end()) break;
1507 
1508  // Made it to start of event record with no heavy jet mother,
1509  // so DO NOT include particle
1510  if (idx == 0) {
1511  workEventJet[i].statusNeg();
1512  break;
1513  }
1514 
1515  // Otherwise next mother and continue
1516  idx = event[idx].mother1();
1517 
1518  // Other jets
1519  } else if (iType == 2) {
1520 
1521  // Only include if originates from other jet
1522  if (typeSet[2].find(idx) != typeSet[2].end()) break;
1523 
1524  // Made it to start of event record with no heavy jet mother,
1525  // so DO NOT include particle
1526  if (idx == 0) {
1527  workEventJet[i].statusNeg();
1528  break;
1529  }
1530 
1531  // Otherwise next mother and continue
1532  idx = event[idx].mother1();
1533 
1534  } // if (iType)
1535  } // while (true)
1536  } // for (i)
1537 }
1538 
1539 //--------------------------------------------------------------------------
1540 
1541 // Step (2b): run jet algorithm and provide common output
1542 // This does nothing, because the jet algorithm is run several times
1543 // in the matching algorithm.
1544 
1545 inline void JetMatchingMadgraph::runJetAlgorithm() {; }
1546 
1547 //--------------------------------------------------------------------------
1548 
1549 // Step (2c): veto decision (returning true vetoes the event)
1550 
1551 inline bool JetMatchingMadgraph::matchPartonsToJets(int iType) {
1552 
1553  // Use different routines for light/heavy/other jets as
1554  // different veto conditions and for clarity
1555  if (iType == 0) {
1556  // Record the jet separations here, also if matchPartonsToJetsLight
1557  // returns preemptively.
1558  setDJR(workEventJet);
1559  set_nMEpartons(origTypeIdx[0].size(), typeIdx[0].size());
1560  // Perform jet matching.
1561  return (matchPartonsToJetsLight() > 0);
1562  } else if (iType == 1) {
1563  return (matchPartonsToJetsHeavy() > 0);
1564  } else {
1565  return (matchPartonsToJetsOther() > 0);
1566  }
1567 
1568 }
1569 
1570 //--------------------------------------------------------------------------
1571 
1572 // Step(2c): light jets
1573 // Return codes are given indicating the reason for a veto.
1574 // Although not currently used, they are a useful debugging tool:
1575 // 0 = no veto
1576 // 1 = veto as number of jets less than number of partons
1577 // 2 = veto as exclusive mode and number of jets greater than
1578 // number of partons
1579 // 3 = veto as inclusive mode and there would be an extra jet
1580 // that is harder than any matched soft jet
1581 // 4 = veto as there is a parton which does not match a jet
1582 
1583 inline int JetMatchingMadgraph::matchPartonsToJetsLight() {
1584 
1585  // Store things that are necessary to perform the kT-MLM veto externally.
1586  workEventJetSave = workEventJet;
1587  // Done if only inputs for an external vetoing procedure should be stored.
1588  if (!performVeto) return false;
1589 
1590  // Count the number of hard partons
1591  int nParton = typeIdx[0].size();
1592 
1593  // Initialize SlowJet with current working event
1594  if (!slowJet->setup(workEventJet) ) {
1595  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1596  "Light: the SlowJet algorithm failed on setup");
1597  return NONE;
1598  }
1599  double localQcutSq = qCutSq;
1600  double dOld = 0.0;
1601  // Cluster in steps to find all hadronic jets at the scale qCut
1602  while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1603  if( slowJet->dNext() > localQcutSq ) break;
1604  dOld = slowJet->dNext();
1605  slowJet->doStep();
1606  }
1607  int nJets = slowJet->sizeJet();
1608  int nClus = slowJet->sizeAll();
1609 
1610  // Debug printout.
1611  if (MATCHINGDEBUG) slowJet->list(true);
1612 
1613  // Count of the number of hadronic jets in SlowJet accounting
1614  int nCLjets = nClus - nJets;
1615  // Get number of partons. Different for MLM and FxFx schemes.
1616  int nRequested = (doFxFx) ? npNLO() : nParton;
1617 
1618  // Veto event if too few hadronic jets
1619  if ( nCLjets < nRequested ) return LESS_JETS;
1620 
1621  // In exclusive mode, do not allow more hadronic jets than partons
1622  if ( exclusive && !doFxFx ) {
1623  if ( nCLjets > nRequested ) return MORE_JETS;
1624  } else {
1625 
1626  // For FxFx, in the non-highest multipicity, all jets need to matched to
1627  // partons. For nCLjets > nRequested, this is not possible. Hence, we can
1628  // veto here already.
1629  if ( doFxFx && nRequested < nJetMax && nCLjets > nRequested )
1630  return MORE_JETS;
1631 
1632  // Now continue in inclusive mode.
1633  // In inclusive mode, there can be more hadronic jets than partons,
1634  // provided that all partons are properly matched to hadronic jets.
1635  // Start by setting up the jet algorithm.
1636  if (!slowJet->setup(workEventJet) ) {
1637  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1638  "Light: the SlowJet algorithm failed on setup");
1639  return NONE;
1640  }
1641 
1642  // For FxFx, continue clustering as long as the jet separation is above
1643  // qCut.
1644  if (doFxFx) {
1645  while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1646  if( slowJet->dNext() > localQcutSq ) break;
1647  slowJet->doStep();
1648  }
1649  // For MLM, cluster into hadronic jets until there are the same number as
1650  // partons.
1651  } else {
1652  while ( slowJet->sizeAll() - slowJet->sizeJet() > nParton )
1653  slowJet->doStep();
1654  }
1655 
1656  // Sort partons in pT. Update local qCut value.
1657  // Hadronic jets are already sorted in pT.
1658  localQcutSq = dOld;
1659  if ( clFact >= 0. && nParton > 0 ) {
1660  vector<double> partonPt;
1661  for (int i = 0; i < nParton; ++i)
1662  partonPt.push_back( eventProcess[typeIdx[0][i]].pT2() );
1663  sort( partonPt.begin(), partonPt.end());
1664  localQcutSq = max( qCutSq, partonPt[0]);
1665  }
1666  nJets = slowJet->sizeJet();
1667  nClus = slowJet->sizeAll();
1668  }
1669  // Update scale if clustering factor is non-zero
1670  if ( clFact != 0. ) localQcutSq *= pow2(clFact);
1671 
1672  Event tempEvent;
1673  tempEvent.init( "(tempEvent)", particleDataPtr);
1674  int nPass = 0;
1675  double pTminEstimate = -1.;
1676  // Construct a master copy of the event containing only the
1677  // hardest nParton hadronic clusters. While constructing the event,
1678  // the parton type (ID_GLUON) and status (98,99) are arbitrary.
1679  for (int i = nJets; i < nClus; ++i) {
1680  tempEvent.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0, slowJet->p(i).px(),
1681  slowJet->p(i).py(), slowJet->p(i).pz(), slowJet->p(i).e() );
1682  ++nPass;
1683  pTminEstimate = max( pTminEstimate, slowJet->pT(i));
1684  if(nPass == nRequested) break;
1685  }
1686 
1687  int tempSize = tempEvent.size();
1688  // This keeps track of which hadronic jets are matched to parton
1689  vector<bool> jetAssigned;
1690  jetAssigned.assign( tempSize, false);
1691 
1692  // This keeps track of which partons are matched to which hadronic
1693  // jets.
1694  vector< vector<bool> > partonMatchesJet;
1695  for (int i=0; i < nParton; ++i )
1696  partonMatchesJet.push_back( vector<bool>(tempEvent.size(),false) );
1697 
1698  // Begin matching.
1699  // Do jet matching for FxFx.
1700  // Make sure that the nPartonsNow hardest hadronic jets are matched to any
1701  // of the nPartonsNow (+1) partons. This matching is done by attaching a jet
1702  // from the list of unmatched hadronic jets, and appending a jet from the
1703  // list of partonic jets, one at a time. The partonic jet will be clustered
1704  // with the hadronic jet or the beam if the distance measure is below the
1705  // cut. The hadronic jet is matched once this happens. Otherwise, another
1706  // partonic jet is tried. When a hadronic jet is matched to a partonic jet,
1707  // it is removed from the list of unmatched hadronic jets. This process
1708  // continues until the nPartonsNow hardest hadronic jets are matched to
1709  // partonic jets, or it is not possible to make a match for a hadronic jet.
1710  int iNow = 0;
1711  int nMatched = 0;
1712  while ( doFxFx && iNow < tempSize ) {
1713 
1714  // Check if this shower jet matches any partonic jet.
1715  Event tempEventJet;
1716  tempEventJet.init("(tempEventJet)", particleDataPtr);
1717  for (int i=0; i < nParton; ++i ) {
1718 
1720  //for (int j=0; j < tempSize; ++j )
1721  // if ( partonMatchesJet[i][j]) continue;
1722 
1723  // Attach a single hadronic jet.
1724  tempEventJet.clear();
1725  tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1726  tempEvent[iNow].px(), tempEvent[iNow].py(),
1727  tempEvent[iNow].pz(), tempEvent[iNow].e() );
1728  // Attach the current parton.
1729  Vec4 pIn = eventProcess[typeIdx[0][i]].p();
1730  tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1731  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1732 
1733  // Setup jet algorithm.
1734  if ( !slowJet->setup(tempEventJet) ) {
1735  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1736  "Light: the SlowJet algorithm failed on setup");
1737  return NONE;
1738  }
1739 
1740  // These are the conditions for the hadronic jet to match the parton
1741  // at the local qCut scale
1742  if ( slowJet->iNext() == tempEventJet.size() - 1
1743  && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1744  jetAssigned[iNow] = true;
1745  partonMatchesJet[i][iNow] = true;
1746  }
1747 
1748  } // End loop over hard partons.
1749 
1750  // Veto if the jet could not be assigned to any parton.
1751  if ( jetAssigned[iNow] ) nMatched++;
1752 
1753  // Continue;
1754  ++iNow;
1755  }
1756 
1757  // Jet matching veto for FxFx
1758  if (doFxFx) {
1759  if ( nRequested < nJetMax && nMatched != nRequested )
1760  return UNMATCHED_PARTON;
1761  if ( nRequested == nJetMax && nMatched < nRequested )
1762  return UNMATCHED_PARTON;
1763  }
1764 
1765  // Do jet matching for MLM.
1766  // Take the list of unmatched hadronic jets and append a parton, one at
1767  // a time. The parton will be clustered with the "closest" hadronic jet
1768  // or the beam if the distance measure is below the cut. When a hadronic
1769  // jet is matched to a parton, it is removed from the list of unmatched
1770  // hadronic jets. This process continues until all hadronic jets are
1771  // matched to partons or it is not possible to make a match.
1772  iNow = 0;
1773  while (!doFxFx && iNow < nParton ) {
1774  Event tempEventJet;
1775  tempEventJet.init("(tempEventJet)", particleDataPtr);
1776  for (int i = 0; i < tempSize; ++i) {
1777  if (jetAssigned[i]) continue;
1778  Vec4 pIn = tempEvent[i].p();
1779  // Append unmatched hadronic jets
1780  tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1781  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1782  }
1783 
1784  Vec4 pIn = eventProcess[typeIdx[0][iNow]].p();
1785  // Append the current parton
1786  tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1787  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1788  if ( !slowJet->setup(tempEventJet) ) {
1789  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1790  "Light: the SlowJet algorithm failed on setup");
1791  return NONE;
1792  }
1793  // These are the conditions for the hadronic jet to match the parton
1794  // at the local qCut scale
1795  if ( slowJet->iNext() == tempEventJet.size() - 1
1796  && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1797  int iKnt = -1;
1798  for (int i = 0; i != tempSize; ++i) {
1799  if (jetAssigned[i]) continue;
1800  ++iKnt;
1801  // Identify the hadronic jet that matches the parton
1802  if (iKnt == slowJet->jNext() ) jetAssigned[i] = true;
1803  }
1804  } else {
1805  return UNMATCHED_PARTON;
1806  }
1807  ++iNow;
1808  }
1809 
1810  // Minimal eT/pT (CellJet/SlowJet) of matched light jets.
1811  // Needed later for heavy jet vetos in inclusive mode.
1812  // This information is not used currently.
1813  if (nParton > 0 && pTminEstimate > 0) eTpTlightMin = pTminEstimate;
1814  else eTpTlightMin = -1.;
1815 
1816  // Record the jet separations.
1817  setDJR(workEventJet);
1818 
1819  // No veto
1820  return NONE;
1821 }
1822 
1823 //--------------------------------------------------------------------------
1824 
1825 // Step(2c): heavy jets
1826 // Return codes are given indicating the reason for a veto.
1827 // Although not currently used, they are a useful debugging tool:
1828 // 0 = no veto as there are no extra jets present
1829 // 1 = veto as in exclusive mode and extra jets present
1830 // 2 = veto as in inclusive mode and extra jets were harder
1831 // than any matched light jet
1832 
1833 inline int JetMatchingMadgraph::matchPartonsToJetsHeavy() {
1834 
1835  // Currently, heavy jets are unmatched
1836  // If there are no extra jets, then accept
1837  // jetMomenta is NEVER used by MadGraph and is always empty.
1838  // This check does nothing.
1839  // Rather, if there is any heavy flavor that is harder than
1840  // what is present at the LHE level, then the event should
1841  // be vetoed.
1842 
1843  // if (jetMomenta.empty()) return NONE;
1844  // Count the number of hard partons
1845  int nParton = typeIdx[1].size();
1846 
1847  Event tempEventJet(workEventJet);
1848 
1849  double scaleF(1.0);
1850  // Rescale the heavy partons that are from the hard process to
1851  // have pT=collider energy. Soft/collinear gluons will cluster
1852  // onto them, leaving a remnant of hard emissions.
1853  for( int i=0; i<nParton; ++i) {
1854  scaleF = eventProcessOrig[0].e()/workEventJet[typeIdx[1][i]].pT();
1855  tempEventJet[typeIdx[1][i]].rescale5(scaleF);
1856  }
1857 
1858  if (!hjSlowJet->setup(tempEventJet) ) {
1859  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1860  "Heavy: the SlowJet algorithm failed on setup");
1861  return NONE;
1862  }
1863 
1864 
1865  while ( hjSlowJet->sizeAll() - hjSlowJet->sizeJet() > 0 ) {
1866  if( hjSlowJet->dNext() > qCutSq ) break;
1867  hjSlowJet->doStep();
1868  }
1869 
1870  int nCLjets(0);
1871  // Count the number of clusters with pT>qCut. This includes the
1872  // original hard partons plus any hard emissions.
1873  for(int idx=0 ; idx< hjSlowJet->sizeAll(); ++idx) {
1874  if( hjSlowJet->pT(idx) > sqrt(qCutSq) ) nCLjets++;
1875  }
1876 
1877  // Debug printout.
1878  if (MATCHINGDEBUG) hjSlowJet->list(true);
1879 
1880  // Count of the number of hadronic jets in SlowJet accounting
1881  // int nCLjets = nClus - nJets;
1882  // Get number of partons. Different for MLM and FxFx schemes.
1883  int nRequested = nParton;
1884 
1885  // Veto event if too few hadronic jets
1886  if ( nCLjets < nRequested ) {
1887  if (MATCHINGDEBUG) cout << "veto : hvy LESS_JETS " << endl;
1888  if (MATCHINGDEBUG) cout << "nCLjets = " << nCLjets << "; nRequest = "
1889  << nRequested << endl;
1890  return LESS_JETS;
1891  }
1892 
1893  // In exclusive mode, do not allow more hadronic jets than partons
1894  if ( exclusive ) {
1895  if ( nCLjets > nRequested ) {
1896  if (MATCHINGDEBUG) cout << "veto : excl hvy MORE_JETS " << endl;
1897  return MORE_JETS;
1898  }
1899  }
1900 
1901  // No extra jets were present so no veto
1902  return NONE;
1903 }
1904 
1905 //--------------------------------------------------------------------------
1906 
1907 // Step(2c): other jets
1908 // Return codes are given indicating the reason for a veto.
1909 // Although not currently used, they are a useful debugging tool:
1910 // 0 = no veto as there are no extra jets present
1911 // 1 = veto as in exclusive mode and extra jets present
1912 // 2 = veto as in inclusive mode and extra jets were harder
1913 // than any matched light jet
1914 
1915 inline int JetMatchingMadgraph::matchPartonsToJetsOther() {
1916 
1917  // Currently, heavy jets are unmatched
1918  // If there are no extra jets, then accept
1919  // jetMomenta is NEVER used by MadGraph and is always empty.
1920  // This check does nothing.
1921  // Rather, if there is any heavy flavor that is harder than
1922  // what is present at the LHE level, then the event should
1923  // be vetoed.
1924 
1925  // if (jetMomenta.empty()) return NONE;
1926  // Count the number of hard partons
1927  int nParton = typeIdx[2].size();
1928 
1929  Event tempEventJet(workEventJet);
1930 
1931  double scaleF(1.0);
1932  // Rescale the heavy partons that are from the hard process to
1933  // have pT=collider energy. Soft/collinear gluons will cluster
1934  // onto them, leaving a remnant of hard emissions.
1935  for( int i=0; i<nParton; ++i) {
1936  scaleF = eventProcessOrig[0].e()/workEventJet[typeIdx[2][i]].pT();
1937  tempEventJet[typeIdx[2][i]].rescale5(scaleF);
1938  }
1939 
1940  if (!hjSlowJet->setup(tempEventJet) ) {
1941  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1942  "Heavy: the SlowJet algorithm failed on setup");
1943  return NONE;
1944  }
1945 
1946 
1947  while ( hjSlowJet->sizeAll() - hjSlowJet->sizeJet() > 0 ) {
1948  if( hjSlowJet->dNext() > qCutSq ) break;
1949  hjSlowJet->doStep();
1950  }
1951 
1952  int nCLjets(0);
1953  // Count the number of clusters with pT>qCut. This includes the
1954  // original hard partons plus any hard emissions.
1955  for(int idx=0 ; idx< hjSlowJet->sizeAll(); ++idx) {
1956  if( hjSlowJet->pT(idx) > sqrt(qCutSq) ) nCLjets++;
1957  }
1958 
1959  // Debug printout.
1960  if (MATCHINGDEBUG) hjSlowJet->list(true);
1961 
1962  // Count of the number of hadronic jets in SlowJet accounting
1963  // int nCLjets = nClus - nJets;
1964  // Get number of partons. Different for MLM and FxFx schemes.
1965  int nRequested = nParton;
1966 
1967  // Veto event if too few hadronic jets
1968  if ( nCLjets < nRequested ) {
1969  if (MATCHINGDEBUG) cout << "veto : other LESS_JETS " << endl;
1970  if (MATCHINGDEBUG) cout << "nCLjets = " << nCLjets << "; nRequest = "
1971  << nRequested << endl;
1972  return LESS_JETS;
1973  }
1974 
1975  // In exclusive mode, do not allow more hadronic jets than partons
1976  if ( exclusive ) {
1977  if ( nCLjets > nRequested ) {
1978  if (MATCHINGDEBUG) cout << "veto : excl other MORE_JETS" << endl;
1979  return MORE_JETS;
1980  }
1981  }
1982 
1983  // No extra jets were present so no veto
1984  return NONE;
1985 }
1986 
1987 //==========================================================================
1988 
1989 } // end namespace Pythia8
1990 
1991 #endif // end Pythia8_JetMatching_H
Definition: AgUStep.h:26