-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathutils_spmv_cuda.h
214 lines (192 loc) · 6.83 KB
/
utils_spmv_cuda.h
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
#ifndef _UTILS_SPMV_CUDA_
#define _UTILS_SPMV_CUDA_
#include "common.h"
#include <cuda_runtime.h>
#include "cusparse.h"
#include "utils_reordering.h"
#define LONGROW_THRESHOLD 2048
#define SHORTROW_THRESHOLD 8
typedef struct SpMV_block
{
int method;
int num_threads;
int num_blocks;
int num_threads_l;
int num_blocks_l;
int m;
int *d_csrRowPtr_l;
int *d_csrColIdx_l;
VALUE_TYPE *d_csrVal_l;
VALUE_TYPE *d_x;
VALUE_TYPE *d_y;
int m_new;
int *d_longrow_idx;
int longrow;
} SpMV_block;
__global__ void spmv_longrow_csr_cuda_executor(const int *d_csrRowPtr,
const int *d_csrColIdx,
const VALUE_TYPE *d_csrVal,
const VALUE_TYPE *d_x,
VALUE_TYPE *d_y,
const int longrow,
const int *d_longrow_idx)
{
const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
const int lane_id = (WARP_SIZE - 1) & threadIdx.x;
for (int i = 0; i < longrow; i++)
{
const int rowid = d_longrow_idx[i];
const int start = d_csrRowPtr[rowid];
const int stop = d_csrRowPtr[rowid + 1];
const int len = stop - start;
VALUE_TYPE sum = 0;
if (global_id < len)
{
sum = d_x[d_csrColIdx[start + global_id]] * d_csrVal[start + global_id];
}
sum = sum_32_shfl(sum);
if (lane_id == 0 && sum != 0)
atomicAdd(&d_y[rowid], sum);
}
}
__global__ void spmv_threadsca_csr_cuda_executor(const int *d_csrRowPtr,
const int *d_csrColIdx,
const VALUE_TYPE *d_csrVal,
const int m,
const VALUE_TYPE *d_x,
VALUE_TYPE *d_y)
{
const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
if (global_id < m)
{
const int rowid = global_id;
const int start = d_csrRowPtr[rowid] - d_csrRowPtr[0];
const int stop = d_csrRowPtr[rowid + 1] - d_csrRowPtr[0];
VALUE_TYPE sum = 0;
if (stop - start <= LONGROW_THRESHOLD)
{
for (int j = start; j < stop; j++)
sum += d_x[d_csrColIdx[j]] * d_csrVal[j];
}
d_y[rowid] = sum;
}
}
__global__ void spmv_threadsca_dcsr_cuda_executor(const int *d_csrRowPtr,
const int *d_csrColIdx,
const VALUE_TYPE *d_csrVal,
const int m,
const VALUE_TYPE *d_x,
VALUE_TYPE *d_y,
const int *d_row_perm)
{
const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
if (global_id < m)
{
const int rowid = global_id;
const int start = d_csrRowPtr[rowid] - d_csrRowPtr[0];
const int stop = d_csrRowPtr[rowid + 1] - d_csrRowPtr[0];
VALUE_TYPE sum = 0;
if (stop - start <= LONGROW_THRESHOLD)
{
for (int j = start; j < stop; j++)
sum += d_x[d_csrColIdx[j]] * d_csrVal[j];
}
d_y[d_row_perm[rowid]] = sum;
}
// const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
// if (!global_id)
// {
// for (int i = 0; i < m; i++)
// {
// const int rowid = i;
// const int start = d_csrRowPtr[rowid] - d_csrRowPtr[0];
// const int stop = d_csrRowPtr[rowid + 1] - d_csrRowPtr[0];
// VALUE_TYPE sum = 0;
// if (stop - start <= LONGROW_THRESHOLD)
// {
// for (int j = start; j < stop; j++)
// sum += d_x[d_csrColIdx[j]] * d_csrVal[j];
// }
// d_y[d_row_perm[rowid]] = sum;
// printf("id = %d sum = %.1lf\n", d_row_perm[rowid], sum);
// }
// }
}
__global__ void spmv_warpvec_csr_cuda_executor(const int *d_csrRowPtr,
const int *d_csrColIdx,
const VALUE_TYPE *d_csrVal,
const int m,
const VALUE_TYPE *d_x,
VALUE_TYPE *d_y)
{
const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
// Initialize
const int rowid = global_id / WARP_SIZE;
if (rowid >= m)
return;
const int lane_id = (WARP_SIZE - 1) & threadIdx.x;
const int start = d_csrRowPtr[rowid] - d_csrRowPtr[0];
const int stop = d_csrRowPtr[rowid + 1] - d_csrRowPtr[0];
if (start == stop)
{
if (!lane_id)
d_y[rowid] = 0;
return;
}
VALUE_TYPE sum = 0;
if (stop - start <= LONGROW_THRESHOLD)
{
for (int j = start + lane_id; j < stop; j += WARP_SIZE)
{
sum += d_x[d_csrColIdx[j]] * d_csrVal[j];
}
sum = sum_32_shfl(sum);
}
//finish
if (!lane_id)
d_y[rowid] = sum;
}
__global__ void spmv_warpvec_dcsr_cuda_executor(const int *d_csrRowPtr,
const int *d_csrColIdx,
const VALUE_TYPE *d_csrVal,
const int m,
const VALUE_TYPE *d_x,
VALUE_TYPE *d_y,
const int *d_row_perm)
{
const int global_id = blockIdx.x * blockDim.x + threadIdx.x;
// Initialize
const int rowid = global_id / WARP_SIZE;
if (rowid >= m)
return;
const int lane_id = (WARP_SIZE - 1) & threadIdx.x;
const int start = d_csrRowPtr[rowid] - d_csrRowPtr[0];
const int stop = d_csrRowPtr[rowid + 1] - d_csrRowPtr[0];
if (start == stop)
{
if (!lane_id)
d_y[rowid] = 0;
return;
}
VALUE_TYPE sum = 0;
if (stop - start <= LONGROW_THRESHOLD)
{
for (int j = start + lane_id; j < stop; j += WARP_SIZE)
{
sum += d_x[d_csrColIdx[j]] * d_csrVal[j];
}
sum = sum_32_shfl(sum);
}
//finish
if (!lane_id)
d_y[d_row_perm[rowid]] = sum;
}
__global__ void subKernel(VALUE_TYPE *b,
VALUE_TYPE *y,
int m)
{
const int global_id = threadIdx.x + blockDim.x * blockIdx.x;
if (global_id < m)
b[global_id] = b[global_id] - y[global_id];
}
#endif