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
20 
21 //
22 // Constructors
23 //
25  ncols_(416),
26  nrows_(160)
27 {
28 }
31  id_(id),
32  ncols_(416),
33  nrows_(160)
34 {
35 }
37 SiPixelDigiModule::SiPixelDigiModule(const uint32_t& id, const int& ncols, const int& nrows) :
38  id_(id),
39  ncols_(ncols),
40  nrows_(nrows)
41 {
42 }
43 //
44 // Destructor
45 //
47 //
48 // Book histograms
49 //
50 void SiPixelDigiModule::book(const edm::ParameterSet& iConfig, int type, bool twoD, bool hiRes, bool reducedSet, bool additInfo) {
51  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
52  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
53  bool isHalfModule = false;
54  if(barrel){
55  isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule();
56  }
57 
58  std::string hid;
59  // Get collection name and instantiate Histo Id builder
60  edm::InputTag src = iConfig.getParameter<edm::InputTag>( "src" );
61 
62 
63  // Get DQM interface
64  DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
65 
66  int nbinx=ncols_/2, nbiny=nrows_/2;
67  std::string twodtitle = "Number of Digis (1bin=four pixels)";
68  std::string pxtitle = "Number of Digis (1bin=two columns)";
69  std::string pytitle = "Number of Digis (1bin=two rows)";
70  std::string twodroctitle = "ROC Occupancy (1bin=one ROC)";
71  std::string twodzeroOccroctitle = "Zero Occupancy ROC Map (1bin=one ROC) for ";
72  if(hiRes){
73  nbinx = ncols_;
74  nbiny = nrows_;
75  twodtitle = "Number of Digis (1bin=one pixel)";
76  pxtitle = "Number of Digis (1bin=one column)";
77  pytitle = "Number of Digis (1bin=one row)";
78  }
79  if(type==0){
80  SiPixelHistogramId* theHistogramId = new SiPixelHistogramId( src.label() );
81  // Number of digis
82  hid = theHistogramId->setHistoId("ndigis",id_);
83  meNDigis_ = theDMBE->book1D(hid,"Number of Digis",25,0.,25.);
84  meNDigis_->setAxisTitle("Number of digis",1);
85  // Charge in ADC counts
86  hid = theHistogramId->setHistoId("adc",id_);
87  meADC_ = theDMBE->book1D(hid,"Digi charge",128,0.,256.);
88  meADC_->setAxisTitle("ADC counts",1);
89  if(!reducedSet)
90  {
91  if(twoD){
92  if(additInfo){
93  // 2D hit map
94  hid = theHistogramId->setHistoId("hitmap",id_);
95  mePixDigis_ = theDMBE->book2D(hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
96  mePixDigis_->setAxisTitle("Columns",1);
97  mePixDigis_->setAxisTitle("Rows",2);
98  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
99  }
100  }
101  else{
102  // projections of 2D hit map
103  hid = theHistogramId->setHistoId("hitmap",id_);
104  mePixDigis_px_ = theDMBE->book1D(hid+"_px",pxtitle,nbinx,0.,float(ncols_));
105  mePixDigis_py_ = theDMBE->book1D(hid+"_py",pytitle,nbiny,0.,float(nrows_));
106  mePixDigis_px_->setAxisTitle("Columns",1);
107  mePixDigis_py_->setAxisTitle("Rows",1);
108  }
109  }
110  delete theHistogramId;
111 
112  }
113 
114  if(type==1 && barrel){
115  uint32_t DBladder = PixelBarrelName(DetId(id_)).ladderName();
116  char sladder[80]; sprintf(sladder,"Ladder_%02i",DBladder);
117  hid = src.label() + "_" + sladder;
118  if(isHalfModule) hid += "H";
119  else hid += "F";
120  // Number of digis
121  meNDigisLad_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
122  meNDigisLad_->setAxisTitle("Number of digis",1);
123  // Charge in ADC counts
124  meADCLad_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
125  meADCLad_->setAxisTitle("ADC counts",1);
126  if(!reducedSet)
127  {
128  if(twoD){
129  // 2D hit map
130  mePixDigisLad_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
131  mePixDigisLad_->setAxisTitle("Columns",1);
132  mePixDigisLad_->setAxisTitle("Rows",2);
133  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
134  }
135  else{
136  // projections of 2D hit map
137  mePixDigisLad_px_ = theDMBE->book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
138  mePixDigisLad_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
139  mePixDigisLad_px_->setAxisTitle("Columns",1);
140  mePixDigisLad_py_->setAxisTitle("Rows",1);
141  }
142  }
143  }
144  if(type==2 && barrel){
145  uint32_t DBlayer = PixelBarrelName(DetId(id_)).layerName();
146  char slayer[80]; sprintf(slayer,"Layer_%i",DBlayer);
147  hid = src.label() + "_" + slayer;
148  if(!additInfo){
149  // Number of digis
150  meNDigisLay_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
151  meNDigisLay_->setAxisTitle("Number of digis",1);
152  // Charge in ADC counts
153  meADCLay_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
154  meADCLay_->setAxisTitle("ADC counts",1);
155  }
156  if(!reducedSet){
157  if(twoD || additInfo){
158  // 2D hit map
159  if(isHalfModule){
160  mePixDigisLay_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),2*nbiny,0.,float(2*nrows_));
161  }
162  else{
163  mePixDigisLay_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
164 
165  }
166  mePixDigisLay_->setAxisTitle("Columns",1);
167  mePixDigisLay_->setAxisTitle("Rows",2);
168 
169  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
170  int yROCbins[3] = {18,30,42};
171  mePixRocsLay_ = theDMBE->book2D("rocmap_"+hid,twodroctitle,32,0.,32.,yROCbins[DBlayer-1],1.5,1.5+float(yROCbins[DBlayer-1]/2));
172  mePixRocsLay_->setAxisTitle("ROCs per Module",1);
173  mePixRocsLay_->setAxisTitle("ROCs per 1/2 Ladder",2);
174  meZeroOccRocsLay_ = theDMBE->book2D("zeroOccROC_map",twodzeroOccroctitle+hid,32,0.,32.,yROCbins[DBlayer-1],1.5,1.5+float(yROCbins[DBlayer-1]/2));
175  meZeroOccRocsLay_->setAxisTitle("ROCs per Module",1);
176  meZeroOccRocsLay_->setAxisTitle("ROCs per 1/2 Ladder",2);
177  }
178  if(!twoD && !additInfo){
179  // projections of 2D hit map
180  mePixDigisLay_px_ = theDMBE->book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
181  if(isHalfModule){
182  mePixDigisLay_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,2*nbiny,0.,float(2*nrows_));
183  }
184  else{
185  mePixDigisLay_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
186  }
187  mePixDigisLay_px_->setAxisTitle("Columns",1);
188  mePixDigisLay_py_->setAxisTitle("Rows",1);
189  }
190  }
191  }
192  if(type==3 && barrel){
193  uint32_t DBmodule = PixelBarrelName(DetId(id_)).moduleName();
194  char smodule[80]; sprintf(smodule,"Ring_%i",DBmodule);
195  hid = src.label() + "_" + smodule;
196  // Number of digis
197  meNDigisPhi_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
198  meNDigisPhi_->setAxisTitle("Number of digis",1);
199  // Charge in ADC counts
200  meADCPhi_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
201  meADCPhi_->setAxisTitle("ADC counts",1);
202  if(!reducedSet)
203  {
204  if(twoD){
205 
206  // 2D hit map
207  if(isHalfModule){
208  mePixDigisPhi_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),2*nbiny,0.,float(2*nrows_));
209  }
210  else {
211  mePixDigisPhi_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
212  }
213  mePixDigisPhi_->setAxisTitle("Columns",1);
214  mePixDigisPhi_->setAxisTitle("Rows",2);
215  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
216  }
217  else{
218  // projections of 2D hit map
219  mePixDigisPhi_px_ = theDMBE->book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
220  if(isHalfModule){
221  mePixDigisPhi_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,2*nbiny,0.,float(2*nrows_));
222  }
223  else{
224  mePixDigisPhi_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
225  }
226  mePixDigisPhi_px_->setAxisTitle("Columns",1);
227  mePixDigisPhi_py_->setAxisTitle("Rows",1);
228  }
229  }
230  }
231  if(type==4 && endcap){
232  uint32_t blade= PixelEndcapName(DetId(id_)).bladeName();
233 
234  char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
235  hid = src.label() + "_" + sblade;
236  // Number of digis
237  meNDigisBlade_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
238  meNDigisBlade_->setAxisTitle("Number of digis",1);
239  // Charge in ADC counts
240  meADCBlade_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
241  meADCBlade_->setAxisTitle("ADC counts",1);
242  }
243  if(type==5 && endcap){
244  uint32_t disk = PixelEndcapName(DetId(id_)).diskName();
245 
246  char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
247  hid = src.label() + "_" + sdisk;
248  if(!additInfo){
249  // Number of digis
250  meNDigisDisk_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
251  meNDigisDisk_->setAxisTitle("Number of digis",1);
252  // Charge in ADC counts
253  meADCDisk_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
254  meADCDisk_->setAxisTitle("ADC counts",1);
255  }
256  if(additInfo){
257  mePixDigisDisk_ = theDMBE->book2D("hitmap_"+hid,twodtitle,260,0.,260.,160,0.,160.);
258  mePixDigisDisk_->setAxisTitle("Columns",1);
259  mePixDigisDisk_->setAxisTitle("Rows",2);
260  //ROC information in disks
261  mePixRocsDisk_ = theDMBE->book2D("rocmap_"+hid,twodroctitle,26,0.,26.,24,1.,13.);
262  mePixRocsDisk_ ->setAxisTitle("ROCs per Module (2 Panels)",1);
263  mePixRocsDisk_ ->setAxisTitle("Blade Number",2);
264  meZeroOccRocsDisk_ = theDMBE->book2D("zeroOccROC_map",twodzeroOccroctitle+hid,26,0.,26.,24,1.,13.);
265  meZeroOccRocsDisk_ ->setAxisTitle("Zero-Occupancy ROCs per Module (2 Panels)",1);
266  meZeroOccRocsDisk_ ->setAxisTitle("Blade Number",2);
267  }
268  }
269  if(type==6 && endcap){
270  uint32_t panel= PixelEndcapName(DetId(id_)).pannelName();
272  char slab[80]; sprintf(slab, "Panel_%i_Ring_%i",panel, module);
273  hid = src.label() + "_" + slab;
274  // Number of digis
275  meNDigisRing_ = theDMBE->book1D("ndigis_"+hid,"Number of Digis",25,0.,25.);
276  meNDigisRing_->setAxisTitle("Number of digis",1);
277  // Charge in ADC counts
278  meADCRing_ = theDMBE->book1D("adc_" + hid,"Digi charge",128,0.,256.);
279  meADCRing_->setAxisTitle("ADC counts",1);
280  if(!reducedSet)
281  {
282  if(twoD){
283  // 2D hit map
284  mePixDigisRing_ = theDMBE->book2D("hitmap_"+hid,twodtitle,nbinx,0.,float(ncols_),nbiny,0.,float(nrows_));
285  mePixDigisRing_->setAxisTitle("Columns",1);
286  mePixDigisRing_->setAxisTitle("Rows",2);
287  //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
288  }
289  else{
290  // projections of 2D hit map
291  mePixDigisRing_px_ = theDMBE->book1D("hitmap_"+hid+"_px",pxtitle,nbinx,0.,float(ncols_));
292  mePixDigisRing_py_ = theDMBE->book1D("hitmap_"+hid+"_py",pytitle,nbiny,0.,float(nrows_));
293  mePixDigisRing_px_->setAxisTitle("Columns",1);
294  mePixDigisRing_py_->setAxisTitle("Rows",1);
295  }
296  }
297  }
298 }
299 
300 
301 //
302 // Fill histograms
303 //
305  bool ladon, bool layon, bool phion,
306  bool bladeon, bool diskon, bool ringon,
307  bool twoD, bool reducedSet, bool twoDimModOn, bool twoDimOnlyLayDisk,
308  int &nDigisA, int &nDigisB) {
309  bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
310  bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
311  bool isHalfModule = false;
312  uint32_t DBladder = 0;
313  if(barrel){
314  isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule();
315  DBladder = PixelBarrelName(DetId(id_)).ladderName();
316  }
317 
318  // Get DQM interface
319  DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
320  //std::cout<<"id_ = "<<id_<<" , dmbe="<<theDMBE->pwd()<<std::endl;
321  //std::cout<<"********************"<<std::endl;
322  edm::DetSetVector<PixelDigi>::const_iterator isearch = input.find(id_); // search digis of detid
323 
324  unsigned int numberOfDigisMod = 0;
325  int numberOfDigis[8]; for(int i=0; i!=8; i++) numberOfDigis[i]=0;
326  nDigisA=0; nDigisB=0;
327  if( isearch != input.end() ) { // Not an empty iterator
328 
329  // Look at digis now
331  for(di = isearch->data.begin(); di != isearch->data.end(); di++) {
332  int adc = di->adc(); // charge
333  int col = di->column(); // column
334  int row = di->row(); // row
335  numberOfDigisMod++;
337  int DBlayer = PixelBarrelName(DetId(id_)).layerName();
338  int DBmodule = PixelBarrelName(DetId(id_)).moduleName();
339  if(barrel){
340  if(isHalfModule){
341  if(DBshell==PixelBarrelName::pI||DBshell==PixelBarrelName::pO){
342  numberOfDigis[0]++; nDigisA++;
343  if(DBlayer==1) numberOfDigis[2]++;
344  if(DBlayer==2) numberOfDigis[3]++;
345  if(DBlayer==3) numberOfDigis[4]++;
346  }
347  if(DBshell==PixelBarrelName::mI||DBshell==PixelBarrelName::mO){
348  numberOfDigis[1]++; nDigisB++;
349  if(DBlayer==1) numberOfDigis[5]++;
350  if(DBlayer==2) numberOfDigis[6]++;
351  if(DBlayer==3) numberOfDigis[7]++;
352  }
353  }else{
354  if(row<80){
355  numberOfDigis[0]++; nDigisA++;
356  if(DBlayer==1) numberOfDigis[2]++;
357  if(DBlayer==2) numberOfDigis[3]++;
358  if(DBlayer==3) numberOfDigis[4]++;
359  }else{
360  numberOfDigis[1]++; nDigisB++;
361  if(DBlayer==1) numberOfDigis[5]++;
362  if(DBlayer==2) numberOfDigis[6]++;
363  if(DBlayer==3) numberOfDigis[7]++;
364  }
365  }
366  }
367  if(modon){
368  if(!reducedSet){
369  if(twoD) {
370  if(twoDimModOn) (mePixDigis_)->Fill((float)col,(float)row);
371  //std::cout << "Col: " << col << ", Row: " << row << ", for DBlayer " << DBlayer << " and isladder " << DBladder << " and module " << PixelBarrelName(DetId(id_)).moduleName() << " and side is " << DBshell << std::endl;
372  //std::cout<<"id_ = "<<id_<<" , dmbe="<<theDMBE->pwd()<<std::endl;
373  }
374  else {
375  (mePixDigis_px_)->Fill((float)col);
376  (mePixDigis_py_)->Fill((float)row);
377  }
378  }
379  (meADC_)->Fill((float)adc);
380  }
381  if(ladon && barrel){
382  (meADCLad_)->Fill((float)adc);
383  if(!reducedSet){
384  if(twoD) (mePixDigisLad_)->Fill((float)col,(float)row);
385  else {
386  (mePixDigisLad_px_)->Fill((float)col);
387  (mePixDigisLad_py_)->Fill((float)row);
388  }
389  }
390  }
391  if((layon || twoDimOnlyLayDisk) && barrel){
392  if(!twoDimOnlyLayDisk) (meADCLay_)->Fill((float)adc);
393  if(!reducedSet){
394  if((layon && twoD) || twoDimOnlyLayDisk){
395  //ROC histos...
396  float rocx = (float)col/52. + 8.0*float(DBmodule-1);
397  float rocy = (float)row/160.+float(DBladder);
398  //Shift 1st ladder (half modules) up by 1 bin
399  if(DBladder==1) rocy = rocy + 0.5;
400  mePixRocsLay_->Fill(rocx,rocy);
401  //Copying full 1/2 module to empty 1/2 module...
402  //if(isHalfModule) mePixRocsLay_->Fill(rocx,rocy+0.5);
403  //end of ROC filling...
404 
405  if(isHalfModule && DBladder==1){
406  (mePixDigisLay_)->Fill((float)col,(float)row+80);
407  }
408  else (mePixDigisLay_)->Fill((float)col,(float)row);
409  }
410  if((layon && !twoD) && !twoDimOnlyLayDisk){
411  (mePixDigisLay_px_)->Fill((float)col);
412  if(isHalfModule && DBladder==1) {
413  (mePixDigisLay_py_)->Fill((float)row+80);
414  }
415  else (mePixDigisLay_py_)->Fill((float)row);
416  }
417  }
418  }
419  if(phion && barrel){
420  (meADCPhi_)->Fill((float)adc);
421  if(!reducedSet)
422  {
423  if(twoD){
424  if(isHalfModule && DBladder==1){
425  (mePixDigisPhi_)->Fill((float)col,(float)row+80);
426  }
427  else (mePixDigisPhi_)->Fill((float)col,(float)row);
428  }
429  else {
430  (mePixDigisPhi_px_)->Fill((float)col);
431  if(isHalfModule && DBladder==1) {
432  (mePixDigisPhi_py_)->Fill((float)row+80);
433  }
434  else (mePixDigisPhi_py_)->Fill((float)row);
435  }
436  }
437  }
438  if(bladeon && endcap){
439  (meADCBlade_)->Fill((float)adc);
440  }
441 
442  if((diskon || twoDimOnlyLayDisk) && endcap){
443  if(!twoDimOnlyLayDisk) (meADCDisk_)->Fill((float)adc);
444  if(twoDimOnlyLayDisk){
445  (mePixDigisDisk_)->Fill((float)col,(float)row);
446  //ROC monitoring
447  int DBpanel= PixelEndcapName(DetId(id_)).pannelName();
448  int DBblade= PixelEndcapName(DetId(id_)).bladeName();
449  float offx = 0.;
450  //This crazy offset takes into account the roc and module fpix configuration
451  for (int i = DBpanel; i < DBmodule; ++i) {offx = offx + float(5+DBpanel-i);}
452  float rocx = (float)col/52. + offx + 14.0*float(DBpanel-1);
453  float rocy = (float)row/160.+float(DBblade);
454  mePixRocsDisk_->Fill(rocx,rocy);
455  //Now handle the 1/2 module cases by cloning those bins and filling...
456  //if (DBpanel==1 && (DBmodule==1||DBmodule==4)){
457  //rocy = rocy + 0.5;
458  //mePixRocsDisk_->Fill(rocx,rocy);}
459  //end ROC monitoring
460  }
461  }
462  if(ringon && endcap){
463  (meADCRing_)->Fill((float)adc);
464  if(!reducedSet)
465  {
466  if(twoD) (mePixDigisRing_)->Fill((float)col,(float)row);
467  else {
468  (mePixDigisRing_px_)->Fill((float)col);
469  (mePixDigisRing_py_)->Fill((float)row);
470  }
471  }
472  }
473  }
474  if(modon) (meNDigis_)->Fill((float)numberOfDigisMod);
475  if(ladon && barrel) (meNDigisLad_)->Fill((float)numberOfDigisMod);
476  if(layon && barrel && !twoDimOnlyLayDisk) (meNDigisLay_)->Fill((float)numberOfDigisMod);
477  if(phion && barrel) (meNDigisPhi_)->Fill((float)numberOfDigisMod);
478  if(bladeon && endcap) (meNDigisBlade_)->Fill((float)numberOfDigisMod);
479  if(diskon && endcap && !twoDimOnlyLayDisk) (meNDigisDisk_)->Fill((float)numberOfDigisMod);
480  if(ringon && endcap) (meNDigisRing_)->Fill((float)numberOfDigisMod);
481  if(barrel){
482  MonitorElement* me=theDMBE->get("Pixel/Barrel/ALLMODS_ndigisCOMB_Barrel");
483  if(me) me->Fill((float)numberOfDigisMod);
484  me=theDMBE->get("Pixel/Barrel/ALLMODS_ndigisCHAN_Barrel");
485  if(me){ if(numberOfDigis[0]>0) me->Fill((float)numberOfDigis[0]); if(numberOfDigis[1]>0) me->Fill((float)numberOfDigis[1]); }
486  me=theDMBE->get("Pixel/Barrel/ALLMODS_ndigisCHAN_BarrelL1");
487  if(me){ if(numberOfDigis[2]>0) me->Fill((float)numberOfDigis[2]); }
488  me=theDMBE->get("Pixel/Barrel/ALLMODS_ndigisCHAN_BarrelL2");
489  if(me){ if(numberOfDigis[3]>0) me->Fill((float)numberOfDigis[3]); }
490  me=theDMBE->get("Pixel/Barrel/ALLMODS_ndigisCHAN_BarrelL3");
491  if(me){ if(numberOfDigis[4]>0) me->Fill((float)numberOfDigis[4]); }
492  }else if(endcap){
493  MonitorElement* me=theDMBE->get("Pixel/Endcap/ALLMODS_ndigisCOMB_Endcap");
494  if(me) me->Fill((float)numberOfDigisMod);
495  }
496  }
497 
498  //std::cout<<"numberOfDigis for this module: "<<numberOfDigis<<std::endl;
499  return numberOfDigisMod;
500 }
int adc(sample_type sample)
get the ADC sample (12 bits)
int plaquetteName() const
plaquetteId (in pannel)
type
Definition: HCALResponse.h:22
T getParameter(std::string const &) const
int fill(const edm::DetSetVector< PixelDigi > &input, 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)
Fill histograms.
int i
Definition: DBlmapReader.cc:9
MonitorElement * meNDigisRing_
MonitorElement * mePixDigisLad_py_
MonitorElement * mePixDigis_
MonitorElement * mePixDigisPhi_py_
iterator find(det_id_type id)
Definition: DetSetVector.h:285
int moduleName() const
module id (index in z)
MonitorElement * meADCRing_
MonitorElement * meNDigis_
MonitorElement * mePixDigisRing_px_
MonitorElement * mePixDigis_py_
~SiPixelDigiModule()
Destructor.
MonitorElement * meADC_
MonitorElement * mePixDigisDisk_
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 * meADCLad_
MonitorElement * meADCBlade_
MonitorElement * meNDigisLay_
SiPixelDigiModule()
Default constructor.
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1468
MonitorElement * mePixDigisLad_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:356
MonitorElement * meADCLay_
std::string setHistoId(std::string variable, uint32_t &rawId)
Set Histogram Id.
Definition: DetId.h:20
int ladderName() const
ladder id (index in phi)
MonitorElement * mePixDigis_px_
int layerName() const
layer id
MonitorElement * mePixDigisLad_px_
Shell shell() const
MonitorElement * mePixDigisLay_
std::string const & label() const
Definition: InputTag.h:25
MonitorElement * meNDigisPhi_
int pannelName() const
pannel id
MonitorElement * meZeroOccRocsDisk_
MonitorElement * mePixDigisRing_
MonitorElement * meADCPhi_
int diskName() const
disk id
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:341
MonitorElement * mePixDigisPhi_
collection_type::const_iterator const_iterator
Definition: DetSet.h:34
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
MonitorElement * meNDigisDisk_
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Definition: vlib.h:209
MonitorElement * meNDigisLad_
MonitorElement * mePixDigisLay_px_
void book(const edm::ParameterSet &iConfig, int type=0, bool twoD=true, bool hiRes=false, bool reducedSet=false, bool additInfo=false)
Book histograms.