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 
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 
22  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  float conversionFactor = (isRun2) ? (id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain) : 1.f;
46  float offset = (isRun2) ? (id[i] < 96 ? VCaltoElectronOffset_L1 : VCaltoElectronOffset) : 0;
47 
48  bool isDeadColumn = false, isNoisyColumn = false;
49 
50  int row = x[i];
51  int col = y[i];
52  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
53  float pedestal = ret.first;
54  float gain = ret.second;
55  // float pedestal = 0; float gain = 1.;
56  if (isDeadColumn | isNoisyColumn) {
57  id[i] = invalidModuleId;
58  adc[i] = 0;
59  printf("bad pixel at %d in %d\n", i, id[i]);
60  } else {
61  float vcal = adc[i] * gain - pedestal * gain;
62  adc[i] = std::max(100, int(vcal * conversionFactor + offset));
63  }
64  }
65  }
66 } // namespace gpuCalibPixel
67 
68 #endif // RecoLocalTracker_SiPixelClusterizer_plugins_gpuCalibPixel_h
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:542
SiPixelGainForHLTonGPU
Definition: SiPixelGainForHLTonGPU.h:28
mps_fire.i
i
Definition: mps_fire.py:428
SiPixelGainForHLTonGPU.h
gpuCalibPixel::VCaltoElectronGain
constexpr float VCaltoElectronGain
Definition: gpuCalibPixel.h:16
gpuClustering::adc
uint16_t *__restrict__ uint16_t const *__restrict__ adc
Definition: gpuClusterChargeCut.h:20
gpuCalibPixel::VCaltoElectronGain_L1
constexpr float VCaltoElectronGain_L1
Definition: gpuCalibPixel.h:17
cuy.col
col
Definition: cuy.py:1010
gpuCalibPixel::calibDigis
__global__ void calibDigis(bool isRun2, 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:21
gpuClustering::numElements
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ int32_t *__restrict__ uint32_t numElements
Definition: gpuClusterChargeCut.h:26
__global__
#define __global__
Definition: cudaCompat.h:19
gpuClustering::moduleStart
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
Definition: gpuClusterChargeCut.h:20
gpuCalibPixel::VCaltoElectronOffset
constexpr float VCaltoElectronOffset
Definition: gpuCalibPixel.h:18
geometryPPS_CMSxz_fromDD_2016_cfi.isRun2
isRun2
Definition: geometryPPS_CMSxz_fromDD_2016_cfi.py:14
gpuClustering::maxNumModules
constexpr uint16_t maxNumModules
Definition: gpuClusteringConstants.h:29
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:112
cms::cudacompat::gridDim
const dim3 gridDim
Definition: cudaCompat.h:33
cms::cudacompat::blockDim
const dim3 blockDim
Definition: cudaCompat.h:30
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
gpuClustering::invalidModuleId
constexpr uint16_t invalidModuleId
Definition: gpuClusteringConstants.h:32
gpuCalibPixel
Definition: gpuCalibPixel.h:11
gpuClusteringConstants.h
EcalCondDBWriter_cfi.pedestal
pedestal
Definition: EcalCondDBWriter_cfi.py:49
cms::cudacompat::threadIdx
const dim3 threadIdx
Definition: cudaCompat.h:29
PedestalClient_cfi.gain
gain
Definition: PedestalClient_cfi.py:37
gpuClustering::nClustersInModule
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
Definition: gpuClusterChargeCut.h:20
gpuCalibPixel::VCaltoElectronOffset_L1
constexpr float VCaltoElectronOffset_L1
Definition: gpuCalibPixel.h:19
cuda_assert.h
hltrates_dqm_sourceclient-live_cfg.offset
offset
Definition: hltrates_dqm_sourceclient-live_cfg.py:82
cms::cudacompat::blockIdx
const dim3 blockIdx
Definition: cudaCompat.h:32