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