CMS 3D CMS Logo

CastorDigiMonitor.cc
Go to the documentation of this file.
5 
6 //****************************************************//
7 //********** CastorDigiMonitor: ******************//
8 //********** Author: Dmytro Volyanskyy *************//
9 //********** Date : 29.08.2008 (first version) ******//
12 //****************************************************//
13 //---- critical revision 26.06.2014 (Vladimir Popov)
14 // add rms check 15.04.2015 (Vladimir Popov)
15 //==================================================================//
16 
18 {
19  fVerbosity = ps.getUntrackedParameter<int>("debug",0);
20  if(fVerbosity) std::cout<<"CastorDigiMonitor Constructor: "<<this<<std::endl;
21 subsystemname_=ps.getUntrackedParameter<std::string>("subSystemFolder","Castor");
22  RatioThresh1 = ps.getUntrackedParameter<double>("ratioThreshold",0.9);
23  Qrms_DEAD = ps.getUntrackedParameter<double>("QrmsDead",0.01); //fC
25  TS_MAX = ps.getUntrackedParameter<double>("qieTSmax",6);
26 }
27 
29 
31  const edm::Run& iRun, const edm::EventSetup& iSetup)
32 {
33  char s[60];
34  if(fVerbosity>0) std::cout << "CastorDigiMonitor::beginRun (start)" << std::endl;
35  char sTileIndex[50];
36  sprintf(sTileIndex,"Tile(=moduleZ*16+sector#phi)");
37 
38  ievt_=0;
39 
40  ibooker.setCurrentFolder(subsystemname_ + "/CastorDigiMonitor");
41 
42  std::string s2 = "CASTOR QIE_capID+er+dv";
43  h2digierr=ibooker.bookProfile2D(s2,s2,14,0.,14., 16,0.,16.,100,0,1.e10,"");
44  h2digierr->getTProfile2D()->GetXaxis()->SetTitle("Module Z");
45  h2digierr->getTProfile2D()->GetYaxis()->SetTitle("Sector #phi");
46  h2digierr->getTProfile2D()->SetMaximum(1.);
47  h2digierr->getTProfile2D()->SetMinimum(QIEerrThreshold);
48  h2digierr->getTProfile2D()->SetOption("colz");
49 
50  sprintf(s,"CASTORreportSummaryMap");
51  h2repsum=ibooker.bookProfile2D(s,s,14,0.,14., 16,0.,16.,100,0,1.e10,"");
52  h2repsum->getTProfile2D()->GetXaxis()->SetTitle("Module Z");
53  h2repsum->getTProfile2D()->GetYaxis()->SetTitle("Sector #phi");
54  h2repsum->getTProfile2D()->SetMaximum(1.);
55  h2repsum->getTProfile2D()->SetMinimum(QIEerrThreshold);
56  h2repsum->getTProfile2D()->SetOption("colz");
57 
58  sprintf(s,"CASTOR DeadChannelsMap");
59  h2status = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
60  h2status->getTH2F()->GetXaxis()->SetTitle("Module Z");
61  h2status->getTH2F()->GetYaxis()->SetTitle("Sector #phi");
62  h2status->getTH2F()->SetOption("colz");
63 
64  sprintf(s,"CASTOR TSmax Significance Map");
65  h2TSratio = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
66  h2TSratio->getTH2F()->GetXaxis()->SetTitle("Module Z");
67  h2TSratio->getTH2F()->GetYaxis()->SetTitle("Sector #phi");
68  h2TSratio->getTH2F()->SetOption("colz");
69 
70  sprintf(s,"CASTOR TSmax Significance All chan");
71  hTSratio = ibooker.book1D(s,s,105,0.,1.05);
72 
73  sprintf(s,"DigiSize");
74  hdigisize = ibooker.book1D(s,s,20,0.,20.);
75  sprintf(s,"ModuleZ(fC)_allTS");
76  hModule = ibooker.book1D(s,s,14,0.,14.);
77  hModule->getTH1F()->GetXaxis()->SetTitle("ModuleZ");
78  hModule->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
79  sprintf(s,"Sector #phi(fC)_allTS");
80  hSector = ibooker.book1D(s,s,16,0.,16.);
81  hSector->getTH1F()->GetXaxis()->SetTitle("Sector #phi");
82  hSector->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
83 
84  sprintf(s,"QfC=f(x=Tile y=TS) (cumulative)");
85  h2QtsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
86  h2QtsvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
87  h2QtsvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
88  h2QtsvsCh->getTH2F()->SetOption("colz");
89 
90  sprintf(s,"QmeanfC=f(Tile TS)");
91  h2QmeantsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
92  h2QmeantsvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
93  h2QmeantsvsCh->getTH2F()->GetYaxis()->SetTitle("Time Slice");
94  h2QmeantsvsCh->getTH2F()->SetOption("colz");
95 
96  sprintf(s,"QrmsfC=f(Tile TS)");
97  h2QrmsTSvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
98  h2QrmsTSvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
99  h2QrmsTSvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
100  h2QrmsTSvsCh->getTH2F()->SetOption("colz");
101 
102  sprintf(s,"CASTOR data quality");
103  h2qualityMap = ibooker.book2D(s,s,14, 0,14, 16, 0,16);
104  h2qualityMap->getTH2F()->GetXaxis()->SetTitle("module Z");
105  h2qualityMap->getTH2F()->GetYaxis()->SetTitle("Sector #phi");
106  h2qualityMap->getTH2F()->SetOption("colz");
107 
108  hReport = ibooker.bookFloat("CASTOR reportSummary");
109 
110  sprintf(s,"QmeanfC_map(allTS)");
111  h2QmeanMap = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
112  h2QmeanMap->getTH2F()->GetXaxis()->SetTitle("Module Z");
113  h2QmeanMap->getTH2F()->GetYaxis()->SetTitle("Sector #phi");
114  h2QmeanMap->getTH2F()->SetOption("textcolz");
115 
116  for(int ts=0; ts<=1; ts++) {
117  sprintf(s,"QIErms_TS=%d",ts);
118  hQIErms[ts] = ibooker.book1D(s,s,1000,0.,100.);
119  hQIErms[ts]->getTH1F()->GetXaxis()->SetTitle("QIErms(fC)");
120  }
121 
122  for(int ind=0; ind<224; ind++) for(int ts=0; ts<10; ts++)
123  QrmsTS[ind][ts] = QmeanTS[ind][ts]= 0.;
124 
125  if(fVerbosity>0)
126 std::cout<<"CastorDigiMonitor::bookingHist(end)"<<std::endl;
127  return;
128 }
129 
130 
132  const CastorDbService& cond) {
133  if(fVerbosity>0) std::cout << "CastorDigiMonitor::processEvent (begin)"<< std::endl;
134 
135  if(castorDigis.empty()) {
136  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++)
137  h2repsum->Fill(mod,sec,0.);
138  if(fVerbosity>0) std::cout<<"CastorPSMonitor::processEvent NO Castor Digis"<<std::endl;
139  return;
140  }
141 
142  for(CastorDigiCollection::const_iterator j=castorDigis.begin();
143  j!=castorDigis.end(); j++)
144  {
145  const CastorDataFrame digi = (const CastorDataFrame)(*j);
146 
147  int capid1 = digi.sample(0).capid();
148  hdigisize->Fill(digi.size());
149  double sum = 0.;
150  int err=0, err2=0;
151  int module = digi.id().module()-1;
152  int sector = digi.id().sector()-1;
153  for (int i=0; i<digi.size(); i++) {
154  int capid = digi.sample(i).capid();
155  int dv = digi.sample(i).dv();
156  int er = digi.sample(i).er();
157  int rawd = digi.sample(i).adc();
158  rawd = rawd&0x7F;
159  err |= (capid != capid1) | er<<1 | (!dv)<<2; // =0
160  err2 += (capid != capid1) | er | (!dv); // =0
161 // if(err !=0) continue;
162  int ind = ModSecToIndex(module,sector);
163  h2QtsvsCh->Fill(ind,i,LedMonAdc2fc[rawd]);
164  float q = LedMonAdc2fc[rawd];
165  sum += q; // sum += LedMonAdc2fc[rawd];
166  QrmsTS[ind][i] += (q*q);
167  QmeanTS[ind][i] += q;
168  if(err != 0 && fVerbosity>0)
169  std::cout<<"event/idigi=" <<ievt_<<"/"<<i<<" cap=cap1_dv_er_err: "<<
170  capid <<"="<< capid1 <<" "<< dv <<" "<< er<<" "<< err << std::endl;
171  if(capid1 < 3) capid1 = capid+1;
172  else capid1 = 0;
173  }
174  h2digierr->Fill(module,sector,err);
175  h2repsum->Fill(module,sector,1.-err2/digi.size());
176 // hBunchOcc->Fill(iBunch,sum);
177  } //end for(CastorDigiCollection::const_iterator ...
178 
179  ievt_++;
180 
181  const float repChanBAD = 0.9;
182  const float repChanWarning = 0.95;
183  if(ievt_ %100 != 0) return;
184  float ModuleSum[14], SectorSum[16];
185  for(int m=0; m<14; m++) ModuleSum[m]=0.;
186  for(int s=0; s<16; s++) SectorSum[s]=0.;
187  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++) {
188  for(int ts=0; ts<=1; ts++) {
189  int ind = ModSecToIndex(mod,sec);
190  double Qmean = QmeanTS[ind][ts]/ievt_;
191  double Qrms = sqrt(QrmsTS[ind][ts]/ievt_ - Qmean*Qmean);
192  hQIErms[ts]->Fill(Qrms);
193  }
194 
195  double sum=0.;
196  for(int ts=1; ts<=TS_MAX; ts++) {
197  int ind = ModSecToIndex(mod,sec) + 1;
198  double a=h2QtsvsCh->getTH2F()->GetBinContent(ind,ts);
199  h2QmeantsvsCh->getTH2F()->SetBinContent(ind,ts,a/double(ievt_));
200  sum += a;
201  double Qmean = QmeanTS[ind-1][ts-1]/ievt_;
202  double Qrms = QrmsTS[ind-1][ts-1]/ievt_ - Qmean*Qmean;
203  h2QrmsTSvsCh->getTH2F()->SetBinContent(ind,ts,sqrt(Qrms));
204  }
205  sum /= double(ievt_);
206  ModuleSum[mod] += sum;
207  SectorSum[sec] += sum;
208  float isum = float(int(sum*10.+0.5))/10.;
209  h2QmeanMap->getTH2F()->SetBinContent(mod+1,sec+1,isum);
210  } // end for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++)
211 
212  for(int mod=0; mod<14; mod++)
213  hModule->getTH1F()->SetBinContent(mod+1,ModuleSum[mod]);
214  for(int sec=0; sec<16; sec++)
215  hSector->getTH1F()->SetBinContent(sec+1,SectorSum[sec]);
216 
217  int nGoodCh = 0;
218  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16;sec++) {
219  int ind = ModSecToIndex(mod,sec);
220  double Qmean = QmeanTS[ind][TSped]/ievt_;
221  double Qrms = QrmsTS[ind][TSped]/ievt_ - Qmean*Qmean;
222  float ChanStatus = 0.;
223  if(Qrms < Qrms_DEAD) ChanStatus = 1.;
224  h2status->getTH2F()->SetBinContent(mod+1,sec+1,ChanStatus);
225 
226  int tsm = 0;
227  float am=0., amin = 10000.;
228  for(int ts=0; ts<TS_MAX; ts++) {
229  float a = h2QmeantsvsCh->getTH2F()->GetBinContent(ind+1,ts+1);
230  if(am < a) {am = a; tsm = ts;}
231  if(a < amin) amin = a;
232  }
233 
234  double sum = 0.;
235  for(int ts=0; ts<TS_MAX; ts++) if(ts != tsm)
236  sum += h2QmeantsvsCh->getTH2F()->GetBinContent(ind+1,ts+1);
237  float r = 1.; // worth case - no peak
238 // if(am > 0.) r = sum/(TS_MAX-1)/am;
239  if(am > amin) r = (sum/(TS_MAX-1)-amin)/(am - amin);
240  h2TSratio->getTH2F()->SetBinContent(mod+1,sec+1,r);
241  hTSratio->Fill(r);
242  float statusTS = 1.0;
243  if(r > RatioThresh1) statusTS = repChanWarning;
244  else if(r > 0.99) statusTS = repChanBAD;
245  float gChanStatus = statusTS;
246  if(ChanStatus > 0.) gChanStatus = repChanBAD; // RMS
247 // if(h2digierr->getTProfile2D()->GetBinContent(mod+1,sec+1)>QIEerrThreshold) gChanStatus = repChanBAD;
248  h2qualityMap->getTH2F()->SetBinContent(mod+1,sec+1,gChanStatus);
249  if(gChanStatus > repChanBAD) ++nGoodCh;
250  }
251  hReport->Fill(float(nGoodCh)/224.);
252  return;
253  }
254 
256  int ind = sector + module*16;
257  if(ind>223) ind=223;
258  return(ind);
259 }
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * h2TSratio
double QrmsTS[224][10]
TProfile2D * getTProfile2D() const
MonitorElement * h2qualityMap
int sector() const
get the sector (1-16)
MonitorElement * hTSratio
TH1F * getTH1F() const
int adc() const
get the ADC sample
Definition: HcalQIESample.h:22
const HcalQIESample & sample(int i) const
access a sample
std::vector< CastorDataFrame >::const_iterator const_iterator
bool dv() const
is the Data Valid bit set?
Definition: HcalQIESample.h:28
int module() const
get the module (1-2 for EM, 1-12 for HAD)
double isum
void Fill(long long x)
MonitorElement * h2status
CastorDigiMonitor(const edm::ParameterSet &ps)
MonitorElement * h2digierr
MonitorElement * hModule
MonitorElement * bookProfile2D(Args &&...args)
Definition: DQMStore.h:166
T sqrt(T t)
Definition: SSEVec.h:18
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
MonitorElement * h2QmeanMap
MonitorElement * h2QtsvsCh
TH2F * getTH2F() const
std::string subsystemname_
const_iterator end() const
MonitorElement * h2QrmsTSvsCh
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
static const float LedMonAdc2fc[128]
int capid() const
get the Capacitor id
Definition: HcalQIESample.h:26
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
int ModSecToIndex(int module, int sector)
MonitorElement * hdigisize
MonitorElement * hSector
Definition: plugin.cc:24
double a
Definition: hdecay.h:121
void processEvent(const CastorDigiCollection &cast, const CastorDbService &cond)
double QmeanTS[224][10]
MonitorElement * h2repsum
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:112
MonitorElement * hQIErms[10]
const HcalCastorDetId & id() const
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
int size() const
total number of samples in the digi
Definition: vlib.h:208
bool er() const
is the error bit set?
Definition: HcalQIESample.h:30
MonitorElement * hReport
const_iterator begin() const
Definition: Run.h:43
MonitorElement * h2QmeantsvsCh