CMS 3D CMS Logo

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