1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
|
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
#include "DataFormats/JetReco/interface/GenJetCollection.h"
#include "DataFormats/Math/interface/deltaPhi.h"
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/Math/interface/LorentzVector.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDFilter.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/EDGetToken.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include <cmath>
#include <cstdlib>
#include <vector>
class VBFGenJetFilter : public edm::global::EDFilter<> {
public:
explicit VBFGenJetFilter(const edm::ParameterSet&);
bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
private:
std::vector<const reco::GenJet*> filterGenJets(const std::vector<reco::GenJet>* jets) const;
std::vector<const reco::GenParticle*> filterGenLeptons(const std::vector<reco::GenParticle>* particles) const;
// Dijet cut
const bool oppositeHemisphere;
const bool leadJetsNoLepMass;
const double ptMin;
const double etaMin;
const double etaMax;
const double minInvMass;
const double maxInvMass;
const double minDeltaPhi;
const double maxDeltaPhi;
const double minDeltaEta;
const double maxDeltaEta;
const double minLeadingJetsInvMass;
const double maxLeadingJetsInvMass;
const double deltaRJetLep;
// Input tags
edm::EDGetTokenT<reco::GenJetCollection> m_inputTag_GenJetCollection;
edm::EDGetTokenT<reco::GenParticleCollection> m_inputTag_GenParticleCollection;
};
using namespace std;
VBFGenJetFilter::VBFGenJetFilter(const edm::ParameterSet& iConfig)
: oppositeHemisphere(iConfig.getUntrackedParameter<bool>("oppositeHemisphere", false)),
leadJetsNoLepMass(iConfig.getUntrackedParameter<bool>("leadJetsNoLepMass", false)),
ptMin(iConfig.getUntrackedParameter<double>("minPt", 20)),
etaMin(iConfig.getUntrackedParameter<double>("minEta", -5.0)),
etaMax(iConfig.getUntrackedParameter<double>("maxEta", 5.0)),
minInvMass(iConfig.getUntrackedParameter<double>("minInvMass", 0.0)),
maxInvMass(iConfig.getUntrackedParameter<double>("maxInvMass", 99999.0)),
minDeltaPhi(iConfig.getUntrackedParameter<double>("minDeltaPhi", -1.0)),
maxDeltaPhi(iConfig.getUntrackedParameter<double>("maxDeltaPhi", 99999.0)),
minDeltaEta(iConfig.getUntrackedParameter<double>("minDeltaEta", -1.0)),
maxDeltaEta(iConfig.getUntrackedParameter<double>("maxDeltaEta", 99999.0)),
minLeadingJetsInvMass(iConfig.getUntrackedParameter<double>("minLeadingJetsInvMass", 0.0)),
maxLeadingJetsInvMass(iConfig.getUntrackedParameter<double>("maxLeadingJetsInvMass", 99999.0)),
deltaRJetLep(iConfig.getUntrackedParameter<double>("deltaRJetLep", 0.3)) {
m_inputTag_GenJetCollection = consumes<reco::GenJetCollection>(
iConfig.getUntrackedParameter<edm::InputTag>("inputTag_GenJetCollection", edm::InputTag("ak5GenJetsNoNu")));
if (leadJetsNoLepMass)
m_inputTag_GenParticleCollection = consumes<reco::GenParticleCollection>(
iConfig.getUntrackedParameter<edm::InputTag>("genParticles", edm::InputTag("genParticles")));
}
vector<const reco::GenParticle*> VBFGenJetFilter::filterGenLeptons(const vector<reco::GenParticle>* particles) const {
vector<const reco::GenParticle*> out;
for (const auto& p : *particles) {
int absPdgId = std::abs(p.pdgId());
if (((absPdgId == 11) || (absPdgId == 13) || (absPdgId == 15)) && p.isHardProcess()) {
out.push_back(&p);
}
}
return out;
}
vector<const reco::GenJet*> VBFGenJetFilter::filterGenJets(const vector<reco::GenJet>* jets) const {
vector<const reco::GenJet*> out;
for (unsigned i = 0; i < jets->size(); i++) {
const reco::GenJet* j = &((*jets)[i]);
if (j->p4().pt() > ptMin && j->p4().eta() > etaMin && j->p4().eta() < etaMax) {
out.push_back(j);
}
}
return out;
}
// ------------ method called to skim the data ------------
bool VBFGenJetFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
edm::Handle<vector<reco::GenJet> > handleGenJets;
iEvent.getByToken(m_inputTag_GenJetCollection, handleGenJets);
const vector<reco::GenJet>* genJets = handleGenJets.product();
// Getting filtered generator jets
vector<const reco::GenJet*> filGenJets = filterGenJets(genJets);
// If we do not find at least 2 jets veto the event
if (filGenJets.size() < 2) {
return false;
}
// Testing dijet mass
if (leadJetsNoLepMass) {
edm::Handle<reco::GenParticleCollection> genParticelesCollection;
iEvent.getByToken(m_inputTag_GenParticleCollection, genParticelesCollection);
const vector<reco::GenParticle>* genParticles = genParticelesCollection.product();
// Getting filtered generator muons
vector<const reco::GenParticle*> filGenLep = filterGenLeptons(genParticles);
// Getting p4 of jet with no lepton
vector<math::XYZTLorentzVector> genJetsWithoutLeptonsP4;
unsigned int jetIdx = 0;
while (genJetsWithoutLeptonsP4.size() < 2 && jetIdx < filGenJets.size()) {
bool jetWhitoutLep = true;
const math::XYZTLorentzVector& p4J = (filGenJets[jetIdx])->p4();
for (unsigned int i = 0; i < filGenLep.size() && jetWhitoutLep; ++i) {
if (reco::deltaR2((filGenLep[i])->p4(), p4J) < deltaRJetLep * deltaRJetLep)
jetWhitoutLep = false;
}
if (jetWhitoutLep)
genJetsWithoutLeptonsP4.push_back(p4J);
++jetIdx;
}
// Checking the invariant mass of the leading jets
if (genJetsWithoutLeptonsP4.size() < 2)
return false;
float invMassLeadingJet = (genJetsWithoutLeptonsP4[0] + genJetsWithoutLeptonsP4[1]).M();
if (invMassLeadingJet > minLeadingJetsInvMass && invMassLeadingJet < maxLeadingJetsInvMass)
return true;
else
return false;
}
for (unsigned a = 0; a < filGenJets.size(); a++) {
for (unsigned b = a + 1; b < filGenJets.size(); b++) {
const reco::GenJet* pA = filGenJets[a];
const reco::GenJet* pB = filGenJets[b];
// Getting the dijet vector
math::XYZTLorentzVector diJet = pA->p4() + pB->p4();
// Testing opposite hemispheres
double dijetProd = pA->p4().eta() * pB->p4().eta();
if (oppositeHemisphere && dijetProd >= 0) {
continue;
}
// Testing dijet mass
double invMass = diJet.mass();
if (invMass <= minInvMass || invMass > maxInvMass) {
continue;
}
// Testing dijet delta eta
double dEta = fabs(pA->p4().eta() - pB->p4().eta());
if (dEta <= minDeltaEta || dEta > maxDeltaEta) {
continue;
}
// Testing dijet delta phi
double dPhi = fabs(reco::deltaPhi(pA->p4().phi(), pB->p4().phi()));
if (dPhi <= minDeltaPhi || dPhi > maxDeltaPhi) {
continue;
}
return true;
}
}
return false;
}
//define this as a plug-in
DEFINE_FWK_MODULE(VBFGenJetFilter);
|