File indexing completed on 2023-03-17 10:48:43
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016 #include <memory>
0017 #include <fstream>
0018
0019
0020 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0021 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
0022 #include "CondFormats/SiPixelObjects/interface/SiPixelFEDChannelContainer.h"
0023 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0024 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
0025 #include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
0026 #include "FWCore/Framework/interface/ESHandle.h"
0027 #include "FWCore/Framework/interface/ESWatcher.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "FWCore/ServiceRegistry/interface/Service.h"
0034
0035
0036
0037
0038
0039 class SiPixelFEDChannelContainerWriteFromASCII : public edm::one::EDAnalyzer<> {
0040 public:
0041 explicit SiPixelFEDChannelContainerWriteFromASCII(const edm::ParameterSet&);
0042 ~SiPixelFEDChannelContainerWriteFromASCII() override;
0043
0044 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0045
0046 private:
0047 void analyze(const edm::Event&, const edm::EventSetup&) override;
0048 void endJob() override;
0049
0050
0051 const std::string m_record;
0052 const std::string m_SnapshotInputs;
0053 const bool printdebug_;
0054 const bool addDefault_;
0055 SiPixelFEDChannelContainer* myQualities;
0056 };
0057
0058
0059
0060
0061 SiPixelFEDChannelContainerWriteFromASCII::SiPixelFEDChannelContainerWriteFromASCII(const edm::ParameterSet& iConfig)
0062 : m_record(iConfig.getParameter<std::string>("record")),
0063 m_SnapshotInputs(iConfig.getParameter<std::string>("snapshots")),
0064 printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", false)),
0065 addDefault_(iConfig.getUntrackedParameter<bool>("addDefault", false)) {
0066
0067 myQualities = new SiPixelFEDChannelContainer();
0068 }
0069
0070 SiPixelFEDChannelContainerWriteFromASCII::~SiPixelFEDChannelContainerWriteFromASCII() { delete myQualities; }
0071
0072
0073
0074
0075
0076
0077 void SiPixelFEDChannelContainerWriteFromASCII::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0078 using namespace edm;
0079
0080 std::ifstream mysnapshots(m_SnapshotInputs);
0081 std::string line;
0082
0083 std::string scenario = "";
0084 unsigned int thedetid(0);
0085
0086 SiPixelFEDChannelContainer::SiPixelFEDChannelCollection theBadFEDChannels;
0087
0088 if (mysnapshots.is_open()) {
0089 while (getline(mysnapshots, line)) {
0090 if (printdebug_) {
0091 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII") << line << std::endl;
0092 }
0093 std::istringstream iss(line);
0094 unsigned int run, ls, detid, fed, link, roc_first, roc_last;
0095 iss >> run >> ls >> detid >> fed >> link >> roc_first >> roc_last;
0096
0097 PixelFEDChannel theBadChannel{fed, link, roc_first, roc_last};
0098
0099 auto newscenario = std::to_string(run) + "_" + std::to_string(ls);
0100 if (newscenario != scenario) {
0101 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII") << "================================" << std::endl;
0102 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII")
0103 << "found a new scenario: " << newscenario << std::endl;
0104 if (!scenario.empty()) {
0105 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII")
0106 << "size of the fed channel vector: " << theBadFEDChannels.size() << std::endl;
0107 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII")
0108 << "================================" << std::endl;
0109 myQualities->setScenario(scenario, theBadFEDChannels);
0110 theBadFEDChannels.clear();
0111 }
0112 scenario = newscenario;
0113 }
0114
0115 if (detid != thedetid) {
0116 if (printdebug_) {
0117 edm::LogVerbatim("SiPixelFEDChannelContainerWriteFromASCII") << "found a new detid!" << detid << std::endl;
0118 }
0119 thedetid = detid;
0120 }
0121 theBadFEDChannels[thedetid].push_back(theBadChannel);
0122 }
0123 }
0124
0125 myQualities->setScenario(scenario, theBadFEDChannels);
0126
0127 if (printdebug_) {
0128 edm::LogInfo("SiPixelFEDChannelContainerWriteFromASCII") << "Content of SiPixelFEDChannelContainer " << std::endl;
0129
0130
0131 myQualities->printAll();
0132 }
0133 }
0134
0135
0136 void SiPixelFEDChannelContainerWriteFromASCII::endJob() {
0137
0138 if (addDefault_) {
0139 SiPixelFEDChannelContainer::SiPixelFEDChannelCollection theBadFEDChannels;
0140 myQualities->setScenario("default", theBadFEDChannels);
0141 }
0142
0143 edm::LogInfo("SiPixelFEDChannelContainerWriteFromASCII")
0144 << "Size of SiPixelFEDChannelContainer object " << myQualities->size() << std::endl
0145 << std::endl;
0146
0147
0148 edm::Service<cond::service::PoolDBOutputService> poolDbService;
0149 if (poolDbService.isAvailable()) {
0150 cond::Time_t valid_time = poolDbService->currentTime();
0151
0152 poolDbService->writeOneIOV(*myQualities, valid_time, m_record);
0153 }
0154 }
0155
0156
0157 void SiPixelFEDChannelContainerWriteFromASCII::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0158 edm::ParameterSetDescription desc;
0159 desc.setComment("Writes payloads of type SiPixelFEDChannelContainer from input ASCII files");
0160 desc.addUntracked<bool>("printDebug", true);
0161 desc.addUntracked<bool>("addDefault", true);
0162 desc.add<std::string>("snapshots", "");
0163 desc.add<std::string>("record", "SiPixelStatusScenariosRcd");
0164 descriptions.add("SiPixelFEDChannelContainerWriteFromASCII", desc);
0165 }
0166
0167
0168 DEFINE_FWK_MODULE(SiPixelFEDChannelContainerWriteFromASCII);