CMS 3D CMS Logo

SiPixelIsAliveCalibration.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelIsAliveCalibrationAnalyzer
4 // Class: SiPixelIsAliveCalibrationAnalyzer
5 //
14 //
15 // Original Author: Freya Blekman
16 // Created: Mon Dec 3 14:07:42 CET 2007
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
25 
27 
32 
33 //
34 // class decleration
35 //
36 
38 public:
40  ~SiPixelIsAliveCalibration() override;
41 
42  void doSetup(const edm::ParameterSet &);
43  bool doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix) override;
44 
45 private:
46  void calibrationSetup(const edm::EventSetup &) override;
47  void calibrationEnd() override;
48  void newDetID(uint32_t detid) override;
49  bool checkCorrectCalibrationType() override;
50  // ----------member data ---------------------------
51 
52  std::map<uint32_t, MonitorElement *> bookkeeper_;
53  std::map<uint32_t, MonitorElement *> summaries_;
54  double mineff_;
55 };
56 
57 //
58 // constants, enums and typedefs
59 //
60 
61 //
62 // static data member definitions
63 //
64 
65 //
66 // constructors and destructor
67 //
70  mineff_(iConfig.getUntrackedParameter<double>("minEfficiencyForAliveDef", 0.8)) {
71  // now do what ever initialization is needed
72 }
73 
75  // do anything here that needs to be done at desctruction time
76  // (e.g. close files, deallocate resources etc.)
77 }
78 
79 //
80 // member functions
81 //
83  setDQMDirectory(detid);
84  std::string tempname = translateDetIdToString(detid);
85  bookkeeper_[detid] = bookDQMHistoPlaquetteSummary2D(detid, "pixelAlive", "pixel alive for " + tempname);
86  int xpix = bookkeeper_[detid]->getNbinsX();
87  int ypix = bookkeeper_[detid]->getNbinsY();
88  int tpix = xpix * ypix;
89  summaries_[detid] = bookDQMHistogram1D(detid,
90  "pixelAliveSummary",
91  bookkeeper_[detid]->getTitle(),
92  calib_->getNTriggers() + 1,
93  -.05,
94  .95 + (1. / (float)calib_->getNTriggers()));
95  summaries_[detid]->setBinContent(1, tpix);
96 }
98  if (calibrationMode_ == "PixelAlive")
99  return true;
100  else if (calibrationMode_ == "unknown") {
101  edm::LogInfo("SiPixelIsAliveCalibration")
102  << "calibration mode is: " << calibrationMode_ << ", continuing anyway...";
103  return true;
104  } else {
105  // edm::LogError("SiPixelIsAliveCalibration")<< "unknown calibration mode
106  // for Pixel ALive, should be \"PixelAlive\" and is \"" <<
107  // calibrationMode_ << "\"";
108  }
109  return false;
110 }
113 bool SiPixelIsAliveCalibration::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix) {
114  edm::LogInfo("SiPixelIsAliveCalibration")
115  << "now looking at DetID " << detid << ", pixel " << ipix->row() << "," << ipix->col() << std::endl;
116 
117  double denom = 0;
118  double nom = 0;
119  for (uint32_t i = 0; i < ipix->getnpoints(); i++) {
120  nom += ipix->getnentries(i);
121  denom += calib_->getNTriggers();
122  if (i > 0)
123  edm::LogWarning("SiPixelIsAliveCalibration::doFits")
124  << " ERROR!!"
125  << " number of vcal points is now " << i << " for detid " << detid << std::endl;
126  }
127  edm::LogInfo("SiPixelIsAliveCalibration")
128  << "DetID/col/row " << detid << "/" << ipix->col() << "/" << ipix->row()
129  << ", now calculating efficiency: " << nom << "/" << denom << "=" << nom / denom << std::endl;
130  double eff = -1;
131  if (denom > 0)
132  eff = nom / denom;
133  if (bookkeeper_[detid]->getBinContent(ipix->col() + 1, ipix->row() + 1) == 0) {
134  bookkeeper_[detid]->Fill(ipix->col(), ipix->row(), eff);
135  summaries_[detid]->Fill(eff);
136  float zerobin = summaries_[detid]->getBinContent(1);
137  summaries_[detid]->setBinContent(1, zerobin - 1);
138  } else
139  bookkeeper_[detid]->setBinContent(ipix->col() + 1, ipix->row() + 1, -2);
140  return true;
141 }
142 
144  // print summary of bad modules:
145  for (std::map<uint32_t, MonitorElement *>::const_iterator idet = bookkeeper_.begin(); idet != bookkeeper_.end();
146  ++idet) {
147  float idead = 0;
148  float iunderthres = 0;
149  float imultiplefill = 0;
150  float itot = 0;
151  uint32_t detid = idet->first;
152 
153  setDQMDirectory(detid);
154  for (int icol = 1; icol <= bookkeeper_[detid]->getNbinsX(); ++icol) {
155  for (int irow = 1; irow <= bookkeeper_[detid]->getNbinsY(); ++irow) {
156  itot++;
157  double val = bookkeeper_[detid]->getBinContent(icol, irow);
158  if (val == 0)
159  idead++;
160  if (val < mineff_)
161  iunderthres++;
162  if (val == -2)
163  imultiplefill++;
164  }
165  }
166  edm::LogInfo("SiPixelIsAliveCalibration")
167  << "summary for " << translateDetIdToString(detid) << "\tfrac dead:" << idead / itot << " frac below "
168  << mineff_ << ":" << iunderthres / itot << " bad " << imultiplefill / itot << std::endl;
169  }
170 }
171 // ---------- method called to for each event ------------
172 
173 // define this as a plug-in
SiPixelIsAliveCalibration(const edm::ParameterSet &)
edm::ESHandle< SiPixelCalibConfiguration > calib_
std::string translateDetIdToString(uint32_t detid)
std::map< uint32_t, MonitorElement * > bookkeeper_
MonitorElement * bookDQMHistoPlaquetteSummary2D(uint32_t detid, std::string name, std::string title)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * bookDQMHistogram1D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX)
bool doFits(uint32_t detid, std::vector< SiPixelCalibDigi >::const_iterator ipix) override
void doSetup(const edm::ParameterSet &)
void newDetID(uint32_t detid) override
void calibrationSetup(const edm::EventSetup &) override
std::map< uint32_t, MonitorElement * > summaries_