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