Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:25

0001 // -*- C++ -*-
0002 //
0003 // Package:    CondFormats/SiPixelObjects
0004 // Class:      SiPixelQualityProbabilitiesTestWriter
0005 //
0006 /**\class SiPixelQualityProbabilitiesTestWriter SiPixelQualityProbabilitiesTestWriter.cc CondFormats/SiPixelObjects/plugins/SiPixelQualityProbabilitiesTestWriter.cc
0007  Description: class to build the SiPixel Quality probabilities 
0008 */
0009 //
0010 // Original Author:  Marco Musich
0011 //         Created:  Wed, 30 Nov 2018 13:22:00 GMT
0012 //
0013 //
0014 
0015 // system include files
0016 #include <memory>
0017 #include <fstream>
0018 
0019 // user include files
0020 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0021 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
0022 #include "CondFormats/SiPixelObjects/interface/SiPixelQualityProbabilities.h"
0023 #include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
0024 #include "FWCore/Framework/interface/ESHandle.h"
0025 #include "FWCore/Framework/interface/ESWatcher.h"
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 
0033 //
0034 // class declaration
0035 //
0036 
0037 class SiPixelQualityProbabilitiesWriteFromASCII : public edm::one::EDAnalyzer<> {
0038 public:
0039   explicit SiPixelQualityProbabilitiesWriteFromASCII(const edm::ParameterSet&);
0040   ~SiPixelQualityProbabilitiesWriteFromASCII() override;
0041 
0042   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0043 
0044 private:
0045   void analyze(const edm::Event&, const edm::EventSetup&) override;
0046   void endJob() override;
0047 
0048   // ----------member data ---------------------------
0049   const std::string m_ProbInputs;
0050   const std::string m_record;
0051   const bool printdebug_;
0052   std::unique_ptr<SiPixelQualityProbabilities> myProbabilities;
0053 };
0054 
0055 //
0056 // constructors and destructor
0057 //
0058 SiPixelQualityProbabilitiesWriteFromASCII::SiPixelQualityProbabilitiesWriteFromASCII(const edm::ParameterSet& iConfig)
0059     : m_ProbInputs(iConfig.getParameter<std::string>("probabilities")),
0060       m_record(iConfig.getParameter<std::string>("record")),
0061       printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", false)) {
0062   //now do what ever initialization is needed
0063   myProbabilities = std::make_unique<SiPixelQualityProbabilities>();
0064 }
0065 
0066 SiPixelQualityProbabilitiesWriteFromASCII::~SiPixelQualityProbabilitiesWriteFromASCII() = default;
0067 
0068 //
0069 // member functions
0070 //
0071 
0072 // ------------ method called for each event  ------------
0073 void SiPixelQualityProbabilitiesWriteFromASCII::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0074   using namespace edm;
0075   std::ifstream myfile(m_ProbInputs);
0076   std::string line;
0077 
0078   SiPixelQualityProbabilities::probabilityVec myProbVector;
0079 
0080   if (myfile.is_open()) {
0081     while (getline(myfile, line)) {
0082       if (printdebug_) {
0083         edm::LogInfo("SiPixelQualityProbabilitiesWriteFromASCII") << line << std::endl;
0084       }
0085       std::istringstream iss(line);
0086       int pileupBinId, nEntries;
0087       iss >> pileupBinId >> nEntries;
0088       edm::LogInfo("SiPixelQualityProbabilitiesWriteFromASCII")
0089           << "PILEUP BIN/ENTRIES:  " << pileupBinId << " " << nEntries << std::endl;
0090       std::vector<std::string> ids(nEntries, "");
0091       std::vector<float> probs(nEntries, 0.0);
0092       for (int i = 0; i < nEntries; ++i) {
0093         iss >> ids.at(i) >> probs.at(i);
0094         if (printdebug_) {
0095           edm::LogInfo("SiPixelQualityProbabilitiesWriteFromASCII") << ids.at(i) << " " << probs.at(i) << std::endl;
0096         }
0097         auto idAndProb = std::make_pair(ids.at(i), probs.at(i));
0098         myProbVector.push_back(idAndProb);
0099       }
0100       if (nEntries > 0) {
0101         myProbabilities->setProbabilities(pileupBinId, myProbVector);
0102       }
0103       myProbVector.clear();
0104     }
0105     myfile.close();
0106   }
0107 
0108   if (printdebug_) {
0109     edm::LogInfo("SiPixelQualityProbabilitiesWriteFromASCII") << "Content of SiPixelQualityProbabilities " << std::endl;
0110     // use buil-in method in the CondFormat
0111     myProbabilities->printAll();
0112   }
0113 }
0114 
0115 // ------------ method called once each job just after ending the event loop  ------------
0116 void SiPixelQualityProbabilitiesWriteFromASCII::endJob() {
0117   edm::LogInfo("SiPixelQualityProbabilitiesWriteFromASCII")
0118       << "Size of SiPixelQualityProbabilities object " << myProbabilities->size() << std::endl
0119       << std::endl;
0120 
0121   // Form the data here
0122   edm::Service<cond::service::PoolDBOutputService> poolDbService;
0123   if (poolDbService.isAvailable()) {
0124     cond::Time_t valid_time = poolDbService->currentTime();
0125     // this writes the payload to begin in current run defined in cfg
0126     poolDbService->writeOneIOV(*myProbabilities, valid_time, m_record);
0127   }
0128 }
0129 
0130 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0131 void SiPixelQualityProbabilitiesWriteFromASCII::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0132   edm::ParameterSetDescription desc;
0133   desc.setComment("Writes payloads of type SiPixelQualityProbabilities");
0134   desc.addUntracked<bool>("printDebug", true);
0135   desc.add<std::string>("record", "SiPixelStatusScenarioProbabilityRcd");
0136   desc.add<std::string>("probabilities", "");
0137   descriptions.add("SiPixelQualityProbabilitiesWriteFromASCII", desc);
0138 }
0139 
0140 //define this as a plug-in
0141 DEFINE_FWK_MODULE(SiPixelQualityProbabilitiesWriteFromASCII);