File indexing completed on 2025-04-24 01:30:28
0001 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
0002 #include "DataFormats/EcalDigi/interface/EcalDigiPhase2HostCollection.h"
0003 #include "DataFormats/EcalDigi/interface/alpaka/EcalDigiPhase2DeviceCollection.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/Utilities/interface/EDGetToken.h"
0006 #include "FWCore/Utilities/interface/EDPutToken.h"
0007 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/global/EDProducer.h"
0008 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/Event.h"
0009 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EventSetup.h"
0010 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0011 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0012
0013 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0014
0015 class EcalPhase2DigiToPortableProducer : public global::EDProducer<> {
0016 public:
0017 explicit EcalPhase2DigiToPortableProducer(edm::ParameterSet const &ps);
0018 ~EcalPhase2DigiToPortableProducer() override = default;
0019 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0020
0021 void produce(edm::StreamID sid, device::Event &event, device::EventSetup const &setup) const override;
0022
0023 private:
0024 const edm::EDGetTokenT<EBDigiCollectionPh2> inputDigiToken_;
0025 const edm::EDPutTokenT<EcalDigiPhase2HostCollection> outputDigiHostToken_;
0026 };
0027
0028 void EcalPhase2DigiToPortableProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0029 edm::ParameterSetDescription desc;
0030
0031 desc.add<edm::InputTag>("BarrelDigis", edm::InputTag("simEcalUnsuppressedDigis", ""));
0032 desc.add<std::string>("digisLabelEB", "ebDigis");
0033
0034 descriptions.addWithDefaultLabel(desc);
0035 }
0036
0037 EcalPhase2DigiToPortableProducer::EcalPhase2DigiToPortableProducer(edm::ParameterSet const &ps)
0038 : EDProducer(ps),
0039 inputDigiToken_{consumes(ps.getParameter<edm::InputTag>("BarrelDigis"))},
0040 outputDigiHostToken_{produces(ps.getParameter<std::string>("digisLabelEB"))} {}
0041
0042 void EcalPhase2DigiToPortableProducer::produce(edm::StreamID sid,
0043 device::Event &event,
0044 device::EventSetup const &setup) const {
0045
0046 const auto &inputDigis = event.get(inputDigiToken_);
0047
0048 const uint32_t size = inputDigis.size();
0049
0050
0051 EcalDigiPhase2HostCollection digisHostColl{static_cast<int32_t>(size), event.queue()};
0052 auto digisHostCollView = digisHostColl.view();
0053
0054
0055 uint32_t i = 0;
0056 for (const auto &inputDigi : inputDigis) {
0057 const unsigned int nSamples = inputDigi.size();
0058
0059 digisHostCollView.id()[i] = inputDigi.id();
0060 if (nSamples > ecalPh2::sampleSize) {
0061 edm::LogError("size_mismatch") << "Number of input samples (" << nSamples
0062 << ") larger than the maximum sample size (" << ecalPh2::sampleSize
0063 << "). Ignoring the excess samples.";
0064 }
0065
0066 for (unsigned int sample = 0; sample < ecalPh2::sampleSize; ++sample) {
0067 if (sample < nSamples) {
0068
0069 EcalLiteDTUSample thisSample = inputDigi[sample];
0070
0071 digisHostCollView.data()[i][sample] = thisSample.raw();
0072 } else {
0073 digisHostCollView.data()[i][sample] = 0;
0074 }
0075 }
0076 ++i;
0077 }
0078 digisHostCollView.size() = i;
0079
0080
0081 event.emplace(outputDigiHostToken_, std::move(digisHostColl));
0082 }
0083 }
0084
0085 DEFINE_FWK_ALPAKA_MODULE(EcalPhase2DigiToPortableProducer);