CMS 3D CMS Logo

gpuCalibPixel.h
Go to the documentation of this file.
1 #ifndef RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
2 #define RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
3 
4 #include <cstdint>
5 #include <cstdio>
6 #include <algorithm>
7 #include <limits>
8 
13 
14 namespace gpuCalibPixel {
15 
17 
18  // calibrationConstants
19  // valid for run2
20  constexpr float VCaltoElectronGain = 47; // L2-4: 47 +- 4.7
21  constexpr float VCaltoElectronGain_L1 = 50; // L1: 49.6 +- 2.6
22  constexpr float VCaltoElectronOffset = -60; // L2-4: -60 +- 130
23  constexpr float VCaltoElectronOffset_L1 = -670; // L1: -670 +- 220
24  constexpr int VCalChargeThreshold = 100;
25  //for phase2
26  constexpr float ElectronPerADCGain = 600;
27  constexpr int8_t Phase2ReadoutMode = 3;
28  constexpr uint16_t Phase2DigiBaseline = 1500;
29  constexpr uint8_t Phase2KinkADC = 8;
30 
31  template <bool isRun2>
32  __global__ void calibDigis(uint16_t* id,
33  uint16_t const* __restrict__ x,
34  uint16_t const* __restrict__ y,
35  uint16_t* adc,
36  SiPixelGainForHLTonGPU const* __restrict__ ped,
37  int numElements,
38  uint32_t* __restrict__ moduleStart, // just to zero first
39  uint32_t* __restrict__ nClustersInModule, // just to zero them
40  uint32_t* __restrict__ clusModuleStart // just to zero first
41  ) {
42  int first = blockDim.x * blockIdx.x + threadIdx.x;
43 
44  // zero for next kernels...
45  if (0 == first)
46  clusModuleStart[0] = moduleStart[0] = 0;
48  nClustersInModule[i] = 0;
49  }
50 
51  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
52  if (invalidModuleId == id[i])
53  continue;
54 
55  bool isDeadColumn = false, isNoisyColumn = false;
56 
57  int row = x[i];
58  int col = y[i];
59  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
60  float pedestal = ret.first;
61  float gain = ret.second;
62  // float pedestal = 0; float gain = 1.;
63  if (isDeadColumn | isNoisyColumn) {
64  printf("bad pixel at %d in %d\n", i, id[i]);
65  id[i] = invalidModuleId;
66  adc[i] = 0;
67  } else {
68  float vcal = float(adc[i]) * gain - pedestal * gain;
69  if constexpr (isRun2) {
70  float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
72  vcal = vcal * conversionFactor + offset;
73  }
74  adc[i] = std::clamp(int(vcal), 100, int(std::numeric_limits<uint16_t>::max()));
75  }
76  }
77  }
78 
79  __global__ void calibDigisPhase2(uint16_t* id,
80  uint16_t* adc,
81  int numElements,
82  uint32_t* __restrict__ moduleStart, // just to zero first
83  uint32_t* __restrict__ nClustersInModule, // just to zero them
84  uint32_t* __restrict__ clusModuleStart // just to zero first
85  ) {
86  int first = blockDim.x * blockIdx.x + threadIdx.x;
87  // zero for next kernels...
88 
89  if (0 == first)
90  clusModuleStart[0] = moduleStart[0] = 0;
92  nClustersInModule[i] = 0;
93  }
94 
95  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
96  if (invalidModuleId == id[i])
97  continue;
98 
99  constexpr int mode = (Phase2ReadoutMode < -1 ? -1 : Phase2ReadoutMode);
100 
101  int adc_int = adc[i];
102 
103  if constexpr (mode < 0)
104  adc_int = int(adc_int * ElectronPerADCGain);
105  else {
106  if (adc_int < Phase2KinkADC)
107  adc_int = int((adc_int - 0.5) * ElectronPerADCGain);
108  else {
109  constexpr int8_t dspp = (Phase2ReadoutMode < 10 ? Phase2ReadoutMode : 10);
110  constexpr int8_t ds = int8_t(dspp <= 1 ? 1 : (dspp - 1) * (dspp - 1));
111 
112  adc_int -= (Phase2KinkADC - 1);
113  adc_int *= ds;
114  adc_int += (Phase2KinkADC - 1);
115 
116  adc_int = ((adc_int + 0.5 * ds) * ElectronPerADCGain);
117  }
118 
119  adc_int += int(Phase2DigiBaseline);
120  }
121 
122  adc[i] = std::min(adc_int, int(std::numeric_limits<uint16_t>::max()));
123  }
124  }
125 
126 } // namespace gpuCalibPixel
127 
128 #endif // RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
const dim3 threadIdx
Definition: cudaCompat.h:29
constexpr float VCaltoElectronOffset_L1
Definition: gpuCalibPixel.h:23
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
const dim3 gridDim
Definition: cudaCompat.h:33
ret
prodAgent to be discontinued
#define __global__
Definition: cudaCompat.h:19
constexpr float VCaltoElectronGain_L1
Definition: gpuCalibPixel.h:21
const dim3 blockDim
Definition: cudaCompat.h:30
constexpr float VCaltoElectronOffset
Definition: gpuCalibPixel.h:22
constexpr int VCalChargeThreshold
Definition: gpuCalibPixel.h:24
constexpr uint32_t numberOfModules
__global__ void calibDigis(uint16_t *id, uint16_t const *__restrict__ x, uint16_t const *__restrict__ y, uint16_t *adc, SiPixelGainForHLTonGPU const *__restrict__ ped, int numElements, uint32_t *__restrict__ moduleStart, uint32_t *__restrict__ nClustersInModule, uint32_t *__restrict__ clusModuleStart)
Definition: gpuCalibPixel.h:32
constexpr float ElectronPerADCGain
Definition: gpuCalibPixel.h:26
constexpr uint16_t Phase2DigiBaseline
Definition: gpuCalibPixel.h:28
const dim3 blockIdx
Definition: cudaCompat.h:32
constexpr uint16_t invalidModuleId
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ int32_t *__restrict__ uint32_t numElements
constexpr uint32_t numberOfModules
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
constexpr int8_t Phase2ReadoutMode
Definition: gpuCalibPixel.h:27
col
Definition: cuy.py:1009
float x
__global__ void calibDigisPhase2(uint16_t *id, uint16_t *adc, int numElements, uint32_t *__restrict__ moduleStart, uint32_t *__restrict__ nClustersInModule, uint32_t *__restrict__ clusModuleStart)
Definition: gpuCalibPixel.h:79
constexpr uint8_t Phase2KinkADC
Definition: gpuCalibPixel.h:29
constexpr float VCaltoElectronGain
Definition: gpuCalibPixel.h:20
uint16_t *__restrict__ uint16_t const *__restrict__ adc