Fix spsv with CSRG handling of descriptors.

repack-nvid
Salvatore Filippone 10 months ago
parent d28ea462d9
commit 62db7c0449

@ -263,42 +263,45 @@ int T_spsvCSRGDevice(T_Cmat *Matrix, TYPE alpha, void *deviceX,
vY=y->v_; vY=y->v_;
CHECK_CUSPARSE( cusparseCreateDnVec(&vecY, cMat->m, vY, CUSPARSE_BASE_TYPE) ); CHECK_CUSPARSE( cusparseCreateDnVec(&vecY, cMat->m, vY, CUSPARSE_BASE_TYPE) );
CHECK_CUSPARSE( cusparseCreateDnVec(&vecX, cMat->n, vX, CUSPARSE_BASE_TYPE) ); CHECK_CUSPARSE( cusparseCreateDnVec(&vecX, cMat->n, vX, CUSPARSE_BASE_TYPE) );
if (T_CSRGIsNullMvDescr(cMat)) {
cMat->spmvDescr = (cusparseSpMatDescr_t *) malloc(sizeof(cusparseSpMatDescr_t *));
}
T_CSRGCreateSpMVDescr(cMat);
// fprintf(stderr,"Entry to SpSVDevice: %d %p\n", // fprintf(stderr,"Entry to SpSVDevice: %d %p\n",
// T_CSRGIsNullSvDescr(cMat),cMat->spsvDescr); // T_CSRGIsNullSvDescr(cMat),cMat->spsvDescr);
if (T_CSRGIsNullSvDescr(cMat)) { if (T_CSRGIsNullSvDescr(cMat)) {
cMat->spsvDescr=(cusparseSpSVDescr_t *) malloc(sizeof(cusparseSpSVDescr_t *)); cMat->spsvDescr=(cusparseSpSVDescr_t *) malloc(sizeof(cusparseSpSVDescr_t *));
cMat->svbsize=0; cMat->svbsize=0;
CHECK_CUSPARSE( cusparseSpSV_createDescr(cMat->spsvDescr) ); CHECK_CUSPARSE( cusparseSpSV_createDescr(cMat->spsvDescr) );
} //fprintf(stderr,"Entry to SpSVDevice: %d %p %d\n",
//fprintf(stderr,"Entry to SpSVDevice: %d %p %d\n", // T_CSRGIsNullSvDescr(cMat),cMat->spsvDescr,cMat->svbsize);
// T_CSRGIsNullSvDescr(cMat),cMat->spsvDescr,cMat->svbsize); CHECK_CUSPARSE(cusparseSpSV_bufferSize(*my_handle,CUSPARSE_OPERATION_NON_TRANSPOSE,
CHECK_CUSPARSE(cusparseSpSV_bufferSize(*my_handle,CUSPARSE_OPERATION_NON_TRANSPOSE, &alpha,*(cMat->spmvDescr),vecX,vecY,
&alpha,*(cMat->spmvDescr),vecX,vecY, CUSPARSE_BASE_TYPE,
CUSPARSE_BASE_TYPE, CUSPARSE_SPSV_ALG_DEFAULT,
CUSPARSE_SPSV_ALG_DEFAULT, *(cMat->spsvDescr),
*(cMat->spsvDescr), &bfsz));
&bfsz)); if (bfsz > cMat->svbsize) {
if (bfsz > cMat->svbsize) { if (cMat->svbuffer != NULL) {
if (cMat->svbuffer != NULL) { CHECK_CUDA(cudaFree(cMat->svbuffer));
CHECK_CUDA(cudaFree(cMat->svbuffer)); cMat->svbuffer = NULL;
cMat->svbuffer = NULL; }
CHECK_CUDA(cudaMalloc((void **) &(cMat->svbuffer), bfsz));
cMat->svbsize=bfsz;
CHECK_CUSPARSE(cusparseSpSV_analysis(*my_handle,
CUSPARSE_OPERATION_NON_TRANSPOSE,
&alpha,
*(cMat->spmvDescr),
vecX, vecY,
CUSPARSE_BASE_TYPE,
CUSPARSE_SPSV_ALG_DEFAULT,
*(cMat->spsvDescr),
cMat->svbuffer));
}
if (T_CSRGIsNullSvBuffer(cMat)) {
fprintf(stderr,"SpSV_SOLVE NULL spsv-buffer\n");
} }
CHECK_CUDA(cudaMalloc((void **) &(cMat->svbuffer), bfsz));
cMat->svbsize=bfsz;
CHECK_CUSPARSE(cusparseSpSV_analysis(*my_handle,
CUSPARSE_OPERATION_NON_TRANSPOSE,
&alpha,
*(cMat->spmvDescr),
vecX, vecY,
CUSPARSE_BASE_TYPE,
CUSPARSE_SPSV_ALG_DEFAULT,
*(cMat->spsvDescr),
cMat->svbuffer));
}
if (T_CSRGIsNullSvBuffer(cMat)) {
fprintf(stderr,"SpSV_SOLVE NULL spsv-buffer\n");
} }
T_CSRGCreateSpMVDescr(cMat);
CHECK_CUSPARSE(cusparseSpSV_solve(*my_handle,CUSPARSE_OPERATION_NON_TRANSPOSE, CHECK_CUSPARSE(cusparseSpSV_solve(*my_handle,CUSPARSE_OPERATION_NON_TRANSPOSE,
&alpha,*(cMat->spmvDescr),vecX,vecY, &alpha,*(cMat->spmvDescr),vecX,vecY,
CUSPARSE_BASE_TYPE, CUSPARSE_BASE_TYPE,

Loading…
Cancel
Save