File indexing completed on 2025-04-30 22:24:39
0001 #include "LCToSCAssociatorByEnergyScoreProducer.h"
0002
0003 #include <memory>
0004
0005 template <typename HIT>
0006 LCToSCAssociatorByEnergyScoreProducer<HIT>::LCToSCAssociatorByEnergyScoreProducer(const edm::ParameterSet &ps)
0007 : hitMap_(consumes<std::unordered_map<DetId, const unsigned int>>(ps.getParameter<edm::InputTag>("hitMapTag"))),
0008 caloGeometry_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
0009 hardScatterOnly_(ps.getParameter<bool>("hardScatterOnly")),
0010 hits_label_(ps.getParameter<std::vector<edm::InputTag>>("hits")) {
0011 for (auto &label : hits_label_) {
0012 if constexpr (std::is_same_v<HIT, HGCRecHit>)
0013 hgcal_hits_token_.push_back(consumes<HGCRecHitCollection>(label));
0014 else
0015 hits_token_.push_back(consumes<std::vector<HIT>>(label));
0016 }
0017
0018 rhtools_ = std::make_shared<hgcal::RecHitTools>();
0019
0020
0021 produces<ticl::LayerClusterToSimClusterAssociator>();
0022 }
0023
0024 template <typename HIT>
0025 LCToSCAssociatorByEnergyScoreProducer<HIT>::~LCToSCAssociatorByEnergyScoreProducer() {}
0026
0027 template <typename HIT>
0028 void LCToSCAssociatorByEnergyScoreProducer<HIT>::produce(edm::StreamID,
0029 edm::Event &iEvent,
0030 const edm::EventSetup &es) const {
0031 edm::ESHandle<CaloGeometry> geom = es.getHandle(caloGeometry_);
0032 rhtools_->setGeometry(*geom);
0033
0034 std::vector<const HIT *> hits;
0035 if constexpr (std::is_same_v<HIT, HGCRecHit>) {
0036 for (auto &token : hgcal_hits_token_) {
0037 edm::Handle<HGCRecHitCollection> hits_handle;
0038 iEvent.getByToken(token, hits_handle);
0039
0040
0041 if (!hits_handle.isValid()) {
0042 edm::LogWarning("LCToSCAssociatorByEnergyScoreProducer")
0043 << "Hit collection not available for token. Skipping this collection.";
0044 continue;
0045 }
0046
0047 for (const auto &hit : *hits_handle) {
0048 hits.push_back(&hit);
0049 }
0050 }
0051 } else {
0052 for (auto &token : hits_token_) {
0053 edm::Handle<std::vector<HIT>> hits_handle;
0054 iEvent.getByToken(token, hits_handle);
0055
0056
0057 if (!hits_handle.isValid()) {
0058 edm::LogWarning("LCToSCAssociatorByEnergyScoreProducer")
0059 << "Hit collection not available for token. Skipping this collection.";
0060 continue;
0061 }
0062
0063 for (const auto &hit : *hits_handle) {
0064 hits.push_back(&hit);
0065 }
0066 }
0067 }
0068 const auto hitMap = &iEvent.get(hitMap_);
0069
0070 auto impl = std::make_unique<LCToSCAssociatorByEnergyScoreImpl<HIT>>(
0071 iEvent.productGetter(), hardScatterOnly_, rhtools_, hitMap, hits);
0072 auto toPut = std::make_unique<ticl::LayerClusterToSimClusterAssociator>(std::move(impl));
0073 iEvent.put(std::move(toPut));
0074 }
0075
0076 template <typename HIT>
0077 void LCToSCAssociatorByEnergyScoreProducer<HIT>::fillDescriptions(edm::ConfigurationDescriptions &cfg) {
0078 edm::ParameterSetDescription desc;
0079 desc.add<bool>("hardScatterOnly", true);
0080 if constexpr (std::is_same_v<HIT, HGCRecHit>) {
0081 desc.add<edm::InputTag>("hitMapTag", edm::InputTag("recHitMapProducer", "hgcalRecHitMap"));
0082 desc.add<std::vector<edm::InputTag>>("hits",
0083 {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
0084 edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
0085 edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
0086 } else {
0087 desc.add<edm::InputTag>("hitMapTag", edm::InputTag("recHitMapProducer", "barrelRecHitMap"));
0088 desc.add<std::vector<edm::InputTag>>("hits",
0089 {edm::InputTag("particleFlowRecHitECAL", ""),
0090 edm::InputTag("particleFlowRecHitHBHE", ""),
0091 edm::InputTag("particleFlowRecHitHO", "")});
0092 }
0093 cfg.addWithDefaultLabel(desc);
0094 }