Write and read tabular data with RNTuple.
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>
ifstream
fin(
gROOT->GetTutorialDir() +
"/io/tree/cernstaff.dat");
auto fldFlag = model->MakeField<
unsigned int>(
"Flag");
auto fldAge = model->MakeField<
int>(
"Age");
auto fldService = model->MakeField<
int>(
"Service");
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");
}
}
void Analyze() {
std::shared_ptr<int>
fldAge = model->MakeField<
int>(
"Age");
std::cout << "The first entry in JSON format:" << std::endl;
auto c =
new TCanvas(
"c",
"", 200, 10, 700, 500);
TH1I h(
"h",
"Age Distribution CERN, 1988", 100, 0, 100);
}
}
Analyze();
}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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())
Throws an exception if the model is null.
1-D histogram with an int per channel (see TH1 documentation)
- Date
- April 2019
- Author
- The ROOT Team
Definition in file ntpl001_staff.C.