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