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