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