CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FWCore/Services/src/Memory.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     Services
00004 // Class  :     Memory
00005 // 
00006 // Implementation:
00007 //
00008 // Original Author:  Jim Kowalkowski
00009 //
00010 // Change Log
00011 //
00012 // 1 - Apr 25, 2008 M. Fischler
00013 //      Collect event summary information and output to XML file and logger
00014 //      at the end of the job.  Involves split-up of updateAndPrint method.
00015 //
00016 // 2 - May 7, 2008 M. Fischler
00017 //      Collect module summary information and output to XML file and logger
00018 //      at the end of the job.
00019 //
00020 // 3 - Jan 14, 2009 Natalia Garcia Nebot
00021 //      Added:  - Average rate of growth in RSS and peak value attained.
00022 //              - Average rate of growth in VSize over time, Peak VSize
00023 //
00024 //
00025 
00026 #include "FWCore/Services/src/Memory.h"
00027 #include "DataFormats/Provenance/interface/ModuleDescription.h"
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 #include "FWCore/MessageLogger/interface/JobReport.h"
00031 #include "FWCore/ServiceRegistry/interface/Service.h"
00032 #include "FWCore/Utilities/interface/Exception.h"
00033 #include "FWCore/Utilities/interface/MallocOpts.h"
00034 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00035 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00036 
00037 #ifdef __linux__
00038 #include <malloc.h>
00039 #endif
00040 #include <sstream>
00041 #include <iostream>
00042 #include <string>
00043 #include <stdio.h>
00044 #include <string.h>
00045 #include <cstring>
00046 
00047 #ifdef __linux__
00048 #define LINUX 1
00049 #endif
00050 
00051 #include <unistd.h>
00052 #include <fcntl.h>
00053 
00054 namespace edm {
00055   namespace service {
00056 
00057     static std::string d2str(double d){
00058       std::ostringstream t;
00059       t << d;
00060       return t.str();
00061     }
00062     
00063     static std::string i2str(int i){
00064       std::ostringstream t;
00065       t << i;
00066       return t.str();
00067     }
00068 
00069     struct linux_proc {
00070       int pid; // %d
00071       char comm[400]; // %s
00072       char state; // %c
00073       int ppid; // %d
00074       int pgrp; // %d
00075       int session; // %d
00076       int tty; // %d
00077       int tpgid; // %d
00078       unsigned int flags; // %u
00079       unsigned int minflt; // %u
00080       unsigned int cminflt; // %u
00081       unsigned int majflt; // %u
00082       unsigned int cmajflt; // %u
00083       int utime; // %d
00084       int stime; // %d
00085       int cutime; // %d
00086       int cstime; // %d
00087       int counter; // %d
00088       int priority; // %d
00089       unsigned int timeout; // %u
00090       unsigned int itrealvalue; // %u
00091       int starttime; // %d
00092       unsigned int vsize; // %u
00093       unsigned int rss; // %u
00094       unsigned int rlim; // %u
00095       unsigned int startcode; // %u
00096       unsigned int endcode; // %u
00097       unsigned int startstack; // %u
00098       unsigned int kstkesp; // %u
00099       unsigned int kstkeip; // %u
00100       int signal; // %d
00101       int blocked; // %d
00102       int sigignore; // %d
00103       int sigcatch; // %d
00104       unsigned int wchan; // %u
00105     };
00106       
00107     procInfo SimpleMemoryCheck::fetch()
00108     {
00109       procInfo ret;
00110       double pr_size=0.0, pr_rssize=0.0;
00111       
00112 #ifdef LINUX
00113       linux_proc pinfo;
00114       int cnt;
00115 
00116       lseek(fd_,0,SEEK_SET);
00117     
00118       if((cnt=read(fd_,buf_,sizeof(buf_)))<0)
00119       {
00120         perror("Read of Proc file failed:");
00121         return procInfo();
00122       }
00123       
00124       if(cnt>0)
00125       {
00126         buf_[cnt]='\0';
00127         
00128         sscanf(buf_,
00129                "%d %s %c %d %d %d %d %d %u %u %u %u %u %d %d %d %d %d %d %u %u %d %u %u %u %u %u %u %u %u %d %d %d %d %u",
00130                &pinfo.pid, // %d
00131                pinfo.comm, // %s
00132                &pinfo.state, // %c
00133                &pinfo.ppid, // %d
00134                &pinfo.pgrp, // %d
00135                &pinfo.session, // %d
00136                &pinfo.tty, // %d
00137                &pinfo.tpgid, // %d
00138                &pinfo.flags, // %u
00139                &pinfo.minflt, // %u
00140                &pinfo.cminflt, // %u
00141                &pinfo.majflt, // %u
00142                &pinfo.cmajflt, // %u
00143                &pinfo.utime, // %d
00144                &pinfo.stime, // %d
00145                &pinfo.cutime, // %d
00146                &pinfo.cstime, // %d
00147                &pinfo.counter, // %d
00148                &pinfo.priority, // %d
00149                &pinfo.timeout, // %u
00150                &pinfo.itrealvalue, // %u
00151                &pinfo.starttime, // %d
00152                &pinfo.vsize, // %u
00153                &pinfo.rss, // %u
00154                &pinfo.rlim, // %u
00155                &pinfo.startcode, // %u
00156                &pinfo.endcode, // %u
00157                &pinfo.startstack, // %u
00158                &pinfo.kstkesp, // %u
00159                &pinfo.kstkeip, // %u
00160                &pinfo.signal, // %d
00161                &pinfo.blocked, // %d
00162                &pinfo.sigignore, // %d
00163                &pinfo.sigcatch, // %d
00164                &pinfo.wchan // %u
00165                );
00166         
00167         // resident set size in pages
00168         pr_size = (double)pinfo.vsize;
00169         pr_rssize = (double)pinfo.rss;
00170         
00171         ret.vsize = pr_size / (1024.0*1024.0);
00172         ret.rss   = (pr_rssize * pg_size_) / (1024.0*1024.0);
00173       }
00174 #else
00175       ret.vsize=0;
00176       ret.rss=0;
00177 #endif
00178       return ret;
00179     }
00180 
00181 
00182     double SimpleMemoryCheck::averageGrowthRate(double current, double past, int count) {
00183       return (current-past)/(double)count;
00184     }
00185     
00186     SimpleMemoryCheck::SimpleMemoryCheck(const ParameterSet& iPS,
00187                                          ActivityRegistry&iReg)
00188     : a_()
00189     , b_()
00190     , current_(&a_)
00191     , previous_(&b_)
00192     , pg_size_(sysconf(_SC_PAGESIZE)) // getpagesize()
00193     , num_to_skip_(iPS.getUntrackedParameter<int>("ignoreTotal"))
00194     , showMallocInfo_(iPS.getUntrackedParameter<bool>("showMallocInfo"))
00195     , oncePerEventMode_
00196         (iPS.getUntrackedParameter<bool>("oncePerEventMode"))
00197     , jobReportOutputOnly_(iPS.getUntrackedParameter<bool>("jobReportOutputOnly"))
00198     , count_()
00199     , growthRateVsize_()
00200     , growthRateRss_()
00201     , moduleSummaryRequested_
00202         (iPS.getUntrackedParameter<bool>("moduleMemorySummary"))
00203                                                                 // changelog 2
00204     {
00205       // pg_size = (double)getpagesize();
00206       std::ostringstream ost;
00207         
00208 #ifdef LINUX
00209       ost << "/proc/" << getpid() << "/stat";
00210       fname_ = ost.str();
00211       
00212       if((fd_=open(ost.str().c_str(),O_RDONLY))<0)
00213       {
00214         throw edm::Exception(errors::Configuration)
00215         << "Memory checker server: Failed to open " << ost.str() << std::endl;
00216       }
00217 #endif
00218       if (!oncePerEventMode_) { // default, prints on increases
00219         iReg.watchPreSourceConstruction(this,
00220              &SimpleMemoryCheck::preSourceConstruction);
00221         iReg.watchPostSourceConstruction(this,
00222              &SimpleMemoryCheck::postSourceConstruction);
00223         iReg.watchPostSource(this,
00224              &SimpleMemoryCheck::postSource);
00225         iReg.watchPostModuleConstruction(this,
00226              &SimpleMemoryCheck::postModuleConstruction);
00227         iReg.watchPostModuleBeginJob(this,
00228              &SimpleMemoryCheck::postModuleBeginJob);
00229         iReg.watchPostProcessEvent(this,
00230              &SimpleMemoryCheck::postEventProcessing);
00231         iReg.watchPostModule(this,
00232              &SimpleMemoryCheck::postModule);
00233         iReg.watchPostBeginJob(this,
00234              &SimpleMemoryCheck::postBeginJob);
00235         iReg.watchPostEndJob(this,
00236              &SimpleMemoryCheck::postEndJob);
00237       } else { 
00238         iReg.watchPostProcessEvent(this,
00239              &SimpleMemoryCheck::postEventProcessing);
00240         iReg.watchPostEndJob(this,
00241              &SimpleMemoryCheck::postEndJob);
00242       }
00243       if (moduleSummaryRequested_) {                            // changelog 2
00244         iReg.watchPreProcessEvent(this,
00245              &SimpleMemoryCheck::preEventProcessing);
00246         iReg.watchPreModule(this,
00247              &SimpleMemoryCheck::preModule);
00248         if (oncePerEventMode_) {
00249         iReg.watchPostModule(this,
00250              &SimpleMemoryCheck::postModule);
00251         }
00252       }
00253        
00254       // The following are not currenty used/implemented below for either
00255       // of the print modes (but are left here for reference)
00256       //  iReg.watchPostBeginJob(this,
00257       //       &SimpleMemoryCheck::postBeginJob);
00258       //  iReg.watchPreProcessEvent(this,
00259       //       &SimpleMemoryCheck::preEventProcessing);
00260       //  iReg.watchPreModule(this,
00261       //       &SimpleMemoryCheck::preModule);
00262 
00263       typedef edm::MallocOpts::opt_type opt_type;
00264       edm::MallocOptionSetter& mopts = edm::getGlobalOptionSetter();
00265       
00266       opt_type 
00267       p_mmap_max=iPS.getUntrackedParameter<int>("M_MMAP_MAX"),
00268       p_trim_thr=iPS.getUntrackedParameter<int>("M_TRIM_THRESHOLD"),
00269       p_top_pad=iPS.getUntrackedParameter<int>("M_TOP_PAD"),
00270       p_mmap_thr=iPS.getUntrackedParameter<int>("M_MMAP_THRESHOLD");
00271           
00272       if(p_mmap_max>=0) mopts.set_mmap_max(p_mmap_max);
00273       if(p_trim_thr>=0) mopts.set_trim_thr(p_trim_thr);
00274       if(p_top_pad>=0) mopts.set_top_pad(p_top_pad);
00275       if(p_mmap_thr>=0) mopts.set_mmap_thr(p_mmap_thr);
00276 
00277       mopts.adjustMallocParams();
00278 
00279       if(mopts.hasErrors())
00280       {
00281         LogWarning("MemoryCheck") 
00282         << "ERROR: Problem with setting malloc options\n"
00283         << mopts.error_message(); 
00284       }
00285       
00286       if(iPS.getUntrackedParameter<bool>("dump")==true)
00287       {
00288         edm::MallocOpts mo = mopts.get();
00289         LogWarning("MemoryCheck") 
00290         << "Malloc options: " << mo << "\n";
00291       }
00292     }
00293 
00294     SimpleMemoryCheck::~SimpleMemoryCheck()
00295     {
00296 #ifdef LINUX
00297       close(fd_);
00298 #endif
00299     }
00300 
00301     void SimpleMemoryCheck::fillDescriptions(edm::ConfigurationDescriptions & descriptions) {
00302       edm::ParameterSetDescription desc;
00303       desc.addUntracked<int>("ignoreTotal",1);
00304       desc.addUntracked<bool>("showMallocInfo",false);
00305       desc.addUntracked<bool>("oncePerEventMode",false);
00306       desc.addUntracked<bool>("jobReportOutputOnly",false);
00307       desc.addUntracked<bool>("moduleMemorySummary",false);
00308       desc.addUntracked<int>("M_MMAP_MAX",-1);
00309       desc.addUntracked<int>("M_TRIM_THRESHOLD",-1);
00310       desc.addUntracked<int>("M_TOP_PAD",-1);
00311       desc.addUntracked<int>("M_MMAP_THRESHOLD",-1);
00312       desc.addUntracked<bool>("dump",false);
00313       descriptions.add("SimpleMemoryCheck", desc);
00314     }
00315 
00316     void SimpleMemoryCheck::postBeginJob()
00317     {
00318         growthRateVsize_ = current_->vsize;
00319         growthRateRss_ = current_->rss;
00320     }
00321  
00322     void SimpleMemoryCheck::preSourceConstruction(const ModuleDescription& md) 
00323     {
00324       updateAndPrint("pre-ctor", md.moduleLabel(), md.moduleName());
00325     }
00326  
00327  
00328     void SimpleMemoryCheck::postSourceConstruction(const ModuleDescription& md)
00329     {
00330       updateAndPrint("ctor", md.moduleLabel(), md.moduleName());
00331     }
00332  
00333     void SimpleMemoryCheck::postSource() 
00334     {
00335       updateAndPrint("module", "source", "source");
00336     }
00337  
00338     void SimpleMemoryCheck::postModuleConstruction(const ModuleDescription& md)
00339     {
00340       updateAndPrint("ctor", md.moduleLabel(), md.moduleName());
00341     }
00342  
00343     void SimpleMemoryCheck::postModuleBeginJob(const ModuleDescription& md) 
00344     {
00345       updateAndPrint("beginJob", md.moduleLabel(), md.moduleName());
00346     }
00347  
00348     void SimpleMemoryCheck::postEndJob() 
00349     {
00350       if(not jobReportOutputOnly_) {
00351         edm::LogAbsolute("MemoryReport")                                // changelog 1
00352         << "MemoryReport> Peak virtual size " << eventT1_.vsize << " Mbytes" 
00353         << "\n"
00354         << " Key events increasing vsize: \n" 
00355         << eventL2_ << "\n"
00356         << eventL1_ << "\n"
00357         << eventM_  << "\n"
00358         << eventR1_ << "\n"
00359         << eventR2_ << "\n"
00360         << eventT3_ << "\n"
00361         << eventT2_ << "\n"
00362         << eventT1_ ;
00363       }
00364       if (moduleSummaryRequested_ and not jobReportOutputOnly_) {                               // changelog 1
00365         edm::LogAbsolute mmr("ModuleMemoryReport"); // at end of if block, mmr
00366                                                     // is destructed, causing
00367                                                     // message to be logged
00368         mmr << "ModuleMemoryReport> Each line has module label and: \n";
00369         mmr << "  (after early ignored events) \n"; 
00370         mmr << 
00371         "    count of times module executed; average increase in vsize \n";
00372         mmr << 
00373         "    maximum increase in vsize; event on which maximum occurred \n";
00374         mmr << "  (during early ignored events) \n";
00375         mmr << "    total and maximum vsize increases \n \n";   
00376         for (SignificantModulesMap::iterator im=modules_.begin(); 
00377              im != modules_.end(); ++im) {
00378           SignificantModule const& m = im->second;
00379           if ( m.totalDeltaVsize == 0 && m.totalEarlyVsize == 0 ) continue;
00380           mmr << im->first << ": ";
00381           mmr << "n = " << m.postEarlyCount;
00382           if ( m.postEarlyCount > 0 ) mmr << " avg = " 
00383             << m.totalDeltaVsize/m.postEarlyCount;
00384           mmr <<  " max = " << m.maxDeltaVsize << " " << m.eventMaxDeltaV;
00385           if ( m.totalEarlyVsize > 0 ) {
00386             mmr << " early total: " << m.totalEarlyVsize;
00387             mmr << " max: " << m.maxEarlyVsize;
00388           }
00389           mmr << "\n";
00390         }
00391       } // end of if; mmr goes out of scope; log message is queued
00392 
00393       Service<JobReport> reportSvc;
00394                                                                 // changelog 1
00395 #define SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
00396 #ifdef  SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
00397 //     std::map<std::string, double> reportData;
00398       std::map<std::string, std::string> reportData;
00399 
00400       if (eventL2_.vsize > 0) 
00401         eventStatOutput("LargeVsizeIncreaseEventL2", eventL2_, reportData);
00402       if (eventL1_.vsize > 0) 
00403         eventStatOutput("LargeVsizeIncreaseEventL1", eventL1_, reportData);
00404       if (eventM_.vsize > 0) 
00405         eventStatOutput("LargestVsizeIncreaseEvent", eventM_,  reportData);
00406       if (eventR1_.vsize > 0) 
00407         eventStatOutput("LargeVsizeIncreaseEventR1", eventR1_, reportData);
00408       if (eventR2_.vsize > 0)
00409         eventStatOutput("LargeVsizeIncreaseEventR2", eventR2_, reportData);
00410       if (eventT3_.vsize > 0) 
00411         eventStatOutput("ThirdLargestVsizeEventT3",  eventT3_, reportData);
00412       if (eventT2_.vsize > 0) 
00413         eventStatOutput("SecondLargestVsizeEventT2", eventT2_, reportData);
00414       if (eventT1_.vsize > 0)
00415         eventStatOutput("LargestVsizeEventT1",       eventT1_, reportData);
00416 
00417       if (eventRssT3_.rss > 0)
00418         eventStatOutput("ThirdLargestRssEvent", eventRssT3_, reportData);
00419       if (eventRssT2_.rss > 0)
00420         eventStatOutput("SecondLargestRssEvent", eventRssT2_, reportData);
00421       if (eventRssT1_.rss > 0)
00422         eventStatOutput("LargestRssEvent", eventRssT1_, reportData);
00423       if (eventDeltaRssT3_.deltaRss > 0)
00424         eventStatOutput("ThirdLargestIncreaseRssEvent", eventDeltaRssT3_, reportData);
00425       if (eventDeltaRssT2_.deltaRss > 0)
00426         eventStatOutput("SecondLargestIncreaseRssEvent", eventDeltaRssT2_, reportData);
00427       if (eventDeltaRssT1_.deltaRss > 0)
00428         eventStatOutput("LargestIncreaseRssEvent", eventDeltaRssT1_, reportData);
00429      
00430 #ifdef __linux__
00431       struct mallinfo minfo = mallinfo();
00432       reportData.insert(
00433         std::make_pair("HEAP_ARENA_SIZE_BYTES", i2str(minfo.arena)));  
00434       reportData.insert(
00435         std::make_pair("HEAP_ARENA_N_UNUSED_CHUNKS", i2str(minfo.ordblks)));  
00436       reportData.insert(
00437         std::make_pair("HEAP_TOP_FREE_BYTES", i2str(minfo.keepcost)));  
00438       reportData.insert(
00439         std::make_pair("HEAP_MAPPED_SIZE_BYTES", i2str(minfo.hblkhd)));  
00440       reportData.insert(
00441         std::make_pair("HEAP_MAPPED_N_CHUNKS", i2str(minfo.hblks)));  
00442       reportData.insert(
00443         std::make_pair("HEAP_USED_BYTES", i2str(minfo.uordblks)));  
00444       reportData.insert(
00445         std::make_pair("HEAP_UNUSED_BYTES", i2str(minfo.fordblks)));  
00446 #endif
00447 
00448       // Report Growth rates for VSize and Rss
00449       reportData.insert(
00450         std::make_pair("AverageGrowthRateVsize", d2str(averageGrowthRate(current_->vsize, growthRateVsize_, count_))));     
00451       reportData.insert(
00452         std::make_pair("PeakValueVsize", d2str(eventT1_.vsize)));
00453       reportData.insert(
00454         std::make_pair("AverageGrowthRateRss", d2str(averageGrowthRate(current_->rss, growthRateRss_, count_))));
00455       reportData.insert(
00456         std::make_pair("PeakValueRss", d2str(eventRssT1_.rss)));
00457 
00458       if (moduleSummaryRequested_) {                            // changelog 2
00459         for (SignificantModulesMap::iterator im=modules_.begin(); 
00460              im != modules_.end(); ++im) {
00461           SignificantModule const& m = im->second;
00462           if ( m.totalDeltaVsize == 0 && m.totalEarlyVsize == 0 ) continue;
00463           std::string label = im->first+":";
00464           reportData.insert(
00465                             std::make_pair(label+"PostEarlyCount", i2str(m.postEarlyCount)));  
00466           if ( m.postEarlyCount > 0 ) {
00467             reportData.insert(
00468                               std::make_pair(label+"AverageDeltaVsize", 
00469                                              d2str(m.totalDeltaVsize/m.postEarlyCount)));  
00470           }
00471           reportData.insert(
00472                             std::make_pair(label+"MaxDeltaVsize", d2str(m.maxDeltaVsize)));  
00473           if ( m.totalEarlyVsize > 0 ) {
00474             reportData.insert(
00475                               std::make_pair(label+"TotalEarlyVsize", d2str(m.totalEarlyVsize)));  
00476             reportData.insert(
00477                               std::make_pair(label+"MaxEarlyDeltaVsize", d2str(m.maxEarlyVsize)));  
00478           }
00479         }
00480       } 
00481 
00482       std::map<std::string, std::string> reportMemoryProperties;
00483 
00484       if (FILE *fmeminfo = fopen("/proc/meminfo", "r")){
00485         char buf[128];
00486         char space[] = " ";
00487         size_t value;
00488         
00489         while (fgets(buf, sizeof(buf), fmeminfo)){
00490           char *token = NULL;
00491           token = strtok(buf, space);
00492           if (token != NULL){
00493             value = atol(strtok(NULL, space));
00494             std::string category = token;
00495             reportMemoryProperties.insert(std::make_pair(category.substr(0,strlen(token)-1), i2str(value)));
00496           }
00497         }
00498         
00499         fclose(fmeminfo);
00500       }
00501         
00502 //      reportSvc->reportMemoryInfo(reportData, reportMemoryProperties);
00503       reportSvc->reportPerformanceSummary("ApplicationMemory", reportData);
00504       reportSvc->reportPerformanceSummary("SystemMemory", reportMemoryProperties);
00505 #endif
00506 
00507 #ifdef SIMPLE_MEMORY_CHECK_DIFFERENT_XML_OUTPUT
00508       std::vector<std::string> reportData;
00509 
00510       if (eventL2_.vsize > 0) reportData.push_back(
00511         eventStatOutput("LargeVsizeIncreaseEventL2", eventL2_));
00512       if (eventL1_.vsize > 0) reportData.push_back(
00513         eventStatOutput("LargeVsizeIncreaseEventL1", eventL1_));
00514       if (eventM_.vsize > 0) reportData.push_back(
00515         eventStatOutput("LargestVsizeIncreaseEvent", eventM_));
00516       if (eventR1_.vsize > 0) reportData.push_back(
00517         eventStatOutput("LargeVsizeIncreaseEventR1", eventR1_));
00518       if (eventR2_.vsize > 0) reportData.push_back(
00519         eventStatOutput("LargeVsizeIncreaseEventR2", eventR2_));
00520       if (eventT3_.vsize > 0) reportData.push_back(
00521         eventStatOutput("ThirdLargestVsizeEventT3", eventT3_));
00522       if (eventT2_.vsize > 0) reportData.push_back(
00523         eventStatOutput("SecondLargestVsizeEventT2", eventT2_));
00524       if (eventT1_.vsize > 0) reportData.push_back(
00525         eventStatOutput("LargestVsizeEventT1", eventT1_));
00526 
00527       if (eventRssT3_.rss > 0)
00528         eventStatOutput("ThirdLargestRssEvent", eventRssT3_, reportData);
00529       if (eventRssT2_.rss > 0)
00530         eventStatOutput("SecondLargestRssEvent", eventRssT2_, reportData);
00531       if (eventRssT1_.rss > 0)
00532         eventStatOutput("LargestRssEvent", eventRssT1_, reportData);
00533       if (eventDeltaRssT3_.deltaRss > 0)
00534         eventStatOutput("ThirdLargestIncreaseRssEvent", eventDeltaRssT3_, reportData);
00535       if (eventDeltaRssT2_.deltaRss > 0)
00536         eventStatOutput("SecondLargestIncreaseRssEvent", eventDeltaRssT2_, reportData);
00537       if (eventDeltaRssT1_.deltaRss > 0)
00538         eventStatOutput("LargestIncreaseRssEvent", eventDeltaRssT1_, reportData);
00539       
00540       struct mallinfo minfo = mallinfo();
00541       reportData.push_back(
00542         mallOutput("HEAP_ARENA_SIZE_BYTES", minfo.arena));  
00543       reportData.push_back(
00544         mallOutput("HEAP_ARENA_N_UNUSED_CHUNKS", minfo.ordblks));  
00545       reportData.push_back(
00546         mallOutput("HEAP_TOP_FREE_BYTES", minfo.keepcost));  
00547       reportData.push_back(
00548         mallOutput("HEAP_MAPPED_SIZE_BYTES", minfo.hblkhd));  
00549       reportData.push_back(
00550         mallOutput("HEAP_MAPPED_N_CHUNKS", minfo.hblks));  
00551       reportData.push_back(
00552         mallOutput("HEAP_USED_BYTES", minfo.uordblks));  
00553       reportData.push_back(
00554         mallOutput("HEAP_UNUSED_BYTES", minfo.fordblks));  
00555         
00556       // Report Growth rates for VSize and Rss
00557       reportData.insert(
00558         std::make_pair("AverageGrowthRateVsize", d2str(averageGrowthRate(current_->vsize, growthRateVsize_, count_))));
00559       reportData.insert(
00560         std::make_pair("PeakValueVsize", d2str(eventT1_.vsize)));
00561       reportData.insert(
00562         std::make_pair("AverageGrowthRateRss", d2str(averageGrowthRate(current_->rss, growthRateRss_, count_))));
00563       reportData.insert(
00564         std::make_pair("PeakValueRss", d2str(eventRssT1_.rss)));
00565 
00566       reportSvc->reportMemoryInfo(reportData);
00567       // This is a form of reportMemoryInfo taking s vector, not a map
00568 #endif
00569     } // postEndJob
00570  
00571     void SimpleMemoryCheck::preEventProcessing(const edm::EventID& iID,
00572                                                const edm::Timestamp& iTime) 
00573     {
00574       currentEventID_ = iID;                                    // changelog 2
00575     }
00576 
00577     void SimpleMemoryCheck::postEventProcessing(const Event& e,
00578                                                 const EventSetup&) 
00579     {
00580       ++count_;
00581       update();
00582       updateEventStats( e.id() );
00583       if (oncePerEventMode_) {
00584         // should probably use be Run:Event or count_ for the label and name
00585         updateMax();
00586         andPrint("event", "", ""); 
00587       } 
00588     }
00589  
00590     void SimpleMemoryCheck::preModule(const ModuleDescription& md) { 
00591       update();                                                 // changelog 2
00592       moduleEntryVsize_ = current_->vsize;
00593     }
00594  
00595     void SimpleMemoryCheck::postModule(const ModuleDescription& md) {
00596       if (!oncePerEventMode_) {
00597         updateAndPrint("module", md.moduleLabel(), md.moduleName());
00598       } else if (moduleSummaryRequested_) {                     // changelog 2
00599         update();
00600       }
00601       if (moduleSummaryRequested_) {                            // changelog 2
00602         double dv = current_->vsize - moduleEntryVsize_;
00603         std::string label =  md.moduleLabel();
00604         updateModuleMemoryStats (modules_[label],dv);
00605       }
00606     }
00607  
00608  
00609     void SimpleMemoryCheck::update() 
00610     {
00611       std::swap(current_,previous_);
00612       *current_ = fetch();
00613     }
00614 
00615     void SimpleMemoryCheck::updateMax() 
00616     {
00617       if ((*current_ > max_) || oncePerEventMode_)
00618         {
00619           if(count_ >= num_to_skip_) {
00620           }
00621           max_ = *current_;
00622         }
00623     }
00624 
00625     void SimpleMemoryCheck::updateEventStats(edm::EventID const & e) {
00626       if (count_ < num_to_skip_) return;
00627       if (count_ == num_to_skip_) {
00628         eventT1_.set(0, 0, e, this);
00629         eventM_.set (0, 0, e, this);
00630         eventRssT1_.set(0, 0, e, this);
00631         eventDeltaRssT1_.set(0, 0, e, this);
00632         return;
00633       }
00634       double vsize = current_->vsize;
00635       double deltaVsize = vsize - eventT1_.vsize;
00636       
00637       // Update significative events for Vsize
00638       if (vsize > eventT1_.vsize) {
00639         double deltaRss = current_->rss - eventT1_.rss;
00640         eventT3_ = eventT2_;
00641         eventT2_ = eventT1_;
00642         eventT1_.set(deltaVsize, deltaRss, e, this);
00643       } else if(vsize > eventT2_.vsize) {
00644         double deltaRss = current_->rss - eventT1_.rss;
00645         eventT3_ = eventT2_;
00646         eventT2_.set(deltaVsize, deltaRss, e, this);
00647       } else if(vsize > eventT3_.vsize) {
00648         double deltaRss = current_->rss - eventT1_.rss;
00649         eventT3_.set(deltaVsize, deltaRss, e, this);
00650       }
00651       
00652       if (deltaVsize > eventM_.deltaVsize) {
00653         double deltaRss = current_->rss - eventM_.rss;
00654         if (eventL1_.deltaVsize >= eventR1_.deltaVsize) {
00655           eventL2_ = eventL1_; 
00656         } else {
00657           eventL2_ = eventR1_; 
00658         }
00659         eventL1_ = eventM_;
00660         eventM_.set(deltaVsize, deltaRss, e, this);
00661         eventR1_ = SignificantEvent();
00662         eventR2_ = SignificantEvent();
00663       } else if (deltaVsize > eventR1_.deltaVsize) {
00664         double deltaRss = current_->rss - eventM_.rss;
00665         eventR2_ = eventR1_;
00666         eventR1_.set(deltaVsize, deltaRss, e, this);
00667       } else if (deltaVsize > eventR2_.deltaVsize) {
00668         double deltaRss = current_->rss - eventR1_.rss;
00669         eventR2_.set(deltaVsize, deltaRss, e, this);
00670       }
00671 
00672       // Update significative events for Rss
00673       double rss = current_->rss;
00674       double deltaRss = rss - eventRssT1_.rss;
00675 
00676       if(rss > eventRssT1_.rss){
00677         eventRssT3_ = eventRssT2_;
00678         eventRssT2_ = eventRssT1_;
00679         eventRssT1_.set(deltaVsize, deltaRss, e, this);
00680       } else if(rss > eventRssT2_.rss) {
00681         eventRssT3_ = eventRssT2_;
00682         eventRssT2_.set(deltaVsize, deltaRss, e, this);
00683       } else if(rss > eventRssT3_.rss) {
00684         eventRssT3_.set(deltaVsize, deltaRss, e, this);
00685       }
00686       if(deltaRss > eventDeltaRssT1_.deltaRss) {
00687         eventDeltaRssT3_ = eventDeltaRssT2_;
00688         eventDeltaRssT2_ = eventDeltaRssT1_;
00689         eventDeltaRssT1_.set(deltaVsize, deltaRss, e, this);
00690       } else if(deltaRss > eventDeltaRssT2_.deltaRss) {
00691         eventDeltaRssT3_ = eventDeltaRssT2_;
00692         eventDeltaRssT2_.set(deltaVsize, deltaRss, e, this);
00693       } else if(deltaRss > eventDeltaRssT3_.deltaRss) {
00694         eventDeltaRssT3_.set(deltaVsize, deltaRss, e, this);
00695       }
00696     }   // updateEventStats
00697       
00698     void SimpleMemoryCheck::andPrint(const std::string& type, 
00699                     const std::string& mdlabel, const std::string& mdname) const
00700     {
00701       if (not jobReportOutputOnly_ && ((*current_ > max_) || oncePerEventMode_)) {
00702         if(count_ >= num_to_skip_) {
00703           double deltaVSIZE = current_->vsize - max_.vsize;
00704           double deltaRSS   = current_->rss - max_.rss;
00705           if (!showMallocInfo_) {  // default
00706             LogWarning("MemoryCheck")
00707             << "MemoryCheck: " << type << " "
00708             << mdname << ":" << mdlabel 
00709             << " VSIZE " << current_->vsize << " " << deltaVSIZE
00710             << " RSS " << current_->rss << " " << deltaRSS
00711             << "\n";
00712           } else {
00713 #ifdef __linux__
00714             struct mallinfo minfo = mallinfo();
00715 #endif
00716             LogWarning("MemoryCheck")
00717             << "MemoryCheck: " << type << " "
00718             << mdname << ":" << mdlabel 
00719             << " VSIZE " << current_->vsize << " " << deltaVSIZE
00720             << " RSS " << current_->rss << " " << deltaRSS
00721 #ifdef __linux__
00722             << " HEAP-ARENA [ SIZE-BYTES " << minfo.arena
00723             << " N-UNUSED-CHUNKS " << minfo.ordblks
00724             << " TOP-FREE-BYTES " << minfo.keepcost << " ]"
00725             << " HEAP-MAPPED [ SIZE-BYTES " << minfo.hblkhd
00726             << " N-CHUNKS " << minfo.hblks << " ]"
00727             << " HEAP-USED-BYTES " << minfo.uordblks
00728             << " HEAP-UNUSED-BYTES " << minfo.fordblks
00729 #endif
00730             << "\n";
00731           }
00732         }
00733       }
00734     }
00735 
00736     void SimpleMemoryCheck::updateAndPrint(const std::string& type, 
00737                     const std::string& mdlabel, const std::string& mdname) 
00738     {
00739       update();
00740       andPrint(type, mdlabel, mdname);
00741       updateMax();
00742     }
00743 
00744 #ifdef SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
00745     void
00746     SimpleMemoryCheck::eventStatOutput(std::string title, 
00747                                        SignificantEvent const& e,
00748                                        std::map<std::string, std::string> &m) const
00749     {
00750       { std::ostringstream os;
00751         os << title << "-a-COUNT";
00752         m.insert(std::make_pair(os.str(), i2str(e.count))); }
00753       { std::ostringstream os;
00754         os << title << "-b-RUN";
00755         m.insert(std::make_pair(os.str(), d2str(static_cast<double>(e.event.run())) )); }
00756       { std::ostringstream os;
00757         os << title << "-c-EVENT";
00758         m.insert(std::make_pair(os.str(), d2str(static_cast<double>(e.event.event())) )); }
00759       { std::ostringstream os;
00760         os << title << "-d-VSIZE";
00761         m.insert(std::make_pair(os.str(), d2str(e.vsize))); }
00762       { std::ostringstream os;
00763         os << title << "-e-DELTV";
00764         m.insert(std::make_pair(os.str(), d2str(e.deltaVsize))); }
00765       { std::ostringstream os;
00766         os << title << "-f-RSS";
00767         m.insert(std::make_pair(os.str(), d2str(e.rss))); }
00768     } // eventStatOutput
00769 #endif
00770 
00771  
00772 #ifdef SIMPLE_MEMORY_CHECK_DIFFERENT_XML_OUTPUT
00773     std::string 
00774     SimpleMemoryCheck::eventStatOutput(std::string title, 
00775                                        SignificantEvent const& e) const
00776     {
00777       std::ostringstream os;
00778       os << "  <" << title << ">\n";
00779       os << "    " << e.count << ": " << e.event;
00780       os << " vsize " << e.vsize-e.deltaVsize << " + " << e.deltaVsize
00781                                               << " = " << e.vsize;
00782       os << "  rss: " << e.rss << "\n";                                       
00783       os << "  </" << title << ">\n";
00784       return os.str();
00785     } // eventStatOutput
00786 
00787     std::string 
00788     SimpleMemoryCheck::mallOutput(std::string title, size_t const& n) const {
00789       std::ostringstream os;
00790       os << "  <" << title << ">\n";
00791       os << "    " << n << "\n";
00792       os << "  </" << title << ">\n";
00793       return os.str();
00794     }
00795 #endif
00796                                                                 // changelog 2
00797     void 
00798     SimpleMemoryCheck::updateModuleMemoryStats(SignificantModule & m, 
00799                                                double dv) {
00800       if(count_ < num_to_skip_) {
00801         m.totalEarlyVsize += dv;
00802         if (dv > m.maxEarlyVsize)  m.maxEarlyVsize = dv;        
00803       } else {
00804         ++m.postEarlyCount;
00805         m.totalDeltaVsize += dv;
00806         if (dv > m.maxDeltaVsize)  {
00807           m.maxDeltaVsize = dv;    
00808           m.eventMaxDeltaV = currentEventID_;
00809         }       
00810       }
00811     } //updateModuleMemoryStats
00812  
00813 
00814 
00815     std::ostream & 
00816     operator<< (std::ostream & os, 
00817                 SimpleMemoryCheck::SignificantEvent const & se) {
00818       os << "[" << se.count << "] "
00819       << se.event << "  vsize = " << se.vsize 
00820       << " deltaVsize = " << se.deltaVsize 
00821       << " rss = " << se.rss << " delta " << se.deltaRss;
00822       return os;
00823     }
00824     
00825     std::ostream & 
00826     operator<< (std::ostream & os, 
00827                 SimpleMemoryCheck::SignificantModule const & sm) {
00828       if ( sm.postEarlyCount > 0 ) {
00829         os << "\nPost Early Events:  TotalDeltaVsize: " << sm.totalDeltaVsize
00830         << " (avg: " << sm.totalDeltaVsize/sm.postEarlyCount 
00831         << "; max: " << sm.maxDeltaVsize  
00832         << " during " << sm.eventMaxDeltaV << ")";
00833       }
00834       if ( sm.totalEarlyVsize > 0 ) {
00835         os << "\n     Early Events:  TotalDeltaVsize: " << sm.totalEarlyVsize
00836         << " (max: " << sm.maxEarlyVsize << ")";
00837       }
00838       
00839       return os;
00840     }
00841 
00842   } // end namespace service
00843 } // end namespace edm
00844