xref: /petsc/src/mat/impls/aij/mpi/mpiviennacl/mpiaijviennacl.cxx (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 #define PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND 1
2 
3 #include <petscconf.h>
4 #include <../src/mat/impls/aij/mpi/mpiaij.h>   /*I "petscmat.h" I*/
5 #include <../src/mat/impls/aij/seq/seqviennacl/viennaclmatimpl.h>
6 
7 PetscErrorCode  MatMPIAIJSetPreallocation_MPIAIJViennaCL(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
8 {
9   Mat_MPIAIJ     *b = (Mat_MPIAIJ*)B->data;
10 
11   PetscFunctionBegin;
12   PetscCall(PetscLayoutSetUp(B->rmap));
13   PetscCall(PetscLayoutSetUp(B->cmap));
14   if (!B->preallocated) {
15     /* Explicitly create the two MATSEQAIJVIENNACL matrices. */
16     PetscCall(MatCreate(PETSC_COMM_SELF,&b->A));
17     PetscCall(MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n));
18     PetscCall(MatSetType(b->A,MATSEQAIJVIENNACL));
19     PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->A));
20     PetscCall(MatCreate(PETSC_COMM_SELF,&b->B));
21     PetscCall(MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N));
22     PetscCall(MatSetType(b->B,MATSEQAIJVIENNACL));
23     PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->B));
24   }
25   PetscCall(MatSeqAIJSetPreallocation(b->A,d_nz,d_nnz));
26   PetscCall(MatSeqAIJSetPreallocation(b->B,o_nz,o_nnz));
27   B->preallocated = PETSC_TRUE;
28   PetscFunctionReturn(0);
29 }
30 
31 PetscErrorCode MatAssemblyEnd_MPIAIJViennaCL(Mat A,MatAssemblyType mode)
32 {
33   Mat_MPIAIJ     *b = (Mat_MPIAIJ*)A->data;
34   PetscBool      v;
35 
36   PetscFunctionBegin;
37   PetscCall(MatAssemblyEnd_MPIAIJ(A,mode));
38   PetscCall(PetscObjectTypeCompare((PetscObject)b->lvec,VECSEQVIENNACL,&v));
39   if (!v) {
40     PetscInt m;
41     PetscCall(VecGetSize(b->lvec,&m));
42     PetscCall(VecDestroy(&b->lvec));
43     PetscCall(VecCreateSeqViennaCL(PETSC_COMM_SELF,m,&b->lvec));
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 PetscErrorCode MatDestroy_MPIAIJViennaCL(Mat A)
49 {
50   PetscFunctionBegin;
51   PetscCall(MatDestroy_MPIAIJ(A));
52   PetscFunctionReturn(0);
53 }
54 
55 PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJViennaCL(Mat A)
56 {
57   PetscFunctionBegin;
58   PetscCall(MatCreate_MPIAIJ(A));
59   A->boundtocpu = PETSC_FALSE;
60   PetscCall(PetscFree(A->defaultvectype));
61   PetscCall(PetscStrallocpy(VECVIENNACL,&A->defaultvectype));
62   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",MatMPIAIJSetPreallocation_MPIAIJViennaCL));
63   A->ops->assemblyend = MatAssemblyEnd_MPIAIJViennaCL;
64   PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATMPIAIJVIENNACL));
65   PetscFunctionReturn(0);
66 }
67 
68 /*@C
69    MatCreateAIJViennaCL - Creates a sparse matrix in AIJ (compressed row) format
70    (the default parallel PETSc format).  This matrix will ultimately be pushed down
71    to GPUs and use the ViennaCL library for calculations. For good matrix
72    assembly performance the user should preallocate the matrix storage by setting
73    the parameter nz (or the array nnz).  By setting these parameters accurately,
74    performance during matrix assembly can be increased substantially.
75 
76    Collective
77 
78    Input Parameters:
79 +  comm - MPI communicator, set to PETSC_COMM_SELF
80 .  m - number of rows
81 .  n - number of columns
82 .  nz - number of nonzeros per row (same for all rows)
83 -  nnz - array containing the number of nonzeros in the various rows
84          (possibly different for each row) or NULL
85 
86    Output Parameter:
87 .  A - the matrix
88 
89    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
90    MatXXXXSetPreallocation() paradigm instead of this routine directly.
91    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
92 
93    Notes:
94    If nnz is given then nz is ignored
95 
96    The AIJ format (also called the Yale sparse matrix format or
97    compressed row storage), is fully compatible with standard Fortran 77
98    storage.  That is, the stored row and column indices can begin at
99    either one (as in Fortran) or zero.  See the users' manual for details.
100 
101    Specify the preallocated storage with either nz or nnz (not both).
102    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
103    allocation.  For large problems you MUST preallocate memory or you
104    will get TERRIBLE performance, see the users' manual chapter on matrices.
105 
106    Level: intermediate
107 
108 .seealso: MatCreate(), MatCreateAIJ(), MatCreateAIJCUSPARSE(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatCreateAIJ(), MATMPIAIJVIENNACL, MATAIJVIENNACL
109 @*/
110 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)
111 {
112   PetscMPIInt    size;
113 
114   PetscFunctionBegin;
115   PetscCall(MatCreate(comm,A));
116   PetscCall(MatSetSizes(*A,m,n,M,N));
117   PetscCallMPI(MPI_Comm_size(comm,&size));
118   if (size > 1) {
119     PetscCall(MatSetType(*A,MATMPIAIJVIENNACL));
120     PetscCall(MatMPIAIJSetPreallocation(*A,d_nz,d_nnz,o_nz,o_nnz));
121   } else {
122     PetscCall(MatSetType(*A,MATSEQAIJVIENNACL));
123     PetscCall(MatSeqAIJSetPreallocation(*A,d_nz,d_nnz));
124   }
125   PetscFunctionReturn(0);
126 }
127 
128 /*MC
129    MATAIJVIENNACL - MATMPIAIJVIENNACL= "aijviennacl" = "mpiaijviennacl" - A matrix type to be used for sparse matrices.
130 
131    A matrix type (CSR format) whose data resides on GPUs.
132    All matrix calculations are performed using the ViennaCL library.
133 
134    This matrix type is identical to MATSEQAIJVIENNACL when constructed with a single process communicator,
135    and MATMPIAIJVIENNACL otherwise.  As a result, for single process communicators,
136    MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
137    for communicators controlling multiple processes.  It is recommended that you call both of
138    the above preallocation routines for simplicity.
139 
140    Options Database Keys:
141 .  -mat_type mpiaijviennacl - sets the matrix type to "mpiaijviennacl" during a call to MatSetFromOptions()
142 
143   Level: beginner
144 
145  .seealso: MatCreateAIJViennaCL(), MATSEQAIJVIENNACL, MatCreateSeqAIJVIENNACL()
146 M*/
147