1 2 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> /*I "petscmat.h" I*/ 3 #include <../src/mat/impls/aij/mpi/mpiaij.h> 4 #include <private/matimpl.h> 5 #include <petscmat.h> 6 7 EXTERN_C_BEGIN 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ" 10 PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 11 { 12 PetscErrorCode ierr; 13 Mat M; 14 Mat_MPIAIJ *mpimat = (Mat_MPIAIJ*)A->data; 15 Mat_SeqAIJ *Aa = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data; 16 PetscInt *d_nnz,*o_nnz; 17 PetscInt i,j,nz; 18 PetscInt m,n,lm,ln; 19 PetscInt rstart,rend; 20 const PetscScalar *vwork; 21 const PetscInt *cwork; 22 23 PetscFunctionBegin; 24 if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); 25 ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 26 ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); 27 ierr = PetscMalloc2(lm*sizeof(PetscInt),PetscInt,&d_nnz,lm*sizeof(PetscInt),PetscInt,&o_nnz);CHKERRQ(ierr); 28 29 for(i=0;i<lm;i++){ 30 d_nnz[i] = Aa->i[i+1] - Aa->diag[i]; 31 o_nnz[i] = Ba->i[i+1] - Ba->i[i]; 32 } 33 34 ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr); 35 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 36 ierr = MatSetType(M,newtype);CHKERRQ(ierr); 37 ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 38 39 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 40 41 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 42 for(i=rstart;i<rend;i++){ 43 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 44 j = 0; 45 while (cwork[j] < i){ j++; nz--;} 46 ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr); 47 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 48 } 49 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51 52 if (A->hermitian){ 53 ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 54 } 55 56 if (reuse == MAT_REUSE_MATRIX) { 57 ierr = MatHeaderReplace(A,M);CHKERRQ(ierr); 58 } else { 59 *newmat = M; 60 } 61 PetscFunctionReturn(0); 62 } 63 EXTERN_C_END 64