Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 #include "DQM/HcalTasks/interface/LEDTask.h"
0003 
0004 using namespace hcaldqm;
0005 using namespace hcaldqm::constants;
0006 using namespace hcaldqm::filter;
0007 
0008 LEDTask::LEDTask(edm::ParameterSet const& ps)
0009     : DQTask(ps), hcalDbServiceToken_(esConsumes<HcalDbService, HcalDbRecord, edm::Transition::BeginRun>()) {
0010   _nevents = ps.getUntrackedParameter<int>("nevents", 2000);
0011   //    tags
0012   _tagQIE11 = ps.getUntrackedParameter<edm::InputTag>("tagHE", edm::InputTag("hcalDigis"));
0013   _tagHO = ps.getUntrackedParameter<edm::InputTag>("tagHO", edm::InputTag("hcalDigis"));
0014   _tagQIE10 = ps.getUntrackedParameter<edm::InputTag>("tagHF", edm::InputTag("hcalDigis"));
0015   _tagTrigger = ps.getUntrackedParameter<edm::InputTag>("tagTrigger", edm::InputTag("tbunpacker"));
0016   _taguMN = ps.getUntrackedParameter<edm::InputTag>("taguMN", edm::InputTag("hcalDigis"));
0017   _tokQIE11 = consumes<QIE11DigiCollection>(_tagQIE11);
0018   _tokHO = consumes<HODigiCollection>(_tagHO);
0019   _tokQIE10 = consumes<QIE10DigiCollection>(_tagQIE10);
0020   _tokTrigger = consumes<HcalTBTriggerData>(_tagTrigger);
0021   _tokuMN = consumes<HcalUMNioDigi>(_taguMN);
0022 
0023   //    constants
0024   _lowHBHE = ps.getUntrackedParameter<double>("lowHBHE", 20);
0025   _lowHO = ps.getUntrackedParameter<double>("lowHO", 20);
0026   _lowHF = ps.getUntrackedParameter<double>("lowHF", 20);
0027 }
0028 
0029 /* virtual */ void LEDTask::bookHistograms(DQMStore::IBooker& ib, edm::Run const& r, edm::EventSetup const& es) {
0030   if (_ptype == fLocal)
0031     if (r.runAuxiliary().run() == 1)
0032       return;
0033 
0034   DQTask::bookHistograms(ib, r, es);
0035 
0036   edm::ESHandle<HcalDbService> dbService = es.getHandle(hcalDbServiceToken_);
0037   _emap = dbService->getHcalMapping();
0038 
0039   // Book LED calibration channels from emap
0040   std::vector<HcalElectronicsId> eids = _emap->allElectronicsId();
0041   for (unsigned i = 0; i < eids.size(); i++) {
0042     HcalElectronicsId eid = eids[i];
0043     DetId id = _emap->lookup(eid);
0044     if (HcalGenericDetId(id.rawId()).isHcalCalibDetId()) {
0045       HcalCalibDetId calibId(id);
0046       if (calibId.calibFlavor() == HcalCalibDetId::CalibrationBox) {
0047         HcalSubdetector this_subdet = HcalEmpty;
0048         switch (calibId.hcalSubdet()) {
0049           case HcalBarrel:
0050             this_subdet = HcalBarrel;
0051             break;
0052           case HcalEndcap:
0053             this_subdet = HcalEndcap;
0054             break;
0055           case HcalOuter:
0056             this_subdet = HcalOuter;
0057             break;
0058           case HcalForward:
0059             this_subdet = HcalForward;
0060             break;
0061           default:
0062             this_subdet = HcalEmpty;
0063             break;
0064         }
0065         _ledCalibrationChannels[this_subdet].push_back(
0066             HcalDetId(HcalOther, calibId.ieta(), calibId.iphi(), calibId.cboxChannel()));
0067       }
0068     }
0069   }
0070 
0071   std::vector<uint32_t> vhashVME;
0072   std::vector<uint32_t> vhashuTCA;
0073   std::vector<uint32_t> vhashC36;
0074   vhashVME.push_back(
0075       HcalElectronicsId(constants::FIBERCH_MIN, constants::FIBER_VME_MIN, SPIGOT_MIN, CRATE_VME_MIN).rawId());
0076   vhashuTCA.push_back(HcalElectronicsId(CRATE_uTCA_MIN, SLOT_uTCA_MIN, FIBER_uTCA_MIN1, FIBERCH_MIN, false).rawId());
0077   _filter_VME.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashVME);
0078   _filter_uTCA.initialize(filter::fFilter, hcaldqm::hashfunctions::fElectronics, vhashuTCA);
0079 
0080   //    INITIALIZE
0081   _cSignalMean_Subdet.initialize(_name,
0082                                  "SignalMean",
0083                                  hcaldqm::hashfunctions::fSubdet,
0084                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0085                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0086                                  0);
0087   _cSignalRMS_Subdet.initialize(_name,
0088                                 "SignalRMS",
0089                                 hcaldqm::hashfunctions::fSubdet,
0090                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_3000),
0091                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0092                                 0);
0093   _cTimingMean_Subdet.initialize(_name,
0094                                  "TimingMean",
0095                                  hcaldqm::hashfunctions::fSubdet,
0096                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0097                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0098                                  0);
0099   _cTimingRMS_Subdet.initialize(_name,
0100                                 "TimingRMS",
0101                                 hcaldqm::hashfunctions::fSubdet,
0102                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0103                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0104                                 0);
0105 
0106   if (_ptype == fLocal) {  // hidefed2crate
0107     _cSignalMean_FEDuTCA.initialize(_name,
0108                                     "SignalMean",
0109                                     hcaldqm::hashfunctions::fFED,
0110                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0111                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0112                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_generic_400000, true),
0113                                     0);
0114     _cSignalRMS_FEDuTCA.initialize(_name,
0115                                    "SignalRMS",
0116                                    hcaldqm::hashfunctions::fFED,
0117                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0118                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0119                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_generic_10000, true),
0120                                    0);
0121     _cTimingMean_FEDuTCA.initialize(_name,
0122                                     "TimingMean",
0123                                     hcaldqm::hashfunctions::fFED,
0124                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0125                                     new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0126                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0127                                     0);
0128     _cTimingRMS_FEDuTCA.initialize(_name,
0129                                    "TimingRMS",
0130                                    hcaldqm::hashfunctions::fFED,
0131                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0132                                    new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0133                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0134                                    0);
0135 
0136     _cShapeCut_FEDSlot.initialize(_name,
0137                                   "Shape",
0138                                   hcaldqm::hashfunctions::fFEDSlot,
0139                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS),
0140                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_generic_400000),
0141                                   0);
0142   }
0143 
0144   _cSignalMean_depth.initialize(_name,
0145                                 "SignalMean",
0146                                 hcaldqm::hashfunctions::fdepth,
0147                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0148                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0149                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_generic_400000, true),
0150                                 0);
0151   _cSignalRMS_depth.initialize(_name,
0152                                "SignalRMS",
0153                                hcaldqm::hashfunctions::fdepth,
0154                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0155                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0156                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::ffC_generic_10000, true),
0157                                0);
0158   _cTimingMean_depth.initialize(_name,
0159                                 "TimingMean",
0160                                 hcaldqm::hashfunctions::fdepth,
0161                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0162                                 new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0163                                 new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0164                                 0);
0165   _cTimingRMS_depth.initialize(_name,
0166                                "TimingRMS",
0167                                hcaldqm::hashfunctions::fdepth,
0168                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0169                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0170                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS200),
0171                                0);
0172 
0173   _cMissing_depth.initialize(_name,
0174                              "Missing",
0175                              hcaldqm::hashfunctions::fdepth,
0176                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0177                              new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0178                              new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0179                              0);
0180   if (_ptype != fOffline) {  // hidefed2crate
0181     _cMissing_FEDuTCA.initialize(_name,
0182                                  "Missing",
0183                                  hcaldqm::hashfunctions::fFED,
0184                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0185                                  new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fFiberuTCAFiberCh),
0186                                  new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0187                                  0);
0188   }
0189 
0190   // Plots for LED in global
0191   if (_ptype == fOnline || _ptype == fLocal) {
0192     _cADCvsTS_SubdetPM.initialize(_name,
0193                                   "ADCvsTS",
0194                                   hcaldqm::hashfunctions::fSubdetPM,
0195                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTiming_TS),
0196                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fQIE10ADC_256),
0197                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0198                                   0);
0199   }
0200   if (_ptype == fOnline) {
0201     _cLowSignal_CrateSlot.initialize(_name,
0202                                      "LowSignal",
0203                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fCrateuTCA),
0204                                      new hcaldqm::quantity::ElectronicsQuantity(hcaldqm::quantity::fSlotuTCA),
0205                                      new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0206                                      0);
0207     _cSumQ_SubdetPM.initialize(_name,
0208                                "SumQ",
0209                                hcaldqm::hashfunctions::fSubdetPM,
0210                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fQIE10fC_400000),
0211                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0212                                0);
0213     _cTDCTime_SubdetPM.initialize(_name,
0214                                   "TDCTime",
0215                                   hcaldqm::hashfunctions::fSubdetPM,
0216                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTime_ns_250),
0217                                   new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN, true),
0218                                   0);
0219     _cTDCTime_depth.initialize(_name,
0220                                "TDCTime",
0221                                hcaldqm::hashfunctions::fdepth,
0222                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fieta),
0223                                new hcaldqm::quantity::DetectorQuantity(hcaldqm::quantity::fiphi),
0224                                new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fTime_ns_250),
0225                                0);
0226     _LED_ADCvsBX_Subdet.initialize(_name,
0227                                    "LED_ADCvsBX",
0228                                    hcaldqm::hashfunctions::fSubdet,
0229                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fBX_36),
0230                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_256_4),
0231                                    new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0232                                    0);
0233   } else if (_ptype == fLocal) {
0234     _LED_ADCvsEvn_Subdet.initialize(_name,
0235                                     "LED_ADCvsEvn",
0236                                     hcaldqm::hashfunctions::fSubdet,
0237                                     new hcaldqm::quantity::EventNumber(_nevents),
0238                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fADC_256_4),
0239                                     new hcaldqm::quantity::ValueQuantity(hcaldqm::quantity::fN),
0240                                     0);
0241   }
0242 
0243   //    initialize compact containers
0244   _xSignalSum.initialize(hcaldqm::hashfunctions::fDChannel);
0245   _xSignalSum2.initialize(hcaldqm::hashfunctions::fDChannel);
0246   _xTimingSum.initialize(hcaldqm::hashfunctions::fDChannel);
0247   _xTimingSum2.initialize(hcaldqm::hashfunctions::fDChannel);
0248   _xEntries.initialize(hcaldqm::hashfunctions::fDChannel);
0249 
0250   //    BOOK
0251   _cSignalMean_Subdet.book(ib, _emap, _subsystem);
0252   _cSignalRMS_Subdet.book(ib, _emap, _subsystem);
0253   _cTimingMean_Subdet.book(ib, _emap, _subsystem);
0254   _cTimingRMS_Subdet.book(ib, _emap, _subsystem);
0255 
0256   _cSignalMean_depth.book(ib, _emap, _subsystem);
0257   _cSignalRMS_depth.book(ib, _emap, _subsystem);
0258   _cTimingMean_depth.book(ib, _emap, _subsystem);
0259   _cTimingRMS_depth.book(ib, _emap, _subsystem);
0260 
0261   _cMissing_depth.book(ib, _emap, _subsystem);
0262   if (_ptype == fLocal) {  // hidefed2crate
0263     _cSignalMean_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0264     _cSignalRMS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0265     _cTimingMean_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0266     _cTimingRMS_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0267     _cShapeCut_FEDSlot.book(ib, _emap, _subsystem);
0268     _cMissing_FEDuTCA.book(ib, _emap, _filter_VME, _subsystem);
0269   }
0270   if (_ptype == fOnline || _ptype == fLocal) {
0271     _cADCvsTS_SubdetPM.book(ib, _emap, _subsystem);
0272   }
0273   if (_ptype == fOnline) {
0274     _cLowSignal_CrateSlot.book(ib, _subsystem);
0275     _cSumQ_SubdetPM.book(ib, _emap, _subsystem);
0276     _cTDCTime_SubdetPM.book(ib, _emap, _subsystem);
0277     _cTDCTime_depth.book(ib, _emap, _subsystem);
0278   }
0279 
0280   if (_ptype == fOnline) {
0281     _LED_ADCvsBX_Subdet.book(ib, _emap, _subsystem);
0282   } else if (_ptype == fLocal) {
0283     _LED_ADCvsEvn_Subdet.book(ib, _emap, _subsystem);
0284   }
0285 
0286   _xSignalSum.book(_emap);
0287   _xSignalSum2.book(_emap);
0288   _xEntries.book(_emap);
0289   _xTimingSum.book(_emap);
0290   _xTimingSum2.book(_emap);
0291 
0292   _ehashmap.initialize(_emap, electronicsmap::fD2EHashMap);
0293 }
0294 
0295 /* virtual */ void LEDTask::_resetMonitors(hcaldqm::UpdateFreq uf) { DQTask::_resetMonitors(uf); }
0296 
0297 /* virtual */ void LEDTask::_dump() {
0298   _cSignalMean_Subdet.reset();
0299   _cSignalRMS_Subdet.reset();
0300   _cTimingMean_Subdet.reset();
0301   _cTimingRMS_Subdet.reset();
0302   _cSignalMean_depth.reset();
0303   _cSignalRMS_depth.reset();
0304   _cTimingMean_depth.reset();
0305   _cTimingRMS_depth.reset();
0306 
0307   if (_ptype == fLocal) {  // hidefed2crate
0308     _cSignalMean_FEDuTCA.reset();
0309     _cSignalRMS_FEDuTCA.reset();
0310     _cTimingMean_FEDuTCA.reset();
0311     _cTimingRMS_FEDuTCA.reset();
0312   }
0313 
0314   std::vector<HcalGenericDetId> dids = _emap->allPrecisionId();
0315   for (std::vector<HcalGenericDetId>::const_iterator it = dids.begin(); it != dids.end(); ++it) {
0316     if (!it->isHcalDetId())
0317       continue;
0318     HcalDetId did = HcalDetId(it->rawId());
0319     HcalElectronicsId eid(_ehashmap.lookup(did));
0320     int n = _xEntries.get(did);
0321     double msig = _xSignalSum.get(did) / n;
0322     double mtim = _xTimingSum.get(did) / n;
0323     double rsig = sqrt(_xSignalSum2.get(did) / n - msig * msig);
0324     double rtim = sqrt(_xTimingSum2.get(did) / n - mtim * mtim);
0325 
0326     //  channels missing or low signal
0327     if (n == 0) {
0328       _cMissing_depth.fill(did);
0329       if (_ptype == fLocal) {  // hidefed2crate
0330         if (!eid.isVMEid())
0331           _cMissing_FEDuTCA.fill(eid);
0332       }
0333       continue;
0334     }
0335     _cSignalMean_Subdet.fill(did, msig);
0336     _cSignalMean_depth.fill(did, msig);
0337     _cSignalRMS_Subdet.fill(did, rsig);
0338     _cSignalRMS_depth.fill(did, rsig);
0339     _cTimingMean_Subdet.fill(did, mtim);
0340     _cTimingMean_depth.fill(did, mtim);
0341     _cTimingRMS_Subdet.fill(did, rtim);
0342     _cTimingRMS_depth.fill(did, rtim);
0343     if (_ptype == fLocal) {  // hidefed2crate
0344       if (!eid.isVMEid()) {
0345         _cSignalMean_FEDuTCA.fill(eid, msig);
0346         _cSignalRMS_FEDuTCA.fill(eid, rsig);
0347         _cTimingMean_FEDuTCA.fill(eid, mtim);
0348         _cTimingRMS_FEDuTCA.fill(eid, rtim);
0349       }
0350     }
0351   }
0352 }
0353 
0354 /* virtual */ void LEDTask::_process(edm::Event const& e, edm::EventSetup const& es) {
0355   edm::Handle<HODigiCollection> c_ho;
0356   edm::Handle<QIE10DigiCollection> c_QIE10;
0357   edm::Handle<QIE11DigiCollection> c_QIE11;
0358 
0359   if (!e.getByToken(_tokHO, c_ho))
0360     _logger.dqmthrow("Collection HODigiCollection isn't available " + _tagHO.label() + " " + _tagHO.instance());
0361   if (!e.getByToken(_tokQIE10, c_QIE10))
0362     _logger.dqmthrow("Collection QIE10DigiCollection isn't available " + _tagQIE10.label() + " " +
0363                      _tagQIE10.instance());
0364   if (!e.getByToken(_tokQIE11, c_QIE11))
0365     _logger.dqmthrow("Collection QIE11DigiCollection isn't available " + _tagQIE11.label() + " " +
0366                      _tagQIE11.instance());
0367 
0368   //    int currentEvent = e.eventAuxiliary().id().event();
0369 
0370   for (QIE11DigiCollection::const_iterator it = c_QIE11->begin(); it != c_QIE11->end(); ++it) {
0371     const QIE11DataFrame digi = static_cast<const QIE11DataFrame>(*it);
0372     HcalDetId const& did = digi.detid();
0373     if ((did.subdet() != HcalBarrel) && (did.subdet() != HcalEndcap)) {
0374       // LED monitoring from calibration channels
0375       if (did.subdet() == HcalOther) {
0376         HcalOtherDetId hodid(digi.detid());
0377         if (hodid.subdet() == HcalCalibration) {
0378           if (std::find(_ledCalibrationChannels[HcalEndcap].begin(), _ledCalibrationChannels[HcalEndcap].end(), did) !=
0379               _ledCalibrationChannels[HcalEndcap].end()) {
0380             for (int i = 0; i < digi.samples(); i++) {
0381               if (_ptype == fOnline) {
0382                 _LED_ADCvsBX_Subdet.fill(HcalDetId(HcalEndcap, 16, 1, 1), e.bunchCrossing(), digi[i].adc());
0383               } else if (_ptype == fLocal) {
0384                 _LED_ADCvsEvn_Subdet.fill(
0385                     HcalDetId(HcalEndcap, 16, 1, 1), e.eventAuxiliary().id().event(), digi[i].adc());
0386               }
0387             }
0388           } else if (std::find(_ledCalibrationChannels[HcalBarrel].begin(),
0389                                _ledCalibrationChannels[HcalBarrel].end(),
0390                                did) != _ledCalibrationChannels[HcalBarrel].end()) {
0391             for (int i = 0; i < digi.samples(); i++) {
0392               if (_ptype == fOnline) {
0393                 _LED_ADCvsBX_Subdet.fill(HcalDetId(HcalBarrel, 1, 1, 1), e.bunchCrossing(), digi[i].adc());
0394               } else if (_ptype == fLocal) {
0395                 _LED_ADCvsEvn_Subdet.fill(
0396                     HcalDetId(HcalBarrel, 1, 1, 1), e.eventAuxiliary().id().event(), digi[i].adc());
0397               }
0398             }
0399           }
0400         }
0401       }
0402       continue;
0403     }
0404     uint32_t rawid = _ehashmap.lookup(did);
0405     if (!rawid) {
0406       std::string unknown_id_string = "Detid " + std::to_string(int(did)) + ", ieta " + std::to_string(did.ieta());
0407       unknown_id_string += ", iphi " + std::to_string(did.iphi()) + ", depth " + std::to_string(did.depth());
0408       unknown_id_string += ", is not in emap. Skipping.";
0409       _logger.warn(unknown_id_string);
0410       continue;
0411     }
0412     HcalElectronicsId const& eid(rawid);
0413 
0414     CaloSamples digi_fC = hcaldqm::utilities::loadADC2fCDB<QIE11DataFrame>(_dbService, did, digi);
0415     //double sumQ = hcaldqm::utilities::sumQ_v10<QIE11DataFrame>(digi, 2.5, 0, digi.samples()-1);
0416     double sumQ = hcaldqm::utilities::sumQDB<QIE11DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0417     if (sumQ >= _lowHBHE) {
0418       //double aveTS = hcaldqm::utilities::aveTS_v10<QIE11DataFrame>(digi, 2.5, 0,digi.samples()-1);
0419       double aveTS = hcaldqm::utilities::aveTSDB<QIE11DataFrame>(_dbService, digi_fC, did, digi, 0, digi.size() - 1);
0420 
0421       _xSignalSum.get(did) += sumQ;
0422       _xSignalSum2.get(did) += sumQ * sumQ;
0423       _xTimingSum.get(did) += aveTS;
0424       _xTimingSum2.get(did) += aveTS * aveTS;
0425       _xEntries.get(did)++;
0426 
0427       if (_ptype == fLocal) {  // hidefed2crate
0428         for (int i = 0; i < digi.samples(); i++) {
0429           //_cShapeCut_FEDSlot.fill(eid, i, digi.sample(i).nominal_fC()-2.5);
0430           _cShapeCut_FEDSlot.fill(
0431               eid, i, hcaldqm::utilities::adc2fCDBMinusPedestal<QIE11DataFrame>(_dbService, digi_fC, did, digi, i));
0432         }
0433       }
0434       if (_ptype == fOnline || _ptype == fLocal) {
0435         for (int iTS = 0; iTS < digi.samples(); ++iTS) {
0436           _cADCvsTS_SubdetPM.fill(did, iTS, digi[iTS].adc());
0437         }
0438       }
0439       if (_ptype == fOnline) {
0440         for (int iTS = 0; iTS < digi.samples(); ++iTS) {
0441           if (digi[iTS].tdc() < 50) {
0442             double time = iTS * 25. + (digi[iTS].tdc() / 2.);
0443             _cTDCTime_SubdetPM.fill(did, time);
0444             _cTDCTime_depth.fill(did, time);
0445           }
0446         }
0447         _cSumQ_SubdetPM.fill(did, sumQ);
0448 
0449         // Low signal in SOI
0450         short soi = -1;
0451         for (int i = 0; i < digi.samples(); i++) {
0452           if (digi[i].soi()) {
0453             soi = i;
0454             break;
0455           }
0456         }
0457         if (digi[soi].adc() < 30) {
0458           _cLowSignal_CrateSlot.fill(eid);
0459         }
0460       }
0461     }
0462   }
0463   for (HODigiCollection::const_iterator it = c_ho->begin(); it != c_ho->end(); ++it) {
0464     const HODataFrame digi = (const HODataFrame)(*it);
0465     HcalDetId did = digi.id();
0466     if (did.subdet() != HcalOuter) {
0467       // LED monitoring from calibration channels
0468       if (did.subdet() == HcalOther) {
0469         HcalOtherDetId hodid(did);
0470         if (hodid.subdet() == HcalCalibration) {
0471           if (std::find(_ledCalibrationChannels[HcalOuter].begin(), _ledCalibrationChannels[HcalOuter].end(), did) !=
0472               _ledCalibrationChannels[HcalOuter].end()) {
0473             for (int i = 0; i < digi.size(); i++) {
0474               if (_ptype == fOnline) {
0475                 _LED_ADCvsBX_Subdet.fill(HcalDetId(HcalOuter, 1, 1, 4), e.bunchCrossing(), digi[i].adc());
0476               } else if (_ptype == fLocal) {
0477                 _LED_ADCvsEvn_Subdet.fill(
0478                     HcalDetId(HcalOuter, 1, 1, 4), e.eventAuxiliary().id().event(), digi[i].adc());
0479               }
0480             }
0481           }
0482         }
0483       }
0484       continue;
0485     }
0486     HcalElectronicsId eid = digi.elecId();
0487     //double sumQ = hcaldqm::utilities::sumQ<HODataFrame>(digi, 8.5, 0, digi.size()-1);
0488     CaloSamples digi_fC = hcaldqm::utilities::loadADC2fCDB<HODataFrame>(_dbService, did, digi);
0489     double sumQ = hcaldqm::utilities::sumQDB<HODataFrame>(_dbService, digi_fC, did, digi, 0, digi.size() - 1);
0490     if (sumQ >= _lowHO) {
0491       //double aveTS = hcaldqm::utilities::aveTS<HODataFrame>(digi, 8.5, 0, digi.size()-1);
0492       double aveTS = hcaldqm::utilities::aveTSDB<HODataFrame>(_dbService, digi_fC, did, digi, 0, digi.size() - 1);
0493 
0494       _xSignalSum.get(did) += sumQ;
0495       _xSignalSum2.get(did) += sumQ * sumQ;
0496       _xTimingSum.get(did) += aveTS;
0497       _xTimingSum2.get(did) += aveTS * aveTS;
0498       _xEntries.get(did)++;
0499 
0500       if (_ptype == fLocal) {  // hidefed2crate
0501         for (int i = 0; i < digi.size(); i++) {
0502           //_cShapeCut_FEDSlot.fill(eid, i, digi.sample(i).nominal_fC()-8.5);
0503           _cShapeCut_FEDSlot.fill(
0504               eid, i, hcaldqm::utilities::adc2fCDBMinusPedestal<HODataFrame>(_dbService, digi_fC, did, digi, i));
0505         }
0506       }
0507       if (_ptype == fOnline || _ptype == fLocal) {
0508         for (int iTS = 0; iTS < digi.size(); ++iTS) {
0509           _cADCvsTS_SubdetPM.fill(did, iTS, digi.sample(iTS).adc());
0510         }
0511       }
0512       if (_ptype == fOnline) {
0513         _cSumQ_SubdetPM.fill(did, sumQ);
0514       }
0515     }
0516   }
0517 
0518   for (QIE10DigiCollection::const_iterator it = c_QIE10->begin(); it != c_QIE10->end(); ++it) {
0519     const QIE10DataFrame digi = static_cast<const QIE10DataFrame>(*it);
0520     HcalDetId did = digi.detid();
0521     if (did.subdet() != HcalForward) {
0522       // LED monitoring from calibration channels
0523       if (did.subdet() == HcalOther) {
0524         HcalOtherDetId hodid(digi.detid());
0525         if (hodid.subdet() == HcalCalibration) {
0526           if (std::find(_ledCalibrationChannels[HcalForward].begin(),
0527                         _ledCalibrationChannels[HcalForward].end(),
0528                         did) != _ledCalibrationChannels[HcalForward].end()) {
0529             for (int i = 0; i < digi.samples(); i++) {
0530               if (_ptype == fOnline) {
0531                 _LED_ADCvsBX_Subdet.fill(HcalDetId(HcalForward, 29, 1, 1), e.bunchCrossing(), digi[i].adc());
0532               } else if (_ptype == fLocal) {
0533                 _LED_ADCvsEvn_Subdet.fill(
0534                     HcalDetId(HcalForward, 29, 1, 1), e.eventAuxiliary().id().event(), digi[i].adc());
0535               }
0536             }
0537           }
0538         }
0539       }
0540       continue;
0541     }
0542     HcalElectronicsId eid = HcalElectronicsId(_ehashmap.lookup(did));
0543     //double sumQ = hcaldqm::utilities::sumQ_v10<QIE10DataFrame>(digi, 2.5, 0, digi.samples()-1);
0544     CaloSamples digi_fC = hcaldqm::utilities::loadADC2fCDB<QIE10DataFrame>(_dbService, did, digi);
0545     double sumQ = hcaldqm::utilities::sumQDB<QIE10DataFrame>(_dbService, digi_fC, did, digi, 0, digi.samples() - 1);
0546     if (sumQ >= _lowHF) {
0547       //double aveTS = hcaldqm::utilities::aveTS_v10<QIE10DataFrame>(digi, 2.5, 0, digi.samples()-1);
0548       double aveTS = hcaldqm::utilities::aveTSDB<QIE10DataFrame>(_dbService, digi_fC, did, digi, 0, digi.size() - 1);
0549 
0550       _xSignalSum.get(did) += sumQ;
0551       _xSignalSum2.get(did) += sumQ * sumQ;
0552       _xTimingSum.get(did) += aveTS;
0553       _xTimingSum2.get(did) += aveTS * aveTS;
0554       _xEntries.get(did)++;
0555 
0556       if (_ptype == fLocal) {  // hidefed2crate
0557         for (int i = 0; i < digi.samples(); ++i) {
0558           // Note: this used to be digi.sample(i).nominal_fC() - 2.5, but this branch doesn't exist in QIE10DataFrame.
0559           // Instead, use lookup table.
0560           //_cShapeCut_FEDSlot.fill(eid, i, constants::adc2fC[digi[i].adc()]);
0561           _cShapeCut_FEDSlot.fill(
0562               eid, i, hcaldqm::utilities::adc2fCDBMinusPedestal<QIE10DataFrame>(_dbService, digi_fC, did, digi, i));
0563         }
0564       }
0565       if (_ptype == fOnline || _ptype == fLocal) {
0566         for (int iTS = 0; iTS < digi.samples(); ++iTS) {
0567           _cADCvsTS_SubdetPM.fill(did, iTS, digi[iTS].adc());
0568         }
0569       }
0570       if (_ptype == fOnline) {
0571         for (int iTS = 0; iTS < digi.samples(); ++iTS) {
0572           if (digi[iTS].le_tdc() < 50) {
0573             double time = iTS * 25. + (digi[iTS].le_tdc() / 2.);
0574             _cTDCTime_SubdetPM.fill(did, time);
0575             _cTDCTime_depth.fill(did, time);
0576           }
0577         }
0578         _cSumQ_SubdetPM.fill(did, sumQ);
0579       }
0580     }
0581   }
0582 
0583   if (_ptype == fOnline && _evsTotal > 0 && _evsTotal % constants::CALIBEVENTS_MIN == 0)
0584     this->_dump();
0585 }
0586 
0587 /* virtual */ bool LEDTask::_isApplicable(edm::Event const& e) {
0588   if (_ptype != fOnline) {
0589     //  local
0590     edm::Handle<HcalTBTriggerData> ctrigger;
0591     if (!e.getByToken(_tokTrigger, ctrigger))
0592       _logger.dqmthrow("Collection HcalTBTriggerData isn't available " + _tagTrigger.label() + " " +
0593                        _tagTrigger.instance());
0594     return ctrigger->wasLEDTrigger();
0595   } else {
0596     //  fOnline mode
0597     edm::Handle<HcalUMNioDigi> cumn;
0598     if (!e.getByToken(_tokuMN, cumn)) {
0599       return false;
0600     }
0601 
0602     return (cumn->eventType() == constants::EVENTTYPE_LED);
0603   }
0604 
0605   return false;
0606 }
0607 
0608 DEFINE_FWK_MODULE(LEDTask);