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#include <ROOT/RNTupleModel.hxx>
15
16#include <TCanvas.h>
17#include <TH1F.h>
18#include <TRandom.h>
19
20#include <cstdint>
21
22// Input and output.
23constexpr char const *kNTupleInputName = "ntpl";
24constexpr char const *kNTupleInputFileName = "ntpl010_input.root";
25constexpr char const *kNTupleOutputName = "ntpl_skim";
26constexpr char const *kNTupleOutputFileName = "ntpl010_skim.root";
27constexpr int kNEvents = 25000;
28
29static void Write()
30{
31 auto model = ROOT::RNTupleModel::Create();
32
33 auto fldVpx = model->MakeField<std::vector<float>>("vpx");
34 auto fldVpy = model->MakeField<std::vector<float>>("vpy");
35 auto fldVpz = model->MakeField<std::vector<float>>("vpz");
36 auto fldN = model->MakeField<float>("n");
37
39
41 for (int i = 0; i < kNEvents; i++) {
42 *fldN = static_cast<int>(gRandom->Rndm(1) * 15);
43
44 fldVpx->clear();
45 fldVpy->clear();
46 fldVpz->clear();
47
48 for (int j = 0; j < *fldN; ++j) {
49 float px, py, pz;
50 gRandom->Rannor(px, py);
51 pz = px * px + py * py;
52
53 fldVpx->emplace_back(px);
54 fldVpy->emplace_back(py);
55 fldVpz->emplace_back(pz);
56 }
57
58 writer->Fill();
59 }
60}
61
62void ntpl010_skim()
63{
64 Write();
65
67
69 // Loop through the top-level fields of the input RNTuple
70 for (const auto &value : reader->GetModel().GetDefaultEntry()) {
71 // Drop "n" field from skimmed dataset
72 if (value.GetField().GetFieldName() == "n")
73 continue;
74
75 // Add a renamed clone of the other fields to the skim model
76 const std::string newName = "skim_" + value.GetField().GetFieldName();
77 skimModel->AddField(value.GetField().Clone(newName));
78 // Connect input and output field
79 skimModel->GetDefaultEntry().BindValue<void>(newName, value.GetPtr<void>());
80 }
81
82 // Add an additional field to the skimmed dataset
83 auto ptrSkip = skimModel->MakeField<std::uint16_t>("skip");
84
86
87 auto hSkip = new TH1F("h", "distribution of skipped entries", 10, 0, 10);
88 auto ptrN = reader->GetModel().GetDefaultEntry().GetPtr<float>("n");
89 for (auto numEntry : *reader) {
90 reader->LoadEntry(numEntry);
91 if (*ptrN <= 7) {
92 (*ptrSkip)++;
93 continue;
94 }
95 writer->Fill();
96 hSkip->Fill(*ptrSkip);
97 *ptrSkip = 0;
98 }
99
100 TCanvas *c1 = new TCanvas("", "Skimming Example", 200, 10, 700, 500);
101 hSkip->DrawCopy();
102}
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 value
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
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
return c1
Definition legend1.C:41