forked from kokkos/kokkos
-
Notifications
You must be signed in to change notification settings - Fork 0
/
TestCuda_InterOp_StreamsMultiGPU.cpp
171 lines (149 loc) · 6.73 KB
/
TestCuda_InterOp_StreamsMultiGPU.cpp
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
//@HEADER
// ************************************************************************
//
// Kokkos v. 4.0
// Copyright (2022) National Technology & Engineering
// Solutions of Sandia, LLC (NTESS).
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// the U.S. Government retains certain rights in this software.
//
// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
// See https://kokkos.org/LICENSE for license information.
// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
//
//@HEADER
#include <TestCuda_Category.hpp>
#include <Test_InterOp_Streams.hpp>
namespace {
struct StreamsAndDevices {
std::array<cudaStream_t, 2> streams;
std::array<int, 2> devices;
StreamsAndDevices() {
int n_devices;
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaGetDeviceCount(&n_devices));
devices = {0, n_devices - 1};
for (int i = 0; i < 2; ++i) {
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaSetDevice(devices[i]));
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaStreamCreate(&streams[i]));
}
}
StreamsAndDevices(const StreamsAndDevices &) = delete;
StreamsAndDevices &operator=(const StreamsAndDevices &) = delete;
~StreamsAndDevices() {
for (int i = 0; i < 2; ++i) {
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaSetDevice(devices[i]));
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaStreamDestroy(streams[i]));
}
}
};
std::array<TEST_EXECSPACE, 2> get_execution_spaces(
const StreamsAndDevices &streams_and_devices) {
TEST_EXECSPACE exec0(streams_and_devices.streams[0]);
TEST_EXECSPACE exec1(streams_and_devices.streams[1]);
// Must return void to use ASSERT_EQ
[&]() {
ASSERT_EQ(exec0.cuda_device(), streams_and_devices.devices[0]);
ASSERT_EQ(exec1.cuda_device(), streams_and_devices.devices[1]);
}();
return {exec0, exec1};
}
// Test Interoperability with Cuda Streams
void test_policies(TEST_EXECSPACE exec0, Kokkos::View<int *, TEST_EXECSPACE> v0,
TEST_EXECSPACE exec, Kokkos::View<int *, TEST_EXECSPACE> v) {
using MemorySpace = typename TEST_EXECSPACE::memory_space;
Kokkos::deep_copy(exec, v, 5);
Kokkos::deep_copy(exec0, v0, 5);
Kokkos::deep_copy(v, v0);
int sum;
int sum0;
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::Range_0",
Kokkos::RangePolicy<TEST_EXECSPACE>(exec0, 0, 100),
Test::FunctorRange<MemorySpace>(v0));
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::Range",
Kokkos::RangePolicy<TEST_EXECSPACE>(exec, 0, 100),
Test::FunctorRange<MemorySpace>(v));
Kokkos::parallel_reduce(
"Test::cuda::raw_cuda_stream::RangeReduce_0",
Kokkos::RangePolicy<TEST_EXECSPACE, Kokkos::LaunchBounds<128, 2>>(exec0,
0, 100),
Test::FunctorRangeReduce<MemorySpace>(v0), sum0);
Kokkos::parallel_reduce(
"Test::cuda::raw_cuda_stream::RangeReduce",
Kokkos::RangePolicy<TEST_EXECSPACE, Kokkos::LaunchBounds<128, 2>>(exec, 0,
100),
Test::FunctorRangeReduce<MemorySpace>(v), sum);
ASSERT_EQ(600, sum0);
ASSERT_EQ(600, sum);
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::MDRange_0",
Kokkos::MDRangePolicy<TEST_EXECSPACE, Kokkos::Rank<2>>(
exec0, {0, 0}, {10, 10}),
Test::FunctorMDRange<MemorySpace>(v0));
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::MDRange",
Kokkos::MDRangePolicy<TEST_EXECSPACE, Kokkos::Rank<2>>(
exec, {0, 0}, {10, 10}),
Test::FunctorMDRange<MemorySpace>(v));
Kokkos::parallel_reduce("Test::cuda::raw_cuda_stream::MDRangeReduce_0",
Kokkos::MDRangePolicy<TEST_EXECSPACE, Kokkos::Rank<2>,
Kokkos::LaunchBounds<128, 2>>(
exec0, {0, 0}, {10, 10}),
Test::FunctorMDRangeReduce<MemorySpace>(v0), sum0);
Kokkos::parallel_reduce("Test::cuda::raw_cuda_stream::MDRangeReduce",
Kokkos::MDRangePolicy<TEST_EXECSPACE, Kokkos::Rank<2>,
Kokkos::LaunchBounds<128, 2>>(
exec, {0, 0}, {10, 10}),
Test::FunctorMDRangeReduce<MemorySpace>(v), sum);
ASSERT_EQ(700, sum0);
ASSERT_EQ(700, sum);
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::Team_0",
Kokkos::TeamPolicy<TEST_EXECSPACE>(exec0, 10, 10),
Test::FunctorTeam<MemorySpace, TEST_EXECSPACE>(v0));
Kokkos::parallel_for("Test::cuda::raw_cuda_stream::Team",
Kokkos::TeamPolicy<TEST_EXECSPACE>(exec, 10, 10),
Test::FunctorTeam<MemorySpace, TEST_EXECSPACE>(v));
Kokkos::parallel_reduce(
"Test::cuda::raw_cuda_stream::Team_0",
Kokkos::TeamPolicy<TEST_EXECSPACE, Kokkos::LaunchBounds<128, 2>>(exec0,
10, 10),
Test::FunctorTeamReduce<MemorySpace, TEST_EXECSPACE>(v0), sum0);
Kokkos::parallel_reduce(
"Test::cuda::raw_cuda_stream::Team",
Kokkos::TeamPolicy<TEST_EXECSPACE, Kokkos::LaunchBounds<128, 2>>(exec, 10,
10),
Test::FunctorTeamReduce<MemorySpace, TEST_EXECSPACE>(v), sum);
ASSERT_EQ(800, sum0);
ASSERT_EQ(800, sum);
}
TEST(cuda_multi_gpu, managed_views) {
StreamsAndDevices streams_and_devices;
{
std::array<TEST_EXECSPACE, 2> execs =
get_execution_spaces(streams_and_devices);
Kokkos::View<int *, TEST_EXECSPACE> view0(
Kokkos::view_alloc("v0", execs[0]), 100);
Kokkos::View<int *, TEST_EXECSPACE> view(Kokkos::view_alloc("v", execs[1]),
100);
test_policies(execs[0], view0, execs[1], view);
}
}
TEST(cuda_multi_gpu, unmanaged_views) {
StreamsAndDevices streams_and_devices;
{
std::array<TEST_EXECSPACE, 2> execs =
get_execution_spaces(streams_and_devices);
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaSetDevice(execs[0].cuda_device()));
int *p0;
KOKKOS_IMPL_CUDA_SAFE_CALL(
cudaMalloc(reinterpret_cast<void **>(&p0), sizeof(int) * 100));
Kokkos::View<int *, TEST_EXECSPACE> view0(p0, 100);
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaSetDevice(execs[1].cuda_device()));
int *p;
KOKKOS_IMPL_CUDA_SAFE_CALL(
cudaMalloc(reinterpret_cast<void **>(&p), sizeof(int) * 100));
Kokkos::View<int *, TEST_EXECSPACE> view(p, 100);
test_policies(execs[0], view0, execs[1], view);
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaFree(p0));
KOKKOS_IMPL_CUDA_SAFE_CALL(cudaFree(p));
}
}
} // namespace