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 PetscBool flg = PETSC_FALSE; 15 16 PetscFunctionBegin; 17 if (reuse != MAT_REUSE_MATRIX) { 18 PetscBool3 sym = A->symmetric, hermitian = A->hermitian, structurally_symmetric = A->structurally_symmetric, spd = A->spd; 19 PetscCall(MatDisAssemble_MPIAIJ(A, PETSC_FALSE)); 20 PetscCall(MatGetSize(A, &m, &n)); 21 PetscCall(MatGetLocalSize(A, &lm, &ln)); 22 PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->A, &d_nnz)); 23 PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->B, &o_nnz)); 24 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M)); 25 PetscCall(MatSetSizes(M, lm, ln, m, n)); 26 PetscCall(MatSetType(M, MATMPIBAIJ)); 27 PetscCall(MatSeqBAIJSetPreallocation(M, bs, 0, d_nnz)); 28 PetscCall(MatMPIBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz)); 29 PetscCall(PetscFree(d_nnz)); 30 PetscCall(PetscFree(o_nnz)); 31 PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 32 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 33 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 34 PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_FALSE)); 35 A->symmetric = sym; 36 A->hermitian = hermitian; 37 A->structurally_symmetric = structurally_symmetric; 38 A->spd = spd; 39 } else { 40 M = *newmat; 41 PetscCall(MatGetOption(M, MAT_NO_OFF_PROC_ENTRIES, &flg)); 42 } 43 44 /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */ 45 /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */ 46 /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */ 47 PetscCall(MatSetOption(M, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 48 PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M)); 49 PetscCall(MatSetOption(M, MAT_NO_OFF_PROC_ENTRIES, flg)); 50 51 if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &M)); 52 else *newmat = M; 53 PetscFunctionReturn(PETSC_SUCCESS); 54 } 55