Skip to content

Commit b569e4a

Browse files
authored
Merge pull request #47035 from makortel/trackingRecHitsSoACollectionAsync
Make implicit copy of `TrackingRecHitsSoACollection<TrackerTraits>` to host asynchronous
2 parents be3a36e + a72babc commit b569e4a

File tree

2 files changed

+17
-8
lines changed

2 files changed

+17
-8
lines changed

DataFormats/TrackingRecHitSoA/interface/alpaka/TrackingRecHitsSoACollection.h

Lines changed: 9 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -51,17 +51,22 @@ namespace cms::alpakatools {
5151
assert(deviceData.nHits() == hostData.nHits());
5252
assert(deviceData.offsetBPIX2() == hostData.offsetBPIX2());
5353
#endif
54-
// Update the contents address of the phiBinner histo container after the copy from device happened
55-
alpaka::wait(queue);
54+
return hostData;
55+
}
56+
57+
// Update the contents address of the phiBinner histo container after the copy from device happened
58+
static void postCopy(TrackingRecHitHost<TrackerTraits>& hostData) {
59+
// Don't bother if zero hits
60+
if (hostData.view().metadata().size() == 0) {
61+
return;
62+
}
5663
typename TrackingRecHitSoA<TrackerTraits>::PhiBinnerView pbv;
5764
pbv.assoc = &(hostData.view().phiBinner());
5865
pbv.offSize = -1;
5966
pbv.offStorage = nullptr;
6067
pbv.contentSize = hostData.nHits();
6168
pbv.contentStorage = hostData.view().phiBinnerStorage();
6269
hostData.view().phiBinner().initStorage(pbv);
63-
64-
return hostData;
6570
}
6671
};
6772
} // namespace cms::alpakatools

DataFormats/TrackingRecHitSoA/test/alpaka/Hits_test.cc

Lines changed: 8 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -51,12 +51,16 @@ int main() {
5151
// requires c++23 to make cms::alpakatools::CopyToHost compile using if constexpr
5252
// see https://www.open-std.org/jtc1/sc22/wg21/docs/papers/2022/p2593r0.html
5353
TrackingRecHitHost<pixelTopology::Phase1> const& host_collection = tkhit;
54+
// wait for the kernel to complete
55+
alpaka::wait(queue);
5456
#else
55-
TrackingRecHitHost<pixelTopology::Phase1> host_collection =
56-
cms::alpakatools::CopyToHost<TrackingRecHitDevice<pixelTopology::Phase1, Device> >::copyAsync(queue, tkhit);
57-
#endif
58-
// wait for the kernel and the potential copy to complete
57+
using CopyT = cms::alpakatools::CopyToHost<TrackingRecHitDevice<pixelTopology::Phase1, Device> >;
58+
TrackingRecHitHost<pixelTopology::Phase1> host_collection = CopyT::copyAsync(queue, tkhit);
59+
// wait for the kernel and the copy to complete
5960
alpaka::wait(queue);
61+
CopyT::postCopy(host_collection);
62+
#endif
63+
6064
assert(tkhit.nHits() == nHits);
6165
assert(tkhit.offsetBPIX2() == 22); // set in the kernel
6266
assert(tkhit.nHits() == host_collection.nHits());

0 commit comments

Comments
 (0)