Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // Author: Benedikt Hegner, Tom Cornelis
0002 // Email:  benedikt.hegner@cern.ch, tom.cornelis@cern.ch
0003 
0004 #include <vector>
0005 #include <memory>
0006 #include <string>
0007 #include <fstream>
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ServiceRegistry/interface/Service.h"
0014 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0015 #include "CondFormats/JetMETObjects/interface/QGLikelihoodObject.h"
0016 #include "FWCore/ParameterSet/interface/FileInPath.h"
0017 
0018 class QGLikelihoodSystematicsDBWriter : public edm::one::EDAnalyzer<> {
0019 public:
0020   QGLikelihoodSystematicsDBWriter(const edm::ParameterSet&);
0021   void beginJob() override;
0022   void analyze(const edm::Event&, const edm::EventSetup&) override {}
0023   void endJob() override {}
0024   ~QGLikelihoodSystematicsDBWriter() override {}
0025 
0026 private:
0027   std::string fileName;
0028   std::string payloadTag;
0029 };
0030 
0031 // Constructor
0032 QGLikelihoodSystematicsDBWriter::QGLikelihoodSystematicsDBWriter(const edm::ParameterSet& pSet) {
0033   fileName = pSet.getParameter<std::string>("src");
0034   payloadTag = pSet.getParameter<std::string>("payload");
0035 }
0036 
0037 // Begin Job
0038 void QGLikelihoodSystematicsDBWriter::beginJob() {
0039   QGLikelihoodSystematicsObject payload;
0040   payload.data.clear();
0041 
0042   std::ifstream database;
0043   database.open(edm::FileInPath(fileName.c_str()).fullPath().c_str(), std::ios::in);
0044   if (!database.is_open()) {
0045     edm::LogError("FileNotFound") << "Could not open file " << fileName << std::endl;
0046     return;
0047   }
0048   std::string line;
0049   while (std::getline(database, line)) {
0050     float ptMin, ptMax, etaMin, etaMax, rhoMin, rhoMax, a_q, b_q, a_g, b_g, lmin, lmax;
0051     char tag[1023], leadchar;
0052     sscanf(line.c_str(), "%c", &leadchar);
0053     if ((leadchar == '#') || (leadchar == '!'))
0054       continue;  //Skip those lines
0055     sscanf(line.c_str(),
0056            "%s %f %f %f %f %f %f %f %f %f %f %f %f",
0057            &tag[0],
0058            &ptMin,
0059            &ptMax,
0060            &rhoMin,
0061            &rhoMax,
0062            &etaMin,
0063            &etaMax,
0064            &a_q,
0065            &b_q,
0066            &a_g,
0067            &b_g,
0068            &lmin,
0069            &lmax);
0070 
0071     QGLikelihoodCategory category;
0072     category.RhoMin = rhoMin;
0073     category.RhoMax = rhoMax;
0074     category.PtMin = ptMin;
0075     category.PtMax = ptMax;
0076     category.EtaMin = etaMin;
0077     category.EtaMax = etaMax;
0078     category.QGIndex = 0;
0079     category.VarIndex = -1;
0080 
0081     //quark entry
0082     QGLikelihoodSystematicsObject::Entry quarkEntry;
0083     quarkEntry.systCategory = category;
0084     quarkEntry.a = a_q;
0085     quarkEntry.b = b_q;
0086     quarkEntry.lmin = lmin;
0087     quarkEntry.lmax = lmax;
0088 
0089     //gluon entry
0090     QGLikelihoodSystematicsObject::Entry gluonEntry = quarkEntry;
0091     gluonEntry.systCategory.QGIndex = 1;
0092     gluonEntry.a = a_g;
0093     gluonEntry.b = b_g;
0094 
0095     payload.data.push_back(quarkEntry);
0096     payload.data.push_back(gluonEntry);
0097   }
0098   database.close();
0099 
0100   // Now write it into the DB
0101   edm::LogInfo("UserOutput") << "Opening PoolDBOutputService" << std::endl;
0102   edm::Service<cond::service::PoolDBOutputService> s;
0103   if (s.isAvailable()) {
0104     edm::LogInfo("UserOutput") << "Setting up payload with " << payload.data.size() << " entries and tag " << payloadTag
0105                                << std::endl;
0106     s->writeOneIOV(payload, s->beginOfTime(), payloadTag);
0107   }
0108   edm::LogInfo("UserOutput") << "Wrote in CondDB QGLikelihoodSystematic payload label: " << payloadTag << std::endl;
0109 }
0110 
0111 DEFINE_FWK_MODULE(QGLikelihoodSystematicsDBWriter);