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 
10 
11 namespace gpuCalibPixel {
12 
14 
15  // valid for run2
16  constexpr float VCaltoElectronGain = 47; // L2-4: 47 +- 4.7
17  constexpr float VCaltoElectronGain_L1 = 50; // L1: 49.6 +- 2.6
18  constexpr float VCaltoElectronOffset = -60; // L2-4: -60 +- 130
19  constexpr float VCaltoElectronOffset_L1 = -670; // L1: -670 +- 220
20 
21  template <bool isRun2>
22  __global__ void calibDigis(uint16_t* id,
23  uint16_t const* __restrict__ x,
24  uint16_t const* __restrict__ y,
25  uint16_t* adc,
26  SiPixelGainForHLTonGPU const* __restrict__ ped,
27  int numElements,
28  uint32_t* __restrict__ moduleStart, // just to zero first
29  uint32_t* __restrict__ nClustersInModule, // just to zero them
30  uint32_t* __restrict__ clusModuleStart // just to zero first
31  ) {
32  int first = blockDim.x * blockIdx.x + threadIdx.x;
33 
34  // zero for next kernels...
35  if (0 == first)
36  clusModuleStart[0] = moduleStart[0] = 0;
37  for (int i = first; i < gpuClustering::maxNumModules; i += gridDim.x * blockDim.x) {
38  nClustersInModule[i] = 0;
39  }
40 
41  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
42  if (invalidModuleId == id[i])
43  continue;
44 
45  bool isDeadColumn = false, isNoisyColumn = false;
46 
47  int row = x[i];
48  int col = y[i];
49  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
50  float pedestal = ret.first;
51  float gain = ret.second;
52  // float pedestal = 0; float gain = 1.;
53  if (isDeadColumn | isNoisyColumn) {
54  printf("bad pixel at %d in %d\n", i, id[i]);
55  id[i] = invalidModuleId;
56  adc[i] = 0;
57  } else {
58  float vcal = float(adc[i]) * gain - pedestal * gain;
59  if constexpr (isRun2) {
60  float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
62  vcal = vcal * conversionFactor + offset;
63  }
64  adc[i] = std::max(100, int(vcal));
65  }
66  }
67  }
68 } // namespace gpuCalibPixel
69 
70 #endif // RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
const dim3 threadIdx
Definition: cudaCompat.h:29
constexpr float VCaltoElectronOffset_L1
Definition: gpuCalibPixel.h:19
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:17
const dim3 blockDim
Definition: cudaCompat.h:30
constexpr float VCaltoElectronOffset
Definition: gpuCalibPixel.h:18
__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:22
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)
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
constexpr uint16_t maxNumModules
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
int col
Definition: cuy.py:1009
constexpr float VCaltoElectronGain
Definition: gpuCalibPixel.h:16
uint16_t *__restrict__ uint16_t const *__restrict__ adc