Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl010_skim.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Example creating a derived RNTuple
5///
6/// \macro_image
7/// \macro_code
8///
9/// \date February 2024
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
23#include <cstdint>
24
25// Import classes from experimental namespace for the time being.
29
30// Input and output.
31constexpr char const *kNTupleInputName = "ntpl";
32constexpr char const *kNTupleInputFileName = "ntpl010_input.root";
33constexpr char const *kNTupleOutputName = "ntpl_skim";
34constexpr char const *kNTupleOutputFileName = "ntpl010_skim.root";
35constexpr int kNEvents = 25000;
36
37static void Write()
38{
39 auto model = RNTupleModel::Create();
40
41 auto fldVpx = model->MakeField<std::vector<float>>("vpx");
42 auto fldVpy = model->MakeField<std::vector<float>>("vpy");
43 auto fldVpz = model->MakeField<std::vector<float>>("vpz");
44 auto fldN = model->MakeField<float>("n");
45
46 auto writer = RNTupleWriter::Recreate(std::move(model), kNTupleInputName, kNTupleInputFileName);
47
49 for (int i = 0; i < kNEvents; i++) {
50 *fldN = static_cast<int>(gRandom->Rndm(1) * 15);
51
52 fldVpx->clear();
53 fldVpy->clear();
54 fldVpz->clear();
55
56 for (int j = 0; j < *fldN; ++j) {
57 float px, py, pz;
58 gRandom->Rannor(px, py);
59 pz = px * px + py * py;
60
61 fldVpx->emplace_back(px);
62 fldVpy->emplace_back(py);
63 fldVpz->emplace_back(pz);
64 }
65
66 writer->Fill();
67 }
68}
69
70void ntpl010_skim()
71{
72 Write();
73
74 auto reader = RNTupleReader::Open(kNTupleInputName, kNTupleInputFileName);
75
76 auto skimModel = RNTupleModel::Create();
77 // Loop through the top-level fields of the input RNTuple
78 for (const auto &value : reader->GetModel().GetDefaultEntry()) {
79 // Drop "n" field from skimmed dataset
80 if (value.GetField().GetFieldName() == "n")
81 continue;
82
83 // Add a renamed clone of the other fields to the skim model
84 const std::string newName = "skim_" + value.GetField().GetFieldName();
85 skimModel->AddField(value.GetField().Clone(newName));
86 // Connect input and output field
87 skimModel->GetDefaultEntry().BindValue<void>(newName, value.GetPtr<void>());
88 }
89
90 // Add an additional field to the skimmed dataset
91 auto ptrSkip = skimModel->MakeField<std::uint16_t>("skip");
92
93 auto writer = RNTupleWriter::Recreate(std::move(skimModel), kNTupleOutputName, kNTupleOutputFileName);
94
95 auto hSkip = new TH1F("h", "distribution of skipped entries", 10, 0, 10);
96 auto ptrN = reader->GetModel().GetDefaultEntry().GetPtr<float>("n");
97 for (auto numEntry : *reader) {
98 reader->LoadEntry(numEntry);
99 if (*ptrN <= 7) {
100 (*ptrSkip)++;
101 continue;
102 }
103 writer->Fill();
104 hSkip->Fill(*ptrSkip);
105 *ptrSkip = 0;
106 }
107
108 TCanvas *c1 = new TCanvas("", "Skimming Example", 200, 10, 700, 500);
109 hSkip->DrawCopy();
110}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
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
return c1
Definition legend1.C:41