File indexing completed on 2023-10-25 10:00:47
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021 #include <iostream>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/MakerMacros.h"
0026 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0027 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0029
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "DataFormats/Common/interface/DetSet.h"
0034 #include "DataFormats/Common/interface/DetSetVector.h"
0035 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0036
0037 #include "DataFormats/SiStripDigi/interface/SiStripProcessedRawDigi.h"
0038 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
0039 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0040 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
0041
0042 #include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
0043 #include "CondFormats/DataRecord/interface/SiStripPedestalsRcd.h"
0044
0045 #include "FWCore/ServiceRegistry/interface/Service.h"
0046 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0047 #include "CommonTools/Utils/interface/TFileDirectory.h"
0048
0049 #include "DataFormats/Provenance/interface/RunLumiEventNumber.h"
0050
0051
0052 #include "TH1F.h"
0053
0054
0055
0056
0057
0058 class SiStripHybridFormatAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0059 public:
0060 explicit SiStripHybridFormatAnalyzer(const edm::ParameterSet&);
0061 ~SiStripHybridFormatAnalyzer() override;
0062 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0063
0064 private:
0065 void analyze(const edm::Event&, const edm::EventSetup&) override;
0066 void beginRun(edm::Run const& iEvent, edm::EventSetup const&) override;
0067 void endRun(edm::Run const& iEvent, edm::EventSetup const&) override;
0068
0069 edm::EDGetTokenT<edm::DetSetVector<SiStripDigi>> srcDigis_;
0070 edm::EDGetTokenT<edm::DetSetVector<SiStripProcessedRawDigi>> srcAPVCM_;
0071 edm::Service<TFileService> fs_;
0072 edm::ESGetToken<SiStripPedestals, SiStripPedestalsRcd> pedestalsToken_;
0073
0074 TH1F* h1Digis_;
0075 TH1F* h1APVCM_;
0076 TH1F* h1BadAPVperEvent_;
0077 TH1F* h1BadAPVperModule_;
0078 TH1F* h1BadAPVperModuleOnlyBadModule_;
0079 TH1F* h1Pedestals_;
0080
0081 TFileDirectory sdDigis_;
0082 TFileDirectory sdMisc_;
0083
0084 uint16_t nModuletoDisplay_;
0085 uint16_t actualModule_;
0086
0087 bool plotAPVCM_;
0088
0089
0090 uint32_t peds_cache_id_;
0091 };
0092
0093 SiStripHybridFormatAnalyzer::SiStripHybridFormatAnalyzer(const edm::ParameterSet& conf) {
0094 usesResource(TFileService::kSharedResource);
0095
0096 srcDigis_ = consumes<edm::DetSetVector<SiStripDigi>>(conf.getParameter<edm::InputTag>("srcDigis"));
0097 srcAPVCM_ = consumes<edm::DetSetVector<SiStripProcessedRawDigi>>(conf.getParameter<edm::InputTag>("srcAPVCM"));
0098 pedestalsToken_ = esConsumes();
0099 nModuletoDisplay_ = conf.getParameter<uint32_t>("nModuletoDisplay");
0100 plotAPVCM_ = conf.getParameter<bool>("plotAPVCM");
0101
0102 sdDigis_ = fs_->mkdir("Digis");
0103 sdMisc_ = fs_->mkdir("Miscellanea");
0104
0105 h1APVCM_ = sdMisc_.make<TH1F>("APV CM", "APV CM", 1601, -100.5, 1500.5);
0106 h1APVCM_->SetXTitle("APV CM [adc]");
0107 h1APVCM_->SetYTitle("Entries");
0108 h1APVCM_->SetLineWidth(2);
0109 h1APVCM_->SetLineStyle(2);
0110
0111 h1BadAPVperEvent_ = sdMisc_.make<TH1F>("BadAPV/Event", "BadAPV/Event", 72786, -0.5, 72785.5);
0112 h1BadAPVperEvent_->SetXTitle("# Bad APVs");
0113 h1BadAPVperEvent_->SetYTitle("Entries");
0114 h1BadAPVperEvent_->SetLineWidth(2);
0115 h1BadAPVperEvent_->SetLineStyle(2);
0116
0117 h1BadAPVperModule_ = sdMisc_.make<TH1F>("BadAPV/Module", "BadAPV/Module", 7, -0.5, 6.5);
0118 h1BadAPVperModule_->SetXTitle("# Bad APVs");
0119 h1BadAPVperModule_->SetYTitle("Entries");
0120 h1BadAPVperModule_->SetLineWidth(2);
0121 h1BadAPVperModule_->SetLineStyle(2);
0122
0123 h1BadAPVperModuleOnlyBadModule_ =
0124 sdMisc_.make<TH1F>("BadAPV/Module Only Bad Modules", "BadAPV/Module Only Bad Modules", 7, -0.5, 6.5);
0125 h1BadAPVperModuleOnlyBadModule_->SetXTitle("# Bad APVs");
0126 h1BadAPVperModuleOnlyBadModule_->SetYTitle("Entries");
0127 h1BadAPVperModuleOnlyBadModule_->SetLineWidth(2);
0128 h1BadAPVperModuleOnlyBadModule_->SetLineStyle(2);
0129
0130 h1Pedestals_ = sdMisc_.make<TH1F>("Pedestals", "Pedestals", 2048, -1023.5, 1023.5);
0131 h1Pedestals_->SetXTitle("Pedestals [adc]");
0132 h1Pedestals_->SetYTitle("Entries");
0133 h1Pedestals_->SetLineWidth(2);
0134 h1Pedestals_->SetLineStyle(2);
0135 }
0136
0137 SiStripHybridFormatAnalyzer::~SiStripHybridFormatAnalyzer() = default;
0138
0139 void SiStripHybridFormatAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0140 edm::ParameterSetDescription desc;
0141 desc.add<edm::InputTag>("srcDigis", edm::InputTag("siStripZeroSuppression", "VirginRaw"));
0142 desc.add<edm::InputTag>("srcAPVCM", edm::InputTag("siStripZeroSuppression", "APVCMVirginRaw"));
0143 desc.add<uint32_t>("nModuletoDisplay", 10000);
0144 desc.add<bool>("plotAPVCM", true);
0145 descriptions.add("siStripHybridFormatAnalyzer", desc);
0146 }
0147
0148 void SiStripHybridFormatAnalyzer::analyze(const edm::Event& e, const edm::EventSetup& es) {
0149 using namespace edm;
0150
0151
0152
0153 if (actualModule_ == 0) {
0154 const auto& pedestalsObj = es.getData(pedestalsToken_);
0155 std::vector<uint32_t> detIdV;
0156 pedestalsObj.getDetIds(detIdV);
0157 std::vector<int> pedestals;
0158 for (const auto det : detIdV) {
0159 pedestals.clear();
0160 SiStripPedestals::Range pedestalsRange = pedestalsObj.getRange(det);
0161 pedestals.resize((pedestalsRange.second - pedestalsRange.first) * 0.8);
0162 pedestalsObj.allPeds(pedestals, pedestalsRange);
0163 for (const int ped : pedestals) {
0164 h1Pedestals_->Fill(ped);
0165 }
0166 }
0167 }
0168
0169
0170
0171
0172 if (plotAPVCM_) {
0173 edm::Handle<edm::DetSetVector<SiStripProcessedRawDigi>> moduleCM;
0174 e.getByToken(srcAPVCM_, moduleCM);
0175
0176 for (const auto& set : *moduleCM) {
0177 for (const auto& itCM : set)
0178 h1APVCM_->Fill(itCM.adc());
0179 }
0180 }
0181
0182
0183
0184 uint32_t nBadAPVevent = 0;
0185
0186 edm::Handle<edm::DetSetVector<SiStripDigi>> moduleDigis;
0187 e.getByToken(srcDigis_, moduleDigis);
0188
0189 edm::DetSetVector<SiStripDigi>::const_iterator itDigiDetSetV = moduleDigis->begin();
0190 for (; itDigiDetSetV != moduleDigis->end(); ++itDigiDetSetV) {
0191 uint32_t detId = itDigiDetSetV->id;
0192 edm::RunNumber_t const run = e.id().run();
0193 edm::EventNumber_t const event = e.id().event();
0194
0195 char detIds[20];
0196 char evs[20];
0197 char runs[20];
0198
0199 if (actualModule_ < nModuletoDisplay_) {
0200 sprintf(detIds, "%ul", detId);
0201 sprintf(evs, "%llu", event);
0202 sprintf(runs, "%u", run);
0203 char* dHistoName = Form("Id_%s_run_%s_ev_%s", detIds, runs, evs);
0204 h1Digis_ = sdDigis_.make<TH1F>(dHistoName, dHistoName, 768, -0.5, 767.5);
0205 h1Digis_->SetXTitle("strip #");
0206 h1Digis_->SetYTitle("adc");
0207 h1Digis_->SetLineWidth(2);
0208 h1Digis_->SetLineStyle(2);
0209 }
0210 uint16_t stripsPerAPV[6] = {0, 0, 0, 0, 0, 0};
0211 for (const auto& itDigi : *itDigiDetSetV) {
0212 uint16_t strip = itDigi.strip();
0213 uint16_t adc = itDigi.adc();
0214 if (actualModule_ < nModuletoDisplay_)
0215 h1Digis_->Fill(strip, adc);
0216 actualModule_++;
0217
0218
0219 stripsPerAPV[strip / 128]++;
0220 }
0221
0222 uint16_t nBadAPVmodule = 0;
0223 for (uint16_t apvN = 0; apvN < 6; apvN++) {
0224 if (stripsPerAPV[apvN] > 64) {
0225 nBadAPVevent++;
0226 nBadAPVmodule++;
0227 }
0228 }
0229 h1BadAPVperModule_->Fill(nBadAPVmodule);
0230 if (nBadAPVmodule)
0231 h1BadAPVperModuleOnlyBadModule_->Fill(nBadAPVmodule);
0232 }
0233 h1BadAPVperEvent_->Fill(nBadAPVevent);
0234 }
0235
0236 void SiStripHybridFormatAnalyzer::beginRun(edm::Run const& iEvent, edm::EventSetup const&) { actualModule_ = 0; }
0237
0238 void SiStripHybridFormatAnalyzer::endRun(edm::Run const& iEvent, edm::EventSetup const&) {}
0239
0240
0241 DEFINE_FWK_MODULE(SiStripHybridFormatAnalyzer);