File indexing completed on 2024-05-29 23:13:06
0001 #include <utility>
0002 #include <vector>
0003
0004 #include <alpaka/alpaka.hpp>
0005
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/global/EDProducer.h"
0011 #include "RecoParticleFlow/PFRecHitProducer/interface/PFRecHitParamsRecord.h"
0012 #include "RecoParticleFlow/PFRecHitProducer/interface/PFRecHitTopologyRecord.h"
0013
0014 #include "CalorimeterDefinitions.h"
0015 #include "PFRecHitProducerKernel.h"
0016
0017 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0018 using namespace particleFlowRecHitProducer;
0019
0020 template <typename CAL>
0021 class PFRecHitSoAProducer : public global::EDProducer<> {
0022 public:
0023 PFRecHitSoAProducer(edm::ParameterSet const& config)
0024 : topologyToken_(esConsumes(config.getParameter<edm::ESInputTag>("topology"))),
0025 pfRecHitsToken_(produces()),
0026 synchronise_(config.getUntrackedParameter<bool>("synchronise")) {
0027 const std::vector<edm::ParameterSet> producers = config.getParameter<std::vector<edm::ParameterSet>>("producers");
0028 recHitsToken_.reserve(producers.size());
0029 for (const edm::ParameterSet& producer : producers) {
0030 recHitsToken_.emplace_back(consumes(producer.getParameter<edm::InputTag>("src")),
0031 esConsumes(producer.getParameter<edm::ESInputTag>("params")));
0032 }
0033 }
0034
0035 void produce(edm::StreamID, device::Event& event, const device::EventSetup& setup) const override {
0036 const typename CAL::TopologyTypeDevice& topology = setup.getData(topologyToken_);
0037
0038 uint32_t num_recHits = 0;
0039 for (const auto& token : recHitsToken_)
0040 num_recHits += event.get(token.first)->metadata().size();
0041
0042 reco::PFRecHitDeviceCollection pfRecHits{(int)num_recHits, event.queue()};
0043
0044 if (num_recHits != 0) {
0045 PFRecHitProducerKernel<CAL> kernel{event.queue(), num_recHits};
0046 for (const auto& token : recHitsToken_)
0047 kernel.processRecHits(
0048 event.queue(), event.get(token.first), setup.getData(token.second), topology, pfRecHits);
0049 kernel.associateTopologyInfo(event.queue(), topology, pfRecHits);
0050 }
0051
0052 if (synchronise_)
0053 alpaka::wait(event.queue());
0054
0055 event.emplace(pfRecHitsToken_, std::move(pfRecHits));
0056 }
0057
0058 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0059 edm::ParameterSetDescription desc;
0060 edm::ParameterSetDescription producers;
0061 producers.add<edm::InputTag>("src", edm::InputTag(""))->setComment("Input CaloRecHitSoA");
0062 producers.add<edm::ESInputTag>("params", edm::ESInputTag(""))->setComment("Quality cut parameters");
0063 std::vector<edm::ParameterSet> producersDefault(1);
0064 producersDefault[0].addParameter<edm::InputTag>("src", edm::InputTag(""));
0065 producersDefault[0].addParameter<edm::ESInputTag>("params", edm::ESInputTag(""));
0066 desc.addVPSet("producers", producers, producersDefault)->setComment("List of inputs and quality cuts");
0067 desc.add<edm::ESInputTag>("topology", edm::ESInputTag(""))->setComment("Topology information");
0068 desc.addUntracked<bool>("synchronise", false)
0069 ->setComment("Add synchronisation point after execution (for benchmarking asynchronous execution)");
0070 descriptions.addWithDefaultLabel(desc);
0071 }
0072
0073 private:
0074 const device::ESGetToken<typename CAL::TopologyTypeDevice, typename CAL::TopologyRecordType> topologyToken_;
0075 std::vector<std::pair<device::EDGetToken<typename CAL::CaloRecHitSoATypeDevice>,
0076 device::ESGetToken<typename CAL::ParameterType, typename CAL::ParameterRecordType>>>
0077 recHitsToken_;
0078 const device::EDPutToken<reco::PFRecHitDeviceCollection> pfRecHitsToken_;
0079 const bool synchronise_;
0080 };
0081
0082 using PFRecHitSoAProducerHCAL = PFRecHitSoAProducer<HCAL>;
0083 using PFRecHitSoAProducerECAL = PFRecHitSoAProducer<ECAL>;
0084 }
0085
0086 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0087 DEFINE_FWK_ALPAKA_MODULE(PFRecHitSoAProducerHCAL);
0088 DEFINE_FWK_ALPAKA_MODULE(PFRecHitSoAProducerECAL);