Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-11 23:37:47

0001 // user include files
0002 #include <unordered_map>
0003 
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/global/EDProducer.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0013 
0014 #include "DataFormats/HGCRecHit/interface/HGCRecHitCollections.h"
0015 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0016 #include "CommonTools/RecoAlgos/interface/MultiVectorManager.h"
0017 
0018 class RecHitMapProducer : public edm::global::EDProducer<> {
0019 public:
0020   RecHitMapProducer(const edm::ParameterSet&);
0021   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0022 
0023   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0024 
0025 private:
0026   std::vector<edm::EDGetTokenT<HGCRecHitCollection>> hgcal_hits_token_;
0027   std::vector<edm::EDGetTokenT<reco::PFRecHitCollection>> barrel_hits_token_;
0028 
0029   bool hgcalOnly_;
0030 };
0031 
0032 DEFINE_FWK_MODULE(RecHitMapProducer);
0033 
0034 using DetIdRecHitMap = std::unordered_map<DetId, const unsigned int>;
0035 
0036 RecHitMapProducer::RecHitMapProducer(const edm::ParameterSet& ps) : hgcalOnly_(ps.getParameter<bool>("hgcalOnly")) {
0037   std::vector<edm::InputTag> tags = ps.getParameter<std::vector<edm::InputTag>>("hits");
0038   for (auto& tag : tags) {
0039     if (tag.label().find("HGCalRecHit") != std::string::npos) {
0040       hgcal_hits_token_.push_back(consumes<HGCRecHitCollection>(tag));
0041     } else {
0042       barrel_hits_token_.push_back(consumes<reco::PFRecHitCollection>(tag));
0043     }
0044   }
0045 
0046   produces<DetIdRecHitMap>("hgcalRecHitMap");
0047   if (!hgcalOnly_)
0048     produces<DetIdRecHitMap>("barrelRecHitMap");
0049 }
0050 
0051 void RecHitMapProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0052   edm::ParameterSetDescription desc;
0053   desc.add<std::vector<edm::InputTag>>("hits",
0054                                        {edm::InputTag("HGCalRecHit", "HGCEERecHits"),
0055                                         edm::InputTag("HGCalRecHit", "HGCHEFRecHits"),
0056                                         edm::InputTag("HGCalRecHit", "HGCHEBRecHits")});
0057   desc.add<bool>("hgcalOnly", true);
0058   descriptions.add("recHitMapProducer", desc);
0059 }
0060 
0061 void RecHitMapProducer::produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const {
0062   auto hitMapHGCal = std::make_unique<DetIdRecHitMap>();
0063 
0064   // Retrieve collections
0065   const auto& ee_hits = evt.getHandle(hgcal_hits_token_[0]);
0066   const auto& fh_hits = evt.getHandle(hgcal_hits_token_[1]);
0067   const auto& bh_hits = evt.getHandle(hgcal_hits_token_[2]);
0068 
0069   // Check validity of all handles
0070   if (!ee_hits.isValid() || !fh_hits.isValid() || !bh_hits.isValid()) {
0071     edm::LogWarning("HGCalRecHitMapProducer") << "One or more hit collections are unavailable. Returning an empty map.";
0072     evt.put(std::move(hitMapHGCal), "hgcalRecHitMap");
0073     return;
0074   }
0075 
0076   // TODO may be worth to avoid dependency on the order
0077   // of the collections, maybe using a map
0078   MultiVectorManager<HGCRecHit> rechitManager;
0079   rechitManager.addVector(*ee_hits);
0080   rechitManager.addVector(*fh_hits);
0081   rechitManager.addVector(*bh_hits);
0082 
0083   for (unsigned int i = 0; i < rechitManager.size(); ++i) {
0084     const auto recHitDetId = rechitManager[i].detid();
0085     hitMapHGCal->emplace(recHitDetId, i);
0086   }
0087 
0088   evt.put(std::move(hitMapHGCal), "hgcalRecHitMap");
0089 
0090   if (!hgcalOnly_) {
0091     auto hitMapBarrel = std::make_unique<DetIdRecHitMap>();
0092     MultiVectorManager<reco::PFRecHit> barrelRechitManager;
0093     barrelRechitManager.addVector(evt.get(barrel_hits_token_[0]));
0094     barrelRechitManager.addVector(evt.get(barrel_hits_token_[1]));
0095     for (unsigned int i = 0; i < barrelRechitManager.size(); ++i) {
0096       const auto recHitDetId = barrelRechitManager[i].detId();
0097       hitMapBarrel->emplace(recHitDetId, i);
0098     }
0099     evt.put(std::move(hitMapBarrel), "barrelRecHitMap");
0100   }
0101 }