Logo ROOT   6.12/07
Reference Guide
TBranch.cxx
Go to the documentation of this file.
1 // @(#)root/tree:$Id$
2 // Author: Rene Brun 12/01/96
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TBranch.h"
13 
14 #include "Compression.h"
15 #include "TBasket.h"
16 #include "TBranchBrowsable.h"
17 #include "TBrowser.h"
18 #include "TBuffer.h"
19 #include "TClass.h"
20 #include "TBufferFile.h"
21 #include "TClonesArray.h"
22 #include "TFile.h"
23 #include "TLeaf.h"
24 #include "TLeafB.h"
25 #include "TLeafC.h"
26 #include "TLeafD.h"
27 #include "TLeafF.h"
28 #include "TLeafI.h"
29 #include "TLeafL.h"
30 #include "TLeafO.h"
31 #include "TLeafObject.h"
32 #include "TLeafS.h"
33 #include "TMessage.h"
34 #include "TROOT.h"
35 #include "TSystem.h"
36 #include "TMath.h"
37 #include "TTree.h"
38 #include "TTreeCache.h"
39 #include "TTreeCacheUnzip.h"
40 #include "TVirtualMutex.h"
41 #include "TVirtualPad.h"
42 
43 #include "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()
207 , fBasketBytes(0)
208 , fBasketEntry(0)
209 , fBasketSeek(0)
210 , fTree(tree)
211 , fMother(0)
212 , fParent(0)
213 , fAddress((char *)address)
215 , fFileName("")
216 , fEntryBuffer(0)
217 , fTransientBuffer(0)
218 , fBrowsables(0)
219 , fSkipZip(kFALSE)
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()
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)
267 , fFileName("")
268 , fEntryBuffer(0)
269 , fTransientBuffer(0)
270 , fBrowsables(0)
271 , fSkipZip(kFALSE)
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 /// Loop on all leaves of this branch to back fill Basket buffer.
600 ///
601 /// Use this routine instead of TBranch::Fill when filling a branch individually
602 /// to catch up with the number of entries already in the TTree.
603 ///
604 /// First it calls TBranch::Fill and then if the number of entries of the branch
605 /// reach one of TTree cluster's boundary, the basket is flushed.
606 ///
607 /// The function returns the number of bytes committed to the memory basket.
608 /// If a write error occurs, the number of bytes returned is -1.
609 /// If no data are written, because e.g. the branch is disabled,
610 /// the number of bytes returned is 0.
611 ///
612 /// To insure that the baskets of each cluster are located close by in the
613 /// file, when back-filling multiple branches make sure to call BackFill
614 /// for the same entry for all the branches consecutively
615 /// ~~~ {.cpp}
616 /// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
617 /// for( auto branch : branchCollection) {
618 /// ... Make change to the data associated with the branch ...
619 /// branch->BackFill();
620 /// }
621 /// }
622 /// // Since we loop over all the branches for each new entry
623 /// // all the baskets for a cluster are consecutive in the file.
624 /// ~~~
625 /// rather than doing all the entries of one branch at a time.
626 /// ~~~ {.cpp}
627 /// // Do NOT do things in the following order, it will lead to
628 /// // poorly clustered files.
629 /// for(auto branch : branchCollection) {
630 /// for( auto e = 0; e < tree->GetEntries(); ++e ) { // loop over entries.
631 /// ... Make change to the data associated with the branch ...
632 /// branch->BackFill();
633 /// }
634 /// }
635 /// // Since we loop over all the entries for one branch
636 /// // all the baskets for that branch are consecutive.
637 /// ~~~
638 
640 
641  // Get the end of the next cluster.
642  auto cluster = GetTree()->GetClusterIterator( GetEntries() );
643  cluster.Next();
644  auto endCluster = cluster.GetNextEntry();
645 
646  auto result = FillImpl(nullptr);
647 
648  if ( result && GetEntries() >= endCluster ) {
649  FlushBaskets();
650  }
651 
652  return result;
653 }
654 
655 ////////////////////////////////////////////////////////////////////////////////
656 /// Browser interface.
657 
659 {
660  if (fNleaves > 1) {
661  fLeaves.Browse(b);
662  } else {
663  // Get the name and strip any extra brackets
664  // in order to get the full arrays.
665  TString name = GetName();
666  Int_t pos = name.First('[');
667  if (pos!=kNPOS) name.Remove(pos);
668 
669  GetTree()->Draw(name, "", b ? b->GetDrawOption() : "");
670  if (gPad) gPad->Update();
671  }
672 }
673 
674  ///////////////////////////////////////////////////////////////////////////////
675  /// Loop on all branch baskets. If the file where branch buffers reside is
676  /// writable, free the disk space associated to the baskets of the branch,
677  /// then call Reset(). If the option contains "all", delete also the baskets
678  /// for the subbranches.
679  /// The branch is reset.
680  ///
681  /// NOTE that this function must be used with extreme care. Deleting branch baskets
682  /// fragments the file and may introduce inefficiencies when adding new entries
683  /// in the Tree or later on when reading the Tree.
684 
686 {
687  TString opt = option;
688  opt.ToLower();
689  TFile *file = GetFile(0);
690 
691  if(fDirectory && (fDirectory != gROOT) && fDirectory->IsWritable()) {
692  for(Int_t i=0; i<fWriteBasket; i++) {
693  if (fBasketSeek[i]) file->MakeFree(fBasketSeek[i],fBasketSeek[i]+fBasketBytes[i]-1);
694  }
695  }
696 
697  // process subbranches
698  if (opt.Contains("all")) {
700  Int_t nb = lb->GetEntriesFast();
701  for (Int_t j = 0; j < nb; j++) {
702  TBranch* branch = (TBranch*) lb->UncheckedAt(j);
703  if (branch) branch->DeleteBaskets("all");
704  }
705  }
706  DropBaskets("all");
707  Reset();
708 }
709 
710 ////////////////////////////////////////////////////////////////////////////////
711 /// Loop on all branch baskets. Drop all baskets from memory except readbasket.
712 /// If the option contains "all", drop all baskets including
713 /// read- and write-baskets (unless they are not stored individually on disk).
714 /// The option "all" also lead to DropBaskets being called on the sub-branches.
715 
717 {
718  Bool_t all = kFALSE;
719  if (options && options[0]) {
720  TString opt = options;
721  opt.ToLower();
722  if (opt.Contains("all")) all = kTRUE;
723  }
724 
725  TBasket *basket;
726  Int_t nbaskets = fBaskets.GetEntriesFast();
727 
728  if ( (fNBaskets>1) || all ) {
729  //slow case
730  for (Int_t i=0;i<nbaskets;i++) {
731  basket = (TBasket*)fBaskets.UncheckedAt(i);
732  if (!basket) continue;
733  if ((i == fReadBasket || i == fWriteBasket) && !all) continue;
734  // if the basket is not yet on file but already has event in it
735  // we must continue to avoid dropping the basket (and thus losing data)
736  if (fBasketBytes[i]==0 && basket->GetNevBuf() > 0) continue;
737  basket->DropBuffers();
738  --fNBaskets;
739  fBaskets.RemoveAt(i);
740  if (basket == fCurrentBasket) {
741  fCurrentBasket = 0;
742  fFirstBasketEntry = -1;
743  fNextBasketEntry = -1;
744  }
745  delete basket;
746  }
747 
748  // process subbranches
749  if (all) {
751  Int_t nb = lb->GetEntriesFast();
752  for (Int_t j = 0; j < nb; j++) {
753  TBranch* branch = (TBranch*) lb->UncheckedAt(j);
754  if (!branch) continue;
755  branch->DropBaskets("all");
756  }
757  }
758  } else {
759  //fast case
760  if (nbaskets > 0) {
761  Int_t i = fBaskets.GetLast();
762  basket = (TBasket*)fBaskets.UncheckedAt(i);
763  if (basket && fBasketBytes[i]!=0) {
764  basket->DropBuffers();
765  if (basket == fCurrentBasket) {
766  fCurrentBasket = 0;
767  fFirstBasketEntry = -1;
768  fNextBasketEntry = -1;
769  }
770  delete basket;
771  fBaskets.AddAt(0,i);
772  fBaskets.SetLast(-1);
773  fNBaskets = 0;
774  }
775  }
776  }
777 
778 }
779 
780 ////////////////////////////////////////////////////////////////////////////////
781 /// Increase BasketEntry buffer of a minimum of 10 locations
782 /// and a maximum of 50 per cent of current size.
783 
785 {
786  Int_t newsize = TMath::Max(10,Int_t(1.5*fMaxBaskets));
789  newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
791  newsize*sizeof(Long64_t),fMaxBaskets*sizeof(Long64_t));
792 
793  fMaxBaskets = newsize;
794 
795  fBaskets.Expand(newsize);
796 
797  for (Int_t i=fWriteBasket;i<fMaxBaskets;i++) {
798  fBasketBytes[i] = 0;
799  fBasketEntry[i] = 0;
800  fBasketSeek[i] = 0;
801  }
802 }
803 
804 ////////////////////////////////////////////////////////////////////////////////
805 /// Loop on all leaves of this branch to fill Basket buffer.
806 ///
807 /// If TBranchIMTHelper is non-null and it is time to WriteBasket, then we will
808 /// use TBB to compress in parallel.
809 ///
810 /// The function returns the number of bytes committed to the memory basket.
811 /// If a write error occurs, the number of bytes returned is -1.
812 /// If no data are written, because e.g. the branch is disabled,
813 /// the number of bytes returned is 0.
814 
816 {
817  if (TestBit(kDoNotProcess)) {
818  return 0;
819  }
820 
821  TBasket* basket = GetBasket(fWriteBasket);
822  if (!basket) {
823  basket = fTree->CreateBasket(this); // create a new basket
824  if (!basket) return 0;
825  ++fNBaskets;
827  }
828  TBuffer* buf = basket->GetBufferRef();
829 
830  // Fill basket buffer.
831 
832  Int_t nsize = 0;
833 
834  if (buf->IsReading()) {
835  basket->SetWriteMode();
836  }
837 
838  buf->ResetMap();
839 
840  Int_t lnew = 0;
841  Int_t nbytes = 0;
842 
843  if (fEntryBuffer) {
844  nbytes = FillEntryBuffer(basket,buf,lnew);
845  } else {
846  Int_t lold = buf->Length();
847  basket->Update(lold);
848  ++fEntries;
849  ++fEntryNumber;
850  (this->*fFillLeaves)(*buf);
851  if (buf->GetMapCount()) {
852  // The map is used.
854  }
855  lnew = buf->Length();
856  nbytes = lnew - lold;
857  }
858 
859  if (fEntryOffsetLen) {
860  Int_t nevbuf = basket->GetNevBuf();
861  // Total size in bytes of EntryOffset table.
862  nsize = nevbuf * sizeof(Int_t);
863  } else {
864  if (!basket->GetNevBufSize()) {
865  basket->SetNevBufSize(nbytes);
866  }
867  }
868 
869  // Should we create a new basket?
870  // fSkipZip force one entry per buffer (old stuff still maintained for CDF)
871  // Transfer full compressed buffer only
872 
873  if ((fSkipZip && (lnew >= TBuffer::kMinimalSize)) || (buf->TestBit(TBufferFile::kNotDecompressed)) || ((lnew + (2 * nsize) + nbytes) >= fBasketSize)) {
875  return nbytes;
876  }
877  Int_t nout = WriteBasketImpl(basket, fWriteBasket, imtHelper);
878  if (nout < 0) Error("TBranch::Fill", "Failed to write out basket.\n");
879  return (nout >= 0) ? nbytes : -1;
880  }
881  return nbytes;
882 }
883 
884 ////////////////////////////////////////////////////////////////////////////////
885 /// Copy the data from fEntryBuffer into the current basket.
886 
888 {
889  Int_t nbytes = 0;
890  Int_t objectStart = 0;
891  Int_t last = 0;
892  Int_t lold = buf->Length();
893 
894  // Handle the special case of fEntryBuffer != 0
895  if (fEntryBuffer->IsA() == TMessage::Class()) {
896  objectStart = 8;
897  }
899  // The buffer given as input has not been decompressed.
900  if (basket->GetNevBuf()) {
901  // If the basket already contains entry we need to close it
902  // out. (This is because we can only transfer full compressed
903  // buffer)
904  WriteBasket(basket,fWriteBasket);
905  // And restart from scratch
906  return Fill();
907  }
908  Int_t startpos = fEntryBuffer->Length();
910  static TBasket toread_fLast;
912  toread_fLast.Streamer(*fEntryBuffer);
914  last = toread_fLast.GetLast();
915  // last now contains the decompressed number of bytes.
916  fEntryBuffer->SetBufferOffset(startpos);
917  buf->SetBufferOffset(0);
919  basket->Update(lold);
920  } else {
921  // We are required to copy starting at the version number (so not
922  // including the class name.
923  // See if byte count is here, if not it class still be a newClass
924  const UInt_t kNewClassTag = 0xFFFFFFFF;
925  const UInt_t kByteCountMask = 0x40000000; // OR the byte count with this
926  UInt_t tag = 0;
927  UInt_t startpos = fEntryBuffer->Length();
928  fEntryBuffer->SetBufferOffset(objectStart);
929  *fEntryBuffer >> tag;
930  if (tag & kByteCountMask) {
931  *fEntryBuffer >> tag;
932  }
933  if (tag == kNewClassTag) {
934  UInt_t maxsize = 256;
935  char* s = new char[maxsize];
936  Int_t name_start = fEntryBuffer->Length();
937  fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end.
938  while (strlen(s) == (maxsize - 1)) {
939  // The classname is too large, try again with a large buffer.
940  fEntryBuffer->SetBufferOffset(name_start);
941  maxsize *= 2;
942  delete[] s;
943  s = new char[maxsize];
944  fEntryBuffer->ReadString(s, maxsize); // Reads at most maxsize - 1 characters, plus null at end
945  }
946  } else {
947  fEntryBuffer->SetBufferOffset(objectStart);
948  }
949  objectStart = fEntryBuffer->Length();
950  fEntryBuffer->SetBufferOffset(startpos);
951  basket->Update(lold, objectStart - fEntryBuffer->GetBufferDisplacement());
952  }
953  fEntries++;
954  fEntryNumber++;
955  UInt_t len = 0;
956  UInt_t startpos = fEntryBuffer->Length();
957  if (startpos > UInt_t(objectStart)) {
958  // We assume this buffer have just been directly filled
959  // the current position in the buffer indicates the end of the object!
960  len = fEntryBuffer->Length() - objectStart;
961  } else {
962  // The buffer have been acquired either via TSocket or via
963  // TBuffer::SetBuffer(newloc,newsize)
964  // Only the actual size of the memory buffer gives us an hint about where
965  // the object ends.
966  len = fEntryBuffer->BufferSize() - objectStart;
967  }
968  buf->WriteBuf(fEntryBuffer->Buffer() + objectStart, len);
970  // The original buffer came pre-compressed and thus the buffer Length
971  // does not really show the really object size
972  // lnew = nbytes = basket->GetLast();
973  nbytes = last;
974  lnew = last;
975  } else {
976  lnew = buf->Length();
977  nbytes = lnew - lold;
978  }
979 
980  return nbytes;
981 }
982 
983 ////////////////////////////////////////////////////////////////////////////////
984 /// Find the immediate sub-branch with passed name.
985 
987 {
988  // We allow the user to pass only the last dotted component of the name.
989  std::string longnm;
990  longnm.reserve(fName.Length()+strlen(name)+3);
991  longnm = fName.Data();
992  if (longnm[longnm.length()-1]==']') {
993  std::size_t dim = longnm.find_first_of("[");
994  if (dim != std::string::npos) {
995  longnm.erase(dim);
996  }
997  }
998  if (longnm[longnm.length()-1] != '.') {
999  longnm += '.';
1000  }
1001  longnm += name;
1002  UInt_t namelen = strlen(name);
1003 
1004  Int_t nbranches = fBranches.GetEntries();
1005  TBranch* branch = 0;
1006  for(Int_t i = 0; i < nbranches; ++i) {
1007  branch = (TBranch*) fBranches.UncheckedAt(i);
1008 
1009  const char *brname = branch->fName.Data();
1010  UInt_t brlen = branch->fName.Length();
1011  if (brname[brlen-1]==']') {
1012  const char *dim = strchr(brname,'[');
1013  if (dim) {
1014  brlen = dim - brname;
1015  }
1016  }
1017  if (namelen == brlen /* same effective size */
1018  && strncmp(name,brname,brlen) == 0) {
1019  return branch;
1020  }
1021  if (brlen == (size_t)longnm.length()
1022  && strncmp(longnm.c_str(),brname,brlen) == 0) {
1023  return branch;
1024  }
1025  }
1026  return 0;
1027 }
1028 
1029 ////////////////////////////////////////////////////////////////////////////////
1030 /// Find the leaf corresponding to the name 'searchname'.
1031 
1032 TLeaf* TBranch::FindLeaf(const char* searchname)
1033 {
1034  TString leafname;
1035  TString leaftitle;
1036  TString longname;
1037  TString longtitle;
1038 
1039  // We allow the user to pass only the last dotted component of the name.
1040  TIter next(GetListOfLeaves());
1041  TLeaf* leaf = 0;
1042  while ((leaf = (TLeaf*) next())) {
1043  leafname = leaf->GetName();
1044  Ssiz_t dim = leafname.First('[');
1045  if (dim >= 0) leafname.Remove(dim);
1046 
1047  if (leafname == searchname) return leaf;
1048 
1049  // The leaf element contains the branch name in its name, let's use the title.
1050  leaftitle = leaf->GetTitle();
1051  dim = leaftitle.First('[');
1052  if (dim >= 0) leaftitle.Remove(dim);
1053 
1054  if (leaftitle == searchname) return leaf;
1055 
1056  TBranch* branch = leaf->GetBranch();
1057  if (branch) {
1058  longname.Form("%s.%s",branch->GetName(),leafname.Data());
1059  dim = longname.First('[');
1060  if (dim>=0) longname.Remove(dim);
1061  if (longname == searchname) return leaf;
1062 
1063  // The leaf element contains the branch name in its name.
1064  longname.Form("%s.%s",branch->GetName(),searchname);
1065  if (longname==leafname) return leaf;
1066 
1067  longtitle.Form("%s.%s",branch->GetName(),leaftitle.Data());
1068  dim = longtitle.First('[');
1069  if (dim>=0) longtitle.Remove(dim);
1070  if (longtitle == searchname) return leaf;
1071 
1072  // The following is for the case where the branch is only
1073  // a sub-branch. Since we do not see it through
1074  // TTree::GetListOfBranches, we need to see it indirectly.
1075  // This is the less sturdy part of this search ... it may
1076  // need refining ...
1077  if (strstr(searchname, ".") && !strcmp(searchname, branch->GetName())) return leaf;
1078  }
1079  }
1080  return 0;
1081 }
1082 
1083 ////////////////////////////////////////////////////////////////////////////////
1084 /// Flush to disk all the baskets of this branch and any of subbranches.
1085 /// Return the number of bytes written or -1 in case of write error.
1086 
1088 {
1089  UInt_t nerror = 0;
1090  Int_t nbytes = 0;
1091 
1092  Int_t maxbasket = fWriteBasket + 1;
1093  // The following protection is not necessary since we should always
1094  // have fWriteBasket < fBasket.GetSize()
1095  //if (fBaskets.GetSize() < maxbasket) {
1096  // maxbasket = fBaskets.GetSize();
1097  //}
1098  for(Int_t i=0; i != maxbasket; ++i) {
1099  if (fBaskets.UncheckedAt(i)) {
1100  Int_t nwrite = FlushOneBasket(i);
1101  if (nwrite<0) {
1102  ++nerror;
1103  } else {
1104  nbytes += nwrite;
1105  }
1106  }
1107  }
1108  Int_t len = fBranches.GetEntriesFast();
1109  for (Int_t i = 0; i < len; ++i) {
1110  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1111  if (!branch) {
1112  continue;
1113  }
1114  Int_t nwrite = branch->FlushBaskets();
1115  if (nwrite<0) {
1116  ++nerror;
1117  } else {
1118  nbytes += nwrite;
1119  }
1120  }
1121  if (nerror) {
1122  return -1;
1123  } else {
1124  return nbytes;
1125  }
1126 }
1127 
1128 ////////////////////////////////////////////////////////////////////////////////
1129 /// If we have a write basket in memory and it contains some entries and
1130 /// has not yet been written to disk, we write it and delete it from memory.
1131 /// Return the number of bytes written;
1132 
1134 {
1135  Int_t nbytes = 0;
1136  if (fDirectory && fBaskets.GetEntries()) {
1137  TBasket *basket = (TBasket*)fBaskets.UncheckedAt(ibasket);
1138 
1139  if (basket) {
1140  if (basket->GetNevBuf()
1141  && fBasketSeek[ibasket]==0) {
1142  // If the basket already contains entry we need to close it out.
1143  // (This is because we can only transfer full compressed buffer)
1144 
1145  if (basket->GetBufferRef()->IsReading()) {
1146  basket->SetWriteMode();
1147  }
1148  nbytes = WriteBasket(basket,ibasket);
1149 
1150  } else {
1151  // If the basket is empty or has already been written.
1152  if ((Int_t)ibasket==fWriteBasket) {
1153  // Nothing to do.
1154  } else {
1155  basket->DropBuffers();
1156  if (basket == fCurrentBasket) {
1157  fCurrentBasket = 0;
1158  fFirstBasketEntry = -1;
1159  fNextBasketEntry = -1;
1160  }
1161  delete basket;
1162  --fNBaskets;
1163  fBaskets[ibasket] = 0;
1164  }
1165  }
1166  }
1167  }
1168  return nbytes;
1169 }
1170 
1171 ////////////////////////////////////////////////////////////////////////////////
1172 /// Return pointer to basket basketnumber in this Branch
1173 
1175 {
1176  // This counter in the sequential case collects errors coming also from
1177  // different files (suppose to have a program reading f1.root, f2.root ...)
1178  // In the mt case, it is made atomic: it safely collects errors from
1179  // different files processed simultaneously.
1180  static std::atomic<Int_t> nerrors(0);
1181 
1182  // reference to an existing basket in memory ?
1183  if (basketnumber <0 || basketnumber > fWriteBasket) return 0;
1184  TBasket *basket = (TBasket*)fBaskets.UncheckedAt(basketnumber);
1185  if (basket) return basket;
1186  if (basketnumber == fWriteBasket) return 0;
1187 
1188  // create/decode basket parameters from buffer
1189  TFile *file = GetFile(0);
1190  if (file == 0) {
1191  return 0;
1192  }
1193  // if cluster pre-fetching or retaining is on, do not re-use existing baskets
1194  // unless a new cluster is used.
1196  basket = GetFreshCluster();
1197  else
1198  basket = GetFreshBasket();
1199 
1200  // fSkipZip is old stuff still maintained for CDF
1202  if (fBasketBytes[basketnumber] == 0) {
1203  fBasketBytes[basketnumber] = basket->ReadBasketBytes(fBasketSeek[basketnumber],file);
1204  }
1205  //add branch to cache (if any)
1206  {
1207  R__LOCKGUARD_IMT2(gROOTMutex); // Lock for parallel TTree I/O
1208  TFileCacheRead *pf = file->GetCacheRead(fTree);
1209  if (pf){
1210  if (pf->IsLearning()) pf->AddBranch(this);
1211  if (fSkipZip) pf->SetSkipZip();
1212  }
1213  }
1214 
1215  //now read basket
1216  Int_t badread = basket->ReadBasketBuffers(fBasketSeek[basketnumber],fBasketBytes[basketnumber],file);
1217  if (R__unlikely(badread || basket->GetSeekKey() != fBasketSeek[basketnumber] || basket->IsZombie())) {
1218  nerrors++;
1219  if (nerrors > 10) return 0;
1220  if (nerrors == 10) {
1221  printf(" file probably overwritten: stopping reporting error messages\n");
1222  if (fBasketSeek[basketnumber] > 2000000000) {
1223  printf("===>File is more than 2 Gigabytes\n");
1224  return 0;
1225  }
1226  if (fBasketSeek[basketnumber] > 1000000000) {
1227  printf("===>Your file is may be bigger than the maximum file size allowed on your system\n");
1228  printf(" Check your AFS maximum file size limit for example\n");
1229  return 0;
1230  }
1231  }
1232  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);
1233  return 0;
1234  }
1235 
1236  ++fNBaskets;
1237  fBaskets.AddAt(basket,basketnumber);
1238  return basket;
1239 }
1240 
1241 ////////////////////////////////////////////////////////////////////////////////
1242 /// Return address of basket in the file
1243 
1245 {
1246  if (basketnumber <0 || basketnumber > fWriteBasket) return 0;
1247  return fBasketSeek[basketnumber];
1248 }
1249 
1250 ////////////////////////////////////////////////////////////////////////////////
1251 /// Returns (and, if 0, creates) browsable objects for this branch
1252 /// See TVirtualBranchBrowsable::FillListOfBrowsables.
1253 
1255  if (fBrowsables) return fBrowsables;
1256  fBrowsables=new TList();
1258  return fBrowsables;
1259 }
1260 
1261 ////////////////////////////////////////////////////////////////////////////////
1262 /// Return the name of the user class whose content is stored in this branch,
1263 /// if any. If this branch was created using the 'leaflist' technique, this
1264 /// function returns an empty string.
1265 
1266 const char * TBranch::GetClassName() const
1267 {
1268  return "";
1269 }
1270 
1271 ////////////////////////////////////////////////////////////////////////////////
1272 /// Return icon name depending on type of branch.
1273 
1274 const char* TBranch::GetIconName() const
1275 {
1276  if (IsFolder())
1277  return "TBranchElement-folder";
1278  else
1279  return "TBranchElement-leaf";
1280 }
1281 
1282 ////////////////////////////////////////////////////////////////////////////////
1283 /// Read all leaves of entry and return total number of bytes read.
1284 ///
1285 /// The input argument "entry" is the entry number in the current tree.
1286 /// In case of a TChain, the entry number in the current Tree must be found
1287 /// before calling this function. For example:
1288 ///
1289 ///~~~ {.cpp}
1290 /// TChain* chain = ...;
1291 /// Long64_t localEntry = chain->LoadTree(entry);
1292 /// branch->GetEntry(localEntry);
1293 ///~~~
1294 ///
1295 /// The function returns the number of bytes read from the input buffer.
1296 /// If entry does not exist, the function returns 0.
1297 /// If an I/O error occurs, the function returns -1.
1298 ///
1299 /// See IMPORTANT REMARKS in TTree::GetEntry.
1300 
1302 {
1303  // Remember which entry we are reading.
1304  fReadEntry = entry;
1305 
1306  Bool_t enabled = !TestBit(kDoNotProcess) || getall;
1307  TBasket *basket; // will be initialized in the if/then clauses.
1308  Long64_t first;
1309  if (R__likely(enabled && fFirstBasketEntry <= entry && entry < fNextBasketEntry)) {
1310  // We have found the basket containing this entry.
1311  // make sure basket buffers are in memory.
1312  basket = fCurrentBasket;
1313  first = fFirstBasketEntry;
1314  } else {
1315  if (!enabled) {
1316  return 0;
1317  }
1318  if ((entry < fFirstEntry) || (entry >= fEntryNumber)) {
1319  return 0;
1320  }
1321  first = fFirstBasketEntry;
1322  Long64_t last = fNextBasketEntry - 1;
1323  // Are we still in the same ReadBasket?
1324  if ((entry < first) || (entry > last)) {
1326  if (fReadBasket < 0) {
1327  fNextBasketEntry = -1;
1328  Error("In the branch %s, no basket contains the entry %d\n", GetName(), entry);
1329  return -1;
1330  }
1331  if (fReadBasket == fWriteBasket) {
1333  } else {
1335  }
1337  }
1338  // We have found the basket containing this entry.
1339  // make sure basket buffers are in memory.
1340  basket = (TBasket*) fBaskets.UncheckedAt(fReadBasket);
1341  if (!basket) {
1342  basket = GetBasket(fReadBasket);
1343  if (!basket) {
1344  fCurrentBasket = 0;
1345  fFirstBasketEntry = -1;
1346  fNextBasketEntry = -1;
1347  return -1;
1348  }
1349  if (fTree->GetClusterPrefetch()) {
1350  TTree::TClusterIterator clusterIterator = fTree->GetClusterIterator(entry);
1351  clusterIterator.Next();
1352  Int_t nextClusterEntry = clusterIterator.GetNextEntry();
1353  for (Int_t i = fReadBasket + 1; i < fMaxBaskets && fBasketEntry[i] < nextClusterEntry; i++) {
1354  GetBasket(i);
1355  }
1356  }
1357  }
1358  fCurrentBasket = basket;
1359  }
1360  basket->PrepareBasket(entry);
1361  TBuffer* buf = basket->GetBufferRef();
1362 
1363  // This test necessary to read very old Root files (NvE).
1364  if (R__unlikely(!buf)) {
1365  TFile* file = GetFile(0);
1366  if (!file) return -1;
1367  basket->ReadBasketBuffers(fBasketSeek[fReadBasket], fBasketBytes[fReadBasket], file);
1368  buf = basket->GetBufferRef();
1369  }
1370  // Set entry offset in buffer.
1371  if (!TestBit(kDoNotUseBufferMap)) {
1372  buf->ResetMap();
1373  }
1374  if (R__unlikely(!buf->IsReading())) {
1375  basket->SetReadMode();
1376  }
1377  Int_t* entryOffset = basket->GetEntryOffset();
1378  Int_t bufbegin = 0;
1379  if (entryOffset) {
1380  bufbegin = entryOffset[entry-first];
1381  buf->SetBufferOffset(bufbegin);
1382  Int_t* displacement = basket->GetDisplacement();
1383  if (R__unlikely(displacement)) {
1384  buf->SetBufferDisplacement(displacement[entry-first]);
1385  }
1386  } else {
1387  bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1388  buf->SetBufferOffset(bufbegin);
1389  }
1390 
1391  // Int_t bufbegin = buf->Length();
1392  (this->*fReadLeaves)(*buf);
1393  return buf->Length() - bufbegin;
1394 }
1395 
1396 ////////////////////////////////////////////////////////////////////////////////
1397 /// Read all leaves of an entry and export buffers to real objects in a TClonesArray list.
1398 ///
1399 /// Returns total number of bytes read.
1400 
1402 {
1403  // Remember which entry we are reading.
1404  fReadEntry = entry;
1405 
1406  if (TestBit(kDoNotProcess)) {
1407  return 0;
1408  }
1409  if ((entry < 0) || (entry >= fEntryNumber)) {
1410  return 0;
1411  }
1412  Int_t nbytes = 0;
1414  Long64_t last = fNextBasketEntry - 1;
1415  // Are we still in the same ReadBasket?
1416  if ((entry < first) || (entry > last)) {
1418  if (fReadBasket < 0) {
1419  fNextBasketEntry = -1;
1420  Error("In the branch %s, no basket contains the entry %d\n", GetName(), entry);
1421  return -1;
1422  }
1423  if (fReadBasket == fWriteBasket) {
1425  } else {
1427  }
1429  }
1430 
1431  // We have found the basket containing this entry.
1432  // Make sure basket buffers are in memory.
1433  TBasket* basket = GetBasket(fReadBasket);
1434  fCurrentBasket = basket;
1435  if (!basket) {
1436  fFirstBasketEntry = -1;
1437  fNextBasketEntry = -1;
1438  return 0;
1439  }
1440  TBuffer* buf = basket->GetBufferRef();
1441  // Set entry offset in buffer and read data from all leaves.
1442  if (!TestBit(kDoNotUseBufferMap)) {
1443  buf->ResetMap();
1444  }
1445  if (R__unlikely(!buf->IsReading())) {
1446  basket->SetReadMode();
1447  }
1448  Int_t* entryOffset = basket->GetEntryOffset();
1449  Int_t bufbegin = 0;
1450  if (entryOffset) {
1451  bufbegin = entryOffset[entry-first];
1452  buf->SetBufferOffset(bufbegin);
1453  Int_t* displacement = basket->GetDisplacement();
1454  if (R__unlikely(displacement)) {
1455  buf->SetBufferDisplacement(displacement[entry-first]);
1456  }
1457  } else {
1458  bufbegin = basket->GetKeylen() + ((entry-first) * basket->GetNevBufSize());
1459  buf->SetBufferOffset(bufbegin);
1460  }
1461  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(0);
1462  leaf->ReadBasketExport(*buf, li, nentries);
1463  nbytes = buf->Length() - bufbegin;
1464  return nbytes;
1465 }
1466 
1467 ////////////////////////////////////////////////////////////////////////////////
1468 /// Fill expectedClass and expectedType with information on the data type of the
1469 /// object/values contained in this branch (and thus the type of pointers
1470 /// expected to be passed to Set[Branch]Address
1471 /// return 0 in case of success and > 0 in case of failure.
1472 
1473 Int_t TBranch::GetExpectedType(TClass *&expectedClass,EDataType &expectedType)
1474 {
1475  expectedClass = 0;
1476  expectedType = kOther_t;
1477  TLeaf* l = (TLeaf*) GetListOfLeaves()->At(0);
1478  if (l) {
1479  expectedType = (EDataType) gROOT->GetType(l->GetTypeName())->GetType();
1480  return 0;
1481  } else {
1482  Error("GetExpectedType", "Did not find any leaves in %s",GetName());
1483  return 1;
1484  }
1485 }
1486 
1487 ////////////////////////////////////////////////////////////////////////////////
1488 /// Return pointer to the file where branch buffers reside, returns 0
1489 /// in case branch buffers reside in the same file as tree header.
1490 /// If mode is 1 the branch buffer file is recreated.
1491 
1493 {
1494  if (fDirectory) return fDirectory->GetFile();
1495 
1496  // check if a file with this name is in the list of Root files
1497  TFile *file = 0;
1498  {
1500  file = (TFile*)gROOT->GetListOfFiles()->FindObject(fFileName.Data());
1501  if (file) {
1502  fDirectory = file;
1503  return file;
1504  }
1505  }
1506 
1507  if (fFileName.Length() == 0) return 0;
1508 
1509  TString bFileName( GetRealFileName() );
1510 
1511  // Open file (new file if mode = 1)
1512  {
1513  TDirectory::TContext ctxt;
1514  if (mode) file = TFile::Open(bFileName, "recreate");
1515  else file = TFile::Open(bFileName);
1516  }
1517  if (!file) return 0;
1518  if (file->IsZombie()) {delete file; return 0;}
1519  fDirectory = (TDirectory*)file;
1520  return file;
1521 }
1522 
1523 ////////////////////////////////////////////////////////////////////////////////
1524 /// Return a fresh basket by either resusing an existing basket that needs
1525 /// to be drop (according to TTree::MemoryFull) or create a new one.
1526 
1528 {
1529  TBasket *basket = 0;
1530  if (GetTree()->MemoryFull(0)) {
1531  if (fNBaskets==1) {
1532  // Steal the existing basket
1533  Int_t oldindex = fBaskets.GetLast();
1534  basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1535  if (!basket) {
1536  fBaskets.SetLast(-2); // For recalculation of Last.
1537  oldindex = fBaskets.GetLast();
1538  basket = (TBasket*)fBaskets.UncheckedAt(oldindex);
1539  }
1540  if (basket && fBasketBytes[oldindex]!=0) {
1541  if (basket == fCurrentBasket) {
1542  fCurrentBasket = 0;
1543  fFirstBasketEntry = -1;
1544  fNextBasketEntry = -1;
1545  }
1546  fBaskets.AddAt(0,oldindex);
1547  fBaskets.SetLast(-1);
1548  fNBaskets = 0;
1549  } else {
1550  basket = fTree->CreateBasket(this);
1551  }
1552  } else if (fNBaskets == 0) {
1553  // There is nothing to drop!
1554  basket = fTree->CreateBasket(this);
1555  } else {
1556  // Memory is full and there is more than one basket,
1557  // Let DropBaskets do it job.
1558  DropBaskets();
1559  basket = fTree->CreateBasket(this);
1560  }
1561  } else {
1562  basket = fTree->CreateBasket(this);
1563  }
1564  return basket;
1565 }
1566 
1567 ////////////////////////////////////////////////////////////////////////////////
1568 /// Drops the cluster two behind the current cluster and returns a fresh basket
1569 /// by either reusing or creating a new one
1570 
1572 {
1573  TBasket *basket = 0;
1574 
1575  // If GetClusterIterator is called with a negative entry then GetStartEntry will be 0
1576  // So we need to check if we reach the zero before we have gone back (1-VirtualSize) clusters
1577  // if this is the case, we want to keep everything in memory so we return a new basket
1579  if (iter.GetStartEntry() == 0) {
1580  return fTree->CreateBasket(this);
1581  }
1582 
1583  // Iterate backwards (1-VirtualSize) clusters to reach cluster to be unloaded from memory,
1584  // skipped if VirtualSize > 0.
1585  for (Int_t j = 0; j < -fTree->GetMaxVirtualSize(); j++) {
1586  if (iter.Previous() == 0) {
1587  return fTree->CreateBasket(this);
1588  }
1589  }
1590 
1591  Int_t entryToUnload = iter.Previous();
1592  // Finds the basket to unload from memory. Since the basket should be close to current
1593  // basket, just iterate backwards until the correct basket is reached. This should
1594  // be fast as long as the number of baskets per cluster is small
1595  Int_t basketToUnload = fReadBasket;
1596  while (fBasketEntry[basketToUnload] != entryToUnload) {
1597  basketToUnload--;
1598  if (basketToUnload < 0) {
1599  return fTree->CreateBasket(this);
1600  }
1601  }
1602 
1603  // Retrieves the basket that is going to be unloaded from memory. If the basket did not
1604  // exist, create a new one
1605  basket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1606  if (basket) {
1607  fBaskets.AddAt(0, basketToUnload);
1608  --fNBaskets;
1609  } else {
1610  basket = fTree->CreateBasket(this);
1611  }
1612  ++basketToUnload;
1613 
1614  // Clear the rest of the baskets. While it would be ideal to reuse these baskets
1615  // for other baskets in the new cluster. It would require the function to go
1616  // beyond its current scope. In the ideal case when each cluster only has 1 basket
1617  // this will perform well
1618  iter.Next();
1619  while (fBasketEntry[basketToUnload] < iter.GetStartEntry()) {
1620  TBasket *oldbasket = (TBasket *)fBaskets.UncheckedAt(basketToUnload);
1621  if (oldbasket) {
1622  oldbasket->DropBuffers();
1623  delete oldbasket;
1624  fBaskets.AddAt(0, basketToUnload);
1625  --fNBaskets;
1626  }
1627  ++basketToUnload;
1628  }
1629  fBaskets.SetLast(-1);
1630  return basket;
1631 }
1632 
1633 ////////////////////////////////////////////////////////////////////////////////
1634 /// Return pointer to the 1st Leaf named name in thisBranch
1635 
1636 TLeaf* TBranch::GetLeaf(const char* name) const
1637 {
1638  Int_t i;
1639  for (i=0;i<fNleaves;i++) {
1640  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
1641  if (!strcmp(leaf->GetName(),name)) return leaf;
1642  }
1643  return 0;
1644 }
1645 
1646 ////////////////////////////////////////////////////////////////////////////////
1647 /// Get real file name
1648 
1650 {
1651  if (fFileName.Length()==0) {
1652  return fFileName;
1653  }
1654  TString bFileName = fFileName;
1655 
1656  // check if branch file name is absolute or a URL (e.g. /castor/...,
1657  // root://host/..., rfio:/path/...)
1658  char *bname = gSystem->ExpandPathName(fFileName.Data());
1659  if (!gSystem->IsAbsoluteFileName(bname) && !strstr(bname, ":/") && fTree && fTree->GetCurrentFile()) {
1660 
1661  // if not, get filename where tree header is stored
1662  const char *tfn = fTree->GetCurrentFile()->GetName();
1663 
1664  // If it is an archive file we need a special treatment
1665  TUrl arc(tfn);
1666  if (strlen(arc.GetAnchor()) > 0) {
1668  bFileName = arc.GetUrl();
1669  } else {
1670  // if this is an absolute path or a URL then prepend this path
1671  // to the branch file name
1672  char *tname = gSystem->ExpandPathName(tfn);
1673  if (gSystem->IsAbsoluteFileName(tname) || strstr(tname, ":/")) {
1674  bFileName = gSystem->DirName(tname);
1675  bFileName += "/";
1676  bFileName += fFileName;
1677  }
1678  delete [] tname;
1679  }
1680  }
1681  delete [] bname;
1682 
1683  return bFileName;
1684 }
1685 
1686 ////////////////////////////////////////////////////////////////////////////////
1687 /// Return all elements of one row unpacked in internal array fValues
1688 /// [Actually just returns 1 (?)]
1689 
1691 {
1692  return 1;
1693 }
1694 
1695 ////////////////////////////////////////////////////////////////////////////////
1696 /// Return whether this branch is in a mode where the object are decomposed
1697 /// or not (Also known as MakeClass mode).
1698 
1700 {
1701  // Regular TBranch and TBrancObject can not be in makeClass mode
1702 
1703  return kFALSE;
1704 }
1705 
1706 ////////////////////////////////////////////////////////////////////////////////
1707 /// Get our top-level parent branch in the tree.
1708 
1710 {
1711  if (fMother) return fMother;
1712 
1713  const TObjArray* array = fTree->GetListOfBranches();
1714  Int_t n = array->GetEntriesFast();
1715  for (Int_t i = 0; i < n; ++i) {
1716  TBranch* branch = (TBranch*) array->UncheckedAt(i);
1717  TBranch* parent = branch->GetSubBranch(this);
1718  if (parent) {
1719  const_cast<TBranch*>(this)->fMother = branch; // We can not yet use the 'mutable' keyword
1720  return branch;
1721  }
1722  }
1723  return 0;
1724 }
1725 
1726 ////////////////////////////////////////////////////////////////////////////////
1727 /// Find the parent branch of child.
1728 /// Return 0 if child is not in this branch hierarchy.
1729 
1731 {
1732  // Handle error condition, if the parameter is us, we cannot find the parent.
1733  if (this == child) {
1734  // Note: We cast away any const-ness of "this".
1735  return (TBranch*) this;
1736  }
1737 
1738  if (child->fParent) {
1739  return child->fParent;
1740  }
1741 
1742  Int_t len = fBranches.GetEntriesFast();
1743  for (Int_t i = 0; i < len; ++i) {
1744  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1745  if (!branch) {
1746  continue;
1747  }
1748  if (branch == child) {
1749  // We are the direct parent of child.
1750  const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
1751  // Note: We cast away any const-ness of "this".
1752  const_cast<TBranch*>(child)->fParent = (TBranch*)this; // We can not yet use the 'mutable' keyword
1753  return (TBranch*) this;
1754  }
1755  // FIXME: This is a tail-recursion!
1756  TBranch* parent = branch->GetSubBranch(child);
1757  if (parent) {
1758  return parent;
1759  }
1760  }
1761  // We failed to find the parent.
1762  return 0;
1763 }
1764 
1765 ////////////////////////////////////////////////////////////////////////////////
1766 /// Return total number of bytes in the branch (including current buffer)
1767 
1769 {
1770  TBufferFile b(TBuffer::kWrite, 10000);
1771  // This intentionally only store the TBranch part and thus slightly
1772  // under-estimate the space used.
1773  // Since the TBranchElement part contains pointers to other branches (branch count),
1774  // doing regular Streaming would end up including those and thus greatly over-estimate
1775  // the size used.
1776  const_cast<TBranch *>(this)->TBranch::Streamer(b);
1777 
1778  Long64_t totbytes = 0;
1779  if (fZipBytes > 0) totbytes = fTotBytes;
1780  return totbytes + b.Length();
1781 }
1782 
1783 ////////////////////////////////////////////////////////////////////////////////
1784 /// Return total number of bytes in the branch (excluding current buffer)
1785 /// if option ="*" includes all sub-branches of this branch too
1786 
1788 {
1789  Long64_t totbytes = fTotBytes;
1790  if (!option) return totbytes;
1791  if (option[0] != '*') return totbytes;
1792  //scan sub-branches
1793  Int_t len = fBranches.GetEntriesFast();
1794  for (Int_t i = 0; i < len; ++i) {
1795  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1796  if (branch) totbytes += branch->GetTotBytes();
1797  }
1798  return totbytes;
1799 }
1800 
1801 ////////////////////////////////////////////////////////////////////////////////
1802 /// Return total number of zip bytes in the branch
1803 /// if option ="*" includes all sub-branches of this branch too
1804 
1806 {
1807  Long64_t zipbytes = fZipBytes;
1808  if (!option) return zipbytes;
1809  if (option[0] != '*') return zipbytes;
1810  //scan sub-branches
1811  Int_t len = fBranches.GetEntriesFast();
1812  for (Int_t i = 0; i < len; ++i) {
1813  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1814  if (branch) zipbytes += branch->GetZipBytes();
1815  }
1816  return zipbytes;
1817 }
1818 
1819 ////////////////////////////////////////////////////////////////////////////////
1820 /// Returns the IO settings currently in use for this branch.
1821 
1823 {
1824  return fIOFeatures;
1825 }
1826 
1827 ////////////////////////////////////////////////////////////////////////////////
1828 /// Return kTRUE if an existing object in a TBranchObject must be deleted.
1829 
1831 {
1832  return TestBit(kAutoDelete);
1833 }
1834 
1835 ////////////////////////////////////////////////////////////////////////////////
1836 /// Return kTRUE if more than one leaf or browsables, kFALSE otherwise.
1837 
1839 {
1840  if (fNleaves > 1) {
1841  return kTRUE;
1842  }
1843  TList* browsables = const_cast<TBranch*>(this)->GetBrowsables();
1844  return browsables && browsables->GetSize();
1845 }
1846 
1847 ////////////////////////////////////////////////////////////////////////////////
1848 /// keep a maximum of fMaxEntries in memory
1849 
1851 {
1852  Int_t dentries = (Int_t) (fEntries - maxEntries);
1853  TBasket* basket = (TBasket*) fBaskets.UncheckedAt(0);
1854  if (basket) basket->MoveEntries(dentries);
1855  fEntries = maxEntries;
1856  fEntryNumber = maxEntries;
1857  //loop on sub branches
1859  for (Int_t i = 0; i < nb; ++i) {
1860  TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
1861  branch->KeepCircular(maxEntries);
1862  }
1863 }
1864 
1865 ////////////////////////////////////////////////////////////////////////////////
1866 /// Baskets associated to this branch are forced to be in memory.
1867 /// You can call TTree::SetMaxVirtualSize(maxmemory) to instruct
1868 /// the system that the total size of the imported baskets does not
1869 /// exceed maxmemory bytes.
1870 ///
1871 /// The function returns the number of baskets that have been put in memory.
1872 /// This method may be called to force all baskets of one or more branches
1873 /// in memory when random access to entries in this branch is required.
1874 /// See also TTree::LoadBaskets to load all baskets of all branches in memory.
1875 
1877 {
1878  Int_t nimported = 0;
1879  Int_t nbaskets = fWriteBasket;
1880  TFile *file = GetFile(0);
1881  if (!file) return 0;
1882  TBasket *basket;
1883  for (Int_t i=0;i<nbaskets;i++) {
1884  basket = (TBasket*)fBaskets.UncheckedAt(i);
1885  if (basket) continue;
1886  basket = GetFreshBasket();
1887  if (fBasketBytes[i] == 0) {
1888  fBasketBytes[i] = basket->ReadBasketBytes(fBasketSeek[i],file);
1889  }
1890  Int_t badread = basket->ReadBasketBuffers(fBasketSeek[i],fBasketBytes[i],file);
1891  if (badread) {
1892  Error("Loadbaskets","Error while reading basket buffer %d of branch %s",i,GetName());
1893  return -1;
1894  }
1895  ++fNBaskets;
1896  fBaskets.AddAt(basket,i);
1897  nimported++;
1898  }
1899  return nimported;
1900 }
1901 
1902 ////////////////////////////////////////////////////////////////////////////////
1903 /// Print TBranch parameters
1904 ///
1905 /// If options contains "basketsInfo" print the entry number, location and size
1906 /// of each baskets.
1907 
1908 void TBranch::Print(Option_t *option) const
1909 {
1910  const int kLINEND = 77;
1911  Float_t cx = 1;
1912 
1913  TString titleContent(GetTitle());
1914  if ( titleContent == GetName() ) {
1915  titleContent.Clear();
1916  }
1917 
1918  if (fLeaves.GetEntries() == 1) {
1919  if (titleContent.Length()>=2 && titleContent[titleContent.Length()-2]=='/' && isalpha(titleContent[titleContent.Length()-1])) {
1920  // The type is already encoded. Nothing to do.
1921  } else {
1922  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(0);
1923  if (titleContent.Length()) {
1924  titleContent.Prepend(" ");
1925  }
1926  // titleContent.Append("type: ");
1927  titleContent.Prepend(leaf->GetTypeName());
1928  }
1929  }
1930  Int_t titleLength = titleContent.Length();
1931 
1932  Int_t aLength = titleLength + strlen(GetName());
1933  aLength += (aLength / 54 + 1) * 80 + 100;
1934  if (aLength < 200) aLength = 200;
1935  char *bline = new char[aLength];
1936 
1937  Long64_t totBytes = GetTotalSize();
1938  if (fZipBytes) cx = (fTotBytes+0.00001)/fZipBytes;
1939  if (titleLength) snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName(),titleContent.Data());
1940  else snprintf(bline,aLength,"*Br%5d :%-9s : %-54s *",fgCount,GetName()," ");
1941  if (strlen(bline) > UInt_t(kLINEND)) {
1942  char *tmp = new char[strlen(bline)+1];
1943  if (titleLength) strlcpy(tmp, titleContent.Data(),strlen(bline)+1);
1944  snprintf(bline,aLength,"*Br%5d :%-9s : ",fgCount,GetName());
1945  int pos = strlen (bline);
1946  int npos = pos;
1947  int beg=0, end;
1948  while (beg < titleLength) {
1949  for (end=beg+1; end < titleLength-1; end ++)
1950  if (tmp[end] == ':') break;
1951  if (npos + end-beg+1 >= 78) {
1952  while (npos < kLINEND) {
1953  bline[pos ++] = ' ';
1954  npos ++;
1955  }
1956  bline[pos ++] = '*';
1957  bline[pos ++] = '\n';
1958  bline[pos ++] = '*';
1959  npos = 1;
1960  for (; npos < 12; npos ++)
1961  bline[pos ++] = ' ';
1962  bline[pos-2] = '|';
1963  }
1964  for (int n = beg; n <= end; n ++)
1965  bline[pos+n-beg] = tmp[n];
1966  pos += end-beg+1;
1967  npos += end-beg+1;
1968  beg = end+1;
1969  }
1970  while (npos < kLINEND) {
1971  bline[pos ++] = ' ';
1972  npos ++;
1973  }
1974  bline[pos ++] = '*';
1975  bline[pos] = '\0';
1976  delete[] tmp;
1977  }
1978  Printf("%s", bline);
1979 
1980  if (fTotBytes > 2000000000) {
1981  Printf("*Entries :%lld : Total Size=%11lld bytes File Size = %lld *",fEntries,totBytes,fZipBytes);
1982  } else {
1983  if (fZipBytes > 0) {
1984  Printf("*Entries :%9lld : Total Size=%11lld bytes File Size = %10lld *",fEntries,totBytes,fZipBytes);
1985  } else {
1986  if (fWriteBasket > 0) {
1987  Printf("*Entries :%9lld : Total Size=%11lld bytes All baskets in memory *",fEntries,totBytes);
1988  } else {
1989  Printf("*Entries :%9lld : Total Size=%11lld bytes One basket in memory *",fEntries,totBytes);
1990  }
1991  }
1992  }
1993  Printf("*Baskets :%9d : Basket Size=%11d bytes Compression= %6.2f *",fWriteBasket,fBasketSize,cx);
1994 
1995  if (strncmp(option,"basketsInfo",strlen("basketsInfo"))==0) {
1996  Int_t nbaskets = fWriteBasket;
1997  for (Int_t i=0;i<nbaskets;i++) {
1998  Printf("*Basket #%4d entry=%6lld pos=%6lld size=%5d",
1999  i, fBasketEntry[i], fBasketSeek[i], fBasketBytes[i]);
2000  }
2001  }
2002 
2003  Printf("*............................................................................*");
2004  delete [] bline;
2005  fgCount++;
2006 }
2007 
2008 ////////////////////////////////////////////////////////////////////////////////
2009 /// Loop on all leaves of this branch to read Basket buffer.
2010 
2012 {
2013  // fLeaves->ReadBasket(basket);
2014 }
2015 
2016 ////////////////////////////////////////////////////////////////////////////////
2017 /// Loop on all leaves of this branch to read Basket buffer.
2018 
2020 {
2021  for (Int_t i = 0; i < fNleaves; ++i) {
2022  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2023  leaf->ReadBasket(b);
2024  }
2025 }
2026 
2027 ////////////////////////////////////////////////////////////////////////////////
2028 /// Read zero leaves without the overhead of a loop.
2029 
2031 {
2032 }
2033 
2034 ////////////////////////////////////////////////////////////////////////////////
2035 /// Read one leaf without the overhead of a loop.
2036 
2038 {
2039  ((TLeaf*) fLeaves.UncheckedAt(0))->ReadBasket(b);
2040 }
2041 
2042 ////////////////////////////////////////////////////////////////////////////////
2043 /// Read two leaves without the overhead of a loop.
2044 
2046 {
2047  ((TLeaf*) fLeaves.UncheckedAt(0))->ReadBasket(b);
2048  ((TLeaf*) fLeaves.UncheckedAt(1))->ReadBasket(b);
2049 }
2050 
2051 ////////////////////////////////////////////////////////////////////////////////
2052 /// Loop on all leaves of this branch to fill Basket buffer.
2053 
2055 {
2056  for (Int_t i = 0; i < fNleaves; ++i) {
2057  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2058  leaf->FillBasket(b);
2059  }
2060 }
2061 
2062 ////////////////////////////////////////////////////////////////////////////////
2063 /// Refresh this branch using new information in b
2064 /// This function is called by TTree::Refresh
2065 
2067 {
2068  if (b==0) return;
2069 
2073  fMaxBaskets = b->fMaxBaskets;
2074  fEntries = b->fEntries;
2075  fTotBytes = b->fTotBytes;
2076  fZipBytes = b->fZipBytes;
2077  fReadBasket = 0;
2078  fReadEntry = -1;
2079  fFirstBasketEntry = -1;
2080  fNextBasketEntry = -1;
2081  fCurrentBasket = 0;
2082  delete [] fBasketBytes;
2083  delete [] fBasketEntry;
2084  delete [] fBasketSeek;
2088  Int_t i;
2089  for (i=0;i<fMaxBaskets;i++) {
2090  fBasketBytes[i] = b->fBasketBytes[i];
2091  fBasketEntry[i] = b->fBasketEntry[i];
2092  fBasketSeek[i] = b->fBasketSeek[i];
2093  }
2094  fBaskets.Delete();
2095  Int_t nbaskets = b->fBaskets.GetSize();
2096  fBaskets.Expand(nbaskets);
2097  // If the current fWritebasket is in memory, take it (just swap)
2098  // from the Tree being read
2100  fBaskets.AddAt(basket,fWriteBasket);
2101  if (basket) {
2102  fNBaskets = 1;
2103  --(b->fNBaskets);
2105  basket->SetBranch(this);
2106  }
2107 }
2108 
2109 ////////////////////////////////////////////////////////////////////////////////
2110 /// Reset a Branch.
2111 ///
2112 /// - Existing buffers are deleted.
2113 /// - Entries, max and min are reset.
2114 
2116 {
2117  fReadBasket = 0;
2118  fReadEntry = -1;
2119  fFirstBasketEntry = -1;
2120  fNextBasketEntry = -1;
2121  fCurrentBasket = 0;
2122  fWriteBasket = 0;
2123  fEntries = 0;
2124  fTotBytes = 0;
2125  fZipBytes = 0;
2126  fEntryNumber = 0;
2127 
2128  if (fBasketBytes) {
2129  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2130  fBasketBytes[i] = 0;
2131  }
2132  }
2133 
2134  if (fBasketEntry) {
2135  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2136  fBasketEntry[i] = 0;
2137  }
2138  }
2139 
2140  if (fBasketSeek) {
2141  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2142  fBasketSeek[i] = 0;
2143  }
2144  }
2145 
2146  fBaskets.Delete();
2147  fNBaskets = 0;
2148 }
2149 
2150 ////////////////////////////////////////////////////////////////////////////////
2151 /// Reset a Branch.
2152 ///
2153 /// - Existing buffers are deleted.
2154 /// - Entries, max and min are reset.
2155 
2157 {
2158  fReadBasket = 0;
2159  fReadEntry = -1;
2160  fFirstBasketEntry = -1;
2161  fNextBasketEntry = -1;
2162  fCurrentBasket = 0;
2163  fWriteBasket = 0;
2164  fEntries = 0;
2165  fTotBytes = 0;
2166  fZipBytes = 0;
2167  fEntryNumber = 0;
2168 
2169  if (fBasketBytes) {
2170  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2171  fBasketBytes[i] = 0;
2172  }
2173  }
2174 
2175  if (fBasketEntry) {
2176  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2177  fBasketEntry[i] = 0;
2178  }
2179  }
2180 
2181  if (fBasketSeek) {
2182  for (Int_t i = 0; i < fMaxBaskets; ++i) {
2183  fBasketSeek[i] = 0;
2184  }
2185  }
2186 
2187  TBasket *reusebasket = (TBasket*)fBaskets[fWriteBasket];
2188  if (reusebasket) {
2189  fBaskets[fWriteBasket] = 0;
2190  } else {
2191  reusebasket = (TBasket*)fBaskets[fReadBasket];
2192  if (reusebasket) {
2193  fBaskets[fReadBasket] = 0;
2194  }
2195  }
2196  fBaskets.Delete();
2197  if (reusebasket) {
2198  fNBaskets = 1;
2199  reusebasket->Reset();
2200  fBaskets[0] = reusebasket;
2201  } else {
2202  fNBaskets = 0;
2203  }
2204 }
2205 
2206 ////////////////////////////////////////////////////////////////////////////////
2207 /// Reset the address of the branch.
2208 
2210 {
2211  fAddress = 0;
2212 
2213  // Reset last read entry number, we have will had new user object now.
2214  fReadEntry = -1;
2215 
2216  for (Int_t i = 0; i < fNleaves; ++i) {
2217  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2218  leaf->SetAddress(0);
2219  }
2220 
2221  Int_t nbranches = fBranches.GetEntriesFast();
2222  for (Int_t i = 0; i < nbranches; ++i) {
2223  TBranch* abranch = (TBranch*) fBranches[i];
2224  // FIXME: This is a tail recursion.
2225  abranch->ResetAddress();
2226  }
2227 }
2228 
2229 ////////////////////////////////////////////////////////////////////////////////
2230 /// Static function resetting fgCount
2231 
2233 {
2234  fgCount = 0;
2235 }
2236 
2237 ////////////////////////////////////////////////////////////////////////////////
2238 /// Set address of this branch.
2239 
2240 void TBranch::SetAddress(void* addr)
2241 {
2242  if (TestBit(kDoNotProcess)) {
2243  return;
2244  }
2245  fReadEntry = -1;
2246  fFirstBasketEntry = -1;
2247  fNextBasketEntry = -1;
2248  fAddress = (char*) addr;
2249  for (Int_t i = 0; i < fNleaves; ++i) {
2250  TLeaf* leaf = (TLeaf*) fLeaves.UncheckedAt(i);
2251  Int_t offset = leaf->GetOffset();
2252  if (TestBit(kIsClone)) {
2253  offset = 0;
2254  }
2255  if (fAddress) leaf->SetAddress(fAddress + offset);
2256  else leaf->SetAddress(0);
2257  }
2258 }
2259 
2260 ////////////////////////////////////////////////////////////////////////////////
2261 /// Set the automatic delete bit.
2262 ///
2263 /// This bit is used by TBranchObject::ReadBasket to decide if an object
2264 /// referenced by a TBranchObject must be deleted or not before reading
2265 /// a new entry.
2266 ///
2267 /// If autodel is kTRUE, this existing object will be deleted, a new object
2268 /// created by the default constructor, then read from disk by the streamer.
2269 ///
2270 /// If autodel is kFALSE, the existing object is not deleted. Root assumes
2271 /// that the user is taking care of deleting any internal object or array
2272 /// (this can be done in the streamer).
2273 
2275 {
2276  if (autodel) {
2277  SetBit(kAutoDelete, 1);
2278  } else {
2279  SetBit(kAutoDelete, 0);
2280  }
2281 }
2282 
2283 ////////////////////////////////////////////////////////////////////////////////
2284 /// Set the basket size
2285 /// The function makes sure that the basket size is greater than fEntryOffsetlen
2286 
2288 {
2289  Int_t minsize = 100 + fName.Length();
2290  if (buffsize < minsize+fEntryOffsetLen) buffsize = minsize+fEntryOffsetLen;
2291  fBasketSize = buffsize;
2292  TBasket *basket = (TBasket*)fBaskets[fWriteBasket];
2293  if (basket) {
2294  basket->AdjustSize(fBasketSize);
2295  }
2296 }
2297 
2298 ////////////////////////////////////////////////////////////////////////////////
2299 /// Set address of this branch directly from a TBuffer to avoid streaming.
2300 ///
2301 /// Note: We do not take ownership of the buffer.
2302 
2304 {
2305  // Check this is possible
2306  if ( (fNleaves != 1)
2307  || (strcmp("TLeafObject",fLeaves.UncheckedAt(0)->ClassName())!=0) ) {
2308  Error("TBranch::SetAddress","Filling from a TBuffer can only be done with a not split object branch. Request ignored.");
2309  } else {
2310  fReadEntry = -1;
2311  fNextBasketEntry = -1;
2312  fFirstBasketEntry = -1;
2313  // Note: We do not take ownership of the buffer.
2314  fEntryBuffer = buf;
2315  }
2316 }
2317 
2318 ////////////////////////////////////////////////////////////////////////////////
2319 /// Set compression algorithm.
2320 
2322 {
2323  if (algorithm < 0 || algorithm >= ROOT::kUndefinedCompressionAlgorithm) algorithm = 0;
2324  if (fCompress < 0) {
2325  fCompress = 100 * algorithm + 1;
2326  } else {
2327  int level = fCompress % 100;
2328  fCompress = 100 * algorithm + level;
2329  }
2330 
2332  for (Int_t i=0;i<nb;i++) {
2333  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2334  branch->SetCompressionAlgorithm(algorithm);
2335  }
2336 }
2337 
2338 ////////////////////////////////////////////////////////////////////////////////
2339 /// Set compression level.
2340 
2342 {
2343  if (level < 0) level = 0;
2344  if (level > 99) level = 99;
2345  if (fCompress < 0) {
2346  fCompress = level;
2347  } else {
2348  int algorithm = fCompress / 100;
2349  if (algorithm >= ROOT::kUndefinedCompressionAlgorithm) algorithm = 0;
2350  fCompress = 100 * algorithm + level;
2351  }
2352 
2354  for (Int_t i=0;i<nb;i++) {
2355  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2356  branch->SetCompressionLevel(level);
2357  }
2358 }
2359 
2360 ////////////////////////////////////////////////////////////////////////////////
2361 /// Set compression settings.
2362 
2364 {
2365  fCompress = settings;
2366 
2368  for (Int_t i=0;i<nb;i++) {
2369  TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
2370  branch->SetCompressionSettings(settings);
2371  }
2372 }
2373 
2374 ////////////////////////////////////////////////////////////////////////////////
2375 /// Update the default value for the branch's fEntryOffsetLen if and only if
2376 /// it was already non zero (and the new value is not zero)
2377 /// If updateExisting is true, also update all the existing branches.
2378 
2379 void TBranch::SetEntryOffsetLen(Int_t newdefault, Bool_t updateExisting)
2380 {
2381  if (fEntryOffsetLen && newdefault) {
2382  fEntryOffsetLen = newdefault;
2383  }
2384  if (updateExisting) {
2385  TIter next( GetListOfBranches() );
2386  TBranch *b;
2387  while ( ( b = (TBranch*)next() ) ) {
2388  b->SetEntryOffsetLen( newdefault, kTRUE );
2389  }
2390  }
2391 }
2392 
2393 ////////////////////////////////////////////////////////////////////////////////
2394 /// Set the number of entries in this branch.
2395 
2397 {
2398  fEntries = entries;
2399  fEntryNumber = entries;
2400 }
2401 
2402 ////////////////////////////////////////////////////////////////////////////////
2403 /// Set file where this branch writes/reads its buffers.
2404 /// By default the branch buffers reside in the file where the
2405 /// Tree was created.
2406 /// If the file name where the tree was created is an absolute
2407 /// path name or an URL (e.g. /castor/... or root://host/...)
2408 /// and if the fname is not an absolute path name or an URL then
2409 /// the path of the tree file is prepended to fname to make the
2410 /// branch file relative to the tree file. In this case one can
2411 /// move the tree + all branch files to a different location in
2412 /// the file system and still access the branch files.
2413 /// The ROOT file will be connected only when necessary.
2414 /// If called by TBranch::Fill (via TBasket::WriteFile), the file
2415 /// will be created with the option "recreate".
2416 /// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2417 /// will be opened in read mode.
2418 /// To open a file in "update" mode or with a certain compression
2419 /// level, use TBranch::SetFile(TFile *file).
2420 
2422 {
2423  if (file == 0) file = fTree->GetCurrentFile();
2424  fDirectory = (TDirectory*)file;
2425  if (file == fTree->GetCurrentFile()) fFileName = "";
2426  else fFileName = file->GetName();
2427 
2428  if (file && fCompress == -1) {
2429  fCompress = file->GetCompressionLevel();
2430  }
2431 
2432  // Apply to all existing baskets.
2433  TIter nextb(GetListOfBaskets());
2434  TBasket *basket;
2435  while ((basket = (TBasket*)nextb())) {
2436  basket->SetParent(file);
2437  }
2438 
2439  // Apply to sub-branches as well.
2440  TIter next(GetListOfBranches());
2441  TBranch *branch;
2442  while ((branch = (TBranch*)next())) {
2443  branch->SetFile(file);
2444  }
2445 }
2446 
2447 ////////////////////////////////////////////////////////////////////////////////
2448 /// Set file where this branch writes/reads its buffers.
2449 /// By default the branch buffers reside in the file where the
2450 /// Tree was created.
2451 /// If the file name where the tree was created is an absolute
2452 /// path name or an URL (e.g. /castor/... or root://host/...)
2453 /// and if the fname is not an absolute path name or an URL then
2454 /// the path of the tree file is prepended to fname to make the
2455 /// branch file relative to the tree file. In this case one can
2456 /// move the tree + all branch files to a different location in
2457 /// the file system and still access the branch files.
2458 /// The ROOT file will be connected only when necessary.
2459 /// If called by TBranch::Fill (via TBasket::WriteFile), the file
2460 /// will be created with the option "recreate".
2461 /// If called by TBranch::GetEntry (via TBranch::GetBasket), the file
2462 /// will be opened in read mode.
2463 /// To open a file in "update" mode or with a certain compression
2464 /// level, use TBranch::SetFile(TFile *file).
2465 
2466 void TBranch::SetFile(const char* fname)
2467 {
2468  fFileName = fname;
2469  fDirectory = 0;
2470 
2471  //apply to sub-branches as well
2472  TIter next(GetListOfBranches());
2473  TBranch *branch;
2474  while ((branch = (TBranch*)next())) {
2475  branch->SetFile(fname);
2476  }
2477 }
2478 
2479 ////////////////////////////////////////////////////////////////////////////////
2480 /// Set the branch in a mode where the object are decomposed
2481 /// (Also known as MakeClass mode).
2482 /// Return whether the setting was possible (it is not possible for
2483 /// TBranch and TBranchObject).
2484 
2486 {
2487  // Regular TBranch and TBrancObject can not be in makeClass mode
2488  return kFALSE;
2489 }
2490 
2491 ////////////////////////////////////////////////////////////////////////////////
2492 /// Set object this branch is pointing to.
2493 
2494 void TBranch::SetObject(void * /* obj */)
2495 {
2496  if (TestBit(kDoNotProcess)) {
2497  return;
2498  }
2499  Warning("SetObject","is not supported in TBranch objects");
2500 }
2501 
2502 ////////////////////////////////////////////////////////////////////////////////
2503 /// Set branch status to Process or DoNotProcess.
2504 
2506 {
2507  if (status) ResetBit(kDoNotProcess);
2508  else SetBit(kDoNotProcess);
2509 }
2510 
2511 ////////////////////////////////////////////////////////////////////////////////
2512 /// Stream a class object
2513 
2514 void TBranch::Streamer(TBuffer& b)
2515 {
2516  if (b.IsReading()) {
2517  UInt_t R__s, R__c;
2518  fTree = 0; // Will be set by TTree::Streamer
2519  fAddress = 0;
2520  gROOT->SetReadingObject(kTRUE);
2521 
2522  // Reset transients.
2524  fCurrentBasket = 0;
2525  fFirstBasketEntry = -1;
2526  fNextBasketEntry = -1;
2527 
2528  Version_t v = b.ReadVersion(&R__s, &R__c);
2529  if (v > 9) {
2530  b.ReadClassBuffer(TBranch::Class(), this, v, R__s, R__c);
2531 
2532  if (fWriteBasket>=fBaskets.GetSize()) {
2534  }
2535  fDirectory = 0;
2537  for (Int_t i=0;i<fNleaves;i++) {
2538  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2539  leaf->SetBranch(this);
2540  }
2541 
2543  for (Int_t j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2544  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2545  if (bk) {
2546  bk->SetBranch(this);
2547  // GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2548  ++n;
2549  }
2550  }
2551  if (fWriteBasket >= fMaxBaskets) {
2552  //old versions may need this fix
2557 
2558  }
2560  gROOT->SetReadingObject(kFALSE);
2561  if (IsA() == TBranch::Class()) {
2562  if (fNleaves == 0) {
2564  } else if (fNleaves == 1) {
2566  } else if (fNleaves == 2) {
2568  } else {
2570  }
2571  }
2572  return;
2573  }
2574  //====process old versions before automatic schema evolution
2575  Int_t n,i,j,ijunk;
2576  if (v > 5) {
2577  Stat_t djunk;
2578  TNamed::Streamer(b);
2579  if (v > 7) TAttFill::Streamer(b);
2580  b >> fCompress;
2581  b >> fBasketSize;
2582  b >> fEntryOffsetLen;
2583  b >> fWriteBasket;
2584  b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2585  b >> fOffset;
2586  b >> fMaxBaskets;
2587  if (v > 6) b >> fSplitLevel;
2588  b >> djunk; fEntries = (Long64_t)djunk;
2589  b >> djunk; fTotBytes = (Long64_t)djunk;
2590  b >> djunk; fZipBytes = (Long64_t)djunk;
2591 
2592  fBranches.Streamer(b);
2593  fLeaves.Streamer(b);
2594  fBaskets.Streamer(b);
2598  Char_t isArray;
2599  b >> isArray;
2600  b.ReadFastArray(fBasketBytes,fMaxBaskets);
2601  b >> isArray;
2602  for (i=0;i<fMaxBaskets;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2603  b >> isArray;
2604  for (i=0;i<fMaxBaskets;i++) {
2605  if (isArray == 2) b >> fBasketSeek[i];
2606  else {Int_t bsize; b >> bsize; fBasketSeek[i] = (Long64_t)bsize;};
2607  }
2608  fFileName.Streamer(b);
2609  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2610  fDirectory = 0;
2612  for (i=0;i<fNleaves;i++) {
2613  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2614  leaf->SetBranch(this);
2615  }
2617  for (j=fWriteBasket,n=0;j>=0 && n<fNBaskets;--j) {
2618  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2619  if (bk) {
2620  bk->SetBranch(this);
2621  //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2622  ++n;
2623  }
2624  }
2625  if (fWriteBasket >= fMaxBaskets) {
2626  //old versions may need this fix
2628  fBasketBytes[fWriteBasket] = fBasketBytes[fWriteBasket-1];
2630  fBasketSeek [fWriteBasket] = fBasketSeek [fWriteBasket-1];
2631 
2632  }
2633  // Check Byte Count is not needed since it was done in ReadBuffer
2634  if (!fSplitLevel && fBranches.GetEntriesFast()) fSplitLevel = 1;
2635  gROOT->SetReadingObject(kFALSE);
2636  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2637  if (IsA() == TBranch::Class()) {
2638  if (fNleaves == 0) {
2640  } else if (fNleaves == 1) {
2642  } else if (fNleaves == 2) {
2644  } else {
2646  }
2647  }
2648  return;
2649  }
2650  //====process very old versions
2651  Stat_t djunk;
2652  TNamed::Streamer(b);
2653  b >> fCompress;
2654  b >> fBasketSize;
2655  b >> fEntryOffsetLen;
2656  b >> fMaxBaskets;
2657  b >> fWriteBasket;
2658  b >> ijunk; fEntryNumber = (Long64_t)ijunk;
2659  b >> djunk; fEntries = (Long64_t)djunk;
2660  b >> djunk; fTotBytes = (Long64_t)djunk;
2661  b >> djunk; fZipBytes = (Long64_t)djunk;
2662  b >> fOffset;
2663  fBranches.Streamer(b);
2664  fLeaves.Streamer(b);
2666  for (i=0;i<fNleaves;i++) {
2667  TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(i);
2668  leaf->SetBranch(this);
2669  }
2670  fBaskets.Streamer(b);
2671  Int_t nbaskets = fBaskets.GetEntries();
2672  for (j=fWriteBasket,n=0;j>0 && n<nbaskets;--j) {
2673  TBasket *bk = (TBasket*)fBaskets.UncheckedAt(j);
2674  if (bk) {
2675  bk->SetBranch(this);
2676  //GetTree()->IncrementTotalBuffers(bk->GetBufferSize());
2677  ++n;
2678  }
2679  }
2681  b >> n;
2682  for (i=0;i<n;i++) {b >> ijunk; fBasketEntry[i] = ijunk;}
2684  if (v > 4) {
2685  n = b.ReadArray(fBasketBytes);
2686  } else {
2687  for (n=0;n<fMaxBaskets;n++) fBasketBytes[n] = 0;
2688  }
2689  if (v < 2) {
2691  for (n=0;n<fWriteBasket;n++) {
2692  TBasket *basket = GetBasket(n);
2693  fBasketSeek[n] = basket ? basket->GetSeekKey() : 0;
2694  }
2695  } else {
2697  b >> n;
2698  for (n=0;n<fMaxBaskets;n++) {
2699  Int_t aseek;
2700  b >> aseek;
2701  fBasketSeek[n] = Long64_t(aseek);
2702  }
2703  }
2704  if (v > 2) {
2705  fFileName.Streamer(b);
2706  }
2707  fDirectory = 0;
2708  if (v < 4) SetAutoDelete(kTRUE);
2710  gROOT->SetReadingObject(kFALSE);
2711  b.CheckByteCount(R__s, R__c, TBranch::IsA());
2712  //====end of old versions
2713  if (IsA() == TBranch::Class()) {
2714  if (fNleaves == 0) {
2716  } else if (fNleaves == 1) {
2718  } else if (fNleaves == 2) {
2720  } else {
2722  }
2723  }
2724  } else {
2725  Int_t maxBaskets = fMaxBaskets;
2727  Int_t lastBasket = fMaxBaskets;
2728  if (fMaxBaskets < 10) fMaxBaskets = 10;
2729 
2730  TBasket **stash = new TBasket *[lastBasket];
2731  for (Int_t i = 0; i < lastBasket; ++i) {
2732  TBasket *ba = (TBasket *)fBaskets.UncheckedAt(i);
2733  if (ba && (fBasketBytes[i] || ba->GetNevBuf()==0)) {
2734  // Already on disk or empty.
2735  stash[i] = ba;
2736  fBaskets[i] = nullptr;
2737  } else {
2738  stash[i] = nullptr;
2739  }
2740  }
2741 
2742  b.WriteClassBuffer(TBranch::Class(), this);
2743 
2744  for (Int_t i = 0; i < lastBasket; ++i) {
2745  if (stash[i]) fBaskets[i] = stash[i];
2746  }
2747 
2748  delete[] stash;
2749  fMaxBaskets = maxBaskets;
2750  }
2751 }
2752 
2753 ////////////////////////////////////////////////////////////////////////////////
2754 /// Write the current basket to disk and return the number of bytes
2755 /// written to the file.
2756 
2758 {
2759  Int_t nevbuf = basket->GetNevBuf();
2760  if (fEntryOffsetLen > 10 && (4*nevbuf) < fEntryOffsetLen ) {
2761  // Make sure that the fEntryOffset array does not stay large unnecessarily.
2762  fEntryOffsetLen = nevbuf < 3 ? 10 : 4*nevbuf; // assume some fluctuations.
2763  } else if (fEntryOffsetLen && nevbuf > fEntryOffsetLen) {
2764  // Increase the array ...
2765  fEntryOffsetLen = 2*nevbuf; // assume some fluctuations.
2766  }
2767 
2768  // Note: captures `basket`, `where`, and `this` by value; modifies the TBranch and basket,
2769  // as we make a copy of the pointer. We cannot capture `basket` by reference as the pointer
2770  // itself might be modified after `WriteBasketImpl` exits.
2771  auto doUpdates = [=]() {
2772  Int_t nout = basket->WriteBuffer(); // Write buffer
2773  if (nout < 0) Error("TBranch::WriteBasketImpl", "basket's WriteBuffer failed.\n");
2774  fBasketBytes[where] = basket->GetNbytes();
2775  fBasketSeek[where] = basket->GetSeekKey();
2776  Int_t addbytes = basket->GetObjlen() + basket->GetKeylen();
2777  TBasket *reusebasket = 0;
2778  if (nout>0) {
2779  // The Basket was written so we can now safely reuse it.
2780  fBaskets[where] = 0;
2781 
2782  reusebasket = basket;
2783  reusebasket->Reset();
2784 
2785  fZipBytes += nout;
2786  fTotBytes += addbytes;
2787  fTree->AddTotBytes(addbytes);
2788  fTree->AddZipBytes(nout);
2789  }
2790 
2791  if (where==fWriteBasket) {
2792  ++fWriteBasket;
2793  if (fWriteBasket >= fMaxBaskets) {
2795  }
2796  if (reusebasket && reusebasket == fCurrentBasket) {
2797  // The 'current' basket has Reset, so if we need it we will need
2798  // to reload it.
2799  fCurrentBasket = 0;
2800  fFirstBasketEntry = -1;
2801  fNextBasketEntry = -1;
2802  }
2803  fBaskets.AddAtAndExpand(reusebasket,fWriteBasket);
2805  } else {
2806  --fNBaskets;
2807  fBaskets[where] = 0;
2808  basket->DropBuffers();
2809  if (basket == fCurrentBasket) {
2810  fCurrentBasket = 0;
2811  fFirstBasketEntry = -1;
2812  fNextBasketEntry = -1;
2813  }
2814  delete basket;
2815  }
2816  return nout;
2817  };
2818  if (imtHelper) {
2819  imtHelper->Run(doUpdates);
2820  return 0;
2821  } else {
2822  return doUpdates();
2823  }
2824 }
2825 
2826 ////////////////////////////////////////////////////////////////////////////////
2827 ///set the first entry number (case of TBranchSTL)
2828 
2830 {
2831  fFirstEntry = entry;
2832  fEntries = 0;
2833  fEntryNumber = entry;
2834  if( fBasketEntry )
2835  fBasketEntry[0] = entry;
2836  for( Int_t i = 0; i < fBranches.GetEntriesFast(); ++i )
2837  ((TBranch*)fBranches[i])->SetFirstEntry( entry );
2838 }
2839 
2840 ////////////////////////////////////////////////////////////////////////////////
2841 /// If the branch address is not set, we set all addresses starting with
2842 /// the top level parent branch.
2843 
2845 {
2846  // Nothing to do for regular branch, the TLeaf already did it.
2847 }
2848 
2849 ////////////////////////////////////////////////////////////////////////////////
2850 /// Refresh the value of fDirectory (i.e. where this branch writes/reads its buffers)
2851 /// with the current value of fTree->GetCurrentFile unless this branch has been
2852 /// redirected to a different file. Also update the sub-branches.
2853 
2855 {
2857  if (fFileName.Length() == 0) {
2858  fDirectory = file;
2859 
2860  // Apply to all existing baskets.
2861  TIter nextb(GetListOfBaskets());
2862  TBasket *basket;
2863  while ((basket = (TBasket*)nextb())) {
2864  basket->SetParent(file);
2865  }
2866  }
2867 
2868  // Apply to sub-branches as well.
2869  TIter next(GetListOfBranches());
2870  TBranch *branch;
2871  while ((branch = (TBranch*)next())) {
2872  branch->UpdateFile();
2873  }
2874 }
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:1266
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:2303
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:815
Long64_t GetNextEntry()
Definition: TTree.h:271
A TLeaf for a 64 bit floating point data type.
Definition: TLeafD.h:26
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
void Update(Int_t newlast)
Definition: TBasket.h:133
virtual void UpdateFile()
Refresh the value of fDirectory (i.e.
Definition: TBranch.cxx:2854
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:2240
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:193
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:1133
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:1473
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:1805
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:887
TBranch * GetSubBranch(const TBranch *br) const
Find the parent branch of child.
Definition: TBranch.cxx:1730
void SetCompressionAlgorithm(Int_t algorithm=0)
Set compression algorithm.
Definition: TBranch.cxx:2321
virtual void DropBaskets(Option_t *option="")
Loop on all branch baskets.
Definition: TBranch.cxx:716
virtual void SetStatus(Bool_t status=1)
Set branch status to Process or DoNotProcess.
Definition: TBranch.cxx:2505
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:2045
virtual void SetFirstEntry(Long64_t entry)
set the first entry number (case of TBranchSTL)
Definition: TBranch.cxx:2829
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:2341
virtual TLeaf * FindLeaf(const char *name)
Find the leaf corresponding to the name &#39;searchname&#39;.
Definition: TBranch.cxx:1032
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:1838
TBasket * fCurrentBasket
! Pointer to the current basket.
Definition: TBranch.h:94
TDirectory * GetDirectory() const
Definition: TBranch.h:173
virtual TBasket * CreateBasket(TBranch *)
Create a basket for this tree and given branch.
Definition: TTree.cxx:3524
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:402
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:1571
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:1649
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:685
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:1690
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:1527
virtual void SetupAddresses()
If the branch address is not set, we set all addresses starting with the top level parent branch...
Definition: TBranch.cxx:2844
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:1254
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:1787
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:658
virtual TLeaf * GetLeaf(const char *name) const
Return pointer to the 1st Leaf named name in thisBranch.
Definition: TBranch.cxx:1636
void ReadLeavesImpl(TBuffer &b)
Loop on all leaves of this branch to read Basket buffer.
Definition: TBranch.cxx:2019
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:2485
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:986
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:5161
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:1830
virtual Int_t GetLenType() const
Definition: TLeaf.h:76
TObjArray * GetListOfBranches()
Definition: TBranch.h:194
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:102
TBasket * GetBasket(Int_t basket)
Return pointer to basket basketnumber in this Branch.
Definition: TBranch.cxx:1174
Int_t fMaxBaskets
Maximum number of Baskets so far.
Definition: TBranch.h:86
virtual void SetAddress(void *add=0)
Definition: TLeaf.h:126
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:2156
Int_t Fill()
Definition: TBranch.h:155
const char * GetIconName() const
Return icon name depending on type of branch.
Definition: TBranch.cxx:1274
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:2054
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:1401
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:2363
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:2396
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:5172
Ssiz_t Length() const
Definition: TString.h:386
virtual void MakeFree(Long64_t first, Long64_t last)
Mark unused bytes on the file.
Definition: TFile.cxx:1410
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:1301
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 ReadFastArray(Bool_t *b, Int_t n)=0
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:1850
virtual Long64_t GetBasketSeek(Int_t basket) const
Return address of basket in the file.
Definition: TBranch.cxx:1244
#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:105
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:2287
TString & Remove(Ssiz_t pos)
Definition: TString.h:619
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
int Ssiz_t
Definition: RtypesCore.h:63
void SetAnchor(const char *anchor)
Definition: TUrl.h:89
Long64_t Next()
Move on to the next cluster and return the starting entry of this next cluster.
Definition: TTree.cxx: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:2209
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:2757
#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:1822
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:195
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:95
Int_t BackFill()
Loop on all leaves of this branch to back fill Basket buffer.
Definition: TBranch.cxx:639
virtual void SetUnsigned()
Definition: TLeaf.h:107
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:1768
virtual void ReadBasketExport(TBuffer &, TClonesArray *, Int_t)
Definition: TLeaf.h:96
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:2379
virtual void SetFile(TFile *file=0)
Set file where this branch writes/reads its buffers.
Definition: TBranch.cxx:2421
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:2066
Long64_t GetEntries() const
Definition: TBranch.h:199
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:2030
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:2115
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:1876
#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:1699
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:200
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
TFileCacheRead * GetCacheRead(TObject *tree=0) const
Return a pointer to the current read cache.
Definition: TFile.cxx:1214
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:2011
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:1087
TObjArray fBranches
-> List of Branches of this branch
Definition: TBranch.h:99
Option_t * GetDrawOption() const
Get option used by the graphics system to draw this object.
Definition: TBrowser.h:104
void SetNevBufSize(Int_t n)
Definition: TBasket.h: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:1492
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:1908
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
Int_t GetCompressionLevel() const
Definition: TFile.h:372
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:1709
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:784
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:2494
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 Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual Int_t ReadArray(Bool_t *&b)=0
virtual TObjArray * GetListOfLeaves()
Definition: TTree.h:408
void ReadLeaves1Impl(TBuffer &b)
Read one leaf without the overhead of a loop.
Definition: TBranch.cxx:2037
virtual void Close(Option_t *option="")
Close a file.
Definition: TFile.cxx:916
static void ResetCount()
Static function resetting fgCount.
Definition: TBranch.cxx:2232
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:2274
char * fAddress
! Address of 1st leaf (variable or object)
Definition: TBranch.h:108