-
Notifications
You must be signed in to change notification settings - Fork 2
/
jacobi3d_cuda_vanilla.cu
149 lines (127 loc) · 4.2 KB
/
jacobi3d_cuda_vanilla.cu
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
#include <iostream>
#include <sstream>
#include <stdexcept>
#include <vector>
#include "eval.h"
#define GET(X, Y, Z) gridOld[((X + dimX) % dimX) + ((Y + dimY) % dimY) * dimX + ((Z + dimZ) % dimZ) * dimX * dimY]
#define SET(X, Y, Z) gridNew[((X + dimX) % dimX) + ((Y + dimY) % dimY) * dimX + ((Z + dimZ) % dimZ) * dimX * dimY]
__global__ void update(double *gridOld, double *gridNew, int dimX, int dimY, int dimZ)
{
int x = threadIdx.x + blockIdx.x * blockDim.x;
int y = threadIdx.y + blockIdx.y * blockDim.y;
int z = threadIdx.z + blockIdx.z * blockDim.z;
SET(x, y, z) = (GET(x, y, z - 1) +
GET(x, y - 1, z ) +
GET(x - 1, y, z ) +
GET(x, y, z ) +
GET(x + 1, y, z ) +
GET(x, y + 1, z ) +
GET(x, y, z + 1)) * (1.0 / 7.0);
}
void init(double *gridNew, int dimX, int dimY, int dimZ)
{
for (int z = 0; z < dimZ; ++z) {
for (int y = 0; y < dimY; ++y) {
for (int x = 0; x < dimX; ++x) {
double value = 0;
if ((x * y * z) == 0) {
value = 1;
}
SET(x, y, z) = value;
}
}
}
}
void print(double *gridOld, int dimX, int dimY, int dimZ)
{
for (int z = 0; z < dimZ; ++z) {
for (int y = 0; y < dimY; ++y) {
for (int x = 0; x < dimX; ++x) {
std::cout << " " << GET(x, y, z);
}
std::cout << "\n";
}
}
}
int divAndRoundUp(int dim, int blockDim)
{
int res = dim / blockDim;
if (dim % blockDim) {
res += 1;
}
return res;
}
void checkForCUDAError()
{
cudaError_t error = cudaGetLastError();
if (error != cudaSuccess) {
std::cerr << "ERROR: " << cudaGetErrorString(error) << "\n";
throw std::runtime_error("CUDA error");
}
}
void benchmark(std::vector<double> *gridOld, std::vector<double> *gridNew, int dimX, int dimY, int dimZ, int repeats, dim3 blockDim)
{
checkForCUDAError();
cudaDeviceSynchronize();
double tStartInit = getUTtime();
int byteSize = dimX * dimY * dimZ * sizeof(double);
double *devGridOld;
double *devGridNew;
cudaMalloc(&devGridOld, byteSize);
cudaMalloc(&devGridNew, byteSize);
cudaMemcpy(devGridOld, &gridOld->front(), byteSize, cudaMemcpyHostToDevice);
cudaMemcpy(devGridNew, &gridNew->front(), byteSize, cudaMemcpyHostToDevice);
dim3 gridDim(divAndRoundUp(dimX, blockDim.x),
divAndRoundUp(dimY, blockDim.y),
divAndRoundUp(dimZ, blockDim.z));
cudaDeviceSynchronize();
double tStartCalc = getUTtime();
for (int t = 0; t < repeats; ++t) {
update<<<gridDim, blockDim>>>(devGridOld, devGridNew, dimX, dimY, dimZ);
std::swap(devGridOld, devGridNew);
}
cudaDeviceSynchronize();
double tEndCalc = getUTtime();
cudaMemcpy(&gridOld->front(), devGridOld, byteSize, cudaMemcpyDeviceToHost);
cudaFree(devGridOld);
cudaFree(devGridNew);
cudaDeviceSynchronize();
double tEnd = getUTtime();
checkForCUDAError();
eval(tStartInit, tStartCalc, tEndCalc, tEnd, dimX, dimY, dimZ, repeats);
}
int main(int argc, char **argv)
{
if ((argc < 6) || (argc > 9)) {
std::cerr << "usage: " << argv[0] << " DIM_X DIM_Y DIM_Z REPEATS CUDA_DEVICE [BLOCK_DIM_X=32] [BLOCK_DIM_Y=32] [BLOCK_DIM_Z=1] \n";
return 1;
}
std::stringstream buf;
for (int i = 1; i < argc; ++i) {
buf << argv[i] << " ";
}
int dimX, dimY, dimZ, repeats, cudaDevice;
buf >> dimX;
buf >> dimY;
buf >> dimZ;
buf >> repeats;
buf >> cudaDevice;
cudaSetDevice(cudaDevice);
dim3 blockDim(32, 32, 1);
if (argc > 5) {
buf >> blockDim.x;
}
if (argc > 6) {
buf >> blockDim.y;
}
if (argc > 7) {
buf >> blockDim.z;
}
int size = dimX * dimY * dimZ;
std::vector<double> gridOld(size);
std::vector<double> gridNew(size);
init(&gridOld[0], dimX, dimY, dimZ);
init(&gridNew[0], dimX, dimY, dimZ);
benchmark(&gridOld, &gridNew, dimX, dimY, dimZ, repeats, blockDim);
print(&gridOld[0], dimX, dimY, dimZ);
}