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