Logo ROOT   6.08/07
Reference Guide
TROOT.cxx
Go to the documentation of this file.
1 // @(#)root/base:$Id$
2 // Author: Rene Brun 08/12/94
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TROOT
13 \ingroup Base
14 
15 ROOT top level object description.
16 
17 The TROOT object is the entry point to the ROOT system.
18 The single instance of TROOT is accessible via the global gROOT.
19 Using the gROOT pointer one has access to basically every object
20 created in a ROOT based program. The TROOT object is essentially a
21 container of several lists pointing to the main ROOT objects.
22 
23 The following lists are accessible from gROOT object:
24 
25 ~~~ {.cpp}
26  gROOT->GetListOfClasses
27  gROOT->GetListOfColors
28  gROOT->GetListOfTypes
29  gROOT->GetListOfGlobals
30  gROOT->GetListOfGlobalFunctions
31  gROOT->GetListOfFiles
32  gROOT->GetListOfMappedFiles
33  gROOT->GetListOfSockets
34  gROOT->GetListOfSecContexts
35  gROOT->GetListOfCanvases
36  gROOT->GetListOfStyles
37  gROOT->GetListOfFunctions
38  gROOT->GetListOfSpecials (for example graphical cuts)
39  gROOT->GetListOfGeometries
40  gROOT->GetListOfBrowsers
41  gROOT->GetListOfCleanups
42  gROOT->GetListOfMessageHandlers
43 ~~~
44 
45 The TROOT class provides also many useful services:
46  - Get pointer to an object in any of the lists above
47  - Time utilities TROOT::Time
48 
49 The ROOT object must be created as a static object. An example
50 of a main program creating an interactive version is shown below:
51 
52 ### Example of a main program
53 
54 ~~~ {.cpp}
55  #include "TRint.h"
56 
57  int main(int argc, char **argv)
58  {
59  TRint *theApp = new TRint("ROOT example", &argc, argv);
60 
61  // Init Intrinsics, build all windows, and enter event loop
62  theApp->Run();
63 
64  return(0);
65  }
66 ~~~
67 */
68 
69 #include "RConfig.h"
70 #include "RConfigure.h"
71 #include "RConfigOptions.h"
72 #include "RVersion.h"
73 #include "RGitCommit.h"
74 
75 #include <string>
76 #include <map>
77 #include <stdlib.h>
78 #ifdef WIN32
79 #include <io.h>
80 #include "Windows4Root.h"
81 #include <Psapi.h>
82 #define RTLD_DEFAULT ((void *)::GetModuleHandle(NULL))
83 #define dlsym(library, function_name) ::GetProcAddress((HMODULE)library, function_name)
84 #define dlopen(library_name, flags) ::LoadLibrary(library_name)
85 #define dlclose(library) ::FreeLibrary((HMODULE)library)
86 char *dlerror() {
87  static char Msg[1000];
88  FormatMessage(FORMAT_MESSAGE_FROM_SYSTEM, NULL, GetLastError(),
89  MAKELANGID(LANG_NEUTRAL, SUBLANG_DEFAULT), Msg,
90  sizeof(Msg), NULL);
91  return Msg;
92 }
93 #else
94 #include <dlfcn.h>
95 #endif
96 
97 #include "Riostream.h"
98 #include "Gtypes.h"
99 #include "TROOT.h"
100 #include "TClass.h"
101 #include "TClassEdit.h"
102 #include "TClassGenerator.h"
103 #include "TDataType.h"
104 #include "TDatime.h"
105 #include "TStyle.h"
106 #include "TObjectTable.h"
107 #include "TClassTable.h"
108 #include "TSystem.h"
109 #include "THashList.h"
110 #include "TObjArray.h"
111 #include "TEnv.h"
112 #include "TError.h"
113 #include "TColor.h"
114 #include "TGlobal.h"
115 #include "TFunction.h"
116 #include "TVirtualPad.h"
117 #include "TBrowser.h"
118 #include "TSystemDirectory.h"
119 #include "TApplication.h"
120 #include "TInterpreter.h"
121 #include "TGuiFactory.h"
122 #include "TMessageHandler.h"
123 #include "TFolder.h"
124 #include "TQObject.h"
125 #include "TProcessUUID.h"
126 #include "TPluginManager.h"
127 #include "TMap.h"
128 #include "TObjString.h"
129 #include "TVirtualMutex.h"
130 #include "TInterpreter.h"
131 #include "TListOfTypes.h"
132 #include "TListOfDataMembers.h"
133 #include "TListOfEnumsWithLock.h"
134 #include "TListOfFunctions.h"
136 #include "TFunctionTemplate.h"
137 #include "ThreadLocalStorage.h"
138 
139 #include <string>
140 namespace std {} using namespace std;
141 
142 #if defined(R__UNIX)
143 #if defined(R__HAS_COCOA)
144 #include "TMacOSXSystem.h"
145 #include "TUrl.h"
146 #else
147 #include "TUnixSystem.h"
148 #endif
149 #elif defined(R__WIN32)
150 #include "TWinNTSystem.h"
151 #endif
152 
153 extern "C" void R__SetZipMode(int);
154 
156 static void *gInterpreterLib = 0;
157 
158 // Mutex for protection of concurrent gROOT access
160 
161 // For accessing TThread::Tsd indirectly.
162 void **(*gThreadTsd)(void*,Int_t) = 0;
163 
164 //-------- Names of next three routines are a small homage to CMZ --------------
165 ////////////////////////////////////////////////////////////////////////////////
166 /// Return version id as an integer, i.e. "2.22/04" -> 22204.
167 
168 static Int_t IVERSQ()
169 {
170  Int_t maj, min, cycle;
171  sscanf(ROOT_RELEASE, "%d.%d/%d", &maj, &min, &cycle);
172  return 10000*maj + 100*min + cycle;
173 }
174 
175 ////////////////////////////////////////////////////////////////////////////////
176 /// Return built date as integer, i.e. "Apr 28 2000" -> 20000428.
177 
178 static Int_t IDATQQ(const char *date)
179 {
180  static const char *months[] = {"Jan","Feb","Mar","Apr","May",
181  "Jun","Jul","Aug","Sep","Oct",
182  "Nov","Dec"};
183 
184  char sm[12];
185  Int_t yy, mm=0, dd;
186  sscanf(date, "%s %d %d", sm, &dd, &yy);
187  for (int i = 0; i < 12; i++)
188  if (!strncmp(sm, months[i], 3)) {
189  mm = i+1;
190  break;
191  }
192  return 10000*yy + 100*mm + dd;
193 }
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 /// Return built time as integer (with min precision), i.e.
197 /// "17:32:37" -> 1732.
198 
199 static Int_t ITIMQQ(const char *time)
200 {
201  Int_t hh, mm, ss;
202  sscanf(time, "%d:%d:%d", &hh, &mm, &ss);
203  return 100*hh + mm;
204 }
205 
206 ////////////////////////////////////////////////////////////////////////////////
207 /// Clean up at program termination before global objects go out of scope.
208 
209 static void CleanUpROOTAtExit()
210 {
211  if (gROOT) {
212  R__LOCKGUARD(gROOTMutex);
213 
214  if (gROOT->GetListOfFiles())
215  gROOT->GetListOfFiles()->Delete("slow");
216  if (gROOT->GetListOfSockets())
217  gROOT->GetListOfSockets()->Delete();
218  if (gROOT->GetListOfMappedFiles())
219  gROOT->GetListOfMappedFiles()->Delete("slow");
220  if (gROOT->GetListOfClosedObjects())
221  gROOT->GetListOfClosedObjects()->Delete("slow");
222  }
223 }
224 
225 ////////////////////////////////////////////////////////////////////////////////
226 /// A module and its headers. Intentionally not a copy:
227 /// If these strings end up in this struct they are
228 /// long lived by definition because they get passed in
229 /// before initialization of TCling.
230 
231 namespace {
232  struct ModuleHeaderInfo_t {
233  ModuleHeaderInfo_t(const char* moduleName,
234  const char** headers,
235  const char** includePaths,
236  const char* payloadCode,
237  const char* fwdDeclCode,
238  void (*triggerFunc)(),
239  const TROOT::FwdDeclArgsToKeepCollection_t& fwdDeclsArgToSkip,
240  const char** classesHeaders):
241  fModuleName(moduleName),
242  fHeaders(headers),
243  fPayloadCode(payloadCode),
244  fFwdDeclCode(fwdDeclCode),
245  fIncludePaths(includePaths),
246  fTriggerFunc(triggerFunc),
247  fClassesHeaders(classesHeaders),
248  fFwdNargsToKeepColl(fwdDeclsArgToSkip){}
249 
250  const char* fModuleName; // module name
251  const char** fHeaders; // 0-terminated array of header files
252  const char* fPayloadCode; // Additional code to be given to cling at library load
253  const char* fFwdDeclCode; // Additional code to let cling know about selected classes and functions
254  const char** fIncludePaths; // 0-terminated array of header files
255  void (*fTriggerFunc)(); // Pointer to the dict initialization used to find the library name
256  const char** fClassesHeaders; // 0-terminated list of classes and related header files
257  const TROOT::FwdDeclArgsToKeepCollection_t fFwdNargsToKeepColl; // Collection of
258  // pairs of template fwd decls and number of
259  };
260 
261  std::vector<ModuleHeaderInfo_t>& GetModuleHeaderInfoBuffer() {
262  static std::vector<ModuleHeaderInfo_t> moduleHeaderInfoBuffer;
263  return moduleHeaderInfoBuffer;
264  }
265 }
266 
270 
271 static void at_exit_of_TROOT() {
274 }
275 
276 // This local static object initializes the ROOT system
277 namespace ROOT {
278 namespace Internal {
279  class TROOTAllocator {
280  // Simple wrapper to separate, time-wise, the call to the
281  // TROOT destructor and the actual free-ing of the memory.
282  //
283  // Since the interpreter implementation (currently TCling) is
284  // loaded via dlopen by libCore, the destruction of its global
285  // variable (i.e. in particular clang's) is scheduled before
286  // those in libCore so we need to schedule the call to the TROOT
287  // destructor before that *but* we want to make sure the memory
288  // stay around until libCore itself is unloaded so that code
289  // using gROOT can 'properly' check for validity.
290  //
291  // The order of loading for is:
292  // libCore.so
293  // libRint.so
294  // ... anything other library hard linked to the executable ...
295  // ... for example libEvent
296  // libCling.so
297  // ... other libraries like libTree for example ....
298  // and the destruction order is (of course) the reverse.
299  // By default the unloading of the dictionary, does use
300  // the service of the interpreter ... which of course
301  // fails if libCling is already unloaded by that information
302  // has not been registered per se.
303  //
304  // To solve this problem, we now schedule the destruction
305  // of the TROOT object to happen _just_ before the
306  // unloading/destruction of libCling so that we can
307  // maximize the amount of clean-up we can do correctly
308  // and we can still allocate the TROOT object's memory
309  // statically.
310  //
311  char fHolder[sizeof(TROOT)];
312  public:
313  TROOTAllocator() {
314  new(&(fHolder[0])) TROOT("root", "The ROOT of EVERYTHING");
315  }
316 
317  ~TROOTAllocator() {
318  if (gROOTLocal) {
319  gROOTLocal->~TROOT();
320  }
321  }
322  };
323 
324  // The global gROOT is defined to be a function (ROOT::GetROOT())
325  // which itself is dereferencing a function pointer.
326 
327  // Initially this function pointer's value is & GetROOT1 whose role is to
328  // create and initialize the TROOT object itself.
329  // At the very end of the TROOT constructor the value of the function pointer
330  // is switch to & GetROOT2 whose role is to initialize the interpreter.
331 
332  // This mechanism was primarily intended to fix the issues with order in which
333  // global TROOT and LLVM globals are initialized. TROOT was initializing
334  // Cling, but Cling could not be used yet due to LLVM globals not being
335  // initialized yet. The solution is to delay initializing the interpreter in
336  // TROOT till after main() when all LLVM globals are initialized.
337 
338  // Technically, the mechanism used actually delay the interpreter
339  // initialization until the first use of gROOT *after* the end of the
340  // TROOT constructor.
341 
342  // So to delay until after the start of main, we also made sure that none
343  // of the ROOT code (mostly the dictionary code) used during library loading
344  // is using gROOT (directly or indirectly).
345 
346  // In practice, the initialization of the interpreter is now delayed until
347  // the first use gROOT (or gInterpreter) after the start of main (but user
348  // could easily break this by using gROOT in their library initialization
349  // code).
350 
351  extern TROOT *gROOTLocal;
352 
354  if (gROOTLocal)
355  return gROOTLocal;
356  static TROOTAllocator alloc;
357  return gROOTLocal;
358  }
359 
361  static Bool_t initInterpreter = kFALSE;
362  if (!initInterpreter) {
363  initInterpreter = kTRUE;
364  gROOTLocal->InitInterpreter();
365  // Load and init threads library
366  gROOTLocal->InitThreads();
367  }
368  return gROOTLocal;
369  }
370  typedef TROOT *(*GetROOTFun_t)();
371 
373 
374  static Func_t GetSymInLibThread(const char *funcname)
375  {
376  const static bool loadSuccess = dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")? false : 0 <= gSystem->Load("libThread");
377  if (loadSuccess) {
378  if (auto sym = gSystem->DynFindSymbol(nullptr, funcname)) {
379  return sym;
380  } else {
381  Error("GetSymInLibThread", "Cannot get symbol %s.", funcname);
382  }
383  }
384  return nullptr;
385  }
386 
387  //////////////////////////////////////////////////////////////////////////////
388  /// Globally enables the parallel branch processing, which is a case of
389  /// implicit multi-threading (IMT) in ROOT, activating the required locks.
390  /// This IMT use case, implemented in TTree::GetEntry, spawns a task for
391  /// each branch of the tree. Therefore, a task takes care of the reading,
392  /// decompression and deserialisation of a given branch.
394  {
395 #ifdef R__USE_IMT
396  if (!IsImplicitMTEnabled())
398  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_EnableParBranchProcessing");
399  if (sym)
400  sym();
401 #else
402  ::Warning("EnableParBranchProcessing", "Cannot enable parallel branch processing, please build ROOT with -Dimt=ON");
403 #endif
404  }
405 
406  //////////////////////////////////////////////////////////////////////////////
407  /// Globally disables the IMT use case of parallel branch processing,
408  /// deactivating the corresponding locks.
410  {
411 #ifdef R__USE_IMT
412  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_DisableParBranchProcessing");
413  if (sym)
414  sym();
415 #else
416  ::Warning("DisableParBranchProcessing", "Cannot disable parallel branch processing, please build ROOT with -Dimt=ON");
417 #endif
418  }
419 
420  //////////////////////////////////////////////////////////////////////////////
421  /// Returns true if parallel branch processing is enabled.
423  {
424 #ifdef R__USE_IMT
425  static Bool_t (*sym)() = (Bool_t(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_IsParBranchProcessingEnabled");
426  if (sym)
427  return sym();
428  else
429  return kFALSE;
430 #else
431  return kFALSE;
432 #endif
433  }
434 
435  ////////////////////////////////////////////////////////////////////////////////
436  /// Globally enables the parallel tree processing, which is a case of
437  /// implicit multi-threading in ROOT, activating the required locks.
438  /// This IMT use case, implemented in TTreeProcessor::Process, receives a user
439  /// function and applies it to subranges of the tree, which correspond to its
440  /// clusters. Hence, for every cluster, a task is spawned to potentially
441  /// process it in parallel with the other clusters.
443  {
444 #ifdef R__USE_IMT
445  if (!IsImplicitMTEnabled())
447  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_EnableParTreeProcessing");
448  if (sym)
449  sym();
450 #else
451  ::Warning("EnableParTreeProcessing", "Cannot enable parallel tree processing, please build ROOT with -Dimt=ON");
452 #endif
453  }
454 
455  //////////////////////////////////////////////////////////////////////////////
456  /// Globally disables the IMT use case of parallel branch processing,
457  /// deactivating the corresponding locks.
459  {
460 #ifdef R__USE_IMT
461  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_DisableParTreeProcessing");
462  if (sym)
463  sym();
464 #else
465  ::Warning("DisableParTreeProcessing", "Cannot disable parallel tree processing, please build ROOT with -Dimt=ON");
466 #endif
467  }
468 
469  ////////////////////////////////////////////////////////////////////////////////
470  /// Returns true if parallel tree processing is enabled.
472  {
473 #ifdef R__USE_IMT
474  static Bool_t (*sym)() = (Bool_t(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_IsParTreeProcessingEnabled");
475  if (sym)
476  return sym();
477  else
478  return kFALSE;
479 #else
480  return kFALSE;
481 #endif
482  }
483 
484 } // end of Internal sub namespace
485 // back to ROOT namespace
486 
488  return (*Internal::gGetROOT)();
489  }
490 
491  TString &GetMacroPath() {
492  static TString macroPath;
493  return macroPath;
494  }
495 
496  ////////////////////////////////////////////////////////////////////////////////
497  /// Enables the global mutex to make ROOT thread safe/aware.
499  {
500  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TThread_Initialize");
501  if (sym)
502  sym();
503  }
504 
505  ////////////////////////////////////////////////////////////////////////////////
506  /// Globally enables the implicit multi-threading in ROOT, activating the
507  /// parallel execution of those methods in ROOT that provide an internal
508  /// parallelisation.
509  /// The 'numthreads' parameter allows to control the number of threads to
510  /// be used by the implicit multi-threading. However, this parameter is just
511  /// a hint for ROOT, which will try to satisfy the request if the execution
512  /// scenario allows it. For example, if ROOT is configured to use an external
513  /// scheduler, setting a value for 'numthreads' might not have any effect.
514  /// @param[in] numthreads Number of threads to use. If not specified or
515  /// set to zero, the number of threads is automatically
516  /// decided by the implementation. Any other value is
517  /// used as a hint.
518  void EnableImplicitMT(UInt_t numthreads)
519  {
520 #ifdef R__USE_IMT
522  static void (*sym)(UInt_t) = (void(*)(UInt_t))Internal::GetSymInLibThread("ROOT_TImplicitMT_EnableImplicitMT");
523  if (sym)
524  sym(numthreads);
525 #else
526  ::Warning("EnableImplicitMT", "Cannot enable implicit multi-threading with %d threads, please build ROOT with -Dimt=ON", numthreads);
527 #endif
528  }
529 
530  ////////////////////////////////////////////////////////////////////////////////
531  /// Disables the implicit multi-threading in ROOT.
533  {
534 #ifdef R__USE_IMT
535  static void (*sym)() = (void(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_DisableImplicitMT");
536  if (sym)
537  sym();
538 #else
539  ::Warning("DisableImplicitMT", "Cannot disable implicit multi-threading, please build ROOT with -Dimt=ON");
540 #endif
541  }
542 
543  ////////////////////////////////////////////////////////////////////////////////
544  /// Returns true if the implicit multi-threading in ROOT is enabled.
546  {
547 #ifdef R__USE_IMT
548  static Bool_t (*sym)() = (Bool_t(*)())Internal::GetSymInLibThread("ROOT_TImplicitMT_IsImplicitMTEnabled");
549  if (sym)
550  return sym();
551  else
552  return kFALSE;
553 #else
554  return kFALSE;
555 #endif
556  }
557 
558 }
559 
561 
562 // Global debug flag (set to > 0 to get debug output).
563 // Can be set either via the interpreter (gDebug is exported to CINT),
564 // via the rootrc resource "Root.Debug", via the shell environment variable
565 // ROOTDEBUG, or via the debugger.
567 
568 
570 
571 ////////////////////////////////////////////////////////////////////////////////
572 /// Default ctor.
573 
575  fLineIsProcessing(0), fVersion(0), fVersionInt(0), fVersionCode(0),
576  fVersionDate(0), fVersionTime(0), fBuiltDate(0), fBuiltTime(0),
577  fTimer(0), fApplication(0), fInterpreter(0), fBatch(kTRUE), fEditHistograms(kTRUE),
578  fFromPopUp(kTRUE),fMustClean(kTRUE),fReadingObject(kFALSE),fForceStyle(kFALSE),
579  fInterrupt(kFALSE),fEscape(kFALSE),fExecutingMacro(kFALSE),fEditorMode(0),
580  fPrimitive(0),fSelectPad(0),fClasses(0),fTypes(0),fGlobals(0),fGlobalFunctions(0),
581  fClosedObjects(0),fFiles(0),fMappedFiles(0),fSockets(0),fCanvases(0),fStyles(0),fFunctions(0),
582  fTasks(0),fColors(0),fGeometries(0),fBrowsers(0),fSpecials(0),fCleanups(0),
583  fMessageHandlers(0),fStreamerInfo(0),fClassGenerators(0),fSecContexts(0),
584  fProofs(0),fClipboard(0),fDataSets(0),fUUIDs(0),fRootFolder(0),fBrowsables(0),
585  fPluginManager(0)
586 {
587 }
588 
589 ////////////////////////////////////////////////////////////////////////////////
590 /// Initialize the ROOT system. The creation of the TROOT object initializes
591 /// the ROOT system. It must be the first ROOT related action that is
592 /// performed by a program. The TROOT object must be created on the stack
593 /// (can not be called via new since "operator new" is protected). The
594 /// TROOT object is either created as a global object (outside the main()
595 /// program), or it is one of the first objects created in main().
596 /// Make sure that the TROOT object stays in scope for as long as ROOT
597 /// related actions are performed. TROOT is a so called singleton so
598 /// only one instance of it can be created. The single TROOT object can
599 /// always be accessed via the global pointer gROOT.
600 /// The name and title arguments can be used to identify the running
601 /// application. The initfunc argument can contain an array of
602 /// function pointers (last element must be 0). These functions are
603 /// executed at the end of the constructor. This way one can easily
604 /// extend the ROOT system without adding permanent dependencies
605 /// (e.g. the graphics system is initialized via such a function).
606 
607 TROOT::TROOT(const char *name, const char *title, VoidFuncPtr_t *initfunc)
608  : TDirectory(), fLineIsProcessing(0), fVersion(0), fVersionInt(0), fVersionCode(0),
609  fVersionDate(0), fVersionTime(0), fBuiltDate(0), fBuiltTime(0),
610  fTimer(0), fApplication(0), fInterpreter(0), fBatch(kTRUE), fEditHistograms(kTRUE),
611  fFromPopUp(kTRUE),fMustClean(kTRUE),fReadingObject(kFALSE),fForceStyle(kFALSE),
612  fInterrupt(kFALSE),fEscape(kFALSE),fExecutingMacro(kFALSE),fEditorMode(0),
613  fPrimitive(0),fSelectPad(0),fClasses(0),fTypes(0),fGlobals(0),fGlobalFunctions(0),
614  fClosedObjects(0),fFiles(0),fMappedFiles(0),fSockets(0),fCanvases(0),fStyles(0),fFunctions(0),
615  fTasks(0),fColors(0),fGeometries(0),fBrowsers(0),fSpecials(0),fCleanups(0),
616  fMessageHandlers(0),fStreamerInfo(0),fClassGenerators(0),fSecContexts(0),
617  fProofs(0),fClipboard(0),fDataSets(0),fUUIDs(0),fRootFolder(0),fBrowsables(0),
618  fPluginManager(0)
619 {
620  if (fgRootInit || ROOT::Internal::gROOTLocal) {
621  //Warning("TROOT", "only one instance of TROOT allowed");
622  return;
623  }
624 
625  R__LOCKGUARD2(gROOTMutex);
626 
627  ROOT::Internal::gROOTLocal = this;
628  gDirectory = 0;
629 
630  // initialize gClassTable is not already done
631  if (!gClassTable) new TClassTable;
632 
633  SetName(name);
634  SetTitle(title);
635 
636  // will be used by global "operator delete" so make sure it is set
637  // before anything is deleted
638  fMappedFiles = 0;
639 
640  // create already here, but only initialize it after gEnv has been created
642 
643  // Initialize Operating System interface
644  InitSystem();
645 
647 
648  // Initialize interface to CINT C++ interpreter
649  fVersionInt = 0; // check in TROOT dtor in case TCling fails
650  fClasses = 0; // might be checked via TCling ctor
651  fEnums = 0;
652 
653  fConfigOptions = R__CONFIGUREOPTIONS;
654  fConfigFeatures = R__CONFIGUREFEATURES;
657  fVersionInt = IVERSQ();
660  fBuiltDate = IDATQQ(__DATE__);
661  fBuiltTime = ITIMQQ(__TIME__);
662 
663  ReadGitInfo();
664 
665  fClasses = new THashTable(800,3);
666  //fIdMap = new IdMap_t;
667  fStreamerInfo = new TObjArray(100);
668  fClassGenerators = new TList;
669 
670  // usedToIdentifyRootClingByDlSym is available when TROOT is part of
671  // rootcling.
672  if (!dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")) {
673  // initialize plugin manager early
675 #if defined(R__MACOSX) && (TARGET_OS_IPHONE || TARGET_IPHONE_SIMULATOR)
676  if (TARGET_OS_IPHONE | TARGET_IPHONE_SIMULATOR) {
677  TEnv plugins(".plugins-ios");
679  }
680 #endif
681  }
682 
683  TSystemDirectory *workdir = new TSystemDirectory("workdir", gSystem->WorkingDirectory());
684 
685  fTimer = 0;
686  fApplication = 0;
687  fColors = new TObjArray(1000); fColors->SetName("ListOfColors");
688  fTypes = 0;
689  fGlobals = 0;
690  fGlobalFunctions = 0;
691  // fList was created in TDirectory::Build but with different sizing.
692  delete fList;
693  fList = new THashList(1000,3);
694  fClosedObjects = new TList; fClosedObjects->SetName("ClosedFiles");
695  fFiles = new TList; fFiles->SetName("Files");
696  fMappedFiles = new TList; fMappedFiles->SetName("MappedFiles");
697  fSockets = new TList; fSockets->SetName("Sockets");
698  fCanvases = new TList; fCanvases->SetName("Canvases");
699  fStyles = new TList; fStyles->SetName("Styles");
700  fFunctions = new TList; fFunctions->SetName("Functions");
701  fTasks = new TList; fTasks->SetName("Tasks");
702  fGeometries = new TList; fGeometries->SetName("Geometries");
703  fBrowsers = new TList; fBrowsers->SetName("Browsers");
704  fSpecials = new TList; fSpecials->SetName("Specials");
705  fBrowsables = new TList; fBrowsables->SetName("Browsables");
706  fCleanups = new THashList; fCleanups->SetName("Cleanups");
707  fMessageHandlers = new TList; fMessageHandlers->SetName("MessageHandlers");
708  fSecContexts = new TList; fSecContexts->SetName("SecContexts");
709  fProofs = new TList; fProofs->SetName("Proofs");
710  fClipboard = new TList; fClipboard->SetName("Clipboard");
711  fDataSets = new TList; fDataSets->SetName("DataSets");
712  fTypes = new TListOfTypes;
713 
715  fUUIDs = new TProcessUUID();
716 
717  fRootFolder = new TFolder();
718  fRootFolder->SetName("root");
719  fRootFolder->SetTitle("root of all folders");
720  fRootFolder->AddFolder("Classes", "List of Active Classes",fClasses);
721  fRootFolder->AddFolder("Colors", "List of Active Colors",fColors);
722  fRootFolder->AddFolder("MapFiles", "List of MapFiles",fMappedFiles);
723  fRootFolder->AddFolder("Sockets", "List of Socket Connections",fSockets);
724  fRootFolder->AddFolder("Canvases", "List of Canvases",fCanvases);
725  fRootFolder->AddFolder("Styles", "List of Styles",fStyles);
726  fRootFolder->AddFolder("Functions", "List of Functions",fFunctions);
727  fRootFolder->AddFolder("Tasks", "List of Tasks",fTasks);
728  fRootFolder->AddFolder("Geometries","List of Geometries",fGeometries);
729  fRootFolder->AddFolder("Browsers", "List of Browsers",fBrowsers);
730  fRootFolder->AddFolder("Specials", "List of Special Objects",fSpecials);
731  fRootFolder->AddFolder("Handlers", "List of Message Handlers",fMessageHandlers);
732  fRootFolder->AddFolder("Cleanups", "List of RecursiveRemove Collections",fCleanups);
733  fRootFolder->AddFolder("StreamerInfo","List of Active StreamerInfo Classes",fStreamerInfo);
734  fRootFolder->AddFolder("SecContexts","List of Security Contexts",fSecContexts);
735  fRootFolder->AddFolder("PROOF Sessions", "List of PROOF sessions",fProofs);
736  fRootFolder->AddFolder("ROOT Memory","List of Objects in the gROOT Directory",fList);
737  fRootFolder->AddFolder("ROOT Files","List of Connected ROOT Files",fFiles);
738 
739  // by default, add the list of files, tasks, canvases and browsers in the Cleanups list
745 
748  fFromPopUp = kFALSE;
750  fInterrupt = kFALSE;
751  fEscape = kFALSE;
752  fMustClean = kTRUE;
753  fPrimitive = 0;
754  fSelectPad = 0;
755  fEditorMode = 0;
756  fDefCanvasName = "c1";
758  fLineIsProcessing = 1; // This prevents WIN32 "Windows" thread to pick ROOT objects with mouse
759  gDirectory = this;
760  gPad = 0;
761 
762  //set name of graphical cut class for the graphics editor
763  //cannot call SetCutClassName at this point because the TClass of TCutG
764  //is not yet build
765  fCutClassName = "TCutG";
766 
767  // Create a default MessageHandler
768  new TMessageHandler((TClass*)0);
769 
770  // Create some styles
771  gStyle = 0;
773  SetStyle(gEnv->GetValue("Canvas.Style", "Modern"));
774 
775  // Setup default (batch) graphics and GUI environment
778  gGXBatch = new TVirtualX("Batch", "ROOT Interface to batch graphics");
780 
781 #if defined(R__WIN32)
782  fBatch = kFALSE;
783 #elif defined(R__HAS_COCOA)
784  fBatch = kFALSE;
785 #else
786  if (gSystem->Getenv("DISPLAY"))
787  fBatch = kFALSE;
788  else
789  fBatch = kTRUE;
790 #endif
791 
792  int i = 0;
793  while (initfunc && initfunc[i]) {
794  (initfunc[i])();
795  fBatch = kFALSE; // put system in graphics mode (backward compatible)
796  i++;
797  }
798 
799  // Set initial/default list of browsable objects
800  fBrowsables->Add(fRootFolder, "root");
801  fBrowsables->Add(fProofs, "PROOF Sessions");
802  fBrowsables->Add(workdir, gSystem->WorkingDirectory());
803  fBrowsables->Add(fFiles, "ROOT Files");
804 
805  atexit(CleanUpROOTAtExit);
806 
808 }
809 
810 ////////////////////////////////////////////////////////////////////////////////
811 /// Clean up and free resources used by ROOT (files, network sockets,
812 /// shared memory segments, etc.).
813 
815 {
816  using namespace ROOT::Internal;
817 
818  if (gROOTLocal == this) {
819 
820  // If the interpreter has not yet been initialized, don't bother
821  gGetROOT = &GetROOT1;
822 
823  // Mark the object as invalid, so that we can veto some actions
824  // (like autoloading) while we are in the destructor.
826 
827  // Turn-off the global mutex to avoid recreating mutexes that have
828  // already been deleted during the destruction phase
829  gGlobalMutex = 0;
830 
831  // Return when error occurred in TCling, i.e. when setup file(s) are
832  // out of date
833  if (!fVersionInt) return;
834 
835  // ATTENTION!!! Order is important!
836 
837 #ifdef R__COMPLETE_MEM_TERMINATION
840  fSpecials->Delete(); SafeDelete(fSpecials); // delete special objects : PostScript, Minuit, Html
841 #endif
842  fClosedObjects->Delete("slow"); // and closed files
843  fFiles->Delete("slow"); // and files
845  fSecContexts->Delete("slow"); SafeDelete(fSecContexts); // and security contexts
846  fSockets->Delete(); SafeDelete(fSockets); // and sockets
847  fMappedFiles->Delete("slow"); // and mapped files
848  delete fUUIDs;
849  TProcessID::Cleanup(); // and list of ProcessIDs
850  TSeqCollection *tl = fMappedFiles; fMappedFiles = 0; delete tl;
851 
853 
854  fFunctions->Delete(); SafeDelete(fFunctions); // etc..
857 #ifdef R__COMPLETE_MEM_TERMINATION
859 #endif
862 
863 #ifdef R__COMPLETE_MEM_TERMINATION
868 #endif
869 
870  // Stop emitting signals
872 
874 
875 #ifdef R__COMPLETE_MEM_TERMINATION
881 
882  fCleanups->Clear();
884  delete gClassTable; gClassTable = 0;
885  delete gEnv; gEnv = 0;
886 
887  if (fTypes) fTypes->Delete();
889  if (fGlobals) fGlobals->Delete();
893  fEnums->Delete();
894  fClasses->Delete(); SafeDelete(fClasses); // TClass'es must be deleted last
895 #endif
896 
897  // Remove shared libraries produced by the TSystem::CompileMacro() call
899 
900  // Cleanup system class
901  delete gSystem;
902 
903  // ROOT-6022:
904  // if (gInterpreterLib) dlclose(gInterpreterLib);
905 #ifdef R__COMPLETE_MEM_TERMINATION
906  // On some 'newer' platform (Fedora Core 17+, Ubuntu 12), the
907  // initialization order is (by default?) is 'wrong' and so we can't
908  // delete the interpreter now .. because any of the static in the
909  // interpreter's library have already been deleted.
910  // On the link line, we must list the most dependent .o file
911  // and end with the least dependent (LLVM libraries), unfortunately,
912  // Fedora Core 17+ or Ubuntu 12 will also execute the initialization
913  // in the same order (hence doing libCore's before LLVM's and
914  // vice et versa for both the destructor. We worked around the
915  // initialization order by delay the TROOT creation until first use.
916  // We can not do the same for destruction as we have no way of knowing
917  // the last access ...
918  // So for now, let's avoid delete TCling except in the special build
919  // checking the completeness of the termination deletion.
920  gDestroyCling(fInterpreter);
921 #endif
922 
923 #ifdef R__COMPLETE_MEM_TERMINATION
925 #endif
926 
927  // Prints memory stats
929 
930  gROOTLocal = 0;
931  fgRootInit = kFALSE;
932  }
933 }
934 
935 ////////////////////////////////////////////////////////////////////////////////
936 /// Add a class to the list and map of classes.
937 /// This routine is deprecated, use TClass::AddClass directly.
938 
940 {
941  TClass::AddClass(cl);
942 }
943 
944 ////////////////////////////////////////////////////////////////////////////////
945 /// Add a class generator. This generator will be called by TClass::GetClass
946 /// in case its does not find a loaded rootcint dictionary to request the
947 /// creation of a TClass object.
948 
950 {
951  if (!generator) return;
952  fClassGenerators->Add(generator);
953 }
954 
955 ////////////////////////////////////////////////////////////////////////////////
956 /// Add browsable objects to TBrowser.
957 
959 {
960  TObject *obj;
961  TIter next(fBrowsables);
962 
963  while ((obj = (TObject *) next())) {
964  const char *opt = next.GetOption();
965  if (opt && strlen(opt))
966  b->Add(obj, opt);
967  else
968  b->Add(obj, obj->GetName());
969  }
970 }
971 
972 ////////////////////////////////////////////////////////////////////////////////
973 /// return class status bit kClassSaved for class cl
974 /// This function is called by the SavePrimitive functions writing
975 /// the C++ code for an object.
976 
978 {
979  if (cl == 0) return kFALSE;
980  if (cl->TestBit(TClass::kClassSaved)) return kTRUE;
982  return kFALSE;
983 }
984 
985 namespace {
986  static void R__ListSlowClose(TList *files)
987  {
988  // Routine to close a list of files using the 'slow' techniques
989  // that also for the deletion ot update the list itself.
990 
991  static TObject harmless;
992  TObjLink *cursor = files->FirstLink();
993  while (cursor) {
994  TDirectory *dir = static_cast<TDirectory*>( cursor->GetObject() );
995  if (dir) {
996  // In order for the iterator to stay valid, we must
997  // prevent the removal of the object (dir) from the list
998  // (which is done in TFile::Close). We can also can not
999  // just move to the next iterator since the Close might
1000  // also (indirectly) remove that file.
1001  // So we SetObject to a harmless value, so that 'dir'
1002  // is not seen as part of the list.
1003  // We will later, remove all the object (see files->Clear()
1004  cursor->SetObject(&harmless); // this must not be zero otherwise things go wrong.
1005  dir->Close();
1006  // Put it back
1007  cursor->SetObject(dir);
1008  }
1009  cursor = cursor->Next();
1010  };
1011  // Now were done, clear the list but do not delete the objects as
1012  // they have been moved to the list of closed objects and must be
1013  // deleted from there in order to avoid a double delete from a
1014  // use objects (on the interpreter stack).
1015  files->Clear("nodelete");
1016  }
1017 }
1018 
1019 ////////////////////////////////////////////////////////////////////////////////
1020 /// Close any files and sockets that gROOT knows about.
1021 /// This can be used to insures that the files and sockets are closed before any library is unloaded!
1022 
1024 {
1025  if (fFiles && fFiles->First()) {
1026  R__ListSlowClose(static_cast<TList*>(fFiles));
1027  }
1028  // and Close TROOT itself.
1029  Close();
1030  // Now sockets.
1031  if (fSockets && fSockets->First()) {
1032  if (0==fCleanups->FindObject(fSockets) ) {
1035  }
1036  CallFunc_t *socketCloser = gInterpreter->CallFunc_Factory();
1037  Long_t offset = 0;
1038  TClass *socketClass = TClass::GetClass("TSocket");
1039  gInterpreter->CallFunc_SetFuncProto(socketCloser, socketClass->GetClassInfo(), "Close", "", &offset);
1040  if (gInterpreter->CallFunc_IsValid(socketCloser)) {
1041  static TObject harmless;
1042  TObjLink *cursor = static_cast<TList*>(fSockets)->FirstLink();
1043  TList notclosed;
1044  while (cursor) {
1045  TObject *socket = cursor->GetObject();
1046  // In order for the iterator to stay valid, we must
1047  // prevent the removal of the object (dir) from the list
1048  // (which is done in TFile::Close). We can also can not
1049  // just move to the next iterator since the Close might
1050  // also (indirectly) remove that file.
1051  // So we SetObject to a harmless value, so that 'dir'
1052  // is not seen as part of the list.
1053  // We will later, remove all the object (see files->Clear()
1054  cursor->SetObject(&harmless); // this must not be zero otherwise things go wrong.
1055 
1056  if (socket->IsA()->InheritsFrom(socketClass)) {
1057  gInterpreter->CallFunc_Exec(socketCloser, ((char*)socket)+offset);
1058  // Put the object in the closed list for later deletion.
1059  socket->SetBit(kMustCleanup);
1060  fClosedObjects->AddLast(socket);
1061  } else {
1062  // Crap ... this is not a socket, likely Proof or something, let's try to find a Close
1063  Long_t other_offset;
1064  CallFunc_t *otherCloser = gInterpreter->CallFunc_Factory();
1065  gInterpreter->CallFunc_SetFuncProto(otherCloser, socket->IsA()->GetClassInfo(), "Close", "", &other_offset);
1066  if (gInterpreter->CallFunc_IsValid(otherCloser)) {
1067  gInterpreter->CallFunc_Exec(otherCloser, ((char*)socket)+other_offset);
1068  // Put the object in the closed list for later deletion.
1069  socket->SetBit(kMustCleanup);
1070  fClosedObjects->AddLast(socket);
1071  } else {
1072  notclosed.AddLast(socket);
1073  }
1074  gInterpreter->CallFunc_Delete(otherCloser);
1075  // Put it back
1076  cursor->SetObject(socket);
1077  }
1078  cursor = cursor->Next();
1079  }
1080  // Now were done, clear the list
1081  fSockets->Clear();
1082  // Read the one we did not close
1083  cursor = notclosed.FirstLink();
1084  while (cursor) {
1085  static_cast<TList*>(fSockets)->AddLast(cursor->GetObject());
1086  cursor = cursor->Next();
1087  }
1088  }
1089  gInterpreter->CallFunc_Delete(socketCloser);
1090  }
1091  if (fMappedFiles && fMappedFiles->First()) {
1092  R__ListSlowClose(static_cast<TList*>(fMappedFiles));
1093  }
1094 
1095 }
1096 
1097 ////////////////////////////////////////////////////////////////////////////////
1098 /// Execute the cleanups necessary at the end of the process, in particular
1099 /// those that must be executed before the library start being unloaded.
1100 
1102 {
1103  CloseFiles();
1104 
1105  if (gInterpreter) {
1106  gInterpreter->ResetGlobals();
1107  }
1108 
1109  // Now a set of simpler things to delete. See the same ordering in
1110  // TROOT::~TROOT
1111  fFunctions->Delete();
1112  fGeometries->Delete();
1113  fBrowsers->Delete();
1114  fCanvases->Delete();
1115  fColors->Delete();
1116  fStyles->Delete();
1117 }
1118 
1119 
1120 ////////////////////////////////////////////////////////////////////////////////
1121 /// Find an object in one Root folder
1122 
1124 {
1125  Error("FindObject","Not yet implemented");
1126  return 0;
1127 }
1128 
1129 ////////////////////////////////////////////////////////////////////////////////
1130 /// Returns address of a ROOT object if it exists
1131 ///
1132 /// If name contains at least one "/" the function calls FindObjectany
1133 /// else
1134 /// This function looks in the following order in the ROOT lists:
1135 /// - List of files
1136 /// - List of memory mapped files
1137 /// - List of functions
1138 /// - List of geometries
1139 /// - List of canvases
1140 /// - List of styles
1141 /// - List of specials
1142 /// - List of materials in current geometry
1143 /// - List of shapes in current geometry
1144 /// - List of matrices in current geometry
1145 /// - List of Nodes in current geometry
1146 /// - Current Directory in memory
1147 /// - Current Directory on file
1148 
1149 TObject *TROOT::FindObject(const char *name) const
1150 {
1151  if (name && strstr(name,"/")) return FindObjectAny(name);
1152 
1153  TObject *temp = 0;
1154 
1155  temp = fFiles->FindObject(name); if (temp) return temp;
1156  temp = fMappedFiles->FindObject(name); if (temp) return temp;
1157  {
1158  R__LOCKGUARD2(gROOTMutex);
1159  temp = fFunctions->FindObject(name);if (temp) return temp;
1160  }
1161  temp = fGeometries->FindObject(name); if (temp) return temp;
1162  temp = fCanvases->FindObject(name); if (temp) return temp;
1163  temp = fStyles->FindObject(name); if (temp) return temp;
1164  {
1165  R__LOCKGUARD2(gROOTMutex);
1166  temp = fSpecials->FindObject(name); if (temp) return temp;
1167  }
1168  TIter next(fGeometries);
1169  TObject *obj;
1170  while ((obj=next())) {
1171  temp = obj->FindObject(name); if (temp) return temp;
1172  }
1173  if (gDirectory) temp = gDirectory->Get(name);
1174  if (temp) return temp;
1175  if (gPad) {
1176  TVirtualPad *canvas = gPad->GetVirtCanvas();
1177  if (fCanvases->FindObject(canvas)) { //this check in case call from TCanvas ctor
1178  temp = canvas->FindObject(name);
1179  if (!temp && canvas != gPad) temp = gPad->FindObject(name);
1180  }
1181  }
1182  return temp;
1183 }
1184 
1185 ////////////////////////////////////////////////////////////////////////////////
1186 /// Returns address and folder of a ROOT object if it exists
1187 ///
1188 /// This function looks in the following order in the ROOT lists:
1189 /// - List of files
1190 /// - List of memory mapped files
1191 /// - List of functions
1192 /// - List of geometries
1193 /// - List of canvases
1194 /// - List of styles
1195 /// - List of specials
1196 /// - List of materials in current geometry
1197 /// - List of shapes in current geometry
1198 /// - List of matrices in current geometry
1199 /// - List of Nodes in current geometry
1200 /// - Current Directory in memory
1201 /// - Current Directory on file
1202 
1203 TObject *TROOT::FindSpecialObject(const char *name, void *&where)
1204 {
1205  TObject *temp = 0;
1206  where = 0;
1207 
1208  if (!temp) {
1209  temp = fFiles->FindObject(name);
1210  where = fFiles;
1211  }
1212  if (!temp) {
1213  temp = fMappedFiles->FindObject(name);
1214  where = fMappedFiles;
1215  }
1216  if (!temp) {
1217  R__LOCKGUARD2(gROOTMutex);
1218  temp = fFunctions->FindObject(name);
1219  where = fFunctions;
1220  }
1221  if (!temp) {
1222  temp = fCanvases->FindObject(name);
1223  where = fCanvases;
1224  }
1225  if (!temp) {
1226  temp = fStyles->FindObject(name);
1227  where = fStyles;
1228  }
1229  if (!temp) {
1230  temp = fSpecials->FindObject(name);
1231  where = fSpecials;
1232  }
1233  if (!temp) {
1234  TObject *glast = fGeometries->Last();
1235  if (glast) {where = glast; temp = glast->FindObject(name);}
1236  }
1237  if (!temp && gDirectory) {
1238  temp = gDirectory->Get(name);
1239  where = gDirectory;
1240  }
1241  if (!temp && gPad) {
1242  TVirtualPad *canvas = gPad->GetVirtCanvas();
1243  if (fCanvases->FindObject(canvas)) { //this check in case call from TCanvas ctor
1244  temp = canvas->FindObject(name);
1245  where = canvas;
1246  if (!temp && canvas != gPad) {
1247  temp = gPad->FindObject(name);
1248  where = gPad;
1249  }
1250  }
1251  }
1252  if (!temp) return 0;
1253  if (temp->TestBit(kNotDeleted)) return temp;
1254  return 0;
1255 }
1256 
1257 ////////////////////////////////////////////////////////////////////////////////
1258 /// Return a pointer to the first object with name starting at //root.
1259 /// This function scans the list of all folders.
1260 /// if no object found in folders, it scans the memory list of all files.
1261 
1262 TObject *TROOT::FindObjectAny(const char *name) const
1263 {
1264  TObject *obj = fRootFolder->FindObjectAny(name);
1265  if (obj) return obj;
1266  return gDirectory->FindObjectAnyFile(name);
1267 }
1268 
1269 ////////////////////////////////////////////////////////////////////////////////
1270 /// Scan the memory lists of all files for an object with name
1271 
1273 {
1274  R__LOCKGUARD(gROOTMutex);
1275  TDirectory *d;
1276  TIter next(GetListOfFiles());
1277  while ((d = (TDirectory*)next())) {
1278  // Call explicitly TDirectory::FindObject to restrict the search to the
1279  // already in memory object.
1280  TObject *obj = d->TDirectory::FindObject(name);
1281  if (obj) return obj;
1282  }
1283  return 0;
1284 }
1285 
1286 ////////////////////////////////////////////////////////////////////////////////
1287 /// Returns class name of a ROOT object including CINT globals.
1288 
1289 const char *TROOT::FindObjectClassName(const char *name) const
1290 {
1291  // Search first in the list of "standard" objects
1292  TObject *obj = FindObject(name);
1293  if (obj) return obj->ClassName();
1294 
1295  // Is it a global variable?
1296  TGlobal *g = GetGlobal(name);
1297  if (g) return g->GetTypeName();
1298 
1299  return 0;
1300 }
1301 
1302 ////////////////////////////////////////////////////////////////////////////////
1303 /// Return path name of obj somewhere in the //root/... path.
1304 /// The function returns the first occurence of the object in the list
1305 /// of folders. The returned string points to a static char array in TROOT.
1306 /// If this function is called in a loop or recursively, it is the
1307 /// user's responsibility to copy this string in their area.
1308 
1309 const char *TROOT::FindObjectPathName(const TObject *) const
1310 {
1311  Error("FindObjectPathName","Not yet implemented");
1312  return "??";
1313 }
1314 
1315 ////////////////////////////////////////////////////////////////////////////////
1316 /// return a TClass object corresponding to 'name' assuming it is an STL container.
1317 /// In particular we looking for possible alternative name (default template
1318 /// parameter, typedefs template arguments, typedefed name).
1319 
1320 TClass *TROOT::FindSTLClass(const char *name, Bool_t load, Bool_t silent) const
1321 {
1322  // Example of inputs are
1323  // vector<int> (*)
1324  // vector<Int_t>
1325  // vector<long long>
1326  // vector<Long_64_t> (*)
1327  // vector<int, allocator<int> >
1328  // vector<Int_t, allocator<int> >
1329  //
1330  // One of the possibly expensive operation is the resolving of the typedef
1331  // which can provoke the parsing of the header files (and/or the loading
1332  // of clang pcms information).
1333 
1335 
1336  // Remove std::, allocator, typedef, add Long64_t, etc. in just one call.
1337  std::string normalized;
1338  TClassEdit::GetNormalizedName(normalized, name);
1339 
1340  TClass *cl = 0;
1341  if (normalized != name) cl = TClass::GetClass(normalized.c_str(),load,silent);
1342 
1343  if (load && cl==0) {
1344  // Create an Emulated class for this container.
1345  cl = gInterpreter->GenerateTClass(normalized.c_str(), kTRUE, silent);
1346  }
1347 
1348  return cl;
1349 }
1350 
1351 ////////////////////////////////////////////////////////////////////////////////
1352 /// Return pointer to class with name. Obsolete, use TClass::GetClass directly
1353 
1354 TClass *TROOT::GetClass(const char *name, Bool_t load, Bool_t silent) const
1355 {
1356  return TClass::GetClass(name,load,silent);
1357 }
1358 
1359 
1360 ////////////////////////////////////////////////////////////////////////////////
1361 /// Return pointer to class from its name. Obsolete, use TClass::GetClass directly
1362 /// See TClass::GetClass
1363 
1364 TClass *TROOT::GetClass(const std::type_info& typeinfo, Bool_t load, Bool_t silent) const
1365 {
1366  return TClass::GetClass(typeinfo,load,silent);
1367 }
1368 
1369 ////////////////////////////////////////////////////////////////////////////////
1370 /// Return address of color with index color.
1371 
1373 {
1375  TObjArray *lcolors = (TObjArray*) GetListOfColors();
1376  if (!lcolors) return 0;
1377  if (color < 0 || color >= lcolors->GetSize()) return 0;
1378  TColor *col = (TColor*)lcolors->At(color);
1379  if (col && col->GetNumber() == color) return col;
1380  TIter next(lcolors);
1381  while ((col = (TColor *) next()))
1382  if (col->GetNumber() == color) return col;
1383 
1384  return 0;
1385 }
1386 
1387 ////////////////////////////////////////////////////////////////////////////////
1388 /// Return a default canvas.
1389 
1391 {
1392  return (TCanvas*)gROOT->ProcessLine("TCanvas::MakeDefCanvas();");
1393 }
1394 
1395 ////////////////////////////////////////////////////////////////////////////////
1396 /// Return pointer to type with name.
1397 
1398 TDataType *TROOT::GetType(const char *name, Bool_t /* load */) const
1399 {
1400  return (TDataType*)gROOT->GetListOfTypes()->FindObject(name);
1401 }
1402 
1403 ////////////////////////////////////////////////////////////////////////////////
1404 /// Return pointer to file with name.
1405 
1406 TFile *TROOT::GetFile(const char *name) const
1407 {
1408  R__LOCKGUARD(gROOTMutex);
1409  return (TFile*)GetListOfFiles()->FindObject(name);
1410 }
1411 
1412 ////////////////////////////////////////////////////////////////////////////////
1413 /// Return pointer to style with name
1414 
1415 TStyle *TROOT::GetStyle(const char *name) const
1416 {
1417  return (TStyle*)GetListOfStyles()->FindObject(name);
1418 }
1419 
1420 ////////////////////////////////////////////////////////////////////////////////
1421 /// Return pointer to function with name.
1422 
1423 TObject *TROOT::GetFunction(const char *name) const
1424 {
1425  if (name == 0 || name[0] == 0) {
1426  return 0;
1427  }
1428 
1429  {
1430  R__LOCKGUARD2(gROOTMutex);
1431  TObject *f1 = fFunctions->FindObject(name);
1432  if (f1) return f1;
1433  }
1434 
1435  gROOT->ProcessLine("TF1::InitStandardFunctions();");
1436 
1437  R__LOCKGUARD2(gROOTMutex);
1438  return fFunctions->FindObject(name);
1439 }
1440 
1441 ////////////////////////////////////////////////////////////////////////////////
1442 
1444 {
1445  if (!gInterpreter) return 0;
1446 
1448 
1449  return (TFunctionTemplate*)fFuncTemplate->FindObject(name);
1450 }
1451 
1452 ////////////////////////////////////////////////////////////////////////////////
1453 /// Return pointer to global variable by name. If load is true force
1454 /// reading of all currently defined globals from CINT (more expensive).
1455 
1456 TGlobal *TROOT::GetGlobal(const char *name, Bool_t load) const
1457 {
1458  return (TGlobal *)gROOT->GetListOfGlobals(load)->FindObject(name);
1459 }
1460 
1461 ////////////////////////////////////////////////////////////////////////////////
1462 /// Return pointer to global variable with address addr.
1463 
1464 TGlobal *TROOT::GetGlobal(const TObject *addr, Bool_t /* load */) const
1465 {
1466  if (addr == 0 || ((Long_t)addr) == -1) return 0;
1467 
1468  TInterpreter::DeclId_t decl = gInterpreter->GetDataMemberAtAddr(addr);
1469  if (decl) {
1470  TListOfDataMembers *globals = ((TListOfDataMembers*)(gROOT->GetListOfGlobals(kFALSE)));
1471  return (TGlobal*)globals->Get(decl);
1472  }
1473  // If we are actually looking for a global that is held by a global
1474  // pointer (for example gRandom), we need to find a pointer with the
1475  // correct value.
1476  decl = gInterpreter->GetDataMemberWithValue(addr);
1477  if (decl) {
1478  TListOfDataMembers *globals = ((TListOfDataMembers*)(gROOT->GetListOfGlobals(kFALSE)));
1479  return (TGlobal*)globals->Get(decl);
1480  }
1481  return 0;
1482 }
1483 
1484 ////////////////////////////////////////////////////////////////////////////////
1485 /// Internal routine returning, and creating if necessary, the list
1486 /// of global function.
1487 
1489 {
1491  return fGlobalFunctions;
1492 }
1493 
1494 ////////////////////////////////////////////////////////////////////////////////
1495 /// Return the collection of functions named "name".
1496 
1498 {
1499  return ((TListOfFunctions*)fFunctions)->GetListForObject(name);
1500 }
1501 
1502 ////////////////////////////////////////////////////////////////////////////////
1503 /// Return pointer to global function by name.
1504 /// If params != 0 it will also resolve overloading other it returns the first
1505 /// name match.
1506 /// If params == 0 and load is true force reading of all currently defined
1507 /// global functions from Cling.
1508 /// The param string must be of the form: "3189,\"aap\",1.3".
1509 
1510 TFunction *TROOT::GetGlobalFunction(const char *function, const char *params,
1511  Bool_t load)
1512 {
1513  if (!params) {
1514  R__LOCKGUARD2(gROOTMutex);
1515  return (TFunction *)GetListOfGlobalFunctions(load)->FindObject(function);
1516  } else {
1517  if (!fInterpreter)
1518  Fatal("GetGlobalFunction", "fInterpreter not initialized");
1519 
1520  R__LOCKGUARD2(gROOTMutex);
1521  TInterpreter::DeclId_t decl = gInterpreter->GetFunctionWithValues(0,
1522  function, params,
1523  false);
1524 
1525  if (!decl) return 0;
1526 
1527  TFunction *f = GetGlobalFunctions()->Get(decl);
1528  if (f) return f;
1529 
1530  Error("GetGlobalFunction",
1531  "\nDid not find matching TFunction <%s> with \"%s\".",
1532  function,params);
1533  return 0;
1534  }
1535 }
1536 
1537 ////////////////////////////////////////////////////////////////////////////////
1538 /// Return pointer to global function by name. If proto != 0
1539 /// it will also resolve overloading. If load is true force reading
1540 /// of all currently defined global functions from CINT (more expensive).
1541 /// The proto string must be of the form: "int, char*, float".
1542 
1544  const char *proto, Bool_t load)
1545 {
1546  if (!proto) {
1547  R__LOCKGUARD2(gROOTMutex);
1548  return (TFunction *)GetListOfGlobalFunctions(load)->FindObject(function);
1549  } else {
1550  if (!fInterpreter)
1551  Fatal("GetGlobalFunctionWithPrototype", "fInterpreter not initialized");
1552 
1553  R__LOCKGUARD2(gROOTMutex);
1554  TInterpreter::DeclId_t decl = gInterpreter->GetFunctionWithPrototype(0,
1555  function, proto);
1556 
1557  if (!decl) return 0;
1558 
1559  TFunction *f = GetGlobalFunctions()->Get(decl);
1560  if (f) return f;
1561 
1562  Error("GetGlobalFunctionWithPrototype",
1563  "\nDid not find matching TFunction <%s> with \"%s\".",
1564  function,proto);
1565  return 0;
1566  }
1567 }
1568 
1569 ////////////////////////////////////////////////////////////////////////////////
1570 /// Return pointer to Geometry with name
1571 
1572 TObject *TROOT::GetGeometry(const char *name) const
1573 {
1574  return GetListOfGeometries()->FindObject(name);
1575 }
1576 
1577 ////////////////////////////////////////////////////////////////////////////////
1578 
1580 {
1581  if(!fEnums.load()) {
1582  R__LOCKGUARD2(gROOTMutex);
1583  // Test again just in case, another thread did the work while we were
1584  // waiting.
1585  if (!fEnums.load()) fEnums = new TListOfEnumsWithLock(0);
1586  }
1587  if (load) {
1588  R__LOCKGUARD2(gROOTMutex);
1589  (*fEnums).Load(); // Refresh the list of enums.
1590  }
1591  return fEnums.load();
1592 }
1593 
1594 ////////////////////////////////////////////////////////////////////////////////
1595 
1597 {
1598  R__LOCKGUARD2(gROOTMutex);
1599  if(!fFuncTemplate) {
1601  }
1602  return fFuncTemplate;
1603 }
1604 
1605 ////////////////////////////////////////////////////////////////////////////////
1606 /// Return list containing the TGlobals currently defined.
1607 /// Since globals are created and deleted during execution of the
1608 /// program, we need to update the list of globals every time we
1609 /// execute this method. However, when calling this function in
1610 /// a (tight) loop where no interpreter symbols will be created
1611 /// you can set load=kFALSE (default).
1612 
1614 {
1615  if (!fGlobals) {
1616  // We add to the list the "funcky-fake" globals.
1617  fGlobals = new TListOfDataMembers(0);
1618  fGlobals->Add(new TGlobalMappedFunction("gROOT", "TROOT*",
1620  fGlobals->Add(new TGlobalMappedFunction("gPad", "TVirtualPad*",
1622  fGlobals->Add(new TGlobalMappedFunction("gInterpreter", "TInterpreter*",
1624  fGlobals->Add(new TGlobalMappedFunction("gVirtualX", "TTVirtualX*",
1626  fGlobals->Add(new TGlobalMappedFunction("gDirectory", "TDirectory*",
1630  }
1631 
1632  if (!fInterpreter)
1633  Fatal("GetListOfGlobals", "fInterpreter not initialized");
1634 
1635  if (load) fGlobals->Load();
1636 
1637  return fGlobals;
1638 }
1639 
1640 ////////////////////////////////////////////////////////////////////////////////
1641 /// Return list containing the TFunctions currently defined.
1642 /// Since functions are created and deleted during execution of the
1643 /// program, we need to update the list of functions every time we
1644 /// execute this method. However, when calling this function in
1645 /// a (tight) loop where no interpreter symbols will be created
1646 /// you can set load=kFALSE (default).
1647 
1649 {
1650  R__LOCKGUARD2(gROOTMutex);
1651 
1652  if (!fGlobalFunctions) {
1654  }
1655 
1656  if (!fInterpreter)
1657  Fatal("GetListOfGlobalFunctions", "fInterpreter not initialized");
1658 
1659  // A thread that calls with load==true and a thread that calls with load==false
1660  // will conflict here (the load==true will be updating the list while the
1661  // other is reading it). To solve the problem, we could use a read-write lock
1662  // inside the list itself.
1663  if (load) fGlobalFunctions->Load();
1664 
1665  return fGlobalFunctions;
1666 }
1667 
1668 ////////////////////////////////////////////////////////////////////////////////
1669 /// Return a dynamic list giving access to all TDataTypes (typedefs)
1670 /// currently defined.
1671 ///
1672 /// The list is populated on demand. Calling
1673 /// ~~~ {.cpp}
1674 /// gROOT->GetListOfTypes()->FindObject(nameoftype);
1675 /// ~~~
1676 /// will return the TDataType corresponding to 'nameoftype'. If the
1677 /// TDataType is not already in the list itself and the type does exist,
1678 /// a new TDataType will be created and added to the list.
1679 ///
1680 /// Calling
1681 /// ~~~ {.cpp}
1682 /// gROOT->GetListOfTypes()->ls(); // or Print()
1683 /// ~~~
1684 /// list only the typedefs that have been previously accessed through the
1685 /// list (plus the builtins types).
1686 
1688 {
1689  if (!fInterpreter)
1690  Fatal("GetListOfTypes", "fInterpreter not initialized");
1691 
1692  return fTypes;
1693 }
1694 
1695 
1696 ////////////////////////////////////////////////////////////////////////////////
1697 /// Execute command when system has been idle for idleTimeInSec seconds.
1698 
1699 void TROOT::Idle(UInt_t idleTimeInSec, const char *command)
1700 {
1701  if (!fApplication.load())
1703 
1704  if (idleTimeInSec <= 0)
1705  (*fApplication).RemoveIdleTimer();
1706  else
1707  (*fApplication).SetIdleTimer(idleTimeInSec, command);
1708 }
1709 
1710 ////////////////////////////////////////////////////////////////////////////////
1711 /// Check whether className is a known class, and only autoload
1712 /// if we can. Helper function for TROOT::IgnoreInclude().
1713 
1714 static TClass* R__GetClassIfKnown(const char* className)
1715 {
1716  // Check whether the class is available for auto-loading first:
1717  const char* libsToLoad = gInterpreter->GetClassSharedLibs(className);
1718  TClass* cla = 0;
1719  if (libsToLoad) {
1720  // trigger autoload, and only create TClass in this case.
1721  return TClass::GetClass(className);
1722  } else if (gROOT->GetListOfClasses()
1723  && (cla = (TClass*)gROOT->GetListOfClasses()->FindObject(className))) {
1724  // cla assigned in if statement
1725  } else if (gClassTable->FindObject(className)) {
1726  return TClass::GetClass(className);
1727  }
1728  return cla;
1729 }
1730 
1731 ////////////////////////////////////////////////////////////////////////////////
1732 /// Return 1 if the name of the given include file corresponds to a class that
1733 /// is known to ROOT, e.g. "TLorentzVector.h" versus TLorentzVector.
1734 
1735 Int_t TROOT::IgnoreInclude(const char *fname, const char * /*expandedfname*/)
1736 {
1737  if (fname == 0) return 0;
1738 
1739  TString stem(fname);
1740  // Remove extension if any, ignore files with extension not being .h*
1741  Int_t where = stem.Last('.');
1742  if (where != kNPOS) {
1743  if (stem.EndsWith(".so") || stem.EndsWith(".sl") ||
1744  stem.EndsWith(".dl") || stem.EndsWith(".a") ||
1745  stem.EndsWith(".dll", TString::kIgnoreCase))
1746  return 0;
1747  stem.Remove(where);
1748  }
1749 
1750  TString className = gSystem->BaseName(stem);
1751  TClass* cla = R__GetClassIfKnown(className);
1752  if (!cla) {
1753  // Try again with modifications to the file name:
1754  className = stem;
1755  className.ReplaceAll("/", "::");
1756  className.ReplaceAll("\\", "::");
1757  if (className.Contains(":::")) {
1758  // "C:\dir" becomes "C:::dir".
1759  // fname corresponds to whatever is stated after #include and
1760  // a full path name usually means that it's not a regular #include
1761  // but e.g. a ".L", so we can assume that this is not a header of
1762  // a class in a namespace (a global-namespace class would have been
1763  // detected already before).
1764  return 0;
1765  }
1766  cla = R__GetClassIfKnown(className);
1767  }
1768 
1769  if (!cla) {
1770  return 0;
1771  }
1772 
1773  // cla is valid, check wether it's actually in the header of the same name:
1774  if (cla->GetDeclFileLine() <= 0) return 0; // to a void an error with VisualC++
1775  TString decfile = gSystem->BaseName(cla->GetDeclFileName());
1776  if (decfile != gSystem->BaseName(fname)) {
1777  return 0;
1778  }
1779  return 1;
1780 }
1781 
1782 ////////////////////////////////////////////////////////////////////////////////
1783 /// Initialize operating system interface.
1784 
1786 {
1787  if (gSystem == 0) {
1788 #if defined(R__UNIX)
1789 #if defined(R__HAS_COCOA)
1790  gSystem = new TMacOSXSystem;
1791 #else
1792  gSystem = new TUnixSystem;
1793 #endif
1794 #elif defined(R__WIN32)
1795  gSystem = new TWinNTSystem;
1796 #else
1797  gSystem = new TSystem;
1798 #endif
1799 
1800  if (gSystem->Init())
1801  fprintf(stderr, "Fatal in <TROOT::InitSystem>: can't init operating system layer\n");
1802 
1803  if (!gSystem->HomeDirectory()) {
1804  fprintf(stderr, "Fatal in <TROOT::InitSystem>: HOME directory not set\n");
1805  fprintf(stderr, "Fix this by defining the HOME shell variable\n");
1806  }
1807 
1808  // read default files
1809  gEnv = new TEnv(".rootrc");
1810 
1811  gDebug = gEnv->GetValue("Root.Debug", 0);
1812 
1813  if (!gEnv->GetValue("Root.ErrorHandlers", 1))
1814  gSystem->ResetSignals();
1815 
1816  // by default the zipmode is 1 (see Bits.h)
1817  Int_t zipmode = gEnv->GetValue("Root.ZipMode", 1);
1818  if (zipmode != 1) R__SetZipMode(zipmode);
1819 
1820  const char *sdeb;
1821  if ((sdeb = gSystem->Getenv("ROOTDEBUG")))
1822  gDebug = atoi(sdeb);
1823 
1824  if (gDebug > 0 && isatty(2))
1825  fprintf(stderr, "Info in <TROOT::InitSystem>: running with gDebug = %d\n", gDebug);
1826 
1827  if (gEnv->GetValue("Root.MemStat", 0))
1829  int msize = gEnv->GetValue("Root.MemStat.size", -1);
1830  int mcnt = gEnv->GetValue("Root.MemStat.cnt", -1);
1831  if (msize != -1 || mcnt != -1)
1832  TStorage::EnableStatistics(msize, mcnt);
1833 
1834  fgMemCheck = gEnv->GetValue("Root.MemCheck", 0);
1835 
1836 #if defined(R__HAS_COCOA)
1837  // create and delete a dummy TUrl so that TObjectStat table does not contain
1838  // objects that are deleted after recording is turned-off (in next line),
1839  // like the TUrl::fgSpecialProtocols list entries which are created in the
1840  // TMacOSXSystem ctor.
1841  { TUrl dummy("/dummy"); }
1842 #endif
1843  TObject::SetObjectStat(gEnv->GetValue("Root.ObjectStat", 0));
1844  }
1845 }
1846 
1847 ////////////////////////////////////////////////////////////////////////////////
1848 /// Load and initialize thread library.
1849 
1851 {
1852  if (gEnv->GetValue("Root.UseThreads", 0) || gEnv->GetValue("Root.EnableThreadSafety", 0)) {
1854  }
1855 }
1856 
1857 ////////////////////////////////////////////////////////////////////////////////
1858 /// Initialize the interpreter. Should be called only after main(),
1859 /// to make sure LLVM/Clang is fully initialized.
1860 
1862 {
1863  // usedToIdentifyRootClingByDlSym is available when TROOT is part of
1864  // rootcling.
1865  if (!dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")
1866  && !dlsym(RTLD_DEFAULT, "usedToIdentifyStaticRoot")) {
1867  // Make sure no llvm symbols are visible before loading libCling. If they
1868  // exist libCling will use those and not ours, causing havoc in the
1869  // interpreter. Look for an extern "C" symbol to avoid mangling; look for a
1870  // symbol from llvm because clang builds on top, so users would likely
1871  // have also their own llvm symbols when providing their own clang.
1872  void *LLVMEnablePrettyStackTraceAddr = 0;
1873  // Can't use gSystem->DynFindSymbol() because that iterates over all *known*
1874  // libraries which is not the same!
1875  LLVMEnablePrettyStackTraceAddr = dlsym(RTLD_DEFAULT, "LLVMEnablePrettyStackTrace");
1876  if (LLVMEnablePrettyStackTraceAddr) {
1877  Error("InitInterpreter()", "LLVM SYMBOLS ARE EXPOSED TO CLING! "
1878  "This will cause problems; please hide them or dlopen() them "
1879  "after the call to TROOT::InitInterpreter()!");
1880  }
1881 
1882  const char *libcling = 0;
1883  char *libclingStorage = 0;
1884 #ifdef ROOTLIBDIR
1885  libcling = ROOTLIBDIR "/libCling."
1886 # ifdef R__WIN32
1887  "dll";
1888 # else
1889  "so";
1890 # endif
1891 #else
1892  libclingStorage = gSystem->DynamicPathName("libCling");
1893  libcling = libclingStorage;
1894 #endif
1895  gInterpreterLib = dlopen(libcling, RTLD_LAZY|RTLD_LOCAL);
1896  delete [] libclingStorage;
1897 
1898  if (!gInterpreterLib) {
1899  TString err = dlerror();
1900  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load library %s\n", err.Data());
1901  exit(1);
1902  }
1903  dlerror(); // reset error message
1904  } else {
1905  gInterpreterLib = RTLD_DEFAULT;
1906  }
1907  CreateInterpreter_t *CreateInterpreter = (CreateInterpreter_t*) dlsym(gInterpreterLib, "CreateInterpreter");
1908  if (!CreateInterpreter) {
1909  TString err = dlerror();
1910  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load symbol %s\n", err.Data());
1911  exit(1);
1912  }
1913  // Schedule the destruction of TROOT.
1914  atexit(at_exit_of_TROOT);
1915 
1916  gDestroyInterpreter = (DestroyInterpreter_t*) dlsym(gInterpreterLib, "DestroyInterpreter");
1917  if (!gDestroyInterpreter) {
1918  TString err = dlerror();
1919  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load symbol %s\n", err.Data());
1920  exit(1);
1921  }
1922 
1923  fInterpreter = CreateInterpreter(gInterpreterLib);
1924 
1927 
1928  fgRootInit = kTRUE;
1929 
1930  // Initialize all registered dictionaries.
1931  for (std::vector<ModuleHeaderInfo_t>::const_iterator
1932  li = GetModuleHeaderInfoBuffer().begin(),
1933  le = GetModuleHeaderInfoBuffer().end(); li != le; ++li) {
1934  // process buffered module registrations
1935  fInterpreter->RegisterModule(li->fModuleName,
1936  li->fHeaders,
1937  li->fIncludePaths,
1938  li->fPayloadCode,
1939  li->fFwdDeclCode,
1940  li->fTriggerFunc,
1941  li->fFwdNargsToKeepColl,
1942  li->fClassesHeaders,
1943  kTRUE /*lateRegistration*/);
1944  }
1945  GetModuleHeaderInfoBuffer().clear();
1946 
1948 
1949  // Read the rules before enabling the auto loading to not inadvertently
1950  // load the libraries for the classes concerned even-though the user is
1951  // *not* using them.
1952  TClass::ReadRules(); // Read the default customization rules ...
1953 
1954  // Enable autoloading
1956 }
1957 
1958 ////////////////////////////////////////////////////////////////////////////////
1959 /// Helper function used by TClass::GetClass().
1960 /// This function attempts to load the dictionary for 'classname'
1961 /// either from the TClassTable or from the list of generator.
1962 /// If silent is 'true', do not warn about missing dictionary for the class.
1963 /// (typically used for class that are used only for transient members)
1964 ///
1965 /// The 'requestedname' is expected to be already normalized.
1966 
1967 TClass *TROOT::LoadClass(const char *requestedname, Bool_t silent) const
1968 {
1969  return TClass::LoadClass(requestedname, silent);
1970 }
1971 
1972 ////////////////////////////////////////////////////////////////////////////////
1973 /// Check if class "classname" is known to the interpreter (in fact,
1974 /// this check is not needed anymore, so classname is ignored). If
1975 /// not it will load library "libname". If the library name does
1976 /// not start with "lib", "lib" will be prepended and a search will
1977 /// be made in the DynamicPath (see .rootrc). If not found a search
1978 /// will be made on libname (without "lib" prepended) and if not found
1979 /// a direct try of libname will be made (in case it contained an
1980 /// absolute path).
1981 /// If check is true it will only check if libname exists and is
1982 /// readable.
1983 /// Returns 0 on successful loading, -1 in case libname does not
1984 /// exist or in case of error and -2 in case of version mismatch.
1985 
1986 Int_t TROOT::LoadClass(const char * /*classname*/, const char *libname,
1987  Bool_t check)
1988 {
1989  Int_t err = -1;
1990 
1991  char *path;
1992  TString lib = libname;
1993  if (!lib.BeginsWith("lib"))
1994  lib = "lib" + lib;
1995  if ((path = gSystem->DynamicPathName(lib, kTRUE))) {
1996  if (check)
1997  err = 0;
1998  else {
1999  err = gSystem->Load(path, 0, kTRUE);
2000  }
2001  delete [] path;
2002  } else {
2003  if (check) {
2004  FileStat_t stat;
2005  if (!gSystem->GetPathInfo(libname, stat)) {
2006  if (R_ISREG(stat.fMode) &&
2008  err = 0;
2009  else
2010  err = -1;
2011  } else
2012  err = -1;
2013  } else {
2014  err = gSystem->Load(libname, 0, kTRUE);
2015  }
2016  }
2017 
2018  if (err == -1) {
2019  //Error("LoadClass", "library %s could not be loaded", libname);
2020  }
2021 
2022  if (err == 1) {
2023  //Error("LoadClass", "library %s already loaded, but class %s unknown",
2024  // libname, classname);
2025  err = 0;
2026  }
2027 
2028  return err;
2029 }
2030 
2031 ////////////////////////////////////////////////////////////////////////////////
2032 /// Return true if the file is local and is (likely) to be a ROOT file
2033 
2034 Bool_t TROOT::IsRootFile(const char *filename) const
2035 {
2036  Bool_t result = kFALSE;
2037  FILE *mayberootfile = fopen(filename,"rb");
2038  if (mayberootfile) {
2039  char header[5];
2040  if (fgets(header,5,mayberootfile)) {
2041  result = strncmp(header,"root",4)==0;
2042  }
2043  fclose(mayberootfile);
2044  }
2045  return result;
2046 }
2047 
2048 ////////////////////////////////////////////////////////////////////////////////
2049 /// To list all objects of the application.
2050 /// Loop on all objects created in the ROOT linked lists.
2051 /// Objects may be files and windows or any other object directly
2052 /// attached to the ROOT linked list.
2053 
2054 void TROOT::ls(Option_t *option) const
2055 {
2056 // TObject::SetDirLevel();
2057 // GetList()->R__FOR_EACH(TObject,ls)(option);
2058  TDirectory::ls(option);
2059 }
2060 
2061 ////////////////////////////////////////////////////////////////////////////////
2062 /// Load a macro in the interpreter's memory. Equivalent to the command line
2063 /// command ".L filename". If the filename has "+" or "++" appended
2064 /// the macro will be compiled by ACLiC. The filename must have the format:
2065 /// [path/]macro.C[+|++[g|O]].
2066 /// The possible error codes are defined by TInterpreter::EErrorCode.
2067 /// If check is true it will only check if filename exists and is
2068 /// readable.
2069 /// Returns 0 on successful loading and -1 in case filename does not
2070 /// exist or in case of error.
2071 
2072 Int_t TROOT::LoadMacro(const char *filename, int *error, Bool_t check)
2073 {
2074  Int_t err = -1;
2075  Int_t lerr, *terr;
2076  if (error)
2077  terr = error;
2078  else
2079  terr = &lerr;
2080 
2081  if (fInterpreter) {
2082  TString aclicMode;
2083  TString arguments;
2084  TString io;
2085  TString fname = gSystem->SplitAclicMode(filename, aclicMode, arguments, io);
2086 
2087  if (arguments.Length()) {
2088  Warning("LoadMacro", "argument(%s) ignored in %s", arguments.Data(), GetMacroPath());
2089  }
2090  char *mac = gSystem->Which(GetMacroPath(), fname, kReadPermission);
2091  if (!mac) {
2092  if (!check)
2093  Error("LoadMacro", "macro %s not found in path %s", fname.Data(), GetMacroPath());
2094  *terr = TInterpreter::kFatal;
2095  } else {
2096  err = 0;
2097  if (!check) {
2098  fname = mac;
2099  fname += aclicMode;
2100  fname += io;
2101  gInterpreter->LoadMacro(fname.Data(), (TInterpreter::EErrorCode*)terr);
2102  if (*terr)
2103  err = -1;
2104  }
2105  }
2106  delete [] mac;
2107  }
2108  return err;
2109 }
2110 
2111 ////////////////////////////////////////////////////////////////////////////////
2112 /// Execute a macro in the interpreter. Equivalent to the command line
2113 /// command ".x filename". If the filename has "+" or "++" appended
2114 /// the macro will be compiled by ACLiC. The filename must have the format:
2115 /// [path/]macro.C[+|++[g|O]][(args)].
2116 /// The possible error codes are defined by TInterpreter::EErrorCode.
2117 /// If padUpdate is true (default) update the current pad.
2118 /// Returns the macro return value.
2119 
2120 Long_t TROOT::Macro(const char *filename, Int_t *error, Bool_t padUpdate)
2121 {
2122  Long_t result = 0;
2123 
2124  if (fInterpreter) {
2125  TString aclicMode;
2126  TString arguments;
2127  TString io;
2128  TString fname = gSystem->SplitAclicMode(filename, aclicMode, arguments, io);
2129 
2130  char *mac = gSystem->Which(GetMacroPath(), fname, kReadPermission);
2131  if (!mac) {
2132  Error("Macro", "macro %s not found in path %s", fname.Data(), GetMacroPath());
2133  if (error)
2134  *error = TInterpreter::kFatal;
2135  } else {
2136  fname = mac;
2137  fname += aclicMode;
2138  fname += arguments;
2139  fname += io;
2140  result = gInterpreter->ExecuteMacro(fname, (TInterpreter::EErrorCode*)error);
2141  }
2142  delete [] mac;
2143 
2144  if (padUpdate && gPad)
2145  gPad->Update();
2146  }
2147 
2148  return result;
2149 }
2150 
2151 ////////////////////////////////////////////////////////////////////////////////
2152 /// Process message id called by obj.
2153 
2154 void TROOT::Message(Int_t id, const TObject *obj)
2155 {
2156  TIter next(fMessageHandlers);
2157  TMessageHandler *mh;
2158  while ((mh = (TMessageHandler*)next())) {
2159  mh->HandleMessage(id,obj);
2160  }
2161 }
2162 
2163 ////////////////////////////////////////////////////////////////////////////////
2164 /// Process interpreter command via TApplication::ProcessLine().
2165 /// On Win32 the line will be processed asynchronously by sending
2166 /// it to the CINT interpreter thread. For explicit synchronous processing
2167 /// use ProcessLineSync(). On non-Win32 platforms there is no difference
2168 /// between ProcessLine() and ProcessLineSync().
2169 /// The possible error codes are defined by TInterpreter::EErrorCode. In
2170 /// particular, error will equal to TInterpreter::kProcessing until the
2171 /// CINT interpreted thread has finished executing the line.
2172 /// Returns the result of the command, cast to a Long_t.
2173 
2174 Long_t TROOT::ProcessLine(const char *line, Int_t *error)
2175 {
2176  TString sline = line;
2177  sline = sline.Strip(TString::kBoth);
2178 
2179  if (!fApplication.load())
2181 
2182  return (*fApplication).ProcessLine(sline, kFALSE, error);
2183 }
2184 
2185 ////////////////////////////////////////////////////////////////////////////////
2186 /// Process interpreter command via TApplication::ProcessLine().
2187 /// On Win32 the line will be processed synchronously (i.e. it will
2188 /// only return when the CINT interpreter thread has finished executing
2189 /// the line). On non-Win32 platforms there is no difference between
2190 /// ProcessLine() and ProcessLineSync().
2191 /// The possible error codes are defined by TInterpreter::EErrorCode.
2192 /// Returns the result of the command, cast to a Long_t.
2193 
2195 {
2196  TString sline = line;
2197  sline = sline.Strip(TString::kBoth);
2198 
2199  if (!fApplication.load())
2201 
2202  return (*fApplication).ProcessLine(sline, kTRUE, error);
2203 }
2204 
2205 ////////////////////////////////////////////////////////////////////////////////
2206 /// Process interpreter command directly via CINT interpreter.
2207 /// Only executable statements are allowed (no variable declarations),
2208 /// In all other cases use TROOT::ProcessLine().
2209 /// The possible error codes are defined by TInterpreter::EErrorCode.
2210 
2212 {
2213  TString sline = line;
2214  sline = sline.Strip(TString::kBoth);
2215 
2216  if (!fApplication.load())
2218 
2219  Long_t result = 0;
2220 
2221  if (fInterpreter) {
2223  result = gInterpreter->Calc(sline, code);
2224  }
2225 
2226  return result;
2227 }
2228 
2229 ////////////////////////////////////////////////////////////////////////////////
2230 /// Read Git commit information and branch name from the
2231 /// etc/gitinfo.txt file.
2232 
2234 {
2235 #ifdef ROOT_GIT_COMMIT
2236  fGitCommit = ROOT_GIT_COMMIT;
2237 #endif
2238 #ifdef ROOT_GIT_BRANCH
2239  fGitBranch = ROOT_GIT_BRANCH;
2240 #endif
2241 
2242  TString gitinfo = "gitinfo.txt";
2243  char *filename = 0;
2244 #ifdef ROOTETCDIR
2245  filename = gSystem->ConcatFileName(ROOTETCDIR, gitinfo);
2246 #else
2247  TString etc = gRootDir;
2248 #ifdef WIN32
2249  etc += "\\etc";
2250 #else
2251  etc += "/etc";
2252 #endif
2253 #if defined(R__MACOSX) && (TARGET_OS_IPHONE || TARGET_IPHONE_SIMULATOR)
2254  // on iOS etc does not exist and gitinfo resides in $ROOTSYS
2255  etc = gRootDir;
2256 #endif
2257  filename = gSystem->ConcatFileName(etc, gitinfo);
2258 #endif
2259 
2260  FILE *fp = fopen(filename, "r");
2261  if (fp) {
2262  TString s;
2263  // read branch name
2264  s.Gets(fp);
2265  fGitBranch = s;
2266  // read commit SHA1
2267  s.Gets(fp);
2268  fGitCommit = s;
2269  // read date/time make was run
2270  s.Gets(fp);
2271  fGitDate = s;
2272  fclose(fp);
2273  }
2274  delete [] filename;
2275 }
2276 
2278  TTHREAD_TLS(Bool_t) fgReadingObject = false;
2279  return fgReadingObject;
2280 }
2281 
2282 ////////////////////////////////////////////////////////////////////////////////
2283 /// Deprecated (will be removed in next release).
2284 
2286 {
2287  return GetReadingObject();
2288 }
2289 
2291 {
2292  GetReadingObject() = flag;
2293 }
2294 
2295 
2296 ////////////////////////////////////////////////////////////////////////////////
2297 /// Return date/time make was run.
2298 
2299 const char *TROOT::GetGitDate()
2300 {
2301  if (fGitDate == "") {
2302  Int_t iday,imonth,iyear, ihour, imin;
2303  static const char *months[] = { "Jan", "Feb", "Mar", "Apr", "May", "Jun",
2304  "Jul", "Aug", "Sep", "Oct", "Nov", "Dec" };
2305  Int_t idate = gROOT->GetBuiltDate();
2306  Int_t itime = gROOT->GetBuiltTime();
2307  iday = idate%100;
2308  imonth = (idate/100)%100;
2309  iyear = idate/10000;
2310  ihour = itime/100;
2311  imin = itime%100;
2312  fGitDate.Form("%s %02d %4d, %02d:%02d:00", months[imonth-1], iday, iyear, ihour, imin);
2313  }
2314  return fGitDate;
2315 }
2316 
2317 ////////////////////////////////////////////////////////////////////////////////
2318 /// Refresh all browsers. Call this method when some command line
2319 /// command or script has changed the browser contents. Not needed
2320 /// for objects that have the kMustCleanup bit set. Most useful to
2321 /// update browsers that show the file system or other objects external
2322 /// to the running ROOT session.
2323 
2325 {
2326  TIter next(GetListOfBrowsers());
2327  TBrowser *b;
2328  while ((b = (TBrowser*) next()))
2329  b->SetRefreshFlag(kTRUE);
2330 }
2331 ////////////////////////////////////////////////////////////////////////////////
2332 /// Insure that the files, canvases and sockets are closed.
2333 
2334 static void CallCloseFiles()
2335 {
2336  if (TROOT::Initialized() && ROOT::Internal::gROOTLocal) {
2337  gROOT->CloseFiles();
2338  }
2339 }
2340 
2341 ////////////////////////////////////////////////////////////////////////////////
2342 /// Called by static dictionary initialization to register clang modules
2343 /// for headers. Calls TCling::RegisterModule() unless gCling
2344 /// is NULL, i.e. during startup, where the information is buffered in
2345 /// the static GetModuleHeaderInfoBuffer().
2346 
2347 void TROOT::RegisterModule(const char* modulename,
2348  const char** headers,
2349  const char** includePaths,
2350  const char* payloadCode,
2351  const char* fwdDeclCode,
2352  void (*triggerFunc)(),
2353  const TInterpreter::FwdDeclArgsToKeepCollection_t& fwdDeclsArgToSkip,
2354  const char** classesHeaders)
2355 {
2356 
2357  // First a side track to insure proper end of process behavior.
2358 
2359  // Register for each loaded dictionary (and thus for each library),
2360  // that we need to Close the ROOT files as soon as this library
2361  // might start being unloaded after main.
2362  //
2363  // By calling atexit here (rather than directly from within the
2364  // library) we make sure that this is not called if the library is
2365  // 'only' dlclosed.
2366 
2367  // On Ubuntu the linker strips the unused libraries. Eventhough
2368  // stressHistogram is explicitly linked against libNet, it is not
2369  // retained and thus is loaded only as needed in the middle part of
2370  // the execution. Concretely this also means that it is loaded
2371  // *after* the construction of the TApplication object and thus
2372  // after the registration (atexit) of the EndOfProcessCleanups
2373  // routine. Consequently, after the end of main, libNet is
2374  // unloaded before EndOfProcessCleanups is called. When
2375  // EndOfProcessCleanups is executed it indirectly needs the TClass
2376  // for TSocket and its search will use resources that have already
2377  // been unloaded (technically the function static in TUnixSystem's
2378  // DynamicPath and the dictionary from libNet).
2379 
2380  // Similarly, the ordering (before this commit) was broken in the
2381  // following case:
2382 
2383  // TApplication creation (EndOfProcessCleanups registration)
2384  // load UserLibrary
2385  // create TFile
2386  // Append UserObject to TFile
2387 
2388  // and after the end of main the order of execution was
2389 
2390  // unload UserLibrary
2391  // call EndOfProcessCleanups
2392  // Write the TFile
2393  // attempt to write the user object.
2394  // ....
2395 
2396  // where what we need is to have the files closen/written before
2397  // the unloading of the library.
2398 
2399  // To solve the problem we now register an atexit function for
2400  // every dictionary thus making sure there is at least one executed
2401  // before the first library tear down after main.
2402 
2403  // If atexit is called directly within a library's code, the
2404  // function will called *either* when the library is 'dlclose'd or
2405  // after then end of main (whichever comes first). We do *not*
2406  // want the files to be closed whenever a library is unloaded via
2407  // dlclose. To avoid this, we add the function (CallCloseFiles)
2408  // from the dictionary indirectly (via ROOT::RegisterModule). In
2409  // this case the function will only only be called either when
2410  // libCore is 'dlclose'd or right after the end of main.
2411 
2412  atexit(CallCloseFiles);
2413 
2414  // Now register with TCling.
2415  if (gCling) {
2416  gCling->RegisterModule(modulename, headers, includePaths, payloadCode, fwdDeclCode,
2417  triggerFunc, fwdDeclsArgToSkip, classesHeaders);
2418  } else {
2419  GetModuleHeaderInfoBuffer()
2420  .push_back(ModuleHeaderInfo_t (modulename, headers, includePaths, payloadCode, fwdDeclCode,
2421  triggerFunc, fwdDeclsArgToSkip,classesHeaders));
2422  }
2423 }
2424 
2425 ////////////////////////////////////////////////////////////////////////////////
2426 /// Remove a class from the list and map of classes.
2427 /// This routine is deprecated, use TClass::RemoveClass directly.
2428 
2430 {
2431  TClass::RemoveClass(oldcl);
2432 }
2433 
2434 ////////////////////////////////////////////////////////////////////////////////
2435 /// Delete all global interpreter objects created since the last call to Reset
2436 ///
2437 /// If option="a" is set reset to startup context (i.e. unload also
2438 /// all loaded files, classes, structs, typedefs, etc.).
2439 ///
2440 /// This function is typically used at the beginning (or end) of an unnamed macro
2441 /// to clean the environment.
2442 ///
2443 /// IMPORTANT WARNING:
2444 /// Do not use this call from within any function (neither compiled nor
2445 /// interpreted. This should only be used from a unnamed macro
2446 /// (which starts with a { (curly braces) ). For example, using TROOT::Reset
2447 /// from within an interpreted function will lead to the unloading of the
2448 /// dictionary and source file, including the one defining the function being
2449 /// executed.
2450 ///
2451 
2452 void TROOT::Reset(Option_t *option)
2453 {
2454  if (IsExecutingMacro()) return; //True when TMacro::Exec runs
2455  if (fInterpreter) {
2456  if (!strncmp(option, "a", 1)) {
2457  fInterpreter->Reset();
2459  } else
2460  gInterpreter->ResetGlobals();
2461 
2462  if (fGlobals) fGlobals->Unload();
2464 
2465  SaveContext();
2466  }
2467 }
2468 
2469 ////////////////////////////////////////////////////////////////////////////////
2470 /// Save the current interpreter context.
2471 
2473 {
2474  if (fInterpreter)
2475  gInterpreter->SaveGlobalsContext();
2476 }
2477 
2478 ////////////////////////////////////////////////////////////////////////////////
2479 /// Set the default graphical cut class name for the graphics editor
2480 /// By default the graphics editor creates an instance of a class TCutG.
2481 /// This function may be called to specify a different class that MUST
2482 /// derive from TCutG
2483 
2484 void TROOT::SetCutClassName(const char *name)
2485 {
2486  if (!name) {
2487  Error("SetCutClassName","Invalid class name");
2488  return;
2489  }
2490  TClass *cl = TClass::GetClass(name);
2491  if (!cl) {
2492  Error("SetCutClassName","Unknown class:%s",name);
2493  return;
2494  }
2495  if (!cl->InheritsFrom("TCutG")) {
2496  Error("SetCutClassName","Class:%s does not derive from TCutG",name);
2497  return;
2498  }
2499  fCutClassName = name;
2500 }
2501 
2502 ////////////////////////////////////////////////////////////////////////////////
2503 /// Set editor mode
2504 
2505 void TROOT::SetEditorMode(const char *mode)
2506 {
2507  fEditorMode = 0;
2508  if (!mode[0]) return;
2509  if (!strcmp(mode,"Arc")) {fEditorMode = kArc; return;}
2510  if (!strcmp(mode,"Line")) {fEditorMode = kLine; return;}
2511  if (!strcmp(mode,"Arrow")) {fEditorMode = kArrow; return;}
2512  if (!strcmp(mode,"Button")) {fEditorMode = kButton; return;}
2513  if (!strcmp(mode,"Diamond")) {fEditorMode = kDiamond; return;}
2514  if (!strcmp(mode,"Ellipse")) {fEditorMode = kEllipse; return;}
2515  if (!strcmp(mode,"Pad")) {fEditorMode = kPad; return;}
2516  if (!strcmp(mode,"Pave")) {fEditorMode = kPave; return;}
2517  if (!strcmp(mode,"PaveLabel")){fEditorMode = kPaveLabel; return;}
2518  if (!strcmp(mode,"PaveText")) {fEditorMode = kPaveText; return;}
2519  if (!strcmp(mode,"PavesText")){fEditorMode = kPavesText; return;}
2520  if (!strcmp(mode,"PolyLine")) {fEditorMode = kPolyLine; return;}
2521  if (!strcmp(mode,"CurlyLine")){fEditorMode = kCurlyLine; return;}
2522  if (!strcmp(mode,"CurlyArc")) {fEditorMode = kCurlyArc; return;}
2523  if (!strcmp(mode,"Text")) {fEditorMode = kText; return;}
2524  if (!strcmp(mode,"Marker")) {fEditorMode = kMarker; return;}
2525  if (!strcmp(mode,"CutG")) {fEditorMode = kCutG; return;}
2526 }
2527 
2528 ////////////////////////////////////////////////////////////////////////////////
2529 /// Change current style to style with name stylename
2530 
2531 void TROOT::SetStyle(const char *stylename)
2532 {
2533  TString style_name = stylename;
2534 
2535  TStyle *style = GetStyle(style_name);
2536  if (style) style->cd();
2537  else Error("SetStyle","Unknown style:%s",style_name.Data());
2538 }
2539 
2540 
2541 //-------- Static Member Functions ---------------------------------------------
2542 
2543 
2544 ////////////////////////////////////////////////////////////////////////////////
2545 /// Decrease the indentation level for ls().
2546 
2548 {
2549  return --fgDirLevel;
2550 }
2551 
2552 ////////////////////////////////////////////////////////////////////////////////
2553 ///return directory level
2554 
2556 {
2557  return fgDirLevel;
2558 }
2559 
2560 ////////////////////////////////////////////////////////////////////////////////
2561 /// Get macro search path. Static utility function.
2562 
2563 const char *TROOT::GetMacroPath()
2564 {
2565  TString &macroPath = ROOT::GetMacroPath();
2566 
2567  if (macroPath.Length() == 0) {
2568  macroPath = gEnv->GetValue("Root.MacroPath", (char*)0);
2569 #if defined(R__WIN32)
2570  macroPath.ReplaceAll("; ", ";");
2571 #else
2572  macroPath.ReplaceAll(": ", ":");
2573 #endif
2574  if (macroPath.Length() == 0)
2575 #if !defined(R__WIN32)
2576  #ifdef ROOTMACRODIR
2577  macroPath = ".:" ROOTMACRODIR;
2578  #else
2579  macroPath = TString(".:") + gRootDir + "/macros";
2580  #endif
2581 #else
2582  #ifdef ROOTMACRODIR
2583  macroPath = ".;" ROOTMACRODIR;
2584  #else
2585  macroPath = TString(".;") + gRootDir + "/macros";
2586  #endif
2587 #endif
2588  }
2589 
2590  return macroPath;
2591 }
2592 
2593 ////////////////////////////////////////////////////////////////////////////////
2594 /// Set or extend the macro search path. Static utility function.
2595 /// If newpath=0 or "" reset to value specified in the rootrc file.
2596 
2597 void TROOT::SetMacroPath(const char *newpath)
2598 {
2599  TString &macroPath = ROOT::GetMacroPath();
2600 
2601  if (!newpath || !*newpath)
2602  macroPath = "";
2603  else
2604  macroPath = newpath;
2605 }
2606 
2607 ////////////////////////////////////////////////////////////////////////////////
2608 /// Increase the indentation level for ls().
2609 
2611 {
2612  return ++fgDirLevel;
2613 }
2614 
2615 ////////////////////////////////////////////////////////////////////////////////
2616 /// Functions used by ls() to indent an object hierarchy.
2617 
2619 {
2620  for (int i = 0; i < fgDirLevel; i++) std::cout.put(' ');
2621 }
2622 
2623 ////////////////////////////////////////////////////////////////////////////////
2624 /// Return kTRUE if the TROOT object has been initialized.
2625 
2627 {
2628  return fgRootInit;
2629 }
2630 
2631 ////////////////////////////////////////////////////////////////////////////////
2632 /// Return kTRUE if the memory leak checker is on.
2633 
2635 {
2636  return fgMemCheck;
2637 }
2638 
2639 ////////////////////////////////////////////////////////////////////////////////
2640 /// Return Indentation level for ls().
2641 
2643 {
2644  fgDirLevel = level;
2645 }
2646 
2647 ////////////////////////////////////////////////////////////////////////////////
2648 /// Convert version code to an integer, i.e. 331527 -> 51507.
2649 
2651 {
2652  return 10000*(code>>16) + 100*((code&65280)>>8) + (code&255);
2653 }
2654 
2655 ////////////////////////////////////////////////////////////////////////////////
2656 /// Convert version as an integer to version code as used in RVersion.h.
2657 
2659 {
2660  int a = v/10000;
2661  int b = (v - a*10000)/100;
2662  int c = v - a*10000 - b*100;
2663  return (a << 16) + (b << 8) + c;
2664 }
2665 
2666 ////////////////////////////////////////////////////////////////////////////////
2667 /// Return ROOT version code as defined in RVersion.h.
2668 
2670 {
2671  return ROOT_VERSION_CODE;
2672 }
2673 
2674 ////////////////////////////////////////////////////////////////////////////////
2675 
2677  static const char** extraInterpArgs = 0;
2678  return extraInterpArgs;
2679 }
2680 
2681 ////////////////////////////////////////////////////////////////////////////////
2682 /// Get the tutorials directory in the installation. Static utility function.
2683 
2685 {
2686 #ifdef ROOTTUTDIR
2687  return ROOTTUTDIR;
2688 #else
2689  static TString tutdir = TString(gRootDir) + "/tutorials";
2690  return tutdir;
2691 #endif
2692 }
TSeqCollection * fStreamerInfo
Definition: TROOT.h:168
void Add(TObject *obj, const char *name=0, Int_t check=-1)
Add object with name to browser.
Definition: TBrowser.cxx:261
const char * FindObjectPathName(const TObject *obj) const
Return path name of obj somewhere in the //root/...
Definition: TROOT.cxx:1309
virtual const char * BaseName(const char *pathname)
Base name of a file name. Base name of /user/root is root.
Definition: TSystem.cxx:929
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1266
Bool_t fExecutingMacro
Definition: TROOT.h:145
void AddClass(TClass *cl)
Add a class to the list and map of classes.
Definition: TROOT.cxx:939
virtual void RegisterModule(const char *, const char **, const char **, const char *, const char *, void(*)(), const FwdDeclArgsToKeepCollection_t &fwdDeclArgsToKeep, const char **classesHeaders, Bool_t lateRegistration=false)=0
virtual void Add(TObject *obj)
TListOfFunctionTemplates * fFuncTemplate
Definition: TROOT.h:151
object has not been deleted
Definition: TObject.h:69
TInterpreter * CreateInterpreter_t(void *shlibHandle)
Definition: TInterpreter.h:513
A TFolder object is a collection of objects and folders.
Definition: TFolder.h:32
TCanvas * MakeDefCanvas() const
Return a default canvas.
Definition: TROOT.cxx:1390
R__EXTERN TGuiFactory * gBatchGuiFactory
Definition: TGuiFactory.h:69
A collection of TDataMember objects designed for fast access given a DeclId_t and for keep track of T...
Semi-Abstract base class defining a generic interface to the underlying, low level, native graphics backend (X11, Win32, MacOS, OpenGL...).
Definition: TVirtualX.h:70
An array of TObjects.
Definition: TObjArray.h:39
virtual void Clear(Option_t *option="")=0
static Int_t DecreaseDirLevel()
Decrease the indentation level for ls().
Definition: TROOT.cxx:2547
static Bool_t BlockAllSignals(Bool_t b)
Block or unblock all signals. Returns the previous block status.
Definition: TQObject.cxx:1325
TFolder * fRootFolder
Definition: TROOT.h:176
TCollection * GetListOfEnums(Bool_t load=kFALSE)
Definition: TROOT.cxx:1579
ROOT top level object description.
Definition: TROOT.h:104
void * DestroyInterpreter_t(TInterpreter *)
Definition: TInterpreter.h:514
This class is a specialized TProcessID managing the list of UUIDs.
Definition: TProcessUUID.h:34
TSeqCollection * fProofs
Definition: TROOT.h:171
const char * GetDeclFileName() const
Definition: TClass.h:386
void RemoveClass(TClass *)
Remove a class from the list and map of classes.
Definition: TROOT.cxx:2429
TCollection * GetListOfGlobalFunctions(Bool_t load=kFALSE)
Return list containing the TFunctions currently defined.
Definition: TROOT.cxx:1648
virtual const char * WorkingDirectory()
Return working directory.
Definition: TSystem.cxx:866
Int_t LoadMacro(const char *filename, Int_t *error=0, Bool_t check=kFALSE)
Load a macro in the interpreter&#39;s memory.
Definition: TROOT.cxx:2072
virtual void Build(TFile *motherFile=0, TDirectory *motherDir=0)
Initialise directory to defaults.
Definition: TDirectory.cxx:200
This namespace contains pre-defined functions to be used in conjuction with TExecutor::Map and TExecu...
Definition: StringConv.hxx:21
static const char **& GetExtraInterpreterArgs()
Definition: TROOT.cxx:2676
TLine * line
TSeqCollection * fGeometries
Definition: TROOT.h:163
static Func_t GetSymInLibThread(const char *funcname)
Definition: TROOT.cxx:374
R__EXTERN TClassTable * gClassTable
Definition: TClassTable.h:97
Handle messages that might be generated by the system.
TDictionary::DeclId_t DeclId_t
Definition: TInterpreter.h:251
return c
const char Option_t
Definition: RtypesCore.h:62
Dictionary for function template This class describes one single function template.
virtual TString SplitAclicMode(const char *filename, TString &mode, TString &args, TString &io) const
This method split a filename of the form: ~~~ {.cpp} [path/]macro.C[+|++[k|f|g|O|c|s|d|v|-]][(args)]...
Definition: TSystem.cxx:4109
void DisableImplicitMT()
Disables the implicit multi-threading in ROOT.
Definition: TROOT.cxx:532
static TVirtualX *& Instance()
Returns gVirtualX global.
Definition: TVirtualX.cxx:57
A collection of TFunction objects designed for fast access given a DeclId_t and for keep track of TFu...
Short_t GetDeclFileLine() const
Definition: TClass.h:387
void EnableImplicitMT(UInt_t numthreads=0)
Globally enables the implicit multi-threading in ROOT, activating the parallel execution of those met...
Definition: TROOT.cxx:518
#define ROOT_RELEASE_TIME
Definition: RVersion.h:19
TList * fList
Definition: TDirectory.h:99
R__EXTERN TVirtualMutex * gInterpreterMutex
Definition: TInterpreter.h:46
This class represents a WWW compatible URL.
Definition: TUrl.h:41
Int_t fVersionDate
Definition: TROOT.h:127
Definition: Buttons.h:33
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
int GetPathInfo(const char *path, Long_t *id, Long_t *size, Long_t *flags, Long_t *modtime)
Get info about a file: id, size, flags, modification time.
Definition: TSystem.cxx:1364
Bool_t fForceStyle
Definition: TROOT.h:142
virtual void HandleMessage(Int_t id, const TObject *obj)
Store message origin, keep statistics and call Notify().
Bool_t ReadingObject() const
Deprecated (will be removed in next release).
Definition: TROOT.cxx:2285
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:157
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
static Int_t RootVersionCode()
Return ROOT version code as defined in RVersion.h.
Definition: TROOT.cxx:2669
static Bool_t MemCheck()
Return kTRUE if the memory leak checker is on.
Definition: TROOT.cxx:2634
TString fVersion
Definition: TROOT.h:124
void R__SetZipMode(int)
virtual const char * HomeDirectory(const char *userName=0)
Return the user&#39;s home directory.
Definition: TSystem.cxx:882
#define ROOT_RELEASE_DATE
Definition: RVersion.h:18
virtual void AddAll(const TCollection *col)
Add all objects from collection col to this collection.
Definition: TCollection.cxx:58
static void SetObjectStat(Bool_t stat)
Turn on/off tracking of objects in the TObjectTable.
Definition: TObject.cxx:1006
This class implements a mutex interface.
Definition: TVirtualMutex.h:34
void DisableParTreeProcessing()
Globally disables the IMT use case of parallel branch processing, deactivating the corresponding lock...
Definition: TROOT.cxx:458
void InitInterpreter()
Initialize the interpreter.
Definition: TROOT.cxx:1861
#define gROOT
Definition: TROOT.h:364
void LoadHandlersFromEnv(TEnv *env)
Load plugin handlers specified in config file, like: Plugin.TFile: ^rfio: TRFIOFile RFI...
TROOT * GetROOT2()
Definition: TROOT.cxx:360
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1819
The TEnv class reads config files, by default named .rootrc.
Definition: TEnv.h:128
void RefreshBrowsers()
Refresh all browsers.
Definition: TROOT.cxx:2324
TString & GetMacroPath()
Definition: TROOT.cxx:491
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
std::vector< std::pair< std::string, int > > FwdDeclArgsToKeepCollection_t
Definition: TInterpreter.h:108
const Bool_t kFALSE
Definition: Rtypes.h:92
#define gInterpreter
Definition: TInterpreter.h:517
virtual char * Which(const char *search, const char *file, EAccessMode mode=kFileExists)
Find location of file in a search path.
Definition: TSystem.cxx:1512
Int_t fVersionCode
Definition: TROOT.h:126
Option_t * GetOption() const
Definition: TCollection.h:160
void Reset(Option_t *option="")
Delete all global interpreter objects created since the last call to Reset.
Definition: TROOT.cxx:2452
Definition: Buttons.h:30
void Load()
Load all the DataMembers known to the interpreter for the scope &#39;fClass&#39; into this collection...
TFunction * Get(DeclId_t id)
Return (after creating it if necessary) the TMethod or TFunction describing the function correspondin...
This class registers for all classes their name, id and dictionary function in a hash table...
Definition: TClassTable.h:40
STL namespace.
#define ROOT_VERSION_CODE
Definition: RVersion.h:21
virtual TObject * FindObjectAny(const char *name) const
Return a pointer to the first object with name starting at this folder.
Definition: TFolder.cxx:352
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
virtual void ResetSignals()
Reset signals handlers to previous behaviour.
Definition: TSystem.cxx:589
std::atomic< TApplication * > fApplication
Definition: TROOT.h:135
void SetStyle(const char *stylename="Default")
Change current style to style with name stylename.
Definition: TROOT.cxx:2531
TDataType * GetType(const char *name, Bool_t load=kFALSE) const
Return pointer to type with name.
Definition: TROOT.cxx:1398
TCollection * GetListOfGlobals(Bool_t load=kFALSE)
Return list containing the TGlobals currently defined.
Definition: TROOT.cxx:1613
TString fGitCommit
Definition: TROOT.h:131
TStyle * GetStyle(const char *name) const
Return pointer to style with name.
Definition: TROOT.cxx:1415
void ls(Option_t *option="") const
To list all objects of the application.
Definition: TROOT.cxx:2054
TROOT * GetROOT1()
Definition: TROOT.cxx:353
Bool_t R_ISREG(Int_t mode)
Definition: TSystem.h:129
virtual void AddLast(TObject *obj)
Add object at the end of the list.
Definition: TList.cxx:137
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
TString fCutClassName
Definition: TROOT.h:179
Bool_t & GetReadingObject()
Definition: TROOT.cxx:2277
static const char * GetTutorialsDir()
Get the tutorials directory in the installation. Static utility function.
Definition: TROOT.cxx:2684
TListOfFunctions * GetGlobalFunctions()
Internal routine returning, and creating if necessary, the list of global function.
Definition: TROOT.cxx:1488
static Int_t fgDirLevel
Definition: TROOT.h:112
Definition: Buttons.h:30
static const char * GetMacroPath()
Get macro search path. Static utility function.
Definition: TROOT.cxx:2563
TInterpreter * fInterpreter
Definition: TROOT.h:136
Int_t fMode
Definition: TSystem.h:138
static TList & GetEarlyRegisteredGlobals()
Definition: TGlobal.cxx:171
const char * GetGitDate()
Return date/time make was run.
Definition: TROOT.cxx:2299
virtual void cd()
Change current style.
Definition: TStyle.cxx:306
AListOfEnums_t fEnums
Definition: TROOT.h:174
TObject * GetFunction(const char *name) const
Return pointer to function with name.
Definition: TROOT.cxx:1423
TSeqCollection * fMappedFiles
Definition: TROOT.h:156
virtual void Delete(Option_t *option="")=0
Delete this object.
#define SafeDelete(p)
Definition: RConfig.h:507
Sequenceable collection abstract base class.
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:188
static void * gInterpreterLib
Definition: TROOT.cxx:156
TVirtualMutex * gROOTMutex
Definition: TROOT.cxx:159
THashTable implements a hash table to store TObject&#39;s.
Definition: THashTable.h:39
static TVirtualPad *& Pad()
Return the current pad for the current thread.
Definition: TVirtualPad.cxx:32
Int_t fVersionInt
Definition: TROOT.h:125
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition: THashList.h:36
Bool_t IsParTreeProcessingEnabled()
Returns true if parallel tree processing is enabled.
Definition: TROOT.cxx:471
void SetEditorMode(const char *mode="")
Set editor mode.
Definition: TROOT.cxx:2505
static Int_t ReadRules()
Read the class.rules files from the default location:.
Definition: TClass.cxx:1732
TSeqCollection * fCleanups
Definition: TROOT.h:166
Int_t fEditorMode
Definition: TROOT.h:146
Bool_t fMustClean
Definition: TROOT.h:140
static Bool_t Initialized()
Return kTRUE if the TROOT object has been initialized.
Definition: TROOT.cxx:2626
void EndOfProcessCleanups()
Execute the cleanups necessary at the end of the process, in particular those that must be executed b...
Definition: TROOT.cxx:1101
static Bool_t fgRootInit
Definition: TROOT.h:113
virtual TObject * FindObjectAny(const char *name) const
Return a pointer to the first object with name starting at //root.
Definition: TROOT.cxx:1262
Long_t ProcessLineSync(const char *line, Int_t *error=0)
Process interpreter command via TApplication::ProcessLine().
Definition: TROOT.cxx:2194
void Message(Int_t id, const TObject *obj)
Process message id called by obj.
Definition: TROOT.cxx:2154
void EnableParTreeProcessing()
Globally enables the parallel tree processing, which is a case of implicit multi-threading in ROOT...
Definition: TROOT.cxx:442
virtual const char * Getenv(const char *env)
Get environment variable.
Definition: TSystem.cxx:1628
TCollection * GetListOfFunctionTemplates()
Definition: TROOT.cxx:1596
static Int_t ITIMQQ(const char *time)
Return built time as integer (with min precision), i.e.
Definition: TROOT.cxx:199
virtual void EnableAutoLoading()=0
void SaveContext()
Save the current interpreter context.
Definition: TROOT.cxx:2472
virtual void Close(Option_t *option="")
Delete all objects from memory and directory structure itself.
Definition: TDirectory.cxx:519
ClassInfo_t * GetClassInfo() const
Definition: TClass.h:391
virtual void Delete(Option_t *option="")
Delete all TFunction object files.
TCollection * GetListOfFunctionOverloads(const char *name) const
Return the collection of functions named "name".
Definition: TROOT.cxx:1497
static void at_exit_of_TROOT()
Definition: TROOT.cxx:271
TSeqCollection * fClipboard
Definition: TROOT.h:172
virtual void Initialize()=0
TFolder * AddFolder(const char *name, const char *title, TCollection *collection=0)
Create a new folder and add it to the list of folders of this folder, return a pointer to the created...
Definition: TFolder.cxx:186
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:59
R__EXTERN TGuiFactory * gGuiFactory
Definition: TGuiFactory.h:68
static void Cleanup()
static function (called by TROOT destructor) to delete all TProcessIDs
Definition: TProcessID.cxx:188
void Error(const char *location, const char *msgfmt,...)
R__EXTERN TPluginManager * gPluginMgr
Definition: Buttons.h:38
Bool_t fEscape
Definition: TROOT.h:144
TSeqCollection * fFiles
Definition: TROOT.h:155
TSeqCollection * GetListOfGeometries() const
Definition: TROOT.h:252
Describes an Operating System directory for the browser.
virtual TObject * FindObject(const char *name) const
Specialize FindObject to do search for the a function just by name or create it if its not already in...
TCollection * fClasses
Definition: TROOT.h:149
A doubly linked list.
Definition: TList.h:47
Long_t Macro(const char *filename, Int_t *error=0, Bool_t padUpdate=kTRUE)
Execute a macro in the interpreter.
Definition: TROOT.cxx:2120
R__EXTERN TVirtualMutex * gGlobalMutex
Definition: TVirtualMutex.h:29
TStyle objects may be created to define special styles.
Definition: TStyle.h:43
R__EXTERN TROOT * gROOTLocal
Definition: TROOT.h:361
Bool_t fEditHistograms
Definition: TROOT.h:138
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:41
virtual void AddLast(TObject *obj)=0
Bool_t ClassSaved(TClass *cl)
return class status bit kClassSaved for class cl This function is called by the SavePrimitive functio...
Definition: TROOT.cxx:977
Int_t fVersionTime
Definition: TROOT.h:128
static void RegisterModule(const char *modulename, const char **headers, const char **includePaths, const char *payLoadCode, const char *fwdDeclCode, void(*triggerFunc)(), const FwdDeclArgsToKeepCollection_t &fwdDeclsArgToSkip, const char **classesHeaders)
Called by static dictionary initialization to register clang modules for headers. ...
Definition: TROOT.cxx:2347
R__EXTERN TVirtualX * gGXBatch
Definition: TVirtualX.h:365
TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE) const
Return pointer to class with name. Obsolete, use TClass::GetClass directly.
Definition: TROOT.cxx:1354
static Bool_t fgMemCheck
Definition: TROOT.h:114
Bool_t IsParBranchProcessingEnabled()
Returns true if parallel branch processing is enabled.
Definition: TROOT.cxx:422
TSeqCollection * fDataSets
Definition: TROOT.h:173
Definition: Buttons.h:33
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:380
TSeqCollection * fMessageHandlers
Definition: TROOT.h:167
void Browse(TBrowser *b)
Add browsable objects to TBrowser.
Definition: TROOT.cxx:958
virtual Bool_t Init()
Initialize the OS interface.
Definition: TSystem.cxx:187
void AddClassGenerator(TClassGenerator *gen)
Add a class generator.
Definition: TROOT.cxx:949
TString fDefCanvasName
Definition: TROOT.h:180
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
void SetCutClassName(const char *name="TCutG")
Set the default graphical cut class name for the graphics editor By default the graphics editor creat...
Definition: TROOT.cxx:2484
TGlobal * GetGlobal(const char *name, Bool_t load=kFALSE) const
Return pointer to global variable by name.
Definition: TROOT.cxx:1456
Definition: Buttons.h:33
SVector< double, 2 > v
Definition: Dict.h:5
if object ctor succeeded but object should not be used
Definition: TObject.h:63
void GetNormalizedName(std::string &norm_name, std::string_view name)
Return the normalized name.
Definition: TClassEdit.cxx:788
Basic data type descriptor (datatype information is obtained from CINT).
Definition: TDataType.h:46
Int_t gDebug
Definition: TROOT.cxx:566
virtual Int_t GetValue(const char *name, Int_t dflt)
Returns the integer value for a resource.
Definition: TEnv.cxx:496
TPluginManager * fPluginManager
Definition: TROOT.h:178
TString fGitBranch
Definition: TROOT.h:132
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:488
TSeqCollection * fSecContexts
Definition: TROOT.h:170
Collection abstract base class.
Definition: TCollection.h:48
void(* VoidFuncPtr_t)()
Definition: Rtypes.h:74
void SetRefreshFlag(Bool_t flag)
Definition: TBrowser.h:101
static Int_t ConvertVersionCode2Int(Int_t code)
Convert version code to an integer, i.e. 331527 -> 51507.
Definition: TROOT.cxx:2650
unsigned int UInt_t
Definition: RtypesCore.h:42
void ReadGitInfo()
Read Git commit information and branch name from the etc/gitinfo.txt file.
Definition: TROOT.cxx:2233
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
virtual TObject * FindObject(const char *name) const
Returns address of a ROOT object if it exists.
Definition: TROOT.cxx:1149
TObject * GetGeometry(const char *name) const
Return pointer to Geometry with name.
Definition: TROOT.cxx:1572
static void EnableStatistics(int size=-1, int ix=-1)
Enable memory usage statistics gathering.
Definition: TStorage.cxx:450
This class implements a plugin library manager.
Int_t LoadClass(const char *classname, const char *libname, Bool_t check=kFALSE)
Check if class "classname" is known to the interpreter (in fact, this check is not needed anymore...
Definition: TROOT.cxx:1986
Objects following this interface can be passed onto the TROOT object to implement a user customized w...
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2618
static TProcessID * AddProcessID()
Static function to add a new TProcessID to the list of PIDs.
Definition: TProcessID.cxx:100
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
TFile * GetFile() const
Definition: TROOT.h:268
virtual void SaveContext()=0
Global variables class (global variables are obtained from CINT).
Definition: TGlobal.h:29
Bool_t InheritsFrom(const char *cl) const
Return kTRUE if this class inherits from a class with name "classname".
Definition: TClass.cxx:4610
A collection of TDataType designed to hold the typedef information and numerical type information...
Definition: TListOfTypes.h:32
void SetName(const char *name)
Definition: TCollection.h:116
TROOT * GetROOT()
Definition: TROOT.cxx:487
void Warning(const char *location, const char *msgfmt,...)
Bool_t IsRootFile(const char *filename) const
Return true if the file is local and is (likely) to be a ROOT file.
Definition: TROOT.cxx:2034
static TInterpreter * Instance()
returns gInterpreter global
void InitSystem()
Initialize operating system interface.
Definition: TROOT.cxx:1785
if object destructor must call RecursiveRemove()
Definition: TObject.h:57
#define gVirtualX
Definition: TVirtualX.h:362
TSeqCollection * fSpecials
Definition: TROOT.h:165
virtual TObjLink * FirstLink() const
Definition: TList.h:101
#define R__LOCKGUARD2(mutex)
static void CreateApplication()
Static function used to create a default application environment.
static void RemoveClass(TClass *cl)
static: Remove a class from the list and map of classes
Definition: TClass.cxx:487
TString fGitDate
Definition: TROOT.h:133
static void SetDirLevel(Int_t level=0)
Return Indentation level for ls().
Definition: TROOT.cxx:2642
long Long_t
Definition: RtypesCore.h:50
The Canvas class.
Definition: TCanvas.h:41
Bool_t IsExecutingMacro() const
Definition: TROOT.h:287
TSeqCollection * GetListOfColors() const
Definition: TROOT.h:240
TListOfFunctions * fGlobalFunctions
Definition: TROOT.h:153
Long_t ProcessLine(const char *line, Int_t *error=0)
Process interpreter command via TApplication::ProcessLine().
Definition: TROOT.cxx:2174
A collection of TEnum objects designed for fast access given a DeclId_t and for keep track of TEnum t...
TSeqCollection * fCanvases
Definition: TROOT.h:158
static void BuildStyles()
Create some standard styles.
Definition: TStyle.cxx:290
virtual void SetName(const char *newname)
Set the name for directory If the directory name is changed after the directory was written once...
virtual const char * GetTypeName() const
Get type of global variable, e,g.
Definition: TGlobal.cxx:111
#define ClassImp(name)
Definition: Rtypes.h:279
double f(double x)
void EnableThreadSafety()
Enables the global mutex to make ROOT thread safe/aware.
Definition: TROOT.cxx:498
TFunction * GetGlobalFunctionWithPrototype(const char *name, const char *proto=0, Bool_t load=kFALSE)
Return pointer to global function by name.
Definition: TROOT.cxx:1543
static DestroyInterpreter_t * gDestroyInterpreter
Definition: TROOT.cxx:155
void SetReadingObject(Bool_t flag=kTRUE)
Definition: TROOT.cxx:2290
virtual Func_t DynFindSymbol(const char *module, const char *entry)
Find specific entry point in specified library.
Definition: TSystem.cxx:1983
char * DynamicPathName(const char *lib, Bool_t quiet=kFALSE)
Find a dynamic library called lib using the system search paths.
Definition: TSystem.cxx:1959
TString fConfigOptions
Definition: TROOT.h:122
void InitThreads()
Load and initialize thread library.
Definition: TROOT.cxx:1850
void Unload()
Mark &#39;all func&#39; as being unloaded.
Describe directory structure in memory.
Definition: TDirectory.h:44
static TDirectory *& CurrentDirectory()
Return the current directory for the current thread.
Definition: TDirectory.cxx:318
R__EXTERN TEnv * gEnv
Definition: TEnv.h:174
TSeqCollection * GetListOfStyles() const
Definition: TROOT.h:249
static void AddClass(TClass *cl)
static: Add a class to the list and map of classes.
Definition: TClass.cxx:461
TSeqCollection * GetListOfFiles() const
Definition: TROOT.h:245
Definition: Buttons.h:31
Bool_t fFromPopUp
Definition: TROOT.h:139
static RooMathCoreReg dummy
Int_t fTimer
Definition: TROOT.h:134
TCanvas * style()
Definition: style.C:1
static void CallCloseFiles()
Insure that the files, canvases and sockets are closed.
Definition: TROOT.cxx:2334
void EnableParBranchProcessing()
Globally enables the parallel branch processing, which is a case of implicit multi-threading (IMT) in...
Definition: TROOT.cxx:393
Bool_t fBatch
Definition: TROOT.h:137
Long_t ProcessLineFast(const char *line, Int_t *error=0)
Process interpreter command directly via CINT interpreter.
Definition: TROOT.cxx:2211
void CloseFiles()
Close any files and sockets that gROOT knows about.
Definition: TROOT.cxx:1023
virtual ~TROOT()
Clean up and free resources used by ROOT (files, network sockets, shared memory segments, etc.).
Definition: TROOT.cxx:814
#define R__LOCKGUARD(mutex)
The color creation and management class.
Definition: TColor.h:23
TSeqCollection * fStyles
Definition: TROOT.h:159
virtual void Add(TObject *obj)=0
virtual void CleanCompiledMacros()
Remove the shared libs produced by the CompileMacro() function.
Definition: TSystem.cxx:4194
Bool_t fReadingObject
Definition: TROOT.h:141
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:2893
virtual void Clear(Option_t *option="")
Remove all objects from the list.
Definition: TList.cxx:349
TCollection * fClassGenerators
Definition: TROOT.h:169
Int_t IgnoreInclude(const char *fname, const char *expandedfname)
Return 1 if the name of the given include file corresponds to a class that is known to ROOT...
Definition: TROOT.cxx:1735
static TClass * LoadClass(const char *requestedname, Bool_t silent)
Helper function used by TClass::GetClass().
Definition: TClass.cxx:5401
virtual void Reset()=0
virtual void Delete(Option_t *option="")
Delete all TDataMember object files.
static Int_t GetDirLevel()
return directory level
Definition: TROOT.cxx:2555
Int_t fBuiltTime
Definition: TROOT.h:130
Mother of all ROOT objects.
Definition: TObject.h:37
Global functions class (global functions are obtained from CINT).
Definition: TFunction.h:30
This ABC is a factory for GUI components.
Definition: TGuiFactory.h:44
TFunction * GetGlobalFunction(const char *name, const char *params=0, Bool_t load=kFALSE)
Return pointer to global function by name.
Definition: TROOT.cxx:1510
std::vector< std::pair< std::string, int > > FwdDeclArgsToKeepCollection_t
Definition: TROOT.h:196
static Int_t ConvertVersionInt2Code(Int_t v)
Convert version as an integer to version code as used in RVersion.h.
Definition: TROOT.cxx:2658
static Int_t IncreaseDirLevel()
Increase the indentation level for ls().
Definition: TROOT.cxx:2610
typedef void((*Func_t)())
TCollection * GetListOfTypes(Bool_t load=kFALSE)
Return a dynamic list giving access to all TDataTypes (typedefs) currently defined.
Definition: TROOT.cxx:1687
Int_t fBuiltDate
Definition: TROOT.h:129
void Unload()
Mark &#39;all func&#39; as being unloaded.
TFunctionTemplate * GetFunctionTemplate(const char *name)
Definition: TROOT.cxx:1443
TColor * GetColor(Int_t color) const
Return address of color with index color.
Definition: TROOT.cxx:1372
void Load()
Load all the functions known to the interpreter for the scope &#39;fClass&#39; into this collection.
TCollection * fFunctions
Definition: TROOT.h:160
const char * FindObjectClassName(const char *name) const
Returns class name of a ROOT object including CINT globals.
Definition: TROOT.cxx:1289
static void SetMacroPath(const char *newpath)
Set or extend the macro search path.
Definition: TROOT.cxx:2597
TROOT *(* GetROOTFun_t)()
Definition: TROOT.cxx:370
static Int_t IDATQQ(const char *date)
Return built date as integer, i.e. "Apr 28 2000" -> 20000428.
Definition: TROOT.cxx:178
virtual void Add(TObject *obj)
Definition: TList.h:81
const Ssiz_t kNPOS
Definition: Rtypes.h:115
TListOfDataMembers * fGlobals
Definition: TROOT.h:152
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:545
R__EXTERN const char * gRootDir
Definition: TSystem.h:233
virtual void ls(Option_t *option="") const
List Directory contents.
TSeqCollection * fClosedObjects
Definition: TROOT.h:154
TClass * FindSTLClass(const char *name, Bool_t load, Bool_t silent=kFALSE) const
return a TClass object corresponding to &#39;name&#39; assuming it is an STL container.
Definition: TROOT.cxx:1320
#define ROOT_RELEASE
Definition: RVersion.h:17
TF1 * f1
Definition: legend1.C:11
static Int_t IVERSQ()
Return version id as an integer, i.e. "2.22/04" -> 22204.
Definition: TROOT.cxx:168
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
Int_t fLineIsProcessing
Definition: TROOT.h:110
TSeqCollection * fTasks
Definition: TROOT.h:161
static void CleanUpROOTAtExit()
Clean up at program termination before global objects go out of scope.
Definition: TROOT.cxx:209
TSeqCollection * fColors
Definition: TROOT.h:162
#define NULL
Definition: Rtypes.h:82
Int_t GetNumber() const
Definition: TColor.h:58
TCollection * fTypes
Definition: TROOT.h:150
#define gPad
Definition: TVirtualPad.h:289
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
const char * proto
Definition: civetweb.c:11652
TSeqCollection * fBrowsers
Definition: TROOT.h:164
static GetROOTFun_t gGetROOT
Definition: TROOT.cxx:372
#define gDirectory
Definition: TDirectory.h:221
static void InitializeColors()
Initialize colors used by the TCanvas based graphics (via TColor objects).
Definition: TColor.cxx:1048
virtual TObject * FindObjectAnyFile(const char *name) const
Scan the memory lists of all files for an object with name.
Definition: TROOT.cxx:1272
double result[121]
TList * fBrowsables
Definition: TROOT.h:177
A collection of TFunction objects designed for fast access given a DeclId_t and for keep track of TFu...
TSeqCollection * fSockets
Definition: TROOT.h:157
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:953
TObject * FindSpecialObject(const char *name, void *&where)
Returns address and folder of a ROOT object if it exists.
Definition: TROOT.cxx:1203
R__EXTERN TInterpreter * gCling
Definition: TInterpreter.h:519
R__DLLEXPORT TInterpreter * CreateInterpreter(void *interpLibHandle)
Definition: TCling.cxx:590
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:416
virtual Int_t GetSize() const
Definition: TCollection.h:95
Bool_t fInterrupt
Definition: TROOT.h:143
Abstract base class defining a generic interface to the underlying Operating System.
Definition: TSystem.h:258
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
void Idle(UInt_t idleTimeInSec, const char *command=0)
Execute command when system has been idle for idleTimeInSec seconds.
Definition: TROOT.cxx:1699
void *(* GlobalFunc_t)()
Definition: TGlobal.h:56
TSeqCollection * GetListOfBrowsers() const
Definition: TROOT.h:253
virtual char * ConcatFileName(const char *dir, const char *name)
Concatenate a directory and a file name. User must delete returned string.
Definition: TSystem.cxx:1045
#define sym(otri1, otri2)
Definition: triangle.c:932
static TClass * R__GetClassIfKnown(const char *className)
Check whether className is a known class, and only autoload if we can.
Definition: TROOT.cxx:1714
char name[80]
Definition: TGX11.cxx:109
TString fConfigFeatures
Definition: TROOT.h:123
TProcessUUID * fUUIDs
Definition: TROOT.h:175
void DisableParBranchProcessing()
Globally disables the IMT use case of parallel branch processing, deactivating the corresponding lock...
Definition: TROOT.cxx:409
static void PrintStatistics()
Print memory usage statistics.
Definition: TStorage.cxx:406
const TObject * fPrimitive
Definition: TROOT.h:147
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:911
TVirtualPad * fSelectPad
Definition: TROOT.h:148
virtual TObject * First() const =0
virtual TObject * Last() const =0
TDictionary * Get(DeclId_t id)
Return (after creating it if necessary) the TDataMember describing the data member corresponding to t...
TROOT()
Default ctor.
Definition: TROOT.cxx:574