File indexing completed on 2024-04-06 12:08:27
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include "FWCore/ServiceRegistry/interface/Service.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/EventSetup.h"
0027
0028
0029 #include "DQMServices/Core/interface/DQMStore.h"
0030 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0031
0032 #include <algorithm>
0033
0034 class SiPixelPhase1EfficiencyExtras : public DQMEDHarvester {
0035 public:
0036 explicit SiPixelPhase1EfficiencyExtras(const edm::ParameterSet& conf);
0037 ~SiPixelPhase1EfficiencyExtras() override;
0038
0039 protected:
0040 void beginRun(edm::Run const& run, edm::EventSetup const& eSetup) override;
0041
0042 void dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) override;
0043
0044 const std::string effFolderName_;
0045 const std::string vtxFolderName_;
0046 const std::string instLumiFolderName_;
0047 };
0048
0049 SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras(const edm::ParameterSet& iConfig)
0050 : effFolderName_(iConfig.getParameter<std::string>("EffFolderName")),
0051 vtxFolderName_(iConfig.getParameter<std::string>("VtxFolderName")),
0052 instLumiFolderName_(iConfig.getParameter<std::string>("InstLumiFolderName")) {
0053 edm::LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::SiPixelPhase1EfficiencyExtras: Hello!";
0054 }
0055
0056 SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras() {
0057 edm::LogInfo("PixelDQM") << "SiPixelPhase1EfficiencyExtras::~SiPixelPhase1EfficiencyExtras: Destructor";
0058 }
0059
0060 void SiPixelPhase1EfficiencyExtras::beginRun(edm::Run const& run, edm::EventSetup const& eSetup) {}
0061
0062
0063
0064
0065 void SiPixelPhase1EfficiencyExtras::dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) {
0066 iBooker.setCurrentFolder(effFolderName_);
0067
0068
0069 MonitorElement* vtx_v_lumi = iGetter.get(vtxFolderName_ + "/NumberOfGoodPVtxVsLS_GenTk");
0070
0071 MonitorElement* scalLumi_v_lumi = iGetter.get(instLumiFolderName_ + "/lumiVsLS");
0072
0073 MonitorElement* eff_v_lumi_forward =
0074 iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXDisk_PXForward");
0075
0076 MonitorElement* eff_v_lumi_barrel =
0077 iGetter.get(effFolderName_ + "/hitefficiency_per_Lumisection_per_PXLayer_PXBarrel");
0078
0079
0080 bool createNvtx = true;
0081 bool createInstLumi = true;
0082
0083
0084 if (!eff_v_lumi_forward) {
0085 edm::LogWarning("SiPixelPhase1EfficiencyExtras")
0086 << "no hitefficiency_per_Lumisection_per_PXDisk_PXForward ME is available in " << effFolderName_;
0087 return;
0088 }
0089 if (!eff_v_lumi_barrel) {
0090 edm::LogWarning("SiPixelPhase1EfficiencyExtras")
0091 << "no hitefficiency_per_Lumisection_per_PXLayer_PXBarrel ME is available in " << effFolderName_;
0092 return;
0093 }
0094 if (!vtx_v_lumi) {
0095 edm::LogWarning("SiPixelPhase1EfficiencyExtras")
0096 << "no NumberOfGoodPVtxVsLS_GenTK ME is available in " << vtxFolderName_;
0097 createNvtx = false;
0098 }
0099 if (!scalLumi_v_lumi) {
0100 edm::LogWarning("SiPixelPhase1EfficiencyExtras") << "no lumiVsLS ME is available in " << instLumiFolderName_;
0101 createInstLumi = false;
0102 }
0103
0104
0105 if (vtx_v_lumi && vtx_v_lumi->getEntries() == 0)
0106 createNvtx = false;
0107 if (scalLumi_v_lumi && scalLumi_v_lumi->getEntries() == 0)
0108 createInstLumi = false;
0109
0110
0111 if (createInstLumi and scalLumi_v_lumi->getTProfile()->GetMaximum() <= 0.)
0112 createInstLumi = false;
0113
0114 double eff = 0.0;
0115
0116
0117 if (createNvtx) {
0118
0119 MonitorElement* eff_v_vtx_barrel =
0120 iBooker.book2D("hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel",
0121 "hitefficiency_per_meanNvtx_per_PXLayer_PXBarrel; meanNvtx; PXLayer",
0122 500,
0123 0,
0124 100,
0125 3,
0126 .5,
0127 3.5);
0128
0129 MonitorElement* eff_v_vtx_forward =
0130 iBooker.book2D("hitefficiency_per_meanNvtx_per_PXDisk_PXForward",
0131 "hitefficiency_per_meanNvtx_per_PXDisk_PXForward; meanNvtx; PXDisk",
0132 500,
0133 0,
0134 100,
0135 7,
0136 -3.5,
0137 3.5);
0138
0139
0140 int numLumiNvtx = int(vtx_v_lumi->getNbinsX());
0141 double nvtx = 0.0;
0142 int binNumVtx = 0;
0143
0144
0145 for (int iLumi = 1; iLumi < numLumiNvtx - 1; iLumi++) {
0146
0147 nvtx = vtx_v_lumi->getBinContent(iLumi);
0148
0149
0150 if (nvtx != 0) {
0151
0152 binNumVtx = eff_v_vtx_barrel->getTH2F()->FindBin(nvtx);
0153
0154
0155 for (int iLayer = 1; iLayer < 8; iLayer++) {
0156
0157 eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
0158
0159
0160 eff_v_vtx_forward->setBinContent(binNumVtx, iLayer, eff);
0161 }
0162
0163
0164 for (int iLayer = 1; iLayer < 5; iLayer++) {
0165
0166 eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
0167
0168
0169 eff_v_vtx_barrel->setBinContent(binNumVtx, iLayer, eff);
0170 }
0171 }
0172 }
0173 }
0174
0175 if (createInstLumi) {
0176
0177 int yMax2 = std::max(1., scalLumi_v_lumi->getTProfile()->GetMaximum());
0178 yMax2 *= 1.1;
0179
0180
0181 MonitorElement* eff_v_scalLumi_barrel =
0182 iBooker.book2D("hitefficiency_per_scalLumi_per_PXLayer_PXBarrel",
0183 "hitefficiency_per_scalLumi_per_PXLayer_PXBarrel; scal inst lumi E30; PXLayer",
0184 500,
0185 0,
0186 yMax2,
0187 3,
0188 .5,
0189 3.5);
0190
0191 MonitorElement* eff_v_scalLumi_forward =
0192 iBooker.book2D("hitefficiency_per_scalLumi_per_PXDisk_PXForward",
0193 "hitefficiency_per_scalLumi_per_PXDisk_PXForward; scal inst lumi E30; PXDisk",
0194 500,
0195 0,
0196 yMax2,
0197 7,
0198 -3.5,
0199 3.5);
0200
0201
0202 int numLumiScal = int(scalLumi_v_lumi->getNbinsX());
0203 double scalLumi = 0.0;
0204 int binNumScal = 0;
0205
0206
0207 for (int iLumi = 1; iLumi < numLumiScal - 1; iLumi++) {
0208
0209 scalLumi = scalLumi_v_lumi->getBinContent(iLumi);
0210
0211
0212 if (scalLumi > 0) {
0213
0214 binNumScal = eff_v_scalLumi_barrel->getTH2F()->FindBin(scalLumi);
0215
0216
0217 for (int iLayer = 1; iLayer < 8; iLayer++) {
0218
0219 eff = eff_v_lumi_forward->getBinContent(iLumi - 1, iLayer);
0220
0221
0222 eff_v_scalLumi_forward->setBinContent(binNumScal, iLayer, eff);
0223 }
0224
0225
0226 for (int iLayer = 1; iLayer < 5; iLayer++) {
0227
0228 eff = eff_v_lumi_barrel->getBinContent(iLumi - 1, iLayer);
0229
0230
0231 eff_v_scalLumi_barrel->setBinContent(binNumScal, iLayer, eff);
0232 }
0233 }
0234 }
0235 } else
0236 return;
0237 }
0238
0239
0240 DEFINE_FWK_MODULE(SiPixelPhase1EfficiencyExtras);