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