CMS 3D CMS Logo

gpuClusterChargeCut.h
Go to the documentation of this file.
1 #ifndef RecoLocalTracker_SiPixelClusterizer_plugins_gpuClusterChargeCut_h
2 #define RecoLocalTracker_SiPixelClusterizer_plugins_gpuClusterChargeCut_h
3 
4 #include <cstdint>
5 #include <cstdio>
6 
11 
12 // local include(s)
14 
15 namespace gpuClustering {
16 
17  template <bool isPhase2>
20  clusterThresholds, // charge cut on cluster in electrons (for layer 1 and for other layers)
21  uint16_t* __restrict__ id, // module id of each pixel (modified if bad cluster)
22  uint16_t const* __restrict__ adc, // charge of each pixel
23  uint32_t const* __restrict__ moduleStart, // index of the first pixel of each module
24  uint32_t* __restrict__ nClustersInModule, // modified: number of clusters found in each module
25  uint32_t const* __restrict__ moduleId, // module id of each module
26  int32_t* __restrict__ clusterId, // modified: cluster id of each pixel
27  uint32_t numElements) {
28  __shared__ int32_t charge[maxNumClustersPerModules];
29  __shared__ uint8_t ok[maxNumClustersPerModules];
30  __shared__ uint16_t newclusId[maxNumClustersPerModules];
31 
33  [[maybe_unused]] constexpr int nMaxModules =
35 
38 
41  for (auto module = firstModule; module < endModule; module += gridDim.x) {
42  auto firstPixel = moduleStart[1 + module];
43  auto thisModuleId = id[firstPixel];
44  while (thisModuleId == invalidModuleId and firstPixel < numElements) {
45  // skip invalid or duplicate pixels
46  ++firstPixel;
47  thisModuleId = id[firstPixel];
48  }
49  if (firstPixel >= numElements) {
50  // reached the end of the input while skipping the invalid pixels, nothing left to do
51  break;
52  }
53  if (thisModuleId != moduleId[module]) {
54  // reached the end of the module while skipping the invalid pixels, skip this module
55  continue;
56  }
57  assert(thisModuleId < nMaxModules);
58 
59  auto nclus = nClustersInModule[thisModuleId];
60  if (nclus == 0)
61  continue;
62 
63  if (threadIdx.x == 0 && nclus > maxNumClustersPerModules)
64  printf("Warning too many clusters in module %d in block %d: %d > %d\n",
65  thisModuleId,
66  blockIdx.x,
67  nclus,
69 
70  auto first = firstPixel + threadIdx.x;
71 
72  if (nclus > maxNumClustersPerModules) {
73  // remove excess FIXME find a way to cut charge first....
74  for (auto i = first; i < numElements; i += blockDim.x) {
75  if (id[i] == invalidModuleId)
76  continue; // not valid
77  if (id[i] != thisModuleId)
78  break; // end of module
80  id[i] = invalidModuleId;
82  }
83  }
85  }
86 
87 #ifdef GPU_DEBUG
88  if (thisModuleId % 100 == 1)
89  if (threadIdx.x == 0)
90  printf("start cluster charge cut for module %d in block %d\n", thisModuleId, blockIdx.x);
91 #endif
92 
94  for (auto i = threadIdx.x; i < nclus; i += blockDim.x) {
95  charge[i] = 0;
96  }
97  __syncthreads();
98 
99  for (auto i = first; i < numElements; i += blockDim.x) {
100  if (id[i] == invalidModuleId)
101  continue; // not valid
102  if (id[i] != thisModuleId)
103  break; // end of module
105  }
106  __syncthreads();
107 
108  auto chargeCut = clusterThresholds.getThresholdForLayerOnCondition(thisModuleId < startBPIX2);
109 
110  bool good = true;
111  for (auto i = threadIdx.x; i < nclus; i += blockDim.x) {
112  newclusId[i] = ok[i] = charge[i] >= chargeCut ? 1 : 0;
113  if (0 == ok[i])
114  good = false;
115  }
116 
117  // if all clusters above threshold do nothing
118  if (__syncthreads_and(good))
119  continue;
120 
121  // renumber
122  __shared__ uint16_t ws[32];
123  cms::cuda::blockPrefixScan(newclusId, nclus, ws);
124 
125  assert(nclus > newclusId[nclus - 1]);
126 
127  nClustersInModule[thisModuleId] = newclusId[nclus - 1];
128 
129  // reassign id
130  for (auto i = first; i < numElements; i += blockDim.x) {
131  if (id[i] == invalidModuleId)
132  continue; // not valid
133  if (id[i] != thisModuleId)
134  break; // end of module
135  if (0 == ok[clusterId[i]])
136  clusterId[i] = id[i] = invalidModuleId;
137  else
138  clusterId[i] = newclusId[clusterId[i]] - 1;
139  }
140 
141  //done
142  __syncthreads();
143  } // loop on modules
144  }
145 
146 } // namespace gpuClustering
147 
148 #endif // RecoLocalTracker_SiPixelClusterizer_plugins_gpuClusterChargeCut_h
const dim3 threadIdx
Definition: cudaCompat.h:29
__shared__ uint8_t ok[maxNumClustersPerModules]
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
constexpr int nMaxModules
const dim3 gridDim
Definition: cudaCompat.h:33
constexpr uint32_t layerStart[numberOfLayers+1]
#define __global__
Definition: cudaCompat.h:19
const dim3 blockDim
Definition: cudaCompat.h:30
float clusterChargeCut(const edm::ParameterSet &conf, const char *name="clusterChargeCut")
constexpr uint32_t numberOfModules
constexpr int32_t getThresholdForLayerOnCondition(bool isLayer1) const noexcept
auto &__restrict__ ws
assert(nMaxModules< maxNumModules)
constexpr int startBPIX2
constexpr int32_t maxNumClustersPerModules
constexpr uint16_t maxNumModules
const dim3 blockIdx
Definition: cudaCompat.h:32
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ int32_t *__restrict__ clusterId
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
auto const good
min quality of good
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ moduleStart
__shared__ uint16_t newclusId[maxNumClustersPerModules]
void __syncthreads()
Definition: cudaCompat.h:132
constexpr uint32_t layerStart[numberOfLayers+1]
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ moduleId
bool __syncthreads_and(bool x)
Definition: cudaCompat.h:135
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61
uint16_t *__restrict__ uint16_t const *__restrict__ adc