CMS 3D CMS Logo

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