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