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