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