File indexing completed on 2025-03-06 03:07:18
0001 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0002 #include "DataFormats/HcalDigi/interface/HcalDigiHostCollection.h"
0003 #include "DataFormats/HcalDigi/interface/alpaka/HcalDigiDeviceCollection.h"
0004
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/Framework/interface/stream/EDProducer.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010
0011 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EDGetToken.h"
0012 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EDPutToken.h"
0013 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/Event.h"
0014 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EventSetup.h"
0015 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/stream/EDProducer.h"
0016
0017 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0018
0019 class HcalDigisSoAProducer : public stream::EDProducer<> {
0020 public:
0021 explicit HcalDigisSoAProducer(edm::ParameterSet const& ps);
0022 ~HcalDigisSoAProducer() override = default;
0023 static void fillDescriptions(edm::ConfigurationDescriptions&);
0024
0025 private:
0026 void produce(device::Event&, device::EventSetup const&) override;
0027
0028 private:
0029
0030 edm::EDGetTokenT<HBHEDigiCollection> hbheDigiToken_;
0031 edm::EDGetTokenT<QIE11DigiCollection> qie11DigiToken_;
0032
0033
0034 using HostCollectionPhase1 = hcal::Phase1DigiHostCollection;
0035 using HostCollectionPhase0 = hcal::Phase0DigiHostCollection;
0036
0037
0038 edm::EDPutTokenT<HostCollectionPhase1> digisF01HEToken_;
0039 edm::EDPutTokenT<HostCollectionPhase0> digisF5HBToken_;
0040 edm::EDPutTokenT<HostCollectionPhase1> digisF3HBToken_;
0041
0042 struct ConfigParameters {
0043 uint32_t maxChannelsF01HE, maxChannelsF5HB, maxChannelsF3HB;
0044 };
0045 ConfigParameters config_;
0046 };
0047
0048 void HcalDigisSoAProducer::fillDescriptions(edm::ConfigurationDescriptions& confDesc) {
0049 edm::ParameterSetDescription desc;
0050
0051 desc.add<edm::InputTag>("hbheDigisLabel", edm::InputTag("hcalDigis"));
0052 desc.add<edm::InputTag>("qie11DigiLabel", edm::InputTag("hcalDigis"));
0053 desc.add<std::string>("digisLabelF01HE", std::string{"f01HEDigis"});
0054 desc.add<std::string>("digisLabelF5HB", std::string{"f5HBDigis"});
0055 desc.add<std::string>("digisLabelF3HB", std::string{"f3HBDigis"});
0056 desc.add<uint32_t>("maxChannelsF01HE", 10000u);
0057 desc.add<uint32_t>("maxChannelsF5HB", 10000u);
0058 desc.add<uint32_t>("maxChannelsF3HB", 10000u);
0059
0060 confDesc.addWithDefaultLabel(desc);
0061 }
0062
0063 HcalDigisSoAProducer::HcalDigisSoAProducer(const edm::ParameterSet& ps)
0064 : EDProducer(ps),
0065 hbheDigiToken_{consumes(ps.getParameter<edm::InputTag>("hbheDigisLabel"))},
0066 qie11DigiToken_{consumes(ps.getParameter<edm::InputTag>("qie11DigiLabel"))},
0067 digisF01HEToken_{produces(ps.getParameter<std::string>("digisLabelF01HE"))},
0068 digisF5HBToken_{produces(ps.getParameter<std::string>("digisLabelF5HB"))},
0069 digisF3HBToken_{produces(ps.getParameter<std::string>("digisLabelF3HB"))} {
0070 config_.maxChannelsF01HE = ps.getParameter<uint32_t>("maxChannelsF01HE");
0071 config_.maxChannelsF5HB = ps.getParameter<uint32_t>("maxChannelsF5HB");
0072 config_.maxChannelsF3HB = ps.getParameter<uint32_t>("maxChannelsF3HB");
0073 }
0074
0075 void HcalDigisSoAProducer::produce(device::Event& event, device::EventSetup const& setup) {
0076 const auto& hbheDigis = event.get(hbheDigiToken_);
0077 const auto& qie11Digis = event.get(qie11DigiToken_);
0078
0079
0080 const int stride = HBHEDataFrame::MAXSAMPLES / 2 + 1;
0081 const int size = hbheDigis.size() * stride;
0082
0083
0084 HostCollectionPhase0 hf5_(size, event.queue());
0085
0086
0087 hf5_.view().stride() = stride;
0088 hf5_.view().size() = hbheDigis.size();
0089
0090 for (unsigned int i = 0; i < hbheDigis.size(); ++i) {
0091 auto const hbhe = hbheDigis[i];
0092 auto const id = hbhe.id().rawId();
0093 auto const presamples = hbhe.presamples();
0094 uint16_t header_word = (1 << 15) | (0x5 << 12) | (0 << 10) | ((hbhe.sample(0).capid() & 0x3) << 8);
0095
0096 auto hf5_vi = hf5_.view()[i];
0097 hf5_vi.ids() = id;
0098 hf5_vi.npresamples() = presamples;
0099 hf5_vi.data()[0] = header_word;
0100
0101 for (unsigned int i = 0; i < hf5_.view().stride() - 1; i++) {
0102 uint16_t s0 = (0 << 7) | (static_cast<uint8_t>(hbhe.sample(2 * i).adc()) & 0x7f);
0103 uint16_t s1 = (0 << 7) | (static_cast<uint8_t>(hbhe.sample(2 * i + 1).adc()) & 0x7f);
0104 uint16_t sample = (s1 << 8) | s0;
0105 hf5_vi.data()[i + 1] = sample;
0106 }
0107 }
0108 event.emplace(digisF5HBToken_, std::move(hf5_));
0109
0110 if (qie11Digis.empty()) {
0111 event.emplace(digisF01HEToken_, 0, event.queue());
0112 event.emplace(digisF3HBToken_, 0, event.queue());
0113
0114 } else {
0115 auto size_f1 = 0;
0116 auto size_f3 = 0;
0117
0118
0119 for (unsigned int i = 0; i < qie11Digis.size(); i++) {
0120 auto const digi = QIE11DataFrame{qie11Digis[i]};
0121
0122 if (digi.flavor() == 0 or digi.flavor() == 1) {
0123 if (digi.detid().subdetId() == HcalEndcap) {
0124 size_f1++;
0125 }
0126 } else if (digi.flavor() == 3) {
0127 if (digi.detid().subdetId() == HcalBarrel) {
0128 size_f3++;
0129 }
0130 }
0131 }
0132 auto const nsamples = qie11Digis.samples();
0133 auto const stride01 = nsamples * QIE11DataFrame::WORDS_PER_SAMPLE + QIE11DataFrame::HEADER_WORDS;
0134
0135
0136 HostCollectionPhase1 hf1_(size_f1, event.queue());
0137 HostCollectionPhase1 hf3_(size_f3, event.queue());
0138
0139
0140 hf1_.view().stride() = stride01;
0141 hf3_.view().stride() = stride01;
0142
0143 unsigned int i_f1 = 0;
0144 unsigned int i_f3 = 0;
0145
0146 for (unsigned int i = 0; i < qie11Digis.size(); i++) {
0147 auto const digi = QIE11DataFrame{qie11Digis[i]};
0148 assert(digi.samples() == qie11Digis.samples() && "collection nsamples must equal per digi samples");
0149
0150 if (digi.flavor() == 0 or digi.flavor() == 1) {
0151 if (digi.detid().subdetId() != HcalEndcap)
0152 continue;
0153 auto hf01_vi = hf1_.view()[i_f1];
0154
0155 hf01_vi.ids() = digi.detid().rawId();
0156 for (int hw = 0; hw < QIE11DataFrame::HEADER_WORDS + digi.samples(); hw++) {
0157 hf01_vi.data()[hw] = (qie11Digis[i][hw]);
0158 }
0159 i_f1++;
0160 } else if (digi.flavor() == 3) {
0161 if (digi.detid().subdetId() != HcalBarrel)
0162 continue;
0163 auto hf03_vi = hf3_.view()[i_f3];
0164
0165 hf03_vi.ids() = digi.detid().rawId();
0166
0167 for (int hw = 0; hw < QIE11DataFrame::HEADER_WORDS + digi.samples(); hw++) {
0168 hf03_vi.data()[hw] = (qie11Digis[i][hw]);
0169 }
0170 i_f3++;
0171 }
0172 }
0173
0174 hf1_.view().size() = size_f1;
0175 hf3_.view().size() = size_f3;
0176
0177 event.emplace(digisF01HEToken_, std::move(hf1_));
0178 event.emplace(digisF3HBToken_, std::move(hf3_));
0179 }
0180 }
0181 }
0182 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0183 DEFINE_FWK_ALPAKA_MODULE(HcalDigisSoAProducer);