1 /*$Id: aijbaij.c,v 1.9 2001/08/07 03:02:55 balay Exp $*/ 2 3 #include "src/mat/impls/baij/seq/baij.h" 4 5 EXTERN_C_BEGIN 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatConvert_SeqBAI_SeqAIJ" 8 int MatConvert_SeqBAIJ_SeqAIJ(Mat A,const MatType newtype,Mat *newmat) { 9 Mat B; 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 PetscScalar *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 = MatCreate(A->comm,A->m,A->n,A->m,A->n,&B);CHKERRQ(ierr); 24 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 25 ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 26 ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 27 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 28 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 29 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 30 31 ierr = PetscMalloc(bs*sizeof(int),&rows);CHKERRQ(ierr); 32 ierr = PetscMalloc(bs*maxlen*sizeof(int),&cols);CHKERRQ(ierr); 33 for (i=0; i<n; i++) { 34 for (j=0; j<bs; j++) { 35 rows[j] = i*bs+j; 36 } 37 ncols = ai[i+1] - ai[i]; 38 for (k=0; k<ncols; k++) { 39 for (j=0; j<bs; j++) { 40 cols[k*bs+j] = bs*(*aj) + j; 41 } 42 aj++; 43 } 44 ierr = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 45 aa += ncols*bs*bs; 46 } 47 ierr = PetscFree(cols);CHKERRQ(ierr); 48 ierr = PetscFree(rows);CHKERRQ(ierr); 49 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51 52 /* Fake support for "inplace" convert. */ 53 if (*newmat == A) { 54 ierr = MatDestroy(A);CHKERRQ(ierr); 55 } 56 *newmat = B; 57 PetscFunctionReturn(0); 58 } 59 EXTERN_C_END 60 61 #include "src/mat/impls/aij/seq/aij.h" 62 63 EXTERN_C_BEGIN 64 #undef __FUNCT__ 65 #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ" 66 int MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,Mat *newmat) { 67 Mat B; 68 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 69 Mat_SeqBAIJ *b; 70 int ierr,*ai=a->i,m=A->M,n=A->N,i,*rowlengths; 71 72 PetscFunctionBegin; 73 if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 74 75 ierr = PetscMalloc(m*sizeof(int),&rowlengths);CHKERRQ(ierr); 76 for (i=0; i<m; i++) { 77 rowlengths[i] = ai[i+1] - ai[i]; 78 } 79 ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 80 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 81 ierr = MatSeqBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr); 82 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 83 84 ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 85 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 86 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 87 88 b = (Mat_SeqBAIJ*)(B->data); 89 90 ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(int));CHKERRQ(ierr); 91 ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(int));CHKERRQ(ierr); 92 ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(int));CHKERRQ(ierr); 93 ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr); 94 95 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97 98 /* Fake support for "inplace" convert. */ 99 if (*newmat == A) { 100 ierr = MatDestroy(A);CHKERRQ(ierr); 101 } 102 *newmat = B; 103 PetscFunctionReturn(0); 104 } 105 EXTERN_C_END 106