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