#ifdef PETSC_RCS_HEADER static char vcid[] = "$Id: baij2.c,v 1.15 1997/07/09 20:55:07 balay Exp bsmith $"; #endif #include "src/mat/impls/baij/seq/baij.h" #include "petsc.h" #include "src/inline/bitarray.h" #undef __FUNC__ #define __FUNC__ "MatIncreaseOverlap_SeqBAIJ" /* ADIC Ignore */ int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; int row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val, ival; int start, end, *ai, *aj,bs,*nidx2; char *table; m = a->mbs; ai = a->i; aj = a->j; bs = a->bs; if (ov < 0) SETERRQ(1,0,"Negative overlap specified"); table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); nidx2 = (int *)PetscMalloc((a->m+1)*sizeof(int)); CHKPTRQ(nidx2); for ( i=0; im) SETERRQ(1,0,"index greater than mat-dim"); if(!BT_LOOKUP(table, ival)) { nidx[isz++] = ival;} } ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); ierr = ISDestroy(is[i]); CHKERRQ(ierr); k = 0; for ( j=0; jdata,*c; int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->nbs,*lens; int row,mat_i,*mat_j,tcol,*mat_ilen; int *irow, *icol, nrows, ncols,*ssmap,bs=a->bs, bs2=a->bs2; int *aj = a->j, *ai = a->i; Scalar *mat_a; Mat C; ierr = ISSorted(iscol,(PetscTruth*)&i); if (!i) SETERRQ(1,0,"IS is not sorted"); ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); ssmap = smap; lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); PetscMemzero(smap,oldcols*sizeof(int)); for ( i=0; iilen[irow[i]]; lens[i] = 0; for ( k=kstart; kdata); if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(1,0,""); if (PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int))) { SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros"); } PetscMemzero(c->ilen,c->mbs*sizeof(int)); C = *B; } else { ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr); } c = (Mat_SeqBAIJ *)(C->data); for (i=0; iilen[row]; mat_i = c->i[i]; mat_j = c->j + mat_i; mat_a = c->a + mat_i*bs2; mat_ilen = c->ilen + i; for ( k=kstart; kj[k]])) { *mat_j++ = tcol - 1; PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(Scalar)); mat_a+=bs2; (*mat_ilen)++; } } } /* Free work space */ ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); PetscFree(smap); PetscFree(lens); ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); *B = C; return 0; } #undef __FUNC__ #define __FUNC__ "MatGetSubMatrix_SeqBAIJ" int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; IS is1,is2; int *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count; ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); /* Verify if the indices corespond to each elementin a block and form the IS with compressed IS */ vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int)); CHKPTRQ(vary); iary = vary + a->mbs; PetscMemzero(vary,(a->mbs)*sizeof(int)); for ( i=0; imbs; i++) { if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:"); if (vary[i]==bs) iary[count++] = i; } ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1); CHKERRQ(ierr); PetscMemzero(vary,(a->mbs)*sizeof(int)); for ( i=0; imbs; i++) { if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:"); if (vary[i]==bs) iary[count++] = i; } ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2); CHKERRQ(ierr); ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); PetscFree(vary); ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B); CHKERRQ(ierr); ISDestroy(is1); ISDestroy(is2); return 0; } extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,MatGetSubMatrixCall,Mat*); #undef __FUNC__ #define __FUNC__ "MatGetSubMatrices_SeqBAIJ" int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, Mat **B) { int ierr,i; if (scall == MAT_INITIAL_MATRIX) { *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); } for ( i=0; i