Actual source code: aijbaij.c
1: #define PETSCMAT_DLL
3: #include ../src/mat/impls/baij/seq/baij.h
8: PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
9: {
10: Mat B;
11: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
13: PetscInt bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k;
14: PetscInt *rowlengths,*rows,*cols,maxlen = 0,ncols;
15: MatScalar *aa = a->a;
18: PetscMalloc(n*bs*sizeof(PetscInt),&rowlengths);
19: for (i=0; i<n; i++) {
20: maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
21: for (j=0; j<bs; j++) {
22: rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
23: }
24: }
25: MatCreate(((PetscObject)A)->comm,&B);
26: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
27: MatSetType(B,newtype);
28: MatSeqAIJSetPreallocation(B,0,rowlengths);
29: MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);
30: PetscFree(rowlengths);
32: PetscMalloc(bs*sizeof(PetscInt),&rows);
33: PetscMalloc(bs*maxlen*sizeof(PetscInt),&cols);
34: for (i=0; i<n; i++) {
35: for (j=0; j<bs; j++) {
36: rows[j] = i*bs+j;
37: }
38: ncols = ai[i+1] - ai[i];
39: for (k=0; k<ncols; k++) {
40: for (j=0; j<bs; j++) {
41: cols[k*bs+j] = bs*(*aj) + j;
42: }
43: aj++;
44: }
45: MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);
46: aa += ncols*bs*bs;
47: }
48: PetscFree(cols);
49: PetscFree(rows);
50: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
51: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
52:
53: B->rmap->bs = A->rmap->bs;
55: if (reuse == MAT_REUSE_MATRIX) {
56: MatHeaderReplace(A,B);
57: } else {
58: *newmat = B;
59: }
60: return(0);
61: }
64: #include ../src/mat/impls/aij/seq/aij.h
69: PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
70: {
71: Mat B;
72: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
73: Mat_SeqBAIJ *b;
75: PetscInt *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths;
78: if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
80: PetscMalloc(m*sizeof(PetscInt),&rowlengths);
81: for (i=0; i<m; i++) {
82: rowlengths[i] = ai[i+1] - ai[i];
83: }
84: MatCreate(((PetscObject)A)->comm,&B);
85: MatSetSizes(B,m,n,m,n);
86: MatSetType(B,newtype);
87: MatSeqBAIJSetPreallocation_SeqBAIJ(B,1,0,rowlengths);
88: PetscFree(rowlengths);
90: MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
91:
92: b = (Mat_SeqBAIJ*)(B->data);
94: PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));
95: PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));
96: PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));
97: PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));
98:
99: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
102: if (reuse == MAT_REUSE_MATRIX) {
103: MatHeaderReplace(A,B);
104: } else {
105: *newmat = B;
106: }
107: return(0);
108: }