Adapted from the cernbuild and cernstaff tree tutorials. Illustrates the type-safe ntuple model interface, which is used to define a data model that is in a second step taken by an ntuple reader or writer.
#include <cassert>
#include <cstdio>
#include <fstream>
#include <iostream>
#include <memory>
#include <string>
#include <sstream>
#include <utility>
constexpr char const* kNTupleFileName = "ntpl001_staff.root";
void Ingest() {
ifstream fin(
gROOT->GetTutorialDir() +
"/io/tree/cernstaff.dat");
assert(fin.is_open());
auto fldCategory = model->MakeField<int>("Category");
auto fldFlag = model->MakeField<unsigned int>("Flag");
auto fldAge = model->MakeField<int>("Age");
auto fldService = model->MakeField<int>("Service");
auto fldChildren = model->MakeField<int>("Children");
auto fldGrade = model->MakeField<int>("Grade");
auto fldStep = model->MakeField<int>("Step");
auto fldHrweek = model->MakeField<int>("Hrweek");
auto fldCost = model->MakeField<int>("Cost");
auto fldDivision = model->MakeField<std::string>("Division");
auto fldNation = model->MakeField<std::string>("Nation");
std::string record;
while (std::getline(fin, record)) {
std::istringstream iss(record);
iss >> *fldCategory >> *fldFlag >> *fldAge >> *fldService >> *fldChildren >> *fldGrade >> *fldStep >> *fldHrweek
>> *fldCost >> *fldDivision >> *fldNation;
}
}
void Analyze() {
std::shared_ptr<int> fldAge = model->MakeField<int>("Age");
reader->PrintInfo();
std::cout << "The first entry in JSON format:" << std::endl;
reader->Show(0);
auto c =
new TCanvas(
"c",
"", 200, 10, 700, 500);
TH1I h(
"h",
"Age Distribution CERN, 1988", 100, 0, 100);
for (auto entryId : *reader) {
reader->LoadEntry(entryId);
}
}
void ntpl001_staff() {
Ingest();
Analyze();
}
static std::unique_ptr< RNTupleModel > Create()
static std::unique_ptr< RNTupleReader > Open(std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleReadOptions &options=ROOT::RNTupleReadOptions())
Open an RNTuple for reading.
static std::unique_ptr< RNTupleWriter > Recreate(std::unique_ptr< ROOT::RNTupleModel > model, std::string_view ntupleName, std::string_view storage, const ROOT::RNTupleWriteOptions &options=ROOT::RNTupleWriteOptions())
Creates an RNTupleWriter backed by storage, overwriting it if one with the same URI exists.
1-D histogram with an int per channel (see TH1 documentation)