#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; 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 == PETSC_BOOL3_TRUE ? 1 : 0; #else const int aconj = 0; #endif PetscFunctionBegin; /* compute rowlengths of newmat */ PetscCall(PetscMalloc2(m,&rowlengths,m+1,&rowstart)); for (i=0; irmap->bs)); } 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; } PetscCheck(bi[mbs] == 2*a->nz - diagcnt,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %" PetscInt_FMT " != 2*a->nz-diagcnt: %" PetscInt_FMT,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; 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) PetscCheck(A->symmetric == PETSC_BOOL3_TRUE,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); #else PetscCheck(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE,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 PetscCheck(n == m,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); PetscCall(PetscMalloc1(m/bs,&rowlengths)); 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) { PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); PetscCall(MatSetSizes(B,m,n,m,n)); PetscCall(MatSetType(B,MATSEQSBAIJ)); PetscCall(MatSeqSBAIJSetPreallocation(B,bs,0,rowlengths)); } 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]; } PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); } else { PetscCall(MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); /* 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 */ PetscCall(MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&B)); } PetscCall(PetscFree(rowlengths)); if (reuse == MAT_INPLACE_MATRIX) { PetscCall(MatHeaderReplace(A,&B)); } 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; 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 == PETSC_BOOL3_TRUE ? 1 : 0; #else const int aconj = 0; #endif PetscFunctionBegin; /* compute browlengths of newmat */ PetscCall(PetscMalloc2(mbs,&browlengths,mbs,&browstart)); 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]; } PetscCheck(bi[mbs] == 2*a->nz - mbs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %" PetscInt_FMT " != 2*a->nz - mbs: %" PetscInt_FMT,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; 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; PetscCheck(A->symmetric,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)"); PetscCheck(n == m,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); PetscCall(MatMissingDiagonal_SeqBAIJ(A,&flg,&dd)); /* check for missing diagonals, then mark diag */ PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %" PetscInt_FMT,dd); PetscCall(PetscMalloc1(mbs,&browlengths)); for (i=0; idiag[i]; } if (reuse != MAT_REUSE_MATRIX) { PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); PetscCall(MatSetSizes(B,m,n,m,n)); PetscCall(MatSetType(B,MATSEQSBAIJ)); PetscCall(MatSeqSBAIJSetPreallocation(B,bs,0,browlengths)); } 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]; } PetscCall(PetscFree(browlengths)); PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); if (reuse == MAT_INPLACE_MATRIX) { PetscCall(MatHeaderReplace(A,&B)); } else *newmat = B; PetscFunctionReturn(0); }