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