Logo ROOT  
Reference Guide
TMemStatMng.cxx
Go to the documentation of this file.
1// @(#)root/memstat:$Id$
2// Author: Anar Manafov (A.Manafov@gsi.de) 2008-03-02
3
4/*************************************************************************
5* Copyright (C) 1995-2010, 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// STD
12#include <cstdlib>
13// ROOT
14#include "TSystem.h"
15#include "TError.h"
16#include "Riostream.h"
17#include "TObject.h"
18#include "TFile.h"
19#include "TTree.h"
20#include "TArrayL64.h"
21#include "TH1.h"
22#include "TMD5.h"
23#include "TMath.h"
24// Memstat
25#include "TMemStatBacktrace.h"
26#include "TMemStatMng.h"
27
28using namespace Memstat;
29
31
32TMemStatMng* TMemStatMng::fgInstance = NULL;
33
34//****************************************************************************//
35//
36//****************************************************************************//
37
38TMemStatMng::TMemStatMng():
39 TObject(),
40#if !defined(__APPLE__)
41 fPreviousMallocHook(TMemStatHook::GetMallocHook()),
42 fPreviousFreeHook(TMemStatHook::GetFreeHook()),
43#endif
44 fDumpFile(NULL),
45 fDumpTree(NULL),
46 fUseGNUBuiltinBacktrace(kFALSE),
47 fBeginTime(0),
48 fPos(0),
49 fTimems(0),
50 fNBytes(0),
51 fBtID(0),
52 fMaxCalls(5000000),
53 fBufferSize(10000),
54 fBufN(0),
55 fBufPos(0),
56 fBufTimems(0),
57 fBufNBytes(0),
58 fBufBtID(0),
59 fIndex(0),
60 fMustWrite(0),
61 fFAddrsList(0),
62 fHbtids(0),
63 fBTCount(0),
64 fBTIDCount(0),
65 fSysInfo(0)
66{
67 // Default constructor
68}
69
70////////////////////////////////////////////////////////////////////////////////
71///Initialize MemStat manager - used only by instance method
72
74{
76
77 fDumpFile = new TFile(Form("memstat_%d.root", gSystem->GetPid()), "recreate");
78 Int_t opt = 200000;
79 if(!fDumpTree) {
80 fDumpTree = new TTree("T", "Memory Statistics");
81 fDumpTree->Branch("pos", &fPos, "pos/l", opt);
82 fDumpTree->Branch("time", &fTimems, "time/I", opt);
83 fDumpTree->Branch("nbytes", &fNBytes, "nbytes/I", opt);
84 fDumpTree->Branch("btid", &fBtID, "btid/I", opt);
85 }
86
87 fBTCount = 0;
88
89 fBTIDCount = 0;
90
91 fFAddrsList = new TObjArray();
93 fFAddrsList->SetName("FAddrsList");
94
95 fHbtids = new TH1I("btids", "table of btids", 10000, 0, 1); //where fHbtids is a member of the manager class
97 // save the histogram and the TObjArray to the tree header
100 // save the system info to a tree header
101 std::string sSysInfo(gSystem->GetBuildNode());
102 sSysInfo += " | ";
103 sSysInfo += gSystem->GetBuildCompilerVersion();
104 sSysInfo += " | ";
105 sSysInfo += gSystem->GetFlagsDebug();
106 sSysInfo += " ";
107 sSysInfo += gSystem->GetFlagsOpt();
108 fSysInfo = new TNamed("SysInfo", sSysInfo.c_str());
109
111 fDumpTree->SetAutoSave(10000000);
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// GetInstance - a static function
116/// Initialize a singleton of MemStat manager
117
119{
120 if(!fgInstance) {
122 fgInstance->Init();
123 }
124 return fgInstance;
125}
126
127////////////////////////////////////////////////////////////////////////////////
128/// Close - a static function
129/// This method stops the manager,
130/// flashes all the buffered data and closes the output tree.
131
133{
134 // TODO: This is a temporary solution until we find a properalgorithm for SaveData
135 //fgInstance->fDumpFile->WriteObject(fgInstance->fFAddrsList, "FAddrsList");
136
137/* std::ofstream f("mem_stat_debug.txt");
138 int *btids = fgInstance->fHbtids->GetArray();
139 if( !btids )
140 return;
141 int btid(1);
142 int count(0);
143 bool bStop(false);
144 int count_empty(0);
145 for (int i = 0; i < fgInstance->fBTChecksums.size(); ++i)
146 {
147 if (bStop)
148 break;
149 count = btids[btid-1];
150 f << "++++++++++++++++++++++++\n";
151 f << "BTID: " << btid << "\n";
152 if ( count <= 0 )
153 ++count_empty;
154 for (int j = btid+1; j <= (btid+count); ++j )
155 {
156 TNamed *nm = (TNamed*)fgInstance->fFAddrsList->At(btids[j]);
157 if( !nm )
158 {
159 f << "Bad ID" << std::endl;
160 bStop = true;
161 }
162 f << "-------> " << nm->GetTitle() << "\n";
163 }
164 btid = btid + count + 1;
165 }
166 f.close();
167 ::Info("TMemStatMng::Close", "btids without a stack %d\n", count_empty);
168*/
169
170 // to be documented
175
176 ::Info("TMemStatMng::Close", "Tree saved to file %s\n", fgInstance->fDumpFile->GetName());
177 ::Info("TMemStatMng::Close", "Tree entries = %d, file size = %g MBytes\n", (Int_t)fgInstance->fDumpTree->GetEntries(),1e-6*Double_t(fgInstance->fDumpFile->GetEND()));
178
179 delete fgInstance->fDumpFile;
180 //fgInstance->fDumpFile->Close();
181 //delete fgInstance->fFAddrsList;
182 //delete fgInstance->fSysInfo;
183
184 delete fgInstance;
185 fgInstance = NULL;
186}
187
188////////////////////////////////////////////////////////////////////////////////
189/// if an instance is destructed - the hooks are reseted to old hooks
190
192{
193 if(this != TMemStatMng::GetInstance())
194 return;
195
196 Info("~TMemStatMng", ">>> All free/malloc calls count: %d", fBTIDCount);
197 Info("~TMemStatMng", ">>> Unique BTIDs count: %zu", fBTChecksums.size());
198
199 Disable();
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Set the maximum number of alloc/free calls to be buffered.
204///if the alloc and free are in the buffer, the corresponding entries
205///are not saved tio the Tree, reducing considerably the Tree output size
206
208{
209 fBufferSize = buffersize;
210 if (fBufferSize < 1) fBufferSize = 1;
211 fBufN = 0;
216 fIndex = new Int_t[fBufferSize];
218}
219
220////////////////////////////////////////////////////////////////////////////////
221/// Set the maximum number of new/delete registered in the output Tree.
222
224{
225 fMaxCalls = maxcalls;
226}
227
228////////////////////////////////////////////////////////////////////////////////
229/// Enable memory hooks
230
232{
233 if(this != GetInstance())
234 return;
235#if defined(__APPLE__)
236 TMemStatHook::trackZoneMalloc(MacAllocHook, MacFreeHook);
237#else
238 // set hook to our functions
241#endif
242}
243
244////////////////////////////////////////////////////////////////////////////////
245/// Disble memory hooks
246
248{
249 //FillTree();
250 if(this != GetInstance())
251 return;
252#if defined(__APPLE__)
253 TMemStatHook::untrackZoneMalloc();
254#else
255 // set hook to our functions
258#endif
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// AllocHook - a static function
263/// a special memory hook for Mac OS X memory zones.
264/// Triggered when memory is allocated.
265
266void TMemStatMng::MacAllocHook(void *ptr, size_t size)
267{
269 // Restore all old hooks
270 instance->Disable();
271
272 // Call our routine
273 instance->AddPointer(ptr, Int_t(size));
274
275 // Restore our own hooks
276 instance->Enable();
277}
278
279////////////////////////////////////////////////////////////////////////////////
280/// AllocHook - a static function
281/// a special memory hook for Mac OS X memory zones.
282/// Triggered when memory is deallocated.
283
285{
287 // Restore all old hooks
288 instance->Disable();
289
290 // Call our routine
291 instance->AddPointer(ptr, -1);
292
293 // Restore our own hooks
294 instance->Enable();
295}
296
297////////////////////////////////////////////////////////////////////////////////
298/// AllocHook - a static function
299/// A glibc memory allocation hook.
300
301void *TMemStatMng::AllocHook(size_t size, const void* /*caller*/)
302{
304 // Restore all old hooks
305 instance->Disable();
306
307 // Call recursively
308 void *result = malloc(size);
309 // Call our routine
310 instance->AddPointer(result, Int_t(size));
311 // TTimer::SingleShot(0, "TYamsMemMng", instance, "SaveData()");
312
313 // Restore our own hooks
314 instance->Enable();
315
316 return result;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// FreeHook - a static function
321/// A glibc memory deallocation hook.
322
323void TMemStatMng::FreeHook(void* ptr, const void* /*caller*/)
324{
326 // Restore all old hooks
327 instance->Disable();
328
329 // Call recursively
330 free(ptr);
331
332 // Call our routine
333 instance->AddPointer(ptr, -1);
334
335 // Restore our own hooks
336 instance->Enable();
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// An internal function, which returns a bitid for a corresponding CRC digest
341/// cache variables
342
344 void **stackPointers)
345{
346 static Int_t old_btid = -1;
347 static SCustomDigest old_digest;
348
349 Int_t ret_val = -1;
350 bool startCheck(false);
351 if(old_btid >= 0) {
352 for(int i = 0; i < g_digestSize; ++i) {
353 if(old_digest.fValue[i] != CRCdigest[i]) {
354 startCheck = true;
355 break;
356 }
357 }
358 ret_val = old_btid;
359 } else {
360 startCheck = true;
361 }
362
363 // return cached value
364 if(!startCheck)
365 return ret_val;
366
367 old_digest = SCustomDigest(CRCdigest);
368 CRCSet_t::const_iterator found = fBTChecksums.find(CRCdigest);
369
370 if(fBTChecksums.end() == found) {
371 // check the size of the BT array container
372 const int nbins = fHbtids->GetNbinsX();
373 //check that the current allocation in fHbtids is enough, otherwise expend it with
374 if(fBTCount + stackEntries + 1 >= nbins) {
375 fHbtids->SetBins(nbins * 2, 0, 1);
376 }
377
378 int *btids = fHbtids->GetArray();
379 // the first value is a number of entries in a given stack
380 btids[fBTCount++] = stackEntries;
381 ret_val = fBTCount;
382 if(stackEntries <= 0) {
383 Warning("AddPointer",
384 "A number of stack entries is equal or less than zero. For btid %d", ret_val);
385 }
386
387 // add new BT's CRC value
388 std::pair<CRCSet_t::iterator, bool> res = fBTChecksums.insert(CRCSet_t::value_type(CRCdigest, ret_val));
389 if(!res.second)
390 Error("AddPointer", "Can't added a new BTID to the container.");
391
392 // save all symbols of this BT
393 for(int i = 0; i < stackEntries; ++i) {
394 ULong_t func_addr = (ULong_t)(stackPointers[i]);
395 Int_t idx = fFAddrs.find(func_addr);
396 // check, whether it's a new symbol
397 if(idx < 0) {
398 TString strFuncAddr;
399 strFuncAddr += func_addr;
400 TString strSymbolInfo;
401 getSymbolFullInfo(stackPointers[i], &strSymbolInfo);
402
403 TNamed *nm = new TNamed(strFuncAddr, strSymbolInfo);
405 idx = fFAddrsList->GetEntriesFast() - 1;
406 // TODO: more detailed error message...
407 if(!fFAddrs.add(func_addr, idx))
408 Error("AddPointer", "Can't add a function return address to the container");
409 }
410
411 // even if we have -1 as an index we add it to the container
412 btids[fBTCount++] = idx;
413 }
414
415 } else {
416 // reuse an existing BT
417 ret_val = found->second;
418 }
419
420 old_btid = ret_val;
421
422 return ret_val;
423}
424
425////////////////////////////////////////////////////////////////////////////////
426/// Add pointer to table.
427/// This method is called every time when any of the hooks are triggered.
428/// The memory de-/allocation information will is recorded.
429
430void TMemStatMng::AddPointer(void *ptr, Int_t size)
431{
432 void *stptr[g_BTStackLevel + 1];
433 const int stackentries = getBacktrace(stptr, g_BTStackLevel, fUseGNUBuiltinBacktrace);
434
435 // save only unique BTs
436 TMD5 md5;
437 md5.Update(reinterpret_cast<UChar_t*>(stptr), sizeof(void*) * stackentries);
438 UChar_t digest[g_digestSize];
439 md5.Final(digest);
440
441 // for Debug. A counter of all (de)allacations.
442 ++fBTIDCount;
443
444 Int_t btid(generateBTID(digest, stackentries, stptr));
445
446 if(btid <= 0)
447 Error("AddPointer", "bad BT id");
448
449 fTimeStamp.Set();
450 Double_t CurTime = fTimeStamp.AsDouble();
451 fBufTimems[fBufN] = Int_t(10000.*(CurTime - fBeginTime));
452 ULong_t ul = (ULong_t)(ptr);
453 fBufPos[fBufN] = (ULong64_t)(ul);
454 fBufNBytes[fBufN] = size;
455 fBufBtID[fBufN] = btid;
456 fBufN++;
457 if (fBufN >= fBufferSize) {
458 FillTree();
459 }
460}
461
462////////////////////////////////////////////////////////////////////////////////
463///loop on all entries in the buffer and fill the output Tree
464///entries with alloc and free in the buffer are eliminated
465
467{
468
469 //eliminate alloc/free pointing to the same location in the current buffer
471 memset(fMustWrite,0,fBufN*sizeof(Bool_t));
472 Int_t i=0,j;
473 while (i<fBufN) {
474 Int_t indi = fIndex[i];
475 Int_t indmin = indi;
476 Int_t indmax = indi;
477 j = i+1;;
478 ULong64_t pos = fBufPos[indi];
479 while (j < fBufN) {
480 Int_t indj = fIndex[j];
481 ULong64_t posj = fBufPos[indj];
482 if (posj != pos) break;
483 if (indmin > indj) indmin = indj;
484 if (indmax < indj) indmax = indj;
485 j++;
486 }
487 if (indmin == indmax) fMustWrite[indmin] = kTRUE;
488 if (fBufNBytes[indmin] == -1) fMustWrite[indmin] = kTRUE;
489 if (fBufNBytes[indmax] > 0) fMustWrite[indmax] = kTRUE;
490 i = j;
491 }
492
493 // now fill the Tree with the remaining allocs/frees
494 for (i=0;i<fBufN;i++) {
495 if (!fMustWrite[i]) continue;
496 fPos = fBufPos[i];
497 fTimems = fBufTimems[i];
498 fNBytes = fBufNBytes[i];
499 fBtID = fBufBtID[i];
500 fDumpTree->Fill();
501 }
502
503 fBufN = 0;
505}
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:43
unsigned char UChar_t
Definition: RtypesCore.h:36
const Bool_t kFALSE
Definition: RtypesCore.h:90
unsigned long ULong_t
Definition: RtypesCore.h:53
unsigned long long ULong64_t
Definition: RtypesCore.h:72
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define ClassImp(name)
Definition: Rtypes.h:361
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
Definition: TSystem.h:556
#define free
Definition: civetweb.c:1539
#define malloc
Definition: civetweb.c:1536
Int_t find(ULong_t addr)
Definition: TMemStatMng.h:41
bool add(ULong_t addr, Int_t idx)
Definition: TMemStatMng.h:36
Bool_t fUseGNUBuiltinBacktrace
Definition: TMemStatMng.h:117
Int_t generateBTID(UChar_t *CRCdigest, Int_t stackEntries, void **stackPointers)
An internal function, which returns a bitid for a corresponding CRC digest cache variables.
static void MacAllocHook(void *ptr, size_t size)
AllocHook - a static function a special memory hook for Mac OS X memory zones.
TMemStatHook::FreeHookFunc_t fPreviousFreeHook
old malloc function
Definition: TMemStatMng.h:98
static void FreeHook(void *ptr, const void *)
FreeHook - a static function A glibc memory deallocation hook.
static void MacFreeHook(void *ptr)
AllocHook - a static function a special memory hook for Mac OS X memory zones.
static TMemStatMng * fgInstance
tree to dump information
Definition: TMemStatMng.h:114
void Disable()
Disble memory hooks.
TObjArray * fFAddrsList
Definition: TMemStatMng.h:136
TMemStatFAddrContainer fFAddrs
Definition: TMemStatMng.h:135
static void * AllocHook(size_t size, const void *)
AllocHook - a static function A glibc memory allocation hook.
TTree * fDumpTree
file to dump current information
Definition: TMemStatMng.h:113
virtual ~TMemStatMng()
if an instance is destructed - the hooks are reseted to old hooks
void Init()
old free function
Definition: TMemStatMng.cxx:73
TTimeStamp fTimeStamp
Definition: TMemStatMng.h:118
ULong64_t * fBufPos
Definition: TMemStatMng.h:127
void Enable()
Enable memory hooks.
void SetBufferSize(Int_t buffersize)
Set the maximum number of alloc/free calls to be buffered.
TMemStatHook::MallocHookFunc_t fPreviousMallocHook
Definition: TMemStatMng.h:97
void AddPointer(void *ptr, Int_t size)
Add pointer to table.
static void Close()
Close - a static function This method stops the manager, flashes all the buffered data and closes the...
void SetMaxCalls(Int_t maxcalls)
Set the maximum number of new/delete registered in the output Tree.
void FillTree()
loop on all entries in the buffer and fill the output Tree entries with alloc and free in the buffer ...
static TMemStatMng * GetInstance()
GetInstance - a static function Initialize a singleton of MemStat manager.
const Int_t * GetArray() const
Definition: TArrayI.h:43
void SetName(const char *name)
Definition: TCollection.h:204
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:53
virtual Long64_t GetEND() const
Definition: TFile.h:222
1-D histogram with an int per channel (see TH1 documentation)}
Definition: TH1.h:530
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8393
virtual Int_t GetNbinsX() const
Definition: TH1.h:292
virtual void SetBins(Int_t nx, Double_t xmin, Double_t xmax)
Redefine x axis parameters.
Definition: TH1.cxx:8223
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:469
This code implements the MD5 message-digest algorithm.
Definition: TMD5.h:44
void Update(const UChar_t *buf, UInt_t len)
Update TMD5 object to reflect the concatenation of another buffer full of bytes.
Definition: TMD5.cxx:108
void Final()
MD5 finalization, ends an MD5 message-digest operation, writing the the message digest and zeroizing ...
Definition: TMD5.cxx:167
static void SetMallocHook(MallocHookFunc_t p)
SetMallocHook - a static function Set pointer to function replacing alloc function.
static void SetFreeHook(FreeHookFunc_t p)
SetFreeHook - a static function Set pointer to function replacing free function.
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
An array of TObjects.
Definition: TObjArray.h:37
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
void Add(TObject *obj)
Definition: TObjArray.h:74
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:877
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:891
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:865
Basic string class.
Definition: TString.h:131
virtual const char * GetBuildNode() const
Return the build node name.
Definition: TSystem.cxx:3875
virtual const char * GetBuildCompilerVersion() const
Return the build compiler version.
Definition: TSystem.cxx:3867
virtual int GetPid()
Get process id.
Definition: TSystem.cxx:705
virtual const char * GetFlagsDebug() const
Return the debug flags.
Definition: TSystem.cxx:3895
virtual const char * GetFlagsOpt() const
Return the optimization flags.
Definition: TSystem.cxx:3903
void Set()
Set Date/Time to current time as reported by the system.
Definition: TTimeStamp.cxx:556
Double_t AsDouble() const
Definition: TTimeStamp.h:138
A TTree represents a columnar dataset.
Definition: TTree.h:78
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4524
virtual void SetAutoSave(Long64_t autos=-300000000)
This function may be called at the start of a program to change the default value for fAutoSave (and ...
Definition: TTree.cxx:8194
virtual Long64_t GetEntries() const
Definition: TTree.h:457
virtual Long64_t AutoSave(Option_t *option="")
AutoSave tree header every fAutoSave bytes.
Definition: TTree.cxx:1484
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition: TTree.h:348
virtual TList * GetUserInfo()
Return a pointer to the list containing user objects associated to this tree.
Definition: TTree.cxx:6259
const UShort_t g_digestSize
Definition: TMemStatMng.h:53
void getSymbolFullInfo(void *_pAddr, TString *_retInfo, const char *const _seporator=" | ")
size_t getBacktrace(void **_trace, size_t _size, Bool_t _bUseGNUBuiltinBacktrace=kFALSE)
Get the backtrace _trace - array of pointers _size - maximal deepness of stack information _bUseGNUBu...
const size_t g_BTStackLevel
Definition: TMemStatDef.h:17
static Roo_reg_AGKInteg1D instance
static constexpr double nm
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Definition: TMathBase.h:362
UChar_t fValue[g_digestSize]
Definition: TMemStatMng.h:62