Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-11 23:28:11

0001 #include <unordered_map>
0002 
0003 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
0004 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0006 #include "DataFormats/ParticleFlowReco/interface/PFRecHitHostCollection.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/stream/EDProducer.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0013 #include "FWCore/Utilities/interface/EDGetToken.h"
0014 #include "FWCore/Utilities/interface/EDPutToken.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0017 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0018 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0019 
0020 class LegacyPFRecHitProducer : public edm::stream::EDProducer<> {
0021 public:
0022   LegacyPFRecHitProducer(edm::ParameterSet const& config)
0023       : alpakaPfRecHitsToken_(consumes(config.getParameter<edm::InputTag>("src"))),
0024         legacyPfRecHitsToken_(produces()),
0025         geomToken_(esConsumes<edm::Transition::BeginRun>()) {}
0026 
0027   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0028     edm::ParameterSetDescription desc;
0029     desc.add<edm::InputTag>("src", edm::InputTag(""));
0030     descriptions.addWithDefaultLabel(desc);
0031   }
0032 
0033 private:
0034   void beginRun(edm::Run const&, edm::EventSetup const&) override;
0035   void produce(edm::Event&, const edm::EventSetup&) override;
0036 
0037   const edm::EDGetTokenT<reco::PFRecHitHostCollection> alpakaPfRecHitsToken_;
0038   const edm::EDPutTokenT<reco::PFRecHitCollection> legacyPfRecHitsToken_;
0039 
0040   edm::ESGetToken<CaloGeometry, CaloGeometryRecord> geomToken_;
0041   std::unordered_map<PFLayer::Layer, const CaloSubdetectorGeometry*> caloGeo_;
0042 };
0043 
0044 void LegacyPFRecHitProducer::beginRun(edm::Run const&, const edm::EventSetup& setup) {
0045   const CaloGeometry& geo = setup.getData(geomToken_);
0046   caloGeo_[PFLayer::HCAL_BARREL1] = geo.getSubdetectorGeometry(DetId::Hcal, HcalSubdetector::HcalBarrel);
0047   caloGeo_[PFLayer::HCAL_ENDCAP] = geo.getSubdetectorGeometry(DetId::Hcal, HcalSubdetector::HcalEndcap);
0048   caloGeo_[PFLayer::ECAL_BARREL] = geo.getSubdetectorGeometry(DetId::Ecal, EcalSubdetector::EcalBarrel);
0049   caloGeo_[PFLayer::ECAL_ENDCAP] = geo.getSubdetectorGeometry(DetId::Ecal, EcalSubdetector::EcalEndcap);
0050 }
0051 
0052 void LegacyPFRecHitProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0053   const reco::PFRecHitHostCollection& pfRecHitsAlpakaSoA = event.get(alpakaPfRecHitsToken_);
0054   const reco::PFRecHitHostCollection::ConstView& alpakaPfRecHits = pfRecHitsAlpakaSoA.const_view();
0055 
0056   reco::PFRecHitCollection out;
0057 
0058   if (alpakaPfRecHits.metadata().size() != 0) {
0059     out.reserve(alpakaPfRecHits.size());
0060     for (size_t i = 0; i < alpakaPfRecHits.size(); i++) {
0061       reco::PFRecHit& pfrh =
0062           out.emplace_back(caloGeo_.at(alpakaPfRecHits[i].layer())->getGeometry(alpakaPfRecHits[i].detId()),
0063                            alpakaPfRecHits[i].detId(),
0064                            alpakaPfRecHits[i].layer(),
0065                            alpakaPfRecHits[i].energy());
0066       pfrh.setTime(alpakaPfRecHits[i].time());
0067       pfrh.setDepth(alpakaPfRecHits[i].depth());
0068 
0069       // order in Alpaka:   N, S, E, W,NE,SW,SE,NW
0070       const short eta[8] = {0, 0, 1, -1, 1, -1, 1, -1};
0071       const short phi[8] = {1, -1, 0, 0, 1, -1, -1, 1};
0072       for (size_t k = 0; k < 8; k++)
0073         if (alpakaPfRecHits[i].neighbours()(k) != -1)
0074           pfrh.addNeighbour(eta[k], phi[k], 0, alpakaPfRecHits[i].neighbours()(k));
0075     }
0076   }
0077 
0078   event.emplace(legacyPfRecHitsToken_, out);
0079 }
0080 
0081 #include "FWCore/Framework/interface/MakerMacros.h"
0082 DEFINE_FWK_MODULE(LegacyPFRecHitProducer);