Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-18 23:17:39

0001 #include "DQM/EcalMonitorClient/interface/LedClient.h"
0002 
0003 #include "DataFormats/EcalDetId/interface/EcalPnDiodeDetId.h"
0004 
0005 #include "CondFormats/EcalObjects/interface/EcalDQMStatusHelper.h"
0006 
0007 #include "DQM/EcalCommon/interface/EcalDQMCommonUtils.h"
0008 #include "DQM/EcalCommon/interface/MESetMulti.h"
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include <cmath>
0013 #include <fstream>
0014 namespace ecaldqm {
0015   LedClient::LedClient()
0016       : DQWorkerClient(),
0017         wlToME_(),
0018         minChannelEntries_(0),
0019         expectedAmplitude_(0),
0020         toleranceAmplitude_(0.),
0021         toleranceAmpRMSRatio_(0.),
0022         expectedTiming_(0),
0023         toleranceTiming_(0.),
0024         toleranceTimRMS_(0.),
0025         expectedPNAmplitude_(0),
0026         tolerancePNAmp_(0.),
0027         tolerancePNRMSRatio_(0.),
0028         forwardFactor_(0.) {}
0029 
0030   void LedClient::setParams(edm::ParameterSet const& _params) {
0031     minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
0032     toleranceAmplitude_ = _params.getUntrackedParameter<double>("toleranceAmplitude");
0033     toleranceAmpRMSRatio_ = _params.getUntrackedParameter<double>("toleranceAmpRMSRatio");
0034     toleranceTiming_ = _params.getUntrackedParameter<double>("toleranceTiming");
0035     toleranceTimRMS_ = _params.getUntrackedParameter<double>("toleranceTimRMS");
0036     tolerancePNAmp_ = _params.getUntrackedParameter<double>("tolerancePNAmp");
0037     tolerancePNRMSRatio_ = _params.getUntrackedParameter<double>("tolerancePNRMSRatio");
0038     forwardFactor_ = _params.getUntrackedParameter<double>("forwardFactor");
0039 
0040     std::vector<int> ledWavelengths(_params.getUntrackedParameter<std::vector<int> >("ledWavelengths"));
0041 
0042     // wavelengths are not necessarily ordered
0043     // create a map wl -> MESet index
0044     // using Amplitude here but any multi-wavelength plot is fine
0045 
0046     MESet::PathReplacements repl;
0047 
0048     MESetMulti const& amplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
0049     unsigned nWL(ledWavelengths.size());
0050     for (unsigned iWL(0); iWL != nWL; ++iWL) {
0051       int wl(ledWavelengths[iWL]);
0052       if (wl != 1 && wl != 2)
0053         throw cms::Exception("InvalidConfiguration") << "Led Wavelength";
0054       repl["wl"] = std::to_string(wl);
0055       wlToME_[wl] = amplitude.getIndex(repl);
0056     }
0057 
0058     expectedAmplitude_.resize(nWL);
0059     expectedTiming_.resize(nWL);
0060     expectedPNAmplitude_.resize(nWL);
0061 
0062     std::vector<double> inExpectedAmplitude(_params.getUntrackedParameter<std::vector<double> >("expectedAmplitude"));
0063     std::vector<double> inExpectedTiming(_params.getUntrackedParameter<std::vector<double> >("expectedTiming"));
0064     std::vector<double> inExpectedPNAmplitude(
0065         _params.getUntrackedParameter<std::vector<double> >("expectedPNAmplitude"));
0066 
0067     for (std::map<int, unsigned>::iterator wlItr(wlToME_.begin()); wlItr != wlToME_.end(); ++wlItr) {
0068       unsigned iME(wlItr->second);
0069       int iWL(wlItr->first - 1);
0070       expectedAmplitude_[iME] = inExpectedAmplitude[iWL];
0071       expectedTiming_[iME] = inExpectedTiming[iWL];
0072       expectedPNAmplitude_[iME] = inExpectedPNAmplitude[iWL];
0073     }
0074 
0075     //Get the list of known problematic Supercrystal ids and store them in the vector SClist_
0076     std::string SClistpath = edm::FileInPath("DQM/EcalMonitorClient/data/LedTowers/SClist.dat").fullPath();
0077     std::ifstream infile;
0078     infile.open((SClistpath).c_str());
0079     uint32_t detid;
0080     int ix, iy, iz;
0081     while (!infile.eof()) {
0082       infile >> ix >> iy >> iz >> detid;
0083       SClist_.push_back(detid);
0084     }
0085 
0086     qualitySummaries_.insert("Quality");
0087     qualitySummaries_.insert("QualitySummary");
0088     qualitySummaries_.insert("PNQualitySummary");
0089   }
0090 
0091   void LedClient::producePlots(ProcessType) {
0092     uint32_t mask(1 << EcalDQMStatusHelper::LED_MEAN_ERROR | 1 << EcalDQMStatusHelper::LED_RMS_ERROR |
0093                   1 << EcalDQMStatusHelper::LED_TIMING_MEAN_ERROR | 1 << EcalDQMStatusHelper::LED_TIMING_RMS_ERROR);
0094 
0095     MESetMulti& meQuality(static_cast<MESetMulti&>(MEs_.at("Quality")));
0096     MESetMulti& meQualitySummary(static_cast<MESetMulti&>(MEs_.at("QualitySummary")));
0097     MESetMulti& meAmplitudeMean(static_cast<MESetMulti&>(MEs_.at("AmplitudeMean")));
0098     MESetMulti& meAmplitudeRMS(static_cast<MESetMulti&>(MEs_.at("AmplitudeRMS")));
0099     MESetMulti& meTimingMean(static_cast<MESetMulti&>(MEs_.at("TimingMean")));
0100     MESetMulti& meTimingRMSMap(static_cast<MESetMulti&>(MEs_.at("TimingRMSMap")));
0101     MESetMulti& mePNQualitySummary(static_cast<MESetMulti&>(MEs_.at("PNQualitySummary")));
0102 
0103     MESetMulti const& sAmplitude(static_cast<MESetMulti const&>(sources_.at("Amplitude")));
0104     MESetMulti const& sTiming(static_cast<MESetMulti const&>(sources_.at("Timing")));
0105     MESetMulti const& sPNAmplitude(static_cast<MESetMulti const&>(sources_.at("PNAmplitude")));
0106     MESet const& sCalibStatus(static_cast<MESet const&>(sources_.at("CalibStatus")));
0107 
0108     for (std::map<int, unsigned>::iterator wlItr(wlToME_.begin()); wlItr != wlToME_.end(); ++wlItr) {
0109       meQuality.use(wlItr->second);
0110       meQualitySummary.use(wlItr->second);
0111       meAmplitudeMean.use(wlItr->second);
0112       meAmplitudeRMS.use(wlItr->second);
0113       meTimingMean.use(wlItr->second);
0114       meTimingRMSMap.use(wlItr->second);
0115       mePNQualitySummary.use(wlItr->second);
0116 
0117       sAmplitude.use(wlItr->second);
0118       sTiming.use(wlItr->second);
0119       sPNAmplitude.use(wlItr->second);
0120 
0121       MESet::iterator qEnd(meQuality.end(GetElectronicsMap()));
0122 
0123       MESet::const_iterator tItr(GetElectronicsMap(), sTiming);
0124       MESet::const_iterator aItr(GetElectronicsMap(), sAmplitude);
0125 
0126       int wl(wlItr->first + 3);
0127       bool enabled(wl < 0 ? false : sCalibStatus.getBinContent(getEcalDQMSetupObjects(), wl) > 0 ? true : false);
0128       for (MESet::iterator qItr(meQuality.beginChannel(GetElectronicsMap())); qItr != qEnd;
0129            qItr.toNextChannel(GetElectronicsMap())) {
0130         DetId id(qItr->getId());
0131 
0132         bool doMask(meQuality.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
0133 
0134         aItr = qItr;
0135 
0136         float aEntries(aItr->getBinEntries());
0137 
0138         if (aEntries < minChannelEntries_) {
0139           qItr->setBinContent(enabled ? (doMask ? kMUnknown : kUnknown) : kMUnknown);
0140           continue;
0141         }
0142 
0143         float aMean(aItr->getBinContent());
0144         float aRms(aItr->getBinError() * sqrt(aEntries));
0145 
0146         meAmplitudeMean.fill(getEcalDQMSetupObjects(), id, aMean);
0147         meAmplitudeRMS.setBinContent(getEcalDQMSetupObjects(), id, aRms);
0148 
0149         tItr = qItr;
0150 
0151         float tEntries(tItr->getBinEntries());
0152 
0153         if (tEntries < minChannelEntries_)
0154           continue;
0155 
0156         float tMean(tItr->getBinContent());
0157         float tRms(tItr->getBinError() * sqrt(tEntries));
0158 
0159         meTimingMean.fill(getEcalDQMSetupObjects(), id, tMean);
0160         meTimingRMSMap.setBinContent(getEcalDQMSetupObjects(), id, tRms);
0161 
0162         //Temporarily disabling all cuts on LED Quality plot.
0163         qItr->setBinContent(doMask ? kMGood : kGood);
0164 
0165         /*
0166         float intensity(aMean / expectedAmplitude_[wlItr->second]);
0167         if (isForward(id))
0168           intensity /= forwardFactor_;
0169 
0170         float aRmsThr(sqrt(pow(aMean * toleranceAmpRMSRatio_, 2) + pow(3., 2)));
0171 
0172         EcalScDetId scid = EEDetId(id).sc();  //Get the Endcap SC id for the given crystal id.
0173 
0174         //For the known bad Supercrystals in the SClist, bad quality flag is only set based on the amplitude RMS
0175         //and everything else is ignored.
0176         if (std::find(SClist_.begin(), SClist_.end(), int(scid)) != SClist_.end()) {
0177           if (aRms > aRmsThr)
0178             qItr->setBinContent(doMask ? kMBad : kBad);
0179           else
0180             qItr->setBinContent(doMask ? kMGood : kGood);
0181         } else {
0182           if (intensity < toleranceAmplitude_ || aRms > aRmsThr ||
0183               std::abs(tMean - expectedTiming_[wlItr->second]) > toleranceTiming_ || tRms > toleranceTimRMS_)
0184             qItr->setBinContent(doMask ? kMBad : kBad);
0185           else
0186             qItr->setBinContent(doMask ? kMGood : kGood);
0187         }*/
0188       }
0189 
0190       towerAverage_(meQualitySummary, meQuality, 0.2);
0191 
0192       for (unsigned iDCC(0); iDCC < nDCC; ++iDCC) {
0193         if (memDCCIndex(iDCC + 1) == unsigned(-1))
0194           continue;
0195         if (iDCC >= kEBmLow && iDCC <= kEBpHigh)
0196           continue;
0197 
0198         for (unsigned iPN(0); iPN < 10; ++iPN) {
0199           EcalPnDiodeDetId id(EcalEndcap, iDCC + 1, iPN + 1);
0200 
0201           bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
0202 
0203           float pEntries(sPNAmplitude.getBinEntries(getEcalDQMSetupObjects(), id));
0204 
0205           if (pEntries < minChannelEntries_) {
0206             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMUnknown : kUnknown);
0207             continue;
0208           }
0209 
0210           float pMean(sPNAmplitude.getBinContent(getEcalDQMSetupObjects(), id));
0211           float pRms(sPNAmplitude.getBinError(getEcalDQMSetupObjects(), id) * sqrt(pEntries));
0212           float intensity(pMean / expectedPNAmplitude_[wlItr->second]);
0213 
0214           if (intensity < tolerancePNAmp_ || pRms > pMean * tolerancePNRMSRatio_)
0215             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMBad : kBad);
0216           else
0217             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMGood : kGood);
0218         }
0219       }
0220     }
0221   }
0222 
0223   DEFINE_ECALDQM_WORKER(LedClient);
0224 }  // namespace ecaldqm