#include <cstdint>
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++) {
fldVpx->clear();
fldVpy->clear();
fldVpz->clear();
for (int j = 0; j < *fldN; ++j) {
float px, py, pz;
pz = px * px + py * py;
fldVpx->emplace_back(px);
fldVpy->emplace_back(py);
fldVpz->emplace_back(pz);
}
}
}
void ntpl010_skim()
{
Write();
auto reader = RNTupleReader::Open(kNTupleInputName, kNTupleInputFileName);
auto skimModel = RNTupleModel::Create();
for (
const auto &
value : reader->GetModel().GetDefaultEntry()) {
if (
value.GetField().GetFieldName() ==
"n")
continue;
const std::string newName =
"skim_" +
value.GetField().GetFieldName();
skimModel->AddField(
value.GetField().Clone(newName));
skimModel->GetDefaultEntry().BindValue<
void>(newName,
value.GetPtr<
void>());
}
auto ptrSkip = skimModel->MakeField<std::uint16_t>("skip");
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;
}
hSkip->Fill(*ptrSkip);
*ptrSkip = 0;
}
hSkip->DrawCopy();
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
R__EXTERN TRandom * gRandom
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.
1-D histogram with a float per channel (see TH1 documentation)
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Double_t Rndm() override
Machine independent random number generator.
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.