1 /*$Id: aijbaij.c,v 1.4 2001/01/15 21:45:50 bsmith Exp balay $*/ 2 3 #include "src/mat/impls/baij/seq/baij.h" 4 5 EXTERN_C_BEGIN 6 #undef __FUNC__ 7 #define __FUNC__ "MatConvert_SeqBAI_SeqAIJ" 8 int MatConvert_SeqBAIJ_SeqAIJ(Mat A,MatType newtype,Mat *B) 9 { 10 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 11 int ierr,bs = a->bs,*ai = a->i,*aj = a->j,n = A->M/bs,i,j,k; 12 int *rowlengths,*rows,*cols,maxlen = 0,ncols; 13 Scalar *aa = a->a; 14 15 PetscFunctionBegin; 16 ierr = PetscMalloc(n*bs*sizeof(int),&rowlengths);CHKERRQ(ierr); 17 for (i=0; i<n; i++) { 18 maxlen = PetscMax(maxlen,(ai[i+1] - ai[i])); 19 for (j=0; j<bs; j++) { 20 rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]); 21 } 22 } 23 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,A->m,A->n,0,rowlengths,B);CHKERRQ(ierr); 24 ierr = MatSetOption(*B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 25 ierr = MatSetOption(*B,MAT_ROWS_SORTED);CHKERRQ(ierr); 26 ierr = MatSetOption(*B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 27 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 28 29 ierr = PetscMalloc(bs*sizeof(int),&rows);CHKERRQ(ierr); 30 ierr = PetscMalloc(bs*maxlen*sizeof(int),&cols);CHKERRQ(ierr); 31 for (i=0; i<n; i++) { 32 for (j=0; j<bs; j++) { 33 rows[j] = i*bs+j; 34 } 35 ncols = ai[i+1] - ai[i]; 36 for (k=0; k<ncols; k++) { 37 for (j=0; j<bs; j++) { 38 cols[k*bs+j] = bs*(*aj) + j; 39 } 40 aj++; 41 } 42 ierr = MatSetValues(*B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES); 43 aa += ncols*bs*bs; 44 } 45 ierr = PetscFree(cols);CHKERRQ(ierr); 46 ierr = PetscFree(rows);CHKERRQ(ierr); 47 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 EXTERN_C_END 52 53