CMS 3D CMS Logo

PatBasicFWLiteJetAnalyzer_Selector.cc
Go to the documentation of this file.
1 /* A macro for making a histogram of Jet Pt with cuts
2 This is a basic way to cut out jets of a certain Pt and Eta using an if statement
3 This example creates a histogram of Jet Pt, using Jets with Pt above 30 and ETA above -2.1 and below 2.1
4 */
5 
6 #include "TFile.h"
7 #include "TH1.h"
8 #include "TH2.h"
9 #include "TCanvas.h"
10 #include "TLegend.h"
11 #include "TSystem.h"
12 #include "TLorentzVector.h"
13 
30 
31 #include <iostream>
32 #include <memory>
33 #include <cmath> //necessary for absolute function fabs()
34 
35 using namespace std;
36 
40 
41 // -*- C++ -*-
42 
43 // CMS includes
49 
50 // Root includes
51 #include "TROOT.h"
52 
53 using namespace std;
54 
56 public:
57  JetIDStudiesSelector(edm::ParameterSet const& caloJetIdParams,
58  edm::ParameterSet const& pfJetIdParams,
60  : jetSel_(new JetIDSelectionFunctor(caloJetIdParams)),
61  pfJetSel_(new PFJetIDSelectionFunctor(pfJetIdParams)),
62  jetSrc_(params.getParameter<edm::InputTag>("jetSrc")),
63  pfJetSrc_(params.getParameter<edm::InputTag>("pfJetSrc")) {
64  bool useCalo = params.getParameter<bool>("useCalo");
65 
66  push_back("Calo Cuts");
67  push_back("Calo Kin Cuts");
68  push_back("Calo Delta Phi");
69  push_back("Calo Jet ID");
70  push_back("PF Cuts");
71  push_back("PF Kin Cuts");
72  push_back("PF Delta Phi");
73  push_back("PF Jet ID");
74 
75  set("Calo Cuts", useCalo);
76  set("Calo Kin Cuts", useCalo);
77  set("Calo Delta Phi", useCalo);
78  set("Calo Jet ID", useCalo);
79  set("PF Cuts", !useCalo);
80  set("PF Kin Cuts", !useCalo);
81  set("PF Delta Phi", !useCalo);
82  set("PF Jet ID", !useCalo);
83 
84  // Indices for fast caching
85  caloCuts_ = index_type(&bits_, std::string("Calo Cuts"));
86  caloKin_ = index_type(&bits_, std::string("Calo Kin Cuts"));
87  caloDeltaPhi_ = index_type(&bits_, std::string("Calo Delta Phi"));
88  caloJetID_ = index_type(&bits_, std::string("Calo Jet ID"));
89  pfCuts_ = index_type(&bits_, std::string("PF Cuts"));
90  pfKin_ = index_type(&bits_, std::string("PF Kin Cuts"));
91  pfDeltaPhi_ = index_type(&bits_, std::string("PF Delta Phi"));
92  pfJetID_ = index_type(&bits_, std::string("PF Jet ID"));
93  }
94 
95  ~JetIDStudiesSelector() override {}
96 
97  bool operator()(edm::EventBase const& event, pat::strbitset& ret) override {
98  pat::strbitset retCaloJet = jetSel_->getBitTemplate();
99  pat::strbitset retPFJet = pfJetSel_->getBitTemplate();
100 
101  if (considerCut(caloCuts_)) {
102  passCut(ret, caloCuts_);
103  event.getByLabel(jetSrc_, h_jets_);
104  // Calo Cuts
105  if (h_jets_->size() >= 2 || ignoreCut(caloKin_)) {
106  passCut(ret, caloKin_);
107  pat::Jet const& jet0 = h_jets_->at(0);
108  pat::Jet const& jet1 = h_jets_->at(1);
109  double dphi = fabs(deltaPhi<double>(jet0.phi(), jet1.phi()));
110 
111  if (fabs(dphi - TMath::Pi()) < 1.0 || ignoreCut(caloDeltaPhi_)) {
112  passCut(ret, caloDeltaPhi_);
113 
114  retCaloJet.set(false);
115  bool pass0 = (*jetSel_)(jet0, retCaloJet);
116  retCaloJet.set(false);
117  bool pass1 = (*jetSel_)(jet1, retCaloJet);
118  if ((pass0 && pass1) || ignoreCut(caloJetID_)) {
119  passCut(ret, caloJetID_);
120  caloJet0_ = edm::Ptr<pat::Jet>(h_jets_, 0);
121  caloJet1_ = edm::Ptr<pat::Jet>(h_jets_, 1);
122 
123  return (bool)ret;
124  } // end if found 2 "loose" jet ID jets
125  } // end if delta phi
126  } // end calo kin cuts
127  } // end if calo cuts
128 
129  if (considerCut(pfCuts_)) {
130  passCut(ret, pfCuts_);
131  event.getByLabel(pfJetSrc_, h_pfjets_);
132  // PF Cuts
133  if (h_pfjets_->size() >= 2 || ignoreCut(pfKin_)) {
134  passCut(ret, pfKin_);
135  pat::Jet const& jet0 = h_pfjets_->at(0);
136  pat::Jet const& jet1 = h_pfjets_->at(1);
137  double dphi = fabs(deltaPhi<double>(jet0.phi(), jet1.phi()));
138 
139  if (fabs(dphi - TMath::Pi()) < 1.0 || ignoreCut(pfDeltaPhi_)) {
140  passCut(ret, pfDeltaPhi_);
141 
142  retPFJet.set(false);
143  bool pass0 = (*pfJetSel_)(jet0, retPFJet);
144  retPFJet.set(false);
145  bool pass1 = (*pfJetSel_)(jet1, retPFJet);
146  if ((pass0 && pass1) || ignoreCut(pfJetID_)) {
147  passCut(ret, pfJetID_);
148  pfJet0_ = edm::Ptr<pat::Jet>(h_pfjets_, 0);
149  pfJet1_ = edm::Ptr<pat::Jet>(h_pfjets_, 1);
150 
151  return (bool)ret;
152  } // end if found 2 "loose" jet ID jets
153  } // end if delta phi
154  } // end pf kin cuts
155  } // end if pf cuts
156 
157  setIgnored(ret);
158 
159  return false;
160  } // end of method
161 
162  std::shared_ptr<JetIDSelectionFunctor> const& jetSel() const { return jetSel_; }
163  std::shared_ptr<PFJetIDSelectionFunctor> const& pfJetSel() const { return pfJetSel_; }
164 
165  vector<pat::Jet> const& allCaloJets() const { return *h_jets_; }
166  vector<pat::Jet> const& allPFJets() const { return *h_pfjets_; }
167 
168  pat::Jet const& caloJet0() const { return *caloJet0_; }
169  pat::Jet const& caloJet1() const { return *caloJet1_; }
170 
171  pat::Jet const& pfJet0() const { return *pfJet0_; }
172  pat::Jet const& pfJet1() const { return *pfJet1_; }
173 
174  // Fast caching indices
175  index_type const& caloCuts() const { return caloCuts_; }
176  index_type const& caloKin() const { return caloKin_; }
177  index_type const& caloDeltaPhi() const { return caloDeltaPhi_; }
178  index_type const& caloJetID() const { return caloJetID_; }
179  index_type const& pfCuts() const { return pfCuts_; }
180  index_type const& pfKin() const { return pfKin_; }
181  index_type const& pfDeltaPhi() const { return pfDeltaPhi_; }
182  index_type const& pfJetID() const { return pfJetID_; }
183 
184 protected:
185  std::shared_ptr<JetIDSelectionFunctor> jetSel_;
186  std::shared_ptr<PFJetIDSelectionFunctor> pfJetSel_;
189 
192 
195 
198 
199  // Fast caching indices
200  index_type caloCuts_;
201  index_type caloKin_;
202  index_type caloDeltaPhi_;
203  index_type caloJetID_;
204  index_type pfCuts_;
205  index_type pfKin_;
206  index_type pfDeltaPhi_;
207  index_type pfJetID_;
208 };
209 
211 // ///////////////////// //
212 // // Main Subroutine // //
213 // ///////////////////// //
215 
216 int main(int argc, char* argv[]) {
217  if (argc < 2) {
218  std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl;
219  return 0;
220  }
221 
222  cout << "Hello from " << argv[0] << "!" << endl;
223 
224  // load framework libraries
225  gSystem->Load("libFWCoreFWLite");
227 
228  cout << "Getting parameters" << endl;
229  // Get the python configuration
230  ProcessDescImpl builder(argv[1]);
231  std::shared_ptr<edm::ProcessDesc> b = builder.processDesc();
232  std::shared_ptr<edm::ParameterSet> parameters = b->getProcessPSet();
233  parameters->registerIt();
234 
235  edm::ParameterSet const& jetStudiesParams = parameters->getParameter<edm::ParameterSet>("jetStudies");
236  edm::ParameterSet const& pfJetStudiesParams = parameters->getParameter<edm::ParameterSet>("pfJetStudies");
237  edm::ParameterSet const& caloJetIDParameters = parameters->getParameter<edm::ParameterSet>("jetIDSelector");
238  edm::ParameterSet const& pfJetIDParameters = parameters->getParameter<edm::ParameterSet>("pfJetIDSelector");
239  edm::ParameterSet const& plotParameters = parameters->getParameter<edm::ParameterSet>("plotParameters");
240  edm::ParameterSet const& inputs = parameters->getParameter<edm::ParameterSet>("inputs");
241  edm::ParameterSet const& outputs = parameters->getParameter<edm::ParameterSet>("outputs");
242 
243  cout << "Making RunLumiSelector" << endl;
244  RunLumiSelector runLumiSel(inputs);
245 
246  cout << "setting up TFileService" << endl;
247  // book a set of histograms
248  fwlite::TFileService fs = fwlite::TFileService(outputs.getParameter<std::string>("outputName"));
249  TFileDirectory theDir = fs.mkdir("histos");
250 
251  cout << "Setting up chain event" << endl;
252  // This object 'event' is used both to get all information from the
253  // event as well as to store histograms, etc.
254  fwlite::ChainEvent ev(inputs.getParameter<std::vector<std::string> >("fileNames"));
255 
256  cout << "Booking histograms" << endl;
257  // Book histograms
258 
259  std::map<std::string, TH1*> hists;
260 
261  hists["hist_nJet"] = theDir.make<TH1D>("hist_nJet", "Number of Calo Jets", 10, 0, 10);
262  hists["hist_nPFJet"] = theDir.make<TH1D>("hist_nPFJet", "Number of PF Jets", 10, 0, 10);
263 
264  hists["hist_jetPt"] = theDir.make<TH1D>("hist_jetPt", "Jet p_{T}", 400, 0, 400);
265  hists["hist_jetEtaVsPhi"] = theDir.make<TH2D>(
266  "hist_jetEtaVsPhi", "Jet #phi versus #eta;#eta;#phi", 50, -5.0, 5.0, 50, -TMath::Pi(), TMath::Pi());
267  hists["hist_jetNTracks"] = theDir.make<TH1D>("hist_jetNTracks", "Jet N_{TRACKS}", 20, 0, 20);
268  hists["hist_jetNTracksVsPt"] = theDir.make<TH2D>(
269  "hist_jetNTracksVsPt", "Number of Tracks versus Jet p_{T};Jet p_{T}(GeV/c) ;N_{Tracks}", 20, 0, 200, 20, 0, 20);
270  hists["hist_jetEMF"] = theDir.make<TH1D>("hist_jetEMF", "Jet EMF", 200, 0, 1);
271  hists["hist_jetPdgID"] = theDir.make<TH1D>("hist_jetPdgID", "PDG Id of Jet Constituents", 10000, 0, 10000);
272  hists["hist_jetGenEmE"] = theDir.make<TH1D>("hist_jetGenEmE", "Gen Jet EM Energy", 200, 0, 200);
273  hists["hist_jetGenHadE"] = theDir.make<TH1D>("hist_jetGenHadE", "Gen Jet HAD Energy", 200, 0, 200);
274  hists["hist_jetGenEMF"] = theDir.make<TH1D>("hist_jetGenEMF", "Gen Jet EMF", 200, 0, 1);
275  hists["hist_jetEoverGenE"] =
276  theDir.make<TH1D>("hist_jetEoverGenE", "Energy of reco Jet / Energy of gen Jet", 200, 0, 2.0);
277  hists["hist_jetCorr"] = theDir.make<TH1D>("hist_jetCorr", "Jet Correction Factor", 200, 0, 1.0);
278  hists["hist_n90Hits"] = theDir.make<TH1D>("hist_n90Hits", "Jet n90Hits", 20, 0, 20);
279  hists["hist_fHPD"] = theDir.make<TH1D>("hist_fHPD", "Jet fHPD", 200, 0, 1);
280  hists["hist_nConstituents"] = theDir.make<TH1D>("hist_nConstituents", "Jet nConstituents", 20, 0, 20);
281  hists["hist_jetCHF"] = theDir.make<TH1D>("hist_jetCHF", "Jet Charged Hadron Fraction", 200, 0, 1.0);
282 
283  hists["hist_good_jetPt"] = theDir.make<TH1D>("hist_good_jetPt", "Jet p_{T}", 400, 0, 400);
284  hists["hist_good_jetEtaVsPhi"] = theDir.make<TH2D>(
285  "hist_good_jetEtaVsPhi", "Jet #phi versus #eta;#eta;#phi", 50, -5.0, 5.0, 50, -TMath::Pi(), TMath::Pi());
286  hists["hist_good_jetNTracks"] = theDir.make<TH1D>("hist_good_jetNTracks", "Jet N_{TRACKS}", 20, 0, 20);
287  hists["hist_good_jetNTracksVsPt"] =
288  theDir.make<TH2D>("hist_good_jetNTracksVsPt",
289  "Number of Tracks versus Jet p_{T};Jet p_{T}(GeV/c) ;N_{Tracks}",
290  20,
291  0,
292  200,
293  20,
294  0,
295  20);
296  hists["hist_good_jetEMF"] = theDir.make<TH1D>("hist_good_jetEMF", "Jet EMF", 200, 0, 1);
297  hists["hist_good_jetPdgID"] = theDir.make<TH1D>("hist_good_jetPdgID", "PDG Id of Jet Constituents", 10000, 0, 10000);
298  hists["hist_good_jetGenEmE"] = theDir.make<TH1D>("hist_good_jetGenEmE", "Gen Jet EM Energy", 200, 0, 200);
299  hists["hist_good_jetGenHadE"] = theDir.make<TH1D>("hist_good_jetGenHadE", "Gen Jet HAD Energy", 200, 0, 200);
300  hists["hist_good_jetGenEMF"] = theDir.make<TH1D>("hist_good_jetGenEMF", "Gen Jet EMF", 200, 0, 1);
301  hists["hist_good_jetEoverGenE"] =
302  theDir.make<TH1D>("hist_good_jetEoverGenE", "Energy of reco Jet / Energy of gen Jet", 200, 0, 2.0);
303  hists["hist_good_jetCorr"] = theDir.make<TH1D>("hist_good_jetCorr", "Jet Correction Factor", 200, 0, 1.0);
304  hists["hist_good_n90Hits"] = theDir.make<TH1D>("hist_good_n90Hits", "Jet n90Hits", 20, 0, 20);
305  hists["hist_good_fHPD"] = theDir.make<TH1D>("hist_good_fHPD", "Jet fHPD", 200, 0, 1);
306  hists["hist_good_nConstituents"] = theDir.make<TH1D>("hist_good_nConstituents", "Jet nConstituents", 20, 0, 20);
307  hists["hist_good_jetCHF"] = theDir.make<TH1D>("hist_good_jetCHF", "Jet Charged Hadron Fraction", 200, 0, 1.0);
308 
309  hists["hist_pf_jetPt"] = theDir.make<TH1D>("hist_pf_jetPt", "PFJet p_{T}", 400, 0, 400);
310  hists["hist_pf_jetEtaVsPhi"] = theDir.make<TH2D>(
311  "hist_pf_jetEtaVsPhi", "PFJet #phi versus #eta;#eta;#phi", 50, -5.0, 5.0, 50, -TMath::Pi(), TMath::Pi());
312  hists["hist_pf_jetNTracks"] = theDir.make<TH1D>("hist_pf_jetNTracks", "PFJet N_{TRACKS}", 20, 0, 20);
313  hists["hist_pf_jetNTracksVsPt"] = theDir.make<TH2D>(
314  "hist_pf_jetNTracksVsPt", "Number of Tracks versus Jet p_{T};Jet p_{T}(GeV/c) ;N_{Tracks}", 20, 0, 200, 20, 0, 20);
315  hists["hist_pf_jetEMF"] = theDir.make<TH1D>("hist_pf_jetEMF", "PFJet EMF", 200, 0, 1);
316  hists["hist_pf_jetCHF"] = theDir.make<TH1D>("hist_pf_jetCHF", "PFJet CHF", 200, 0, 1);
317  hists["hist_pf_jetNHF"] = theDir.make<TH1D>("hist_pf_jetNHF", "PFJet NHF", 200, 0, 1);
318  hists["hist_pf_jetCEF"] = theDir.make<TH1D>("hist_pf_jetCEF", "PFJet CEF", 200, 0, 1);
319  hists["hist_pf_jetNEF"] = theDir.make<TH1D>("hist_pf_jetNEF", "PFJet NEF", 200, 0, 1);
320  hists["hist_pf_jetPdgID"] = theDir.make<TH1D>("hist_pf_jetPdgID", "PDG Id of Jet Constituents", 10000, 0, 10000);
321  hists["hist_pf_jetGenEmE"] = theDir.make<TH1D>("hist_pf_jetGenEmE", "Gen Jet EM Energy", 200, 0, 200);
322  hists["hist_pf_jetGenHadE"] = theDir.make<TH1D>("hist_pf_jetGenHadE", "Gen Jet HAD Energy", 200, 0, 200);
323  hists["hist_pf_jetGenEMF"] = theDir.make<TH1D>("hist_pf_jetGenEMF", "Gen Jet EMF", 200, 0, 1);
324  hists["hist_pf_jetEoverGenE"] =
325  theDir.make<TH1D>("hist_pf_jetEoverGenE", "Energy of reco Jet / Energy of gen Jet", 200, 0, 2.0);
326  hists["hist_pf_jetCorr"] = theDir.make<TH1D>("hist_pf_jetCorr", "PFJet Correction Factor", 200, 0, 1.0);
327  hists["hist_pf_nConstituents"] = theDir.make<TH1D>("hist_pf_nConstituents", "PFJet nConstituents", 20, 0, 20);
328 
329  hists["hist_pf_good_jetPt"] = theDir.make<TH1D>("hist_pf_good_jetPt", "PFJet p_{T}", 400, 0, 400);
330  hists["hist_pf_good_jetEtaVsPhi"] = theDir.make<TH2D>(
331  "hist_pf_good_jetEtaVsPhi", "PFJet #phi versus #eta;#eta;#phi", 50, -5.0, 5.0, 50, -TMath::Pi(), TMath::Pi());
332  hists["hist_pf_good_jetNTracks"] = theDir.make<TH1D>("hist_pf_good_jetNTracks", "PFJet N_{TRACKS}", 20, 0, 20);
333  hists["hist_pf_good_jetNTracksVsPt"] =
334  theDir.make<TH2D>("hist_pf_good_jetNTracksVsPt",
335  "Number of Tracks versus Jet p_{T};Jet p_{T}(GeV/c) ;N_{Tracks}",
336  20,
337  0,
338  200,
339  20,
340  0,
341  20);
342  hists["hist_pf_good_jetEMF"] = theDir.make<TH1D>("hist_pf_good_jetEMF", "PFJet EMF", 200, 0, 1);
343  hists["hist_pf_good_jetCHF"] = theDir.make<TH1D>("hist_pf_good_jetCHF", "PFJet CHF", 200, 0, 1);
344  hists["hist_pf_good_jetNHF"] = theDir.make<TH1D>("hist_pf_good_jetNHF", "PFJet NHF", 200, 0, 1);
345  hists["hist_pf_good_jetCEF"] = theDir.make<TH1D>("hist_pf_good_jetCEF", "PFJet CEF", 200, 0, 1);
346  hists["hist_pf_good_jetNEF"] = theDir.make<TH1D>("hist_pf_good_jetNEF", "PFJet NEF", 200, 0, 1);
347  hists["hist_pf_good_jetPdgID"] =
348  theDir.make<TH1D>("hist_pf_good_jetPdgID", "PDG Id of Jet Constituents", 10000, 0, 10000);
349  hists["hist_pf_good_jetGenEmE"] = theDir.make<TH1D>("hist_pf_good_jetGenEmE", "Gen Jet EM Energy", 200, 0, 200);
350  hists["hist_pf_good_jetGenHadE"] = theDir.make<TH1D>("hist_pf_good_jetGenHadE", "Gen Jet HAD Energy", 200, 0, 200);
351  hists["hist_pf_good_jetGenEMF"] = theDir.make<TH1D>("hist_pf_good_jetGenEMF", "Gen Jet EMF", 200, 0, 1);
352  hists["hist_pf_good_jetEoverGenE"] =
353  theDir.make<TH1D>("hist_pf_good_jetEoverGenE", "Energy of reco Jet / Energy of gen Jet", 200, 0, 2.0);
354  hists["hist_pf_good_jetCorr"] = theDir.make<TH1D>("hist_pf_good_jetCorr", "PFJet Correction Factor", 200, 0, 1.0);
355  hists["hist_pf_good_nConstituents"] =
356  theDir.make<TH1D>("hist_pf_good_nConstituents", "PFJet nConstituents", 20, 0, 20);
357 
358  hists["hist_mjj"] = theDir.make<TH1D>("hist_mjj", "Dijet mass", 300, 0, 300);
359  hists["hist_dR_jj"] = theDir.make<TH1D>("hist_dR_jj", "#Delta R_{JJ}", 200, 0, TMath::TwoPi());
360  hists["hist_imbalance_jj"] = theDir.make<TH1D>("hist_imbalance_jj", "Dijet imbalance", 200, -1.0, 1.0);
361 
362  hists["hist_pf_mjj"] = theDir.make<TH1D>("hist_pf_mjj", "Dijet mass", 300, 0, 300);
363  hists["hist_pf_dR_jj"] = theDir.make<TH1D>("hist_pf_dR_jj", "#Delta R_{JJ}", 200, 0, TMath::TwoPi());
364  hists["hist_pf_imbalance_jj"] = theDir.make<TH1D>("hist_pf_imbalance_jj", "Dijet imbalance", 200, -1.0, 1.0);
365 
366  cout << "Making functors" << endl;
367  JetIDStudiesSelector caloSelector(caloJetIDParameters, pfJetIDParameters, jetStudiesParams);
368 
369  JetIDStudiesSelector pfSelector(caloJetIDParameters, pfJetIDParameters, pfJetStudiesParams);
370 
371  bool doTracks = plotParameters.getParameter<bool>("doTracks");
372  bool useMC = plotParameters.getParameter<bool>("useMC");
373 
374  cout << "About to begin looping" << endl;
375 
376  int nev = 0;
377  //loop through each event
378  for (ev.toBegin(); !ev.atEnd(); ++ev, ++nev) {
379  edm::EventBase const& event = ev;
380 
381  if (runLumiSel(ev) == false)
382  continue;
383 
384  if (nev % 100 == 0)
385  cout << "Processing run " << event.id().run() << ", lumi " << event.id().luminosityBlock() << ", event "
386  << event.id().event() << endl;
387 
388  pat::strbitset retCalo = caloSelector.getBitTemplate();
389  caloSelector(event, retCalo);
390 
391  pat::strbitset retPF = pfSelector.getBitTemplate();
392  pfSelector(event, retPF);
393 
397  if (retCalo.test(caloSelector.caloKin())) {
398  if (retCalo.test(caloSelector.caloDeltaPhi())) {
399  vector<pat::Jet> const& allCaloJets = caloSelector.allCaloJets();
400 
401  for (std::vector<pat::Jet>::const_iterator jetBegin = allCaloJets.begin(),
402  jetEnd = jetBegin + 2,
403  ijet = jetBegin;
404  ijet != jetEnd;
405  ++ijet) {
406  const pat::Jet& jet = *ijet;
407 
408  double pt = jet.pt();
409 
410  const reco::TrackRefVector& jetTracks = jet.associatedTracks();
411 
412  hists["hist_jetPt"]->Fill(pt);
413  hists["hist_jetEtaVsPhi"]->Fill(jet.eta(), jet.phi());
414  hists["hist_jetNTracks"]->Fill(jetTracks.size());
415  hists["hist_jetNTracksVsPt"]->Fill(pt, jetTracks.size());
416  hists["hist_jetEMF"]->Fill(jet.emEnergyFraction());
417  hists["hist_jetCorr"]->Fill(jet.jecFactor("Uncorrected"));
418  hists["hist_n90Hits"]->Fill(static_cast<int>(jet.jetID().n90Hits));
419  hists["hist_fHPD"]->Fill(jet.jetID().fHPD);
420  hists["hist_nConstituents"]->Fill(jet.nConstituents());
421 
422  if (useMC && jet.genJet() != nullptr) {
423  hists["hist_jetGenEmE"]->Fill(jet.genJet()->emEnergy());
424  hists["hist_jetGenHadE"]->Fill(jet.genJet()->hadEnergy());
425  hists["hist_jetEoverGenE"]->Fill(jet.energy() / jet.genJet()->energy());
426  hists["hist_jetGenEMF"]->Fill(jet.genJet()->emEnergy() / jet.genJet()->energy());
427  }
428  if (doTracks) {
429  TLorentzVector p4_tracks(0, 0, 0, 0);
430  for (reco::TrackRefVector::const_iterator itrk = jetTracks.begin(), itrkEnd = jetTracks.end();
431  itrk != itrkEnd;
432  ++itrk) {
433  TLorentzVector p4_trk;
434  double M_PION = 0.140;
435  p4_trk.SetPtEtaPhiM((*itrk)->pt(), (*itrk)->eta(), (*itrk)->phi(), M_PION);
436  p4_tracks += p4_trk;
437  }
438  hists["hist_jetCHF"]->Fill(p4_tracks.Energy() / jet.energy());
439  }
440 
441  } // end loop over jets
442 
443  if (retCalo.test(caloSelector.caloJetID())) {
444  pat::Jet const& jet0 = caloSelector.caloJet0();
445  pat::Jet const& jet1 = caloSelector.caloJet1();
446 
447  TLorentzVector p4_j0(jet0.px(), jet0.py(), jet0.pz(), jet0.energy());
448  TLorentzVector p4_j1(jet1.px(), jet1.py(), jet1.pz(), jet1.energy());
449 
450  TLorentzVector p4_jj = p4_j0 + p4_j1;
451 
452  hists["hist_mjj"]->Fill(p4_jj.M());
453  hists["hist_dR_jj"]->Fill(p4_j0.DeltaR(p4_j1));
454  hists["hist_imbalance_jj"]->Fill((p4_j0.Perp() - p4_j1.Perp()) / (p4_j0.Perp() + p4_j1.Perp()));
455 
456  hists["hist_good_jetPt"]->Fill(jet0.pt());
457  hists["hist_good_jetEtaVsPhi"]->Fill(jet0.eta(), jet0.phi());
458  hists["hist_good_jetNTracks"]->Fill(jet0.associatedTracks().size());
459  hists["hist_good_jetNTracksVsPt"]->Fill(jet0.pt(), jet0.associatedTracks().size());
460  hists["hist_good_jetEMF"]->Fill(jet0.emEnergyFraction());
461  hists["hist_good_jetCorr"]->Fill(jet0.jecFactor("Uncorrected"));
462  hists["hist_good_n90Hits"]->Fill(static_cast<int>(jet0.jetID().n90Hits));
463  hists["hist_good_fHPD"]->Fill(jet0.jetID().fHPD);
464  hists["hist_good_nConstituents"]->Fill(jet0.nConstituents());
465 
466  hists["hist_good_jetPt"]->Fill(jet1.pt());
467  hists["hist_good_jetEtaVsPhi"]->Fill(jet1.eta(), jet1.phi());
468  hists["hist_good_jetNTracks"]->Fill(jet1.associatedTracks().size());
469  hists["hist_good_jetNTracksVsPt"]->Fill(jet1.pt(), jet1.associatedTracks().size());
470  hists["hist_good_jetEMF"]->Fill(jet1.emEnergyFraction());
471  hists["hist_good_jetCorr"]->Fill(jet1.jecFactor("Uncorrected"));
472  hists["hist_good_n90Hits"]->Fill(static_cast<int>(jet1.jetID().n90Hits));
473  hists["hist_good_fHPD"]->Fill(jet1.jetID().fHPD);
474  hists["hist_good_nConstituents"]->Fill(jet1.nConstituents());
475 
476  } // end if passed calo jet id
477  } // end if passed dphi cuts
478  } // end if passed kin cuts
479 
483  if (retPF.test(pfSelector.pfDeltaPhi())) {
484  vector<pat::Jet> const& allPFJets = pfSelector.allPFJets();
485 
486  for (std::vector<pat::Jet>::const_iterator jetBegin = allPFJets.begin(), jetEnd = jetBegin + 2, ijet = jetBegin;
487  ijet != jetEnd;
488  ++ijet) {
489  const pat::Jet& jet = *ijet;
490 
491  double pt = jet.pt();
492 
493  hists["hist_pf_jetPt"]->Fill(pt);
494  hists["hist_pf_jetEtaVsPhi"]->Fill(jet.eta(), jet.phi());
495  hists["hist_pf_nConstituents"]->Fill(jet.nConstituents());
496  hists["hist_pf_jetCEF"]->Fill(jet.chargedEmEnergyFraction());
497  hists["hist_pf_jetNEF"]->Fill(jet.neutralEmEnergyFraction());
498  hists["hist_pf_jetCHF"]->Fill(jet.chargedHadronEnergyFraction());
499  hists["hist_pf_jetNHF"]->Fill(jet.neutralHadronEnergyFraction());
500 
501  if (useMC && jet.genJet() != nullptr) {
502  hists["hist_pf_jetGenEmE"]->Fill(jet.genJet()->emEnergy());
503  hists["hist_pf_jetGenHadE"]->Fill(jet.genJet()->hadEnergy());
504  hists["hist_pf_jetEoverGenE"]->Fill(jet.energy() / jet.genJet()->energy());
505 
506  hists["hist_pf_jetGenEMF"]->Fill(jet.genJet()->emEnergy() / jet.genJet()->energy());
507  }
508 
509  } // end loop over jets
510 
511  if (retPF.test(pfSelector.pfJetID())) {
512  pat::Jet const& jet0 = pfSelector.pfJet0();
513  pat::Jet const& jet1 = pfSelector.pfJet1();
514 
515  TLorentzVector p4_j0(jet0.px(), jet0.py(), jet0.pz(), jet0.energy());
516  TLorentzVector p4_j1(jet1.px(), jet1.py(), jet1.pz(), jet1.energy());
517 
518  TLorentzVector p4_jj = p4_j0 + p4_j1;
519 
520  hists["hist_pf_mjj"]->Fill(p4_jj.M());
521  hists["hist_pf_dR_jj"]->Fill(p4_j0.DeltaR(p4_j1));
522  hists["hist_pf_imbalance_jj"]->Fill((p4_j0.Perp() - p4_j1.Perp()) / (p4_j0.Perp() + p4_j1.Perp()));
523 
524  hists["hist_pf_good_jetPt"]->Fill(jet0.pt());
525  hists["hist_pf_good_jetEtaVsPhi"]->Fill(jet0.eta(), jet0.phi());
526  hists["hist_pf_good_nConstituents"]->Fill(jet0.nConstituents());
527  hists["hist_pf_good_jetCEF"]->Fill(jet0.chargedEmEnergyFraction());
528  hists["hist_pf_good_jetNEF"]->Fill(jet0.neutralEmEnergyFraction());
529  hists["hist_pf_good_jetCHF"]->Fill(jet0.chargedHadronEnergyFraction());
530  hists["hist_pf_good_jetNHF"]->Fill(jet0.neutralHadronEnergyFraction());
531 
532  hists["hist_pf_good_jetPt"]->Fill(jet1.pt());
533  hists["hist_pf_good_jetEtaVsPhi"]->Fill(jet1.eta(), jet1.phi());
534  hists["hist_pf_good_nConstituents"]->Fill(jet1.nConstituents());
535  hists["hist_pf_good_jetCEF"]->Fill(jet1.chargedEmEnergyFraction());
536  hists["hist_pf_good_jetNEF"]->Fill(jet1.neutralEmEnergyFraction());
537  hists["hist_pf_good_jetCHF"]->Fill(jet1.chargedHadronEnergyFraction());
538  hists["hist_pf_good_jetNHF"]->Fill(jet1.neutralHadronEnergyFraction());
539 
540  } // end if 2 good PF jets
541 
542  } // end if delta phi pf cuts
543 
544  } // end loop over events
545 
546  cout << "Calo jet selection" << endl;
547  caloSelector.print(std::cout);
548  cout << "PF jet selection" << endl;
549  pfSelector.print(std::cout);
550 
551  return 0;
552 }
JetIDStudiesSelector::pfJetID
index_type const & pfJetID() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:182
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:355
TFileService.h
TFileDirectory::mkdir
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileDirectory.cc:112
JetIDStudiesSelector::pfCuts_
index_type pfCuts_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:204
Handle.h
cmsBatch.argv
argv
Definition: cmsBatch.py:279
TwoPi
const double TwoPi
Definition: CosmicMuonParameters.h:19
JetIDStudiesSelector::caloJet0_
edm::Ptr< pat::Jet > caloJet0_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:193
PFJetIDSelectionFunctor.h
JetIDStudiesSelector::pfKin_
index_type pfKin_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:205
dir2webdir.argc
argc
Definition: dir2webdir.py:39
JetIDStudiesSelector::~JetIDStudiesSelector
~JetIDStudiesSelector() override
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:95
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
JetIDSelectionFunctor.h
TFileDirectory::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileDirectory.h:53
edm
HLT enums.
Definition: AlignableModifier.h:19
Muon.h
pat::strbitset::set
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:126
Photon.h
JetIDStudiesSelector::caloDeltaPhi_
index_type caloDeltaPhi_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:202
gather_cfg.cout
cout
Definition: gather_cfg.py:144
PatBasicFWLiteJetAnalyzer_Selector_cfg.outputs
outputs
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:48
pat::Jet::neutralHadronEnergyFraction
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:404
JetIDStudiesSelector::jetSrc_
edm::InputTag jetSrc_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:187
HLT_2018_cff.jetTracks
jetTracks
Definition: HLT_2018_cff.py:49366
JetIDStudiesSelector::pfDeltaPhi_
index_type pfDeltaPhi_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:206
PatBasicFWLiteJetAnalyzer_Selector_cfg.plotParameters
plotParameters
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:34
TFileDirectory
Definition: TFileDirectory.h:24
edm::RefVector< TrackCollection >
compare.hists
hists
Definition: compare.py:319
JetIDStudiesSelector::caloJetID
index_type const & caloJetID() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:178
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition: LeafCandidate.h:146
reco::JetID::fHPD
float fHPD
Definition: JetID.h:41
main
int main(int argc, char *argv[])
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:216
align_cfg.TFileService
TFileService
Definition: align_cfg.py:143
PatBasicFWLiteJetAnalyzer_Selector_cfg.useMC
useMC
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:36
edm::Handle
Definition: AssociativeIterator.h:50
Tau.h
JetIDStudiesSelector::caloKin_
index_type caloKin_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:201
PyBind11ProcessDesc::processDesc
std::unique_ptr< edm::ProcessDesc > processDesc() const
Definition: PyBind11ProcessDesc.cc:103
JetIDStudiesSelector::pfJetID_
index_type pfJetID_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:207
PFJetIDSelectionFunctor
PF Jet selector for pat::Jets.
Definition: PFJetIDSelectionFunctor.h:25
EventSelector
Selector< edm::EventBase > EventSelector
Definition: EventSelector.h:18
parameters
parameters
Definition: BeamSpot_PayloadInspector.cc:14
JetIDStudiesSelector::caloCuts
index_type const & caloCuts() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:175
pat::Jet
Analysis-level calorimeter jet class.
Definition: Jet.h:77
pat::Jet::associatedTracks
const reco::TrackRefVector & associatedTracks() const
method to return a vector of refs to the tracks associated to this jet
JetIDStudiesSelector::pfKin
index_type const & pfKin() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:180
JetIDStudiesSelector::pfDeltaPhi
index_type const & pfDeltaPhi() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:181
reco::LeafCandidate::py
double py() const final
y coordinate of momentum vector
Definition: LeafCandidate.h:142
ChainEvent.h
JetIDSelectionFunctor
Jet selector for pat::Jets and for CaloJets.
Definition: JetIDSelectionFunctor.h:28
JetIDStudiesSelector::JetIDStudiesSelector
JetIDStudiesSelector(edm::ParameterSet const &caloJetIdParams, edm::ParameterSet const &pfJetIdParams, edm::ParameterSet const &params)
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:57
JetIDStudiesSelector::jetSel_
std::shared_ptr< JetIDSelectionFunctor > jetSel_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:185
reco::Jet::nConstituents
virtual int nConstituents() const
Definition: Jet.h:65
JetIDStudiesSelector
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:55
JetIDStudiesSelector::caloCuts_
index_type caloCuts_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:200
ProcessDesc.h
JetIDStudiesSelector::caloKin
index_type const & caloKin() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:176
distMuonMETValueMapProducer_cff.useCalo
useCalo
Definition: distMuonMETValueMapProducer_cff.py:25
JetIDStudiesSelector::caloDeltaPhi
index_type const & caloDeltaPhi() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:177
b
double b
Definition: hdecay.h:118
RunLumiSelector.h
JetIDStudiesSelector::caloJetID_
index_type caloJetID_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:203
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
EventSelector.h
PyBind11ProcessDesc
Definition: PyBind11ProcessDesc.h:16
JetIDStudiesSelector::pfJetSel
std::shared_ptr< PFJetIDSelectionFunctor > const & pfJetSel() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:163
JetIDStudiesSelector::caloJet1_
edm::Ptr< pat::Jet > caloJet1_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:194
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
FWLiteEnabler::enable
static void enable()
enable automatic library loading
Definition: FWLiteEnabler.cc:46
edm::ParameterSet
Definition: ParameterSet.h:36
JetIDStudiesSelector::pfJet1_
edm::Ptr< pat::Jet > pfJet1_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:197
reco::LeafCandidate::eta
double eta() const final
momentum pseudorapidity
Definition: LeafCandidate.h:152
pat::Jet::neutralEmEnergyFraction
float neutralEmEnergyFraction() const
neutralEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:412
JetIDStudiesSelector::pfJet0
pat::Jet const & pfJet0() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:171
JetIDStudiesSelector::pfJet0_
edm::Ptr< pat::Jet > pfJet0_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:196
Ptr.h
PixelMapPlotter.inputs
inputs
Definition: PixelMapPlotter.py:490
PatBasicFWLiteJetAnalyzer_Selector_cfg.doTracks
doTracks
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:35
pat::strbitset::test
bool test(std::string s) const
test
Definition: strbitset.h:287
MET.h
Jet.h
JetIDStudiesSelector::operator()
bool operator()(edm::EventBase const &event, pat::strbitset &ret) override
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:97
InputTag.h
JetIDStudiesSelector::pfJetSrc_
edm::InputTag pfJetSrc_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:188
edm::Ptr< pat::Jet >
pat::Jet::chargedEmEnergyFraction
float chargedEmEnergyFraction() const
chargedEmEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:408
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
pat::strbitset
Definition: strbitset.h:23
RunLumiSelector
Definition: RunLumiSelector.h:17
JetIDStudiesSelector::allPFJets
vector< pat::Jet > const & allPFJets() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:166
JetIDStudiesSelector::h_pfjets_
edm::Handle< vector< pat::Jet > > h_pfjets_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:191
std
Definition: JetResolutionObject.h:76
pat::Jet::emEnergyFraction
float emEnergyFraction() const
returns the jet electromagnetic energy fraction
Definition: Jet.h:326
pat::Jet::chargedHadronEnergyFraction
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:400
JetIDStudiesSelector::caloJet1
pat::Jet const & caloJet1() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:169
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition: LeafCandidate.h:148
JetIDStudiesSelector::pfJetSel_
std::shared_ptr< PFJetIDSelectionFunctor > pfJetSel_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:186
strbitset.h
metsig::jet
Definition: SignAlgoResolutions.h:47
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
edm::EventBase
Definition: EventBase.h:46
edm::RefVectorIterator
Definition: EDProductfwd.h:33
ProcessDescImpl.h
pat::Jet::jetID
reco::JetID const & jetID() const
accessing Jet ID information
Definition: Jet.h:506
fwlite::TFileService
Definition: TFileService.h:12
reco::JetID::n90Hits
short n90Hits
Definition: JetID.h:43
reco::LeafCandidate::energy
double energy() const final
energy
Definition: LeafCandidate.h:125
Electron.h
JetIDStudiesSelector::allCaloJets
vector< pat::Jet > const & allCaloJets() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:165
JetIDStudiesSelector::jetSel
std::shared_ptr< JetIDSelectionFunctor > const & jetSel() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:162
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
JetIDStudiesSelector::pfJet1
pat::Jet const & pfJet1() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:172
fwlite::ChainEvent
Definition: ChainEvent.h:46
JetIDStudiesSelector::h_jets_
edm::Handle< vector< pat::Jet > > h_jets_
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:190
event
Definition: event.py:1
edm::RefVector::size
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
Selector.h
reco::LeafCandidate::px
double px() const final
x coordinate of momentum vector
Definition: LeafCandidate.h:140
reco::LeafCandidate::pz
double pz() const final
z coordinate of momentum vector
Definition: LeafCandidate.h:144
edm::InputTag
Definition: InputTag.h:15
FWLiteEnabler.h
pat::Jet::jecFactor
float jecFactor(const std::string &level, const std::string &flavor="none", const std::string &set="") const
JetIDStudiesSelector::caloJet0
pat::Jet const & caloJet0() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:168
JetIDStudiesSelector::pfCuts
index_type const & pfCuts() const
Definition: PatBasicFWLiteJetAnalyzer_Selector.cc:179