Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 <ROOT/RConfig.hxx>
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 FARPROC dlsym(void *library, const char *function_name)
94 {
95  HMODULE hMods[1024];
96  DWORD cbNeeded;
97  FARPROC address = NULL;
98  unsigned int i;
99  if (library == RTLD_DEFAULT) {
100  if (EnumProcessModules(::GetCurrentProcess(), hMods, sizeof(hMods), &cbNeeded)) {
101  for (i = 0; i < (cbNeeded / sizeof(HMODULE)); i++) {
102  address = ::GetProcAddress((HMODULE)hMods[i], function_name);
103  if (address)
104  return address;
105  }
106  }
107  return address;
108  } else {
109  return ::GetProcAddress((HMODULE)library, function_name);
110  }
111 }
112 #else
113 #include <dlfcn.h>
114 #endif
115 
116 #include "Riostream.h"
117 #include "ROOT/FoundationUtils.hxx"
118 #include "TROOT.h"
119 #include "TClass.h"
120 #include "TClassEdit.h"
121 #include "TClassGenerator.h"
122 #include "TDataType.h"
123 #include "TDatime.h"
124 #include "TStyle.h"
125 #include "TObjectTable.h"
126 #include "TClassTable.h"
127 #include "TSystem.h"
128 #include "THashList.h"
129 #include "TObjArray.h"
130 #include "TEnv.h"
131 #include "TError.h"
132 #include "TColor.h"
133 #include "TGlobal.h"
134 #include "TFunction.h"
135 #include "TVirtualPad.h"
136 #include "TBrowser.h"
137 #include "TSystemDirectory.h"
138 #include "TApplication.h"
139 #include "TInterpreter.h"
140 #include "TGuiFactory.h"
141 #include "TMessageHandler.h"
142 #include "TFolder.h"
143 #include "TQObject.h"
144 #include "TProcessUUID.h"
145 #include "TPluginManager.h"
146 #include "TMap.h"
147 #include "TObjString.h"
148 #include "TVirtualMutex.h"
149 #include "TInterpreter.h"
150 #include "TListOfTypes.h"
151 #include "TListOfDataMembers.h"
152 #include "TListOfEnumsWithLock.h"
153 #include "TListOfFunctions.h"
155 #include "TFunctionTemplate.h"
156 #include "ThreadLocalStorage.h"
157 #include "TVirtualRWMutex.h"
158 
159 #include <string>
160 namespace std {} using namespace std;
161 
162 #if defined(R__UNIX)
163 #if defined(R__HAS_COCOA)
164 #include "TMacOSXSystem.h"
165 #include "TUrl.h"
166 #else
167 #include "TUnixSystem.h"
168 #endif
169 #elif defined(R__WIN32)
170 #include "TWinNTSystem.h"
171 #endif
172 
173 extern "C" void R__SetZipMode(int);
174 
175 static DestroyInterpreter_t *gDestroyInterpreter = 0;
176 static void *gInterpreterLib = 0;
177 
178 // Mutex for protection of concurrent gROOT access
179 TVirtualMutex* gROOTMutex = 0;
180 ROOT::TVirtualRWMutex *ROOT::gCoreMutex = nullptr;
181 
182 // For accessing TThread::Tsd indirectly.
183 void **(*gThreadTsd)(void*,Int_t) = 0;
184 
185 //-------- Names of next three routines are a small homage to CMZ --------------
186 ////////////////////////////////////////////////////////////////////////////////
187 /// Return version id as an integer, i.e. "2.22/04" -> 22204.
188 
189 static Int_t IVERSQ()
190 {
191  Int_t maj, min, cycle;
192  sscanf(ROOT_RELEASE, "%d.%d/%d", &maj, &min, &cycle);
193  return 10000*maj + 100*min + cycle;
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// Return built date as integer, i.e. "Apr 28 2000" -> 20000428.
198 
199 static Int_t IDATQQ(const char *date)
200 {
201  static const char *months[] = {"Jan","Feb","Mar","Apr","May",
202  "Jun","Jul","Aug","Sep","Oct",
203  "Nov","Dec"};
204 
205  char sm[12];
206  Int_t yy, mm=0, dd;
207  sscanf(date, "%s %d %d", sm, &dd, &yy);
208  for (int i = 0; i < 12; i++)
209  if (!strncmp(sm, months[i], 3)) {
210  mm = i+1;
211  break;
212  }
213  return 10000*yy + 100*mm + dd;
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Return built time as integer (with min precision), i.e.
218 /// "17:32:37" -> 1732.
219 
220 static Int_t ITIMQQ(const char *time)
221 {
222  Int_t hh, mm, ss;
223  sscanf(time, "%d:%d:%d", &hh, &mm, &ss);
224  return 100*hh + mm;
225 }
226 
227 ////////////////////////////////////////////////////////////////////////////////
228 /// Clean up at program termination before global objects go out of scope.
229 
230 static void CleanUpROOTAtExit()
231 {
232  if (gROOT) {
233  R__LOCKGUARD(gROOTMutex);
234 
235  if (gROOT->GetListOfFiles())
236  gROOT->GetListOfFiles()->Delete("slow");
237  if (gROOT->GetListOfSockets())
238  gROOT->GetListOfSockets()->Delete();
239  if (gROOT->GetListOfMappedFiles())
240  gROOT->GetListOfMappedFiles()->Delete("slow");
241  if (gROOT->GetListOfClosedObjects())
242  gROOT->GetListOfClosedObjects()->Delete("slow");
243  }
244 }
245 
246 ////////////////////////////////////////////////////////////////////////////////
247 /// A module and its headers. Intentionally not a copy:
248 /// If these strings end up in this struct they are
249 /// long lived by definition because they get passed in
250 /// before initialization of TCling.
251 
252 namespace {
253  struct ModuleHeaderInfo_t {
254  ModuleHeaderInfo_t(const char* moduleName,
255  const char** headers,
256  const char** includePaths,
257  const char* payloadCode,
258  const char* fwdDeclCode,
259  void (*triggerFunc)(),
260  const TROOT::FwdDeclArgsToKeepCollection_t& fwdDeclsArgToSkip,
261  const char **classesHeaders,
262  bool hasCxxModule):
263  fModuleName(moduleName),
264  fHeaders(headers),
265  fPayloadCode(payloadCode),
266  fFwdDeclCode(fwdDeclCode),
267  fIncludePaths(includePaths),
268  fTriggerFunc(triggerFunc),
269  fClassesHeaders(classesHeaders),
270  fFwdNargsToKeepColl(fwdDeclsArgToSkip),
271  fHasCxxModule(hasCxxModule) {}
272 
273  const char* fModuleName; // module name
274  const char** fHeaders; // 0-terminated array of header files
275  const char* fPayloadCode; // Additional code to be given to cling at library load
276  const char* fFwdDeclCode; // Additional code to let cling know about selected classes and functions
277  const char** fIncludePaths; // 0-terminated array of header files
278  void (*fTriggerFunc)(); // Pointer to the dict initialization used to find the library name
279  const char** fClassesHeaders; // 0-terminated list of classes and related header files
280  const TROOT::FwdDeclArgsToKeepCollection_t fFwdNargsToKeepColl; // Collection of
281  // pairs of template fwd decls and number of
282  bool fHasCxxModule; // Whether this module has a C++ module alongside it.
283  };
284 
285  std::vector<ModuleHeaderInfo_t>& GetModuleHeaderInfoBuffer() {
286  static std::vector<ModuleHeaderInfo_t> moduleHeaderInfoBuffer;
287  return moduleHeaderInfoBuffer;
288  }
289 }
290 
291 Int_t TROOT::fgDirLevel = 0;
292 Bool_t TROOT::fgRootInit = kFALSE;
293 Bool_t TROOT::fgMemCheck = kFALSE;
294 
295 static void at_exit_of_TROOT() {
296  if (ROOT::Internal::gROOTLocal)
297  ROOT::Internal::gROOTLocal->~TROOT();
298 }
299 
300 // This local static object initializes the ROOT system
301 namespace ROOT {
302 namespace Internal {
303  class TROOTAllocator {
304  // Simple wrapper to separate, time-wise, the call to the
305  // TROOT destructor and the actual free-ing of the memory.
306  //
307  // Since the interpreter implementation (currently TCling) is
308  // loaded via dlopen by libCore, the destruction of its global
309  // variable (i.e. in particular clang's) is scheduled before
310  // those in libCore so we need to schedule the call to the TROOT
311  // destructor before that *but* we want to make sure the memory
312  // stay around until libCore itself is unloaded so that code
313  // using gROOT can 'properly' check for validity.
314  //
315  // The order of loading for is:
316  // libCore.so
317  // libRint.so
318  // ... anything other library hard linked to the executable ...
319  // ... for example libEvent
320  // libCling.so
321  // ... other libraries like libTree for example ....
322  // and the destruction order is (of course) the reverse.
323  // By default the unloading of the dictionary, does use
324  // the service of the interpreter ... which of course
325  // fails if libCling is already unloaded by that information
326  // has not been registered per se.
327  //
328  // To solve this problem, we now schedule the destruction
329  // of the TROOT object to happen _just_ before the
330  // unloading/destruction of libCling so that we can
331  // maximize the amount of clean-up we can do correctly
332  // and we can still allocate the TROOT object's memory
333  // statically.
334  //
335  char fHolder[sizeof(TROOT)];
336  public:
337  TROOTAllocator() {
338  new(&(fHolder[0])) TROOT("root", "The ROOT of EVERYTHING");
339  }
340 
341  ~TROOTAllocator() {
342  if (gROOTLocal) {
343  gROOTLocal->~TROOT();
344  }
345  }
346  };
347 
348  // The global gROOT is defined to be a function (ROOT::GetROOT())
349  // which itself is dereferencing a function pointer.
350 
351  // Initially this function pointer's value is & GetROOT1 whose role is to
352  // create and initialize the TROOT object itself.
353  // At the very end of the TROOT constructor the value of the function pointer
354  // is switch to & GetROOT2 whose role is to initialize the interpreter.
355 
356  // This mechanism was primarily intended to fix the issues with order in which
357  // global TROOT and LLVM globals are initialized. TROOT was initializing
358  // Cling, but Cling could not be used yet due to LLVM globals not being
359  // Initialized yet. The solution is to delay initializing the interpreter in
360  // TROOT till after main() when all LLVM globals are initialized.
361 
362  // Technically, the mechanism used actually delay the interpreter
363  // initialization until the first use of gROOT *after* the end of the
364  // TROOT constructor.
365 
366  // So to delay until after the start of main, we also made sure that none
367  // of the ROOT code (mostly the dictionary code) used during library loading
368  // is using gROOT (directly or indirectly).
369 
370  // In practice, the initialization of the interpreter is now delayed until
371  // the first use gROOT (or gInterpreter) after the start of main (but user
372  // could easily break this by using gROOT in their library initialization
373  // code).
374 
375  extern TROOT *gROOTLocal;
376 
377  TROOT *GetROOT1() {
378  if (gROOTLocal)
379  return gROOTLocal;
380  static TROOTAllocator alloc;
381  return gROOTLocal;
382  }
383 
384  TROOT *GetROOT2() {
385  static Bool_t initInterpreter = kFALSE;
386  if (!initInterpreter) {
387  initInterpreter = kTRUE;
388  gROOTLocal->InitInterpreter();
389  // Load and init threads library
390  gROOTLocal->InitThreads();
391  }
392  return gROOTLocal;
393  }
394  typedef TROOT *(*GetROOTFun_t)();
395 
396  static GetROOTFun_t gGetROOT = &GetROOT1;
397 
398  static Func_t GetSymInLibImt(const char *funcname)
399  {
400  const static bool loadSuccess = dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")? false : 0 <= gSystem->Load("libImt");
401  if (loadSuccess) {
402  if (auto sym = gSystem->DynFindSymbol(nullptr, funcname)) {
403  return sym;
404  } else {
405  Error("GetSymInLibImt", "Cannot get symbol %s.", funcname);
406  }
407  }
408  return nullptr;
409  }
410 
411  //////////////////////////////////////////////////////////////////////////////
412  /// Globally enables the parallel branch processing, which is a case of
413  /// implicit multi-threading (IMT) in ROOT, activating the required locks.
414  /// This IMT use case, implemented in TTree::GetEntry, spawns a task for
415  /// each branch of the tree. Therefore, a task takes care of the reading,
416  /// decompression and deserialisation of a given branch.
417  void EnableParBranchProcessing()
418  {
419 #ifdef R__USE_IMT
420  if (!IsImplicitMTEnabled())
421  EnableImplicitMT();
422  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_EnableParBranchProcessing");
423  if (sym)
424  sym();
425 #else
426  ::Warning("EnableParBranchProcessing", "Cannot enable parallel branch processing, please build ROOT with -Dimt=ON");
427 #endif
428  }
429 
430  //////////////////////////////////////////////////////////////////////////////
431  /// Globally disables the IMT use case of parallel branch processing,
432  /// deactivating the corresponding locks.
433  void DisableParBranchProcessing()
434  {
435 #ifdef R__USE_IMT
436  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_DisableParBranchProcessing");
437  if (sym)
438  sym();
439 #else
440  ::Warning("DisableParBranchProcessing", "Cannot disable parallel branch processing, please build ROOT with -Dimt=ON");
441 #endif
442  }
443 
444  //////////////////////////////////////////////////////////////////////////////
445  /// Returns true if parallel branch processing is enabled.
446  Bool_t IsParBranchProcessingEnabled()
447  {
448 #ifdef R__USE_IMT
449  static Bool_t (*sym)() = (Bool_t(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_IsParBranchProcessingEnabled");
450  if (sym)
451  return sym();
452  else
453  return kFALSE;
454 #else
455  return kFALSE;
456 #endif
457  }
458 
459  ////////////////////////////////////////////////////////////////////////////////
460  /// Globally enables the parallel tree processing, which is a case of
461  /// implicit multi-threading in ROOT, activating the required locks.
462  /// This IMT use case, implemented in TTreeProcessor::Process, receives a user
463  /// function and applies it to subranges of the tree, which correspond to its
464  /// clusters. Hence, for every cluster, a task is spawned to potentially
465  /// process it in parallel with the other clusters.
466  void EnableParTreeProcessing()
467  {
468 #ifdef R__USE_IMT
469  if (!IsImplicitMTEnabled())
470  EnableImplicitMT();
471  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_EnableParTreeProcessing");
472  if (sym)
473  sym();
474 #else
475  ::Warning("EnableParTreeProcessing", "Cannot enable parallel tree processing, please build ROOT with -Dimt=ON");
476 #endif
477  }
478 
479  //////////////////////////////////////////////////////////////////////////////
480  /// Globally disables the IMT use case of parallel branch processing,
481  /// deactivating the corresponding locks.
482  void DisableParTreeProcessing()
483  {
484 #ifdef R__USE_IMT
485  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_DisableParTreeProcessing");
486  if (sym)
487  sym();
488 #else
489  ::Warning("DisableParTreeProcessing", "Cannot disable parallel tree processing, please build ROOT with -Dimt=ON");
490 #endif
491  }
492 
493  ////////////////////////////////////////////////////////////////////////////////
494  /// Returns true if parallel tree processing is enabled.
495  Bool_t IsParTreeProcessingEnabled()
496  {
497 #ifdef R__USE_IMT
498  static Bool_t (*sym)() = (Bool_t(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_IsParTreeProcessingEnabled");
499  if (sym)
500  return sym();
501  else
502  return kFALSE;
503 #else
504  return kFALSE;
505 #endif
506  }
507 
508  ////////////////////////////////////////////////////////////////////////////////
509  /// Keeps track of the status of ImplicitMT w/o resorting to the load of
510  /// libImt
511  static Bool_t &IsImplicitMTEnabledImpl()
512  {
513  static Bool_t isImplicitMTEnabled = kFALSE;
514  return isImplicitMTEnabled;
515  }
516 
517 } // end of Internal sub namespace
518 // back to ROOT namespace
519 
520  TROOT *GetROOT() {
521  return (*Internal::gGetROOT)();
522  }
523 
524  TString &GetMacroPath() {
525  static TString macroPath;
526  return macroPath;
527  }
528 
529  // clang-format off
530  ////////////////////////////////////////////////////////////////////////////////
531  /// Enables the global mutex to make ROOT thread safe/aware.
532  ///
533  /// The following becomes safe:
534  /// - concurrent construction and destruction of TObjects, including the ones registered in ROOT's global lists (e.g. gROOT->GetListOfCleanups(), gROOT->GetListOfFiles())
535  /// - concurrent usage of _different_ ROOT objects from different threads, including ones with global state (e.g. TFile, TTree, TChain) with the exception of graphics classes (e.g. TCanvas)
536  /// - concurrent calls to ROOT's type system classes, e.g. TClass and TEnum
537  /// - concurrent calls to the interpreter through gInterpreter
538  /// - concurrent loading of ROOT plug-ins
539  ///
540  /// In addition, gDirectory, gFile and gPad become a thread-local variable.
541  /// In all threads, gDirectory defaults to gROOT, a singleton which supports thread-safe insertion and deletion of contents.
542  /// gFile and gPad default to nullptr, as it is for single-thread programs.
543  ///
544  /// The ROOT graphics subsystem is not made thread-safe by this method. In particular drawing or printing different
545  /// canvases from different threads (and analogous operations such as invoking `Draw` on a `TObject`) is not thread-safe.
546  ///
547  /// Note that there is no `DisableThreadSafety()`. ROOT's thread-safety features cannot be disabled once activated.
548  // clang-format on
549  void EnableThreadSafety()
550  {
551  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TThread_Initialize");
552  if (sym)
553  sym();
554  }
555 
556  ////////////////////////////////////////////////////////////////////////////////
557  /// @param[in] numthreads Number of threads to use. If not specified or
558  /// set to zero, the number of threads is automatically
559  /// decided by the implementation. Any other value is
560  /// used as a hint.
561  ///
562  /// ROOT must be built with the compilation flag `imt=ON` for this feature to be available.
563  /// The following objects and methods automatically take advantage of
564  /// multi-threading if a call to `EnableImplicitMT` has been made before usage:
565  ///
566  /// - RDataFrame internally runs the event-loop by parallelizing over clusters of entries
567  /// - TTree::GetEntry reads multiple branches in parallel
568  /// - TTree::FlushBaskets writes multiple baskets to disk in parallel
569  /// - TTreeCacheUnzip decompresses the baskets contained in a TTreeCache in parallel
570  /// - THx::Fit performs in parallel the evaluation of the objective function over the data
571  /// - TMVA::DNN trains the deep neural networks in parallel
572  /// - TMVA::BDT trains the classifier in parallel and multiclass BDTs are evaluated in parallel
573  ///
574  /// EnableImplicitMT calls in turn EnableThreadSafety.
575  /// The 'numthreads' parameter allows to control the number of threads to
576  /// be used by the implicit multi-threading. However, this parameter is just
577  /// a hint for ROOT: it will try to satisfy the request if the execution
578  /// scenario allows it. For example, if ROOT is configured to use an external
579  /// scheduler, setting a value for 'numthreads' might not have any effect.
580  void EnableImplicitMT(UInt_t numthreads)
581  {
582 #ifdef R__USE_IMT
583  if (ROOT::Internal::IsImplicitMTEnabledImpl())
584  return;
585  EnableThreadSafety();
586  static void (*sym)(UInt_t) = (void(*)(UInt_t))Internal::GetSymInLibImt("ROOT_TImplicitMT_EnableImplicitMT");
587  if (sym)
588  sym(numthreads);
589  ROOT::Internal::IsImplicitMTEnabledImpl() = true;
590 #else
591  ::Warning("EnableImplicitMT", "Cannot enable implicit multi-threading with %d threads, please build ROOT with -Dimt=ON", numthreads);
592 #endif
593  }
594 
595  ////////////////////////////////////////////////////////////////////////////////
596  /// Disables the implicit multi-threading in ROOT (see EnableImplicitMT).
597  void DisableImplicitMT()
598  {
599 #ifdef R__USE_IMT
600  static void (*sym)() = (void(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_DisableImplicitMT");
601  if (sym)
602  sym();
603  ROOT::Internal::IsImplicitMTEnabledImpl() = kFALSE;
604 #else
605  ::Warning("DisableImplicitMT", "Cannot disable implicit multi-threading, please build ROOT with -Dimt=ON");
606 #endif
607  }
608 
609  ////////////////////////////////////////////////////////////////////////////////
610  /// Returns true if the implicit multi-threading in ROOT is enabled.
611  Bool_t IsImplicitMTEnabled()
612  {
613  return ROOT::Internal::IsImplicitMTEnabledImpl();
614  }
615 
616  ////////////////////////////////////////////////////////////////////////////////
617  /// Returns the size of the pool used for implicit multi-threading.
618  UInt_t GetImplicitMTPoolSize()
619  {
620 #ifdef R__USE_IMT
621  static UInt_t (*sym)() = (UInt_t(*)())Internal::GetSymInLibImt("ROOT_TImplicitMT_GetImplicitMTPoolSize");
622  if (sym)
623  return sym();
624  else
625  return 0;
626 #else
627  return 0;
628 #endif
629  }
630 
631 }
632 
633 TROOT *ROOT::Internal::gROOTLocal = ROOT::GetROOT();
634 
635 // Global debug flag (set to > 0 to get debug output).
636 // Can be set either via the interpreter (gDebug is exported to CINT),
637 // via the rootrc resource "Root.Debug", via the shell environment variable
638 // ROOTDEBUG, or via the debugger.
639 Int_t gDebug;
640 
641 
642 ClassImp(TROOT);
643 
644 ////////////////////////////////////////////////////////////////////////////////
645 /// Default ctor.
646 
647 TROOT::TROOT() : TDirectory(),
648  fLineIsProcessing(0), fVersion(0), fVersionInt(0), fVersionCode(0),
649  fVersionDate(0), fVersionTime(0), fBuiltDate(0), fBuiltTime(0),
650  fTimer(0), fApplication(0), fInterpreter(0), fBatch(kTRUE),
651  fIsWebDisplay(kFALSE), fIsWebDisplayBatch(kFALSE), fEditHistograms(kTRUE),
652  fFromPopUp(kTRUE),fMustClean(kTRUE),fReadingObject(kFALSE),fForceStyle(kFALSE),
653  fInterrupt(kFALSE),fEscape(kFALSE),fExecutingMacro(kFALSE),fEditorMode(0),
654  fPrimitive(0),fSelectPad(0),fClasses(0),fTypes(0),fGlobals(0),fGlobalFunctions(0),
655  fClosedObjects(0),fFiles(0),fMappedFiles(0),fSockets(0),fCanvases(0),fStyles(0),fFunctions(0),
656  fTasks(0),fColors(0),fGeometries(0),fBrowsers(0),fSpecials(0),fCleanups(0),
657  fMessageHandlers(0),fStreamerInfo(0),fClassGenerators(0),fSecContexts(0),
658  fProofs(0),fClipboard(0),fDataSets(0),fUUIDs(0),fRootFolder(0),fBrowsables(0),
659  fPluginManager(0)
660 {
661 }
662 
663 ////////////////////////////////////////////////////////////////////////////////
664 /// Initialize the ROOT system. The creation of the TROOT object initializes
665 /// the ROOT system. It must be the first ROOT related action that is
666 /// performed by a program. The TROOT object must be created on the stack
667 /// (can not be called via new since "operator new" is protected). The
668 /// TROOT object is either created as a global object (outside the main()
669 /// program), or it is one of the first objects created in main().
670 /// Make sure that the TROOT object stays in scope for as long as ROOT
671 /// related actions are performed. TROOT is a so called singleton so
672 /// only one instance of it can be created. The single TROOT object can
673 /// always be accessed via the global pointer gROOT.
674 /// The name and title arguments can be used to identify the running
675 /// application. The initfunc argument can contain an array of
676 /// function pointers (last element must be 0). These functions are
677 /// executed at the end of the constructor. This way one can easily
678 /// extend the ROOT system without adding permanent dependencies
679 /// (e.g. the graphics system is initialized via such a function).
680 
681 TROOT::TROOT(const char *name, const char *title, VoidFuncPtr_t *initfunc)
682  : TDirectory(), fLineIsProcessing(0), fVersion(0), fVersionInt(0), fVersionCode(0),
683  fVersionDate(0), fVersionTime(0), fBuiltDate(0), fBuiltTime(0),
684  fTimer(0), fApplication(0), fInterpreter(0), fBatch(kTRUE),
685  fIsWebDisplay(kFALSE), fIsWebDisplayBatch(kFALSE), fEditHistograms(kTRUE),
686  fFromPopUp(kTRUE),fMustClean(kTRUE),fReadingObject(kFALSE),fForceStyle(kFALSE),
687  fInterrupt(kFALSE),fEscape(kFALSE),fExecutingMacro(kFALSE),fEditorMode(0),
688  fPrimitive(0),fSelectPad(0),fClasses(0),fTypes(0),fGlobals(0),fGlobalFunctions(0),
689  fClosedObjects(0),fFiles(0),fMappedFiles(0),fSockets(0),fCanvases(0),fStyles(0),fFunctions(0),
690  fTasks(0),fColors(0),fGeometries(0),fBrowsers(0),fSpecials(0),fCleanups(0),
691  fMessageHandlers(0),fStreamerInfo(0),fClassGenerators(0),fSecContexts(0),
692  fProofs(0),fClipboard(0),fDataSets(0),fUUIDs(0),fRootFolder(0),fBrowsables(0),
693  fPluginManager(0)
694 {
695  if (fgRootInit || ROOT::Internal::gROOTLocal) {
696  //Warning("TROOT", "only one instance of TROOT allowed");
697  return;
698  }
699 
700  R__LOCKGUARD(gROOTMutex);
701 
702  ROOT::Internal::gROOTLocal = this;
703  gDirectory = 0;
704 
705  SetName(name);
706  SetTitle(title);
707 
708  // will be used by global "operator delete" so make sure it is set
709  // before anything is deleted
710  fMappedFiles = 0;
711 
712  // create already here, but only initialize it after gEnv has been created
713  gPluginMgr = fPluginManager = new TPluginManager;
714 
715  // Initialize Operating System interface
716  InitSystem();
717 
718  // Initialize static directory functions
719  GetRootSys();
720  GetBinDir();
721  GetLibDir();
722  GetIncludeDir();
723  GetEtcDir();
724  GetDataDir();
725  GetDocDir();
726  GetMacroDir();
727  GetTutorialDir();
728  GetSourceDir();
729  GetIconPath();
730  GetTTFFontDir();
731 
732  gRootDir = GetRootSys().Data();
733 
734  TDirectory::BuildDirectory(nullptr, nullptr);
735 
736  // Initialize interface to CINT C++ interpreter
737  fVersionInt = 0; // check in TROOT dtor in case TCling fails
738  fClasses = 0; // might be checked via TCling ctor
739  fEnums = 0;
740 
741  fConfigOptions = R__CONFIGUREOPTIONS;
742  fConfigFeatures = R__CONFIGUREFEATURES;
743  fVersion = ROOT_RELEASE;
744  fVersionCode = ROOT_VERSION_CODE;
745  fVersionInt = IVERSQ();
746  fVersionDate = IDATQQ(ROOT_RELEASE_DATE);
747  fVersionTime = ITIMQQ(ROOT_RELEASE_TIME);
748  fBuiltDate = IDATQQ(__DATE__);
749  fBuiltTime = ITIMQQ(__TIME__);
750 
751  ReadGitInfo();
752 
753  fClasses = new THashTable(800,3); fClasses->UseRWLock();
754  //fIdMap = new IdMap_t;
755  fStreamerInfo = new TObjArray(100); fStreamerInfo->UseRWLock();
756  fClassGenerators = new TList;
757 
758  // usedToIdentifyRootClingByDlSym is available when TROOT is part of
759  // rootcling.
760  if (!dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")) {
761  // initialize plugin manager early
762  fPluginManager->LoadHandlersFromEnv(gEnv);
763  }
764 
765  TSystemDirectory *workdir = new TSystemDirectory("workdir", gSystem->WorkingDirectory());
766 
767  auto setNameLocked = [](TSeqCollection *l, const char *collection_name) {
768  l->SetName(collection_name);
769  l->UseRWLock();
770  return l;
771  };
772 
773  fTimer = 0;
774  fApplication = 0;
775  fColors = setNameLocked(new TObjArray(1000), "ListOfColors");
776  fTypes = 0;
777  fGlobals = 0;
778  fGlobalFunctions = 0;
779  // fList was created in TDirectory::Build but with different sizing.
780  delete fList;
781  fList = new THashList(1000,3); fList->UseRWLock();
782  fClosedObjects = setNameLocked(new TList, "ClosedFiles");
783  fFiles = setNameLocked(new TList, "Files");
784  fMappedFiles = setNameLocked(new TList, "MappedFiles");
785  fSockets = setNameLocked(new TList, "Sockets");
786  fCanvases = setNameLocked(new TList, "Canvases");
787  fStyles = setNameLocked(new TList, "Styles");
788  fFunctions = setNameLocked(new TList, "Functions");
789  fTasks = setNameLocked(new TList, "Tasks");
790  fGeometries = setNameLocked(new TList, "Geometries");
791  fBrowsers = setNameLocked(new TList, "Browsers");
792  fSpecials = setNameLocked(new TList, "Specials");
793  fBrowsables = (TList*)setNameLocked(new TList, "Browsables");
794  fCleanups = setNameLocked(new THashList, "Cleanups");
795  fMessageHandlers = setNameLocked(new TList, "MessageHandlers");
796  fSecContexts = setNameLocked(new TList, "SecContexts");
797  fProofs = setNameLocked(new TList, "Proofs");
798  fClipboard = setNameLocked(new TList, "Clipboard");
799  fDataSets = setNameLocked(new TList, "DataSets");
800  fTypes = new TListOfTypes; fTypes->UseRWLock();
801 
802  TProcessID::AddProcessID();
803  fUUIDs = new TProcessUUID();
804 
805  fRootFolder = new TFolder();
806  fRootFolder->SetName("root");
807  fRootFolder->SetTitle("root of all folders");
808  fRootFolder->AddFolder("Classes", "List of Active Classes",fClasses);
809  fRootFolder->AddFolder("Colors", "List of Active Colors",fColors);
810  fRootFolder->AddFolder("MapFiles", "List of MapFiles",fMappedFiles);
811  fRootFolder->AddFolder("Sockets", "List of Socket Connections",fSockets);
812  fRootFolder->AddFolder("Canvases", "List of Canvases",fCanvases);
813  fRootFolder->AddFolder("Styles", "List of Styles",fStyles);
814  fRootFolder->AddFolder("Functions", "List of Functions",fFunctions);
815  fRootFolder->AddFolder("Tasks", "List of Tasks",fTasks);
816  fRootFolder->AddFolder("Geometries","List of Geometries",fGeometries);
817  fRootFolder->AddFolder("Browsers", "List of Browsers",fBrowsers);
818  fRootFolder->AddFolder("Specials", "List of Special Objects",fSpecials);
819  fRootFolder->AddFolder("Handlers", "List of Message Handlers",fMessageHandlers);
820  fRootFolder->AddFolder("Cleanups", "List of RecursiveRemove Collections",fCleanups);
821  fRootFolder->AddFolder("StreamerInfo","List of Active StreamerInfo Classes",fStreamerInfo);
822  fRootFolder->AddFolder("SecContexts","List of Security Contexts",fSecContexts);
823  fRootFolder->AddFolder("PROOF Sessions", "List of PROOF sessions",fProofs);
824  fRootFolder->AddFolder("ROOT Memory","List of Objects in the gROOT Directory",fList);
825  fRootFolder->AddFolder("ROOT Files","List of Connected ROOT Files",fFiles);
826 
827  // by default, add the list of files, tasks, canvases and browsers in the Cleanups list
828  fCleanups->Add(fCanvases); fCanvases->SetBit(kMustCleanup);
829  fCleanups->Add(fBrowsers); fBrowsers->SetBit(kMustCleanup);
830  fCleanups->Add(fTasks); fTasks->SetBit(kMustCleanup);
831  fCleanups->Add(fFiles); fFiles->SetBit(kMustCleanup);
832  fCleanups->Add(fClosedObjects); fClosedObjects->SetBit(kMustCleanup);
833  // And add TROOT's TDirectory personality
834  fCleanups->Add(fList);
835 
836  fExecutingMacro= kFALSE;
837  fForceStyle = kFALSE;
838  fFromPopUp = kFALSE;
839  fReadingObject = kFALSE;
840  fInterrupt = kFALSE;
841  fEscape = kFALSE;
842  fMustClean = kTRUE;
843  fPrimitive = 0;
844  fSelectPad = 0;
845  fEditorMode = 0;
846  fDefCanvasName = "c1";
847  fEditHistograms= kFALSE;
848  fLineIsProcessing = 1; // This prevents WIN32 "Windows" thread to pick ROOT objects with mouse
849  gDirectory = this;
850  gPad = 0;
851 
852  //set name of graphical cut class for the graphics editor
853  //cannot call SetCutClassName at this point because the TClass of TCutG
854  //is not yet build
855  fCutClassName = "TCutG";
856 
857  // Create a default MessageHandler
858  new TMessageHandler((TClass*)0);
859 
860  // Create some styles
861  gStyle = 0;
862  TStyle::BuildStyles();
863  SetStyle(gEnv->GetValue("Canvas.Style", "Modern"));
864 
865  // Setup default (batch) graphics and GUI environment
866  gBatchGuiFactory = new TGuiFactory;
867  gGuiFactory = gBatchGuiFactory;
868  gGXBatch = new TVirtualX("Batch", "ROOT Interface to batch graphics");
869  gVirtualX = gGXBatch;
870 
871 #if defined(R__WIN32)
872  fBatch = kFALSE;
873 #elif defined(R__HAS_COCOA)
874  fBatch = kFALSE;
875 #else
876  if (gSystem->Getenv("DISPLAY"))
877  fBatch = kFALSE;
878  else
879  fBatch = kTRUE;
880 #endif
881 
882  int i = 0;
883  while (initfunc && initfunc[i]) {
884  (initfunc[i])();
885  fBatch = kFALSE; // put system in graphics mode (backward compatible)
886  i++;
887  }
888 
889  // Set initial/default list of browsable objects
890  fBrowsables->Add(fRootFolder, "root");
891  fBrowsables->Add(fProofs, "PROOF Sessions");
892  fBrowsables->Add(workdir, gSystem->WorkingDirectory());
893  fBrowsables->Add(fFiles, "ROOT Files");
894 
895  atexit(CleanUpROOTAtExit);
896 
897  ROOT::Internal::gGetROOT = &ROOT::Internal::GetROOT2;
898 }
899 
900 ////////////////////////////////////////////////////////////////////////////////
901 /// Clean up and free resources used by ROOT (files, network sockets,
902 /// shared memory segments, etc.).
903 
904 TROOT::~TROOT()
905 {
906  using namespace ROOT::Internal;
907 
908  if (gROOTLocal == this) {
909 
910  // If the interpreter has not yet been initialized, don't bother
911  gGetROOT = &GetROOT1;
912 
913  // Mark the object as invalid, so that we can veto some actions
914  // (like autoloading) while we are in the destructor.
915  SetBit(TObject::kInvalidObject);
916 
917  // Turn-off the global mutex to avoid recreating mutexes that have
918  // already been deleted during the destruction phase
919  gGlobalMutex = 0;
920 
921  // Return when error occurred in TCling, i.e. when setup file(s) are
922  // out of date
923  if (!fVersionInt) return;
924 
925  // ATTENTION!!! Order is important!
926 
927  SafeDelete(fBrowsables);
928 
929  // FIXME: Causes rootcling to deadlock, debug and uncomment
930  // SafeDelete(fRootFolder);
931 
932 #ifdef R__COMPLETE_MEM_TERMINATION
933  fSpecials->Delete(); SafeDelete(fSpecials); // delete special objects : PostScript, Minuit, Html
934 #endif
935 
936  fClosedObjects->Delete("slow"); // and closed files
937  fFiles->Delete("slow"); // and files
938  SafeDelete(fFiles);
939  fSecContexts->Delete("slow"); SafeDelete(fSecContexts); // and security contexts
940  fSockets->Delete(); SafeDelete(fSockets); // and sockets
941  fMappedFiles->Delete("slow"); // and mapped files
942  TSeqCollection *tl = fMappedFiles; fMappedFiles = 0; delete tl;
943 
944  SafeDelete(fClosedObjects);
945 
946  delete fUUIDs;
947  TProcessID::Cleanup(); // and list of ProcessIDs
948 
949  fFunctions->Delete(); SafeDelete(fFunctions); // etc..
950  fGeometries->Delete(); SafeDelete(fGeometries);
951  fBrowsers->Delete(); SafeDelete(fBrowsers);
952  SafeDelete(fCanvases);
953  fColors->Delete(); SafeDelete(fColors);
954  fStyles->Delete(); SafeDelete(fStyles);
955 
956 #ifdef R__COMPLETE_MEM_TERMINATION
957  if (gGuiFactory != gBatchGuiFactory) SafeDelete(gGuiFactory);
958  SafeDelete(gBatchGuiFactory);
959  if (gGXBatch != gVirtualX) SafeDelete(gGXBatch);
960  SafeDelete(gVirtualX);
961 #endif
962 
963  // Stop emitting signals
964  TQObject::BlockAllSignals(kTRUE);
965 
966  fMessageHandlers->Delete(); SafeDelete(fMessageHandlers);
967 
968 #ifdef R__COMPLETE_MEM_TERMINATION
969  SafeDelete(fCanvases);
970  SafeDelete(fTasks);
971  SafeDelete(fProofs);
972  SafeDelete(fDataSets);
973  SafeDelete(fClipboard);
974 
975  fCleanups->Clear();
976  delete fPluginManager; gPluginMgr = fPluginManager = 0;
977  delete gClassTable; gClassTable = 0;
978  delete gEnv; gEnv = 0;
979 
980  if (fTypes) fTypes->Delete();
981  SafeDelete(fTypes);
982  if (fGlobals) fGlobals->Delete();
983  SafeDelete(fGlobals);
984  if (fGlobalFunctions) fGlobalFunctions->Delete();
985  SafeDelete(fGlobalFunctions);
986  fEnums.load()->Delete();
987 
988  // FIXME: Causes segfault in rootcling, debug and uncomment
989  // fClasses->Delete(); SafeDelete(fClasses); // TClass'es must be deleted last
990 #endif
991 
992  // Remove shared libraries produced by the TSystem::CompileMacro() call
993  gSystem->CleanCompiledMacros();
994 
995  // Cleanup system class
996  delete gSystem;
997 
998  // ROOT-6022:
999  // if (gInterpreterLib) dlclose(gInterpreterLib);
1000 #ifdef R__COMPLETE_MEM_TERMINATION
1001  // On some 'newer' platform (Fedora Core 17+, Ubuntu 12), the
1002  // initialization order is (by default?) is 'wrong' and so we can't
1003  // delete the interpreter now .. because any of the static in the
1004  // interpreter's library have already been deleted.
1005  // On the link line, we must list the most dependent .o file
1006  // and end with the least dependent (LLVM libraries), unfortunately,
1007  // Fedora Core 17+ or Ubuntu 12 will also execute the initialization
1008  // in the same order (hence doing libCore's before LLVM's and
1009  // vice et versa for both the destructor. We worked around the
1010  // initialization order by delay the TROOT creation until first use.
1011  // We can not do the same for destruction as we have no way of knowing
1012  // the last access ...
1013  // So for now, let's avoid delete TCling except in the special build
1014  // checking the completeness of the termination deletion.
1015 
1016  // TODO: Should we do more cleanup here than just call delete?
1017  // Segfaults rootcling in some cases, debug and uncomment:
1018  //
1019  // delete fInterpreter;
1020 
1021  // We cannot delete fCleanups because of the logic in atexit which needs it.
1022  SafeDelete(fCleanups);
1023 #endif
1024 
1025 #ifndef _MSC_VER
1026  // deleting the interpreter makes things crash at exit in some cases
1027  delete fInterpreter;
1028 #endif
1029 
1030  // Prints memory stats
1031  TStorage::PrintStatistics();
1032 
1033  gROOTLocal = 0;
1034  fgRootInit = kFALSE;
1035  }
1036 }
1037 
1038 ////////////////////////////////////////////////////////////////////////////////
1039 /// Add a class to the list and map of classes.
1040 /// This routine is deprecated, use TClass::AddClass directly.
1041 
1042 void TROOT::AddClass(TClass *cl)
1043 {
1044  TClass::AddClass(cl);
1045 }
1046 
1047 ////////////////////////////////////////////////////////////////////////////////
1048 /// Add a class generator. This generator will be called by TClass::GetClass
1049 /// in case its does not find a loaded rootcint dictionary to request the
1050 /// creation of a TClass object.
1051 
1052 void TROOT::AddClassGenerator(TClassGenerator *generator)
1053 {
1054  if (!generator) return;
1055  fClassGenerators->Add(generator);
1056 }
1057 
1058 ////////////////////////////////////////////////////////////////////////////////
1059 /// Append object to this directory.
1060 ///
1061 /// If replace is true:
1062 /// remove any existing objects with the same same (if the name is not "")
1063 
1064 void TROOT::Append(TObject *obj, Bool_t replace /* = kFALSE */)
1065 {
1066  R__LOCKGUARD(gROOTMutex);
1067  TDirectory::Append(obj,replace);
1068 }
1069 
1070 ////////////////////////////////////////////////////////////////////////////////
1071 /// Add browsable objects to TBrowser.
1072 
1073 void TROOT::Browse(TBrowser *b)
1074 {
1075  TObject *obj;
1076  TIter next(fBrowsables);
1077 
1078  while ((obj = (TObject *) next())) {
1079  const char *opt = next.GetOption();
1080  if (opt && strlen(opt))
1081  b->Add(obj, opt);
1082  else
1083  b->Add(obj, obj->GetName());
1084  }
1085 }
1086 
1087 ////////////////////////////////////////////////////////////////////////////////
1088 /// return class status bit kClassSaved for class cl
1089 /// This function is called by the SavePrimitive functions writing
1090 /// the C++ code for an object.
1091 
1092 Bool_t TROOT::ClassSaved(TClass *cl)
1093 {
1094  if (cl == 0) return kFALSE;
1095  if (cl->TestBit(TClass::kClassSaved)) return kTRUE;
1096  cl->SetBit(TClass::kClassSaved);
1097  return kFALSE;
1098 }
1099 
1100 namespace {
1101  static void R__ListSlowClose(TList *files)
1102  {
1103  // Routine to close a list of files using the 'slow' techniques
1104  // that also for the deletion ot update the list itself.
1105 
1106  static TObject harmless;
1107  TObjLink *cursor = files->FirstLink();
1108  while (cursor) {
1109  TDirectory *dir = static_cast<TDirectory*>( cursor->GetObject() );
1110  if (dir) {
1111  // In order for the iterator to stay valid, we must
1112  // prevent the removal of the object (dir) from the list
1113  // (which is done in TFile::Close). We can also can not
1114  // just move to the next iterator since the Close might
1115  // also (indirectly) remove that file.
1116  // So we SetObject to a harmless value, so that 'dir'
1117  // is not seen as part of the list.
1118  // We will later, remove all the object (see files->Clear()
1119  cursor->SetObject(&harmless); // this must not be zero otherwise things go wrong.
1120  // See related comment at the files->Clear("nodelete");
1121  dir->Close("nodelete");
1122  // Put it back
1123  cursor->SetObject(dir);
1124  }
1125  cursor = cursor->Next();
1126  };
1127  // Now were done, clear the list but do not delete the objects as
1128  // they have been moved to the list of closed objects and must be
1129  // deleted from there in order to avoid a double delete from a
1130  // use objects (on the interpreter stack).
1131  files->Clear("nodelete");
1132  }
1133 
1134  static void R__ListSlowDeleteContent(TList *files)
1135  {
1136  // Routine to delete the content of list of files using the 'slow' techniques
1137 
1138  static TObject harmless;
1139  TObjLink *cursor = files->FirstLink();
1140  while (cursor) {
1141  TDirectory *dir = dynamic_cast<TDirectory*>( cursor->GetObject() );
1142  if (dir) {
1143  // In order for the iterator to stay valid, we must
1144  // prevent the removal of the object (dir) from the list
1145  // (which is done in TFile::Close). We can also can not
1146  // just move to the next iterator since the Close might
1147  // also (indirectly) remove that file.
1148  // So we SetObject to a harmless value, so that 'dir'
1149  // is not seen as part of the list.
1150  // We will later, remove all the object (see files->Clear()
1151  cursor->SetObject(&harmless); // this must not be zero otherwise things go wrong.
1152  // See related comment at the files->Clear("nodelete");
1153  dir->GetList()->Delete("slow");
1154  // Put it back
1155  cursor->SetObject(dir);
1156  }
1157  cursor = cursor->Next();
1158  };
1159  }
1160 }
1161 
1162 ////////////////////////////////////////////////////////////////////////////////
1163 /// Close any files and sockets that gROOT knows about.
1164 /// This can be used to insures that the files and sockets are closed before any library is unloaded!
1165 
1166 void TROOT::CloseFiles()
1167 {
1168  if (fFiles && fFiles->First()) {
1169  R__ListSlowClose(static_cast<TList*>(fFiles));
1170  }
1171  // and Close TROOT itself.
1172  Close("slow");
1173  // Now sockets.
1174  if (fSockets && fSockets->First()) {
1175  if (0==fCleanups->FindObject(fSockets) ) {
1176  fCleanups->Add(fSockets);
1177  fSockets->SetBit(kMustCleanup);
1178  }
1179  CallFunc_t *socketCloser = gInterpreter->CallFunc_Factory();
1180  Long_t offset = 0;
1181  TClass *socketClass = TClass::GetClass("TSocket");
1182  gInterpreter->CallFunc_SetFuncProto(socketCloser, socketClass->GetClassInfo(), "Close", "", &offset);
1183  if (gInterpreter->CallFunc_IsValid(socketCloser)) {
1184  static TObject harmless;
1185  TObjLink *cursor = static_cast<TList*>(fSockets)->FirstLink();
1186  TList notclosed;
1187  while (cursor) {
1188  TObject *socket = cursor->GetObject();
1189  // In order for the iterator to stay valid, we must
1190  // prevent the removal of the object (dir) from the list
1191  // (which is done in TFile::Close). We can also can not
1192  // just move to the next iterator since the Close might
1193  // also (indirectly) remove that file.
1194  // So we SetObject to a harmless value, so that 'dir'
1195  // is not seen as part of the list.
1196  // We will later, remove all the object (see files->Clear()
1197  cursor->SetObject(&harmless); // this must not be zero otherwise things go wrong.
1198 
1199  if (socket->IsA()->InheritsFrom(socketClass)) {
1200  gInterpreter->CallFunc_Exec(socketCloser, ((char*)socket)+offset);
1201  // Put the object in the closed list for later deletion.
1202  socket->SetBit(kMustCleanup);
1203  fClosedObjects->AddLast(socket);
1204  } else {
1205  // Crap ... this is not a socket, likely Proof or something, let's try to find a Close
1206  Long_t other_offset;
1207  CallFunc_t *otherCloser = gInterpreter->CallFunc_Factory();
1208  gInterpreter->CallFunc_SetFuncProto(otherCloser, socket->IsA()->GetClassInfo(), "Close", "", &other_offset);
1209  if (gInterpreter->CallFunc_IsValid(otherCloser)) {
1210  gInterpreter->CallFunc_Exec(otherCloser, ((char*)socket)+other_offset);
1211  // Put the object in the closed list for later deletion.
1212  socket->SetBit(kMustCleanup);
1213  fClosedObjects->AddLast(socket);
1214  } else {
1215  notclosed.AddLast(socket);
1216  }
1217  gInterpreter->CallFunc_Delete(otherCloser);
1218  // Put it back
1219  cursor->SetObject(socket);
1220  }
1221  cursor = cursor->Next();
1222  }
1223  // Now were done, clear the list
1224  fSockets->Clear();
1225  // Read the one we did not close
1226  cursor = notclosed.FirstLink();
1227  while (cursor) {
1228  static_cast<TList*>(fSockets)->AddLast(cursor->GetObject());
1229  cursor = cursor->Next();
1230  }
1231  }
1232  gInterpreter->CallFunc_Delete(socketCloser);
1233  }
1234  if (fMappedFiles && fMappedFiles->First()) {
1235  R__ListSlowClose(static_cast<TList*>(fMappedFiles));
1236  }
1237 
1238 }
1239 
1240 ////////////////////////////////////////////////////////////////////////////////
1241 /// Execute the cleanups necessary at the end of the process, in particular
1242 /// those that must be executed before the library start being unloaded.
1243 
1244 void TROOT::EndOfProcessCleanups()
1245 {
1246  // This will not delete the objects 'held' by the TFiles so that
1247  // they can still be 'reacheable' when ResetGlobals is run.
1248  CloseFiles();
1249 
1250  if (gInterpreter) {
1251  gInterpreter->ResetGlobals();
1252  }
1253 
1254  // Now delete the objects 'held' by the TFiles so that it
1255  // is done before the tear down of the libraries.
1256  if (fClosedObjects && fClosedObjects->First()) {
1257  R__ListSlowDeleteContent(static_cast<TList*>(fClosedObjects));
1258  }
1259 
1260  // Now a set of simpler things to delete. See the same ordering in
1261  // TROOT::~TROOT
1262  fFunctions->Delete();
1263  fGeometries->Delete();
1264  fBrowsers->Delete();
1265  fCanvases->Delete("slow");
1266  fColors->Delete();
1267  fStyles->Delete();
1268 
1269  TQObject::BlockAllSignals(kTRUE);
1270 
1271  if (gInterpreter) {
1272  gInterpreter->ShutDown();
1273  }
1274 }
1275 
1276 
1277 ////////////////////////////////////////////////////////////////////////////////
1278 /// Find an object in one Root folder
1279 
1280 TObject *TROOT::FindObject(const TObject *) const
1281 {
1282  Error("FindObject","Not yet implemented");
1283  return 0;
1284 }
1285 
1286 ////////////////////////////////////////////////////////////////////////////////
1287 /// Returns address of a ROOT object if it exists
1288 ///
1289 /// If name contains at least one "/" the function calls FindObjectany
1290 /// else
1291 /// This function looks in the following order in the ROOT lists:
1292 /// - List of files
1293 /// - List of memory mapped files
1294 /// - List of functions
1295 /// - List of geometries
1296 /// - List of canvases
1297 /// - List of styles
1298 /// - List of specials
1299 /// - List of materials in current geometry
1300 /// - List of shapes in current geometry
1301 /// - List of matrices in current geometry
1302 /// - List of Nodes in current geometry
1303 /// - Current Directory in memory
1304 /// - Current Directory on file
1305 
1306 TObject *TROOT::FindObject(const char *name) const
1307 {
1308  if (name && strstr(name,"/")) return FindObjectAny(name);
1309 
1310  TObject *temp = 0;
1311 
1312  temp = fFiles->FindObject(name); if (temp) return temp;
1313  temp = fMappedFiles->FindObject(name); if (temp) return temp;
1314  {
1315  R__LOCKGUARD(gROOTMutex);
1316  temp = fFunctions->FindObject(name);if (temp) return temp;
1317  }
1318  temp = fGeometries->FindObject(name); if (temp) return temp;
1319  temp = fCanvases->FindObject(name); if (temp) return temp;
1320  temp = fStyles->FindObject(name); if (temp) return temp;
1321  {
1322  R__LOCKGUARD(gROOTMutex);
1323  temp = fSpecials->FindObject(name); if (temp) return temp;
1324  }
1325  TIter next(fGeometries);
1326  TObject *obj;
1327  while ((obj=next())) {
1328  temp = obj->FindObject(name); if (temp) return temp;
1329  }
1330  if (gDirectory) temp = gDirectory->Get(name);
1331  if (temp) return temp;
1332  if (gPad) {
1333  TVirtualPad *canvas = gPad->GetVirtCanvas();
1334  if (fCanvases->FindObject(canvas)) { //this check in case call from TCanvas ctor
1335  temp = canvas->FindObject(name);
1336  if (!temp && canvas != gPad) temp = gPad->FindObject(name);
1337  }
1338  }
1339  return temp;
1340 }
1341 
1342 ////////////////////////////////////////////////////////////////////////////////
1343 /// Returns address and folder of a ROOT object if it exists
1344 ///
1345 /// This function looks in the following order in the ROOT lists:
1346 /// - List of files
1347 /// - List of memory mapped files
1348 /// - List of functions
1349 /// - List of geometries
1350 /// - List of canvases
1351 /// - List of styles
1352 /// - List of specials
1353 /// - List of materials in current geometry
1354 /// - List of shapes in current geometry
1355 /// - List of matrices in current geometry
1356 /// - List of Nodes in current geometry
1357 /// - Current Directory in memory
1358 /// - Current Directory on file
1359 
1360 TObject *TROOT::FindSpecialObject(const char *name, void *&where)
1361 {
1362  TObject *temp = 0;
1363  where = 0;
1364 
1365  if (!temp) {
1366  temp = fFiles->FindObject(name);
1367  where = fFiles;
1368  }
1369  if (!temp) {
1370  temp = fMappedFiles->FindObject(name);
1371  where = fMappedFiles;
1372  }
1373  if (!temp) {
1374  R__LOCKGUARD(gROOTMutex);
1375  temp = fFunctions->FindObject(name);
1376  where = fFunctions;
1377  }
1378  if (!temp) {
1379  temp = fCanvases->FindObject(name);
1380  where = fCanvases;
1381  }
1382  if (!temp) {
1383  temp = fStyles->FindObject(name);
1384  where = fStyles;
1385  }
1386  if (!temp) {
1387  temp = fSpecials->FindObject(name);
1388  where = fSpecials;
1389  }
1390  if (!temp) {
1391  TObject *glast = fGeometries->Last();
1392  if (glast) {where = glast; temp = glast->FindObject(name);}
1393  }
1394  if (!temp && gDirectory) {
1395  temp = gDirectory->Get(name);
1396  where = gDirectory;
1397  }
1398  if (!temp && gPad) {
1399  TVirtualPad *canvas = gPad->GetVirtCanvas();
1400  if (fCanvases->FindObject(canvas)) { //this check in case call from TCanvas ctor
1401  temp = canvas->FindObject(name);
1402  where = canvas;
1403  if (!temp && canvas != gPad) {
1404  temp = gPad->FindObject(name);
1405  where = gPad;
1406  }
1407  }
1408  }
1409  if (!temp) return 0;
1410  if (temp->TestBit(kNotDeleted)) return temp;
1411  return 0;
1412 }
1413 
1414 ////////////////////////////////////////////////////////////////////////////////
1415 /// Return a pointer to the first object with name starting at //root.
1416 /// This function scans the list of all folders.
1417 /// if no object found in folders, it scans the memory list of all files.
1418 
1419 TObject *TROOT::FindObjectAny(const char *name) const
1420 {
1421  TObject *obj = fRootFolder->FindObjectAny(name);
1422  if (obj) return obj;
1423  return gDirectory->FindObjectAnyFile(name);
1424 }
1425 
1426 ////////////////////////////////////////////////////////////////////////////////
1427 /// Scan the memory lists of all files for an object with name
1428 
1429 TObject *TROOT::FindObjectAnyFile(const char *name) const
1430 {
1431  R__LOCKGUARD(gROOTMutex);
1432  TDirectory *d;
1433  TIter next(GetListOfFiles());
1434  while ((d = (TDirectory*)next())) {
1435  // Call explicitly TDirectory::FindObject to restrict the search to the
1436  // already in memory object.
1437  TObject *obj = d->TDirectory::FindObject(name);
1438  if (obj) return obj;
1439  }
1440  return 0;
1441 }
1442 
1443 ////////////////////////////////////////////////////////////////////////////////
1444 /// Returns class name of a ROOT object including CINT globals.
1445 
1446 const char *TROOT::FindObjectClassName(const char *name) const
1447 {
1448  // Search first in the list of "standard" objects
1449  TObject *obj = FindObject(name);
1450  if (obj) return obj->ClassName();
1451 
1452  // Is it a global variable?
1453  TGlobal *g = GetGlobal(name);
1454  if (g) return g->GetTypeName();
1455 
1456  return 0;
1457 }
1458 
1459 ////////////////////////////////////////////////////////////////////////////////
1460 /// Return path name of obj somewhere in the //root/... path.
1461 /// The function returns the first occurence of the object in the list
1462 /// of folders. The returned string points to a static char array in TROOT.
1463 /// If this function is called in a loop or recursively, it is the
1464 /// user's responsibility to copy this string in their area.
1465 
1466 const char *TROOT::FindObjectPathName(const TObject *) const
1467 {
1468  Error("FindObjectPathName","Not yet implemented");
1469  return "??";
1470 }
1471 
1472 ////////////////////////////////////////////////////////////////////////////////
1473 /// return a TClass object corresponding to 'name' assuming it is an STL container.
1474 /// In particular we looking for possible alternative name (default template
1475 /// parameter, typedefs template arguments, typedefed name).
1476 
1477 TClass *TROOT::FindSTLClass(const char *name, Bool_t load, Bool_t silent) const
1478 {
1479  // Example of inputs are
1480  // vector<int> (*)
1481  // vector<Int_t>
1482  // vector<long long>
1483  // vector<Long_64_t> (*)
1484  // vector<int, allocator<int> >
1485  // vector<Int_t, allocator<int> >
1486  //
1487  // One of the possibly expensive operation is the resolving of the typedef
1488  // which can provoke the parsing of the header files (and/or the loading
1489  // of clang pcms information).
1490 
1491  R__LOCKGUARD(gInterpreterMutex);
1492 
1493  // Remove std::, allocator, typedef, add Long64_t, etc. in just one call.
1494  std::string normalized;
1495  TClassEdit::GetNormalizedName(normalized, name);
1496 
1497  TClass *cl = 0;
1498  if (normalized != name) cl = TClass::GetClass(normalized.c_str(),load,silent);
1499 
1500  if (load && cl==0) {
1501  // Create an Emulated class for this container.
1502  cl = gInterpreter->GenerateTClass(normalized.c_str(), kTRUE, silent);
1503  }
1504 
1505  return cl;
1506 }
1507 
1508 ////////////////////////////////////////////////////////////////////////////////
1509 /// Return pointer to class with name. Obsolete, use TClass::GetClass directly
1510 
1511 TClass *TROOT::GetClass(const char *name, Bool_t load, Bool_t silent) const
1512 {
1513  return TClass::GetClass(name,load,silent);
1514 }
1515 
1516 
1517 ////////////////////////////////////////////////////////////////////////////////
1518 /// Return pointer to class from its name. Obsolete, use TClass::GetClass directly
1519 /// See TClass::GetClass
1520 
1521 TClass *TROOT::GetClass(const std::type_info& typeinfo, Bool_t load, Bool_t silent) const
1522 {
1523  return TClass::GetClass(typeinfo,load,silent);
1524 }
1525 
1526 ////////////////////////////////////////////////////////////////////////////////
1527 /// Return address of color with index color.
1528 
1529 TColor *TROOT::GetColor(Int_t color) const
1530 {
1531  TColor::InitializeColors();
1532  TObjArray *lcolors = (TObjArray*) GetListOfColors();
1533  if (!lcolors) return 0;
1534  if (color < 0 || color >= lcolors->GetSize()) return 0;
1535  TColor *col = (TColor*)lcolors->At(color);
1536  if (col && col->GetNumber() == color) return col;
1537  TIter next(lcolors);
1538  while ((col = (TColor *) next()))
1539  if (col->GetNumber() == color) return col;
1540 
1541  return 0;
1542 }
1543 
1544 ////////////////////////////////////////////////////////////////////////////////
1545 /// Return a default canvas.
1546 
1547 TCanvas *TROOT::MakeDefCanvas() const
1548 {
1549  return (TCanvas*)gROOT->ProcessLine("TCanvas::MakeDefCanvas();");
1550 }
1551 
1552 ////////////////////////////////////////////////////////////////////////////////
1553 /// Return pointer to type with name.
1554 
1555 TDataType *TROOT::GetType(const char *name, Bool_t /* load */) const
1556 {
1557  return (TDataType*)gROOT->GetListOfTypes()->FindObject(name);
1558 }
1559 
1560 ////////////////////////////////////////////////////////////////////////////////
1561 /// Return pointer to file with name.
1562 
1563 TFile *TROOT::GetFile(const char *name) const
1564 {
1565  R__LOCKGUARD(gROOTMutex);
1566  return (TFile*)GetListOfFiles()->FindObject(name);
1567 }
1568 
1569 ////////////////////////////////////////////////////////////////////////////////
1570 /// Return pointer to style with name
1571 
1572 TStyle *TROOT::GetStyle(const char *name) const
1573 {
1574  return (TStyle*)GetListOfStyles()->FindObject(name);
1575 }
1576 
1577 ////////////////////////////////////////////////////////////////////////////////
1578 /// Return pointer to function with name.
1579 
1580 TObject *TROOT::GetFunction(const char *name) const
1581 {
1582  if (name == 0 || name[0] == 0) {
1583  return 0;
1584  }
1585 
1586  {
1587  R__LOCKGUARD(gROOTMutex);
1588  TObject *f1 = fFunctions->FindObject(name);
1589  if (f1) return f1;
1590  }
1591 
1592  gROOT->ProcessLine("TF1::InitStandardFunctions();");
1593 
1594  R__LOCKGUARD(gROOTMutex);
1595  return fFunctions->FindObject(name);
1596 }
1597 
1598 ////////////////////////////////////////////////////////////////////////////////
1599 
1600 TFunctionTemplate *TROOT::GetFunctionTemplate(const char *name)
1601 {
1602  if (!gInterpreter) return 0;
1603 
1604  if (!fFuncTemplate) fFuncTemplate = new TListOfFunctionTemplates(0);
1605 
1606  return (TFunctionTemplate*)fFuncTemplate->FindObject(name);
1607 }
1608 
1609 ////////////////////////////////////////////////////////////////////////////////
1610 /// Return pointer to global variable by name. If load is true force
1611 /// reading of all currently defined globals from CINT (more expensive).
1612 
1613 TGlobal *TROOT::GetGlobal(const char *name, Bool_t load) const
1614 {
1615  return (TGlobal *)gROOT->GetListOfGlobals(load)->FindObject(name);
1616 }
1617 
1618 ////////////////////////////////////////////////////////////////////////////////
1619 /// Return pointer to global variable with address addr.
1620 
1621 TGlobal *TROOT::GetGlobal(const TObject *addr, Bool_t /* load */) const
1622 {
1623  if (addr == 0 || ((Long_t)addr) == -1) return 0;
1624 
1625  TInterpreter::DeclId_t decl = gInterpreter->GetDataMemberAtAddr(addr);
1626  if (decl) {
1627  TListOfDataMembers *globals = ((TListOfDataMembers*)(gROOT->GetListOfGlobals(kFALSE)));
1628  return (TGlobal*)globals->Get(decl);
1629  }
1630  // If we are actually looking for a global that is held by a global
1631  // pointer (for example gRandom), we need to find a pointer with the
1632  // correct value.
1633  decl = gInterpreter->GetDataMemberWithValue(addr);
1634  if (decl) {
1635  TListOfDataMembers *globals = ((TListOfDataMembers*)(gROOT->GetListOfGlobals(kFALSE)));
1636  return (TGlobal*)globals->Get(decl);
1637  }
1638  return 0;
1639 }
1640 
1641 ////////////////////////////////////////////////////////////////////////////////
1642 /// Internal routine returning, and creating if necessary, the list
1643 /// of global function.
1644 
1645 TListOfFunctions *TROOT::GetGlobalFunctions()
1646 {
1647  if (!fGlobalFunctions) fGlobalFunctions = new TListOfFunctions(0);
1648  return fGlobalFunctions;
1649 }
1650 
1651 ////////////////////////////////////////////////////////////////////////////////
1652 /// Return the collection of functions named "name".
1653 
1654 TCollection *TROOT::GetListOfFunctionOverloads(const char* name) const
1655 {
1656  return ((TListOfFunctions*)fGlobalFunctions)->GetListForObject(name);
1657 }
1658 
1659 ////////////////////////////////////////////////////////////////////////////////
1660 /// Return pointer to global function by name.
1661 /// If params != 0 it will also resolve overloading other it returns the first
1662 /// name match.
1663 /// If params == 0 and load is true force reading of all currently defined
1664 /// global functions from Cling.
1665 /// The param string must be of the form: "3189,\"aap\",1.3".
1666 
1667 TFunction *TROOT::GetGlobalFunction(const char *function, const char *params,
1668  Bool_t load)
1669 {
1670  if (!params) {
1671  R__LOCKGUARD(gROOTMutex);
1672  return (TFunction *)GetListOfGlobalFunctions(load)->FindObject(function);
1673  } else {
1674  if (!fInterpreter)
1675  Fatal("GetGlobalFunction", "fInterpreter not initialized");
1676 
1677  R__LOCKGUARD(gROOTMutex);
1678  TInterpreter::DeclId_t decl = gInterpreter->GetFunctionWithValues(0,
1679  function, params,
1680  false);
1681 
1682  if (!decl) return 0;
1683 
1684  TFunction *f = GetGlobalFunctions()->Get(decl);
1685  if (f) return f;
1686 
1687  Error("GetGlobalFunction",
1688  "\nDid not find matching TFunction <%s> with \"%s\".",
1689  function,params);
1690  return 0;
1691  }
1692 }
1693 
1694 ////////////////////////////////////////////////////////////////////////////////
1695 /// Return pointer to global function by name. If proto != 0
1696 /// it will also resolve overloading. If load is true force reading
1697 /// of all currently defined global functions from CINT (more expensive).
1698 /// The proto string must be of the form: "int, char*, float".
1699 
1700 TFunction *TROOT::GetGlobalFunctionWithPrototype(const char *function,
1701  const char *proto, Bool_t load)
1702 {
1703  if (!proto) {
1704  R__LOCKGUARD(gROOTMutex);
1705  return (TFunction *)GetListOfGlobalFunctions(load)->FindObject(function);
1706  } else {
1707  if (!fInterpreter)
1708  Fatal("GetGlobalFunctionWithPrototype", "fInterpreter not initialized");
1709 
1710  R__LOCKGUARD(gROOTMutex);
1711  TInterpreter::DeclId_t decl = gInterpreter->GetFunctionWithPrototype(0,
1712  function, proto);
1713 
1714  if (!decl) return 0;
1715 
1716  TFunction *f = GetGlobalFunctions()->Get(decl);
1717  if (f) return f;
1718 
1719  Error("GetGlobalFunctionWithPrototype",
1720  "\nDid not find matching TFunction <%s> with \"%s\".",
1721  function,proto);
1722  return 0;
1723  }
1724 }
1725 
1726 ////////////////////////////////////////////////////////////////////////////////
1727 /// Return pointer to Geometry with name
1728 
1729 TObject *TROOT::GetGeometry(const char *name) const
1730 {
1731  return GetListOfGeometries()->FindObject(name);
1732 }
1733 
1734 ////////////////////////////////////////////////////////////////////////////////
1735 
1736 TCollection *TROOT::GetListOfEnums(Bool_t load /* = kTRUE */)
1737 {
1738  if(!fEnums.load()) {
1739  R__LOCKGUARD(gROOTMutex);
1740  // Test again just in case, another thread did the work while we were
1741  // waiting.
1742  if (!fEnums.load()) fEnums = new TListOfEnumsWithLock(0);
1743  }
1744  if (load) {
1745  R__LOCKGUARD(gROOTMutex);
1746  (*fEnums).Load(); // Refresh the list of enums.
1747  }
1748  return fEnums.load();
1749 }
1750 
1751 ////////////////////////////////////////////////////////////////////////////////
1752 
1753 TCollection *TROOT::GetListOfFunctionTemplates()
1754 {
1755  R__LOCKGUARD(gROOTMutex);
1756  if(!fFuncTemplate) {
1757  fFuncTemplate = new TListOfFunctionTemplates(0);
1758  }
1759  return fFuncTemplate;
1760 }
1761 
1762 ////////////////////////////////////////////////////////////////////////////////
1763 /// Return list containing the TGlobals currently defined.
1764 /// Since globals are created and deleted during execution of the
1765 /// program, we need to update the list of globals every time we
1766 /// execute this method. However, when calling this function in
1767 /// a (tight) loop where no interpreter symbols will be created
1768 /// you can set load=kFALSE (default).
1769 
1770 TCollection *TROOT::GetListOfGlobals(Bool_t load)
1771 {
1772  if (!fGlobals) {
1773  fGlobals = new TListOfDataMembers(0);
1774  // We add to the list the "funcky-fake" globals.
1775 
1776  // provide special functor for gROOT, while ROOT::GetROOT() does not return reference
1777  TGlobalMappedFunction::MakeFunctor("gROOT", "TROOT*", ROOT::GetROOT, [] {
1778  ROOT::GetROOT();
1779  return (void *)&ROOT::Internal::gROOTLocal;
1780  });
1781 
1782  TGlobalMappedFunction::MakeFunctor("gPad", "TVirtualPad*", TVirtualPad::Pad);
1783  TGlobalMappedFunction::MakeFunctor("gVirtualX", "TVirtualX*", TVirtualX::Instance);
1784  TGlobalMappedFunction::MakeFunctor("gDirectory", "TDirectory*", TDirectory::CurrentDirectory);
1785 
1786  // Don't let TGlobalMappedFunction delete our globals, now that we take them.
1787  fGlobals->AddAll(&TGlobalMappedFunction::GetEarlyRegisteredGlobals());
1788  TGlobalMappedFunction::GetEarlyRegisteredGlobals().SetOwner(kFALSE);
1789  TGlobalMappedFunction::GetEarlyRegisteredGlobals().Clear();
1790  }
1791 
1792  if (!fInterpreter)
1793  Fatal("GetListOfGlobals", "fInterpreter not initialized");
1794 
1795  if (load) fGlobals->Load();
1796 
1797  return fGlobals;
1798 }
1799 
1800 ////////////////////////////////////////////////////////////////////////////////
1801 /// Return list containing the TFunctions currently defined.
1802 /// Since functions are created and deleted during execution of the
1803 /// program, we need to update the list of functions every time we
1804 /// execute this method. However, when calling this function in
1805 /// a (tight) loop where no interpreter symbols will be created
1806 /// you can set load=kFALSE (default).
1807 
1808 TCollection *TROOT::GetListOfGlobalFunctions(Bool_t load)
1809 {
1810  R__LOCKGUARD(gROOTMutex);
1811 
1812  if (!fGlobalFunctions) {
1813  fGlobalFunctions = new TListOfFunctions(0);
1814  }
1815 
1816  if (!fInterpreter)
1817  Fatal("GetListOfGlobalFunctions", "fInterpreter not initialized");
1818 
1819  // A thread that calls with load==true and a thread that calls with load==false
1820  // will conflict here (the load==true will be updating the list while the
1821  // other is reading it). To solve the problem, we could use a read-write lock
1822  // inside the list itself.
1823  if (load) fGlobalFunctions->Load();
1824 
1825  return fGlobalFunctions;
1826 }
1827 
1828 ////////////////////////////////////////////////////////////////////////////////
1829 /// Return a dynamic list giving access to all TDataTypes (typedefs)
1830 /// currently defined.
1831 ///
1832 /// The list is populated on demand. Calling
1833 /// ~~~ {.cpp}
1834 /// gROOT->GetListOfTypes()->FindObject(nameoftype);
1835 /// ~~~
1836 /// will return the TDataType corresponding to 'nameoftype'. If the
1837 /// TDataType is not already in the list itself and the type does exist,
1838 /// a new TDataType will be created and added to the list.
1839 ///
1840 /// Calling
1841 /// ~~~ {.cpp}
1842 /// gROOT->GetListOfTypes()->ls(); // or Print()
1843 /// ~~~
1844 /// list only the typedefs that have been previously accessed through the
1845 /// list (plus the builtins types).
1846 
1847 TCollection *TROOT::GetListOfTypes(Bool_t /* load */)
1848 {
1849  if (!fInterpreter)
1850  Fatal("GetListOfTypes", "fInterpreter not initialized");
1851 
1852  return fTypes;
1853 }
1854 
1855 
1856 ////////////////////////////////////////////////////////////////////////////////
1857 /// Execute command when system has been idle for idleTimeInSec seconds.
1858 
1859 void TROOT::Idle(UInt_t idleTimeInSec, const char *command)
1860 {
1861  if (!fApplication.load())
1862  TApplication::CreateApplication();
1863 
1864  if (idleTimeInSec <= 0)
1865  (*fApplication).RemoveIdleTimer();
1866  else
1867  (*fApplication).SetIdleTimer(idleTimeInSec, command);
1868 }
1869 
1870 ////////////////////////////////////////////////////////////////////////////////
1871 /// Check whether className is a known class, and only autoload
1872 /// if we can. Helper function for TROOT::IgnoreInclude().
1873 
1874 static TClass* R__GetClassIfKnown(const char* className)
1875 {
1876  // Check whether the class is available for auto-loading first:
1877  const char* libsToLoad = gInterpreter->GetClassSharedLibs(className);
1878  TClass* cla = 0;
1879  if (libsToLoad) {
1880  // trigger autoload, and only create TClass in this case.
1881  return TClass::GetClass(className);
1882  } else if (gROOT->GetListOfClasses()
1883  && (cla = (TClass*)gROOT->GetListOfClasses()->FindObject(className))) {
1884  // cla assigned in if statement
1885  } else if (gClassTable->FindObject(className)) {
1886  return TClass::GetClass(className);
1887  }
1888  return cla;
1889 }
1890 
1891 ////////////////////////////////////////////////////////////////////////////////
1892 /// Return 1 if the name of the given include file corresponds to a class that
1893 /// is known to ROOT, e.g. "TLorentzVector.h" versus TLorentzVector.
1894 
1895 Int_t TROOT::IgnoreInclude(const char *fname, const char * /*expandedfname*/)
1896 {
1897  if (fname == 0) return 0;
1898 
1899  TString stem(fname);
1900  // Remove extension if any, ignore files with extension not being .h*
1901  Int_t where = stem.Last('.');
1902  if (where != kNPOS) {
1903  if (stem.EndsWith(".so") || stem.EndsWith(".sl") ||
1904  stem.EndsWith(".dl") || stem.EndsWith(".a") ||
1905  stem.EndsWith(".dll", TString::kIgnoreCase))
1906  return 0;
1907  stem.Remove(where);
1908  }
1909 
1910  TString className = gSystem->BaseName(stem);
1911  TClass* cla = R__GetClassIfKnown(className);
1912  if (!cla) {
1913  // Try again with modifications to the file name:
1914  className = stem;
1915  className.ReplaceAll("/", "::");
1916  className.ReplaceAll("\\", "::");
1917  if (className.Contains(":::")) {
1918  // "C:\dir" becomes "C:::dir".
1919  // fname corresponds to whatever is stated after #include and
1920  // a full path name usually means that it's not a regular #include
1921  // but e.g. a ".L", so we can assume that this is not a header of
1922  // a class in a namespace (a global-namespace class would have been
1923  // detected already before).
1924  return 0;
1925  }
1926  cla = R__GetClassIfKnown(className);
1927  }
1928 
1929  if (!cla) {
1930  return 0;
1931  }
1932 
1933  // cla is valid, check wether it's actually in the header of the same name:
1934  if (cla->GetDeclFileLine() <= 0) return 0; // to a void an error with VisualC++
1935  TString decfile = gSystem->BaseName(cla->GetDeclFileName());
1936  if (decfile != gSystem->BaseName(fname)) {
1937  return 0;
1938  }
1939  return 1;
1940 }
1941 
1942 ////////////////////////////////////////////////////////////////////////////////
1943 /// Initialize operating system interface.
1944 
1945 void TROOT::InitSystem()
1946 {
1947  if (gSystem == 0) {
1948 #if defined(R__UNIX)
1949 #if defined(R__HAS_COCOA)
1950  gSystem = new TMacOSXSystem;
1951 #else
1952  gSystem = new TUnixSystem;
1953 #endif
1954 #elif defined(R__WIN32)
1955  gSystem = new TWinNTSystem;
1956 #else
1957  gSystem = new TSystem;
1958 #endif
1959 
1960  if (gSystem->Init())
1961  fprintf(stderr, "Fatal in <TROOT::InitSystem>: can't init operating system layer\n");
1962 
1963  if (!gSystem->HomeDirectory()) {
1964  fprintf(stderr, "Fatal in <TROOT::InitSystem>: HOME directory not set\n");
1965  fprintf(stderr, "Fix this by defining the HOME shell variable\n");
1966  }
1967 
1968  // read default files
1969  gEnv = new TEnv(".rootrc");
1970 
1971  gDebug = gEnv->GetValue("Root.Debug", 0);
1972 
1973  if (!gEnv->GetValue("Root.ErrorHandlers", 1))
1974  gSystem->ResetSignals();
1975 
1976  // The old "Root.ZipMode" had a discrepancy between documentation vs actual meaning.
1977  // Also, a value with the meaning "default" wasn't available. To solved this,
1978  // "Root.ZipMode" was replaced by "Root.CompressionAlgorithm". Warn about usage of
1979  // the old value, if it's set to 0, but silently translate the setting to
1980  // "Root.CompressionAlgorithm" for values > 1.
1981  Int_t oldzipmode = gEnv->GetValue("Root.ZipMode", -1);
1982  if (oldzipmode == 0) {
1983  fprintf(stderr, "Warning in <TROOT::InitSystem>: ignoring old rootrc entry \"Root.ZipMode = 0\"!\n");
1984  } else {
1985  if (oldzipmode == -1 || oldzipmode == 1) {
1986  // Not set or default value, use "default" for "Root.CompressionAlgorithm":
1987  oldzipmode = 0;
1988  }
1989  // else keep the old zipmode (e.g. "3") as "Root.CompressionAlgorithm"
1990  // if "Root.CompressionAlgorithm" isn't set; see below.
1991  }
1992 
1993  Int_t zipmode = gEnv->GetValue("Root.CompressionAlgorithm", oldzipmode);
1994  if (zipmode != 0) R__SetZipMode(zipmode);
1995 
1996  const char *sdeb;
1997  if ((sdeb = gSystem->Getenv("ROOTDEBUG")))
1998  gDebug = atoi(sdeb);
1999 
2000  if (gDebug > 0 && isatty(2))
2001  fprintf(stderr, "Info in <TROOT::InitSystem>: running with gDebug = %d\n", gDebug);
2002 
2003  if (gEnv->GetValue("Root.MemStat", 0))
2004  TStorage::EnableStatistics();
2005  int msize = gEnv->GetValue("Root.MemStat.size", -1);
2006  int mcnt = gEnv->GetValue("Root.MemStat.cnt", -1);
2007  if (msize != -1 || mcnt != -1)
2008  TStorage::EnableStatistics(msize, mcnt);
2009 
2010  fgMemCheck = gEnv->GetValue("Root.MemCheck", 0);
2011 
2012 #if defined(R__HAS_COCOA)
2013  // create and delete a dummy TUrl so that TObjectStat table does not contain
2014  // objects that are deleted after recording is turned-off (in next line),
2015  // like the TUrl::fgSpecialProtocols list entries which are created in the
2016  // TMacOSXSystem ctor.
2017  { TUrl dummy("/dummy"); }
2018 #endif
2019  TObject::SetObjectStat(gEnv->GetValue("Root.ObjectStat", 0));
2020  }
2021 }
2022 
2023 ////////////////////////////////////////////////////////////////////////////////
2024 /// Load and initialize thread library.
2025 
2026 void TROOT::InitThreads()
2027 {
2028  if (gEnv->GetValue("Root.UseThreads", 0) || gEnv->GetValue("Root.EnableThreadSafety", 0)) {
2029  ROOT::EnableThreadSafety();
2030  }
2031 }
2032 
2033 ////////////////////////////////////////////////////////////////////////////////
2034 /// Initialize the interpreter. Should be called only after main(),
2035 /// to make sure LLVM/Clang is fully initialized.
2036 
2037 void TROOT::InitInterpreter()
2038 {
2039  // usedToIdentifyRootClingByDlSym is available when TROOT is part of
2040  // rootcling.
2041  if (!dlsym(RTLD_DEFAULT, "usedToIdentifyRootClingByDlSym")
2042  && !dlsym(RTLD_DEFAULT, "usedToIdentifyStaticRoot")) {
2043  // Make sure no llvm symbols are visible before loading libCling. If they
2044  // exist libCling will use those and not ours, causing havoc in the
2045  // interpreter. Look for an extern "C" symbol to avoid mangling; look for a
2046  // symbol from llvm because clang builds on top, so users would likely
2047  // have also their own llvm symbols when providing their own clang.
2048  void *LLVMEnablePrettyStackTraceAddr = 0;
2049  // Can't use gSystem->DynFindSymbol() because that iterates over all *known*
2050  // libraries which is not the same!
2051  LLVMEnablePrettyStackTraceAddr = dlsym(RTLD_DEFAULT, "LLVMEnablePrettyStackTrace");
2052  // FIXME: When we configure with -Dclingtest=On we intentionally export the symbols. Silence this error.
2053  if (LLVMEnablePrettyStackTraceAddr) {
2054  Error("InitInterpreter()", "LLVM SYMBOLS ARE EXPOSED TO CLING! "
2055  "This will cause problems; please hide them or dlopen() them "
2056  "after the call to TROOT::InitInterpreter()!");
2057  }
2058 
2059  char *libRIO = gSystem->DynamicPathName("libRIO");
2060  void *libRIOHandle = dlopen(libRIO, RTLD_NOW|RTLD_GLOBAL);
2061  delete [] libRIO;
2062  if (!libRIOHandle) {
2063  TString err = dlerror();
2064  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load library %s\n", err.Data());
2065  exit(1);
2066  }
2067 
2068  char *libcling = gSystem->DynamicPathName("libCling");
2069  gInterpreterLib = dlopen(libcling, RTLD_LAZY|RTLD_LOCAL);
2070  delete [] libcling;
2071 
2072  if (!gInterpreterLib) {
2073  TString err = dlerror();
2074  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load library %s\n", err.Data());
2075  exit(1);
2076  }
2077  dlerror(); // reset error message
2078  } else {
2079  gInterpreterLib = RTLD_DEFAULT;
2080  }
2081  CreateInterpreter_t *CreateInterpreter = (CreateInterpreter_t*) dlsym(gInterpreterLib, "CreateInterpreter");
2082  if (!CreateInterpreter) {
2083  TString err = dlerror();
2084  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load symbol %s\n", err.Data());
2085  exit(1);
2086  }
2087  // Schedule the destruction of TROOT.
2088  atexit(at_exit_of_TROOT);
2089 
2090  gDestroyInterpreter = (DestroyInterpreter_t*) dlsym(gInterpreterLib, "DestroyInterpreter");
2091  if (!gDestroyInterpreter) {
2092  TString err = dlerror();
2093  fprintf(stderr, "Fatal in <TROOT::InitInterpreter>: cannot load symbol %s\n", err.Data());
2094  exit(1);
2095  }
2096 
2097  const char *interpArgs[] = {
2098 #ifdef NDEBUG
2099  "-DNDEBUG",
2100 #else
2101  "-UNDEBUG",
2102 #endif
2103 #ifdef DEBUG
2104  "-DDEBUG",
2105 #else
2106  "-UDEBUG",
2107 #endif
2108 #ifdef _DEBUG
2109  "-D_DEBUG",
2110 #else
2111  "-U_DEBUG",
2112 #endif
2113  nullptr};
2114 
2115  fInterpreter = CreateInterpreter(gInterpreterLib, interpArgs);
2116 
2117  fCleanups->Add(fInterpreter);
2118  fInterpreter->SetBit(kMustCleanup);
2119 
2120  fgRootInit = kTRUE;
2121 
2122  // initialize gClassTable is not already done
2123  if (!gClassTable)
2124  new TClassTable;
2125 
2126  // Initialize all registered dictionaries.
2127  for (std::vector<ModuleHeaderInfo_t>::const_iterator
2128  li = GetModuleHeaderInfoBuffer().begin(),
2129  le = GetModuleHeaderInfoBuffer().end(); li != le; ++li) {
2130  // process buffered module registrations
2131  fInterpreter->RegisterModule(li->fModuleName,
2132  li->fHeaders,
2133  li->fIncludePaths,
2134  li->fPayloadCode,
2135  li->fFwdDeclCode,
2136  li->fTriggerFunc,
2137  li->fFwdNargsToKeepColl,
2138  li->fClassesHeaders,
2139  kTRUE /*lateRegistration*/,
2140  li->fHasCxxModule);
2141  }
2142  GetModuleHeaderInfoBuffer().clear();
2143 
2144  fInterpreter->Initialize();
2145 }
2146 
2147 ////////////////////////////////////////////////////////////////////////////////
2148 /// Helper function used by TClass::GetClass().
2149 /// This function attempts to load the dictionary for 'classname'
2150 /// either from the TClassTable or from the list of generator.
2151 /// If silent is 'true', do not warn about missing dictionary for the class.
2152 /// (typically used for class that are used only for transient members)
2153 ///
2154 /// The 'requestedname' is expected to be already normalized.
2155 
2156 TClass *TROOT::LoadClass(const char *requestedname, Bool_t silent) const
2157 {
2158  return TClass::LoadClass(requestedname, silent);
2159 }
2160 
2161 ////////////////////////////////////////////////////////////////////////////////
2162 /// Check if class "classname" is known to the interpreter (in fact,
2163 /// this check is not needed anymore, so classname is ignored). If
2164 /// not it will load library "libname". If the library couldn't be found with original
2165 /// libname and if the name was not prefixed with lib, try to prefix with "lib" and search again.
2166 /// If DynamicPathName still couldn't find the library, return -1.
2167 /// If check is true it will only check if libname exists and is
2168 /// readable.
2169 /// Returns 0 on successful loading, -1 in case libname does not
2170 /// exist or in case of error and -2 in case of version mismatch.
2171 
2172 Int_t TROOT::LoadClass(const char * /*classname*/, const char *libname,
2173  Bool_t check)
2174 {
2175  TString lib(libname);
2176 
2177  // Check if libname exists in path or not
2178  if (char *path = gSystem->DynamicPathName(lib, kTRUE)) {
2179  // If check == true, only check if it exists and if it's readable
2180  if (check) {
2181  delete [] path;
2182  return 0;
2183  }
2184 
2185  // If check == false, try to load the library
2186  else {
2187  int err = gSystem->Load(path, 0, kTRUE);
2188  delete [] path;
2189 
2190  // TSystem::Load returns 1 when the library was already loaded, return success in this case.
2191  if (err == 1)
2192  err = 0;
2193  return err;
2194  }
2195  } else {
2196  // This is the branch where libname didn't exist
2197  if (check) {
2198  FileStat_t stat;
2199  if (!gSystem->GetPathInfo(libname, stat) && (R_ISREG(stat.fMode) &&
2200  !gSystem->AccessPathName(libname, kReadPermission)))
2201  return 0;
2202  }
2203 
2204  // Take care of user who didn't write the whole name
2205  if (!lib.BeginsWith("lib")) {
2206  lib = "lib" + lib;
2207  return LoadClass("", lib.Data(), check);
2208  }
2209  }
2210 
2211  // Execution reaches here when library was prefixed with lib, check is false and couldn't find
2212  // the library name.
2213  return -1;
2214 }
2215 
2216 ////////////////////////////////////////////////////////////////////////////////
2217 /// Return true if the file is local and is (likely) to be a ROOT file
2218 
2219 Bool_t TROOT::IsRootFile(const char *filename) const
2220 {
2221  Bool_t result = kFALSE;
2222  FILE *mayberootfile = fopen(filename,"rb");
2223  if (mayberootfile) {
2224  char header[5];
2225  if (fgets(header,5,mayberootfile)) {
2226  result = strncmp(header,"root",4)==0;
2227  }
2228  fclose(mayberootfile);
2229  }
2230  return result;
2231 }
2232 
2233 ////////////////////////////////////////////////////////////////////////////////
2234 /// To list all objects of the application.
2235 /// Loop on all objects created in the ROOT linked lists.
2236 /// Objects may be files and windows or any other object directly
2237 /// attached to the ROOT linked list.
2238 
2239 void TROOT::ls(Option_t *option) const
2240 {
2241 // TObject::SetDirLevel();
2242 // GetList()->R__FOR_EACH(TObject,ls)(option);
2243  TDirectory::ls(option);
2244 }
2245 
2246 ////////////////////////////////////////////////////////////////////////////////
2247 /// Load a macro in the interpreter's memory. Equivalent to the command line
2248 /// command ".L filename". If the filename has "+" or "++" appended
2249 /// the macro will be compiled by ACLiC. The filename must have the format:
2250 /// [path/]macro.C[+|++[g|O]].
2251 /// The possible error codes are defined by TInterpreter::EErrorCode.
2252 /// If check is true it will only check if filename exists and is
2253 /// readable.
2254 /// Returns 0 on successful loading and -1 in case filename does not
2255 /// exist or in case of error.
2256 
2257 Int_t TROOT::LoadMacro(const char *filename, int *error, Bool_t check)
2258 {
2259  Int_t err = -1;
2260  Int_t lerr, *terr;
2261  if (error)
2262  terr = error;
2263  else
2264  terr = &lerr;
2265 
2266  if (fInterpreter) {
2267  TString aclicMode;
2268  TString arguments;
2269  TString io;
2270  TString fname = gSystem->SplitAclicMode(filename, aclicMode, arguments, io);
2271 
2272  if (arguments.Length()) {
2273  Warning("LoadMacro", "argument(%s) ignored in %s", arguments.Data(), GetMacroPath());
2274  }
2275  char *mac = gSystem->Which(GetMacroPath(), fname, kReadPermission);
2276  if (!mac) {
2277  if (!check)
2278  Error("LoadMacro", "macro %s not found in path %s", fname.Data(), GetMacroPath());
2279  *terr = TInterpreter::kFatal;
2280  } else {
2281  err = 0;
2282  if (!check) {
2283  fname = mac;
2284  fname += aclicMode;
2285  fname += io;
2286  gInterpreter->LoadMacro(fname.Data(), (TInterpreter::EErrorCode*)terr);
2287  if (*terr)
2288  err = -1;
2289  }
2290  }
2291  delete [] mac;
2292  }
2293  return err;
2294 }
2295 
2296 ////////////////////////////////////////////////////////////////////////////////
2297 /// Execute a macro in the interpreter. Equivalent to the command line
2298 /// command ".x filename". If the filename has "+" or "++" appended
2299 /// the macro will be compiled by ACLiC. The filename must have the format:
2300 /// [path/]macro.C[+|++[g|O]][(args)].
2301 /// The possible error codes are defined by TInterpreter::EErrorCode.
2302 /// If padUpdate is true (default) update the current pad.
2303 /// Returns the macro return value.
2304 
2305 Long_t TROOT::Macro(const char *filename, Int_t *error, Bool_t padUpdate)
2306 {
2307  Long_t result = 0;
2308 
2309  if (fInterpreter) {
2310  TString aclicMode;
2311  TString arguments;
2312  TString io;
2313  TString fname = gSystem->SplitAclicMode(filename, aclicMode, arguments, io);
2314 
2315  char *mac = gSystem->Which(GetMacroPath(), fname, kReadPermission);
2316  if (!mac) {
2317  Error("Macro", "macro %s not found in path %s", fname.Data(), GetMacroPath());
2318  if (error)
2319  *error = TInterpreter::kFatal;
2320  } else {
2321  fname = mac;
2322  fname += aclicMode;
2323  fname += arguments;
2324  fname += io;
2325  result = gInterpreter->ExecuteMacro(fname, (TInterpreter::EErrorCode*)error);
2326  }
2327  delete [] mac;
2328 
2329  if (padUpdate && gPad)
2330  gPad->Update();
2331  }
2332 
2333  return result;
2334 }
2335 
2336 ////////////////////////////////////////////////////////////////////////////////
2337 /// Process message id called by obj.
2338 
2339 void TROOT::Message(Int_t id, const TObject *obj)
2340 {
2341  TIter next(fMessageHandlers);
2342  TMessageHandler *mh;
2343  while ((mh = (TMessageHandler*)next())) {
2344  mh->HandleMessage(id,obj);
2345  }
2346 }
2347 
2348 ////////////////////////////////////////////////////////////////////////////////
2349 /// Process interpreter command via TApplication::ProcessLine().
2350 /// On Win32 the line will be processed asynchronously by sending
2351 /// it to the CINT interpreter thread. For explicit synchronous processing
2352 /// use ProcessLineSync(). On non-Win32 platforms there is no difference
2353 /// between ProcessLine() and ProcessLineSync().
2354 /// The possible error codes are defined by TInterpreter::EErrorCode. In
2355 /// particular, error will equal to TInterpreter::kProcessing until the
2356 /// CINT interpreted thread has finished executing the line.
2357 /// Returns the result of the command, cast to a Long_t.
2358 
2359 Long_t TROOT::ProcessLine(const char *line, Int_t *error)
2360 {
2361  TString sline = line;
2362  sline = sline.Strip(TString::kBoth);
2363 
2364  if (!fApplication.load())
2365  TApplication::CreateApplication();
2366 
2367  return (*fApplication).ProcessLine(sline, kFALSE, error);
2368 }
2369 
2370 ////////////////////////////////////////////////////////////////////////////////
2371 /// Process interpreter command via TApplication::ProcessLine().
2372 /// On Win32 the line will be processed synchronously (i.e. it will
2373 /// only return when the CINT interpreter thread has finished executing
2374 /// the line). On non-Win32 platforms there is no difference between
2375 /// ProcessLine() and ProcessLineSync().
2376 /// The possible error codes are defined by TInterpreter::EErrorCode.
2377 /// Returns the result of the command, cast to a Long_t.
2378 
2379 Long_t TROOT::ProcessLineSync(const char *line, Int_t *error)
2380 {
2381  TString sline = line;
2382  sline = sline.Strip(TString::kBoth);
2383 
2384  if (!fApplication.load())
2385  TApplication::CreateApplication();
2386 
2387  return (*fApplication).ProcessLine(sline, kTRUE, error);
2388 }
2389 
2390 ////////////////////////////////////////////////////////////////////////////////
2391 /// Process interpreter command directly via CINT interpreter.
2392 /// Only executable statements are allowed (no variable declarations),
2393 /// In all other cases use TROOT::ProcessLine().
2394 /// The possible error codes are defined by TInterpreter::EErrorCode.
2395 
2396 Long_t TROOT::ProcessLineFast(const char *line, Int_t *error)
2397 {
2398  TString sline = line;
2399  sline = sline.Strip(TString::kBoth);
2400 
2401  if (!fApplication.load())
2402  TApplication::CreateApplication();
2403 
2404  Long_t result = 0;
2405 
2406  if (fInterpreter) {
2407  TInterpreter::EErrorCode *code = (TInterpreter::EErrorCode*)error;
2408  result = gInterpreter->Calc(sline, code);
2409  }
2410 
2411  return result;
2412 }
2413 
2414 ////////////////////////////////////////////////////////////////////////////////
2415 /// Read Git commit information and branch name from the
2416 /// etc/gitinfo.txt file.
2417 
2418 void TROOT::ReadGitInfo()
2419 {
2420 #ifdef ROOT_GIT_COMMIT
2421  fGitCommit = ROOT_GIT_COMMIT;
2422 #endif
2423 #ifdef ROOT_GIT_BRANCH
2424  fGitBranch = ROOT_GIT_BRANCH;
2425 #endif
2426 
2427  TString gitinfo = "gitinfo.txt";
2428  char *filename = gSystem->ConcatFileName(TROOT::GetEtcDir(), gitinfo);
2429 
2430  FILE *fp = fopen(filename, "r");
2431  if (fp) {
2432  TString s;
2433  // read branch name
2434  s.Gets(fp);
2435  fGitBranch = s;
2436  // read commit SHA1
2437  s.Gets(fp);
2438  fGitCommit = s;
2439  // read date/time make was run
2440  s.Gets(fp);
2441  fGitDate = s;
2442  fclose(fp);
2443  }
2444  delete [] filename;
2445 }
2446 
2447 Bool_t &GetReadingObject() {
2448  TTHREAD_TLS(Bool_t) fgReadingObject = false;
2449  return fgReadingObject;
2450 }
2451 
2452 ////////////////////////////////////////////////////////////////////////////////
2453 /// Deprecated (will be removed in next release).
2454 
2455 Bool_t TROOT::ReadingObject() const
2456 {
2457  return GetReadingObject();
2458 }
2459 
2460 void TROOT::SetReadingObject(Bool_t flag)
2461 {
2462  GetReadingObject() = flag;
2463 }
2464 
2465 
2466 ////////////////////////////////////////////////////////////////////////////////
2467 /// Return date/time make was run.
2468 
2469 const char *TROOT::GetGitDate()
2470 {
2471  if (fGitDate == "") {
2472  Int_t iday,imonth,iyear, ihour, imin;
2473  static const char *months[] = { "Jan", "Feb", "Mar", "Apr", "May", "Jun",
2474  "Jul", "Aug", "Sep", "Oct", "Nov", "Dec" };
2475  Int_t idate = gROOT->GetBuiltDate();
2476  Int_t itime = gROOT->GetBuiltTime();
2477  iday = idate%100;
2478  imonth = (idate/100)%100;
2479  iyear = idate/10000;
2480  ihour = itime/100;
2481  imin = itime%100;
2482  fGitDate.Form("%s %02d %4d, %02d:%02d:00", months[imonth-1], iday, iyear, ihour, imin);
2483  }
2484  return fGitDate;
2485 }
2486 
2487 ////////////////////////////////////////////////////////////////////////////////
2488 /// Recursively remove this object from the list of Cleanups.
2489 /// Typically RecursiveRemove is implemented by classes that can contain
2490 /// mulitple references to a same object or shared ownership of the object
2491 /// with others.
2492 
2493 void TROOT::RecursiveRemove(TObject *obj)
2494 {
2495  R__READ_LOCKGUARD(ROOT::gCoreMutex);
2496 
2497  fCleanups->RecursiveRemove(obj);
2498 }
2499 
2500 ////////////////////////////////////////////////////////////////////////////////
2501 /// Refresh all browsers. Call this method when some command line
2502 /// command or script has changed the browser contents. Not needed
2503 /// for objects that have the kMustCleanup bit set. Most useful to
2504 /// update browsers that show the file system or other objects external
2505 /// to the running ROOT session.
2506 
2507 void TROOT::RefreshBrowsers()
2508 {
2509  TIter next(GetListOfBrowsers());
2510  TBrowser *b;
2511  while ((b = (TBrowser*) next()))
2512  b->SetRefreshFlag(kTRUE);
2513 }
2514 ////////////////////////////////////////////////////////////////////////////////
2515 /// Insure that the files, canvases and sockets are closed.
2516 
2517 static void CallCloseFiles()
2518 {
2519  if (TROOT::Initialized() && ROOT::Internal::gROOTLocal) {
2520  gROOT->CloseFiles();
2521  }
2522 }
2523 
2524 ////////////////////////////////////////////////////////////////////////////////
2525 /// Called by static dictionary initialization to register clang modules
2526 /// for headers. Calls TCling::RegisterModule() unless gCling
2527 /// is NULL, i.e. during startup, where the information is buffered in
2528 /// the static GetModuleHeaderInfoBuffer().
2529 
2530 void TROOT::RegisterModule(const char* modulename,
2531  const char** headers,
2532  const char** includePaths,
2533  const char* payloadCode,
2534  const char* fwdDeclCode,
2535  void (*triggerFunc)(),
2536  const TInterpreter::FwdDeclArgsToKeepCollection_t& fwdDeclsArgToSkip,
2537  const char** classesHeaders,
2538  bool hasCxxModule)
2539 {
2540 
2541  // First a side track to insure proper end of process behavior.
2542 
2543  // Register for each loaded dictionary (and thus for each library),
2544  // that we need to Close the ROOT files as soon as this library
2545  // might start being unloaded after main.
2546  //
2547  // By calling atexit here (rather than directly from within the
2548  // library) we make sure that this is not called if the library is
2549  // 'only' dlclosed.
2550 
2551  // On Ubuntu the linker strips the unused libraries. Eventhough
2552  // stressHistogram is explicitly linked against libNet, it is not
2553  // retained and thus is loaded only as needed in the middle part of
2554  // the execution. Concretely this also means that it is loaded
2555  // *after* the construction of the TApplication object and thus
2556  // after the registration (atexit) of the EndOfProcessCleanups
2557  // routine. Consequently, after the end of main, libNet is
2558  // unloaded before EndOfProcessCleanups is called. When
2559  // EndOfProcessCleanups is executed it indirectly needs the TClass
2560  // for TSocket and its search will use resources that have already
2561  // been unloaded (technically the function static in TUnixSystem's
2562  // DynamicPath and the dictionary from libNet).
2563 
2564  // Similarly, the ordering (before this commit) was broken in the
2565  // following case:
2566 
2567  // TApplication creation (EndOfProcessCleanups registration)
2568  // load UserLibrary
2569  // create TFile
2570  // Append UserObject to TFile
2571 
2572  // and after the end of main the order of execution was
2573 
2574  // unload UserLibrary
2575  // call EndOfProcessCleanups
2576  // Write the TFile
2577  // attempt to write the user object.
2578  // ....
2579 
2580  // where what we need is to have the files closen/written before
2581  // the unloading of the library.
2582 
2583  // To solve the problem we now register an atexit function for
2584  // every dictionary thus making sure there is at least one executed
2585  // before the first library tear down after main.
2586 
2587  // If atexit is called directly within a library's code, the
2588  // function will called *either* when the library is 'dlclose'd or
2589  // after then end of main (whichever comes first). We do *not*
2590  // want the files to be closed whenever a library is unloaded via
2591  // dlclose. To avoid this, we add the function (CallCloseFiles)
2592  // from the dictionary indirectly (via ROOT::RegisterModule). In
2593  // this case the function will only only be called either when
2594  // libCore is 'dlclose'd or right after the end of main.
2595 
2596  atexit(CallCloseFiles);
2597 
2598  // Now register with TCling.
2599  if (TROOT::Initialized()) {
2600  gCling->RegisterModule(modulename, headers, includePaths, payloadCode, fwdDeclCode, triggerFunc,
2601  fwdDeclsArgToSkip, classesHeaders, false, hasCxxModule);
2602  } else {
2603  GetModuleHeaderInfoBuffer().push_back(ModuleHeaderInfo_t(modulename, headers, includePaths, payloadCode,
2604  fwdDeclCode, triggerFunc, fwdDeclsArgToSkip,
2605  classesHeaders, hasCxxModule));
2606  }
2607 }
2608 
2609 ////////////////////////////////////////////////////////////////////////////////
2610 /// Remove an object from the in-memory list.
2611 /// Since TROOT is global resource, this is lock protected.
2612 
2613 TObject *TROOT::Remove(TObject* obj)
2614 {
2615  R__LOCKGUARD(gROOTMutex);
2616  return TDirectory::Remove(obj);
2617 }
2618 
2619 ////////////////////////////////////////////////////////////////////////////////
2620 /// Remove a class from the list and map of classes.
2621 /// This routine is deprecated, use TClass::RemoveClass directly.
2622 
2623 void TROOT::RemoveClass(TClass *oldcl)
2624 {
2625  TClass::RemoveClass(oldcl);
2626 }
2627 
2628 ////////////////////////////////////////////////////////////////////////////////
2629 /// Delete all global interpreter objects created since the last call to Reset
2630 ///
2631 /// If option="a" is set reset to startup context (i.e. unload also
2632 /// all loaded files, classes, structs, typedefs, etc.).
2633 ///
2634 /// This function is typically used at the beginning (or end) of an unnamed macro
2635 /// to clean the environment.
2636 ///
2637 /// IMPORTANT WARNING:
2638 /// Do not use this call from within any function (neither compiled nor
2639 /// interpreted. This should only be used from a unnamed macro
2640 /// (which starts with a { (curly braces) ). For example, using TROOT::Reset
2641 /// from within an interpreted function will lead to the unloading of the
2642 /// dictionary and source file, including the one defining the function being
2643 /// executed.
2644 ///
2645 
2646 void TROOT::Reset(Option_t *option)
2647 {
2648  if (IsExecutingMacro()) return; //True when TMacro::Exec runs
2649  if (fInterpreter) {
2650  if (!strncmp(option, "a", 1)) {
2651  fInterpreter->Reset();
2652  fInterpreter->SaveContext();
2653  } else
2654  gInterpreter->ResetGlobals();
2655 
2656  if (fGlobals) fGlobals->Unload();
2657  if (fGlobalFunctions) fGlobalFunctions->Unload();
2658 
2659  SaveContext();
2660  }
2661 }
2662 
2663 ////////////////////////////////////////////////////////////////////////////////
2664 /// Save the current interpreter context.
2665 
2666 void TROOT::SaveContext()
2667 {
2668  if (fInterpreter)
2669  gInterpreter->SaveGlobalsContext();
2670 }
2671 
2672 ////////////////////////////////////////////////////////////////////////////////
2673 /// Set the default graphical cut class name for the graphics editor
2674 /// By default the graphics editor creates an instance of a class TCutG.
2675 /// This function may be called to specify a different class that MUST
2676 /// derive from TCutG
2677 
2678 void TROOT::SetCutClassName(const char *name)
2679 {
2680  if (!name) {
2681  Error("SetCutClassName","Invalid class name");
2682  return;
2683  }
2684  TClass *cl = TClass::GetClass(name);
2685  if (!cl) {
2686  Error("SetCutClassName","Unknown class:%s",name);
2687  return;
2688  }
2689  if (!cl->InheritsFrom("TCutG")) {
2690  Error("SetCutClassName","Class:%s does not derive from TCutG",name);
2691  return;
2692  }
2693  fCutClassName = name;
2694 }
2695 
2696 ////////////////////////////////////////////////////////////////////////////////
2697 /// Set editor mode
2698 
2699 void TROOT::SetEditorMode(const char *mode)
2700 {
2701  fEditorMode = 0;
2702  if (!mode[0]) return;
2703  if (!strcmp(mode,"Arc")) {fEditorMode = kArc; return;}
2704  if (!strcmp(mode,"Line")) {fEditorMode = kLine; return;}
2705  if (!strcmp(mode,"Arrow")) {fEditorMode = kArrow; return;}
2706  if (!strcmp(mode,"Button")) {fEditorMode = kButton; return;}
2707  if (!strcmp(mode,"Diamond")) {fEditorMode = kDiamond; return;}
2708  if (!strcmp(mode,"Ellipse")) {fEditorMode = kEllipse; return;}
2709  if (!strcmp(mode,"Pad")) {fEditorMode = kPad; return;}
2710  if (!strcmp(mode,"Pave")) {fEditorMode = kPave; return;}
2711  if (!strcmp(mode,"PaveLabel")){fEditorMode = kPaveLabel; return;}
2712  if (!strcmp(mode,"PaveText")) {fEditorMode = kPaveText; return;}
2713  if (!strcmp(mode,"PavesText")){fEditorMode = kPavesText; return;}
2714  if (!strcmp(mode,"PolyLine")) {fEditorMode = kPolyLine; return;}
2715  if (!strcmp(mode,"CurlyLine")){fEditorMode = kCurlyLine; return;}
2716  if (!strcmp(mode,"CurlyArc")) {fEditorMode = kCurlyArc; return;}
2717  if (!strcmp(mode,"Text")) {fEditorMode = kText; return;}
2718  if (!strcmp(mode,"Marker")) {fEditorMode = kMarker; return;}
2719  if (!strcmp(mode,"CutG")) {fEditorMode = kCutG; return;}
2720 }
2721 
2722 ////////////////////////////////////////////////////////////////////////////////
2723 /// Change current style to style with name stylename
2724 
2725 void TROOT::SetStyle(const char *stylename)
2726 {
2727  TString style_name = stylename;
2728 
2729  TStyle *style = GetStyle(style_name);
2730  if (style) style->cd();
2731  else Error("SetStyle","Unknown style:%s",style_name.Data());
2732 }
2733 
2734 
2735 //-------- Static Member Functions ---------------------------------------------
2736 
2737 
2738 ////////////////////////////////////////////////////////////////////////////////
2739 /// Decrease the indentation level for ls().
2740 
2741 Int_t TROOT::DecreaseDirLevel()
2742 {
2743  return --fgDirLevel;
2744 }
2745 
2746 ////////////////////////////////////////////////////////////////////////////////
2747 ///return directory level
2748 
2749 Int_t TROOT::GetDirLevel()
2750 {
2751  return fgDirLevel;
2752 }
2753 
2754 ////////////////////////////////////////////////////////////////////////////////
2755 /// Get macro search path. Static utility function.
2756 
2757 const char *TROOT::GetMacroPath()
2758 {
2759  TString &macroPath = ROOT::GetMacroPath();
2760 
2761  if (macroPath.Length() == 0) {
2762  macroPath = gEnv->GetValue("Root.MacroPath", (char*)0);
2763 #if defined(R__WIN32)
2764  macroPath.ReplaceAll("; ", ";");
2765 #else
2766  macroPath.ReplaceAll(": ", ":");
2767 #endif
2768  if (macroPath.Length() == 0)
2769 #if !defined(R__WIN32)
2770  macroPath = ".:" + TROOT::GetMacroDir();
2771 #else
2772  macroPath = ".;" + TROOT::GetMacroDir();
2773 #endif
2774  }
2775 
2776  return macroPath;
2777 }
2778 
2779 ////////////////////////////////////////////////////////////////////////////////
2780 /// Set or extend the macro search path. Static utility function.
2781 /// If newpath=0 or "" reset to value specified in the rootrc file.
2782 
2783 void TROOT::SetMacroPath(const char *newpath)
2784 {
2785  TString &macroPath = ROOT::GetMacroPath();
2786 
2787  if (!newpath || !*newpath)
2788  macroPath = "";
2789  else
2790  macroPath = newpath;
2791 }
2792 
2793 ////////////////////////////////////////////////////////////////////////////////
2794 /// \brief Specify where web graphics shall be rendered
2795 ///
2796 /// The input parameter `webdisplay` defines where web graphics is rendered.
2797 /// `webdisplay` parameter may contain:
2798 ///
2799 /// - "off": turns off the web display and comes back to normal graphics in
2800 /// interactive mode.
2801 /// - "batch": turns the web display in batch mode. It can be prepended with
2802 /// another string which is considered as the new current web display.
2803 /// - "nobatch": turns the web display in interactive mode. It can be
2804 /// prepended with another string which is considered as the new current web display.
2805 ///
2806 /// If the option "off" is not set, this method turns the normal graphics to
2807 /// "Batch" to avoid the loading of local graphics libraries.
2808 
2809 void TROOT::SetWebDisplay(const char *webdisplay)
2810 {
2811  const char *wd = webdisplay;
2812  if (!wd)
2813  wd = "";
2814 
2815  if (!strcmp(wd, "off")) {
2816  fIsWebDisplay = kFALSE;
2817  fIsWebDisplayBatch = kFALSE;
2818  fWebDisplay = "";
2819  } else {
2820  fIsWebDisplay = kTRUE;
2821  if (!strncmp(wd, "batch", 5)) {
2822  fIsWebDisplayBatch = kTRUE;
2823  wd += 5;
2824  } else if (!strncmp(wd, "nobatch", 7)) {
2825  fIsWebDisplayBatch = kFALSE;
2826  wd += 7;
2827  } else {
2828  fIsWebDisplayBatch = kFALSE;
2829  }
2830  fWebDisplay = wd;
2831  }
2832 }
2833 
2834 ////////////////////////////////////////////////////////////////////////////////
2835 /// Increase the indentation level for ls().
2836 
2837 Int_t TROOT::IncreaseDirLevel()
2838 {
2839  return ++fgDirLevel;
2840 }
2841 
2842 ////////////////////////////////////////////////////////////////////////////////
2843 /// Functions used by ls() to indent an object hierarchy.
2844 
2845 void TROOT::IndentLevel()
2846 {
2847  for (int i = 0; i < fgDirLevel; i++) std::cout.put(' ');
2848 }
2849 
2850 ////////////////////////////////////////////////////////////////////////////////
2851 /// Initialize ROOT explicitly.
2852 
2853 void TROOT::Initialize() {
2854  (void) gROOT;
2855 }
2856 
2857 ////////////////////////////////////////////////////////////////////////////////
2858 /// Return kTRUE if the TROOT object has been initialized.
2859 
2860 Bool_t TROOT::Initialized()
2861 {
2862  return fgRootInit;
2863 }
2864 
2865 ////////////////////////////////////////////////////////////////////////////////
2866 /// Return kTRUE if the memory leak checker is on.
2867 
2868 Bool_t TROOT::MemCheck()
2869 {
2870  return fgMemCheck;
2871 }
2872 
2873 ////////////////////////////////////////////////////////////////////////////////
2874 /// Return Indentation level for ls().
2875 
2876 void TROOT::SetDirLevel(Int_t level)
2877 {
2878  fgDirLevel = level;
2879 }
2880 
2881 ////////////////////////////////////////////////////////////////////////////////
2882 /// Convert version code to an integer, i.e. 331527 -> 51507.
2883 
2884 Int_t TROOT::ConvertVersionCode2Int(Int_t code)
2885 {
2886  return 10000*(code>>16) + 100*((code&65280)>>8) + (code&255);
2887 }
2888 
2889 ////////////////////////////////////////////////////////////////////////////////
2890 /// Convert version as an integer to version code as used in RVersion.h.
2891 
2892 Int_t TROOT::ConvertVersionInt2Code(Int_t v)
2893 {
2894  int a = v/10000;
2895  int b = (v - a*10000)/100;
2896  int c = v - a*10000 - b*100;
2897  return (a << 16) + (b << 8) + c;
2898 }
2899 
2900 ////////////////////////////////////////////////////////////////////////////////
2901 /// Return ROOT version code as defined in RVersion.h.
2902 
2903 Int_t TROOT::RootVersionCode()
2904 {
2905  return ROOT_VERSION_CODE;
2906 }
2907 ////////////////////////////////////////////////////////////////////////////////
2908 /// Provide command line arguments to the interpreter construction.
2909 /// These arguments are added to the existing flags (e.g. `-DNDEBUG`).
2910 /// They are evaluated once per process, at the time where TROOT (and thus
2911 /// TInterpreter) is constructed.
2912 /// Returns the new flags.
2913 
2914 const std::vector<std::string> &TROOT::AddExtraInterpreterArgs(const std::vector<std::string> &args) {
2915  static std::vector<std::string> sArgs = {};
2916  sArgs.insert(sArgs.begin(), args.begin(), args.end());
2917  return sArgs;
2918 }
2919 
2920 ////////////////////////////////////////////////////////////////////////////////
2921 /// INTERNAL function!
2922 /// Used by rootcling to inject interpreter arguments through a C-interface layer.
2923 
2924 const char**& TROOT::GetExtraInterpreterArgs() {
2925  static const char** extraInterpArgs = 0;
2926  return extraInterpArgs;
2927 }
2928 
2929 ////////////////////////////////////////////////////////////////////////////////
2930 
2931 #ifdef ROOTPREFIX
2932 static Bool_t IgnorePrefix() {
2933  static Bool_t ignorePrefix = gSystem->Getenv("ROOTIGNOREPREFIX");
2934  return ignorePrefix;
2935 }
2936 #endif
2937 
2938 ////////////////////////////////////////////////////////////////////////////////
2939 /// Get the rootsys directory in the installation. Static utility function.
2940 
2941 const TString& TROOT::GetRootSys() {
2942  // Avoid returning a reference to a temporary because of the conversion
2943  // between std::string and TString.
2944  const static TString rootsys = ROOT::FoundationUtils::GetRootSys();
2945  return rootsys;
2946 }
2947 
2948 ////////////////////////////////////////////////////////////////////////////////
2949 /// Get the binary directory in the installation. Static utility function.
2950 
2951 const TString& TROOT::GetBinDir() {
2952 #ifdef ROOTBINDIR
2953  if (IgnorePrefix()) {
2954 #endif
2955  static TString rootbindir;
2956  if (rootbindir.IsNull()) {
2957  rootbindir = "bin";
2958  gSystem->PrependPathName(GetRootSys(), rootbindir);
2959  }
2960  return rootbindir;
2961 #ifdef ROOTBINDIR
2962  } else {
2963  const static TString rootbindir = ROOTBINDIR;
2964  return rootbindir;
2965  }
2966 #endif
2967 }
2968 
2969 ////////////////////////////////////////////////////////////////////////////////
2970 /// Get the library directory in the installation. Static utility function.
2971 
2972 const TString& TROOT::GetLibDir() {
2973 #ifdef ROOTLIBDIR
2974  if (IgnorePrefix()) {
2975 #endif
2976  static TString rootlibdir;
2977  if (rootlibdir.IsNull()) {
2978  rootlibdir = "lib";
2979  gSystem->PrependPathName(GetRootSys(), rootlibdir);
2980  }
2981  return rootlibdir;
2982 #ifdef ROOTLIBDIR
2983  } else {
2984  const static TString rootlibdir = ROOTLIBDIR;
2985  return rootlibdir;
2986  }
2987 #endif
2988 }
2989 
2990 ////////////////////////////////////////////////////////////////////////////////
2991 /// Get the include directory in the installation. Static utility function.
2992 
2993 const TString& TROOT::GetIncludeDir() {
2994  // Avoid returning a reference to a temporary because of the conversion
2995  // between std::string and TString.
2996  const static TString includedir = ROOT::FoundationUtils::GetIncludeDir();
2997  return includedir;
2998 }
2999 
3000 ////////////////////////////////////////////////////////////////////////////////
3001 /// Get the sysconfig directory in the installation. Static utility function.
3002 
3003 const TString& TROOT::GetEtcDir() {
3004  // Avoid returning a reference to a temporary because of the conversion
3005  // between std::string and TString.
3006  const static TString etcdir = ROOT::FoundationUtils::GetEtcDir();
3007  return etcdir;
3008 }
3009 
3010 ////////////////////////////////////////////////////////////////////////////////
3011 /// Get the data directory in the installation. Static utility function.
3012 
3013 const TString& TROOT::GetDataDir() {
3014 #ifdef ROOTDATADIR
3015  if (IgnorePrefix()) {
3016 #endif
3017  return GetRootSys();
3018 #ifdef ROOTDATADIR
3019  } else {
3020  const static TString rootdatadir = ROOTDATADIR;
3021  return rootdatadir;
3022  }
3023 #endif
3024 }
3025 
3026 ////////////////////////////////////////////////////////////////////////////////
3027 /// Get the documentation directory in the installation. Static utility function.
3028 
3029 const TString& TROOT::GetDocDir() {
3030 #ifdef ROOTDOCDIR
3031  if (IgnorePrefix()) {
3032 #endif
3033  return GetRootSys();
3034 #ifdef ROOTDOCDIR
3035  } else {
3036  const static TString rootdocdir = ROOTDOCDIR;
3037  return rootdocdir;
3038  }
3039 #endif
3040 }
3041 
3042 ////////////////////////////////////////////////////////////////////////////////
3043 /// Get the macro directory in the installation. Static utility function.
3044 
3045 const TString& TROOT::GetMacroDir() {
3046 #ifdef ROOTMACRODIR
3047  if (IgnorePrefix()) {
3048 #endif
3049  static TString rootmacrodir;
3050  if (rootmacrodir.IsNull()) {
3051  rootmacrodir = "macros";
3052  gSystem->PrependPathName(GetRootSys(), rootmacrodir);
3053  }
3054  return rootmacrodir;
3055 #ifdef ROOTMACRODIR
3056  } else {
3057  const static TString rootmacrodir = ROOTMACRODIR;
3058  return rootmacrodir;
3059  }
3060 #endif
3061 }
3062 
3063 ////////////////////////////////////////////////////////////////////////////////
3064 /// Get the tutorials directory in the installation. Static utility function.
3065 
3066 const TString& TROOT::GetTutorialDir() {
3067 #ifdef ROOTTUTDIR
3068  if (IgnorePrefix()) {
3069 #endif
3070  static TString roottutdir;
3071  if (roottutdir.IsNull()) {
3072  roottutdir = "tutorials";
3073  gSystem->PrependPathName(GetRootSys(), roottutdir);
3074  }
3075  return roottutdir;
3076 #ifdef ROOTTUTDIR
3077  } else {
3078  const static TString roottutdir = ROOTTUTDIR;
3079  return roottutdir;
3080  }
3081 #endif
3082 }
3083 
3084 ////////////////////////////////////////////////////////////////////////////////
3085 /// Shut down ROOT.
3086 
3087 void TROOT::ShutDown()
3088 {
3089  if (gROOT)
3090  gROOT->EndOfProcessCleanups();
3091  else if (gInterpreter)
3092  gInterpreter->ShutDown();
3093 }
3094 
3095 ////////////////////////////////////////////////////////////////////////////////
3096 /// Get the source directory in the installation. Static utility function.
3097 
3098 const TString& TROOT::GetSourceDir() {
3099 #ifdef ROOTSRCDIR
3100  if (IgnorePrefix()) {
3101 #endif
3102  static TString rootsrcdir;
3103  if (rootsrcdir.IsNull()) {
3104  rootsrcdir = "src";
3105  gSystem->PrependPathName(GetRootSys(), rootsrcdir);
3106  }
3107  return rootsrcdir;
3108 #ifdef ROOTSRCDIR
3109  } else {
3110  const static TString rootsrcdir = ROOTSRCDIR;
3111  return rootsrcdir;
3112  }
3113 #endif
3114 }
3115 
3116 ////////////////////////////////////////////////////////////////////////////////
3117 /// Get the icon path in the installation. Static utility function.
3118 
3119 const TString& TROOT::GetIconPath() {
3120 #ifdef ROOTICONPATH
3121  if (IgnorePrefix()) {
3122 #endif
3123  static TString rooticonpath;
3124  if (rooticonpath.IsNull()) {
3125  rooticonpath = "icons";
3126  gSystem->PrependPathName(GetRootSys(), rooticonpath);
3127  }
3128  return rooticonpath;
3129 #ifdef ROOTICONPATH
3130  } else {
3131  const static TString rooticonpath = ROOTICONPATH;
3132  return rooticonpath;
3133  }
3134 #endif
3135 }
3136 
3137 ////////////////////////////////////////////////////////////////////////////////
3138 /// Get the fonts directory in the installation. Static utility function.
3139 
3140 const TString& TROOT::GetTTFFontDir() {
3141 #ifdef TTFFONTDIR
3142  if (IgnorePrefix()) {
3143 #endif
3144  static TString ttffontdir;
3145  if (ttffontdir.IsNull()) {
3146  ttffontdir = "fonts";
3147  gSystem->PrependPathName(GetRootSys(), ttffontdir);
3148  }
3149  return ttffontdir;
3150 #ifdef TTFFONTDIR
3151  } else {
3152  const static TString ttffontdir = TTFFONTDIR;
3153  return ttffontdir;
3154  }
3155 #endif
3156 }
3157 
3158 ////////////////////////////////////////////////////////////////////////////////
3159 /// Get the tutorials directory in the installation. Static utility function.
3160 /// Backward compatibility function - do not use for new code
3161 
3162 const char *TROOT::GetTutorialsDir() {
3163  return GetTutorialDir();
3164 }