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).
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 if (ondisk) {
576 fBasketBytes[where] = basket->GetNbytes(); // not for in mem
577 fBasketSeek[where] = basket->GetSeekKey(); // not for in mem
579 ++fWriteBasket;
580 } else {
581 ++fNBaskets;
584 }
585
586 fEntries += basket->GetNevBuf();
587 fEntryNumber += basket->GetNevBuf();
588 if (ondisk) {
589 fTotBytes += basket->GetObjlen() + basket->GetKeylen() ;
590 fZipBytes += basket->GetNbytes();
591 fTree->AddTotBytes(basket->GetObjlen() + basket->GetKeylen());
592 fTree->AddZipBytes(basket->GetNbytes());
593 }
594}
595
596////////////////////////////////////////////////////////////////////////////////
597/// Add the start entry of the write basket (not yet created)
598
600{
601 if (fWriteBasket >= fMaxBaskets) {
603 }
604 Int_t where = fWriteBasket;
605
606 if (where && startEntry < fBasketEntry[where-1]) {
607 // Need to find the right location and move the possible baskets
608
609 Fatal("AddBasket","The last basket must have the highest entry number (%s/%lld/%d).",GetName(),startEntry,fWriteBasket);
610
611 }
612 fBasketEntry[where] = startEntry;
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Loop on all leaves of this branch to back fill Basket buffer.
618///
619/// Use this routine instead of TBranch::Fill when filling a branch individually
620/// to catch up with the number of entries already in the TTree.
621///
622/// First it calls TBranch::Fill and then if the number of entries of the branch
623/// reach one of TTree cluster's boundary, the basket is flushed.
624///
625/// The function returns the number of bytes committed to the memory basket.
626/// If a write error occurs, the number of bytes returned is -1.
627/// If no data are written, because e.g. the branch is disabled,
628/// the number of bytes returned is 0.
629///
630/// To insure that the baskets of each cluster are located close by in the
631/// file, when back-filling multiple branches make sure to call BackFill
632/// for the same entry for all the branches consecutively
633/// ~~~ {.cpp}
634/// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
635/// for( auto branch : branchCollection) {
636/// ... Make change to the data associated with the branch ...
637/// branch->BackFill();
638/// }
639/// }
640/// // Since we loop over all the branches for each new entry
641/// // all the baskets for a cluster are consecutive in the file.
642/// ~~~
643/// rather than doing all the entries of one branch at a time.
644/// ~~~ {.cpp}
645/// // Do NOT do things in the following order, it will lead to
646/// // poorly clustered files.
647/// for(auto branch : branchCollection) {
648/// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
649/// ... Make change to the data associated with the branch ...
650/// branch->BackFill();
651/// }
652/// }
653/// // Since we loop over all the entries for one branch
654/// // all the baskets for that branch are consecutive.
655/// ~~~
656
658
659 // Get the end of the next cluster.
660 auto cluster = GetTree()->GetClusterIterator( GetEntries() );
661 cluster.Next();
662 auto endCluster = cluster.GetNextEntry();
663
664 auto result = FillImpl(nullptr);
665
666 if ( result && GetEntries() >= endCluster ) {
667 FlushBaskets();
668 }
669
670 return result;
671}
672
673////////////////////////////////////////////////////////////////////////////////
674/// Browser interface.
675
677{
678 if (fNleaves > 1) {
680 } else {
681 // Get the name and strip any extra brackets
682 // in order to get the full arrays.
683 TString name = GetName();
684 Int_t pos = name.First('[');
685 if (pos!=kNPOS) name.Remove(pos);
686
687 GetTree()->Draw(name, "", b ? b->GetDrawOption() : "");
688 if (gPad) gPad->Update();
689 }
690}
691
692 ///////////////////////////////////////////////////////////////////////////////
693 /// Loop on all branch baskets. If the file where branch buffers reside is
694 /// writable, free the disk space associated to the baskets of the branch,
695 /// then call Reset(). If the option contains "all", delete also the baskets
696 /// for the subbranches.
697 /// The branch is reset.
698 ///
699 /// NOTE that this function must be used with extreme care. Deleting branch baskets
700 /// fragments the file and may introduce inefficiencies when adding new entries
701 /// in the Tree or later on when reading the Tree.
702
704{
705 TString opt = option;
706 opt.ToLower();
707 TFile *file = GetFile(0);
708
710 for(Int_t i=0; i<fWriteBasket; i++) {
711 if (fBasketSeek[i]) file->MakeFree(fBasketSeek[i],fBasketSeek[i]+fBasketBytes[i]-1);
712 }
713 }
714
715 // process subbranches
716 if (opt.Contains("all")) {
718 Int_t nb = lb->GetEntriesFast();
719 for (Int_t j = 0; j < nb; j++) {
720 TBranch* branch = (TBranch*) lb->UncheckedAt(j);
721 if (branch) branch->DeleteBaskets("all");
722 }
723 }
724 DropBaskets("all");
725 Reset();
726}
727
728////////////////////////////////////////////////////////////////////////////////
729/// Loop on all branch baskets. Drop all baskets from memory except readbasket.
730/// If the option contains "all", drop all baskets including
731/// read- and write-baskets (unless they are not stored individually on disk).
732/// The option "all" also lead to DropBaskets being called on the sub-branches.
733
735{
736 Bool_t all = kFALSE;
737 if (options && options[0]) {
738 TString opt = options;
739 opt.ToLower();
740 if (opt.Contains("all")) all = kTRUE;
741 }
742
743 TBasket *basket;
744 Int_t nbaskets = fBaskets.GetEntriesFast();
745
746 if ( (fNBaskets>1) || all ) {
747 //slow case
748 for (Int_t i=0;i<nbaskets;i++) {
749 basket = (TBasket*)fBaskets.UncheckedAt(i);
750 if (!basket) continue;
751 if ((i == fReadBasket || i == fWriteBasket) && !all) continue;
752 // if the basket is not yet on file but already has event in it
753 // we must continue to avoid dropping the basket (and thus losing data)
754 if (fBasketBytes[i]==0 && basket->GetNevBuf() > 0) continue;
755 basket->DropBuffers();
756 --fNBaskets;
758 if (basket == fCurrentBasket) {
759 fCurrentBasket = 0;
761 fNextBasketEntry = -1;
762 }
763 delete basket;
764 }
765
766 // process subbranches
767 if (all) {
769 Int_t nb = lb->GetEntriesFast();
770 for (Int_t j = 0; j < nb; j++) {
771 TBranch* branch = (TBranch*) lb->UncheckedAt(j);
772 if (!branch) continue;
773 branch->DropBaskets("all");
774 }
775 }
776 } else {
777 //fast case
778 if (nbaskets > 0) {
779 Int_t i = fBaskets.GetLast();
780 basket = (TBasket*)fBaskets.UncheckedAt(i);
781 if (basket && fBasketBytes[i]!=0) {
782 basket->DropBuffers();
783 if (basket == fCurrentBasket) {
784 fCurrentBasket = 0;
786 fNextBasketEntry = -1;
787 }
788 delete basket;
789 fBaskets.AddAt(0,i);
790 fBaskets.SetLast(-1);
791 fNBaskets = 0;
792 }
793 }
794 }
795
796}
797
798////////////////////////////////////////////////////////////////////////////////
799/// Increase BasketEntry buffer of a minimum of 10 locations
800/// and a maximum of 50 per cent of current size.
801
803{
804 Int_t newsize = TMath::Max(10,Int_t(1.5*fMaxBaskets));
807 newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
809 newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
810
811 fMaxBaskets = newsize;
812
813 fBaskets.Expand(newsize);
814
815 for (Int_t i=fWriteBasket;i<fMaxBaskets;i++) {
816 fBasketBytes[i] = 0;
817 fBasketEntry[i] = 0;
818 fBasketSeek[i] = 0;
819 }
820}
821
822////////////////////////////////////////////////////////////////////////////////
823/// Loop on all leaves of this branch to fill Basket buffer.
824///
825/// If TBranchIMTHelper is non-null and it is time to WriteBasket, then we will
826/// use TBB to compress in parallel.
827///
828/// The function returns the number of bytes committed to the memory basket.
829/// If a write error occurs, the number of bytes returned is -1.
830/// If no data are written, because e.g. the branch is disabled,
831/// the number of bytes returned is 0.
832
834{
835 if (TestBit(kDoNotProcess)) {
836 return 0;
837 }
838
840 if (!basket) {
841 basket = fTree->CreateBasket(this); // create a new basket
842 if (!basket) return 0;
843 ++fNBaskets;
845 }
846 TBuffer* buf = basket->GetBufferRef();
847
848 // Fill basket buffer.
849
850 Int_t nsize = 0;
851
852 if (buf->IsReading()) {
853 basket->SetWriteMode();
854 }
855
857 buf->ResetMap();
858 }
859
860 Int_t lnew = 0;
861 Int_t nbytes = 0;
862
863 if (fEntryBuffer) {
864 nbytes = FillEntryBuffer(basket,buf,lnew);
865 } else {
866 Int_t lold = buf->Length();
867 basket->Update(lold);
868 ++fEntries;
869 ++fEntryNumber;
870 (this->*fFillLeaves)(*buf);
871 if (buf->GetMapCount()) {
872 // The map is used.
874 }
875 lnew = buf->Length();
876 nbytes = lnew - lold;
877 }
878
879 if (fEntryOffsetLen) {
880 Int_t nevbuf = basket->GetNevBuf();
881 // Total size in bytes of EntryOffset table.
882 nsize = nevbuf * sizeof(Int_t);
883 } else {
884 if (!basket->GetNevBufSize()) {
885 basket->SetNevBufSize(nbytes);
886 }
887 }
888
889 // Should we create a new basket?
890 // fSkipZip force one entry per buffer (old stuff still maintained for CDF)
891 // Transfer full compressed buffer only
892
893 // If GetAutoFlush() is less than zero, then we are determining the end of the autocluster
894 // based upon the number of bytes already flushed. This is incompatible with one-basket-per-cluster
895 // (since we will grow the basket indefinitely and never flush!). Hence, we wait until the
896 // first event cluster is written out and *then* enable one-basket-per-cluster mode.
897 bool noFlushAtCluster = !fTree->TestBit(TTree::kOnlyFlushAtCluster) || (fTree->GetAutoFlush() < 0);
898
899 if (noFlushAtCluster && !fTree->TestBit(TTree::kCircular) &&
901 ((lnew + (2 * nsize) + nbytes) >= fBasketSize))) {
902 Int_t nout = WriteBasketImpl(basket, fWriteBasket, imtHelper);
903 if (nout < 0) Error("TBranch::Fill", "Failed to write out basket.\n");
904 return (nout >= 0) ? nbytes : -1;
905 }
906 return nbytes;
907}
908
909////////////////////////////////////////////////////////////////////////////////
910/// Copy the data from fEntryBuffer into the current basket.
911
913{
914 Int_t nbytes = 0;
915 Int_t objectStart = 0;
916 Int_t last = 0;
917 Int_t lold = buf->Length();
918
919 // Handle the special case of fEntryBuffer != 0
920 if (fEntryBuffer->IsA() == TMessage::Class()) {
921 objectStart = 8;
922 }
924 // The buffer given as input has not been decompressed.
925 if (basket->GetNevBuf()) {
926 // If the basket already contains entry we need to close it
927 // out. (This is because we can only transfer full compressed
928 // buffer)
930 // And restart from scratch
931 return Fill();
932 }
933 Int_t startpos = fEntryBuffer->Length();
935 static TBasket toread_fLast;
937 toread_fLast.Streamer(*fEntryBuffer);
939 last = toread_fLast.GetLast();
940 // last now contains the decompressed number of bytes.
941 fEntryBuffer->SetBufferOffset(startpos);
942 buf->SetBufferOffset(0);
944 basket->Update(lold);
945 } else {
946 // We are required to copy starting at the version number (so not
947 // including the class name.
948 // See if byte count is here, if not it class still be a newClass
949 const UInt_t kNewClassTag = 0xFFFFFFFF;
950 const UInt_t kByteCountMask = 0x40000000; // OR the byte count with this
951 UInt_t tag = 0;
952 UInt_t startpos = fEntryBuffer->Length();
953 fEntryBuffer->SetBufferOffset(objectStart);
954 *fEntryBuffer >> tag;
955 if (tag & kByteCountMask) {
956 *fEntryBuffer >> tag;
957 }
958 if (tag == kNewClassTag) {
959 UInt_t maxsize = 256;
960 char* s = new char[maxsize];
961 Int_t name_start = fEntryBuffer->Length();
962 fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end.
963 while (strlen(s) == (maxsize - 1)) {
964 // The classname is too large, try again with a large buffer.
965 fEntryBuffer->SetBufferOffset(name_start);
966 maxsize *= 2;
967 delete[] s;
968 s = new char[maxsize];
969 fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end
970 }
971 } else {
972 fEntryBuffer->SetBufferOffset(objectStart);
973 }
974 objectStart = fEntryBuffer->Length();
975 fEntryBuffer->SetBufferOffset(startpos);
976 basket->Update(lold, objectStart - fEntryBuffer->GetBufferDisplacement());
977 }
978 fEntries++;
979 fEntryNumber++;
980 UInt_t len = 0;
981 UInt_t startpos = fEntryBuffer->Length();
982 if (startpos > UInt_t(objectStart)) {
983 // We assume this buffer have just been directly filled
984 // the current position in the buffer indicates the end of the object!
985 len = fEntryBuffer->Length() - objectStart;
986 } else {
987 // The buffer have been acquired either via TSocket or via
988 // TBuffer::SetBuffer(newloc,newsize)
989 // Only the actual size of the memory buffer gives us an hint about where
990 // the object ends.
991 len = fEntryBuffer->BufferSize() - objectStart;
992 }
993 buf->WriteBuf(fEntryBuffer->Buffer() + objectStart, len);
995 // The original buffer came pre-compressed and thus the buffer Length
996 // does not really show the really object size
997 // lnew = nbytes = basket->GetLast();
998 nbytes = last;
999 lnew = last;
1000 } else {
1001 lnew = buf->Length();
1002 nbytes = lnew - lold;
1003 }
1004
1005 return nbytes;
1006}
1007
1008////////////////////////////////////////////////////////////////////////////////
1009/// Find the immediate sub-branch with passed name.
1010
1012{
1013 // We allow the user to pass only the last dotted component of the name.
1014 std::string longnm;
1015 longnm.reserve(fName.Length()+strlen(name)+3);
1016 longnm = fName.Data();
1017 if (longnm[longnm.length()-1]==']') {
1018 std::size_t dim = longnm.find_first_of("[");
1019 if (dim != std::string::npos) {
1020 longnm.erase(dim);
1021 }
1022 }
1023 if (longnm[longnm.length()-1] != '.') {
1024 longnm += '.';
1025 }
1026 longnm += name;
1027 UInt_t namelen = strlen(name);
1028
1029 Int_t nbranches = fBranches.GetEntries();
1030 TBranch* branch = 0;
1031 for(Int_t i = 0; i < nbranches; ++i) {
1032 branch = (TBranch*) fBranches.UncheckedAt(i);
1033
1034 const char *brname = branch->fName.Data();
1035 UInt_t brlen = branch->fName.Length();
1036 if (brname[brlen-1]==']') {
1037 const char *dim = strchr(brname,'[');
1038 if (dim) {
1039 brlen = dim - brname;
1040 }
1041 }
1042 if (namelen == brlen /* same effective size */
1043 && strncmp(name,brname,brlen) == 0) {
1044 return branch;
1045 }
1046 if (brlen == (size_t)longnm.length()
1047 && strncmp(longnm.c_str(),brname,brlen) == 0) {
1048 return branch;
1049 }
1050 }
1051 return 0;
1052}
1053
1054////////////////////////////////////////////////////////////////////////////////
1055/// Find the leaf corresponding to the name 'searchname'.
1056
1057TLeaf* TBranch::FindLeaf(const char* searchname)
1058{
1059 TString leafname;
1060 TString leaftitle;
1061 TString longname;
1062 TString longtitle;
1063
1064 // We allow the user to pass only the last dotted component of the name.
1065 TIter next(GetListOfLeaves());
1066 TLeaf* leaf = 0;
1067 while ((leaf = (TLeaf*) next())) {
1068 leafname = leaf->GetName();
1069 Ssiz_t dim = leafname.First('[');
1070 if (dim >= 0) leafname.Remove(dim);
1071
1072 if (leafname == searchname) return leaf;
1073
1074 // The leaf element contains the branch name in its name, let's use the title.
1075 leaftitle = leaf->GetTitle();
1076 dim = leaftitle.First('[');
1077 if (dim >= 0) leaftitle.Remove(dim);
1078
1079 if (leaftitle == searchname) return leaf;
1080
1081 TBranch* branch = leaf->GetBranch();
1082 if (branch) {
1083 longname.Form("%s.%s",branch->GetName(),leafname.Data());
1084 dim = longname.First('[');
1085 if (dim>=0) longname.Remove(dim);
1086 if (longname == searchname) return leaf;
1087
1088 // The leaf element contains the branch name in its name.
1089 longname.Form("%s.%s",branch->GetName(),searchname);
1090 if (longname==leafname) return leaf;
1091
1092 longtitle.Form("%s.%s",branch->GetName(),leaftitle.Data());
1093 dim = longtitle.First('[');
1094 if (dim>=0) longtitle.Remove(dim);
1095 if (longtitle == searchname) return leaf;
1096
1097 // The following is for the case where the branch is only
1098 // a sub-branch. Since we do not see it through
1099 // TTree::GetListOfBranches, we need to see it indirectly.
1100 // This is the less sturdy part of this search ... it may
1101 // need refining ...
1102 if (strstr(searchname, ".") && !strcmp(searchname, branch->GetName())) return leaf;
1103 }
1104 }
1105 return 0;
1106}
1107
1108////////////////////////////////////////////////////////////////////////////////
1109/// Flush to disk all the baskets of this branch and any of subbranches.
1110/// Return the number of bytes written or -1 in case of write error.
1111
1113{
1114 UInt_t nerror = 0;
1115 Int_t nbytes = 0;
1116
1117 Int_t maxbasket = fWriteBasket + 1;
1118 // The following protection is not necessary since we should always
1119 // have fWriteBasket < fBasket.GetSize()
1120 //if (fBaskets.GetSize() < maxbasket) {
1121 // maxbasket = fBaskets.GetSize();
1122 //}
1123 for(Int_t i=0; i != maxbasket; ++i) {
1124 if (fBaskets.UncheckedAt(i)) {
1125 Int_t nwrite = FlushOneBasket(i);
1126 if (nwrite<0) {
1127 ++nerror;
1128 } else {
1129 nbytes += nwrite;
1130 }
1131 }
1132 }
1134 for (Int_t i = 0; i < len; ++i) {
1135 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1136 if (!branch) {
1137 continue;
1138 }
1139 Int_t nwrite = branch->FlushBaskets();
1140 if (nwrite<0) {
1141 ++nerror;
1142 } else {
1143 nbytes += nwrite;
1144 }
1145 }
1146 if (nerror) {
1147 return -1;
1148 } else {
1149 return nbytes;
1150 }
1151}
1152
1153////////////////////////////////////////////////////////////////////////////////
1154/// If we have a write basket in memory and it contains some entries and
1155/// has not yet been written to disk, we write it and delete it from memory.
1156/// Return the number of bytes written;
1157
1159{
1160 Int_t nbytes = 0;
1161 if (fDirectory && fBaskets.GetEntries()) {
1162 TBasket *basket = (TBasket*)fBaskets.UncheckedAt(ibasket);
1163
1164 if (basket) {
1165 if (basket->GetNevBuf()
1166 && fBasketSeek[ibasket]==0) {
1167 // If the basket already contains entry we need to close it out.
1168 // (This is because we can only transfer full compressed buffer)
1169
1170 if (basket->GetBufferRef()->IsReading()) {
1171 basket->SetWriteMode();
1172 }
1173 nbytes = WriteBasket(basket,ibasket);
1174
1175 } else {
1176 // If the basket is empty or has already been written.
1177 if ((Int_t)ibasket==fWriteBasket) {
1178 // Nothing to do.
1179 } else {
1180 basket->DropBuffers();
1181 if (basket == fCurrentBasket) {
1182 fCurrentBasket = 0;
1183 fFirstBasketEntry = -1;
1184 fNextBasketEntry = -1;
1185 }
1186 delete basket;
1187 --fNBaskets;
1188 fBaskets[ibasket] = 0;
1189 }
1190 }
1191 }
1192 }
1193 return nbytes;
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Return pointer to basket basketnumber in this Branch
1198///
1199/// If a new buffer must be created and the user_buffer argument is non-null,
1200/// then the memory in the user_bufer will be shared with the returned TBasket.
1201
1202TBasket* TBranch::GetBasketImpl(Int_t basketnumber, TBuffer *user_buffer)
1203{
1204 // This counter in the sequential case collects errors coming also from
1205 // different files (suppose to have a program reading f1.root, f2.root ...)
1206 // In the mt case, it is made atomic: it safely collects errors from
1207 // different files processed simultaneously.
1208 static std::atomic<Int_t> nerrors(0);
1209
1210 // reference to an existing basket in memory ?
1211 if (basketnumber <0 || basketnumber > fWriteBasket) return 0;
1212 TBasket *basket = (TBasket*)fBaskets.UncheckedAt(basketnumber);
1213 if (basket) return basket;
1214 if (basketnumber == fWriteBasket) return 0;
1215
1216 // create/decode basket parameters from buffer
1217 TFile *file = GetFile(0);
1218 if (file == 0) {
1219 return 0;
1220 }
1221 // if cluster pre-fetching or retaining is on, do not re-use existing baskets
1222 // unless a new cluster is used.
1224 basket = GetFreshCluster();
1225 else
1226 basket = GetFreshBasket(basketnumber, user_buffer);
1227
1228 // fSkipZip is old stuff still maintained for CDF
1230 if (fBasketBytes[basketnumber] == 0) {
1231 fBasketBytes[basketnumber] = basket->ReadBasketBytes(fBasketSeek[basketnumber],file);
1232 }
1233 //add branch to cache (if any)
1234 {
1235 R__LOCKGUARD_IMT(gROOTMutex); // Lock for parallel TTree I/O
1237 if (pf){
1238 if (pf->IsLearning()) pf->LearnBranch(this, kFALSE);
1239 if (fSkipZip) pf->SetSkipZip();
1240 }
1241 }
1242
1243 //now read basket
1244 Int_t badread = basket->ReadBasketBuffers(fBasketSeek[basketnumber],fBasketBytes[basketnumber],file);
1245 if (R__unlikely(badread || basket->GetSeekKey() != fBasketSeek[basketnumber] || basket->IsZombie())) {
1246 nerrors++;
1247 if (nerrors > 10) return 0;
1248 if (nerrors == 10) {
1249 printf(" file probably overwritten: stopping reporting error messages\n");
1250 if (fBasketSeek[basketnumber] > 2000000000) {
1251 printf("===>File is more than 2 Gigabytes\n");
1252 return 0;
1253 }
1254 if (fBasketSeek[basketnumber] > 1000000000) {
1255 printf("===>Your file is may be bigger than the maximum file size allowed on your system\n");
1256 printf(" Check your AFS maximum file size limit for example\n");
1257 return 0;
1258 }
1259 }
1260 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);
1261 return 0;
1262 }
1263
1264 ++fNBaskets;
1265
1266 fCacheInfo.SetUsed(basketnumber);
1267 auto perfStats = GetTree()->GetPerfStats();
1268 if (perfStats)
1269 perfStats->SetUsed(this, basketnumber);
1270
1271 fBaskets.AddAt(basket,basketnumber);
1272 return basket;
1273}
1274
1275////////////////////////////////////////////////////////////////////////////////
1276/// Return address of basket in the file
1277
1279{
1280 if (basketnumber <0 || basketnumber > fWriteBasket) return 0;
1281 return fBasketSeek[basketnumber];
1282}
1283
1284////////////////////////////////////////////////////////////////////////////////
1285/// Returns (and, if 0, creates) browsable objects for this branch
1286/// See TVirtualBranchBrowsable::FillListOfBrowsables.
1287
1289 if (fBrowsables) return fBrowsables;
1290 fBrowsables=new TList();
1292 return fBrowsables;
1293}
1294
1295////////////////////////////////////////////////////////////////////////////////
1296/// Return the name of the user class whose content is stored in this branch,
1297/// if any. If this branch was created using the 'leaflist' technique, this
1298/// function returns an empty string.
1299
1300const char * TBranch::GetClassName() const
1301{
1302 return "";
1303}
1304
1305////////////////////////////////////////////////////////////////////////////////
1306/// Return icon name depending on type of branch.
1307
1308const char* TBranch::GetIconName() const
1309{
1310 if (IsFolder())
1311 return "TBranchElement-folder";
1312 else
1313 return "TBranchElement-leaf";
1314}
1315
1316////////////////////////////////////////////////////////////////////////////////
1317/// A helper function to locate the correct basket - and its first entry.
1318/// Extracted to a common private function because it is needed by both GetEntry
1319/// and GetBulkEntries. It should not be called directly.
1320///
1321/// If a new basket must be constructed and the user_buffer is provided, then
1322/// the user_buffer will back the memory of the newly-constructed basket.
1323///
1324/// Assumes that this branch is enabled.
1326 TBuffer *user_buffer)
1327{
1328 Long64_t updatedNext = fNextBasketEntry;
1329 Long64_t entry = fReadEntry;
1330 if (R__likely(fCurrentBasket && fFirstBasketEntry <= entry && entry < fNextBasketEntry)) {
1331 // We have found the basket containing this entry.
1332 // make sure basket buffers are in memory.
1333 basket = fCurrentBasket;
1335 } else {
1336 if ((entry < fFirstEntry) || (entry >= fEntryNumber)) {
1337 return 0;
1338 }
1340 Long64_t last = fNextBasketEntry - 1;
1341 // Are we still in the same ReadBasket?
1342 if ((entry < first) || (entry > last)) {
1344 if (fReadBasket < 0) {
1345 fNextBasketEntry = -1;
1346 Error("GetBasketAndFirst", "In the branch %s, no basket contains the entry %lld\n", GetName(), entry);
1347 return -1;
1348 }
1349 if (fReadBasket == fWriteBasket) {
1351 } else {
1353 }
1354 updatedNext = fNextBasketEntry;
1356 }
1357 // We have found the basket containing this entry.
1358 // make sure basket buffers are in memory.
1360 if (!basket) {
1361 basket = GetBasketImpl(fReadBasket, user_buffer);
1362 if (!basket) {
1363 fCurrentBasket = 0;
1364 fFirstBasketEntry = -1;
1365 fNextBasketEntry = -1;
1366 return -1;
1367 }
1368 if (fTree->GetClusterPrefetch()) {
1369 TTree::TClusterIterator clusterIterator = fTree->GetClusterIterator(entry);
1370 clusterIterator.Next();
1371 Int_t nextClusterEntry = clusterIterator.GetNextEntry();
1372 for (Int_t i = fReadBasket + 1; i < fMaxBaskets && fBasketEntry[i] < nextClusterEntry; i++) {
1373 GetBasket(i);
1374 }
1375 }
1376 // Getting the next basket might reset the current one and
1377 // cause a reset of the first / next basket entries back to -1.
1379 fNextBasketEntry = updatedNext;
1380 }
1381 if (user_buffer) {
1382 // Disassociate basket from memory buffer for bulk IO
1383 // When the user provides a memory buffer (i.e., for bulk IO), we should
1384 // make sure to drop all references to that buffer in the TTree afterward.
1385 fCurrentBasket = nullptr;
1386 fBaskets[fReadBasket] = nullptr;
1387 } else {
1388 fCurrentBasket = basket;
1389 }
1390 }
1391 return 1;
1392}
1393
1394////////////////////////////////////////////////////////////////////////////////
1395/// Returns true if this branch supports bulk IO, false otherwise.
1396///
1397/// This will return true if all the various preconditions necessary hold true
1398/// to perform bulk IO (reasonable type, single TLeaf, etc); the bulk IO may
1399/// still fail, depending on the contents of the individual TBaskets loaded.
1401 return (fNleaves == 1) &&
1402 (static_cast<TLeaf*>(fLeaves.UncheckedAt(0))->GetDeserializeType() != TLeaf::DeserializeType::kDestructive);
1403}
1404
1405////////////////////////////////////////////////////////////////////////////////
1406/// Read as many events as possible into the given buffer, using zero-copy
1407/// mechanisms.
1408///
1409/// Returns -1 in case of a failure. On success, returns a (non-zero) number of
1410/// events of the type held by this branch currently in the buffer.
1411///
1412/// On success, the caller should be able to access the contents of buf as
1413///
1414/// static_cast<T*>(buf.GetCurrent())
1415///
1416/// where T is the type stored on this branch. The array's length is the return
1417/// value of this function.
1418///
1419/// NOTES:
1420/// - This interface is meant to be used by higher-level, type-safe wrappers, not
1421/// by end-users.
1422/// - This only returns events
1423///
1424
1426{
1427 // TODO: eventually support multiple leaves.
1428 if (R__unlikely(fNleaves != 1)) return -1;
1429 TLeaf *leaf = static_cast<TLeaf*>(fLeaves.UncheckedAt(0));
1431
1432 // Remember which entry we are reading.
1433 fReadEntry = entry;
1434
1435 Bool_t enabled = !TestBit(kDoNotProcess);
1436 if (R__unlikely(!enabled)) return -1;
1437 TBasket *basket = nullptr;
1439 Int_t result = GetBasketAndFirst(basket, first, &user_buf);
1440 if (R__unlikely(result <= 0)) return -1;
1441 // Only support reading from full clusters.
1442 if (R__unlikely(entry != first)) {
1443 //printf("Failed to read from full cluster; first entry is %ld; requested entry is %ld.\n", first, entry);
1444 return -1;
1445 }
1446
1447 basket->PrepareBasket(entry);
1448 TBuffer* buf = basket->GetBufferRef();
1449
1450 // Test for very old ROOT files.
1451 if (R__unlikely(!buf)) {
1452 Error("GetBulkEntries", "Failed to get a new buffer.\n");
1453 return -1;
1454 }
1455 // Test for displacements, which aren't supported in fast mode.
1456 if (R__unlikely(basket->GetDisplacement())) {
1457 Error("GetBulkEntries", "Basket has displacement.\n");
1458 return -1;
1459 }
1460
1461 Int_t bufbegin = basket->GetKeylen();
1462 buf->SetBufferOffset(bufbegin);
1463
1465 //printf("Requesting %d events; fNextBasketEntry=%lld; first=%lld.\n", N, fNextBasketEntry, first);
1466 if (R__unlikely(!leaf->ReadBasketFast(*buf, N))) {
1467 Error("GetBulkEntries", "Leaf failed to read.\n");
1468 return -1;
1469 }
1470 user_buf.SetBufferOffset(bufbegin);
1471
1472 fCurrentBasket = nullptr;
1473 fBaskets[fReadBasket] = nullptr;
1474 R__ASSERT(fExtraBasket == nullptr && "fExtraBasket should have been set to nullptr by GetFreshBasket");
1475 fExtraBasket = basket;
1476 basket->DisownBuffer();
1477
1478 return N;
1479}
1480
1481// TODO: Template this and the call above; only difference is the TLeaf function (ReadBasketFast vs
1482// ReadBasketSerialized
1484{
1485 // TODO: eventually support multiple leaves.
1486 if (R__unlikely(fNleaves != 1)) { return -1; }
1487 TLeaf *leaf = static_cast<TLeaf*>(fLeaves.UncheckedAt(0));
1489 Error("GetEntriesSerialized", "Encountered a branch with destructive deserialization; failing.\n");
1490 return -1;
1491 }
1492
1493 // Remember which entry we are reading.
1494 fReadEntry = entry;
1495
1496 Bool_t enabled = !TestBit(kDoNotProcess);
1497 if (R__unlikely(!enabled)) { return -1; }
1498 TBasket *basket = nullptr;
1500 Int_t result = GetBasketAndFirst(basket, first, &user_buf);
1501 if (R__unlikely(result <= 0)) { return -1; }
1502 // Only support reading from full clusters.
1503 if (R__unlikely(entry != first)) {
1504 Error("GetEntriesSerialized", "Failed to read from full cluster; first entry is %lld; requested entry is %lld.\n", first, entry);
1505 return -1;
1506 }
1507
1508 basket->PrepareBasket(entry);
1509 TBuffer* buf = basket->GetBufferRef();
1510
1511 // Test for very old ROOT files.
1512 if (R__unlikely(!buf)) {
1513 Error("GetEntriesSerialized", "Failed to get a new buffer.\n");
1514 return -1;
1515 }
1516 // Test for displacements, which aren't supported in fast mode.
1517 if (R__unlikely(basket->GetDisplacement())) {
1518 Error("GetEntriesSerialized", "Basket has displacement.\n");
1519 return -1;
1520 }
1521
1522 Int_t bufbegin = basket->GetKeylen();
1523 buf->SetBufferOffset(bufbegin);
1524
1526 //Info("GetEntriesSerialized", "Requesting %d events; fNextBasketEntry=%lld; first=%lld.\n", N, fNextBasketEntry, first);
1527
1528 if (R__unlikely(!leaf->ReadBasketSerialized(*buf, N))) {
1529 Error("GetEntriesSerialized", "Leaf failed to read.\n");
1530 return -1;
1531 }
1532 user_buf.SetBufferOffset(bufbegin);
1533
1534 if (count_buf) {
1535 TLeaf *count_leaf = leaf->GetLeafCount();
1536 if (count_leaf) {
1537 //printf("Getting leaf count entries.\n");
1538 TBranch *count_branch = count_leaf->GetBranch();
1539 if (R__unlikely(count_branch->GetEntriesSerialized(entry, *count_buf) < 0)) {
1540 Error("GetEntriesSerialized", "Failed to read count leaf.\n");
1541 return -1;
1542 }
1543 } else {
1544 // TODO: if you ask for a count on a fixed-size branch, maybe we should
1545 // just fail?
1546 Int_t entry_count_serialized;
1547 char *tmp_ptr = reinterpret_cast<char*>(&entry_count_serialized);
1548 tobuf(tmp_ptr, leaf->GetLenType() * leaf->GetNdata());
1549 Int_t cur_offset = count_buf->GetCurrent() - count_buf->Buffer();
1550 for (int idx=0; idx<N; idx++) {
1551 *count_buf << entry_count_serialized;
1552 }
1553 count_buf->SetBufferOffset(cur_offset);
1554 }
1555 }
1556
1557 return N;
1558}
1559
1560////////////////////////////////////////////////////////////////////////////////
1561/// Read all leaves of entry and return total number of bytes read.
1562///
1563/// The input argument "entry" is the entry number in the current tree.
1564/// In case of a TChain, the entry number in the current Tree must be found
1565/// before calling this function. For example:
1566///
1567///~~~ {.cpp}
1568/// TChain* chain = ...;
1569/// Long64_t localEntry = chain->LoadTree(entry);
1570/// branch->GetEntry(localEntry);
1571///~~~
1572///
1573/// The function returns the number of bytes read from the input buffer.
1574/// If entry does not exist, the function returns 0.
1575/// If an I/O error occurs, the function returns -1.
1576///
1577/// See IMPORTANT REMARKS in TTree::GetEntry.
1578
1580{
1581 // Remember which entry we are reading.
1582 fReadEntry = entry;
1583
1584 if (R__unlikely(TestBit(kDoNotProcess) && !getall)) { return 0; }
1585
1586 TBasket *basket; // will be initialized in the if/then clauses.
1588
1589 Int_t result = GetBasketAndFirst(basket, first, nullptr);
1590 if (R__unlikely(result <= 0)) { return result; }
1591
1592 basket->PrepareBasket(entry);
1593 TBuffer* buf = basket->GetBufferRef();
1594
1595 // This test necessary to read very old Root files (NvE).
1596 if (R__unlikely(!buf)) {
1597 TFile* file = GetFile(0);
1598 if (!file) return -1;
1600 buf = basket->GetBufferRef();
1601 }
1602
1603 // Set entry offset in buffer.
1605 buf->ResetMap();
1606 }
1607 if (R__unlikely(!buf->IsReading())) {
1608 basket->SetReadMode();
1609 }
1610
1611 Int_t* entryOffset = basket->GetEntryOffset();
1612 Int_t bufbegin = 0;
1613 if (entryOffset) {
1614 bufbegin = entryOffset[entry-first];
1615 buf->SetBufferOffset(bufbegin);
1616 Int_t* displacement = basket->GetDisplacement();
1617 if (R__unlikely(displacement)) {
1618 buf->SetBufferDisplacement(displacement[entry-first]);
1619 }
1620 } else {
1621 bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1622 buf->SetBufferOffset(bufbegin);
1623 }
1624
1625 // Int_t bufbegin = buf->Length();
1626 (this->*fReadLeaves)(*buf);
1627 return buf->Length() - bufbegin;
1628}
1629
1630////////////////////////////////////////////////////////////////////////////////
1631/// Read all leaves of an entry and export buffers to real objects in a TClonesArray list.
1632///
1633/// Returns total number of bytes read.
1634
1636{
1637 // Remember which entry we are reading.
1638 fReadEntry = entry;
1639
1640 if (TestBit(kDoNotProcess)) {
1641 return 0;
1642 }
1643 if ((entry < 0) || (entry >= fEntryNumber)) {
1644 return 0;
1645 }
1646 Int_t nbytes = 0;
1648 Long64_t last = fNextBasketEntry - 1;
1649 // Are we still in the same ReadBasket?
1650 if ((entry < first) || (entry > last)) {
1652 if (fReadBasket < 0) {
1653 fNextBasketEntry = -1;
1654 Error("In the branch %s, no basket contains the entry %d\n", GetName(), entry);
1655 return -1;
1656 }
1657 if (fReadBasket == fWriteBasket) {
1659 } else {
1661 }
1663 }
1664
1665 // We have found the basket containing this entry.
1666 // Make sure basket buffers are in memory.
1667 TBasket* basket = GetBasketImpl(fReadBasket, nullptr);
1668 fCurrentBasket = basket;
1669 if (!basket) {
1670 fFirstBasketEntry = -1;
1671 fNextBasketEntry = -1;
1672 return 0;
1673 }
1674 TBuffer* buf = basket->GetBufferRef();
1675 // Set entry offset in buffer and read data from all leaves.
1677 buf->ResetMap();
1678 }
1679 if (R__unlikely(!buf->IsReading())) {
1680 basket->SetReadMode();
1681 }
1682 Int_t* entryOffset = basket->GetEntryOffset();
1683 Int_t bufbegin = 0;
1684 if (entryOffset) {
1685 bufbegin = entryOffset[entry-first];
1686 buf->SetBufferOffset(bufbegin);
1687 Int_t* displacement = basket->GetDisplacement();
1688 if (R__unlikely(displacement)) {
1689 buf->SetBufferDisplacement(displacement[entry-first]);
1690 }
1691 } else {
1692 bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1693 buf->SetBufferOffset(bufbegin);
1694 }
1695 TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(0);
1696 leaf->ReadBasketExport(*buf, li, nentries);
1697 nbytes = buf->Length() - bufbegin;
1698 return nbytes;
1699}
1700
1701////////////////////////////////////////////////////////////////////////////////
1702/// Fill expectedClass and expectedType with information on the data type of the
1703/// object/values contained in this branch (and thus the type of pointers
1704/// expected to be passed to Set[Branch]Address
1705/// return 0 in case of success and > 0 in case of failure.
1706
1707Int_t TBranch::GetExpectedType(TClass *&expectedClass,EDataType &expectedType)
1708{
1709 expectedClass = 0;
1710 expectedType = kOther_t;
1711 TLeaf* l = (TLeaf*) GetListOfLeaves()->At(0);
1712 if (l) {
1713 expectedType = (EDataType) gROOT->GetType(l->GetTypeName())->GetType();
1714 return 0;
1715 } else {
1716 Error("GetExpectedType", "Did not find any leaves in %s",GetName());
1717 return 1;
1718 }
1719}
1720
1721////////////////////////////////////////////////////////////////////////////////
1722/// Return pointer to the file where branch buffers reside, returns 0
1723/// in case branch buffers reside in the same file as tree header.
1724/// If mode is 1 the branch buffer file is recreated.
1725
1727{
1728 if (fDirectory) return fDirectory->GetFile();
1729
1730 // check if a file with this name is in the list of Root files
1731 TFile *file = 0;
1732 {
1734 file = (TFile*)gROOT->GetListOfFiles()->FindObject(fFileName.Data());
1735 if (file) {
1736 fDirectory = file;
1737 return file;
1738 }
1739 }
1740
1741 if (fFileName.Length() == 0) return 0;
1742
1743 TString bFileName( GetRealFileName() );
1744
1745 // Open file (new file if mode = 1)
1746 {
1748 if (mode) file = TFile::Open(bFileName, "recreate");
1749 else file = TFile::Open(bFileName);
1750 }
1751 if (!file) return 0;
1752 if (file->IsZombie()) {delete file; return 0;}
1754 return file;
1755}
1756
1757////////////////////////////////////////////////////////////////////////////////
1758/// Return a fresh basket by either resusing an existing basket that needs
1759/// to be drop (according to TTree::MemoryFull) or create a new one.
1760///
1761/// If the user_buffer argument is non-null, then the memory in the
1762/// user-provided buffer will be utilized by the underlying basket.
1763///
1764/// The basket number is used to estimate the required buffer size
1765/// and try to optimize memory usage and number of memory allocation.
1766
1767TBasket* TBranch::GetFreshBasket(Int_t basketnumber, TBuffer* user_buffer)
1768{
1769 TBasket *basket = 0;
1770 if (user_buffer && fExtraBasket) {
1771 basket = fExtraBasket;
1772 fExtraBasket = nullptr;
1773 basket->AdoptBuffer(user_buffer);
1774 } else {
1775 if (GetTree()->MemoryFull(0)) {
1776 if (fNBaskets==1) {
1777 // Steal the existing basket
1778 Int_t oldindex = fBaskets.GetLast();
1779 basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1780 if (!basket) {
1781 fBaskets.SetLast(-2); // For recalculation of Last.
1782 oldindex = fBaskets.GetLast();
1783 if (oldindex != fBaskets.LowerBound()-1) {
1784 basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1785 }
1786 }
1787 if (basket && fBasketBytes[oldindex]!=0) {
1788 if (basket == fCurrentBasket) {
1789 fCurrentBasket = 0;
1790 fFirstBasketEntry = -1;
1791 fNextBasketEntry = -1;
1792 }
1793 fBaskets.AddAt(0,oldindex);
1794 fBaskets.SetLast(-1);
1795 fNBaskets = 0;
1796 basket->ReadResetBuffer(basketnumber);
1797#ifdef R__TRACK_BASKET_ALLOC_TIME
1798 fTree->AddAllocationTime(basket->GetResetAllocationTime());
1799#endif
1801 } else {
1802 basket = fTree->CreateBasket(this);
1803 }
1804 } else if (fNBaskets == 0) {
1805 // There is nothing to drop!
1806 basket = fTree->CreateBasket(this);
1807 } else {
1808 // Memory is full and there is more than one basket,
1809 // Let DropBaskets do it job.
1810 DropBaskets();
1811 basket = fTree->CreateBasket(this);
1812 }
1813 } else {
1814 basket = fTree->CreateBasket(this);
1815 }
1816 if (user_buffer)
1817 basket->AdoptBuffer(user_buffer);
1818 }
1819 return basket;
1820}
1821
1822////////////////////////////////////////////////////////////////////////////////
1823/// Drops the cluster two behind the current cluster and returns a fresh basket
1824/// by either reusing or creating a new one
1825
1827{
1828 TBasket *basket = 0;
1829
1830 // If GetClusterIterator is called with a negative entry then GetStartEntry will be 0
1831 // So we need to check if we reach the zero before we have gone back (1-VirtualSize) clusters
1832 // if this is the case, we want to keep everything in memory so we return a new basket
1834 if (iter.GetStartEntry() == 0) {
1835 return fTree->CreateBasket(this);
1836 }
1837
1838 // Iterate backwards (1-VirtualSize) clusters to reach cluster to be unloaded from memory,
1839 // skipped if VirtualSize > 0.
1840 for (Int_t j = 0; j < -fTree->GetMaxVirtualSize(); j++) {
1841 if (iter.Previous() == 0) {
1842 return fTree->CreateBasket(this);
1843 }
1844 }
1845
1846 Int_t entryToUnload = iter.Previous();
1847 // Finds the basket to unload from memory. Since the basket should be close to current
1848 // basket, just iterate backwards until the correct basket is reached. This should
1849 // be fast as long as the number of baskets per cluster is small
1850 Int_t basketToUnload = fReadBasket;
1851 while (fBasketEntry[basketToUnload] != entryToUnload) {
1852 basketToUnload--;
1853 if (basketToUnload < 0) {
1854 return fTree->CreateBasket(this);
1855 }
1856 }
1857
1858 // Retrieves the basket that is going to be unloaded from memory. If the basket did not
1859 // exist, create a new one
1860 basket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1861 if (basket) {
1862 fBaskets.AddAt(0, basketToUnload);
1863 --fNBaskets;
1864 } else {
1865 basket = fTree->CreateBasket(this);
1866 }
1867 ++basketToUnload;
1868
1869 // Clear the rest of the baskets. While it would be ideal to reuse these baskets
1870 // for other baskets in the new cluster. It would require the function to go
1871 // beyond its current scope. In the ideal case when each cluster only has 1 basket
1872 // this will perform well
1873 iter.Next();
1874 while (fBasketEntry[basketToUnload] < iter.GetStartEntry()) {
1875 TBasket *oldbasket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1876 if (oldbasket) {
1877 oldbasket->DropBuffers();
1878 delete oldbasket;
1879 fBaskets.AddAt(0, basketToUnload);
1880 --fNBaskets;
1881 }
1882 ++basketToUnload;
1883 }
1884 fBaskets.SetLast(-1);
1885 return basket;
1886}
1887
1888////////////////////////////////////////////////////////////////////////////////
1889/// Return the 'full' name of the branch. In particular prefix the mother's name
1890/// when it does not end in a trailing dot and thus is not part of the branch name
1892{
1893 TBranch* mother = GetMother();
1894 if (!mother || mother==this) {
1895 return fName;
1896 }
1897 TString motherName(mother->GetName());
1898 if (motherName.Length() && (motherName[motherName.Length()-1] == '.')) {
1899 return fName;
1900 }
1901 return motherName + "." + fName;
1902}
1903
1904////////////////////////////////////////////////////////////////////////////////
1905/// Return pointer to the 1st Leaf named name in thisBranch
1906
1907TLeaf* TBranch::GetLeaf(const char* name) const
1908{
1909 Int_t i;
1910 for (i=0;i<fNleaves;i++) {
1911 TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
1912 if (!strcmp(leaf->GetName(),name)) return leaf;
1913 }
1914 return 0;
1915}
1916
1917////////////////////////////////////////////////////////////////////////////////
1918/// Get real file name
1919
1921{
1922 if (fFileName.Length()==0) {
1923 return fFileName;
1924 }
1925 TString bFileName = fFileName;
1926
1927 // check if branch file name is absolute or a URL (e.g. root://host/...)
1928 char *bname = gSystem->ExpandPathName(fFileName.Data());
1929 if (!gSystem->IsAbsoluteFileName(bname) && !strstr(bname, ":/") && fTree && fTree->GetCurrentFile()) {
1930
1931 // if not, get filename where tree header is stored
1932 const char *tfn = fTree->GetCurrentFile()->GetName();
1933
1934 // If it is an archive file we need a special treatment
1935 TUrl arc(tfn);
1936 if (strlen(arc.GetAnchor()) > 0) {
1938 bFileName = arc.GetUrl();
1939 } else {
1940 // if this is an absolute path or a URL then prepend this path
1941 // to the branch file name
1942 char *tname = gSystem->ExpandPathName(tfn);
1943 if (gSystem->IsAbsoluteFileName(tname) || strstr(tname, ":/")) {
1944 bFileName = gSystem->DirName(tname);
1945 bFileName += "/";
1946 bFileName += fFileName;
1947 }
1948 delete [] tname;
1949 }
1950 }
1951 delete [] bname;
1952
1953 return bFileName;
1954}
1955
1956////////////////////////////////////////////////////////////////////////////////
1957/// Return all elements of one row unpacked in internal array fValues
1958/// [Actually just returns 1 (?)]
1959
1961{
1962 return 1;
1963}
1964
1965////////////////////////////////////////////////////////////////////////////////
1966/// Return whether this branch is in a mode where the object are decomposed
1967/// or not (Also known as MakeClass mode).
1968
1970{
1971 // Regular TBranch and TBrancObject can not be in makeClass mode
1972
1973 return kFALSE;
1974}
1975
1976////////////////////////////////////////////////////////////////////////////////
1977/// Get our top-level parent branch in the tree.
1978
1980{
1981 if (fMother) return fMother;
1982
1983 const TObjArray* array = fTree->GetListOfBranches();
1984 Int_t n = array->GetEntriesFast();
1985 for (Int_t i = 0; i < n; ++i) {
1986 TBranch* branch = (TBranch*) array->UncheckedAt(i);
1987 TBranch* parent = branch->GetSubBranch(this);
1988 if (parent) {
1989 const_cast<TBranch*>(this)->fMother = branch; // We can not yet use the 'mutable' keyword
1990 return branch;
1991 }
1992 }
1993 return 0;
1994}
1995
1996////////////////////////////////////////////////////////////////////////////////
1997/// Find the parent branch of child.
1998/// Return 0 if child is not in this branch hierarchy.
1999
2001{
2002 // Handle error condition, if the parameter is us, we cannot find the parent.
2003 if (this == child) {
2004 // Note: We cast away any const-ness of "this".
2005 return (TBranch*) this;
2006 }
2007
2008 if (child->fParent) {
2009 return child->fParent;
2010 }
2011
2013 for (Int_t i = 0; i < len; ++i) {
2014 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2015 if (!branch) {
2016 continue;
2017 }
2018 if (branch == child) {
2019 // We are the direct parent of child.
2020 const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
2021 // Note: We cast away any const-ness of "this".
2022 const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
2023 return (TBranch*) this;
2024 }
2025 // FIXME: This is a tail-recursion!
2026 TBranch* parent = branch->GetSubBranch(child);
2027 if (parent) {
2028 return parent;
2029 }
2030 }
2031 // We failed to find the parent.
2032 return 0;
2033}
2034
2035////////////////////////////////////////////////////////////////////////////////
2036/// Return total number of bytes in the branch (including current buffer)
2037
2039{
2041 // This intentionally only store the TBranch part and thus slightly
2042 // under-estimate the space used.
2043 // Since the TBranchElement part contains pointers to other branches (branch count),
2044 // doing regular Streaming would end up including those and thus greatly over-estimate
2045 // the size used.
2046 const_cast<TBranch *>(this)->TBranch::Streamer(b);
2047
2048 Long64_t totbytes = 0;
2049 if (fZipBytes > 0) totbytes = fTotBytes;
2050 return totbytes + b.Length();
2051}
2052
2053////////////////////////////////////////////////////////////////////////////////
2054/// Return total number of bytes in the branch (excluding current buffer)
2055/// if option ="*" includes all sub-branches of this branch too
2056
2058{
2059 Long64_t totbytes = fTotBytes;
2060 if (!option) return totbytes;
2061 if (option[0] != '*') return totbytes;
2062 //scan sub-branches
2064 for (Int_t i = 0; i < len; ++i) {
2065 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2066 if (branch) totbytes += branch->GetTotBytes(option);
2067 }
2068 return totbytes;
2069}
2070
2071////////////////////////////////////////////////////////////////////////////////
2072/// Return total number of zip bytes in the branch
2073/// if option ="*" includes all sub-branches of this branch too
2074
2076{
2077 Long64_t zipbytes = fZipBytes;
2078 if (!option) return zipbytes;
2079 if (option[0] != '*') return zipbytes;
2080 //scan sub-branches
2082 for (Int_t i = 0; i < len; ++i) {
2083 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2084 if (branch) zipbytes += branch->GetZipBytes(option);
2085 }
2086 return zipbytes;
2087}
2088
2089////////////////////////////////////////////////////////////////////////////////
2090/// Returns the IO settings currently in use for this branch.
2091
2093{
2094 return fIOFeatures;
2095}
2096
2097////////////////////////////////////////////////////////////////////////////////
2098/// Return kTRUE if an existing object in a TBranchObject must be deleted.
2099
2101{
2102 return TestBit(kAutoDelete);
2103}
2104
2105////////////////////////////////////////////////////////////////////////////////
2106/// Return kTRUE if more than one leaf or browsables, kFALSE otherwise.
2107
2109{
2110 if (fNleaves > 1) {
2111 return kTRUE;
2112 }
2113 TList* browsables = const_cast<TBranch*>(this)->GetBrowsables();
2114 return browsables && browsables->GetSize();
2115}
2116
2117////////////////////////////////////////////////////////////////////////////////
2118/// keep a maximum of fMaxEntries in memory
2119
2121{
2122 Int_t dentries = (Int_t) (fEntries - maxEntries);
2123 TBasket* basket = (TBasket*) fBaskets.UncheckedAt(0);
2124 if (basket) basket->MoveEntries(dentries);
2125 fEntries = maxEntries;
2126 fEntryNumber = maxEntries;
2127 //loop on sub branches
2129 for (Int_t i = 0; i < nb; ++i) {
2130 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2131 branch->KeepCircular(maxEntries);
2132 }
2133}
2134
2135////////////////////////////////////////////////////////////////////////////////
2136/// Baskets associated to this branch are forced to be in memory.
2137/// You can call TTree::SetMaxVirtualSize(maxmemory) to instruct
2138/// the system that the total size of the imported baskets does not
2139/// exceed maxmemory bytes.
2140///
2141/// The function returns the number of baskets that have been put in memory.
2142/// This method may be called to force all baskets of one or more branches
2143/// in memory when random access to entries in this branch is required.
2144/// See also TTree::LoadBaskets to load all baskets of all branches in memory.
2145
2147{
2148 Int_t nimported = 0;
2149 Int_t nbaskets = fWriteBasket;
2150 TFile *file = GetFile(0);
2151 if (!file) return 0;
2152 TBasket *basket;
2153 for (Int_t i=0;i<nbaskets;i++) {
2154 basket = (TBasket*)fBaskets.UncheckedAt(i);
2155 if (basket) continue;
2156 basket = GetFreshBasket(i, nullptr);
2157 if (fBasketBytes[i] == 0) {
2159 }
2160 Int_t badread = basket->ReadBasketBuffers(fBasketSeek[i],fBasketBytes[i],file);
2161 if (badread) {
2162 Error("Loadbaskets","Error while reading basket buffer %d of branch %s",i,GetName());
2163 return -1;
2164 }
2165 ++fNBaskets;
2166 fBaskets.AddAt(basket,i);
2167 nimported++;
2168 }
2169 return nimported;
2170}
2171
2172////////////////////////////////////////////////////////////////////////////////
2173/// Print TBranch parameters
2174///
2175/// If options contains "basketsInfo" print the entry number, location and size
2176/// of each baskets.
2177
2178void TBranch::Print(Option_t *option) const
2179{
2180 const int kLINEND = 77;
2181 Float_t cx = 1;
2182
2183 TString titleContent(GetTitle());
2184 if ( titleContent == GetName() ) {
2185 titleContent.Clear();
2186 }
2187
2188 if (fLeaves.GetEntries() == 1) {
2189 if (titleContent.Length()>=2 && titleContent[titleContent.Length()-2]=='/' && isalpha(titleContent[titleContent.Length()-1])) {
2190 // The type is already encoded. Nothing to do.
2191 } else {
2192 TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(0);
2193 if (titleContent.Length()) {
2194 titleContent.Prepend(" ");
2195 }
2196 // titleContent.Append("type: ");
2197 titleContent.Prepend(leaf->GetTypeName());
2198 }
2199 }
2200 Int_t titleLength = titleContent.Length();
2201
2202 Int_t aLength = titleLength + strlen(GetName());
2203 aLength += (aLength / 54 + 1) * 80 + 100;
2204 if (aLength < 200) aLength = 200;
2205 char *bline = new char[aLength];
2206
2207 Long64_t totBytes = GetTotalSize();
2208 if (fZipBytes) cx = (fTotBytes+0.00001)/fZipBytes;
2209 if (titleLength) snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName(),titleContent.Data());
2210 else snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName()," ");
2211 if (strlen(bline) > UInt_t(kLINEND)) {
2212 char *tmp = new char[strlen(bline)+1];
2213 if (titleLength) strlcpy(tmp, titleContent.Data(),strlen(bline)+1);
2214 snprintf(bline,aLength,"*Br%5d :%-9s : ",fgCount,GetName());
2215 int pos = strlen (bline);
2216 int npos = pos;
2217 int beg=0, end;
2218 while (beg < titleLength) {
2219 for (end=beg+1; end < titleLength-1; end ++)
2220 if (tmp[end] == ':') break;
2221 if (npos + end-beg+1 >= 78) {
2222 while (npos < kLINEND) {
2223 bline[pos ++] = ' ';
2224 npos ++;
2225 }
2226 bline[pos ++] = '*';
2227 bline[pos ++] = '\n';
2228 bline[pos ++] = '*';
2229 npos = 1;
2230 for (; npos < 12; npos ++)
2231 bline[pos ++] = ' ';
2232 bline[pos-2] = '|';
2233 }
2234 for (int n = beg; n <= end; n ++)
2235 bline[pos+n-beg] = tmp[n];
2236 pos += end-beg+1;
2237 npos += end-beg+1;
2238 beg = end+1;
2239 }
2240 while (npos < kLINEND) {
2241 bline[pos ++] = ' ';
2242 npos ++;
2243 }
2244 bline[pos ++] = '*';
2245 bline[pos] = '\0';
2246 delete[] tmp;
2247 }
2248 Printf("%s", bline);
2249
2250 if (fTotBytes > 2000000000) {
2251 Printf("*Entries :%lld : Total Size=%11lld bytes File Size = %lld *",fEntries,totBytes,fZipBytes);
2252 } else {
2253 if (fZipBytes > 0) {
2254 Printf("*Entries :%9lld : Total Size=%11lld bytes File Size = %10lld *",fEntries,totBytes,fZipBytes);
2255 } else {
2256 if (fWriteBasket > 0) {
2257 Printf("*Entries :%9lld : Total Size=%11lld bytes All baskets in memory *",fEntries,totBytes);
2258 } else {
2259 Printf("*Entries :%9lld : Total Size=%11lld bytes One basket in memory *",fEntries,totBytes);
2260 }
2261 }
2262 }
2263 Printf("*Baskets :%9d : Basket Size=%11d bytes Compression= %6.2f *",fWriteBasket,fBasketSize,cx);
2264
2265 if (strncmp(option,"basketsInfo",strlen("basketsInfo"))==0) {
2266 Int_t nbaskets = fWriteBasket;
2267 for (Int_t i=0;i<nbaskets;i++) {
2268 Printf("*Basket #%4d entry=%6lld pos=%6lld size=%5d",
2269 i, fBasketEntry[i], fBasketSeek[i], fBasketBytes[i]);
2270 }
2271 }
2272
2273 Printf("*............................................................................*");
2274 delete [] bline;
2275 fgCount++;
2276}
2277
2278////////////////////////////////////////////////////////////////////////////////
2279/// Print the information we have about which basket is currently cached and
2280/// whether they have been 'used'/'read' from the cache.
2281
2283{
2285}
2286
2287////////////////////////////////////////////////////////////////////////////////
2288/// Loop on all leaves of this branch to read Basket buffer.
2289
2291{
2292 // fLeaves->ReadBasket(basket);
2293}
2294
2295////////////////////////////////////////////////////////////////////////////////
2296/// Loop on all leaves of this branch to read Basket buffer.
2297
2299{
2300 for (Int_t i = 0; i < fNleaves; ++i) {
2301 TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2302 leaf->ReadBasket(b);
2303 }
2304}
2305
2306////////////////////////////////////////////////////////////////////////////////
2307/// Read zero leaves without the overhead of a loop.
2308
2310{
2311}
2312
2313////////////////////////////////////////////////////////////////////////////////
2314/// Read one leaf without the overhead of a loop.
2315
2317{
2319}
2320
2321////////////////////////////////////////////////////////////////////////////////
2322/// Read two leaves without the overhead of a loop.
2323
2325{
2328}
2329
2330////////////////////////////////////////////////////////////////////////////////
2331/// Loop on all leaves of this branch to fill Basket buffer.
2332
2334{
2335 for (Int_t i = 0; i < fNleaves; ++i) {
2336 TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2337 leaf->FillBasket(b);
2338 }
2339}
2340
2341////////////////////////////////////////////////////////////////////////////////
2342/// Refresh this branch using new information in b
2343/// This function is called by TTree::Refresh
2344
2346{
2347 if (b==0) return;
2348
2349 fEntryOffsetLen = b->fEntryOffsetLen;
2350 fWriteBasket = b->fWriteBasket;
2351 fEntryNumber = b->fEntryNumber;
2352 fMaxBaskets = b->fMaxBaskets;
2353 fEntries = b->fEntries;
2354 fTotBytes = b->fTotBytes;
2355 fZipBytes = b->fZipBytes;
2356 fReadBasket = 0;
2357 fReadEntry = -1;
2358 fFirstBasketEntry = -1;
2359 fNextBasketEntry = -1;
2360 fCurrentBasket = 0;
2361 delete [] fBasketBytes;
2362 delete [] fBasketEntry;
2363 delete [] fBasketSeek;
2367 Int_t i;
2368 for (i=0;i<fMaxBaskets;i++) {
2369 fBasketBytes[i] = b->fBasketBytes[i];
2370 fBasketEntry[i] = b->fBasketEntry[i];
2371 fBasketSeek[i] = b->fBasketSeek[i];
2372 }
2373 fBaskets.Delete();
2374 Int_t nbaskets = b->fBaskets.GetSize();
2375 fBaskets.Expand(nbaskets);
2376 // If the current fWritebasket is in memory, take it (just swap)
2377 // from the Tree being read
2378 TBasket *basket = (TBasket*)b->fBaskets.UncheckedAt(fWriteBasket);
2379 fBaskets.AddAt(basket,fWriteBasket);
2380 if (basket) {
2381 fNBaskets = 1;
2382 --(b->fNBaskets);
2383 b->fBaskets.RemoveAt(fWriteBasket);
2384 basket->SetBranch(this);
2385 }
2386}
2387
2388////////////////////////////////////////////////////////////////////////////////
2389/// Reset a Branch.
2390///
2391/// - Existing buffers are deleted.
2392/// - Entries, max and min are reset.
2393
2395{
2396 fReadBasket = 0;
2397 fReadEntry = -1;
2398 fFirstBasketEntry = -1;
2399 fNextBasketEntry = -1;
2400 fCurrentBasket = 0;
2401 fWriteBasket = 0;
2402 fEntries = 0;
2403 fTotBytes = 0;
2404 fZipBytes = 0;
2405 fEntryNumber = 0;
2406
2407 if (fBasketBytes) {
2408 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2409 fBasketBytes[i] = 0;
2410 }
2411 }
2412
2413 if (fBasketEntry) {
2414 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2415 fBasketEntry[i] = 0;
2416 }
2417 }
2418
2419 if (fBasketSeek) {
2420 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2421 fBasketSeek[i] = 0;
2422 }
2423 }
2424
2425 fBaskets.Delete();
2426 fNBaskets = 0;
2427}
2428
2429////////////////////////////////////////////////////////////////////////////////
2430/// Reset a Branch.
2431///
2432/// - Existing buffers are deleted.
2433/// - Entries, max and min are reset.
2434
2436{
2437 fReadBasket = 0;
2438 fReadEntry = -1;
2439 fFirstBasketEntry = -1;
2440 fNextBasketEntry = -1;
2441 fCurrentBasket = 0;
2442 fWriteBasket = 0;
2443 fEntries = 0;
2444 fTotBytes = 0;
2445 fZipBytes = 0;
2446 fEntryNumber = 0;
2447
2448 if (fBasketBytes) {
2449 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2450 fBasketBytes[i] = 0;
2451 }
2452 }
2453
2454 if (fBasketEntry) {
2455 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2456 fBasketEntry[i] = 0;
2457 }
2458 }
2459
2460 if (fBasketSeek) {
2461 for (Int_t i = 0; i < fMaxBaskets; ++i) {
2462 fBasketSeek[i] = 0;
2463 }
2464 }
2465
2466 TBasket *reusebasket = (TBasket*)fBaskets[fWriteBasket];
2467 if (reusebasket) {
2469 } else {
2470 reusebasket = (TBasket*)fBaskets[fReadBasket];
2471 if (reusebasket) {
2472 fBaskets[fReadBasket] = 0;
2473 }
2474 }
2475 fBaskets.Delete();
2476 if (reusebasket) {
2477 fNBaskets = 1;
2478 reusebasket->WriteReset();
2479 fBaskets[0] = reusebasket;
2480 } else {
2481 fNBaskets = 0;
2482 }
2483}
2484
2485////////////////////////////////////////////////////////////////////////////////
2486/// Reset the address of the branch.
2487
2489{
2490 fAddress = 0;
2491
2492 // Reset last read entry number, we have will had new user object now.
2493 fReadEntry = -1;
2494
2495 for (Int_t i = 0; i < fNleaves; ++i) {
2496 TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2497 leaf->SetAddress(0);
2498 }
2499
2500 Int_t nbranches = fBranches.GetEntriesFast();
2501 for (Int_t i = 0; i < nbranches; ++i) {
2502 TBranch* abranch = (TBranch*) fBranches[i];
2503 // FIXME: This is a tail recursion.
2504 abranch->ResetAddress();
2505 }
2506}
2507
2508////////////////////////////////////////////////////////////////////////////////
2509/// Static function resetting fgCount
2510
2512{
2513 fgCount = 0;
2514}
2515
2516////////////////////////////////////////////////////////////////////////////////
2517/// Set address of this branch.
2518
2519void TBranch::SetAddress(void* addr)
2520{
2521 if (TestBit(kDoNotProcess)) {
2522 return;
2523 }
2524 fReadEntry = -1;
2525 fFirstBasketEntry = -1;
2526 fNextBasketEntry = -1;
2527 fAddress = (char*) addr;
2528 for (Int_t i = 0; i < fNleaves; ++i) {
2529 TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2530 Int_t offset = leaf->GetOffset();
2531 if (TestBit(kIsClone)) {
2532 offset = 0;
2533 }
2534 if (fAddress) leaf->SetAddress(fAddress + offset);
2535 else leaf->SetAddress(0);
2536 }
2537}
2538
2539////////////////////////////////////////////////////////////////////////////////
2540/// Set the automatic delete bit.
2541///
2542/// This bit is used by TBranchObject::ReadBasket to decide if an object
2543/// referenced by a TBranchObject must be deleted or not before reading
2544/// a new entry.
2545///
2546/// If autodel is kTRUE, this existing object will be deleted, a new object
2547/// created by the default constructor, then read from disk by the streamer.
2548///
2549/// If autodel is kFALSE, the existing object is not deleted. Root assumes
2550/// that the user is taking care of deleting any internal object or array
2551/// (this can be done in the streamer).
2552
2554{
2555 if (autodel) {
2556 SetBit(kAutoDelete, 1);
2557 } else {
2558 SetBit(kAutoDelete, 0);
2559 }
2560}
2561
2562////////////////////////////////////////////////////////////////////////////////
2563/// Set the basket size
2564/// The function makes sure that the basket size is greater than fEntryOffsetlen
2565
2567{
2568 Int_t minsize = 100 + fName.Length();
2569 if (buffsize < minsize+fEntryOffsetLen) buffsize = minsize+fEntryOffsetLen;
2570 fBasketSize = buffsize;
2571 TBasket *basket = (TBasket*)fBaskets[fWriteBasket];
2572 if (basket) {
2573 basket->AdjustSize(fBasketSize);
2574 }
2575}
2576
2577////////////////////////////////////////////////////////////////////////////////
2578/// Set address of this branch directly from a TBuffer to avoid streaming.
2579///
2580/// Note: We do not take ownership of the buffer.
2581
2583{
2584 // Check this is possible
2585 if ( (fNleaves != 1)
2586 || (strcmp("TLeafObject",fLeaves.UncheckedAt(0)->ClassName())!=0) ) {
2587 Error("TBranch::SetAddress","Filling from a TBuffer can only be done with a not split object branch. Request ignored.");
2588 } else {
2589 fReadEntry = -1;
2590 fNextBasketEntry = -1;
2591 fFirstBasketEntry = -1;
2592 // Note: We do not take ownership of the buffer.
2593 fEntryBuffer = buf;
2594 }
2595}
2596
2597////////////////////////////////////////////////////////////////////////////////
2598/// Set compression algorithm.
2599
2601{
2602 if (algorithm < 0 || algorithm >= ROOT::RCompressionSetting::EAlgorithm::kUndefined) algorithm = 0;
2603 if (fCompress < 0) {
2605 } else {
2606 int level = fCompress % 100;
2607 fCompress = 100 * algorithm + level;
2608 }
2609
2611 for (Int_t i=0;i<nb;i++) {
2612 TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2613 branch->SetCompressionAlgorithm(algorithm);
2614 }
2615}
2616
2617////////////////////////////////////////////////////////////////////////////////
2618/// Set compression level.
2619
2621{
2622 if (level < 0) level = 0;
2623 if (level > 99) level = 99;
2624 if (fCompress < 0) {
2625 fCompress = level;
2626 } else {
2627 int algorithm = fCompress / 100;
2628 if (algorithm >= ROOT::RCompressionSetting::EAlgorithm::kUndefined) algorithm = 0;
2629 fCompress = 100 * algorithm + level;
2630 }
2631
2633 for (Int_t i=0;i<nb;i++) {
2634 TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2635 branch->SetCompressionLevel(level);
2636 }
2637}
2638
2639////////////////////////////////////////////////////////////////////////////////
2640/// Set compression settings.
2641
2643{
2644 fCompress = settings;
2645
2647 for (Int_t i=0;i<nb;i++) {
2648 TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2649 branch->SetCompressionSettings(settings);
2650 }
2651}
2652
2653////////////////////////////////////////////////////////////////////////////////
2654/// Update the default value for the branch's fEntryOffsetLen if and only if
2655/// it was already non zero (and the new value is not zero)
2656/// If updateExisting is true, also update all the existing branches.
2657
2658void TBranch::SetEntryOffsetLen(Int_t newdefault, Bool_t updateExisting)
2659{
2660 if (fEntryOffsetLen && newdefault) {
2661 fEntryOffsetLen = newdefault;
2662 }
2663 if (updateExisting) {
2664 TIter next( GetListOfBranches() );
2665 TBranch *b;
2666 while ( ( b = (TBranch*)next() ) ) {
2667 b->SetEntryOffsetLen( newdefault, kTRUE );
2668 }
2669 }
2670}
2671
2672////////////////////////////////////////////////////////////////////////////////
2673/// Set the number of entries in this branch.
2674
2676{
2677 fEntries = entries;
2678 fEntryNumber = entries;
2679}
2680
2681////////////////////////////////////////////////////////////////////////////////
2682/// Set file where this branch writes/reads its buffers.
2683/// By default the branch buffers reside in the file where the
2684/// Tree was created.
2685/// If the file name where the tree was created is an absolute
2686/// path name or an URL (e.g. or root://host/...)
2687/// and if the fname is not an absolute path name or an URL then
2688/// the path of the tree file is prepended to fname to make the
2689/// branch file relative to the tree file. In this case one can
2690/// move the tree + all branch files to a different location in
2691/// the file system and still access the branch files.
2692/// The ROOT file will be connected only when necessary.
2693/// If called by TBranch::Fill (via TBasket::WriteFile), the file
2694/// will be created with the option "recreate".
2695/// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2696/// will be opened in read mode.
2697/// To open a file in "update" mode or with a certain compression
2698/// level, use TBranch::SetFile(TFile *file).
2699
2701{
2702 if (file == 0) file = fTree->GetCurrentFile();
2704 if (file == fTree->GetCurrentFile()) fFileName = "";
2705 else fFileName = file->GetName();
2706
2707 if (file && fCompress == -1) {
2708 fCompress = file->GetCompressionLevel();
2709 }
2710
2711 // Apply to all existing baskets.
2712 TIter nextb(GetListOfBaskets());
2713 TBasket *basket;
2714 while ((basket = (TBasket*)nextb())) {
2715 basket->SetParent(file);
2716 }
2717
2718 // Apply to sub-branches as well.
2719 TIter next(GetListOfBranches());
2720 TBranch *branch;
2721 while ((branch = (TBranch*)next())) {
2722 branch->SetFile(file);
2723 }
2724}
2725
2726////////////////////////////////////////////////////////////////////////////////
2727/// Set file where this branch writes/reads its buffers.
2728/// By default the branch buffers reside in the file where the
2729/// Tree was created.
2730/// If the file name where the tree was created is an absolute
2731/// path name or an URL (e.g. root://host/...)
2732/// and if the fname is not an absolute path name or an URL then
2733/// the path of the tree file is prepended to fname to make the
2734/// branch file relative to the tree file. In this case one can
2735/// move the tree + all branch files to a different location in
2736/// the file system and still access the branch files.
2737/// The ROOT file will be connected only when necessary.
2738/// If called by TBranch::Fill (via TBasket::WriteFile), the file
2739/// will be created with the option "recreate".
2740/// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2741/// will be opened in read mode.
2742/// To open a file in "update" mode or with a certain compression
2743/// level, use TBranch::SetFile(TFile *file).
2744
2745void TBranch::SetFile(const char* fname)
2746{
2747 fFileName = fname;
2748 fDirectory = 0;
2749
2750 //apply to sub-branches as well
2751 TIter next(GetListOfBranches());
2752 TBranch *branch;
2753 while ((branch = (TBranch*)next())) {
2754 branch->SetFile(fname);
2755 }
2756}
2757
2758////////////////////////////////////////////////////////////////////////////////
2759/// Set the branch in a mode where the object are decomposed
2760/// (Also known as MakeClass mode).
2761/// Return whether the setting was possible (it is not possible for
2762/// TBranch and TBranchObject).
2763
2765{
2766 // Regular TBranch and TBrancObject can not be in makeClass mode
2767 return kFALSE;
2768}
2769
2770////////////////////////////////////////////////////////////////////////////////
2771/// Set object this branch is pointing to.
2772
2773void TBranch::SetObject(void * /* obj */)
2774{
2775 if (TestBit(kDoNotProcess)) {
2776 return;
2777 }
2778 Warning("SetObject","is not supported in TBranch objects");
2779}
2780
2781////////////////////////////////////////////////////////////////////////////////
2782/// Set branch status to Process or DoNotProcess.
2783
2785{
2786 if (status) ResetBit(kDoNotProcess);
2787 else SetBit(kDoNotProcess);
2788}
2789
2790////////////////////////////////////////////////////////////////////////////////
2791/// Stream a class object
2792
2793void TBranch::Streamer(TBuffer& b)
2794{
2795 if (b.IsReading()) {
2796 UInt_t R__s, R__c;
2797 fTree = 0; // Will be set by TTree::Streamer
2798 fAddress = 0;
2799 gROOT->SetReadingObject(kTRUE);
2800
2801 // Reset transients.
2803 fCurrentBasket = 0;
2804 fFirstBasketEntry = -1;
2805 fNextBasketEntry = -1;
2806
2807 Version_t v = b.ReadVersion(&R__s, &R__c);
2808 if (v > 9) {
2809 b.ReadClassBuffer(TBranch::Class(), this, v, R__s, R__c);
2810
2811 if (fWriteBasket>=fBaskets.GetSize()) {
2813 }
2814 fDirectory = 0;
2816 for (Int_t i=0;i<fNleaves;i++) {
2817 TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2818 leaf->SetBranch(this);
2819 }
2820
2822 for (Int_t j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2824 if (bk) {
2825 bk->SetBranch(this);
2826 // GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2827 ++n;
2828 }
2829 }
2830 if (fWriteBasket >= fMaxBaskets) {
2831 //old versions may need this fix
2836
2837 }
2839 gROOT->SetReadingObject(kFALSE);
2840 if (IsA() == TBranch::Class()) {
2841 if (fNleaves == 0) {
2843 } else if (fNleaves == 1) {
2845 } else if (fNleaves == 2) {
2847 } else {
2849 }
2850 }
2851 return;
2852 }
2853 //====process old versions before automatic schema evolution
2854 Int_t n,i,j,ijunk;
2855 if (v > 5) {
2856 Stat_t djunk;
2857 TNamed::Streamer(b);
2858 if (v > 7) TAttFill::Streamer(b);
2859 b >> fCompress;
2860 b >> fBasketSize;
2861 b >> fEntryOffsetLen;
2862 b >> fWriteBasket;
2863 b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2864 b >> fOffset;
2865 b >> fMaxBaskets;
2866 if (v > 6) b >> fSplitLevel;
2867 b >> djunk; fEntries = (Long64_t)djunk;
2868 b >> djunk; fTotBytes = (Long64_t)djunk;
2869 b >> djunk; fZipBytes = (Long64_t)djunk;
2870
2871 fBranches.Streamer(b);
2872 fLeaves.Streamer(b);
2873 fBaskets.Streamer(b);
2877 Char_t isArray;
2878 b >> isArray;
2879 b.ReadFastArray(fBasketBytes,fMaxBaskets);
2880 b >> isArray;
2881 for (i=0;i<fMaxBaskets;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2882 b >> isArray;
2883 for (i=0;i<fMaxBaskets;i++) {
2884 if (isArray == 2) b >> fBasketSeek[i];
2885 else {Int_t bsize; b >> bsize; fBasketSeek[i] = (Long64_t)bsize;};
2886 }
2887 fFileName.Streamer(b);
2888 b.CheckByteCount(R__s, R__c, TBranch::IsA());
2889 fDirectory = 0;
2891 for (i=0;i<fNleaves;i++) {
2892 TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2893 leaf->SetBranch(this);
2894 }
2896 for (j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2898 if (bk) {
2899 bk->SetBranch(this);
2900 //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2901 ++n;
2902 }
2903 }
2904 if (fWriteBasket >= fMaxBaskets) {
2905 //old versions may need this fix
2910
2911 }
2912 // Check Byte Count is not needed since it was done in ReadBuffer
2914 gROOT->SetReadingObject(kFALSE);
2915 b.CheckByteCount(R__s, R__c, TBranch::IsA());
2916 if (IsA() == TBranch::Class()) {
2917 if (fNleaves == 0) {
2919 } else if (fNleaves == 1) {
2921 } else if (fNleaves == 2) {
2923 } else {
2925 }
2926 }
2927 return;
2928 }
2929 //====process very old versions
2930 Stat_t djunk;
2931 TNamed::Streamer(b);
2932 b >> fCompress;
2933 b >> fBasketSize;
2934 b >> fEntryOffsetLen;
2935 b >> fMaxBaskets;
2936 b >> fWriteBasket;
2937 b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2938 b >> djunk; fEntries = (Long64_t)djunk;
2939 b >> djunk; fTotBytes = (Long64_t)djunk;
2940 b >> djunk; fZipBytes = (Long64_t)djunk;
2941 b >> fOffset;
2942 fBranches.Streamer(b);
2943 fLeaves.Streamer(b);
2945 for (i=0;i<fNleaves;i++) {
2946 TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2947 leaf->SetBranch(this);
2948 }
2949 fBaskets.Streamer(b);
2950 Int_t nbaskets = fBaskets.GetEntries();
2951 for (j=fWriteBasket,n=0;j>0 && n<nbaskets;--j) {
2953 if (bk) {
2954 bk->SetBranch(this);
2955 //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2956 ++n;
2957 }
2958 }
2960 b >> n;
2961 for (i=0;i<n;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2963 if (v > 4) {
2964 n = b.ReadArray(fBasketBytes);
2965 } else {
2966 for (n=0;n<fMaxBaskets;n++) fBasketBytes[n] = 0;
2967 }
2968 if (v < 2) {
2970 for (n=0;n<fWriteBasket;n++) {
2971 TBasket *basket = GetBasketImpl(n, nullptr);
2972 fBasketSeek[n] = basket ? basket->GetSeekKey() : 0;
2973 }
2974 } else {
2976 b >> n;
2977 for (n=0;n<fMaxBaskets;n++) {
2978 Int_t aseek;
2979 b >> aseek;
2980 fBasketSeek[n] = Long64_t(aseek);
2981 }
2982 }
2983 if (v > 2) {
2984 fFileName.Streamer(b);
2985 }
2986 fDirectory = 0;
2987 if (v < 4) SetAutoDelete(kTRUE);
2989 gROOT->SetReadingObject(kFALSE);
2990 b.CheckByteCount(R__s, R__c, TBranch::IsA());
2991 //====end of old versions
2992 if (IsA() == TBranch::Class()) {
2993 if (fNleaves == 0) {
2995 } else if (fNleaves == 1) {
2997 } else if (fNleaves == 2) {
2999 } else {
3001 }
3002 }
3003 } else {
3004 Int_t maxBaskets = fMaxBaskets;
3006 Int_t lastBasket = fMaxBaskets;
3007 if (fMaxBaskets < 10) fMaxBaskets = 10;
3008
3009 TBasket **stash = new TBasket *[lastBasket];
3010 for (Int_t i = 0; i < lastBasket; ++i) {
3011 TBasket *ba = (TBasket *)fBaskets.UncheckedAt(i);
3012 if (ba && (fBasketBytes[i] || ba->GetNevBuf()==0)) {
3013 // Already on disk or empty.
3014 stash[i] = ba;
3015 fBaskets[i] = nullptr;
3016 } else {
3017 stash[i] = nullptr;
3018 }
3019 }
3020
3021 b.WriteClassBuffer(TBranch::Class(), this);
3022
3023 for (Int_t i = 0; i < lastBasket; ++i) {
3024 if (stash[i]) fBaskets[i] = stash[i];
3025 }
3026
3027 delete[] stash;
3028 fMaxBaskets = maxBaskets;
3029 }
3030}
3031
3032////////////////////////////////////////////////////////////////////////////////
3033/// Write the current basket to disk and return the number of bytes
3034/// written to the file.
3035
3037{
3038 Int_t nevbuf = basket->GetNevBuf();
3039 if (fEntryOffsetLen > 10 && (4*nevbuf) < fEntryOffsetLen ) {
3040 // Make sure that the fEntryOffset array does not stay large unnecessarily.
3041 fEntryOffsetLen = nevbuf < 3 ? 10 : 4*nevbuf; // assume some fluctuations.
3042 } else if (fEntryOffsetLen && nevbuf > fEntryOffsetLen) {
3043 // Increase the array ...
3044 fEntryOffsetLen = 2*nevbuf; // assume some fluctuations.
3045 }
3046
3047 // Note: captures `basket`, `where`, and `this` by value; modifies the TBranch and basket,
3048 // as we make a copy of the pointer. We cannot capture `basket` by reference as the pointer
3049 // itself might be modified after `WriteBasketImpl` exits.
3050 auto doUpdates = [=]() {
3051 Int_t nout = basket->WriteBuffer(); // Write buffer
3052 if (nout < 0) Error("TBranch::WriteBasketImpl", "basket's WriteBuffer failed.\n");
3053 fBasketBytes[where] = basket->GetNbytes();
3054 fBasketSeek[where] = basket->GetSeekKey();
3055 Int_t addbytes = basket->GetObjlen() + basket->GetKeylen();
3056 TBasket *reusebasket = 0;
3057 if (nout>0) {
3058 // The Basket was written so we can now safely reuse it.
3059 fBaskets[where] = 0;
3060
3061 reusebasket = basket;
3062 reusebasket->WriteReset();
3063
3064 fZipBytes += nout;
3065 fTotBytes += addbytes;
3066 fTree->AddTotBytes(addbytes);
3067 fTree->AddZipBytes(nout);
3068#ifdef R__TRACK_BASKET_ALLOC_TIME
3069 fTree->AddAllocationTime(reusebasket->GetResetAllocationTime());
3070#endif
3072 }
3073
3074 if (where==fWriteBasket) {
3075 ++fWriteBasket;
3076 if (fWriteBasket >= fMaxBaskets) {
3078 }
3079 if (reusebasket && reusebasket == fCurrentBasket) {
3080 // The 'current' basket has Reset, so if we need it we will need
3081 // to reload it.
3082 fCurrentBasket = 0;
3083 fFirstBasketEntry = -1;
3084 fNextBasketEntry = -1;
3085 }
3088 } else {
3089 --fNBaskets;
3090 fBaskets[where] = 0;
3091 basket->DropBuffers();
3092 if (basket == fCurrentBasket) {
3093 fCurrentBasket = 0;
3094 fFirstBasketEntry = -1;
3095 fNextBasketEntry = -1;
3096 }
3097 delete basket;
3098 }
3099 return nout;
3100 };
3101 if (imtHelper) {
3102 imtHelper->Run(doUpdates);
3103 return 0;
3104 } else {
3105 return doUpdates();
3106 }
3107}
3108
3109////////////////////////////////////////////////////////////////////////////////
3110///set the first entry number (case of TBranchSTL)
3111
3113{
3114 fFirstEntry = entry;
3115 fEntries = 0;
3116 fEntryNumber = entry;
3117 if( fBasketEntry )
3118 fBasketEntry[0] = entry;
3119 for( Int_t i = 0; i < fBranches.GetEntriesFast(); ++i )
3120 ((TBranch*)fBranches[i])->SetFirstEntry( entry );
3121}
3122
3123////////////////////////////////////////////////////////////////////////////////
3124/// If the branch address is not set, we set all addresses starting with
3125/// the top level parent branch.
3126
3128{
3129 // Nothing to do for regular branch, the TLeaf already did it.
3130}
3131
3132////////////////////////////////////////////////////////////////////////////////
3133/// Refresh the value of fDirectory (i.e. where this branch writes/reads its buffers)
3134/// with the current value of fTree->GetCurrentFile unless this branch has been
3135/// redirected to a different file. Also update the sub-branches.
3136
3138{
3140 if (fFileName.Length() == 0) {
3141 fDirectory = file;
3142
3143 // Apply to all existing baskets.
3144 TIter nextb(GetListOfBaskets());
3145 TBasket *basket;
3146 while ((basket = (TBasket*)nextb())) {
3147 basket->SetParent(file);
3148 }
3149 }
3150
3151 // Apply to sub-branches as well.
3152 TIter next(GetListOfBranches());
3153 TBranch *branch;
3154 while ((branch = (TBranch*)next())) {
3155 branch->UpdateFile();
3156 }
3157}
void tobuf(char *&buf, Bool_t x)
Definition: Bytes.h:57
void Class()
Definition: Class.C:29
#define R__likely(expr)
Definition: RConfig.hxx:612
#define R__unlikely(expr)
Definition: RConfig.hxx:611
#define b(i)
Definition: RSha256.hxx:100
const Ssiz_t kNPOS
Definition: RtypesCore.h:111
int Int_t
Definition: RtypesCore.h:41
short Version_t
Definition: RtypesCore.h:61
int Ssiz_t
Definition: RtypesCore.h:63
char Char_t
Definition: RtypesCore.h:29
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Stat_t
Definition: RtypesCore.h:73
long long Long64_t
Definition: RtypesCore.h:69
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
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:415
void Printf(const char *fmt,...)
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
#define R__LOCKGUARD_IMT(mutex)
#define R__LOCKGUARD(mutex)
#define gPad
Definition: TVirtualPad.h:286
#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:62
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:717
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:709
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:460
virtual Int_t DropBuffers()
Drop buffers of this basket if it is not the current basket.
Definition: TBasket.cxx:169
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:306
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:123
virtual void SetReadMode()
Set read mode of basket.
Definition: TBasket.cxx:921
virtual void SetWriteMode()
Set write mode of basket.
Definition: TBasket.cxx:930
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:694
virtual void ReadResetBuffer(Int_t basketnumber)
Reset the read basket TBuffer memory allocation if needed.
Definition: TBasket.cxx:729
virtual Int_t WriteBuffer()
Write buffer of this basket on the current file.
Definition: TBasket.cxx:1127
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:802
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:1907
virtual void SetupAddresses()
If the branch address is not set, we set all addresses starting with the top level parent branch.
Definition: TBranch.cxx:3127
virtual void ResetAddress()
Reset the address of the branch.
Definition: TBranch.cxx:2488
virtual void SetAutoDelete(Bool_t autodel=kTRUE)
Set the automatic delete bit.
Definition: TBranch.cxx:2553
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:1300
const char * GetIconName() const
Return icon name depending on type of branch.
Definition: TBranch.cxx:1308
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:1767
TBasket * GetBasketImpl(Int_t basket, TBuffer *user_buffer)
Return pointer to basket basketnumber in this Branch.
Definition: TBranch.cxx:1202
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:703
virtual Long64_t GetBasketSeek(Int_t basket) const
Return address of basket in the file.
Definition: TBranch.cxx:1278
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:2642
Int_t BackFill()
Loop on all leaves of this branch to back fill Basket buffer.
Definition: TBranch.cxx:657
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:2290
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:1891
@ 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:1325
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:734
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:1288
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:2298
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:2658
void ExpandBasketArrays()
Increase BasketEntry buffer of a minimum of 10 locations and a maximum of 50 per cent of current size...
Definition: TBranch.cxx:802
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:1579
TIOFeatures GetIOFeatures() const
Returns the IO settings currently in use for this branch.
Definition: TBranch.cxx:2092
void FillLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to fill Basket buffer.
Definition: TBranch.cxx:2333
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:2511
TBranch * GetSubBranch(const TBranch *br) const
Find the parent branch of child.
Definition: TBranch.cxx:2000
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:2773
Int_t FlushBaskets()
Flush to disk all the baskets of this branch and any of subbranches.
Definition: TBranch.cxx:1112
void ReadLeaves2Impl(TBuffer &b)
Read two leaves without the overhead of a loop.
Definition: TBranch.cxx:2324
TBasket * GetFreshCluster()
Drops the cluster two behind the current cluster and returns a fresh basket by either reusing or crea...
Definition: TBranch.cxx:1826
virtual void SetAddress(void *add)
Set address of this branch.
Definition: TBranch.cxx:2519
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:599
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:2316
Int_t GetBulkEntries(Long64_t, TBuffer &)
Read as many events as possible into the given buffer, using zero-copy mechanisms.
Definition: TBranch.cxx:1425
Bool_t IsFolder() const
Return kTRUE if more than one leaf or browsables, kFALSE otherwise.
Definition: TBranch.cxx:2108
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:1635
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:2600
Bool_t IsAutoDelete() const
Return kTRUE if an existing object in a TBranchObject must be deleted.
Definition: TBranch.cxx:2100
virtual TLeaf * FindLeaf(const char *name)
Find the leaf corresponding to the name 'searchname'.
Definition: TBranch.cxx:1057
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:2582
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:3036
virtual void UpdateFile()
Refresh the value of fDirectory (i.e.
Definition: TBranch.cxx:3137
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:2178
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:912
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:1726
virtual void Browse(TBrowser *b)
Browser interface.
Definition: TBranch.cxx:676
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:2120
virtual void ResetAfterMerge(TFileMergeInfo *)
Reset a Branch.
Definition: TBranch.cxx:2435
void ReadLeaves0Impl(TBuffer &b)
Read zero leaves without the overhead of a loop.
Definition: TBranch.cxx:2309
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:1969
TString GetRealFileName() const
Get real file name.
Definition: TBranch.cxx:1920
virtual TBranch * FindBranch(const char *name)
Find the immediate sub-branch with passed name.
Definition: TBranch.cxx:1011
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:2282
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:2146
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:2057
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:2764
virtual void SetFile(TFile *file=0)
Set file where this branch writes/reads its buffers.
Definition: TBranch.cxx:2700
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:1158
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:1707
virtual void SetFirstEntry(Long64_t entry)
set the first entry number (case of TBranchSTL)
Definition: TBranch.cxx:3112
Long64_t GetTotalSize(Option_t *option="") const
Return total number of bytes in the branch (including current buffer)
Definition: TBranch.cxx:2038
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:2075
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:2566
virtual void Refresh(TBranch *b)
Refresh this branch using new information in b This function is called by TTree::Refresh.
Definition: TBranch.cxx:2345
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:2784
TObjArray * GetListOfLeaves()
Definition: TBranch.h:245
virtual void SetEntries(Long64_t entries)
Set the number of entries in this branch.
Definition: TBranch.cxx:2675
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:1960
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:2394
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:1979
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:833
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:2620
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:1400
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:75
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:41
Describe directory structure in memory.
Definition: TDirectory.h:34
virtual TFile * GetFile() const
Definition: TDirectory.h:157
virtual Bool_t IsWritable() const
Definition: TDirectory.h:173
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:48
Int_t GetCompressionSettings() const
Definition: TFile.h:394
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:3923
virtual Long64_t GetSeekKey() const
Definition: TKey.h:86
Int_t GetKeylen() const
Definition: TKey.h:81
Int_t GetObjlen() const
Definition: TKey.h:84
Int_t GetNbytes() const
Definition: TKey.h:83
virtual void SetParent(const TObject *parent)
Set parent in key buffer.
Definition: TKey.cxx:1273
TBuffer * GetBufferRef() const
Definition: TKey.h:76
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:379
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:157
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:386
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:234
void Add(TObject *obj)
Definition: TObjArray.h:74
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
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:355
void SetLast(Int_t last)
Set index of last object in array, effectively truncating the array.
Definition: TObjArray.cxx:774
Int_t GetLast() const
Return index of last object in array.
Definition: TObjArray.cxx:576
Int_t LowerBound() const
Definition: TObjArray.h:97
virtual void AddAt(TObject *obj, Int_t idx)
Add object at position ids.
Definition: TObjArray.cxx:253
virtual TObject * RemoveAt(Int_t idx)
Remove object at index idx.
Definition: TObjArray.cxx:693
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
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:866
R__ALWAYS_INLINE Bool_t IsZombie() const
Definition: TObject.h:134
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:880
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:908
void MakeZombie()
Definition: TObject.h:49
void ResetBit(UInt_t f)
Definition: TObject.h:171
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:1265
virtual const char * DirName(const char *pathname)
Return the directory name in pathname.
Definition: TSystem.cxx:1014
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition: TSystem.cxx:942
virtual Bool_t IsAbsoluteFileName(const char *dir)
Return true if dir is an absolute pathname.
Definition: TSystem.cxx:959
Helper class to iterate over cluster of baskets.
Definition: TTree.h:258
Long64_t Previous()
Move on to the previous cluster and return the starting entry of this previous cluster.
Definition: TTree.cxx:680
Long64_t GetStartEntry()
Definition: TTree.h:290
Long64_t Next()
Move on to the next cluster and return the starting entry of this next cluster.
Definition: TTree.cxx:636
Long64_t GetNextEntry()
Definition: TTree.h:295
A TTree represents a columnar dataset.
Definition: TTree.h:72
virtual TVirtualPerfStats * GetPerfStats() const
Definition: TTree.h:493
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition: TTree.cxx:5326
void AddAllocationCount(UInt_t count)
Definition: TTree.h:325
virtual TObjArray * GetListOfLeaves()
Definition: TTree.h:476
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition: TTree.cxx:5338
virtual void IncrementTotalBuffers(Int_t nbytes)
Definition: TTree.h:533
TDirectory * GetDirectory() const
Definition: TTree.h:449
TTreeCache * GetReadCache(TFile *file) const
Find and return the TTreeCache registered with the file and which may contain branches for us.
Definition: TTree.cxx:6174
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:475
virtual void AddZipBytes(Int_t zip)
Definition: TTree.h:320
virtual TBasket * CreateBasket(TBranch *)
Create a basket for this tree and given branch.
Definition: TTree.cxx:3647
@ kOnlyFlushAtCluster
If set, the branch's buffers will grow until an event cluster boundary is hit, guaranteeing a basket ...
Definition: TTree.h:244
@ kCircular
Definition: TTree.h:240
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:419
virtual Bool_t GetClusterPrefetch() const
Definition: TTree.h:444
virtual void AddTotBytes(Int_t tot)
Definition: TTree.h:319
virtual Long64_t GetAutoFlush() const
Definition: TTree.h:434
virtual Long64_t GetMaxVirtualSize() const
Definition: TTree.h:487
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:385
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:95
auto * l
Definition: textangle.C:4