xref: /petsc/src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu (revision 3d96e9964ff330fd2a9eee374bcd4376da7efe60)
1 /*
2   Defines the basic matrix operations for the AIJ (compressed row)
3   matrix storage format using the CUSPARSE library,
4 */
5 #define PETSC_SKIP_SPINLOCK
6 
7 #include <petscconf.h>
8 #include <../src/mat/impls/aij/seq/aij.h>          /*I "petscmat.h" I*/
9 #include <../src/mat/impls/sbaij/seq/sbaij.h>
10 #include <../src/vec/vec/impls/dvecimpl.h>
11 #include <petsc/private/vecimpl.h>
12 #undef VecType
13 #include <../src/mat/impls/aij/seq/seqcusparse/cusparsematimpl.h>
14 
15 const char *const MatCUSPARSEStorageFormats[] = {"CSR","ELL","HYB","MatCUSPARSEStorageFormat","MAT_CUSPARSE_",0};
16 
17 static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
18 static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,const MatFactorInfo*);
19 static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
20 
21 static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
22 static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat,Mat,IS,IS,const MatFactorInfo*);
23 static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat,Mat,const MatFactorInfo*);
24 
25 static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat,Vec,Vec);
26 static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
27 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
28 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat,Vec,Vec);
29 static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat);
30 static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat,Vec,Vec);
31 static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
32 static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat,Vec,Vec);
33 static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat,Vec,Vec,Vec);
34 
35 static PetscErrorCode CsrMatrix_Destroy(CsrMatrix**);
36 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct**);
37 static PetscErrorCode Mat_SeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct**,MatCUSPARSEStorageFormat);
38 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors**);
39 static PetscErrorCode Mat_SeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE**);
40 
41 PetscErrorCode MatCUSPARSESetStream(Mat A,const cudaStream_t stream)
42 {
43   cusparseStatus_t   stat;
44   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
45 
46   PetscFunctionBegin;
47   cusparsestruct->stream = stream;
48   stat = cusparseSetStream(cusparsestruct->handle,cusparsestruct->stream);CHKERRCUDA(stat);
49   PetscFunctionReturn(0);
50 }
51 
52 PetscErrorCode MatCUSPARSESetHandle(Mat A,const cusparseHandle_t handle)
53 {
54   cusparseStatus_t   stat;
55   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
56 
57   PetscFunctionBegin;
58   if (cusparsestruct->handle != handle) {
59     if (cusparsestruct->handle) {
60       stat = cusparseDestroy(cusparsestruct->handle);CHKERRCUDA(stat);
61     }
62     cusparsestruct->handle = handle;
63   }
64   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
65   PetscFunctionReturn(0);
66 }
67 
68 PetscErrorCode MatCUSPARSEClearHandle(Mat A)
69 {
70   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
71   PetscFunctionBegin;
72   if (cusparsestruct->handle)
73     cusparsestruct->handle = 0;
74   PetscFunctionReturn(0);
75 }
76 
77 PetscErrorCode MatFactorGetSolverPackage_seqaij_cusparse(Mat A,const MatSolverPackage *type)
78 {
79   PetscFunctionBegin;
80   *type = MATSOLVERCUSPARSE;
81   PetscFunctionReturn(0);
82 }
83 
84 /*MC
85   MATSOLVERCUSPARSE = "cusparse" - A matrix type providing triangular solvers for seq matrices
86   on a single GPU of type, seqaijcusparse, aijcusparse, or seqaijcusp, aijcusp. Currently supported
87   algorithms are ILU(k) and ICC(k). Typically, deeper factorizations (larger k) results in poorer
88   performance in the triangular solves. Full LU, and Cholesky decompositions can be solved through the
89   CUSPARSE triangular solve algorithm. However, the performance can be quite poor and thus these
90   algorithms are not recommended. This class does NOT support direct solver operations.
91 
92   Level: beginner
93 
94 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage, MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
95 M*/
96 
97 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat A,MatFactorType ftype,Mat *B)
98 {
99   PetscErrorCode ierr;
100   PetscInt       n = A->rmap->n;
101 
102   PetscFunctionBegin;
103   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
104   (*B)->factortype = ftype;
105   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
106   ierr = MatSetType(*B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
107 
108   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
109     ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr);
110     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJCUSPARSE;
111     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJCUSPARSE;
112   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
113     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJCUSPARSE;
114     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJCUSPARSE;
115   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported for CUSPARSE Matrix Types");
116 
117   ierr = MatSeqAIJSetPreallocation(*B,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
118   ierr = PetscObjectComposeFunction((PetscObject)(*B),"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_cusparse);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 PETSC_INTERN PetscErrorCode MatCUSPARSESetFormat_SeqAIJCUSPARSE(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
123 {
124   Mat_SeqAIJCUSPARSE *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
125 
126   PetscFunctionBegin;
127 #if CUDA_VERSION>=4020
128   switch (op) {
129   case MAT_CUSPARSE_MULT:
130     cusparsestruct->format = format;
131     break;
132   case MAT_CUSPARSE_ALL:
133     cusparsestruct->format = format;
134     break;
135   default:
136     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPARSEFormatOperation. MAT_CUSPARSE_MULT and MAT_CUSPARSE_ALL are currently supported.",op);
137   }
138 #else
139   if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ELL (Ellpack) and HYB (Hybrid) storage format require CUDA 4.2 or later.");
140 #endif
141   PetscFunctionReturn(0);
142 }
143 
144 /*@
145    MatCUSPARSESetFormat - Sets the storage format of CUSPARSE matrices for a particular
146    operation. Only the MatMult operation can use different GPU storage formats
147    for MPIAIJCUSPARSE matrices.
148    Not Collective
149 
150    Input Parameters:
151 +  A - Matrix of type SEQAIJCUSPARSE
152 .  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.
153 -  format - MatCUSPARSEStorageFormat (one of MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB. The latter two require CUDA 4.2)
154 
155    Output Parameter:
156 
157    Level: intermediate
158 
159 .seealso: MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
160 @*/
161 PetscErrorCode MatCUSPARSESetFormat(Mat A,MatCUSPARSEFormatOperation op,MatCUSPARSEStorageFormat format)
162 {
163   PetscErrorCode ierr;
164 
165   PetscFunctionBegin;
166   PetscValidHeaderSpecific(A, MAT_CLASSID,1);
167   ierr = PetscTryMethod(A, "MatCUSPARSESetFormat_C",(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat),(A,op,format));CHKERRQ(ierr);
168   PetscFunctionReturn(0);
169 }
170 
171 static PetscErrorCode MatSetFromOptions_SeqAIJCUSPARSE(PetscOptionItems *PetscOptionsObject,Mat A)
172 {
173   PetscErrorCode           ierr;
174   MatCUSPARSEStorageFormat format;
175   PetscBool                flg;
176   Mat_SeqAIJCUSPARSE       *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
177 
178   PetscFunctionBegin;
179   ierr = PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSPARSE options");CHKERRQ(ierr);
180   ierr = PetscObjectOptionsBegin((PetscObject)A);
181   if (A->factortype==MAT_FACTOR_NONE) {
182     ierr = PetscOptionsEnum("-mat_cusparse_mult_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV",
183                             "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
184     if (flg) {
185       ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_MULT,format);CHKERRQ(ierr);
186     }
187   }
188   ierr = PetscOptionsEnum("-mat_cusparse_storage_format","sets storage format of (seq)aijcusparse gpu matrices for SpMV and TriSolve",
189                           "MatCUSPARSESetFormat",MatCUSPARSEStorageFormats,(PetscEnum)cusparsestruct->format,(PetscEnum*)&format,&flg);CHKERRQ(ierr);
190   if (flg) {
191     ierr = MatCUSPARSESetFormat(A,MAT_CUSPARSE_ALL,format);CHKERRQ(ierr);
192   }
193   ierr = PetscOptionsEnd();CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 
196 }
197 
198 static PetscErrorCode MatILUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
199 {
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   ierr = MatILUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
204   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
205   PetscFunctionReturn(0);
206 }
207 
208 static PetscErrorCode MatLUFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
209 {
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   ierr = MatLUFactorSymbolic_SeqAIJ(B,A,isrow,iscol,info);CHKERRQ(ierr);
214   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJCUSPARSE;
215   PetscFunctionReturn(0);
216 }
217 
218 static PetscErrorCode MatICCFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
219 {
220   PetscErrorCode ierr;
221 
222   PetscFunctionBegin;
223   ierr = MatICCFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
224   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
225   PetscFunctionReturn(0);
226 }
227 
228 static PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJCUSPARSE(Mat B,Mat A,IS perm,const MatFactorInfo *info)
229 {
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   ierr = MatCholeskyFactorSymbolic_SeqAIJ(B,A,perm,info);CHKERRQ(ierr);
234   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJCUSPARSE;
235   PetscFunctionReturn(0);
236 }
237 
238 static PetscErrorCode MatSeqAIJCUSPARSEBuildILULowerTriMatrix(Mat A)
239 {
240   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
241   PetscInt                          n = A->rmap->n;
242   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
243   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
244   cusparseStatus_t                  stat;
245   const PetscInt                    *ai = a->i,*aj = a->j,*vi;
246   const MatScalar                   *aa = a->a,*v;
247   PetscInt                          *AiLo, *AjLo;
248   PetscScalar                       *AALo;
249   PetscInt                          i,nz, nzLower, offset, rowOffset;
250   PetscErrorCode                    ierr;
251 
252   PetscFunctionBegin;
253   if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) {
254     try {
255       /* first figure out the number of nonzeros in the lower triangular matrix including 1's on the diagonal. */
256       nzLower=n+ai[n]-ai[1];
257 
258       /* Allocate Space for the lower triangular matrix */
259       ierr = cudaMallocHost((void**) &AiLo, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
260       ierr = cudaMallocHost((void**) &AjLo, nzLower*sizeof(PetscInt));CHKERRCUDA(ierr);
261       ierr = cudaMallocHost((void**) &AALo, nzLower*sizeof(PetscScalar));CHKERRCUDA(ierr);
262 
263       /* Fill the lower triangular matrix */
264       AiLo[0]  = (PetscInt) 0;
265       AiLo[n]  = nzLower;
266       AjLo[0]  = (PetscInt) 0;
267       AALo[0]  = (MatScalar) 1.0;
268       v        = aa;
269       vi       = aj;
270       offset   = 1;
271       rowOffset= 1;
272       for (i=1; i<n; i++) {
273         nz = ai[i+1] - ai[i];
274         /* additional 1 for the term on the diagonal */
275         AiLo[i]    = rowOffset;
276         rowOffset += nz+1;
277 
278         ierr = PetscMemcpy(&(AjLo[offset]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr);
279         ierr = PetscMemcpy(&(AALo[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
280 
281         offset      += nz;
282         AjLo[offset] = (PetscInt) i;
283         AALo[offset] = (MatScalar) 1.0;
284         offset      += 1;
285 
286         v  += nz;
287         vi += nz;
288       }
289 
290       /* allocate space for the triangular factor information */
291       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
292 
293       /* Create the matrix description */
294       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
295       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
296       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
297       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_LOWER);CHKERRCUDA(stat);
298       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
299 
300       /* Create the solve analysis information */
301       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
302 
303       /* set the operation */
304       loTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
305 
306       /* set the matrix */
307       loTriFactor->csrMat = new CsrMatrix;
308       loTriFactor->csrMat->num_rows = n;
309       loTriFactor->csrMat->num_cols = n;
310       loTriFactor->csrMat->num_entries = nzLower;
311 
312       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
313       loTriFactor->csrMat->row_offsets->assign(AiLo, AiLo+n+1);
314 
315       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzLower);
316       loTriFactor->csrMat->column_indices->assign(AjLo, AjLo+nzLower);
317 
318       loTriFactor->csrMat->values = new THRUSTARRAY(nzLower);
319       loTriFactor->csrMat->values->assign(AALo, AALo+nzLower);
320 
321       /* perform the solve analysis */
322       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
323                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
324                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
325                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
326 
327       /* assign the pointer. Is this really necessary? */
328       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
329 
330       ierr = cudaFreeHost(AiLo);CHKERRCUDA(ierr);
331       ierr = cudaFreeHost(AjLo);CHKERRCUDA(ierr);
332       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
333     } catch(char *ex) {
334       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
335     }
336   }
337   PetscFunctionReturn(0);
338 }
339 
340 static PetscErrorCode MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(Mat A)
341 {
342   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
343   PetscInt                          n = A->rmap->n;
344   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
345   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
346   cusparseStatus_t                  stat;
347   const PetscInt                    *aj = a->j,*adiag = a->diag,*vi;
348   const MatScalar                   *aa = a->a,*v;
349   PetscInt                          *AiUp, *AjUp;
350   PetscScalar                       *AAUp;
351   PetscInt                          i,nz, nzUpper, offset;
352   PetscErrorCode                    ierr;
353 
354   PetscFunctionBegin;
355   if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) {
356     try {
357       /* next, figure out the number of nonzeros in the upper triangular matrix. */
358       nzUpper = adiag[0]-adiag[n];
359 
360       /* Allocate Space for the upper triangular matrix */
361       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
362       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
363       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
364 
365       /* Fill the upper triangular matrix */
366       AiUp[0]=(PetscInt) 0;
367       AiUp[n]=nzUpper;
368       offset = nzUpper;
369       for (i=n-1; i>=0; i--) {
370         v  = aa + adiag[i+1] + 1;
371         vi = aj + adiag[i+1] + 1;
372 
373         /* number of elements NOT on the diagonal */
374         nz = adiag[i] - adiag[i+1]-1;
375 
376         /* decrement the offset */
377         offset -= (nz+1);
378 
379         /* first, set the diagonal elements */
380         AjUp[offset] = (PetscInt) i;
381         AAUp[offset] = (MatScalar)1./v[nz];
382         AiUp[i]      = AiUp[i+1] - (nz+1);
383 
384         ierr = PetscMemcpy(&(AjUp[offset+1]), vi, nz*sizeof(PetscInt));CHKERRQ(ierr);
385         ierr = PetscMemcpy(&(AAUp[offset+1]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
386       }
387 
388       /* allocate space for the triangular factor information */
389       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
390 
391       /* Create the matrix description */
392       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
393       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
394       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
395       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
396       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
397 
398       /* Create the solve analysis information */
399       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
400 
401       /* set the operation */
402       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
403 
404       /* set the matrix */
405       upTriFactor->csrMat = new CsrMatrix;
406       upTriFactor->csrMat->num_rows = n;
407       upTriFactor->csrMat->num_cols = n;
408       upTriFactor->csrMat->num_entries = nzUpper;
409 
410       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(n+1);
411       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+n+1);
412 
413       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(nzUpper);
414       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+nzUpper);
415 
416       upTriFactor->csrMat->values = new THRUSTARRAY(nzUpper);
417       upTriFactor->csrMat->values->assign(AAUp, AAUp+nzUpper);
418 
419       /* perform the solve analysis */
420       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
421                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
422                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
423                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
424 
425       /* assign the pointer. Is this really necessary? */
426       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
427 
428       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
429       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
430       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
431     } catch(char *ex) {
432       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
433     }
434   }
435   PetscFunctionReturn(0);
436 }
437 
438 static PetscErrorCode MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(Mat A)
439 {
440   PetscErrorCode               ierr;
441   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
442   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
443   IS                           isrow = a->row,iscol = a->icol;
444   PetscBool                    row_identity,col_identity;
445   const PetscInt               *r,*c;
446   PetscInt                     n = A->rmap->n;
447 
448   PetscFunctionBegin;
449   ierr = MatSeqAIJCUSPARSEBuildILULowerTriMatrix(A);CHKERRQ(ierr);
450   ierr = MatSeqAIJCUSPARSEBuildILUUpperTriMatrix(A);CHKERRQ(ierr);
451 
452   cusparseTriFactors->workVector = new THRUSTARRAY;
453   cusparseTriFactors->workVector->resize(n);
454   cusparseTriFactors->nnz=a->nz;
455 
456   A->valid_GPU_matrix = PETSC_CUDA_BOTH;
457   /*lower triangular indices */
458   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
459   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
460   if (!row_identity) {
461     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
462     cusparseTriFactors->rpermIndices->assign(r, r+n);
463   }
464   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
465 
466   /*upper triangular indices */
467   ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
468   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
469   if (!col_identity) {
470     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
471     cusparseTriFactors->cpermIndices->assign(c, c+n);
472   }
473   ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
474   PetscFunctionReturn(0);
475 }
476 
477 static PetscErrorCode MatSeqAIJCUSPARSEBuildICCTriMatrices(Mat A)
478 {
479   Mat_SeqAIJ                        *a = (Mat_SeqAIJ*)A->data;
480   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
481   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
482   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
483   cusparseStatus_t                  stat;
484   PetscErrorCode                    ierr;
485   PetscInt                          *AiUp, *AjUp;
486   PetscScalar                       *AAUp;
487   PetscScalar                       *AALo;
488   PetscInt                          nzUpper = a->nz,n = A->rmap->n,i,offset,nz,j;
489   Mat_SeqSBAIJ                      *b = (Mat_SeqSBAIJ*)A->data;
490   const PetscInt                    *ai = b->i,*aj = b->j,*vj;
491   const MatScalar                   *aa = b->a,*v;
492 
493   PetscFunctionBegin;
494   if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) {
495     try {
496       /* Allocate Space for the upper triangular matrix */
497       ierr = cudaMallocHost((void**) &AiUp, (n+1)*sizeof(PetscInt));CHKERRCUDA(ierr);
498       ierr = cudaMallocHost((void**) &AjUp, nzUpper*sizeof(PetscInt));CHKERRCUDA(ierr);
499       ierr = cudaMallocHost((void**) &AAUp, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
500       ierr = cudaMallocHost((void**) &AALo, nzUpper*sizeof(PetscScalar));CHKERRCUDA(ierr);
501 
502       /* Fill the upper triangular matrix */
503       AiUp[0]=(PetscInt) 0;
504       AiUp[n]=nzUpper;
505       offset = 0;
506       for (i=0; i<n; i++) {
507         /* set the pointers */
508         v  = aa + ai[i];
509         vj = aj + ai[i];
510         nz = ai[i+1] - ai[i] - 1; /* exclude diag[i] */
511 
512         /* first, set the diagonal elements */
513         AjUp[offset] = (PetscInt) i;
514         AAUp[offset] = (MatScalar)1.0/v[nz];
515         AiUp[i]      = offset;
516         AALo[offset] = (MatScalar)1.0/v[nz];
517 
518         offset+=1;
519         if (nz>0) {
520           ierr = PetscMemcpy(&(AjUp[offset]), vj, nz*sizeof(PetscInt));CHKERRQ(ierr);
521           ierr = PetscMemcpy(&(AAUp[offset]), v, nz*sizeof(PetscScalar));CHKERRQ(ierr);
522           for (j=offset; j<offset+nz; j++) {
523             AAUp[j] = -AAUp[j];
524             AALo[j] = AAUp[j]/v[nz];
525           }
526           offset+=nz;
527         }
528       }
529 
530       /* allocate space for the triangular factor information */
531       upTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
532 
533       /* Create the matrix description */
534       stat = cusparseCreateMatDescr(&upTriFactor->descr);CHKERRCUDA(stat);
535       stat = cusparseSetMatIndexBase(upTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
536       stat = cusparseSetMatType(upTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
537       stat = cusparseSetMatFillMode(upTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
538       stat = cusparseSetMatDiagType(upTriFactor->descr, CUSPARSE_DIAG_TYPE_UNIT);CHKERRCUDA(stat);
539 
540       /* Create the solve analysis information */
541       stat = cusparseCreateSolveAnalysisInfo(&upTriFactor->solveInfo);CHKERRCUDA(stat);
542 
543       /* set the operation */
544       upTriFactor->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
545 
546       /* set the matrix */
547       upTriFactor->csrMat = new CsrMatrix;
548       upTriFactor->csrMat->num_rows = A->rmap->n;
549       upTriFactor->csrMat->num_cols = A->cmap->n;
550       upTriFactor->csrMat->num_entries = a->nz;
551 
552       upTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
553       upTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
554 
555       upTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
556       upTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
557 
558       upTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
559       upTriFactor->csrMat->values->assign(AAUp, AAUp+a->nz);
560 
561       /* perform the solve analysis */
562       stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactor->solveOp,
563                                upTriFactor->csrMat->num_rows, upTriFactor->csrMat->num_entries, upTriFactor->descr,
564                                upTriFactor->csrMat->values->data().get(), upTriFactor->csrMat->row_offsets->data().get(),
565                                upTriFactor->csrMat->column_indices->data().get(), upTriFactor->solveInfo);CHKERRCUDA(stat);
566 
567       /* assign the pointer. Is this really necessary? */
568       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtr = upTriFactor;
569 
570       /* allocate space for the triangular factor information */
571       loTriFactor = new Mat_SeqAIJCUSPARSETriFactorStruct;
572 
573       /* Create the matrix description */
574       stat = cusparseCreateMatDescr(&loTriFactor->descr);CHKERRCUDA(stat);
575       stat = cusparseSetMatIndexBase(loTriFactor->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
576       stat = cusparseSetMatType(loTriFactor->descr, CUSPARSE_MATRIX_TYPE_TRIANGULAR);CHKERRCUDA(stat);
577       stat = cusparseSetMatFillMode(loTriFactor->descr, CUSPARSE_FILL_MODE_UPPER);CHKERRCUDA(stat);
578       stat = cusparseSetMatDiagType(loTriFactor->descr, CUSPARSE_DIAG_TYPE_NON_UNIT);CHKERRCUDA(stat);
579 
580       /* Create the solve analysis information */
581       stat = cusparseCreateSolveAnalysisInfo(&loTriFactor->solveInfo);CHKERRCUDA(stat);
582 
583       /* set the operation */
584       loTriFactor->solveOp = CUSPARSE_OPERATION_TRANSPOSE;
585 
586       /* set the matrix */
587       loTriFactor->csrMat = new CsrMatrix;
588       loTriFactor->csrMat->num_rows = A->rmap->n;
589       loTriFactor->csrMat->num_cols = A->cmap->n;
590       loTriFactor->csrMat->num_entries = a->nz;
591 
592       loTriFactor->csrMat->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
593       loTriFactor->csrMat->row_offsets->assign(AiUp, AiUp+A->rmap->n+1);
594 
595       loTriFactor->csrMat->column_indices = new THRUSTINTARRAY32(a->nz);
596       loTriFactor->csrMat->column_indices->assign(AjUp, AjUp+a->nz);
597 
598       loTriFactor->csrMat->values = new THRUSTARRAY(a->nz);
599       loTriFactor->csrMat->values->assign(AALo, AALo+a->nz);
600 
601       /* perform the solve analysis */
602       stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactor->solveOp,
603                                loTriFactor->csrMat->num_rows, loTriFactor->csrMat->num_entries, loTriFactor->descr,
604                                loTriFactor->csrMat->values->data().get(), loTriFactor->csrMat->row_offsets->data().get(),
605                                loTriFactor->csrMat->column_indices->data().get(), loTriFactor->solveInfo);CHKERRCUDA(stat);
606 
607       /* assign the pointer. Is this really necessary? */
608       ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtr = loTriFactor;
609 
610       A->valid_GPU_matrix = PETSC_CUDA_BOTH;
611       ierr = cudaFreeHost(AiUp);CHKERRCUDA(ierr);
612       ierr = cudaFreeHost(AjUp);CHKERRCUDA(ierr);
613       ierr = cudaFreeHost(AAUp);CHKERRCUDA(ierr);
614       ierr = cudaFreeHost(AALo);CHKERRCUDA(ierr);
615     } catch(char *ex) {
616       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
617     }
618   }
619   PetscFunctionReturn(0);
620 }
621 
622 static PetscErrorCode MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(Mat A)
623 {
624   PetscErrorCode               ierr;
625   Mat_SeqAIJ                   *a                  = (Mat_SeqAIJ*)A->data;
626   Mat_SeqAIJCUSPARSETriFactors *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
627   IS                           ip = a->row;
628   const PetscInt               *rip;
629   PetscBool                    perm_identity;
630   PetscInt                     n = A->rmap->n;
631 
632   PetscFunctionBegin;
633   ierr = MatSeqAIJCUSPARSEBuildICCTriMatrices(A);CHKERRQ(ierr);
634   cusparseTriFactors->workVector = new THRUSTARRAY;
635   cusparseTriFactors->workVector->resize(n);
636   cusparseTriFactors->nnz=(a->nz-n)*2 + n;
637 
638   /*lower triangular indices */
639   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
640   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
641   if (!perm_identity) {
642     cusparseTriFactors->rpermIndices = new THRUSTINTARRAY(n);
643     cusparseTriFactors->rpermIndices->assign(rip, rip+n);
644     cusparseTriFactors->cpermIndices = new THRUSTINTARRAY(n);
645     cusparseTriFactors->cpermIndices->assign(rip, rip+n);
646   }
647   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 static PetscErrorCode MatLUFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
652 {
653   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
654   IS             isrow = b->row,iscol = b->col;
655   PetscBool      row_identity,col_identity;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   ierr = MatLUFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
660   /* determine which version of MatSolve needs to be used. */
661   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
662   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
663   if (row_identity && col_identity) {
664     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
665     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
666   } else {
667     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
668     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
669   }
670 
671   /* get the triangular factors */
672   ierr = MatSeqAIJCUSPARSEILUAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
673   PetscFunctionReturn(0);
674 }
675 
676 static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJCUSPARSE(Mat B,Mat A,const MatFactorInfo *info)
677 {
678   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
679   IS             ip = b->row;
680   PetscBool      perm_identity;
681   PetscErrorCode ierr;
682 
683   PetscFunctionBegin;
684   ierr = MatCholeskyFactorNumeric_SeqAIJ(B,A,info);CHKERRQ(ierr);
685 
686   /* determine which version of MatSolve needs to be used. */
687   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
688   if (perm_identity) {
689     B->ops->solve = MatSolve_SeqAIJCUSPARSE_NaturalOrdering;
690     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering;
691   } else {
692     B->ops->solve = MatSolve_SeqAIJCUSPARSE;
693     B->ops->solvetranspose = MatSolveTranspose_SeqAIJCUSPARSE;
694   }
695 
696   /* get the triangular factors */
697   ierr = MatSeqAIJCUSPARSEICCAnalysisAndCopyToGPU(B);CHKERRQ(ierr);
698   PetscFunctionReturn(0);
699 }
700 
701 static PetscErrorCode MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(Mat A)
702 {
703   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
704   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
705   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
706   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
707   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
708   cusparseStatus_t                  stat;
709   cusparseIndexBase_t               indexBase;
710   cusparseMatrixType_t              matrixType;
711   cusparseFillMode_t                fillMode;
712   cusparseDiagType_t                diagType;
713 
714   PetscFunctionBegin;
715 
716   /*********************************************/
717   /* Now the Transpose of the Lower Tri Factor */
718   /*********************************************/
719 
720   /* allocate space for the transpose of the lower triangular factor */
721   loTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
722 
723   /* set the matrix descriptors of the lower triangular factor */
724   matrixType = cusparseGetMatType(loTriFactor->descr);
725   indexBase = cusparseGetMatIndexBase(loTriFactor->descr);
726   fillMode = cusparseGetMatFillMode(loTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
727     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
728   diagType = cusparseGetMatDiagType(loTriFactor->descr);
729 
730   /* Create the matrix description */
731   stat = cusparseCreateMatDescr(&loTriFactorT->descr);CHKERRCUDA(stat);
732   stat = cusparseSetMatIndexBase(loTriFactorT->descr, indexBase);CHKERRCUDA(stat);
733   stat = cusparseSetMatType(loTriFactorT->descr, matrixType);CHKERRCUDA(stat);
734   stat = cusparseSetMatFillMode(loTriFactorT->descr, fillMode);CHKERRCUDA(stat);
735   stat = cusparseSetMatDiagType(loTriFactorT->descr, diagType);CHKERRCUDA(stat);
736 
737   /* Create the solve analysis information */
738   stat = cusparseCreateSolveAnalysisInfo(&loTriFactorT->solveInfo);CHKERRCUDA(stat);
739 
740   /* set the operation */
741   loTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
742 
743   /* allocate GPU space for the CSC of the lower triangular factor*/
744   loTriFactorT->csrMat = new CsrMatrix;
745   loTriFactorT->csrMat->num_rows = loTriFactor->csrMat->num_rows;
746   loTriFactorT->csrMat->num_cols = loTriFactor->csrMat->num_cols;
747   loTriFactorT->csrMat->num_entries = loTriFactor->csrMat->num_entries;
748   loTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(loTriFactor->csrMat->num_rows+1);
749   loTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(loTriFactor->csrMat->num_entries);
750   loTriFactorT->csrMat->values = new THRUSTARRAY(loTriFactor->csrMat->num_entries);
751 
752   /* compute the transpose of the lower triangular factor, i.e. the CSC */
753   stat = cusparse_csr2csc(cusparseTriFactors->handle, loTriFactor->csrMat->num_rows,
754                           loTriFactor->csrMat->num_cols, loTriFactor->csrMat->num_entries,
755                           loTriFactor->csrMat->values->data().get(),
756                           loTriFactor->csrMat->row_offsets->data().get(),
757                           loTriFactor->csrMat->column_indices->data().get(),
758                           loTriFactorT->csrMat->values->data().get(),
759                           loTriFactorT->csrMat->column_indices->data().get(),
760                           loTriFactorT->csrMat->row_offsets->data().get(),
761                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
762 
763   /* perform the solve analysis on the transposed matrix */
764   stat = cusparse_analysis(cusparseTriFactors->handle, loTriFactorT->solveOp,
765                            loTriFactorT->csrMat->num_rows, loTriFactorT->csrMat->num_entries,
766                            loTriFactorT->descr, loTriFactorT->csrMat->values->data().get(),
767                            loTriFactorT->csrMat->row_offsets->data().get(), loTriFactorT->csrMat->column_indices->data().get(),
768                            loTriFactorT->solveInfo);CHKERRCUDA(stat);
769 
770   /* assign the pointer. Is this really necessary? */
771   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->loTriFactorPtrTranspose = loTriFactorT;
772 
773   /*********************************************/
774   /* Now the Transpose of the Upper Tri Factor */
775   /*********************************************/
776 
777   /* allocate space for the transpose of the upper triangular factor */
778   upTriFactorT = new Mat_SeqAIJCUSPARSETriFactorStruct;
779 
780   /* set the matrix descriptors of the upper triangular factor */
781   matrixType = cusparseGetMatType(upTriFactor->descr);
782   indexBase = cusparseGetMatIndexBase(upTriFactor->descr);
783   fillMode = cusparseGetMatFillMode(upTriFactor->descr)==CUSPARSE_FILL_MODE_UPPER ?
784     CUSPARSE_FILL_MODE_LOWER : CUSPARSE_FILL_MODE_UPPER;
785   diagType = cusparseGetMatDiagType(upTriFactor->descr);
786 
787   /* Create the matrix description */
788   stat = cusparseCreateMatDescr(&upTriFactorT->descr);CHKERRCUDA(stat);
789   stat = cusparseSetMatIndexBase(upTriFactorT->descr, indexBase);CHKERRCUDA(stat);
790   stat = cusparseSetMatType(upTriFactorT->descr, matrixType);CHKERRCUDA(stat);
791   stat = cusparseSetMatFillMode(upTriFactorT->descr, fillMode);CHKERRCUDA(stat);
792   stat = cusparseSetMatDiagType(upTriFactorT->descr, diagType);CHKERRCUDA(stat);
793 
794   /* Create the solve analysis information */
795   stat = cusparseCreateSolveAnalysisInfo(&upTriFactorT->solveInfo);CHKERRCUDA(stat);
796 
797   /* set the operation */
798   upTriFactorT->solveOp = CUSPARSE_OPERATION_NON_TRANSPOSE;
799 
800   /* allocate GPU space for the CSC of the upper triangular factor*/
801   upTriFactorT->csrMat = new CsrMatrix;
802   upTriFactorT->csrMat->num_rows = upTriFactor->csrMat->num_rows;
803   upTriFactorT->csrMat->num_cols = upTriFactor->csrMat->num_cols;
804   upTriFactorT->csrMat->num_entries = upTriFactor->csrMat->num_entries;
805   upTriFactorT->csrMat->row_offsets = new THRUSTINTARRAY32(upTriFactor->csrMat->num_rows+1);
806   upTriFactorT->csrMat->column_indices = new THRUSTINTARRAY32(upTriFactor->csrMat->num_entries);
807   upTriFactorT->csrMat->values = new THRUSTARRAY(upTriFactor->csrMat->num_entries);
808 
809   /* compute the transpose of the upper triangular factor, i.e. the CSC */
810   stat = cusparse_csr2csc(cusparseTriFactors->handle, upTriFactor->csrMat->num_rows,
811                           upTriFactor->csrMat->num_cols, upTriFactor->csrMat->num_entries,
812                           upTriFactor->csrMat->values->data().get(),
813                           upTriFactor->csrMat->row_offsets->data().get(),
814                           upTriFactor->csrMat->column_indices->data().get(),
815                           upTriFactorT->csrMat->values->data().get(),
816                           upTriFactorT->csrMat->column_indices->data().get(),
817                           upTriFactorT->csrMat->row_offsets->data().get(),
818                           CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
819 
820   /* perform the solve analysis on the transposed matrix */
821   stat = cusparse_analysis(cusparseTriFactors->handle, upTriFactorT->solveOp,
822                            upTriFactorT->csrMat->num_rows, upTriFactorT->csrMat->num_entries,
823                            upTriFactorT->descr, upTriFactorT->csrMat->values->data().get(),
824                            upTriFactorT->csrMat->row_offsets->data().get(), upTriFactorT->csrMat->column_indices->data().get(),
825                            upTriFactorT->solveInfo);CHKERRCUDA(stat);
826 
827   /* assign the pointer. Is this really necessary? */
828   ((Mat_SeqAIJCUSPARSETriFactors*)A->spptr)->upTriFactorPtrTranspose = upTriFactorT;
829   PetscFunctionReturn(0);
830 }
831 
832 static PetscErrorCode MatSeqAIJCUSPARSEGenerateTransposeForMult(Mat A)
833 {
834   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
835   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
836   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
837   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
838   cusparseStatus_t             stat;
839   cusparseIndexBase_t          indexBase;
840   cudaError_t                  err;
841 
842   PetscFunctionBegin;
843 
844   /* allocate space for the triangular factor information */
845   matstructT = new Mat_SeqAIJCUSPARSEMultStruct;
846   stat = cusparseCreateMatDescr(&matstructT->descr);CHKERRCUDA(stat);
847   indexBase = cusparseGetMatIndexBase(matstruct->descr);
848   stat = cusparseSetMatIndexBase(matstructT->descr, indexBase);CHKERRCUDA(stat);
849   stat = cusparseSetMatType(matstructT->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
850 
851   /* set alpha and beta */
852   err = cudaMalloc((void **)&(matstructT->alpha),sizeof(PetscScalar));CHKERRCUDA(err);
853   err = cudaMemcpy(matstructT->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
854   err = cudaMalloc((void **)&(matstructT->beta),sizeof(PetscScalar));CHKERRCUDA(err);
855   err = cudaMemcpy(matstructT->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
856   stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
857 
858   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
859     CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
860     CsrMatrix *matrixT= new CsrMatrix;
861     matrixT->num_rows = A->rmap->n;
862     matrixT->num_cols = A->cmap->n;
863     matrixT->num_entries = a->nz;
864     matrixT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
865     matrixT->column_indices = new THRUSTINTARRAY32(a->nz);
866     matrixT->values = new THRUSTARRAY(a->nz);
867 
868     /* compute the transpose of the upper triangular factor, i.e. the CSC */
869     indexBase = cusparseGetMatIndexBase(matstruct->descr);
870     stat = cusparse_csr2csc(cusparsestruct->handle, matrix->num_rows,
871                             matrix->num_cols, matrix->num_entries,
872                             matrix->values->data().get(),
873                             matrix->row_offsets->data().get(),
874                             matrix->column_indices->data().get(),
875                             matrixT->values->data().get(),
876                             matrixT->column_indices->data().get(),
877                             matrixT->row_offsets->data().get(),
878                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
879 
880     /* assign the pointer */
881     matstructT->mat = matrixT;
882 
883   } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
884 #if CUDA_VERSION>=5000
885     /* First convert HYB to CSR */
886     CsrMatrix *temp= new CsrMatrix;
887     temp->num_rows = A->rmap->n;
888     temp->num_cols = A->cmap->n;
889     temp->num_entries = a->nz;
890     temp->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
891     temp->column_indices = new THRUSTINTARRAY32(a->nz);
892     temp->values = new THRUSTARRAY(a->nz);
893 
894 
895     stat = cusparse_hyb2csr(cusparsestruct->handle,
896                             matstruct->descr, (cusparseHybMat_t)matstruct->mat,
897                             temp->values->data().get(),
898                             temp->row_offsets->data().get(),
899                             temp->column_indices->data().get());CHKERRCUDA(stat);
900 
901     /* Next, convert CSR to CSC (i.e. the matrix transpose) */
902     CsrMatrix *tempT= new CsrMatrix;
903     tempT->num_rows = A->rmap->n;
904     tempT->num_cols = A->cmap->n;
905     tempT->num_entries = a->nz;
906     tempT->row_offsets = new THRUSTINTARRAY32(A->rmap->n+1);
907     tempT->column_indices = new THRUSTINTARRAY32(a->nz);
908     tempT->values = new THRUSTARRAY(a->nz);
909 
910     stat = cusparse_csr2csc(cusparsestruct->handle, temp->num_rows,
911                             temp->num_cols, temp->num_entries,
912                             temp->values->data().get(),
913                             temp->row_offsets->data().get(),
914                             temp->column_indices->data().get(),
915                             tempT->values->data().get(),
916                             tempT->column_indices->data().get(),
917                             tempT->row_offsets->data().get(),
918                             CUSPARSE_ACTION_NUMERIC, indexBase);CHKERRCUDA(stat);
919 
920     /* Last, convert CSC to HYB */
921     cusparseHybMat_t hybMat;
922     stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
923     cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
924       CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
925     stat = cusparse_csr2hyb(cusparsestruct->handle, A->rmap->n, A->cmap->n,
926                             matstructT->descr, tempT->values->data().get(),
927                             tempT->row_offsets->data().get(),
928                             tempT->column_indices->data().get(),
929                             hybMat, 0, partition);CHKERRCUDA(stat);
930 
931     /* assign the pointer */
932     matstructT->mat = hybMat;
933 
934     /* delete temporaries */
935     if (tempT) {
936       if (tempT->values) delete (THRUSTARRAY*) tempT->values;
937       if (tempT->column_indices) delete (THRUSTINTARRAY32*) tempT->column_indices;
938       if (tempT->row_offsets) delete (THRUSTINTARRAY32*) tempT->row_offsets;
939       delete (CsrMatrix*) tempT;
940     }
941     if (temp) {
942       if (temp->values) delete (THRUSTARRAY*) temp->values;
943       if (temp->column_indices) delete (THRUSTINTARRAY32*) temp->column_indices;
944       if (temp->row_offsets) delete (THRUSTINTARRAY32*) temp->row_offsets;
945       delete (CsrMatrix*) temp;
946     }
947 #else
948     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.");
949 #endif
950   }
951   /* assign the compressed row indices */
952   matstructT->cprowIndices = new THRUSTINTARRAY;
953 
954   /* assign the pointer */
955   ((Mat_SeqAIJCUSPARSE*)A->spptr)->matTranspose = matstructT;
956   PetscFunctionReturn(0);
957 }
958 
959 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
960 {
961   PetscInt                              n = xx->map->n;
962   const PetscScalar                     *barray;
963   PetscScalar                           *xarray;
964   thrust::device_ptr<const PetscScalar> bGPU;
965   thrust::device_ptr<PetscScalar>       xGPU;
966   cusparseStatus_t                      stat;
967   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
968   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
969   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
970   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
971   PetscErrorCode                        ierr;
972 
973   PetscFunctionBegin;
974   /* Analyze the matrix and create the transpose ... on the fly */
975   if (!loTriFactorT && !upTriFactorT) {
976     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
977     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
978     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
979   }
980 
981   /* Get the GPU pointers */
982   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
983   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
984   xGPU = thrust::device_pointer_cast(xarray);
985   bGPU = thrust::device_pointer_cast(barray);
986 
987   /* First, reorder with the row permutation */
988   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
989                thrust::make_permutation_iterator(bGPU+n, cusparseTriFactors->rpermIndices->end()),
990                xGPU);
991 
992   /* First, solve U */
993   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
994                         upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr,
995                         upTriFactorT->csrMat->values->data().get(),
996                         upTriFactorT->csrMat->row_offsets->data().get(),
997                         upTriFactorT->csrMat->column_indices->data().get(),
998                         upTriFactorT->solveInfo,
999                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1000 
1001   /* Then, solve L */
1002   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
1003                         loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr,
1004                         loTriFactorT->csrMat->values->data().get(),
1005                         loTriFactorT->csrMat->row_offsets->data().get(),
1006                         loTriFactorT->csrMat->column_indices->data().get(),
1007                         loTriFactorT->solveInfo,
1008                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1009 
1010   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1011   thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()),
1012                thrust::make_permutation_iterator(xGPU+n, cusparseTriFactors->cpermIndices->end()),
1013                tempGPU->begin());
1014 
1015   /* Copy the temporary to the full solution. */
1016   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU);
1017 
1018   /* restore */
1019   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1020   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1021   ierr = WaitForGPU();CHKERRCUDA(ierr);
1022 
1023   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 static PetscErrorCode MatSolveTranspose_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
1028 {
1029   const PetscScalar                 *barray;
1030   PetscScalar                       *xarray;
1031   cusparseStatus_t                  stat;
1032   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1033   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1034   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorT = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1035   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1036   PetscErrorCode                    ierr;
1037 
1038   PetscFunctionBegin;
1039   /* Analyze the matrix and create the transpose ... on the fly */
1040   if (!loTriFactorT && !upTriFactorT) {
1041     ierr = MatSeqAIJCUSPARSEAnalyzeTransposeForSolve(A);CHKERRQ(ierr);
1042     loTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtrTranspose;
1043     upTriFactorT       = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtrTranspose;
1044   }
1045 
1046   /* Get the GPU pointers */
1047   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1048   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1049 
1050   /* First, solve U */
1051   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactorT->solveOp,
1052                         upTriFactorT->csrMat->num_rows, &ALPHA, upTriFactorT->descr,
1053                         upTriFactorT->csrMat->values->data().get(),
1054                         upTriFactorT->csrMat->row_offsets->data().get(),
1055                         upTriFactorT->csrMat->column_indices->data().get(),
1056                         upTriFactorT->solveInfo,
1057                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1058 
1059   /* Then, solve L */
1060   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactorT->solveOp,
1061                         loTriFactorT->csrMat->num_rows, &ALPHA, loTriFactorT->descr,
1062                         loTriFactorT->csrMat->values->data().get(),
1063                         loTriFactorT->csrMat->row_offsets->data().get(),
1064                         loTriFactorT->csrMat->column_indices->data().get(),
1065                         loTriFactorT->solveInfo,
1066                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1067 
1068   /* restore */
1069   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1070   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1071   ierr = WaitForGPU();CHKERRCUDA(ierr);
1072   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 static PetscErrorCode MatSolve_SeqAIJCUSPARSE(Mat A,Vec bb,Vec xx)
1077 {
1078   const PetscScalar                     *barray;
1079   PetscScalar                           *xarray;
1080   thrust::device_ptr<const PetscScalar> bGPU;
1081   thrust::device_ptr<PetscScalar>       xGPU;
1082   cusparseStatus_t                      stat;
1083   Mat_SeqAIJCUSPARSETriFactors          *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1084   Mat_SeqAIJCUSPARSETriFactorStruct     *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1085   Mat_SeqAIJCUSPARSETriFactorStruct     *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1086   THRUSTARRAY                           *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1087   PetscErrorCode                        ierr;
1088 
1089   PetscFunctionBegin;
1090 
1091   /* Get the GPU pointers */
1092   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1093   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1094   xGPU = thrust::device_pointer_cast(xarray);
1095   bGPU = thrust::device_pointer_cast(barray);
1096 
1097   /* First, reorder with the row permutation */
1098   thrust::copy(thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->begin()),
1099                thrust::make_permutation_iterator(bGPU, cusparseTriFactors->rpermIndices->end()),
1100                xGPU);
1101 
1102   /* Next, solve L */
1103   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
1104                         loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr,
1105                         loTriFactor->csrMat->values->data().get(),
1106                         loTriFactor->csrMat->row_offsets->data().get(),
1107                         loTriFactor->csrMat->column_indices->data().get(),
1108                         loTriFactor->solveInfo,
1109                         xarray, tempGPU->data().get());CHKERRCUDA(stat);
1110 
1111   /* Then, solve U */
1112   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
1113                         upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr,
1114                         upTriFactor->csrMat->values->data().get(),
1115                         upTriFactor->csrMat->row_offsets->data().get(),
1116                         upTriFactor->csrMat->column_indices->data().get(),
1117                         upTriFactor->solveInfo,
1118                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1119 
1120   /* Last, copy the solution, xGPU, into a temporary with the column permutation ... can't be done in place. */
1121   thrust::copy(thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->begin()),
1122                thrust::make_permutation_iterator(xGPU, cusparseTriFactors->cpermIndices->end()),
1123                tempGPU->begin());
1124 
1125   /* Copy the temporary to the full solution. */
1126   thrust::copy(tempGPU->begin(), tempGPU->end(), xGPU);
1127 
1128   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1129   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1130   ierr = WaitForGPU();CHKERRCUDA(ierr);
1131   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1132   PetscFunctionReturn(0);
1133 }
1134 
1135 static PetscErrorCode MatSolve_SeqAIJCUSPARSE_NaturalOrdering(Mat A,Vec bb,Vec xx)
1136 {
1137   const PetscScalar                 *barray;
1138   PetscScalar                       *xarray;
1139   cusparseStatus_t                  stat;
1140   Mat_SeqAIJCUSPARSETriFactors      *cusparseTriFactors = (Mat_SeqAIJCUSPARSETriFactors*)A->spptr;
1141   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->loTriFactorPtr;
1142   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactor = (Mat_SeqAIJCUSPARSETriFactorStruct*)cusparseTriFactors->upTriFactorPtr;
1143   THRUSTARRAY                       *tempGPU = (THRUSTARRAY*)cusparseTriFactors->workVector;
1144   PetscErrorCode                    ierr;
1145 
1146   PetscFunctionBegin;
1147   /* Get the GPU pointers */
1148   ierr = VecCUDAGetArrayWrite(xx,&xarray);CHKERRQ(ierr);
1149   ierr = VecCUDAGetArrayRead(bb,&barray);CHKERRQ(ierr);
1150 
1151   /* First, solve L */
1152   stat = cusparse_solve(cusparseTriFactors->handle, loTriFactor->solveOp,
1153                         loTriFactor->csrMat->num_rows, &ALPHA, loTriFactor->descr,
1154                         loTriFactor->csrMat->values->data().get(),
1155                         loTriFactor->csrMat->row_offsets->data().get(),
1156                         loTriFactor->csrMat->column_indices->data().get(),
1157                         loTriFactor->solveInfo,
1158                         barray, tempGPU->data().get());CHKERRCUDA(stat);
1159 
1160   /* Next, solve U */
1161   stat = cusparse_solve(cusparseTriFactors->handle, upTriFactor->solveOp,
1162                         upTriFactor->csrMat->num_rows, &ALPHA, upTriFactor->descr,
1163                         upTriFactor->csrMat->values->data().get(),
1164                         upTriFactor->csrMat->row_offsets->data().get(),
1165                         upTriFactor->csrMat->column_indices->data().get(),
1166                         upTriFactor->solveInfo,
1167                         tempGPU->data().get(), xarray);CHKERRCUDA(stat);
1168 
1169   ierr = VecCUDARestoreArrayRead(bb,&barray);CHKERRQ(ierr);
1170   ierr = VecCUDARestoreArrayWrite(xx,&xarray);CHKERRQ(ierr);
1171   ierr = WaitForGPU();CHKERRCUDA(ierr);
1172   ierr = PetscLogFlops(2.0*cusparseTriFactors->nnz - A->cmap->n);CHKERRQ(ierr);
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 static PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat A)
1177 {
1178 
1179   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1180   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1181   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1182   PetscInt                     m = A->rmap->n,*ii,*ridx;
1183   PetscErrorCode               ierr;
1184   cusparseStatus_t             stat;
1185   cudaError_t                  err;
1186 
1187   PetscFunctionBegin;
1188   if (A->valid_GPU_matrix == PETSC_CUDA_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUDA_CPU) {
1189     ierr = PetscLogEventBegin(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
1190     if (A->assembled && A->nonzerostate == cusparsestruct->nonzerostate && cusparsestruct->format == MAT_CUSPARSE_CSR) {
1191       CsrMatrix *matrix = (CsrMatrix*)matstruct->mat;
1192       /* copy values only */
1193       matrix->values->assign(a->a, a->a+a->nz);
1194     } else {
1195       Mat_SeqAIJCUSPARSEMultStruct_Destroy(&matstruct,cusparsestruct->format);
1196       try {
1197         cusparsestruct->nonzerorow=0;
1198         for (int j = 0; j<m; j++) cusparsestruct->nonzerorow += ((a->i[j+1]-a->i[j])>0);
1199 
1200         if (a->compressedrow.use) {
1201           m    = a->compressedrow.nrows;
1202           ii   = a->compressedrow.i;
1203           ridx = a->compressedrow.rindex;
1204         } else {
1205           /* Forcing compressed row on the GPU */
1206           int k=0;
1207           ierr = PetscMalloc1(cusparsestruct->nonzerorow+1, &ii);CHKERRQ(ierr);
1208           ierr = PetscMalloc1(cusparsestruct->nonzerorow, &ridx);CHKERRQ(ierr);
1209           ii[0]=0;
1210           for (int j = 0; j<m; j++) {
1211             if ((a->i[j+1]-a->i[j])>0) {
1212               ii[k]  = a->i[j];
1213               ridx[k]= j;
1214               k++;
1215             }
1216           }
1217           ii[cusparsestruct->nonzerorow] = a->nz;
1218           m = cusparsestruct->nonzerorow;
1219         }
1220 
1221         /* allocate space for the triangular factor information */
1222         matstruct = new Mat_SeqAIJCUSPARSEMultStruct;
1223         stat = cusparseCreateMatDescr(&matstruct->descr);CHKERRCUDA(stat);
1224         stat = cusparseSetMatIndexBase(matstruct->descr, CUSPARSE_INDEX_BASE_ZERO);CHKERRCUDA(stat);
1225         stat = cusparseSetMatType(matstruct->descr, CUSPARSE_MATRIX_TYPE_GENERAL);CHKERRCUDA(stat);
1226 
1227         err = cudaMalloc((void **)&(matstruct->alpha),sizeof(PetscScalar));CHKERRCUDA(err);
1228         err = cudaMemcpy(matstruct->alpha,&ALPHA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
1229         err = cudaMalloc((void **)&(matstruct->beta),sizeof(PetscScalar));CHKERRCUDA(err);
1230         err = cudaMemcpy(matstruct->beta,&BETA,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(err);
1231         stat = cusparseSetPointerMode(cusparsestruct->handle, CUSPARSE_POINTER_MODE_DEVICE);CHKERRCUDA(stat);
1232 
1233         /* Build a hybrid/ellpack matrix if this option is chosen for the storage */
1234         if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1235           /* set the matrix */
1236           CsrMatrix *matrix= new CsrMatrix;
1237           matrix->num_rows = m;
1238           matrix->num_cols = A->cmap->n;
1239           matrix->num_entries = a->nz;
1240           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1241           matrix->row_offsets->assign(ii, ii + m+1);
1242 
1243           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1244           matrix->column_indices->assign(a->j, a->j+a->nz);
1245 
1246           matrix->values = new THRUSTARRAY(a->nz);
1247           matrix->values->assign(a->a, a->a+a->nz);
1248 
1249           /* assign the pointer */
1250           matstruct->mat = matrix;
1251 
1252         } else if (cusparsestruct->format==MAT_CUSPARSE_ELL || cusparsestruct->format==MAT_CUSPARSE_HYB) {
1253 #if CUDA_VERSION>=4020
1254           CsrMatrix *matrix= new CsrMatrix;
1255           matrix->num_rows = m;
1256           matrix->num_cols = A->cmap->n;
1257           matrix->num_entries = a->nz;
1258           matrix->row_offsets = new THRUSTINTARRAY32(m+1);
1259           matrix->row_offsets->assign(ii, ii + m+1);
1260 
1261           matrix->column_indices = new THRUSTINTARRAY32(a->nz);
1262           matrix->column_indices->assign(a->j, a->j+a->nz);
1263 
1264           matrix->values = new THRUSTARRAY(a->nz);
1265           matrix->values->assign(a->a, a->a+a->nz);
1266 
1267           cusparseHybMat_t hybMat;
1268           stat = cusparseCreateHybMat(&hybMat);CHKERRCUDA(stat);
1269           cusparseHybPartition_t partition = cusparsestruct->format==MAT_CUSPARSE_ELL ?
1270             CUSPARSE_HYB_PARTITION_MAX : CUSPARSE_HYB_PARTITION_AUTO;
1271           stat = cusparse_csr2hyb(cusparsestruct->handle, matrix->num_rows, matrix->num_cols,
1272               matstruct->descr, matrix->values->data().get(),
1273               matrix->row_offsets->data().get(),
1274               matrix->column_indices->data().get(),
1275               hybMat, 0, partition);CHKERRCUDA(stat);
1276           /* assign the pointer */
1277           matstruct->mat = hybMat;
1278 
1279           if (matrix) {
1280             if (matrix->values) delete (THRUSTARRAY*)matrix->values;
1281             if (matrix->column_indices) delete (THRUSTINTARRAY32*)matrix->column_indices;
1282             if (matrix->row_offsets) delete (THRUSTINTARRAY32*)matrix->row_offsets;
1283             delete (CsrMatrix*)matrix;
1284           }
1285 #endif
1286         }
1287 
1288         /* assign the compressed row indices */
1289         matstruct->cprowIndices = new THRUSTINTARRAY(m);
1290         matstruct->cprowIndices->assign(ridx,ridx+m);
1291 
1292         /* assign the pointer */
1293         cusparsestruct->mat = matstruct;
1294 
1295         if (!a->compressedrow.use) {
1296           ierr = PetscFree(ii);CHKERRQ(ierr);
1297           ierr = PetscFree(ridx);CHKERRQ(ierr);
1298         }
1299         cusparsestruct->workVector = new THRUSTARRAY;
1300         cusparsestruct->workVector->resize(m);
1301       } catch(char *ex) {
1302         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1303       }
1304       cusparsestruct->nonzerostate = A->nonzerostate;
1305     }
1306     ierr = WaitForGPU();CHKERRCUDA(ierr);
1307     A->valid_GPU_matrix = PETSC_CUDA_BOTH;
1308     ierr = PetscLogEventEnd(MAT_CUSPARSECopyToGPU,A,0,0,0);CHKERRQ(ierr);
1309   }
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 static PetscErrorCode MatCreateVecs_SeqAIJCUSPARSE(Mat mat, Vec *right, Vec *left)
1314 {
1315   PetscErrorCode ierr;
1316   PetscInt rbs,cbs;
1317 
1318   PetscFunctionBegin;
1319   ierr = MatGetBlockSizes(mat,&rbs,&cbs);CHKERRQ(ierr);
1320   if (right) {
1321     ierr = VecCreate(PetscObjectComm((PetscObject)mat),right);CHKERRQ(ierr);
1322     ierr = VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);CHKERRQ(ierr);
1323     ierr = VecSetBlockSize(*right,cbs);CHKERRQ(ierr);
1324     ierr = VecSetType(*right,VECSEQCUDA);CHKERRQ(ierr);
1325     ierr = PetscLayoutReference(mat->cmap,&(*right)->map);CHKERRQ(ierr);
1326   }
1327   if (left) {
1328     ierr = VecCreate(PetscObjectComm((PetscObject)mat),left);CHKERRQ(ierr);
1329     ierr = VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);CHKERRQ(ierr);
1330     ierr = VecSetBlockSize(*left,rbs);CHKERRQ(ierr);
1331     ierr = VecSetType(*left,VECSEQCUDA);CHKERRQ(ierr);
1332     ierr = PetscLayoutReference(mat->rmap,&(*left)->map);CHKERRQ(ierr);
1333   }
1334   PetscFunctionReturn(0);
1335 }
1336 
1337 struct VecCUDAPlusEquals
1338 {
1339   template <typename Tuple>
1340   __host__ __device__
1341   void operator()(Tuple t)
1342   {
1343     thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t);
1344   }
1345 };
1346 
1347 static PetscErrorCode MatMult_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
1348 {
1349   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1350   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1351   Mat_SeqAIJCUSPARSEMultStruct *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1352   const PetscScalar            *xarray;
1353   PetscScalar                  *yarray;
1354   PetscErrorCode               ierr;
1355   cusparseStatus_t             stat;
1356 
1357   PetscFunctionBegin;
1358   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
1359   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1360   ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1361   ierr = VecSet(yy,0);CHKERRQ(ierr);
1362   ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1363   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1364     CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1365     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1366                              mat->num_rows, mat->num_cols, mat->num_entries,
1367                              matstruct->alpha, matstruct->descr, mat->values->data().get(), mat->row_offsets->data().get(),
1368                              mat->column_indices->data().get(), xarray, matstruct->beta,
1369                              yarray);CHKERRCUDA(stat);
1370   } else {
1371 #if CUDA_VERSION>=4020
1372     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1373     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1374                              matstruct->alpha, matstruct->descr, hybMat,
1375                              xarray, matstruct->beta,
1376                              yarray);CHKERRCUDA(stat);
1377 #endif
1378   }
1379   ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1380   ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1381   if (!cusparsestruct->stream) {
1382     ierr = WaitForGPU();CHKERRCUDA(ierr);
1383   }
1384   ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
1385   PetscFunctionReturn(0);
1386 }
1387 
1388 static PetscErrorCode MatMultTranspose_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy)
1389 {
1390   Mat_SeqAIJ                   *a = (Mat_SeqAIJ*)A->data;
1391   Mat_SeqAIJCUSPARSE           *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1392   Mat_SeqAIJCUSPARSEMultStruct *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1393   const PetscScalar            *xarray;
1394   PetscScalar                  *yarray;
1395   PetscErrorCode               ierr;
1396   cusparseStatus_t             stat;
1397 
1398   PetscFunctionBegin;
1399   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
1400   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1401   if (!matstructT) {
1402     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1403     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1404   }
1405   ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1406   ierr = VecSet(yy,0);CHKERRQ(ierr);
1407   ierr = VecCUDAGetArrayWrite(yy,&yarray);CHKERRQ(ierr);
1408 
1409   if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1410     CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1411     stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1412                              mat->num_rows, mat->num_cols,
1413                              mat->num_entries, matstructT->alpha, matstructT->descr,
1414                              mat->values->data().get(), mat->row_offsets->data().get(),
1415                              mat->column_indices->data().get(), xarray, matstructT->beta,
1416                              yarray);CHKERRCUDA(stat);
1417   } else {
1418 #if CUDA_VERSION>=4020
1419     cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1420     stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1421                              matstructT->alpha, matstructT->descr, hybMat,
1422                              xarray, matstructT->beta,
1423                              yarray);CHKERRCUDA(stat);
1424 #endif
1425   }
1426   ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1427   ierr = VecCUDARestoreArrayWrite(yy,&yarray);CHKERRQ(ierr);
1428   if (!cusparsestruct->stream) {
1429     ierr = WaitForGPU();CHKERRCUDA(ierr);
1430   }
1431   ierr = PetscLogFlops(2.0*a->nz - cusparsestruct->nonzerorow);CHKERRQ(ierr);
1432   PetscFunctionReturn(0);
1433 }
1434 
1435 
1436 static PetscErrorCode MatMultAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
1437 {
1438   Mat_SeqAIJ                      *a = (Mat_SeqAIJ*)A->data;
1439   Mat_SeqAIJCUSPARSE              *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1440   Mat_SeqAIJCUSPARSEMultStruct    *matstruct = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->mat;
1441   thrust::device_ptr<PetscScalar> zptr;
1442   const PetscScalar               *xarray;
1443   PetscScalar                     *zarray;
1444   PetscErrorCode                  ierr;
1445   cusparseStatus_t                stat;
1446 
1447   PetscFunctionBegin;
1448   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
1449   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1450   try {
1451     ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
1452     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1453     ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr);
1454     zptr = thrust::device_pointer_cast(zarray);
1455 
1456     /* multiply add */
1457     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1458       CsrMatrix *mat = (CsrMatrix*)matstruct->mat;
1459     /* here we need to be careful to set the number of rows in the multiply to the
1460        number of compressed rows in the matrix ... which is equivalent to the
1461        size of the workVector */
1462       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1463                                mat->num_rows, mat->num_cols,
1464                                mat->num_entries, matstruct->alpha, matstruct->descr,
1465                                mat->values->data().get(), mat->row_offsets->data().get(),
1466                                mat->column_indices->data().get(), xarray, matstruct->beta,
1467                                cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1468     } else {
1469 #if CUDA_VERSION>=4020
1470       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstruct->mat;
1471       if (cusparsestruct->workVector->size()) {
1472         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1473             matstruct->alpha, matstruct->descr, hybMat,
1474             xarray, matstruct->beta,
1475             cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1476       }
1477 #endif
1478     }
1479 
1480     /* scatter the data from the temporary into the full vector with a += operation */
1481     thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))),
1482         thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstruct->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1483         VecCUDAPlusEquals());
1484     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1485     ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr);
1486 
1487   } catch(char *ex) {
1488     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1489   }
1490   ierr = WaitForGPU();CHKERRCUDA(ierr);
1491   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1492   PetscFunctionReturn(0);
1493 }
1494 
1495 static PetscErrorCode MatMultTransposeAdd_SeqAIJCUSPARSE(Mat A,Vec xx,Vec yy,Vec zz)
1496 {
1497   Mat_SeqAIJ                      *a = (Mat_SeqAIJ*)A->data;
1498   Mat_SeqAIJCUSPARSE              *cusparsestruct = (Mat_SeqAIJCUSPARSE*)A->spptr;
1499   Mat_SeqAIJCUSPARSEMultStruct    *matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1500   thrust::device_ptr<PetscScalar> zptr;
1501   const PetscScalar               *xarray;
1502   PetscScalar                     *zarray;
1503   PetscErrorCode                  ierr;
1504   cusparseStatus_t                stat;
1505 
1506   PetscFunctionBegin;
1507   /* The line below is necessary due to the operations that modify the matrix on the CPU (axpy, scale, etc) */
1508   ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1509   if (!matstructT) {
1510     ierr = MatSeqAIJCUSPARSEGenerateTransposeForMult(A);CHKERRQ(ierr);
1511     matstructT = (Mat_SeqAIJCUSPARSEMultStruct*)cusparsestruct->matTranspose;
1512   }
1513 
1514   try {
1515     ierr = VecCopy_SeqCUDA(yy,zz);CHKERRQ(ierr);
1516     ierr = VecCUDAGetArrayRead(xx,&xarray);CHKERRQ(ierr);
1517     ierr = VecCUDAGetArrayReadWrite(zz,&zarray);CHKERRQ(ierr);
1518     zptr = thrust::device_pointer_cast(zarray);
1519 
1520     /* multiply add with matrix transpose */
1521     if (cusparsestruct->format==MAT_CUSPARSE_CSR) {
1522       CsrMatrix *mat = (CsrMatrix*)matstructT->mat;
1523       /* here we need to be careful to set the number of rows in the multiply to the
1524          number of compressed rows in the matrix ... which is equivalent to the
1525          size of the workVector */
1526       stat = cusparse_csr_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1527                                mat->num_rows, mat->num_cols,
1528                                mat->num_entries, matstructT->alpha, matstructT->descr,
1529                                mat->values->data().get(), mat->row_offsets->data().get(),
1530                                mat->column_indices->data().get(), xarray, matstructT->beta,
1531                                cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1532     } else {
1533 #if CUDA_VERSION>=4020
1534       cusparseHybMat_t hybMat = (cusparseHybMat_t)matstructT->mat;
1535       if (cusparsestruct->workVector->size()) {
1536         stat = cusparse_hyb_spmv(cusparsestruct->handle, CUSPARSE_OPERATION_NON_TRANSPOSE,
1537             matstructT->alpha, matstructT->descr, hybMat,
1538             xarray, matstructT->beta,
1539             cusparsestruct->workVector->data().get());CHKERRCUDA(stat);
1540       }
1541 #endif
1542     }
1543 
1544     /* scatter the data from the temporary into the full vector with a += operation */
1545     thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))),
1546         thrust::make_zip_iterator(thrust::make_tuple(cusparsestruct->workVector->begin(), thrust::make_permutation_iterator(zptr, matstructT->cprowIndices->begin()))) + cusparsestruct->workVector->size(),
1547         VecCUDAPlusEquals());
1548 
1549     ierr = VecCUDARestoreArrayRead(xx,&xarray);CHKERRQ(ierr);
1550     ierr = VecCUDARestoreArrayReadWrite(zz,&zarray);CHKERRQ(ierr);
1551 
1552   } catch(char *ex) {
1553     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSPARSE error: %s", ex);
1554   }
1555   ierr = WaitForGPU();CHKERRCUDA(ierr);
1556   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 static PetscErrorCode MatAssemblyEnd_SeqAIJCUSPARSE(Mat A,MatAssemblyType mode)
1561 {
1562   PetscErrorCode ierr;
1563 
1564   PetscFunctionBegin;
1565   ierr = MatAssemblyEnd_SeqAIJ(A,mode);CHKERRQ(ierr);
1566   if (A->factortype==MAT_FACTOR_NONE) {
1567     ierr = MatSeqAIJCUSPARSECopyToGPU(A);CHKERRQ(ierr);
1568   }
1569   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
1570   A->ops->mult             = MatMult_SeqAIJCUSPARSE;
1571   A->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1572   A->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1573   A->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
1574   PetscFunctionReturn(0);
1575 }
1576 
1577 /* --------------------------------------------------------------------------------*/
1578 /*@
1579    MatCreateSeqAIJCUSPARSE - Creates a sparse matrix in AIJ (compressed row) format
1580    (the default parallel PETSc format). This matrix will ultimately pushed down
1581    to NVidia GPUs and use the CUSPARSE library for calculations. For good matrix
1582    assembly performance the user should preallocate the matrix storage by setting
1583    the parameter nz (or the array nnz).  By setting these parameters accurately,
1584    performance during matrix assembly can be increased by more than a factor of 50.
1585 
1586    Collective on MPI_Comm
1587 
1588    Input Parameters:
1589 +  comm - MPI communicator, set to PETSC_COMM_SELF
1590 .  m - number of rows
1591 .  n - number of columns
1592 .  nz - number of nonzeros per row (same for all rows)
1593 -  nnz - array containing the number of nonzeros in the various rows
1594          (possibly different for each row) or NULL
1595 
1596    Output Parameter:
1597 .  A - the matrix
1598 
1599    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
1600    MatXXXXSetPreallocation() paradgm instead of this routine directly.
1601    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
1602 
1603    Notes:
1604    If nnz is given then nz is ignored
1605 
1606    The AIJ format (also called the Yale sparse matrix format or
1607    compressed row storage), is fully compatible with standard Fortran 77
1608    storage.  That is, the stored row and column indices can begin at
1609    either one (as in Fortran) or zero.  See the users' manual for details.
1610 
1611    Specify the preallocated storage with either nz or nnz (not both).
1612    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
1613    allocation.  For large problems you MUST preallocate memory or you
1614    will get TERRIBLE performance, see the users' manual chapter on matrices.
1615 
1616    By default, this format uses inodes (identical nodes) when possible, to
1617    improve numerical efficiency of matrix-vector products and solves. We
1618    search for consecutive rows with the same nonzero structure, thereby
1619    reusing matrix information to achieve increased efficiency.
1620 
1621    Level: intermediate
1622 
1623 .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATSEQAIJCUSPARSE, MATAIJCUSPARSE
1624 @*/
1625 PetscErrorCode  MatCreateSeqAIJCUSPARSE(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
1626 {
1627   PetscErrorCode ierr;
1628 
1629   PetscFunctionBegin;
1630   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1631   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1632   ierr = MatSetType(*A,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
1633   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr);
1634   PetscFunctionReturn(0);
1635 }
1636 
1637 static PetscErrorCode MatDestroy_SeqAIJCUSPARSE(Mat A)
1638 {
1639   PetscErrorCode   ierr;
1640 
1641   PetscFunctionBegin;
1642   if (A->factortype==MAT_FACTOR_NONE) {
1643     if (A->valid_GPU_matrix != PETSC_CUDA_UNALLOCATED) {
1644       ierr = Mat_SeqAIJCUSPARSE_Destroy((Mat_SeqAIJCUSPARSE**)&A->spptr);CHKERRQ(ierr);
1645     }
1646   } else {
1647     ierr = Mat_SeqAIJCUSPARSETriFactors_Destroy((Mat_SeqAIJCUSPARSETriFactors**)&A->spptr);CHKERRQ(ierr);
1648   }
1649   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1650   PetscFunctionReturn(0);
1651 }
1652 
1653 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat B)
1654 {
1655   PetscErrorCode ierr;
1656   cusparseStatus_t stat;
1657   cusparseHandle_t handle=0;
1658 
1659   PetscFunctionBegin;
1660   ierr = MatCreate_SeqAIJ(B);CHKERRQ(ierr);
1661   if (B->factortype==MAT_FACTOR_NONE) {
1662     /* you cannot check the inode.use flag here since the matrix was just created.
1663        now build a GPU matrix data structure */
1664     B->spptr = new Mat_SeqAIJCUSPARSE;
1665     ((Mat_SeqAIJCUSPARSE*)B->spptr)->mat          = 0;
1666     ((Mat_SeqAIJCUSPARSE*)B->spptr)->matTranspose = 0;
1667     ((Mat_SeqAIJCUSPARSE*)B->spptr)->workVector   = 0;
1668     ((Mat_SeqAIJCUSPARSE*)B->spptr)->format       = MAT_CUSPARSE_CSR;
1669     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream       = 0;
1670     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle       = 0;
1671     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1672     ((Mat_SeqAIJCUSPARSE*)B->spptr)->handle       = handle;
1673     ((Mat_SeqAIJCUSPARSE*)B->spptr)->stream       = 0;
1674   } else {
1675     /* NEXT, set the pointers to the triangular factors */
1676     B->spptr = new Mat_SeqAIJCUSPARSETriFactors;
1677     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtr          = 0;
1678     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtr          = 0;
1679     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->loTriFactorPtrTranspose = 0;
1680     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->upTriFactorPtrTranspose = 0;
1681     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->rpermIndices            = 0;
1682     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->cpermIndices            = 0;
1683     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->workVector              = 0;
1684     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = 0;
1685     stat = cusparseCreate(&handle);CHKERRCUDA(stat);
1686     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->handle                  = handle;
1687     ((Mat_SeqAIJCUSPARSETriFactors*)B->spptr)->nnz                     = 0;
1688   }
1689 
1690   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJCUSPARSE;
1691   B->ops->destroy          = MatDestroy_SeqAIJCUSPARSE;
1692   B->ops->getvecs          = MatCreateVecs_SeqAIJCUSPARSE;
1693   B->ops->setfromoptions   = MatSetFromOptions_SeqAIJCUSPARSE;
1694   B->ops->mult             = MatMult_SeqAIJCUSPARSE;
1695   B->ops->multadd          = MatMultAdd_SeqAIJCUSPARSE;
1696   B->ops->multtranspose    = MatMultTranspose_SeqAIJCUSPARSE;
1697   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJCUSPARSE;
1698 
1699   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSPARSE);CHKERRQ(ierr);
1700 
1701   B->valid_GPU_matrix = PETSC_CUDA_UNALLOCATED;
1702 
1703   ierr = PetscObjectComposeFunction((PetscObject)B, "MatCUSPARSESetFormat_C", MatCUSPARSESetFormat_SeqAIJCUSPARSE);CHKERRQ(ierr);
1704   PetscFunctionReturn(0);
1705 }
1706 
1707 /*M
1708    MATSEQAIJCUSPARSE - MATAIJCUSPARSE = "(seq)aijcusparse" - A matrix type to be used for sparse matrices.
1709 
1710    A matrix type type whose data resides on Nvidia GPUs. These matrices can be in either
1711    CSR, ELL, or Hybrid format. The ELL and HYB formats require CUDA 4.2 or later.
1712    All matrix calculations are performed on Nvidia GPUs using the CUSPARSE library.
1713 
1714    Options Database Keys:
1715 +  -mat_type aijcusparse - sets the matrix type to "seqaijcusparse" during a call to MatSetFromOptions()
1716 .  -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).
1717 .  -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).
1718 
1719   Level: beginner
1720 
1721 .seealso: MatCreateSeqAIJCUSPARSE(), MATAIJCUSPARSE, MatCreateAIJCUSPARSE(), MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat, MatCUSPARSEFormatOperation
1722 M*/
1723 
1724 PETSC_EXTERN PetscErrorCode MatGetFactor_seqaijcusparse_cusparse(Mat,MatFactorType,Mat*);
1725 
1726 
1727 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_CUSPARSE(void)
1728 {
1729   PetscErrorCode ierr;
1730 
1731   PetscFunctionBegin;
1732   ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_LU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
1733   ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_CHOLESKY,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
1734   ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ILU,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
1735   ierr = MatSolverPackageRegister(MATSOLVERCUSPARSE,MATSEQAIJCUSPARSE,MAT_FACTOR_ICC,MatGetFactor_seqaijcusparse_cusparse);CHKERRQ(ierr);
1736   PetscFunctionReturn(0);
1737 }
1738 
1739 
1740 static PetscErrorCode Mat_SeqAIJCUSPARSE_Destroy(Mat_SeqAIJCUSPARSE **cusparsestruct)
1741 {
1742   cusparseStatus_t stat;
1743   cusparseHandle_t handle;
1744 
1745   PetscFunctionBegin;
1746   if (*cusparsestruct) {
1747     Mat_SeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->mat,(*cusparsestruct)->format);
1748     Mat_SeqAIJCUSPARSEMultStruct_Destroy(&(*cusparsestruct)->matTranspose,(*cusparsestruct)->format);
1749     delete (*cusparsestruct)->workVector;
1750     if (handle = (*cusparsestruct)->handle) {
1751       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
1752     }
1753     delete *cusparsestruct;
1754     *cusparsestruct = 0;
1755   }
1756   PetscFunctionReturn(0);
1757 }
1758 
1759 static PetscErrorCode CsrMatrix_Destroy(CsrMatrix **mat)
1760 {
1761   PetscFunctionBegin;
1762   if (*mat) {
1763     delete (*mat)->values;
1764     delete (*mat)->column_indices;
1765     delete (*mat)->row_offsets;
1766     delete *mat;
1767     *mat = 0;
1768   }
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(Mat_SeqAIJCUSPARSETriFactorStruct **trifactor)
1773 {
1774   cusparseStatus_t stat;
1775   PetscErrorCode   ierr;
1776 
1777   PetscFunctionBegin;
1778   if (*trifactor) {
1779     if ((*trifactor)->descr) { stat = cusparseDestroyMatDescr((*trifactor)->descr);CHKERRCUDA(stat); }
1780     if ((*trifactor)->solveInfo) { stat = cusparseDestroySolveAnalysisInfo((*trifactor)->solveInfo);CHKERRCUDA(stat); }
1781     ierr = CsrMatrix_Destroy(&(*trifactor)->csrMat);CHKERRQ(ierr);
1782     delete *trifactor;
1783     *trifactor = 0;
1784   }
1785   PetscFunctionReturn(0);
1786 }
1787 
1788 static PetscErrorCode Mat_SeqAIJCUSPARSEMultStruct_Destroy(Mat_SeqAIJCUSPARSEMultStruct **matstruct,MatCUSPARSEStorageFormat format)
1789 {
1790   CsrMatrix        *mat;
1791   cusparseStatus_t stat;
1792   cudaError_t      err;
1793 
1794   PetscFunctionBegin;
1795   if (*matstruct) {
1796     if ((*matstruct)->mat) {
1797       if (format==MAT_CUSPARSE_ELL || format==MAT_CUSPARSE_HYB) {
1798         cusparseHybMat_t hybMat = (cusparseHybMat_t)(*matstruct)->mat;
1799         stat = cusparseDestroyHybMat(hybMat);CHKERRCUDA(stat);
1800       } else {
1801         mat = (CsrMatrix*)(*matstruct)->mat;
1802         CsrMatrix_Destroy(&mat);
1803       }
1804     }
1805     if ((*matstruct)->descr) { stat = cusparseDestroyMatDescr((*matstruct)->descr);CHKERRCUDA(stat); }
1806     delete (*matstruct)->cprowIndices;
1807     if ((*matstruct)->alpha) { err=cudaFree((*matstruct)->alpha);CHKERRCUDA(err); }
1808     if ((*matstruct)->beta) { err=cudaFree((*matstruct)->beta);CHKERRCUDA(err); }
1809     delete *matstruct;
1810     *matstruct = 0;
1811   }
1812   PetscFunctionReturn(0);
1813 }
1814 
1815 static PetscErrorCode Mat_SeqAIJCUSPARSETriFactors_Destroy(Mat_SeqAIJCUSPARSETriFactors** trifactors)
1816 {
1817   cusparseHandle_t handle;
1818   cusparseStatus_t stat;
1819 
1820   PetscFunctionBegin;
1821   if (*trifactors) {
1822     Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->loTriFactorPtr);
1823     Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->upTriFactorPtr);
1824     Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->loTriFactorPtrTranspose);
1825     Mat_SeqAIJCUSPARSETriFactorStruct_Destroy(&(*trifactors)->upTriFactorPtrTranspose);
1826     delete (*trifactors)->rpermIndices;
1827     delete (*trifactors)->cpermIndices;
1828     delete (*trifactors)->workVector;
1829     if (handle = (*trifactors)->handle) {
1830       stat = cusparseDestroy(handle);CHKERRCUDA(stat);
1831     }
1832     delete *trifactors;
1833     *trifactors = 0;
1834   }
1835   PetscFunctionReturn(0);
1836 }
1837 
1838