CMS 3D CMS Logo

Functions | Variables
gpuCalibPixel Namespace Reference

Functions

__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)
 

Variables

constexpr float VCaltoElectronGain = 47
 
constexpr float VCaltoElectronGain_L1 = 50
 
constexpr float VCaltoElectronOffset = -60
 
constexpr float VCaltoElectronOffset_L1 = -670
 

Function Documentation

◆ calibDigis()

__global__ void gpuCalibPixel::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 at line 21 of file gpuCalibPixel.h.

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  printf("bad pixel at %d in %d\n", i, id[i]);
58  id[i] = invalidModuleId;
59  adc[i] = 0;
60  } else {
61  float vcal = adc[i] * gain - pedestal * gain;
62  adc[i] = std::max(100, int(vcal * conversionFactor + offset));
63  }
64  }
65  }

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

Variable Documentation

◆ VCaltoElectronGain

constexpr float gpuCalibPixel::VCaltoElectronGain = 47
constexpr

Definition at line 16 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronGain_L1

constexpr float gpuCalibPixel::VCaltoElectronGain_L1 = 50
constexpr

Definition at line 17 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronOffset

constexpr float gpuCalibPixel::VCaltoElectronOffset = -60
constexpr

Definition at line 18 of file gpuCalibPixel.h.

Referenced by calibDigis().

◆ VCaltoElectronOffset_L1

constexpr float gpuCalibPixel::VCaltoElectronOffset_L1 = -670
constexpr

Definition at line 19 of file gpuCalibPixel.h.

Referenced by calibDigis().

runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:543
mps_fire.i
i
Definition: mps_fire.py:428
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:1009
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
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:18
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:125
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:20
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
genVertex_cff.x
x
Definition: genVertex_cff.py:13
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
gpuCalibPixel::VCaltoElectronOffset_L1
constexpr float VCaltoElectronOffset_L1
Definition: gpuCalibPixel.h:19
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