File indexing completed on 2024-04-06 12:07:23
0001 #include <iostream>
0002 #include <fstream>
0003 #include <vector>
0004
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 #include "DQM/EcalPreshowerMonitorModule/interface/ESFEDIntegrityTask.h"
0013
0014 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0015 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0016 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0017 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0018 #include "DataFormats/EcalRawData/interface/ESDCCHeaderBlock.h"
0019 #include "DataFormats/EcalRawData/interface/ESKCHIPBlock.h"
0020
0021 using namespace cms;
0022 using namespace edm;
0023 using namespace std;
0024
0025 ESFEDIntegrityTask::ESFEDIntegrityTask(const ParameterSet& ps) {
0026 prefixME_ = ps.getUntrackedParameter<string>("prefixME", "");
0027 fedDirName_ = ps.getUntrackedParameter<string>("FEDDirName", "FEDIntegrity");
0028 debug_ = ps.getUntrackedParameter<bool>("debug", false);
0029
0030 dccCollections_ = consumes<ESRawDataCollection>(ps.getParameter<InputTag>("ESDCCCollections"));
0031 FEDRawDataCollection_ = consumes<FEDRawDataCollection>(ps.getParameter<edm::InputTag>("FEDRawDataCollection"));
0032
0033 meESFedsEntries_ = nullptr;
0034 meESFedsFatal_ = nullptr;
0035 meESFedsNonFatal_ = nullptr;
0036
0037 ievt_ = 0;
0038 }
0039
0040 void ESFEDIntegrityTask::bookHistograms(DQMStore::IBooker& iBooker, Run const&, EventSetup const&) {
0041 char histo[200];
0042
0043 iBooker.setCurrentFolder(prefixME_ + "/" + fedDirName_);
0044
0045 sprintf(histo, "FEDEntries");
0046 meESFedsEntries_ = iBooker.book1D(histo, histo, 56, 520, 576);
0047
0048 sprintf(histo, "FEDFatal");
0049 meESFedsFatal_ = iBooker.book1D(histo, histo, 56, 520, 576);
0050
0051 sprintf(histo, "FEDNonFatal");
0052 meESFedsNonFatal_ = iBooker.book1D(histo, histo, 56, 520, 576);
0053 }
0054
0055 void ESFEDIntegrityTask::analyze(const Event& e, const EventSetup& c) {
0056 ievt_++;
0057
0058 int gt_L1A = 0;
0059
0060 int esDCC_L1A_MostFreqCounts = 0;
0061 int esDCC_BX_MostFreqCounts = 0;
0062 int esDCC_OrbitNumber_MostFreqCounts = 0;
0063 int gtFedDataSize = 0;
0064
0065 Handle<ESRawDataCollection> dccs;
0066 Handle<FEDRawDataCollection> allFedRawData;
0067
0068 if (e.getByToken(FEDRawDataCollection_, allFedRawData)) {
0069
0070 for (int esFED = 520; esFED <= 575; ++esFED) {
0071 const FEDRawData& fedData = allFedRawData->FEDData(esFED);
0072 int length = fedData.size() / sizeof(uint64_t);
0073
0074 if (length > 0)
0075 if (meESFedsEntries_)
0076 meESFedsEntries_->Fill(esFED);
0077 }
0078
0079
0080 const FEDRawData& gtFedData = allFedRawData->FEDData(812);
0081
0082 gtFedDataSize = gtFedData.size() / sizeof(uint64_t);
0083
0084 if (gtFedDataSize > 0) {
0085 FEDHeader header(gtFedData.data());
0086
0087 gt_L1A = header.lvl1ID();
0088
0089
0090 } else {
0091 map<int, int> esDCC_L1A_FreqMap;
0092 map<int, int> esDCC_BX_FreqMap;
0093 map<int, int> esDCC_OrbitNumber_FreqMap;
0094
0095 if (e.getByToken(dccCollections_, dccs)) {
0096 for (ESRawDataCollection::const_iterator dccItr = dccs->begin(); dccItr != dccs->end(); ++dccItr) {
0097 const ESDCCHeaderBlock& esdcc = (*dccItr);
0098
0099 esDCC_L1A_FreqMap[esdcc.getLV1()]++;
0100 esDCC_BX_FreqMap[esdcc.getBX()]++;
0101 esDCC_OrbitNumber_FreqMap[esdcc.getOrbitNumber()]++;
0102
0103 if (esDCC_L1A_FreqMap[esdcc.getLV1()] > esDCC_L1A_MostFreqCounts) {
0104 esDCC_L1A_MostFreqCounts = esDCC_L1A_FreqMap[esdcc.getLV1()];
0105 gt_L1A = esdcc.getLV1();
0106 }
0107
0108 if (esDCC_BX_FreqMap[esdcc.getBX()] > esDCC_BX_MostFreqCounts) {
0109 esDCC_BX_MostFreqCounts = esDCC_BX_FreqMap[esdcc.getBX()];
0110
0111 }
0112
0113 if (esDCC_OrbitNumber_FreqMap[esdcc.getOrbitNumber()] > esDCC_OrbitNumber_MostFreqCounts) {
0114 esDCC_OrbitNumber_MostFreqCounts = esDCC_OrbitNumber_FreqMap[esdcc.getOrbitNumber()];
0115
0116 }
0117 }
0118 } else {
0119 LogWarning("ESFEDIntegrityTask") << "dccCollections not available";
0120 }
0121 }
0122
0123 } else {
0124 LogWarning("ESFEDIntegrityTask") << "FEDRawDataCollection not available";
0125 }
0126
0127 vector<int> fiberStatus;
0128 if (e.getByToken(dccCollections_, dccs)) {
0129 for (ESRawDataCollection::const_iterator dccItr = dccs->begin(); dccItr != dccs->end(); ++dccItr) {
0130 const ESDCCHeaderBlock& dcc = (*dccItr);
0131
0132 if (dcc.getDCCErrors() > 0) {
0133 if (meESFedsFatal_)
0134 meESFedsFatal_->Fill(dcc.fedId());
0135
0136 } else {
0137 if (debug_)
0138 cout << dcc.fedId() << " " << dcc.getOptoRX0() << " " << dcc.getOptoRX1() << " " << dcc.getOptoRX2() << endl;
0139 fiberStatus = dcc.getFEChannelStatus();
0140
0141 if (dcc.getOptoRX0() == 128) {
0142 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 3.);
0143 } else if (dcc.getOptoRX0() == 129) {
0144 for (unsigned int i = 0; i < 12; ++i) {
0145 if (fiberStatus[i] == 8 || fiberStatus[i] == 10 || fiberStatus[i] == 11 || fiberStatus[i] == 12)
0146 if (meESFedsNonFatal_)
0147 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 12.);
0148 }
0149 }
0150 if (dcc.getOptoRX1() == 128) {
0151 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 3.);
0152 } else if (dcc.getOptoRX1() == 129) {
0153 for (unsigned int i = 12; i < 24; ++i) {
0154 if (fiberStatus[i] == 8 || fiberStatus[i] == 10 || fiberStatus[i] == 11 || fiberStatus[i] == 12)
0155 if (meESFedsNonFatal_)
0156 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 12.);
0157 }
0158 }
0159 if (dcc.getOptoRX2() == 128) {
0160 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 3.);
0161 } else if (dcc.getOptoRX2() == 129) {
0162 for (unsigned int i = 24; i < 36; ++i) {
0163 if (fiberStatus[i] == 8 || fiberStatus[i] == 10 || fiberStatus[i] == 11 || fiberStatus[i] == 12)
0164 if (meESFedsNonFatal_)
0165 meESFedsNonFatal_->Fill(dcc.fedId(), 1. / 12.);
0166 }
0167 }
0168 }
0169
0170 if (dcc.getLV1() != gt_L1A)
0171 meESFedsNonFatal_->Fill(dcc.fedId());
0172
0173
0174 }
0175 }
0176
0177
0178
0179
0180
0181
0182
0183
0184
0185
0186
0187
0188 }
0189
0190 DEFINE_FWK_MODULE(ESFEDIntegrityTask);