Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
ntpl002_vector.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Write and read STL vectors with RNTuple.

Adapted from the hvector tree tutorial.

#include <TCanvas.h>
#include <TH1F.h>
#include <TRandom.h>
#include <TSystem.h>
#include <cstdio>
#include <iostream>
#include <memory>
#include <vector>
#include <utility>
// Where to store the ntuple of this example
constexpr char const* kNTupleFileName = "ntpl002_vector.root";
// Update the histogram GUI every so many fills
constexpr int kUpdateGuiFreq = 1000;
// Number of events to generate
constexpr int kNEvents = 25000;
// Generate kNEvents with vectors in kNTupleFileName
void Write()
{
// We create a unique pointer to an empty data model
// Creating fields of std::vector is the same as creating fields of simple types. As a result, we get
// shared pointers of the given type
std::shared_ptr<std::vector<float>> fldVpx = model->MakeField<std::vector<float>>("vpx");
auto fldVpy = model->MakeField<std::vector<float>>("vpy");
auto fldVpz = model->MakeField<std::vector<float>>("vpz");
auto fldVrand = model->MakeField<std::vector<float>>("vrand");
// We hand-over the data model to a newly created ntuple of name "F", stored in kNTupleFileName
// In return, we get a unique pointer to an ntuple that we can fill
auto writer = ROOT::RNTupleWriter::Recreate(std::move(model), "F", kNTupleFileName);
TH1F hpx("hpx", "This is the px distribution", 100, -4, 4);
hpx.SetFillColor(48);
auto c1 = new TCanvas("c1", "Dynamic Filling Example", 200, 10, 700, 500);
for (int i = 0; i < kNEvents; i++) {
int npx = static_cast<int>(gRandom->Rndm(1) * 15);
fldVpx->clear();
fldVpy->clear();
fldVpz->clear();
fldVrand->clear();
// Set the field data for the current event
for (int j = 0; j < npx; ++j) {
float px, py, pz;
gRandom->Rannor(px, py);
pz = px*px + py*py;
auto random = gRandom->Rndm(1);
hpx.Fill(px);
fldVpx->emplace_back(px);
fldVpy->emplace_back(py);
fldVpz->emplace_back(pz);
fldVrand->emplace_back(random);
}
// Gui updates
if (i && (i % kUpdateGuiFreq) == 0) {
if (i == kUpdateGuiFreq) hpx.Draw();
c1->Modified();
c1->Update();
break;
}
writer->Fill();
}
hpx.DrawCopy();
// The ntuple unique pointer goes out of scope here. On destruction, the ntuple flushes unwritten data to disk
// and closes the attached ROOT file.
}
// For all of the events, histogram only one of the written vectors
void Read()
{
// Get a unique pointer to an empty RNTuple model
// We only define the fields that are needed for reading
auto fldVpx = model->MakeField<std::vector<float>>("vpx");
// Create an ntuple without imposing a specific data model. We could generate the data model from the ntuple
// but here we prefer the view because we only want to access a single field
auto reader = ROOT::RNTupleReader::Open(std::move(model), "F", kNTupleFileName);
// Quick overview of the ntuple's key meta-data
reader->PrintInfo();
std::cout << "Entry number 42 in JSON format:" << std::endl;
reader->Show(41);
TCanvas *c2 = new TCanvas("c2", "Dynamic Filling Example", 200, 10, 700, 500);
TH1F h("h", "This is the px distribution", 100, -4, 4);
h.SetFillColor(48);
// Iterate through all the events using i as event number and as an index for accessing the view
for (auto entryId : *reader) {
reader->LoadEntry(entryId);
for (auto px : *fldVpx) {
h.Fill(px);
}
if (entryId && (entryId % kUpdateGuiFreq) == 0) {
if (entryId == kUpdateGuiFreq) h.Draw();
c2->Modified();
c2->Update();
break;
}
}
// Prevent the histogram from disappearing
h.DrawCopy();
}
{
Write();
Read();
}
#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:650
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
Date
April 2019
Author
The ROOT Team

Definition in file ntpl002_vector.C.