File indexing completed on 2022-12-01 00:47:52
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 "HeterogeneousCore/CUDAUtilities/interface/cuda_assert.h"
0012 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0013
0014 namespace gpuCalibPixel {
0015
0016 using gpuClustering::invalidModuleId;
0017
0018
0019
0020 constexpr float VCaltoElectronGain = 47;
0021 constexpr float VCaltoElectronGain_L1 = 50;
0022 constexpr float VCaltoElectronOffset = -60;
0023 constexpr float VCaltoElectronOffset_L1 = -670;
0024 constexpr int VCalChargeThreshold = 100;
0025
0026 constexpr float ElectronPerADCGain = 1500;
0027 constexpr int8_t Phase2ReadoutMode = 3;
0028 constexpr uint16_t Phase2DigiBaseline = 1000;
0029 constexpr uint8_t Phase2KinkADC = 8;
0030
0031 template <bool isRun2>
0032 __global__ void calibDigis(uint16_t* id,
0033 uint16_t const* __restrict__ x,
0034 uint16_t const* __restrict__ y,
0035 uint16_t* adc,
0036 SiPixelGainForHLTonGPU const* __restrict__ ped,
0037 int numElements,
0038 uint32_t* __restrict__ moduleStart,
0039 uint32_t* __restrict__ nClustersInModule,
0040 uint32_t* __restrict__ clusModuleStart
0041 ) {
0042 int first = blockDim.x * blockIdx.x + threadIdx.x;
0043
0044
0045 if (0 == first)
0046 clusModuleStart[0] = moduleStart[0] = 0;
0047 for (int i = first; i < phase1PixelTopology::numberOfModules; i += gridDim.x * blockDim.x) {
0048 nClustersInModule[i] = 0;
0049 }
0050
0051 for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
0052 if (invalidModuleId == id[i])
0053 continue;
0054
0055 bool isDeadColumn = false, isNoisyColumn = false;
0056
0057 int row = x[i];
0058 int col = y[i];
0059 auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
0060 float pedestal = ret.first;
0061 float gain = ret.second;
0062
0063 if (isDeadColumn | isNoisyColumn) {
0064 printf("bad pixel at %d in %d\n", i, id[i]);
0065 id[i] = invalidModuleId;
0066 adc[i] = 0;
0067 } else {
0068 float vcal = float(adc[i]) * gain - pedestal * gain;
0069 if constexpr (isRun2) {
0070 float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
0071 float offset = id[i] < 96 ? VCaltoElectronOffset_L1 : VCaltoElectronOffset;
0072 vcal = vcal * conversionFactor + offset;
0073 }
0074 adc[i] = std::clamp(int(vcal), 100, int(std::numeric_limits<uint16_t>::max()));
0075 }
0076 }
0077 }
0078
0079 __global__ void calibDigisPhase2(uint16_t* id,
0080 uint16_t* adc,
0081 int numElements,
0082 uint32_t* __restrict__ moduleStart,
0083 uint32_t* __restrict__ nClustersInModule,
0084 uint32_t* __restrict__ clusModuleStart
0085 ) {
0086 int first = blockDim.x * blockIdx.x + threadIdx.x;
0087
0088
0089 if (0 == first)
0090 clusModuleStart[0] = moduleStart[0] = 0;
0091 for (int i = first; i < phase2PixelTopology::numberOfModules; i += gridDim.x * blockDim.x) {
0092 nClustersInModule[i] = 0;
0093 }
0094
0095 for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
0096 if (invalidModuleId == id[i])
0097 continue;
0098
0099 constexpr int mode = (Phase2ReadoutMode < -1 ? -1 : Phase2ReadoutMode);
0100
0101 int adc_int = adc[i];
0102
0103 if constexpr (mode < 0)
0104 adc_int = int(adc_int * ElectronPerADCGain);
0105 else {
0106 if (adc_int < Phase2KinkADC)
0107 adc_int = int((adc_int + 0.5) * ElectronPerADCGain);
0108 else {
0109 constexpr int8_t dspp = (Phase2ReadoutMode < 10 ? Phase2ReadoutMode : 10);
0110 constexpr int8_t ds = int8_t(dspp <= 1 ? 1 : (dspp - 1) * (dspp - 1));
0111
0112 adc_int -= Phase2KinkADC;
0113 adc_int *= ds;
0114 adc_int += Phase2KinkADC;
0115
0116 adc_int = ((adc_int + 0.5 * ds) * ElectronPerADCGain);
0117 }
0118
0119 adc_int += int(Phase2DigiBaseline);
0120 }
0121 adc[i] = std::min(adc_int, int(std::numeric_limits<uint16_t>::max()));
0122 }
0123 }
0124
0125 }
0126
0127 #endif