Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:58

0001 #include "DQMOffline/Trigger/plugins/METDQM.h"
0002 
0003 METDQM::METDQM() = default;
0004 
0005 METDQM::~METDQM() = default;
0006 
0007 void METDQM::initialise(const edm::ParameterSet& iConfig) {
0008   met_variable_binning_ =
0009       iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("metBinning");
0010   met_binning_ =
0011       getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("metPSet"));
0012   phi_binning_ =
0013       getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("phiPSet"));
0014   ls_binning_ =
0015       getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lsPSet"));
0016 }
0017 
0018 void METDQM::bookHistograms(DQMStore::IBooker& ibooker) {
0019   std::string histname, histtitle;
0020 
0021   histname = "met";
0022   histtitle = "PFMET";
0023   bookME(ibooker, metME_, histname, histtitle, met_binning_.nbins, met_binning_.xmin, met_binning_.xmax);
0024   setMETitle(metME_, "PF MET [GeV]", "events / [GeV]");
0025 
0026   histname = "met_variable";
0027   histtitle = "PFMET";
0028   bookME(ibooker, metME_variableBinning_, histname, histtitle, met_variable_binning_);
0029   setMETitle(metME_variableBinning_, "PF MET [GeV]", "events / [GeV]");
0030 
0031   histname = "metVsLS";
0032   histtitle = "PFMET vs LS";
0033   bookME(ibooker,
0034          metVsLS_,
0035          histname,
0036          histtitle,
0037          ls_binning_.nbins,
0038          ls_binning_.xmin,
0039          ls_binning_.xmax,
0040          met_binning_.xmin,
0041          met_binning_.xmax);
0042   setMETitle(metVsLS_, "LS", "PF MET [GeV]");
0043 
0044   histname = "metPhi";
0045   histtitle = "PFMET phi";
0046   bookME(ibooker, metPhiME_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0047   setMETitle(metPhiME_, "PF MET #phi", "events / 0.1 rad");
0048 }
0049 
0050 void METDQM::fillHistograms(const double& met, const double& phi, const int& ls, const bool passCond) {
0051   // filling histograms (denominator)
0052   metME_.denominator->Fill(met);
0053   metME_variableBinning_.denominator->Fill(met);
0054   metPhiME_.denominator->Fill(phi);
0055 
0056   metVsLS_.denominator->Fill(ls, met);
0057 
0058   // applying selection for numerator
0059   if (passCond) {
0060     // filling histograms (num_genTriggerEventFlag_)
0061     metME_.numerator->Fill(met);
0062     metME_variableBinning_.numerator->Fill(met);
0063     metPhiME_.numerator->Fill(phi);
0064     metVsLS_.numerator->Fill(ls, met);
0065   }
0066 }
0067 
0068 void METDQM::fillMetDescription(edm::ParameterSetDescription& histoPSet) {
0069   edm::ParameterSetDescription metPSet;
0070   fillHistoPSetDescription(metPSet);
0071   histoPSet.add<edm::ParameterSetDescription>("metPSet", metPSet);
0072 
0073   std::vector<double> bins = {0.,   20.,  40.,  60.,  80.,  90.,  100., 110., 120., 130., 140., 150., 160.,
0074                               170., 180., 190., 200., 220., 240., 260., 280., 300., 350., 400., 450., 1000.};
0075   histoPSet.add<std::vector<double> >("metBinning", bins);
0076 
0077   edm::ParameterSetDescription phiPSet;
0078   fillHistoPSetDescription(phiPSet);
0079   histoPSet.add<edm::ParameterSetDescription>("phiPSet", phiPSet);
0080 
0081   edm::ParameterSetDescription lsPSet;
0082   fillHistoLSPSetDescription(lsPSet);
0083   histoPSet.add<edm::ParameterSetDescription>("lsPSet", lsPSet);
0084 }