File indexing completed on 2024-06-22 02:24:07
0001 #include <utility>
0002
0003 #include <alpaka/alpaka.hpp>
0004
0005 #include "DataFormats/Common/interface/SortedCollection.h"
0006 #include "DataFormats/ParticleFlowReco/interface/CaloRecHitHostCollection.h"
0007 #include "DataFormats/ParticleFlowReco/interface/alpaka/CaloRecHitDeviceCollection.h"
0008 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "FWCore/Utilities/interface/EDGetToken.h"
0012 #include "FWCore/Utilities/interface/InputTag.h"
0013 #include "FWCore/Utilities/interface/StreamID.h"
0014 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/global/EDProducer.h"
0015 #include "CalorimeterDefinitions.h"
0016
0017 #include "DataFormats/HcalRecHit/interface/HcalRecHitHostCollection.h"
0018 #include "DataFormats/HcalRecHit/interface/alpaka/HcalRecHitDeviceCollection.h"
0019
0020 #define DEBUG false
0021
0022 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0023 using namespace particleFlowRecHitProducer;
0024
0025 template <typename CAL>
0026 class CaloRecHitSoAProducer : public global::EDProducer<> {
0027 public:
0028 CaloRecHitSoAProducer(edm::ParameterSet const& config)
0029 : recHitsToken_(consumes(config.getParameter<edm::InputTag>("src"))),
0030 deviceToken_(produces()),
0031 synchronise_(config.getUntrackedParameter<bool>("synchronise")) {}
0032
0033 void produce(edm::StreamID sid, device::Event& event, device::EventSetup const&) const override {
0034 const edm::SortedCollection<typename CAL::CaloRecHitType>& recHits = event.get(recHitsToken_);
0035 const int32_t num_recHits = recHits.size();
0036 if (DEBUG)
0037 printf("Found %d recHits\n", num_recHits);
0038
0039 hcal::RecHitHostCollection hostProduct{num_recHits, event.queue()};
0040 auto& view = hostProduct.view();
0041
0042 for (int i = 0; i < num_recHits; i++) {
0043 convertRecHit(view[i], recHits[i]);
0044
0045 if (DEBUG && i < 10)
0046 printf("recHit %4d %u %f %f\n", i, view.detId(i), view.energy(i), view.timeM0(i));
0047 }
0048
0049 hcal::RecHitDeviceCollection deviceProduct{num_recHits, event.queue()};
0050 alpaka::memcpy(event.queue(), deviceProduct.buffer(), hostProduct.buffer());
0051 if (synchronise_)
0052 alpaka::wait(event.queue());
0053 event.emplace(deviceToken_, std::move(deviceProduct));
0054 }
0055
0056 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0057 edm::ParameterSetDescription desc;
0058 desc.add<edm::InputTag>("src", edm::InputTag(""))->setComment("Input calorimeter rec hit collection");
0059 desc.addUntracked<bool>("synchronise", false)
0060 ->setComment("Add synchronisation point after execution (for benchmarking asynchronous execution)");
0061 descriptions.addWithDefaultLabel(desc);
0062 }
0063
0064 private:
0065 const edm::EDGetTokenT<edm::SortedCollection<typename CAL::CaloRecHitType>> recHitsToken_;
0066 const device::EDPutToken<hcal::RecHitDeviceCollection> deviceToken_;
0067 const bool synchronise_;
0068
0069 static void convertRecHit(hcal::RecHitHostCollection::View::element to, const typename CAL::CaloRecHitType& from);
0070 };
0071
0072 template <>
0073 void CaloRecHitSoAProducer<HCAL>::convertRecHit(hcal::RecHitHostCollection::View::element to,
0074 const HCAL::CaloRecHitType& from) {
0075
0076 to.detId() = from.id().rawId();
0077 to.energy() = from.energy();
0078 to.timeM0() = from.time();
0079 }
0080
0081
0082
0083
0084
0085
0086
0087
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098
0099
0100 using HCALRecHitSoAProducer = CaloRecHitSoAProducer<HCAL>;
0101
0102
0103
0104 }
0105
0106 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0107 DEFINE_FWK_ALPAKA_MODULE(HCALRecHitSoAProducer);
0108