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