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