19ae82921SPaul Mullowney /* 29ae82921SPaul Mullowney Defines the basic matrix operations for the AIJ (compressed row) 3bc3f50f2SPaul Mullowney matrix storage format using the CUSPARSE library, 49ae82921SPaul Mullowney */ 5dced61a5SBarry Smith #define PETSC_SKIP_SPINLOCK 69ae82921SPaul Mullowney 73d13b8fdSMatthew G. Knepley #include <petscconf.h> 83d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 9087f3262SPaul Mullowney #include <../src/mat/impls/sbaij/seq/sbaij.h> 103d13b8fdSMatthew G. Knepley #include <../src/vec/vec/impls/dvecimpl.h> 11af0996ceSBarry Smith #include <petsc/private/vecimpl.h> 129ae82921SPaul Mullowney #undef VecType 133d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h> 14bc3f50f2SPaul Mullowney 15e057df02SPaul Mullowney const char *const MatCUSPARSEStorageFormats[] = {"CSR","ELL","HYB","MatCUSPARSEStorageFormat","MAT_CUSPARSE_",0}; 169ae82921SPaul Mullowney 17087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*); 18087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*); 19087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*); 20087f3262SPaul Mullowney 216fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*); 226fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*); 236fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*); 24087f3262SPaul Mullowney 256fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat,Vec,Vec); 266fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec); 276fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec); 286fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec); 294416b707SBarry Smith static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat); 306fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat,Vec,Vec); 316fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec); 326fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec); 336fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec); 349ae82921SPaul Mullowney 357f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix**); 36470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct**); 37470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct**,MatCUSPARSEStorageFormat); 38470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors**); 39470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE**); 407f756511SDominic Meiser 41b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetStream(Mat A,const cudaStream_t stream) 42b06137fdSPaul Mullowney { 43b06137fdSPaul Mullowney cusparseStatus_t stat; 44b06137fdSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 45b06137fdSPaul Mullowney 46b06137fdSPaul Mullowney PetscFunctionBegin; 47b06137fdSPaul Mullowney cusparsestruct->stream = stream; 48c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetStream(cusparsestruct->handle,cusparsestruct->stream);CHKERRCUDA(stat); 49b06137fdSPaul Mullowney PetscFunctionReturn(0); 50b06137fdSPaul Mullowney } 51b06137fdSPaul Mullowney 52b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetHandle(Mat A,const cusparseHandle_t handle) 53b06137fdSPaul Mullowney { 54b06137fdSPaul Mullowney cusparseStatus_t stat; 55b06137fdSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 56b06137fdSPaul Mullowney 57b06137fdSPaul Mullowney PetscFunctionBegin; 586b1cf21dSAlejandro Lamas Daviña if (cusparsestruct->handle != handle) { 5916a2e217SAlejandro Lamas Daviña if (cusparsestruct->handle) { 60c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUDA(stat); 6116a2e217SAlejandro Lamas Daviña } 62b06137fdSPaul Mullowney cusparsestruct->handle = handle; 636b1cf21dSAlejandro Lamas Daviña } 64c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 65b06137fdSPaul Mullowney PetscFunctionReturn(0); 66b06137fdSPaul Mullowney } 67b06137fdSPaul Mullowney 68b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSEClearHandle(Mat A) 69b06137fdSPaul Mullowney { 70b06137fdSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 71b06137fdSPaul Mullowney PetscFunctionBegin; 72b06137fdSPaul Mullowney if (cusparsestruct->handle) 73b06137fdSPaul Mullowney cusparsestruct->handle = 0; 74b06137fdSPaul Mullowney PetscFunctionReturn(0); 75b06137fdSPaul Mullowney } 76b06137fdSPaul Mullowney 77ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_cusparse(Mat A,MatSolverType *type) 789ae82921SPaul Mullowney { 799ae82921SPaul Mullowney PetscFunctionBegin; 809ae82921SPaul Mullowney *type = MATSOLVERCUSPARSE; 819ae82921SPaul Mullowney PetscFunctionReturn(0); 829ae82921SPaul Mullowney } 839ae82921SPaul Mullowney 84c708e6cdSJed Brown /*MC 85087f3262SPaul Mullowney MATSOLVERCUSPARSE = "cusparse" - A matrix type providing triangular solvers for seq matrices 86087f3262SPaul Mullowney on a single GPU of type, seqaijcusparse, aijcusparse, or seqaijcusp, aijcusp. Currently supported 87087f3262SPaul Mullowney algorithms are ILU(k) and ICC(k). Typically, deeper factorizations (larger k) results in poorer 88087f3262SPaul Mullowney performance in the triangular solves. Full LU, and Cholesky decompositions can be solved through the 89087f3262SPaul Mullowney CUSPARSE triangular solve algorithm. However, the performance can be quite poor and thus these 90087f3262SPaul Mullowney algorithms are not recommended. This class does NOT support direct solver operations. 91c708e6cdSJed Brown 929ae82921SPaul Mullowney Level: beginner 93c708e6cdSJed Brown 943ca39a21SBarry Smith .seealso: PCFactorSetMatSolverType(), MatSolverType, MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 95c708e6cdSJed Brown M*/ 969ae82921SPaul Mullowney 9742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat A,MatFactorType ftype,Mat *B) 989ae82921SPaul Mullowney { 999ae82921SPaul Mullowney PetscErrorCode ierr; 100bc3f50f2SPaul Mullowney PetscInt n = A->rmap->n; 1019ae82921SPaul Mullowney 1029ae82921SPaul Mullowney PetscFunctionBegin; 103bc3f50f2SPaul Mullowney ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 104404133a2SPaul Mullowney (*B)->factortype = ftype; 105bc3f50f2SPaul Mullowney ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 1069ae82921SPaul Mullowney ierr = MatSetType(*B,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 1072205254eSKarl Rupp 108087f3262SPaul Mullowney if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 10933d57670SJed Brown ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr); 1109ae82921SPaul Mullowney (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJCUSPARSE; 1119ae82921SPaul Mullowney (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJCUSPARSE; 112087f3262SPaul Mullowney } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 113087f3262SPaul Mullowney (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJCUSPARSE; 114087f3262SPaul Mullowney (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJCUSPARSE; 1159ae82921SPaul Mullowney } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported for CUSPARSE Matrix Types"); 116bc3f50f2SPaul Mullowney 117fa03d054SJed Brown ierr = MatSeqAIJSetPreallocation(*B,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 1183ca39a21SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)(*B),"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cusparse);CHKERRQ(ierr); 1199ae82921SPaul Mullowney PetscFunctionReturn(0); 1209ae82921SPaul Mullowney } 1219ae82921SPaul Mullowney 122bc3f50f2SPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetFormat_SeqAIJCUSPARSE(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format) 123ca45077fSPaul Mullowney { 124aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1256e111a19SKarl Rupp 126ca45077fSPaul Mullowney PetscFunctionBegin; 1272692e278SPaul Mullowney #if CUDA_VERSION>=4020 128ca45077fSPaul Mullowney switch (op) { 129e057df02SPaul Mullowney case MAT_CUSPARSE_MULT: 130aa372e3fSPaul Mullowney cusparsestruct->format = format; 131ca45077fSPaul Mullowney break; 132e057df02SPaul Mullowney case MAT_CUSPARSE_ALL: 133aa372e3fSPaul Mullowney cusparsestruct->format = format; 134ca45077fSPaul Mullowney break; 135ca45077fSPaul Mullowney default: 13636d62e41SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op); 137ca45077fSPaul Mullowney } 1382692e278SPaul Mullowney #else 1396c4ed002SBarry Smith 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."); 1402692e278SPaul Mullowney #endif 141ca45077fSPaul Mullowney PetscFunctionReturn(0); 142ca45077fSPaul Mullowney } 1439ae82921SPaul Mullowney 144e057df02SPaul Mullowney /*@ 145e057df02SPaul Mullowney MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular 146e057df02SPaul Mullowney operation. Only the MatMult operation can use different GPU storage formats 147aa372e3fSPaul Mullowney for MPIAIJCUSPARSE matrices. 148e057df02SPaul Mullowney Not Collective 149e057df02SPaul Mullowney 150e057df02SPaul Mullowney Input Parameters: 1518468deeeSKarl Rupp + A - Matrix of type SEQAIJCUSPARSE 15236d62e41SPaul Mullowney . 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. 1532692e278SPaul Mullowney - format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2) 154e057df02SPaul Mullowney 155e057df02SPaul Mullowney Output Parameter: 156e057df02SPaul Mullowney 157e057df02SPaul Mullowney Level: intermediate 158e057df02SPaul Mullowney 1598468deeeSKarl Rupp .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 160e057df02SPaul Mullowney @*/ 161e057df02SPaul Mullowney PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format) 162e057df02SPaul Mullowney { 163e057df02SPaul Mullowney PetscErrorCode ierr; 1646e111a19SKarl Rupp 165e057df02SPaul Mullowney PetscFunctionBegin; 166e057df02SPaul Mullowney PetscValidHeaderSpecific(A, MAT_CLASSID,1); 167e057df02SPaul Mullowney ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr); 168e057df02SPaul Mullowney PetscFunctionReturn(0); 169e057df02SPaul Mullowney } 170e057df02SPaul Mullowney 1714416b707SBarry Smith static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat A) 1729ae82921SPaul Mullowney { 1739ae82921SPaul Mullowney PetscErrorCode ierr; 174e057df02SPaul Mullowney MatCUSPARSEStorageFormat format; 1759ae82921SPaul Mullowney PetscBool flg; 176a183c035SDominic Meiser Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1776e111a19SKarl Rupp 1789ae82921SPaul Mullowney PetscFunctionBegin; 179e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSPARSE options");CHKERRQ(ierr); 1809ae82921SPaul Mullowney if (A->factortype==MAT_FACTOR_NONE) { 181e057df02SPaul Mullowney ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV", 182a183c035SDominic Meiser "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr); 183e057df02SPaul Mullowney if (flg) { 184e057df02SPaul Mullowney ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr); 185045c96e1SPaul Mullowney } 1869ae82921SPaul Mullowney } 1874c87dfd4SPaul Mullowney ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve", 188a183c035SDominic Meiser "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr); 1894c87dfd4SPaul Mullowney if (flg) { 1904c87dfd4SPaul Mullowney ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr); 1914c87dfd4SPaul Mullowney } 192*0af67c1bSStefano Zampini ierr = PetscOptionsTail();CHKERRQ(ierr); 1939ae82921SPaul Mullowney PetscFunctionReturn(0); 1949ae82921SPaul Mullowney 1959ae82921SPaul Mullowney } 1969ae82921SPaul Mullowney 1976fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 1989ae82921SPaul Mullowney { 1999ae82921SPaul Mullowney PetscErrorCode ierr; 2009ae82921SPaul Mullowney 2019ae82921SPaul Mullowney PetscFunctionBegin; 2029ae82921SPaul Mullowney ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr); 2039ae82921SPaul Mullowney B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE; 2049ae82921SPaul Mullowney PetscFunctionReturn(0); 2059ae82921SPaul Mullowney } 2069ae82921SPaul Mullowney 2076fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 2089ae82921SPaul Mullowney { 2099ae82921SPaul Mullowney PetscErrorCode ierr; 2109ae82921SPaul Mullowney 2119ae82921SPaul Mullowney PetscFunctionBegin; 2129ae82921SPaul Mullowney ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr); 2139ae82921SPaul Mullowney B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE; 2149ae82921SPaul Mullowney PetscFunctionReturn(0); 2159ae82921SPaul Mullowney } 2169ae82921SPaul Mullowney 217087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info) 218087f3262SPaul Mullowney { 219087f3262SPaul Mullowney PetscErrorCode ierr; 220087f3262SPaul Mullowney 221087f3262SPaul Mullowney PetscFunctionBegin; 222087f3262SPaul Mullowney ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr); 223087f3262SPaul Mullowney B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE; 224087f3262SPaul Mullowney PetscFunctionReturn(0); 225087f3262SPaul Mullowney } 226087f3262SPaul Mullowney 227087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info) 228087f3262SPaul Mullowney { 229087f3262SPaul Mullowney PetscErrorCode ierr; 230087f3262SPaul Mullowney 231087f3262SPaul Mullowney PetscFunctionBegin; 232087f3262SPaul Mullowney ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr); 233087f3262SPaul Mullowney B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE; 234087f3262SPaul Mullowney PetscFunctionReturn(0); 235087f3262SPaul Mullowney } 236087f3262SPaul Mullowney 237087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A) 2389ae82921SPaul Mullowney { 2399ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2409ae82921SPaul Mullowney PetscInt n = A->rmap->n; 2419ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 242aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 2439ae82921SPaul Mullowney cusparseStatus_t stat; 2449ae82921SPaul Mullowney const PetscInt *ai = a->i,*aj = a->j,*vi; 2459ae82921SPaul Mullowney const MatScalar *aa = a->a,*v; 2469ae82921SPaul Mullowney PetscInt *AiLo, *AjLo; 2479ae82921SPaul Mullowney PetscScalar *AALo; 2489ae82921SPaul Mullowney PetscInt i,nz, nzLower, offset, rowOffset; 249b175d8bbSPaul Mullowney PetscErrorCode ierr; 2509ae82921SPaul Mullowney 2519ae82921SPaul Mullowney PetscFunctionBegin; 252b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 2539ae82921SPaul Mullowney try { 2549ae82921SPaul Mullowney /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */ 2559ae82921SPaul Mullowney nzLower=n+ai[n]-ai[1]; 2569ae82921SPaul Mullowney 2579ae82921SPaul Mullowney /* Allocate Space for the lower triangular matrix */ 258c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 259c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUDA(ierr); 260c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUDA(ierr); 2619ae82921SPaul Mullowney 2629ae82921SPaul Mullowney /* Fill the lower triangular matrix */ 2639ae82921SPaul Mullowney AiLo[0] = (PetscInt) 0; 2649ae82921SPaul Mullowney AiLo[n] = nzLower; 2659ae82921SPaul Mullowney AjLo[0] = (PetscInt) 0; 2669ae82921SPaul Mullowney AALo[0] = (MatScalar) 1.0; 2679ae82921SPaul Mullowney v = aa; 2689ae82921SPaul Mullowney vi = aj; 2699ae82921SPaul Mullowney offset = 1; 2709ae82921SPaul Mullowney rowOffset= 1; 2719ae82921SPaul Mullowney for (i=1; i<n; i++) { 2729ae82921SPaul Mullowney nz = ai[i+1] - ai[i]; 273e057df02SPaul Mullowney /* additional 1 for the term on the diagonal */ 2749ae82921SPaul Mullowney AiLo[i] = rowOffset; 2759ae82921SPaul Mullowney rowOffset += nz+1; 2769ae82921SPaul Mullowney 2779ae82921SPaul Mullowney ierr = PetscMemcpy(&(AjLo[offset]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr); 2789ae82921SPaul Mullowney ierr = PetscMemcpy(&(AALo[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 2799ae82921SPaul Mullowney 2809ae82921SPaul Mullowney offset += nz; 2819ae82921SPaul Mullowney AjLo[offset] = (PetscInt) i; 2829ae82921SPaul Mullowney AALo[offset] = (MatScalar) 1.0; 2839ae82921SPaul Mullowney offset += 1; 2849ae82921SPaul Mullowney 2859ae82921SPaul Mullowney v += nz; 2869ae82921SPaul Mullowney vi += nz; 2879ae82921SPaul Mullowney } 2882205254eSKarl Rupp 289aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 290aa372e3fSPaul Mullowney loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 2912205254eSKarl Rupp 292aa372e3fSPaul Mullowney /* Create the matrix description */ 293c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat); 294c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 295c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 296c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUDA(stat); 297c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat); 298aa372e3fSPaul Mullowney 299aa372e3fSPaul Mullowney /* Create the solve analysis information */ 300c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat); 301aa372e3fSPaul Mullowney 302aa372e3fSPaul Mullowney /* set the operation */ 303aa372e3fSPaul Mullowney loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 304aa372e3fSPaul Mullowney 305aa372e3fSPaul Mullowney /* set the matrix */ 306aa372e3fSPaul Mullowney loTriFactor->csrMat = new CsrMatrix; 307aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows = n; 308aa372e3fSPaul Mullowney loTriFactor->csrMat->num_cols = n; 309aa372e3fSPaul Mullowney loTriFactor->csrMat->num_entries = nzLower; 310aa372e3fSPaul Mullowney 311aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1); 312aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1); 313aa372e3fSPaul Mullowney 314aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower); 315aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower); 316aa372e3fSPaul Mullowney 317aa372e3fSPaul Mullowney loTriFactor->csrMat->values = new THRUSTARRAY(nzLower); 318aa372e3fSPaul Mullowney loTriFactor->csrMat->values->assign(AALo, AALo+nzLower); 319aa372e3fSPaul Mullowney 320aa372e3fSPaul Mullowney /* perform the solve analysis */ 321aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp, 322aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr, 323aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(), 324c41cb2e2SAlejandro Lamas Daviña loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat); 325aa372e3fSPaul Mullowney 326aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 327aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor; 3282205254eSKarl Rupp 329c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AiLo);CHKERRCUDA(ierr); 330c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AjLo);CHKERRCUDA(ierr); 331c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr); 3329ae82921SPaul Mullowney } catch(char *ex) { 3339ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 3349ae82921SPaul Mullowney } 3359ae82921SPaul Mullowney } 3369ae82921SPaul Mullowney PetscFunctionReturn(0); 3379ae82921SPaul Mullowney } 3389ae82921SPaul Mullowney 339087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A) 3409ae82921SPaul Mullowney { 3419ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 3429ae82921SPaul Mullowney PetscInt n = A->rmap->n; 3439ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 344aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 3459ae82921SPaul Mullowney cusparseStatus_t stat; 3469ae82921SPaul Mullowney const PetscInt *aj = a->j,*adiag = a->diag,*vi; 3479ae82921SPaul Mullowney const MatScalar *aa = a->a,*v; 3489ae82921SPaul Mullowney PetscInt *AiUp, *AjUp; 3499ae82921SPaul Mullowney PetscScalar *AAUp; 3509ae82921SPaul Mullowney PetscInt i,nz, nzUpper, offset; 3519ae82921SPaul Mullowney PetscErrorCode ierr; 3529ae82921SPaul Mullowney 3539ae82921SPaul Mullowney PetscFunctionBegin; 354b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 3559ae82921SPaul Mullowney try { 3569ae82921SPaul Mullowney /* next, figure out the number of nonzeros in the upper triangular matrix. */ 3579ae82921SPaul Mullowney nzUpper = adiag[0]-adiag[n]; 3589ae82921SPaul Mullowney 3599ae82921SPaul Mullowney /* Allocate Space for the upper triangular matrix */ 360c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 361c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 362c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 3639ae82921SPaul Mullowney 3649ae82921SPaul Mullowney /* Fill the upper triangular matrix */ 3659ae82921SPaul Mullowney AiUp[0]=(PetscInt) 0; 3669ae82921SPaul Mullowney AiUp[n]=nzUpper; 3679ae82921SPaul Mullowney offset = nzUpper; 3689ae82921SPaul Mullowney for (i=n-1; i>=0; i--) { 3699ae82921SPaul Mullowney v = aa + adiag[i+1] + 1; 3709ae82921SPaul Mullowney vi = aj + adiag[i+1] + 1; 3719ae82921SPaul Mullowney 372e057df02SPaul Mullowney /* number of elements NOT on the diagonal */ 3739ae82921SPaul Mullowney nz = adiag[i] - adiag[i+1]-1; 3749ae82921SPaul Mullowney 375e057df02SPaul Mullowney /* decrement the offset */ 3769ae82921SPaul Mullowney offset -= (nz+1); 3779ae82921SPaul Mullowney 378e057df02SPaul Mullowney /* first, set the diagonal elements */ 3799ae82921SPaul Mullowney AjUp[offset] = (PetscInt) i; 38009f51544SAlejandro Lamas Daviña AAUp[offset] = (MatScalar)1./v[nz]; 3819ae82921SPaul Mullowney AiUp[i] = AiUp[i+1] - (nz+1); 3829ae82921SPaul Mullowney 3839ae82921SPaul Mullowney ierr = PetscMemcpy(&(AjUp[offset+1]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr); 3849ae82921SPaul Mullowney ierr = PetscMemcpy(&(AAUp[offset+1]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 3859ae82921SPaul Mullowney } 3862205254eSKarl Rupp 387aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 388aa372e3fSPaul Mullowney upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 3892205254eSKarl Rupp 390aa372e3fSPaul Mullowney /* Create the matrix description */ 391c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 392c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 393c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 394c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 395c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 396aa372e3fSPaul Mullowney 397aa372e3fSPaul Mullowney /* Create the solve analysis information */ 398c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 399aa372e3fSPaul Mullowney 400aa372e3fSPaul Mullowney /* set the operation */ 401aa372e3fSPaul Mullowney upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 402aa372e3fSPaul Mullowney 403aa372e3fSPaul Mullowney /* set the matrix */ 404aa372e3fSPaul Mullowney upTriFactor->csrMat = new CsrMatrix; 405aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows = n; 406aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols = n; 407aa372e3fSPaul Mullowney upTriFactor->csrMat->num_entries = nzUpper; 408aa372e3fSPaul Mullowney 409aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1); 410aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1); 411aa372e3fSPaul Mullowney 412aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper); 413aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper); 414aa372e3fSPaul Mullowney 415aa372e3fSPaul Mullowney upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper); 416aa372e3fSPaul Mullowney upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper); 417aa372e3fSPaul Mullowney 418aa372e3fSPaul Mullowney /* perform the solve analysis */ 419aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 420aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 421aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 422c41cb2e2SAlejandro Lamas Daviña upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 423aa372e3fSPaul Mullowney 424aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 425aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 4262205254eSKarl Rupp 427c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 428c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 429c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 4309ae82921SPaul Mullowney } catch(char *ex) { 4319ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 4329ae82921SPaul Mullowney } 4339ae82921SPaul Mullowney } 4349ae82921SPaul Mullowney PetscFunctionReturn(0); 4359ae82921SPaul Mullowney } 4369ae82921SPaul Mullowney 437087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A) 4389ae82921SPaul Mullowney { 4399ae82921SPaul Mullowney PetscErrorCode ierr; 4409ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 4419ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 4429ae82921SPaul Mullowney IS isrow = a->row,iscol = a->icol; 4439ae82921SPaul Mullowney PetscBool row_identity,col_identity; 4449ae82921SPaul Mullowney const PetscInt *r,*c; 4459ae82921SPaul Mullowney PetscInt n = A->rmap->n; 4469ae82921SPaul Mullowney 4479ae82921SPaul Mullowney PetscFunctionBegin; 448087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr); 449087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr); 4502205254eSKarl Rupp 451e65717acSKarl Rupp cusparseTriFactors->workVector = new THRUSTARRAY(n); 452aa372e3fSPaul Mullowney cusparseTriFactors->nnz=a->nz; 4539ae82921SPaul Mullowney 454b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 455e057df02SPaul Mullowney /*lower triangular indices */ 4569ae82921SPaul Mullowney ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 4579ae82921SPaul Mullowney ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 4582205254eSKarl Rupp if (!row_identity) { 459aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 460aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices->assign(r, r+n); 4612205254eSKarl Rupp } 4629ae82921SPaul Mullowney ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 4639ae82921SPaul Mullowney 464e057df02SPaul Mullowney /*upper triangular indices */ 4659ae82921SPaul Mullowney ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 4669ae82921SPaul Mullowney ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 4672205254eSKarl Rupp if (!col_identity) { 468aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 469aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices->assign(c, c+n); 4702205254eSKarl Rupp } 4719ae82921SPaul Mullowney ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 4729ae82921SPaul Mullowney PetscFunctionReturn(0); 4739ae82921SPaul Mullowney } 4749ae82921SPaul Mullowney 475087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A) 476087f3262SPaul Mullowney { 477087f3262SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 478087f3262SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 479aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 480aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 481087f3262SPaul Mullowney cusparseStatus_t stat; 482087f3262SPaul Mullowney PetscErrorCode ierr; 483087f3262SPaul Mullowney PetscInt *AiUp, *AjUp; 484087f3262SPaul Mullowney PetscScalar *AAUp; 485087f3262SPaul Mullowney PetscScalar *AALo; 486087f3262SPaul Mullowney PetscInt nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j; 487087f3262SPaul Mullowney Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)A->data; 488087f3262SPaul Mullowney const PetscInt *ai = b->i,*aj = b->j,*vj; 489087f3262SPaul Mullowney const MatScalar *aa = b->a,*v; 490087f3262SPaul Mullowney 491087f3262SPaul Mullowney PetscFunctionBegin; 492b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 493087f3262SPaul Mullowney try { 494087f3262SPaul Mullowney /* Allocate Space for the upper triangular matrix */ 495c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 496c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 497c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 498c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 499087f3262SPaul Mullowney 500087f3262SPaul Mullowney /* Fill the upper triangular matrix */ 501087f3262SPaul Mullowney AiUp[0]=(PetscInt) 0; 502087f3262SPaul Mullowney AiUp[n]=nzUpper; 503087f3262SPaul Mullowney offset = 0; 504087f3262SPaul Mullowney for (i=0; i<n; i++) { 505087f3262SPaul Mullowney /* set the pointers */ 506087f3262SPaul Mullowney v = aa + ai[i]; 507087f3262SPaul Mullowney vj = aj + ai[i]; 508087f3262SPaul Mullowney nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */ 509087f3262SPaul Mullowney 510087f3262SPaul Mullowney /* first, set the diagonal elements */ 511087f3262SPaul Mullowney AjUp[offset] = (PetscInt) i; 51209f51544SAlejandro Lamas Daviña AAUp[offset] = (MatScalar)1.0/v[nz]; 513087f3262SPaul Mullowney AiUp[i] = offset; 51409f51544SAlejandro Lamas Daviña AALo[offset] = (MatScalar)1.0/v[nz]; 515087f3262SPaul Mullowney 516087f3262SPaul Mullowney offset+=1; 517087f3262SPaul Mullowney if (nz>0) { 518087f3262SPaul Mullowney ierr = PetscMemcpy(&(AjUp[offset]), vj, nz*sizeof(PetscInt));CHKERRQ(ierr); 519087f3262SPaul Mullowney ierr = PetscMemcpy(&(AAUp[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr); 520087f3262SPaul Mullowney for (j=offset; j<offset+nz; j++) { 521087f3262SPaul Mullowney AAUp[j] = -AAUp[j]; 522087f3262SPaul Mullowney AALo[j] = AAUp[j]/v[nz]; 523087f3262SPaul Mullowney } 524087f3262SPaul Mullowney offset+=nz; 525087f3262SPaul Mullowney } 526087f3262SPaul Mullowney } 527087f3262SPaul Mullowney 528aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 529aa372e3fSPaul Mullowney upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 530087f3262SPaul Mullowney 531aa372e3fSPaul Mullowney /* Create the matrix description */ 532c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 533c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 534c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 535c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 536c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat); 537087f3262SPaul Mullowney 538aa372e3fSPaul Mullowney /* Create the solve analysis information */ 539c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 540aa372e3fSPaul Mullowney 541aa372e3fSPaul Mullowney /* set the operation */ 542aa372e3fSPaul Mullowney upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 543aa372e3fSPaul Mullowney 544aa372e3fSPaul Mullowney /* set the matrix */ 545aa372e3fSPaul Mullowney upTriFactor->csrMat = new CsrMatrix; 546aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows = A->rmap->n; 547aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols = A->cmap->n; 548aa372e3fSPaul Mullowney upTriFactor->csrMat->num_entries = a->nz; 549aa372e3fSPaul Mullowney 550aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 551aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 552aa372e3fSPaul Mullowney 553aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 554aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 555aa372e3fSPaul Mullowney 556aa372e3fSPaul Mullowney upTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 557aa372e3fSPaul Mullowney upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz); 558aa372e3fSPaul Mullowney 559aa372e3fSPaul Mullowney /* perform the solve analysis */ 560aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 561aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 562aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 563c41cb2e2SAlejandro Lamas Daviña upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 564aa372e3fSPaul Mullowney 565aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 566aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 567aa372e3fSPaul Mullowney 568aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 569aa372e3fSPaul Mullowney loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 570aa372e3fSPaul Mullowney 571aa372e3fSPaul Mullowney /* Create the matrix description */ 572c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat); 573c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 574c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 575c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 576c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 577aa372e3fSPaul Mullowney 578aa372e3fSPaul Mullowney /* Create the solve analysis information */ 579c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat); 580aa372e3fSPaul Mullowney 581aa372e3fSPaul Mullowney /* set the operation */ 582aa372e3fSPaul Mullowney loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE; 583aa372e3fSPaul Mullowney 584aa372e3fSPaul Mullowney /* set the matrix */ 585aa372e3fSPaul Mullowney loTriFactor->csrMat = new CsrMatrix; 586aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows = A->rmap->n; 587aa372e3fSPaul Mullowney loTriFactor->csrMat->num_cols = A->cmap->n; 588aa372e3fSPaul Mullowney loTriFactor->csrMat->num_entries = a->nz; 589aa372e3fSPaul Mullowney 590aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 591aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 592aa372e3fSPaul Mullowney 593aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 594aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 595aa372e3fSPaul Mullowney 596aa372e3fSPaul Mullowney loTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 597aa372e3fSPaul Mullowney loTriFactor->csrMat->values->assign(AALo, AALo+a->nz); 598aa372e3fSPaul Mullowney 599aa372e3fSPaul Mullowney /* perform the solve analysis */ 600aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp, 601aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr, 602aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(), 603c41cb2e2SAlejandro Lamas Daviña loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat); 604aa372e3fSPaul Mullowney 605aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 606aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor; 607087f3262SPaul Mullowney 608b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 609c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 610c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 611c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 612c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr); 613087f3262SPaul Mullowney } catch(char *ex) { 614087f3262SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 615087f3262SPaul Mullowney } 616087f3262SPaul Mullowney } 617087f3262SPaul Mullowney PetscFunctionReturn(0); 618087f3262SPaul Mullowney } 619087f3262SPaul Mullowney 620087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A) 6219ae82921SPaul Mullowney { 6229ae82921SPaul Mullowney PetscErrorCode ierr; 623087f3262SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 624087f3262SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 625087f3262SPaul Mullowney IS ip = a->row; 626087f3262SPaul Mullowney const PetscInt *rip; 627087f3262SPaul Mullowney PetscBool perm_identity; 628087f3262SPaul Mullowney PetscInt n = A->rmap->n; 629087f3262SPaul Mullowney 630087f3262SPaul Mullowney PetscFunctionBegin; 631087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr); 632e65717acSKarl Rupp cusparseTriFactors->workVector = new THRUSTARRAY(n); 633aa372e3fSPaul Mullowney cusparseTriFactors->nnz=(a->nz-n)*2 + n; 634aa372e3fSPaul Mullowney 635087f3262SPaul Mullowney /*lower triangular indices */ 636087f3262SPaul Mullowney ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 637087f3262SPaul Mullowney ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 638087f3262SPaul Mullowney if (!perm_identity) { 639aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 640aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices->assign(rip, rip+n); 641aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 642aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices->assign(rip, rip+n); 643087f3262SPaul Mullowney } 644087f3262SPaul Mullowney ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 645087f3262SPaul Mullowney PetscFunctionReturn(0); 646087f3262SPaul Mullowney } 647087f3262SPaul Mullowney 6486fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 6499ae82921SPaul Mullowney { 6509ae82921SPaul Mullowney Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 6519ae82921SPaul Mullowney IS isrow = b->row,iscol = b->col; 6529ae82921SPaul Mullowney PetscBool row_identity,col_identity; 653b175d8bbSPaul Mullowney PetscErrorCode ierr; 6549ae82921SPaul Mullowney 6559ae82921SPaul Mullowney PetscFunctionBegin; 6569ae82921SPaul Mullowney ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 657e057df02SPaul Mullowney /* determine which version of MatSolve needs to be used. */ 6589ae82921SPaul Mullowney ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 6599ae82921SPaul Mullowney ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 660bda325fcSPaul Mullowney if (row_identity && col_identity) { 661bda325fcSPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 662bda325fcSPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 663bda325fcSPaul Mullowney } else { 664bda325fcSPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE; 665bda325fcSPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 666bda325fcSPaul Mullowney } 6678dc1d2a3SPaul Mullowney 668e057df02SPaul Mullowney /* get the triangular factors */ 669087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 6709ae82921SPaul Mullowney PetscFunctionReturn(0); 6719ae82921SPaul Mullowney } 6729ae82921SPaul Mullowney 673087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 674087f3262SPaul Mullowney { 675087f3262SPaul Mullowney Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 676087f3262SPaul Mullowney IS ip = b->row; 677087f3262SPaul Mullowney PetscBool perm_identity; 678b175d8bbSPaul Mullowney PetscErrorCode ierr; 679087f3262SPaul Mullowney 680087f3262SPaul Mullowney PetscFunctionBegin; 681087f3262SPaul Mullowney ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 682087f3262SPaul Mullowney 683087f3262SPaul Mullowney /* determine which version of MatSolve needs to be used. */ 684087f3262SPaul Mullowney ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 685087f3262SPaul Mullowney if (perm_identity) { 686087f3262SPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 687087f3262SPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 688087f3262SPaul Mullowney } else { 689087f3262SPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE; 690087f3262SPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 691087f3262SPaul Mullowney } 692087f3262SPaul Mullowney 693087f3262SPaul Mullowney /* get the triangular factors */ 694087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 695087f3262SPaul Mullowney PetscFunctionReturn(0); 696087f3262SPaul Mullowney } 6979ae82921SPaul Mullowney 698b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A) 699bda325fcSPaul Mullowney { 700bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 701aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 702aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 703aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 704aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 705bda325fcSPaul Mullowney cusparseStatus_t stat; 706aa372e3fSPaul Mullowney cusparseIndexBase_t indexBase; 707aa372e3fSPaul Mullowney cusparseMatrixType_t matrixType; 708aa372e3fSPaul Mullowney cusparseFillMode_t fillMode; 709aa372e3fSPaul Mullowney cusparseDiagType_t diagType; 710b175d8bbSPaul Mullowney 711bda325fcSPaul Mullowney PetscFunctionBegin; 712bda325fcSPaul Mullowney 713aa372e3fSPaul Mullowney /*********************************************/ 714aa372e3fSPaul Mullowney /* Now the Transpose of the Lower Tri Factor */ 715aa372e3fSPaul Mullowney /*********************************************/ 716aa372e3fSPaul Mullowney 717aa372e3fSPaul Mullowney /* allocate space for the transpose of the lower triangular factor */ 718aa372e3fSPaul Mullowney loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 719aa372e3fSPaul Mullowney 720aa372e3fSPaul Mullowney /* set the matrix descriptors of the lower triangular factor */ 721aa372e3fSPaul Mullowney matrixType = cusparseGetMatType(loTriFactor->descr); 722aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(loTriFactor->descr); 723aa372e3fSPaul Mullowney fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 724aa372e3fSPaul Mullowney CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 725aa372e3fSPaul Mullowney diagType = cusparseGetMatDiagType(loTriFactor->descr); 726aa372e3fSPaul Mullowney 727aa372e3fSPaul Mullowney /* Create the matrix description */ 728c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat); 729c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat); 730c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat); 731c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat); 732c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat); 733aa372e3fSPaul Mullowney 734aa372e3fSPaul Mullowney /* Create the solve analysis information */ 735c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat); 736aa372e3fSPaul Mullowney 737aa372e3fSPaul Mullowney /* set the operation */ 738aa372e3fSPaul Mullowney loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 739aa372e3fSPaul Mullowney 740aa372e3fSPaul Mullowney /* allocate GPU space for the CSC of the lower triangular factor*/ 741aa372e3fSPaul Mullowney loTriFactorT->csrMat = new CsrMatrix; 742aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows; 743aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols; 744aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries; 745aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1); 746aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries); 747aa372e3fSPaul Mullowney loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries); 748aa372e3fSPaul Mullowney 749aa372e3fSPaul Mullowney /* compute the transpose of the lower triangular factor, i.e. the CSC */ 750aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows, 751aa372e3fSPaul Mullowney loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries, 752aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 753aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 754aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 755aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 756aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 757aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 758c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 759aa372e3fSPaul Mullowney 760aa372e3fSPaul Mullowney /* perform the solve analysis on the transposed matrix */ 761aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp, 762aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries, 763aa372e3fSPaul Mullowney loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(), 764aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(), 765c41cb2e2SAlejandro Lamas Daviña loTriFactorT->solveInfo);CHKERRCUDA(stat); 766aa372e3fSPaul Mullowney 767aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 768aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT; 769aa372e3fSPaul Mullowney 770aa372e3fSPaul Mullowney /*********************************************/ 771aa372e3fSPaul Mullowney /* Now the Transpose of the Upper Tri Factor */ 772aa372e3fSPaul Mullowney /*********************************************/ 773aa372e3fSPaul Mullowney 774aa372e3fSPaul Mullowney /* allocate space for the transpose of the upper triangular factor */ 775aa372e3fSPaul Mullowney upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 776aa372e3fSPaul Mullowney 777aa372e3fSPaul Mullowney /* set the matrix descriptors of the upper triangular factor */ 778aa372e3fSPaul Mullowney matrixType = cusparseGetMatType(upTriFactor->descr); 779aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(upTriFactor->descr); 780aa372e3fSPaul Mullowney fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 781aa372e3fSPaul Mullowney CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 782aa372e3fSPaul Mullowney diagType = cusparseGetMatDiagType(upTriFactor->descr); 783aa372e3fSPaul Mullowney 784aa372e3fSPaul Mullowney /* Create the matrix description */ 785c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat); 786c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat); 787c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat); 788c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat); 789c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat); 790aa372e3fSPaul Mullowney 791aa372e3fSPaul Mullowney /* Create the solve analysis information */ 792c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat); 793aa372e3fSPaul Mullowney 794aa372e3fSPaul Mullowney /* set the operation */ 795aa372e3fSPaul Mullowney upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 796aa372e3fSPaul Mullowney 797aa372e3fSPaul Mullowney /* allocate GPU space for the CSC of the upper triangular factor*/ 798aa372e3fSPaul Mullowney upTriFactorT->csrMat = new CsrMatrix; 799aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows; 800aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols; 801aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries; 802aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1); 803aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries); 804aa372e3fSPaul Mullowney upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries); 805aa372e3fSPaul Mullowney 806aa372e3fSPaul Mullowney /* compute the transpose of the upper triangular factor, i.e. the CSC */ 807aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows, 808aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries, 809aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 810aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 811aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 812aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 813aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 814aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 815c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 816aa372e3fSPaul Mullowney 817aa372e3fSPaul Mullowney /* perform the solve analysis on the transposed matrix */ 818aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp, 819aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries, 820aa372e3fSPaul Mullowney upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(), 821aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(), 822c41cb2e2SAlejandro Lamas Daviña upTriFactorT->solveInfo);CHKERRCUDA(stat); 823aa372e3fSPaul Mullowney 824aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 825aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT; 826bda325fcSPaul Mullowney PetscFunctionReturn(0); 827bda325fcSPaul Mullowney } 828bda325fcSPaul Mullowney 829b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A) 830bda325fcSPaul Mullowney { 831aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 832aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 833aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 834bda325fcSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 835bda325fcSPaul Mullowney cusparseStatus_t stat; 836aa372e3fSPaul Mullowney cusparseIndexBase_t indexBase; 837b06137fdSPaul Mullowney cudaError_t err; 838b175d8bbSPaul Mullowney 839bda325fcSPaul Mullowney PetscFunctionBegin; 840aa372e3fSPaul Mullowney 841aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 842aa372e3fSPaul Mullowney matstructT = new Mat_SeqAIJCUSPARSEMultStruct; 843c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat); 844aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(matstruct->descr); 845c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat); 846c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 847aa372e3fSPaul Mullowney 848b06137fdSPaul Mullowney /* set alpha and beta */ 849c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstructT->alpha),sizeof(PetscScalar));CHKERRCUDA(err); 850c41cb2e2SAlejandro Lamas Daviña err = cudaMemcpy(matstructT->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 851c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstructT->beta),sizeof(PetscScalar));CHKERRCUDA(err); 852c41cb2e2SAlejandro Lamas Daviña err = cudaMemcpy(matstructT->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 853c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 854b06137fdSPaul Mullowney 855aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 856aa372e3fSPaul Mullowney CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 857aa372e3fSPaul Mullowney CsrMatrix *matrixT= new CsrMatrix; 858554b8892SKarl Rupp matrixT->num_rows = A->cmap->n; 859554b8892SKarl Rupp matrixT->num_cols = A->rmap->n; 860aa372e3fSPaul Mullowney matrixT->num_entries = a->nz; 861aa372e3fSPaul Mullowney matrixT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 862aa372e3fSPaul Mullowney matrixT->column_indices = new THRUSTINTARRAY32(a->nz); 863aa372e3fSPaul Mullowney matrixT->values = new THRUSTARRAY(a->nz); 864aa372e3fSPaul Mullowney 865aa372e3fSPaul Mullowney /* compute the transpose of the upper triangular factor, i.e. the CSC */ 866aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(matstruct->descr); 867aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows, 868aa372e3fSPaul Mullowney matrix->num_cols, matrix->num_entries, 869aa372e3fSPaul Mullowney matrix->values->data().get(), 870aa372e3fSPaul Mullowney matrix->row_offsets->data().get(), 871aa372e3fSPaul Mullowney matrix->column_indices->data().get(), 872aa372e3fSPaul Mullowney matrixT->values->data().get(), 873aa372e3fSPaul Mullowney matrixT->column_indices->data().get(), 874aa372e3fSPaul Mullowney matrixT->row_offsets->data().get(), 875c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 876aa372e3fSPaul Mullowney 877aa372e3fSPaul Mullowney /* assign the pointer */ 878aa372e3fSPaul Mullowney matstructT->mat = matrixT; 879aa372e3fSPaul Mullowney 880aa372e3fSPaul Mullowney } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 8812692e278SPaul Mullowney #if CUDA_VERSION>=5000 882aa372e3fSPaul Mullowney /* First convert HYB to CSR */ 883aa372e3fSPaul Mullowney CsrMatrix *temp= new CsrMatrix; 884aa372e3fSPaul Mullowney temp->num_rows = A->rmap->n; 885aa372e3fSPaul Mullowney temp->num_cols = A->cmap->n; 886aa372e3fSPaul Mullowney temp->num_entries = a->nz; 887aa372e3fSPaul Mullowney temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 888aa372e3fSPaul Mullowney temp->column_indices = new THRUSTINTARRAY32(a->nz); 889aa372e3fSPaul Mullowney temp->values = new THRUSTARRAY(a->nz); 890aa372e3fSPaul Mullowney 8912692e278SPaul Mullowney 892aa372e3fSPaul Mullowney stat = cusparse_hyb2csr(cusparsestruct->handle, 893aa372e3fSPaul Mullowney matstruct->descr, (cusparseHybMat_t)matstruct->mat, 894aa372e3fSPaul Mullowney temp->values->data().get(), 895aa372e3fSPaul Mullowney temp->row_offsets->data().get(), 896c41cb2e2SAlejandro Lamas Daviña temp->column_indices->data().get());CHKERRCUDA(stat); 897aa372e3fSPaul Mullowney 898aa372e3fSPaul Mullowney /* Next, convert CSR to CSC (i.e. the matrix transpose) */ 899aa372e3fSPaul Mullowney CsrMatrix *tempT= new CsrMatrix; 900aa372e3fSPaul Mullowney tempT->num_rows = A->rmap->n; 901aa372e3fSPaul Mullowney tempT->num_cols = A->cmap->n; 902aa372e3fSPaul Mullowney tempT->num_entries = a->nz; 903aa372e3fSPaul Mullowney tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 904aa372e3fSPaul Mullowney tempT->column_indices = new THRUSTINTARRAY32(a->nz); 905aa372e3fSPaul Mullowney tempT->values = new THRUSTARRAY(a->nz); 906aa372e3fSPaul Mullowney 907aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows, 908aa372e3fSPaul Mullowney temp->num_cols, temp->num_entries, 909aa372e3fSPaul Mullowney temp->values->data().get(), 910aa372e3fSPaul Mullowney temp->row_offsets->data().get(), 911aa372e3fSPaul Mullowney temp->column_indices->data().get(), 912aa372e3fSPaul Mullowney tempT->values->data().get(), 913aa372e3fSPaul Mullowney tempT->column_indices->data().get(), 914aa372e3fSPaul Mullowney tempT->row_offsets->data().get(), 915c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 916aa372e3fSPaul Mullowney 917aa372e3fSPaul Mullowney /* Last, convert CSC to HYB */ 918aa372e3fSPaul Mullowney cusparseHybMat_t hybMat; 919c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 920aa372e3fSPaul Mullowney cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 921aa372e3fSPaul Mullowney CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 922aa372e3fSPaul Mullowney stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n, 923aa372e3fSPaul Mullowney matstructT->descr, tempT->values->data().get(), 924aa372e3fSPaul Mullowney tempT->row_offsets->data().get(), 925aa372e3fSPaul Mullowney tempT->column_indices->data().get(), 926c41cb2e2SAlejandro Lamas Daviña hybMat, 0, partition);CHKERRCUDA(stat); 927aa372e3fSPaul Mullowney 928aa372e3fSPaul Mullowney /* assign the pointer */ 929aa372e3fSPaul Mullowney matstructT->mat = hybMat; 930aa372e3fSPaul Mullowney 931aa372e3fSPaul Mullowney /* delete temporaries */ 932aa372e3fSPaul Mullowney if (tempT) { 933aa372e3fSPaul Mullowney if (tempT->values) delete (THRUSTARRAY*) tempT->values; 934aa372e3fSPaul Mullowney if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices; 935aa372e3fSPaul Mullowney if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets; 936aa372e3fSPaul Mullowney delete (CsrMatrix*) tempT; 937087f3262SPaul Mullowney } 938aa372e3fSPaul Mullowney if (temp) { 939aa372e3fSPaul Mullowney if (temp->values) delete (THRUSTARRAY*) temp->values; 940aa372e3fSPaul Mullowney if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices; 941aa372e3fSPaul Mullowney if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets; 942aa372e3fSPaul Mullowney delete (CsrMatrix*) temp; 943aa372e3fSPaul Mullowney } 9442692e278SPaul Mullowney #else 9452692e278SPaul Mullowney 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."); 9462692e278SPaul Mullowney #endif 947aa372e3fSPaul Mullowney } 948aa372e3fSPaul Mullowney /* assign the compressed row indices */ 949aa372e3fSPaul Mullowney matstructT->cprowIndices = new THRUSTINTARRAY; 950554b8892SKarl Rupp matstructT->cprowIndices->resize(A->cmap->n); 951554b8892SKarl Rupp thrust::sequence(matstructT->cprowIndices->begin(), matstructT->cprowIndices->end()); 952aa372e3fSPaul Mullowney 953aa372e3fSPaul Mullowney /* assign the pointer */ 954aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT; 955bda325fcSPaul Mullowney PetscFunctionReturn(0); 956bda325fcSPaul Mullowney } 957bda325fcSPaul Mullowney 9586fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 959bda325fcSPaul Mullowney { 960c41cb2e2SAlejandro Lamas Daviña PetscInt n = xx->map->n; 961465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 962465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 963465f34aeSAlejandro Lamas Daviña thrust::device_ptr<const PetscScalar> bGPU; 964465f34aeSAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> xGPU; 965bda325fcSPaul Mullowney cusparseStatus_t stat; 966bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 967aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 968aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 969aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 970b175d8bbSPaul Mullowney PetscErrorCode ierr; 971bda325fcSPaul Mullowney 972bda325fcSPaul Mullowney PetscFunctionBegin; 973aa372e3fSPaul Mullowney /* Analyze the matrix and create the transpose ... on the fly */ 974aa372e3fSPaul Mullowney if (!loTriFactorT && !upTriFactorT) { 975bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 976aa372e3fSPaul Mullowney loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 977aa372e3fSPaul Mullowney upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 978bda325fcSPaul Mullowney } 979bda325fcSPaul Mullowney 980bda325fcSPaul Mullowney /* Get the GPU pointers */ 981c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 982c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 983c41cb2e2SAlejandro Lamas Daviña xGPU = thrust::device_pointer_cast(xarray); 984c41cb2e2SAlejandro Lamas Daviña bGPU = thrust::device_pointer_cast(barray); 985bda325fcSPaul Mullowney 986aa372e3fSPaul Mullowney /* First, reorder with the row permutation */ 987c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 988c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()), 989c41cb2e2SAlejandro Lamas Daviña xGPU); 990aa372e3fSPaul Mullowney 991aa372e3fSPaul Mullowney /* First, solve U */ 992aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 993aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr, 994aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 995aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 996aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 997aa372e3fSPaul Mullowney upTriFactorT->solveInfo, 998c41cb2e2SAlejandro Lamas Daviña xarray, tempGPU->data().get());CHKERRCUDA(stat); 999aa372e3fSPaul Mullowney 1000aa372e3fSPaul Mullowney /* Then, solve L */ 1001aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 1002aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr, 1003aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 1004aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 1005aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 1006aa372e3fSPaul Mullowney loTriFactorT->solveInfo, 1007c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1008aa372e3fSPaul Mullowney 1009aa372e3fSPaul Mullowney /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1010c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1011c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()), 1012aa372e3fSPaul Mullowney tempGPU->begin()); 1013aa372e3fSPaul Mullowney 1014aa372e3fSPaul Mullowney /* Copy the temporary to the full solution. */ 1015c41cb2e2SAlejandro Lamas Daviña thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 1016bda325fcSPaul Mullowney 1017bda325fcSPaul Mullowney /* restore */ 1018c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1019c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1020c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1021087f3262SPaul Mullowney 1022aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1023bda325fcSPaul Mullowney PetscFunctionReturn(0); 1024bda325fcSPaul Mullowney } 1025bda325fcSPaul Mullowney 10266fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 1027bda325fcSPaul Mullowney { 1028465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1029465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 1030bda325fcSPaul Mullowney cusparseStatus_t stat; 1031bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1032aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1033aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1034aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1035b175d8bbSPaul Mullowney PetscErrorCode ierr; 1036bda325fcSPaul Mullowney 1037bda325fcSPaul Mullowney PetscFunctionBegin; 1038aa372e3fSPaul Mullowney /* Analyze the matrix and create the transpose ... on the fly */ 1039aa372e3fSPaul Mullowney if (!loTriFactorT && !upTriFactorT) { 1040bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 1041aa372e3fSPaul Mullowney loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1042aa372e3fSPaul Mullowney upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1043bda325fcSPaul Mullowney } 1044bda325fcSPaul Mullowney 1045bda325fcSPaul Mullowney /* Get the GPU pointers */ 1046c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1047c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1048bda325fcSPaul Mullowney 1049aa372e3fSPaul Mullowney /* First, solve U */ 1050aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 1051aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr, 1052aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 1053aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 1054aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 1055aa372e3fSPaul Mullowney upTriFactorT->solveInfo, 1056c41cb2e2SAlejandro Lamas Daviña barray, tempGPU->data().get());CHKERRCUDA(stat); 1057aa372e3fSPaul Mullowney 1058aa372e3fSPaul Mullowney /* Then, solve L */ 1059aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 1060aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr, 1061aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 1062aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 1063aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 1064aa372e3fSPaul Mullowney loTriFactorT->solveInfo, 1065c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1066bda325fcSPaul Mullowney 1067bda325fcSPaul Mullowney /* restore */ 1068c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1069c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1070c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1071aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1072bda325fcSPaul Mullowney PetscFunctionReturn(0); 1073bda325fcSPaul Mullowney } 1074bda325fcSPaul Mullowney 10756fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 10769ae82921SPaul Mullowney { 1077465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1078465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 1079465f34aeSAlejandro Lamas Daviña thrust::device_ptr<const PetscScalar> bGPU; 1080465f34aeSAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> xGPU; 10819ae82921SPaul Mullowney cusparseStatus_t stat; 10829ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1083aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1084aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1085aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1086b175d8bbSPaul Mullowney PetscErrorCode ierr; 10879ae82921SPaul Mullowney 10889ae82921SPaul Mullowney PetscFunctionBegin; 1089ebc8f436SDominic Meiser 1090e057df02SPaul Mullowney /* Get the GPU pointers */ 1091c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1092c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1093c41cb2e2SAlejandro Lamas Daviña xGPU = thrust::device_pointer_cast(xarray); 1094c41cb2e2SAlejandro Lamas Daviña bGPU = thrust::device_pointer_cast(barray); 10959ae82921SPaul Mullowney 1096aa372e3fSPaul Mullowney /* First, reorder with the row permutation */ 1097c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 1098c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()), 1099c41cb2e2SAlejandro Lamas Daviña xGPU); 1100aa372e3fSPaul Mullowney 1101aa372e3fSPaul Mullowney /* Next, solve L */ 1102aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 1103aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr, 1104aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 1105aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 1106aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 1107aa372e3fSPaul Mullowney loTriFactor->solveInfo, 1108c41cb2e2SAlejandro Lamas Daviña xarray, tempGPU->data().get());CHKERRCUDA(stat); 1109aa372e3fSPaul Mullowney 1110aa372e3fSPaul Mullowney /* Then, solve U */ 1111aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 1112aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr, 1113aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 1114aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 1115aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 1116aa372e3fSPaul Mullowney upTriFactor->solveInfo, 1117c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1118aa372e3fSPaul Mullowney 1119aa372e3fSPaul Mullowney /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1120c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1121c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->end()), 1122aa372e3fSPaul Mullowney tempGPU->begin()); 1123aa372e3fSPaul Mullowney 1124aa372e3fSPaul Mullowney /* Copy the temporary to the full solution. */ 1125c41cb2e2SAlejandro Lamas Daviña thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 11269ae82921SPaul Mullowney 1127c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1128c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1129c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1130aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 11319ae82921SPaul Mullowney PetscFunctionReturn(0); 11329ae82921SPaul Mullowney } 11339ae82921SPaul Mullowney 11346fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 11359ae82921SPaul Mullowney { 1136465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1137465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 11389ae82921SPaul Mullowney cusparseStatus_t stat; 11399ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1140aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1141aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1142aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1143b175d8bbSPaul Mullowney PetscErrorCode ierr; 11449ae82921SPaul Mullowney 11459ae82921SPaul Mullowney PetscFunctionBegin; 1146e057df02SPaul Mullowney /* Get the GPU pointers */ 1147c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1148c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 11499ae82921SPaul Mullowney 1150aa372e3fSPaul Mullowney /* First, solve L */ 1151aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 1152aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr, 1153aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 1154aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 1155aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 1156aa372e3fSPaul Mullowney loTriFactor->solveInfo, 1157c41cb2e2SAlejandro Lamas Daviña barray, tempGPU->data().get());CHKERRCUDA(stat); 1158aa372e3fSPaul Mullowney 1159aa372e3fSPaul Mullowney /* Next, solve U */ 1160aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 1161aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr, 1162aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 1163aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 1164aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 1165aa372e3fSPaul Mullowney upTriFactor->solveInfo, 1166c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 11679ae82921SPaul Mullowney 1168c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1169c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1170c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1171aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 11729ae82921SPaul Mullowney PetscFunctionReturn(0); 11739ae82921SPaul Mullowney } 11749ae82921SPaul Mullowney 11756fa9248bSJed Brown static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A) 11769ae82921SPaul Mullowney { 11779ae82921SPaul Mullowney 1178aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1179aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 11809ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 11819ae82921SPaul Mullowney PetscInt m = A->rmap->n,*ii,*ridx; 11829ae82921SPaul Mullowney PetscErrorCode ierr; 1183aa372e3fSPaul Mullowney cusparseStatus_t stat; 1184b06137fdSPaul Mullowney cudaError_t err; 11859ae82921SPaul Mullowney 11869ae82921SPaul Mullowney PetscFunctionBegin; 1187b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 11889ae82921SPaul Mullowney ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 118934d6c7a5SJose E. Roman if (A->assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) { 119034d6c7a5SJose E. Roman CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 119134d6c7a5SJose E. Roman /* copy values only */ 119234d6c7a5SJose E. Roman matrix->values->assign(a->a, a->a+a->nz); 119334d6c7a5SJose E. Roman } else { 1194470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&matstruct,cusparsestruct->format); 11959ae82921SPaul Mullowney try { 1196aa372e3fSPaul Mullowney cusparsestruct->nonzerorow=0; 1197aa372e3fSPaul Mullowney for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0); 11989ae82921SPaul Mullowney 11999ae82921SPaul Mullowney if (a->compressedrow.use) { 12009ae82921SPaul Mullowney m = a->compressedrow.nrows; 12019ae82921SPaul Mullowney ii = a->compressedrow.i; 12029ae82921SPaul Mullowney ridx = a->compressedrow.rindex; 12039ae82921SPaul Mullowney } else { 1204b06137fdSPaul Mullowney /* Forcing compressed row on the GPU */ 12059ae82921SPaul Mullowney int k=0; 1206854ce69bSBarry Smith ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr); 1207854ce69bSBarry Smith ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr); 12089ae82921SPaul Mullowney ii[0]=0; 12099ae82921SPaul Mullowney for (int j = 0; j<m; j++) { 12109ae82921SPaul Mullowney if ((a->i[j+1]-a->i[j])>0) { 12119ae82921SPaul Mullowney ii[k] = a->i[j]; 12129ae82921SPaul Mullowney ridx[k]= j; 12139ae82921SPaul Mullowney k++; 12149ae82921SPaul Mullowney } 12159ae82921SPaul Mullowney } 1216aa372e3fSPaul Mullowney ii[cusparsestruct->nonzerorow] = a->nz; 1217aa372e3fSPaul Mullowney m = cusparsestruct->nonzerorow; 12189ae82921SPaul Mullowney } 12199ae82921SPaul Mullowney 1220aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 1221aa372e3fSPaul Mullowney matstruct = new Mat_SeqAIJCUSPARSEMultStruct; 1222c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat); 1223c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 1224c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 12259ae82921SPaul Mullowney 1226c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstruct->alpha),sizeof(PetscScalar));CHKERRCUDA(err); 1227c41cb2e2SAlejandro Lamas Daviña err = cudaMemcpy(matstruct->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 1228c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstruct->beta),sizeof(PetscScalar));CHKERRCUDA(err); 1229c41cb2e2SAlejandro Lamas Daviña err = cudaMemcpy(matstruct->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 1230c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 1231b06137fdSPaul Mullowney 1232aa372e3fSPaul Mullowney /* Build a hybrid/ellpack matrix if this option is chosen for the storage */ 1233aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1234aa372e3fSPaul Mullowney /* set the matrix */ 1235aa372e3fSPaul Mullowney CsrMatrix *matrix= new CsrMatrix; 1236a65300a6SPaul Mullowney matrix->num_rows = m; 1237aa372e3fSPaul Mullowney matrix->num_cols = A->cmap->n; 1238aa372e3fSPaul Mullowney matrix->num_entries = a->nz; 1239a65300a6SPaul Mullowney matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1240a65300a6SPaul Mullowney matrix->row_offsets->assign(ii, ii + m+1); 12419ae82921SPaul Mullowney 1242aa372e3fSPaul Mullowney matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1243aa372e3fSPaul Mullowney matrix->column_indices->assign(a->j, a->j+a->nz); 1244aa372e3fSPaul Mullowney 1245aa372e3fSPaul Mullowney matrix->values = new THRUSTARRAY(a->nz); 1246aa372e3fSPaul Mullowney matrix->values->assign(a->a, a->a+a->nz); 1247aa372e3fSPaul Mullowney 1248aa372e3fSPaul Mullowney /* assign the pointer */ 1249aa372e3fSPaul Mullowney matstruct->mat = matrix; 1250aa372e3fSPaul Mullowney 1251aa372e3fSPaul Mullowney } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 12522692e278SPaul Mullowney #if CUDA_VERSION>=4020 1253aa372e3fSPaul Mullowney CsrMatrix *matrix= new CsrMatrix; 1254a65300a6SPaul Mullowney matrix->num_rows = m; 1255aa372e3fSPaul Mullowney matrix->num_cols = A->cmap->n; 1256aa372e3fSPaul Mullowney matrix->num_entries = a->nz; 1257a65300a6SPaul Mullowney matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1258a65300a6SPaul Mullowney matrix->row_offsets->assign(ii, ii + m+1); 1259aa372e3fSPaul Mullowney 1260aa372e3fSPaul Mullowney matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1261aa372e3fSPaul Mullowney matrix->column_indices->assign(a->j, a->j+a->nz); 1262aa372e3fSPaul Mullowney 1263aa372e3fSPaul Mullowney matrix->values = new THRUSTARRAY(a->nz); 1264aa372e3fSPaul Mullowney matrix->values->assign(a->a, a->a+a->nz); 1265aa372e3fSPaul Mullowney 1266aa372e3fSPaul Mullowney cusparseHybMat_t hybMat; 1267c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 1268aa372e3fSPaul Mullowney cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 1269aa372e3fSPaul Mullowney CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 1270a65300a6SPaul Mullowney stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols, 1271aa372e3fSPaul Mullowney matstruct->descr, matrix->values->data().get(), 1272aa372e3fSPaul Mullowney matrix->row_offsets->data().get(), 1273aa372e3fSPaul Mullowney matrix->column_indices->data().get(), 1274c41cb2e2SAlejandro Lamas Daviña hybMat, 0, partition);CHKERRCUDA(stat); 1275aa372e3fSPaul Mullowney /* assign the pointer */ 1276aa372e3fSPaul Mullowney matstruct->mat = hybMat; 1277aa372e3fSPaul Mullowney 1278aa372e3fSPaul Mullowney if (matrix) { 1279aa372e3fSPaul Mullowney if (matrix->values) delete (THRUSTARRAY*)matrix->values; 1280aa372e3fSPaul Mullowney if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices; 1281aa372e3fSPaul Mullowney if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets; 1282aa372e3fSPaul Mullowney delete (CsrMatrix*)matrix; 1283087f3262SPaul Mullowney } 12842692e278SPaul Mullowney #endif 1285087f3262SPaul Mullowney } 1286ca45077fSPaul Mullowney 1287aa372e3fSPaul Mullowney /* assign the compressed row indices */ 1288aa372e3fSPaul Mullowney matstruct->cprowIndices = new THRUSTINTARRAY(m); 1289aa372e3fSPaul Mullowney matstruct->cprowIndices->assign(ridx,ridx+m); 1290aa372e3fSPaul Mullowney 1291aa372e3fSPaul Mullowney /* assign the pointer */ 1292aa372e3fSPaul Mullowney cusparsestruct->mat = matstruct; 1293aa372e3fSPaul Mullowney 12949ae82921SPaul Mullowney if (!a->compressedrow.use) { 12959ae82921SPaul Mullowney ierr = PetscFree(ii);CHKERRQ(ierr); 12969ae82921SPaul Mullowney ierr = PetscFree(ridx);CHKERRQ(ierr); 12979ae82921SPaul Mullowney } 1298e65717acSKarl Rupp cusparsestruct->workVector = new THRUSTARRAY(m); 12999ae82921SPaul Mullowney } catch(char *ex) { 13009ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 13019ae82921SPaul Mullowney } 130234d6c7a5SJose E. Roman cusparsestruct->nonzerostate = A->nonzerostate; 130334d6c7a5SJose E. Roman } 13047d0e52d8SJose E. Roman ierr = WaitForGPU();CHKERRCUDA(ierr); 1305b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 13069ae82921SPaul Mullowney ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 13079ae82921SPaul Mullowney } 13089ae82921SPaul Mullowney PetscFunctionReturn(0); 13099ae82921SPaul Mullowney } 13109ae82921SPaul Mullowney 1311c41cb2e2SAlejandro Lamas Daviña struct VecCUDAPlusEquals 1312aa372e3fSPaul Mullowney { 1313aa372e3fSPaul Mullowney template <typename Tuple> 1314aa372e3fSPaul Mullowney __host__ __device__ 1315aa372e3fSPaul Mullowney void operator()(Tuple t) 1316aa372e3fSPaul Mullowney { 1317aa372e3fSPaul Mullowney thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t); 1318aa372e3fSPaul Mullowney } 1319aa372e3fSPaul Mullowney }; 1320aa372e3fSPaul Mullowney 13216fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 13229ae82921SPaul Mullowney { 13239ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1324aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 13259ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstruct; /* Do not initialize yet, because GPU data may not be available yet */ 1326465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1327465f34aeSAlejandro Lamas Daviña PetscScalar *yarray; 1328b175d8bbSPaul Mullowney PetscErrorCode ierr; 1329aa372e3fSPaul Mullowney cusparseStatus_t stat; 13309ae82921SPaul Mullowney 13319ae82921SPaul Mullowney PetscFunctionBegin; 133234d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 133334d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1334c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1335ca6ae6e6SAlejandro Lamas Daviña ierr = VecSet(yy,0);CHKERRQ(ierr); 1336c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr); 13379ff858a8SKarl Rupp matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 1338aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1339aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstruct->mat; 1340aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1341aa372e3fSPaul Mullowney mat->num_rows, mat->num_cols, mat->num_entries, 1342b06137fdSPaul Mullowney matstruct->alpha, matstruct->descr, mat->values->data().get(), mat->row_offsets->data().get(), 1343c41cb2e2SAlejandro Lamas Daviña mat->column_indices->data().get(), xarray, matstruct->beta, 1344c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 1345aa372e3fSPaul Mullowney } else { 13462692e278SPaul Mullowney #if CUDA_VERSION>=4020 1347aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat; 1348aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1349b06137fdSPaul Mullowney matstruct->alpha, matstruct->descr, hybMat, 1350c41cb2e2SAlejandro Lamas Daviña xarray, matstruct->beta, 1351c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 13522692e278SPaul Mullowney #endif 13539ae82921SPaul Mullowney } 1354c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1355c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr); 1356aa372e3fSPaul Mullowney if (!cusparsestruct->stream) { 1357c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1358ca45077fSPaul Mullowney } 1359aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr); 13609ae82921SPaul Mullowney PetscFunctionReturn(0); 13619ae82921SPaul Mullowney } 13629ae82921SPaul Mullowney 13636fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 1364ca45077fSPaul Mullowney { 1365ca45077fSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1366aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 13679ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstructT; 1368465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1369465f34aeSAlejandro Lamas Daviña PetscScalar *yarray; 1370b175d8bbSPaul Mullowney PetscErrorCode ierr; 1371aa372e3fSPaul Mullowney cusparseStatus_t stat; 1372ca45077fSPaul Mullowney 1373ca45077fSPaul Mullowney PetscFunctionBegin; 137434d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 137534d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 13769ff858a8SKarl Rupp matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1377aa372e3fSPaul Mullowney if (!matstructT) { 1378bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1379aa372e3fSPaul Mullowney matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1380bda325fcSPaul Mullowney } 1381c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1382ca6ae6e6SAlejandro Lamas Daviña ierr = VecSet(yy,0);CHKERRQ(ierr); 1383c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr); 1384aa372e3fSPaul Mullowney 1385aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1386aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1387aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1388aa372e3fSPaul Mullowney mat->num_rows, mat->num_cols, 1389b06137fdSPaul Mullowney mat->num_entries, matstructT->alpha, matstructT->descr, 1390aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 1391c41cb2e2SAlejandro Lamas Daviña mat->column_indices->data().get(), xarray, matstructT->beta, 1392c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 1393aa372e3fSPaul Mullowney } else { 13942692e278SPaul Mullowney #if CUDA_VERSION>=4020 1395aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1396aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1397b06137fdSPaul Mullowney matstructT->alpha, matstructT->descr, hybMat, 1398c41cb2e2SAlejandro Lamas Daviña xarray, matstructT->beta, 1399c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 14002692e278SPaul Mullowney #endif 1401ca45077fSPaul Mullowney } 1402c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1403c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr); 1404aa372e3fSPaul Mullowney if (!cusparsestruct->stream) { 1405c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1406ca45077fSPaul Mullowney } 1407aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr); 1408ca45077fSPaul Mullowney PetscFunctionReturn(0); 1409ca45077fSPaul Mullowney } 1410ca45077fSPaul Mullowney 1411aa372e3fSPaul Mullowney 14126fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 14139ae82921SPaul Mullowney { 14149ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1415aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 14169ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstruct; 1417c41cb2e2SAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> zptr; 1418465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1419465f34aeSAlejandro Lamas Daviña PetscScalar *zarray; 1420b175d8bbSPaul Mullowney PetscErrorCode ierr; 1421aa372e3fSPaul Mullowney cusparseStatus_t stat; 14226e111a19SKarl Rupp 14239ae82921SPaul Mullowney PetscFunctionBegin; 142434d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 142534d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 14269ff858a8SKarl Rupp matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 14279ae82921SPaul Mullowney try { 1428c41cb2e2SAlejandro Lamas Daviña ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 1429c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 143034d6c7a5SJose E. Roman ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1431c41cb2e2SAlejandro Lamas Daviña zptr = thrust::device_pointer_cast(zarray); 14329ae82921SPaul Mullowney 1433e057df02SPaul Mullowney /* multiply add */ 1434aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1435aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstruct->mat; 1436b06137fdSPaul Mullowney /* here we need to be careful to set the number of rows in the multiply to the 1437b06137fdSPaul Mullowney number of compressed rows in the matrix ... which is equivalent to the 1438b06137fdSPaul Mullowney size of the workVector */ 1439aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1440a65300a6SPaul Mullowney mat->num_rows, mat->num_cols, 1441b06137fdSPaul Mullowney mat->num_entries, matstruct->alpha, matstruct->descr, 1442aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 1443c41cb2e2SAlejandro Lamas Daviña mat->column_indices->data().get(), xarray, matstruct->beta, 1444c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1445aa372e3fSPaul Mullowney } else { 14462692e278SPaul Mullowney #if CUDA_VERSION>=4020 1447aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat; 1448a65300a6SPaul Mullowney if (cusparsestruct->workVector->size()) { 1449aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1450b06137fdSPaul Mullowney matstruct->alpha, matstruct->descr, hybMat, 1451c41cb2e2SAlejandro Lamas Daviña xarray, matstruct->beta, 1452c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1453a65300a6SPaul Mullowney } 14542692e278SPaul Mullowney #endif 1455aa372e3fSPaul Mullowney } 1456aa372e3fSPaul Mullowney 1457aa372e3fSPaul Mullowney /* scatter the data from the temporary into the full vector with a += operation */ 1458c41cb2e2SAlejandro Lamas Daviña thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))), 1459c41cb2e2SAlejandro Lamas Daviña thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(), 1460c41cb2e2SAlejandro Lamas Daviña VecCUDAPlusEquals()); 1461c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 146234d6c7a5SJose E. Roman ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 14639ae82921SPaul Mullowney 14649ae82921SPaul Mullowney } catch(char *ex) { 14659ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 14669ae82921SPaul Mullowney } 1467c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 14689ae82921SPaul Mullowney ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 14699ae82921SPaul Mullowney PetscFunctionReturn(0); 14709ae82921SPaul Mullowney } 14719ae82921SPaul Mullowney 14726fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 1473ca45077fSPaul Mullowney { 1474ca45077fSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1475aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 14769ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstructT; 1477c41cb2e2SAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> zptr; 1478465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1479465f34aeSAlejandro Lamas Daviña PetscScalar *zarray; 1480b175d8bbSPaul Mullowney PetscErrorCode ierr; 1481aa372e3fSPaul Mullowney cusparseStatus_t stat; 14826e111a19SKarl Rupp 1483ca45077fSPaul Mullowney PetscFunctionBegin; 148434d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 148534d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 14869ff858a8SKarl Rupp matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1487aa372e3fSPaul Mullowney if (!matstructT) { 1488bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1489aa372e3fSPaul Mullowney matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1490bda325fcSPaul Mullowney } 1491aa372e3fSPaul Mullowney 1492ca45077fSPaul Mullowney try { 1493c41cb2e2SAlejandro Lamas Daviña ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 1494c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 149534d6c7a5SJose E. Roman ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1496c41cb2e2SAlejandro Lamas Daviña zptr = thrust::device_pointer_cast(zarray); 1497ca45077fSPaul Mullowney 1498e057df02SPaul Mullowney /* multiply add with matrix transpose */ 1499aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1500aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1501b06137fdSPaul Mullowney /* here we need to be careful to set the number of rows in the multiply to the 1502b06137fdSPaul Mullowney number of compressed rows in the matrix ... which is equivalent to the 1503b06137fdSPaul Mullowney size of the workVector */ 1504aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1505a65300a6SPaul Mullowney mat->num_rows, mat->num_cols, 1506b06137fdSPaul Mullowney mat->num_entries, matstructT->alpha, matstructT->descr, 1507aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 1508c41cb2e2SAlejandro Lamas Daviña mat->column_indices->data().get(), xarray, matstructT->beta, 1509c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1510aa372e3fSPaul Mullowney } else { 15112692e278SPaul Mullowney #if CUDA_VERSION>=4020 1512aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1513a65300a6SPaul Mullowney if (cusparsestruct->workVector->size()) { 1514aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1515b06137fdSPaul Mullowney matstructT->alpha, matstructT->descr, hybMat, 1516c41cb2e2SAlejandro Lamas Daviña xarray, matstructT->beta, 1517c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1518a65300a6SPaul Mullowney } 15192692e278SPaul Mullowney #endif 1520aa372e3fSPaul Mullowney } 1521aa372e3fSPaul Mullowney 1522aa372e3fSPaul Mullowney /* scatter the data from the temporary into the full vector with a += operation */ 1523c41cb2e2SAlejandro Lamas Daviña thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))), 1524554b8892SKarl Rupp thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))) + A->cmap->n, 1525c41cb2e2SAlejandro Lamas Daviña VecCUDAPlusEquals()); 1526ca45077fSPaul Mullowney 1527c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 152834d6c7a5SJose E. Roman ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr); 1529ca45077fSPaul Mullowney 1530ca45077fSPaul Mullowney } catch(char *ex) { 1531ca45077fSPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 1532ca45077fSPaul Mullowney } 1533c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1534ca45077fSPaul Mullowney ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1535ca45077fSPaul Mullowney PetscFunctionReturn(0); 1536ca45077fSPaul Mullowney } 1537ca45077fSPaul Mullowney 15386fa9248bSJed Brown static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode) 15399ae82921SPaul Mullowney { 15409ae82921SPaul Mullowney PetscErrorCode ierr; 15416e111a19SKarl Rupp 15429ae82921SPaul Mullowney PetscFunctionBegin; 15439ae82921SPaul Mullowney ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr); 1544bc3f50f2SPaul Mullowney if (A->factortype==MAT_FACTOR_NONE) { 1545e057df02SPaul Mullowney ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1546bc3f50f2SPaul Mullowney } 15479ae82921SPaul Mullowney if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1548bbf3fe20SPaul Mullowney A->ops->mult = MatMult_SeqAIJCUSPARSE; 1549bbf3fe20SPaul Mullowney A->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1550bbf3fe20SPaul Mullowney A->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1551bbf3fe20SPaul Mullowney A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 15529ae82921SPaul Mullowney PetscFunctionReturn(0); 15539ae82921SPaul Mullowney } 15549ae82921SPaul Mullowney 15559ae82921SPaul Mullowney /* --------------------------------------------------------------------------------*/ 1556e057df02SPaul Mullowney /*@ 15579ae82921SPaul Mullowney MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format 1558e057df02SPaul Mullowney (the default parallel PETSc format). This matrix will ultimately pushed down 1559e057df02SPaul Mullowney to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix 1560e057df02SPaul Mullowney assembly performance the user should preallocate the matrix storage by setting 1561e057df02SPaul Mullowney the parameter nz (or the array nnz). By setting these parameters accurately, 1562e057df02SPaul Mullowney performance during matrix assembly can be increased by more than a factor of 50. 15639ae82921SPaul Mullowney 15649ae82921SPaul Mullowney Collective on MPI_Comm 15659ae82921SPaul Mullowney 15669ae82921SPaul Mullowney Input Parameters: 15679ae82921SPaul Mullowney + comm - MPI communicator, set to PETSC_COMM_SELF 15689ae82921SPaul Mullowney . m - number of rows 15699ae82921SPaul Mullowney . n - number of columns 15709ae82921SPaul Mullowney . nz - number of nonzeros per row (same for all rows) 15719ae82921SPaul Mullowney - nnz - array containing the number of nonzeros in the various rows 15720298fd71SBarry Smith (possibly different for each row) or NULL 15739ae82921SPaul Mullowney 15749ae82921SPaul Mullowney Output Parameter: 15759ae82921SPaul Mullowney . A - the matrix 15769ae82921SPaul Mullowney 15779ae82921SPaul Mullowney It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 15789ae82921SPaul Mullowney MatXXXXSetPreallocation() paradgm instead of this routine directly. 15799ae82921SPaul Mullowney [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 15809ae82921SPaul Mullowney 15819ae82921SPaul Mullowney Notes: 15829ae82921SPaul Mullowney If nnz is given then nz is ignored 15839ae82921SPaul Mullowney 15849ae82921SPaul Mullowney The AIJ format (also called the Yale sparse matrix format or 15859ae82921SPaul Mullowney compressed row storage), is fully compatible with standard Fortran 77 15869ae82921SPaul Mullowney storage. That is, the stored row and column indices can begin at 15879ae82921SPaul Mullowney either one (as in Fortran) or zero. See the users' manual for details. 15889ae82921SPaul Mullowney 15899ae82921SPaul Mullowney Specify the preallocated storage with either nz or nnz (not both). 15900298fd71SBarry Smith Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 15919ae82921SPaul Mullowney allocation. For large problems you MUST preallocate memory or you 15929ae82921SPaul Mullowney will get TERRIBLE performance, see the users' manual chapter on matrices. 15939ae82921SPaul Mullowney 15949ae82921SPaul Mullowney By default, this format uses inodes (identical nodes) when possible, to 15959ae82921SPaul Mullowney improve numerical efficiency of matrix-vector products and solves. We 15969ae82921SPaul Mullowney search for consecutive rows with the same nonzero structure, thereby 15979ae82921SPaul Mullowney reusing matrix information to achieve increased efficiency. 15989ae82921SPaul Mullowney 15999ae82921SPaul Mullowney Level: intermediate 16009ae82921SPaul Mullowney 1601e057df02SPaul Mullowney .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE 16029ae82921SPaul Mullowney @*/ 16039ae82921SPaul Mullowney PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 16049ae82921SPaul Mullowney { 16059ae82921SPaul Mullowney PetscErrorCode ierr; 16069ae82921SPaul Mullowney 16079ae82921SPaul Mullowney PetscFunctionBegin; 16089ae82921SPaul Mullowney ierr = MatCreate(comm,A);CHKERRQ(ierr); 16099ae82921SPaul Mullowney ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 16109ae82921SPaul Mullowney ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 16119ae82921SPaul Mullowney ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr); 16129ae82921SPaul Mullowney PetscFunctionReturn(0); 16139ae82921SPaul Mullowney } 16149ae82921SPaul Mullowney 16156fa9248bSJed Brown static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A) 16169ae82921SPaul Mullowney { 16179ae82921SPaul Mullowney PetscErrorCode ierr; 1618ab25e6cbSDominic Meiser 16199ae82921SPaul Mullowney PetscFunctionBegin; 16209ae82921SPaul Mullowney if (A->factortype==MAT_FACTOR_NONE) { 1621b8ced49eSKarl Rupp if (A->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 1622470880abSPatrick Sanan ierr = MatSeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr); 16239ae82921SPaul Mullowney } 16249ae82921SPaul Mullowney } else { 1625470880abSPatrick Sanan ierr = MatSeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr); 1626aa372e3fSPaul Mullowney } 16279ae82921SPaul Mullowney ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 16289ae82921SPaul Mullowney PetscFunctionReturn(0); 16299ae82921SPaul Mullowney } 16309ae82921SPaul Mullowney 16319ff858a8SKarl Rupp static PetscErrorCode MatDuplicate_SeqAIJCUSPARSE(Mat A,MatDuplicateOption cpvalues,Mat *B) 16329ff858a8SKarl Rupp { 16339ff858a8SKarl Rupp PetscErrorCode ierr; 16349ff858a8SKarl Rupp Mat C; 16359ff858a8SKarl Rupp cusparseStatus_t stat; 16369ff858a8SKarl Rupp cusparseHandle_t handle=0; 16379ff858a8SKarl Rupp 16389ff858a8SKarl Rupp PetscFunctionBegin; 16399ff858a8SKarl Rupp ierr = MatDuplicate_SeqAIJ(A,cpvalues,B);CHKERRQ(ierr); 16409ff858a8SKarl Rupp C = *B; 164134136279SStefano Zampini ierr = PetscFree(C->defaultvectype);CHKERRQ(ierr); 164234136279SStefano Zampini ierr = PetscStrallocpy(VECCUDA,&C->defaultvectype);CHKERRQ(ierr); 164334136279SStefano Zampini 16449ff858a8SKarl Rupp /* inject CUSPARSE-specific stuff */ 16459ff858a8SKarl Rupp if (C->factortype==MAT_FACTOR_NONE) { 16469ff858a8SKarl Rupp /* you cannot check the inode.use flag here since the matrix was just created. 16479ff858a8SKarl Rupp now build a GPU matrix data structure */ 16489ff858a8SKarl Rupp C->spptr = new Mat_SeqAIJCUSPARSE; 16499ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->mat = 0; 16509ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->matTranspose = 0; 16519ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->workVector = 0; 16529ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->format = MAT_CUSPARSE_CSR; 16539ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream = 0; 16549ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle = 0; 16559ff858a8SKarl Rupp stat = cusparseCreate(&handle);CHKERRCUDA(stat); 16569ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle = handle; 16579ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream = 0; 16589ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->nonzerostate = 0; 16599ff858a8SKarl Rupp } else { 16609ff858a8SKarl Rupp /* NEXT, set the pointers to the triangular factors */ 16619ff858a8SKarl Rupp C->spptr = new Mat_SeqAIJCUSPARSETriFactors; 16629ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtr = 0; 16639ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtr = 0; 16649ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtrTranspose = 0; 16659ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtrTranspose = 0; 16669ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->rpermIndices = 0; 16679ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->cpermIndices = 0; 16689ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->workVector = 0; 16699ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle = 0; 16709ff858a8SKarl Rupp stat = cusparseCreate(&handle);CHKERRCUDA(stat); 16719ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle = handle; 16729ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->nnz = 0; 16739ff858a8SKarl Rupp } 16749ff858a8SKarl Rupp 16759ff858a8SKarl Rupp C->ops->assemblyend = MatAssemblyEnd_SeqAIJCUSPARSE; 16769ff858a8SKarl Rupp C->ops->destroy = MatDestroy_SeqAIJCUSPARSE; 16779ff858a8SKarl Rupp C->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSPARSE; 16789ff858a8SKarl Rupp C->ops->mult = MatMult_SeqAIJCUSPARSE; 16799ff858a8SKarl Rupp C->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 16809ff858a8SKarl Rupp C->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 16819ff858a8SKarl Rupp C->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 16829ff858a8SKarl Rupp C->ops->duplicate = MatDuplicate_SeqAIJCUSPARSE; 16839ff858a8SKarl Rupp 16849ff858a8SKarl Rupp ierr = PetscObjectChangeTypeName((PetscObject)C,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 16859ff858a8SKarl Rupp 1686b8ced49eSKarl Rupp C->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED; 16879ff858a8SKarl Rupp 16889ff858a8SKarl Rupp ierr = PetscObjectComposeFunction((PetscObject)C, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr); 16899ff858a8SKarl Rupp PetscFunctionReturn(0); 16909ff858a8SKarl Rupp } 16919ff858a8SKarl Rupp 16928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B) 16939ae82921SPaul Mullowney { 16949ae82921SPaul Mullowney PetscErrorCode ierr; 1695aa372e3fSPaul Mullowney cusparseStatus_t stat; 1696aa372e3fSPaul Mullowney cusparseHandle_t handle=0; 16979ae82921SPaul Mullowney 16989ae82921SPaul Mullowney PetscFunctionBegin; 16999ae82921SPaul Mullowney ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr); 170034136279SStefano Zampini ierr = PetscFree(B->defaultvectype);CHKERRQ(ierr); 170134136279SStefano Zampini ierr = PetscStrallocpy(VECCUDA,&B->defaultvectype);CHKERRQ(ierr); 170234136279SStefano Zampini 17039ae82921SPaul Mullowney if (B->factortype==MAT_FACTOR_NONE) { 1704e057df02SPaul Mullowney /* you cannot check the inode.use flag here since the matrix was just created. 1705e057df02SPaul Mullowney now build a GPU matrix data structure */ 17069ae82921SPaul Mullowney B->spptr = new Mat_SeqAIJCUSPARSE; 17079ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat = 0; 1708aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0; 1709aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector = 0; 1710e057df02SPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->format = MAT_CUSPARSE_CSR; 1711aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 1712aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = 0; 1713c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1714aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = handle; 1715aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 17169ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)B->spptr)->nonzerostate = 0; 17179ae82921SPaul Mullowney } else { 17189ae82921SPaul Mullowney /* NEXT, set the pointers to the triangular factors */ 1719debe9ee2SPaul Mullowney B->spptr = new Mat_SeqAIJCUSPARSETriFactors; 17209ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr = 0; 17219ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr = 0; 1722aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0; 1723aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0; 1724aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices = 0; 1725aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices = 0; 1726aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector = 0; 1727aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = 0; 1728c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1729aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = handle; 1730aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz = 0; 17319ae82921SPaul Mullowney } 1732aa372e3fSPaul Mullowney 17339ae82921SPaul Mullowney B->ops->assemblyend = MatAssemblyEnd_SeqAIJCUSPARSE; 17349ae82921SPaul Mullowney B->ops->destroy = MatDestroy_SeqAIJCUSPARSE; 17359ae82921SPaul Mullowney B->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSPARSE; 1736ca45077fSPaul Mullowney B->ops->mult = MatMult_SeqAIJCUSPARSE; 1737ca45077fSPaul Mullowney B->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1738ca45077fSPaul Mullowney B->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1739ca45077fSPaul Mullowney B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 17409ff858a8SKarl Rupp B->ops->duplicate = MatDuplicate_SeqAIJCUSPARSE; 17412205254eSKarl Rupp 17429ae82921SPaul Mullowney ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 17432205254eSKarl Rupp 1744b8ced49eSKarl Rupp B->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED; 17452205254eSKarl Rupp 1746bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr); 17479ae82921SPaul Mullowney PetscFunctionReturn(0); 17489ae82921SPaul Mullowney } 17499ae82921SPaul Mullowney 17503ca39a21SBarry Smith /*MC 1751e057df02SPaul Mullowney MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices. 1752e057df02SPaul Mullowney 1753e057df02SPaul Mullowney A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either 17542692e278SPaul Mullowney CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later. 17552692e278SPaul Mullowney All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library. 1756e057df02SPaul Mullowney 1757e057df02SPaul Mullowney Options Database Keys: 1758e057df02SPaul Mullowney + -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions() 1759aa372e3fSPaul Mullowney . -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). 1760aa372e3fSPaul Mullowney . -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). 1761e057df02SPaul Mullowney 1762e057df02SPaul Mullowney Level: beginner 1763e057df02SPaul Mullowney 17648468deeeSKarl Rupp .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 1765e057df02SPaul Mullowney M*/ 17667f756511SDominic Meiser 176742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*); 176842c9c57cSBarry Smith 17690f39cd5aSBarry Smith 17703ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_CUSPARSE(void) 177142c9c57cSBarry Smith { 177242c9c57cSBarry Smith PetscErrorCode ierr; 177342c9c57cSBarry Smith 177442c9c57cSBarry Smith PetscFunctionBegin; 17753ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17763ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17773ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17783ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 177942c9c57cSBarry Smith PetscFunctionReturn(0); 178042c9c57cSBarry Smith } 178129b38603SBarry Smith 178281e08676SBarry Smith 1783470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct) 17847f756511SDominic Meiser { 17857f756511SDominic Meiser cusparseStatus_t stat; 17867f756511SDominic Meiser cusparseHandle_t handle; 17877f756511SDominic Meiser 17887f756511SDominic Meiser PetscFunctionBegin; 17897f756511SDominic Meiser if (*cusparsestruct) { 1790470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format); 1791470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format); 17927f756511SDominic Meiser delete (*cusparsestruct)->workVector; 17937f756511SDominic Meiser if (handle = (*cusparsestruct)->handle) { 1794c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroy(handle);CHKERRCUDA(stat); 17957f756511SDominic Meiser } 17967f756511SDominic Meiser delete *cusparsestruct; 17977f756511SDominic Meiser *cusparsestruct = 0; 17987f756511SDominic Meiser } 17997f756511SDominic Meiser PetscFunctionReturn(0); 18007f756511SDominic Meiser } 18017f756511SDominic Meiser 18027f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat) 18037f756511SDominic Meiser { 18047f756511SDominic Meiser PetscFunctionBegin; 18057f756511SDominic Meiser if (*mat) { 18067f756511SDominic Meiser delete (*mat)->values; 18077f756511SDominic Meiser delete (*mat)->column_indices; 18087f756511SDominic Meiser delete (*mat)->row_offsets; 18097f756511SDominic Meiser delete *mat; 18107f756511SDominic Meiser *mat = 0; 18117f756511SDominic Meiser } 18127f756511SDominic Meiser PetscFunctionReturn(0); 18137f756511SDominic Meiser } 18147f756511SDominic Meiser 1815470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor) 18167f756511SDominic Meiser { 18177f756511SDominic Meiser cusparseStatus_t stat; 18187f756511SDominic Meiser PetscErrorCode ierr; 18197f756511SDominic Meiser 18207f756511SDominic Meiser PetscFunctionBegin; 18217f756511SDominic Meiser if (*trifactor) { 1822c41cb2e2SAlejandro Lamas Daviña if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); } 1823c41cb2e2SAlejandro Lamas Daviña if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); } 18247f756511SDominic Meiser ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr); 18257f756511SDominic Meiser delete *trifactor; 18267f756511SDominic Meiser *trifactor = 0; 18277f756511SDominic Meiser } 18287f756511SDominic Meiser PetscFunctionReturn(0); 18297f756511SDominic Meiser } 18307f756511SDominic Meiser 1831470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format) 18327f756511SDominic Meiser { 18337f756511SDominic Meiser CsrMatrix *mat; 18347f756511SDominic Meiser cusparseStatus_t stat; 18357f756511SDominic Meiser cudaError_t err; 18367f756511SDominic Meiser 18377f756511SDominic Meiser PetscFunctionBegin; 18387f756511SDominic Meiser if (*matstruct) { 18397f756511SDominic Meiser if ((*matstruct)->mat) { 18407f756511SDominic Meiser if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) { 18417f756511SDominic Meiser cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat; 1842c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat); 18437f756511SDominic Meiser } else { 18447f756511SDominic Meiser mat = (CsrMatrix*)(*matstruct)->mat; 18457f756511SDominic Meiser CsrMatrix_Destroy(&mat); 18467f756511SDominic Meiser } 18477f756511SDominic Meiser } 1848c41cb2e2SAlejandro Lamas Daviña if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); } 18497f756511SDominic Meiser delete (*matstruct)->cprowIndices; 1850c41cb2e2SAlejandro Lamas Daviña if ((*matstruct)->alpha) { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); } 1851c41cb2e2SAlejandro Lamas Daviña if ((*matstruct)->beta) { err=cudaFree((*matstruct)->beta);CHKERRCUDA(err); } 18527f756511SDominic Meiser delete *matstruct; 18537f756511SDominic Meiser *matstruct = 0; 18547f756511SDominic Meiser } 18557f756511SDominic Meiser PetscFunctionReturn(0); 18567f756511SDominic Meiser } 18577f756511SDominic Meiser 1858470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors) 18597f756511SDominic Meiser { 18607f756511SDominic Meiser cusparseHandle_t handle; 18617f756511SDominic Meiser cusparseStatus_t stat; 18627f756511SDominic Meiser 18637f756511SDominic Meiser PetscFunctionBegin; 18647f756511SDominic Meiser if (*trifactors) { 1865470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtr); 1866470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtr); 1867470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose); 1868470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose); 18697f756511SDominic Meiser delete (*trifactors)->rpermIndices; 18707f756511SDominic Meiser delete (*trifactors)->cpermIndices; 18717f756511SDominic Meiser delete (*trifactors)->workVector; 18727f756511SDominic Meiser if (handle = (*trifactors)->handle) { 1873c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroy(handle);CHKERRCUDA(stat); 18747f756511SDominic Meiser } 18757f756511SDominic Meiser delete *trifactors; 18767f756511SDominic Meiser *trifactors = 0; 18777f756511SDominic Meiser } 18787f756511SDominic Meiser PetscFunctionReturn(0); 18797f756511SDominic Meiser } 18807f756511SDominic Meiser 1881