#include "FWCore/Utilities/interface/InputTag.h"0029 #include "DataFormats/Provenance/interface/EventID.h"
0030 #include "FWCore/MessageLogger/interface/ErrorSummaryEntry.h"
0031 #include "FWCore/MessageLogger/interface/LoggedErrorsSummary.h"
0032
0033
0034
0035
0036
0037 class HLTLogMonitorFilter : public edm::EDFilter {
0038 public:
0039 explicit HLTLogMonitorFilter(const edm::ParameterSet &);
0040 ~HLTLogMonitorFilter() override;
0041
0042 struct CategoryEntry {
0043 uint32_t threshold;
0044 uint32_t prescale;
0045 uint64_t counter;
0046 uint64_t accepted;
0047 bool done;
0048
0049 CategoryEntry(uint32_t t = 0) :
0050 threshold(t),
0051 prescale(1),
0052 counter(0),
0053 accepted(0),
0054 done(false)
0055 { }
0056
0057 bool accept() {
0058 if (not done) {
0059 done = true;
0060 ++counter;
0061 }
0062
0063
0064 if ((threshold == 0) or (counter % prescale))
0065 return false;
0066
0067
0068 if (counter == prescale * threshold)
0069 prescale *= threshold;
0070
0071 ++accepted;
0072 return true;
0073 }
0074
0075 };
0076
0077 typedef std::map<std::string, CategoryEntry> CategoryMap;
0078
0079
0080
0081
0082 bool filter(edm::Event&, const edm::EventSetup &) override;
0083
0084
0085 void beginJob() override;
0086
0087
0088 void endJob() override;
0089
0090
0091 bool knownCategory(const std::string & category);
0092
0093
0094 CategoryEntry & addCategory(const std::string & category, uint32_t threshold);
0095
0096
0097 CategoryEntry & getCategory(const std::string & category);
0098
0099
0100 void summary();
0101
0102
0103
0104 uint32_t m_prescale;
0105 CategoryMap m_data;
0106 };
0107
0108
0109
0110 #include <sstream>
0111 #include <iomanip>
0112 #include <memory>
0113 #include <boost/range.hpp>
0114 #include <boost/algorithm/string.hpp>
0115
0116
0117 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0118 #include "FWCore/MessageLogger/interface/MessageSender.h"
0119
0120
0121
0122
0123 HLTLogMonitorFilter::HLTLogMonitorFilter(const edm::ParameterSet & config) :
0124 m_prescale(),
0125 m_data()
0126 {
0127 m_prescale = config.getParameter<uint32_t>("default_threshold");
0128
0129 typedef std::vector<edm::ParameterSet> VPSet;
0130 const VPSet & categories = config.getParameter<VPSet>("categories");
0131 for (auto const & categorie : categories) {
0132 const std::string & name = categorie.getParameter<std::string>("name");
0133 uint32_t threshold = categorie.getParameter<uint32_t>("threshold");
0134 addCategory(name, threshold);
0135 }
0136
0137 produces<std::vector<edm::ErrorSummaryEntry> >();
0138 }
0139
0140 HLTLogMonitorFilter::~HLTLogMonitorFilter() = default;
0141
0142
0143
0144
0145
0146
0147 bool HLTLogMonitorFilter::filter(edm::Event & event, const edm::EventSetup & setup) {
0148
0149 if (not edm::FreshErrorsExist(event.streamID().value()))
0150 return false;
0151
0152
0153 for(auto& entry : m_data)
0154 entry.second.done = false;
0155
0156
0157 bool accept = false;
0158 std::string category;
0159
0160 std::vector<edm::ErrorSummaryEntry> errorSummary{ edm::LoggedErrorsSummary(event.streamID().value()) };
0161 for( auto const& entry : errorSummary ) {
0162
0163 typedef boost::split_iterator<std::string::const_iterator> splitter;
0164 for (splitter i = boost::make_split_iterator(entry.category, boost::first_finder("|", boost::is_equal()));
0165 i != splitter();
0166 ++i)
0167 {
0168
0169
0170 category.assign(i->begin(), i->end());
0171
0172
0173 if (getCategory(category).accept())
0174 accept = true;
0175 }
0176 }
0177
0178
0179 std::unique_ptr<std::vector<edm::ErrorSummaryEntry> > errors(new std::vector<edm::ErrorSummaryEntry>());
0180 if (accept) {
0181 errors->swap(errorSummary);
0182 }
0183 event.put(std::move(errors));
0184
0185 return accept;
0186 }
0187
0188
0189 void HLTLogMonitorFilter::beginJob() {
0190 edm::EnableLoggedErrorsSummary();
0191 }
0192
0193 void HLTLogMonitorFilter::endJob() {
0194 edm::DisableLoggedErrorsSummary();
0195 summary();
0196 }
0197
0198
0199 bool HLTLogMonitorFilter::knownCategory(const std::string & category) {
0200 return (m_data.find( category ) != m_data.end());
0201 }
0202
0203
0204 HLTLogMonitorFilter::CategoryEntry & HLTLogMonitorFilter::addCategory(const std::string & category, uint32_t threshold) {
0205
0206 std::pair<CategoryMap::iterator, bool> result = m_data.insert( std::make_pair(category, CategoryEntry(threshold)) );
0207 if (not result.second)
0208 throw cms::Exception("Configuration") << "Duplicate entry for category " << category;
0209 return result.first->second;
0210 }
0211
0212
0213 HLTLogMonitorFilter::CategoryEntry & HLTLogMonitorFilter::getCategory(const std::string & category) {
0214
0215 auto i = m_data.find(category);
0216 if (i != m_data.end())
0217 return i->second;
0218 else
0219 return m_data.insert( std::make_pair(category, CategoryEntry(m_prescale)) ).first->second;
0220 }
0221
0222
0223 void HLTLogMonitorFilter::summary() {
0224 std::stringstream out;
0225 out << "Log-Report ---------- HLTLogMonitorFilter Summary ------------\n"
0226 << "Log-Report Threshold Prescale Issued Accepted Rejected Category\n";
0227 for(auto const& entry : m_data) {
0228 out << "Log-Report "
0229 << std::right
0230 << std::setw(10) << entry.second.threshold << ' '
0231 << std::setw(10) << entry.second.prescale << ' '
0232 << std::setw(10) << entry.second.counter << ' '
0233 << std::setw(10) << entry.second.accepted << ' '
0234 << std::setw(10) << (entry.second.counter - entry.second.accepted) << ' '
0235 << std::left << entry.first << '\n';
0236 }
0237 edm::LogVerbatim("Report") << out.str();
0238 }
0239
0240
0241 #include "FWCore/Framework/interface/MakerMacros.h"
0242 DEFINE_FWK_MODULE(HLTLogMonitorFilter);