Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl002_vector.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Write and read STL vectors with RNTuple. Adapted from the hvector tree tutorial.
5///
6/// \macro_image
7/// \macro_code
8///
9/// \date April 2019
10/// \author The ROOT Team
11
12// NOTE: The RNTuple classes are experimental at this point.
13// Functionality, interface, and data format is still subject to changes.
14// Do not use for real data!
15
16// Until C++ runtime modules are universally used, we explicitly load the ntuple library. Otherwise
17// triggering autoloading from the use of templated types would require an exhaustive enumeration
18// of "all" template instances in the LinkDef file.
19R__LOAD_LIBRARY(ROOTNTuple)
20
21#include <ROOT/RNTuple.hxx>
22#include <ROOT/RNTupleModel.hxx>
23
24#include <TCanvas.h>
25#include <TH1F.h>
26#include <TRandom.h>
27#include <TSystem.h>
28
29#include <cstdio>
30#include <iostream>
31#include <memory>
32#include <vector>
33#include <utility>
34
35// Import classes from experimental namespace for the time being
36using RNTupleModel = ROOT::Experimental::RNTupleModel;
37using RNTupleReader = ROOT::Experimental::RNTupleReader;
38using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
39
40// Where to store the ntuple of this example
41constexpr char const* kNTupleFileName = "ntpl002_vector.root";
42
43// Update the histogram GUI every so many fills
44constexpr int kUpdateGuiFreq = 1000;
45
46// Number of events to generate
47constexpr int kNEvents = 25000;
48
49// Generate kNEvents with vectors in kNTupleFileName
50void Write()
51{
52 // We create a unique pointer to an empty data model
53 auto model = RNTupleModel::Create();
54
55 // Creating fields of std::vector is the same as creating fields of simple types. As a result, we get
56 // shared pointers of the given type
57 std::shared_ptr<std::vector<float>> fldVpx = model->MakeField<std::vector<float>>("vpx");
58 auto fldVpy = model->MakeField<std::vector<float>>("vpy");
59 auto fldVpz = model->MakeField<std::vector<float>>("vpz");
60 auto fldVrand = model->MakeField<std::vector<float>>("vrand");
61
62 // We hand-over the data model to a newly created ntuple of name "F", stored in kNTupleFileName
63 // In return, we get a unique pointer to an ntuple that we can fill
64 auto ntuple = RNTupleWriter::Recreate(std::move(model), "F", kNTupleFileName);
65
66 TH1F hpx("hpx", "This is the px distribution", 100, -4, 4);
67 hpx.SetFillColor(48);
68
69 auto c1 = new TCanvas("c1", "Dynamic Filling Example", 200, 10, 700, 500);
70
72 for (int i = 0; i < kNEvents; i++) {
73 int npx = static_cast<int>(gRandom->Rndm(1) * 15);
74
75 fldVpx->clear();
76 fldVpy->clear();
77 fldVpz->clear();
78 fldVrand->clear();
79
80 // Set the field data for the current event
81 for (int j = 0; j < npx; ++j) {
82 float px, py, pz;
83 gRandom->Rannor(px, py);
84 pz = px*px + py*py;
85 auto random = gRandom->Rndm(1);
86
87 hpx.Fill(px);
88
89 fldVpx->emplace_back(px);
90 fldVpy->emplace_back(py);
91 fldVpz->emplace_back(pz);
92 fldVrand->emplace_back(random);
93 }
94
95 // Gui updates
96 if (i && (i % kUpdateGuiFreq) == 0) {
97 if (i == kUpdateGuiFreq) hpx.Draw();
98 c1->Modified();
99 c1->Update();
100 if (gSystem->ProcessEvents())
101 break;
102 }
103
104 ntuple->Fill();
105 }
106
107 hpx.DrawCopy();
108
109 // The ntuple unique pointer goes out of scope here. On destruction, the ntuple flushes unwritten data to disk
110 // and closes the attached ROOT file.
111}
112
113
114// For all of the events, histogram only one of the written vectors
115void Read()
116{
117 // Get a unique pointer to an empty RNTuple model
118 auto model = RNTupleModel::Create();
119
120 // We only define the fields that are needed for reading
121 auto fldVpx = model->MakeField<std::vector<float>>("vpx");
122
123 // Create an ntuple without imposing a specific data model. We could generate the data model from the ntuple
124 // but here we prefer the view because we only want to access a single field
125 auto ntuple = RNTupleReader::Open(std::move(model), "F", kNTupleFileName);
126
127 // Quick overview of the ntuple's key meta-data
128 ntuple->PrintInfo();
129
130 std::cout << "Entry number 42 in JSON format:" << std::endl;
131 ntuple->Show(41);
132 // In a future version of RNTuple, there will be support for ntuple->Scan()
133
134 TCanvas *c2 = new TCanvas("c2", "Dynamic Filling Example", 200, 10, 700, 500);
135 TH1F h("h", "This is the px distribution", 100, -4, 4);
136 h.SetFillColor(48);
137
138 // Iterate through all the events using i as event number and as an index for accessing the view
139 for (auto entryId : *ntuple) {
140 ntuple->LoadEntry(entryId);
141
142 for (auto px : *fldVpx) {
143 h.Fill(px);
144 }
145
146 if (entryId && (entryId % kUpdateGuiFreq) == 0) {
147 if (entryId == kUpdateGuiFreq) h.Draw();
148 c2->Modified();
149 c2->Update();
150 if (gSystem->ProcessEvents())
151 break;
152 }
153 }
154
155 // Prevent the histogram from disappearing
156 h.DrawCopy();
157}
158
159
160void ntpl002_vector()
161{
162 Write();
163 Read();
164}
#define h(i)
Definition RSha256.hxx:106
#define R__LOAD_LIBRARY(LIBRARY)
Definition Rtypes.h:491
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
R__EXTERN TSystem * gSystem
Definition TSystem.h:560
The RNTupleModel encapulates the schema of an ntuple.
An RNTuple that is used to read data from storage.
Definition RNTuple.hxx:110
An RNTuple that gets filled with entries (data) and writes them to storage.
Definition RNTuple.hxx:368
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:608
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:552
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition TRandom.cxx:500
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:419
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14