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 } 1920af67c1bSStefano 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 277580bdb30SBarry Smith ierr = PetscArraycpy(&(AjLo[offset]), vi, nz);CHKERRQ(ierr); 278580bdb30SBarry Smith ierr = PetscArraycpy(&(AALo[offset]), v, nz);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); 332*4863603aSSatish Balay ierr = PetscLogCpuToGpu((n+1+nzLower)*sizeof(int)+nzLower*sizeof(PetscScalar));CHKERRQ(ierr); 3339ae82921SPaul Mullowney } catch(char *ex) { 3349ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 3359ae82921SPaul Mullowney } 3369ae82921SPaul Mullowney } 3379ae82921SPaul Mullowney PetscFunctionReturn(0); 3389ae82921SPaul Mullowney } 3399ae82921SPaul Mullowney 340087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A) 3419ae82921SPaul Mullowney { 3429ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 3439ae82921SPaul Mullowney PetscInt n = A->rmap->n; 3449ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 345aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 3469ae82921SPaul Mullowney cusparseStatus_t stat; 3479ae82921SPaul Mullowney const PetscInt *aj = a->j,*adiag = a->diag,*vi; 3489ae82921SPaul Mullowney const MatScalar *aa = a->a,*v; 3499ae82921SPaul Mullowney PetscInt *AiUp, *AjUp; 3509ae82921SPaul Mullowney PetscScalar *AAUp; 3519ae82921SPaul Mullowney PetscInt i,nz, nzUpper, offset; 3529ae82921SPaul Mullowney PetscErrorCode ierr; 3539ae82921SPaul Mullowney 3549ae82921SPaul Mullowney PetscFunctionBegin; 355b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 3569ae82921SPaul Mullowney try { 3579ae82921SPaul Mullowney /* next, figure out the number of nonzeros in the upper triangular matrix. */ 3589ae82921SPaul Mullowney nzUpper = adiag[0]-adiag[n]; 3599ae82921SPaul Mullowney 3609ae82921SPaul Mullowney /* Allocate Space for the upper triangular matrix */ 361c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 362c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 363c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 3649ae82921SPaul Mullowney 3659ae82921SPaul Mullowney /* Fill the upper triangular matrix */ 3669ae82921SPaul Mullowney AiUp[0]=(PetscInt) 0; 3679ae82921SPaul Mullowney AiUp[n]=nzUpper; 3689ae82921SPaul Mullowney offset = nzUpper; 3699ae82921SPaul Mullowney for (i=n-1; i>=0; i--) { 3709ae82921SPaul Mullowney v = aa + adiag[i+1] + 1; 3719ae82921SPaul Mullowney vi = aj + adiag[i+1] + 1; 3729ae82921SPaul Mullowney 373e057df02SPaul Mullowney /* number of elements NOT on the diagonal */ 3749ae82921SPaul Mullowney nz = adiag[i] - adiag[i+1]-1; 3759ae82921SPaul Mullowney 376e057df02SPaul Mullowney /* decrement the offset */ 3779ae82921SPaul Mullowney offset -= (nz+1); 3789ae82921SPaul Mullowney 379e057df02SPaul Mullowney /* first, set the diagonal elements */ 3809ae82921SPaul Mullowney AjUp[offset] = (PetscInt) i; 38109f51544SAlejandro Lamas Daviña AAUp[offset] = (MatScalar)1./v[nz]; 3829ae82921SPaul Mullowney AiUp[i] = AiUp[i+1] - (nz+1); 3839ae82921SPaul Mullowney 384580bdb30SBarry Smith ierr = PetscArraycpy(&(AjUp[offset+1]), vi, nz);CHKERRQ(ierr); 385580bdb30SBarry Smith ierr = PetscArraycpy(&(AAUp[offset+1]), v, nz);CHKERRQ(ierr); 3869ae82921SPaul Mullowney } 3872205254eSKarl Rupp 388aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 389aa372e3fSPaul Mullowney upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 3902205254eSKarl Rupp 391aa372e3fSPaul Mullowney /* Create the matrix description */ 392c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 393c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 394c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 395c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 396c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 397aa372e3fSPaul Mullowney 398aa372e3fSPaul Mullowney /* Create the solve analysis information */ 399c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 400aa372e3fSPaul Mullowney 401aa372e3fSPaul Mullowney /* set the operation */ 402aa372e3fSPaul Mullowney upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 403aa372e3fSPaul Mullowney 404aa372e3fSPaul Mullowney /* set the matrix */ 405aa372e3fSPaul Mullowney upTriFactor->csrMat = new CsrMatrix; 406aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows = n; 407aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols = n; 408aa372e3fSPaul Mullowney upTriFactor->csrMat->num_entries = nzUpper; 409aa372e3fSPaul Mullowney 410aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1); 411aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1); 412aa372e3fSPaul Mullowney 413aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper); 414aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper); 415aa372e3fSPaul Mullowney 416aa372e3fSPaul Mullowney upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper); 417aa372e3fSPaul Mullowney upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper); 418aa372e3fSPaul Mullowney 419aa372e3fSPaul Mullowney /* perform the solve analysis */ 420aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 421aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 422aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 423c41cb2e2SAlejandro Lamas Daviña upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 424aa372e3fSPaul Mullowney 425aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 426aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 4272205254eSKarl Rupp 428c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 429c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 430c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 431*4863603aSSatish Balay ierr = PetscLogCpuToGpu((n+1+nzUpper)*sizeof(int)+nzUpper*sizeof(PetscScalar));CHKERRQ(ierr); 4329ae82921SPaul Mullowney } catch(char *ex) { 4339ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 4349ae82921SPaul Mullowney } 4359ae82921SPaul Mullowney } 4369ae82921SPaul Mullowney PetscFunctionReturn(0); 4379ae82921SPaul Mullowney } 4389ae82921SPaul Mullowney 439087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A) 4409ae82921SPaul Mullowney { 4419ae82921SPaul Mullowney PetscErrorCode ierr; 4429ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 4439ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 4449ae82921SPaul Mullowney IS isrow = a->row,iscol = a->icol; 4459ae82921SPaul Mullowney PetscBool row_identity,col_identity; 4469ae82921SPaul Mullowney const PetscInt *r,*c; 4479ae82921SPaul Mullowney PetscInt n = A->rmap->n; 4489ae82921SPaul Mullowney 4499ae82921SPaul Mullowney PetscFunctionBegin; 450087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr); 451087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr); 4522205254eSKarl Rupp 453e65717acSKarl Rupp cusparseTriFactors->workVector = new THRUSTARRAY(n); 454aa372e3fSPaul Mullowney cusparseTriFactors->nnz=a->nz; 4559ae82921SPaul Mullowney 456b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 457e057df02SPaul Mullowney /*lower triangular indices */ 4589ae82921SPaul Mullowney ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 4599ae82921SPaul Mullowney ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 4602205254eSKarl Rupp if (!row_identity) { 461aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 462aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices->assign(r, r+n); 4632205254eSKarl Rupp } 4649ae82921SPaul Mullowney ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 4659ae82921SPaul Mullowney 466e057df02SPaul Mullowney /*upper triangular indices */ 4679ae82921SPaul Mullowney ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr); 4689ae82921SPaul Mullowney ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 4692205254eSKarl Rupp if (!col_identity) { 470aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 471aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices->assign(c, c+n); 4722205254eSKarl Rupp } 473*4863603aSSatish Balay 474*4863603aSSatish Balay if(!row_identity && !col_identity) { 475*4863603aSSatish Balay ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr); 476*4863603aSSatish Balay } else if(!row_identity) { 477*4863603aSSatish Balay ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr); 478*4863603aSSatish Balay } else if(!col_identity) { 479*4863603aSSatish Balay ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr); 480*4863603aSSatish Balay } 481*4863603aSSatish Balay 4829ae82921SPaul Mullowney ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 4839ae82921SPaul Mullowney PetscFunctionReturn(0); 4849ae82921SPaul Mullowney } 4859ae82921SPaul Mullowney 486087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A) 487087f3262SPaul Mullowney { 488087f3262SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 489087f3262SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 490aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 491aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 492087f3262SPaul Mullowney cusparseStatus_t stat; 493087f3262SPaul Mullowney PetscErrorCode ierr; 494087f3262SPaul Mullowney PetscInt *AiUp, *AjUp; 495087f3262SPaul Mullowney PetscScalar *AAUp; 496087f3262SPaul Mullowney PetscScalar *AALo; 497087f3262SPaul Mullowney PetscInt nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j; 498087f3262SPaul Mullowney Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)A->data; 499087f3262SPaul Mullowney const PetscInt *ai = b->i,*aj = b->j,*vj; 500087f3262SPaul Mullowney const MatScalar *aa = b->a,*v; 501087f3262SPaul Mullowney 502087f3262SPaul Mullowney PetscFunctionBegin; 503b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 504087f3262SPaul Mullowney try { 505087f3262SPaul Mullowney /* Allocate Space for the upper triangular matrix */ 506c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr); 507c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr); 508c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 509c41cb2e2SAlejandro Lamas Daviña ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr); 510087f3262SPaul Mullowney 511087f3262SPaul Mullowney /* Fill the upper triangular matrix */ 512087f3262SPaul Mullowney AiUp[0]=(PetscInt) 0; 513087f3262SPaul Mullowney AiUp[n]=nzUpper; 514087f3262SPaul Mullowney offset = 0; 515087f3262SPaul Mullowney for (i=0; i<n; i++) { 516087f3262SPaul Mullowney /* set the pointers */ 517087f3262SPaul Mullowney v = aa + ai[i]; 518087f3262SPaul Mullowney vj = aj + ai[i]; 519087f3262SPaul Mullowney nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */ 520087f3262SPaul Mullowney 521087f3262SPaul Mullowney /* first, set the diagonal elements */ 522087f3262SPaul Mullowney AjUp[offset] = (PetscInt) i; 52309f51544SAlejandro Lamas Daviña AAUp[offset] = (MatScalar)1.0/v[nz]; 524087f3262SPaul Mullowney AiUp[i] = offset; 52509f51544SAlejandro Lamas Daviña AALo[offset] = (MatScalar)1.0/v[nz]; 526087f3262SPaul Mullowney 527087f3262SPaul Mullowney offset+=1; 528087f3262SPaul Mullowney if (nz>0) { 529f22e0265SBarry Smith ierr = PetscArraycpy(&(AjUp[offset]), vj, nz);CHKERRQ(ierr); 530580bdb30SBarry Smith ierr = PetscArraycpy(&(AAUp[offset]), v, nz);CHKERRQ(ierr); 531087f3262SPaul Mullowney for (j=offset; j<offset+nz; j++) { 532087f3262SPaul Mullowney AAUp[j] = -AAUp[j]; 533087f3262SPaul Mullowney AALo[j] = AAUp[j]/v[nz]; 534087f3262SPaul Mullowney } 535087f3262SPaul Mullowney offset+=nz; 536087f3262SPaul Mullowney } 537087f3262SPaul Mullowney } 538087f3262SPaul Mullowney 539aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 540aa372e3fSPaul Mullowney upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 541087f3262SPaul Mullowney 542aa372e3fSPaul Mullowney /* Create the matrix description */ 543c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat); 544c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 545c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 546c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 547c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat); 548087f3262SPaul Mullowney 549aa372e3fSPaul Mullowney /* Create the solve analysis information */ 550c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat); 551aa372e3fSPaul Mullowney 552aa372e3fSPaul Mullowney /* set the operation */ 553aa372e3fSPaul Mullowney upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 554aa372e3fSPaul Mullowney 555aa372e3fSPaul Mullowney /* set the matrix */ 556aa372e3fSPaul Mullowney upTriFactor->csrMat = new CsrMatrix; 557aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows = A->rmap->n; 558aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols = A->cmap->n; 559aa372e3fSPaul Mullowney upTriFactor->csrMat->num_entries = a->nz; 560aa372e3fSPaul Mullowney 561aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 562aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 563aa372e3fSPaul Mullowney 564aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 565aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 566aa372e3fSPaul Mullowney 567aa372e3fSPaul Mullowney upTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 568aa372e3fSPaul Mullowney upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz); 569aa372e3fSPaul Mullowney 570aa372e3fSPaul Mullowney /* perform the solve analysis */ 571aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp, 572aa372e3fSPaul Mullowney upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr, 573aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(), 574c41cb2e2SAlejandro Lamas Daviña upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat); 575aa372e3fSPaul Mullowney 576aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 577aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor; 578aa372e3fSPaul Mullowney 579aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 580aa372e3fSPaul Mullowney loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct; 581aa372e3fSPaul Mullowney 582aa372e3fSPaul Mullowney /* Create the matrix description */ 583c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat); 584c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 585c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat); 586c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat); 587c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat); 588aa372e3fSPaul Mullowney 589aa372e3fSPaul Mullowney /* Create the solve analysis information */ 590c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat); 591aa372e3fSPaul Mullowney 592aa372e3fSPaul Mullowney /* set the operation */ 593aa372e3fSPaul Mullowney loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE; 594aa372e3fSPaul Mullowney 595aa372e3fSPaul Mullowney /* set the matrix */ 596aa372e3fSPaul Mullowney loTriFactor->csrMat = new CsrMatrix; 597aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows = A->rmap->n; 598aa372e3fSPaul Mullowney loTriFactor->csrMat->num_cols = A->cmap->n; 599aa372e3fSPaul Mullowney loTriFactor->csrMat->num_entries = a->nz; 600aa372e3fSPaul Mullowney 601aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 602aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1); 603aa372e3fSPaul Mullowney 604aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz); 605aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz); 606aa372e3fSPaul Mullowney 607aa372e3fSPaul Mullowney loTriFactor->csrMat->values = new THRUSTARRAY(a->nz); 608aa372e3fSPaul Mullowney loTriFactor->csrMat->values->assign(AALo, AALo+a->nz); 609*4863603aSSatish Balay ierr = PetscLogCpuToGpu(2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)));CHKERRQ(ierr); 610aa372e3fSPaul Mullowney 611aa372e3fSPaul Mullowney /* perform the solve analysis */ 612aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp, 613aa372e3fSPaul Mullowney loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr, 614aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(), 615c41cb2e2SAlejandro Lamas Daviña loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat); 616aa372e3fSPaul Mullowney 617aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 618aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor; 619087f3262SPaul Mullowney 620b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 621c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr); 622c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr); 623c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr); 624c41cb2e2SAlejandro Lamas Daviña ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr); 625087f3262SPaul Mullowney } catch(char *ex) { 626087f3262SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 627087f3262SPaul Mullowney } 628087f3262SPaul Mullowney } 629087f3262SPaul Mullowney PetscFunctionReturn(0); 630087f3262SPaul Mullowney } 631087f3262SPaul Mullowney 632087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A) 6339ae82921SPaul Mullowney { 6349ae82921SPaul Mullowney PetscErrorCode ierr; 635087f3262SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 636087f3262SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 637087f3262SPaul Mullowney IS ip = a->row; 638087f3262SPaul Mullowney const PetscInt *rip; 639087f3262SPaul Mullowney PetscBool perm_identity; 640087f3262SPaul Mullowney PetscInt n = A->rmap->n; 641087f3262SPaul Mullowney 642087f3262SPaul Mullowney PetscFunctionBegin; 643087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr); 644e65717acSKarl Rupp cusparseTriFactors->workVector = new THRUSTARRAY(n); 645aa372e3fSPaul Mullowney cusparseTriFactors->nnz=(a->nz-n)*2 + n; 646aa372e3fSPaul Mullowney 647087f3262SPaul Mullowney /*lower triangular indices */ 648087f3262SPaul Mullowney ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 649087f3262SPaul Mullowney ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 650087f3262SPaul Mullowney if (!perm_identity) { 651aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n); 652aa372e3fSPaul Mullowney cusparseTriFactors->rpermIndices->assign(rip, rip+n); 653aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n); 654aa372e3fSPaul Mullowney cusparseTriFactors->cpermIndices->assign(rip, rip+n); 655*4863603aSSatish Balay ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr); 656087f3262SPaul Mullowney } 657087f3262SPaul Mullowney ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 658087f3262SPaul Mullowney PetscFunctionReturn(0); 659087f3262SPaul Mullowney } 660087f3262SPaul Mullowney 6616fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 6629ae82921SPaul Mullowney { 6639ae82921SPaul Mullowney Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 6649ae82921SPaul Mullowney IS isrow = b->row,iscol = b->col; 6659ae82921SPaul Mullowney PetscBool row_identity,col_identity; 666b175d8bbSPaul Mullowney PetscErrorCode ierr; 6679ae82921SPaul Mullowney 6689ae82921SPaul Mullowney PetscFunctionBegin; 6699ae82921SPaul Mullowney ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 670e057df02SPaul Mullowney /* determine which version of MatSolve needs to be used. */ 6719ae82921SPaul Mullowney ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 6729ae82921SPaul Mullowney ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 673bda325fcSPaul Mullowney if (row_identity && col_identity) { 674bda325fcSPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 675bda325fcSPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 676bda325fcSPaul Mullowney } else { 677bda325fcSPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE; 678bda325fcSPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 679bda325fcSPaul Mullowney } 6808dc1d2a3SPaul Mullowney 681e057df02SPaul Mullowney /* get the triangular factors */ 682087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 6839ae82921SPaul Mullowney PetscFunctionReturn(0); 6849ae82921SPaul Mullowney } 6859ae82921SPaul Mullowney 686087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info) 687087f3262SPaul Mullowney { 688087f3262SPaul Mullowney Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 689087f3262SPaul Mullowney IS ip = b->row; 690087f3262SPaul Mullowney PetscBool perm_identity; 691b175d8bbSPaul Mullowney PetscErrorCode ierr; 692087f3262SPaul Mullowney 693087f3262SPaul Mullowney PetscFunctionBegin; 694087f3262SPaul Mullowney ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr); 695087f3262SPaul Mullowney 696087f3262SPaul Mullowney /* determine which version of MatSolve needs to be used. */ 697087f3262SPaul Mullowney ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 698087f3262SPaul Mullowney if (perm_identity) { 699087f3262SPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering; 700087f3262SPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering; 701087f3262SPaul Mullowney } else { 702087f3262SPaul Mullowney B->ops->solve = MatSolve_SeqAIJCUSPARSE; 703087f3262SPaul Mullowney B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE; 704087f3262SPaul Mullowney } 705087f3262SPaul Mullowney 706087f3262SPaul Mullowney /* get the triangular factors */ 707087f3262SPaul Mullowney ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr); 708087f3262SPaul Mullowney PetscFunctionReturn(0); 709087f3262SPaul Mullowney } 7109ae82921SPaul Mullowney 711b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A) 712bda325fcSPaul Mullowney { 713bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 714aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 715aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 716aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 717aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 718bda325fcSPaul Mullowney cusparseStatus_t stat; 719aa372e3fSPaul Mullowney cusparseIndexBase_t indexBase; 720aa372e3fSPaul Mullowney cusparseMatrixType_t matrixType; 721aa372e3fSPaul Mullowney cusparseFillMode_t fillMode; 722aa372e3fSPaul Mullowney cusparseDiagType_t diagType; 723b175d8bbSPaul Mullowney 724bda325fcSPaul Mullowney PetscFunctionBegin; 725bda325fcSPaul Mullowney 726aa372e3fSPaul Mullowney /*********************************************/ 727aa372e3fSPaul Mullowney /* Now the Transpose of the Lower Tri Factor */ 728aa372e3fSPaul Mullowney /*********************************************/ 729aa372e3fSPaul Mullowney 730aa372e3fSPaul Mullowney /* allocate space for the transpose of the lower triangular factor */ 731aa372e3fSPaul Mullowney loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 732aa372e3fSPaul Mullowney 733aa372e3fSPaul Mullowney /* set the matrix descriptors of the lower triangular factor */ 734aa372e3fSPaul Mullowney matrixType = cusparseGetMatType(loTriFactor->descr); 735aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(loTriFactor->descr); 736aa372e3fSPaul Mullowney fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 737aa372e3fSPaul Mullowney CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 738aa372e3fSPaul Mullowney diagType = cusparseGetMatDiagType(loTriFactor->descr); 739aa372e3fSPaul Mullowney 740aa372e3fSPaul Mullowney /* Create the matrix description */ 741c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat); 742c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat); 743c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat); 744c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat); 745c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat); 746aa372e3fSPaul Mullowney 747aa372e3fSPaul Mullowney /* Create the solve analysis information */ 748c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat); 749aa372e3fSPaul Mullowney 750aa372e3fSPaul Mullowney /* set the operation */ 751aa372e3fSPaul Mullowney loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 752aa372e3fSPaul Mullowney 753aa372e3fSPaul Mullowney /* allocate GPU space for the CSC of the lower triangular factor*/ 754aa372e3fSPaul Mullowney loTriFactorT->csrMat = new CsrMatrix; 755aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows; 756aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols; 757aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries; 758aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1); 759aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries); 760aa372e3fSPaul Mullowney loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries); 761aa372e3fSPaul Mullowney 762aa372e3fSPaul Mullowney /* compute the transpose of the lower triangular factor, i.e. the CSC */ 763aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows, 764aa372e3fSPaul Mullowney loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries, 765aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 766aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 767aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 768aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 769aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 770aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 771c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 772aa372e3fSPaul Mullowney 773aa372e3fSPaul Mullowney /* perform the solve analysis on the transposed matrix */ 774aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp, 775aa372e3fSPaul Mullowney loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries, 776aa372e3fSPaul Mullowney loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(), 777aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(), 778c41cb2e2SAlejandro Lamas Daviña loTriFactorT->solveInfo);CHKERRCUDA(stat); 779aa372e3fSPaul Mullowney 780aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 781aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT; 782aa372e3fSPaul Mullowney 783aa372e3fSPaul Mullowney /*********************************************/ 784aa372e3fSPaul Mullowney /* Now the Transpose of the Upper Tri Factor */ 785aa372e3fSPaul Mullowney /*********************************************/ 786aa372e3fSPaul Mullowney 787aa372e3fSPaul Mullowney /* allocate space for the transpose of the upper triangular factor */ 788aa372e3fSPaul Mullowney upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct; 789aa372e3fSPaul Mullowney 790aa372e3fSPaul Mullowney /* set the matrix descriptors of the upper triangular factor */ 791aa372e3fSPaul Mullowney matrixType = cusparseGetMatType(upTriFactor->descr); 792aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(upTriFactor->descr); 793aa372e3fSPaul Mullowney fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ? 794aa372e3fSPaul Mullowney CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER; 795aa372e3fSPaul Mullowney diagType = cusparseGetMatDiagType(upTriFactor->descr); 796aa372e3fSPaul Mullowney 797aa372e3fSPaul Mullowney /* Create the matrix description */ 798c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat); 799c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat); 800c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat); 801c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat); 802c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat); 803aa372e3fSPaul Mullowney 804aa372e3fSPaul Mullowney /* Create the solve analysis information */ 805c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat); 806aa372e3fSPaul Mullowney 807aa372e3fSPaul Mullowney /* set the operation */ 808aa372e3fSPaul Mullowney upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE; 809aa372e3fSPaul Mullowney 810aa372e3fSPaul Mullowney /* allocate GPU space for the CSC of the upper triangular factor*/ 811aa372e3fSPaul Mullowney upTriFactorT->csrMat = new CsrMatrix; 812aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows; 813aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols; 814aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries; 815aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1); 816aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries); 817aa372e3fSPaul Mullowney upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries); 818aa372e3fSPaul Mullowney 819aa372e3fSPaul Mullowney /* compute the transpose of the upper triangular factor, i.e. the CSC */ 820aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows, 821aa372e3fSPaul Mullowney upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries, 822aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 823aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 824aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 825aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 826aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 827aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 828c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 829aa372e3fSPaul Mullowney 830aa372e3fSPaul Mullowney /* perform the solve analysis on the transposed matrix */ 831aa372e3fSPaul Mullowney stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp, 832aa372e3fSPaul Mullowney upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries, 833aa372e3fSPaul Mullowney upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(), 834aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(), 835c41cb2e2SAlejandro Lamas Daviña upTriFactorT->solveInfo);CHKERRCUDA(stat); 836aa372e3fSPaul Mullowney 837aa372e3fSPaul Mullowney /* assign the pointer. Is this really necessary? */ 838aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT; 839bda325fcSPaul Mullowney PetscFunctionReturn(0); 840bda325fcSPaul Mullowney } 841bda325fcSPaul Mullowney 842b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A) 843bda325fcSPaul Mullowney { 844aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 845aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 846aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 847bda325fcSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 848bda325fcSPaul Mullowney cusparseStatus_t stat; 849aa372e3fSPaul Mullowney cusparseIndexBase_t indexBase; 850b06137fdSPaul Mullowney cudaError_t err; 851*4863603aSSatish Balay PetscErrorCode ierr; 852b175d8bbSPaul Mullowney 853bda325fcSPaul Mullowney PetscFunctionBegin; 854aa372e3fSPaul Mullowney 855aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 856aa372e3fSPaul Mullowney matstructT = new Mat_SeqAIJCUSPARSEMultStruct; 857c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat); 858aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(matstruct->descr); 859c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat); 860c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 861aa372e3fSPaul Mullowney 862b06137fdSPaul Mullowney /* set alpha and beta */ 863c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstructT->alpha), sizeof(PetscScalar));CHKERRCUDA(err); 8647656d835SStefano Zampini err = cudaMalloc((void **)&(matstructT->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err); 8657656d835SStefano Zampini err = cudaMalloc((void **)&(matstructT->beta_one), sizeof(PetscScalar));CHKERRCUDA(err); 8667656d835SStefano Zampini err = cudaMemcpy(matstructT->alpha, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 8677656d835SStefano Zampini err = cudaMemcpy(matstructT->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 8687656d835SStefano Zampini err = cudaMemcpy(matstructT->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 869c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 870b06137fdSPaul Mullowney 871aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 872aa372e3fSPaul Mullowney CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 873aa372e3fSPaul Mullowney CsrMatrix *matrixT= new CsrMatrix; 874554b8892SKarl Rupp matrixT->num_rows = A->cmap->n; 875554b8892SKarl Rupp matrixT->num_cols = A->rmap->n; 876aa372e3fSPaul Mullowney matrixT->num_entries = a->nz; 877aa372e3fSPaul Mullowney matrixT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 878aa372e3fSPaul Mullowney matrixT->column_indices = new THRUSTINTARRAY32(a->nz); 879aa372e3fSPaul Mullowney matrixT->values = new THRUSTARRAY(a->nz); 880aa372e3fSPaul Mullowney 881aa372e3fSPaul Mullowney /* compute the transpose of the upper triangular factor, i.e. the CSC */ 882aa372e3fSPaul Mullowney indexBase = cusparseGetMatIndexBase(matstruct->descr); 883aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows, 884aa372e3fSPaul Mullowney matrix->num_cols, matrix->num_entries, 885aa372e3fSPaul Mullowney matrix->values->data().get(), 886aa372e3fSPaul Mullowney matrix->row_offsets->data().get(), 887aa372e3fSPaul Mullowney matrix->column_indices->data().get(), 888aa372e3fSPaul Mullowney matrixT->values->data().get(), 889aa372e3fSPaul Mullowney matrixT->column_indices->data().get(), 890aa372e3fSPaul Mullowney matrixT->row_offsets->data().get(), 891c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 892aa372e3fSPaul Mullowney 893aa372e3fSPaul Mullowney /* assign the pointer */ 894aa372e3fSPaul Mullowney matstructT->mat = matrixT; 895*4863603aSSatish Balay ierr = PetscLogCpuToGpu(((A->rmap->n+1)+(a->nz))*sizeof(int)+(3+a->nz)*sizeof(PetscScalar));CHKERRQ(ierr); 896aa372e3fSPaul Mullowney } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 8972692e278SPaul Mullowney #if CUDA_VERSION>=5000 898aa372e3fSPaul Mullowney /* First convert HYB to CSR */ 899aa372e3fSPaul Mullowney CsrMatrix *temp= new CsrMatrix; 900aa372e3fSPaul Mullowney temp->num_rows = A->rmap->n; 901aa372e3fSPaul Mullowney temp->num_cols = A->cmap->n; 902aa372e3fSPaul Mullowney temp->num_entries = a->nz; 903aa372e3fSPaul Mullowney temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 904aa372e3fSPaul Mullowney temp->column_indices = new THRUSTINTARRAY32(a->nz); 905aa372e3fSPaul Mullowney temp->values = new THRUSTARRAY(a->nz); 906aa372e3fSPaul Mullowney 9072692e278SPaul Mullowney 908aa372e3fSPaul Mullowney stat = cusparse_hyb2csr(cusparsestruct->handle, 909aa372e3fSPaul Mullowney matstruct->descr, (cusparseHybMat_t)matstruct->mat, 910aa372e3fSPaul Mullowney temp->values->data().get(), 911aa372e3fSPaul Mullowney temp->row_offsets->data().get(), 912c41cb2e2SAlejandro Lamas Daviña temp->column_indices->data().get());CHKERRCUDA(stat); 913aa372e3fSPaul Mullowney 914aa372e3fSPaul Mullowney /* Next, convert CSR to CSC (i.e. the matrix transpose) */ 915aa372e3fSPaul Mullowney CsrMatrix *tempT= new CsrMatrix; 916aa372e3fSPaul Mullowney tempT->num_rows = A->rmap->n; 917aa372e3fSPaul Mullowney tempT->num_cols = A->cmap->n; 918aa372e3fSPaul Mullowney tempT->num_entries = a->nz; 919aa372e3fSPaul Mullowney tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1); 920aa372e3fSPaul Mullowney tempT->column_indices = new THRUSTINTARRAY32(a->nz); 921aa372e3fSPaul Mullowney tempT->values = new THRUSTARRAY(a->nz); 922aa372e3fSPaul Mullowney 923aa372e3fSPaul Mullowney stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows, 924aa372e3fSPaul Mullowney temp->num_cols, temp->num_entries, 925aa372e3fSPaul Mullowney temp->values->data().get(), 926aa372e3fSPaul Mullowney temp->row_offsets->data().get(), 927aa372e3fSPaul Mullowney temp->column_indices->data().get(), 928aa372e3fSPaul Mullowney tempT->values->data().get(), 929aa372e3fSPaul Mullowney tempT->column_indices->data().get(), 930aa372e3fSPaul Mullowney tempT->row_offsets->data().get(), 931c41cb2e2SAlejandro Lamas Daviña CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat); 932aa372e3fSPaul Mullowney 933aa372e3fSPaul Mullowney /* Last, convert CSC to HYB */ 934aa372e3fSPaul Mullowney cusparseHybMat_t hybMat; 935c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 936aa372e3fSPaul Mullowney cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 937aa372e3fSPaul Mullowney CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 938aa372e3fSPaul Mullowney stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n, 939aa372e3fSPaul Mullowney matstructT->descr, tempT->values->data().get(), 940aa372e3fSPaul Mullowney tempT->row_offsets->data().get(), 941aa372e3fSPaul Mullowney tempT->column_indices->data().get(), 942c41cb2e2SAlejandro Lamas Daviña hybMat, 0, partition);CHKERRCUDA(stat); 943aa372e3fSPaul Mullowney 944aa372e3fSPaul Mullowney /* assign the pointer */ 945aa372e3fSPaul Mullowney matstructT->mat = hybMat; 946*4863603aSSatish Balay ierr = PetscLogCpuToGpu((2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)))+3*sizeof(PetscScalar));CHKERRQ(ierr); 947aa372e3fSPaul Mullowney 948aa372e3fSPaul Mullowney /* delete temporaries */ 949aa372e3fSPaul Mullowney if (tempT) { 950aa372e3fSPaul Mullowney if (tempT->values) delete (THRUSTARRAY*) tempT->values; 951aa372e3fSPaul Mullowney if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices; 952aa372e3fSPaul Mullowney if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets; 953aa372e3fSPaul Mullowney delete (CsrMatrix*) tempT; 954087f3262SPaul Mullowney } 955aa372e3fSPaul Mullowney if (temp) { 956aa372e3fSPaul Mullowney if (temp->values) delete (THRUSTARRAY*) temp->values; 957aa372e3fSPaul Mullowney if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices; 958aa372e3fSPaul Mullowney if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets; 959aa372e3fSPaul Mullowney delete (CsrMatrix*) temp; 960aa372e3fSPaul Mullowney } 9612692e278SPaul Mullowney #else 9622692e278SPaul 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."); 9632692e278SPaul Mullowney #endif 964aa372e3fSPaul Mullowney } 965aa372e3fSPaul Mullowney /* assign the compressed row indices */ 966aa372e3fSPaul Mullowney matstructT->cprowIndices = new THRUSTINTARRAY; 967554b8892SKarl Rupp matstructT->cprowIndices->resize(A->cmap->n); 968554b8892SKarl Rupp thrust::sequence(matstructT->cprowIndices->begin(), matstructT->cprowIndices->end()); 969aa372e3fSPaul Mullowney /* assign the pointer */ 970aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT; 971bda325fcSPaul Mullowney PetscFunctionReturn(0); 972bda325fcSPaul Mullowney } 973bda325fcSPaul Mullowney 9746fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 975bda325fcSPaul Mullowney { 976c41cb2e2SAlejandro Lamas Daviña PetscInt n = xx->map->n; 977465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 978465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 979465f34aeSAlejandro Lamas Daviña thrust::device_ptr<const PetscScalar> bGPU; 980465f34aeSAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> xGPU; 981bda325fcSPaul Mullowney cusparseStatus_t stat; 982bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 983aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 984aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 985aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 986b175d8bbSPaul Mullowney PetscErrorCode ierr; 987bda325fcSPaul Mullowney 988bda325fcSPaul Mullowney PetscFunctionBegin; 989aa372e3fSPaul Mullowney /* Analyze the matrix and create the transpose ... on the fly */ 990aa372e3fSPaul Mullowney if (!loTriFactorT && !upTriFactorT) { 991bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 992aa372e3fSPaul Mullowney loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 993aa372e3fSPaul Mullowney upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 994bda325fcSPaul Mullowney } 995bda325fcSPaul Mullowney 996bda325fcSPaul Mullowney /* Get the GPU pointers */ 997c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 998c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 999c41cb2e2SAlejandro Lamas Daviña xGPU = thrust::device_pointer_cast(xarray); 1000c41cb2e2SAlejandro Lamas Daviña bGPU = thrust::device_pointer_cast(barray); 1001bda325fcSPaul Mullowney 1002aa372e3fSPaul Mullowney /* First, reorder with the row permutation */ 1003c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 1004c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()), 1005c41cb2e2SAlejandro Lamas Daviña xGPU); 1006aa372e3fSPaul Mullowney 1007aa372e3fSPaul Mullowney /* First, solve U */ 1008aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 10097656d835SStefano Zampini upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr, 1010aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 1011aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 1012aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 1013aa372e3fSPaul Mullowney upTriFactorT->solveInfo, 1014c41cb2e2SAlejandro Lamas Daviña xarray, tempGPU->data().get());CHKERRCUDA(stat); 1015aa372e3fSPaul Mullowney 1016aa372e3fSPaul Mullowney /* Then, solve L */ 1017aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 10187656d835SStefano Zampini loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr, 1019aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 1020aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 1021aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 1022aa372e3fSPaul Mullowney loTriFactorT->solveInfo, 1023c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1024aa372e3fSPaul Mullowney 1025aa372e3fSPaul Mullowney /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1026c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1027c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()), 1028aa372e3fSPaul Mullowney tempGPU->begin()); 1029aa372e3fSPaul Mullowney 1030aa372e3fSPaul Mullowney /* Copy the temporary to the full solution. */ 1031c41cb2e2SAlejandro Lamas Daviña thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 1032bda325fcSPaul Mullowney 1033bda325fcSPaul Mullowney /* restore */ 1034c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1035c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1036c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1037087f3262SPaul Mullowney 1038aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1039bda325fcSPaul Mullowney PetscFunctionReturn(0); 1040bda325fcSPaul Mullowney } 1041bda325fcSPaul Mullowney 10426fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 1043bda325fcSPaul Mullowney { 1044465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1045465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 1046bda325fcSPaul Mullowney cusparseStatus_t stat; 1047bda325fcSPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1048aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1049aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1050aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1051b175d8bbSPaul Mullowney PetscErrorCode ierr; 1052bda325fcSPaul Mullowney 1053bda325fcSPaul Mullowney PetscFunctionBegin; 1054aa372e3fSPaul Mullowney /* Analyze the matrix and create the transpose ... on the fly */ 1055aa372e3fSPaul Mullowney if (!loTriFactorT && !upTriFactorT) { 1056bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr); 1057aa372e3fSPaul Mullowney loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose; 1058aa372e3fSPaul Mullowney upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose; 1059bda325fcSPaul Mullowney } 1060bda325fcSPaul Mullowney 1061bda325fcSPaul Mullowney /* Get the GPU pointers */ 1062c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1063c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1064bda325fcSPaul Mullowney 1065aa372e3fSPaul Mullowney /* First, solve U */ 1066aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp, 10677656d835SStefano Zampini upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr, 1068aa372e3fSPaul Mullowney upTriFactorT->csrMat->values->data().get(), 1069aa372e3fSPaul Mullowney upTriFactorT->csrMat->row_offsets->data().get(), 1070aa372e3fSPaul Mullowney upTriFactorT->csrMat->column_indices->data().get(), 1071aa372e3fSPaul Mullowney upTriFactorT->solveInfo, 1072c41cb2e2SAlejandro Lamas Daviña barray, tempGPU->data().get());CHKERRCUDA(stat); 1073aa372e3fSPaul Mullowney 1074aa372e3fSPaul Mullowney /* Then, solve L */ 1075aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp, 10767656d835SStefano Zampini loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr, 1077aa372e3fSPaul Mullowney loTriFactorT->csrMat->values->data().get(), 1078aa372e3fSPaul Mullowney loTriFactorT->csrMat->row_offsets->data().get(), 1079aa372e3fSPaul Mullowney loTriFactorT->csrMat->column_indices->data().get(), 1080aa372e3fSPaul Mullowney loTriFactorT->solveInfo, 1081c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1082bda325fcSPaul Mullowney 1083bda325fcSPaul Mullowney /* restore */ 1084c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1085c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1086c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1087aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 1088bda325fcSPaul Mullowney PetscFunctionReturn(0); 1089bda325fcSPaul Mullowney } 1090bda325fcSPaul Mullowney 10916fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx) 10929ae82921SPaul Mullowney { 1093465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1094465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 1095465f34aeSAlejandro Lamas Daviña thrust::device_ptr<const PetscScalar> bGPU; 1096465f34aeSAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> xGPU; 10979ae82921SPaul Mullowney cusparseStatus_t stat; 10989ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1099aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1100aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1101aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1102b175d8bbSPaul Mullowney PetscErrorCode ierr; 11039ae82921SPaul Mullowney 11049ae82921SPaul Mullowney PetscFunctionBegin; 1105ebc8f436SDominic Meiser 1106e057df02SPaul Mullowney /* Get the GPU pointers */ 1107c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1108c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 1109c41cb2e2SAlejandro Lamas Daviña xGPU = thrust::device_pointer_cast(xarray); 1110c41cb2e2SAlejandro Lamas Daviña bGPU = thrust::device_pointer_cast(barray); 11119ae82921SPaul Mullowney 1112aa372e3fSPaul Mullowney /* First, reorder with the row permutation */ 1113c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()), 1114c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()), 1115c41cb2e2SAlejandro Lamas Daviña xGPU); 1116aa372e3fSPaul Mullowney 1117aa372e3fSPaul Mullowney /* Next, solve L */ 1118aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 11197656d835SStefano Zampini loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr, 1120aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 1121aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 1122aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 1123aa372e3fSPaul Mullowney loTriFactor->solveInfo, 1124c41cb2e2SAlejandro Lamas Daviña xarray, tempGPU->data().get());CHKERRCUDA(stat); 1125aa372e3fSPaul Mullowney 1126aa372e3fSPaul Mullowney /* Then, solve U */ 1127aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 11287656d835SStefano Zampini upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr, 1129aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 1130aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 1131aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 1132aa372e3fSPaul Mullowney upTriFactor->solveInfo, 1133c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 1134aa372e3fSPaul Mullowney 1135aa372e3fSPaul Mullowney /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */ 1136c41cb2e2SAlejandro Lamas Daviña thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()), 1137c41cb2e2SAlejandro Lamas Daviña thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->end()), 1138aa372e3fSPaul Mullowney tempGPU->begin()); 1139aa372e3fSPaul Mullowney 1140aa372e3fSPaul Mullowney /* Copy the temporary to the full solution. */ 1141c41cb2e2SAlejandro Lamas Daviña thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU); 11429ae82921SPaul Mullowney 1143c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1144c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1145c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1146aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 11479ae82921SPaul Mullowney PetscFunctionReturn(0); 11489ae82921SPaul Mullowney } 11499ae82921SPaul Mullowney 11506fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx) 11519ae82921SPaul Mullowney { 1152465f34aeSAlejandro Lamas Daviña const PetscScalar *barray; 1153465f34aeSAlejandro Lamas Daviña PetscScalar *xarray; 11549ae82921SPaul Mullowney cusparseStatus_t stat; 11559ae82921SPaul Mullowney Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr; 1156aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr; 1157aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr; 1158aa372e3fSPaul Mullowney THRUSTARRAY *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector; 1159b175d8bbSPaul Mullowney PetscErrorCode ierr; 11609ae82921SPaul Mullowney 11619ae82921SPaul Mullowney PetscFunctionBegin; 1162e057df02SPaul Mullowney /* Get the GPU pointers */ 1163c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr); 1164c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr); 11659ae82921SPaul Mullowney 1166aa372e3fSPaul Mullowney /* First, solve L */ 1167aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp, 11687656d835SStefano Zampini loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr, 1169aa372e3fSPaul Mullowney loTriFactor->csrMat->values->data().get(), 1170aa372e3fSPaul Mullowney loTriFactor->csrMat->row_offsets->data().get(), 1171aa372e3fSPaul Mullowney loTriFactor->csrMat->column_indices->data().get(), 1172aa372e3fSPaul Mullowney loTriFactor->solveInfo, 1173c41cb2e2SAlejandro Lamas Daviña barray, tempGPU->data().get());CHKERRCUDA(stat); 1174aa372e3fSPaul Mullowney 1175aa372e3fSPaul Mullowney /* Next, solve U */ 1176aa372e3fSPaul Mullowney stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp, 11777656d835SStefano Zampini upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr, 1178aa372e3fSPaul Mullowney upTriFactor->csrMat->values->data().get(), 1179aa372e3fSPaul Mullowney upTriFactor->csrMat->row_offsets->data().get(), 1180aa372e3fSPaul Mullowney upTriFactor->csrMat->column_indices->data().get(), 1181aa372e3fSPaul Mullowney upTriFactor->solveInfo, 1182c41cb2e2SAlejandro Lamas Daviña tempGPU->data().get(), xarray);CHKERRCUDA(stat); 11839ae82921SPaul Mullowney 1184c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr); 1185c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr); 1186c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1187aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr); 11889ae82921SPaul Mullowney PetscFunctionReturn(0); 11899ae82921SPaul Mullowney } 11909ae82921SPaul Mullowney 11916fa9248bSJed Brown static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A) 11929ae82921SPaul Mullowney { 11939ae82921SPaul Mullowney 1194aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 1195aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 11969ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 11979ae82921SPaul Mullowney PetscInt m = A->rmap->n,*ii,*ridx; 11989ae82921SPaul Mullowney PetscErrorCode ierr; 1199aa372e3fSPaul Mullowney cusparseStatus_t stat; 1200b06137fdSPaul Mullowney cudaError_t err; 12019ae82921SPaul Mullowney 12029ae82921SPaul Mullowney PetscFunctionBegin; 1203b8ced49eSKarl Rupp if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) { 12049ae82921SPaul Mullowney ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 120534d6c7a5SJose E. Roman if (A->assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) { 120634d6c7a5SJose E. Roman CsrMatrix *matrix = (CsrMatrix*)matstruct->mat; 120734d6c7a5SJose E. Roman /* copy values only */ 120834d6c7a5SJose E. Roman matrix->values->assign(a->a, a->a+a->nz); 1209*4863603aSSatish Balay ierr = PetscLogCpuToGpu((a->nz)*sizeof(PetscScalar));CHKERRQ(ierr); 121034d6c7a5SJose E. Roman } else { 1211470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&matstruct,cusparsestruct->format); 12129ae82921SPaul Mullowney try { 1213aa372e3fSPaul Mullowney cusparsestruct->nonzerorow=0; 1214aa372e3fSPaul Mullowney for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0); 12159ae82921SPaul Mullowney 12169ae82921SPaul Mullowney if (a->compressedrow.use) { 12179ae82921SPaul Mullowney m = a->compressedrow.nrows; 12189ae82921SPaul Mullowney ii = a->compressedrow.i; 12199ae82921SPaul Mullowney ridx = a->compressedrow.rindex; 12209ae82921SPaul Mullowney } else { 1221b06137fdSPaul Mullowney /* Forcing compressed row on the GPU */ 12229ae82921SPaul Mullowney int k=0; 1223854ce69bSBarry Smith ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr); 1224854ce69bSBarry Smith ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr); 12259ae82921SPaul Mullowney ii[0]=0; 12269ae82921SPaul Mullowney for (int j = 0; j<m; j++) { 12279ae82921SPaul Mullowney if ((a->i[j+1]-a->i[j])>0) { 12289ae82921SPaul Mullowney ii[k] = a->i[j]; 12299ae82921SPaul Mullowney ridx[k]= j; 12309ae82921SPaul Mullowney k++; 12319ae82921SPaul Mullowney } 12329ae82921SPaul Mullowney } 1233aa372e3fSPaul Mullowney ii[cusparsestruct->nonzerorow] = a->nz; 1234aa372e3fSPaul Mullowney m = cusparsestruct->nonzerorow; 12359ae82921SPaul Mullowney } 12369ae82921SPaul Mullowney 1237aa372e3fSPaul Mullowney /* allocate space for the triangular factor information */ 1238aa372e3fSPaul Mullowney matstruct = new Mat_SeqAIJCUSPARSEMultStruct; 1239c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat); 1240c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat); 1241c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat); 12429ae82921SPaul Mullowney 1243c41cb2e2SAlejandro Lamas Daviña err = cudaMalloc((void **)&(matstruct->alpha), sizeof(PetscScalar));CHKERRCUDA(err); 12447656d835SStefano Zampini err = cudaMalloc((void **)&(matstruct->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err); 12457656d835SStefano Zampini err = cudaMalloc((void **)&(matstruct->beta_one), sizeof(PetscScalar));CHKERRCUDA(err); 12467656d835SStefano Zampini err = cudaMemcpy(matstruct->alpha, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 12477656d835SStefano Zampini err = cudaMemcpy(matstruct->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 12487656d835SStefano Zampini err = cudaMemcpy(matstruct->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err); 1249c41cb2e2SAlejandro Lamas Daviña stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat); 1250b06137fdSPaul Mullowney 1251aa372e3fSPaul Mullowney /* Build a hybrid/ellpack matrix if this option is chosen for the storage */ 1252aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1253aa372e3fSPaul Mullowney /* set the matrix */ 1254aa372e3fSPaul Mullowney CsrMatrix *matrix= new CsrMatrix; 1255a65300a6SPaul Mullowney matrix->num_rows = m; 1256aa372e3fSPaul Mullowney matrix->num_cols = A->cmap->n; 1257aa372e3fSPaul Mullowney matrix->num_entries = a->nz; 1258a65300a6SPaul Mullowney matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1259a65300a6SPaul Mullowney matrix->row_offsets->assign(ii, ii + m+1); 12609ae82921SPaul Mullowney 1261aa372e3fSPaul Mullowney matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1262aa372e3fSPaul Mullowney matrix->column_indices->assign(a->j, a->j+a->nz); 1263aa372e3fSPaul Mullowney 1264aa372e3fSPaul Mullowney matrix->values = new THRUSTARRAY(a->nz); 1265aa372e3fSPaul Mullowney matrix->values->assign(a->a, a->a+a->nz); 1266aa372e3fSPaul Mullowney 1267aa372e3fSPaul Mullowney /* assign the pointer */ 1268aa372e3fSPaul Mullowney matstruct->mat = matrix; 1269aa372e3fSPaul Mullowney 1270aa372e3fSPaul Mullowney } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) { 12712692e278SPaul Mullowney #if CUDA_VERSION>=4020 1272aa372e3fSPaul Mullowney CsrMatrix *matrix= new CsrMatrix; 1273a65300a6SPaul Mullowney matrix->num_rows = m; 1274aa372e3fSPaul Mullowney matrix->num_cols = A->cmap->n; 1275aa372e3fSPaul Mullowney matrix->num_entries = a->nz; 1276a65300a6SPaul Mullowney matrix->row_offsets = new THRUSTINTARRAY32(m+1); 1277a65300a6SPaul Mullowney matrix->row_offsets->assign(ii, ii + m+1); 1278aa372e3fSPaul Mullowney 1279aa372e3fSPaul Mullowney matrix->column_indices = new THRUSTINTARRAY32(a->nz); 1280aa372e3fSPaul Mullowney matrix->column_indices->assign(a->j, a->j+a->nz); 1281aa372e3fSPaul Mullowney 1282aa372e3fSPaul Mullowney matrix->values = new THRUSTARRAY(a->nz); 1283aa372e3fSPaul Mullowney matrix->values->assign(a->a, a->a+a->nz); 1284aa372e3fSPaul Mullowney 1285aa372e3fSPaul Mullowney cusparseHybMat_t hybMat; 1286c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat); 1287aa372e3fSPaul Mullowney cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ? 1288aa372e3fSPaul Mullowney CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO; 1289a65300a6SPaul Mullowney stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols, 1290aa372e3fSPaul Mullowney matstruct->descr, matrix->values->data().get(), 1291aa372e3fSPaul Mullowney matrix->row_offsets->data().get(), 1292aa372e3fSPaul Mullowney matrix->column_indices->data().get(), 1293c41cb2e2SAlejandro Lamas Daviña hybMat, 0, partition);CHKERRCUDA(stat); 1294aa372e3fSPaul Mullowney /* assign the pointer */ 1295aa372e3fSPaul Mullowney matstruct->mat = hybMat; 1296aa372e3fSPaul Mullowney 1297aa372e3fSPaul Mullowney if (matrix) { 1298aa372e3fSPaul Mullowney if (matrix->values) delete (THRUSTARRAY*)matrix->values; 1299aa372e3fSPaul Mullowney if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices; 1300aa372e3fSPaul Mullowney if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets; 1301aa372e3fSPaul Mullowney delete (CsrMatrix*)matrix; 1302087f3262SPaul Mullowney } 13032692e278SPaul Mullowney #endif 1304087f3262SPaul Mullowney } 1305ca45077fSPaul Mullowney 1306aa372e3fSPaul Mullowney /* assign the compressed row indices */ 1307aa372e3fSPaul Mullowney matstruct->cprowIndices = new THRUSTINTARRAY(m); 1308aa372e3fSPaul Mullowney matstruct->cprowIndices->assign(ridx,ridx+m); 1309*4863603aSSatish Balay ierr = PetscLogCpuToGpu(((m+1)+(a->nz))*sizeof(int)+m*sizeof(PetscInt)+(3+(a->nz))*sizeof(PetscScalar));CHKERRQ(ierr); 1310aa372e3fSPaul Mullowney 1311aa372e3fSPaul Mullowney /* assign the pointer */ 1312aa372e3fSPaul Mullowney cusparsestruct->mat = matstruct; 1313aa372e3fSPaul Mullowney 13149ae82921SPaul Mullowney if (!a->compressedrow.use) { 13159ae82921SPaul Mullowney ierr = PetscFree(ii);CHKERRQ(ierr); 13169ae82921SPaul Mullowney ierr = PetscFree(ridx);CHKERRQ(ierr); 13179ae82921SPaul Mullowney } 1318e65717acSKarl Rupp cusparsestruct->workVector = new THRUSTARRAY(m); 13199ae82921SPaul Mullowney } catch(char *ex) { 13209ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 13219ae82921SPaul Mullowney } 132234d6c7a5SJose E. Roman cusparsestruct->nonzerostate = A->nonzerostate; 132334d6c7a5SJose E. Roman } 13247d0e52d8SJose E. Roman ierr = WaitForGPU();CHKERRCUDA(ierr); 1325b8ced49eSKarl Rupp A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH; 13269ae82921SPaul Mullowney ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr); 13279ae82921SPaul Mullowney } 13289ae82921SPaul Mullowney PetscFunctionReturn(0); 13299ae82921SPaul Mullowney } 13309ae82921SPaul Mullowney 1331c41cb2e2SAlejandro Lamas Daviña struct VecCUDAPlusEquals 1332aa372e3fSPaul Mullowney { 1333aa372e3fSPaul Mullowney template <typename Tuple> 1334aa372e3fSPaul Mullowney __host__ __device__ 1335aa372e3fSPaul Mullowney void operator()(Tuple t) 1336aa372e3fSPaul Mullowney { 1337aa372e3fSPaul Mullowney thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t); 1338aa372e3fSPaul Mullowney } 1339aa372e3fSPaul Mullowney }; 1340aa372e3fSPaul Mullowney 13416fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 13429ae82921SPaul Mullowney { 1343b175d8bbSPaul Mullowney PetscErrorCode ierr; 13449ae82921SPaul Mullowney 13459ae82921SPaul Mullowney PetscFunctionBegin; 13467656d835SStefano Zampini ierr = MatMultAdd_SeqAIJCUSPARSE(A,xx,NULL,yy);CHKERRQ(ierr); 13479ae82921SPaul Mullowney PetscFunctionReturn(0); 13489ae82921SPaul Mullowney } 13499ae82921SPaul Mullowney 13506fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy) 1351ca45077fSPaul Mullowney { 1352ca45077fSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1353aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 13549ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstructT; 1355465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1356465f34aeSAlejandro Lamas Daviña PetscScalar *yarray; 1357b175d8bbSPaul Mullowney PetscErrorCode ierr; 1358aa372e3fSPaul Mullowney cusparseStatus_t stat; 1359ca45077fSPaul Mullowney 1360ca45077fSPaul Mullowney PetscFunctionBegin; 136134d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 136234d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 13639ff858a8SKarl Rupp matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1364aa372e3fSPaul Mullowney if (!matstructT) { 1365bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1366aa372e3fSPaul Mullowney matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1367bda325fcSPaul Mullowney } 1368c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1369ca6ae6e6SAlejandro Lamas Daviña ierr = VecSet(yy,0);CHKERRQ(ierr); 1370c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr); 1371aa372e3fSPaul Mullowney 1372aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1373aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1374aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1375aa372e3fSPaul Mullowney mat->num_rows, mat->num_cols, 1376b06137fdSPaul Mullowney mat->num_entries, matstructT->alpha, matstructT->descr, 1377aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 13787656d835SStefano Zampini mat->column_indices->data().get(), xarray, matstructT->beta_zero, 1379c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 1380aa372e3fSPaul Mullowney } else { 13812692e278SPaul Mullowney #if CUDA_VERSION>=4020 1382aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1383aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1384b06137fdSPaul Mullowney matstructT->alpha, matstructT->descr, hybMat, 13857656d835SStefano Zampini xarray, matstructT->beta_zero, 1386c41cb2e2SAlejandro Lamas Daviña yarray);CHKERRCUDA(stat); 13872692e278SPaul Mullowney #endif 1388ca45077fSPaul Mullowney } 1389c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1390c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr); 1391aa372e3fSPaul Mullowney if (!cusparsestruct->stream) { 1392c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1393ca45077fSPaul Mullowney } 1394aa372e3fSPaul Mullowney ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr); 1395ca45077fSPaul Mullowney PetscFunctionReturn(0); 1396ca45077fSPaul Mullowney } 1397ca45077fSPaul Mullowney 1398aa372e3fSPaul Mullowney 13996fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 14009ae82921SPaul Mullowney { 14019ae82921SPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1402aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 14039ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstruct; 1404465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 14057656d835SStefano Zampini PetscScalar *zarray,*dptr,*beta; 1406b175d8bbSPaul Mullowney PetscErrorCode ierr; 1407aa372e3fSPaul Mullowney cusparseStatus_t stat; 14086e111a19SKarl Rupp 14099ae82921SPaul Mullowney PetscFunctionBegin; 141034d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 141134d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 14129ff858a8SKarl Rupp matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat; 14139ae82921SPaul Mullowney try { 1414c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1415f2d70e9dSBarry Smith ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr); 14167656d835SStefano Zampini dptr = cusparsestruct->workVector->size() == (thrust::detail::vector_base<PetscScalar, thrust::device_malloc_allocator<PetscScalar> >::size_type)(A->rmap->n) ? zarray : cusparsestruct->workVector->data().get(); 14177656d835SStefano Zampini beta = (yy == zz && dptr == zarray) ? matstruct->beta_one : matstruct->beta_zero; 14189ae82921SPaul Mullowney 14197656d835SStefano Zampini /* csr_spmv is multiply add */ 1420aa372e3fSPaul Mullowney if (cusparsestruct->format == MAT_CUSPARSE_CSR) { 1421b06137fdSPaul Mullowney /* here we need to be careful to set the number of rows in the multiply to the 1422b06137fdSPaul Mullowney number of compressed rows in the matrix ... which is equivalent to the 1423b06137fdSPaul Mullowney size of the workVector */ 14247656d835SStefano Zampini CsrMatrix *mat = (CsrMatrix*)matstruct->mat; 1425aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1426a65300a6SPaul Mullowney mat->num_rows, mat->num_cols, 1427b06137fdSPaul Mullowney mat->num_entries, matstruct->alpha, matstruct->descr, 1428aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 14297656d835SStefano Zampini mat->column_indices->data().get(), xarray, beta, 14307656d835SStefano Zampini dptr);CHKERRCUDA(stat); 1431aa372e3fSPaul Mullowney } else { 14322692e278SPaul Mullowney #if CUDA_VERSION>=4020 1433aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat; 1434a65300a6SPaul Mullowney if (cusparsestruct->workVector->size()) { 1435aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1436b06137fdSPaul Mullowney matstruct->alpha, matstruct->descr, hybMat, 14377656d835SStefano Zampini xarray, beta, 14387656d835SStefano Zampini dptr);CHKERRCUDA(stat); 1439a65300a6SPaul Mullowney } 14402692e278SPaul Mullowney #endif 1441aa372e3fSPaul Mullowney } 1442aa372e3fSPaul Mullowney 14437656d835SStefano Zampini if (yy) { 14447656d835SStefano Zampini if (dptr != zarray) { 14457656d835SStefano Zampini ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 14467656d835SStefano Zampini } else if (zz != yy) { 14477656d835SStefano Zampini ierr = VecAXPY_SeqCUDA(zz,1.0,yy);CHKERRQ(ierr); 14487656d835SStefano Zampini } 14497656d835SStefano Zampini } else if (dptr != zarray) { 14507656d835SStefano Zampini ierr = VecSet(zz,0);CHKERRQ(ierr); 14517656d835SStefano Zampini } 1452aa372e3fSPaul Mullowney /* scatter the data from the temporary into the full vector with a += operation */ 14537656d835SStefano Zampini if (dptr != zarray) { 14547656d835SStefano Zampini thrust::device_ptr<PetscScalar> zptr; 14557656d835SStefano Zampini 14567656d835SStefano Zampini zptr = thrust::device_pointer_cast(zarray); 1457c41cb2e2SAlejandro Lamas Daviña thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))), 1458c41cb2e2SAlejandro Lamas Daviña thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(), 1459c41cb2e2SAlejandro Lamas Daviña VecCUDAPlusEquals()); 14607656d835SStefano Zampini } 1461c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1462f2d70e9dSBarry Smith ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr); 14639ae82921SPaul Mullowney } catch(char *ex) { 14649ae82921SPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 14659ae82921SPaul Mullowney } 14667656d835SStefano Zampini if (!yy) { /* MatMult */ 14677656d835SStefano Zampini if (!cusparsestruct->stream) { 1468c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 14697656d835SStefano Zampini } 14707656d835SStefano Zampini } 14719ae82921SPaul Mullowney ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 14729ae82921SPaul Mullowney PetscFunctionReturn(0); 14739ae82921SPaul Mullowney } 14749ae82921SPaul Mullowney 14756fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz) 1476ca45077fSPaul Mullowney { 1477ca45077fSPaul Mullowney Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1478aa372e3fSPaul Mullowney Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr; 14799ff858a8SKarl Rupp Mat_SeqAIJCUSPARSEMultStruct *matstructT; 1480c41cb2e2SAlejandro Lamas Daviña thrust::device_ptr<PetscScalar> zptr; 1481465f34aeSAlejandro Lamas Daviña const PetscScalar *xarray; 1482465f34aeSAlejandro Lamas Daviña PetscScalar *zarray; 1483b175d8bbSPaul Mullowney PetscErrorCode ierr; 1484aa372e3fSPaul Mullowney cusparseStatus_t stat; 14856e111a19SKarl Rupp 1486ca45077fSPaul Mullowney PetscFunctionBegin; 148734d6c7a5SJose E. Roman /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */ 148834d6c7a5SJose E. Roman ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 14899ff858a8SKarl Rupp matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1490aa372e3fSPaul Mullowney if (!matstructT) { 1491bda325fcSPaul Mullowney ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr); 1492aa372e3fSPaul Mullowney matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose; 1493bda325fcSPaul Mullowney } 1494aa372e3fSPaul Mullowney 1495ca45077fSPaul Mullowney try { 1496c41cb2e2SAlejandro Lamas Daviña ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr); 1497c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr); 1498f2d70e9dSBarry Smith ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr); 1499c41cb2e2SAlejandro Lamas Daviña zptr = thrust::device_pointer_cast(zarray); 1500ca45077fSPaul Mullowney 1501e057df02SPaul Mullowney /* multiply add with matrix transpose */ 1502aa372e3fSPaul Mullowney if (cusparsestruct->format==MAT_CUSPARSE_CSR) { 1503aa372e3fSPaul Mullowney CsrMatrix *mat = (CsrMatrix*)matstructT->mat; 1504b06137fdSPaul Mullowney /* here we need to be careful to set the number of rows in the multiply to the 1505b06137fdSPaul Mullowney number of compressed rows in the matrix ... which is equivalent to the 1506b06137fdSPaul Mullowney size of the workVector */ 1507aa372e3fSPaul Mullowney stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1508a65300a6SPaul Mullowney mat->num_rows, mat->num_cols, 1509b06137fdSPaul Mullowney mat->num_entries, matstructT->alpha, matstructT->descr, 1510aa372e3fSPaul Mullowney mat->values->data().get(), mat->row_offsets->data().get(), 15117656d835SStefano Zampini mat->column_indices->data().get(), xarray, matstructT->beta_zero, 1512c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1513aa372e3fSPaul Mullowney } else { 15142692e278SPaul Mullowney #if CUDA_VERSION>=4020 1515aa372e3fSPaul Mullowney cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat; 1516a65300a6SPaul Mullowney if (cusparsestruct->workVector->size()) { 1517aa372e3fSPaul Mullowney stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE, 1518b06137fdSPaul Mullowney matstructT->alpha, matstructT->descr, hybMat, 15197656d835SStefano Zampini xarray, matstructT->beta_zero, 1520c41cb2e2SAlejandro Lamas Daviña cusparsestruct->workVector->data().get());CHKERRCUDA(stat); 1521a65300a6SPaul Mullowney } 15222692e278SPaul Mullowney #endif 1523aa372e3fSPaul Mullowney } 1524aa372e3fSPaul Mullowney 1525aa372e3fSPaul Mullowney /* scatter the data from the temporary into the full vector with a += operation */ 1526c41cb2e2SAlejandro Lamas Daviña thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))), 1527554b8892SKarl Rupp thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))) + A->cmap->n, 1528c41cb2e2SAlejandro Lamas Daviña VecCUDAPlusEquals()); 1529ca45077fSPaul Mullowney 1530c41cb2e2SAlejandro Lamas Daviña ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr); 1531f2d70e9dSBarry Smith ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr); 1532ca45077fSPaul Mullowney 1533ca45077fSPaul Mullowney } catch(char *ex) { 1534ca45077fSPaul Mullowney SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex); 1535ca45077fSPaul Mullowney } 1536c41cb2e2SAlejandro Lamas Daviña ierr = WaitForGPU();CHKERRCUDA(ierr); 1537ca45077fSPaul Mullowney ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1538ca45077fSPaul Mullowney PetscFunctionReturn(0); 1539ca45077fSPaul Mullowney } 1540ca45077fSPaul Mullowney 15416fa9248bSJed Brown static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode) 15429ae82921SPaul Mullowney { 15439ae82921SPaul Mullowney PetscErrorCode ierr; 15446e111a19SKarl Rupp 15459ae82921SPaul Mullowney PetscFunctionBegin; 15469ae82921SPaul Mullowney ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr); 1547bc3f50f2SPaul Mullowney if (A->factortype==MAT_FACTOR_NONE) { 1548e057df02SPaul Mullowney ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); 1549bc3f50f2SPaul Mullowney } 15509ae82921SPaul Mullowney if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1551bbf3fe20SPaul Mullowney A->ops->mult = MatMult_SeqAIJCUSPARSE; 1552bbf3fe20SPaul Mullowney A->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1553bbf3fe20SPaul Mullowney A->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1554bbf3fe20SPaul Mullowney A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 15559ae82921SPaul Mullowney PetscFunctionReturn(0); 15569ae82921SPaul Mullowney } 15579ae82921SPaul Mullowney 15589ae82921SPaul Mullowney /* --------------------------------------------------------------------------------*/ 1559e057df02SPaul Mullowney /*@ 15609ae82921SPaul Mullowney MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format 1561e057df02SPaul Mullowney (the default parallel PETSc format). This matrix will ultimately pushed down 1562e057df02SPaul Mullowney to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix 1563e057df02SPaul Mullowney assembly performance the user should preallocate the matrix storage by setting 1564e057df02SPaul Mullowney the parameter nz (or the array nnz). By setting these parameters accurately, 1565e057df02SPaul Mullowney performance during matrix assembly can be increased by more than a factor of 50. 15669ae82921SPaul Mullowney 1567d083f849SBarry Smith Collective 15689ae82921SPaul Mullowney 15699ae82921SPaul Mullowney Input Parameters: 15709ae82921SPaul Mullowney + comm - MPI communicator, set to PETSC_COMM_SELF 15719ae82921SPaul Mullowney . m - number of rows 15729ae82921SPaul Mullowney . n - number of columns 15739ae82921SPaul Mullowney . nz - number of nonzeros per row (same for all rows) 15749ae82921SPaul Mullowney - nnz - array containing the number of nonzeros in the various rows 15750298fd71SBarry Smith (possibly different for each row) or NULL 15769ae82921SPaul Mullowney 15779ae82921SPaul Mullowney Output Parameter: 15789ae82921SPaul Mullowney . A - the matrix 15799ae82921SPaul Mullowney 15809ae82921SPaul Mullowney It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 15819ae82921SPaul Mullowney MatXXXXSetPreallocation() paradgm instead of this routine directly. 15829ae82921SPaul Mullowney [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 15839ae82921SPaul Mullowney 15849ae82921SPaul Mullowney Notes: 15859ae82921SPaul Mullowney If nnz is given then nz is ignored 15869ae82921SPaul Mullowney 15879ae82921SPaul Mullowney The AIJ format (also called the Yale sparse matrix format or 15889ae82921SPaul Mullowney compressed row storage), is fully compatible with standard Fortran 77 15899ae82921SPaul Mullowney storage. That is, the stored row and column indices can begin at 15909ae82921SPaul Mullowney either one (as in Fortran) or zero. See the users' manual for details. 15919ae82921SPaul Mullowney 15929ae82921SPaul Mullowney Specify the preallocated storage with either nz or nnz (not both). 15930298fd71SBarry Smith Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 15949ae82921SPaul Mullowney allocation. For large problems you MUST preallocate memory or you 15959ae82921SPaul Mullowney will get TERRIBLE performance, see the users' manual chapter on matrices. 15969ae82921SPaul Mullowney 15979ae82921SPaul Mullowney By default, this format uses inodes (identical nodes) when possible, to 15989ae82921SPaul Mullowney improve numerical efficiency of matrix-vector products and solves. We 15999ae82921SPaul Mullowney search for consecutive rows with the same nonzero structure, thereby 16009ae82921SPaul Mullowney reusing matrix information to achieve increased efficiency. 16019ae82921SPaul Mullowney 16029ae82921SPaul Mullowney Level: intermediate 16039ae82921SPaul Mullowney 1604e057df02SPaul Mullowney .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE 16059ae82921SPaul Mullowney @*/ 16069ae82921SPaul Mullowney PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 16079ae82921SPaul Mullowney { 16089ae82921SPaul Mullowney PetscErrorCode ierr; 16099ae82921SPaul Mullowney 16109ae82921SPaul Mullowney PetscFunctionBegin; 16119ae82921SPaul Mullowney ierr = MatCreate(comm,A);CHKERRQ(ierr); 16129ae82921SPaul Mullowney ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 16139ae82921SPaul Mullowney ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 16149ae82921SPaul Mullowney ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr); 16159ae82921SPaul Mullowney PetscFunctionReturn(0); 16169ae82921SPaul Mullowney } 16179ae82921SPaul Mullowney 16186fa9248bSJed Brown static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A) 16199ae82921SPaul Mullowney { 16209ae82921SPaul Mullowney PetscErrorCode ierr; 1621ab25e6cbSDominic Meiser 16229ae82921SPaul Mullowney PetscFunctionBegin; 16239ae82921SPaul Mullowney if (A->factortype==MAT_FACTOR_NONE) { 1624b8ced49eSKarl Rupp if (A->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) { 1625470880abSPatrick Sanan ierr = MatSeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr); 16269ae82921SPaul Mullowney } 16279ae82921SPaul Mullowney } else { 1628470880abSPatrick Sanan ierr = MatSeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr); 1629aa372e3fSPaul Mullowney } 16309ae82921SPaul Mullowney ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 16319ae82921SPaul Mullowney PetscFunctionReturn(0); 16329ae82921SPaul Mullowney } 16339ae82921SPaul Mullowney 16349ff858a8SKarl Rupp static PetscErrorCode MatDuplicate_SeqAIJCUSPARSE(Mat A,MatDuplicateOption cpvalues,Mat *B) 16359ff858a8SKarl Rupp { 16369ff858a8SKarl Rupp PetscErrorCode ierr; 16379ff858a8SKarl Rupp Mat C; 16389ff858a8SKarl Rupp cusparseStatus_t stat; 16399ff858a8SKarl Rupp cusparseHandle_t handle=0; 16409ff858a8SKarl Rupp 16419ff858a8SKarl Rupp PetscFunctionBegin; 16429ff858a8SKarl Rupp ierr = MatDuplicate_SeqAIJ(A,cpvalues,B);CHKERRQ(ierr); 16439ff858a8SKarl Rupp C = *B; 164434136279SStefano Zampini ierr = PetscFree(C->defaultvectype);CHKERRQ(ierr); 164534136279SStefano Zampini ierr = PetscStrallocpy(VECCUDA,&C->defaultvectype);CHKERRQ(ierr); 164634136279SStefano Zampini 16479ff858a8SKarl Rupp /* inject CUSPARSE-specific stuff */ 16489ff858a8SKarl Rupp if (C->factortype==MAT_FACTOR_NONE) { 16499ff858a8SKarl Rupp /* you cannot check the inode.use flag here since the matrix was just created. 16509ff858a8SKarl Rupp now build a GPU matrix data structure */ 16519ff858a8SKarl Rupp C->spptr = new Mat_SeqAIJCUSPARSE; 16529ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->mat = 0; 16539ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->matTranspose = 0; 16549ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->workVector = 0; 16559ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->format = MAT_CUSPARSE_CSR; 16569ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream = 0; 16579ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle = 0; 16589ff858a8SKarl Rupp stat = cusparseCreate(&handle);CHKERRCUDA(stat); 16599ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle = handle; 16609ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream = 0; 16619ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)C->spptr)->nonzerostate = 0; 16629ff858a8SKarl Rupp } else { 16639ff858a8SKarl Rupp /* NEXT, set the pointers to the triangular factors */ 16649ff858a8SKarl Rupp C->spptr = new Mat_SeqAIJCUSPARSETriFactors; 16659ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtr = 0; 16669ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtr = 0; 16679ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtrTranspose = 0; 16689ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtrTranspose = 0; 16699ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->rpermIndices = 0; 16709ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->cpermIndices = 0; 16719ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->workVector = 0; 16729ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle = 0; 16739ff858a8SKarl Rupp stat = cusparseCreate(&handle);CHKERRCUDA(stat); 16749ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle = handle; 16759ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->nnz = 0; 16769ff858a8SKarl Rupp } 16779ff858a8SKarl Rupp 16789ff858a8SKarl Rupp C->ops->assemblyend = MatAssemblyEnd_SeqAIJCUSPARSE; 16799ff858a8SKarl Rupp C->ops->destroy = MatDestroy_SeqAIJCUSPARSE; 16809ff858a8SKarl Rupp C->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSPARSE; 16819ff858a8SKarl Rupp C->ops->mult = MatMult_SeqAIJCUSPARSE; 16829ff858a8SKarl Rupp C->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 16839ff858a8SKarl Rupp C->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 16849ff858a8SKarl Rupp C->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 16859ff858a8SKarl Rupp C->ops->duplicate = MatDuplicate_SeqAIJCUSPARSE; 16869ff858a8SKarl Rupp 16879ff858a8SKarl Rupp ierr = PetscObjectChangeTypeName((PetscObject)C,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 16889ff858a8SKarl Rupp 1689b8ced49eSKarl Rupp C->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED; 16909ff858a8SKarl Rupp 16919ff858a8SKarl Rupp ierr = PetscObjectComposeFunction((PetscObject)C, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr); 16929ff858a8SKarl Rupp PetscFunctionReturn(0); 16939ff858a8SKarl Rupp } 16949ff858a8SKarl Rupp 169502fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat B) 16969ae82921SPaul Mullowney { 16979ae82921SPaul Mullowney PetscErrorCode ierr; 1698aa372e3fSPaul Mullowney cusparseStatus_t stat; 1699aa372e3fSPaul Mullowney cusparseHandle_t handle=0; 17009ae82921SPaul Mullowney 17019ae82921SPaul Mullowney PetscFunctionBegin; 170234136279SStefano Zampini ierr = PetscFree(B->defaultvectype);CHKERRQ(ierr); 170334136279SStefano Zampini ierr = PetscStrallocpy(VECCUDA,&B->defaultvectype);CHKERRQ(ierr); 170434136279SStefano Zampini 17059ae82921SPaul Mullowney if (B->factortype==MAT_FACTOR_NONE) { 1706e057df02SPaul Mullowney /* you cannot check the inode.use flag here since the matrix was just created. 1707e057df02SPaul Mullowney now build a GPU matrix data structure */ 17089ae82921SPaul Mullowney B->spptr = new Mat_SeqAIJCUSPARSE; 17099ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat = 0; 1710aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0; 1711aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector = 0; 1712e057df02SPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->format = MAT_CUSPARSE_CSR; 1713aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 1714aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = 0; 1715c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1716aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle = handle; 1717aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream = 0; 17189ff858a8SKarl Rupp ((Mat_SeqAIJCUSPARSE*)B->spptr)->nonzerostate = 0; 17199ae82921SPaul Mullowney } else { 17209ae82921SPaul Mullowney /* NEXT, set the pointers to the triangular factors */ 1721debe9ee2SPaul Mullowney B->spptr = new Mat_SeqAIJCUSPARSETriFactors; 17229ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr = 0; 17239ae82921SPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr = 0; 1724aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0; 1725aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0; 1726aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices = 0; 1727aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices = 0; 1728aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector = 0; 1729aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = 0; 1730c41cb2e2SAlejandro Lamas Daviña stat = cusparseCreate(&handle);CHKERRCUDA(stat); 1731aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle = handle; 1732aa372e3fSPaul Mullowney ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz = 0; 17339ae82921SPaul Mullowney } 1734aa372e3fSPaul Mullowney 17359ae82921SPaul Mullowney B->ops->assemblyend = MatAssemblyEnd_SeqAIJCUSPARSE; 17369ae82921SPaul Mullowney B->ops->destroy = MatDestroy_SeqAIJCUSPARSE; 17379ae82921SPaul Mullowney B->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSPARSE; 1738ca45077fSPaul Mullowney B->ops->mult = MatMult_SeqAIJCUSPARSE; 1739ca45077fSPaul Mullowney B->ops->multadd = MatMultAdd_SeqAIJCUSPARSE; 1740ca45077fSPaul Mullowney B->ops->multtranspose = MatMultTranspose_SeqAIJCUSPARSE; 1741ca45077fSPaul Mullowney B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE; 17429ff858a8SKarl Rupp B->ops->duplicate = MatDuplicate_SeqAIJCUSPARSE; 17432205254eSKarl Rupp 17449ae82921SPaul Mullowney ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr); 17452205254eSKarl Rupp 1746b8ced49eSKarl Rupp B->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED; 17472205254eSKarl Rupp 1748bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr); 17499ae82921SPaul Mullowney PetscFunctionReturn(0); 17509ae82921SPaul Mullowney } 17519ae82921SPaul Mullowney 175202fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B) 175302fe1965SBarry Smith { 175402fe1965SBarry Smith PetscErrorCode ierr; 175502fe1965SBarry Smith 175602fe1965SBarry Smith PetscFunctionBegin; 175702fe1965SBarry Smith ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr); 175802fe1965SBarry Smith ierr = MatConvert_SeqAIJ_SeqAIJCUSPARSE(B); 175902fe1965SBarry Smith PetscFunctionReturn(0); 176002fe1965SBarry Smith } 176102fe1965SBarry Smith 17623ca39a21SBarry Smith /*MC 1763e057df02SPaul Mullowney MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices. 1764e057df02SPaul Mullowney 1765e057df02SPaul Mullowney A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either 17662692e278SPaul Mullowney CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later. 17672692e278SPaul Mullowney All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library. 1768e057df02SPaul Mullowney 1769e057df02SPaul Mullowney Options Database Keys: 1770e057df02SPaul Mullowney + -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions() 1771aa372e3fSPaul 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). 1772aa372e3fSPaul 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). 1773e057df02SPaul Mullowney 1774e057df02SPaul Mullowney Level: beginner 1775e057df02SPaul Mullowney 17768468deeeSKarl Rupp .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation 1777e057df02SPaul Mullowney M*/ 17787f756511SDominic Meiser 177942c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*); 178042c9c57cSBarry Smith 17810f39cd5aSBarry Smith 17823ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_CUSPARSE(void) 178342c9c57cSBarry Smith { 178442c9c57cSBarry Smith PetscErrorCode ierr; 178542c9c57cSBarry Smith 178642c9c57cSBarry Smith PetscFunctionBegin; 17873ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17883ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17893ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 17903ca39a21SBarry Smith ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr); 179142c9c57cSBarry Smith PetscFunctionReturn(0); 179242c9c57cSBarry Smith } 179329b38603SBarry Smith 179481e08676SBarry Smith 1795470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct) 17967f756511SDominic Meiser { 17977f756511SDominic Meiser cusparseStatus_t stat; 17987f756511SDominic Meiser cusparseHandle_t handle; 17997f756511SDominic Meiser 18007f756511SDominic Meiser PetscFunctionBegin; 18017f756511SDominic Meiser if (*cusparsestruct) { 1802470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format); 1803470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format); 18047f756511SDominic Meiser delete (*cusparsestruct)->workVector; 18057f756511SDominic Meiser if (handle = (*cusparsestruct)->handle) { 1806c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroy(handle);CHKERRCUDA(stat); 18077f756511SDominic Meiser } 18087f756511SDominic Meiser delete *cusparsestruct; 18097f756511SDominic Meiser *cusparsestruct = 0; 18107f756511SDominic Meiser } 18117f756511SDominic Meiser PetscFunctionReturn(0); 18127f756511SDominic Meiser } 18137f756511SDominic Meiser 18147f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat) 18157f756511SDominic Meiser { 18167f756511SDominic Meiser PetscFunctionBegin; 18177f756511SDominic Meiser if (*mat) { 18187f756511SDominic Meiser delete (*mat)->values; 18197f756511SDominic Meiser delete (*mat)->column_indices; 18207f756511SDominic Meiser delete (*mat)->row_offsets; 18217f756511SDominic Meiser delete *mat; 18227f756511SDominic Meiser *mat = 0; 18237f756511SDominic Meiser } 18247f756511SDominic Meiser PetscFunctionReturn(0); 18257f756511SDominic Meiser } 18267f756511SDominic Meiser 1827470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor) 18287f756511SDominic Meiser { 18297f756511SDominic Meiser cusparseStatus_t stat; 18307f756511SDominic Meiser PetscErrorCode ierr; 18317f756511SDominic Meiser 18327f756511SDominic Meiser PetscFunctionBegin; 18337f756511SDominic Meiser if (*trifactor) { 1834c41cb2e2SAlejandro Lamas Daviña if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); } 1835c41cb2e2SAlejandro Lamas Daviña if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); } 18367f756511SDominic Meiser ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr); 18377f756511SDominic Meiser delete *trifactor; 18387f756511SDominic Meiser *trifactor = 0; 18397f756511SDominic Meiser } 18407f756511SDominic Meiser PetscFunctionReturn(0); 18417f756511SDominic Meiser } 18427f756511SDominic Meiser 1843470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format) 18447f756511SDominic Meiser { 18457f756511SDominic Meiser CsrMatrix *mat; 18467f756511SDominic Meiser cusparseStatus_t stat; 18477f756511SDominic Meiser cudaError_t err; 18487f756511SDominic Meiser 18497f756511SDominic Meiser PetscFunctionBegin; 18507f756511SDominic Meiser if (*matstruct) { 18517f756511SDominic Meiser if ((*matstruct)->mat) { 18527f756511SDominic Meiser if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) { 18537f756511SDominic Meiser cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat; 1854c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat); 18557f756511SDominic Meiser } else { 18567f756511SDominic Meiser mat = (CsrMatrix*)(*matstruct)->mat; 18577f756511SDominic Meiser CsrMatrix_Destroy(&mat); 18587f756511SDominic Meiser } 18597f756511SDominic Meiser } 1860c41cb2e2SAlejandro Lamas Daviña if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); } 18617f756511SDominic Meiser delete (*matstruct)->cprowIndices; 1862c41cb2e2SAlejandro Lamas Daviña if ((*matstruct)->alpha) { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); } 18637656d835SStefano Zampini if ((*matstruct)->beta_zero) { err=cudaFree((*matstruct)->beta_zero);CHKERRCUDA(err); } 18647656d835SStefano Zampini if ((*matstruct)->beta_one) { err=cudaFree((*matstruct)->beta_one);CHKERRCUDA(err); } 18657f756511SDominic Meiser delete *matstruct; 18667f756511SDominic Meiser *matstruct = 0; 18677f756511SDominic Meiser } 18687f756511SDominic Meiser PetscFunctionReturn(0); 18697f756511SDominic Meiser } 18707f756511SDominic Meiser 1871470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors) 18727f756511SDominic Meiser { 18737f756511SDominic Meiser cusparseHandle_t handle; 18747f756511SDominic Meiser cusparseStatus_t stat; 18757f756511SDominic Meiser 18767f756511SDominic Meiser PetscFunctionBegin; 18777f756511SDominic Meiser if (*trifactors) { 1878470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtr); 1879470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtr); 1880470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose); 1881470880abSPatrick Sanan MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose); 18827f756511SDominic Meiser delete (*trifactors)->rpermIndices; 18837f756511SDominic Meiser delete (*trifactors)->cpermIndices; 18847f756511SDominic Meiser delete (*trifactors)->workVector; 18857f756511SDominic Meiser if (handle = (*trifactors)->handle) { 1886c41cb2e2SAlejandro Lamas Daviña stat = cusparseDestroy(handle);CHKERRCUDA(stat); 18877f756511SDominic Meiser } 18887f756511SDominic Meiser delete *trifactors; 18897f756511SDominic Meiser *trifactors = 0; 18907f756511SDominic Meiser } 18917f756511SDominic Meiser PetscFunctionReturn(0); 18927f756511SDominic Meiser } 18937f756511SDominic Meiser 1894