This tutorial illustrates data conversion for a simple, tabular data model. For reading, the tutorial shows the use of an ntuple View, which selectively accesses specific fields. If a view is used for reading, there is no need to define the data model as an RNTupleModel first. The advantage of a view is that it directly accesses RNTuple's data buffers without making an additional memory copy.
#include <ROOT/RField.hxx>
#include <cassert>
#include <memory>
#include <vector>
constexpr char const* kTreeFileName = "http://root.cern.ch/files/LHCb/lhcb_B2HHH_MagnetUp.root";
constexpr char const* kNTupleFileName = "ntpl003_lhcbOpenData.root";
void Convert() {
assert(
f && !
f->IsZombie());
auto model = RNTupleModel::CreateBare();
TLeaf *
l =
static_cast<TLeaf*
>(
b->GetListOfLeaves()->First());
auto field = RFieldBase::Create(
l->
GetName(),
l->GetTypeName()).Unwrap();
std::cout <<
"Convert leaf " <<
l->
GetName() <<
" [" <<
l->GetTypeName() <<
"]"
<<
" --> " <<
"field " << field->
GetName() <<
" [" << field->GetType() <<
"]" << std::endl;
model->AddField(std::move(field));
}
auto ntuple = RNTupleWriter::Recreate(std::move(model), "DecayTree", kNTupleFileName);
auto entry = ntuple->GetModel()->CreateEntry();
auto l =
static_cast<TLeaf *
>(
b->GetListOfLeaves()->First());
tree->SetBranchAddress(
b->GetName(), fieldDataPtr);
}
auto nEntries =
tree->GetEntries();
for (decltype(nEntries) i = 0; i < nEntries; ++i) {
ntuple->Fill(*entry);
if (i && i % 100000 == 0)
std::cout << "Wrote " << i << " entries" << std::endl;
}
}
void ntpl003_lhcbOpenData()
{
Convert();
auto ntuple = RNTupleReader::Open("DecayTree", kNTupleFileName);
auto viewFlightDistance = ntuple->GetView<double>("B_FlightDistance");
auto c =
new TCanvas(
"c",
"B Flight Distance", 200, 10, 700, 500);
TH1F h(
"h",
"B Flight Distance", 200, 0, 140);
for (auto i : ntuple->GetEntryRange()) {
h.Fill(viewFlightDistance(i));
}
}
#define R__LOAD_LIBRARY(LIBRARY)
A field translates read and write calls from/to underlying columns to/from tree values.
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
An RNTuple that gets filled with entries (data) and writes them to storage.
A TTree is a list of TBranches.
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
1-D histogram with a float per channel (see TH1 documentation)}
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
virtual Double_t GetValue(Int_t i=0) const
virtual const char * GetName() const
Returns name of object.
A TTree represents a columnar dataset.