File indexing completed on 2023-03-17 10:58:54
0001 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0002
0003 #include "FWCore/Version/interface/GetReleaseVersion.h"
0004
0005 class DQMHarvestingMetadata : public DQMEDHarvester {
0006 public:
0007
0008 DQMHarvestingMetadata(const edm::ParameterSet& ps);
0009
0010
0011 ~DQMHarvestingMetadata() override;
0012
0013 protected:
0014
0015 void dqmEndRun(DQMStore::IBooker& ibooker,
0016 DQMStore::IGetter& igetter,
0017 edm::Run const& iRun,
0018 edm::EventSetup const& ) override;
0019 void dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0020 DQMStore::IGetter& igetter,
0021 edm::LuminosityBlock const& iLumi,
0022 edm::EventSetup const& ) override;
0023
0024 void dqmEndJob(DQMStore::IBooker&, DQMStore::IGetter&) override{};
0025
0026 private:
0027 std::string eventInfoFolder_;
0028 std::string subsystemname_;
0029
0030 MonitorElement* runId_;
0031 MonitorElement* runList_;
0032 MonitorElement* runStartTimeStamp_;
0033 MonitorElement* lumisecId_;
0034 MonitorElement* firstLumisecId_;
0035 MonitorElement* lastLumisecId_;
0036
0037 MonitorElement* processTimeStamp_;
0038 MonitorElement* hostName_;
0039 MonitorElement* workingDir_;
0040 MonitorElement* cmsswVer_;
0041
0042 double currentTime_;
0043 };
0044
0045 static inline double stampToReal(edm::Timestamp time) {
0046 return (time.value() >> 32) + 1e-6 * (time.value() & 0xffffffff);
0047 }
0048
0049 static inline double stampToReal(const timeval& time) { return time.tv_sec + 1e-6 * time.tv_usec; }
0050
0051 DQMHarvestingMetadata::DQMHarvestingMetadata(const edm::ParameterSet& ps) {
0052 struct timeval now;
0053 gettimeofday(&now, nullptr);
0054 currentTime_ = stampToReal(now);
0055
0056
0057 std::string folder = ps.getUntrackedParameter<std::string>("eventInfoFolder", "EventInfo");
0058 subsystemname_ = ps.getUntrackedParameter<std::string>("subSystemFolder", "YourSubsystem");
0059
0060 eventInfoFolder_ = subsystemname_ + "/" + folder;
0061 }
0062
0063 DQMHarvestingMetadata::~DQMHarvestingMetadata() = default;
0064
0065 void DQMHarvestingMetadata::dqmEndRun(DQMStore::IBooker& ibooker,
0066 DQMStore::IGetter& igetter,
0067 edm::Run const& iRun,
0068 edm::EventSetup const& ) {
0069 ibooker.setCurrentFolder(eventInfoFolder_);
0070
0071 runList_ = ibooker.bookString("Run", "");
0072 runId_ = ibooker.bookInt("iRun");
0073 runStartTimeStamp_ = ibooker.bookFloat("runStartTimeStamp");
0074
0075 if (runList_->getStringValue().empty()) {
0076 std::string run = std::to_string(iRun.id().run());
0077 runList_->Fill(run);
0078 runId_->Fill(iRun.id().run());
0079
0080 runStartTimeStamp_->Fill(stampToReal(iRun.beginTime()));
0081 } else {
0082 std::string run = runList_->getStringValue() + "," + std::to_string(iRun.id().run());
0083 runList_->Fill(run);
0084
0085 runId_->Fill(999999);
0086 }
0087
0088 processTimeStamp_ = ibooker.bookFloat("processTimeStamp");
0089 processTimeStamp_->Fill(currentTime_);
0090
0091 char hostname[65];
0092 gethostname(hostname, 64);
0093 hostname[64] = 0;
0094 hostName_ = ibooker.bookString("hostName", hostname);
0095 char* pwd = getcwd(nullptr, 0);
0096 workingDir_ = ibooker.bookString("workingDir", pwd);
0097 free(pwd);
0098 cmsswVer_ = ibooker.bookString("CMSSW_Version", edm::getReleaseVersion());
0099 }
0100
0101 void DQMHarvestingMetadata::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0102 DQMStore::IGetter& igetter,
0103 edm::LuminosityBlock const& iLumi,
0104 edm::EventSetup const& ) {
0105 int lumi = iLumi.luminosityBlock();
0106
0107 ibooker.setCurrentFolder(eventInfoFolder_);
0108 firstLumisecId_ = ibooker.bookInt("firstLumiSection");
0109 lastLumisecId_ = ibooker.bookInt("lastLumiSection");
0110 lumisecId_ = ibooker.bookInt("iLumiSection");
0111 lumisecId_->Fill(lumi);
0112
0113 if (firstLumisecId_->getIntValue() == 0 || firstLumisecId_->getIntValue() > lumi) {
0114 firstLumisecId_->Fill(lumi);
0115 }
0116 if (lastLumisecId_->getIntValue() == 0 || lastLumisecId_->getIntValue() > lumi) {
0117 lastLumisecId_->Fill(lumi);
0118 }
0119 }
0120
0121 #include "FWCore/Framework/interface/MakerMacros.h"
0122 DEFINE_FWK_MODULE(DQMHarvestingMetadata);