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