xref: /petsc/src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu (revision ebc8f436579ae4907d6cb7dbf060cef5ed1ab4cb)
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 */
59ae82921SPaul Mullowney 
63d13b8fdSMatthew G. Knepley #include <petscconf.h>
73d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/aij.h>          /*I "petscmat.h" I*/
8087f3262SPaul Mullowney #include <../src/mat/impls/sbaij/seq/sbaij.h>
93d13b8fdSMatthew G. Knepley #include <../src/vec/vec/impls/dvecimpl.h>
103d13b8fdSMatthew G. Knepley #include <petsc-private/vecimpl.h>
119ae82921SPaul Mullowney #undef VecType
123d13b8fdSMatthew G. Knepley #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h>
13bc3f50f2SPaul Mullowney 
14e057df02SPaul Mullowney const char *const MatCUSPARSEStorageFormats[] = {"CSR","ELL","HYB","MatCUSPARSEStorageFormat","MAT_CUSPARSE_",0};
159ae82921SPaul Mullowney 
16087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
17087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
18087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
19087f3262SPaul Mullowney 
206fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
216fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
226fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
23087f3262SPaul Mullowney 
246fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat,Vec,Vec);
256fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
266fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
276fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
286fa9248bSJed Brown static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(Mat);
296fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat,Vec,Vec);
306fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
316fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
326fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
339ae82921SPaul Mullowney 
349ae82921SPaul Mullowney #undef __FUNCT__
35b06137fdSPaul Mullowney #define __FUNCT__ "MatCUSPARSESetStream"
36b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetStream(Mat A,const cudaStream_t stream)
37b06137fdSPaul Mullowney {
38b06137fdSPaul Mullowney   cusparseStatus_t   stat;
39b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
40b06137fdSPaul Mullowney 
41b06137fdSPaul Mullowney   PetscFunctionBegin;
42b06137fdSPaul Mullowney   cusparsestruct->stream = stream;
43b06137fdSPaul Mullowney   stat = cusparseSetStream(cusparsestruct->handle,cusparsestruct->stream);CHKERRCUSP(stat);
44b06137fdSPaul Mullowney   PetscFunctionReturn(0);
45b06137fdSPaul Mullowney }
46b06137fdSPaul Mullowney 
47b06137fdSPaul Mullowney #undef __FUNCT__
48b06137fdSPaul Mullowney #define __FUNCT__ "MatCUSPARSESetHandle"
49b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSESetHandle(Mat A,const cusparseHandle_t handle)
50b06137fdSPaul Mullowney {
51b06137fdSPaul Mullowney   cusparseStatus_t   stat;
52b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
53b06137fdSPaul Mullowney 
54b06137fdSPaul Mullowney   PetscFunctionBegin;
55b06137fdSPaul Mullowney   if (cusparsestruct->handle)
56b06137fdSPaul Mullowney     stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUSP(stat);
57b06137fdSPaul Mullowney   cusparsestruct->handle = handle;
58b06137fdSPaul Mullowney   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUSP(stat);
59b06137fdSPaul Mullowney   PetscFunctionReturn(0);
60b06137fdSPaul Mullowney }
61b06137fdSPaul Mullowney 
62b06137fdSPaul Mullowney #undef __FUNCT__
63b06137fdSPaul Mullowney #define __FUNCT__ "MatCUSPARSEClearHandle"
64b06137fdSPaul Mullowney PetscErrorCode MatCUSPARSEClearHandle(Mat A)
65b06137fdSPaul Mullowney {
66b06137fdSPaul Mullowney   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
67b06137fdSPaul Mullowney   PetscFunctionBegin;
68b06137fdSPaul Mullowney   if (cusparsestruct->handle)
69b06137fdSPaul Mullowney     cusparsestruct->handle = 0;
70b06137fdSPaul Mullowney   PetscFunctionReturn(0);
71b06137fdSPaul Mullowney }
72b06137fdSPaul Mullowney 
73b06137fdSPaul Mullowney #undef __FUNCT__
749ae82921SPaul Mullowney #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_cusparse"
759ae82921SPaul Mullowney PetscErrorCode MatFactorGetSolverPackage_seqaij_cusparse(Mat A,const MatSolverPackage *type)
769ae82921SPaul Mullowney {
779ae82921SPaul Mullowney   PetscFunctionBegin;
789ae82921SPaul Mullowney   *type = MATSOLVERCUSPARSE;
799ae82921SPaul Mullowney   PetscFunctionReturn(0);
809ae82921SPaul Mullowney }
819ae82921SPaul Mullowney 
82c708e6cdSJed Brown /*MC
83087f3262SPaul Mullowney   MATSOLVERCUSPARSE = "cusparse" - A matrix type providing triangular solvers for seq matrices
84087f3262SPaul Mullowney   on a single GPU of type, seqaijcusparse, aijcusparse, or seqaijcusp, aijcusp. Currently supported
85087f3262SPaul Mullowney   algorithms are ILU(k) and ICC(k). Typically, deeper factorizations (larger k) results in poorer
86087f3262SPaul Mullowney   performance in the triangular solves. Full LU, and Cholesky decompositions can be solved through the
87087f3262SPaul Mullowney   CUSPARSE triangular solve algorithm. However, the performance can be quite poor and thus these
88087f3262SPaul Mullowney   algorithms are not recommended. This class does NOT support direct solver operations.
89c708e6cdSJed Brown 
909ae82921SPaul Mullowney   Level: beginner
91c708e6cdSJed Brown 
92c708e6cdSJed Brown .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
93c708e6cdSJed Brown M*/
949ae82921SPaul Mullowney 
959ae82921SPaul Mullowney #undef __FUNCT__
969ae82921SPaul Mullowney #define __FUNCT__ "MatGetFactor_seqaij_cusparse"
973c3a0fd4SSatish Balay PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_cusparse(Mat A,MatFactorType ftype,Mat *B)
989ae82921SPaul Mullowney {
999ae82921SPaul Mullowney   PetscErrorCode ierr;
100bc3f50f2SPaul Mullowney   PetscInt       n = A->rmap->n;
1019ae82921SPaul Mullowney 
1029ae82921SPaul Mullowney   PetscFunctionBegin;
103bc3f50f2SPaul Mullowney   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
104404133a2SPaul Mullowney   (*B)->factortype = ftype;
105bc3f50f2SPaul Mullowney   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
1069ae82921SPaul Mullowney   ierr = MatSetType(*B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
1072205254eSKarl Rupp 
108087f3262SPaul Mullowney   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
10933d57670SJed Brown     ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr);
1109ae82921SPaul Mullowney     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJCUSPARSE;
1119ae82921SPaul Mullowney     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJCUSPARSE;
112087f3262SPaul Mullowney   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
113087f3262SPaul Mullowney     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJCUSPARSE;
114087f3262SPaul Mullowney     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJCUSPARSE;
1159ae82921SPaul Mullowney   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported for CUSPARSE Matrix Types");
116bc3f50f2SPaul Mullowney 
117fa03d054SJed Brown   ierr = MatSeqAIJSetPreallocation(*B,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
11862a20339SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)(*B),"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_cusparse);CHKERRQ(ierr);
1199ae82921SPaul Mullowney   PetscFunctionReturn(0);
1209ae82921SPaul Mullowney }
1219ae82921SPaul Mullowney 
1229ae82921SPaul Mullowney #undef __FUNCT__
123e057df02SPaul Mullowney #define __FUNCT__ "MatCUSPARSESetFormat_SeqAIJCUSPARSE"
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;
1292692e278SPaul Mullowney #if CUDA_VERSION>=4020
130ca45077fSPaul Mullowney   switch (op) {
131e057df02SPaul Mullowney   case MAT_CUSPARSE_MULT:
132aa372e3fSPaul Mullowney     cusparsestruct->format = format;
133ca45077fSPaul Mullowney     break;
134e057df02SPaul Mullowney   case MAT_CUSPARSE_ALL:
135aa372e3fSPaul Mullowney     cusparsestruct->format = format;
136ca45077fSPaul Mullowney     break;
137ca45077fSPaul Mullowney   default:
13836d62e41SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op);
139ca45077fSPaul Mullowney   }
1402692e278SPaul Mullowney #else
1412692e278SPaul Mullowney   if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB)
1422692e278SPaul Mullowney     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ELL (Ellpack) and HYB (Hybrid) storage format require CUDA 4.2 or later.");
1432692e278SPaul Mullowney #endif
144ca45077fSPaul Mullowney   PetscFunctionReturn(0);
145ca45077fSPaul Mullowney }
1469ae82921SPaul Mullowney 
147e057df02SPaul Mullowney /*@
148e057df02SPaul Mullowney    MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular
149e057df02SPaul Mullowney    operation. Only the MatMult operation can use different GPU storage formats
150aa372e3fSPaul Mullowney    for MPIAIJCUSPARSE matrices.
151e057df02SPaul Mullowney    Not Collective
152e057df02SPaul Mullowney 
153e057df02SPaul Mullowney    Input Parameters:
1548468deeeSKarl Rupp +  A - Matrix of type SEQAIJCUSPARSE
15536d62e41SPaul 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.
1562692e278SPaul Mullowney -  format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2)
157e057df02SPaul Mullowney 
158e057df02SPaul Mullowney    Output Parameter:
159e057df02SPaul Mullowney 
160e057df02SPaul Mullowney    Level: intermediate
161e057df02SPaul Mullowney 
1628468deeeSKarl Rupp .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
163e057df02SPaul Mullowney @*/
164e057df02SPaul Mullowney #undef __FUNCT__
165e057df02SPaul Mullowney #define __FUNCT__ "MatCUSPARSESetFormat"
166e057df02SPaul Mullowney PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
167e057df02SPaul Mullowney {
168e057df02SPaul Mullowney   PetscErrorCode ierr;
1696e111a19SKarl Rupp 
170e057df02SPaul Mullowney   PetscFunctionBegin;
171e057df02SPaul Mullowney   PetscValidHeaderSpecific(A, MAT_CLASSID,1);
172e057df02SPaul Mullowney   ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr);
173e057df02SPaul Mullowney   PetscFunctionReturn(0);
174e057df02SPaul Mullowney }
175e057df02SPaul Mullowney 
1769ae82921SPaul Mullowney #undef __FUNCT__
1779ae82921SPaul Mullowney #define __FUNCT__ "MatSetFromOptions_SeqAIJCUSPARSE"
1786fa9248bSJed Brown static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(Mat A)
1799ae82921SPaul Mullowney {
1809ae82921SPaul Mullowney   PetscErrorCode           ierr;
181e057df02SPaul Mullowney   MatCUSPARSEStorageFormat format;
1829ae82921SPaul Mullowney   PetscBool                flg;
1836e111a19SKarl Rupp 
1849ae82921SPaul Mullowney   PetscFunctionBegin;
185e057df02SPaul Mullowney   ierr = PetscOptionsHead("SeqAIJCUSPARSE options");CHKERRQ(ierr);
186e057df02SPaul Mullowney   ierr = PetscObjectOptionsBegin((PetscObject)A);
1879ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
188e057df02SPaul Mullowney     ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV",
1897083f85cSSatish Balay                             "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)MAT_CUSPARSE_CSR,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
190e057df02SPaul Mullowney     if (flg) {
191e057df02SPaul Mullowney       ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr);
192045c96e1SPaul Mullowney     }
1939ae82921SPaul Mullowney   }
1944c87dfd4SPaul Mullowney   ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve",
1957083f85cSSatish Balay                           "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)MAT_CUSPARSE_CSR,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
1964c87dfd4SPaul Mullowney   if (flg) {
1974c87dfd4SPaul Mullowney     ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr);
1984c87dfd4SPaul Mullowney   }
1999ae82921SPaul Mullowney   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2009ae82921SPaul Mullowney   PetscFunctionReturn(0);
2019ae82921SPaul Mullowney 
2029ae82921SPaul Mullowney }
2039ae82921SPaul Mullowney 
2049ae82921SPaul Mullowney #undef __FUNCT__
2059ae82921SPaul Mullowney #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJCUSPARSE"
2066fa9248bSJed Brown static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2079ae82921SPaul Mullowney {
2089ae82921SPaul Mullowney   PetscErrorCode ierr;
2099ae82921SPaul Mullowney 
2109ae82921SPaul Mullowney   PetscFunctionBegin;
2119ae82921SPaul Mullowney   ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
2129ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2139ae82921SPaul Mullowney   PetscFunctionReturn(0);
2149ae82921SPaul Mullowney }
2159ae82921SPaul Mullowney 
2169ae82921SPaul Mullowney #undef __FUNCT__
2179ae82921SPaul Mullowney #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJCUSPARSE"
2186fa9248bSJed Brown static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2199ae82921SPaul Mullowney {
2209ae82921SPaul Mullowney   PetscErrorCode ierr;
2219ae82921SPaul Mullowney 
2229ae82921SPaul Mullowney   PetscFunctionBegin;
2239ae82921SPaul Mullowney   ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
2249ae82921SPaul Mullowney   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
2259ae82921SPaul Mullowney   PetscFunctionReturn(0);
2269ae82921SPaul Mullowney }
2279ae82921SPaul Mullowney 
2289ae82921SPaul Mullowney #undef __FUNCT__
229087f3262SPaul Mullowney #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJCUSPARSE"
230087f3262SPaul Mullowney static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
231087f3262SPaul Mullowney {
232087f3262SPaul Mullowney   PetscErrorCode ierr;
233087f3262SPaul Mullowney 
234087f3262SPaul Mullowney   PetscFunctionBegin;
235087f3262SPaul Mullowney   ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
236087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
237087f3262SPaul Mullowney   PetscFunctionReturn(0);
238087f3262SPaul Mullowney }
239087f3262SPaul Mullowney 
240087f3262SPaul Mullowney #undef __FUNCT__
241087f3262SPaul Mullowney #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJCUSPARSE"
242087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
243087f3262SPaul Mullowney {
244087f3262SPaul Mullowney   PetscErrorCode ierr;
245087f3262SPaul Mullowney 
246087f3262SPaul Mullowney   PetscFunctionBegin;
247087f3262SPaul Mullowney   ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
248087f3262SPaul Mullowney   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
249087f3262SPaul Mullowney   PetscFunctionReturn(0);
250087f3262SPaul Mullowney }
251087f3262SPaul Mullowney 
252087f3262SPaul Mullowney #undef __FUNCT__
253087f3262SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEBuildILULowerTriMatrix"
254087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A)
2559ae82921SPaul Mullowney {
2569ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
2579ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
2589ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
259aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
2609ae82921SPaul Mullowney   cusparseStatus_t                  stat;
2619ae82921SPaul Mullowney   const PetscInt                    *ai = a->i,*aj = a->j,*vi;
2629ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
2639ae82921SPaul Mullowney   PetscInt                          *AiLo, *AjLo;
2649ae82921SPaul Mullowney   PetscScalar                       *AALo;
2659ae82921SPaul Mullowney   PetscInt                          i,nz, nzLower, offset, rowOffset;
266b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
2679ae82921SPaul Mullowney 
2689ae82921SPaul Mullowney   PetscFunctionBegin;
2699ae82921SPaul Mullowney   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUSP_CPU) {
2709ae82921SPaul Mullowney     try {
2719ae82921SPaul Mullowney       /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */
2729ae82921SPaul Mullowney       nzLower=n+ai[n]-ai[1];
2739ae82921SPaul Mullowney 
2749ae82921SPaul Mullowney       /* Allocate Space for the lower triangular matrix */
2759ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUSP(ierr);
2769ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUSP(ierr);
2779ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUSP(ierr);
2789ae82921SPaul Mullowney 
2799ae82921SPaul Mullowney       /* Fill the lower triangular matrix */
2809ae82921SPaul Mullowney       AiLo[0]  = (PetscInt) 0;
2819ae82921SPaul Mullowney       AiLo[n]  = nzLower;
2829ae82921SPaul Mullowney       AjLo[0]  = (PetscInt) 0;
2839ae82921SPaul Mullowney       AALo[0]  = (MatScalar) 1.0;
2849ae82921SPaul Mullowney       v        = aa;
2859ae82921SPaul Mullowney       vi       = aj;
2869ae82921SPaul Mullowney       offset   = 1;
2879ae82921SPaul Mullowney       rowOffset= 1;
2889ae82921SPaul Mullowney       for (i=1; i<n; i++) {
2899ae82921SPaul Mullowney         nz = ai[i+1] - ai[i];
290e057df02SPaul Mullowney         /* additional 1 for the term on the diagonal */
2919ae82921SPaul Mullowney         AiLo[i]    = rowOffset;
2929ae82921SPaul Mullowney         rowOffset += nz+1;
2939ae82921SPaul Mullowney 
2949ae82921SPaul Mullowney         ierr = PetscMemcpy(&(AjLo[offset]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr);
2959ae82921SPaul Mullowney         ierr = PetscMemcpy(&(AALo[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
2969ae82921SPaul Mullowney 
2979ae82921SPaul Mullowney         offset      += nz;
2989ae82921SPaul Mullowney         AjLo[offset] = (PetscInt) i;
2999ae82921SPaul Mullowney         AALo[offset] = (MatScalar) 1.0;
3009ae82921SPaul Mullowney         offset      += 1;
3019ae82921SPaul Mullowney 
3029ae82921SPaul Mullowney         v  += nz;
3039ae82921SPaul Mullowney         vi += nz;
3049ae82921SPaul Mullowney       }
3052205254eSKarl Rupp 
306aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
307aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
3082205254eSKarl Rupp 
309aa372e3fSPaul Mullowney       /* Create the matrix description */
310aa372e3fSPaul Mullowney       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUSP(stat);
311aa372e3fSPaul Mullowney       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUSP(stat);
312aa372e3fSPaul Mullowney       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUSP(stat);
313aa372e3fSPaul Mullowney       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUSP(stat);
314aa372e3fSPaul Mullowney       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUSP(stat);
315aa372e3fSPaul Mullowney 
316aa372e3fSPaul Mullowney       /* Create the solve analysis information */
317aa372e3fSPaul Mullowney       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUSP(stat);
318aa372e3fSPaul Mullowney 
319aa372e3fSPaul Mullowney       /* set the operation */
320aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
321aa372e3fSPaul Mullowney 
322aa372e3fSPaul Mullowney       /* set the matrix */
323aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
324aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = n;
325aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = n;
326aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = nzLower;
327aa372e3fSPaul Mullowney 
328aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
329aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1);
330aa372e3fSPaul Mullowney 
331aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower);
332aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower);
333aa372e3fSPaul Mullowney 
334aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(nzLower);
335aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+nzLower);
336aa372e3fSPaul Mullowney 
337aa372e3fSPaul Mullowney       /* perform the solve analysis */
338aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
339aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
340aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
341aa372e3fSPaul Mullowney                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUSP(stat);
342aa372e3fSPaul Mullowney 
343aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
344aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
3452205254eSKarl Rupp 
3469ae82921SPaul Mullowney       ierr = cudaFreeHost(AiLo);CHKERRCUSP(ierr);
3479ae82921SPaul Mullowney       ierr = cudaFreeHost(AjLo);CHKERRCUSP(ierr);
3489ae82921SPaul Mullowney       ierr = cudaFreeHost(AALo);CHKERRCUSP(ierr);
3499ae82921SPaul Mullowney     } catch(char *ex) {
3509ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
3519ae82921SPaul Mullowney     }
3529ae82921SPaul Mullowney   }
3539ae82921SPaul Mullowney   PetscFunctionReturn(0);
3549ae82921SPaul Mullowney }
3559ae82921SPaul Mullowney 
3569ae82921SPaul Mullowney #undef __FUNCT__
357087f3262SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEBuildILUUpperTriMatrix"
358087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A)
3599ae82921SPaul Mullowney {
3609ae82921SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
3619ae82921SPaul Mullowney   PetscInt                          n = A->rmap->n;
3629ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
363aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
3649ae82921SPaul Mullowney   cusparseStatus_t                  stat;
3659ae82921SPaul Mullowney   const PetscInt                    *aj = a->j,*adiag = a->diag,*vi;
3669ae82921SPaul Mullowney   const MatScalar                   *aa = a->a,*v;
3679ae82921SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
3689ae82921SPaul Mullowney   PetscScalar                       *AAUp;
3699ae82921SPaul Mullowney   PetscInt                          i,nz, nzUpper, offset;
3709ae82921SPaul Mullowney   PetscErrorCode                    ierr;
3719ae82921SPaul Mullowney 
3729ae82921SPaul Mullowney   PetscFunctionBegin;
3739ae82921SPaul Mullowney   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUSP_CPU) {
3749ae82921SPaul Mullowney     try {
3759ae82921SPaul Mullowney       /* next, figure out the number of nonzeros in the upper triangular matrix. */
3769ae82921SPaul Mullowney       nzUpper = adiag[0]-adiag[n];
3779ae82921SPaul Mullowney 
3789ae82921SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
3799ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUSP(ierr);
3809ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUSP(ierr);
3819ae82921SPaul Mullowney       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUSP(ierr);
3829ae82921SPaul Mullowney 
3839ae82921SPaul Mullowney       /* Fill the upper triangular matrix */
3849ae82921SPaul Mullowney       AiUp[0]=(PetscInt) 0;
3859ae82921SPaul Mullowney       AiUp[n]=nzUpper;
3869ae82921SPaul Mullowney       offset = nzUpper;
3879ae82921SPaul Mullowney       for (i=n-1; i>=0; i--) {
3889ae82921SPaul Mullowney         v  = aa + adiag[i+1] + 1;
3899ae82921SPaul Mullowney         vi = aj + adiag[i+1] + 1;
3909ae82921SPaul Mullowney 
391e057df02SPaul Mullowney         /* number of elements NOT on the diagonal */
3929ae82921SPaul Mullowney         nz = adiag[i] - adiag[i+1]-1;
3939ae82921SPaul Mullowney 
394e057df02SPaul Mullowney         /* decrement the offset */
3959ae82921SPaul Mullowney         offset -= (nz+1);
3969ae82921SPaul Mullowney 
397e057df02SPaul Mullowney         /* first, set the diagonal elements */
3989ae82921SPaul Mullowney         AjUp[offset] = (PetscInt) i;
3999ae82921SPaul Mullowney         AAUp[offset] = 1./v[nz];
4009ae82921SPaul Mullowney         AiUp[i]      = AiUp[i+1] - (nz+1);
4019ae82921SPaul Mullowney 
4029ae82921SPaul Mullowney         ierr = PetscMemcpy(&(AjUp[offset+1]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr);
4039ae82921SPaul Mullowney         ierr = PetscMemcpy(&(AAUp[offset+1]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
4049ae82921SPaul Mullowney       }
4052205254eSKarl Rupp 
406aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
407aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
4082205254eSKarl Rupp 
409aa372e3fSPaul Mullowney       /* Create the matrix description */
410aa372e3fSPaul Mullowney       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUSP(stat);
411aa372e3fSPaul Mullowney       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUSP(stat);
412aa372e3fSPaul Mullowney       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUSP(stat);
413aa372e3fSPaul Mullowney       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUSP(stat);
414aa372e3fSPaul Mullowney       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUSP(stat);
415aa372e3fSPaul Mullowney 
416aa372e3fSPaul Mullowney       /* Create the solve analysis information */
417aa372e3fSPaul Mullowney       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUSP(stat);
418aa372e3fSPaul Mullowney 
419aa372e3fSPaul Mullowney       /* set the operation */
420aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
421aa372e3fSPaul Mullowney 
422aa372e3fSPaul Mullowney       /* set the matrix */
423aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
424aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = n;
425aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = n;
426aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = nzUpper;
427aa372e3fSPaul Mullowney 
428aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
429aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1);
430aa372e3fSPaul Mullowney 
431aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper);
432aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper);
433aa372e3fSPaul Mullowney 
434aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper);
435aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper);
436aa372e3fSPaul Mullowney 
437aa372e3fSPaul Mullowney       /* perform the solve analysis */
438aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
439aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
440aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
441aa372e3fSPaul Mullowney                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUSP(stat);
442aa372e3fSPaul Mullowney 
443aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
444aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
4452205254eSKarl Rupp 
4469ae82921SPaul Mullowney       ierr = cudaFreeHost(AiUp);CHKERRCUSP(ierr);
4479ae82921SPaul Mullowney       ierr = cudaFreeHost(AjUp);CHKERRCUSP(ierr);
4489ae82921SPaul Mullowney       ierr = cudaFreeHost(AAUp);CHKERRCUSP(ierr);
4499ae82921SPaul Mullowney     } catch(char *ex) {
4509ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
4519ae82921SPaul Mullowney     }
4529ae82921SPaul Mullowney   }
4539ae82921SPaul Mullowney   PetscFunctionReturn(0);
4549ae82921SPaul Mullowney }
4559ae82921SPaul Mullowney 
4569ae82921SPaul Mullowney #undef __FUNCT__
457087f3262SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU"
458087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A)
4599ae82921SPaul Mullowney {
4609ae82921SPaul Mullowney   PetscErrorCode               ierr;
4619ae82921SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
4629ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
4639ae82921SPaul Mullowney   IS                           isrow = a->row,iscol = a->icol;
4649ae82921SPaul Mullowney   PetscBool                    row_identity,col_identity;
4659ae82921SPaul Mullowney   const PetscInt               *r,*c;
4669ae82921SPaul Mullowney   PetscInt                     n = A->rmap->n;
4679ae82921SPaul Mullowney 
4689ae82921SPaul Mullowney   PetscFunctionBegin;
469087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr);
470087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr);
4712205254eSKarl Rupp 
472aa372e3fSPaul Mullowney   cusparseTriFactors->workVector = new THRUSTARRAY;
473aa372e3fSPaul Mullowney   cusparseTriFactors->workVector->resize(n);
474aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=a->nz;
4759ae82921SPaul Mullowney 
4769ae82921SPaul Mullowney   A->valid_GPU_matrix = PETSC_CUSP_BOTH;
477e057df02SPaul Mullowney   /*lower triangular indices */
4789ae82921SPaul Mullowney   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
4799ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
4802205254eSKarl Rupp   if (!row_identity) {
481aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
482aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(r, r+n);
4832205254eSKarl Rupp   }
4849ae82921SPaul Mullowney   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
4859ae82921SPaul Mullowney 
486e057df02SPaul Mullowney   /*upper triangular indices */
4879ae82921SPaul Mullowney   ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
4889ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4892205254eSKarl Rupp   if (!col_identity) {
490aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
491aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices->assign(c, c+n);
4922205254eSKarl Rupp   }
4939ae82921SPaul Mullowney   ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
4949ae82921SPaul Mullowney   PetscFunctionReturn(0);
4959ae82921SPaul Mullowney }
4969ae82921SPaul Mullowney 
4979ae82921SPaul Mullowney #undef __FUNCT__
498087f3262SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEBuildICCTriMatrices"
499087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A)
500087f3262SPaul Mullowney {
501087f3262SPaul Mullowney   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
502087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
503aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
504aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
505087f3262SPaul Mullowney   cusparseStatus_t                  stat;
506087f3262SPaul Mullowney   PetscErrorCode                    ierr;
507087f3262SPaul Mullowney   PetscInt                          *AiUp, *AjUp;
508087f3262SPaul Mullowney   PetscScalar                       *AAUp;
509087f3262SPaul Mullowney   PetscScalar                       *AALo;
510087f3262SPaul Mullowney   PetscInt                          nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j;
511087f3262SPaul Mullowney   Mat_SeqSBAIJ                      *b = (Mat_SeqSBAIJ*)A->data;
512087f3262SPaul Mullowney   const PetscInt                    *ai = b->i,*aj = b->j,*vj;
513087f3262SPaul Mullowney   const MatScalar                   *aa = b->a,*v;
514087f3262SPaul Mullowney 
515087f3262SPaul Mullowney   PetscFunctionBegin;
516087f3262SPaul Mullowney   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUSP_CPU) {
517087f3262SPaul Mullowney     try {
518087f3262SPaul Mullowney       /* Allocate Space for the upper triangular matrix */
519087f3262SPaul Mullowney       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUSP(ierr);
520087f3262SPaul Mullowney       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUSP(ierr);
521087f3262SPaul Mullowney       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUSP(ierr);
522087f3262SPaul Mullowney       ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUSP(ierr);
523087f3262SPaul Mullowney 
524087f3262SPaul Mullowney       /* Fill the upper triangular matrix */
525087f3262SPaul Mullowney       AiUp[0]=(PetscInt) 0;
526087f3262SPaul Mullowney       AiUp[n]=nzUpper;
527087f3262SPaul Mullowney       offset = 0;
528087f3262SPaul Mullowney       for (i=0; i<n; i++) {
529087f3262SPaul Mullowney         /* set the pointers */
530087f3262SPaul Mullowney         v  = aa + ai[i];
531087f3262SPaul Mullowney         vj = aj + ai[i];
532087f3262SPaul Mullowney         nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */
533087f3262SPaul Mullowney 
534087f3262SPaul Mullowney         /* first, set the diagonal elements */
535087f3262SPaul Mullowney         AjUp[offset] = (PetscInt) i;
536087f3262SPaul Mullowney         AAUp[offset] = 1.0/v[nz];
537087f3262SPaul Mullowney         AiUp[i]      = offset;
538087f3262SPaul Mullowney         AALo[offset] = 1.0/v[nz];
539087f3262SPaul Mullowney 
540087f3262SPaul Mullowney         offset+=1;
541087f3262SPaul Mullowney         if (nz>0) {
542087f3262SPaul Mullowney           ierr = PetscMemcpy(&(AjUp[offset]), vj, nz*sizeof(PetscInt));CHKERRQ(ierr);
543087f3262SPaul Mullowney           ierr = PetscMemcpy(&(AAUp[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
544087f3262SPaul Mullowney           for (j=offset; j<offset+nz; j++) {
545087f3262SPaul Mullowney             AAUp[j] = -AAUp[j];
546087f3262SPaul Mullowney             AALo[j] = AAUp[j]/v[nz];
547087f3262SPaul Mullowney           }
548087f3262SPaul Mullowney           offset+=nz;
549087f3262SPaul Mullowney         }
550087f3262SPaul Mullowney       }
551087f3262SPaul Mullowney 
552aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
553aa372e3fSPaul Mullowney       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
554087f3262SPaul Mullowney 
555aa372e3fSPaul Mullowney       /* Create the matrix description */
556aa372e3fSPaul Mullowney       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUSP(stat);
557aa372e3fSPaul Mullowney       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUSP(stat);
558aa372e3fSPaul Mullowney       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUSP(stat);
559aa372e3fSPaul Mullowney       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUSP(stat);
560aa372e3fSPaul Mullowney       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUSP(stat);
561087f3262SPaul Mullowney 
562aa372e3fSPaul Mullowney       /* Create the solve analysis information */
563aa372e3fSPaul Mullowney       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUSP(stat);
564aa372e3fSPaul Mullowney 
565aa372e3fSPaul Mullowney       /* set the operation */
566aa372e3fSPaul Mullowney       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
567aa372e3fSPaul Mullowney 
568aa372e3fSPaul Mullowney       /* set the matrix */
569aa372e3fSPaul Mullowney       upTriFactor->csrMat = new CsrMatrix;
570aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_rows = A->rmap->n;
571aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_cols = A->cmap->n;
572aa372e3fSPaul Mullowney       upTriFactor->csrMat->num_entries = a->nz;
573aa372e3fSPaul Mullowney 
574aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
575aa372e3fSPaul Mullowney       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
576aa372e3fSPaul Mullowney 
577aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
578aa372e3fSPaul Mullowney       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
579aa372e3fSPaul Mullowney 
580aa372e3fSPaul Mullowney       upTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
581aa372e3fSPaul Mullowney       upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz);
582aa372e3fSPaul Mullowney 
583aa372e3fSPaul Mullowney       /* perform the solve analysis */
584aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
585aa372e3fSPaul Mullowney                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
586aa372e3fSPaul Mullowney                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
587aa372e3fSPaul Mullowney                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUSP(stat);
588aa372e3fSPaul Mullowney 
589aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
590aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
591aa372e3fSPaul Mullowney 
592aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
593aa372e3fSPaul Mullowney       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
594aa372e3fSPaul Mullowney 
595aa372e3fSPaul Mullowney       /* Create the matrix description */
596aa372e3fSPaul Mullowney       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUSP(stat);
597aa372e3fSPaul Mullowney       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUSP(stat);
598aa372e3fSPaul Mullowney       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUSP(stat);
599aa372e3fSPaul Mullowney       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUSP(stat);
600aa372e3fSPaul Mullowney       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUSP(stat);
601aa372e3fSPaul Mullowney 
602aa372e3fSPaul Mullowney       /* Create the solve analysis information */
603aa372e3fSPaul Mullowney       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUSP(stat);
604aa372e3fSPaul Mullowney 
605aa372e3fSPaul Mullowney       /* set the operation */
606aa372e3fSPaul Mullowney       loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE;
607aa372e3fSPaul Mullowney 
608aa372e3fSPaul Mullowney       /* set the matrix */
609aa372e3fSPaul Mullowney       loTriFactor->csrMat = new CsrMatrix;
610aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_rows = A->rmap->n;
611aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_cols = A->cmap->n;
612aa372e3fSPaul Mullowney       loTriFactor->csrMat->num_entries = a->nz;
613aa372e3fSPaul Mullowney 
614aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
615aa372e3fSPaul Mullowney       loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
616aa372e3fSPaul Mullowney 
617aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
618aa372e3fSPaul Mullowney       loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
619aa372e3fSPaul Mullowney 
620aa372e3fSPaul Mullowney       loTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
621aa372e3fSPaul Mullowney       loTriFactor->csrMat->values->assign(AALo, AALo+a->nz);
622aa372e3fSPaul Mullowney 
623aa372e3fSPaul Mullowney       /* perform the solve analysis */
624aa372e3fSPaul Mullowney       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
625aa372e3fSPaul Mullowney                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
626aa372e3fSPaul Mullowney                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
627aa372e3fSPaul Mullowney                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUSP(stat);
628aa372e3fSPaul Mullowney 
629aa372e3fSPaul Mullowney       /* assign the pointer. Is this really necessary? */
630aa372e3fSPaul Mullowney       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
631087f3262SPaul Mullowney 
632087f3262SPaul Mullowney       A->valid_GPU_matrix = PETSC_CUSP_BOTH;
633087f3262SPaul Mullowney       ierr = cudaFreeHost(AiUp);CHKERRCUSP(ierr);
634087f3262SPaul Mullowney       ierr = cudaFreeHost(AjUp);CHKERRCUSP(ierr);
635087f3262SPaul Mullowney       ierr = cudaFreeHost(AAUp);CHKERRCUSP(ierr);
636087f3262SPaul Mullowney       ierr = cudaFreeHost(AALo);CHKERRCUSP(ierr);
637087f3262SPaul Mullowney     } catch(char *ex) {
638087f3262SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
639087f3262SPaul Mullowney     }
640087f3262SPaul Mullowney   }
641087f3262SPaul Mullowney   PetscFunctionReturn(0);
642087f3262SPaul Mullowney }
643087f3262SPaul Mullowney 
644087f3262SPaul Mullowney #undef __FUNCT__
645087f3262SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU"
646087f3262SPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A)
6479ae82921SPaul Mullowney {
6489ae82921SPaul Mullowney   PetscErrorCode               ierr;
649087f3262SPaul Mullowney   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
650087f3262SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
651087f3262SPaul Mullowney   IS                           ip = a->row;
652087f3262SPaul Mullowney   const PetscInt               *rip;
653087f3262SPaul Mullowney   PetscBool                    perm_identity;
654087f3262SPaul Mullowney   PetscInt                     n = A->rmap->n;
655087f3262SPaul Mullowney 
656087f3262SPaul Mullowney   PetscFunctionBegin;
657087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr);
658aa372e3fSPaul Mullowney   cusparseTriFactors->workVector = new THRUSTARRAY;
659aa372e3fSPaul Mullowney   cusparseTriFactors->workVector->resize(n);
660aa372e3fSPaul Mullowney   cusparseTriFactors->nnz=(a->nz-n)*2 + n;
661aa372e3fSPaul Mullowney 
662087f3262SPaul Mullowney   /*lower triangular indices */
663087f3262SPaul Mullowney   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
664087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
665087f3262SPaul Mullowney   if (!perm_identity) {
666aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
667aa372e3fSPaul Mullowney     cusparseTriFactors->rpermIndices->assign(rip, rip+n);
668aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
669aa372e3fSPaul Mullowney     cusparseTriFactors->cpermIndices->assign(rip, rip+n);
670087f3262SPaul Mullowney   }
671087f3262SPaul Mullowney   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
672087f3262SPaul Mullowney   PetscFunctionReturn(0);
673087f3262SPaul Mullowney }
674087f3262SPaul Mullowney 
675087f3262SPaul Mullowney #undef __FUNCT__
6769ae82921SPaul Mullowney #define __FUNCT__ "MatLUFactorNumeric_SeqAIJCUSPARSE"
6776fa9248bSJed Brown static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
6789ae82921SPaul Mullowney {
6799ae82921SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
6809ae82921SPaul Mullowney   IS             isrow = b->row,iscol = b->col;
6819ae82921SPaul Mullowney   PetscBool      row_identity,col_identity;
682b175d8bbSPaul Mullowney   PetscErrorCode ierr;
6839ae82921SPaul Mullowney 
6849ae82921SPaul Mullowney   PetscFunctionBegin;
6859ae82921SPaul Mullowney   ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
686e057df02SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
6879ae82921SPaul Mullowney   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
6889ae82921SPaul Mullowney   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
689bda325fcSPaul Mullowney   if (row_identity && col_identity) {
690bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
691bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
692bda325fcSPaul Mullowney   } else {
693bda325fcSPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
694bda325fcSPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
695bda325fcSPaul Mullowney   }
6968dc1d2a3SPaul Mullowney 
697e057df02SPaul Mullowney   /* get the triangular factors */
698087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
6999ae82921SPaul Mullowney   PetscFunctionReturn(0);
7009ae82921SPaul Mullowney }
7019ae82921SPaul Mullowney 
702087f3262SPaul Mullowney #undef __FUNCT__
703087f3262SPaul Mullowney #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJCUSPARSE"
704087f3262SPaul Mullowney static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
705087f3262SPaul Mullowney {
706087f3262SPaul Mullowney   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
707087f3262SPaul Mullowney   IS             ip = b->row;
708087f3262SPaul Mullowney   PetscBool      perm_identity;
709b175d8bbSPaul Mullowney   PetscErrorCode ierr;
710087f3262SPaul Mullowney 
711087f3262SPaul Mullowney   PetscFunctionBegin;
712087f3262SPaul Mullowney   ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
713087f3262SPaul Mullowney 
714087f3262SPaul Mullowney   /* determine which version of MatSolve needs to be used. */
715087f3262SPaul Mullowney   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
716087f3262SPaul Mullowney   if (perm_identity) {
717087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
718087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
719087f3262SPaul Mullowney   } else {
720087f3262SPaul Mullowney     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
721087f3262SPaul Mullowney     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
722087f3262SPaul Mullowney   }
723087f3262SPaul Mullowney 
724087f3262SPaul Mullowney   /* get the triangular factors */
725087f3262SPaul Mullowney   ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
726087f3262SPaul Mullowney   PetscFunctionReturn(0);
727087f3262SPaul Mullowney }
7289ae82921SPaul Mullowney 
729bda325fcSPaul Mullowney #undef __FUNCT__
730bda325fcSPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEAnalyzeTransposeForSolve"
731b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A)
732bda325fcSPaul Mullowney {
733bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
734aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
735aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
736aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
737aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
738bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
739aa372e3fSPaul Mullowney   cusparseIndexBase_t               indexBase;
740aa372e3fSPaul Mullowney   cusparseMatrixType_t              matrixType;
741aa372e3fSPaul Mullowney   cusparseFillMode_t                fillMode;
742aa372e3fSPaul Mullowney   cusparseDiagType_t                diagType;
743b175d8bbSPaul Mullowney 
744bda325fcSPaul Mullowney   PetscFunctionBegin;
745bda325fcSPaul Mullowney 
746aa372e3fSPaul Mullowney   /*********************************************/
747aa372e3fSPaul Mullowney   /* Now the Transpose of the Lower Tri Factor */
748aa372e3fSPaul Mullowney   /*********************************************/
749aa372e3fSPaul Mullowney 
750aa372e3fSPaul Mullowney   /* allocate space for the transpose of the lower triangular factor */
751aa372e3fSPaul Mullowney   loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
752aa372e3fSPaul Mullowney 
753aa372e3fSPaul Mullowney   /* set the matrix descriptors of the lower triangular factor */
754aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(loTriFactor->descr);
755aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(loTriFactor->descr);
756aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
757aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
758aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(loTriFactor->descr);
759aa372e3fSPaul Mullowney 
760aa372e3fSPaul Mullowney   /* Create the matrix description */
761aa372e3fSPaul Mullowney   stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUSP(stat);
762aa372e3fSPaul Mullowney   stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUSP(stat);
763aa372e3fSPaul Mullowney   stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUSP(stat);
764aa372e3fSPaul Mullowney   stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUSP(stat);
765aa372e3fSPaul Mullowney   stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUSP(stat);
766aa372e3fSPaul Mullowney 
767aa372e3fSPaul Mullowney   /* Create the solve analysis information */
768aa372e3fSPaul Mullowney   stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUSP(stat);
769aa372e3fSPaul Mullowney 
770aa372e3fSPaul Mullowney   /* set the operation */
771aa372e3fSPaul Mullowney   loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
772aa372e3fSPaul Mullowney 
773aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the lower triangular factor*/
774aa372e3fSPaul Mullowney   loTriFactorT->csrMat = new CsrMatrix;
775aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows;
776aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols;
777aa372e3fSPaul Mullowney   loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries;
778aa372e3fSPaul Mullowney   loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1);
779aa372e3fSPaul Mullowney   loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries);
780aa372e3fSPaul Mullowney   loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries);
781aa372e3fSPaul Mullowney 
782aa372e3fSPaul Mullowney   /* compute the transpose of the lower triangular factor, i.e. the CSC */
783aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows,
784aa372e3fSPaul Mullowney                           loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries,
785aa372e3fSPaul Mullowney                           loTriFactor->csrMat->values->data().get(),
786aa372e3fSPaul Mullowney                           loTriFactor->csrMat->row_offsets->data().get(),
787aa372e3fSPaul Mullowney                           loTriFactor->csrMat->column_indices->data().get(),
788aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->values->data().get(),
789aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->column_indices->data().get(),
790aa372e3fSPaul Mullowney                           loTriFactorT->csrMat->row_offsets->data().get(),
791aa372e3fSPaul Mullowney                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUSP(stat);
792aa372e3fSPaul Mullowney 
793aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
794aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp,
795aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries,
796aa372e3fSPaul Mullowney                            loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(),
797aa372e3fSPaul Mullowney                            loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(),
798aa372e3fSPaul Mullowney                            loTriFactorT->solveInfo);CHKERRCUSP(stat);
799aa372e3fSPaul Mullowney 
800aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
801aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT;
802aa372e3fSPaul Mullowney 
803aa372e3fSPaul Mullowney   /*********************************************/
804aa372e3fSPaul Mullowney   /* Now the Transpose of the Upper Tri Factor */
805aa372e3fSPaul Mullowney   /*********************************************/
806aa372e3fSPaul Mullowney 
807aa372e3fSPaul Mullowney   /* allocate space for the transpose of the upper triangular factor */
808aa372e3fSPaul Mullowney   upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
809aa372e3fSPaul Mullowney 
810aa372e3fSPaul Mullowney   /* set the matrix descriptors of the upper triangular factor */
811aa372e3fSPaul Mullowney   matrixType = cusparseGetMatType(upTriFactor->descr);
812aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(upTriFactor->descr);
813aa372e3fSPaul Mullowney   fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
814aa372e3fSPaul Mullowney     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
815aa372e3fSPaul Mullowney   diagType = cusparseGetMatDiagType(upTriFactor->descr);
816aa372e3fSPaul Mullowney 
817aa372e3fSPaul Mullowney   /* Create the matrix description */
818aa372e3fSPaul Mullowney   stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUSP(stat);
819aa372e3fSPaul Mullowney   stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUSP(stat);
820aa372e3fSPaul Mullowney   stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUSP(stat);
821aa372e3fSPaul Mullowney   stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUSP(stat);
822aa372e3fSPaul Mullowney   stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUSP(stat);
823aa372e3fSPaul Mullowney 
824aa372e3fSPaul Mullowney   /* Create the solve analysis information */
825aa372e3fSPaul Mullowney   stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUSP(stat);
826aa372e3fSPaul Mullowney 
827aa372e3fSPaul Mullowney   /* set the operation */
828aa372e3fSPaul Mullowney   upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
829aa372e3fSPaul Mullowney 
830aa372e3fSPaul Mullowney   /* allocate GPU space for the CSC of the upper triangular factor*/
831aa372e3fSPaul Mullowney   upTriFactorT->csrMat = new CsrMatrix;
832aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows;
833aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols;
834aa372e3fSPaul Mullowney   upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries;
835aa372e3fSPaul Mullowney   upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1);
836aa372e3fSPaul Mullowney   upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries);
837aa372e3fSPaul Mullowney   upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries);
838aa372e3fSPaul Mullowney 
839aa372e3fSPaul Mullowney   /* compute the transpose of the upper triangular factor, i.e. the CSC */
840aa372e3fSPaul Mullowney   stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows,
841aa372e3fSPaul Mullowney                           upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries,
842aa372e3fSPaul Mullowney                           upTriFactor->csrMat->values->data().get(),
843aa372e3fSPaul Mullowney                           upTriFactor->csrMat->row_offsets->data().get(),
844aa372e3fSPaul Mullowney                           upTriFactor->csrMat->column_indices->data().get(),
845aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->values->data().get(),
846aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->column_indices->data().get(),
847aa372e3fSPaul Mullowney                           upTriFactorT->csrMat->row_offsets->data().get(),
848aa372e3fSPaul Mullowney                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUSP(stat);
849aa372e3fSPaul Mullowney 
850aa372e3fSPaul Mullowney   /* perform the solve analysis on the transposed matrix */
851aa372e3fSPaul Mullowney   stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp,
852aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries,
853aa372e3fSPaul Mullowney                            upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(),
854aa372e3fSPaul Mullowney                            upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(),
855aa372e3fSPaul Mullowney                            upTriFactorT->solveInfo);CHKERRCUSP(stat);
856aa372e3fSPaul Mullowney 
857aa372e3fSPaul Mullowney   /* assign the pointer. Is this really necessary? */
858aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT;
859bda325fcSPaul Mullowney   PetscFunctionReturn(0);
860bda325fcSPaul Mullowney }
861bda325fcSPaul Mullowney 
862bda325fcSPaul Mullowney #undef __FUNCT__
863bda325fcSPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSEGenerateTransposeForMult"
864b175d8bbSPaul Mullowney static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A)
865bda325fcSPaul Mullowney {
866aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
867aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
868aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
869bda325fcSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
870bda325fcSPaul Mullowney   cusparseStatus_t             stat;
871aa372e3fSPaul Mullowney   cusparseIndexBase_t          indexBase;
872b06137fdSPaul Mullowney   cudaError_t                  err;
873b175d8bbSPaul Mullowney 
874bda325fcSPaul Mullowney   PetscFunctionBegin;
875aa372e3fSPaul Mullowney 
876aa372e3fSPaul Mullowney   /* allocate space for the triangular factor information */
877aa372e3fSPaul Mullowney   matstructT = new Mat_SeqAIJCUSPARSEMultStruct;
878aa372e3fSPaul Mullowney   stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUSP(stat);
879aa372e3fSPaul Mullowney   indexBase = cusparseGetMatIndexBase(matstruct->descr);
880aa372e3fSPaul Mullowney   stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUSP(stat);
881aa372e3fSPaul Mullowney   stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUSP(stat);
882aa372e3fSPaul Mullowney 
883b06137fdSPaul Mullowney   /* set alpha and beta */
884b06137fdSPaul Mullowney   err = cudaMalloc((void **)&(matstructT->alpha),sizeof(PetscScalar));CHKERRCUSP(err);
885b06137fdSPaul Mullowney   err = cudaMemcpy(matstructT->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUSP(err);
886b06137fdSPaul Mullowney   err = cudaMalloc((void **)&(matstructT->beta),sizeof(PetscScalar));CHKERRCUSP(err);
887b06137fdSPaul Mullowney   err = cudaMemcpy(matstructT->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUSP(err);
888b06137fdSPaul Mullowney   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUSP(stat);
889b06137fdSPaul Mullowney 
890aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
891aa372e3fSPaul Mullowney     CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
892aa372e3fSPaul Mullowney     CsrMatrix *matrixT= new CsrMatrix;
893aa372e3fSPaul Mullowney     matrixT->num_rows = A->rmap->n;
894aa372e3fSPaul Mullowney     matrixT->num_cols = A->cmap->n;
895aa372e3fSPaul Mullowney     matrixT->num_entries = a->nz;
896aa372e3fSPaul Mullowney     matrixT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
897aa372e3fSPaul Mullowney     matrixT->column_indices = new THRUSTINTARRAY32(a->nz);
898aa372e3fSPaul Mullowney     matrixT->values = new THRUSTARRAY(a->nz);
899aa372e3fSPaul Mullowney 
900aa372e3fSPaul Mullowney     /* compute the transpose of the upper triangular factor, i.e. the CSC */
901aa372e3fSPaul Mullowney     indexBase = cusparseGetMatIndexBase(matstruct->descr);
902aa372e3fSPaul Mullowney     stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows,
903aa372e3fSPaul Mullowney                             matrix->num_cols, matrix->num_entries,
904aa372e3fSPaul Mullowney                             matrix->values->data().get(),
905aa372e3fSPaul Mullowney                             matrix->row_offsets->data().get(),
906aa372e3fSPaul Mullowney                             matrix->column_indices->data().get(),
907aa372e3fSPaul Mullowney                             matrixT->values->data().get(),
908aa372e3fSPaul Mullowney                             matrixT->column_indices->data().get(),
909aa372e3fSPaul Mullowney                             matrixT->row_offsets->data().get(),
910aa372e3fSPaul Mullowney                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUSP(stat);
911aa372e3fSPaul Mullowney 
912aa372e3fSPaul Mullowney     /* assign the pointer */
913aa372e3fSPaul Mullowney     matstructT->mat = matrixT;
914aa372e3fSPaul Mullowney 
915aa372e3fSPaul Mullowney   } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
9162692e278SPaul Mullowney #if CUDA_VERSION>=5000
917aa372e3fSPaul Mullowney     /* First convert HYB to CSR */
918aa372e3fSPaul Mullowney     CsrMatrix *temp= new CsrMatrix;
919aa372e3fSPaul Mullowney     temp->num_rows = A->rmap->n;
920aa372e3fSPaul Mullowney     temp->num_cols = A->cmap->n;
921aa372e3fSPaul Mullowney     temp->num_entries = a->nz;
922aa372e3fSPaul Mullowney     temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
923aa372e3fSPaul Mullowney     temp->column_indices = new THRUSTINTARRAY32(a->nz);
924aa372e3fSPaul Mullowney     temp->values = new THRUSTARRAY(a->nz);
925aa372e3fSPaul Mullowney 
9262692e278SPaul Mullowney 
927aa372e3fSPaul Mullowney     stat = cusparse_hyb2csr(cusparsestruct->handle,
928aa372e3fSPaul Mullowney                             matstruct->descr, (cusparseHybMat_t)matstruct->mat,
929aa372e3fSPaul Mullowney                             temp->values->data().get(),
930aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
931aa372e3fSPaul Mullowney                             temp->column_indices->data().get());CHKERRCUSP(stat);
932aa372e3fSPaul Mullowney 
933aa372e3fSPaul Mullowney     /* Next, convert CSR to CSC (i.e. the matrix transpose) */
934aa372e3fSPaul Mullowney     CsrMatrix *tempT= new CsrMatrix;
935aa372e3fSPaul Mullowney     tempT->num_rows = A->rmap->n;
936aa372e3fSPaul Mullowney     tempT->num_cols = A->cmap->n;
937aa372e3fSPaul Mullowney     tempT->num_entries = a->nz;
938aa372e3fSPaul Mullowney     tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
939aa372e3fSPaul Mullowney     tempT->column_indices = new THRUSTINTARRAY32(a->nz);
940aa372e3fSPaul Mullowney     tempT->values = new THRUSTARRAY(a->nz);
941aa372e3fSPaul Mullowney 
942aa372e3fSPaul Mullowney     stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows,
943aa372e3fSPaul Mullowney                             temp->num_cols, temp->num_entries,
944aa372e3fSPaul Mullowney                             temp->values->data().get(),
945aa372e3fSPaul Mullowney                             temp->row_offsets->data().get(),
946aa372e3fSPaul Mullowney                             temp->column_indices->data().get(),
947aa372e3fSPaul Mullowney                             tempT->values->data().get(),
948aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
949aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
950aa372e3fSPaul Mullowney                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUSP(stat);
951aa372e3fSPaul Mullowney 
952aa372e3fSPaul Mullowney     /* Last, convert CSC to HYB */
953aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat;
954aa372e3fSPaul Mullowney     stat = cusparseCreateHybMat(&hybMat);CHKERRCUSP(stat);
955aa372e3fSPaul Mullowney     cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
956aa372e3fSPaul Mullowney       CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
957aa372e3fSPaul Mullowney     stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n,
958aa372e3fSPaul Mullowney                             matstructT->descr, tempT->values->data().get(),
959aa372e3fSPaul Mullowney                             tempT->row_offsets->data().get(),
960aa372e3fSPaul Mullowney                             tempT->column_indices->data().get(),
961aa372e3fSPaul Mullowney                             hybMat, 0, partition);CHKERRCUSP(stat);
962aa372e3fSPaul Mullowney 
963aa372e3fSPaul Mullowney     /* assign the pointer */
964aa372e3fSPaul Mullowney     matstructT->mat = hybMat;
965aa372e3fSPaul Mullowney 
966aa372e3fSPaul Mullowney     /* delete temporaries */
967aa372e3fSPaul Mullowney     if (tempT) {
968aa372e3fSPaul Mullowney       if (tempT->values) delete (THRUSTARRAY*) tempT->values;
969aa372e3fSPaul Mullowney       if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices;
970aa372e3fSPaul Mullowney       if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets;
971aa372e3fSPaul Mullowney       delete (CsrMatrix*) tempT;
972087f3262SPaul Mullowney     }
973aa372e3fSPaul Mullowney     if (temp) {
974aa372e3fSPaul Mullowney       if (temp->values) delete (THRUSTARRAY*) temp->values;
975aa372e3fSPaul Mullowney       if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices;
976aa372e3fSPaul Mullowney       if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets;
977aa372e3fSPaul Mullowney       delete (CsrMatrix*) temp;
978aa372e3fSPaul Mullowney     }
9792692e278SPaul Mullowney #else
9802692e278SPaul Mullowney     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ELL (Ellpack) and HYB (Hybrid) storage format for the Matrix Transpose (in MatMultTranspose) require CUDA 5.0 or later.");
9812692e278SPaul Mullowney #endif
982aa372e3fSPaul Mullowney   }
983aa372e3fSPaul Mullowney   /* assign the compressed row indices */
984aa372e3fSPaul Mullowney   matstructT->cprowIndices = new THRUSTINTARRAY;
985aa372e3fSPaul Mullowney 
986aa372e3fSPaul Mullowney   /* assign the pointer */
987aa372e3fSPaul Mullowney   ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT;
988bda325fcSPaul Mullowney   PetscFunctionReturn(0);
989bda325fcSPaul Mullowney }
990bda325fcSPaul Mullowney 
991bda325fcSPaul Mullowney #undef __FUNCT__
992bda325fcSPaul Mullowney #define __FUNCT__ "MatSolveTranspose_SeqAIJCUSPARSE"
9936fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
994bda325fcSPaul Mullowney {
995bda325fcSPaul Mullowney   CUSPARRAY                         *xGPU, *bGPU;
996bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
997bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
998aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
999aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1000aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1001b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
1002bda325fcSPaul Mullowney 
1003bda325fcSPaul Mullowney   PetscFunctionBegin;
1004aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
1005aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
1006bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1007aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1008aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1009bda325fcSPaul Mullowney   }
1010bda325fcSPaul Mullowney 
1011bda325fcSPaul Mullowney   /* Get the GPU pointers */
1012bda325fcSPaul Mullowney   ierr = VecCUSPGetArrayWrite(xx,&xGPU);CHKERRQ(ierr);
1013bda325fcSPaul Mullowney   ierr = VecCUSPGetArrayRead(bb,&bGPU);CHKERRQ(ierr);
1014bda325fcSPaul Mullowney 
1015aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
1016aa372e3fSPaul Mullowney   thrust::copy(thrust::make_permutation_iterator(bGPU->begin(), cusparseTriFactors->rpermIndices->begin()),
1017aa372e3fSPaul Mullowney                thrust::make_permutation_iterator(bGPU->end(), cusparseTriFactors->rpermIndices->end()),
1018aa372e3fSPaul Mullowney                xGPU->begin());
1019aa372e3fSPaul Mullowney 
1020aa372e3fSPaul Mullowney   /* First, solve U */
1021aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
1022aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr,
1023aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1024aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1025aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1026aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1027aa372e3fSPaul Mullowney                         xGPU->data().get(), tempGPU->data().get());CHKERRCUSP(stat);
1028aa372e3fSPaul Mullowney 
1029aa372e3fSPaul Mullowney   /* Then, solve L */
1030aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
1031aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr,
1032aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1033aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1034aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1035aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1036aa372e3fSPaul Mullowney                         tempGPU->data().get(), xGPU->data().get());CHKERRCUSP(stat);
1037aa372e3fSPaul Mullowney 
1038aa372e3fSPaul Mullowney   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1039aa372e3fSPaul Mullowney   thrust::copy(thrust::make_permutation_iterator(xGPU->begin(), cusparseTriFactors->cpermIndices->begin()),
1040aa372e3fSPaul Mullowney                thrust::make_permutation_iterator(xGPU->end(), cusparseTriFactors->cpermIndices->end()),
1041aa372e3fSPaul Mullowney                tempGPU->begin());
1042aa372e3fSPaul Mullowney 
1043aa372e3fSPaul Mullowney   /* Copy the temporary to the full solution. */
1044aa372e3fSPaul Mullowney   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU->begin());
1045bda325fcSPaul Mullowney 
1046bda325fcSPaul Mullowney   /* restore */
1047bda325fcSPaul Mullowney   ierr = VecCUSPRestoreArrayRead(bb,&bGPU);CHKERRQ(ierr);
1048bda325fcSPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(xx,&xGPU);CHKERRQ(ierr);
1049bda325fcSPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
1050087f3262SPaul Mullowney 
1051aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1052bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1053bda325fcSPaul Mullowney }
1054bda325fcSPaul Mullowney 
1055bda325fcSPaul Mullowney #undef __FUNCT__
1056bda325fcSPaul Mullowney #define __FUNCT__ "MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering"
10576fa9248bSJed Brown static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
1058bda325fcSPaul Mullowney {
1059bda325fcSPaul Mullowney   CUSPARRAY                         *xGPU,*bGPU;
1060bda325fcSPaul Mullowney   cusparseStatus_t                  stat;
1061bda325fcSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1062aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1063aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1064aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1065b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
1066bda325fcSPaul Mullowney 
1067bda325fcSPaul Mullowney   PetscFunctionBegin;
1068aa372e3fSPaul Mullowney   /* Analyze the matrix and create the transpose ... on the fly */
1069aa372e3fSPaul Mullowney   if (!loTriFactorT && !upTriFactorT) {
1070bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1071aa372e3fSPaul Mullowney     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1072aa372e3fSPaul Mullowney     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1073bda325fcSPaul Mullowney   }
1074bda325fcSPaul Mullowney 
1075bda325fcSPaul Mullowney   /* Get the GPU pointers */
1076bda325fcSPaul Mullowney   ierr = VecCUSPGetArrayWrite(xx,&xGPU);CHKERRQ(ierr);
1077bda325fcSPaul Mullowney   ierr = VecCUSPGetArrayRead(bb,&bGPU);CHKERRQ(ierr);
1078bda325fcSPaul Mullowney 
1079aa372e3fSPaul Mullowney   /* First, solve U */
1080aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
1081aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr,
1082aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->values->data().get(),
1083aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->row_offsets->data().get(),
1084aa372e3fSPaul Mullowney                         upTriFactorT->csrMat->column_indices->data().get(),
1085aa372e3fSPaul Mullowney                         upTriFactorT->solveInfo,
1086aa372e3fSPaul Mullowney                         bGPU->data().get(), tempGPU->data().get());CHKERRCUSP(stat);
1087aa372e3fSPaul Mullowney 
1088aa372e3fSPaul Mullowney   /* Then, solve L */
1089aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
1090aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr,
1091aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->values->data().get(),
1092aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->row_offsets->data().get(),
1093aa372e3fSPaul Mullowney                         loTriFactorT->csrMat->column_indices->data().get(),
1094aa372e3fSPaul Mullowney                         loTriFactorT->solveInfo,
1095aa372e3fSPaul Mullowney                         tempGPU->data().get(), xGPU->data().get());CHKERRCUSP(stat);
1096bda325fcSPaul Mullowney 
1097bda325fcSPaul Mullowney   /* restore */
1098bda325fcSPaul Mullowney   ierr = VecCUSPRestoreArrayRead(bb,&bGPU);CHKERRQ(ierr);
1099bda325fcSPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(xx,&xGPU);CHKERRQ(ierr);
1100bda325fcSPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
1101aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1102bda325fcSPaul Mullowney   PetscFunctionReturn(0);
1103bda325fcSPaul Mullowney }
1104bda325fcSPaul Mullowney 
11059ae82921SPaul Mullowney #undef __FUNCT__
11069ae82921SPaul Mullowney #define __FUNCT__ "MatSolve_SeqAIJCUSPARSE"
11076fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
11089ae82921SPaul Mullowney {
1109bda325fcSPaul Mullowney   CUSPARRAY                         *xGPU,*bGPU;
11109ae82921SPaul Mullowney   cusparseStatus_t                  stat;
11119ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1112aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1113aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1114aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1115b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
1116*ebc8f436SDominic Meiser   VecType                           t;
1117*ebc8f436SDominic Meiser   PetscBool                         flg;
11189ae82921SPaul Mullowney 
11199ae82921SPaul Mullowney   PetscFunctionBegin;
1120*ebc8f436SDominic Meiser   ierr = VecGetType(bb,&t);CHKERRQ(ierr);
1121*ebc8f436SDominic Meiser   ierr = PetscStrcmp(t,VECSEQCUSP,&flg);CHKERRQ(ierr);
1122*ebc8f436SDominic Meiser   if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector of type %s passed into MatSolve_SeqAIJCUSPARSE (Arg #2). Can only deal with %s\n.",t,VECSEQCUSP);
1123*ebc8f436SDominic Meiser   ierr = VecGetType(xx,&t);CHKERRQ(ierr);
1124*ebc8f436SDominic Meiser   ierr = PetscStrcmp(t,VECSEQCUSP,&flg);CHKERRQ(ierr);
1125*ebc8f436SDominic Meiser   if (!flg) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector of type %s passed into MatSolve_SeqAIJCUSPARSE (Arg #3). Can only deal with %s\n.",t,VECSEQCUSP);
1126*ebc8f436SDominic Meiser 
1127e057df02SPaul Mullowney   /* Get the GPU pointers */
11289ae82921SPaul Mullowney   ierr = VecCUSPGetArrayWrite(xx,&xGPU);CHKERRQ(ierr);
11299ae82921SPaul Mullowney   ierr = VecCUSPGetArrayRead(bb,&bGPU);CHKERRQ(ierr);
11309ae82921SPaul Mullowney 
1131aa372e3fSPaul Mullowney   /* First, reorder with the row permutation */
1132aa372e3fSPaul Mullowney   thrust::copy(thrust::make_permutation_iterator(bGPU->begin(), cusparseTriFactors->rpermIndices->begin()),
1133aa372e3fSPaul Mullowney                thrust::make_permutation_iterator(bGPU->end(), cusparseTriFactors->rpermIndices->end()),
1134aa372e3fSPaul Mullowney                xGPU->begin());
1135aa372e3fSPaul Mullowney 
1136aa372e3fSPaul Mullowney   /* Next, solve L */
1137aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
1138aa372e3fSPaul Mullowney                         loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr,
1139aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1140aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1141aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1142aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
1143aa372e3fSPaul Mullowney                         xGPU->data().get(), tempGPU->data().get());CHKERRCUSP(stat);
1144aa372e3fSPaul Mullowney 
1145aa372e3fSPaul Mullowney   /* Then, solve U */
1146aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
1147aa372e3fSPaul Mullowney                         upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr,
1148aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1149aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1150aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1151aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
1152aa372e3fSPaul Mullowney                         tempGPU->data().get(), xGPU->data().get());CHKERRCUSP(stat);
1153aa372e3fSPaul Mullowney 
1154aa372e3fSPaul Mullowney   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1155aa372e3fSPaul Mullowney   thrust::copy(thrust::make_permutation_iterator(xGPU->begin(), cusparseTriFactors->cpermIndices->begin()),
1156aa372e3fSPaul Mullowney                thrust::make_permutation_iterator(xGPU->end(), cusparseTriFactors->cpermIndices->end()),
1157aa372e3fSPaul Mullowney                tempGPU->begin());
1158aa372e3fSPaul Mullowney 
1159aa372e3fSPaul Mullowney   /* Copy the temporary to the full solution. */
1160aa372e3fSPaul Mullowney   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU->begin());
11619ae82921SPaul Mullowney 
11629ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayRead(bb,&bGPU);CHKERRQ(ierr);
11639ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(xx,&xGPU);CHKERRQ(ierr);
11649ae82921SPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
1165aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
11669ae82921SPaul Mullowney   PetscFunctionReturn(0);
11679ae82921SPaul Mullowney }
11689ae82921SPaul Mullowney 
11699ae82921SPaul Mullowney #undef __FUNCT__
11709ae82921SPaul Mullowney #define __FUNCT__ "MatSolve_SeqAIJCUSPARSE_NaturalOrdering"
11716fa9248bSJed Brown static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
11729ae82921SPaul Mullowney {
1173bda325fcSPaul Mullowney   CUSPARRAY                         *xGPU,*bGPU;
11749ae82921SPaul Mullowney   cusparseStatus_t                  stat;
11759ae82921SPaul Mullowney   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1176aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1177aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1178aa372e3fSPaul Mullowney   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1179b175d8bbSPaul Mullowney   PetscErrorCode                    ierr;
11809ae82921SPaul Mullowney 
11819ae82921SPaul Mullowney   PetscFunctionBegin;
1182e057df02SPaul Mullowney   /* Get the GPU pointers */
11839ae82921SPaul Mullowney   ierr = VecCUSPGetArrayWrite(xx,&xGPU);CHKERRQ(ierr);
11849ae82921SPaul Mullowney   ierr = VecCUSPGetArrayRead(bb,&bGPU);CHKERRQ(ierr);
11859ae82921SPaul Mullowney 
1186aa372e3fSPaul Mullowney   /* First, solve L */
1187aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
1188aa372e3fSPaul Mullowney                         loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr,
1189aa372e3fSPaul Mullowney                         loTriFactor->csrMat->values->data().get(),
1190aa372e3fSPaul Mullowney                         loTriFactor->csrMat->row_offsets->data().get(),
1191aa372e3fSPaul Mullowney                         loTriFactor->csrMat->column_indices->data().get(),
1192aa372e3fSPaul Mullowney                         loTriFactor->solveInfo,
1193aa372e3fSPaul Mullowney                         bGPU->data().get(), tempGPU->data().get());CHKERRCUSP(stat);
1194aa372e3fSPaul Mullowney 
1195aa372e3fSPaul Mullowney   /* Next, solve U */
1196aa372e3fSPaul Mullowney   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
1197aa372e3fSPaul Mullowney                         upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr,
1198aa372e3fSPaul Mullowney                         upTriFactor->csrMat->values->data().get(),
1199aa372e3fSPaul Mullowney                         upTriFactor->csrMat->row_offsets->data().get(),
1200aa372e3fSPaul Mullowney                         upTriFactor->csrMat->column_indices->data().get(),
1201aa372e3fSPaul Mullowney                         upTriFactor->solveInfo,
1202aa372e3fSPaul Mullowney                         tempGPU->data().get(), xGPU->data().get());CHKERRCUSP(stat);
12039ae82921SPaul Mullowney 
12049ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayRead(bb,&bGPU);CHKERRQ(ierr);
12059ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(xx,&xGPU);CHKERRQ(ierr);
12069ae82921SPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
1207aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
12089ae82921SPaul Mullowney   PetscFunctionReturn(0);
12099ae82921SPaul Mullowney }
12109ae82921SPaul Mullowney 
12119ae82921SPaul Mullowney #undef __FUNCT__
1212e057df02SPaul Mullowney #define __FUNCT__ "MatSeqAIJCUSPARSECopyToGPU"
12136fa9248bSJed Brown static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A)
12149ae82921SPaul Mullowney {
12159ae82921SPaul Mullowney 
1216aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1217aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
12189ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
12199ae82921SPaul Mullowney   PetscInt                     m = A->rmap->n,*ii,*ridx;
12209ae82921SPaul Mullowney   PetscErrorCode               ierr;
1221aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1222b06137fdSPaul Mullowney   cudaError_t                  err;
12239ae82921SPaul Mullowney 
12249ae82921SPaul Mullowney   PetscFunctionBegin;
12259ae82921SPaul Mullowney   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUSP_CPU) {
12269ae82921SPaul Mullowney     ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
1227aa372e3fSPaul Mullowney     if (matstruct) {
1228aa372e3fSPaul Mullowney       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Data structure should not be initialized.");
12299ae82921SPaul Mullowney     }
12309ae82921SPaul Mullowney     try {
1231aa372e3fSPaul Mullowney       cusparsestruct->nonzerorow=0;
1232aa372e3fSPaul Mullowney       for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0);
12339ae82921SPaul Mullowney 
12349ae82921SPaul Mullowney       if (a->compressedrow.use) {
12359ae82921SPaul Mullowney         m    = a->compressedrow.nrows;
12369ae82921SPaul Mullowney         ii   = a->compressedrow.i;
12379ae82921SPaul Mullowney         ridx = a->compressedrow.rindex;
12389ae82921SPaul Mullowney       } else {
1239b06137fdSPaul Mullowney         /* Forcing compressed row on the GPU */
12409ae82921SPaul Mullowney         int k=0;
1241785e854fSJed Brown         ierr = PetscMalloc1((cusparsestruct->nonzerorow+1), &ii);CHKERRQ(ierr);
1242785e854fSJed Brown         ierr = PetscMalloc1((cusparsestruct->nonzerorow), &ridx);CHKERRQ(ierr);
12439ae82921SPaul Mullowney         ii[0]=0;
12449ae82921SPaul Mullowney         for (int j = 0; j<m; j++) {
12459ae82921SPaul Mullowney           if ((a->i[j+1]-a->i[j])>0) {
12469ae82921SPaul Mullowney             ii[k]  = a->i[j];
12479ae82921SPaul Mullowney             ridx[k]= j;
12489ae82921SPaul Mullowney             k++;
12499ae82921SPaul Mullowney           }
12509ae82921SPaul Mullowney         }
1251aa372e3fSPaul Mullowney         ii[cusparsestruct->nonzerorow] = a->nz;
1252aa372e3fSPaul Mullowney         m = cusparsestruct->nonzerorow;
12539ae82921SPaul Mullowney       }
12549ae82921SPaul Mullowney 
1255aa372e3fSPaul Mullowney       /* allocate space for the triangular factor information */
1256aa372e3fSPaul Mullowney       matstruct = new Mat_SeqAIJCUSPARSEMultStruct;
1257aa372e3fSPaul Mullowney       stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUSP(stat);
1258aa372e3fSPaul Mullowney       stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUSP(stat);
1259aa372e3fSPaul Mullowney       stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUSP(stat);
12609ae82921SPaul Mullowney 
1261b06137fdSPaul Mullowney       err = cudaMalloc((void **)&(matstruct->alpha),sizeof(PetscScalar));CHKERRCUSP(err);
1262b06137fdSPaul Mullowney       err = cudaMemcpy(matstruct->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUSP(err);
1263b06137fdSPaul Mullowney       err = cudaMalloc((void **)&(matstruct->beta),sizeof(PetscScalar));CHKERRCUSP(err);
1264b06137fdSPaul Mullowney       err = cudaMemcpy(matstruct->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUSP(err);
1265b06137fdSPaul Mullowney       stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUSP(stat);
1266b06137fdSPaul Mullowney 
1267aa372e3fSPaul Mullowney       /* Build a hybrid/ellpack matrix if this option is chosen for the storage */
1268aa372e3fSPaul Mullowney       if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1269aa372e3fSPaul Mullowney         /* set the matrix */
1270aa372e3fSPaul Mullowney 	CsrMatrix *matrix= new CsrMatrix;
1271a65300a6SPaul Mullowney 	matrix->num_rows = m;
1272aa372e3fSPaul Mullowney 	matrix->num_cols = A->cmap->n;
1273aa372e3fSPaul Mullowney 	matrix->num_entries = a->nz;
1274a65300a6SPaul Mullowney 	matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1275a65300a6SPaul Mullowney 	matrix->row_offsets->assign(ii, ii + m+1);
12769ae82921SPaul Mullowney 
1277aa372e3fSPaul Mullowney 	matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1278aa372e3fSPaul Mullowney 	matrix->column_indices->assign(a->j, a->j+a->nz);
1279aa372e3fSPaul Mullowney 
1280aa372e3fSPaul Mullowney 	matrix->values = new THRUSTARRAY(a->nz);
1281aa372e3fSPaul Mullowney 	matrix->values->assign(a->a, a->a+a->nz);
1282aa372e3fSPaul Mullowney 
1283aa372e3fSPaul Mullowney         /* assign the pointer */
1284aa372e3fSPaul Mullowney         matstruct->mat = matrix;
1285aa372e3fSPaul Mullowney 
1286aa372e3fSPaul Mullowney       } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
12872692e278SPaul Mullowney #if CUDA_VERSION>=4020
1288aa372e3fSPaul Mullowney 	CsrMatrix *matrix= new CsrMatrix;
1289a65300a6SPaul Mullowney 	matrix->num_rows = m;
1290aa372e3fSPaul Mullowney 	matrix->num_cols = A->cmap->n;
1291aa372e3fSPaul Mullowney 	matrix->num_entries = a->nz;
1292a65300a6SPaul Mullowney 	matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1293a65300a6SPaul Mullowney 	matrix->row_offsets->assign(ii, ii + m+1);
1294aa372e3fSPaul Mullowney 
1295aa372e3fSPaul Mullowney 	matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1296aa372e3fSPaul Mullowney 	matrix->column_indices->assign(a->j, a->j+a->nz);
1297aa372e3fSPaul Mullowney 
1298aa372e3fSPaul Mullowney 	matrix->values = new THRUSTARRAY(a->nz);
1299aa372e3fSPaul Mullowney 	matrix->values->assign(a->a, a->a+a->nz);
1300aa372e3fSPaul Mullowney 
1301aa372e3fSPaul Mullowney         cusparseHybMat_t hybMat;
1302aa372e3fSPaul Mullowney         stat = cusparseCreateHybMat(&hybMat);CHKERRCUSP(stat);
1303aa372e3fSPaul Mullowney         cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
1304aa372e3fSPaul Mullowney           CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
1305a65300a6SPaul Mullowney         stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols,
1306aa372e3fSPaul Mullowney                                 matstruct->descr, matrix->values->data().get(),
1307aa372e3fSPaul Mullowney                                 matrix->row_offsets->data().get(),
1308aa372e3fSPaul Mullowney                                 matrix->column_indices->data().get(),
1309aa372e3fSPaul Mullowney                                 hybMat, 0, partition);CHKERRCUSP(stat);
1310aa372e3fSPaul Mullowney         /* assign the pointer */
1311aa372e3fSPaul Mullowney         matstruct->mat = hybMat;
1312aa372e3fSPaul Mullowney 
1313aa372e3fSPaul Mullowney 	if (matrix) {
1314aa372e3fSPaul Mullowney 	  if (matrix->values) delete (THRUSTARRAY*)matrix->values;
1315aa372e3fSPaul Mullowney 	  if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices;
1316aa372e3fSPaul Mullowney 	  if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets;
1317aa372e3fSPaul Mullowney 	  delete (CsrMatrix*)matrix;
1318087f3262SPaul Mullowney 	}
13192692e278SPaul Mullowney #endif
1320087f3262SPaul Mullowney       }
1321ca45077fSPaul Mullowney 
1322aa372e3fSPaul Mullowney       /* assign the compressed row indices */
1323aa372e3fSPaul Mullowney       matstruct->cprowIndices = new THRUSTINTARRAY(m);
1324aa372e3fSPaul Mullowney       matstruct->cprowIndices->assign(ridx,ridx+m);
1325aa372e3fSPaul Mullowney 
1326aa372e3fSPaul Mullowney       /* assign the pointer */
1327aa372e3fSPaul Mullowney       cusparsestruct->mat = matstruct;
1328aa372e3fSPaul Mullowney 
13299ae82921SPaul Mullowney       if (!a->compressedrow.use) {
13309ae82921SPaul Mullowney         ierr = PetscFree(ii);CHKERRQ(ierr);
13319ae82921SPaul Mullowney         ierr = PetscFree(ridx);CHKERRQ(ierr);
13329ae82921SPaul Mullowney       }
1333aa372e3fSPaul Mullowney       cusparsestruct->workVector = new THRUSTARRAY;
1334aa372e3fSPaul Mullowney       cusparsestruct->workVector->resize(m);
13359ae82921SPaul Mullowney     } catch(char *ex) {
13369ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
13379ae82921SPaul Mullowney     }
13389ae82921SPaul Mullowney     ierr = WaitForGPU();CHKERRCUSP(ierr);
13392205254eSKarl Rupp 
13409ae82921SPaul Mullowney     A->valid_GPU_matrix = PETSC_CUSP_BOTH;
13412205254eSKarl Rupp 
13429ae82921SPaul Mullowney     ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
13439ae82921SPaul Mullowney   }
13449ae82921SPaul Mullowney   PetscFunctionReturn(0);
13459ae82921SPaul Mullowney }
13469ae82921SPaul Mullowney 
13479ae82921SPaul Mullowney #undef __FUNCT__
13489ae82921SPaul Mullowney #define __FUNCT__ "MatGetVecs_SeqAIJCUSPARSE"
13496fa9248bSJed Brown static PetscErrorCode MatGetVecs_SeqAIJCUSPARSE(Mat mat, Vec *right, Vec *left)
13509ae82921SPaul Mullowney {
13519ae82921SPaul Mullowney   PetscErrorCode ierr;
135233d57670SJed Brown   PetscInt rbs,cbs;
13539ae82921SPaul Mullowney 
13549ae82921SPaul Mullowney   PetscFunctionBegin;
135533d57670SJed Brown   ierr = MatGetBlockSizes(mat,&rbs,&cbs);CHKERRQ(ierr);
13569ae82921SPaul Mullowney   if (right) {
1357ce94432eSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)mat),right);CHKERRQ(ierr);
13589ae82921SPaul Mullowney     ierr = VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);CHKERRQ(ierr);
135933d57670SJed Brown     ierr = VecSetBlockSize(*right,cbs);CHKERRQ(ierr);
13609ae82921SPaul Mullowney     ierr = VecSetType(*right,VECSEQCUSP);CHKERRQ(ierr);
13619ae82921SPaul Mullowney     ierr = PetscLayoutReference(mat->cmap,&(*right)->map);CHKERRQ(ierr);
13629ae82921SPaul Mullowney   }
13639ae82921SPaul Mullowney   if (left) {
1364ce94432eSBarry Smith     ierr = VecCreate(PetscObjectComm((PetscObject)mat),left);CHKERRQ(ierr);
13659ae82921SPaul Mullowney     ierr = VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);CHKERRQ(ierr);
136633d57670SJed Brown     ierr = VecSetBlockSize(*left,rbs);CHKERRQ(ierr);
13679ae82921SPaul Mullowney     ierr = VecSetType(*left,VECSEQCUSP);CHKERRQ(ierr);
13689ae82921SPaul Mullowney     ierr = PetscLayoutReference(mat->rmap,&(*left)->map);CHKERRQ(ierr);
13699ae82921SPaul Mullowney   }
13709ae82921SPaul Mullowney   PetscFunctionReturn(0);
13719ae82921SPaul Mullowney }
13729ae82921SPaul Mullowney 
1373aa372e3fSPaul Mullowney struct VecCUSPPlusEquals
1374aa372e3fSPaul Mullowney {
1375aa372e3fSPaul Mullowney   template <typename Tuple>
1376aa372e3fSPaul Mullowney   __host__ __device__
1377aa372e3fSPaul Mullowney   void operator()(Tuple t)
1378aa372e3fSPaul Mullowney   {
1379aa372e3fSPaul Mullowney     thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t);
1380aa372e3fSPaul Mullowney   }
1381aa372e3fSPaul Mullowney };
1382aa372e3fSPaul Mullowney 
13839ae82921SPaul Mullowney #undef __FUNCT__
13849ae82921SPaul Mullowney #define __FUNCT__ "MatMult_SeqAIJCUSPARSE"
13856fa9248bSJed Brown static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
13869ae82921SPaul Mullowney {
13879ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1388aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1389aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1390bda325fcSPaul Mullowney   CUSPARRAY                    *xarray,*yarray;
1391b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1392aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
13939ae82921SPaul Mullowney 
13949ae82921SPaul Mullowney   PetscFunctionBegin;
1395e057df02SPaul Mullowney   /* The line below should not be necessary as it has been moved to MatAssemblyEnd_SeqAIJCUSPARSE
1396e057df02SPaul Mullowney      ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); */
13979ae82921SPaul Mullowney   ierr = VecCUSPGetArrayRead(xx,&xarray);CHKERRQ(ierr);
13989ae82921SPaul Mullowney   ierr = VecCUSPGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1399aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1400aa372e3fSPaul Mullowney     CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1401aa372e3fSPaul Mullowney     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1402aa372e3fSPaul Mullowney                              mat->num_rows, mat->num_cols, mat->num_entries,
1403b06137fdSPaul Mullowney                              matstruct->alpha, matstruct->descr, mat->values->data().get(), mat->row_offsets->data().get(),
1404b06137fdSPaul Mullowney                              mat->column_indices->data().get(), xarray->data().get(), matstruct->beta,
1405aa372e3fSPaul Mullowney                              yarray->data().get());CHKERRCUSP(stat);
1406aa372e3fSPaul Mullowney   } else {
14072692e278SPaul Mullowney #if CUDA_VERSION>=4020
1408aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1409aa372e3fSPaul Mullowney     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1410b06137fdSPaul Mullowney                              matstruct->alpha, matstruct->descr, hybMat,
1411b06137fdSPaul Mullowney                              xarray->data().get(), matstruct->beta,
1412aa372e3fSPaul Mullowney                              yarray->data().get());CHKERRCUSP(stat);
14132692e278SPaul Mullowney #endif
14149ae82921SPaul Mullowney   }
14159ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
14169ae82921SPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1417aa372e3fSPaul Mullowney   if (!cusparsestruct->stream) {
14189ae82921SPaul Mullowney     ierr = WaitForGPU();CHKERRCUSP(ierr);
1419ca45077fSPaul Mullowney   }
1420aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
14219ae82921SPaul Mullowney   PetscFunctionReturn(0);
14229ae82921SPaul Mullowney }
14239ae82921SPaul Mullowney 
14249ae82921SPaul Mullowney #undef __FUNCT__
1425ca45077fSPaul Mullowney #define __FUNCT__ "MatMultTranspose_SeqAIJCUSPARSE"
14266fa9248bSJed Brown static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
1427ca45077fSPaul Mullowney {
1428ca45077fSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1429aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1430aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1431bda325fcSPaul Mullowney   CUSPARRAY                    *xarray,*yarray;
1432b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1433aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
1434ca45077fSPaul Mullowney 
1435ca45077fSPaul Mullowney   PetscFunctionBegin;
1436e057df02SPaul Mullowney   /* The line below should not be necessary as it has been moved to MatAssemblyEnd_SeqAIJCUSPARSE
1437e057df02SPaul Mullowney      ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); */
1438aa372e3fSPaul Mullowney   if (!matstructT) {
1439bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1440aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1441bda325fcSPaul Mullowney   }
1442ca45077fSPaul Mullowney   ierr = VecCUSPGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1443ca45077fSPaul Mullowney   ierr = VecCUSPGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1444aa372e3fSPaul Mullowney 
1445aa372e3fSPaul Mullowney   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1446aa372e3fSPaul Mullowney     CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1447aa372e3fSPaul Mullowney     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1448aa372e3fSPaul Mullowney                              mat->num_rows, mat->num_cols,
1449b06137fdSPaul Mullowney                              mat->num_entries, matstructT->alpha, matstructT->descr,
1450aa372e3fSPaul Mullowney                              mat->values->data().get(), mat->row_offsets->data().get(),
1451b06137fdSPaul Mullowney                              mat->column_indices->data().get(), xarray->data().get(), matstructT->beta,
1452aa372e3fSPaul Mullowney                              yarray->data().get());CHKERRCUSP(stat);
1453aa372e3fSPaul Mullowney   } else {
14542692e278SPaul Mullowney #if CUDA_VERSION>=4020
1455aa372e3fSPaul Mullowney     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1456aa372e3fSPaul Mullowney     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1457b06137fdSPaul Mullowney                              matstructT->alpha, matstructT->descr, hybMat,
1458b06137fdSPaul Mullowney                              xarray->data().get(), matstructT->beta,
1459aa372e3fSPaul Mullowney                              yarray->data().get());CHKERRCUSP(stat);
14602692e278SPaul Mullowney #endif
1461ca45077fSPaul Mullowney   }
1462ca45077fSPaul Mullowney   ierr = VecCUSPRestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1463ca45077fSPaul Mullowney   ierr = VecCUSPRestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1464aa372e3fSPaul Mullowney   if (!cusparsestruct->stream) {
1465ca45077fSPaul Mullowney     ierr = WaitForGPU();CHKERRCUSP(ierr);
1466ca45077fSPaul Mullowney   }
1467aa372e3fSPaul Mullowney   ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
1468ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1469ca45077fSPaul Mullowney }
1470ca45077fSPaul Mullowney 
1471aa372e3fSPaul Mullowney 
1472ca45077fSPaul Mullowney #undef __FUNCT__
14739ae82921SPaul Mullowney #define __FUNCT__ "MatMultAdd_SeqAIJCUSPARSE"
14746fa9248bSJed Brown static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
14759ae82921SPaul Mullowney {
14769ae82921SPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1477aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1478aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1479bda325fcSPaul Mullowney   CUSPARRAY                    *xarray,*yarray,*zarray;
1480b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1481aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
14826e111a19SKarl Rupp 
14839ae82921SPaul Mullowney   PetscFunctionBegin;
1484e057df02SPaul Mullowney   /* The line below should not be necessary as it has been moved to MatAssemblyEnd_SeqAIJCUSPARSE
1485e057df02SPaul Mullowney      ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); */
14869ae82921SPaul Mullowney   try {
14879ae82921SPaul Mullowney     ierr = VecCopy_SeqCUSP(yy,zz);CHKERRQ(ierr);
14889ae82921SPaul Mullowney     ierr = VecCUSPGetArrayRead(xx,&xarray);CHKERRQ(ierr);
14899ae82921SPaul Mullowney     ierr = VecCUSPGetArrayRead(yy,&yarray);CHKERRQ(ierr);
14909ae82921SPaul Mullowney     ierr = VecCUSPGetArrayWrite(zz,&zarray);CHKERRQ(ierr);
14919ae82921SPaul Mullowney 
1492e057df02SPaul Mullowney     /* multiply add */
1493aa372e3fSPaul Mullowney     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1494aa372e3fSPaul Mullowney       CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1495b06137fdSPaul Mullowney       /* here we need to be careful to set the number of rows in the multiply to the
1496b06137fdSPaul Mullowney 	 number of compressed rows in the matrix ... which is equivalent to the
1497b06137fdSPaul Mullowney 	 size of the workVector */
1498aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1499a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1500b06137fdSPaul Mullowney                                mat->num_entries, matstruct->alpha, matstruct->descr,
1501aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
1502b06137fdSPaul Mullowney                                mat->column_indices->data().get(), xarray->data().get(), matstruct->beta,
1503aa372e3fSPaul Mullowney                                cusparsestruct->workVector->data().get());CHKERRCUSP(stat);
1504aa372e3fSPaul Mullowney     } else {
15052692e278SPaul Mullowney #if CUDA_VERSION>=4020
1506aa372e3fSPaul Mullowney       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1507a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1508aa372e3fSPaul Mullowney 	stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1509b06137fdSPaul Mullowney 				 matstruct->alpha, matstruct->descr, hybMat,
1510b06137fdSPaul Mullowney 				 xarray->data().get(), matstruct->beta,
1511aa372e3fSPaul Mullowney 				 cusparsestruct->workVector->data().get());CHKERRCUSP(stat);
1512a65300a6SPaul Mullowney       }
15132692e278SPaul Mullowney #endif
1514aa372e3fSPaul Mullowney     }
1515aa372e3fSPaul Mullowney 
1516aa372e3fSPaul Mullowney     /* scatter the data from the temporary into the full vector with a += operation */
1517aa372e3fSPaul Mullowney     thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zarray->begin(), matstruct->cprowIndices->begin()))),
1518aa372e3fSPaul Mullowney 		     thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zarray->begin(), matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1519aa372e3fSPaul Mullowney 		     VecCUSPPlusEquals());
15209ae82921SPaul Mullowney     ierr = VecCUSPRestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
15219ae82921SPaul Mullowney     ierr = VecCUSPRestoreArrayRead(yy,&yarray);CHKERRQ(ierr);
15229ae82921SPaul Mullowney     ierr = VecCUSPRestoreArrayWrite(zz,&zarray);CHKERRQ(ierr);
15239ae82921SPaul Mullowney 
15249ae82921SPaul Mullowney   } catch(char *ex) {
15259ae82921SPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
15269ae82921SPaul Mullowney   }
15279ae82921SPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
15289ae82921SPaul Mullowney   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
15299ae82921SPaul Mullowney   PetscFunctionReturn(0);
15309ae82921SPaul Mullowney }
15319ae82921SPaul Mullowney 
15329ae82921SPaul Mullowney #undef __FUNCT__
1533b175d8bbSPaul Mullowney #define __FUNCT__ "MatMultTransposeAdd_SeqAIJCUSPARSE"
15346fa9248bSJed Brown static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
1535ca45077fSPaul Mullowney {
1536ca45077fSPaul Mullowney   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1537aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1538aa372e3fSPaul Mullowney   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1539ca45077fSPaul Mullowney   CUSPARRAY                    *xarray,*yarray,*zarray;
1540b175d8bbSPaul Mullowney   PetscErrorCode               ierr;
1541aa372e3fSPaul Mullowney   cusparseStatus_t             stat;
15426e111a19SKarl Rupp 
1543ca45077fSPaul Mullowney   PetscFunctionBegin;
1544e057df02SPaul Mullowney   /* The line below should not be necessary as it has been moved to MatAssemblyEnd_SeqAIJCUSPARSE
1545e057df02SPaul Mullowney      ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr); */
1546aa372e3fSPaul Mullowney   if (!matstructT) {
1547bda325fcSPaul Mullowney     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1548aa372e3fSPaul Mullowney     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1549bda325fcSPaul Mullowney   }
1550aa372e3fSPaul Mullowney 
1551ca45077fSPaul Mullowney   try {
1552ca45077fSPaul Mullowney     ierr = VecCopy_SeqCUSP(yy,zz);CHKERRQ(ierr);
1553ca45077fSPaul Mullowney     ierr = VecCUSPGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1554ca45077fSPaul Mullowney     ierr = VecCUSPGetArrayRead(yy,&yarray);CHKERRQ(ierr);
1555ca45077fSPaul Mullowney     ierr = VecCUSPGetArrayWrite(zz,&zarray);CHKERRQ(ierr);
1556ca45077fSPaul Mullowney 
1557e057df02SPaul Mullowney     /* multiply add with matrix transpose */
1558aa372e3fSPaul Mullowney     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1559aa372e3fSPaul Mullowney       CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1560b06137fdSPaul Mullowney       /* here we need to be careful to set the number of rows in the multiply to the
1561b06137fdSPaul Mullowney 	 number of compressed rows in the matrix ... which is equivalent to the
1562b06137fdSPaul Mullowney 	 size of the workVector */
1563aa372e3fSPaul Mullowney       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1564a65300a6SPaul Mullowney                                mat->num_rows, mat->num_cols,
1565b06137fdSPaul Mullowney 			       mat->num_entries, matstructT->alpha, matstructT->descr,
1566aa372e3fSPaul Mullowney                                mat->values->data().get(), mat->row_offsets->data().get(),
1567b06137fdSPaul Mullowney                                mat->column_indices->data().get(), xarray->data().get(), matstructT->beta,
1568aa372e3fSPaul Mullowney                                cusparsestruct->workVector->data().get());CHKERRCUSP(stat);
1569aa372e3fSPaul Mullowney     } else {
15702692e278SPaul Mullowney #if CUDA_VERSION>=4020
1571aa372e3fSPaul Mullowney       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1572a65300a6SPaul Mullowney       if (cusparsestruct->workVector->size()) {
1573aa372e3fSPaul Mullowney 	stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1574b06137fdSPaul Mullowney 				 matstructT->alpha, matstructT->descr, hybMat,
1575b06137fdSPaul Mullowney 				 xarray->data().get(), matstructT->beta,
1576aa372e3fSPaul Mullowney 				 cusparsestruct->workVector->data().get());CHKERRCUSP(stat);
1577a65300a6SPaul Mullowney       }
15782692e278SPaul Mullowney #endif
1579aa372e3fSPaul Mullowney     }
1580aa372e3fSPaul Mullowney 
1581aa372e3fSPaul Mullowney     /* scatter the data from the temporary into the full vector with a += operation */
1582aa372e3fSPaul Mullowney     thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zarray->begin(), matstructT->cprowIndices->begin()))),
1583aa372e3fSPaul Mullowney 		     thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zarray->begin(), matstructT->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1584aa372e3fSPaul Mullowney 		     VecCUSPPlusEquals());
1585ca45077fSPaul Mullowney 
1586ca45077fSPaul Mullowney     ierr = VecCUSPRestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1587ca45077fSPaul Mullowney     ierr = VecCUSPRestoreArrayRead(yy,&yarray);CHKERRQ(ierr);
1588ca45077fSPaul Mullowney     ierr = VecCUSPRestoreArrayWrite(zz,&zarray);CHKERRQ(ierr);
1589ca45077fSPaul Mullowney 
1590ca45077fSPaul Mullowney   } catch(char *ex) {
1591ca45077fSPaul Mullowney     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1592ca45077fSPaul Mullowney   }
1593ca45077fSPaul Mullowney   ierr = WaitForGPU();CHKERRCUSP(ierr);
1594ca45077fSPaul Mullowney   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1595ca45077fSPaul Mullowney   PetscFunctionReturn(0);
1596ca45077fSPaul Mullowney }
1597ca45077fSPaul Mullowney 
1598ca45077fSPaul Mullowney #undef __FUNCT__
15999ae82921SPaul Mullowney #define __FUNCT__ "MatAssemblyEnd_SeqAIJCUSPARSE"
16006fa9248bSJed Brown static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode)
16019ae82921SPaul Mullowney {
16029ae82921SPaul Mullowney   PetscErrorCode ierr;
16036e111a19SKarl Rupp 
16049ae82921SPaul Mullowney   PetscFunctionBegin;
16059ae82921SPaul Mullowney   ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr);
1606bc3f50f2SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
1607e057df02SPaul Mullowney     ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1608bc3f50f2SPaul Mullowney   }
16099ae82921SPaul Mullowney   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
1610bbf3fe20SPaul Mullowney   A->ops->mult             = MatMult_SeqAIJCUSPARSE;
1611bbf3fe20SPaul Mullowney   A->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1612bbf3fe20SPaul Mullowney   A->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1613bbf3fe20SPaul Mullowney   A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
16149ae82921SPaul Mullowney   PetscFunctionReturn(0);
16159ae82921SPaul Mullowney }
16169ae82921SPaul Mullowney 
16179ae82921SPaul Mullowney /* --------------------------------------------------------------------------------*/
16189ae82921SPaul Mullowney #undef __FUNCT__
16199ae82921SPaul Mullowney #define __FUNCT__ "MatCreateSeqAIJCUSPARSE"
1620e057df02SPaul Mullowney /*@
16219ae82921SPaul Mullowney    MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format
1622e057df02SPaul Mullowney    (the default parallel PETSc format). This matrix will ultimately pushed down
1623e057df02SPaul Mullowney    to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix
1624e057df02SPaul Mullowney    assembly performance the user should preallocate the matrix storage by setting
1625e057df02SPaul Mullowney    the parameter nz (or the array nnz).  By setting these parameters accurately,
1626e057df02SPaul Mullowney    performance during matrix assembly can be increased by more than a factor of 50.
16279ae82921SPaul Mullowney 
16289ae82921SPaul Mullowney    Collective on MPI_Comm
16299ae82921SPaul Mullowney 
16309ae82921SPaul Mullowney    Input Parameters:
16319ae82921SPaul Mullowney +  comm - MPI communicator, set to PETSC_COMM_SELF
16329ae82921SPaul Mullowney .  m - number of rows
16339ae82921SPaul Mullowney .  n - number of columns
16349ae82921SPaul Mullowney .  nz - number of nonzeros per row (same for all rows)
16359ae82921SPaul Mullowney -  nnz - array containing the number of nonzeros in the various rows
16360298fd71SBarry Smith          (possibly different for each row) or NULL
16379ae82921SPaul Mullowney 
16389ae82921SPaul Mullowney    Output Parameter:
16399ae82921SPaul Mullowney .  A - the matrix
16409ae82921SPaul Mullowney 
16419ae82921SPaul Mullowney    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
16429ae82921SPaul Mullowney    MatXXXXSetPreallocation() paradgm instead of this routine directly.
16439ae82921SPaul Mullowney    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
16449ae82921SPaul Mullowney 
16459ae82921SPaul Mullowney    Notes:
16469ae82921SPaul Mullowney    If nnz is given then nz is ignored
16479ae82921SPaul Mullowney 
16489ae82921SPaul Mullowney    The AIJ format (also called the Yale sparse matrix format or
16499ae82921SPaul Mullowney    compressed row storage), is fully compatible with standard Fortran 77
16509ae82921SPaul Mullowney    storage.  That is, the stored row and column indices can begin at
16519ae82921SPaul Mullowney    either one (as in Fortran) or zero.  See the users' manual for details.
16529ae82921SPaul Mullowney 
16539ae82921SPaul Mullowney    Specify the preallocated storage with either nz or nnz (not both).
16540298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
16559ae82921SPaul Mullowney    allocation.  For large problems you MUST preallocate memory or you
16569ae82921SPaul Mullowney    will get TERRIBLE performance, see the users' manual chapter on matrices.
16579ae82921SPaul Mullowney 
16589ae82921SPaul Mullowney    By default, this format uses inodes (identical nodes) when possible, to
16599ae82921SPaul Mullowney    improve numerical efficiency of matrix-vector products and solves. We
16609ae82921SPaul Mullowney    search for consecutive rows with the same nonzero structure, thereby
16619ae82921SPaul Mullowney    reusing matrix information to achieve increased efficiency.
16629ae82921SPaul Mullowney 
16639ae82921SPaul Mullowney    Level: intermediate
16649ae82921SPaul Mullowney 
1665e057df02SPaul Mullowney .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE
16669ae82921SPaul Mullowney @*/
16679ae82921SPaul Mullowney PetscErrorCode  MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
16689ae82921SPaul Mullowney {
16699ae82921SPaul Mullowney   PetscErrorCode ierr;
16709ae82921SPaul Mullowney 
16719ae82921SPaul Mullowney   PetscFunctionBegin;
16729ae82921SPaul Mullowney   ierr = MatCreate(comm,A);CHKERRQ(ierr);
16739ae82921SPaul Mullowney   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
16749ae82921SPaul Mullowney   ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
16759ae82921SPaul Mullowney   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
16769ae82921SPaul Mullowney   PetscFunctionReturn(0);
16779ae82921SPaul Mullowney }
16789ae82921SPaul Mullowney 
16799ae82921SPaul Mullowney #undef __FUNCT__
16809ae82921SPaul Mullowney #define __FUNCT__ "MatDestroy_SeqAIJCUSPARSE"
16816fa9248bSJed Brown static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A)
16829ae82921SPaul Mullowney {
16839ae82921SPaul Mullowney   PetscErrorCode   ierr;
1684aa372e3fSPaul Mullowney   cusparseStatus_t stat;
1685b06137fdSPaul Mullowney   cudaError_t      err;
16869ae82921SPaul Mullowney   PetscFunctionBegin;
16879ae82921SPaul Mullowney   if (A->factortype==MAT_FACTOR_NONE) {
16889ae82921SPaul Mullowney     try {
16899ae82921SPaul Mullowney       if (A->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
1690aa372e3fSPaul Mullowney         Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1691aa372e3fSPaul Mullowney         Mat_SeqAIJCUSPARSEMultStruct *matstruct  = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1692aa372e3fSPaul Mullowney         Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1693aa372e3fSPaul Mullowney 
1694aa372e3fSPaul Mullowney         /* delete all the members in each of the data structures */
1695aa372e3fSPaul Mullowney         if (matstruct) {
1696aa372e3fSPaul Mullowney           if (matstruct->descr) { stat = cusparseDestroyMatDescr(matstruct->descr);CHKERRCUSP(stat); }
1697aa372e3fSPaul Mullowney           if (matstruct->mat) {
1698aa372e3fSPaul Mullowney             if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
16992692e278SPaul Mullowney #if CUDA_VERSION>=4020
1700aa372e3fSPaul Mullowney               cusparseHybMat_t hybMat = (cusparseHybMat_t) matstruct->mat;
1701aa372e3fSPaul Mullowney               stat = cusparseDestroyHybMat(hybMat);CHKERRCUSP(stat);
17022692e278SPaul Mullowney #endif
1703aa372e3fSPaul Mullowney             } else {
1704aa372e3fSPaul Mullowney 	      CsrMatrix* mat = (CsrMatrix*)matstruct->mat;
1705aa372e3fSPaul Mullowney 	      if (mat->values) delete (THRUSTARRAY*)mat->values;
1706aa372e3fSPaul Mullowney 	      if (mat->column_indices) delete (THRUSTINTARRAY32*)mat->column_indices;
1707aa372e3fSPaul Mullowney 	      if (mat->row_offsets) delete (THRUSTINTARRAY32*)mat->row_offsets;
1708aa372e3fSPaul Mullowney               delete (CsrMatrix*)matstruct->mat;
17099ae82921SPaul Mullowney             }
1710aa372e3fSPaul Mullowney           }
1711b06137fdSPaul Mullowney 	  if (matstruct->alpha) { err=cudaFree(matstruct->alpha);CHKERRCUSP(err); }
1712b06137fdSPaul Mullowney 	  if (matstruct->beta) { err=cudaFree(matstruct->beta);CHKERRCUSP(err); }
1713aa372e3fSPaul Mullowney           if (matstruct->cprowIndices) delete (THRUSTINTARRAY*)matstruct->cprowIndices;
1714aa372e3fSPaul Mullowney           if (matstruct) delete (Mat_SeqAIJCUSPARSEMultStruct*)matstruct;
1715aa372e3fSPaul Mullowney         }
1716aa372e3fSPaul Mullowney         if (matstructT) {
1717aa372e3fSPaul Mullowney           if (matstructT->descr) { stat = cusparseDestroyMatDescr(matstructT->descr);CHKERRCUSP(stat); }
1718aa372e3fSPaul Mullowney           if (matstructT->mat) {
1719aa372e3fSPaul Mullowney             if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
17202692e278SPaul Mullowney #if CUDA_VERSION>=4020
1721aa372e3fSPaul Mullowney               cusparseHybMat_t hybMat = (cusparseHybMat_t) matstructT->mat;
1722aa372e3fSPaul Mullowney               stat = cusparseDestroyHybMat(hybMat);CHKERRCUSP(stat);
17232692e278SPaul Mullowney #endif
1724aa372e3fSPaul Mullowney             } else {
1725aa372e3fSPaul Mullowney 	      CsrMatrix* matT = (CsrMatrix*)matstructT->mat;
1726aa372e3fSPaul Mullowney 	      if (matT->values) delete (THRUSTARRAY*)matT->values;
1727aa372e3fSPaul Mullowney 	      if (matT->column_indices) delete (THRUSTINTARRAY32*)matT->column_indices;
1728aa372e3fSPaul Mullowney 	      if (matT->row_offsets) delete (THRUSTINTARRAY32*)matT->row_offsets;
1729aa372e3fSPaul Mullowney               delete (CsrMatrix*)matstructT->mat;
1730aa372e3fSPaul Mullowney             }
1731aa372e3fSPaul Mullowney           }
1732b06137fdSPaul Mullowney 	  if (matstructT->alpha) { err=cudaFree(matstructT->alpha);CHKERRCUSP(err); }
1733b06137fdSPaul Mullowney 	  if (matstructT->beta) { err=cudaFree(matstructT->beta);CHKERRCUSP(err); }
1734aa372e3fSPaul Mullowney           if (matstructT->cprowIndices) delete (THRUSTINTARRAY*)matstructT->cprowIndices;
1735aa372e3fSPaul Mullowney           if (matstructT) delete (Mat_SeqAIJCUSPARSEMultStruct*)matstructT;
1736aa372e3fSPaul Mullowney         }
1737aa372e3fSPaul Mullowney         if (cusparsestruct->workVector) delete (THRUSTARRAY*)cusparsestruct->workVector;
1738aa372e3fSPaul Mullowney         if (cusparsestruct->handle) { stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUSP(stat); }
1739aa372e3fSPaul Mullowney         if (cusparsestruct) delete (Mat_SeqAIJCUSPARSE*)cusparsestruct;
17409ae82921SPaul Mullowney         A->valid_GPU_matrix = PETSC_CUSP_UNALLOCATED;
1741aa372e3fSPaul Mullowney       } else {
1742aa372e3fSPaul Mullowney         Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1743aa372e3fSPaul Mullowney         if (cusparsestruct->handle) { stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUSP(stat); }
1744aa372e3fSPaul Mullowney         if (cusparsestruct) delete (Mat_SeqAIJCUSPARSE*)cusparsestruct;
1745aa372e3fSPaul Mullowney       }
17469ae82921SPaul Mullowney     } catch(char *ex) {
17479ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
17489ae82921SPaul Mullowney     }
17499ae82921SPaul Mullowney   } else {
1750e057df02SPaul Mullowney     /* The triangular factors */
17519ae82921SPaul Mullowney     try {
17529ae82921SPaul Mullowney       Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1753aa372e3fSPaul Mullowney       Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1754aa372e3fSPaul Mullowney       Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1755aa372e3fSPaul Mullowney       Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1756aa372e3fSPaul Mullowney       Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1757aa372e3fSPaul Mullowney 
1758aa372e3fSPaul Mullowney       /* delete all the members in each of the data structures */
1759aa372e3fSPaul Mullowney       if (loTriFactor) {
1760aa372e3fSPaul Mullowney         if (loTriFactor->descr) { stat = cusparseDestroyMatDescr(loTriFactor->descr);CHKERRCUSP(stat); }
1761aa372e3fSPaul Mullowney         if (loTriFactor->solveInfo) { stat = cusparseDestroySolveAnalysisInfo(loTriFactor->solveInfo);CHKERRCUSP(stat); }
1762aa372e3fSPaul Mullowney         if (loTriFactor->csrMat) {
1763aa372e3fSPaul Mullowney 	  if (loTriFactor->csrMat->values) delete (THRUSTARRAY*)loTriFactor->csrMat->values;
1764aa372e3fSPaul Mullowney 	  if (loTriFactor->csrMat->column_indices) delete (THRUSTINTARRAY32*)loTriFactor->csrMat->column_indices;
1765aa372e3fSPaul Mullowney 	  if (loTriFactor->csrMat->row_offsets) delete (THRUSTINTARRAY32*)loTriFactor->csrMat->row_offsets;
1766aa372e3fSPaul Mullowney 	  delete (CsrMatrix*)loTriFactor->csrMat;
1767aa372e3fSPaul Mullowney 	}
1768aa372e3fSPaul Mullowney         if (loTriFactor) delete (Mat_SeqAIJCUSPARSETriFactorStruct*)loTriFactor;
1769aa372e3fSPaul Mullowney       }
1770aa372e3fSPaul Mullowney 
1771aa372e3fSPaul Mullowney       if (upTriFactor) {
1772aa372e3fSPaul Mullowney         if (upTriFactor->descr) { stat = cusparseDestroyMatDescr(upTriFactor->descr);CHKERRCUSP(stat); }
1773aa372e3fSPaul Mullowney         if (upTriFactor->solveInfo) { stat = cusparseDestroySolveAnalysisInfo(upTriFactor->solveInfo);CHKERRCUSP(stat); }
1774aa372e3fSPaul Mullowney         if (upTriFactor->csrMat) {
1775aa372e3fSPaul Mullowney 	  if (upTriFactor->csrMat->values) delete (THRUSTARRAY*)upTriFactor->csrMat->values;
1776aa372e3fSPaul Mullowney 	  if (upTriFactor->csrMat->column_indices) delete (THRUSTINTARRAY32*)upTriFactor->csrMat->column_indices;
1777aa372e3fSPaul Mullowney 	  if (upTriFactor->csrMat->row_offsets) delete (THRUSTINTARRAY32*)upTriFactor->csrMat->row_offsets;
1778aa372e3fSPaul Mullowney 	  delete (CsrMatrix*)upTriFactor->csrMat;
1779aa372e3fSPaul Mullowney 	}
1780aa372e3fSPaul Mullowney         if (upTriFactor) delete (Mat_SeqAIJCUSPARSETriFactorStruct*)upTriFactor;
1781aa372e3fSPaul Mullowney       }
1782aa372e3fSPaul Mullowney 
1783aa372e3fSPaul Mullowney       if (loTriFactorT) {
1784aa372e3fSPaul Mullowney         if (loTriFactorT->descr) { stat = cusparseDestroyMatDescr(loTriFactorT->descr);CHKERRCUSP(stat); }
1785aa372e3fSPaul Mullowney         if (loTriFactorT->solveInfo) { stat = cusparseDestroySolveAnalysisInfo(loTriFactorT->solveInfo);CHKERRCUSP(stat); }
1786aa372e3fSPaul Mullowney         if (loTriFactorT->csrMat) {
1787aa372e3fSPaul Mullowney 	  if (loTriFactorT->csrMat->values) delete (THRUSTARRAY*)loTriFactorT->csrMat->values;
1788aa372e3fSPaul Mullowney 	  if (loTriFactorT->csrMat->column_indices) delete (THRUSTINTARRAY32*)loTriFactorT->csrMat->column_indices;
1789aa372e3fSPaul Mullowney 	  if (loTriFactorT->csrMat->row_offsets) delete (THRUSTINTARRAY32*)loTriFactorT->csrMat->row_offsets;
1790aa372e3fSPaul Mullowney 	  delete (CsrMatrix*)loTriFactorT->csrMat;
1791aa372e3fSPaul Mullowney 	}
1792aa372e3fSPaul Mullowney         if (loTriFactorT) delete (Mat_SeqAIJCUSPARSETriFactorStruct*)loTriFactorT;
1793aa372e3fSPaul Mullowney       }
1794aa372e3fSPaul Mullowney 
1795aa372e3fSPaul Mullowney       if (upTriFactorT) {
1796aa372e3fSPaul Mullowney         if (upTriFactorT->descr) { stat = cusparseDestroyMatDescr(upTriFactorT->descr);CHKERRCUSP(stat); }
1797aa372e3fSPaul Mullowney         if (upTriFactorT->solveInfo) { stat = cusparseDestroySolveAnalysisInfo(upTriFactorT->solveInfo);CHKERRCUSP(stat); }
1798aa372e3fSPaul Mullowney         if (upTriFactorT->csrMat) {
1799aa372e3fSPaul Mullowney 	  if (upTriFactorT->csrMat->values) delete (THRUSTARRAY*)upTriFactorT->csrMat->values;
1800aa372e3fSPaul Mullowney 	  if (upTriFactorT->csrMat->column_indices) delete (THRUSTINTARRAY32*)upTriFactorT->csrMat->column_indices;
1801aa372e3fSPaul Mullowney 	  if (upTriFactorT->csrMat->row_offsets) delete (THRUSTINTARRAY32*)upTriFactorT->csrMat->row_offsets;
1802aa372e3fSPaul Mullowney 	  delete (CsrMatrix*)upTriFactorT->csrMat;
1803aa372e3fSPaul Mullowney 	}
1804aa372e3fSPaul Mullowney         if (upTriFactorT) delete (Mat_SeqAIJCUSPARSETriFactorStruct*)upTriFactorT;
1805aa372e3fSPaul Mullowney       }
1806aa372e3fSPaul Mullowney       if (cusparseTriFactors->rpermIndices) delete (THRUSTINTARRAY*)cusparseTriFactors->rpermIndices;
1807aa372e3fSPaul Mullowney       if (cusparseTriFactors->cpermIndices) delete (THRUSTINTARRAY*)cusparseTriFactors->cpermIndices;
1808aa372e3fSPaul Mullowney       if (cusparseTriFactors->workVector) delete (THRUSTARRAY*)cusparseTriFactors->workVector;
1809aa372e3fSPaul Mullowney       if (cusparseTriFactors->handle) { stat = cusparseDestroy(cusparseTriFactors->handle);CHKERRCUSP(stat); }
1810aa372e3fSPaul Mullowney       if (cusparseTriFactors) delete (Mat_SeqAIJCUSPARSETriFactors*)cusparseTriFactors;
18119ae82921SPaul Mullowney     } catch(char *ex) {
18129ae82921SPaul Mullowney       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
18139ae82921SPaul Mullowney     }
18149ae82921SPaul Mullowney   }
18159ae82921SPaul Mullowney   /*this next line is because MatDestroy tries to PetscFree spptr if it is not zero, and PetscFree only works if the memory was allocated with PetscNew or PetscMalloc, which don't call the constructor */
18169ae82921SPaul Mullowney   A->spptr = 0;
18179ae82921SPaul Mullowney 
18189ae82921SPaul Mullowney   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
18199ae82921SPaul Mullowney   PetscFunctionReturn(0);
18209ae82921SPaul Mullowney }
18219ae82921SPaul Mullowney 
18229ae82921SPaul Mullowney #undef __FUNCT__
18239ae82921SPaul Mullowney #define __FUNCT__ "MatCreate_SeqAIJCUSPARSE"
18248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B)
18259ae82921SPaul Mullowney {
18269ae82921SPaul Mullowney   PetscErrorCode ierr;
1827aa372e3fSPaul Mullowney   cusparseStatus_t stat;
1828aa372e3fSPaul Mullowney   cusparseHandle_t handle=0;
18299ae82921SPaul Mullowney 
18309ae82921SPaul Mullowney   PetscFunctionBegin;
18319ae82921SPaul Mullowney   ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr);
18329ae82921SPaul Mullowney   if (B->factortype==MAT_FACTOR_NONE) {
1833e057df02SPaul Mullowney     /* you cannot check the inode.use flag here since the matrix was just created.
1834e057df02SPaul Mullowney        now build a GPU matrix data structure */
18359ae82921SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSE;
18369ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat          = 0;
1837aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0;
1838aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector   = 0;
1839e057df02SPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->format       = MAT_CUSPARSE_CSR;
1840aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream       = 0;
1841aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle       = 0;
1842aa372e3fSPaul Mullowney     stat = cusparseCreate(&handle);CHKERRCUSP(stat);
1843aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle       = handle;
1844aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream       = 0;
18459ae82921SPaul Mullowney   } else {
18469ae82921SPaul Mullowney     /* NEXT, set the pointers to the triangular factors */
1847debe9ee2SPaul Mullowney     B->spptr = new Mat_SeqAIJCUSPARSETriFactors;
18489ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr          = 0;
18499ae82921SPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr          = 0;
1850aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0;
1851aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0;
1852aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices            = 0;
1853aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices            = 0;
1854aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector              = 0;
1855aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = 0;
1856aa372e3fSPaul Mullowney     stat = cusparseCreate(&handle);CHKERRCUSP(stat);
1857aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = handle;
1858aa372e3fSPaul Mullowney     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz                     = 0;
18599ae82921SPaul Mullowney   }
1860aa372e3fSPaul Mullowney 
1861e057df02SPaul Mullowney   /* Here we overload MatGetFactor_petsc_C which enables -mat_type aijcusparse to use the
1862e057df02SPaul Mullowney      default cusparse tri solve. Note the difference with the implementation in
1863e057df02SPaul Mullowney      MatCreate_SeqAIJCUSP in ../seqcusp/aijcusp.cu */
1864bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqaij_cusparse);CHKERRQ(ierr);
18652205254eSKarl Rupp 
18669ae82921SPaul Mullowney   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
18679ae82921SPaul Mullowney   B->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
18689ae82921SPaul Mullowney   B->ops->getvecs          = MatGetVecs_SeqAIJCUSPARSE;
18699ae82921SPaul Mullowney   B->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
1870ca45077fSPaul Mullowney   B->ops->mult             = MatMult_SeqAIJCUSPARSE;
1871ca45077fSPaul Mullowney   B->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1872ca45077fSPaul Mullowney   B->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1873ca45077fSPaul Mullowney   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
18742205254eSKarl Rupp 
18759ae82921SPaul Mullowney   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
18762205254eSKarl Rupp 
18779ae82921SPaul Mullowney   B->valid_GPU_matrix = PETSC_CUSP_UNALLOCATED;
18782205254eSKarl Rupp 
1879bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
18809ae82921SPaul Mullowney   PetscFunctionReturn(0);
18819ae82921SPaul Mullowney }
18829ae82921SPaul Mullowney 
1883e057df02SPaul Mullowney /*M
1884e057df02SPaul Mullowney    MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices.
1885e057df02SPaul Mullowney 
1886e057df02SPaul Mullowney    A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either
18872692e278SPaul Mullowney    CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later.
18882692e278SPaul Mullowney    All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library.
1889e057df02SPaul Mullowney 
1890e057df02SPaul Mullowney    Options Database Keys:
1891e057df02SPaul Mullowney +  -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions()
1892aa372e3fSPaul 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).
1893aa372e3fSPaul Mullowney .  -mat_cusparse_mult_storage_format csr - sets the storage format of matrices (for MatMult) during a call to MatSetFromOptions(). Other options include ell (ellpack) or hyb (hybrid).
1894e057df02SPaul Mullowney 
1895e057df02SPaul Mullowney   Level: beginner
1896e057df02SPaul Mullowney 
18978468deeeSKarl Rupp .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
1898e057df02SPaul Mullowney M*/
1899