Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl005_introspection.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Write and read an RNTuple from a user-defined class. Adapted from tv3.C
5/// Illustrates various RNTuple introspection methods.
6///
7/// \macro_image
8/// \macro_code
9///
10/// \date April 2020
11/// \author The ROOT Team
12
13// NOTE: The RNTuple classes are experimental at this point.
14// Functionality and interface are still subject to changes.
15
16#include <ROOT/RNTupleModel.hxx>
20
21#include <Compression.h>
22#include <TCanvas.h>
23#include <TH1.h>
24#include <TRandom.h>
25#include <TSystem.h>
26
27#include <cassert>
28
29// 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
35constexpr char const* kNTupleFileName = "ntpl005_introspection.root";
36
37// Store entries of type Vector3 in the ntuple
38class Vector3 {
39private:
40 double fX = 0;
41 double fY = 0;
42 double fZ = 0;
43
44public:
45 double x() const { return fX; }
46 double y() const { return fY; }
47 double z() const { return fZ; }
48
49 void SetXYZ(double x, double y, double z) {
50 fX = x;
51 fY = y;
52 fZ = z;
53 }
54};
55
56
57void Generate()
58{
59 auto model = RNTupleModel::Create();
60 auto fldVector3 = model->MakeField<Vector3>("v3");
61
62 // Explicitly enforce a certain compression algorithm
65
66 auto ntuple = RNTupleWriter::Recreate(std::move(model), "Vector3", kNTupleFileName, options);
67 TRandom r;
68 for (unsigned int i = 0; i < 500000; ++i) {
69 fldVector3->SetXYZ(r.Gaus(0,1), r.Landau(0,1), r.Gaus(100,10));
70 ntuple->Fill();
71 }
72}
73
74
76 Generate();
77
78 auto ntuple = RNTupleReader::Open("Vector3", kNTupleFileName);
79
80 // Display the schema of the ntuple
81 ntuple->PrintInfo();
82
83 // Display information about the storage layout of the data
84 ntuple->PrintInfo(ENTupleInfo::kStorageDetails);
85
86 // Display the first entry
87 ntuple->Show(0);
88
89 // Collect I/O runtime counters when processing the data set.
90 // Maintaining the counters comes with a small performance overhead, so it has to be explicitly enabled
91 ntuple->EnableMetrics();
92
93 // Plot the y components of vector3
94 TCanvas *c1 = new TCanvas("c1","RNTuple Demo", 10, 10, 600, 800);
95 c1->Divide(1,2);
96 c1->cd(1);
97 TH1F h1("x", "x component of Vector3", 100, -3, 3);
98 {
99 /// We enclose viewX in a scope in order to indicate to the RNTuple when we are not
100 /// anymore interested in v3.fX
101 auto viewX = ntuple->GetView<double>("v3.fX");
102 for (auto i : ntuple->GetEntryRange()) {
103 h1.Fill(viewX(i));
104 }
105 }
106 h1.DrawCopy();
107
108 c1->cd(2);
109 TH1F h2("y", "y component of Vector3", 100, -5, 20);
110 auto viewY = ntuple->GetView<double>("v3.fY");
111 for (auto i : ntuple->GetEntryRange()) {
112 h2.Fill(viewY(i));
113 }
114 h2.DrawCopy();
115
116 // Display the I/O operation statistics performed by the RNTuple reader
117 ntuple->PrintInfo(ENTupleInfo::kMetrics);
118
119 // We read 2 out of the 3 Vector3 members and thus should have requested approximately 2/3 of the file
122 assert(retval == 0);
123 float fileSize = static_cast<float>(fileStat.fSize);
124 float nbytesRead = ntuple->GetMetrics().GetCounter("RNTupleReader.RPageSourceFile.szReadPayload")->GetValueAsInt() +
125 ntuple->GetMetrics().GetCounter("RNTupleReader.RPageSourceFile.szReadOverhead")->GetValueAsInt();
126
127 std::cout << "File size: " << fileSize / 1024. / 1024. << " MiB" << std::endl;
128 std::cout << "Read from file: " << nbytesRead / 1024. / 1024. << " MiB" << std::endl;
129 std::cout << "Ratio: " << nbytesRead / fileSize << std::endl;
130}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
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.
Common user-tunable settings for storing ntuples.
void SetCompression(std::uint32_t val)
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:647
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3315
virtual TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const
Copy this histogram and Draw in the current pad.
Definition TH1.cxx:3084
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
int GetPathInfo(const char *path, Long_t *id, Long_t *size, Long_t *flags, Long_t *modtime)
Get info about a file: id, size, flags, modification time.
Definition TSystem.cxx:1398
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
ENTupleInfo
Listing of the different options that can be printed by RNTupleReader::GetInfo()
@ kUseGeneralPurpose
Use the new recommended general-purpose setting; it is a best trade-off between compression ratio/dec...
Definition Compression.h:58