forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathgpuPixelRecHits.h
240 lines (202 loc) · 10.7 KB
/
gpuPixelRecHits.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
#ifndef RecoLocalTracker_SiPixelRecHits_plugins_alpaka_gpuPixelRecHits_h
#define RecoLocalTracker_SiPixelRecHits_plugins_alpaka_gpuPixelRecHits_h
#include <algorithm>
#include <cstdint>
#include <cstdio>
#include <limits>
#include <alpaka/alpaka.hpp>
#include "HeterogeneousCore/AlpakaInterface/interface/config.h"
#include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
#include "DataFormats/BeamSpotAlpaka/interface/alpaka/BeamSpotAlpaka.h"
#include "DataFormats/Math/interface/approx_atan2.h"
#include "DataFormats/SiPixelClusterSoA/interface/gpuClusteringConstants.h"
#include "DataFormats/SiPixelDigiSoA/interface/alpaka/SiPixelDigisDevice.h"
#include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsLayout.h"
#include "DataFormats/PixelCPEFastParams/interface/PixelCPEFastParams.h"
#include "DataFormats/TrackerCommon/interface/SimplePixelTopology.h"
#include "RecoLocalTracker/SiPixelRecHits/interface/pixelCPEforDevice.h"
//#define GPU_DEBUG 1
namespace ALPAKA_ACCELERATOR_NAMESPACE {
namespace gpuPixelRecHits {
template <typename TrackerTraits>
class getHits {
public:
template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
ALPAKA_FN_ACC void operator()(const TAcc& acc,
pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* __restrict__ cpeParams,
BeamSpotPOD const* __restrict__ bs,
SiPixelDigisLayoutSoAConstView digis,
uint32_t numElements,
SiPixelClustersLayoutSoAConstView clusters,
TrackingRecHitAlpakaSoAView<TrackerTraits> hits) const {
// FIXME
// the compiler seems NOT to optimize loads from views (even in a simple test case)
// The whole gimnastic here of copying or not is a pure heuristic exercise that seems to produce the fastest code with the above signature
// not using views (passing a gazzilion of array pointers) seems to produce the fastest code (but it is harder to mantain)
ALPAKA_ASSERT_OFFLOAD(cpeParams);
const uint32_t blockIdx(alpaka::getIdx<alpaka::Grid, alpaka::Blocks>(acc)[0u]);
const uint32_t threadIdxLocal(alpaka::getIdx<alpaka::Block, alpaka::Threads>(acc)[0u]);
// copy average geometry corrected by beamspot . FIXME (move it somewhere else???)
if (0 == blockIdx) {
auto& agc = hits.averageGeometry();
auto const& ag = cpeParams->averageGeometry();
auto nLadders = TrackerTraits::numberOfLaddersInBarrel;
cms::alpakatools::for_each_element_in_block_strided(acc, nLadders, [&](uint32_t il) {
agc.ladderZ[il] = ag.ladderZ[il] - bs->z;
agc.ladderX[il] = ag.ladderX[il] - bs->x;
agc.ladderY[il] = ag.ladderY[il] - bs->y;
agc.ladderR[il] = sqrt(agc.ladderX[il] * agc.ladderX[il] + agc.ladderY[il] * agc.ladderY[il]);
agc.ladderMinZ[il] = ag.ladderMinZ[il] - bs->z;
agc.ladderMaxZ[il] = ag.ladderMaxZ[il] - bs->z;
});
if (0 == threadIdxLocal) {
agc.endCapZ[0] = ag.endCapZ[0] - bs->z;
agc.endCapZ[1] = ag.endCapZ[1] - bs->z;
}
}
// to be moved in common namespace...
using gpuClustering::invalidModuleId;
constexpr int32_t MaxHitsInIter = pixelCPEforDevice::MaxHitsInIter;
using ClusParams = pixelCPEforDevice::ClusParams;
// as usual one block per module
auto& clusParams = alpaka::declareSharedVar<ClusParams, __COUNTER__>(acc);
auto me = clusters[blockIdx].moduleId();
int nclus = clusters[me].clusInModule();
if (0 == nclus)
return;
#ifdef GPU_DEBUG
if (threadIdx == 0) {
auto k = clusters[1 + blockIdx.x].moduleStart();
while (digis[k].moduleId() == invalidModuleId)
++k;
ALPAKA_ASSERT_OFFLOAD(digis[k].moduleId() == me);
}
if (me % 100 == 1)
if (threadIdx == 0)
printf(
"hitbuilder: %d clusters in module %d. will write at %d\n", nclus, me, clusters[me].clusModuleStart());
#endif
for (int startClus = 0, endClus = nclus; startClus < endClus; startClus += MaxHitsInIter) {
auto first = clusters[me].clusModuleStart() + startClus;
int nClusInIter = alpaka::math::min(acc, MaxHitsInIter, endClus - startClus);
int lastClus = startClus + nClusInIter;
assert(nClusInIter <= nclus);
assert(nClusInIter > 0);
assert(lastClus <= nclus);
assert(nclus > MaxHitsInIter || (0 == startClus && nClusInIter == nclus && lastClus == nclus));
// init
cms::alpakatools::for_each_element_in_block_strided(acc, nClusInIter, [&](uint32_t ic) {
clusParams.minRow[ic] = std::numeric_limits<uint32_t>::max();
clusParams.maxRow[ic] = 0;
clusParams.minCol[ic] = std::numeric_limits<uint32_t>::max();
clusParams.maxCol[ic] = 0;
clusParams.charge[ic] = 0;
clusParams.q_f_X[ic] = 0;
clusParams.q_l_X[ic] = 0;
clusParams.q_f_Y[ic] = 0;
clusParams.q_l_Y[ic] = 0;
});
alpaka::syncBlockThreads(acc);
// one thread per "digi"
const uint32_t blockDimension(alpaka::getWorkDiv<alpaka::Block, alpaka::Elems>(acc)[0u]);
const auto& [firstElementIdxNoStride, endElementIdxNoStride] =
cms::alpakatools::element_index_range_in_block(acc, first);
uint32_t rowsColsFirstElementIdx = firstElementIdxNoStride;
uint32_t rowsColsEndElementIdx = endElementIdxNoStride;
for (uint32_t i = rowsColsFirstElementIdx; i < numElements; ++i) {
if (not cms::alpakatools::next_valid_element_index_strided(
i, rowsColsFirstElementIdx, rowsColsEndElementIdx, blockDimension, numElements))
break;
auto id = digis[i].moduleId();
if (id == invalidModuleId)
continue; // not valid
if (id != me)
break; // end of module
auto cl = digis[i].clus();
if (cl < startClus || cl >= lastClus)
continue;
cl -= startClus;
ALPAKA_ASSERT_OFFLOAD(cl >= 0);
ALPAKA_ASSERT_OFFLOAD(cl < MaxHitsInIter);
auto x = digis[i].xx();
auto y = digis[i].yy();
alpaka::atomicMin(acc, &clusParams.minRow[cl], (uint32_t)x, alpaka::hierarchy::Threads{});
alpaka::atomicMax(acc, &clusParams.maxRow[cl], (uint32_t)x, alpaka::hierarchy::Threads{});
alpaka::atomicMin(acc, &clusParams.minCol[cl], (uint32_t)y, alpaka::hierarchy::Threads{});
alpaka::atomicMax(acc, &clusParams.maxCol[cl], (uint32_t)y, alpaka::hierarchy::Threads{});
}
alpaka::syncBlockThreads(acc);
auto pixmx = cpeParams->detParams(me).pixmx;
uint32_t chargeFirstElementIdx = firstElementIdxNoStride;
uint32_t chargeEndElementIdx = endElementIdxNoStride;
for (uint32_t i = chargeFirstElementIdx; i < numElements; ++i) {
if (not cms::alpakatools::next_valid_element_index_strided(
i, chargeFirstElementIdx, chargeEndElementIdx, blockDimension, numElements))
break;
auto id = digis[i].moduleId();
if (id == invalidModuleId)
continue; // not valid
if (id != me)
break; // end of module
auto cl = digis[i].clus();
if (cl < startClus || cl >= lastClus)
continue;
cl -= startClus;
ALPAKA_ASSERT_OFFLOAD(cl >= 0);
ALPAKA_ASSERT_OFFLOAD(cl < MaxHitsInIter);
auto x = digis[i].xx();
auto y = digis[i].yy();
auto ch = digis[i].adc();
alpaka::atomicAdd(acc, &clusParams.charge[cl], (int32_t)ch, alpaka::hierarchy::Threads{});
ch = alpaka::math::min(acc, ch, pixmx);
if (clusParams.minRow[cl] == x)
alpaka::atomicAdd(acc, &clusParams.q_f_X[cl], (int32_t)ch, alpaka::hierarchy::Threads{});
if (clusParams.maxRow[cl] == x)
alpaka::atomicAdd(acc, &clusParams.q_l_X[cl], (int32_t)ch, alpaka::hierarchy::Threads{});
if (clusParams.minCol[cl] == y)
alpaka::atomicAdd(acc, &clusParams.q_f_Y[cl], (int32_t)ch, alpaka::hierarchy::Threads{});
if (clusParams.maxCol[cl] == y)
alpaka::atomicAdd(acc, &clusParams.q_l_Y[cl], (int32_t)ch, alpaka::hierarchy::Threads{});
}
alpaka::syncBlockThreads(acc);
// next one cluster per thread...
cms::alpakatools::for_each_element_in_block_strided(acc, nClusInIter, [&](uint32_t ic) {
auto h = first + ic; // output index in global memory
assert(h < hits.nHits());
assert(h < clusters[me + 1].clusModuleStart());
pixelCPEforDevice::position<TrackerTraits>(
cpeParams->commonParams(), cpeParams->detParams(me), clusParams, ic);
pixelCPEforDevice::errorFromDB<TrackerTraits>(
cpeParams->commonParams(), cpeParams->detParams(me), clusParams, ic);
// store it
hits[h].chargeAndStatus().charge = clusParams.charge[ic];
hits[h].chargeAndStatus().status = clusParams.status[ic];
hits[h].detectorIndex() = me;
float xl, yl;
hits[h].xLocal() = xl = clusParams.xpos[ic];
hits[h].yLocal() = yl = clusParams.ypos[ic];
hits[h].clusterSizeX() = clusParams.xsize[ic];
hits[h].clusterSizeY() = clusParams.ysize[ic];
hits[h].xerrLocal() = clusParams.xerr[ic] * clusParams.xerr[ic] + cpeParams->detParams(me).apeXX;
hits[h].yerrLocal() = clusParams.yerr[ic] * clusParams.yerr[ic] + cpeParams->detParams(me).apeYY;
// keep it local for computations
float xg, yg, zg;
// to global and compute phi...
cpeParams->detParams(me).frame.toGlobal(xl, yl, xg, yg, zg);
// here correct for the beamspot...
xg -= bs->x;
yg -= bs->y;
zg -= bs->z;
hits[h].xGlobal() = xg;
hits[h].yGlobal() = yg;
hits[h].zGlobal() = zg;
hits[h].rGlobal() = alpaka::math::sqrt(acc, xg * xg + yg * yg);
hits[h].iphi() = unsafe_atan2s<7>(yg, xg);
});
alpaka::syncBlockThreads(acc);
} // end loop on batches
}
};
} // namespace gpuPixelRecHits
} // namespace ALPAKA_ACCELERATOR_NAMESPACE
#endif // RecoLocalTracker_SiPixelRecHits_plugins_gpuPixelRecHits_h