File indexing completed on 2024-04-06 12:26:19
0001 #ifndef RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
0002 #define RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
0003
0004 #include <cstdint>
0005 #include <cstdio>
0006 #include <algorithm>
0007 #include <limits>
0008
0009 #include "CUDADataFormats/SiPixelCluster/interface/gpuClusteringConstants.h"
0010 #include "CondFormats/SiPixelObjects/interface/SiPixelGainForHLTonGPU.h"
0011 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0012 #include "HeterogeneousCore/CUDAUtilities/interface/cuda_assert.h"
0013 #include "RecoLocalTracker/SiPixelClusterizer/interface/SiPixelClusterThresholds.h"
0014
0015 namespace gpuCalibPixel {
0016
0017 using gpuClustering::invalidModuleId;
0018
0019
0020 __global__ void calibDigis(SiPixelClusterThresholds clusterThresholds,
0021 uint16_t* id,
0022 uint16_t const* __restrict__ x,
0023 uint16_t const* __restrict__ y,
0024 uint16_t* adc,
0025 SiPixelGainForHLTonGPU const* __restrict__ ped,
0026 int numElements,
0027 uint32_t* __restrict__ moduleStart,
0028 uint32_t* __restrict__ nClustersInModule,
0029 uint32_t* __restrict__ clusModuleStart
0030 ) {
0031 int first = blockDim.x * blockIdx.x + threadIdx.x;
0032
0033 const float VCaltoElectronGain = clusterThresholds.vCaltoElectronGain;
0034 const float VCaltoElectronGain_L1 = clusterThresholds.vCaltoElectronGain_L1;
0035 const float VCaltoElectronOffset = clusterThresholds.vCaltoElectronOffset;
0036 const float VCaltoElectronOffset_L1 = clusterThresholds.vCaltoElectronOffset_L1;
0037
0038
0039 if (0 == first)
0040 clusModuleStart[0] = moduleStart[0] = 0;
0041 for (int i = first; i < phase1PixelTopology::numberOfModules; i += gridDim.x * blockDim.x) {
0042 nClustersInModule[i] = 0;
0043 }
0044
0045 for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
0046 if (invalidModuleId == id[i])
0047 continue;
0048
0049 bool isDeadColumn = false, isNoisyColumn = false;
0050
0051 int row = x[i];
0052 int col = y[i];
0053
0054 auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
0055 float pedestal = ret.first;
0056 float gain = ret.second;
0057
0058 if (isDeadColumn | isNoisyColumn) {
0059 printf("bad pixel at %d in %d\n", i, id[i]);
0060 id[i] = invalidModuleId;
0061 adc[i] = 0;
0062 } else {
0063 float vcal = float(adc[i]) * gain - pedestal * gain;
0064
0065 float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
0066 float offset = id[i] < 96 ? VCaltoElectronOffset_L1 : VCaltoElectronOffset;
0067 vcal = vcal * conversionFactor + offset;
0068
0069 adc[i] = std::clamp(int(vcal), 100, int(std::numeric_limits<uint16_t>::max()));
0070 }
0071 }
0072 }
0073
0074 __global__ void calibDigisPhase2(SiPixelClusterThresholds clusterThresholds,
0075 uint16_t* id,
0076 uint16_t* adc,
0077 int numElements,
0078 uint32_t* __restrict__ moduleStart,
0079 uint32_t* __restrict__ nClustersInModule,
0080 uint32_t* __restrict__ clusModuleStart
0081 ) {
0082 int first = blockDim.x * blockIdx.x + threadIdx.x;
0083
0084
0085 const float ElectronPerADCGain = clusterThresholds.electronPerADCGain;
0086 const int8_t Phase2ReadoutMode = clusterThresholds.phase2ReadoutMode;
0087 const uint16_t Phase2DigiBaseline = clusterThresholds.phase2DigiBaseline;
0088 const uint8_t Phase2KinkADC = clusterThresholds.phase2KinkADC;
0089
0090 if (0 == first)
0091 clusModuleStart[0] = moduleStart[0] = 0;
0092 for (int i = first; i < phase2PixelTopology::numberOfModules; i += gridDim.x * blockDim.x) {
0093 nClustersInModule[i] = 0;
0094 }
0095
0096 for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
0097 if (invalidModuleId == id[i])
0098 continue;
0099
0100 const int mode = (Phase2ReadoutMode < -1 ? -1 : Phase2ReadoutMode);
0101
0102 int adc_int = adc[i];
0103
0104 if (mode < 0)
0105 adc_int = int(adc_int * ElectronPerADCGain);
0106 else {
0107 if (adc_int < Phase2KinkADC)
0108 adc_int = int((adc_int + 0.5) * ElectronPerADCGain);
0109 else {
0110 const int8_t dspp = (Phase2ReadoutMode < 10 ? Phase2ReadoutMode : 10);
0111 const int8_t ds = int8_t(dspp <= 1 ? 1 : (dspp - 1) * (dspp - 1));
0112
0113 adc_int -= Phase2KinkADC;
0114 adc_int *= ds;
0115 adc_int += Phase2KinkADC;
0116
0117 adc_int = ((adc_int + 0.5 * ds) * ElectronPerADCGain);
0118 }
0119
0120 adc_int += int(Phase2DigiBaseline);
0121 }
0122 adc[i] = std::min(adc_int, int(std::numeric_limits<uint16_t>::max()));
0123 }
0124 }
0125
0126 }
0127
0128 #endif