Actual source code: mpiaijsbaij.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  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>

  9: PETSC_EXTERN PetscErrorCode MatConvert_MPIAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
 10: {
 11:   PetscErrorCode    ierr;
 12:   Mat               M;
 13:   Mat_MPIAIJ        *mpimat = (Mat_MPIAIJ*)A->data;
 14:   Mat_SeqAIJ        *Aa     = (Mat_SeqAIJ*)mpimat->A->data,*Ba = (Mat_SeqAIJ*)mpimat->B->data;
 15:   PetscInt          *d_nnz,*o_nnz;
 16:   PetscInt          i,j,nz;
 17:   PetscInt          m,n,lm,ln;
 18:   PetscInt          rstart,rend;
 19:   const PetscScalar *vwork;
 20:   const PetscInt    *cwork;

 23:   if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
 24:   MatGetSize(A,&m,&n);
 25:   MatGetLocalSize(A,&lm,&ln);
 26:   PetscMalloc2(lm,&d_nnz,lm,&o_nnz);

 28:   MatMarkDiagonal_SeqAIJ(mpimat->A);
 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:   }

 34:   MatCreate(PetscObjectComm((PetscObject)A),&M);
 35:   MatSetSizes(M,lm,ln,m,n);
 36:   MatSetType(M,MATMPISBAIJ);
 37:   MatSeqSBAIJSetPreallocation(M,1,0,d_nnz);
 38:   MatMPISBAIJSetPreallocation(M,1,0,d_nnz,0,o_nnz);

 40:   PetscFree2(d_nnz,o_nnz);

 42:   MatGetOwnershipRange(A,&rstart,&rend);
 43:   for (i=rstart; i<rend; i++) {
 44:     MatGetRow(A,i,&nz,&cwork,&vwork);
 45:     j    = 0;
 46:     while (cwork[j] < i) {
 47:       j++; nz--;
 48:     }
 49:     MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
 50:     MatRestoreRow(A,i,&nz,&cwork,&vwork);
 51:   }
 52:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 53:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

 55:   if (reuse == MAT_REUSE_MATRIX) {
 56:     MatHeaderReplace(A,M);
 57:   } else {
 58:     *newmat = M;
 59:   }
 60:   return(0);
 61: }
 62: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
 65: PETSC_EXTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
 66: {
 67:   PetscErrorCode    ierr;
 68:   Mat               M;
 69:   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ*)A->data;
 70:   Mat_SeqBAIJ       *Aa     = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
 71:   PetscInt          *d_nnz,*o_nnz;
 72:   PetscInt          i,j,nz;
 73:   PetscInt          m,n,lm,ln;
 74:   PetscInt          rstart,rend;
 75:   const PetscScalar *vwork;
 76:   const PetscInt    *cwork;
 77:   PetscInt          bs = A->rmap->bs;

 80:   MatGetSize(A,&m,&n);
 81:   MatGetLocalSize(A,&lm,&ln);
 82:   PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);

 84:   MatMarkDiagonal_SeqBAIJ(mpimat->A);
 85:   for (i=0; i<lm/bs; i++) {
 86:     d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
 87:     o_nnz[i] = Ba->i[i+1] - Ba->i[i];
 88:   }

 90:   MatCreate(PetscObjectComm((PetscObject)A),&M);
 91:   MatSetSizes(M,lm,ln,m,n);
 92:   MatSetType(M,MATMPISBAIJ);
 93:   MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
 94:   MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);

 96:   PetscFree2(d_nnz,o_nnz);

 98:   MatGetOwnershipRange(A,&rstart,&rend);
 99:   MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
100:   for (i=rstart; i<rend; i++) {
101:     MatGetRow(A,i,&nz,&cwork,&vwork);
102:     j    = 0;
103:     MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
104:     MatRestoreRow(A,i,&nz,&cwork,&vwork);
105:   }
106:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
107:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

109:   if (reuse == MAT_REUSE_MATRIX) {
110:     MatHeaderReplace(A,M);
111:   } else {
112:     *newmat = M;
113:   }
114:   return(0);
115: }