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