CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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)
 

Variables

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

Function Documentation

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 22 of file gpuCalibPixel.h.

References cms::cudacompat::blockDim, cms::cudacompat::blockIdx, cuy::col, first, cms::cudacompat::gridDim, mps_fire::i, gpuClustering::invalidModuleId, HLT_FULL_cff::isRun2, SiStripPI::max, gpuClustering::maxNumModules, gpuClustering::numElements, hltrates_dqm_sourceclient-live_cfg::offset, EcalCondDBWriter_cfi::pedestal, gpuVertexFinder::printf(), runTheMatrix::ret, cms::cudacompat::threadIdx, VCaltoElectronGain, VCaltoElectronGain_L1, VCaltoElectronOffset, and VCaltoElectronOffset_L1.

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  bool isDeadColumn = false, isNoisyColumn = false;
46 
47  int row = x[i];
48  int col = y[i];
49  auto ret = ped->getPedAndGain(id[i], col, row, isDeadColumn, isNoisyColumn);
50  float pedestal = ret.first;
51  float gain = ret.second;
52  // float pedestal = 0; float gain = 1.;
53  if (isDeadColumn | isNoisyColumn) {
54  printf("bad pixel at %d in %d\n", i, id[i]);
55  id[i] = invalidModuleId;
56  adc[i] = 0;
57  } else {
58  float vcal = float(adc[i]) * gain - pedestal * gain;
59  if constexpr (isRun2) {
60  float conversionFactor = id[i] < 96 ? VCaltoElectronGain_L1 : VCaltoElectronGain;
62  vcal = vcal * conversionFactor + offset;
63  }
64  adc[i] = std::max(100, int(vcal));
65  }
66  }
67  }
const dim3 threadIdx
Definition: cudaCompat.h:29
tuple ret
prodAgent to be discontinued
tuple VCaltoElectronOffset
const dim3 gridDim
Definition: cudaCompat.h:33
const dim3 blockDim
Definition: cudaCompat.h:30
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
uint16_t const *__restrict__ x
Definition: gpuClustering.h:39
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
constexpr uint16_t maxNumModules
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ int32_t *__restrict__ uint32_t numElements
const dim3 blockIdx
Definition: cudaCompat.h:32
constexpr uint16_t invalidModuleId
tuple VCaltoElectronGain
tuple VCaltoElectronOffset_L1
tuple VCaltoElectronGain_L1
int col
Definition: cuy.py:1009
uint16_t *__restrict__ uint16_t const *__restrict__ adc

Variable Documentation

constexpr float gpuCalibPixel::VCaltoElectronGain = 47

Definition at line 16 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronGain_L1 = 50

Definition at line 17 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronOffset = -60

Definition at line 18 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronOffset_L1 = -670

Definition at line 19 of file gpuCalibPixel.h.

Referenced by calibDigis().