CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: SiPixelIsAliveCalibration.cc,v 1.20 2008/04/17 11:47:05 fblekman Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
26 
28 
33 
34 //
35 // class decleration
36 //
37 
39  public:
42 
43  void doSetup(const edm::ParameterSet &);
44  virtual bool doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix);
45 
46 
47  private:
48 
49  virtual void calibrationSetup(const edm::EventSetup&) ;
50  virtual void calibrationEnd() ;
51  virtual void newDetID(uint32_t detid);
52  virtual bool checkCorrectCalibrationType();
53  // ----------member data ---------------------------
54 
55  std::map<uint32_t,MonitorElement *> bookkeeper_;
56  std::map<uint32_t,MonitorElement *> summaries_;
57  double mineff_;
58 };
59 
60 //
61 // constants, enums and typedefs
62 //
63 
64 //
65 // static data member definitions
66 //
67 
68 //
69 // constructors and destructor
70 //
73  mineff_(iConfig.getUntrackedParameter<double>("minEfficiencyForAliveDef",0.8))
74 {
75  //now do what ever initialization is needed
76 
77 }
78 
79 
81 {
82 
83  // do anything here that needs to be done at desctruction time
84  // (e.g. close files, deallocate resources etc.)
85 
86 }
87 
88 
89 //
90 // member functions
91 //
92 void
94  setDQMDirectory(detid);
95  std::string tempname=translateDetIdToString(detid);
96  bookkeeper_[detid]= bookDQMHistoPlaquetteSummary2D(detid,"pixelAlive","pixel alive for "+tempname);
97  int xpix = bookkeeper_[detid]->getNbinsX();
98  int ypix = bookkeeper_[detid]->getNbinsY();
99  int tpix = xpix*ypix;
100  summaries_[detid]= bookDQMHistogram1D(detid,"pixelAliveSummary",bookkeeper_[detid]->getTitle(),calib_->getNTriggers()+1,-.05,.95+(1./(float)calib_->getNTriggers()));
101  summaries_[detid]->setBinContent(1, tpix);
102 }
103 bool
105  if(calibrationMode_=="PixelAlive")
106  return true;
107  else if(calibrationMode_=="unknown"){
108  edm::LogInfo("SiPixelIsAliveCalibration") << "calibration mode is: " << calibrationMode_ << ", continuing anyway..." ;
109  return true;
110  }
111  else{
112  // edm::LogError("SiPixelIsAliveCalibration")<< "unknown calibration mode for Pixel ALive, should be \"PixelAlive\" and is \"" << calibrationMode_ << "\"";
113  }
114  return false;
115 
116 }
118 
119 }
120 void
122 
123 }
124 bool
125 SiPixelIsAliveCalibration::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix){
126  edm::LogInfo("SiPixelIsAliveCalibration") << "now looking at DetID " << detid << ", pixel " << ipix->row() << "," << ipix->col() << std::endl;
127 
128  double denom=0;
129  double nom=0;
130  for(uint32_t i=0; i<ipix->getnpoints(); i++){
131  nom+=ipix->getnentries(i);
132  denom+=calib_->getNTriggers();
133  if(i>0)
134  edm::LogWarning("SiPixelIsAliveCalibration::doFits") << " ERROR!!" << " number of vcal points is now " << i << " for detid " << detid << std::endl;
135  }
136  edm::LogInfo("SiPixelIsAliveCalibration") << "DetID/col/row " << detid << "/"<< ipix->col() << "/" << ipix->row() << ", now calculating efficiency: " << nom << "/" << denom <<"=" << nom/denom << std::endl;
137  double eff = -1;
138  if(denom>0)
139  eff = nom/denom;
140  if(bookkeeper_[detid]->getBinContent(ipix->col()+1,ipix->row()+1)==0){
141  bookkeeper_[detid]->Fill(ipix->col(), ipix->row(), eff);
142  summaries_[detid]->Fill(eff);
143  float zerobin = summaries_[detid]->getBinContent(1);
144  summaries_[detid]->setBinContent(1, zerobin-1);
145  }
146  else
147  bookkeeper_[detid]->setBinContent(ipix->col()+1,ipix->row()+1,-2);
148  return true;
149 }
150 
151 void
153  // print summary of bad modules:
154  for(std::map<uint32_t, MonitorElement *>::const_iterator idet= bookkeeper_.begin(); idet!=bookkeeper_.end();++idet){
155  float idead = 0;
156  float iunderthres=0;
157  float imultiplefill=0;
158  float itot=0;
159  uint32_t detid=idet->first;
160 
161  setDQMDirectory(detid);
162  for(int icol=1; icol <= bookkeeper_[detid]->getNbinsX(); ++icol){
163  for(int irow=1; irow <= bookkeeper_[detid]->getNbinsY(); ++irow){
164  itot++;
165  double val = bookkeeper_[detid]->getBinContent(icol,irow);
166  if(val==0)
167  idead++;
168  if(val<mineff_)
169  iunderthres++;
170  if(val==-2)
171  imultiplefill++;
172 
173  }
174  }
175  edm::LogInfo("SiPixelIsAliveCalibration") << "summary for " << translateDetIdToString(detid) << "\tfrac dead:" << idead/itot << " frac below " << mineff_ << ":" << iunderthres/itot << " bad " << imultiplefill/itot << std::endl;
176  }
177 
178 }
179 // ---------- method called to for each event ------------
180 
181 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
SiPixelIsAliveCalibration(const edm::ParameterSet &)
edm::ESHandle< SiPixelCalibConfiguration > calib_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void newDetID(uint32_t detid)
std::string translateDetIdToString(uint32_t detid)
MonitorElement * bookDQMHistoPlaquetteSummary2D(uint32_t detid, std::string name, std::string title)
std::map< uint32_t, MonitorElement * > bookkeeper_
virtual bool doFits(uint32_t detid, std::vector< SiPixelCalibDigi >::const_iterator ipix)
MonitorElement * bookDQMHistogram1D(uint32_t detid, std::string name, std::string title, int nchX, double lowX, double highX)
std::map< uint32_t, MonitorElement * > summaries_
virtual void calibrationSetup(const edm::EventSetup &)
void doSetup(const edm::ParameterSet &)