Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
ntpl001_staff.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_ntuple
3 /// \notebook
4 /// Write and read tabular data with RNTuple. Adapted from the cernbuild and cernstaff tree tutorials.
5 /// Illustrates the type-safe ntuple model interface, which is used to define a data model that is in a second step
6 /// taken by an ntuple reader or writer.
7 ///
8 /// \macro_image
9 /// \macro_code
10 ///
11 /// \date April 2019
12 /// \author The ROOT Team
13 
14 // NOTE: The RNTuple classes are experimental at this point.
15 // Functionality, interface, and data format is still subject to changes.
16 // Do not use for real data!
17 
18 // The following line should disappear in a future version of RNTuple, when
19 // the common template specializations of RField are part of the LinkDef.h
20 R__LOAD_LIBRARY(ROOTNTuple)
21 
22 #include <ROOT/RNTuple.hxx>
23 #include <ROOT/RNTupleModel.hxx>
24 
25 #include <TCanvas.h>
26 #include <TH1I.h>
27 #include <TROOT.h>
28 #include <TString.h>
29 
30 #include <cassert>
31 #include <cstdio>
32 #include <fstream>
33 #include <iostream>
34 #include <memory>
35 #include <string>
36 #include <sstream>
37 #include <utility>
38 
39 // Import classes from experimental namespace for the time being
40 using RNTupleModel = ROOT::Experimental::RNTupleModel;
41 using RNTupleReader = ROOT::Experimental::RNTupleReader;
42 using RNTupleWriter = ROOT::Experimental::RNTupleWriter;
43 
44 constexpr char const* kNTupleFileName = "ntpl001_staff.root";
45 
46 void Ingest() {
47  // The input file cernstaff.dat is a copy of the CERN staff data base from 1988
48  ifstream fin(gROOT->GetTutorialDir() + "/tree/cernstaff.dat");
49  assert(fin.is_open());
50 
51  // We create a unique pointer to an empty data model
52  auto model = RNTupleModel::Create();
53 
54  // To define the data model, we create fields with a given C++ type and name. Fields are roughly TTree branches.
55  // MakeField returns a shared pointer to a memory location that we can populate to fill the ntuple with data
56  auto fldCategory = model->MakeField<int>("Category");
57  auto fldFlag = model->MakeField<unsigned int>("Flag");
58  auto fldAge = model->MakeField<int>("Age");
59  auto fldService = model->MakeField<int>("Service");
60  auto fldChildren = model->MakeField<int>("Children");
61  auto fldGrade = model->MakeField<int>("Grade");
62  auto fldStep = model->MakeField<int>("Step");
63  auto fldHrweek = model->MakeField<int>("Hrweek");
64  auto fldCost = model->MakeField<int>("Cost");
65  auto fldDivision = model->MakeField<std::string>("Division");
66  auto fldNation = model->MakeField<std::string>("Nation");
67 
68  // We hand-over the data model to a newly created ntuple of name "Staff", stored in kNTupleFileName
69  // In return, we get a unique pointer to an ntuple that we can fill
70  auto ntuple = RNTupleWriter::Recreate(std::move(model), "Staff", kNTupleFileName);
71 
72  std::string record;
73  while (std::getline(fin, record)) {
74  std::istringstream iss(record);
75  iss >> *fldCategory >> *fldFlag >> *fldAge >> *fldService >> *fldChildren >> *fldGrade >> *fldStep >> *fldHrweek
76  >> *fldCost >> *fldDivision >> *fldNation;
77  ntuple->Fill();
78  }
79 
80  // The ntuple unique pointer goes out of scope here. On destruction, the ntuple flushes unwritten data to disk
81  // and closes the attached ROOT file.
82 }
83 
84 void Analyze() {
85  // Get a unique pointer to an empty RNTuple model
86  auto model = RNTupleModel::Create();
87 
88  // We only define the fields that are needed for reading
89  std::shared_ptr<int> fldAge = model->MakeField<int>("Age");
90 
91  // Create an ntuple and attach the read model to it
92  auto ntuple = RNTupleReader::Open(std::move(model), "Staff", kNTupleFileName);
93 
94  // Quick overview of the ntuple and list of fields.
95  ntuple->PrintInfo();
96  // In a future version of RNTuple, there will be support for ntuple->Show() and ntuple->Scan()
97 
98  auto c = new TCanvas("c", "", 200, 10, 700, 500);
99  TH1I h("h", "Age Distribution CERN, 1988", 100, 0, 100);
100  h.SetFillColor(48);
101 
102  for (auto entryId : *ntuple) {
103  // Populate fldAge
104  ntuple->LoadEntry(entryId);
105  h.Fill(*fldAge);
106  }
107 
108  h.DrawCopy();
109 }
110 
111 void ntpl001_staff() {
112  Ingest();
113  Analyze();
114 }