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

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 30 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, phase1PixelTopology::numberOfModules, 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.

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

References cms::cudacompat::blockDim, cms::cudacompat::blockIdx, ElectronPerADCGain, first, cms::cudacompat::gridDim, mps_fire::i, gpuClustering::invalidModuleId, universalConfigTemplate::mode, phase2PixelTopology::numberOfModules, gpuClustering::numElements, Phase2DigiBaseline, Phase2KinkADC, Phase2ReadoutMode, and cms::cudacompat::threadIdx.

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

Variable Documentation

constexpr float gpuCalibPixel::ElectronPerADCGain = 600

Definition at line 24 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

constexpr uint16_t gpuCalibPixel::Phase2DigiBaseline = 1500

Definition at line 26 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

constexpr uint8_t gpuCalibPixel::Phase2KinkADC = 8

Definition at line 27 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

constexpr int8_t gpuCalibPixel::Phase2ReadoutMode = 3

Definition at line 25 of file gpuCalibPixel.h.

Referenced by calibDigisPhase2().

constexpr int gpuCalibPixel::VCalChargeThreshold = 100

Definition at line 22 of file gpuCalibPixel.h.

constexpr float gpuCalibPixel::VCaltoElectronGain = 47

Definition at line 18 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronGain_L1 = 50

Definition at line 19 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronOffset = -60

Definition at line 20 of file gpuCalibPixel.h.

Referenced by calibDigis().

constexpr float gpuCalibPixel::VCaltoElectronOffset_L1 = -670

Definition at line 21 of file gpuCalibPixel.h.

Referenced by calibDigis().