#define PETSCMAT_DLL #include "../src/mat/impls/aij/seq/aij.h" #include "../src/mat/impls/baij/seq/baij.h" #include "../src/mat/impls/sbaij/seq/sbaij.h" EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) { Mat B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; Mat_SeqAIJ *b; PetscErrorCode ierr; PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp; PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0; MatScalar *av,*bv; PetscFunctionBegin; /* compute rowlengths of newmat */ ierr = PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);CHKERRQ(ierr); for (i=0; ij; k = 0; for (i=0; icomm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,newtype);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); B->rmap->bs = A->rmap->bs; b = (Mat_SeqAIJ*)(B->data); bi = b->i; bj = b->j; bv = b->a; /* set b->i */ bi[0] = 0; rowstart[0] = 0; for (i=0; iilen[i*bs+j] = rowlengths[i*bs]; rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs]; } bi[i+1] = bi[i] + rowlengths[i*bs]/bs; } if (bi[mbs] != 2*a->nz - diagcnt) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-diagcnt: %D\n",bi[mbs],2*a->nz - diagcnt); /* set b->j and b->a */ aj = a->j; av = a->a; for (i=0; idata; Mat_SeqSBAIJ *b; PetscErrorCode ierr; PetscInt *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths; MatScalar *av,*bv; PetscFunctionBegin; if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); for (i=0; idiag[i]; } ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,newtype);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); b = (Mat_SeqSBAIJ*)(B->data); bi = b->i; bj = b->j; bv = b->a; bi[0] = 0; for (i=0; ij + a->diag[i]; av = a->a + a->diag[i]; for (j=0; jilen[i] = rowlengths[i]; } ierr = PetscFree(rowlengths);CHKERRQ(ierr); ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); if (A->hermitian){ ierr = MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); } if (reuse == MAT_REUSE_MATRIX) { ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); } else { *newmat = B; } PetscFunctionReturn(0); } EXTERN_C_END EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) { Mat B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; Mat_SeqBAIJ *b; PetscErrorCode ierr; PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp; PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs; MatScalar *av,*bv; PetscFunctionBegin; /* compute browlengths of newmat */ ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr); for (i=0; ij; for (i=0; icomm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,newtype);CHKERRQ(ierr); ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); b = (Mat_SeqBAIJ*)(B->data); bi = b->i; bj = b->j; bv = b->a; /* set b->i */ bi[0] = 0; for (i=0; iilen[i] = browlengths[i]; bi[i+1] = bi[i] + browlengths[i]; browstart[i] = bi[i]; } if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs); /* set b->j and b->a */ aj = a->j; av = a->a; for (i=0; idata; Mat_SeqSBAIJ *b; PetscErrorCode ierr; PetscInt *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths; PetscInt bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd; MatScalar *av,*bv; PetscTruth flg; PetscFunctionBegin; if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ if (flg) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd); ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); for (i=0; idiag[i]; } ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,newtype);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); b = (Mat_SeqSBAIJ*)(B->data); bi = b->i; bj = b->j; bv = b->a; bi[0] = 0; for (i=0; ij + a->diag[i]; av = a->a + (a->diag[i])*bs2; for (j=0; jilen[i] = browlengths[i]; } ierr = PetscFree(browlengths);CHKERRQ(ierr); ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); if (reuse == MAT_REUSE_MATRIX) { ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); } else { *newmat = B; } PetscFunctionReturn(0); } EXTERN_C_END