File indexing completed on 2024-04-06 12:09:55
0001 #include "DQMOffline/Trigger/interface/HLTTauDQMOfflineSource.h"
0002
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004
0005 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0006
0007 using namespace std;
0008 using namespace edm;
0009 using namespace reco;
0010 using namespace trigger;
0011
0012
0013
0014
0015 HLTTauDQMOfflineSource::HLTTauDQMOfflineSource(const edm::ParameterSet& ps)
0016 : hltProcessName_(ps.getUntrackedParameter<std::string>("HLTProcessName", "HLT")),
0017 triggerResultsSrc_(ps.getUntrackedParameter<edm::InputTag>("TriggerResultsSrc")),
0018 triggerResultsToken_(consumes<edm::TriggerResults>(triggerResultsSrc_)),
0019 triggerEventSrc_(ps.getUntrackedParameter<edm::InputTag>("TriggerEventSrc")),
0020 triggerEventToken_(consumes<trigger::TriggerEvent>(triggerEventSrc_)),
0021 pathRegex_(ps.getUntrackedParameter<std::string>("Paths")),
0022 nPtBins_(ps.getUntrackedParameter<int>("PtHistoBins", 20)),
0023 nEtaBins_(ps.getUntrackedParameter<int>("EtaHistoBins", 12)),
0024 nPhiBins_(ps.getUntrackedParameter<int>("PhiHistoBins", 18)),
0025 ptMax_(ps.getUntrackedParameter<double>("PtHistoMax", 200)),
0026 highPtMax_(ps.getUntrackedParameter<double>("HighPtHistoMax", 1000)),
0027 l1MatchDr_(ps.getUntrackedParameter<double>("L1MatchDeltaR", 0.5)),
0028 hltMatchDr_(ps.getUntrackedParameter<double>("HLTMatchDeltaR", 0.5)),
0029 dqmBaseFolder_(ps.getUntrackedParameter<std::string>("DQMBaseFolder")),
0030 counterEvt_(0),
0031 prescaleEvt_(ps.getUntrackedParameter<int>("prescaleEvt", -1)) {
0032 edm::ParameterSet matching = ps.getParameter<edm::ParameterSet>("Matching");
0033 doRefAnalysis_ = matching.getUntrackedParameter<bool>("doMatching");
0034
0035 iWrapper = new HistoWrapper(ps);
0036
0037 if (ps.exists("L1Plotter") && !ps.exists("TagAndProbe")) {
0038 l1Plotter_ = std::make_unique<HLTTauDQML1Plotter>(ps.getUntrackedParameter<edm::ParameterSet>("L1Plotter"),
0039 consumesCollector(),
0040 nPhiBins_,
0041 ptMax_,
0042 highPtMax_,
0043 doRefAnalysis_,
0044 l1MatchDr_,
0045 dqmBaseFolder_);
0046 }
0047 if (ps.exists("PathSummaryPlotter")) {
0048 pathSummaryPlotter_ = std::make_unique<HLTTauDQMPathSummaryPlotter>(
0049 ps.getUntrackedParameter<edm::ParameterSet>("PathSummaryPlotter"), doRefAnalysis_, dqmBaseFolder_, hltMatchDr_);
0050 }
0051 tagAndProbe_ = false;
0052 if (ps.exists("TagAndProbe")) {
0053 tagAndProbePaths = ps.getUntrackedParameter<std::vector<edm::ParameterSet> >("TagAndProbe");
0054 tagAndProbe_ = true;
0055 }
0056
0057 if (doRefAnalysis_) {
0058 using VPSet = std::vector<edm::ParameterSet>;
0059 VPSet matchObjects = matching.getUntrackedParameter<VPSet>("matchFilters");
0060 for (const edm::ParameterSet& pset : matchObjects) {
0061 refObjects_.push_back(RefObject{pset.getUntrackedParameter<int>("matchObjectID"),
0062 consumes<LVColl>(pset.getUntrackedParameter<edm::InputTag>("FilterName"))});
0063 }
0064 }
0065 }
0066
0067 HLTTauDQMOfflineSource::~HLTTauDQMOfflineSource() = default;
0068
0069
0070 void HLTTauDQMOfflineSource::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0071
0072 bool hltMenuChanged = false;
0073 if (HLTCP_.init(iRun, iSetup, hltProcessName_, hltMenuChanged)) {
0074 LogDebug("HLTTauDQMOffline") << "dqmBeginRun(), hltMenuChanged " << hltMenuChanged;
0075 if (hltMenuChanged) {
0076 if (!tagAndProbe_) {
0077
0078 std::vector<std::string> foundPaths;
0079 std::smatch what;
0080 LogDebug("HLTTauDQMOffline") << "Looking for paths with regex " << pathRegex_;
0081 const std::regex pathRegex(pathRegex_);
0082 for (const std::string& pathName : HLTCP_.triggerNames()) {
0083 if (std::regex_search(pathName, what, pathRegex)) {
0084 LogDebug("HLTTauDQMOffline") << "Found path " << pathName;
0085 foundPaths.emplace_back(pathName);
0086 }
0087 }
0088 std::sort(foundPaths.begin(), foundPaths.end());
0089
0090
0091 std::vector<const HLTTauDQMPath*> pathObjects;
0092 pathPlotters_.reserve(foundPaths.size());
0093 pathObjects.reserve(foundPaths.size());
0094 for (const std::string& pathName : foundPaths) {
0095 pathPlotters_.emplace_back(pathName,
0096 HLTCP_,
0097 doRefAnalysis_,
0098 dqmBaseFolder_,
0099 hltProcessName_,
0100 nPtBins_,
0101 nEtaBins_,
0102 nPhiBins_,
0103 ptMax_,
0104 highPtMax_,
0105 l1MatchDr_,
0106 hltMatchDr_);
0107 if (pathPlotters_.back().isValid()) {
0108 pathObjects.push_back(pathPlotters_.back().getPathObject());
0109 }
0110 }
0111
0112
0113 if (pathSummaryPlotter_) {
0114 pathSummaryPlotter_->setPathObjects(pathObjects);
0115 }
0116 } else {
0117
0118 std::vector<std::string> foundPaths;
0119 std::smatch what;
0120
0121 for (const edm::ParameterSet& tpset : tagAndProbePaths) {
0122 std::vector<std::string> moduleLabels;
0123 edm::ParameterSet denpset = tpset.getParameter<edm::ParameterSet>("denominator");
0124 std::vector<std::string> denominators = denpset.getParameter<std::vector<std::string> >("hltPaths");
0125 std::vector<std::string> updatedDenominators;
0126 for (size_t i = 0; i < denominators.size(); ++i) {
0127 const std::regex denRegex_(denominators[i]);
0128 for (const std::string& pathName : HLTCP_.triggerNames()) {
0129 if (std::regex_search(pathName, what, denRegex_)) {
0130 updatedDenominators.push_back(pathName);
0131 moduleLabels = HLTCP_.moduleLabels(pathName);
0132 }
0133 }
0134 }
0135 denpset.addParameter<std::vector<std::string> >("hltPaths", updatedDenominators);
0136
0137 edm::ParameterSet numpset = tpset.getParameter<edm::ParameterSet>("numerator");
0138 std::vector<std::string> numerators = numpset.getParameter<std::vector<std::string> >("hltPaths");
0139
0140 const std::regex numRegex_(numerators[0]);
0141 for (const std::string& pathName : HLTCP_.triggerNames()) {
0142 if (std::regex_search(pathName, what, numRegex_)) {
0143 edm::ParameterSet new_tpset = tpset;
0144 new_tpset.addParameter<std::string>("name", pathName);
0145 std::vector<std::string> updatedHltPaths;
0146 updatedHltPaths.push_back(pathName);
0147 numpset.addParameter<std::vector<std::string> >("hltPaths", updatedHltPaths);
0148 new_tpset.addParameter<edm::ParameterSet>("numerator", numpset);
0149 new_tpset.addParameter<edm::ParameterSet>("denominator", denpset);
0150
0151 tagandprobePlotters_.emplace_back(
0152 new HLTTauDQMTagAndProbePlotter(new_tpset, moduleLabels, dqmBaseFolder_));
0153 }
0154 }
0155 }
0156 }
0157 }
0158 } else {
0159 edm::LogWarning("HLTTauDQMOffline") << "HLT config extraction failure with process name '" << hltProcessName_
0160 << "'";
0161 }
0162 }
0163
0164
0165 void HLTTauDQMOfflineSource::bookHistograms(DQMStore::IBooker& iBooker,
0166 const edm::Run& iRun,
0167 const EventSetup& iSetup) {
0168 if (l1Plotter_) {
0169 l1Plotter_->bookHistograms(*iWrapper, iBooker);
0170 }
0171 for (auto& pathPlotter : pathPlotters_) {
0172 pathPlotter.bookHistograms(*iWrapper, iBooker);
0173 }
0174 for (auto& tpPlotter : tagandprobePlotters_) {
0175 tpPlotter->bookHistograms(*iWrapper, iBooker, iRun, iSetup);
0176 }
0177 if (pathSummaryPlotter_) {
0178 pathSummaryPlotter_->bookHistograms(*iWrapper, iBooker);
0179 }
0180 }
0181
0182
0183 void HLTTauDQMOfflineSource::analyze(const Event& iEvent, const EventSetup& iSetup) {
0184
0185 if (counterEvt_ > prescaleEvt_) {
0186
0187 counterEvt_ = 0;
0188
0189 edm::Handle<edm::TriggerResults> triggerResultsHandle;
0190 iEvent.getByToken(triggerResultsToken_, triggerResultsHandle);
0191 if (!triggerResultsHandle.isValid()) {
0192 edm::LogWarning("HLTTauDQMOffline") << "Unable to read edm::TriggerResults with label " << triggerResultsSrc_;
0193 return;
0194 }
0195
0196 edm::Handle<trigger::TriggerEvent> triggerEventHandle;
0197 iEvent.getByToken(triggerEventToken_, triggerEventHandle);
0198 if (!triggerEventHandle.isValid()) {
0199 edm::LogWarning("HLTTauDQMOffline") << "Unable to read trigger::TriggerEvent with label " << triggerEventSrc_;
0200 return;
0201 }
0202
0203
0204 HLTTauDQMOfflineObjects refC;
0205 if (doRefAnalysis_) {
0206 for (RefObject& refObj : refObjects_) {
0207 edm::Handle<LVColl> collHandle;
0208 iEvent.getByToken(refObj.token, collHandle);
0209 if (!collHandle.isValid())
0210 continue;
0211
0212 if (refObj.objID == 11) {
0213 refC.electrons.insert(refC.electrons.end(), collHandle->begin(), collHandle->end());
0214 } else if (refObj.objID == 13) {
0215 refC.muons.insert(refC.muons.end(), collHandle->begin(), collHandle->end());
0216 } else if (refObj.objID == 15) {
0217 refC.taus.insert(refC.taus.end(), collHandle->begin(), collHandle->end());
0218 } else if (refObj.objID == 0) {
0219 refC.met.insert(refC.met.end(), collHandle->begin(), collHandle->end());
0220 }
0221 }
0222 }
0223
0224
0225 for (auto& pathPlotter : pathPlotters_) {
0226 if (pathPlotter.isValid())
0227 pathPlotter.analyze(*triggerResultsHandle, *triggerEventHandle, refC);
0228 }
0229
0230 if (pathSummaryPlotter_ && pathSummaryPlotter_->isValid()) {
0231 pathSummaryPlotter_->analyze(*triggerResultsHandle, *triggerEventHandle, refC);
0232 }
0233
0234
0235 if (l1Plotter_ && l1Plotter_->isValid()) {
0236 l1Plotter_->analyze(iEvent, iSetup, refC);
0237 }
0238
0239
0240 for (auto& tpPlotter : tagandprobePlotters_) {
0241 if (tpPlotter->isValid())
0242 tpPlotter->analyze(iEvent, *triggerResultsHandle, *triggerEventHandle, refC);
0243 }
0244
0245 } else {
0246 counterEvt_++;
0247 }
0248 }
0249
0250 #include "FWCore/Framework/interface/MakerMacros.h"
0251 DEFINE_FWK_MODULE(HLTTauDQMOfflineSource);