Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:09:37

0001 /** \class HLTriggerJSONMonitoring
0002  *
0003  *  Description: This class outputs JSON files with HLT monitoring information.
0004  *
0005  */
0006 
0007 #include <atomic>
0008 #include <fstream>
0009 #include <sstream>
0010 #include <vector>
0011 #include <string>
0012 #include <memory>
0013 #include <algorithm>
0014 
0015 #include <fmt/printf.h>
0016 
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "DataFormats/Common/interface/TriggerResults.h"
0019 #include "EventFilter/Utilities/interface/EvFDaqDirector.h"
0020 #include "EventFilter/Utilities/interface/FastMonitor.h"
0021 #include "EventFilter/Utilities/interface/FastMonitoringService.h"
0022 #include "EventFilter/Utilities/interface/JSONSerializer.h"
0023 #include "EventFilter/Utilities/interface/JsonMonitorable.h"
0024 #include "FWCore/Framework/interface/ESHandle.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/LuminosityBlock.h"
0027 #include "FWCore/Framework/interface/Run.h"
0028 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0030 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "FWCore/ServiceRegistry/interface/Service.h"
0033 #include "FWCore/Utilities/interface/Adler32Calculator.h"
0034 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0035 
0036 //note this was updated 20/10/22 to change the logic such that instead
0037 //of having it passing the L1Seed, it is now number passing pre prescale
0038 //this is indentical for "standard" paths and more meaningful for
0039 //the special paths which are affected
0040 //a standard path logic goes "trigger type -> l1 seed -> prescale -> other selection"
0041 
0042 struct HLTriggerJSONMonitoringData {
0043   // variables accumulated event by event in each stream
0044   struct stream {
0045     unsigned int processed;               // number of events processed
0046     std::vector<unsigned int> hltWasRun;  // number of events where each path was run
0047     std::vector<unsigned int> hltPrePS;   // number of events where each path made it to the prescale module
0048     std::vector<unsigned int> hltPostPS;  // number of events where each path passed the prescale
0049     std::vector<unsigned int> hltAccept;  // number of events accepted by each path
0050     std::vector<unsigned int> hltReject;  // number of events rejected by each path
0051     std::vector<unsigned int> hltErrors;  // number of events with errors in each path
0052     std::vector<unsigned int> datasets;   // number of events accepted by each dataset
0053   };
0054 
0055   // variables initialised for each run
0056   struct run {
0057     std::string streamDestination;
0058     std::string streamMergeType;
0059     std::string baseRunDir;   // base directory from EvFDaqDirector
0060     std::string jsdFileName;  // definition file name for JSON with rates
0061 
0062     HLTConfigProvider hltConfig;  // HLT configuration for the current run
0063     std::vector<int> posPre;      // position of last HLT prescale filter in each path, or -1 if not present
0064     std::vector<std::vector<unsigned int>> datasets;  // list of paths in each dataset
0065     std::vector<unsigned int> indicesOfTriggerPaths;  // indices of triggers (without DatasetPaths) in TriggerNames
0066   };
0067 
0068   // variables accumulated over the whole lumisection
0069   struct lumisection {
0070     jsoncollector::HistoJ<unsigned int> processed;  // number of events processed
0071     jsoncollector::HistoJ<unsigned int> hltWasRun;  // number of events where each path was run
0072     jsoncollector::HistoJ<unsigned int> hltPrePS;   // number of events where each path made it to the prescale module
0073     jsoncollector::HistoJ<unsigned int> hltPostPS;  // number of events where each path passed the prescale
0074     jsoncollector::HistoJ<unsigned int> hltAccept;  // number of events accepted by each path
0075     jsoncollector::HistoJ<unsigned int> hltReject;  // number of events rejected by each path
0076     jsoncollector::HistoJ<unsigned int> hltErrors;  // number of events with errors in each path
0077     jsoncollector::HistoJ<unsigned int> datasets;   // number of events accepted by each dataset
0078   };
0079 };
0080 
0081 class HLTriggerJSONMonitoring : public edm::global::EDAnalyzer<
0082                                     // per-stream information
0083                                     edm::StreamCache<HLTriggerJSONMonitoringData::stream>,
0084                                     // per-run accounting
0085                                     edm::RunCache<HLTriggerJSONMonitoringData::run>,
0086                                     // accumulate per-lumisection statistics
0087                                     edm::LuminosityBlockSummaryCache<HLTriggerJSONMonitoringData::lumisection>> {
0088 public:
0089   // constructor
0090   explicit HLTriggerJSONMonitoring(const edm::ParameterSet&);
0091 
0092   // destructor
0093   ~HLTriggerJSONMonitoring() override = default;
0094 
0095   // validate the configuration and optionally fill the default values
0096   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0097 
0098   // called for each Event
0099   void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;
0100 
0101   // -- inherited from edm::StreamCache<HLTriggerJSONMonitoringData::stream>
0102 
0103   // called once for each Stream being used in the job to create the cache object that will be used for that particular Stream
0104   std::unique_ptr<HLTriggerJSONMonitoringData::stream> beginStream(edm::StreamID) const override;
0105 
0106   // called when the Stream is switching from one LuminosityBlock to a new LuminosityBlock.
0107   void streamBeginLuminosityBlock(edm::StreamID, edm::LuminosityBlock const&, edm::EventSetup const&) const override;
0108 
0109   // -- inherited from edm::RunCache<HLTriggerJSONMonitoringData::run>
0110 
0111   // called each time the Source sees a new Run, and guaranteed to finish before any Stream calls streamBeginRun for that same Run
0112   std::shared_ptr<HLTriggerJSONMonitoringData::run> globalBeginRun(edm::Run const&,
0113                                                                    edm::EventSetup const&) const override;
0114 
0115   // called after all Streams have finished processing a given Run (i.e. streamEndRun for all Streams have completed)
0116   void globalEndRun(edm::Run const&, edm::EventSetup const&) const override;
0117 
0118   // -- inherited from edm::LuminosityBlockSummaryCache<HLTriggerJSONMonitoringData::lumisection>
0119 
0120   // called each time the Source sees a new LuminosityBlock
0121   std::shared_ptr<HLTriggerJSONMonitoringData::lumisection> globalBeginLuminosityBlockSummary(
0122       edm::LuminosityBlock const&, edm::EventSetup const&) const override;
0123 
0124   // called when a Stream has finished processing a LuminosityBlock, after streamEndLuminosityBlock
0125   void streamEndLuminosityBlockSummary(edm::StreamID,
0126                                        edm::LuminosityBlock const&,
0127                                        edm::EventSetup const&,
0128                                        HLTriggerJSONMonitoringData::lumisection*) const override;
0129 
0130   // called after the streamEndLuminosityBlockSummary method for all Streams have finished processing a given LuminosityBlock
0131   void globalEndLuminosityBlockSummary(edm::LuminosityBlock const&,
0132                                        edm::EventSetup const&,
0133                                        HLTriggerJSONMonitoringData::lumisection*) const override;
0134 
0135 private:
0136   static constexpr const char* streamName_ = "streamHLTRates";
0137 
0138   static constexpr const char* datasetPathNamePrefix_ = "Dataset_";
0139 
0140   static void writeJsdFile(HLTriggerJSONMonitoringData::run const&);
0141   static void writeIniFile(HLTriggerJSONMonitoringData::run const&, unsigned int);
0142 
0143   // configuration
0144   const edm::InputTag triggerResults_;                               // InputTag for TriggerResults
0145   const edm::EDGetTokenT<edm::TriggerResults> triggerResultsToken_;  // Token for TriggerResults
0146 };
0147 
0148 // constructor
0149 HLTriggerJSONMonitoring::HLTriggerJSONMonitoring(edm::ParameterSet const& config)
0150     : triggerResults_(config.getParameter<edm::InputTag>("triggerResults")),
0151       triggerResultsToken_(consumes(triggerResults_)) {
0152   if (edm::Service<evf::EvFDaqDirector>().isAvailable()) {
0153     //output initemp file. This lets hltd know number of streams early
0154     std::string initFileName = edm::Service<evf::EvFDaqDirector>()->getInitTempFilePath("streamHLTRates");
0155     std::ofstream file(initFileName);
0156     if (!file)
0157       throw cms::Exception("HLTriggerJsonMonitoring")
0158           << "Cannot create INITEMP file: " << initFileName << " error: " << strerror(errno);
0159     file.close();
0160   }
0161 }
0162 
0163 // validate the configuration and optionally fill the default values
0164 void HLTriggerJSONMonitoring::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0165   edm::ParameterSetDescription desc;
0166   desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults", "", "@currentProcess"));
0167   descriptions.add("HLTriggerJSONMonitoring", desc);
0168 }
0169 
0170 // called once for each Stream being used in the job to create the cache object that will be used for that particular Stream
0171 std::unique_ptr<HLTriggerJSONMonitoringData::stream> HLTriggerJSONMonitoring::beginStream(edm::StreamID) const {
0172   return std::make_unique<HLTriggerJSONMonitoringData::stream>();
0173 }
0174 
0175 // called each time the Source sees a new Run, and guaranteed to finish before any Stream calls streamBeginRun for that same Run
0176 std::shared_ptr<HLTriggerJSONMonitoringData::run> HLTriggerJSONMonitoring::globalBeginRun(
0177     edm::Run const& run, edm::EventSetup const& setup) const {
0178   auto rundata = std::make_shared<HLTriggerJSONMonitoringData::run>();
0179 
0180   // set the DAQ parameters
0181   if (edm::Service<evf::EvFDaqDirector>().isAvailable()) {
0182     rundata->streamDestination = edm::Service<evf::EvFDaqDirector>()->getStreamDestinations(streamName_);
0183     rundata->streamMergeType =
0184         edm::Service<evf::EvFDaqDirector>()->getStreamMergeType(streamName_, evf::MergeTypeJSNDATA);
0185     rundata->baseRunDir = edm::Service<evf::EvFDaqDirector>()->baseRunDir();
0186   } else {
0187     rundata->streamDestination = "";
0188     rundata->streamMergeType = "";
0189     rundata->baseRunDir = ".";
0190   }
0191 
0192   // initialize HLTConfigProvider
0193   bool changed = true;
0194   if (not rundata->hltConfig.init(run, setup, triggerResults_.process(), changed)) {
0195     edm::LogError("HLTriggerJSONMonitoring") << "HLTConfigProvider initialization failed!";
0196   } else if (changed) {
0197     // triggerNames from TriggerResults (includes DatasetPaths)
0198     auto const& triggerNames = rundata->hltConfig.triggerNames();
0199     auto const triggerNamesSize = triggerNames.size();
0200 
0201     // update the list of indices of the HLT Paths (without DatasetPaths) in the TriggerNames list
0202     rundata->indicesOfTriggerPaths.clear();
0203     rundata->indicesOfTriggerPaths.reserve(triggerNamesSize);
0204     for (auto triggerNameIdx = 0u; triggerNameIdx < triggerNamesSize; ++triggerNameIdx) {
0205       // skip DatasetPaths
0206       if (triggerNames[triggerNameIdx].find(datasetPathNamePrefix_) != 0) {
0207         rundata->indicesOfTriggerPaths.emplace_back(triggerNameIdx);
0208       }
0209     }
0210     auto const triggersSize = rundata->indicesOfTriggerPaths.size();
0211 
0212     // update the list of paths in each dataset
0213     auto const& datasets = rundata->hltConfig.datasetContents();
0214     auto const& datasetNames = rundata->hltConfig.datasetNames();
0215     auto const datasetsSize = datasetNames.size();
0216     rundata->datasets.resize(datasetsSize);
0217     for (auto ds = 0u; ds < datasetsSize; ++ds) {
0218       auto& dataset = rundata->datasets[ds];
0219       // check if TriggerNames include the DatasetPath corresponding to this Dataset
0220       //  - DatasetPaths are normal cms.Path objects
0221       //  - in Run-3 HLT menus, DatasetPaths are used to define PrimaryDatasets
0222       auto const datasetPathName = datasetPathNamePrefix_ + datasetNames[ds];
0223       auto const datasetPathExists =
0224           std::find(triggerNames.begin(), triggerNames.end(), datasetPathName) != triggerNames.end();
0225       if (datasetPathExists) {
0226         // if a DatasetPath exists, only that Path is assigned to the Dataset
0227         //  - this way, the counts of the Dataset properly include prescales on the DatasetPath
0228         //    and smart-Prescales applied by the DatasetPath to its triggers
0229         dataset.reserve(1);
0230         auto const index = rundata->hltConfig.triggerIndex(datasetPathName);
0231         if (index < triggerNamesSize)
0232           dataset.push_back(index);
0233       } else {
0234         auto const paths = datasets[ds].size();
0235         dataset.reserve(paths);
0236         for (auto p = 0u; p < paths; p++) {
0237           auto const index = rundata->hltConfig.triggerIndex(datasets[ds][p]);
0238           if (index < triggerNamesSize)
0239             dataset.push_back(index);
0240         }
0241       }
0242     }
0243 
0244     // find the positions of the prescale filters
0245     rundata->posPre.resize(triggersSize);
0246     for (auto i = 0u; i < triggersSize; ++i) {
0247       rundata->posPre[i] = -1;
0248       auto const trigNameIndx = rundata->indicesOfTriggerPaths[i];
0249       auto const& moduleLabels = rundata->hltConfig.moduleLabels(trigNameIndx);
0250       for (auto j = 0u; j < moduleLabels.size(); ++j) {
0251         auto const& moduleType = rundata->hltConfig.moduleType(moduleLabels[j]);
0252         if (moduleType == "HLTPrescaler") {
0253           rundata->posPre[i] = j;
0254           break;
0255         }
0256       }
0257     }
0258   }
0259 
0260   // write the per-run .jsd file
0261   rundata->jsdFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.jsd", run.run(), getpid());
0262   writeJsdFile(*rundata);
0263 
0264   // write the per-run .ini file
0265   //iniFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.ini", run.run(), getpid());
0266   writeIniFile(*rundata, run.run());
0267 
0268   return rundata;
0269 }
0270 
0271 // called after all Streams have finished processing a given Run (i.e. streamEndRun for all Streams have completed)
0272 void HLTriggerJSONMonitoring::globalEndRun(edm::Run const&, edm::EventSetup const&) const {}
0273 
0274 // called for each Event
0275 void HLTriggerJSONMonitoring::analyze(edm::StreamID sid, edm::Event const& event, edm::EventSetup const&) const {
0276   auto& stream = *streamCache(sid);
0277   auto const& rundata = *runCache(event.getRun().index());
0278 
0279   ++stream.processed;
0280 
0281   // check that the HLTConfigProvider for the current run has been successfully initialised
0282   if (not rundata.hltConfig.inited())
0283     return;
0284 
0285   // get hold of TriggerResults
0286   edm::Handle<edm::TriggerResults> handle;
0287   if (not event.getByToken(triggerResultsToken_, handle) or not handle.isValid()) {
0288     edm::LogError("HLTriggerJSONMonitoring")
0289         << "TriggerResults with label [" + triggerResults_.encode() + "] not present or invalid";
0290     return;
0291   }
0292   edm::TriggerResults const& results = *handle;
0293   assert(results.size() == rundata.hltConfig.triggerNames().size());
0294 
0295   // check the results for each HLT path
0296   for (auto idx = 0u; idx < rundata.indicesOfTriggerPaths.size(); ++idx) {
0297     auto const triggerPathIdx = rundata.indicesOfTriggerPaths[idx];
0298     auto const& status = results[triggerPathIdx];
0299     if (status.wasrun()) {
0300       ++stream.hltWasRun[idx];
0301       if (status.accept()) {
0302         ++stream.hltPrePS[idx];
0303         ++stream.hltPostPS[idx];
0304         ++stream.hltAccept[idx];
0305       } else {
0306         int const index = (int)status.index();
0307         if (index >= rundata.posPre[idx])
0308           ++stream.hltPrePS[idx];
0309         if (index > rundata.posPre[idx])
0310           ++stream.hltPostPS[idx];
0311         if (status.error())
0312           ++stream.hltErrors[idx];
0313         else
0314           ++stream.hltReject[idx];
0315       }
0316     }
0317   }
0318 
0319   // check the decision for each dataset
0320   // FIXME this ignores the prescales, "smart" prescales, and event selection applied in the OutputModule itself
0321   for (auto i = 0u; i < rundata.datasets.size(); ++i)
0322     if (std::any_of(rundata.datasets[i].begin(), rundata.datasets[i].end(), [&](unsigned int path) {
0323           return results.accept(path);
0324         }))
0325       ++stream.datasets[i];
0326 }
0327 
0328 // called each time the Source sees a new LuminosityBlock
0329 std::shared_ptr<HLTriggerJSONMonitoringData::lumisection> HLTriggerJSONMonitoring::globalBeginLuminosityBlockSummary(
0330     edm::LuminosityBlock const& lumi, edm::EventSetup const&) const {
0331   unsigned int triggers = 0;
0332   unsigned int datasets = 0;
0333   auto const& rundata = *runCache(lumi.getRun().index());
0334   if (rundata.hltConfig.inited()) {
0335     triggers = rundata.indicesOfTriggerPaths.size();
0336     datasets = rundata.hltConfig.datasetNames().size();
0337   };
0338 
0339   // the API of jsoncollector::HistoJ does not really match our use case,
0340   // but it is the only vector-like object available in JsonMonitorable.h
0341   auto lumidata = std::make_shared<HLTriggerJSONMonitoringData::lumisection>(HLTriggerJSONMonitoringData::lumisection{
0342       jsoncollector::HistoJ<unsigned int>(1),         // processed
0343       jsoncollector::HistoJ<unsigned int>(triggers),  // hltWasRun
0344       jsoncollector::HistoJ<unsigned int>(triggers),  // hltPrePS
0345       jsoncollector::HistoJ<unsigned int>(triggers),  // hltPostPS
0346       jsoncollector::HistoJ<unsigned int>(triggers),  // hltAccept
0347       jsoncollector::HistoJ<unsigned int>(triggers),  // hltReject
0348       jsoncollector::HistoJ<unsigned int>(triggers),  // hltErrors
0349       jsoncollector::HistoJ<unsigned int>(datasets)   // datasets
0350   });
0351   // repeated calls to `update` necessary to set the internal element counter
0352   lumidata->processed.update(0);
0353   for (unsigned int i = 0; i < triggers; ++i)
0354     lumidata->hltWasRun.update(0);
0355   for (unsigned int i = 0; i < triggers; ++i)
0356     lumidata->hltPrePS.update(0);
0357   for (unsigned int i = 0; i < triggers; ++i)
0358     lumidata->hltPostPS.update(0);
0359   for (unsigned int i = 0; i < triggers; ++i)
0360     lumidata->hltAccept.update(0);
0361   for (unsigned int i = 0; i < triggers; ++i)
0362     lumidata->hltReject.update(0);
0363   for (unsigned int i = 0; i < triggers; ++i)
0364     lumidata->hltErrors.update(0);
0365   for (unsigned int i = 0; i < datasets; ++i)
0366     lumidata->datasets.update(0);
0367 
0368   return lumidata;
0369 }
0370 
0371 // called when the Stream is switching from one LuminosityBlock to a new LuminosityBlock.
0372 void HLTriggerJSONMonitoring::streamBeginLuminosityBlock(edm::StreamID sid,
0373                                                          edm::LuminosityBlock const& lumi,
0374                                                          edm::EventSetup const&) const {
0375   auto& stream = *streamCache(sid);
0376 
0377   unsigned int triggers = 0;
0378   unsigned int datasets = 0;
0379   auto const& rundata = *runCache(lumi.getRun().index());
0380   if (rundata.hltConfig.inited()) {
0381     triggers = rundata.indicesOfTriggerPaths.size();
0382     datasets = rundata.hltConfig.datasetNames().size();
0383   };
0384 
0385   // reset the stream counters
0386   stream.processed = 0;
0387   stream.hltWasRun.assign(triggers, 0);
0388   stream.hltPrePS.assign(triggers, 0);
0389   stream.hltPostPS.assign(triggers, 0);
0390   stream.hltAccept.assign(triggers, 0);
0391   stream.hltReject.assign(triggers, 0);
0392   stream.hltErrors.assign(triggers, 0);
0393   stream.datasets.assign(datasets, 0);
0394 }
0395 
0396 // called when a Stream has finished processing a LuminosityBlock, after streamEndLuminosityBlock
0397 void HLTriggerJSONMonitoring::streamEndLuminosityBlockSummary(edm::StreamID sid,
0398                                                               edm::LuminosityBlock const& lumi,
0399                                                               edm::EventSetup const&,
0400                                                               HLTriggerJSONMonitoringData::lumisection* lumidata) const {
0401   auto const& stream = *streamCache(sid);
0402   auto const& rundata = *runCache(lumi.getRun().index());
0403   lumidata->processed.value()[0] += stream.processed;
0404 
0405   // check that the HLTConfigProvider for the current run has been successfully initialised
0406   if (not rundata.hltConfig.inited())
0407     return;
0408 
0409   auto const triggers = rundata.indicesOfTriggerPaths.size();
0410   for (auto i = 0u; i < triggers; ++i) {
0411     lumidata->hltWasRun.value()[i] += stream.hltWasRun[i];
0412     lumidata->hltPrePS.value()[i] += stream.hltPrePS[i];
0413     lumidata->hltPostPS.value()[i] += stream.hltPostPS[i];
0414     lumidata->hltAccept.value()[i] += stream.hltAccept[i];
0415     lumidata->hltReject.value()[i] += stream.hltReject[i];
0416     lumidata->hltErrors.value()[i] += stream.hltErrors[i];
0417   }
0418   auto const datasets = rundata.hltConfig.datasetNames().size();
0419   for (auto i = 0u; i < datasets; ++i)
0420     lumidata->datasets.value()[i] += stream.datasets[i];
0421 }
0422 
0423 // called after the streamEndLuminosityBlockSummary method for all Streams have finished processing a given LuminosityBlock
0424 void HLTriggerJSONMonitoring::globalEndLuminosityBlockSummary(edm::LuminosityBlock const& lumi,
0425                                                               edm::EventSetup const&,
0426                                                               HLTriggerJSONMonitoringData::lumisection* lumidata) const {
0427   unsigned int ls = lumi.luminosityBlock();
0428   unsigned int run = lumi.run();
0429 
0430   bool writeFiles = true;
0431   if (edm::Service<evf::MicroStateService>().isAvailable()) {
0432     evf::FastMonitoringService* fms =
0433         (evf::FastMonitoringService*)(edm::Service<evf::MicroStateService>().operator->());
0434     if (fms)
0435       writeFiles = fms->shouldWriteFiles(ls);
0436   }
0437   if (not writeFiles)
0438     return;
0439 
0440   unsigned int processed = lumidata->processed.value().at(0);
0441   auto const& rundata = *runCache(lumi.getRun().index());
0442   Json::StyledWriter writer;
0443 
0444   // [SIC]
0445   char hostname[33];
0446   gethostname(hostname, 32);
0447   std::string sourceHost(hostname);
0448 
0449   // [SIC]
0450   std::stringstream sOutDef;
0451   sOutDef << rundata.baseRunDir << "/"
0452           << "output_" << getpid() << ".jsd";
0453 
0454   std::string jsndataFileList = "";
0455   unsigned int jsndataSize = 0;
0456   unsigned int jsndataAdler32 = 1;  // adler32 checksum for an empty file
0457 
0458   if (processed) {
0459     // write the .jsndata files which contain the actual rates
0460     Json::Value jsndata;
0461     jsndata[jsoncollector::DataPoint::SOURCE] = sourceHost;
0462     jsndata[jsoncollector::DataPoint::DEFINITION] = rundata.jsdFileName;
0463     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->processed.toJsonValue());
0464     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltWasRun.toJsonValue());
0465     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltPrePS.toJsonValue());
0466     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltPostPS.toJsonValue());
0467     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltAccept.toJsonValue());
0468     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltReject.toJsonValue());
0469     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltErrors.toJsonValue());
0470     jsndata[jsoncollector::DataPoint::DATA].append(lumidata->datasets.toJsonValue());
0471 
0472     auto jsndataFileName = fmt::sprintf("run%06d_ls%04d_streamHLTRates_pid%05d.jsndata", run, ls, getpid());
0473 
0474     std::string result = writer.write(jsndata);
0475     std::ofstream jsndataFile(rundata.baseRunDir + "/" + jsndataFileName);
0476     jsndataFile << result;
0477     jsndataFile.close();
0478 
0479     jsndataFileList = jsndataFileName;
0480     jsndataSize = result.size();
0481     jsndataAdler32 = cms::Adler32(result.c_str(), result.size());
0482   }
0483 
0484   // create a metadata json file for the "HLT rates" pseudo-stream
0485   unsigned int jsnProcessed = processed;
0486   unsigned int jsnAccepted = processed;
0487   unsigned int jsnErrorEvents = 0;
0488   unsigned int jsnRetCodeMask = 0;
0489   std::string jsnInputFiles = "";
0490   unsigned int jsnHLTErrorEvents = 0;
0491 
0492   Json::Value jsn;
0493   jsn[jsoncollector::DataPoint::SOURCE] = sourceHost;
0494   jsn[jsoncollector::DataPoint::DEFINITION] = sOutDef.str();
0495   jsn[jsoncollector::DataPoint::DATA].append(jsnProcessed);
0496   jsn[jsoncollector::DataPoint::DATA].append(jsnAccepted);
0497   jsn[jsoncollector::DataPoint::DATA].append(jsnErrorEvents);
0498   jsn[jsoncollector::DataPoint::DATA].append(jsnRetCodeMask);
0499   jsn[jsoncollector::DataPoint::DATA].append(jsndataFileList);
0500   jsn[jsoncollector::DataPoint::DATA].append(jsndataSize);
0501   jsn[jsoncollector::DataPoint::DATA].append(jsnInputFiles);
0502   jsn[jsoncollector::DataPoint::DATA].append(jsndataAdler32);
0503   jsn[jsoncollector::DataPoint::DATA].append(rundata.streamDestination);
0504   jsn[jsoncollector::DataPoint::DATA].append(rundata.streamMergeType);
0505   jsn[jsoncollector::DataPoint::DATA].append(jsnHLTErrorEvents);
0506 
0507   auto jsnFileName = fmt::sprintf("run%06d_ls%04d_streamHLTRates_pid%05d.jsn", run, ls, getpid());
0508   std::ofstream jsnFile(rundata.baseRunDir + "/" + jsnFileName);
0509   jsnFile << writer.write(jsn);
0510   jsnFile.close();
0511 }
0512 
0513 void HLTriggerJSONMonitoring::writeJsdFile(HLTriggerJSONMonitoringData::run const& rundata) {
0514   std::ofstream file(rundata.baseRunDir + "/" + rundata.jsdFileName);
0515   file << R"""({
0516    "data" : [
0517       { "name" : "Processed", "type" : "integer", "operation" : "histo"},
0518       { "name" : "Path-WasRun", "type" : "integer", "operation" : "histo"},
0519       { "name" : "Path-AfterL1Seed", "type" : "integer", "operation" : "histo"},
0520       { "name" : "Path-AfterPrescale", "type" : "integer", "operation" : "histo"},
0521       { "name" : "Path-Accepted", "type" : "integer", "operation" : "histo"},
0522       { "name" : "Path-Rejected", "type" : "integer", "operation" : "histo"},
0523       { "name" : "Path-Errors", "type" : "integer", "operation" : "histo"},
0524       { "name" : "Dataset-Accepted", "type" : "integer", "operation" : "histo"}
0525    ]
0526 }
0527 )""";
0528   file.close();
0529 }
0530 
0531 void HLTriggerJSONMonitoring::writeIniFile(HLTriggerJSONMonitoringData::run const& rundata, unsigned int run) {
0532   Json::Value content;
0533 
0534   Json::Value triggerNames(Json::arrayValue);
0535   for (auto idx : rundata.indicesOfTriggerPaths)
0536     triggerNames.append(rundata.hltConfig.triggerNames()[idx]);
0537   content["Path-Names"] = triggerNames;
0538 
0539   Json::Value datasetNames(Json::arrayValue);
0540   for (auto const& name : rundata.hltConfig.datasetNames())
0541     datasetNames.append(name);
0542   content["Dataset-Names"] = datasetNames;
0543 
0544   std::string iniFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.ini", run, getpid());
0545   std::ofstream file(rundata.baseRunDir + "/" + iniFileName);
0546   Json::StyledWriter writer;
0547   file << writer.write(content);
0548   file.close();
0549 }
0550 
0551 // declare as a framework plugin
0552 #include "FWCore/ServiceRegistry/interface/ServiceMaker.h"
0553 #include "FWCore/Framework/interface/MakerMacros.h"
0554 DEFINE_FWK_MODULE(HLTriggerJSONMonitoring);