Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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// Until C++ runtime modules are universally used, we explicitly load the ntuple library. Otherwise
22// triggering autoloading from the use of templated types would require an exhaustive enumeration
23// of "all" template instances in the LinkDef file.
24R__LOAD_LIBRARY(ROOTNTuple)
25
26#include <ROOT/RField.hxx>
27#include <ROOT/RNTuple.hxx>
28#include <ROOT/RNTupleModel.hxx>
29
30#include <TBranch.h>
31#include <TCanvas.h>
32#include <TFile.h>
33#include <TH1F.h>
34#include <TLeaf.h>
35#include <TTree.h>
36
37#include <cassert>
38#include <memory>
39#include <vector>
40
41// Import classes from experimental namespace for the time being
42using RNTupleModel = ROOT::Experimental::RNTupleModel;
44using RNTupleReader = ROOT::Experimental::RNTupleReader;
45using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
46
47constexpr char const* kTreeFileName = "http://root.cern.ch/files/LHCb/lhcb_B2HHH_MagnetUp.root";
48constexpr char const* kNTupleFileName = "ntpl003_lhcbOpenData.root";
49
50
51void Convert() {
52 std::unique_ptr<TFile> f(TFile::Open(kTreeFileName));
53 assert(f && ! f->IsZombie());
54
55 // Get a unique pointer to an empty RNTuple model
56 auto model = RNTupleModel::Create();
57
58 // We create RNTuple fields based on the types found in the TTree
59 // This simple approach only works for trees with simple branches and only one leaf per branch
60 auto tree = f->Get<TTree>("DecayTree");
61 for (auto b : TRangeDynCast<TBranch>(*tree->GetListOfBranches())) {
62 // The dynamic cast to TBranch should never fail for GetListOfBranches()
63 assert(b);
64
65 // We assume every branch has a single leaf
66 TLeaf *l = static_cast<TLeaf*>(b->GetListOfLeaves()->First());
67
68 // Create an ntuple field with the same name and type than the tree branch
69 auto field = RFieldBase::Create(l->GetName(), l->GetTypeName()).Unwrap();
70 std::cout << "Convert leaf " << l->GetName() << " [" << l->GetTypeName() << "]"
71 << " --> " << "field " << field->GetName() << " [" << field->GetType() << "]" << std::endl;
72
73 // Hand over ownership of the field to the ntuple model. This will also create a memory location attached
74 // to the model's default entry, that will be used to place the data supposed to be written
75 model->AddField(std::move(field));
76
77 // We connect the model's default entry's memory location for the new field to the branch, so that we can
78 // fill the ntuple with the data read from the TTree
79 void *fieldDataPtr = model->GetDefaultEntry()->GetValue(l->GetName()).GetRawPtr();
80 tree->SetBranchAddress(b->GetName(), fieldDataPtr);
81 }
82
83 // The new ntuple takes ownership of the model
84 auto ntuple = RNTupleWriter::Recreate(std::move(model), "DecayTree", kNTupleFileName);
85
86 auto nEntries = tree->GetEntries();
87 for (decltype(nEntries) i = 0; i < nEntries; ++i) {
88 tree->GetEntry(i);
89 ntuple->Fill();
90
91 if (i && i % 100000 == 0)
92 std::cout << "Wrote " << i << " entries" << std::endl;
93 }
94}
95
96
97void ntpl003_lhcbOpenData()
98{
99 Convert();
100
101 // Create histogram of the flight distance
102
103 // We open the ntuple without specifiying an explicit model first, but instead use a view on the field we are
104 // interested in
105 auto ntuple = RNTupleReader::Open("DecayTree", kNTupleFileName);
106
107 // The view wraps a read-only double value and accesses directly the ntuple's data buffers
108 auto viewFlightDistance = ntuple->GetView<double>("B_FlightDistance");
109
110 auto c = new TCanvas("c", "B Flight Distance", 200, 10, 700, 500);
111 TH1F h("h", "B Flight Distance", 200, 0, 140);
112 h.SetFillColor(48);
113
114 for (auto i : ntuple->GetEntryRange()) {
115 // Note that we do not load an entry in this loop, i.e. we avoid the memory copy of loading the data into
116 // the memory location given by the entry
117 h.Fill(viewFlightDistance(i));
118 }
119
120 h.DrawCopy();
121}
#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
#define R__LOAD_LIBRARY(LIBRARY)
Definition Rtypes.h:472
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.
Definition RNTuple.hxx:90
An RNTuple that gets filled with entries (data) and writes them to storage.
Definition RNTuple.hxx:213
A TTree is a list of TBranches.
Definition TBranch.h:89
The Canvas class.
Definition TCanvas.h:23
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.
Definition TFile.cxx:3997
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition TLeaf.h:57
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
TRangeDynCast is an adaptater class that allows the typed iteration through a TCollection.
A TTree represents a columnar dataset.
Definition TTree.h:79
constexpr const char * kNTupleFileName
Definition tree.py:1
auto * l
Definition textangle.C:4