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