1 #define PETSCMAT_DLL 2 3 #include "mpisbaij.h" /*I "petscmat.h" I*/ 4 #include "../src/mat/impls/aij/mpi/mpiaij.h" 5 #include "private/matimpl.h" 6 #include "petscmat.h" 7 8 EXTERN_C_BEGIN 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatConvert_MPIAIJ_MPISBAIJ" 11 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 12 { 13 PetscErrorCode ierr; 14 Mat M; 15 PetscInt *d_nnz,*o_nnz; 16 PetscInt i,j,k,nz; 17 PetscInt m,n,lm,ln; 18 PetscInt rstart,rend; 19 const PetscScalar *vwork; 20 const PetscInt *cwork; 21 22 PetscFunctionBegin; 23 ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 24 ierr = MatGetLocalSize(A,&lm,&ln);CHKERRQ(ierr); 25 ierr = PetscMalloc(lm*sizeof(PetscInt),&d_nnz);CHKERRQ(ierr); 26 ierr = PetscMalloc(lm*sizeof(PetscInt),&o_nnz);CHKERRQ(ierr); 27 28 k = 0; 29 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 30 for(i=rstart;i<rend;i++){ 31 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 32 d_nnz[k] = o_nnz[k] = 0; 33 for(j=0;j<nz;j++){ 34 if(cwork[j] >= rstart && cwork[j] <= rend-1) /* diagonal portion */ 35 d_nnz[k] += 1; 36 else{ 37 /* insert values only in upper triangular portion */ 38 if(cwork[j] >= rend) 39 o_nnz[k] += 1; /* off-diagonal portion */ 40 } 41 } 42 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 43 k++; 44 } 45 46 ierr = MatCreate(((PetscObject)A)->comm,&M);CHKERRQ(ierr); 47 ierr = MatSetSizes(M,lm,ln,m,n);CHKERRQ(ierr); 48 ierr = MatSetType(M,newtype);CHKERRQ(ierr); 49 ierr = MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);CHKERRQ(ierr); 50 51 ierr = PetscFree(d_nnz);CHKERRQ(ierr); 52 ierr = PetscFree(o_nnz);CHKERRQ(ierr); 53 54 for(i=rstart;i<rend;i++){ 55 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 56 for(j=0;j<nz;j++){ 57 if(cwork[j] >= rstart && cwork[j] <= rend-1){ /* diagonal portion */ 58 ierr = MatSetValues(M,1,&i,1,&cwork[j],&vwork[j],INSERT_VALUES);CHKERRQ(ierr); 59 }else{ 60 /* insert values only in upper triangular portion */ 61 if(cwork[j] >= rend){ 62 ierr = MatSetValues(M,1,&i,1,&cwork[j],&vwork[j],INSERT_VALUES);CHKERRQ(ierr); 63 } 64 } 65 } 66 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 67 k++; 68 } 69 ierr = MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70 ierr = MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71 72 if (A->hermitian){ 73 ierr = MatSetOption(M,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 74 } 75 76 if (reuse == MAT_REUSE_MATRIX) { 77 ierr = MatHeaderReplace(A,M);CHKERRQ(ierr); 78 } else 79 *newmat = M; 80 81 PetscFunctionReturn(0); 82 } 83 EXTERN_C_END 84