-
Notifications
You must be signed in to change notification settings - Fork 1.8k
/
cuSolverSp_LowlevelQR.cpp
449 lines (375 loc) · 15.6 KB
/
cuSolverSp_LowlevelQR.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
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
/* Copyright (c) 2022, NVIDIA CORPORATION. All rights reserved.
*
* Redistribution and use in source and binary forms, with or without
* modification, are permitted provided that the following conditions
* are met:
* * Redistributions of source code must retain the above copyright
* notice, this list of conditions and the following disclaimer.
* * Redistributions in binary form must reproduce the above copyright
* notice, this list of conditions and the following disclaimer in the
* documentation and/or other materials provided with the distribution.
* * Neither the name of NVIDIA CORPORATION nor the names of its
* contributors may be used to endorse or promote products derived
* from this software without specific prior written permission.
*
* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS ``AS IS'' AND ANY
* EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
* IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
* PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
* CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
* EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
* PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
* PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY
* OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
* (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
* OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*/
#include <assert.h>
#include <ctype.h>
#include <cuda_runtime.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "cusolverSp.h"
#include "cusolverSp_LOWLEVEL_PREVIEW.h"
#include "helper_cuda.h"
#include "helper_cusolver.h"
template <typename T_ELEM>
int loadMMSparseMatrix(char *filename, char elem_type, bool csrFormat, int *m,
int *n, int *nnz, T_ELEM **aVal, int **aRowInd,
int **aColInd, int extendSymMatrix);
void UsageSP(void) {
printf("<options>\n");
printf("-h : display this help\n");
printf("-file=<filename> : filename containing a matrix in MM format\n");
printf("-device=<device_id> : <device_id> if want to run on specific GPU\n");
exit(0);
}
void parseCommandLineArguments(int argc, char *argv[], struct testOpts &opts) {
memset(&opts, 0, sizeof(opts));
if (checkCmdLineFlag(argc, (const char **)argv, "-h")) {
UsageSP();
}
if (checkCmdLineFlag(argc, (const char **)argv, "file")) {
char *fileName = 0;
getCmdLineArgumentString(argc, (const char **)argv, "file", &fileName);
if (fileName) {
opts.sparse_mat_filename = fileName;
} else {
printf("\nIncorrect filename passed to -file \n ");
UsageSP();
}
}
}
int main(int argc, char *argv[]) {
struct testOpts opts;
cusolverSpHandle_t cusolverSpH =
NULL; // reordering, permutation and 1st LU factorization
cusparseHandle_t cusparseH = NULL; // residual evaluation
cudaStream_t stream = NULL;
cusparseMatDescr_t descrA = NULL; // A is a base-0 general matrix
csrqrInfoHost_t h_info =
NULL; // opaque info structure for LU with parital pivoting
csrqrInfo_t d_info =
NULL; // opaque info structure for LU with parital pivoting
int rowsA = 0; // number of rows of A
int colsA = 0; // number of columns of A
int nnzA = 0; // number of nonzeros of A
int baseA = 0; // base index in CSR format
// CSR(A) from I/O
int *h_csrRowPtrA = NULL; // <int> n+1
int *h_csrColIndA = NULL; // <int> nnzA
double *h_csrValA = NULL; // <double> nnzA
double *h_x = NULL; // <double> n, x = A \ b
double *h_b = NULL; // <double> n, b = ones(m,1)
double *h_bcopy = NULL; // <double> n, b = ones(m,1)
double *h_r = NULL; // <double> n, r = b - A*x
size_t size_internal = 0;
size_t size_chol = 0; // size of working space for csrlu
void *buffer_cpu = NULL; // working space for Cholesky
void *buffer_gpu = NULL; // working space for Cholesky
int *d_csrRowPtrA = NULL; // <int> n+1
int *d_csrColIndA = NULL; // <int> nnzA
double *d_csrValA = NULL; // <double> nnzA
double *d_x = NULL; // <double> n, x = A \ b
double *d_b = NULL; // <double> n, a copy of h_b
double *d_r = NULL; // <double> n, r = b - A*x
// the constants used in residual evaluation, r = b - A*x
const double minus_one = -1.0;
const double one = 1.0;
const double zero = 0.0;
// the constant used in cusolverSp
// singularity is -1 if A is invertible under tol
// tol determines the condition of singularity
int singularity = 0;
const double tol = 1.e-14;
double x_inf = 0.0; // |x|
double r_inf = 0.0; // |r|
double A_inf = 0.0; // |A|
parseCommandLineArguments(argc, argv, opts);
findCudaDevice(argc, (const char **)argv);
if (opts.sparse_mat_filename == NULL) {
opts.sparse_mat_filename = sdkFindFilePath("lap2D_5pt_n32.mtx", argv[0]);
if (opts.sparse_mat_filename != NULL)
printf("Using default input file [%s]\n", opts.sparse_mat_filename);
else
printf("Could not find lap2D_5pt_n32.mtx\n");
} else {
printf("Using input file [%s]\n", opts.sparse_mat_filename);
}
printf("step 1: read matrix market format\n");
if (opts.sparse_mat_filename) {
if (loadMMSparseMatrix<double>(opts.sparse_mat_filename, 'd', true, &rowsA,
&colsA, &nnzA, &h_csrValA, &h_csrRowPtrA,
&h_csrColIndA, true)) {
return 1;
}
baseA = h_csrRowPtrA[0]; // baseA = {0,1}
} else {
fprintf(stderr, "Error: input matrix is not provided\n");
return 1;
}
if (rowsA != colsA) {
fprintf(stderr, "Error: only support square matrix\n");
return 1;
}
printf("sparse matrix A is %d x %d with %d nonzeros, base=%d\n", rowsA, colsA,
nnzA, baseA);
checkCudaErrors(cusolverSpCreate(&cusolverSpH));
checkCudaErrors(cusparseCreate(&cusparseH));
checkCudaErrors(cudaStreamCreate(&stream));
checkCudaErrors(cusolverSpSetStream(cusolverSpH, stream));
checkCudaErrors(cusparseSetStream(cusparseH, stream));
checkCudaErrors(cusparseCreateMatDescr(&descrA));
checkCudaErrors(cusparseSetMatType(descrA, CUSPARSE_MATRIX_TYPE_GENERAL));
if (baseA) {
checkCudaErrors(cusparseSetMatIndexBase(descrA, CUSPARSE_INDEX_BASE_ONE));
} else {
checkCudaErrors(cusparseSetMatIndexBase(descrA, CUSPARSE_INDEX_BASE_ZERO));
}
h_x = (double *)malloc(sizeof(double) * colsA);
h_b = (double *)malloc(sizeof(double) * rowsA);
h_bcopy = (double *)malloc(sizeof(double) * rowsA);
h_r = (double *)malloc(sizeof(double) * rowsA);
assert(NULL != h_x);
assert(NULL != h_b);
assert(NULL != h_bcopy);
assert(NULL != h_r);
checkCudaErrors(
cudaMalloc((void **)&d_csrRowPtrA, sizeof(int) * (rowsA + 1)));
checkCudaErrors(cudaMalloc((void **)&d_csrColIndA, sizeof(int) * nnzA));
checkCudaErrors(cudaMalloc((void **)&d_csrValA, sizeof(double) * nnzA));
checkCudaErrors(cudaMalloc((void **)&d_x, sizeof(double) * colsA));
checkCudaErrors(cudaMalloc((void **)&d_b, sizeof(double) * rowsA));
checkCudaErrors(cudaMalloc((void **)&d_r, sizeof(double) * rowsA));
for (int row = 0; row < rowsA; row++) {
h_b[row] = 1.0;
}
memcpy(h_bcopy, h_b, sizeof(double) * rowsA);
printf("step 2: create opaque info structure\n");
checkCudaErrors(cusolverSpCreateCsrqrInfoHost(&h_info));
printf("step 3: analyze qr(A) to know structure of L\n");
checkCudaErrors(cusolverSpXcsrqrAnalysisHost(cusolverSpH, rowsA, colsA, nnzA,
descrA, h_csrRowPtrA,
h_csrColIndA, h_info));
printf("step 4: workspace for qr(A)\n");
checkCudaErrors(cusolverSpDcsrqrBufferInfoHost(
cusolverSpH, rowsA, colsA, nnzA, descrA, h_csrValA, h_csrRowPtrA,
h_csrColIndA, h_info, &size_internal, &size_chol));
if (buffer_cpu) {
free(buffer_cpu);
}
buffer_cpu = (void *)malloc(sizeof(char) * size_chol);
assert(NULL != buffer_cpu);
printf("step 5: compute A = L*L^T \n");
checkCudaErrors(cusolverSpDcsrqrSetupHost(cusolverSpH, rowsA, colsA, nnzA,
descrA, h_csrValA, h_csrRowPtrA,
h_csrColIndA, zero, h_info));
checkCudaErrors(cusolverSpDcsrqrFactorHost(cusolverSpH, rowsA, colsA, nnzA,
NULL, NULL, h_info, buffer_cpu));
printf("step 6: check if the matrix is singular \n");
checkCudaErrors(
cusolverSpDcsrqrZeroPivotHost(cusolverSpH, h_info, tol, &singularity));
if (0 <= singularity) {
fprintf(stderr, "Error: A is not invertible, singularity=%d\n",
singularity);
return 1;
}
printf("step 7: solve A*x = b \n");
checkCudaErrors(cusolverSpDcsrqrSolveHost(cusolverSpH, rowsA, colsA, h_b, h_x,
h_info, buffer_cpu));
printf("step 8: evaluate residual r = b - A*x (result on CPU)\n");
// use GPU gemv to compute r = b - A*x
checkCudaErrors(cudaMemcpy(d_csrRowPtrA, h_csrRowPtrA,
sizeof(int) * (rowsA + 1),
cudaMemcpyHostToDevice));
checkCudaErrors(cudaMemcpy(d_csrColIndA, h_csrColIndA, sizeof(int) * nnzA,
cudaMemcpyHostToDevice));
checkCudaErrors(cudaMemcpy(d_csrValA, h_csrValA, sizeof(double) * nnzA,
cudaMemcpyHostToDevice));
checkCudaErrors(
cudaMemcpy(d_r, h_bcopy, sizeof(double) * rowsA, cudaMemcpyHostToDevice));
checkCudaErrors(
cudaMemcpy(d_x, h_x, sizeof(double) * colsA, cudaMemcpyHostToDevice));
/* Wrap raw data into cuSPARSE generic API objects */
cusparseSpMatDescr_t matA = NULL;
if (baseA) {
checkCudaErrors(cusparseCreateCsr(&matA, rowsA, colsA, nnzA, d_csrRowPtrA,
d_csrColIndA, d_csrValA,
CUSPARSE_INDEX_32I, CUSPARSE_INDEX_32I,
CUSPARSE_INDEX_BASE_ONE, CUDA_R_64F));
} else {
checkCudaErrors(cusparseCreateCsr(&matA, rowsA, colsA, nnzA, d_csrRowPtrA,
d_csrColIndA, d_csrValA,
CUSPARSE_INDEX_32I, CUSPARSE_INDEX_32I,
CUSPARSE_INDEX_BASE_ZERO, CUDA_R_64F));
}
cusparseDnVecDescr_t vecx = NULL;
checkCudaErrors(cusparseCreateDnVec(&vecx, colsA, d_x, CUDA_R_64F));
cusparseDnVecDescr_t vecAx = NULL;
checkCudaErrors(cusparseCreateDnVec(&vecAx, rowsA, d_r, CUDA_R_64F));
/* Allocate workspace for cuSPARSE */
size_t bufferSize = 0;
checkCudaErrors(cusparseSpMV_bufferSize(
cusparseH, CUSPARSE_OPERATION_NON_TRANSPOSE, &minus_one, matA, vecx, &one,
vecAx, CUDA_R_64F, CUSPARSE_SPMV_ALG_DEFAULT, &bufferSize));
void *buffer = NULL;
checkCudaErrors(cudaMalloc(&buffer, bufferSize));
checkCudaErrors(cusparseSpMV(cusparseH, CUSPARSE_OPERATION_NON_TRANSPOSE,
&minus_one, matA, vecx, &one, vecAx, CUDA_R_64F,
CUSPARSE_SPMV_ALG_DEFAULT, buffer));
checkCudaErrors(
cudaMemcpy(h_r, d_r, sizeof(double) * rowsA, cudaMemcpyDeviceToHost));
x_inf = vec_norminf(colsA, h_x);
r_inf = vec_norminf(rowsA, h_r);
A_inf = csr_mat_norminf(rowsA, colsA, nnzA, descrA, h_csrValA, h_csrRowPtrA,
h_csrColIndA);
printf("(CPU) |b - A*x| = %E \n", r_inf);
printf("(CPU) |A| = %E \n", A_inf);
printf("(CPU) |x| = %E \n", x_inf);
printf("(CPU) |b - A*x|/(|A|*|x|) = %E \n", r_inf / (A_inf * x_inf));
printf("step 9: create opaque info structure\n");
checkCudaErrors(cusolverSpCreateCsrqrInfo(&d_info));
checkCudaErrors(cudaMemcpy(d_csrRowPtrA, h_csrRowPtrA,
sizeof(int) * (rowsA + 1),
cudaMemcpyHostToDevice));
checkCudaErrors(cudaMemcpy(d_csrColIndA, h_csrColIndA, sizeof(int) * nnzA,
cudaMemcpyHostToDevice));
checkCudaErrors(cudaMemcpy(d_csrValA, h_csrValA, sizeof(double) * nnzA,
cudaMemcpyHostToDevice));
checkCudaErrors(
cudaMemcpy(d_b, h_bcopy, sizeof(double) * rowsA, cudaMemcpyHostToDevice));
printf("step 10: analyze qr(A) to know structure of L\n");
checkCudaErrors(cusolverSpXcsrqrAnalysis(cusolverSpH, rowsA, colsA, nnzA,
descrA, d_csrRowPtrA, d_csrColIndA,
d_info));
printf("step 11: workspace for qr(A)\n");
checkCudaErrors(cusolverSpDcsrqrBufferInfo(
cusolverSpH, rowsA, colsA, nnzA, descrA, d_csrValA, d_csrRowPtrA,
d_csrColIndA, d_info, &size_internal, &size_chol));
printf("GPU buffer size = %lld bytes\n", (signed long long)size_chol);
if (buffer_gpu) {
checkCudaErrors(cudaFree(buffer_gpu));
}
checkCudaErrors(cudaMalloc(&buffer_gpu, sizeof(char) * size_chol));
printf("step 12: compute A = L*L^T \n");
checkCudaErrors(cusolverSpDcsrqrSetup(cusolverSpH, rowsA, colsA, nnzA, descrA,
d_csrValA, d_csrRowPtrA, d_csrColIndA,
zero, d_info));
checkCudaErrors(cusolverSpDcsrqrFactor(cusolverSpH, rowsA, colsA, nnzA, NULL,
NULL, d_info, buffer_gpu));
printf("step 13: check if the matrix is singular \n");
checkCudaErrors(
cusolverSpDcsrqrZeroPivot(cusolverSpH, d_info, tol, &singularity));
if (0 <= singularity) {
fprintf(stderr, "Error: A is not invertible, singularity=%d\n",
singularity);
return 1;
}
printf("step 14: solve A*x = b \n");
checkCudaErrors(cusolverSpDcsrqrSolve(cusolverSpH, rowsA, colsA, d_b, d_x,
d_info, buffer_gpu));
checkCudaErrors(
cudaMemcpy(d_r, h_bcopy, sizeof(double) * rowsA, cudaMemcpyHostToDevice));
checkCudaErrors(cusparseSpMV(cusparseH, CUSPARSE_OPERATION_NON_TRANSPOSE,
&minus_one, matA, vecx, &one, vecAx, CUDA_R_64F,
CUSPARSE_SPMV_ALG_DEFAULT, buffer));
checkCudaErrors(
cudaMemcpy(h_r, d_r, sizeof(double) * rowsA, cudaMemcpyDeviceToHost));
r_inf = vec_norminf(rowsA, h_r);
printf("(GPU) |b - A*x| = %E \n", r_inf);
printf("(GPU) |b - A*x|/(|A|*|x|) = %E \n", r_inf / (A_inf * x_inf));
if (cusolverSpH) {
checkCudaErrors(cusolverSpDestroy(cusolverSpH));
}
if (cusparseH) {
checkCudaErrors(cusparseDestroy(cusparseH));
}
if (stream) {
checkCudaErrors(cudaStreamDestroy(stream));
}
if (descrA) {
checkCudaErrors(cusparseDestroyMatDescr(descrA));
}
if (h_info) {
checkCudaErrors(cusolverSpDestroyCsrqrInfoHost(h_info));
}
if (d_info) {
checkCudaErrors(cusolverSpDestroyCsrqrInfo(d_info));
}
if (matA) {
checkCudaErrors(cusparseDestroySpMat(matA));
}
if (vecx) {
checkCudaErrors(cusparseDestroyDnVec(vecx));
}
if (vecAx) {
checkCudaErrors(cusparseDestroyDnVec(vecAx));
}
if (h_csrValA) {
free(h_csrValA);
}
if (h_csrRowPtrA) {
free(h_csrRowPtrA);
}
if (h_csrColIndA) {
free(h_csrColIndA);
}
if (h_x) {
free(h_x);
}
if (h_b) {
free(h_b);
}
if (h_bcopy) {
free(h_bcopy);
}
if (h_r) {
free(h_r);
}
if (buffer_cpu) {
free(buffer_cpu);
}
if (buffer_gpu) {
checkCudaErrors(cudaFree(buffer_gpu));
}
if (d_csrValA) {
checkCudaErrors(cudaFree(d_csrValA));
}
if (d_csrRowPtrA) {
checkCudaErrors(cudaFree(d_csrRowPtrA));
}
if (d_csrColIndA) {
checkCudaErrors(cudaFree(d_csrColIndA));
}
if (d_x) {
checkCudaErrors(cudaFree(d_x));
}
if (d_b) {
checkCudaErrors(cudaFree(d_b));
}
if (d_r) {
checkCudaErrors(cudaFree(d_r));
}
return 0;
}