Logo ROOT   6.10/09
Reference Guide
TChain.cxx
Go to the documentation of this file.
1 // @(#)root/tree:$Id$
2 // Author: Rene Brun 03/02/97
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 /** \class TChain
13 \ingroup tree
14 
15 A chain is a collection of files containing TTree objects.
16 When the chain is created, the first parameter is the default name
17 for the Tree to be processed later on.
18 
19 Enter a new element in the chain via the TChain::Add function.
20 Once a chain is defined, one can use the normal TTree functions
21 to Draw,Scan,etc.
22 
23 Use TChain::SetBranchStatus to activate one or more branches for all
24 the trees in the chain.
25 */
26 
27 #include "TChain.h"
28 
29 #include "TBranch.h"
30 #include "TBrowser.h"
31 #include "TChainElement.h"
32 #include "TClass.h"
33 #include "TCut.h"
34 #include "TError.h"
35 #include "TMath.h"
36 #include "TFile.h"
37 #include "TFileInfo.h"
38 #include "TFriendElement.h"
39 #include "TLeaf.h"
40 #include "TList.h"
41 #include "TObjString.h"
42 #include "TPluginManager.h"
43 #include "TROOT.h"
44 #include "TRegexp.h"
45 #include "TSelector.h"
46 #include "TSystem.h"
47 #include "TTree.h"
48 #include "TTreeCache.h"
49 #include "TUrl.h"
50 #include "TVirtualIndex.h"
51 #include "TEventList.h"
52 #include "TEntryList.h"
53 #include "TEntryListFromFile.h"
54 #include "TFileStager.h"
55 #include "TFilePrefetch.h"
56 #include "TVirtualMutex.h"
57 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Default constructor.
62 
64 : TTree()
65 , fTreeOffsetLen(100)
66 , fNtrees(0)
67 , fTreeNumber(-1)
68 , fTreeOffset(0)
69 , fCanDeleteRefs(kFALSE)
70 , fTree(0)
71 , fFile(0)
72 , fFiles(0)
73 , fStatus(0)
74 , fProofChain(0)
75 {
76  fTreeOffset = new Long64_t[fTreeOffsetLen];
77  fFiles = new TObjArray(fTreeOffsetLen);
78  fStatus = new TList();
79  fTreeOffset[0] = 0;
80  if (gDirectory) gDirectory->Remove(this);
81  gROOT->GetListOfSpecials()->Add(this);
82  fFile = 0;
83  fDirectory = 0;
84 
85  // Reset PROOF-related bits
86  ResetBit(kProofUptodate);
87  ResetBit(kProofLite);
88 
89  // Add to the global list
90  gROOT->GetListOfDataSets()->Add(this);
91 
92  // Make sure we are informed if the TFile is deleted.
94  gROOT->GetListOfCleanups()->Add(this);
95 }
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// Create a chain.
99 ///
100 /// A TChain is a collection of TFile objects.
101 /// the first parameter "name" is the name of the TTree object
102 /// in the files added with Add.
103 /// Use TChain::Add to add a new element to this chain.
104 ///
105 /// In case the Tree is in a subdirectory, do, eg:
106 /// ~~~ {.cpp}
107 /// TChain ch("subdir/treename");
108 /// ~~~
109 /// Example:
110 /// Suppose we have 3 files f1.root, f2.root and f3.root. Each file
111 /// contains a TTree object named "T".
112 /// ~~~ {.cpp}
113 /// TChain ch("T"); creates a chain to process a Tree called "T"
114 /// ch.Add("f1.root");
115 /// ch.Add("f2.root");
116 /// ch.Add("f3.root");
117 /// ch.Draw("x");
118 /// ~~~
119 /// The Draw function above will process the variable "x" in Tree "T"
120 /// reading sequentially the 3 files in the chain ch.
121 ///
122 /// The TChain data structure:
123 ///
124 /// Each TChainElement has a name equal to the tree name of this TChain
125 /// and a title equal to the file name. So, to loop over the
126 /// TFiles that have been added to this chain:
127 /// ~~~ {.cpp}
128 /// TObjArray *fileElements=chain->GetListOfFiles();
129 /// TIter next(fileElements);
130 /// TChainElement *chEl=0;
131 /// while (( chEl=(TChainElement*)next() )) {
132 /// TFile f(chEl->GetTitle());
133 /// ... do something with f ...
134 /// }
135 /// ~~~
136 
137 TChain::TChain(const char* name, const char* title)
138 :TTree(name, title)
139 , fTreeOffsetLen(100)
140 , fNtrees(0)
141 , fTreeNumber(-1)
142 , fTreeOffset(0)
143 , fCanDeleteRefs(kFALSE)
144 , fTree(0)
145 , fFile(0)
146 , fFiles(0)
147 , fStatus(0)
148 , fProofChain(0)
149 {
150  //
151  //*-*
152 
155  fStatus = new TList();
156  fTreeOffset[0] = 0;
157  if (gDirectory) gDirectory->Remove(this);
158  gROOT->GetListOfSpecials()->Add(this);
159  fFile = 0;
160  fDirectory = 0;
161 
162  // Reset PROOF-related bits
165 
166  // Add to the global list
167  gROOT->GetListOfDataSets()->Add(this);
168 
169  // Make sure we are informed if the TFile is deleted.
171  gROOT->GetListOfCleanups()->Add(this);
172 }
173 
174 ////////////////////////////////////////////////////////////////////////////////
175 /// Destructor.
176 
178 {
179  bool rootAlive = gROOT && !gROOT->TestBit(TObject::kInvalidObject);
180 
181  if (rootAlive) {
183  gROOT->GetListOfCleanups()->Remove(this);
184  }
185 
187  fStatus->Delete();
188  delete fStatus;
189  fStatus = 0;
190  fFiles->Delete();
191  delete fFiles;
192  fFiles = 0;
193 
194  //first delete cache if exists
195  if (fFile && fFile->GetCacheRead(fTree)) {
196  delete fFile->GetCacheRead(fTree);
197  fFile->SetCacheRead(0, fTree);
198  }
199 
200  delete fFile;
201  fFile = 0;
202  // Note: We do *not* own the tree.
203  fTree = 0;
204  delete[] fTreeOffset;
205  fTreeOffset = 0;
206 
207  if (rootAlive) gROOT->GetListOfSpecials()->Remove(this);
208 
209  // Remove from the global list
210  if (rootAlive) gROOT->GetListOfDataSets()->Remove(this);
211 
212  // This is the same as fFile, don't delete it a second time.
213  fDirectory = 0;
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Add all files referenced by the passed chain to this chain.
218 /// The function returns the total number of files connected.
219 
221 {
222  if (!chain) return 0;
223 
224  // Check for enough space in fTreeOffset.
225  if ((fNtrees + chain->GetNtrees()) >= fTreeOffsetLen) {
226  fTreeOffsetLen += 2 * chain->GetNtrees();
227  Long64_t* trees = new Long64_t[fTreeOffsetLen];
228  for (Int_t i = 0; i <= fNtrees; i++) {
229  trees[i] = fTreeOffset[i];
230  }
231  delete[] fTreeOffset;
232  fTreeOffset = trees;
233  }
234  chain->GetEntries(); //to force the computation of nentries
235  TIter next(chain->GetListOfFiles());
236  Int_t nf = 0;
237  TChainElement* element = 0;
238  while ((element = (TChainElement*) next())) {
239  Long64_t nentries = element->GetEntries();
242  } else {
244  }
245  fNtrees++;
246  fEntries += nentries;
247  TChainElement* newelement = new TChainElement(element->GetName(), element->GetTitle());
248  newelement->SetPacketSize(element->GetPacketSize());
249  newelement->SetNumberEntries(nentries);
250  fFiles->Add(newelement);
251  nf++;
252  }
253  if (fProofChain)
254  // This updates the proxy chain when we will really use PROOF
256 
257  return nf;
258 }
259 
260 ////////////////////////////////////////////////////////////////////////////////
261 /// Add a new file to this chain.
262 ///
263 /// Argument name may have either of two formats. The first:
264 /// ~~~ {.cpp}
265 /// [//machine]/path/file_name.root[/tree_name]
266 /// ~~~
267 /// If tree_name is missing the chain name will be assumed.
268 /// Wildcard treatment is triggered by the any of the special characters []*?
269 /// which may be used in the file name, eg. specifying "xxx*.root" adds
270 /// all files starting with xxx in the current file system directory.
271 ///
272 /// Alternatively name may have the format of a url, eg.
273 /// ~~~ {.cpp}
274 /// root://machine/path/file_name.root
275 /// or root://machine/path/file_name.root/tree_name
276 /// or root://machine/path/file_name.root/tree_name?query
277 /// ~~~
278 /// where "query" is to be interpreted by the remote server. Wildcards may be
279 /// supported in urls, depending on the protocol plugin and the remote server.
280 /// http or https urls can contain a query identifier without tree_name, but
281 /// generally urls can not be written with them because of ambiguity with the
282 /// wildcard character. (Also see the documentation for TChain::AddFile,
283 /// which does not support wildcards but allows the url to contain query)
284 ///
285 /// NB. To add all the files of a TChain to a chain, use Add(TChain *chain).
286 ///
287 /// A. if nentries <= 0, the file is connected and the tree header read
288 /// in memory to get the number of entries.
289 ///
290 /// B. if (nentries > 0, the file is not connected, nentries is assumed to be
291 /// the number of entries in the file. In this case, no check is made that
292 /// the file exists and the Tree existing in the file. This second mode
293 /// is interesting in case the number of entries in the file is already stored
294 /// in a run data base for example.
295 ///
296 /// C. if (nentries == TTree::kMaxEntries) (default), the file is not connected.
297 /// the number of entries in each file will be read only when the file
298 /// will need to be connected to read an entry.
299 /// This option is the default and very efficient if one process
300 /// the chain sequentially. Note that in case TChain::GetEntry(entry)
301 /// is called and entry refers to an entry in the 3rd file, for example,
302 /// this forces the Tree headers in the first and second file
303 /// to be read to find the number of entries in these files.
304 /// Note that if one calls TChain::GetEntriesFast() after having created
305 /// a chain with this default, GetEntriesFast will return TTree::kMaxEntries!
306 /// TChain::GetEntries will force of the Tree headers in the chain to be
307 /// read to read the number of entries in each Tree.
308 ///
309 /// D. The TChain data structure
310 /// Each TChainElement has a name equal to the tree name of this TChain
311 /// and a title equal to the file name. So, to loop over the
312 /// TFiles that have been added to this chain:
313 /// ~~~ {.cpp}
314 /// TObjArray *fileElements=chain->GetListOfFiles();
315 /// TIter next(fileElements);
316 /// TChainElement *chEl=0;
317 /// while (( chEl=(TChainElement*)next() )) {
318 /// TFile f(chEl->GetTitle());
319 /// ... do something with f ...
320 /// }
321 /// ~~~
322 /// Return value:
323 ///
324 /// - If nentries>0 (including the default of TTree::kMaxEntries) and no
325 /// wildcarding is used, ALWAYS returns 1 without regard to whether
326 /// the file exists or contains the correct tree.
327 ///
328 /// - If wildcarding is used, regardless of the value of nentries,
329 /// returns the number of files matching the name without regard to
330 /// whether they contain the correct tree.
331 ///
332 /// - If nentries<=0 and wildcarding is not used, return 1 if the file
333 /// exists and contains the correct tree and 0 otherwise.
334 
335 Int_t TChain::Add(const char* name, Long64_t nentries /* = TTree::kMaxEntries */)
336 {
337  TString basename, treename, query, suffix;
338  ParseTreeFilename(name, basename, treename, query, suffix, kTRUE);
339 
340  // case with one single file
341  if (!basename.MaybeWildcard()) {
342  return AddFile(name, nentries);
343  }
344 
345  // wildcarding used in name
346  Int_t nf = 0;
347 
348  Int_t slashpos = basename.Last('/');
349  TString directory;
350  if (slashpos>=0) {
351  directory = basename(0,slashpos); // Copy the directory name
352  basename.Remove(0,slashpos+1); // and remove it from basename
353  } else {
354  directory = gSystem->UnixPathName(gSystem->WorkingDirectory());
355  }
356 
357  const char *file;
358  const char *epath = gSystem->ExpandPathName(directory.Data());
359  void *dir = gSystem->OpenDirectory(epath);
360  delete [] epath;
361  if (dir) {
362  //create a TList to store the file names (not yet sorted)
363  TList l;
364  TRegexp re(basename,kTRUE);
365  while ((file = gSystem->GetDirEntry(dir))) {
366  if (!strcmp(file,".") || !strcmp(file,"..")) continue;
367  TString s = file;
368  if ( (basename!=file) && s.Index(re) == kNPOS) continue;
369  l.Add(new TObjString(file));
370  }
371  gSystem->FreeDirectory(dir);
372  //sort the files in alphanumeric order
373  l.Sort();
374  TIter next(&l);
375  TObjString *obj;
376  while ((obj = (TObjString*)next())) {
377  file = obj->GetName();
378  nf += AddFile(TString::Format("%s/%s%s",directory.Data(),file,suffix.Data()),nentries);
379  }
380  l.Delete();
381  }
382  if (fProofChain)
383  // This updates the proxy chain when we will really use PROOF
385 
386  return nf;
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Add a new file to this chain.
391 ///
392 /// Filename formats are similar to TChain::Add. Wildcards are not
393 /// applied. urls may also contain query and fragment identifiers
394 /// where the tree name can be specified in the url fragment.
395 ///
396 /// eg.
397 /// ~~~ {.cpp}
398 /// root://machine/path/file_name.root?query#tree_name
399 /// ~~~
400 /// If tree_name is given as a part of the file name it is used to
401 /// as the name of the tree to load from the file. Otherwise if tname
402 /// argument is specified the chain will load the tree named tname from
403 /// the file, otherwise the original treename specified in the TChain
404 /// constructor will be used.
405 ///
406 /// A. If nentries <= 0, the file is opened and the tree header read
407 /// into memory to get the number of entries.
408 ///
409 /// B. If nentries > 0, the file is not opened, and nentries is assumed
410 /// to be the number of entries in the file. In this case, no check
411 /// is made that the file exists nor that the tree exists in the file.
412 /// This second mode is interesting in case the number of entries in
413 /// the file is already stored in a run database for example.
414 ///
415 /// C. If nentries == TTree::kMaxEntries (default), the file is not opened.
416 /// The number of entries in each file will be read only when the file
417 /// is opened to read an entry. This option is the default and very
418 /// efficient if one processes the chain sequentially. Note that in
419 /// case GetEntry(entry) is called and entry refers to an entry in the
420 /// third file, for example, this forces the tree headers in the first
421 /// and second file to be read to find the number of entries in those
422 /// files. Note that if one calls GetEntriesFast() after having created
423 /// a chain with this default, GetEntriesFast() will return TTree::kMaxEntries!
424 /// Using the GetEntries() function instead will force all of the tree
425 /// headers in the chain to be read to read the number of entries in
426 /// each tree.
427 ///
428 /// D. The TChain data structure
429 /// Each TChainElement has a name equal to the tree name of this TChain
430 /// and a title equal to the file name. So, to loop over the
431 /// TFiles that have been added to this chain:
432 /// ~~~ {.cpp}
433 /// TObjArray *fileElements=chain->GetListOfFiles();
434 /// TIter next(fileElements);
435 /// TChainElement *chEl=0;
436 /// while (( chEl=(TChainElement*)next() )) {
437 /// TFile f(chEl->GetTitle());
438 /// ... do something with f ...
439 /// }
440 /// ~~~
441 /// The function returns 1 if the file is successfully connected, 0 otherwise.
442 
443 Int_t TChain::AddFile(const char* name, Long64_t nentries /* = TTree::kMaxEntries */, const char* tname /* = "" */)
444 {
445  if(name==0 || name[0]=='\0') {
446  Error("AddFile", "No file name; no files connected");
447  return 0;
448  }
449 
450  const char *treename = GetName();
451  if (tname && strlen(tname) > 0) treename = tname;
452 
453  TString basename, tn, query, suffix;
454  ParseTreeFilename(name, basename, tn, query, suffix, kFALSE);
455 
456  if (!tn.IsNull()) {
457  treename = tn.Data();
458  }
459 
460  Int_t nch = basename.Length() + query.Length();
461  char *filename = new char[nch+1];
462  strlcpy(filename,basename.Data(),nch+1);
463  strlcat(filename,query.Data(),nch+1);
464 
465  //Check enough space in fTreeOffset
466  if (fNtrees+1 >= fTreeOffsetLen) {
467  fTreeOffsetLen *= 2;
468  Long64_t *trees = new Long64_t[fTreeOffsetLen];
469  for (Int_t i=0;i<=fNtrees;i++) trees[i] = fTreeOffset[i];
470  delete [] fTreeOffset;
471  fTreeOffset = trees;
472  }
473 
474  // Open the file to get the number of entries.
475  Int_t pksize = 0;
476  if (nentries <= 0) {
477  TFile* file;
478  {
480  file = TFile::Open(filename);
481  }
482  if (!file || file->IsZombie()) {
483  delete file;
484  file = 0;
485  delete[] filename;
486  filename = 0;
487  return 0;
488  }
489 
490  // Check that tree with the right name exists in the file.
491  // Note: We are not the owner of obj, the file is!
492  TObject* obj = file->Get(treename);
493  if (!obj || !obj->InheritsFrom(TTree::Class())) {
494  Error("AddFile", "cannot find tree with name %s in file %s", treename, filename);
495  delete file;
496  file = 0;
497  delete[] filename;
498  filename = 0;
499  return 0;
500  }
501  TTree* tree = (TTree*) obj;
502  nentries = tree->GetEntries();
503  pksize = tree->GetPacketSize();
504  // Note: This deletes the tree we fetched.
505  delete file;
506  file = 0;
507  }
508 
509  if (nentries > 0) {
510  if (nentries != TTree::kMaxEntries) {
512  fEntries += nentries;
513  } else {
516  }
517  fNtrees++;
518 
519  TChainElement* element = new TChainElement(treename, filename);
520  element->SetPacketSize(pksize);
521  element->SetNumberEntries(nentries);
522  fFiles->Add(element);
523  } else {
524  Warning("AddFile", "Adding tree with no entries from file: %s", filename);
525  }
526 
527  delete [] filename;
528  if (fProofChain)
529  // This updates the proxy chain when we will really use PROOF
531 
532  return 1;
533 }
534 
535 ////////////////////////////////////////////////////////////////////////////////
536 /// Add all files referenced in the list to the chain. The object type in the
537 /// list must be either TFileInfo or TObjString or TUrl .
538 /// The function return 1 if successful, 0 otherwise.
539 
540 Int_t TChain::AddFileInfoList(TCollection* filelist, Long64_t nfiles /* = TTree::kMaxEntries */)
541 {
542  if (!filelist)
543  return 0;
544  TIter next(filelist);
545 
546  TObject *o = 0;
547  Long64_t cnt=0;
548  while ((o = next())) {
549  // Get the url
550  TString cn = o->ClassName();
551  const char *url = 0;
552  if (cn == "TFileInfo") {
553  TFileInfo *fi = (TFileInfo *)o;
554  url = (fi->GetCurrentUrl()) ? fi->GetCurrentUrl()->GetUrl() : 0;
555  if (!url) {
556  Warning("AddFileInfoList", "found TFileInfo with empty Url - ignoring");
557  continue;
558  }
559  } else if (cn == "TUrl") {
560  url = ((TUrl*)o)->GetUrl();
561  } else if (cn == "TObjString") {
562  url = ((TObjString*)o)->GetName();
563  }
564  if (!url) {
565  Warning("AddFileInfoList", "object is of type %s : expecting TFileInfo, TUrl"
566  " or TObjString - ignoring", o->ClassName());
567  continue;
568  }
569  // Good entry
570  cnt++;
571  AddFile(url);
572  if (cnt >= nfiles)
573  break;
574  }
575  if (fProofChain) {
576  // This updates the proxy chain when we will really use PROOF
578  }
579 
580  return 1;
581 }
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// Add a TFriendElement to the list of friends of this chain.
585 ///
586 /// A TChain has a list of friends similar to a tree (see TTree::AddFriend).
587 /// You can add a friend to a chain with the TChain::AddFriend method, and you
588 /// can retrieve the list of friends with TChain::GetListOfFriends.
589 /// This example has four chains each has 20 ROOT trees from 20 ROOT files.
590 /// ~~~ {.cpp}
591 /// TChain ch("t"); // a chain with 20 trees from 20 files
592 /// TChain ch1("t1");
593 /// TChain ch2("t2");
594 /// TChain ch3("t3");
595 /// ~~~
596 /// Now we can add the friends to the first chain.
597 /// ~~~ {.cpp}
598 /// ch.AddFriend("t1")
599 /// ch.AddFriend("t2")
600 /// ch.AddFriend("t3")
601 /// ~~~
602 /// \image html tchain_friend.png
603 ///
604 ///
605 /// The parameter is the name of friend chain (the name of a chain is always
606 /// the name of the tree from which it was created).
607 /// The original chain has access to all variable in its friends.
608 /// We can use the TChain::Draw method as if the values in the friends were
609 /// in the original chain.
610 /// To specify the chain to use in the Draw method, use the syntax:
611 /// ~~~ {.cpp}
612 /// <chainname>.<branchname>.<varname>
613 /// ~~~
614 /// If the variable name is enough to uniquely identify the variable, you can
615 /// leave out the chain and/or branch name.
616 /// For example, this generates a 3-d scatter plot of variable "var" in the
617 /// TChain ch versus variable v1 in TChain t1 versus variable v2 in TChain t2.
618 /// ~~~ {.cpp}
619 /// ch.Draw("var:t1.v1:t2.v2");
620 /// ~~~
621 /// When a TChain::Draw is executed, an automatic call to TTree::AddFriend
622 /// connects the trees in the chain. When a chain is deleted, its friend
623 /// elements are also deleted.
624 ///
625 /// The number of entries in the friend must be equal or greater to the number
626 /// of entries of the original chain. If the friend has fewer entries a warning
627 /// is given and the resulting histogram will have missing entries.
628 /// For additional information see TTree::AddFriend.
629 
630 TFriendElement* TChain::AddFriend(const char* chain, const char* dummy /* = "" */)
631 {
632  if (!fFriends) {
633  fFriends = new TList();
634  }
635  TFriendElement* fe = new TFriendElement(this, chain, dummy);
636 
637  R__ASSERT(fe); // There used to be a "if (fe)" test ... Keep this assert until we are sure that fe is never null
638 
639  fFriends->Add(fe);
640 
641  if (fProofChain)
642  // This updates the proxy chain when we will really use PROOF
644 
645  // We need to invalidate the loading of the current tree because its list
646  // of real friends is now obsolete. It is repairable only from LoadTree.
648 
649  TTree* tree = fe->GetTree();
650  if (!tree) {
651  Warning("AddFriend", "Unknown TChain %s", chain);
652  }
653  return fe;
654 }
655 
656 ////////////////////////////////////////////////////////////////////////////////
657 /// Add the whole chain or tree as a friend of this chain.
658 
660 {
661  if (!fFriends) fFriends = new TList();
662  TFriendElement *fe = new TFriendElement(this,chain,dummy);
663 
664  R__ASSERT(fe); // There used to be a "if (fe)" test ... Keep this assert until we are sure that fe is never null
665 
666  fFriends->Add(fe);
667 
668  if (fProofChain)
669  // This updates the proxy chain when we will really use PROOF
671 
672  // We need to invalidate the loading of the current tree because its list
673  // of real friend is now obsolete. It is repairable only from LoadTree
675 
676  TTree *t = fe->GetTree();
677  if (!t) {
678  Warning("AddFriend","Unknown TChain %s",chain);
679  }
680  return fe;
681 }
682 
683 ////////////////////////////////////////////////////////////////////////////////
684 /// Add the whole chain or tree as a friend of this chain.
685 
686 TFriendElement* TChain::AddFriend(TTree* chain, const char* alias, Bool_t /* warn = kFALSE */)
687 {
688  if (!fFriends) fFriends = new TList();
689  TFriendElement *fe = new TFriendElement(this,chain,alias);
690  R__ASSERT(fe);
691 
692  fFriends->Add(fe);
693 
694  if (fProofChain)
695  // This updates the proxy chain when we will really use PROOF
697 
698  // We need to invalidate the loading of the current tree because its list
699  // of real friend is now obsolete. It is repairable only from LoadTree
701 
702  TTree *t = fe->GetTree();
703  if (!t) {
704  Warning("AddFriend","Unknown TChain %s",chain->GetName());
705  }
706  return fe;
707 }
708 
709 ////////////////////////////////////////////////////////////////////////////////
710 /// Browse the contents of the chain.
711 
713 {
714  TTree::Browse(b);
715 }
716 
717 ////////////////////////////////////////////////////////////////////////////////
718 /// When closing a file during the chain processing, the file
719 /// may be closed with option "R" if flag is set to kTRUE.
720 /// by default flag is kTRUE.
721 /// When closing a file with option "R", all TProcessIDs referenced by this
722 /// file are deleted.
723 /// Calling TFile::Close("R") might be necessary in case one reads a long list
724 /// of files having TRef, writing some of the referenced objects or TRef
725 /// to a new file. If the TRef or referenced objects of the file being closed
726 /// will not be referenced again, it is possible to minimize the size
727 /// of the TProcessID data structures in memory by forcing a delete of
728 /// the unused TProcessID.
729 
730 void TChain::CanDeleteRefs(Bool_t flag /* = kTRUE */)
731 {
732  fCanDeleteRefs = flag;
733 }
734 
735 ////////////////////////////////////////////////////////////////////////////////
736 /// Initialize the packet descriptor string.
737 
739 {
740  TIter next(fFiles);
741  TChainElement* element = 0;
742  while ((element = (TChainElement*) next())) {
743  element->CreatePackets();
744  }
745 }
746 
747 ////////////////////////////////////////////////////////////////////////////////
748 /// Override the TTree::DirectoryAutoAdd behavior:
749 /// we never auto add.
750 
752 {
753 }
754 
755 ////////////////////////////////////////////////////////////////////////////////
756 /// Draw expression varexp for selected entries.
757 /// Returns -1 in case of error or number of selected events in case of success.
758 ///
759 /// This function accepts TCut objects as arguments.
760 /// Useful to use the string operator +, example:
761 /// ntuple.Draw("x",cut1+cut2+cut3);
762 ///
763 
764 Long64_t TChain::Draw(const char* varexp, const TCut& selection,
765  Option_t* option, Long64_t nentries, Long64_t firstentry)
766 {
767  if (fProofChain) {
768  // Make sure the element list is uptodate
769  if (!TestBit(kProofUptodate))
770  SetProof(kTRUE, kTRUE);
773  return fProofChain->Draw(varexp, selection, option, nentries, firstentry);
774  }
775 
776  return TChain::Draw(varexp, selection.GetTitle(), option, nentries, firstentry);
777 }
778 
779 ////////////////////////////////////////////////////////////////////////////////
780 /// Process all entries in this chain and draw histogram corresponding to
781 /// expression varexp.
782 /// Returns -1 in case of error or number of selected events in case of success.
783 
784 Long64_t TChain::Draw(const char* varexp, const char* selection,
785  Option_t* option,Long64_t nentries, Long64_t firstentry)
786 {
787  if (fProofChain) {
788  // Make sure the element list is uptodate
789  if (!TestBit(kProofUptodate))
790  SetProof(kTRUE, kTRUE);
793  return fProofChain->Draw(varexp, selection, option, nentries, firstentry);
794  }
795  GetPlayer();
796  if (LoadTree(firstentry) < 0) return 0;
797  return TTree::Draw(varexp,selection,option,nentries,firstentry);
798 }
799 
800 ////////////////////////////////////////////////////////////////////////////////
801 /// See TTree::GetReadEntry().
802 
803 TBranch* TChain::FindBranch(const char* branchname)
804 {
806  // Make sure the element list is uptodate
807  if (!TestBit(kProofUptodate))
808  SetProof(kTRUE, kTRUE);
809  return fProofChain->FindBranch(branchname);
810  }
811  if (fTree) {
812  return fTree->FindBranch(branchname);
813  }
814  LoadTree(0);
815  if (fTree) {
816  return fTree->FindBranch(branchname);
817  }
818  return 0;
819 }
820 
821 ////////////////////////////////////////////////////////////////////////////////
822 /// See TTree::GetReadEntry().
823 
824 TLeaf* TChain::FindLeaf(const char* searchname)
825 {
827  // Make sure the element list is uptodate
828  if (!TestBit(kProofUptodate))
829  SetProof(kTRUE, kTRUE);
830  return fProofChain->FindLeaf(searchname);
831  }
832  if (fTree) {
833  return fTree->FindLeaf(searchname);
834  }
835  LoadTree(0);
836  if (fTree) {
837  return fTree->FindLeaf(searchname);
838  }
839  return 0;
840 }
841 
842 ////////////////////////////////////////////////////////////////////////////////
843 /// Returns the expanded value of the alias. Search in the friends if any.
844 
845 const char* TChain::GetAlias(const char* aliasName) const
846 {
847  const char* alias = TTree::GetAlias(aliasName);
848  if (alias) {
849  return alias;
850  }
851  if (fTree) {
852  return fTree->GetAlias(aliasName);
853  }
854  const_cast<TChain*>(this)->LoadTree(0);
855  if (fTree) {
856  return fTree->GetAlias(aliasName);
857  }
858  return 0;
859 }
860 
861 ////////////////////////////////////////////////////////////////////////////////
862 /// Return pointer to the branch name in the current tree.
863 
865 {
867  // Make sure the element list is uptodate
868  if (!TestBit(kProofUptodate))
869  SetProof(kTRUE, kTRUE);
870  return fProofChain->GetBranch(name);
871  }
872  if (fTree) {
873  return fTree->GetBranch(name);
874  }
875  LoadTree(0);
876  if (fTree) {
877  return fTree->GetBranch(name);
878  }
879  return 0;
880 }
881 
882 ////////////////////////////////////////////////////////////////////////////////
883 /// See TTree::GetReadEntry().
884 
885 Bool_t TChain::GetBranchStatus(const char* branchname) const
886 {
888  // Make sure the element list is uptodate
889  if (!TestBit(kProofUptodate))
890  Warning("GetBranchStatus", "PROOF proxy not up-to-date:"
891  " run TChain::SetProof(kTRUE, kTRUE) first");
892  return fProofChain->GetBranchStatus(branchname);
893  }
894  return TTree::GetBranchStatus(branchname);
895 }
896 
897 ////////////////////////////////////////////////////////////////////////////////
898 /// Return an iterator over the cluster of baskets starting at firstentry.
899 ///
900 /// This iterator is not yet supported for TChain object.
901 
903 {
904  Fatal("GetClusterIterator","Not support for TChain object");
905  return TTree::GetClusterIterator(-1);
906 }
907 
908 ////////////////////////////////////////////////////////////////////////////////
909 /// Return absolute entry number in the chain.
910 /// The input parameter entry is the entry number in
911 /// the current tree of this chain.
912 
914 {
915  return entry + fTreeOffset[fTreeNumber];
916 }
917 
918 ////////////////////////////////////////////////////////////////////////////////
919 /// Return the total number of entries in the chain.
920 /// In case the number of entries in each tree is not yet known,
921 /// the offset table is computed.
922 
924 {
926  // Make sure the element list is uptodate
927  if (!TestBit(kProofUptodate))
928  Warning("GetEntries", "PROOF proxy not up-to-date:"
929  " run TChain::SetProof(kTRUE, kTRUE) first");
930  return fProofChain->GetEntries();
931  }
932  if (fEntries == TTree::kMaxEntries) {
933  const_cast<TChain*>(this)->LoadTree(TTree::kMaxEntries-1);
934  }
935  return fEntries;
936 }
937 
938 ////////////////////////////////////////////////////////////////////////////////
939 /// Get entry from the file to memory.
940 ///
941 /// - getall = 0 : get only active branches
942 /// - getall = 1 : get all branches
943 ///
944 /// Return the total number of bytes read,
945 /// 0 bytes read indicates a failure.
946 
948 {
949  Long64_t treeReadEntry = LoadTree(entry);
950  if (treeReadEntry < 0) {
951  return 0;
952  }
953  if (!fTree) {
954  return 0;
955  }
956  return fTree->GetEntry(treeReadEntry, getall);
957 }
958 
959 ////////////////////////////////////////////////////////////////////////////////
960 /// Return entry number corresponding to entry.
961 ///
962 /// if no TEntryList set returns entry
963 /// else returns entry \#entry from this entry list and
964 /// also computes the global entry number (loads all tree headers)
965 
967 {
968 
969  if (fEntryList){
970  Int_t treenum = 0;
971  Long64_t localentry = fEntryList->GetEntryAndTree(entry, treenum);
972  //find the global entry number
973  //same const_cast as in the GetEntries() function
974  if (localentry<0) return -1;
975  if (treenum != fTreeNumber){
976  if (fTreeOffset[treenum]==TTree::kMaxEntries){
977  for (Int_t i=0; i<=treenum; i++){
979  (const_cast<TChain*>(this))->LoadTree(fTreeOffset[i-1]);
980  }
981  }
982  //(const_cast<TChain*>(this))->LoadTree(fTreeOffset[treenum]);
983  }
984  Long64_t globalentry = fTreeOffset[treenum] + localentry;
985  return globalentry;
986  }
987  return entry;
988 }
989 
990 ////////////////////////////////////////////////////////////////////////////////
991 /// Return entry corresponding to major and minor number.
992 ///
993 /// The function returns the total number of bytes read.
994 /// If the Tree has friend trees, the corresponding entry with
995 /// the index values (major,minor) is read. Note that the master Tree
996 /// and its friend may have different entry serial numbers corresponding
997 /// to (major,minor).
998 
1000 {
1001  Long64_t serial = GetEntryNumberWithIndex(major, minor);
1002  if (serial < 0) return -1;
1003  return GetEntry(serial);
1004 }
1005 
1006 ////////////////////////////////////////////////////////////////////////////////
1007 /// Return a pointer to the current file.
1008 /// If no file is connected, the first file is automatically loaded.
1009 
1011 {
1012  if (fFile) {
1013  return fFile;
1014  }
1015  // Force opening the first file in the chain.
1016  const_cast<TChain*>(this)->LoadTree(0);
1017  return fFile;
1018 }
1019 
1020 ////////////////////////////////////////////////////////////////////////////////
1021 /// Return a pointer to the leaf name in the current tree.
1022 
1023 TLeaf* TChain::GetLeaf(const char* branchname, const char *leafname)
1024 {
1025  if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
1026  // Make sure the element list is uptodate
1027  if (!TestBit(kProofUptodate))
1028  SetProof(kTRUE, kTRUE);
1029  return fProofChain->GetLeaf(branchname, leafname);
1030  }
1031  if (fTree) {
1032  return fTree->GetLeaf(branchname, leafname);
1033  }
1034  LoadTree(0);
1035  if (fTree) {
1036  return fTree->GetLeaf(branchname, leafname);
1037  }
1038  return 0;
1039 }
1040 
1041 ////////////////////////////////////////////////////////////////////////////////
1042 /// Return a pointer to the leaf name in the current tree.
1043 
1045 {
1046  if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
1047  // Make sure the element list is uptodate
1048  if (!TestBit(kProofUptodate))
1049  SetProof(kTRUE, kTRUE);
1050  return fProofChain->GetLeaf(name);
1051  }
1052  if (fTree) {
1053  return fTree->GetLeaf(name);
1054  }
1055  LoadTree(0);
1056  if (fTree) {
1057  return fTree->GetLeaf(name);
1058  }
1059  return 0;
1060 }
1061 
1062 ////////////////////////////////////////////////////////////////////////////////
1063 /// Return a pointer to the list of branches of the current tree.
1064 ///
1065 /// Warning: If there is no current TTree yet, this routine will open the
1066 /// first in the chain.
1067 ///
1068 /// Returns 0 on failure.
1069 
1071 {
1072  if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
1073  // Make sure the element list is uptodate
1074  if (!TestBit(kProofUptodate))
1075  SetProof(kTRUE, kTRUE);
1076  return fProofChain->GetListOfBranches();
1077  }
1078  if (fTree) {
1079  return fTree->GetListOfBranches();
1080  }
1081  LoadTree(0);
1082  if (fTree) {
1083  return fTree->GetListOfBranches();
1084  }
1085  return 0;
1086 }
1087 
1088 ////////////////////////////////////////////////////////////////////////////////
1089 /// Return a pointer to the list of leaves of the current tree.
1090 ///
1091 /// Warning: May set the current tree!
1092 
1094 {
1095  if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
1096  // Make sure the element list is uptodate
1097  if (!TestBit(kProofUptodate))
1098  SetProof(kTRUE, kTRUE);
1099  return fProofChain->GetListOfLeaves();
1100  }
1101  if (fTree) {
1102  return fTree->GetListOfLeaves();
1103  }
1104  LoadTree(0);
1105  if (fTree) {
1106  return fTree->GetListOfLeaves();
1107  }
1108  return 0;
1109 }
1110 
1111 ////////////////////////////////////////////////////////////////////////////////
1112 /// Return maximum of column with name columname.
1113 
1114 Double_t TChain::GetMaximum(const char* columname)
1115 {
1116  Double_t theMax = -DBL_MAX;
1117  for (Int_t file = 0; file < fNtrees; file++) {
1119  LoadTree(first);
1120  Double_t curmax = fTree->GetMaximum(columname);
1121  if (curmax > theMax) {
1122  theMax = curmax;
1123  }
1124  }
1125  return theMax;
1126 }
1127 
1128 ////////////////////////////////////////////////////////////////////////////////
1129 /// Return minimum of column with name columname.
1130 
1131 Double_t TChain::GetMinimum(const char* columname)
1132 {
1133  Double_t theMin = DBL_MAX;
1134  for (Int_t file = 0; file < fNtrees; file++) {
1136  LoadTree(first);
1137  Double_t curmin = fTree->GetMinimum(columname);
1138  if (curmin < theMin) {
1139  theMin = curmin;
1140  }
1141  }
1142  return theMin;
1143 }
1144 
1145 ////////////////////////////////////////////////////////////////////////////////
1146 /// Return the number of branches of the current tree.
1147 ///
1148 /// Warning: May set the current tree!
1149 
1151 {
1152  if (fTree) {
1153  return fTree->GetNbranches();
1154  }
1155  LoadTree(0);
1156  if (fTree) {
1157  return fTree->GetNbranches();
1158  }
1159  return 0;
1160 }
1161 
1162 ////////////////////////////////////////////////////////////////////////////////
1163 /// See TTree::GetReadEntry().
1164 
1166 {
1167  if (fProofChain && !(fProofChain->TestBit(kProofLite))) {
1168  // Make sure the element list is uptodate
1169  if (!TestBit(kProofUptodate))
1170  Warning("GetBranchStatus", "PROOF proxy not up-to-date:"
1171  " run TChain::SetProof(kTRUE, kTRUE) first");
1172  return fProofChain->GetReadEntry();
1173  }
1174  return TTree::GetReadEntry();
1175 }
1176 
1177 ////////////////////////////////////////////////////////////////////////////////
1178 /// Return the chain weight.
1179 ///
1180 /// By default the weight is the weight of the current tree.
1181 /// However, if the weight has been set in TChain::SetWeight()
1182 /// with the option "global", then that weight will be returned.
1183 ///
1184 /// Warning: May set the current tree!
1185 
1187 {
1188  if (TestBit(kGlobalWeight)) {
1189  return fWeight;
1190  } else {
1191  if (fTree) {
1192  return fTree->GetWeight();
1193  }
1194  const_cast<TChain*>(this)->LoadTree(0);
1195  if (fTree) {
1196  return fTree->GetWeight();
1197  }
1198  return 0;
1199  }
1200 }
1201 
1202 ////////////////////////////////////////////////////////////////////////////////
1203 /// Set the TTree to be reloaded as soon as possible. In particular this
1204 /// is needed when adding a Friend.
1205 ///
1206 /// If the tree has clones, copy them into the chain
1207 /// clone list so we can change their branch addresses
1208 /// when necessary.
1209 ///
1210 /// This is to support the syntax:
1211 /// ~~~ {.cpp}
1212 /// TTree* clone = chain->GetTree()->CloneTree(0);
1213 /// ~~~
1214 
1216 {
1217  if (fTree && fTree->GetListOfClones()) {
1218  for (TObjLink* lnk = fTree->GetListOfClones()->FirstLink(); lnk; lnk = lnk->Next()) {
1219  TTree* clone = (TTree*) lnk->GetObject();
1220  AddClone(clone);
1221  }
1222  }
1223  fTreeNumber = -1;
1224  fTree = 0;
1225 }
1226 
1227 ////////////////////////////////////////////////////////////////////////////////
1228 /// Dummy function.
1229 /// It could be implemented and load all baskets of all trees in the chain.
1230 /// For the time being use TChain::Merge and TTree::LoadBasket
1231 /// on the resulting tree.
1232 
1234 {
1235  Error("LoadBaskets", "Function not yet implemented for TChain.");
1236  return 0;
1237 }
1238 
1239 ////////////////////////////////////////////////////////////////////////////////
1240 /// Find the tree which contains entry, and set it as the current tree.
1241 ///
1242 /// Returns the entry number in that tree.
1243 ///
1244 /// The input argument entry is the entry serial number in the whole chain.
1245 ///
1246 /// In case of error, LoadTree returns a negative number:
1247 /// * -1: The chain is empty.
1248 /// * -2: The requested entry number is less than zero or too large for the chain.
1249 /// or too large for the large TTree.
1250 /// * -3: The file corresponding to the entry could not be correctly open
1251 /// * -4: The TChainElement corresponding to the entry is missing or
1252 /// the TTree is missing from the file.
1253 /// * -5: Internal error, please report the circumstance when this happen
1254 /// as a ROOT issue.
1255 ///
1256 /// Note: This is the only routine which sets the value of fTree to
1257 /// a non-zero pointer.
1258 
1260 {
1261  // We already have been visited while recursively looking
1262  // through the friends tree, let's return.
1263  if (kLoadTree & fFriendLockStatus) {
1264  return 0;
1265  }
1266 
1267  if (!fNtrees) {
1268  // -- The chain is empty.
1269  return -1;
1270  }
1271 
1272  if ((entry < 0) || ((entry > 0) && (entry >= fEntries && entry!=(TTree::kMaxEntries-1) ))) {
1273  // -- Invalid entry number.
1274  if (fTree) fTree->LoadTree(-1);
1275  fReadEntry = -1;
1276  return -2;
1277  }
1278 
1279  // Find out which tree in the chain contains the passed entry.
1280  Int_t treenum = fTreeNumber;
1281  if ((fTreeNumber == -1) || (entry < fTreeOffset[fTreeNumber]) || (entry >= fTreeOffset[fTreeNumber+1]) || (entry==TTree::kMaxEntries-1)) {
1282  // -- Entry is *not* in the chain's current tree.
1283  // Do a linear search of the tree offset array.
1284  // FIXME: We could be smarter by starting at the
1285  // current tree number and going forwards,
1286  // then wrapping around at the end.
1287  for (treenum = 0; treenum < fNtrees; treenum++) {
1288  if (entry < fTreeOffset[treenum+1]) {
1289  break;
1290  }
1291  }
1292  }
1293 
1294  // Calculate the entry number relative to the found tree.
1295  Long64_t treeReadEntry = entry - fTreeOffset[treenum];
1296  fReadEntry = entry;
1297 
1298  // If entry belongs to the current tree return entry.
1299  if (fTree && treenum == fTreeNumber) {
1300  // First set the entry the tree on its owns friends
1301  // (the friends of the chain will be updated in the
1302  // next loop).
1303  fTree->LoadTree(treeReadEntry);
1304  if (fFriends) {
1305  // The current tree has not changed but some of its friends might.
1306  //
1307  // An alternative would move this code to each of
1308  // the functions calling LoadTree (and to overload a few more).
1309  TIter next(fFriends);
1310  TFriendLock lock(this, kLoadTree);
1311  TFriendElement* fe = 0;
1312  TFriendElement* fetree = 0;
1313  Bool_t needUpdate = kFALSE;
1314  while ((fe = (TFriendElement*) next())) {
1315  TObjLink* lnk = 0;
1316  if (fTree->GetListOfFriends()) {
1317  lnk = fTree->GetListOfFriends()->FirstLink();
1318  }
1319  fetree = 0;
1320  while (lnk) {
1321  TObject* obj = lnk->GetObject();
1322  if (obj->TestBit(TFriendElement::kFromChain) && obj->GetName() && !strcmp(fe->GetName(), obj->GetName())) {
1323  fetree = (TFriendElement*) obj;
1324  break;
1325  }
1326  lnk = lnk->Next();
1327  }
1328  TTree* at = fe->GetTree();
1329  if (at->InheritsFrom(TChain::Class())) {
1330  Int_t oldNumber = ((TChain*) at)->GetTreeNumber();
1331  TTree* old = at->GetTree();
1332  TTree* oldintree = fetree ? fetree->GetTree() : 0;
1333  at->LoadTreeFriend(entry, this);
1334  Int_t newNumber = ((TChain*) at)->GetTreeNumber();
1335  if ((oldNumber != newNumber) || (old != at->GetTree()) || (oldintree && (oldintree != at->GetTree()))) {
1336  // We can not compare just the tree pointers because
1337  // they could be reused. So we compare the tree
1338  // number instead.
1339  needUpdate = kTRUE;
1340  fTree->RemoveFriend(oldintree);
1342  }
1343  } else {
1344  // else we assume it is a simple tree If the tree is a
1345  // direct friend of the chain, it should be scanned
1346  // used the chain entry number and NOT the tree entry
1347  // number (treeReadEntry) hence we redo:
1348  at->LoadTreeFriend(entry, this);
1349  }
1350  }
1351  if (needUpdate) {
1352  // Update the branch/leaf addresses and
1353  // thelist of leaves in all TTreeFormula of the TTreePlayer (if any).
1354 
1355  // Set the branch statuses for the newly opened file.
1356  TChainElement *frelement;
1357  TIter fnext(fStatus);
1358  while ((frelement = (TChainElement*) fnext())) {
1359  Int_t status = frelement->GetStatus();
1360  fTree->SetBranchStatus(frelement->GetName(), status);
1361  }
1362 
1363  // Set the branch addresses for the newly opened file.
1364  fnext.Reset();
1365  while ((frelement = (TChainElement*) fnext())) {
1366  void* addr = frelement->GetBaddress();
1367  if (addr) {
1368  TBranch* br = fTree->GetBranch(frelement->GetName());
1369  TBranch** pp = frelement->GetBranchPtr();
1370  if (pp) {
1371  // FIXME: What if br is zero here?
1372  *pp = br;
1373  }
1374  if (br) {
1375  // FIXME: We may have to tell the branch it should
1376  // not be an owner of the object pointed at.
1377  br->SetAddress(addr);
1378  if (TestBit(kAutoDelete)) {
1379  br->SetAutoDelete(kTRUE);
1380  }
1381  }
1382  }
1383  }
1384  if (fPlayer) {
1386  }
1387  // Notify user if requested.
1388  if (fNotify) {
1389  fNotify->Notify();
1390  }
1391  }
1392  }
1393  return treeReadEntry;
1394  }
1395 
1396  // Delete the current tree and open the new tree.
1397 
1398  TTreeCache* tpf = 0;
1399  // Delete file unless the file owns this chain!
1400  // FIXME: The "unless" case here causes us to leak memory.
1401  if (fFile) {
1402  if (!fDirectory->GetList()->FindObject(this)) {
1403  if (fTree) {
1404  // (fFile != 0 && fTree == 0) can happen when
1405  // InvalidateCurrentTree is called (for example from
1406  // AddFriend). Having fTree === 0 is necessary in that
1407  // case because in some cases GetTree is used as a check
1408  // to see if a TTree is already loaded.
1409  // However, this prevent using the following to reuse
1410  // the TTreeCache object.
1411  tpf = (TTreeCache*) fFile->GetCacheRead(fTree);
1412  if (tpf) {
1413  tpf->ResetCache();
1414  }
1415 
1416  fFile->SetCacheRead(0, fTree);
1417  // If the tree has clones, copy them into the chain
1418  // clone list so we can change their branch addresses
1419  // when necessary.
1420  //
1421  // This is to support the syntax:
1422  //
1423  // TTree* clone = chain->GetTree()->CloneTree(0);
1424  //
1425  // We need to call the invalidate exactly here, since
1426  // we no longer need the value of fTree and it is
1427  // about to be deleted.
1429  }
1430 
1431  if (fCanDeleteRefs) {
1432  fFile->Close("R");
1433  }
1434  delete fFile;
1435  fFile = 0;
1436  } else {
1437  // If the tree has clones, copy them into the chain
1438  // clone list so we can change their branch addresses
1439  // when necessary.
1440  //
1441  // This is to support the syntax:
1442  //
1443  // TTree* clone = chain->GetTree()->CloneTree(0);
1444  //
1446  }
1447  }
1448 
1449  TChainElement* element = (TChainElement*) fFiles->At(treenum);
1450  if (!element) {
1451  if (treeReadEntry) {
1452  return -4;
1453  }
1454  // Last attempt, just in case all trees in the chain have 0 entries.
1455  element = (TChainElement*) fFiles->At(0);
1456  if (!element) {
1457  return -4;
1458  }
1459  }
1460 
1461  // FIXME: We leak memory here, we've just lost the open file
1462  // if we did not delete it above.
1463  {
1464  TDirectory::TContext ctxt;
1465  fFile = TFile::Open(element->GetTitle());
1466  if (fFile) fFile->SetBit(kMustCleanup);
1467  }
1468 
1469  // ----- Begin of modifications by MvL
1470  Int_t returnCode = 0;
1471  if (!fFile || fFile->IsZombie()) {
1472  if (fFile) {
1473  delete fFile;
1474  fFile = 0;
1475  }
1476  // Note: We do *not* own fTree.
1477  fTree = 0;
1478  returnCode = -3;
1479  } else {
1480  // Note: We do *not* own fTree after this, the file does!
1481  fTree = (TTree*) fFile->Get(element->GetName());
1482  if (!fTree) {
1483  // Now that we do not check during the addition, we need to check here!
1484  Error("LoadTree", "Cannot find tree with name %s in file %s", element->GetName(), element->GetTitle());
1485  delete fFile;
1486  fFile = 0;
1487  // We do not return yet so that 'fEntries' can be updated with the
1488  // sum of the entries of all the other trees.
1489  returnCode = -4;
1490  }
1491  }
1492 
1493  fTreeNumber = treenum;
1494  // FIXME: We own fFile, we must be careful giving away a pointer to it!
1495  // FIXME: We may set fDirectory to zero here!
1496  fDirectory = fFile;
1497 
1498  // Reuse cache from previous file (if any).
1499  if (tpf) {
1500  if (fFile) {
1501  tpf->ResetCache();
1502  fFile->SetCacheRead(tpf, fTree);
1503  // FIXME: fTree may be zero here.
1504  tpf->UpdateBranches(fTree);
1505  } else {
1506  // FIXME: One of the file in the chain is missing
1507  // we have no place to hold the pointer to the
1508  // TTreeCache.
1509  delete tpf;
1510  tpf = 0;
1511  }
1512  } else {
1513  if (fCacheUserSet) {
1514  this->SetCacheSize(fCacheSize);
1515  }
1516  }
1517 
1518  // Check if fTreeOffset has really been set.
1519  Long64_t nentries = 0;
1520  if (fTree) {
1521  nentries = fTree->GetEntries();
1522  }
1523 
1524  if (fTreeOffset[fTreeNumber+1] != (fTreeOffset[fTreeNumber] + nentries)) {
1525  fTreeOffset[fTreeNumber+1] = fTreeOffset[fTreeNumber] + nentries;
1526  fEntries = fTreeOffset[fNtrees];
1527  element->SetNumberEntries(nentries);
1528  // Below we must test >= in case the tree has no entries.
1529  if (entry >= fTreeOffset[fTreeNumber+1]) {
1530  if ((fTreeNumber < (fNtrees - 1)) && (entry < fTreeOffset[fTreeNumber+2])) {
1531  // The request entry is not in the tree 'fTreeNumber' we will need
1532  // to look further.
1533 
1534  // Before moving on, let's record the result.
1535  element->SetLoadResult(returnCode);
1536 
1537  // Before trying to read the file file/tree, notify the user
1538  // that we have switched trees if requested; the user might need
1539  // to properly account for the number of files/trees even if they
1540  // have no entries.
1541  if (fNotify) {
1542  fNotify->Notify();
1543  }
1544 
1545  // Load the next TTree.
1546  return LoadTree(entry);
1547  } else {
1548  treeReadEntry = fReadEntry = -2;
1549  }
1550  }
1551  }
1552 
1553 
1554  if (!fTree) {
1555  // The Error message already issued. However if we reach here
1556  // we need to make sure that we do not use fTree.
1557  //
1558  // Force a reload of the tree next time.
1559  fTreeNumber = -1;
1560 
1561  element->SetLoadResult(returnCode);
1562  return returnCode;
1563  }
1564  // ----- End of modifications by MvL
1565 
1566  // Copy the chain's clone list into the new tree's
1567  // clone list so that branch addresses stay synchronized.
1568  if (fClones) {
1569  for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
1570  TTree* clone = (TTree*) lnk->GetObject();
1571  ((TChain*) fTree)->TTree::AddClone(clone);
1572  }
1573  }
1574 
1575  // Since some of the friends of this chain might simple trees
1576  // (i.e., not really chains at all), we need to execute this
1577  // before calling LoadTree(entry) on the friends (so that
1578  // they use the correct read entry number).
1579 
1580  // Change the new current tree to the new entry.
1581  Long64_t loadResult = fTree->LoadTree(treeReadEntry);
1582  if (loadResult == treeReadEntry) {
1583  element->SetLoadResult(0);
1584  } else {
1585  // This is likely to be an internal error, if treeReadEntry was not in range
1586  // (or intentionally -2 for TChain::GetEntries) then something happened
1587  // that is very odd/surprising.
1588  element->SetLoadResult(-5);
1589  }
1590 
1591 
1592  // Change the chain friends to the new entry.
1593  if (fFriends) {
1594  // An alternative would move this code to each of the function
1595  // calling LoadTree (and to overload a few more).
1596  TIter next(fFriends);
1597  TFriendLock lock(this, kLoadTree);
1598  TFriendElement* fe = 0;
1599  while ((fe = (TFriendElement*) next())) {
1600  TTree* t = fe->GetTree();
1601  if (!t) continue;
1602  if (t->GetTreeIndex()) {
1604  }
1605  if (t->GetTree() && t->GetTree()->GetTreeIndex()) {
1607  }
1608  t->LoadTreeFriend(entry, this);
1609  TTree* friend_t = t->GetTree();
1610  if (friend_t) {
1612  }
1613  }
1614  }
1615 
1618 
1619  SetChainOffset(fTreeOffset[fTreeNumber]);
1620 
1621  // Set the branch statuses for the newly opened file.
1622  TIter next(fStatus);
1623  while ((element = (TChainElement*) next())) {
1624  Int_t status = element->GetStatus();
1625  fTree->SetBranchStatus(element->GetName(), status);
1626  }
1627 
1628  // Set the branch addresses for the newly opened file.
1629  next.Reset();
1630  while ((element = (TChainElement*) next())) {
1631  void* addr = element->GetBaddress();
1632  if (addr) {
1633  TBranch* br = fTree->GetBranch(element->GetName());
1634  TBranch** pp = element->GetBranchPtr();
1635  if (pp) {
1636  // FIXME: What if br is zero here?
1637  *pp = br;
1638  }
1639  if (br) {
1640  // FIXME: We may have to tell the branch it should
1641  // not be an owner of the object pointed at.
1642  br->SetAddress(addr);
1643  if (TestBit(kAutoDelete)) {
1644  br->SetAutoDelete(kTRUE);
1645  }
1646  }
1647  }
1648  }
1649 
1650  // Update the addresses of the chain's cloned trees, if any.
1651  if (fClones) {
1652  for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
1653  TTree* clone = (TTree*) lnk->GetObject();
1654  CopyAddresses(clone);
1655  }
1656  }
1657 
1658  // Update list of leaves in all TTreeFormula's of the TTreePlayer (if any).
1659  if (fPlayer) {
1661  }
1662 
1663  // Notify user we have switched trees if requested.
1664  if (fNotify) {
1665  fNotify->Notify();
1666  }
1667 
1668  // Return the new local entry number.
1669  return treeReadEntry;
1670 }
1671 
1672 ////////////////////////////////////////////////////////////////////////////////
1673 /// Check / locate the files in the chain.
1674 /// By default only the files not yet looked up are checked.
1675 /// Use force = kTRUE to check / re-check every file.
1676 
1678 {
1679  TIter next(fFiles);
1680  TChainElement* element = 0;
1681  Int_t nelements = fFiles->GetEntries();
1682  printf("\n");
1683  printf("TChain::Lookup - Looking up %d files .... \n", nelements);
1684  Int_t nlook = 0;
1685  TFileStager *stg = 0;
1686  while ((element = (TChainElement*) next())) {
1687  // Do not do it more than needed
1688  if (element->HasBeenLookedUp() && !force) continue;
1689  // Count
1690  nlook++;
1691  // Get the Url
1692  TUrl elemurl(element->GetTitle(), kTRUE);
1693  // Save current options and anchor
1694  TString anchor = elemurl.GetAnchor();
1695  TString options = elemurl.GetOptions();
1696  // Reset options and anchor
1697  elemurl.SetOptions("");
1698  elemurl.SetAnchor("");
1699  // Locate the file
1700  TString eurl(elemurl.GetUrl());
1701  if (!stg || !stg->Matches(eurl)) {
1702  SafeDelete(stg);
1703  {
1704  TDirectory::TContext ctxt;
1705  stg = TFileStager::Open(eurl);
1706  }
1707  if (!stg) {
1708  Error("Lookup", "TFileStager instance cannot be instantiated");
1709  break;
1710  }
1711  }
1712  Int_t n1 = (nelements > 100) ? (Int_t) nelements / 100 : 1;
1713  if (stg->Locate(eurl.Data(), eurl) == 0) {
1714  if (nlook > 0 && !(nlook % n1)) {
1715  printf("Lookup | %3d %% finished\r", 100 * nlook / nelements);
1716  fflush(stdout);
1717  }
1718  // Get the effective end-point Url
1719  elemurl.SetUrl(eurl);
1720  // Restore original options and anchor, if any
1721  elemurl.SetOptions(options);
1722  elemurl.SetAnchor(anchor);
1723  // Save it into the element
1724  element->SetTitle(elemurl.GetUrl());
1725  // Remember
1726  element->SetLookedUp();
1727  } else {
1728  // Failure: remove
1729  fFiles->Remove(element);
1730  if (gSystem->AccessPathName(eurl))
1731  Error("Lookup", "file %s does not exist\n", eurl.Data());
1732  else
1733  Error("Lookup", "file %s cannot be read\n", eurl.Data());
1734  }
1735  }
1736  if (nelements > 0)
1737  printf("Lookup | %3d %% finished\n", 100 * nlook / nelements);
1738  else
1739  printf("\n");
1740  fflush(stdout);
1741  SafeDelete(stg);
1742 }
1743 
1744 ////////////////////////////////////////////////////////////////////////////////
1745 /// Loop on nentries of this chain starting at firstentry. (NOT IMPLEMENTED)
1746 
1747 void TChain::Loop(Option_t* option, Long64_t nentries, Long64_t firstentry)
1748 {
1749  Error("Loop", "Function not yet implemented");
1750 
1751  if (option || nentries || firstentry) { } // keep warnings away
1752 
1753 #if 0
1754  if (LoadTree(firstentry) < 0) return;
1755 
1756  if (firstentry < 0) firstentry = 0;
1757  Long64_t lastentry = firstentry + nentries -1;
1758  if (lastentry > fEntries-1) {
1759  lastentry = fEntries -1;
1760  }
1761 
1762  GetPlayer();
1763  GetSelector();
1764  fSelector->Start(option);
1765 
1766  Long64_t entry = firstentry;
1767  Int_t tree,e0,en;
1768  for (tree=0;tree<fNtrees;tree++) {
1769  e0 = fTreeOffset[tree];
1770  en = fTreeOffset[tree+1] - 1;
1771  if (en > lastentry) en = lastentry;
1772  if (entry > en) continue;
1773 
1774  LoadTree(entry);
1775  fSelector->BeginFile();
1776 
1777  while (entry <= en) {
1778  fSelector->Execute(fTree, entry - e0);
1779  entry++;
1780  }
1781  fSelector->EndFile();
1782  }
1783 
1784  fSelector->Finish(option);
1785 #endif
1786 }
1787 
1788 ////////////////////////////////////////////////////////////////////////////////
1789 /// List the chain.
1790 
1791 void TChain::ls(Option_t* option) const
1792 {
1793  TObject::ls(option);
1794  TIter next(fFiles);
1795  TChainElement* file = 0;
1797  while ((file = (TChainElement*)next())) {
1798  file->ls(option);
1799  }
1801 }
1802 
1803 ////////////////////////////////////////////////////////////////////////////////
1804 /// Merge all the entries in the chain into a new tree in a new file.
1805 ///
1806 /// See important note in the following function Merge().
1807 ///
1808 /// If the chain is expecting the input tree inside a directory,
1809 /// this directory is NOT created by this routine.
1810 ///
1811 /// So in a case where we have:
1812 /// ~~~ {.cpp}
1813 /// TChain ch("mydir/mytree");
1814 /// ch.Merge("newfile.root");
1815 /// ~~~
1816 /// The resulting file will have not subdirectory. To recreate
1817 /// the directory structure do:
1818 /// ~~~ {.cpp}
1819 /// TFile* file = TFile::Open("newfile.root", "RECREATE");
1820 /// file->mkdir("mydir")->cd();
1821 /// ch.Merge(file);
1822 /// ~~~
1823 
1824 Long64_t TChain::Merge(const char* name, Option_t* option)
1825 {
1826  TFile *file = TFile::Open(name, "recreate", "chain files", 1);
1827  return Merge(file, 0, option);
1828 }
1829 
1830 ////////////////////////////////////////////////////////////////////////////////
1831 /// Merge all chains in the collection. (NOT IMPLEMENTED)
1832 
1833 Long64_t TChain::Merge(TCollection* /* list */, Option_t* /* option */ )
1834 {
1835  Error("Merge", "not implemented");
1836  return -1;
1837 }
1838 
1839 ////////////////////////////////////////////////////////////////////////////////
1840 /// Merge all chains in the collection. (NOT IMPLEMENTED)
1841 
1843 {
1844  Error("Merge", "not implemented");
1845  return -1;
1846 }
1847 
1848 ////////////////////////////////////////////////////////////////////////////////
1849 /// Merge all the entries in the chain into a new tree in the current file.
1850 ///
1851 /// Note: The "file" parameter is *not* the file where the new
1852 /// tree will be inserted. The new tree is inserted into
1853 /// gDirectory, which is usually the most recently opened
1854 /// file, or the directory most recently cd()'d to.
1855 ///
1856 /// If option = "C" is given, the compression level for all branches
1857 /// in the new Tree is set to the file compression level. By default,
1858 /// the compression level of all branches is the original compression
1859 /// level in the old trees.
1860 ///
1861 /// If basketsize > 1000, the basket size for all branches of the
1862 /// new tree will be set to basketsize.
1863 ///
1864 /// Example using the file generated in $ROOTSYS/test/Event
1865 /// merge two copies of Event.root
1866 /// ~~~ {.cpp}
1867 /// gSystem.Load("libEvent");
1868 /// TChain ch("T");
1869 /// ch.Add("Event1.root");
1870 /// ch.Add("Event2.root");
1871 /// ch.Merge("all.root");
1872 /// ~~~
1873 /// If the chain is expecting the input tree inside a directory,
1874 /// this directory is NOT created by this routine.
1875 ///
1876 /// So if you do:
1877 /// ~~~ {.cpp}
1878 /// TChain ch("mydir/mytree");
1879 /// ch.Merge("newfile.root");
1880 /// ~~~
1881 /// The resulting file will not have subdirectories. In order to
1882 /// preserve the directory structure do the following instead:
1883 /// ~~~ {.cpp}
1884 /// TFile* file = TFile::Open("newfile.root", "RECREATE");
1885 /// file->mkdir("mydir")->cd();
1886 /// ch.Merge(file);
1887 /// ~~~
1888 /// If 'option' contains the word 'fast' the merge will be done without
1889 /// unzipping or unstreaming the baskets (i.e., a direct copy of the raw
1890 /// bytes on disk).
1891 ///
1892 /// When 'fast' is specified, 'option' can also contains a
1893 /// sorting order for the baskets in the output file.
1894 ///
1895 /// There is currently 3 supported sorting order:
1896 /// ~~~ {.cpp}
1897 /// SortBasketsByOffset (the default)
1898 /// SortBasketsByBranch
1899 /// SortBasketsByEntry
1900 /// ~~~
1901 /// When using SortBasketsByOffset the baskets are written in
1902 /// the output file in the same order as in the original file
1903 /// (i.e. the basket are sorted on their offset in the original
1904 /// file; Usually this also means that the baskets are sorted
1905 /// on the index/number of the _last_ entry they contain)
1906 ///
1907 /// When using SortBasketsByBranch all the baskets of each
1908 /// individual branches are stored contiguously. This tends to
1909 /// optimize reading speed when reading a small number (1->5) of
1910 /// branches, since all their baskets will be clustered together
1911 /// instead of being spread across the file. However it might
1912 /// decrease the performance when reading more branches (or the full
1913 /// entry).
1914 ///
1915 /// When using SortBasketsByEntry the baskets with the lowest
1916 /// starting entry are written first. (i.e. the baskets are
1917 /// sorted on the index/number of the first entry they contain).
1918 /// This means that on the file the baskets will be in the order
1919 /// in which they will be needed when reading the whole tree
1920 /// sequentially.
1921 ///
1922 /// ## IMPORTANT Note 1: AUTOMATIC FILE OVERFLOW
1923 ///
1924 /// When merging many files, it may happen that the resulting file
1925 /// reaches a size > TTree::fgMaxTreeSize (default = 100 GBytes).
1926 /// In this case the current file is automatically closed and a new
1927 /// file started. If the name of the merged file was "merged.root",
1928 /// the subsequent files will be named "merged_1.root", "merged_2.root",
1929 /// etc. fgMaxTreeSize may be modified via the static function
1930 /// TTree::SetMaxTreeSize.
1931 /// When in fast mode, the check and switch is only done in between each
1932 /// input file.
1933 ///
1934 /// ## IMPORTANT Note 2: The output file is automatically closed and deleted.
1935 ///
1936 /// This is required because in general the automatic file overflow described
1937 /// above may happen during the merge.
1938 /// If only the current file is produced (the file passed as first argument),
1939 /// one can instruct Merge to not close and delete the file by specifying
1940 /// the option "keep".
1941 ///
1942 /// The function returns the total number of files produced.
1943 /// To check that all files have been merged use something like:
1944 /// ~~~ {.cpp}
1945 /// if (newchain->GetEntries()!=oldchain->GetEntries()) {
1946 /// ... not all the file have been copied ...
1947 /// }
1948 /// ~~~
1949 
1951 {
1952  // We must have been passed a file, we will use it
1953  // later to reset the compression level of the branches.
1954  if (!file) {
1955  // FIXME: We need an error message here.
1956  return 0;
1957  }
1958 
1959  // Options
1960  Bool_t fastClone = kFALSE;
1961  TString opt = option;
1962  opt.ToLower();
1963  if (opt.Contains("fast")) {
1964  fastClone = kTRUE;
1965  }
1966 
1967  // The chain tree must have a list of branches
1968  // because we may try to change their basket
1969  // size later.
1970  TObjArray* lbranches = GetListOfBranches();
1971  if (!lbranches) {
1972  // FIXME: We need an error message here.
1973  return 0;
1974  }
1975 
1976  // The chain must have a current tree because
1977  // that is the one we will clone.
1978  if (!fTree) {
1979  // -- LoadTree() has not yet been called, no current tree.
1980  // FIXME: We need an error message here.
1981  return 0;
1982  }
1983 
1984  // Copy the chain's current tree without
1985  // copying any entries, we will do that later.
1986  TTree* newTree = CloneTree(0);
1987  if (!newTree) {
1988  // FIXME: We need an error message here.
1989  return 0;
1990  }
1991 
1992  // Strip out the (potential) directory name.
1993  // FIXME: The merged chain may or may not have the
1994  // same name as the original chain. This is
1995  // bad because the chain name determines the
1996  // names of the trees in the chain by default.
1997  newTree->SetName(gSystem->BaseName(GetName()));
1998 
1999  // FIXME: Why do we do this?
2000  newTree->SetAutoSave(2000000000);
2001 
2002  // Circularity is incompatible with merging, it may
2003  // force us to throw away entries, which is not what
2004  // we are supposed to do.
2005  newTree->SetCircular(0);
2006 
2007  // Reset the compression level of the branches.
2008  if (opt.Contains("c")) {
2009  TBranch* branch = 0;
2010  TIter nextb(newTree->GetListOfBranches());
2011  while ((branch = (TBranch*) nextb())) {
2013  }
2014  }
2015 
2016  // Reset the basket size of the branches.
2017  if (basketsize > 1000) {
2018  TBranch* branch = 0;
2019  TIter nextb(newTree->GetListOfBranches());
2020  while ((branch = (TBranch*) nextb())) {
2021  branch->SetBasketSize(basketsize);
2022  }
2023  }
2024 
2025  // Copy the entries.
2026  if (fastClone) {
2027  if ( newTree->CopyEntries( this, -1, option ) < 0 ) {
2028  // There was a problem!
2029  Error("Merge", "TTree has not been cloned\n");
2030  }
2031  } else {
2032  newTree->CopyEntries( this, -1, option );
2033  }
2034 
2035  // Write the new tree header.
2036  newTree->Write();
2037 
2038  // Get our return value.
2039  Int_t nfiles = newTree->GetFileNumber() + 1;
2040 
2041  // Close and delete the current file of the new tree.
2042  if (!opt.Contains("keep")) {
2043  // Delete the currentFile and the TTree object.
2044  delete newTree->GetCurrentFile();
2045  }
2046  return nfiles;
2047 }
2048 
2049 ////////////////////////////////////////////////////////////////////////////////
2050 /// Get the tree url or filename and other information from the name
2051 ///
2052 /// A treename and a url's query section is split off from name. The
2053 /// splitting depends on whether the resulting filename is to be
2054 /// subsequently treated for wildcards or not, since the question mark is
2055 /// both the url query identifier and a wildcard. Wildcard matching is not
2056 /// done in this method itself.
2057 /// ~~~ {.cpp}
2058 /// [xxx://host]/a/path/file.root[/treename]
2059 /// [xxx://host]/a/path/file.root[/treename][?query]
2060 /// [xxx://host]/a/path/file.root[?query[#treename]]
2061 /// ~~~
2062 ///
2063 /// Note that in a case like this
2064 /// ~~~ {.cpp}
2065 /// [xxx://host]/a/path/file#treename
2066 /// ~~~
2067 /// i.e. anchor but no options (query), the filename will be the full path, as the
2068 /// ancho may be the internal file name of an archive.
2069 ///
2070 /// \param[in] name is the original name
2071 /// \param[in] wildcards indicates if the resulting filename will be treated for
2072 /// wildcards. For backwards compatibility, with most protocols
2073 /// this flag suppresses the search for the url fragment
2074 /// identifier and limits the query identifier search to cases
2075 /// where the tree name is given as a trailing slash-separated
2076 /// string at the end of the file name.
2077 /// \param[out] filename the url or filename to be opened or matched
2078 /// \param[out] treename the treename, which may be found as a trailing part of the
2079 /// name or in a url fragment section. If not found this will
2080 /// be empty.
2081 /// \param[out] query is the url query section, including the leading question
2082 /// mark. If not found or the query section is only followed by
2083 /// a fragment this will be empty.
2084 /// \param[out] suffix the portion of name which was removed to form filename.
2085 
2086 void TChain::ParseTreeFilename(const char *name, TString &filename, TString &treename, TString &query, TString &suffix,
2087  Bool_t) const
2088 {
2089  Ssiz_t pIdx = kNPOS;
2090  filename = name;
2091  treename.Clear();
2092  query.Clear();
2093  suffix.Clear();
2094 
2095  // General case
2096  TUrl url(name, kTRUE);
2097 
2098  TString fn = url.GetFile();
2099  // Extract query (and treename, if any)
2100  if (url.GetOptions() && (strlen(url.GetOptions()) > 0)) {
2101  query.Form("?%s", url.GetOptions());
2102  // The treename can be passed as anchor in this case
2103  treename = url.GetAnchor();
2104  // Suffix
2105  suffix = url.GetFileAndOptions();
2106  suffix.Replace(suffix.Index(fn), fn.Length(), "");
2107  // Remove it from the file name
2108  filename.Remove(filename.Index(fn) + fn.Length());
2109  }
2110 
2111  // Special case: [...]file.root/treename
2112  static const char *dotr = ".root/";
2113  static Ssiz_t dotrl = strlen(dotr);
2114  // Find the last one
2115  Ssiz_t js = filename.Index(dotr);
2116  while (js != kNPOS) {
2117  pIdx = js;
2118  js = filename.Index(dotr, js + 1);
2119  }
2120  if (pIdx != kNPOS) {
2121  TString tn = filename(pIdx + dotrl, filename.Length());
2122  if (!tn.EndsWith(".root")) {
2123  // Good treename
2124  treename = tn;
2125  filename.Remove(pIdx + dotrl - 1);
2126  suffix.Insert(0, TString::Format("/%s", treename.Data()));
2127  }
2128  }
2129 }
2130 
2131 ////////////////////////////////////////////////////////////////////////////////
2132 /// Print the header information of each tree in the chain.
2133 /// See TTree::Print for a list of options.
2134 
2135 void TChain::Print(Option_t *option) const
2136 {
2137  TIter next(fFiles);
2138  TChainElement *element;
2139  while ((element = (TChainElement*)next())) {
2140  Printf("******************************************************************************");
2141  Printf("*Chain :%-10s: %-54s *", GetName(), element->GetTitle());
2142  Printf("******************************************************************************");
2143  TFile *file = TFile::Open(element->GetTitle());
2144  if (file && !file->IsZombie()) {
2145  TTree *tree = (TTree*)file->Get(element->GetName());
2146  if (tree) tree->Print(option);
2147  }
2148  delete file;
2149  }
2150 }
2151 
2152 ////////////////////////////////////////////////////////////////////////////////
2153 /// Process all entries in this chain, calling functions in filename.
2154 /// The return value is -1 in case of error and TSelector::GetStatus() in
2155 /// in case of success.
2156 /// See TTree::Process.
2157 
2158 Long64_t TChain::Process(const char *filename, Option_t *option, Long64_t nentries, Long64_t firstentry)
2159 {
2160  if (fProofChain) {
2161  // Make sure the element list is uptodate
2162  if (!TestBit(kProofUptodate))
2163  SetProof(kTRUE, kTRUE);
2166  return fProofChain->Process(filename, option, nentries, firstentry);
2167  }
2168 
2169  if (LoadTree(firstentry) < 0) {
2170  return 0;
2171  }
2172  return TTree::Process(filename, option, nentries, firstentry);
2173 }
2174 
2175 ////////////////////////////////////////////////////////////////////////////////
2176 /// Process this chain executing the code in selector.
2177 /// The return value is -1 in case of error and TSelector::GetStatus() in
2178 /// in case of success.
2179 
2181 {
2182  if (fProofChain) {
2183  // Make sure the element list is uptodate
2184  if (!TestBit(kProofUptodate))
2185  SetProof(kTRUE, kTRUE);
2188  return fProofChain->Process(selector, option, nentries, firstentry);
2189  }
2190 
2191  return TTree::Process(selector, option, nentries, firstentry);
2192 }
2193 
2194 ////////////////////////////////////////////////////////////////////////////////
2195 /// Make sure that obj (which is being deleted or will soon be) is no
2196 /// longer referenced by this TTree.
2197 
2199 {
2200  if (fFile == obj) {
2201  fFile = 0;
2202  fDirectory = 0;
2203  fTree = 0;
2204  }
2205  if (fDirectory == obj) {
2206  fDirectory = 0;
2207  fTree = 0;
2208  }
2209  if (fTree == obj) {
2210  fTree = 0;
2211  }
2212 }
2213 
2214 ////////////////////////////////////////////////////////////////////////////////
2215 /// Remove a friend from the list of friends.
2216 
2217 void TChain::RemoveFriend(TTree* oldFriend)
2218 {
2219  // We already have been visited while recursively looking
2220  // through the friends tree, let return
2221 
2222  if (!fFriends) {
2223  return;
2224  }
2225 
2226  TTree::RemoveFriend(oldFriend);
2227 
2228  if (fProofChain)
2229  // This updates the proxy chain when we will really use PROOF
2231 
2232  // We need to invalidate the loading of the current tree because its list
2233  // of real friends is now obsolete. It is repairable only from LoadTree.
2235 }
2236 
2237 ////////////////////////////////////////////////////////////////////////////////
2238 /// Resets the state of this chain.
2239 
2241 {
2242  delete fFile;
2243  fFile = 0;
2244  fNtrees = 0;
2245  fTreeNumber = -1;
2246  fTree = 0;
2247  fFile = 0;
2248  fFiles->Delete();
2249  fStatus->Delete();
2250  fTreeOffset[0] = 0;
2251  TChainElement* element = new TChainElement("*", "");
2252  fStatus->Add(element);
2253  fDirectory = 0;
2254 
2255  TTree::Reset();
2256 }
2257 
2258 ////////////////////////////////////////////////////////////////////////////////
2259 /// Resets the state of this chain after a merge (keep the customization but
2260 /// forget the data).
2261 
2263 {
2264  fNtrees = 0;
2265  fTreeNumber = -1;
2266  fTree = 0;
2267  fFile = 0;
2268  fFiles->Delete();
2269  fTreeOffset[0] = 0;
2270 
2271  TTree::ResetAfterMerge(info);
2272 }
2273 
2274 ////////////////////////////////////////////////////////////////////////////////
2275 /// Loop on tree and print entries passing selection.
2276 /// - If varexp is 0 (or "") then print only first 8 columns.
2277 /// - If varexp = "*" print all columns.
2278 /// - Otherwise a columns selection can be made using "var1:var2:var3".
2279 /// See TTreePlayer::Scan for more information.
2280 
2281 Long64_t TChain::Scan(const char* varexp, const char* selection, Option_t* option, Long64_t nentries, Long64_t firstentry)
2282 {
2283  if (LoadTree(firstentry) < 0) {
2284  return 0;
2285  }
2286  return TTree::Scan(varexp, selection, option, nentries, firstentry);
2287 }
2288 
2289 ////////////////////////////////////////////////////////////////////////////////
2290 /// Set the global branch kAutoDelete bit.
2291 ///
2292 /// When LoadTree loads a new Tree, the branches for which
2293 /// the address is set will have the option AutoDelete set
2294 /// For more details on AutoDelete, see TBranch::SetAutoDelete.
2295 
2297 {
2298  if (autodelete) {
2299  SetBit(kAutoDelete, 1);
2300  } else {
2301  SetBit(kAutoDelete, 0);
2302  }
2303 }
2304 
2306 {
2307  // Set the cache size of the underlying TTree,
2308  // See TTree::SetCacheSize.
2309  // Returns 0 cache state ok (exists or not, as appropriate)
2310  // -1 on error
2311 
2312  Int_t res = 0;
2313 
2314  // remember user has requested this cache setting
2315  fCacheUserSet = kTRUE;
2316 
2317  if (fTree) {
2318  res = fTree->SetCacheSize(cacheSize);
2319  } else {
2320  // If we don't have a TTree yet only record the cache size wanted
2321  res = 0;
2322  }
2323  fCacheSize = cacheSize; // Record requested size.
2324  return res;
2325 }
2326 
2327 ////////////////////////////////////////////////////////////////////////////////
2328 /// Reset the addresses of the branch.
2329 
2331 {
2332  TChainElement* element = (TChainElement*) fStatus->FindObject(branch->GetName());
2333  if (element) {
2334  element->SetBaddress(0);
2335  }
2336  if (fTree) {
2337  fTree->ResetBranchAddress(branch);
2338  }
2339 }
2340 
2341 ////////////////////////////////////////////////////////////////////////////////
2342 /// Reset the addresses of the branches.
2343 
2345 {
2346  TIter next(fStatus);
2347  TChainElement* element = 0;
2348  while ((element = (TChainElement*) next())) {
2349  element->SetBaddress(0);
2350  }
2351  if (fTree) {
2353  }
2354 }
2355 
2356 ////////////////////////////////////////////////////////////////////////////////
2357 /// Set branch address.
2358 ///
2359 /// \param[in] bname is the name of a branch.
2360 /// \param[in] add is the address of the branch.
2361 /// \param[in] ptr
2362 ///
2363 /// Note: See the comments in TBranchElement::SetAddress() for a more
2364 /// detailed discussion of the meaning of the add parameter.
2365 ///
2366 /// IMPORTANT REMARK:
2367 ///
2368 /// In case TChain::SetBranchStatus is called, it must be called
2369 /// BEFORE calling this function.
2370 ///
2371 /// See TTree::CheckBranchAddressType for the semantic of the return value.
2372 
2373 Int_t TChain::SetBranchAddress(const char *bname, void* add, TBranch** ptr)
2374 {
2375  Int_t res = kNoCheck;
2376 
2377  // Check if bname is already in the status list.
2378  // If not, create a TChainElement object and set its address.
2379  TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
2380  if (!element) {
2381  element = new TChainElement(bname, "");
2382  fStatus->Add(element);
2383  }
2384  element->SetBaddress(add);
2385  element->SetBranchPtr(ptr);
2386  // Also set address in current tree.
2387  // FIXME: What about the chain clones?
2388  if (fTreeNumber >= 0) {
2389  TBranch* branch = fTree->GetBranch(bname);
2390  if (ptr) {
2391  *ptr = branch;
2392  }
2393  if (branch) {
2394  res = CheckBranchAddressType(branch, TClass::GetClass(element->GetBaddressClassName()), (EDataType) element->GetBaddressType(), element->GetBaddressIsPtr());
2395  if (fClones) {
2396  void* oldAdd = branch->GetAddress();
2397  for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
2398  TTree* clone = (TTree*) lnk->GetObject();
2399  TBranch* cloneBr = clone->GetBranch(bname);
2400  if (cloneBr && (cloneBr->GetAddress() == oldAdd)) {
2401  // the clone's branch is still pointing to us
2402  cloneBr->SetAddress(add);
2403  }
2404  }
2405  }
2406  branch->SetAddress(add);
2407  } else {
2408  Error("SetBranchAddress", "unknown branch -> %s", bname);
2409  return kMissingBranch;
2410  }
2411  } else {
2412  if (ptr) {
2413  *ptr = 0;
2414  }
2415  }
2416  return res;
2417 }
2418 
2419 ////////////////////////////////////////////////////////////////////////////////
2420 /// Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
2421 /// See TTree::CheckBranchAddressType for the semantic of the return value.
2422 ///
2423 /// Note: See the comments in TBranchElement::SetAddress() for a more
2424 /// detailed discussion of the meaning of the add parameter.
2425 
2426 Int_t TChain::SetBranchAddress(const char* bname, void* add, TClass* realClass, EDataType datatype, Bool_t isptr)
2427 {
2428  return SetBranchAddress(bname, add, 0, realClass, datatype, isptr);
2429 }
2430 
2431 ////////////////////////////////////////////////////////////////////////////////
2432 /// Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
2433 /// See TTree::CheckBranchAddressType for the semantic of the return value.
2434 ///
2435 /// Note: See the comments in TBranchElement::SetAddress() for a more
2436 /// detailed discussion of the meaning of the add parameter.
2437 
2438 Int_t TChain::SetBranchAddress(const char* bname, void* add, TBranch** ptr, TClass* realClass, EDataType datatype, Bool_t isptr)
2439 {
2440  TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
2441  if (!element) {
2442  element = new TChainElement(bname, "");
2443  fStatus->Add(element);
2444  }
2445  if (realClass) {
2446  element->SetBaddressClassName(realClass->GetName());
2447  }
2448  element->SetBaddressType((UInt_t) datatype);
2449  element->SetBaddressIsPtr(isptr);
2450  element->SetBranchPtr(ptr);
2451  return SetBranchAddress(bname, add, ptr);
2452 }
2453 
2454 ////////////////////////////////////////////////////////////////////////////////
2455 /// Set branch status to Process or DoNotProcess
2456 ///
2457 /// \param[in] bname is the name of a branch. if bname="*", apply to all branches.
2458 /// \param[in] status = 1 branch will be processed,
2459 /// = 0 branch will not be processed
2460 /// \param[out] found
2461 ///
2462 /// See IMPORTANT REMARKS in TTree::SetBranchStatus and TChain::SetBranchAddress
2463 ///
2464 /// If found is not 0, the number of branch(es) found matching the regular
2465 /// expression is returned in *found AND the error message 'unknown branch'
2466 /// is suppressed.
2467 
2468 void TChain::SetBranchStatus(const char* bname, Bool_t status, UInt_t* found)
2469 {
2470  // FIXME: We never explicitly set found to zero!
2471 
2472  // Check if bname is already in the status list,
2473  // if not create a TChainElement object and set its status.
2474  TChainElement* element = (TChainElement*) fStatus->FindObject(bname);
2475  if (element) {
2476  fStatus->Remove(element);
2477  } else {
2478  element = new TChainElement(bname, "");
2479  }
2480  fStatus->Add(element);
2481  element->SetStatus(status);
2482  // Also set status in current tree.
2483  if (fTreeNumber >= 0) {
2484  fTree->SetBranchStatus(bname, status, found);
2485  } else if (found) {
2486  *found = 1;
2487  }
2488 }
2489 
2490 ////////////////////////////////////////////////////////////////////////////////
2491 /// Remove reference to this chain from current directory and add
2492 /// reference to new directory dir. dir can be 0 in which case the chain
2493 /// does not belong to any directory.
2494 
2496 {
2497  if (fDirectory == dir) return;
2498  if (fDirectory) fDirectory->Remove(this);
2499  fDirectory = dir;
2500  if (fDirectory) {
2501  fDirectory->Append(this);
2502  fFile = fDirectory->GetFile();
2503  } else {
2504  fFile = 0;
2505  }
2506 }
2507 
2508 ////////////////////////////////////////////////////////////////////////////////
2509 /// Set the input entry list (processing the entries of the chain will then be
2510 /// limited to the entries in the list)
2511 /// This function finds correspondance between the sub-lists of the TEntryList
2512 /// and the trees of the TChain
2513 /// By default (opt=""), both the file names of the chain elements and
2514 /// the file names of the TEntryList sublists are expanded to full path name.
2515 /// If opt = "ne", the file names are taken as they are and not expanded
2516 
2518 {
2519  if (fEntryList){
2520  //check, if the chain is the owner of the previous entry list
2521  //(it happens, if the previous entry list was created from a user-defined
2522  //TEventList in SetEventList() function)
2523  if (fEntryList->TestBit(kCanDelete)) {
2524  TEntryList *tmp = fEntryList;
2525  fEntryList = 0; // Avoid problem with RecursiveRemove.
2526  delete tmp;
2527  } else {
2528  fEntryList = 0;
2529  }
2530  }
2531  if (!elist){
2532  fEntryList = 0;
2533  fEventList = 0;
2534  return;
2535  }
2536  if (!elist->TestBit(kCanDelete)){
2537  //this is a direct call to SetEntryList, not via SetEventList
2538  fEventList = 0;
2539  }
2540  if (elist->GetN() == 0){
2541  fEntryList = elist;
2542  return;
2543  }
2544  if (fProofChain){
2545  //for processing on proof, event list and entry list can't be
2546  //set at the same time.
2547  fEventList = 0;
2548  fEntryList = elist;
2549  return;
2550  }
2551 
2552  Int_t ne = fFiles->GetEntries();
2553  Int_t listfound=0;
2554  TString treename, filename;
2555 
2556  TEntryList *templist = 0;
2557  for (Int_t ie = 0; ie<ne; ie++){
2558  auto chainElement = (TChainElement*)fFiles->UncheckedAt(ie);
2559  treename = chainElement->GetName();
2560  filename = chainElement->GetTitle();
2561  templist = elist->GetEntryList(treename, filename, opt);
2562  if (templist) {
2563  listfound++;
2564  templist->SetTreeNumber(ie);
2565  }
2566  }
2567 
2568  if (listfound == 0){
2569  Error("SetEntryList", "No list found for the trees in this chain");
2570  fEntryList = 0;
2571  return;
2572  }
2573  fEntryList = elist;
2574  TList *elists = elist->GetLists();
2575  Bool_t shift = kFALSE;
2576  TIter next(elists);
2577 
2578  //check, if there are sub-lists in the entry list, that don't
2579  //correspond to any trees in the chain
2580  while((templist = (TEntryList*)next())){
2581  if (templist->GetTreeNumber() < 0){
2582  shift = kTRUE;
2583  break;
2584  }
2585  }
2586  fEntryList->SetShift(shift);
2587 
2588 }
2589 
2590 ////////////////////////////////////////////////////////////////////////////////
2591 /// Set the input entry list (processing the entries of the chain will then be
2592 /// limited to the entries in the list). This function creates a special kind
2593 /// of entry list (TEntryListFromFile object) that loads lists, corresponding
2594 /// to the chain elements, one by one, so that only one list is in memory at a time.
2595 ///
2596 /// If there is an error opening one of the files, this file is skipped and the
2597 /// next file is loaded
2598 ///
2599 /// File naming convention:
2600 ///
2601 /// - by default, filename_elist.root is used, where filename is the
2602 /// name of the chain element
2603 /// - xxx$xxx.root - $ sign is replaced by the name of the chain element
2604 ///
2605 /// If the list name is not specified (by passing filename_elist.root/listname to
2606 /// the TChain::SetEntryList() function, the first object of class TEntryList
2607 /// in the file is taken.
2608 ///
2609 /// It is assumed, that there are as many list files, as there are elements in
2610 /// the chain and they are in the same order
2611 
2612 void TChain::SetEntryListFile(const char *filename, Option_t * /*opt*/)
2613 {
2614 
2615  if (fEntryList){
2616  //check, if the chain is the owner of the previous entry list
2617  //(it happens, if the previous entry list was created from a user-defined
2618  //TEventList in SetEventList() function)
2619  if (fEntryList->TestBit(kCanDelete)) {
2620  TEntryList *tmp = fEntryList;
2621  fEntryList = 0; // Avoid problem with RecursiveRemove.
2622  delete tmp;
2623  } else {
2624  fEntryList = 0;
2625  }
2626  }
2627 
2628  fEventList = 0;
2629 
2630  TString basename(filename);
2631 
2632  Int_t dotslashpos = basename.Index(".root/");
2633  TString behind_dot_root = "";
2634  if (dotslashpos>=0) {
2635  // Copy the list name specification
2636  behind_dot_root = basename(dotslashpos+6,basename.Length()-dotslashpos+6);
2637  // and remove it from basename
2638  basename.Remove(dotslashpos+5);
2639  }
2640  fEntryList = new TEntryListFromFile(basename.Data(), behind_dot_root.Data(), fNtrees);
2643  ((TEntryListFromFile*)fEntryList)->SetFileNames(fFiles);
2644 }
2645 
2646 ////////////////////////////////////////////////////////////////////////////////
2647 /// This function transfroms the given TEventList into a TEntryList
2648 ///
2649 /// NOTE, that this function loads all tree headers, because the entry numbers
2650 /// in the TEventList are global and have to be recomputed, taking into account
2651 /// the number of entries in each tree.
2652 ///
2653 /// The new TEntryList is owned by the TChain and gets deleted when the chain
2654 /// is deleted. This TEntryList is returned by GetEntryList() function, and after
2655 /// GetEntryList() function is called, the TEntryList is not owned by the chain
2656 /// any more and will not be deleted with it.
2657 
2659 {
2660  fEventList = evlist;
2661  if (fEntryList) {
2662  if (fEntryList->TestBit(kCanDelete)) {
2663  TEntryList *tmp = fEntryList;
2664  fEntryList = 0; // Avoid problem with RecursiveRemove.
2665  delete tmp;
2666  } else {
2667  fEntryList = 0;
2668  }
2669  }
2670 
2671  if (!evlist) {
2672  fEntryList = 0;
2673  fEventList = 0;
2674  return;
2675  }
2676 
2677  if(fProofChain) {
2678  //on proof, fEventList and fEntryList shouldn't be set at the same time
2679  if (fEntryList){
2680  //check, if the chain is the owner of the previous entry list
2681  //(it happens, if the previous entry list was created from a user-defined
2682  //TEventList in SetEventList() function)
2683  if (fEntryList->TestBit(kCanDelete)){
2684  TEntryList *tmp = fEntryList;
2685  fEntryList = 0; // Avoid problem with RecursiveRemove.
2686  delete tmp;
2687  } else {
2688  fEntryList = 0;
2689  }
2690  }
2691  return;
2692  }
2693 
2694  char enlistname[100];
2695  snprintf(enlistname,100, "%s_%s", evlist->GetName(), "entrylist");
2696  TEntryList *enlist = new TEntryList(enlistname, evlist->GetTitle());
2697  enlist->SetDirectory(0);
2698 
2699  Int_t nsel = evlist->GetN();
2700  Long64_t globalentry, localentry;
2701  const char *treename;
2702  const char *filename;
2704  //Load all the tree headers if the tree offsets are not known
2705  //It is assumed here, that loading the last tree will load all
2706  //previous ones
2707  printf("loading trees\n");
2708  (const_cast<TChain*>(this))->LoadTree(evlist->GetEntry(evlist->GetN()-1));
2709  }
2710  for (Int_t i=0; i<nsel; i++){
2711  globalentry = evlist->GetEntry(i);
2712  //add some protection from globalentry<0 here
2713  Int_t treenum = 0;
2714  while (globalentry>=fTreeOffset[treenum])
2715  treenum++;
2716  treenum--;
2717  localentry = globalentry - fTreeOffset[treenum];
2718  // printf("globalentry=%lld, treeoffset=%lld, localentry=%lld\n", globalentry, fTreeOffset[treenum], localentry);
2719  treename = ((TNamed*)fFiles->At(treenum))->GetName();
2720  filename = ((TNamed*)fFiles->At(treenum))->GetTitle();
2721  //printf("entering for tree %s %s\n", treename, filename);
2722  enlist->SetTree(treename, filename);
2723  enlist->Enter(localentry);
2724  }
2725  enlist->SetBit(kCanDelete, kTRUE);
2726  enlist->SetReapplyCut(evlist->GetReapplyCut());
2727  SetEntryList(enlist);
2728 }
2729 
2730 ////////////////////////////////////////////////////////////////////////////////
2731 /// Set number of entries per packet for parallel root.
2732 
2734 {
2735  fPacketSize = size;
2736  TIter next(fFiles);
2737  TChainElement *element;
2738  while ((element = (TChainElement*)next())) {
2739  element->SetPacketSize(size);
2740  }
2741 }
2742 
2743 ////////////////////////////////////////////////////////////////////////////////
2744 /// Enable/Disable PROOF processing on the current default Proof (gProof).
2745 ///
2746 /// "Draw" and "Processed" commands will be handled by PROOF.
2747 /// The refresh and gettreeheader are meaningful only if on == kTRUE.
2748 /// If refresh is kTRUE the underlying fProofChain (chain proxy) is always
2749 /// rebuilt (even if already existing).
2750 /// If gettreeheader is kTRUE the header of the tree will be read from the
2751 /// PROOF cluster: this is only needed for browsing and should be used with
2752 /// care because it may take a long time to execute.
2753 
2754 void TChain::SetProof(Bool_t on, Bool_t refresh, Bool_t gettreeheader)
2755 {
2756  if (!on) {
2757  // Disable
2759  // Reset related bit
2761  } else {
2762  if (fProofChain && !refresh &&
2763  (!gettreeheader || (gettreeheader && fProofChain->GetTree()))) {
2764  return;
2765  }
2768 
2769  // Make instance of TChainProof via the plugin manager
2770  TPluginHandler *h;
2771  if ((h = gROOT->GetPluginManager()->FindHandler("TChain", "proof"))) {
2772  if (h->LoadPlugin() == -1)
2773  return;
2774  if (!(fProofChain = reinterpret_cast<TChain *>(h->ExecPlugin(2, this, gettreeheader))))
2775  Error("SetProof", "creation of TProofChain failed");
2776  // Set related bits
2778  }
2779  }
2780 }
2781 
2782 ////////////////////////////////////////////////////////////////////////////////
2783 /// Set chain weight.
2784 ///
2785 /// The weight is used by TTree::Draw to automatically weight each
2786 /// selected entry in the resulting histogram.
2787 /// For example the equivalent of
2788 /// ~~~ {.cpp}
2789 /// chain.Draw("x","w")
2790 /// ~~~
2791 /// is
2792 /// ~~~ {.cpp}
2793 /// chain.SetWeight(w,"global");
2794 /// chain.Draw("x");
2795 /// ~~~
2796 /// By default the weight used will be the weight
2797 /// of each Tree in the TChain. However, one can force the individual
2798 /// weights to be ignored by specifying the option "global".
2799 /// In this case, the TChain global weight will be used for all Trees.
2800 
2802 {
2803  fWeight = w;
2804  TString opt = option;
2805  opt.ToLower();
2807  if (opt.Contains("global")) {
2809  }
2810 }
2811 
2812 ////////////////////////////////////////////////////////////////////////////////
2813 /// Stream a class object.
2814 
2815 void TChain::Streamer(TBuffer& b)
2816 {
2817  if (b.IsReading()) {
2818  // Remove using the 'old' name.
2819  {
2821  gROOT->GetListOfCleanups()->Remove(this);
2822  }
2823 
2824  UInt_t R__s, R__c;
2825  Version_t R__v = b.ReadVersion(&R__s, &R__c);
2826  if (R__v > 2) {
2827  b.ReadClassBuffer(TChain::Class(), this, R__v, R__s, R__c);
2828  } else {
2829  //====process old versions before automatic schema evolution
2830  TTree::Streamer(b);
2831  b >> fTreeOffsetLen;
2832  b >> fNtrees;
2833  fFiles->Streamer(b);
2834  if (R__v > 1) {
2835  fStatus->Streamer(b);
2837  b.ReadFastArray(fTreeOffset,fTreeOffsetLen);
2838  }
2839  b.CheckByteCount(R__s, R__c, TChain::IsA());
2840  //====end of old versions
2841  }
2842  // Re-add using the new name.
2843  {
2845  gROOT->GetListOfCleanups()->Add(this);
2846  }
2847 
2848  } else {
2849  b.WriteClassBuffer(TChain::Class(),this);
2850  }
2851 }
2852 
2853 ////////////////////////////////////////////////////////////////////////////////
2854 /// Dummy function kept for back compatibility.
2855 /// The cache is now activated automatically when processing TTrees/TChain.
2856 
2857 void TChain::UseCache(Int_t /* maxCacheSize */, Int_t /* pageSize */)
2858 {
2859 }
virtual Int_t LoadBaskets(Long64_t maxmemory)
Dummy function.
Definition: TChain.cxx:1233
virtual void SetBaddressIsPtr(Bool_t isptr)
Definition: TChainElement.h:67
virtual Bool_t GetReapplyCut() const
Definition: TEventList.h:57
virtual void SetAutoDelete(Bool_t autodel=kTRUE)
Set the global branch kAutoDelete bit.
Definition: TChain.cxx:2296
virtual TBranch * FindBranch(const char *name)
Return the branch that correspond to the path &#39;branchname&#39;, which can include the name of the tree or...
Definition: TTree.cxx:4589
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition: TSystem.cxx:931
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Double_t fWeight
Tree weight (see TTree::SetWeight)
Definition: TTree.h:87
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1272
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:32
Bool_t IsReading() const
Definition: TBuffer.h:81
virtual TList * GetListOfClones()
Definition: TTree.h:404
virtual UInt_t GetBaddressType() const
Definition: TChainElement.h:56
An array of TObjects.
Definition: TObjArray.h:37
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition: TROOT.cxx:2569
virtual void Browse(TBrowser *)
Browse the contents of the chain.
Definition: TChain.cxx:712
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual TObjArray * GetListOfLeaves()
Return a pointer to the list of leaves of the current tree.
Definition: TChain.cxx:1093
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:409
virtual TTree * GetTree()
Return pointer to friend TTree.
virtual void SetAddress(void *add)
Set address of this branch.
Definition: TBranch.cxx:2148
long long Long64_t
Definition: RtypesCore.h:69
virtual Long64_t GetEntryNumber(Long64_t entry) const
Return entry number corresponding to entry.
Definition: TChain.cxx:966
virtual Long64_t GetN() const
Definition: TEntryList.h:75
virtual void ls(Option_t *option="") const
List the chain.
Definition: TChain.cxx:1791
virtual Bool_t GetBranchStatus(const char *branchname) const
See TTree::GetReadEntry().
Definition: TChain.cxx:885
virtual const char * WorkingDirectory()
Return working directory.
Definition: TSystem.cxx:868
virtual TLeaf * GetLeaf(const char *branchname, const char *leafname)
Return a pointer to the leaf name in the current tree.
Definition: TChain.cxx:1023
virtual void SetBranchStatus(const char *bname, Bool_t status=1, UInt_t *found=0)
Set branch status to Process or DoNotProcess.
Definition: TTree.cxx:8019
virtual void Reset(Option_t *option="")
Resets the state of this chain.
Definition: TChain.cxx:2240
short Version_t
Definition: RtypesCore.h:61
Manages entry lists from different files, when they are not loaded in memory at the same time...
Collectable string class.
Definition: TObjString.h:28
Int_t fNtrees
Number of trees.
Definition: TChain.h:37
const char Option_t
Definition: RtypesCore.h:62
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:329
virtual Int_t GetTreeNumber() const
Definition: TEntryList.h:78
virtual void Print(Option_t *option="") const
Print a summary of the tree contents.
Definition: TTree.cxx:6841
virtual Int_t GetStatus() const
Definition: TChainElement.h:62
Bool_t fCanDeleteRefs
! If true, TProcessIDs are deleted when closing a file
Definition: TChain.h:40
virtual void SetCacheRead(TFileCacheRead *cache, TObject *tree=0, ECacheAction action=kDisconnect)
Set a pointer to the read cache.
Definition: TFile.cxx:2201
virtual void * GetBaddress() const
Definition: TChainElement.h:53
const Ssiz_t kNPOS
Definition: RtypesCore.h:115
This class represents a WWW compatible URL.
Definition: TUrl.h:35
TList * fFriends
pointer to list of friend elements
Definition: TTree.h:123
TUrl * GetCurrentUrl() const
Return the current url.
Definition: TFileInfo.cxx:248
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
virtual Int_t SetCacheSize(Long64_t cacheSize=-1)
Set maximum size of the file cache .
Definition: TChain.cxx:2305
virtual void ResetBranchAddress(TBranch *)
Reset the addresses of the branch.
Definition: TChain.cxx:2330
virtual Int_t GetPacketSize() const
Definition: TTree.h:422
TH1 * h
Definition: legend2.C:5
virtual TList * GetLists() const
Definition: TEntryList.h:73
virtual void SetWeight(Double_t w=1, Option_t *option="")
Set chain weight.
Definition: TChain.cxx:2801
A specialized TFileCacheRead object for a TTree.
Definition: TTreeCache.h:30
virtual TLeaf * GetLeaf(const char *branchname, const char *leafname)
Return pointer to the 1st Leaf named name in any Branch of this Tree or any branch in the list of fri...
Definition: TTree.cxx:5854
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:46
virtual void SetDirectory(TDirectory *dir)
Remove reference to this chain from current directory and add reference to new directory dir...
Definition: TChain.cxx:2495
virtual TList * GetListOfFriends() const
Definition: TTree.h:407
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
void SetLoadResult(Int_t result)
Definition: TChainElement.h:70
Regular expression class.
Definition: TRegexp.h:31
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
TDirectory * fDirectory
! Pointer to directory holding this tree
Definition: TTree.h:114
const char * GetFileAndOptions() const
Return the file and its options (the string specified behind the ?).
Definition: TUrl.cxx:501
virtual Double_t GetMaximum(const char *columname)
Return maximum of column with name columname.
Definition: TChain.cxx:1114
Int_t fMakeClass
! not zero when processing code generated by MakeClass
Definition: TTree.h:111
#define R__ASSERT(e)
Definition: TError.h:96
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:375
virtual TObject * Remove(TObject *obj)
Remove object from array.
Definition: TObjArray.cxx:653
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:587
virtual void SetAutoSave(Long64_t autos=-300000000)
This function may be called at the start of a program to change the default value for fAutoSave (and ...
Definition: TTree.cxx:7830
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
Definition: TTree.cxx:5321
Int_t LoadPlugin()
Load the plugin library for this handler.
Basic string class.
Definition: TString.h:129
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual void UpdateFormulaLeaves(const TTree *parent)=0
R__EXTERN TVirtualMutex * gROOTMutex
Definition: TROOT.h:57
virtual void Browse(TBrowser *)
Browse content of the TTree.
Definition: TTree.cxx:2520
virtual void SetShift(Bool_t shift)
Definition: TEntryList.h:99
const char * GetOptions() const
Definition: TUrl.h:74
virtual void SetLookedUp(Bool_t y=kTRUE)
Set/Reset the looked-up bit.
virtual TTree * CloneTree(Long64_t nentries=-1, Option_t *option="")
Create a clone of this tree and copy nentries.
Definition: TTree.cxx:3012
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
virtual Long64_t Scan(const char *varexp="", const char *selection="", Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Loop on tree and print entries passing selection.
Definition: TChain.cxx:2281
A TChainElement describes a component of a TChain.
Definition: TChainElement.h:28
virtual void Loop(Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Loop on nentries of this chain starting at firstentry. (NOT IMPLEMENTED)
Definition: TChain.cxx:1747
virtual void RecursiveRemove(TObject *obj)
Make sure that obj (which is being deleted or will soon be) is no longer referenced by this TTree...
Definition: TChain.cxx:2198
TVirtualTreePlayer * GetPlayer()
Load the TTreePlayer (if not already done).
Definition: TTree.cxx:5968
virtual void SetMaxVirtualSize(Long64_t size=0)
Definition: TTree.h:542
Bool_t MaybeWildcard() const
Returns true if string contains one of the wildcard characters "[]*?".
Definition: TString.cxx:908
TTree * fTree
! Pointer to current tree (Note: We do not own this tree.)
Definition: TChain.h:41
TString & Insert(Ssiz_t pos, const char *s)
Definition: TString.h:597
virtual void DirectoryAutoAdd(TDirectory *)
Override the TTree::DirectoryAutoAdd behavior: we never auto add.
Definition: TChain.cxx:751
virtual const char * GetAlias(const char *aliasName) const
Returns the expanded value of the alias. Search in the friends if any.
Definition: TChain.cxx:845
virtual Int_t GetFileNumber() const
Definition: TTree.h:394
void Reset()
Definition: TCollection.h:156
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
const char * GetUrl(Bool_t withDeflt=kFALSE) const
Return full URL.
Definition: TUrl.cxx:387
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:501
virtual const char * UnixPathName(const char *unixpathname)
Convert from a Unix pathname to a local pathname.
Definition: TSystem.cxx:1043
if object in a list can be deleted
Definition: TObject.h:58
virtual Double_t GetMinimum(const char *columname)
Return minimum of column with name columname.
Definition: TTree.cxx:5939
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition: TString.h:630
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3909
virtual Bool_t HasBeenLookedUp()
Definition: TChainElement.h:63
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
virtual void Sort(Bool_t order=kSortAscending)
Sort linked list.
Definition: TList.cxx:775
virtual Double_t GetMinimum(const char *columname)
Return minimum of column with name columname.
Definition: TChain.cxx:1131
virtual Double_t GetWeight() const
Definition: TTree.h:458
virtual void CanDeleteRefs(Bool_t flag=kTRUE)
When closing a file during the chain processing, the file may be closed with option "R" if flag is se...
Definition: TChain.cxx:730
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition: TChain.cxx:902
virtual Int_t AddFileInfoList(TCollection *list, Long64_t nfiles=TTree::kMaxEntries)
Add all files referenced in the list to the chain.
Definition: TChain.cxx:540
const char * GetFile() const
Definition: TUrl.h:72
virtual const char * GetDirEntry(void *dirp)
Get a directory entry. Returns 0 if no more entries.
Definition: TSystem.cxx:850
#define SafeDelete(p)
Definition: RConfig.h:499
virtual Long64_t Merge(const char *name, Option_t *option="")
Merge all the entries in the chain into a new tree in a new file.
Definition: TChain.cxx:1824
TObjArray * GetListOfFiles() const
Definition: TChain.h:107
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:405
Helper class to iterate over cluster of baskets.
Definition: TTree.h:241
TVirtualTreePlayer * fPlayer
! Pointer to current Tree player
Definition: TTree.h:126
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:135
virtual Long64_t Draw(const char *varexp, const TCut &selection, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Draw expression varexp for selected entries.
Definition: TChain.cxx:764
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2345
virtual void SetBaddressClassName(const char *clname)
Definition: TChainElement.h:66
void Class()
Definition: Class.C:29
virtual Long64_t CopyEntries(TTree *tree, Long64_t nentries=-1, Option_t *option="")
Copy nentries from given tree to this tree.
Definition: TTree.cxx:3374
virtual void SetPacketSize(Int_t size=100)
Set number of entries per packet for parallel root.
virtual void SetEntryListFile(const char *filename="", Option_t *opt="")
Set the input entry list (processing the entries of the chain will then be limited to the entries in ...
Definition: TChain.cxx:2612
virtual TBranch * FindBranch(const char *name)
See TTree::GetReadEntry().
Definition: TChain.cxx:803
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual Bool_t Notify()
This method must be overridden to handle object notification.
Definition: TObject.cxx:499
virtual Long64_t GetEntries() const
Definition: TChainElement.h:58
virtual Long64_t GetReadEntry() const
Definition: TTree.h:424
virtual const char * GetAlias(const char *aliasName) const
Returns the expanded value of the alias. Search in the friends if any.
Definition: TTree.cxx:4938
virtual Long64_t LoadTree(Long64_t entry)
Set current entry.
Definition: TTree.cxx:6121
Bool_t fCacheUserSet
! true if the cache setting was explicitly given by user
Definition: TTree.h:132
virtual void SetEntryList(TEntryList *elist, Option_t *opt="")
Set the input entry list (processing the entries of the chain will then be limited to the entries in ...
Definition: TChain.cxx:2517
virtual Long64_t GetEntry(Int_t index) const
Return value of entry at index in the list.
Definition: TEventList.cxx:222
virtual Int_t GetNbranches()
Definition: TTree.h:419
void Clear()
Clear string without changing its capacity.
Definition: TString.cxx:1150
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition: TTree.cxx:5152
void InvalidateCurrentTree()
Set the TTree to be reloaded as soon as possible.
Definition: TChain.cxx:1215
virtual Int_t GetN() const
Definition: TEventList.h:56
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Return true if string ends with the specified string.
Definition: TString.cxx:2231
UInt_t fFriendLockStatus
! Record which method is locking the friend recursion
Definition: TTree.h:129
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch name in the current tree.
Definition: TChain.cxx:864
virtual void ResetCache()
This will simply clear the cache.
virtual void SetBranchPtr(TBranch **ptr)
Definition: TChainElement.h:69
virtual TFriendElement * AddFriend(const char *treename, const char *filename="")
Add a TFriendElement to the list of friends.
Definition: TTree.cxx:1193
virtual Long64_t Scan(const char *varexp="", const char *selection="", Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Loop over tree entries and print entries passing selection.
Definition: TTree.cxx:7648
const char * GetAnchor() const
Definition: TUrl.h:73
virtual Long64_t LoadTree(Long64_t entry)
Find the tree which contains entry, and set it as the current tree.
Definition: TChain.cxx:1259
virtual TList * GetList() const
Definition: TDirectory.h:147
virtual Long64_t GetEntries() const
Return the total number of entries in the chain.
Definition: TChain.cxx:923
virtual Double_t GetMaximum(const char *columname)
Return maximum of column with name columname.
Definition: TTree.cxx:5900
virtual void RemoveFriend(TTree *)
Remove a friend from the list of friends.
Definition: TTree.cxx:7531
virtual TFile * GetFile() const
Definition: TDirectory.h:145
virtual void ls(Option_t *option="") const
The ls function lists the contents of a class on stdout.
Definition: TObject.cxx:485
Long64_t * fTreeOffset
[fTreeOffsetLen] Array of variables
Definition: TChain.h:39
virtual void SetStatus(Int_t status)
Definition: TChainElement.h:74
Int_t fTreeOffsetLen
Current size of fTreeOffset array.
Definition: TChain.h:36
virtual TTree * GetTree() const
Definition: TTree.h:432
virtual Int_t GetTreeNumber() const
Definition: TTree.h:434
A specialized string object used for TTree selections.
Definition: TCut.h:25
A doubly linked list.
Definition: TList.h:43
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:4979
Int_t fTreeNumber
! Current Tree number in fTreeOffset table
Definition: TChain.h:38
const char * GetName() const
Returns name of object.
Definition: TObjString.h:39
virtual void RemoveFriend(TTree *)
Remove a friend from the list of friends.
Definition: TChain.cxx:2217
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
TEntryList * fEntryList
! Pointer to event selection list (if one)
Definition: TTree.h:119
virtual Long64_t GetReadEntry() const
See TTree::GetReadEntry().
Definition: TChain.cxx:1165
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TTree.cxx:9163
virtual Bool_t GetBaddressIsPtr() const
Definition: TChainElement.h:55
virtual void UpdateFormulaLeaves()=0
void SetCompressionSettings(Int_t settings=1)
Set compression settings.
Definition: TBranch.cxx:2271
R__EXTERN TSystem * gSystem
Definition: TSystem.h:539
virtual Int_t GetNbranches()
Return the number of branches of the current tree.
Definition: TChain.cxx:1150
if object ctor succeeded but object should not be used
Definition: TObject.h:65
virtual Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor=0) const
Return entry number corresponding to major and minor number.
Definition: TTree.cxx:5584
Long_t ExecPlugin(int nargs, const T &... params)
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:679
Long64_t fReadEntry
! Number of the entry being processed
Definition: TTree.h:104
virtual TVirtualIndex * GetTreeIndex() const
Definition: TTree.h:433
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
Collection abstract base class.
Definition: TCollection.h:42
virtual TBranch ** GetBranchPtr() const
Definition: TChainElement.h:57
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2332
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void SetProof(Bool_t on=kTRUE, Bool_t refresh=kFALSE, Bool_t gettreeheader=kFALSE)
Enable/Disable PROOF processing on the current default Proof (gProof).
Definition: TChain.cxx:2754
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition: TTree.cxx:5163
Ssiz_t Length() const
Definition: TString.h:388
virtual void Append(TObject *obj, Bool_t replace=kFALSE)
Append object to this directory.
Definition: TDirectory.cxx:153
void ParseTreeFilename(const char *name, TString &filename, TString &treename, TString &query, TString &suffix, Bool_t wildcards) const
Get the tree url or filename and other information from the name.
Definition: TChain.cxx:2086
A TEventList object is a list of selected events (entries) in a TTree.
Definition: TEventList.h:31
virtual TLeaf * FindLeaf(const char *name)
Find leaf..
Definition: TTree.cxx:4661
TLine * l
Definition: textangle.C:4
virtual void SetNumberEntries(Long64_t n)
Definition: TChainElement.h:72
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:71
virtual void ResetBranchAddress(TBranch *)
Tell all of our branches to set their addresses to zero.
Definition: TTree.cxx:7619
TEventList * fEventList
! Pointer to event selection list (if one)
Definition: TTree.h:118
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
virtual Bool_t GetBranchStatus(const char *branchname) const
Return status of branch with name branchname.
Definition: TTree.cxx:5075
if object destructor must call RecursiveRemove()
Definition: TObject.h:59
virtual void FreeDirectory(void *dirp)
Free a directory.
Definition: TSystem.cxx:842
void AddClone(TTree *)
Add a cloned tree to our list of trees to be notified whenever we change our branch addresses or when...
Definition: TTree.cxx:1105
virtual TObjLink * FirstLink() const
Definition: TList.h:97
#define Printf
Definition: TGeoToOCC.h:18
virtual Int_t AddFile(const char *name, Long64_t nentries=TTree::kMaxEntries, const char *tname="")
Add a new file to this chain.
Definition: TChain.cxx:443
#define R__LOCKGUARD2(mutex)
const Bool_t kFALSE
Definition: RtypesCore.h:92
virtual void CreatePackets()
Initialize the packet descriptor string.
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:2195
virtual Int_t CheckBranchAddressType(TBranch *branch, TClass *ptrClass, EDataType datatype, Bool_t ptr)
Check whether or not the address described by the last 3 parameters matches the content of the branch...
Definition: TTree.cxx:2748
virtual Int_t GetEntryWithIndex(Int_t major, Int_t minor=0)
Return entry corresponding to major and minor number.
Definition: TChain.cxx:999
virtual Long64_t Process(const char *filename, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Process this tree executing the TSelector code in the specified filename.
Definition: TTree.cxx:7035
TString & Remove(Ssiz_t pos)
Definition: TString.h:621
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
int Ssiz_t
Definition: RtypesCore.h:63
virtual void SetMakeClass(Int_t make)
Set all the branches in this TTree to be in decomposed object mode (also known as MakeClass mode)...
Definition: TTree.cxx:8616
TFile * GetFile() const
Return a pointer to the current file.
Definition: TChain.cxx:1010
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Get entry from the file to memory.
Definition: TChain.cxx:947
virtual Long64_t LoadTreeFriend(Long64_t entry, TTree *T)
Load entry on behalf of our master tree, we may use an index.
Definition: TTree.cxx:6213
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:89
Long64_t fMaxVirtualSize
Maximum total size of buffers kept in memory.
Definition: TTree.h:96
virtual TObject * Remove(TObject *)
Remove an object from the in-memory list.
#define ClassImp(name)
Definition: Rtypes.h:336
Bool_t IsZombie() const
Definition: TObject.h:122
virtual TFriendElement * AddFriend(const char *chainname, const char *dummy="")
Add a TFriendElement to the list of friends of this chain.
Definition: TChain.cxx:630
virtual Long64_t GetEntryAndTree(Int_t index, Int_t &treenum)
Return the index of "index"-th non-zero entry in the TTree or TChain and the # of the corresponding t...
Definition: TEntryList.cxx:733
double Double_t
Definition: RtypesCore.h:55
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:875
Describe directory structure in memory.
Definition: TDirectory.h:34
void Lookup(Bool_t force=kFALSE)
Check / locate the files in the chain.
Definition: TChain.cxx:1677
virtual void ResetBranchAddresses()
Tell all of our branches to drop their current objects and allocate new ones.
Definition: TTree.cxx:7629
static RooMathCoreReg dummy
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:355
int nentries
Definition: THbookFile.cxx:89
virtual void UseCache(Int_t maxCacheSize=10, Int_t pageSize=0)
Dummy function kept for back compatibility.
Definition: TChain.cxx:2857
EDataType
Definition: TDataType.h:28
virtual void SetBranchStatus(const char *bname, Bool_t status=1, UInt_t *found=0)
Set branch status to Process or DoNotProcess.
Definition: TChain.cxx:2468
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
static TFileStager * Open(const char *stager)
Open a stager, after having loaded the relevant plug-in.
virtual Double_t GetWeight() const
Return the chain weight.
Definition: TChain.cxx:1186
virtual void SetBaddressType(UInt_t type)
Definition: TChainElement.h:68
Helper class to prevent infinite recursion in the usage of TTree Friends.
Definition: TTree.h:172
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2885
virtual Int_t SetCacheSize(Long64_t cachesize=-1)
Set maximum size of the file cache .
Definition: TTree.cxx:8161
virtual void Print(Option_t *option="") const
Print the header information of each tree in the chain.
Definition: TChain.cxx:2135
virtual void ResetBranchAddresses()
Reset the addresses of the branches.
Definition: TChain.cxx:2344
virtual Long64_t GetEntries() const
Definition: TTree.h:381
Bool_t IsNull() const
Definition: TString.h:385
Int_t GetNtrees() const
Definition: TChain.h:95
virtual void SetBaddress(void *add)
Definition: TChainElement.h:65
virtual Int_t GetPacketSize() const
Definition: TChainElement.h:61
virtual void ResetAfterMerge(TFileMergeInfo *)
Resets the state of this TTree after a merge (keep the customization but forget the data)...
Definition: TTree.cxx:7588
Mother of all ROOT objects.
Definition: TObject.h:37
virtual TTree * GetTree() const
Definition: TChain.h:115
static Int_t IncreaseDirLevel()
Increase the indentation level for ls().
Definition: TROOT.cxx:2624
Long64_t fEntries
Number of entries.
Definition: TTree.h:81
virtual const char * GetBaddressClassName() const
Definition: TChainElement.h:54
TList * fClones
! List of cloned trees which share our addresses
Definition: TTree.h:127
virtual Long64_t Process(const char *filename, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Process all entries in this chain, calling functions in filename.
Definition: TChain.cxx:2158
virtual void Add(TObject *obj)
Definition: TList.h:77
TChain()
Default constructor.
Definition: TChain.cxx:63
Definition: file.py:1
A TFriendElement TF describes a TTree object TF in a file.
A chain is a collection of files containing TTree objects.
Definition: TChain.h:33
virtual void ResetAfterMerge(TFileMergeInfo *)
Resets the state of this chain after a merge (keep the customization but forget the data)...
Definition: TChain.cxx:2262
virtual void SetPacketSize(Int_t size=100)
Set number of entries per packet for parallel root.
Definition: TChain.cxx:2733
virtual TEntryList * GetEntryList(const char *treename, const char *filename, Option_t *opt="")
Return the entry list, corresponding to treename and filename By default, the filename is first tried...
Definition: TEntryList.cxx:781
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 SetTreeNumber(Int_t index)
Definition: TEntryList.h:104
TObjArray * fFiles
-> List of file names containing the trees (TChainElement, owned)
Definition: TChain.h:43
#define snprintf
Definition: civetweb.c:822
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Set branch address.
Definition: TChain.cxx:2373
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:494
virtual void * OpenDirectory(const char *name)
Open a directory. Returns 0 if directory does not exist.
Definition: TSystem.cxx:833
virtual void SetCircular(Long64_t maxEntries)
Enable/Disable circularity for this tree.
Definition: TTree.cxx:8365
virtual void CreatePackets()
Initialize the packet descriptor string.
Definition: TChain.cxx:738
virtual void Reset(Option_t *option="")
Reset baskets, buffers and entries count in all branches and leaves.
Definition: TTree.cxx:7557
virtual Bool_t Matches(const char *s)
Definition: TFileStager.h:46
Definition: tree.py:1
TFileCacheRead * GetCacheRead(TObject *tree=0) const
Return a pointer to the current read cache.
Definition: TFile.cxx:1202
virtual char * GetAddress() const
Definition: TBranch.h:151
TObject * fNotify
! Object to be notified when loading a Tree
Definition: TTree.h:113
void Add(TObject *obj)
Definition: TObjArray.h:73
A TTree object has a header with a name and a title.
Definition: TTree.h:78
TFile * fFile
! Pointer to current file (We own the file).
Definition: TChain.h:42
TChain * fProofChain
! chain proxy when going to be processed by PROOF
Definition: TChain.h:45
#define gDirectory
Definition: TDirectory.h:211
Class describing a generic file including meta information.
Definition: TFileInfo.h:38
void ResetBit(UInt_t f)
Definition: TObject.h:158
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition: TSystem.cxx:1250
Definition: first.py:1
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:901
virtual Int_t Locate(const char *u, TString &f)
Just check if the file exists locally.
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
A TTree is a list of TBranches.
Definition: TBranch.h:57
virtual TLeaf * FindLeaf(const char *name)
See TTree::GetReadEntry().
Definition: TChain.cxx:824
A TSelector object is used by the TTree::Draw, TTree::Scan, TTree::Process to navigate in a TTree and...
Definition: TSelector.h:33
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:25
virtual void UpdateBranches(TTree *tree)
Update pointer to current Tree and recompute pointers to the branches in the cache.
const Bool_t kTRUE
Definition: RtypesCore.h:91
Int_t fPacketSize
! Number of entries in one packet for parallel root
Definition: TTree.h:106
virtual ~TChain()
Destructor.
Definition: TChain.cxx:177
virtual TObjArray * GetListOfBranches()
Return a pointer to the list of branches of the current tree.
Definition: TChain.cxx:1070
static constexpr Long64_t kMaxEntries
Definition: TTree.h:213
virtual Long64_t GetChainEntryNumber(Long64_t entry) const
Return absolute entry number in the chain.
Definition: TChain.cxx:913
virtual void SetEventList(TEventList *evlist)
This function transfroms the given TEventList into a TEntryList.
Definition: TChain.cxx:2658
const char * cnt
Definition: TXMLSetup.cxx:75
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859
TList * fStatus
-> List of active/inactive branches (TChainElement, owned)
Definition: TChain.h:44
virtual void SetName(const char *name)
Change the name of this tree.
Definition: TTree.cxx:8644
virtual Int_t Add(TChain *chain)
Add all files referenced by the passed chain to this chain.
Definition: TChain.cxx:220
Int_t GetCompressionSettings() const
Definition: TFile.h:374
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 void CopyAddresses(TTree *, Bool_t undo=kFALSE)
Set branch addresses of passed tree equal to ours.
Definition: TTree.cxx:3167
virtual TObjArray * GetListOfLeaves()
Definition: TTree.h:406
virtual void Close(Option_t *option="")
Close a file.
Definition: TFile.cxx:904
virtual void ls(Option_t *option="") const
List files in the chain.
Long64_t fCacheSize
! Maximum size of file buffers
Definition: TTree.h:102
virtual void SetChainOffset(Long64_t offset=0)
Definition: TTree.h:528
const char * Data() const
Definition: TString.h:347
virtual void SetAutoDelete(Bool_t autodel=kTRUE)
Set the automatic delete bit.
Definition: TBranch.cxx:2182