46 static constexpr std::uint64_t gUpdateFrequencyBytes = 50 * 1000 * 1000;
47 std::uint64_t fNbytesNext = gUpdateFrequencyBytes;
50 ~RDefaultProgressCallback()
override {}
51 void Call(std::uint64_t nbytesWritten, std::uint64_t neventsWritten)
final
54 if (nbytesWritten < fNbytesNext)
56 std::cout <<
"Wrote " << nbytesWritten / 1000 / 1000 <<
"MB, " << neventsWritten <<
" entries" << std::endl;
57 fNbytesNext += gUpdateFrequencyBytes;
58 if (nbytesWritten > fNbytesNext) {
60 fNbytesNext = nbytesWritten + gUpdateFrequencyBytes;
64 void Finish(std::uint64_t nbytesWritten, std::uint64_t neventsWritten)
final
66 std::cout <<
"Done, wrote " << nbytesWritten / 1000 / 1000 <<
"MB, " << neventsWritten <<
" entries" << std::endl;
89std::unique_ptr<ROOT::Experimental::RNTupleImporter>
91 std::string_view destFileName)
93 auto importer = std::unique_ptr<RNTupleImporter>(
new RNTupleImporter());
94 importer->fNTupleName = treeName;
95 importer->fSourceFile = std::unique_ptr<TFile>(
TFile::Open(std::string(sourceFileName).c_str()));
96 if (!importer->fSourceFile || importer->fSourceFile->IsZombie()) {
97 throw RException(
R__FAIL(
"cannot open source file " + std::string(sourceFileName)));
100 importer->fSourceTree = importer->fSourceFile->Get<
TTree>(std::string(treeName).c_str());
101 if (!importer->fSourceTree) {
102 throw RException(
R__FAIL(
"cannot read TTree " + std::string(treeName) +
" from " + std::string(sourceFileName)));
106 importer->fSourceTree->SetImplicitMT(
false);
107 auto result = importer->InitDestination(destFileName);
115std::unique_ptr<ROOT::Experimental::RNTupleImporter>
118 auto importer = std::unique_ptr<RNTupleImporter>(
new RNTupleImporter());
125 importer->fNTupleName = sourceTree->
GetName();
128 importer->fSourceTree = sourceTree;
132 auto result = importer->InitDestination(destFileName);
154 std::cout <<
"Importing '" <<
f.fField->GetFieldName() <<
"' [" <<
f.fField->GetTypeName() <<
']' << std::endl;
177 const auto firstLeaf =
static_cast<TLeaf *
>(
b->GetListOfLeaves()->First());
180 const bool isLeafList =
b->GetNleaves() > 1;
181 const bool isCountLeaf = firstLeaf->IsRange();
183 if (isLeafList && isClass)
184 return R__FAIL(
"unsupported: classes in leaf list, branch " + std::string(
b->GetName()));
185 if (isLeafList && isCountLeaf)
186 return R__FAIL(
"unsupported: count leaf arrays in leaf list, branch " + std::string(
b->GetName()));
191 Int_t firstLeafCountval;
192 const bool isCString = !isLeafList && (firstLeaf->IsA() ==
TLeafC::Class()) &&
193 (!firstLeaf->GetLeafCounter(firstLeafCountval)) && (firstLeafCountval == 1);
202 c.fMaxLength = firstLeaf->GetMaximum();
203 c.fCountVal = std::make_unique<Int_t>();
210 std::size_t branchBufferSize = 0;
212 std::vector<std::unique_ptr<RFieldBase>> recordItems;
215 return R__FAIL(
"unsupported: TObject branches, branch: " + std::string(
b->GetName()));
220 Int_t countval =
l->GetLenStatic();
221 auto *countleaf =
l->GetLeafCount();
222 const bool isLeafCountArray = (countleaf !=
nullptr);
223 const bool isFixedSizeArray = !isCString && (countleaf ==
nullptr) && (countval > 1);
228 std::string fieldName =
b->GetName();
229 std::string fieldType =
l->GetTypeName();
235 fieldType =
"std::string";
238 fieldType =
b->GetClassName();
240 if (isFixedSizeArray)
241 fieldType =
"std::array<" + fieldType +
"," + std::to_string(countval) +
">";
245 std::replace(fieldName.begin(), fieldName.end(),
'.',
'_');
253 auto field = fieldOrError.Unwrap();
255 branchBufferSize =
l->GetMaximum();
256 f.fValue = std::make_unique<RFieldBase::RValue>(field->CreateValue());
257 f.fFieldBuffer =
f.fValue->GetPtr<
void>().get();
264 branchBufferSize =
sizeof(
void *) * countval;
265 }
else if (isLeafCountArray) {
268 branchBufferSize =
l->GetOffset() + field->GetValueSize();
271 f.fField = field.get();
274 recordItems.emplace_back(std::move(field));
275 }
else if (isLeafCountArray) {
276 f.fValue = std::make_unique<RFieldBase::RValue>(field->CreateValue());
277 f.fFieldBuffer =
f.fValue->GetPtr<
void>().get();
278 f.fIsInUntypedCollection =
true;
279 const std::string countleafName = countleaf->GetName();
286 fModel->AddField(std::move(field));
290 if (!recordItems.empty()) {
291 auto recordField = std::make_unique<RRecordField>(
b->GetName(), std::move(recordItems));
295 fModel->AddField(std::move(recordField));
300 ib.
fBranchBuffer = std::make_unique<unsigned char[]>(branchBufferSize);
304 return R__FAIL(
"unable to load class " + std::string(
b->GetClassName()) +
" for branch " +
305 std::string(
b->GetName()));
307 auto ptrBuf =
reinterpret_cast<void **
>(ib.
fBranchBuffer.get());
322 int iLeafCountCollection = 0;
326 auto &countLeafName =
p.first;
328 c.fCollectionModel->Freeze();
329 c.fCollectionEntry =
c.fCollectionModel->CreateBareEntry();
330 for (
auto idx :
c.fImportFieldIndexes) {
333 c.fCollectionEntry->BindRawPtr(
name, buffer);
335 c.fFieldName =
"_collection" + std::to_string(iLeafCountCollection);
336 c.fCollectionWriter =
fModel->MakeCollection(
c.fFieldName, std::move(
c.fCollectionModel));
338 for (
auto idx :
c.fImportFieldIndexes) {
340 auto projectedField =
343 fModel->AddProjectedField(std::move(projectedField), [&
name, &
c](
const std::string &fieldName) {
344 if (fieldName ==
name)
347 return c.fFieldName +
"." +
name;
351 auto projectedField =
352 RFieldBase::Create(countLeafName,
"ROOT::Experimental::RNTupleCardinality<std::uint32_t>").Unwrap();
353 fModel->AddProjectedField(std::move(projectedField), [&
c](
const std::string &) {
return c.fFieldName; });
354 iLeafCountCollection++;
360 if (
f.fIsInUntypedCollection)
362 fEntry->BindRawPtr(
f.fField->GetFieldName(),
f.fFieldBuffer);
365 fEntry->BindRawPtr<
void>(
c.fFieldName,
c.fCollectionWriter->GetOffsetPtr());
381 std::unique_ptr<Internal::RPageSink> sink =
383 sink->GetMetrics().Enable();
384 auto ctrZippedBytes = sink->GetMetrics().GetCounter(
"RPageSinkFile.szWritePayload");
387 sink = std::make_unique<Internal::RPageSinkBuf>(std::move(sink));
402 for (
decltype(nEntries) i = 0; i < nEntries; ++i) {
406 for (
Int_t l = 0;
l < *
c.fCountVal; ++
l) {
407 for (
auto &t :
c.fTransformations) {
412 c.fCollectionWriter->Fill(*
c.fCollectionEntry);
414 for (
auto &t :
c.fTransformations)
425 ntplWriter->Fill(*
fEntry);
#define R__FORWARD_ERROR(res)
Short-hand to return an RResult<T> in an error state (i.e. after checking)
#define R__FAIL(msg)
Short-hand to return an RResult<T> in an error state; the RError is implicitly converted into RResult...
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Base class for all ROOT issued exceptions.
static RResult< std::unique_ptr< RFieldBase > > Create(const std::string &fieldName, const std::string &canonicalType, const std::string &typeAlias, bool fContinueOnError=false)
Factory method to resurrect a field from the stored on-disk type information.
virtual size_t GetValueSize() const =0
The number of bytes taken by a value of the appropriate type.
Used to report every ~50MB (compressed), and at the end about the status of the import.
virtual void Finish(std::uint64_t nbytesWritten, std::uint64_t neventsWritten)=0
virtual void Call(std::uint64_t nbytesWritten, std::uint64_t neventsWritten)=0
bool fConvertDotsInBranchNames
Whether or not dot characters in branch names should be converted to underscores.
std::unique_ptr< TFile > fDestFile
RNTupleWriteOptions fWriteOptions
std::int64_t fMaxEntries
The maximum number of entries to import. When this value is -1 (default), import all entries.
std::map< std::string, RImportLeafCountCollection > fLeafCountCollections
Maps the count leaf to the information about the corresponding untyped collection.
RNTupleImporter()=default
std::vector< RImportBranch > fImportBranches
static std::unique_ptr< RNTupleImporter > Create(std::string_view sourceFileName, std::string_view treeName, std::string_view destFileName)
Opens the input file for reading and the output file for writing (update).
std::unique_ptr< RProgressCallback > fProgressCallback
RResult< void > PrepareSchema()
Sets up the connection from TTree branches to RNTuple fields, including initialization of the memory ...
ROOT::Experimental::RResult< void > InitDestination(std::string_view destFileName)
void Import()
Import works in two steps:
std::unique_ptr< REntry > fEntry
std::string fDestFileName
bool fIsQuiet
No standard output, conversely if set to false, schema information and progress is printed.
std::vector< RImportField > fImportFields
std::unique_ptr< RNTupleModel > fModel
std::vector< std::unique_ptr< RImportTransformation > > fImportTransformations
The list of transformations to be performed for every entry.
static std::unique_ptr< RNTupleModel > CreateBare()
A bare model has no default entry.
bool GetUseBufferedWrite() const
The type-erased field for a RVec<Type>
The class is used as a return type for operations that can fail; wraps a value of type T or an RError...
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
TClass * IsA() const override
const char * GetName() const override
Returns name of object.
virtual const char * GetName() const
Returns name of object.
A TTree represents a columnar dataset.
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=nullptr)
Change branch address, dealing with clone trees properly.
virtual Long64_t GetEntries() const
virtual void SetImplicitMT(bool enabled)
virtual TObjArray * GetListOfBranches()
virtual TTree * GetTree() const
virtual Long64_t LoadTree(Long64_t entry)
Set current entry.
TClass * IsA() const override
std::unique_ptr< RNTupleWriter > CreateRNTupleWriter(std::unique_ptr< RNTupleModel > model, std::unique_ptr< Internal::RPageSink > sink)
std::string fBranchName
Top-level branch name from the input TTree.
std::unique_ptr< unsigned char[]> fBranchBuffer
The destination of SetBranchAddress() for fBranchName
void * fFieldBuffer
Usually points to the corresponding RImportBranch::fBranchBuffer but not always.
RFieldBase * fField
The field is kept during schema preparation and transferred to the fModel before the writing starts.
bool fIsClass
Field imported from a branch with stramer info (e.g., STL, user-defined class)
When the schema is set up and the import started, it needs to be reset before the next Import() call ...
Leaf count arrays require special treatment.
std::unique_ptr< RNTupleModel > fCollectionModel
The model for the collection itself.