File indexing completed on 2025-03-10 00:37:00
0001 #include <memory>
0002
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/stream/EDProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Framework/interface/ESWatcher.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/Utilities/interface/StreamID.h"
0010
0011 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0012 #include "DataFormats/HGCalDigi/interface/HGCalElectronicsId.h"
0013 #include "DataFormats/HGCalDigi/interface/HGCalDigiHost.h"
0014 #include "DataFormats/HGCalDigi/interface/HGCalECONDPacketInfoHost.h"
0015 #include "DataFormats/HGCalDigi/interface/HGCalFEDPacketInfoHost.h"
0016 #include "DataFormats/HGCalDigi/interface/HGCalRawDataDefinitions.h"
0017
0018 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0019 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0020
0021 #include "CondFormats/DataRecord/interface/HGCalModuleConfigurationRcd.h"
0022 #include "CondFormats/HGCalObjects/interface/HGCalConfiguration.h"
0023
0024 #include "oneapi/tbb/task_arena.h"
0025 #include "oneapi/tbb.h"
0026
0027 #include "EventFilter/HGCalRawToDigi/interface/HGCalUnpacker.h"
0028 class HGCalRawToDigi : public edm::stream::EDProducer<> {
0029 public:
0030 explicit HGCalRawToDigi(const edm::ParameterSet&);
0031 uint16_t callUnpacker(unsigned fedId,
0032 const FEDRawData& fed_data,
0033 const HGCalMappingModuleIndexer& moduleIndexer,
0034 const HGCalConfiguration& config,
0035 hgcaldigi::HGCalDigiHost& digis,
0036 hgcaldigi::HGCalECONDPacketInfoHost& econdPacketInfo);
0037 static void fillDescriptions(edm::ConfigurationDescriptions&);
0038
0039 private:
0040 void produce(edm::Event&, const edm::EventSetup&) override;
0041 void beginRun(edm::Run const&, edm::EventSetup const&) override;
0042
0043
0044 const edm::EDGetTokenT<FEDRawDataCollection> fedRawToken_;
0045
0046
0047 const edm::EDPutTokenT<hgcaldigi::HGCalDigiHost> digisToken_;
0048 const edm::EDPutTokenT<hgcaldigi::HGCalECONDPacketInfoHost> econdPacketInfoToken_;
0049 const edm::EDPutTokenT<hgcaldigi::HGCalFEDPacketInfoHost> fedPacketInfoToken_;
0050
0051
0052
0053
0054
0055
0056
0057 edm::ESGetToken<HGCalMappingModuleIndexer, HGCalElectronicsMappingRcd> moduleIndexToken_;
0058 edm::ESGetToken<HGCalConfiguration, HGCalModuleConfigurationRcd> configToken_;
0059
0060
0061
0062
0063
0064
0065 HGCalUnpacker unpacker_;
0066
0067 const bool doSerial_;
0068 bool headersOnly_;
0069 };
0070
0071 HGCalRawToDigi::HGCalRawToDigi(const edm::ParameterSet& iConfig)
0072 : fedRawToken_(consumes<FEDRawDataCollection>(iConfig.getParameter<edm::InputTag>("src"))),
0073 digisToken_(produces<hgcaldigi::HGCalDigiHost>()),
0074 econdPacketInfoToken_(produces<hgcaldigi::HGCalECONDPacketInfoHost>()),
0075 fedPacketInfoToken_(produces<hgcaldigi::HGCalFEDPacketInfoHost>()),
0076
0077 moduleIndexToken_(esConsumes()),
0078 configToken_(esConsumes()),
0079 doSerial_(iConfig.getParameter<bool>("doSerial")),
0080 headersOnly_(iConfig.getParameter<bool>("headersOnly")) {}
0081
0082 void HGCalRawToDigi::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
0083
0084
0085 }
0086
0087 void HGCalRawToDigi::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0088
0089 const auto& moduleIndexer = iSetup.getData(moduleIndexToken_);
0090
0091 const auto& config = iSetup.getData(configToken_);
0092
0093 hgcaldigi::HGCalDigiHost digis(moduleIndexer.getMaxDataSize(), cms::alpakatools::host());
0094 hgcaldigi::HGCalECONDPacketInfoHost econdPacketInfo(moduleIndexer.getMaxModuleSize(), cms::alpakatools::host());
0095 hgcaldigi::HGCalFEDPacketInfoHost fedPacketInfo(moduleIndexer.fedCount(), cms::alpakatools::host());
0096
0097
0098 const auto& raw_data = iEvent.get(fedRawToken_);
0099
0100 for (int32_t i = 0; i < digis.view().metadata().size(); i++) {
0101 digis.view()[i].flags() = hgcal::DIGI_FLAG::NotAvailable;
0102 }
0103
0104
0105 if (doSerial_) {
0106 for (unsigned fedId = 0; fedId < moduleIndexer.fedCount(); ++fedId) {
0107 const auto& fed_data = raw_data.FEDData(fedId);
0108 fedPacketInfo.view()[fedId].FEDPayload() = fed_data.size();
0109 if (fed_data.size() == 0)
0110 continue;
0111 fedPacketInfo.view()[fedId].FEDUnpackingFlag() =
0112 callUnpacker(fedId, fed_data, moduleIndexer, config, digis, econdPacketInfo);
0113 }
0114 }
0115
0116 else {
0117 oneapi::tbb::this_task_arena::isolate([&]() {
0118 oneapi::tbb::parallel_for(0U, moduleIndexer.fedCount(), [&](unsigned fedId) {
0119 const auto& fed_data = raw_data.FEDData(fedId);
0120 fedPacketInfo.view()[fedId].FEDPayload() = fed_data.size();
0121 if (fed_data.size() == 0)
0122 return;
0123 fedPacketInfo.view()[fedId].FEDUnpackingFlag() =
0124 callUnpacker(fedId, fed_data, moduleIndexer, config, digis, econdPacketInfo);
0125 return;
0126 });
0127 });
0128 }
0129
0130
0131 iEvent.emplace(digisToken_, std::move(digis));
0132 iEvent.emplace(econdPacketInfoToken_, std::move(econdPacketInfo));
0133 iEvent.emplace(fedPacketInfoToken_, std::move(fedPacketInfo));
0134 }
0135
0136
0137 uint16_t HGCalRawToDigi::callUnpacker(unsigned fedId,
0138 const FEDRawData& fed_data,
0139 const HGCalMappingModuleIndexer& moduleIndexer,
0140 const HGCalConfiguration& config,
0141 hgcaldigi::HGCalDigiHost& digis,
0142 hgcaldigi::HGCalECONDPacketInfoHost& econdPacketInfo) {
0143 uint16_t status =
0144 unpacker_.parseFEDData(fedId, fed_data, moduleIndexer, config, digis, econdPacketInfo, headersOnly_);
0145 return status;
0146 }
0147
0148
0149 void HGCalRawToDigi::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0150 edm::ParameterSetDescription desc;
0151 desc.add<edm::InputTag>("src", edm::InputTag("rawDataCollector"));
0152 desc.add<std::vector<unsigned int> >("fedIds", {});
0153 desc.add<bool>("doSerial", true)->setComment("do not attempt to paralleize unpacking of different FEDs");
0154 desc.add<bool>("headersOnly", false)->setComment("unpack only headers");
0155 descriptions.add("hgcalDigis", desc);
0156 }
0157
0158
0159 DEFINE_FWK_MODULE(HGCalRawToDigi);