CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
B2GDQM.cc
Go to the documentation of this file.
2 
3 #include <memory>
4 
5 // DQM
7 
8 // Framework
17 
18 // Candidate handling
25 
26 // Vertex utilities
29 
30 // Other
34 
35 // Math
43 
44 // vertexing
45 
46 // Transient tracks
54 
55 // JetCorrection
57 
58 // ROOT
59 #include "TLorentzVector.h"
60 
61 // STDLIB
62 #include <iostream>
63 #include <iomanip>
64 #include <cstdio>
65 #include <string>
66 #include <sstream>
67 #include <cmath>
68 
69 using namespace edm;
70 using namespace std;
71 using namespace reco;
72 using namespace trigger;
73 
74 //
75 // -- Constructor
76 //
78  edm::LogInfo("B2GDQM") << " Starting B2GDQM "
79  << "\n";
80 
81  typedef std::vector<edm::InputTag> vtag;
82 
83  // Get parameters from configuration file
84  // Trigger
85  theTriggerResultsCollection = ps.getParameter<InputTag>("triggerResultsCollection");
86  triggerToken_ = consumes<edm::TriggerResults>(theTriggerResultsCollection);
87 
88  // Jets
89  jetLabels_ = ps.getParameter<std::vector<edm::InputTag> >("jetLabels");
90  for (std::vector<edm::InputTag>::const_iterator jetlabel = jetLabels_.begin(), jetlabelEnd = jetLabels_.end();
91  jetlabel != jetlabelEnd;
92  ++jetlabel) {
93  jetTokens_.push_back(consumes<edm::View<reco::Jet> >(*jetlabel));
94  }
95  sdjetLabel_ = ps.getParameter<edm::InputTag>("sdjetLabel");
96  sdjetToken_ = consumes<edm::View<reco::BasicJet> >(sdjetLabel_);
97 
98  muonToken_ = consumes<edm::View<reco::Muon> >(ps.getParameter<edm::InputTag>("muonSrc"));
99  electronToken_ = consumes<edm::View<reco::GsfElectron> >(ps.getParameter<edm::InputTag>("electronSrc"));
100 
101  jetPtMins_ = ps.getParameter<std::vector<double> >("jetPtMins");
102  allHadPtCut_ = ps.getParameter<double>("allHadPtCut");
103  allHadRapidityCut_ = ps.getParameter<double>("allHadRapidityCut");
104  allHadDeltaPhiCut_ = ps.getParameter<double>("allHadDeltaPhiCut");
105 
106  semiMu_HadJetPtCut_ = ps.getParameter<double>("semiMu_HadJetPtCut");
107  semiMu_LepJetPtCut_ = ps.getParameter<double>("semiMu_LepJetPtCut");
108  semiMu_dphiHadCut_ = ps.getParameter<double>("semiMu_dphiHadCut");
109  semiMu_dRMin_ = ps.getParameter<double>("semiMu_dRMin");
110  semiMu_ptRel_ = ps.getParameter<double>("semiMu_ptRel");
111  muonSelect_ = std::make_shared<StringCutObjectSelector<reco::Muon> >(
112 
113  ps.getParameter<std::string>("muonSelect"));
114 
115  semiE_HadJetPtCut_ = ps.getParameter<double>("semiE_HadJetPtCut");
116  semiE_LepJetPtCut_ = ps.getParameter<double>("semiE_LepJetPtCut");
117  semiE_dphiHadCut_ = ps.getParameter<double>("semiE_dphiHadCut");
118  semiE_dRMin_ = ps.getParameter<double>("semiE_dRMin");
119  semiE_ptRel_ = ps.getParameter<double>("semiE_ptRel");
120  elecSelect_ = std::make_shared<StringCutObjectSelector<reco::GsfElectron> >(
121 
122  ps.getParameter<std::string>("elecSelect"));
123 
124  PFJetCorService_ = ps.getParameter<std::string>("PFJetCorService");
125 
126  // MET
127  PFMETLabel_ = ps.getParameter<InputTag>("pfMETCollection");
128  PFMETToken_ = consumes<std::vector<reco::PFMET> >(PFMETLabel_);
129 }
130 
131 //
132 // -- Destructor
133 //
135  edm::LogInfo("B2GDQM") << " Deleting B2GDQM "
136  << "\n";
137 }
138 
139 //
140 // -- Book histograms
141 //
143  bei.setCurrentFolder("Physics/B2G");
144 
145  //--- Jets
146 
147  for (unsigned int icoll = 0; icoll < jetLabels_.size(); ++icoll) {
148  std::stringstream ss;
149  ss << "Physics/B2G/" << jetLabels_[icoll].label();
150  bei.setCurrentFolder(ss.str());
151  pfJet_pt.push_back(bei.book1D("pfJet_pt", "Pt of PFJet (GeV)", 50, 0.0, 1000));
152  pfJet_y.push_back(bei.book1D("pfJet_y", "Rapidity of PFJet", 60, -6.0, 6.0));
153  pfJet_phi.push_back(bei.book1D("pfJet_phi", "#phi of PFJet (radians)", 60, -3.14159, 3.14159));
154  pfJet_m.push_back(bei.book1D("pfJet_m", "Mass of PFJet (GeV)", 50, 0.0, 500));
155  pfJet_chef.push_back(bei.book1D("pfJet_pfchef", "PFJetID CHEF", 50, 0.0, 1.0));
156  pfJet_nhef.push_back(bei.book1D("pfJet_pfnhef", "PFJetID NHEF", 50, 0.0, 1.0));
157  pfJet_cemf.push_back(bei.book1D("pfJet_pfcemf", "PFJetID CEMF", 50, 0.0, 1.0));
158  pfJet_nemf.push_back(bei.book1D("pfJet_pfnemf", "PFJetID NEMF", 50, 0.0, 1.0));
159 
160  boostedJet_subjetPt.push_back(bei.book1D("boostedJet_subjetPt", "Pt of subjets (GeV)", 50, 0.0, 500));
161  boostedJet_subjetY.push_back(bei.book1D("boostedJet_subjetY", "Rapidity of subjets", 60, -6.0, 6.0));
162  boostedJet_subjetPhi.push_back(
163  bei.book1D("boostedJet_subjetPhi", "#phi of subjets (radians)", 60, -3.14159, 3.14159));
164  boostedJet_subjetM.push_back(bei.book1D("boostedJet_subjetM", "Mass of subjets (GeV)", 50, 0.0, 250.));
165  boostedJet_subjetN.push_back(bei.book1D("boostedJet_subjetN", "Number of subjets", 10, 0, 10));
166  boostedJet_massDrop.push_back(bei.book1D("boostedJet_massDrop", "Mass drop for W-like jets", 50, 0.0, 1.0));
167  boostedJet_wMass.push_back(bei.book1D("boostedJet_wMass", "W Mass for top-like jets", 50, 0.0, 250.0));
168  }
169 
170  bei.setCurrentFolder("Physics/B2G/MET");
171  pfMet_pt = bei.book1D("pfMet_pt", "Pf Missing p_{T}; GeV", 50, 0.0, 500);
172  pfMet_phi = bei.book1D("pfMet_phi", "Pf Missing p_{T} #phi;#phi (radians)", 35, -3.5, 3.5);
173 
174  //--- Mu+Jets
175  bei.setCurrentFolder("Physics/B2G/SemiMu");
176  semiMu_muPt = bei.book1D("semiMu_muPt", "Pt of Muon in #mu+Jets Channel (GeV)", 50, 0.0, 1000);
177  semiMu_muEta = bei.book1D("semiMu_muEta", "#eta of Muon in #mu+Jets Channel", 60, -6.0, 6.0);
178  semiMu_muPhi = bei.book1D("semiMu_muPhi", "#phi of Muon in #mu+Jets Channel (radians)", 60, -3.14159, 3.14159);
179  semiMu_muDRMin = bei.book1D("semiMu_muDRMin", "#Delta R(E,nearest jet) in #mu+Jets Channel", 50, 0, 10.0);
180  semiMu_muPtRel = bei.book1D("semiMu_muPtRel", "p_{T}^{REL} in #mu+Jets Channel", 60, 0, 300.);
181  semiMu_hadJetDR = bei.book1D("semiMu_hadJetDR", "#Delta R(E,had jet) in #mu+Jets Channel", 50, 0, 10.0);
182  semiMu_hadJetPt =
183  bei.book1D("semiMu_hadJetPt", "Pt of Leading Hadronic Jet in #mu+Jets Channel (GeV)", 50, 0.0, 1000);
184  semiMu_hadJetY = bei.book1D("semiMu_hadJetY", "Rapidity of Leading Hadronic Jet in #mu+Jets Channel", 60, -6.0, 6.0);
185  semiMu_hadJetPhi = bei.book1D(
186  "semiMu_hadJetPhi", "#phi of Leading Hadronic Jet in #mu+Jets Channel (radians)", 60, -3.14159, 3.14159);
187  semiMu_hadJetMass =
188  bei.book1D("semiMu_hadJetMass", "Mass of Leading Hadronic Jet in #mu+Jets Channel (GeV)", 50, 0.0, 500);
189  semiMu_hadJetWMass =
190  bei.book1D("semiMu_hadJetwMass", "W Mass for Leading Hadronic Jet in #mu+Jets Channel (GeV)", 50, 0.0, 250.0);
191  semiMu_mttbar = bei.book1D("semiMu_mttbar", "Mass of #mu+Jets ttbar Candidate", 100, 0., 5000.);
192 
193  //--- E+Jets
194  bei.setCurrentFolder("Physics/B2G/SemiE");
195  semiE_ePt = bei.book1D("semiE_ePt", "Pt of Electron in e+Jets Channel (GeV)", 50, 0.0, 1000);
196  semiE_eEta = bei.book1D("semiE_eEta", "#eta of Electron in e+Jets Channel", 60, -6.0, 6.0);
197  semiE_ePhi = bei.book1D("semiE_ePhi", "#phi of Electron in e+Jets Channel (radians)", 60, -3.14159, 3.14159);
198  semiE_eDRMin = bei.book1D("semiE_eDRMin", "#Delta R(E,nearest jet) in e+Jets Channel", 50, 0, 10.0);
199  semiE_ePtRel = bei.book1D("semiE_ePtRel", "p_{T}^{REL} in e+Jets Channel", 60, 0, 300.);
200  semiE_hadJetDR = bei.book1D("semiE_hadJetDR", "#Delta R(E,had jet) in e+Jets Channel", 50, 0, 10.0);
201  semiE_hadJetPt = bei.book1D("semiE_hadJetPt", "Pt of Leading Hadronic Jet in e+Jets Channel (GeV)", 50, 0.0, 1000);
202  semiE_hadJetY = bei.book1D("semiE_hadJetY", "Rapidity of Leading Hadronic Jet in e+Jets Channel", 60, -6.0, 6.0);
203  semiE_hadJetPhi =
204  bei.book1D("semiE_hadJetPhi", "#phi of Leading Hadronic Jet in e+Jets Channel (radians)", 60, -3.14159, 3.14159);
205  semiE_hadJetMass =
206  bei.book1D("semiE_hadJetMass", "Mass of Leading Hadronic Jet in e+Jets Channel (GeV)", 50, 0.0, 500);
207  semiE_hadJetWMass =
208  bei.book1D("semiE_hadJetwMass", "W Mass for Leading Hadronic Jet in e+Jets Channel (GeV)", 50, 0.0, 250.0);
209  semiE_mttbar = bei.book1D("semiE_mttbar", "Mass of e+Jets ttbar Candidate", 100, 0., 5000.);
210 
211  //--- All-hadronic
212  bei.setCurrentFolder("Physics/B2G/AllHad");
213  allHad_pt0 = bei.book1D("allHad_pt0", "Pt of Leading All-Hadronic PFJet (GeV)", 50, 0.0, 1000);
214  allHad_y0 = bei.book1D("allHad_y0", "Rapidity of Leading All-Hadronic PFJet", 60, -6.0, 6.0);
215  allHad_phi0 = bei.book1D("allHad_phi0", "#phi of Leading All-Hadronic PFJet (radians)", 60, -3.14159, 3.14159);
216  allHad_mass0 = bei.book1D("allHad_mass0", "Mass of Leading All-Hadronic PFJet (GeV)", 50, 0.0, 500);
217  allHad_wMass0 = bei.book1D("allHad_wMass0", "W Mass for Leading All-Hadronic PFJet (GeV)", 50, 0.0, 250.0);
218  allHad_pt1 = bei.book1D("allHad_pt1", "Pt of Subleading All-Hadronic PFJet (GeV)", 50, 0.0, 1000);
219  allHad_y1 = bei.book1D("allHad_y1", "Rapidity of Subleading All-Hadronic PFJet", 60, -6.0, 6.0);
220  allHad_phi1 = bei.book1D("allHad_phi1", "#phi of Subleading All-Hadronic PFJet (radians)", 60, -3.14159, 3.14159);
221  allHad_mass1 = bei.book1D("allHad_mass1", "Mass of Subleading All-Hadronic PFJet (GeV)", 50, 0.0, 500);
222  allHad_wMass1 = bei.book1D("allHad_wMass1", "W Mass for Subleading All-Hadronic PFJet (GeV)", 50, 0.0, 250.0);
223  allHad_mttbar = bei.book1D("allHad_mttbar", "Mass of All-Hadronic ttbar Candidate", 100, 0., 5000.);
224 }
225 
226 //
227 // -- Analyze
228 //
229 void B2GDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
230  analyzeJets(iEvent, iSetup);
231  analyzeSemiMu(iEvent, iSetup);
232  analyzeSemiE(iEvent, iSetup);
233  analyzeAllHad(iEvent, iSetup);
234 }
235 
236 void B2GDQM::analyzeJets(const Event& iEvent, const edm::EventSetup& iSetup) {
237  // Loop over the different types of jets,
238  // Loop over the jets in that collection,
239  // fill PF jet information as well as substructure
240  // information for boosted jets.
241  // Utilizes the CMS top-tagging algorithm and the "mass drop" W-tagger.
242  for (unsigned int icoll = 0; icoll < jetLabels_.size(); ++icoll) {
244  bool ValidPFJets = iEvent.getByToken(jetTokens_[icoll], pfJetCollection);
245  if (!ValidPFJets)
246  continue;
247  edm::View<reco::Jet> const& pfjets = *pfJetCollection;
248 
249  // Jet Correction
250  int countJet = 0;
251  // const JetCorrector* pfcorrector =
252  // JetCorrector::getJetCorrector(PFJetCorService_,iSetup);
253 
254  for (edm::View<reco::Jet>::const_iterator jet = pfjets.begin(), jetEnd = pfjets.end(); jet != jetEnd; ++jet) {
255  if (jet->pt() < jetPtMins_[icoll])
256  continue;
257  pfJet_pt[icoll]->Fill(jet->pt());
258  pfJet_y[icoll]->Fill(jet->rapidity());
259  pfJet_phi[icoll]->Fill(jet->phi());
260  pfJet_m[icoll]->Fill(jet->mass());
261 
262  // Dynamic cast the base class (reco::Jet) to the derived class (PFJet)
263  // to access the PFJet information
264  reco::PFJet const* pfjet = dynamic_cast<reco::PFJet const*>(&*jet);
265 
266  if (pfjet != nullptr) {
267  pfJet_chef[icoll]->Fill(pfjet->chargedHadronEnergyFraction());
268  pfJet_nhef[icoll]->Fill(pfjet->neutralHadronEnergyFraction());
269  pfJet_cemf[icoll]->Fill(pfjet->chargedEmEnergyFraction());
270  pfJet_nemf[icoll]->Fill(pfjet->neutralEmEnergyFraction());
271  }
272 
273  // Dynamic cast the base class (reco::Jet) to the derived class (BasicJet)
274  // to access the substructure information
275  reco::BasicJet const* basicjet = dynamic_cast<reco::BasicJet const*>(&*jet);
276 
277  if (basicjet != nullptr) {
278  boostedJet_subjetN[icoll]->Fill(jet->numberOfDaughters());
279 
280  for (unsigned int ida = 0; ida < jet->numberOfDaughters(); ++ida) {
281  reco::Candidate const* subjet = jet->daughter(ida);
282  boostedJet_subjetPt[icoll]->Fill(subjet->pt());
283  boostedJet_subjetY[icoll]->Fill(subjet->rapidity());
284  boostedJet_subjetPhi[icoll]->Fill(subjet->phi());
285  boostedJet_subjetM[icoll]->Fill(subjet->mass());
286  }
287  // Check the various tagging algorithms
288  if ((jetLabels_[icoll].label() == "ak8PFJetsPuppiSoftdrop")) {
289  if (jet->numberOfDaughters() > 1) {
290  reco::Candidate const* da0 = jet->daughter(0);
291  reco::Candidate const* da1 = jet->daughter(1);
292  if (da0->mass() > da1->mass()) {
293  boostedJet_wMass[icoll]->Fill(da0->mass());
294  boostedJet_massDrop[icoll]->Fill(da0->mass() / jet->mass());
295  } else {
296  boostedJet_wMass[icoll]->Fill(da1->mass());
297  boostedJet_massDrop[icoll]->Fill(da1->mass() / jet->mass());
298  }
299 
300  } else {
301  boostedJet_massDrop[icoll]->Fill(-1.0);
302  }
303 
304  } // end if collection is AK8 PFJets Puppi soft-drop
305 
306  } // end if basic jet != 0
307  countJet++;
308  }
309  }
310 
311  // PFMETs
312  edm::Handle<std::vector<reco::PFMET> > pfMETCollection;
313  bool ValidPFMET = iEvent.getByToken(PFMETToken_, pfMETCollection);
314  if (!ValidPFMET)
315  return;
316 
317  pfMet_pt->Fill((*pfMETCollection)[0].pt());
318  pfMet_phi->Fill((*pfMETCollection)[0].phi());
319 }
320 
321 void B2GDQM::analyzeAllHad(const Event& iEvent, const edm::EventSetup& iSetup) {
323  bool validJets = iEvent.getByToken(sdjetToken_, jetCollection);
324  if (!validJets)
325  return;
326 
327  // Require two back-to-back jets at high pt with |delta y| < 1.0
328  if (jetCollection->size() < 2)
329  return;
330  edm::Ptr<reco::BasicJet> jet0 = jetCollection->ptrAt(0);
331  edm::Ptr<reco::BasicJet> jet1 = jetCollection->ptrAt(1);
332  if (jet0.isAvailable() == false || jet1.isAvailable() == false)
333  return;
334  if (jet0->pt() < allHadPtCut_ || jet1->pt() < allHadPtCut_)
335  return;
336  if (std::abs(jet0->rapidity() - jet1->rapidity()) > allHadRapidityCut_)
337  return;
338  if (std::abs(reco::deltaPhi(jet0->phi(), jet1->phi())) < M_PI * 0.5)
339  return;
340 
341  allHad_pt0->Fill(jet0->pt());
342  allHad_y0->Fill(jet0->rapidity());
343  allHad_phi0->Fill(jet0->phi());
344  allHad_mass0->Fill(jet0->mass());
345  if (jet0->numberOfDaughters() > 2) {
346  double wMass =
347  jet0->daughter(0)->mass() >= jet0->daughter(1)->mass() ? jet0->daughter(0)->mass() : jet0->daughter(1)->mass();
348  allHad_wMass0->Fill(wMass);
349  } else {
350  allHad_wMass0->Fill(-1.0);
351  }
352 
353  allHad_pt1->Fill(jet1->pt());
354  allHad_y1->Fill(jet1->rapidity());
355  allHad_phi1->Fill(jet1->phi());
356  allHad_mass1->Fill(jet1->mass());
357  if (jet1->numberOfDaughters() > 2) {
358  double wMass =
359  jet1->daughter(0)->mass() >= jet1->daughter(1)->mass() ? jet1->daughter(0)->mass() : jet1->daughter(1)->mass();
360  allHad_wMass1->Fill(wMass);
361  } else {
362  allHad_wMass1->Fill(-1.0);
363  }
364 
365  auto p4cand = (jet0->p4() + jet1->p4());
366  allHad_mttbar->Fill(p4cand.mass());
367 }
368 
369 void B2GDQM::analyzeSemiMu(const Event& iEvent, const edm::EventSetup& iSetup) {
370  edm::Handle<edm::View<reco::Muon> > muonCollection;
371  bool validMuons = iEvent.getByToken(muonToken_, muonCollection);
372 
373  if (!validMuons)
374  return;
375  if (muonCollection->empty())
376  return;
377  reco::Muon const& muon = muonCollection->at(0);
378  if (!(*muonSelect_)(muon))
379  return;
380 
382  bool validJets = iEvent.getByToken(sdjetToken_, jetCollection);
383  if (!validJets)
384  return;
385  if (jetCollection->size() < 2)
386  return;
387 
388  double pt0 = -1.0;
389  double dRMin = 999.0;
390  edm::Ptr<reco::BasicJet> hadJet; // highest pt jet with dphi(lep,jet) > pi/2
391  edm::Ptr<reco::BasicJet> lepJet; // closest jet to lepton with pt > ptMin
392 
393  for (auto ijet = jetCollection->begin(), ijetBegin = ijet, ijetEnd = jetCollection->end(); ijet != ijetEnd; ++ijet) {
394  // Hadronic jets
395  if (std::abs(reco::deltaPhi(muon, *ijet)) > M_PI * 0.5) {
396  if (ijet->pt() > pt0 && ijet->p() > semiMu_HadJetPtCut_) {
397  hadJet = jetCollection->ptrAt(ijet - ijetBegin);
398  pt0 = hadJet->pt();
399  }
400  }
401  // Leptonic jets
402  else if (ijet->pt() > semiMu_LepJetPtCut_) {
403  auto idRMin = reco::deltaR(muon, *ijet);
404  if (idRMin < dRMin) {
405  lepJet = jetCollection->ptrAt(ijet - ijetBegin);
406  dRMin = idRMin;
407  }
408  }
409  }
410  if (hadJet.isAvailable() == false || lepJet.isAvailable() == false)
411  return;
412 
413  auto lepJetP4 = lepJet->p4();
414  const auto& muonP4 = muon.p4();
415 
416  double tot = lepJetP4.mag2();
417  double ss = muonP4.Dot(lepJet->p4());
418  double per = muonP4.mag2();
419  if (tot > 0.0)
420  per -= ss * ss / tot;
421  if (per < 0)
422  per = 0;
423  double ptRel = per;
424  bool pass2D = dRMin > semiMu_dRMin_ || ptRel > semiMu_ptRel_;
425 
426  if (!pass2D)
427  return;
428 
429  semiMu_muPt->Fill(muon.pt());
430  semiMu_muEta->Fill(muon.eta());
431  semiMu_muPhi->Fill(muon.phi());
432  semiMu_muDRMin->Fill(dRMin);
433  semiMu_muPtRel->Fill(ptRel);
434 
435  semiMu_hadJetDR->Fill(reco::deltaR(muon, *hadJet));
436  semiMu_mttbar->Fill(0.0);
437 
438  semiMu_hadJetPt->Fill(hadJet->pt());
439  semiMu_hadJetY->Fill(hadJet->rapidity());
440  semiMu_hadJetPhi->Fill(hadJet->phi());
441  semiMu_hadJetMass->Fill(hadJet->mass());
442  if (hadJet->numberOfDaughters() > 2) {
443  double wMass = hadJet->daughter(0)->mass() >= hadJet->daughter(1)->mass() ? hadJet->daughter(0)->mass()
444  : hadJet->daughter(1)->mass();
445  semiMu_hadJetWMass->Fill(wMass);
446  } else {
447  semiMu_hadJetWMass->Fill(-1.0);
448  }
449 }
450 
451 void B2GDQM::analyzeSemiE(const Event& iEvent, const edm::EventSetup& iSetup) {
452  edm::Handle<edm::View<reco::GsfElectron> > electronCollection;
453  bool validElectrons = iEvent.getByToken(electronToken_, electronCollection);
454 
455  if (!validElectrons)
456  return;
457  if (electronCollection->empty())
458  return;
459  reco::GsfElectron const& electron = electronCollection->at(0);
460  if (!(*elecSelect_)(electron))
461  return;
462 
464  bool validJets = iEvent.getByToken(sdjetToken_, jetCollection);
465  if (!validJets)
466  return;
467  if (jetCollection->size() < 2)
468  return;
469 
470  double pt0 = -1.0;
471  double dRMin = 999.0;
472  edm::Ptr<reco::BasicJet> hadJet; // highest pt jet with dphi(lep,jet) > pi/2
473  edm::Ptr<reco::BasicJet> lepJet; // closest jet to lepton with pt > ptMin
474 
475  for (auto ijet = jetCollection->begin(), ijetBegin = ijet, ijetEnd = jetCollection->end(); ijet != ijetEnd; ++ijet) {
476  // Hadronic jets
477  if (std::abs(reco::deltaPhi(electron, *ijet)) > M_PI * 0.5) {
478  if (ijet->pt() > pt0 && ijet->p() > semiE_HadJetPtCut_) {
479  hadJet = jetCollection->ptrAt(ijet - ijetBegin);
480  pt0 = hadJet->pt();
481  }
482  }
483  // Leptonic jets
484  else if (ijet->pt() > semiE_LepJetPtCut_) {
485  auto idRMin = reco::deltaR(electron, *ijet);
486  if (idRMin < dRMin) {
487  lepJet = jetCollection->ptrAt(ijet - ijetBegin);
488  dRMin = idRMin;
489  }
490  }
491  }
492  if (hadJet.isAvailable() == false || lepJet.isAvailable() == false)
493  return;
494 
495  auto lepJetP4 = lepJet->p4();
496  const auto& electronP4 = electron.p4();
497 
498  double tot = lepJetP4.mag2();
499  double ss = electronP4.Dot(lepJet->p4());
500  double per = electronP4.mag2();
501  if (tot > 0.0)
502  per -= ss * ss / tot;
503  if (per < 0)
504  per = 0;
505  double ptRel = per;
506  bool pass2D = dRMin > semiE_dRMin_ || ptRel > semiE_ptRel_;
507 
508  if (!pass2D)
509  return;
510 
511  semiE_ePt->Fill(electron.pt());
512  semiE_eEta->Fill(electron.eta());
513  semiE_ePhi->Fill(electron.phi());
514  semiE_eDRMin->Fill(dRMin);
515  semiE_ePtRel->Fill(ptRel);
516 
517  semiE_hadJetDR->Fill(reco::deltaR(electron, *hadJet));
518  semiE_mttbar->Fill(0.0);
519 
520  semiE_hadJetPt->Fill(hadJet->pt());
521  semiE_hadJetY->Fill(hadJet->rapidity());
522  semiE_hadJetPhi->Fill(hadJet->phi());
523  semiE_hadJetMass->Fill(hadJet->mass());
524  if (hadJet->numberOfDaughters() > 2) {
525  double wMass = hadJet->daughter(0)->mass() >= hadJet->daughter(1)->mass() ? hadJet->daughter(0)->mass()
526  : hadJet->daughter(1)->mass();
527  semiE_hadJetWMass->Fill(wMass);
528  } else {
529  semiE_hadJetWMass->Fill(-1.0);
530  }
531 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
virtual void analyzeSemiMu(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:369
~B2GDQM() override
Definition: B2GDQM.cc:134
virtual void analyzeAllHad(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:321
double pt() const final
transverse momentum
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:217
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
virtual double pt() const =0
transverse momentum
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:97
virtual double mass() const =0
mass
bool isAvailable() const
Definition: Ptr.h:230
virtual double rapidity() const =0
rapidity
void bookHistograms(DQMStore::IBooker &bei, edm::Run const &, edm::EventSetup const &) override
Definition: B2GDQM.cc:142
Jets made from CaloTowers.
Definition: BasicJet.h:19
Jets made from PFObjects.
Definition: PFJet.h:20
tuple pfJetCollection
char const * label
float neutralEmEnergyFraction() const
neutralEmEnergyFraction
Definition: PFJet.h:149
int iEvent
Definition: GenABIO.cc:224
const_iterator begin() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:101
virtual void analyzeSemiE(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:451
#define M_PI
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Log< level::Info, false > LogInfo
float chargedEmEnergyFraction() const
chargedEmEnergyFraction
Definition: PFJet.h:141
B2GDQM(const edm::ParameterSet &ps)
Definition: B2GDQM.cc:77
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
Definition: B2GDQM.cc:229
virtual void analyzeJets(edm::Event const &e, edm::EventSetup const &eSetup)
Definition: B2GDQM.cc:236
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const_iterator end() const
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
double phi() const final
momentum azimuthal angle
Definition: Run.h:45
virtual double phi() const =0
momentum azimuthal angle
double eta() const final
momentum pseudorapidity
int icoll
Definition: AMPTWrapper.h:146