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