CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
JetDQM Class Reference

#include <JetDQM.h>

Inheritance diagram for JetDQM:
TriggerDQMBase

Public Member Functions

void bookHistograms (DQMStore::IBooker &)
 
void fillHistograms (const std::vector< reco::PFJet > &jets, const reco::PFMET &pfmet, const int &ls, const bool passCond)
 
void initialise (const edm::ParameterSet &iConfig)
 
 JetDQM ()
 
 ~JetDQM () override
 
- Public Member Functions inherited from TriggerDQMBase
 TriggerDQMBase ()=default
 
virtual ~TriggerDQMBase ()=default
 

Static Public Member Functions

static void fillJetDescription (edm::ParameterSetDescription &histoPSet)
 
- Static Public Member Functions inherited from TriggerDQMBase
static void fillHistoLSPSetDescription (edm::ParameterSetDescription &pset)
 
static void fillHistoPSetDescription (edm::ParameterSetDescription &pset)
 
static MEbinning getHistoLSPSet (const edm::ParameterSet &pset)
 
static MEbinning getHistoPSet (const edm::ParameterSet &pset)
 

Private Attributes

ObjME cjetetaME_
 
ObjME cjetptME_
 
MEbinning detajj_binning_
 
ObjME detajjME_
 
MEbinning dphijj_binning_
 
ObjME dphijjME_
 
ObjME fjetetaME_
 
ObjME fjetptME_
 
ObjME jet1etaME_
 
ObjME jet1etaVsLS_
 
std::vector< double > jet1pt_variable_binning_
 
ObjME jet1ptME_
 
ObjME jet2etaME_
 
std::vector< double > jet2pt_variable_binning_
 
ObjME jet2ptME_
 
MEbinning jeteta_binning_
 
std::vector< double > jetpt_variable_binning_
 
MEbinning ls_binning_
 
MEbinning mindphijmet_binning_
 
ObjME mindphijmetME_
 
ObjME mindphijmetVsLS_
 
std::vector< double > mjj_variable_binning_
 
ObjME mjjME_
 
ObjME mjjVsLS_
 

Additional Inherited Members

- Protected Member Functions inherited from TriggerDQMBase
void bookME (DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, unsigned nbins, double xmin, double xmax)
 
void bookME (DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, const std::vector< double > &binningX)
 
void bookME (DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, unsigned nbinsX, double xmin, double xmax, double ymin, double ymax)
 
void bookME (DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, unsigned nbinsX, double xmin, double xmax, unsigned nbinsY, double ymin, double ymax)
 
void bookME (DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, const std::vector< double > &binningX, const std::vector< double > &binningY)
 
void setMETitle (ObjME &me, const std::string &titleX, const std::string &titleY)
 

Detailed Description

Definition at line 13 of file JetDQM.h.

Constructor & Destructor Documentation

JetDQM::JetDQM ( )
default
JetDQM::~JetDQM ( )
overridedefault

Member Function Documentation

void JetDQM::bookHistograms ( DQMStore::IBooker ibooker)

Definition at line 23 of file JetDQM.cc.

References TriggerDQMBase::bookME(), cjetetaME_, cjetptME_, detajj_binning_, detajjME_, dphijj_binning_, dphijjME_, fjetetaME_, fjetptME_, jet1etaME_, jet1etaVsLS_, jet1pt_variable_binning_, jet1ptME_, jet2etaME_, jet2pt_variable_binning_, jet2ptME_, jeteta_binning_, jetpt_variable_binning_, ls_binning_, mindphijmet_binning_, mindphijmetME_, mindphijmetVsLS_, mjj_variable_binning_, mjjME_, mjjVsLS_, TriggerDQMBase::MEbinning::nbins, TriggerDQMBase::setMETitle(), AlCaHLTBitMon_QueryRunRegistry::string, TriggerDQMBase::MEbinning::xmax, and TriggerDQMBase::MEbinning::xmin.

Referenced by ObjMonitor::bookHistograms().

24 {
25 
26  std::string histname, histtitle;
27 
28  histname = "jet1pt"; histtitle = "PFJet1 pT";
29  bookME(ibooker,jet1ptME_,histname,histtitle,jet1pt_variable_binning_);
30  setMETitle(jet1ptME_,"PFJet1 p_{T} [GeV]","events / [GeV]");
31 
32  histname = "jet2pt"; histtitle = "PFJet2 pT";
33  bookME(ibooker,jet2ptME_,histname,histtitle,jet2pt_variable_binning_);
34  setMETitle(jet2ptME_,"PFJet2 p_{T} [GeV]","events / [GeV]");
35 
36  histname = "jet1eta"; histtitle = "PFJet1 eta";
38  setMETitle(jet1etaME_,"PFJet1 #eta","events / [rad]");
39 
40  histname = "jet2eta"; histtitle = "PFJet2 eta";
42  setMETitle(jet2etaME_,"PFJet2 #eta","events / [rad]");
43 
44  histname = "cjetpt"; histtitle = "central PFJet pT";
45  bookME(ibooker,cjetptME_,histname,histtitle,jetpt_variable_binning_);
46  setMETitle(cjetptME_,"central PFJet p_{T} [GeV]","events / [GeV]");
47 
48  histname = "fjetpt"; histtitle = "forward PFJet pT";
49  bookME(ibooker,fjetptME_,histname,histtitle,jetpt_variable_binning_);
50  setMETitle(fjetptME_,"forward PFJet p_{T} [GeV]","events / [GeV]");
51 
52  histname = "cjeteta"; histtitle = "central PFJet eta";
54  setMETitle(cjetetaME_,"central PFJet #eta","events / [rad]");
55 
56  histname = "fjeteta"; histtitle = "forward PFJet eta";
58  setMETitle(fjetetaME_,"forward PFJet #eta","events / [rad]");
59 
60  histname = "mjj"; histtitle = "PFDiJet M";
61  bookME(ibooker,mjjME_,histname,histtitle,mjj_variable_binning_);
62  setMETitle(mjjME_,"PFDiJet M [GeV]","events / [GeV]");
63 
64  histname = "detajj"; histtitle = "PFDiJet DeltaEta";
66  setMETitle(detajjME_,"PFDiJet #Delta#eta","events / [rad]");
67 
68  histname = "dphijj"; histtitle = "PFDiJet DeltaPhi";
70  setMETitle(dphijjME_,"PFDiJet #Delta#phi","events / [rad]");
71 
72  histname = "mindphijmet"; histtitle = "minDeltaPhi(PFJets,MET)";
74  setMETitle(mindphijmetME_,"min#Delta#phi(jets,MET)","events / [rad]");
75 
76  histname = "jet1etaVsLS"; histtitle = "PFJet1 eta vs LS";
78  setMETitle(jet1etaVsLS_,"LS","PF Jet1 #eta");
79 
80  histname = "mjjVsLS"; histtitle = "PFDiJet M vs LS";
81  bookME(ibooker,mjjVsLS_,histname,histtitle,ls_binning_.nbins, ls_binning_.xmin, ls_binning_.xmax,0,14000);
82  setMETitle(mjjVsLS_,"LS","PFDiJet M [GeV]");
83 
84  histname = "mindphijmetVsLS"; histtitle = "minDeltaPhi(PFJets,MET) vs LS";
86  setMETitle(mindphijmetVsLS_,"LS","min#Delta#phi(jets,MET)");
87 
88 }
ObjME dphijjME_
Definition: JetDQM.h:53
ObjME cjetptME_
Definition: JetDQM.h:48
ObjME detajjME_
Definition: JetDQM.h:52
ObjME fjetetaME_
Definition: JetDQM.h:47
ObjME mindphijmetME_
Definition: JetDQM.h:55
ObjME jet1etaME_
Definition: JetDQM.h:43
MEbinning ls_binning_
Definition: JetDQM.h:38
void setMETitle(ObjME &me, const std::string &titleX, const std::string &titleY)
std::vector< double > jet1pt_variable_binning_
Definition: JetDQM.h:30
ObjME jet2etaME_
Definition: JetDQM.h:44
MEbinning mindphijmet_binning_
Definition: JetDQM.h:36
ObjME jet2ptME_
Definition: JetDQM.h:42
ObjME jet1ptME_
Definition: JetDQM.h:41
std::vector< double > mjj_variable_binning_
Definition: JetDQM.h:32
ObjME mjjVsLS_
Definition: JetDQM.h:58
MEbinning jeteta_binning_
Definition: JetDQM.h:33
ObjME mindphijmetVsLS_
Definition: JetDQM.h:59
std::vector< double > jetpt_variable_binning_
Definition: JetDQM.h:29
ObjME cjetetaME_
Definition: JetDQM.h:46
std::vector< double > jet2pt_variable_binning_
Definition: JetDQM.h:31
MEbinning dphijj_binning_
Definition: JetDQM.h:35
MEbinning detajj_binning_
Definition: JetDQM.h:34
ObjME fjetptME_
Definition: JetDQM.h:49
ObjME jet1etaVsLS_
Definition: JetDQM.h:57
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:51
void JetDQM::fillHistograms ( const std::vector< reco::PFJet > &  jets,
const reco::PFMET pfmet,
const int &  ls,
const bool  passCond 
)

Definition at line 90 of file JetDQM.cc.

References cjetetaME_, cjetptME_, reco::deltaPhi(), TriggerDQMBase::ObjME::denominator, detajjME_, dphijjME_, HcalObjRepresent::Fill(), fjetetaME_, fjetptME_, jet1etaME_, jet1etaVsLS_, jet1ptME_, jet2etaME_, jet2ptME_, ResonanceBuilder::mass, mindphijmetME_, mindphijmetVsLS_, mjjME_, mjjVsLS_, TriggerDQMBase::ObjME::numerator, phi, reco::LeafCandidate::phi(), and EnergyCorrector::pt.

Referenced by ObjMonitor::analyze().

93  {
94 
95  // filling histograms (denominator)
96  if (!jets.empty()){
97  double eta1 = jets[0].eta();
98  jet1ptME_.denominator -> Fill(jets[0].pt());
99  jet1etaME_.denominator -> Fill(eta1);
100  jet1etaVsLS_.denominator -> Fill(ls, eta1);
101  if (jets.size()>1){
102  double eta2 = jets[1].eta();
103  jet2ptME_.denominator -> Fill(jets[1].pt());
104  jet2etaME_.denominator -> Fill(eta2);
105  if (fabs(eta1)<fabs(eta2)){
106  cjetetaME_.denominator -> Fill(eta1);
107  fjetetaME_.denominator -> Fill(eta2);
108  cjetptME_.denominator -> Fill(jets[0].pt());
109  fjetptME_.denominator -> Fill(jets[1].pt());
110  } else {
111  cjetetaME_.denominator -> Fill(eta2);
112  fjetetaME_.denominator -> Fill(eta1);
113  cjetptME_.denominator -> Fill(jets[1].pt());
114  fjetptME_.denominator -> Fill(jets[0].pt());
115  }
116  double mass = (jets[0].p4()+jets[1].p4()).mass();
117  mjjME_.denominator -> Fill(mass);
118  mjjVsLS_.denominator -> Fill(ls, mass);
119  detajjME_.denominator -> Fill(fabs(eta1-eta2));
120  dphijjME_.denominator -> Fill(fabs(reco::deltaPhi(jets[0].phi(),jets[1].phi())));
121 
122  double mindphi = fabs(reco::deltaPhi(jets[0].phi(),pfmet.phi()));
123  for (unsigned ij(0); ij<jets.size(); ++ij){
124  if (ij>4) break;
125  double dphi = fabs(reco::deltaPhi(jets[ij].phi(),pfmet.phi()));
126  if (dphi<mindphi) mindphi=dphi;
127  }
128 
129  mindphijmetME_.denominator -> Fill(mindphi);
130  mindphijmetVsLS_.denominator -> Fill(ls, mindphi);
131  }
132  }//at least 1 jet
133 
134 
135  // applying selection for numerator
136  if (passCond){
137  // filling histograms (num_genTriggerEventFlag_)
138  if (!jets.empty()){
139  double eta1 = jets[0].eta();
140  jet1ptME_.numerator -> Fill(jets[0].pt());
141  jet1etaME_.numerator -> Fill(eta1);
142  jet1etaVsLS_.numerator -> Fill(ls, eta1);
143  if (jets.size()>1){
144  double eta2 = jets[1].eta();
145  jet2ptME_.numerator -> Fill(jets[1].pt());
146  jet2etaME_.numerator -> Fill(eta2);
147  if (fabs(eta1)<fabs(eta2)){
148  cjetetaME_.numerator -> Fill(eta1);
149  fjetetaME_.numerator -> Fill(eta2);
150  cjetptME_.numerator -> Fill(jets[0].pt());
151  fjetptME_.numerator -> Fill(jets[1].pt());
152  } else {
153  cjetetaME_.numerator -> Fill(eta2);
154  fjetetaME_.numerator -> Fill(eta1);
155  cjetptME_.numerator -> Fill(jets[1].pt());
156  fjetptME_.numerator -> Fill(jets[0].pt());
157  }
158  double mass = (jets[0].p4()+jets[1].p4()).mass();
159  mjjME_.numerator -> Fill(mass);
160  mjjVsLS_.numerator -> Fill(ls, mass);
161  detajjME_.numerator -> Fill(fabs(eta1-eta2));
162  dphijjME_.numerator -> Fill(fabs(reco::deltaPhi(jets[0].phi(),jets[1].phi())));
163 
164  double mindphi = fabs(reco::deltaPhi(jets[0].phi(),pfmet.phi()));
165  for (unsigned ij(0); ij<jets.size(); ++ij){
166  if (ij>4) break;
167  double dphi = fabs(reco::deltaPhi(jets[ij].phi(),pfmet.phi()));
168  if (dphi<mindphi) mindphi=dphi;
169  }
170 
171  mindphijmetME_.numerator -> Fill(mindphi);
172  mindphijmetVsLS_.numerator -> Fill(ls, mindphi);
173  }
174  }//at least 1 jet
175 
176  }
177 
178 }
ObjME dphijjME_
Definition: JetDQM.h:53
MonitorElement * numerator
ObjME cjetptME_
Definition: JetDQM.h:48
ObjME detajjME_
Definition: JetDQM.h:52
ObjME fjetetaME_
Definition: JetDQM.h:47
ObjME mindphijmetME_
Definition: JetDQM.h:55
ObjME jet1etaME_
Definition: JetDQM.h:43
ObjME jet2etaME_
Definition: JetDQM.h:44
ObjME jet2ptME_
Definition: JetDQM.h:42
MonitorElement * denominator
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
vector< PseudoJet > jets
ObjME jet1ptME_
Definition: JetDQM.h:41
ObjME mjjVsLS_
Definition: JetDQM.h:58
ObjME mindphijmetVsLS_
Definition: JetDQM.h:59
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
def ls(path, rec=False)
Definition: eostools.py:348
ObjME cjetetaME_
Definition: JetDQM.h:46
ObjME fjetptME_
Definition: JetDQM.h:49
double phi() const final
momentum azimuthal angle
ObjME jet1etaVsLS_
Definition: JetDQM.h:57
ObjME mjjME_
Definition: JetDQM.h:51
void JetDQM::fillJetDescription ( edm::ParameterSetDescription histoPSet)
static

Definition at line 180 of file JetDQM.cc.

References edm::ParameterSetDescription::add(), create_public_pileup_plots::bins, TriggerDQMBase::fillHistoLSPSetDescription(), and TriggerDQMBase::fillHistoPSetDescription().

Referenced by ObjMonitor::fillDescriptions().

180  {
181 
182  edm::ParameterSetDescription jetetaPSet;
183  fillHistoPSetDescription(jetetaPSet);
184  histoPSet.add<edm::ParameterSetDescription>("jetetaPSet", jetetaPSet);
185 
186  edm::ParameterSetDescription detajjPSet;
187  fillHistoPSetDescription(detajjPSet);
188  histoPSet.add<edm::ParameterSetDescription>("detajjPSet", detajjPSet);
189 
190  edm::ParameterSetDescription dphijjPSet;
191  fillHistoPSetDescription(dphijjPSet);
192  histoPSet.add<edm::ParameterSetDescription>("dphijjPSet", dphijjPSet);
193 
194  edm::ParameterSetDescription mindphijmetPSet;
195  fillHistoPSetDescription(mindphijmetPSet);
196  histoPSet.add<edm::ParameterSetDescription>("mindphijmetPSet", mindphijmetPSet);
197 
198  std::vector<double> bins = {0.,20.,40.,60.,80.,100.,120.,140.,160.,180.,200.,220.,240.,260.,280.,300.,350.,400.,450.,500,750,1000.,1500.};
199  histoPSet.add<std::vector<double> >("jetptBinning", bins);
200 
201  std::vector<double> pt1bins = {0.,20.,40.,60.,80.,90.,100.,110.,120.,130.,140.,150.,160.,180.,210.,240.,270.,300.,330.,360.,400.,450.,500.,750.,1000.,1500.};
202  histoPSet.add<std::vector<double> >("jet1ptBinning", pt1bins);
203 
204  std::vector<double> pt2bins = {0.,20.,40.,45.,50.,55.,60.,65.,70.,80.,90.,100.,110.,120.,150.,180.,210.,240.,270.,300.,350.,400.,1000.};
205  histoPSet.add<std::vector<double> >("jet2ptBinning", pt2bins);
206 
207  std::vector<double> mjjbins = {0.,200,400,600,620,640,660,680,700,720,740,760,780,800,850,900,950,1000,1200,1400,1600,1800,2000,2500,3000,4000,6000};
208  histoPSet.add<std::vector<double> >("mjjBinning", mjjbins);
209 
212  histoPSet.add<edm::ParameterSetDescription>("jetlsPSet", lsPSet);
213 
214 }
static void fillHistoLSPSetDescription(edm::ParameterSetDescription &pset)
static void fillHistoPSetDescription(edm::ParameterSetDescription &pset)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void JetDQM::initialise ( const edm::ParameterSet iConfig)

Definition at line 8 of file JetDQM.cc.

References detajj_binning_, dphijj_binning_, TriggerDQMBase::getHistoPSet(), edm::ParameterSet::getParameter(), jet1pt_variable_binning_, jet2pt_variable_binning_, jeteta_binning_, jetpt_variable_binning_, ls_binning_, mindphijmet_binning_, and mjj_variable_binning_.

Referenced by ObjMonitor::ObjMonitor().

8  {
9 
10  jetpt_variable_binning_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetptBinning");
11  jet1pt_variable_binning_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jet1ptBinning");
12  jet2pt_variable_binning_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jet2ptBinning");
13  mjj_variable_binning_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("mjjBinning");
14 
20 
21 }
T getParameter(std::string const &) const
MEbinning ls_binning_
Definition: JetDQM.h:38
std::vector< double > jet1pt_variable_binning_
Definition: JetDQM.h:30
MEbinning mindphijmet_binning_
Definition: JetDQM.h:36
static MEbinning getHistoPSet(const edm::ParameterSet &pset)
std::vector< double > mjj_variable_binning_
Definition: JetDQM.h:32
MEbinning jeteta_binning_
Definition: JetDQM.h:33
std::vector< double > jetpt_variable_binning_
Definition: JetDQM.h:29
std::vector< double > jet2pt_variable_binning_
Definition: JetDQM.h:31
MEbinning dphijj_binning_
Definition: JetDQM.h:35
MEbinning detajj_binning_
Definition: JetDQM.h:34

Member Data Documentation

ObjME JetDQM::cjetetaME_
private

Definition at line 46 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::cjetptME_
private

Definition at line 48 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

MEbinning JetDQM::detajj_binning_
private

Definition at line 34 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::detajjME_
private

Definition at line 52 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

MEbinning JetDQM::dphijj_binning_
private

Definition at line 35 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::dphijjME_
private

Definition at line 53 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::fjetetaME_
private

Definition at line 47 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::fjetptME_
private

Definition at line 49 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::jet1etaME_
private

Definition at line 43 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::jet1etaVsLS_
private

Definition at line 57 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

std::vector<double> JetDQM::jet1pt_variable_binning_
private

Definition at line 30 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::jet1ptME_
private

Definition at line 41 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::jet2etaME_
private

Definition at line 44 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

std::vector<double> JetDQM::jet2pt_variable_binning_
private

Definition at line 31 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::jet2ptME_
private

Definition at line 42 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

MEbinning JetDQM::jeteta_binning_
private

Definition at line 33 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

std::vector<double> JetDQM::jetpt_variable_binning_
private

Definition at line 29 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

MEbinning JetDQM::ls_binning_
private

Definition at line 38 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

MEbinning JetDQM::mindphijmet_binning_
private

Definition at line 36 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::mindphijmetME_
private

Definition at line 55 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::mindphijmetVsLS_
private

Definition at line 59 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

std::vector<double> JetDQM::mjj_variable_binning_
private

Definition at line 32 of file JetDQM.h.

Referenced by bookHistograms(), and initialise().

ObjME JetDQM::mjjME_
private

Definition at line 51 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().

ObjME JetDQM::mjjVsLS_
private

Definition at line 58 of file JetDQM.h.

Referenced by bookHistograms(), and fillHistograms().