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 <petsc-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,PetscInt,&d_nnz,lm,PetscInt,&o_nnz);CHKERRQ(ierr); 28 29 ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr); 30 for (i=0;i<lm;i++) { 31 d_nnz[i] = Aa->i[i+1] - Aa->diag[i]; 32 o_nnz[i] = Ba->i[i+1] - Ba->i[i]; 33 } 34 35 ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr); 36 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 37 ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); 38 ierr = MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);CHKERRQ(ierr); 39 ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 40 41 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 42 43 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 44 for (i=rstart;i<rend;i++) { 45 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 46 j = 0; 47 while (cwork[j] < i) { j++; nz--;} 48 ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr); 49 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 50 } 51 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53 54 if (reuse == MAT_REUSE_MATRIX) { 55 ierr = MatHeaderReplace(A,M);CHKERRQ(ierr); 56 } else { 57 *newmat = M; 58 } 59 PetscFunctionReturn(0); 60 } 61 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */ 62 #undef __FUNCT__ 63 #define __FUNCT__ "MatConvert_MPIBAIJ_MPISBAIJ" 64 PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 65 { 66 PetscErrorCode ierr; 67 Mat M; 68 Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data; 69 Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data; 70 PetscInt *d_nnz,*o_nnz; 71 PetscInt i,j,nz; 72 PetscInt m,n,lm,ln; 73 PetscInt rstart,rend; 74 const PetscScalar *vwork; 75 const PetscInt *cwork; 76 PetscInt bs = A->rmap->bs; 77 78 PetscFunctionBegin; 79 ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 80 ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); 81 ierr = PetscMalloc2(lm/bs,PetscInt,&d_nnz,lm/bs,PetscInt,&o_nnz);CHKERRQ(ierr); 82 83 ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr); 84 for (i=0;i<lm/bs;i++) { 85 d_nnz[i] = Aa->i[i+1] - Aa->diag[i]; 86 o_nnz[i] = Ba->i[i+1] - Ba->i[i]; 87 } 88 89 ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr); 90 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 91 ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); 92 ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr); 93 ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 94 95 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 96 97 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 98 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 99 for (i=rstart;i<rend;i++) { 100 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 101 j = 0; 102 ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr); 103 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 104 } 105 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107 108 if (reuse == MAT_REUSE_MATRIX) { 109 ierr = MatHeaderReplace(A,M);CHKERRQ(ierr); 110 } else { 111 *newmat = M; 112 } 113 PetscFunctionReturn(0); 114 } 115 EXTERN_C_END 116