xref: /petsc/src/mat/impls/baij/mpi/mpiaijbaij.c (revision 7e1a0bbe36d2be40a00a95404ece00db4857f70d)
1 #include <../src/mat/impls/baij/mpi/mpibaij.h> /*I "petscmat.h" I*/
2 #include <../src/mat/impls/aij/mpi/mpiaij.h>
3 #include <petsc/private/matimpl.h>
4 #include <petscmat.h>
5 
6 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat, PetscInt **);
7 
8 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPIBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
9 {
10   Mat         M;
11   Mat_MPIAIJ *mpimat = (Mat_MPIAIJ *)A->data;
12   PetscInt   *d_nnz, *o_nnz;
13   PetscInt    m, n, lm, ln, bs = A->rmap->bs;
14 
15   PetscFunctionBegin;
16   if (reuse != MAT_REUSE_MATRIX) {
17     PetscBool3 sym = A->symmetric, hermitian = A->hermitian, structurally_symmetric = A->structurally_symmetric, spd = A->spd;
18     PetscCall(MatDisAssemble_MPIAIJ(A, PETSC_FALSE));
19     PetscCall(MatGetSize(A, &m, &n));
20     PetscCall(MatGetLocalSize(A, &lm, &ln));
21     PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->A, &d_nnz));
22     PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->B, &o_nnz));
23     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
24     PetscCall(MatSetSizes(M, lm, ln, m, n));
25     PetscCall(MatSetType(M, MATMPIBAIJ));
26     PetscCall(MatSeqBAIJSetPreallocation(M, bs, 0, d_nnz));
27     PetscCall(MatMPIBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
28     PetscCall(PetscFree(d_nnz));
29     PetscCall(PetscFree(o_nnz));
30     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
31     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
32     A->symmetric              = sym;
33     A->hermitian              = hermitian;
34     A->structurally_symmetric = structurally_symmetric;
35     A->spd                    = spd;
36   } else M = *newmat;
37 
38   /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
39   /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
40   /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
41   PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M));
42 
43   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &M));
44   else *newmat = M;
45   PetscFunctionReturn(PETSC_SUCCESS);
46 }
47