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