Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:34

0001 #include "DQM/HcalTasks/interface/PedestalTask.h"
0002 #include "FWCore/Framework/interface/Run.h"
0003 
0004 using namespace hcaldqm;
0005 using namespace hcaldqm::constants;
0006 PedestalTask::PedestalTask(edm::ParameterSet const& ps)
0007     : DQTask(ps), hcalDbServiceToken_(esConsumes<HcalDbService, HcalDbRecord, edm::Transition::BeginRun>()) {
0008   //    tags
0009   _tagQIE11 = ps.getUntrackedParameter<edm::InputTag>("tagHE", edm::InputTag("hcalDigis"));
0010   _tagHO = ps.getUntrackedParameter<edm::InputTag>("tagHO", edm::InputTag("hcalDigis"));
0011   _tagQIE10 = ps.getUntrackedParameter<edm::InputTag>("tagHF", edm::InputTag("hcalDigis"));
0012   _tagTrigger = ps.getUntrackedParameter<edm::InputTag>("tagTrigger", edm::InputTag("tbunpacker"));
0013   _taguMN = ps.getUntrackedParameter<edm::InputTag>("taguMN", edm::InputTag("hcalDigis"));
0014   _tokQIE11 = consumes<QIE11DigiCollection>(_tagQIE11);
0015   _tokHO = consumes<HODigiCollection>(_tagHO);
0016   _tokQIE10 = consumes<QIE10DigiCollection>(_tagQIE10);
0017   _tokTrigger = consumes<HcalTBTriggerData>(_tagTrigger);
0018   _tokuMN = consumes<HcalUMNioDigi>(_taguMN);
0019 
0020   _vflags.resize(2);
0021   _vflags[fMsn] = hcaldqm::flag::Flag("Msn");
0022   _vflags[fBadM] = hcaldqm::flag::Flag("BadM");
0023   //_vflags[fBadR]=hcaldqm::flag::Flag("BadR");
0024 
0025   _thresh_mean = ps.getUntrackedParameter<double>("thresh_mean", 0.25);
0026   _thresh_rms = ps.getUntrackedParameter<double>("thresh_mean", 0.25);
0027   _thresh_badm = ps.getUntrackedParameter<double>("thresh_badm", 0.1);
0028   _thresh_badr = ps.getUntrackedParameter<double>("thresh_badr", 0.1);
0029   _thresh_missing_high = ps.getUntrackedParameter<double>("thresh_missing_high", 0.2);
0030   _thresh_missing_low = ps.getUntrackedParameter<double>("thresh_missing_low", 0.05);
0031 }
0032 
0033 /* virtual */ void PedestalTask::bookHistograms(DQMStore::IBooker& ib, edm::Run const& r, edm::EventSetup const& es) {
0034   if (_ptype == fLocal)
0035     if (r.runAuxiliary().run() == 1)
0036       return;
0037   DQTask::bookHistograms(ib, r, es);
0038 
0039   edm::ESHandle<HcalDbService> dbs = es.getHandle(hcalDbServiceToken_);
0040   _emap = dbs->getHcalMapping();
0041   std::vector<uint32_t> vhashVME;
0042   std::vector<uint32_t> vhashuTCA;
0043   std::vector<uint32_t> vhashC38;
0044   vhashVME.push_back(
0045       HcalElectronicsId(constants::FIBERCH_MIN, constants::FIBER_VME_MIN, SPIGOT_MIN, CRATE_VME_MIN).rawId());
0046   vhashuTCA.push_back(HcalElectronicsId(CRATE_uTCA_MIN, SLOT_uTCA_MIN, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0047   vhashC38.push_back(HcalElectronicsId(38, SLOT_uTCA_MIN, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0048   _filter_VME.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashVME);
0049   _filter_uTCA.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashuTCA);
0050   _filter_C38.initialize(filter::fFilter, hcaldqm::hashfunctions::fCrate, vhashC38);
0051 
0052   //    Containers XXX
0053   _xPedSum1LS.initialize(hcaldqm::hashfunctions::fDChannel);
0054   _xPedSum21LS.initialize(hcaldqm::hashfunctions::fDChannel);
0055   _xPedEntries1LS.initialize(hcaldqm::hashfunctions::fDChannel);
0056   _xPedSumTotal.initialize(hcaldqm::hashfunctions::fDChannel);
0057   _xPedSum2Total.initialize(hcaldqm::hashfunctions::fDChannel);
0058   _xPedEntriesTotal.initialize(hcaldqm::hashfunctions::fDChannel);
0059 
0060 #ifndef HIDE_PEDESTAL_CONDITIONS
0061   _xPedRefMean.initialize(hcaldqm::hashfunctions::fDChannel);
0062   _xPedRefRMS.initialize(hcaldqm::hashfunctions::fDChannel);
0063 #endif
0064 
0065   //    Containers
0066   _cMean1LS_Subdet.initialize(_name,
0067                               "Mean1LS",
0068                               hcaldqm::hashfunctions::fSubdet,
0069                               new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0070                               new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0071                               0);
0072   _cRMS1LS_Subdet.initialize(_name,
0073                              "RMS1LS",
0074                              hcaldqm::hashfunctions::fSubdet,
0075                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0076                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0077                              0);
0078   _cMean1LS_depth.initialize(_name,
0079                              "Mean1LS",
0080                              hcaldqm::hashfunctions::fdepth,
0081                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0082                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0083                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0084                              0);
0085   _cRMS1LS_depth.initialize(_name,
0086                             "RMS1LS",
0087                             hcaldqm::hashfunctions::fdepth,
0088                             new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0089                             new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0090                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0091                             0);
0092   if (_ptype != fOffline) {  // hidefed2crate
0093     _cMean1LS_FEDuTCA.initialize(_name,
0094                                  "Mean1LS",
0095                                  hcaldqm::hashfunctions::fFED,
0096                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0097                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0098                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0099                                  0);
0100     _cRMS1LS_FEDuTCA.initialize(_name,
0101                                 "RMS1LS",
0102                                 hcaldqm::hashfunctions::fFED,
0103                                 new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0104                                 new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0105                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0106                                 0);
0107   }
0108 
0109   _cMeanTotal_Subdet.initialize(_name,
0110                                 "Mean",
0111                                 hcaldqm::hashfunctions::fSubdet,
0112                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0113                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0114                                 0);
0115   _cRMSTotal_Subdet.initialize(_name,
0116                                "RMS",
0117                                hcaldqm::hashfunctions::fSubdet,
0118                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0119                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0120                                0);
0121   _cMeanTotal_depth.initialize(_name,
0122                                "Mean",
0123                                hcaldqm::hashfunctions::fdepth,
0124                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0125                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0126                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0127                                0);
0128   _cRMSTotal_depth.initialize(_name,
0129                               "RMS",
0130                               hcaldqm::hashfunctions::fdepth,
0131                               new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0132                               new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0133                               new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0134                               0);
0135 
0136   _cMeanDBRef1LS_Subdet.initialize(_name,
0137                                    "MeanDBRef1LS",
0138                                    hcaldqm::hashfunctions::fSubdet,
0139                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0140                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0141                                    0);
0142   _cRMSDBRef1LS_Subdet.initialize(_name,
0143                                   "RMSDBRef1LS",
0144                                   hcaldqm::hashfunctions::fSubdet,
0145                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0146                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0147                                   0);
0148   _cMeanDBRef1LS_depth.initialize(_name,
0149                                   "MeanDBRef1LS",
0150                                   hcaldqm::hashfunctions::fdepth,
0151                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0152                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0153                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0154                                   0);
0155   _cRMSDBRef1LS_depth.initialize(_name,
0156                                  "RMSDBRef1LS",
0157                                  hcaldqm::hashfunctions::fdepth,
0158                                  new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0159                                  new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0160                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0161                                  0);
0162 
0163   _cMeanDBRefTotal_Subdet.initialize(_name,
0164                                      "MeanDBRef",
0165                                      hcaldqm::hashfunctions::fSubdet,
0166                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0167                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0168                                      0);
0169   _cRMSDBRefTotal_Subdet.initialize(_name,
0170                                     "RMSDBRef",
0171                                     hcaldqm::hashfunctions::fSubdet,
0172                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0173                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0174                                     0);
0175   _cMeanDBRefTotal_depth.initialize(_name,
0176                                     "MeanDBRef",
0177                                     hcaldqm::hashfunctions::fdepth,
0178                                     new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0179                                     new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0180                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0181                                     0);
0182   _cRMSDBRefTotal_depth.initialize(_name,
0183                                    "RMSDBRef",
0184                                    hcaldqm::hashfunctions::fdepth,
0185                                    new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0186                                    new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0187                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fAroundZero),
0188                                    0);
0189 
0190   _cMissingvsLS_Subdet.initialize(_name,
0191                                   "MissingvsLS",
0192                                   hcaldqm::hashfunctions::fSubdet,
0193                                   new hcaldqm::quantity::LumiSection(_maxLS),
0194                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0195                                   0);
0196   _cOccupancyvsLS_Subdet.initialize(_name,
0197                                     "OccupancyvsLS",
0198                                     hcaldqm::hashfunctions::fSubdet,
0199                                     new hcaldqm::quantity::LumiSection(_maxLS),
0200                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0201                                     0);
0202   _cOccupancyEAvsLS_Subdet.initialize(_name,
0203                                       "OccupancyEAvsLS",
0204                                       hcaldqm::hashfunctions::fSubdet,
0205                                       new hcaldqm::quantity::LumiSection(_maxLS),
0206                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN_to8000),
0207                                       0);
0208   _cNBadMeanvsLS_Subdet.initialize(_name,
0209                                    "NBadMeanvsLS",
0210                                    hcaldqm::hashfunctions::fSubdet,
0211                                    new hcaldqm::quantity::LumiSection(_maxLS),
0212                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0213                                    0);
0214   _cNBadRMSvsLS_Subdet.initialize(_name,
0215                                   "NBadRMSvsLS",
0216                                   hcaldqm::hashfunctions::fSubdet,
0217                                   new hcaldqm::quantity::LumiSection(_maxLS),
0218                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0219                                   0);
0220 
0221   _cMissing1LS_depth.initialize(_name,
0222                                 "Missing1LS",
0223                                 hcaldqm::hashfunctions::fdepth,
0224                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0225                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0226                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0227                                 0);
0228   _cMeanBad1LS_depth.initialize(_name,
0229                                 "MeanBad1LS",
0230                                 hcaldqm::hashfunctions::fdepth,
0231                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0232                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0233                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0234                                 0);
0235   _cRMSBad1LS_depth.initialize(_name,
0236                                "RMSBad1LS",
0237                                hcaldqm::hashfunctions::fdepth,
0238                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0239                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0240                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0241                                0);
0242 
0243   _cMissingTotal_depth.initialize(_name,
0244                                   "Missing",
0245                                   hcaldqm::hashfunctions::fdepth,
0246                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0247                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0248                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0249                                   0);
0250   _cMeanBadTotal_depth.initialize(_name,
0251                                   "MeanBad",
0252                                   hcaldqm::hashfunctions::fdepth,
0253                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0254                                   new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0255                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0256                                   0);
0257   _cRMSBadTotal_depth.initialize(_name,
0258                                  "RMSBad",
0259                                  hcaldqm::hashfunctions::fdepth,
0260                                  new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0261                                  new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0262                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0263                                  0);
0264 
0265   _cADC_SubdetPM.initialize(_name,
0266                             "ADC",
0267                             hcaldqm::hashfunctions::fSubdetPM,
0268                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_256),
0269                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0270                             0);
0271 
0272   if (_ptype != fOffline) {  // hidefed2crate
0273     std::vector<int> vFEDs = hcaldqm::utilities::getFEDList(_emap);
0274     std::vector<int> vFEDsVME = hcaldqm::utilities::getFEDVMEList(_emap);
0275     std::vector<int> vFEDsuTCA = hcaldqm::utilities::getFEDuTCAList(_emap);
0276     for (std::vector<int>::const_iterator it = vFEDsVME.begin(); it != vFEDsVME.end(); ++it)
0277       _vhashFEDs.push_back(
0278           HcalElectronicsId(constants::FIBERCH_MIN, FIBER_VME_MIN, SPIGOT_MIN, (*it) - FED_VME_MIN).rawId());
0279     for (std::vector<int>::const_iterator it = vFEDsuTCA.begin(); it != vFEDsuTCA.end(); ++it) {
0280       std::pair<uint16_t, uint16_t> cspair = utilities::fed2crate(*it);
0281       _vhashFEDs.push_back(HcalElectronicsId(cspair.first, cspair.second, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0282     }
0283     _xNChs.initialize(hcaldqm::hashfunctions::fFED);
0284     _xNMsn1LS.initialize(hcaldqm::hashfunctions::fFED);
0285     _xNBadMean1LS.initialize(hcaldqm::hashfunctions::fFED);
0286     _xNBadRMS1LS.initialize(hcaldqm::hashfunctions::fFED);
0287     _cMeanTotal_FEDuTCA.initialize(_name,
0288                                    "Mean",
0289                                    hcaldqm::hashfunctions::fFED,
0290                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0291                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0292                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_15),
0293                                    0);
0294     _cRMSTotal_FEDuTCA.initialize(_name,
0295                                   "RMS",
0296                                   hcaldqm::hashfunctions::fFED,
0297                                   new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0298                                   new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0299                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0300                                   0);
0301     _cMeanDBRef1LS_FEDuTCA.initialize(_name,
0302                                       "MeanDBRef1LS",
0303                                       hcaldqm::hashfunctions::fFED,
0304                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0305                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0306                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0307                                       0);
0308     _cRMSDBRef1LS_FEDuTCA.initialize(_name,
0309                                      "RMSDBRef1LS",
0310                                      hcaldqm::hashfunctions::fFED,
0311                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0312                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0313                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0314                                      0);
0315     _cMeanDBRefTotal_FEDuTCA.initialize(
0316         _name,
0317         "MeanDBRef",
0318         hcaldqm::hashfunctions::fFED,
0319         new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0320         new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0321         new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0322         0);
0323     _cRMSDBRefTotal_FEDuTCA.initialize(_name,
0324                                        "RMSDBRef",
0325                                        hcaldqm::hashfunctions::fFED,
0326                                        new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0327                                        new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0328                                        new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_5),
0329                                        0);
0330     _cMissing1LS_FEDuTCA.initialize(_name,
0331                                     "Missing1LS",
0332                                     hcaldqm::hashfunctions::fFED,
0333                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0334                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0335                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0336                                     0);
0337     _cMeanBad1LS_FEDuTCA.initialize(_name,
0338                                     "MeanBad1LS",
0339                                     hcaldqm::hashfunctions::fFED,
0340                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0341                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0342                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0343                                     0);
0344     _cRMSBad1LS_FEDuTCA.initialize(_name,
0345                                    "RMSBad1LS",
0346                                    hcaldqm::hashfunctions::fFED,
0347                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0348                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0349                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0350                                    0);
0351     _cMissingTotal_FEDuTCA.initialize(_name,
0352                                       "Missing",
0353                                       hcaldqm::hashfunctions::fFED,
0354                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0355                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0356                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0357                                       0);
0358     _cMeanBadTotal_FEDuTCA.initialize(_name,
0359                                       "MeanBad",
0360                                       hcaldqm::hashfunctions::fFED,
0361                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0362                                       new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0363                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0364                                       0);
0365     _cRMSBadTotal_FEDuTCA.initialize(_name,
0366                                      "RMSBad",
0367                                      hcaldqm::hashfunctions::fFED,
0368                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0369                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0370                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0371                                      0);
0372     _cSummaryvsLS_FED.initialize(_name,
0373                                  "SummaryvsLS",
0374                                  hcaldqm::hashfunctions::fFED,
0375                                  new hcaldqm::quantity::LumiSection(_maxLS),
0376                                  new hcaldqm::quantity::FlagQuantity(_vflags),
0377                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fState),
0378                                  0);
0379     _cSummaryvsLS.initialize(_name,
0380                              "SummaryvsLS",
0381                              new hcaldqm::quantity::LumiSection(_maxLS),
0382                              new hcaldqm::quantity::FEDQuantity(vFEDs),
0383                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fState),
0384                              0);
0385   }
0386 
0387   //    book plots
0388   _cADC_SubdetPM.book(ib, _emap, _subsystem);
0389   _cMean1LS_Subdet.book(ib, _emap, _subsystem);
0390   _cRMS1LS_Subdet.book(ib, _emap, _subsystem);
0391   _cMean1LS_depth.book(ib, _emap, _subsystem);
0392   _cRMS1LS_depth.book(ib, _emap, _subsystem);
0393   _cMeanDBRef1LS_Subdet.book(ib, _emap, _subsystem);
0394   _cRMSDBRef1LS_Subdet.book(ib, _emap, _subsystem);
0395   _cMeanDBRef1LS_depth.book(ib, _emap, _subsystem);
0396   _cRMSDBRef1LS_depth.book(ib, _emap, _subsystem);
0397   _cMissing1LS_depth.book(ib, _emap, _subsystem);
0398   _cMeanBad1LS_depth.book(ib, _emap, _subsystem);
0399   _cRMSBad1LS_depth.book(ib, _emap, _subsystem);
0400 
0401   _cMeanTotal_Subdet.book(ib, _emap, _subsystem);
0402   _cRMSTotal_Subdet.book(ib, _emap, _subsystem);
0403   _cMeanTotal_depth.book(ib, _emap, _subsystem);
0404   _cRMSTotal_depth.book(ib, _emap, _subsystem);
0405   _cMeanDBRefTotal_Subdet.book(ib, _emap, _subsystem);
0406   _cRMSDBRefTotal_Subdet.book(ib, _emap, _subsystem);
0407   _cMeanDBRefTotal_depth.book(ib, _emap, _subsystem);
0408   _cRMSDBRefTotal_depth.book(ib, _emap, _subsystem);
0409   _cMissingTotal_depth.book(ib, _emap, _subsystem);
0410   _cMeanBadTotal_depth.book(ib, _emap, _subsystem);
0411   _cRMSBadTotal_depth.book(ib, _emap, _subsystem);
0412 
0413   if (_ptype != fOffline) {  // hidefed2crate
0414     _cMean1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0415     _cRMS1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0416     _cMeanDBRef1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0417     _cRMSDBRef1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0418     _cMissing1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0419     _cRMSBad1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0420     _cMeanBad1LS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0421 
0422     _cMeanTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0423     _cRMSTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0424     _cMeanDBRefTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0425     _cRMSDBRefTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0426     _cMissingTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0427     _cRMSBadTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0428     _cMeanBadTotal_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0429   }
0430 
0431   _cMissingvsLS_Subdet.book(ib, _emap, _subsystem);
0432   _cOccupancyvsLS_Subdet.book(ib, _emap, _subsystem);
0433   _cOccupancyEAvsLS_Subdet.book(ib, _emap, _subsystem);
0434   _cNBadMeanvsLS_Subdet.book(ib, _emap, _subsystem);
0435   _cNBadRMSvsLS_Subdet.book(ib, _emap, _subsystem);
0436   if (_ptype != fOffline) {  // hidefed2crate
0437     _cSummaryvsLS_FED.book(ib, _emap, _subsystem);
0438     _cSummaryvsLS.book(ib, _subsystem);
0439   }
0440 
0441   //    book compact containers
0442   _xPedSum1LS.book(_emap);
0443   _xPedSum21LS.book(_emap);
0444   _xPedEntries1LS.book(_emap);
0445   _xPedSumTotal.book(_emap);
0446   _xPedSum2Total.book(_emap);
0447   _xPedEntriesTotal.book(_emap);
0448 
0449 #ifndef HIDE_PEDESTAL_CONDITIONS
0450   _xPedRefMean.book(_emap);
0451   _xPedRefRMS.book(_emap);
0452 #endif
0453 
0454   if (_ptype != fOffline) {  // hidefed2crate
0455     _xNChs.book(_emap);
0456     _xNMsn1LS.book(_emap);
0457     _xNBadMean1LS.book(_emap);
0458     _xNBadRMS1LS.book(_emap);
0459   }
0460 
0461   _ehashmap.initialize(_emap, electronicsmap::fD2EHashMap);
0462 
0463   //    load conditions pedestals
0464   std::vector<HcalGenericDetId> dids = _emap->allPrecisionId();
0465   for (std::vector<HcalGenericDetId>::const_iterator it = dids.begin(); it != dids.end(); ++it) {
0466     //  skip if calib or whatever
0467     if (!it->isHcalDetId())
0468       continue;
0469     //  skip Crate 38
0470     if (_filter_C38.filter(HcalElectronicsId(_ehashmap.lookup(*it))))
0471       continue;
0472 #ifndef HIDE_PEDESTAL_CONDITIONS
0473     HcalDetId did = HcalDetId(it->rawId());
0474 
0475     HcalPedestal const* peds = dbs->getPedestal(did);
0476     float const* means = peds->getValues();
0477     float const* rmss = peds->getWidths();
0478     double msum = 0;
0479     double rsum = 0;
0480     for (uint32_t i = 0; i < 4; i++) {
0481       msum += means[i];
0482       rsum += rmss[i];
0483     }
0484     msum /= 4;
0485     rsum /= 4;
0486     _xPedRefMean.set(did, msum);
0487     _xPedRefRMS.set(did, rsum);
0488 #endif
0489   }
0490 }
0491 
0492 /* virtual */ void PedestalTask::_resetMonitors(hcaldqm::UpdateFreq uf) {
0493   DQTask::_resetMonitors(uf);
0494 
0495   switch (uf) {
0496     case hcaldqm::f50LS:
0497       _cADC_SubdetPM.reset();
0498       break;
0499     default:
0500       break;
0501   }
0502 }
0503 
0504 /* virtual */ void PedestalTask::_dump() {
0505   //    reset what's needed
0506 
0507   //    Mean/RMS actual values
0508   _cMean1LS_Subdet.reset();
0509   _cRMS1LS_Subdet.reset();
0510   _cMean1LS_depth.reset();
0511   _cRMS1LS_depth.reset();
0512   if (_ptype != fOffline) {  // hidefed2crate
0513     _cMean1LS_FEDuTCA.reset();
0514     _cRMS1LS_FEDuTCA.reset();
0515   }
0516 
0517   _cMeanTotal_Subdet.reset();
0518   _cRMSTotal_Subdet.reset();
0519   _cMeanTotal_depth.reset();
0520   _cRMSTotal_depth.reset();
0521   if (_ptype != fOffline) {  // hidefed2crate
0522     _cMeanTotal_FEDuTCA.reset();
0523     _cRMSTotal_FEDuTCA.reset();
0524   }
0525 
0526   //    DB Conditions Comparison
0527   _cMeanDBRef1LS_Subdet.reset();
0528   _cMeanDBRef1LS_depth.reset();
0529   _cRMSDBRef1LS_Subdet.reset();
0530   _cRMSDBRef1LS_depth.reset();
0531 
0532   _cMeanDBRefTotal_Subdet.reset();
0533   _cMeanDBRefTotal_depth.reset();
0534   _cRMSDBRefTotal_Subdet.reset();
0535   _cRMSDBRefTotal_depth.reset();
0536 
0537   if (_ptype != fOffline) {  // hidefed2crate
0538     _cMeanDBRef1LS_FEDuTCA.reset();
0539     _cRMSDBRef1LS_FEDuTCA.reset();
0540 
0541     _cMeanDBRefTotal_FEDuTCA.reset();
0542     _cRMSDBRefTotal_FEDuTCA.reset();
0543   }
0544 
0545   //    missing channels
0546   _cMissing1LS_depth.reset();
0547   _cMeanBad1LS_depth.reset();
0548   _cRMSBad1LS_depth.reset();
0549 
0550   _cMissingTotal_depth.reset();
0551   _cMeanBadTotal_depth.reset();
0552   _cRMSBadTotal_depth.reset();
0553 
0554   //    Missing or Bad
0555   if (_ptype != fOffline) {  // hidefed2crate
0556     _cMissing1LS_FEDuTCA.reset();
0557     _cMeanBad1LS_FEDuTCA.reset();
0558     _cRMSBad1LS_FEDuTCA.reset();
0559 
0560     _cMissingTotal_FEDuTCA.reset();
0561     _cMeanBadTotal_FEDuTCA.reset();
0562     _cRMSBadTotal_FEDuTCA.reset();
0563 
0564     //  reset some XXX containers
0565     _xNChs.reset();
0566     _xNMsn1LS.reset();
0567     _xNBadMean1LS.reset();
0568     _xNBadRMS1LS.reset();
0569   }
0570   // - ITERATE OVER ALL TEH CHANNELS
0571   // - FIND THE ONES THAT ARE MISSING
0572   // - FIND THE ONES WITH BAD PEDESTAL MEANs
0573   // - FIND THE ONES WITH BAD PEDESTAL RMSs
0574   std::vector<HcalGenericDetId> dids = _emap->allPrecisionId();
0575   for (std::vector<HcalGenericDetId>::const_iterator it = dids.begin(); it != dids.end(); ++it) {
0576     if (!it->isHcalDetId())
0577       continue;
0578     HcalElectronicsId eid(_ehashmap.lookup(*it));
0579     if (_filter_C38.filter(eid))
0580       continue;
0581 
0582     //  filter out channels with bad quality
0583     if (_xQuality.exists(HcalDetId(*it))) {
0584       HcalChannelStatus cs(it->rawId(), _xQuality.get(HcalDetId(*it)));
0585       if (cs.isBitSet(HcalChannelStatus::HcalCellMask) || cs.isBitSet(HcalChannelStatus::HcalCellDead))
0586         continue;
0587     }
0588 
0589     HcalDetId did = HcalDetId(it->rawId());
0590     double sum1LS = _xPedSum1LS.get(did);
0591 #ifndef HIDE_PEDESTAL_CONDITIONS
0592     double refm = _xPedRefMean.get(did);
0593 #endif
0594     double sum21LS = _xPedSum21LS.get(did);
0595 #ifndef HIDE_PEDESTAL_CONDITIONS
0596     double refr = _xPedRefRMS.get(did);
0597 #endif
0598     double n1LS = _xPedEntries1LS.get(did);
0599 
0600     double sumTotal = _xPedSumTotal.get(did);
0601     double sum2Total = _xPedSum2Total.get(did);
0602     double nTotal = _xPedEntriesTotal.get(did);
0603 
0604     if (_ptype != fOffline) {  // hidefed2crate
0605       _xNChs.get(eid)++;
0606     }
0607     // IF A CHANNEL IS MISSING FOR THIS LS
0608     if (n1LS == 0) {
0609       _cMissing1LS_depth.fill(did);
0610       _cMissingvsLS_Subdet.fill(did, _currentLS);
0611       if (_ptype != fOffline) {  // hidefed2crate
0612         if (!eid.isVMEid())
0613           _cMissing1LS_FEDuTCA.fill(eid);
0614         _xNMsn1LS.get(eid)++;
0615       }
0616       //    ALSO CHECK
0617       //    IF A CHANNEL HAS BEEN MISSING FOR ALL LSs SO FAR
0618       if (nTotal == 0) {
0619         _cMissingTotal_depth.fill(did);
0620         if (_ptype != fOffline) {  // hidefed2crate
0621           if (!eid.isVMEid())
0622             _cMissingTotal_FEDuTCA.fill(eid);
0623         }
0624       }
0625       continue;
0626     }
0627 
0628     //  if not missing, fill the occupancy...
0629     _cOccupancyvsLS_Subdet.fill(did, _currentLS);
0630 
0631     //  compute the means and diffs for this LS
0632     sum1LS /= n1LS;
0633     double rms1LS = sqrt(sum21LS / n1LS - sum1LS * sum1LS);
0634 #ifndef HIDE_PEDESTAL_CONDITIONS
0635     double diffm1LS = sum1LS - refm;
0636     double diffr1LS = rms1LS - refr;
0637 #endif
0638     //  compute the means and diffs for the whole Run
0639     sumTotal /= nTotal;
0640     double rmsTotal = sqrt(sum2Total / nTotal - sumTotal * sumTotal);
0641 #ifndef HIDE_PEDESTAL_CONDITIONS
0642     double diffmTotal = sumTotal - refm;
0643     double diffrTotal = rmsTotal - refr;
0644 #endif
0645     //  FILL ACTUAL MEANs AND RMSs FOR THIS LS
0646     _cMean1LS_Subdet.fill(did, sum1LS);
0647     _cMean1LS_depth.fill(did, sum1LS);
0648     _cRMS1LS_Subdet.fill(did, rms1LS);
0649     _cRMS1LS_depth.fill(did, rms1LS);
0650 
0651     //  FILL THE DIFFERENCES FOR THIS LS
0652 #ifndef HIDE_PEDESTAL_CONDITIONS
0653     _cMeanDBRef1LS_Subdet.fill(did, diffm1LS);
0654     _cMeanDBRef1LS_depth.fill(did, diffm1LS);
0655     _cRMSDBRef1LS_Subdet.fill(did, diffr1LS);
0656     _cRMSDBRef1LS_depth.fill(did, diffr1LS);
0657 #endif
0658     //  FILL ACTUAL MEANs AND RMSs FOR THE WHOLE RUN
0659     _cMeanTotal_Subdet.fill(did, sumTotal);
0660     _cMeanTotal_depth.fill(did, sumTotal);
0661     _cRMSTotal_Subdet.fill(did, rmsTotal);
0662     _cRMSTotal_depth.fill(did, rmsTotal);
0663 
0664     //  FILL THE DIFFERENCES FOR THE WHOLE RUN
0665 #ifndef HIDE_PEDESTAL_CONDITIONS
0666     _cMeanDBRefTotal_Subdet.fill(did, diffmTotal);
0667     _cMeanDBRefTotal_depth.fill(did, diffmTotal);
0668     _cRMSDBRefTotal_Subdet.fill(did, diffrTotal);
0669     _cRMSDBRefTotal_depth.fill(did, diffrTotal);
0670 #endif
0671     //  FOR THIS LS
0672     if (_ptype != fOffline) {  // hidefed2crate
0673       if (!eid.isVMEid()) {
0674         _cMean1LS_FEDuTCA.fill(eid, sum1LS);
0675         _cRMS1LS_FEDuTCA.fill(eid, rms1LS);
0676         _cMeanDBRef1LS_FEDuTCA.fill(eid, diffm1LS);
0677         _cRMSDBRef1LS_FEDuTCA.fill(eid, diffr1LS);
0678       }
0679 
0680       //    FOR THE WHOLE RUN
0681       if (!eid.isVMEid()) {
0682         _cMeanTotal_FEDuTCA.fill(eid, sumTotal);
0683         _cRMSTotal_FEDuTCA.fill(eid, rmsTotal);
0684         _cMeanDBRefTotal_FEDuTCA.fill(eid, diffmTotal);
0685         _cRMSDBRefTotal_FEDuTCA.fill(eid, diffrTotal);
0686       }
0687     }
0688 
0689     //  FOR THE CURRENT LS COMPARE MEANS AND RMSS
0690 #ifndef HIDE_PEDESTAL_CONDITIONS
0691     if (fabs(diffm1LS) > _thresh_mean) {
0692       _cMeanBad1LS_depth.fill(did);
0693       _cNBadMeanvsLS_Subdet.fill(did, _currentLS);
0694       if (_ptype != fOffline) {  // hidefed2crate
0695         if (!eid.isVMEid())
0696           _cMeanBad1LS_FEDuTCA.fill(eid);
0697         _xNBadMean1LS.get(eid)++;
0698       }
0699     }
0700     if (fabs(diffr1LS) > _thresh_rms) {
0701       _cRMSBad1LS_depth.fill(did);
0702       _cNBadRMSvsLS_Subdet.fill(did, _currentLS);
0703       if (_ptype != fOffline) {  // hidefed2crate
0704         if (!eid.isVMEid())
0705           _cRMSBad1LS_FEDuTCA.fill(eid);
0706         _xNBadRMS1LS.get(eid)++;
0707       }
0708     }
0709 
0710     //  FOR THIS RUN
0711     if (fabs(diffmTotal) > _thresh_mean) {
0712       _cMeanBadTotal_depth.fill(did);
0713       if (_ptype != fOffline) {  // hidefed2crate
0714         if (!eid.isVMEid())
0715           _cMeanBadTotal_FEDuTCA.fill(eid);
0716       }
0717     }
0718     if (fabs(diffrTotal) > _thresh_rms) {
0719       _cRMSBadTotal_depth.fill(did);
0720       if (_ptype != fOffline) {  // hidefed2crate
0721         if (!eid.isVMEid())
0722           _cRMSBadTotal_FEDuTCA.fill(eid);
0723       }
0724     }
0725 #endif
0726   }
0727 
0728   //    SET THE FLAGS FOR THIS LS
0729   if (_ptype != fOffline) {  // hidefed2crate
0730     for (std::vector<uint32_t>::const_iterator it = _vhashFEDs.begin(); it != _vhashFEDs.end(); ++it) {
0731       hcaldqm::flag::Flag fSum("PED");
0732       HcalElectronicsId eid = HcalElectronicsId(*it);
0733 
0734       std::vector<uint32_t>::const_iterator jt = std::find(_vcdaqEids.begin(), _vcdaqEids.end(), (*it));
0735       if (jt == _vcdaqEids.end()) {
0736         //  not @cDAQ
0737         for (uint32_t iflag = 0; iflag < _vflags.size(); iflag++)
0738           _cSummaryvsLS_FED.setBinContent(eid, _currentLS, int(iflag), int(hcaldqm::flag::fNCDAQ));
0739         _cSummaryvsLS.setBinContent(eid, _currentLS, int(hcaldqm::flag::fNCDAQ));
0740         continue;
0741       }
0742 
0743       //    @cDAQ
0744       if (hcaldqm::utilities::isFEDHBHE(eid) || hcaldqm::utilities::isFEDHO(eid) || hcaldqm::utilities::isFEDHF(eid)) {
0745         double frmissing = (_xNChs.get(eid) == 0) ? 0 : double(_xNMsn1LS.get(eid)) / double(_xNChs.get(eid));
0746         double frbadm = (_xNChs.get(eid) == 0) ? 0 : _xNBadMean1LS.get(eid) / _xNChs.get(eid);
0747         //double frbadr = _xNBadRMS1LS.get(eid)/_xNChs.get(eid);
0748 
0749         if (frmissing >= _thresh_missing_high)
0750           _vflags[fMsn]._state = hcaldqm::flag::fBAD;
0751         else if (frmissing >= _thresh_missing_low)
0752           _vflags[fMsn]._state = hcaldqm::flag::fPROBLEMATIC;
0753         else
0754           _vflags[fMsn]._state = hcaldqm::flag::fGOOD;
0755         if (frbadm >= _thresh_badm)
0756           _vflags[fBadM]._state = hcaldqm::flag::fBAD;
0757         else
0758           _vflags[fBadM]._state = hcaldqm::flag::fGOOD;
0759         // BadR removed May 9, 2018 - the pedestal RMS isn't stable enough to monitor right now.
0760         //if (frbadr>=_thresh_badr)
0761         //  _vflags[fBadR]._state = hcaldqm::flag::fBAD;
0762         //else
0763         //  _vflags[fBadR]._state = hcaldqm::flag::fGOOD;
0764       }
0765 
0766       int iflag = 0;
0767       for (std::vector<hcaldqm::flag::Flag>::iterator ft = _vflags.begin(); ft != _vflags.end(); ++ft) {
0768         _cSummaryvsLS_FED.setBinContent(eid, _currentLS, iflag, int(ft->_state));
0769         fSum += (*ft);
0770         iflag++;
0771         ft->reset();
0772       }
0773       _cSummaryvsLS.setBinContent(eid, _currentLS, fSum._state);
0774     }
0775   }
0776 
0777   //    reset the pedestal containers instead of writting reset... func
0778   _xPedSum1LS.reset();
0779   _xPedSum21LS.reset();
0780   _xPedEntries1LS.reset();
0781 }
0782 
0783 std::shared_ptr<hcaldqm::Cache> PedestalTask::globalBeginLuminosityBlock(edm::LuminosityBlock const& lb,
0784                                                                          edm::EventSetup const& es) const {
0785   return DQTask::globalBeginLuminosityBlock(lb, es);
0786 }
0787 
0788 /* virtual */ void PedestalTask::dqmEndRun(edm::Run const& r, edm::EventSetup const&) {
0789   if (_ptype == fLocal) {
0790     if (r.runAuxiliary().run() == 1)
0791       return;
0792     else
0793       this->_dump();
0794   } else if (_ptype == fOnline)
0795     return;
0796 }
0797 
0798 /* virtual */ void PedestalTask::globalEndLuminosityBlock(edm::LuminosityBlock const& lb, edm::EventSetup const& es) {
0799   auto lumiCache = luminosityBlockCache(lb.index());
0800   _currentLS = lumiCache->currentLS;
0801   _xQuality.reset();
0802   _xQuality = lumiCache->xQuality;
0803 
0804   if (_ptype == fLocal)
0805     return;
0806   this->_dump();
0807 
0808   DQTask::globalEndLuminosityBlock(lb, es);
0809 }
0810 
0811 /* virtual */ void PedestalTask::_process(edm::Event const& e, edm::EventSetup const& es) {
0812   edm::Handle<HODigiCollection> c_ho;
0813   edm::Handle<QIE10DigiCollection> c_QIE10;
0814   edm::Handle<QIE11DigiCollection> c_QIE11;
0815 
0816   if (!e.getByToken(_tokHO, c_ho))
0817     _logger.dqmthrow("Collection HODigiCollection isn't available" + _tagHO.label() + " " + _tagHO.instance());
0818   if (!e.getByToken(_tokQIE10, c_QIE10))
0819     _logger.dqmthrow("Collection QIE10DigiCollection isn't available" + _tagQIE10.label() + " " + _tagQIE10.instance());
0820   if (!e.getByToken(_tokQIE11, c_QIE11))
0821     _logger.dqmthrow("Collection QIE11DigiCollection isn't available" + _tagQIE11.label() + " " + _tagQIE11.instance());
0822 
0823   auto lumiCache = luminosityBlockCache(e.getLuminosityBlock().index());
0824   _currentLS = lumiCache->currentLS;
0825 
0826   int nHB(0), nHE(0), nHO(0), nHF(0);
0827   for (QIE11DigiCollection::const_iterator it = c_QIE11->begin(); it != c_QIE11->end(); ++it) {
0828     const QIE11DataFrame digi = static_cast<const QIE11DataFrame>(*it);
0829     HcalDetId const& did = digi.detid();
0830     // Require barrel or endcap. As of 2017, some calibration channels are ending up in this collection.
0831     if ((did.subdet() != HcalEndcap) && (did.subdet() != HcalBarrel)) {
0832       continue;
0833     }
0834     int digiSizeToUse = floor(digi.samples() / constants::CAPS_NUM) * constants::CAPS_NUM - 1;
0835     did.subdet() == HcalBarrel ? nHB++ : nHE++;
0836 
0837     for (int i = 0; i < digiSizeToUse; i++) {
0838       _cADC_SubdetPM.fill(did, digi[i].adc());
0839 
0840       _xPedSum1LS.get(did) += digi[i].adc();
0841       _xPedSum21LS.get(did) += digi[i].adc() * digi[i].adc();
0842       _xPedEntries1LS.get(did)++;
0843 
0844       _xPedSumTotal.get(did) += digi[i].adc();
0845       _xPedSum2Total.get(did) += digi[i].adc() * digi[i].adc();
0846       _xPedEntriesTotal.get(did)++;
0847     }
0848   }
0849 
0850   _cOccupancyEAvsLS_Subdet.fill(HcalDetId(HcalBarrel, 1, 1, 1), _currentLS, nHB);
0851   _cOccupancyEAvsLS_Subdet.fill(HcalDetId(HcalEndcap, 1, 1, 1), _currentLS, nHE);
0852 
0853   for (HODigiCollection::const_iterator it = c_ho->begin(); it != c_ho->end(); ++it) {
0854     const HODataFrame digi = (const HODataFrame)(*it);
0855     HcalDetId did = digi.id();
0856     int digiSizeToUse = floor(digi.size() / constants::CAPS_NUM) * constants::CAPS_NUM - 1;
0857     nHO++;
0858     for (int i = 0; i < digiSizeToUse; i++) {
0859       _cADC_SubdetPM.fill(did, it->sample(i).adc());
0860 
0861       _xPedSum1LS.get(did) += it->sample(i).adc();
0862       _xPedSum21LS.get(did) += it->sample(i).adc() * it->sample(i).adc();
0863       _xPedEntries1LS.get(did)++;
0864 
0865       _xPedSumTotal.get(did) += it->sample(i).adc();
0866       _xPedSum2Total.get(did) += it->sample(i).adc() * it->sample(i).adc();
0867       _xPedEntriesTotal.get(did)++;
0868     }
0869   }
0870   _cOccupancyEAvsLS_Subdet.fill(HcalDetId(HcalOuter, 1, 1, 1), _currentLS, nHO);
0871 
0872   for (QIE10DigiCollection::const_iterator it = c_QIE10->begin(); it != c_QIE10->end(); ++it) {
0873     const QIE10DataFrame digi = static_cast<const QIE10DataFrame>(*it);
0874     HcalDetId did = digi.detid();
0875     if (did.subdet() != HcalForward) {
0876       continue;
0877     }
0878     // HF has 3 samples in global, so impossible to make divisible by 4
0879     int digiSizeToUse =
0880         (digi.samples() >= 4 ? floor(digi.samples() / constants::CAPS_NUM) * constants::CAPS_NUM - 1 : digi.samples());
0881     nHF++;
0882     for (int i = 0; i < digiSizeToUse; i++) {
0883       _cADC_SubdetPM.fill(did, digi[i].adc());
0884 
0885       _xPedSum1LS.get(did) += digi[i].adc();
0886       _xPedSum21LS.get(did) += digi[i].adc() * digi[i].adc();
0887       _xPedEntries1LS.get(did)++;
0888 
0889       _xPedSumTotal.get(did) += digi[i].adc();
0890       _xPedSum2Total.get(did) += digi[i].adc() * digi[i].adc();
0891       _xPedEntriesTotal.get(did)++;
0892     }
0893   }
0894   _cOccupancyEAvsLS_Subdet.fill(HcalDetId(HcalForward, 1, 1, 1), _currentLS, nHF);
0895 }
0896 
0897 /* virtual */ bool PedestalTask::_isApplicable(edm::Event const& e) {
0898   if (_ptype == fOnline) {
0899     edm::Handle<HcalUMNioDigi> cumn;
0900     if (!e.getByToken(_tokuMN, cumn))
0901       return false;
0902 
0903     //  for online just check the event type not the user Word
0904     uint8_t eventType = cumn->eventType();
0905     if (eventType == constants::EVENTTYPE_PEDESTAL)
0906       return true;
0907   } else {
0908     //  local
0909     edm::Handle<HcalTBTriggerData> ctrigger;
0910     if (!e.getByToken(_tokTrigger, ctrigger))
0911       _logger.dqmthrow("Collection HcalTBTriggerData isn't available" + _tagTrigger.label() + " " +
0912                        _tagTrigger.instance());
0913     return ctrigger->wasSpillIgnorantPedestalTrigger();
0914   }
0915 
0916   return false;
0917 }
0918 
0919 DEFINE_FWK_MODULE(PedestalTask);