Actual source code: mpiaijsbaij.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
  3:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  4:  #include <petsc/private/matimpl.h>
  5:  #include <petscmat.h>

  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;

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

 27:     MatMarkDiagonal_SeqAIJ(mpimat->A);
 28:     for (i=0; i<lm/bs; i++) {
 29:       d_nnz[i] = (Aa->i[i*bs+1] - Aa->diag[i*bs])/bs;
 30:       o_nnz[i] = (Ba->i[i*bs+1] - Ba->i[i*bs])/bs;
 31:     }

 33:     MatCreate(PetscObjectComm((PetscObject)A),&M);
 34:     MatSetSizes(M,lm,ln,m,n);
 35:     MatSetType(M,MATMPISBAIJ);
 36:     MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
 37:     MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);
 38:     PetscFree2(d_nnz,o_nnz);
 39:   } else M = *newmat;

 41:   if (bs == 1) {
 42:     MatGetOwnershipRange(A,&rstart,&rend);
 43:     for (i=rstart; i<rend; i++) {
 44:       MatGetRow(A,i,&nz,&cwork,&vwork);
 45:       if (nz) {
 46:         j = 0;
 47:         while (cwork[j] < i) {
 48:           j++; nz--;
 49:         }
 50:         MatSetValues(M,1,&i,nz,cwork+j,vwork+j,INSERT_VALUES);
 51:       }
 52:       MatRestoreRow(A,i,&nz,&cwork,&vwork);
 53:     }
 54:     MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 55:     MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
 56:   } else {
 57:     MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
 58:     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
 59:     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
 60:     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
 61:     MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&M);
 62:   }

 64:   if (reuse == MAT_INPLACE_MATRIX) {
 65:     MatHeaderReplace(A,&M);
 66:   } else *newmat = M;
 67:   return(0);
 68: }
 69: /* contributed by Dahai Guo <dhguo@ncsa.uiuc.edu> April 2011 */
 70: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
 71: {
 72:   PetscErrorCode    ierr;
 73:   Mat               M;
 74:   Mat_MPIBAIJ       *mpimat = (Mat_MPIBAIJ*)A->data;
 75:   Mat_SeqBAIJ       *Aa     = (Mat_SeqBAIJ*)mpimat->A->data,*Ba = (Mat_SeqBAIJ*)mpimat->B->data;
 76:   PetscInt          *d_nnz,*o_nnz;
 77:   PetscInt          i,nz;
 78:   PetscInt          m,n,lm,ln;
 79:   PetscInt          rstart,rend;
 80:   const PetscScalar *vwork;
 81:   const PetscInt    *cwork;
 82:   PetscInt          bs = A->rmap->bs;

 85:   if (reuse != MAT_REUSE_MATRIX) {
 86:     MatGetSize(A,&m,&n);
 87:     MatGetLocalSize(A,&lm,&ln);
 88:     PetscMalloc2(lm/bs,&d_nnz,lm/bs,&o_nnz);

 90:     MatMarkDiagonal_SeqBAIJ(mpimat->A);
 91:     for (i=0; i<lm/bs; i++) {
 92:       d_nnz[i] = Aa->i[i+1] - Aa->diag[i];
 93:       o_nnz[i] = Ba->i[i+1] - Ba->i[i];
 94:     }

 96:     MatCreate(PetscObjectComm((PetscObject)A),&M);
 97:     MatSetSizes(M,lm,ln,m,n);
 98:     MatSetType(M,MATMPISBAIJ);
 99:     MatSeqSBAIJSetPreallocation(M,bs,0,d_nnz);
100:     MatMPISBAIJSetPreallocation(M,bs,0,d_nnz,0,o_nnz);

102:     PetscFree2(d_nnz,o_nnz);
103:   } else M = *newmat;

105:   MatGetOwnershipRange(A,&rstart,&rend);
106:   MatSetOption(M,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
107:   for (i=rstart; i<rend; i++) {
108:     MatGetRow(A,i,&nz,&cwork,&vwork);
109:     MatSetValues(M,1,&i,nz,cwork,vwork,INSERT_VALUES);
110:     MatRestoreRow(A,i,&nz,&cwork,&vwork);
111:   }
112:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
113:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);

115:   if (reuse == MAT_INPLACE_MATRIX) {
116:     MatHeaderReplace(A,&M);
117:   } else *newmat = M;
118:   return(0);
119: }