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