Actual source code: mpiaijbaij.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  2: #include <../src/mat/impls/baij/mpi/mpibaij.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_MPIAIJ_MPIBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
  8: {
  9:   PetscErrorCode    ierr;
 10:   Mat               M;
 11:   Mat_MPIAIJ        *mpimat = (Mat_MPIAIJ*)A->data;
 12:   Mat_SeqAIJ        *Aa     = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data;
 13:   PetscInt          *d_nnz,*o_nnz;
 14:   PetscInt          i,m,n,lm,ln,bs=PetscAbs(A->rmap->bs);

 17:   if (reuse != MAT_REUSE_MATRIX) {
 18:     MatGetSize(A,&m,&n);
 19:     MatGetLocalSize(A,&lm,&ln);
 20:     PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);

 22:     for (i=0; i<lm/bs; i++) {
 23:       d_nnz[i] = (Aa->i[i*bs+1] - Aa->i[i*bs])/bs;
 24:       o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs;
 25:     }

 27:     MatCreate(PetscObjectComm((PetscObject)A),&M);
 28:     MatSetSizes(M,lm,ln,m,n);
 29:     MatSetType(M,MATMPIBAIJ);
 30:     MatSeqBAIJSetPreallocation(M,bs,0,d_nnz);
 31:     MatMPIBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
 32:     PetscFree2(d_nnz,o_nnz);
 33:   } else M = *newmat;

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

 40:   if (reuse == MAT_INPLACE_MATRIX) {
 41:     MatHeaderReplace(A,&M);
 42:   } else *newmat = M;
 43:   return(0);
 44: }