Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TBranch.cxx
Go to the documentation of this file.
1// @(#)root/tree:$Id$
2// Author: Rene Brun 12/01/96
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "TBranchCacheInfo.h"
13
14#include "TBranch.h"
15
16#include "Bytes.h"
17#include "Compression.h"
18#include "TBasket.h"
19#include "TBranchBrowsable.h"
20#include "TBrowser.h"
21#include "TBuffer.h"
22#include "TClass.h"
23#include "TBufferFile.h"
24#include "TClonesArray.h"
25#include "TFile.h"
26#include "TLeaf.h"
27#include "TLeafB.h"
28#include "TLeafC.h"
29#include "TLeafD.h"
30#include "TLeafD32.h"
31#include "TLeafF.h"
32#include "TLeafF16.h"
33#include "TLeafI.h"
34#include "TLeafL.h"
35#include "TLeafG.h"
36#include "TLeafO.h"
37#include "TLeafObject.h"
38#include "TLeafS.h"
39#include "TMessage.h"
40#include "TROOT.h"
41#include "TSystem.h"
42#include "TMath.h"
43#include "TTree.h"
44#include "TTreeCache.h"
45#include "TTreeCacheUnzip.h"
46#include "TVirtualMutex.h"
47#include "TVirtualPad.h"
48#include "TVirtualPerfStats.h"
49#include "strlcpy.h"
50#include "snprintf.h"
51
52#include "TBranchIMTHelper.h"
53
54#include "ROOT/TIOFeatures.hxx"
55
56#include <atomic>
57#include <cstddef>
58#include <cstring>
59#include <cstdio>
60
61
63
64/** \class TBranch
65\ingroup tree
66
67A TTree is a list of TBranches
68
69A TBranch supports:
70 - The list of TLeaf describing this branch.
71 - The list of TBasket (branch buffers).
72
73See TBranch structure in TTree.
74
75See also specialized branches:
76 - TBranchObject in case the branch is one object
77 - TBranchClones in case the branch is an array of clone objects
78*/
79
80
81
82
83////////////////////////////////////////////////////////////////////////////////
84/// Default constructor. Used for I/O by default.
85
87: TNamed()
88, TAttFill(0, 1001)
89, fCompress(0)
90, fBasketSize(32000)
91, fEntryOffsetLen(1000)
92, fWriteBasket(0)
93, fEntryNumber(0)
94, fExtraBasket(nullptr)
95, fOffset(0)
96, fMaxBaskets(10)
97, fNBaskets(0)
98, fSplitLevel(0)
99, fNleaves(0)
100, fReadBasket(0)
101, fReadEntry(-1)
102, fFirstBasketEntry(-1)
103, fNextBasketEntry(-1)
104, fCurrentBasket(nullptr)
105, fEntries(0)
106, fFirstEntry(0)
107, fTotBytes(0)
108, fZipBytes(0)
109, fBranches()
110, fLeaves()
111, fBaskets(fMaxBaskets)
112, fBasketBytes(nullptr)
113, fBasketEntry(nullptr)
114, fBasketSeek(nullptr)
115, fTree(nullptr)
116, fMother(nullptr)
117, fParent(nullptr)
118, fAddress(nullptr)
119, fDirectory(nullptr)
120, fFileName("")
121, fEntryBuffer(nullptr)
122, fTransientBuffer(nullptr)
123, fBrowsables(nullptr)
124, fBulk(*this)
125, fSkipZip(false)
126, fReadLeaves(&TBranch::ReadLeavesImpl)
127, fFillLeaves(&TBranch::FillLeavesImpl)
128{
130}
131
132////////////////////////////////////////////////////////////////////////////////
133/// Create a Branch as a child of a Tree
134///
135/// * address is the address of the first item of a structure
136/// or the address of a pointer to an object (see example in TTree.cxx).
137/// * leaflist is the concatenation of all the variable names and types
138/// separated by a colon character :
139/// The variable name and the variable type are separated by a
140/// slash (/). The variable type must be 1 character. (Characters
141/// after the first are legal and will be appended to the visible
142/// name of the leaf, but have no effect.) If no type is given, the
143/// type of the variable is assumed to be the same as the previous
144/// variable. If the first variable does not have a type, it is
145/// assumed of type F by default. The list of currently supported
146/// types is given below:
147/// - `C` : a character string terminated by the 0 character
148/// - `B` : an 8 bit signed integer (`Char_t`); Treated as a character when in an array.
149/// - `b` : an 8 bit unsigned integer (`UChar_t`)
150/// - `S` : a 16 bit signed integer (`Short_t`)
151/// - `s` : a 16 bit unsigned integer (`UShort_t`)
152/// - `I` : a 32 bit signed integer (`Int_t`)
153/// - `i` : a 32 bit unsigned integer (`UInt_t`)
154/// - `F` : a 32 bit floating point (`Float_t`)
155/// - `f` : a 24 bit floating point with truncated mantissa (`Float16_t`)
156/// - `D` : a 64 bit floating point (`Double_t`)
157/// - `d` : a 24 bit truncated floating point (`Double32_t`)
158/// - `L` : a 64 bit signed integer (`Long64_t`)
159/// - `l` : a 64 bit unsigned integer (`ULong64_t`)
160/// - `G` : a long signed integer (`Long_t`, which `sizeof` is platform dependent), stored as a 64 bit integer but usually held in memory as a 64 bit integer on 64 bit machines and 32 bit on 32 bit machines. Due to this difference, this data type is **not cross-platform**.
161/// - `g` : a long unsigned integer (`ULong_t`, which `sizeof` is platform dependent), stored as a 64 bit unsigned integer but held in memory usually as a 64 bit integer on 64 bit machines and 32 bit on 32 bit machines. Due to this difference, this data type is **not cross-platform**.
162/// - `O` : [the letter `o`, not a zero] a boolean (`bool`)
163///
164/// Arrays of values are supported with the following syntax:
165/// - If leaf name has the form var[nelem], where nelem is alphanumeric, then
166/// if nelem is a leaf name, it is used as the variable size of the array,
167/// otherwise return 0.
168/// The leaf referred to by nelem **MUST** be an int (/I),
169/// - If leaf name has the form var[nelem], where nelem is a non-negative integers, then
170/// it is used as the fixed size of the array.
171/// - If leaf name has the form of a multi dimension array (e.g. var[nelem][nelem2])
172/// where nelem and nelem2 are non-negative integers) then
173/// it is used as a 2 dimensional array of fixed size.
174/// - In case of the truncated floating point types (Float16_t and Double32_t) you can
175/// furthermore specify the range in the style [xmin,xmax] or [xmin,xmax,nbits] after
176/// the type character. See `TStreamerElement::GetRange()` for further information.
177/// - Any of other form is not supported.
178///
179/// Note that the TTree will assume that all the item are contiguous in memory.
180/// On some platform, this is not always true of the member of a struct or a class,
181/// due to padding and alignment. Sorting your data member in order of decreasing
182/// sizeof usually leads to their being contiguous in memory.
183///
184/// * bufsize is the buffer size in bytes for this branch
185/// The default value is 32000 bytes and should be ok for most cases.
186/// You can specify a larger value (e.g. 256000) if your Tree is not split
187/// and each entry is large (Megabytes)
188/// A small value for bufsize is optimum if you intend to access
189/// the entries in the Tree randomly and your Tree is in split mode.
190///
191/// See an example of a Branch definition in the TTree constructor.
192///
193/// Note that in case the data type is an object, this branch can contain
194/// only this object.
195///
196/// Note that this function is invoked by TTree::Branch
197
198TBranch::TBranch(TTree *tree, const char *name, void *address, const char *leaflist, Int_t basketsize, Int_t compress)
200, TAttFill(0, 1001)
201, fCompress(compress)
202, fBasketSize((basketsize < 100) ? 100 : basketsize)
203, fEntryOffsetLen(0)
204, fWriteBasket(0)
205, fEntryNumber(0)
206, fExtraBasket(nullptr)
207, fIOFeatures(tree ? tree->GetIOFeatures().GetFeatures() : 0)
208, fOffset(0)
209, fMaxBaskets(10)
210, fNBaskets(0)
211, fSplitLevel(0)
212, fNleaves(0)
213, fReadBasket(0)
214, fReadEntry(-1)
215, fFirstBasketEntry(-1)
216, fNextBasketEntry(-1)
217, fCurrentBasket(nullptr)
218, fEntries(0)
219, fFirstEntry(0)
220, fTotBytes(0)
221, fZipBytes(0)
222, fBranches()
223, fLeaves()
224, fBaskets(fMaxBaskets)
225, fBasketBytes(nullptr)
226, fBasketEntry(nullptr)
227, fBasketSeek(nullptr)
228, fTree(tree)
229, fMother(nullptr)
230, fParent(nullptr)
231, fAddress((char *)address)
232, fDirectory(fTree->GetDirectory())
233, fFileName("")
234, fEntryBuffer(nullptr)
235, fTransientBuffer(nullptr)
236, fBrowsables(nullptr)
237, fBulk(*this)
238, fSkipZip(false)
239, fReadLeaves(&TBranch::ReadLeavesImpl)
240, fFillLeaves(&TBranch::FillLeavesImpl)
241{
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Create a Branch as a child of another Branch
247///
248/// See documentation for
249/// TBranch::TBranch(TTree *, const char *, void *, const char *, Int_t, Int_t)
250
251TBranch::TBranch(TBranch *parent, const char *name, void *address, const char *leaflist, Int_t basketsize,
254, TAttFill(0, 1001)
255, fCompress(compress)
256, fBasketSize((basketsize < 100) ? 100 : basketsize)
257, fEntryOffsetLen(0)
258, fWriteBasket(0)
259, fEntryNumber(0)
260, fExtraBasket(nullptr)
261, fIOFeatures(parent->fIOFeatures)
262, fOffset(0)
263, fMaxBaskets(10)
264, fNBaskets(0)
265, fSplitLevel(0)
266, fNleaves(0)
267, fReadBasket(0)
268, fReadEntry(-1)
269, fFirstBasketEntry(-1)
270, fNextBasketEntry(-1)
271, fCurrentBasket(nullptr)
272, fEntries(0)
273, fFirstEntry(0)
274, fTotBytes(0)
275, fZipBytes(0)
276, fBranches()
277, fLeaves()
278, fBaskets(fMaxBaskets)
279, fBasketBytes(nullptr)
280, fBasketEntry(nullptr)
281, fBasketSeek(nullptr)
282, fTree(parent ? parent->GetTree() : nullptr)
283, fMother(parent ? parent->GetMother() : nullptr)
284, fParent(parent)
285, fAddress((char *)address)
286, fDirectory(fTree ? fTree->GetDirectory() : nullptr)
287, fFileName("")
288, fEntryBuffer(nullptr)
289, fTransientBuffer(nullptr)
290, fBrowsables(nullptr)
291, fBulk(*this)
292, fSkipZip(false)
293, fReadLeaves(&TBranch::ReadLeavesImpl)
294, fFillLeaves(&TBranch::FillLeavesImpl)
295{
297}
298
299void TBranch::Init(const char* name, const char* leaflist, Int_t compress)
300{
301 // Initialization routine called from the constructor. This should NOT be made virtual.
302
304 if ((compress == -1) && fTree->GetDirectory()) {
306 if (bfile) {
307 fCompress = bfile->GetCompressionSettings();
308 }
309 }
310
314
315 for (Int_t i = 0; i < fMaxBaskets; ++i) {
316 fBasketBytes[i] = 0;
317 fBasketEntry[i] = 0;
318 fBasketSeek[i] = 0;
319 }
320
321 //
322 // Decode the leaflist (search for : as separator).
323 //
324
325 char* nameBegin = const_cast<char*>(leaflist);
326 Int_t offset = 0;
327 auto len = strlen(leaflist);
328 // FIXME: Make these string streams instead.
329 char* leafname = new char[len + 1];
330 char* leaftype = new char[320];
331 // Note: The default leaf type is a float.
332 strlcpy(leaftype, "F",320);
333 char* pos = const_cast<char*>(leaflist);
334 const char* leaflistEnd = leaflist + len;
335 for (; pos <= leaflistEnd; ++pos) {
336 // -- Scan leaf specification and create leaves.
337 if ((*pos == ':') || (*pos == 0)) {
338 // -- Reached end of a leaf spec, create a leaf.
339 Int_t lenName = pos - nameBegin;
340 char* ctype = nullptr;
341 if (lenName) {
343 leafname[lenName] = 0;
344 ctype = strstr(leafname, "/");
345 if (ctype) {
346 *ctype = 0;
347 strlcpy(leaftype, ctype + 1,320);
348 }
349 }
350 if (lenName == 0 || ctype == leafname) {
351 Warning("TBranch","No name was given to the leaf number '%d' in the leaflist of the branch '%s'.",fNleaves,name);
352 snprintf(leafname, len + 1, "__noname%d", fNleaves);
353 }
354 TLeaf* leaf = nullptr;
355 if (leaftype[1] == '[' && !strchr(leaftype, ',')) {
356 Warning("TBranch", "Array size for branch '%s' must be specified after leaf name, not after the type name!", name);
357 // and continue for backward compatibility?
358 } else if (leaftype[1] && !strchr(leaftype, ',')) {
359 Warning("TBranch", "Extra characters after type tag '%s' for branch '%s'; must be one character.", leaftype, name);
360 // and continue for backward compatibility?
361 }
362 if (*leaftype == 'C') {
363 leaf = new TLeafC(this, leafname, leaftype);
364 } else if (*leaftype == 'O') {
365 leaf = new TLeafO(this, leafname, leaftype);
366 } else if (*leaftype == 'B') {
367 leaf = new TLeafB(this, leafname, leaftype);
368 } else if (*leaftype == 'b') {
369 leaf = new TLeafB(this, leafname, leaftype);
370 leaf->SetUnsigned();
371 } else if (*leaftype == 'S') {
372 leaf = new TLeafS(this, leafname, leaftype);
373 } else if (*leaftype == 's') {
374 leaf = new TLeafS(this, leafname, leaftype);
375 leaf->SetUnsigned();
376 } else if (*leaftype == 'I') {
377 leaf = new TLeafI(this, leafname, leaftype);
378 } else if (*leaftype == 'i') {
379 leaf = new TLeafI(this, leafname, leaftype);
380 leaf->SetUnsigned();
381 } else if (*leaftype == 'F') {
382 leaf = new TLeafF(this, leafname, leaftype);
383 } else if (*leaftype == 'f') {
384 leaf = new TLeafF16(this, leafname, leaftype);
385 } else if (*leaftype == 'L') {
386 leaf = new TLeafL(this, leafname, leaftype);
387 } else if (*leaftype == 'l') {
388 leaf = new TLeafL(this, leafname, leaftype);
389 leaf->SetUnsigned();
390 } else if (*leaftype == 'D') {
391 leaf = new TLeafD(this, leafname, leaftype);
392 } else if (*leaftype == 'd') {
393 leaf = new TLeafD32(this, leafname, leaftype);
394 } else if (*leaftype == 'G') {
395 leaf = new TLeafG(this, leafname, leaftype);
396 } else if (*leaftype == 'g') {
397 leaf = new TLeafG(this, leafname, leaftype);
398 leaf->SetUnsigned();
399 }
400 if (!leaf) {
401 Error("TLeaf", "Illegal data type for %s/%s", name, leaflist);
402 delete[] leaftype;
403 delete [] leafname;
404 MakeZombie();
405 return;
406 }
407 if (leaf->IsZombie()) {
408 delete leaf;
409 leaf = nullptr;
410 auto msg = "Illegal leaf: %s/%s. If this is a variable size C array it's possible that the branch holding the size is not available.";
411 Error("TBranch", msg, name, leaflist);
412 delete [] leafname;
413 delete[] leaftype;
414 MakeZombie();
415 return;
416 }
417 leaf->SetBranch(this);
418 leaf->SetAddress((char*) (fAddress + offset));
419 leaf->SetOffset(offset);
420 if (leaf->GetLeafCount()) {
421 // -- Leaf is a varying length array, we need an offset array.
422 fEntryOffsetLen = 1000;
423 }
424 if (leaf->InheritsFrom(TLeafC::Class())) {
425 // -- Leaf is a character string, we need an offset array.
426 fEntryOffsetLen = 1000;
427 }
428 ++fNleaves;
430 fTree->GetListOfLeaves()->Add(leaf);
431 if (*pos == 0) {
432 // -- We reached the end of the leaf specification.
433 break;
434 }
435 nameBegin = pos + 1;
436 offset += leaf->GetLenType() * leaf->GetLen();
437 }
438 }
439 delete[] leafname;
440 leafname = nullptr;
441 delete[] leaftype;
442 leaftype = nullptr;
443
444}
445
446////////////////////////////////////////////////////////////////////////////////
447/// Destructor.
448
450{
451 delete fBrowsables;
452 fBrowsables = nullptr;
453
454 // Note: We do *not* have ownership of the buffer.
455 fEntryBuffer = nullptr;
456
457 delete [] fBasketSeek;
458 fBasketSeek = nullptr;
459
460 delete [] fBasketEntry;
461 fBasketEntry = nullptr;
462
463 delete [] fBasketBytes;
464 fBasketBytes = nullptr;
465
467 delete fExtraBasket;
469 fNBaskets = 0;
470 fCurrentBasket = nullptr;
472 fNextBasketEntry = -1;
473
474 // Remove our leaves from our tree's list of leaves.
475 if (fTree) {
477 if (lst && lst->GetLast()!=-1) {
478 lst->RemoveAll(&fLeaves);
479 }
480 }
481 // And delete our leaves.
482 fLeaves.Delete();
483
485
486 // If we are in a directory and that directory is not the same
487 // directory that our tree is in, then try to find an open file
488 // with the name fFileName. If we find one, delete that file.
489 // We are attempting to close any alternate file which we have
490 // been directed to write our baskets to.
491 // FIXME: We make no attempt to check if someone else might be
492 // using this file. This is very user hostile. A violation
493 // of the principle of least surprises.
494 //
495 // Warning. Must use FindObject by name instead of fDirectory->GetFile()
496 // because two branches may point to the same file and the file
497 // may have already been deleted in the previous branch.
498 if (fDirectory && (!fTree || fDirectory != fTree->GetDirectory())) {
500
502 TFile* file = (TFile*)gROOT->GetListOfFiles()->FindObject(bFileName);
503 if (file){
504 file->Close();
505 delete file;
506 file = nullptr;
507 }
508 }
509
510 fTree = nullptr;
511 fDirectory = nullptr;
512
513 if (fTransientBuffer) {
514 delete fTransientBuffer;
515 fTransientBuffer = nullptr;
516 }
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// Returns the transient buffer currently used by this TBranch for reading/writing baskets.
521
533
534////////////////////////////////////////////////////////////////////////////////
535/// Add the basket to this branch.
536///
537/// Warning: if the basket are not 'flushed/copied' in the same
538/// order as they were created, this will induce a slow down in
539/// the insert (since we'll need to move all the record that are
540/// entere 'too early').
541/// Warning we also assume that the __current__ write basket is
542/// not present (aka has been removed) or is empty (no entries).
543
545{
546 TBasket *basket = &b;
547
548 basket->SetBranch(this);
549
550 if (fWriteBasket >= fMaxBaskets) {
552 }
554
555 if (where && startEntry < fBasketEntry[where-1]) {
556 // Need to find the right location and move the possible baskets
557
558 if (!ondisk) {
559 Warning("AddBasket","The assumption that out-of-order basket only comes from disk based ntuple is false.");
560 }
561
562 if (startEntry < fBasketEntry[0]) {
563 where = 0;
564 } else {
565 for(Int_t i=fWriteBasket-1; i>=0; --i) {
566 if (fBasketEntry[i] < startEntry) {
567 where = i+1;
568 break;
569 } else if (fBasketEntry[i] == startEntry) {
570 Error("AddBasket","An out-of-order basket matches the entry number of an existing basket.");
571 }
572 }
573 }
574
575 if (where < fWriteBasket) {
576 // We shall move the content of the array
577 for (Int_t j=fWriteBasket; j > where; --j) {
581 }
582 }
583 }
585
587 if (existing && existing->GetNevBuf()) {
588 Fatal("AddBasket", "Dropping non-empty 'write' basket in %s %s",
589 GetTree()->GetName(), GetName());
590 }
591 delete existing;
592 if (ondisk) {
593 fBasketBytes[where] = basket->GetNbytes(); // not for in mem
594 fBasketSeek[where] = basket->GetSeekKey(); // not for in mem
596 ++fWriteBasket;
597 } else {
598 ++fNBaskets;
599 // The basket we are adding becomes the new 'write' basket.
601 fTree->IncrementTotalBuffers(basket->GetBufferSize());
602 }
603
604 fEntries += basket->GetNevBuf();
605 fEntryNumber += basket->GetNevBuf();
606 if (ondisk) {
607 fTotBytes += basket->GetObjlen() + basket->GetKeylen() ;
608 fZipBytes += basket->GetNbytes();
609 fTree->AddTotBytes(basket->GetObjlen() + basket->GetKeylen());
610 fTree->AddZipBytes(basket->GetNbytes());
611 }
612}
613
614////////////////////////////////////////////////////////////////////////////////
615/// Add the start entry of the write basket (not yet created)
616
618{
619 if (fWriteBasket >= fMaxBaskets) {
621 }
623
624 if (where && startEntry < fBasketEntry[where-1]) {
625 // Need to find the right location and move the possible baskets
626
627 Fatal("AddBasket","The last basket must have the highest entry number (%s/%lld/%d).",GetName(),startEntry,fWriteBasket);
628
629 }
630 // The first basket (should) always start at zero. If we are asked to update
631 // it, this likely to be from merging 'empty' branches (base class node and the likes)
632 if (where) {
635 }
636}
637
638////////////////////////////////////////////////////////////////////////////////
639/// Loop on all leaves of this branch to back fill Basket buffer.
640///
641/// Use this routine instead of TBranch::Fill when filling a branch individually
642/// to catch up with the number of entries already in the TTree.
643///
644/// First it calls TBranch::Fill and then if the number of entries of the branch
645/// reach one of TTree cluster's boundary, the basket is flushed.
646///
647/// The function returns the number of bytes committed to the memory basket.
648/// If a write error occurs, the number of bytes returned is -1.
649/// If no data are written, because e.g. the branch is disabled,
650/// the number of bytes returned is 0.
651///
652/// To insure that the baskets of each cluster are located close by in the
653/// file, when back-filling multiple branches make sure to call BackFill
654/// for the same entry for all the branches consecutively
655/// ~~~ {.cpp}
656/// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
657/// for( auto branch : branchCollection) {
658/// ... Make change to the data associated with the branch ...
659/// branch->BackFill();
660/// }
661/// }
662/// // Since we loop over all the branches for each new entry
663/// // all the baskets for a cluster are consecutive in the file.
664/// ~~~
665/// rather than doing all the entries of one branch at a time.
666/// ~~~ {.cpp}
667/// // Do NOT do things in the following order, it will lead to
668/// // poorly clustered files.
669/// for(auto branch : branchCollection) {
670/// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
671/// ... Make change to the data associated with the branch ...
672/// branch->BackFill();
673/// }
674/// }
675/// // Since we loop over all the entries for one branch
676/// // all the baskets for that branch are consecutive.
677/// ~~~
678
680
681 // Get the end of the next cluster.
683 cluster.Next();
684 auto endCluster = cluster.GetNextEntry();
685
686 auto result = FillImpl(nullptr);
687
688 if ( result && GetEntries() >= endCluster ) {
689 FlushBaskets();
690 }
691
692 return result;
693}
694
695////////////////////////////////////////////////////////////////////////////////
696/// Browser interface.
697
699{
700 if (fNleaves > 1) {
702 } else {
703 // Get the name and strip any extra brackets
704 // in order to get the full arrays.
705 TString name = GetName();
706 Int_t pos = name.First('[');
707 if (pos!=kNPOS) name.Remove(pos);
708
709 GetTree()->Draw(name, "", b ? b->GetDrawOption() : "");
710 if (gPad) gPad->Update();
711 }
712}
713
714 ///////////////////////////////////////////////////////////////////////////////
715 /// Loop on all branch baskets. If the file where branch buffers reside is
716 /// writable, free the disk space associated to the baskets of the branch,
717 /// then call Reset(). If the option contains "all", delete also the baskets
718 /// for the subbranches.
719 /// The branch is reset.
720 ///
721 /// NOTE that this function must be used with extreme care. Deleting branch baskets
722 /// fragments the file and may introduce inefficiencies when adding new entries
723 /// in the Tree or later on when reading the Tree.
724
726{
727 TString opt = option;
728 opt.ToLower();
729 TFile *file = GetFile(0);
730
732 for(Int_t i=0; i<fWriteBasket; i++) {
733 if (fBasketSeek[i]) file->MakeFree(fBasketSeek[i],fBasketSeek[i]+fBasketBytes[i]-1);
734 }
735 }
736
737 // process subbranches
738 if (opt.Contains("all")) {
740 Int_t nb = lb->GetEntriesFast();
741 for (Int_t j = 0; j < nb; j++) {
742 TBranch* branch = (TBranch*) lb->UncheckedAt(j);
743 if (branch) branch->DeleteBaskets("all");
744 }
745 }
746 DropBaskets("all");
747 Reset();
748}
749
750////////////////////////////////////////////////////////////////////////////////
751/// Loop on all branch baskets. Drop all baskets from memory except readbasket.
752/// If the option contains "all", drop all baskets including
753/// read- and write-baskets (unless they are not stored individually on disk).
754/// The option "all" also lead to DropBaskets being called on the sub-branches.
755
757{
758 bool all = false;
759 if (options && options[0]) {
760 TString opt = options;
761 opt.ToLower();
762 if (opt.Contains("all")) all = true;
763 }
764
767
768 if ( (fNBaskets>1) || all ) {
769 //slow case
770 for (Int_t i=0;i<nbaskets;i++) {
772 if (!basket) continue;
773 if ((i == fReadBasket || i == fWriteBasket) && !all) continue;
774 // if the basket is not yet on file but already has event in it
775 // we must continue to avoid dropping the basket (and thus losing data)
776 if (fBasketBytes[i]==0 && basket->GetNevBuf() > 0) continue;
777 basket->DropBuffers();
778 --fNBaskets;
780 if (basket == fCurrentBasket) {
781 fCurrentBasket = nullptr;
783 fNextBasketEntry = -1;
784 }
785 delete basket;
786 }
787
788 // process subbranches
789 if (all) {
791 Int_t nb = lb->GetEntriesFast();
792 for (Int_t j = 0; j < nb; j++) {
793 TBranch* branch = (TBranch*) lb->UncheckedAt(j);
794 if (!branch) continue;
795 branch->DropBaskets("all");
796 }
797 }
798 } else {
799 //fast case
800 if (nbaskets > 0) {
801 Int_t i = fBaskets.GetLast();
803 if (basket && fBasketBytes[i]!=0) {
804 basket->DropBuffers();
805 if (basket == fCurrentBasket) {
806 fCurrentBasket = nullptr;
808 fNextBasketEntry = -1;
809 }
810 delete basket;
811 fBaskets.AddAt(nullptr,i);
812 fBaskets.SetLast(-1);
813 fNBaskets = 0;
814 }
815 }
816 }
817
818}
819
820////////////////////////////////////////////////////////////////////////////////
821/// Increase BasketEntry buffer of a minimum of 10 locations
822/// and a maximum of 50 per cent of current size.
823
843
844////////////////////////////////////////////////////////////////////////////////
845/// Loop on all leaves of this branch to fill Basket buffer.
846///
847/// If TBranchIMTHelper is non-null and it is time to WriteBasket, then we will
848/// use TBB to compress in parallel.
849///
850/// The function returns the number of bytes committed to the memory basket.
851/// If a write error occurs, the number of bytes returned is -1.
852/// If no data are written, because e.g. the branch is disabled,
853/// the number of bytes returned is 0.
854
856{
857 if (TestBit(kDoNotProcess)) {
858 return 0;
859 }
860
862 if (!basket) {
863 basket = fTree->CreateBasket(this); // create a new basket
864 if (!basket) return 0;
865 ++fNBaskets;
867 }
868 TBuffer* buf = basket->GetBufferRef();
869
870 // Fill basket buffer.
871
872 Int_t nsize = 0;
873
874 if (buf->IsReading()) {
875 basket->SetWriteMode();
876 }
877
879 buf->ResetMap();
880 }
881
882 Int_t lnew = 0;
883 Int_t nbytes = 0;
884
885 if (fEntryBuffer) {
887 } else {
888 Int_t lold = buf->Length();
889 basket->Update(lold);
890 ++fEntries;
891 ++fEntryNumber;
892 (this->*fFillLeaves)(*buf);
893 if (buf->GetMapCount()) {
894 // The map is used.
896 }
897 lnew = buf->Length();
898 nbytes = lnew - lold;
899 }
900
901 if (fEntryOffsetLen) {
902 Int_t nevbuf = basket->GetNevBuf();
903 // Total size in bytes of EntryOffset table.
904 nsize = nevbuf * sizeof(Int_t);
905 } else {
906 if (!basket->GetNevBufSize()) {
907 basket->SetNevBufSize(nbytes);
908 }
909 }
910
911 // Should we create a new basket?
912 // fSkipZip force one entry per buffer (old stuff still maintained for CDF)
913 // Transfer full compressed buffer only
914
915 // If GetAutoFlush() is less than zero, then we are determining the end of the autocluster
916 // based upon the number of bytes already flushed. This is incompatible with one-basket-per-cluster
917 // (since we will grow the basket indefinitely and never flush!). Hence, we wait until the
918 // first event cluster is written out and *then* enable one-basket-per-cluster mode.
920
923 ((lnew + (2 * nsize) + nbytes) >= fBasketSize))) {
925 if (nout < 0) Error("TBranch::Fill", "Failed to write out basket.\n");
926 return (nout >= 0) ? nbytes : -1;
927 }
928 return nbytes;
929}
930
931////////////////////////////////////////////////////////////////////////////////
932/// Copy the data from fEntryBuffer into the current basket.
933
935{
936 Int_t nbytes = 0;
937 Int_t objectStart = 0;
938 Int_t last = 0;
939 Int_t lold = buf->Length();
940
941 // Handle the special case of fEntryBuffer != 0
942 if (fEntryBuffer->IsA() == TMessage::Class()) {
943 objectStart = 8;
944 }
946 // The buffer given as input has not been decompressed.
947 if (basket->GetNevBuf()) {
948 // If the basket already contains entry we need to close it
949 // out. (This is because we can only transfer full compressed
950 // buffer)
952 // And restart from scratch
953 return Fill();
954 }
957 static TBasket toread_fLast;
959 toread_fLast.Streamer(*fEntryBuffer);
961 last = toread_fLast.GetLast();
962 // last now contains the decompressed number of bytes.
964 buf->SetBufferOffset(0);
966 basket->Update(lold);
967 } else {
968 // We are required to copy starting at the version number (so not
969 // including the class name.
970 // See if byte count is here, if not it class still be a newClass
971 const UInt_t kNewClassTag = 0xFFFFFFFF;
972 const UInt_t kByteCountMask = 0x40000000; // OR the byte count with this
973 UInt_t tag = 0;
976 *fEntryBuffer >> tag;
977 if (tag & kByteCountMask) {
978 *fEntryBuffer >> tag;
979 }
980 if (tag == kNewClassTag) {
981 UInt_t maxsize = 256;
982 char* s = new char[maxsize];
984 fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end.
985 while (strlen(s) == (maxsize - 1)) {
986 // The classname is too large, try again with a large buffer.
988 maxsize *= 2;
989 delete[] s;
990 s = new char[maxsize];
991 fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end
992 }
993 delete[] s;
994 } else {
996 }
1000 }
1001 fEntries++;
1002 fEntryNumber++;
1003 UInt_t len = 0;
1005 if (startpos > UInt_t(objectStart)) {
1006 // We assume this buffer have just been directly filled
1007 // the current position in the buffer indicates the end of the object!
1009 } else {
1010 // The buffer have been acquired either via TSocket or via
1011 // TBuffer::SetBuffer(newloc,newsize)
1012 // Only the actual size of the memory buffer gives us an hint about where
1013 // the object ends.
1015 }
1018 // The original buffer came pre-compressed and thus the buffer Length
1019 // does not really show the really object size
1020 // lnew = nbytes = basket->GetLast();
1021 nbytes = last;
1022 lnew = last;
1023 } else {
1024 lnew = buf->Length();
1025 nbytes = lnew - lold;
1026 }
1027
1028 return nbytes;
1029}
1030
1031////////////////////////////////////////////////////////////////////////////////
1032/// Find the immediate sub-branch with passed name.
1033
1035{
1036 // We allow the user to pass only the last dotted component of the name.
1037 std::string longnm;
1038 longnm.reserve(fName.Length()+strlen(name)+3);
1039 longnm = fName.Data();
1040 if (longnm[longnm.length()-1]==']') {
1041 std::size_t dim = longnm.find_first_of('[');
1042 if (dim != std::string::npos) {
1043 longnm.erase(dim);
1044 }
1045 }
1046 if (longnm[longnm.length()-1] != '.') {
1047 longnm += '.';
1048 }
1049 longnm += name;
1050 UInt_t namelen = strlen(name);
1051
1053 TBranch* branch = nullptr;
1054 for(Int_t i = 0; i < nbranches; ++i) {
1056
1057 const char *brname = branch->fName.Data();
1058 UInt_t brlen = branch->fName.Length();
1059 if (brname[brlen-1]==']') {
1060 const char *dim = strchr(brname,'[');
1061 if (dim) {
1062 brlen = dim - brname;
1063 }
1064 }
1065 if (namelen == brlen /* same effective size */
1066 && strncmp(name,brname,brlen) == 0) {
1067 return branch;
1068 }
1069 if (brlen == (size_t)longnm.length()
1070 && strncmp(longnm.c_str(),brname,brlen) == 0) {
1071 return branch;
1072 }
1073 }
1074 return nullptr;
1075}
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// Find the leaf corresponding to the name 'searchname'.
1079
1081{
1086
1087 // We allow the user to pass only the last dotted component of the name.
1088 TIter next(GetListOfLeaves());
1089 TLeaf* leaf = nullptr;
1090 while ((leaf = (TLeaf*) next())) {
1091 leafname = leaf->GetName();
1092 Ssiz_t dim = leafname.First('[');
1093 if (dim >= 0) leafname.Remove(dim);
1094
1095 if (leafname == searchname) return leaf;
1096
1097 // The leaf element contains the branch name in its name, let's use the title.
1098 leaftitle = leaf->GetTitle();
1099 dim = leaftitle.First('[');
1100 if (dim >= 0) leaftitle.Remove(dim);
1101
1102 if (leaftitle == searchname) return leaf;
1103
1104 TBranch* branch = leaf->GetBranch();
1105 if (branch) {
1106 longname.Form("%s.%s",branch->GetName(),leafname.Data());
1107 dim = longname.First('[');
1108 if (dim>=0) longname.Remove(dim);
1109 if (longname == searchname) return leaf;
1110
1111 // The leaf element contains the branch name in its name.
1112 longname.Form("%s.%s",branch->GetName(),searchname);
1113 if (longname==leafname) return leaf;
1114
1115 longtitle.Form("%s.%s",branch->GetName(),leaftitle.Data());
1116 dim = longtitle.First('[');
1117 if (dim>=0) longtitle.Remove(dim);
1118 if (longtitle == searchname) return leaf;
1119
1120 // The following is for the case where the branch is only
1121 // a sub-branch. Since we do not see it through
1122 // TTree::GetListOfBranches, we need to see it indirectly.
1123 // This is the less sturdy part of this search ... it may
1124 // need refining ...
1125 if (strstr(searchname, ".") && !strcmp(searchname, branch->GetName())) return leaf;
1126 }
1127 }
1128 return nullptr;
1129}
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// Flush to disk all the baskets of this branch and any of subbranches.
1133/// Return the number of bytes written or -1 in case of write error.
1134
1136{
1137 UInt_t nerror = 0;
1138 Int_t nbytes = 0;
1139
1141 // The following protection is not necessary since we should always
1142 // have fWriteBasket < fBasket.GetSize()
1143 //if (fBaskets.GetSize() < maxbasket) {
1144 // maxbasket = fBaskets.GetSize();
1145 //}
1146 for(Int_t i=0; i != maxbasket; ++i) {
1147 if (fBaskets.UncheckedAt(i)) {
1149 if (nwrite<0) {
1150 ++nerror;
1151 } else {
1152 nbytes += nwrite;
1153 }
1154 }
1155 }
1157 for (Int_t i = 0; i < len; ++i) {
1159 if (!branch) {
1160 continue;
1161 }
1162 Int_t nwrite = branch->FlushBaskets();
1163 if (nwrite<0) {
1164 ++nerror;
1165 } else {
1166 nbytes += nwrite;
1167 }
1168 }
1169 if (nerror) {
1170 return -1;
1171 } else {
1172 return nbytes;
1173 }
1174}
1175
1176////////////////////////////////////////////////////////////////////////////////
1177/// If we have a write basket in memory and it contains some entries and
1178/// has not yet been written to disk, we write it and delete it from memory.
1179/// Return the number of bytes written;
1180
1182{
1183 Int_t nbytes = 0;
1186
1187 if (basket) {
1188 if (basket->GetNevBuf()
1189 && fBasketSeek[ibasket]==0) {
1190 // If the basket already contains entry we need to close it out.
1191 // (This is because we can only transfer full compressed buffer)
1192
1193 if (basket->GetBufferRef()->IsReading()) {
1194 basket->SetWriteMode();
1195 }
1197
1198 } else {
1199 // If the basket is empty or has already been written.
1200 if ((Int_t)ibasket==fWriteBasket) {
1201 // Nothing to do.
1202 } else {
1203 basket->DropBuffers();
1204 if (basket == fCurrentBasket) {
1205 fCurrentBasket = nullptr;
1206 fFirstBasketEntry = -1;
1207 fNextBasketEntry = -1;
1208 }
1209 delete basket;
1210 --fNBaskets;
1211 fBaskets[ibasket] = nullptr;
1212 }
1213 }
1214 }
1215 }
1216 return nbytes;
1217}
1218
1219////////////////////////////////////////////////////////////////////////////////
1220/// Return pointer to basket basketnumber in this Branch
1221///
1222/// If a new buffer must be created and the user_buffer argument is non-null,
1223/// then the memory in the user_buffer will be shared with the returned TBasket.
1224
1226{
1227 // This counter in the sequential case collects errors coming also from
1228 // different files (suppose to have a program reading f1.root, f2.root ...)
1229 // In the mt case, it is made atomic: it safely collects errors from
1230 // different files processed simultaneously.
1231 static std::atomic<Int_t> nerrors(0);
1232
1233 // reference to an existing basket in memory ?
1236 if (basket) return basket;
1237 if (basketnumber == fWriteBasket) return nullptr;
1238
1239 // create/decode basket parameters from buffer
1240 TFile *file = GetFile(0);
1241 if (file == nullptr) {
1242 return nullptr;
1243 }
1244 // if cluster pre-fetching or retaining is on, do not re-use existing baskets
1245 // unless a new cluster is used.
1248 else
1250
1251 // fSkipZip is old stuff still maintained for CDF
1253 if (fBasketBytes[basketnumber] == 0) {
1254 fBasketBytes[basketnumber] = basket->ReadBasketBytes(fBasketSeek[basketnumber],file);
1255 }
1256 //add branch to cache (if any)
1257 {
1258 R__LOCKGUARD_IMT(gROOTMutex); // Lock for parallel TTree I/O
1260 if (pf){
1261 if (pf->IsLearning()) pf->LearnBranch(this, false);
1262 if (fSkipZip) pf->SetSkipZip();
1263 }
1264 }
1265
1266 //now read basket
1268 if (R__unlikely(badread || basket->GetSeekKey() != fBasketSeek[basketnumber] || basket->IsZombie())) {
1269 nerrors++;
1270 if (nerrors > 10) return nullptr;
1271 if (nerrors == 10) {
1272 printf(" file probably overwritten: stopping reporting error messages\n");
1273 if (fBasketSeek[basketnumber] > 2000000000) {
1274 printf("===>File is more than 2 Gigabytes\n");
1275 return nullptr;
1276 }
1277 if (fBasketSeek[basketnumber] > 1000000000) {
1278 printf("===>Your file is may be bigger than the maximum file size allowed on your system\n");
1279 printf(" Check your AFS maximum file size limit for example\n");
1280 return nullptr;
1281 }
1282 }
1283 Error("GetBasket","File: %s at byte:%lld, branch:%s, entry:%lld, badread=%d, nerrors=%d, basketnumber=%d",file->GetName(),basket->GetSeekKey(),GetName(),fReadEntry,badread,nerrors.load(),basketnumber);
1284 return nullptr;
1285 }
1286
1287 ++fNBaskets;
1288
1290 auto perfStats = GetTree()->GetPerfStats();
1291 if (perfStats)
1292 perfStats->SetUsed(this, basketnumber);
1293
1295 return basket;
1296}
1297
1298////////////////////////////////////////////////////////////////////////////////
1299/// Return address of basket in the file
1300
1306
1307////////////////////////////////////////////////////////////////////////////////
1308/// Returns (and, if 0, creates) browsable objects for this branch
1309/// See TVirtualBranchBrowsable::FillListOfBrowsables.
1310
1317
1318////////////////////////////////////////////////////////////////////////////////
1319/// Return the name of the user class whose content is stored in this branch,
1320/// if any. If this branch was created using the 'leaflist' technique, this
1321/// function returns an empty string.
1322
1323const char * TBranch::GetClassName() const
1324{
1325 return "";
1326}
1327
1328////////////////////////////////////////////////////////////////////////////////
1329/// Return icon name depending on type of branch.
1330
1331const char* TBranch::GetIconName() const
1332{
1333 if (IsFolder())
1334 return "TBranchElement-folder";
1335 else
1336 return "TBranchElement-leaf";
1337}
1338
1339////////////////////////////////////////////////////////////////////////////////
1340/// A helper function to locate the correct basket - and its first entry.
1341/// Extracted to a common private function because it is needed by both GetEntry
1342/// and GetBulkEntries. It should not be called directly.
1343///
1344/// If a new basket must be constructed and the user_buffer is provided, then
1345/// the user_buffer will back the memory of the newly-constructed basket.
1346///
1347/// Assumes that this branch is enabled.
1348///
1349/// Returns -1 if the entry does not exist
1350/// Returns -2 in case of error
1351/// Returns the index of the basket in case of success.
1354{
1358 // We have found the basket containing this entry.
1359 // make sure basket buffers are in memory.
1361 first = fFirstBasketEntry;
1362 return fReadBasket;
1363 } else {
1364 if ((entry < fFirstEntry) || (entry >= fEntryNumber)) {
1365 return -1;
1366 }
1367 first = fFirstBasketEntry;
1368 Long64_t last = fNextBasketEntry - 1;
1369 // Are we still in the same ReadBasket?
1370 if ((entry < first) || (entry > last)) {
1372 if (fReadBasket < 0) {
1373 fNextBasketEntry = -1;
1374 Error("GetBasketAndFirst", "In the branch %s, no basket contains the entry %lld\n", GetName(), entry);
1375 return -2;
1376 }
1377 if (fReadBasket == fWriteBasket) {
1379 } else {
1381 }
1384 }
1385 // We have found the basket containing this entry.
1386 // make sure basket buffers are in memory.
1388 if (!basket) {
1390 if (!basket) {
1391 fCurrentBasket = nullptr;
1392 fFirstBasketEntry = -1;
1393 fNextBasketEntry = -1;
1394 return -2;
1395 }
1396 if (fTree->GetClusterPrefetch()) {
1398 clusterIterator.Next();
1399 Int_t nextClusterEntry = clusterIterator.GetNextEntry();
1400 for (Int_t i = fReadBasket + 1; i < fMaxBaskets && fBasketEntry[i] < nextClusterEntry; i++) {
1401 GetBasket(i);
1402 }
1403 }
1404 // Getting the next basket might reset the current one and
1405 // cause a reset of the first / next basket entries back to -1.
1406 fFirstBasketEntry = first;
1408 if (user_buffer) {
1409 // Disassociate basket from memory buffer for bulk IO
1410 // When the user provides a memory buffer (i.e., for bulk IO), we should
1411 // make sure to drop all references to that buffer in the TTree afterward.
1412 fCurrentBasket = nullptr;
1413 fBaskets[fReadBasket] = nullptr;
1414 } else {
1416 }
1417 } else {
1419 }
1420 return fReadBasket;
1421 }
1422}
1423
1424////////////////////////////////////////////////////////////////////////////////
1425/// Returns true if this branch supports bulk IO, false otherwise.
1426///
1427/// This will return true if all the various preconditions necessary hold true
1428/// to perform bulk IO (reasonable type, single TLeaf, etc); the bulk IO may
1429/// still fail, depending on the contents of the individual TBaskets loaded.
1431 return (fNleaves == 1) &&
1432 (static_cast<TLeaf*>(fLeaves.UncheckedAt(0))->GetDeserializeType() != TLeaf::DeserializeType::kExternal);
1433}
1434
1435////////////////////////////////////////////////////////////////////////////////
1436/// \brief Read a basket of events into the given buffer with byte swapping.
1437///
1438/// \return On success, the number of events of the type held by this branch
1439/// that have been read into the buffer. -1 on failure.
1440///
1441/// On success, the caller should be able to access the contents of buf as they
1442/// are with:
1443///
1444/// ~~~{.cpp}
1445/// static_cast<T*>(buf.GetCurrent())
1446/// ~~~
1447///
1448/// where T is the type stored on this branch.
1449///
1450/// When `count_buf` points to a valid TBuffer and the branch has a branch count,
1451/// `count_buf` will be filled (via a call to GetEntriesSerialized) with the data
1452/// from the branchCount. After deserialization those value can be used to calculate
1453/// the number of elements corresponding to each entries.
1454///
1455/// For each entry the number of elements is the multiplication of
1456///
1457/// ~~~{.cpp}
1458/// TLeaf *leaf = static_cast<TLeaf*>(branch->GetListOfLeaves()->At(0));
1459/// auto len = leaf->GetLen();
1460/// ~~~
1461///
1462/// and the value in the BranchCount corresponding to that entry (can be obtained
1463/// from `branch->GetBranchCount()`).
1464///
1465/// \note This interface is not meant to be exposed to end users, but rather it should
1466/// be wrapped by higher-level interfaces.
1467///
1468/// \note See TBranch::GetEntriesSerialized() for an alternative that does not
1469/// perform byte swapping (useful to save one pass over data in some cases).
1470///
1472{
1473 // TODO: eventually support multiple leaves.
1474 if (R__unlikely(fNleaves != 1)) return -1;
1475 TLeaf *leaf = static_cast<TLeaf*>(fLeaves.UncheckedAt(0));
1476 if (R__unlikely(leaf->GetDeserializeType() == TLeaf::DeserializeType::kExternal)) {
1477 return -1;
1478 }
1479
1480 // Remember which entry we are reading.
1481 fReadEntry = entry;
1482
1483 bool enabled = !TestBit(kDoNotProcess);
1484 if (R__unlikely(!enabled)) return -1;
1485 TBasket *basket = nullptr;
1486 Long64_t first;
1488 if (R__unlikely(result < 0)) return -1;
1489 // Only support reading from full clusters.
1490 if (R__unlikely(entry != first)) {
1491 //printf("Failed to read from full cluster; first entry is %ld; requested entry is %ld.\n", first, entry);
1492 return -1;
1493 }
1494
1495 basket->PrepareBasket(entry);
1496 TBuffer* buf = basket->GetBufferRef();
1497
1498 // Test for very old ROOT files.
1499 if (R__unlikely(!buf)) {
1500 Error("GetBulkEntries", "Failed to get a new buffer.\n");
1501 return -1;
1502 }
1503 // Test for displacements, which aren't supported in fast mode.
1504 if (R__unlikely(basket->GetDisplacement())) {
1505 Error("GetBulkEntries", "Basket has displacement.\n");
1506 return -1;
1507 }
1508
1509 if (&user_buf != buf) {
1510 // The basket was already in memory and might (and might not) be backed by persistent
1511 // storage.
1513 if (fBasketSeek[fReadBasket]) {
1514 // It is backed, so we can be destructive
1515 user_buf.SetBuffer(buf->Buffer(), buf->BufferSize());
1517 fCurrentBasket = nullptr;
1518 fBaskets[fReadBasket] = nullptr;
1519 } else {
1520 // This is the only copy, we can't return it as is to the user, just make a copy.
1521 if (user_buf.BufferSize() < buf->BufferSize()) {
1522 user_buf.AutoExpand(buf->BufferSize());
1523 }
1524 memcpy(user_buf.Buffer(), buf->Buffer(), buf->BufferSize());
1525 }
1526 }
1527
1528 Int_t bufbegin = basket->GetKeylen();
1529 user_buf.SetBufferOffset(bufbegin);
1530
1532 //printf("Requesting %d events; fNextBasketEntry=%lld; first=%lld.\n", N, fNextBasketEntry, first);
1533 if (R__unlikely(!leaf->ReadBasketFast(user_buf, N))) {
1534 Error("GetBulkEntries", "Leaf failed to read.\n");
1535 return -1;
1536 }
1537 user_buf.SetBufferOffset(bufbegin);
1538
1539 if (fCurrentBasket == nullptr) {
1540 R__ASSERT(fExtraBasket == nullptr && "fExtraBasket should have been set to nullptr by GetFreshBasket");
1542 basket->DisownBuffer();
1543 }
1544
1545 return N;
1546}
1547
1548////////////////////////////////////////////////////////////////////////////////
1549/// \brief Read a basket of events into the given buffer without byte swapping.
1550///
1551/// \return On success, the number of events of the type held by this branch
1552/// that have been read into the buffer. -1 on failure.
1553///
1554/// On success, the caller still need to deserialize the content. For example for
1555/// a scalar branch and `N` the return value (i.e. number of entries)
1556///
1557/// ~~~{.cpp}
1558/// rawdata = static_cast<char*>(buf.GetCurrent());
1559/// for (std::size_t i = 0u; i < N; ++i, ++target)
1560/// frombuf(rawdata, target); // `frombuf` also advances the `rawdata` pointer
1561/// ~~~
1562///
1563/// where target is a pointer or array to the type stored on this branch.
1564///
1565/// When `count_buf` points to a valid TBuffer and the branch has a branch count,
1566/// `count_buf` will be filled (via a call to GetEntriesSerialized()) with the data
1567/// from the branchCount. After deserialization those value can be used to calculate
1568/// the number of elements corresponding to each entries.
1569///
1570/// For each entry the number of elements is the multiplication of
1571///
1572/// ~~~{.cpp}
1573/// TLeaf *leaf = dynamic_cast<TLeaf*>(branch->GetListOfLeaves()->At(0));
1574/// auto len = leaf->GetLen();
1575/// ~~~
1576///
1577/// and the value in the BranchCount corresponding to that entry (can be obtained
1578/// from `branch->GetBranchCount()`).
1579///
1580/// \note This interface is not meant to be exposed to end users, but rather it should
1581/// be wrapped by higher-level interfaces.
1582///
1583/// \note See TBranch::GetBulkEntries() for an alternative that also performs byte swapping.
1584///
1586{
1587 // TODO: Template this and TBranch::GetBulkEntries; only difference is the TLeaf function (ReadBasketFast vs
1588 // ReadBasketSerialized
1589
1590 // TODO: eventually support multiple leaves.
1591 if (R__unlikely(fNleaves != 1)) { return -1; }
1592 TLeaf *leaf = static_cast<TLeaf*>(fLeaves.UncheckedAt(0));
1593 if (R__unlikely(leaf->GetDeserializeType() == TLeaf::DeserializeType::kDestructive)) {
1594 Error("GetEntriesSerialized", "Encountered a branch with destructive deserialization; failing.");
1595 return -1;
1596 }
1597
1598 // Remember which entry we are reading.
1599 fReadEntry = entry;
1600
1601 bool enabled = !TestBit(kDoNotProcess);
1602 if (R__unlikely(!enabled)) { return -1; }
1603 TBasket *basket = nullptr;
1604 Long64_t first;
1606 if (R__unlikely(result < 0)) { return -1; }
1607 // Only support reading from full clusters.
1608 if (R__unlikely(entry != first)) {
1609 Error("GetEntriesSerialized", "Failed to read from full cluster; first entry is %lld; requested entry is %lld.\n", first, entry);
1610 return -1;
1611 }
1612
1613 basket->PrepareBasket(entry);
1614 TBuffer* buf = basket->GetBufferRef();
1615
1616 // Test for very old ROOT files.
1617 if (R__unlikely(!buf)) {
1618 Error("GetEntriesSerialized", "Failed to get a new buffer.\n");
1619 return -1;
1620 }
1621 // Test for displacements, which aren't supported in fast mode.
1622 if (R__unlikely(basket->GetDisplacement())) {
1623 Error("GetEntriesSerialized", "Basket has displacement.\n");
1624 return -1;
1625 }
1626
1627 if (&user_buf != buf) {
1628 // The basket was already in memory and might (and might not) be backed by persistent
1629 // storage.
1631 if (fBasketSeek[fReadBasket]) {
1632 // It is backed, so we can be destructive
1633 user_buf.SetBuffer(buf->Buffer(), buf->BufferSize());
1635 fCurrentBasket = nullptr;
1636 fBaskets[fReadBasket] = nullptr;
1637 } else {
1638 // This is the only copy, we can't return it as is to the user, just make a copy.
1639 if (user_buf.BufferSize() < buf->BufferSize()) {
1640 user_buf.AutoExpand(buf->BufferSize());
1641 }
1642 memcpy(user_buf.Buffer(), buf->Buffer(), buf->BufferSize());
1643 }
1644 }
1645
1646 Int_t bufbegin = basket->GetKeylen();
1647 user_buf.SetBufferOffset(bufbegin);
1648
1650 //Info("GetEntriesSerialized", "Requesting %d events; fNextBasketEntry=%lld; first=%lld.\n", N, fNextBasketEntry, first);
1651
1652 user_buf.SetBufferOffset(bufbegin);
1653
1654 if (count_buf) {
1655 TLeaf *count_leaf = leaf->GetLeafCount();
1656 if (count_leaf) {
1657 //printf("Getting leaf count entries.\n");
1658 TBranch *count_branch = count_leaf->GetBranch();
1659 if (R__unlikely(count_branch->GetEntriesSerialized(entry, *count_buf) < 0)) {
1660 Error("GetEntriesSerialized", "Failed to read count leaf.\n");
1661 return -1;
1662 }
1663 } else {
1664 // TODO: if you ask for a count on a fixed-size branch, maybe we should
1665 // just fail?
1667 char *tmp_ptr = reinterpret_cast<char*>(&entry_count_serialized);
1668 tobuf(tmp_ptr, leaf->GetLenType() * leaf->GetNdata());
1669 Int_t cur_offset = count_buf->GetCurrent() - count_buf->Buffer();
1670 for (int idx=0; idx<N; idx++) {
1672 }
1673 count_buf->SetBufferOffset(cur_offset);
1674 }
1675 }
1676
1677 if (fCurrentBasket == nullptr) {
1678 R__ASSERT(fExtraBasket == nullptr && "fExtraBasket should have been set to nullptr by GetFreshBasket");
1680 basket->DisownBuffer();
1681 }
1682
1683 return N;
1684}
1685
1686////////////////////////////////////////////////////////////////////////////////
1687/// Read all leaves of entry and return total number of bytes read.
1688///
1689/// The input argument "entry" is the entry number in the current tree.
1690/// In case of a TChain, the entry number in the current Tree must be found
1691/// before calling this function. For example:
1692///
1693///~~~ {.cpp}
1694/// TChain* chain = ...;
1695/// Long64_t localEntry = chain->LoadTree(entry);
1696/// branch->GetEntry(localEntry);
1697///~~~
1698///
1699/// The function returns the number of bytes read from the input buffer.
1700/// If entry does not exist, the function returns 0.
1701/// If an I/O error occurs, the function returns -1.
1702///
1703/// See IMPORTANT REMARKS in TTree::GetEntry.
1704
1706{
1707 // Remember which entry we are reading.
1708 fReadEntry = entry;
1709
1710 if (R__unlikely(TestBit(kDoNotProcess) && !getall)) { return 0; }
1711
1712 TBasket *basket; // will be initialized in the if/then clauses.
1713 Long64_t first;
1714
1715 Int_t result = GetBasketAndFirst(basket, first, nullptr);
1716 if (R__unlikely(result < 0)) { return result + 1; }
1717
1718 basket->PrepareBasket(entry);
1719 TBuffer* buf = basket->GetBufferRef();
1720
1721 // This test necessary to read very old Root files (NvE).
1722 if (R__unlikely(!buf)) {
1723 TFile* file = GetFile(0);
1724 if (!file) return -1;
1725 basket->ReadBasketBuffers(fBasketSeek[fReadBasket], fBasketBytes[fReadBasket], file);
1726 buf = basket->GetBufferRef();
1727 }
1728
1729 // Set entry offset in buffer.
1731 buf->ResetMap();
1732 }
1733 if (R__unlikely(!buf->IsReading())) {
1734 basket->SetReadMode();
1735 }
1736
1737 Int_t* entryOffset = basket->GetEntryOffset();
1738 Int_t bufbegin = 0;
1739 if (entryOffset) {
1740 bufbegin = entryOffset[entry-first];
1742 Int_t* displacement = basket->GetDisplacement();
1745 }
1746 } else {
1747 bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1749 }
1750
1751 // Int_t bufbegin = buf->Length();
1752 (this->*fReadLeaves)(*buf);
1753 return buf->Length() - bufbegin;
1754}
1755
1756////////////////////////////////////////////////////////////////////////////////
1757/// Read all leaves of an entry and export buffers to real objects in a TClonesArray list.
1758///
1759/// Returns total number of bytes read.
1760
1762{
1763 // Remember which entry we are reading.
1764 fReadEntry = entry;
1765
1766 if (TestBit(kDoNotProcess)) {
1767 return 0;
1768 }
1769 if ((entry < 0) || (entry >= fEntryNumber)) {
1770 return 0;
1771 }
1772 Int_t nbytes = 0;
1774 Long64_t last = fNextBasketEntry - 1;
1775 // Are we still in the same ReadBasket?
1776 if ((entry < first) || (entry > last)) {
1778 if (fReadBasket < 0) {
1779 fNextBasketEntry = -1;
1780 Error("In the branch %s, no basket contains the entry %d\n", GetName(), entry);
1781 return -1;
1782 }
1783 if (fReadBasket == fWriteBasket) {
1785 } else {
1787 }
1789 }
1790
1791 // We have found the basket containing this entry.
1792 // Make sure basket buffers are in memory.
1795 if (!basket) {
1796 fFirstBasketEntry = -1;
1797 fNextBasketEntry = -1;
1798 return 0;
1799 }
1800 TBuffer* buf = basket->GetBufferRef();
1801 // Set entry offset in buffer and read data from all leaves.
1803 buf->ResetMap();
1804 }
1805 if (R__unlikely(!buf->IsReading())) {
1806 basket->SetReadMode();
1807 }
1808 Int_t* entryOffset = basket->GetEntryOffset();
1809 Int_t bufbegin = 0;
1810 if (entryOffset) {
1811 bufbegin = entryOffset[entry-first];
1813 Int_t* displacement = basket->GetDisplacement();
1816 }
1817 } else {
1818 bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1820 }
1822 leaf->ReadBasketExport(*buf, li, nentries);
1823 nbytes = buf->Length() - bufbegin;
1824 return nbytes;
1825}
1826
1827////////////////////////////////////////////////////////////////////////////////
1828/// Fill expectedClass and expectedType with information on the data type of the
1829/// object/values contained in this branch (and thus the type of pointers
1830/// expected to be passed to Set[Branch]Address
1831/// return 0 in case of success and > 0 in case of failure.
1832
1834{
1835 expectedClass = nullptr;
1837 TLeaf* l = (TLeaf*) GetListOfLeaves()->At(0);
1838 if (l) {
1839 expectedType = (EDataType) gROOT->GetType(l->GetTypeName())->GetType();
1840 return 0;
1841 } else {
1842 Error("GetExpectedType", "Did not find any leaves in %s",GetName());
1843 return 1;
1844 }
1845}
1846
1847////////////////////////////////////////////////////////////////////////////////
1848/// Return pointer to the file where branch buffers reside, returns 0
1849/// in case branch buffers reside in the same file as tree header.
1850/// If mode is 1 the branch buffer file is recreated.
1851
1853{
1854 if (fDirectory) return fDirectory->GetFile();
1855
1856 // check if a file with this name is in the list of Root files
1857 TFile *file = nullptr;
1858 {
1860 file = (TFile*)gROOT->GetListOfFiles()->FindObject(fFileName.Data());
1861 if (file) {
1862 fDirectory = file;
1863 return file;
1864 }
1865 }
1866
1867 if (fFileName.Length() == 0) return nullptr;
1868
1870
1871 // Open file (new file if mode = 1)
1872 {
1874 if (mode) file = TFile::Open(bFileName, "recreate");
1875 else file = TFile::Open(bFileName);
1876 }
1877 if (!file) return nullptr;
1878 if (file->IsZombie()) {delete file; return nullptr;}
1879 fDirectory = (TDirectory*)file;
1880 return file;
1881}
1882
1883////////////////////////////////////////////////////////////////////////////////
1884/// Return a fresh basket by either reusing an existing basket that needs
1885/// to be drop (according to TTree::MemoryFull) or create a new one.
1886///
1887/// If the user_buffer argument is non-null, then the memory in the
1888/// user-provided buffer will be utilized by the underlying basket.
1889///
1890/// The basket number is used to estimate the required buffer size
1891/// and try to optimize memory usage and number of memory allocation.
1892
1894{
1895 TBasket *basket = nullptr;
1896 if (user_buffer && fExtraBasket) {
1898 fExtraBasket = nullptr;
1899 basket->AdoptBuffer(user_buffer);
1900 } else {
1901 if (GetTree()->MemoryFull(0)) {
1902 if (fNBaskets==1) {
1903 // Steal the existing basket
1906 if (!basket) {
1907 fBaskets.SetLast(-2); // For recalculation of Last.
1909 if (oldindex != fBaskets.LowerBound()-1) {
1911 }
1912 }
1913 if (basket && fBasketBytes[oldindex]!=0) {
1914 if (basket == fCurrentBasket) {
1915 fCurrentBasket = nullptr;
1916 fFirstBasketEntry = -1;
1917 fNextBasketEntry = -1;
1918 }
1919 fBaskets.AddAt(nullptr,oldindex);
1920 fBaskets.SetLast(-1);
1921 fNBaskets = 0;
1922 basket->ReadResetBuffer(basketnumber);
1923#ifdef R__TRACK_BASKET_ALLOC_TIME
1924 fTree->AddAllocationTime(basket->GetResetAllocationTime());
1925#endif
1926 fTree->AddAllocationCount(basket->GetResetAllocationCount());
1927 } else {
1928 basket = fTree->CreateBasket(this);
1929 }
1930 } else if (fNBaskets == 0) {
1931 // There is nothing to drop!
1932 basket = fTree->CreateBasket(this);
1933 } else {
1934 // Memory is full and there is more than one basket,
1935 // Let DropBaskets do it job.
1936 DropBaskets();
1937 basket = fTree->CreateBasket(this);
1938 }
1939 } else {
1940 basket = fTree->CreateBasket(this);
1941 }
1942 if (user_buffer)
1943 basket->AdoptBuffer(user_buffer);
1944 }
1945 return basket;
1946}
1947
1948////////////////////////////////////////////////////////////////////////////////
1949/// Drops the cluster two behind the current cluster and returns a fresh basket
1950/// by either reusing or creating a new one
1951
1953{
1954 TBasket *basket = nullptr;
1955
1956 auto CreateOrReuseBasket = [this, user_buffer]() -> TBasket* {
1957 TBasket *newbasket = nullptr;
1958 if (fExtraBasket) {
1960 fExtraBasket = nullptr;
1961 } else {
1962 newbasket = fTree->CreateBasket(this);
1963 }
1964 if (user_buffer)
1965 newbasket->AdoptBuffer(user_buffer);
1966 return newbasket;
1967 };
1968
1969 // If GetClusterIterator is called with a negative entry then GetStartEntry will be 0
1970 // So we need to check if we reach the zero before we have gone back (1-VirtualSize) clusters
1971 // if this is the case, we want to keep everything in memory so we return a new basket
1973 if (iter.GetStartEntry() == 0) {
1974 return CreateOrReuseBasket();
1975 }
1976
1977 // Iterate backwards (1-VirtualSize) clusters to reach cluster to be unloaded from memory,
1978 // skipped if VirtualSize > 0.
1979 for (Int_t j = 0; j < -fTree->GetMaxVirtualSize(); j++) {
1980 if (iter.Previous() == 0) {
1981 return CreateOrReuseBasket();
1982 }
1983 }
1984
1985 Int_t entryToUnload = iter.Previous();
1986 // Finds the basket to unload from memory. Since the basket should be close to current
1987 // basket, just iterate backwards until the correct basket is reached. This should
1988 // be fast as long as the number of baskets per cluster is small
1992 if (basketToUnload < 0) {
1993 return CreateOrReuseBasket();
1994 }
1995 }
1996
1997 // Retrieves the basket that is going to be unloaded from memory. If the basket did not
1998 // exist, create a new one
2000 if (basket) {
2001 fBaskets.AddAt(nullptr, basketToUnload);
2002 --fNBaskets;
2003 } else {
2005 }
2007
2008 // Clear the rest of the baskets. While it would be ideal to reuse these baskets
2009 // for other baskets in the new cluster. It would require the function to go
2010 // beyond its current scope. In the ideal case when each cluster only has 1 basket
2011 // this will perform well
2012 iter.Next();
2013 while (fBasketEntry[basketToUnload] < iter.GetStartEntry()) {
2015 if (oldbasket) {
2016 oldbasket->DropBuffers();
2017 delete oldbasket;
2018 fBaskets.AddAt(nullptr, basketToUnload);
2019 --fNBaskets;
2020 }
2022 }
2023 fBaskets.SetLast(-1);
2024 return basket;
2025}
2026
2027////////////////////////////////////////////////////////////////////////////////
2028/// Return the 'full' name of the branch. In particular prefix the mother's name
2029/// when it does not end in a trailing dot and thus is not part of the branch name
2031{
2033 if (!mother || mother==this) {
2034 return fName;
2035 }
2036
2037 const auto motherName = mother->GetName();
2038 const auto len = strlen(motherName);
2039 if (len > 0 && (motherName[len-1] == '.')) {
2040 return fName;
2041 }
2042
2043 // Reserve the final size to avoid allocations
2044 TString result{static_cast<Ssiz_t>(len + 1 + fName.Length())};
2046 result += ".";
2047 result += fName;
2048 return result;
2049}
2050
2051////////////////////////////////////////////////////////////////////////////////
2052/// Return pointer to the 1st Leaf named name in thisBranch
2053
2054TLeaf* TBranch::GetLeaf(const char* name) const
2055{
2056 Int_t i;
2057 for (i=0;i<fNleaves;i++) {
2059 if (!strcmp(leaf->GetName(),name)) return leaf;
2060 }
2061 return nullptr;
2062}
2063
2064////////////////////////////////////////////////////////////////////////////////
2065/// Get real file name
2066
2068{
2069 if (fFileName.Length()==0) {
2070 return fFileName;
2071 }
2073
2074 // check if branch file name is absolute or a URL (e.g. root://host/...)
2075 char *bname = gSystem->ExpandPathName(fFileName.Data());
2076 if (!gSystem->IsAbsoluteFileName(bname) && !strstr(bname, ":/") && fTree && fTree->GetCurrentFile()) {
2077
2078 // if not, get filename where tree header is stored
2079 const char *tfn = fTree->GetCurrentFile()->GetName();
2080
2081 // If it is an archive file we need a special treatment
2082 TUrl arc(tfn);
2083 if (strlen(arc.GetAnchor()) > 0) {
2084 arc.SetAnchor(gSystem->BaseName(fFileName));
2085 bFileName = arc.GetUrl();
2086 } else {
2087 // if this is an absolute path or a URL then prepend this path
2088 // to the branch file name
2089 char *tname = gSystem->ExpandPathName(tfn);
2090 if (gSystem->IsAbsoluteFileName(tname) || strstr(tname, ":/")) {
2092 bFileName += "/";
2094 }
2095 delete [] tname;
2096 }
2097 }
2098 delete [] bname;
2099
2100 return bFileName;
2101}
2102
2103////////////////////////////////////////////////////////////////////////////////
2104/// Return all elements of one row unpacked in internal array fValues
2105/// [Actually just returns 1 (?)]
2106
2108{
2109 return 1;
2110}
2111
2112////////////////////////////////////////////////////////////////////////////////
2113/// Return whether this branch is in a mode where the object are decomposed
2114/// or not (Also known as MakeClass mode).
2115
2117{
2118 // Regular TBranch and TBrancObject can not be in makeClass mode
2119
2120 return false;
2121}
2122
2123////////////////////////////////////////////////////////////////////////////////
2124/// Get our top-level parent branch in the tree.
2125
2127{
2128 if (fMother) return fMother;
2129
2130 {
2131 TBranch *parent = fParent;
2132 while(parent) {
2133 if (parent->fMother) {
2134 const_cast<TBranch*>(this)->fMother = parent->fMother; // We can not yet use the 'mutable' keyword
2135 return fMother;
2136 }
2137 if (!parent->fParent) {
2138 // This is the top node
2139 const_cast<TBranch*>(this)->fMother = parent; // We can not yet use the 'mutable' keyword
2140 return fMother;
2141 }
2142 parent = parent->fParent;
2143 }
2144 }
2145
2146 const TObjArray* array = fTree->GetListOfBranches();
2147 Int_t n = array->GetEntriesFast();
2148 for (Int_t i = 0; i < n; ++i) {
2149 TBranch* branch = (TBranch*) array->UncheckedAt(i);
2150 TBranch* parent = branch->GetSubBranch(this);
2151 if (parent) {
2152 const_cast<TBranch*>(this)->fMother = branch; // We can not yet use the 'mutable' keyword
2153 return branch;
2154 }
2155 }
2156 return nullptr;
2157}
2158
2159////////////////////////////////////////////////////////////////////////////////
2160/// Find the parent branch of child.
2161/// Return 0 if child is not in this branch hierarchy.
2162
2164{
2165 // Handle error condition, if the parameter is us, we cannot find the parent.
2166 if (this == child) {
2167 // Note: We cast away any const-ness of "this".
2168 return (TBranch*) this;
2169 }
2170
2171 if (child->fParent) {
2172 return child->fParent;
2173 }
2174
2176 for (Int_t i = 0; i < len; ++i) {
2178 if (!branch) {
2179 continue;
2180 }
2181 if (branch == child) {
2182 // We are the direct parent of child.
2183 // Note: We cast away any const-ness of "this".
2184 const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
2185 return (TBranch*) this;
2186 }
2187 // FIXME: This is a tail-recursion!
2188 TBranch* parent = branch->GetSubBranch(child);
2189 if (parent) {
2190 return parent;
2191 }
2192 }
2193 // We failed to find the parent.
2194 return nullptr;
2195}
2196
2197////////////////////////////////////////////////////////////////////////////////
2198/// Return total number of bytes in the branch (including current buffer)
2199
2201{
2203 // This intentionally only store the TBranch part and thus slightly
2204 // under-estimate the space used.
2205 // Since the TBranchElement part contains pointers to other branches (branch count),
2206 // doing regular Streaming would end up including those and thus greatly over-estimate
2207 // the size used.
2208 const_cast<TBranch *>(this)->TBranch::Streamer(b);
2209
2210 Long64_t totbytes = 0;
2211 if (fZipBytes > 0) totbytes = fTotBytes;
2212 return totbytes + b.Length();
2213}
2214
2215////////////////////////////////////////////////////////////////////////////////
2216/// Return total number of bytes in the branch (excluding current buffer)
2217/// if option ="*" includes all sub-branches of this branch too
2218
2220{
2222 if (!option) return totbytes;
2223 if (option[0] != '*') return totbytes;
2224 //scan sub-branches
2226 for (Int_t i = 0; i < len; ++i) {
2228 if (branch) totbytes += branch->GetTotBytes(option);
2229 }
2230 return totbytes;
2231}
2232
2233////////////////////////////////////////////////////////////////////////////////
2234/// Return total number of zip bytes in the branch
2235/// if option ="*" includes all sub-branches of this branch too
2236
2238{
2240 if (!option) return zipbytes;
2241 if (option[0] != '*') return zipbytes;
2242 //scan sub-branches
2244 for (Int_t i = 0; i < len; ++i) {
2246 if (branch) zipbytes += branch->GetZipBytes(option);
2247 }
2248 return zipbytes;
2249}
2250
2251////////////////////////////////////////////////////////////////////////////////
2252/// Returns the IO settings currently in use for this branch.
2253
2258
2259////////////////////////////////////////////////////////////////////////////////
2260/// Return true if an existing object in a TBranchObject must be deleted.
2261
2263{
2264 return TestBit(kAutoDelete);
2265}
2266
2267////////////////////////////////////////////////////////////////////////////////
2268/// Return true if more than one leaf or browsables, false otherwise.
2269
2271{
2272 if (fNleaves > 1) {
2273 return true;
2274 }
2275 TList* browsables = const_cast<TBranch*>(this)->GetBrowsables();
2276 return browsables && browsables->GetSize();
2277}
2278
2279////////////////////////////////////////////////////////////////////////////////
2280/// keep a maximum of fMaxEntries in memory
2281
2283{
2286 if (basket) basket->MoveEntries(dentries);
2289 //loop on sub branches
2291 for (Int_t i = 0; i < nb; ++i) {
2293 branch->KeepCircular(maxEntries);
2294 }
2295}
2296
2297////////////////////////////////////////////////////////////////////////////////
2298/// Baskets associated to this branch are forced to be in memory.
2299/// You can call TTree::SetMaxVirtualSize(maxmemory) to instruct
2300/// the system that the total size of the imported baskets does not
2301/// exceed maxmemory bytes.
2302///
2303/// The function returns the number of baskets that have been put in memory.
2304/// This method may be called to force all baskets of one or more branches
2305/// in memory when random access to entries in this branch is required.
2306/// See also TTree::LoadBaskets to load all baskets of all branches in memory.
2307
2309{
2310 Int_t nimported = 0;
2312 TFile *file = GetFile(0);
2313 if (!file) return 0;
2314 TBasket *basket;
2315 for (Int_t i=0;i<nbaskets;i++) {
2317 if (basket) continue;
2318 basket = GetFreshBasket(i, nullptr);
2319 if (fBasketBytes[i] == 0) {
2320 fBasketBytes[i] = basket->ReadBasketBytes(fBasketSeek[i],file);
2321 }
2322 Int_t badread = basket->ReadBasketBuffers(fBasketSeek[i],fBasketBytes[i],file);
2323 if (badread) {
2324 Error("Loadbaskets","Error while reading basket buffer %d of branch %s",i,GetName());
2325 return -1;
2326 }
2327 ++fNBaskets;
2329 nimported++;
2330 }
2331 return nimported;
2332}
2333
2334////////////////////////////////////////////////////////////////////////////////
2335/// Print TBranch parameters
2336///
2337/// If options contains "basketsInfo" print the entry number, location and size
2338/// of each baskets.
2339
2341{
2342 const int kLINEND = 77;
2343 Float_t cx = 1;
2344
2346 if ( titleContent == GetName() ) {
2347 titleContent.Clear();
2348 }
2349
2350 if (fLeaves.GetEntries() == 1) {
2351 if (titleContent.Length()>=2 && titleContent[titleContent.Length()-2]=='/' && isalpha(titleContent[titleContent.Length()-1])) {
2352 // The type is already encoded. Nothing to do.
2353 } else {
2355 if (titleContent.Length()) {
2356 titleContent.Prepend(" ");
2357 }
2358 // titleContent.Append("type: ");
2359 titleContent.Prepend(leaf->GetTypeName());
2360 }
2361 }
2362 Int_t titleLength = titleContent.Length();
2363
2365 aLength += (aLength / 54 + 1) * 80 + 100;
2366 if (aLength < 200) aLength = 200;
2367 char *bline = new char[aLength];
2368
2370 if (fZipBytes) cx = (fTotBytes+0.00001)/fZipBytes;
2371 if (titleLength) snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName(),titleContent.Data());
2372 else snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName()," ");
2373 if (strlen(bline) > UInt_t(kLINEND)) {
2374 char *tmp = new char[strlen(bline)+1];
2375 if (titleLength) strlcpy(tmp, titleContent.Data(),strlen(bline)+1);
2376 snprintf(bline,aLength,"*Br%5d :%-9s : ",fgCount,GetName());
2377 int pos = strlen (bline);
2378 int npos = pos;
2379 int beg=0, end;
2380 while (beg < titleLength) {
2381 for (end=beg+1; end < titleLength-1; end ++)
2382 if (tmp[end] == ':') break;
2383 if (npos + end-beg+1 >= 78) {
2384 while (npos < kLINEND) {
2385 bline[pos ++] = ' ';
2386 npos ++;
2387 }
2388 bline[pos ++] = '*';
2389 bline[pos ++] = '\n';
2390 bline[pos ++] = '*';
2391 npos = 1;
2392 for (; npos < 12; npos ++)
2393 bline[pos ++] = ' ';
2394 bline[pos-2] = '|';
2395 }
2396 for (int n = beg; n <= end; n ++)
2397 bline[pos+n-beg] = tmp[n];
2398 pos += end-beg+1;
2399 npos += end-beg+1;
2400 beg = end+1;
2401 }
2402 while (npos < kLINEND) {
2403 bline[pos ++] = ' ';
2404 npos ++;
2405 }
2406 bline[pos ++] = '*';
2407 bline[pos] = '\0';
2408 delete[] tmp;
2409 }
2410 Printf("%s", bline);
2411
2412 if (fTotBytes > 2000000000) {
2413 Printf("*Entries :%lld : Total Size=%11lld bytes File Size = %lld *",fEntries,totBytes,fZipBytes);
2414 } else {
2415 if (fZipBytes > 0) {
2416 Printf("*Entries :%9lld : Total Size=%11lld bytes File Size = %10lld *",fEntries,totBytes,fZipBytes);
2417 } else {
2418 if (fWriteBasket > 0) {
2419 Printf("*Entries :%9lld : Total Size=%11lld bytes All baskets in memory *",fEntries,totBytes);
2420 } else {
2421 Printf("*Entries :%9lld : Total Size=%11lld bytes One basket in memory *",fEntries,totBytes);
2422 }
2423 }
2424 }
2425 Printf("*Baskets :%9d : Basket Size=%11d bytes Compression= %6.2f *",fWriteBasket,fBasketSize,cx);
2426
2427 if (strncmp(option,"basketsInfo",std::char_traits<char>::length("basketsInfo"))==0) {
2429 for (Int_t i=0;i<nbaskets;i++) {
2430 Printf("*Basket #%4d entry=%6lld pos=%6lld size=%5d",
2431 i, fBasketEntry[i], fBasketSeek[i], fBasketBytes[i]);
2432 }
2433 }
2434
2435 Printf("*............................................................................*");
2436 delete [] bline;
2437 fgCount++;
2438}
2439
2440////////////////////////////////////////////////////////////////////////////////
2441/// Print the information we have about which basket is currently cached and
2442/// whether they have been 'used'/'read' from the cache.
2443
2448
2449////////////////////////////////////////////////////////////////////////////////
2450/// Loop on all leaves of this branch to read Basket buffer.
2451
2453{
2454 // fLeaves->ReadBasket(basket);
2455}
2456
2457////////////////////////////////////////////////////////////////////////////////
2458/// Loop on all leaves of this branch to read Basket buffer.
2459
2461{
2462 for (Int_t i = 0; i < fNleaves; ++i) {
2464 leaf->ReadBasket(b);
2465 }
2466}
2467
2468////////////////////////////////////////////////////////////////////////////////
2469/// Read zero leaves without the overhead of a loop.
2470
2474
2475////////////////////////////////////////////////////////////////////////////////
2476/// Read one leaf without the overhead of a loop.
2477
2482
2483////////////////////////////////////////////////////////////////////////////////
2484/// Read two leaves without the overhead of a loop.
2485
2491
2492////////////////////////////////////////////////////////////////////////////////
2493/// Loop on all leaves of this branch to fill Basket buffer.
2494
2496{
2497 for (Int_t i = 0; i < fNleaves; ++i) {
2499 leaf->FillBasket(b);
2500 }
2501}
2502
2503////////////////////////////////////////////////////////////////////////////////
2504/// Refresh this branch using new information in b
2505/// This function is called by TTree::Refresh
2506
2508{
2509 if (b==nullptr) return;
2510
2511 fEntryOffsetLen = b->fEntryOffsetLen;
2512 fWriteBasket = b->fWriteBasket;
2513 fEntryNumber = b->fEntryNumber;
2514 fMaxBaskets = b->fMaxBaskets;
2515 fEntries = b->fEntries;
2516 fTotBytes = b->fTotBytes;
2517 fZipBytes = b->fZipBytes;
2518 fReadBasket = 0;
2519 fReadEntry = -1;
2520 fFirstBasketEntry = -1;
2521 fNextBasketEntry = -1;
2522 fCurrentBasket = nullptr;
2523 delete [] fBasketBytes;
2524 delete [] fBasketEntry;
2525 delete [] fBasketSeek;
2529 Int_t i;
2530 for (i=0;i<fMaxBaskets;i++) {
2531 fBasketBytes[i] = b->fBasketBytes[i];
2532 fBasketEntry[i] = b->fBasketEntry[i];
2533 fBasketSeek[i] = b->fBasketSeek[i];
2534 }
2535 fBaskets.Delete();
2536 Int_t nbaskets = b->fBaskets.GetSize();
2538 // If the current fWritebasket is in memory, take it (just swap)
2539 // from the Tree being read
2540 TBasket *basket = (TBasket*)b->fBaskets.UncheckedAt(fWriteBasket);
2542 if (basket) {
2543 fNBaskets = 1;
2544 --(b->fNBaskets);
2545 b->fBaskets.RemoveAt(fWriteBasket);
2546 basket->SetBranch(this);
2547 }
2548}
2549
2550////////////////////////////////////////////////////////////////////////////////
2551/// Reset a Branch.
2552///
2553/// - Existing buffers are deleted.
2554/// - Entries, max and min are reset.
2555
2557{
2558 fReadBasket = 0;
2559 fReadEntry = -1;
2560 fFirstBasketEntry = -1;
2561 fNextBasketEntry = -1;
2562 fCurrentBasket = nullptr;
2563 fWriteBasket = 0;
2564 fEntries = 0;
2565 fTotBytes = 0;
2566 fZipBytes = 0;
2567 fEntryNumber = 0;
2568
2569 if (fBasketBytes) {
2570 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2571 fBasketBytes[i] = 0;
2572 }
2573 }
2574
2575 if (fBasketEntry) {
2576 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2577 fBasketEntry[i] = 0;
2578 }
2579 }
2580
2581 if (fBasketSeek) {
2582 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2583 fBasketSeek[i] = 0;
2584 }
2585 }
2586
2587 fBaskets.Delete();
2588 fNBaskets = 0;
2589}
2590
2591////////////////////////////////////////////////////////////////////////////////
2592/// Reset a Branch.
2593///
2594/// - Existing buffers are deleted.
2595/// - Entries, max and min are reset.
2596
2598{
2599 fReadBasket = 0;
2600 fReadEntry = -1;
2601 fFirstBasketEntry = -1;
2602 fNextBasketEntry = -1;
2603 fCurrentBasket = nullptr;
2604 fWriteBasket = 0;
2605 fEntries = 0;
2606 fTotBytes = 0;
2607 fZipBytes = 0;
2608 fEntryNumber = 0;
2609
2610 if (fBasketBytes) {
2611 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2612 fBasketBytes[i] = 0;
2613 }
2614 }
2615
2616 if (fBasketEntry) {
2617 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2618 fBasketEntry[i] = 0;
2619 }
2620 }
2621
2622 if (fBasketSeek) {
2623 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2624 fBasketSeek[i] = 0;
2625 }
2626 }
2627
2629 if (reusebasket) {
2630 fBaskets[fWriteBasket] = nullptr;
2631 } else {
2633 if (reusebasket) {
2634 fBaskets[fReadBasket] = nullptr;
2635 }
2636 }
2637 fBaskets.Delete();
2638 if (reusebasket) {
2639 fNBaskets = 1;
2640 reusebasket->WriteReset();
2641 fBaskets[0] = reusebasket;
2642 } else {
2643 fNBaskets = 0;
2644 }
2645}
2646
2647////////////////////////////////////////////////////////////////////////////////
2648/// Reset the address of the branch.
2649
2651{
2652 fAddress = nullptr;
2653
2654 // Reset last read entry number, we have will had new user object now.
2655 fReadEntry = -1;
2656
2657 for (Int_t i = 0; i < fNleaves; ++i) {
2659 leaf->SetAddress(nullptr);
2660 }
2661
2663 for (Int_t i = 0; i < nbranches; ++i) {
2665 // FIXME: This is a tail recursion.
2666 abranch->ResetAddress();
2667 }
2668}
2669
2670////////////////////////////////////////////////////////////////////////////////
2671/// Static function resetting fgCount
2672
2674{
2675 fgCount = 0;
2676}
2677
2678////////////////////////////////////////////////////////////////////////////////
2679/// Set address of this branch.
2680/// @see TBranchElement::SetAddress
2681/// @note TBranch::SetAddress is a lower level interface and has less ability
2682/// to check for incorrect setup than TTree::SetBranchAddress. Without
2683/// TTree::SetMakeClass, if the branch is within an object, the input of
2684/// SetAddress is expected to be the start of the object (and thus the offset
2685/// of the data member is added to the provided address). The explicit purpose
2686/// of TTree::SetMakeClass is to disable this addition of the offset. Note that
2687/// TTree::SetBranchAddress will detect this case and automatically call
2688/// SetMakeClass for a data member within a class.
2689/// For example, the tutorial https://root.cern/doc/master/tree108__tree_8C.html
2690/// generates a ROOT file with a TTree. To read the temperature values,
2691/// you need either `tree->SetBranchAddress("fTemperature", &temp);` or
2692/// `tree->SetMakeClass(1); tree->GetBranch("fTemperature")->SetAddress(&temp);`.
2693
2695{
2696 if (TestBit(kDoNotProcess)) {
2697 return;
2698 }
2699 fReadEntry = -1;
2700 fFirstBasketEntry = -1;
2701 fNextBasketEntry = -1;
2702 fAddress = (char*) addr;
2703 for (Int_t i = 0; i < fNleaves; ++i) {
2705 Int_t offset = leaf->GetOffset();
2706 if (TestBit(kIsClone)) {
2707 offset = 0;
2708 }
2709 if (fAddress) leaf->SetAddress(fAddress + offset);
2710 else leaf->SetAddress(nullptr);
2711 }
2712}
2713
2714////////////////////////////////////////////////////////////////////////////////
2715/// Set the automatic delete bit.
2716///
2717/// This bit is used by TBranchObject::ReadBasket to decide if an object
2718/// referenced by a TBranchObject must be deleted or not before reading
2719/// a new entry.
2720///
2721/// If autodel is true, this existing object will be deleted, a new object
2722/// created by the default constructor, then read from disk by the streamer.
2723///
2724/// If autodel is false, the existing object is not deleted. Root assumes
2725/// that the user is taking care of deleting any internal object or array
2726/// (this can be done in the streamer).
2727
2729{
2730 if (autodel) {
2731 SetBit(kAutoDelete, true);
2732 } else {
2733 SetBit(kAutoDelete, false);
2734 }
2735}
2736
2737////////////////////////////////////////////////////////////////////////////////
2738/// Set the basket size
2739/// The function makes sure that the basket size is greater than fEntryOffsetlen
2740
2751
2752////////////////////////////////////////////////////////////////////////////////
2753/// Set address of this branch directly from a TBuffer to avoid streaming.
2754///
2755/// Note: We do not take ownership of the buffer.
2756
2758{
2759 // Check this is possible
2760 if ( (fNleaves != 1)
2761 || (strcmp("TLeafObject",fLeaves.UncheckedAt(0)->ClassName())!=0) ) {
2762 Error("TBranch::SetAddress","Filling from a TBuffer can only be done with a not split object branch. Request ignored.");
2763 } else {
2764 fReadEntry = -1;
2765 fNextBasketEntry = -1;
2766 fFirstBasketEntry = -1;
2767 // Note: We do not take ownership of the buffer.
2768 fEntryBuffer = buf;
2769 }
2770}
2771
2772////////////////////////////////////////////////////////////////////////////////
2773/// Set compression algorithm.
2774
2776{
2778 if (fCompress < 0) {
2780 } else {
2781 int level = fCompress % 100;
2782 fCompress = 100 * algorithm + level;
2783 }
2784
2786 for (Int_t i=0;i<nb;i++) {
2788 branch->SetCompressionAlgorithm(algorithm);
2789 }
2790}
2791
2792////////////////////////////////////////////////////////////////////////////////
2793/// Set compression level.
2794
2796{
2797 if (level < 0) level = 0;
2798 if (level > 99) level = 99;
2799 if (fCompress < 0) {
2800 fCompress = level;
2801 } else {
2802 int algorithm = fCompress / 100;
2804 fCompress = 100 * algorithm + level;
2805 }
2806
2808 for (Int_t i=0;i<nb;i++) {
2810 branch->SetCompressionLevel(level);
2811 }
2812}
2813
2814////////////////////////////////////////////////////////////////////////////////
2815/// Set compression settings.
2816
2818{
2820
2822 for (Int_t i=0;i<nb;i++) {
2824 branch->SetCompressionSettings(settings);
2825 }
2826}
2827
2828////////////////////////////////////////////////////////////////////////////////
2829/// Update the default value for the branch's fEntryOffsetLen if and only if
2830/// it was already non zero (and the new value is not zero)
2831/// If updateExisting is true, also update all the existing branches.
2832
2834{
2835 if (fEntryOffsetLen && newdefault) {
2837 }
2838 if (updateExisting) {
2839 TIter next( GetListOfBranches() );
2840 TBranch *b;
2841 while ( ( b = (TBranch*)next() ) ) {
2842 b->SetEntryOffsetLen( newdefault, true );
2843 }
2844 }
2845}
2846
2847////////////////////////////////////////////////////////////////////////////////
2848/// Set the number of entries in this branch.
2849
2851{
2852 fEntries = entries;
2853 fEntryNumber = entries;
2854}
2855
2856////////////////////////////////////////////////////////////////////////////////
2857/// Set file where this branch writes/reads its buffers.
2858/// By default the branch buffers reside in the file where the
2859/// Tree was created.
2860/// If the file name where the tree was created is an absolute
2861/// path name or an URL (e.g. or root://host/...)
2862/// and if the fname is not an absolute path name or an URL then
2863/// the path of the tree file is prepended to fname to make the
2864/// branch file relative to the tree file. In this case one can
2865/// move the tree + all branch files to a different location in
2866/// the file system and still access the branch files.
2867/// The ROOT file will be connected only when necessary.
2868/// If called by TBranch::Fill (via TBasket::WriteFile), the file
2869/// will be created with the option "recreate".
2870/// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2871/// will be opened in read mode.
2872/// To open a file in "update" mode or with a certain compression
2873/// level, use TBranch::SetFile(TFile *file).
2874
2876{
2877 if (file == nullptr) file = fTree->GetCurrentFile();
2878 fDirectory = (TDirectory*)file;
2879 if (file == fTree->GetCurrentFile()) fFileName = "";
2880 else fFileName = file->GetName();
2881
2882 if (file && fCompress == -1) {
2884 }
2885
2886 // Apply to all existing baskets.
2888 TBasket *basket;
2889 while ((basket = (TBasket*)nextb())) {
2890 basket->SetParent(file);
2891 }
2892
2893 // Apply to sub-branches as well.
2894 TIter next(GetListOfBranches());
2895 TBranch *branch;
2896 while ((branch = (TBranch*)next())) {
2897 branch->SetFile(file);
2898 }
2899}
2900
2901////////////////////////////////////////////////////////////////////////////////
2902/// Set file where this branch writes/reads its buffers.
2903/// By default the branch buffers reside in the file where the
2904/// Tree was created.
2905/// If the file name where the tree was created is an absolute
2906/// path name or an URL (e.g. root://host/...)
2907/// and if the fname is not an absolute path name or an URL then
2908/// the path of the tree file is prepended to fname to make the
2909/// branch file relative to the tree file. In this case one can
2910/// move the tree + all branch files to a different location in
2911/// the file system and still access the branch files.
2912/// The ROOT file will be connected only when necessary.
2913/// If called by TBranch::Fill (via TBasket::WriteFile), the file
2914/// will be created with the option "recreate".
2915/// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2916/// will be opened in read mode.
2917/// To open a file in "update" mode or with a certain compression
2918/// level, use TBranch::SetFile(TFile *file).
2919
2920void TBranch::SetFile(const char* fname)
2921{
2922 fFileName = fname;
2923 fDirectory = nullptr;
2924
2925 //apply to sub-branches as well
2926 TIter next(GetListOfBranches());
2927 TBranch *branch;
2928 while ((branch = (TBranch*)next())) {
2929 branch->SetFile(fname);
2930 }
2931}
2932
2933////////////////////////////////////////////////////////////////////////////////
2934/// Set the branch in a mode where the object are decomposed
2935/// (Also known as MakeClass mode).
2936/// Return whether the setting was possible (it is not possible for
2937/// TBranch and TBranchObject).
2938
2939bool TBranch::SetMakeClass(bool /* decomposeObj */)
2940{
2941 // Regular TBranch and TBrancObject can not be in makeClass mode
2942 return false;
2943}
2944
2945////////////////////////////////////////////////////////////////////////////////
2946/// Set object this branch is pointing to.
2947
2948void TBranch::SetObject(void * /* obj */)
2949{
2950 if (TestBit(kDoNotProcess)) {
2951 return;
2952 }
2953 Warning("SetObject","is not supported in TBranch objects");
2954}
2955
2956////////////////////////////////////////////////////////////////////////////////
2957/// Set branch status to Process or DoNotProcess.
2958
2959void TBranch::SetStatus(bool status)
2960{
2961 if (status) ResetBit(kDoNotProcess);
2962 else SetBit(kDoNotProcess);
2963}
2964
2965////////////////////////////////////////////////////////////////////////////////
2966/// Stream a class object
2967
2969{
2970 if (b.IsReading()) {
2971 UInt_t R__s, R__c;
2972 fTree = nullptr; // Will be set by TTree::Streamer
2973 fAddress = nullptr;
2974 gROOT->SetReadingObject(true);
2975
2976 // Reset transients.
2978 fCurrentBasket = nullptr;
2979 fFirstBasketEntry = -1;
2980 fNextBasketEntry = -1;
2981
2982 Version_t v = b.ReadVersion(&R__s, &R__c);
2983 if (v > 9) {
2984 b.ReadClassBuffer(TBranch::Class(), this, v, R__s, R__c);
2985
2986 if (fWriteBasket>=fBaskets.GetSize()) {
2988 }
2989 fDirectory = nullptr;
2991 for (Int_t i=0;i<fNleaves;i++) {
2993 leaf->SetBranch(this);
2994 }
2996 for (Int_t i=0;i<nbranches;i++) {
2998 br->fParent = this;
2999 }
3000
3001 fNBaskets = 0;
3002 for (Int_t j = fWriteBasket; j>=0; --j) {
3004 if (bk) {
3005 bk->SetBranch(this);
3006 // GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
3007 ++fNBaskets;
3008 }
3009 }
3010 if (fWriteBasket >= fMaxBaskets) {
3011 //old versions may need this fix
3016
3017 }
3019 gROOT->SetReadingObject(false);
3020 if (IsA() == TBranch::Class()) {
3021 if (fNleaves == 0) {
3023 } else if (fNleaves == 1) {
3025 } else if (fNleaves == 2) {
3027 } else {
3029 }
3030 }
3031 return;
3032 }
3033 //====process old versions before automatic schema evolution
3034 Int_t n,i,j,ijunk;
3035 if (v > 5) {
3036 Stat_t djunk;
3038 if (v > 7) TAttFill::Streamer(b);
3039 b >> fCompress;
3040 b >> fBasketSize;
3041 b >> fEntryOffsetLen;
3042 b >> fWriteBasket;
3044 b >> fOffset;
3045 b >> fMaxBaskets;
3046 if (v > 6) b >> fSplitLevel;
3047 b >> djunk; fEntries = (Long64_t)djunk;
3050
3058 b >> isArray;
3059 b.ReadFastArray(fBasketBytes,fMaxBaskets);
3060 b >> isArray;
3061 for (i=0;i<fMaxBaskets;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
3062 b >> isArray;
3063 for (i=0;i<fMaxBaskets;i++) {
3064 if (isArray == 2) b >> fBasketSeek[i];
3065 else {Int_t bsize; b >> bsize; fBasketSeek[i] = (Long64_t)bsize;};
3066 }
3068 b.CheckByteCount(R__s, R__c, TBranch::IsA());
3069 fDirectory = nullptr;
3071 for (i=0;i<fNleaves;i++) {
3073 leaf->SetBranch(this);
3074 }
3075 fNBaskets = 0;
3076 for (j = fWriteBasket; j >= 0; --j) {
3078 if (bk) {
3079 bk->SetBranch(this);
3080 //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
3081 ++fNBaskets;
3082 }
3083 }
3084 if (fWriteBasket >= fMaxBaskets) {
3085 //old versions may need this fix
3090
3091 }
3092 // Check Byte Count is not needed since it was done in ReadBuffer
3094 gROOT->SetReadingObject(false);
3095 b.CheckByteCount(R__s, R__c, TBranch::IsA());
3096 if (IsA() == TBranch::Class()) {
3097 if (fNleaves == 0) {
3099 } else if (fNleaves == 1) {
3101 } else if (fNleaves == 2) {
3103 } else {
3105 }
3106 }
3107 return;
3108 }
3109 //====process very old versions
3110 Stat_t djunk;
3112 b >> fCompress;
3113 b >> fBasketSize;
3114 b >> fEntryOffsetLen;
3115 b >> fMaxBaskets;
3116 b >> fWriteBasket;
3118 b >> djunk; fEntries = (Long64_t)djunk;
3121 b >> fOffset;
3125 for (i=0;i<fNleaves;i++) {
3127 leaf->SetBranch(this);
3128 }
3130 for (j = fWriteBasket; j > 0; --j) {
3132 if (bk) {
3133 bk->SetBranch(this);
3134 //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
3135 }
3136 }
3138 b >> n;
3139 for (i=0;i<n;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
3141 if (v > 4) {
3142 n = b.ReadArray(fBasketBytes);
3143 } else {
3144 for (n=0;n<fMaxBaskets;n++) fBasketBytes[n] = 0;
3145 }
3146 if (v < 2) {
3148 for (n=0;n<fWriteBasket;n++) {
3149 TBasket *basket = GetBasketImpl(n, nullptr);
3150 fBasketSeek[n] = basket ? basket->GetSeekKey() : 0;
3151 }
3152 } else {
3154 b >> n;
3155 for (n=0;n<fMaxBaskets;n++) {
3156 Int_t aseek;
3157 b >> aseek;
3159 }
3160 }
3161 if (v > 2) {
3163 }
3164 fDirectory = nullptr;
3165 if (v < 4) SetAutoDelete(true);
3167 gROOT->SetReadingObject(false);
3168 b.CheckByteCount(R__s, R__c, TBranch::IsA());
3169 //====end of old versions
3170 if (IsA() == TBranch::Class()) {
3171 if (fNleaves == 0) {
3173 } else if (fNleaves == 1) {
3175 } else if (fNleaves == 2) {
3177 } else {
3179 }
3180 }
3181 } else {
3185 if (fMaxBaskets < 10) fMaxBaskets = 10;
3186
3187 TBasket **stash = new TBasket *[lastBasket];
3188 for (Int_t i = 0; i < lastBasket; ++i) {
3190 if (ba && (fBasketBytes[i] || ba->GetNevBuf()==0)) {
3191 // Already on disk or empty.
3192 stash[i] = ba;
3193 fBaskets[i] = nullptr;
3194 } else {
3195 stash[i] = nullptr;
3196 }
3197 }
3198
3199 b.WriteClassBuffer(TBranch::Class(), this);
3200
3201 for (Int_t i = 0; i < lastBasket; ++i) {
3202 if (stash[i]) fBaskets[i] = stash[i];
3203 }
3204
3205 delete[] stash;
3207 }
3208}
3209
3210////////////////////////////////////////////////////////////////////////////////
3211/// Write the current basket to disk and return the number of bytes
3212/// written to the file.
3213
3215{
3216 Int_t nevbuf = basket->GetNevBuf();
3217 if (fEntryOffsetLen > 10 && (4*nevbuf) < fEntryOffsetLen ) {
3218 // Make sure that the fEntryOffset array does not stay large unnecessarily.
3219 fEntryOffsetLen = nevbuf < 3 ? 10 : 4*nevbuf; // assume some fluctuations.
3220 } else if (fEntryOffsetLen && nevbuf > fEntryOffsetLen) {
3221 // Increase the array ...
3222 fEntryOffsetLen = 2*nevbuf; // assume some fluctuations.
3223 }
3224
3225 // Note: captures `basket`, `where`, and `this` by value; modifies the TBranch and basket,
3226 // as we make a copy of the pointer. We cannot capture `basket` by reference as the pointer
3227 // itself might be modified after `WriteBasketImpl` exits.
3228 auto doUpdates = [this, basket, where]() {
3229 Int_t nout = basket->WriteBuffer(); // Write buffer
3230 if (nout < 0)
3231 Error("WriteBasketImpl", "basket's WriteBuffer failed.");
3232 fBasketBytes[where] = basket->GetNbytes();
3233 fBasketSeek[where] = basket->GetSeekKey();
3234 Int_t addbytes = basket->GetObjlen() + basket->GetKeylen();
3235 TBasket *reusebasket = nullptr;
3236 if (nout>0) {
3237 // The Basket was written so we can now safely reuse it.
3238 fBaskets[where] = nullptr;
3239
3241 reusebasket->WriteReset();
3242
3243 fZipBytes += nout;
3247#ifdef R__TRACK_BASKET_ALLOC_TIME
3248 fTree->AddAllocationTime(reusebasket->GetResetAllocationTime());
3249#endif
3250 fTree->AddAllocationCount(reusebasket->GetResetAllocationCount());
3251 }
3252
3253 if (where==fWriteBasket) {
3254 ++fWriteBasket;
3255 if (fWriteBasket >= fMaxBaskets) {
3257 }
3259 // The 'current' basket has Reset, so if we need it we will need
3260 // to reload it.
3261 fCurrentBasket = nullptr;
3262 fFirstBasketEntry = -1;
3263 fNextBasketEntry = -1;
3264 }
3267 } else {
3268 --fNBaskets;
3269 fBaskets[where] = nullptr;
3270 basket->DropBuffers();
3271 if (basket == fCurrentBasket) {
3272 fCurrentBasket = nullptr;
3273 fFirstBasketEntry = -1;
3274 fNextBasketEntry = -1;
3275 }
3276 delete basket;
3277 }
3278 return nout;
3279 };
3280 if (imtHelper) {
3281 imtHelper->Run(doUpdates);
3282 return 0;
3283 } else {
3284 return doUpdates();
3285 }
3286}
3287
3288////////////////////////////////////////////////////////////////////////////////
3289///set the first entry number (case of TBranchSTL)
3290
3292{
3294 fEntries = 0;
3296 if( fBasketEntry )
3297 fBasketEntry[0] = entry;
3298 for( Int_t i = 0; i < fBranches.GetEntriesFast(); ++i )
3299 ((TBranch*)fBranches[i])->SetFirstEntry( entry );
3300}
3301
3302////////////////////////////////////////////////////////////////////////////////
3303/// If the branch address is not set, we set all addresses starting with
3304/// the top level parent branch.
3305
3307{
3308 SetAddress(nullptr); // in some cases, this triggers setting of the address
3309}
3310
3311////////////////////////////////////////////////////////////////////////////////
3312/// Refresh the value of fDirectory (i.e. where this branch writes/reads its buffers)
3313/// with the current value of fTree->GetCurrentFile unless this branch has been
3314/// redirected to a different file. Also update the sub-branches.
3315
3317{
3318 TFile *file = fTree->GetCurrentFile();
3319 if (fFileName.Length() == 0) {
3320 fDirectory = file;
3321
3322 // Apply to all existing baskets.
3324 TBasket *basket;
3325 while ((basket = (TBasket*)nextb())) {
3326 basket->SetParent(file);
3327 }
3328 }
3329
3330 // Apply to sub-branches as well.
3331 TIter next(GetListOfBranches());
3332 TBranch *branch;
3333 while ((branch = (TBranch*)next())) {
3334 branch->UpdateFile();
3335 }
3336}
void tobuf(char *&buf, Bool_t x)
Definition Bytes.h:55
#define R__likely(expr)
Definition RConfig.hxx:595
#define R__unlikely(expr)
Definition RConfig.hxx:594
#define b(i)
Definition RSha256.hxx:100
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
short Version_t
Class version identifier (short)
Definition RtypesCore.h:79
char Char_t
Character 1 byte (char)
Definition RtypesCore.h:51
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int)
Definition RtypesCore.h:60
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Ssiz_t kNPOS
The equivalent of std::string::npos for the ROOT class TString.
Definition RtypesCore.h:131
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
const UInt_t kNewClassTag
const UInt_t kByteCountMask
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
EDataType
Definition TDataType.h:28
@ kOther_t
Definition TDataType.h:32
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
#define N
Option_t Option_t option
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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
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 Int_t Int_t Window_t child
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Option_t Option_t TPoint TPoint const char mode
char name[80]
Definition TGX11.cxx:110
int nentries
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define gROOT
Definition TROOT.h:411
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
Definition TString.cxx:2509
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define R__LOCKGUARD_IMT(mutex)
#define R__LOCKGUARD(mutex)
#define gPad
#define snprintf
Definition civetweb.c:1579
void SetUsed(Int_t basketNumber)
Mark if the basket has been marked as 'used'.
void Print(const char *owner, Long64_t *entries) const
Print the info we have for the baskets.
A helper class for managing IMT work during TTree:Fill operations.
TIOFeatures provides the end-user with the ability to change the IO behavior of data written via a TT...
Fill Area Attributes class.
Definition TAttFill.h:20
virtual void Streamer(TBuffer &)
Manages buffers for branches of a Tree.
Definition TBasket.h:34
A TTree is a list of TBranches.
Definition TBranch.h:93
virtual TLeaf * GetLeaf(const char *name) const
Return pointer to the 1st Leaf named name in thisBranch.
Definition TBranch.cxx:2054
virtual bool GetMakeClass() const
Return whether this branch is in a mode where the object are decomposed or not (Also known as MakeCla...
Definition TBranch.cxx:2116
virtual void SetupAddresses()
If the branch address is not set, we set all addresses starting with the top level parent branch.
Definition TBranch.cxx:3306
virtual void ResetAddress()
Reset the address of the branch.
Definition TBranch.cxx:2650
TString fFileName
Name of file where buffers are stored ("" if in same file as Tree header)
Definition TBranch.h:149
virtual const char * GetClassName() const
Return the name of the user class whose content is stored in this branch, if any.
Definition TBranch.cxx:1323
TBasket * GetFreshBasket(Int_t basketnumber, TBuffer *user_buffer)
Return a fresh basket by either reusing an existing basket that needs to be drop (according to TTree:...
Definition TBranch.cxx:1893
TBasket * GetBasketImpl(Int_t basket, TBuffer *user_buffer)
Return pointer to basket basketnumber in this Branch.
Definition TBranch.cxx:1225
Int_t fEntryOffsetLen
Initial Length of fEntryOffset table in the basket buffers.
Definition TBranch.h:119
virtual void DeleteBaskets(Option_t *option="")
Loop on all branch baskets.
Definition TBranch.cxx:725
virtual Long64_t GetBasketSeek(Int_t basket) const
Return address of basket in the file.
Definition TBranch.cxx:1301
TBranch()
Default constructor. Used for I/O by default.
Definition TBranch.cxx:86
void SetCompressionSettings(Int_t settings=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault)
Set compression settings.
Definition TBranch.cxx:2817
const char * GetIconName() const override
Return icon name depending on type of branch.
Definition TBranch.cxx:1331
Int_t BackFill()
Loop on all leaves of this branch to back fill Basket buffer.
Definition TBranch.cxx:679
~TBranch() override
Destructor.
Definition TBranch.cxx:449
Int_t fMaxBaskets
Maximum number of Baskets so far.
Definition TBranch.h:125
Long64_t fTotBytes
Total number of bytes in all leaves before compression.
Definition TBranch.h:136
TBuffer * fTransientBuffer
! Pointer to the current transient buffer.
Definition TBranch.h:151
virtual void ReadBasket(TBuffer &b)
Loop on all leaves of this branch to read Basket buffer.
Definition TBranch.cxx:2452
TTree * GetTree() const
Definition TBranch.h:252
static TClass * Class()
FillLeaves_t fFillLeaves
! Pointer to the FillLeaves implementation to use.
Definition TBranch.h:163
virtual TString GetFullName() const
Return the 'full' name of the branch.
Definition TBranch.cxx:2030
@ kAutoDelete
Definition TBranch.h:110
@ kDoNotUseBufferMap
If set, at least one of the entry in the branch will use the buffer's map of classname and objects.
Definition TBranch.h:112
@ kIsClone
To indicate a TBranchClones.
Definition TBranch.h:106
@ kDoNotProcess
Active bit for branches.
Definition TBranch.h:105
TObjArray fLeaves
-> List of leaves of this branch
Definition TBranch.h:139
Int_t GetBasketAndFirst(TBasket *&basket, Long64_t &first, TBuffer *user_buffer)
A helper function to locate the correct basket - and its first entry.
Definition TBranch.cxx:1352
char * fAddress
! Address of 1st leaf (variable or object)
Definition TBranch.h:147
virtual void DropBaskets(Option_t *option="")
Loop on all branch baskets.
Definition TBranch.cxx:756
TObjArray * GetListOfBranches()
Definition TBranch.h:246
virtual TList * GetBrowsables()
Returns (and, if 0, creates) browsable objects for this branch See TVirtualBranchBrowsable::FillListO...
Definition TBranch.cxx:1311
TList * fBrowsables
! List of TVirtualBranchBrowsables used for Browse()
Definition TBranch.h:152
void ReadLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to read Basket buffer.
Definition TBranch.cxx:2460
Int_t fOffset
Offset of this branch.
Definition TBranch.h:124
Long64_t * fBasketEntry
[fMaxBaskets] Table of first entry in each basket
Definition TBranch.h:142
Int_t GetEntriesSerialized(Long64_t N, TBuffer &user_buf)
Definition TBranch.h:185
void ExpandBasketArrays()
Increase BasketEntry buffer of a minimum of 10 locations and a maximum of 50 per cent of current size...
Definition TBranch.cxx:824
void Init(const char *name, const char *leaflist, Int_t compress)
Definition TBranch.cxx:299
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all leaves of entry and return total number of bytes read.
Definition TBranch.cxx:1705
TIOFeatures GetIOFeatures() const
Returns the IO settings currently in use for this branch.
Definition TBranch.cxx:2254
TClass * IsA() const override
Definition TBranch.h:295
void FillLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to fill Basket buffer.
Definition TBranch.cxx:2495
Long64_t fReadEntry
! Current entry number when reading
Definition TBranch.h:130
void Print(Option_t *option="") const override
Print TBranch parameters.
Definition TBranch.cxx:2340
static void ResetCount()
Static function resetting fgCount.
Definition TBranch.cxx:2673
TBranch * GetSubBranch(const TBranch *br) const
Find the parent branch of child.
Definition TBranch.cxx:2163
ReadLeaves_t fReadLeaves
! Pointer to the ReadLeaves implementation to use.
Definition TBranch.h:161
virtual void SetObject(void *objadd)
Set object this branch is pointing to.
Definition TBranch.cxx:2948
Int_t FlushBaskets()
Flush to disk all the baskets of this branch and any of subbranches.
Definition TBranch.cxx:1135
void ReadLeaves2Impl(TBuffer &b)
Read two leaves without the overhead of a loop.
Definition TBranch.cxx:2486
virtual void AddBasket(TBasket &b, bool ondisk, Long64_t startEntry)
Add the basket to this branch.
Definition TBranch.cxx:544
virtual void SetAddress(void *add)
Set address of this branch.
Definition TBranch.cxx:2694
static Int_t fgCount
! branch counter
Definition TBranch.h:116
virtual void AddLastBasket(Long64_t startEntry)
Add the start entry of the write basket (not yet created)
Definition TBranch.cxx:617
TBasket * GetBasket(Int_t basket)
Definition TBranch.h:213
Int_t fNBaskets
! Number of baskets in memory
Definition TBranch.h:126
void ReadLeaves1Impl(TBuffer &b)
Read one leaf without the overhead of a loop.
Definition TBranch.cxx:2478
Int_t GetBulkEntries(Long64_t, TBuffer &)
Read a basket of events into the given buffer with byte swapping.
Definition TBranch.cxx:1471
virtual void SetFile(TFile *file=nullptr)
Set file where this branch writes/reads its buffers.
Definition TBranch.cxx:2875
virtual Int_t GetEntryExport(Long64_t entry, Int_t getall, TClonesArray *list, Int_t n)
Read all leaves of an entry and export buffers to real objects in a TClonesArray list.
Definition TBranch.cxx:1761
virtual void SetAutoDelete(bool autodel=true)
Set the automatic delete bit.
Definition TBranch.cxx:2728
Long64_t fZipBytes
Total number of bytes in all leaves after compression.
Definition TBranch.h:137
TIOFeatures fIOFeatures
IO features for newly-created baskets.
Definition TBranch.h:123
void Browse(TBrowser *b) override
Browser interface.
Definition TBranch.cxx:698
void SetCompressionAlgorithm(Int_t algorithm=ROOT::RCompressionSetting::EAlgorithm::kUseGlobal)
Set compression algorithm.
Definition TBranch.cxx:2775
virtual void SetEntryOffsetLen(Int_t len, bool updateSubBranches=false)
Update the default value for the branch's fEntryOffsetLen if and only if it was already non zero (and...
Definition TBranch.cxx:2833
virtual TLeaf * FindLeaf(const char *name)
Find the leaf corresponding to the name 'searchname'.
Definition TBranch.cxx:1080
CacheInfo_t fCacheInfo
! Hold info about which basket are in the cache and if they have been retrieved from the cache.
Definition TBranch.h:158
TObjArray * GetListOfBaskets()
Definition TBranch.h:245
virtual void SetBufferAddress(TBuffer *entryBuffer)
Set address of this branch directly from a TBuffer to avoid streaming.
Definition TBranch.cxx:2757
Long64_t GetEntries() const
Definition TBranch.h:251
Int_t fNleaves
! Number of leaves
Definition TBranch.h:128
Int_t fSplitLevel
Branch split level.
Definition TBranch.h:127
Int_t WriteBasketImpl(TBasket *basket, Int_t where, ROOT::Internal::TBranchIMTHelper *)
Write the current basket to disk and return the number of bytes written to the file.
Definition TBranch.cxx:3214
virtual void UpdateFile()
Refresh the value of fDirectory (i.e.
Definition TBranch.cxx:3316
Int_t * fBasketBytes
[fMaxBaskets] Length of baskets on file
Definition TBranch.h:141
Long64_t fNextBasketEntry
! Next entry that will requires us to go to the next basket
Definition TBranch.h:132
bool IsAutoDelete() const
Return true if an existing object in a TBranchObject must be deleted.
Definition TBranch.cxx:2262
Int_t FillEntryBuffer(TBasket *basket, TBuffer *buf, Int_t &lnew)
Copy the data from fEntryBuffer into the current basket.
Definition TBranch.cxx:934
virtual TFile * GetFile(Int_t mode=0)
Return pointer to the file where branch buffers reside, returns 0 in case branch buffers reside in th...
Definition TBranch.cxx:1852
TObjArray fBranches
-> List of Branches of this branch
Definition TBranch.h:138
virtual void KeepCircular(Long64_t maxEntries)
keep a maximum of fMaxEntries in memory
Definition TBranch.cxx:2282
virtual void SetStatus(bool status=true)
Set branch status to Process or DoNotProcess.
Definition TBranch.cxx:2959
virtual void ResetAfterMerge(TFileMergeInfo *)
Reset a Branch.
Definition TBranch.cxx:2597
void ReadLeaves0Impl(TBuffer &b)
Read zero leaves without the overhead of a loop.
Definition TBranch.cxx:2471
bool SupportsBulkRead() const
Returns true if this branch supports bulk IO, false otherwise.
Definition TBranch.cxx:1430
TString GetRealFileName() const
Get real file name.
Definition TBranch.cxx:2067
virtual TBranch * FindBranch(const char *name)
Find the immediate sub-branch with passed name.
Definition TBranch.cxx:1034
TDirectory * fDirectory
! Pointer to directory where this branch buffers are stored
Definition TBranch.h:148
void PrintCacheInfo() const
Print the information we have about which basket is currently cached and whether they have been 'used...
Definition TBranch.cxx:2444
TObjArray fBaskets
-> List of baskets of this branch
Definition TBranch.h:140
virtual Int_t LoadBaskets()
Baskets associated to this branch are forced to be in memory.
Definition TBranch.cxx:2308
TBranch * fMother
! Pointer to top-level parent branch in the tree.
Definition TBranch.h:145
Long64_t GetTotBytes(Option_t *option="") const
Return total number of bytes in the branch (excluding current buffer) if option ="*" includes all sub...
Definition TBranch.cxx:2219
TBranch * fParent
! Pointer to parent branch.
Definition TBranch.h:146
Int_t WriteBasket(TBasket *basket, Int_t where)
Definition TBranch.h:175
Int_t FlushOneBasket(UInt_t which)
If we have a write basket in memory and it contains some entries and has not yet been written to disk...
Definition TBranch.cxx:1181
bool fSkipZip
! After being read, the buffer will not be unzipped.
Definition TBranch.h:155
virtual Int_t GetExpectedType(TClass *&clptr, EDataType &type)
Fill expectedClass and expectedType with information on the data type of the object/values contained ...
Definition TBranch.cxx:1833
bool IsFolder() const override
Return true if more than one leaf or browsables, false otherwise.
Definition TBranch.cxx:2270
virtual void SetFirstEntry(Long64_t entry)
set the first entry number (case of TBranchSTL)
Definition TBranch.cxx:3291
Long64_t GetTotalSize(Option_t *option="") const
Return total number of bytes in the branch (including current buffer)
Definition TBranch.cxx:2200
Long64_t GetZipBytes(Option_t *option="") const
Return total number of zip bytes in the branch if option ="*" includes all sub-branches of this branc...
Definition TBranch.cxx:2237
virtual void Refresh(TBranch *b)
Refresh this branch using new information in b This function is called by TTree::Refresh.
Definition TBranch.cxx:2507
virtual bool SetMakeClass(bool decomposeObj=true)
Set the branch in a mode where the object are decomposed (Also known as MakeClass mode).
Definition TBranch.cxx:2939
Int_t fWriteBasket
Last basket number written.
Definition TBranch.h:120
Long64_t * fBasketSeek
[fMaxBaskets] Addresses of baskets on file
Definition TBranch.h:143
TObjArray * GetListOfLeaves()
Definition TBranch.h:247
virtual void SetEntries(Long64_t entries)
Set the number of entries in this branch.
Definition TBranch.cxx:2850
Int_t fReadBasket
! Current basket number when reading
Definition TBranch.h:129
Long64_t fFirstEntry
Number of the first entry in this branch.
Definition TBranch.h:135
TBasket * fExtraBasket
! Allocated basket not currently holding any data.
Definition TBranch.h:122
virtual Int_t GetRow(Int_t row)
Return all elements of one row unpacked in internal array fValues [Actually just returns 1 (?...
Definition TBranch.cxx:2107
Int_t fBasketSize
Initial Size of Basket Buffer.
Definition TBranch.h:118
Int_t Fill()
Definition TBranch.h:205
virtual void Reset(Option_t *option="")
Reset a Branch.
Definition TBranch.cxx:2556
virtual void SetBasketSize(Int_t bufsize)
Set the basket size The function makes sure that the basket size is greater than fEntryOffsetlen.
Definition TBranch.cxx:2741
Long64_t fEntryNumber
Current entry number (last one filled in this branch)
Definition TBranch.h:121
TBranch * GetMother() const
Get our top-level parent branch in the tree.
Definition TBranch.cxx:2126
Int_t fCompress
Compression level and algorithm.
Definition TBranch.h:117
TBuffer * GetTransientBuffer(Int_t size)
Returns the transient buffer currently used by this TBranch for reading/writing baskets.
Definition TBranch.cxx:522
virtual Int_t FillImpl(ROOT::Internal::TBranchIMTHelper *)
Loop on all leaves of this branch to fill Basket buffer.
Definition TBranch.cxx:855
TBuffer * fEntryBuffer
! Buffer used to directly pass the content without streaming
Definition TBranch.h:150
TBasket * fCurrentBasket
! Pointer to the current basket.
Definition TBranch.h:133
Long64_t fFirstBasketEntry
! First entry in the current basket.
Definition TBranch.h:131
void SetCompressionLevel(Int_t level=ROOT::RCompressionSetting::ELevel::kUseMin)
Set compression level.
Definition TBranch.cxx:2795
void Streamer(TBuffer &) override
Stream a class object.
Definition TBranch.cxx:2968
Long64_t fEntries
Number of entries.
Definition TBranch.h:134
TBasket * GetFreshCluster(TBuffer *user_buffer)
Drops the cluster two behind the current cluster and returns a fresh basket by either reusing or crea...
Definition TBranch.cxx:1952
TTree * fTree
! Pointer to Tree header
Definition TBranch.h:144
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
The concrete implementation of TBuffer for writing/reading to/from a ROOT file or socket.
Definition TBufferFile.h:47
@ kNotDecompressed
Definition TBufferIO.h:66
Buffer base class used for serializing objects.
Definition TBuffer.h:43
void SetWriteMode()
Set buffer in write mode.
Definition TBuffer.cxx:314
void Expand(Int_t newsize, Bool_t copy=kTRUE)
Expand (or shrink) the I/O buffer to newsize bytes.
Definition TBuffer.cxx:222
virtual Int_t GetBufferDisplacement() const =0
Int_t BufferSize() const
Definition TBuffer.h:98
@ kIsOwner
Definition TBuffer.h:75
@ kWrite
Definition TBuffer.h:73
@ kRead
Definition TBuffer.h:73
@ kMinimalSize
Definition TBuffer.h:78
Bool_t IsReading() const
Definition TBuffer.h:86
virtual void ResetMap()=0
void SetBufferOffset(Int_t offset=0)
Definition TBuffer.h:93
virtual Int_t GetMapCount() const =0
void SetReadMode()
Set buffer in read mode.
Definition TBuffer.cxx:301
virtual void WriteBuf(const void *buf, Int_t max)=0
TClass * IsA() const override
Definition TBuffer.h:340
virtual void SetBufferDisplacement()=0
virtual char * ReadString(char *s, Int_t max)=0
Int_t Length() const
Definition TBuffer.h:100
char * Buffer() const
Definition TBuffer.h:96
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
An array of clone (identical) objects.
void Browse(TBrowser *b) override
Browse this collection (called by TBrowser).
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
Describe directory structure in memory.
Definition TDirectory.h:45
virtual TFile * GetFile() const
Definition TDirectory.h:221
TObject * FindObject(const char *name) const override
Find object by name in the list of memory objects.
virtual Bool_t IsWritable() const
Definition TDirectory.h:238
A cache when reading files over the network.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
Int_t GetCompressionLevel() const
Definition TFile.h:474
virtual void MakeFree(Long64_t first, Long64_t last)
Mark unused bytes on the file.
Definition TFile.cxx:1473
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.
Definition TFile.cxx:3764
void Close(Option_t *option="") override
Close a file.
Definition TFile.cxx:958
A TLeaf for an 8 bit Integer data type.
Definition TLeafB.h:26
A TLeaf for a variable length string.
Definition TLeafC.h:26
static TClass * Class()
A TLeaf for a 24 bit truncated floating point data type.
Definition TLeafD32.h:28
A TLeaf for a 64 bit floating point data type.
Definition TLeafD.h:26
A TLeaf for a 24 bit truncated floating point data type.
Definition TLeafF16.h:27
A TLeaf for a 32 bit floating point data type.
Definition TLeafF.h:26
A TLeaf for a long integer data type (Long_t, non-portable size).
Definition TLeafG.h:27
A TLeaf for an Integer data type.
Definition TLeafI.h:27
A TLeaf for a 64 bit Integer data type.
Definition TLeafL.h:27
A TLeaf for a bool data type.
Definition TLeafO.h:26
A TLeaf for a 16 bit Integer data type.
Definition TLeafS.h:26
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition TLeaf.h:57
A doubly linked list.
Definition TList.h:38
static TClass * Class()
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
void Streamer(TBuffer &) override
Stream an object of class TObject.
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
TString fName
Definition TNamed.h:32
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
virtual void Expand(Int_t newSize)
Expand or shrink the array to newSize elements.
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
void Streamer(TBuffer &) override
Stream all objects in the array to or from the I/O buffer.
Int_t GetEntries() const override
Return the number of objects in array (i.e.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
TObject * Remove(TObject *obj) override
Remove object from array.
void SetLast(Int_t last)
Set index of last object in array, effectively truncating the array.
TObject * RemoveAt(Int_t idx) override
Remove object at index idx.
Int_t GetLast() const override
Return index of last object in array.
Int_t LowerBound() const
Definition TObjArray.h:91
void Add(TObject *obj) override
Definition TObjArray.h:68
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1057
R__ALWAYS_INLINE Bool_t IsZombie() const
Definition TObject.h:159
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1099
void MakeZombie()
Definition TObject.h:53
void ResetBit(UInt_t f)
Definition TObject.h:201
static Int_t * ReAllocInt(Int_t *vp, size_t size, size_t oldsize)
Reallocate (i.e.
Definition TStorage.cxx:257
static void * ReAlloc(void *vp, size_t size, size_t oldsize)
Reallocate (i.e.
Definition TStorage.cxx:182
Basic string class.
Definition TString.h:138
Ssiz_t Length() const
Definition TString.h:425
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
const char * Data() const
Definition TString.h:384
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1418
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition TSystem.cxx:1285
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition TSystem.cxx:944
virtual Bool_t IsAbsoluteFileName(const char *dir)
Return true if dir is an absolute pathname.
Definition TSystem.cxx:961
virtual TString GetDirName(const char *pathname)
Return the directory name in pathname.
Definition TSystem.cxx:1042
Helper class to iterate over cluster of baskets.
Definition TTree.h:306
Long64_t Previous()
Move on to the previous cluster and return the starting entry of this previous cluster.
Definition TTree.cxx:721
Long64_t GetStartEntry()
Definition TTree.h:338
Long64_t Next()
Move on to the next cluster and return the starting entry of this next cluster.
Definition TTree.cxx:677
A TTree represents a columnar dataset.
Definition TTree.h:89
virtual TVirtualPerfStats * GetPerfStats() const
Definition TTree.h:585
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition TTree.cxx:5543
void AddAllocationCount(UInt_t count)
Definition TTree.h:373
virtual TObjArray * GetListOfLeaves()
Definition TTree.h:568
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition TTree.cxx:5555
void Draw(Option_t *opt) override
Default Draw method for all objects.
Definition TTree.h:470
virtual void IncrementTotalBuffers(Int_t nbytes)
Definition TTree.h:625
TDirectory * GetDirectory() const
Definition TTree.h:501
TTreeCache * GetReadCache(TFile *file) const
Find and return the TTreeCache registered with the file and which may contain branches for us.
Definition TTree.cxx:6399
virtual bool GetClusterPrefetch() const
Definition TTree.h:496
virtual TObjArray * GetListOfBranches()
Definition TTree.h:567
virtual void AddZipBytes(Int_t zip)
Definition TTree.h:368
virtual TBasket * CreateBasket(TBranch *)
Create a basket for this tree and given branch.
Definition TTree.cxx:3771
@ kOnlyFlushAtCluster
If set, the branch's buffers will grow until an event cluster boundary is hit, guaranteeing a basket ...
Definition TTree.h:292
@ kCircular
Definition TTree.h:288
virtual void AddTotBytes(Int_t tot)
Definition TTree.h:367
virtual Long64_t GetAutoFlush() const
Definition TTree.h:486
virtual Long64_t GetMaxVirtualSize() const
Definition TTree.h:579
This class represents a WWW compatible URL.
Definition TUrl.h:33
static Int_t FillListOfBrowsables(TList &list, const TBranch *branch, const TVirtualBranchBrowsable *parent=nullptr)
Askes all registered generators to fill their browsables into the list.
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:348
@ kUndefined
Undefined compression algorithm (must be kept the last of the list in case a new algorithm is added).
@ kUseMin
Compression level reserved when we are not sure what to use (1 is for the fastest compression)
Definition Compression.h:72
TLine l
Definition textangle.C:4