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#include <ROOT/RNTupleModel.hxx>
19
20#include <TCanvas.h>
21#include <TH1F.h>
22#include <TRandom.h>
23#include <TSystem.h>
24
25#include <cstdio>
26#include <iostream>
27#include <memory>
28#include <vector>
29#include <utility>
30
31// Import classes from experimental namespace for the time being
32using RNTupleModel = ROOT::Experimental::RNTupleModel;
33using RNTupleReader = ROOT::Experimental::RNTupleReader;
34using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
35
36// Where to store the ntuple of this example
37constexpr char const* kNTupleFileName = "ntpl002_vector.root";
38
39// Update the histogram GUI every so many fills
40constexpr int kUpdateGuiFreq = 1000;
41
42// Number of events to generate
43constexpr int kNEvents = 25000;
44
45// Generate kNEvents with vectors in kNTupleFileName
46void Write()
47{
48 // We create a unique pointer to an empty data model
49 auto model = RNTupleModel::Create();
50
51 // Creating fields of std::vector is the same as creating fields of simple types. As a result, we get
52 // shared pointers of the given type
53 std::shared_ptr<std::vector<float>> fldVpx = model->MakeField<std::vector<float>>("vpx");
54 auto fldVpy = model->MakeField<std::vector<float>>("vpy");
55 auto fldVpz = model->MakeField<std::vector<float>>("vpz");
56 auto fldVrand = model->MakeField<std::vector<float>>("vrand");
57
58 // We hand-over the data model to a newly created ntuple of name "F", stored in kNTupleFileName
59 // In return, we get a unique pointer to an ntuple that we can fill
60 auto ntuple = RNTupleWriter::Recreate(std::move(model), "F", kNTupleFileName);
61
62 TH1F hpx("hpx", "This is the px distribution", 100, -4, 4);
63 hpx.SetFillColor(48);
64
65 auto c1 = new TCanvas("c1", "Dynamic Filling Example", 200, 10, 700, 500);
66
68 for (int i = 0; i < kNEvents; i++) {
69 int npx = static_cast<int>(gRandom->Rndm(1) * 15);
70
71 fldVpx->clear();
72 fldVpy->clear();
73 fldVpz->clear();
74 fldVrand->clear();
75
76 // Set the field data for the current event
77 for (int j = 0; j < npx; ++j) {
78 float px, py, pz;
79 gRandom->Rannor(px, py);
80 pz = px*px + py*py;
81 auto random = gRandom->Rndm(1);
82
83 hpx.Fill(px);
84
85 fldVpx->emplace_back(px);
86 fldVpy->emplace_back(py);
87 fldVpz->emplace_back(pz);
88 fldVrand->emplace_back(random);
89 }
90
91 // Gui updates
92 if (i && (i % kUpdateGuiFreq) == 0) {
93 if (i == kUpdateGuiFreq) hpx.Draw();
94 c1->Modified();
95 c1->Update();
97 break;
98 }
99
100 ntuple->Fill();
101 }
102
103 hpx.DrawCopy();
104
105 // The ntuple unique pointer goes out of scope here. On destruction, the ntuple flushes unwritten data to disk
106 // and closes the attached ROOT file.
107}
108
109
110// For all of the events, histogram only one of the written vectors
111void Read()
112{
113 // Get a unique pointer to an empty RNTuple model
114 auto model = RNTupleModel::Create();
115
116 // We only define the fields that are needed for reading
117 auto fldVpx = model->MakeField<std::vector<float>>("vpx");
118
119 // Create an ntuple without imposing a specific data model. We could generate the data model from the ntuple
120 // but here we prefer the view because we only want to access a single field
121 auto ntuple = RNTupleReader::Open(std::move(model), "F", kNTupleFileName);
122
123 // Quick overview of the ntuple's key meta-data
124 ntuple->PrintInfo();
125
126 std::cout << "Entry number 42 in JSON format:" << std::endl;
127 ntuple->Show(41);
128 // In a future version of RNTuple, there will be support for ntuple->Scan()
129
130 TCanvas *c2 = new TCanvas("c2", "Dynamic Filling Example", 200, 10, 700, 500);
131 TH1F h("h", "This is the px distribution", 100, -4, 4);
132 h.SetFillColor(48);
133
134 // Iterate through all the events using i as event number and as an index for accessing the view
135 for (auto entryId : *ntuple) {
136 ntuple->LoadEntry(entryId);
137
138 for (auto px : *fldVpx) {
139 h.Fill(px);
140 }
141
142 if (entryId && (entryId % kUpdateGuiFreq) == 0) {
143 if (entryId == kUpdateGuiFreq) h.Draw();
144 c2->Modified();
145 c2->Update();
146 if (gSystem->ProcessEvents())
147 break;
148 }
149 }
150
151 // Prevent the histogram from disappearing
152 h.DrawCopy();
153}
154
155
156void ntpl002_vector()
157{
158 Write();
159 Read();
160}
#define h(i)
Definition RSha256.hxx:106
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
R__EXTERN TSystem * gSystem
Definition TSystem.h:555
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.
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:559
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:507
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:416
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14