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) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, 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 and Stefan Prestel (implementation of shower-kT
11 // MLM-style matching and flavour treatment for Madgraph input, and FxFx NLO
12 // 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 "GeneratorInput.h"
24 using namespace Pythia8;
25 
26 //==========================================================================
27 
28 // Declaration of main JetMatching class to perform MLM matching.
29 // Note that it is defined with virtual inheritance, so that it can
30 // be combined with other UserHooks classes, see e.g. main33.cc.
31 
32 class JetMatching : virtual public UserHooks {
33 
34 public:
35 
36  // Constructor and destructor
37  JetMatching() : cellJet(NULL), slowJet(NULL), slowJetHard(NULL) {}
38  ~JetMatching() {
39  if (cellJet) delete cellJet;
40  if (slowJet) delete slowJet;
41  if (slowJetHard) delete slowJetHard;
42  }
43 
44  // Initialisation
45  virtual bool initAfterBeams() = 0;
46 
47  // Process level vetos
48  bool canVetoProcessLevel() { return doMerge; }
49  bool doVetoProcessLevel(Event& process) {
50  eventProcessOrig = process;
51  return false;
52  }
53 
54  // Parton level vetos (before beam remnants and resonance decays)
55  bool canVetoPartonLevelEarly() { return doMerge; }
56  bool doVetoPartonLevelEarly(const Event& event);
57 
58  // Shower step vetoes (after the first emission, for Shower-kT scheme)
59  int numberVetoStep() {return 1;}
60  bool canVetoStep() { return false; }
61  bool doVetoStep(int, int, int, const Event& ) { return false; }
62 
63 protected:
64 
65  // Constants to be changed for debug printout or extra checks.
66  static const bool MATCHINGDEBUG, MATCHINGCHECK;
67 
68  // Different steps of the matching algorithm.
69  virtual void sortIncomingProcess(const Event &)=0;
70  virtual void jetAlgorithmInput(const Event &, int)=0;
71  virtual void runJetAlgorithm()=0;
72  virtual bool matchPartonsToJets(int)=0;
73  virtual int matchPartonsToJetsLight()=0;
74  virtual int matchPartonsToJetsHeavy()=0;
75 
76  enum vetoStatus { NONE, LESS_JETS, MORE_JETS, HARD_JET, UNMATCHED_PARTON };
77  enum partonTypes { ID_CHARM=4, ID_BOT=5, ID_TOP=6, ID_LEPMIN=11,
78  ID_LEPMAX=16, ID_GLUON=21, ID_PHOTON=22 };
79 
80  // Master switch for merging
81  bool doMerge;
82  // Switch for merging in the shower-kT scheme. Needed here because
83  // the scheme uses different UserHooks functionality.
84  bool doShowerKt;
85 
86  // Maximum and current number of jets
87  int nJetMax, nJet;
88 
89  // Jet algorithm parameters
90  int jetAlgorithm;
91  double eTjetMin, coneRadius, etaJetMax, etaJetMaxAlgo;
92 
93  // Internal jet algorithms
94  CellJet* cellJet;
95  SlowJet* slowJet;
96  SlowJet* slowJetHard;
97 
98  // SlowJet specific
99  int slowJetPower;
100 
101  // Event records to store original incoming process, final-state of the
102  // incoming process and what will be passed to the jet algorithm.
103  // Not completely necessary to store all steps, but makes tracking the
104  // steps of the algorithm a lot easier.
105  Event eventProcessOrig, eventProcess, workEventJet;
106 
107  // Sort final-state of incoming process into light/heavy jets and 'other'
108  vector<int> typeIdx[3];
109  set<int> typeSet[3];
110 
111  // Momenta output of jet algorithm (to provide same output regardless of
112  // the selected jet algorithm)
113  vector<Vec4> jetMomenta;
114 
115  // CellJet specific
116  int nEta, nPhi;
117  double eTseed, eTthreshold;
118 
119  // Merging procedure parameters
120  int jetAllow, jetMatch, exclusiveMode;
121  double coneMatchLight, coneRadiusHeavy, coneMatchHeavy;
122  bool exclusive;
123 
124  // Store the minimum eT/pT of matched light jets
125  double eTpTlightMin;
126 
127 };
128 
129 //==========================================================================
130 
131 // Declaration of main UserHooks class to perform Alpgen matching.
132 
133 class JetMatchingAlpgen : virtual public JetMatching {
134 
135 public:
136 
137  // Constructor and destructor
138  JetMatchingAlpgen() { }
139  ~JetMatchingAlpgen() { }
140 
141  // Initialisation
142  bool initAfterBeams();
143 
144 private:
145 
146  // Different steps of the matching algorithm.
147  void sortIncomingProcess(const Event &);
148  void jetAlgorithmInput(const Event &, int);
149  void runJetAlgorithm();
150  bool matchPartonsToJets(int);
151  int matchPartonsToJetsLight();
152  int matchPartonsToJetsHeavy();
153 
154  // Sorting utility
155  void sortTypeIdx(vector < int > &vecIn);
156 
157  // Constants
158  static const double GHOSTENERGY, ZEROTHRESHOLD;
159 
160 };
161 
162 //==========================================================================
163 
164 // Declaration of main UserHooks class to perform Madgraph matching.
165 
166 class JetMatchingMadgraph : virtual public JetMatching {
167 
168 public:
169 
170  // Constructor and destructor
171  JetMatchingMadgraph() { }
172  ~JetMatchingMadgraph() { }
173 
174  // Initialisation
175  bool initAfterBeams();
176 
177  // Shower step vetoes (after the first emission, for Shower-kT scheme)
178  int numberVetoStep() {return 1;}
179  bool canVetoStep() { return doShowerKt; }
180  bool doVetoStep(int, int, int, const Event& );
181 
182 protected:
183 
184  // Different steps of the matching algorithm.
185  void sortIncomingProcess(const Event &);
186  void jetAlgorithmInput(const Event &, int);
187  void runJetAlgorithm();
188  bool matchPartonsToJets(int);
189  int matchPartonsToJetsLight();
190  int matchPartonsToJetsHeavy();
191  bool doShowerKtVeto(double pTfirst);
192 
193  // Variables.
194  vector<int> origTypeIdx[3];
195  int nQmatch;
196  double qCut, qCutSq, clFact;
197  bool doFxFx;
198  int nPartonsNow;
199  double qCutME, qCutMESq;
200 
201 };
202 
203 //==========================================================================
204 
205 // Main implementation of JetMatching class.
206 // This may be split out to a separate C++ file if desired,
207 // but currently included here for ease of use.
208 
209 //--------------------------------------------------------------------------
210 
211 // Constants to be changed for debug printout or extra checks.
212 const bool JetMatching::MATCHINGDEBUG = false;
213 const bool JetMatching::MATCHINGCHECK = false;
214 
215 //--------------------------------------------------------------------------
216 
217 // Early parton level veto (before beam remnants and resonance showers)
218 
219 bool JetMatching::doVetoPartonLevelEarly(const Event& event) {
220 
221  // 1) Sort the original incoming process. After this step is performed,
222  // the following assignments have been made:
223  // eventProcessOrig - the original incoming process
224  // eventProcess - the final-state of the incoming process with
225  // resonance decays removed (and resonances
226  // themselves now with positive status code)
227  // typeIdx[0/1/2] - Indices into 'eventProcess' of
228  // light jets/heavy jets/other
229  // typeSet[0/1/2] - Indices into 'event' of light jets/heavy jets/other
230  // workEvent - partons from the hardest subsystem + ISR + FSR only
231  sortIncomingProcess(event);
232 
233  // For the shower-kT scheme, do not perform any veto here, as any vetoing
234  // will already have taken place in doVetoStep.
235  if ( doShowerKt ) return false;
236 
237  // Debug printout.
238  if (MATCHINGDEBUG) {
239  // Begin
240  cout << endl << "-------- Begin Madgraph Debug --------" << endl;
241  // Original incoming process
242  cout << endl << "Original incoming process:";
243  eventProcessOrig.list();
244  // Final-state of original incoming process
245  cout << endl << "Final-state incoming process:";
246  eventProcess.list();
247  // List categories of sorted particles
248  for (size_t i = 0; i < typeIdx[0].size(); i++)
249  cout << ((i == 0) ? "Light jets: " : ", ") << setw(3) << typeIdx[0][i];
250  if( typeIdx[0].size()== 0 )
251  cout << "Light jets: None";
252 
253  for (size_t i = 0; i < typeIdx[1].size(); i++)
254  cout << ((i == 0) ? "\nHeavy jets: " : ", ") << setw(3) << typeIdx[1][i];
255  for (size_t i = 0; i < typeIdx[2].size(); i++)
256  cout << ((i == 0) ? "\nOther: " : ", ") << setw(3) << typeIdx[2][i];
257  // Full event at this stage
258  cout << endl << endl << "Event:";
259  event.list();
260  // Work event (partons from hardest subsystem + ISR + FSR)
261  cout << endl << "Work event:";
262  workEvent.list();
263  }
264 
265  // 2) Light/heavy jets: iType = 0 (light jets), 1 (heavy jets)
266  int iTypeEnd = (typeIdx[1].empty()) ? 1 : 2;
267  for (int iType = 0; iType < iTypeEnd; iType++) {
268 
269  // 2a) Find particles which will be passed from the jet algorithm.
270  // Input from 'workEvent' and output in 'workEventJet'.
271  jetAlgorithmInput(event, iType);
272 
273  // Debug printout.
274  if (MATCHINGDEBUG) {
275  // Jet algorithm event
276  cout << endl << "Jet algorithm event (iType = " << iType << "):";
277  workEventJet.list();
278  }
279 
280  // 2b) Run jet algorithm on 'workEventJet'.
281  // Output is stored in jetMomenta.
282  runJetAlgorithm();
283 
284  // 2c) Match partons to jets and decide if veto is necessary
285  if (matchPartonsToJets(iType) == true) {
286  // Debug printout.
287  if (MATCHINGDEBUG) {
288  cout << endl << "Event vetoed" << endl
289  << "---------- End MLM Debug ----------" << endl;
290  }
291  return true;
292  }
293  }
294 
295  // Debug printout.
296  if (MATCHINGDEBUG) {
297  cout << endl << "Event accepted" << endl
298  << "---------- End MLM Debug ----------" << endl;
299  }
300 
301  // If we reached here, then no veto
302  return false;
303 
304 }
305 
306 //==========================================================================
307 
308 // Main implementation of Alpgen UserHooks class.
309 // This may be split out to a separate C++ file if desired,
310 // but currently included here for ease of use.
311 
312 //--------------------------------------------------------------------------
313 
314 // Constants: could be changed here if desired, but normally should not.
315 // These are of technical nature, as described for each.
316 
317 // The energy of ghost particles. For technical reasons, this cannot be
318 // set arbitrarily low, see 'Particle::TINY' in 'Event.cc' for details.
319 const double JetMatchingAlpgen::GHOSTENERGY = 1e-15;
320 
321 // A zero threshold value for double comparisons.
322 const double JetMatchingAlpgen::ZEROTHRESHOLD = 1e-10;
323 
324 //--------------------------------------------------------------------------
325 
326 // Function to sort typeIdx vectors into descending eT/pT order.
327 // Uses a selection sort, as number of partons generally small
328 // and so efficiency not a worry.
329 
330 void JetMatchingAlpgen::sortTypeIdx(vector < int > &vecIn) {
331  for (size_t i = 0; i < vecIn.size(); i++) {
332  size_t jMax = i;
333  double vMax = (jetAlgorithm == 1) ?
334  eventProcess[vecIn[i]].eT() :
335  eventProcess[vecIn[i]].pT();
336  for (size_t j = i + 1; j < vecIn.size(); j++) {
337  double vNow = (jetAlgorithm == 1)
338  ? eventProcess[vecIn[j]].eT() : eventProcess[vecIn[j]].pT();
339  if (vNow > vMax) {
340  vMax = vNow;
341  jMax = j;
342  }
343  }
344  if (jMax != i) swap(vecIn[i], vecIn[jMax]);
345  }
346 }
347 
348 //--------------------------------------------------------------------------
349 
350 // Initialisation routine automatically called from Pythia::init().
351 // Setup all parts needed for the merging.
352 
353 bool JetMatchingAlpgen::initAfterBeams() {
354 
355  // Read in parameters
356  doMerge = settingsPtr->flag("JetMatching:merge");
357  jetAlgorithm = settingsPtr->mode("JetMatching:jetAlgorithm");
358  nJet = settingsPtr->mode("JetMatching:nJet");
359  nJetMax = settingsPtr->mode("JetMatching:nJetMax");
360  eTjetMin = settingsPtr->parm("JetMatching:eTjetMin");
361  coneRadius = settingsPtr->parm("JetMatching:coneRadius");
362  etaJetMax = settingsPtr->parm("JetMatching:etaJetMax");
363  doShowerKt = settingsPtr->flag("JetMatching:doShowerKt");
364 
365  // Use etaJetMax + coneRadius in input to jet algorithms
366  etaJetMaxAlgo = etaJetMax + coneRadius;
367 
368  // CellJet specific
369  nEta = settingsPtr->mode("JetMatching:nEta");
370  nPhi = settingsPtr->mode("JetMatching:nPhi");
371  eTseed = settingsPtr->parm("JetMatching:eTseed");
372  eTthreshold = settingsPtr->parm("JetMatching:eTthreshold");
373 
374  // SlowJet specific
375  slowJetPower = settingsPtr->mode("JetMatching:slowJetPower");
376  coneMatchLight = settingsPtr->parm("JetMatching:coneMatchLight");
377  coneRadiusHeavy = settingsPtr->parm("JetMatching:coneRadiusHeavy");
378  if (coneRadiusHeavy < 0.) coneRadiusHeavy = coneRadius;
379  coneMatchHeavy = settingsPtr->parm("JetMatching:coneMatchHeavy");
380 
381  // Matching procedure
382  jetAllow = settingsPtr->mode("JetMatching:jetAllow");
383  jetMatch = settingsPtr->mode("JetMatching:jetMatch");
384  exclusiveMode = settingsPtr->mode("JetMatching:exclusive");
385 
386  // If not merging, then done
387  if (!doMerge) return true;
388 
389  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
390  if (exclusiveMode == 2) {
391 
392  // No nJet or nJetMax, so default to exclusive mode
393  if (nJet < 0 || nJetMax < 0) {
394  infoPtr->errorMsg("Warning in JetMatchingAlpgen:init: "
395  "missing jet multiplicity information; running in exclusive mode");
396  exclusive = true;
397 
398  // Inclusive if nJet == nJetMax, exclusive otherwise
399  } else {
400  exclusive = (nJet == nJetMax) ? false : true;
401  }
402 
403  // Otherwise, just set as given
404  } else {
405  exclusive = (exclusiveMode == 0) ? false : true;
406  }
407 
408  // Initialise chosen jet algorithm. CellJet.
409  if (jetAlgorithm == 1) {
410 
411  // Extra options for CellJet. nSel = 1 means that all final-state
412  // particles are taken and we retain control of what to select.
413  // smear/resolution/upperCut are not used and are set to default values.
414  int nSel = 2, smear = 0;
415  double resolution = 0.5, upperCut = 2.;
416  cellJet = new CellJet(etaJetMaxAlgo, nEta, nPhi, nSel,
417  smear, resolution, upperCut, eTthreshold);
418 
419  // SlowJet
420  } else if (jetAlgorithm == 2) {
421  slowJet = new SlowJet(slowJetPower, coneRadius, eTjetMin, etaJetMaxAlgo);
422  }
423 
424  // Check the jetMatch parameter; option 2 only works with SlowJet
425  if (jetAlgorithm == 1 && jetMatch == 2) {
426  infoPtr->errorMsg("Warning in JetMatchingAlpgen:init: "
427  "jetMatch = 2 only valid with SlowJet algorithm. "
428  "Reverting to jetMatch = 1");
429  jetMatch = 1;
430  }
431 
432  // Setup local event records
433  eventProcessOrig.init("(eventProcessOrig)", particleDataPtr);
434  eventProcess.init("(eventProcess)", particleDataPtr);
435  workEventJet.init("(workEventJet)", particleDataPtr);
436 
437  // Print information
438  string jetStr = (jetAlgorithm == 1) ? "CellJet" :
439  (slowJetPower == -1) ? "anti-kT" :
440  (slowJetPower == 0) ? "C/A" :
441  (slowJetPower == 1) ? "kT" : "unknown";
442  string modeStr = (exclusive) ? "exclusive" : "inclusive";
443  stringstream nJetStr, nJetMaxStr;
444  if (nJet >= 0) nJetStr << nJet; else nJetStr << "unknown";
445  if (nJetMax >= 0) nJetMaxStr << nJetMax; else nJetMaxStr << "unknown";
446  cout << endl
447  << " *------- MLM matching parameters -------*" << endl
448  << " | nJet | " << setw(14)
449  << nJetStr.str() << " |" << endl
450  << " | nJetMax | " << setw(14)
451  << nJetMaxStr.str() << " |" << endl
452  << " | Jet algorithm | " << setw(14)
453  << jetStr << " |" << endl
454  << " | eTjetMin | " << setw(14)
455  << eTjetMin << " |" << endl
456  << " | coneRadius | " << setw(14)
457  << coneRadius << " |" << endl
458  << " | etaJetMax | " << setw(14)
459  << etaJetMax << " |" << endl
460  << " | jetAllow | " << setw(14)
461  << jetAllow << " |" << endl
462  << " | jetMatch | " << setw(14)
463  << jetMatch << " |" << endl
464  << " | coneMatchLight | " << setw(14)
465  << coneMatchLight << " |" << endl
466  << " | coneRadiusHeavy | " << setw(14)
467  << coneRadiusHeavy << " |" << endl
468  << " | coneMatchHeavy | " << setw(14)
469  << coneMatchHeavy << " |" << endl
470  << " | Mode | " << setw(14)
471  << modeStr << " |" << endl
472  << " *-----------------------------------------*" << endl;
473 
474  return true;
475 }
476 
477 //--------------------------------------------------------------------------
478 
479 // Step (1): sort the incoming particles
480 
481 void JetMatchingAlpgen::sortIncomingProcess(const Event &event) {
482 
483  // Remove resonance decays from original process and keep only final
484  // state. Resonances will have positive status code after this step.
485  omitResonanceDecays(eventProcessOrig, true);
486  eventProcess = workEvent;
487 
488  // Sort original process final state into light/heavy jets and 'other'.
489  // Criteria:
490  // 1 <= ID <= 5 and massless, or ID == 21 --> light jet (typeIdx[0])
491  // 4 <= ID <= 6 and massive --> heavy jet (typeIdx[1])
492  // All else --> other (typeIdx[2])
493  // Note that 'typeIdx' stores indices into 'eventProcess' (after resonance
494  // decays are omitted), while 'typeSet' stores indices into the original
495  // process record, 'eventProcessOrig', but these indices are also valid
496  // in 'event'.
497  for (int i = 0; i < 3; i++) {
498  typeIdx[i].clear();
499  typeSet[i].clear();
500  }
501  for (int i = 0; i < eventProcess.size(); i++) {
502  // Ignore nonfinal and default to 'other'
503  if (!eventProcess[i].isFinal()) continue;
504  int idx = 2;
505 
506  // Light jets
507  if (eventProcess[i].id() == ID_GLUON
508  || (eventProcess[i].idAbs() <= ID_BOT
509  && abs(eventProcess[i].m()) < ZEROTHRESHOLD)) idx = 0;
510 
511  // Heavy jets
512  else if (eventProcess[i].idAbs() >= ID_CHARM
513  && eventProcess[i].idAbs() <= ID_TOP) idx = 1;
514 
515  // Store
516  typeIdx[idx].push_back(i);
517  typeSet[idx].insert(eventProcess[i].daughter1());
518  }
519 
520  // Extract partons from hardest subsystem + ISR + FSR only into
521  // workEvent. Note no resonance showers or MPIs.
522  subEvent(event);
523 }
524 
525 //--------------------------------------------------------------------------
526 
527 // Step (2a): pick which particles to pass to the jet algorithm
528 
529 void JetMatchingAlpgen::jetAlgorithmInput(const Event &event, int iType) {
530 
531  // Take input from 'workEvent' and put output in 'workEventJet'
532  workEventJet = workEvent;
533 
534  // Loop over particles and decide what to pass to the jet algorithm
535  for (int i = 0; i < workEventJet.size(); ++i) {
536  if (!workEventJet[i].isFinal()) continue;
537 
538  // jetAllow option to disallow certain particle types
539  if (jetAllow == 1) {
540 
541  // Original AG+Py6 algorithm explicitly excludes tops,
542  // leptons and photons.
543  int id = workEventJet[i].idAbs();
544  if ( (id >= ID_LEPMIN && id <= ID_LEPMAX) || id == ID_TOP
545  || id == ID_PHOTON) {
546  workEventJet[i].statusNeg();
547  continue;
548  }
549  }
550 
551  // Get the index of this particle in original event
552  int idx = workEventJet[i].daughter1();
553 
554  // Start with particle idx, and afterwards track mothers
555  while (true) {
556 
557  // Light jets
558  if (iType == 0) {
559 
560  // Do not include if originates from heavy jet or 'other'
561  if (typeSet[1].find(idx) != typeSet[1].end() ||
562  typeSet[2].find(idx) != typeSet[2].end()) {
563  workEventJet[i].statusNeg();
564  break;
565  }
566 
567  // Made it to start of event record so done
568  if (idx == 0) break;
569  // Otherwise next mother and continue
570  idx = event[idx].mother1();
571 
572  // Heavy jets
573  } else if (iType == 1) {
574 
575  // Only include if originates from heavy jet
576  if (typeSet[1].find(idx) != typeSet[1].end()) break;
577 
578  // Made it to start of event record with no heavy jet mother,
579  // so DO NOT include particle
580  if (idx == 0) {
581  workEventJet[i].statusNeg();
582  break;
583  }
584 
585  // Otherwise next mother and continue
586  idx = event[idx].mother1();
587 
588  } // if (iType)
589  } // while (true)
590  } // for (i)
591 
592  // For jetMatch = 2, insert ghost particles corresponding to
593  // each hard parton in the original process
594  if (jetMatch == 2) {
595  for (int i = 0; i < int(typeIdx[iType].size()); i++) {
596  // Get y/phi of the parton
597  Vec4 pIn = eventProcess[typeIdx[iType][i]].p();
598  double y = pIn.rap();
599  double phi = pIn.phi();
600 
601  // Create a ghost particle and add to the workEventJet
602  double e = GHOSTENERGY;
603  double e2y = exp(2. * y);
604  double pz = e * (e2y - 1.) / (e2y + 1.);
605  double pt = sqrt(e*e - pz*pz);
606  double px = pt * cos(phi);
607  double py = pt * sin(phi);
608  workEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0, px, py, pz, e);
609 
610  // Extra check on reconstructed y/phi values. If many warnings
611  // of this type, GHOSTENERGY may be set too low.
612  if (MATCHINGCHECK) {
613  int lastIdx = workEventJet.size() - 1;
614  if (abs(y - workEventJet[lastIdx].y()) > ZEROTHRESHOLD ||
615  abs(phi - workEventJet[lastIdx].phi()) > ZEROTHRESHOLD)
616  infoPtr->errorMsg("Warning in JetMatchingAlpgen:jetAlgorithmInput: "
617  "ghost particle y/phi mismatch");
618  }
619 
620  } // for (i)
621  } // if (jetMatch == 2)
622 }
623 
624 //--------------------------------------------------------------------------
625 
626 // Step (2b): run jet algorithm and provide common output
627 
628 void JetMatchingAlpgen::runJetAlgorithm() {
629 
630  // Run the jet clustering algorithm
631  if (jetAlgorithm == 1)
632  cellJet->analyze(workEventJet, eTjetMin, coneRadius, eTseed);
633  else
634  slowJet->analyze(workEventJet);
635 
636  // Extract four-momenta of jets with |eta| < etaJetMax and
637  // put into jetMomenta. Note that this is done backwards as
638  // jets are removed with SlowJet.
639  jetMomenta.clear();
640  int iJet = (jetAlgorithm == 1) ? cellJet->size() - 1:
641  slowJet->sizeJet() - 1;
642  for (int i = iJet; i > -1; i--) {
643  Vec4 jetMom = (jetAlgorithm == 1) ? cellJet->pMassive(i) :
644  slowJet->p(i);
645  double eta = jetMom.eta();
646 
647  if (abs(eta) > etaJetMax) {
648  if (jetAlgorithm == 2) slowJet->removeJet(i);
649  continue;
650  }
651  jetMomenta.push_back(jetMom);
652  }
653 
654  // Reverse jetMomenta to restore eT/pT ordering
655  reverse(jetMomenta.begin(), jetMomenta.end());
656 }
657 
658 //--------------------------------------------------------------------------
659 
660 // Step (2c): veto decision (returning true vetoes the event)
661 
662 bool JetMatchingAlpgen::matchPartonsToJets(int iType) {
663 
664  // Use two different routines for light/heavy jets as
665  // different veto conditions and for clarity
666  if (iType == 0) return (matchPartonsToJetsLight() > 0);
667  else return (matchPartonsToJetsHeavy() > 0);
668 }
669 
670 //--------------------------------------------------------------------------
671 
672 // Step(2c): light jets
673 // Return codes are given indicating the reason for a veto.
674 // Although not currently used, they are a useful debugging tool:
675 // 0 = no veto
676 // 1 = veto as number of jets less than number of partons
677 // 2 = veto as exclusive mode and number of jets greater than
678 // number of partons
679 // 3 = veto as inclusive mode and there would be an extra jet
680 // that is harder than any matched soft jet
681 // 4 = veto as there is a parton which does not match a jet
682 
683 int JetMatchingAlpgen::matchPartonsToJetsLight() {
684 
685  // Always veto if number of jets is less than original number of jets
686  if (jetMomenta.size() < typeIdx[0].size()) return LESS_JETS;
687  // Veto if in exclusive mode and number of jets bigger than original
688  if (exclusive && jetMomenta.size() > typeIdx[0].size()) return MORE_JETS;
689 
690  // Sort partons by eT/pT
691  sortTypeIdx(typeIdx[0]);
692 
693  // Number of hard partons
694  int nParton = typeIdx[0].size();
695 
696  // Keep track of which jets have been assigned a hard parton
697  vector < bool > jetAssigned;
698  jetAssigned.assign(jetMomenta.size(), false);
699 
700  // Jet matching procedure: (1) deltaR between partons and jets
701  if (jetMatch == 1) {
702 
703  // Loop over light hard partons and get 4-momentum
704  for (int i = 0; i < nParton; i++) {
705  Vec4 p1 = eventProcess[typeIdx[0][i]].p();
706 
707  // Track which jet has the minimal dR measure with this parton
708  int jMin = -1;
709  double dRmin = 0.;
710 
711  // Loop over all jets (skipping those already assigned).
712  for (int j = 0; j < int(jetMomenta.size()); j++) {
713  if (jetAssigned[j]) continue;
714 
715  // DeltaR between parton/jet and store if minimum
716  double dR = (jetAlgorithm == 1)
717  ? REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
718  if (jMin < 0 || dR < dRmin) {
719  dRmin = dR;
720  jMin = j;
721  }
722  } // for (j)
723 
724  // Check for jet-parton match
725  if (jMin >= 0 && dRmin < coneRadius * coneMatchLight) {
726 
727  // If the matched jet is not one of the nParton hardest jets,
728  // the extra left over jet would be harder than some of the
729  // matched jets. This is disallowed, so veto.
730  if (jMin >= nParton) return HARD_JET;
731 
732  // Mark jet as assigned.
733  jetAssigned[jMin] = true;
734 
735  // If no match, then event will be vetoed in all cases
736  } else return UNMATCHED_PARTON;
737 
738  } // for (i)
739 
740  // Jet matching procedure: (2) ghost particles in SlowJet
741  } else {
742 
743  // Loop over added 'ghost' particles and find if assigned to a jet
744  for (int i = workEventJet.size() - nParton;
745  i < workEventJet.size(); i++) {
746  int jMin = slowJet->jetAssignment(i);
747 
748  // Veto if:
749  // 1) not one of nParton hardest jets
750  // 2) not assigned to a jet
751  // 3) jet has already been assigned
752  if (jMin >= nParton) return HARD_JET;
753  if (jMin < 0 || jetAssigned[jMin]) return UNMATCHED_PARTON;
754 
755  // Mark jet as assigned
756  jetAssigned[jMin] = true;
757 
758  } // for (i)
759  } // if (jetMatch)
760 
761  // Minimal eT/pT (CellJet/SlowJet) of matched light jets. Needed
762  // later for heavy jet vetos in inclusive mode.
763  if (nParton > 0)
764  eTpTlightMin = (jetAlgorithm == 1) ? jetMomenta[nParton - 1].eT()
765  : jetMomenta[nParton - 1].pT();
766  else
767  eTpTlightMin = -1.;
768 
769  // No veto
770  return NONE;
771 }
772 
773 //--------------------------------------------------------------------------
774 
775 // Step(2c): heavy jets
776 // Return codes are given indicating the reason for a veto.
777 // Although not currently used, they are a useful debugging tool:
778 // 0 = no veto as there are no extra jets present
779 // 1 = veto as in exclusive mode and extra jets present
780 // 2 = veto as in inclusive mode and extra jets were harder
781 // than any matched light jet
782 
783 int JetMatchingAlpgen::matchPartonsToJetsHeavy() {
784 
785  // If there are no extra jets, then accept
786  if (jetMomenta.empty()) return NONE;
787 
788  // Number of hard partons
789  int nParton = typeIdx[1].size();
790 
791  // Remove jets that are close to heavy quarks
792  set < int > removeJets;
793 
794  // Jet matching procedure: (1) deltaR between partons and jets
795  if (jetMatch == 1) {
796 
797  // Loop over heavy hard partons and get 4-momentum
798  for (int i = 0; i < nParton; i++) {
799  Vec4 p1 = eventProcess[typeIdx[1][i]].p();
800 
801  // Loop over all jets, find dR and mark for removal if match
802  for (int j = 0; j < int(jetMomenta.size()); j++) {
803  double dR = (jetAlgorithm == 1) ?
804  REtaPhi(p1, jetMomenta[j]) : RRapPhi(p1, jetMomenta[j]);
805  if (dR < coneRadiusHeavy * coneMatchHeavy)
806  removeJets.insert(j);
807 
808  } // for (j)
809  } // for (i)
810 
811  // Jet matching procedure: (2) ghost particles in SlowJet
812  } else {
813 
814  // Loop over added 'ghost' particles and if assigned to a jet
815  // then mark this jet for removal
816  for (int i = workEventJet.size() - nParton;
817  i < workEventJet.size(); i++) {
818  int jMin = slowJet->jetAssignment(i);
819  if (jMin >= 0) removeJets.insert(jMin);
820  }
821 
822  }
823 
824  // Remove jets (backwards order to not disturb indices)
825  for (set < int >::reverse_iterator it = removeJets.rbegin();
826  it != removeJets.rend(); it++)
827  jetMomenta.erase(jetMomenta.begin() + *it);
828 
829  // Handle case if there are still extra jets
830  if (!jetMomenta.empty()) {
831 
832  // Exclusive mode, so immediate veto
833  if (exclusive) return MORE_JETS;
834 
835  // Inclusive mode; extra jets must be softer than any matched light jet
836  else if (eTpTlightMin >= 0.)
837  for (size_t j = 0; j < jetMomenta.size(); j++) {
838  // CellJet uses eT, SlowJet uses pT
839  if ( (jetAlgorithm == 1 && jetMomenta[j].eT() > eTpTlightMin) ||
840  (jetAlgorithm == 2 && jetMomenta[j].pT() > eTpTlightMin) )
841  return HARD_JET;
842  }
843 
844  } // if (!jetMomenta.empty())
845 
846  // No extra jets were present so no veto
847  return NONE;
848 }
849 
850 //==========================================================================
851 
852 // Main implementation of Madgraph UserHooks class.
853 // This may be split out to a separate C++ file if desired,
854 // but currently included here for ease of use.
855 
856 //--------------------------------------------------------------------------
857 
858 // Initialisation routine automatically called from Pythia::init().
859 // Setup all parts needed for the merging.
860 
861 bool JetMatchingMadgraph::initAfterBeams() {
862 
863  // Read in Madgraph specific configuration variables
864  bool setMad = settingsPtr->flag("JetMatching:setMad");
865 
866  // If Madgraph parameters are present, then parse in MadgraphPar object
867  MadgraphPar par(infoPtr);
868  string parStr = infoPtr->header("MGRunCard");
869  if (!parStr.empty()) {
870  par.parse(parStr);
871  par.printParams();
872  }
873 
874  // Set Madgraph merging parameters from the file if requested
875  if (setMad) {
876  if ( par.haveParam("xqcut") && par.haveParam("maxjetflavor")
877  && par.haveParam("alpsfact") && par.haveParam("ickkw") ) {
878  settingsPtr->flag("JetMatching:merge", par.getParam("ickkw"));
879  settingsPtr->parm("JetMatching:qCut", par.getParam("xqcut"));
880  settingsPtr->mode("JetMatching:nQmatch",
881  par.getParamAsInt("maxjetflavor"));
882  settingsPtr->parm("JetMatching:clFact",
883  clFact = par.getParam("alpsfact"));
884  if (par.getParamAsInt("ickkw") == 0)
885  infoPtr->errorMsg("Error in JetMatchingMadgraph:init: "
886  "Madgraph file parameters are not set for merging");
887 
888  // Warn if setMad requested, but one or more parameters not present
889  } else {
890  infoPtr->errorMsg("Warning in JetMatchingMadgraph:init: "
891  "Madgraph merging parameters not found");
892  if (!par.haveParam("xqcut")) infoPtr->errorMsg("Warning in "
893  "JetMatchingMadgraph:init: No xqcut");
894  if (!par.haveParam("ickkw")) infoPtr->errorMsg("Warning in "
895  "JetMatchingMadgraph:init: No ickkw");
896  if (!par.haveParam("maxjetflavor")) infoPtr->errorMsg("Warning in "
897  "JetMatchingMadgraph:init: No maxjetflavor");
898  if (!par.haveParam("alpsfact")) infoPtr->errorMsg("Warning in "
899  "JetMatchingMadgraph:init: No alpsfact");
900  }
901  }
902 
903  // Read in FxFx matching parameters
904  doFxFx = settingsPtr->flag("JetMatching:doFxFx");
905  nPartonsNow = settingsPtr->mode("JetMatching:nPartonsNow");
906  qCutME = settingsPtr->parm("JetMatching:qCutME");
907  qCutMESq = pow(qCutME,2);
908 
909  // Read in Madgraph merging parameters
910  doMerge = settingsPtr->flag("JetMatching:merge");
911  doShowerKt = settingsPtr->flag("JetMatching:doShowerKt");
912  qCut = settingsPtr->parm("JetMatching:qCut");
913  nQmatch = settingsPtr->mode("JetMatching:nQmatch");
914  clFact = settingsPtr->parm("JetMatching:clFact");
915 
916  // Read in jet algorithm parameters
917  jetAlgorithm = settingsPtr->mode("JetMatching:jetAlgorithm");
918  nJetMax = settingsPtr->mode("JetMatching:nJetMax");
919  eTjetMin = settingsPtr->parm("JetMatching:eTjetMin");
920  coneRadius = settingsPtr->parm("JetMatching:coneRadius");
921  etaJetMax = settingsPtr->parm("JetMatching:etaJetMax");
922  slowJetPower = settingsPtr->mode("JetMatching:slowJetPower");
923 
924  // Matching procedure
925  jetAllow = settingsPtr->mode("JetMatching:jetAllow");
926  exclusiveMode = settingsPtr->mode("JetMatching:exclusive");
927  qCutSq = pow(qCut,2);
928  etaJetMaxAlgo = etaJetMax;
929 
930  // If not merging, then done
931  if (!doMerge) return true;
932 
933  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
934  if (exclusiveMode == 2) {
935 
936  // No nJet or nJetMax, so default to exclusive mode
937  if (nJetMax < 0) {
938  infoPtr->errorMsg("Warning in JetMatchingMadgraph:init: "
939  "missing jet multiplicity information; running in exclusive mode");
940  exclusiveMode = 1;
941  }
942  }
943 
944  // Initialise chosen jet algorithm.
945  // Currently, this only supports the kT-algorithm in SlowJet.
946  // Use the QCD distance measure by default.
947  jetAlgorithm = 2;
948  slowJetPower = 1;
949  slowJet = new SlowJet(slowJetPower, coneRadius, eTjetMin,
950  etaJetMaxAlgo, 2, 2, NULL, false);
951 
952  // For FxFx, also initialise jet algorithm to define matrix element jets.
953  // Currently, this only supports the kT-algorithm in SlowJet.
954  // Use the QCD distance measure by default.
955  slowJetHard = new SlowJet(slowJetPower, coneRadius, qCutME,
956  etaJetMaxAlgo, 2, 2, NULL, false);
957 
958  // Setup local event records
959  eventProcessOrig.init("(eventProcessOrig)", particleDataPtr);
960  eventProcess.init("(eventProcess)", particleDataPtr);
961  workEventJet.init("(workEventJet)", particleDataPtr);
962 
963  // Print information
964  string jetStr = (jetAlgorithm == 1) ? "CellJet" :
965  (slowJetPower == -1) ? "anti-kT" :
966  (slowJetPower == 0) ? "C/A" :
967  (slowJetPower == 1) ? "kT" : "unknown";
968  string modeStr = (exclusiveMode) ? "exclusive" : "inclusive";
969  cout << endl
970  << " *----- Madgraph matching parameters -----*" << endl
971  << " | qCut | " << setw(14)
972  << qCut << " |" << endl
973  << " | nQmatch | " << setw(14)
974  << nQmatch << " |" << endl
975  << " | clFact | " << setw(14)
976  << clFact << " |" << endl
977  << " | Jet algorithm | " << setw(14)
978  << jetStr << " |" << endl
979  << " | eTjetMin | " << setw(14)
980  << eTjetMin << " |" << endl
981  << " | etaJetMax | " << setw(14)
982  << etaJetMax << " |" << endl
983  << " | jetAllow | " << setw(14)
984  << jetAllow << " |" << endl
985  << " | Mode | " << setw(14)
986  << modeStr << " |" << endl
987  << " *-----------------------------------------*" << endl;
988 
989  return true;
990 }
991 
992 //--------------------------------------------------------------------------
993 
994 bool JetMatchingMadgraph::doVetoStep(int iPos, int nISR, int nFSR,
995  const Event& event) {
996 
997  // Do not perform any veto if not in the Shower-kT scheme.
998  if ( !doShowerKt ) return false;
999 
1000  // Default to no veto in case the hard input matrix element already has too
1001  // many partons (same as in Pythia6).
1002  if ( int(typeIdx[0].size()) > nJetMax ) return false;
1003 
1004  // Do nothing for emissions after the first one.
1005  if ( nISR + nFSR > 1 ) return false;
1006 
1007  // Do nothing in resonance decay showers.
1008  if (iPos == 5) return false;
1009 
1010  // Clear the event of MPI systems and resonace decay products. Store trimmed
1011  // event in workEvent.
1012  sortIncomingProcess(event);
1013 
1014  // Get (kinematical) pT of first emission
1015  double pTfirst = 0.;
1016  for (int i = 0; i < workEvent.size(); i++){
1017  // Since this event only contains one emission, this test is enough
1018  // to isolate this emission.
1019  if ( workEvent[i].isFinal()
1020  && (workEvent[i].statusAbs()==43 || workEvent[i].statusAbs()==51)) {
1021  // Only check partons originating from QCD splittings.
1022  int jPos = 1;
1023  bool QCDemission = true;
1024  while ( workEvent[jPos].statusAbs() > 23 ) {
1025  if ( workEvent[jPos].id() == 22 || workEvent[jPos].id() == 23
1026  || workEvent[jPos].idAbs() == 24){
1027  QCDemission = false;
1028  break;
1029  }
1030  jPos = workEvent[jPos].mother1();
1031  }
1032  // Get kinematical pT.
1033  if (QCDemission) {
1034  pTfirst = workEvent[i].pT();
1035  break;
1036  }
1037  }
1038  }
1039 
1040  // Check veto.
1041  if ( doShowerKtVeto(pTfirst) ) return true;
1042 
1043  // No veto if come this far.
1044  return false;
1045 
1046 }
1047 
1048 //--------------------------------------------------------------------------
1049 
1050 bool JetMatchingMadgraph::doShowerKtVeto(double pTfirst) {
1051 
1052  // Only check veto in the shower-kT scheme.
1053  if ( !doShowerKt ) return false;
1054 
1055  // Reset veto code
1056  bool doVeto = false;
1057 
1058  // Find the (kinematical) pT of the softest (light) parton in the hard
1059  // process.
1060  int nParton = typeIdx[0].size();
1061  double pTminME=1e10;
1062  for ( int i = 0; i < nParton; ++i)
1063  pTminME = min(pTminME,eventProcess[typeIdx[0][i]].pT());
1064 
1065  // Veto if the softest hard process parton is below Qcut.
1066  if ( nParton > 0 && pow(pTminME,2) < qCutSq ) doVeto = true;
1067 
1068  // For non-highest multiplicity, veto if the hardest emission is harder
1069  // than Qcut.
1070  if ( exclusive && pow(pTfirst,2) > qCutSq ) {
1071  doVeto = true;
1072  // For highest multiplicity sample, veto if the hardest emission is harder
1073  // than the hard process parton.
1074  } else if ( !exclusive && nParton > 0 && pTfirst > pTminME ) {
1075  doVeto = true;
1076  }
1077 
1078  // Return veto
1079  return doVeto;
1080 
1081 }
1082 
1083 //--------------------------------------------------------------------------
1084 
1085 // Step (1): sort the incoming particles
1086 
1087 void JetMatchingMadgraph::sortIncomingProcess(const Event &event) {
1088 
1089  // Remove resonance decays from original process and keep only final
1090  // state. Resonances will have positive status code after this step.
1091  omitResonanceDecays(eventProcessOrig, true);
1092 
1093  // For FxFx, pre-cluster partons in the event into jets.
1094  if (doFxFx) {
1095 
1096  // Get final state partons
1097  eventProcess.clear();
1098  workEventJet.clear();
1099  for( int i=0; i < workEvent.size(); ++i) {
1100  // Original AG+Py6 algorithm explicitly excludes tops,
1101  // leptons and photons.
1102  int id = workEvent[i].idAbs();
1103  if ((id >= ID_LEPMIN && id <= ID_LEPMAX) || id == ID_TOP
1104  || id == ID_PHOTON || id == 23 || id == 24 || id == 25) {
1105  eventProcess.append(workEvent[i]);
1106  } else {
1107  workEventJet.append(workEvent[i]);
1108  }
1109  }
1110 
1111  // Initialize SlowJetHard jet algorithm with current working event
1112  if (!slowJetHard->setup(workEventJet) ) {
1113  infoPtr->errorMsg("Warning in JetMatchingMadgraph:sortIncomingProcess"
1114  ": the SlowJet algorithm failed on setup");
1115  return;
1116  }
1117 
1118  // Get matrix element cut scale.
1119  double localQcutSq = qCutMESq;
1120  // Cluster in steps to find all hadronic jets at the scale qCutME
1121  while ( slowJetHard->sizeAll() - slowJetHard->sizeJet() > 0 ) {
1122  // Done if next step is above qCut
1123  if( slowJetHard->dNext() > localQcutSq ) break;
1124  // Done if we're at or below the number of partons in the Born state.
1125  if( slowJetHard->sizeAll()-slowJetHard->sizeJet() <= nPartonsNow) break;
1126  slowJetHard->doStep();
1127  }
1128 
1129  // Construct a master copy of the event containing only the
1130  // hardest nPartonsNow hadronic clusters. While constructing the event,
1131  // the parton type (ID_GLUON) and status (98,99) are arbitrary.
1132  int nJets = slowJetHard->sizeJet();
1133  int nClus = slowJetHard->sizeAll();
1134  int nNow = 0;
1135  for (int i = nJets; i < nClus; ++i) {
1136  vector<int> parts;
1137  if (i < nClus-nJets) parts = slowJetHard->clusConstituents(i);
1138  else parts = slowJetHard->constituents(nClus-nJets-i);
1139  int flavour = ID_GLUON;
1140  for(int j=0; j < int(parts.size()); ++j)
1141  if (workEventJet[parts[j]].id() == ID_BOT)
1142  flavour = ID_BOT;
1143  eventProcess.append( flavour, 98,
1144  workEventJet[parts.back()].mother1(),
1145  workEventJet[parts.back()].mother2(),
1146  workEventJet[parts.back()].daughter1(),
1147  workEventJet[parts.back()].daughter2(),
1148  0, 0, slowJetHard->p(i).px(), slowJetHard->p(i).py(),
1149  slowJetHard->p(i).pz(), slowJetHard->p(i).e() );
1150  nNow++;
1151  }
1152 
1153  // Done. Clean-up
1154  workEventJet.clear();
1155 
1156  // For MLM matching, simply take hard process state from workEvent,
1157  // without any preclustering.
1158  } else {
1159  eventProcess = workEvent;
1160  }
1161 
1162  // Sort original process final state into light/heavy jets and 'other'.
1163  // Criteria:
1164  // 1 <= ID <= nQmatch, or ID == 21 --> light jet (typeIdx[0])
1165  // nQMatch < ID --> heavy jet (typeIdx[1])
1166  // All else --> other (typeIdx[2])
1167  // Note that 'typeIdx' stores indices into 'eventProcess' (after resonance
1168  // decays are omitted), while 'typeSet' stores indices into the original
1169  // process record, 'eventProcessOrig', but these indices are also valid
1170  // in 'event'.
1171  for (int i = 0; i < 3; i++) {
1172  typeIdx[i].clear();
1173  typeSet[i].clear();
1174  origTypeIdx[i].clear();
1175  }
1176  for (int i = 0; i < eventProcess.size(); i++) {
1177  // Ignore non-final state and default to 'other'
1178  if (!eventProcess[i].isFinal()) continue;
1179  int idx = 2;
1180  int orig_idx = 2;
1181 
1182  // Light jets: all gluons and quarks with id less than or equal to nQmatch
1183  if (eventProcess[i].id() == ID_GLUON
1184  || (eventProcess[i].idAbs() <= nQmatch) ) {
1185  orig_idx = 0;
1186  // Crucial point: MG puts the scale of a non-QCD particle to eCM. For
1187  // such particles, we should keep the default "2"
1188  if ( eventProcess[i].scale() < 1.999*sqrt(infoPtr->eA()*infoPtr->eB()) )
1189  idx = 0;
1190  }
1191 
1192  // Heavy jets: all quarks with id greater than nQmatch
1193  else if (eventProcess[i].idAbs() > nQmatch
1194  && eventProcess[i].idAbs() <= ID_TOP) {
1195  idx = 1;
1196  orig_idx = 1;
1197 
1198  } else {
1199  idx = 2;
1200  orig_idx = 2;
1201  }
1202 
1203  // Store
1204  typeIdx[idx].push_back(i);
1205  typeSet[idx].insert(eventProcess[i].daughter1());
1206  origTypeIdx[orig_idx].push_back(i);
1207  }
1208 
1209  // Exclusive mode; if set to 2, then set based on nJet/nJetMax
1210  if (exclusiveMode == 2) {
1211 
1212  // Inclusive if nJet == nJetMax, exclusive otherwise
1213  int nParton = origTypeIdx[0].size();
1214  exclusive = (nParton == nJetMax) ? false : true;
1215 
1216  // Otherwise, just set as given
1217  } else {
1218  exclusive = (exclusiveMode == 0) ? false : true;
1219  }
1220 
1221  // Extract partons from hardest subsystem + ISR + FSR only into
1222  // workEvent. Note no resonance showers or MPIs.
1223  subEvent(event);
1224 }
1225 
1226 //--------------------------------------------------------------------------
1227 
1228 // Step (2a): pick which particles to pass to the jet algorithm
1229 
1230 void JetMatchingMadgraph::jetAlgorithmInput(const Event &event, int iType) {
1231 
1232  // Take input from 'workEvent' and put output in 'workEventJet'
1233  workEventJet = workEvent;
1234 
1235  // Loop over particles and decide what to pass to the jet algorithm
1236  for (int i = 0; i < workEventJet.size(); ++i) {
1237  if (!workEventJet[i].isFinal()) continue;
1238 
1239  // jetAllow option to disallow certain particle types
1240  if (jetAllow == 1) {
1241 
1242  // Original AG+Py6 algorithm explicitly excludes tops,
1243  // leptons and photons.
1244  int id = workEventJet[i].idAbs();
1245  if ((id >= ID_LEPMIN && id <= ID_LEPMAX) || id == ID_TOP
1246  || id == ID_PHOTON || (id > nQmatch && id!=21)) {
1247  workEventJet[i].statusNeg();
1248  continue;
1249  }
1250  }
1251 
1252  // Get the index of this particle in original event
1253  int idx = workEventJet[i].daughter1();
1254 
1255  // Start with particle idx, and afterwards track mothers
1256  while (true) {
1257 
1258  // Light jets
1259  if (iType == 0) {
1260 
1261  // Do not include if originates from heavy jet or 'other'
1262  if (typeSet[1].find(idx) != typeSet[1].end() ||
1263  typeSet[2].find(idx) != typeSet[2].end()) {
1264  workEventJet[i].statusNeg();
1265  break;
1266  }
1267 
1268  // Made it to start of event record so done
1269  if (idx == 0) break;
1270  // Otherwise next mother and continue
1271  idx = event[idx].mother1();
1272 
1273  // Heavy jets
1274  } else if (iType == 1) {
1275 
1276  // Only include if originates from heavy jet
1277  if (typeSet[1].find(idx) != typeSet[1].end()) break;
1278 
1279  // Made it to start of event record with no heavy jet mother,
1280  // so DO NOT include particle
1281  if (idx == 0) {
1282  workEventJet[i].statusNeg();
1283  break;
1284  }
1285 
1286  // Otherwise next mother and continue
1287  idx = event[idx].mother1();
1288 
1289  } // if (iType)
1290  } // while (true)
1291  } // for (i)
1292 }
1293 
1294 //--------------------------------------------------------------------------
1295 
1296 // Step (2b): run jet algorithm and provide common output
1297 // This does nothing, because the jet algorithm is run several times
1298 // in the matching algorithm.
1299 
1300 void JetMatchingMadgraph::runJetAlgorithm() {; }
1301 
1302 //--------------------------------------------------------------------------
1303 
1304 // Step (2c): veto decision (returning true vetoes the event)
1305 
1306 bool JetMatchingMadgraph::matchPartonsToJets(int iType) {
1307 
1308  // Use two different routines for light/heavy jets as
1309  // different veto conditions and for clarity
1310  if (iType == 0) return (matchPartonsToJetsLight() > 0);
1311  else return (matchPartonsToJetsHeavy() > 0);
1312 }
1313 
1314 //--------------------------------------------------------------------------
1315 
1316 // Step(2c): light jets
1317 // Return codes are given indicating the reason for a veto.
1318 // Although not currently used, they are a useful debugging tool:
1319 // 0 = no veto
1320 // 1 = veto as number of jets less than number of partons
1321 // 2 = veto as exclusive mode and number of jets greater than
1322 // number of partons
1323 // 3 = veto as inclusive mode and there would be an extra jet
1324 // that is harder than any matched soft jet
1325 // 4 = veto as there is a parton which does not match a jet
1326 
1327 int JetMatchingMadgraph::matchPartonsToJetsLight() {
1328 
1329  // Count the number of hard partons
1330  //BUG!! int nParton = origTypeIdx[0].size();
1331  int nParton = typeIdx[0].size();
1332 
1333  // Initialize SlowJet with current working event
1334  if (!slowJet->setup(workEventJet) ) {
1335  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1336  "Light: the SlowJet algorithm failed on setup");
1337  return NONE;
1338  }
1339  double localQcutSq = qCutSq;
1340  double dOld = 0.0;
1341  // Cluster in steps to find all hadronic jets at the scale qCut
1342  while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1343  if( slowJet->dNext() > localQcutSq ) break;
1344  dOld = slowJet->dNext();
1345  slowJet->doStep();
1346  }
1347  int nJets = slowJet->sizeJet();
1348  int nClus = slowJet->sizeAll();
1349 
1350  // Debug printout.
1351  if (MATCHINGDEBUG) slowJet->list(true);
1352 
1353  // Count of the number of hadronic jets in SlowJet accounting
1354  int nCLjets = nClus - nJets;
1355  // Get number of partons. Different for MLM and FxFx schemes.
1356  int nRequested = (doFxFx) ? nPartonsNow : nParton;
1357 
1358  // Veto event if too few hadronic jets
1359  if ( nCLjets < nRequested ) return LESS_JETS;
1360 
1361  // In exclusive mode, do not allow more hadronic jets than partons
1362  if ( exclusive && !doFxFx ) {
1363  if ( nCLjets > nRequested ) return MORE_JETS;
1364  } else {
1365 
1366  // For FxFx, in the non-highest multipicity, all jets need to matched to
1367  // partons. For nCLjets > nRequested, this is not possible. Hence, we can
1368  // veto here already.
1369  if ( doFxFx && nRequested < nJetMax && nCLjets > nRequested )
1370  return MORE_JETS;
1371 
1372  // Now continue in inclusive mode.
1373  // In inclusive mode, there can be more hadronic jets than partons,
1374  // provided that all partons are properly matched to hadronic jets.
1375  // Start by setting up the jet algorithm.
1376  if (!slowJet->setup(workEventJet) ) {
1377  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1378  "Light: the SlowJet algorithm failed on setup");
1379  return NONE;
1380  }
1381 
1382  // For FxFx, continue clustering as long as the jet separation is above
1383  // qCut.
1384  if (doFxFx) {
1385  while ( slowJet->sizeAll() - slowJet->sizeJet() > 0 ) {
1386  if( slowJet->dNext() > localQcutSq ) break;
1387  slowJet->doStep();
1388  }
1389  // For MLM, cluster into hadronic jets until there are the same number as
1390  // partons.
1391  } else {
1392  while ( slowJet->sizeAll() - slowJet->sizeJet() > nParton )
1393  slowJet->doStep();
1394  }
1395 
1396  // Sort partons in pT. Update local qCut value.
1397  // Hadronic jets are already sorted in pT.
1398  localQcutSq = dOld;
1399  if ( clFact >= 0. && nParton > 0 ) {
1400  vector<double> partonPt;
1401  for (int i = 0; i < nParton; ++i)
1402  partonPt.push_back( eventProcess[typeIdx[0][i]].pT2() );
1403  sort( partonPt.begin(), partonPt.end());
1404  localQcutSq = max( qCutSq, partonPt[0]);
1405  }
1406  nJets = slowJet->sizeJet();
1407  nClus = slowJet->sizeAll();
1408  }
1409  // Update scale if clustering factor is non-zero
1410  if ( clFact != 0. ) localQcutSq *= pow2(clFact);
1411 
1412  Event tempEvent;
1413  tempEvent.init( "(tempEvent)", particleDataPtr);
1414  int nPass = 0;
1415  double pTminEstimate = -1.;
1416  // Construct a master copy of the event containing only the
1417  // hardest nParton hadronic clusters. While constructing the event,
1418  // the parton type (ID_GLUON) and status (98,99) are arbitrary.
1419  for (int i = nJets; i < nClus; ++i) {
1420  tempEvent.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0, slowJet->p(i).px(),
1421  slowJet->p(i).py(), slowJet->p(i).pz(), slowJet->p(i).e() );
1422  ++nPass;
1423  pTminEstimate = max( pTminEstimate, slowJet->pT(i));
1424  if(nPass == nRequested) break;
1425  }
1426 
1427  int tempSize = tempEvent.size();
1428  // This keeps track of which hadronic jets are matched to parton
1429  vector<bool> jetAssigned;
1430  jetAssigned.assign( tempSize, false);
1431 
1432  // This keeps track of which partons are matched to which hadronic
1433  // jets.
1434  vector< vector<bool> > partonMatchesJet;
1435  for (int i=0; i < nParton; ++i )
1436  partonMatchesJet.push_back( vector<bool>(tempEvent.size(),false) );
1437 
1438  // Begin matching.
1439  // Do jet matching for FxFx.
1440  // Make sure that the nPartonsNow hardest hadronic jets are matched to any
1441  // of the nPartonsNow (+1) partons. This matching is done by attaching a jet
1442  // from the list of unmatched hadronic jets, and appending a jet from the
1443  // list of partonic jets, one at a time. The partonic jet will be clustered
1444  // with the hadronic jet or the beam if the distance measure is below the
1445  // cut. The hadronic jet is matched once this happens. Otherwise, another
1446  // partonic jet is tried. When a hadronic jet is matched to a partonic jet,
1447  // it is removed from the list of unmatched hadronic jets. This process
1448  // continues until the nPartonsNow hardest hadronic jets are matched to
1449  // partonic jets, or it is not possible to make a match for a hadronic jet.
1450  int iNow = 0;
1451  while ( doFxFx && iNow < tempSize ) {
1452 
1453  // Check if this shower jet matches any partonic jet.
1454  Event tempEventJet;
1455  tempEventJet.init("(tempEventJet)", particleDataPtr);
1456  for (int i=0; i < nParton; ++i ) {
1457 
1459  //for (int j=0; j < tempSize; ++j )
1460  // if ( partonMatchesJet[i][j]) continue;
1461 
1462  // Attach a single hadronic jet.
1463  tempEventJet.clear();
1464  tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1465  tempEvent[iNow].px(), tempEvent[iNow].py(),
1466  tempEvent[iNow].pz(), tempEvent[iNow].e() );
1467  // Attach the current parton.
1468  Vec4 pIn = eventProcess[typeIdx[0][i]].p();
1469  tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1470  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1471 
1472  // Setup jet algorithm.
1473  if ( !slowJet->setup(tempEventJet) ) {
1474  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1475  "Light: the SlowJet algorithm failed on setup");
1476  return NONE;
1477  }
1478 
1479  // These are the conditions for the hadronic jet to match the parton
1480  // at the local qCut scale
1481  if ( slowJet->iNext() == tempEventJet.size() - 1
1482  && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1483  jetAssigned[iNow] = true;
1484  partonMatchesJet[i][iNow] = true;
1485  }
1486 
1487  } // End loop over hard partons.
1488 
1489  // Veto if the jet could not be assigned to any parton.
1490  if ( !jetAssigned[iNow] ) return UNMATCHED_PARTON;
1491 
1492  // Continue;
1493  ++iNow;
1494  }
1495 
1496  // Do jet matching for MLM.
1497  // Take the list of unmatched hadronic jets and append a parton, one at
1498  // a time. The parton will be clustered with the "closest" hadronic jet
1499  // or the beam if the distance measure is below the cut. When a hadronic
1500  // jet is matched to a parton, it is removed from the list of unmatched
1501  // hadronic jets. This process continues until all hadronic jets are
1502  // matched to partons or it is not possible to make a match.
1503  iNow = 0;
1504  while (!doFxFx && iNow < nParton ) {
1505  Event tempEventJet;
1506  tempEventJet.init("(tempEventJet)", particleDataPtr);
1507  for (int i = 0; i < tempSize; ++i) {
1508  if (jetAssigned[i]) continue;
1509  Vec4 pIn = tempEvent[i].p();
1510  // Append unmatched hadronic jets
1511  tempEventJet.append( ID_GLUON, 98, 0, 0, 0, 0, 0, 0,
1512  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1513  }
1514 
1515  Vec4 pIn = eventProcess[typeIdx[0][iNow]].p();
1516  // Append the current parton
1517  tempEventJet.append( ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
1518  pIn.px(), pIn.py(), pIn.pz(), pIn.e() );
1519  if ( !slowJet->setup(tempEventJet) ) {
1520  infoPtr->errorMsg("Warning in JetMatchingMadgraph:matchPartonsToJets"
1521  "Light: the SlowJet algorithm failed on setup");
1522  return NONE;
1523  }
1524  // These are the conditions for the hadronic jet to match the parton
1525  // at the local qCut scale
1526  if ( slowJet->iNext() == tempEventJet.size() - 1
1527  && slowJet->jNext() > -1 && slowJet->dNext() < localQcutSq ) {
1528  int iKnt = -1;
1529  for (int i = 0; i != tempSize; ++i) {
1530  if (jetAssigned[i]) continue;
1531  ++iKnt;
1532  // Identify the hadronic jet that matches the parton
1533  if (iKnt == slowJet->jNext() ) jetAssigned[i] = true;
1534  }
1535  } else {
1536  return UNMATCHED_PARTON;
1537  }
1538  ++iNow;
1539  }
1540 
1541  // Minimal eT/pT (CellJet/SlowJet) of matched light jets.
1542  // Needed later for heavy jet vetos in inclusive mode.
1543  // This information is not used currently.
1544  if (nParton > 0 && pTminEstimate > 0) eTpTlightMin = pTminEstimate;
1545  else eTpTlightMin = -1.;
1546 
1547  // No veto
1548  return NONE;
1549 }
1550 
1551 //--------------------------------------------------------------------------
1552 
1553 // Step(2c): heavy jets
1554 // Return codes are given indicating the reason for a veto.
1555 // Although not currently used, they are a useful debugging tool:
1556 // 0 = no veto as there are no extra jets present
1557 // 1 = veto as in exclusive mode and extra jets present
1558 // 2 = veto as in inclusive mode and extra jets were harder
1559 // than any matched light jet
1560 
1561 int JetMatchingMadgraph::matchPartonsToJetsHeavy() {
1562 
1563  // Currently, heavy jets are unmatched
1564  // If there are no extra jets, then accept
1565  if (jetMomenta.empty()) return NONE;
1566 
1567  // No extra jets were present so no veto
1568  return NONE;
1569 }
1570 
1571 //==========================================================================
1572 
1573 #endif // end Pythia8_JetMatching_H