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 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(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,j,nz; 15 PetscInt m,n,lm,ln; 16 PetscInt rstart,rend,bs=PetscAbs(A->rmap->bs); 17 const PetscScalar *vwork; 18 const PetscInt *cwork; 19 20 PetscFunctionBegin; 21 if (!A->symmetric && !A->hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric or hermitian. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE) or MatSetOption(mat,MAT_HERMITIAN,PETSC_TRUE)"); 22 if (reuse != MAT_REUSE_MATRIX) { 23 ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 24 ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); 25 ierr = PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr); 26 27 ierr = MatMarkDiagonal_SeqAIJ(mpimat->A);CHKERRQ(ierr); 28 for (i=0; i<lm/bs; i++) { 29 if (Aa->i[i*bs+1] == Aa->diag[i*bs]) { /* misses diagonal entry */ 30 d_nnz[i] = (Aa->i[i*bs+1] - Aa->i[i*bs])/bs; 31 } else { 32 d_nnz[i] = (Aa->i[i*bs+1] - Aa->diag[i*bs])/bs; 33 } 34 o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs; 35 } 36 37 ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr); 38 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 39 ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); 40 ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr); 41 ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 42 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 43 } else M = *newmat; 44 45 if (bs == 1) { 46 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 47 for (i=rstart; i<rend; i++) { 48 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 49 if (nz) { 50 j = 0; 51 while (cwork[j] < i) { 52 j++; nz--; 53 } 54 ierr = MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);CHKERRQ(ierr); 55 } 56 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 57 } 58 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60 } else { 61 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 62 /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */ 63 /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */ 64 /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */ 65 ierr = MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&M);CHKERRQ(ierr); 66 } 67 68 if (reuse == MAT_INPLACE_MATRIX) { 69 ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr); 70 } else *newmat = M; 71 PetscFunctionReturn(0); 72 } 73 /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */ 74 PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 75 { 76 PetscErrorCode ierr; 77 Mat M; 78 Mat_MPIBAIJ *mpimat = (Mat_MPIBAIJ*)A->data; 79 Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data; 80 PetscInt *d_nnz,*o_nnz; 81 PetscInt i,nz; 82 PetscInt m,n,lm,ln; 83 PetscInt rstart,rend; 84 const PetscScalar *vwork; 85 const PetscInt *cwork; 86 PetscInt bs = A->rmap->bs; 87 88 PetscFunctionBegin; 89 if (reuse != MAT_REUSE_MATRIX) { 90 ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 91 ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); 92 ierr = PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);CHKERRQ(ierr); 93 94 ierr = MatMarkDiagonal_SeqBAIJ(mpimat->A);CHKERRQ(ierr); 95 for (i=0; i<lm/bs; i++) { 96 d_nnz[i] = Aa->i[i+1] - Aa->diag[i]; 97 o_nnz[i] = Ba->i[i+1] - Ba->i[i]; 98 } 99 100 ierr = MatCreate(PetscObjectComm((PetscObject)A),&M);CHKERRQ(ierr); 101 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 102 ierr = MatSetType(M,MATMPISBAIJ);CHKERRQ(ierr); 103 ierr = MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);CHKERRQ(ierr); 104 ierr = MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 105 106 ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 107 } else M = *newmat; 108 109 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 110 ierr = MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); 111 for (i=rstart; i<rend; i++) { 112 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 113 ierr = MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 114 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 115 } 116 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 117 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118 119 if (reuse == MAT_INPLACE_MATRIX) { 120 ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr); 121 } else *newmat = M; 122 PetscFunctionReturn(0); 123 } 124