CMS 3D CMS Logo

Functions
gpuCalibPixel Namespace Reference

Functions

__global__ void calibDigis (SiPixelClusterThresholds clusterThresholds, 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)
 
__global__ void calibDigisPhase2 (SiPixelClusterThresholds clusterThresholds, uint16_t *id, uint16_t *adc, int numElements, uint32_t *__restrict__ moduleStart, uint32_t *__restrict__ nClustersInModule, uint32_t *__restrict__ clusModuleStart)
 

Function Documentation

◆ calibDigis()

__global__ void gpuCalibPixel::calibDigis ( SiPixelClusterThresholds  clusterThresholds,
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 at line 22 of file gpuCalibPixel.h.

References gpuClustering::adc, cms::cudacompat::blockDim, cms::cudacompat::blockIdx, cuy::col, dqmdumpme::first, nano_mu_digi_cff::float, PedestalClient_cfi::gain, cms::cudacompat::gridDim, mps_fire::i, gpuClustering::invalidModuleId, SiStripPI::max, gpuClustering::moduleStart, gpuClustering::nClustersInModule, phase1PixelTopology::numberOfModules, gpuClustering::numElements, hltrates_dqm_sourceclient-live_cfg::offset, EcalCondDBWriter_cfi::pedestal, runTheMatrix::ret, cms::cudacompat::threadIdx, SiPixelClusterThresholds::vCaltoElectronGain, HLT_2023v12_cff::VCaltoElectronGain, SiPixelClusterThresholds::vCaltoElectronGain_L1, HLT_2023v12_cff::VCaltoElectronGain_L1, SiPixelClusterThresholds::vCaltoElectronOffset, HLT_2023v12_cff::VCaltoElectronOffset, SiPixelClusterThresholds::vCaltoElectronOffset_L1, HLT_2023v12_cff::VCaltoElectronOffset_L1, and x.

32  {
33  int first = blockDim.x * blockIdx.x + threadIdx.x;
34 
35  const float VCaltoElectronGain = clusterThresholds.vCaltoElectronGain;
36  const float VCaltoElectronGain_L1 = clusterThresholds.vCaltoElectronGain_L1;
37  const float VCaltoElectronOffset = clusterThresholds.vCaltoElectronOffset;
38  const float VCaltoElectronOffset_L1 = clusterThresholds.vCaltoElectronOffset_L1;
39 
40  // zero for next kernels...
41  if (0 == first)
42  clusModuleStart[0] = moduleStart[0] = 0;
44  nClustersInModule[i] = 0;
45  }
46 
47  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
48  if (invalidModuleId == id[i])
49  continue;
50 
51  bool isDeadColumn = false, isNoisyColumn = false;
52 
53  int row = x[i];
54  int col = y[i];
55  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
56  float pedestal = ret.first;
57  float gain = ret.second;
58  // float pedestal = 0; float gain = 1.;
59  if (isDeadColumn | isNoisyColumn) {
60  printf("bad pixel at %d in %d\n", i, id[i]);
61  id[i] = invalidModuleId;
62  adc[i] = 0;
63  } else {
64  float vcal = float(adc[i]) * gain - pedestal * gain;
65 
66  float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
68  vcal = vcal * conversionFactor + offset;
69 
70  adc[i] = std::clamp(int(vcal), 100, int(std::numeric_limits<uint16_t>::max()));
71  }
72  }
73  }
const dim3 threadIdx
Definition: cudaCompat.h:29
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
const dim3 blockDim
Definition: cudaCompat.h:30
constexpr uint16_t numberOfModules
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
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
col
Definition: cuy.py:1009
float x
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ calibDigisPhase2()

__global__ void gpuCalibPixel::calibDigisPhase2 ( SiPixelClusterThresholds  clusterThresholds,
uint16_t *  id,
uint16_t *  adc,
int  numElements,
uint32_t *__restrict__  moduleStart,
uint32_t *__restrict__  nClustersInModule,
uint32_t *__restrict__  clusModuleStart 
)

Definition at line 75 of file gpuCalibPixel.h.

References gpuClustering::adc, cms::cudacompat::blockDim, cms::cudacompat::blockIdx, SiPixelClusterThresholds::electronPerADCGain, HLT_2023v12_cff::ElectronPerADCGain, dqmdumpme::first, cms::cudacompat::gridDim, mps_fire::i, createfilelist::int, gpuClustering::invalidModuleId, SiStripPI::max, SiStripPI::min, ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, gpuClustering::moduleStart, gpuClustering::nClustersInModule, phase2PixelTopology::numberOfModules, gpuClustering::numElements, SiPixelClusterThresholds::phase2DigiBaseline, HLT_2023v12_cff::Phase2DigiBaseline, SiPixelClusterThresholds::phase2KinkADC, HLT_2023v12_cff::Phase2KinkADC, SiPixelClusterThresholds::phase2ReadoutMode, HLT_2023v12_cff::Phase2ReadoutMode, and cms::cudacompat::threadIdx.

82  {
83  int first = blockDim.x * blockIdx.x + threadIdx.x;
84  // zero for next kernels...
85 
86  const float ElectronPerADCGain = clusterThresholds.electronPerADCGain;
87  const int8_t Phase2ReadoutMode = clusterThresholds.phase2ReadoutMode;
88  const uint16_t Phase2DigiBaseline = clusterThresholds.phase2DigiBaseline;
89  const uint8_t Phase2KinkADC = clusterThresholds.phase2KinkADC;
90 
91  if (0 == first)
92  clusModuleStart[0] = moduleStart[0] = 0;
94  nClustersInModule[i] = 0;
95  }
96 
97  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
98  if (invalidModuleId == id[i])
99  continue;
100 
101  const int mode = (Phase2ReadoutMode < -1 ? -1 : Phase2ReadoutMode);
102 
103  int adc_int = adc[i];
104 
105  if (mode < 0)
106  adc_int = int(adc_int * ElectronPerADCGain);
107  else {
108  if (adc_int < Phase2KinkADC)
109  adc_int = int((adc_int + 0.5) * ElectronPerADCGain);
110  else {
111  const int8_t dspp = (Phase2ReadoutMode < 10 ? Phase2ReadoutMode : 10);
112  const int8_t ds = int8_t(dspp <= 1 ? 1 : (dspp - 1) * (dspp - 1));
113 
114  adc_int -= Phase2KinkADC;
115  adc_int *= ds;
116  adc_int += Phase2KinkADC;
117 
118  adc_int = ((adc_int + 0.5 * ds) * ElectronPerADCGain);
119  }
120 
121  adc_int += int(Phase2DigiBaseline);
122  }
123  adc[i] = std::min(adc_int, int(std::numeric_limits<uint16_t>::max()));
124  }
125  }
const dim3 threadIdx
Definition: cudaCompat.h:29
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
const dim3 gridDim
Definition: cudaCompat.h:33
const dim3 blockDim
Definition: cudaCompat.h:30
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
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
constexpr uint16_t numberOfModules
uint16_t *__restrict__ uint16_t const *__restrict__ adc