1 #ifndef RecoLocalTracker_SiPixelClusterizer_plugins_gpuClustering_h
2 #define RecoLocalTracker_SiPixelClusterizer_plugins_gpuClustering_h
12 namespace gpuClustering {
18 __global__ void countModules(uint16_t
const* __restrict__
id,
30 if (
j < 0
or id[
j] !=
id[i]) {
33 moduleStart[loc + 1] =
i;
38 __global__ void findClus(uint16_t
const* __restrict__
id,
39 uint16_t
const* __restrict__
x,
40 uint16_t
const* __restrict__
y,
41 uint32_t
const* __restrict__ moduleStart,
44 int32_t* __restrict__ clusterId,
51 auto firstPixel = moduleStart[1 +
module];
52 auto thisModuleId =
id[firstPixel];
56 if (thisModuleId % 100 == 1)
58 printf(
"start clusterizer for module %d in block %d\n", thisModuleId,
blockIdx.x);
71 if (
id[i] != thisModuleId) {
79 constexpr uint32_t maxPixInModule = 6000;
89 assert((msize == numElements)
or ((msize < numElements) and (
id[msize] != thisModuleId)));
93 if (msize - firstPixel > maxPixInModule) {
94 printf(
"too many pixels in module %d: %d > %d\n", thisModuleId, msize - firstPixel, maxPixInModule);
95 msize = maxPixInModule + firstPixel;
100 assert(msize - firstPixel <= maxPixInModule);
103 __shared__ uint32_t totGood;
109 for (
int i = first; i < msize; i +=
blockDim.x) {
124 assert(hist.size() == totGood);
125 if (thisModuleId % 100 == 1)
127 printf(
"histo size %d\n", hist.size());
129 for (
int i = first; i < msize; i +=
blockDim.x) {
132 hist.fill(y[i], i - firstPixel);
139 printf(
"THIS IS NOT SUPPOSED TO HAPPEN too many hits in module %d: %d for block size %d\n",
144 auto maxiter = hist.size();
147 constexpr
int maxNeighbours = 10;
159 __shared__ uint32_t n40, n60;
163 if (hist.size(
j) > 60)
165 if (hist.size(
j) > 40)
171 printf(
"columns with more than 60 px %d in %d\n", n60, thisModuleId);
173 printf(
"columns with more than 40 px %d in %d\n", n40, thisModuleId);
181 auto p = hist.begin() +
j;
182 auto i = *
p + firstPixel;
184 assert(
id[i] == thisModuleId);
186 auto e = hist.end(be);
190 auto m = (*p) + firstPixel;
192 assert(
int(y[
m]) -
int(y[i]) >= 0);
193 assert(
int(y[m]) -
int(y[i]) <= 1);
194 if (
std::abs(
int(x[m]) -
int(x[i])) > 1)
209 if (1 == nloops % 2) {
211 auto p = hist.begin() +
j;
212 auto i = *
p + firstPixel;
213 auto m = clusterId[
i];
214 while (
m != clusterId[
m])
221 auto p = hist.begin() +
j;
222 auto i = *
p + firstPixel;
223 for (
int kk = 0;
kk < nnn[
k]; ++
kk) {
225 auto m =
l + firstPixel;
229 if (old != clusterId[i]) {
248 if (thisModuleId % 100 == 1)
250 printf(
"# loops %d\n", nloops);
260 for (
int i = first; i < msize; i +=
blockDim.x) {
263 if (clusterId[i] == i) {
264 auto old =
atomicInc(&foundClusters, 0xffffffff);
265 clusterId[
i] = -(old + 1);
271 for (
int i = first; i < msize; i +=
blockDim.x) {
274 if (clusterId[i] >= 0) {
276 clusterId[
i] = clusterId[clusterId[
i]];
282 for (
int i = first; i < msize; i +=
blockDim.x) {
287 clusterId[
i] = -clusterId[
i] - 1;
293 moduleId[
module] = thisModuleId;
295 if (foundClusters > gMaxHit) {
297 if (foundClusters > 8)
298 printf(
"max hit %d in %d\n", foundClusters, thisModuleId);
302 if (thisModuleId % 100 == 1)
303 printf(
"%d clusters in module %d\n", foundClusters, thisModuleId);
310 #endif // RecoLocalTracker_SiPixelClusterizer_plugins_gpuClustering_h
__shared__ uint8_t ok[maxNumClustersPerModules]
bool __syncthreads_or(bool x)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
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 numColsInModule
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ nClustersInModule
std::function< unsigned int(align::ID)> Counter
cms::cuda::HistoContainer< uint8_t, 256, 16000, 8, uint16_t > Hist
T1 atomicInc(T1 *a, T2 b)
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
Abs< T >::type abs(const T &t)
uint16_t const *__restrict__ x
static const MaxIter maxiter
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
constexpr uint16_t invalidModuleId
uint16_t *__restrict__ uint16_t const *__restrict__ uint32_t const *__restrict__ uint32_t *__restrict__ uint32_t const *__restrict__ moduleId
T1 atomicMin_block(T1 *a, T2 b)
constexpr int invalidClusterId
bool __syncthreads_and(bool x)
T1 atomicMin(T1 *a, T2 b)
uint16_t const *__restrict__ uint16_t const *__restrict__ y
__shared__ unsigned int foundClusters
T1 atomicAdd(T1 *a, T2 b)