xref: /petsc/src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu (revision a8bd5306f2fe3d0bceaf87eb6e700b41c78cecdd)
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;
127ca45077fSPaul Mullowney   switch (op) {
128e057df02SPaul Mullowney   case MAT_CUSPARSE_MULT:
129aa372e3fSPaul Mullowney     cusparsestruct->format = format;
130ca45077fSPaul Mullowney     break;
131e057df02SPaul Mullowney   case MAT_CUSPARSE_ALL:
132aa372e3fSPaul Mullowney     cusparsestruct->format = format;
133ca45077fSPaul Mullowney     break;
134ca45077fSPaul Mullowney   default:
13536d62e41SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op);
136ca45077fSPaul Mullowney   }
137ca45077fSPaul Mullowney   PetscFunctionReturn(0);
138ca45077fSPaul Mullowney }
1399ae82921SPaul Mullowney 
140e057df02SPaul Mullowney /*@
141e057df02SPaul Mullowney    MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular
142e057df02SPaul Mullowney    operation. Only the MatMult operation can use different GPU storage formats
143aa372e3fSPaul Mullowney    for MPIAIJCUSPARSE matrices.
144e057df02SPaul Mullowney    Not Collective
145e057df02SPaul Mullowney 
146e057df02SPaul Mullowney    Input Parameters:
1478468deeeSKarl Rupp +  A - Matrix of type SEQAIJCUSPARSE
14836d62e41SPaul 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.
1492692e278SPaul Mullowney -  format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2)
150e057df02SPaul Mullowney 
151e057df02SPaul Mullowney    Output Parameter:
152e057df02SPaul Mullowney 
153e057df02SPaul Mullowney    Level: intermediate
154e057df02SPaul Mullowney 
1558468deeeSKarl Rupp .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
156e057df02SPaul Mullowney @*/
157e057df02SPaul Mullowney PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
158e057df02SPaul Mullowney {
159e057df02SPaul Mullowney   PetscErrorCode ierr;
1606e111a19SKarl Rupp 
161e057df02SPaul Mullowney   PetscFunctionBegin;
162e057df02SPaul Mullowney   PetscValidHeaderSpecific(A, MAT_CLASSID,1);
163e057df02SPaul Mullowney   ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr);
164e057df02SPaul Mullowney   PetscFunctionReturn(0);
165e057df02SPaul Mullowney }
166e057df02SPaul Mullowney 
1674416b707SBarry Smith static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat A)
1689ae82921SPaul Mullowney {
1699ae82921SPaul Mullowney   PetscErrorCode           ierr;
170e057df02SPaul Mullowney   MatCUSPARSEStorageFormat format;
1719ae82921SPaul Mullowney   PetscBool                flg;
172a183c035SDominic Meiser   Mat_SeqAIJCUSPARSE       *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1736e111a19SKarl Rupp 
1749ae82921SPaul Mullowney   PetscFunctionBegin;
175e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSPARSE options");CHKERRQ(ierr);
1769ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
177e057df02SPaul Mullowney     ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV",
178a183c035SDominic Meiser                             "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
179e057df02SPaul Mullowney     if (flg) {
180e057df02SPaul Mullowney       ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr);
181045c96e1SPaul Mullowney     }
1829ae82921SPaul Mullowney   }
1834c87dfd4SPaul Mullowney   ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve",
184a183c035SDominic Meiser                           "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
1854c87dfd4SPaul Mullowney   if (flg) {
1864c87dfd4SPaul Mullowney     ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr);
1874c87dfd4SPaul Mullowney   }
1880af67c1bSStefano Zampini   ierr = PetscOptionsTail();CHKERRQ(ierr);
1899ae82921SPaul Mullowney   PetscFunctionReturn(0);
1909ae82921SPaul Mullowney 
1919ae82921SPaul Mullowney }
1929ae82921SPaul Mullowney 
1936fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1949ae82921SPaul Mullowney {
1959ae82921SPaul Mullowney   PetscErrorCode ierr;
1969ae82921SPaul Mullowney 
1979ae82921SPaul Mullowney   PetscFunctionBegin;
1989ae82921SPaul Mullowney   ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
1999ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2009ae82921SPaul Mullowney   PetscFunctionReturn(0);
2019ae82921SPaul Mullowney }
2029ae82921SPaul Mullowney 
2036fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2049ae82921SPaul Mullowney {
2059ae82921SPaul Mullowney   PetscErrorCode ierr;
2069ae82921SPaul Mullowney 
2079ae82921SPaul Mullowney   PetscFunctionBegin;
2089ae82921SPaul Mullowney   ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
2099ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2109ae82921SPaul Mullowney   PetscFunctionReturn(0);
2119ae82921SPaul Mullowney }
2129ae82921SPaul Mullowney 
213087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
214087f3262SPaul Mullowney {
215087f3262SPaul Mullowney   PetscErrorCode ierr;
216087f3262SPaul Mullowney 
217087f3262SPaul Mullowney   PetscFunctionBegin;
218087f3262SPaul Mullowney   ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
219087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
220087f3262SPaul Mullowney   PetscFunctionReturn(0);
221087f3262SPaul Mullowney }
222087f3262SPaul Mullowney 
223087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
224087f3262SPaul Mullowney {
225087f3262SPaul Mullowney   PetscErrorCode ierr;
226087f3262SPaul Mullowney 
227087f3262SPaul Mullowney   PetscFunctionBegin;
228087f3262SPaul Mullowney   ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
229087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
230087f3262SPaul Mullowney   PetscFunctionReturn(0);
231087f3262SPaul Mullowney }
232087f3262SPaul Mullowney 
233087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A)
2349ae82921SPaul Mullowney {
2359ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
2369ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
2379ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
238aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
2399ae82921SPaul Mullowney   cusparseStatus_t                  stat;
2409ae82921SPaul Mullowney   const PetscInt                    *ai = a->i,*aj = a->j,*vi;
2419ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
2429ae82921SPaul Mullowney   PetscInt                          *AiLo, *AjLo;
2439ae82921SPaul Mullowney   PetscScalar                       *AALo;
2449ae82921SPaul Mullowney   PetscInt                          i,nz, nzLower, offset, rowOffset;
245b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
2469ae82921SPaul Mullowney 
2479ae82921SPaul Mullowney   PetscFunctionBegin;
248b8ced49eSKarl Rupp   if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) {
2499ae82921SPaul Mullowney     try {
2509ae82921SPaul Mullowney       /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */
2519ae82921SPaul Mullowney       nzLower=n+ai[n]-ai[1];
2529ae82921SPaul Mullowney 
2539ae82921SPaul Mullowney       /* Allocate Space for the lower triangular matrix */
254c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
255c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUDA(ierr);
256c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUDA(ierr);
2579ae82921SPaul Mullowney 
2589ae82921SPaul Mullowney       /* Fill the lower triangular matrix */
2599ae82921SPaul Mullowney       AiLo[0]  = (PetscInt) 0;
2609ae82921SPaul Mullowney       AiLo[n]  = nzLower;
2619ae82921SPaul Mullowney       AjLo[0]  = (PetscInt) 0;
2629ae82921SPaul Mullowney       AALo[0]  = (MatScalar) 1.0;
2639ae82921SPaul Mullowney       v        = aa;
2649ae82921SPaul Mullowney       vi       = aj;
2659ae82921SPaul Mullowney       offset   = 1;
2669ae82921SPaul Mullowney       rowOffset= 1;
2679ae82921SPaul Mullowney       for (i=1; i<n; i++) {
2689ae82921SPaul Mullowney         nz = ai[i+1] - ai[i];
269e057df02SPaul Mullowney         /* additional 1 for the term on the diagonal */
2709ae82921SPaul Mullowney         AiLo[i]    = rowOffset;
2719ae82921SPaul Mullowney         rowOffset += nz+1;
2729ae82921SPaul Mullowney 
273580bdb30SBarry Smith         ierr = PetscArraycpy(&(AjLo[offset]), vi, nz);CHKERRQ(ierr);
274580bdb30SBarry Smith         ierr = PetscArraycpy(&(AALo[offset]), v, nz);CHKERRQ(ierr);
2759ae82921SPaul Mullowney 
2769ae82921SPaul Mullowney         offset      += nz;
2779ae82921SPaul Mullowney         AjLo[offset] = (PetscInt) i;
2789ae82921SPaul Mullowney         AALo[offset] = (MatScalar) 1.0;
2799ae82921SPaul Mullowney         offset      += 1;
2809ae82921SPaul Mullowney 
2819ae82921SPaul Mullowney         v  += nz;
2829ae82921SPaul Mullowney         vi += nz;
2839ae82921SPaul Mullowney       }
2842205254eSKarl Rupp 
285aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
286aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
2872205254eSKarl Rupp 
288aa372e3fSPaul Mullowney       /* Create the matrix description */
289c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
290c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
291c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
292c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUDA(stat);
293c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
294aa372e3fSPaul Mullowney 
295aa372e3fSPaul Mullowney       /* Create the solve analysis information */
296c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
297aa372e3fSPaul Mullowney 
298aa372e3fSPaul Mullowney       /* set the operation */
299aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
300aa372e3fSPaul Mullowney 
301aa372e3fSPaul Mullowney       /* set the matrix */
302aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
303aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = n;
304aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = n;
305aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = nzLower;
306aa372e3fSPaul Mullowney 
307aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
308aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1);
309aa372e3fSPaul Mullowney 
310aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower);
311aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower);
312aa372e3fSPaul Mullowney 
313aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(nzLower);
314aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+nzLower);
315aa372e3fSPaul Mullowney 
316aa372e3fSPaul Mullowney       /* perform the solve analysis */
317aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
318aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
319aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
320c41cb2e2SAlejandro Lamas Daviña                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
321aa372e3fSPaul Mullowney 
322aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
323aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
3242205254eSKarl Rupp 
325c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiLo);CHKERRCUDA(ierr);
326c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjLo);CHKERRCUDA(ierr);
327c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
3284863603aSSatish Balay       ierr = PetscLogCpuToGpu((n+1+nzLower)*sizeof(int)+nzLower*sizeof(PetscScalar));CHKERRQ(ierr);
3299ae82921SPaul Mullowney     } catch(char *ex) {
3309ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
3319ae82921SPaul Mullowney     }
3329ae82921SPaul Mullowney   }
3339ae82921SPaul Mullowney   PetscFunctionReturn(0);
3349ae82921SPaul Mullowney }
3359ae82921SPaul Mullowney 
336087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A)
3379ae82921SPaul Mullowney {
3389ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
3399ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
3409ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
341aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
3429ae82921SPaul Mullowney   cusparseStatus_t                  stat;
3439ae82921SPaul Mullowney   const PetscInt                    *aj = a->j,*adiag = a->diag,*vi;
3449ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
3459ae82921SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
3469ae82921SPaul Mullowney   PetscScalar                       *AAUp;
3479ae82921SPaul Mullowney   PetscInt                          i,nz, nzUpper, offset;
3489ae82921SPaul Mullowney   PetscErrorCode                    ierr;
3499ae82921SPaul Mullowney 
3509ae82921SPaul Mullowney   PetscFunctionBegin;
351b8ced49eSKarl Rupp   if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) {
3529ae82921SPaul Mullowney     try {
3539ae82921SPaul Mullowney       /* next, figure out the number of nonzeros in the upper triangular matrix. */
3549ae82921SPaul Mullowney       nzUpper = adiag[0]-adiag[n];
3559ae82921SPaul Mullowney 
3569ae82921SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
357c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
358c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
359c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
3609ae82921SPaul Mullowney 
3619ae82921SPaul Mullowney       /* Fill the upper triangular matrix */
3629ae82921SPaul Mullowney       AiUp[0]=(PetscInt) 0;
3639ae82921SPaul Mullowney       AiUp[n]=nzUpper;
3649ae82921SPaul Mullowney       offset = nzUpper;
3659ae82921SPaul Mullowney       for (i=n-1; i>=0; i--) {
3669ae82921SPaul Mullowney         v  = aa + adiag[i+1] + 1;
3679ae82921SPaul Mullowney         vi = aj + adiag[i+1] + 1;
3689ae82921SPaul Mullowney 
369e057df02SPaul Mullowney         /* number of elements NOT on the diagonal */
3709ae82921SPaul Mullowney         nz = adiag[i] - adiag[i+1]-1;
3719ae82921SPaul Mullowney 
372e057df02SPaul Mullowney         /* decrement the offset */
3739ae82921SPaul Mullowney         offset -= (nz+1);
3749ae82921SPaul Mullowney 
375e057df02SPaul Mullowney         /* first, set the diagonal elements */
3769ae82921SPaul Mullowney         AjUp[offset] = (PetscInt) i;
37709f51544SAlejandro Lamas Daviña         AAUp[offset] = (MatScalar)1./v[nz];
3789ae82921SPaul Mullowney         AiUp[i]      = AiUp[i+1] - (nz+1);
3799ae82921SPaul Mullowney 
380580bdb30SBarry Smith         ierr = PetscArraycpy(&(AjUp[offset+1]), vi, nz);CHKERRQ(ierr);
381580bdb30SBarry Smith         ierr = PetscArraycpy(&(AAUp[offset+1]), v, nz);CHKERRQ(ierr);
3829ae82921SPaul Mullowney       }
3832205254eSKarl Rupp 
384aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
385aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
3862205254eSKarl Rupp 
387aa372e3fSPaul Mullowney       /* Create the matrix description */
388c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
389c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
390c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
391c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
392c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
393aa372e3fSPaul Mullowney 
394aa372e3fSPaul Mullowney       /* Create the solve analysis information */
395c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
396aa372e3fSPaul Mullowney 
397aa372e3fSPaul Mullowney       /* set the operation */
398aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
399aa372e3fSPaul Mullowney 
400aa372e3fSPaul Mullowney       /* set the matrix */
401aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
402aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = n;
403aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = n;
404aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = nzUpper;
405aa372e3fSPaul Mullowney 
406aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
407aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1);
408aa372e3fSPaul Mullowney 
409aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper);
410aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper);
411aa372e3fSPaul Mullowney 
412aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper);
413aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper);
414aa372e3fSPaul Mullowney 
415aa372e3fSPaul Mullowney       /* perform the solve analysis */
416aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
417aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
418aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
419c41cb2e2SAlejandro Lamas Daviña                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
420aa372e3fSPaul Mullowney 
421aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
422aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
4232205254eSKarl Rupp 
424c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
425c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
426c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
4274863603aSSatish Balay       ierr = PetscLogCpuToGpu((n+1+nzUpper)*sizeof(int)+nzUpper*sizeof(PetscScalar));CHKERRQ(ierr);
4289ae82921SPaul Mullowney     } catch(char *ex) {
4299ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
4309ae82921SPaul Mullowney     }
4319ae82921SPaul Mullowney   }
4329ae82921SPaul Mullowney   PetscFunctionReturn(0);
4339ae82921SPaul Mullowney }
4349ae82921SPaul Mullowney 
435087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A)
4369ae82921SPaul Mullowney {
4379ae82921SPaul Mullowney   PetscErrorCode               ierr;
4389ae82921SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
4399ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
4409ae82921SPaul Mullowney   IS                           isrow = a->row,iscol = a->icol;
4419ae82921SPaul Mullowney   PetscBool                    row_identity,col_identity;
4429ae82921SPaul Mullowney   const PetscInt               *r,*c;
4439ae82921SPaul Mullowney   PetscInt                     n = A->rmap->n;
4449ae82921SPaul Mullowney 
4459ae82921SPaul Mullowney   PetscFunctionBegin;
446087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr);
447087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr);
4482205254eSKarl Rupp 
449e65717acSKarl Rupp   cusparseTriFactors->workVector = new THRUSTARRAY(n);
450aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=a->nz;
4519ae82921SPaul Mullowney 
452b8ced49eSKarl Rupp   A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH;
453e057df02SPaul Mullowney   /*lower triangular indices */
4549ae82921SPaul Mullowney   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
4559ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
4562205254eSKarl Rupp   if (!row_identity) {
457aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
458aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(r, r+n);
4592205254eSKarl Rupp   }
4609ae82921SPaul Mullowney   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
4619ae82921SPaul Mullowney 
462e057df02SPaul Mullowney   /*upper triangular indices */
4639ae82921SPaul Mullowney   ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
4649ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4652205254eSKarl Rupp   if (!col_identity) {
466aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
467aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices->assign(c, c+n);
4682205254eSKarl Rupp   }
4694863603aSSatish Balay 
4704863603aSSatish Balay   if(!row_identity && !col_identity) {
4714863603aSSatish Balay     ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr);
4724863603aSSatish Balay   } else if(!row_identity) {
4734863603aSSatish Balay     ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr);
4744863603aSSatish Balay   } else if(!col_identity) {
4754863603aSSatish Balay     ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr);
4764863603aSSatish Balay   }
4774863603aSSatish Balay 
4789ae82921SPaul Mullowney   ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
4799ae82921SPaul Mullowney   PetscFunctionReturn(0);
4809ae82921SPaul Mullowney }
4819ae82921SPaul Mullowney 
482087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A)
483087f3262SPaul Mullowney {
484087f3262SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
485087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
486aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
487aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
488087f3262SPaul Mullowney   cusparseStatus_t                  stat;
489087f3262SPaul Mullowney   PetscErrorCode                    ierr;
490087f3262SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
491087f3262SPaul Mullowney   PetscScalar                       *AAUp;
492087f3262SPaul Mullowney   PetscScalar                       *AALo;
493087f3262SPaul Mullowney   PetscInt                          nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j;
494087f3262SPaul Mullowney   Mat_SeqSBAIJ                      *b = (Mat_SeqSBAIJ*)A->data;
495087f3262SPaul Mullowney   const PetscInt                    *ai = b->i,*aj = b->j,*vj;
496087f3262SPaul Mullowney   const MatScalar                   *aa = b->a,*v;
497087f3262SPaul Mullowney 
498087f3262SPaul Mullowney   PetscFunctionBegin;
499b8ced49eSKarl Rupp   if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) {
500087f3262SPaul Mullowney     try {
501087f3262SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
502c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
503c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
504c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
505c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
506087f3262SPaul Mullowney 
507087f3262SPaul Mullowney       /* Fill the upper triangular matrix */
508087f3262SPaul Mullowney       AiUp[0]=(PetscInt) 0;
509087f3262SPaul Mullowney       AiUp[n]=nzUpper;
510087f3262SPaul Mullowney       offset = 0;
511087f3262SPaul Mullowney       for (i=0; i<n; i++) {
512087f3262SPaul Mullowney         /* set the pointers */
513087f3262SPaul Mullowney         v  = aa + ai[i];
514087f3262SPaul Mullowney         vj = aj + ai[i];
515087f3262SPaul Mullowney         nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */
516087f3262SPaul Mullowney 
517087f3262SPaul Mullowney         /* first, set the diagonal elements */
518087f3262SPaul Mullowney         AjUp[offset] = (PetscInt) i;
51909f51544SAlejandro Lamas Daviña         AAUp[offset] = (MatScalar)1.0/v[nz];
520087f3262SPaul Mullowney         AiUp[i]      = offset;
52109f51544SAlejandro Lamas Daviña         AALo[offset] = (MatScalar)1.0/v[nz];
522087f3262SPaul Mullowney 
523087f3262SPaul Mullowney         offset+=1;
524087f3262SPaul Mullowney         if (nz>0) {
525f22e0265SBarry Smith           ierr = PetscArraycpy(&(AjUp[offset]), vj, nz);CHKERRQ(ierr);
526580bdb30SBarry Smith           ierr = PetscArraycpy(&(AAUp[offset]), v, nz);CHKERRQ(ierr);
527087f3262SPaul Mullowney           for (j=offset; j<offset+nz; j++) {
528087f3262SPaul Mullowney             AAUp[j] = -AAUp[j];
529087f3262SPaul Mullowney             AALo[j] = AAUp[j]/v[nz];
530087f3262SPaul Mullowney           }
531087f3262SPaul Mullowney           offset+=nz;
532087f3262SPaul Mullowney         }
533087f3262SPaul Mullowney       }
534087f3262SPaul Mullowney 
535aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
536aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
537087f3262SPaul Mullowney 
538aa372e3fSPaul Mullowney       /* Create the matrix description */
539c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
540c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
541c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
542c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
543c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
544087f3262SPaul Mullowney 
545aa372e3fSPaul Mullowney       /* Create the solve analysis information */
546c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
547aa372e3fSPaul Mullowney 
548aa372e3fSPaul Mullowney       /* set the operation */
549aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
550aa372e3fSPaul Mullowney 
551aa372e3fSPaul Mullowney       /* set the matrix */
552aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
553aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = A->rmap->n;
554aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = A->cmap->n;
555aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = a->nz;
556aa372e3fSPaul Mullowney 
557aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
558aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
559aa372e3fSPaul Mullowney 
560aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
561aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
562aa372e3fSPaul Mullowney 
563aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
564aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz);
565aa372e3fSPaul Mullowney 
566aa372e3fSPaul Mullowney       /* perform the solve analysis */
567aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
568aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
569aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
570c41cb2e2SAlejandro Lamas Daviña                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
571aa372e3fSPaul Mullowney 
572aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
573aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
574aa372e3fSPaul Mullowney 
575aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
576aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
577aa372e3fSPaul Mullowney 
578aa372e3fSPaul Mullowney       /* Create the matrix description */
579c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
580c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
581c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
582c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
583c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
584aa372e3fSPaul Mullowney 
585aa372e3fSPaul Mullowney       /* Create the solve analysis information */
586c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
587aa372e3fSPaul Mullowney 
588aa372e3fSPaul Mullowney       /* set the operation */
589aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE;
590aa372e3fSPaul Mullowney 
591aa372e3fSPaul Mullowney       /* set the matrix */
592aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
593aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = A->rmap->n;
594aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = A->cmap->n;
595aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = a->nz;
596aa372e3fSPaul Mullowney 
597aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
598aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
599aa372e3fSPaul Mullowney 
600aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
601aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
602aa372e3fSPaul Mullowney 
603aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
604aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+a->nz);
6054863603aSSatish Balay       ierr = PetscLogCpuToGpu(2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)));CHKERRQ(ierr);
606aa372e3fSPaul Mullowney 
607aa372e3fSPaul Mullowney       /* perform the solve analysis */
608aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
609aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
610aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
611c41cb2e2SAlejandro Lamas Daviña                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
612aa372e3fSPaul Mullowney 
613aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
614aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
615087f3262SPaul Mullowney 
616b8ced49eSKarl Rupp       A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH;
617c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
618c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
619c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
620c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
621087f3262SPaul Mullowney     } catch(char *ex) {
622087f3262SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
623087f3262SPaul Mullowney     }
624087f3262SPaul Mullowney   }
625087f3262SPaul Mullowney   PetscFunctionReturn(0);
626087f3262SPaul Mullowney }
627087f3262SPaul Mullowney 
628087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A)
6299ae82921SPaul Mullowney {
6309ae82921SPaul Mullowney   PetscErrorCode               ierr;
631087f3262SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
632087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
633087f3262SPaul Mullowney   IS                           ip = a->row;
634087f3262SPaul Mullowney   const PetscInt               *rip;
635087f3262SPaul Mullowney   PetscBool                    perm_identity;
636087f3262SPaul Mullowney   PetscInt                     n = A->rmap->n;
637087f3262SPaul Mullowney 
638087f3262SPaul Mullowney   PetscFunctionBegin;
639087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr);
640e65717acSKarl Rupp   cusparseTriFactors->workVector = new THRUSTARRAY(n);
641aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=(a->nz-n)*2 + n;
642aa372e3fSPaul Mullowney 
643087f3262SPaul Mullowney   /*lower triangular indices */
644087f3262SPaul Mullowney   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
645087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
646087f3262SPaul Mullowney   if (!perm_identity) {
647aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
648aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(rip, rip+n);
649aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
650aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices->assign(rip, rip+n);
6514863603aSSatish Balay     ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr);
652087f3262SPaul Mullowney   }
653087f3262SPaul Mullowney   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
654087f3262SPaul Mullowney   PetscFunctionReturn(0);
655087f3262SPaul Mullowney }
656087f3262SPaul Mullowney 
6576fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
6589ae82921SPaul Mullowney {
6599ae82921SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
6609ae82921SPaul Mullowney   IS             isrow = b->row,iscol = b->col;
6619ae82921SPaul Mullowney   PetscBool      row_identity,col_identity;
662b175d8bbSPaul Mullowney   PetscErrorCode ierr;
6639ae82921SPaul Mullowney 
6649ae82921SPaul Mullowney   PetscFunctionBegin;
6659ae82921SPaul Mullowney   ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
666e057df02SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
6679ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
6689ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
669bda325fcSPaul Mullowney   if (row_identity && col_identity) {
670bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
671bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
672bda325fcSPaul Mullowney   } else {
673bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
674bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
675bda325fcSPaul Mullowney   }
6768dc1d2a3SPaul Mullowney 
677e057df02SPaul Mullowney   /* get the triangular factors */
678087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
6799ae82921SPaul Mullowney   PetscFunctionReturn(0);
6809ae82921SPaul Mullowney }
6819ae82921SPaul Mullowney 
682087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
683087f3262SPaul Mullowney {
684087f3262SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
685087f3262SPaul Mullowney   IS             ip = b->row;
686087f3262SPaul Mullowney   PetscBool      perm_identity;
687b175d8bbSPaul Mullowney   PetscErrorCode ierr;
688087f3262SPaul Mullowney 
689087f3262SPaul Mullowney   PetscFunctionBegin;
690087f3262SPaul Mullowney   ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
691087f3262SPaul Mullowney 
692087f3262SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
693087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
694087f3262SPaul Mullowney   if (perm_identity) {
695087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
696087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
697087f3262SPaul Mullowney   } else {
698087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
699087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
700087f3262SPaul Mullowney   }
701087f3262SPaul Mullowney 
702087f3262SPaul Mullowney   /* get the triangular factors */
703087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
704087f3262SPaul Mullowney   PetscFunctionReturn(0);
705087f3262SPaul Mullowney }
7069ae82921SPaul Mullowney 
707b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A)
708bda325fcSPaul Mullowney {
709bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
710aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
711aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
712aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
713aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
714bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
715aa372e3fSPaul Mullowney   cusparseIndexBase_t               indexBase;
716aa372e3fSPaul Mullowney   cusparseMatrixType_t              matrixType;
717aa372e3fSPaul Mullowney   cusparseFillMode_t                fillMode;
718aa372e3fSPaul Mullowney   cusparseDiagType_t                diagType;
719b175d8bbSPaul Mullowney 
720bda325fcSPaul Mullowney   PetscFunctionBegin;
721bda325fcSPaul Mullowney 
722aa372e3fSPaul Mullowney   /*********************************************/
723aa372e3fSPaul Mullowney   /* Now the Transpose of the Lower Tri Factor */
724aa372e3fSPaul Mullowney   /*********************************************/
725aa372e3fSPaul Mullowney 
726aa372e3fSPaul Mullowney   /* allocate space for the transpose of the lower triangular factor */
727aa372e3fSPaul Mullowney   loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
728aa372e3fSPaul Mullowney 
729aa372e3fSPaul Mullowney   /* set the matrix descriptors of the lower triangular factor */
730aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(loTriFactor->descr);
731aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(loTriFactor->descr);
732aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
733aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
734aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(loTriFactor->descr);
735aa372e3fSPaul Mullowney 
736aa372e3fSPaul Mullowney   /* Create the matrix description */
737c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat);
738c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat);
739c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat);
740c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat);
741c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat);
742aa372e3fSPaul Mullowney 
743aa372e3fSPaul Mullowney   /* Create the solve analysis information */
744c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat);
745aa372e3fSPaul Mullowney 
746aa372e3fSPaul Mullowney   /* set the operation */
747aa372e3fSPaul Mullowney   loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
748aa372e3fSPaul Mullowney 
749aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the lower triangular factor*/
750aa372e3fSPaul Mullowney   loTriFactorT->csrMat = new CsrMatrix;
751aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows;
752aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols;
753aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries;
754aa372e3fSPaul Mullowney   loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1);
755aa372e3fSPaul Mullowney   loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries);
756aa372e3fSPaul Mullowney   loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries);
757aa372e3fSPaul Mullowney 
758aa372e3fSPaul Mullowney   /* compute the transpose of the lower triangular factor, i.e. the CSC */
759aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows,
760aa372e3fSPaul Mullowney                           loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries,
761aa372e3fSPaul Mullowney                           loTriFactor->csrMat->values->data().get(),
762aa372e3fSPaul Mullowney                           loTriFactor->csrMat->row_offsets->data().get(),
763aa372e3fSPaul Mullowney                           loTriFactor->csrMat->column_indices->data().get(),
764aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->values->data().get(),
765aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->column_indices->data().get(),
766aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->row_offsets->data().get(),
767c41cb2e2SAlejandro Lamas Daviña                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
768aa372e3fSPaul Mullowney 
769aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
770aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp,
771aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries,
772aa372e3fSPaul Mullowney                            loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(),
773aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(),
774c41cb2e2SAlejandro Lamas Daviña                            loTriFactorT->solveInfo);CHKERRCUDA(stat);
775aa372e3fSPaul Mullowney 
776aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
777aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT;
778aa372e3fSPaul Mullowney 
779aa372e3fSPaul Mullowney   /*********************************************/
780aa372e3fSPaul Mullowney   /* Now the Transpose of the Upper Tri Factor */
781aa372e3fSPaul Mullowney   /*********************************************/
782aa372e3fSPaul Mullowney 
783aa372e3fSPaul Mullowney   /* allocate space for the transpose of the upper triangular factor */
784aa372e3fSPaul Mullowney   upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
785aa372e3fSPaul Mullowney 
786aa372e3fSPaul Mullowney   /* set the matrix descriptors of the upper triangular factor */
787aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(upTriFactor->descr);
788aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(upTriFactor->descr);
789aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
790aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
791aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(upTriFactor->descr);
792aa372e3fSPaul Mullowney 
793aa372e3fSPaul Mullowney   /* Create the matrix description */
794c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat);
795c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat);
796c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat);
797c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat);
798c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat);
799aa372e3fSPaul Mullowney 
800aa372e3fSPaul Mullowney   /* Create the solve analysis information */
801c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat);
802aa372e3fSPaul Mullowney 
803aa372e3fSPaul Mullowney   /* set the operation */
804aa372e3fSPaul Mullowney   upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
805aa372e3fSPaul Mullowney 
806aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the upper triangular factor*/
807aa372e3fSPaul Mullowney   upTriFactorT->csrMat = new CsrMatrix;
808aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows;
809aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols;
810aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries;
811aa372e3fSPaul Mullowney   upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1);
812aa372e3fSPaul Mullowney   upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries);
813aa372e3fSPaul Mullowney   upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries);
814aa372e3fSPaul Mullowney 
815aa372e3fSPaul Mullowney   /* compute the transpose of the upper triangular factor, i.e. the CSC */
816aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows,
817aa372e3fSPaul Mullowney                           upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries,
818aa372e3fSPaul Mullowney                           upTriFactor->csrMat->values->data().get(),
819aa372e3fSPaul Mullowney                           upTriFactor->csrMat->row_offsets->data().get(),
820aa372e3fSPaul Mullowney                           upTriFactor->csrMat->column_indices->data().get(),
821aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->values->data().get(),
822aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->column_indices->data().get(),
823aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->row_offsets->data().get(),
824c41cb2e2SAlejandro Lamas Daviña                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
825aa372e3fSPaul Mullowney 
826aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
827aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp,
828aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries,
829aa372e3fSPaul Mullowney                            upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(),
830aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(),
831c41cb2e2SAlejandro Lamas Daviña                            upTriFactorT->solveInfo);CHKERRCUDA(stat);
832aa372e3fSPaul Mullowney 
833aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
834aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT;
835bda325fcSPaul Mullowney   PetscFunctionReturn(0);
836bda325fcSPaul Mullowney }
837bda325fcSPaul Mullowney 
838b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A)
839bda325fcSPaul Mullowney {
840aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
841aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
842aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
843bda325fcSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
844bda325fcSPaul Mullowney   cusparseStatus_t             stat;
845aa372e3fSPaul Mullowney   cusparseIndexBase_t          indexBase;
846b06137fdSPaul Mullowney   cudaError_t                  err;
8474863603aSSatish Balay   PetscErrorCode               ierr;
848b175d8bbSPaul Mullowney 
849bda325fcSPaul Mullowney   PetscFunctionBegin;
850aa372e3fSPaul Mullowney 
851aa372e3fSPaul Mullowney   /* allocate space for the triangular factor information */
852aa372e3fSPaul Mullowney   matstructT = new Mat_SeqAIJCUSPARSEMultStruct;
853c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat);
854aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(matstruct->descr);
855c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat);
856c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
857aa372e3fSPaul Mullowney 
858b06137fdSPaul Mullowney   /* set alpha and beta */
859c41cb2e2SAlejandro Lamas Daviña   err = cudaMalloc((void **)&(matstructT->alpha),    sizeof(PetscScalar));CHKERRCUDA(err);
8607656d835SStefano Zampini   err = cudaMalloc((void **)&(matstructT->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err);
8617656d835SStefano Zampini   err = cudaMalloc((void **)&(matstructT->beta_one), sizeof(PetscScalar));CHKERRCUDA(err);
8627656d835SStefano Zampini   err = cudaMemcpy(matstructT->alpha,    &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
8637656d835SStefano Zampini   err = cudaMemcpy(matstructT->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
8647656d835SStefano Zampini   err = cudaMemcpy(matstructT->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
865c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
866b06137fdSPaul Mullowney 
867aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
868aa372e3fSPaul Mullowney     CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
869aa372e3fSPaul Mullowney     CsrMatrix *matrixT= new CsrMatrix;
870554b8892SKarl Rupp     matrixT->num_rows = A->cmap->n;
871554b8892SKarl Rupp     matrixT->num_cols = A->rmap->n;
872aa372e3fSPaul Mullowney     matrixT->num_entries = a->nz;
873*a8bd5306SMark Adams     matrixT->row_offsets = new THRUSTINTARRAY32(matrixT->num_rows+1);
874aa372e3fSPaul Mullowney     matrixT->column_indices = new THRUSTINTARRAY32(a->nz);
875aa372e3fSPaul Mullowney     matrixT->values = new THRUSTARRAY(a->nz);
876aa372e3fSPaul Mullowney 
877aa372e3fSPaul Mullowney     /* compute the transpose of the upper triangular factor, i.e. the CSC */
878aa372e3fSPaul Mullowney     indexBase = cusparseGetMatIndexBase(matstruct->descr);
879aa372e3fSPaul Mullowney     stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows,
880aa372e3fSPaul Mullowney                             matrix->num_cols, matrix->num_entries,
881aa372e3fSPaul Mullowney                             matrix->values->data().get(),
882aa372e3fSPaul Mullowney                             matrix->row_offsets->data().get(),
883aa372e3fSPaul Mullowney                             matrix->column_indices->data().get(),
884aa372e3fSPaul Mullowney                             matrixT->values->data().get(),
885aa372e3fSPaul Mullowney                             matrixT->column_indices->data().get(),
886aa372e3fSPaul Mullowney                             matrixT->row_offsets->data().get(),
887c41cb2e2SAlejandro Lamas Daviña                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
888aa372e3fSPaul Mullowney 
889aa372e3fSPaul Mullowney     /* assign the pointer */
890aa372e3fSPaul Mullowney     matstructT->mat = matrixT;
8914863603aSSatish Balay     ierr = PetscLogCpuToGpu(((A->rmap->n+1)+(a->nz))*sizeof(int)+(3+a->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
892aa372e3fSPaul Mullowney   } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
893aa372e3fSPaul Mullowney     /* First convert HYB to CSR */
894aa372e3fSPaul Mullowney     CsrMatrix *temp= new CsrMatrix;
895aa372e3fSPaul Mullowney     temp->num_rows = A->rmap->n;
896aa372e3fSPaul Mullowney     temp->num_cols = A->cmap->n;
897aa372e3fSPaul Mullowney     temp->num_entries = a->nz;
898aa372e3fSPaul Mullowney     temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
899aa372e3fSPaul Mullowney     temp->column_indices = new THRUSTINTARRAY32(a->nz);
900aa372e3fSPaul Mullowney     temp->values = new THRUSTARRAY(a->nz);
901aa372e3fSPaul Mullowney 
9022692e278SPaul Mullowney 
903aa372e3fSPaul Mullowney     stat = cusparse_hyb2csr(cusparsestruct->handle,
904aa372e3fSPaul Mullowney                             matstruct->descr, (cusparseHybMat_t)matstruct->mat,
905aa372e3fSPaul Mullowney                             temp->values->data().get(),
906aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
907c41cb2e2SAlejandro Lamas Daviña                             temp->column_indices->data().get());CHKERRCUDA(stat);
908aa372e3fSPaul Mullowney 
909aa372e3fSPaul Mullowney     /* Next, convert CSR to CSC (i.e. the matrix transpose) */
910aa372e3fSPaul Mullowney     CsrMatrix *tempT= new CsrMatrix;
911aa372e3fSPaul Mullowney     tempT->num_rows = A->rmap->n;
912aa372e3fSPaul Mullowney     tempT->num_cols = A->cmap->n;
913aa372e3fSPaul Mullowney     tempT->num_entries = a->nz;
914aa372e3fSPaul Mullowney     tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
915aa372e3fSPaul Mullowney     tempT->column_indices = new THRUSTINTARRAY32(a->nz);
916aa372e3fSPaul Mullowney     tempT->values = new THRUSTARRAY(a->nz);
917aa372e3fSPaul Mullowney 
918aa372e3fSPaul Mullowney     stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows,
919aa372e3fSPaul Mullowney                             temp->num_cols, temp->num_entries,
920aa372e3fSPaul Mullowney                             temp->values->data().get(),
921aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
922aa372e3fSPaul Mullowney                             temp->column_indices->data().get(),
923aa372e3fSPaul Mullowney                             tempT->values->data().get(),
924aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
925aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
926c41cb2e2SAlejandro Lamas Daviña                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
927aa372e3fSPaul Mullowney 
928aa372e3fSPaul Mullowney     /* Last, convert CSC to HYB */
929aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat;
930c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
931aa372e3fSPaul Mullowney     cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
932aa372e3fSPaul Mullowney       CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
933aa372e3fSPaul Mullowney     stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n,
934aa372e3fSPaul Mullowney                             matstructT->descr, tempT->values->data().get(),
935aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
936aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
937c41cb2e2SAlejandro Lamas Daviña                             hybMat, 0, partition);CHKERRCUDA(stat);
938aa372e3fSPaul Mullowney 
939aa372e3fSPaul Mullowney     /* assign the pointer */
940aa372e3fSPaul Mullowney     matstructT->mat = hybMat;
9414863603aSSatish Balay     ierr = PetscLogCpuToGpu((2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)))+3*sizeof(PetscScalar));CHKERRQ(ierr);
942aa372e3fSPaul Mullowney 
943aa372e3fSPaul Mullowney     /* delete temporaries */
944aa372e3fSPaul Mullowney     if (tempT) {
945aa372e3fSPaul Mullowney       if (tempT->values) delete (THRUSTARRAY*) tempT->values;
946aa372e3fSPaul Mullowney       if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices;
947aa372e3fSPaul Mullowney       if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets;
948aa372e3fSPaul Mullowney       delete (CsrMatrix*) tempT;
949087f3262SPaul Mullowney     }
950aa372e3fSPaul Mullowney     if (temp) {
951aa372e3fSPaul Mullowney       if (temp->values) delete (THRUSTARRAY*) temp->values;
952aa372e3fSPaul Mullowney       if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices;
953aa372e3fSPaul Mullowney       if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets;
954aa372e3fSPaul Mullowney       delete (CsrMatrix*) temp;
955aa372e3fSPaul Mullowney     }
956aa372e3fSPaul Mullowney   }
957aa372e3fSPaul Mullowney   /* assign the compressed row indices */
958aa372e3fSPaul Mullowney   matstructT->cprowIndices = new THRUSTINTARRAY;
959554b8892SKarl Rupp   matstructT->cprowIndices->resize(A->cmap->n);
960554b8892SKarl Rupp   thrust::sequence(matstructT->cprowIndices->begin(), matstructT->cprowIndices->end());
961aa372e3fSPaul Mullowney   /* assign the pointer */
962aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT;
963bda325fcSPaul Mullowney   PetscFunctionReturn(0);
964bda325fcSPaul Mullowney }
965bda325fcSPaul Mullowney 
9666fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
967bda325fcSPaul Mullowney {
968c41cb2e2SAlejandro Lamas Daviña   PetscInt                              n = xx->map->n;
969465f34aeSAlejandro Lamas Daviña   const PetscScalar                     *barray;
970465f34aeSAlejandro Lamas Daviña   PetscScalar                           *xarray;
971465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<const PetscScalar> bGPU;
972465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<PetscScalar>       xGPU;
973bda325fcSPaul Mullowney   cusparseStatus_t                      stat;
974bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
975aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
976aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
977aa372e3fSPaul Mullowney   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
978b175d8bbSPaul Mullowney   PetscErrorCode                        ierr;
979bda325fcSPaul Mullowney 
980bda325fcSPaul Mullowney   PetscFunctionBegin;
981aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
982aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
983bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
984aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
985aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
986bda325fcSPaul Mullowney   }
987bda325fcSPaul Mullowney 
988bda325fcSPaul Mullowney   /* Get the GPU pointers */
989c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
990c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
991c41cb2e2SAlejandro Lamas Daviña   xGPU = thrust::device_pointer_cast(xarray);
992c41cb2e2SAlejandro Lamas Daviña   bGPU = thrust::device_pointer_cast(barray);
993bda325fcSPaul Mullowney 
9947a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
995aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
996c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
997c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()),
998c41cb2e2SAlejandro Lamas Daviña                xGPU);
999aa372e3fSPaul Mullowney 
1000aa372e3fSPaul Mullowney   /* First, solve U */
1001aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
10027656d835SStefano Zampini                         upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr,
1003aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1004aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1005aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1006aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1007c41cb2e2SAlejandro Lamas Daviña                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1008aa372e3fSPaul Mullowney 
1009aa372e3fSPaul Mullowney   /* Then, solve L */
1010aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
10117656d835SStefano Zampini                         loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr,
1012aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1013aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1014aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1015aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1016c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1017aa372e3fSPaul Mullowney 
1018aa372e3fSPaul Mullowney   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1019c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()),
1020c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()),
1021aa372e3fSPaul Mullowney                tempGPU->begin());
1022aa372e3fSPaul Mullowney 
1023aa372e3fSPaul Mullowney   /* Copy the temporary to the full solution. */
1024c41cb2e2SAlejandro Lamas Daviña   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU);
1025bda325fcSPaul Mullowney 
1026bda325fcSPaul Mullowney   /* restore */
1027c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1028c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1029c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1030661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1031958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1032bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1033bda325fcSPaul Mullowney }
1034bda325fcSPaul Mullowney 
10356fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
1036bda325fcSPaul Mullowney {
1037465f34aeSAlejandro Lamas Daviña   const PetscScalar                 *barray;
1038465f34aeSAlejandro Lamas Daviña   PetscScalar                       *xarray;
1039bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
1040bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1041aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1042aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1043aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1044b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
1045bda325fcSPaul Mullowney 
1046bda325fcSPaul Mullowney   PetscFunctionBegin;
1047aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
1048aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
1049bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1050aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1051aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1052bda325fcSPaul Mullowney   }
1053bda325fcSPaul Mullowney 
1054bda325fcSPaul Mullowney   /* Get the GPU pointers */
1055c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1056c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1057bda325fcSPaul Mullowney 
10587a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1059aa372e3fSPaul Mullowney   /* First, solve U */
1060aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
10617656d835SStefano Zampini                         upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr,
1062aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1063aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1064aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1065aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1066c41cb2e2SAlejandro Lamas Daviña                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1067aa372e3fSPaul Mullowney 
1068aa372e3fSPaul Mullowney   /* Then, solve L */
1069aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
10707656d835SStefano Zampini                         loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr,
1071aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1072aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1073aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1074aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1075c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1076bda325fcSPaul Mullowney 
1077bda325fcSPaul Mullowney   /* restore */
1078c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1079c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1080c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1081661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1082958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1083bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1084bda325fcSPaul Mullowney }
1085bda325fcSPaul Mullowney 
10866fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
10879ae82921SPaul Mullowney {
1088465f34aeSAlejandro Lamas Daviña   const PetscScalar                     *barray;
1089465f34aeSAlejandro Lamas Daviña   PetscScalar                           *xarray;
1090465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<const PetscScalar> bGPU;
1091465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<PetscScalar>       xGPU;
10929ae82921SPaul Mullowney   cusparseStatus_t                      stat;
10939ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1094aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1095aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1096aa372e3fSPaul Mullowney   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1097b175d8bbSPaul Mullowney   PetscErrorCode                        ierr;
10989ae82921SPaul Mullowney 
10999ae82921SPaul Mullowney   PetscFunctionBegin;
1100ebc8f436SDominic Meiser 
1101e057df02SPaul Mullowney   /* Get the GPU pointers */
1102c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1103c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1104c41cb2e2SAlejandro Lamas Daviña   xGPU = thrust::device_pointer_cast(xarray);
1105c41cb2e2SAlejandro Lamas Daviña   bGPU = thrust::device_pointer_cast(barray);
11069ae82921SPaul Mullowney 
11077a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1108aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
1109c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
1110c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()),
1111c41cb2e2SAlejandro Lamas Daviña                xGPU);
1112aa372e3fSPaul Mullowney 
1113aa372e3fSPaul Mullowney   /* Next, solve L */
1114aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
11157656d835SStefano Zampini                         loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr,
1116aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1117aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1118aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1119aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
1120c41cb2e2SAlejandro Lamas Daviña                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1121aa372e3fSPaul Mullowney 
1122aa372e3fSPaul Mullowney   /* Then, solve U */
1123aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
11247656d835SStefano Zampini                         upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr,
1125aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1126aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1127aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1128aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
1129c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1130aa372e3fSPaul Mullowney 
1131aa372e3fSPaul Mullowney   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1132c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()),
1133c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->end()),
1134aa372e3fSPaul Mullowney                tempGPU->begin());
1135aa372e3fSPaul Mullowney 
1136aa372e3fSPaul Mullowney   /* Copy the temporary to the full solution. */
1137c41cb2e2SAlejandro Lamas Daviña   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU);
11389ae82921SPaul Mullowney 
1139c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1140c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1141c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1142661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1143958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
11449ae82921SPaul Mullowney   PetscFunctionReturn(0);
11459ae82921SPaul Mullowney }
11469ae82921SPaul Mullowney 
11476fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
11489ae82921SPaul Mullowney {
1149465f34aeSAlejandro Lamas Daviña   const PetscScalar                 *barray;
1150465f34aeSAlejandro Lamas Daviña   PetscScalar                       *xarray;
11519ae82921SPaul Mullowney   cusparseStatus_t                  stat;
11529ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1153aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1154aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1155aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1156b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
11579ae82921SPaul Mullowney 
11589ae82921SPaul Mullowney   PetscFunctionBegin;
1159e057df02SPaul Mullowney   /* Get the GPU pointers */
1160c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1161c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
11629ae82921SPaul Mullowney 
11637a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1164aa372e3fSPaul Mullowney   /* First, solve L */
1165aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
11667656d835SStefano Zampini                         loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr,
1167aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1168aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1169aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1170aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
1171c41cb2e2SAlejandro Lamas Daviña                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1172aa372e3fSPaul Mullowney 
1173aa372e3fSPaul Mullowney   /* Next, solve U */
1174aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
11757656d835SStefano Zampini                         upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr,
1176aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1177aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1178aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1179aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
1180c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
11819ae82921SPaul Mullowney 
1182c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1183c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1184c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1185661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1186958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
11879ae82921SPaul Mullowney   PetscFunctionReturn(0);
11889ae82921SPaul Mullowney }
11899ae82921SPaul Mullowney 
11906fa9248bSJed Brown static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A)
11919ae82921SPaul Mullowney {
1192aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1193aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
11949ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
11959ae82921SPaul Mullowney   PetscInt                     m = A->rmap->n,*ii,*ridx;
11969ae82921SPaul Mullowney   PetscErrorCode               ierr;
1197aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1198b06137fdSPaul Mullowney   cudaError_t                  err;
11999ae82921SPaul Mullowney 
12009ae82921SPaul Mullowney   PetscFunctionBegin;
1201fdc842d1SBarry Smith   if (A->pinnedtocpu) PetscFunctionReturn(0);
1202b8ced49eSKarl Rupp   if (A->valid_GPU_matrix == PETSC_OFFLOAD_UNALLOCATED || A->valid_GPU_matrix == PETSC_OFFLOAD_CPU) {
12039ae82921SPaul Mullowney     ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
120434d6c7a5SJose E. Roman     if (A->assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) {
120534d6c7a5SJose E. Roman       CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
120634d6c7a5SJose E. Roman       /* copy values only */
120734d6c7a5SJose E. Roman       matrix->values->assign(a->a, a->a+a->nz);
12084863603aSSatish Balay       ierr = PetscLogCpuToGpu((a->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
120934d6c7a5SJose E. Roman     } else {
1210470880abSPatrick Sanan       MatSeqAIJCUSPARSEMultStruct_Destroy(&matstruct,cusparsestruct->format);
12119ae82921SPaul Mullowney       try {
1212aa372e3fSPaul Mullowney         cusparsestruct->nonzerorow=0;
1213aa372e3fSPaul Mullowney         for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0);
12149ae82921SPaul Mullowney 
12159ae82921SPaul Mullowney         if (a->compressedrow.use) {
12169ae82921SPaul Mullowney           m    = a->compressedrow.nrows;
12179ae82921SPaul Mullowney           ii   = a->compressedrow.i;
12189ae82921SPaul Mullowney           ridx = a->compressedrow.rindex;
12199ae82921SPaul Mullowney         } else {
1220b06137fdSPaul Mullowney           /* Forcing compressed row on the GPU */
12219ae82921SPaul Mullowney           int k=0;
1222854ce69bSBarry Smith           ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr);
1223854ce69bSBarry Smith           ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr);
12249ae82921SPaul Mullowney           ii[0]=0;
12259ae82921SPaul Mullowney           for (int j = 0; j<m; j++) {
12269ae82921SPaul Mullowney             if ((a->i[j+1]-a->i[j])>0) {
12279ae82921SPaul Mullowney               ii[k]  = a->i[j];
12289ae82921SPaul Mullowney               ridx[k]= j;
12299ae82921SPaul Mullowney               k++;
12309ae82921SPaul Mullowney             }
12319ae82921SPaul Mullowney           }
1232aa372e3fSPaul Mullowney           ii[cusparsestruct->nonzerorow] = a->nz;
1233aa372e3fSPaul Mullowney           m = cusparsestruct->nonzerorow;
12349ae82921SPaul Mullowney         }
12359ae82921SPaul Mullowney 
1236aa372e3fSPaul Mullowney         /* allocate space for the triangular factor information */
1237aa372e3fSPaul Mullowney         matstruct = new Mat_SeqAIJCUSPARSEMultStruct;
1238c41cb2e2SAlejandro Lamas Daviña         stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat);
1239c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
1240c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
12419ae82921SPaul Mullowney 
1242c41cb2e2SAlejandro Lamas Daviña         err = cudaMalloc((void **)&(matstruct->alpha),    sizeof(PetscScalar));CHKERRCUDA(err);
12437656d835SStefano Zampini         err = cudaMalloc((void **)&(matstruct->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err);
12447656d835SStefano Zampini         err = cudaMalloc((void **)&(matstruct->beta_one), sizeof(PetscScalar));CHKERRCUDA(err);
12457656d835SStefano Zampini         err = cudaMemcpy(matstruct->alpha,    &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
12467656d835SStefano Zampini         err = cudaMemcpy(matstruct->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
12477656d835SStefano Zampini         err = cudaMemcpy(matstruct->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
1248c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
1249b06137fdSPaul Mullowney 
1250aa372e3fSPaul Mullowney         /* Build a hybrid/ellpack matrix if this option is chosen for the storage */
1251aa372e3fSPaul Mullowney         if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1252aa372e3fSPaul Mullowney           /* set the matrix */
1253aa372e3fSPaul Mullowney           CsrMatrix *matrix= new CsrMatrix;
1254a65300a6SPaul Mullowney           matrix->num_rows = m;
1255aa372e3fSPaul Mullowney           matrix->num_cols = A->cmap->n;
1256aa372e3fSPaul Mullowney           matrix->num_entries = a->nz;
1257a65300a6SPaul Mullowney           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1258a65300a6SPaul Mullowney           matrix->row_offsets->assign(ii, ii + m+1);
12599ae82921SPaul Mullowney 
1260aa372e3fSPaul Mullowney           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1261aa372e3fSPaul Mullowney           matrix->column_indices->assign(a->j, a->j+a->nz);
1262aa372e3fSPaul Mullowney 
1263aa372e3fSPaul Mullowney           matrix->values = new THRUSTARRAY(a->nz);
1264aa372e3fSPaul Mullowney           matrix->values->assign(a->a, a->a+a->nz);
1265aa372e3fSPaul Mullowney 
1266aa372e3fSPaul Mullowney           /* assign the pointer */
1267aa372e3fSPaul Mullowney           matstruct->mat = matrix;
1268aa372e3fSPaul Mullowney 
1269aa372e3fSPaul Mullowney         } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
1270aa372e3fSPaul Mullowney           CsrMatrix *matrix= new CsrMatrix;
1271a65300a6SPaul Mullowney           matrix->num_rows = m;
1272aa372e3fSPaul Mullowney           matrix->num_cols = A->cmap->n;
1273aa372e3fSPaul Mullowney           matrix->num_entries = a->nz;
1274a65300a6SPaul Mullowney           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1275a65300a6SPaul Mullowney           matrix->row_offsets->assign(ii, ii + m+1);
1276aa372e3fSPaul Mullowney 
1277aa372e3fSPaul Mullowney           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1278aa372e3fSPaul Mullowney           matrix->column_indices->assign(a->j, a->j+a->nz);
1279aa372e3fSPaul Mullowney 
1280aa372e3fSPaul Mullowney           matrix->values = new THRUSTARRAY(a->nz);
1281aa372e3fSPaul Mullowney           matrix->values->assign(a->a, a->a+a->nz);
1282aa372e3fSPaul Mullowney 
1283aa372e3fSPaul Mullowney           cusparseHybMat_t hybMat;
1284c41cb2e2SAlejandro Lamas Daviña           stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
1285aa372e3fSPaul Mullowney           cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
1286aa372e3fSPaul Mullowney             CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
1287a65300a6SPaul Mullowney           stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols,
1288aa372e3fSPaul Mullowney               matstruct->descr, matrix->values->data().get(),
1289aa372e3fSPaul Mullowney               matrix->row_offsets->data().get(),
1290aa372e3fSPaul Mullowney               matrix->column_indices->data().get(),
1291c41cb2e2SAlejandro Lamas Daviña               hybMat, 0, partition);CHKERRCUDA(stat);
1292aa372e3fSPaul Mullowney           /* assign the pointer */
1293aa372e3fSPaul Mullowney           matstruct->mat = hybMat;
1294aa372e3fSPaul Mullowney 
1295aa372e3fSPaul Mullowney           if (matrix) {
1296aa372e3fSPaul Mullowney             if (matrix->values) delete (THRUSTARRAY*)matrix->values;
1297aa372e3fSPaul Mullowney             if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices;
1298aa372e3fSPaul Mullowney             if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets;
1299aa372e3fSPaul Mullowney             delete (CsrMatrix*)matrix;
1300087f3262SPaul Mullowney           }
1301087f3262SPaul Mullowney         }
1302ca45077fSPaul Mullowney 
1303aa372e3fSPaul Mullowney         /* assign the compressed row indices */
1304aa372e3fSPaul Mullowney         matstruct->cprowIndices = new THRUSTINTARRAY(m);
1305aa372e3fSPaul Mullowney         matstruct->cprowIndices->assign(ridx,ridx+m);
13064863603aSSatish Balay         ierr = PetscLogCpuToGpu(((m+1)+(a->nz))*sizeof(int)+m*sizeof(PetscInt)+(3+(a->nz))*sizeof(PetscScalar));CHKERRQ(ierr);
1307aa372e3fSPaul Mullowney 
1308aa372e3fSPaul Mullowney         /* assign the pointer */
1309aa372e3fSPaul Mullowney         cusparsestruct->mat = matstruct;
1310aa372e3fSPaul Mullowney 
13119ae82921SPaul Mullowney         if (!a->compressedrow.use) {
13129ae82921SPaul Mullowney           ierr = PetscFree(ii);CHKERRQ(ierr);
13139ae82921SPaul Mullowney           ierr = PetscFree(ridx);CHKERRQ(ierr);
13149ae82921SPaul Mullowney         }
1315e65717acSKarl Rupp         cusparsestruct->workVector = new THRUSTARRAY(m);
13169ae82921SPaul Mullowney       } catch(char *ex) {
13179ae82921SPaul Mullowney         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
13189ae82921SPaul Mullowney       }
131934d6c7a5SJose E. Roman       cusparsestruct->nonzerostate = A->nonzerostate;
132034d6c7a5SJose E. Roman     }
13217d0e52d8SJose E. Roman     ierr = WaitForGPU();CHKERRCUDA(ierr);
1322b8ced49eSKarl Rupp     A->valid_GPU_matrix = PETSC_OFFLOAD_BOTH;
13239ae82921SPaul Mullowney     ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
13249ae82921SPaul Mullowney   }
13259ae82921SPaul Mullowney   PetscFunctionReturn(0);
13269ae82921SPaul Mullowney }
13279ae82921SPaul Mullowney 
1328c41cb2e2SAlejandro Lamas Daviña struct VecCUDAPlusEquals
1329aa372e3fSPaul Mullowney {
1330aa372e3fSPaul Mullowney   template <typename Tuple>
1331aa372e3fSPaul Mullowney   __host__ __device__
1332aa372e3fSPaul Mullowney   void operator()(Tuple t)
1333aa372e3fSPaul Mullowney   {
1334aa372e3fSPaul Mullowney     thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t);
1335aa372e3fSPaul Mullowney   }
1336aa372e3fSPaul Mullowney };
1337aa372e3fSPaul Mullowney 
13386fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
13399ae82921SPaul Mullowney {
1340b175d8bbSPaul Mullowney   PetscErrorCode ierr;
13419ae82921SPaul Mullowney 
13429ae82921SPaul Mullowney   PetscFunctionBegin;
13437656d835SStefano Zampini   ierr = MatMultAdd_SeqAIJCUSPARSE(A,xx,NULL,yy);CHKERRQ(ierr);
13449ae82921SPaul Mullowney   PetscFunctionReturn(0);
13459ae82921SPaul Mullowney }
13469ae82921SPaul Mullowney 
13476fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
1348ca45077fSPaul Mullowney {
1349ca45077fSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1350aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
13519ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct *matstructT;
1352465f34aeSAlejandro Lamas Daviña   const PetscScalar            *xarray;
1353465f34aeSAlejandro Lamas Daviña   PetscScalar                  *yarray;
1354b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1355aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1356ca45077fSPaul Mullowney 
1357ca45077fSPaul Mullowney   PetscFunctionBegin;
135834d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
135934d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
13609ff858a8SKarl Rupp   matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1361aa372e3fSPaul Mullowney   if (!matstructT) {
1362bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1363aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1364bda325fcSPaul Mullowney   }
1365c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1366ca6ae6e6SAlejandro Lamas Daviña   ierr = VecSet(yy,0);CHKERRQ(ierr);
1367c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1368aa372e3fSPaul Mullowney 
13697a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1370aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1371aa372e3fSPaul Mullowney     CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1372aa372e3fSPaul Mullowney     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1373aa372e3fSPaul Mullowney                              mat->num_rows, mat->num_cols,
1374b06137fdSPaul Mullowney                              mat->num_entries, matstructT->alpha, matstructT->descr,
1375aa372e3fSPaul Mullowney                              mat->values->data().get(), mat->row_offsets->data().get(),
13767656d835SStefano Zampini                              mat->column_indices->data().get(), xarray, matstructT->beta_zero,
1377c41cb2e2SAlejandro Lamas Daviña                              yarray);CHKERRCUDA(stat);
1378aa372e3fSPaul Mullowney   } else {
1379aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1380aa372e3fSPaul Mullowney     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1381b06137fdSPaul Mullowney                              matstructT->alpha, matstructT->descr, hybMat,
13827656d835SStefano Zampini                              xarray, matstructT->beta_zero,
1383c41cb2e2SAlejandro Lamas Daviña                              yarray);CHKERRCUDA(stat);
1384ca45077fSPaul Mullowney   }
1385c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1386c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1387aa372e3fSPaul Mullowney   if (!cusparsestruct->stream) {
1388c41cb2e2SAlejandro Lamas Daviña     ierr = WaitForGPU();CHKERRCUDA(ierr);
1389ca45077fSPaul Mullowney   }
1390661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1391958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
1392ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1393ca45077fSPaul Mullowney }
1394ca45077fSPaul Mullowney 
1395aa372e3fSPaul Mullowney 
13966fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
13979ae82921SPaul Mullowney {
13989ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1399aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
14009ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct *matstruct;
1401465f34aeSAlejandro Lamas Daviña   const PetscScalar            *xarray;
14027656d835SStefano Zampini   PetscScalar                  *zarray,*dptr,*beta;
1403b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1404aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1405c1fb3f03SStefano Zampini   PetscBool                    cmpr; /* if the matrix has been compressed (zero rows) */
14066e111a19SKarl Rupp 
14079ae82921SPaul Mullowney   PetscFunctionBegin;
140834d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
140934d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
14109ff858a8SKarl Rupp   matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
14119ae82921SPaul Mullowney   try {
1412c1fb3f03SStefano Zampini     cmpr = (PetscBool)(cusparsestruct->workVector->size() == (thrust::detail::vector_base<PetscScalar, thrust::device_malloc_allocator<PetscScalar> >::size_type)(A->rmap->n));
1413c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1414c1fb3f03SStefano Zampini     if (yy && !cmpr) { /* MatMultAdd with noncompressed storage -> need uptodate zz vector */
1415f2d70e9dSBarry Smith       ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr);
1416c1fb3f03SStefano Zampini     } else {
1417c1fb3f03SStefano Zampini       ierr = VecCUDAGetArrayWrite(zz,&zarray);CHKERRQ(ierr);
1418c1fb3f03SStefano Zampini     }
1419c1fb3f03SStefano Zampini     dptr = cmpr ? zarray : cusparsestruct->workVector->data().get();
14207656d835SStefano Zampini     beta = (yy == zz && dptr == zarray) ? matstruct->beta_one : matstruct->beta_zero;
14219ae82921SPaul Mullowney 
14227a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
14237656d835SStefano Zampini     /* csr_spmv is multiply add */
1424aa372e3fSPaul Mullowney     if (cusparsestruct->format == MAT_CUSPARSE_CSR) {
1425b06137fdSPaul Mullowney       /* here we need to be careful to set the number of rows in the multiply to the
1426b06137fdSPaul Mullowney          number of compressed rows in the matrix ... which is equivalent to the
1427b06137fdSPaul Mullowney          size of the workVector */
14287656d835SStefano Zampini       CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1429aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1430a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1431b06137fdSPaul Mullowney                                mat->num_entries, matstruct->alpha, matstruct->descr,
1432aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
14337656d835SStefano Zampini                                mat->column_indices->data().get(), xarray, beta,
14347656d835SStefano Zampini                                dptr);CHKERRCUDA(stat);
1435aa372e3fSPaul Mullowney     } else {
1436a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1437301298b4SMark Adams 	cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1438aa372e3fSPaul Mullowney         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1439b06137fdSPaul Mullowney                                  matstruct->alpha, matstruct->descr, hybMat,
14407656d835SStefano Zampini                                  xarray, beta,
14417656d835SStefano Zampini                                  dptr);CHKERRCUDA(stat);
1442a65300a6SPaul Mullowney       }
1443aa372e3fSPaul Mullowney     }
1444958c4211Shannah_mairs     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1445aa372e3fSPaul Mullowney 
14467656d835SStefano Zampini     if (yy) {
14477656d835SStefano Zampini       if (dptr != zarray) {
14487656d835SStefano Zampini         ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
14497656d835SStefano Zampini       } else if (zz != yy) {
14507656d835SStefano Zampini         ierr = VecAXPY_SeqCUDA(zz,1.0,yy);CHKERRQ(ierr);
14517656d835SStefano Zampini       }
14527656d835SStefano Zampini     } else if (dptr != zarray) {
1453c1fb3f03SStefano Zampini       ierr = VecSet_SeqCUDA(zz,0);CHKERRQ(ierr);
14547656d835SStefano Zampini     }
1455aa372e3fSPaul Mullowney     /* scatter the data from the temporary into the full vector with a += operation */
14567a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
14577656d835SStefano Zampini     if (dptr != zarray) {
14587656d835SStefano Zampini       thrust::device_ptr<PetscScalar> zptr;
14597656d835SStefano Zampini 
14607656d835SStefano Zampini       zptr = thrust::device_pointer_cast(zarray);
1461c41cb2e2SAlejandro Lamas Daviña       thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))),
1462c41cb2e2SAlejandro Lamas Daviña                        thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1463c41cb2e2SAlejandro Lamas Daviña                        VecCUDAPlusEquals());
14647656d835SStefano Zampini     }
1465958c4211Shannah_mairs     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1466c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1467c1fb3f03SStefano Zampini     if (yy && !cmpr) {
1468f2d70e9dSBarry Smith       ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr);
1469c1fb3f03SStefano Zampini     } else {
1470c1fb3f03SStefano Zampini       ierr = VecCUDARestoreArrayWrite(zz,&zarray);CHKERRQ(ierr);
1471c1fb3f03SStefano Zampini     }
14729ae82921SPaul Mullowney   } catch(char *ex) {
14739ae82921SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
14749ae82921SPaul Mullowney   }
14757656d835SStefano Zampini   if (!yy) { /* MatMult */
14767656d835SStefano Zampini     if (!cusparsestruct->stream) {
1477c41cb2e2SAlejandro Lamas Daviña       ierr = WaitForGPU();CHKERRCUDA(ierr);
14787656d835SStefano Zampini     }
14797656d835SStefano Zampini   }
1480958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz);CHKERRQ(ierr);
14819ae82921SPaul Mullowney   PetscFunctionReturn(0);
14829ae82921SPaul Mullowney }
14839ae82921SPaul Mullowney 
14846fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
1485ca45077fSPaul Mullowney {
1486ca45077fSPaul Mullowney   Mat_SeqAIJ                      *a = (Mat_SeqAIJ*)A->data;
1487aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE              *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
14889ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct    *matstructT;
1489c41cb2e2SAlejandro Lamas Daviña   thrust::device_ptr<PetscScalar> zptr;
1490465f34aeSAlejandro Lamas Daviña   const PetscScalar               *xarray;
1491465f34aeSAlejandro Lamas Daviña   PetscScalar                     *zarray;
1492b175d8bbSPaul Mullowney   PetscErrorCode                  ierr;
1493aa372e3fSPaul Mullowney   cusparseStatus_t                stat;
14946e111a19SKarl Rupp 
1495ca45077fSPaul Mullowney   PetscFunctionBegin;
149634d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
149734d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
14989ff858a8SKarl Rupp   matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1499aa372e3fSPaul Mullowney   if (!matstructT) {
1500bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1501aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1502661c2d29Shannah_mairs     ierr = WaitForGPU();CHKERRCUDA(ierr);
1503bda325fcSPaul Mullowney   }
1504aa372e3fSPaul Mullowney 
1505ca45077fSPaul Mullowney   try {
1506c41cb2e2SAlejandro Lamas Daviña     ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
1507c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1508f2d70e9dSBarry Smith     ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr);
1509c41cb2e2SAlejandro Lamas Daviña     zptr = thrust::device_pointer_cast(zarray);
1510ca45077fSPaul Mullowney 
15117a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1512e057df02SPaul Mullowney     /* multiply add with matrix transpose */
1513aa372e3fSPaul Mullowney     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1514aa372e3fSPaul Mullowney       CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1515b06137fdSPaul Mullowney       /* here we need to be careful to set the number of rows in the multiply to the
1516b06137fdSPaul Mullowney          number of compressed rows in the matrix ... which is equivalent to the
1517b06137fdSPaul Mullowney          size of the workVector */
1518aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1519a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1520b06137fdSPaul Mullowney                                mat->num_entries, matstructT->alpha, matstructT->descr,
1521aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
15227656d835SStefano Zampini                                mat->column_indices->data().get(), xarray, matstructT->beta_zero,
1523c41cb2e2SAlejandro Lamas Daviña                                cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1524aa372e3fSPaul Mullowney     } else {
1525aa372e3fSPaul Mullowney       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1526a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1527aa372e3fSPaul Mullowney         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1528b06137fdSPaul Mullowney             matstructT->alpha, matstructT->descr, hybMat,
15297656d835SStefano Zampini             xarray, matstructT->beta_zero,
1530c41cb2e2SAlejandro Lamas Daviña             cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1531a65300a6SPaul Mullowney       }
1532aa372e3fSPaul Mullowney     }
1533aa372e3fSPaul Mullowney 
1534aa372e3fSPaul Mullowney     /* scatter the data from the temporary into the full vector with a += operation */
1535c41cb2e2SAlejandro Lamas Daviña     thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))),
1536554b8892SKarl Rupp         thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))) + A->cmap->n,
1537c41cb2e2SAlejandro Lamas Daviña         VecCUDAPlusEquals());
1538c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1539f2d70e9dSBarry Smith     ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr);
1540661c2d29Shannah_mairs     ierr = WaitForGPU();CHKERRCUDA(ierr);
1541661c2d29Shannah_mairs     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1542ca45077fSPaul Mullowney 
1543ca45077fSPaul Mullowney   } catch(char *ex) {
1544ca45077fSPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1545c41cb2e2SAlejandro Lamas Daviña     ierr = WaitForGPU();CHKERRCUDA(ierr);
1546661c2d29Shannah_mairs   }
1547958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz);CHKERRQ(ierr);
1548ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1549ca45077fSPaul Mullowney }
1550ca45077fSPaul Mullowney 
15516fa9248bSJed Brown static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode)
15529ae82921SPaul Mullowney {
15539ae82921SPaul Mullowney   PetscErrorCode ierr;
15546e111a19SKarl Rupp 
15559ae82921SPaul Mullowney   PetscFunctionBegin;
15569ae82921SPaul Mullowney   ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr);
1557489de41dSStefano Zampini   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
1558bc3f50f2SPaul Mullowney   if (A->factortype == MAT_FACTOR_NONE) {
1559e057df02SPaul Mullowney     ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1560bc3f50f2SPaul Mullowney   }
1561bbf3fe20SPaul Mullowney   A->ops->mult             = MatMult_SeqAIJCUSPARSE;
1562bbf3fe20SPaul Mullowney   A->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1563bbf3fe20SPaul Mullowney   A->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1564bbf3fe20SPaul Mullowney   A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
15659ae82921SPaul Mullowney   PetscFunctionReturn(0);
15669ae82921SPaul Mullowney }
15679ae82921SPaul Mullowney 
15689ae82921SPaul Mullowney /* --------------------------------------------------------------------------------*/
1569e057df02SPaul Mullowney /*@
15709ae82921SPaul Mullowney    MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format
1571e057df02SPaul Mullowney    (the default parallel PETSc format). This matrix will ultimately pushed down
1572e057df02SPaul Mullowney    to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix
1573e057df02SPaul Mullowney    assembly performance the user should preallocate the matrix storage by setting
1574e057df02SPaul Mullowney    the parameter nz (or the array nnz).  By setting these parameters accurately,
1575e057df02SPaul Mullowney    performance during matrix assembly can be increased by more than a factor of 50.
15769ae82921SPaul Mullowney 
1577d083f849SBarry Smith    Collective
15789ae82921SPaul Mullowney 
15799ae82921SPaul Mullowney    Input Parameters:
15809ae82921SPaul Mullowney +  comm - MPI communicator, set to PETSC_COMM_SELF
15819ae82921SPaul Mullowney .  m - number of rows
15829ae82921SPaul Mullowney .  n - number of columns
15839ae82921SPaul Mullowney .  nz - number of nonzeros per row (same for all rows)
15849ae82921SPaul Mullowney -  nnz - array containing the number of nonzeros in the various rows
15850298fd71SBarry Smith          (possibly different for each row) or NULL
15869ae82921SPaul Mullowney 
15879ae82921SPaul Mullowney    Output Parameter:
15889ae82921SPaul Mullowney .  A - the matrix
15899ae82921SPaul Mullowney 
15909ae82921SPaul Mullowney    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
15919ae82921SPaul Mullowney    MatXXXXSetPreallocation() paradgm instead of this routine directly.
15929ae82921SPaul Mullowney    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
15939ae82921SPaul Mullowney 
15949ae82921SPaul Mullowney    Notes:
15959ae82921SPaul Mullowney    If nnz is given then nz is ignored
15969ae82921SPaul Mullowney 
15979ae82921SPaul Mullowney    The AIJ format (also called the Yale sparse matrix format or
15989ae82921SPaul Mullowney    compressed row storage), is fully compatible with standard Fortran 77
15999ae82921SPaul Mullowney    storage.  That is, the stored row and column indices can begin at
16009ae82921SPaul Mullowney    either one (as in Fortran) or zero.  See the users' manual for details.
16019ae82921SPaul Mullowney 
16029ae82921SPaul Mullowney    Specify the preallocated storage with either nz or nnz (not both).
16030298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
16049ae82921SPaul Mullowney    allocation.  For large problems you MUST preallocate memory or you
16059ae82921SPaul Mullowney    will get TERRIBLE performance, see the users' manual chapter on matrices.
16069ae82921SPaul Mullowney 
16079ae82921SPaul Mullowney    By default, this format uses inodes (identical nodes) when possible, to
16089ae82921SPaul Mullowney    improve numerical efficiency of matrix-vector products and solves. We
16099ae82921SPaul Mullowney    search for consecutive rows with the same nonzero structure, thereby
16109ae82921SPaul Mullowney    reusing matrix information to achieve increased efficiency.
16119ae82921SPaul Mullowney 
16129ae82921SPaul Mullowney    Level: intermediate
16139ae82921SPaul Mullowney 
1614e057df02SPaul Mullowney .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE
16159ae82921SPaul Mullowney @*/
16169ae82921SPaul Mullowney PetscErrorCode  MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
16179ae82921SPaul Mullowney {
16189ae82921SPaul Mullowney   PetscErrorCode ierr;
16199ae82921SPaul Mullowney 
16209ae82921SPaul Mullowney   PetscFunctionBegin;
16219ae82921SPaul Mullowney   ierr = MatCreate(comm,A);CHKERRQ(ierr);
16229ae82921SPaul Mullowney   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
16239ae82921SPaul Mullowney   ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
16249ae82921SPaul Mullowney   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
16259ae82921SPaul Mullowney   PetscFunctionReturn(0);
16269ae82921SPaul Mullowney }
16279ae82921SPaul Mullowney 
16286fa9248bSJed Brown static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A)
16299ae82921SPaul Mullowney {
16309ae82921SPaul Mullowney   PetscErrorCode   ierr;
1631ab25e6cbSDominic Meiser 
16329ae82921SPaul Mullowney   PetscFunctionBegin;
16339ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
1634b8ced49eSKarl Rupp     if (A->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1635470880abSPatrick Sanan       ierr = MatSeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr);
16369ae82921SPaul Mullowney     }
16379ae82921SPaul Mullowney   } else {
1638470880abSPatrick Sanan     ierr = MatSeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr);
1639aa372e3fSPaul Mullowney   }
16409ae82921SPaul Mullowney   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
16419ae82921SPaul Mullowney   PetscFunctionReturn(0);
16429ae82921SPaul Mullowney }
16439ae82921SPaul Mullowney 
16449ff858a8SKarl Rupp static PetscErrorCode MatDuplicate_SeqAIJCUSPARSE(Mat A,MatDuplicateOption cpvalues,Mat *B)
16459ff858a8SKarl Rupp {
16469ff858a8SKarl Rupp   PetscErrorCode ierr;
16479ff858a8SKarl Rupp   Mat C;
16489ff858a8SKarl Rupp   cusparseStatus_t stat;
16499ff858a8SKarl Rupp   cusparseHandle_t handle=0;
16509ff858a8SKarl Rupp 
16519ff858a8SKarl Rupp   PetscFunctionBegin;
16529ff858a8SKarl Rupp   ierr = MatDuplicate_SeqAIJ(A,cpvalues,B);CHKERRQ(ierr);
16539ff858a8SKarl Rupp   C    = *B;
165434136279SStefano Zampini   ierr = PetscFree(C->defaultvectype);CHKERRQ(ierr);
165534136279SStefano Zampini   ierr = PetscStrallocpy(VECCUDA,&C->defaultvectype);CHKERRQ(ierr);
165634136279SStefano Zampini 
16579ff858a8SKarl Rupp   /* inject CUSPARSE-specific stuff */
16589ff858a8SKarl Rupp   if (C->factortype==MAT_FACTOR_NONE) {
16599ff858a8SKarl Rupp     /* you cannot check the inode.use flag here since the matrix was just created.
16609ff858a8SKarl Rupp        now build a GPU matrix data structure */
16619ff858a8SKarl Rupp     C->spptr = new Mat_SeqAIJCUSPARSE;
16629ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->mat          = 0;
16639ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->matTranspose = 0;
16649ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->workVector   = 0;
16659ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->format       = MAT_CUSPARSE_CSR;
16669ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream       = 0;
16679ff858a8SKarl Rupp     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
16689ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle       = handle;
16699ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->nonzerostate = 0;
16709ff858a8SKarl Rupp   } else {
16719ff858a8SKarl Rupp     /* NEXT, set the pointers to the triangular factors */
16729ff858a8SKarl Rupp     C->spptr = new Mat_SeqAIJCUSPARSETriFactors;
16739ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtr          = 0;
16749ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtr          = 0;
16759ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtrTranspose = 0;
16769ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtrTranspose = 0;
16779ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->rpermIndices            = 0;
16789ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->cpermIndices            = 0;
16799ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->workVector              = 0;
16809ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle                  = 0;
16819ff858a8SKarl Rupp     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
16829ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle                  = handle;
16839ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->nnz                     = 0;
16849ff858a8SKarl Rupp   }
16859ff858a8SKarl Rupp 
16869ff858a8SKarl Rupp   C->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
16879ff858a8SKarl Rupp   C->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
16889ff858a8SKarl Rupp   C->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
16899ff858a8SKarl Rupp   C->ops->mult             = MatMult_SeqAIJCUSPARSE;
16909ff858a8SKarl Rupp   C->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
16919ff858a8SKarl Rupp   C->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
16929ff858a8SKarl Rupp   C->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
16939ff858a8SKarl Rupp   C->ops->duplicate        = MatDuplicate_SeqAIJCUSPARSE;
16949ff858a8SKarl Rupp 
16959ff858a8SKarl Rupp   ierr = PetscObjectChangeTypeName((PetscObject)C,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
16969ff858a8SKarl Rupp 
1697b8ced49eSKarl Rupp   C->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED;
16989ff858a8SKarl Rupp 
16999ff858a8SKarl Rupp   ierr = PetscObjectComposeFunction((PetscObject)C, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
17009ff858a8SKarl Rupp   PetscFunctionReturn(0);
17019ff858a8SKarl Rupp }
17029ff858a8SKarl Rupp 
170302fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat B)
17049ae82921SPaul Mullowney {
17059ae82921SPaul Mullowney   PetscErrorCode ierr;
1706aa372e3fSPaul Mullowney   cusparseStatus_t stat;
1707aa372e3fSPaul Mullowney   cusparseHandle_t handle=0;
17089ae82921SPaul Mullowney 
17099ae82921SPaul Mullowney   PetscFunctionBegin;
171034136279SStefano Zampini   ierr = PetscFree(B->defaultvectype);CHKERRQ(ierr);
171134136279SStefano Zampini   ierr = PetscStrallocpy(VECCUDA,&B->defaultvectype);CHKERRQ(ierr);
171234136279SStefano Zampini 
17139ae82921SPaul Mullowney   if (B->factortype==MAT_FACTOR_NONE) {
1714e057df02SPaul Mullowney     /* you cannot check the inode.use flag here since the matrix was just created.
1715e057df02SPaul Mullowney        now build a GPU matrix data structure */
17169ae82921SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSE;
17179ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat          = 0;
1718aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0;
1719aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector   = 0;
1720e057df02SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->format       = MAT_CUSPARSE_CSR;
1721aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream       = 0;
1722c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1723aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle       = handle;
17249ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)B->spptr)->nonzerostate = 0;
17259ae82921SPaul Mullowney   } else {
17269ae82921SPaul Mullowney     /* NEXT, set the pointers to the triangular factors */
1727debe9ee2SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSETriFactors;
17289ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr          = 0;
17299ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr          = 0;
1730aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0;
1731aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0;
1732aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices            = 0;
1733aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices            = 0;
1734aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector              = 0;
1735c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1736aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = handle;
1737aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz                     = 0;
17389ae82921SPaul Mullowney   }
1739aa372e3fSPaul Mullowney 
17409ae82921SPaul Mullowney   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
17419ae82921SPaul Mullowney   B->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
17429ae82921SPaul Mullowney   B->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
1743ca45077fSPaul Mullowney   B->ops->mult             = MatMult_SeqAIJCUSPARSE;
1744ca45077fSPaul Mullowney   B->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1745ca45077fSPaul Mullowney   B->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1746ca45077fSPaul Mullowney   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
17479ff858a8SKarl Rupp   B->ops->duplicate        = MatDuplicate_SeqAIJCUSPARSE;
17482205254eSKarl Rupp 
17499ae82921SPaul Mullowney   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
17502205254eSKarl Rupp 
1751b8ced49eSKarl Rupp   B->valid_GPU_matrix = PETSC_OFFLOAD_UNALLOCATED;
17522205254eSKarl Rupp 
1753bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
17549ae82921SPaul Mullowney   PetscFunctionReturn(0);
17559ae82921SPaul Mullowney }
17569ae82921SPaul Mullowney 
175702fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B)
175802fe1965SBarry Smith {
175902fe1965SBarry Smith   PetscErrorCode ierr;
176002fe1965SBarry Smith 
176102fe1965SBarry Smith   PetscFunctionBegin;
176202fe1965SBarry Smith   ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr);
1763489de41dSStefano Zampini   ierr = MatConvert_SeqAIJ_SeqAIJCUSPARSE(B);CHKERRQ(ierr);
176402fe1965SBarry Smith   PetscFunctionReturn(0);
176502fe1965SBarry Smith }
176602fe1965SBarry Smith 
17673ca39a21SBarry Smith /*MC
1768e057df02SPaul Mullowney    MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices.
1769e057df02SPaul Mullowney 
1770e057df02SPaul Mullowney    A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either
17712692e278SPaul Mullowney    CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later.
17722692e278SPaul Mullowney    All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library.
1773e057df02SPaul Mullowney 
1774e057df02SPaul Mullowney    Options Database Keys:
1775e057df02SPaul Mullowney +  -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions()
1776aa372e3fSPaul 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).
1777a2b725a8SWilliam Gropp -  -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).
1778e057df02SPaul Mullowney 
1779e057df02SPaul Mullowney   Level: beginner
1780e057df02SPaul Mullowney 
17818468deeeSKarl Rupp .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
1782e057df02SPaul Mullowney M*/
17837f756511SDominic Meiser 
178442c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*);
178542c9c57cSBarry Smith 
17860f39cd5aSBarry Smith 
17873ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_CUSPARSE(void)
178842c9c57cSBarry Smith {
178942c9c57cSBarry Smith   PetscErrorCode ierr;
179042c9c57cSBarry Smith 
179142c9c57cSBarry Smith   PetscFunctionBegin;
17923ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
17933ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
17943ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
17953ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
179642c9c57cSBarry Smith   PetscFunctionReturn(0);
179742c9c57cSBarry Smith }
179829b38603SBarry Smith 
179981e08676SBarry Smith 
1800470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct)
18017f756511SDominic Meiser {
18027f756511SDominic Meiser   cusparseStatus_t stat;
18037f756511SDominic Meiser   cusparseHandle_t handle;
18047f756511SDominic Meiser 
18057f756511SDominic Meiser   PetscFunctionBegin;
18067f756511SDominic Meiser   if (*cusparsestruct) {
1807470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format);
1808470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format);
18097f756511SDominic Meiser     delete (*cusparsestruct)->workVector;
18107f756511SDominic Meiser     if (handle = (*cusparsestruct)->handle) {
1811c41cb2e2SAlejandro Lamas Daviña       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
18127f756511SDominic Meiser     }
18137f756511SDominic Meiser     delete *cusparsestruct;
18147f756511SDominic Meiser     *cusparsestruct = 0;
18157f756511SDominic Meiser   }
18167f756511SDominic Meiser   PetscFunctionReturn(0);
18177f756511SDominic Meiser }
18187f756511SDominic Meiser 
18197f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat)
18207f756511SDominic Meiser {
18217f756511SDominic Meiser   PetscFunctionBegin;
18227f756511SDominic Meiser   if (*mat) {
18237f756511SDominic Meiser     delete (*mat)->values;
18247f756511SDominic Meiser     delete (*mat)->column_indices;
18257f756511SDominic Meiser     delete (*mat)->row_offsets;
18267f756511SDominic Meiser     delete *mat;
18277f756511SDominic Meiser     *mat = 0;
18287f756511SDominic Meiser   }
18297f756511SDominic Meiser   PetscFunctionReturn(0);
18307f756511SDominic Meiser }
18317f756511SDominic Meiser 
1832470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor)
18337f756511SDominic Meiser {
18347f756511SDominic Meiser   cusparseStatus_t stat;
18357f756511SDominic Meiser   PetscErrorCode   ierr;
18367f756511SDominic Meiser 
18377f756511SDominic Meiser   PetscFunctionBegin;
18387f756511SDominic Meiser   if (*trifactor) {
1839c41cb2e2SAlejandro Lamas Daviña     if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); }
1840c41cb2e2SAlejandro Lamas Daviña     if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); }
18417f756511SDominic Meiser     ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr);
18427f756511SDominic Meiser     delete *trifactor;
18437f756511SDominic Meiser     *trifactor = 0;
18447f756511SDominic Meiser   }
18457f756511SDominic Meiser   PetscFunctionReturn(0);
18467f756511SDominic Meiser }
18477f756511SDominic Meiser 
1848470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format)
18497f756511SDominic Meiser {
18507f756511SDominic Meiser   CsrMatrix        *mat;
18517f756511SDominic Meiser   cusparseStatus_t stat;
18527f756511SDominic Meiser   cudaError_t      err;
18537f756511SDominic Meiser 
18547f756511SDominic Meiser   PetscFunctionBegin;
18557f756511SDominic Meiser   if (*matstruct) {
18567f756511SDominic Meiser     if ((*matstruct)->mat) {
18577f756511SDominic Meiser       if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) {
18587f756511SDominic Meiser         cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat;
1859c41cb2e2SAlejandro Lamas Daviña         stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat);
18607f756511SDominic Meiser       } else {
18617f756511SDominic Meiser         mat = (CsrMatrix*)(*matstruct)->mat;
18627f756511SDominic Meiser         CsrMatrix_Destroy(&mat);
18637f756511SDominic Meiser       }
18647f756511SDominic Meiser     }
1865c41cb2e2SAlejandro Lamas Daviña     if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); }
18667f756511SDominic Meiser     delete (*matstruct)->cprowIndices;
1867c41cb2e2SAlejandro Lamas Daviña     if ((*matstruct)->alpha)     { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); }
18687656d835SStefano Zampini     if ((*matstruct)->beta_zero) { err=cudaFree((*matstruct)->beta_zero);CHKERRCUDA(err); }
18697656d835SStefano Zampini     if ((*matstruct)->beta_one)  { err=cudaFree((*matstruct)->beta_one);CHKERRCUDA(err); }
18707f756511SDominic Meiser     delete *matstruct;
18717f756511SDominic Meiser     *matstruct = 0;
18727f756511SDominic Meiser   }
18737f756511SDominic Meiser   PetscFunctionReturn(0);
18747f756511SDominic Meiser }
18757f756511SDominic Meiser 
1876470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors)
18777f756511SDominic Meiser {
18787f756511SDominic Meiser   cusparseHandle_t handle;
18797f756511SDominic Meiser   cusparseStatus_t stat;
18807f756511SDominic Meiser 
18817f756511SDominic Meiser   PetscFunctionBegin;
18827f756511SDominic Meiser   if (*trifactors) {
1883470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtr);
1884470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtr);
1885470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose);
1886470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose);
18877f756511SDominic Meiser     delete (*trifactors)->rpermIndices;
18887f756511SDominic Meiser     delete (*trifactors)->cpermIndices;
18897f756511SDominic Meiser     delete (*trifactors)->workVector;
18907f756511SDominic Meiser     if (handle = (*trifactors)->handle) {
1891c41cb2e2SAlejandro Lamas Daviña       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
18927f756511SDominic Meiser     }
18937f756511SDominic Meiser     delete *trifactors;
18947f756511SDominic Meiser     *trifactors = 0;
18957f756511SDominic Meiser   }
18967f756511SDominic Meiser   PetscFunctionReturn(0);
18977f756511SDominic Meiser }
18987f756511SDominic Meiser 
1899