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