File indexing completed on 2024-04-06 12:02:15
0001 #include <cmath>
0002
0003 #include "CondFormats/HcalObjects/interface/HcalConvertedPedestalWidthsGPU.h"
0004 #include "FWCore/Utilities/interface/typelookup.h"
0005 #include "HeterogeneousCore/CUDAUtilities/interface/copyAsync.h"
0006
0007 namespace {
0008 float convert(
0009 float const value, float const width, int const i, HcalQIECoder const& coder, HcalQIEShape const& shape) {
0010 float const y = value;
0011 float const x = width;
0012 unsigned const x1 = static_cast<unsigned>(std::floor(y));
0013 unsigned const x2 = static_cast<unsigned>(std::floor(y + 1.));
0014 unsigned iun = static_cast<unsigned>(i);
0015 float const y1 = coder.charge(shape, x1, iun);
0016 float const y2 = coder.charge(shape, x2, iun);
0017 return (y2 - y1) * x;
0018 }
0019 }
0020
0021
0022 HcalConvertedPedestalWidthsGPU::HcalConvertedPedestalWidthsGPU(HcalPedestals const& pedestals,
0023 HcalPedestalWidths const& pedestalWidths,
0024 HcalQIEData const& qieData,
0025 HcalQIETypes const& qieTypes)
0026 : totalChannels_{pedestals.getAllContainers()[0].second.size() + pedestals.getAllContainers()[1].second.size()},
0027 values_(totalChannels_ * 4) {
0028 #ifdef HCAL_MAHI_CPUDEBUG
0029 std::cout << "hello from converted pedestal widths" << std::endl;
0030 std::cout << "pedestals HB values = " << pedestals.getAllContainers()[0].second.size()
0031 << " HE values = " << pedestals.getAllContainers()[1].second.size() << std::endl;
0032 std::cout << "qiedata HB values = " << qieData.getAllContainers()[0].second.size()
0033 << " HE values = " << qieData.getAllContainers()[1].second.size() << std::endl;
0034 #endif
0035
0036
0037 auto const pedestalsAll = pedestals.getAllContainers();
0038 auto const pedestalWidthsAll = pedestalWidths.getAllContainers();
0039 auto const qieDataAll = qieData.getAllContainers();
0040 auto const qieTypesAll = qieTypes.getAllContainers();
0041
0042
0043 auto const unitIsADC = pedestals.isADC();
0044
0045
0046 auto const& pedestalBarrelValues = pedestalsAll[0].second;
0047 auto const& pedestalWidthBarrelValues = pedestalWidthsAll[0].second;
0048 auto const& qieDataBarrelValues = qieDataAll[0].second;
0049 auto const& qieTypesBarrelValues = qieTypesAll[0].second;
0050
0051 #ifdef HCAL_MAHI_CPUDEBUG
0052 assert(pedestalWidthBarrelValues.size() == pedestalBarrelValues.size());
0053 assert(pedestalBarrelValues.size() == qieDataBarrelValues.size());
0054 assert(pedestalBarrelValues.size() == qieTypesBarrelValues.size());
0055 #endif
0056
0057 for (uint64_t i = 0; i < pedestalBarrelValues.size(); ++i) {
0058 auto const& qieCoder = qieDataBarrelValues[i];
0059 auto const qieType = qieTypesBarrelValues[i].getValue() > 1 ? 1 : 0;
0060 auto const& qieShape = qieData.getShape(qieType);
0061
0062 values_[i * 4] =
0063 unitIsADC
0064 ? convert(
0065 pedestalBarrelValues[i].getValue(0), pedestalWidthBarrelValues[i].getWidth(0), 0, qieCoder, qieShape)
0066 : pedestalWidthBarrelValues[i].getWidth(0);
0067 values_[i * 4 + 1] =
0068 unitIsADC
0069 ? convert(
0070 pedestalBarrelValues[i].getValue(1), pedestalWidthBarrelValues[i].getWidth(1), 1, qieCoder, qieShape)
0071 : pedestalWidthBarrelValues[i].getWidth(1);
0072 values_[i * 4 + 2] =
0073 unitIsADC
0074 ? convert(
0075 pedestalBarrelValues[i].getValue(2), pedestalWidthBarrelValues[i].getWidth(2), 2, qieCoder, qieShape)
0076 : pedestalWidthBarrelValues[i].getWidth(2);
0077 values_[i * 4 + 3] =
0078 unitIsADC
0079 ? convert(
0080 pedestalBarrelValues[i].getValue(3), pedestalWidthBarrelValues[i].getWidth(3), 3, qieCoder, qieShape)
0081 : pedestalWidthBarrelValues[i].getWidth(3);
0082 }
0083
0084
0085 auto const& pedestalEndcapValues = pedestalsAll[1].second;
0086 auto const& pedestalWidthEndcapValues = pedestalWidthsAll[1].second;
0087 auto const& qieDataEndcapValues = qieDataAll[1].second;
0088 auto const& qieTypesEndcapValues = qieTypesAll[1].second;
0089
0090 #ifdef HCAL_MAHI_CPUDEBUG
0091 assert(pedestalWidthEndcapValues.size() == pedestalEndcapValues.size());
0092 assert(pedestalEndcapValues.size() == qieDataEndcapValues.size());
0093 assert(pedestalEndcapValues.size() == qieTypesEndcapValues.size());
0094 #endif
0095
0096 auto const offset = pedestalWidthBarrelValues.size();
0097 for (uint64_t i = 0; i < pedestalEndcapValues.size(); ++i) {
0098 auto const& qieCoder = qieDataEndcapValues[i];
0099 auto const qieType = qieTypesEndcapValues[i].getValue() > 1 ? 1 : 0;
0100 auto const& qieShape = qieData.getShape(qieType);
0101 auto const off = offset + i;
0102
0103 values_[off * 4] =
0104 unitIsADC
0105 ? convert(
0106 pedestalEndcapValues[i].getValue(0), pedestalWidthEndcapValues[i].getWidth(0), 0, qieCoder, qieShape)
0107 : pedestalWidthEndcapValues[i].getWidth(0);
0108 values_[off * 4 + 1] =
0109 unitIsADC
0110 ? convert(
0111 pedestalEndcapValues[i].getValue(1), pedestalWidthEndcapValues[i].getWidth(1), 1, qieCoder, qieShape)
0112 : pedestalWidthEndcapValues[i].getWidth(1);
0113 values_[off * 4 + 2] =
0114 unitIsADC
0115 ? convert(
0116 pedestalEndcapValues[i].getValue(2), pedestalWidthEndcapValues[i].getWidth(2), 2, qieCoder, qieShape)
0117 : pedestalWidthEndcapValues[i].getWidth(2);
0118 values_[off * 4 + 3] =
0119 unitIsADC
0120 ? convert(
0121 pedestalEndcapValues[i].getValue(3), pedestalWidthEndcapValues[i].getWidth(3), 3, qieCoder, qieShape)
0122 : pedestalWidthEndcapValues[i].getWidth(3);
0123
0124 #ifdef HCAL_MAHI_CPUDEBUG
0125 if (pedestalEndcapValues[i].rawId() == DETID_TO_DEBUG) {
0126 for (int i = 0; i < 4; i++)
0127 printf("pedestalWidth(%d) = %f original pedestalWidth(%d) = %f\n",
0128 i,
0129 values_[off * 4 + i],
0130 i,
0131 pedestalWidthEndcapValues[i].getWidth(3));
0132 }
0133 #endif
0134 }
0135 }
0136
0137 HcalConvertedPedestalWidthsGPU::Product const& HcalConvertedPedestalWidthsGPU::getProduct(cudaStream_t stream) const {
0138 auto const& product = product_.dataForCurrentDeviceAsync(
0139 stream, [this](HcalConvertedPedestalWidthsGPU::Product& product, cudaStream_t stream) {
0140
0141 product.values = cms::cuda::make_device_unique<float[]>(values_.size(), stream);
0142
0143
0144 cms::cuda::copyAsync(product.values, values_, stream);
0145 });
0146
0147 return product;
0148 }
0149
0150 TYPELOOKUP_DATA_REG(HcalConvertedPedestalWidthsGPU);