CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelDigiModule.cc
Go to the documentation of this file.
6 // STL
7 #include <vector>
8 #include <memory>
9 #include <string>
10 #include <iostream>
11 #include <stdlib.h>
12 #include <sstream>
13 #include <cstdio>
14 
15 // Data Formats
22 
23 //
24 // Constructors
25 //
27  ncols_(416),
28  nrows_(160)
29 {
30 }
33  id_(id),
34  ncols_(416),
35  nrows_(160)
36 {
37 }
39 SiPixelDigiModule::SiPixelDigiModule(const uint32_t& id, const int& ncols, const int& nrows) :
40  id_(id),
41  ncols_(ncols),
42  nrows_(nrows)
43 {
44 }
45 //
46 // Destructor
47 //
49 //
50 // Book histograms
51 //
52 void SiPixelDigiModule::book(const edm::ParameterSet& iConfig, DQMStore::IBooker & iBooker, int type, bool twoD, bool hiRes, bool reducedSet, bool additInfo, bool isUpgrade) {
53 
54  //isUpgrade = iConfig.getUntrackedParameter<bool>("isUpgrade");
55 
56  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
57  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
58  bool isHalfModule = false;
59  if(barrel){
60  if (!isUpgrade) {
61  isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule();
62  } else if (isUpgrade) {
63  isHalfModule = PixelBarrelNameUpgrade(DetId(id_)).isHalfModule();
64  }
65  }
66 
67  std::string hid;
68  // Get collection name and instantiate Histo Id builder
69  edm::InputTag src = iConfig.getParameter<edm::InputTag>( "src" );
70 
71  int nbinx=ncols_/2, nbiny=nrows_/2;
72  std::string twodtitle = "Number of Digis (1bin=four pixels)";
73  std::string pxtitle = "Number of Digis (1bin=two columns)";
74  std::string pytitle = "Number of Digis (1bin=two rows)";
75  std::string twodroctitle = "ROC Occupancy (1bin=one ROC)";
76  std::string twodzeroOccroctitle = "Zero Occupancy ROC Map (1bin=one ROC) for ";
77  if(hiRes){
78  nbinx = ncols_;
79  nbiny = nrows_;
80  twodtitle = "Number of Digis (1bin=one pixel)";
81  pxtitle = "Number of Digis (1bin=one column)";
82  pytitle = "Number of Digis (1bin=one row)";
83  }
84  if(type==0){
85  SiPixelHistogramId* theHistogramId = new SiPixelHistogramId( src.label() );
86  // Number of digis
87  hid = theHistogramId->setHistoId("ndigis",id_);
88  meNDigis_ = iBooker.book1D(hid,"Number of Digis",25,0.,25.);
89  meNDigis_->setAxisTitle("Number of digis",1);
90  // Charge in ADC counts
91  hid = theHistogramId->setHistoId("adc",id_);
92  meADC_ = iBooker.book1D(hid,"Digi charge",128,0.,256.);
93  meADC_->setAxisTitle("ADC counts",1);
94  if(!reducedSet)
95  {
96  if(twoD){
97  if(additInfo){
98  // 2D hit map
99  hid = theHistogramId->setHistoId("hitmap",id_);
100  mePixDigis_ = iBooker.book2D(hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
101  mePixDigis_->setAxisTitle("Columns",1);
102  mePixDigis_->setAxisTitle("Rows",2);
103  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
104  }
105  }
106  else{
107  // projections of 2D hit map
108  hid = theHistogramId->setHistoId("hitmap",id_);
109  mePixDigis_px_ = iBooker.book1D(hid+"_px",pxtitle,nbinx,0.,float(ncols_));
110  mePixDigis_py_ = iBooker.book1D(hid+"_py",pytitle,nbiny,0.,float(nrows_));
111  mePixDigis_px_->setAxisTitle("Columns",1);
112  mePixDigis_py_->setAxisTitle("Rows",1);
113  }
114  }
115  delete theHistogramId;
116 
117  }
118 
119  if(type==1 && barrel){
120  uint32_t DBladder;
121  if (!isUpgrade) { DBladder = PixelBarrelName(DetId(id_)).ladderName();}
122  else { DBladder = PixelBarrelNameUpgrade(DetId(id_)).ladderName();}
123  char sladder[80]; sprintf(sladder,"Ladder_%02i",DBladder);
124  hid = src.label() + "_" + sladder;
125  if(isHalfModule) hid += "H";
126  else hid += "F";
127  // Number of digis
128  meNDigisLad_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
129  meNDigisLad_->setAxisTitle("Number of digis",1);
130  // Charge in ADC counts
131  meADCLad_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
132  meADCLad_->setAxisTitle("ADC counts",1);
133  if(!reducedSet)
134  {
135  if(twoD){
136  // 2D hit map
137  mePixDigisLad_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
138  mePixDigisLad_->setAxisTitle("Columns",1);
139  mePixDigisLad_->setAxisTitle("Rows",2);
140  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
141  }
142  else{
143  // projections of 2D hit map
144  mePixDigisLad_px_ = iBooker.book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
145  mePixDigisLad_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
146  mePixDigisLad_px_->setAxisTitle("Columns",1);
147  mePixDigisLad_py_->setAxisTitle("Rows",1);
148  }
149  }
150  }
151  if(type==2 && barrel){
152  uint32_t DBlayer;
153  if (!isUpgrade) { DBlayer = PixelBarrelName(DetId(id_)).layerName(); }
154  else { DBlayer = PixelBarrelNameUpgrade(DetId(id_)).layerName(); }
155  char slayer[80]; sprintf(slayer,"Layer_%i",DBlayer);
156  hid = src.label() + "_" + slayer;
157  if(!additInfo){
158  // Number of digis
159  meNDigisLay_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
160  meNDigisLay_->setAxisTitle("Number of digis",1);
161  // Charge in ADC counts
162  meADCLay_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
163  meADCLay_->setAxisTitle("ADC counts",1);
164  }
165  if(!reducedSet){
166  if(twoD || additInfo){
167  // 2D hit map
168  if(isHalfModule){
169  mePixDigisLay_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),2*nbiny,0.,float(2*nrows_));
170  }
171  else{
172  mePixDigisLay_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
173 
174  }
175  mePixDigisLay_->setAxisTitle("Columns",1);
176  mePixDigisLay_->setAxisTitle("Rows",2);
177 
178  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
179  int yROCbins[3] = {18,30,42};
180  mePixRocsLay_ = iBooker.book2D("rocmap_"+hid,twodroctitle,32,0.,32.,yROCbins[DBlayer-1],1.5,1.5+float(yROCbins[DBlayer-1]/2));
181  mePixRocsLay_->setAxisTitle("ROCs per Module",1);
182  mePixRocsLay_->setAxisTitle("ROCs per 1/2 Ladder",2);
183  meZeroOccRocsLay_ = iBooker.book2D("zeroOccROC_map",twodzeroOccroctitle+hid,32,0.,32.,yROCbins[DBlayer-1],1.5,1.5+float(yROCbins[DBlayer-1]/2));
184  meZeroOccRocsLay_->setAxisTitle("ROCs per Module",1);
185  meZeroOccRocsLay_->setAxisTitle("ROCs per 1/2 Ladder",2);
186  }
187  if(!twoD && !additInfo){
188  // projections of 2D hit map
189  mePixDigisLay_px_ = iBooker.book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
190  if(isHalfModule){
191  mePixDigisLay_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,2*nbiny,0.,float(2*nrows_));
192  }
193  else{
194  mePixDigisLay_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
195  }
196  mePixDigisLay_px_->setAxisTitle("Columns",1);
197  mePixDigisLay_py_->setAxisTitle("Rows",1);
198  }
199  }
200  }
201  if(type==3 && barrel){
202  uint32_t DBmodule;
203  if (!isUpgrade) { DBmodule = PixelBarrelName(DetId(id_)).moduleName(); }
204  else { DBmodule = PixelBarrelNameUpgrade(DetId(id_)).moduleName(); }
205  char smodule[80]; sprintf(smodule,"Ring_%i",DBmodule);
206  hid = src.label() + "_" + smodule;
207  // Number of digis
208  meNDigisPhi_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
209  meNDigisPhi_->setAxisTitle("Number of digis",1);
210  // Charge in ADC counts
211  meADCPhi_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
212  meADCPhi_->setAxisTitle("ADC counts",1);
213  if(!reducedSet)
214  {
215  if(twoD){
216 
217  // 2D hit map
218  if(isHalfModule){
219  mePixDigisPhi_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),2*nbiny,0.,float(2*nrows_));
220  }
221  else {
222  mePixDigisPhi_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
223  }
224  mePixDigisPhi_->setAxisTitle("Columns",1);
225  mePixDigisPhi_->setAxisTitle("Rows",2);
226  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
227  }
228  else{
229  // projections of 2D hit map
230  mePixDigisPhi_px_ = iBooker.book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
231  if(isHalfModule){
232  mePixDigisPhi_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,2*nbiny,0.,float(2*nrows_));
233  }
234  else{
235  mePixDigisPhi_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
236  }
237  mePixDigisPhi_px_->setAxisTitle("Columns",1);
238  mePixDigisPhi_py_->setAxisTitle("Rows",1);
239  }
240  }
241  }
242  if(type==4 && endcap){
243  uint32_t blade;
244  if (!isUpgrade) { blade= PixelEndcapName(DetId(id_)).bladeName(); }
245  else { blade= PixelEndcapNameUpgrade(DetId(id_)).bladeName(); }
246 
247  char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
248  hid = src.label() + "_" + sblade;
249  // Number of digis
250  meNDigisBlade_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
251  meNDigisBlade_->setAxisTitle("Number of digis",1);
252  // Charge in ADC counts
253  meADCBlade_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
254  meADCBlade_->setAxisTitle("ADC counts",1);
255  }
256  if(type==5 && endcap){
257  uint32_t disk;
258  if (!isUpgrade) { disk = PixelEndcapName(DetId(id_)).diskName(); }
259  else { disk = PixelEndcapNameUpgrade(DetId(id_)).diskName(); }
260 
261  char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
262  hid = src.label() + "_" + sdisk;
263  if(!additInfo){
264  // Number of digis
265  meNDigisDisk_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
266  meNDigisDisk_->setAxisTitle("Number of digis",1);
267  // Charge in ADC counts
268  meADCDisk_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
269  meADCDisk_->setAxisTitle("ADC counts",1);
270  }
271  if(additInfo){
272  mePixDigisDisk_ = iBooker.book2D("hitmap_"+hid,twodtitle,260,0.,260.,160,0.,160.);
273  mePixDigisDisk_->setAxisTitle("Columns",1);
274  mePixDigisDisk_->setAxisTitle("Rows",2);
275  //ROC information in disks
276  mePixRocsDisk_ = iBooker.book2D("rocmap_"+hid,twodroctitle,26,0.,26.,24,1.,13.);
277  mePixRocsDisk_ ->setAxisTitle("ROCs per Module (2 Panels)",1);
278  mePixRocsDisk_ ->setAxisTitle("Blade Number",2);
279  meZeroOccRocsDisk_ = iBooker.book2D("zeroOccROC_map",twodzeroOccroctitle+hid,26,0.,26.,24,1.,13.);
280  meZeroOccRocsDisk_ ->setAxisTitle("Zero-Occupancy ROCs per Module (2 Panels)",1);
281  meZeroOccRocsDisk_ ->setAxisTitle("Blade Number",2);
282  }
283  }
284  if(type==6 && endcap){
285  uint32_t panel;
286  uint32_t module;
287  if (!isUpgrade) {
288  panel= PixelEndcapName(DetId(id_)).pannelName();
290  } else {
293  }
294 
295  char slab[80]; sprintf(slab, "Panel_%i_Ring_%i",panel, module);
296  hid = src.label() + "_" + slab;
297  // Number of digis
298  meNDigisRing_ = iBooker.book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
299  meNDigisRing_->setAxisTitle("Number of digis",1);
300  // Charge in ADC counts
301  meADCRing_ = iBooker.book1D("adc_" + hid,"Digi charge",128,0.,256.);
302  meADCRing_->setAxisTitle("ADC counts",1);
303  if(!reducedSet)
304  {
305  if(twoD){
306  // 2D hit map
307  mePixDigisRing_ = iBooker.book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
308  mePixDigisRing_->setAxisTitle("Columns",1);
309  mePixDigisRing_->setAxisTitle("Rows",2);
310  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
311  }
312  else{
313  // projections of 2D hit map
314  mePixDigisRing_px_ = iBooker.book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
315  mePixDigisRing_py_ = iBooker.book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
316  mePixDigisRing_px_->setAxisTitle("Columns",1);
317  mePixDigisRing_py_->setAxisTitle("Rows",1);
318  }
319  }
320  }
321 }
322 
323 
324 //
325 // Fill histograms
326 //
328  MonitorElement* combBarrel, MonitorElement* chanBarrel, MonitorElement* chanBarrelL1, MonitorElement* chanBarrelL2, MonitorElement* chanBarrelL3, MonitorElement* chanBarrelL4, MonitorElement* combEndcap,
329  bool modon, bool ladon, bool layon, bool phion,
330  bool bladeon, bool diskon, bool ringon,
331  bool twoD, bool reducedSet, bool twoDimModOn, bool twoDimOnlyLayDisk,
332  int &nDigisA, int &nDigisB, bool isUpgrade) {
333  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
334  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
335  bool isHalfModule = false;
336  uint32_t DBladder = 0;
337  if(barrel && !isUpgrade){
338  isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule();
339  DBladder = PixelBarrelName(DetId(id_)).ladderName();
340  } else if (barrel && isUpgrade) {
341  isHalfModule = PixelBarrelNameUpgrade(DetId(id_)).isHalfModule();
343  }
344 
345  edm::DetSetVector<PixelDigi>::const_iterator isearch = input.find(id_); // search digis of detid
346 
347  unsigned int numberOfDigisMod = 0;
348  int msize;
349  if (isUpgrade) {msize=10;} else {msize=8;}
350  int numberOfDigis[msize]; for(int i=0; i!=msize; i++) numberOfDigis[i]=0;
351  nDigisA=0; nDigisB=0;
352  if( isearch != input.end() ) { // Not an empty iterator
353 
354  // Look at digis now
356  for(di = isearch->data.begin(); di != isearch->data.end(); di++) {
357  int adc = di->adc(); // charge
358  int col = di->column(); // column
359  int row = di->row(); // row
360  numberOfDigisMod++;
361 
362  int DBlayer = 0;
363  int DBmodule =0;
364 
365  if (!isUpgrade) {
367  DBlayer = PixelBarrelName(DetId(id_)).layerName();
368  DBmodule = PixelBarrelName(DetId(id_)).moduleName();
369  if(barrel){
370  if(isHalfModule){
371  if(DBshell==PixelBarrelName::pI||DBshell==PixelBarrelName::pO){
372  numberOfDigis[0]++; nDigisA++;
373  if(DBlayer==1) numberOfDigis[2]++;
374  if(DBlayer==2) numberOfDigis[3]++;
375  if(DBlayer==3) numberOfDigis[4]++;
376  }
377  if(DBshell==PixelBarrelName::mI||DBshell==PixelBarrelName::mO){
378  numberOfDigis[1]++; nDigisB++;
379  if(DBlayer==1) numberOfDigis[5]++;
380  if(DBlayer==2) numberOfDigis[6]++;
381  if(DBlayer==3) numberOfDigis[7]++;
382  }
383  }else{
384  if(row<80){
385  numberOfDigis[0]++; nDigisA++;
386  if(DBlayer==1) numberOfDigis[2]++;
387  if(DBlayer==2) numberOfDigis[3]++;
388  if(DBlayer==3) numberOfDigis[4]++;
389  }else{
390  numberOfDigis[1]++; nDigisB++;
391  if(DBlayer==1) numberOfDigis[5]++;
392  if(DBlayer==2) numberOfDigis[6]++;
393  if(DBlayer==3) numberOfDigis[7]++;
394  }
395  }
396  }
397  } else if (isUpgrade) {
398  //PixelBarrelNameUpgrade::Shell DBshell = PixelBarrelNameUpgrade(DetId(id_)).shell();
401  if(barrel){
402  if(row<80){
403  numberOfDigis[0]++; nDigisA++;
404  if(DBlayer==1) numberOfDigis[2]++;
405  if(DBlayer==2) numberOfDigis[3]++;
406  if(DBlayer==3) numberOfDigis[4]++;
407  if(DBlayer==4) numberOfDigis[5]++;
408  }else{
409  numberOfDigis[1]++; nDigisB++;
410  if(DBlayer==1) numberOfDigis[6]++;
411  if(DBlayer==2) numberOfDigis[7]++;
412  if(DBlayer==3) numberOfDigis[8]++;
413  if(DBlayer==4) numberOfDigis[9]++;
414  }
415  }
416  }
417 
418  if(modon){
419  if(!reducedSet){
420  if(twoD) {
421  if(twoDimModOn) (mePixDigis_)->Fill((float)col,(float)row);
422  }
423  else {
424  (mePixDigis_px_)->Fill((float)col);
425  (mePixDigis_py_)->Fill((float)row);
426  }
427  }
428  (meADC_)->Fill((float)adc);
429  }
430  if(ladon && barrel){
431  (meADCLad_)->Fill((float)adc);
432  if(!reducedSet){
433  if(twoD) (mePixDigisLad_)->Fill((float)col,(float)row);
434  else {
435  (mePixDigisLad_px_)->Fill((float)col);
436  (mePixDigisLad_py_)->Fill((float)row);
437  }
438  }
439  }
440  if((layon || twoDimOnlyLayDisk) && barrel){
441  if(!twoDimOnlyLayDisk) (meADCLay_)->Fill((float)adc);
442  if(!reducedSet){
443  if((layon && twoD) || twoDimOnlyLayDisk){
444  //ROC histos...
445  float rocx = (float)col/52. + 8.0*float(DBmodule-1);
446  float rocy = (float)row/160.+float(DBladder);
447  //Shift 1st ladder (half modules) up by 1 bin
448  if(DBladder==1) rocy = rocy + 0.5;
449  mePixRocsLay_->Fill(rocx,rocy);
450 
451  if(isHalfModule && DBladder==1){
452  (mePixDigisLay_)->Fill((float)col,(float)row+80);
453  }
454  else (mePixDigisLay_)->Fill((float)col,(float)row);
455  }
456  if((layon && !twoD) && !twoDimOnlyLayDisk){
457  (mePixDigisLay_px_)->Fill((float)col);
458  if(isHalfModule && DBladder==1) {
459  (mePixDigisLay_py_)->Fill((float)row+80);
460  }
461  else (mePixDigisLay_py_)->Fill((float)row);
462  }
463  }
464  }
465  if(phion && barrel){
466  (meADCPhi_)->Fill((float)adc);
467  if(!reducedSet)
468  {
469  if(twoD){
470  if(isHalfModule && DBladder==1){
471  (mePixDigisPhi_)->Fill((float)col,(float)row+80);
472  }
473  else (mePixDigisPhi_)->Fill((float)col,(float)row);
474  }
475  else {
476  (mePixDigisPhi_px_)->Fill((float)col);
477  if(isHalfModule && DBladder==1) {
478  (mePixDigisPhi_py_)->Fill((float)row+80);
479  }
480  else (mePixDigisPhi_py_)->Fill((float)row);
481  }
482  }
483  }
484  if(bladeon && endcap){
485  (meADCBlade_)->Fill((float)adc);
486  }
487 
488  if((diskon || twoDimOnlyLayDisk) && endcap){
489  if(!twoDimOnlyLayDisk) (meADCDisk_)->Fill((float)adc);
490  if(twoDimOnlyLayDisk){
491  (mePixDigisDisk_)->Fill((float)col,(float)row);
492  //ROC monitoring
493  int DBpanel;
494  int DBblade;
495  if (!isUpgrade) {
496  DBpanel= PixelEndcapName(DetId(id_)).pannelName();
497  DBblade= PixelEndcapName(DetId(id_)).bladeName();
498  } else {
501  }
502  float offx = 0.;
503  //This crazy offset takes into account the roc and module fpix configuration
504  for (int i = DBpanel; i < DBmodule; ++i) {offx = offx + float(5+DBpanel-i);}
505  float rocx = (float)col/52. + offx + 14.0*float(DBpanel-1);
506  float rocy = (float)row/160.+float(DBblade);
507  mePixRocsDisk_->Fill(rocx,rocy);
508  }
509  }
510  if(ringon && endcap){
511  (meADCRing_)->Fill((float)adc);
512  if(!reducedSet)
513  {
514  if(twoD) (mePixDigisRing_)->Fill((float)col,(float)row);
515  else {
516  (mePixDigisRing_px_)->Fill((float)col);
517  (mePixDigisRing_py_)->Fill((float)row);
518  }
519  }
520  }
521  }
522  if(modon) (meNDigis_)->Fill((float)numberOfDigisMod);
523  if(ladon && barrel) (meNDigisLad_)->Fill((float)numberOfDigisMod);
524  if(layon && barrel && !twoDimOnlyLayDisk) (meNDigisLay_)->Fill((float)numberOfDigisMod);
525  if(phion && barrel) (meNDigisPhi_)->Fill((float)numberOfDigisMod);
526  if(bladeon && endcap) (meNDigisBlade_)->Fill((float)numberOfDigisMod);
527  if(diskon && endcap && !twoDimOnlyLayDisk) (meNDigisDisk_)->Fill((float)numberOfDigisMod);
528  if(ringon && endcap) (meNDigisRing_)->Fill((float)numberOfDigisMod);
529  if(barrel){
530  if(combBarrel) combBarrel->Fill((float)numberOfDigisMod);
531  if(chanBarrel){ if(numberOfDigis[0]>0) chanBarrel->Fill((float)numberOfDigis[0]); if(numberOfDigis[1]>0) chanBarrel->Fill((float)numberOfDigis[1]); }
532  if(chanBarrelL1){ if(numberOfDigis[2]>0) chanBarrelL1->Fill((float)numberOfDigis[2]); }
533  if(chanBarrelL2){ if(numberOfDigis[3]>0) chanBarrelL2->Fill((float)numberOfDigis[3]); }
534  if(chanBarrelL3){ if(numberOfDigis[4]>0) chanBarrelL3->Fill((float)numberOfDigis[4]); }
535  if(chanBarrelL4){ if(numberOfDigis[5]>0) chanBarrelL4->Fill((float)numberOfDigis[5]); }
536  }else if(endcap){
537  if(combEndcap) combEndcap->Fill((float)numberOfDigisMod);
538  }
539  }
540 
541  //std::cout<<"numberOfDigis for this module: "<<numberOfDigis<<std::endl;
542  return numberOfDigisMod;
543 }
544 
545 // This was done in the Source file, but is moved to the Module for thread safety reasons. Using ME that is booked here.
549 }
550 
551 //Moved from source. Gets the zero and low eff ROCs from each module. Called in source for each module.
553  int nZeroROC = 0;
554  int nLoEffROC = 0;
555  float SF = 1.0;
558  for (int i = 1; i < mePixRocsDisk_->getNbinsX(); ++i){
559  for (int j = 1; j < mePixRocsDisk_->getNbinsY(); ++j){
560  float localX = float(i) - 0.5;
561  float localY = float(j)/2.0 + 0.75;
562  if (mePixRocsDisk_->getBinContent(i,j) < 1 ) {nZeroROC++; meZeroOccRocsDisk_->Fill(localX,localY);}
563  if (mePixRocsDisk_->getBinContent(i,j)*SF < 0.25){nLoEffROC++;}
564  }
565  }
566  return std::pair<int,int>(nZeroROC,nLoEffROC);
567  }
570  for (int i = 1; i < mePixRocsLay_->getNbinsX(); ++i){
571  for (int j = 1; j < mePixRocsLay_->getNbinsY(); ++j){
572  float localX = float(i) - 0.5;
573  float localY = float(j)/2.0 + 1.25;
574  if (mePixRocsLay_->getBinContent(i,j) < 1 ) {nZeroROC++; meZeroOccRocsLay_->Fill(localX,localY);}
575  if (mePixRocsLay_->getBinContent(i,j)*SF < 0.25){nLoEffROC++;}
576  }
577  }
578  return std::pair<int,int>(nZeroROC,nLoEffROC);
579  }
580  return std::pair<int,int>(0,0);
581 }
int adc(sample_type sample)
get the ADC sample (12 bits)
int plaquetteName() const
plaquetteId (in pannel)
int fill(const edm::DetSetVector< PixelDigi > &input, MonitorElement *combBarrel, MonitorElement *chanBarrel, MonitorElement *chanBarrelL1, MonitorElement *chanBarrelL2, MonitorElement *chanBarrelL3, MonitorElement *chanBarrelL4, MonitorElement *combEndcap, const bool modon, const bool ladon, const bool layon, const bool phion, const bool bladeon, const bool diskon, const bool ringon, const bool twoD, const bool reducedSet, const bool twoDimModOn, const bool twoDimOnlyLayDisk, int &nDigisA, int &nDigisB, bool isUpgrade)
Fill histograms.
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
MonitorElement * meNDigisRing_
MonitorElement * mePixDigisLad_py_
MonitorElement * mePixDigis_
MonitorElement * mePixDigisPhi_py_
void book(const edm::ParameterSet &iConfig, DQMStore::IBooker &iBooker, int type=0, bool twoD=true, bool hiRes=false, bool reducedSet=false, bool additInfo=false, bool isUpgrade=false)
Book histograms.
iterator find(det_id_type id)
Definition: DetSetVector.h:292
int moduleName() const
module id (index in z)
MonitorElement * meADCRing_
MonitorElement * meNDigis_
MonitorElement * mePixDigisRing_px_
MonitorElement * mePixDigis_py_
bool isHalfModule() const
full or half module
std::pair< int, int > getZeroLoEffROCs()
~SiPixelDigiModule()
Destructor.
MonitorElement * meADC_
double getEntries(void) const
get # of entries
MonitorElement * mePixDigisDisk_
int plaquetteName() const
plaquetteId (in pannel)
static std::string const input
Definition: EdmProvDump.cc:44
int getNbinsY(void) const
get # of bins in Y-axis
void Fill(long long x)
MonitorElement * meNDigisBlade_
int bladeName() const
blade id
MonitorElement * meADCDisk_
MonitorElement * mePixDigisPhi_px_
bool isHalfModule() const
full or half module
MonitorElement * mePixDigisRing_py_
MonitorElement * mePixRocsDisk_
MonitorElement * mePixRocsLay_
MonitorElement * meZeroOccRocsLay_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * mePixDigisLay_py_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
MonitorElement * meADCLad_
int j
Definition: DBlmapReader.cc:9
MonitorElement * meADCBlade_
MonitorElement * meNDigisLay_
SiPixelDigiModule()
Default constructor.
int bladeName() const
blade id
MonitorElement * mePixDigisLad_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:363
MonitorElement * meADCLay_
int diskName() const
disk id
std::string setHistoId(std::string variable, uint32_t &rawId)
Set Histogram Id.
Definition: DetId.h:18
int ladderName() const
ladder id (index in phi)
int ladderName() const
ladder id (index in phi)
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:131
MonitorElement * mePixDigis_px_
int moduleName() const
module id (index in z)
int layerName() const
layer id
MonitorElement * mePixDigisLad_px_
Shell shell() const
MonitorElement * mePixDigisLay_
std::string const & label() const
Definition: InputTag.h:42
MonitorElement * meNDigisPhi_
double getBinContent(int binx) const
get content of bin (1-D)
int pannelName() const
pannel id
MonitorElement * meZeroOccRocsDisk_
int getNbinsX(void) const
get # of bins in X-axis
MonitorElement * mePixDigisRing_
MonitorElement * meADCPhi_
int diskName() const
disk id
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:348
MonitorElement * mePixDigisPhi_
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
int col
Definition: cuy.py:1008
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
MonitorElement * meNDigisDisk_
int pannelName() const
pannel id
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void Reset(void)
reset ME (ie. contents, errors, etc)
int layerName() const
layer id
MonitorElement * meNDigisLad_
MonitorElement * mePixDigisLay_px_