Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:13:05

0001 // -*- C++ -*-
0002 //
0003 // Package:     Services
0004 // Class  :     Memory
0005 //
0006 // Implementation:
0007 //
0008 // Original Author:  Jim Kowalkowski
0009 //
0010 // Change Log
0011 //
0012 // 1 - Apr 25, 2008 M. Fischler
0013 //        Collect event summary information and output to XML file and logger
0014 //        at the end of the job.  Involves split-up of updateAndPrint method.
0015 //
0016 // 2 - May 7, 2008 M. Fischler
0017 //      Collect module summary information and output to XML file and logger
0018 //        at the end of the job.
0019 //
0020 // 3 - Jan 14, 2009 Natalia Garcia Nebot
0021 //        Added:        - Average rate of growth in RSS and peak value attained.
0022 //                - Average rate of growth in VSize over time, Peak VSize
0023 //
0024 //
0025 #include "FWCore/ServiceRegistry/interface/ServiceMaker.h"
0026 
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/ServiceRegistry/interface/ActivityRegistry.h"
0029 #include "DataFormats/Provenance/interface/EventID.h"
0030 #include "FWCore/Services/plugins/ProcInfoFetcher.h"
0031 
0032 #include "DataFormats/Provenance/interface/ModuleDescription.h"
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/MessageLogger/interface/JobReport.h"
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0038 #include "FWCore/ServiceRegistry/interface/Service.h"
0039 #include "FWCore/ServiceRegistry/interface/StreamContext.h"
0040 #include "FWCore/ServiceRegistry/interface/ModuleCallingContext.h"
0041 #include "FWCore/Utilities/interface/Exception.h"
0042 #include "FWCore/Utilities/interface/get_underlying_safe.h"
0043 
0044 #include <cstring>
0045 #include <iostream>
0046 #include <memory>
0047 #ifdef __linux__
0048 #include <malloc.h>
0049 #endif
0050 #include <sstream>
0051 //#include <stdio.h>
0052 #include <string>
0053 //#include <string.h>
0054 
0055 #include <cstdio>
0056 #include <atomic>
0057 
0058 namespace edm {
0059   class EventID;
0060   class Timestamp;
0061 
0062   namespace service {
0063     struct smapsInfo {
0064       smapsInfo() : private_(), pss_() {}
0065       smapsInfo(double private_sz, double pss_sz) : private_(private_sz), pss_(pss_sz) {}
0066 
0067       bool operator==(const smapsInfo& p) const { return private_ == p.private_ && pss_ == p.pss_; }
0068 
0069       bool operator>(const smapsInfo& p) const { return private_ > p.private_ || pss_ > p.pss_; }
0070 
0071       double private_;  // in MB
0072       double pss_;      // in MB
0073     };
0074 
0075     class SimpleMemoryCheck {
0076     public:
0077       SimpleMemoryCheck(const ParameterSet&, ActivityRegistry&);
0078       ~SimpleMemoryCheck();
0079 
0080       static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0081 
0082       void preSourceConstruction(const ModuleDescription&);
0083       void postSourceConstruction(const ModuleDescription&);
0084       void postSourceEvent(StreamID);
0085 
0086       void postBeginJob();
0087 
0088       void postEvent(StreamContext const&);
0089 
0090       void postModuleBeginJob(const ModuleDescription&);
0091       void postModuleConstruction(const ModuleDescription&);
0092 
0093       void preModule(StreamContext const&, ModuleCallingContext const&);
0094       void postModule(StreamContext const&, ModuleCallingContext const&);
0095 
0096       void postEndJob();
0097 
0098     private:
0099       ProcInfo fetch();
0100       smapsInfo fetchSmaps();
0101       double pageSize() const { return pg_size_; }
0102       double averageGrowthRate(double current, double past, int count);
0103       void update();
0104       void updateMax();
0105       void andPrint(const std::string& type, const std::string& mdlabel, const std::string& mdname) const;
0106       void updateAndPrint(const std::string& type, const std::string& mdlabel, const std::string& mdname);
0107       void openFiles();
0108 
0109       char const* smapsLineBuffer() const { return get_underlying_safe(smapsLineBuffer_); }
0110       char*& smapsLineBuffer() { return get_underlying_safe(smapsLineBuffer_); }
0111 
0112       ProcInfo a_;
0113       ProcInfo b_;
0114       ProcInfo max_;
0115       edm::propagate_const<ProcInfo*> current_;
0116       edm::propagate_const<ProcInfo*> previous_;
0117 
0118       smapsInfo currentSmaps_;
0119 
0120       ProcInfoFetcher piFetcher_;
0121       double pg_size_;
0122       int num_to_skip_;
0123       //options
0124       bool showMallocInfo_;
0125       bool oncePerEventMode_;
0126       bool jobReportOutputOnly_;
0127       bool monitorPssAndPrivate_;
0128       std::atomic<int> count_;
0129 
0130       //smaps
0131       edm::propagate_const<FILE*> smapsFile_;
0132       edm::propagate_const<char*> smapsLineBuffer_;
0133       size_t smapsLineBufferLen_;
0134 
0135       //Rates of growth
0136       double growthRateVsize_;
0137       double growthRateRss_;
0138 
0139       // Event summary statistics               changeLog 1
0140       struct SignificantEvent {
0141         int count;
0142         double vsize;
0143         double deltaVsize;
0144         double rss;
0145         double deltaRss;
0146         bool monitorPssAndPrivate;
0147         double privateSize;
0148         double pss;
0149         edm::EventID event;
0150         SignificantEvent()
0151             : count(0),
0152               vsize(0),
0153               deltaVsize(0),
0154               rss(0),
0155               deltaRss(0),
0156               monitorPssAndPrivate(false),
0157               privateSize(0),
0158               pss(0),
0159               event() {}
0160         void set(double deltaV, double deltaR, edm::EventID const& e, SimpleMemoryCheck* t) {
0161           count = t->count_;
0162           vsize = t->current_->vsize;
0163           deltaVsize = deltaV;
0164           rss = t->current_->rss;
0165           deltaRss = deltaR;
0166           monitorPssAndPrivate = t->monitorPssAndPrivate_;
0167           if (monitorPssAndPrivate) {
0168             privateSize = t->currentSmaps_.private_;
0169             pss = t->currentSmaps_.pss_;
0170           }
0171           event = e;
0172         }
0173       };  // SignificantEvent
0174       friend struct SignificantEvent;
0175       friend std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantEvent const& se);
0176 
0177       /*
0178        Significative events for deltaVsize:
0179        - eventM_: Event which makes the biggest value for deltaVsize
0180        - eventL1_: Event which makes the second biggest value for deltaVsize
0181        - eventL2_: Event which makes the third biggest value for deltaVsize
0182        - eventR1_: Event which makes the second biggest value for deltaVsize
0183        - eventR2_: Event which makes the third biggest value for deltaVsize
0184        M>L1>L2 and M>R1>R2
0185        Unknown relation between Ls and Rs events ???????
0186        Significative events for vsize:
0187        - eventT1_: Event whith the biggest value for vsize
0188        - eventT2_: Event whith the second biggest value for vsize
0189        - eventT3_: Event whith the third biggest value for vsize
0190        T1>T2>T3
0191        */
0192       SignificantEvent eventM_;
0193       SignificantEvent eventL1_;
0194       SignificantEvent eventL2_;
0195       SignificantEvent eventR1_;
0196       SignificantEvent eventR2_;
0197       SignificantEvent eventT1_;
0198       SignificantEvent eventT2_;
0199       SignificantEvent eventT3_;
0200 
0201       /*
0202        Significative event for deltaRss:
0203        - eventRssT1_: Event whith the biggest value for rss
0204        - eventRssT2_: Event whith the second biggest value for rss
0205        - eventRssT3_: Event whith the third biggest value for rss
0206        T1>T2>T3
0207        Significative events for deltaRss:
0208        - eventDeltaRssT1_: Event whith the biggest value for deltaRss
0209        - eventDeltaRssT2_: Event whith the second biggest value for deltaRss
0210        - eventDeltaRssT3_: Event whith the third biggest value for deltaRss
0211        T1>T2>T3
0212        */
0213       SignificantEvent eventRssT1_;
0214       SignificantEvent eventRssT2_;
0215       SignificantEvent eventRssT3_;
0216       SignificantEvent eventDeltaRssT1_;
0217       SignificantEvent eventDeltaRssT2_;
0218       SignificantEvent eventDeltaRssT3_;
0219 
0220       void updateEventStats(edm::EventID const& e);
0221       std::string eventStatOutput(std::string title, SignificantEvent const& e) const;
0222       void eventStatOutput(std::string title, SignificantEvent const& e, std::map<std::string, std::string>& m) const;
0223       std::string mallOutput(std::string title, size_t const& n) const;
0224 
0225       // Module summary statistices
0226       struct SignificantModule {
0227         int postEarlyCount;
0228         double totalDeltaVsize;
0229         double maxDeltaVsize;
0230         edm::EventID eventMaxDeltaV;
0231         double totalEarlyVsize;
0232         double maxEarlyVsize;
0233         SignificantModule()
0234             : postEarlyCount(0),
0235               totalDeltaVsize(0),
0236               maxDeltaVsize(0),
0237               eventMaxDeltaV(),
0238               totalEarlyVsize(0),
0239               maxEarlyVsize(0) {}
0240         void set(double deltaV, bool early);
0241       };  // SignificantModule
0242       friend struct SignificantModule;
0243       friend std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantModule const& se);
0244       bool moduleSummaryRequested_;
0245       typedef std::map<std::string, SignificantModule> SignificantModulesMap;
0246       SignificantModulesMap modules_;
0247       double moduleEntryVsize_;
0248       void updateModuleMemoryStats(SignificantModule& m, double dv, edm::EventID const&);
0249 
0250       //Used to guarantee we only do one measurement at a time
0251       std::atomic<bool> measurementUnderway_;
0252       std::atomic<bool> moduleMeasurementUnderway_;
0253       std::atomic<unsigned int> moduleStreamID_;
0254       std::atomic<unsigned int> moduleID_;
0255 
0256     };  // SimpleMemoryCheck
0257 
0258     std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantEvent const& se);
0259 
0260     std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantModule const& se);
0261 
0262   }  // namespace service
0263 }  // namespace edm
0264 
0265 #ifdef __linux__
0266 #define LINUX 1
0267 #endif
0268 
0269 #include <fcntl.h>
0270 #include <unistd.h>
0271 
0272 namespace edm {
0273   namespace service {
0274 
0275     static std::string d2str(double d) {
0276       std::ostringstream t;
0277       t << d;
0278       return t.str();
0279     }
0280 
0281     static std::string i2str(int i) {
0282       std::ostringstream t;
0283       t << i;
0284       return t.str();
0285     }
0286 
0287     ProcInfo SimpleMemoryCheck::fetch() { return piFetcher_.fetch(); }
0288 
0289     smapsInfo SimpleMemoryCheck::fetchSmaps() {
0290       smapsInfo ret;
0291       ret.private_ = 0;
0292       ret.pss_ = 0;
0293 #ifdef LINUX
0294       fseek(smapsFile_, 0, SEEK_SET);
0295       ssize_t read;
0296 
0297       /*
0298        The format of the report is
0299        Private_Clean:        0 kB
0300        Private_Dirty:       72 kB
0301        Swap:                 0 kB
0302        Pss:                 72 kB
0303        */
0304 
0305       while ((read = getline(&smapsLineBuffer(), &smapsLineBufferLen_, smapsFile_)) != -1) {
0306         if (read > 14) {
0307           //Private
0308           if (0 == strncmp("Private_", smapsLineBuffer_, 8)) {
0309             unsigned int value = atoi(smapsLineBuffer_ + 14);
0310             //Convert from kB to MB
0311             ret.private_ += static_cast<double>(value) / 1024.;
0312           } else if (0 == strncmp("Pss:", smapsLineBuffer_, 4)) {
0313             unsigned int value = atoi(smapsLineBuffer_ + 4);
0314             //Convert from kB to MB
0315             ret.pss_ += static_cast<double>(value) / 1024.;
0316           }
0317         }
0318       }
0319 #endif
0320       return ret;
0321     }
0322 
0323     double SimpleMemoryCheck::averageGrowthRate(double current, double past, int count) {
0324       return (current - past) / (double)count;
0325     }
0326 
0327     SimpleMemoryCheck::SimpleMemoryCheck(ParameterSet const& iPS, ActivityRegistry& iReg)
0328         : a_(),
0329           b_(),
0330           current_(&a_),
0331           previous_(&b_),
0332           pg_size_(sysconf(_SC_PAGESIZE))  // getpagesize()
0333           ,
0334           num_to_skip_(iPS.getUntrackedParameter<int>("ignoreTotal")),
0335           showMallocInfo_(iPS.getUntrackedParameter<bool>("showMallocInfo")),
0336           oncePerEventMode_(iPS.getUntrackedParameter<bool>("oncePerEventMode")),
0337           jobReportOutputOnly_(iPS.getUntrackedParameter<bool>("jobReportOutputOnly")),
0338           monitorPssAndPrivate_(iPS.getUntrackedParameter<bool>("monitorPssAndPrivate")),
0339           count_(),
0340           smapsFile_(nullptr),
0341           smapsLineBuffer_(nullptr),
0342           smapsLineBufferLen_(0),
0343           growthRateVsize_(),
0344           growthRateRss_(),
0345           moduleSummaryRequested_(iPS.getUntrackedParameter<bool>("moduleMemorySummary")),
0346           measurementUnderway_(false) {
0347       // changelog 2
0348       // pg_size = (double)getpagesize();
0349       std::ostringstream ost;
0350 
0351       openFiles();
0352 
0353       if (!oncePerEventMode_) {  // default, prints on increases
0354         iReg.watchPreSourceConstruction(this, &SimpleMemoryCheck::preSourceConstruction);
0355         iReg.watchPostSourceConstruction(this, &SimpleMemoryCheck::postSourceConstruction);
0356         iReg.watchPostSourceEvent(this, &SimpleMemoryCheck::postSourceEvent);
0357         iReg.watchPostModuleConstruction(this, &SimpleMemoryCheck::postModuleConstruction);
0358         iReg.watchPostModuleBeginJob(this, &SimpleMemoryCheck::postModuleBeginJob);
0359         iReg.watchPostEvent(this, &SimpleMemoryCheck::postEvent);
0360         iReg.watchPostModuleEvent(this, &SimpleMemoryCheck::postModule);
0361         iReg.watchPostBeginJob(this, &SimpleMemoryCheck::postBeginJob);
0362         iReg.watchPostEndJob(this, &SimpleMemoryCheck::postEndJob);
0363       } else {
0364         iReg.watchPostEvent(this, &SimpleMemoryCheck::postEvent);
0365         iReg.watchPostEndJob(this, &SimpleMemoryCheck::postEndJob);
0366       }
0367       if (moduleSummaryRequested_) {  // changelog 2
0368         iReg.watchPreModuleEvent(this, &SimpleMemoryCheck::preModule);
0369         if (oncePerEventMode_) {
0370           iReg.watchPostModuleEvent(this, &SimpleMemoryCheck::postModule);
0371         }
0372       }
0373 
0374       // The following are not currenty used/implemented below for either
0375       // of the print modes (but are left here for reference)
0376       //  iReg.watchPostBeginJob(this,
0377       //       &SimpleMemoryCheck::postBeginJob);
0378       //  iReg.watchPreProcessEvent(this,
0379       //       &SimpleMemoryCheck::preEventProcessing);
0380       //  iReg.watchPreModule(this,
0381       //       &SimpleMemoryCheck::preModule);
0382     }
0383 
0384     SimpleMemoryCheck::~SimpleMemoryCheck() {
0385 #ifdef LINUX
0386       if (nullptr != smapsFile_) {
0387         fclose(smapsFile_);
0388       }
0389 #endif
0390       if (smapsLineBuffer_) {
0391         //getline will create the memory using malloc
0392         free(smapsLineBuffer_);
0393       }
0394     }
0395 
0396     void SimpleMemoryCheck::fillDescriptions(ConfigurationDescriptions& descriptions) {
0397       ParameterSetDescription desc;
0398       desc.addUntracked<int>("ignoreTotal", 1);
0399       desc.addUntracked<bool>("showMallocInfo", false);
0400       desc.addUntracked<bool>("oncePerEventMode", false);
0401       desc.addUntracked<bool>("jobReportOutputOnly", false);
0402       desc.addUntracked<bool>("monitorPssAndPrivate", false);
0403       desc.addUntracked<bool>("moduleMemorySummary", false);
0404       desc.addUntracked<bool>("dump", false);
0405       descriptions.add("SimpleMemoryCheck", desc);
0406     }
0407 
0408     void SimpleMemoryCheck::openFiles() {
0409 #ifdef LINUX
0410       if (monitorPssAndPrivate_) {
0411         std::ostringstream smapsNameOst;
0412         smapsNameOst << "/proc/" << getpid() << "/smaps";
0413         if ((smapsFile_ = fopen(smapsNameOst.str().c_str(), "r")) == nullptr) {
0414           throw Exception(errors::Configuration) << "Failed to open smaps file " << smapsNameOst.str() << std::endl;
0415         }
0416       }
0417 #endif
0418     }
0419 
0420     void SimpleMemoryCheck::postBeginJob() {
0421       growthRateVsize_ = current_->vsize;
0422       growthRateRss_ = current_->rss;
0423     }
0424 
0425     void SimpleMemoryCheck::preSourceConstruction(ModuleDescription const& md) {
0426       bool expected = false;
0427       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0428         std::shared_ptr<void> guard(
0429             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0430         updateAndPrint("pre-ctor", md.moduleLabel(), md.moduleName());
0431       }
0432     }
0433 
0434     void SimpleMemoryCheck::postSourceConstruction(ModuleDescription const& md) {
0435       bool expected = false;
0436       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0437         std::shared_ptr<void> guard(
0438             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0439         updateAndPrint("ctor", md.moduleLabel(), md.moduleName());
0440       }
0441     }
0442 
0443     void SimpleMemoryCheck::postSourceEvent(StreamID sid) {
0444       bool expected = false;
0445       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0446         std::shared_ptr<void> guard(
0447             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0448         updateAndPrint("module", "source", "source");
0449       }
0450     }
0451 
0452     void SimpleMemoryCheck::postModuleConstruction(ModuleDescription const& md) {
0453       bool expected = false;
0454       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0455         std::shared_ptr<void> guard(
0456             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0457         updateAndPrint("ctor", md.moduleLabel(), md.moduleName());
0458       }
0459     }
0460 
0461     void SimpleMemoryCheck::postModuleBeginJob(ModuleDescription const& md) {
0462       bool expected = false;
0463       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0464         std::shared_ptr<void> guard(
0465             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0466         updateAndPrint("beginJob", md.moduleLabel(), md.moduleName());
0467       }
0468     }
0469 
0470     void SimpleMemoryCheck::postEndJob() {
0471       if (not jobReportOutputOnly_) {
0472         LogAbsolute("MemoryReport")  // changelog 1
0473             << "MemoryReport> Peak virtual size " << eventT1_.vsize << " Mbytes"
0474             << "\n"
0475             << " Key events increasing vsize: \n"
0476             << eventL2_ << "\n"
0477             << eventL1_ << "\n"
0478             << eventM_ << "\n"
0479             << eventR1_ << "\n"
0480             << eventR2_ << "\n"
0481             << eventT3_ << "\n"
0482             << eventT2_ << "\n"
0483             << eventT1_ << "\nMemoryReport> Peak rss size " << eventRssT1_.rss
0484             << " Mbytes"
0485                "\n Key events increasing rss:\n"
0486             << eventRssT3_ << "\n"
0487             << eventRssT2_ << "\n"
0488             << eventRssT1_ << "\n"
0489             << eventDeltaRssT3_ << "\n"
0490             << eventDeltaRssT2_ << "\n"
0491             << eventDeltaRssT1_;
0492         ;
0493       }
0494       if (moduleSummaryRequested_ and not jobReportOutputOnly_) {  // changelog 1
0495         LogAbsolute mmr("ModuleMemoryReport");                     // at end of if block, mmr
0496                                                                    // is destructed, causing
0497                                                                    // message to be logged
0498         mmr << "ModuleMemoryReport> Each line has module label and: \n";
0499         mmr << "  (after early ignored events) \n";
0500         mmr << "    count of times module executed; average increase in vsize \n";
0501         mmr << "    maximum increase in vsize; event on which maximum occurred \n";
0502         mmr << "  (during early ignored events) \n";
0503         mmr << "    total and maximum vsize increases \n \n";
0504         for (SignificantModulesMap::iterator im = modules_.begin(); im != modules_.end(); ++im) {
0505           SignificantModule const& m = im->second;
0506           if (m.totalDeltaVsize == 0 && m.totalEarlyVsize == 0)
0507             continue;
0508           mmr << im->first << ": ";
0509           mmr << "n = " << m.postEarlyCount;
0510           if (m.postEarlyCount > 0) {
0511             mmr << " avg = " << m.totalDeltaVsize / m.postEarlyCount;
0512           }
0513           mmr << " max = " << m.maxDeltaVsize << " " << m.eventMaxDeltaV;
0514           if (m.totalEarlyVsize > 0) {
0515             mmr << " early total: " << m.totalEarlyVsize;
0516             mmr << " max: " << m.maxEarlyVsize;
0517           }
0518           mmr << "\n";
0519         }
0520       }  // end of if; mmr goes out of scope; log message is queued
0521 
0522       Service<JobReport> reportSvc;
0523       // changelog 1
0524 #define SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
0525 #ifdef SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
0526       //     std::map<std::string, double> reportData;
0527       std::map<std::string, std::string> reportData;
0528 
0529       if (eventL2_.vsize > 0)
0530         eventStatOutput("LargeVsizeIncreaseEventL2", eventL2_, reportData);
0531       if (eventL1_.vsize > 0)
0532         eventStatOutput("LargeVsizeIncreaseEventL1", eventL1_, reportData);
0533       if (eventM_.vsize > 0)
0534         eventStatOutput("LargestVsizeIncreaseEvent", eventM_, reportData);
0535       if (eventR1_.vsize > 0)
0536         eventStatOutput("LargeVsizeIncreaseEventR1", eventR1_, reportData);
0537       if (eventR2_.vsize > 0)
0538         eventStatOutput("LargeVsizeIncreaseEventR2", eventR2_, reportData);
0539       if (eventT3_.vsize > 0)
0540         eventStatOutput("ThirdLargestVsizeEventT3", eventT3_, reportData);
0541       if (eventT2_.vsize > 0)
0542         eventStatOutput("SecondLargestVsizeEventT2", eventT2_, reportData);
0543       if (eventT1_.vsize > 0)
0544         eventStatOutput("LargestVsizeEventT1", eventT1_, reportData);
0545 
0546       if (eventRssT3_.rss > 0)
0547         eventStatOutput("ThirdLargestRssEvent", eventRssT3_, reportData);
0548       if (eventRssT2_.rss > 0)
0549         eventStatOutput("SecondLargestRssEvent", eventRssT2_, reportData);
0550       if (eventRssT1_.rss > 0)
0551         eventStatOutput("LargestRssEvent", eventRssT1_, reportData);
0552       if (eventDeltaRssT3_.deltaRss > 0)
0553         eventStatOutput("ThirdLargestIncreaseRssEvent", eventDeltaRssT3_, reportData);
0554       if (eventDeltaRssT2_.deltaRss > 0)
0555         eventStatOutput("SecondLargestIncreaseRssEvent", eventDeltaRssT2_, reportData);
0556       if (eventDeltaRssT1_.deltaRss > 0)
0557         eventStatOutput("LargestIncreaseRssEvent", eventDeltaRssT1_, reportData);
0558 
0559 #ifdef __linux__
0560 #if (__GLIBC__ > 2) || (__GLIBC__ == 2 && __GLIBC_MINOR__ >= 33)
0561       struct mallinfo2 minfo = mallinfo2();
0562 #else
0563       struct mallinfo minfo = mallinfo();
0564 #endif
0565       reportData.insert(std::make_pair("HEAP_ARENA_SIZE_BYTES", std::to_string(minfo.arena)));
0566       reportData.insert(std::make_pair("HEAP_ARENA_N_UNUSED_CHUNKS", std::to_string(minfo.ordblks)));
0567       reportData.insert(std::make_pair("HEAP_TOP_FREE_BYTES", std::to_string(minfo.keepcost)));
0568       reportData.insert(std::make_pair("HEAP_MAPPED_SIZE_BYTES", std::to_string(minfo.hblkhd)));
0569       reportData.insert(std::make_pair("HEAP_MAPPED_N_CHUNKS", std::to_string(minfo.hblks)));
0570       reportData.insert(std::make_pair("HEAP_USED_BYTES", std::to_string(minfo.uordblks)));
0571       reportData.insert(std::make_pair("HEAP_UNUSED_BYTES", std::to_string(minfo.fordblks)));
0572 #endif
0573 
0574       // Report Growth rates for VSize and Rss
0575       reportData.insert(std::make_pair("AverageGrowthRateVsize",
0576                                        d2str(averageGrowthRate(current_->vsize, growthRateVsize_, count_))));
0577       reportData.insert(std::make_pair("PeakValueVsize", d2str(eventT1_.vsize)));
0578       reportData.insert(
0579           std::make_pair("AverageGrowthRateRss", d2str(averageGrowthRate(current_->rss, growthRateRss_, count_))));
0580       reportData.insert(std::make_pair("PeakValueRss", d2str(eventRssT1_.rss)));
0581 
0582       if (moduleSummaryRequested_) {  // changelog 2
0583         for (SignificantModulesMap::iterator im = modules_.begin(); im != modules_.end(); ++im) {
0584           SignificantModule const& m = im->second;
0585           if (m.totalDeltaVsize == 0 && m.totalEarlyVsize == 0)
0586             continue;
0587           std::string label = im->first + ":";
0588           reportData.insert(std::make_pair(label + "PostEarlyCount", i2str(m.postEarlyCount)));
0589           if (m.postEarlyCount > 0) {
0590             reportData.insert(std::make_pair(label + "AverageDeltaVsize", d2str(m.totalDeltaVsize / m.postEarlyCount)));
0591           }
0592           reportData.insert(std::make_pair(label + "MaxDeltaVsize", d2str(m.maxDeltaVsize)));
0593           if (m.totalEarlyVsize > 0) {
0594             reportData.insert(std::make_pair(label + "TotalEarlyVsize", d2str(m.totalEarlyVsize)));
0595             reportData.insert(std::make_pair(label + "MaxEarlyDeltaVsize", d2str(m.maxEarlyVsize)));
0596           }
0597         }
0598       }
0599 
0600       std::map<std::string, std::string> reportMemoryProperties;
0601 
0602       if (FILE* fmeminfo = fopen("/proc/meminfo", "r")) {
0603         char buf[128];
0604         char space[] = " ";
0605         size_t value;
0606         while (fgets(buf, sizeof(buf), fmeminfo)) {
0607           char* saveptr;
0608           char* token = nullptr;
0609           token = strtok_r(buf, space, &saveptr);
0610           if (token != nullptr) {
0611             value = atol(strtok_r(nullptr, space, &saveptr));
0612             std::string category = token;
0613             reportMemoryProperties.insert(std::make_pair(category.substr(0, strlen(token) - 1), i2str(value)));
0614           }
0615         }
0616 
0617         fclose(fmeminfo);
0618       }
0619 
0620       //      reportSvc->reportMemoryInfo(reportData, reportMemoryProperties);
0621       reportSvc->reportPerformanceSummary("ApplicationMemory", reportData);
0622       reportSvc->reportPerformanceSummary("SystemMemory", reportMemoryProperties);
0623 #endif
0624 
0625 #ifdef SIMPLE_MEMORY_CHECK_DIFFERENT_XML_OUTPUT
0626       std::vector<std::string> reportData;
0627 
0628       if (eventL2_.vsize > 0)
0629         reportData.push_back(eventStatOutput("LargeVsizeIncreaseEventL2", eventL2_));
0630       if (eventL1_.vsize > 0)
0631         reportData.push_back(eventStatOutput("LargeVsizeIncreaseEventL1", eventL1_));
0632       if (eventM_.vsize > 0)
0633         reportData.push_back(eventStatOutput("LargestVsizeIncreaseEvent", eventM_));
0634       if (eventR1_.vsize > 0)
0635         reportData.push_back(eventStatOutput("LargeVsizeIncreaseEventR1", eventR1_));
0636       if (eventR2_.vsize > 0)
0637         reportData.push_back(eventStatOutput("LargeVsizeIncreaseEventR2", eventR2_));
0638       if (eventT3_.vsize > 0)
0639         reportData.push_back(eventStatOutput("ThirdLargestVsizeEventT3", eventT3_));
0640       if (eventT2_.vsize > 0)
0641         reportData.push_back(eventStatOutput("SecondLargestVsizeEventT2", eventT2_));
0642       if (eventT1_.vsize > 0)
0643         reportData.push_back(eventStatOutput("LargestVsizeEventT1", eventT1_));
0644 
0645       if (eventRssT3_.rss > 0)
0646         eventStatOutput("ThirdLargestRssEvent", eventRssT3_, reportData);
0647       if (eventRssT2_.rss > 0)
0648         eventStatOutput("SecondLargestRssEvent", eventRssT2_, reportData);
0649       if (eventRssT1_.rss > 0)
0650         eventStatOutput("LargestRssEvent", eventRssT1_, reportData);
0651       if (eventDeltaRssT3_.deltaRss > 0)
0652         eventStatOutput("ThirdLargestIncreaseRssEvent", eventDeltaRssT3_, reportData);
0653       if (eventDeltaRssT2_.deltaRss > 0)
0654         eventStatOutput("SecondLargestIncreaseRssEvent", eventDeltaRssT2_, reportData);
0655       if (eventDeltaRssT1_.deltaRss > 0)
0656         eventStatOutput("LargestIncreaseRssEvent", eventDeltaRssT1_, reportData);
0657 
0658 #if (__GLIBC__ > 2) || (__GLIBC__ == 2 && __GLIBC_MINOR__ >= 33)
0659       struct mallinfo2 minfo = mallinfo2();
0660 #else
0661       struct mallinfo minfo = mallinfo();
0662 #endif
0663       reportData.push_back(mallOutput("HEAP_ARENA_SIZE_BYTES", minfo.arena));
0664       reportData.push_back(mallOutput("HEAP_ARENA_N_UNUSED_CHUNKS", minfo.ordblks));
0665       reportData.push_back(mallOutput("HEAP_TOP_FREE_BYTES", minfo.keepcost));
0666       reportData.push_back(mallOutput("HEAP_MAPPED_SIZE_BYTES", minfo.hblkhd));
0667       reportData.push_back(mallOutput("HEAP_MAPPED_N_CHUNKS", minfo.hblks));
0668       reportData.push_back(mallOutput("HEAP_USED_BYTES", minfo.uordblks));
0669       reportData.push_back(mallOutput("HEAP_UNUSED_BYTES", minfo.fordblks));
0670 
0671       // Report Growth rates for VSize and Rss
0672       reportData.insert(std::make_pair("AverageGrowthRateVsize",
0673                                        d2str(averageGrowthRate(current_->vsize, growthRateVsize_, count_))));
0674       reportData.insert(std::make_pair("PeakValueVsize", d2str(eventT1_.vsize)));
0675       reportData.insert(
0676           std::make_pair("AverageGrowthRateRss", d2str(averageGrowthRate(current_->rss, growthRateRss_, count_))));
0677       reportData.insert(std::make_pair("PeakValueRss", d2str(eventRssT1_.rss)));
0678 
0679       reportSvc->reportMemoryInfo(reportData);
0680       // This is a form of reportMemoryInfo taking s vector, not a map
0681 #endif
0682     }  // postEndJob
0683 
0684     void SimpleMemoryCheck::postEvent(StreamContext const& iContext) {
0685       ++count_;
0686       bool expected = false;
0687       if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0688         std::shared_ptr<void> guard(
0689             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0690         update();
0691         if (monitorPssAndPrivate_) {
0692           currentSmaps_ = fetchSmaps();
0693         }
0694         updateEventStats(iContext.eventID());
0695         if (oncePerEventMode_) {
0696           // should probably use be Run:Event or count_ for the label and name
0697           updateMax();
0698           andPrint("event", "", "");
0699         }
0700       }
0701     }
0702 
0703     void SimpleMemoryCheck::preModule(StreamContext const& iStreamContext, ModuleCallingContext const& iModuleContext) {
0704       bool expectedMeasurementUnderway = false;
0705       if (measurementUnderway_.compare_exchange_strong(expectedMeasurementUnderway, true, std::memory_order_acq_rel)) {
0706         std::shared_ptr<void> guard(
0707             nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0708         bool expectedModuleMeasurementUnderway = false;
0709         if (moduleMeasurementUnderway_.compare_exchange_strong(expectedModuleMeasurementUnderway, true)) {
0710           update();
0711           // changelog 2
0712           moduleEntryVsize_ = current_->vsize;
0713           moduleStreamID_.store(iStreamContext.streamID().value(), std::memory_order_release);
0714           moduleID_.store(iModuleContext.moduleDescription()->id(), std::memory_order_release);
0715         }
0716       }
0717     }
0718 
0719     void SimpleMemoryCheck::postModule(StreamContext const& iStreamContext,
0720                                        ModuleCallingContext const& iModuleContext) {
0721       if (!oncePerEventMode_) {
0722         bool expected = false;
0723         if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0724           std::shared_ptr<void> guard(
0725               nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0726           auto const md = iModuleContext.moduleDescription();
0727           updateAndPrint("module", md->moduleLabel(), md->moduleName());
0728         }
0729       }
0730 
0731       if (moduleSummaryRequested_) {
0732         //is this the module instance we are measuring?
0733         if (moduleMeasurementUnderway_.load(std::memory_order_acquire) and
0734             (iStreamContext.streamID().value() == moduleStreamID_.load(std::memory_order_acquire)) and
0735             (iModuleContext.moduleDescription()->id() == moduleID_.load(std::memory_order_acquire))) {
0736           //Need to release our module measurement lock
0737           std::shared_ptr<void> guardModuleMeasurementUnderway(
0738               nullptr, [this](void const*) { moduleMeasurementUnderway_.store(false, std::memory_order_release); });
0739           bool expected = false;
0740           if (measurementUnderway_.compare_exchange_strong(expected, true, std::memory_order_acq_rel)) {
0741             std::shared_ptr<void> guardMeasurementUnderway(
0742                 nullptr, [this](void const*) { measurementUnderway_.store(false, std::memory_order_release); });
0743             if (oncePerEventMode_) {
0744               update();
0745             }
0746             // changelog 2
0747             double dv = current_->vsize - moduleEntryVsize_;
0748             std::string label = iModuleContext.moduleDescription()->moduleLabel();
0749             updateModuleMemoryStats(modules_[label], dv, iStreamContext.eventID());
0750           }
0751         }
0752       }
0753     }
0754 
0755     void SimpleMemoryCheck::update() {
0756       std::swap(current_, previous_);
0757       *current_ = fetch();
0758     }
0759 
0760     void SimpleMemoryCheck::updateMax() {
0761       auto v = *current_;
0762       if ((v > max_) || oncePerEventMode_) {
0763         if (max_.vsize < v.vsize) {
0764           max_.vsize = v.vsize;
0765         }
0766         if (max_.rss < v.rss) {
0767           max_.rss = v.rss;
0768         }
0769       }
0770     }
0771 
0772     void SimpleMemoryCheck::updateEventStats(EventID const& e) {
0773       if (count_ < num_to_skip_)
0774         return;
0775       if (count_ == num_to_skip_) {
0776         eventT1_.set(0, 0, e, this);
0777         eventM_.set(0, 0, e, this);
0778         eventRssT1_.set(0, 0, e, this);
0779         eventDeltaRssT1_.set(0, 0, e, this);
0780         return;
0781       }
0782       double vsize = current_->vsize;
0783       double deltaVsize = vsize - eventT1_.vsize;
0784 
0785       // Update significative events for Vsize
0786       if (vsize > eventT1_.vsize) {
0787         double deltaRss = current_->rss - eventT1_.rss;
0788         eventT3_ = eventT2_;
0789         eventT2_ = eventT1_;
0790         eventT1_.set(deltaVsize, deltaRss, e, this);
0791       } else if (vsize > eventT2_.vsize) {
0792         double deltaRss = current_->rss - eventT1_.rss;
0793         eventT3_ = eventT2_;
0794         eventT2_.set(deltaVsize, deltaRss, e, this);
0795       } else if (vsize > eventT3_.vsize) {
0796         double deltaRss = current_->rss - eventT1_.rss;
0797         eventT3_.set(deltaVsize, deltaRss, e, this);
0798       }
0799 
0800       if (deltaVsize > eventM_.deltaVsize) {
0801         double deltaRss = current_->rss - eventM_.rss;
0802         if (eventL1_.deltaVsize >= eventR1_.deltaVsize) {
0803           eventL2_ = eventL1_;
0804         } else {
0805           eventL2_ = eventR1_;
0806         }
0807         eventL1_ = eventM_;
0808         eventM_.set(deltaVsize, deltaRss, e, this);
0809         eventR1_ = SignificantEvent();
0810         eventR2_ = SignificantEvent();
0811       } else if (deltaVsize > eventR1_.deltaVsize) {
0812         double deltaRss = current_->rss - eventM_.rss;
0813         eventR2_ = eventR1_;
0814         eventR1_.set(deltaVsize, deltaRss, e, this);
0815       } else if (deltaVsize > eventR2_.deltaVsize) {
0816         double deltaRss = current_->rss - eventR1_.rss;
0817         eventR2_.set(deltaVsize, deltaRss, e, this);
0818       }
0819 
0820       // Update significative events for Rss
0821       double rss = current_->rss;
0822       double deltaRss = rss - eventRssT1_.rss;
0823 
0824       if (rss > eventRssT1_.rss) {
0825         eventRssT3_ = eventRssT2_;
0826         eventRssT2_ = eventRssT1_;
0827         eventRssT1_.set(deltaVsize, deltaRss, e, this);
0828       } else if (rss > eventRssT2_.rss) {
0829         eventRssT3_ = eventRssT2_;
0830         eventRssT2_.set(deltaVsize, deltaRss, e, this);
0831       } else if (rss > eventRssT3_.rss) {
0832         eventRssT3_.set(deltaVsize, deltaRss, e, this);
0833       }
0834       if (deltaRss > eventDeltaRssT1_.deltaRss) {
0835         eventDeltaRssT3_ = eventDeltaRssT2_;
0836         eventDeltaRssT2_ = eventDeltaRssT1_;
0837         eventDeltaRssT1_.set(deltaVsize, deltaRss, e, this);
0838       } else if (deltaRss > eventDeltaRssT2_.deltaRss) {
0839         eventDeltaRssT3_ = eventDeltaRssT2_;
0840         eventDeltaRssT2_.set(deltaVsize, deltaRss, e, this);
0841       } else if (deltaRss > eventDeltaRssT3_.deltaRss) {
0842         eventDeltaRssT3_.set(deltaVsize, deltaRss, e, this);
0843       }
0844     }  // updateEventStats
0845 
0846     void SimpleMemoryCheck::andPrint(std::string const& type,
0847                                      std::string const& mdlabel,
0848                                      std::string const& mdname) const {
0849       if (not jobReportOutputOnly_ && ((*current_ > max_) || oncePerEventMode_)) {
0850         if (count_ >= num_to_skip_) {
0851           double deltaVSIZE = current_->vsize - max_.vsize;
0852           double deltaRSS = current_->rss - max_.rss;
0853           if (!showMallocInfo_) {  // default
0854             LogWarning("MemoryCheck") << "MemoryCheck: " << type << " " << mdname << ":" << mdlabel << " VSIZE "
0855                                       << current_->vsize << " " << deltaVSIZE << " RSS " << current_->rss << " "
0856                                       << deltaRSS;
0857           } else {
0858 #ifdef __linux__
0859 #if (__GLIBC__ > 2) || (__GLIBC__ == 2 && __GLIBC_MINOR__ >= 33)
0860             struct mallinfo2 minfo = mallinfo2();
0861 #else
0862             struct mallinfo minfo = mallinfo();
0863 #endif
0864 #endif
0865             LogWarning("MemoryCheck") << "MemoryCheck: " << type << " " << mdname << ":" << mdlabel << " VSIZE "
0866                                       << current_->vsize << " " << deltaVSIZE << " RSS " << current_->rss << " "
0867                                       << deltaRSS
0868 #ifdef __linux__
0869                                       << " HEAP-ARENA [ SIZE-BYTES " << minfo.arena << " N-UNUSED-CHUNKS "
0870                                       << minfo.ordblks << " TOP-FREE-BYTES " << minfo.keepcost << " ]"
0871                                       << " HEAP-MAPPED [ SIZE-BYTES " << minfo.hblkhd << " N-CHUNKS " << minfo.hblks
0872                                       << " ]"
0873                                       << " HEAP-USED-BYTES " << minfo.uordblks << " HEAP-UNUSED-BYTES "
0874                                       << minfo.fordblks
0875 #endif
0876                 ;
0877           }
0878         }
0879       }
0880     }
0881 
0882     void SimpleMemoryCheck::updateAndPrint(std::string const& type,
0883                                            std::string const& mdlabel,
0884                                            std::string const& mdname) {
0885       update();
0886       andPrint(type, mdlabel, mdname);
0887       updateMax();
0888     }
0889 
0890 #ifdef SIMPLE_MEMORY_CHECK_ORIGINAL_XML_OUTPUT
0891     void SimpleMemoryCheck::eventStatOutput(std::string title,
0892                                             SignificantEvent const& e,
0893                                             std::map<std::string, std::string>& m) const {
0894       {
0895         std::ostringstream os;
0896         os << title << "-a-COUNT";
0897         m.insert(std::make_pair(os.str(), i2str(e.count)));
0898       }
0899       {
0900         std::ostringstream os;
0901         os << title << "-b-RUN";
0902         m.insert(std::make_pair(os.str(), d2str(static_cast<double>(e.event.run()))));
0903       }
0904       {
0905         std::ostringstream os;
0906         os << title << "-c-EVENT";
0907         m.insert(std::make_pair(os.str(), d2str(static_cast<double>(e.event.event()))));
0908       }
0909       {
0910         std::ostringstream os;
0911         os << title << "-d-VSIZE";
0912         m.insert(std::make_pair(os.str(), d2str(e.vsize)));
0913       }
0914       {
0915         std::ostringstream os;
0916         os << title << "-e-DELTV";
0917         m.insert(std::make_pair(os.str(), d2str(e.deltaVsize)));
0918       }
0919       {
0920         std::ostringstream os;
0921         os << title << "-f-RSS";
0922         m.insert(std::make_pair(os.str(), d2str(e.rss)));
0923       }
0924       if (monitorPssAndPrivate_) {
0925         {
0926           std::ostringstream os;
0927           os << title << "-g-PRIVATE";
0928           m.insert(std::make_pair(os.str(), d2str(e.privateSize)));
0929         }
0930         {
0931           std::ostringstream os;
0932           os << title << "-h-PSS";
0933           m.insert(std::make_pair(os.str(), d2str(e.pss)));
0934         }
0935       }
0936     }  // eventStatOutput
0937 #endif
0938 
0939 #ifdef SIMPLE_MEMORY_CHECK_DIFFERENT_XML_OUTPUT
0940     std::string SimpleMemoryCheck::eventStatOutput(std::string title, SignificantEvent const& e) const {
0941       std::ostringstream os;
0942       os << "  <" << title << ">\n";
0943       os << "    " << e.count << ": " << e.event;
0944       os << " vsize " << e.vsize - e.deltaVsize << " + " << e.deltaVsize << " = " << e.vsize;
0945       os << "  rss: " << e.rss << "\n";
0946       os << "  </" << title << ">\n";
0947       return os.str();
0948     }  // eventStatOutput
0949 
0950     std::string SimpleMemoryCheck::mallOutput(std::string title, size_t const& n) const {
0951       std::ostringstream os;
0952       os << "  <" << title << ">\n";
0953       os << "    " << n << "\n";
0954       os << "  </" << title << ">\n";
0955       return os.str();
0956     }
0957 #endif
0958     // changelog 2
0959     void SimpleMemoryCheck::updateModuleMemoryStats(SignificantModule& m,
0960                                                     double dv,
0961                                                     edm::EventID const& currentEventID) {
0962       if (count_ < num_to_skip_) {
0963         m.totalEarlyVsize += dv;
0964         if (dv > m.maxEarlyVsize)
0965           m.maxEarlyVsize = dv;
0966       } else {
0967         ++m.postEarlyCount;
0968         m.totalDeltaVsize += dv;
0969         if (dv > m.maxDeltaVsize) {
0970           m.maxDeltaVsize = dv;
0971           m.eventMaxDeltaV = currentEventID;
0972         }
0973       }
0974     }  //updateModuleMemoryStats
0975 
0976     std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantEvent const& se) {
0977       os << "[" << se.count << "] " << se.event << "  vsize = " << se.vsize << " deltaVsize = " << se.deltaVsize
0978          << " rss = " << se.rss << " delta = " << se.deltaRss;
0979 
0980       if (se.monitorPssAndPrivate) {
0981         os << " private = " << se.privateSize << " pss = " << se.pss;
0982       }
0983       return os;
0984     }
0985 
0986     std::ostream& operator<<(std::ostream& os, SimpleMemoryCheck::SignificantModule const& sm) {
0987       if (sm.postEarlyCount > 0) {
0988         os << "\nPost Early Events:  TotalDeltaVsize: " << sm.totalDeltaVsize
0989            << " (avg: " << sm.totalDeltaVsize / sm.postEarlyCount << "; max: " << sm.maxDeltaVsize << " during "
0990            << sm.eventMaxDeltaV << ")";
0991       }
0992       if (sm.totalEarlyVsize > 0) {
0993         os << "\n     Early Events:  TotalDeltaVsize: " << sm.totalEarlyVsize << " (max: " << sm.maxEarlyVsize << ")";
0994       }
0995 
0996       return os;
0997     }
0998 
0999   }  // end namespace service
1000 }  // end namespace edm
1001 
1002 #if defined(__linux__)
1003 using edm::service::SimpleMemoryCheck;
1004 DEFINE_FWK_SERVICE(SimpleMemoryCheck);
1005 #endif