5 #include "fastjet/ClusterSequence.hh" 15 #include "fastjet/JetDefinition.hh" 41 typedef fastjet::JetDefinition
JetDef;
53 minLeptonPt_(
pset.getParameter<double>(
"minLeptonPt")),
54 maxLeptonEta_(
pset.getParameter<double>(
"maxLeptonEta")),
55 minJetPt_(
pset.getParameter<double>(
"minJetPt")),
56 maxJetEta_(
pset.getParameter<double>(
"maxJetEta")),
57 wMass_(
pset.getParameter<double>(
"wMass")),
58 tMass_(
pset.getParameter<double>(
"tMass")),
59 minLeptonPtDilepton_(
pset.getParameter<double>(
"minLeptonPtDilepton")),
60 maxLeptonEtaDilepton_(
pset.getParameter<double>(
"maxLeptonEtaDilepton")),
61 minDileptonMassDilepton_(
pset.getParameter<double>(
"minDileptonMassDilepton")),
62 minLeptonPtSemilepton_(
pset.getParameter<double>(
"minLeptonPtSemilepton")),
63 maxLeptonEtaSemilepton_(
pset.getParameter<double>(
"maxLeptonEtaSemilepton")),
64 minVetoLeptonPtSemilepton_(
pset.getParameter<double>(
"minVetoLeptonPtSemilepton")),
65 maxVetoLeptonEtaSemilepton_(
pset.getParameter<double>(
"maxVetoLeptonEtaSemilepton")),
66 minMETSemiLepton_(
pset.getParameter<double>(
"minMETSemiLepton")),
67 minMtWSemiLepton_(
pset.getParameter<double>(
"minMtWSemiLepton")) {
75 produces<reco::GenParticleCollection>(
"neutrinos");
76 produces<reco::GenJetCollection>(
"leptons");
77 produces<reco::GenJetCollection>(
"jets");
79 produces<reco::GenParticleCollection>();
100 std::set<int> bHadronIdxs;
101 for (
size_t i = 0,
n = genParticleHandle->size();
i <
n; ++
i) {
109 bHadronIdxs.insert(-
i - 1);
114 std::vector<size_t> leptonIdxs;
115 for (
size_t i = 0,
n = finalStateHandle->size();
i <
n; ++
i) {
117 const int absPdgId =
abs(
p.pdgId());
122 if (
p.numberOfMothers() == 0)
124 if (
p.mother()->status() == 4)
132 leptonIdxs.push_back(
i);
147 std::vector<fastjet::PseudoJet> fjLepInputs;
148 fjLepInputs.reserve(leptonIdxs.size());
149 for (
auto index : leptonIdxs) {
154 fjLepInputs.push_back(fastjet::PseudoJet(
p.px(),
p.py(),
p.pz(),
p.energy()));
155 fjLepInputs.back().set_user_index(
index);
159 fastjet::ClusterSequence fjLepClusterSeq(fjLepInputs, *
fjLepDef_);
160 std::vector<fastjet::PseudoJet> fjLepJets = fastjet::sorted_by_pt(fjLepClusterSeq.inclusive_jets(
minLeptonPt_));
163 leptons->reserve(fjLepJets.size());
164 std::set<size_t> lepDauIdxs;
165 for (
auto& fjJet : fjLepJets) {
170 const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
172 std::vector<reco::CandidatePtr> constituents;
174 for (
auto& fjConstituent : fjConstituents) {
175 const size_t index = fjConstituent.user_index();
177 const int absPdgId =
abs(
cand->pdgId());
178 if (absPdgId == 11
or absPdgId == 13) {
183 constituents.push_back(
cand);
188 const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
192 lepJet.setPdgId(lepCand->
pdgId());
193 lepJet.setCharge(lepCand->
charge());
195 const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
196 lepJet.setJetArea(jetArea);
201 for (
auto& fjConstituent : fjConstituents) {
202 lepDauIdxs.insert(fjConstituent.user_index());
209 std::vector<fastjet::PseudoJet> fjJetInputs;
210 fjJetInputs.reserve(nStables);
211 for (
size_t i = 0,
n = finalStateHandle->size();
i <
n; ++
i) {
219 if (absId == 12
or absId == 14
or absId == 16)
221 if (lepDauIdxs.find(
i) != lepDauIdxs.end())
224 fjJetInputs.push_back(fastjet::PseudoJet(
p.px(),
p.py(),
p.pz(),
p.energy()));
225 fjJetInputs.back().set_user_index(
i);
228 for (
auto index : bHadronIdxs) {
233 const double scale = 1
e-20 /
p.p();
235 fjJetInputs.back().set_user_index(
index);
239 fastjet::ClusterSequence fjJetClusterSeq(fjJetInputs, *
fjJetDef_);
240 std::vector<fastjet::PseudoJet> fjJets = fastjet::sorted_by_pt(fjJetClusterSeq.inclusive_jets(
minJetPt_));
243 jets->reserve(fjJets.size());
244 std::vector<size_t> bjetIdxs, ljetIdxs;
245 for (
auto& fjJet : fjJets) {
250 const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
252 std::vector<reco::CandidatePtr> constituents;
253 bool hasBHadron =
false;
254 for (
size_t j = 0,
m = fjConstituents.size();
j <
m; ++
j) {
255 const int index = fjConstituents[
j].user_index();
256 if (bHadronIdxs.find(
index) != bHadronIdxs.end()) {
261 constituents.push_back(
cand);
264 const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
268 const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
269 genJet.setJetArea(jetArea);
272 bjetIdxs.push_back(
jets->size());
274 ljetIdxs.push_back(
jets->size());
277 jets->push_back(genJet);
283 if (bjetIdxs.size() < 2)
289 const int q1 =
leptons->at(0).charge();
290 const int q2 =
leptons->at(1).charge();
305 int selNu1 = -1, selNu2 = -1;
306 for (
int i = 0;
i < 2; ++
i) {
308 for (
int j = 0;
j < 2; ++
j) {
312 const double newDm = dm1 + dm2;
326 const auto w1LVec = lepton1.p4() + nu1.p4();
327 const auto w2LVec = lepton2.p4() + nu2.p4();
331 int selB1 = -1, selB2 = -1;
332 for (
auto i : bjetIdxs) {
334 for (
auto j : bjetIdxs) {
338 const double newDm = dm1 + dm2;
352 const auto t1LVec = w1LVec +
bJet1.p4();
353 const auto t2LVec = w2LVec +
bJet2.p4();
356 const int lepQ1 = lepton1.charge();
357 const int lepQ2 = lepton2.charge();
384 }
else if (!
leptons->empty() and !
neutrinos->empty() and ljetIdxs.size() >= 2) {
386 const auto& lepton =
leptons->at(0);
391 bool hasVetoLepton =
false;
392 for (
auto vetoLeptonItr =
next(
leptons->begin()); vetoLeptonItr !=
leptons->end(); ++vetoLeptonItr) {
395 hasVetoLepton =
true;
402 double metX = 0, metY = 0;
403 for (
const auto& neutrino : *
neutrinos) {
404 metX += neutrino.px();
405 metY += neutrino.py();
411 const double mtW =
std::sqrt(2 * lepton.pt() * metPt *
cos(lepton.phi() - atan2(metX, metY)));
421 const double lpz = lepton.pz(), le = lepton.energy(), lpt = lepton.pt();
422 const double cf = (
wMass_ *
wMass_) / 2 + lepton.px() * metX + lepton.py() * metY;
423 const double det = cf * cf * lpz * lpz - lpt * lpt * (le * le * metPt * metPt - cf * cf);
424 const double pz = (cf * lpz + (cf < 0 ? -
sqrt(det) :
sqrt(det))) / lpt / lpt;
426 const auto w1LVec = lepton.p4() + nu1P4;
431 for (
auto b1Itr = bjetIdxs.begin(); b1Itr != bjetIdxs.end(); ++b1Itr) {
432 const double dR =
deltaR(
jets->at(*b1Itr).p4(), w1LVec);
441 const auto t1LVec = w1LVec +
bJet1.p4();
444 const auto& wJet1 =
jets->at(ljetIdxs[0]);
445 const auto& wJet2 =
jets->at(ljetIdxs[1]);
446 const auto w2LVec = wJet1.p4() + wJet2.p4();
451 for (
auto i : bjetIdxs) {
464 const auto t2LVec = w2LVec +
bJet2.p4();
466 const int q = lepton.charge();
530 for (
size_t i = 0,
n =
p->numberOfMothers();
i <
n; ++
i) {
547 const unsigned int absPdgId =
abs(
p->pdgId());
553 for (
int i = 0,
n =
p->numberOfDaughters();
i <
n; ++
i) {
565 if (absPdgId >= 1000000000)
571 const int nq3 = (absPdgId / 10) % 10;
572 const int nq2 = (absPdgId / 100) % 10;
573 const int nq1 = (absPdgId / 1000) % 10;
577 if (nq1 == 0 and nq2 == 5)
const edm::EDGetTokenT< edm::View< reco::Candidate > > genParticleToken_
PseudoTopProducer(const edm::ParameterSet &pset)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
reco::Particle::Point genVertex_
MPlex< T, D1, D2, N > hypot(const MPlex< T, D1, D2, N > &a, const MPlex< T, D1, D2, N > &b)
const double maxLeptonEtaSemilepton_
common ppss p3p6s2 common epss epspn46 common const1 w2
const double maxLeptonEta_
const double minVetoLeptonPtSemilepton_
fastjet::JetDefinition JetDef
const double maxLeptonEtaDilepton_
virtual double pt() const =0
transverse momentum
std::shared_ptr< JetDef > fjJetDef_
virtual int status() const =0
status word
std::vector< GenJet > GenJetCollection
collection of GenJet objects
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, CaloGeometry const &geometry, HcalTopology const &topology)
std::shared_ptr< JetDef > fjLepDef_
const double minMtWSemiLepton_
const double minMETSemiLepton_
reco::Particle::LorentzVector LorentzVector
const double minLeptonPtSemilepton_
bool isNull() const
Checks for null.
bJet1
do the razor with one or two b jets (medium CSV)
bool isNonnull() const
Checks for non-null.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Cos< T >::type cos(const T &t)
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
math::XYZPoint Point
point in the space
Abs< T >::type abs(const T &t)
Jets made from MC generator particles.
#define DEFINE_FWK_MODULE(type)
virtual int charge() const =0
electric charge
virtual int pdgId() const =0
PDG identifier.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
bool isFromHadron(const reco::Candidate *p) const
const double maxVetoLeptonEtaSemilepton_
const double minLeptonPt_
const double minDileptonMassDilepton_
const edm::EDGetTokenT< edm::View< reco::Candidate > > finalStateToken_
const double minLeptonPtDilepton_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
static constexpr float b1
bool isBHadron(const reco::Candidate *p) const