From 536f97234ebec974a2ae8875fe77d59ec85c2d10 Mon Sep 17 00:00:00 2001 From: Jonas Hahnfeld Date: Tue, 7 Jan 2025 16:59:38 +0100 Subject: [PATCH] [ntuple] Forbid RNTupleParallelWriter with a non-binary TFile --- tree/ntuple/v7/src/RNTupleParallelWriter.cxx | 12 ++++++++++++ tree/ntuple/v7/test/ntuple_parallel_writer.cxx | 10 ++++++++++ 2 files changed, 22 insertions(+) diff --git a/tree/ntuple/v7/src/RNTupleParallelWriter.cxx b/tree/ntuple/v7/src/RNTupleParallelWriter.cxx index 143ef581357a4..98344e4b5a965 100644 --- a/tree/ntuple/v7/src/RNTupleParallelWriter.cxx +++ b/tree/ntuple/v7/src/RNTupleParallelWriter.cxx @@ -21,7 +21,9 @@ #include #include +#include #include +#include namespace { @@ -168,6 +170,16 @@ std::unique_ptr ROOT::Experimental::RNTupleParallelWriter::Append(std::unique_ptr model, std::string_view ntupleName, TDirectory &fileOrDirectory, const RNTupleWriteOptions &options) { + auto file = fileOrDirectory.GetFile(); + if (!file) { + throw RException( + R__FAIL("RNTupleParallelWriter only supports writing to a ROOT file. Cannot write into a directory " + "that is not backed by a file")); + } + if (!file->IsBinary()) { + throw RException(R__FAIL("RNTupleParallelWriter only supports writing to a ROOT file. Cannot write into " + + std::string(file->GetName()))); + } if (!options.GetUseBufferedWrite()) { throw RException(R__FAIL("parallel writing requires buffering")); } diff --git a/tree/ntuple/v7/test/ntuple_parallel_writer.cxx b/tree/ntuple/v7/test/ntuple_parallel_writer.cxx index 99b0329245726..b59233344db58 100644 --- a/tree/ntuple/v7/test/ntuple_parallel_writer.cxx +++ b/tree/ntuple/v7/test/ntuple_parallel_writer.cxx @@ -293,6 +293,16 @@ TEST(RNTupleParallelWriter, ForbidModelWithSubfields) } } +TEST(RNTupleParallelWriter, ForbidNonRootTFiles) +{ + FileRaii fileGuard("test_ntuple_parallel_forbid_xml.xml"); + + auto model = RNTupleModel::Create(); + auto file = std::unique_ptr(TFile::Open(fileGuard.GetPath().c_str(), "RECREATE")); + // Opening an XML TFile should fail + EXPECT_THROW(RNTupleParallelWriter::Append(std::move(model), "ntpl", *file), ROOT::RException); +} + TEST(RNTupleFillContext, FlushColumns) { FileRaii fileGuard("test_ntuple_context_flush.root");