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#include <ROOT/RField.hxx>
22#include <ROOT/RNTupleModel.hxx>
25
26#include <TBranch.h>
27#include <TCanvas.h>
28#include <TFile.h>
29#include <TH1F.h>
30#include <TLeaf.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;
39using RFieldBase = ROOT::Experimental::RFieldBase;
40using RNTupleReader = ROOT::Experimental::RNTupleReader;
41using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
42
43constexpr char const* kTreeFileName = "http://root.cern/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 without a default entry
52 auto model = RNTupleModel::CreateBare();
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()).Unwrap();
66 std::cout << "Convert leaf " << l->GetName() << " [" << l->GetTypeName() << "]"
67 << " --> "
68 << "field " << field->GetFieldName() << " [" << field->GetTypeName() << "]" << std::endl;
69
70 // Hand over ownership of the field to the ntuple model. This will also create a memory location attached
71 // to the model's default entry, that will be used to place the data supposed to be written
72 model->AddField(std::move(field));
73 }
74
75 // The new ntuple takes ownership of the model
76 auto ntuple = RNTupleWriter::Recreate(std::move(model), "DecayTree", kNTupleFileName);
77
78 auto entry = ntuple->GetModel().CreateEntry();
79 for (auto b : TRangeDynCast<TBranch>(*tree->GetListOfBranches())) {
80 auto l = static_cast<TLeaf *>(b->GetListOfLeaves()->First());
81 // We connect the model's default entry's memory location for the new field to the branch, so that we can
82 // fill the ntuple with the data read from the TTree
83 std::shared_ptr<void> fieldDataPtr = entry->GetPtr<void>(l->GetName());
84 tree->SetBranchAddress(b->GetName(), fieldDataPtr.get());
85 }
86
87 auto nEntries = tree->GetEntries();
88 for (decltype(nEntries) i = 0; i < nEntries; ++i) {
89 tree->GetEntry(i);
90 ntuple->Fill(*entry);
91
92 if (i && i % 100000 == 0)
93 std::cout << "Wrote " << i << " entries" << std::endl;
94 }
95}
96
97
98void ntpl003_lhcbOpenData()
99{
100 Convert();
101
102 // Create histogram of the flight distance
103
104 // We open the ntuple without specifying an explicit model first, but instead use a view on the field we are
105 // interested in
106 auto ntuple = RNTupleReader::Open("DecayTree", kNTupleFileName);
107
108 // The view wraps a read-only double value and accesses directly the ntuple's data buffers
109 auto viewFlightDistance = ntuple->GetView<double>("B_FlightDistance");
110
111 auto c = new TCanvas("c", "B Flight Distance", 200, 10, 700, 500);
112 TH1F h("h", "B Flight Distance", 200, 0, 140);
113 h.SetFillColor(48);
114
115 for (auto i : ntuple->GetEntryRange()) {
116 // Note that we do not load an entry in this loop, i.e. we avoid the memory copy of loading the data into
117 // the memory location given by the entry
118 h.Fill(viewFlightDistance(i));
119 }
120
121 h.DrawCopy();
122}
#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
A field translates read and write calls from/to underlying columns to/from tree values.
Definition RField.hxx:95
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.
Definition TBranch.h:93
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:4082
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
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 TObject.cxx:439
A TTree represents a columnar dataset.
Definition TTree.h:79
TLine l
Definition textangle.C:4