forked from csc-training/high-level-gpu-programming
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpingpong.cpp
234 lines (178 loc) · 6.88 KB
/
pingpong.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
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
#include <cstdio>
#include <cstdlib>
#include <mpi.h>
#include <unistd.h>
#include <cuda_runtime_api.h>
#include "error_checks.h"
// Kernel call wrapper function prototype
void call_kernel(double *data, int N, int blocksize, int tib);
/*
This routine can be used to inspect the properties of a node
Return arguments:
nodeRank (int *) -- My rank in the node communicator
nodeProcs (int *) -- Total number of processes in this node
devCount (int *) -- Number of CUDA devices available in the node
*/
void getNodeInfo(int *nodeRank, int *nodeProcs, int *devCount)
{
MPI_Comm intranodecomm;
MPI_Comm_split_type(MPI_COMM_WORLD, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &intranodecomm);
MPI_Comm_rank(intranodecomm, nodeRank);
MPI_Comm_size(intranodecomm, nodeProcs);
MPI_Comm_free(&intranodecomm);
CUDA_CHECK( cudaGetDeviceCount(devCount) );
}
/* Test routine for CPU-to-CPU copy */
void CPUtoCPUtest(int rank, double *data, int N, double &timer)
{
double start, stop;
start = MPI_Wtime();
if (rank == 0) {
MPI_Send(data, N, MPI_DOUBLE, 1, 11, MPI_COMM_WORLD);
MPI_Recv(data, N, MPI_DOUBLE, 1, 12, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
} else {
MPI_Recv(data, N, MPI_DOUBLE, 0, 11, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
/* Add one*/
for (int i = 0; i < N; ++i)
data[i] += 1.0;
MPI_Send(data, N, MPI_DOUBLE, 0, 12, MPI_COMM_WORLD);
}
stop = MPI_Wtime();
timer = stop - start;
}
/* Test routine for GPU-CPU-to-CPU-GPU copy */
void GPUtoGPUtestManual(int rank, double *hA, double *dA, int N, double &timer)
{
double start, stop;
start = MPI_Wtime();
//Transfer that uses manual copies to host, and MPI on
//host.
if (rank == 0) { //Sender process
CUDA_CHECK( cudaMemcpy(hA, dA, sizeof(double)*N,
cudaMemcpyDeviceToHost) );
/* Send data to rank 1 for addition */
MPI_Send(hA, N, MPI_DOUBLE, 1, 11, MPI_COMM_WORLD);
/* Receive the added data back */
MPI_Recv(hA, N, MPI_DOUBLE, 1, 12, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
CUDA_CHECK( cudaMemcpy(dA, hA, sizeof(double)*N,
cudaMemcpyHostToDevice) );
} else { // Adder process
MPI_Recv(hA, N, MPI_DOUBLE, 0, 11, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
CUDA_CHECK( cudaMemcpy(dA, hA, sizeof(double)*N,
cudaMemcpyHostToDevice) );
int tib = 128;
int grid = (N + tib - 1) / tib;
call_kernel(dA, N, grid, tib);
CUDA_CHECK( cudaMemcpy(hA, dA, sizeof(double)*N,
cudaMemcpyDeviceToHost) );
MPI_Send(hA, N, MPI_DOUBLE, 0, 12, MPI_COMM_WORLD);
}
stop = MPI_Wtime();
timer = stop - start;
}
/* Test routine for GPU-CPU-to-CPU-GPU copy */
void GPUtoGPUtestCudaAware(int rank, double *dA, int N, double &timer)
{
double start, stop;
start = MPI_Wtime();
//Transfer that uses CUDA-aware MPI to transfer data
if (rank == 0) { //Sender process
/* Send data to rank 1 for addition */
MPI_Send(dA, N, MPI_DOUBLE, 1, 11, MPI_COMM_WORLD);
/* Receive the added data back */
MPI_Recv(dA, N, MPI_DOUBLE, 1, 12, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
} else { // Adder process
int tib = 128;
int grid = (N + tib - 1) / tib;
MPI_Recv(dA, N, MPI_DOUBLE, 0, 11, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
call_kernel(dA, N, grid, tib);
CUDA_CHECK(cudaStreamSynchronize(0));
MPI_Send(dA, N, MPI_DOUBLE, 0, 12, MPI_COMM_WORLD);
}
stop = MPI_Wtime();
timer = stop - start;
}
/* Simple ping-pong main program */
int main(int argc, char *argv[])
{
int rank, nprocs, noderank, nodenprocs, devcount;
double GPUtime, CPUtime;
double *dA, *hA;
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
if (argc < 2) {
printf("Need the vector length as argument\n");
exit(EXIT_FAILURE);
}
int N = atoi(argv[1]);
getNodeInfo(&noderank, &nodenprocs, &devcount);
/* Due to the test, we need exactly two processes with one GPU for
each */
if (nprocs != 2) {
printf("Need exactly two processes!\n");
exit(EXIT_FAILURE);
}
if (devcount == 0) {
printf("Could now find any CUDA devices.\n");
exit(EXIT_FAILURE);
}
if (nodenprocs > devcount) {
printf("Not enough GPUs for all processes in the node.\n");
exit(EXIT_FAILURE);
}
//Select the device according to the node rank
CUDA_CHECK( cudaSetDevice(noderank) );
//allocate device memories
CUDA_CHECK( cudaMallocHost((void **)&hA, sizeof(double) * N) );
CUDA_CHECK( cudaMalloc((void **)&dA, sizeof(double) * N) );
// Dummy transfer to remove the overhead of the first communication
CPUtoCPUtest(rank, hA, N, CPUtime);
/* Re-initialize and copy the data to the device memory to prepare for
* MPI test */
for (int i = 0; i < N; ++i)
hA[i] = 1.0;
/* CPU-to-CPU test */
CPUtoCPUtest(rank, hA, N, CPUtime);
if (rank == 0) {
double errorsum = 0;
for (int i = 0; i < N; ++i)
errorsum += hA[i] - 2.0;
printf("CPU-CPU time %f, errorsum %f\n", CPUtime, errorsum);
}
// Dummy transfer to remove the overhead of the first communication
GPUtoGPUtestCudaAware(rank, dA, N, GPUtime);
/* Re-initialize and copy the data to the device memory */
for (int i = 0; i < N; ++i)
hA[i] = 1.0;
CUDA_CHECK( cudaMemcpy(dA, hA, sizeof(double)*N, cudaMemcpyHostToDevice) );
/* GPU-to-GPU test, Cuda-aware */
GPUtoGPUtestCudaAware(rank, dA, N, GPUtime);
/*Check results, copy device array back to Host*/
CUDA_CHECK( cudaMemcpy(hA, dA, sizeof(double)*N, cudaMemcpyDeviceToHost) );
if (rank == 0) {
double errorsum = 0;
for (int i = 0; i < N; ++i)
errorsum += hA[i] - 2.0;
printf("GPU-GPU cuda-aware time %f, errorsum %f\n", GPUtime, errorsum);
}
// Dummy transfer to remove the overhead of the first communication
GPUtoGPUtestManual(rank, hA, dA, N, GPUtime);
/* Re-initialize and copy the data to the device memory to prepare for
* MPI test */
for (int i = 0; i < N; ++i)
hA[i] = 1.0;
CUDA_CHECK( cudaMemcpy(dA, hA, sizeof(double)*N, cudaMemcpyHostToDevice) );
/* GPU-to-GPU test, Manual option*/
GPUtoGPUtestManual(rank, hA, dA, N, GPUtime);
/*Check results, copy device array back to Host*/
CUDA_CHECK( cudaMemcpy(hA, dA, sizeof(double)*N, cudaMemcpyDeviceToHost) );
if (rank == 0) {
double errorsum = 0;
for (int i = 0; i < N; ++i)
errorsum += hA[i] - 2.0;
printf("GPU-GPU manual time %f, errorsum %f\n", GPUtime, errorsum);
}
MPI_Finalize();
return 0;
}