File indexing completed on 2025-02-13 02:58:34
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012 #include <memory>
0013 #include <chrono>
0014 #include <ctime>
0015
0016
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/Utilities/interface/InputTag.h"
0022
0023
0024 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0025
0026
0027 #include "DataFormats/METReco/interface/GenMETCollection.h"
0028 #include "DataFormats/METReco/interface/GenMET.h"
0029
0030
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0033
0034
0035 #include "DQMOffline/Trigger/interface/FunctionDefs.h"
0036
0037
0038 #include "Validation/HLTrigger/interface/HLTGenResHistColl.h"
0039
0040
0041 #include "DQMServices/Core/interface/DQMStore.h"
0042 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0043
0044
0045 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0046 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0047
0048
0049 #include "Validation/HLTrigger/interface/HLTGenValObject.h"
0050 #include "Validation/HLTrigger/interface/HLTGenValObjectMgr.h"
0051
0052 class HLTGenResSource : public DQMEDAnalyzer {
0053 public:
0054 explicit HLTGenResSource(const edm::ParameterSet&);
0055 ~HLTGenResSource() override = default;
0056 HLTGenResSource(const HLTGenResSource&) = delete;
0057 HLTGenResSource& operator=(const HLTGenResSource&) = delete;
0058
0059 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0060
0061 private:
0062 void analyze(const edm::Event&, const edm::EventSetup&) override;
0063 void bookHistograms(DQMStore::IBooker&, edm::Run const& run, edm::EventSetup const& c) override;
0064 void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0065 void initCfgs(const edm::Run&, const edm::EventSetup&);
0066
0067
0068 std::vector<HLTGenValObject> getObjectCollection(const edm::Event&);
0069 std::vector<HLTGenValObject> getGenParticles(const edm::Event&);
0070 reco::GenParticle getLastCopyPreFSR(reco::GenParticle part);
0071 reco::GenParticle getLastCopy(reco::GenParticle part);
0072 bool passGenJetID(const reco::GenJet& jet);
0073
0074
0075
0076 HLTGenValObjectMgr genObjMgr_;
0077 const edm::EDGetTokenT<trigger::TriggerEvent> trigEventToken_;
0078
0079 bool initalised_;
0080 bool booked_;
0081
0082 std::string objType_;
0083 std::string dirName_;
0084 std::vector<edm::ParameterSet> histConfigs_;
0085 std::vector<edm::ParameterSet> histConfigs2D_;
0086 std::vector<edm::ParameterSet> binnings_;
0087 std::string hltProcessName_;
0088 HLTConfigProvider hltConfig_;
0089
0090
0091
0092 std::string infoString_ = "{";
0093
0094
0095 std::vector<HLTGenResHistColl> collections_;
0096 };
0097
0098 HLTGenResSource::HLTGenResSource(const edm::ParameterSet& iConfig)
0099 : genObjMgr_(iConfig.getParameter<edm::ParameterSet>("genConfig"), consumesCollector()),
0100
0101 trigEventToken_(consumes<trigger::TriggerEvent>(iConfig.getParameter<edm::InputTag>("trigEvent"))),
0102 initalised_(false),
0103 booked_(false)
0104
0105 {
0106
0107 dirName_ = iConfig.getParameter<std::string>("dqmDirName");
0108 hltProcessName_ = iConfig.getParameter<std::string>("hltProcessName");
0109
0110 auto resCollections = iConfig.getParameter<std::vector<edm::ParameterSet>>("resCollConfigs");
0111 for (const auto& resCollConfig : resCollections) {
0112 collections_.emplace_back(HLTGenResHistColl(resCollConfig, hltProcessName_));
0113 }
0114 }
0115
0116 void HLTGenResSource::initCfgs(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0117
0118 auto t = std::time(nullptr);
0119 auto tm = *std::localtime(&t);
0120
0121
0122 std::ostringstream timeStringStream;
0123 timeStringStream << std::put_time(&tm, "%d-%m-%Y %H-%M-%S");
0124 auto timeString = timeStringStream.str();
0125 infoString_ += "\"date & time\":\"" + timeString + "\",";
0126
0127
0128 std::string cmsswVersion = std::getenv("CMSSW_VERSION");
0129 infoString_ += std::string("\"CMSSW release\":\"") + cmsswVersion + "\",";
0130
0131
0132 bool changedConfig;
0133 if (!hltConfig_.init(iRun, iSetup, hltProcessName_, changedConfig)) {
0134 edm::LogError("HLTGenResSource") << "Initialization of HLTConfigProvider failed!";
0135 return;
0136 }
0137
0138
0139 infoString_ += std::string("\"global tag\":\"") + hltConfig_.globalTag() + "\",";
0140
0141
0142 infoString_ += std::string("\"HLT ConfDB table\":\"") + hltConfig_.tableName() + "\"}";
0143
0144 initalised_ = true;
0145 }
0146
0147 void HLTGenResSource::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0148 if (!initalised_)
0149 initCfgs(iRun, iSetup);
0150 }
0151
0152
0153 void HLTGenResSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0154
0155
0156
0157 edm::Handle<trigger::TriggerEvent> triggerEvent;
0158 iEvent.getByToken(trigEventToken_, triggerEvent);
0159
0160
0161 for (auto& collection : collections_) {
0162 const std::vector<HLTGenValObject> objects = genObjMgr_.getGenValObjects(iEvent, collection.objType());
0163 for (const auto& object : objects) {
0164 collection.fillHists(object, triggerEvent);
0165 }
0166 }
0167 }
0168
0169
0170 void HLTGenResSource::bookHistograms(DQMStore::IBooker& iBooker, const edm::Run& run, const edm::EventSetup& setup) {
0171 if (booked_)
0172 return;
0173 iBooker.setCurrentFolder(dirName_);
0174
0175 iBooker.bookString("HLTGenResInfo", infoString_);
0176
0177 for (long unsigned int collNr = 0; collNr < collections_.size(); collNr++) {
0178 collections_[collNr].bookHists(iBooker);
0179 }
0180 booked_ = true;
0181 }
0182
0183
0184 void HLTGenResSource::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0185 edm::ParameterSetDescription desc;
0186
0187 desc.add<std::string>("dqmDirName", "HLTGenVal");
0188 desc.add<std::string>("hltProcessName", "HLT");
0189 desc.add<edm::InputTag>("trigEvent", edm::InputTag("hltTriggerSummaryAOD"));
0190 desc.add<edm::ParameterSetDescription>("genConfig", HLTGenValObjectMgr::makePSetDescription());
0191 desc.addVPSet("resCollConfigs", HLTGenResHistColl::makePSetDescription(), std::vector<edm::ParameterSet>());
0192
0193 descriptions.addDefault(desc);
0194 }
0195
0196
0197 DEFINE_FWK_MODULE(HLTGenResSource);