Logo ROOT   master
Reference Guide
TBranch.cxx
Go to the documentation of this file.
1 // @(#)root/tree:$Id$
2 // Author: Rene Brun 12/01/96
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TBranchCacheInfo.h"
13 
14 #include "TBranch.h"
15 
16 #include "Bytes.h"
17 #include "Compression.h"
18 #include "TBasket.h"
19 #include "TBranchBrowsable.h"
20 #include "TBrowser.h"
21 #include "TBuffer.h"
22 #include "TClass.h"
23 #include "TBufferFile.h"
24 #include "TClonesArray.h"
25 #include "TFile.h"
26 #include "TLeaf.h"
27 #include "TLeafB.h"
28 #include "TLeafC.h"
29 #include "TLeafD.h"
30 #include "TLeafD32.h"
31 #include "TLeafF.h"
32 #include "TLeafF16.h"
33 #include "TLeafI.h"
34 #include "TLeafL.h"
35 #include "TLeafO.h"
36 #include "TLeafObject.h"
37 #include "TLeafS.h"
38 #include "TMessage.h"
39 #include "TROOT.h"
40 #include "TSystem.h"
41 #include "TMath.h"
42 #include "TTree.h"
43 #include "TTreeCache.h"
44 #include "TTreeCacheUnzip.h"
45 #include "TVirtualMutex.h"
46 #include "TVirtualPad.h"
47 #include "TVirtualPerfStats.h"
48 
49 #include "TBranchIMTHelper.h"
50 
51 #include "ROOT/TIOFeatures.hxx"
52 
53 #include <atomic>
54 #include <cstddef>
55 #include <string.h>
56 #include <stdio.h>
57 
58 
60 
61 /** \class TBranch
62 \ingroup tree
63 
64 A TTree is a list of TBranches
65 
66 A TBranch supports:
67  - The list of TLeaf describing this branch.
68  - The list of TBasket (branch buffers).
69 
70 See TBranch structure in TTree.
71 
72 See also specialized branches:
73  - TBranchObject in case the branch is one object
74  - TBranchClones in case the branch is an array of clone objects
75 */
76 
78 
79 
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 /// Default constructor. Used for I/O by default.
83 
85 : TNamed()
86 , TAttFill(0, 1001)
87 , fCompress(0)
88 , fBasketSize(32000)
89 , fEntryOffsetLen(1000)
90 , fWriteBasket(0)
91 , fEntryNumber(0)
92 , fExtraBasket(nullptr)
93 , fOffset(0)
94 , fMaxBaskets(10)
95 , fNBaskets(0)
96 , fSplitLevel(0)
97 , fNleaves(0)
98 , fReadBasket(0)
99 , fReadEntry(-1)
100 , fFirstBasketEntry(-1)
101 , fNextBasketEntry(-1)
102 , fCurrentBasket(0)
103 , fEntries(0)
104 , fFirstEntry(0)
105 , fTotBytes(0)
106 , fZipBytes(0)
107 , fBranches()
108 , fLeaves()
109 , fBaskets(fMaxBaskets)
110 , fBasketBytes(0)
111 , fBasketEntry(0)
112 , fBasketSeek(0)
113 , fTree(0)
114 , fMother(0)
115 , fParent(0)
116 , fAddress(0)
117 , fDirectory(0)
118 , fFileName("")
119 , fEntryBuffer(0)
120 , fTransientBuffer(0)
121 , fBrowsables(0)
122 , fBulk(*this)
123 , fSkipZip(kFALSE)
124 , fReadLeaves(&TBranch::ReadLeavesImpl)
125 , fFillLeaves(&TBranch::FillLeavesImpl)
126 {
128 }
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// Create a Branch as a child of a Tree
132 ///
133 /// * address is the address of the first item of a structure
134 /// or the address of a pointer to an object (see example in TTree.cxx).
135 /// * leaflist is the concatenation of all the variable names and types
136 /// separated by a colon character :
137 /// The variable name and the variable type are separated by a
138 /// slash (/). The variable type must be 1 character. (Characters
139 /// after the first are legal and will be appended to the visible
140 /// name of the leaf, but have no effect.) If no type is given, the
141 /// type of the variable is assumed to be the same as the previous
142 /// variable. If the first variable does not have a type, it is
143 /// assumed of type F by default. The list of currently supported
144 /// types is given below:
145 /// - `C` : a character string terminated by the 0 character
146 /// - `B` : an 8 bit signed integer (`Char_t`)
147 /// - `b` : an 8 bit unsigned integer (`UChar_t`)
148 /// - `S` : a 16 bit signed integer (`Short_t`)
149 /// - `s` : a 16 bit unsigned integer (`UShort_t`)
150 /// - `I` : a 32 bit signed integer (`Int_t`)
151 /// - `i` : a 32 bit unsigned integer (`UInt_t`)
152 /// - `F` : a 32 bit floating point (`Float_t`)
153 /// - `f` : a 24 bit floating point with truncated mantissa (`Float16_t`)
154 /// - `D` : a 64 bit floating point (`Double_t`)
155 /// - `d` : a 24 bit truncated floating point (`Double32_t`)
156 /// - `L` : a 64 bit signed integer (`Long64_t`)
157 /// - `l` : a 64 bit unsigned integer (`ULong64_t`)
158 /// - `O` : [the letter `o`, not a zero] a boolean (`Bool_t`)
159 ///
160 /// Arrays of values are supported with the following syntax:
161 /// - If leaf name has the form var[nelem], where nelem is alphanumeric, then
162 /// if nelem is a leaf name, it is used as the variable size of the array,
163 /// otherwise return 0.
164 /// The leaf referred to by nelem **MUST** be an int (/I),
165 /// - If leaf name has the form var[nelem], where nelem is a non-negative integers, then
166 /// it is used as the fixed size of the array.
167 /// - If leaf name has the form of a multi dimension array (e.g. var[nelem][nelem2])
168 /// where nelem and nelem2 are non-negative integers) then
169 /// it is used as a 2 dimensional array of fixed size.
170 /// - In case of the truncated floating point types (Float16_t and Double32_t) you can
171 /// furthermore specify the range in the style [xmin,xmax] or [xmin,xmax,nbits] after
172 /// the type character. See `TStreamerElement::GetRange()` for further information.
173 /// - Any of other form is not supported.
174 ///
175 /// Note that the TTree will assume that all the item are contiguous in memory.
176 /// On some platform, this is not always true of the member of a struct or a class,
177 /// due to padding and alignment. Sorting your data member in order of decreasing
178 /// sizeof usually leads to their being contiguous in memory.
179 ///
180 /// * bufsize is the buffer size in bytes for this branch
181 /// The default value is 32000 bytes and should be ok for most cases.
182 /// You can specify a larger value (e.g. 256000) if your Tree is not split
183 /// and each entry is large (Megabytes)
184 /// A small value for bufsize is optimum if you intend to access
185 /// the entries in the Tree randomly and your Tree is in split mode.
186 ///
187 /// See an example of a Branch definition in the TTree constructor.
188 ///
189 /// Note that in case the data type is an object, this branch can contain
190 /// only this object.
191 ///
192 /// Note that this function is invoked by TTree::Branch
193 
194 TBranch::TBranch(TTree *tree, const char *name, void *address, const char *leaflist, Int_t basketsize, Int_t compress)
195  : TNamed(name, leaflist)
196 , TAttFill(0, 1001)
197 , fCompress(compress)
198 , fBasketSize((basketsize < 100) ? 100 : basketsize)
199 , fEntryOffsetLen(0)
200 , fWriteBasket(0)
201 , fEntryNumber(0)
202 , fExtraBasket(nullptr)
203 , fIOFeatures(tree ? tree->GetIOFeatures().GetFeatures() : 0)
204 , fOffset(0)
205 , fMaxBaskets(10)
206 , fNBaskets(0)
207 , fSplitLevel(0)
208 , fNleaves(0)
209 , fReadBasket(0)
210 , fReadEntry(-1)
211 , fFirstBasketEntry(-1)
212 , fNextBasketEntry(-1)
213 , fCurrentBasket(0)
214 , fEntries(0)
215 , fFirstEntry(0)
216 , fTotBytes(0)
217 , fZipBytes(0)
218 , fBranches()
219 , fLeaves()
220 , fBaskets(fMaxBaskets)
221 , fBasketBytes(0)
222 , fBasketEntry(0)
223 , fBasketSeek(0)
224 , fTree(tree)
225 , fMother(0)
226 , fParent(0)
227 , fAddress((char *)address)
228 , fDirectory(fTree->GetDirectory())
229 , fFileName("")
230 , fEntryBuffer(0)
231 , fTransientBuffer(0)
232 , fBrowsables(0)
233 , fBulk(*this)
234 , fSkipZip(kFALSE)
235 , fReadLeaves(&TBranch::ReadLeavesImpl)
236 , fFillLeaves(&TBranch::FillLeavesImpl)
237 {
238  Init(name,leaflist,compress);
239 }
240 
241 ////////////////////////////////////////////////////////////////////////////////
242 /// Create a Branch as a child of another Branch
243 ///
244 /// See documentation for
245 /// TBranch::TBranch(TTree *, const char *, void *, const char *, Int_t, Int_t)
246 
247 TBranch::TBranch(TBranch *parent, const char *name, void *address, const char *leaflist, Int_t basketsize,
248  Int_t compress)
249 : TNamed(name, leaflist)
250 , TAttFill(0, 1001)
251 , fCompress(compress)
252 , fBasketSize((basketsize < 100) ? 100 : basketsize)
253 , fEntryOffsetLen(0)
254 , fWriteBasket(0)
255 , fEntryNumber(0)
256 , fExtraBasket(nullptr)
257 , fIOFeatures(parent->fIOFeatures)
258 , fOffset(0)
259 , fMaxBaskets(10)
260 , fNBaskets(0)
261 , fSplitLevel(0)
262 , fNleaves(0)
263 , fReadBasket(0)
264 , fReadEntry(-1)
265 , fFirstBasketEntry(-1)
266 , fNextBasketEntry(-1)
267 , fCurrentBasket(0)
268 , fEntries(0)
269 , fFirstEntry(0)
270 , fTotBytes(0)
271 , fZipBytes(0)
272 , fBranches()
273 , fLeaves()
274 , fBaskets(fMaxBaskets)
275 , fBasketBytes(0)
276 , fBasketEntry(0)
277 , fBasketSeek(0)
278 , fTree(parent ? parent->GetTree() : 0)
279 , fMother(parent ? parent->GetMother() : 0)
280 , fParent(parent)
281 , fAddress((char *)address)
282 , fDirectory(fTree ? fTree->GetDirectory() : 0)
283 , fFileName("")
284 , fEntryBuffer(0)
285 , fTransientBuffer(0)
286 , fBrowsables(0)
287 , fBulk(*this)
288 , fSkipZip(kFALSE)
289 , fReadLeaves(&TBranch::ReadLeavesImpl)
290 , fFillLeaves(&TBranch::FillLeavesImpl)
291 {
292  Init(name,leaflist,compress);
293 }
294 
295 void TBranch::Init(const char* name, const char* leaflist, Int_t compress)
296 {
297  // Initialization routine called from the constructor. This should NOT be made virtual.
298 
300  if ((compress == -1) && fTree->GetDirectory()) {
301  TFile* bfile = fTree->GetDirectory()->GetFile();
302  if (bfile) {
304  }
305  }
306 
310 
311  for (Int_t i = 0; i < fMaxBaskets; ++i) {
312  fBasketBytes[i] = 0;
313  fBasketEntry[i] = 0;
314  fBasketSeek[i] = 0;
315  }
316 
317  //
318  // Decode the leaflist (search for : as separator).
319  //
320 
321  char* nameBegin = const_cast<char*>(leaflist);
322  Int_t offset = 0;
323  auto len = strlen(leaflist);
324  // FIXME: Make these string streams instead.
325  char* leafname = new char[len + 1];
326  char* leaftype = new char[320];
327  // Note: The default leaf type is a float.
328  strlcpy(leaftype, "F",320);
329  char* pos = const_cast<char*>(leaflist);
330  const char* leaflistEnd = leaflist + len;
331  for (; pos <= leaflistEnd; ++pos) {
332  // -- Scan leaf specification and create leaves.
333  if ((*pos == ':') || (*pos == 0)) {
334  // -- Reached end of a leaf spec, create a leaf.
335  Int_t lenName = pos - nameBegin;
336  char* ctype = 0;
337  if (lenName) {
338  strncpy(leafname, nameBegin, lenName);
339  leafname[lenName] = 0;
340  ctype = strstr(leafname, "/");
341  if (ctype) {
342  *ctype = 0;
343  strlcpy(leaftype, ctype + 1,320);
344  }
345  }
346  if (lenName == 0 || ctype == leafname) {
347  Warning("TBranch","No name was given to the leaf number '%d' in the leaflist of the branch '%s'.",fNleaves,name);
348  snprintf(leafname,640,"__noname%d",fNleaves);
349  }
350  TLeaf* leaf = 0;
351  if (leaftype[1] == '[' && !strchr(leaftype, ',')) {
352  Warning("TBranch", "Array size for branch '%s' must be specified after leaf name, not after the type name!", name);
353  // and continue for backward compatibility?
354  } else if (leaftype[1] && !strchr(leaftype, ',')) {
355  Warning("TBranch", "Extra characters after type tag '%s' for branch '%s'; must be one character.", leaftype, name);
356  // and continue for backward compatibility?
357  }
358  if (*leaftype == 'C') {
359  leaf = new TLeafC(this, leafname, leaftype);
360  } else if (*leaftype == 'O') {
361  leaf = new TLeafO(this, leafname, leaftype);
362  } else if (*leaftype == 'B') {
363  leaf = new TLeafB(this, leafname, leaftype);
364  } else if (*leaftype == 'b') {
365  leaf = new TLeafB(this, leafname, leaftype);
366  leaf->SetUnsigned();
367  } else if (*leaftype == 'S') {
368  leaf = new TLeafS(this, leafname, leaftype);
369  } else if (*leaftype == 's') {
370  leaf = new TLeafS(this, leafname, leaftype);
371  leaf->SetUnsigned();
372  } else if (*leaftype == 'I') {
373  leaf = new TLeafI(this, leafname, leaftype);
374  } else if (*leaftype == 'i') {
375  leaf = new TLeafI(this, leafname, leaftype);
376  leaf->SetUnsigned();
377  } else if (*leaftype == 'F') {
378  leaf = new TLeafF(this, leafname, leaftype);
379  } else if (*leaftype == 'f') {
380  leaf = new TLeafF16(this, leafname, leaftype);
381  } else if (*leaftype == 'L') {
382  leaf = new TLeafL(this, leafname, leaftype);
383  } else if (*leaftype == 'l') {
384  leaf = new TLeafL(this, leafname, leaftype);
385  leaf->SetUnsigned();
386  } else if (*leaftype == 'D') {
387  leaf = new TLeafD(this, leafname, leaftype);
388  } else if (*leaftype == 'd') {
389  leaf = new TLeafD32(this, leafname, leaftype);
390  }
391  if (!leaf) {
392  Error("TLeaf", "Illegal data type for %s/%s", name, leaflist);
393  delete[] leaftype;
394  delete [] leafname;
395  MakeZombie();
396  return;
397  }
398  if (leaf->IsZombie()) {
399  delete leaf;
400  leaf = 0;
401  auto msg = "Illegal leaf: %s/%s. If this is a variable size C array it's possible that the branch holding the size is not available.";
402  Error("TBranch", msg, name, leaflist);
403  delete [] leafname;
404  delete[] leaftype;
405  MakeZombie();
406  return;
407  }
408  leaf->SetBranch(this);
409  leaf->SetAddress((char*) (fAddress + offset));
410  leaf->SetOffset(offset);
411  if (leaf->GetLeafCount()) {
412  // -- Leaf is a varying length array, we need an offset array.
413  fEntryOffsetLen = 1000;
414  }
415  if (leaf->InheritsFrom(TLeafC::Class())) {
416  // -- Leaf is a character string, we need an offset array.
417  fEntryOffsetLen = 1000;
418  }
419  ++fNleaves;
420  fLeaves.Add(leaf);
421  fTree->GetListOfLeaves()->Add(leaf);
422  if (*pos == 0) {
423  // -- We reached the end of the leaf specification.
424  break;
425  }
426  nameBegin = pos + 1;
427  offset += leaf->GetLenType() * leaf->GetLen();
428  }
429  }
430  delete[] leafname;
431  leafname = 0;
432  delete[] leaftype;
433  leaftype = 0;
434 
435 }
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// Destructor.
439 
441 {
442  delete fBrowsables;
443  fBrowsables = 0;
444 
445  // Note: We do *not* have ownership of the buffer.
446  fEntryBuffer = 0;
447 
448  delete [] fBasketSeek;
449  fBasketSeek = 0;
450 
451  delete [] fBasketEntry;
452  fBasketEntry = 0;
453 
454  delete [] fBasketBytes;
455  fBasketBytes = 0;
456 
457  fBaskets.Delete();
458  fNBaskets = 0;
459  fCurrentBasket = 0;
460  fFirstBasketEntry = -1;
461  fNextBasketEntry = -1;
462 
463  // Remove our leaves from our tree's list of leaves.
464  if (fTree) {
465  TObjArray* lst = fTree->GetListOfLeaves();
466  if (lst && lst->GetLast()!=-1) {
467  lst->RemoveAll(&fLeaves);
468  }
469  }
470  // And delete our leaves.
471  fLeaves.Delete();
472 
473  fBranches.Delete();
474 
475  // If we are in a directory and that directory is not the same
476  // directory that our tree is in, then try to find an open file
477  // with the name fFileName. If we find one, delete that file.
478  // We are attempting to close any alternate file which we have
479  // been directed to write our baskets to.
480  // FIXME: We make no attempt to check if someone else might be
481  // using this file. This is very user hostile. A violation
482  // of the principle of least surprises.
483  //
484  // Warning. Must use FindObject by name instead of fDirectory->GetFile()
485  // because two branches may point to the same file and the file
486  // may have already been deleted in the previous branch.
487  if (fDirectory && (!fTree || fDirectory != fTree->GetDirectory())) {
488  TString bFileName( GetRealFileName() );
489 
491  TFile* file = (TFile*)gROOT->GetListOfFiles()->FindObject(bFileName);
492  if (file){
493  file->Close();
494  delete file;
495  file = 0;
496  }
497  }
498 
499  fTree = 0;
500  fDirectory = 0;
501 
502  if (fTransientBuffer) {
503  delete fTransientBuffer;
504  fTransientBuffer = 0;
505  }
506 }
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Returns the transient buffer currently used by this TBranch for reading/writing baskets.
510 
512 {
513  if (fTransientBuffer) {
514  if (fTransientBuffer->BufferSize() < size) {
515  fTransientBuffer->Expand(size);
516  }
517  return fTransientBuffer;
518  }
520  return fTransientBuffer;
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// Add the basket to this branch.
525 ///
526 /// Warning: if the basket are not 'flushed/copied' in the same
527 /// order as they were created, this will induce a slow down in
528 /// the insert (since we'll need to move all the record that are
529 /// entere 'too early').
530 /// Warning we also assume that the __current__ write basket is
531 /// not present (aka has been removed).
532 
533 void TBranch::AddBasket(TBasket& b, Bool_t ondisk, Long64_t startEntry)
534 {
535  TBasket *basket = &b;
536 
537  basket->SetBranch(this);
538 
539  if (fWriteBasket >= fMaxBaskets) {
541  }
542  Int_t where = fWriteBasket;
543 
544  if (where && startEntry < fBasketEntry[where-1]) {
545  // Need to find the right location and move the possible baskets
546 
547  if (!ondisk) {
548  Warning("AddBasket","The assumption that out-of-order basket only comes from disk based ntuple is false.");
549  }
550 
551  if (startEntry < fBasketEntry[0]) {
552  where = 0;
553  } else {
554  for(Int_t i=fWriteBasket-1; i>=0; --i) {
555  if (fBasketEntry[i] < startEntry) {
556  where = i+1;
557  break;
558  } else if (fBasketEntry[i] == startEntry) {
559  Error("AddBasket","An out-of-order basket matches the entry number of an existing basket.");
560  }
561  }
562  }
563 
564  if (where < fWriteBasket) {
565  // We shall move the content of the array
566  for (Int_t j=fWriteBasket; j > where; --j) {
567  fBasketEntry[j] = fBasketEntry[j-1];
568  fBasketBytes[j] = fBasketBytes[j-1];
569  fBasketSeek[j] = fBasketSeek[j-1];
570  }
571  }
572  }
573  fBasketEntry[where] = startEntry;
574 
575  if (ondisk) {
576  fBasketBytes[where] = basket->GetNbytes(); // not for in mem
577  fBasketSeek[where] = basket->GetSeekKey(); // not for in mem
579  ++fWriteBasket;
580  } else {
581  ++fNBaskets;
584  }
585 
586  fEntries += basket->GetNevBuf();
587  fEntryNumber += basket->GetNevBuf();
588  if (ondisk) {
589  fTotBytes += basket->GetObjlen() + basket->GetKeylen() ;
590  fZipBytes += basket->GetNbytes();
591  fTree->AddTotBytes(basket->GetObjlen() + basket->GetKeylen());
592  fTree->AddZipBytes(basket->GetNbytes());
593  }
594 }
595 
596 ////////////////////////////////////////////////////////////////////////////////
597 /// Add the start entry of the write basket (not yet created)
598 
600 {
601  if (fWriteBasket >= fMaxBaskets) {
603  }
604  Int_t where = fWriteBasket;
605 
606  if (where && startEntry < fBasketEntry[where-1]) {
607  // Need to find the right location and move the possible baskets
608 
609  Fatal("AddBasket","The last basket must have the highest entry number (%s/%lld/%d).",GetName(),startEntry,fWriteBasket);
610 
611  }
612  fBasketEntry[where] = startEntry;
614 }
615 
616 ////////////////////////////////////////////////////////////////////////////////
617 /// Loop on all leaves of this branch to back fill Basket buffer.
618 ///
619 /// Use this routine instead of TBranch::Fill when filling a branch individually
620 /// to catch up with the number of entries already in the TTree.
621 ///
622 /// First it calls TBranch::Fill and then if the number of entries of the branch
623 /// reach one of TTree cluster's boundary, the basket is flushed.
624 ///
625 /// The function returns the number of bytes committed to the memory basket.
626 /// If a write error occurs, the number of bytes returned is -1.
627 /// If no data are written, because e.g. the branch is disabled,
628 /// the number of bytes returned is 0.
629 ///
630 /// To insure that the baskets of each cluster are located close by in the
631 /// file, when back-filling multiple branches make sure to call BackFill
632 /// for the same entry for all the branches consecutively
633 /// ~~~ {.cpp}
634 /// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
635 /// for( auto branch : branchCollection) {
636 /// ... Make change to the data associated with the branch ...
637 /// branch->BackFill();
638 /// }
639 /// }
640 /// // Since we loop over all the branches for each new entry
641 /// // all the baskets for a cluster are consecutive in the file.
642 /// ~~~
643 /// rather than doing all the entries of one branch at a time.
644 /// ~~~ {.cpp}
645 /// // Do NOT do things in the following order, it will lead to
646 /// // poorly clustered files.
647 /// for(auto branch : branchCollection) {
648 /// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
649 /// ... Make change to the data associated with the branch ...
650 /// branch->BackFill();
651 /// }
652 /// }
653 /// // Since we loop over all the entries for one branch
654 /// // all the baskets for that branch are consecutive.
655 /// ~~~
656 
658 
659  // Get the end of the next cluster.
660  auto cluster = GetTree()->GetClusterIterator( GetEntries() );
661  cluster.Next();
662  auto endCluster = cluster.GetNextEntry();
663 
664  auto result = FillImpl(nullptr);
665 
666  if ( result && GetEntries() >= endCluster ) {
667  FlushBaskets();
668  }
669 
670  return result;
671 }
672 
673 ////////////////////////////////////////////////////////////////////////////////
674 /// Browser interface.
675 
677 {
678  if (fNleaves > 1) {
679  fLeaves.Browse(b);
680  } else {
681  // Get the name and strip any extra brackets
682  // in order to get the full arrays.
683  TString name = GetName();
684  Int_t pos = name.First('[');
685  if (pos!=kNPOS) name.Remove(pos);
686 
687  GetTree()->Draw(name, "", b ? b->GetDrawOption() : "");
688  if (gPad) gPad->Update();
689  }
690 }
691 
692  ///////////////////////////////////////////////////////////////////////////////
693  /// Loop on all branch baskets. If the file where branch buffers reside is
694  /// writable, free the disk space associated to the baskets of the branch,
695  /// then call Reset(). If the option contains "all", delete also the baskets
696  /// for the subbranches.
697  /// The branch is reset.
698  ///
699  /// NOTE that this function must be used with extreme care. Deleting branch baskets
700  /// fragments the file and may introduce inefficiencies when adding new entries
701  /// in the Tree or later on when reading the Tree.
702 
704 {
705  TString opt = option;
706  opt.ToLower();
707  TFile *file = GetFile(0);
708 
709  if(fDirectory && (fDirectory != gROOT) && fDirectory->IsWritable()) {
710  for(Int_t i=0; i<fWriteBasket; i++) {
711  if (fBasketSeek[i]) file->MakeFree(fBasketSeek[i],fBasketSeek[i]+fBasketBytes[i]-1);
712  }
713  }
714 
715  // process subbranches
716  if (opt.Contains("all")) {
718  Int_t nb = lb->GetEntriesFast();
719  for (Int_t j = 0; j < nb; j++) {
720  TBranch* branch = (TBranch*) lb->UncheckedAt(j);
721  if (branch) branch->DeleteBaskets("all");
722  }
723  }
724  DropBaskets("all");
725  Reset();
726 }
727 
728 ////////////////////////////////////////////////////////////////////////////////
729 /// Loop on all branch baskets. Drop all baskets from memory except readbasket.
730 /// If the option contains "all", drop all baskets including
731 /// read- and write-baskets (unless they are not stored individually on disk).
732 /// The option "all" also lead to DropBaskets being called on the sub-branches.
733 
735 {
736  Bool_t all = kFALSE;
737  if (options && options[0]) {
738  TString opt = options;
739  opt.ToLower();
740  if (opt.Contains("all")) all = kTRUE;
741  }
742 
743  TBasket *basket;
744  Int_t nbaskets = fBaskets.GetEntriesFast();
745 
746  if ( (fNBaskets>1) || all ) {
747  //slow case
748  for (Int_t i=0;i<nbaskets;i++) {
749  basket = (TBasket*)fBaskets.UncheckedAt(i);
750  if (!basket) continue;
751  if ((i == fReadBasket || i == fWriteBasket) && !all) continue;
752  // if the basket is not yet on file but already has event in it
753  // we must continue to avoid dropping the basket (and thus losing data)
754  if (fBasketBytes[i]==0 && basket->GetNevBuf() > 0) continue;
755  basket->DropBuffers();
756  --fNBaskets;
757  fBaskets.RemoveAt(i);
758  if (basket == fCurrentBasket) {
759  fCurrentBasket = 0;
760  fFirstBasketEntry = -1;
761  fNextBasketEntry = -1;
762  }
763  delete basket;
764  }
765 
766  // process subbranches
767  if (all) {
769  Int_t nb = lb->GetEntriesFast();
770  for (Int_t j = 0; j < nb; j++) {
771  TBranch* branch = (TBranch*) lb->UncheckedAt(j);
772  if (!branch) continue;
773  branch->DropBaskets("all");
774  }
775  }
776  } else {
777  //fast case
778  if (nbaskets > 0) {
779  Int_t i = fBaskets.GetLast();
780  basket = (TBasket*)fBaskets.UncheckedAt(i);
781  if (basket && fBasketBytes[i]!=0) {
782  basket->DropBuffers();
783  if (basket == fCurrentBasket) {
784  fCurrentBasket = 0;
785  fFirstBasketEntry = -1;
786  fNextBasketEntry = -1;
787  }
788  delete basket;
789  fBaskets.AddAt(0,i);
790  fBaskets.SetLast(-1);
791  fNBaskets = 0;
792  }
793  }
794  }
795 
796 }
797 
798 ////////////////////////////////////////////////////////////////////////////////
799 /// Increase BasketEntry buffer of a minimum of 10 locations
800 /// and a maximum of 50 per cent of current size.
801 
803 {
804  Int_t newsize = TMath::Max(10,Int_t(1.5*fMaxBaskets));
807  newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
809  newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
810 
811  fMaxBaskets = newsize;
812 
813  fBaskets.Expand(newsize);
814 
815  for (Int_t i=fWriteBasket;i<fMaxBaskets;i++) {
816  fBasketBytes[i] = 0;
817  fBasketEntry[i] = 0;
818  fBasketSeek[i] = 0;
819  }
820 }
821 
822 ////////////////////////////////////////////////////////////////////////////////
823 /// Loop on all leaves of this branch to fill Basket buffer.
824 ///
825 /// If TBranchIMTHelper is non-null and it is time to WriteBasket, then we will
826 /// use TBB to compress in parallel.
827 ///
828 /// The function returns the number of bytes committed to the memory basket.
829 /// If a write error occurs, the number of bytes returned is -1.
830 /// If no data are written, because e.g. the branch is disabled,
831 /// the number of bytes returned is 0.
832 
834 {
835  if (TestBit(kDoNotProcess)) {
836  return 0;
837  }
838 
840  if (!basket) {
841  basket = fTree->CreateBasket(this); // create a new basket
842  if (!basket) return 0;
843  ++fNBaskets;
845  }
846  TBuffer* buf = basket->GetBufferRef();
847 
848  // Fill basket buffer.
849 
850  Int_t nsize = 0;
851 
852  if (buf->IsReading()) {
853  basket->SetWriteMode();
854  }
855 
856  if (!TestBit(kDoNotUseBufferMap)) {
857  buf->ResetMap();
858  }
859 
860  Int_t lnew = 0;
861  Int_t nbytes = 0;
862 
863  if (fEntryBuffer) {
864  nbytes = FillEntryBuffer(basket,buf,lnew);
865  } else {
866  Int_t lold = buf->Length();
867  basket->Update(lold);
868  ++fEntries;
869  ++fEntryNumber;
870  (this->*fFillLeaves)(*buf);
871  if (buf->GetMapCount()) {
872  // The map is used.
874  }
875  lnew = buf->Length();
876  nbytes = lnew - lold;
877  }
878 
879  if (fEntryOffsetLen) {
880  Int_t nevbuf = basket->GetNevBuf();
881  // Total size in bytes of EntryOffset table.
882  nsize = nevbuf * sizeof(Int_t);
883  } else {
884  if (!basket->GetNevBufSize()) {
885  basket->SetNevBufSize(nbytes);
886  }
887  }
888 
889  // Should we create a new basket?
890  // fSkipZip force one entry per buffer (old stuff still maintained for CDF)
891  // Transfer full compressed buffer only
892 
893  // If GetAutoFlush() is less than zero, then we are determining the end of the autocluster
894  // based upon the number of bytes already flushed. This is incompatible with one-basket-per-cluster
895  // (since we will grow the basket indefinitely and never flush!). Hence, we wait until the
896  // first event cluster is written out and *then* enable one-basket-per-cluster mode.
897  bool noFlushAtCluster = !fTree->TestBit(TTree::kOnlyFlushAtCluster) || (fTree->GetAutoFlush() < 0);
898 
899  if (noFlushAtCluster && !fTree->TestBit(TTree::kCircular) &&
901  ((lnew + (2 * nsize) + nbytes) >= fBasketSize))) {
902  Int_t nout = WriteBasketImpl(basket, fWriteBasket, imtHelper);
903  if (nout < 0) Error("TBranch::Fill", "Failed to write out basket.\n");
904  return (nout >= 0) ? nbytes : -1;
905  }
906  return nbytes;
907 }
908 
909 ////////////////////////////////////////////////////////////////////////////////
910 /// Copy the data from fEntryBuffer into the current basket.
911 
913 {
914  Int_t nbytes = 0;
915  Int_t objectStart = 0;
916  Int_t last = 0;
917  Int_t lold = buf->Length();
918 
919  // Handle the special case of fEntryBuffer != 0
920  if (fEntryBuffer->IsA() == TMessage::Class()) {
921  objectStart = 8;
922  }
924  // The buffer given as input has not been decompressed.
925  if (basket->GetNevBuf()) {
926  // If the basket already contains entry we need to close it
927  // out. (This is because we can only transfer full compressed
928  // buffer)
929  WriteBasket(basket,fWriteBasket);
930  // And restart from scratch
931  return Fill();
932  }
933  Int_t startpos = fEntryBuffer->Length();
935  static TBasket toread_fLast;
937  toread_fLast.Streamer(*fEntryBuffer);
939  last = toread_fLast.GetLast();
940  // last now contains the decompressed number of bytes.
941  fEntryBuffer->SetBufferOffset(startpos);
942  buf->SetBufferOffset(0);
944  basket->Update(lold);
945  } else {
946  // We are required to copy starting at the version number (so not
947  // including the class name.
948  // See if byte count is here, if not it class still be a newClass
949  const UInt_t kNewClassTag = 0xFFFFFFFF;
950  const UInt_t kByteCountMask = 0x40000000; // OR the byte count with this
951  UInt_t tag = 0;
952  UInt_t startpos = fEntryBuffer->Length();
953  fEntryBuffer->SetBufferOffset(objectStart);
954  *fEntryBuffer >> tag;
955  if (tag & kByteCountMask) {
956  *fEntryBuffer >> tag;
957  }
958  if (tag == kNewClassTag) {
959  UInt_t maxsize = 256;
960  char* s = new char[maxsize];
961  Int_t name_start = fEntryBuffer->Length();
962  fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end.
963  while (strlen(s) == (maxsize - 1)) {
964  // The classname is too large, try again with a large buffer.
965  fEntryBuffer->SetBufferOffset(name_start);
966  maxsize *= 2;
967  delete[] s;
968  s = new char[maxsize];
969  fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end
970  }
971  delete[] s;
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 
1058 TLeaf* 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  }
1134  Int_t len = fBranches.GetEntriesFast();
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 
1203 TBasket* 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(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 
1301 const char * TBranch::GetClassName() const
1302 {
1303  return "";
1304 }
1305 
1306 ////////////////////////////////////////////////////////////////////////////////
1307 /// Return icon name depending on type of branch.
1308 
1309 const 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.
1360  basket = (TBasket*) fBaskets.UncheckedAt(fReadBasket);
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;
1439  Long64_t first;
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;
1500  Long64_t first;
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.
1588  Long64_t first;
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.
1605  if (!TestBit(kDoNotUseBufferMap)) {
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.
1677  if (!TestBit(kDoNotUseBufferMap)) {
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 
1708 Int_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  {
1748  TDirectory::TContext ctxt;
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 
1766 {
1767  TBasket *basket = 0;
1768  if (user_buffer && fExtraBasket) {
1769  basket = fExtraBasket;
1770  fExtraBasket = nullptr;
1771  basket->AdoptBuffer(user_buffer);
1772  } else {
1773  if (GetTree()->MemoryFull(0)) {
1774  if (fNBaskets==1) {
1775  // Steal the existing basket
1776  Int_t oldindex = fBaskets.GetLast();
1777  basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1778  if (!basket) {
1779  fBaskets.SetLast(-2); // For recalculation of Last.
1780  oldindex = fBaskets.GetLast();
1781  if (oldindex != fBaskets.LowerBound()-1) {
1782  basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1783  }
1784  }
1785  if (basket && fBasketBytes[oldindex]!=0) {
1786  if (basket == fCurrentBasket) {
1787  fCurrentBasket = 0;
1788  fFirstBasketEntry = -1;
1789  fNextBasketEntry = -1;
1790  }
1791  fBaskets.AddAt(0,oldindex);
1792  fBaskets.SetLast(-1);
1793  fNBaskets = 0;
1794  } else {
1795  basket = fTree->CreateBasket(this);
1796  }
1797  } else if (fNBaskets == 0) {
1798  // There is nothing to drop!
1799  basket = fTree->CreateBasket(this);
1800  } else {
1801  // Memory is full and there is more than one basket,
1802  // Let DropBaskets do it job.
1803  DropBaskets();
1804  basket = fTree->CreateBasket(this);
1805  }
1806  } else {
1807  basket = fTree->CreateBasket(this);
1808  }
1809  if (user_buffer)
1810  basket->AdoptBuffer(user_buffer);
1811  }
1812  return basket;
1813 }
1814 
1815 ////////////////////////////////////////////////////////////////////////////////
1816 /// Drops the cluster two behind the current cluster and returns a fresh basket
1817 /// by either reusing or creating a new one
1818 
1820 {
1821  TBasket *basket = 0;
1822 
1823  // If GetClusterIterator is called with a negative entry then GetStartEntry will be 0
1824  // So we need to check if we reach the zero before we have gone back (1-VirtualSize) clusters
1825  // if this is the case, we want to keep everything in memory so we return a new basket
1827  if (iter.GetStartEntry() == 0) {
1828  return fTree->CreateBasket(this);
1829  }
1830 
1831  // Iterate backwards (1-VirtualSize) clusters to reach cluster to be unloaded from memory,
1832  // skipped if VirtualSize > 0.
1833  for (Int_t j = 0; j < -fTree->GetMaxVirtualSize(); j++) {
1834  if (iter.Previous() == 0) {
1835  return fTree->CreateBasket(this);
1836  }
1837  }
1838 
1839  Int_t entryToUnload = iter.Previous();
1840  // Finds the basket to unload from memory. Since the basket should be close to current
1841  // basket, just iterate backwards until the correct basket is reached. This should
1842  // be fast as long as the number of baskets per cluster is small
1843  Int_t basketToUnload = fReadBasket;
1844  while (fBasketEntry[basketToUnload] != entryToUnload) {
1845  basketToUnload--;
1846  if (basketToUnload < 0) {
1847  return fTree->CreateBasket(this);
1848  }
1849  }
1850 
1851  // Retrieves the basket that is going to be unloaded from memory. If the basket did not
1852  // exist, create a new one
1853  basket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1854  if (basket) {
1855  fBaskets.AddAt(0, basketToUnload);
1856  --fNBaskets;
1857  } else {
1858  basket = fTree->CreateBasket(this);
1859  }
1860  ++basketToUnload;
1861 
1862  // Clear the rest of the baskets. While it would be ideal to reuse these baskets
1863  // for other baskets in the new cluster. It would require the function to go
1864  // beyond its current scope. In the ideal case when each cluster only has 1 basket
1865  // this will perform well
1866  iter.Next();
1867  while (fBasketEntry[basketToUnload] < iter.GetStartEntry()) {
1868  TBasket *oldbasket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1869  if (oldbasket) {
1870  oldbasket->DropBuffers();
1871  delete oldbasket;
1872  fBaskets.AddAt(0, basketToUnload);
1873  --fNBaskets;
1874  }
1875  ++basketToUnload;
1876  }
1877  fBaskets.SetLast(-1);
1878  return basket;
1879 }
1880 
1881 ////////////////////////////////////////////////////////////////////////////////
1882 /// Return pointer to the 1st Leaf named name in thisBranch
1883 
1884 TLeaf* TBranch::GetLeaf(const char* name) const
1885 {
1886  Int_t i;
1887  for (i=0;i<fNleaves;i++) {
1888  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
1889  if (!strcmp(leaf->GetName(),name)) return leaf;
1890  }
1891  return 0;
1892 }
1893 
1894 ////////////////////////////////////////////////////////////////////////////////
1895 /// Get real file name
1896 
1898 {
1899  if (fFileName.Length()==0) {
1900  return fFileName;
1901  }
1902  TString bFileName = fFileName;
1903 
1904  // check if branch file name is absolute or a URL (e.g. root://host/...)
1905  char *bname = gSystem->ExpandPathName(fFileName.Data());
1906  if (!gSystem->IsAbsoluteFileName(bname) && !strstr(bname, ":/") && fTree && fTree->GetCurrentFile()) {
1907 
1908  // if not, get filename where tree header is stored
1909  const char *tfn = fTree->GetCurrentFile()->GetName();
1910 
1911  // If it is an archive file we need a special treatment
1912  TUrl arc(tfn);
1913  if (strlen(arc.GetAnchor()) > 0) {
1915  bFileName = arc.GetUrl();
1916  } else {
1917  // if this is an absolute path or a URL then prepend this path
1918  // to the branch file name
1919  char *tname = gSystem->ExpandPathName(tfn);
1920  if (gSystem->IsAbsoluteFileName(tname) || strstr(tname, ":/")) {
1921  bFileName = gSystem->GetDirName(tname);
1922  bFileName += "/";
1923  bFileName += fFileName;
1924  }
1925  delete [] tname;
1926  }
1927  }
1928  delete [] bname;
1929 
1930  return bFileName;
1931 }
1932 
1933 ////////////////////////////////////////////////////////////////////////////////
1934 /// Return all elements of one row unpacked in internal array fValues
1935 /// [Actually just returns 1 (?)]
1936 
1938 {
1939  return 1;
1940 }
1941 
1942 ////////////////////////////////////////////////////////////////////////////////
1943 /// Return whether this branch is in a mode where the object are decomposed
1944 /// or not (Also known as MakeClass mode).
1945 
1947 {
1948  // Regular TBranch and TBrancObject can not be in makeClass mode
1949 
1950  return kFALSE;
1951 }
1952 
1953 ////////////////////////////////////////////////////////////////////////////////
1954 /// Get our top-level parent branch in the tree.
1955 
1957 {
1958  if (fMother) return fMother;
1959 
1960  const TObjArray* array = fTree->GetListOfBranches();
1961  Int_t n = array->GetEntriesFast();
1962  for (Int_t i = 0; i < n; ++i) {
1963  TBranch* branch = (TBranch*) array->UncheckedAt(i);
1964  TBranch* parent = branch->GetSubBranch(this);
1965  if (parent) {
1966  const_cast<TBranch*>(this)->fMother = branch; // We can not yet use the 'mutable' keyword
1967  return branch;
1968  }
1969  }
1970  return 0;
1971 }
1972 
1973 ////////////////////////////////////////////////////////////////////////////////
1974 /// Find the parent branch of child.
1975 /// Return 0 if child is not in this branch hierarchy.
1976 
1978 {
1979  // Handle error condition, if the parameter is us, we cannot find the parent.
1980  if (this == child) {
1981  // Note: We cast away any const-ness of "this".
1982  return (TBranch*) this;
1983  }
1984 
1985  if (child->fParent) {
1986  return child->fParent;
1987  }
1988 
1989  Int_t len = fBranches.GetEntriesFast();
1990  for (Int_t i = 0; i < len; ++i) {
1991  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1992  if (!branch) {
1993  continue;
1994  }
1995  if (branch == child) {
1996  // We are the direct parent of child.
1997  const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
1998  // Note: We cast away any const-ness of "this".
1999  const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
2000  return (TBranch*) this;
2001  }
2002  // FIXME: This is a tail-recursion!
2003  TBranch* parent = branch->GetSubBranch(child);
2004  if (parent) {
2005  return parent;
2006  }
2007  }
2008  // We failed to find the parent.
2009  return 0;
2010 }
2011 
2012 ////////////////////////////////////////////////////////////////////////////////
2013 /// Return total number of bytes in the branch (including current buffer)
2014 
2016 {
2017  TBufferFile b(TBuffer::kWrite, 10000);
2018  // This intentionally only store the TBranch part and thus slightly
2019  // under-estimate the space used.
2020  // Since the TBranchElement part contains pointers to other branches (branch count),
2021  // doing regular Streaming would end up including those and thus greatly over-estimate
2022  // the size used.
2023  const_cast<TBranch *>(this)->TBranch::Streamer(b);
2024 
2025  Long64_t totbytes = 0;
2026  if (fZipBytes > 0) totbytes = fTotBytes;
2027  return totbytes + b.Length();
2028 }
2029 
2030 ////////////////////////////////////////////////////////////////////////////////
2031 /// Return total number of bytes in the branch (excluding current buffer)
2032 /// if option ="*" includes all sub-branches of this branch too
2033 
2035 {
2036  Long64_t totbytes = fTotBytes;
2037  if (!option) return totbytes;
2038  if (option[0] != '*') return totbytes;
2039  //scan sub-branches
2040  Int_t len = fBranches.GetEntriesFast();
2041  for (Int_t i = 0; i < len; ++i) {
2042  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2043  if (branch) totbytes += branch->GetTotBytes(option);
2044  }
2045  return totbytes;
2046 }
2047 
2048 ////////////////////////////////////////////////////////////////////////////////
2049 /// Return total number of zip bytes in the branch
2050 /// if option ="*" includes all sub-branches of this branch too
2051 
2053 {
2054  Long64_t zipbytes = fZipBytes;
2055  if (!option) return zipbytes;
2056  if (option[0] != '*') return zipbytes;
2057  //scan sub-branches
2058  Int_t len = fBranches.GetEntriesFast();
2059  for (Int_t i = 0; i < len; ++i) {
2060  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2061  if (branch) zipbytes += branch->GetZipBytes(option);
2062  }
2063  return zipbytes;
2064 }
2065 
2066 ////////////////////////////////////////////////////////////////////////////////
2067 /// Returns the IO settings currently in use for this branch.
2068 
2070 {
2071  return fIOFeatures;
2072 }
2073 
2074 ////////////////////////////////////////////////////////////////////////////////
2075 /// Return kTRUE if an existing object in a TBranchObject must be deleted.
2076 
2078 {
2079  return TestBit(kAutoDelete);
2080 }
2081 
2082 ////////////////////////////////////////////////////////////////////////////////
2083 /// Return kTRUE if more than one leaf or browsables, kFALSE otherwise.
2084 
2086 {
2087  if (fNleaves > 1) {
2088  return kTRUE;
2089  }
2090  TList* browsables = const_cast<TBranch*>(this)->GetBrowsables();
2091  return browsables && browsables->GetSize();
2092 }
2093 
2094 ////////////////////////////////////////////////////////////////////////////////
2095 /// keep a maximum of fMaxEntries in memory
2096 
2098 {
2099  Int_t dentries = (Int_t) (fEntries - maxEntries);
2100  TBasket* basket = (TBasket*) fBaskets.UncheckedAt(0);
2101  if (basket) basket->MoveEntries(dentries);
2102  fEntries = maxEntries;
2103  fEntryNumber = maxEntries;
2104  //loop on sub branches
2106  for (Int_t i = 0; i < nb; ++i) {
2107  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
2108  branch->KeepCircular(maxEntries);
2109  }
2110 }
2111 
2112 ////////////////////////////////////////////////////////////////////////////////
2113 /// Baskets associated to this branch are forced to be in memory.
2114 /// You can call TTree::SetMaxVirtualSize(maxmemory) to instruct
2115 /// the system that the total size of the imported baskets does not
2116 /// exceed maxmemory bytes.
2117 ///
2118 /// The function returns the number of baskets that have been put in memory.
2119 /// This method may be called to force all baskets of one or more branches
2120 /// in memory when random access to entries in this branch is required.
2121 /// See also TTree::LoadBaskets to load all baskets of all branches in memory.
2122 
2124 {
2125  Int_t nimported = 0;
2126  Int_t nbaskets = fWriteBasket;
2127  TFile *file = GetFile(0);
2128  if (!file) return 0;
2129  TBasket *basket;
2130  for (Int_t i=0;i<nbaskets;i++) {
2131  basket = (TBasket*)fBaskets.UncheckedAt(i);
2132  if (basket) continue;
2133  basket = GetFreshBasket(nullptr);
2134  if (fBasketBytes[i] == 0) {
2135  fBasketBytes[i] = basket->ReadBasketBytes(fBasketSeek[i],file);
2136  }
2137  Int_t badread = basket->ReadBasketBuffers(fBasketSeek[i],fBasketBytes[i],file);
2138  if (badread) {
2139  Error("Loadbaskets","Error while reading basket buffer %d of branch %s",i,GetName());
2140  return -1;
2141  }
2142  ++fNBaskets;
2143  fBaskets.AddAt(basket,i);
2144  nimported++;
2145  }
2146  return nimported;
2147 }
2148 
2149 ////////////////////////////////////////////////////////////////////////////////
2150 /// Print TBranch parameters
2151 ///
2152 /// If options contains "basketsInfo" print the entry number, location and size
2153 /// of each baskets.
2154 
2155 void TBranch::Print(Option_t *option) const
2156 {
2157  const int kLINEND = 77;
2158  Float_t cx = 1;
2159 
2160  TString titleContent(GetTitle());
2161  if ( titleContent == GetName() ) {
2162  titleContent.Clear();
2163  }
2164 
2165  if (fLeaves.GetEntries() == 1) {
2166  if (titleContent.Length()>=2 && titleContent[titleContent.Length()-2]=='/' && isalpha(titleContent[titleContent.Length()-1])) {
2167  // The type is already encoded. Nothing to do.
2168  } else {
2169  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(0);
2170  if (titleContent.Length()) {
2171  titleContent.Prepend(" ");
2172  }
2173  // titleContent.Append("type: ");
2174  titleContent.Prepend(leaf->GetTypeName());
2175  }
2176  }
2177  Int_t titleLength = titleContent.Length();
2178 
2179  Int_t aLength = titleLength + strlen(GetName());
2180  aLength += (aLength / 54 + 1) * 80 + 100;
2181  if (aLength < 200) aLength = 200;
2182  char *bline = new char[aLength];
2183 
2184  Long64_t totBytes = GetTotalSize();
2185  if (fZipBytes) cx = (fTotBytes+0.00001)/fZipBytes;
2186  if (titleLength) snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName(),titleContent.Data());
2187  else snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName()," ");
2188  if (strlen(bline) > UInt_t(kLINEND)) {
2189  char *tmp = new char[strlen(bline)+1];
2190  if (titleLength) strlcpy(tmp, titleContent.Data(),strlen(bline)+1);
2191  snprintf(bline,aLength,"*Br%5d :%-9s : ",fgCount,GetName());
2192  int pos = strlen (bline);
2193  int npos = pos;
2194  int beg=0, end;
2195  while (beg < titleLength) {
2196  for (end=beg+1; end < titleLength-1; end ++)
2197  if (tmp[end] == ':') break;
2198  if (npos + end-beg+1 >= 78) {
2199  while (npos < kLINEND) {
2200  bline[pos ++] = ' ';
2201  npos ++;
2202  }
2203  bline[pos ++] = '*';
2204  bline[pos ++] = '\n';
2205  bline[pos ++] = '*';
2206  npos = 1;
2207  for (; npos < 12; npos ++)
2208  bline[pos ++] = ' ';
2209  bline[pos-2] = '|';
2210  }
2211  for (int n = beg; n <= end; n ++)
2212  bline[pos+n-beg] = tmp[n];
2213  pos += end-beg+1;
2214  npos += end-beg+1;
2215  beg = end+1;
2216  }
2217  while (npos < kLINEND) {
2218  bline[pos ++] = ' ';
2219  npos ++;
2220  }
2221  bline[pos ++] = '*';
2222  bline[pos] = '\0';
2223  delete[] tmp;
2224  }
2225  Printf("%s", bline);
2226 
2227  if (fTotBytes > 2000000000) {
2228  Printf("*Entries :%lld : Total Size=%11lld bytes File Size = %lld *",fEntries,totBytes,fZipBytes);
2229  } else {
2230  if (fZipBytes > 0) {
2231  Printf("*Entries :%9lld : Total Size=%11lld bytes File Size = %10lld *",fEntries,totBytes,fZipBytes);
2232  } else {
2233  if (fWriteBasket > 0) {
2234  Printf("*Entries :%9lld : Total Size=%11lld bytes All baskets in memory *",fEntries,totBytes);
2235  } else {
2236  Printf("*Entries :%9lld : Total Size=%11lld bytes One basket in memory *",fEntries,totBytes);
2237  }
2238  }
2239  }
2240  Printf("*Baskets :%9d : Basket Size=%11d bytes Compression= %6.2f *",fWriteBasket,fBasketSize,cx);
2241 
2242  if (strncmp(option,"basketsInfo",strlen("basketsInfo"))==0) {
2243  Int_t nbaskets = fWriteBasket;
2244  for (Int_t i=0;i<nbaskets;i++) {
2245  Printf("*Basket #%4d entry=%6lld pos=%6lld size=%5d",
2246  i, fBasketEntry[i], fBasketSeek[i], fBasketBytes[i]);
2247  }
2248  }
2249 
2250  Printf("*............................................................................*");
2251  delete [] bline;
2252  fgCount++;
2253 }
2254 
2255 ////////////////////////////////////////////////////////////////////////////////
2256 /// Print the information we have about which basket is currently cached and
2257 /// whether they have been 'used'/'read' from the cache.
2258 
2260 {
2262 }
2263 
2264 ////////////////////////////////////////////////////////////////////////////////
2265 /// Loop on all leaves of this branch to read Basket buffer.
2266 
2268 {
2269  // fLeaves->ReadBasket(basket);
2270 }
2271 
2272 ////////////////////////////////////////////////////////////////////////////////
2273 /// Loop on all leaves of this branch to read Basket buffer.
2274 
2276 {
2277  for (Int_t i = 0; i < fNleaves; ++i) {
2278  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2279  leaf->ReadBasket(b);
2280  }
2281 }
2282 
2283 ////////////////////////////////////////////////////////////////////////////////
2284 /// Read zero leaves without the overhead of a loop.
2285 
2287 {
2288 }
2289 
2290 ////////////////////////////////////////////////////////////////////////////////
2291 /// Read one leaf without the overhead of a loop.
2292 
2294 {
2295  ((TLeaf*) fLeaves.UncheckedAt(0))->ReadBasket(b);
2296 }
2297 
2298 ////////////////////////////////////////////////////////////////////////////////
2299 /// Read two leaves without the overhead of a loop.
2300 
2302 {
2303  ((TLeaf*) fLeaves.UncheckedAt(0))->ReadBasket(b);
2304  ((TLeaf*) fLeaves.UncheckedAt(1))->ReadBasket(b);
2305 }
2306 
2307 ////////////////////////////////////////////////////////////////////////////////
2308 /// Loop on all leaves of this branch to fill Basket buffer.
2309 
2311 {
2312  for (Int_t i = 0; i < fNleaves; ++i) {
2313  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2314  leaf->FillBasket(b);
2315  }
2316 }
2317 
2318 ////////////////////////////////////////////////////////////////////////////////
2319 /// Refresh this branch using new information in b
2320 /// This function is called by TTree::Refresh
2321 
2323 {
2324  if (b==0) return;
2325 
2326  fEntryOffsetLen = b->fEntryOffsetLen;
2327  fWriteBasket = b->fWriteBasket;
2328  fEntryNumber = b->fEntryNumber;
2329  fMaxBaskets = b->fMaxBaskets;
2330  fEntries = b->fEntries;
2331  fTotBytes = b->fTotBytes;
2332  fZipBytes = b->fZipBytes;
2333  fReadBasket = 0;
2334  fReadEntry = -1;
2335  fFirstBasketEntry = -1;
2336  fNextBasketEntry = -1;
2337  fCurrentBasket = 0;
2338  delete [] fBasketBytes;
2339  delete [] fBasketEntry;
2340  delete [] fBasketSeek;
2344  Int_t i;
2345  for (i=0;i<fMaxBaskets;i++) {
2346  fBasketBytes[i] = b->fBasketBytes[i];
2347  fBasketEntry[i] = b->fBasketEntry[i];
2348  fBasketSeek[i] = b->fBasketSeek[i];
2349  }
2350  fBaskets.Delete();
2351  Int_t nbaskets = b->fBaskets.GetSize();
2352  fBaskets.Expand(nbaskets);
2353  // If the current fWritebasket is in memory, take it (just swap)
2354  // from the Tree being read
2355  TBasket *basket = (TBasket*)b->fBaskets.UncheckedAt(fWriteBasket);
2356  fBaskets.AddAt(basket,fWriteBasket);
2357  if (basket) {
2358  fNBaskets = 1;
2359  --(b->fNBaskets);
2360  b->fBaskets.RemoveAt(fWriteBasket);
2361  basket->SetBranch(this);
2362  }
2363 }
2364 
2365 ////////////////////////////////////////////////////////////////////////////////
2366 /// Reset a Branch.
2367 ///
2368 /// - Existing buffers are deleted.
2369 /// - Entries, max and min are reset.
2370 
2372 {
2373  fReadBasket = 0;
2374  fReadEntry = -1;
2375  fFirstBasketEntry = -1;
2376  fNextBasketEntry = -1;
2377  fCurrentBasket = 0;
2378  fWriteBasket = 0;
2379  fEntries = 0;
2380  fTotBytes = 0;
2381  fZipBytes = 0;
2382  fEntryNumber = 0;
2383 
2384  if (fBasketBytes) {
2385  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2386  fBasketBytes[i] = 0;
2387  }
2388  }
2389 
2390  if (fBasketEntry) {
2391  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2392  fBasketEntry[i] = 0;
2393  }
2394  }
2395 
2396  if (fBasketSeek) {
2397  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2398  fBasketSeek[i] = 0;
2399  }
2400  }
2401 
2402  fBaskets.Delete();
2403  fNBaskets = 0;
2404 }
2405 
2406 ////////////////////////////////////////////////////////////////////////////////
2407 /// Reset a Branch.
2408 ///
2409 /// - Existing buffers are deleted.
2410 /// - Entries, max and min are reset.
2411 
2413 {
2414  fReadBasket = 0;
2415  fReadEntry = -1;
2416  fFirstBasketEntry = -1;
2417  fNextBasketEntry = -1;
2418  fCurrentBasket = 0;
2419  fWriteBasket = 0;
2420  fEntries = 0;
2421  fTotBytes = 0;
2422  fZipBytes = 0;
2423  fEntryNumber = 0;
2424 
2425  if (fBasketBytes) {
2426  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2427  fBasketBytes[i] = 0;
2428  }
2429  }
2430 
2431  if (fBasketEntry) {
2432  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2433  fBasketEntry[i] = 0;
2434  }
2435  }
2436 
2437  if (fBasketSeek) {
2438  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2439  fBasketSeek[i] = 0;
2440  }
2441  }
2442 
2443  TBasket *reusebasket = (TBasket*)fBaskets[fWriteBasket];
2444  if (reusebasket) {
2445  fBaskets[fWriteBasket] = 0;
2446  } else {
2447  reusebasket = (TBasket*)fBaskets[fReadBasket];
2448  if (reusebasket) {
2449  fBaskets[fReadBasket] = 0;
2450  }
2451  }
2452  fBaskets.Delete();
2453  if (reusebasket) {
2454  fNBaskets = 1;
2455  reusebasket->Reset();
2456  fBaskets[0] = reusebasket;
2457  } else {
2458  fNBaskets = 0;
2459  }
2460 }
2461 
2462 ////////////////////////////////////////////////////////////////////////////////
2463 /// Reset the address of the branch.
2464 
2466 {
2467  fAddress = 0;
2468 
2469  // Reset last read entry number, we have will had new user object now.
2470  fReadEntry = -1;
2471 
2472  for (Int_t i = 0; i < fNleaves; ++i) {
2473  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2474  leaf->SetAddress(0);
2475  }
2476 
2477  Int_t nbranches = fBranches.GetEntriesFast();
2478  for (Int_t i = 0; i < nbranches; ++i) {
2479  TBranch* abranch = (TBranch*) fBranches[i];
2480  // FIXME: This is a tail recursion.
2481  abranch->ResetAddress();
2482  }
2483 }
2484 
2485 ////////////////////////////////////////////////////////////////////////////////
2486 /// Static function resetting fgCount
2487 
2489 {
2490  fgCount = 0;
2491 }
2492 
2493 ////////////////////////////////////////////////////////////////////////////////
2494 /// Set address of this branch.
2495 
2496 void TBranch::SetAddress(void* addr)
2497 {
2498  if (TestBit(kDoNotProcess)) {
2499  return;
2500  }
2501  fReadEntry = -1;
2502  fFirstBasketEntry = -1;
2503  fNextBasketEntry = -1;
2504  fAddress = (char*) addr;
2505  for (Int_t i = 0; i < fNleaves; ++i) {
2506  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2507  Int_t offset = leaf->GetOffset();
2508  if (TestBit(kIsClone)) {
2509  offset = 0;
2510  }
2511  if (fAddress) leaf->SetAddress(fAddress + offset);
2512  else leaf->SetAddress(0);
2513  }
2514 }
2515 
2516 ////////////////////////////////////////////////////////////////////////////////
2517 /// Set the automatic delete bit.
2518 ///
2519 /// This bit is used by TBranchObject::ReadBasket to decide if an object
2520 /// referenced by a TBranchObject must be deleted or not before reading
2521 /// a new entry.
2522 ///
2523 /// If autodel is kTRUE, this existing object will be deleted, a new object
2524 /// created by the default constructor, then read from disk by the streamer.
2525 ///
2526 /// If autodel is kFALSE, the existing object is not deleted. Root assumes
2527 /// that the user is taking care of deleting any internal object or array
2528 /// (this can be done in the streamer).
2529 
2531 {
2532  if (autodel) {
2533  SetBit(kAutoDelete, 1);
2534  } else {
2535  SetBit(kAutoDelete, 0);
2536  }
2537 }
2538 
2539 ////////////////////////////////////////////////////////////////////////////////
2540 /// Set the basket size
2541 /// The function makes sure that the basket size is greater than fEntryOffsetlen
2542 
2544 {
2545  Int_t minsize = 100 + fName.Length();
2546  if (buffsize < minsize+fEntryOffsetLen) buffsize = minsize+fEntryOffsetLen;
2547  fBasketSize = buffsize;
2548  TBasket *basket = (TBasket*)fBaskets[fWriteBasket];
2549  if (basket) {
2550  basket->AdjustSize(fBasketSize);
2551  }
2552 }
2553 
2554 ////////////////////////////////////////////////////////////////////////////////
2555 /// Set address of this branch directly from a TBuffer to avoid streaming.
2556 ///
2557 /// Note: We do not take ownership of the buffer.
2558 
2560 {
2561  // Check this is possible
2562  if ( (fNleaves != 1)
2563  || (strcmp("TLeafObject",fLeaves.UncheckedAt(0)->ClassName())!=0) ) {
2564  Error("TBranch::SetAddress","Filling from a TBuffer can only be done with a not split object branch. Request ignored.");
2565  } else {
2566  fReadEntry = -1;
2567  fNextBasketEntry = -1;
2568  fFirstBasketEntry = -1;
2569  // Note: We do not take ownership of the buffer.
2570  fEntryBuffer = buf;
2571  }
2572 }
2573 
2574 ////////////////////////////////////////////////////////////////////////////////
2575 /// Set compression algorithm.
2576 
2578 {
2579  if (algorithm < 0 || algorithm >= ROOT::RCompressionSetting::EAlgorithm::kUndefined) algorithm = 0;
2580  if (fCompress < 0) {
2582  } else {
2583  int level = fCompress % 100;
2584  fCompress = 100 * algorithm + level;
2585  }
2586 
2588  for (Int_t i=0;i<nb;i++) {
2589  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2590  branch->SetCompressionAlgorithm(algorithm);
2591  }
2592 }
2593 
2594 ////////////////////////////////////////////////////////////////////////////////
2595 /// Set compression level.
2596 
2598 {
2599  if (level < 0) level = 0;
2600  if (level > 99) level = 99;
2601  if (fCompress < 0) {
2602  fCompress = level;
2603  } else {
2604  int algorithm = fCompress / 100;
2605  if (algorithm >= ROOT::RCompressionSetting::EAlgorithm::kUndefined) algorithm = 0;
2606  fCompress = 100 * algorithm + level;
2607  }
2608 
2610  for (Int_t i=0;i<nb;i++) {
2611  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2612  branch->SetCompressionLevel(level);
2613  }
2614 }
2615 
2616 ////////////////////////////////////////////////////////////////////////////////
2617 /// Set compression settings.
2618 
2620 {
2621  fCompress = settings;
2622 
2624  for (Int_t i=0;i<nb;i++) {
2625  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2626  branch->SetCompressionSettings(settings);
2627  }
2628 }
2629 
2630 ////////////////////////////////////////////////////////////////////////////////
2631 /// Update the default value for the branch's fEntryOffsetLen if and only if
2632 /// it was already non zero (and the new value is not zero)
2633 /// If updateExisting is true, also update all the existing branches.
2634 
2635 void TBranch::SetEntryOffsetLen(Int_t newdefault, Bool_t updateExisting)
2636 {
2637  if (fEntryOffsetLen && newdefault) {
2638  fEntryOffsetLen = newdefault;
2639  }
2640  if (updateExisting) {
2641  TIter next( GetListOfBranches() );
2642  TBranch *b;
2643  while ( ( b = (TBranch*)next() ) ) {
2644  b->SetEntryOffsetLen( newdefault, kTRUE );
2645  }
2646  }
2647 }
2648 
2649 ////////////////////////////////////////////////////////////////////////////////
2650 /// Set the number of entries in this branch.
2651 
2653 {
2654  fEntries = entries;
2655  fEntryNumber = entries;
2656 }
2657 
2658 ////////////////////////////////////////////////////////////////////////////////
2659 /// Set file where this branch writes/reads its buffers.
2660 /// By default the branch buffers reside in the file where the
2661 /// Tree was created.
2662 /// If the file name where the tree was created is an absolute
2663 /// path name or an URL (e.g. or root://host/...)
2664 /// and if the fname is not an absolute path name or an URL then
2665 /// the path of the tree file is prepended to fname to make the
2666 /// branch file relative to the tree file. In this case one can
2667 /// move the tree + all branch files to a different location in
2668 /// the file system and still access the branch files.
2669 /// The ROOT file will be connected only when necessary.
2670 /// If called by TBranch::Fill (via TBasket::WriteFile), the file
2671 /// will be created with the option "recreate".
2672 /// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2673 /// will be opened in read mode.
2674 /// To open a file in "update" mode or with a certain compression
2675 /// level, use TBranch::SetFile(TFile *file).
2676 
2678 {
2679  if (file == 0) file = fTree->GetCurrentFile();
2681  if (file == fTree->GetCurrentFile()) fFileName = "";
2682  else fFileName = file->GetName();
2683 
2684  if (file && fCompress == -1) {
2685  fCompress = file->GetCompressionLevel();
2686  }
2687 
2688  // Apply to all existing baskets.
2689  TIter nextb(GetListOfBaskets());
2690  TBasket *basket;
2691  while ((basket = (TBasket*)nextb())) {
2692  basket->SetParent(file);
2693  }
2694 
2695  // Apply to sub-branches as well.
2696  TIter next(GetListOfBranches());
2697  TBranch *branch;
2698  while ((branch = (TBranch*)next())) {
2699  branch->SetFile(file);
2700  }
2701 }
2702 
2703 ////////////////////////////////////////////////////////////////////////////////
2704 /// Set file where this branch writes/reads its buffers.
2705 /// By default the branch buffers reside in the file where the
2706 /// Tree was created.
2707 /// If the file name where the tree was created is an absolute
2708 /// path name or an URL (e.g. root://host/...)
2709 /// and if the fname is not an absolute path name or an URL then
2710 /// the path of the tree file is prepended to fname to make the
2711 /// branch file relative to the tree file. In this case one can
2712 /// move the tree + all branch files to a different location in
2713 /// the file system and still access the branch files.
2714 /// The ROOT file will be connected only when necessary.
2715 /// If called by TBranch::Fill (via TBasket::WriteFile), the file
2716 /// will be created with the option "recreate".
2717 /// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2718 /// will be opened in read mode.
2719 /// To open a file in "update" mode or with a certain compression
2720 /// level, use TBranch::SetFile(TFile *file).
2721 
2722 void TBranch::SetFile(const char* fname)
2723 {
2724  fFileName = fname;
2725  fDirectory = 0;
2726 
2727  //apply to sub-branches as well
2728  TIter next(GetListOfBranches());
2729  TBranch *branch;
2730  while ((branch = (TBranch*)next())) {
2731  branch->SetFile(fname);
2732  }
2733 }
2734 
2735 ////////////////////////////////////////////////////////////////////////////////
2736 /// Set the branch in a mode where the object are decomposed
2737 /// (Also known as MakeClass mode).
2738 /// Return whether the setting was possible (it is not possible for
2739 /// TBranch and TBranchObject).
2740 
2742 {
2743  // Regular TBranch and TBrancObject can not be in makeClass mode
2744  return kFALSE;
2745 }
2746 
2747 ////////////////////////////////////////////////////////////////////////////////
2748 /// Set object this branch is pointing to.
2749 
2750 void TBranch::SetObject(void * /* obj */)
2751 {
2752  if (TestBit(kDoNotProcess)) {
2753  return;
2754  }
2755  Warning("SetObject","is not supported in TBranch objects");
2756 }
2757 
2758 ////////////////////////////////////////////////////////////////////////////////
2759 /// Set branch status to Process or DoNotProcess.
2760 
2762 {
2763  if (status) ResetBit(kDoNotProcess);
2764  else SetBit(kDoNotProcess);
2765 }
2766 
2767 ////////////////////////////////////////////////////////////////////////////////
2768 /// Stream a class object
2769 
2770 void TBranch::Streamer(TBuffer& b)
2771 {
2772  if (b.IsReading()) {
2773  UInt_t R__s, R__c;
2774  fTree = 0; // Will be set by TTree::Streamer
2775  fAddress = 0;
2776  gROOT->SetReadingObject(kTRUE);
2777 
2778  // Reset transients.
2780  fCurrentBasket = 0;
2781  fFirstBasketEntry = -1;
2782  fNextBasketEntry = -1;
2783 
2784  Version_t v = b.ReadVersion(&R__s, &R__c);
2785  if (v > 9) {
2786  b.ReadClassBuffer(TBranch::Class(), this, v, R__s, R__c);
2787 
2788  if (fWriteBasket>=fBaskets.GetSize()) {
2790  }
2791  fDirectory = 0;
2793  for (Int_t i=0;i<fNleaves;i++) {
2794  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2795  leaf->SetBranch(this);
2796  }
2797 
2799  for (Int_t j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2800  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2801  if (bk) {
2802  bk->SetBranch(this);
2803  // GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2804  ++n;
2805  }
2806  }
2807  if (fWriteBasket >= fMaxBaskets) {
2808  //old versions may need this fix
2813 
2814  }
2816  gROOT->SetReadingObject(kFALSE);
2817  if (IsA() == TBranch::Class()) {
2818  if (fNleaves == 0) {
2820  } else if (fNleaves == 1) {
2822  } else if (fNleaves == 2) {
2824  } else {
2826  }
2827  }
2828  return;
2829  }
2830  //====process old versions before automatic schema evolution
2831  Int_t n,i,j,ijunk;
2832  if (v > 5) {
2833  Stat_t djunk;
2834  TNamed::Streamer(b);
2835  if (v > 7) TAttFill::Streamer(b);
2836  b >> fCompress;
2837  b >> fBasketSize;
2838  b >> fEntryOffsetLen;
2839  b >> fWriteBasket;
2840  b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2841  b >> fOffset;
2842  b >> fMaxBaskets;
2843  if (v > 6) b >> fSplitLevel;
2844  b >> djunk; fEntries = (Long64_t)djunk;
2845  b >> djunk; fTotBytes = (Long64_t)djunk;
2846  b >> djunk; fZipBytes = (Long64_t)djunk;
2847 
2848  fBranches.Streamer(b);
2849  fLeaves.Streamer(b);
2850  fBaskets.Streamer(b);
2854  Char_t isArray;
2855  b >> isArray;
2856  b.ReadFastArray(fBasketBytes,fMaxBaskets);
2857  b >> isArray;
2858  for (i=0;i<fMaxBaskets;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2859  b >> isArray;
2860  for (i=0;i<fMaxBaskets;i++) {
2861  if (isArray == 2) b >> fBasketSeek[i];
2862  else {Int_t bsize; b >> bsize; fBasketSeek[i] = (Long64_t)bsize;};
2863  }
2864  fFileName.Streamer(b);
2865  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2866  fDirectory = 0;
2868  for (i=0;i<fNleaves;i++) {
2869  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2870  leaf->SetBranch(this);
2871  }
2873  for (j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2874  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2875  if (bk) {
2876  bk->SetBranch(this);
2877  //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2878  ++n;
2879  }
2880  }
2881  if (fWriteBasket >= fMaxBaskets) {
2882  //old versions may need this fix
2887 
2888  }
2889  // Check Byte Count is not needed since it was done in ReadBuffer
2891  gROOT->SetReadingObject(kFALSE);
2892  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2893  if (IsA() == TBranch::Class()) {
2894  if (fNleaves == 0) {
2896  } else if (fNleaves == 1) {
2898  } else if (fNleaves == 2) {
2900  } else {
2902  }
2903  }
2904  return;
2905  }
2906  //====process very old versions
2907  Stat_t djunk;
2908  TNamed::Streamer(b);
2909  b >> fCompress;
2910  b >> fBasketSize;
2911  b >> fEntryOffsetLen;
2912  b >> fMaxBaskets;
2913  b >> fWriteBasket;
2914  b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2915  b >> djunk; fEntries = (Long64_t)djunk;
2916  b >> djunk; fTotBytes = (Long64_t)djunk;
2917  b >> djunk; fZipBytes = (Long64_t)djunk;
2918  b >> fOffset;
2919  fBranches.Streamer(b);
2920  fLeaves.Streamer(b);
2922  for (i=0;i<fNleaves;i++) {
2923  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2924  leaf->SetBranch(this);
2925  }
2926  fBaskets.Streamer(b);
2927  Int_t nbaskets = fBaskets.GetEntries();
2928  for (j=fWriteBasket,n=0;j>0 && n<nbaskets;--j) {
2929  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2930  if (bk) {
2931  bk->SetBranch(this);
2932  //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2933  ++n;
2934  }
2935  }
2937  b >> n;
2938  for (i=0;i<n;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2940  if (v > 4) {
2941  n = b.ReadArray(fBasketBytes);
2942  } else {
2943  for (n=0;n<fMaxBaskets;n++) fBasketBytes[n] = 0;
2944  }
2945  if (v < 2) {
2947  for (n=0;n<fWriteBasket;n++) {
2948  TBasket *basket = GetBasketImpl(n, nullptr);
2949  fBasketSeek[n] = basket ? basket->GetSeekKey() : 0;
2950  }
2951  } else {
2953  b >> n;
2954  for (n=0;n<fMaxBaskets;n++) {
2955  Int_t aseek;
2956  b >> aseek;
2957  fBasketSeek[n] = Long64_t(aseek);
2958  }
2959  }
2960  if (v > 2) {
2961  fFileName.Streamer(b);
2962  }
2963  fDirectory = 0;
2964  if (v < 4) SetAutoDelete(kTRUE);
2966  gROOT->SetReadingObject(kFALSE);
2967  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2968  //====end of old versions
2969  if (IsA() == TBranch::Class()) {
2970  if (fNleaves == 0) {
2972  } else if (fNleaves == 1) {
2974  } else if (fNleaves == 2) {
2976  } else {
2978  }
2979  }
2980  } else {
2981  Int_t maxBaskets = fMaxBaskets;
2983  Int_t lastBasket = fMaxBaskets;
2984  if (fMaxBaskets < 10) fMaxBaskets = 10;
2985 
2986  TBasket **stash = new TBasket *[lastBasket];
2987  for (Int_t i = 0; i < lastBasket; ++i) {
2988  TBasket *ba = (TBasket *)fBaskets.UncheckedAt(i);
2989  if (ba && (fBasketBytes[i] || ba->GetNevBuf()==0)) {
2990  // Already on disk or empty.
2991  stash[i] = ba;
2992  fBaskets[i] = nullptr;
2993  } else {
2994  stash[i] = nullptr;
2995  }
2996  }
2997 
2998  b.WriteClassBuffer(TBranch::Class(), this);
2999 
3000  for (Int_t i = 0; i < lastBasket; ++i) {
3001  if (stash[i]) fBaskets[i] = stash[i];
3002  }
3003 
3004  delete[] stash;
3005  fMaxBaskets = maxBaskets;
3006  }
3007 }
3008 
3009 ////////////////////////////////////////////////////////////////////////////////
3010 /// Write the current basket to disk and return the number of bytes
3011 /// written to the file.
3012 
3014 {
3015  Int_t nevbuf = basket->GetNevBuf();
3016  if (fEntryOffsetLen > 10 && (4*nevbuf) < fEntryOffsetLen ) {
3017  // Make sure that the fEntryOffset array does not stay large unnecessarily.
3018  fEntryOffsetLen = nevbuf < 3 ? 10 : 4*nevbuf; // assume some fluctuations.
3019  } else if (fEntryOffsetLen && nevbuf > fEntryOffsetLen) {
3020  // Increase the array ...
3021  fEntryOffsetLen = 2*nevbuf; // assume some fluctuations.
3022  }
3023 
3024  // Note: captures `basket`, `where`, and `this` by value; modifies the TBranch and basket,
3025  // as we make a copy of the pointer. We cannot capture `basket` by reference as the pointer
3026  // itself might be modified after `WriteBasketImpl` exits.
3027  auto doUpdates = [=]() {
3028  Int_t nout = basket->WriteBuffer(); // Write buffer
3029  if (nout < 0) Error("TBranch::WriteBasketImpl", "basket's WriteBuffer failed.\n");
3030  fBasketBytes[where] = basket->GetNbytes();
3031  fBasketSeek[where] = basket->GetSeekKey();
3032  Int_t addbytes = basket->GetObjlen() + basket->GetKeylen();
3033  TBasket *reusebasket = 0;
3034  if (nout>0) {
3035  // The Basket was written so we can now safely reuse it.
3036  fBaskets[where] = 0;
3037 
3038  reusebasket = basket;
3039  reusebasket->Reset();
3040 
3041  fZipBytes += nout;
3042  fTotBytes += addbytes;
3043  fTree->AddTotBytes(addbytes);
3044  fTree->AddZipBytes(nout);
3045 #ifdef R__TRACK_BASKET_ALLOC_TIME
3046  fTree->AddAllocationTime(reusebasket->GetResetAllocationTime());
3047 #endif
3049  }
3050 
3051  if (where==fWriteBasket) {
3052  ++fWriteBasket;
3053  if (fWriteBasket >= fMaxBaskets) {
3055  }
3056  if (reusebasket && reusebasket == fCurrentBasket) {
3057  // The 'current' basket has Reset, so if we need it we will need
3058  // to reload it.
3059  fCurrentBasket = 0;
3060  fFirstBasketEntry = -1;
3061  fNextBasketEntry = -1;
3062  }
3063  fBaskets.AddAtAndExpand(reusebasket,fWriteBasket);
3065  } else {
3066  --fNBaskets;
3067  fBaskets[where] = 0;
3068  basket->DropBuffers();
3069  if (basket == fCurrentBasket) {
3070  fCurrentBasket = 0;
3071  fFirstBasketEntry = -1;
3072  fNextBasketEntry = -1;
3073  }
3074  delete basket;
3075  }
3076  return nout;
3077  };
3078  if (imtHelper) {
3079  imtHelper->Run(doUpdates);
3080  return 0;
3081  } else {
3082  return doUpdates();
3083  }
3084 }
3085 
3086 ////////////////////////////////////////////////////////////////////////////////
3087 ///set the first entry number (case of TBranchSTL)
3088 
3090 {
3091  fFirstEntry = entry;
3092  fEntries = 0;
3093  fEntryNumber = entry;
3094  if( fBasketEntry )
3095  fBasketEntry[0] = entry;
3096  for( Int_t i = 0; i < fBranches.GetEntriesFast(); ++i )
3097  ((TBranch*)fBranches[i])->SetFirstEntry( entry );
3098 }
3099 
3100 ////////////////////////////////////////////////////////////////////////////////
3101 /// If the branch address is not set, we set all addresses starting with
3102 /// the top level parent branch.
3103 
3105 {
3106  // Nothing to do for regular branch, the TLeaf already did it.
3107 }
3108 
3109 ////////////////////////////////////////////////////////////////////////////////
3110 /// Refresh the value of fDirectory (i.e. where this branch writes/reads its buffers)
3111 /// with the current value of fTree->GetCurrentFile unless this branch has been
3112 /// redirected to a different file. Also update the sub-branches.
3113 
3115 {
3117  if (fFileName.Length() == 0) {
3118  fDirectory = file;
3119 
3120  // Apply to all existing baskets.
3121  TIter nextb(GetListOfBaskets());
3122  TBasket *basket;
3123  while ((basket = (TBasket*)nextb())) {
3124  basket->SetParent(file);
3125  }
3126  }
3127 
3128  // Apply to sub-branches as well.
3129  TIter next(GetListOfBranches());
3130  TBranch *branch;
3131  while ((branch = (TBranch*)next())) {
3132  branch->UpdateFile();
3133  }
3134 }
void Init(const char *name, const char *leaflist, Int_t compress)
Definition: TBranch.cxx:295
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition: TSystem.cxx:926
Int_t * GetEntryOffset()
Definition: TBasket.h:122
Int_t ReadBasketBytes(Long64_t pos, TFile *file)
Read basket buffers in memory and cleanup.
Definition: TBasket.cxx:694
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Int_t fNBaskets
! Number of baskets in memory
Definition: TBranch.h:122
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
void SetBufferOffset(Int_t offset=0)
Definition: TBuffer.h:92
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:49
Long64_t Previous()
Move on to the previous cluster and return the starting entry of this previous cluster.
Definition: TTree.cxx:677
Bool_t IsReading() const
Definition: TBuffer.h:85
TTreeCache * GetReadCache(TFile *file) const
Find and return the TTreeCache registered with the file and which may contain branches for us...
Definition: TTree.cxx:6156
virtual void AddTotBytes(Int_t tot)
Definition: TTree.h:317
A TLeaf for an 8 bit Integer data type.
Definition: TLeafB.h:26
virtual Bool_t IsAbsoluteFileName(const char *dir)
Return true if dir is an absolute pathname.
Definition: TSystem.cxx:943
Compression level reserved when we are not sure what to use (1 is for the fastest compression) ...
Definition: Compression.h:68
virtual void SetBufferAddress(TBuffer *entryBuffer)
Set address of this branch directly from a TBuffer to avoid streaming.
Definition: TBranch.cxx:2559
A TLeaf for a 24 bit truncated floating point data type.
Definition: TLeafF16.h:26
An array of TObjects.
Definition: TObjArray.h:37
virtual void SetReadMode()
Set read mode of basket.
Definition: TBasket.cxx:848
virtual Int_t FillImpl(ROOT::Internal::TBranchIMTHelper *)
Loop on all leaves of this branch to fill Basket buffer.
Definition: TBranch.cxx:833
Long64_t GetNextEntry()
Definition: TTree.h:293
A TLeaf for a 64 bit floating point data type.
Definition: TLeafD.h:26
void Update(Int_t newlast)
Definition: TBasket.h:150
virtual void UpdateFile()
Refresh the value of fDirectory (i.e.
Definition: TBranch.cxx:3114
virtual void SetUsed(TBranch *b, size_t basketNumber)=0
The concrete implementation of TBuffer for writing/reading to/from a ROOT file or socket...
Definition: TBufferFile.h:46
virtual void SetAddress(void *add)
Set address of this branch.
Definition: TBranch.cxx:2496
Bool_t SupportsBulkRead() const
Returns true if this branch supports bulk IO, false otherwise.
Definition: TBranch.cxx:1401
Int_t ReadBasketBuffers(Long64_t pos, Int_t len, TFile *file)
Read basket buffers in memory and cleanup.
Definition: TBasket.cxx:460
Int_t fOffset
Offset of this branch.
Definition: TBranch.h:120
long long Long64_t
Definition: RtypesCore.h:69
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3925
virtual void IncrementTotalBuffers(Int_t nbytes)
Definition: TTree.h:531
short Version_t
Definition: RtypesCore.h:61
char * GetCurrent() const
Definition: TBuffer.h:96
TObjArray * GetListOfBaskets()
Definition: TBranch.h:238
virtual ~TBranch()
Destructor.
Definition: TBranch.cxx:440
virtual TString GetDirName(const char *pathname)
Return the directory name in pathname.
Definition: TSystem.cxx:1023
Int_t GetNevBufSize() const
Definition: TBasket.h:128
virtual void AddZipBytes(Int_t zip)
Definition: TTree.h:318
Long64_t fEntries
Number of entries.
Definition: TBranch.h:130
float Float_t
Definition: RtypesCore.h:53
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
A cache when reading files over the network.
const char Option_t
Definition: RtypesCore.h:62
virtual Bool_t IsLearning() const
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:355
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:2052
virtual Bool_t GetClusterPrefetch() const
Definition: TTree.h:442
Long64_t fZipBytes
Total number of bytes in all leaves after compression.
Definition: TBranch.h:133
const Ssiz_t kNPOS
Definition: RtypesCore.h:111
This class represents a WWW compatible URL.
Definition: TUrl.h:35
const UInt_t kByteCountMask
Definition: TBufferFile.cxx:51
ReadLeaves_t fReadLeaves
! Pointer to the ReadLeaves implementation to use.
Definition: TBranch.h:157
If set, the branch&#39;s buffers will grow until an event cluster boundary is hit, guaranteeing a basket ...
Definition: TTree.h:242
Int_t FillEntryBuffer(TBasket *basket, TBuffer *buf, Int_t &lnew)
Copy the data from fEntryBuffer into the current basket.
Definition: TBranch.cxx:912
TBranch * GetSubBranch(const TBranch *br) const
Find the parent branch of child.
Definition: TBranch.cxx:1977
virtual void DropBaskets(Option_t *option="")
Loop on all branch baskets.
Definition: TBranch.cxx:734
virtual void SetStatus(Bool_t status=1)
Set branch status to Process or DoNotProcess.
Definition: TBranch.cxx:2761
virtual const char * GetTypeName() const
Definition: TLeaf.h:129
void ReadLeaves2Impl(TBuffer &b)
Read two leaves without the overhead of a loop.
Definition: TBranch.cxx:2301
virtual void SetFirstEntry(Long64_t entry)
set the first entry number (case of TBranchSTL)
Definition: TBranch.cxx:3089
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:48
virtual TLeaf * FindLeaf(const char *name)
Find the leaf corresponding to the name &#39;searchname&#39;.
Definition: TBranch.cxx:1058
virtual Long64_t GetAutoFlush() const
Definition: TTree.h:432
#define R__unlikely(expr)
Definition: RConfig.hxx:604
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
Bool_t IsFolder() const
Return kTRUE if more than one leaf or browsables, kFALSE otherwise.
Definition: TBranch.cxx:2085
TBasket * fCurrentBasket
! Pointer to the current basket.
Definition: TBranch.h:129
#define R__ASSERT(e)
Definition: TError.h:96
virtual TBasket * CreateBasket(TBranch *)
Create a basket for this tree and given branch.
Definition: TTree.cxx:3639
#define N
#define gROOT
Definition: TROOT.h:415
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
Int_t GetBulkEntries(Long64_t, TBuffer &)
Read as many events as possible into the given buffer, using zero-copy mechanisms.
Definition: TBranch.cxx:1426
virtual Int_t GetOffset() const
Definition: TLeaf.h:127
virtual Int_t LearnBranch(TBranch *, Bool_t=kFALSE)
TObjArray fBaskets
-> List of baskets of this branch
Definition: TBranch.h:136
Basic string class.
Definition: TString.h:131
#define R__likely(expr)
Definition: RConfig.hxx:605
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
int Int_t
Definition: RtypesCore.h:41
A TLeaf for a bool data type.
Definition: TLeafO.h:26
R__EXTERN TVirtualMutex * gROOTMutex
Definition: TROOT.h:59
TBasket * GetFreshCluster()
Drops the cluster two behind the current cluster and returns a fresh basket by either reusing or crea...
Definition: TBranch.cxx:1819
Int_t fNleaves
! Number of leaves
Definition: TBranch.h:124
TIOFeatures provides the end-user with the ability to change the IO behavior of data written via a TT...
Definition: TIOFeatures.hxx:62
const UInt_t kNewClassTag
Definition: TBufferFile.cxx:49
TString GetRealFileName() const
Get real file name.
Definition: TBranch.cxx:1897
Type GetType(const std::string &Name)
Definition: Systematics.cxx:34
TObjArray fLeaves
-> List of leaves of this branch
Definition: TBranch.h:135
TString & Prepend(const char *cs)
Definition: TString.h:656
Long64_t * fBasketSeek
[fMaxBaskets] Addresses of baskets on file
Definition: TBranch.h:139
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
virtual void DeleteBaskets(Option_t *option="")
Loop on all branch baskets.
Definition: TBranch.cxx:703
A TLeaf for a variable length string.
Definition: TLeafC.h:26
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:1937
TBuffer * GetBufferRef() const
Definition: TKey.h:76
virtual void SetupAddresses()
If the branch address is not set, we set all addresses starting with the top level parent branch...
Definition: TBranch.cxx:3104
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:693
const char * GetUrl(Bool_t withDeflt=kFALSE) const
Return full URL.
Definition: TUrl.cxx:387
Long64_t fFirstBasketEntry
! First entry in the current basket.
Definition: TBranch.h:127
void SetBranch(TBranch *branch)
Definition: TBasket.h:146
Int_t * fBasketBytes
[fMaxBaskets] Length of baskets on file
Definition: TBranch.h:137
virtual TList * GetBrowsables()
Returns (and, if 0, creates) browsable objects for this branch See TVirtualBranchBrowsable::FillListO...
Definition: TBranch.cxx:1289
Long64_t fEntryNumber
Current entry number (last one filled in this branch)
Definition: TBranch.h:117
Int_t Length() const
Definition: TBuffer.h:99
#define R__LOCKGUARD_IMT(mutex)
TBasket * GetFreshBasket(TBuffer *user_buffer)
Return a fresh basket by either resusing an existing basket that needs to be drop (according to TTree...
Definition: TBranch.cxx:1765
static Int_t FillListOfBrowsables(TList &list, const TBranch *branch, const TVirtualBranchBrowsable *parent=0)
Askes all registered generators to fill their browsables into the list.
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:2034
Int_t fWriteBasket
Last basket number written.
Definition: TBranch.h:116
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:473
Helper class to iterate over cluster of baskets.
Definition: TTree.h:255
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:127
Long64_t fTotBytes
Total number of bytes in all leaves before compression.
Definition: TBranch.h:132
Fill Area Attributes class.
Definition: TAttFill.h:19
Bool_t GetResetAllocationCount() const
Definition: TBasket.h:141
void Class()
Definition: Class.C:29
void SetLast(Int_t last)
Set index of last object in array, effectively truncating the array.
Definition: TObjArray.cxx:774
virtual Bool_t IsWritable() const
Definition: TDirectory.h:173
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual void Browse(TBrowser *b)
Browser interface.
Definition: TBranch.cxx:676
virtual TLeaf * GetLeaf(const char *name) const
Return pointer to the 1st Leaf named name in thisBranch.
Definition: TBranch.cxx:1884
A TLeaf for a 24 bit truncated floating point data type.
Definition: TLeafD32.h:26
void ReadLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to read Basket buffer.
Definition: TBranch.cxx:2275
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:2741
virtual TVirtualPerfStats * GetPerfStats() const
Definition: TTree.h:491
char * Buffer() const
Definition: TBuffer.h:95
void Clear()
Clear string without changing its capacity.
Definition: TString.cxx:1176
virtual TBranch * FindBranch(const char *name)
Find the immediate sub-branch with passed name.
Definition: TBranch.cxx:1012
virtual void Reset()
Reset the basket to the starting state.
Definition: TBasket.cxx:729
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition: TTree.cxx:5311
virtual Int_t WriteBuffer()
Write buffer of this basket on the current file.
Definition: TBasket.cxx:1054
Int_t GetEntriesSerialized(Long64_t N, TBuffer &user_buf)
Definition: TBranch.h:179
Bool_t IsAutoDelete() const
Return kTRUE if an existing object in a TBranchObject must be deleted.
Definition: TBranch.cxx:2077
virtual bool ReadBasketFast(TBuffer &, Long64_t)
Definition: TLeaf.h:144
virtual Int_t GetLenType() const
Definition: TLeaf.h:123
TObjArray * GetListOfBranches()
Definition: TBranch.h:239
void SetCompressionSettings(Int_t settings=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault)
Set compression settings.
Definition: TBranch.cxx:2619
const char * GetAnchor() const
Definition: TUrl.h:72
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition: TString.cxx:499
virtual void SetBranch(TBranch *branch)
Definition: TLeaf.h:151
void tobuf(char *&buf, Bool_t x)
Definition: Bytes.h:57
Int_t fMaxBaskets
Maximum number of Baskets so far.
Definition: TBranch.h:121
static constexpr double s
virtual DeserializeType GetDeserializeType() const
Definition: TLeaf.h:108
virtual void SetAddress(void *add=0)
Definition: TLeaf.h:175
virtual void FillBasket(TBuffer &b)
Pack leaf elements in Basket output buffer.
Definition: TLeaf.cxx:157
virtual TFile * GetFile() const
Definition: TDirectory.h:157
Int_t fSplitLevel
Branch split level.
Definition: TBranch.h:123
R__ALWAYS_INLINE Bool_t IsZombie() const
Definition: TObject.h:134
TBranch()
Default constructor. Used for I/O by default.
Definition: TBranch.cxx:84
void Expand(Int_t newsize, Bool_t copy=kTRUE)
Expand (or shrink) the I/O buffer to newsize bytes.
Definition: TBuffer.cxx:223
A doubly linked list.
Definition: TList.h:44
virtual void ResetAfterMerge(TFileMergeInfo *)
Reset a Branch.
Definition: TBranch.cxx:2412
Int_t Fill()
Definition: TBranch.h:199
const char * GetIconName() const
Return icon name depending on type of branch.
Definition: TBranch.cxx:1309
void PrintCacheInfo() const
Print the information we have about which basket is currently cached and whether they have been &#39;used...
Definition: TBranch.cxx:2259
Int_t GetObjlen() const
Definition: TKey.h:84
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
Int_t GetLast() const
Return index of last object in array.
Definition: TObjArray.cxx:576
virtual void AddLastBasket(Long64_t startEntry)
Add the start entry of the write basket (not yet created)
Definition: TBranch.cxx:599
virtual Int_t GetLen() const
Return the number of effective elements of this leaf, for the current entry.
Definition: TLeaf.cxx:368
Int_t * GetDisplacement() const
Definition: TBasket.h:121
Int_t WriteBasket(TBasket *basket, Int_t where)
Definition: TBranch.h:171
TBuffer * fEntryBuffer
! Buffer used to directly pass the content without streaming
Definition: TBranch.h:146
virtual char * ReadString(char *s, Int_t max)=0
static Int_t fgCount
! branch counter
Definition: TBranch.h:112
void FillLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to fill Basket buffer.
Definition: TBranch.cxx:2310
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
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:234
Int_t fBasketSize
Initial Size of Basket Buffer.
Definition: TBranch.h:114
CacheInfo_t fCacheInfo
! Hold info about which basket are in the cache and if they have been retrieved from the cache...
Definition: TBranch.h:154
R__EXTERN TSystem * gSystem
Definition: TSystem.h:557
virtual bool ReadBasketSerialized(TBuffer &, Long64_t)
Definition: TLeaf.h:145
virtual void SetEntries(Long64_t entries)
Set the number of entries in this branch.
Definition: TBranch.cxx:2652
TList * fBrowsables
! List of TVirtualBranchBrowsables used for Browse()
Definition: TBranch.h:148
virtual void AdjustSize(Int_t newsize)
Increase the size of the current fBuffer up to newsize.
Definition: TBasket.cxx:123
void SetCompressionLevel(Int_t level=ROOT::RCompressionSetting::ELevel::kUseMin)
Set compression level.
Definition: TBranch.cxx:2597
virtual TObject * RemoveAt(Int_t idx)
Remove object at index idx.
Definition: TObjArray.cxx:693
Int_t GetLast() const
Definition: TBasket.h:129
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:442
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2289
FillLeaves_t fFillLeaves
! Pointer to the FillLeaves implementation to use.
Definition: TBranch.h:159
unsigned int UInt_t
Definition: RtypesCore.h:42
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:887
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition: TTree.cxx:5323
Ssiz_t Length() const
Definition: TString.h:405
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:111
Long64_t fNextBasketEntry
! Next entry that will requires us to go to the next basket
Definition: TBranch.h:128
Manages buffers for branches of a Tree.
Definition: TBasket.h:34
void SetReadMode()
Set buffer in read mode.
Definition: TBuffer.cxx:302
void SetCompressionAlgorithm(Int_t algorithm=ROOT::RCompressionSetting::EAlgorithm::kUseGlobal)
Set compression algorithm.
Definition: TBranch.cxx:2577
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
TBasket * fExtraBasket
! Allocated basket not currently holding any data.
Definition: TBranch.h:118
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition: TClass.h:75
A TLeaf for a 32 bit floating point data type.
Definition: TLeafF.h:26
void AddAllocationCount(UInt_t count)
Definition: TTree.h:323
TBuffer * GetTransientBuffer(Int_t size)
Returns the transient buffer currently used by this TBranch for reading/writing baskets.
Definition: TBranch.cxx:511
TString fName
Definition: TNamed.h:32
virtual void AddAt(TObject *obj, Int_t idx)
Add object at position ids.
Definition: TObjArray.cxx:253
virtual void KeepCircular(Long64_t maxEntries)
keep a maximum of fMaxEntries in memory
Definition: TBranch.cxx:2097
virtual Long64_t GetBasketSeek(Int_t basket) const
Return address of basket in the file.
Definition: TBranch.cxx:1279
virtual void SetParent(const TObject *parent)
Set parent in key buffer.
Definition: TKey.cxx:1292
virtual void SetOffset(Int_t offset=0)
Definition: TLeaf.h:154
const Bool_t kFALSE
Definition: RtypesCore.h:88
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:2543
void SetUsed(Int_t basketNumber)
TString & Remove(Ssiz_t pos)
Definition: TString.h:668
void SetAnchor(const char *anchor)
Definition: TUrl.h:88
Long64_t Next()
Move on to the next cluster and return the starting entry of this next cluster.
Definition: TTree.cxx:633
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:90
Bool_t fSkipZip
! After being read, the buffer will not be unzipped.
Definition: TBranch.h:151
#define ClassImp(name)
Definition: Rtypes.h:365
Int_t fReadBasket
! Current basket number when reading
Definition: TBranch.h:125
virtual void ResetAddress()
Reset the address of the branch.
Definition: TBranch.cxx:2465
void AdoptBuffer(TBuffer *user_buffer)
Adopt a buffer from an external entity.
Definition: TBasket.cxx:717
virtual Long64_t GetSeekKey() const
Definition: TKey.h:86
void Printf(const char *fmt,...)
Long64_t GetStartEntry()
Definition: TTree.h:288
Describe directory structure in memory.
Definition: TDirectory.h:34
TDirectory * GetDirectory() const
Definition: TTree.h:447
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:3013
void DisownBuffer()
Disown all references to the internal buffer - some other object likely now owns it.
Definition: TBasket.cxx:709
Int_t GetNbytes() const
Definition: TKey.h:83
Int_t LowerBound() const
Definition: TObjArray.h:97
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:417
int nentries
Definition: THbookFile.cxx:89
virtual void Expand(Int_t newSize)
Expand or shrink the array to newSize elements.
Definition: TObjArray.cxx:386
TIOFeatures GetIOFeatures() const
Returns the IO settings currently in use for this branch.
Definition: TBranch.cxx:2069
EDataType
Definition: TDataType.h:28
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
#define R__LOCKGUARD(mutex)
TTree * fTree
! Pointer to Tree header
Definition: TBranch.h:140
void Browse(TBrowser *b)
Browse this collection (called by TBrowser).
TObjArray * GetListOfLeaves()
Definition: TBranch.h:240
Int_t BufferSize() const
Definition: TBuffer.h:97
TDirectory * fDirectory
! Pointer to directory where this branch buffers are stored
Definition: TBranch.h:144
virtual void ReadBasket(TBuffer &)
Definition: TLeaf.h:142
Int_t BackFill()
Loop on all leaves of this branch to back fill Basket buffer.
Definition: TBranch.cxx:657
virtual void SetUnsigned()
Definition: TLeaf.h:156
TBasket * GetBasket(Int_t basket)
Definition: TBranch.h:207
Int_t GetKeylen() const
Definition: TKey.h:81
Long64_t GetTotalSize(Option_t *option="") const
Return total number of bytes in the branch (including current buffer)
Definition: TBranch.cxx:2015
virtual void ReadBasketExport(TBuffer &, TClonesArray *, Int_t)
Definition: TLeaf.h:143
static void * ReAlloc(void *vp, size_t size)
Reallocate (i.e.
Definition: TStorage.cxx:183
Small helper to keep current directory context.
Definition: TDirectory.h:41
char Char_t
Definition: RtypesCore.h:29
virtual void SetEntryOffsetLen(Int_t len, Bool_t updateSubBranches=kFALSE)
Update the default value for the branch&#39;s fEntryOffsetLen if and only if it was already non zero (and...
Definition: TBranch.cxx:2635
virtual void SetFile(TFile *file=0)
Set file where this branch writes/reads its buffers.
Definition: TBranch.cxx:2677
A TLeaf for a 16 bit Integer data type.
Definition: TLeafS.h:26
virtual void Refresh(TBranch *b)
Refresh this branch using new information in b This function is called by TTree::Refresh.
Definition: TBranch.cxx:2322
Long64_t GetEntries() const
Definition: TBranch.h:244
Int_t GetNevBuf() const
Definition: TBasket.h:127
An array of clone (identical) objects.
Definition: TClonesArray.h:32
void ReadLeaves0Impl(TBuffer &b)
Read zero leaves without the overhead of a loop.
Definition: TBranch.cxx:2286
Long64_t * fBasketEntry
[fMaxBaskets] Table of first entry in each basket
Definition: TBranch.h:138
virtual Int_t DropBuffers()
Drop buffers of this basket if it is not the current basket.
Definition: TBasket.cxx:169
virtual void SetSkipZip(Bool_t=kTRUE)
auto * l
Definition: textangle.C:4
virtual void PrepareBasket(Long64_t)
Definition: TBasket.h:131
Definition: file.py:1
virtual void RemoveAll(TCollection *col)
Remove all objects in collection col from this collection.
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
virtual void Reset(Option_t *option="")
Reset a Branch.
Definition: TBranch.cxx:2371
void MakeZombie()
Definition: TObject.h:49
TIOFeatures fIOFeatures
IO features for newly-created baskets.
Definition: TBranch.h:119
Long64_t fReadEntry
! Current entry number when reading
Definition: TBranch.h:126
void Print(const char *owner, Long64_t *entries) const
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
virtual void AddBasket(TBasket &b, Bool_t ondisk, Long64_t startEntry)
Add the basket to this branch.
Definition: TBranch.cxx:533
TBranch * fMother
! Pointer to top-level parent branch in the tree.
Definition: TBranch.h:141
virtual Int_t LoadBaskets()
Baskets associated to this branch are forced to be in memory.
Definition: TBranch.cxx:2123
#define snprintf
Definition: civetweb.c:1540
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:1946
virtual void WriteBuf(const void *buf, Int_t max)=0
A TLeaf for a 64 bit Integer data type.
Definition: TLeafL.h:27
TTree * GetTree() const
Definition: TBranch.h:245
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
#define gPad
Definition: TVirtualPad.h:287
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
Int_t GetBufferSize() const
Definition: TBasket.h:120
virtual Int_t GetMapCount() const =0
Definition: tree.py:1
void Add(TObject *obj)
Definition: TObjArray.h:74
A TTree represents a columnar dataset.
Definition: TTree.h:72
Int_t fEntryOffsetLen
Initial Length of fEntryOffset table in the basket buffers.
Definition: TBranch.h:115
virtual void ResetMap()=0
void ResetBit(UInt_t f)
Definition: TObject.h:171
Undefined compression algorithm (must be kept the last of the list in case a new algorithm is added)...
Definition: Compression.h:100
virtual void ReadBasket(TBuffer &b)
Loop on all leaves of this branch to read Basket buffer.
Definition: TBranch.cxx:2267
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition: TSystem.cxx:1265
Definition: first.py:1
Int_t FlushBaskets()
Flush to disk all the baskets of this branch and any of subbranches.
Definition: TBranch.cxx:1113
TObjArray fBranches
-> List of Branches of this branch
Definition: TBranch.h:134
void SetNevBufSize(Int_t n)
Definition: TBasket.h:147
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:915
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 SetBufferDisplacement()=0
TBranch * GetBranch() const
Definition: TLeaf.h:107
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
Definition: TCollection.h:182
A TTree is a list of TBranches.
Definition: TBranch.h:90
TBuffer * fTransientBuffer
! Pointer to the current transient buffer.
Definition: TBranch.h:147
Int_t fCompress
Compression level and algorithm.
Definition: TBranch.h:113
virtual void Print(Option_t *option="") const
Print TBranch parameters.
Definition: TBranch.cxx:2155
TBasket * GetBasketImpl(Int_t basket, TBuffer *user_buffer)
Return pointer to basket basketnumber in this Branch.
Definition: TBranch.cxx:1203
TString fFileName
Name of file where buffers are stored ("" if in same file as Tree header)
Definition: TBranch.h:145
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual Int_t GetBufferDisplacement() const =0
TBranch * fParent
! Pointer to parent branch.
Definition: TBranch.h:142
const Int_t n
Definition: legend1.C:16
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Definition: TMathBase.h:278
Long64_t fFirstEntry
Number of the first entry in this branch.
Definition: TBranch.h:131
virtual void SetWriteMode()
Set write mode of basket.
Definition: TBasket.cxx:857
TBranch * GetMother() const
Get our top-level parent branch in the tree.
Definition: TBranch.cxx:1956
char name[80]
Definition: TGX11.cxx:109
void ExpandBasketArrays()
Increase BasketEntry buffer of a minimum of 10 locations and a maximum of 50 per cent of current size...
Definition: TBranch.cxx:802
A TLeaf for an Integer data type.
Definition: TLeafI.h:27
virtual void SetObject(void *objadd)
Set object this branch is pointing to.
Definition: TBranch.cxx:2750
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:873
static Int_t * ReAllocInt(Int_t *vp, size_t size, size_t oldsize)
Reallocate (i.e.
Definition: TStorage.cxx:295
void SetWriteMode()
Set buffer in write mode.
Definition: TBuffer.cxx:315
Int_t GetCompressionSettings() const
Definition: TFile.h:394
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual TObjArray * GetListOfLeaves()
Definition: TTree.h:474
void ReadLeaves1Impl(TBuffer &b)
Read one leaf without the overhead of a loop.
Definition: TBranch.cxx:2293
static void ResetCount()
Static function resetting fgCount.
Definition: TBranch.cxx:2488
virtual Long64_t GetMaxVirtualSize() const
Definition: TTree.h:485
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
const char * Data() const
Definition: TString.h:364
virtual Int_t GetNdata() const
Definition: TLeaf.h:126
virtual void SetAutoDelete(Bool_t autodel=kTRUE)
Set the automatic delete bit.
Definition: TBranch.cxx:2530
char * fAddress
! Address of 1st leaf (variable or object)
Definition: TBranch.h:143