1 2 #include <../src/mat/impls/baij/seq/baij.h> 3 4 EXTERN_C_BEGIN 5 #undef __FUNCT__ 6 #define __FUNCT__ "MatConvert_SeqBAIJ_SeqAIJ" 7 PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 8 { 9 Mat B; 10 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 11 PetscErrorCode ierr; 12 PetscInt bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k; 13 PetscInt *rowlengths,*rows,*cols,maxlen = 0,ncols; 14 MatScalar *aa = a->a; 15 16 PetscFunctionBegin; 17 ierr = PetscMalloc(n*bs*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 18 for (i=0; i<n; i++) { 19 maxlen = PetscMax(maxlen,(ai[i+1] - ai[i])); 20 for (j=0; j<bs; j++) { 21 rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]); 22 } 23 } 24 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 25 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 26 ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 27 ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 28 ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 29 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 30 31 ierr = PetscMalloc(bs*sizeof(PetscInt),&rows);CHKERRQ(ierr); 32 ierr = PetscMalloc(bs*maxlen*sizeof(PetscInt),&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 B->rmap->bs = A->rmap->bs; 53 54 if (reuse == MAT_REUSE_MATRIX) { 55 ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 56 } else { 57 *newmat = B; 58 } 59 PetscFunctionReturn(0); 60 } 61 EXTERN_C_END 62 63 #include <../src/mat/impls/aij/seq/aij.h> 64 65 EXTERN_C_BEGIN 66 #undef __FUNCT__ 67 #define __FUNCT__ "MatConvert_SeqAIJ_SeqBAIJ" 68 PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 69 { 70 Mat B; 71 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 72 Mat_SeqBAIJ *b; 73 PetscErrorCode ierr; 74 PetscInt *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths; 75 76 PetscFunctionBegin; 77 if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square"); 78 79 ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 80 for (i=0; i<m; i++) { 81 rowlengths[i] = ai[i+1] - ai[i]; 82 } 83 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 84 ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 85 ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr); 86 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,1,0,rowlengths);CHKERRQ(ierr); 87 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 88 89 ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 90 91 b = (Mat_SeqBAIJ*)(B->data); 92 93 ierr = PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 94 ierr = PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));CHKERRQ(ierr); 95 ierr = PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));CHKERRQ(ierr); 96 ierr = PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));CHKERRQ(ierr); 97 98 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100 101 if (reuse == MAT_REUSE_MATRIX) { 102 ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 103 } else { 104 *newmat = B; 105 } 106 PetscFunctionReturn(0); 107 } 108 EXTERN_C_END 109