Actual source code: aijbaij.c
petsc-3.8.0 2017-09-26
2: #include <../src/mat/impls/baij/seq/baij.h>
4: PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
5: {
6: Mat B;
7: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
9: PetscInt bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k;
10: PetscInt *rowlengths,*rows,*cols,maxlen = 0,ncols;
11: MatScalar *aa = a->a;
14: PetscMalloc1(n*bs,&rowlengths);
15: for (i=0; i<n; i++) {
16: maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
17: for (j=0; j<bs; j++) {
18: rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
19: }
20: }
21: MatCreate(PetscObjectComm((PetscObject)A),&B);
22: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
23: MatSetType(B,MATSEQAIJ);
24: MatSeqAIJSetPreallocation(B,0,rowlengths);
25: MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);
26: PetscFree(rowlengths);
28: PetscMalloc1(bs,&rows);
29: PetscMalloc1(bs*maxlen,&cols);
30: for (i=0; i<n; i++) {
31: for (j=0; j<bs; j++) {
32: rows[j] = i*bs+j;
33: }
34: ncols = ai[i+1] - ai[i];
35: for (k=0; k<ncols; k++) {
36: for (j=0; j<bs; j++) {
37: cols[k*bs+j] = bs*(*aj) + j;
38: }
39: aj++;
40: }
41: MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);
42: aa += ncols*bs*bs;
43: }
44: PetscFree(cols);
45: PetscFree(rows);
46: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
49: B->rmap->bs = A->rmap->bs;
51: if (reuse == MAT_INPLACE_MATRIX) {
52: MatHeaderReplace(A,&B);
53: } else {
54: *newmat = B;
55: }
56: return(0);
57: }
59: #include <../src/mat/impls/aij/seq/aij.h>
61: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
62: {
63: Mat B;
64: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
65: Mat_SeqBAIJ *b;
67: PetscInt *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths;
70: if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
71: if (A->rmap->bs > 1) {
72: MatConvert_Basic(A,newtype,reuse,newmat);
73: return(0);
74: }
76: PetscMalloc1(m,&rowlengths);
77: for (i=0; i<m; i++) {
78: rowlengths[i] = ai[i+1] - ai[i];
79: }
80: MatCreate(PetscObjectComm((PetscObject)A),&B);
81: MatSetSizes(B,m,n,m,n);
82: MatSetType(B,MATSEQBAIJ);
83: MatSeqBAIJSetPreallocation(B,1,0,rowlengths);
84: PetscFree(rowlengths);
86: MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
88: b = (Mat_SeqBAIJ*)(B->data);
90: PetscMemcpy(b->i,a->i,(m+1)*sizeof(PetscInt));
91: PetscMemcpy(b->ilen,a->ilen,m*sizeof(PetscInt));
92: PetscMemcpy(b->j,a->j,a->nz*sizeof(PetscInt));
93: PetscMemcpy(b->a,a->a,a->nz*sizeof(MatScalar));
95: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
96: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
98: if (reuse == MAT_INPLACE_MATRIX) {
99: MatHeaderReplace(A,&B);
100: } else {
101: *newmat = B;
102: }
103: return(0);
104: }