TMemStatShow.cxx

Go to the documentation of this file.
00001 // @(#)root/treeviewer:$Id: TMemStatShow.cxx 37335 2010-12-06 14:51:05Z brun $
00002 // Author: Rene Brun   21/09/2010
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2010, 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 // Utility class post-processing the file generated by TMemStat (default memstat.root)
00014 //
00015 // TMemStat records all the calls to malloc and free and write a TTree
00016 // with the position where the memory is allocated/freed , as well as
00017 // the number of bytes.
00018 //
00019 // To use the class TMemStat, add the following statement at the beginning
00020 // of your script or program
00021 //     TMemStat mm("gnubuiltin");
00022 // or in an interactive session do something like:
00023 //    root > TMemStat mm("gnubuiltin");
00024 //    root > .x somescript.C
00025 //    root > .q
00026 //
00027 // another (may be more practical way) is to modify $ROOTSYS/etc/system.rootrc
00028 // and activate the variable
00029 //    Root.TMemStat:           1
00030 //
00031 // The file collected by TMemStat is named memstat_ProcessID and can be analyzed and results shown
00032 // by executing the static function Show.
00033 // When TMemStat is active it recors every call to malloc/free in a ROOT Tree.
00034 // You must be careful when running jobs with many millions (or more) of calls
00035 // to malloc/free because the generated Tree may become very large.
00036 // The TMemStat constructor TMemStat(const char* system, Int_t buffersize, Int_t maxcalls)
00037 // has its 3 arguments optional:
00038 //   -system refers to the internal algorithm to compute the back traces.
00039 //    the recommended value is "gnubuiltin"
00040 //   -buffersize is the number of calls to malloc or free that can be stored in one memory buffer.
00041 //    when the buffer is full, the calls to malloc/free pointing to the same location
00042 //    are eliminated and not written to the final Tree. The default value 100000
00043 //    is such that between 50 and 90% of the calls are eliminated depending on the application.
00044 //    You can set buffersize <=1 to keep every single call to malloc/free.
00045 //   -maxcalls can set a limit for the maximum number of calls to be registered in the Tree.
00046 //    The default value is 5000000.
00047 // The 3 arguments can be set  in $ROOTSYS/etc/system.rootrc
00048 //    Root.TMemStat.system      gnubuiltin
00049 //    Root.TMemStat.buffersize  100000
00050 //    Root.TMemStat.maxcalls    5000000
00051 //
00052 // TMemStat::Show creates 3 canvases.
00053 // -In canvas1 it displays a dynamic histogram showing for pages (10 kbytes by default)
00054 //  the percentage of the page used.
00055 //  A summary pave shows the total memory still in use when the TMemStat object
00056 //  goes out of scope and the average occupancy of the pages.
00057 //  The average occupancy gives a good indication of the memory fragmentation.
00058 //  When moving the mouse on this canvas, a tooltip shows the backtrace for the allocations
00059 //  at the address at the mouse position.
00060 //
00061 // -In canvas2 it displays the histogram of memory leaks in decreasing order.
00062 //  when moving the mouse on this canvas, a tooltip shows the backtrace for the leak
00063 //  in the bin below the mouse.
00064 //
00065 // -In canvas3 it displays the histogram of the nbigleaks largest leaks (default is 20)
00066 //    for each leak, the number of allocs and average alloc size is shown.
00067 //
00068 //
00069 // Simply do:
00070 //   root > TMemStat::Show()
00071 // or specifying arguments
00072 //   root > TMemStat::Show(0.1,20,"mydir/mymemstat.root");
00073 //
00074 // The first argument to Show is the percentage of the time of the original job
00075 // that produced the file after which the display is updated. By default update=0.1,
00076 // ie 10 time intervals will be shown.
00077 // The second argument is nbigleaks. if <=0 canvas2 and canvas3 are not shown
00078 // The third argument is the imput file name (result of TMemStat).
00079 // If this argument is omitted, Show will take the most recent file
00080 // generated by TMemStat.
00081 //
00082 // You can restrict the address range to be analyzed via TMemStatShow::SetAddressRange
00083 // You can restrict the entry range to be analyzed via TMemStatShow::SetEntryRange
00084 //
00085 //Author: Rene Brun 7 July 2010
00086 //___________________________________________________________________________
00087       
00088 #include "TMemStatShow.h"
00089 #include "TMath.h"
00090 #include "TFile.h"
00091 #include "TTree.h"
00092 #include "TCanvas.h"
00093 #include "TStyle.h"
00094 #include "TH1.h"
00095 #include "TPaveText.h"
00096 #include "TPaveLabel.h"
00097 #include "TSystem.h"
00098 #include "TGClient.h"
00099 #include "TGToolTip.h"
00100 #include "TRootCanvas.h"
00101 
00102    static MemInfo_t minfo;
00103 
00104    TTree     *TMemStatShow::fgT = 0;         //TMemStat Tree
00105    TH1D      *TMemStatShow::fgHalloc = 0;    //histogram with allocations
00106    TH1D      *TMemStatShow::fgHfree = 0;     //histogram with frees
00107    TH1D      *TMemStatShow::fgH = 0;         //histogram with allocations - frees
00108    TH1I      *TMemStatShow::fgHleaks = 0;    //histogram with leaks
00109    TH1I      *TMemStatShow::fgHentry = 0;    //histogram with entry numbers in the TObjArray
00110    TH1I      *TMemStatShow::fgHdiff = 0;     //histogram with diff of entry number between alloc/free
00111    
00112    TGToolTip *TMemStatShow::fgTip1 = 0;      //pointer to tool tip for canvas 1
00113    TGToolTip *TMemStatShow::fgTip2 = 0;      //pointer to tool tip for canvas 2
00114    TObjArray *TMemStatShow::fgBtidlist = 0;  //list of back trace ids
00115    Double_t  *TMemStatShow::fgV1 = 0;        //pointer to V1 array of TTree::Draw (pos)
00116    Double_t  *TMemStatShow::fgV2 = 0;        //pointer to V2 array of TTree::Draw (nbytes)
00117    Double_t  *TMemStatShow::fgV3 = 0;        //pointer to V3 array of TTree::Draw (time)
00118    Double_t  *TMemStatShow::fgV4 = 0;        //pointer to V4 array of TTree::Draw (btid)
00119    TCanvas   *TMemStatShow::fgC1 = 0;        //pointer to canvas showing allocs/deallocs vs time
00120    TCanvas   *TMemStatShow::fgC2 = 0;        //pointer to canvas with leaks in decreasing order
00121    TCanvas   *TMemStatShow::fgC3 = 0;        //pointer to canvas showing the main leaks
00122 
00123    Long64_t TMemStatShow::fgEntryFirst   = 0; //first address to process
00124    Long64_t TMemStatShow::fgEntryN       = 0; //number of addresses in bytes to process
00125    Long64_t TMemStatShow::fgAddressFirst = 0; //first entry to process
00126    Long64_t TMemStatShow::fgAddressN     = 0; //number of entries to process
00127    
00128 //___________________________________________________________________________
00129 void TMemStatShow::SetAddressRange(Long64_t nbytes, Long64_t first)
00130 {
00131    //specify a memory address range to process (static function).
00132    // This function can be used to restrict the range of memory addresses 
00133    // to be analyzed. For example whem TmemStat is run on a 64 bits machine and
00134    // the results visualized on a 32 bits machine, it might be necessary to
00135    // restrict the analysis range to the addresses below 2 Gigabytes, eg
00136    //   TMemStatShow::SetMemoryRange(500000000,0); //analyse only the first 500 MBytes
00137    // -first : first address to process (default is 0)
00138    // -nbytes : number of addresses in bytes to process starting at first
00139    //             if 0 (default), then all addresses are processed
00140    
00141    fgAddressFirst = first;
00142    fgAddressN     = nbytes;
00143 }
00144    
00145 //___________________________________________________________________________
00146 void TMemStatShow::SetEntryRange(Long64_t nentries, Long64_t first)
00147 {
00148    //specify a range of entries to process (static function)
00149    // -first : first entry to process (default is 0)
00150    // -nentries : number of entries to process starting at first
00151    //             if 0 (default), then all entries are processed
00152    // call this function when the amount of data collected in the Tree is large
00153    // and therefore making the analysis slow.
00154    
00155    fgEntryFirst = first;
00156    fgEntryN     = nentries;
00157 }
00158 
00159 //___________________________________________________________________________
00160 void TMemStatShow::Show(double update, int nbigleaks, const char* fname) 
00161 {
00162    // function called by TMemStat::Show
00163    // Open the memstat data file, then call TTree::Draw to precompute
00164    // the arrays of positions and nbytes per entry.
00165    // update is the time interval in the data file  in seconds after which
00166    // the display is updated. For example is the job producing the memstat.root file
00167    // took 100s to execute, an update of 0.1s will generate 1000 time views of
00168    // the memory use.
00169    // the histogram hbigleaks will contain the nbigleaks largest leaks
00170    // if fname=="*" (default), the most recent file memstat*.root will be taken.
00171    
00172    
00173    TString s;
00174    if (!fname || strlen(fname) <5 || strstr(fname,"*")) {
00175       //take the most recent file memstat*.root
00176       s = gSystem->GetFromPipe("ls -lrt memstat*.root");
00177       Int_t ns = s.Length();
00178       fname = strstr(s.Data()+ns-25,"memstat");
00179    }
00180    printf("Analyzing file: %s\n",fname);
00181    TFile *f = TFile::Open(fname);
00182    if (!f) {
00183       printf("Cannot open file %s\n",fname);
00184       return;
00185    }
00186    fgT = (TTree*)f->Get("T");
00187    if (!fgT) {
00188       printf("cannot find the TMemStat TTree named T in file %s\n",fname);
00189       return;
00190    }
00191    if (update <= 0) {
00192       printf("Illegal update value %g, changed to 0.01\n",update);
00193       update = 0.01;
00194    }
00195    if (update < 0.001) printf("Warning update parameter is very small, processing may be slow\n");
00196 
00197    //autorestrict the amount of data to analyze 
00198    gSystem->GetMemInfo(&minfo);
00199    Int_t nfree = minfo.fMemTotal - minfo.fMemUsed;  //in Mbytes
00200    printf("TMemStat::Show info: you are running on a machine with %d free MBytes of memory\n",nfree);
00201    Long64_t nfreebytes = 200000*Long64_t(nfree); //use only 20% of the memory available
00202    if (fgAddressN <=0) fgAddressN = nfreebytes;
00203    Long64_t nentries = fgT->GetEntries();
00204    if (fgEntryN > 0 && nentries > fgEntryN) nentries = fgEntryN;
00205    if (2*8*nentries > 4*nfreebytes) {
00206       nentries = 4*nfreebytes/16;
00207       printf("not enough memory, restricting analysis to %lld entries\n",nentries);
00208    }
00209    fgT->SetEstimate(nentries);
00210    Long64_t nsel = fgT->Draw("pos","pos>0","goff",nentries);
00211    fgV1 = fgT->GetV1();
00212    Long64_t ivmin = (Long64_t)TMath::MinElement(nsel,fgV1);
00213    Long64_t ivmax = (Long64_t)TMath::MaxElement(nsel,fgV1);
00214    if (ivmax-ivmin > fgAddressN) ivmax = ivmin+fgAddressN;
00215    printf("TMemStatShow::Show will analyze only %lld bytes in its first pass\n",ivmax);
00216    
00217    
00218    //initialize statics
00219    fgTip1 = 0;
00220    fgTip2 = 0;
00221    fgBtidlist = 0;
00222       
00223    Long64_t ne = nfreebytes/32;
00224    if (ne < nentries) nentries = ne;
00225    fgT->SetEstimate(nentries+10);
00226    printf("sel: ivmin=%lld, ivmax=%lld, nentries=%lld\n",ivmin,ivmax,nentries);
00227    nsel = fgT->Draw("pos:nbytes:time:btid",
00228       TString::Format("pos>%g && pos<%g",Double_t(ivmin),Double_t(ivmax)),
00229       "goff",nentries,fgEntryFirst);
00230 
00231    //now we compute the best binning for the histogram
00232    Int_t nbytes;
00233    Double_t pos;
00234    fgV1 = fgT->GetV1();
00235    fgV2 = fgT->GetV2();
00236    fgV3 = fgT->GetV3();
00237    fgV4 = fgT->GetV4();
00238    ivmin = (Long64_t)TMath::MinElement(nsel,fgV1);
00239    ivmax = (Long64_t)TMath::MaxElement(nsel,fgV1);
00240    Long64_t bw = 1000;
00241    Double_t dvv = (Double_t(ivmax) - Double_t(ivmin))/Double_t(bw);
00242    Long64_t nbins = Long64_t(dvv);
00243    ivmin = ivmin -ivmin%bw;
00244    ivmax = ivmin+bw*nbins;
00245    Long64_t nvm = Long64_t(ivmax-ivmin+1);
00246    printf("==>The data Tree contains %lld entries with addresses in range[%lld,%lld]\n",nsel,ivmin,ivmax);
00247    //ne = (1000000*nfree-nvm*12)/32;
00248    ne = 1000000*nfree/32;
00249    if (ne < 0) return;    
00250    if (ne < nentries) {
00251       //we take only the first side of the allocations
00252       //we are mostly interested by the small allocations, so we select
00253       //only values in the first gigabyte
00254       nsel = fgT->Draw("pos:nbytes:time:btid",
00255          TString::Format("pos>=%g && pos<%g",Double_t(ivmin),Double_t(ivmax)),"goff",ne,fgEntryFirst);
00256       fgV1 = fgT->GetV1();
00257       fgV2 = fgT->GetV2();
00258       fgV3 = fgT->GetV3();
00259       fgV4 = fgT->GetV4();
00260       ivmin = (Long64_t)TMath::MinElement(nsel,fgV1);
00261       ivmax = (Long64_t)TMath::MaxElement(nsel,fgV1);
00262       bw = 10000;
00263       dvv = (Double_t(ivmax) - Double_t(ivmin))/Double_t(bw);
00264       nbins = Long64_t(dvv+0.5);
00265       ivmin = ivmin -ivmin%bw;
00266       ivmax = ivmin+bw*nbins;
00267       printf("==>Address range or/and Entry range is too large\n");
00268       printf("==>restricting the analysis range to [%lld,%lld] and %lld entries\n",ivmin,ivmax,ne);
00269       printf("==>you can restrict the address range with TMemStatShow::SetAddressRange\n");
00270       printf("==>you can restrict the entries range with TMemStatShow::SetEntryRange\n");
00271    } 
00272    update *= 0.0001*fgV3[nsel-1]; //convert time per cent in seconds
00273    nvm = Long64_t(ivmax-ivmin);
00274    Long64_t *nbold = new Long64_t[nvm];
00275    Int_t *ientry  = new Int_t[nvm];
00276    if (!nbold || !ientry) {
00277       printf("you do not have enough memory to run, %lld bytes needed\n",12*nvm);
00278       return;
00279    }
00280    memset(nbold,0,nvm*8);
00281    memset(ientry,0,nvm*4);
00282    Double_t dv = (ivmax-ivmin)/nbins;
00283    TH1D *h = new TH1D("h",Form("%s;pos;per cent of pages used",fname),nbins,ivmin,ivmax);
00284    fgH = h;
00285    TAxis *axis = h->GetXaxis();
00286    gStyle->SetOptStat("ie");
00287    h->SetFillColor(kRed);
00288    h->SetMinimum(0);
00289    h->SetMaximum(100);
00290    fgHalloc = new TH1D("fgHalloc",Form("%s;pos;number of mallocs",fname),nbins,ivmin,ivmax);
00291    fgHfree  = new TH1D("fgHfree", Form("%s;pos;number of frees",fname),nbins,ivmin,ivmax);
00292    fgHdiff = new TH1I("fgHdiff","",1000,0,1e5);
00293    //open a canvas and draw the empty histogram
00294    fgC1 = new TCanvas("fgC1","c1",1200,600);
00295    fgC1->SetFrameFillColor(kYellow-3);
00296    fgC1->SetGridx();
00297    fgC1->SetGridy();
00298    h->Draw();
00299    //create a TPaveText to show the summary results
00300    TPaveText *pvt = new TPaveText(.5,.9,.75,.99,"brNDC");
00301    pvt->Draw();
00302    //create a TPaveLabel to show the time
00303    TPaveLabel *ptime = new TPaveLabel(.905,.7,.995,.76,"time","brNDC");
00304    ptime->SetFillColor(kYellow-3);
00305    ptime->Draw();
00306    //draw producer identifier
00307    TNamed *named = (TNamed*)fgT->GetUserInfo()->FindObject("SysInfo");
00308    TText tmachine;
00309    tmachine.SetTextSize(0.02);
00310    tmachine.SetNDC();
00311    if (named) tmachine.DrawText(0.01,0.01,named->GetTitle());
00312 
00313    //start loop on selected rows
00314    Int_t bin,nb=0,j;
00315    Long64_t ipos;
00316    Double_t dbin,rest,time;
00317    Double_t updateLast = 0;
00318    Int_t nleaks = 0;
00319    Int_t i;
00320    for (i=0;i<nsel;i++) {
00321       pos    = fgV1[i];
00322       ipos = (Long64_t)(pos-ivmin);
00323       nbytes = (Int_t)fgV2[i];
00324       time = 0.0001*fgV3[i];
00325       bin = axis->FindBin(pos);
00326       if (bin<1 || bin>nbins) continue;
00327       dbin = axis->GetBinUpEdge(bin)-pos;
00328       if (nbytes > 0) {
00329          ientry[ipos] = i;
00330          fgHalloc->Fill(pos);
00331          if (dbin > nbytes) dbin = nbytes;
00332          //fill bytes in the first page
00333          h->AddBinContent(bin,100*dbin/dv);
00334          //fill bytes in full following pages
00335          nb = Int_t((nbytes-dbin)/dv);
00336          if (bin+nb >nbins) nb = nbins-bin;
00337          for (j=1;j<=nb;j++) h->AddBinContent(bin+j,100);
00338          //fill the bytes remaining in last page
00339          rest = nbytes-nb*dv-dbin;
00340          if (rest > 0) h->AddBinContent(bin+nb+1,100*rest/dv);
00341          //we save nbytes at pos. This info will be used when we free this slot
00342          //if (nbold[ipos] > 0) printf("reallocating %d bytes (was %lld) at %lld, entry=%d\n",nbytes,nbold[ipos],ipos,i);
00343          if (nbold[ipos] == 0) {
00344             nleaks++;            
00345             //save the Tree entry number where we made this allocation
00346             ientry[ipos] = i;
00347          }
00348          nbold[ipos] = nbytes;
00349       } else {
00350          fgHfree->Fill(pos);
00351          nbytes = nbold[ipos];
00352          if (bin+nb >nbins) nb = nbins-bin;
00353          nbold[ipos] = 0; nleaks--;
00354          fgHdiff->Fill(i-ientry[ipos]);
00355          if (nbytes <= 0) continue;
00356          //fill bytes free in the first page
00357          if (dbin > nbytes) dbin = nbytes;
00358          h->AddBinContent(bin,-100*dbin/dv);
00359          //fill bytes free in full following pages
00360          nb = Int_t((nbytes-dbin)/dv);
00361          if (bin+nb >nbins) nb = nbins-bin;
00362          for (j=1;j<=nb;j++) h->AddBinContent(bin+j,-100);
00363          //fill the bytes free in  in last page
00364          rest = nbytes-nb*dv-dbin;
00365          if (rest > 0) h->AddBinContent(bin+nb+1,-100*rest/dv);
00366 
00367       }
00368       if (time -updateLast > update) {
00369          //update canvas at regular intervals
00370          updateLast = time;
00371          h->SetEntries(i);
00372          fgC1->Modified();
00373          pvt->GetListOfLines()->Delete();
00374          Double_t mbytes = 0;
00375          Int_t nonEmpty = 0;
00376          Double_t w;
00377          for (Int_t k=1;k<nbins;k++) {
00378             w = h->GetBinContent(k);
00379             if (w > 0) {
00380                nonEmpty++;
00381                mbytes += 0.01*w*dv;
00382             }
00383          }
00384          Double_t occupancy = mbytes/(nonEmpty*0.01*dv);
00385          pvt->AddText(Form("memory used = %g Mbytes",mbytes*1e-6));
00386          pvt->AddText(Form("page occupancy = %f per cent",occupancy));
00387          pvt->AddText("(for non empty pages only)");
00388          ptime->SetLabel(Form("%g sec",time));
00389          
00390          fgC1->Update();
00391          gSystem->ProcessEvents();
00392       }
00393    }
00394    h->SetEntries(nsel);
00395    if (nleaks < 0) nleaks=0;
00396    Int_t nlmax = nleaks;
00397    nleaks += 1000;
00398    Int_t *lindex  = new Int_t[nleaks];
00399    Int_t *entry   = new Int_t[nleaks];
00400    Int_t *ileaks  = new Int_t[nleaks];
00401 
00402    nleaks =0;
00403    for (Int_t ii=0;ii<nvm;ii++) {
00404       if (nbold[ii] > 0) {
00405          ileaks[nleaks] = (Int_t)nbold[ii];
00406          entry[nleaks]  = ientry[ii];
00407          nleaks++;
00408          if (nleaks > nlmax) break;
00409       }
00410    }
00411    TMath::Sort(nleaks,ileaks,lindex);
00412    fgHentry = new TH1I("fgHentry","leak entry index",nleaks,0,nleaks);
00413    fgHleaks = new TH1I("fgHleaks","leaks;leak number;nbytes in leak",nleaks,0,nleaks);
00414    for (Int_t k=0;k<nleaks;k++) {
00415       Int_t kk = lindex[k];
00416       i = entry[kk];
00417       fgHentry->SetBinContent(k+1,i);
00418       fgHleaks->SetBinContent(k+1,ileaks[kk]);
00419    }
00420    delete [] ileaks;
00421    delete [] entry;
00422    delete [] lindex;
00423    delete [] nbold;
00424    delete [] ientry;
00425    fgHentry->SetEntries(nleaks);
00426    fgHleaks->SetEntries(nleaks);
00427    
00428    
00429    //construct the first tooltip
00430    fgC1->Modified();
00431    fgC1->Update();
00432    TRootCanvas *rc1 = (TRootCanvas *)fgC1->GetCanvasImp();
00433    TGMainFrame *frm1 = dynamic_cast<TGMainFrame *>(rc1);
00434    // create the tooltip with a timeout of 250 ms
00435    if (!fgTip1) fgTip1 = new TGToolTip(gClient->GetDefaultRoot(), frm1, "", 250);
00436    fgC1->Connect("ProcessedEvent(Int_t, Int_t, Int_t, TObject*)",
00437                 "TMemStatShow", 0, "EventInfo1(Int_t, Int_t, Int_t, TObject*)");
00438    if (nbigleaks <= 0) return;
00439 
00440    //---------------------------------------------------------------------------
00441    //open a second canvas and draw the histogram with leaks in decreasing order
00442    fgC2 = new TCanvas("fgC2","c2",1200,600);
00443    fgC2->SetFrameFillColor(kCyan-6);
00444    fgC2->SetGridx();
00445    fgC2->SetGridy();
00446    fgC2->SetLogy();
00447    fgHleaks->SetFillColor(kRed-3);
00448    if (nleaks > 1000) fgHleaks->GetXaxis()->SetRange(1,1000);
00449    fgHleaks->Draw();
00450    //draw producer identifier
00451    if (named) tmachine.DrawText(0.01,0.01,named->GetTitle());
00452    
00453    //construct the second tooltip
00454    TRootCanvas *rc2 = (TRootCanvas *)fgC2->GetCanvasImp();
00455    TGMainFrame *frm2 = dynamic_cast<TGMainFrame *>(rc2);
00456    // create the tooltip with a timeout of 250 ms
00457    if (!fgTip2) fgTip2 = new TGToolTip(gClient->GetDefaultRoot(), frm2, "", 250);
00458    fgC2->Connect("ProcessedEvent(Int_t, Int_t, Int_t, TObject*)",
00459                "TMemStatShow", 0, "EventInfo2(Int_t, Int_t, Int_t, TObject*)");
00460    
00461    //---------------------------------------------------------------------------
00462    //open a third canvas and draw the histogram with the nbigleaks largest leaks
00463    fgC3 = new TCanvas("fgC3","c3",1200,600);
00464    fgC3->SetFrameFillColor(kCyan-6);
00465    fgC3->SetGridx();
00466    fgC3->SetGridy();
00467    fgC3->SetLogx();
00468    fgC3->SetLeftMargin(0.05);
00469    fgC3->SetRightMargin(0.7);
00470    
00471    //fill histogram htotleaks accumulating in the same bin all leaks
00472    //from btids having identical nchar first characters   
00473    TH1I *htotleaks = new TH1I("htotleaks","main leaks sorted by btids",100,0,0);
00474    Int_t l;
00475    for (l=1;l<=nleaks;l++) {
00476       TString btstring = "";
00477       TMemStatShow::FillBTString(l,1,btstring);
00478       htotleaks->Fill(btstring.Data()+2,fgHleaks->GetBinContent(l));
00479    }
00480    Double_t tsize = 0.03;
00481    if (nbigleaks > 30) tsize = 0.02;
00482    htotleaks->LabelsOption(">");
00483    htotleaks->GetXaxis()->SetRange(1,nbigleaks); 
00484    htotleaks->GetXaxis()->SetLabelSize(tsize);
00485    htotleaks->GetYaxis()->SetLabelSize(tsize);
00486    htotleaks->SetFillColor(kBlue-3);
00487    htotleaks->Draw("hbar2 y+");
00488    
00489    //now loop on all the sorted bins and count the number of leaks
00490    Double_t xr = 0.96*fgC3->GetLeftMargin();
00491    Double_t xr2 = 1.04*fgC3->GetLeftMargin();
00492    Double_t ytop = 1-fgC3->GetTopMargin();
00493    Double_t ylow = fgC3->GetBottomMargin();
00494    Double_t dy = (ytop-ylow)/nbigleaks;
00495    TString btstring;
00496    TText tnl;
00497    tnl.SetNDC();
00498    tnl.SetTextSize(tsize);
00499    tnl.SetTextAlign(32);
00500    TText tnl2;
00501    tnl2.SetNDC();
00502    tnl2.SetTextSize(tsize);
00503    tnl2.SetTextAlign(12);
00504    tnl2.SetTextColor(kYellow);
00505    for (Int_t lb=1;lb<=nbigleaks;lb++) {
00506       if (htotleaks->GetBinContent(lb) <= 0) continue;
00507       const char *label = htotleaks->GetXaxis()->GetBinLabel(lb);
00508       Int_t nchlabel = strlen(label);
00509       if (nchlabel == 0) htotleaks->GetXaxis()->SetBinLabel(lb,"???");
00510       Int_t nl =0;
00511       for (l=1;l<=nleaks;l++) {
00512          btstring = "";
00513          TMemStatShow::FillBTString(l,1,btstring);
00514          if (nchlabel > 0) {     
00515             if (!strncmp(btstring.Data()+2,label,nchlabel)) nl++;
00516          } else {
00517             if (btstring.Length() == 0) nl++;
00518          }
00519       }
00520       Double_t yr = ylow +(lb-0.5)*dy;
00521       tnl.DrawText(xr,yr,Form("%d",nl));
00522       Int_t nbmean = Int_t(htotleaks->GetBinContent(lb)/nl);
00523       if (lb == 1) tnl2.DrawText(xr2,yr,Form("%d bytes/alloc",nbmean));
00524       else         tnl2.DrawText(xr2,yr,Form("%d",nbmean));
00525    }
00526    tnl.DrawText(xr,ytop+0.015,"nallocs");
00527    tnl.DrawText(1-fgC3->GetRightMargin(),0.5*ylow,"nbytes");
00528    //draw producer identifier
00529    if (named) tmachine.DrawText(0.01,0.01,named->GetTitle());
00530    
00531 }
00532 
00533 //______________________________________________________________________
00534 void TMemStatShow::EventInfo1(Int_t event, Int_t px, Int_t , TObject *selected)
00535 {
00536    // static: draw the tooltip showing the backtrace for the allocatios histogram
00537    if (!fgTip1) return;
00538    fgTip1->Hide();
00539    if (event == kMouseLeave)
00540       return;
00541    Double_t xpx  = fgC1->AbsPixeltoX(px);
00542    Double_t xpx1 = fgC1->AbsPixeltoX(px+1);
00543    Int_t bin  = fgH->GetXaxis()->FindBin(xpx);
00544    Int_t bin1 = fgH->GetXaxis()->FindBin(xpx1);
00545    //to take into account consecutive bins on the same pixel
00546    while (bin <= bin1) {
00547       if (fgH->GetBinContent(bin) > 0) break;
00548       bin++;
00549    }
00550    if (fgH->GetBinContent(bin) <= 0) return;
00551    if (bin <=0 || bin > fgH->GetXaxis()->GetNbins()) return;
00552    Double_t posmin = fgH->GetXaxis()->GetBinLowEdge(bin);
00553    Double_t posmax = fgH->GetXaxis()->GetBinUpEdge(bin);
00554    Int_t nsel   = (Int_t)fgT->GetSelectedRows();
00555    Int_t entry  = 0;
00556    Int_t nhits  = 0;
00557    Int_t nbytes = 0;
00558    //search for all allocations in this bin and select last one only
00559    for (Int_t i=0;i<nsel;i++) {
00560       if (fgV2[i] < 0) continue;
00561       if (fgV1[i] < posmax && fgV1[i]+fgV2[i] >posmin) {
00562          entry = i;
00563          nbytes = (Int_t)fgV2[i];
00564          nhits++;
00565       }
00566    }
00567    if (!nhits) return;
00568       
00569    Double_t time = 0.0001*fgV3[entry];
00570    TString ttip;
00571    TMemStatShow::FillBTString(entry,0,ttip);
00572    
00573    if (selected) {
00574       TString form1 = TString::Format("  Alloc(%d) at %lld of %d bytes, time=%gseconds\n\n",nhits,Long64_t(fgV1[entry]),nbytes,time);
00575       fgTip1->SetText(TString::Format("%s%s",form1.Data(),ttip.Data() ));
00576       fgTip1->SetPosition(px+15, 100);
00577       fgTip1->Reset();
00578    }
00579 }
00580 
00581 //______________________________________________________________________
00582 void TMemStatShow::EventInfo2(Int_t event, Int_t px, Int_t , TObject *selected)
00583 {
00584    // static: draw the tooltip showing the backtrace for the histogram of leaks
00585    if (!fgTip2) return;
00586    fgTip2->Hide();
00587    if (event == kMouseLeave)
00588       return;
00589    Double_t xpx  = fgC2->AbsPixeltoX(px);
00590    Int_t bin = fgHleaks->GetXaxis()->FindBin(xpx);
00591    if (bin <=0 || bin > fgHleaks->GetXaxis()->GetNbins()) return;
00592    Int_t nbytes  = (Int_t)fgHleaks->GetBinContent(bin);
00593    Int_t entry   = (Int_t)fgHentry->GetBinContent(bin);
00594    Double_t time = 0.0001*fgV3[entry];
00595    TString ttip;
00596    TMemStatShow::FillBTString(entry,0,ttip);
00597       
00598    if (selected) {
00599       TString form1 = TString::Format("  Leak number=%d, leaking %d bytes at entry=%d    time=%gseconds\n\n",bin,nbytes,entry,time);
00600       fgTip2->SetText(TString::Format("%s%s",form1.Data(),ttip.Data() ));
00601       fgTip2->SetPosition(px+15, 100);
00602       fgTip2->Reset();
00603    }
00604 }
00605 
00606 //______________________________________________________________________
00607 void TMemStatShow::FillBTString(Int_t entry,Int_t mode,TString &btstring)
00608 {
00609    // static: fill btstring with the traceback corresponding to entry in T
00610    //          btstring must be initialized in calling function
00611    
00612    Int_t btid   = (Int_t)fgV4[entry];
00613    TH1I *hbtids = (TH1I*)fgT->GetUserInfo()->FindObject("btids");
00614    if (!hbtids) return;
00615    if (!fgBtidlist) fgBtidlist = (TObjArray*)fgT->GetUserInfo()->FindObject("FAddrsList");
00616    if (!fgBtidlist) fgBtidlist = (TObjArray*)gFile->Get("FAddrsList"); //old memstat files
00617    if (!fgBtidlist) return;
00618    Int_t nbt = (Int_t)hbtids->GetBinContent(btid-1);
00619    for (Int_t i=0;i<nbt;i++) {
00620       Int_t j = (Int_t)hbtids->GetBinContent(btid+i);
00621       TNamed *nm = (TNamed*)fgBtidlist->At(j);
00622       if (nm==0) break;
00623       char *title = (char*)nm->GetTitle();
00624       Int_t nch = strlen(title);
00625       if (nch < 10) continue;
00626       if (strstr(title,"malloc")) continue;
00627       if (strstr(title,"memstat")) continue;
00628       if (strstr(title,"TMemStatHook")) continue;
00629       char *bar = strchr(title+5,'|');
00630       if (!bar) bar = title;
00631       
00632       if (strstr(bar,"operator new")) continue;
00633       if (strstr(bar,"libMemStat")) continue;
00634       if (strstr(bar,"G__Exception")) continue;
00635       if (mode) {
00636          btstring += TString::Format("%s ",bar);
00637          if (btstring.Length() > 80) return;
00638       } else {
00639          btstring += TString::Format("%2d %s\n",i,bar+1);
00640       }
00641    }
00642 }

Generated on Tue Jul 5 15:43:44 2011 for ROOT_528-00b_version by  doxygen 1.5.1