CMS 3D CMS Logo

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