CMS 3D CMS Logo

LeptonInJetProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
7 
10 
13 
16 
21 
23 #include "fastjet/PseudoJet.hh"
24 #include <fastjet/JetDefinition.hh>
25 #include <TLorentzVector.h>
26 #include <TMath.h>
27 
29 
30 template <typename T>
32 public:
33  explicit LeptonInJetProducer(const edm::ParameterSet &iConfig)
34  : srcJet_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("src"))),
35  srcEle_(consumes<edm::View<pat::Electron>>(iConfig.getParameter<edm::InputTag>("srcEle"))),
36  srcMu_(consumes<edm::View<pat::Muon>>(iConfig.getParameter<edm::InputTag>("srcMu"))) {
37  produces<edm::ValueMap<float>>("lsf3");
38  produces<edm::ValueMap<int>>("muIdx3SJ");
39  produces<edm::ValueMap<int>>("eleIdx3SJ");
40  }
41  ~LeptonInJetProducer() override{};
42 
43  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
44 
45 private:
46  void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override;
47 
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,
51  float ilPt,
52  float ilEta,
53  float ilPhi,
54  int ilId,
55  double dr,
56  int nsj) const;
57 
61 };
62 
63 // ------------ method called to produce the data ------------
64 template <typename T>
66  // needs jet collection (srcJet), leptons collection
67  auto srcJet = iEvent.getHandle(srcJet_);
68  const auto &eleProd = iEvent.get(srcEle_);
69  const auto &muProd = iEvent.get(srcMu_);
70 
71  unsigned int nJet = srcJet->size();
72  unsigned int nEle = eleProd.size();
73  unsigned int nMu = muProd.size();
74 
75  std::vector<float> vlsf3;
76  std::vector<int> vmuIdx3SJ;
77  std::vector<int> veleIdx3SJ;
78 
79  // Find leptons in jets
80  for (unsigned int ij = 0; ij < nJet; ij++) {
81  const pat::Jet &itJet = (*srcJet)[ij];
82  if (itJet.pt() <= 10)
83  continue;
84  std::vector<fastjet::PseudoJet> lClusterParticles;
85  float lepPt(-1), lepEta(-1), lepPhi(-1);
86  int lepId(-1);
87  for (auto const &d : itJet.daughterPtrVector()) {
88  fastjet::PseudoJet p(d->px(), d->py(), d->pz(), d->energy());
89  lClusterParticles.emplace_back(p);
90  }
91 
92  int ele_pfmatch_index = -1;
93  int mu_pfmatch_index = -1;
94 
95  // match to leading and closest electron or muon
96  double dRmin(0.8), dRele(999), dRmu(999), dRtmp(999);
97  for (unsigned int il(0); il < nEle; il++) {
98  const auto &lep = eleProd.at(il);
99  if (matchByCommonSourceCandidatePtr(lep, itJet)) {
100  dRtmp = reco::deltaR(itJet.eta(), itJet.phi(), lep.eta(), lep.phi());
101  if (dRtmp < dRmin && dRtmp < dRele && lep.pt() > lepPt) {
102  lepPt = lep.pt();
103  lepEta = lep.eta();
104  lepPhi = lep.phi();
105  lepId = 11;
106  ele_pfmatch_index = il;
107  dRele = dRtmp;
108  break;
109  }
110  }
111  }
112  for (unsigned int il(0); il < nMu; il++) {
113  const auto &lep = muProd.at(il);
114  if (matchByCommonSourceCandidatePtr(lep, itJet)) {
115  dRtmp = reco::deltaR(itJet.eta(), itJet.phi(), lep.eta(), lep.phi());
116  if (dRtmp < dRmin && dRtmp < dRele && dRtmp < dRmu && lep.pt() > lepPt) {
117  lepPt = lep.pt();
118  lepEta = lep.eta();
119  lepPhi = lep.phi();
120  lepId = 13;
121  ele_pfmatch_index = -1;
122  mu_pfmatch_index = il;
123  dRmu = dRtmp;
124  break;
125  }
126  }
127  }
128 
129  std::vector<fastjet::PseudoJet> psub_3;
130  std::sort(lClusterParticles.begin(), lClusterParticles.end(), orderPseudoJet);
131  auto lsf_3 = calculateLSF(lClusterParticles, psub_3, lepPt, lepEta, lepPhi, lepId, 2.0, 3);
132  vlsf3.push_back(std::get<0>(lsf_3));
133  veleIdx3SJ.push_back(ele_pfmatch_index);
134  vmuIdx3SJ.push_back(mu_pfmatch_index);
135  }
136 
137  // Filling table
138  auto lsf3V = std::make_unique<edm::ValueMap<float>>();
139  edm::ValueMap<float>::Filler fillerlsf3(*lsf3V);
140  fillerlsf3.insert(srcJet, vlsf3.begin(), vlsf3.end());
141  fillerlsf3.fill();
142  iEvent.put(std::move(lsf3V), "lsf3");
143 
144  auto muIdx3SJV = std::make_unique<edm::ValueMap<int>>();
145  edm::ValueMap<int>::Filler fillermuIdx3SJ(*muIdx3SJV);
146  fillermuIdx3SJ.insert(srcJet, vmuIdx3SJ.begin(), vmuIdx3SJ.end());
147  fillermuIdx3SJ.fill();
148  iEvent.put(std::move(muIdx3SJV), "muIdx3SJ");
149 
150  auto eleIdx3SJV = std::make_unique<edm::ValueMap<int>>();
151  edm::ValueMap<int>::Filler fillereleIdx3SJ(*eleIdx3SJV);
152  fillereleIdx3SJ.insert(srcJet, veleIdx3SJ.begin(), veleIdx3SJ.end());
153  fillereleIdx3SJ.fill();
154  iEvent.put(std::move(eleIdx3SJV), "eleIdx3SJ");
155 }
156 
157 template <typename T>
158 bool LeptonInJetProducer<T>::orderPseudoJet(fastjet::PseudoJet j1, fastjet::PseudoJet j2) {
159  return j1.perp2() > j2.perp2();
160 }
161 
162 template <typename T>
163 std::tuple<float, float> LeptonInJetProducer<T>::calculateLSF(std::vector<fastjet::PseudoJet> iCParticles,
164  std::vector<fastjet::PseudoJet> &lsubjets,
165  float ilPt,
166  float ilEta,
167  float ilPhi,
168  int ilId,
169  double dr,
170  int nsj) const {
171  float lsf(-1), lmd(-1);
172  if (ilPt > 0 && (ilId == 11 || ilId == 13)) {
173  TLorentzVector ilep;
174  if (ilId == 11)
175  ilep.SetPtEtaPhiM(ilPt, ilEta, ilPhi, 0.000511);
176  if (ilId == 13)
177  ilep.SetPtEtaPhiM(ilPt, ilEta, ilPhi, 0.105658);
178  fastjet::JetDefinition lCJet_def(fastjet::kt_algorithm, dr);
179  fastjet::ClusterSequence lCClust_seq(iCParticles, lCJet_def);
180  if (dr > 0.5) {
181  lsubjets = sorted_by_pt(lCClust_seq.exclusive_jets_up_to(nsj));
182  } else {
183  lsubjets = sorted_by_pt(lCClust_seq.inclusive_jets());
184  }
185  int lId(-1);
186  double dRmin = 999.;
187  for (unsigned int i0 = 0; i0 < lsubjets.size(); i0++) {
188  double dR = reco::deltaR(lsubjets[i0].eta(), lsubjets[i0].phi(), ilep.Eta(), ilep.Phi());
189  if (dR < dRmin) {
190  dRmin = dR;
191  lId = i0;
192  }
193  }
194  if (lId != -1) {
195  TLorentzVector pVec;
196  pVec.SetPtEtaPhiM(lsubjets[lId].pt(), lsubjets[lId].eta(), lsubjets[lId].phi(), lsubjets[lId].m());
197  lsf = ilep.Pt() / pVec.Pt();
198  lmd = (ilep - pVec).M() / pVec.M();
199  }
200  }
201  return std::tuple<float, float>(lsf, lmd);
202 }
203 
204 template <typename T>
207  desc.add<edm::InputTag>("src")->setComment("jet input collection");
208  desc.add<edm::InputTag>("srcEle")->setComment("electron input collection");
209  desc.add<edm::InputTag>("srcMu")->setComment("muon input collection");
210  descriptions.addWithDefaultLabel(desc);
211 }
212 
214 
215 //define this as a plug-in
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)
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
bool matchByCommonSourceCandidatePtr(const C1 &c1, const C2 &c2)
Definition: MatchingUtils.h:9
Definition: Muon.py:1
Definition: Jet.py:1
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)
Definition: MakerMacros.h:16
LeptonInJetProducer< pat::Jet > LepInJetProducer
d
Definition: ztail.py:151
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
LeptonInJetProducer(const edm::ParameterSet &iConfig)
Analysis-level calorimeter jet class.
Definition: Jet.h:77
const reco::CompositePtrCandidate::daughters & daughterPtrVector() const override
references to daughtes
HLT enums.
static bool orderPseudoJet(fastjet::PseudoJet j1, fastjet::PseudoJet j2)
edm::EDGetTokenT< edm::View< pat::Muon > > srcMu_
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511
double eta() const final
momentum pseudorapidity