xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 2475b7ca256cea2a4b7cbf2d8babcda14e5fa36e)
1 
2 #include <../src/mat/impls/baij/seq/baij.h>
3 
4 PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
5 {
6   Mat            B;
7   Mat_SeqAIJ     *b;
8   PetscBool      roworiented;
9   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
10   PetscErrorCode ierr;
11   PetscInt       bs = A->rmap->bs,*ai = a->i,*aj = a->j,n = A->rmap->N/bs,i,j,k;
12   PetscInt       *rowlengths,*rows,*cols,maxlen = 0,ncols;
13   MatScalar      *aa = a->a;
14 
15   PetscFunctionBegin;
16   if (reuse == MAT_REUSE_MATRIX) {
17     B = *newmat;
18     for (i=0; i<n; i++) {
19       maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
20     }
21   } else {
22     ierr = PetscMalloc1(n*bs,&rowlengths);CHKERRQ(ierr);
23     for (i=0; i<n; i++) {
24       maxlen = PetscMax(maxlen,(ai[i+1] - ai[i]));
25       for (j=0; j<bs; j++) {
26         rowlengths[i*bs+j] = bs*(ai[i+1] - ai[i]);
27       }
28     }
29     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
30     ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
31     ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
32     ierr = MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr);
33     ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
34     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
35   }
36   b = (Mat_SeqAIJ*)B->data;
37   roworiented = b->roworiented;
38 
39   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
40   ierr = PetscMalloc1(bs,&rows);CHKERRQ(ierr);
41   ierr = PetscMalloc1(bs*maxlen,&cols);CHKERRQ(ierr);
42   for (i=0; i<n; i++) {
43     for (j=0; j<bs; j++) {
44       rows[j] = i*bs+j;
45     }
46     ncols = ai[i+1] - ai[i];
47     for (k=0; k<ncols; k++) {
48       for (j=0; j<bs; j++) {
49         cols[k*bs+j] = bs*(*aj) + j;
50       }
51       aj++;
52     }
53     ierr = MatSetValues(B,bs,rows,bs*ncols,cols,aa,INSERT_VALUES);CHKERRQ(ierr);
54     aa  += ncols*bs*bs;
55   }
56   ierr = PetscFree(cols);CHKERRQ(ierr);
57   ierr = PetscFree(rows);CHKERRQ(ierr);
58   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
60   ierr = MatSetOption(B,MAT_ROW_ORIENTED,roworiented);CHKERRQ(ierr);
61 
62   if (reuse == MAT_INPLACE_MATRIX) {
63     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
64   } else {
65     *newmat = B;
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 #include <../src/mat/impls/aij/seq/aij.h>
71 
72 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
73 {
74   Mat            B;
75   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
76   Mat_SeqBAIJ    *b;
77   PetscErrorCode ierr;
78   PetscInt       *ai=a->i,m=A->rmap->N,n=A->cmap->N,i,*rowlengths;
79 
80   PetscFunctionBegin;
81   if (A->rmap->bs > 1) {
82     ierr = MatConvert_Basic(A,newtype,reuse,newmat);CHKERRQ(ierr);
83     PetscFunctionReturn(0);
84   }
85   if (reuse != MAT_REUSE_MATRIX) {
86     ierr = PetscMalloc1(m,&rowlengths);CHKERRQ(ierr);
87     for (i=0; i<m; i++) {
88       rowlengths[i] = ai[i+1] - ai[i];
89     }
90     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
91     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
92     ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
93     ierr = MatSeqBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr);
94     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
95   } else B = *newmat;
96 
97   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
98 
99   b = (Mat_SeqBAIJ*)(B->data);
100 
101   ierr = PetscArraycpy(b->i,a->i,m+1);CHKERRQ(ierr);
102   ierr = PetscArraycpy(b->ilen,a->ilen,m);CHKERRQ(ierr);
103   ierr = PetscArraycpy(b->j,a->j,a->nz);CHKERRQ(ierr);
104   ierr = PetscArraycpy(b->a,a->a,a->nz);CHKERRQ(ierr);
105 
106   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108 
109   if (reuse == MAT_INPLACE_MATRIX) {
110     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
111   } else {
112     *newmat = B;
113   }
114   PetscFunctionReturn(0);
115 }
116