Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTPrescaleRecorder
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Martin Grunewald
0007  *
0008  */
0009 
0010 #include "HLTrigger/HLTcore/interface/HLTPrescaleRecorder.h"
0011 #include "DataFormats/Provenance/interface/ProcessHistory.h"
0012 #include "DataFormats/Provenance/interface/Timestamp.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 
0016 #include <sys/time.h>
0017 #include <string>
0018 #include <ostream>
0019 
0020 using namespace std;
0021 using namespace edm;
0022 using namespace trigger;
0023 
0024 //
0025 // constructors and destructor
0026 //
0027 HLTPrescaleRecorder::HLTPrescaleRecorder(const edm::ParameterSet& ps)
0028     : src_(ps.getParameter<int>("src")),
0029       run_(ps.getParameter<bool>("run")),
0030       lumi_(ps.getParameter<bool>("lumi")),
0031       event_(ps.getParameter<bool>("event")),
0032       condDB_(ps.getParameter<bool>("condDB")),
0033       psetName_(ps.getParameter<string>("psetName")),
0034       hltInputTag_(ps.getParameter<InputTag>("hltInputTag")),
0035       hltInputToken_(),
0036       hltDBTag_(ps.getParameter<string>("hltDBTag")),
0037       hltPrescaleTableCondToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", hltDBTag_))),
0038       ps_(nullptr),
0039       db_(nullptr),
0040       hltHandle_(),
0041       hlt_() {
0042   if (src_ == 1) {
0043     // Run
0044     hltInputToken_ = consumes<trigger::HLTPrescaleTable, edm::InRun>(hltInputTag_);
0045   } else if (src_ == 2) {
0046     // Lumi
0047     hltInputToken_ = consumes<trigger::HLTPrescaleTable, edm::InLumi>(hltInputTag_);
0048   } else if (src_ == 3) {
0049     // Event
0050     hltInputToken_ = consumes<trigger::HLTPrescaleTable>(hltInputTag_);
0051   }
0052 
0053   LogInfo("HLTPrescaleRecorder") << "src:run-lumi-event-condDB+psetName+tags: " << src_ << ":" << run_ << "-" << lumi_
0054                                  << "-" << event_ << "-" << condDB_ << "+" << psetName_ << "+" << hltInputTag_.encode()
0055                                  << "+" << hltDBTag_;
0056 
0057   if (edm::Service<edm::service::PrescaleService>().isAvailable()) {
0058     ps_ = edm::Service<edm::service::PrescaleService>().operator->();
0059   } else if (src_ == 0) {
0060     LogError("HLTPrescaleRecorder") << "PrescaleService requested as source but unavailable!";
0061   }
0062 
0063   if (edm::Service<cond::service::PoolDBOutputService>().isAvailable()) {
0064     db_ = edm::Service<cond::service::PoolDBOutputService>().operator->();
0065   } else if (condDB_) {
0066     LogError("HLTPrescaleRecorder") << "PoolDBOutputService requested as destination but unavailable!";
0067   }
0068 
0069   if (run_)
0070     produces<HLTPrescaleTable, edm::Transition::EndRun>("Run");
0071   if (lumi_)
0072     produces<HLTPrescaleTable, edm::Transition::EndLuminosityBlock>("Lumi");
0073   if (event_)
0074     produces<HLTPrescaleTable, edm::Transition::Event>("Event");
0075 }
0076 
0077 HLTPrescaleRecorder::~HLTPrescaleRecorder() = default;
0078 
0079 //
0080 // member functions
0081 //
0082 
0083 void HLTPrescaleRecorder::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0084   edm::ParameterSetDescription desc;
0085   //  # (single) source:
0086   //  # -1:PrescaleServicePSet, 0:PrescaleService,
0087   //  #  1:Run, 2:Lumi, 3:Event, 4:CondDB
0088   desc.add<int>("src", 0);
0089   //  # (multiple) destinations
0090   desc.add<bool>("run", true);
0091   desc.add<bool>("lumi", true);
0092   desc.add<bool>("event", true);
0093   desc.add<bool>("condDB", true);
0094   //  # src=-1
0095   desc.add<std::string>("psetName", "");
0096   //  # src= 1,2,3
0097   desc.add<edm::InputTag>("hltInputTag", edm::InputTag("", "", ""));
0098   //  # src= 4
0099   desc.add<std::string>("hltDBTag", "");
0100   //
0101   descriptions.add("hltPrescaleRecorder", desc);
0102 }
0103 
0104 void HLTPrescaleRecorder::beginRun(edm::Run const& iRun, const edm::EventSetup& iSetup) {
0105   hlt_ = HLTPrescaleTable();
0106 
0107   if (src_ == -1) {
0108     /// From PrescaleTable tracked PSet
0109     ParameterSet const& pPSet(getProcessParameterSetContainingModule(moduleDescription()));
0110     ParameterSet iPS(pPSet.getParameter<ParameterSet>(psetName_));
0111 
0112     string defaultLabel(iPS.getParameter<std::string>("lvl1DefaultLabel"));
0113     vector<string> labels(iPS.getParameter<std::vector<std::string> >("lvl1Labels"));
0114     vector<ParameterSet> vpTable(iPS.getParameter<std::vector<ParameterSet> >("prescaleTable"));
0115 
0116     unsigned int set(0);
0117     const unsigned int n(labels.size());
0118     for (unsigned int i = 0; i != n; ++i) {
0119       if (labels[i] == defaultLabel)
0120         set = i;
0121     }
0122 
0123     map<string, vector<unsigned int> > table;
0124     const unsigned int m(vpTable.size());
0125     for (unsigned int i = 0; i != m; ++i) {
0126       table[vpTable[i].getParameter<std::string>("pathName")] =
0127           vpTable[i].getParameter<std::vector<unsigned int> >("prescales");
0128     }
0129     hlt_ = HLTPrescaleTable(set, labels, table);
0130 
0131   } else if (src_ == 0) {
0132     /// From PrescaleService
0133     /// default index updated at lumi block boundaries
0134     if (ps_ != nullptr) {
0135       hlt_ = HLTPrescaleTable(ps_->getLvl1IndexDefault(), ps_->getLvl1Labels(), ps_->getPrescaleTable());
0136     } else {
0137       hlt_ = HLTPrescaleTable();
0138       LogError("HLTPrescaleRecorder") << "PrescaleService not found!";
0139     }
0140   } else if (src_ == 1) {
0141     /// From Run Block
0142     if (iRun.getByToken(hltInputToken_, hltHandle_)) {
0143       hlt_ = *hltHandle_;
0144     } else {
0145       LogError("HLTPrescaleRecorder") << "HLTPrescaleTable not found in Run!";
0146     }
0147   } else if (src_ == 4) {
0148     /// From CondDB (needs ESProducer module as well)
0149     auto const& hltPrescaleTableCond = iSetup.getData(hltPrescaleTableCondToken_);
0150     hlt_ = hltPrescaleTableCond.hltPrescaleTable();
0151   }
0152 
0153   return;
0154 }
0155 
0156 void HLTPrescaleRecorder::beginLuminosityBlock(edm::LuminosityBlock const& iLumi, const edm::EventSetup& iSetup) {
0157   if (src_ == 0) {
0158     /// From PrescaleService
0159     /// default index updated at lumi block boundaries
0160     if (ps_ != nullptr) {
0161       hlt_ = HLTPrescaleTable(ps_->getLvl1IndexDefault(), ps_->getLvl1Labels(), ps_->getPrescaleTable());
0162     } else {
0163       hlt_ = HLTPrescaleTable();
0164       LogError("HLTPrescaleRecorder") << "PrescaleService not found!";
0165     }
0166   } else if (src_ == 2) {
0167     /// From Lumi Block
0168     if (iLumi.getByToken(hltInputToken_, hltHandle_)) {
0169       hlt_ = *hltHandle_;
0170     } else {
0171       hlt_ = HLTPrescaleTable();
0172       LogError("HLTPrescaleRecorder") << "HLTPrescaleTable not found in LumiBlock!";
0173     }
0174   }
0175 
0176   return;
0177 }
0178 
0179 void HLTPrescaleRecorder::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0180   if (src_ == 3) {
0181     /// From Event Block
0182     if (iEvent.getByToken(hltInputToken_, hltHandle_)) {
0183       hlt_ = *hltHandle_;
0184     } else {
0185       hlt_ = HLTPrescaleTable();
0186       LogError("HLTPrescaleRecorder") << "HLTPrescaleTable not found in Event!";
0187     }
0188   }
0189 
0190   if (event_) {
0191     /// Writing to Event
0192     unique_ptr<HLTPrescaleTable> product(new HLTPrescaleTable(hlt_));
0193     iEvent.put(std::move(product), "Event");
0194   }
0195 
0196   return;
0197 }
0198 void HLTPrescaleRecorder::endLuminosityBlock(edm::LuminosityBlock const& iLumi, const edm::EventSetup& iSetup) {}
0199 
0200 void HLTPrescaleRecorder::endLuminosityBlockProduce(edm::LuminosityBlock& iLumi, const edm::EventSetup& iSetup) {
0201   if (lumi_) {
0202     /// Writing to Lumi Block
0203     unique_ptr<HLTPrescaleTable> product(new HLTPrescaleTable(hlt_));
0204     iLumi.put(std::move(product), "Lumi");
0205   }
0206   return;
0207 }
0208 
0209 void HLTPrescaleRecorder::endRun(edm::Run const& iRun, const edm::EventSetup& iSetup) {
0210   /// Dump to logfile
0211   ostringstream oss;
0212   const unsigned int n(hlt_.size());
0213   oss << "PrescaleTable: # of labels = " << n << endl;
0214   const vector<string>& labels(hlt_.labels());
0215   for (unsigned int i = 0; i != n; ++i) {
0216     oss << " " << i << "/'" << labels.at(i) << "'";
0217   }
0218   oss << endl;
0219   const map<string, vector<unsigned int> >& table(hlt_.table());
0220   oss << "PrescaleTable: # of paths = " << table.size() << endl;
0221   const map<string, vector<unsigned int> >::const_iterator tb(table.begin());
0222   const map<string, vector<unsigned int> >::const_iterator te(table.end());
0223   for (map<string, vector<unsigned int> >::const_iterator ti = tb; ti != te; ++ti) {
0224     for (unsigned int i = 0; i != n; ++i) {
0225       oss << " " << ti->second.at(i);
0226     }
0227     oss << " " << ti->first << endl;
0228   }
0229   LogVerbatim("HLTPrescaleRecorder") << oss.str();
0230 
0231   if (condDB_) {
0232     /// Writing to CondDB (needs PoolDBOutputService)
0233     if (db_ != nullptr) {
0234       HLTPrescaleTableCond product(hlt_);
0235       ::timeval tv;
0236       gettimeofday(&tv, nullptr);
0237       edm::Timestamp tstamp((unsigned long long)tv.tv_sec);
0238       db_->writeOneIOV(product, tstamp.value(), "HLTPrescaleTableRcd");
0239     } else {
0240       LogError("HLTPrescaleRecorder") << "PoolDBOutputService not available!";
0241     }
0242   }
0243 
0244   return;
0245 }
0246 
0247 void HLTPrescaleRecorder::endRunProduce(edm::Run& iRun, const edm::EventSetup& iSetup) {
0248   if (run_) {
0249     /// Writing to Run Block
0250     unique_ptr<HLTPrescaleTable> product(new HLTPrescaleTable(hlt_));
0251     iRun.put(std::move(product), "Run");
0252   }
0253 }