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