Actual source code: aijcusp.cu

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: /*
  2:   Defines the basic matrix operations for the AIJ (compressed row)
  3:   matrix storage format.
  4: */
  5: #define PETSC_SKIP_COMPLEX

  7: #include <petscconf.h>
  8: PETSC_CUDA_EXTERN_C_BEGIN
  9: #include <../src/mat/impls/aij/seq/aij.h>          /*I "petscmat.h" I*/
 10: #include <petscbt.h>
 11: #include <../src/vec/vec/impls/dvecimpl.h>
 12: #include <petsc/private/vecimpl.h>
 13: PETSC_CUDA_EXTERN_C_END
 14: #undef VecType
 15: #include <../src/mat/impls/aij/seq/seqcusp/cuspmatimpl.h>

 17: const char *const MatCUSPStorageFormats[] = {"CSR","DIA","ELL","MatCUSPStorageFormat","MAT_CUSP_",0};

 21: PetscErrorCode MatCUSPSetStream(Mat A,const cudaStream_t stream)
 22: {
 23:   Mat_SeqAIJCUSP *cuspstruct = (Mat_SeqAIJCUSP*)A->spptr;

 26:   cuspstruct->stream = stream;
 27:   return(0);
 28: }

 32: PetscErrorCode MatCUSPSetFormat_SeqAIJCUSP(Mat A,MatCUSPFormatOperation op,MatCUSPStorageFormat format)
 33: {
 34:   Mat_SeqAIJCUSP *cuspMat = (Mat_SeqAIJCUSP*)A->spptr;

 37:   switch (op) {
 38:   case MAT_CUSP_MULT:
 39:     cuspMat->format = format;
 40:     break;
 41:   case MAT_CUSP_ALL:
 42:     cuspMat->format = format;
 43:     break;
 44:   default:
 45:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported operation %d for MatCUSPFormatOperation. Only MAT_CUSP_MULT and MAT_CUSP_ALL are currently supported.",op);
 46:   }
 47:   return(0);
 48: }

 50: /*@
 51:    MatCUSPSetFormat - Sets the storage format of CUSP matrices for a particular
 52:    operation. Only the MatMult operation can use different GPU storage formats
 53:    for AIJCUSP matrices.

 55:    Not Collective

 57:    Input Parameters:
 58: +  A - Matrix of type SEQAIJCUSP
 59: .  op - MatCUSPFormatOperation. SEQAIJCUSP matrices support MAT_CUSP_MULT and MAT_CUSP_ALL. MPIAIJCUSP matrices support MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, and MAT_CUSP_ALL.
 60: -  format - MatCUSPStorageFormat (one of MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL)

 62:    Output Parameter:

 64:    Level: intermediate

 66: .seealso: MatCUSPStorageFormat, MatCUSPFormatOperation
 67: @*/
 70: PetscErrorCode MatCUSPSetFormat(Mat A,MatCUSPFormatOperation op,MatCUSPStorageFormat format)
 71: {

 76:   PetscTryMethod(A, "MatCUSPSetFormat_C",(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat),(A,op,format));
 77:   return(0);
 78: }

 82: PetscErrorCode MatSetFromOptions_SeqAIJCUSP(PetscOptions *PetscOptionsObject,Mat A)
 83: {
 84:   Mat_SeqAIJCUSP       *cuspMat = (Mat_SeqAIJCUSP*)A->spptr;
 85:   PetscErrorCode       ierr;
 86:   MatCUSPStorageFormat format;
 87:   PetscBool            flg;

 90:   PetscOptionsHead(PetscOptionsObject,"SeqAIJCUSP options");
 91:   PetscObjectOptionsBegin((PetscObject)A);
 92:   PetscOptionsEnum("-mat_cusp_mult_storage_format","sets storage format of (seq)aijcusp gpu matrices for SpMV",
 93:                           "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)cuspMat->format,(PetscEnum*)&format,&flg);
 94:   if (flg) {
 95:     MatCUSPSetFormat(A,MAT_CUSP_MULT,format);
 96:   }
 97:   PetscOptionsEnum("-mat_cusp_storage_format","sets storage format of (seq)aijcusp gpu matrices for SpMV",
 98:                           "MatCUSPSetFormat",MatCUSPStorageFormats,(PetscEnum)cuspMat->format,(PetscEnum*)&format,&flg);
 99:   if (flg) {
100:     MatCUSPSetFormat(A,MAT_CUSP_ALL,format);
101:   }
102:   PetscOptionsEnd();
103:   return(0);

105: }


110: PetscErrorCode MatCUSPCopyToGPU(Mat A)
111: {

113:   Mat_SeqAIJCUSP *cuspstruct = (Mat_SeqAIJCUSP*)A->spptr;
114:   Mat_SeqAIJ     *a          = (Mat_SeqAIJ*)A->data;
115:   PetscInt       m           = A->rmap->n,*ii,*ridx;
116:   CUSPMATRIX     *mat;

120:   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED || A->valid_GPU_matrix == PETSC_CUSP_CPU) {
121:     PetscLogEventBegin(MAT_CUSPCopyToGPU,A,0,0,0);
122:     /*
123:       It may be possible to reuse nonzero structure with new matrix values but
124:       for simplicity and insured correctness we delete and build a new matrix on
125:       the GPU. Likely a very small performance hit.
126:     */
127:     if (cuspstruct->mat) {
128:       try {
129:         if (cuspstruct->format==MAT_CUSP_ELL)
130:           delete (CUSPMATRIXELL *) cuspstruct->mat;
131:         else if (cuspstruct->format==MAT_CUSP_DIA)
132:           delete (CUSPMATRIXDIA *) cuspstruct->mat;
133:         else
134:           delete (CUSPMATRIX *) cuspstruct->mat;
135:       } catch(char *ex) {
136:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
137:       }
138:     }
139:     try {
140:       cuspstruct->nonzerorow=0;
141:       for (int j = 0; j<m; j++) cuspstruct->nonzerorow += ((a->i[j+1]-a->i[j])>0);
142:       if (a->compressedrow.use) {
143:         m    = a->compressedrow.nrows;
144:         ii   = a->compressedrow.i;
145:         ridx = a->compressedrow.rindex;
146:       } else {
147:         /* Forcing compressed row on the GPU */
148:         int k=0;
149:         PetscMalloc1(cuspstruct->nonzerorow+1, &ii);
150:         PetscMalloc1(cuspstruct->nonzerorow, &ridx);
151:         ii[0]=0;
152:         for (int j = 0; j<m; j++) {
153:           if ((a->i[j+1]-a->i[j])>0) {
154:             ii[k]  = a->i[j];
155:             ridx[k]= j;
156:             k++;
157:           }
158:         }
159:         ii[cuspstruct->nonzerorow] = a->nz;
160:         m = cuspstruct->nonzerorow;
161:       }

163:       /* now build matrix */
164:       mat = new CUSPMATRIX;
165:       mat->resize(m,A->cmap->n,a->nz);
166:       mat->row_offsets.assign(ii,ii+m+1);
167:       mat->column_indices.assign(a->j,a->j+a->nz);
168:       mat->values.assign(a->a,a->a+a->nz);

170:       /* convert to other formats if selected */
171:       if (a->compressedrow.use || cuspstruct->format==MAT_CUSP_CSR) {
172:         cuspstruct->mat = mat;
173:         cuspstruct->format = MAT_CUSP_CSR;
174:       } else {
175:         if (cuspstruct->format==MAT_CUSP_ELL) {
176:           CUSPMATRIXELL *ellMat = new CUSPMATRIXELL(*mat);
177:           cuspstruct->mat = ellMat;
178:         } else {
179:           CUSPMATRIXDIA *diaMat = new CUSPMATRIXDIA(*mat);
180:           cuspstruct->mat = diaMat;
181:         }
182:         delete (CUSPMATRIX*) mat;
183:       }

185:       /* assign the compressed row indices */
186:       if (cuspstruct->indices) delete (CUSPINTARRAYGPU*)cuspstruct->indices;
187:       cuspstruct->indices = new CUSPINTARRAYGPU;
188:       cuspstruct->indices->assign(ridx,ridx+m);

190:       /* free the temporaries */
191:       if (!a->compressedrow.use) {
192:         PetscFree(ii);
193:         PetscFree(ridx);
194:       }
195:       if (cuspstruct->tempvec) delete (CUSPARRAY*)cuspstruct->tempvec;
196:       cuspstruct->tempvec = new CUSPARRAY;
197:       cuspstruct->tempvec->resize(m);
198:     } catch(char *ex) {
199:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
200:     }
201:     WaitForGPU();CHKERRCUSP(ierr);

203:     A->valid_GPU_matrix = PETSC_CUSP_BOTH;

205:     PetscLogEventEnd(MAT_CUSPCopyToGPU,A,0,0,0);
206:   }
207:   return(0);
208: }

212: PetscErrorCode MatCUSPCopyFromGPU(Mat A, CUSPMATRIX *Agpu)
213: {
214:   Mat_SeqAIJCUSP *cuspstruct = (Mat_SeqAIJCUSP*) A->spptr;
215:   Mat_SeqAIJ     *a          = (Mat_SeqAIJ*) A->data;
216:   PetscInt       m           = A->rmap->n;
217:   CUSPMATRIX *mat;

221:   /* if the data is stored in non-CSR format, create a temporary */
222:   if (cuspstruct->format==MAT_CUSP_ELL) {
223:     mat = new CUSPMATRIX(*((CUSPMATRIXELL*)cuspstruct->mat));
224:   } else if (cuspstruct->format==MAT_CUSP_DIA) {
225:     mat = new CUSPMATRIX(*((CUSPMATRIXDIA*)cuspstruct->mat));
226:   } else {
227:     mat = (CUSPMATRIX*) cuspstruct->mat;
228:   }

230:   if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED) {
231:     if (A->valid_GPU_matrix == PETSC_CUSP_UNALLOCATED) {
232:       try {
233:         mat = Agpu;
234:         if (a->compressedrow.use) {
235:           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Cannot handle row compression for GPU matrices");
236:         } else {
237:           PetscInt i;

239:           if (m+1 != (PetscInt) mat->row_offsets.size()) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "GPU matrix has %d rows, should be %d", mat->row_offsets.size()-1, m);
240:           a->nz           = mat->values.size();
241:           a->maxnz        = a->nz; /* Since we allocate exactly the right amount */
242:           A->preallocated = PETSC_TRUE;
243:           if (a->singlemalloc) {
244:             if (a->a) {PetscFree3(a->a,a->j,a->i);}
245:           } else {
246:             if (a->i) {PetscFree(a->i);}
247:             if (a->j) {PetscFree(a->j);}
248:             if (a->a) {PetscFree(a->a);}
249:           }
250:           PetscMalloc3(a->nz,&a->a,a->nz,&a->j,m+1,&a->i);
251:           PetscLogObjectMemory((PetscObject)A, a->nz*(sizeof(PetscScalar)+sizeof(PetscInt))+(m+1)*sizeof(PetscInt));

253:           a->singlemalloc = PETSC_TRUE;
254:           thrust::copy(mat->row_offsets.begin(), mat->row_offsets.end(), a->i);
255:           thrust::copy(mat->column_indices.begin(), mat->column_indices.end(), a->j);
256:           thrust::copy(mat->values.begin(), mat->values.end(), a->a);
257:           /* Setup row lengths */
258:           if (a->imax) {PetscFree2(a->imax,a->ilen);}
259:           PetscMalloc2(m,&a->imax,m,&a->ilen);
260:           PetscLogObjectMemory((PetscObject)A, 2*m*sizeof(PetscInt));
261:           for (i = 0; i < m; ++i) a->imax[i] = a->ilen[i] = a->i[i+1] - a->i[i];
262:           /* a->diag?*/
263:         }
264:         cuspstruct->tempvec = new CUSPARRAY;
265:         cuspstruct->tempvec->resize(m);
266:       } catch(char *ex) {
267:         SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "CUSP error: %s", ex);
268:       }
269:     }
270:     /* This assembly prevents resetting the flag to PETSC_CUSP_CPU and recopying */
271:     MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
272:     MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

274:     A->valid_GPU_matrix = PETSC_CUSP_BOTH;

276:     /* delete the temporary */
277:     if (cuspstruct->format==MAT_CUSP_ELL || cuspstruct->format==MAT_CUSP_DIA)
278:       delete (CUSPMATRIX*) mat;
279:   } else {
280:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Only valid for unallocated GPU matrices");
281:   }
282:   return(0);
283: }

287: PetscErrorCode MatCreateVecs_SeqAIJCUSP(Mat mat, Vec *right, Vec *left)
288: {
290:   PetscInt rbs,cbs;

293:   MatGetBlockSizes(mat,&rbs,&cbs);
294:   if (right) {
295:     VecCreate(PetscObjectComm((PetscObject)mat),right);
296:     VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);
297:     VecSetBlockSize(*right,cbs);
298:     VecSetType(*right,VECSEQCUSP);
299:     PetscLayoutReference(mat->cmap,&(*right)->map);
300:   }
301:   if (left) {
302:     VecCreate(PetscObjectComm((PetscObject)mat),left);
303:     VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);
304:     VecSetBlockSize(*left,rbs);
305:     VecSetType(*left,VECSEQCUSP);
306:     PetscLayoutReference(mat->rmap,&(*left)->map);
307:   }
308:   return(0);
309: }

313: PetscErrorCode MatMult_SeqAIJCUSP(Mat A,Vec xx,Vec yy)
314: {
315:   Mat_SeqAIJ       *a = (Mat_SeqAIJ*)A->data;
316:   PetscErrorCode   ierr;
317:   Mat_SeqAIJCUSP   *cuspstruct = (Mat_SeqAIJCUSP*)A->spptr;
318:   PetscBool        usecprow = a->compressedrow.use;
319:   CUSPARRAY        *xarray=NULL,*yarray=NULL;
320:   static PetscBool cite = PETSC_FALSE;

323:   /*
324:     DM 1/28/2014: As of version 0.4.0 cusp does not handle the case of
325:     zero matrices well.  It produces segfaults on some platforms.
326:     Therefore we manually check for the case of a zero matrix here.
327:   */
328:   if (a->nz == 0) {
329:     return(0);
330:   }
331:   PetscCitationsRegister("@incollection{msk2013,\n  author = {Victor Minden and Barry F. Smith and Matthew G. Knepley},\n  title = {Preliminary Implementation of {PETSc} Using {GPUs}},\n  booktitle = {GPU Solutions to Multi-scale Problems in Science and Engineering},\n  series = {Lecture Notes in Earth System Sciences},\n  editor = {David A. Yuen and Long Wang and Xuebin Chi and Lennart Johnsson and Wei Ge and Yaolin Shi},\n  publisher = {Springer Berlin Heidelberg},\n  pages = {131--140},\n  year = {2013},\n}\n",&cite);
332:   VecCUSPGetArrayRead(xx,&xarray);
333:   VecCUSPGetArrayWrite(yy,&yarray);
334:   try {
335:     if (usecprow) {
336:       /* use compressed row format */
337:       CUSPMATRIX *mat = (CUSPMATRIX*)cuspstruct->mat;
338:       cusp::multiply(*mat,*xarray,*cuspstruct->tempvec);
339:       VecSet_SeqCUSP(yy,0.0);
340:       thrust::copy(cuspstruct->tempvec->begin(),cuspstruct->tempvec->end(),thrust::make_permutation_iterator(yarray->begin(),cuspstruct->indices->begin()));
341:     } else {
342:       /* do not use compressed row format */
343:       if (cuspstruct->format==MAT_CUSP_ELL) {
344:         CUSPMATRIXELL *mat = (CUSPMATRIXELL*)cuspstruct->mat;
345:         cusp::multiply(*mat,*xarray,*yarray);
346:       } else if (cuspstruct->format==MAT_CUSP_DIA) {
347:         CUSPMATRIXDIA *mat = (CUSPMATRIXDIA*)cuspstruct->mat;
348:         cusp::multiply(*mat,*xarray,*yarray);
349:       } else {
350:         CUSPMATRIX *mat = (CUSPMATRIX*)cuspstruct->mat;
351:         cusp::multiply(*mat,*xarray,*yarray);
352:       }
353:     }

355:   } catch (char *ex) {
356:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
357:   }
358:   VecCUSPRestoreArrayRead(xx,&xarray);
359:   VecCUSPRestoreArrayWrite(yy,&yarray);
360:   if (!cuspstruct->stream) {
361:     WaitForGPU();CHKERRCUSP(ierr);
362:   }
363:   PetscLogFlops(2.0*a->nz - cuspstruct->nonzerorow);
364:   return(0);
365: }


368: struct VecCUSPPlusEquals
369: {
370:   template <typename Tuple>
371:   __host__ __device__
372:   void operator()(Tuple t)
373:   {
374:     thrust::get<1>(t) = thrust::get<1>(t) + thrust::get<0>(t);
375:   }
376: };

380: PetscErrorCode MatMultAdd_SeqAIJCUSP(Mat A,Vec xx,Vec yy,Vec zz)
381: {
382:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
384:   PetscBool      usecprow = a->compressedrow.use;
385:   Mat_SeqAIJCUSP *cuspstruct = (Mat_SeqAIJCUSP*)A->spptr;
386:   CUSPARRAY      *xarray     = NULL,*yarray=NULL,*zarray=NULL;

389:   /*
390:     DM 1/28/2014: As of version 0.4.0 cusp does not handle the case of
391:     zero matrices well.  It produces segfaults on some platforms.
392:     Therefore we manually check for the case of a zero matrix here.
393:   */
394:   if (a->nz == 0) {
395:     return(0);
396:   }
397:   try {
398:     VecCopy_SeqCUSP(yy,zz);
399:     VecCUSPGetArrayRead(xx,&xarray);
400:     VecCUSPGetArrayRead(yy,&yarray);
401:     VecCUSPGetArrayWrite(zz,&zarray);

403:     if (usecprow) {
404:       /* use compressed row format */
405:       CUSPMATRIX *mat = (CUSPMATRIX*)cuspstruct->mat;
406:       cusp::multiply(*mat,*xarray,*cuspstruct->tempvec);
407:       thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->begin(),
408:                                                                     thrust::make_permutation_iterator(zarray->begin(), cuspstruct->indices->begin()))),
409:                        thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->end(),
410:                                                                     thrust::make_permutation_iterator(zarray->end(),cuspstruct->indices->end()))),
411:                        VecCUSPPlusEquals());
412:     } else {

414:       if (cuspstruct->format==MAT_CUSP_ELL) {
415:         CUSPMATRIXELL *mat = (CUSPMATRIXELL*)cuspstruct->mat;
416:         cusp::multiply(*mat,*xarray,*cuspstruct->tempvec);
417:       } else if (cuspstruct->format==MAT_CUSP_DIA) {
418:         CUSPMATRIXDIA *mat = (CUSPMATRIXDIA*)cuspstruct->mat;
419:         cusp::multiply(*mat,*xarray,*cuspstruct->tempvec);
420:       } else {
421:         CUSPMATRIX *mat = (CUSPMATRIX*)cuspstruct->mat;
422:         cusp::multiply(*mat,*xarray,*cuspstruct->tempvec);
423:       }

425:       if (zarray->size() == cuspstruct->indices->size()) {
426:         thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->begin(),zarray->begin())),
427:             thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->end(),zarray->end())),
428:             VecCUSPPlusEquals());
429:       } else {
430:         thrust::for_each(thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->begin(),
431:                 thrust::make_permutation_iterator(zarray->begin(), cuspstruct->indices->begin()))),
432:             thrust::make_zip_iterator(thrust::make_tuple(cuspstruct->tempvec->end(),
433:                 thrust::make_permutation_iterator(zarray->end(),cuspstruct->indices->end()))),
434:             VecCUSPPlusEquals());
435:       }
436:     }
437:     VecCUSPRestoreArrayRead(xx,&xarray);
438:     VecCUSPRestoreArrayRead(yy,&yarray);
439:     VecCUSPRestoreArrayWrite(zz,&zarray);

441:   } catch(char *ex) {
442:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
443:   }
444:   WaitForGPU();CHKERRCUSP(ierr);
445:   PetscLogFlops(2.0*a->nz);
446:   return(0);
447: }

451: PetscErrorCode MatAssemblyEnd_SeqAIJCUSP(Mat A,MatAssemblyType mode)
452: {

456:   MatAssemblyEnd_SeqAIJ(A,mode);
457:   MatCUSPCopyToGPU(A);
458:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);
459:   A->ops->mult    = MatMult_SeqAIJCUSP;
460:   A->ops->multadd = MatMultAdd_SeqAIJCUSP;
461:   return(0);
462: }

464: /* --------------------------------------------------------------------------------*/
467: /*@
468:    MatCreateSeqAIJCUSP - Creates a sparse matrix in AIJ (compressed row) format
469:    (the default parallel PETSc format).  This matrix will ultimately pushed down
470:    to NVidia GPUs and use the CUSP library for calculations. For good matrix
471:    assembly performance the user should preallocate the matrix storage by setting
472:    the parameter nz (or the array nnz).  By setting these parameters accurately,
473:    performance during matrix assembly can be increased by more than a factor of 50.


476:    Collective on MPI_Comm

478:    Input Parameters:
479: +  comm - MPI communicator, set to PETSC_COMM_SELF
480: .  m - number of rows
481: .  n - number of columns
482: .  nz - number of nonzeros per row (same for all rows)
483: -  nnz - array containing the number of nonzeros in the various rows
484:          (possibly different for each row) or NULL

486:    Output Parameter:
487: .  A - the matrix

489:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
490:    MatXXXXSetPreallocation() paradigm instead of this routine directly.
491:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

493:    Notes:
494:    If nnz is given then nz is ignored

496:    The AIJ format (also called the Yale sparse matrix format or
497:    compressed row storage), is fully compatible with standard Fortran 77
498:    storage.  That is, the stored row and column indices can begin at
499:    either one (as in Fortran) or zero.  See the users' manual for details.

501:    Specify the preallocated storage with either nz or nnz (not both).
502:    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
503:    allocation.  For large problems you MUST preallocate memory or you
504:    will get TERRIBLE performance, see the users' manual chapter on matrices.

506:    By default, this format uses inodes (identical nodes) when possible, to
507:    improve numerical efficiency of matrix-vector products and solves. We
508:    search for consecutive rows with the same nonzero structure, thereby
509:    reusing matrix information to achieve increased efficiency.

511:    Level: intermediate

513: .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ()

515: @*/
516: PetscErrorCode  MatCreateSeqAIJCUSP(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
517: {

521:   MatCreate(comm,A);
522:   MatSetSizes(*A,m,n,m,n);
523:   MatSetType(*A,MATSEQAIJCUSP);
524:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
525:   return(0);
526: }

530: PetscErrorCode MatDestroy_SeqAIJCUSP(Mat A)
531: {
533:   Mat_SeqAIJCUSP *cuspcontainer = (Mat_SeqAIJCUSP*)A->spptr;

536:   try {
537:     if (A->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) {
538:       if (cuspcontainer->format==MAT_CUSP_ELL) {
539:         delete (CUSPMATRIXELL*)(cuspcontainer->mat);
540:       } else if (cuspcontainer->format==MAT_CUSP_DIA) {
541:         delete (CUSPMATRIXDIA*)(cuspcontainer->mat);
542:       } else {
543:         delete (CUSPMATRIX*)(cuspcontainer->mat);
544:       }
545:       if (cuspcontainer->indices) delete (CUSPINTARRAYGPU*)cuspcontainer->indices;
546:       if (cuspcontainer->tempvec) delete (CUSPARRAY*)cuspcontainer->tempvec;
547:     }
548:     delete cuspcontainer;
549:     A->valid_GPU_matrix = PETSC_CUSP_UNALLOCATED;
550:   } catch(char *ex) {
551:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
552:   }
553:   /*this next line is because MatDestroy tries to PetscFree spptr if it is not zero, and PetscFree only works if the memory was allocated with PetscNew or PetscMalloc, which don't call the constructor */
554:   A->spptr = 0;
555:   MatDestroy_SeqAIJ(A);
556:   return(0);
557: }

559: extern PetscErrorCode MatSetValuesBatch_SeqAIJCUSP(Mat, PetscInt, PetscInt, PetscInt*,const PetscScalar*);

561: PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage_seqaij_cusparse(Mat,const MatSolverPackage*);

565: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSP(Mat B)
566: {
568:   Mat_SeqAIJ     *aij;

571:   MatCreate_SeqAIJ(B);
572:   aij             = (Mat_SeqAIJ*)B->data;
573:   aij->inode.use  = PETSC_FALSE;
574:   B->ops->mult    = MatMult_SeqAIJCUSP;
575:   B->ops->multadd = MatMultAdd_SeqAIJCUSP;
576:   B->spptr        = new Mat_SeqAIJCUSP;

578:   if (B->factortype==MAT_FACTOR_NONE) {
579:     ((Mat_SeqAIJCUSP*)B->spptr)->mat     = 0;
580:     ((Mat_SeqAIJCUSP*)B->spptr)->tempvec = 0;
581:     ((Mat_SeqAIJCUSP*)B->spptr)->indices = 0;
582:     ((Mat_SeqAIJCUSP*)B->spptr)->nonzerorow = 0;
583:     ((Mat_SeqAIJCUSP*)B->spptr)->stream = 0;
584:     ((Mat_SeqAIJCUSP*)B->spptr)->format = MAT_CUSP_CSR;
585:   }

587:   B->ops->assemblyend    = MatAssemblyEnd_SeqAIJCUSP;
588:   B->ops->destroy        = MatDestroy_SeqAIJCUSP;
589:   B->ops->getvecs        = MatCreateVecs_SeqAIJCUSP;
590:   B->ops->setvaluesbatch = MatSetValuesBatch_SeqAIJCUSP;
591:   B->ops->setfromoptions = MatSetFromOptions_SeqAIJCUSP;

593:   PetscObjectComposeFunction((PetscObject)B,"MatCUSPSetFormat_C", MatCUSPSetFormat_SeqAIJCUSP);
594:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJCUSP);

596:   B->valid_GPU_matrix = PETSC_CUSP_UNALLOCATED;
597:   return(0);
598: }

600: /*M
601:    MATSEQAIJCUSP - MATAIJCUSP = "aijcusp" = "seqaijcusp" - A matrix type to be used for sparse matrices.

603:    A matrix type type whose data resides on Nvidia GPUs. These matrices are in CSR format by
604:    default. All matrix calculations are performed using the CUSP library. DIA and ELL formats are 
605:    also available.

607:    Options Database Keys:
608: +  -mat_type aijcusp - sets the matrix type to "seqaijcusp" during a call to MatSetFromOptions()
609: .  -mat_cusp_storage_format csr - sets the storage format of matrices for MatMult during a call to MatSetFromOptions(). Other options include dia (diagonal) or ell (ellpack).
610: -  -mat_cusp_mult_storage_format csr - sets the storage format of matrices for MatMult during a call to MatSetFromOptions(). Other options include dia (diagonal) or ell (ellpack).

612:   Level: beginner

614: .seealso: MatCreateSeqAIJCUSP(), MATAIJCUSP, MatCreateAIJCUSP(), MatCUSPSetFormat(), MatCUSPStorageFormat, MatCUSPFormatOperation
615: M*/