23 #include "fastjet/PseudoJet.hh" 24 #include <fastjet/JetDefinition.hh> 25 #include <TLorentzVector.h> 37 produces<edm::ValueMap<float>>(
"lsf3");
38 produces<edm::ValueMap<int>>(
"muIdx3SJ");
39 produces<edm::ValueMap<int>>(
"eleIdx3SJ");
48 static bool orderPseudoJet(fastjet::PseudoJet j1, fastjet::PseudoJet j2);
49 std::tuple<float, float>
calculateLSF(std::vector<fastjet::PseudoJet> iCParticles,
50 std::vector<fastjet::PseudoJet> &ljets,
68 const auto &eleProd =
iEvent.get(srcEle_);
69 const auto &muProd =
iEvent.get(srcMu_);
71 unsigned int nJet =
srcJet->size();
72 unsigned int nEle = eleProd.size();
73 unsigned int nMu = muProd.size();
75 std::vector<float> vlsf3;
76 std::vector<int> vmuIdx3SJ;
77 std::vector<int> veleIdx3SJ;
80 for (
unsigned int ij = 0; ij < nJet; ij++) {
81 const pat::Jet &itJet = (*srcJet)[ij];
84 std::vector<fastjet::PseudoJet> lClusterParticles;
85 float lepPt(-1), lepEta(-1), lepPhi(-1);
88 fastjet::PseudoJet
p(
d->px(),
d->py(),
d->pz(),
d->energy());
89 lClusterParticles.emplace_back(
p);
91 int ele_pfmatch_index = -1;
92 int mu_pfmatch_index = -1;
95 double dRmin(0.8), dRele(999), dRmu(999), dRtmp(999);
96 for (
unsigned int il(0); il < nEle; il++) {
97 const auto &lep = eleProd.at(il);
100 if (dRtmp <
dRmin && dRtmp < dRele && lep.pt() > lepPt) {
105 ele_pfmatch_index = il;
111 for (
unsigned int il(0); il <
nMu; il++) {
112 const auto &lep = muProd.at(il);
115 if (dRtmp <
dRmin && dRtmp < dRele && dRtmp < dRmu && lep.pt() > lepPt) {
120 ele_pfmatch_index = -1;
121 mu_pfmatch_index = il;
128 std::vector<fastjet::PseudoJet> psub_3;
129 std::sort(lClusterParticles.begin(), lClusterParticles.end(), orderPseudoJet);
130 auto lsf_3 = calculateLSF(lClusterParticles, psub_3, lepPt, lepEta, lepPhi, lepId, 2.0, 3);
131 vlsf3.push_back(std::get<0>(lsf_3));
132 veleIdx3SJ.push_back(ele_pfmatch_index);
133 vmuIdx3SJ.push_back(mu_pfmatch_index);
137 auto lsf3V = std::make_unique<edm::ValueMap<float>>();
139 fillerlsf3.insert(
srcJet, vlsf3.begin(), vlsf3.end());
143 auto muIdx3SJV = std::make_unique<edm::ValueMap<int>>();
145 fillermuIdx3SJ.insert(
srcJet, vmuIdx3SJ.begin(), vmuIdx3SJ.end());
146 fillermuIdx3SJ.fill();
149 auto eleIdx3SJV = std::make_unique<edm::ValueMap<int>>();
151 fillereleIdx3SJ.insert(
srcJet, veleIdx3SJ.begin(), veleIdx3SJ.end());
152 fillereleIdx3SJ.fill();
156 template <
typename T>
158 return j1.perp2() > j2.perp2();
161 template <
typename T>
163 std::vector<fastjet::PseudoJet> &lsubjets,
170 float lsf(-1), lmd(-1);
171 if (ilPt > 0 && (ilId == 11 || ilId == 13)) {
174 ilep.SetPtEtaPhiM(ilPt, ilEta, ilPhi, 0.000511);
176 ilep.SetPtEtaPhiM(ilPt, ilEta, ilPhi, 0.105658);
177 fastjet::JetDefinition lCJet_def(fastjet::kt_algorithm,
dr);
178 fastjet::ClusterSequence lCClust_seq(iCParticles, lCJet_def);
180 lsubjets = sorted_by_pt(lCClust_seq.exclusive_jets_up_to(nsj));
182 lsubjets = sorted_by_pt(lCClust_seq.inclusive_jets());
186 for (
unsigned int i0 = 0; i0 < lsubjets.size(); i0++) {
195 pVec.SetPtEtaPhiM(lsubjets[lId].
pt(), lsubjets[lId].
eta(), lsubjets[lId].
phi(), lsubjets[lId].
m());
196 lsf = ilep.Pt() / pVec.Pt();
197 lmd = (ilep - pVec).M() / pVec.M();
200 return std::tuple<float, float>(lsf, lmd);
203 template <
typename T>
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
double pt() const final
transverse momentum
edm::EDGetTokenT< edm::View< pat::Electron > > srcEle_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
bool matchByCommonSourceCandidatePtr(const C1 &c1, const C2 &c2)
~LeptonInJetProducer() override
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
std::tuple< float, float > calculateLSF(std::vector< fastjet::PseudoJet > iCParticles, std::vector< fastjet::PseudoJet > &ljets, float ilPt, float ilEta, float ilPhi, int ilId, double dr, int nsj) const
#define DEFINE_FWK_MODULE(type)
LeptonInJetProducer< pat::Jet > LepInJetProducer
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
LeptonInJetProducer(const edm::ParameterSet &iConfig)
Analysis-level calorimeter jet class.
const reco::CompositePtrCandidate::daughters & daughterPtrVector() const override
references to daughtes
static bool orderPseudoJet(fastjet::PseudoJet j1, fastjet::PseudoJet j2)
edm::EDGetTokenT< edm::View< pat::Muon > > srcMu_
double phi() const final
momentum azimuthal angle
double eta() const final
momentum pseudorapidity