File indexing completed on 2025-05-27 01:56:29
0001
0002 #include "FWCore/Framework/interface/ESWatcher.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "FWCore/Utilities/interface/StreamID.h"
0010 #include <iomanip> // for std::setw
0011 #include <future>
0012 #include <chrono>
0013
0014
0015 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/stream/EDProducer.h"
0016 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EDPutToken.h"
0017 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0018 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/Event.h"
0019 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0020 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToDevice.h"
0021
0022
0023 #include "DataFormats/HGCalDigi/interface/HGCalDigiHost.h"
0024 #include "DataFormats/HGCalDigi/interface/alpaka/HGCalDigiDevice.h"
0025 #include "DataFormats/HGCalRecHit/interface/HGCalRecHitHost.h"
0026 #include "DataFormats/HGCalRecHit/interface/alpaka/HGCalRecHitDevice.h"
0027
0028
0029 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0030 #include "CondFormats/DataRecord/interface/HGCalModuleConfigurationRcd.h"
0031 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0032 #include "CondFormats/HGCalObjects/interface/HGCalCalibrationParameterHost.h"
0033 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalCalibrationParameterDevice.h"
0034 #include "RecoLocalCalo/HGCalRecAlgos/interface/alpaka/HGCalRecHitCalibrationAlgorithms.h"
0035 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalMappingParameterDevice.h"
0036 #include "CondFormats/DataRecord/interface/HGCalDenseIndexInfoRcd.h"
0037
0038
0039
0040
0041 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0042
0043 using namespace cms::alpakatools;
0044
0045 class HGCalRecHitsProducer : public stream::EDProducer<> {
0046 public:
0047 explicit HGCalRecHitsProducer(const edm::ParameterSet&);
0048 static void fillDescriptions(edm::ConfigurationDescriptions&);
0049
0050 private:
0051 void produce(device::Event&, device::EventSetup const&) override;
0052 edm::ESWatcher<HGCalElectronicsMappingRcd> calibWatcher_;
0053 const edm::EDGetTokenT<hgcaldigi::HGCalDigiHost> digisToken_;
0054 edm::ESGetToken<hgcalrechit::HGCalCalibParamHost, HGCalModuleConfigurationRcd> calibToken_;
0055 device::ESGetToken<hgcal::HGCalMappingCellParamDevice, HGCalElectronicsMappingRcd> mappingToken_;
0056 device::ESGetToken<hgcal::HGCalDenseIndexInfoDevice, HGCalDenseIndexInfoRcd> indexingToken_;
0057 const device::EDPutToken<hgcalrechit::HGCalRecHitDevice> recHitsToken_;
0058 const HGCalRecHitCalibrationAlgorithms calibrator_;
0059 const int n_hits_scale;
0060 };
0061
0062 HGCalRecHitsProducer::HGCalRecHitsProducer(const edm::ParameterSet& iConfig)
0063 : EDProducer(iConfig),
0064 digisToken_{consumes<hgcaldigi::HGCalDigiHost>(iConfig.getParameter<edm::InputTag>("digis"))},
0065 calibToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("calibSource"))},
0066 mappingToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("mappingSource"))},
0067 indexingToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("indexingSource"))},
0068 recHitsToken_{produces()},
0069 calibrator_{iConfig.getParameter<int>("n_blocks"), iConfig.getParameter<int>("n_threads")},
0070 n_hits_scale{iConfig.getParameter<int>("n_hits_scale")} {
0071 #ifndef HGCAL_PERF_TEST
0072 if (n_hits_scale > 1) {
0073 throw cms::Exception("RuntimeError") << "Build with `HGCAL_PERF_TEST` flag to activate `n_hits_scale`.";
0074 }
0075 #endif
0076 }
0077
0078 void HGCalRecHitsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0079 edm::ParameterSetDescription desc;
0080 desc.add<edm::InputTag>("digis", edm::InputTag("hgcalDigis", "DIGI", "TEST"));
0081 desc.add("calibSource", edm::ESInputTag{})->setComment("Label for calibration parameters");
0082 desc.add("mappingSource", edm::ESInputTag{})->setComment("Label for cell mapping parameters");
0083 desc.add("indexingSource", edm::ESInputTag{})->setComment("Label for cell dense indexer");
0084 desc.add<int>("n_blocks", -1);
0085 desc.add<int>("n_threads", -1);
0086 desc.add<int>("n_hits_scale", -1);
0087 descriptions.addWithDefaultLabel(desc);
0088 }
0089
0090 void HGCalRecHitsProducer::produce(device::Event& iEvent, device::EventSetup const& iSetup) {
0091 auto& queue = iEvent.queue();
0092
0093
0094 auto const& hostCalibParamProvider = iSetup.getData(calibToken_);
0095 auto const& deviceMappingCellParamProvider = iSetup.getData(mappingToken_);
0096 auto const& deviceIndexingParamProvider = iSetup.getData(indexingToken_);
0097 auto const& hostDigisIn = iEvent.get(digisToken_);
0098
0099
0100 LogDebug("HGCalCalibrationParameter").log([&](auto& log) {
0101 if (calibWatcher_.check(iSetup)) {
0102 for (int i = 0; i < hostCalibParamProvider.view().metadata().size(); i++) {
0103 log << "idx = " << i << ", "
0104 << "ADC_ped = " << hostCalibParamProvider.view()[i].ADC_ped() << ", "
0105 << "CM_slope = " << hostCalibParamProvider.view()[i].CM_slope() << ", "
0106 << "CM_ped = " << hostCalibParamProvider.view()[i].CM_ped() << ", "
0107 << "BXm1_slope = " << hostCalibParamProvider.view()[i].BXm1_slope() << ", ";
0108 }
0109 }
0110 });
0111
0112 #ifdef HGCAL_PERF_TEST
0113 uint32_t oldSize = hostDigisIn.view().metadata().size();
0114 uint32_t newSize = oldSize * (n_hits_scale > 0 ? (unsigned)n_hits_scale : 1);
0115 auto hostDigis = HGCalDigiHost(newSize, queue);
0116 auto hostCalibParam = HGCalCalibParamHost(newSize, queue);
0117
0118 for (uint32_t i = 0; i < newSize; i++) {
0119 hostDigis.view()[i].tctp() = hostDigisIn.view()[i % oldSize].tctp();
0120 hostDigis.view()[i].adcm1() = hostDigisIn.view()[i % oldSize].adcm1();
0121 hostDigis.view()[i].adc() = hostDigisIn.view()[i % oldSize].adc();
0122 hostDigis.view()[i].tot() = hostDigisIn.view()[i % oldSize].tot();
0123 hostDigis.view()[i].toa() = hostDigisIn.view()[i % oldSize].toa();
0124 hostDigis.view()[i].cm() = hostDigisIn.view()[i % oldSize].cm();
0125 hostDigis.view()[i].flags() = hostDigisIn.view()[i % oldSize].flags();
0126
0127 hostCalibParam.view()[i].ADC_ped() = hostCalibParamProvider.view()[i % oldSize].ADC_ped();
0128 hostCalibParam.view()[i].Noise() = hostCalibParamProvider.view()[i % oldSize].Noise();
0129 hostCalibParam.view()[i].CM_slope() = hostCalibParamProvider.view()[i % oldSize].CM_slope();
0130 hostCalibParam.view()[i].CM_ped() = hostCalibParamProvider.view()[i % oldSize].CM_ped();
0131 hostCalibParam.view()[i].BXm1_slope() = hostCalibParamProvider.view()[i % oldSize].BXm1_slope();
0132 hostCalibParam.view()[i].TOTtoADC() = hostCalibParamProvider.view()[i % oldSize].TOTtoADC();
0133 hostCalibParam.view()[i].TOT_ped() = hostCalibParamProvider.view()[i % oldSize].TOT_ped();
0134 hostCalibParam.view()[i].TOT_lin() = hostCalibParamProvider.view()[i % oldSize].TOT_lin();
0135 hostCalibParam.view()[i].TOT_P0() = hostCalibParamProvider.view()[i % oldSize].TOT_P0();
0136 hostCalibParam.view()[i].TOT_P1() = hostCalibParamProvider.view()[i % oldSize].TOT_P1();
0137 hostCalibParam.view()[i].TOT_P2() = hostCalibParamProvider.view()[i % oldSize].TOT_P2();
0138 hostCalibParam.view()[i].TOAtops() = hostCalibParamProvider.view()[i % oldSize].TOAtops();
0139 hostCalibParam.view()[i].MIPS_scale() = hostCalibParamProvider.view()[i % oldSize].MIPS_scale();
0140 hostCalibParam.view()[i].valid() = hostCalibParamProvider.view()[i % oldSize].valid();
0141 }
0142 #else
0143 const auto& hostDigis = hostDigisIn;
0144 const auto& hostCalibParam = hostCalibParamProvider;
0145 #endif
0146
0147 LogDebug("HGCalRecHitsProducer") << "Loaded host digis: " << hostDigis.view().metadata().size();
0148
0149 LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- calling calibrate method";
0150
0151 #ifdef EDM_ML_DEBUG
0152 alpaka::wait(queue);
0153 auto start = std::chrono::steady_clock::now();
0154 #endif
0155
0156 LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- Copying the calib to the device\n\n" << std::endl;
0157 HGCalCalibParamDevice deviceCalibParam(hostCalibParam.view().metadata().size(), queue);
0158 alpaka::memcpy(queue, deviceCalibParam.buffer(), hostCalibParam.const_buffer());
0159
0160 #ifdef HGCAL_PERF_TEST
0161 auto tmpRecHits = calibrator_.calibrate(queue, hostDigis, deviceCalibParam);
0162 HGCalRecHitDevice recHits(oldSize, queue);
0163 alpaka::memcpy(queue, recHits.buffer(), tmpRecHits.const_buffer(), oldSize);
0164 #else
0165 auto recHits = calibrator_.calibrate(
0166 queue, hostDigis, deviceCalibParam, deviceMappingCellParamProvider, deviceIndexingParamProvider);
0167 #endif
0168
0169 #ifdef EDM_ML_DEBUG
0170 alpaka::wait(queue);
0171 auto stop = std::chrono::steady_clock::now();
0172 std::chrono::duration<float> elapsed = stop - start;
0173 LogDebug("HGCalRecHitsProducer") << "Time spent calibrating " << hostDigis.view().metadata().size()
0174 << " digis: " << elapsed.count();
0175 #endif
0176
0177 LogDebug("HGCalRecHitsProducer") << "\n\nINFO -- storing rec hits in the event";
0178 iEvent.emplace(recHitsToken_, std::move(recHits));
0179 }
0180
0181 }
0182
0183
0184 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0185 DEFINE_FWK_ALPAKA_MODULE(HGCalRecHitsProducer);