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