Actual source code: mpiaijsbaij.c


  2: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
  3: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  4: #include <petsc/private/matimpl.h>
  5: #include <petscmat.h>

  7: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat, PetscInt **);
  8: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat, PetscInt **);

 10: /* The code is virtually identical to MatConvert_MPIAIJ_MPIBAIJ() */
 11: PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
 12: {
 13:   Mat         M;
 14:   Mat_MPIAIJ *mpimat = (Mat_MPIAIJ *)A->data;
 15:   PetscInt   *d_nnz, *o_nnz;
 16:   PetscInt    m, n, lm, ln, bs = PetscAbs(A->rmap->bs);

 18:   PetscFunctionBegin;
 19:   if (reuse != MAT_REUSE_MATRIX) {
 20:     PetscCall(MatDisAssemble_MPIAIJ(A));
 21:     PetscCall(MatGetSize(A, &m, &n));
 22:     PetscCall(MatGetLocalSize(A, &lm, &ln));
 23:     PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(mpimat->A, &d_nnz));
 24:     PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(mpimat->B, &o_nnz));
 25:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
 26:     PetscCall(MatSetSizes(M, lm, ln, m, n));
 27:     PetscCall(MatSetType(M, MATMPISBAIJ));
 28:     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
 29:     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));
 30:     PetscCall(PetscFree(d_nnz));
 31:     PetscCall(PetscFree(o_nnz));
 32:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 33:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
 34:   } else M = *newmat;

 36:   /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
 37:   /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
 38:   /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
 39:   PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &M));

 41:   if (reuse == MAT_INPLACE_MATRIX) {
 42:     PetscCall(MatHeaderReplace(A, &M));
 43:   } else *newmat = M;
 44:   PetscFunctionReturn(PETSC_SUCCESS);
 45: }

 47: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
 48: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
 49: {
 50:   Mat                M;
 51:   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ *)A->data;
 52:   Mat_SeqBAIJ       *Aa = (Mat_SeqBAIJ *)mpimat->A->data, *Ba = (Mat_SeqBAIJ *)mpimat->B->data;
 53:   PetscInt          *d_nnz, *o_nnz;
 54:   PetscInt           i, nz;
 55:   PetscInt           m, n, lm, ln;
 56:   PetscInt           rstart, rend;
 57:   const PetscScalar *vwork;
 58:   const PetscInt    *cwork;
 59:   PetscInt           bs = A->rmap->bs;

 61:   PetscFunctionBegin;
 62:   if (reuse != MAT_REUSE_MATRIX) {
 63:     PetscCall(MatGetSize(A, &m, &n));
 64:     PetscCall(MatGetLocalSize(A, &lm, &ln));
 65:     PetscCall(PetscMalloc2(lm / bs, &d_nnz, lm / bs, &o_nnz));

 67:     PetscCall(MatMarkDiagonal_SeqBAIJ(mpimat->A));
 68:     for (i = 0; i < lm / bs; i++) {
 69:       d_nnz[i] = Aa->i[i + 1] - Aa->diag[i];
 70:       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
 71:     }

 73:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &M));
 74:     PetscCall(MatSetSizes(M, lm, ln, m, n));
 75:     PetscCall(MatSetType(M, MATMPISBAIJ));
 76:     PetscCall(MatSeqSBAIJSetPreallocation(M, bs, 0, d_nnz));
 77:     PetscCall(MatMPISBAIJSetPreallocation(M, bs, 0, d_nnz, 0, o_nnz));

 79:     PetscCall(PetscFree2(d_nnz, o_nnz));
 80:   } else M = *newmat;

 82:   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
 83:   PetscCall(MatSetOption(M, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
 84:   for (i = rstart; i < rend; i++) {
 85:     PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork));
 86:     PetscCall(MatSetValues(M, 1, &i, nz, cwork, vwork, INSERT_VALUES));
 87:     PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork));
 88:   }
 89:   PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
 90:   PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));

 92:   if (reuse == MAT_INPLACE_MATRIX) {
 93:     PetscCall(MatHeaderReplace(A, &M));
 94:   } else *newmat = M;
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }