MemCheck.cxx

Go to the documentation of this file.
00001 // @(#)root/new:$Id: MemCheck.cxx 35583 2010-09-22 14:00:00Z rdm $
00002 // Author: D.Bertini and M.Ivanov   10/08/2000
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2001, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //****************************************************************************//
00013 //
00014 //
00015 // MemCheck is used to check the memory in ROOT based applications.
00016 //
00017 // Principe:
00018 //  A memory leak often arises whenever memory allocated through
00019 //  (new, new[]) is never returned via a corresponding (delete, delete[]).
00020 //  Redefining a special version of (operator new, operator delete) will
00021 //  allow the bookkeeping of created pointers to chunks of memory and, at
00022 //  the same time, pointers to the current function (and all of its callers
00023 //  scanning the stack) in order to trace back where the memory is not freed.
00024 //  This specific bookkeeping of pointers will be done by a kind of
00025 //  "memory info" container class TMemHashTable that will be used via
00026 //  the ROOT memory management defined inside the NewDelete.cxx file.
00027 //
00028 //  To activate the memory checker you have to set in the .rootrc file
00029 //  the resource Root.MemCheck to 1 (e.g.: Root.MemCheck: 1) and you
00030 //  have to link with libNew.so (e.g. use root-config --new --libs) or
00031 //  use rootn.exe. When all this is the case you will find at the end
00032 //  of the program execution a file "memcheck.out" in the directory
00033 //  where you started your ROOT program. Alternatively you can set
00034 //  the resource Root.MemCheckFile to the name of a file to which
00035 //  the leak information will be written. The contents of this
00036 //  "memcheck.out" file can be analyzed and transformed into printable
00037 //  text via the memprobe program (in $ROOTSYS/bin).
00038 //
00039 // (c) 2000 : Gesellschaft fuer Schwerionenforschung GmbH
00040 //            Planckstrasse, 1
00041 //            D-64291 Darmstadt
00042 //            Germany
00043 //
00044 // Created 10/08/2000 by: D.Bertini and M.Ivanov.
00045 // Based on ideas from LeakTracer by Erwin Andreasen.
00046 //
00047 // - Updated:
00048 //    Date: 12/02/2001 Adapt script to new GDB 5.0, new glibc2.2.x and gcc 2.96.
00049 //    Date: 23/10/2000 (hash mechanism speeding up the bookkeeping)
00050 //
00051 // - Documentation:
00052 //
00053 //    http://www-hades.gsi.de/~dbertini/mem.html
00054 //
00055 //****************************************************************************//
00056 
00057 #include <stdio.h>
00058 #include <string.h>
00059 #include <signal.h>
00060 #include <stdlib.h>
00061 #include "MemCheck.h"
00062 #include "TSystem.h"
00063 #include "TEnv.h"
00064 #include "TError.h"
00065 
00066 #define stack_history_size 20
00067 
00068 
00069 static TMemHashTable gMemHashTable;
00070 
00071 
00072 //****************************************************************************//
00073 //                                 Storage of Stack information
00074 //****************************************************************************//
00075 
00076 //______________________________________________________________________________
00077 void TStackInfo::Init(int stacksize, void **stackptrs)
00078 {
00079    //Initialize the stack
00080    fSize = stacksize;
00081    memcpy(&(this[1]), stackptrs, stacksize * sizeof(void *));
00082    fTotalAllocCount = fTotalAllocSize = fAllocCount = fAllocSize = 0;
00083 }
00084 
00085 //______________________________________________________________________________
00086 ULong_t TStackInfo::HashStack(unsigned int size, void **ptr)
00087 {
00088    // Hash stack information.
00089 
00090    ULong_t hash = 0;
00091    for (unsigned int i = 0; i < size; i++)
00092       hash ^= TString::Hash(&ptr[i], sizeof(void*));
00093    return hash;
00094 }
00095 
00096 //______________________________________________________________________________
00097 int TStackInfo::IsEqual(unsigned int size, void **ptr)
00098 {
00099    // Return 0 if stack information not equal otherwise return 1.
00100 
00101    if (size != fSize)
00102       return 0;
00103    void **stptr = (void **) &(this[1]);
00104    for (unsigned int i = 0; i < size; i++)
00105       if (ptr[i] != stptr[i])
00106          return 0;
00107    return 1;
00108 }
00109 
00110 
00111 //****************************************************************************//
00112 //                                   Global Stack Table
00113 //****************************************************************************//
00114 
00115 //______________________________________________________________________________
00116 void TStackTable::Init()
00117 {
00118    //Initialize table.
00119 
00120    fSize = 65536;
00121    fCount = 0;
00122    fTable = (char *) malloc(fSize);
00123    if (!fTable)
00124       _exit(1);
00125    memset(fTable, 0, fSize);
00126    fNext = fTable;
00127    //initialize hash table
00128    fHashSize = 65536;
00129    fHashTable = (TStackInfo **) malloc(sizeof(TStackInfo *) * fHashSize);
00130    memset(fHashTable, 0, sizeof(TStackInfo *) * fHashSize);
00131 }
00132 
00133 //______________________________________________________________________________
00134 void TStackTable::Expand(int newsize)
00135 {
00136    // Expand stack buffer to the new size.
00137 
00138    char *tableold = fTable;
00139    fTable = (char *) realloc(fTable, newsize);
00140    fSize = newsize;
00141    int nextindex = (char *) fNext - tableold;
00142    memset(&fTable[nextindex], 0, fSize - nextindex);
00143    fNext = (char *) (&fTable[nextindex]);
00144    //
00145    //update list
00146    TStackInfo *info = (TStackInfo *) fTable;
00147    while (((char *) info->Next() - fTable) <= nextindex) {
00148       if (info->fNextHash != 0)
00149          info->fNextHash = (TStackInfo *)
00150              & fTable[(char *) info->fNextHash - tableold];
00151       info = info->Next();
00152    }
00153    //
00154    //update hash table
00155    for (int i = 0; i < fHashSize; i++)
00156       if (fHashTable[i] != 0)
00157          fHashTable[i] =
00158              (TStackInfo *) & fTable[((char *) fHashTable[i]) - tableold];
00159    //  printf("new table %p\n",fTable);
00160 }
00161 
00162 //______________________________________________________________________________
00163 TStackInfo *TStackTable::AddInfo(int size, void **stackptrs)
00164 {
00165    // Add stack information to table.
00166 
00167    // add next stack to table
00168    TStackInfo *info = (TStackInfo *) fNext;
00169    if (((char *) info + size * sizeof(void *)
00170         + sizeof(TStackInfo) - fTable) > fSize) {
00171       //need expand
00172       Expand(2 * fSize);
00173       info = (TStackInfo *) fNext;
00174    }
00175    info->Init(size, stackptrs);
00176    info->fNextHash = 0;
00177    fNext = (char *) info->Next();
00178 
00179    //add info to hash table
00180    int hash = int(info->Hash() % fHashSize);
00181    TStackInfo *info2 = fHashTable[hash];
00182    if (info2 == 0) {
00183       fHashTable[hash] = info;
00184    } else {
00185       while (info2->fNextHash)
00186          info2 = info2->fNextHash;
00187       info2->fNextHash = info;
00188    }
00189    fCount++;
00190    return info;
00191 }
00192 
00193 //______________________________________________________________________________
00194 TStackInfo *TStackTable::FindInfo(int size, void **stackptrs)
00195 {
00196    // Try to find stack info in hash table if doesn't find it will add it.
00197 
00198    int hash = int(TStackInfo::HashStack(size, (void **) stackptrs) % fHashSize);
00199    TStackInfo *info = fHashTable[hash];
00200    if (info == 0) {
00201       info = AddInfo(size, stackptrs);
00202       //printf("f0 %p    - %d\n",info,(char*)info-fTable);
00203       return info;
00204    }
00205    while (info->IsEqual(size, stackptrs) == 0) {
00206       if (info->fNextHash == 0) {
00207          info = AddInfo(size, stackptrs);
00208          //  printf("f1 %p    - %d\n",info,(char*)info-fTable);
00209          return info;
00210       } else
00211          info = info->fNextHash;
00212    }
00213    //printf("f2  %p   - %d\n",info,(char*)info-fTable);
00214    return info;
00215 };
00216 
00217 //______________________________________________________________________________
00218 int TStackTable::GetIndex(TStackInfo * info)
00219 {
00220    //return index of info
00221    return (char *) info - fTable;
00222 }
00223 
00224 //______________________________________________________________________________
00225 TStackInfo *TStackTable::GetInfo(int index)
00226 {
00227    //return TStackInfo class corresponding to index
00228    return (TStackInfo *) & fTable[index];
00229 }
00230 
00231 
00232 Int_t        TMemHashTable::fgSize = 0;
00233 Int_t        TMemHashTable::fgAllocCount = 0;
00234 TMemTable  **TMemHashTable::fgLeak = 0;
00235 TDeleteTable TMemHashTable::fgMultDeleteTable;
00236 TStackTable  TMemHashTable::fgStackTable;
00237 
00238 
00239 static void *get_stack_pointer(int level);
00240 
00241 //______________________________________________________________________________
00242 void TMemHashTable::Init()
00243 {
00244    //Initialize the hash table
00245    fgStackTable.Init();
00246    fgSize = 65536;
00247    fgAllocCount = 0;
00248    fgLeak = (TMemTable **) malloc(sizeof(void *) * fgSize);
00249    fgMultDeleteTable.fLeaks = 0;
00250    fgMultDeleteTable.fAllocCount = 0;
00251    fgMultDeleteTable.fTableSize = 0;
00252 
00253    for (int i = 0; i < fgSize; i++) {
00254       fgLeak[i] = (TMemTable *) malloc(sizeof(TMemTable));
00255       fgLeak[i]->fAllocCount = 0;
00256       fgLeak[i]->fMemSize = 0;
00257       fgLeak[i]->fFirstFreeSpot = 0;
00258       fgLeak[i]->fTableSize = 0;
00259       fgLeak[i]->fLeaks = 0;
00260    }
00261 }
00262 
00263 //______________________________________________________________________________
00264 void TMemHashTable::RehashLeak(int newSize)
00265 {
00266    // Rehash leak pointers.
00267 
00268    if (newSize <= fgSize)
00269       return;
00270    TMemTable **newLeak = (TMemTable **) malloc(sizeof(void *) * newSize);
00271    for (int i = 0; i < newSize; i++) {
00272       //build new branches
00273       newLeak[i] = (TMemTable *) malloc(sizeof(TMemTable));
00274       newLeak[i]->fAllocCount = 0;
00275       newLeak[i]->fMemSize = 0;
00276       newLeak[i]->fFirstFreeSpot = 0;
00277       newLeak[i]->fTableSize = 0;
00278       newLeak[i]->fLeaks = 0;
00279    }
00280    for (int ib = 0; ib < fgSize; ib++) {
00281       TMemTable *branch = fgLeak[ib];
00282       for (int i = 0; i < branch->fTableSize; i++)
00283          if (branch->fLeaks[i].fAddress != 0) {
00284             int hash = int(TString::Hash(&branch->fLeaks[i].fAddress, sizeof(void*)) % newSize);
00285             TMemTable *newbranch = newLeak[hash];
00286             if (newbranch->fAllocCount >= newbranch->fTableSize) {
00287                int newTableSize =
00288                    newbranch->fTableSize ==
00289                    0 ? 16 : newbranch->fTableSize * 2;
00290                newbranch->fLeaks =
00291                    (TMemInfo *) realloc(newbranch->fLeaks,
00292                                         sizeof(TMemInfo) * newTableSize);
00293                if (!newbranch->fLeaks) {
00294                   Error("TMemHashTable::AddPointer", "realloc failure");
00295                   _exit(1);
00296                }
00297                memset(newbranch->fLeaks + newbranch->fTableSize, 0,
00298                       sizeof(TMemInfo) * (newTableSize -
00299                                           newbranch->fTableSize));
00300                newbranch->fTableSize = newTableSize;
00301             }
00302             memcpy(&newbranch->fLeaks[newbranch->fAllocCount],
00303                    &branch->fLeaks[i], sizeof(TMemInfo));
00304             newbranch->fAllocCount++;
00305             newbranch->fMemSize += branch->fLeaks[i].fSize;
00306          }
00307       free(branch->fLeaks);
00308       free(branch);
00309    }                 //loop over all old branches and rehash information
00310    free(fgLeak);
00311    fgLeak = newLeak;
00312    fgSize = newSize;
00313 }
00314 
00315 //______________________________________________________________________________
00316 void *TMemHashTable::AddPointer(size_t size, void *ptr)
00317 {
00318    // Add pointer to table.
00319 
00320    void *p = 0;
00321 
00322    if (ptr == 0) {
00323       p = malloc(size);
00324       if (!p) {
00325          Error("TMemHashTable::AddPointer", "malloc failure");
00326          _exit(1);
00327       }
00328    } else {
00329       p = realloc((char *) ptr, size);
00330       if (!p) {
00331          Error("TMemHashTable::AddPointer", "realloc failure");
00332          _exit(1);
00333       }
00334       return p;
00335    }
00336 
00337    if (!fgSize)
00338       Init();
00339    fgAllocCount++;
00340    if ((fgAllocCount / fgSize) > 128)
00341       RehashLeak(fgSize * 2);
00342    int hash = int(TString::Hash(&p, sizeof(void*)) % fgSize);
00343    TMemTable *branch = fgLeak[hash];
00344    branch->fAllocCount++;
00345    branch->fMemSize += size;
00346    for (;;) {
00347       for (int i = branch->fFirstFreeSpot; i < branch->fTableSize; i++)
00348          if (branch->fLeaks[i].fAddress == 0) {
00349             branch->fLeaks[i].fAddress = p;
00350             branch->fLeaks[i].fSize = size;
00351             void *sp = 0;
00352             int j = 0;
00353             void *stptr[stack_history_size + 1];
00354             for (j = 0; (j < stack_history_size); j++) {
00355                sp = get_stack_pointer(j + 1);
00356                if (sp == 0)
00357                   break;
00358                stptr[j] = sp;
00359             }
00360             TStackInfo *info = fgStackTable.FindInfo(j, stptr);
00361             info->Inc(size);
00362             branch->fLeaks[i].fStackIndex = fgStackTable.GetIndex(info);
00363             branch->fFirstFreeSpot = i + 1;
00364             return p;
00365          }
00366 
00367       int newTableSize =
00368           branch->fTableSize == 0 ? 16 : branch->fTableSize * 2;
00369       branch->fLeaks =
00370           (TMemInfo *) realloc(branch->fLeaks,
00371                                sizeof(TMemInfo) * newTableSize);
00372       if (!branch->fLeaks) {
00373          Error("TMemHashTable::AddPointer", "realloc failure (2)");
00374          _exit(1);
00375       }
00376       memset(branch->fLeaks + branch->fTableSize, 0, sizeof(TMemInfo) *
00377              (newTableSize - branch->fTableSize));
00378       branch->fTableSize = newTableSize;
00379    }
00380 }
00381 
00382 //______________________________________________________________________________
00383 void TMemHashTable::FreePointer(void *p)
00384 {
00385    // Free pointer.
00386 
00387    if (p == 0)
00388       return;
00389    int hash = int(TString::Hash(&p, sizeof(void*)) % fgSize);
00390    fgAllocCount--;
00391    TMemTable *branch = fgLeak[hash];
00392    for (int i = 0; i < branch->fTableSize; i++) {
00393       if (branch->fLeaks[i].fAddress == p) {
00394          branch->fLeaks[i].fAddress = 0;
00395          branch->fMemSize -= branch->fLeaks[i].fSize;
00396          if (i < branch->fFirstFreeSpot)
00397             branch->fFirstFreeSpot = i;
00398          free(p);
00399          TStackInfo *info =
00400              fgStackTable.GetInfo(branch->fLeaks[i].fStackIndex);
00401          info->Dec(branch->fLeaks[i].fSize);
00402          branch->fAllocCount--;
00403          return;
00404       }
00405    }
00406    //
00407    //if try to delete non existing pointer
00408    //printf("***TMemHashTable::FreePointer: Multiple deletion %8p ** ?  \n",p);
00409    //  printf("-+-+%8p  \n",p);
00410    //free(p);
00411    if (fgMultDeleteTable.fTableSize + 1 > fgMultDeleteTable.fAllocCount) {
00412       int newTableSize =
00413           fgMultDeleteTable.fTableSize ==
00414           0 ? 16 : fgMultDeleteTable.fTableSize * 2;
00415       fgMultDeleteTable.fLeaks =
00416           (TMemInfo *) realloc(fgMultDeleteTable.fLeaks,
00417                                sizeof(TMemInfo) * newTableSize);
00418       fgMultDeleteTable.fAllocCount = newTableSize;
00419    }
00420 
00421    fgMultDeleteTable.fLeaks[fgMultDeleteTable.fTableSize].fAddress = 0;
00422    void *sp = 0;
00423    void *stptr[stack_history_size + 1];
00424    int j;
00425    for (j = 0; (j < stack_history_size); j++) {
00426       sp = get_stack_pointer(j + 1);
00427       if (sp == 0)
00428          break;
00429       stptr[j] = sp;
00430    }
00431    TStackInfo *info = fgStackTable.FindInfo(j, stptr);
00432    info->Dec(0);
00433    fgMultDeleteTable.fLeaks[fgMultDeleteTable.fTableSize].fStackIndex =
00434        fgStackTable.GetIndex(info);
00435    fgMultDeleteTable.fTableSize++;
00436 }
00437 
00438 //______________________________________________________________________________
00439 void TMemHashTable::Dump()
00440 {
00441    // Print memory check information.
00442 
00443    const char *filename;
00444    if (gEnv)
00445       filename = gEnv->GetValue("Root.MemCheckFile", "memcheck.out");
00446    else
00447       filename = "memcheck.out";
00448 
00449    char *fn = 0;
00450    if (gSystem)
00451       fn = gSystem->ExpandPathName(filename);
00452 
00453    FILE *fp;
00454    if (!(fp = fn ? fopen(fn, "w") : fopen(filename, "w")))
00455       Error("TMenHashTable::Dump", "could not open %s", filename);
00456    else {
00457       /*
00458          for (int i = 0; i <  fgMultDeleteTable.fTableSize; i++){
00459          fprintf(fp, "size %9ld  ",(long)0);
00460          fprintf(fp, "stack:");
00461          TStackInfo *info = fgStackTable.GetInfo(fgMultDeleteTable.fLeaks[i].fStackIndex);
00462          for (int j=0; info->StackAt(j);j++)
00463          fprintf(fp, "%8p  ", info->StackAt(j));
00464          fprintf(fp, "\n");
00465          }
00466        */
00467       TStackInfo *info = fgStackTable.First();
00468       while (info->fSize) {
00469          fprintf(fp, "size %d:%d:%d:%d  ",
00470                  info->fTotalAllocCount, info->fTotalAllocSize,
00471                  info->fAllocCount, info->fAllocSize);
00472          fprintf(fp, "stack:");
00473          for (int j = 0; info->StackAt(j); j++)
00474             fprintf(fp, "%8p  ", info->StackAt(j));
00475          fprintf(fp, "\n");
00476          info = info->Next();
00477       }
00478       fclose(fp);
00479    }
00480    delete [] fn;
00481 }
00482 
00483 //______________________________________________________________________________
00484 static void *get_stack_pointer(int level)
00485 {
00486    // These special __builtin calls are supported by gcc "only".
00487    // For other compiler one will need to implement this again !
00488 
00489    void *p = 0;
00490 #if defined(R__GNU) && (defined(R__LINUX) || defined(R__HURD)) && \
00491    !defined(__alpha__)
00492    switch (level) {
00493    case 0:
00494       if (__builtin_frame_address(1))
00495          p = __builtin_return_address(1);
00496       break;
00497    case 1:
00498       if (__builtin_frame_address(2))
00499          p = __builtin_return_address(2);
00500       break;
00501    case 2:
00502       if (__builtin_frame_address(3))
00503          p = __builtin_return_address(3);
00504       break;
00505    case 3:
00506       if (__builtin_frame_address(4))
00507          p = __builtin_return_address(4);
00508       break;
00509    case 4:
00510       if (__builtin_frame_address(5))
00511          p = __builtin_return_address(5);
00512       break;
00513    case 5:
00514       if (__builtin_frame_address(6))
00515          p = __builtin_return_address(6);
00516       break;
00517    case 6:
00518       if (__builtin_frame_address(7))
00519          p = __builtin_return_address(7);
00520       break;
00521    case 7:
00522       if (__builtin_frame_address(8))
00523          p = __builtin_return_address(8);
00524       break;
00525    case 8:
00526       if (__builtin_frame_address(9))
00527          p = __builtin_return_address(9);
00528       break;
00529    case 9:
00530       if (__builtin_frame_address(10))
00531          p = __builtin_return_address(10);
00532       break;
00533    case 10:
00534       if (__builtin_frame_address(11))
00535          p = __builtin_return_address(11);
00536       break;
00537    case 11:
00538       if (__builtin_frame_address(12))
00539          p = __builtin_return_address(12);
00540       break;
00541    case 12:
00542       if (__builtin_frame_address(13))
00543          p = __builtin_return_address(13);
00544       break;
00545    case 13:
00546       if (__builtin_frame_address(14))
00547          p = __builtin_return_address(14);
00548       break;
00549    case 14:
00550       if (__builtin_frame_address(15))
00551          p = __builtin_return_address(15);
00552       break;
00553    case 15:
00554       if (__builtin_frame_address(16))
00555          p = __builtin_return_address(16);
00556       break;
00557    case 16:
00558       if (__builtin_frame_address(17))
00559          p = __builtin_return_address(17);
00560       break;
00561    case 17:
00562       if (__builtin_frame_address(18))
00563          p = __builtin_return_address(18);
00564       break;
00565    case 18:
00566       if (__builtin_frame_address(19))
00567          p = __builtin_return_address(19);
00568       break;
00569    case 19:
00570       if (__builtin_frame_address(20))
00571          p = __builtin_return_address(20);
00572       break;
00573 
00574    default:
00575       p = 0;
00576    }
00577 #else
00578    if (level) { }
00579 #endif
00580    return p;
00581 }

Generated on Tue Jul 5 14:11:56 2011 for ROOT_528-00b_version by  doxygen 1.5.1