Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:17:47

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