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
711TFriendElement* TChain::AddFriend(const char* chain, TFile* dummy)
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();
742 TFriendElement *fe = new TFriendElement(this,chain,alias);
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;
992 const auto readEntry = fReadEntry;
993 auto *thisChain = const_cast<TChain *>(this);
994 thisChain->LoadTree(TTree::kMaxEntries - 1);
995 thisChain->InvalidateCurrentTree();
996 if (readEntry >= 0)
997 thisChain->LoadTree(readEntry);
998 else
999 thisChain->fReadEntry = readEntry;
1000 }
1001 return fEntries;
1002}
1003
1004////////////////////////////////////////////////////////////////////////////////
1005/// Get entry from the file to memory.
1006///
1007/// - getall = 0 : get only active branches
1008/// - getall = 1 : get all branches
1009///
1010/// Return the total number of bytes read,
1011/// 0 bytes read indicates a failure.
1012
1014{
1016 if (treeReadEntry < 0) {
1017 return 0;
1018 }
1019 if (!fTree) {
1020 return 0;
1021 }
1023}
1024
1025////////////////////////////////////////////////////////////////////////////////
1026/// Return entry number corresponding to entry.
1027///
1028/// if no TEntryList set returns entry
1029/// else returns entry \#entry from this entry list and
1030/// also computes the global entry number (loads all tree headers)
1031
1033{
1034
1035 if (fEntryList){
1036 Int_t treenum = 0;
1038 //find the global entry number
1039 //same const_cast as in the GetEntries() function
1040 if (localentry<0) return -1;
1041 if (treenum != fTreeNumber){
1043 for (Int_t i=0; i<=treenum; i++){
1045 (const_cast<TChain*>(this))->LoadTree(fTreeOffset[i-1]);
1046 }
1047 }
1048 //(const_cast<TChain*>(this))->LoadTree(fTreeOffset[treenum]);
1049 }
1051 return globalentry;
1052 }
1053 return entry;
1054}
1055
1056////////////////////////////////////////////////////////////////////////////////
1057/// Return entry corresponding to major and minor number.
1058///
1059/// The function returns the total number of bytes read.
1060/// If the Tree has friend trees, the corresponding entry with
1061/// the index values (major,minor) is read. Note that the master Tree
1062/// and its friend may have different entry serial numbers corresponding
1063/// to (major,minor).
1064
1066{
1068 if (serial < 0) return -1;
1069 return GetEntry(serial);
1070}
1071
1072////////////////////////////////////////////////////////////////////////////////
1073/// Return a pointer to the current file.
1074/// If no file is connected, the first file is automatically loaded.
1075
1077{
1078 if (fFile) {
1079 return fFile;
1080 }
1081 // Force opening the first file in the chain.
1082 const_cast<TChain*>(this)->LoadTree(0);
1083 return fFile;
1084}
1085
1086////////////////////////////////////////////////////////////////////////////////
1087/// Return a pointer to the leaf name in the current tree.
1088
1089TLeaf* TChain::GetLeaf(const char* branchname, const char *leafname)
1090{
1092 // Make sure the element list is up to date
1093 if (!TestBit(kProofUptodate))
1094 SetProof(true, true);
1096 }
1097 if (fTree) {
1098 return fTree->GetLeaf(branchname, leafname);
1099 }
1100 LoadTree(0);
1101 if (fTree) {
1102 return fTree->GetLeaf(branchname, leafname);
1103 }
1104 return nullptr;
1105}
1106
1107////////////////////////////////////////////////////////////////////////////////
1108/// Return a pointer to the leaf name in the current tree.
1109
1111{
1113 // Make sure the element list is up to date
1114 if (!TestBit(kProofUptodate))
1115 SetProof(true, true);
1116 return fProofChain->GetLeaf(name);
1117 }
1118 if (fTree) {
1119 return fTree->GetLeaf(name);
1120 }
1121 LoadTree(0);
1122 if (fTree) {
1123 return fTree->GetLeaf(name);
1124 }
1125 return nullptr;
1126}
1127
1128////////////////////////////////////////////////////////////////////////////////
1129/// Return a pointer to the list of branches of the current tree.
1130///
1131/// Warning: If there is no current TTree yet, this routine will open the
1132/// first in the chain.
1133///
1134/// Returns 0 on failure.
1135
1137{
1139 // Make sure the element list is up to date
1140 if (!TestBit(kProofUptodate))
1141 SetProof(true, true);
1143 }
1144 if (fTree) {
1145 return fTree->GetListOfBranches();
1146 }
1147 LoadTree(0);
1148 if (fTree) {
1149 return fTree->GetListOfBranches();
1150 }
1151 return nullptr;
1152}
1153
1154////////////////////////////////////////////////////////////////////////////////
1155/// Return a pointer to the list of leaves of the current tree.
1156///
1157/// Warning: May set the current tree!
1158
1160{
1162 // Make sure the element list is up to date
1163 if (!TestBit(kProofUptodate))
1164 SetProof(true, true);
1165 return fProofChain->GetListOfLeaves();
1166 }
1167 if (fTree) {
1168 return fTree->GetListOfLeaves();
1169 }
1170 LoadTree(0);
1171 if (fTree) {
1172 return fTree->GetListOfLeaves();
1173 }
1174 return nullptr;
1175}
1176
1177////////////////////////////////////////////////////////////////////////////////
1178/// Return maximum of column with name columname.
1179
1181{
1183 TLeaf *leaf = nullptr;
1184 TBranch *branch = nullptr;
1185 Int_t treenumber = -1;
1186 for (Long64_t i = 0; i < fEntries; ++i) {
1188 if (entryNumber < 0)
1189 break;
1190 Long64_t localEntryNumber = this->LoadTree(entryNumber);
1191 if (localEntryNumber < 0)
1192 break;
1193 if (treenumber != this->GetTreeNumber()) {
1194 leaf = this->GetLeaf(columname);
1195 if (leaf)
1196 branch = leaf->GetBranch();
1197 }
1198 treenumber = this->GetTreeNumber();
1199 if (!branch)
1200 continue;
1201 branch->GetEntry(localEntryNumber);
1202 for (Int_t j = 0; j < leaf->GetLen(); ++j) {
1203 Double_t val = leaf->GetValue(j);
1204 if (val > cmax) {
1205 cmax = val;
1206 }
1207 }
1208 }
1209 return cmax;
1210}
1211
1212////////////////////////////////////////////////////////////////////////////////
1213/// Return minimum of column with name columname.
1214
1216{
1218 TLeaf *leaf = nullptr;
1219 TBranch *branch = nullptr;
1220 Int_t treenumber = -1;
1221 for (Long64_t i = 0; i < fEntries; ++i) {
1223 if (entryNumber < 0)
1224 break;
1225 Long64_t localEntryNumber = this->LoadTree(entryNumber);
1226 if (localEntryNumber < 0)
1227 break;
1228 if (treenumber != this->GetTreeNumber()) {
1229 leaf = this->GetLeaf(columname);
1230 if (leaf)
1231 branch = leaf->GetBranch();
1232 }
1233 treenumber = this->GetTreeNumber();
1234 if (!branch)
1235 continue;
1236 branch->GetEntry(localEntryNumber);
1237 for (Int_t j = 0; j < leaf->GetLen(); ++j) {
1238 Double_t val = leaf->GetValue(j);
1239 if (val < cmin) {
1240 cmin = val;
1241 }
1242 }
1243 }
1244 return cmin;
1245}
1246
1247////////////////////////////////////////////////////////////////////////////////
1248/// Return the number of branches of the current tree.
1249///
1250/// Warning: May set the current tree!
1251
1253{
1254 if (fTree) {
1255 return fTree->GetNbranches();
1256 }
1257 LoadTree(0);
1258 if (fTree) {
1259 return fTree->GetNbranches();
1260 }
1261 return 0;
1262}
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// See TTree::GetReadEntry().
1266
1268{
1270 // Make sure the element list is up to date
1271 if (!TestBit(kProofUptodate))
1272 Warning("GetBranchStatus", "PROOF proxy not up-to-date:"
1273 " run TChain::SetProof(true, true) first");
1274 return fProofChain->GetReadEntry();
1275 }
1276 return TTree::GetReadEntry();
1277}
1278
1279////////////////////////////////////////////////////////////////////////////////
1280/// Return the chain weight.
1281///
1282/// By default the weight is the weight of the current tree.
1283/// However, if the weight has been set in TChain::SetWeight()
1284/// with the option "global", then that weight will be returned.
1285///
1286/// Warning: May set the current tree!
1287
1289{
1290 if (TestBit(kGlobalWeight)) {
1291 return fWeight;
1292 } else {
1293 if (fTree) {
1294 return fTree->GetWeight();
1295 }
1296 const_cast<TChain*>(this)->LoadTree(0);
1297 if (fTree) {
1298 return fTree->GetWeight();
1299 }
1300 return 0;
1301 }
1302}
1303
1304////////////////////////////////////////////////////////////////////////////////
1305/// Move content to a new file. (NOT IMPLEMENTED for TChain)
1306bool TChain::InPlaceClone(TDirectory * /* new directory */, const char * /* options */)
1307{
1308 Error("InPlaceClone", "not implemented");
1309 return false;
1310}
1311
1312////////////////////////////////////////////////////////////////////////////////
1313/// Set the TTree to be reloaded as soon as possible. In particular this
1314/// is needed when adding a Friend.
1315///
1316/// If the tree has clones, copy them into the chain
1317/// clone list so we can change their branch addresses
1318/// when necessary.
1319///
1320/// This is to support the syntax:
1321/// ~~~ {.cpp}
1322/// TTree* clone = chain->GetTree()->CloneTree(0);
1323/// ~~~
1324
1326{
1327 if (fTree && fTree->GetListOfClones()) {
1328 for (TObjLink* lnk = fTree->GetListOfClones()->FirstLink(); lnk; lnk = lnk->Next()) {
1329 TTree* clone = (TTree*) lnk->GetObject();
1330 AddClone(clone);
1331 }
1332 }
1333 fTreeNumber = -1;
1334 fTree = nullptr;
1335}
1336
1337////////////////////////////////////////////////////////////////////////////////
1338/// Dummy function.
1339/// It could be implemented and load all baskets of all trees in the chain.
1340/// For the time being use TChain::Merge and TTree::LoadBasket
1341/// on the resulting tree.
1342
1344{
1345 Error("LoadBaskets", "Function not yet implemented for TChain.");
1346 return 0;
1347}
1348
1349////////////////////////////////////////////////////////////////////////////////
1350/// Find the tree which contains entry, and set it as the current tree.
1351///
1352/// Returns the entry number in that tree.
1353///
1354/// The input argument entry is the entry serial number in the whole chain.
1355///
1356/// In case of error, LoadTree returns a negative number:
1357/// * -1: The chain is empty.
1358/// * -2: The requested entry number is less than zero or too large for the chain.
1359/// * -3: The file corresponding to the entry could not be correctly opened
1360/// * -4: The TChainElement corresponding to the entry is missing or
1361/// the TTree is missing from the file.
1362/// * -5: Internal error, please report the circumstance when this happen
1363/// as a ROOT issue.
1364/// * -6: An error occurred within the notify callback.
1365///
1366/// Calls fNotify->Notify() (if fNotify is not null) when starting the processing of a new sub-tree.
1367/// See TNotifyLink for more information on the notification mechanism.
1368///
1369/// \note This is the only routine which sets the value of fTree to a non-zero pointer.
1370///
1372{
1373 // We already have been visited while recursively looking
1374 // through the friends tree, let's return.
1376 return 0;
1377 }
1378
1379 if (!fNtrees) {
1380 // -- The chain is empty.
1381 return -1;
1382 }
1383
1384 if ((entry < 0) || ((entry > 0) && (entry >= fEntries && entry!=(TTree::kMaxEntries-1) ))) {
1385 // -- Invalid entry number.
1386 if (fTree) fTree->LoadTree(-1);
1387 fReadEntry = -1;
1388 return -2;
1389 }
1390
1391 // Find out which tree in the chain contains the passed entry.
1394 // -- Entry is *not* in the chain's current tree.
1395 // Do a linear search of the tree offset array.
1396 // FIXME: We could be smarter by starting at the
1397 // current tree number and going forwards,
1398 // then wrapping around at the end.
1399 for (treenum = 0; treenum < fNtrees; treenum++) {
1400 if (entry < fTreeOffset[treenum+1]) {
1401 break;
1402 }
1403 }
1404 }
1405
1406 // Calculate the entry number relative to the found tree.
1408 fReadEntry = entry;
1409
1410 // If entry belongs to the current tree return entry.
1411 if (fTree && treenum == fTreeNumber) {
1412 // First set the entry the tree on its owns friends
1413 // (the friends of the chain will be updated in the
1414 // next loop).
1416
1417 if (fFriends) {
1418 // The current tree has not changed but some of its friends might.
1419 //
1420 TIter next(fFriends);
1421 TFriendLock lock(this, kLoadTree);
1422 TFriendElement* fe = nullptr;
1423 while ((fe = (TFriendElement*) next())) {
1424 TTree* at = fe->GetTree();
1425 // If the tree is a
1426 // direct friend of the chain, it should be scanned
1427 // used the chain entry number and NOT the tree entry
1428 // number (treeReadEntry) hence we do:
1429 at->LoadTreeFriend(entry, this);
1430 }
1431 bool needUpdate = false;
1432 if (fTree->GetListOfFriends()) {
1434 if (fetree->IsUpdated()) {
1435 needUpdate = true;
1436 fetree->ResetUpdated();
1437 }
1438 }
1439 }
1440 if (needUpdate) {
1441 // Update the branch/leaf addresses and
1442 // the list of leaves in all TTreeFormula of the TTreePlayer (if any).
1443
1444 // Set the branch statuses for the newly opened file.
1447 while ((frelement = (TChainElement*) fnext())) {
1448 Int_t status = frelement->GetStatus();
1449 fTree->SetBranchStatus(frelement->GetName(), status);
1450 }
1451
1452 // Set the branch addresses for the newly opened file.
1453 fnext.Reset();
1454 while ((frelement = (TChainElement*) fnext())) {
1455 void* addr = frelement->GetBaddress();
1456 if (addr) {
1457 TBranch* br = fTree->GetBranch(frelement->GetName());
1458 TBranch** pp = frelement->GetBranchPtr();
1459 if (pp) {
1460 // FIXME: What if br is zero here?
1461 *pp = br;
1462 }
1463 if (br) {
1464 if (!frelement->GetCheckedType()) {
1465 Int_t res = CheckBranchAddressType(br, TClass::GetClass(frelement->GetBaddressClassName()),
1466 (EDataType) frelement->GetBaddressType(), frelement->GetBaddressIsPtr());
1467 if ((res & kNeedEnableDecomposedObj) && !br->GetMakeClass()) {
1468 br->SetMakeClass(true);
1469 }
1470 frelement->SetDecomposedObj(br->GetMakeClass());
1471 frelement->SetCheckedType(true);
1472 }
1473 // FIXME: We may have to tell the branch it should
1474 // not be an owner of the object pointed at.
1475 br->SetAddress(addr);
1476 if (TestBit(kAutoDelete)) {
1477 br->SetAutoDelete(true);
1478 }
1479 }
1480 }
1481 }
1482 if (fPlayer) {
1484 }
1485 // Notify user if requested.
1486 if (fNotify) {
1487 if(!fNotify->Notify()) return -6;
1488 }
1489 }
1490 }
1491 return treeReadEntry;
1492 }
1493
1494 if (fExternalFriends) {
1496 external_fe->MarkUpdated();
1497 }
1498 }
1499
1500 // Delete the current tree and open the new tree.
1501 TTreeCache* tpf = nullptr;
1502 // Delete file unless the file owns this chain!
1503 // FIXME: The "unless" case here causes us to leak memory.
1504 if (fFile) {
1505 if (!fDirectory->GetList()->FindObject(this)) {
1506 if (fTree) {
1507 // (fFile != 0 && fTree == 0) can happen when
1508 // InvalidateCurrentTree is called (for example from
1509 // AddFriend). Having fTree === 0 is necessary in that
1510 // case because in some cases GetTree is used as a check
1511 // to see if a TTree is already loaded.
1512 // However, this prevent using the following to reuse
1513 // the TTreeCache object.
1515 if (tpf) {
1516 tpf->ResetCache();
1517 }
1518
1519 fFile->SetCacheRead(nullptr, fTree);
1520 // If the tree has clones, copy them into the chain
1521 // clone list so we can change their branch addresses
1522 // when necessary.
1523 //
1524 // This is to support the syntax:
1525 //
1526 // TTree* clone = chain->GetTree()->CloneTree(0);
1527 //
1528 // We need to call the invalidate exactly here, since
1529 // we no longer need the value of fTree and it is
1530 // about to be deleted.
1532 }
1533
1534 if (fCanDeleteRefs) {
1535 fFile->Close("R");
1536 }
1537 delete fFile;
1538 fFile = nullptr;
1539 } else {
1540 // If the tree has clones, copy them into the chain
1541 // clone list so we can change their branch addresses
1542 // when necessary.
1543 //
1544 // This is to support the syntax:
1545 //
1546 // TTree* clone = chain->GetTree()->CloneTree(0);
1547 //
1549 }
1550 }
1551
1553 if (!element) {
1554 if (treeReadEntry) {
1555 return -4;
1556 }
1557 // Last attempt, just in case all trees in the chain have 0 entries.
1558 element = (TChainElement*) fFiles->At(0);
1559 if (!element) {
1560 return -4;
1561 }
1562 }
1563
1564 // FIXME: We leak memory here, we've just lost the open file
1565 // if we did not delete it above.
1566 {
1568 const char *option = fGlobalRegistration ? "READ" : "READ_WITHOUT_GLOBALREGISTRATION";
1569 fFile = TFile::Open(element->GetTitle(), option);
1572 }
1573
1574 // ----- Begin of modifications by MvL
1575 Int_t returnCode = 0;
1576 if (!fFile || fFile->IsZombie()) {
1577 if (fFile) {
1578 delete fFile;
1579 fFile = nullptr;
1580 }
1581 // Note: We do *not* own fTree.
1582 fTree = nullptr;
1583 returnCode = -3;
1584 } else {
1585 if (fPerfStats)
1587
1588 // Note: We do *not* own fTree after this, the file does!
1589 fTree = dynamic_cast<TTree*>(fFile->Get(element->GetName()));
1590 if (!fTree) {
1591 // Now that we do not check during the addition, we need to check here!
1592 Error("LoadTree", "Cannot find tree with name %s in file %s", element->GetName(), element->GetTitle());
1593 delete fFile;
1594 fFile = nullptr;
1595 // We do not return yet so that 'fEntries' can be updated with the
1596 // sum of the entries of all the other trees.
1597 returnCode = -4;
1598 } else if (!fGlobalRegistration) {
1600 }
1601 // Propagate the IMT settings
1602 if (fTree) {
1604 }
1605 }
1606
1608 // FIXME: We own fFile, we must be careful giving away a pointer to it!
1609 // FIXME: We may set fDirectory to zero here!
1610 fDirectory = fFile;
1611
1612 // Reuse cache from previous file (if any).
1613 if (tpf) {
1614 if (fFile) {
1615 // FIXME: fTree may be zero here.
1616 tpf->UpdateBranches(fTree);
1617 tpf->ResetCache();
1619 } else {
1620 // FIXME: One of the file in the chain is missing
1621 // we have no place to hold the pointer to the
1622 // TTreeCache.
1623 delete tpf;
1624 tpf = nullptr;
1625 }
1626 } else {
1627 if (fCacheUserSet) {
1628 this->SetCacheSize(fCacheSize);
1629 }
1630 }
1631
1632 // Check if fTreeOffset has really been set.
1633 Long64_t nentries = 0;
1634 if (fTree) {
1636 }
1637
1641 element->SetNumberEntries(nentries);
1642 // Below we must test >= in case the tree has no entries.
1643 if (entry >= fTreeOffset[fTreeNumber+1]) {
1644 if ((fTreeNumber < (fNtrees - 1)) && (entry < fTreeOffset[fTreeNumber+2])) {
1645 // The request entry is not in the tree 'fTreeNumber' we will need
1646 // to look further.
1647
1648 // Before moving on, let's record the result.
1649 element->SetLoadResult(returnCode);
1650
1651 // Before trying to read the file file/tree, notify the user
1652 // that we have switched trees if requested; the user might need
1653 // to properly account for the number of files/trees even if they
1654 // have no entries.
1655 if (fNotify) {
1656 if(!fNotify->Notify()) return -6;
1657 }
1658
1659 // Load the next TTree.
1660 return LoadTree(entry);
1661 } else {
1663 }
1664 }
1665 }
1666
1667
1668 if (!fTree) {
1669 // The Error message already issued. However if we reach here
1670 // we need to make sure that we do not use fTree.
1671 //
1672 // Force a reload of the tree next time.
1673 fTreeNumber = -1;
1674
1675 element->SetLoadResult(returnCode);
1676 return returnCode;
1677 }
1678 // ----- End of modifications by MvL
1679
1680 // Copy the chain's clone list into the new tree's
1681 // clone list so that branch addresses stay synchronized.
1682 if (fClones) {
1683 for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
1684 TTree* clone = (TTree*) lnk->GetObject();
1685 ((TChain*) fTree)->TTree::AddClone(clone);
1686 }
1687 }
1688
1689 // Since some of the friends of this chain might simple trees
1690 // (i.e., not really chains at all), we need to execute this
1691 // before calling LoadTree(entry) on the friends (so that
1692 // they use the correct read entry number).
1693
1694 // Change the new current tree to the new entry.
1696 if (loadResult == treeReadEntry) {
1697 element->SetLoadResult(0);
1698 } else {
1699 // This is likely to be an internal error, if treeReadEntry was not in range
1700 // (or intentionally -2 for TChain::GetEntries) then something happened
1701 // that is very odd/surprising.
1702 element->SetLoadResult(-5);
1703 }
1704
1705
1706 // Change the chain friends to the new entry.
1707 if (fFriends) {
1708 // An alternative would move this code to each of the function
1709 // calling LoadTree (and to overload a few more).
1710 TIter next(fFriends);
1711 TFriendLock lock(this, kLoadTree);
1712 TFriendElement* fe = nullptr;
1713 while ((fe = (TFriendElement*) next())) {
1714 TTree* t = fe->GetTree();
1715 if (!t) continue;
1716 if (t->GetTreeIndex()) {
1717 t->GetTreeIndex()->UpdateFormulaLeaves(GetTree());
1718 }
1719 if (t->GetTree() && t->GetTree()->GetTreeIndex()) {
1720 t->GetTree()->GetTreeIndex()->UpdateFormulaLeaves(GetTree());
1721 }
1722 if (treeReadEntry == -2) {
1723 // an entry after the end of the chain was requested (it usually happens when GetEntries is called)
1724 t->LoadTree(entry);
1725 } else {
1726 t->LoadTreeFriend(entry, this);
1727 }
1728 TTree* friend_t = t->GetTree();
1729 if (friend_t) {
1730 auto localfe = fTree->AddFriend(t, fe->GetName());
1732 }
1733 }
1734 }
1735
1738
1741
1742 // Set the branch statuses for the newly opened file.
1743 TIter next(fStatus);
1744 while ((element = (TChainElement*) next())) {
1745 Int_t status = element->GetStatus();
1746 fTree->SetBranchStatus(element->GetName(), status);
1747 }
1748
1749 // Set the branch addresses for the newly opened file.
1750 next.Reset();
1751 while ((element = (TChainElement*) next())) {
1752 void* addr = element->GetBaddress();
1753 if (addr) {
1754 TBranch* br = fTree->GetBranch(element->GetName());
1755 TBranch** pp = element->GetBranchPtr();
1756 if (pp) {
1757 // FIXME: What if br is zero here?
1758 *pp = br;
1759 }
1760 if (br) {
1761 if (!element->GetCheckedType()) {
1762 Int_t res = CheckBranchAddressType(br, TClass::GetClass(element->GetBaddressClassName()),
1763 (EDataType) element->GetBaddressType(), element->GetBaddressIsPtr());
1764 if ((res & kNeedEnableDecomposedObj) && !br->GetMakeClass()) {
1765 br->SetMakeClass(true);
1766 }
1767 element->SetDecomposedObj(br->GetMakeClass());
1768 element->SetCheckedType(true);
1769 }
1770 // FIXME: We may have to tell the branch it should
1771 // not be an owner of the object pointed at.
1772 br->SetAddress(addr);
1773 if (TestBit(kAutoDelete)) {
1774 br->SetAutoDelete(true);
1775 }
1776 }
1777 }
1778 }
1779
1780 // Update the addresses of the chain's cloned trees, if any.
1781 if (fClones) {
1782 for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
1783 TTree* clone = (TTree*) lnk->GetObject();
1784 CopyAddresses(clone);
1785 }
1786 }
1787
1788 // Update list of leaves in all TTreeFormula's of the TTreePlayer (if any).
1789 if (fPlayer) {
1791 }
1792
1793 // Notify user we have switched trees if requested.
1794 if (fNotify) {
1795 if(!fNotify->Notify()) return -6;
1796 }
1797
1798 // Return the new local entry number.
1799 return treeReadEntry;
1800}
1801
1802////////////////////////////////////////////////////////////////////////////////
1803/// Check / locate the files in the chain.
1804/// By default only the files not yet looked up are checked.
1805/// Use force = true to check / re-check every file.
1806
1808{
1809 TIter next(fFiles);
1810 TChainElement* element = nullptr;
1812 printf("\n");
1813 printf("TChain::Lookup - Looking up %d files .... \n", nelements);
1814 Int_t nlook = 0;
1815 TFileStager *stg = nullptr;
1816 while ((element = (TChainElement*) next())) {
1817 // Do not do it more than needed
1818 if (element->HasBeenLookedUp() && !force) continue;
1819 // Count
1820 nlook++;
1821 // Get the Url
1822 TUrl elemurl(element->GetTitle(), true);
1823 // Save current options and anchor
1824 TString anchor = elemurl.GetAnchor();
1825 TString options = elemurl.GetOptions();
1826 // Reset options and anchor
1827 elemurl.SetOptions("");
1828 elemurl.SetAnchor("");
1829 // Locate the file
1830 TString eurl(elemurl.GetUrl());
1831 if (!stg || !stg->Matches(eurl)) {
1832 SafeDelete(stg);
1833 {
1836 }
1837 if (!stg) {
1838 Error("Lookup", "TFileStager instance cannot be instantiated");
1839 break;
1840 }
1841 }
1842 Int_t n1 = (nelements > 100) ? (Int_t) nelements / 100 : 1;
1843 if (stg->Locate(eurl.Data(), eurl) == 0) {
1844 if (nlook > 0 && !(nlook % n1)) {
1845 printf("Lookup | %3d %% finished\r", 100 * nlook / nelements);
1846 fflush(stdout);
1847 }
1848 // Get the effective end-point Url
1849 elemurl.SetUrl(eurl);
1850 // Restore original options and anchor, if any
1851 elemurl.SetOptions(options);
1852 elemurl.SetAnchor(anchor);
1853 // Save it into the element
1854 element->SetTitle(elemurl.GetUrl());
1855 // Remember
1856 element->SetLookedUp();
1857 } else {
1858 // Failure: remove
1861 Error("Lookup", "file %s does not exist\n", eurl.Data());
1862 else
1863 Error("Lookup", "file %s cannot be read\n", eurl.Data());
1864 }
1865 }
1866 if (nelements > 0)
1867 printf("Lookup | %3d %% finished\n", 100 * nlook / nelements);
1868 else
1869 printf("\n");
1870 fflush(stdout);
1871 SafeDelete(stg);
1872}
1873
1874////////////////////////////////////////////////////////////////////////////////
1875/// Loop on nentries of this chain starting at firstentry. (NOT IMPLEMENTED)
1876
1878{
1879 Error("Loop", "Function not yet implemented");
1880
1881 if (option || nentries || firstentry) { } // keep warnings away
1882
1883#if 0
1884 if (LoadTree(firstentry) < 0) return;
1885
1886 if (firstentry < 0) firstentry = 0;
1888 if (lastentry > fEntries-1) {
1889 lastentry = fEntries -1;
1890 }
1891
1892 GetPlayer();
1893 GetSelector();
1894 fSelector->Start(option);
1895
1897 Int_t tree,e0,en;
1898 for (tree=0;tree<fNtrees;tree++) {
1899 e0 = fTreeOffset[tree];
1900 en = fTreeOffset[tree+1] - 1;
1901 if (en > lastentry) en = lastentry;
1902 if (entry > en) continue;
1903
1904 LoadTree(entry);
1905 fSelector->BeginFile();
1906
1907 while (entry <= en) {
1909 entry++;
1910 }
1911 fSelector->EndFile();
1912 }
1913
1914 fSelector->Finish(option);
1915#endif
1916}
1917
1918////////////////////////////////////////////////////////////////////////////////
1919/// List the chain.
1920
1922{
1924 TIter next(fFiles);
1925 TChainElement* file = nullptr;
1927 while ((file = (TChainElement*)next())) {
1928 file->ls(option);
1929 }
1931}
1932
1933////////////////////////////////////////////////////////////////////////////////
1934/// Merge all the entries in the chain into a new tree in a new file.
1935///
1936/// See important note in the following function Merge().
1937///
1938/// If the chain is expecting the input tree inside a directory,
1939/// this directory is NOT created by this routine.
1940///
1941/// So in a case where we have:
1942/// ~~~ {.cpp}
1943/// TChain ch("mydir/mytree");
1944/// ch.Merge("newfile.root");
1945/// ~~~
1946/// The resulting file will have not subdirectory. To recreate
1947/// the directory structure do:
1948/// ~~~ {.cpp}
1949/// TFile* file = TFile::Open("newfile.root", "RECREATE");
1950/// file->mkdir("mydir")->cd();
1951/// ch.Merge(file, 0);
1952/// ~~~
1953
1955{
1956 TFile *file = TFile::Open(name, "recreate", "chain files", 1);
1957 return Merge(file, 0, option);
1958}
1959
1960////////////////////////////////////////////////////////////////////////////////
1961/// Merge all chains in the collection. (NOT IMPLEMENTED)
1962
1963Long64_t TChain::Merge(TCollection* /* list */, Option_t* /* option */ )
1964{
1965 Error("Merge", "not implemented");
1966 return -1;
1967}
1968
1969////////////////////////////////////////////////////////////////////////////////
1970/// Merge all chains in the collection. (NOT IMPLEMENTED)
1971
1973{
1974 Error("Merge", "not implemented");
1975 return -1;
1976}
1977
1978////////////////////////////////////////////////////////////////////////////////
1979/// Merge all the entries in the chain into a new tree in the current file.
1980///
1981/// Note: The "file" parameter is *not* the file where the new
1982/// tree will be inserted. The new tree is inserted into
1983/// gDirectory, which is usually the most recently opened
1984/// file, or the directory most recently cd()'d to.
1985///
1986/// If option = "C" is given, the compression level for all branches
1987/// in the new Tree is set to the file compression level. By default,
1988/// the compression level of all branches is the original compression
1989/// level in the old trees.
1990///
1991/// If basketsize > 1000, the basket size for all branches of the
1992/// new tree will be set to basketsize.
1993///
1994/// Example using the file generated in $ROOTSYS/test/Event
1995/// merge two copies of Event.root
1996/// ~~~ {.cpp}
1997/// gSystem.Load("libEvent");
1998/// TChain ch("T");
1999/// ch.Add("Event1.root");
2000/// ch.Add("Event2.root");
2001/// ch.Merge("all.root");
2002/// ~~~
2003/// If the chain is expecting the input tree inside a directory,
2004/// this directory is NOT created by this routine.
2005///
2006/// So if you do:
2007/// ~~~ {.cpp}
2008/// TChain ch("mydir/mytree");
2009/// ch.Merge("newfile.root");
2010/// ~~~
2011/// The resulting file will not have subdirectories. In order to
2012/// preserve the directory structure do the following instead:
2013/// ~~~ {.cpp}
2014/// TFile* file = TFile::Open("newfile.root", "RECREATE");
2015/// file->mkdir("mydir")->cd();
2016/// ch.Merge(file, 0);
2017/// ~~~
2018/// If 'option' contains the word 'fast' the merge will be done without
2019/// unzipping or unstreaming the baskets (i.e., a direct copy of the raw
2020/// bytes on disk).
2021///
2022/// When 'fast' is specified, 'option' can also contains a
2023/// sorting order for the baskets in the output file.
2024///
2025/// There is currently 3 supported sorting order:
2026/// ~~~ {.cpp}
2027/// SortBasketsByOffset (the default)
2028/// SortBasketsByBranch
2029/// SortBasketsByEntry
2030/// ~~~
2031/// When using SortBasketsByOffset the baskets are written in
2032/// the output file in the same order as in the original file
2033/// (i.e. the basket are sorted on their offset in the original
2034/// file; Usually this also means that the baskets are sorted
2035/// on the index/number of the _last_ entry they contain)
2036///
2037/// When using SortBasketsByBranch all the baskets of each
2038/// individual branches are stored contiguously. This tends to
2039/// optimize reading speed when reading a small number (1->5) of
2040/// branches, since all their baskets will be clustered together
2041/// instead of being spread across the file. However it might
2042/// decrease the performance when reading more branches (or the full
2043/// entry).
2044///
2045/// When using SortBasketsByEntry the baskets with the lowest
2046/// starting entry are written first. (i.e. the baskets are
2047/// sorted on the index/number of the first entry they contain).
2048/// This means that on the file the baskets will be in the order
2049/// in which they will be needed when reading the whole tree
2050/// sequentially.
2051///
2052/// ## IMPORTANT Note 1: AUTOMATIC FILE OVERFLOW
2053///
2054/// When merging many files, it may happen that the resulting file
2055/// reaches a size > TTree::fgMaxTreeSize (default = 100 GBytes).
2056/// In this case the current file is automatically closed and a new
2057/// file started. If the name of the merged file was "merged.root",
2058/// the subsequent files will be named "merged_1.root", "merged_2.root",
2059/// etc. fgMaxTreeSize may be modified via the static function
2060/// TTree::SetMaxTreeSize.
2061/// When in fast mode, the check and switch is only done in between each
2062/// input file.
2063///
2064/// ## IMPORTANT Note 2: The output file is automatically closed and deleted.
2065///
2066/// This is required because in general the automatic file overflow described
2067/// above may happen during the merge.
2068/// If only the current file is produced (the file passed as first argument),
2069/// one can instruct Merge to not close and delete the file by specifying
2070/// the option "keep".
2071///
2072/// The function returns the total number of files produced.
2073/// To check that all files have been merged use something like:
2074/// ~~~ {.cpp}
2075/// if (newchain->GetEntries()!=oldchain->GetEntries()) {
2076/// ... not all the file have been copied ...
2077/// }
2078/// ~~~
2079
2081{
2082 // We must have been passed a file, we will use it
2083 // later to reset the compression level of the branches.
2084 if (!file) {
2085 // FIXME: We need an error message here.
2086 return 0;
2087 }
2088
2089 // Options
2090 bool fastClone = false;
2091 TString opt = option;
2092 opt.ToLower();
2093 if (opt.Contains("fast")) {
2094 fastClone = true;
2095 }
2096
2097 // The chain tree must have a list of branches
2098 // because we may try to change their basket
2099 // size later.
2101 if (!lbranches) {
2102 // FIXME: We need an error message here.
2103 return 0;
2104 }
2105
2106 // The chain must have a current tree because
2107 // that is the one we will clone.
2108 if (!fTree) {
2109 // -- LoadTree() has not yet been called, no current tree.
2110 // FIXME: We need an error message here.
2111 return 0;
2112 }
2113
2114 // Copy the chain's current tree without
2115 // copying any entries, we will do that later.
2116 TTree* newTree = CloneTree(0);
2117 if (!newTree) {
2118 // FIXME: We need an error message here.
2119 return 0;
2120 }
2121
2122 // Strip out the (potential) directory name.
2123 // FIXME: The merged chain may or may not have the
2124 // same name as the original chain. This is
2125 // bad because the chain name determines the
2126 // names of the trees in the chain by default.
2127 newTree->SetName(gSystem->BaseName(GetName()));
2128
2129 // FIXME: Why do we do this?
2130 newTree->SetAutoSave(2000000000);
2131
2132 // Circularity is incompatible with merging, it may
2133 // force us to throw away entries, which is not what
2134 // we are supposed to do.
2135 newTree->SetCircular(0);
2136
2137 // Reset the compression level of the branches.
2138 if (opt.Contains("c")) {
2139 TBranch* branch = nullptr;
2140 TIter nextb(newTree->GetListOfBranches());
2141 while ((branch = (TBranch*) nextb())) {
2142 branch->SetCompressionSettings(file->GetCompressionSettings());
2143 }
2144 }
2145
2146 // Reset the basket size of the branches.
2147 if (basketsize > 1000) {
2148 TBranch* branch = nullptr;
2149 TIter nextb(newTree->GetListOfBranches());
2150 while ((branch = (TBranch*) nextb())) {
2151 branch->SetBasketSize(basketsize);
2152 }
2153 }
2154
2155 // Copy the entries.
2156 if (fastClone) {
2157 if ( newTree->CopyEntries( this, -1, option ) < 0 ) {
2158 // There was a problem!
2159 Error("Merge", "TTree has not been cloned\n");
2160 }
2161 } else {
2162 newTree->CopyEntries( this, -1, option );
2163 }
2164
2165 // Write the new tree header.
2166 newTree->Write();
2167
2168 // Get our return value.
2169 Int_t nfiles = newTree->GetFileNumber() + 1;
2170
2171 // Close and delete the current file of the new tree.
2172 if (!opt.Contains("keep")) {
2173 // Delete the currentFile and the TTree object.
2174 delete newTree->GetCurrentFile();
2175 }
2176 return nfiles;
2177}
2178
2179////////////////////////////////////////////////////////////////////////////////
2180/// Get the tree url or filename and other information from the name
2181///
2182/// A treename and a url's query section is split off from name. The
2183/// splitting depends on whether the resulting filename is to be
2184/// subsequently treated for wildcards or not, since the question mark is
2185/// both the url query identifier and a wildcard. Wildcard matching is not
2186/// done in this method itself.
2187/// ~~~ {.cpp}
2188/// [xxx://host]/a/path/file_name[?query[#treename]]
2189/// ~~~
2190///
2191/// The following way to specify the treename is still supported with the
2192/// constrain that the file name contains the sub-string '.root'.
2193/// This is now deprecated and will be removed in future versions.
2194/// ~~~ {.cpp}
2195/// [xxx://host]/a/path/file.root[.oext][/treename]
2196/// [xxx://host]/a/path/file.root[.oext][/treename][?query]
2197/// ~~~
2198///
2199/// Note that in a case like this
2200/// ~~~ {.cpp}
2201/// [xxx://host]/a/path/file#treename
2202/// ~~~
2203/// i.e. anchor but no options (query), the filename will be the full path, as
2204/// the anchor may be the internal file name of an archive. Use '?#%treename' to
2205/// pass the treename if the query field is empty.
2206///
2207/// \param[in] name is the original name
2208/// \param[out] filename the url or filename to be opened or matched
2209/// \param[out] treename the treename, which may be found in a url fragment section
2210/// as a trailing part of the name (deprecated).
2211/// If not found this will be empty.
2212/// Exception: a fragment containing the '=' character is _not_
2213/// interpreted as a treename
2214/// \param[out] query is the url query section, including the leading question
2215/// mark. If not found or the query section is only followed by
2216/// a fragment this will be empty.
2217/// \param[out] suffix the portion of name which was removed to from filename.
2218
2220 TString &suffix) const
2221{
2222 Ssiz_t pIdx = kNPOS;
2223 filename.Clear();
2224 treename.Clear();
2225 query.Clear();
2226 suffix.Clear();
2227
2228 // General case
2229 TUrl url(name, true);
2230 filename = (strcmp(url.GetProtocol(), "file")) ? url.GetUrl() : url.GetFileAndOptions();
2231
2232 TString fn = url.GetFile();
2233 // Extract query, if any
2234 if (url.GetOptions() && (strlen(url.GetOptions()) > 0))
2235 query.Form("?%s", url.GetOptions());
2236 // The treename can be passed as anchor
2237 const char *anchor = url.GetAnchor();
2238 if (anchor && anchor[0] != '\0') {
2239 // Support "?#tree_name" and "?query#tree_name"
2240 // "#tree_name" (no '?' is for tar archives)
2241 // If the treename would contain a '=', treat the anchor as part of the query instead. This makes sure
2242 // that Davix parameters are passed.
2243 if (!query.IsNull() || strstr(name, "?#")) {
2244 if (strstr(anchor, "=")) {
2245 query.Append("#");
2246 query.Append(anchor);
2247 } else {
2248 treename = anchor;
2249 }
2250 } else {
2251 // The anchor is part of the file name
2252 fn = url.GetFileAndOptions();
2253 }
2254 }
2255 // Suffix
2256 suffix = url.GetFileAndOptions();
2257 // Get options from suffix by removing the file name
2258 suffix.Replace(suffix.Index(fn), fn.Length(), "");
2259 // Remove the options suffix from the original file name
2260 filename.Replace(filename.Index(suffix), suffix.Length(), "");
2261
2262 // Special case: [...]file.root/treename
2263 static const char *dotr = ".root";
2264 static Ssiz_t dotrl = strlen(dotr);
2265 // Find the last one
2266 Ssiz_t js = filename.Index(dotr);
2267 while (js != kNPOS) {
2268 pIdx = js;
2269 js = filename.Index(dotr, js + 1);
2270 }
2271 if (pIdx != kNPOS) {
2272 static const char *slash = "/";
2273 static Ssiz_t slashl = strlen(slash);
2274 // Find the last one
2275 Ssiz_t ppIdx = filename.Index(slash, pIdx + dotrl);
2276 if (ppIdx != kNPOS) {
2277 // Good treename with the old recipe
2278 treename = filename(ppIdx + slashl, filename.Length());
2279 filename.Remove(ppIdx + slashl - 1);
2280 suffix.Insert(0, TString::Format("/%s", treename.Data()));
2281 }
2282 }
2283}
2284
2285////////////////////////////////////////////////////////////////////////////////
2286/// Print the header information of each tree in the chain.
2287/// See TTree::Print for a list of options.
2288
2290{
2291 TIter next(fFiles);
2293 while ((element = (TChainElement*)next())) {
2294 Printf("******************************************************************************");
2295 Printf("*Chain :%-10s: %-54s *", GetName(), element->GetTitle());
2296 Printf("******************************************************************************");
2297 TFile *file = TFile::Open(element->GetTitle());
2298 if (file && !file->IsZombie()) {
2299 TTree *tree = (TTree*)file->Get(element->GetName());
2300 if (tree) tree->Print(option);
2301 }
2302 delete file;
2303 }
2304}
2305
2306////////////////////////////////////////////////////////////////////////////////
2307/// Process all entries in this chain, calling functions in filename.
2308/// The return value is -1 in case of error and TSelector::GetStatus() in
2309/// in case of success.
2310/// See TTree::Process.
2311
2313{
2314 if (fProofChain) {
2315 // Make sure the element list is up to date
2316 if (!TestBit(kProofUptodate))
2317 SetProof(true, true);
2321 }
2322
2323 if (LoadTree(firstentry) < 0) {
2324 return 0;
2325 }
2327}
2328
2329////////////////////////////////////////////////////////////////////////////////
2330/// Process this chain executing the code in selector.
2331/// The return value is -1 in case of error and TSelector::GetStatus() in
2332/// in case of success.
2333
2335{
2336 if (fProofChain) {
2337 // Make sure the element list is up to date
2338 if (!TestBit(kProofUptodate))
2339 SetProof(true, true);
2342 return fProofChain->Process(selector, option, nentries, firstentry);
2343 }
2344
2345 return TTree::Process(selector, option, nentries, firstentry);
2346}
2347
2348////////////////////////////////////////////////////////////////////////////////
2349/// Make sure that obj (which is being deleted or will soon be) is no
2350/// longer referenced by this TTree.
2351
2353{
2354 if (fFile == obj) {
2355 fFile = nullptr;
2356 fDirectory = nullptr;
2357 fTree = nullptr;
2358 }
2359 if (fDirectory == obj) {
2360 fDirectory = nullptr;
2361 fTree = nullptr;
2362 }
2363 if (fTree == obj) {
2364 fTree = nullptr;
2365 }
2366}
2367
2368////////////////////////////////////////////////////////////////////////////////
2369/// Remove a friend from the list of friends.
2370
2372{
2373 // We already have been visited while recursively looking
2374 // through the friends tree, let return
2375
2376 if (!fFriends) {
2377 return;
2378 }
2379
2381
2382 if (fProofChain)
2383 // This updates the proxy chain when we will really use PROOF
2385
2386 // We need to invalidate the loading of the current tree because its list
2387 // of real friends is now obsolete. It is repairable only from LoadTree.
2389}
2390
2391////////////////////////////////////////////////////////////////////////////////
2392/// Resets the state of this chain.
2393
2395{
2396 delete fFile;
2397 fFile = nullptr;
2398 fNtrees = 0;
2399 fTreeNumber = -1;
2400 fTree = nullptr;
2401 fFile = nullptr;
2402 fFiles->Delete();
2403 fStatus->Delete();
2404 fTreeOffset[0] = 0;
2405 TChainElement* element = new TChainElement("*", "");
2407 fDirectory = nullptr;
2408
2409 TTree::Reset();
2410}
2411
2412////////////////////////////////////////////////////////////////////////////////
2413/// Resets the state of this chain after a merge (keep the customization but
2414/// forget the data).
2415
2417{
2418 fNtrees = 0;
2419 fTreeNumber = -1;
2420 fTree = nullptr;
2421 fFile = nullptr;
2422 fFiles->Delete();
2423 fTreeOffset[0] = 0;
2424
2426}
2427
2428////////////////////////////////////////////////////////////////////////////////
2429/// Save TChain as a C++ statements on output stream out.
2430/// With the option "friend" save the description of all the
2431/// TChain's friend trees or chains as well.
2432
2433void TChain::SavePrimitive(std::ostream &out, Option_t *option)
2434{
2435 static Int_t chCounter = 0;
2436
2437 TString chName = gInterpreter->MapCppName(GetName());
2438 if (chName.IsNull())
2439 chName = "_chain";
2440 ++chCounter;
2441 chName += chCounter;
2442
2443 TString opt = option;
2444 opt.ToLower();
2445
2446 out << " TChain *" << chName.Data() << " = new TChain(\"" << GetName() << "\");" << std::endl;
2447
2448 if (opt.Contains("friend")) {
2449 opt.ReplaceAll("friend", "");
2450 for (TObject *frel : *fFriends) {
2451 TTree *frtree = ((TFriendElement *)frel)->GetTree();
2452 if (dynamic_cast<TChain *>(frtree)) {
2453 if (strcmp(frtree->GetName(), GetName()) != 0)
2454 --chCounter; // make friends get the same chain counter
2455 frtree->SavePrimitive(out, opt.Data());
2456 out << " " << chName.Data() << "->AddFriend(\"" << frtree->GetName() << "\");" << std::endl;
2457 } else { // ordinary friend TTree
2458 TDirectory *file = frtree->GetDirectory();
2459 if (file && dynamic_cast<TFile *>(file))
2460 out << " " << chName.Data() << "->AddFriend(\"" << frtree->GetName() << "\", \"" << file->GetName()
2461 << "\");" << std::endl;
2462 }
2463 }
2464 }
2465 out << std::endl;
2466
2467 for (TObject *el : *fFiles) {
2469 // Save tree file if it is really loaded to the chain
2470 if (chel->GetLoadResult() == 0 && chel->GetEntries() != 0) {
2471 if (chel->GetEntries() == TTree::kMaxEntries) // tree number of entries is not yet known
2472 out << " " << chName.Data() << "->AddFile(\"" << chel->GetTitle() << "\");" << std::endl;
2473 else
2474 out << " " << chName.Data() << "->AddFile(\"" << chel->GetTitle() << "\"," << chel->GetEntries() << ");"
2475 << std::endl;
2476 }
2477 }
2478 out << std::endl;
2479
2480 SaveMarkerAttributes(out, chName.Data(), 1, 1, 1);
2481}
2482
2483////////////////////////////////////////////////////////////////////////////////
2484/// Loop on tree and print entries passing selection.
2485/// - If varexp is 0 (or "") then print only first 8 columns.
2486/// - If varexp = "*" print all columns.
2487/// - Otherwise a columns selection can be made using "var1:var2:var3".
2488/// See TTreePlayer::Scan for more information.
2489
2491{
2492 if (LoadTree(firstentry) < 0) {
2493 return 0;
2494 }
2496}
2497
2498////////////////////////////////////////////////////////////////////////////////
2499/// Set the global branch kAutoDelete bit.
2500///
2501/// When LoadTree loads a new Tree, the branches for which
2502/// the address is set will have the option AutoDelete set
2503/// For more details on AutoDelete, see TBranch::SetAutoDelete.
2504
2506{
2507 if (autodelete) {
2508 SetBit(kAutoDelete, true);
2509 } else {
2510 SetBit(kAutoDelete, false);
2511 }
2512}
2513
2515{
2516 // Set the cache size of the underlying TTree,
2517 // See TTree::SetCacheSize.
2518 // Returns 0 cache state ok (exists or not, as appropriate)
2519 // -1 on error
2520
2521 Int_t res = 0;
2522
2523 // remember user has requested this cache setting
2524 fCacheUserSet = true;
2525
2526 if (fTree) {
2527 res = fTree->SetCacheSize(cacheSize);
2528 } else {
2529 // If we don't have a TTree yet only record the cache size wanted
2530 res = 0;
2531 }
2532 fCacheSize = cacheSize; // Record requested size.
2533 return res;
2534}
2535
2536////////////////////////////////////////////////////////////////////////////////
2537/// Reset the addresses of the branch.
2538
2540{
2542 if (element) {
2543 element->SetBaddress(nullptr);
2544 }
2545 if (fTree) {
2547 }
2548}
2549
2550////////////////////////////////////////////////////////////////////////////////
2551/// Reset the addresses of the branches.
2552
2554{
2555 // We already have been visited while recursively looking
2556 // through the friends tree, let return
2558 return;
2559 }
2560 TIter next(fStatus);
2561 TChainElement* element = nullptr;
2562 while ((element = (TChainElement*) next())) {
2563 element->SetBaddress(nullptr);
2564 }
2565 if (fTree) {
2567 }
2568 if (fFriends) {
2571 auto *frTree = frEl->GetTree();
2572 if (frTree) {
2573 frTree->ResetBranchAddresses();
2574 }
2575 }
2576 }
2577}
2578
2579////////////////////////////////////////////////////////////////////////////////
2580/// Set branch address.
2581///
2582/// \param[in] bname is the name of a branch.
2583/// \param[in] add is the address of the branch.
2584/// \param[in] ptr
2585///
2586/// Note: See the comments in TBranchElement::SetAddress() for a more
2587/// detailed discussion of the meaning of the add parameter.
2588///
2589/// IMPORTANT REMARK:
2590///
2591/// In case TChain::SetBranchStatus is called, it must be called
2592/// BEFORE calling this function.
2593///
2594/// See TTree::CheckBranchAddressType for the semantic of the return value.
2595
2596Int_t TChain::SetBranchAddress(const char *bname, void* add, TBranch** ptr)
2597{
2598 Int_t res = kNoCheck;
2599
2600 // Check if bname is already in the status list.
2601 // If not, create a TChainElement object and set its address.
2603 if (!element) {
2604 element = new TChainElement(bname, "");
2606 }
2607 element->SetBaddress(add);
2608 element->SetBranchPtr(ptr);
2609 // Also set address in current tree.
2610 // FIXME: What about the chain clones?
2611 if (fTreeNumber >= 0) {
2612 TBranch* branch = fTree->GetBranch(bname);
2613 if (ptr) {
2614 *ptr = branch;
2615 }
2616 if (branch) {
2617 res = CheckBranchAddressType(branch, TClass::GetClass(element->GetBaddressClassName()), (EDataType) element->GetBaddressType(), element->GetBaddressIsPtr());
2618 if ((res & kNeedEnableDecomposedObj) && !branch->GetMakeClass()) {
2619 branch->SetMakeClass(true);
2620 }
2621 element->SetDecomposedObj(branch->GetMakeClass());
2622 element->SetCheckedType(true);
2623 if (fClones) {
2624 void* oldAdd = branch->GetAddress();
2625 for (TObjLink* lnk = fClones->FirstLink(); lnk; lnk = lnk->Next()) {
2626 TTree* clone = (TTree*) lnk->GetObject();
2627 TBranch* cloneBr = clone->GetBranch(bname);
2628 if (cloneBr && (cloneBr->GetAddress() == oldAdd)) {
2629 // the clone's branch is still pointing to us
2630 cloneBr->SetAddress(add);
2631 if ((res & kNeedEnableDecomposedObj) && !cloneBr->GetMakeClass()) {
2632 cloneBr->SetMakeClass(true);
2633 }
2634 }
2635 }
2636 }
2637
2638 branch->SetAddress(add);
2639 } else {
2640 Error("SetBranchAddress", "unknown branch -> %s", bname);
2641 return kMissingBranch;
2642 }
2643 } else {
2644 if (ptr) {
2645 *ptr = nullptr;
2646 }
2647 }
2648 return res;
2649}
2650
2651////////////////////////////////////////////////////////////////////////////////
2652/// Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
2653/// See TTree::CheckBranchAddressType for the semantic of the return value.
2654///
2655/// Note: See the comments in TBranchElement::SetAddress() for a more
2656/// detailed discussion of the meaning of the add parameter.
2657
2659{
2660 return SetBranchAddress(bname, add, nullptr, realClass, datatype, isptr);
2661}
2662
2663////////////////////////////////////////////////////////////////////////////////
2664/// Check if bname is already in the status list, and if not, create a TChainElement object and set its address.
2665/// See TTree::CheckBranchAddressType for the semantic of the return value.
2666///
2667/// Note: See the comments in TBranchElement::SetAddress() for a more
2668/// detailed discussion of the meaning of the add parameter.
2669
2670Int_t TChain::SetBranchAddress(const char* bname, void* add, TBranch** ptr, TClass* realClass, EDataType datatype, bool isptr)
2671{
2673 if (!element) {
2674 element = new TChainElement(bname, "");
2676 }
2677 if (realClass) {
2678 element->SetBaddressClassName(realClass->GetName());
2679 }
2680 element->SetBaddressType((UInt_t) datatype);
2681 element->SetBaddressIsPtr(isptr);
2682 element->SetBranchPtr(ptr);
2683 return SetBranchAddress(bname, add, ptr);
2684}
2685
2686////////////////////////////////////////////////////////////////////////////////
2687/// Set branch status to Process or DoNotProcess
2688///
2689/// \param[in] bname is the name of a branch. if bname="*", apply to all branches.
2690/// \param[in] status = 1 branch will be processed,
2691/// = 0 branch will not be processed
2692/// \param[out] found
2693///
2694/// See IMPORTANT REMARKS in TTree::SetBranchStatus and TChain::SetBranchAddress
2695///
2696/// If found is not 0, the number of branch(es) found matching the regular
2697/// expression is returned in *found AND the error message 'unknown branch'
2698/// is suppressed.
2699
2700void TChain::SetBranchStatus(const char* bname, bool status, UInt_t* found)
2701{
2702 // FIXME: We never explicitly set found to zero!
2703
2704 // Check if bname is already in the status list,
2705 // if not create a TChainElement object and set its status.
2707 if (element) {
2709 } else {
2710 element = new TChainElement(bname, "");
2711 }
2713 element->SetStatus(status);
2714 // Also set status in current tree.
2715 if (fTreeNumber >= 0) {
2716 fTree->SetBranchStatus(bname, status, found);
2717 } else if (found) {
2718 *found = 1;
2719 }
2720}
2721
2722////////////////////////////////////////////////////////////////////////////////
2723/// Remove reference to this chain from current directory and add
2724/// reference to new directory dir. dir can be 0 in which case the chain
2725/// does not belong to any directory.
2726
2728{
2729 if (fDirectory == dir) return;
2730 if (fDirectory) fDirectory->Remove(this);
2731 fDirectory = dir;
2732 if (fDirectory) {
2733 fDirectory->Append(this);
2735 } else {
2736 fFile = nullptr;
2737 }
2738}
2739
2740////////////////////////////////////////////////////////////////////////////////
2741/// \brief Set the input entry list (processing the entries of the chain will
2742/// then be limited to the entries in the list).
2743///
2744/// \param[in] elist The entry list to be assigned to this chain.
2745/// \param[in] opt An option string. Possible values are:
2746/// - "" (default): both the file names of the chain elements and the file
2747/// names of the TEntryList sublists are expanded to full path name.
2748/// - "ne": the file names are taken as they are and not expanded
2749/// - "sync": the TChain will go through the TEntryList in lockstep with the
2750/// trees in the chain rather than performing a lookup based on
2751/// treename and filename. This is mostly useful when the TEntryList
2752/// has multiple sublists for the same tree and filename.
2753/// \throws std::runtime_error If option "sync" was chosen and either:
2754/// - \p elist doesn't have sub entry lists.
2755/// - the number of sub entry lists in \p elist is different than the
2756/// number of trees in the chain.
2757/// - any of the sub entry lists in \p elist doesn't correspond to the
2758/// tree of the chain with the same index (i.e. it doesn't share the
2759/// same tree name and file name).
2760///
2761/// This function finds correspondence between the sub-lists of the TEntryList
2762/// and the trees of the TChain.
2763
2765{
2766 if (fEntryList){
2767 //check, if the chain is the owner of the previous entry list
2768 //(it happens, if the previous entry list was created from a user-defined
2769 //TEventList in SetEventList() function)
2771 TEntryList *tmp = fEntryList;
2772 fEntryList = nullptr; // Avoid problem with RecursiveRemove.
2773 delete tmp;
2774 } else {
2775 fEntryList = nullptr;
2776 }
2777 }
2778 if (!elist){
2779 fEntryList = nullptr;
2780 fEventList = nullptr;
2781 return;
2782 }
2783 if (!elist->TestBit(kCanDelete)){
2784 //this is a direct call to SetEntryList, not via SetEventList
2785 fEventList = nullptr;
2786 }
2787 if (elist->GetN() == 0){
2788 fEntryList = elist;
2789 return;
2790 }
2791 if (fProofChain){
2792 //for processing on proof, event list and entry list can't be
2793 //set at the same time.
2794 fEventList = nullptr;
2795 fEntryList = elist;
2796 return;
2797 }
2798
2800 Int_t listfound=0;
2802
2803 TEntryList *templist = nullptr;
2804
2805 const auto *subentrylists = elist->GetLists();
2806 if(strcmp(opt, "sync") == 0){
2807 if(!subentrylists){
2808 std::string msg{"In 'TChain::SetEntryList': "};
2809 msg += "the input TEntryList doesn't have sub entry lists. Please make sure too add them through ";
2810 msg += "TEntryList::AddSubList";
2811 throw std::runtime_error(msg);
2812 }
2813 const auto nsubelists = subentrylists->GetEntries();
2814 if(nsubelists != ne){
2815 std::string msg{"In 'TChain::SetEntryList': "};
2816 msg += "the number of sub entry lists in the input TEntryList (";
2817 msg += std::to_string(nsubelists);
2818 msg += ") is not equal to the number of files in the chain (";
2819 msg += std::to_string(ne);
2820 msg += ")";
2821 throw std::runtime_error(msg);
2822 }
2823 }
2824
2825 for (Int_t ie = 0; ie<ne; ie++){
2827 treename = chainElement->GetName();
2828 filename = chainElement->GetTitle();
2829
2830 if(strcmp(opt, "sync") == 0){
2831 // If the user asked for "sync" option, there should be a 1:1 mapping
2832 // between trees in the chain and sub entry lists in the argument elist
2833 // We have already checked that the input TEntryList has a number of
2834 // sub entry lists equal to the number of files in the chain.
2835 templist = static_cast<TEntryList*>(subentrylists->At(ie));
2836 auto elisttreename = templist->GetTreeName();
2837 auto elistfilename = templist->GetFileName();
2838
2840 std::string msg{"In 'TChain::SetEntryList': "};
2841 msg += "the sub entry list at index ";
2842 msg += std::to_string(ie);
2843 msg += " doesn't correspond to treename '";
2844 msg += treename;
2845 msg += "' and filename '";
2846 msg += filename;
2847 msg += "': it has treename '";
2848 msg += elisttreename;
2849 msg += "' and filename '";
2850 msg += elistfilename;
2851 msg += "'";
2852 throw std::runtime_error(msg);
2853 }
2854
2855 }else{
2856 templist = elist->GetEntryList(treename, filename, opt);
2857 }
2858
2859 if (templist) {
2860 listfound++;
2861 templist->SetTreeNumber(ie);
2862 }
2863 }
2864
2865 if (listfound == 0){
2866 Error("SetEntryList", "No list found for the trees in this chain");
2867 fEntryList = nullptr;
2868 return;
2869 }
2870 fEntryList = elist;
2871 TList *elists = elist->GetLists();
2872 bool shift = false;
2873 TIter next(elists);
2874
2875 //check, if there are sub-lists in the entry list, that don't
2876 //correspond to any trees in the chain
2877 while((templist = (TEntryList*)next())){
2878 if (templist->GetTreeNumber() < 0){
2879 shift = true;
2880 break;
2881 }
2882 }
2883 fEntryList->SetShift(shift);
2884
2885}
2886
2887////////////////////////////////////////////////////////////////////////////////
2888/// Set the input entry list (processing the entries of the chain will then be
2889/// limited to the entries in the list). This function creates a special kind
2890/// of entry list (TEntryListFromFile object) that loads lists, corresponding
2891/// to the chain elements, one by one, so that only one list is in memory at a time.
2892///
2893/// If there is an error opening one of the files, this file is skipped and the
2894/// next file is loaded
2895///
2896/// File naming convention:
2897///
2898/// - by default, filename_elist.root is used, where filename is the
2899/// name of the chain element
2900/// - xxx$xxx.root - $ sign is replaced by the name of the chain element
2901///
2902/// If the list name is not specified (by passing filename_elist.root/listname to
2903/// the TChain::SetEntryList() function, the first object of class TEntryList
2904/// in the file is taken.
2905///
2906/// It is assumed, that there are as many list files, as there are elements in
2907/// the chain and they are in the same order
2908
2909void TChain::SetEntryListFile(const char *filename, Option_t * /*opt*/)
2910{
2911
2912 if (fEntryList){
2913 //check, if the chain is the owner of the previous entry list
2914 //(it happens, if the previous entry list was created from a user-defined
2915 //TEventList in SetEventList() function)
2917 TEntryList *tmp = fEntryList;
2918 fEntryList = nullptr; // Avoid problem with RecursiveRemove.
2919 delete tmp;
2920 } else {
2921 fEntryList = nullptr;
2922 }
2923 }
2924
2925 fEventList = nullptr;
2926
2928
2929 Int_t dotslashpos = basename.Index(".root/");
2931 if (dotslashpos>=0) {
2932 // Copy the list name specification
2934 // and remove it from basename
2935 basename.Remove(dotslashpos+5);
2936 }
2939 fEntryList->SetDirectory(nullptr);
2940 ((TEntryListFromFile*)fEntryList)->SetFileNames(fFiles);
2941}
2942
2943////////////////////////////////////////////////////////////////////////////////
2944/// This function transfroms the given TEventList into a TEntryList
2945///
2946/// NOTE, that this function loads all tree headers, because the entry numbers
2947/// in the TEventList are global and have to be recomputed, taking into account
2948/// the number of entries in each tree.
2949///
2950/// The new TEntryList is owned by the TChain and gets deleted when the chain
2951/// is deleted. This TEntryList is returned by GetEntryList() function, and after
2952/// GetEntryList() function is called, the TEntryList is not owned by the chain
2953/// any more and will not be deleted with it.
2954
2956{
2958 if (fEntryList) {
2960 TEntryList *tmp = fEntryList;
2961 fEntryList = nullptr; // Avoid problem with RecursiveRemove.
2962 delete tmp;
2963 } else {
2964 fEntryList = nullptr;
2965 }
2966 }
2967
2968 if (!evlist) {
2969 fEntryList = nullptr;
2970 fEventList = nullptr;
2971 return;
2972 }
2973
2974 if(fProofChain) {
2975 //on proof, fEventList and fEntryList shouldn't be set at the same time
2976 if (fEntryList){
2977 //check, if the chain is the owner of the previous entry list
2978 //(it happens, if the previous entry list was created from a user-defined
2979 //TEventList in SetEventList() function)
2981 TEntryList *tmp = fEntryList;
2982 fEntryList = nullptr; // Avoid problem with RecursiveRemove.
2983 delete tmp;
2984 } else {
2985 fEntryList = nullptr;
2986 }
2987 }
2988 return;
2989 }
2990
2991 char enlistname[100];
2992 snprintf(enlistname,100, "%s_%s", evlist->GetName(), "entrylist");
2993 TEntryList *enlist = new TEntryList(enlistname, evlist->GetTitle());
2994 enlist->SetDirectory(nullptr);
2995
2996 Int_t nsel = evlist->GetN();
2998 const char *treename;
2999 const char *filename;
3001 //Load all the tree headers if the tree offsets are not known
3002 //It is assumed here, that loading the last tree will load all
3003 //previous ones
3004 printf("loading trees\n");
3005 (const_cast<TChain*>(this))->LoadTree(evlist->GetEntry(evlist->GetN()-1));
3006 }
3007 for (Int_t i=0; i<nsel; i++){
3008 globalentry = evlist->GetEntry(i);
3009 //add some protection from globalentry<0 here
3010 Int_t treenum = 0;
3012 treenum++;
3013 treenum--;
3015 // printf("globalentry=%lld, treeoffset=%lld, localentry=%lld\n", globalentry, fTreeOffset[treenum], localentry);
3018 //printf("entering for tree %s %s\n", treename, filename);
3019 enlist->SetTree(treename, filename);
3020 enlist->Enter(localentry);
3021 }
3022 enlist->SetBit(kCanDelete, true);
3023 enlist->SetReapplyCut(evlist->GetReapplyCut());
3025}
3026
3027////////////////////////////////////////////////////////////////////////////////
3028/// Change the name of this TChain.
3029
3030void TChain::SetName(const char* name)
3031{
3032 if (fGlobalRegistration) {
3033 // Should this be extended to include the call to TTree::SetName?
3034 R__WRITE_LOCKGUARD(ROOT::gCoreMutex); // Take the lock once rather than 3 times.
3035 gROOT->GetListOfCleanups()->Remove(this);
3036 gROOT->GetListOfSpecials()->Remove(this);
3037 gROOT->GetListOfDataSets()->Remove(this);
3038 }
3040 if (fGlobalRegistration) {
3041 // Should this be extended to include the call to TTree::SetName?
3042 R__WRITE_LOCKGUARD(ROOT::gCoreMutex); // Take the lock once rather than 3 times.
3043 gROOT->GetListOfCleanups()->Add(this);
3044 gROOT->GetListOfSpecials()->Add(this);
3045 gROOT->GetListOfDataSets()->Add(this);
3046 }
3047}
3048
3049////////////////////////////////////////////////////////////////////////////////
3050/// Set number of entries per packet for parallel root.
3051
3053{
3054 fPacketSize = size;
3055 TIter next(fFiles);
3057 while ((element = (TChainElement*)next())) {
3058 element->SetPacketSize(size);
3059 }
3060}
3061
3062////////////////////////////////////////////////////////////////////////////////
3063/// Enable/Disable PROOF processing on the current default Proof (gProof).
3064///
3065/// "Draw" and "Processed" commands will be handled by PROOF.
3066/// The refresh and gettreeheader are meaningful only if on == true.
3067/// If refresh is true the underlying fProofChain (chain proxy) is always
3068/// rebuilt (even if already existing).
3069/// If gettreeheader is true the header of the tree will be read from the
3070/// PROOF cluster: this is only needed for browsing and should be used with
3071/// care because it may take a long time to execute.
3072
3074{
3075 if (!on) {
3076 // Disable
3078 // Reset related bit
3080 } else {
3081 if (fProofChain && !refresh &&
3083 return;
3084 }
3087
3088 // Make instance of TChainProof via the plugin manager
3090 if ((h = gROOT->GetPluginManager()->FindHandler("TChain", "proof"))) {
3091 if (h->LoadPlugin() == -1)
3092 return;
3093 if (!(fProofChain = reinterpret_cast<TChain *>(h->ExecPlugin(2, this, gettreeheader))))
3094 Error("SetProof", "creation of TProofChain failed");
3095 // Set related bits
3097 }
3098 }
3099}
3100
3101////////////////////////////////////////////////////////////////////////////////
3102/// Set chain weight.
3103///
3104/// The weight is used by TTree::Draw to automatically weight each
3105/// selected entry in the resulting histogram.
3106/// For example the equivalent of
3107/// ~~~ {.cpp}
3108/// chain.Draw("x","w")
3109/// ~~~
3110/// is
3111/// ~~~ {.cpp}
3112/// chain.SetWeight(w,"global");
3113/// chain.Draw("x");
3114/// ~~~
3115/// By default the weight used will be the weight
3116/// of each Tree in the TChain. However, one can force the individual
3117/// weights to be ignored by specifying the option "global".
3118/// In this case, the TChain global weight will be used for all Trees.
3119
3121{
3122 fWeight = w;
3123 TString opt = option;
3124 opt.ToLower();
3126 if (opt.Contains("global")) {
3128 }
3129}
3130
3131////////////////////////////////////////////////////////////////////////////////
3132/// Stream a class object.
3133
3135{
3136 if (b.IsReading()) {
3137 // Remove using the 'old' name.
3138 {
3140 gROOT->GetListOfCleanups()->Remove(this);
3141 }
3142
3143 UInt_t R__s, R__c;
3144 Version_t R__v = b.ReadVersion(&R__s, &R__c);
3145 if (R__v > 2) {
3146 b.ReadClassBuffer(TChain::Class(), this, R__v, R__s, R__c);
3147 } else {
3148 //====process old versions before automatic schema evolution
3150 b >> fTreeOffsetLen;
3151 b >> fNtrees;
3152 fFiles->Streamer(b);
3153 if (R__v > 1) {
3154 fStatus->Streamer(b);
3156 b.ReadFastArray(fTreeOffset,fTreeOffsetLen);
3157 }
3158 b.CheckByteCount(R__s, R__c, TChain::IsA());
3159 //====end of old versions
3160 }
3161 // Re-add using the new name.
3162 {
3164 gROOT->GetListOfCleanups()->Add(this);
3165 }
3166
3167 } else {
3168 b.WriteClassBuffer(TChain::Class(),this);
3169 }
3170}
3171
3172////////////////////////////////////////////////////////////////////////////////
3173/// Dummy function kept for back compatibility.
3174/// The cache is now activated automatically when processing TTrees/TChain.
3175
3176void TChain::UseCache(Int_t /* maxCacheSize */, Int_t /* pageSize */)
3177{
3178}
#define SafeDelete(p)
Definition RConfig.hxx:533
#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:406
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:2514
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:1136
void Streamer(TBuffer &) override
Stream a class object.
Definition TChain.cxx:3134
Long64_t GetEntryNumber(Long64_t entry) const override
Return entry number corresponding to entry.
Definition TChain.cxx:1032
Double_t GetWeight() const override
Return the chain weight.
Definition TChain.cxx:1288
virtual void SetAutoDelete(bool autodel=true)
Set the global branch kAutoDelete bit.
Definition TChain.cxx:2505
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:2764
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:1343
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:2289
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:2352
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:3073
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:2394
void ResetAfterMerge(TFileMergeInfo *) override
Resets the state of this chain after a merge (keep the customization but forget the data).
Definition TChain.cxx:2416
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:2553
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save TChain as a C++ statements on output stream out.
Definition TChain.cxx:2433
void SetEventList(TEventList *evlist) override
This function transfroms the given TEventList into a TEntryList.
Definition TChain.cxx:2955
TObjArray * GetListOfLeaves() override
Return a pointer to the list of leaves of the current tree.
Definition TChain.cxx:1159
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:1065
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:2490
Double_t GetMaximum(const char *columname) override
Return maximum of column with name columname.
Definition TChain.cxx:1180
void SetWeight(Double_t w=1, Option_t *option="") override
Set chain weight.
Definition TChain.cxx:3120
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:1877
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:2596
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:2312
void RemoveFriend(TTree *) override
Remove a friend from the list of friends.
Definition TChain.cxx:2371
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:2539
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:2727
Double_t GetMinimum(const char *columname) override
Return minimum of column with name columname.
Definition TChain.cxx:1215
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:1954
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:2219
bool InPlaceClone(TDirectory *newdirectory, const char *options="") override
Move content to a new file. (NOT IMPLEMENTED for TChain)
Definition TChain.cxx:1306
Int_t GetNbranches() override
Return the number of branches of the current tree.
Definition TChain.cxx:1252
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:2909
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:3030
void ls(Option_t *option="") const override
List the chain.
Definition TChain.cxx:1921
void SetBranchStatus(const char *bname, bool status=true, UInt_t *found=nullptr) override
Set branch status to Process or DoNotProcess.
Definition TChain.cxx:2700
Long64_t LoadTree(Long64_t entry) override
Find the tree which contains entry, and set it as the current tree.
Definition TChain.cxx:1371
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:1267
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:3176
void InvalidateCurrentTree()
Set the TTree to be reloaded as soon as possible.
Definition TChain.cxx:1325
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:1076
Int_t GetEntry(Long64_t entry=0, Int_t getall=0) override
Get entry from the file to memory.
Definition TChain.cxx:1013
virtual void SetPacketSize(Int_t size=100)
Set number of entries per packet for parallel root.
Definition TChain.cxx:3052
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:1807
TLeaf * GetLeaf(const char *branchname, const char *leafname) override
Return a pointer to the leaf name in the current tree.
Definition TChain.cxx:1089
@ 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:3073
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:2886
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition TROOT.cxx:2746
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:271
Helper class to prevent infinite recursion in the usage of TTree Friends.
Definition TTree.h:188
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual TFriendElement * AddFriend(const char *treename, const char *filename="")
Add a TFriendElement to the list of friends.
Definition TTree.cxx: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:8632
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:130
bool fIMTEnabled
! true if implicit multi-threading is enabled for this tree
Definition TTree.h:142
virtual bool GetBranchStatus(const char *branchname) const
Return status of branch with name branchname.
Definition TTree.cxx:5385
UInt_t fFriendLockStatus
! Record which method is locking the friend recursion
Definition TTree.h:137
TEventList * fEventList
! Pointer to event selection list (if one)
Definition TTree.h:125
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:5645
virtual TClusterIterator GetClusterIterator(Long64_t firstentry)
Return an iterator over the cluster of baskets starting at firstentry.
Definition TTree.cxx:5472
virtual void ResetBranchAddress(TBranch *)
Tell a branch to set its address to zero.
Definition TTree.cxx:8149
bool fCacheUserSet
! true if the cache setting was explicitly given by user
Definition TTree.h:141
virtual Long64_t GetEntryNumberWithIndex(Long64_t major, Long64_t minor=0) const
Return entry number corresponding to major and minor number.
Definition TTree.cxx:5917
virtual TObjArray * GetListOfLeaves()
Definition TTree.h:530
void Streamer(TBuffer &) override
Stream a class object.
Definition TTree.cxx:9654
TVirtualTreePlayer * GetPlayer()
Load the TTreePlayer (if not already done).
Definition TTree.cxx:6312
virtual Double_t GetWeight() const
Definition TTree.h:585
void Draw(Option_t *opt) override
Default Draw method for all objects.
Definition TTree.h:432
TVirtualTreePlayer * fPlayer
! Pointer to current Tree player
Definition TTree.h:134
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:9286
TTreeCache * GetReadCache(TFile *file) const
Find and return the TTreeCache registered with the file and which may contain branches for us.
Definition TTree.cxx:6325
Long64_t fEntries
Number of entries.
Definition TTree.h:84
TEntryList * fEntryList
! Pointer to event selection list (if one)
Definition TTree.h:126
virtual TVirtualIndex * GetTreeIndex() const
Definition TTree.h:559
TList * fExternalFriends
! List of TFriendsElement pointing to us and need to be notified of LoadTree. Content not owned.
Definition TTree.h:131
virtual void SetMaxVirtualSize(Long64_t size=0)
Definition TTree.h:666
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:7534
virtual void ResetAfterMerge(TFileMergeInfo *)
Resets the state of this TTree after a merge (keep the customization but forget the data).
Definition TTree.cxx:8118
virtual Long64_t GetEntries() const
Definition TTree.h:464
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:6202
virtual void Reset(Option_t *option="")
Reset baskets, buffers and entries count in all branches and leaves.
Definition TTree.cxx:8087
virtual void SetImplicitMT(bool enabled)
Definition TTree.h:662
Long64_t fMaxVirtualSize
Maximum total size of buffers kept in memory.
Definition TTree.h:99
TVirtualPerfStats * fPerfStats
! pointer to the current perf stats object
Definition TTree.h:132
Double_t fWeight
Tree weight (see TTree::SetWeight)
Definition TTree.h:90
virtual Long64_t GetReadEntry() const
Definition TTree.h:550
virtual TObjArray * GetListOfBranches()
Definition TTree.h:529
virtual TTree * GetTree() const
Definition TTree.h:558
virtual Long64_t LoadTree(Long64_t entry)
Set current entry.
Definition TTree.cxx:6480
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:8061
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:6564
TObject * fNotify
Object to be notified when loading a Tree.
Definition TTree.h:120
virtual TList * GetListOfClones()
Definition TTree.h:528
Long64_t fCacheSize
! Maximum size of file buffers
Definition TTree.h:105
TList * fClones
! List of cloned trees which share our addresses
Definition TTree.h:135
static TClass * Class()
@ kResetBranchAddresses
Definition TTree.h:225
@ kLoadTree
Definition TTree.h:221
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:531
Long64_t fReadEntry
! Number of the entry being processed
Definition TTree.h:107
virtual Int_t GetNbranches()
Definition TTree.h:543
virtual TLeaf * FindLeaf(const char *name)
Find leaf..
Definition TTree.cxx:4921
TDirectory * fDirectory
! Pointer to directory holding this tree
Definition TTree.h:121
@ kNeedEnableDecomposedObj
Definition TTree.h:245
@ kNoCheck
Definition TTree.h:244
@ kMissingBranch
Definition TTree.h:234
virtual void ResetBranchAddresses()
Tell all of our branches to drop their current objects and allocate new ones.
Definition TTree.cxx:8159
void SetName(const char *name) override
Change the name of this tree.
Definition TTree.cxx:9314
virtual Int_t SetCacheSize(Long64_t cachesize=-1)
Set maximum size of the file cache .
Definition TTree.cxx:8781
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:650
Int_t fPacketSize
! Number of entries in one packet for parallel root
Definition TTree.h:109
virtual Long64_t GetChainOffset() const
Definition TTree.h:457
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:8197
Int_t fMakeClass
! not zero when processing code generated by MakeClass
Definition TTree.h:115
static constexpr Long64_t kMaxEntries
Definition TTree.h:230
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