Logo ROOT   6.18/05
Reference Guide
ntpl003_lhcbOpenData.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Convert LHCb run 1 open data from a TTree to RNTuple.
5/// This tutorial illustrates data conversion for a simple, tabular data model.
6/// For reading, the tutorial shows the use of an ntuple View, which selectively accesses specific fields.
7/// If a view is used for reading, there is no need to define the data model as an RNTupleModel first.
8/// The advantage of a view is that it directly accesses RNTuple's data buffers without making an additional
9/// memory copy.
10///
11/// \macro_image
12/// \macro_code
13///
14/// \date April 2019
15/// \author The ROOT Team
16
17// NOTE: The RNTuple classes are experimental at this point.
18// Functionality, interface, and data format is still subject to changes.
19// Do not use for real data!
20
21#include <ROOT/RField.hxx>
22#include <ROOT/RNTuple.hxx>
23#include <ROOT/RNTupleModel.hxx>
24
25#include <TBranch.h>
26#include <TCanvas.h>
27#include <TFile.h>
28#include <TH1F.h>
29#include <TLeaf.h>
30#include <TSystem.h>
31#include <TTree.h>
32
33#include <cassert>
34#include <memory>
35#include <vector>
36
37// Import classes from experimental namespace for the time being
38using RNTupleModel = ROOT::Experimental::RNTupleModel;
40using RNTupleReader = ROOT::Experimental::RNTupleReader;
41using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
42
43constexpr char const* kTreeFileName = "http://root.cern.ch/files/LHCb/lhcb_B2HHH_MagnetUp.root";
44constexpr char const* kNTupleFileName = "ntpl003_lhcbOpenData.root";
45
46
47void Convert() {
48 std::unique_ptr<TFile> f(TFile::Open(kTreeFileName));
49 assert(f && ! f->IsZombie());
50
51 // Get a unique pointer to an empty RNTuple model
52 auto model = RNTupleModel::Create();
53
54 // We create RNTuple fields based on the types found in the TTree
55 // This simple approach only works for trees with simple branches and only one leaf per branch
56 auto tree = f->Get<TTree>("DecayTree");
57 for (auto b : TRangeDynCast<TBranch>(*tree->GetListOfBranches())) {
58 // The dynamic cast to TBranch should never fail for GetListOfBranches()
59 assert(b);
60
61 // We assume every branch has a single leaf
62 TLeaf *l = static_cast<TLeaf*>(b->GetListOfLeaves()->First());
63
64 // Create an ntuple field with the same name and type than the tree branch
65 auto field = RFieldBase::Create(l->GetName(), l->GetTypeName());
66 std::cout << "Convert leaf " << l->GetName() << " [" << l->GetTypeName() << "]"
67 << " --> " << "field " << field->GetName() << " [" << field->GetType() << "]" << std::endl;
68
69 // Hand over ownership of the field to the ntuple model. This will also create a memory location attached
70 // to the model's default entry, that will be used to place the data supposed to be written
71 model->AddField(std::unique_ptr<RFieldBase>(field));
72
73 // We connect the model's default entry's memory location for the new field to the branch, so that we can
74 // fill the ntuple with the data read from the TTree
75 void *fieldDataPtr = model->GetDefaultEntry()->GetValue(l->GetName()).GetRawPtr();
76 tree->SetBranchAddress(b->GetName(), fieldDataPtr);
77 }
78
79 // The new ntuple takes ownership of the model
80 auto ntuple = RNTupleWriter::Recreate(std::move(model), "DecayTree", kNTupleFileName);
81
82 auto nEntries = tree->GetEntries();
83 for (decltype(nEntries) i = 0; i < nEntries; ++i) {
84 tree->GetEntry(i);
85 ntuple->Fill();
86
87 if (i && i % 100000 == 0)
88 std::cout << "Wrote " << i << " entries" << std::endl;
89 }
90}
91
92
93void ntpl003_lhcbOpenData()
94{
95 if (gSystem->AccessPathName(kNTupleFileName))
96 Convert();
97
98 // Create histogram of the flight distance
99
100 // We open the ntuple without specifiying an explicit model first, but instead use a view on the field we are
101 // interested in
102 auto ntuple = RNTupleReader::Open("DecayTree", kNTupleFileName);
103
104 // The view wraps a read-only double value and accesses directly the ntuple's data buffers
105 auto viewFlightDistance = ntuple->GetView<double>("B_FlightDistance");
106
107 auto c = new TCanvas("c", "B Flight Distance", 200, 10, 700, 500);
108 TH1F h("h", "B Flight Distance", 200, 0, 140);
109 h.SetFillColor(48);
110
111 for (auto i : ntuple->GetViewRange()) {
112 // Note that we do not load an entry in this loop, i.e. we avoid the memory copy of loading the data into
113 // the memory location given by the entry
114 h.Fill(viewFlightDistance(i));
115 }
116
117 h.DrawCopy();
118}
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define c(i)
Definition: RSha256.hxx:101
#define h(i)
Definition: RSha256.hxx:106
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
A field translates read and write calls from/to underlying columns to/from tree values.
Definition: RField.hxx:53
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
Definition: RNTuple.hxx:94
An RNTuple that gets filled with entries (data) and writes them to storage.
Definition: RNTuple.hxx:170
The Canvas class.
Definition: TCanvas.h:31
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3980
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:49
TRangeDynCast is an adaptater class that allows the typed iteration through a TCollection.
Definition: TCollection.h:411
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1286
A TTree represents a columnar dataset.
Definition: TTree.h:71
Definition: tree.py:1
auto * l
Definition: textangle.C:4