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