1 /* 2 Defines the basic matrix operations for the AIJ (compressed row) 3 matrix storage format using the CUSPARSE library, 4 */ 5 #define PETSC_SKIP_SPINLOCK 6 7 #include <petscconf.h> 8 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 9 #include <../src/mat/impls/sbaij/seq/sbaij.h> 10 #include <../src/vec/vec/impls/dvecimpl.h> 11 #include <petsc/private/vecimpl.h> 12 #undef VecType 13 #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h> 14 15 const char *const MatCUSPARSEStorageFormats[] = {"CSR","ELL","HYB","MatCUSPARSEStorageFormat","MAT_CUSPARSE_",0}; 16 17 static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*); 18 static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*); 19 static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*); 20 21 static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*); 22 static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*); 23 static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*); 24 25 static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat,Vec,Vec); 26 static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec); 27 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec); 28 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec); 29 static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat); 30 static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat,Vec,Vec); 31 static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec); 32 static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec); 33 static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec); 34 35 static PetscErrorCode CsrMatrix_Destroy(CsrMatrix**); 36 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct**); 37 static PetscErrorCode Mat_SeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct**,MatCUSPARSEStorageFormat); 38 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors**); 39 static PetscErrorCode Mat_SeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE**); 40 41 PetscErrorCode MatCUSPARSESetStream(Mat A,const cudaStream_t stream) 42 { 43 cusparseStatus_t stat; 44 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 45 46 PetscFunctionBegin; 47 cusparsestruct->stream = stream; 48 stat = cusparseSetStream(cusparsestruct->handle,cusparsestruct->stream);CHKERRCUDA(stat); 49 PetscFunctionReturn(0); 50 } 51 52 PetscErrorCode MatCUSPARSESetHandle(Mat A,const cusparseHandle_t handle) 53 { 54 cusparseStatus_t stat; 55 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 56 57 PetscFunctionBegin; 58 if (cusparsestruct->handle != handle) { 59 if (cusparsestruct->handle) { 60 stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUDA(stat); 61 } 62 cusparsestruct->handle = handle; 63 } 64 stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 65 PetscFunctionReturn(0); 66 } 67 68 PetscErrorCode MatCUSPARSEClearHandle(Mat A) 69 { 70 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 71 PetscFunctionBegin; 72 if (cusparsestruct->handle) 73 cusparsestruct->handle = 0; 74 PetscFunctionReturn(0); 75 } 76 77 PetscErrorCode MatFactorGetSolverPackage_seqaij_cusparse(Mat A,const MatSolverPackage *type) 78 { 79 PetscFunctionBegin; 80 *type = MATSOLVERCUSPARSE; 81 PetscFunctionReturn(0); 82 } 83 84 /*MC 85 MATSOLVERCUSPARSE = "cusparse" - A matrix type providing triangular solvers for seq matrices 86 on a single GPU of type, seqaijcusparse, aijcusparse, or seqaijcusp, aijcusp. Currently supported 87 algorithms are ILU(k) and ICC(k). Typically, deeper factorizations (larger k) results in poorer 88 performance in the triangular solves. Full LU, and Cholesky decompositions can be solved through the 89 CUSPARSE triangular solve algorithm. However, the performance can be quite poor and thus these 90 algorithms are not recommended. This class does NOT support direct solver operations. 91 92 Level: beginner 93 94 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 95 M*/ 96 97 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat A,MatFactorType ftype,Mat *B) 98 { 99 PetscErrorCode ierr; 100 PetscInt n = A->rmap->n; 101 102 PetscFunctionBegin; 103 ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 104 (*B)->factortype = ftype; 105 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 106 ierr = MatSetType(*B,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 107 108 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 109 ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr); 110 (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJCUSPARSE; 111 (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJCUSPARSE; 112 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 113 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJCUSPARSE; 114 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJCUSPARSE; 115 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported for CUSPARSE Matrix Types"); 116 117 ierr = MatSeqAIJSetPreallocation(*B,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 118 ierr = PetscObjectComposeFunction((PetscObject)(*B),"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_cusparse);CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 PETSC_INTERN PetscErrorCode MatCUSPARSESetFormat_SeqAIJCUSPARSE(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format) 123 { 124 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 125 126 PetscFunctionBegin; 127 #if CUDA_VERSION>=4020 128 switch (op) { 129 case MAT_CUSPARSE_MULT: 130 cusparsestruct->format = format; 131 break; 132 case MAT_CUSPARSE_ALL: 133 cusparsestruct->format = format; 134 break; 135 default: 136 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op); 137 } 138 #else 139 if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ELL (Ellpack) and HYB (Hybrid) storage format require CUDA 4.2 or later."); 140 #endif 141 PetscFunctionReturn(0); 142 } 143 144 /*@ 145 MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular 146 operation. Only the MatMult operation can use different GPU storage formats 147 for MPIAIJCUSPARSE matrices. 148 Not Collective 149 150 Input Parameters: 151 + A - Matrix of type SEQAIJCUSPARSE 152 . op - MatCUSPARSEFormatOperation. SEQAIJCUSPARSE matrices support MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL. MPIAIJCUSPARSE matrices support MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, and MAT_CUSPARSE_ALL. 153 - format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2) 154 155 Output Parameter: 156 157 Level: intermediate 158 159 .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 160 @*/ 161 PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format) 162 { 163 PetscErrorCode ierr; 164 165 PetscFunctionBegin; 166 PetscValidHeaderSpecific(A, MAT_CLASSID,1); 167 ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr); 168 PetscFunctionReturn(0); 169 } 170 171 static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat A) 172 { 173 PetscErrorCode ierr; 174 MatCUSPARSEStorageFormat format; 175 PetscBool flg; 176 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 177 178 PetscFunctionBegin; 179 ierr = PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSPARSE options");CHKERRQ(ierr); 180 ierr = PetscObjectOptionsBegin((PetscObject)A); 181 if (A->factortype==MAT_FACTOR_NONE) { 182 ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV", 183 "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr); 184 if (flg) { 185 ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr); 186 } 187 } 188 ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve", 189 "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr); 190 if (flg) { 191 ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr); 192 } 193 ierr = PetscOptionsEnd();CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 196 } 197 198 static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 199 { 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr); 204 B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE; 205 PetscFunctionReturn(0); 206 } 207 208 static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 209 { 210 PetscErrorCode ierr; 211 212 PetscFunctionBegin; 213 ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr); 214 B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE; 215 PetscFunctionReturn(0); 216 } 217 218 static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info) 219 { 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr); 224 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE; 225 PetscFunctionReturn(0); 226 } 227 228 static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info) 229 { 230 PetscErrorCode ierr; 231 232 PetscFunctionBegin; 233 ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr); 234 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE; 235 PetscFunctionReturn(0); 236 } 237 238 static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A) 239 { 240 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 241 PetscInt n = A->rmap->n; 242 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 243 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 244 cusparseStatus_t stat; 245 const PetscInt *ai = a->i,*aj = a->j,*vi; 246 const MatScalar *aa = a->a,*v; 247 PetscInt *AiLo, *AjLo; 248 PetscScalar *AALo; 249 PetscInt i,nz, nzLower, offset, rowOffset; 250 PetscErrorCode ierr; 251 252 PetscFunctionBegin; 253 if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) { 254 try { 255 /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */ 256 nzLower=n+ai[n]-ai[1]; 257 258 /* Allocate Space for the lower triangular matrix */ 259 ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 260 ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUDA(ierr); 261 ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUDA(ierr); 262 263 /* Fill the lower triangular matrix */ 264 AiLo[0] = (PetscInt) 0; 265 AiLo[n] = nzLower; 266 AjLo[0] = (PetscInt) 0; 267 AALo[0] = (MatScalar) 1.0; 268 v = aa; 269 vi = aj; 270 offset = 1; 271 rowOffset= 1; 272 for (i=1; i<n; i++) { 273 nz = ai[i+1] - ai[i]; 274 /* additional 1 for the term on the diagonal */ 275 AiLo[i] = rowOffset; 276 rowOffset += nz+1; 277 278 ierr = PetscMemcpy(&(AjLo[offset]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr); 279 ierr = PetscMemcpy(&(AALo[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 280 281 offset += nz; 282 AjLo[offset] = (PetscInt) i; 283 AALo[offset] = (MatScalar) 1.0; 284 offset += 1; 285 286 v += nz; 287 vi += nz; 288 } 289 290 /* allocate space for the triangular factor information */ 291 loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 292 293 /* Create the matrix description */ 294 stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat); 295 stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 296 stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 297 stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUDA(stat); 298 stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat); 299 300 /* Create the solve analysis information */ 301 stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat); 302 303 /* set the operation */ 304 loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 305 306 /* set the matrix */ 307 loTriFactor->csrMat = new CsrMatrix; 308 loTriFactor->csrMat->num_rows = n; 309 loTriFactor->csrMat->num_cols = n; 310 loTriFactor->csrMat->num_entries = nzLower; 311 312 loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1); 313 loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1); 314 315 loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower); 316 loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower); 317 318 loTriFactor->csrMat->values = new THRUSTARRAY(nzLower); 319 loTriFactor->csrMat->values->assign(AALo, AALo+nzLower); 320 321 /* perform the solve analysis */ 322 stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp, 323 loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr, 324 loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(), 325 loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat); 326 327 /* assign the pointer. Is this really necessary? */ 328 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor; 329 330 ierr = cudaFreeHost(AiLo);CHKERRCUDA(ierr); 331 ierr = cudaFreeHost(AjLo);CHKERRCUDA(ierr); 332 ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr); 333 } catch(char *ex) { 334 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 335 } 336 } 337 PetscFunctionReturn(0); 338 } 339 340 static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A) 341 { 342 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 343 PetscInt n = A->rmap->n; 344 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 345 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 346 cusparseStatus_t stat; 347 const PetscInt *aj = a->j,*adiag = a->diag,*vi; 348 const MatScalar *aa = a->a,*v; 349 PetscInt *AiUp, *AjUp; 350 PetscScalar *AAUp; 351 PetscInt i,nz, nzUpper, offset; 352 PetscErrorCode ierr; 353 354 PetscFunctionBegin; 355 if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) { 356 try { 357 /* next, figure out the number of nonzeros in the upper triangular matrix. */ 358 nzUpper = adiag[0]-adiag[n]; 359 360 /* Allocate Space for the upper triangular matrix */ 361 ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 362 ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 363 ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 364 365 /* Fill the upper triangular matrix */ 366 AiUp[0]=(PetscInt) 0; 367 AiUp[n]=nzUpper; 368 offset = nzUpper; 369 for (i=n-1; i>=0; i--) { 370 v = aa + adiag[i+1] + 1; 371 vi = aj + adiag[i+1] + 1; 372 373 /* number of elements NOT on the diagonal */ 374 nz = adiag[i] - adiag[i+1]-1; 375 376 /* decrement the offset */ 377 offset -= (nz+1); 378 379 /* first, set the diagonal elements */ 380 AjUp[offset] = (PetscInt) i; 381 AAUp[offset] = (MatScalar)1./v[nz]; 382 AiUp[i] = AiUp[i+1] - (nz+1); 383 384 ierr = PetscMemcpy(&(AjUp[offset+1]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr); 385 ierr = PetscMemcpy(&(AAUp[offset+1]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 386 } 387 388 /* allocate space for the triangular factor information */ 389 upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 390 391 /* Create the matrix description */ 392 stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 393 stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 394 stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 395 stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 396 stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 397 398 /* Create the solve analysis information */ 399 stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 400 401 /* set the operation */ 402 upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 403 404 /* set the matrix */ 405 upTriFactor->csrMat = new CsrMatrix; 406 upTriFactor->csrMat->num_rows = n; 407 upTriFactor->csrMat->num_cols = n; 408 upTriFactor->csrMat->num_entries = nzUpper; 409 410 upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1); 411 upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1); 412 413 upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper); 414 upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper); 415 416 upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper); 417 upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper); 418 419 /* perform the solve analysis */ 420 stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 421 upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 422 upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 423 upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 424 425 /* assign the pointer. Is this really necessary? */ 426 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 427 428 ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 429 ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 430 ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 431 } catch(char *ex) { 432 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 433 } 434 } 435 PetscFunctionReturn(0); 436 } 437 438 static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A) 439 { 440 PetscErrorCode ierr; 441 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 442 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 443 IS isrow = a->row,iscol = a->icol; 444 PetscBool row_identity,col_identity; 445 const PetscInt *r,*c; 446 PetscInt n = A->rmap->n; 447 448 PetscFunctionBegin; 449 ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr); 450 ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr); 451 452 cusparseTriFactors->workVector = new THRUSTARRAY; 453 cusparseTriFactors->workVector->resize(n); 454 cusparseTriFactors->nnz=a->nz; 455 456 A->valid_GPU_matrix = PETSC_CUDA_BOTH; 457 /*lower triangular indices */ 458 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 459 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 460 if (!row_identity) { 461 cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 462 cusparseTriFactors->rpermIndices->assign(r, r+n); 463 } 464 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 465 466 /*upper triangular indices */ 467 ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 468 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 469 if (!col_identity) { 470 cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 471 cusparseTriFactors->cpermIndices->assign(c, c+n); 472 } 473 ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 474 PetscFunctionReturn(0); 475 } 476 477 static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A) 478 { 479 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 480 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 481 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 482 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 483 cusparseStatus_t stat; 484 PetscErrorCode ierr; 485 PetscInt *AiUp, *AjUp; 486 PetscScalar *AAUp; 487 PetscScalar *AALo; 488 PetscInt nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j; 489 Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)A->data; 490 const PetscInt *ai = b->i,*aj = b->j,*vj; 491 const MatScalar *aa = b->a,*v; 492 493 PetscFunctionBegin; 494 if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) { 495 try { 496 /* Allocate Space for the upper triangular matrix */ 497 ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 498 ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 499 ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 500 ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 501 502 /* Fill the upper triangular matrix */ 503 AiUp[0]=(PetscInt) 0; 504 AiUp[n]=nzUpper; 505 offset = 0; 506 for (i=0; i<n; i++) { 507 /* set the pointers */ 508 v = aa + ai[i]; 509 vj = aj + ai[i]; 510 nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */ 511 512 /* first, set the diagonal elements */ 513 AjUp[offset] = (PetscInt) i; 514 AAUp[offset] = (MatScalar)1.0/v[nz]; 515 AiUp[i] = offset; 516 AALo[offset] = (MatScalar)1.0/v[nz]; 517 518 offset+=1; 519 if (nz>0) { 520 ierr = PetscMemcpy(&(AjUp[offset]), vj, nz*sizeof(PetscInt));CHKERRQ(ierr); 521 ierr = PetscMemcpy(&(AAUp[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 522 for (j=offset; j<offset+nz; j++) { 523 AAUp[j] = -AAUp[j]; 524 AALo[j] = AAUp[j]/v[nz]; 525 } 526 offset+=nz; 527 } 528 } 529 530 /* allocate space for the triangular factor information */ 531 upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 532 533 /* Create the matrix description */ 534 stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 535 stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 536 stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 537 stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 538 stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat); 539 540 /* Create the solve analysis information */ 541 stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 542 543 /* set the operation */ 544 upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 545 546 /* set the matrix */ 547 upTriFactor->csrMat = new CsrMatrix; 548 upTriFactor->csrMat->num_rows = A->rmap->n; 549 upTriFactor->csrMat->num_cols = A->cmap->n; 550 upTriFactor->csrMat->num_entries = a->nz; 551 552 upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 553 upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 554 555 upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 556 upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 557 558 upTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 559 upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz); 560 561 /* perform the solve analysis */ 562 stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 563 upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 564 upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 565 upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 566 567 /* assign the pointer. Is this really necessary? */ 568 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 569 570 /* allocate space for the triangular factor information */ 571 loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 572 573 /* Create the matrix description */ 574 stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat); 575 stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 576 stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 577 stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 578 stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 579 580 /* Create the solve analysis information */ 581 stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat); 582 583 /* set the operation */ 584 loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE; 585 586 /* set the matrix */ 587 loTriFactor->csrMat = new CsrMatrix; 588 loTriFactor->csrMat->num_rows = A->rmap->n; 589 loTriFactor->csrMat->num_cols = A->cmap->n; 590 loTriFactor->csrMat->num_entries = a->nz; 591 592 loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 593 loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 594 595 loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 596 loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 597 598 loTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 599 loTriFactor->csrMat->values->assign(AALo, AALo+a->nz); 600 601 /* perform the solve analysis */ 602 stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp, 603 loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr, 604 loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(), 605 loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat); 606 607 /* assign the pointer. Is this really necessary? */ 608 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor; 609 610 A->valid_GPU_matrix = PETSC_CUDA_BOTH; 611 ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 612 ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 613 ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 614 ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr); 615 } catch(char *ex) { 616 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 617 } 618 } 619 PetscFunctionReturn(0); 620 } 621 622 static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A) 623 { 624 PetscErrorCode ierr; 625 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 626 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 627 IS ip = a->row; 628 const PetscInt *rip; 629 PetscBool perm_identity; 630 PetscInt n = A->rmap->n; 631 632 PetscFunctionBegin; 633 ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr); 634 cusparseTriFactors->workVector = new THRUSTARRAY; 635 cusparseTriFactors->workVector->resize(n); 636 cusparseTriFactors->nnz=(a->nz-n)*2 + n; 637 638 /*lower triangular indices */ 639 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 640 ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 641 if (!perm_identity) { 642 cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 643 cusparseTriFactors->rpermIndices->assign(rip, rip+n); 644 cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 645 cusparseTriFactors->cpermIndices->assign(rip, rip+n); 646 } 647 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 648 PetscFunctionReturn(0); 649 } 650 651 static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 652 { 653 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 654 IS isrow = b->row,iscol = b->col; 655 PetscBool row_identity,col_identity; 656 PetscErrorCode ierr; 657 658 PetscFunctionBegin; 659 ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 660 /* determine which version of MatSolve needs to be used. */ 661 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 662 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 663 if (row_identity && col_identity) { 664 B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 665 B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 666 } else { 667 B->ops->solve = MatSolve_SeqAIJCUSPARSE; 668 B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 669 } 670 671 /* get the triangular factors */ 672 ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 673 PetscFunctionReturn(0); 674 } 675 676 static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 677 { 678 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 679 IS ip = b->row; 680 PetscBool perm_identity; 681 PetscErrorCode ierr; 682 683 PetscFunctionBegin; 684 ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 685 686 /* determine which version of MatSolve needs to be used. */ 687 ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 688 if (perm_identity) { 689 B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 690 B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 691 } else { 692 B->ops->solve = MatSolve_SeqAIJCUSPARSE; 693 B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 694 } 695 696 /* get the triangular factors */ 697 ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A) 702 { 703 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 704 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 705 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 706 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 707 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 708 cusparseStatus_t stat; 709 cusparseIndexBase_t indexBase; 710 cusparseMatrixType_t matrixType; 711 cusparseFillMode_t fillMode; 712 cusparseDiagType_t diagType; 713 714 PetscFunctionBegin; 715 716 /*********************************************/ 717 /* Now the Transpose of the Lower Tri Factor */ 718 /*********************************************/ 719 720 /* allocate space for the transpose of the lower triangular factor */ 721 loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 722 723 /* set the matrix descriptors of the lower triangular factor */ 724 matrixType = cusparseGetMatType(loTriFactor->descr); 725 indexBase = cusparseGetMatIndexBase(loTriFactor->descr); 726 fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 727 CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 728 diagType = cusparseGetMatDiagType(loTriFactor->descr); 729 730 /* Create the matrix description */ 731 stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat); 732 stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat); 733 stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat); 734 stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat); 735 stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat); 736 737 /* Create the solve analysis information */ 738 stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat); 739 740 /* set the operation */ 741 loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 742 743 /* allocate GPU space for the CSC of the lower triangular factor*/ 744 loTriFactorT->csrMat = new CsrMatrix; 745 loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows; 746 loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols; 747 loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries; 748 loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1); 749 loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries); 750 loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries); 751 752 /* compute the transpose of the lower triangular factor, i.e. the CSC */ 753 stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows, 754 loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries, 755 loTriFactor->csrMat->values->data().get(), 756 loTriFactor->csrMat->row_offsets->data().get(), 757 loTriFactor->csrMat->column_indices->data().get(), 758 loTriFactorT->csrMat->values->data().get(), 759 loTriFactorT->csrMat->column_indices->data().get(), 760 loTriFactorT->csrMat->row_offsets->data().get(), 761 CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 762 763 /* perform the solve analysis on the transposed matrix */ 764 stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp, 765 loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries, 766 loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(), 767 loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(), 768 loTriFactorT->solveInfo);CHKERRCUDA(stat); 769 770 /* assign the pointer. Is this really necessary? */ 771 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT; 772 773 /*********************************************/ 774 /* Now the Transpose of the Upper Tri Factor */ 775 /*********************************************/ 776 777 /* allocate space for the transpose of the upper triangular factor */ 778 upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 779 780 /* set the matrix descriptors of the upper triangular factor */ 781 matrixType = cusparseGetMatType(upTriFactor->descr); 782 indexBase = cusparseGetMatIndexBase(upTriFactor->descr); 783 fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 784 CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 785 diagType = cusparseGetMatDiagType(upTriFactor->descr); 786 787 /* Create the matrix description */ 788 stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat); 789 stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat); 790 stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat); 791 stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat); 792 stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat); 793 794 /* Create the solve analysis information */ 795 stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat); 796 797 /* set the operation */ 798 upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 799 800 /* allocate GPU space for the CSC of the upper triangular factor*/ 801 upTriFactorT->csrMat = new CsrMatrix; 802 upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows; 803 upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols; 804 upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries; 805 upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1); 806 upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries); 807 upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries); 808 809 /* compute the transpose of the upper triangular factor, i.e. the CSC */ 810 stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows, 811 upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries, 812 upTriFactor->csrMat->values->data().get(), 813 upTriFactor->csrMat->row_offsets->data().get(), 814 upTriFactor->csrMat->column_indices->data().get(), 815 upTriFactorT->csrMat->values->data().get(), 816 upTriFactorT->csrMat->column_indices->data().get(), 817 upTriFactorT->csrMat->row_offsets->data().get(), 818 CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 819 820 /* perform the solve analysis on the transposed matrix */ 821 stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp, 822 upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries, 823 upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(), 824 upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(), 825 upTriFactorT->solveInfo);CHKERRCUDA(stat); 826 827 /* assign the pointer. Is this really necessary? */ 828 ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT; 829 PetscFunctionReturn(0); 830 } 831 832 static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A) 833 { 834 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 835 Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 836 Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 837 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 838 cusparseStatus_t stat; 839 cusparseIndexBase_t indexBase; 840 cudaError_t err; 841 842 PetscFunctionBegin; 843 844 /* allocate space for the triangular factor information */ 845 matstructT = new Mat_SeqAIJCUSPARSEMultStruct; 846 stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat); 847 indexBase = cusparseGetMatIndexBase(matstruct->descr); 848 stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat); 849 stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 850 851 /* set alpha and beta */ 852 err = cudaMalloc((void **)&(matstructT->alpha),sizeof(PetscScalar));CHKERRCUDA(err); 853 err = cudaMemcpy(matstructT->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 854 err = cudaMalloc((void **)&(matstructT->beta),sizeof(PetscScalar));CHKERRCUDA(err); 855 err = cudaMemcpy(matstructT->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 856 stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 857 858 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 859 CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 860 CsrMatrix *matrixT= new CsrMatrix; 861 matrixT->num_rows = A->rmap->n; 862 matrixT->num_cols = A->cmap->n; 863 matrixT->num_entries = a->nz; 864 matrixT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 865 matrixT->column_indices = new THRUSTINTARRAY32(a->nz); 866 matrixT->values = new THRUSTARRAY(a->nz); 867 868 /* compute the transpose of the upper triangular factor, i.e. the CSC */ 869 indexBase = cusparseGetMatIndexBase(matstruct->descr); 870 stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows, 871 matrix->num_cols, matrix->num_entries, 872 matrix->values->data().get(), 873 matrix->row_offsets->data().get(), 874 matrix->column_indices->data().get(), 875 matrixT->values->data().get(), 876 matrixT->column_indices->data().get(), 877 matrixT->row_offsets->data().get(), 878 CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 879 880 /* assign the pointer */ 881 matstructT->mat = matrixT; 882 883 } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 884 #if CUDA_VERSION>=5000 885 /* First convert HYB to CSR */ 886 CsrMatrix *temp= new CsrMatrix; 887 temp->num_rows = A->rmap->n; 888 temp->num_cols = A->cmap->n; 889 temp->num_entries = a->nz; 890 temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 891 temp->column_indices = new THRUSTINTARRAY32(a->nz); 892 temp->values = new THRUSTARRAY(a->nz); 893 894 895 stat = cusparse_hyb2csr(cusparsestruct->handle, 896 matstruct->descr, (cusparseHybMat_t)matstruct->mat, 897 temp->values->data().get(), 898 temp->row_offsets->data().get(), 899 temp->column_indices->data().get());CHKERRCUDA(stat); 900 901 /* Next, convert CSR to CSC (i.e. the matrix transpose) */ 902 CsrMatrix *tempT= new CsrMatrix; 903 tempT->num_rows = A->rmap->n; 904 tempT->num_cols = A->cmap->n; 905 tempT->num_entries = a->nz; 906 tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 907 tempT->column_indices = new THRUSTINTARRAY32(a->nz); 908 tempT->values = new THRUSTARRAY(a->nz); 909 910 stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows, 911 temp->num_cols, temp->num_entries, 912 temp->values->data().get(), 913 temp->row_offsets->data().get(), 914 temp->column_indices->data().get(), 915 tempT->values->data().get(), 916 tempT->column_indices->data().get(), 917 tempT->row_offsets->data().get(), 918 CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 919 920 /* Last, convert CSC to HYB */ 921 cusparseHybMat_t hybMat; 922 stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 923 cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 924 CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 925 stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n, 926 matstructT->descr, tempT->values->data().get(), 927 tempT->row_offsets->data().get(), 928 tempT->column_indices->data().get(), 929 hybMat, 0, partition);CHKERRCUDA(stat); 930 931 /* assign the pointer */ 932 matstructT->mat = hybMat; 933 934 /* delete temporaries */ 935 if (tempT) { 936 if (tempT->values) delete (THRUSTARRAY*) tempT->values; 937 if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices; 938 if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets; 939 delete (CsrMatrix*) tempT; 940 } 941 if (temp) { 942 if (temp->values) delete (THRUSTARRAY*) temp->values; 943 if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices; 944 if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets; 945 delete (CsrMatrix*) temp; 946 } 947 #else 948 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ELL (Ellpack) and HYB (Hybrid) storage format for the Matrix Transpose (in MatMultTranspose) require CUDA 5.0 or later."); 949 #endif 950 } 951 /* assign the compressed row indices */ 952 matstructT->cprowIndices = new THRUSTINTARRAY; 953 954 /* assign the pointer */ 955 ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT; 956 PetscFunctionReturn(0); 957 } 958 959 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 960 { 961 PetscInt n = xx->map->n; 962 const PetscScalar *barray; 963 PetscScalar *xarray; 964 thrust::device_ptr<const PetscScalar> bGPU; 965 thrust::device_ptr<PetscScalar> xGPU; 966 cusparseStatus_t stat; 967 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 968 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 969 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 970 THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 971 PetscErrorCode ierr; 972 973 PetscFunctionBegin; 974 /* Analyze the matrix and create the transpose ... on the fly */ 975 if (!loTriFactorT && !upTriFactorT) { 976 ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 977 loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 978 upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 979 } 980 981 /* Get the GPU pointers */ 982 ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 983 ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 984 xGPU = thrust::device_pointer_cast(xarray); 985 bGPU = thrust::device_pointer_cast(barray); 986 987 /* First, reorder with the row permutation */ 988 thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 989 thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()), 990 xGPU); 991 992 /* First, solve U */ 993 stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 994 upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr, 995 upTriFactorT->csrMat->values->data().get(), 996 upTriFactorT->csrMat->row_offsets->data().get(), 997 upTriFactorT->csrMat->column_indices->data().get(), 998 upTriFactorT->solveInfo, 999 xarray, tempGPU->data().get());CHKERRCUDA(stat); 1000 1001 /* Then, solve L */ 1002 stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 1003 loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr, 1004 loTriFactorT->csrMat->values->data().get(), 1005 loTriFactorT->csrMat->row_offsets->data().get(), 1006 loTriFactorT->csrMat->column_indices->data().get(), 1007 loTriFactorT->solveInfo, 1008 tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1009 1010 /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1011 thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1012 thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()), 1013 tempGPU->begin()); 1014 1015 /* Copy the temporary to the full solution. */ 1016 thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 1017 1018 /* restore */ 1019 ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1020 ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1021 ierr = WaitForGPU();CHKERRCUDA(ierr); 1022 1023 ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1024 PetscFunctionReturn(0); 1025 } 1026 1027 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 1028 { 1029 const PetscScalar *barray; 1030 PetscScalar *xarray; 1031 cusparseStatus_t stat; 1032 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1033 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1034 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1035 THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1036 PetscErrorCode ierr; 1037 1038 PetscFunctionBegin; 1039 /* Analyze the matrix and create the transpose ... on the fly */ 1040 if (!loTriFactorT && !upTriFactorT) { 1041 ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 1042 loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1043 upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1044 } 1045 1046 /* Get the GPU pointers */ 1047 ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1048 ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1049 1050 /* First, solve U */ 1051 stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 1052 upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr, 1053 upTriFactorT->csrMat->values->data().get(), 1054 upTriFactorT->csrMat->row_offsets->data().get(), 1055 upTriFactorT->csrMat->column_indices->data().get(), 1056 upTriFactorT->solveInfo, 1057 barray, tempGPU->data().get());CHKERRCUDA(stat); 1058 1059 /* Then, solve L */ 1060 stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 1061 loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr, 1062 loTriFactorT->csrMat->values->data().get(), 1063 loTriFactorT->csrMat->row_offsets->data().get(), 1064 loTriFactorT->csrMat->column_indices->data().get(), 1065 loTriFactorT->solveInfo, 1066 tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1067 1068 /* restore */ 1069 ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1070 ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1071 ierr = WaitForGPU();CHKERRCUDA(ierr); 1072 ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1073 PetscFunctionReturn(0); 1074 } 1075 1076 static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 1077 { 1078 const PetscScalar *barray; 1079 PetscScalar *xarray; 1080 thrust::device_ptr<const PetscScalar> bGPU; 1081 thrust::device_ptr<PetscScalar> xGPU; 1082 cusparseStatus_t stat; 1083 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1084 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1085 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1086 THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1087 PetscErrorCode ierr; 1088 VecType t; 1089 PetscBool flg; 1090 1091 PetscFunctionBegin; 1092 ierr = VecGetType(bb,&t);CHKERRQ(ierr); 1093 ierr = PetscStrcmp(t,VECSEQCUDA,&flg);CHKERRQ(ierr); 1094 if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector of type %s passed into MatSolve_SeqAIJCUSPARSE (Arg #2). Can only deal with %s\n.",t,VECSEQCUDA); 1095 ierr = VecGetType(xx,&t);CHKERRQ(ierr); 1096 ierr = PetscStrcmp(t,VECSEQCUDA,&flg);CHKERRQ(ierr); 1097 if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector of type %s passed into MatSolve_SeqAIJCUSPARSE (Arg #3). Can only deal with %s\n.",t,VECSEQCUDA); 1098 1099 /* Get the GPU pointers */ 1100 ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1101 ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1102 xGPU = thrust::device_pointer_cast(xarray); 1103 bGPU = thrust::device_pointer_cast(barray); 1104 1105 /* First, reorder with the row permutation */ 1106 thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 1107 thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()), 1108 xGPU); 1109 1110 /* Next, solve L */ 1111 stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 1112 loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr, 1113 loTriFactor->csrMat->values->data().get(), 1114 loTriFactor->csrMat->row_offsets->data().get(), 1115 loTriFactor->csrMat->column_indices->data().get(), 1116 loTriFactor->solveInfo, 1117 xarray, tempGPU->data().get());CHKERRCUDA(stat); 1118 1119 /* Then, solve U */ 1120 stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 1121 upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr, 1122 upTriFactor->csrMat->values->data().get(), 1123 upTriFactor->csrMat->row_offsets->data().get(), 1124 upTriFactor->csrMat->column_indices->data().get(), 1125 upTriFactor->solveInfo, 1126 tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1127 1128 /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1129 thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1130 thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->end()), 1131 tempGPU->begin()); 1132 1133 /* Copy the temporary to the full solution. */ 1134 thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 1135 1136 ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1137 ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1138 ierr = WaitForGPU();CHKERRCUDA(ierr); 1139 ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1140 PetscFunctionReturn(0); 1141 } 1142 1143 static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 1144 { 1145 const PetscScalar *barray; 1146 PetscScalar *xarray; 1147 cusparseStatus_t stat; 1148 Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1149 Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1150 Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1151 THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1152 PetscErrorCode ierr; 1153 1154 PetscFunctionBegin; 1155 /* Get the GPU pointers */ 1156 ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1157 ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1158 1159 /* First, solve L */ 1160 stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 1161 loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr, 1162 loTriFactor->csrMat->values->data().get(), 1163 loTriFactor->csrMat->row_offsets->data().get(), 1164 loTriFactor->csrMat->column_indices->data().get(), 1165 loTriFactor->solveInfo, 1166 barray, tempGPU->data().get());CHKERRCUDA(stat); 1167 1168 /* Next, solve U */ 1169 stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 1170 upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr, 1171 upTriFactor->csrMat->values->data().get(), 1172 upTriFactor->csrMat->row_offsets->data().get(), 1173 upTriFactor->csrMat->column_indices->data().get(), 1174 upTriFactor->solveInfo, 1175 tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1176 1177 ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1178 ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1179 ierr = WaitForGPU();CHKERRCUDA(ierr); 1180 ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1181 PetscFunctionReturn(0); 1182 } 1183 1184 static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A) 1185 { 1186 1187 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1188 Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 1189 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1190 PetscInt m = A->rmap->n,*ii,*ridx; 1191 PetscErrorCode ierr; 1192 cusparseStatus_t stat; 1193 cudaError_t err; 1194 1195 PetscFunctionBegin; 1196 if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) { 1197 ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 1198 if (A->assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) { 1199 CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 1200 /* copy values only */ 1201 matrix->values->assign(a->a, a->a+a->nz); 1202 } else { 1203 Mat_SeqAIJCUSPARSEMultStruct_Destroy(&matstruct,cusparsestruct->format); 1204 try { 1205 cusparsestruct->nonzerorow=0; 1206 for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0); 1207 1208 if (a->compressedrow.use) { 1209 m = a->compressedrow.nrows; 1210 ii = a->compressedrow.i; 1211 ridx = a->compressedrow.rindex; 1212 } else { 1213 /* Forcing compressed row on the GPU */ 1214 int k=0; 1215 ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr); 1216 ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr); 1217 ii[0]=0; 1218 for (int j = 0; j<m; j++) { 1219 if ((a->i[j+1]-a->i[j])>0) { 1220 ii[k] = a->i[j]; 1221 ridx[k]= j; 1222 k++; 1223 } 1224 } 1225 ii[cusparsestruct->nonzerorow] = a->nz; 1226 m = cusparsestruct->nonzerorow; 1227 } 1228 1229 /* allocate space for the triangular factor information */ 1230 matstruct = new Mat_SeqAIJCUSPARSEMultStruct; 1231 stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat); 1232 stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 1233 stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 1234 1235 err = cudaMalloc((void **)&(matstruct->alpha),sizeof(PetscScalar));CHKERRCUDA(err); 1236 err = cudaMemcpy(matstruct->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 1237 err = cudaMalloc((void **)&(matstruct->beta),sizeof(PetscScalar));CHKERRCUDA(err); 1238 err = cudaMemcpy(matstruct->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 1239 stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 1240 1241 /* Build a hybrid/ellpack matrix if this option is chosen for the storage */ 1242 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1243 /* set the matrix */ 1244 CsrMatrix *matrix= new CsrMatrix; 1245 matrix->num_rows = m; 1246 matrix->num_cols = A->cmap->n; 1247 matrix->num_entries = a->nz; 1248 matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1249 matrix->row_offsets->assign(ii, ii + m+1); 1250 1251 matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1252 matrix->column_indices->assign(a->j, a->j+a->nz); 1253 1254 matrix->values = new THRUSTARRAY(a->nz); 1255 matrix->values->assign(a->a, a->a+a->nz); 1256 1257 /* assign the pointer */ 1258 matstruct->mat = matrix; 1259 1260 } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 1261 #if CUDA_VERSION>=4020 1262 CsrMatrix *matrix= new CsrMatrix; 1263 matrix->num_rows = m; 1264 matrix->num_cols = A->cmap->n; 1265 matrix->num_entries = a->nz; 1266 matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1267 matrix->row_offsets->assign(ii, ii + m+1); 1268 1269 matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1270 matrix->column_indices->assign(a->j, a->j+a->nz); 1271 1272 matrix->values = new THRUSTARRAY(a->nz); 1273 matrix->values->assign(a->a, a->a+a->nz); 1274 1275 cusparseHybMat_t hybMat; 1276 stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 1277 cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 1278 CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 1279 stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols, 1280 matstruct->descr, matrix->values->data().get(), 1281 matrix->row_offsets->data().get(), 1282 matrix->column_indices->data().get(), 1283 hybMat, 0, partition);CHKERRCUDA(stat); 1284 /* assign the pointer */ 1285 matstruct->mat = hybMat; 1286 1287 if (matrix) { 1288 if (matrix->values) delete (THRUSTARRAY*)matrix->values; 1289 if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices; 1290 if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets; 1291 delete (CsrMatrix*)matrix; 1292 } 1293 #endif 1294 } 1295 1296 /* assign the compressed row indices */ 1297 matstruct->cprowIndices = new THRUSTINTARRAY(m); 1298 matstruct->cprowIndices->assign(ridx,ridx+m); 1299 1300 /* assign the pointer */ 1301 cusparsestruct->mat = matstruct; 1302 1303 if (!a->compressedrow.use) { 1304 ierr = PetscFree(ii);CHKERRQ(ierr); 1305 ierr = PetscFree(ridx);CHKERRQ(ierr); 1306 } 1307 cusparsestruct->workVector = new THRUSTARRAY; 1308 cusparsestruct->workVector->resize(m); 1309 } catch(char *ex) { 1310 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 1311 } 1312 cusparsestruct->nonzerostate = A->nonzerostate; 1313 } 1314 ierr = WaitForGPU();CHKERRCUDA(ierr); 1315 A->valid_GPU_matrix = PETSC_CUDA_BOTH; 1316 ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 1317 } 1318 PetscFunctionReturn(0); 1319 } 1320 1321 static PetscErrorCode MatCreateVecs_SeqAIJCUSPARSE(Mat mat, Vec *right, Vec *left) 1322 { 1323 PetscErrorCode ierr; 1324 PetscInt rbs,cbs; 1325 1326 PetscFunctionBegin; 1327 ierr = MatGetBlockSizes(mat,&rbs,&cbs);CHKERRQ(ierr); 1328 if (right) { 1329 ierr = VecCreate(PetscObjectComm((PetscObject)mat),right);CHKERRQ(ierr); 1330 ierr = VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);CHKERRQ(ierr); 1331 ierr = VecSetBlockSize(*right,cbs);CHKERRQ(ierr); 1332 ierr = VecSetType(*right,VECSEQCUDA);CHKERRQ(ierr); 1333 ierr = PetscLayoutReference(mat->cmap,&(*right)->map);CHKERRQ(ierr); 1334 } 1335 if (left) { 1336 ierr = VecCreate(PetscObjectComm((PetscObject)mat),left);CHKERRQ(ierr); 1337 ierr = VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);CHKERRQ(ierr); 1338 ierr = VecSetBlockSize(*left,rbs);CHKERRQ(ierr); 1339 ierr = VecSetType(*left,VECSEQCUDA);CHKERRQ(ierr); 1340 ierr = PetscLayoutReference(mat->rmap,&(*left)->map);CHKERRQ(ierr); 1341 } 1342 PetscFunctionReturn(0); 1343 } 1344 1345 struct VecCUDAPlusEquals 1346 { 1347 template <typename Tuple> 1348 __host__ __device__ 1349 void operator()(Tuple t) 1350 { 1351 thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t); 1352 } 1353 }; 1354 1355 static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 1356 { 1357 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1358 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1359 Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 1360 const PetscScalar *xarray; 1361 PetscScalar *yarray; 1362 PetscErrorCode ierr; 1363 cusparseStatus_t stat; 1364 1365 PetscFunctionBegin; 1366 /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 1367 ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1368 ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1369 ierr = VecSet(yy,0);CHKERRQ(ierr); 1370 ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr); 1371 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1372 CsrMatrix *mat = (CsrMatrix*)matstruct->mat; 1373 stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1374 mat->num_rows, mat->num_cols, mat->num_entries, 1375 matstruct->alpha, matstruct->descr, mat->values->data().get(), mat->row_offsets->data().get(), 1376 mat->column_indices->data().get(), xarray, matstruct->beta, 1377 yarray);CHKERRCUDA(stat); 1378 } else { 1379 #if CUDA_VERSION>=4020 1380 cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat; 1381 stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1382 matstruct->alpha, matstruct->descr, hybMat, 1383 xarray, matstruct->beta, 1384 yarray);CHKERRCUDA(stat); 1385 #endif 1386 } 1387 ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1388 ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr); 1389 if (!cusparsestruct->stream) { 1390 ierr = WaitForGPU();CHKERRCUDA(ierr); 1391 } 1392 ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr); 1393 PetscFunctionReturn(0); 1394 } 1395 1396 static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 1397 { 1398 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1399 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1400 Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1401 const PetscScalar *xarray; 1402 PetscScalar *yarray; 1403 PetscErrorCode ierr; 1404 cusparseStatus_t stat; 1405 1406 PetscFunctionBegin; 1407 /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 1408 ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1409 if (!matstructT) { 1410 ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1411 matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1412 } 1413 ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1414 ierr = VecSet(yy,0);CHKERRQ(ierr); 1415 ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr); 1416 1417 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1418 CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1419 stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1420 mat->num_rows, mat->num_cols, 1421 mat->num_entries, matstructT->alpha, matstructT->descr, 1422 mat->values->data().get(), mat->row_offsets->data().get(), 1423 mat->column_indices->data().get(), xarray, matstructT->beta, 1424 yarray);CHKERRCUDA(stat); 1425 } else { 1426 #if CUDA_VERSION>=4020 1427 cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1428 stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1429 matstructT->alpha, matstructT->descr, hybMat, 1430 xarray, matstructT->beta, 1431 yarray);CHKERRCUDA(stat); 1432 #endif 1433 } 1434 ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1435 ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr); 1436 if (!cusparsestruct->stream) { 1437 ierr = WaitForGPU();CHKERRCUDA(ierr); 1438 } 1439 ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr); 1440 PetscFunctionReturn(0); 1441 } 1442 1443 1444 static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 1445 { 1446 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1447 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1448 Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 1449 thrust::device_ptr<PetscScalar> zptr; 1450 const PetscScalar *xarray; 1451 PetscScalar *zarray; 1452 PetscErrorCode ierr; 1453 cusparseStatus_t stat; 1454 1455 PetscFunctionBegin; 1456 /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 1457 ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1458 try { 1459 ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 1460 ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1461 ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1462 zptr = thrust::device_pointer_cast(zarray); 1463 1464 /* multiply add */ 1465 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1466 CsrMatrix *mat = (CsrMatrix*)matstruct->mat; 1467 /* here we need to be careful to set the number of rows in the multiply to the 1468 number of compressed rows in the matrix ... which is equivalent to the 1469 size of the workVector */ 1470 stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1471 mat->num_rows, mat->num_cols, 1472 mat->num_entries, matstruct->alpha, matstruct->descr, 1473 mat->values->data().get(), mat->row_offsets->data().get(), 1474 mat->column_indices->data().get(), xarray, matstruct->beta, 1475 cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1476 } else { 1477 #if CUDA_VERSION>=4020 1478 cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat; 1479 if (cusparsestruct->workVector->size()) { 1480 stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1481 matstruct->alpha, matstruct->descr, hybMat, 1482 xarray, matstruct->beta, 1483 cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1484 } 1485 #endif 1486 } 1487 1488 /* scatter the data from the temporary into the full vector with a += operation */ 1489 thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))), 1490 thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(), 1491 VecCUDAPlusEquals()); 1492 ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1493 ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1494 1495 } catch(char *ex) { 1496 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 1497 } 1498 ierr = WaitForGPU();CHKERRCUDA(ierr); 1499 ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1500 PetscFunctionReturn(0); 1501 } 1502 1503 static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 1504 { 1505 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1506 Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1507 Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1508 thrust::device_ptr<PetscScalar> zptr; 1509 const PetscScalar *xarray; 1510 PetscScalar *zarray; 1511 PetscErrorCode ierr; 1512 cusparseStatus_t stat; 1513 1514 PetscFunctionBegin; 1515 /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 1516 ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1517 if (!matstructT) { 1518 ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1519 matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1520 } 1521 1522 try { 1523 ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 1524 ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1525 ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1526 zptr = thrust::device_pointer_cast(zarray); 1527 1528 /* multiply add with matrix transpose */ 1529 if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1530 CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1531 /* here we need to be careful to set the number of rows in the multiply to the 1532 number of compressed rows in the matrix ... which is equivalent to the 1533 size of the workVector */ 1534 stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1535 mat->num_rows, mat->num_cols, 1536 mat->num_entries, matstructT->alpha, matstructT->descr, 1537 mat->values->data().get(), mat->row_offsets->data().get(), 1538 mat->column_indices->data().get(), xarray, matstructT->beta, 1539 cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1540 } else { 1541 #if CUDA_VERSION>=4020 1542 cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1543 if (cusparsestruct->workVector->size()) { 1544 stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1545 matstructT->alpha, matstructT->descr, hybMat, 1546 xarray, matstructT->beta, 1547 cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1548 } 1549 #endif 1550 } 1551 1552 /* scatter the data from the temporary into the full vector with a += operation */ 1553 thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))), 1554 thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))) + cusparsestruct->workVector->size(), 1555 VecCUDAPlusEquals()); 1556 1557 ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1558 ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1559 1560 } catch(char *ex) { 1561 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 1562 } 1563 ierr = WaitForGPU();CHKERRCUDA(ierr); 1564 ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1565 PetscFunctionReturn(0); 1566 } 1567 1568 static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode) 1569 { 1570 PetscErrorCode ierr; 1571 1572 PetscFunctionBegin; 1573 ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr); 1574 if (A->factortype==MAT_FACTOR_NONE) { 1575 ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1576 } 1577 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1578 A->ops->mult = MatMult_SeqAIJCUSPARSE; 1579 A->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1580 A->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1581 A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 1582 PetscFunctionReturn(0); 1583 } 1584 1585 /* --------------------------------------------------------------------------------*/ 1586 /*@ 1587 MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format 1588 (the default parallel PETSc format). This matrix will ultimately pushed down 1589 to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix 1590 assembly performance the user should preallocate the matrix storage by setting 1591 the parameter nz (or the array nnz). By setting these parameters accurately, 1592 performance during matrix assembly can be increased by more than a factor of 50. 1593 1594 Collective on MPI_Comm 1595 1596 Input Parameters: 1597 + comm - MPI communicator, set to PETSC_COMM_SELF 1598 . m - number of rows 1599 . n - number of columns 1600 . nz - number of nonzeros per row (same for all rows) 1601 - nnz - array containing the number of nonzeros in the various rows 1602 (possibly different for each row) or NULL 1603 1604 Output Parameter: 1605 . A - the matrix 1606 1607 It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 1608 MatXXXXSetPreallocation() paradgm instead of this routine directly. 1609 [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 1610 1611 Notes: 1612 If nnz is given then nz is ignored 1613 1614 The AIJ format (also called the Yale sparse matrix format or 1615 compressed row storage), is fully compatible with standard Fortran 77 1616 storage. That is, the stored row and column indices can begin at 1617 either one (as in Fortran) or zero. See the users' manual for details. 1618 1619 Specify the preallocated storage with either nz or nnz (not both). 1620 Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 1621 allocation. For large problems you MUST preallocate memory or you 1622 will get TERRIBLE performance, see the users' manual chapter on matrices. 1623 1624 By default, this format uses inodes (identical nodes) when possible, to 1625 improve numerical efficiency of matrix-vector products and solves. We 1626 search for consecutive rows with the same nonzero structure, thereby 1627 reusing matrix information to achieve increased efficiency. 1628 1629 Level: intermediate 1630 1631 .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE 1632 @*/ 1633 PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1634 { 1635 PetscErrorCode ierr; 1636 1637 PetscFunctionBegin; 1638 ierr = MatCreate(comm,A);CHKERRQ(ierr); 1639 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1640 ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 1641 ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr); 1642 PetscFunctionReturn(0); 1643 } 1644 1645 static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A) 1646 { 1647 PetscErrorCode ierr; 1648 1649 PetscFunctionBegin; 1650 if (A->factortype==MAT_FACTOR_NONE) { 1651 if (A->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) { 1652 ierr = Mat_SeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr); 1653 } 1654 } else { 1655 ierr = Mat_SeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr); 1656 } 1657 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 1658 PetscFunctionReturn(0); 1659 } 1660 1661 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B) 1662 { 1663 PetscErrorCode ierr; 1664 cusparseStatus_t stat; 1665 cusparseHandle_t handle=0; 1666 1667 PetscFunctionBegin; 1668 ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr); 1669 if (B->factortype==MAT_FACTOR_NONE) { 1670 /* you cannot check the inode.use flag here since the matrix was just created. 1671 now build a GPU matrix data structure */ 1672 B->spptr = new Mat_SeqAIJCUSPARSE; 1673 ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat = 0; 1674 ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0; 1675 ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector = 0; 1676 ((Mat_SeqAIJCUSPARSE*)B->spptr)->format = MAT_CUSPARSE_CSR; 1677 ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 1678 ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = 0; 1679 stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1680 ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = handle; 1681 ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 1682 } else { 1683 /* NEXT, set the pointers to the triangular factors */ 1684 B->spptr = new Mat_SeqAIJCUSPARSETriFactors; 1685 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr = 0; 1686 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr = 0; 1687 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0; 1688 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0; 1689 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices = 0; 1690 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices = 0; 1691 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector = 0; 1692 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = 0; 1693 stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1694 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = handle; 1695 ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz = 0; 1696 } 1697 1698 B->ops->assemblyend = MatAssemblyEnd_SeqAIJCUSPARSE; 1699 B->ops->destroy = MatDestroy_SeqAIJCUSPARSE; 1700 B->ops->getvecs = MatCreateVecs_SeqAIJCUSPARSE; 1701 B->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSPARSE; 1702 B->ops->mult = MatMult_SeqAIJCUSPARSE; 1703 B->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1704 B->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1705 B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 1706 1707 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 1708 1709 B->valid_GPU_matrix = PETSC_CUDA_UNALLOCATED; 1710 1711 ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr); 1712 PetscFunctionReturn(0); 1713 } 1714 1715 /*M 1716 MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices. 1717 1718 A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either 1719 CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later. 1720 All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library. 1721 1722 Options Database Keys: 1723 + -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions() 1724 . -mat_cusparse_storage_format csr - sets the storage format of matrices (for MatMult and factors in MatSolve) during a call to MatSetFromOptions(). Other options include ell (ellpack) or hyb (hybrid). 1725 . -mat_cusparse_mult_storage_format csr - sets the storage format of matrices (for MatMult) during a call to MatSetFromOptions(). Other options include ell (ellpack) or hyb (hybrid). 1726 1727 Level: beginner 1728 1729 .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 1730 M*/ 1731 1732 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*); 1733 1734 1735 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_CUSPARSE(void) 1736 { 1737 PetscErrorCode ierr; 1738 1739 PetscFunctionBegin; 1740 ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 1741 ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 1742 ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 1743 ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 1744 PetscFunctionReturn(0); 1745 } 1746 1747 1748 static PetscErrorCode Mat_SeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct) 1749 { 1750 cusparseStatus_t stat; 1751 cusparseHandle_t handle; 1752 1753 PetscFunctionBegin; 1754 if (*cusparsestruct) { 1755 Mat_SeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format); 1756 Mat_SeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format); 1757 delete (*cusparsestruct)->workVector; 1758 if (handle = (*cusparsestruct)->handle) { 1759 stat = cusparseDestroy(handle);CHKERRCUDA(stat); 1760 } 1761 delete *cusparsestruct; 1762 *cusparsestruct = 0; 1763 } 1764 PetscFunctionReturn(0); 1765 } 1766 1767 static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat) 1768 { 1769 PetscFunctionBegin; 1770 if (*mat) { 1771 delete (*mat)->values; 1772 delete (*mat)->column_indices; 1773 delete (*mat)->row_offsets; 1774 delete *mat; 1775 *mat = 0; 1776 } 1777 PetscFunctionReturn(0); 1778 } 1779 1780 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor) 1781 { 1782 cusparseStatus_t stat; 1783 PetscErrorCode ierr; 1784 1785 PetscFunctionBegin; 1786 if (*trifactor) { 1787 if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); } 1788 if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); } 1789 ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr); 1790 delete *trifactor; 1791 *trifactor = 0; 1792 } 1793 PetscFunctionReturn(0); 1794 } 1795 1796 static PetscErrorCode Mat_SeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format) 1797 { 1798 CsrMatrix *mat; 1799 cusparseStatus_t stat; 1800 cudaError_t err; 1801 1802 PetscFunctionBegin; 1803 if (*matstruct) { 1804 if ((*matstruct)->mat) { 1805 if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) { 1806 cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat; 1807 stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat); 1808 } else { 1809 mat = (CsrMatrix*)(*matstruct)->mat; 1810 CsrMatrix_Destroy(&mat); 1811 } 1812 } 1813 if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); } 1814 delete (*matstruct)->cprowIndices; 1815 if ((*matstruct)->alpha) { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); } 1816 if ((*matstruct)->beta) { err=cudaFree((*matstruct)->beta);CHKERRCUDA(err); } 1817 delete *matstruct; 1818 *matstruct = 0; 1819 } 1820 PetscFunctionReturn(0); 1821 } 1822 1823 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors) 1824 { 1825 cusparseHandle_t handle; 1826 cusparseStatus_t stat; 1827 1828 PetscFunctionBegin; 1829 if (*trifactors) { 1830 Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->loTriFactorPtr); 1831 Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->upTriFactorPtr); 1832 Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose); 1833 Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose); 1834 delete (*trifactors)->rpermIndices; 1835 delete (*trifactors)->cpermIndices; 1836 delete (*trifactors)->workVector; 1837 if (handle = (*trifactors)->handle) { 1838 stat = cusparseDestroy(handle);CHKERRCUDA(stat); 1839 } 1840 delete *trifactors; 1841 *trifactors = 0; 1842 } 1843 PetscFunctionReturn(0); 1844 } 1845 1846