File indexing completed on 2023-03-17 10:56:15
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include <iostream>
0011
0012
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016
0017
0018 #include "DQMServices/Core/interface/DQMStore.h"
0019
0020
0021 #include "DataFormats/Common/interface/DetSetVector.h"
0022 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0023
0024
0025 #include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
0026
0027 namespace {
0028
0029 class SiPixelPhase1Digis final : public SiPixelPhase1Base {
0030
0031 enum {
0032 ADC,
0033 NDIGIS,
0034 NDIGISINCLUSIVE,
0035 NDIGIS_FED,
0036 NDIGIS_FEDtrend,
0037 EVENT,
0038 MAP,
0039 OCCUPANCY,
0040
0041 MAX_HIST
0042 };
0043
0044 public:
0045 explicit SiPixelPhase1Digis(const edm::ParameterSet& conf);
0046
0047 void analyze(const edm::Event&, const edm::EventSetup&) override;
0048
0049 private:
0050 edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> srcToken_;
0051 };
0052
0053 SiPixelPhase1Digis::SiPixelPhase1Digis(const edm::ParameterSet& iConfig) : SiPixelPhase1Base(iConfig) {
0054 srcToken_ = consumes<edm::DetSetVector<PixelDigi>>(iConfig.getParameter<edm::InputTag>("src"));
0055 }
0056
0057 void SiPixelPhase1Digis::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0058 if (!checktrigger(iEvent, iSetup, DCS))
0059 return;
0060
0061 edm::Handle<edm::DetSetVector<PixelDigi>> input;
0062 iEvent.getByToken(srcToken_, input);
0063 if (!input.isValid())
0064 return;
0065 bool hasDigis = false;
0066
0067 edm::DetSetVector<PixelDigi>::const_iterator it;
0068 for (it = input->begin(); it != input->end(); ++it) {
0069 for (PixelDigi const& digi : *it) {
0070 hasDigis = true;
0071 histo[ADC].fill((double)digi.adc(), DetId(it->detId()), &iEvent, digi.column(), digi.row());
0072 histo[MAP].fill(DetId(it->detId()), &iEvent, digi.column(), digi.row());
0073 histo[OCCUPANCY].fill(DetId(it->detId()), &iEvent, digi.column(), digi.row());
0074 histo[NDIGIS].fill(DetId(it->detId()), &iEvent);
0075 histo[NDIGISINCLUSIVE].fill(DetId(it->detId()), &iEvent);
0076 histo[NDIGIS_FED].fill(DetId(it->detId()), &iEvent);
0077 histo[NDIGIS_FEDtrend].fill(DetId(it->detId()), &iEvent);
0078 }
0079 }
0080 if (hasDigis)
0081 histo[EVENT].fill(DetId(0), &iEvent);
0082 histo[NDIGIS].executePerEventHarvesting(&iEvent);
0083 histo[NDIGISINCLUSIVE].executePerEventHarvesting(&iEvent);
0084 histo[NDIGIS_FED].executePerEventHarvesting(&iEvent);
0085 histo[NDIGIS_FEDtrend].executePerEventHarvesting(&iEvent);
0086 }
0087
0088 }
0089
0090 DEFINE_FWK_MODULE(SiPixelPhase1Digis);