Actual source code: mpiaijviennacl.cxx

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #include <petscconf.h>
  2:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  3:  #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h>

  5: PetscErrorCode  MatMPIAIJSetPreallocation_MPIAIJViennaCL(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
  6: {
  7:   Mat_MPIAIJ     *b = (Mat_MPIAIJ*)B->data;

 11:   PetscLayoutSetUp(B->rmap);
 12:   PetscLayoutSetUp(B->cmap);
 13:   if (!B->preallocated) {
 14:     /* Explicitly create the two MATSEQAIJVIENNACL matrices. */
 15:     MatCreate(PETSC_COMM_SELF,&b->A);
 16:     MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
 17:     MatSetType(b->A,MATSEQAIJVIENNACL);
 18:     PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
 19:     MatCreate(PETSC_COMM_SELF,&b->B);
 20:     MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
 21:     MatSetType(b->B,MATSEQAIJVIENNACL);
 22:     PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
 23:   }
 24:   MatSeqAIJSetPreallocation(b->A,d_nz,d_nnz);
 25:   MatSeqAIJSetPreallocation(b->B,o_nz,o_nnz);
 26:   B->preallocated = PETSC_TRUE;
 27:   return(0);
 28: }

 30: PetscErrorCode MatAssemblyEnd_MPIAIJViennaCL(Mat A,MatAssemblyType mode)
 31: {
 32:   Mat_MPIAIJ     *b = (Mat_MPIAIJ*)A->data;
 34:   PetscBool      v;

 37:   MatAssemblyEnd_MPIAIJ(A,mode);
 38:   PetscObjectTypeCompare((PetscObject)b->lvec,VECSEQVIENNACL,&v);
 39:   if (!v) {
 40:     PetscInt m;
 41:     VecGetSize(b->lvec,&m);
 42:     VecDestroy(&b->lvec);
 43:     VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&b->lvec);
 44:   }
 45:   return(0);
 46: }

 48: PetscErrorCode MatDestroy_MPIAIJViennaCL(Mat A)
 49: {

 53:   MatDestroy_MPIAIJ(A);
 54:   return(0);
 55: }

 57: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJViennaCL(Mat A)
 58: {

 62:   MatCreate_MPIAIJ(A);
 63:   PetscFree(A->defaultvectype);
 64:   PetscStrallocpy(VECVIENNACL,&A->defaultvectype);
 65:   PetscObjectComposeFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",MatMPIAIJSetPreallocation_MPIAIJViennaCL);
 66:   A->ops->assemblyend = MatAssemblyEnd_MPIAIJViennaCL;
 67:   PetscObjectChangeTypeName((PetscObject)A,MATMPIAIJVIENNACL);
 68:   return(0);
 69: }


 72: /*@
 73:    MatCreateAIJViennaCL - Creates a sparse matrix in AIJ (compressed row) format
 74:    (the default parallel PETSc format).  This matrix will ultimately be pushed down
 75:    to GPUs and use the ViennaCL library for calculations. For good matrix
 76:    assembly performance the user should preallocate the matrix storage by setting
 77:    the parameter nz (or the array nnz).  By setting these parameters accurately,
 78:    performance during matrix assembly can be increased substantially.


 81:    Collective on MPI_Comm

 83:    Input Parameters:
 84: +  comm - MPI communicator, set to PETSC_COMM_SELF
 85: .  m - number of rows
 86: .  n - number of columns
 87: .  nz - number of nonzeros per row (same for all rows)
 88: -  nnz - array containing the number of nonzeros in the various rows
 89:          (possibly different for each row) or NULL

 91:    Output Parameter:
 92: .  A - the matrix

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

 98:    Notes:
 99:    If nnz is given then nz is ignored

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

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

111:    Level: intermediate

113: .seealso: MatCreate(), MatCreateAIJ(), MatCreateAIJCUSPARSE(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATMPIAIJVIENNACL, MATAIJVIENNACL
114: @*/
115: PetscErrorCode  MatCreateAIJViennaCL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
116: {
118:   PetscMPIInt    size;

121:   MatCreate(comm,A);
122:   MatSetSizes(*A,m,n,M,N);
123:   MPI_Comm_size(comm,&size);
124:   if (size > 1) {
125:     MatSetType(*A,MATMPIAIJVIENNACL);
126:     MatMPIAIJSetPreallocation(*A,d_nz,d_nnz,o_nz,o_nnz);
127:   } else {
128:     MatSetType(*A,MATSEQAIJVIENNACL);
129:     MatSeqAIJSetPreallocation(*A,d_nz,d_nnz);
130:   }
131:   return(0);
132: }

134: /*MC
135:    MATAIJVIENNACL - MATMPIAIJVIENNACL= "aijviennacl" = "mpiaijviennacl" - A matrix type to be used for sparse matrices.

137:    A matrix type (CSR format) whose data resides on GPUs.
138:    All matrix calculations are performed using the ViennaCL library.

140:    This matrix type is identical to MATSEQAIJVIENNACL when constructed with a single process communicator,
141:    and MATMPIAIJVIENNACL otherwise.  As a result, for single process communicators,
142:    MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
143:    for communicators controlling multiple processes.  It is recommended that you call both of
144:    the above preallocation routines for simplicity.

146:    Options Database Keys:
147: +  -mat_type mpiaijviennacl - sets the matrix type to "mpiaijviennacl" during a call to MatSetFromOptions()

149:   Level: beginner

151:  .seealso: MatCreateAIJViennaCL(), MATSEQAIJVIENNACL, MatCreateSeqAIJVIENNACL()
152: M*/