Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:22

0001 #include "DQM/SiPixelMonitorDigi/interface/SiPixelDigiModule.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "DQM/SiPixelCommon/interface/SiPixelHistogramId.h"
0004 /// Framework
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 // STL
0007 #include <vector>
0008 #include <memory>
0009 #include <string>
0010 #include <iostream>
0011 #include <cstdlib>
0012 #include <sstream>
0013 #include <cstdio>
0014 
0015 // Data Formats
0016 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0017 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0018 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0019 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0020 #include "DataFormats/DetId/interface/DetId.h"
0021 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0022 
0023 //
0024 // Constructors
0025 //
0026 SiPixelDigiModule::SiPixelDigiModule() : id_(0), ncols_(416), nrows_(160) {}
0027 ///
0028 SiPixelDigiModule::SiPixelDigiModule(const uint32_t& id) : id_(id), ncols_(416), nrows_(160) {}
0029 ///
0030 SiPixelDigiModule::SiPixelDigiModule(const uint32_t& id, const int& ncols, const int& nrows)
0031     : id_(id), ncols_(ncols), nrows_(nrows) {}
0032 //
0033 // Destructor
0034 //
0035 SiPixelDigiModule::~SiPixelDigiModule() {}
0036 //
0037 // Book histograms
0038 //
0039 void SiPixelDigiModule::book(const edm::ParameterSet& iConfig,
0040                              const TrackerTopology* pTT,
0041                              DQMStore::IBooker& iBooker,
0042                              int type,
0043                              bool twoD,
0044                              bool hiRes,
0045                              bool reducedSet,
0046                              bool additInfo,
0047                              bool isUpgrade) {
0048   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0049   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0050   bool isHalfModule = false;
0051   if (barrel) {
0052     isHalfModule = PixelBarrelName(DetId(id_), pTT, isUpgrade).isHalfModule();
0053   }
0054 
0055   std::string hid;
0056   // Get collection name and instantiate Histo Id builder
0057   edm::InputTag src = iConfig.getParameter<edm::InputTag>("src");
0058 
0059   int nbinx = ncols_ / 2, nbiny = nrows_ / 2;
0060   std::string twodtitle = "Number of Digis (1bin=four pixels)";
0061   std::string pxtitle = "Number of Digis (1bin=two columns)";
0062   std::string pytitle = "Number of Digis (1bin=two rows)";
0063   std::string twodroctitle = "ROC Occupancy (1bin=one ROC)";
0064   std::string twodzeroOccroctitle = "Zero Occupancy ROC Map (1bin=one ROC) for ";
0065   if (hiRes) {
0066     nbinx = ncols_;
0067     nbiny = nrows_;
0068     twodtitle = "Number of Digis (1bin=one pixel)";
0069     pxtitle = "Number of Digis (1bin=one column)";
0070     pytitle = "Number of Digis (1bin=one row)";
0071   }
0072   if (type == 0) {
0073     SiPixelHistogramId* theHistogramId = new SiPixelHistogramId(src.label());
0074     // Number of digis
0075     hid = theHistogramId->setHistoId("ndigis", id_);
0076     meNDigis_ = iBooker.book1D(hid, "Number of Digis", 25, 0., 25.);
0077     meNDigis_->setAxisTitle("Number of digis", 1);
0078     // Charge in ADC counts
0079     hid = theHistogramId->setHistoId("adc", id_);
0080     meADC_ = iBooker.book1D(hid, "Digi charge", 128, 0., 256.);
0081     meADC_->setAxisTitle("ADC counts", 1);
0082     if (!reducedSet) {
0083       if (twoD) {
0084         if (additInfo) {
0085           // 2D hit map
0086           hid = theHistogramId->setHistoId("hitmap", id_);
0087           mePixDigis_ = iBooker.book2D(hid, twodtitle, nbinx, 0., float(ncols_), nbiny, 0., float(nrows_));
0088           mePixDigis_->setAxisTitle("Columns", 1);
0089           mePixDigis_->setAxisTitle("Rows", 2);
0090           //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
0091         }
0092       } else {
0093         // projections of 2D hit map
0094         hid = theHistogramId->setHistoId("hitmap", id_);
0095         mePixDigis_px_ = iBooker.book1D(hid + "_px", pxtitle, nbinx, 0., float(ncols_));
0096         mePixDigis_py_ = iBooker.book1D(hid + "_py", pytitle, nbiny, 0., float(nrows_));
0097         mePixDigis_px_->setAxisTitle("Columns", 1);
0098         mePixDigis_py_->setAxisTitle("Rows", 1);
0099       }
0100     }
0101     delete theHistogramId;
0102   }
0103 
0104   if (type == 1 && barrel) {
0105     uint32_t DBladder;
0106     DBladder = PixelBarrelName(DetId(id_), pTT, isUpgrade).ladderName();
0107     char sladder[80];
0108     sprintf(sladder, "Ladder_%02i", DBladder);
0109     hid = src.label() + "_" + sladder;
0110     if (isHalfModule)
0111       hid += "H";
0112     else
0113       hid += "F";
0114     // Number of digis
0115     meNDigisLad_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0116     meNDigisLad_->setAxisTitle("Number of digis", 1);
0117     // Charge in ADC counts
0118     meADCLad_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0119     meADCLad_->setAxisTitle("ADC counts", 1);
0120     if (!reducedSet) {
0121       if (twoD) {
0122         // 2D hit map
0123         mePixDigisLad_ = iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), nbiny, 0., float(nrows_));
0124         mePixDigisLad_->setAxisTitle("Columns", 1);
0125         mePixDigisLad_->setAxisTitle("Rows", 2);
0126         //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
0127       } else {
0128         // projections of 2D hit map
0129         mePixDigisLad_px_ = iBooker.book1D("hitmap_" + hid + "_px", pxtitle, nbinx, 0., float(ncols_));
0130         mePixDigisLad_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, nbiny, 0., float(nrows_));
0131         mePixDigisLad_px_->setAxisTitle("Columns", 1);
0132         mePixDigisLad_py_->setAxisTitle("Rows", 1);
0133       }
0134     }
0135   }
0136   if (type == 2 && barrel) {
0137     uint32_t DBlayer;
0138     DBlayer = PixelBarrelName(DetId(id_), pTT, isUpgrade).layerName();
0139     char slayer[80];
0140     sprintf(slayer, "Layer_%i", DBlayer);
0141     hid = src.label() + "_" + slayer;
0142     if (!additInfo) {
0143       // Number of digis
0144       meNDigisLay_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0145       meNDigisLay_->setAxisTitle("Number of digis", 1);
0146       // Charge in ADC counts
0147       meADCLay_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0148       meADCLay_->setAxisTitle("ADC counts", 1);
0149     }
0150     if (!reducedSet) {
0151       if (twoD || additInfo) {
0152         // 2D hit map
0153         if (isHalfModule) {
0154           mePixDigisLay_ =
0155               iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), 2 * nbiny, 0., float(2 * nrows_));
0156         } else {
0157           mePixDigisLay_ =
0158               iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), nbiny, 0., float(nrows_));
0159         }
0160         mePixDigisLay_->setAxisTitle("Columns", 1);
0161         mePixDigisLay_->setAxisTitle("Rows", 2);
0162 
0163         //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
0164         int yROCbins[3] = {18, 30, 42};
0165         mePixRocsLay_ = iBooker.book2D("rocmap_" + hid,
0166                                        twodroctitle,
0167                                        32,
0168                                        0.,
0169                                        32.,
0170                                        yROCbins[DBlayer - 1],
0171                                        1.5,
0172                                        1.5 + float(yROCbins[DBlayer - 1] / 2));
0173         mePixRocsLay_->setAxisTitle("ROCs per Module", 1);
0174         mePixRocsLay_->setAxisTitle("ROCs per 1/2 Ladder", 2);
0175         meZeroOccRocsLay_ = iBooker.book2D("zeroOccROC_map",
0176                                            twodzeroOccroctitle + hid,
0177                                            32,
0178                                            0.,
0179                                            32.,
0180                                            yROCbins[DBlayer - 1],
0181                                            1.5,
0182                                            1.5 + float(yROCbins[DBlayer - 1] / 2));
0183         meZeroOccRocsLay_->setAxisTitle("ROCs per Module", 1);
0184         meZeroOccRocsLay_->setAxisTitle("ROCs per 1/2 Ladder", 2);
0185       }
0186       if (!twoD && !additInfo) {
0187         // projections of 2D hit map
0188         mePixDigisLay_px_ = iBooker.book1D("hitmap_" + hid + "_px", pxtitle, nbinx, 0., float(ncols_));
0189         if (isHalfModule) {
0190           mePixDigisLay_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, 2 * nbiny, 0., float(2 * nrows_));
0191         } else {
0192           mePixDigisLay_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, nbiny, 0., float(nrows_));
0193         }
0194         mePixDigisLay_px_->setAxisTitle("Columns", 1);
0195         mePixDigisLay_py_->setAxisTitle("Rows", 1);
0196       }
0197     }
0198   }
0199   if (type == 3 && barrel) {
0200     uint32_t DBmodule;
0201     DBmodule = PixelBarrelName(DetId(id_), pTT, isUpgrade).moduleName();
0202     char smodule[80];
0203     sprintf(smodule, "Ring_%i", DBmodule);
0204     hid = src.label() + "_" + smodule;
0205     // Number of digis
0206     meNDigisPhi_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0207     meNDigisPhi_->setAxisTitle("Number of digis", 1);
0208     // Charge in ADC counts
0209     meADCPhi_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0210     meADCPhi_->setAxisTitle("ADC counts", 1);
0211     if (!reducedSet) {
0212       if (twoD) {
0213         // 2D hit map
0214         if (isHalfModule) {
0215           mePixDigisPhi_ =
0216               iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), 2 * nbiny, 0., float(2 * nrows_));
0217         } else {
0218           mePixDigisPhi_ =
0219               iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), nbiny, 0., float(nrows_));
0220         }
0221         mePixDigisPhi_->setAxisTitle("Columns", 1);
0222         mePixDigisPhi_->setAxisTitle("Rows", 2);
0223         //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
0224       } else {
0225         // projections of 2D hit map
0226         mePixDigisPhi_px_ = iBooker.book1D("hitmap_" + hid + "_px", pxtitle, nbinx, 0., float(ncols_));
0227         if (isHalfModule) {
0228           mePixDigisPhi_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, 2 * nbiny, 0., float(2 * nrows_));
0229         } else {
0230           mePixDigisPhi_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, nbiny, 0., float(nrows_));
0231         }
0232         mePixDigisPhi_px_->setAxisTitle("Columns", 1);
0233         mePixDigisPhi_py_->setAxisTitle("Rows", 1);
0234       }
0235     }
0236   }
0237   if (type == 4 && endcap) {
0238     uint32_t blade;
0239     blade = PixelEndcapName(DetId(id_), pTT, isUpgrade).bladeName();
0240 
0241     char sblade[80];
0242     sprintf(sblade, "Blade_%02i", blade);
0243     hid = src.label() + "_" + sblade;
0244     // Number of digis
0245     meNDigisBlade_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0246     meNDigisBlade_->setAxisTitle("Number of digis", 1);
0247     // Charge in ADC counts
0248     meADCBlade_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0249     meADCBlade_->setAxisTitle("ADC counts", 1);
0250   }
0251   if (type == 5 && endcap) {
0252     uint32_t disk;
0253     disk = PixelEndcapName(DetId(id_), pTT, isUpgrade).diskName();
0254 
0255     char sdisk[80];
0256     sprintf(sdisk, "Disk_%i", disk);
0257     hid = src.label() + "_" + sdisk;
0258     if (!additInfo) {
0259       // Number of digis
0260       meNDigisDisk_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0261       meNDigisDisk_->setAxisTitle("Number of digis", 1);
0262       // Charge in ADC counts
0263       meADCDisk_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0264       meADCDisk_->setAxisTitle("ADC counts", 1);
0265     }
0266     if (additInfo) {
0267       mePixDigisDisk_ = iBooker.book2D("hitmap_" + hid, twodtitle, 260, 0., 260., 160, 0., 160.);
0268       mePixDigisDisk_->setAxisTitle("Columns", 1);
0269       mePixDigisDisk_->setAxisTitle("Rows", 2);
0270       //ROC information in disks
0271       mePixRocsDisk_ = iBooker.book2D("rocmap_" + hid, twodroctitle, 26, 0., 26., 24, 1., 13.);
0272       mePixRocsDisk_->setAxisTitle("ROCs per Module (2 Panels)", 1);
0273       mePixRocsDisk_->setAxisTitle("Blade Number", 2);
0274       meZeroOccRocsDisk_ = iBooker.book2D("zeroOccROC_map", twodzeroOccroctitle + hid, 26, 0., 26., 24, 1., 13.);
0275       meZeroOccRocsDisk_->setAxisTitle("Zero-Occupancy ROCs per Module (2 Panels)", 1);
0276       meZeroOccRocsDisk_->setAxisTitle("Blade Number", 2);
0277     }
0278   }
0279   if (type == 6 && endcap) {
0280     uint32_t panel;
0281     uint32_t module;
0282     panel = PixelEndcapName(DetId(id_), pTT, isUpgrade).pannelName();
0283     module = PixelEndcapName(DetId(id_), pTT, isUpgrade).plaquetteName();
0284 
0285     char slab[80];
0286     sprintf(slab, "Panel_%i_Ring_%i", panel, module);
0287     hid = src.label() + "_" + slab;
0288     // Number of digis
0289     meNDigisRing_ = iBooker.book1D("ndigis_" + hid, "Number of Digis", 25, 0., 25.);
0290     meNDigisRing_->setAxisTitle("Number of digis", 1);
0291     // Charge in ADC counts
0292     meADCRing_ = iBooker.book1D("adc_" + hid, "Digi charge", 128, 0., 256.);
0293     meADCRing_->setAxisTitle("ADC counts", 1);
0294     if (!reducedSet) {
0295       if (twoD) {
0296         // 2D hit map
0297         mePixDigisRing_ =
0298             iBooker.book2D("hitmap_" + hid, twodtitle, nbinx, 0., float(ncols_), nbiny, 0., float(nrows_));
0299         mePixDigisRing_->setAxisTitle("Columns", 1);
0300         mePixDigisRing_->setAxisTitle("Rows", 2);
0301         //std::cout << "During booking: type is "<< type << ", ID is "<< id_ << ", pwd for booking is " << theDMBE->pwd() << ", Plot name: " << hid << std::endl;
0302       } else {
0303         // projections of 2D hit map
0304         mePixDigisRing_px_ = iBooker.book1D("hitmap_" + hid + "_px", pxtitle, nbinx, 0., float(ncols_));
0305         mePixDigisRing_py_ = iBooker.book1D("hitmap_" + hid + "_py", pytitle, nbiny, 0., float(nrows_));
0306         mePixDigisRing_px_->setAxisTitle("Columns", 1);
0307         mePixDigisRing_py_->setAxisTitle("Rows", 1);
0308       }
0309     }
0310   }
0311 }
0312 
0313 //
0314 // Fill histograms
0315 //
0316 int SiPixelDigiModule::fill(const edm::DetSetVector<PixelDigi>& input,
0317                             const TrackerTopology* pTT,
0318                             MonitorElement* combBarrel,
0319                             MonitorElement* chanBarrel,
0320                             std::vector<MonitorElement*>& chanBarrelL,
0321                             MonitorElement* combEndcap,
0322                             bool modon,
0323                             bool ladon,
0324                             bool layon,
0325                             bool phion,
0326                             bool bladeon,
0327                             bool diskon,
0328                             bool ringon,
0329                             bool twoD,
0330                             bool reducedSet,
0331                             bool twoDimModOn,
0332                             bool twoDimOnlyLayDisk,
0333                             int& nDigisA,
0334                             int& nDigisB,
0335                             bool isUpgrade) {
0336   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0337   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0338   bool isHalfModule = false;
0339   uint32_t DBladder = 0;
0340   if (barrel) {
0341     isHalfModule = PixelBarrelName(DetId(id_), pTT, isUpgrade).isHalfModule();
0342     DBladder = PixelBarrelName(DetId(id_), pTT, isUpgrade).ladderName();
0343   }
0344 
0345   edm::DetSetVector<PixelDigi>::const_iterator isearch = input.find(id_);  // search  digis of detid
0346 
0347   unsigned int numberOfDigisMod = 0;
0348   int msize;
0349   if (isUpgrade) {
0350     msize = 10;
0351   } else {
0352     msize = 8;
0353   }
0354   int numberOfDigis[msize];
0355   for (int i = 0; i != msize; i++)
0356     numberOfDigis[i] = 0;
0357   nDigisA = 0;
0358   nDigisB = 0;
0359   if (isearch != input.end()) {  // Not an empty iterator
0360 
0361     // Look at digis now
0362     edm::DetSet<PixelDigi>::const_iterator di;
0363     for (di = isearch->data.begin(); di != isearch->data.end(); di++) {
0364       int adc = di->adc();     // charge
0365       int col = di->column();  // column
0366       int row = di->row();     // row
0367       numberOfDigisMod++;
0368 
0369       int DBlayer = 0;
0370       int DBmodule = 0;
0371 
0372       if (!isUpgrade) {
0373         PixelBarrelName::Shell DBshell = PixelBarrelName(DetId(id_), pTT, isUpgrade).shell();
0374         DBlayer = PixelBarrelName(DetId(id_), pTT, isUpgrade).layerName();
0375         DBmodule = PixelBarrelName(DetId(id_), pTT, isUpgrade).moduleName();
0376         if (barrel) {
0377           if (isHalfModule) {
0378             if (DBshell == PixelBarrelName::pI || DBshell == PixelBarrelName::pO) {
0379               numberOfDigis[0]++;
0380               nDigisA++;
0381               if (DBlayer == 1)
0382                 numberOfDigis[2]++;
0383               if (DBlayer == 2)
0384                 numberOfDigis[3]++;
0385               if (DBlayer == 3)
0386                 numberOfDigis[4]++;
0387             }
0388             if (DBshell == PixelBarrelName::mI || DBshell == PixelBarrelName::mO) {
0389               numberOfDigis[1]++;
0390               nDigisB++;
0391               if (DBlayer == 1)
0392                 numberOfDigis[5]++;
0393               if (DBlayer == 2)
0394                 numberOfDigis[6]++;
0395               if (DBlayer == 3)
0396                 numberOfDigis[7]++;
0397             }
0398           } else {
0399             if (row < 80) {
0400               numberOfDigis[0]++;
0401               nDigisA++;
0402               if (DBlayer == 1)
0403                 numberOfDigis[2]++;
0404               if (DBlayer == 2)
0405                 numberOfDigis[3]++;
0406               if (DBlayer == 3)
0407                 numberOfDigis[4]++;
0408             } else {
0409               numberOfDigis[1]++;
0410               nDigisB++;
0411               if (DBlayer == 1)
0412                 numberOfDigis[5]++;
0413               if (DBlayer == 2)
0414                 numberOfDigis[6]++;
0415               if (DBlayer == 3)
0416                 numberOfDigis[7]++;
0417             }
0418           }
0419         }
0420       } else if (isUpgrade) {
0421         DBlayer = PixelBarrelName(DetId(id_), pTT, isUpgrade).layerName();
0422         DBmodule = PixelBarrelName(DetId(id_), pTT, isUpgrade).moduleName();
0423         if (barrel) {
0424           if (row < 80) {
0425             numberOfDigis[0]++;
0426             nDigisA++;
0427             if (DBlayer == 1)
0428               numberOfDigis[2]++;
0429             if (DBlayer == 2)
0430               numberOfDigis[3]++;
0431             if (DBlayer == 3)
0432               numberOfDigis[4]++;
0433             if (DBlayer == 4)
0434               numberOfDigis[5]++;
0435           } else {
0436             numberOfDigis[1]++;
0437             nDigisB++;
0438             if (DBlayer == 1)
0439               numberOfDigis[6]++;
0440             if (DBlayer == 2)
0441               numberOfDigis[7]++;
0442             if (DBlayer == 3)
0443               numberOfDigis[8]++;
0444             if (DBlayer == 4)
0445               numberOfDigis[9]++;
0446           }
0447         }
0448       }
0449 
0450       if (modon) {
0451         if (!reducedSet) {
0452           if (twoD) {
0453             if (twoDimModOn)
0454               (mePixDigis_)->Fill((float)col, (float)row);
0455           } else {
0456             (mePixDigis_px_)->Fill((float)col);
0457             (mePixDigis_py_)->Fill((float)row);
0458           }
0459         }
0460         (meADC_)->Fill((float)adc);
0461       }
0462       if (ladon && barrel) {
0463         (meADCLad_)->Fill((float)adc);
0464         if (!reducedSet) {
0465           if (twoD)
0466             (mePixDigisLad_)->Fill((float)col, (float)row);
0467           else {
0468             (mePixDigisLad_px_)->Fill((float)col);
0469             (mePixDigisLad_py_)->Fill((float)row);
0470           }
0471         }
0472       }
0473       if ((layon || twoDimOnlyLayDisk) && barrel) {
0474         if (!twoDimOnlyLayDisk)
0475           (meADCLay_)->Fill((float)adc);
0476         if (!reducedSet) {
0477           if ((layon && twoD) || twoDimOnlyLayDisk) {
0478             //ROC histos...
0479             float rocx = (float)col / 52. + 8.0 * float(DBmodule - 1);
0480             float rocy = (float)row / 160. + float(DBladder);
0481             //Shift 1st ladder (half modules) up by 1 bin
0482             if (DBladder == 1)
0483               rocy = rocy + 0.5;
0484             mePixRocsLay_->Fill(rocx, rocy);
0485 
0486             if (isHalfModule && DBladder == 1) {
0487               (mePixDigisLay_)->Fill((float)col, (float)row + 80);
0488             } else
0489               (mePixDigisLay_)->Fill((float)col, (float)row);
0490           }
0491           if ((layon && !twoD) && !twoDimOnlyLayDisk) {
0492             (mePixDigisLay_px_)->Fill((float)col);
0493             if (isHalfModule && DBladder == 1) {
0494               (mePixDigisLay_py_)->Fill((float)row + 80);
0495             } else
0496               (mePixDigisLay_py_)->Fill((float)row);
0497           }
0498         }
0499       }
0500       if (phion && barrel) {
0501         (meADCPhi_)->Fill((float)adc);
0502         if (!reducedSet) {
0503           if (twoD) {
0504             if (isHalfModule && DBladder == 1) {
0505               (mePixDigisPhi_)->Fill((float)col, (float)row + 80);
0506             } else
0507               (mePixDigisPhi_)->Fill((float)col, (float)row);
0508           } else {
0509             (mePixDigisPhi_px_)->Fill((float)col);
0510             if (isHalfModule && DBladder == 1) {
0511               (mePixDigisPhi_py_)->Fill((float)row + 80);
0512             } else
0513               (mePixDigisPhi_py_)->Fill((float)row);
0514           }
0515         }
0516       }
0517       if (bladeon && endcap) {
0518         (meADCBlade_)->Fill((float)adc);
0519       }
0520 
0521       if ((diskon || twoDimOnlyLayDisk) && endcap) {
0522         if (!twoDimOnlyLayDisk)
0523           (meADCDisk_)->Fill((float)adc);
0524         if (twoDimOnlyLayDisk) {
0525           (mePixDigisDisk_)->Fill((float)col, (float)row);
0526           //ROC monitoring
0527           int DBpanel;
0528           int DBblade;
0529           DBpanel = PixelEndcapName(DetId(id_), pTT, isUpgrade).pannelName();
0530           DBblade = PixelEndcapName(DetId(id_), pTT, isUpgrade).bladeName();
0531           float offx = 0.;
0532           //This crazy offset takes into account the roc and module fpix configuration
0533           for (int i = DBpanel; i < DBmodule; ++i) {
0534             offx = offx + float(5 + DBpanel - i);
0535           }
0536           float rocx = (float)col / 52. + offx + 14.0 * float(DBpanel - 1);
0537           float rocy = (float)row / 160. + float(DBblade);
0538           mePixRocsDisk_->Fill(rocx, rocy);
0539         }
0540       }
0541       if (ringon && endcap) {
0542         (meADCRing_)->Fill((float)adc);
0543         if (!reducedSet) {
0544           if (twoD)
0545             (mePixDigisRing_)->Fill((float)col, (float)row);
0546           else {
0547             (mePixDigisRing_px_)->Fill((float)col);
0548             (mePixDigisRing_py_)->Fill((float)row);
0549           }
0550         }
0551       }
0552     }
0553     if (modon)
0554       (meNDigis_)->Fill((float)numberOfDigisMod);
0555     if (ladon && barrel)
0556       (meNDigisLad_)->Fill((float)numberOfDigisMod);
0557     if (layon && barrel && !twoDimOnlyLayDisk)
0558       (meNDigisLay_)->Fill((float)numberOfDigisMod);
0559     if (phion && barrel)
0560       (meNDigisPhi_)->Fill((float)numberOfDigisMod);
0561     if (bladeon && endcap)
0562       (meNDigisBlade_)->Fill((float)numberOfDigisMod);
0563     if (diskon && endcap && !twoDimOnlyLayDisk)
0564       (meNDigisDisk_)->Fill((float)numberOfDigisMod);
0565     if (ringon && endcap)
0566       (meNDigisRing_)->Fill((float)numberOfDigisMod);
0567     if (barrel) {
0568       if (combBarrel)
0569         combBarrel->Fill((float)numberOfDigisMod);
0570       if (chanBarrel) {
0571         if (numberOfDigis[0] > 0)
0572           chanBarrel->Fill((float)numberOfDigis[0]);
0573         if (numberOfDigis[1] > 0)
0574           chanBarrel->Fill((float)numberOfDigis[1]);
0575       }
0576       int j = 2;
0577       for (std::vector<MonitorElement*>::iterator i = chanBarrelL.begin(); i != chanBarrelL.end(); i++) {
0578         if (numberOfDigis[j] > 0)
0579           (*i)->Fill((float)numberOfDigis[j]);
0580         j++;
0581       }
0582     } else if (endcap) {
0583       if (combEndcap)
0584         combEndcap->Fill((float)numberOfDigisMod);
0585     }
0586   }
0587 
0588   //std::cout<<"numberOfDigis for this module: "<<numberOfDigis<<std::endl;
0589   return numberOfDigisMod;
0590 }
0591 
0592 // This was done in the Source file, but is moved to the Module for thread safety reasons. Using ME that is booked here.
0593 void SiPixelDigiModule::resetRocMap() {
0594   if (mePixRocsDisk_)
0595     mePixRocsDisk_->Reset();
0596   if (mePixRocsLay_)
0597     mePixRocsLay_->Reset();
0598 }
0599 
0600 //Moved from source. Gets the zero and low eff ROCs from each module. Called in source for each module.
0601 std::pair<int, int> SiPixelDigiModule::getZeroLoEffROCs() {
0602   int nZeroROC = 0;
0603   int nLoEffROC = 0;
0604   float SF = 1.0;
0605   if (mePixRocsDisk_ && meZeroOccRocsDisk_) {
0606     if (mePixRocsDisk_->getEntries() > 0)
0607       SF = float(mePixRocsDisk_->getNbinsX() * mePixRocsDisk_->getNbinsY() / mePixRocsDisk_->getEntries());
0608     for (int i = 1; i < mePixRocsDisk_->getNbinsX() + 1; ++i) {
0609       for (int j = 1; j < mePixRocsDisk_->getNbinsY() + 1; ++j) {
0610         float localX = float(i) - 0.5;
0611         float localY = float(j) / 2.0 + 0.75;
0612         if (mePixRocsDisk_->getBinContent(i, j) < 1) {
0613           nZeroROC++;
0614           meZeroOccRocsDisk_->Fill(localX, localY);
0615         }
0616         if (mePixRocsDisk_->getBinContent(i, j) * SF < 0.25) {
0617           nLoEffROC++;
0618         }
0619       }
0620     }
0621     return std::pair<int, int>(nZeroROC, nLoEffROC);
0622   }
0623   if (mePixRocsLay_ && meZeroOccRocsLay_) {
0624     if (mePixRocsLay_->getEntries() > 0)
0625       SF = float(mePixRocsLay_->getNbinsX() * mePixRocsLay_->getNbinsY() / mePixRocsLay_->getEntries());
0626     for (int i = 1; i < mePixRocsLay_->getNbinsX() + 1; ++i) {
0627       for (int j = 1; j < mePixRocsLay_->getNbinsY() + 1; ++j) {
0628         float localX = float(i) - 0.5;
0629         float localY = float(j) / 2.0 + 1.25;
0630         if (mePixRocsLay_->getBinContent(i, j) < 8) {
0631           nZeroROC++;
0632           meZeroOccRocsLay_->Fill(localX, localY);
0633         }  //in some regions of pixel there are modules with no HV but enabled ROCs that sometime give a fake hit, so the dead rocs have to be counted to have less than 8 hits in 10 LS
0634         if (mePixRocsLay_->getBinContent(i, j) * SF < 0.25) {
0635           nLoEffROC++;
0636         }
0637       }
0638     }
0639     return std::pair<int, int>(nZeroROC, nLoEffROC);
0640   }
0641   return std::pair<int, int>(0, 0);
0642 }