Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-29 23:12:54

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