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