File indexing completed on 2025-01-18 03:42:14
0001 #include "FWCore/ParameterSet/interface/FileInPath.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003
0004 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESGetToken.h"
0005 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ESProducer.h"
0006 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/ModuleFactory.h"
0007 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0008 #include "HeterogeneousCore/AlpakaInterface/interface/host.h"
0009 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0010
0011 #include "DataFormats/HGCalDigi/interface/HGCalElectronicsId.h"
0012 #include "CondFormats/HGCalObjects/interface/HGCalConfiguration.h"
0013 #include "CondFormats/HGCalObjects/interface/HGCalMappingModuleIndexer.h"
0014 #include "CondFormats/DataRecord/interface/HGCalElectronicsMappingRcd.h"
0015 #include "CondFormats/DataRecord/interface/HGCalModuleConfigurationRcd.h" // depends on HGCalElectronicsMappingRcd
0016 #include "CondFormats/HGCalObjects/interface/HGCalCalibrationParameterHost.h" // for HGCalConfigParamHost
0017 #include "CondFormats/HGCalObjects/interface/alpaka/HGCalCalibrationParameterDevice.h"
0018
0019 #include <string>
0020 #include <iostream>
0021 #include <iomanip> // for std::setw
0022 #include <sstream>
0023 #include <fstream> // needed to read json file with std::ifstream
0024 #include <nlohmann/json.hpp>
0025 using json = nlohmann::json;
0026
0027 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0028
0029 namespace hgcalrechit {
0030
0031 class HGCalCalibrationESProducer : public ESProducer {
0032 public:
0033 HGCalCalibrationESProducer(const edm::ParameterSet& iConfig)
0034 : ESProducer(iConfig), filename_(iConfig.getParameter<edm::FileInPath>("filename")) {
0035 auto cc = setWhatProduced(this);
0036 indexToken_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("indexSource"));
0037 configToken_ = cc.consumes(iConfig.getParameter<edm::ESInputTag>("configSource"));
0038 }
0039
0040 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0041 edm::ParameterSetDescription desc;
0042 desc.add<edm::FileInPath>("filename")->setComment("Path to JSON file with calibration parameters");
0043 desc.add<edm::ESInputTag>("indexSource", edm::ESInputTag(""))
0044 ->setComment("Label for module indexer to set SoA size");
0045 desc.add<edm::ESInputTag>("configSource", edm::ESInputTag(""))
0046 ->setComment("Label for ROC configuration parameters");
0047 descriptions.addWithDefaultLabel(desc);
0048 }
0049
0050 template <typename T>
0051 static void fill_SoA_column(
0052 T* column_SoA, const std::vector<T>& values, const int offset, const int nrows, int arr_offset = 0) {
0053
0054 const int nrows_vals = values.size();
0055 if (arr_offset < 0) {
0056 arr_offset = 0;
0057 if (nrows_vals != arr_offset + nrows) {
0058 edm::LogWarning("HGCalCalibrationESProducer")
0059 << " Expected " << nrows << " rows, but got " << nrows_vals << "!";
0060 }
0061 } else if (nrows_vals < arr_offset + nrows) {
0062 edm::LogWarning("HGCalCalibrationESProducer")
0063 << " Tried to copy " << nrows << " rows to SoA with offset " << arr_offset << ", but only have "
0064 << nrows_vals << " values in JSON!";
0065 }
0066 auto begin = values.begin() + arr_offset;
0067 auto end = (begin + nrows > values.end()) ? values.end() : begin + nrows;
0068 std::copy(begin, end, &column_SoA[offset]);
0069 }
0070
0071 std::optional<hgcalrechit::HGCalCalibParamHost> produce(const HGCalModuleConfigurationRcd& iRecord) {
0072 auto const& moduleMap = iRecord.get(indexToken_);
0073 auto const& config = iRecord.get(configToken_);
0074
0075
0076 const uint32_t nchans = moduleMap.getMaxDataSize();
0077 hgcalrechit::HGCalCalibParamHost product(nchans, cms::alpakatools::host());
0078
0079
0080 std::ifstream infile(filename_.fullPath().c_str());
0081 json calib_data = json::parse(infile);
0082 for (const auto& it : calib_data.items()) {
0083 const std::string& module = it.key();
0084 if (module == "Metadata")
0085 continue;
0086 const auto& [ifed, imod] = moduleMap.getIndexForFedAndModule(module);
0087 const uint32_t offset =
0088 moduleMap.getIndexForModuleData(module);
0089 const uint32_t nrows =
0090 calib_data[module]["Channel"].size();
0091 const uint32_t nrocs =
0092 config.feds[ifed].econds[imod].rocs.size();
0093
0094
0095 uint32_t nchans = (nrows % 39 == 0 ? 39 : 37);
0096 if (nrows % 37 != 0 and nrows % 39 != 0) {
0097 edm::LogWarning("HGCalCalibrationESProducer")
0098 << " nchannels%nchannels_per_erX!=0 nchannels=" << nrows << ", nchannels_per_erX=37 or 39!";
0099 }
0100
0101
0102 for (std::size_t iroc = 0; iroc < nrocs; ++iroc) {
0103 const uint32_t offset_arr = iroc * nchans;
0104 const uint32_t offset_soa = offset + offset_arr;
0105 if (offset_arr + nchans > nrows) {
0106 edm::LogWarning("HGCalCalibrationESProducer")
0107 << " offset + nchannels_per_eRx = " << offset_arr << " + " << nchans << " = " << offset_arr + nchans
0108 << " > " << nrows << " = nchannels ";
0109 }
0110 fill_SoA_column<float>(
0111 product.view().ADC_ped(), calib_data[module]["ADC_ped"], offset_soa, nchans, offset_arr);
0112 fill_SoA_column<float>(product.view().Noise(), calib_data[module]["Noise"], offset_soa, nchans, offset_arr);
0113 fill_SoA_column<float>(
0114 product.view().CM_slope(), calib_data[module]["CM_slope"], offset_soa, nchans, offset_arr);
0115 fill_SoA_column<float>(
0116 product.view().CM_ped(), calib_data[module]["CM_ped"], offset_soa, nchans, offset_arr);
0117 fill_SoA_column<float>(
0118 product.view().BXm1_slope(), calib_data[module]["BXm1_slope"], offset_soa, nchans, offset_arr);
0119 }
0120
0121 fill_SoA_column<float>(product.view().TOTtoADC(), calib_data[module]["TOTtoADC"], offset, nrows);
0122 fill_SoA_column<float>(product.view().TOT_ped(), calib_data[module]["TOT_ped"], offset, nrows);
0123 fill_SoA_column<float>(product.view().TOT_lin(), calib_data[module]["TOT_lin"], offset, nrows);
0124 fill_SoA_column<float>(product.view().TOT_P0(), calib_data[module]["TOT_P0"], offset, nrows);
0125 fill_SoA_column<float>(product.view().TOT_P1(), calib_data[module]["TOT_P1"], offset, nrows);
0126 fill_SoA_column<float>(product.view().TOT_P2(), calib_data[module]["TOT_P2"], offset, nrows);
0127 fill_SoA_column<float>(product.view().TOAtops(), calib_data[module]["TOAtops"], offset, nrows);
0128 fill_SoA_column<float>(product.view().MIPS_scale(), calib_data[module]["MIPS_scale"], offset, nrows);
0129 fill_SoA_column<unsigned char>(product.view().valid(), calib_data[module]["Valid"], offset, nrows);
0130 }
0131
0132 return product;
0133 }
0134
0135 private:
0136 edm::ESGetToken<HGCalMappingModuleIndexer, HGCalElectronicsMappingRcd> indexToken_;
0137 edm::ESGetToken<HGCalConfiguration, HGCalModuleConfigurationRcd> configToken_;
0138 const edm::FileInPath filename_;
0139 };
0140
0141 }
0142
0143 }
0144
0145 DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(hgcalrechit::HGCalCalibrationESProducer);