#include <../src/mat/impls/aij/seq/aij.h> #include <../src/mat/impls/baij/seq/baij.h> #include <../src/mat/impls/sbaij/seq/sbaij.h> PETSC_INTERN PetscErrorCode 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; #if defined(PETSC_USE_COMPLEX) const int aconj = A->hermitian ? 1 : 0; #else const int aconj = 0; #endif PetscFunctionBegin; /* compute rowlengths of newmat */ ierr = PetscMalloc2(m,&rowlengths,m+1,&rowstart);CHKERRQ(ierr); for (i=0; irmap->bs);CHKERRQ(ierr); } else B = *newmat; 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_COMM_SELF,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,bs=PetscAbs(A->rmap->bs); MatScalar *av,*bv; PetscBool miss = PETSC_FALSE; PetscFunctionBegin; #if !defined(PETSC_USE_COMPLEX) if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); #else if (!A->symmetric && !A->hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be either symmetric or hermitian. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE) and/or MatSetOption(mat,MAT_HERMITIAN,PETSC_TRUE)"); #endif if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); ierr = PetscMalloc1(m/bs,&rowlengths);CHKERRQ(ierr); for (i=0; idiag[i*bs] == ai[i*bs+1]) { /* missing diagonal */ rowlengths[i] = (ai[i*bs+1] - ai[i*bs])/bs; /* allocate some extra space */ miss = PETSC_TRUE; } else { rowlengths[i] = (ai[i*bs+1] - a->diag[i*bs])/bs; } } if (reuse != MAT_REUSE_MATRIX) { ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation(B,bs,0,rowlengths);CHKERRQ(ierr); } else B = *newmat; if (bs == 1 && !miss) { 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 = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } else { ierr = MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr); /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */ /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before */ /* MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below */ ierr = MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr); } ierr = PetscFree(rowlengths);CHKERRQ(ierr); if (reuse == MAT_INPLACE_MATRIX) { ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); } else *newmat = B; PetscFunctionReturn(0); } PETSC_INTERN PetscErrorCode 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,col,row; MatScalar *av,*bv; #if defined(PETSC_USE_COMPLEX) const int aconj = A->hermitian ? 1 : 0; #else const int aconj = 0; #endif PetscFunctionBegin; /* compute browlengths of newmat */ ierr = PetscMalloc2(mbs,&browlengths,mbs,&browstart);CHKERRQ(ierr); for (i=0; idata); 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_COMM_SELF,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; PetscBool flg; PetscFunctionBegin; if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); if (n != m) SETERRQ(PETSC_COMM_SELF,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_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd); ierr = PetscMalloc1(mbs,&browlengths);CHKERRQ(ierr); for (i=0; idiag[i]; } if (reuse != MAT_REUSE_MATRIX) { ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); } else B = *newmat; 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_INPLACE_MATRIX) { ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); } else *newmat = B; PetscFunctionReturn(0); }