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