File indexing completed on 2023-10-25 09:58:02
0001 #include "NanoAODRNTuples.h"
0002
0003 #include "DataFormats/NanoAOD/interface/MergeableCounterTable.h"
0004 #include "FWCore/Framework/interface/RunForOutput.h"
0005
0006 #include <ROOT/RNTuple.hxx>
0007 #include <ROOT/RNTupleModel.hxx>
0008 #include <ROOT/RNTupleOptions.hxx>
0009 #include <ROOT/RPageStorageFile.hxx>
0010 using ROOT::Experimental::RNTupleModel;
0011 using ROOT::Experimental::RNTupleWriteOptions;
0012 using ROOT::Experimental::RNTupleWriter;
0013 using ROOT::Experimental::Detail::RPageSinkFile;
0014
0015 #include "RNTupleFieldPtr.h"
0016 #include "SummaryTableOutputFields.h"
0017
0018 void LumiNTuple::createFields(const edm::LuminosityBlockID& id, TFile& file) {
0019 auto model = RNTupleModel::Create();
0020 m_run = RNTupleFieldPtr<UInt_t>("run", "", *model);
0021 m_luminosityBlock = RNTupleFieldPtr<UInt_t>("luminosityBlock", "", *model);
0022
0023
0024 RNTupleWriteOptions options;
0025 options.SetCompression(file.GetCompressionSettings());
0026 m_ntuple = std::make_unique<RNTupleWriter>(std::move(model),
0027 std::make_unique<RPageSinkFile>("LuminosityBlocks", file, options));
0028 }
0029
0030 void LumiNTuple::fill(const edm::LuminosityBlockID& id, TFile& file) {
0031 if (!m_ntuple) {
0032 createFields(id, file);
0033 }
0034 m_run.fill(id.run());
0035 m_luminosityBlock.fill(id.value());
0036 m_ntuple->Fill();
0037 }
0038
0039 void LumiNTuple::finalizeWrite() { m_ntuple.reset(); }
0040
0041 void RunNTuple::registerToken(const edm::EDGetToken& token) { m_tokens.push_back(token); }
0042
0043 void RunNTuple::createFields(const edm::RunForOutput& iRun, TFile& file) {
0044 auto model = RNTupleModel::Create();
0045 m_run = RNTupleFieldPtr<UInt_t>("run", "", *model);
0046
0047 edm::Handle<nanoaod::MergeableCounterTable> handle;
0048 for (const auto& token : m_tokens) {
0049 iRun.getByToken(token, handle);
0050 const nanoaod::MergeableCounterTable& tab = *handle;
0051 m_tables.push_back(SummaryTableOutputFields(tab, *model));
0052 }
0053
0054
0055 RNTupleWriteOptions options;
0056 options.SetCompression(file.GetCompressionSettings());
0057 m_ntuple = std::make_unique<RNTupleWriter>(std::move(model), std::make_unique<RPageSinkFile>("Runs", file, options));
0058 }
0059
0060 void RunNTuple::fill(const edm::RunForOutput& iRun, TFile& file) {
0061 if (!m_ntuple) {
0062 createFields(iRun, file);
0063 }
0064 m_run.fill(iRun.id().run());
0065 edm::Handle<nanoaod::MergeableCounterTable> handle;
0066 for (std::size_t i = 0; i < m_tokens.size(); i++) {
0067 iRun.getByToken(m_tokens.at(i), handle);
0068 const nanoaod::MergeableCounterTable& tab = *handle;
0069 m_tables.at(i).fill(tab);
0070 }
0071 m_ntuple->Fill();
0072 }
0073
0074 void RunNTuple::finalizeWrite() { m_ntuple.reset(); }
0075
0076 void PSetNTuple::createFields(TFile& file) {
0077
0078 auto pairModel = RNTupleModel::Create();
0079 m_psetId = RNTupleFieldPtr<std::string>("first", "", *pairModel);
0080 m_psetBlob = RNTupleFieldPtr<std::string>("second", "", *pairModel);
0081 auto model = RNTupleModel::Create();
0082 m_collection = model->MakeCollection(edm::poolNames::idToParameterSetBlobsBranchName(), std::move(pairModel));
0083
0084 RNTupleWriteOptions options;
0085 options.SetCompression(file.GetCompressionSettings());
0086 m_ntuple = std::make_unique<RNTupleWriter>(
0087 std::move(model), std::make_unique<RPageSinkFile>(edm::poolNames::parameterSetsTreeName(), file, options));
0088 }
0089
0090 void PSetNTuple::fill(edm::pset::Registry* pset, TFile& file) {
0091 if (!pset) {
0092 throw cms::Exception("LogicError", "null edm::pset::Registry::Instance pointer");
0093 }
0094 if (!m_ntuple) {
0095 createFields(file);
0096 }
0097 for (const auto& ps : *pset) {
0098 std::ostringstream oss;
0099 oss << ps.first;
0100 m_psetId.fill(oss.str());
0101 m_psetBlob.fill(ps.second.toString());
0102 m_collection->Fill();
0103 m_ntuple->Fill();
0104 }
0105 }
0106
0107 void PSetNTuple::finalizeWrite() { m_ntuple.reset(); }
0108
0109
0110 void MetadataNTuple::createFields(TFile& file) {
0111 auto procHistModel = RNTupleModel::Create();
0112
0113 m_phId = RNTupleFieldPtr<std::string>("transients_phid_", "", *procHistModel);
0114 auto model = RNTupleModel::Create();
0115 m_procHist = model->MakeCollection(edm::poolNames::processHistoryBranchName(), std::move(procHistModel));
0116 RNTupleWriteOptions options;
0117 options.SetCompression(file.GetCompressionSettings());
0118 m_ntuple = std::make_unique<RNTupleWriter>(
0119 std::move(model), std::make_unique<RPageSinkFile>(edm::poolNames::metaDataTreeName(), file, options));
0120 }
0121
0122 void MetadataNTuple::fill(const edm::ProcessHistoryRegistry& procHist, TFile& file) {
0123 if (!m_ntuple) {
0124 createFields(file);
0125 }
0126 for (const auto& ph : procHist) {
0127 std::string phid;
0128 ph.second.id().toString(phid);
0129 m_phId.fill(phid);
0130 m_procHist->Fill();
0131 }
0132 m_ntuple->Fill();
0133 }
0134
0135 void MetadataNTuple::finalizeWrite() { m_ntuple.reset(); }