CMS 3D CMS Logo

Functions | Variables
gpuCalibPixel Namespace Reference

Functions

template<bool isRun2>
__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)
 
__global__ void calibDigisPhase2 (uint16_t *id, uint16_t *adc, int numElements, uint32_t *__restrict__ moduleStart, uint32_t *__restrict__ nClustersInModule, uint32_t *__restrict__ clusModuleStart)
 

Variables

constexpr float ElectronPerADCGain = 600
 
constexpr uint16_t Phase2DigiBaseline = 1500
 
constexpr uint8_t Phase2KinkADC = 8
 
constexpr int8_t Phase2ReadoutMode = 3
 
constexpr int VCalChargeThreshold = 100
 
constexpr float VCaltoElectronGain = 47
 
constexpr float VCaltoElectronGain_L1 = 50
 
constexpr float VCaltoElectronOffset = -60
 
constexpr float VCaltoElectronOffset_L1 = -670
 

Function Documentation

◆ calibDigis()

template<bool isRun2>
__global__ void gpuCalibPixel::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 at line 32 of file gpuCalibPixel.h.

References gpuClustering::adc, cms::cudacompat::blockDim, cms::cudacompat::blockIdx, cuy::col, first, dqmMemoryStats::float, PedestalClient_cfi::gain, cms::cudacompat::gridDim, mps_fire::i, gpuClustering::invalidModuleId, geometryPPS_CMSxz_fromDD_2016_cfi::isRun2, SiStripPI::max, gpuClustering::moduleStart, gpuClustering::nClustersInModule, phase1PixelTopology::numberOfModules, gpuClustering::numElements, hltrates_dqm_sourceclient-live_cfg::offset, EcalCondDBWriter_cfi::pedestal, runTheMatrix::ret, cms::cudacompat::threadIdx, VCaltoElectronGain, VCaltoElectronGain_L1, VCaltoElectronOffset, VCaltoElectronOffset_L1, and x.

41  {
42  int first = blockDim.x * blockIdx.x + threadIdx.x;
43 
44  // zero for next kernels...
45  if (0 == first)
46  clusModuleStart[0] = moduleStart[0] = 0;
48  nClustersInModule[i] = 0;
49  }
50 
51  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
52  if (invalidModuleId == id[i])
53  continue;
54 
55  bool isDeadColumn = false, isNoisyColumn = false;
56 
57  int row = x[i];
58  int col = y[i];
59  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
60  float pedestal = ret.first;
61  float gain = ret.second;
62  // float pedestal = 0; float gain = 1.;
63  if (isDeadColumn | isNoisyColumn) {
64  printf("bad pixel at %d in %d\n", i, id[i]);
65  id[i] = invalidModuleId;
66  adc[i] = 0;
67  } else {
68  float vcal = float(adc[i]) * gain - pedestal * gain;
69  if constexpr (isRun2) {
70  float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
72  vcal = vcal * conversionFactor + offset;
73  }
74  adc[i] = std::clamp(int(vcal), 100, int(std::numeric_limits<uint16_t>::max()));
75  }
76  }
77  }
const dim3 threadIdx
Definition: cudaCompat.h:29
constexpr float VCaltoElectronOffset_L1
Definition: gpuCalibPixel.h:23
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
constexpr float VCaltoElectronGain_L1
Definition: gpuCalibPixel.h:21
const dim3 blockDim
Definition: cudaCompat.h:30
constexpr float VCaltoElectronOffset
Definition: gpuCalibPixel.h:22
constexpr uint32_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
constexpr float VCaltoElectronGain
Definition: gpuCalibPixel.h:20
uint16_t *__restrict__ uint16_t const *__restrict__ adc

◆ calibDigisPhase2()

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

Definition at line 79 of file gpuCalibPixel.h.

References gpuClustering::adc, cms::cudacompat::blockDim, cms::cudacompat::blockIdx, ElectronPerADCGain, 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, Phase2DigiBaseline, Phase2KinkADC, Phase2ReadoutMode, and cms::cudacompat::threadIdx.

85  {
86  int first = blockDim.x * blockIdx.x + threadIdx.x;
87  // zero for next kernels...
88 
89  if (0 == first)
90  clusModuleStart[0] = moduleStart[0] = 0;
92  nClustersInModule[i] = 0;
93  }
94 
95  for (int i = first; i < numElements; i += gridDim.x * blockDim.x) {
96  if (invalidModuleId == id[i])
97  continue;
98 
99  constexpr int mode = (Phase2ReadoutMode < -1 ? -1 : Phase2ReadoutMode);
100 
101  int adc_int = adc[i];
102 
103  if constexpr (mode < 0)
104  adc_int = int(adc_int * ElectronPerADCGain);
105  else {
106  if (adc_int < Phase2KinkADC)
107  adc_int = int((adc_int - 0.5) * ElectronPerADCGain);
108  else {
109  constexpr int8_t dspp = (Phase2ReadoutMode < 10 ? Phase2ReadoutMode : 10);
110  constexpr int8_t ds = int8_t(dspp <= 1 ? 1 : (dspp - 1) * (dspp - 1));
111 
112  adc_int -= (Phase2KinkADC - 1);
113  adc_int *= ds;
114  adc_int += (Phase2KinkADC - 1);
115 
116  adc_int = ((adc_int + 0.5 * ds) * ElectronPerADCGain);
117  }
118 
119  adc_int += int(Phase2DigiBaseline);
120  }
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
constexpr float ElectronPerADCGain
Definition: gpuCalibPixel.h:26
constexpr uint16_t Phase2DigiBaseline
Definition: gpuCalibPixel.h:28
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
constexpr uint32_t numberOfModules
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
constexpr int8_t Phase2ReadoutMode
Definition: gpuCalibPixel.h:27
constexpr uint8_t Phase2KinkADC
Definition: gpuCalibPixel.h:29
uint16_t *__restrict__ uint16_t const *__restrict__ adc

Variable Documentation

◆ ElectronPerADCGain

constexpr float gpuCalibPixel::ElectronPerADCGain = 600

Definition at line 26 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

◆ Phase2DigiBaseline

constexpr uint16_t gpuCalibPixel::Phase2DigiBaseline = 1500

Definition at line 28 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

◆ Phase2KinkADC

constexpr uint8_t gpuCalibPixel::Phase2KinkADC = 8

Definition at line 29 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

◆ Phase2ReadoutMode

constexpr int8_t gpuCalibPixel::Phase2ReadoutMode = 3

Definition at line 27 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

◆ VCalChargeThreshold

constexpr int gpuCalibPixel::VCalChargeThreshold = 100

Definition at line 24 of file gpuCalibPixel.h.

◆ VCaltoElectronGain

constexpr float gpuCalibPixel::VCaltoElectronGain = 47

Definition at line 20 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronGain_L1

constexpr float gpuCalibPixel::VCaltoElectronGain_L1 = 50

Definition at line 21 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronOffset

constexpr float gpuCalibPixel::VCaltoElectronOffset = -60

Definition at line 22 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronOffset_L1

constexpr float gpuCalibPixel::VCaltoElectronOffset_L1 = -670

Definition at line 23 of file gpuCalibPixel.h.

Referenced by calibDigis().