CMS 3D CMS Logo

JetMatchingEWKFxFx.h
Go to the documentation of this file.
1 #ifndef GeneratorInterface_Pythia8Interface_JetMatchingEWKFxFx_h
2 #define GeneratorInterface_Pythia8Interface_JetMatchingEWKFxFx_h
3 
4 // Class declaration for the new JetMatching plugin
5 // Author: Carlos Vico (U. Oviedo)
6 // taken from: https://amcatnlo.web.cern.ch/amcatnlo/JetMatching.h
7 
8 // Includes
11 #include "Pythia8/Pythia.h"
12 #include "Pythia8Plugins/JetMatching.h"
13 #include "Pythia8Plugins/GeneratorInput.h"
14 #include <memory>
16 
17 // The plugin must inherit from the original Pythia8::JetMatching
18 // class
19 
20 class JetMatchingEWKFxFx : public Pythia8::JetMatching {
21 public:
22  // Constructor and destructor
23  JetMatchingEWKFxFx(const edm::ParameterSet& iConfig);
24  ~JetMatchingEWKFxFx() override {}
25 
26  // Method declaration
27  bool initAfterBeams() override;
28 
29  bool canVetoPartonLevelEarly() override { return true; }
30  bool doVetoPartonLevelEarly(const Pythia8::Event& event) override;
31 
32  bool canVetoProcessLevel() override { return true; }
33  bool doVetoProcessLevel(Pythia8::Event& event) override;
34 
35  // Shower step vetoes (after the first emission, for Shower-kT scheme)
36  int numberVetoStep() override { return 1; }
37  bool canVetoStep() override { return doShowerKt; }
38  bool doVetoStep(int, int, int, const Pythia8::Event&) override;
39 
40  Pythia8::SlowJet* slowJetDJR;
41 
42  Pythia8::vector<double> getDJR() { return DJR; }
43  Pythia8::pair<int, int> nMEpartons() { return nMEpartonsSave; }
44 
45  Pythia8::Event getWorkEventJet() { return workEventJetSave; }
46  Pythia8::Event getProcessSubset() { return processSubsetSave; }
47  bool getExclusive() { return exclusive; }
48  double getPTfirst() { return pTfirstSave; }
49 
50  // Different steps of the matching algorithm.
51  void sortIncomingProcess(const Pythia8::Event&) override;
52 
53  void jetAlgorithmInput(const Pythia8::Event&, int) override;
54  void runJetAlgorithm() override;
55  bool matchPartonsToJets(int) override;
56  int matchPartonsToJetsLight() override;
57  int matchPartonsToJetsHeavy() override;
59  bool doShowerKtVeto(double pTfirst);
60 
61  // Functions to clear and set the jet clustering scales.
62  void clearDJR() { DJR.resize(0); }
63  void setDJR(const Pythia8::Event& event);
64  // Functions to clear and set the jet clustering scales.
65  void clear_nMEpartons() { nMEpartonsSave.first = nMEpartonsSave.second = -1; }
66  void set_nMEpartons(const int nOrig, const int nMatch) {
68  nMEpartonsSave.first = nOrig;
69  nMEpartonsSave.second = nMatch;
70  };
71 
72  // Function to get the current number of partons in the Born state, as
73  // read from LHE.
74  int npNLO();
75 
76 private:
77  Pythia8::Event processSubsetSave;
78  Pythia8::Event workEventJetSave;
79  double pTfirstSave;
80 
82 
83  Pythia8::vector<int> origTypeIdx[3];
84  int nQmatch;
85  double qCut, qCutSq, clFact;
86  bool doFxFx;
88  double qCutME, qCutMESq;
89 
90  Pythia8::vector<double> DJR;
91 
92  Pythia8::pair<int, int> nMEpartonsSave;
93 };
94 
95 // Register in the UserHook factory
97 #endif
bool matchPartonsToJets(int) override
void set_nMEpartons(const int nOrig, const int nMatch)
void jetAlgorithmInput(const Pythia8::Event &, int) override
bool doVetoPartonLevelEarly(const Pythia8::Event &event) override
Pythia8::Event getProcessSubset()
int matchPartonsToJetsLight() override
Pythia8::vector< double > getDJR()
Pythia8::vector< int > origTypeIdx[3]
Pythia8::Event processSubsetSave
Pythia8::pair< int, int > nMEpartons()
bool doVetoStep(int, int, int, const Pythia8::Event &) override
int matchPartonsToJetsHeavy() override
void runJetAlgorithm() override
bool doShowerKtVeto(double pTfirst)
int numberVetoStep() override
bool canVetoPartonLevelEarly() override
void sortIncomingProcess(const Pythia8::Event &) override
bool initAfterBeams() override
Pythia8::SlowJet * slowJetDJR
void setDJR(const Pythia8::Event &event)
Pythia8::pair< int, int > nMEpartonsSave
JetMatchingEWKFxFx(const edm::ParameterSet &iConfig)
bool canVetoStep() override
Pythia8::Event workEventJetSave
~JetMatchingEWKFxFx() override
#define REGISTER_USERHOOK(type)
Definition: CustomHook.h:8
Pythia8::Event getWorkEventJet()
bool doVetoProcessLevel(Pythia8::Event &event) override
Definition: event.py:1
Pythia8::vector< double > DJR
bool canVetoProcessLevel() override