xref: /petsc/src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu (revision 20291eb5db195094b15008b4c3f060d35e3ffbed)
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
653800007SKarl Rupp #define PETSC_SKIP_CXX_COMPLEX_FIX
799acd6aaSStefano Zampini #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1
89ae82921SPaul Mullowney 
93d13b8fdSMatthew G. Knepley #include <petscconf.h>
103d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/aij.h>          /*I "petscmat.h" I*/
11087f3262SPaul Mullowney #include <../src/mat/impls/sbaij/seq/sbaij.h>
123d13b8fdSMatthew G. Knepley #include <../src/vec/vec/impls/dvecimpl.h>
13af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
149ae82921SPaul Mullowney #undef VecType
153d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h>
16bc3f50f2SPaul Mullowney 
17e057df02SPaul Mullowney const char *const MatCUSPARSEStorageFormats[] = {"CSR","ELL","HYB","MatCUSPARSEStorageFormat","MAT_CUSPARSE_",0};
189ae82921SPaul Mullowney 
19087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
20087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
21087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
22087f3262SPaul Mullowney 
236fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
246fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
256fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
26087f3262SPaul Mullowney 
276fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat,Vec,Vec);
286fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
296fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
306fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
314416b707SBarry Smith static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat);
326fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat,Vec,Vec);
336fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
346fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
356fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
369ae82921SPaul Mullowney 
377f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix**);
38470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct**);
39470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct**,MatCUSPARSEStorageFormat);
40470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors**);
41470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE**);
427f756511SDominic Meiser 
43b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetStream(Mat A,const cudaStream_t stream)
44b06137fdSPaul Mullowney {
45b06137fdSPaul Mullowney   cusparseStatus_t   stat;
46b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
47b06137fdSPaul Mullowney 
48b06137fdSPaul Mullowney   PetscFunctionBegin;
49b06137fdSPaul Mullowney   cusparsestruct->stream = stream;
50c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetStream(cusparsestruct->handle,cusparsestruct->stream);CHKERRCUDA(stat);
51b06137fdSPaul Mullowney   PetscFunctionReturn(0);
52b06137fdSPaul Mullowney }
53b06137fdSPaul Mullowney 
54b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetHandle(Mat A,const cusparseHandle_t handle)
55b06137fdSPaul Mullowney {
56b06137fdSPaul Mullowney   cusparseStatus_t   stat;
57b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
58b06137fdSPaul Mullowney 
59b06137fdSPaul Mullowney   PetscFunctionBegin;
606b1cf21dSAlejandro Lamas Daviña   if (cusparsestruct->handle != handle) {
6116a2e217SAlejandro Lamas Daviña     if (cusparsestruct->handle) {
62c41cb2e2SAlejandro Lamas Daviña       stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUDA(stat);
6316a2e217SAlejandro Lamas Daviña     }
64b06137fdSPaul Mullowney     cusparsestruct->handle = handle;
656b1cf21dSAlejandro Lamas Daviña   }
66c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
67b06137fdSPaul Mullowney   PetscFunctionReturn(0);
68b06137fdSPaul Mullowney }
69b06137fdSPaul Mullowney 
70b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSEClearHandle(Mat A)
71b06137fdSPaul Mullowney {
72b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
73b06137fdSPaul Mullowney   PetscFunctionBegin;
74b06137fdSPaul Mullowney   if (cusparsestruct->handle)
75b06137fdSPaul Mullowney     cusparsestruct->handle = 0;
76b06137fdSPaul Mullowney   PetscFunctionReturn(0);
77b06137fdSPaul Mullowney }
78b06137fdSPaul Mullowney 
79ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_cusparse(Mat A,MatSolverType *type)
809ae82921SPaul Mullowney {
819ae82921SPaul Mullowney   PetscFunctionBegin;
829ae82921SPaul Mullowney   *type = MATSOLVERCUSPARSE;
839ae82921SPaul Mullowney   PetscFunctionReturn(0);
849ae82921SPaul Mullowney }
859ae82921SPaul Mullowney 
86c708e6cdSJed Brown /*MC
87087f3262SPaul Mullowney   MATSOLVERCUSPARSE = "cusparse" - A matrix type providing triangular solvers for seq matrices
88087f3262SPaul Mullowney   on a single GPU of type, seqaijcusparse, aijcusparse, or seqaijcusp, aijcusp. Currently supported
89087f3262SPaul Mullowney   algorithms are ILU(k) and ICC(k). Typically, deeper factorizations (larger k) results in poorer
90087f3262SPaul Mullowney   performance in the triangular solves. Full LU, and Cholesky decompositions can be solved through the
91087f3262SPaul Mullowney   CUSPARSE triangular solve algorithm. However, the performance can be quite poor and thus these
92087f3262SPaul Mullowney   algorithms are not recommended. This class does NOT support direct solver operations.
93c708e6cdSJed Brown 
949ae82921SPaul Mullowney   Level: beginner
95c708e6cdSJed Brown 
963ca39a21SBarry Smith .seealso: PCFactorSetMatSolverType(), MatSolverType, MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
97c708e6cdSJed Brown M*/
989ae82921SPaul Mullowney 
9942c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat A,MatFactorType ftype,Mat *B)
1009ae82921SPaul Mullowney {
1019ae82921SPaul Mullowney   PetscErrorCode ierr;
102bc3f50f2SPaul Mullowney   PetscInt       n = A->rmap->n;
1039ae82921SPaul Mullowney 
1049ae82921SPaul Mullowney   PetscFunctionBegin;
105bc3f50f2SPaul Mullowney   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
106404133a2SPaul Mullowney   (*B)->factortype = ftype;
107bc3f50f2SPaul Mullowney   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
1089ae82921SPaul Mullowney   ierr = MatSetType(*B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
1092205254eSKarl Rupp 
110087f3262SPaul Mullowney   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
11133d57670SJed Brown     ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr);
1129ae82921SPaul Mullowney     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJCUSPARSE;
1139ae82921SPaul Mullowney     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJCUSPARSE;
114087f3262SPaul Mullowney   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
115087f3262SPaul Mullowney     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJCUSPARSE;
116087f3262SPaul Mullowney     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJCUSPARSE;
1179ae82921SPaul Mullowney   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported for CUSPARSE Matrix Types");
118bc3f50f2SPaul Mullowney 
119fa03d054SJed Brown   ierr = MatSeqAIJSetPreallocation(*B,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
1203ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)(*B),"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cusparse);CHKERRQ(ierr);
1219ae82921SPaul Mullowney   PetscFunctionReturn(0);
1229ae82921SPaul Mullowney }
1239ae82921SPaul Mullowney 
124bc3f50f2SPaul Mullowney PETSC_INTERN PetscErrorCode MatCUSPARSESetFormat_SeqAIJCUSPARSE(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
125ca45077fSPaul Mullowney {
126aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1276e111a19SKarl Rupp 
128ca45077fSPaul Mullowney   PetscFunctionBegin;
129ca45077fSPaul Mullowney   switch (op) {
130e057df02SPaul Mullowney   case MAT_CUSPARSE_MULT:
131aa372e3fSPaul Mullowney     cusparsestruct->format = format;
132ca45077fSPaul Mullowney     break;
133e057df02SPaul Mullowney   case MAT_CUSPARSE_ALL:
134aa372e3fSPaul Mullowney     cusparsestruct->format = format;
135ca45077fSPaul Mullowney     break;
136ca45077fSPaul Mullowney   default:
13736d62e41SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op);
138ca45077fSPaul Mullowney   }
139ca45077fSPaul Mullowney   PetscFunctionReturn(0);
140ca45077fSPaul Mullowney }
1419ae82921SPaul Mullowney 
142e057df02SPaul Mullowney /*@
143e057df02SPaul Mullowney    MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular
144e057df02SPaul Mullowney    operation. Only the MatMult operation can use different GPU storage formats
145aa372e3fSPaul Mullowney    for MPIAIJCUSPARSE matrices.
146e057df02SPaul Mullowney    Not Collective
147e057df02SPaul Mullowney 
148e057df02SPaul Mullowney    Input Parameters:
1498468deeeSKarl Rupp +  A - Matrix of type SEQAIJCUSPARSE
15036d62e41SPaul 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.
1512692e278SPaul Mullowney -  format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2)
152e057df02SPaul Mullowney 
153e057df02SPaul Mullowney    Output Parameter:
154e057df02SPaul Mullowney 
155e057df02SPaul Mullowney    Level: intermediate
156e057df02SPaul Mullowney 
1578468deeeSKarl Rupp .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
158e057df02SPaul Mullowney @*/
159e057df02SPaul Mullowney PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
160e057df02SPaul Mullowney {
161e057df02SPaul Mullowney   PetscErrorCode ierr;
1626e111a19SKarl Rupp 
163e057df02SPaul Mullowney   PetscFunctionBegin;
164e057df02SPaul Mullowney   PetscValidHeaderSpecific(A, MAT_CLASSID,1);
165e057df02SPaul Mullowney   ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr);
166e057df02SPaul Mullowney   PetscFunctionReturn(0);
167e057df02SPaul Mullowney }
168e057df02SPaul Mullowney 
1694416b707SBarry Smith static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat A)
1709ae82921SPaul Mullowney {
1719ae82921SPaul Mullowney   PetscErrorCode           ierr;
172e057df02SPaul Mullowney   MatCUSPARSEStorageFormat format;
1739ae82921SPaul Mullowney   PetscBool                flg;
174a183c035SDominic Meiser   Mat_SeqAIJCUSPARSE       *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1756e111a19SKarl Rupp 
1769ae82921SPaul Mullowney   PetscFunctionBegin;
177e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSPARSE options");CHKERRQ(ierr);
1789ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
179e057df02SPaul Mullowney     ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV",
180a183c035SDominic Meiser                             "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
181e057df02SPaul Mullowney     if (flg) {
182e057df02SPaul Mullowney       ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr);
183045c96e1SPaul Mullowney     }
1849ae82921SPaul Mullowney   }
1854c87dfd4SPaul Mullowney   ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve",
186a183c035SDominic Meiser                           "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
1874c87dfd4SPaul Mullowney   if (flg) {
1884c87dfd4SPaul Mullowney     ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr);
1894c87dfd4SPaul Mullowney   }
1900af67c1bSStefano Zampini   ierr = PetscOptionsTail();CHKERRQ(ierr);
1919ae82921SPaul Mullowney   PetscFunctionReturn(0);
1929ae82921SPaul Mullowney 
1939ae82921SPaul Mullowney }
1949ae82921SPaul Mullowney 
1956fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1969ae82921SPaul Mullowney {
1979ae82921SPaul Mullowney   PetscErrorCode ierr;
1989ae82921SPaul Mullowney 
1999ae82921SPaul Mullowney   PetscFunctionBegin;
2009ae82921SPaul Mullowney   ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
2019ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2029ae82921SPaul Mullowney   PetscFunctionReturn(0);
2039ae82921SPaul Mullowney }
2049ae82921SPaul Mullowney 
2056fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2069ae82921SPaul Mullowney {
2079ae82921SPaul Mullowney   PetscErrorCode ierr;
2089ae82921SPaul Mullowney 
2099ae82921SPaul Mullowney   PetscFunctionBegin;
2109ae82921SPaul Mullowney   ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
2119ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2129ae82921SPaul Mullowney   PetscFunctionReturn(0);
2139ae82921SPaul Mullowney }
2149ae82921SPaul Mullowney 
215087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
216087f3262SPaul Mullowney {
217087f3262SPaul Mullowney   PetscErrorCode ierr;
218087f3262SPaul Mullowney 
219087f3262SPaul Mullowney   PetscFunctionBegin;
220087f3262SPaul Mullowney   ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
221087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
222087f3262SPaul Mullowney   PetscFunctionReturn(0);
223087f3262SPaul Mullowney }
224087f3262SPaul Mullowney 
225087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
226087f3262SPaul Mullowney {
227087f3262SPaul Mullowney   PetscErrorCode ierr;
228087f3262SPaul Mullowney 
229087f3262SPaul Mullowney   PetscFunctionBegin;
230087f3262SPaul Mullowney   ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
231087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
232087f3262SPaul Mullowney   PetscFunctionReturn(0);
233087f3262SPaul Mullowney }
234087f3262SPaul Mullowney 
235087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A)
2369ae82921SPaul Mullowney {
2379ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
2389ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
2399ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
240aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
2419ae82921SPaul Mullowney   cusparseStatus_t                  stat;
2429ae82921SPaul Mullowney   const PetscInt                    *ai = a->i,*aj = a->j,*vi;
2439ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
2449ae82921SPaul Mullowney   PetscInt                          *AiLo, *AjLo;
2459ae82921SPaul Mullowney   PetscScalar                       *AALo;
2469ae82921SPaul Mullowney   PetscInt                          i,nz, nzLower, offset, rowOffset;
247b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
2489ae82921SPaul Mullowney 
2499ae82921SPaul Mullowney   PetscFunctionBegin;
250cf00fe3bSKarl Rupp   if (!n) PetscFunctionReturn(0);
251c70f7ee4SJunchao Zhang   if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) {
2529ae82921SPaul Mullowney     try {
2539ae82921SPaul Mullowney       /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */
2549ae82921SPaul Mullowney       nzLower=n+ai[n]-ai[1];
2559ae82921SPaul Mullowney 
2569ae82921SPaul Mullowney       /* Allocate Space for the lower triangular matrix */
257c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
258c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUDA(ierr);
259c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUDA(ierr);
2609ae82921SPaul Mullowney 
2619ae82921SPaul Mullowney       /* Fill the lower triangular matrix */
2629ae82921SPaul Mullowney       AiLo[0]  = (PetscInt) 0;
2639ae82921SPaul Mullowney       AiLo[n]  = nzLower;
2649ae82921SPaul Mullowney       AjLo[0]  = (PetscInt) 0;
2659ae82921SPaul Mullowney       AALo[0]  = (MatScalar) 1.0;
2669ae82921SPaul Mullowney       v        = aa;
2679ae82921SPaul Mullowney       vi       = aj;
2689ae82921SPaul Mullowney       offset   = 1;
2699ae82921SPaul Mullowney       rowOffset= 1;
2709ae82921SPaul Mullowney       for (i=1; i<n; i++) {
2719ae82921SPaul Mullowney         nz = ai[i+1] - ai[i];
272e057df02SPaul Mullowney         /* additional 1 for the term on the diagonal */
2739ae82921SPaul Mullowney         AiLo[i]    = rowOffset;
2749ae82921SPaul Mullowney         rowOffset += nz+1;
2759ae82921SPaul Mullowney 
276580bdb30SBarry Smith         ierr = PetscArraycpy(&(AjLo[offset]), vi, nz);CHKERRQ(ierr);
277580bdb30SBarry Smith         ierr = PetscArraycpy(&(AALo[offset]), v, nz);CHKERRQ(ierr);
2789ae82921SPaul Mullowney 
2799ae82921SPaul Mullowney         offset      += nz;
2809ae82921SPaul Mullowney         AjLo[offset] = (PetscInt) i;
2819ae82921SPaul Mullowney         AALo[offset] = (MatScalar) 1.0;
2829ae82921SPaul Mullowney         offset      += 1;
2839ae82921SPaul Mullowney 
2849ae82921SPaul Mullowney         v  += nz;
2859ae82921SPaul Mullowney         vi += nz;
2869ae82921SPaul Mullowney       }
2872205254eSKarl Rupp 
288aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
289aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
2902205254eSKarl Rupp 
291aa372e3fSPaul Mullowney       /* Create the matrix description */
292c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
293c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
294c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
295c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUDA(stat);
296c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
297aa372e3fSPaul Mullowney 
298aa372e3fSPaul Mullowney       /* Create the solve analysis information */
299c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
300aa372e3fSPaul Mullowney 
301aa372e3fSPaul Mullowney       /* set the operation */
302aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
303aa372e3fSPaul Mullowney 
304aa372e3fSPaul Mullowney       /* set the matrix */
305aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
306aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = n;
307aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = n;
308aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = nzLower;
309aa372e3fSPaul Mullowney 
310aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
311aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1);
312aa372e3fSPaul Mullowney 
313aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower);
314aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower);
315aa372e3fSPaul Mullowney 
316aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(nzLower);
317aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+nzLower);
318aa372e3fSPaul Mullowney 
319aa372e3fSPaul Mullowney       /* perform the solve analysis */
320aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
321aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
322aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
323c41cb2e2SAlejandro Lamas Daviña                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
324aa372e3fSPaul Mullowney 
325aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
326aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
3272205254eSKarl Rupp 
328c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiLo);CHKERRCUDA(ierr);
329c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjLo);CHKERRCUDA(ierr);
330c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
3314863603aSSatish Balay       ierr = PetscLogCpuToGpu((n+1+nzLower)*sizeof(int)+nzLower*sizeof(PetscScalar));CHKERRQ(ierr);
3329ae82921SPaul Mullowney     } catch(char *ex) {
3339ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
3349ae82921SPaul Mullowney     }
3359ae82921SPaul Mullowney   }
3369ae82921SPaul Mullowney   PetscFunctionReturn(0);
3379ae82921SPaul Mullowney }
3389ae82921SPaul Mullowney 
339087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A)
3409ae82921SPaul Mullowney {
3419ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
3429ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
3439ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
344aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
3459ae82921SPaul Mullowney   cusparseStatus_t                  stat;
3469ae82921SPaul Mullowney   const PetscInt                    *aj = a->j,*adiag = a->diag,*vi;
3479ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
3489ae82921SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
3499ae82921SPaul Mullowney   PetscScalar                       *AAUp;
3509ae82921SPaul Mullowney   PetscInt                          i,nz, nzUpper, offset;
3519ae82921SPaul Mullowney   PetscErrorCode                    ierr;
3529ae82921SPaul Mullowney 
3539ae82921SPaul Mullowney   PetscFunctionBegin;
354cf00fe3bSKarl Rupp   if (!n) PetscFunctionReturn(0);
355c70f7ee4SJunchao Zhang   if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) {
3569ae82921SPaul Mullowney     try {
3579ae82921SPaul Mullowney       /* next, figure out the number of nonzeros in the upper triangular matrix. */
3589ae82921SPaul Mullowney       nzUpper = adiag[0]-adiag[n];
3599ae82921SPaul Mullowney 
3609ae82921SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
361c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
362c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
363c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
3649ae82921SPaul Mullowney 
3659ae82921SPaul Mullowney       /* Fill the upper triangular matrix */
3669ae82921SPaul Mullowney       AiUp[0]=(PetscInt) 0;
3679ae82921SPaul Mullowney       AiUp[n]=nzUpper;
3689ae82921SPaul Mullowney       offset = nzUpper;
3699ae82921SPaul Mullowney       for (i=n-1; i>=0; i--) {
3709ae82921SPaul Mullowney         v  = aa + adiag[i+1] + 1;
3719ae82921SPaul Mullowney         vi = aj + adiag[i+1] + 1;
3729ae82921SPaul Mullowney 
373e057df02SPaul Mullowney         /* number of elements NOT on the diagonal */
3749ae82921SPaul Mullowney         nz = adiag[i] - adiag[i+1]-1;
3759ae82921SPaul Mullowney 
376e057df02SPaul Mullowney         /* decrement the offset */
3779ae82921SPaul Mullowney         offset -= (nz+1);
3789ae82921SPaul Mullowney 
379e057df02SPaul Mullowney         /* first, set the diagonal elements */
3809ae82921SPaul Mullowney         AjUp[offset] = (PetscInt) i;
38109f51544SAlejandro Lamas Daviña         AAUp[offset] = (MatScalar)1./v[nz];
3829ae82921SPaul Mullowney         AiUp[i]      = AiUp[i+1] - (nz+1);
3839ae82921SPaul Mullowney 
384580bdb30SBarry Smith         ierr = PetscArraycpy(&(AjUp[offset+1]), vi, nz);CHKERRQ(ierr);
385580bdb30SBarry Smith         ierr = PetscArraycpy(&(AAUp[offset+1]), v, nz);CHKERRQ(ierr);
3869ae82921SPaul Mullowney       }
3872205254eSKarl Rupp 
388aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
389aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
3902205254eSKarl Rupp 
391aa372e3fSPaul Mullowney       /* Create the matrix description */
392c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
393c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
394c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
395c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
396c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
397aa372e3fSPaul Mullowney 
398aa372e3fSPaul Mullowney       /* Create the solve analysis information */
399c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
400aa372e3fSPaul Mullowney 
401aa372e3fSPaul Mullowney       /* set the operation */
402aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
403aa372e3fSPaul Mullowney 
404aa372e3fSPaul Mullowney       /* set the matrix */
405aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
406aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = n;
407aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = n;
408aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = nzUpper;
409aa372e3fSPaul Mullowney 
410aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
411aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1);
412aa372e3fSPaul Mullowney 
413aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper);
414aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper);
415aa372e3fSPaul Mullowney 
416aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper);
417aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper);
418aa372e3fSPaul Mullowney 
419aa372e3fSPaul Mullowney       /* perform the solve analysis */
420aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
421aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
422aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
423c41cb2e2SAlejandro Lamas Daviña                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
424aa372e3fSPaul Mullowney 
425aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
426aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
4272205254eSKarl Rupp 
428c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
429c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
430c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
4314863603aSSatish Balay       ierr = PetscLogCpuToGpu((n+1+nzUpper)*sizeof(int)+nzUpper*sizeof(PetscScalar));CHKERRQ(ierr);
4329ae82921SPaul Mullowney     } catch(char *ex) {
4339ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
4349ae82921SPaul Mullowney     }
4359ae82921SPaul Mullowney   }
4369ae82921SPaul Mullowney   PetscFunctionReturn(0);
4379ae82921SPaul Mullowney }
4389ae82921SPaul Mullowney 
439087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A)
4409ae82921SPaul Mullowney {
4419ae82921SPaul Mullowney   PetscErrorCode               ierr;
4429ae82921SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
4439ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
4449ae82921SPaul Mullowney   IS                           isrow = a->row,iscol = a->icol;
4459ae82921SPaul Mullowney   PetscBool                    row_identity,col_identity;
4469ae82921SPaul Mullowney   const PetscInt               *r,*c;
4479ae82921SPaul Mullowney   PetscInt                     n = A->rmap->n;
4489ae82921SPaul Mullowney 
4499ae82921SPaul Mullowney   PetscFunctionBegin;
450087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr);
451087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr);
4522205254eSKarl Rupp 
453e65717acSKarl Rupp   cusparseTriFactors->workVector = new THRUSTARRAY(n);
454aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=a->nz;
4559ae82921SPaul Mullowney 
456c70f7ee4SJunchao Zhang   A->offloadmask = PETSC_OFFLOAD_BOTH;
457e057df02SPaul Mullowney   /* lower triangular indices */
4589ae82921SPaul Mullowney   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
4599ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
4602205254eSKarl Rupp   if (!row_identity) {
461aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
462aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(r, r+n);
4632205254eSKarl Rupp   }
4649ae82921SPaul Mullowney   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
4659ae82921SPaul Mullowney 
466e057df02SPaul Mullowney   /* upper triangular indices */
4679ae82921SPaul Mullowney   ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
4689ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4692205254eSKarl Rupp   if (!col_identity) {
470aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
471aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices->assign(c, c+n);
4722205254eSKarl Rupp   }
4734863603aSSatish Balay 
4744863603aSSatish Balay   if (!row_identity && !col_identity) {
4754863603aSSatish Balay     ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr);
4764863603aSSatish Balay   } else if(!row_identity) {
4774863603aSSatish Balay     ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr);
4784863603aSSatish Balay   } else if(!col_identity) {
4794863603aSSatish Balay     ierr = PetscLogCpuToGpu(n*sizeof(PetscInt));CHKERRQ(ierr);
4804863603aSSatish Balay   }
4814863603aSSatish Balay 
4829ae82921SPaul Mullowney   ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
4839ae82921SPaul Mullowney   PetscFunctionReturn(0);
4849ae82921SPaul Mullowney }
4859ae82921SPaul Mullowney 
486087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A)
487087f3262SPaul Mullowney {
488087f3262SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
489087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
490aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
491aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
492087f3262SPaul Mullowney   cusparseStatus_t                  stat;
493087f3262SPaul Mullowney   PetscErrorCode                    ierr;
494087f3262SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
495087f3262SPaul Mullowney   PetscScalar                       *AAUp;
496087f3262SPaul Mullowney   PetscScalar                       *AALo;
497087f3262SPaul Mullowney   PetscInt                          nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j;
498087f3262SPaul Mullowney   Mat_SeqSBAIJ                      *b = (Mat_SeqSBAIJ*)A->data;
499087f3262SPaul Mullowney   const PetscInt                    *ai = b->i,*aj = b->j,*vj;
500087f3262SPaul Mullowney   const MatScalar                   *aa = b->a,*v;
501087f3262SPaul Mullowney 
502087f3262SPaul Mullowney   PetscFunctionBegin;
503cf00fe3bSKarl Rupp   if (!n) PetscFunctionReturn(0);
504c70f7ee4SJunchao Zhang   if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) {
505087f3262SPaul Mullowney     try {
506087f3262SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
507c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
508c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
509c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
510c41cb2e2SAlejandro Lamas Daviña       ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
511087f3262SPaul Mullowney 
512087f3262SPaul Mullowney       /* Fill the upper triangular matrix */
513087f3262SPaul Mullowney       AiUp[0]=(PetscInt) 0;
514087f3262SPaul Mullowney       AiUp[n]=nzUpper;
515087f3262SPaul Mullowney       offset = 0;
516087f3262SPaul Mullowney       for (i=0; i<n; i++) {
517087f3262SPaul Mullowney         /* set the pointers */
518087f3262SPaul Mullowney         v  = aa + ai[i];
519087f3262SPaul Mullowney         vj = aj + ai[i];
520087f3262SPaul Mullowney         nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */
521087f3262SPaul Mullowney 
522087f3262SPaul Mullowney         /* first, set the diagonal elements */
523087f3262SPaul Mullowney         AjUp[offset] = (PetscInt) i;
52409f51544SAlejandro Lamas Daviña         AAUp[offset] = (MatScalar)1.0/v[nz];
525087f3262SPaul Mullowney         AiUp[i]      = offset;
52609f51544SAlejandro Lamas Daviña         AALo[offset] = (MatScalar)1.0/v[nz];
527087f3262SPaul Mullowney 
528087f3262SPaul Mullowney         offset+=1;
529087f3262SPaul Mullowney         if (nz>0) {
530f22e0265SBarry Smith           ierr = PetscArraycpy(&(AjUp[offset]), vj, nz);CHKERRQ(ierr);
531580bdb30SBarry Smith           ierr = PetscArraycpy(&(AAUp[offset]), v, nz);CHKERRQ(ierr);
532087f3262SPaul Mullowney           for (j=offset; j<offset+nz; j++) {
533087f3262SPaul Mullowney             AAUp[j] = -AAUp[j];
534087f3262SPaul Mullowney             AALo[j] = AAUp[j]/v[nz];
535087f3262SPaul Mullowney           }
536087f3262SPaul Mullowney           offset+=nz;
537087f3262SPaul Mullowney         }
538087f3262SPaul Mullowney       }
539087f3262SPaul Mullowney 
540aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
541aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
542087f3262SPaul Mullowney 
543aa372e3fSPaul Mullowney       /* Create the matrix description */
544c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
545c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
546c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
547c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
548c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
549087f3262SPaul Mullowney 
550aa372e3fSPaul Mullowney       /* Create the solve analysis information */
551c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
552aa372e3fSPaul Mullowney 
553aa372e3fSPaul Mullowney       /* set the operation */
554aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
555aa372e3fSPaul Mullowney 
556aa372e3fSPaul Mullowney       /* set the matrix */
557aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
558aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = A->rmap->n;
559aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = A->cmap->n;
560aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = a->nz;
561aa372e3fSPaul Mullowney 
562aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
563aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
564aa372e3fSPaul Mullowney 
565aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
566aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
567aa372e3fSPaul Mullowney 
568aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
569aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz);
570aa372e3fSPaul Mullowney 
571aa372e3fSPaul Mullowney       /* perform the solve analysis */
572aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
573aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
574aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
575c41cb2e2SAlejandro Lamas Daviña                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
576aa372e3fSPaul Mullowney 
577aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
578aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
579aa372e3fSPaul Mullowney 
580aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
581aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
582aa372e3fSPaul Mullowney 
583aa372e3fSPaul Mullowney       /* Create the matrix description */
584c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
585c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
586c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
587c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
588c41cb2e2SAlejandro Lamas Daviña       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
589aa372e3fSPaul Mullowney 
590aa372e3fSPaul Mullowney       /* Create the solve analysis information */
591c41cb2e2SAlejandro Lamas Daviña       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
592aa372e3fSPaul Mullowney 
593aa372e3fSPaul Mullowney       /* set the operation */
594aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE;
595aa372e3fSPaul Mullowney 
596aa372e3fSPaul Mullowney       /* set the matrix */
597aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
598aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = A->rmap->n;
599aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = A->cmap->n;
600aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = a->nz;
601aa372e3fSPaul Mullowney 
602aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
603aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
604aa372e3fSPaul Mullowney 
605aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
606aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
607aa372e3fSPaul Mullowney 
608aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
609aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+a->nz);
6104863603aSSatish Balay       ierr = PetscLogCpuToGpu(2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)));CHKERRQ(ierr);
611aa372e3fSPaul Mullowney 
612aa372e3fSPaul Mullowney       /* perform the solve analysis */
613aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
614aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
615aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
616c41cb2e2SAlejandro Lamas Daviña                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
617aa372e3fSPaul Mullowney 
618aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
619aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
620087f3262SPaul Mullowney 
621c70f7ee4SJunchao Zhang       A->offloadmask = PETSC_OFFLOAD_BOTH;
622c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
623c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
624c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
625c41cb2e2SAlejandro Lamas Daviña       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
626087f3262SPaul Mullowney     } catch(char *ex) {
627087f3262SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
628087f3262SPaul Mullowney     }
629087f3262SPaul Mullowney   }
630087f3262SPaul Mullowney   PetscFunctionReturn(0);
631087f3262SPaul Mullowney }
632087f3262SPaul Mullowney 
633087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A)
6349ae82921SPaul Mullowney {
6359ae82921SPaul Mullowney   PetscErrorCode               ierr;
636087f3262SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
637087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
638087f3262SPaul Mullowney   IS                           ip = a->row;
639087f3262SPaul Mullowney   const PetscInt               *rip;
640087f3262SPaul Mullowney   PetscBool                    perm_identity;
641087f3262SPaul Mullowney   PetscInt                     n = A->rmap->n;
642087f3262SPaul Mullowney 
643087f3262SPaul Mullowney   PetscFunctionBegin;
644087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr);
645e65717acSKarl Rupp   cusparseTriFactors->workVector = new THRUSTARRAY(n);
646aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=(a->nz-n)*2 + n;
647aa372e3fSPaul Mullowney 
648087f3262SPaul Mullowney   /* lower triangular indices */
649087f3262SPaul Mullowney   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
650087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
651087f3262SPaul Mullowney   if (!perm_identity) {
6524e4bbfaaSStefano Zampini     IS             iip;
6534e4bbfaaSStefano Zampini     const PetscInt *irip;
6544e4bbfaaSStefano Zampini 
6554e4bbfaaSStefano Zampini     ierr = ISInvertPermutation(ip,PETSC_DECIDE,&iip);CHKERRQ(ierr);
6564e4bbfaaSStefano Zampini     ierr = ISGetIndices(iip,&irip);CHKERRQ(ierr);
657aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
658aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(rip, rip+n);
659aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
6604e4bbfaaSStefano Zampini     cusparseTriFactors->cpermIndices->assign(irip, irip+n);
6614e4bbfaaSStefano Zampini     ierr = ISRestoreIndices(iip,&irip);CHKERRQ(ierr);
6624e4bbfaaSStefano Zampini     ierr = ISDestroy(&iip);CHKERRQ(ierr);
6634863603aSSatish Balay     ierr = PetscLogCpuToGpu(2*n*sizeof(PetscInt));CHKERRQ(ierr);
664087f3262SPaul Mullowney   }
665087f3262SPaul Mullowney   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
666087f3262SPaul Mullowney   PetscFunctionReturn(0);
667087f3262SPaul Mullowney }
668087f3262SPaul Mullowney 
6696fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
6709ae82921SPaul Mullowney {
6719ae82921SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
6729ae82921SPaul Mullowney   IS             isrow = b->row,iscol = b->col;
6739ae82921SPaul Mullowney   PetscBool      row_identity,col_identity;
674b175d8bbSPaul Mullowney   PetscErrorCode ierr;
6759ae82921SPaul Mullowney 
6769ae82921SPaul Mullowney   PetscFunctionBegin;
6779ae82921SPaul Mullowney   ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
678e057df02SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
6799ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
6809ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
681bda325fcSPaul Mullowney   if (row_identity && col_identity) {
682bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
683bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
6844e4bbfaaSStefano Zampini     B->ops->matsolve = NULL;
6854e4bbfaaSStefano Zampini     B->ops->matsolvetranspose = NULL;
686bda325fcSPaul Mullowney   } else {
687bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
688bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
6894e4bbfaaSStefano Zampini     B->ops->matsolve = NULL;
6904e4bbfaaSStefano Zampini     B->ops->matsolvetranspose = NULL;
691bda325fcSPaul Mullowney   }
6928dc1d2a3SPaul Mullowney 
693e057df02SPaul Mullowney   /* get the triangular factors */
694087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
6959ae82921SPaul Mullowney   PetscFunctionReturn(0);
6969ae82921SPaul Mullowney }
6979ae82921SPaul Mullowney 
698087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
699087f3262SPaul Mullowney {
700087f3262SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
701087f3262SPaul Mullowney   IS             ip = b->row;
702087f3262SPaul Mullowney   PetscBool      perm_identity;
703b175d8bbSPaul Mullowney   PetscErrorCode ierr;
704087f3262SPaul Mullowney 
705087f3262SPaul Mullowney   PetscFunctionBegin;
706087f3262SPaul Mullowney   ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
707087f3262SPaul Mullowney 
708087f3262SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
709087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
710087f3262SPaul Mullowney   if (perm_identity) {
711087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
712087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
7134e4bbfaaSStefano Zampini     B->ops->matsolve = NULL;
7144e4bbfaaSStefano Zampini     B->ops->matsolvetranspose = NULL;
715087f3262SPaul Mullowney   } else {
716087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
717087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
7184e4bbfaaSStefano Zampini     B->ops->matsolve = NULL;
7194e4bbfaaSStefano Zampini     B->ops->matsolvetranspose = NULL;
720087f3262SPaul Mullowney   }
721087f3262SPaul Mullowney 
722087f3262SPaul Mullowney   /* get the triangular factors */
723087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
724087f3262SPaul Mullowney   PetscFunctionReturn(0);
725087f3262SPaul Mullowney }
7269ae82921SPaul Mullowney 
727b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A)
728bda325fcSPaul Mullowney {
729bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
730aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
731aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
732aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
733aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
734bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
735aa372e3fSPaul Mullowney   cusparseIndexBase_t               indexBase;
736aa372e3fSPaul Mullowney   cusparseMatrixType_t              matrixType;
737aa372e3fSPaul Mullowney   cusparseFillMode_t                fillMode;
738aa372e3fSPaul Mullowney   cusparseDiagType_t                diagType;
739b175d8bbSPaul Mullowney 
740bda325fcSPaul Mullowney   PetscFunctionBegin;
741bda325fcSPaul Mullowney 
742aa372e3fSPaul Mullowney   /*********************************************/
743aa372e3fSPaul Mullowney   /* Now the Transpose of the Lower Tri Factor */
744aa372e3fSPaul Mullowney   /*********************************************/
745aa372e3fSPaul Mullowney 
746aa372e3fSPaul Mullowney   /* allocate space for the transpose of the lower triangular factor */
747aa372e3fSPaul Mullowney   loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
748aa372e3fSPaul Mullowney 
749aa372e3fSPaul Mullowney   /* set the matrix descriptors of the lower triangular factor */
750aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(loTriFactor->descr);
751aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(loTriFactor->descr);
752aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
753aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
754aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(loTriFactor->descr);
755aa372e3fSPaul Mullowney 
756aa372e3fSPaul Mullowney   /* Create the matrix description */
757c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat);
758c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat);
759c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat);
760c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat);
761c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat);
762aa372e3fSPaul Mullowney 
763aa372e3fSPaul Mullowney   /* Create the solve analysis information */
764c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat);
765aa372e3fSPaul Mullowney 
766aa372e3fSPaul Mullowney   /* set the operation */
767aa372e3fSPaul Mullowney   loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
768aa372e3fSPaul Mullowney 
769aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the lower triangular factor*/
770aa372e3fSPaul Mullowney   loTriFactorT->csrMat = new CsrMatrix;
771aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows;
772aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols;
773aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries;
774aa372e3fSPaul Mullowney   loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1);
775aa372e3fSPaul Mullowney   loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries);
776aa372e3fSPaul Mullowney   loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries);
777aa372e3fSPaul Mullowney 
778aa372e3fSPaul Mullowney   /* compute the transpose of the lower triangular factor, i.e. the CSC */
779aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows,
780aa372e3fSPaul Mullowney                           loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries,
781aa372e3fSPaul Mullowney                           loTriFactor->csrMat->values->data().get(),
782aa372e3fSPaul Mullowney                           loTriFactor->csrMat->row_offsets->data().get(),
783aa372e3fSPaul Mullowney                           loTriFactor->csrMat->column_indices->data().get(),
784aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->values->data().get(),
785aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->column_indices->data().get(),
786aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->row_offsets->data().get(),
787c41cb2e2SAlejandro Lamas Daviña                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
788aa372e3fSPaul Mullowney 
789aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
790aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp,
791aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries,
792aa372e3fSPaul Mullowney                            loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(),
793aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(),
794c41cb2e2SAlejandro Lamas Daviña                            loTriFactorT->solveInfo);CHKERRCUDA(stat);
795aa372e3fSPaul Mullowney 
796aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
797aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT;
798aa372e3fSPaul Mullowney 
799aa372e3fSPaul Mullowney   /*********************************************/
800aa372e3fSPaul Mullowney   /* Now the Transpose of the Upper Tri Factor */
801aa372e3fSPaul Mullowney   /*********************************************/
802aa372e3fSPaul Mullowney 
803aa372e3fSPaul Mullowney   /* allocate space for the transpose of the upper triangular factor */
804aa372e3fSPaul Mullowney   upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
805aa372e3fSPaul Mullowney 
806aa372e3fSPaul Mullowney   /* set the matrix descriptors of the upper triangular factor */
807aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(upTriFactor->descr);
808aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(upTriFactor->descr);
809aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
810aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
811aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(upTriFactor->descr);
812aa372e3fSPaul Mullowney 
813aa372e3fSPaul Mullowney   /* Create the matrix description */
814c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat);
815c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat);
816c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat);
817c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat);
818c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat);
819aa372e3fSPaul Mullowney 
820aa372e3fSPaul Mullowney   /* Create the solve analysis information */
821c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat);
822aa372e3fSPaul Mullowney 
823aa372e3fSPaul Mullowney   /* set the operation */
824aa372e3fSPaul Mullowney   upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
825aa372e3fSPaul Mullowney 
826aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the upper triangular factor*/
827aa372e3fSPaul Mullowney   upTriFactorT->csrMat = new CsrMatrix;
828aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows;
829aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols;
830aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries;
831aa372e3fSPaul Mullowney   upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1);
832aa372e3fSPaul Mullowney   upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries);
833aa372e3fSPaul Mullowney   upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries);
834aa372e3fSPaul Mullowney 
835aa372e3fSPaul Mullowney   /* compute the transpose of the upper triangular factor, i.e. the CSC */
836aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows,
837aa372e3fSPaul Mullowney                           upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries,
838aa372e3fSPaul Mullowney                           upTriFactor->csrMat->values->data().get(),
839aa372e3fSPaul Mullowney                           upTriFactor->csrMat->row_offsets->data().get(),
840aa372e3fSPaul Mullowney                           upTriFactor->csrMat->column_indices->data().get(),
841aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->values->data().get(),
842aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->column_indices->data().get(),
843aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->row_offsets->data().get(),
844c41cb2e2SAlejandro Lamas Daviña                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
845aa372e3fSPaul Mullowney 
846aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
847aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp,
848aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries,
849aa372e3fSPaul Mullowney                            upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(),
850aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(),
851c41cb2e2SAlejandro Lamas Daviña                            upTriFactorT->solveInfo);CHKERRCUDA(stat);
852aa372e3fSPaul Mullowney 
853aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
854aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT;
855bda325fcSPaul Mullowney   PetscFunctionReturn(0);
856bda325fcSPaul Mullowney }
857bda325fcSPaul Mullowney 
858b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A)
859bda325fcSPaul Mullowney {
860aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
861aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
862aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
863bda325fcSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
864bda325fcSPaul Mullowney   cusparseStatus_t             stat;
865aa372e3fSPaul Mullowney   cusparseIndexBase_t          indexBase;
866b06137fdSPaul Mullowney   cudaError_t                  err;
8674863603aSSatish Balay   PetscErrorCode               ierr;
868b175d8bbSPaul Mullowney 
869bda325fcSPaul Mullowney   PetscFunctionBegin;
870aa372e3fSPaul Mullowney 
871aa372e3fSPaul Mullowney   /* allocate space for the triangular factor information */
872aa372e3fSPaul Mullowney   matstructT = new Mat_SeqAIJCUSPARSEMultStruct;
873c41cb2e2SAlejandro Lamas Daviña   stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat);
874aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(matstruct->descr);
875c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat);
876c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
877aa372e3fSPaul Mullowney 
878b06137fdSPaul Mullowney   /* set alpha and beta */
879c41cb2e2SAlejandro Lamas Daviña   err = cudaMalloc((void **)&(matstructT->alpha),    sizeof(PetscScalar));CHKERRCUDA(err);
8807656d835SStefano Zampini   err = cudaMalloc((void **)&(matstructT->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err);
8817656d835SStefano Zampini   err = cudaMalloc((void **)&(matstructT->beta_one), sizeof(PetscScalar));CHKERRCUDA(err);
8827656d835SStefano Zampini   err = cudaMemcpy(matstructT->alpha,    &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
8837656d835SStefano Zampini   err = cudaMemcpy(matstructT->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
8847656d835SStefano Zampini   err = cudaMemcpy(matstructT->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
885c41cb2e2SAlejandro Lamas Daviña   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
886b06137fdSPaul Mullowney 
887aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
888aa372e3fSPaul Mullowney     CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
889aa372e3fSPaul Mullowney     CsrMatrix *matrixT= new CsrMatrix;
890554b8892SKarl Rupp     matrixT->num_rows = A->cmap->n;
891554b8892SKarl Rupp     matrixT->num_cols = A->rmap->n;
892aa372e3fSPaul Mullowney     matrixT->num_entries = a->nz;
893a8bd5306SMark Adams     matrixT->row_offsets = new THRUSTINTARRAY32(matrixT->num_rows+1);
894aa372e3fSPaul Mullowney     matrixT->column_indices = new THRUSTINTARRAY32(a->nz);
895aa372e3fSPaul Mullowney     matrixT->values = new THRUSTARRAY(a->nz);
896a3fdcf43SKarl Rupp 
89781902715SJunchao Zhang     cusparsestruct->rowoffsets_gpu = new THRUSTINTARRAY32(A->rmap->n+1);
89881902715SJunchao Zhang     cusparsestruct->rowoffsets_gpu->assign(a->i,a->i+A->rmap->n+1);
89981902715SJunchao Zhang     /* compute the transpose, i.e. the CSC */
900aa372e3fSPaul Mullowney     indexBase = cusparseGetMatIndexBase(matstruct->descr);
901a3fdcf43SKarl Rupp     stat = cusparse_csr2csc(cusparsestruct->handle, A->rmap->n,
902a3fdcf43SKarl Rupp                             A->cmap->n, matrix->num_entries,
903aa372e3fSPaul Mullowney                             matrix->values->data().get(),
90481902715SJunchao Zhang                             cusparsestruct->rowoffsets_gpu->data().get(),
905aa372e3fSPaul Mullowney                             matrix->column_indices->data().get(),
906aa372e3fSPaul Mullowney                             matrixT->values->data().get(),
907aa372e3fSPaul Mullowney                             matrixT->column_indices->data().get(),
908aa372e3fSPaul Mullowney                             matrixT->row_offsets->data().get(),
909c41cb2e2SAlejandro Lamas Daviña                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
910aa372e3fSPaul Mullowney     /* assign the pointer */
911aa372e3fSPaul Mullowney     matstructT->mat = matrixT;
9124863603aSSatish Balay     ierr = PetscLogCpuToGpu(((A->rmap->n+1)+(a->nz))*sizeof(int)+(3+a->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
913aa372e3fSPaul Mullowney   } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
914aa372e3fSPaul Mullowney     /* First convert HYB to CSR */
915aa372e3fSPaul Mullowney     CsrMatrix *temp= new CsrMatrix;
916aa372e3fSPaul Mullowney     temp->num_rows = A->rmap->n;
917aa372e3fSPaul Mullowney     temp->num_cols = A->cmap->n;
918aa372e3fSPaul Mullowney     temp->num_entries = a->nz;
919aa372e3fSPaul Mullowney     temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
920aa372e3fSPaul Mullowney     temp->column_indices = new THRUSTINTARRAY32(a->nz);
921aa372e3fSPaul Mullowney     temp->values = new THRUSTARRAY(a->nz);
922aa372e3fSPaul Mullowney 
9232692e278SPaul Mullowney 
924aa372e3fSPaul Mullowney     stat = cusparse_hyb2csr(cusparsestruct->handle,
925aa372e3fSPaul Mullowney                             matstruct->descr, (cusparseHybMat_t)matstruct->mat,
926aa372e3fSPaul Mullowney                             temp->values->data().get(),
927aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
928c41cb2e2SAlejandro Lamas Daviña                             temp->column_indices->data().get());CHKERRCUDA(stat);
929aa372e3fSPaul Mullowney 
930aa372e3fSPaul Mullowney     /* Next, convert CSR to CSC (i.e. the matrix transpose) */
931aa372e3fSPaul Mullowney     CsrMatrix *tempT= new CsrMatrix;
932aa372e3fSPaul Mullowney     tempT->num_rows = A->rmap->n;
933aa372e3fSPaul Mullowney     tempT->num_cols = A->cmap->n;
934aa372e3fSPaul Mullowney     tempT->num_entries = a->nz;
935aa372e3fSPaul Mullowney     tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
936aa372e3fSPaul Mullowney     tempT->column_indices = new THRUSTINTARRAY32(a->nz);
937aa372e3fSPaul Mullowney     tempT->values = new THRUSTARRAY(a->nz);
938aa372e3fSPaul Mullowney 
939aa372e3fSPaul Mullowney     stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows,
940aa372e3fSPaul Mullowney                             temp->num_cols, temp->num_entries,
941aa372e3fSPaul Mullowney                             temp->values->data().get(),
942aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
943aa372e3fSPaul Mullowney                             temp->column_indices->data().get(),
944aa372e3fSPaul Mullowney                             tempT->values->data().get(),
945aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
946aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
947c41cb2e2SAlejandro Lamas Daviña                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
948aa372e3fSPaul Mullowney 
949aa372e3fSPaul Mullowney     /* Last, convert CSC to HYB */
950aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat;
951c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
952aa372e3fSPaul Mullowney     cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
953aa372e3fSPaul Mullowney       CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
954aa372e3fSPaul Mullowney     stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n,
955aa372e3fSPaul Mullowney                             matstructT->descr, tempT->values->data().get(),
956aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
957aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
958c41cb2e2SAlejandro Lamas Daviña                             hybMat, 0, partition);CHKERRCUDA(stat);
959aa372e3fSPaul Mullowney 
960aa372e3fSPaul Mullowney     /* assign the pointer */
961aa372e3fSPaul Mullowney     matstructT->mat = hybMat;
9624863603aSSatish Balay     ierr = PetscLogCpuToGpu((2*(((A->rmap->n+1)+(a->nz))*sizeof(int)+(a->nz)*sizeof(PetscScalar)))+3*sizeof(PetscScalar));CHKERRQ(ierr);
963aa372e3fSPaul Mullowney 
964aa372e3fSPaul Mullowney     /* delete temporaries */
965aa372e3fSPaul Mullowney     if (tempT) {
966aa372e3fSPaul Mullowney       if (tempT->values) delete (THRUSTARRAY*) tempT->values;
967aa372e3fSPaul Mullowney       if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices;
968aa372e3fSPaul Mullowney       if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets;
969aa372e3fSPaul Mullowney       delete (CsrMatrix*) tempT;
970087f3262SPaul Mullowney     }
971aa372e3fSPaul Mullowney     if (temp) {
972aa372e3fSPaul Mullowney       if (temp->values) delete (THRUSTARRAY*) temp->values;
973aa372e3fSPaul Mullowney       if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices;
974aa372e3fSPaul Mullowney       if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets;
975aa372e3fSPaul Mullowney       delete (CsrMatrix*) temp;
976aa372e3fSPaul Mullowney     }
977aa372e3fSPaul Mullowney   }
978aa372e3fSPaul Mullowney   /* assign the compressed row indices */
979aa372e3fSPaul Mullowney   matstructT->cprowIndices = new THRUSTINTARRAY;
980554b8892SKarl Rupp   matstructT->cprowIndices->resize(A->cmap->n);
981554b8892SKarl Rupp   thrust::sequence(matstructT->cprowIndices->begin(), matstructT->cprowIndices->end());
982aa372e3fSPaul Mullowney   /* assign the pointer */
983aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT;
984bda325fcSPaul Mullowney   PetscFunctionReturn(0);
985bda325fcSPaul Mullowney }
986bda325fcSPaul Mullowney 
9874e4bbfaaSStefano Zampini /* Why do we need to analyze the tranposed matrix again? Can't we just use op(A) = CUSPARSE_OPERATION_TRANSPOSE in MatSolve_SeqAIJCUSPARSE? */
9886fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
989bda325fcSPaul Mullowney {
990c41cb2e2SAlejandro Lamas Daviña   PetscInt                              n = xx->map->n;
991465f34aeSAlejandro Lamas Daviña   const PetscScalar                     *barray;
992465f34aeSAlejandro Lamas Daviña   PetscScalar                           *xarray;
993465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<const PetscScalar> bGPU;
994465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<PetscScalar>       xGPU;
995bda325fcSPaul Mullowney   cusparseStatus_t                      stat;
996bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
997aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
998aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
999aa372e3fSPaul Mullowney   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1000b175d8bbSPaul Mullowney   PetscErrorCode                        ierr;
1001bda325fcSPaul Mullowney 
1002bda325fcSPaul Mullowney   PetscFunctionBegin;
1003aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
1004aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
1005bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1006aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1007aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1008bda325fcSPaul Mullowney   }
1009bda325fcSPaul Mullowney 
1010bda325fcSPaul Mullowney   /* Get the GPU pointers */
1011c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1012c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1013c41cb2e2SAlejandro Lamas Daviña   xGPU = thrust::device_pointer_cast(xarray);
1014c41cb2e2SAlejandro Lamas Daviña   bGPU = thrust::device_pointer_cast(barray);
1015bda325fcSPaul Mullowney 
10167a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1017aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
1018c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
1019c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()),
1020c41cb2e2SAlejandro Lamas Daviña                xGPU);
1021aa372e3fSPaul Mullowney 
1022aa372e3fSPaul Mullowney   /* First, solve U */
1023aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
10247656d835SStefano Zampini                         upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr,
1025aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1026aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1027aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1028aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1029c41cb2e2SAlejandro Lamas Daviña                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1030aa372e3fSPaul Mullowney 
1031aa372e3fSPaul Mullowney   /* Then, solve L */
1032aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
10337656d835SStefano Zampini                         loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr,
1034aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1035aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1036aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1037aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1038c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1039aa372e3fSPaul Mullowney 
1040aa372e3fSPaul Mullowney   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1041c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()),
1042c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()),
1043aa372e3fSPaul Mullowney                tempGPU->begin());
1044aa372e3fSPaul Mullowney 
1045aa372e3fSPaul Mullowney   /* Copy the temporary to the full solution. */
1046c41cb2e2SAlejandro Lamas Daviña   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU);
1047bda325fcSPaul Mullowney 
1048bda325fcSPaul Mullowney   /* restore */
1049c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1050c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1051c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1052661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1053958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1054bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1055bda325fcSPaul Mullowney }
1056bda325fcSPaul Mullowney 
10576fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
1058bda325fcSPaul Mullowney {
1059465f34aeSAlejandro Lamas Daviña   const PetscScalar                 *barray;
1060465f34aeSAlejandro Lamas Daviña   PetscScalar                       *xarray;
1061bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
1062bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1063aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1064aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1065aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1066b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
1067bda325fcSPaul Mullowney 
1068bda325fcSPaul Mullowney   PetscFunctionBegin;
1069aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
1070aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
1071bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1072aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1073aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1074bda325fcSPaul Mullowney   }
1075bda325fcSPaul Mullowney 
1076bda325fcSPaul Mullowney   /* Get the GPU pointers */
1077c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1078c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1079bda325fcSPaul Mullowney 
10807a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1081aa372e3fSPaul Mullowney   /* First, solve U */
1082aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
10837656d835SStefano Zampini                         upTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactorT->descr,
1084aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1085aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1086aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1087aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1088c41cb2e2SAlejandro Lamas Daviña                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1089aa372e3fSPaul Mullowney 
1090aa372e3fSPaul Mullowney   /* Then, solve L */
1091aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
10927656d835SStefano Zampini                         loTriFactorT->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactorT->descr,
1093aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1094aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1095aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1096aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1097c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1098bda325fcSPaul Mullowney 
1099bda325fcSPaul Mullowney   /* restore */
1100c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1101c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1102c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1103661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1104958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1105bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1106bda325fcSPaul Mullowney }
1107bda325fcSPaul Mullowney 
11086fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
11099ae82921SPaul Mullowney {
1110465f34aeSAlejandro Lamas Daviña   const PetscScalar                     *barray;
1111465f34aeSAlejandro Lamas Daviña   PetscScalar                           *xarray;
1112465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<const PetscScalar> bGPU;
1113465f34aeSAlejandro Lamas Daviña   thrust::device_ptr<PetscScalar>       xGPU;
11149ae82921SPaul Mullowney   cusparseStatus_t                      stat;
11159ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1116aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1117aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1118aa372e3fSPaul Mullowney   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1119b175d8bbSPaul Mullowney   PetscErrorCode                        ierr;
11209ae82921SPaul Mullowney 
11219ae82921SPaul Mullowney   PetscFunctionBegin;
1122ebc8f436SDominic Meiser 
1123e057df02SPaul Mullowney   /* Get the GPU pointers */
1124c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1125c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1126c41cb2e2SAlejandro Lamas Daviña   xGPU = thrust::device_pointer_cast(xarray);
1127c41cb2e2SAlejandro Lamas Daviña   bGPU = thrust::device_pointer_cast(barray);
11289ae82921SPaul Mullowney 
11297a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1130aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
1131c41cb2e2SAlejandro Lamas Daviña   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
1132c41cb2e2SAlejandro Lamas Daviña                thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()),
11334e4bbfaaSStefano Zampini                tempGPU->begin());
1134aa372e3fSPaul Mullowney 
1135aa372e3fSPaul Mullowney   /* Next, solve L */
1136aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
11377656d835SStefano Zampini                         loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr,
1138aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1139aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1140aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1141aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
11424e4bbfaaSStefano Zampini                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1143aa372e3fSPaul Mullowney 
1144aa372e3fSPaul Mullowney   /* Then, solve U */
1145aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
11467656d835SStefano Zampini                         upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr,
1147aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1148aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1149aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1150aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
11514e4bbfaaSStefano Zampini                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1152aa372e3fSPaul Mullowney 
11534e4bbfaaSStefano Zampini   /* Last, reorder with the column permutation */
11544e4bbfaaSStefano Zampini   thrust::copy(thrust::make_permutation_iterator(tempGPU->begin(), cusparseTriFactors->cpermIndices->begin()),
11554e4bbfaaSStefano Zampini                thrust::make_permutation_iterator(tempGPU->begin(), cusparseTriFactors->cpermIndices->end()),
11564e4bbfaaSStefano Zampini                xGPU);
11579ae82921SPaul Mullowney 
1158c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1159c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1160c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1161661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1162958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
11639ae82921SPaul Mullowney   PetscFunctionReturn(0);
11649ae82921SPaul Mullowney }
11659ae82921SPaul Mullowney 
11666fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
11679ae82921SPaul Mullowney {
1168465f34aeSAlejandro Lamas Daviña   const PetscScalar                 *barray;
1169465f34aeSAlejandro Lamas Daviña   PetscScalar                       *xarray;
11709ae82921SPaul Mullowney   cusparseStatus_t                  stat;
11719ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1172aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1173aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1174aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1175b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
11769ae82921SPaul Mullowney 
11779ae82921SPaul Mullowney   PetscFunctionBegin;
1178e057df02SPaul Mullowney   /* Get the GPU pointers */
1179c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1180c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
11819ae82921SPaul Mullowney 
11827a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1183aa372e3fSPaul Mullowney   /* First, solve L */
1184aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
11857656d835SStefano Zampini                         loTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, loTriFactor->descr,
1186aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1187aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1188aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1189aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
1190c41cb2e2SAlejandro Lamas Daviña                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1191aa372e3fSPaul Mullowney 
1192aa372e3fSPaul Mullowney   /* Next, solve U */
1193aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
11947656d835SStefano Zampini                         upTriFactor->csrMat->num_rows, &PETSC_CUSPARSE_ONE, upTriFactor->descr,
1195aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1196aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1197aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1198aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
1199c41cb2e2SAlejandro Lamas Daviña                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
12009ae82921SPaul Mullowney 
1201c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1202c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1203c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1204661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1205958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
12069ae82921SPaul Mullowney   PetscFunctionReturn(0);
12079ae82921SPaul Mullowney }
12089ae82921SPaul Mullowney 
12096fa9248bSJed Brown static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A)
12109ae82921SPaul Mullowney {
1211aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
12127c700b8dSJunchao Zhang   Mat_SeqAIJCUSPARSEMultStruct *matstruct = cusparsestruct->mat;
12139ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
12149ae82921SPaul Mullowney   PetscInt                     m = A->rmap->n,*ii,*ridx;
12159ae82921SPaul Mullowney   PetscErrorCode               ierr;
1216aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1217b06137fdSPaul Mullowney   cudaError_t                  err;
12189ae82921SPaul Mullowney 
12199ae82921SPaul Mullowney   PetscFunctionBegin;
1220c70f7ee4SJunchao Zhang   if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) {
12219ae82921SPaul Mullowney     ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
122281902715SJunchao Zhang     if (A->was_assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) {
122381902715SJunchao Zhang       /* Copy values only */
122481902715SJunchao Zhang       CsrMatrix *mat,*matT;
122581902715SJunchao Zhang       mat  = (CsrMatrix*)cusparsestruct->mat->mat;
122681902715SJunchao Zhang       mat->values->assign(a->a, a->a+a->nz);
12274863603aSSatish Balay       ierr = PetscLogCpuToGpu((a->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
122881902715SJunchao Zhang 
122981902715SJunchao Zhang       /* Update matT when it was built before */
123081902715SJunchao Zhang       if (cusparsestruct->matTranspose) {
123181902715SJunchao Zhang         cusparseIndexBase_t indexBase = cusparseGetMatIndexBase(cusparsestruct->mat->descr);
123281902715SJunchao Zhang         matT = (CsrMatrix*)cusparsestruct->matTranspose->mat;
123381902715SJunchao Zhang         stat = cusparse_csr2csc(cusparsestruct->handle, A->rmap->n,
123481902715SJunchao Zhang                                  A->cmap->n, mat->num_entries,
123581902715SJunchao Zhang                                  mat->values->data().get(),
123681902715SJunchao Zhang                                  cusparsestruct->rowoffsets_gpu->data().get(),
123781902715SJunchao Zhang                                  mat->column_indices->data().get(),
123881902715SJunchao Zhang                                  matT->values->data().get(),
123981902715SJunchao Zhang                                  matT->column_indices->data().get(),
124081902715SJunchao Zhang                                  matT->row_offsets->data().get(),
124181902715SJunchao Zhang                                  CUSPARSE_ACTION_NUMERIC,indexBase);CHKERRCUDA(stat);
124281902715SJunchao Zhang       }
124334d6c7a5SJose E. Roman     } else {
12447c700b8dSJunchao Zhang       ierr = MatSeqAIJCUSPARSEMultStruct_Destroy(&cusparsestruct->mat,cusparsestruct->format);CHKERRQ(ierr);
12457c700b8dSJunchao Zhang       ierr = MatSeqAIJCUSPARSEMultStruct_Destroy(&cusparsestruct->matTranspose,cusparsestruct->format);CHKERRQ(ierr);
12467c700b8dSJunchao Zhang       delete cusparsestruct->workVector;
124781902715SJunchao Zhang       delete cusparsestruct->rowoffsets_gpu;
12489ae82921SPaul Mullowney       try {
1249aa372e3fSPaul Mullowney         cusparsestruct->nonzerorow=0;
1250aa372e3fSPaul Mullowney         for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0);
12519ae82921SPaul Mullowney 
12529ae82921SPaul Mullowney         if (a->compressedrow.use) {
12539ae82921SPaul Mullowney           m    = a->compressedrow.nrows;
12549ae82921SPaul Mullowney           ii   = a->compressedrow.i;
12559ae82921SPaul Mullowney           ridx = a->compressedrow.rindex;
12569ae82921SPaul Mullowney         } else {
1257b06137fdSPaul Mullowney           /* Forcing compressed row on the GPU */
12589ae82921SPaul Mullowney           int k=0;
1259854ce69bSBarry Smith           ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr);
1260854ce69bSBarry Smith           ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr);
12619ae82921SPaul Mullowney           ii[0]=0;
12629ae82921SPaul Mullowney           for (int j = 0; j<m; j++) {
12639ae82921SPaul Mullowney             if ((a->i[j+1]-a->i[j])>0) {
12649ae82921SPaul Mullowney               ii[k]  = a->i[j];
12659ae82921SPaul Mullowney               ridx[k]= j;
12669ae82921SPaul Mullowney               k++;
12679ae82921SPaul Mullowney             }
12689ae82921SPaul Mullowney           }
1269aa372e3fSPaul Mullowney           ii[cusparsestruct->nonzerorow] = a->nz;
1270aa372e3fSPaul Mullowney           m = cusparsestruct->nonzerorow;
12719ae82921SPaul Mullowney         }
12729ae82921SPaul Mullowney 
1273aa372e3fSPaul Mullowney         /* allocate space for the triangular factor information */
1274aa372e3fSPaul Mullowney         matstruct = new Mat_SeqAIJCUSPARSEMultStruct;
1275c41cb2e2SAlejandro Lamas Daviña         stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat);
1276c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
1277c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
12789ae82921SPaul Mullowney 
1279c41cb2e2SAlejandro Lamas Daviña         err = cudaMalloc((void **)&(matstruct->alpha),    sizeof(PetscScalar));CHKERRCUDA(err);
12807656d835SStefano Zampini         err = cudaMalloc((void **)&(matstruct->beta_zero),sizeof(PetscScalar));CHKERRCUDA(err);
12817656d835SStefano Zampini         err = cudaMalloc((void **)&(matstruct->beta_one), sizeof(PetscScalar));CHKERRCUDA(err);
12827656d835SStefano Zampini         err = cudaMemcpy(matstruct->alpha,    &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
12837656d835SStefano Zampini         err = cudaMemcpy(matstruct->beta_zero,&PETSC_CUSPARSE_ZERO,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
12847656d835SStefano Zampini         err = cudaMemcpy(matstruct->beta_one, &PETSC_CUSPARSE_ONE, sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
1285c41cb2e2SAlejandro Lamas Daviña         stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
1286b06137fdSPaul Mullowney 
1287aa372e3fSPaul Mullowney         /* Build a hybrid/ellpack matrix if this option is chosen for the storage */
1288aa372e3fSPaul Mullowney         if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1289aa372e3fSPaul Mullowney           /* set the matrix */
1290aa372e3fSPaul Mullowney           CsrMatrix *matrix= new CsrMatrix;
1291a65300a6SPaul Mullowney           matrix->num_rows = m;
1292aa372e3fSPaul Mullowney           matrix->num_cols = A->cmap->n;
1293aa372e3fSPaul Mullowney           matrix->num_entries = a->nz;
1294a65300a6SPaul Mullowney           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1295a65300a6SPaul Mullowney           matrix->row_offsets->assign(ii, ii + m+1);
12969ae82921SPaul Mullowney 
1297aa372e3fSPaul Mullowney           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1298aa372e3fSPaul Mullowney           matrix->column_indices->assign(a->j, a->j+a->nz);
1299aa372e3fSPaul Mullowney 
1300aa372e3fSPaul Mullowney           matrix->values = new THRUSTARRAY(a->nz);
1301aa372e3fSPaul Mullowney           matrix->values->assign(a->a, a->a+a->nz);
1302aa372e3fSPaul Mullowney 
1303aa372e3fSPaul Mullowney           /* assign the pointer */
1304aa372e3fSPaul Mullowney           matstruct->mat = matrix;
1305aa372e3fSPaul Mullowney 
1306aa372e3fSPaul Mullowney         } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
1307aa372e3fSPaul Mullowney           CsrMatrix *matrix= new CsrMatrix;
1308a65300a6SPaul Mullowney           matrix->num_rows = m;
1309aa372e3fSPaul Mullowney           matrix->num_cols = A->cmap->n;
1310aa372e3fSPaul Mullowney           matrix->num_entries = a->nz;
1311a65300a6SPaul Mullowney           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1312a65300a6SPaul Mullowney           matrix->row_offsets->assign(ii, ii + m+1);
1313aa372e3fSPaul Mullowney 
1314aa372e3fSPaul Mullowney           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1315aa372e3fSPaul Mullowney           matrix->column_indices->assign(a->j, a->j+a->nz);
1316aa372e3fSPaul Mullowney 
1317aa372e3fSPaul Mullowney           matrix->values = new THRUSTARRAY(a->nz);
1318aa372e3fSPaul Mullowney           matrix->values->assign(a->a, a->a+a->nz);
1319aa372e3fSPaul Mullowney 
1320aa372e3fSPaul Mullowney           cusparseHybMat_t hybMat;
1321c41cb2e2SAlejandro Lamas Daviña           stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
1322aa372e3fSPaul Mullowney           cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
1323aa372e3fSPaul Mullowney             CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
1324a65300a6SPaul Mullowney           stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols,
1325aa372e3fSPaul Mullowney               matstruct->descr, matrix->values->data().get(),
1326aa372e3fSPaul Mullowney               matrix->row_offsets->data().get(),
1327aa372e3fSPaul Mullowney               matrix->column_indices->data().get(),
1328c41cb2e2SAlejandro Lamas Daviña               hybMat, 0, partition);CHKERRCUDA(stat);
1329aa372e3fSPaul Mullowney           /* assign the pointer */
1330aa372e3fSPaul Mullowney           matstruct->mat = hybMat;
1331aa372e3fSPaul Mullowney 
1332aa372e3fSPaul Mullowney           if (matrix) {
1333aa372e3fSPaul Mullowney             if (matrix->values) delete (THRUSTARRAY*)matrix->values;
1334aa372e3fSPaul Mullowney             if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices;
1335aa372e3fSPaul Mullowney             if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets;
1336aa372e3fSPaul Mullowney             delete (CsrMatrix*)matrix;
1337087f3262SPaul Mullowney           }
1338087f3262SPaul Mullowney         }
1339ca45077fSPaul Mullowney 
1340aa372e3fSPaul Mullowney         /* assign the compressed row indices */
1341aa372e3fSPaul Mullowney         matstruct->cprowIndices = new THRUSTINTARRAY(m);
1342aa372e3fSPaul Mullowney         matstruct->cprowIndices->assign(ridx,ridx+m);
13434863603aSSatish Balay         ierr = PetscLogCpuToGpu(((m+1)+(a->nz))*sizeof(int)+m*sizeof(PetscInt)+(3+(a->nz))*sizeof(PetscScalar));CHKERRQ(ierr);
1344aa372e3fSPaul Mullowney 
1345aa372e3fSPaul Mullowney         /* assign the pointer */
1346aa372e3fSPaul Mullowney         cusparsestruct->mat = matstruct;
1347aa372e3fSPaul Mullowney 
13489ae82921SPaul Mullowney         if (!a->compressedrow.use) {
13499ae82921SPaul Mullowney           ierr = PetscFree(ii);CHKERRQ(ierr);
13509ae82921SPaul Mullowney           ierr = PetscFree(ridx);CHKERRQ(ierr);
13519ae82921SPaul Mullowney         }
1352e65717acSKarl Rupp         cusparsestruct->workVector = new THRUSTARRAY(m);
13539ae82921SPaul Mullowney       } catch(char *ex) {
13549ae82921SPaul Mullowney         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
13559ae82921SPaul Mullowney       }
135634d6c7a5SJose E. Roman       cusparsestruct->nonzerostate = A->nonzerostate;
135734d6c7a5SJose E. Roman     }
13587d0e52d8SJose E. Roman     ierr = WaitForGPU();CHKERRCUDA(ierr);
1359c70f7ee4SJunchao Zhang     A->offloadmask = PETSC_OFFLOAD_BOTH;
13609ae82921SPaul Mullowney     ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
13619ae82921SPaul Mullowney   }
13629ae82921SPaul Mullowney   PetscFunctionReturn(0);
13639ae82921SPaul Mullowney }
13649ae82921SPaul Mullowney 
1365c41cb2e2SAlejandro Lamas Daviña struct VecCUDAPlusEquals
1366aa372e3fSPaul Mullowney {
1367aa372e3fSPaul Mullowney   template <typename Tuple>
1368aa372e3fSPaul Mullowney   __host__ __device__
1369aa372e3fSPaul Mullowney   void operator()(Tuple t)
1370aa372e3fSPaul Mullowney   {
1371aa372e3fSPaul Mullowney     thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t);
1372aa372e3fSPaul Mullowney   }
1373aa372e3fSPaul Mullowney };
1374aa372e3fSPaul Mullowney 
13756fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
13769ae82921SPaul Mullowney {
1377b175d8bbSPaul Mullowney   PetscErrorCode ierr;
13789ae82921SPaul Mullowney 
13799ae82921SPaul Mullowney   PetscFunctionBegin;
13807656d835SStefano Zampini   ierr = MatMultAdd_SeqAIJCUSPARSE(A,xx,NULL,yy);CHKERRQ(ierr);
13819ae82921SPaul Mullowney   PetscFunctionReturn(0);
13829ae82921SPaul Mullowney }
13839ae82921SPaul Mullowney 
13846fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
1385ca45077fSPaul Mullowney {
1386ca45077fSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1387aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
13889ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct *matstructT;
1389465f34aeSAlejandro Lamas Daviña   const PetscScalar            *xarray;
1390465f34aeSAlejandro Lamas Daviña   PetscScalar                  *yarray;
1391b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1392aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1393ca45077fSPaul Mullowney 
1394ca45077fSPaul Mullowney   PetscFunctionBegin;
139534d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
139634d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
13979ff858a8SKarl Rupp   matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1398aa372e3fSPaul Mullowney   if (!matstructT) {
1399bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1400aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1401bda325fcSPaul Mullowney   }
1402c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1403c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1404f6ae8131SMark   if (yy->map->n) {
1405f6ae8131SMark     PetscInt                     n = yy->map->n;
1406f6ae8131SMark     cudaError_t                  err;
1407f6ae8131SMark     err = cudaMemset(yarray,0,n*sizeof(PetscScalar));CHKERRCUDA(err); /* hack to fix numerical errors from reading output vector yy, apparently */
1408f6ae8131SMark   }
1409aa372e3fSPaul Mullowney 
14107a052e47Shannah_mairs   ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1411aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1412aa372e3fSPaul Mullowney     CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1413aa372e3fSPaul Mullowney     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1414aa372e3fSPaul Mullowney                              mat->num_rows, mat->num_cols,
1415b06137fdSPaul Mullowney                              mat->num_entries, matstructT->alpha, matstructT->descr,
1416aa372e3fSPaul Mullowney                              mat->values->data().get(), mat->row_offsets->data().get(),
14177656d835SStefano Zampini                              mat->column_indices->data().get(), xarray, matstructT->beta_zero,
1418c41cb2e2SAlejandro Lamas Daviña                              yarray);CHKERRCUDA(stat);
1419aa372e3fSPaul Mullowney   } else {
1420aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1421aa372e3fSPaul Mullowney     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1422b06137fdSPaul Mullowney                              matstructT->alpha, matstructT->descr, hybMat,
14237656d835SStefano Zampini                              xarray, matstructT->beta_zero,
1424c41cb2e2SAlejandro Lamas Daviña                              yarray);CHKERRCUDA(stat);
1425ca45077fSPaul Mullowney   }
1426c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1427c41cb2e2SAlejandro Lamas Daviña   ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1428c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1429661c2d29Shannah_mairs   ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1430958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
1431ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1432ca45077fSPaul Mullowney }
1433ca45077fSPaul Mullowney 
1434aa372e3fSPaul Mullowney 
14356fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
14369ae82921SPaul Mullowney {
14379ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1438aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
14399ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct *matstruct;
1440465f34aeSAlejandro Lamas Daviña   const PetscScalar            *xarray;
14417656d835SStefano Zampini   PetscScalar                  *zarray,*dptr,*beta;
1442b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1443aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1444c1fb3f03SStefano Zampini   PetscBool                    cmpr; /* if the matrix has been compressed (zero rows) */
14456e111a19SKarl Rupp 
14469ae82921SPaul Mullowney   PetscFunctionBegin;
144734d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
144834d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
14499ff858a8SKarl Rupp   matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
14509ae82921SPaul Mullowney   try {
1451c1fb3f03SStefano Zampini     cmpr = (PetscBool)(cusparsestruct->workVector->size() == (thrust::detail::vector_base<PetscScalar, thrust::device_malloc_allocator<PetscScalar> >::size_type)(A->rmap->n));
1452c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1453c1fb3f03SStefano Zampini     if (yy && !cmpr) { /* MatMultAdd with noncompressed storage -> need uptodate zz vector */
1454f2d70e9dSBarry Smith       ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr);
1455c1fb3f03SStefano Zampini     } else {
1456c1fb3f03SStefano Zampini       ierr = VecCUDAGetArrayWrite(zz,&zarray);CHKERRQ(ierr);
1457c1fb3f03SStefano Zampini     }
1458c1fb3f03SStefano Zampini     dptr = cmpr ? zarray : cusparsestruct->workVector->data().get();
14597656d835SStefano Zampini     beta = (yy == zz && dptr == zarray) ? matstruct->beta_one : matstruct->beta_zero;
14609ae82921SPaul Mullowney 
14617a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
14627656d835SStefano Zampini     /* csr_spmv is multiply add */
1463aa372e3fSPaul Mullowney     if (cusparsestruct->format == MAT_CUSPARSE_CSR) {
1464b06137fdSPaul Mullowney       /* here we need to be careful to set the number of rows in the multiply to the
1465b06137fdSPaul Mullowney          number of compressed rows in the matrix ... which is equivalent to the
1466b06137fdSPaul Mullowney          size of the workVector */
14677656d835SStefano Zampini       CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1468aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1469a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1470b06137fdSPaul Mullowney                                mat->num_entries, matstruct->alpha, matstruct->descr,
1471aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
14727656d835SStefano Zampini                                mat->column_indices->data().get(), xarray, beta,
14737656d835SStefano Zampini                                dptr);CHKERRCUDA(stat);
1474aa372e3fSPaul Mullowney     } else {
1475a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1476301298b4SMark Adams         cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1477aa372e3fSPaul Mullowney         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1478b06137fdSPaul Mullowney                                  matstruct->alpha, matstruct->descr, hybMat,
14797656d835SStefano Zampini                                  xarray, beta,
14807656d835SStefano Zampini                                  dptr);CHKERRCUDA(stat);
1481a65300a6SPaul Mullowney       }
1482aa372e3fSPaul Mullowney     }
1483958c4211Shannah_mairs     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1484aa372e3fSPaul Mullowney 
14857656d835SStefano Zampini     if (yy) {
14867656d835SStefano Zampini       if (dptr != zarray) {
14877656d835SStefano Zampini         ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
14887656d835SStefano Zampini       } else if (zz != yy) {
14897656d835SStefano Zampini         ierr = VecAXPY_SeqCUDA(zz,1.0,yy);CHKERRQ(ierr);
14907656d835SStefano Zampini       }
14917656d835SStefano Zampini     } else if (dptr != zarray) {
1492c1fb3f03SStefano Zampini       ierr = VecSet_SeqCUDA(zz,0);CHKERRQ(ierr);
14937656d835SStefano Zampini     }
1494aa372e3fSPaul Mullowney     /* scatter the data from the temporary into the full vector with a += operation */
14957a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
14967656d835SStefano Zampini     if (dptr != zarray) {
14977656d835SStefano Zampini       thrust::device_ptr<PetscScalar> zptr;
14987656d835SStefano Zampini 
14997656d835SStefano Zampini       zptr = thrust::device_pointer_cast(zarray);
1500c41cb2e2SAlejandro Lamas Daviña       thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))),
1501c41cb2e2SAlejandro Lamas Daviña                        thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1502c41cb2e2SAlejandro Lamas Daviña                        VecCUDAPlusEquals());
15037656d835SStefano Zampini     }
1504958c4211Shannah_mairs     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1505c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1506c1fb3f03SStefano Zampini     if (yy && !cmpr) {
1507f2d70e9dSBarry Smith       ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr);
1508c1fb3f03SStefano Zampini     } else {
1509c1fb3f03SStefano Zampini       ierr = VecCUDARestoreArrayWrite(zz,&zarray);CHKERRQ(ierr);
1510c1fb3f03SStefano Zampini     }
15119ae82921SPaul Mullowney   } catch(char *ex) {
15129ae82921SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
15139ae82921SPaul Mullowney   }
1514c41cb2e2SAlejandro Lamas Daviña   ierr = WaitForGPU();CHKERRCUDA(ierr);
1515958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz);CHKERRQ(ierr);
15169ae82921SPaul Mullowney   PetscFunctionReturn(0);
15179ae82921SPaul Mullowney }
15189ae82921SPaul Mullowney 
15196fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
1520ca45077fSPaul Mullowney {
1521ca45077fSPaul Mullowney   Mat_SeqAIJ                      *a = (Mat_SeqAIJ*)A->data;
1522aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE              *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
15239ff858a8SKarl Rupp   Mat_SeqAIJCUSPARSEMultStruct    *matstructT;
1524465f34aeSAlejandro Lamas Daviña   const PetscScalar               *xarray;
1525*20291eb5SJunchao Zhang   PetscScalar                     *zarray,*beta;
1526b175d8bbSPaul Mullowney   PetscErrorCode                  ierr;
1527aa372e3fSPaul Mullowney   cusparseStatus_t                stat;
15286e111a19SKarl Rupp 
1529ca45077fSPaul Mullowney   PetscFunctionBegin;
153034d6c7a5SJose E. Roman   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
153134d6c7a5SJose E. Roman   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
15329ff858a8SKarl Rupp   matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1533aa372e3fSPaul Mullowney   if (!matstructT) {
1534bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1535aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1536bda325fcSPaul Mullowney   }
1537aa372e3fSPaul Mullowney 
1538*20291eb5SJunchao Zhang   /* Note unlike Mat, MatTranspose uses non-compressed row storage */
1539ca45077fSPaul Mullowney   try {
1540c41cb2e2SAlejandro Lamas Daviña     ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
1541c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1542f2d70e9dSBarry Smith     ierr = VecCUDAGetArray(zz,&zarray);CHKERRQ(ierr);
1543*20291eb5SJunchao Zhang     beta = (yy == zz) ? matstructT->beta_one : matstructT->beta_zero;
1544ca45077fSPaul Mullowney 
15457a052e47Shannah_mairs     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
1546e057df02SPaul Mullowney     /* multiply add with matrix transpose */
1547aa372e3fSPaul Mullowney     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1548aa372e3fSPaul Mullowney       CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1549aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1550a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1551b06137fdSPaul Mullowney                                mat->num_entries, matstructT->alpha, matstructT->descr,
1552aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
1553a3fdcf43SKarl Rupp                                mat->column_indices->data().get(), xarray, beta,
1554*20291eb5SJunchao Zhang                                zarray);CHKERRCUDA(stat);
1555aa372e3fSPaul Mullowney     } else {
1556aa372e3fSPaul Mullowney       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1557a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1558aa372e3fSPaul Mullowney         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1559b06137fdSPaul Mullowney                                  matstructT->alpha, matstructT->descr, hybMat,
1560a3fdcf43SKarl Rupp                                  xarray, beta,
1561*20291eb5SJunchao Zhang                                  zarray);CHKERRCUDA(stat);
1562a65300a6SPaul Mullowney       }
1563aa372e3fSPaul Mullowney     }
1564a3fdcf43SKarl Rupp     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
1565a3fdcf43SKarl Rupp 
1566*20291eb5SJunchao Zhang     if (zz != yy) {ierr = VecAXPY_SeqCUDA(zz,1.0,yy);CHKERRQ(ierr);}
1567c41cb2e2SAlejandro Lamas Daviña     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1568f2d70e9dSBarry Smith     ierr = VecCUDARestoreArray(zz,&zarray);CHKERRQ(ierr);
1569ca45077fSPaul Mullowney   } catch(char *ex) {
1570ca45077fSPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1571661c2d29Shannah_mairs   }
1572*20291eb5SJunchao Zhang   ierr = WaitForGPU();CHKERRCUDA(ierr);
1573958c4211Shannah_mairs   ierr = PetscLogGpuFlops(2.0*a->nz);CHKERRQ(ierr);
1574ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1575ca45077fSPaul Mullowney }
1576ca45077fSPaul Mullowney 
15776fa9248bSJed Brown static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode)
15789ae82921SPaul Mullowney {
15799ae82921SPaul Mullowney   PetscErrorCode ierr;
15806e111a19SKarl Rupp 
15819ae82921SPaul Mullowney   PetscFunctionBegin;
15829ae82921SPaul Mullowney   ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr);
1583489de41dSStefano Zampini   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
1584bc3f50f2SPaul Mullowney   if (A->factortype == MAT_FACTOR_NONE) {
1585e057df02SPaul Mullowney     ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1586bc3f50f2SPaul Mullowney   }
1587bbf3fe20SPaul Mullowney   A->ops->mult             = MatMult_SeqAIJCUSPARSE;
1588bbf3fe20SPaul Mullowney   A->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1589bbf3fe20SPaul Mullowney   A->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1590bbf3fe20SPaul Mullowney   A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
15919ae82921SPaul Mullowney   PetscFunctionReturn(0);
15929ae82921SPaul Mullowney }
15939ae82921SPaul Mullowney 
15949ae82921SPaul Mullowney /* --------------------------------------------------------------------------------*/
1595e057df02SPaul Mullowney /*@
15969ae82921SPaul Mullowney    MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format
1597e057df02SPaul Mullowney    (the default parallel PETSc format). This matrix will ultimately pushed down
1598e057df02SPaul Mullowney    to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix
1599e057df02SPaul Mullowney    assembly performance the user should preallocate the matrix storage by setting
1600e057df02SPaul Mullowney    the parameter nz (or the array nnz).  By setting these parameters accurately,
1601e057df02SPaul Mullowney    performance during matrix assembly can be increased by more than a factor of 50.
16029ae82921SPaul Mullowney 
1603d083f849SBarry Smith    Collective
16049ae82921SPaul Mullowney 
16059ae82921SPaul Mullowney    Input Parameters:
16069ae82921SPaul Mullowney +  comm - MPI communicator, set to PETSC_COMM_SELF
16079ae82921SPaul Mullowney .  m - number of rows
16089ae82921SPaul Mullowney .  n - number of columns
16099ae82921SPaul Mullowney .  nz - number of nonzeros per row (same for all rows)
16109ae82921SPaul Mullowney -  nnz - array containing the number of nonzeros in the various rows
16110298fd71SBarry Smith          (possibly different for each row) or NULL
16129ae82921SPaul Mullowney 
16139ae82921SPaul Mullowney    Output Parameter:
16149ae82921SPaul Mullowney .  A - the matrix
16159ae82921SPaul Mullowney 
16169ae82921SPaul Mullowney    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
16179ae82921SPaul Mullowney    MatXXXXSetPreallocation() paradgm instead of this routine directly.
16189ae82921SPaul Mullowney    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
16199ae82921SPaul Mullowney 
16209ae82921SPaul Mullowney    Notes:
16219ae82921SPaul Mullowney    If nnz is given then nz is ignored
16229ae82921SPaul Mullowney 
16239ae82921SPaul Mullowney    The AIJ format (also called the Yale sparse matrix format or
16249ae82921SPaul Mullowney    compressed row storage), is fully compatible with standard Fortran 77
16259ae82921SPaul Mullowney    storage.  That is, the stored row and column indices can begin at
16269ae82921SPaul Mullowney    either one (as in Fortran) or zero.  See the users' manual for details.
16279ae82921SPaul Mullowney 
16289ae82921SPaul Mullowney    Specify the preallocated storage with either nz or nnz (not both).
16290298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
16309ae82921SPaul Mullowney    allocation.  For large problems you MUST preallocate memory or you
16319ae82921SPaul Mullowney    will get TERRIBLE performance, see the users' manual chapter on matrices.
16329ae82921SPaul Mullowney 
16339ae82921SPaul Mullowney    By default, this format uses inodes (identical nodes) when possible, to
16349ae82921SPaul Mullowney    improve numerical efficiency of matrix-vector products and solves. We
16359ae82921SPaul Mullowney    search for consecutive rows with the same nonzero structure, thereby
16369ae82921SPaul Mullowney    reusing matrix information to achieve increased efficiency.
16379ae82921SPaul Mullowney 
16389ae82921SPaul Mullowney    Level: intermediate
16399ae82921SPaul Mullowney 
1640e057df02SPaul Mullowney .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE
16419ae82921SPaul Mullowney @*/
16429ae82921SPaul Mullowney PetscErrorCode  MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
16439ae82921SPaul Mullowney {
16449ae82921SPaul Mullowney   PetscErrorCode ierr;
16459ae82921SPaul Mullowney 
16469ae82921SPaul Mullowney   PetscFunctionBegin;
16479ae82921SPaul Mullowney   ierr = MatCreate(comm,A);CHKERRQ(ierr);
16489ae82921SPaul Mullowney   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
16499ae82921SPaul Mullowney   ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
16509ae82921SPaul Mullowney   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
16519ae82921SPaul Mullowney   PetscFunctionReturn(0);
16529ae82921SPaul Mullowney }
16539ae82921SPaul Mullowney 
16546fa9248bSJed Brown static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A)
16559ae82921SPaul Mullowney {
16569ae82921SPaul Mullowney   PetscErrorCode   ierr;
1657ab25e6cbSDominic Meiser 
16589ae82921SPaul Mullowney   PetscFunctionBegin;
16599ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
1660c70f7ee4SJunchao Zhang     if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) {
1661470880abSPatrick Sanan       ierr = MatSeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr);
16629ae82921SPaul Mullowney     }
16639ae82921SPaul Mullowney   } else {
1664470880abSPatrick Sanan     ierr = MatSeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr);
1665aa372e3fSPaul Mullowney   }
16669ae82921SPaul Mullowney   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
16679ae82921SPaul Mullowney   PetscFunctionReturn(0);
16689ae82921SPaul Mullowney }
16699ae82921SPaul Mullowney 
16709ff858a8SKarl Rupp static PetscErrorCode MatDuplicate_SeqAIJCUSPARSE(Mat A,MatDuplicateOption cpvalues,Mat *B)
16719ff858a8SKarl Rupp {
16729ff858a8SKarl Rupp   PetscErrorCode ierr;
16739ff858a8SKarl Rupp   Mat C;
16749ff858a8SKarl Rupp   cusparseStatus_t stat;
16759ff858a8SKarl Rupp   cusparseHandle_t handle=0;
16769ff858a8SKarl Rupp 
16779ff858a8SKarl Rupp   PetscFunctionBegin;
16789ff858a8SKarl Rupp   ierr = MatDuplicate_SeqAIJ(A,cpvalues,B);CHKERRQ(ierr);
16799ff858a8SKarl Rupp   C    = *B;
168034136279SStefano Zampini   ierr = PetscFree(C->defaultvectype);CHKERRQ(ierr);
168134136279SStefano Zampini   ierr = PetscStrallocpy(VECCUDA,&C->defaultvectype);CHKERRQ(ierr);
168234136279SStefano Zampini 
16839ff858a8SKarl Rupp   /* inject CUSPARSE-specific stuff */
16849ff858a8SKarl Rupp   if (C->factortype==MAT_FACTOR_NONE) {
16859ff858a8SKarl Rupp     /* you cannot check the inode.use flag here since the matrix was just created.
16869ff858a8SKarl Rupp        now build a GPU matrix data structure */
16879ff858a8SKarl Rupp     C->spptr = new Mat_SeqAIJCUSPARSE;
16889ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->mat            = 0;
16899ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->matTranspose   = 0;
16909ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->workVector     = 0;
169181902715SJunchao Zhang     ((Mat_SeqAIJCUSPARSE*)C->spptr)->rowoffsets_gpu = 0;
16929ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->format         = MAT_CUSPARSE_CSR;
16939ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->stream         = 0;
16949ff858a8SKarl Rupp     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
16959ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->handle         = handle;
16969ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)C->spptr)->nonzerostate   = 0;
16979ff858a8SKarl Rupp   } else {
16989ff858a8SKarl Rupp     /* NEXT, set the pointers to the triangular factors */
16999ff858a8SKarl Rupp     C->spptr = new Mat_SeqAIJCUSPARSETriFactors;
17009ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtr          = 0;
17019ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtr          = 0;
17029ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->loTriFactorPtrTranspose = 0;
17039ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->upTriFactorPtrTranspose = 0;
17049ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->rpermIndices            = 0;
17059ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->cpermIndices            = 0;
17069ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->workVector              = 0;
17079ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle                  = 0;
17089ff858a8SKarl Rupp     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
17099ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->handle                  = handle;
17109ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSETriFactors*)C->spptr)->nnz                     = 0;
17119ff858a8SKarl Rupp   }
17129ff858a8SKarl Rupp 
17139ff858a8SKarl Rupp   C->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
17149ff858a8SKarl Rupp   C->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
17159ff858a8SKarl Rupp   C->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
17169ff858a8SKarl Rupp   C->ops->mult             = MatMult_SeqAIJCUSPARSE;
17179ff858a8SKarl Rupp   C->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
17189ff858a8SKarl Rupp   C->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
17199ff858a8SKarl Rupp   C->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
17209ff858a8SKarl Rupp   C->ops->duplicate        = MatDuplicate_SeqAIJCUSPARSE;
17219ff858a8SKarl Rupp 
17229ff858a8SKarl Rupp   ierr = PetscObjectChangeTypeName((PetscObject)C,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
17239ff858a8SKarl Rupp 
1724c70f7ee4SJunchao Zhang   C->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
17259ff858a8SKarl Rupp 
17269ff858a8SKarl Rupp   ierr = PetscObjectComposeFunction((PetscObject)C, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
17279ff858a8SKarl Rupp   PetscFunctionReturn(0);
17289ff858a8SKarl Rupp }
17299ff858a8SKarl Rupp 
173002fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat B)
17319ae82921SPaul Mullowney {
17329ae82921SPaul Mullowney   PetscErrorCode ierr;
1733aa372e3fSPaul Mullowney   cusparseStatus_t stat;
1734aa372e3fSPaul Mullowney   cusparseHandle_t handle=0;
17359ae82921SPaul Mullowney 
17369ae82921SPaul Mullowney   PetscFunctionBegin;
173734136279SStefano Zampini   ierr = PetscFree(B->defaultvectype);CHKERRQ(ierr);
173834136279SStefano Zampini   ierr = PetscStrallocpy(VECCUDA,&B->defaultvectype);CHKERRQ(ierr);
173934136279SStefano Zampini 
17409ae82921SPaul Mullowney   if (B->factortype==MAT_FACTOR_NONE) {
1741e057df02SPaul Mullowney     /* you cannot check the inode.use flag here since the matrix was just created.
1742e057df02SPaul Mullowney        now build a GPU matrix data structure */
17439ae82921SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSE;
17449ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat            = 0;
1745aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose   = 0;
1746aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector     = 0;
174781902715SJunchao Zhang     ((Mat_SeqAIJCUSPARSE*)B->spptr)->rowoffsets_gpu = 0;
1748e057df02SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->format         = MAT_CUSPARSE_CSR;
1749aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream         = 0;
1750c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1751aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle         = handle;
17529ff858a8SKarl Rupp     ((Mat_SeqAIJCUSPARSE*)B->spptr)->nonzerostate   = 0;
17539ae82921SPaul Mullowney   } else {
17549ae82921SPaul Mullowney     /* NEXT, set the pointers to the triangular factors */
1755debe9ee2SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSETriFactors;
17569ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr          = 0;
17579ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr          = 0;
1758aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0;
1759aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0;
1760aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices            = 0;
1761aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices            = 0;
1762aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector              = 0;
1763c41cb2e2SAlejandro Lamas Daviña     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1764aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = handle;
1765aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz                     = 0;
17669ae82921SPaul Mullowney   }
1767aa372e3fSPaul Mullowney 
17689ae82921SPaul Mullowney   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
17699ae82921SPaul Mullowney   B->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
17709ae82921SPaul Mullowney   B->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
1771ca45077fSPaul Mullowney   B->ops->mult             = MatMult_SeqAIJCUSPARSE;
1772ca45077fSPaul Mullowney   B->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1773ca45077fSPaul Mullowney   B->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1774ca45077fSPaul Mullowney   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
17759ff858a8SKarl Rupp   B->ops->duplicate        = MatDuplicate_SeqAIJCUSPARSE;
17762205254eSKarl Rupp 
17779ae82921SPaul Mullowney   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
17782205254eSKarl Rupp 
1779c70f7ee4SJunchao Zhang   B->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
17802205254eSKarl Rupp 
1781bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
17829ae82921SPaul Mullowney   PetscFunctionReturn(0);
17839ae82921SPaul Mullowney }
17849ae82921SPaul Mullowney 
178502fe1965SBarry Smith PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B)
178602fe1965SBarry Smith {
178702fe1965SBarry Smith   PetscErrorCode ierr;
178802fe1965SBarry Smith 
178902fe1965SBarry Smith   PetscFunctionBegin;
179002fe1965SBarry Smith   ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr);
1791489de41dSStefano Zampini   ierr = MatConvert_SeqAIJ_SeqAIJCUSPARSE(B);CHKERRQ(ierr);
179202fe1965SBarry Smith   PetscFunctionReturn(0);
179302fe1965SBarry Smith }
179402fe1965SBarry Smith 
17953ca39a21SBarry Smith /*MC
1796e057df02SPaul Mullowney    MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices.
1797e057df02SPaul Mullowney 
1798e057df02SPaul Mullowney    A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either
17992692e278SPaul Mullowney    CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later.
18002692e278SPaul Mullowney    All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library.
1801e057df02SPaul Mullowney 
1802e057df02SPaul Mullowney    Options Database Keys:
1803e057df02SPaul Mullowney +  -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions()
1804aa372e3fSPaul 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).
1805a2b725a8SWilliam 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).
1806e057df02SPaul Mullowney 
1807e057df02SPaul Mullowney   Level: beginner
1808e057df02SPaul Mullowney 
18098468deeeSKarl Rupp .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
1810e057df02SPaul Mullowney M*/
18117f756511SDominic Meiser 
181242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*);
181342c9c57cSBarry Smith 
18140f39cd5aSBarry Smith 
18153ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_CUSPARSE(void)
181642c9c57cSBarry Smith {
181742c9c57cSBarry Smith   PetscErrorCode ierr;
181842c9c57cSBarry Smith 
181942c9c57cSBarry Smith   PetscFunctionBegin;
18203ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
18213ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
18223ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
18233ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
182442c9c57cSBarry Smith   PetscFunctionReturn(0);
182542c9c57cSBarry Smith }
182629b38603SBarry Smith 
182781e08676SBarry Smith 
1828470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct)
18297f756511SDominic Meiser {
18307f756511SDominic Meiser   cusparseStatus_t stat;
18317f756511SDominic Meiser   cusparseHandle_t handle;
18327f756511SDominic Meiser 
18337f756511SDominic Meiser   PetscFunctionBegin;
18347f756511SDominic Meiser   if (*cusparsestruct) {
1835470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format);
1836470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format);
18377f756511SDominic Meiser     delete (*cusparsestruct)->workVector;
183881902715SJunchao Zhang     delete (*cusparsestruct)->rowoffsets_gpu;
18397f756511SDominic Meiser     if (handle = (*cusparsestruct)->handle) {
1840c41cb2e2SAlejandro Lamas Daviña       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
18417f756511SDominic Meiser     }
18427f756511SDominic Meiser     delete *cusparsestruct;
18437f756511SDominic Meiser     *cusparsestruct = 0;
18447f756511SDominic Meiser   }
18457f756511SDominic Meiser   PetscFunctionReturn(0);
18467f756511SDominic Meiser }
18477f756511SDominic Meiser 
18487f756511SDominic Meiser static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat)
18497f756511SDominic Meiser {
18507f756511SDominic Meiser   PetscFunctionBegin;
18517f756511SDominic Meiser   if (*mat) {
18527f756511SDominic Meiser     delete (*mat)->values;
18537f756511SDominic Meiser     delete (*mat)->column_indices;
18547f756511SDominic Meiser     delete (*mat)->row_offsets;
18557f756511SDominic Meiser     delete *mat;
18567f756511SDominic Meiser     *mat = 0;
18577f756511SDominic Meiser   }
18587f756511SDominic Meiser   PetscFunctionReturn(0);
18597f756511SDominic Meiser }
18607f756511SDominic Meiser 
1861470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor)
18627f756511SDominic Meiser {
18637f756511SDominic Meiser   cusparseStatus_t stat;
18647f756511SDominic Meiser   PetscErrorCode   ierr;
18657f756511SDominic Meiser 
18667f756511SDominic Meiser   PetscFunctionBegin;
18677f756511SDominic Meiser   if (*trifactor) {
1868c41cb2e2SAlejandro Lamas Daviña     if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); }
1869c41cb2e2SAlejandro Lamas Daviña     if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); }
18707f756511SDominic Meiser     ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr);
18717f756511SDominic Meiser     delete *trifactor;
18727f756511SDominic Meiser     *trifactor = 0;
18737f756511SDominic Meiser   }
18747f756511SDominic Meiser   PetscFunctionReturn(0);
18757f756511SDominic Meiser }
18767f756511SDominic Meiser 
1877470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format)
18787f756511SDominic Meiser {
18797f756511SDominic Meiser   CsrMatrix        *mat;
18807f756511SDominic Meiser   cusparseStatus_t stat;
18817f756511SDominic Meiser   cudaError_t      err;
18827f756511SDominic Meiser 
18837f756511SDominic Meiser   PetscFunctionBegin;
18847f756511SDominic Meiser   if (*matstruct) {
18857f756511SDominic Meiser     if ((*matstruct)->mat) {
18867f756511SDominic Meiser       if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) {
18877f756511SDominic Meiser         cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat;
1888c41cb2e2SAlejandro Lamas Daviña         stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat);
18897f756511SDominic Meiser       } else {
18907f756511SDominic Meiser         mat = (CsrMatrix*)(*matstruct)->mat;
18917f756511SDominic Meiser         CsrMatrix_Destroy(&mat);
18927f756511SDominic Meiser       }
18937f756511SDominic Meiser     }
1894c41cb2e2SAlejandro Lamas Daviña     if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); }
18957f756511SDominic Meiser     delete (*matstruct)->cprowIndices;
1896c41cb2e2SAlejandro Lamas Daviña     if ((*matstruct)->alpha)     { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); }
18977656d835SStefano Zampini     if ((*matstruct)->beta_zero) { err=cudaFree((*matstruct)->beta_zero);CHKERRCUDA(err); }
18987656d835SStefano Zampini     if ((*matstruct)->beta_one)  { err=cudaFree((*matstruct)->beta_one);CHKERRCUDA(err); }
18997f756511SDominic Meiser     delete *matstruct;
19007f756511SDominic Meiser     *matstruct = 0;
19017f756511SDominic Meiser   }
19027f756511SDominic Meiser   PetscFunctionReturn(0);
19037f756511SDominic Meiser }
19047f756511SDominic Meiser 
1905470880abSPatrick Sanan static PetscErrorCode MatSeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors)
19067f756511SDominic Meiser {
19077f756511SDominic Meiser   cusparseHandle_t handle;
19087f756511SDominic Meiser   cusparseStatus_t stat;
19097f756511SDominic Meiser 
19107f756511SDominic Meiser   PetscFunctionBegin;
19117f756511SDominic Meiser   if (*trifactors) {
1912470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtr);
1913470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtr);
1914470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose);
1915470880abSPatrick Sanan     MatSeqAIJCUSPARSEMultStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose);
19167f756511SDominic Meiser     delete (*trifactors)->rpermIndices;
19177f756511SDominic Meiser     delete (*trifactors)->cpermIndices;
19187f756511SDominic Meiser     delete (*trifactors)->workVector;
19197f756511SDominic Meiser     if (handle = (*trifactors)->handle) {
1920c41cb2e2SAlejandro Lamas Daviña       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
19217f756511SDominic Meiser     }
19227f756511SDominic Meiser     delete *trifactors;
19237f756511SDominic Meiser     *trifactors = 0;
19247f756511SDominic Meiser   }
19257f756511SDominic Meiser   PetscFunctionReturn(0);
19267f756511SDominic Meiser }
19277f756511SDominic Meiser 
1928