Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQM/EcalMonitorTasks/interface/PresampleTask.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 
0004 #include "DQM/EcalCommon/interface/EcalDQMCommonUtils.h"
0005 
0006 #include "DataFormats/EcalRawData/interface/EcalDCCHeaderBlock.h"
0007 
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 
0010 namespace ecaldqm {
0011   PresampleTask::PresampleTask()
0012       : DQWorkerTask(), doPulseMaxCheck_(true), pulseMaxPosition_(0), nSamples_(0), mePedestalByLS(nullptr) {}
0013 
0014   void PresampleTask::setParams(edm::ParameterSet const& _params) {
0015     doPulseMaxCheck_ = _params.getUntrackedParameter<bool>("doPulseMaxCheck");
0016     pulseMaxPosition_ = _params.getUntrackedParameter<int>("pulseMaxPosition");
0017     nSamples_ = _params.getUntrackedParameter<int>("nSamples");
0018   }
0019   void PresampleTask::setTokens(edm::ConsumesCollector& _collector) { Pedtoken_ = _collector.esConsumes(); }
0020 
0021   bool PresampleTask::filterRunType(short const* _runType) {
0022     for (int iFED(0); iFED < nDCC; iFED++) {
0023       if (_runType[iFED] == EcalDCCHeaderBlock::COSMIC || _runType[iFED] == EcalDCCHeaderBlock::MTCC ||
0024           _runType[iFED] == EcalDCCHeaderBlock::COSMICS_GLOBAL ||
0025           _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_GLOBAL || _runType[iFED] == EcalDCCHeaderBlock::COSMICS_LOCAL ||
0026           _runType[iFED] == EcalDCCHeaderBlock::PHYSICS_LOCAL)
0027         return true;
0028     }
0029 
0030     return false;
0031   }
0032 
0033   void PresampleTask::beginRun(edm::Run const&, edm::EventSetup const& _es) { FillPedestal = true; }
0034 
0035   void PresampleTask::beginEvent(edm::Event const& _evt,
0036                                  edm::EventSetup const& _es,
0037                                  bool const& ByLumiResetSwitch,
0038                                  bool&) {
0039     if (ByLumiResetSwitch) {
0040       // Fill separate MEs with only 10 LSs worth of stats
0041       // Used to correctly fill Presample Trend plots:
0042       // 1 pt:10 LS in Trend plots
0043       mePedestalByLS = &MEs_.at("PedestalByLS");
0044       if (timestamp_.iLumi % 10 == 0)
0045         mePedestalByLS->reset(GetElectronicsMap());
0046     }
0047 
0048     MESet& mePedestalProjEtaG1(MEs_.at("PedestalProjEtaG1"));
0049     MESet& mePedestalProjEtaG6(MEs_.at("PedestalProjEtaG6"));
0050     MESet& mePedestalProjEtaG12(MEs_.at("PedestalProjEtaG12"));
0051 
0052     if (FillPedestal) {
0053       const EcalPedestals* myped = &_es.getData(Pedtoken_);
0054 
0055       for (int i = 0; i < EBDetId::kSizeForDenseIndexing; i++) {
0056         if (!EBDetId::validDenseIndex(i))
0057           continue;
0058         EBDetId ebid(EBDetId::unhashIndex(i));
0059         EcalPedestals::const_iterator it = myped->find(ebid.rawId());
0060         if (it != myped->end()) {
0061           mePedestalProjEtaG1.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x1);
0062           mePedestalProjEtaG6.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x6);
0063           mePedestalProjEtaG12.fill(getEcalDQMSetupObjects(), ebid, (*it).rms_x12);
0064         }
0065       }
0066       for (int i = 0; i < EEDetId::kSizeForDenseIndexing; i++) {
0067         if (!EEDetId::validDenseIndex(i))
0068           continue;
0069         EEDetId eeid(EEDetId::unhashIndex(i));
0070         EcalPedestals::const_iterator it = myped->find(eeid.rawId());
0071         if (it != myped->end()) {
0072           mePedestalProjEtaG1.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x1);
0073           mePedestalProjEtaG6.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x6);
0074           mePedestalProjEtaG12.fill(getEcalDQMSetupObjects(), eeid, (*it).rms_x12);
0075         }
0076       }
0077 
0078       FillPedestal = false;
0079     }
0080   }
0081 
0082   template <typename DigiCollection>
0083   void PresampleTask::runOnDigis(DigiCollection const& _digis) {
0084     MESet& mePedestal(MEs_.at("Pedestal"));  // contains cumulative run stats => not suitable for Trend plots
0085 
0086     for (typename DigiCollection::const_iterator digiItr(_digis.begin()); digiItr != _digis.end(); ++digiItr) {
0087       DetId id(digiItr->id());
0088 
0089       // EcalDataFrame is not a derived class of edm::DataFrame, but can take edm::DataFrame in the constructor
0090       EcalDataFrame dataFrame(*digiItr);
0091 
0092       // Check that the digi pulse maximum occurs on the 6th sample
0093       // For cosmics: disable this check to preserve statistics
0094       if (doPulseMaxCheck_) {
0095         bool gainSwitch(false);
0096         int iMax(-1);
0097         int maxADC(0);
0098         for (int iSample(0); iSample < EcalDataFrame::MAXSAMPLES; ++iSample) {
0099           int adc(dataFrame.sample(iSample).adc());
0100           if (adc > maxADC) {
0101             iMax = iSample;
0102             maxADC = adc;
0103           }
0104           if (iSample < nSamples_ && dataFrame.sample(iSample).gainId() != 1) {
0105             gainSwitch = true;
0106             break;
0107           }
0108         }  // iSample
0109         if (iMax != pulseMaxPosition_ || gainSwitch)
0110           continue;
0111       }  // PulseMaxCheck
0112 
0113       for (int iSample(0); iSample < nSamples_; ++iSample) {
0114         mePedestal.fill(getEcalDQMSetupObjects(), id, double(dataFrame.sample(iSample).adc()));
0115         mePedestalByLS->fill(getEcalDQMSetupObjects(), id, double(dataFrame.sample(iSample).adc()));
0116       }
0117 
0118     }  // _digis loop
0119   }    // runOnDigis
0120 
0121   DEFINE_ECALDQM_WORKER(PresampleTask);
0122 }  // namespace ecaldqm