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 Mat M; 9 Mat_MPIAIJ *mpimat = (Mat_MPIAIJ *)A->data; 10 Mat_SeqAIJ *Aa = (Mat_SeqAIJ *)mpimat->A->data, *Ba = (Mat_SeqAIJ *)mpimat->B->data; 11 PetscInt *d_nnz, *o_nnz; 12 PetscInt i, m, n, lm, ln, bs = PetscAbs(A->rmap->bs); 13 14 PetscFunctionBegin; 15 if (reuse != MAT_REUSE_MATRIX) { 16 PetscCall(MatGetSize(A, &m, &n)); 17 PetscCall(MatGetLocalSize(A, &lm, &ln)); 18 PetscCall(PetscMalloc2(lm / bs, &d_nnz, lm / bs, &o_nnz)); 19 20 for (i = 0; i < lm / bs; i++) { 21 d_nnz[i] = (Aa->i[i * bs + 1] - Aa->i[i * bs]) / bs; 22 o_nnz[i] = (Ba->i[i * bs + 1] - Ba->i[i * bs]) / bs; 23 } 24 25 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M)); 26 PetscCall(MatSetSizes(M, lm, ln, m, n)); 27 PetscCall(MatSetType(M, MATMPIBAIJ)); 28 PetscCall(MatSeqBAIJSetPreallocation(M, bs, 0, d_nnz)); 29 PetscCall(MatMPIBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz)); 30 PetscCall(PetscFree2(d_nnz, o_nnz)); 31 } else M = *newmat; 32 33 /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */ 34 /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */ 35 /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */ 36 PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M)); 37 38 if (reuse == MAT_INPLACE_MATRIX) { 39 PetscCall(MatHeaderReplace(A, &M)); 40 } else *newmat = M; 41 PetscFunctionReturn(0); 42 } 43