File indexing completed on 2024-04-06 12:25:48
0001 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalRecoParamsWithPulseShapesGPU.h"
0002
0003 #include "CondFormats/HcalObjects/interface/HcalRecoParams.h"
0004 #include "CalibCalorimetry/HcalAlgos/interface/HcalPulseShapes.h"
0005 #include "RecoLocalCalo/HcalRecAlgos/interface/PulseShapeFunctor.h"
0006
0007 #include "FWCore/Utilities/interface/typelookup.h"
0008 #include "HeterogeneousCore/CUDAUtilities/interface/cudaCheck.h"
0009
0010 #include <unordered_map>
0011
0012
0013 HcalRecoParamsWithPulseShapesGPU::HcalRecoParamsWithPulseShapesGPU(HcalRecoParams const& recoParams)
0014 : totalChannels_{recoParams.getAllContainers()[0].second.size() + recoParams.getAllContainers()[1].second.size()},
0015 param1_(totalChannels_),
0016 param2_(totalChannels_),
0017 ids_(totalChannels_) {
0018 #ifdef HCAL_MAHI_CPUDEBUG
0019 printf("hello from a reco params with pulse shapes\n");
0020 #endif
0021
0022 auto const containers = recoParams.getAllContainers();
0023
0024 HcalPulseShapes pulseShapes;
0025 std::unordered_map<unsigned int, uint32_t> idCache;
0026
0027
0028 auto const& barrelValues = containers[0].second;
0029 for (uint64_t i = 0; i < barrelValues.size(); ++i) {
0030 param1_[i] = barrelValues[i].param1();
0031 param2_[i] = barrelValues[i].param2();
0032
0033 auto const pulseShapeId = barrelValues[i].pulseShapeID();
0034
0035
0036
0037 if (pulseShapeId == 0) {
0038 ids_[i] = 0;
0039 continue;
0040 }
0041 if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
0042
0043 auto const newId = idCache.size();
0044 idCache[pulseShapeId] = newId;
0045
0046 ids_[i] = newId;
0047
0048
0049 acc25nsVec_.resize(acc25nsVec_.size() + hcal::constants::maxPSshapeBin);
0050 diff25nsItvlVec_.resize(diff25nsItvlVec_.size() + hcal::constants::maxPSshapeBin);
0051 accVarLenIdxMinusOneVec_.resize(accVarLenIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
0052 diffVarItvlIdxMinusOneVec_.resize(diffVarItvlIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
0053 accVarLenIdxZEROVec_.resize(accVarLenIdxZEROVec_.size() + hcal::constants::nsPerBX);
0054 diffVarItvlIdxZEROVec_.resize(diffVarItvlIdxZEROVec_.size() + hcal::constants::nsPerBX);
0055
0056
0057 auto const& pulseShape = pulseShapes.getShape(pulseShapeId);
0058 FitterFuncs::PulseShapeFunctor functor{pulseShape, false, false, false, 1, 0, 0, hcal::constants::maxSamples};
0059 auto const offset256 = newId * hcal::constants::maxPSshapeBin;
0060 auto const offset25 = newId * hcal::constants::nsPerBX;
0061 auto const numShapes = newId;
0062 for (int i = 0; i < hcal::constants::maxPSshapeBin; i++) {
0063 acc25nsVec_[offset256 * numShapes + i] = functor.acc25nsVec()[i];
0064 diff25nsItvlVec_[offset256 * numShapes + i] = functor.diff25nsItvlVec()[i];
0065 }
0066
0067 for (int i = 0; i < hcal::constants::nsPerBX; i++) {
0068 accVarLenIdxMinusOneVec_[offset25 * numShapes + i] = functor.accVarLenIdxMinusOneVec()[i];
0069 diffVarItvlIdxMinusOneVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxMinusOneVec()[i];
0070 accVarLenIdxZEROVec_[offset25 * numShapes + i] = functor.accVarLenIdxZEROVec()[i];
0071 diffVarItvlIdxZEROVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxZEROVec()[i];
0072 }
0073 } else {
0074
0075 ids_[i] = iter->second;
0076 }
0077 #ifdef HCAL_MAHI_CPUDEBUG
0078 if (barrelValues[i].rawId() == DETID_TO_DEBUG) {
0079 printf("recoShapeId = %u myid = %u\n", pulseShapeId, ids_[i]);
0080 }
0081 #endif
0082 }
0083
0084
0085 auto const& endcapValues = containers[1].second;
0086 auto const offset = barrelValues.size();
0087 for (uint64_t i = 0; i < endcapValues.size(); ++i) {
0088 param1_[i + offset] = endcapValues[i].param1();
0089 param2_[i + offset] = endcapValues[i].param2();
0090
0091 auto const pulseShapeId = endcapValues[i].pulseShapeID();
0092
0093
0094
0095 if (pulseShapeId == 0) {
0096 ids_[i + offset] = 0;
0097 continue;
0098 }
0099 if (auto const iter = idCache.find(pulseShapeId); iter == idCache.end()) {
0100
0101 auto const newId = idCache.size();
0102 idCache[pulseShapeId] = newId;
0103
0104 ids_[i + offset] = newId;
0105
0106
0107 acc25nsVec_.resize(acc25nsVec_.size() + hcal::constants::maxPSshapeBin);
0108 diff25nsItvlVec_.resize(diff25nsItvlVec_.size() + hcal::constants::maxPSshapeBin);
0109 accVarLenIdxMinusOneVec_.resize(accVarLenIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
0110 diffVarItvlIdxMinusOneVec_.resize(diffVarItvlIdxMinusOneVec_.size() + hcal::constants::nsPerBX);
0111 accVarLenIdxZEROVec_.resize(accVarLenIdxZEROVec_.size() + hcal::constants::nsPerBX);
0112 diffVarItvlIdxZEROVec_.resize(diffVarItvlIdxZEROVec_.size() + hcal::constants::nsPerBX);
0113
0114
0115 auto const& pulseShape = pulseShapes.getShape(pulseShapeId);
0116 FitterFuncs::PulseShapeFunctor functor{pulseShape, false, false, false, 1, 0, 0, hcal::constants::maxSamples};
0117 auto const offset256 = newId * hcal::constants::maxPSshapeBin;
0118 auto const offset25 = newId * hcal::constants::nsPerBX;
0119 auto const numShapes = newId;
0120 for (int i = 0; i < hcal::constants::maxPSshapeBin; i++) {
0121 acc25nsVec_[offset256 * numShapes + i] = functor.acc25nsVec()[i];
0122 diff25nsItvlVec_[offset256 * numShapes + i] = functor.diff25nsItvlVec()[i];
0123 }
0124
0125 for (int i = 0; i < hcal::constants::nsPerBX; i++) {
0126 accVarLenIdxMinusOneVec_[offset25 * numShapes + i] = functor.accVarLenIdxMinusOneVec()[i];
0127 diffVarItvlIdxMinusOneVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxMinusOneVec()[i];
0128 accVarLenIdxZEROVec_[offset25 * numShapes + i] = functor.accVarLenIdxZEROVec()[i];
0129 diffVarItvlIdxZEROVec_[offset25 * numShapes + i] = functor.diffVarItvlIdxZEROVec()[i];
0130 }
0131 } else {
0132
0133 ids_[i + offset] = iter->second;
0134 }
0135 }
0136
0137 #ifdef HCAL_MAHI_CPUDEBUG
0138 for (auto const& p : idCache)
0139 printf("recoPulseShapeId = %u id = %u\n", p.first, p.second);
0140 #endif
0141 }
0142
0143 HcalRecoParamsWithPulseShapesGPU::Product::~Product() {
0144
0145 cudaCheck(cudaFree(param1));
0146 cudaCheck(cudaFree(param2));
0147 cudaCheck(cudaFree(ids));
0148 cudaCheck(cudaFree(acc25nsVec));
0149 cudaCheck(cudaFree(diff25nsItvlVec));
0150 cudaCheck(cudaFree(accVarLenIdxMinusOneVec));
0151 cudaCheck(cudaFree(diffVarItvlIdxMinusOneVec));
0152 cudaCheck(cudaFree(accVarLenIdxZEROVec));
0153 cudaCheck(cudaFree(diffVarItvlIdxZEROVec));
0154 }
0155
0156 HcalRecoParamsWithPulseShapesGPU::Product const& HcalRecoParamsWithPulseShapesGPU::getProduct(
0157 cudaStream_t cudaStream) const {
0158 auto const& product = product_.dataForCurrentDeviceAsync(
0159 cudaStream, [this](HcalRecoParamsWithPulseShapesGPU::Product& product, cudaStream_t cudaStream) {
0160
0161 cudaCheck(cudaMalloc((void**)&product.param1, this->param1_.size() * sizeof(uint32_t)));
0162 cudaCheck(cudaMalloc((void**)&product.param2, this->param2_.size() * sizeof(uint32_t)));
0163 cudaCheck(cudaMalloc((void**)&product.ids, this->ids_.size() * sizeof(uint32_t)));
0164 cudaCheck(cudaMalloc((void**)&product.acc25nsVec, this->acc25nsVec_.size() * sizeof(float)));
0165 cudaCheck(cudaMalloc((void**)&product.diff25nsItvlVec, this->diff25nsItvlVec_.size() * sizeof(float)));
0166 cudaCheck(cudaMalloc((void**)&product.accVarLenIdxMinusOneVec,
0167 this->accVarLenIdxMinusOneVec_.size() * sizeof(float)));
0168 cudaCheck(cudaMalloc((void**)&product.diffVarItvlIdxMinusOneVec,
0169 this->diffVarItvlIdxMinusOneVec_.size() * sizeof(float)));
0170 cudaCheck(cudaMalloc((void**)&product.accVarLenIdxZEROVec, this->accVarLenIdxZEROVec_.size() * sizeof(float)));
0171 cudaCheck(
0172 cudaMalloc((void**)&product.diffVarItvlIdxZEROVec, this->diffVarItvlIdxZEROVec_.size() * sizeof(float)));
0173
0174
0175 cudaCheck(cudaMemcpyAsync(product.param1,
0176 this->param1_.data(),
0177 this->param1_.size() * sizeof(uint32_t),
0178 cudaMemcpyHostToDevice,
0179 cudaStream));
0180 cudaCheck(cudaMemcpyAsync(product.param2,
0181 this->param2_.data(),
0182 this->param2_.size() * sizeof(uint32_t),
0183 cudaMemcpyHostToDevice,
0184 cudaStream));
0185 cudaCheck(cudaMemcpyAsync(
0186 product.ids, this->ids_.data(), this->ids_.size() * sizeof(uint32_t), cudaMemcpyHostToDevice, cudaStream));
0187 cudaCheck(cudaMemcpyAsync(product.acc25nsVec,
0188 this->acc25nsVec_.data(),
0189 this->acc25nsVec_.size() * sizeof(float),
0190 cudaMemcpyHostToDevice,
0191 cudaStream));
0192 cudaCheck(cudaMemcpyAsync(product.diff25nsItvlVec,
0193 this->diff25nsItvlVec_.data(),
0194 this->diff25nsItvlVec_.size() * sizeof(float),
0195 cudaMemcpyHostToDevice,
0196 cudaStream));
0197 cudaCheck(cudaMemcpyAsync(product.accVarLenIdxMinusOneVec,
0198 this->accVarLenIdxMinusOneVec_.data(),
0199 this->accVarLenIdxMinusOneVec_.size() * sizeof(float),
0200 cudaMemcpyHostToDevice,
0201 cudaStream));
0202 cudaCheck(cudaMemcpyAsync(product.diffVarItvlIdxMinusOneVec,
0203 this->diffVarItvlIdxMinusOneVec_.data(),
0204 this->diffVarItvlIdxMinusOneVec_.size() * sizeof(float),
0205 cudaMemcpyHostToDevice,
0206 cudaStream));
0207 cudaCheck(cudaMemcpyAsync(product.accVarLenIdxZEROVec,
0208 this->accVarLenIdxZEROVec_.data(),
0209 this->accVarLenIdxZEROVec_.size() * sizeof(float),
0210 cudaMemcpyHostToDevice,
0211 cudaStream));
0212 cudaCheck(cudaMemcpyAsync(product.diffVarItvlIdxZEROVec,
0213 this->diffVarItvlIdxZEROVec_.data(),
0214 this->diffVarItvlIdxZEROVec_.size() * sizeof(float),
0215 cudaMemcpyHostToDevice,
0216 cudaStream));
0217 });
0218
0219 return product;
0220 }
0221
0222 TYPELOOKUP_DATA_REG(HcalRecoParamsWithPulseShapesGPU);