File indexing completed on 2024-10-16 05:06:39
0001
0002
0003
0004 #include "HitToLayerClusterAssociatorProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Framework/interface/ESHandle.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/Utilities/interface/EDGetToken.h"
0011 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0012 #include "DataFormats/HGCRecHit/interface/HGCRecHitCollections.h"
0013 #include "SimDataFormats/Associations/interface/TICLAssociationMap.h"
0014 #include "CommonTools/RecoAlgos/interface/MultiVectorManager.h"
0015
0016 HitToLayerClusterAssociatorProducer::HitToLayerClusterAssociatorProducer(const edm::ParameterSet &pset)
0017 : LCCollectionToken_(consumes<std::vector<reco::CaloCluster>>(pset.getParameter<edm::InputTag>("layer_clusters"))),
0018 hitMapToken_(consumes<std::unordered_map<DetId, unsigned int>>(pset.getParameter<edm::InputTag>("hitMap"))) {
0019 auto hitsTags = pset.getParameter<std::vector<edm::InputTag>>("hits");
0020 for (const auto &tag : hitsTags) {
0021 hitsTokens_.push_back(consumes<HGCRecHitCollection>(tag));
0022 }
0023 produces<ticl::AssociationMap<ticl::mapWithFraction>>("hitToLayerClusterMap");
0024 }
0025
0026 HitToLayerClusterAssociatorProducer::~HitToLayerClusterAssociatorProducer() {}
0027
0028 void HitToLayerClusterAssociatorProducer::produce(edm::StreamID,
0029 edm::Event &iEvent,
0030 const edm::EventSetup &iSetup) const {
0031 using namespace edm;
0032
0033 Handle<std::vector<reco::CaloCluster>> layer_clusters;
0034 iEvent.getByToken(LCCollectionToken_, layer_clusters);
0035
0036 Handle<std::unordered_map<DetId, unsigned int>> hitMap;
0037 iEvent.getByToken(hitMapToken_, hitMap);
0038
0039 MultiVectorManager<HGCRecHit> rechitManager;
0040 for (const auto &token : hitsTokens_) {
0041 Handle<HGCRecHitCollection> hitsHandle;
0042 iEvent.getByToken(token, hitsHandle);
0043 rechitManager.addVector(*hitsHandle);
0044 }
0045
0046
0047 auto hitToLayerClusterMap = std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(rechitManager.size());
0048
0049
0050 for (unsigned int lcId = 0; lcId < layer_clusters->size(); ++lcId) {
0051 const auto &layer_cluster = (*layer_clusters)[lcId];
0052
0053 for (const auto &hitAndFraction : layer_cluster.hitsAndFractions()) {
0054 auto hitMapIter = hitMap->find(hitAndFraction.first);
0055 if (hitMapIter != hitMap->end()) {
0056 unsigned int rechitIndex = hitMapIter->second;
0057 float fraction = hitAndFraction.second;
0058 hitToLayerClusterMap->insert(rechitIndex, lcId, fraction);
0059 }
0060 }
0061 }
0062 iEvent.put(std::move(hitToLayerClusterMap), "hitToLayerClusterMap");
0063 }
0064
0065 void HitToLayerClusterAssociatorProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0066 edm::ParameterSetDescription desc;
0067 desc.add<edm::InputTag>("layer_clusters", edm::InputTag("hgcalMergeLayerClusters"));
0068 desc.add<edm::InputTag>("hitMap", edm::InputTag("recHitMapProducer", "hgcalRecHitMap"));
0069 desc.add<std::vector<edm::InputTag>>("hits",
0070 {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
0071 edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
0072 edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
0073 descriptions.add("hitToLayerClusterAssociator", desc);
0074 }
0075
0076
0077 DEFINE_FWK_MODULE(HitToLayerClusterAssociatorProducer);