CMS 3D CMS Logo

JetDQM.cc
Go to the documentation of this file.
2 //#include "DataFormats/Math/interface/deltaR.h"
3 
4 JetDQM::JetDQM() = default;
5 
6 JetDQM::~JetDQM() = default;
7 
8 void JetDQM::initialise(const edm::ParameterSet& iConfig) {
10  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetptBinning");
12  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jet1ptBinning");
14  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jet2ptBinning");
16  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("mjjBinning");
18  getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("jetetaPSet"));
20  getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("detajjPSet"));
22  getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("dphijjPSet"));
24  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("mindphijmetPSet"));
25  ls_binning_ =
26  getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("jetlsPSet"));
27 }
28 
30  std::string histname, histtitle;
31 
32  histname = "jet1pt";
33  histtitle = "PFJet1 pT";
34  bookME(ibooker, jet1ptME_, histname, histtitle, jet1pt_variable_binning_);
35  setMETitle(jet1ptME_, "PFJet1 p_{T} [GeV]", "events / [GeV]");
36 
37  histname = "jet2pt";
38  histtitle = "PFJet2 pT";
39  bookME(ibooker, jet2ptME_, histname, histtitle, jet2pt_variable_binning_);
40  setMETitle(jet2ptME_, "PFJet2 p_{T} [GeV]", "events / [GeV]");
41 
42  histname = "jet1eta";
43  histtitle = "PFJet1 eta";
45  setMETitle(jet1etaME_, "PFJet1 #eta", "events / [rad]");
46 
47  histname = "jet2eta";
48  histtitle = "PFJet2 eta";
50  setMETitle(jet2etaME_, "PFJet2 #eta", "events / [rad]");
51 
52  histname = "cjetpt";
53  histtitle = "central PFJet pT";
54  bookME(ibooker, cjetptME_, histname, histtitle, jetpt_variable_binning_);
55  setMETitle(cjetptME_, "central PFJet p_{T} [GeV]", "events / [GeV]");
56 
57  histname = "fjetpt";
58  histtitle = "forward PFJet pT";
59  bookME(ibooker, fjetptME_, histname, histtitle, jetpt_variable_binning_);
60  setMETitle(fjetptME_, "forward PFJet p_{T} [GeV]", "events / [GeV]");
61 
62  histname = "cjeteta";
63  histtitle = "central PFJet eta";
65  setMETitle(cjetetaME_, "central PFJet #eta", "events / [rad]");
66 
67  histname = "fjeteta";
68  histtitle = "forward PFJet eta";
70  setMETitle(fjetetaME_, "forward PFJet #eta", "events / [rad]");
71 
72  histname = "mjj";
73  histtitle = "PFDiJet M";
74  bookME(ibooker, mjjME_, histname, histtitle, mjj_variable_binning_);
75  setMETitle(mjjME_, "PFDiJet M [GeV]", "events / [GeV]");
76 
77  histname = "detajj";
78  histtitle = "PFDiJet DeltaEta";
80  setMETitle(detajjME_, "PFDiJet #Delta#eta", "events / [rad]");
81 
82  histname = "dphijj";
83  histtitle = "PFDiJet DeltaPhi";
85  setMETitle(dphijjME_, "PFDiJet #Delta#phi", "events / [rad]");
86 
87  histname = "mindphijmet";
88  histtitle = "minDeltaPhi(PFJets,MET)";
89  bookME(ibooker,
91  histname,
92  histtitle,
96  setMETitle(mindphijmetME_, "min#Delta#phi(jets,MET)", "events / [rad]");
97 
98  histname = "jet1etaVsLS";
99  histtitle = "PFJet1 eta vs LS";
100  bookME(ibooker,
101  jet1etaVsLS_,
102  histname,
103  histtitle,
109  setMETitle(jet1etaVsLS_, "LS", "PF Jet1 #eta");
110 
111  histname = "mjjVsLS";
112  histtitle = "PFDiJet M vs LS";
113  bookME(ibooker, mjjVsLS_, histname, histtitle, ls_binning_.nbins, ls_binning_.xmin, ls_binning_.xmax, 0, 14000);
114  setMETitle(mjjVsLS_, "LS", "PFDiJet M [GeV]");
115 
116  histname = "mindphijmetVsLS";
117  histtitle = "minDeltaPhi(PFJets,MET) vs LS";
118  bookME(ibooker,
120  histname,
121  histtitle,
127  setMETitle(mindphijmetVsLS_, "LS", "min#Delta#phi(jets,MET)");
128 }
129 
130 void JetDQM::fillHistograms(const std::vector<reco::PFJet>& jets,
131  const reco::PFMET& pfmet,
132  const int ls,
133  const bool passCond) {
134  // filling histograms (denominator)
135  if (!jets.empty()) {
136  double eta1 = jets[0].eta();
140  if (jets.size() > 1) {
141  double eta2 = jets[1].eta();
144  if (fabs(eta1) < fabs(eta2)) {
149  } else {
154  }
155  double mass = (jets[0].p4() + jets[1].p4()).mass();
158  detajjME_.denominator->Fill(fabs(eta1 - eta2));
159  dphijjME_.denominator->Fill(fabs(reco::deltaPhi(jets[0].phi(), jets[1].phi())));
160 
161  double mindphi = fabs(reco::deltaPhi(jets[0].phi(), pfmet.phi()));
162  for (unsigned ij(0); ij < jets.size(); ++ij) {
163  if (ij > 4)
164  break;
165  double dphi = fabs(reco::deltaPhi(jets[ij].phi(), pfmet.phi()));
166  if (dphi < mindphi)
167  mindphi = dphi;
168  }
169 
170  mindphijmetME_.denominator->Fill(mindphi);
171  mindphijmetVsLS_.denominator->Fill(ls, mindphi);
172  }
173  } //at least 1 jet
174 
175  // applying selection for numerator
176  if (passCond) {
177  // filling histograms (num_genTriggerEventFlag_)
178  if (!jets.empty()) {
179  double eta1 = jets[0].eta();
180  jet1ptME_.numerator->Fill(jets[0].pt());
183  if (jets.size() > 1) {
184  double eta2 = jets[1].eta();
185  jet2ptME_.numerator->Fill(jets[1].pt());
187  if (fabs(eta1) < fabs(eta2)) {
190  cjetptME_.numerator->Fill(jets[0].pt());
191  fjetptME_.numerator->Fill(jets[1].pt());
192  } else {
195  cjetptME_.numerator->Fill(jets[1].pt());
196  fjetptME_.numerator->Fill(jets[0].pt());
197  }
198  double mass = (jets[0].p4() + jets[1].p4()).mass();
201  detajjME_.numerator->Fill(fabs(eta1 - eta2));
202  dphijjME_.numerator->Fill(fabs(reco::deltaPhi(jets[0].phi(), jets[1].phi())));
203 
204  double mindphi = fabs(reco::deltaPhi(jets[0].phi(), pfmet.phi()));
205  for (unsigned ij(0); ij < jets.size(); ++ij) {
206  if (ij > 4)
207  break;
208  double dphi = fabs(reco::deltaPhi(jets[ij].phi(), pfmet.phi()));
209  if (dphi < mindphi)
210  mindphi = dphi;
211  }
212 
213  mindphijmetME_.numerator->Fill(mindphi);
214  mindphijmetVsLS_.numerator->Fill(ls, mindphi);
215  }
216  } //at least 1 jet
217  }
218 }
219 
224 
228 
232 
235  histoPSet.add<edm::ParameterSetDescription>("mindphijmetPSet", mindphijmetPSet);
236 
237  std::vector<double> bins = {0., 20., 40., 60., 80., 100., 120., 140., 160., 180., 200., 220.,
238  240., 260., 280., 300., 350., 400., 450., 500, 750, 1000., 1500.};
239  histoPSet.add<std::vector<double> >("jetptBinning", bins);
240 
241  std::vector<double> pt1bins = {0., 20., 40., 60., 80., 90., 100., 110., 120., 130., 140., 150., 160.,
242  180., 210., 240., 270., 300., 330., 360., 400., 450., 500., 750., 1000., 1500.};
243  histoPSet.add<std::vector<double> >("jet1ptBinning", pt1bins);
244 
245  std::vector<double> pt2bins = {0., 20., 40., 45., 50., 55., 60., 65., 70., 80., 90., 100.,
246  110., 120., 150., 180., 210., 240., 270., 300., 350., 400., 1000.};
247  histoPSet.add<std::vector<double> >("jet2ptBinning", pt2bins);
248 
249  std::vector<double> mjjbins = {0., 200, 400, 600, 620, 640, 660, 680, 700, 720, 740, 760, 780, 800,
250  850, 900, 950, 1000, 1200, 1400, 1600, 1800, 2000, 2500, 3000, 4000, 6000};
251  histoPSet.add<std::vector<double> >("mjjBinning", mjjbins);
252 
255  histoPSet.add<edm::ParameterSetDescription>("jetlsPSet", lsPSet);
256 }
ObjME dphijjME_
Definition: JetDQM.h:50
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void fillHistograms(const std::vector< reco::PFJet > &jets, const reco::PFMET &pfmet, const int ls, const bool passCond)
Definition: JetDQM.cc:130
void initialise(const edm::ParameterSet &iConfig)
Definition: JetDQM.cc:8
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
~JetDQM() override
ObjME cjetptME_
Definition: JetDQM.h:44
static void fillHistoLSPSetDescription(edm::ParameterSetDescription &pset)
ObjME detajjME_
Definition: JetDQM.h:49
ObjME fjetetaME_
Definition: JetDQM.h:43
static void fillJetDescription(edm::ParameterSetDescription &histoPSet)
Definition: JetDQM.cc:220
static void fillHistoPSetDescription(edm::ParameterSetDescription &pset)
ObjME mindphijmetME_
Definition: JetDQM.h:53
ObjME jet1etaME_
Definition: JetDQM.h:38
MEbinning ls_binning_
Definition: JetDQM.h:33
void setMETitle(ObjME &me, const std::string &titleX, const std::string &titleY)
std::vector< double > jet1pt_variable_binning_
Definition: JetDQM.h:25
ObjME jet2etaME_
Definition: JetDQM.h:39
void Fill(long long x)
MEbinning mindphijmet_binning_
Definition: JetDQM.h:32
static MEbinning getHistoPSet(const edm::ParameterSet &pset)
ObjME jet2ptME_
Definition: JetDQM.h:37
void bookHistograms(DQMStore::IBooker &)
Definition: JetDQM.cc:29
ObjME jet1ptME_
Definition: JetDQM.h:36
MonitorElement * denominator
MonitorElement * numerator
std::vector< double > mjj_variable_binning_
Definition: JetDQM.h:27
ObjME mjjVsLS_
Definition: JetDQM.h:56
MEbinning jeteta_binning_
Definition: JetDQM.h:29
ObjME mindphijmetVsLS_
Definition: JetDQM.h:57
def ls(path, rec=False)
Definition: eostools.py:349
std::vector< double > jetpt_variable_binning_
Definition: JetDQM.h:24
void bookME(DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, const uint nbins, const double xmin, const double xmax, const bool bookDen=true)
ObjME cjetetaME_
Definition: JetDQM.h:42
std::vector< double > jet2pt_variable_binning_
Definition: JetDQM.h:26
MEbinning dphijj_binning_
Definition: JetDQM.h:31
MEbinning detajj_binning_
Definition: JetDQM.h:30
double phi() const final
momentum azimuthal angle
ObjME fjetptME_
Definition: JetDQM.h:45
ObjME jet1etaVsLS_
Definition: JetDQM.h:55
ObjME mjjME_
Definition: JetDQM.h:48