Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
clingwrapper.cxx
Go to the documentation of this file.
1// Bindings
2#include "capi.h"
3#include "cpp_cppyy.h"
4#include "callcontext.h"
5
6// ROOT
7#include "TBaseClass.h"
8#include "TClass.h"
9#include "TClassRef.h"
10#include "TClassTable.h"
11#include "TClassEdit.h"
12#include "TCollection.h"
13#include "TDataMember.h"
14#include "TDataType.h"
15#include "TEnum.h"
16#include "TEnumConstant.h"
17#include "TEnv.h"
18#include "TError.h"
19#include "TException.h"
20#include "TFunction.h"
21#include "TFunctionTemplate.h"
22#include "TGlobal.h"
23#include "THashList.h"
24#include "TInterpreter.h"
25#include "TList.h"
26#include "TListOfDataMembers.h"
27#include "TListOfEnums.h"
28#include "TMethod.h"
29#include "TMethodArg.h"
30#include "TROOT.h"
31#include "TSystem.h"
32
33// Standard
34#include <assert.h>
35#include <algorithm> // for std::count, std::remove
36#include <stdexcept>
37#include <map>
38#include <new>
39#include <set>
40#include <sstream>
41#include <signal.h>
42#include <stdlib.h> // for getenv
43#include <string.h>
44#include <typeinfo>
45
46// temp
47#include <iostream>
49// --temp
50
51
52// small number that allows use of stack for argument passing
53const int SMALL_ARGS_N = 8;
54
55// data for life time management ---------------------------------------------
56typedef std::vector<TClassRef> ClassRefs_t;
58static const ClassRefs_t::size_type GLOBAL_HANDLE = 1;
59static const ClassRefs_t::size_type STD_HANDLE = GLOBAL_HANDLE + 1;
60
61typedef std::map<std::string, ClassRefs_t::size_type> Name2ClassRefIndex_t;
63
64namespace {
65
66static inline Cppyy::TCppType_t find_memoized(const std::string& name)
67{
68 auto icr = g_name2classrefidx.find(name);
69 if (icr != g_name2classrefidx.end())
70 return (Cppyy::TCppType_t)icr->second;
71 return (Cppyy::TCppType_t)0;
72}
73
74class CallWrapper {
75public:
76 typedef const void* DeclId_t;
77
78public:
79 CallWrapper(TFunction* f) : fDecl(f->GetDeclId()), fName(f->GetName()), fTF(nullptr) {}
80 CallWrapper(DeclId_t fid, const std::string& n) : fDecl(fid), fName(n), fTF(nullptr) {}
81 ~CallWrapper() {
82 if (fTF && fDecl == fTF->GetDeclId())
83 delete fTF;
84 }
85
86public:
88 DeclId_t fDecl;
89 std::string fName;
90 TFunction* fTF;
91};
92
93}
94
95static std::vector<CallWrapper*> gWrapperHolder;
96static inline CallWrapper* new_CallWrapper(TFunction* f)
97{
98 CallWrapper* wrap = new CallWrapper(f);
99 gWrapperHolder.push_back(wrap);
100 return wrap;
101}
102
103static inline CallWrapper* new_CallWrapper(CallWrapper::DeclId_t fid, const std::string& n)
104{
105 CallWrapper* wrap = new CallWrapper(fid, n);
106 gWrapperHolder.push_back(wrap);
107 return wrap;
108}
109
110typedef std::vector<TGlobal*> GlobalVars_t;
112
113static std::set<std::string> gSTLNames;
114
115
116// data ----------------------------------------------------------------------
118
119// builtin types (including a few common STL templates as long as they live in
120// the global namespace b/c of choices upstream)
121static std::set<std::string> g_builtins =
122 {"bool", "char", "signed char", "unsigned char", "wchar_t", "short", "unsigned short",
123 "int", "unsigned int", "long", "unsigned long", "long long", "unsigned long long",
124 "float", "double", "long double", "void",
125 "allocator", "array", "basic_string", "complex", "initializer_list", "less", "list",
126 "map", "pair", "set", "vector"};
127
128// smart pointer types
129static std::set<std::string> gSmartPtrTypes =
130 {"auto_ptr", "std::auto_ptr", "shared_ptr", "std::shared_ptr",
131 "unique_ptr", "std::unique_ptr", "weak_ptr", "std::weak_ptr"};
132
133// to filter out ROOT names
134static std::set<std::string> gInitialNames;
135static std::set<std::string> gRootSOs;
136
137// configuration
138static bool gEnableFastPath = true;
139
140
141// global initialization -----------------------------------------------------
142namespace {
143
144// names copied from TUnixSystem
145#ifdef WIN32
146const int SIGBUS = 0; // simple placeholders for ones that don't exist
147const int SIGSYS = 0;
148const int SIGPIPE = 0;
149const int SIGQUIT = 0;
150const int SIGWINCH = 0;
151const int SIGALRM = 0;
152const int SIGCHLD = 0;
153const int SIGURG = 0;
154const int SIGUSR1 = 0;
155const int SIGUSR2 = 0;
156#endif
157
158static struct Signalmap_t {
159 int fCode;
160 const char *fSigName;
161} gSignalMap[kMAXSIGNALS] = { // the order of the signals should be identical
162 { SIGBUS, "bus error" }, // to the one in TSysEvtHandler.h
163 { SIGSEGV, "segmentation violation" },
164 { SIGSYS, "bad argument to system call" },
165 { SIGPIPE, "write on a pipe with no one to read it" },
166 { SIGILL, "illegal instruction" },
167 { SIGABRT, "abort" },
168 { SIGQUIT, "quit" },
169 { SIGINT, "interrupt" },
170 { SIGWINCH, "window size change" },
171 { SIGALRM, "alarm clock" },
172 { SIGCHLD, "death of a child" },
173 { SIGURG, "urgent data arrived on an I/O channel" },
174 { SIGFPE, "floating point exception" },
175 { SIGTERM, "termination signal" },
176 { SIGUSR1, "user-defined signal 1" },
177 { SIGUSR2, "user-defined signal 2" }
178};
179
180static void inline do_trace(int sig) {
181 std::cerr << " *** Break *** " << (sig < kMAXSIGNALS ? gSignalMap[sig].fSigName : "") << std::endl;
183}
184
185class TExceptionHandlerImp : public TExceptionHandler {
186public:
187 virtual void HandleException(Int_t sig) {
188 if (TROOT::Initialized()) {
189 if (gException) {
190 gInterpreter->RewindDictionary();
191 gInterpreter->ClearFileBusy();
192 }
193
194 if (!getenv("CPPYY_CRASH_QUIET"))
195 do_trace(sig);
196
197 // jump back, if catch point set
198 Throw(sig);
199 }
200
201 do_trace(sig);
202 gSystem->Exit(128 + sig);
203 }
204};
205
206class ApplicationStarter {
207public:
208 ApplicationStarter() {
209 // initialize ROOT early to guarantee proper order of shutdown later on (gROOT is a
210 // macro that resolves to the ROOT::GetROOT() function call)
211 (void)gROOT;
212
213 // setup dummy holders for global and std namespaces
214 assert(g_classrefs.size() == GLOBAL_HANDLE);
216 g_classrefs.push_back(TClassRef(""));
217
218 // aliases for std (setup already in pythonify)
220 g_name2classrefidx["::std"] = g_name2classrefidx["std"];
221 g_classrefs.push_back(TClassRef("std"));
222
223 // add a dummy global to refer to as null at index 0
224 g_globalvars.push_back(nullptr);
225
226 // disable fast path if requested
227 if (getenv("CPPYY_DISABLE_FASTPATH")) gEnableFastPath = false;
228
229 // fill the set of STL names
230 const char* stl_names[] = {"allocator", "auto_ptr", "bad_alloc", "bad_cast",
231 "bad_exception", "bad_typeid", "basic_filebuf", "basic_fstream", "basic_ifstream",
232 "basic_ios", "basic_iostream", "basic_istream", "basic_istringstream",
233 "basic_ofstream", "basic_ostream", "basic_ostringstream", "basic_streambuf",
234 "basic_string", "basic_stringbuf", "basic_stringstream", "binary_function",
235 "binary_negate", "bitset", "byte", "char_traits", "codecvt_byname", "codecvt", "collate",
236 "collate_byname", "compare", "complex", "ctype_byname", "ctype", "default_delete",
237 "deque", "divides", "domain_error", "equal_to", "exception", "forward_list", "fpos",
238 "function", "greater_equal", "greater", "gslice_array", "gslice", "hash", "indirect_array",
239 "integer_sequence", "invalid_argument", "ios_base", "istream_iterator", "istreambuf_iterator",
240 "istrstream", "iterator_traits", "iterator", "length_error", "less_equal", "less",
241 "list", "locale", "localedef utility", "locale utility", "logic_error", "logical_and",
242 "logical_not", "logical_or", "map", "mask_array", "mem_fun", "mem_fun_ref", "messages",
243 "messages_byname", "minus", "modulus", "money_get", "money_put", "moneypunct",
244 "moneypunct_byname", "multimap", "multiplies", "multiset", "negate", "not_equal_to",
245 "num_get", "num_put", "numeric_limits", "numpunct", "numpunct_byname",
246 "ostream_iterator", "ostreambuf_iterator", "ostrstream", "out_of_range",
247 "overflow_error", "pair", "plus", "pointer_to_binary_function",
248 "pointer_to_unary_function", "priority_queue", "queue", "range_error",
249 "raw_storage_iterator", "reverse_iterator", "runtime_error", "set", "shared_ptr",
250 "slice_array", "slice", "stack", "string", "strstream", "strstreambuf",
251 "time_get_byname", "time_get", "time_put_byname", "time_put", "unary_function",
252 "unary_negate", "unique_ptr", "underflow_error", "unordered_map", "unordered_multimap",
253 "unordered_multiset", "unordered_set", "valarray", "vector", "weak_ptr", "wstring"};
254 for (auto& name : stl_names)
255 gSTLNames.insert(name);
256
257 // set opt level (default to 2 if not given; Cling itself defaults to 0)
258 int optLevel = 2;
259 if (getenv("CPPYY_OPT_LEVEL")) optLevel = atoi(getenv("CPPYY_OPT_LEVEL"));
260 if (optLevel != 0) {
261 std::ostringstream s;
262 s << "#pragma cling optimize " << optLevel;
263 gInterpreter->ProcessLine(s.str().c_str());
264 }
265
266 // load frequently used headers
267 const char* code =
268 "#include <iostream>\n"
269 "#include <string>\n"
270 "#include <DllImport.h>\n" // defines R__EXTERN
271 "#include <vector>\n"
272 "#include <utility>";
273 gInterpreter->ProcessLine(code);
274
275 // make sure we run in batch mode as far as ROOT graphics is concerned
276 if (!getenv("ROOTSYS"))
277 gROOT->SetBatch(kTRUE);
278
279 // create helpers for comparing thingies
280 gInterpreter->Declare(
281 "namespace __cppyy_internal { template<class C1, class C2>"
282 " bool is_equal(const C1& c1, const C2& c2) { return (bool)(c1 == c2); } }");
283 gInterpreter->Declare(
284 "namespace __cppyy_internal { template<class C1, class C2>"
285 " bool is_not_equal(const C1& c1, const C2& c2) { return (bool)(c1 != c2); } }");
286
287 // retrieve all initial (ROOT) C++ names in the global scope to allow filtering later
288 if (!getenv("CPPYY_NO_ROOT_FILTER")) {
289 gROOT->GetListOfGlobals(true); // force initialize
290 gROOT->GetListOfGlobalFunctions(true); // id.
291 std::set<std::string> initial;
293 gInitialNames = initial;
294
295#ifndef WIN32
296 gRootSOs.insert("libCore.so ");
297 gRootSOs.insert("libRIO.so ");
298 gRootSOs.insert("libThread.so ");
299 gRootSOs.insert("libMathCore.so ");
300#else
301 gRootSOs.insert("libCore.dll ");
302 gRootSOs.insert("libRIO.dll ");
303 gRootSOs.insert("libThread.dll ");
304 gRootSOs.insert("libMathCore.dll ");
305#endif
306 }
307
308 // start off with a reasonable size placeholder for wrappers
309 gWrapperHolder.reserve(1024);
310
311 // create an exception handler to process signals
312 gExceptionHandler = new TExceptionHandlerImp{};
313 }
314
315 ~ApplicationStarter() {
316 for (auto wrap : gWrapperHolder)
317 delete wrap;
318 delete gExceptionHandler; gExceptionHandler = nullptr;
319 }
320} _applicationStarter;
321
322} // unnamed namespace
323
324
325// local helpers -------------------------------------------------------------
326static inline
328{
329 assert((ClassRefs_t::size_type)scope < g_classrefs.size());
330 return g_classrefs[(ClassRefs_t::size_type)scope];
331}
332
333static inline
335 CallWrapper* wrap = ((CallWrapper*)method);
336 if (!wrap->fTF || wrap->fTF->GetDeclId() != wrap->fDecl) {
337 MethodInfo_t* mi = gInterpreter->MethodInfo_Factory(wrap->fDecl);
338 wrap->fTF = new TFunction(mi);
339 }
340 return wrap->fTF;
341}
342
343static inline
344char* cppstring_to_cstring(const std::string& cppstr)
345{
346 char* cstr = (char*)malloc(cppstr.size()+1);
347 memcpy(cstr, cppstr.c_str(), cppstr.size()+1);
348 return cstr;
349}
350
351static inline
352bool match_name(const std::string& tname, const std::string fname)
353{
354// either match exactly, or match the name as template
355 if (fname.rfind(tname, 0) == 0) {
356 if ((tname.size() == fname.size()) ||
357 (tname.size() < fname.size() && fname[tname.size()] == '<'))
358 return true;
359 }
360 return false;
361}
362
363static inline
364bool is_missclassified_stl(const std::string& name)
365{
366 std::string::size_type pos = name.find('<');
367 if (pos != std::string::npos)
368 return gSTLNames.find(name.substr(0, pos)) != gSTLNames.end();
369 return gSTLNames.find(name) != gSTLNames.end();
370}
371
372
373// direct interpreter access -------------------------------------------------
374bool Cppyy::Compile(const std::string& code)
375{
376 return gInterpreter->Declare(code.c_str());
377}
378
379
380// name to opaque C++ scope representation -----------------------------------
381std::string Cppyy::ResolveName(const std::string& cppitem_name)
382{
383// Fully resolve the given name to the final type name.
384
385// try memoized type cache, in case seen before
386 TCppType_t klass = find_memoized(cppitem_name);
387 if (klass) return GetScopedFinalName(klass);
388
389// remove global scope '::' if present
390 std::string tclean = cppitem_name.compare(0, 2, "::") == 0 ?
391 cppitem_name.substr(2, std::string::npos) : cppitem_name;
392
393// classes (most common)
394 tclean = TClassEdit::CleanType(tclean.c_str());
395 if (tclean.empty() /* unknown, eg. an operator */) return cppitem_name;
396
397// reduce [N] to []
398 if (tclean[tclean.size()-1] == ']')
399 tclean = tclean.substr(0, tclean.rfind('[')) + "[]";
400
401 if (tclean.rfind("byte", 0) == 0 || tclean.rfind("std::byte", 0) == 0)
402 return tclean;
403
404// check data types list (accept only builtins as typedefs will
405// otherwise not be resolved)
406 TDataType* dt = gROOT->GetType(tclean.c_str());
407 if (dt && dt->GetType() != kOther_t) return dt->GetFullTypeName();
408
409// special case for enums
410 if (IsEnum(cppitem_name))
411 return ResolveEnum(cppitem_name);
412
413// special case for clang's builtin __type_pack_element (which does not resolve)
414 if (cppitem_name.rfind("__type_pack_element", 0) != std::string::npos) {
415 // shape is "__type_pack_element<index,type1,type2,...,typeN>cpd": extract
416 // first the index, and from there the indexed type; finally, restore the
417 // qualifiers
418 const char* str = cppitem_name.c_str();
419 char* endptr = nullptr;
420 unsigned long index = strtoul(str+20, &endptr, 0);
421
422 std::string tmplvars{endptr};
423 auto start = tmplvars.find(',') + 1;
424 auto end = tmplvars.find(',', start);
425 while (index != 0) {
426 start = end+1;
427 end = tmplvars.find(',', start);
428 if (end == std::string::npos) end = tmplvars.rfind('>');
429 --index;
430 }
431
432 std::string resolved = tmplvars.substr(start, end-start);
433 auto cpd = tmplvars.rfind('>');
434 if (cpd != std::string::npos && cpd+1 != tmplvars.size())
435 return resolved + tmplvars.substr(cpd+1, std::string::npos);
436 return resolved;
437 }
438
439// typedefs
440 return TClassEdit::ResolveTypedef(tclean.c_str(), true);
441}
442
443static std::map<std::string, std::string> resolved_enum_types;
444std::string Cppyy::ResolveEnum(const std::string& enum_type)
445{
446// The underlying type of a an enum may be any kind of integer.
447// Resolve that type via a workaround (note: this function assumes
448// that the enum_type name is a valid enum type name)
449 auto res = resolved_enum_types.find(enum_type);
450 if (res != resolved_enum_types.end())
451 return res->second;
452
453// desugar the type before resolving
454 std::string et_short = TClassEdit::ShortType(enum_type.c_str(), 1);
455 if (et_short.find("(anonymous") == std::string::npos) {
456 std::ostringstream decl;
457 // TODO: now presumed fixed with https://sft.its.cern.ch/jira/browse/ROOT-6988
458 for (auto& itype : {"unsigned int"}) {
459 decl << "std::is_same<"
460 << itype
461 << ", std::underlying_type<"
462 << et_short
463 << ">::type>::value;";
464 if (gInterpreter->ProcessLine(decl.str().c_str())) {
465 // TODO: "re-sugaring" like this is brittle, but the top
466 // should be re-translated into AST-based code anyway
467 std::string resugared;
468 if (et_short.size() != enum_type.size()) {
469 auto pos = enum_type.find(et_short);
470 if (pos != std::string::npos) {
471 resugared = enum_type.substr(0, pos) + itype;
472 if (pos+et_short.size() < enum_type.size())
473 resugared += enum_type.substr(pos+et_short.size(), std::string::npos);
474 }
475 }
476 if (resugared.empty()) resugared = itype;
477 resolved_enum_types[enum_type] = resugared;
478 return resugared;
479 }
480 }
481 }
482
483// failed or anonymous ... signal upstream to special case this
484 int ipos = (int)enum_type.size()-1;
485 for (; 0 <= ipos; --ipos) {
486 char c = enum_type[ipos];
487 if (isspace(c)) continue;
488 if (isalnum(c) || c == '_' || c == '>' || c == ')') break;
489 }
490 bool isConst = enum_type.find("const ", 6) != std::string::npos;
491 std::string restype = isConst ? "const " : "";
492 restype += "internal_enum_type_t"+enum_type.substr((std::string::size_type)ipos+1, std::string::npos);
493 resolved_enum_types[enum_type] = restype;
494 return restype; // should default to some int variant
495}
496
497Cppyy::TCppScope_t Cppyy::GetScope(const std::string& sname)
498{
499// First, try cache
500 TCppType_t result = find_memoized(sname);
501 if (result) return result;
502
503// Second, skip builtins before going through the more expensive steps of resolving
504// typedefs and looking up TClass
505 if (g_builtins.find(sname) != g_builtins.end())
506 return (TCppScope_t)0;
507
508// TODO: scope_name should always be final already?
509// Resolve name fully before lookup to make sure all aliases point to the same scope
510 std::string scope_name = ResolveName(sname);
511 bool bHasAlias = sname != scope_name;
512 if (bHasAlias) {
513 result = find_memoized(scope_name);
514 if (result) return result;
515 }
516
517// both failed, but may be STL name that's missing 'std::' now, but didn't before
518 bool b_scope_name_missclassified = is_missclassified_stl(scope_name);
519 if (b_scope_name_missclassified) {
520 result = find_memoized("std::"+scope_name);
521 if (result) g_name2classrefidx["std::"+scope_name] = (ClassRefs_t::size_type)result;
522 }
523 bool b_sname_missclassified = bHasAlias ? is_missclassified_stl(sname) : false;
524 if (b_sname_missclassified) {
525 if (!result) result = find_memoized("std::"+sname);
526 if (result) g_name2classrefidx["std::"+sname] = (ClassRefs_t::size_type)result;
527 }
528
529 if (result) return result;
530
531// use TClass directly, to enable auto-loading; class may be stubbed (eg. for
532// function returns) or forward declared, leading to a non-null TClass that is
533// otherwise invalid/unusable
534 TClassRef cr(TClass::GetClass(scope_name.c_str(), true /* load */, true /* silent */));
535 if (!cr.GetClass())
536 return (TCppScope_t)0;
537
538// memoize found/created TClass
539 ClassRefs_t::size_type sz = g_classrefs.size();
540 g_name2classrefidx[scope_name] = sz;
541 if (bHasAlias) g_name2classrefidx[sname] = sz;
542 g_classrefs.push_back(TClassRef(scope_name.c_str()));
543
544// TODO: make ROOT/meta NOT remove std :/
545 if (b_scope_name_missclassified)
546 g_name2classrefidx["std::"+scope_name] = sz;
547 if (b_sname_missclassified)
548 g_name2classrefidx["std::"+sname] = sz;
549
550 return (TCppScope_t)sz;
551}
552
553bool Cppyy::IsTemplate(const std::string& template_name)
554{
555 return (bool)gInterpreter->CheckClassTemplate(template_name.c_str());
556}
557
558namespace {
559 class AutoCastRTTI {
560 public:
561 virtual ~AutoCastRTTI() {}
562 };
563}
564
566{
567 TClassRef& cr = type_from_handle(klass);
568 if (!cr.GetClass() || !obj) return klass;
569
570#ifdef _WIN64
571// Cling does not provide a consistent ImageBase address for calculating relative addresses
572// as used in Windows 64b RTTI. So, check for our own RTTI extension instead. If that fails,
573// see whether the unmangled raw_name is available (e.g. if this is an MSVC compiled rather
574// than JITed class) and pass on if it is.
575 volatile const char* raw = nullptr; // to prevent too aggressive reordering
576 try {
577 // this will filter those objects that do not have RTTI to begin with (throws)
578 AutoCastRTTI* pcst = (AutoCastRTTI*)obj;
579 raw = typeid(*pcst).raw_name();
580
581 // check the signature id (0 == absolute, 1 == relative, 2 == ours)
582 void* vfptr = *(void**)((intptr_t)obj);
583 void* meta = (void*)((intptr_t)*((void**)((intptr_t)vfptr-sizeof(void*))));
584 if (*(intptr_t*)meta == 2) {
585 // access the extra data item which is an absolute pointer to the RTTI
586 void* ptdescr = (void*)((intptr_t)meta + 4*sizeof(unsigned long)+sizeof(void*));
587 if (ptdescr && *(void**)ptdescr) {
588 auto rtti = *(std::type_info**)ptdescr;
589 raw = rtti->raw_name();
590 if (raw && raw[0] != '\0') // likely unnecessary
591 return (TCppType_t)GetScope(rtti->name());
592 }
593
594 return klass; // do not fall through if no RTTI info available
595 }
596
597 // if the raw name is the empty string (no guarantees that this is so as truly, the
598 // address is corrupt, but it is common to be empty), then there is no accessible RTTI
599 // and getting the unmangled name will crash ...
600 if (!raw || raw[0] == '\0')
601 return klass;
602 } catch (std::bad_typeid) {
603 return klass; // can't risk passing to ROOT/meta as it may do RTTI
604 }
605#endif
606
607 TClass* clActual = cr->GetActualClass((void*)obj);
608 if (clActual && clActual != cr.GetClass()) {
609 auto itt = g_name2classrefidx.find(clActual->GetName());
610 if (itt != g_name2classrefidx.end())
611 return (TCppType_t)itt->second;
612 return (TCppType_t)GetScope(clActual->GetName());
613 }
614
615 return klass;
616}
617
619{
620 TClassRef& cr = type_from_handle(klass);
621 if (cr.GetClass() && cr->GetClassInfo())
622 return (size_t)gInterpreter->ClassInfo_Size(cr->GetClassInfo());
623 return (size_t)0;
624}
625
626size_t Cppyy::SizeOf(const std::string& type_name)
627{
628 TDataType* dt = gROOT->GetType(type_name.c_str());
629 if (dt) return dt->Size();
630 return SizeOf(GetScope(type_name));
631}
632
633bool Cppyy::IsBuiltin(const std::string& type_name)
634{
635 TDataType* dt = gROOT->GetType(TClassEdit::CleanType(type_name.c_str(), 1).c_str());
636 if (dt) return dt->GetType() != kOther_t;
637 return false;
638}
639
640bool Cppyy::IsComplete(const std::string& type_name)
641{
642// verify whether the dictionary of this class is fully available
643 bool b = false;
644
645 int oldEIL = gErrorIgnoreLevel;
646 gErrorIgnoreLevel = 3000;
647 TClass* klass = TClass::GetClass(TClassEdit::ShortType(type_name.c_str(), 1).c_str());
648 if (klass && klass->GetClassInfo()) // works for normal case w/ dict
649 b = gInterpreter->ClassInfo_IsLoaded(klass->GetClassInfo());
650 else { // special case for forward declared classes
651 ClassInfo_t* ci = gInterpreter->ClassInfo_Factory(type_name.c_str());
652 if (ci) {
653 b = gInterpreter->ClassInfo_IsLoaded(ci);
654 gInterpreter->ClassInfo_Delete(ci); // we own the fresh class info
655 }
656 }
657 gErrorIgnoreLevel = oldEIL;
658 return b;
659}
660
661// memory management ---------------------------------------------------------
663{
665 return (TCppObject_t)malloc(gInterpreter->ClassInfo_Size(cr->GetClassInfo()));
666}
667
668void Cppyy::Deallocate(TCppType_t /* type */, TCppObject_t instance)
669{
670 ::operator delete(instance);
671}
672
674{
676 return (TCppObject_t)cr->New();
677}
678
679static std::map<Cppyy::TCppType_t, bool> sHasOperatorDelete;
681{
684 cr->Destructor((void*)instance);
685 else {
686 ROOT::DelFunc_t fdel = cr->GetDelete();
687 if (fdel) fdel((void*)instance);
688 else {
689 auto ib = sHasOperatorDelete.find(type);
690 if (ib == sHasOperatorDelete.end()) {
692 ib = sHasOperatorDelete.find(type);
693 }
694 ib->second ? cr->Destructor((void*)instance) : free((void*)instance);
695 }
696 }
697}
698
699
700// method/function dispatching -----------------------------------------------
702{
703// TODO: method should be a callfunc, so that no mapping would be needed.
704 CallWrapper* wrap = (CallWrapper*)method;
705
706 CallFunc_t* callf = gInterpreter->CallFunc_Factory();
707 MethodInfo_t* meth = gInterpreter->MethodInfo_Factory(wrap->fDecl);
708 gInterpreter->CallFunc_SetFunc(callf, meth);
709 gInterpreter->MethodInfo_Delete(meth);
710
711 if (!(callf && gInterpreter->CallFunc_IsValid(callf))) {
712 // TODO: propagate this error to caller w/o use of Python C-API
713 /*
714 PyErr_Format(PyExc_RuntimeError, "could not resolve %s::%s(%s)",
715 const_cast<TClassRef&>(klass).GetClassName(),
716 wrap.fName, callString.c_str()); */
717 std::cerr << "TODO: report unresolved function error to Python\n";
718 if (callf) gInterpreter->CallFunc_Delete(callf);
720 }
721
722// generate the wrapper and JIT it; ignore wrapper generation errors (will simply
723// result in a nullptr that is reported upstream if necessary; often, however,
724// there is a different overload available that will do)
725 auto oldErrLvl = gErrorIgnoreLevel;
727 wrap->fFaceptr = gInterpreter->CallFunc_IFacePtr(callf);
728 gErrorIgnoreLevel = oldErrLvl;
729
730 gInterpreter->CallFunc_Delete(callf); // does not touch IFacePtr
731 return wrap->fFaceptr;
732}
733
734static inline
735bool copy_args(Parameter* args, size_t nargs, void** vargs)
736{
737 bool runRelease = false;
738 for (size_t i = 0; i < nargs; ++i) {
739 switch (args[i].fTypeCode) {
740 case 'X': /* (void*)type& with free */
741 runRelease = true;
742 case 'V': /* (void*)type& */
743 vargs[i] = args[i].fValue.fVoidp;
744 break;
745 case 'r': /* const type& */
746 vargs[i] = args[i].fRef;
747 break;
748 default: /* all other types in union */
749 vargs[i] = (void*)&args[i].fValue.fVoidp;
750 break;
751 }
752 }
753 return runRelease;
754}
755
756static inline
757void release_args(Parameter* args, size_t nargs) {
758 for (size_t i = 0; i < nargs; ++i) {
759 if (args[i].fTypeCode == 'X')
760 free(args[i].fValue.fVoidp);
761 }
762}
763
764static inline bool WrapperCall(Cppyy::TCppMethod_t method, size_t nargs, void* args_, void* self, void* result)
765{
766 Parameter* args = (Parameter*)args_;
767
768 CallWrapper* wrap = (CallWrapper*)method;
769 const TInterpreter::CallFuncIFacePtr_t& faceptr = wrap->fFaceptr.fGeneric ? wrap->fFaceptr : GetCallFunc(method);
770 if (!faceptr.fGeneric)
771 return false; // happens with compilation error
772
774 bool runRelease = false;
775 if (nargs <= SMALL_ARGS_N) {
776 void* smallbuf[SMALL_ARGS_N];
777 if (nargs) runRelease = copy_args(args, nargs, smallbuf);
778 faceptr.fGeneric(self, (int)nargs, smallbuf, result);
779 } else {
780 std::vector<void*> buf(nargs);
781 runRelease = copy_args(args, nargs, buf.data());
782 faceptr.fGeneric(self, (int)nargs, buf.data(), result);
783 }
784 if (runRelease) release_args(args, nargs);
785 return true;
786 }
787
789 bool runRelease = false;
790 if (nargs <= SMALL_ARGS_N) {
791 void* smallbuf[SMALL_ARGS_N];
792 if (nargs) runRelease = copy_args(args, nargs, (void**)smallbuf);
793 faceptr.fCtor((void**)smallbuf, result, (unsigned long)nargs);
794 } else {
795 std::vector<void*> buf(nargs);
796 runRelease = copy_args(args, nargs, buf.data());
797 faceptr.fCtor(buf.data(), result, (unsigned long)nargs);
798 }
799 if (runRelease) release_args(args, nargs);
800 return true;
801 }
802
804 std::cerr << " DESTRUCTOR NOT IMPLEMENTED YET! " << std::endl;
805 return false;
806 }
807
808 return false;
809}
810
811template<typename T>
812static inline
813T CallT(Cppyy::TCppMethod_t method, Cppyy::TCppObject_t self, size_t nargs, void* args)
814{
815 T t{};
816 if (WrapperCall(method, nargs, args, (void*)self, &t))
817 return t;
818 return (T)-1;
819}
820
821#define CPPYY_IMP_CALL(typecode, rtype) \
822rtype Cppyy::Call##typecode(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)\
823{ \
824 return CallT<rtype>(method, self, nargs, args); \
825}
826
827void Cppyy::CallV(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)
828{
829 if (!WrapperCall(method, nargs, args, (void*)self, nullptr))
830 return /* TODO ... report error */;
831}
832
833CPPYY_IMP_CALL(B, unsigned char)
834CPPYY_IMP_CALL(C, char )
835CPPYY_IMP_CALL(H, short )
836CPPYY_IMP_CALL(I, int )
837CPPYY_IMP_CALL(L, long )
839CPPYY_IMP_CALL(F, float )
840CPPYY_IMP_CALL(D, double )
842
843void* Cppyy::CallR(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)
844{
845 void* r = nullptr;
846 if (WrapperCall(method, nargs, args, (void*)self, &r))
847 return r;
848 return nullptr;
849}
850
852 TCppMethod_t method, TCppObject_t self, size_t nargs, void* args, size_t* length)
853{
854 char* cstr = nullptr;
855 TClassRef cr("std::string");
856 std::string* cppresult = (std::string*)malloc(sizeof(std::string));
857 if (WrapperCall(method, nargs, args, self, (void*)cppresult)) {
858 cstr = cppstring_to_cstring(*cppresult);
859 *length = cppresult->size();
860 cppresult->std::string::~basic_string();
861 } else
862 *length = 0;
863 free((void*)cppresult);
864 return cstr;
865}
866
868 TCppMethod_t method, TCppType_t /* klass */, size_t nargs, void* args)
869{
870 void* obj = nullptr;
871 if (WrapperCall(method, nargs, args, nullptr, &obj))
872 return (TCppObject_t)obj;
873 return (TCppObject_t)0;
874}
875
877{
879 cr->Destructor((void*)self, true);
880}
881
883 TCppObject_t self, size_t nargs, void* args, TCppType_t result_type)
884{
885 TClassRef& cr = type_from_handle(result_type);
886 void* obj = ::operator new(gInterpreter->ClassInfo_Size(cr->GetClassInfo()));
887 if (WrapperCall(method, nargs, args, self, obj))
888 return (TCppObject_t)obj;
889 return (TCppObject_t)0;
890}
891
893{
894 if (check_enabled && !gEnableFastPath) return (TCppFuncAddr_t)nullptr;
895 TFunction* f = m2f(method);
896 return (TCppFuncAddr_t)gInterpreter->FindSym(f->GetMangledName());
897}
898
899
900// handling of function argument buffer --------------------------------------
902{
903 return new Parameter[nargs];
904}
905
907{
908 delete [] (Parameter*)args;
909}
910
912{
913 return sizeof(Parameter);
914}
915
917{
918 return offsetof(Parameter, fTypeCode);
919}
920
921
922// scope reflection information ----------------------------------------------
924{
925// Test if this scope represents a namespace.
926 if (scope == GLOBAL_HANDLE)
927 return true;
928 TClassRef& cr = type_from_handle(scope);
929 if (cr.GetClass())
930 return cr->Property() & kIsNamespace;
931 return false;
932}
933
935{
936// Test if this type may not be instantiated.
937 TClassRef& cr = type_from_handle(klass);
938 if (cr.GetClass())
939 return cr->Property() & kIsAbstract;
940 return false;
941}
942
943bool Cppyy::IsEnum(const std::string& type_name)
944{
945 if (type_name.empty()) return false;
946 std::string tn_short = TClassEdit::ShortType(type_name.c_str(), 1);
947 if (tn_short.empty()) return false;
948 return gInterpreter->ClassInfo_IsEnum(tn_short.c_str());
949}
950
951// helpers for stripping scope names
952static
953std::string outer_with_template(const std::string& name)
954{
955// Cut down to the outer-most scope from <name>, taking proper care of templates.
956 int tpl_open = 0;
957 for (std::string::size_type pos = 0; pos < name.size(); ++pos) {
958 std::string::value_type c = name[pos];
959
960 // count '<' and '>' to be able to skip template contents
961 if (c == '<')
962 ++tpl_open;
963 else if (c == '>')
964 --tpl_open;
965
966 // collect name up to "::"
967 else if (tpl_open == 0 && \
968 c == ':' && pos+1 < name.size() && name[pos+1] == ':') {
969 // found the extend of the scope ... done
970 return name.substr(0, pos-1);
971 }
972 }
973
974// whole name is apparently a single scope
975 return name;
976}
977
978static
979std::string outer_no_template(const std::string& name)
980{
981// Cut down to the outer-most scope from <name>, drop templates
982 std::string::size_type first_scope = name.find(':');
983 if (first_scope == std::string::npos)
984 return name.substr(0, name.find('<'));
985 std::string::size_type first_templ = name.find('<');
986 if (first_templ == std::string::npos)
987 return name.substr(0, first_scope);
988 return name.substr(0, std::min(first_templ, first_scope));
989}
990
991#define FILL_COLL(type, filter) { \
992 TIter itr{coll}; \
993 type* obj = nullptr; \
994 while ((obj = (type*)itr.Next())) { \
995 const char* nm = obj->GetName(); \
996 if (nm && nm[0] != '_' && !(obj->Property() & (filter))) { \
997 if (gInitialNames.find(nm) == gInitialNames.end()) \
998 cppnames.insert(nm); \
999 }}}
1000
1001static inline
1002void cond_add(Cppyy::TCppScope_t scope, const std::string& ns_scope,
1003 std::set<std::string>& cppnames, const char* name, bool nofilter = false)
1004{
1005 if (!name || name[0] == '_' || strstr(name, ".h") != 0 || strncmp(name, "operator", 8) == 0)
1006 return;
1007
1008 if (scope == GLOBAL_HANDLE) {
1009 std::string to_add = outer_no_template(name);
1010 if ((nofilter || gInitialNames.find(to_add) == gInitialNames.end()) && !is_missclassified_stl(name))
1011 cppnames.insert(outer_no_template(name));
1012 } else if (scope == STD_HANDLE) {
1013 if (strncmp(name, "std::", 5) == 0) {
1014 name += 5;
1015#ifdef __APPLE__
1016 if (strncmp(name, "__1::", 5) == 0) name += 5;
1017#endif
1018 } else if (!is_missclassified_stl(name))
1019 return;
1020 cppnames.insert(outer_no_template(name));
1021 } else {
1022 if (strncmp(name, ns_scope.c_str(), ns_scope.size()) == 0)
1023 cppnames.insert(outer_with_template(name + ns_scope.size()));
1024 }
1025}
1026
1027void Cppyy::GetAllCppNames(TCppScope_t scope, std::set<std::string>& cppnames)
1028{
1029// Collect all known names of C++ entities under scope. This is useful for IDEs
1030// employing tab-completion, for example. Note that functions names need not be
1031// unique as they can be overloaded.
1032 TClassRef& cr = type_from_handle(scope);
1033 if (scope != GLOBAL_HANDLE && !(cr.GetClass() && cr->Property()))
1034 return;
1035
1036 std::string ns_scope = GetFinalName(scope);
1037 if (scope != GLOBAL_HANDLE) ns_scope += "::";
1038
1039// add existing values from read rootmap files if within this scope
1040 TCollection* coll = gInterpreter->GetMapfile()->GetTable();
1041 {
1042 TIter itr{coll};
1043 TEnvRec* ev = nullptr;
1044 while ((ev = (TEnvRec*)itr.Next())) {
1045 // TEnv contains rootmap entries and user-side rootmap files may be already
1046 // loaded on startup. Thus, filter on file name rather than load time.
1047 if (gRootSOs.find(ev->GetValue()) == gRootSOs.end())
1048 cond_add(scope, ns_scope, cppnames, ev->GetName(), true);
1049 }
1050 }
1051
1052// do we care about the class table or are the rootmap and list of types enough?
1053/*
1054 gClassTable->Init();
1055 const int N = gClassTable->Classes();
1056 for (int i = 0; i < N; ++i)
1057 cond_add(scope, ns_scope, cppnames, gClassTable->Next());
1058*/
1059
1060// any other types (e.g. that may have come from parsing headers)
1061 coll = gROOT->GetListOfTypes();
1062 {
1063 TIter itr{coll};
1064 TDataType* dt = nullptr;
1065 while ((dt = (TDataType*)itr.Next())) {
1066 if (!(dt->Property() & kIsFundamental)) {
1067 cond_add(scope, ns_scope, cppnames, dt->GetName());
1068 }
1069 }
1070 }
1071
1072// add functions
1073 coll = (scope == GLOBAL_HANDLE) ?
1074 gROOT->GetListOfGlobalFunctions() : cr->GetListOfMethods();
1075 {
1076 TIter itr{coll};
1077 TFunction* obj = nullptr;
1078 while ((obj = (TFunction*)itr.Next())) {
1079 const char* nm = obj->GetName();
1080 // skip templated functions, adding only the un-instantiated ones
1081 if (nm && nm[0] != '_' && strstr(nm, "<") == 0 && strncmp(nm, "operator", 8) != 0) {
1082 if (gInitialNames.find(nm) == gInitialNames.end())
1083 cppnames.insert(nm);
1084 }
1085 }
1086 }
1087
1088// add uninstantiated templates
1089 coll = (scope == GLOBAL_HANDLE) ?
1090 gROOT->GetListOfFunctionTemplates() : cr->GetListOfFunctionTemplates();
1092
1093// add (global) data members
1094 if (scope == GLOBAL_HANDLE) {
1095 coll = gROOT->GetListOfGlobals();
1097 } else {
1098 coll = cr->GetListOfDataMembers();
1100 coll = cr->GetListOfUsingDataMembers();
1102 }
1103
1104// add enums values only for user classes/namespaces
1105 if (scope != GLOBAL_HANDLE && scope != STD_HANDLE) {
1106 coll = cr->GetListOfEnums();
1108 }
1109
1110#ifdef __APPLE__
1111// special case for Apple, add version namespace '__1' entries to std
1112 if (scope == STD_HANDLE)
1113 GetAllCppNames(GetScope("std::__1"), cppnames);
1114#endif
1115}
1116
1117
1118// class reflection information ----------------------------------------------
1119std::vector<Cppyy::TCppScope_t> Cppyy::GetUsingNamespaces(TCppScope_t scope)
1120{
1121 std::vector<Cppyy::TCppScope_t> res;
1122 if (!IsNamespace(scope))
1123 return res;
1124
1125#ifdef __APPLE__
1126 if (scope == STD_HANDLE) {
1127 res.push_back(GetScope("__1"));
1128 return res;
1129 }
1130#endif
1131
1132 TClassRef& cr = type_from_handle(scope);
1133 if (!cr.GetClass() || !cr->GetClassInfo())
1134 return res;
1135
1136 const std::vector<std::string>& v = gInterpreter->GetUsingNamespaces(cr->GetClassInfo());
1137 res.reserve(v.size());
1138 for (auto uid : v) {
1139 Cppyy::TCppScope_t uscope = GetScope(uid);
1140 if (uscope) res.push_back(uscope);
1141 }
1142
1143 return res;
1144}
1145
1146
1147// class reflection information ----------------------------------------------
1149{
1150 if (klass == GLOBAL_HANDLE)
1151 return "";
1152 TClassRef& cr = type_from_handle(klass);
1153 std::string clName = cr->GetName();
1154// TODO: why is this template splitting needed?
1155 std::string::size_type pos = clName.substr(0, clName.find('<')).rfind("::");
1156 if (pos != std::string::npos)
1157 return clName.substr(pos+2, std::string::npos);
1158 return clName;
1159}
1160
1162{
1163 if (klass == GLOBAL_HANDLE)
1164 return "";
1165 TClassRef& cr = type_from_handle(klass);
1166 if (cr.GetClass()) {
1167 std::string name = cr->GetName();
1169 return std::string("std::")+cr->GetName();
1170 return cr->GetName();
1171 }
1172 return "";
1173}
1174
1176{
1177 TClassRef& cr = type_from_handle(klass);
1178 if (!cr.GetClass())
1179 return false;
1180
1181 TFunction* f = cr->GetMethod(("~"+GetFinalName(klass)).c_str(), "");
1182 if (f && (f->Property() & kIsVirtual))
1183 return true;
1184
1185 return false;
1186}
1187
1189{
1190 int is_complex = 1;
1191 size_t nbases = 0;
1192
1193 TClassRef& cr = type_from_handle(klass);
1194 if (cr.GetClass() && cr->GetListOfBases() != 0)
1195 nbases = GetNumBases(klass);
1196
1197 if (1 < nbases)
1198 is_complex = 1;
1199 else if (nbases == 0)
1200 is_complex = 0;
1201 else { // one base class only
1202 TBaseClass* base = (TBaseClass*)cr->GetListOfBases()->At(0);
1203 if (base->Property() & kIsVirtualBase)
1204 is_complex = 1; // TODO: verify; can be complex, need not be.
1205 else
1206 is_complex = HasComplexHierarchy(GetScope(base->GetName()));
1207 }
1208
1209 return is_complex;
1210}
1211
1213{
1214// Get the total number of base classes that this class has.
1215 TClassRef& cr = type_from_handle(klass);
1216 if (cr.GetClass() && cr->GetListOfBases() != 0)
1217 return (TCppIndex_t)cr->GetListOfBases()->GetSize();
1218 return (TCppIndex_t)0;
1219}
1220
1222{
1223 TClassRef& cr = type_from_handle(klass);
1224 return ((TBaseClass*)cr->GetListOfBases()->At((int)ibase))->GetName();
1225}
1226
1228{
1229 if (derived == base)
1230 return true;
1231 TClassRef& derived_type = type_from_handle(derived);
1232 TClassRef& base_type = type_from_handle(base);
1233 return derived_type->GetBaseClass(base_type) != 0;
1234}
1235
1237{
1238 TClassRef& cr = type_from_handle(klass);
1239 const std::string& tn = cr->GetName();
1240 if (gSmartPtrTypes.find(tn.substr(0, tn.find("<"))) != gSmartPtrTypes.end())
1241 return true;
1242 return false;
1243}
1244
1246 const std::string& tname, TCppType_t* raw, TCppMethod_t* deref)
1247{
1248 const std::string& rn = ResolveName(tname);
1249 if (gSmartPtrTypes.find(rn.substr(0, rn.find("<"))) != gSmartPtrTypes.end()) {
1250 if (!raw && !deref) return true;
1251
1252 TClassRef& cr = type_from_handle(GetScope(tname));
1253 if (cr.GetClass()) {
1254 TFunction* func = cr->GetMethod("operator->", "");
1255 if (!func) {
1256 gInterpreter->UpdateListOfMethods(cr.GetClass());
1257 func = cr->GetMethod("operator->", "");
1258 }
1259 if (func) {
1260 if (deref) *deref = (TCppMethod_t)new_CallWrapper(func);
1261 if (raw) *raw = GetScope(TClassEdit::ShortType(
1262 func->GetReturnTypeNormalizedName().c_str(), 1));
1263 return (!deref || *deref) && (!raw || *raw);
1264 }
1265 }
1266 }
1267
1268 return false;
1269}
1270
1271void Cppyy::AddSmartPtrType(const std::string& type_name)
1272{
1273 gSmartPtrTypes.insert(ResolveName(type_name));
1274}
1275
1276
1277// type offsets --------------------------------------------------------------
1279 TCppObject_t address, int direction, bool rerror)
1280{
1281// calculate offsets between declared and actual type, up-cast: direction > 0; down-cast: direction < 0
1282 if (derived == base || !(base && derived))
1283 return (ptrdiff_t)0;
1284
1285 TClassRef& cd = type_from_handle(derived);
1286 TClassRef& cb = type_from_handle(base);
1287
1288 if (!cd.GetClass() || !cb.GetClass())
1289 return (ptrdiff_t)0;
1290
1291 ptrdiff_t offset = -1;
1292 if (!(cd->GetClassInfo() && cb->GetClassInfo())) { // gInterpreter requirement
1293 // would like to warn, but can't quite determine error from intentional
1294 // hiding by developers, so only cover the case where we really should have
1295 // had a class info, but apparently don't:
1296 if (cd->IsLoaded()) {
1297 // warn to allow diagnostics
1298 std::ostringstream msg;
1299 msg << "failed offset calculation between " << cb->GetName() << " and " << cd->GetName();
1300 // TODO: propagate this warning to caller w/o use of Python C-API
1301 // PyErr_Warn(PyExc_RuntimeWarning, const_cast<char*>(msg.str().c_str()));
1302 std::cerr << "Warning: " << msg.str() << '\n';
1303 }
1304
1305 // return -1 to signal caller NOT to apply offset
1306 return rerror ? (ptrdiff_t)offset : 0;
1307 }
1308
1309 offset = gInterpreter->ClassInfo_GetBaseOffset(
1310 cd->GetClassInfo(), cb->GetClassInfo(), (void*)address, direction > 0);
1311 if (offset == -1) // Cling error, treat silently
1312 return rerror ? (ptrdiff_t)offset : 0;
1313
1314 return (ptrdiff_t)(direction < 0 ? -offset : offset);
1315}
1316
1317
1318// method/function reflection information ------------------------------------
1320{
1321 if (IsNamespace(scope))
1322 return (TCppIndex_t)0; // enforce lazy
1323
1324 TClassRef& cr = type_from_handle(scope);
1325 if (cr.GetClass() && cr->GetListOfMethods(true)) {
1326 Cppyy::TCppIndex_t nMethods = (TCppIndex_t)cr->GetListOfMethods(false)->GetSize();
1327 if (nMethods == (TCppIndex_t)0) {
1328 std::string clName = GetScopedFinalName(scope);
1329 if (clName.find('<') != std::string::npos) {
1330 // chicken-and-egg problem: TClass does not know about methods until
1331 // instantiation, so force it
1332 if (clName.find("std::", 0, 5) == std::string::npos && \
1333 is_missclassified_stl(clName)) {
1334 // TODO: this is too simplistic for template arguments missing std::
1335 clName = "std::" + clName;
1336 }
1337 std::ostringstream stmt;
1338 stmt << "template class " << clName << ";";
1339 gInterpreter->Declare(stmt.str().c_str());
1340
1341 // now reload the methods
1342 return (TCppIndex_t)cr->GetListOfMethods(true)->GetSize();
1343 }
1344 }
1345 return nMethods;
1346 }
1347
1348 return (TCppIndex_t)0; // unknown class?
1349}
1350
1351std::vector<Cppyy::TCppIndex_t> Cppyy::GetMethodIndicesFromName(
1352 TCppScope_t scope, const std::string& name)
1353{
1354 std::vector<TCppIndex_t> indices;
1355 TClassRef& cr = type_from_handle(scope);
1356 if (cr.GetClass()) {
1357 gInterpreter->UpdateListOfMethods(cr.GetClass());
1358 int imeth = 0;
1359 TFunction* func = nullptr;
1360 TIter next(cr->GetListOfMethods());
1361 while ((func = (TFunction*)next())) {
1362 if (match_name(name, func->GetName())) {
1363 if (func->Property() & kIsPublic)
1364 indices.push_back((TCppIndex_t)imeth);
1365 }
1366 ++imeth;
1367 }
1368 } else if (scope == GLOBAL_HANDLE) {
1369 TCollection* funcs = gROOT->GetListOfGlobalFunctions(true);
1370
1371 // tickle deserialization
1372 if (!funcs->FindObject(name.c_str()))
1373 return indices;
1374
1375 TFunction* func = nullptr;
1376 TIter ifunc(funcs);
1377 while ((func = (TFunction*)ifunc.Next())) {
1378 if (match_name(name, func->GetName()))
1379 indices.push_back((TCppIndex_t)new_CallWrapper(func));
1380 }
1381 }
1382
1383 return indices;
1384}
1385
1387{
1388 TClassRef& cr = type_from_handle(scope);
1389 if (cr.GetClass()) {
1390 TFunction* f = (TFunction*)cr->GetListOfMethods(false)->At((int)idx);
1391 if (f) return (Cppyy::TCppMethod_t)new_CallWrapper(f);
1392 return (Cppyy::TCppMethod_t)nullptr;
1393 }
1394
1395 assert(klass == (Cppyy::TCppType_t)GLOBAL_HANDLE);
1396 return (Cppyy::TCppMethod_t)idx;
1397}
1398
1400{
1401 if (method) {
1402 const std::string& name = ((CallWrapper*)method)->fName;
1403
1404 if (name.compare(0, 8, "operator") != 0)
1405 // strip template instantiation part, if any
1406 return name.substr(0, name.find('<'));
1407 return name;
1408 }
1409 return "<unknown>";
1410}
1411
1413{
1414 if (method) {
1415 std::string name = ((CallWrapper*)method)->fName;
1416 name.erase(std::remove(name.begin(), name.end(), ' '), name.end());
1417 return name;
1418 }
1419 return "<unknown>";
1420}
1421
1423{
1424 if (method)
1425 return m2f(method)->GetMangledName();
1426 return "<unknown>";
1427}
1428
1430{
1431 if (method) {
1432 TFunction* f = m2f(method);
1433 if (f->ExtraProperty() & kIsConstructor)
1434 return "constructor";
1435 std::string restype = f->GetReturnTypeName();
1436 // TODO: this is ugly, but we can't use GetReturnTypeName() for ostreams
1437 // and maybe others, whereas GetReturnTypeNormalizedName() has proven to
1438 // be save in all cases (Note: 'int8_t' covers 'int8_t' and 'uint8_t')
1439 if (restype.find("int8_t") != std::string::npos)
1440 return restype;
1441 restype = f->GetReturnTypeNormalizedName();
1442 if (restype == "(lambda)") {
1443 std::ostringstream s;
1444 // TODO: what if there are parameters to the lambda?
1445 s << "__cling_internal::FT<decltype("
1446 << GetMethodFullName(method) << "(";
1447 for (Cppyy::TCppIndex_t i = 0; i < Cppyy::GetMethodNumArgs(method); ++i) {
1448 if (i != 0) s << ", ";
1449 s << Cppyy::GetMethodArgType(method, i) << "{}";
1450 }
1451 s << "))>::F";
1452 TClass* cl = TClass::GetClass(s.str().c_str());
1453 if (cl) return cl->GetName();
1454 // TODO: signal some type of error (or should that be upstream?
1455 }
1456 return restype;
1457 }
1458 return "<unknown>";
1459}
1460
1462{
1463 if (method)
1464 return m2f(method)->GetNargs();
1465 return 0;
1466}
1467
1469{
1470 if (method) {
1471 TFunction* f = m2f(method);
1472 return (TCppIndex_t)(f->GetNargs() - f->GetNargsOpt());
1473 }
1474 return (TCppIndex_t)0;
1475}
1476
1478{
1479 if (method) {
1480 TFunction* f = m2f(method);
1481 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1482 return arg->GetName();
1483 }
1484 return "<unknown>";
1485}
1486
1488{
1489 if (method) {
1490 TFunction* f = m2f(method);
1491 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1492 return arg->GetTypeNormalizedName();
1493 }
1494 return "<unknown>";
1495}
1496
1498{
1499 if (method) {
1500 TFunction* f = m2f(method);
1501 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1502 const char* def = arg->GetDefault();
1503 if (def)
1504 return def;
1505 }
1506
1507 return "";
1508}
1509
1510std::string Cppyy::GetMethodSignature(TCppMethod_t method, bool show_formalargs, TCppIndex_t maxargs)
1511{
1512 TFunction* f = m2f(method);
1513 if (f) {
1514 std::ostringstream sig;
1515 sig << "(";
1516 int nArgs = f->GetNargs();
1517 if (maxargs != (TCppIndex_t)-1) nArgs = std::min(nArgs, (int)maxargs);
1518 for (int iarg = 0; iarg < nArgs; ++iarg) {
1519 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At(iarg);
1520 sig << arg->GetFullTypeName();
1521 if (show_formalargs) {
1522 const char* argname = arg->GetName();
1523 if (argname && argname[0] != '\0') sig << " " << argname;
1524 const char* defvalue = arg->GetDefault();
1525 if (defvalue && defvalue[0] != '\0') sig << " = " << defvalue;
1526 }
1527 if (iarg != nArgs-1) sig << (show_formalargs ? ", " : ",");
1528 }
1529 sig << ")";
1530 return sig.str();
1531 }
1532 return "<unknown>";
1533}
1534
1535std::string Cppyy::GetMethodPrototype(TCppScope_t scope, TCppMethod_t method, bool show_formalargs)
1536{
1537 std::string scName = GetScopedFinalName(scope);
1538 TFunction* f = m2f(method);
1539 if (f) {
1540 std::ostringstream sig;
1541 sig << f->GetReturnTypeName() << " "
1542 << scName << "::" << f->GetName();
1543 sig << GetMethodSignature(method, show_formalargs);
1544 return sig.str();
1545 }
1546 return "<unknown>";
1547}
1548
1550{
1551 if (method) {
1552 TFunction* f = m2f(method);
1553 return f->Property() & kIsConstMethod;
1554 }
1555 return false;
1556}
1557
1559{
1560 if (scope == (TCppScope_t)GLOBAL_HANDLE) {
1561 TCollection* coll = gROOT->GetListOfFunctionTemplates();
1562 if (coll) return (TCppIndex_t)coll->GetSize();
1563 } else {
1564 TClassRef& cr = type_from_handle(scope);
1565 if (cr.GetClass()) {
1566 TCollection* coll = cr->GetListOfFunctionTemplates(true);
1567 if (coll) return (TCppIndex_t)coll->GetSize();
1568 }
1569 }
1570
1571// failure ...
1572 return (TCppIndex_t)0;
1573}
1574
1576{
1577 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1578 return ((THashList*)gROOT->GetListOfFunctionTemplates())->At((int)imeth)->GetName();
1579 else {
1580 TClassRef& cr = type_from_handle(scope);
1581 if (cr.GetClass())
1582 return cr->GetListOfFunctionTemplates(false)->At((int)imeth)->GetName();
1583 }
1584
1585// failure ...
1586 assert(!"should not be called unless GetNumTemplatedMethods() succeeded");
1587 return "";
1588}
1589
1591{
1592 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1593 return false;
1594
1595 TClassRef& cr = type_from_handle(scope);
1596 if (cr.GetClass()) {
1598 return f->ExtraProperty() & kIsConstructor;
1599 }
1600
1601 return false;
1602}
1603
1604bool Cppyy::ExistsMethodTemplate(TCppScope_t scope, const std::string& name)
1605{
1606 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1607 return (bool)gROOT->GetFunctionTemplate(name.c_str());
1608 else {
1609 TClassRef& cr = type_from_handle(scope);
1610 if (cr.GetClass())
1611 return (bool)cr->GetFunctionTemplate(name.c_str());
1612 }
1613
1614// failure ...
1615 return false;
1616}
1617
1619{
1620 TClassRef& cr = type_from_handle(scope);
1621 if (cr.GetClass()) {
1622 TFunction* f = (TFunction*)cr->GetListOfMethods(false)->At((int)idx);
1623 if (f && strstr(f->GetName(), "<")) return true;
1624 return false;
1625 }
1626
1627 assert(scope == (Cppyy::TCppType_t)GLOBAL_HANDLE);
1628 if (((CallWrapper*)idx)->fName.find('<') != std::string::npos) return true;
1629 return false;
1630}
1631
1632// helpers for Cppyy::GetMethodTemplate()
1633static std::map<TDictionary::DeclId_t, CallWrapper*> gMethodTemplates;
1634
1636 TCppScope_t scope, const std::string& name, const std::string& proto)
1637{
1638// There is currently no clean way of extracting a templated method out of ROOT/meta
1639// for a variety of reasons, none of them fundamental. The game played below is to
1640// first get any pre-existing functions already managed by ROOT/meta, but if that fails,
1641// to do an explicit lookup that ignores the prototype (i.e. the full name should be
1642// enough), and finally to ignore the template arguments part of the name as this fails
1643// in cling if there are default parameters.
1644// It would be possible to get the prototype from the created functions and use that to
1645// do a new lookup, after which ROOT/meta will manage the function. However, neither
1646// TFunction::GetPrototype() nor TFunction::GetSignature() is of the proper form, so
1647// we'll/ manage the new TFunctions instead and will assume that they are cached on the
1648// calling side to prevent multiple creations.
1649 TFunction* func = nullptr; ClassInfo_t* cl = nullptr;
1650 if (scope == (cppyy_scope_t)GLOBAL_HANDLE) {
1651 func = gROOT->GetGlobalFunctionWithPrototype(name.c_str(), proto.c_str());
1652 if (func && name.back() == '>' && name != func->GetName())
1653 func = nullptr; // happens if implicit conversion matches the overload
1654 } else {
1655 TClassRef& cr = type_from_handle(scope);
1656 if (cr.GetClass()) {
1657 func = cr->GetMethodWithPrototype(name.c_str(), proto.c_str());
1658 if (!func) {
1659 cl = cr->GetClassInfo();
1660 // try base classes to cover a common 'using' case (TODO: this is stupid and misses
1661 // out on base classes; fix that with improved access to Cling)
1662 TCppIndex_t nbases = GetNumBases(scope);
1663 for (TCppIndex_t i = 0; i < nbases; ++i) {
1664 TClassRef& base = type_from_handle(GetScope(GetBaseName(scope, i)));
1665 if (base.GetClass()) {
1666 func = base->GetMethodWithPrototype(name.c_str(), proto.c_str());
1667 if (func) break;
1668 }
1669 }
1670 }
1671 }
1672 }
1673
1674 if (!func && name.back() == '>' && (cl || scope == (cppyy_scope_t)GLOBAL_HANDLE)) {
1675 // try again, ignoring proto in case full name is complete template
1676 auto declid = gInterpreter->GetFunction(cl, name.c_str());
1677 if (declid) {
1678 auto existing = gMethodTemplates.find(declid);
1679 if (existing == gMethodTemplates.end()) {
1680 auto cw = new_CallWrapper(declid, name);
1681 existing = gMethodTemplates.insert(std::make_pair(declid, cw)).first;
1682 }
1683 return (TCppMethod_t)existing->second;
1684 }
1685 }
1686
1687 if (func) {
1688 // make sure we didn't match a non-templated overload
1689 if (func->ExtraProperty() & kIsTemplateSpec)
1690 return (TCppMethod_t)new_CallWrapper(func);
1691
1692 // disregard this non-templated method as it will be considered when appropriate
1693 return (TCppMethod_t)nullptr;
1694 }
1695
1696// try again with template arguments removed from name, if applicable
1697 if (name.back() == '>') {
1698 auto pos = name.find('<');
1699 if (pos != std::string::npos) {
1700 TCppMethod_t cppmeth = GetMethodTemplate(scope, name.substr(0, pos), proto);
1701 if (cppmeth) {
1702 // allow if requested template names match up to the result
1703 const std::string& alt = GetMethodFullName(cppmeth);
1704 if (name.size() < alt.size() && alt.find('<') == pos) {
1705 const std::string& partial = name.substr(pos, name.size()-1-pos);
1706 if (strncmp(partial.c_str(), alt.substr(pos, alt.size()-1-pos).c_str(), partial.size()) == 0)
1707 return cppmeth;
1708 }
1709 }
1710 }
1711 }
1712
1713// failure ...
1714 return (TCppMethod_t)nullptr;
1715}
1716
1717static inline
1718std::string type_remap(const std::string& n1, const std::string& n2)
1719{
1720// Operator lookups of (C++ string, Python str) should succeeded, for the combos of
1721// string/str, wstring/str, string/unicode and wstring/unicode; since C++ does not have a
1722// operator+(std::string, std::wstring), we'll have to look up the same type and rely on
1723// the converters in CPyCppyy/_cppyy.
1724 if (n1 == "str") {
1725 if (n2 == "std::basic_string<wchar_t,std::char_traits<wchar_t>,std::allocator<wchar_t> >")
1726 return n2; // match like for like
1727 return "std::string"; // probably best bet
1728 } else if (n1 == "float")
1729 return "double"; // debatable, but probably intended
1730 return n1;
1731}
1732
1734 TCppType_t scope, const std::string& lc, const std::string& rc, const std::string& opname)
1735{
1736// Find a global operator function with a matching signature; prefer by-ref, but
1737// fall back on by-value if that fails.
1738 std::string lcname1 = TClassEdit::CleanType(lc.c_str());
1739 const std::string& rcname = rc.empty() ? rc : type_remap(TClassEdit::CleanType(rc.c_str()), lcname1);
1740 const std::string& lcname = type_remap(lcname1, rcname);
1741
1742 std::string proto = lcname + "&" + (rc.empty() ? rc : (", " + rcname + "&"));
1743 if (scope == (cppyy_scope_t)GLOBAL_HANDLE) {
1744 TFunction* func = gROOT->GetGlobalFunctionWithPrototype(opname.c_str(), proto.c_str());
1745 if (func) return (TCppIndex_t)new_CallWrapper(func);
1746 proto = lcname + (rc.empty() ? rc : (", " + rcname));
1747 func = gROOT->GetGlobalFunctionWithPrototype(opname.c_str(), proto.c_str());
1748 if (func) return (TCppIndex_t)new_CallWrapper(func);
1749 } else {
1750 TClassRef& cr = type_from_handle(scope);
1751 if (cr.GetClass()) {
1752 TFunction* func = cr->GetMethodWithPrototype(opname.c_str(), proto.c_str());
1753 if (func) return (TCppIndex_t)cr->GetListOfMethods()->IndexOf(func);
1754 proto = lcname + (rc.empty() ? rc : (", " + rcname));
1755 func = cr->GetMethodWithPrototype(opname.c_str(), proto.c_str());
1756 if (func) return (TCppIndex_t)cr->GetListOfMethods()->IndexOf(func);
1757 }
1758 }
1759
1760// failure ...
1761 return (TCppIndex_t)-1;
1762}
1763
1764// method properties ---------------------------------------------------------
1766{
1767 if (method) {
1768 TFunction* f = m2f(method);
1769 return f->Property() & kIsPublic;
1770 }
1771 return false;
1772}
1773
1775{
1776 if (method) {
1777 TFunction* f = m2f(method);
1778 return f->Property() & kIsProtected;
1779 }
1780 return false;
1781}
1782
1784{
1785 if (method) {
1786 TFunction* f = m2f(method);
1787 return f->ExtraProperty() & kIsConstructor;
1788 }
1789 return false;
1790}
1791
1793{
1794 if (method) {
1795 TFunction* f = m2f(method);
1796 return f->ExtraProperty() & kIsDestructor;
1797 }
1798 return false;
1799}
1800
1802{
1803 if (method) {
1804 TFunction* f = m2f(method);
1805 return f->Property() & kIsStatic;
1806 }
1807 return false;
1808}
1809
1810// data member reflection information ----------------------------------------
1812{
1813 if (IsNamespace(scope))
1814 return (TCppIndex_t)0; // enforce lazy
1815
1816 TClassRef& cr = type_from_handle(scope);
1817 if (cr.GetClass()) {
1819 if (cr->GetListOfDataMembers())
1820 sum = cr->GetListOfDataMembers()->GetSize();
1821 if (cr->GetListOfUsingDataMembers())
1823 return sum;
1824 }
1825
1826 return (TCppIndex_t)0; // unknown class?
1827}
1828
1830{
1831 if (!cr.GetClass() || !cr->GetListOfDataMembers())
1832 return nullptr;
1833
1834 int numDMs = cr->GetListOfDataMembers()->GetSize();
1835 if ((int)idata < numDMs)
1836 return (TDataMember*)cr->GetListOfDataMembers()->At((int)idata);
1837 return (TDataMember*)cr->GetListOfUsingDataMembers()->At((int)idata - numDMs);
1838}
1839
1841{
1842 TClassRef& cr = type_from_handle(scope);
1843 if (cr.GetClass()) {
1844 TDataMember *m = GetDataMemberByIndex(cr, (int)idata);
1845 return m->GetName();
1846 }
1847 assert(scope == GLOBAL_HANDLE);
1848 TGlobal* gbl = g_globalvars[idata];
1849 return gbl->GetName();
1850}
1851
1853{
1854 if (scope == GLOBAL_HANDLE) {
1855 TGlobal* gbl = g_globalvars[idata];
1856 std::string fullType = gbl->GetFullTypeName();
1857
1858 if ((int)gbl->GetArrayDim() > 1)
1859 fullType.append("*");
1860 else if ((int)gbl->GetArrayDim() == 1) {
1861 std::ostringstream s;
1862 s << '[' << gbl->GetMaxIndex(0) << ']' << std::ends;
1863 fullType.append(s.str());
1864 }
1865 return fullType;
1866 }
1867
1868 TClassRef& cr = type_from_handle(scope);
1869 if (cr.GetClass()) {
1870 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
1871 // TODO: fix this upstream. Usually, we want m->GetFullTypeName(), because it does
1872 // not resolve typedefs, but it looses scopes for inner classes/structs, so in that
1873 // case m->GetTrueTypeName() should be used (this also cleans up the cases where
1874 // the "full type" retains spurious "struct" or "union" in the name).
1875 std::string fullType = m->GetFullTypeName();
1876 if (fullType != m->GetTrueTypeName()) {
1877 const std::string& trueName = m->GetTrueTypeName();
1878 if (fullType.find("::") == std::string::npos && trueName.find("::") != std::string::npos)
1879 fullType = trueName;
1880 }
1881
1882 if ((int)m->GetArrayDim() > 1 || (!m->IsBasic() && m->IsaPointer()))
1883 fullType.append("*");
1884 else if ((int)m->GetArrayDim() == 1) {
1885 std::ostringstream s;
1886 s << '[' << m->GetMaxIndex(0) << ']' << std::ends;
1887 fullType.append(s.str());
1888 }
1889 return fullType;
1890 }
1891
1892 return "<unknown>";
1893}
1894
1896{
1897 if (scope == GLOBAL_HANDLE) {
1898 TGlobal* gbl = g_globalvars[idata];
1899 if (!gbl->GetAddress() || gbl->GetAddress() == (void*)-1) {
1900 // CLING WORKAROUND: make sure variable is loaded
1901 intptr_t addr = (intptr_t)gInterpreter->ProcessLine((std::string("&")+gbl->GetName()+";").c_str());
1902 if (gbl->GetAddress() && gbl->GetAddress() != (void*)-1)
1903 return (intptr_t)gbl->GetAddress(); // now loaded!
1904 return addr; // last resort ...
1905 }
1906 return (intptr_t)gbl->GetAddress();
1907 }
1908
1909 TClassRef& cr = type_from_handle(scope);
1910 if (cr.GetClass()) {
1911 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
1912 // CLING WORKAROUND: the following causes templates to be instantiated first within the proper
1913 // scope, making the lookup succeed and preventing spurious duplicate instantiations later. Also,
1914 // if the variable is not yet loaded, pull it in through gInterpreter.
1915 if (m->Property() & kIsStatic) {
1916 if (strchr(cr->GetName(), '<'))
1917 gInterpreter->ProcessLine(((std::string)cr->GetName()+"::"+m->GetName()+";").c_str());
1918 if ((intptr_t)m->GetOffsetCint() == (intptr_t)-1)
1919 return (intptr_t)gInterpreter->ProcessLine((std::string("&")+cr->GetName()+"::"+m->GetName()+";").c_str());
1920 }
1921 return (intptr_t)m->GetOffsetCint(); // yes, CINT (GetOffset() is both wrong
1922 // and caches that wrong result!
1923 }
1924
1925 return (intptr_t)-1;
1926}
1927
1929{
1930 if (scope == GLOBAL_HANDLE) {
1931 TGlobal* gb = (TGlobal*)gROOT->GetListOfGlobals(false /* load */)->FindObject(name.c_str());
1932 if (!gb) gb = (TGlobal*)gROOT->GetListOfGlobals(true /* load */)->FindObject(name.c_str());
1933 if (!gb) {
1934 // some enums are not loaded as they are not considered part of
1935 // the global scope, but of the enum scope; get them w/o checking
1936 TDictionary::DeclId_t did = gInterpreter->GetDataMember(nullptr, name.c_str());
1937 if (did) {
1938 DataMemberInfo_t* t = gInterpreter->DataMemberInfo_Factory(did, nullptr);
1939 ((TListOfDataMembers*)gROOT->GetListOfGlobals())->Get(t, true);
1940 gb = (TGlobal*)gROOT->GetListOfGlobals(false /* load */)->FindObject(name.c_str());
1941 }
1942 }
1943
1944 if (gb && strcmp(gb->GetFullTypeName(), "(lambda)") == 0) {
1945 // lambdas use a compiler internal closure type, so we wrap
1946 // them, then return the wrapper's type
1947 // TODO: this current leaks the std::function; also, if possible,
1948 // should instantiate through TClass rather then ProcessLine
1949 std::ostringstream s;
1950 s << "auto __cppyy_internal_wrap_" << name << " = "
1951 "new __cling_internal::FT<decltype(" << name << ")>::F"
1952 "{" << name << "};";
1953 gInterpreter->ProcessLine(s.str().c_str());
1954 TGlobal* wrap = (TGlobal*)gROOT->GetListOfGlobals(true)->FindObject(
1955 ("__cppyy_internal_wrap_"+name).c_str());
1956 if (wrap && wrap->GetAddress()) gb = wrap;
1957 }
1958
1959 if (gb) {
1960 // TODO: do we ever need a reverse lookup?
1961 g_globalvars.push_back(gb);
1962 return TCppIndex_t(g_globalvars.size() - 1);
1963 }
1964
1965 } else {
1966 TClassRef& cr = type_from_handle(scope);
1967 if (cr.GetClass()) {
1968 TDataMember* dm =
1970 // TODO: turning this into an index is silly ...
1971 if (dm) return (TCppIndex_t)cr->GetListOfDataMembers()->IndexOf(dm);
1973 if (dm)
1974 return (TCppIndex_t)cr->GetListOfDataMembers()->IndexOf(dm)
1975 + cr->GetListOfDataMembers()->GetSize();
1976 }
1977 }
1978
1979 return (TCppIndex_t)-1;
1980}
1981
1982
1983// data member properties ----------------------------------------------------
1985{
1986 if (scope == GLOBAL_HANDLE)
1987 return true;
1988 TClassRef& cr = type_from_handle(scope);
1989 if (cr->Property() & kIsNamespace)
1990 return true;
1991 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
1992 return m->Property() & kIsPublic;
1993}
1994
1996{
1997 if (scope == GLOBAL_HANDLE)
1998 return true;
1999 TClassRef& cr = type_from_handle(scope);
2000 if (cr->Property() & kIsNamespace)
2001 return true;
2002 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2003 return m->Property() & kIsProtected;
2004}
2005
2007{
2008 if (scope == GLOBAL_HANDLE)
2009 return true;
2010 TClassRef& cr = type_from_handle(scope);
2011 if (cr->Property() & kIsNamespace)
2012 return true;
2013 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2014 return m->Property() & kIsStatic;
2015}
2016
2018{
2019 if (scope == GLOBAL_HANDLE) {
2020 TGlobal* gbl = g_globalvars[idata];
2021 return gbl->Property() & kIsConstant;
2022 }
2023 TClassRef& cr = type_from_handle(scope);
2024 if (cr.GetClass()) {
2025 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2026 return m->Property() & kIsConstant;
2027 }
2028 return false;
2029}
2030
2032{
2033// TODO: currently, ROOT/meta does not properly distinguish between variables of enum
2034// type, and values of enums. The latter are supposed to be const. This code relies on
2035// odd features (bugs?) to figure out the difference, but this should really be fixed
2036// upstream and/or deserves a new API.
2037
2038 if (scope == GLOBAL_HANDLE) {
2039 TGlobal* gbl = g_globalvars[idata];
2040
2041 // make use of an oddity: enum global variables do not have their kIsStatic bit
2042 // set, whereas enum global values do
2043 return (gbl->Property() & kIsEnum) && (gbl->Property() & kIsStatic);
2044 }
2045
2046 TClassRef& cr = type_from_handle(scope);
2047 if (cr.GetClass()) {
2048 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2049 std::string ti = m->GetTypeName();
2050
2051 // can't check anonymous enums by type name, so just accept them as enums
2052 if (ti.rfind("(anonymous)") != std::string::npos)
2053 return m->Property() & kIsEnum;
2054
2055 // since there seems to be no distinction between data of enum type and enum values,
2056 // check the list of constants for the type to see if there's a match
2057 if (ti.rfind(cr->GetName(), 0) != std::string::npos) {
2058 std::string::size_type s = strlen(cr->GetName())+2;
2059 if (s < ti.size()) {
2060 TEnum* ee = ((TListOfEnums*)cr->GetListOfEnums())->GetObject(ti.substr(s, std::string::npos).c_str());
2061 if (ee) return ee->GetConstant(m->GetName());
2062 }
2063 }
2064 }
2065
2066// this default return only means that the data will be writable, not that it will
2067// be unreadable or otherwise misrepresented
2068 return false;
2069}
2070
2071int Cppyy::GetDimensionSize(TCppScope_t scope, TCppIndex_t idata, int dimension)
2072{
2073 if (scope == GLOBAL_HANDLE) {
2074 TGlobal* gbl = g_globalvars[idata];
2075 return gbl->GetMaxIndex(dimension);
2076 }
2077 TClassRef& cr = type_from_handle(scope);
2078 if (cr.GetClass()) {
2079 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2080 return m->GetMaxIndex(dimension);
2081 }
2082 return -1;
2083}
2084
2085
2086// enum properties -----------------------------------------------------------
2087Cppyy::TCppEnum_t Cppyy::GetEnum(TCppScope_t scope, const std::string& enum_name)
2088{
2089 if (scope == GLOBAL_HANDLE)
2090 return (TCppEnum_t)gROOT->GetListOfEnums(kTRUE)->FindObject(enum_name.c_str());
2091
2092 TClassRef& cr = type_from_handle(scope);
2093 if (cr.GetClass())
2094 return (TCppEnum_t)cr->GetListOfEnums(kTRUE)->FindObject(enum_name.c_str());
2095
2096 return (TCppEnum_t)0;
2097}
2098
2100{
2101 return (TCppIndex_t)((TEnum*)etype)->GetConstants()->GetSize();
2102}
2103
2105{
2106 return ((TEnumConstant*)((TEnum*)etype)->GetConstants()->At(idata))->GetName();
2107}
2108
2110{
2111 TEnumConstant* ecst = (TEnumConstant*)((TEnum*)etype)->GetConstants()->At(idata);
2112 return (long long)ecst->GetValue();
2113}
2114
2115
2116//- C-linkage wrappers -------------------------------------------------------
2117
2118extern "C" {
2119/* direct interpreter access ---------------------------------------------- */
2120int cppyy_compile(const char* code) {
2121 return Cppyy::Compile(code);
2122}
2123
2124
2125/* name to opaque C++ scope representation -------------------------------- */
2126char* cppyy_resolve_name(const char* cppitem_name) {
2127 return cppstring_to_cstring(Cppyy::ResolveName(cppitem_name));
2128}
2129
2130char* cppyy_resolve_enum(const char* enum_type) {
2131 return cppstring_to_cstring(Cppyy::ResolveEnum(enum_type));
2132}
2133
2134cppyy_scope_t cppyy_get_scope(const char* scope_name) {
2135 return cppyy_scope_t(Cppyy::GetScope(scope_name));
2136}
2137
2139 return cppyy_type_t(Cppyy::GetActualClass(klass, (void*)obj));
2140}
2141
2143 return Cppyy::SizeOf(klass);
2144}
2145
2146size_t cppyy_size_of_type(const char* type_name) {
2147 return Cppyy::SizeOf(type_name);
2148}
2149
2150
2151/* memory management ------------------------------------------------------ */
2154}
2155
2157 Cppyy::Deallocate(type, (void*)self);
2158}
2159
2162}
2163
2165 Cppyy::Destruct(type, (void*)self);
2166}
2167
2168
2169/* method/function dispatching -------------------------------------------- */
2170/* Exception types:
2171 1: default (unknown exception)
2172 2: standard exception
2173*/
2174#define CPPYY_HANDLE_EXCEPTION \
2175 catch (std::exception& e) { \
2176 cppyy_exctype_t* etype = (cppyy_exctype_t*)((Parameter*)args+nargs); \
2177 *etype = (cppyy_exctype_t)2; \
2178 *((char**)(etype+1)) = cppstring_to_cstring(e.what()); \
2179 } \
2180 catch (...) { \
2181 cppyy_exctype_t* etype = (cppyy_exctype_t*)((Parameter*)args+nargs); \
2182 *etype = (cppyy_exctype_t)1; \
2183 *((char**)(etype+1)) = \
2184 cppstring_to_cstring("unhandled, unknown C++ exception"); \
2185 }
2186
2187void cppyy_call_v(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2188 try {
2189 Cppyy::CallV(method, (void*)self, nargs, args);
2191}
2192
2193unsigned char cppyy_call_b(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2194 try {
2195 return (unsigned char)Cppyy::CallB(method, (void*)self, nargs, args);
2197 return (unsigned char)-1;
2198}
2199
2200char cppyy_call_c(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2201 try {
2202 return (char)Cppyy::CallC(method, (void*)self, nargs, args);
2204 return (char)-1;
2205}
2206
2207short cppyy_call_h(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2208 try {
2209 return (short)Cppyy::CallH(method, (void*)self, nargs, args);
2211 return (short)-1;
2212}
2213
2214int cppyy_call_i(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2215 try {
2216 return (int)Cppyy::CallI(method, (void*)self, nargs, args);
2218 return (int)-1;
2219}
2220
2221long cppyy_call_l(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2222 try {
2223 return (long)Cppyy::CallL(method, (void*)self, nargs, args);
2225 return (long)-1;
2226}
2227
2228long long cppyy_call_ll(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2229 try {
2230 return (long long)Cppyy::CallLL(method, (void*)self, nargs, args);
2232 return (long long)-1;
2233}
2234
2235float cppyy_call_f(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2236 try {
2237 return (float)Cppyy::CallF(method, (void*)self, nargs, args);
2239 return (float)-1;
2240}
2241
2242double cppyy_call_d(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2243 try {
2244 return (double)Cppyy::CallD(method, (void*)self, nargs, args);
2246 return (double)-1;
2247}
2248
2249long double cppyy_call_ld(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2250 try {
2251 return (long double)Cppyy::CallLD(method, (void*)self, nargs, args);
2253 return (long double)-1;
2254}
2255
2256double cppyy_call_nld(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2257 return (double)cppyy_call_ld(method, self, nargs, args);
2258}
2259
2260void* cppyy_call_r(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2261 try {
2262 return (void*)Cppyy::CallR(method, (void*)self, nargs, args);
2264 return (void*)nullptr;
2265}
2266
2268 cppyy_method_t method, cppyy_object_t self, int nargs, void* args, size_t* lsz) {
2269 try {
2270 return Cppyy::CallS(method, (void*)self, nargs, args, lsz);
2272 return (char*)nullptr;
2273}
2274
2276 cppyy_method_t method, cppyy_type_t klass, int nargs, void* args) {
2277 try {
2278 return cppyy_object_t(Cppyy::CallConstructor(method, klass, nargs, args));
2280 return (cppyy_object_t)0;
2281}
2282
2284 Cppyy::CallDestructor(klass, self);
2285}
2286
2288 int nargs, void* args, cppyy_type_t result_type) {
2289 try {
2290 return cppyy_object_t(Cppyy::CallO(method, (void*)self, nargs, args, result_type));
2292 return (cppyy_object_t)0;
2293}
2294
2296 return cppyy_funcaddr_t(Cppyy::GetFunctionAddress(method, true));
2297}
2298
2299
2300/* handling of function argument buffer ----------------------------------- */
2302// for calls through C interface, require extra space for reporting exceptions
2303 return malloc(nargs*sizeof(Parameter)+sizeof(cppyy_exctype_t)+sizeof(char**));
2304}
2305
2307 free(args);
2308}
2309
2311 return (size_t)Cppyy::GetFunctionArgSizeof();
2312}
2313
2315 return (size_t)Cppyy::GetFunctionArgTypeoffset();
2316}
2317
2318
2319/* scope reflection information ------------------------------------------- */
2321 return (int)Cppyy::IsNamespace(scope);
2322}
2323
2324int cppyy_is_template(const char* template_name) {
2325 return (int)Cppyy::IsTemplate(template_name);
2326}
2327
2329 return (int)Cppyy::IsAbstract(type);
2330}
2331
2332int cppyy_is_enum(const char* type_name) {
2333 return (int)Cppyy::IsEnum(type_name);
2334}
2335
2336const char** cppyy_get_all_cpp_names(cppyy_scope_t scope, size_t* count) {
2337 std::set<std::string> cppnames;
2338 Cppyy::GetAllCppNames(scope, cppnames);
2339 const char** c_cppnames = (const char**)malloc(cppnames.size()*sizeof(const char*));
2340 int i = 0;
2341 for (const auto& name : cppnames) {
2342 c_cppnames[i] = cppstring_to_cstring(name);
2343 ++i;
2344 }
2345 *count = cppnames.size();
2346 return c_cppnames;
2347}
2348
2349
2350/* namespace reflection information --------------------------------------- */
2352 const std::vector<Cppyy::TCppScope_t>& uv = Cppyy::GetUsingNamespaces((Cppyy::TCppScope_t)scope);
2353
2354 if (uv.empty())
2355 return (cppyy_index_t*)nullptr;
2356
2357 cppyy_scope_t* llresult = (cppyy_scope_t*)malloc(sizeof(cppyy_scope_t)*(uv.size()+1));
2358 for (int i = 0; i < (int)uv.size(); ++i) llresult[i] = uv[i];
2359 llresult[uv.size()] = (cppyy_scope_t)0;
2360 return llresult;
2361}
2362
2363
2364/* class reflection information ------------------------------------------- */
2367}
2368
2371}
2372
2374 return (int)Cppyy::HasVirtualDestructor(type);
2375}
2376
2378 return (int)Cppyy::HasComplexHierarchy(type);
2379}
2380
2382 return (int)Cppyy::GetNumBases(type);
2383}
2384
2385char* cppyy_base_name(cppyy_type_t type, int base_index) {
2386 return cppstring_to_cstring(Cppyy::GetBaseName (type, base_index));
2387}
2388
2390 return (int)Cppyy::IsSubtype(derived, base);
2391}
2392
2394 return (int)Cppyy::IsSmartPtr(type);
2395}
2396
2397int cppyy_smartptr_info(const char* name, cppyy_type_t* raw, cppyy_method_t* deref) {
2398 return (int)Cppyy::GetSmartPtrInfo(name, raw, deref);
2399}
2400
2401void cppyy_add_smartptr_type(const char* type_name) {
2402 Cppyy::AddSmartPtrType(type_name);
2403}
2404
2405
2406/* calculate offsets between declared and actual type, up-cast: direction > 0; down-cast: direction < 0 */
2407ptrdiff_t cppyy_base_offset(cppyy_type_t derived, cppyy_type_t base, cppyy_object_t address, int direction) {
2408 return (ptrdiff_t)Cppyy::GetBaseOffset(derived, base, (void*)address, direction, 0);
2409}
2410
2411
2412/* method/function reflection information --------------------------------- */
2414 return (int)Cppyy::GetNumMethods(scope);
2415}
2416
2418{
2419 std::vector<cppyy_index_t> result = Cppyy::GetMethodIndicesFromName(scope, name);
2420
2421 if (result.empty())
2422 return (cppyy_index_t*)nullptr;
2423
2424 cppyy_index_t* llresult = (cppyy_index_t*)malloc(sizeof(cppyy_index_t)*(result.size()+1));
2425 for (int i = 0; i < (int)result.size(); ++i) llresult[i] = result[i];
2426 llresult[result.size()] = -1;
2427 return llresult;
2428}
2429
2431 return cppyy_method_t(Cppyy::GetMethod(scope, idx));
2432}
2433
2436}
2437
2440}
2441
2444}
2445
2448}
2449
2451 return (int)Cppyy::GetMethodNumArgs((Cppyy::TCppMethod_t)method);
2452}
2453
2455 return (int)Cppyy::GetMethodReqArgs((Cppyy::TCppMethod_t)method);
2456}
2457
2458char* cppyy_method_arg_name(cppyy_method_t method, int arg_index) {
2460}
2461
2462char* cppyy_method_arg_type(cppyy_method_t method, int arg_index) {
2464}
2465
2466char* cppyy_method_arg_default(cppyy_method_t method, int arg_index) {
2468}
2469
2470char* cppyy_method_signature(cppyy_method_t method, int show_formalargs) {
2471 return cppstring_to_cstring(Cppyy::GetMethodSignature((Cppyy::TCppMethod_t)method, (bool)show_formalargs));
2472}
2473
2474char* cppyy_method_signature_max(cppyy_method_t method, int show_formalargs, int maxargs) {
2475 return cppstring_to_cstring(Cppyy::GetMethodSignature((Cppyy::TCppMethod_t)method, (bool)show_formalargs, (Cppyy::TCppIndex_t)maxargs));
2476}
2477
2478char* cppyy_method_prototype(cppyy_scope_t scope, cppyy_method_t method, int show_formalargs) {
2480 (Cppyy::TCppScope_t)scope, (Cppyy::TCppMethod_t)method, (bool)show_formalargs));
2481}
2482
2484 return (int)Cppyy::IsConstMethod(method);
2485}
2486
2488 return (int)Cppyy::GetNumTemplatedMethods(scope);
2489}
2490
2493}
2494
2497}
2498
2500 return (int)Cppyy::ExistsMethodTemplate(scope, name);
2501}
2502
2504 return (int)Cppyy::IsMethodTemplate(scope, idx);
2505}
2506
2509}
2510
2513}
2514
2515
2516/* method properties ------------------------------------------------------ */
2518 return (int)Cppyy::IsPublicMethod((Cppyy::TCppMethod_t)method);
2519}
2520
2522 return (int)Cppyy::IsProtectedMethod((Cppyy::TCppMethod_t)method);
2523}
2524
2526 return (int)Cppyy::IsConstructor((Cppyy::TCppMethod_t)method);
2527}
2528
2530 return (int)Cppyy::IsDestructor((Cppyy::TCppMethod_t)method);
2531}
2532
2534 return (int)Cppyy::IsStaticMethod((Cppyy::TCppMethod_t)method);
2535}
2536
2537
2538/* data member reflection information ------------------------------------- */
2540 return (int)Cppyy::GetNumDatamembers(scope);
2541}
2542
2543char* cppyy_datamember_name(cppyy_scope_t scope, int datamember_index) {
2544 return cppstring_to_cstring(Cppyy::GetDatamemberName(scope, datamember_index));
2545}
2546
2547char* cppyy_datamember_type(cppyy_scope_t scope, int datamember_index) {
2548 return cppstring_to_cstring(Cppyy::GetDatamemberType(scope, datamember_index));
2549}
2550
2551intptr_t cppyy_datamember_offset(cppyy_scope_t scope, int datamember_index) {
2552 return intptr_t(Cppyy::GetDatamemberOffset(scope, datamember_index));
2553}
2554
2556 return (int)Cppyy::GetDatamemberIndex(scope, name);
2557}
2558
2559
2560
2561/* data member properties ------------------------------------------------- */
2563 return (int)Cppyy::IsPublicData(type, datamember_index);
2564}
2565
2567 return (int)Cppyy::IsProtectedData(type, datamember_index);
2568}
2569
2571 return (int)Cppyy::IsStaticData(type, datamember_index);
2572}
2573
2575 return (int)Cppyy::IsConstData(scope, idata);
2576}
2577
2579 return (int)Cppyy::IsEnumData(scope, idata);
2580}
2581
2582int cppyy_get_dimension_size(cppyy_scope_t scope, cppyy_index_t idata, int dimension) {
2583 return Cppyy::GetDimensionSize(scope, idata, dimension);
2584}
2585
2586
2587/* misc helpers ----------------------------------------------------------- */
2589void* cppyy_load_dictionary(const char* lib_name) {
2590 int result = gSystem->Load(lib_name);
2591 return (void*)(result == 0 /* success */ || result == 1 /* already loaded */);
2592}
2593
2594#if defined(_MSC_VER)
2595long long cppyy_strtoll(const char* str) {
2596 return _strtoi64(str, NULL, 0);
2597}
2598
2599extern "C" {
2600unsigned long long cppyy_strtoull(const char* str) {
2601 return _strtoui64(str, NULL, 0);
2602}
2603}
2604#else
2605long long cppyy_strtoll(const char* str) {
2606 return strtoll(str, NULL, 0);
2607}
2608
2609extern "C" {
2610unsigned long long cppyy_strtoull(const char* str) {
2611 return strtoull(str, NULL, 0);
2612}
2613}
2614#endif
2615
2616void cppyy_free(void* ptr) {
2617 free(ptr);
2618}
2619
2620cppyy_object_t cppyy_charp2stdstring(const char* str, size_t sz) {
2621 return (cppyy_object_t)new std::string(str, sz);
2622}
2623
2624const char* cppyy_stdstring2charp(cppyy_object_t ptr, size_t* lsz) {
2625 *lsz = ((std::string*)ptr)->size();
2626 return ((std::string*)ptr)->data();
2627}
2628
2630 return (cppyy_object_t)new std::string(*(std::string*)ptr);
2631}
2632
2634 return (double)*(long double*)p;
2635}
2636
2637void cppyy_double2longdouble(double d, void* p) {
2638 *(long double*)p = d;
2639}
2640
2642 return (int)(*(std::vector<bool>*)ptr)[idx];
2643}
2644
2645void cppyy_vectorbool_setitem(cppyy_object_t ptr, int idx, int value) {
2646 (*(std::vector<bool>*)ptr)[idx] = (bool)value;
2647}
2648
2649} // end C-linkage wrappers
long
ROOT::R::TRInterface & r
Definition Object.C:4
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
int Int_t
Definition RtypesCore.h:45
long double LongDouble_t
Definition RtypesCore.h:61
long long Long64_t
Definition RtypesCore.h:73
const Bool_t kTRUE
Definition RtypesCore.h:91
@ kMAXSIGNALS
Definition Rtypes.h:59
@ kOther_t
Definition TDataType.h:32
@ kIsDestructor
@ kIsTemplateSpec
@ kIsConstructor
@ kClassHasExplicitDtor
@ kClassHasImplicitDtor
@ kIsPublic
Definition TDictionary.h:75
@ kIsConstant
Definition TDictionary.h:88
@ kIsConstMethod
Definition TDictionary.h:96
@ kIsEnum
Definition TDictionary.h:68
@ kIsPrivate
Definition TDictionary.h:77
@ kIsFundamental
Definition TDictionary.h:70
@ kIsAbstract
Definition TDictionary.h:71
@ kIsStatic
Definition TDictionary.h:80
@ kIsProtected
Definition TDictionary.h:76
@ kIsVirtual
Definition TDictionary.h:72
@ kIsNamespace
Definition TDictionary.h:95
@ kIsVirtualBase
Definition TDictionary.h:89
const Int_t kFatal
Definition TError.h:51
R__EXTERN Int_t gErrorIgnoreLevel
Definition TError.h:129
R__EXTERN TExceptionHandler * gExceptionHandler
Definition TException.h:84
R__EXTERN ExceptionContext_t * gException
Definition TException.h:74
R__EXTERN void Throw(int code)
If an exception context has been set (using the TRY and RETRY macros) jump back to where it was set.
char name[80]
Definition TGX11.cxx:110
int type
Definition TGX11.cxx:121
#define gInterpreter
#define gROOT
Definition TROOT.h:406
typedef void((*Func_t)())
R__EXTERN TSystem * gSystem
Definition TSystem.h:559
static struct Signalmap_t gSignalMap[kMAXSIGNALS]
size_t cppyy_scope_t
Definition capi.h:12
cppyy_scope_t cppyy_type_t
Definition capi.h:13
intptr_t cppyy_method_t
Definition capi.h:15
size_t cppyy_index_t
Definition capi.h:17
void * cppyy_object_t
Definition capi.h:14
unsigned long cppyy_exctype_t
Definition capi.h:20
void * cppyy_funcaddr_t
Definition capi.h:18
const char * proto
Definition civetweb.c:16604
#define free
Definition civetweb.c:1539
#define malloc
Definition civetweb.c:1536
Each class (see TClass) has a linked list of its base class(es).
Definition TBaseClass.h:33
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
TClassRef is used to implement a permanent reference to a TClass object.
Definition TClassRef.h:28
TClass * GetClass() const
Definition TClassRef.h:70
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:80
TList * GetListOfUsingDataMembers(Bool_t load=kTRUE)
Return list containing the TDataMembers of using declarations of a class.
Definition TClass.cxx:3763
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:4955
TMethod * GetMethod(const char *method, const char *params, Bool_t objectIsConst=kFALSE)
Find the best method (if there is one) matching the parameters.
Definition TClass.cxx:4388
TMethod * GetMethodWithPrototype(const char *method, const char *proto, Bool_t objectIsConst=kFALSE, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Find the method with a given prototype.
Definition TClass.cxx:4433
const TList * GetListOfAllPublicMethods(Bool_t load=kTRUE)
Returns a list of all public methods of this class and its base classes.
Definition TClass.cxx:3822
void Destructor(void *obj, Bool_t dtorOnly=kFALSE)
Explicitly call destructor for object.
Definition TClass.cxx:5377
TList * GetListOfFunctionTemplates(Bool_t load=kTRUE)
Return TListOfFunctionTemplates for a class.
Definition TClass.cxx:3775
TList * GetListOfEnums(Bool_t load=kTRUE)
Return a list containing the TEnums of a class.
Definition TClass.cxx:3663
TList * GetListOfMethods(Bool_t load=kTRUE)
Return list containing the TMethods of a class.
Definition TClass.cxx:3789
TClass * GetBaseClass(const char *classname)
Return pointer to the base class "classname".
Definition TClass.cxx:2644
TList * GetListOfDataMembers(Bool_t load=kTRUE)
Return list containing the TDataMembers of a class.
Definition TClass.cxx:3747
TList * GetListOfBases()
Return list containing the TBaseClass(es) of a class.
Definition TClass.cxx:3613
Bool_t IsLoaded() const
Return true if the shared library of this class is currently in the a process's memory.
Definition TClass.cxx:5889
ClassInfo_t * GetClassInfo() const
Definition TClass.h:430
Long_t Property() const
Returns the properties of the TClass as a bit field stored as a Long_t value.
Definition TClass.cxx:6063
Long_t ClassProperty() const
Return the C++ property of this class, eg.
Definition TClass.cxx:2386
ROOT::DelFunc_t GetDelete() const
Return the wrapper around delete ThiObject.
Definition TClass.cxx:7437
TClass * GetActualClass(const void *object) const
Return a pointer the the real class of the object.
Definition TClass.cxx:2597
TFunctionTemplate * GetFunctionTemplate(const char *name)
Definition TClass.cxx:3584
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:2957
Collection abstract base class.
Definition TCollection.h:63
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
All ROOT classes may have RTTI (run time type identification) support added.
Definition TDataMember.h:31
Basic data type descriptor (datatype information is obtained from CINT).
Definition TDataType.h:44
Int_t GetType() const
Definition TDataType.h:68
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
const char * GetFullTypeName() const
Get full type description of typedef, e,g.: "class TDirectory*".
Int_t Size() const
Get size of basic typedef'ed type.
const void * DeclId_t
The TEnumConstant class implements the constants of the enum type.
Long64_t GetValue() const
The TEnum class implements the enum type.
Definition TEnum.h:33
const TEnumConstant * GetConstant(const char *name) const
Definition TEnum.h:61
Definition TEnv.h:87
const char * GetValue() const
Definition TEnv.h:111
const char * GetName() const
Returns name of object.
Definition TEnv.h:110
virtual void HandleException(int sig)=0
Dictionary for function template This class describes one single function template.
Global functions class (global functions are obtained from CINT).
Definition TFunction.h:30
virtual const char * GetMangledName() const
Returns the mangled name as defined by CINT, or 0 in case of error.
Long_t Property() const
Get property description word. For meaning of bits see EProperty.
Int_t GetNargs() const
Number of function arguments.
Long_t ExtraProperty() const
Get property description word. For meaning of bits see EProperty.
std::string GetReturnTypeNormalizedName() const
Get the normalized name of the return type.
Global variables class (global variables are obtained from CINT).
Definition TGlobal.h:28
virtual Int_t GetArrayDim() const
Return number of array dimensions.
Definition TGlobal.cxx:85
virtual Long_t Property() const
Get property description word. For meaning of bits see EProperty.
Definition TGlobal.cxx:148
virtual Int_t GetMaxIndex(Int_t dim) const
Return maximum index for array dimension "dim".
Definition TGlobal.cxx:101
virtual void * GetAddress() const
Return address of global.
Definition TGlobal.cxx:77
virtual const char * GetFullTypeName() const
Get full type description of global variable, e,g.: "class TDirectory*".
Definition TGlobal.cxx:120
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition THashList.h:34
TObject * Next()
A collection of TDataMember objects designed for fast access given a DeclId_t and for keep track of T...
A collection of TEnum objects designed for fast access given a DeclId_t and for keep track of TEnum t...
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition TList.cxx:578
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
Each ROOT method (see TMethod) has a linked list of its arguments.
Definition TMethodArg.h:36
const char * GetFullTypeName() const
Get full type description of method argument, e.g.: "class TDirectory*".
const char * GetDefault() const
Get default value of method argument.
std::string GetTypeNormalizedName() const
Get the normalized name of the return type.
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:359
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:323
static Bool_t Initialized()
Return kTRUE if the TROOT object has been initialized.
Definition TROOT.cxx:2826
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition TSystem.cxx:1853
virtual void Exit(int code, Bool_t mode=kTRUE)
Exit the application.
Definition TSystem.cxx:717
virtual void StackTrace()
Print a stack trace.
Definition TSystem.cxx:733
char * cppyy_method_mangled_name(cppyy_method_t method)
int cppyy_is_staticdata(cppyy_type_t type, cppyy_index_t datamember_index)
int cppyy_vectorbool_getitem(cppyy_object_t ptr, int idx)
int cppyy_has_virtual_destructor(cppyy_type_t type)
int cppyy_is_publicmethod(cppyy_method_t method)
void * cppyy_call_r(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
cppyy_object_t cppyy_constructor(cppyy_method_t method, cppyy_type_t klass, int nargs, void *args)
static void cond_add(Cppyy::TCppScope_t scope, const std::string &ns_scope, std::set< std::string > &cppnames, const char *name, bool nofilter=false)
size_t cppyy_function_arg_typeoffset()
int cppyy_exists_method_template(cppyy_scope_t scope, const char *name)
int cppyy_get_num_templated_methods(cppyy_scope_t scope)
char * cppyy_call_s(cppyy_method_t method, cppyy_object_t self, int nargs, void *args, size_t *lsz)
char * cppyy_resolve_enum(const char *enum_type)
cppyy_object_t cppyy_charp2stdstring(const char *str, size_t sz)
static Name2ClassRefIndex_t g_name2classrefidx
int cppyy_is_constructor(cppyy_method_t method)
int cppyy_is_enum_data(cppyy_scope_t scope, cppyy_index_t idata)
#define FILL_COLL(type, filter)
char * cppyy_scoped_final_name(cppyy_type_t type)
int cppyy_num_bases(cppyy_type_t type)
const int SMALL_ARGS_N
static std::string outer_with_template(const std::string &name)
cppyy_method_t cppyy_get_method_template(cppyy_scope_t scope, const char *name, const char *proto)
static std::string outer_no_template(const std::string &name)
int cppyy_is_smartptr(cppyy_type_t type)
static std::string type_remap(const std::string &n1, const std::string &n2)
static std::map< Cppyy::TCppType_t, bool > sHasOperatorDelete
long long cppyy_strtoll(const char *str)
static std::set< std::string > gRootSOs
int cppyy_is_subtype(cppyy_type_t derived, cppyy_type_t base)
static TClassRef & type_from_handle(Cppyy::TCppScope_t scope)
double cppyy_call_nld(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
void cppyy_free(void *ptr)
int cppyy_smartptr_info(const char *name, cppyy_type_t *raw, cppyy_method_t *deref)
unsigned long long cppyy_strtoull(const char *str)
char * cppyy_method_arg_type(cppyy_method_t method, int arg_index)
cppyy_object_t cppyy_stdstring2stdstring(cppyy_object_t ptr)
char * cppyy_base_name(cppyy_type_t type, int base_index)
cppyy_scope_t * cppyy_get_using_namespaces(cppyy_scope_t scope)
static bool copy_args(Parameter *args, size_t nargs, void **vargs)
static TInterpreter::CallFuncIFacePtr_t GetCallFunc(Cppyy::TCppMethod_t method)
int cppyy_is_publicdata(cppyy_type_t type, cppyy_index_t datamember_index)
CPyCppyy::Parameter Parameter
intptr_t cppyy_datamember_offset(cppyy_scope_t scope, int datamember_index)
char * cppyy_datamember_type(cppyy_scope_t scope, int datamember_index)
static bool is_missclassified_stl(const std::string &name)
std::map< std::string, ClassRefs_t::size_type > Name2ClassRefIndex_t
void cppyy_destructor(cppyy_type_t klass, cppyy_object_t self)
static const ClassRefs_t::size_type STD_HANDLE
char cppyy_call_c(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
ptrdiff_t cppyy_base_offset(cppyy_type_t derived, cppyy_type_t base, cppyy_object_t address, int direction)
long double cppyy_call_ld(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
int cppyy_is_const_data(cppyy_scope_t scope, cppyy_index_t idata)
const char ** cppyy_get_all_cpp_names(cppyy_scope_t scope, size_t *count)
void cppyy_call_v(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
static std::set< std::string > gSmartPtrTypes
double cppyy_longdouble2double(void *p)
int cppyy_method_req_args(cppyy_method_t method)
int cppyy_is_staticmethod(cppyy_method_t method)
static GlobalVars_t g_globalvars
cppyy_index_t * cppyy_method_indices_from_name(cppyy_scope_t scope, const char *name)
static char * cppstring_to_cstring(const std::string &cppstr)
char * cppyy_method_signature(cppyy_method_t method, int show_formalargs)
int cppyy_is_abstract(cppyy_type_t type)
static void release_args(Parameter *args, size_t nargs)
int cppyy_call_i(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
static std::set< std::string > g_builtins
unsigned char cppyy_call_b(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
int cppyy_is_template(const char *template_name)
char * cppyy_method_result_type(cppyy_method_t method)
cppyy_type_t cppyy_actual_class(cppyy_type_t klass, cppyy_object_t obj)
int cppyy_is_templated_constructor(cppyy_scope_t scope, cppyy_index_t imeth)
char * cppyy_method_name(cppyy_method_t method)
int cppyy_is_protectedmethod(cppyy_method_t method)
char * cppyy_final_name(cppyy_type_t type)
char * cppyy_get_templated_method_name(cppyy_scope_t scope, cppyy_index_t imeth)
char * cppyy_method_full_name(cppyy_method_t method)
long cppyy_call_l(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
void cppyy_double2longdouble(double d, void *p)
void cppyy_deallocate_function_args(void *args)
int cppyy_method_is_template(cppyy_scope_t scope, cppyy_index_t idx)
static std::vector< CallWrapper * > gWrapperHolder
static bool gEnableFastPath
void cppyy_destruct(cppyy_type_t type, cppyy_object_t self)
int cppyy_is_enum(const char *type_name)
size_t cppyy_size_of_type(const char *type_name)
int cppyy_is_destructor(cppyy_method_t method)
cppyy_object_t cppyy_allocate(cppyy_type_t type)
char * cppyy_datamember_name(cppyy_scope_t scope, int datamember_index)
cppyy_object_t cppyy_construct(cppyy_type_t type)
static ClassRefs_t g_classrefs(1)
static bool WrapperCall(Cppyy::TCppMethod_t method, size_t nargs, void *args_, void *self, void *result)
static std::map< std::string, std::string > resolved_enum_types
int cppyy_is_namespace(cppyy_scope_t scope)
cppyy_object_t cppyy_call_o(cppyy_method_t method, cppyy_object_t self, int nargs, void *args, cppyy_type_t result_type)
long long cppyy_call_ll(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
cppyy_scope_t cppyy_get_scope(const char *scope_name)
void cppyy_vectorbool_setitem(cppyy_object_t ptr, int idx, int value)
static CallWrapper * new_CallWrapper(TFunction *f)
void cppyy_add_smartptr_type(const char *type_name)
void * cppyy_allocate_function_args(int nargs)
float cppyy_call_f(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
short cppyy_call_h(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
int cppyy_is_protecteddata(cppyy_type_t type, cppyy_index_t datamember_index)
static std::set< std::string > gSTLNames
static const ClassRefs_t::size_type GLOBAL_HANDLE
int cppyy_get_dimension_size(cppyy_scope_t scope, cppyy_index_t idata, int dimension)
char * cppyy_method_signature_max(cppyy_method_t method, int show_formalargs, int maxargs)
double cppyy_call_d(cppyy_method_t method, cppyy_object_t self, int nargs, void *args)
int cppyy_num_methods(cppyy_scope_t scope)
char * cppyy_method_prototype(cppyy_scope_t scope, cppyy_method_t method, int show_formalargs)
std::vector< TGlobal * > GlobalVars_t
char * cppyy_resolve_name(const char *cppitem_name)
cppyy_method_t cppyy_get_method(cppyy_scope_t scope, cppyy_index_t idx)
std::vector< TClassRef > ClassRefs_t
size_t cppyy_function_arg_sizeof()
int cppyy_is_const_method(cppyy_method_t method)
static std::map< TDictionary::DeclId_t, CallWrapper * > gMethodTemplates
int cppyy_num_datamembers(cppyy_scope_t scope)
#define CPPYY_HANDLE_EXCEPTION
size_t cppyy_size_of_klass(cppyy_type_t klass)
static bool match_name(const std::string &tname, const std::string fname)
int cppyy_datamember_index(cppyy_scope_t scope, const char *name)
cppyy_index_t cppyy_get_global_operator(cppyy_scope_t scope, cppyy_scope_t lc, cppyy_scope_t rc, const char *op)
static TFunction * m2f(Cppyy::TCppMethod_t method)
void cppyy_deallocate(cppyy_type_t type, cppyy_object_t self)
const char * cppyy_stdstring2charp(cppyy_object_t ptr, size_t *lsz)
int cppyy_compile(const char *code)
static TDataMember * GetDataMemberByIndex(TClassRef cr, int idata)
#define CPPYY_IMP_CALL(typecode, rtype)
int cppyy_method_num_args(cppyy_method_t method)
static std::set< std::string > gInitialNames
char * cppyy_method_arg_default(cppyy_method_t method, int arg_index)
RPY_EXTERN void * cppyy_load_dictionary(const char *lib_name)
static T CallT(Cppyy::TCppMethod_t method, Cppyy::TCppObject_t self, size_t nargs, void *args)
cppyy_funcaddr_t cppyy_function_address(cppyy_method_t method)
int cppyy_has_complex_hierarchy(cppyy_type_t type)
char * cppyy_method_arg_name(cppyy_method_t method, int arg_index)
const Int_t n
Definition legend1.C:16
#define F(x, y, z)
#define I(x, y, z)
#define H(x, y, z)
RPY_EXPORTED short CallH(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
size_t TCppIndex_t
Definition cpp_cppyy.h:24
RPY_EXPORTED std::vector< TCppIndex_t > GetMethodIndicesFromName(TCppScope_t scope, const std::string &name)
RPY_EXPORTED bool Compile(const std::string &code)
RPY_EXPORTED TCppObject_t Allocate(TCppType_t type)
RPY_EXPORTED TCppType_t GetActualClass(TCppType_t klass, TCppObject_t obj)
RPY_EXPORTED TCppIndex_t GetNumEnumData(TCppEnum_t)
RPY_EXPORTED ptrdiff_t GetBaseOffset(TCppType_t derived, TCppType_t base, TCppObject_t address, int direction, bool rerror=false)
intptr_t TCppMethod_t
Definition cpp_cppyy.h:22
RPY_EXPORTED std::string GetMethodArgDefault(TCppMethod_t, TCppIndex_t iarg)
RPY_EXPORTED void AddSmartPtrType(const std::string &)
RPY_EXPORTED bool IsBuiltin(const std::string &type_name)
RPY_EXPORTED bool GetSmartPtrInfo(const std::string &, TCppType_t *raw, TCppMethod_t *deref)
RPY_EXPORTED size_t GetFunctionArgTypeoffset()
RPY_EXPORTED void * AllocateFunctionArgs(size_t nargs)
RPY_EXPORTED bool IsProtectedMethod(TCppMethod_t method)
RPY_EXPORTED std::string GetMethodSignature(TCppMethod_t, bool show_formalargs, TCppIndex_t maxargs=(TCppIndex_t) -1)
RPY_EXPORTED std::string GetMethodPrototype(TCppScope_t scope, TCppMethod_t, bool show_formalargs)
RPY_EXPORTED unsigned char CallB(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetTemplatedMethodName(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED std::vector< TCppScope_t > GetUsingNamespaces(TCppScope_t)
RPY_EXPORTED TCppIndex_t GetGlobalOperator(TCppType_t scope, const std::string &lc, const std::string &rc, const std::string &op)
void * TCppObject_t
Definition cpp_cppyy.h:21
RPY_EXPORTED Long64_t CallLL(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED bool IsPublicData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool HasComplexHierarchy(TCppType_t type)
RPY_EXPORTED TCppObject_t Construct(TCppType_t type)
RPY_EXPORTED std::string GetEnumDataName(TCppEnum_t, TCppIndex_t idata)
RPY_EXPORTED int GetDimensionSize(TCppScope_t scope, TCppIndex_t idata, int dimension)
RPY_EXPORTED float CallF(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED intptr_t GetDatamemberOffset(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool IsTemplate(const std::string &template_name)
TCppScope_t TCppType_t
Definition cpp_cppyy.h:19
RPY_EXPORTED bool IsComplete(const std::string &type_name)
RPY_EXPORTED bool IsProtectedData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool IsEnum(const std::string &type_name)
RPY_EXPORTED long CallL(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED void * CallR(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetMethodMangledName(TCppMethod_t)
RPY_EXPORTED std::string GetDatamemberName(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool IsEnumData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool IsMethodTemplate(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED TCppIndex_t GetNumTemplatedMethods(TCppScope_t scope)
RPY_EXPORTED size_t GetFunctionArgSizeof()
void * TCppEnum_t
Definition cpp_cppyy.h:20
RPY_EXPORTED bool IsConstructor(TCppMethod_t method)
RPY_EXPORTED TCppIndex_t GetDatamemberIndex(TCppScope_t scope, const std::string &name)
RPY_EXPORTED bool IsDestructor(TCppMethod_t method)
RPY_EXPORTED void CallDestructor(TCppType_t type, TCppObject_t self)
RPY_EXPORTED bool IsStaticData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED TCppMethod_t GetMethodTemplate(TCppScope_t scope, const std::string &name, const std::string &proto)
RPY_EXPORTED TCppMethod_t GetMethod(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED std::string GetFinalName(TCppType_t type)
RPY_EXPORTED bool IsPublicMethod(TCppMethod_t method)
RPY_EXPORTED void CallV(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED char * CallS(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args, size_t *length)
RPY_EXPORTED TCppIndex_t GetMethodReqArgs(TCppMethod_t)
RPY_EXPORTED bool IsStaticMethod(TCppMethod_t method)
RPY_EXPORTED std::string GetMethodName(TCppMethod_t)
RPY_EXPORTED long long GetEnumDataValue(TCppEnum_t, TCppIndex_t idata)
RPY_EXPORTED bool IsConstData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED TCppIndex_t GetNumBases(TCppType_t type)
RPY_EXPORTED bool IsTemplatedConstructor(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED TCppObject_t CallO(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args, TCppType_t result_type)
RPY_EXPORTED char CallC(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetMethodArgName(TCppMethod_t, TCppIndex_t iarg)
RPY_EXPORTED double CallD(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED bool IsSubtype(TCppType_t derived, TCppType_t base)
RPY_EXPORTED bool IsConstMethod(TCppMethod_t)
size_t TCppScope_t
Definition cpp_cppyy.h:18
RPY_EXPORTED LongDouble_t CallLD(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetMethodFullName(TCppMethod_t)
RPY_EXPORTED std::string ResolveEnum(const std::string &enum_type)
RPY_EXPORTED std::string ResolveName(const std::string &cppitem_name)
RPY_EXPORTED TCppFuncAddr_t GetFunctionAddress(TCppMethod_t method, bool check_enabled=true)
RPY_EXPORTED TCppObject_t CallConstructor(TCppMethod_t method, TCppType_t type, size_t nargs, void *args)
RPY_EXPORTED size_t SizeOf(TCppType_t klass)
RPY_EXPORTED std::string GetBaseName(TCppType_t type, TCppIndex_t ibase)
RPY_EXPORTED TCppIndex_t GetMethodNumArgs(TCppMethod_t)
RPY_EXPORTED std::string GetScopedFinalName(TCppType_t type)
RPY_EXPORTED bool IsNamespace(TCppScope_t scope)
RPY_EXPORTED void GetAllCppNames(TCppScope_t scope, std::set< std::string > &cppnames)
RPY_EXPORTED TCppIndex_t GetNumDatamembers(TCppScope_t scope)
RPY_EXPORTED bool IsSmartPtr(TCppType_t type)
RPY_EXPORTED TCppScope_t GetScope(const std::string &scope_name)
RPY_EXPORTED std::string GetDatamemberType(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED void Destruct(TCppType_t type, TCppObject_t instance)
RPY_EXPORTED std::string GetMethodResultType(TCppMethod_t)
RPY_EXPORTED bool HasVirtualDestructor(TCppType_t type)
RPY_EXPORTED bool ExistsMethodTemplate(TCppScope_t scope, const std::string &name)
RPY_EXPORTED TCppScope_t gGlobalScope
Definition cpp_cppyy.h:51
void * TCppFuncAddr_t
Definition cpp_cppyy.h:25
RPY_EXPORTED std::string GetMethodArgType(TCppMethod_t, TCppIndex_t iarg)
RPY_EXPORTED void DeallocateFunctionArgs(void *args)
RPY_EXPORTED int CallI(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED TCppIndex_t GetNumMethods(TCppScope_t scope)
RPY_EXPORTED void Deallocate(TCppType_t type, TCppObject_t instance)
RPY_EXPORTED TCppEnum_t GetEnum(TCppScope_t scope, const std::string &enum_name)
RPY_EXPORTED bool IsAbstract(TCppType_t type)
void(* DelFunc_t)(void *)
Definition Rtypes.h:110
std::string ResolveTypedef(const char *tname, bool resolveAll=false)
std::string CleanType(const char *typeDesc, int mode=0, const char **tail=0)
Cleanup type description, redundant blanks removed and redundant tail ignored return *tail = pointer ...
std::string ShortType(const char *typeDesc, int mode)
Return the absolute type of typeDesc.
#define RPY_EXTERN
union CPyCppyy::Parameter::Value fValue
auto * m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345