Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-04-30 22:24:39

0001 // Author: Felice Pantaleo, felice.pantaleo@cern.ch 06/2024
0002 #include "HitToSimClusterCaloParticleAssociatorProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/Utilities/interface/EDGetToken.h"
0009 #include "DataFormats/HGCalReco/interface/Trackster.h"
0010 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0011 #include "SimDataFormats/Associations/interface/TICLAssociationMap.h"
0012 #include "DataFormats/Provenance/interface/ProductID.h"
0013 #include "DataFormats/HGCRecHit/interface/HGCRecHitCollections.h"
0014 #include "CommonTools/RecoAlgos/interface/MultiVectorManager.h"
0015 #include "SimDataFormats/CaloAnalysis/interface/CaloParticle.h"
0016 #include "SimDataFormats/CaloAnalysis/interface/SimCluster.h"
0017 
0018 HitToSimClusterCaloParticleAssociatorProducer::HitToSimClusterCaloParticleAssociatorProducer(
0019     const edm::ParameterSet &pset)
0020     : simClusterToken_(consumes<std::vector<SimCluster>>(pset.getParameter<edm::InputTag>("simClusters"))),
0021       caloParticleToken_(consumes<std::vector<CaloParticle>>(pset.getParameter<edm::InputTag>("caloParticles"))),
0022       hitMapToken_(consumes<std::unordered_map<DetId, const unsigned int>>(pset.getParameter<edm::InputTag>("hitMap"))),
0023       hitsTags_(pset.getParameter<std::vector<edm::InputTag>>("hits")) {
0024   for (const auto &tag : hitsTags_) {
0025     hitsTokens_.push_back(consumes<HGCRecHitCollection>(tag));
0026   }
0027   produces<ticl::AssociationMap<ticl::mapWithFraction>>("hitToSimClusterMap");
0028   produces<ticl::AssociationMap<ticl::mapWithFraction>>("hitToCaloParticleMap");
0029 }
0030 
0031 void HitToSimClusterCaloParticleAssociatorProducer::produce(edm::StreamID,
0032                                                             edm::Event &iEvent,
0033                                                             const edm::EventSetup &iSetup) const {
0034   using namespace edm;
0035 
0036   Handle<std::vector<CaloParticle>> caloParticlesHandle;
0037   iEvent.getByToken(caloParticleToken_, caloParticlesHandle);
0038   const auto &caloParticles = *caloParticlesHandle;
0039 
0040   Handle<std::vector<SimCluster>> simClustersHandle;
0041   iEvent.getByToken(simClusterToken_, simClustersHandle);
0042   Handle<std::unordered_map<DetId, const unsigned int>> hitMap;
0043   iEvent.getByToken(hitMapToken_, hitMap);
0044 
0045   MultiVectorManager<HGCRecHit> rechitManager;
0046   // Loop over tokens with index to access corresponding InputTag
0047   for (size_t i = 0; i < hitsTokens_.size(); ++i) {
0048     const auto &token = hitsTokens_[i];
0049     Handle<HGCRecHitCollection> hitsHandle;
0050     iEvent.getByToken(token, hitsHandle);
0051 
0052     // Error handling with tag name
0053     if (!hitsHandle.isValid()) {
0054       edm::LogWarning("HitToSimClusterCaloParticleAssociatorProducer")
0055           << "Missing HGCRecHitCollection for tag: " << hitsTags_[i].encode();
0056       continue;
0057     }
0058     rechitManager.addVector(*hitsHandle);
0059   }
0060 
0061   // Check if rechitManager is empty after processing hitsTokens_
0062   if (rechitManager.size() == 0) {
0063     edm::LogWarning("HitToSimClusterCaloParticleAssociatorProducer")
0064         << "No valid HGCRecHitCollections found. Association maps will be empty.";
0065     // Store empty maps in the event
0066     iEvent.put(std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(), "hitToSimClusterMap");
0067     iEvent.put(std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(), "hitToCaloParticleMap");
0068     return;
0069   }
0070 
0071   // Create association maps
0072   auto hitToSimClusterMap = std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(rechitManager.size());
0073   auto hitToCaloParticleMap = std::make_unique<ticl::AssociationMap<ticl::mapWithFraction>>(rechitManager.size());
0074 
0075   // Loop over caloParticles
0076   for (unsigned int cpId = 0; cpId < caloParticles.size(); ++cpId) {
0077     const auto &caloParticle = caloParticles[cpId];
0078     // Loop over simClusters in caloParticle
0079     for (const auto &simCluster : caloParticle.simClusters()) {
0080       // Loop over hits in simCluster
0081       for (const auto &hitAndFraction : simCluster->hits_and_fractions()) {
0082         auto hitMapIter = hitMap->find(hitAndFraction.first);
0083         if (hitMapIter != hitMap->end()) {
0084           unsigned int rechitIndex = hitMapIter->second;
0085           float fraction = hitAndFraction.second;
0086           hitToSimClusterMap->insert(rechitIndex, simCluster.key(), fraction);
0087           hitToCaloParticleMap->insert(rechitIndex, cpId, fraction);
0088         }
0089       }
0090     }
0091   }
0092   iEvent.put(std::move(hitToSimClusterMap), "hitToSimClusterMap");
0093   iEvent.put(std::move(hitToCaloParticleMap), "hitToCaloParticleMap");
0094 }
0095 
0096 void HitToSimClusterCaloParticleAssociatorProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0097   edm::ParameterSetDescription desc;
0098   desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
0099   desc.add<edm::InputTag>("simClusters", edm::InputTag("mix", "MergedCaloTruth"));
0100 
0101   desc.add<edm::InputTag>("hitMap", edm::InputTag("recHitMapProducer", "hgcalRecHitMap"));
0102   desc.add<std::vector<edm::InputTag>>("hits",
0103                                        {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
0104                                         edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
0105                                         edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
0106   descriptions.add("hitToSimClusterCaloParticleAssociator", desc);
0107 }
0108 
0109 // Define this as a plug-in
0110 DEFINE_FWK_MODULE(HitToSimClusterCaloParticleAssociatorProducer);