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