Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl010_skim.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Example creating a derived RNTuple

// NOTE: The RNTuple classes are experimental at this point.
// Functionality, interface, and data format is still subject to changes.
// Do not use for real data!
#include <TCanvas.h>
#include <TH1F.h>
#include <TRandom.h>
#include <cstdint>
// Import classes from experimental namespace for the time being.
// Input and output.
constexpr char const *kNTupleInputName = "ntpl";
constexpr char const *kNTupleInputFileName = "ntpl010_input.root";
constexpr char const *kNTupleOutputName = "ntpl_skim";
constexpr char const *kNTupleOutputFileName = "ntpl010_skim.root";
constexpr int kNEvents = 25000;
static void Write()
{
auto model = RNTupleModel::Create();
auto fldVpx = model->MakeField<std::vector<float>>("vpx");
auto fldVpy = model->MakeField<std::vector<float>>("vpy");
auto fldVpz = model->MakeField<std::vector<float>>("vpz");
auto fldN = model->MakeField<float>("n");
auto writer = RNTupleWriter::Recreate(std::move(model), kNTupleInputName, kNTupleInputFileName);
for (int i = 0; i < kNEvents; i++) {
*fldN = static_cast<int>(gRandom->Rndm(1) * 15);
fldVpx->clear();
fldVpy->clear();
fldVpz->clear();
for (int j = 0; j < *fldN; ++j) {
float px, py, pz;
gRandom->Rannor(px, py);
pz = px * px + py * py;
fldVpx->emplace_back(px);
fldVpy->emplace_back(py);
fldVpz->emplace_back(pz);
}
writer->Fill();
}
}
void ntpl010_skim()
{
Write();
auto reader = RNTupleReader::Open(kNTupleInputName, kNTupleInputFileName);
auto skimModel = RNTupleModel::Create();
// Loop through the top-level fields of the input RNTuple
for (const auto &value : reader->GetModel().GetDefaultEntry()) {
// Drop "n" field from skimmed dataset
if (value.GetField().GetFieldName() == "n")
continue;
// Add a renamed clone of the other fields to the skim model
const std::string newName = "skim_" + value.GetField().GetFieldName();
skimModel->AddField(value.GetField().Clone(newName));
// Connect input and output field
skimModel->GetDefaultEntry().BindValue<void>(newName, value.GetPtr<void>());
}
// Add an additional field to the skimmed dataset
auto ptrSkip = skimModel->MakeField<std::uint16_t>("skip", 0);
auto writer = RNTupleWriter::Recreate(std::move(skimModel), kNTupleOutputName, kNTupleOutputFileName);
auto hSkip = new TH1F("h", "distribution of skipped entries", 10, 0, 10);
auto ptrN = reader->GetModel().GetDefaultEntry().GetPtr<float>("n");
for (auto numEntry : *reader) {
reader->LoadEntry(numEntry);
if (*ptrN <= 7) {
(*ptrSkip)++;
continue;
}
writer->Fill();
hSkip->Fill(*ptrSkip);
*ptrSkip = 0;
}
TCanvas *c1 = new TCanvas("", "Skimming Example", 200, 10, 700, 500);
hSkip->DrawCopy();
}
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.
Representation of an RNTuple data set in a ROOT file.
Definition RNTuple.hxx:61
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:622
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
Date
February 2024
Author
The ROOT Team

Definition in file ntpl010_skim.C.