-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathopt_32.cuh
135 lines (109 loc) · 5.4 KB
/
opt_32.cuh
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
template <int kernel, int block_y>
__launch_bounds__(32 * block_y)
__global__ void opt1wl79_32x32x32(float *in) {
int idx = threadIdx.x;
int idy = threadIdx.y;
const int planes = block_y;
const int block_size = 32 * 32 * 32;
const int snx = 33;
const int sny = 32;
const int snxy = snx * sny;
__shared__ float smem[snxy * planes];
__shared__ float smem2[snxy * planes];
size_t block_idx =
block_size *
(blockIdx.x + gridDim.x * (blockIdx.y + gridDim.y * blockIdx.z));
const int num_batches_z = 32 / planes;
for (int batch_z = 0; batch_z < num_batches_z; ++batch_z) {
// Load all (x,y) planes into shared memory
for (int z = 0; z < planes; ++z) {
// Process an entire 32 x 32 plane
for (int tile_y = 0; tile_y < 32 / block_y; ++tile_y) {
size_t sptr = idx + snx * (tile_y * block_y + idy) + snxy * z;
size_t gptr = idx + 32 * (tile_y * block_y + idy) + 1024 * z;
smem[sptr] = in[batch_z * planes * 1024 + gptr + block_idx];
}
}
__syncthreads();
// Apply wavelet transform line by line in the x-direction
for (int z = 0; z < planes / block_y; ++z) {
if (kernel == 0) {
ds79_compute_shared<32>(
&smem[snx * idx + snxy * (idy + z * block_y)],
&smem2[snx * idx + snxy * (idy + z * block_y)], 1);
} else {
us79_compute_shared<32>(
&smem[snx * idx + snxy * (idy + z * block_y)],
&smem2[snx * idx + snxy * (idy + z * block_y)], 1);
}
}
__syncthreads();
// Apply wavelet transform line by line in the y-direction
for (int z = 0; z < planes / block_y; ++z) {
if (kernel == 0) {
ds79_compute_shared<32>(
&smem[idx + snxy * (idy + z * block_y)],
&smem2[idx + snxy * (idy + z * block_y)], snx);
} else {
us79_compute_shared<32>(
&smem[idx + snxy * (idy + z * block_y)],
&smem2[idx + snxy * (idy + z * block_y)], snx);
}
}
__syncthreads();
// Write result to global memory
// Write all (x,y) planes back to global memory
for (int z = 0; z < planes; ++z) {
// Process an entire 32 x 32 plane
for (int tile_y = 0; tile_y < 32 / block_y; ++tile_y) {
size_t sptr = idx + snx * (tile_y * block_y + idy) + snxy * z;
size_t gptr = idx + 32 * (tile_y * block_y + idy) + 1024 * z;
in[batch_z * planes * 1024 + gptr + block_idx] = smem[sptr];
}
}
__syncthreads();
}
const int num_batches_y = 32 / planes;
for (int batch_y = 0; batch_y < num_batches_y; ++batch_y) {
// Load all (x,z) planes into shared memory
for (int y = 0; y < planes; ++y) {
// Process an entire 32 x 32 plane
for (int tile_z = 0; tile_z < 32 / block_y; ++tile_z) {
size_t sptr = idx + snx * (tile_z * block_y + idy) + snxy * y;
size_t gptr = idx + 1024 * (tile_z * block_y + idy) + 32 * y;
smem[sptr] = in[batch_y * planes * 32 + gptr + block_idx];
}
}
__syncthreads();
// Apply wavelet transform line by line in the z-direction
for (int y = 0; y < planes / block_y; ++y) {
if (kernel == 0) {
ds79_compute_shared<32>(
&smem[idx + snxy * (idy + y * block_y)],
&smem2[idx + snxy * (idy + y * block_y)], snx);
} else {
us79_compute_shared<32>(
&smem[idx + snxy * (idy + y * block_y)],
&smem2[idx + snxy * (idy + y * block_y)], snx);
}
}
__syncthreads();
// Write all (x,z) planes back to global memory
for (int y = 0; y < planes; ++y) {
// Process an entire 32 x 32 plane
for (int tile_z = 0; tile_z < 32 / block_y; ++tile_z) {
size_t sptr = idx + snx * (tile_z * block_y + idy) + snxy * y;
size_t gptr = idx + 1024 * (tile_z * block_y + idy) + 32 * y;
in[batch_y * planes * 32 + gptr + block_idx] = smem[sptr];
}
}
}
}
template <int mode>
void opt1wl79_32x32x32_h(float *in, const int bx, const int by, const int bz) {
const int block_y = 4;
dim3 threads(32, block_y, 1);
dim3 blocks(bx, by, bz);
opt1wl79_32x32x32<mode, block_y><<<blocks, threads>>>(in);
cudaErrCheck(cudaPeekAtLastError());
}