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