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 20 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, QIE10Task_cfi::ped, EcalCondDBWriter_cfi::pedestal, runTheMatrix::ret, cms::cudacompat::threadIdx, SiPixelClusterThresholds::vCaltoElectronGain, HLT_2024v14_cff::VCaltoElectronGain, SiPixelClusterThresholds::vCaltoElectronGain_L1, HLT_2024v14_cff::VCaltoElectronGain_L1, SiPixelClusterThresholds::vCaltoElectronOffset, HLT_2024v14_cff::VCaltoElectronOffset, SiPixelClusterThresholds::vCaltoElectronOffset_L1, HLT_2024v14_cff::VCaltoElectronOffset_L1, and x.

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

References gpuClustering::adc, cms::cudacompat::blockDim, cms::cudacompat::blockIdx, SiPixelClusterThresholds::electronPerADCGain, HLT_2024v14_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_2024v14_cff::Phase2DigiBaseline, SiPixelClusterThresholds::phase2KinkADC, HLT_2024v14_cff::Phase2KinkADC, SiPixelClusterThresholds::phase2ReadoutMode, HLT_2024v14_cff::Phase2ReadoutMode, and cms::cudacompat::threadIdx.

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