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 void HandleException(Int_t sig) override {
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 // create helpers for comparing thingies
276 gInterpreter->Declare(
277 "namespace __cppyy_internal { template<class C1, class C2>"
278 " bool is_equal(const C1& c1, const C2& c2) { return (bool)(c1 == c2); } }");
279 gInterpreter->Declare(
280 "namespace __cppyy_internal { template<class C1, class C2>"
281 " bool is_not_equal(const C1& c1, const C2& c2) { return (bool)(c1 != c2); } }");
282
283 // retrieve all initial (ROOT) C++ names in the global scope to allow filtering later
284 if (!getenv("CPPYY_NO_ROOT_FILTER")) {
285 gROOT->GetListOfGlobals(true); // force initialize
286 gROOT->GetListOfGlobalFunctions(true); // id.
287 std::set<std::string> initial;
289 gInitialNames = initial;
290
291#ifndef WIN32
292 gRootSOs.insert("libCore.so ");
293 gRootSOs.insert("libRIO.so ");
294 gRootSOs.insert("libThread.so ");
295 gRootSOs.insert("libMathCore.so ");
296#else
297 gRootSOs.insert("libCore.dll ");
298 gRootSOs.insert("libRIO.dll ");
299 gRootSOs.insert("libThread.dll ");
300 gRootSOs.insert("libMathCore.dll ");
301#endif
302 }
303
304 // start off with a reasonable size placeholder for wrappers
305 gWrapperHolder.reserve(1024);
306
307 // create an exception handler to process signals
308 gExceptionHandler = new TExceptionHandlerImp{};
309 }
310
311 ~ApplicationStarter() {
312 for (auto wrap : gWrapperHolder)
313 delete wrap;
314 delete gExceptionHandler; gExceptionHandler = nullptr;
315 }
316} _applicationStarter;
317
318} // unnamed namespace
319
320
321// local helpers -------------------------------------------------------------
322static inline
324{
325 assert((ClassRefs_t::size_type)scope < g_classrefs.size());
326 return g_classrefs[(ClassRefs_t::size_type)scope];
327}
328
329static inline
331 CallWrapper* wrap = ((CallWrapper*)method);
332 if (!wrap->fTF || wrap->fTF->GetDeclId() != wrap->fDecl) {
333 MethodInfo_t* mi = gInterpreter->MethodInfo_Factory(wrap->fDecl);
334 wrap->fTF = new TFunction(mi);
335 }
336 return wrap->fTF;
337}
338
339static inline
340char* cppstring_to_cstring(const std::string& cppstr)
341{
342 char* cstr = (char*)malloc(cppstr.size()+1);
343 memcpy(cstr, cppstr.c_str(), cppstr.size()+1);
344 return cstr;
345}
346
347static inline
348bool match_name(const std::string& tname, const std::string fname)
349{
350// either match exactly, or match the name as template
351 if (fname.rfind(tname, 0) == 0) {
352 if ((tname.size() == fname.size()) ||
353 (tname.size() < fname.size() && fname[tname.size()] == '<'))
354 return true;
355 }
356 return false;
357}
358
359static inline
360bool is_missclassified_stl(const std::string& name)
361{
362 std::string::size_type pos = name.find('<');
363 if (pos != std::string::npos)
364 return gSTLNames.find(name.substr(0, pos)) != gSTLNames.end();
365 return gSTLNames.find(name) != gSTLNames.end();
366}
367
368
369// direct interpreter access -------------------------------------------------
370bool Cppyy::Compile(const std::string& code)
371{
372 return gInterpreter->Declare(code.c_str());
373}
374
375
376// name to opaque C++ scope representation -----------------------------------
377std::string Cppyy::ResolveName(const std::string& cppitem_name)
378{
379// Fully resolve the given name to the final type name.
380
381// try memoized type cache, in case seen before
382 TCppType_t klass = find_memoized(cppitem_name);
383 if (klass) return GetScopedFinalName(klass);
384
385// remove global scope '::' if present
386 std::string tclean = cppitem_name.compare(0, 2, "::") == 0 ?
387 cppitem_name.substr(2, std::string::npos) : cppitem_name;
388
389// classes (most common)
390 tclean = TClassEdit::CleanType(tclean.c_str());
391 if (tclean.empty() /* unknown, eg. an operator */) return cppitem_name;
392
393// reduce [N] to []
394 if (tclean[tclean.size()-1] == ']')
395 tclean = tclean.substr(0, tclean.rfind('[')) + "[]";
396
397 if (tclean.rfind("byte", 0) == 0 || tclean.rfind("std::byte", 0) == 0)
398 return tclean;
399
400// check data types list (accept only builtins as typedefs will
401// otherwise not be resolved)
402 TDataType* dt = gROOT->GetType(tclean.c_str());
403 if (dt && dt->GetType() != kOther_t) return dt->GetFullTypeName();
404
405// special case for enums
406 if (IsEnum(cppitem_name))
407 return ResolveEnum(cppitem_name);
408
409// special case for clang's builtin __type_pack_element (which does not resolve)
410 if (cppitem_name.rfind("__type_pack_element", 0) != std::string::npos) {
411 // shape is "__type_pack_element<index,type1,type2,...,typeN>cpd": extract
412 // first the index, and from there the indexed type; finally, restore the
413 // qualifiers
414 const char* str = cppitem_name.c_str();
415 char* endptr = nullptr;
416 unsigned long index = strtoul(str+20, &endptr, 0);
417
418 std::string tmplvars{endptr};
419 auto start = tmplvars.find(',') + 1;
420 auto end = tmplvars.find(',', start);
421 while (index != 0) {
422 start = end+1;
423 end = tmplvars.find(',', start);
424 if (end == std::string::npos) end = tmplvars.rfind('>');
425 --index;
426 }
427
428 std::string resolved = tmplvars.substr(start, end-start);
429 auto cpd = tmplvars.rfind('>');
430 if (cpd != std::string::npos && cpd+1 != tmplvars.size())
431 return resolved + tmplvars.substr(cpd+1, std::string::npos);
432 return resolved;
433 }
434
435// typedefs
436 return TClassEdit::ResolveTypedef(tclean.c_str(), true);
437}
438
439static std::map<std::string, std::string> resolved_enum_types;
440std::string Cppyy::ResolveEnum(const std::string& enum_type)
441{
442// The underlying type of a an enum may be any kind of integer.
443// Resolve that type via a workaround (note: this function assumes
444// that the enum_type name is a valid enum type name)
445 auto res = resolved_enum_types.find(enum_type);
446 if (res != resolved_enum_types.end())
447 return res->second;
448
449// desugar the type before resolving
450 std::string et_short = TClassEdit::ShortType(enum_type.c_str(), 1);
451 if (et_short.find("(unnamed") == std::string::npos) {
452 std::ostringstream decl;
453 // TODO: now presumed fixed with https://sft.its.cern.ch/jira/browse/ROOT-6988
454 for (auto& itype : {"unsigned int"}) {
455 decl << "std::is_same<"
456 << itype
457 << ", std::underlying_type<"
458 << et_short
459 << ">::type>::value;";
460 if (gInterpreter->ProcessLine(decl.str().c_str())) {
461 // TODO: "re-sugaring" like this is brittle, but the top
462 // should be re-translated into AST-based code anyway
463 std::string resugared;
464 if (et_short.size() != enum_type.size()) {
465 auto pos = enum_type.find(et_short);
466 if (pos != std::string::npos) {
467 resugared = enum_type.substr(0, pos) + itype;
468 if (pos+et_short.size() < enum_type.size())
469 resugared += enum_type.substr(pos+et_short.size(), std::string::npos);
470 }
471 }
472 if (resugared.empty()) resugared = itype;
473 resolved_enum_types[enum_type] = resugared;
474 return resugared;
475 }
476 }
477 }
478
479// failed or anonymous ... signal upstream to special case this
480 int ipos = (int)enum_type.size()-1;
481 for (; 0 <= ipos; --ipos) {
482 char c = enum_type[ipos];
483 if (isspace(c)) continue;
484 if (isalnum(c) || c == '_' || c == '>' || c == ')') break;
485 }
486 bool isConst = enum_type.find("const ", 6) != std::string::npos;
487 std::string restype = isConst ? "const " : "";
488 restype += "internal_enum_type_t"+enum_type.substr((std::string::size_type)ipos+1, std::string::npos);
489 resolved_enum_types[enum_type] = restype;
490 return restype; // should default to some int variant
491}
492
493Cppyy::TCppScope_t Cppyy::GetScope(const std::string& sname)
494{
495// First, try cache
496 TCppType_t result = find_memoized(sname);
497 if (result) return result;
498
499// Second, skip builtins before going through the more expensive steps of resolving
500// typedefs and looking up TClass
501 if (g_builtins.find(sname) != g_builtins.end())
502 return (TCppScope_t)0;
503
504// TODO: scope_name should always be final already?
505// Resolve name fully before lookup to make sure all aliases point to the same scope
506 std::string scope_name = ResolveName(sname);
507 bool bHasAlias = sname != scope_name;
508 if (bHasAlias) {
509 result = find_memoized(scope_name);
510 if (result) return result;
511 }
512
513// both failed, but may be STL name that's missing 'std::' now, but didn't before
514 bool b_scope_name_missclassified = is_missclassified_stl(scope_name);
515 if (b_scope_name_missclassified) {
516 result = find_memoized("std::"+scope_name);
517 if (result) g_name2classrefidx["std::"+scope_name] = (ClassRefs_t::size_type)result;
518 }
519 bool b_sname_missclassified = bHasAlias ? is_missclassified_stl(sname) : false;
520 if (b_sname_missclassified) {
521 if (!result) result = find_memoized("std::"+sname);
522 if (result) g_name2classrefidx["std::"+sname] = (ClassRefs_t::size_type)result;
523 }
524
525 if (result) return result;
526
527// use TClass directly, to enable auto-loading; class may be stubbed (eg. for
528// function returns) or forward declared, leading to a non-null TClass that is
529// otherwise invalid/unusable
530 TClassRef cr(TClass::GetClass(scope_name.c_str(), true /* load */, true /* silent */));
531 if (!cr.GetClass())
532 return (TCppScope_t)0;
533
534// memoize found/created TClass
535 ClassRefs_t::size_type sz = g_classrefs.size();
536 g_name2classrefidx[scope_name] = sz;
537 if (bHasAlias) g_name2classrefidx[sname] = sz;
538 g_classrefs.push_back(TClassRef(scope_name.c_str()));
539
540// TODO: make ROOT/meta NOT remove std :/
541 if (b_scope_name_missclassified)
542 g_name2classrefidx["std::"+scope_name] = sz;
543 if (b_sname_missclassified)
544 g_name2classrefidx["std::"+sname] = sz;
545
546 return (TCppScope_t)sz;
547}
548
549bool Cppyy::IsTemplate(const std::string& template_name)
550{
551 return (bool)gInterpreter->CheckClassTemplate(template_name.c_str());
552}
553
554namespace {
555 class AutoCastRTTI {
556 public:
557 virtual ~AutoCastRTTI() {}
558 };
559}
560
562{
563 TClassRef& cr = type_from_handle(klass);
564 if (!cr.GetClass() || !obj) return klass;
565
566#ifdef _WIN64
567// Cling does not provide a consistent ImageBase address for calculating relative addresses
568// as used in Windows 64b RTTI. So, check for our own RTTI extension instead. If that fails,
569// see whether the unmangled raw_name is available (e.g. if this is an MSVC compiled rather
570// than JITed class) and pass on if it is.
571 volatile const char* raw = nullptr; // to prevent too aggressive reordering
572 try {
573 // this will filter those objects that do not have RTTI to begin with (throws)
574 AutoCastRTTI* pcst = (AutoCastRTTI*)obj;
575 raw = typeid(*pcst).raw_name();
576
577 // check the signature id (0 == absolute, 1 == relative, 2 == ours)
578 void* vfptr = *(void**)((intptr_t)obj);
579 void* meta = (void*)((intptr_t)*((void**)((intptr_t)vfptr-sizeof(void*))));
580 if (*(intptr_t*)meta == 2) {
581 // access the extra data item which is an absolute pointer to the RTTI
582 void* ptdescr = (void*)((intptr_t)meta + 4*sizeof(unsigned long)+sizeof(void*));
583 if (ptdescr && *(void**)ptdescr) {
584 auto rtti = *(std::type_info**)ptdescr;
585 raw = rtti->raw_name();
586 if (raw && raw[0] != '\0') // likely unnecessary
587 return (TCppType_t)GetScope(rtti->name());
588 }
589
590 return klass; // do not fall through if no RTTI info available
591 }
592
593 // if the raw name is the empty string (no guarantees that this is so as truly, the
594 // address is corrupt, but it is common to be empty), then there is no accessible RTTI
595 // and getting the unmangled name will crash ...
596 if (!raw)
597 return klass;
598 } catch (std::bad_typeid) {
599 return klass; // can't risk passing to ROOT/meta as it may do RTTI
600 }
601#endif
602
603 TClass* clActual = cr->GetActualClass((void*)obj);
604 if (clActual && clActual != cr.GetClass()) {
605 auto itt = g_name2classrefidx.find(clActual->GetName());
606 if (itt != g_name2classrefidx.end())
607 return (TCppType_t)itt->second;
608 return (TCppType_t)GetScope(clActual->GetName());
609 }
610
611 return klass;
612}
613
615{
616 TClassRef& cr = type_from_handle(klass);
617 if (cr.GetClass() && cr->GetClassInfo())
618 return (size_t)gInterpreter->ClassInfo_Size(cr->GetClassInfo());
619 return (size_t)0;
620}
621
622size_t Cppyy::SizeOf(const std::string& type_name)
623{
624 TDataType* dt = gROOT->GetType(type_name.c_str());
625 if (dt) return dt->Size();
626 return SizeOf(GetScope(type_name));
627}
628
629bool Cppyy::IsBuiltin(const std::string& type_name)
630{
631 TDataType* dt = gROOT->GetType(TClassEdit::CleanType(type_name.c_str(), 1).c_str());
632 if (dt) return dt->GetType() != kOther_t;
633 return false;
634}
635
636bool Cppyy::IsComplete(const std::string& type_name)
637{
638// verify whether the dictionary of this class is fully available
639 bool b = false;
640
641 int oldEIL = gErrorIgnoreLevel;
642 gErrorIgnoreLevel = 3000;
643 TClass* klass = TClass::GetClass(TClassEdit::ShortType(type_name.c_str(), 1).c_str());
644 if (klass && klass->GetClassInfo()) // works for normal case w/ dict
645 b = gInterpreter->ClassInfo_IsLoaded(klass->GetClassInfo());
646 else { // special case for forward declared classes
647 ClassInfo_t* ci = gInterpreter->ClassInfo_Factory(type_name.c_str());
648 if (ci) {
649 b = gInterpreter->ClassInfo_IsLoaded(ci);
650 gInterpreter->ClassInfo_Delete(ci); // we own the fresh class info
651 }
652 }
653 gErrorIgnoreLevel = oldEIL;
654 return b;
655}
656
657// memory management ---------------------------------------------------------
659{
661 return (TCppObject_t)malloc(gInterpreter->ClassInfo_Size(cr->GetClassInfo()));
662}
663
664void Cppyy::Deallocate(TCppType_t /* type */, TCppObject_t instance)
665{
666 ::operator delete(instance);
667}
668
670{
672 return (TCppObject_t)cr->New();
673}
674
675static std::map<Cppyy::TCppType_t, bool> sHasOperatorDelete;
677{
680 cr->Destructor((void*)instance);
681 else {
682 ROOT::DelFunc_t fdel = cr->GetDelete();
683 if (fdel) fdel((void*)instance);
684 else {
685 auto ib = sHasOperatorDelete.find(type);
686 if (ib == sHasOperatorDelete.end()) {
688 ib = sHasOperatorDelete.find(type);
689 }
690 ib->second ? cr->Destructor((void*)instance) : free((void*)instance);
691 }
692 }
693}
694
695
696// method/function dispatching -----------------------------------------------
698{
699// TODO: method should be a callfunc, so that no mapping would be needed.
700 CallWrapper* wrap = (CallWrapper*)method;
701
702 CallFunc_t* callf = gInterpreter->CallFunc_Factory();
703 MethodInfo_t* meth = gInterpreter->MethodInfo_Factory(wrap->fDecl);
704 gInterpreter->CallFunc_SetFunc(callf, meth);
705 gInterpreter->MethodInfo_Delete(meth);
706
707 if (!(callf && gInterpreter->CallFunc_IsValid(callf))) {
708 // TODO: propagate this error to caller w/o use of Python C-API
709 /*
710 PyErr_Format(PyExc_RuntimeError, "could not resolve %s::%s(%s)",
711 const_cast<TClassRef&>(klass).GetClassName(),
712 wrap.fName, callString.c_str()); */
713 std::cerr << "TODO: report unresolved function error to Python\n";
714 if (callf) gInterpreter->CallFunc_Delete(callf);
716 }
717
718// generate the wrapper and JIT it; ignore wrapper generation errors (will simply
719// result in a nullptr that is reported upstream if necessary; often, however,
720// there is a different overload available that will do)
721 auto oldErrLvl = gErrorIgnoreLevel;
723 wrap->fFaceptr = gInterpreter->CallFunc_IFacePtr(callf);
724 gErrorIgnoreLevel = oldErrLvl;
725
726 gInterpreter->CallFunc_Delete(callf); // does not touch IFacePtr
727 return wrap->fFaceptr;
728}
729
730static inline
731bool copy_args(Parameter* args, size_t nargs, void** vargs)
732{
733 bool runRelease = false;
734 for (size_t i = 0; i < nargs; ++i) {
735 switch (args[i].fTypeCode) {
736 case 'X': /* (void*)type& with free */
737 runRelease = true;
738 case 'V': /* (void*)type& */
739 vargs[i] = args[i].fValue.fVoidp;
740 break;
741 case 'r': /* const type& */
742 vargs[i] = args[i].fRef;
743 break;
744 default: /* all other types in union */
745 vargs[i] = (void*)&args[i].fValue.fVoidp;
746 break;
747 }
748 }
749 return runRelease;
750}
751
752static inline
753void release_args(Parameter* args, size_t nargs) {
754 for (size_t i = 0; i < nargs; ++i) {
755 if (args[i].fTypeCode == 'X')
756 free(args[i].fValue.fVoidp);
757 }
758}
759
760static inline bool WrapperCall(Cppyy::TCppMethod_t method, size_t nargs, void* args_, void* self, void* result)
761{
762 Parameter* args = (Parameter*)args_;
763
764 CallWrapper* wrap = (CallWrapper*)method;
765 const TInterpreter::CallFuncIFacePtr_t& faceptr = wrap->fFaceptr.fGeneric ? wrap->fFaceptr : GetCallFunc(method);
766 if (!faceptr.fGeneric)
767 return false; // happens with compilation error
768
770 bool runRelease = false;
771 if (nargs <= SMALL_ARGS_N) {
772 void* smallbuf[SMALL_ARGS_N];
773 if (nargs) runRelease = copy_args(args, nargs, smallbuf);
774 faceptr.fGeneric(self, (int)nargs, smallbuf, result);
775 } else {
776 std::vector<void*> buf(nargs);
777 runRelease = copy_args(args, nargs, buf.data());
778 faceptr.fGeneric(self, (int)nargs, buf.data(), result);
779 }
780 if (runRelease) release_args(args, nargs);
781 return true;
782 }
783
785 bool runRelease = false;
786 if (nargs <= SMALL_ARGS_N) {
787 void* smallbuf[SMALL_ARGS_N];
788 if (nargs) runRelease = copy_args(args, nargs, (void**)smallbuf);
789 faceptr.fCtor((void**)smallbuf, result, (unsigned long)nargs);
790 } else {
791 std::vector<void*> buf(nargs);
792 runRelease = copy_args(args, nargs, buf.data());
793 faceptr.fCtor(buf.data(), result, (unsigned long)nargs);
794 }
795 if (runRelease) release_args(args, nargs);
796 return true;
797 }
798
800 std::cerr << " DESTRUCTOR NOT IMPLEMENTED YET! " << std::endl;
801 return false;
802 }
803
804 return false;
805}
806
807template<typename T>
808static inline
809T CallT(Cppyy::TCppMethod_t method, Cppyy::TCppObject_t self, size_t nargs, void* args)
810{
811 T t{};
812 if (WrapperCall(method, nargs, args, (void*)self, &t))
813 return t;
814 return (T)-1;
815}
816
817#define CPPYY_IMP_CALL(typecode, rtype) \
818rtype Cppyy::Call##typecode(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)\
819{ \
820 return CallT<rtype>(method, self, nargs, args); \
821}
822
823void Cppyy::CallV(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)
824{
825 if (!WrapperCall(method, nargs, args, (void*)self, nullptr))
826 return /* TODO ... report error */;
827}
828
829CPPYY_IMP_CALL(B, unsigned char)
836CPPYY_IMP_CALL(D, double )
838
839void* Cppyy::CallR(TCppMethod_t method, TCppObject_t self, size_t nargs, void* args)
840{
841 void* r = nullptr;
842 if (WrapperCall(method, nargs, args, (void*)self, &r))
843 return r;
844 return nullptr;
845}
846
848 TCppMethod_t method, TCppObject_t self, size_t nargs, void* args, size_t* length)
849{
850 char* cstr = nullptr;
851 TClassRef cr("std::string");
852 std::string* cppresult = (std::string*)malloc(sizeof(std::string));
853 if (WrapperCall(method, nargs, args, self, (void*)cppresult)) {
854 cstr = cppstring_to_cstring(*cppresult);
855 *length = cppresult->size();
856 cppresult->std::string::~basic_string();
857 } else
858 *length = 0;
859 free((void*)cppresult);
860 return cstr;
861}
862
864 TCppMethod_t method, TCppType_t /* klass */, size_t nargs, void* args)
865{
866 void* obj = nullptr;
867 if (WrapperCall(method, nargs, args, nullptr, &obj))
868 return (TCppObject_t)obj;
869 return (TCppObject_t)0;
870}
871
873{
875 cr->Destructor((void*)self, true);
876}
877
879 TCppObject_t self, size_t nargs, void* args, TCppType_t result_type)
880{
881 TClassRef& cr = type_from_handle(result_type);
882 void* obj = ::operator new(gInterpreter->ClassInfo_Size(cr->GetClassInfo()));
883 if (WrapperCall(method, nargs, args, self, obj))
884 return (TCppObject_t)obj;
885 ::operator delete(obj);
886 return (TCppObject_t)0;
887}
888
890{
891 if (check_enabled && !gEnableFastPath) return (TCppFuncAddr_t)nullptr;
892 TFunction* f = m2f(method);
893 return (TCppFuncAddr_t)gInterpreter->FindSym(f->GetMangledName());
894}
895
896
897// handling of function argument buffer --------------------------------------
899{
900 return new Parameter[nargs];
901}
902
904{
905 delete [] (Parameter*)args;
906}
907
909{
910 return sizeof(Parameter);
911}
912
914{
915 return offsetof(Parameter, fTypeCode);
916}
917
918
919// scope reflection information ----------------------------------------------
921{
922// Test if this scope represents a namespace.
923 if (scope == GLOBAL_HANDLE)
924 return true;
925 TClassRef& cr = type_from_handle(scope);
926 if (cr.GetClass())
927 return cr->Property() & kIsNamespace;
928 return false;
929}
930
932{
933// Test if this type may not be instantiated.
934 TClassRef& cr = type_from_handle(klass);
935 if (cr.GetClass())
936 return cr->Property() & kIsAbstract;
937 return false;
938}
939
940bool Cppyy::IsEnum(const std::string& type_name)
941{
942 if (type_name.empty()) return false;
943 std::string tn_short = TClassEdit::ShortType(type_name.c_str(), 1);
944 if (tn_short.empty()) return false;
945 return gInterpreter->ClassInfo_IsEnum(tn_short.c_str());
946}
947
948// helpers for stripping scope names
949static
950std::string outer_with_template(const std::string& name)
951{
952// Cut down to the outer-most scope from <name>, taking proper care of templates.
953 int tpl_open = 0;
954 for (std::string::size_type pos = 0; pos < name.size(); ++pos) {
955 std::string::value_type c = name[pos];
956
957 // count '<' and '>' to be able to skip template contents
958 if (c == '<')
959 ++tpl_open;
960 else if (c == '>')
961 --tpl_open;
962
963 // collect name up to "::"
964 else if (tpl_open == 0 && \
965 c == ':' && pos+1 < name.size() && name[pos+1] == ':') {
966 // found the extend of the scope ... done
967 return name.substr(0, pos-1);
968 }
969 }
970
971// whole name is apparently a single scope
972 return name;
973}
974
975static
976std::string outer_no_template(const std::string& name)
977{
978// Cut down to the outer-most scope from <name>, drop templates
979 std::string::size_type first_scope = name.find(':');
980 if (first_scope == std::string::npos)
981 return name.substr(0, name.find('<'));
982 std::string::size_type first_templ = name.find('<');
983 if (first_templ == std::string::npos)
984 return name.substr(0, first_scope);
985 return name.substr(0, std::min(first_templ, first_scope));
986}
987
988#define FILL_COLL(type, filter) { \
989 TIter itr{coll}; \
990 type* obj = nullptr; \
991 while ((obj = (type*)itr.Next())) { \
992 const char* nm = obj->GetName(); \
993 if (nm && nm[0] != '_' && !(obj->Property() & (filter))) { \
994 if (gInitialNames.find(nm) == gInitialNames.end()) \
995 cppnames.insert(nm); \
996 }}}
997
998static inline
999void cond_add(Cppyy::TCppScope_t scope, const std::string& ns_scope,
1000 std::set<std::string>& cppnames, const char* name, bool nofilter = false)
1001{
1002 if (!name || name[0] == '_' || strstr(name, ".h") != 0 || strncmp(name, "operator", 8) == 0)
1003 return;
1004
1005 if (scope == GLOBAL_HANDLE) {
1006 std::string to_add = outer_no_template(name);
1007 if ((nofilter || gInitialNames.find(to_add) == gInitialNames.end()) && !is_missclassified_stl(name))
1008 cppnames.insert(outer_no_template(name));
1009 } else if (scope == STD_HANDLE) {
1010 if (strncmp(name, "std::", 5) == 0) {
1011 name += 5;
1012#ifdef __APPLE__
1013 if (strncmp(name, "__1::", 5) == 0) name += 5;
1014#endif
1015 } else if (!is_missclassified_stl(name))
1016 return;
1017 cppnames.insert(outer_no_template(name));
1018 } else {
1019 if (strncmp(name, ns_scope.c_str(), ns_scope.size()) == 0)
1020 cppnames.insert(outer_with_template(name + ns_scope.size()));
1021 }
1022}
1023
1024void Cppyy::GetAllCppNames(TCppScope_t scope, std::set<std::string>& cppnames)
1025{
1026// Collect all known names of C++ entities under scope. This is useful for IDEs
1027// employing tab-completion, for example. Note that functions names need not be
1028// unique as they can be overloaded.
1029 TClassRef& cr = type_from_handle(scope);
1030 if (scope != GLOBAL_HANDLE && !(cr.GetClass() && cr->Property()))
1031 return;
1032
1033 std::string ns_scope = GetFinalName(scope);
1034 if (scope != GLOBAL_HANDLE) ns_scope += "::";
1035
1036// add existing values from read rootmap files if within this scope
1037 TCollection* coll = gInterpreter->GetMapfile()->GetTable();
1038 {
1039 TIter itr{coll};
1040 TEnvRec* ev = nullptr;
1041 while ((ev = (TEnvRec*)itr.Next())) {
1042 // TEnv contains rootmap entries and user-side rootmap files may be already
1043 // loaded on startup. Thus, filter on file name rather than load time.
1044 if (gRootSOs.find(ev->GetValue()) == gRootSOs.end())
1045 cond_add(scope, ns_scope, cppnames, ev->GetName(), true);
1046 }
1047 }
1048
1049// do we care about the class table or are the rootmap and list of types enough?
1050/*
1051 gClassTable->Init();
1052 const int N = gClassTable->Classes();
1053 for (int i = 0; i < N; ++i)
1054 cond_add(scope, ns_scope, cppnames, gClassTable->Next());
1055*/
1056
1057// any other types (e.g. that may have come from parsing headers)
1058 coll = gROOT->GetListOfTypes();
1059 {
1060 TIter itr{coll};
1061 TDataType* dt = nullptr;
1062 while ((dt = (TDataType*)itr.Next())) {
1063 if (!(dt->Property() & kIsFundamental)) {
1064 cond_add(scope, ns_scope, cppnames, dt->GetName());
1065 }
1066 }
1067 }
1068
1069// add functions
1070 coll = (scope == GLOBAL_HANDLE) ?
1071 gROOT->GetListOfGlobalFunctions() : cr->GetListOfMethods();
1072 {
1073 TIter itr{coll};
1074 TFunction* obj = nullptr;
1075 while ((obj = (TFunction*)itr.Next())) {
1076 const char* nm = obj->GetName();
1077 // skip templated functions, adding only the un-instantiated ones
1078 if (nm && nm[0] != '_' && strstr(nm, "<") == 0 && strncmp(nm, "operator", 8) != 0) {
1079 if (gInitialNames.find(nm) == gInitialNames.end())
1080 cppnames.insert(nm);
1081 }
1082 }
1083 }
1084
1085// add uninstantiated templates
1086 coll = (scope == GLOBAL_HANDLE) ?
1087 gROOT->GetListOfFunctionTemplates() : cr->GetListOfFunctionTemplates();
1089
1090// add (global) data members
1091 if (scope == GLOBAL_HANDLE) {
1092 coll = gROOT->GetListOfGlobals();
1094 } else {
1095 coll = cr->GetListOfDataMembers();
1097 coll = cr->GetListOfUsingDataMembers();
1099 }
1100
1101// add enums values only for user classes/namespaces
1102 if (scope != GLOBAL_HANDLE && scope != STD_HANDLE) {
1103 coll = cr->GetListOfEnums();
1105 }
1106
1107#ifdef __APPLE__
1108// special case for Apple, add version namespace '__1' entries to std
1109 if (scope == STD_HANDLE)
1110 GetAllCppNames(GetScope("std::__1"), cppnames);
1111#endif
1112}
1113
1114
1115// class reflection information ----------------------------------------------
1116std::vector<Cppyy::TCppScope_t> Cppyy::GetUsingNamespaces(TCppScope_t scope)
1117{
1118 std::vector<Cppyy::TCppScope_t> res;
1119 if (!IsNamespace(scope))
1120 return res;
1121
1122#ifdef __APPLE__
1123 if (scope == STD_HANDLE) {
1124 res.push_back(GetScope("__1"));
1125 return res;
1126 }
1127#endif
1128
1129 TClassRef& cr = type_from_handle(scope);
1130 if (!cr.GetClass() || !cr->GetClassInfo())
1131 return res;
1132
1133 const std::vector<std::string>& v = gInterpreter->GetUsingNamespaces(cr->GetClassInfo());
1134 res.reserve(v.size());
1135 for (const auto& uid : v) {
1136 Cppyy::TCppScope_t uscope = GetScope(uid);
1137 if (uscope) res.push_back(uscope);
1138 }
1139
1140 return res;
1141}
1142
1143
1144// class reflection information ----------------------------------------------
1146{
1147 if (klass == GLOBAL_HANDLE)
1148 return "";
1149 TClassRef& cr = type_from_handle(klass);
1150 std::string clName = cr->GetName();
1151// TODO: why is this template splitting needed?
1152 std::string::size_type pos = clName.substr(0, clName.find('<')).rfind("::");
1153 if (pos != std::string::npos)
1154 return clName.substr(pos+2, std::string::npos);
1155 return clName;
1156}
1157
1159{
1160 if (klass == GLOBAL_HANDLE)
1161 return "";
1162 TClassRef& cr = type_from_handle(klass);
1163 if (cr.GetClass()) {
1164 std::string name = cr->GetName();
1166 return std::string("std::")+cr->GetName();
1167 return cr->GetName();
1168 }
1169 return "";
1170}
1171
1173{
1174 TClassRef& cr = type_from_handle(klass);
1175 if (!cr.GetClass())
1176 return false;
1177
1178 TFunction* f = cr->GetMethod(("~"+GetFinalName(klass)).c_str(), "");
1179 if (f && (f->Property() & kIsVirtual))
1180 return true;
1181
1182 return false;
1183}
1184
1186{
1187 int is_complex = 1;
1188 size_t nbases = 0;
1189
1190 TClassRef& cr = type_from_handle(klass);
1191 if (cr.GetClass() && cr->GetListOfBases() != 0)
1192 nbases = GetNumBases(klass);
1193
1194 if (1 < nbases)
1195 is_complex = 1;
1196 else if (nbases == 0)
1197 is_complex = 0;
1198 else { // one base class only
1199 TBaseClass* base = (TBaseClass*)cr->GetListOfBases()->At(0);
1200 if (base->Property() & kIsVirtualBase)
1201 is_complex = 1; // TODO: verify; can be complex, need not be.
1202 else
1203 is_complex = HasComplexHierarchy(GetScope(base->GetName()));
1204 }
1205
1206 return is_complex;
1207}
1208
1210{
1211// Get the total number of base classes that this class has.
1212 TClassRef& cr = type_from_handle(klass);
1213 if (cr.GetClass() && cr->GetListOfBases() != 0)
1214 return (TCppIndex_t)cr->GetListOfBases()->GetSize();
1215 return (TCppIndex_t)0;
1216}
1217
1218////////////////////////////////////////////////////////////////////////////////
1219/// \fn Cppyy::TCppIndex_t GetLongestInheritancePath(TClass *klass)
1220/// \brief Retrieve number of base classes in the longest branch of the
1221/// inheritance tree of the input class.
1222/// \param[in] klass The class to start the retrieval process from.
1223///
1224/// This is a helper function for Cppyy::GetNumBasesLongestBranch.
1225/// Given an inheritance tree, the function assigns weight 1 to each class that
1226/// has at least one base. Starting from the input class, the function is
1227/// called recursively on all the bases. For each base the return value is one
1228/// (the weight of the base itself) plus the maximum value retrieved for their
1229/// bases in turn. For example, given the following inheritance tree:
1230///
1231/// ~~~{.cpp}
1232/// class A {}; class B: public A {};
1233/// class X {}; class Y: public X {}; class Z: public Y {};
1234/// class C: public B, Z {};
1235/// ~~~
1236///
1237/// calling this function on an instance of `C` will return 3, the steps
1238/// required to go from C to X.
1240{
1241
1242 auto directbases = klass->GetListOfBases();
1243 if (!directbases) {
1244 // This is a leaf with no bases
1245 return 0;
1246 }
1247 auto ndirectbases = directbases->GetSize();
1248 if (ndirectbases == 0) {
1249 // This is a leaf with no bases
1250 return 0;
1251 } else {
1252 // If there is at least one direct base
1253 std::vector<Cppyy::TCppIndex_t> nbases_branches;
1254 nbases_branches.reserve(ndirectbases);
1255
1256 // Traverse all direct bases of the current class and call the function
1257 // recursively
1258 for (auto baseclass : TRangeDynCast<TBaseClass>(directbases)) {
1259 if (!baseclass)
1260 continue;
1261 if (auto baseclass_tclass = baseclass->GetClassPointer()) {
1262 nbases_branches.emplace_back(GetLongestInheritancePath(baseclass_tclass));
1263 }
1264 }
1265
1266 // Get longest path among the direct bases of the current class
1267 auto longestbranch = std::max_element(std::begin(nbases_branches), std::end(nbases_branches));
1268
1269 // Add 1 to include the current class in the count
1270 return 1 + *longestbranch;
1271 }
1272}
1273
1274////////////////////////////////////////////////////////////////////////////////
1275/// \fn Cppyy::TCppIndex_t Cppyy::GetNumBasesLongestBranch(TCppType_t klass)
1276/// \brief Retrieve number of base classes in the longest branch of the
1277/// inheritance tree.
1278/// \param[in] klass The class to start the retrieval process from.
1279///
1280/// The function converts the input class to a `TClass *` and calls
1281/// GetLongestInheritancePath.
1283{
1284
1285 const auto &cr = type_from_handle(klass);
1286
1287 if (auto klass_tclass = cr.GetClass()) {
1288 return GetLongestInheritancePath(klass_tclass);
1289 }
1290
1291 // In any other case, return zero
1292 return 0;
1293}
1294
1296{
1297 TClassRef& cr = type_from_handle(klass);
1298 return ((TBaseClass*)cr->GetListOfBases()->At((int)ibase))->GetName();
1299}
1300
1302{
1303 if (derived == base)
1304 return true;
1305 TClassRef& derived_type = type_from_handle(derived);
1306 TClassRef& base_type = type_from_handle(base);
1307 return derived_type->GetBaseClass(base_type) != 0;
1308}
1309
1311{
1312 TClassRef& cr = type_from_handle(klass);
1313 const std::string& tn = cr->GetName();
1314 if (gSmartPtrTypes.find(tn.substr(0, tn.find("<"))) != gSmartPtrTypes.end())
1315 return true;
1316 return false;
1317}
1318
1320 const std::string& tname, TCppType_t* raw, TCppMethod_t* deref)
1321{
1322 const std::string& rn = ResolveName(tname);
1323 if (gSmartPtrTypes.find(rn.substr(0, rn.find("<"))) != gSmartPtrTypes.end()) {
1324 if (!raw && !deref) return true;
1325
1326 TClassRef& cr = type_from_handle(GetScope(tname));
1327 if (cr.GetClass()) {
1328 TFunction* func = cr->GetMethod("operator->", "");
1329 if (!func) {
1330 gInterpreter->UpdateListOfMethods(cr.GetClass());
1331 func = cr->GetMethod("operator->", "");
1332 }
1333 if (func) {
1334 if (deref) *deref = (TCppMethod_t)new_CallWrapper(func);
1335 if (raw) *raw = GetScope(TClassEdit::ShortType(
1336 func->GetReturnTypeNormalizedName().c_str(), 1));
1337 return (!deref || *deref) && (!raw || *raw);
1338 }
1339 }
1340 }
1341
1342 return false;
1343}
1344
1345void Cppyy::AddSmartPtrType(const std::string& type_name)
1346{
1347 gSmartPtrTypes.insert(ResolveName(type_name));
1348}
1349
1350
1351// type offsets --------------------------------------------------------------
1353 TCppObject_t address, int direction, bool rerror)
1354{
1355// calculate offsets between declared and actual type, up-cast: direction > 0; down-cast: direction < 0
1356 if (derived == base || !(base && derived))
1357 return (ptrdiff_t)0;
1358
1359 TClassRef& cd = type_from_handle(derived);
1360 TClassRef& cb = type_from_handle(base);
1361
1362 if (!cd.GetClass() || !cb.GetClass())
1363 return (ptrdiff_t)0;
1364
1365 ptrdiff_t offset = -1;
1366 if (!(cd->GetClassInfo() && cb->GetClassInfo())) { // gInterpreter requirement
1367 // would like to warn, but can't quite determine error from intentional
1368 // hiding by developers, so only cover the case where we really should have
1369 // had a class info, but apparently don't:
1370 if (cd->IsLoaded()) {
1371 // warn to allow diagnostics
1372 std::ostringstream msg;
1373 msg << "failed offset calculation between " << cb->GetName() << " and " << cd->GetName();
1374 // TODO: propagate this warning to caller w/o use of Python C-API
1375 // PyErr_Warn(PyExc_RuntimeWarning, const_cast<char*>(msg.str().c_str()));
1376 std::cerr << "Warning: " << msg.str() << '\n';
1377 }
1378
1379 // return -1 to signal caller NOT to apply offset
1380 return rerror ? (ptrdiff_t)offset : 0;
1381 }
1382
1383 offset = gInterpreter->ClassInfo_GetBaseOffset(
1384 cd->GetClassInfo(), cb->GetClassInfo(), (void*)address, direction > 0);
1385 if (offset == -1) // Cling error, treat silently
1386 return rerror ? (ptrdiff_t)offset : 0;
1387
1388 return (ptrdiff_t)(direction < 0 ? -offset : offset);
1389}
1390
1391
1392// method/function reflection information ------------------------------------
1394{
1395 if (IsNamespace(scope))
1396 return (TCppIndex_t)0; // enforce lazy
1397
1398 TClassRef& cr = type_from_handle(scope);
1399 if (cr.GetClass() && cr->GetListOfMethods(true)) {
1400 Cppyy::TCppIndex_t nMethods = (TCppIndex_t)cr->GetListOfMethods(false)->GetSize();
1401 if (nMethods == (TCppIndex_t)0) {
1402 std::string clName = GetScopedFinalName(scope);
1403 if (clName.find('<') != std::string::npos) {
1404 // chicken-and-egg problem: TClass does not know about methods until
1405 // instantiation, so force it
1406 if (clName.find("std::", 0, 5) == std::string::npos && \
1407 is_missclassified_stl(clName)) {
1408 // TODO: this is too simplistic for template arguments missing std::
1409 clName = "std::" + clName;
1410 }
1411 std::ostringstream stmt;
1412 stmt << "template class " << clName << ";";
1413 gInterpreter->Declare(stmt.str().c_str());
1414
1415 // now reload the methods
1416 return (TCppIndex_t)cr->GetListOfMethods(true)->GetSize();
1417 }
1418 }
1419 return nMethods;
1420 }
1421
1422 return (TCppIndex_t)0; // unknown class?
1423}
1424
1425std::vector<Cppyy::TCppIndex_t> Cppyy::GetMethodIndicesFromName(
1426 TCppScope_t scope, const std::string& name)
1427{
1428 std::vector<TCppIndex_t> indices;
1429 TClassRef& cr = type_from_handle(scope);
1430 if (cr.GetClass()) {
1431 gInterpreter->UpdateListOfMethods(cr.GetClass());
1432 int imeth = 0;
1433 TFunction* func = nullptr;
1434 TIter next(cr->GetListOfMethods());
1435 while ((func = (TFunction*)next())) {
1436 if (match_name(name, func->GetName())) {
1437 if (func->Property() & kIsPublic)
1438 indices.push_back((TCppIndex_t)imeth);
1439 }
1440 ++imeth;
1441 }
1442 } else if (scope == GLOBAL_HANDLE) {
1443 TCollection* funcs = gROOT->GetListOfGlobalFunctions(true);
1444
1445 // tickle deserialization
1446 if (!funcs->FindObject(name.c_str()))
1447 return indices;
1448
1449 TFunction* func = nullptr;
1450 TIter ifunc(funcs);
1451 while ((func = (TFunction*)ifunc.Next())) {
1452 if (match_name(name, func->GetName()))
1453 indices.push_back((TCppIndex_t)new_CallWrapper(func));
1454 }
1455 }
1456
1457 return indices;
1458}
1459
1461{
1462 TClassRef& cr = type_from_handle(scope);
1463 if (cr.GetClass()) {
1464 TFunction* f = (TFunction*)cr->GetListOfMethods(false)->At((int)idx);
1465 if (f) return (Cppyy::TCppMethod_t)new_CallWrapper(f);
1466 return (Cppyy::TCppMethod_t)nullptr;
1467 }
1468
1469 assert(klass == (Cppyy::TCppType_t)GLOBAL_HANDLE);
1470 return (Cppyy::TCppMethod_t)idx;
1471}
1472
1474{
1475 if (method) {
1476 const std::string& name = ((CallWrapper*)method)->fName;
1477
1478 if (name.compare(0, 8, "operator") != 0)
1479 // strip template instantiation part, if any
1480 return name.substr(0, name.find('<'));
1481 return name;
1482 }
1483 return "<unknown>";
1484}
1485
1487{
1488 if (method) {
1489 std::string name = ((CallWrapper*)method)->fName;
1490 name.erase(std::remove(name.begin(), name.end(), ' '), name.end());
1491 return name;
1492 }
1493 return "<unknown>";
1494}
1495
1497{
1498 if (method)
1499 return m2f(method)->GetMangledName();
1500 return "<unknown>";
1501}
1502
1504{
1505 if (method) {
1506 TFunction* f = m2f(method);
1507 if (f->ExtraProperty() & kIsConstructor)
1508 return "constructor";
1509 std::string restype = f->GetReturnTypeName();
1510 // TODO: this is ugly, but we can't use GetReturnTypeName() for ostreams
1511 // and maybe others, whereas GetReturnTypeNormalizedName() has proven to
1512 // be save in all cases (Note: 'int8_t' covers 'int8_t' and 'uint8_t')
1513 if (restype.find("int8_t") != std::string::npos)
1514 return restype;
1515 restype = f->GetReturnTypeNormalizedName();
1516 if (restype == "(lambda)") {
1517 std::ostringstream s;
1518 // TODO: what if there are parameters to the lambda?
1519 s << "__cling_internal::FT<decltype("
1520 << GetMethodFullName(method) << "(";
1521 for (Cppyy::TCppIndex_t i = 0; i < Cppyy::GetMethodNumArgs(method); ++i) {
1522 if (i != 0) s << ", ";
1523 s << Cppyy::GetMethodArgType(method, i) << "{}";
1524 }
1525 s << "))>::F";
1526 TClass* cl = TClass::GetClass(s.str().c_str());
1527 if (cl) return cl->GetName();
1528 // TODO: signal some type of error (or should that be upstream?
1529 }
1530 return restype;
1531 }
1532 return "<unknown>";
1533}
1534
1536{
1537 if (method)
1538 return m2f(method)->GetNargs();
1539 return 0;
1540}
1541
1543{
1544 if (method) {
1545 TFunction* f = m2f(method);
1546 return (TCppIndex_t)(f->GetNargs() - f->GetNargsOpt());
1547 }
1548 return (TCppIndex_t)0;
1549}
1550
1552{
1553 if (method) {
1554 TFunction* f = m2f(method);
1555 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1556 return arg->GetName();
1557 }
1558 return "<unknown>";
1559}
1560
1562{
1563 if (method) {
1564 TFunction* f = m2f(method);
1565 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1566 return arg->GetTypeNormalizedName();
1567 }
1568 return "<unknown>";
1569}
1570
1572{
1573 if (method) {
1574 TFunction* f = m2f(method);
1575 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At((int)iarg);
1576 const char* def = arg->GetDefault();
1577 if (def)
1578 return def;
1579 }
1580
1581 return "";
1582}
1583
1584std::string Cppyy::GetMethodSignature(TCppMethod_t method, bool show_formalargs, TCppIndex_t maxargs)
1585{
1586 TFunction* f = m2f(method);
1587 if (f) {
1588 std::ostringstream sig;
1589 sig << "(";
1590 int nArgs = f->GetNargs();
1591 if (maxargs != (TCppIndex_t)-1) nArgs = std::min(nArgs, (int)maxargs);
1592 for (int iarg = 0; iarg < nArgs; ++iarg) {
1593 TMethodArg* arg = (TMethodArg*)f->GetListOfMethodArgs()->At(iarg);
1594 sig << arg->GetFullTypeName();
1595 if (show_formalargs) {
1596 const char* argname = arg->GetName();
1597 if (argname && argname[0] != '\0') sig << " " << argname;
1598 const char* defvalue = arg->GetDefault();
1599 if (defvalue && defvalue[0] != '\0') sig << " = " << defvalue;
1600 }
1601 if (iarg != nArgs-1) sig << (show_formalargs ? ", " : ",");
1602 }
1603 sig << ")";
1604 return sig.str();
1605 }
1606 return "<unknown>";
1607}
1608
1609std::string Cppyy::GetMethodPrototype(TCppScope_t scope, TCppMethod_t method, bool show_formalargs)
1610{
1611 std::string scName = GetScopedFinalName(scope);
1612 TFunction* f = m2f(method);
1613 if (f) {
1614 std::ostringstream sig;
1615 sig << f->GetReturnTypeName() << " "
1616 << scName << "::" << f->GetName();
1617 sig << GetMethodSignature(method, show_formalargs);
1618 return sig.str();
1619 }
1620 return "<unknown>";
1621}
1622
1624{
1625 if (method) {
1626 TFunction* f = m2f(method);
1627 return f->Property() & kIsConstMethod;
1628 }
1629 return false;
1630}
1631
1633{
1634 if (scope == (TCppScope_t)GLOBAL_HANDLE) {
1635 TCollection* coll = gROOT->GetListOfFunctionTemplates();
1636 if (coll) return (TCppIndex_t)coll->GetSize();
1637 } else {
1638 TClassRef& cr = type_from_handle(scope);
1639 if (cr.GetClass()) {
1640 TCollection* coll = cr->GetListOfFunctionTemplates(true);
1641 if (coll) return (TCppIndex_t)coll->GetSize();
1642 }
1643 }
1644
1645// failure ...
1646 return (TCppIndex_t)0;
1647}
1648
1650{
1651 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1652 return ((THashList*)gROOT->GetListOfFunctionTemplates())->At((int)imeth)->GetName();
1653 else {
1654 TClassRef& cr = type_from_handle(scope);
1655 if (cr.GetClass())
1656 return cr->GetListOfFunctionTemplates(false)->At((int)imeth)->GetName();
1657 }
1658
1659// failure ...
1660 assert(!"should not be called unless GetNumTemplatedMethods() succeeded");
1661 return "";
1662}
1663
1665{
1666 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1667 return false;
1668
1669 TClassRef& cr = type_from_handle(scope);
1670 if (cr.GetClass()) {
1672 return f->ExtraProperty() & kIsConstructor;
1673 }
1674
1675 return false;
1676}
1677
1678bool Cppyy::ExistsMethodTemplate(TCppScope_t scope, const std::string& name)
1679{
1680 if (scope == (TCppScope_t)GLOBAL_HANDLE)
1681 return (bool)gROOT->GetFunctionTemplate(name.c_str());
1682 else {
1683 TClassRef& cr = type_from_handle(scope);
1684 if (cr.GetClass())
1685 return (bool)cr->GetFunctionTemplate(name.c_str());
1686 }
1687
1688// failure ...
1689 return false;
1690}
1691
1693{
1694 TClassRef& cr = type_from_handle(scope);
1695 if (cr.GetClass()) {
1696 TFunction* f = (TFunction*)cr->GetListOfMethods(false)->At((int)idx);
1697 if (f && strstr(f->GetName(), "<")) return true;
1698 return false;
1699 }
1700
1701 assert(scope == (Cppyy::TCppType_t)GLOBAL_HANDLE);
1702 if (((CallWrapper*)idx)->fName.find('<') != std::string::npos) return true;
1703 return false;
1704}
1705
1706// helpers for Cppyy::GetMethodTemplate()
1707static std::map<TDictionary::DeclId_t, CallWrapper*> gMethodTemplates;
1708
1710 TCppScope_t scope, const std::string& name, const std::string& proto)
1711{
1712// There is currently no clean way of extracting a templated method out of ROOT/meta
1713// for a variety of reasons, none of them fundamental. The game played below is to
1714// first get any pre-existing functions already managed by ROOT/meta, but if that fails,
1715// to do an explicit lookup that ignores the prototype (i.e. the full name should be
1716// enough), and finally to ignore the template arguments part of the name as this fails
1717// in cling if there are default parameters.
1718// It would be possible to get the prototype from the created functions and use that to
1719// do a new lookup, after which ROOT/meta will manage the function. However, neither
1720// TFunction::GetPrototype() nor TFunction::GetSignature() is of the proper form, so
1721// we'll/ manage the new TFunctions instead and will assume that they are cached on the
1722// calling side to prevent multiple creations.
1723 TFunction* func = nullptr; ClassInfo_t* cl = nullptr;
1724 if (scope == (cppyy_scope_t)GLOBAL_HANDLE) {
1725 func = gROOT->GetGlobalFunctionWithPrototype(name.c_str(), proto.c_str());
1726 if (func && name.back() == '>' && name != func->GetName())
1727 func = nullptr; // happens if implicit conversion matches the overload
1728 } else {
1729 TClassRef& cr = type_from_handle(scope);
1730 if (cr.GetClass()) {
1731 func = cr->GetMethodWithPrototype(name.c_str(), proto.c_str());
1732 if (!func) {
1733 cl = cr->GetClassInfo();
1734 // try base classes to cover a common 'using' case (TODO: this is stupid and misses
1735 // out on base classes; fix that with improved access to Cling)
1736 TCppIndex_t nbases = GetNumBases(scope);
1737 for (TCppIndex_t i = 0; i < nbases; ++i) {
1738 TClassRef& base = type_from_handle(GetScope(GetBaseName(scope, i)));
1739 if (base.GetClass()) {
1740 func = base->GetMethodWithPrototype(name.c_str(), proto.c_str());
1741 if (func) break;
1742 }
1743 }
1744 }
1745 }
1746 }
1747
1748 if (!func && name.back() == '>' && (cl || scope == (cppyy_scope_t)GLOBAL_HANDLE)) {
1749 // try again, ignoring proto in case full name is complete template
1750 auto declid = gInterpreter->GetFunction(cl, name.c_str());
1751 if (declid) {
1752 auto existing = gMethodTemplates.find(declid);
1753 if (existing == gMethodTemplates.end()) {
1754 auto cw = new_CallWrapper(declid, name);
1755 existing = gMethodTemplates.insert(std::make_pair(declid, cw)).first;
1756 }
1757 return (TCppMethod_t)existing->second;
1758 }
1759 }
1760
1761 if (func) {
1762 // make sure we didn't match a non-templated overload
1763 if (func->ExtraProperty() & kIsTemplateSpec)
1764 return (TCppMethod_t)new_CallWrapper(func);
1765
1766 // disregard this non-templated method as it will be considered when appropriate
1767 return (TCppMethod_t)nullptr;
1768 }
1769
1770// try again with template arguments removed from name, if applicable
1771 if (name.back() == '>') {
1772 auto pos = name.find('<');
1773 if (pos != std::string::npos) {
1774 TCppMethod_t cppmeth = GetMethodTemplate(scope, name.substr(0, pos), proto);
1775 if (cppmeth) {
1776 // allow if requested template names match up to the result
1777 const std::string& alt = GetMethodFullName(cppmeth);
1778 if (name.size() < alt.size() && alt.find('<') == pos) {
1779 const std::string& partial = name.substr(pos, name.size()-1-pos);
1780 if (strncmp(partial.c_str(), alt.substr(pos, alt.size()-1-pos).c_str(), partial.size()) == 0)
1781 return cppmeth;
1782 }
1783 }
1784 }
1785 }
1786
1787// failure ...
1788 return (TCppMethod_t)nullptr;
1789}
1790
1791static inline
1792std::string type_remap(const std::string& n1, const std::string& n2)
1793{
1794// Operator lookups of (C++ string, Python str) should succeeded, for the combos of
1795// string/str, wstring/str, string/unicode and wstring/unicode; since C++ does not have a
1796// operator+(std::string, std::wstring), we'll have to look up the same type and rely on
1797// the converters in CPyCppyy/_cppyy.
1798 if (n1 == "str") {
1799 if (n2 == "std::basic_string<wchar_t,std::char_traits<wchar_t>,std::allocator<wchar_t> >")
1800 return n2; // match like for like
1801 return "std::string"; // probably best bet
1802 } else if (n1 == "float")
1803 return "double"; // debatable, but probably intended
1804 return n1;
1805}
1806
1808 TCppType_t scope, const std::string& lc, const std::string& rc, const std::string& opname)
1809{
1810// Find a global operator function with a matching signature; prefer by-ref, but
1811// fall back on by-value if that fails.
1812 std::string lcname1 = TClassEdit::CleanType(lc.c_str());
1813 const std::string& rcname = rc.empty() ? rc : type_remap(TClassEdit::CleanType(rc.c_str()), lcname1);
1814 const std::string& lcname = type_remap(lcname1, rcname);
1815
1816 std::string proto = lcname + "&" + (rc.empty() ? rc : (", " + rcname + "&"));
1817 if (scope == (cppyy_scope_t)GLOBAL_HANDLE) {
1818 TFunction* func = gROOT->GetGlobalFunctionWithPrototype(opname.c_str(), proto.c_str());
1819 if (func) return (TCppIndex_t)new_CallWrapper(func);
1820 proto = lcname + (rc.empty() ? rc : (", " + rcname));
1821 func = gROOT->GetGlobalFunctionWithPrototype(opname.c_str(), proto.c_str());
1822 if (func) return (TCppIndex_t)new_CallWrapper(func);
1823 } else {
1824 TClassRef& cr = type_from_handle(scope);
1825 if (cr.GetClass()) {
1826 TFunction* func = cr->GetMethodWithPrototype(opname.c_str(), proto.c_str());
1827 if (func) return (TCppIndex_t)cr->GetListOfMethods()->IndexOf(func);
1828 proto = lcname + (rc.empty() ? rc : (", " + rcname));
1829 func = cr->GetMethodWithPrototype(opname.c_str(), proto.c_str());
1830 if (func) return (TCppIndex_t)cr->GetListOfMethods()->IndexOf(func);
1831 }
1832 }
1833
1834// failure ...
1835 return (TCppIndex_t)-1;
1836}
1837
1838// method properties ---------------------------------------------------------
1840{
1841 if (method) {
1842 TFunction* f = m2f(method);
1843 return f->Property() & kIsPublic;
1844 }
1845 return false;
1846}
1847
1849{
1850 if (method) {
1851 TFunction* f = m2f(method);
1852 return f->Property() & kIsProtected;
1853 }
1854 return false;
1855}
1856
1858{
1859 if (method) {
1860 TFunction* f = m2f(method);
1861 return f->ExtraProperty() & kIsConstructor;
1862 }
1863 return false;
1864}
1865
1867{
1868 if (method) {
1869 TFunction* f = m2f(method);
1870 return f->ExtraProperty() & kIsDestructor;
1871 }
1872 return false;
1873}
1874
1876{
1877 if (method) {
1878 TFunction* f = m2f(method);
1879 return f->Property() & kIsStatic;
1880 }
1881 return false;
1882}
1883
1884// data member reflection information ----------------------------------------
1886{
1887 if (IsNamespace(scope))
1888 return (TCppIndex_t)0; // enforce lazy
1889
1890 TClassRef& cr = type_from_handle(scope);
1891 if (cr.GetClass()) {
1893 if (cr->GetListOfDataMembers())
1894 sum = cr->GetListOfDataMembers()->GetSize();
1895 if (cr->GetListOfUsingDataMembers())
1897 return sum;
1898 }
1899
1900 return (TCppIndex_t)0; // unknown class?
1901}
1902
1904{
1905 if (!cr.GetClass() || !cr->GetListOfDataMembers())
1906 return nullptr;
1907
1908 int numDMs = cr->GetListOfDataMembers()->GetSize();
1909 if ((int)idata < numDMs)
1910 return (TDataMember*)cr->GetListOfDataMembers()->At((int)idata);
1911 return (TDataMember*)cr->GetListOfUsingDataMembers()->At((int)idata - numDMs);
1912}
1913
1915{
1916 TClassRef& cr = type_from_handle(scope);
1917 if (cr.GetClass()) {
1918 TDataMember *m = GetDataMemberByIndex(cr, (int)idata);
1919 return m->GetName();
1920 }
1921 assert(scope == GLOBAL_HANDLE);
1922 TGlobal* gbl = g_globalvars[idata];
1923 return gbl->GetName();
1924}
1925
1927{
1928 if (scope == GLOBAL_HANDLE) {
1929 TGlobal* gbl = g_globalvars[idata];
1930 std::string fullType = gbl->GetFullTypeName();
1931
1932 if ((int)gbl->GetArrayDim() > 1)
1933 fullType.append("*");
1934 else if ((int)gbl->GetArrayDim() == 1) {
1935 std::ostringstream s;
1936 s << '[' << gbl->GetMaxIndex(0) << ']' << std::ends;
1937 fullType.append(s.str());
1938 }
1939 return fullType;
1940 }
1941
1942 TClassRef& cr = type_from_handle(scope);
1943 if (cr.GetClass()) {
1944 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
1945 // TODO: fix this upstream. Usually, we want m->GetFullTypeName(), because it does
1946 // not resolve typedefs, but it looses scopes for inner classes/structs, so in that
1947 // case m->GetTrueTypeName() should be used (this also cleans up the cases where
1948 // the "full type" retains spurious "struct" or "union" in the name).
1949 std::string fullType = m->GetFullTypeName();
1950 if (fullType != m->GetTrueTypeName()) {
1951 const std::string& trueName = m->GetTrueTypeName();
1952 if (fullType.find("::") == std::string::npos && trueName.find("::") != std::string::npos)
1953 fullType = trueName;
1954 }
1955
1956 if ((int)m->GetArrayDim() > 1 || (!m->IsBasic() && m->IsaPointer()))
1957 fullType.append("*");
1958 else if ((int)m->GetArrayDim() == 1) {
1959 std::ostringstream s;
1960 s << '[' << m->GetMaxIndex(0) << ']' << std::ends;
1961 fullType.append(s.str());
1962 }
1963 return fullType;
1964 }
1965
1966 return "<unknown>";
1967}
1968
1970{
1971 if (scope == GLOBAL_HANDLE) {
1972 TGlobal* gbl = g_globalvars[idata];
1973 if (!gbl->GetAddress() || gbl->GetAddress() == (void*)-1) {
1974 // CLING WORKAROUND: make sure variable is loaded
1975 intptr_t addr = (intptr_t)gInterpreter->ProcessLine((std::string("&")+gbl->GetName()+";").c_str());
1976 if (gbl->GetAddress() && gbl->GetAddress() != (void*)-1)
1977 return (intptr_t)gbl->GetAddress(); // now loaded!
1978 return addr; // last resort ...
1979 }
1980 return (intptr_t)gbl->GetAddress();
1981 }
1982
1983 TClassRef& cr = type_from_handle(scope);
1984 if (cr.GetClass()) {
1985 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
1986 // CLING WORKAROUND: the following causes templates to be instantiated first within the proper
1987 // scope, making the lookup succeed and preventing spurious duplicate instantiations later. Also,
1988 // if the variable is not yet loaded, pull it in through gInterpreter.
1989 if (m->Property() & kIsStatic) {
1990 if (strchr(cr->GetName(), '<'))
1991 gInterpreter->ProcessLine(((std::string)cr->GetName()+"::"+m->GetName()+";").c_str());
1992 if ((intptr_t)m->GetOffsetCint() == (intptr_t)-1)
1993 return (intptr_t)gInterpreter->ProcessLine((std::string("&")+cr->GetName()+"::"+m->GetName()+";").c_str());
1994 }
1995 return (intptr_t)m->GetOffsetCint(); // yes, CINT (GetOffset() is both wrong
1996 // and caches that wrong result!
1997 }
1998
1999 return (intptr_t)-1;
2000}
2001
2003{
2004 if (scope == GLOBAL_HANDLE) {
2005 TGlobal* gb = (TGlobal*)gROOT->GetListOfGlobals(false /* load */)->FindObject(name.c_str());
2006 if (!gb) gb = (TGlobal*)gROOT->GetListOfGlobals(true /* load */)->FindObject(name.c_str());
2007 if (!gb) {
2008 // some enums are not loaded as they are not considered part of
2009 // the global scope, but of the enum scope; get them w/o checking
2010 TDictionary::DeclId_t did = gInterpreter->GetDataMember(nullptr, name.c_str());
2011 if (did) {
2012 DataMemberInfo_t* t = gInterpreter->DataMemberInfo_Factory(did, nullptr);
2013 ((TListOfDataMembers*)gROOT->GetListOfGlobals())->Get(t, true);
2014 gb = (TGlobal*)gROOT->GetListOfGlobals(false /* load */)->FindObject(name.c_str());
2015 }
2016 }
2017
2018 if (gb && strcmp(gb->GetFullTypeName(), "(lambda)") == 0) {
2019 // lambdas use a compiler internal closure type, so we wrap
2020 // them, then return the wrapper's type
2021 // TODO: this current leaks the std::function; also, if possible,
2022 // should instantiate through TClass rather then ProcessLine
2023 std::ostringstream s;
2024 s << "auto __cppyy_internal_wrap_" << name << " = "
2025 "new __cling_internal::FT<decltype(" << name << ")>::F"
2026 "{" << name << "};";
2027 gInterpreter->ProcessLine(s.str().c_str());
2028 TGlobal* wrap = (TGlobal*)gROOT->GetListOfGlobals(true)->FindObject(
2029 ("__cppyy_internal_wrap_"+name).c_str());
2030 if (wrap && wrap->GetAddress()) gb = wrap;
2031 }
2032
2033 if (gb) {
2034 // TODO: do we ever need a reverse lookup?
2035 g_globalvars.push_back(gb);
2036 return TCppIndex_t(g_globalvars.size() - 1);
2037 }
2038
2039 } else {
2040 TClassRef& cr = type_from_handle(scope);
2041 if (cr.GetClass()) {
2042 TDataMember* dm =
2044 // TODO: turning this into an index is silly ...
2045 if (dm) return (TCppIndex_t)cr->GetListOfDataMembers()->IndexOf(dm);
2047 if (dm)
2048 return (TCppIndex_t)cr->GetListOfDataMembers()->IndexOf(dm)
2049 + cr->GetListOfDataMembers()->GetSize();
2050 }
2051 }
2052
2053 return (TCppIndex_t)-1;
2054}
2055
2056
2057// data member properties ----------------------------------------------------
2059{
2060 if (scope == GLOBAL_HANDLE)
2061 return true;
2062 TClassRef& cr = type_from_handle(scope);
2063 if (cr->Property() & kIsNamespace)
2064 return true;
2065 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2066 return m->Property() & kIsPublic;
2067}
2068
2070{
2071 if (scope == GLOBAL_HANDLE)
2072 return true;
2073 TClassRef& cr = type_from_handle(scope);
2074 if (cr->Property() & kIsNamespace)
2075 return true;
2076 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2077 return m->Property() & kIsProtected;
2078}
2079
2081{
2082 if (scope == GLOBAL_HANDLE)
2083 return true;
2084 TClassRef& cr = type_from_handle(scope);
2085 if (cr->Property() & kIsNamespace)
2086 return true;
2087 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2088 return m->Property() & kIsStatic;
2089}
2090
2092{
2093 if (scope == GLOBAL_HANDLE) {
2094 TGlobal* gbl = g_globalvars[idata];
2095 return gbl->Property() & kIsConstant;
2096 }
2097 TClassRef& cr = type_from_handle(scope);
2098 if (cr.GetClass()) {
2099 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2100 return m->Property() & kIsConstant;
2101 }
2102 return false;
2103}
2104
2106{
2107// TODO: currently, ROOT/meta does not properly distinguish between variables of enum
2108// type, and values of enums. The latter are supposed to be const. This code relies on
2109// odd features (bugs?) to figure out the difference, but this should really be fixed
2110// upstream and/or deserves a new API.
2111
2112 if (scope == GLOBAL_HANDLE) {
2113 TGlobal* gbl = g_globalvars[idata];
2114
2115 // make use of an oddity: enum global variables do not have their kIsStatic bit
2116 // set, whereas enum global values do
2117 return (gbl->Property() & kIsEnum) && (gbl->Property() & kIsStatic);
2118 }
2119
2120 TClassRef& cr = type_from_handle(scope);
2121 if (cr.GetClass()) {
2122 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2123 std::string ti = m->GetTypeName();
2124
2125 // can't check anonymous enums by type name, so just accept them as enums
2126 if (ti.rfind("(unnamed)") != std::string::npos)
2127 return m->Property() & kIsEnum;
2128
2129 // since there seems to be no distinction between data of enum type and enum values,
2130 // check the list of constants for the type to see if there's a match
2131 if (ti.rfind(cr->GetName(), 0) != std::string::npos) {
2132 std::string::size_type s = strlen(cr->GetName())+2;
2133 if (s < ti.size()) {
2134 TEnum* ee = ((TListOfEnums*)cr->GetListOfEnums())->GetObject(ti.substr(s, std::string::npos).c_str());
2135 if (ee) return ee->GetConstant(m->GetName());
2136 }
2137 }
2138 }
2139
2140// this default return only means that the data will be writable, not that it will
2141// be unreadable or otherwise misrepresented
2142 return false;
2143}
2144
2145int Cppyy::GetDimensionSize(TCppScope_t scope, TCppIndex_t idata, int dimension)
2146{
2147 if (scope == GLOBAL_HANDLE) {
2148 TGlobal* gbl = g_globalvars[idata];
2149 return gbl->GetMaxIndex(dimension);
2150 }
2151 TClassRef& cr = type_from_handle(scope);
2152 if (cr.GetClass()) {
2153 TDataMember* m = GetDataMemberByIndex(cr, (int)idata);
2154 return m->GetMaxIndex(dimension);
2155 }
2156 return -1;
2157}
2158
2159
2160// enum properties -----------------------------------------------------------
2161Cppyy::TCppEnum_t Cppyy::GetEnum(TCppScope_t scope, const std::string& enum_name)
2162{
2163 if (scope == GLOBAL_HANDLE)
2164 return (TCppEnum_t)gROOT->GetListOfEnums(kTRUE)->FindObject(enum_name.c_str());
2165
2166 TClassRef& cr = type_from_handle(scope);
2167 if (cr.GetClass())
2168 return (TCppEnum_t)cr->GetListOfEnums(kTRUE)->FindObject(enum_name.c_str());
2169
2170 return (TCppEnum_t)0;
2171}
2172
2174{
2175 return (TCppIndex_t)((TEnum*)etype)->GetConstants()->GetSize();
2176}
2177
2179{
2180 return ((TEnumConstant*)((TEnum*)etype)->GetConstants()->At(idata))->GetName();
2181}
2182
2184{
2185 TEnumConstant* ecst = (TEnumConstant*)((TEnum*)etype)->GetConstants()->At(idata);
2186 return (long long)ecst->GetValue();
2187}
2188
2189
2190//- C-linkage wrappers -------------------------------------------------------
2191
2192extern "C" {
2193/* direct interpreter access ---------------------------------------------- */
2194int cppyy_compile(const char* code) {
2195 return Cppyy::Compile(code);
2196}
2197
2198
2199/* name to opaque C++ scope representation -------------------------------- */
2200char* cppyy_resolve_name(const char* cppitem_name) {
2201 return cppstring_to_cstring(Cppyy::ResolveName(cppitem_name));
2202}
2203
2204char* cppyy_resolve_enum(const char* enum_type) {
2205 return cppstring_to_cstring(Cppyy::ResolveEnum(enum_type));
2206}
2207
2208cppyy_scope_t cppyy_get_scope(const char* scope_name) {
2209 return cppyy_scope_t(Cppyy::GetScope(scope_name));
2210}
2211
2213 return cppyy_type_t(Cppyy::GetActualClass(klass, (void*)obj));
2214}
2215
2217 return Cppyy::SizeOf(klass);
2218}
2219
2220size_t cppyy_size_of_type(const char* type_name) {
2221 return Cppyy::SizeOf(type_name);
2222}
2223
2224
2225/* memory management ------------------------------------------------------ */
2228}
2229
2231 Cppyy::Deallocate(type, (void*)self);
2232}
2233
2236}
2237
2239 Cppyy::Destruct(type, (void*)self);
2240}
2241
2242
2243/* method/function dispatching -------------------------------------------- */
2244/* Exception types:
2245 1: default (unknown exception)
2246 2: standard exception
2247*/
2248#define CPPYY_HANDLE_EXCEPTION \
2249 catch (std::exception& e) { \
2250 cppyy_exctype_t* etype = (cppyy_exctype_t*)((Parameter*)args+nargs); \
2251 *etype = (cppyy_exctype_t)2; \
2252 *((char**)(etype+1)) = cppstring_to_cstring(e.what()); \
2253 } \
2254 catch (...) { \
2255 cppyy_exctype_t* etype = (cppyy_exctype_t*)((Parameter*)args+nargs); \
2256 *etype = (cppyy_exctype_t)1; \
2257 *((char**)(etype+1)) = \
2258 cppstring_to_cstring("unhandled, unknown C++ exception"); \
2259 }
2260
2261void cppyy_call_v(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2262 try {
2263 Cppyy::CallV(method, (void*)self, nargs, args);
2265}
2266
2267unsigned char cppyy_call_b(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2268 try {
2269 return (unsigned char)Cppyy::CallB(method, (void*)self, nargs, args);
2271 return (unsigned char)-1;
2272}
2273
2274char cppyy_call_c(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2275 try {
2276 return (char)Cppyy::CallC(method, (void*)self, nargs, args);
2278 return (char)-1;
2279}
2280
2281short cppyy_call_h(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2282 try {
2283 return (short)Cppyy::CallH(method, (void*)self, nargs, args);
2285 return (short)-1;
2286}
2287
2288int cppyy_call_i(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2289 try {
2290 return (int)Cppyy::CallI(method, (void*)self, nargs, args);
2292 return (int)-1;
2293}
2294
2295long cppyy_call_l(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2296 try {
2297 return (long)Cppyy::CallL(method, (void*)self, nargs, args);
2299 return (long)-1;
2300}
2301
2302long long cppyy_call_ll(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2303 try {
2304 return (long long)Cppyy::CallLL(method, (void*)self, nargs, args);
2306 return (long long)-1;
2307}
2308
2309float cppyy_call_f(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2310 try {
2311 return (float)Cppyy::CallF(method, (void*)self, nargs, args);
2313 return (float)-1;
2314}
2315
2316double cppyy_call_d(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2317 try {
2318 return (double)Cppyy::CallD(method, (void*)self, nargs, args);
2320 return (double)-1;
2321}
2322
2323long double cppyy_call_ld(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2324 try {
2325 return (long double)Cppyy::CallLD(method, (void*)self, nargs, args);
2327 return (long double)-1;
2328}
2329
2330double cppyy_call_nld(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2331 return (double)cppyy_call_ld(method, self, nargs, args);
2332}
2333
2334void* cppyy_call_r(cppyy_method_t method, cppyy_object_t self, int nargs, void* args) {
2335 try {
2336 return (void*)Cppyy::CallR(method, (void*)self, nargs, args);
2338 return (void*)nullptr;
2339}
2340
2342 cppyy_method_t method, cppyy_object_t self, int nargs, void* args, size_t* lsz) {
2343 try {
2344 return Cppyy::CallS(method, (void*)self, nargs, args, lsz);
2346 return (char*)nullptr;
2347}
2348
2350 cppyy_method_t method, cppyy_type_t klass, int nargs, void* args) {
2351 try {
2352 return cppyy_object_t(Cppyy::CallConstructor(method, klass, nargs, args));
2354 return (cppyy_object_t)0;
2355}
2356
2358 Cppyy::CallDestructor(klass, self);
2359}
2360
2362 int nargs, void* args, cppyy_type_t result_type) {
2363 try {
2364 return cppyy_object_t(Cppyy::CallO(method, (void*)self, nargs, args, result_type));
2366 return (cppyy_object_t)0;
2367}
2368
2370 return cppyy_funcaddr_t(Cppyy::GetFunctionAddress(method, true));
2371}
2372
2373
2374/* handling of function argument buffer ----------------------------------- */
2376// for calls through C interface, require extra space for reporting exceptions
2377 return malloc(nargs*sizeof(Parameter)+sizeof(cppyy_exctype_t)+sizeof(char**));
2378}
2379
2381 free(args);
2382}
2383
2385 return (size_t)Cppyy::GetFunctionArgSizeof();
2386}
2387
2389 return (size_t)Cppyy::GetFunctionArgTypeoffset();
2390}
2391
2392
2393/* scope reflection information ------------------------------------------- */
2395 return (int)Cppyy::IsNamespace(scope);
2396}
2397
2398int cppyy_is_template(const char* template_name) {
2399 return (int)Cppyy::IsTemplate(template_name);
2400}
2401
2403 return (int)Cppyy::IsAbstract(type);
2404}
2405
2406int cppyy_is_enum(const char* type_name) {
2407 return (int)Cppyy::IsEnum(type_name);
2408}
2409
2410const char** cppyy_get_all_cpp_names(cppyy_scope_t scope, size_t* count) {
2411 std::set<std::string> cppnames;
2412 Cppyy::GetAllCppNames(scope, cppnames);
2413 const char** c_cppnames = (const char**)malloc(cppnames.size()*sizeof(const char*));
2414 int i = 0;
2415 for (const auto& name : cppnames) {
2416 c_cppnames[i] = cppstring_to_cstring(name);
2417 ++i;
2418 }
2419 *count = cppnames.size();
2420 return c_cppnames;
2421}
2422
2423
2424/* namespace reflection information --------------------------------------- */
2426 const std::vector<Cppyy::TCppScope_t>& uv = Cppyy::GetUsingNamespaces((Cppyy::TCppScope_t)scope);
2427
2428 if (uv.empty())
2429 return (cppyy_index_t*)nullptr;
2430
2431 cppyy_scope_t* llresult = (cppyy_scope_t*)malloc(sizeof(cppyy_scope_t)*(uv.size()+1));
2432 for (int i = 0; i < (int)uv.size(); ++i) llresult[i] = uv[i];
2433 llresult[uv.size()] = (cppyy_scope_t)0;
2434 return llresult;
2435}
2436
2437
2438/* class reflection information ------------------------------------------- */
2441}
2442
2445}
2446
2448 return (int)Cppyy::HasVirtualDestructor(type);
2449}
2450
2452 return (int)Cppyy::HasComplexHierarchy(type);
2453}
2454
2456 return (int)Cppyy::GetNumBases(type);
2457}
2458
2461}
2462
2463char* cppyy_base_name(cppyy_type_t type, int base_index) {
2464 return cppstring_to_cstring(Cppyy::GetBaseName (type, base_index));
2465}
2466
2468 return (int)Cppyy::IsSubtype(derived, base);
2469}
2470
2472 return (int)Cppyy::IsSmartPtr(type);
2473}
2474
2475int cppyy_smartptr_info(const char* name, cppyy_type_t* raw, cppyy_method_t* deref) {
2476 return (int)Cppyy::GetSmartPtrInfo(name, raw, deref);
2477}
2478
2479void cppyy_add_smartptr_type(const char* type_name) {
2480 Cppyy::AddSmartPtrType(type_name);
2481}
2482
2483
2484/* calculate offsets between declared and actual type, up-cast: direction > 0; down-cast: direction < 0 */
2485ptrdiff_t cppyy_base_offset(cppyy_type_t derived, cppyy_type_t base, cppyy_object_t address, int direction) {
2486 return (ptrdiff_t)Cppyy::GetBaseOffset(derived, base, (void*)address, direction, 0);
2487}
2488
2489
2490/* method/function reflection information --------------------------------- */
2492 return (int)Cppyy::GetNumMethods(scope);
2493}
2494
2496{
2497 std::vector<cppyy_index_t> result = Cppyy::GetMethodIndicesFromName(scope, name);
2498
2499 if (result.empty())
2500 return (cppyy_index_t*)nullptr;
2501
2502 cppyy_index_t* llresult = (cppyy_index_t*)malloc(sizeof(cppyy_index_t)*(result.size()+1));
2503 for (int i = 0; i < (int)result.size(); ++i) llresult[i] = result[i];
2504 llresult[result.size()] = -1;
2505 return llresult;
2506}
2507
2509 return cppyy_method_t(Cppyy::GetMethod(scope, idx));
2510}
2511
2514}
2515
2518}
2519
2522}
2523
2526}
2527
2529 return (int)Cppyy::GetMethodNumArgs((Cppyy::TCppMethod_t)method);
2530}
2531
2533 return (int)Cppyy::GetMethodReqArgs((Cppyy::TCppMethod_t)method);
2534}
2535
2536char* cppyy_method_arg_name(cppyy_method_t method, int arg_index) {
2538}
2539
2540char* cppyy_method_arg_type(cppyy_method_t method, int arg_index) {
2542}
2543
2544char* cppyy_method_arg_default(cppyy_method_t method, int arg_index) {
2546}
2547
2548char* cppyy_method_signature(cppyy_method_t method, int show_formalargs) {
2549 return cppstring_to_cstring(Cppyy::GetMethodSignature((Cppyy::TCppMethod_t)method, (bool)show_formalargs));
2550}
2551
2552char* cppyy_method_signature_max(cppyy_method_t method, int show_formalargs, int maxargs) {
2553 return cppstring_to_cstring(Cppyy::GetMethodSignature((Cppyy::TCppMethod_t)method, (bool)show_formalargs, (Cppyy::TCppIndex_t)maxargs));
2554}
2555
2556char* cppyy_method_prototype(cppyy_scope_t scope, cppyy_method_t method, int show_formalargs) {
2558 (Cppyy::TCppScope_t)scope, (Cppyy::TCppMethod_t)method, (bool)show_formalargs));
2559}
2560
2562 return (int)Cppyy::IsConstMethod(method);
2563}
2564
2566 return (int)Cppyy::GetNumTemplatedMethods(scope);
2567}
2568
2571}
2572
2575}
2576
2578 return (int)Cppyy::ExistsMethodTemplate(scope, name);
2579}
2580
2582 return (int)Cppyy::IsMethodTemplate(scope, idx);
2583}
2584
2587}
2588
2591}
2592
2593
2594/* method properties ------------------------------------------------------ */
2596 return (int)Cppyy::IsPublicMethod((Cppyy::TCppMethod_t)method);
2597}
2598
2600 return (int)Cppyy::IsProtectedMethod((Cppyy::TCppMethod_t)method);
2601}
2602
2604 return (int)Cppyy::IsConstructor((Cppyy::TCppMethod_t)method);
2605}
2606
2608 return (int)Cppyy::IsDestructor((Cppyy::TCppMethod_t)method);
2609}
2610
2612 return (int)Cppyy::IsStaticMethod((Cppyy::TCppMethod_t)method);
2613}
2614
2615
2616/* data member reflection information ------------------------------------- */
2618 return (int)Cppyy::GetNumDatamembers(scope);
2619}
2620
2621char* cppyy_datamember_name(cppyy_scope_t scope, int datamember_index) {
2622 return cppstring_to_cstring(Cppyy::GetDatamemberName(scope, datamember_index));
2623}
2624
2625char* cppyy_datamember_type(cppyy_scope_t scope, int datamember_index) {
2626 return cppstring_to_cstring(Cppyy::GetDatamemberType(scope, datamember_index));
2627}
2628
2629intptr_t cppyy_datamember_offset(cppyy_scope_t scope, int datamember_index) {
2630 return intptr_t(Cppyy::GetDatamemberOffset(scope, datamember_index));
2631}
2632
2634 return (int)Cppyy::GetDatamemberIndex(scope, name);
2635}
2636
2637
2638
2639/* data member properties ------------------------------------------------- */
2641 return (int)Cppyy::IsPublicData(type, datamember_index);
2642}
2643
2645 return (int)Cppyy::IsProtectedData(type, datamember_index);
2646}
2647
2649 return (int)Cppyy::IsStaticData(type, datamember_index);
2650}
2651
2653 return (int)Cppyy::IsConstData(scope, idata);
2654}
2655
2657 return (int)Cppyy::IsEnumData(scope, idata);
2658}
2659
2660int cppyy_get_dimension_size(cppyy_scope_t scope, cppyy_index_t idata, int dimension) {
2661 return Cppyy::GetDimensionSize(scope, idata, dimension);
2662}
2663
2664
2665/* misc helpers ----------------------------------------------------------- */
2667void* cppyy_load_dictionary(const char* lib_name) {
2668 int result = gSystem->Load(lib_name);
2669 return (void*)(result == 0 /* success */ || result == 1 /* already loaded */);
2670}
2671
2672#if defined(_MSC_VER)
2673long long cppyy_strtoll(const char* str) {
2674 return _strtoi64(str, NULL, 0);
2675}
2676
2677extern "C" {
2678unsigned long long cppyy_strtoull(const char* str) {
2679 return _strtoui64(str, NULL, 0);
2680}
2681}
2682#else
2683long long cppyy_strtoll(const char* str) {
2684 return strtoll(str, NULL, 0);
2685}
2686
2687extern "C" {
2688unsigned long long cppyy_strtoull(const char* str) {
2689 return strtoull(str, NULL, 0);
2690}
2691}
2692#endif
2693
2694void cppyy_free(void* ptr) {
2695 free(ptr);
2696}
2697
2698cppyy_object_t cppyy_charp2stdstring(const char* str, size_t sz) {
2699 return (cppyy_object_t)new std::string(str, sz);
2700}
2701
2702const char* cppyy_stdstring2charp(cppyy_object_t ptr, size_t* lsz) {
2703 *lsz = ((std::string*)ptr)->size();
2704 return ((std::string*)ptr)->data();
2705}
2706
2708 return (cppyy_object_t)new std::string(*(std::string*)ptr);
2709}
2710
2712 return (double)*(long double*)p;
2713}
2714
2715void cppyy_double2longdouble(double d, void* p) {
2716 *(long double*)p = d;
2717}
2718
2720 return (int)(*(std::vector<bool>*)ptr)[idx];
2721}
2722
2724 (*(std::vector<bool>*)ptr)[idx] = (bool)value;
2725}
2726
2727} // end C-linkage wrappers
#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:80
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
@ 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.cxx:41
Int_t gErrorIgnoreLevel
Definition TError.cxx:43
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.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h length
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void funcs
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
#define gInterpreter
#define gROOT
Definition TROOT.h:405
R__EXTERN TSystem * gSystem
Definition TSystem.h:560
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:17502
#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 override
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:81
TList * GetListOfUsingDataMembers(Bool_t load=kTRUE)
Return list containing the TDataMembers of using declarations of a class.
Definition TClass.cxx:3786
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:4978
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:4411
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:4456
const TList * GetListOfAllPublicMethods(Bool_t load=kTRUE)
Returns a list of all public methods of this class and its base classes.
Definition TClass.cxx:3845
void Destructor(void *obj, Bool_t dtorOnly=kFALSE)
Explicitly call destructor for object.
Definition TClass.cxx:5400
TList * GetListOfFunctionTemplates(Bool_t load=kTRUE)
Return TListOfFunctionTemplates for a class.
Definition TClass.cxx:3798
TList * GetListOfEnums(Bool_t load=kTRUE)
Return a list containing the TEnums of a class.
Definition TClass.cxx:3686
TList * GetListOfMethods(Bool_t load=kTRUE)
Return list containing the TMethods of a class.
Definition TClass.cxx:3812
TClass * GetBaseClass(const char *classname)
Return pointer to the base class "classname".
Definition TClass.cxx:2655
TList * GetListOfDataMembers(Bool_t load=kTRUE)
Return list containing the TDataMembers of a class.
Definition TClass.cxx:3770
TList * GetListOfBases()
Return list containing the TBaseClass(es) of a class.
Definition TClass.cxx:3636
Bool_t IsLoaded() const
Return true if the shared library of this class is currently in the a process's memory.
Definition TClass.cxx:5912
ClassInfo_t * GetClassInfo() const
Definition TClass.h:431
Long_t ClassProperty() const
Return the C++ property of this class, eg.
Definition TClass.cxx:2396
Long_t Property() const override
Returns the properties of the TClass as a bit field stored as a Long_t value.
Definition TClass.cxx:6086
ROOT::DelFunc_t GetDelete() const
Return the wrapper around delete ThiObject.
Definition TClass.cxx:7463
TClass * GetActualClass(const void *object) const
Return a pointer to the real class of the object.
Definition TClass.cxx:2607
TFunctionTemplate * GetFunctionTemplate(const char *name)
Definition TClass.cxx:3607
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:2968
Collection abstract base class.
Definition TCollection.h:65
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
const char * GetFullTypeName() const
Get full type description of typedef, e,g.: "class TDirectory*".
Long_t Property() const override
Get property description word. For meaning of bits see EProperty.
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:64
Definition TEnv.h:86
const char * GetValue() const
Definition TEnv.h:110
const char * GetName() const override
Returns name of object.
Definition TEnv.h:109
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 override
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:89
virtual Int_t GetMaxIndex(Int_t dim) const
Return maximum index for array dimension "dim".
Definition TGlobal.cxx:105
Long_t Property() const override
Get property description word. For meaning of bits see EProperty.
Definition TGlobal.cxx:152
virtual void * GetAddress() const
Return address of global.
Definition TGlobal.cxx:81
virtual const char * GetFullTypeName() const
Get full type description of global variable, e,g.: "class TDirectory*".
Definition TGlobal.cxx:124
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...
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:578
TObject * At(Int_t idx) const override
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.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:403
static Bool_t Initialized()
Return kTRUE if the TROOT object has been initialized.
Definition TROOT.cxx:2849
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:1858
virtual void Exit(int code, Bool_t mode=kTRUE)
Exit the application.
Definition TSystem.cxx:719
virtual void StackTrace()
Print a stack trace.
Definition TSystem.cxx:737
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)
void * cppyy_load_dictionary(const char *lib_name)
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)
int cppyy_num_bases_longest_branch(cppyy_type_t type)
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)
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)
Cppyy::TCppIndex_t GetLongestInheritancePath(TClass *klass)
Retrieve number of base classes in the longest branch of the inheritance tree of the input class.
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)
size_t TCppIndex_t
Definition cpp_cppyy.h:24
RPY_EXPORTED TCppIndex_t GetNumMethods(TCppScope_t scope)
RPY_EXPORTED std::string GetMethodMangledName(TCppMethod_t)
RPY_EXPORTED TCppObject_t CallO(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args, TCppType_t result_type)
RPY_EXPORTED int CallI(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED ptrdiff_t GetBaseOffset(TCppType_t derived, TCppType_t base, TCppObject_t address, int direction, bool rerror=false)
RPY_EXPORTED void DeallocateFunctionArgs(void *args)
RPY_EXPORTED bool IsEnumData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED unsigned char CallB(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED bool IsAbstract(TCppType_t type)
RPY_EXPORTED size_t SizeOf(TCppType_t klass)
RPY_EXPORTED TCppObject_t CallConstructor(TCppMethod_t method, TCppType_t type, size_t nargs, void *args)
intptr_t TCppMethod_t
Definition cpp_cppyy.h:22
RPY_EXPORTED void * AllocateFunctionArgs(size_t nargs)
RPY_EXPORTED bool IsTemplate(const std::string &template_name)
RPY_EXPORTED TCppIndex_t GetMethodReqArgs(TCppMethod_t)
RPY_EXPORTED bool IsEnum(const std::string &type_name)
RPY_EXPORTED TCppObject_t Construct(TCppType_t type)
RPY_EXPORTED std::vector< TCppIndex_t > GetMethodIndicesFromName(TCppScope_t scope, const std::string &name)
RPY_EXPORTED bool ExistsMethodTemplate(TCppScope_t scope, const std::string &name)
RPY_EXPORTED std::string GetMethodName(TCppMethod_t)
RPY_EXPORTED bool IsConstData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED void AddSmartPtrType(const std::string &)
RPY_EXPORTED void CallDestructor(TCppType_t type, TCppObject_t self)
RPY_EXPORTED TCppScope_t gGlobalScope
Definition cpp_cppyy.h:51
RPY_EXPORTED bool IsProtectedData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED std::string GetMethodSignature(TCppMethod_t, bool show_formalargs, TCppIndex_t maxargs=(TCppIndex_t) -1)
RPY_EXPORTED int GetDimensionSize(TCppScope_t scope, TCppIndex_t idata, int dimension)
RPY_EXPORTED bool IsSubtype(TCppType_t derived, TCppType_t base)
RPY_EXPORTED TCppMethod_t GetMethodTemplate(TCppScope_t scope, const std::string &name, const std::string &proto)
void * TCppObject_t
Definition cpp_cppyy.h:21
RPY_EXPORTED bool IsConstructor(TCppMethod_t method)
RPY_EXPORTED bool GetSmartPtrInfo(const std::string &, TCppType_t *raw, TCppMethod_t *deref)
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 size_t GetFunctionArgTypeoffset()
RPY_EXPORTED TCppObject_t Allocate(TCppType_t type)
RPY_EXPORTED void Destruct(TCppType_t type, TCppObject_t instance)
RPY_EXPORTED std::string ResolveName(const std::string &cppitem_name)
TCppScope_t TCppType_t
Definition cpp_cppyy.h:19
RPY_EXPORTED std::string ResolveEnum(const std::string &enum_type)
RPY_EXPORTED long long GetEnumDataValue(TCppEnum_t, TCppIndex_t idata)
RPY_EXPORTED TCppIndex_t GetMethodNumArgs(TCppMethod_t)
RPY_EXPORTED TCppType_t GetActualClass(TCppType_t klass, TCppObject_t obj)
RPY_EXPORTED std::string GetBaseName(TCppType_t type, TCppIndex_t ibase)
RPY_EXPORTED bool IsNamespace(TCppScope_t scope)
void * TCppEnum_t
Definition cpp_cppyy.h:20
RPY_EXPORTED float CallF(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetScopedFinalName(TCppType_t type)
RPY_EXPORTED void Deallocate(TCppType_t type, TCppObject_t instance)
RPY_EXPORTED bool IsPublicData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED std::string GetMethodArgType(TCppMethod_t, TCppIndex_t iarg)
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 GetEnumDataName(TCppEnum_t, TCppIndex_t idata)
RPY_EXPORTED void GetAllCppNames(TCppScope_t scope, std::set< std::string > &cppnames)
RPY_EXPORTED bool IsComplete(const std::string &type_name)
RPY_EXPORTED bool IsBuiltin(const std::string &type_name)
RPY_EXPORTED bool IsStaticMethod(TCppMethod_t method)
RPY_EXPORTED TCppIndex_t GetDatamemberIndex(TCppScope_t scope, const std::string &name)
RPY_EXPORTED void CallV(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED bool IsStaticData(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED std::string GetDatamemberType(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED TCppMethod_t GetMethod(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED bool IsDestructor(TCppMethod_t method)
RPY_EXPORTED bool IsSmartPtr(TCppType_t type)
RPY_EXPORTED LongDouble_t CallLD(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED std::string GetTemplatedMethodName(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED size_t GetFunctionArgSizeof()
RPY_EXPORTED TCppScope_t GetScope(const std::string &scope_name)
RPY_EXPORTED bool HasVirtualDestructor(TCppType_t type)
RPY_EXPORTED bool IsConstMethod(TCppMethod_t)
RPY_EXPORTED bool HasComplexHierarchy(TCppType_t type)
RPY_EXPORTED std::vector< TCppScope_t > GetUsingNamespaces(TCppScope_t)
size_t TCppScope_t
Definition cpp_cppyy.h:18
RPY_EXPORTED bool IsTemplatedConstructor(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED TCppIndex_t GetNumDatamembers(TCppScope_t scope)
RPY_EXPORTED TCppIndex_t GetGlobalOperator(TCppType_t scope, const std::string &lc, const std::string &rc, const std::string &op)
RPY_EXPORTED TCppFuncAddr_t GetFunctionAddress(TCppMethod_t method, bool check_enabled=true)
RPY_EXPORTED TCppIndex_t GetNumEnumData(TCppEnum_t)
RPY_EXPORTED TCppIndex_t GetNumBases(TCppType_t type)
RPY_EXPORTED TCppIndex_t GetNumBasesLongestBranch(TCppType_t type)
Retrieve number of base classes in the longest branch of the inheritance tree.
RPY_EXPORTED TCppIndex_t GetNumTemplatedMethods(TCppScope_t scope)
RPY_EXPORTED std::string GetMethodPrototype(TCppScope_t scope, TCppMethod_t, bool show_formalargs)
RPY_EXPORTED std::string GetMethodResultType(TCppMethod_t)
RPY_EXPORTED std::string GetFinalName(TCppType_t type)
RPY_EXPORTED char * CallS(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args, size_t *length)
RPY_EXPORTED std::string GetMethodArgDefault(TCppMethod_t, TCppIndex_t iarg)
RPY_EXPORTED bool IsMethodTemplate(TCppScope_t scope, TCppIndex_t imeth)
RPY_EXPORTED std::string GetDatamemberName(TCppScope_t scope, TCppIndex_t idata)
RPY_EXPORTED bool IsPublicMethod(TCppMethod_t method)
RPY_EXPORTED intptr_t GetDatamemberOffset(TCppScope_t scope, TCppIndex_t idata)
void * TCppFuncAddr_t
Definition cpp_cppyy.h:25
RPY_EXPORTED std::string GetMethodFullName(TCppMethod_t)
RPY_EXPORTED short CallH(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED bool IsProtectedMethod(TCppMethod_t method)
RPY_EXPORTED bool Compile(const std::string &code)
RPY_EXPORTED Long64_t CallLL(TCppMethod_t method, TCppObject_t self, size_t nargs, void *args)
RPY_EXPORTED TCppEnum_t GetEnum(TCppScope_t scope, const std::string &enum_name)
void(* DelFunc_t)(void *)
Definition Rtypes.h:111
std::string ResolveTypedef(const char *tname, bool resolveAll=false)
std::string CleanType(const char *typeDesc, int mode=0, const char **tail=nullptr)
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
const char * fSigName
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345