Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQM/HcalTasks/interface/LaserTask.h"
0002 
0003 using namespace hcaldqm;
0004 using namespace hcaldqm::constants;
0005 using namespace hcaldqm::filter;
0006 LaserTask::LaserTask(edm::ParameterSet const& ps)
0007     : DQTask(ps), hcalDbServiceToken_(esConsumes<HcalDbService, HcalDbRecord, edm::Transition::BeginRun>()) {
0008   _nevents = ps.getUntrackedParameter<int>("nevents", 2000);
0009 
0010   //    tags
0011   _tagQIE11 = ps.getUntrackedParameter<edm::InputTag>("tagHE", edm::InputTag("hcalDigis"));
0012   _tagHO = ps.getUntrackedParameter<edm::InputTag>("tagHO", edm::InputTag("hcalDigis"));
0013   _tagQIE10 = ps.getUntrackedParameter<edm::InputTag>("tagHF", edm::InputTag("hcalDigis"));
0014   _taguMN = ps.getUntrackedParameter<edm::InputTag>("taguMN", edm::InputTag("hcalDigis"));
0015   _tagLaserMon = ps.getUntrackedParameter<edm::InputTag>("tagLaserMon", edm::InputTag("LASERMON"));
0016   _tokQIE11 = consumes<QIE11DigiCollection>(_tagQIE11);
0017   _tokHO = consumes<HODigiCollection>(_tagHO);
0018   _tokQIE10 = consumes<QIE10DigiCollection>(_tagQIE10);
0019   _tokuMN = consumes<HcalUMNioDigi>(_taguMN);
0020   _tokLaserMon = consumes<QIE10DigiCollection>(_tagLaserMon);
0021 
0022   _vflags.resize(nLaserFlag);
0023   _vflags[fBadTiming] = hcaldqm::flag::Flag("BadTiming");
0024   _vflags[fMissingLaserMon] = hcaldqm::flag::Flag("MissingLaserMon");
0025 
0026   //    constants
0027   _lowHBHE = ps.getUntrackedParameter<double>("lowHBHE", 50);
0028   _lowHO = ps.getUntrackedParameter<double>("lowHO", 100);
0029   _lowHF = ps.getUntrackedParameter<double>("lowHF", 50);
0030   _laserType = (uint32_t)ps.getUntrackedParameter<uint32_t>("laserType");
0031 
0032   // Laser mon digi ordering list
0033   _vLaserMonIPhi = ps.getUntrackedParameter<std::vector<int>>("vLaserMonIPhi");
0034   _laserMonIEta = ps.getUntrackedParameter<int>("laserMonIEta");
0035   _laserMonCBox = ps.getUntrackedParameter<int>("laserMonCBox");
0036   _laserMonDigiOverlap = ps.getUntrackedParameter<int>("laserMonDigiOverlap");
0037   _laserMonTS0 = ps.getUntrackedParameter<int>("laserMonTS0");
0038   _laserMonThreshold = ps.getUntrackedParameter<double>("laserMonThreshold", 1.e5);
0039   _thresh_frac_timingreflm = ps.getUntrackedParameter<double>("_thresh_frac_timingreflm", 0.01);
0040   _thresh_min_lmsumq = ps.getUntrackedParameter<double>("thresh_min_lmsumq", 50000.);
0041 
0042   std::vector<double> vTimingRangeHB =
0043       ps.getUntrackedParameter<std::vector<double>>("thresh_timingreflm_HB", std::vector<double>({-70, -10.}));
0044   std::vector<double> vTimingRangeHE =
0045       ps.getUntrackedParameter<std::vector<double>>("thresh_timingreflm_HE", std::vector<double>({-60., 0.}));
0046   std::vector<double> vTimingRangeHO =
0047       ps.getUntrackedParameter<std::vector<double>>("thresh_timingreflm_HO", std::vector<double>({-50., 20.}));
0048   std::vector<double> vTimingRangeHF =
0049       ps.getUntrackedParameter<std::vector<double>>("thresh_timingreflm_HF", std::vector<double>({-50., 20.}));
0050   _thresh_timingreflm[HcalBarrel] = std::make_pair(vTimingRangeHB[0], vTimingRangeHB[1]);
0051   _thresh_timingreflm[HcalEndcap] = std::make_pair(vTimingRangeHE[0], vTimingRangeHE[1]);
0052   _thresh_timingreflm[HcalOuter] = std::make_pair(vTimingRangeHO[0], vTimingRangeHO[1]);
0053   _thresh_timingreflm[HcalForward] = std::make_pair(vTimingRangeHF[0], vTimingRangeHF[1]);
0054 }
0055 
0056 /* virtual */ void LaserTask::bookHistograms(DQMStore::IBooker& ib, edm::Run const& r, edm::EventSetup const& es) {
0057   if (_ptype == fLocal)
0058     if (r.runAuxiliary().run() == 1)
0059       return;
0060 
0061   DQTask::bookHistograms(ib, r, es);
0062 
0063   edm::ESHandle<HcalDbService> dbService = es.getHandle(hcalDbServiceToken_);
0064   _emap = dbService->getHcalMapping();
0065 
0066   std::vector<uint32_t> vhashVME;
0067   std::vector<uint32_t> vhashuTCA;
0068   std::vector<uint32_t> vhashC36;
0069   vhashVME.push_back(
0070       HcalElectronicsId(constants::FIBERCH_MIN, constants::FIBER_VME_MIN, SPIGOT_MIN, CRATE_VME_MIN).rawId());
0071   vhashuTCA.push_back(HcalElectronicsId(CRATE_uTCA_MIN, SLOT_uTCA_MIN, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0072   _filter_VME.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashVME);
0073   _filter_uTCA.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashuTCA);
0074 
0075   //    INITIALIZE
0076   _cSignalMean_Subdet.initialize(_name,
0077                                  "SignalMean",
0078                                  hcaldqm::hashfunctions::fSubdet,
0079                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fQIE10fC_100000Coarse),
0080                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0081                                  0);
0082   _cSignalRMS_Subdet.initialize(_name,
0083                                 "SignalRMS",
0084                                 hcaldqm::hashfunctions::fSubdet,
0085                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_1000),
0086                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0087                                 0);
0088   _cTimingMean_Subdet.initialize(_name,
0089                                  "TimingMean",
0090                                  hcaldqm::hashfunctions::fSubdet,
0091                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0092                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0093                                  0);
0094   _cTimingRMS_Subdet.initialize(_name,
0095                                 "TimingRMS",
0096                                 hcaldqm::hashfunctions::fSubdet,
0097                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0098                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0099                                 0);
0100 
0101   _cADC_SubdetPM.initialize(_name,
0102                             "ADC",
0103                             hcaldqm::hashfunctions::fSubdetPM,
0104                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_128),
0105                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0106                             0);
0107 
0108   if (_ptype == fLocal) {  // hidefed2crate
0109     _cSignalMean_FEDuTCA.initialize(_name,
0110                                     "SignalMean",
0111                                     hcaldqm::hashfunctions::fFED,
0112                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0113                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0114                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fQIE10fC_100000Coarse),
0115                                     0);
0116     _cSignalRMS_FEDuTCA.initialize(_name,
0117                                    "SignalRMS",
0118                                    hcaldqm::hashfunctions::fFED,
0119                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0120                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0121                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0122                                    0);
0123     _cTimingMean_FEDuTCA.initialize(_name,
0124                                     "TimingMean",
0125                                     hcaldqm::hashfunctions::fFED,
0126                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0127                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0128                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0129                                     0);
0130     _cTimingRMS_FEDuTCA.initialize(_name,
0131                                    "TimingRMS",
0132                                    hcaldqm::hashfunctions::fFED,
0133                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0134                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0135                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0136                                    0);
0137     _cMissing_FEDuTCA.initialize(_name,
0138                                  "Missing",
0139                                  hcaldqm::hashfunctions::fFED,
0140                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0141                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0142                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0143                                  0);
0144 
0145     _cShapeCut_FEDSlot.initialize(_name,
0146                                   "Shape",
0147                                   hcaldqm::hashfunctions::fFEDSlot,
0148                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS),
0149                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0150                                   0);
0151   }
0152   _cTimingvsEvent_SubdetPM.initialize(_name,
0153                                       "TimingvsEvent",
0154                                       hcaldqm::hashfunctions::fSubdetPM,
0155                                       new hcaldqm::quantity::EventNumber(_nevents),
0156                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0157                                       0);
0158   _cSignalvsEvent_SubdetPM.initialize(_name,
0159                                       "SignalvsEvent",
0160                                       hcaldqm::hashfunctions::fSubdetPM,
0161                                       new hcaldqm::quantity::EventNumber(_nevents),
0162                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0163                                       0);
0164   _cTimingvsLS_SubdetPM.initialize(_name,
0165                                    "TimingvsLS",
0166                                    hcaldqm::hashfunctions::fSubdetPM,
0167                                    new hcaldqm::quantity::LumiSection(_maxLS),
0168                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0169                                    0);
0170   _cSignalvsLS_SubdetPM.initialize(_name,
0171                                    "SignalvsLS",
0172                                    hcaldqm::hashfunctions::fSubdetPM,
0173                                    new hcaldqm::quantity::LumiSection(_maxLS),
0174                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0175                                    0);
0176   _cTimingvsBX_SubdetPM.initialize(_name,
0177                                    "TimingvsBX",
0178                                    hcaldqm::hashfunctions::fSubdetPM,
0179                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fBX),
0180                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0181                                    0);
0182   _cSignalvsBX_SubdetPM.initialize(_name,
0183                                    "SignalvsBX",
0184                                    hcaldqm::hashfunctions::fSubdetPM,
0185                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fBX),
0186                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0187                                    0);
0188 
0189   _cSignalMean_depth.initialize(_name,
0190                                 "SignalMean",
0191                                 hcaldqm::hashfunctions::fdepth,
0192                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0193                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0194                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fQIE10fC_100000Coarse),
0195                                 0);
0196   _cSignalRMS_depth.initialize(_name,
0197                                "SignalRMS",
0198                                hcaldqm::hashfunctions::fdepth,
0199                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0200                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0201                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_1000),
0202                                0);
0203   _cTimingMean_depth.initialize(_name,
0204                                 "TimingMean",
0205                                 hcaldqm::hashfunctions::fdepth,
0206                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0207                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0208                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0209                                 0);
0210   _cTimingRMS_depth.initialize(_name,
0211                                "TimingRMS",
0212                                hcaldqm::hashfunctions::fdepth,
0213                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0214                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0215                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0216                                0);
0217 
0218   _cMissing_depth.initialize(_name,
0219                              "Missing",
0220                              hcaldqm::hashfunctions::fdepth,
0221                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0222                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0223                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0224                              0);
0225 
0226   //    initialize compact containers
0227   _xSignalSum.initialize(hcaldqm::hashfunctions::fDChannel);
0228   _xSignalSum2.initialize(hcaldqm::hashfunctions::fDChannel);
0229   _xTimingSum.initialize(hcaldqm::hashfunctions::fDChannel);
0230   _xTimingSum2.initialize(hcaldqm::hashfunctions::fDChannel);
0231   _xEntries.initialize(hcaldqm::hashfunctions::fDChannel);
0232 
0233   // LaserMon containers
0234   if (_ptype == fOnline || _ptype == fLocal) {
0235     _cLaserMonSumQ.initialize(_name,
0236                               "LaserMonSumQ",
0237                               new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_1000000),
0238                               new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0239                               0);
0240     _cLaserMonSumQ.showOverflowX(true);
0241 
0242     _cLaserMonTiming.initialize(_name,
0243                                 "LaserMonTiming",
0244                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTime_ns_250),
0245                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0246                                 0);
0247     _cLaserMonTiming.showOverflowX(true);
0248 
0249     if (_ptype == fOnline) {
0250       _cLaserMonSumQ_LS.initialize(_name,
0251                                    "LaserMonSumQ_LS",
0252                                    new hcaldqm::quantity::LumiSectionCoarse(_maxLS, 10),
0253                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_1000000));
0254       _cLaserMonTiming_LS.initialize(_name,
0255                                      "LaserMonTiming_LS",
0256                                      new hcaldqm::quantity::LumiSectionCoarse(_maxLS, 10),
0257                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_100TS));
0258       _cTimingDiffLS_SubdetPM.initialize(_name,
0259                                          "TimingDiff_DigiMinusLaserMon",
0260                                          hcaldqm::hashfunctions::fSubdetPM,
0261                                          new hcaldqm::quantity::LumiSectionCoarse(_maxLS, 10),
0262                                          new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fRBX),
0263                                          new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTimingDiff_ns),
0264                                          0);
0265       _cTimingDiffLS_SubdetPM.showOverflowY(true);
0266     } else if (_ptype == fLocal) {
0267       _cLaserMonSumQ_Event.initialize(_name,
0268                                       "LaserMonSumQ_Event",
0269                                       new hcaldqm::quantity::EventNumber(_nevents),
0270                                       new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_1000000));
0271       _cLaserMonTiming_Event.initialize(_name,
0272                                         "LaserMonTiming_Event",
0273                                         new hcaldqm::quantity::EventNumber(_nevents),
0274                                         new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_100TS));
0275       _cTimingDiffEvent_SubdetPM.initialize(_name,
0276                                             "TimingDiff_DigiMinusLaserMon",
0277                                             hcaldqm::hashfunctions::fSubdetPM,
0278                                             new hcaldqm::quantity::EventNumber(_nevents),
0279                                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fRBX),
0280                                             new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTimingDiff_ns),
0281                                             0);
0282       _cTimingDiffEvent_SubdetPM.showOverflowY(true);
0283     }
0284     _cTiming_DigivsLaserMon_SubdetPM.initialize(
0285         _name,
0286         "Timing_DigivsLaserMon",
0287         hcaldqm::hashfunctions::fSubdetPM,
0288         new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTime_ns_250_coarse),
0289         new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTime_ns_250_coarse),
0290         new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0291         0);
0292     _cTiming_DigivsLaserMon_SubdetPM.showOverflowX(true);
0293     _cTiming_DigivsLaserMon_SubdetPM.showOverflowY(true);
0294 
0295     _xTimingRefLMSum.initialize(hcaldqm::hashfunctions::fDChannel);
0296     _xTimingRefLMSum2.initialize(hcaldqm::hashfunctions::fDChannel);
0297     _xNBadTimingRefLM.initialize(hcaldqm::hashfunctions::fFED);
0298     _xNChs.initialize(hcaldqm::hashfunctions::fFED);
0299     _xMissingLaserMon = 0;
0300 
0301     std::vector<int> vFEDs = hcaldqm::utilities::getFEDList(_emap);
0302     std::vector<int> vFEDsVME = hcaldqm::utilities::getFEDVMEList(_emap);
0303     std::vector<int> vFEDsuTCA = hcaldqm::utilities::getFEDuTCAList(_emap);
0304     for (std::vector<int>::const_iterator it = vFEDsVME.begin(); it != vFEDsVME.end(); ++it)
0305       _vhashFEDs.push_back(
0306           HcalElectronicsId(constants::FIBERCH_MIN, FIBER_VME_MIN, SPIGOT_MIN, (*it) - FED_VME_MIN).rawId());
0307     for (std::vector<int>::const_iterator it = vFEDsuTCA.begin(); it != vFEDsuTCA.end(); ++it) {
0308       std::pair<uint16_t, uint16_t> cspair = utilities::fed2crate(*it);
0309       _vhashFEDs.push_back(HcalElectronicsId(cspair.first, cspair.second, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0310     }
0311 
0312     _cSummaryvsLS_FED.initialize(_name,
0313                                  "SummaryvsLS",
0314                                  hcaldqm::hashfunctions::fFED,
0315                                  new hcaldqm::quantity::LumiSection(_maxLS),
0316                                  new hcaldqm::quantity::FlagQuantity(_vflags),
0317                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fState),
0318                                  0);
0319     _cSummaryvsLS.initialize(_name,
0320                              "SummaryvsLS",
0321                              new hcaldqm::quantity::LumiSection(_maxLS),
0322                              new hcaldqm::quantity::FEDQuantity(vFEDs),
0323                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fState),
0324                              0);
0325   }  // End if (_ptype == fOnline || _ptype == fLocal) {
0326 
0327   //    BOOK
0328   _cSignalMean_Subdet.book(ib, _emap, _subsystem);
0329   _cSignalRMS_Subdet.book(ib, _emap, _subsystem);
0330   _cTimingMean_Subdet.book(ib, _emap, _subsystem);
0331   _cTimingRMS_Subdet.book(ib, _emap, _subsystem);
0332 
0333   _cSignalMean_depth.book(ib, _emap, _subsystem);
0334   _cSignalRMS_depth.book(ib, _emap, _subsystem);
0335   _cTimingMean_depth.book(ib, _emap, _subsystem);
0336   _cTimingRMS_depth.book(ib, _emap, _subsystem);
0337 
0338   if (_ptype == fLocal) {
0339     _cTimingvsEvent_SubdetPM.book(ib, _emap, _subsystem);
0340     _cSignalvsEvent_SubdetPM.book(ib, _emap, _subsystem);
0341   } else {
0342     _cTimingvsLS_SubdetPM.book(ib, _emap, _subsystem);
0343     _cSignalvsLS_SubdetPM.book(ib, _emap, _subsystem);
0344     _cTimingvsBX_SubdetPM.book(ib, _emap, _subsystem);
0345     _cSignalvsBX_SubdetPM.book(ib, _emap, _subsystem);
0346   }
0347 
0348   if (_ptype == fLocal) {  // hidefed2crate
0349     _cSignalMean_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0350     _cSignalRMS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0351     _cTimingMean_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0352     _cTimingRMS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0353     _cMissing_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0354     _cShapeCut_FEDSlot.book(ib, _emap, _subsystem);
0355   }
0356   _cADC_SubdetPM.book(ib, _emap, _subsystem);
0357 
0358   _cMissing_depth.book(ib, _emap, _subsystem);
0359 
0360   _xSignalSum.book(_emap);
0361   _xSignalSum2.book(_emap);
0362   _xEntries.book(_emap);
0363   _xTimingSum.book(_emap);
0364   _xTimingSum2.book(_emap);
0365 
0366   if (_ptype == fOnline || _ptype == fLocal) {
0367     _cLaserMonSumQ.book(ib, _subsystem);
0368     _cLaserMonTiming.book(ib, _subsystem);
0369     if (_ptype == fOnline) {
0370       _cLaserMonSumQ_LS.book(ib, _subsystem);
0371       _cLaserMonTiming_LS.book(ib, _subsystem);
0372       _cTimingDiffLS_SubdetPM.book(ib, _emap, _subsystem);
0373     } else if (_ptype == fLocal) {
0374       _cLaserMonSumQ_Event.book(ib, _subsystem);
0375       _cLaserMonTiming_Event.book(ib, _subsystem);
0376       _cTimingDiffEvent_SubdetPM.book(ib, _emap, _subsystem);
0377     }
0378     _cTiming_DigivsLaserMon_SubdetPM.book(ib, _emap, _subsystem);
0379     _xTimingRefLMSum.book(_emap);
0380     _xTimingRefLMSum2.book(_emap);
0381     _xNBadTimingRefLM.book(_emap);
0382     _xNChs.book(_emap);
0383     _cSummaryvsLS_FED.book(ib, _emap, _subsystem);
0384     _cSummaryvsLS.book(ib, _subsystem);
0385   }
0386 
0387   _ehashmap.initialize(_emap, electronicsmap::fD2EHashMap);
0388 }
0389 
0390 /* virtual */ void LaserTask::_resetMonitors(hcaldqm::UpdateFreq uf) { DQTask::_resetMonitors(uf); }
0391 
0392 /* virtual */ void LaserTask::_dump() {
0393   _cSignalMean_Subdet.reset();
0394   _cSignalRMS_Subdet.reset();
0395   _cTimingMean_Subdet.reset();
0396   _cTimingRMS_Subdet.reset();
0397   _cSignalMean_depth.reset();
0398   _cSignalRMS_depth.reset();
0399   _cTimingMean_depth.reset();
0400   _cTimingRMS_depth.reset();
0401 
0402   if (_ptype == fLocal) {  // hidefed2crate
0403     _cSignalMean_FEDuTCA.reset();
0404     _cSignalRMS_FEDuTCA.reset();
0405     _cTimingMean_FEDuTCA.reset();
0406     _cTimingRMS_FEDuTCA.reset();
0407   }
0408   if (_ptype != fOffline) {
0409     _xNChs.reset();
0410   }
0411 
0412   std::vector<HcalGenericDetId> dids = _emap->allPrecisionId();
0413   for (std::vector<HcalGenericDetId>::const_iterator it = dids.begin(); it != dids.end(); ++it) {
0414     if (!it->isHcalDetId())
0415       continue;
0416     HcalDetId did = HcalDetId(it->rawId());
0417     HcalElectronicsId eid(_ehashmap.lookup(*it));
0418     int n = _xEntries.get(did);
0419     //  channels missing or low signal
0420     if (n == 0) {
0421       _cMissing_depth.fill(did);
0422       if (_ptype == fLocal) {  // hidefed2crate
0423         if (!eid.isVMEid())
0424           _cMissing_FEDuTCA.fill(eid);
0425       }
0426       continue;
0427     }
0428 
0429     _xNChs.get(eid)++;
0430 
0431     double msig = _xSignalSum.get(did) / n;
0432     double mtim = _xTimingSum.get(did) / n;
0433     double rsig = sqrt(_xSignalSum2.get(did) / n - msig * msig);
0434     double rtim = sqrt(_xTimingSum2.get(did) / n - mtim * mtim);
0435 
0436     _cSignalMean_Subdet.fill(did, msig);
0437     _cSignalMean_depth.fill(did, msig);
0438     _cSignalRMS_Subdet.fill(did, rsig);
0439     _cSignalRMS_depth.fill(did, rsig);
0440     _cTimingMean_Subdet.fill(did, mtim);
0441     _cTimingMean_depth.fill(did, mtim);
0442     _cTimingRMS_Subdet.fill(did, rtim);
0443     _cTimingRMS_depth.fill(did, rtim);
0444     if (_ptype == fLocal) {  // hidefed2crate
0445       if (!eid.isVMEid()) {
0446         _cSignalMean_FEDuTCA.fill(eid, msig);
0447         _cSignalRMS_FEDuTCA.fill(eid, rsig);
0448         _cTimingMean_FEDuTCA.fill(eid, mtim);
0449         _cTimingRMS_FEDuTCA.fill(eid, rtim);
0450       }
0451     }
0452 
0453     // Bad timing
0454 
0455     double timingreflm_mean = _xTimingRefLMSum.get(did) / n;
0456     //double timingreflm_rms = sqrt(_xTimingRefLMSum2.get(did) / n - timingreflm_mean * timingreflm_mean);
0457 
0458     if ((timingreflm_mean < _thresh_timingreflm[did.subdet()].first) ||
0459         (timingreflm_mean > _thresh_timingreflm[did.subdet()].second)) {
0460       _xNBadTimingRefLM.get(eid)++;
0461     }
0462   }
0463   if (_ptype != fOffline) {  // hidefed2crate
0464     for (std::vector<uint32_t>::const_iterator it = _vhashFEDs.begin(); it != _vhashFEDs.end(); ++it) {
0465       hcaldqm::flag::Flag fSum("LASER");
0466       HcalElectronicsId eid = HcalElectronicsId(*it);
0467       std::vector<uint32_t>::const_iterator jt = std::find(_vcdaqEids.begin(), _vcdaqEids.end(), (*it));
0468       if (jt == _vcdaqEids.end()) {
0469         //  not @cDAQ
0470         for (uint32_t iflag = 0; iflag < _vflags.size(); iflag++)
0471           _cSummaryvsLS_FED.setBinContent(eid, _currentLS, int(iflag), int(hcaldqm::flag::fNCDAQ));
0472         _cSummaryvsLS.setBinContent(eid, _currentLS, int(hcaldqm::flag::fNCDAQ));
0473         continue;
0474       }
0475       //    @cDAQ
0476       if (hcaldqm::utilities::isFEDHBHE(eid) || hcaldqm::utilities::isFEDHO(eid) || hcaldqm::utilities::isFEDHF(eid)) {
0477         int min_nchs = 10;
0478         if (_xNChs.get(eid) < min_nchs) {
0479           _vflags[fBadTiming]._state = hcaldqm::flag::fNA;
0480         } else {
0481           double frbadtimingreflm = double(_xNBadTimingRefLM.get(eid)) / double(_xNChs.get(eid));
0482           if (frbadtimingreflm > _thresh_frac_timingreflm) {
0483             _vflags[fBadTiming]._state = hcaldqm::flag::fBAD;
0484           } else {
0485             _vflags[fBadTiming]._state = hcaldqm::flag::fGOOD;
0486           }
0487         }
0488         if (_xMissingLaserMon) {
0489           _vflags[fMissingLaserMon]._state = hcaldqm::flag::fBAD;
0490         } else {
0491           _vflags[fMissingLaserMon]._state = hcaldqm::flag::fGOOD;
0492         }
0493       }
0494 
0495       // Set SummaryVsLS bins
0496       int iflag = 0;
0497       for (std::vector<hcaldqm::flag::Flag>::iterator ft = _vflags.begin(); ft != _vflags.end(); ++ft) {
0498         _cSummaryvsLS_FED.setBinContent(eid, _currentLS, iflag, int(ft->_state));
0499         fSum += (*ft);
0500         iflag++;
0501         ft->reset();
0502       }
0503       _cSummaryvsLS.setBinContent(eid, _currentLS, fSum._state);
0504     }  // End loop over FEDs
0505 
0506     // Reset per-LS containers
0507     _xNBadTimingRefLM.reset();
0508     _xMissingLaserMon = 0;
0509   }  // End if _ptype != fOffline
0510 }
0511 
0512 /* virtual */ void LaserTask::_process(edm::Event const& e, edm::EventSetup const& es) {
0513   edm::Handle<QIE11DigiCollection> c_QIE11;
0514   edm::Handle<HODigiCollection> c_ho;
0515   edm::Handle<QIE10DigiCollection> c_QIE10;
0516 
0517   if (!e.getByToken(_tokQIE11, c_QIE11))
0518     _logger.dqmthrow("Collection QIE11DigiCollection isn't available " + _tagQIE11.label() + " " +
0519                      _tagQIE11.instance());
0520   if (!e.getByToken(_tokHO, c_ho))
0521     _logger.dqmthrow("Collection HODigiCollection isn't available " + _tagHO.label() + " " + _tagHO.instance());
0522   if (!e.getByToken(_tokQIE10, c_QIE10))
0523     _logger.dqmthrow("Collection QIE10DigiCollection isn't available " + _tagQIE10.label() + " " +
0524                      _tagQIE10.instance());
0525 
0526   //    int currentEvent = e.eventAuxiliary().id().event();
0527   int bx = e.bunchCrossing();
0528 
0529   // LASERMON
0530   edm::Handle<QIE10DigiCollection> cLaserMon;
0531   if (!e.getByToken(_tokLaserMon, cLaserMon)) {
0532     _logger.dqmthrow("QIE10DigiCollection for laserMonDigis isn't available.");
0533   }
0534   std::vector<int> laserMonADC;
0535   processLaserMon(cLaserMon, laserMonADC);
0536 
0537   // SumQ = peak +/- 3 TSes
0538   // Timing = fC-weighted average (TS-TS0) * 25 ns, also in peak +/- 3 TSes
0539   int peakTS = -1;
0540   double peakLaserMonADC = -1;
0541   for (unsigned int iTS = 0; iTS < laserMonADC.size(); ++iTS) {
0542     if (laserMonADC[iTS] > peakLaserMonADC) {
0543       peakLaserMonADC = laserMonADC[iTS];
0544       peakTS = iTS;
0545     }
0546   }
0547 
0548   double laserMonSumQ = 0;
0549   double laserMonTiming = 0.;
0550 
0551   if (peakTS >= 0) {
0552     int minTS = std::max(0, peakTS - 3);
0553     int maxTS = std::min(int(laserMonADC.size() - 1), peakTS + 3);
0554     for (int iTS = minTS; iTS <= maxTS; ++iTS) {
0555       double this_fC = hcaldqm::constants::adc2fC[laserMonADC[iTS]];
0556       laserMonSumQ += this_fC;
0557       laserMonTiming += 25. * (iTS - _laserMonTS0) * this_fC;
0558     }
0559   }
0560   if (laserMonSumQ > 0.) {
0561     laserMonTiming = laserMonTiming / laserMonSumQ;
0562   } else {
0563     ++_xMissingLaserMon;
0564   }
0565 
0566   if (laserMonSumQ > _laserMonThreshold) {
0567     _cLaserMonSumQ.fill(laserMonSumQ);
0568     _cLaserMonTiming.fill(laserMonTiming);
0569     if (_ptype == fOnline) {
0570       _cLaserMonSumQ_LS.fill(_currentLS, laserMonSumQ);
0571       _cLaserMonTiming_LS.fill(_currentLS, laserMonTiming);
0572     } else if (_ptype == fLocal) {
0573       int currentEvent = e.eventAuxiliary().id().event();
0574       _cLaserMonSumQ_Event.fill(currentEvent, laserMonSumQ);
0575       _cLaserMonTiming_Event.fill(currentEvent, laserMonTiming);
0576     }
0577   }
0578 
0579   for (QIE11DigiCollection::const_iterator it = c_QIE11->begin(); it != c_QIE11->end(); ++it) {
0580     const QIE11DataFrame digi = static_cast<const QIE11DataFrame>(*it);
0581     HcalDetId const& did = digi.detid();
0582     if ((did.subdet() != HcalBarrel) && (did.subdet() != HcalEndcap)) {
0583       continue;
0584     }
0585     uint32_t rawid = _ehashmap.lookup(did);
0586     HcalElectronicsId const& eid(rawid);
0587 
0588     CaloSamples digi_fC = hcaldqm::utilities::loadADC2fCDB<QIE11DataFrame>(_dbService, did, digi);
0589     //double sumQ = hcaldqm::utilities::sumQ_v10<QIE11DataFrame>(digi, 2.5, 0, digi.samples()-1);
0590     double sumQ = hcaldqm::utilities::sumQDB<QIE11DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0591     if (sumQ < _lowHBHE)
0592       continue;
0593 
0594     //double aveTS = hcaldqm::utilities::aveTS_v10<QIE11DataFrame>(digi, 2.5, 0,digi.samples()-1);
0595     double aveTS = hcaldqm::utilities::aveTSDB<QIE11DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0596     _xSignalSum.get(did) += sumQ;
0597     _xSignalSum2.get(did) += sumQ * sumQ;
0598     _xTimingSum.get(did) += aveTS;
0599     _xTimingSum2.get(did) += aveTS * aveTS;
0600     _xEntries.get(did)++;
0601 
0602     for (int i = 0; i < digi.samples(); i++) {
0603       if (_ptype == fLocal) {
0604         _cShapeCut_FEDSlot.fill(
0605             eid, i, hcaldqm::utilities::adc2fCDBMinusPedestal<QIE11DataFrame>(_dbService, digi_fC, did, digi, i));
0606       }
0607       _cADC_SubdetPM.fill(did, digi[i].adc());
0608     }
0609 
0610     //  select based on local global
0611     double digiTimingSOI = (aveTS - digi.presamples()) * 25.;
0612     double deltaTiming = digiTimingSOI - laserMonTiming;
0613     _cTiming_DigivsLaserMon_SubdetPM.fill(did, laserMonTiming, digiTimingSOI);
0614     _xTimingRefLMSum.get(did) += deltaTiming;
0615     _xTimingRefLMSum2.get(did) += deltaTiming * deltaTiming;
0616     if (_ptype == fLocal) {
0617       int currentEvent = e.eventAuxiliary().id().event();
0618       _cTimingvsEvent_SubdetPM.fill(did, currentEvent, aveTS);
0619       _cSignalvsEvent_SubdetPM.fill(did, currentEvent, sumQ);
0620       _cTimingDiffEvent_SubdetPM.fill(did, currentEvent, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0621     } else {
0622       _cTimingvsLS_SubdetPM.fill(did, _currentLS, aveTS);
0623       _cSignalvsLS_SubdetPM.fill(did, _currentLS, sumQ);
0624       _cTimingvsBX_SubdetPM.fill(did, bx, aveTS);
0625       _cSignalvsBX_SubdetPM.fill(did, bx, sumQ);
0626       _cTimingDiffLS_SubdetPM.fill(did, _currentLS, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0627     }
0628   }
0629   for (HODigiCollection::const_iterator it = c_ho->begin(); it != c_ho->end(); ++it) {
0630     const HODataFrame digi = (const HODataFrame)(*it);
0631     double sumQ = hcaldqm::utilities::sumQ<HODataFrame>(digi, 8.5, 0, digi.size() - 1);
0632     if (sumQ < _lowHO)
0633       continue;
0634     HcalDetId did = digi.id();
0635     HcalElectronicsId eid = digi.elecId();
0636 
0637     double aveTS = hcaldqm::utilities::aveTS<HODataFrame>(digi, 8.5, 0, digi.size() - 1);
0638     _xSignalSum.get(did) += sumQ;
0639     _xSignalSum2.get(did) += sumQ * sumQ;
0640     _xTimingSum.get(did) += aveTS;
0641     _xTimingSum2.get(did) += aveTS * aveTS;
0642     _xEntries.get(did)++;
0643 
0644     for (int i = 0; i < digi.size(); i++) {
0645       if (_ptype == fLocal) {  // hidefed2crate
0646         _cShapeCut_FEDSlot.fill(eid, i, digi.sample(i).nominal_fC() - 8.5);
0647       }
0648       _cADC_SubdetPM.fill(did, digi.sample(i).adc());
0649     }
0650 
0651     //  select based on local global
0652     double digiTimingSOI = (aveTS - digi.presamples()) * 25.;
0653     double deltaTiming = digiTimingSOI - laserMonTiming;
0654     _cTiming_DigivsLaserMon_SubdetPM.fill(did, laserMonTiming, digiTimingSOI);
0655     _xTimingRefLMSum.get(did) += deltaTiming;
0656     _xTimingRefLMSum2.get(did) += deltaTiming * deltaTiming;
0657     if (_ptype == fLocal) {
0658       int currentEvent = e.eventAuxiliary().id().event();
0659       _cTimingvsEvent_SubdetPM.fill(did, currentEvent, aveTS);
0660       _cSignalvsEvent_SubdetPM.fill(did, currentEvent, sumQ);
0661       _cTimingDiffEvent_SubdetPM.fill(did, currentEvent, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0662     } else {
0663       _cTimingvsLS_SubdetPM.fill(did, _currentLS, aveTS);
0664       _cSignalvsLS_SubdetPM.fill(did, _currentLS, sumQ);
0665       _cTimingvsBX_SubdetPM.fill(did, bx, aveTS);
0666       _cSignalvsBX_SubdetPM.fill(did, bx, sumQ);
0667       _cTimingDiffLS_SubdetPM.fill(did, _currentLS, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0668     }
0669   }
0670   for (QIE10DigiCollection::const_iterator it = c_QIE10->begin(); it != c_QIE10->end(); ++it) {
0671     const QIE10DataFrame digi = (const QIE10DataFrame)(*it);
0672     HcalDetId did = digi.detid();
0673     if (did.subdet() != HcalForward) {
0674       continue;
0675     }
0676     HcalElectronicsId eid = HcalElectronicsId(_ehashmap.lookup(did));
0677 
0678     CaloSamples digi_fC = hcaldqm::utilities::loadADC2fCDB<QIE10DataFrame>(_dbService, did, digi);
0679     double sumQ = hcaldqm::utilities::sumQDB<QIE10DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0680     //double sumQ = hcaldqm::utilities::sumQ_v10<QIE10DataFrame>(digi, 2.5, 0, digi.samples()-1);
0681     if (sumQ < _lowHF)
0682       continue;
0683 
0684     //double aveTS = hcaldqm::utilities::aveTS_v10<QIE10DataFrame>(digi, 2.5, 0, digi.samples()-1);
0685     double aveTS = hcaldqm::utilities::aveTSDB<QIE10DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0686 
0687     _xSignalSum.get(did) += sumQ;
0688     _xSignalSum2.get(did) += sumQ * sumQ;
0689     _xTimingSum.get(did) += aveTS;
0690     _xTimingSum2.get(did) += aveTS * aveTS;
0691     _xEntries.get(did)++;
0692 
0693     for (int i = 0; i < digi.samples(); i++) {
0694       if (_ptype == fLocal) {  // hidefed2crate
0695         _cShapeCut_FEDSlot.fill(
0696             eid, (int)i, hcaldqm::utilities::adc2fCDBMinusPedestal<QIE10DataFrame>(_dbService, digi_fC, did, digi, i));
0697       }
0698       _cADC_SubdetPM.fill(did, digi[i].adc());
0699     }
0700 
0701     //  select based on local global
0702     double digiTimingSOI = (aveTS - digi.presamples()) * 25.;
0703     double deltaTiming = digiTimingSOI - laserMonTiming;
0704     _cTiming_DigivsLaserMon_SubdetPM.fill(did, laserMonTiming, digiTimingSOI);
0705     _xTimingRefLMSum.get(did) += deltaTiming;
0706     _xTimingRefLMSum2.get(did) += deltaTiming * deltaTiming;
0707     if (_ptype == fLocal) {
0708       int currentEvent = e.eventAuxiliary().id().event();
0709       _cTimingvsEvent_SubdetPM.fill(did, currentEvent, aveTS);
0710       _cSignalvsEvent_SubdetPM.fill(did, currentEvent, sumQ);
0711       _cTimingDiffEvent_SubdetPM.fill(did, currentEvent, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0712     } else {
0713       _cTimingvsLS_SubdetPM.fill(did, _currentLS, aveTS);
0714       _cSignalvsLS_SubdetPM.fill(did, _currentLS, sumQ);
0715       _cTimingvsBX_SubdetPM.fill(did, bx, aveTS);
0716       _cSignalvsBX_SubdetPM.fill(did, bx, sumQ);
0717       _cTimingDiffLS_SubdetPM.fill(did, _currentLS, hcaldqm::utilities::getRBX(did.iphi()), deltaTiming);
0718     }
0719   }
0720 }
0721 
0722 void LaserTask::processLaserMon(edm::Handle<QIE10DigiCollection>& col, std::vector<int>& iLaserMonADC) {
0723   for (QIE10DigiCollection::const_iterator it = col->begin(); it != col->end(); ++it) {
0724     const QIE10DataFrame digi = (const QIE10DataFrame)(*it);
0725     HcalCalibDetId hcdid(digi.id());
0726 
0727     if ((hcdid.ieta() != _laserMonIEta) || (hcdid.cboxChannel() != _laserMonCBox)) {
0728       continue;
0729     }
0730 
0731     unsigned int digiIndex =
0732         std::find(_vLaserMonIPhi.begin(), _vLaserMonIPhi.end(), hcdid.iphi()) - _vLaserMonIPhi.begin();
0733     if (digiIndex == _vLaserMonIPhi.size()) {
0734       continue;
0735     }
0736 
0737     // First digi: initialize the vectors to -1 (need to know the length of the digi)
0738     if (iLaserMonADC.empty()) {
0739       int totalNSamples = (digi.samples() - _laserMonDigiOverlap) * _vLaserMonIPhi.size();
0740       for (int i = 0; i < totalNSamples; ++i) {
0741         iLaserMonADC.push_back(-1);
0742       }
0743     }
0744 
0745     for (int subindex = 0; subindex < digi.samples() - _laserMonDigiOverlap; ++subindex) {
0746       int totalIndex = (digi.samples() - _laserMonDigiOverlap) * digiIndex + subindex;
0747       iLaserMonADC[totalIndex] = (digi[subindex].ok() ? digi[subindex].adc() : -1);
0748     }
0749   }
0750 }
0751 
0752 /* virtual */ void LaserTask::globalEndLuminosityBlock(edm::LuminosityBlock const& lb, edm::EventSetup const& es) {
0753   auto lumiCache = luminosityBlockCache(lb.index());
0754   _currentLS = lumiCache->currentLS;
0755 
0756   if (_ptype == fLocal)
0757     return;
0758   this->_dump();
0759 
0760   DQTask::globalEndLuminosityBlock(lb, es);
0761 }
0762 
0763 /* virtual */ bool LaserTask::_isApplicable(edm::Event const& e) {
0764   if (_ptype != fOnline)
0765     return true;
0766   else {
0767     //  fOnline mode
0768     edm::Handle<HcalUMNioDigi> cumn;
0769     if (!e.getByToken(_tokuMN, cumn))
0770       return false;
0771 
0772     //  event type check first
0773     uint8_t eventType = cumn->eventType();
0774     if (eventType != constants::EVENTTYPE_LASER)
0775       return false;
0776 
0777     //  check if this analysis task is of the right laser type
0778     uint32_t laserType = cumn->valueUserWord(0);
0779     if (laserType == _laserType)
0780       return true;
0781   }
0782 
0783   return false;
0784 }
0785 
0786 DEFINE_FWK_MODULE(LaserTask);