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