SiPixelQualityProbabilitiesTestWriter

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151
// -*- C++ -*-
//
// Package:    CondFormats/SiPixelObjects
// Class:      SiPixelQualityProbabilitiesTestWriter
//
/**\class SiPixelQualityProbabilitiesTestWriter SiPixelQualityProbabilitiesTestWriter.cc CondFormats/SiPixelObjects/plugins/SiPixelQualityProbabilitiesTestWriter.cc
 Description: class to build the SiPixel Quality probabilities 
*/
//
// Original Author:  Marco Musich
//         Created:  Wed, 30 Nov 2018 13:22:00 GMT
//
//

// system include files
#include <memory>
#include <fstream>

// user include files
#include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
#include "CondFormats/SiPixelObjects/interface/SiPixelQualityProbabilities.h"
#include "CondFormats/DataRecord/interface/SiPixelQualityFromDbRcd.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/ESWatcher.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

//
// class declaration
//

class SiPixelQualityProbabilitiesTestWriter : public edm::one::EDAnalyzer<> {
public:
  explicit SiPixelQualityProbabilitiesTestWriter(const edm::ParameterSet&);
  ~SiPixelQualityProbabilitiesTestWriter() override;

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
  void analyze(const edm::Event&, const edm::EventSetup&) override;
  void endJob() override;

  // ----------member data ---------------------------
  const std::string m_ProbInputs;
  const std::string m_SnapshotInputs;
  const std::string m_record;
  const bool printdebug_;
  std::unique_ptr<SiPixelQualityProbabilities> myProbabilities;
};

//
// constructors and destructor
//
SiPixelQualityProbabilitiesTestWriter::SiPixelQualityProbabilitiesTestWriter(const edm::ParameterSet& iConfig)
    : m_ProbInputs(iConfig.getParameter<std::string>("probabilities")),
      m_SnapshotInputs(iConfig.getParameter<std::string>("snapshots")),
      m_record(iConfig.getParameter<std::string>("record")),
      printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", false)) {
  //now do what ever initialization is needed
  myProbabilities = std::make_unique<SiPixelQualityProbabilities>();
}

SiPixelQualityProbabilitiesTestWriter::~SiPixelQualityProbabilitiesTestWriter() = default;

//
// member functions
//

// ------------ method called for each event  ------------
void SiPixelQualityProbabilitiesTestWriter::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
  using namespace edm;
  std::ifstream myfile(m_ProbInputs);
  std::ifstream mysnapshots(m_SnapshotInputs);
  std::string line1, line2;
  std::map<int, std::string> snapshotIdToString;

  if (mysnapshots.is_open()) {
    while (getline(mysnapshots, line1)) {
      //edm::LogInfo("SiPixelQualityProbabilitiesTestWriter") << line1 << std::endl;
      std::istringstream iss(line1);
      int id, run, ls;
      iss >> id >> run >> ls;
      snapshotIdToString[id] = std::to_string(run) + "_" + std::to_string(ls);
    }
  }

  SiPixelQualityProbabilities::probabilityVec myProbVector;

  if (myfile.is_open()) {
    while (getline(myfile, line2)) {
      edm::LogInfo("SiPixelQualityProbabilitiesTestWriter") << line2 << std::endl;
      std::istringstream iss(line2);
      int pileupBinId, nEntries;
      iss >> pileupBinId >> nEntries;
      edm::LogInfo("SiPixelQualityProbabilitiesTestWriter")
          << "PILEUP BIN/ENTRIES:  " << pileupBinId << " " << nEntries << std::endl;
      std::vector<int> ids(nEntries, 0);
      std::vector<float> probs(nEntries, 0.0);
      for (int i = 0; i < nEntries; ++i) {
        iss >> ids.at(i) >> probs.at(i);
        //edm::LogInfo("SiPixelQualityProbabilitiesTestWriter") << ids.at(i) << " " << probs.at(i)<< std::endl;
        auto idAndProb = std::make_pair(snapshotIdToString.at(ids.at(i)), probs.at(i));
        myProbVector.push_back(idAndProb);
      }
      if (nEntries > 0)
        myProbabilities->setProbabilities(pileupBinId, myProbVector);
      myProbVector.clear();
    }
    myfile.close();
  }

  if (printdebug_) {
    edm::LogInfo("SiPixelQualityProbabilitiesTestWriter") << "Content of SiPixelQualityProbabilities " << std::endl;
    // use buil-in method in the CondFormat
    myProbabilities->printAll();
  }
}

// ------------ method called once each job just after ending the event loop  ------------
void SiPixelQualityProbabilitiesTestWriter::endJob() {
  edm::LogInfo("SiPixelQualityProbabilitiesTestWriter")
      << "Size of SiPixelQualityProbabilities object " << myProbabilities->size() << std::endl
      << std::endl;

  // Form the data here
  edm::Service<cond::service::PoolDBOutputService> poolDbService;
  if (poolDbService.isAvailable()) {
    cond::Time_t valid_time = poolDbService->currentTime();
    // this writes the payload to begin in current run defined in cfg
    poolDbService->writeOneIOV(*myProbabilities, valid_time, m_record);
  }
}

// ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
void SiPixelQualityProbabilitiesTestWriter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.setComment("Writes payloads of type SiPixelQualityProbabilities");
  desc.addUntracked<bool>("printDebug", true);
  desc.add<std::string>("record", "SiPixelStatusScenarioProbabilityRcd");
  desc.add<std::string>("snapshots", "");
  desc.add<std::string>("probabilities", "");
  descriptions.add("SiPixelQualityProbabilitiesTestWriter", desc);
}

//define this as a plug-in
DEFINE_FWK_MODULE(SiPixelQualityProbabilitiesTestWriter);