Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQM/EcalMonitorClient/interface/PedestalClient.h"
0002 
0003 #include "DQM/EcalCommon/interface/MESetMulti.h"
0004 #include "DQM/EcalCommon/interface/EcalDQMCommonUtils.h"
0005 
0006 #include "DataFormats/EcalDetId/interface/EcalPnDiodeDetId.h"
0007 
0008 #include "CondFormats/EcalObjects/interface/EcalDQMStatusHelper.h"
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include <iomanip>
0013 
0014 namespace ecaldqm {
0015   PedestalClient::PedestalClient()
0016       : DQWorkerClient(),
0017         gainToME_(),
0018         pnGainToME_(),
0019         minChannelEntries_(0),
0020         expectedMean_(0.),
0021         toleranceMean_(0.),
0022         toleranceRMSEB_(0),
0023         toleranceRMSEE_(0),
0024         expectedPNMean_(0.),
0025         tolerancePNMean_(0.),
0026         tolerancePNRMS_(0) {}
0027 
0028   void PedestalClient::setParams(edm::ParameterSet const& _params) {
0029     minChannelEntries_ = _params.getUntrackedParameter<int>("minChannelEntries");
0030     expectedMean_ = _params.getUntrackedParameter<double>("expectedMean");
0031     toleranceMean_ = _params.getUntrackedParameter<double>("toleranceMean");
0032     expectedPNMean_ = _params.getUntrackedParameter<double>("expectedPNMean");
0033     tolerancePNMean_ = _params.getUntrackedParameter<double>("tolerancePNMean");
0034 
0035     std::vector<int> MGPAGains(_params.getUntrackedParameter<std::vector<int> >("MGPAGains"));
0036     std::vector<int> MGPAGainsPN(_params.getUntrackedParameter<std::vector<int> >("MGPAGainsPN"));
0037 
0038     MESet::PathReplacements repl;
0039 
0040     MESetMulti const& pedestal(static_cast<MESetMulti const&>(sources_.at("Pedestal")));
0041     unsigned nG(MGPAGains.size());
0042     for (unsigned iG(0); iG != nG; ++iG) {
0043       int gain(MGPAGains[iG]);
0044       if (gain != 1 && gain != 6 && gain != 12)
0045         throw cms::Exception("InvalidConfiguration") << "MGPA gain";
0046       repl["gain"] = std::to_string(gain);
0047       gainToME_[gain] = pedestal.getIndex(repl);
0048     }
0049 
0050     repl.clear();
0051 
0052     MESetMulti const& pnPedestal(static_cast<MESetMulti const&>(sources_.at("PNPedestal")));
0053     unsigned nGPN(MGPAGainsPN.size());
0054     for (unsigned iG(0); iG != nGPN; ++iG) {
0055       int gain(MGPAGainsPN[iG]);
0056       if (gain != 1 && gain != 16)
0057         throw cms::Exception("InvalidConfiguration") << "PN MGPA gain";
0058       repl["pngain"] = std::to_string(gain);
0059       pnGainToME_[gain] = pnPedestal.getIndex(repl);
0060     }
0061 
0062     toleranceRMSEB_.resize(nG);
0063     toleranceRMSEE_.resize(nG);
0064 
0065     std::vector<double> inToleranceRMSEB(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEB"));
0066     std::vector<double> inToleranceRMSEE(_params.getUntrackedParameter<std::vector<double> >("toleranceRMSEE"));
0067 
0068     for (std::map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
0069       unsigned iME(gainItr->second);
0070       unsigned iGain(0);
0071       switch (gainItr->first) {
0072         case 1:
0073           iGain = 0;
0074           break;
0075         case 6:
0076           iGain = 1;
0077           break;
0078         case 12:
0079           iGain = 2;
0080           break;
0081       }
0082 
0083       toleranceRMSEB_[iME] = inToleranceRMSEB[iGain];
0084       toleranceRMSEE_[iME] = inToleranceRMSEE[iGain];
0085     }
0086 
0087     tolerancePNRMS_.resize(nGPN);
0088 
0089     std::vector<double> inTolerancePNRMS(_params.getUntrackedParameter<std::vector<double> >("tolerancePNRMS"));
0090 
0091     for (std::map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
0092       unsigned iME(gainItr->second);
0093       unsigned iGain(0);
0094       switch (gainItr->first) {
0095         case 1:
0096           iGain = 0;
0097           break;
0098         case 16:
0099           iGain = 1;
0100           break;
0101       }
0102 
0103       tolerancePNRMS_[iME] = inTolerancePNRMS[iGain];
0104     }
0105 
0106     qualitySummaries_.insert("Quality");
0107     qualitySummaries_.insert("QualitySummary");
0108     qualitySummaries_.insert("PNQualitySummary");
0109   }
0110 
0111   void PedestalClient::producePlots(ProcessType) {
0112     using namespace std;
0113 
0114     MESetMulti& meQuality(static_cast<MESetMulti&>(MEs_.at("Quality")));
0115     MESetMulti& meQualitySummary(static_cast<MESetMulti&>(MEs_.at("QualitySummary")));
0116     MESetMulti& meMean(static_cast<MESetMulti&>(MEs_.at("Mean")));
0117     MESetMulti& meRMS(static_cast<MESetMulti&>(MEs_.at("RMS")));
0118     MESetMulti& mePNQualitySummary(static_cast<MESetMulti&>(MEs_.at("PNQualitySummary")));
0119     MESetMulti& mePNRMS(static_cast<MESetMulti&>(MEs_.at("PNRMS")));
0120 
0121     MESetMulti const& sPedestal(static_cast<MESetMulti const&>(sources_.at("Pedestal")));
0122     MESetMulti const& sPNPedestal(static_cast<MESetMulti const&>(sources_.at("PNPedestal")));
0123 
0124     for (map<int, unsigned>::iterator gainItr(gainToME_.begin()); gainItr != gainToME_.end(); ++gainItr) {
0125       meQuality.use(gainItr->second);
0126       meQualitySummary.use(gainItr->second);
0127       meMean.use(gainItr->second);
0128       meRMS.use(gainItr->second);
0129 
0130       sPedestal.use(gainItr->second);
0131 
0132       uint32_t mask(0);
0133       switch (gainItr->first) {
0134         case 1:
0135           mask |= (1 << EcalDQMStatusHelper::PEDESTAL_LOW_GAIN_MEAN_ERROR |
0136                    1 << EcalDQMStatusHelper::PEDESTAL_LOW_GAIN_RMS_ERROR);
0137           break;
0138         case 6:
0139           mask |= (1 << EcalDQMStatusHelper::PEDESTAL_MIDDLE_GAIN_MEAN_ERROR |
0140                    1 << EcalDQMStatusHelper::PEDESTAL_MIDDLE_GAIN_RMS_ERROR);
0141           break;
0142         case 12:
0143           mask |= (1 << EcalDQMStatusHelper::PEDESTAL_HIGH_GAIN_MEAN_ERROR |
0144                    1 << EcalDQMStatusHelper::PEDESTAL_HIGH_GAIN_RMS_ERROR);
0145           break;
0146         default:
0147           break;
0148       }
0149 
0150       MESet::iterator qEnd(meQuality.end(GetElectronicsMap()));
0151       MESet::const_iterator pItr(GetElectronicsMap(), sPedestal);
0152       for (MESet::iterator qItr(meQuality.beginChannel(GetElectronicsMap())); qItr != qEnd;
0153            qItr.toNextChannel(GetElectronicsMap())) {
0154         DetId id(qItr->getId());
0155 
0156         bool doMask(meQuality.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
0157 
0158         pItr = qItr;
0159 
0160         float entries(pItr->getBinEntries());
0161 
0162         if (entries < minChannelEntries_) {
0163           qItr->setBinContent(doMask ? kMUnknown : kUnknown);
0164           continue;
0165         }
0166 
0167         float mean(pItr->getBinContent());
0168         float rms(pItr->getBinError() * sqrt(entries));
0169 
0170         meMean.fill(getEcalDQMSetupObjects(), id, mean);
0171         meRMS.fill(getEcalDQMSetupObjects(), id, rms);
0172 
0173         float toleranceRMS_ =
0174             (id.subdetId() == EcalBarrel) ? toleranceRMSEB_[gainItr->second] : toleranceRMSEE_[gainItr->second];
0175 
0176         if (abs(mean - expectedMean_) > toleranceMean_ || rms > toleranceRMS_)
0177           qItr->setBinContent(doMask ? kMBad : kBad);
0178         else
0179           qItr->setBinContent(doMask ? kMGood : kGood);
0180       }
0181 
0182       towerAverage_(meQualitySummary, meQuality, 0.2);
0183     }
0184 
0185     for (map<int, unsigned>::iterator gainItr(pnGainToME_.begin()); gainItr != pnGainToME_.end(); ++gainItr) {
0186       mePNQualitySummary.use(gainItr->second);
0187       mePNRMS.use(gainItr->second);
0188 
0189       sPNPedestal.use(gainItr->second);
0190 
0191       uint32_t mask(0);
0192       switch (gainItr->first) {
0193         case 1:
0194           mask |= (1 << EcalDQMStatusHelper::PEDESTAL_LOW_GAIN_MEAN_ERROR |
0195                    1 << EcalDQMStatusHelper::PEDESTAL_LOW_GAIN_RMS_ERROR);
0196           break;
0197         case 16:
0198           mask |= (1 << EcalDQMStatusHelper::PEDESTAL_HIGH_GAIN_MEAN_ERROR |
0199                    1 << EcalDQMStatusHelper::PEDESTAL_HIGH_GAIN_RMS_ERROR);
0200           break;
0201         default:
0202           break;
0203       }
0204 
0205       for (unsigned iDCC(0); iDCC < nDCC; ++iDCC) {
0206         if (memDCCIndex(iDCC + 1) == unsigned(-1))
0207           continue;
0208 
0209         for (unsigned iPN(0); iPN < 10; ++iPN) {
0210           int subdet(0);
0211           if (iDCC >= kEBmLow && iDCC <= kEBpHigh)
0212             subdet = EcalBarrel;
0213           else
0214             subdet = EcalEndcap;
0215 
0216           EcalPnDiodeDetId id(subdet, iDCC + 1, iPN + 1);
0217 
0218           bool doMask(mePNQualitySummary.maskMatches(id, mask, statusManager_, GetTrigTowerMap()));
0219 
0220           float entries(sPNPedestal.getBinEntries(getEcalDQMSetupObjects(), id));
0221 
0222           if (entries < minChannelEntries_) {
0223             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMUnknown : kUnknown);
0224             continue;
0225           }
0226 
0227           float mean(sPNPedestal.getBinContent(getEcalDQMSetupObjects(), id));
0228           float rms(sPNPedestal.getBinError(getEcalDQMSetupObjects(), id) * sqrt(entries));
0229 
0230           mePNRMS.fill(getEcalDQMSetupObjects(), id, rms);
0231 
0232           if (abs(mean - expectedPNMean_) > tolerancePNMean_ || rms > tolerancePNRMS_[gainItr->second])
0233             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMBad : kBad);
0234           else
0235             mePNQualitySummary.setBinContent(getEcalDQMSetupObjects(), id, doMask ? kMGood : kGood);
0236         }
0237       }
0238     }
0239   }
0240 
0241   DEFINE_ECALDQM_WORKER(PedestalClient);
0242 }  // namespace ecaldqm