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