#include <../src/mat/impls/sbaij/mpi/mpisbaij.h> /*I "petscmat.h" I*/ #include <../src/mat/impls/aij/mpi/mpiaij.h> #include #include #undef __FUNCT__ #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ" PETSC_EXTERN PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) { PetscErrorCode ierr; Mat M; Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)A->data; Mat_SeqAIJ *Aa = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data; PetscInt *d_nnz,*o_nnz; PetscInt i,j,nz; PetscInt m,n,lm,ln; PetscInt rstart,rend; const PetscScalar *vwork; const PetscInt *cwork; PetscFunctionBegin; if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); ierr = PetscMalloc2(lm,&d_nnz,lm,&o_nnz);CHKERRQ(ierr); ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr); for (i=0; ii[i+1] - Aa->diag[i]; o_nnz[i] = Ba->i[i+1] - Ba->i[i]; } ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr); ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);CHKERRQ(ierr); ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr); ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); for (i=rstart; i April 2011 */ #undef __FUNCT__ #define __FUNCT__ "MatConvert_MPIBAIJ_MPISBAIJ" PETSC_EXTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) { PetscErrorCode ierr; Mat M; Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data; Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data; PetscInt *d_nnz,*o_nnz; PetscInt i,j,nz; PetscInt m,n,lm,ln; PetscInt rstart,rend; const PetscScalar *vwork; const PetscInt *cwork; PetscInt bs = A->rmap->bs; PetscFunctionBegin; ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); ierr = PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr); ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr); for (i=0; ii[i+1] - Aa->diag[i]; o_nnz[i] = Ba->i[i+1] - Ba->i[i]; } ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr); ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr); ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr); ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); for (i=rstart; i