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