Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TTreeGeneratorBase.cxx
Go to the documentation of this file.
1// @(#)root/treeplayer:$Id$
2// Author: Akos Hajdu 13/08/2015
3
4/*************************************************************************
5 * Copyright (C) 1995-2015, Rene Brun and Fons Rademakers and al. *
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#include "TTreeGeneratorBase.h"
13
14#include "TBranchElement.h"
15#include "TClass.h"
16#include "TClassEdit.h"
17#include "TClonesArray.h"
18#include "TDirectory.h"
19#include "TFile.h"
20#include "TStreamerElement.h"
21#include "TStreamerInfo.h"
22#include "TTree.h"
23#include "TError.h"
26
27/** \class ROOT::Internal::TTreeGeneratorBase
28Base class for code generators like TTreeProxyGenerator and TTreeReaderGenerator
29*/
30
31namespace ROOT {
32namespace Internal {
33
34 ////////////////////////////////////////////////////////////////////////////////
35 /// Constructor.
36
37 TTreeGeneratorBase::TTreeGeneratorBase(TTree *tree, const char *option) : fTree(tree), fOptionStr(option) { }
38
39 ////////////////////////////////////////////////////////////////////////////////
40 /// Add a header inclusion request. If the header is already included it will
41 /// not be included again.
42
44 {
45 if (cl==nullptr) return;
46
47 // Check if already included
49 if (obj) return;
50
51 TString directive;
52
53 // Extract inner class from collection
56 }
57
58 // Construct directive
59 Int_t stlType;
60 if (0 == strcmp(cl->GetName(), "string")) { // Check if it is a string
61 directive = "#include <string>\n";
62 } else if (cl->GetCollectionProxy() && (stlType = cl->GetCollectionType())) { // Check if it is an STL container
63 const char *what = "";
64 switch(stlType) {
65 case ROOT::kSTLvector: what = "vector"; break;
66 case ROOT::kSTLlist: what = "list"; break;
67 case ROOT::kSTLforwardlist: what = "forward_list"; break;
68 case -ROOT::kSTLdeque: // same as positive
69 case ROOT::kSTLdeque: what = "deque"; break;
70 case -ROOT::kSTLmap: // same as positive
71 case ROOT::kSTLmap: what = "map"; break;
72 case -ROOT::kSTLmultimap: // same as positive
73 case ROOT::kSTLmultimap: what = "map"; break;
74 case -ROOT::kSTLset: // same as positive
75 case ROOT::kSTLset: what = "set"; break;
76 case -ROOT::kSTLmultiset: // same as positive
77 case ROOT::kSTLmultiset: what = "set"; break;
78 case -ROOT::kSTLunorderedset: // same as positive
79 case ROOT::kSTLunorderedset: what = "unordered_set"; break;
80 case -ROOT::kSTLunorderedmultiset: // same as positive
81 case ROOT::kSTLunorderedmultiset: what = "unordered_multiset"; break;
82 case -ROOT::kSTLunorderedmap: // same as positive
83 case ROOT::kSTLunorderedmap: what = "unordered_map"; break;
84 case -ROOT::kSTLunorderedmultimap: // same as positive
85 case ROOT::kSTLunorderedmultimap: what = "unordered_multimap"; break;
86 case -ROOT::kROOTRVec: // same as positive
87 case ROOT::kROOTRVec: what = "ROOT/RVec.hxx"; break;
88 }
89 if (what[0]) {
90 directive = "#include <";
91 directive.Append(what);
92 directive.Append(">\n");
93 }
94 } else if (TClassEdit::IsStdPair(cl->GetName())) {
96 // 4 elements expected: "pair", "first type name", "second type name", "trailing stars"
97 // However legacy code had a test for 3, we will leave it here until
98 // a test is developed (or found :) ) that exercise these lines of code.
99 if (split.fElements.size() == 3 || split.fElements.size() == 4) {
100 for (int arg = 1; arg < 3; ++arg) {
101 TClass* clArg = TClass::GetClass(split.fElements[arg].c_str());
102 if (clArg) AddHeader(clArg);
103 }
104 }
105 } else if (cl->GetDeclFileName() && strlen(cl->GetDeclFileName()) ) { // Custom file
106 const char *filename = cl->GetDeclFileName();
107
108 if (!filename) return;
109
110#ifdef R__WIN32
111 TString inclPath("include;prec_stl"); // GetHtml()->GetIncludePath());
112#else
113 TString inclPath("include:prec_stl"); // GetHtml()->GetIncludePath());
114#endif
115 Ssiz_t posDelim = 0;
116 TString inclDir;
117 TString sIncl(filename);
118#ifdef R__WIN32
119 const char* pdelim = ";";
120 static const char ddelim = '\\';
121#else
122 const char* pdelim = ":";
123 static const char ddelim = '/';
124#endif
125 while (inclPath.Tokenize(inclDir, posDelim, pdelim))
126 {
127 if (sIncl.BeginsWith(inclDir)) {
128 filename += inclDir.Length();
129 if (filename[0] == ddelim || filename[0] == '/') {
130 ++filename;
131 }
132 break;
133 }
134 }
135 directive = Form("#include \"%s\"\n",filename);
136 }
137 // Add directive (if it is not added already)
138 if (directive.Length()) {
139 TIter i( &fListOfHeaders );
140 for(TNamed *n = (TNamed *)i(); n; n = (TNamed*)i()) {
141 if (directive == n->GetTitle()) {
142 return;
143 }
144 }
145 fListOfHeaders.Add(new TNamed(cl->GetName(), directive.Data()));
146 }
147 }
148
149 ////////////////////////////////////////////////////////////////////////////////
150 /// Add a header inclusion request. If the header is already included it will
151 /// not be included again.
152
153 void TTreeGeneratorBase::AddHeader(const char *classname)
154 {
155 AddHeader(TClass::GetClass(classname));
156 }
157
158 ////////////////////////////////////////////////////////////////////////////////
159 /// Get name of class inside a container.
160
162 {
163 TString cname = branch->GetClonesName();
164 if (cname.Length()==0) {
165 // We may have any unsplit clones array
166 Long64_t i = branch->GetTree()->GetReadEntry();
167 if (i<0) i = 0;
168 branch->GetEntry(i);
169 char *obj = branch->GetObject();
170
171 TBranchElement *parent = (TBranchElement*)branch->GetMother()->GetSubBranch(branch);
172 const char *pclname = parent->GetClassName();
173
174 TClass *clparent = TClass::GetClass(pclname);
175 // TClass *clm = TClass::GetClass(GetClassName());
176 Int_t lOffset = 0; // offset in the local streamerInfo.
177 if (clparent) {
178 const char *ename = nullptr;
179 if (element) {
180 ename = element->GetName();
181 lOffset = clparent->GetStreamerInfo()->GetOffset(ename);
182 } else {
183 lOffset = 0;
184 }
185 }
186 else Error("AnalyzeBranch", "Missing parent for %s.", branch->GetName());
187
188 TClonesArray *arr;
189 if (ispointer) {
190 arr = (TClonesArray*)*(void**)(obj+lOffset);
191 } else {
192 arr = (TClonesArray*)(obj+lOffset);
193 }
194 cname = arr->GetClass()->GetName();
195 }
196 if (cname.Length()==0) {
197 Error("AnalyzeBranch",
198 "Introspection of TClonesArray in older file not implemented yet.");
199 }
200 return cname;
201 }
202
203 ////////////////////////////////////////////////////////////////////////////////
204 /// Check if element is a base class and if yes, return the base class.
205
207 {
208 TStreamerBase *base = dynamic_cast<TStreamerBase*>(element);
209 if (base) {
211 if (info) return info;
212 }
213 return nullptr;
214 }
215
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Return the correct TStreamerInfo of class 'cl' in the list of branches
218 /// (current) [Assuming these branches correspond to a flattened version of
219 /// the class.]
220
222 {
223 TVirtualStreamerInfo *objInfo = nullptr;
224 TBranchElement *b = nullptr;
225 TString cname = cl->GetName();
226
227 while( ( b = (TBranchElement*)current() ) ) {
228 if ( cname == b->GetInfo()->GetName() ) {
229 objInfo = b->GetInfo();
230 break;
231 }
232 }
233 if (objInfo == nullptr && branch->GetTree()->GetDirectory()->GetFile()) {
234 const TList *infolist = branch->GetTree()->GetDirectory()->GetFile()->GetStreamerInfoCache();
235 if (infolist) {
237 if (i) {
238 // NOTE: Is this correct for Foreigh classes?
240 }
241 }
242 }
243 if (objInfo == nullptr) {
244 // We still haven't found it ... this is likely to be an STL collection .. anyway, use the current StreamerInfo.
245 objInfo = cl->GetStreamerInfo();
246 }
247 return objInfo;
248 }
249
250} // namespace Internal
251} // namespace ROOT
#define b(i)
Definition RSha256.hxx:100
long long Long64_t
Definition RtypesCore.h:80
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
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 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 cname
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
TList fListOfHeaders
List of included headers.
void AddHeader(TClass *cl)
Add a header inclusion request.
TVirtualStreamerInfo * GetStreamerInfo(TBranch *branch, TIter current, TClass *cl)
Return the correct TStreamerInfo of class 'cl' in the list of branches (current) [Assuming these bran...
TVirtualStreamerInfo * GetBaseClass(TStreamerElement *element)
Check if element is a base class and if yes, return the base class.
TString GetContainedClassName(TBranchElement *branch, TStreamerElement *element, bool ispointer)
Get name of class inside a container.
TTreeGeneratorBase(TTree *tree, const char *option)
Constructor.
A Branch for the case of an object.
const char * GetClassName() const override
Return the name of the user class whose content is stored in this branch, if any.
Int_t GetEntry(Long64_t entry=0, Int_t getall=0) override
Read all branches of a BranchElement and return total number of bytes.
char * GetObject() const
Return a pointer to our object.
virtual const char * GetClonesName() const
A TTree is a list of TBranches.
Definition TBranch.h:93
TTree * GetTree() const
Definition TBranch.h:252
TBranch * GetSubBranch(const TBranch *br) const
Find the parent branch of child.
Definition TBranch.cxx:2164
TBranch * GetMother() const
Get our top-level parent branch in the tree.
Definition TBranch.cxx:2127
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
ROOT::ESTLType GetCollectionType() const
Return the 'type' of the STL the TClass is representing.
Definition TClass.cxx:2886
TVirtualStreamerInfo * GetStreamerInfo(Int_t version=0, Bool_t isTransient=kFALSE) const
returns a pointer to the TVirtualStreamerInfo object for version If the object does not exist,...
Definition TClass.cxx:4599
TVirtualCollectionProxy * GetCollectionProxy() const
Return the proxy describing the collection (if any).
Definition TClass.cxx:2897
const char * GetDeclFileName() const
Return name of the file containing the declaration of this class.
Definition TClass.cxx:3463
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:2968
An array of clone (identical) objects.
TClass * GetClass() const
virtual TFile * GetFile() const
Definition TDirectory.h:220
const TList * GetStreamerInfoCache()
Returns the cached list of StreamerInfos used in this file.
Definition TFile.cxx:1366
A doubly linked list.
Definition TList.h:38
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:83
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
Mother of all ROOT objects.
Definition TObject.h:41
TVirtualStreamerInfo * GetBaseStreamerInfo() const
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:417
const char * Data() const
Definition TString.h:376
TObjArray * Tokenize(const TString &delim) const
This function is used to isolate sequential tokens in a TString.
Definition TString.cxx:2264
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
Definition TString.h:623
TString & Append(const char *cs)
Definition TString.h:572
A TTree represents a columnar dataset.
Definition TTree.h:79
TDirectory * GetDirectory() const
Definition TTree.h:462
virtual Long64_t GetReadEntry() const
Definition TTree.h:509
virtual TClass * GetValueClass() const =0
If the value type is a user-defined class, return a pointer to the TClass representing the value type...
Abstract Interface class describing Streamer information for one class.
virtual Int_t GetOffset(const char *) const =0
virtual Int_t GetClassVersion() const =0
const Int_t n
Definition legend1.C:16
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
@ kSTLmap
Definition ESTLType.h:33
@ kSTLunorderedmultiset
Definition ESTLType.h:43
@ kROOTRVec
Definition ESTLType.h:46
@ kSTLset
Definition ESTLType.h:35
@ kSTLmultiset
Definition ESTLType.h:36
@ kSTLdeque
Definition ESTLType.h:32
@ kSTLvector
Definition ESTLType.h:30
@ kSTLunorderedmultimap
Definition ESTLType.h:45
@ kSTLunorderedset
Definition ESTLType.h:42
@ kSTLlist
Definition ESTLType.h:31
@ kSTLforwardlist
Definition ESTLType.h:41
@ kSTLunorderedmap
Definition ESTLType.h:44
@ kSTLmultimap
Definition ESTLType.h:34
bool IsStdPair(std::string_view name)
Definition TClassEdit.h:184
static const char * what
Definition stlLoader.cc:5
std::vector< std::string > fElements
Definition TClassEdit.h:141