HLTriggerJSONMonitoring

HLTriggerJSONMonitoringData

lumisection

run

stream

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551
/** \class HLTriggerJSONMonitoring
 *
 *  Description: This class outputs JSON files with HLT monitoring information.
 *
 */

#include <atomic>
#include <fstream>
#include <sstream>
#include <vector>
#include <string>
#include <memory>
#include <algorithm>

#include <fmt/printf.h>

#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/TriggerResults.h"
#include "EventFilter/Utilities/interface/EvFDaqDirector.h"
#include "EventFilter/Utilities/interface/FastMonitor.h"
#include "EventFilter/Utilities/interface/FastMonitoringService.h"
#include "EventFilter/Utilities/interface/JSONSerializer.h"
#include "EventFilter/Utilities/interface/JsonMonitorable.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/LuminosityBlock.h"
#include "FWCore/Framework/interface/Run.h"
#include "FWCore/Framework/interface/global/EDAnalyzer.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Utilities/interface/Adler32Calculator.h"
#include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"

//note this was updated 20/10/22 to change the logic such that instead
//of having it passing the L1Seed, it is now number passing pre prescale
//this is indentical for "standard" paths and more meaningful for
//the special paths which are affected
//a standard path logic goes "trigger type -> l1 seed -> prescale -> other selection"

using namespace jsoncollector;

struct HLTriggerJSONMonitoringData {
  // variables accumulated event by event in each stream
  struct stream {
    unsigned int processed;               // number of events processed
    std::vector<unsigned int> hltWasRun;  // number of events where each path was run
    std::vector<unsigned int> hltPrePS;   // number of events where each path made it to the prescale module
    std::vector<unsigned int> hltPostPS;  // number of events where each path passed the prescale
    std::vector<unsigned int> hltAccept;  // number of events accepted by each path
    std::vector<unsigned int> hltReject;  // number of events rejected by each path
    std::vector<unsigned int> hltErrors;  // number of events with errors in each path
    std::vector<unsigned int> datasets;   // number of events accepted by each dataset
  };

  // variables initialised for each run
  struct run {
    std::string streamDestination;
    std::string streamMergeType;
    std::string baseRunDir;   // base directory from EvFDaqDirector
    std::string jsdFileName;  // definition file name for JSON with rates

    HLTConfigProvider hltConfig;  // HLT configuration for the current run
    std::vector<int> posPre;      // position of last HLT prescale filter in each path, or -1 if not present
    std::vector<std::vector<unsigned int>> datasets;  // list of paths in each dataset
    std::vector<unsigned int> indicesOfTriggerPaths;  // indices of triggers (without DatasetPaths) in TriggerNames
  };

  // variables accumulated over the whole lumisection
  struct lumisection {
    jsoncollector::HistoJ<unsigned int> processed;  // number of events processed
    jsoncollector::HistoJ<unsigned int> hltWasRun;  // number of events where each path was run
    jsoncollector::HistoJ<unsigned int> hltPrePS;   // number of events where each path made it to the prescale module
    jsoncollector::HistoJ<unsigned int> hltPostPS;  // number of events where each path passed the prescale
    jsoncollector::HistoJ<unsigned int> hltAccept;  // number of events accepted by each path
    jsoncollector::HistoJ<unsigned int> hltReject;  // number of events rejected by each path
    jsoncollector::HistoJ<unsigned int> hltErrors;  // number of events with errors in each path
    jsoncollector::HistoJ<unsigned int> datasets;   // number of events accepted by each dataset
  };
};

class HLTriggerJSONMonitoring : public edm::global::EDAnalyzer<
                                    // per-stream information
                                    edm::StreamCache<HLTriggerJSONMonitoringData::stream>,
                                    // per-run accounting
                                    edm::RunCache<HLTriggerJSONMonitoringData::run>,
                                    // accumulate per-lumisection statistics
                                    edm::LuminosityBlockSummaryCache<HLTriggerJSONMonitoringData::lumisection>> {
public:
  // constructor
  explicit HLTriggerJSONMonitoring(const edm::ParameterSet&);

  // destructor
  ~HLTriggerJSONMonitoring() override = default;

  // validate the configuration and optionally fill the default values
  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

  // called for each Event
  void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;

  // -- inherited from edm::StreamCache<HLTriggerJSONMonitoringData::stream>

  // called once for each Stream being used in the job to create the cache object that will be used for that particular Stream
  std::unique_ptr<HLTriggerJSONMonitoringData::stream> beginStream(edm::StreamID) const override;

  // called when the Stream is switching from one LuminosityBlock to a new LuminosityBlock.
  void streamBeginLuminosityBlock(edm::StreamID, edm::LuminosityBlock const&, edm::EventSetup const&) const override;

  // -- inherited from edm::RunCache<HLTriggerJSONMonitoringData::run>

  // called each time the Source sees a new Run, and guaranteed to finish before any Stream calls streamBeginRun for that same Run
  std::shared_ptr<HLTriggerJSONMonitoringData::run> globalBeginRun(edm::Run const&,
                                                                   edm::EventSetup const&) const override;

  // called after all Streams have finished processing a given Run (i.e. streamEndRun for all Streams have completed)
  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override;

  // -- inherited from edm::LuminosityBlockSummaryCache<HLTriggerJSONMonitoringData::lumisection>

  // called each time the Source sees a new LuminosityBlock
  std::shared_ptr<HLTriggerJSONMonitoringData::lumisection> globalBeginLuminosityBlockSummary(
      edm::LuminosityBlock const&, edm::EventSetup const&) const override;

  // called when a Stream has finished processing a LuminosityBlock, after streamEndLuminosityBlock
  void streamEndLuminosityBlockSummary(edm::StreamID,
                                       edm::LuminosityBlock const&,
                                       edm::EventSetup const&,
                                       HLTriggerJSONMonitoringData::lumisection*) const override;

  // called after the streamEndLuminosityBlockSummary method for all Streams have finished processing a given LuminosityBlock
  void globalEndLuminosityBlockSummary(edm::LuminosityBlock const&,
                                       edm::EventSetup const&,
                                       HLTriggerJSONMonitoringData::lumisection*) const override;

private:
  static constexpr const char* streamName_ = "streamHLTRates";

  static constexpr const char* datasetPathNamePrefix_ = "Dataset_";

  static void writeJsdFile(HLTriggerJSONMonitoringData::run const&);
  static void writeIniFile(HLTriggerJSONMonitoringData::run const&, unsigned int);

  // configuration
  const edm::InputTag triggerResults_;                               // InputTag for TriggerResults
  const edm::EDGetTokenT<edm::TriggerResults> triggerResultsToken_;  // Token for TriggerResults
};

// constructor
HLTriggerJSONMonitoring::HLTriggerJSONMonitoring(edm::ParameterSet const& config)
    : triggerResults_(config.getParameter<edm::InputTag>("triggerResults")),
      triggerResultsToken_(consumes(triggerResults_)) {
  if (edm::Service<evf::EvFDaqDirector>().isAvailable()) {
    //output initemp file. This lets hltd know number of streams early
    std::string initFileName = edm::Service<evf::EvFDaqDirector>()->getInitTempFilePath("streamHLTRates");
    std::ofstream file(initFileName);
    if (!file)
      throw cms::Exception("HLTriggerJsonMonitoring")
          << "Cannot create INITEMP file: " << initFileName << " error: " << strerror(errno);
    file.close();
  }
}

// validate the configuration and optionally fill the default values
void HLTriggerJSONMonitoring::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults", "", "@currentProcess"));
  descriptions.add("HLTriggerJSONMonitoring", desc);
}

// called once for each Stream being used in the job to create the cache object that will be used for that particular Stream
std::unique_ptr<HLTriggerJSONMonitoringData::stream> HLTriggerJSONMonitoring::beginStream(edm::StreamID) const {
  return std::make_unique<HLTriggerJSONMonitoringData::stream>();
}

// called each time the Source sees a new Run, and guaranteed to finish before any Stream calls streamBeginRun for that same Run
std::shared_ptr<HLTriggerJSONMonitoringData::run> HLTriggerJSONMonitoring::globalBeginRun(
    edm::Run const& run, edm::EventSetup const& setup) const {
  auto rundata = std::make_shared<HLTriggerJSONMonitoringData::run>();

  // set the DAQ parameters
  if (edm::Service<evf::EvFDaqDirector>().isAvailable()) {
    rundata->streamDestination = edm::Service<evf::EvFDaqDirector>()->getStreamDestinations(streamName_);
    rundata->streamMergeType =
        edm::Service<evf::EvFDaqDirector>()->getStreamMergeType(streamName_, evf::MergeTypeJSNDATA);
    rundata->baseRunDir = edm::Service<evf::EvFDaqDirector>()->baseRunDir();
  } else {
    rundata->streamDestination = "";
    rundata->streamMergeType = "";
    rundata->baseRunDir = ".";
  }

  // initialize HLTConfigProvider
  bool changed = true;
  if (not rundata->hltConfig.init(run, setup, triggerResults_.process(), changed)) {
    edm::LogError("HLTriggerJSONMonitoring") << "HLTConfigProvider initialization failed!";
  } else if (changed) {
    // triggerNames from TriggerResults (includes DatasetPaths)
    auto const& triggerNames = rundata->hltConfig.triggerNames();
    auto const triggerNamesSize = triggerNames.size();

    // update the list of indices of the HLT Paths (without DatasetPaths) in the TriggerNames list
    rundata->indicesOfTriggerPaths.clear();
    rundata->indicesOfTriggerPaths.reserve(triggerNamesSize);
    for (auto triggerNameIdx = 0u; triggerNameIdx < triggerNamesSize; ++triggerNameIdx) {
      // skip DatasetPaths
      if (triggerNames[triggerNameIdx].find(datasetPathNamePrefix_) != 0) {
        rundata->indicesOfTriggerPaths.emplace_back(triggerNameIdx);
      }
    }
    auto const triggersSize = rundata->indicesOfTriggerPaths.size();

    // update the list of paths in each dataset
    auto const& datasets = rundata->hltConfig.datasetContents();
    auto const& datasetNames = rundata->hltConfig.datasetNames();
    auto const datasetsSize = datasetNames.size();
    rundata->datasets.resize(datasetsSize);
    for (auto ds = 0u; ds < datasetsSize; ++ds) {
      auto& dataset = rundata->datasets[ds];
      // check if TriggerNames include the DatasetPath corresponding to this Dataset
      //  - DatasetPaths are normal cms.Path objects
      //  - in Run-3 HLT menus, DatasetPaths are used to define PrimaryDatasets
      auto const datasetPathName = datasetPathNamePrefix_ + datasetNames[ds];
      auto const datasetPathExists =
          std::find(triggerNames.begin(), triggerNames.end(), datasetPathName) != triggerNames.end();
      if (datasetPathExists) {
        // if a DatasetPath exists, only that Path is assigned to the Dataset
        //  - this way, the counts of the Dataset properly include prescales on the DatasetPath
        //    and smart-Prescales applied by the DatasetPath to its triggers
        dataset.reserve(1);
        auto const index = rundata->hltConfig.triggerIndex(datasetPathName);
        if (index < triggerNamesSize)
          dataset.push_back(index);
      } else {
        auto const paths = datasets[ds].size();
        dataset.reserve(paths);
        for (auto p = 0u; p < paths; p++) {
          auto const index = rundata->hltConfig.triggerIndex(datasets[ds][p]);
          if (index < triggerNamesSize)
            dataset.push_back(index);
        }
      }
    }

    // find the positions of the prescale filters
    rundata->posPre.resize(triggersSize);
    for (auto i = 0u; i < triggersSize; ++i) {
      rundata->posPre[i] = -1;
      auto const trigNameIndx = rundata->indicesOfTriggerPaths[i];
      auto const& moduleLabels = rundata->hltConfig.moduleLabels(trigNameIndx);
      for (auto j = 0u; j < moduleLabels.size(); ++j) {
        auto const& moduleType = rundata->hltConfig.moduleType(moduleLabels[j]);
        if (moduleType == "HLTPrescaler") {
          rundata->posPre[i] = j;
          break;
        }
      }
    }
  }

  // write the per-run .jsd file
  rundata->jsdFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.jsd", run.run(), getpid());
  writeJsdFile(*rundata);

  // write the per-run .ini file
  //iniFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.ini", run.run(), getpid());
  writeIniFile(*rundata, run.run());

  return rundata;
}

// called after all Streams have finished processing a given Run (i.e. streamEndRun for all Streams have completed)
void HLTriggerJSONMonitoring::globalEndRun(edm::Run const&, edm::EventSetup const&) const {}

// called for each Event
void HLTriggerJSONMonitoring::analyze(edm::StreamID sid, edm::Event const& event, edm::EventSetup const&) const {
  auto& stream = *streamCache(sid);
  auto const& rundata = *runCache(event.getRun().index());

  ++stream.processed;

  // check that the HLTConfigProvider for the current run has been successfully initialised
  if (not rundata.hltConfig.inited())
    return;

  // get hold of TriggerResults
  edm::Handle<edm::TriggerResults> handle;
  if (not event.getByToken(triggerResultsToken_, handle) or not handle.isValid()) {
    edm::LogError("HLTriggerJSONMonitoring")
        << "TriggerResults with label [" + triggerResults_.encode() + "] not present or invalid";
    return;
  }
  edm::TriggerResults const& results = *handle;
  assert(results.size() == rundata.hltConfig.triggerNames().size());

  // check the results for each HLT path
  for (auto idx = 0u; idx < rundata.indicesOfTriggerPaths.size(); ++idx) {
    auto const triggerPathIdx = rundata.indicesOfTriggerPaths[idx];
    auto const& status = results[triggerPathIdx];
    if (status.wasrun()) {
      ++stream.hltWasRun[idx];
      if (status.accept()) {
        ++stream.hltPrePS[idx];
        ++stream.hltPostPS[idx];
        ++stream.hltAccept[idx];
      } else {
        int const index = (int)status.index();
        if (index >= rundata.posPre[idx])
          ++stream.hltPrePS[idx];
        if (index > rundata.posPre[idx])
          ++stream.hltPostPS[idx];
        if (status.error())
          ++stream.hltErrors[idx];
        else
          ++stream.hltReject[idx];
      }
    }
  }

  // check the decision for each dataset
  // FIXME this ignores the prescales, "smart" prescales, and event selection applied in the OutputModule itself
  for (auto i = 0u; i < rundata.datasets.size(); ++i)
    if (std::any_of(rundata.datasets[i].begin(), rundata.datasets[i].end(), [&](unsigned int path) {
          return results.accept(path);
        }))
      ++stream.datasets[i];
}

// called each time the Source sees a new LuminosityBlock
std::shared_ptr<HLTriggerJSONMonitoringData::lumisection> HLTriggerJSONMonitoring::globalBeginLuminosityBlockSummary(
    edm::LuminosityBlock const& lumi, edm::EventSetup const&) const {
  unsigned int triggers = 0;
  unsigned int datasets = 0;
  auto const& rundata = *runCache(lumi.getRun().index());
  if (rundata.hltConfig.inited()) {
    triggers = rundata.indicesOfTriggerPaths.size();
    datasets = rundata.hltConfig.datasetNames().size();
  };

  // the API of jsoncollector::HistoJ does not really match our use case,
  // but it is the only vector-like object available in JsonMonitorable.h
  auto lumidata = std::make_shared<HLTriggerJSONMonitoringData::lumisection>(HLTriggerJSONMonitoringData::lumisection{
      jsoncollector::HistoJ<unsigned int>(1),         // processed
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltWasRun
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltPrePS
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltPostPS
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltAccept
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltReject
      jsoncollector::HistoJ<unsigned int>(triggers),  // hltErrors
      jsoncollector::HistoJ<unsigned int>(datasets)   // datasets
  });
  // repeated calls to `update` necessary to set the internal element counter
  lumidata->processed.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltWasRun.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltPrePS.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltPostPS.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltAccept.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltReject.update(0);
  for (unsigned int i = 0; i < triggers; ++i)
    lumidata->hltErrors.update(0);
  for (unsigned int i = 0; i < datasets; ++i)
    lumidata->datasets.update(0);

  return lumidata;
}

// called when the Stream is switching from one LuminosityBlock to a new LuminosityBlock.
void HLTriggerJSONMonitoring::streamBeginLuminosityBlock(edm::StreamID sid,
                                                         edm::LuminosityBlock const& lumi,
                                                         edm::EventSetup const&) const {
  auto& stream = *streamCache(sid);

  unsigned int triggers = 0;
  unsigned int datasets = 0;
  auto const& rundata = *runCache(lumi.getRun().index());
  if (rundata.hltConfig.inited()) {
    triggers = rundata.indicesOfTriggerPaths.size();
    datasets = rundata.hltConfig.datasetNames().size();
  };

  // reset the stream counters
  stream.processed = 0;
  stream.hltWasRun.assign(triggers, 0);
  stream.hltPrePS.assign(triggers, 0);
  stream.hltPostPS.assign(triggers, 0);
  stream.hltAccept.assign(triggers, 0);
  stream.hltReject.assign(triggers, 0);
  stream.hltErrors.assign(triggers, 0);
  stream.datasets.assign(datasets, 0);
}

// called when a Stream has finished processing a LuminosityBlock, after streamEndLuminosityBlock
void HLTriggerJSONMonitoring::streamEndLuminosityBlockSummary(edm::StreamID sid,
                                                              edm::LuminosityBlock const& lumi,
                                                              edm::EventSetup const&,
                                                              HLTriggerJSONMonitoringData::lumisection* lumidata) const {
  auto const& stream = *streamCache(sid);
  auto const& rundata = *runCache(lumi.getRun().index());
  lumidata->processed.value()[0] += stream.processed;

  // check that the HLTConfigProvider for the current run has been successfully initialised
  if (not rundata.hltConfig.inited())
    return;

  auto const triggers = rundata.indicesOfTriggerPaths.size();
  for (auto i = 0u; i < triggers; ++i) {
    lumidata->hltWasRun.value()[i] += stream.hltWasRun[i];
    lumidata->hltPrePS.value()[i] += stream.hltPrePS[i];
    lumidata->hltPostPS.value()[i] += stream.hltPostPS[i];
    lumidata->hltAccept.value()[i] += stream.hltAccept[i];
    lumidata->hltReject.value()[i] += stream.hltReject[i];
    lumidata->hltErrors.value()[i] += stream.hltErrors[i];
  }
  auto const datasets = rundata.hltConfig.datasetNames().size();
  for (auto i = 0u; i < datasets; ++i)
    lumidata->datasets.value()[i] += stream.datasets[i];
}

// called after the streamEndLuminosityBlockSummary method for all Streams have finished processing a given LuminosityBlock
void HLTriggerJSONMonitoring::globalEndLuminosityBlockSummary(edm::LuminosityBlock const& lumi,
                                                              edm::EventSetup const&,
                                                              HLTriggerJSONMonitoringData::lumisection* lumidata) const {
  unsigned int ls = lumi.luminosityBlock();
  unsigned int run = lumi.run();

  if (edm::Service<evf::FastMonitoringService>().isAvailable()) {
    if (!edm::Service<evf::FastMonitoringService>()->shouldWriteFiles(ls))
      return;
  }

  unsigned int processed = lumidata->processed.value().at(0);
  auto const& rundata = *runCache(lumi.getRun().index());
  Json::StyledWriter writer;

  // [SIC]
  char hostname[33];
  gethostname(hostname, 32);
  std::string sourceHost(hostname);

  // [SIC]
  std::stringstream sOutDef;
  sOutDef << rundata.baseRunDir << "/"
          << "output_" << getpid() << ".jsd";

  std::string jsndataFileList = "";
  unsigned int jsndataSize = 0;
  unsigned int jsndataAdler32 = 1;  // adler32 checksum for an empty file

  if (processed) {
    // write the .jsndata files which contain the actual rates
    Json::Value jsndata;
    jsndata[jsoncollector::DataPoint::SOURCE] = sourceHost;
    jsndata[jsoncollector::DataPoint::DEFINITION] = rundata.jsdFileName;
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->processed.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltWasRun.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltPrePS.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltPostPS.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltAccept.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltReject.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->hltErrors.toJsonValue());
    jsndata[jsoncollector::DataPoint::DATA].append(lumidata->datasets.toJsonValue());

    auto jsndataFileName = fmt::sprintf("run%06d_ls%04d_streamHLTRates_pid%05d.jsndata", run, ls, getpid());

    std::string result = writer.write(jsndata);
    std::ofstream jsndataFile(rundata.baseRunDir + "/" + jsndataFileName);
    jsndataFile << result;
    jsndataFile.close();

    jsndataFileList = jsndataFileName;
    jsndataSize = result.size();
    jsndataAdler32 = cms::Adler32(result.c_str(), result.size());
  }

  // create a metadata json file for the "HLT rates" pseudo-stream
  unsigned int jsnProcessed = processed;
  unsigned int jsnAccepted = processed;
  unsigned int jsnErrorEvents = 0;
  unsigned int jsnRetCodeMask = 0;
  std::string jsnInputFiles = "";
  unsigned int jsnHLTErrorEvents = 0;

  Json::Value jsn;
  jsn[jsoncollector::DataPoint::SOURCE] = sourceHost;
  jsn[jsoncollector::DataPoint::DEFINITION] = sOutDef.str();
  jsn[jsoncollector::DataPoint::DATA].append(jsnProcessed);
  jsn[jsoncollector::DataPoint::DATA].append(jsnAccepted);
  jsn[jsoncollector::DataPoint::DATA].append(jsnErrorEvents);
  jsn[jsoncollector::DataPoint::DATA].append(jsnRetCodeMask);
  jsn[jsoncollector::DataPoint::DATA].append(jsndataFileList);
  jsn[jsoncollector::DataPoint::DATA].append(jsndataSize);
  jsn[jsoncollector::DataPoint::DATA].append(jsnInputFiles);
  jsn[jsoncollector::DataPoint::DATA].append(jsndataAdler32);
  jsn[jsoncollector::DataPoint::DATA].append(rundata.streamDestination);
  jsn[jsoncollector::DataPoint::DATA].append(rundata.streamMergeType);
  jsn[jsoncollector::DataPoint::DATA].append(jsnHLTErrorEvents);

  auto jsnFileName = fmt::sprintf("run%06d_ls%04d_streamHLTRates_pid%05d.jsn", run, ls, getpid());
  std::ofstream jsnFile(rundata.baseRunDir + "/" + jsnFileName);
  jsnFile << writer.write(jsn);
  jsnFile.close();
}

void HLTriggerJSONMonitoring::writeJsdFile(HLTriggerJSONMonitoringData::run const& rundata) {
  std::ofstream file(rundata.baseRunDir + "/" + rundata.jsdFileName);
  file << R"""({
   "data" : [
      { "name" : "Processed", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-WasRun", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-AfterL1Seed", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-AfterPrescale", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-Accepted", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-Rejected", "type" : "integer", "operation" : "histo"},
      { "name" : "Path-Errors", "type" : "integer", "operation" : "histo"},
      { "name" : "Dataset-Accepted", "type" : "integer", "operation" : "histo"}
   ]
}
)""";
  file.close();
}

void HLTriggerJSONMonitoring::writeIniFile(HLTriggerJSONMonitoringData::run const& rundata, unsigned int run) {
  Json::Value content;

  Json::Value triggerNames(Json::arrayValue);
  for (auto idx : rundata.indicesOfTriggerPaths)
    triggerNames.append(rundata.hltConfig.triggerNames()[idx]);
  content["Path-Names"] = triggerNames;

  Json::Value datasetNames(Json::arrayValue);
  for (auto const& name : rundata.hltConfig.datasetNames())
    datasetNames.append(name);
  content["Dataset-Names"] = datasetNames;

  std::string iniFileName = fmt::sprintf("run%06d_ls0000_streamHLTRates_pid%05d.ini", run, getpid());
  std::ofstream file(rundata.baseRunDir + "/" + iniFileName);
  Json::StyledWriter writer;
  file << writer.write(content);
  file.close();
}

// declare as a framework plugin
#include "FWCore/ServiceRegistry/interface/ServiceMaker.h"
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(HLTriggerJSONMonitoring);