File indexing completed on 2024-04-06 11:59:34
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <memory>
0022
0023
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025
0026 #include "CalibTracker/SiPixelTools/interface/SiPixelOfflineCalibAnalysisBase.h"
0027 #include "DQMServices/Core/interface/DQMStore.h"
0028 #include "DataFormats/SiPixelDigi/interface/SiPixelCalibDigi.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030
0031
0032
0033
0034
0035 class SiPixelIsAliveCalibration : public SiPixelOfflineCalibAnalysisBase {
0036 public:
0037 explicit SiPixelIsAliveCalibration(const edm::ParameterSet &);
0038 ~SiPixelIsAliveCalibration() override;
0039
0040 void doSetup(const edm::ParameterSet &);
0041 bool doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix) override;
0042
0043 private:
0044 void calibrationSetup(const edm::EventSetup &) override;
0045 void calibrationEnd() override;
0046 void newDetID(uint32_t detid) override;
0047 bool checkCorrectCalibrationType() override;
0048
0049
0050 std::map<uint32_t, MonitorElement *> bookkeeper_;
0051 std::map<uint32_t, MonitorElement *> summaries_;
0052 double mineff_;
0053 };
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066 SiPixelIsAliveCalibration::SiPixelIsAliveCalibration(const edm::ParameterSet &iConfig)
0067 : SiPixelOfflineCalibAnalysisBase(iConfig),
0068 mineff_(iConfig.getUntrackedParameter<double>("minEfficiencyForAliveDef", 0.8)) {
0069
0070 }
0071
0072 SiPixelIsAliveCalibration::~SiPixelIsAliveCalibration() {
0073
0074
0075 }
0076
0077
0078
0079
0080 void SiPixelIsAliveCalibration::newDetID(uint32_t detid) {
0081 setDQMDirectory(detid);
0082 std::string tempname = translateDetIdToString(detid);
0083 bookkeeper_[detid] = bookDQMHistoPlaquetteSummary2D(detid, "pixelAlive", "pixel alive for " + tempname);
0084 int xpix = bookkeeper_[detid]->getNbinsX();
0085 int ypix = bookkeeper_[detid]->getNbinsY();
0086 int tpix = xpix * ypix;
0087 summaries_[detid] = bookDQMHistogram1D(detid,
0088 "pixelAliveSummary",
0089 bookkeeper_[detid]->getTitle(),
0090 calib_->getNTriggers() + 1,
0091 -.05,
0092 .95 + (1. / (float)calib_->getNTriggers()));
0093 summaries_[detid]->setBinContent(1, tpix);
0094 }
0095 bool SiPixelIsAliveCalibration::checkCorrectCalibrationType() {
0096 if (calibrationMode_ == "PixelAlive")
0097 return true;
0098 else if (calibrationMode_ == "unknown") {
0099 edm::LogInfo("SiPixelIsAliveCalibration")
0100 << "calibration mode is: " << calibrationMode_ << ", continuing anyway...";
0101 return true;
0102 } else {
0103
0104
0105
0106 }
0107 return false;
0108 }
0109 void SiPixelIsAliveCalibration::doSetup(const edm::ParameterSet &iConf) {}
0110 void SiPixelIsAliveCalibration::calibrationSetup(const edm::EventSetup &iSetup) {}
0111 bool SiPixelIsAliveCalibration::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix) {
0112 edm::LogInfo("SiPixelIsAliveCalibration")
0113 << "now looking at DetID " << detid << ", pixel " << ipix->row() << "," << ipix->col() << std::endl;
0114
0115 double denom = 0;
0116 double nom = 0;
0117 for (uint32_t i = 0; i < ipix->getnpoints(); i++) {
0118 nom += ipix->getnentries(i);
0119 denom += calib_->getNTriggers();
0120 if (i > 0)
0121 edm::LogWarning("SiPixelIsAliveCalibration::doFits")
0122 << " ERROR!!"
0123 << " number of vcal points is now " << i << " for detid " << detid << std::endl;
0124 }
0125 edm::LogInfo("SiPixelIsAliveCalibration")
0126 << "DetID/col/row " << detid << "/" << ipix->col() << "/" << ipix->row()
0127 << ", now calculating efficiency: " << nom << "/" << denom << "=" << nom / denom << std::endl;
0128 double eff = -1;
0129 if (denom > 0)
0130 eff = nom / denom;
0131 if (bookkeeper_[detid]->getBinContent(ipix->col() + 1, ipix->row() + 1) == 0) {
0132 bookkeeper_[detid]->Fill(ipix->col(), ipix->row(), eff);
0133 summaries_[detid]->Fill(eff);
0134 float zerobin = summaries_[detid]->getBinContent(1);
0135 summaries_[detid]->setBinContent(1, zerobin - 1);
0136 } else
0137 bookkeeper_[detid]->setBinContent(ipix->col() + 1, ipix->row() + 1, -2);
0138 return true;
0139 }
0140
0141 void SiPixelIsAliveCalibration::calibrationEnd() {
0142
0143 for (std::map<uint32_t, MonitorElement *>::const_iterator idet = bookkeeper_.begin(); idet != bookkeeper_.end();
0144 ++idet) {
0145 float idead = 0;
0146 float iunderthres = 0;
0147 float imultiplefill = 0;
0148 float itot = 0;
0149 uint32_t detid = idet->first;
0150
0151 setDQMDirectory(detid);
0152 for (int icol = 1; icol <= bookkeeper_[detid]->getNbinsX(); ++icol) {
0153 for (int irow = 1; irow <= bookkeeper_[detid]->getNbinsY(); ++irow) {
0154 itot++;
0155 double val = bookkeeper_[detid]->getBinContent(icol, irow);
0156 if (val == 0)
0157 idead++;
0158 if (val < mineff_)
0159 iunderthres++;
0160 if (val == -2)
0161 imultiplefill++;
0162 }
0163 }
0164 edm::LogInfo("SiPixelIsAliveCalibration")
0165 << "summary for " << translateDetIdToString(detid) << "\tfrac dead:" << idead / itot << " frac below "
0166 << mineff_ << ":" << iunderthres / itot << " bad " << imultiplefill / itot << std::endl;
0167 }
0168 }
0169
0170
0171
0172 DEFINE_FWK_MODULE(SiPixelIsAliveCalibration);