xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
1 
2 #include "../src/mat/impls/aij/seq/aij.h"
3 #include "../src/mat/impls/baij/seq/baij.h"
4 #include "../src/mat/impls/sbaij/seq/sbaij.h"
5 
6 EXTERN_C_BEGIN
7 #undef __FUNCT__
8 #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqAIJ"
9 PetscErrorCode  MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
10 {
11   Mat            B;
12   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
13   Mat_SeqAIJ     *b;
14   PetscErrorCode ierr;
15   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp;
16   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs,diagcnt=0;
17   MatScalar      *av,*bv;
18 
19   PetscFunctionBegin;
20   /* compute rowlengths of newmat */
21   ierr = PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);CHKERRQ(ierr);
22 
23   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
24   aj = a->j;
25   k = 0;
26   for (i=0; i<mbs; i++) {
27     nz = ai[i+1] - ai[i];
28     if (nz) {
29       rowlengths[k] += nz;   /* no. of upper triangular blocks */
30       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
31       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
32 	rowlengths[(*aj)*bs]++; aj++;
33       }
34     }
35     rowlengths[k] *= bs;
36     for (j=1; j<bs; j++) {
37       rowlengths[k+j] = rowlengths[k];
38     }
39     k += bs;
40     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
41   }
42 
43   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
44   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
45   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
46   ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr);
47   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
48   B->rmap->bs = A->rmap->bs;
49 
50   b  = (Mat_SeqAIJ*)(B->data);
51   bi = b->i;
52   bj = b->j;
53   bv = b->a;
54 
55   /* set b->i */
56   bi[0] = 0; rowstart[0] = 0;
57   for (i=0; i<mbs; i++){
58     for (j=0; j<bs; j++){
59       b->ilen[i*bs+j]    = rowlengths[i*bs];
60       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
61     }
62     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
63   }
64   if (bi[mbs] != 2*a->nz - diagcnt) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-diagcnt: %D\n",bi[mbs],2*a->nz - diagcnt);
65 
66   /* set b->j and b->a */
67   aj = a->j; av = a->a;
68   for (i=0; i<mbs; i++) {
69     nz = ai[i+1] - ai[i];
70     /* diagonal block */
71     if (nz && *aj == i) {
72       nz--;
73       for (j=0; j<bs; j++){   /* row i*bs+j */
74 	itmp = i*bs+j;
75 	for (k=0; k<bs; k++){ /* col i*bs+k */
76 	  *(bj + rowstart[itmp]) = (*aj)*bs+k;
77 	  *(bv + rowstart[itmp]) = *(av+k*bs+j);
78 	  rowstart[itmp]++;
79 	}
80       }
81       aj++; av += bs2;
82     }
83 
84     while (nz--){
85       /* lower triangular blocks */
86       for (j=0; j<bs; j++){   /* row (*aj)*bs+j */
87         itmp = (*aj)*bs+j;
88         for (k=0; k<bs; k++){ /* col i*bs+k */
89           *(bj + rowstart[itmp]) = i*bs+k;
90           *(bv + rowstart[itmp]) = *(av+j*bs+k);
91           rowstart[itmp]++;
92         }
93       }
94       /* upper triangular blocks */
95       for (j=0; j<bs; j++){   /* row i*bs+j */
96         itmp = i*bs+j;
97         for (k=0; k<bs; k++){ /* col (*aj)*bs+k */
98           *(bj + rowstart[itmp]) = (*aj)*bs+k;
99           *(bv + rowstart[itmp]) = *(av+k*bs+j);
100           rowstart[itmp]++;
101         }
102       }
103       aj++; av += bs2;
104     }
105   }
106   ierr = PetscFree2(rowlengths,rowstart);CHKERRQ(ierr);
107   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109 
110   if (reuse == MAT_REUSE_MATRIX) {
111     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
112   } else {
113     *newmat = B;
114   }
115   PetscFunctionReturn(0);
116 }
117 EXTERN_C_END
118 
119 EXTERN_C_BEGIN
120 #undef __FUNCT__
121 #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ"
122 PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) {
123   Mat            B;
124   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
125   Mat_SeqSBAIJ   *b;
126   PetscErrorCode ierr;
127   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
128   MatScalar      *av,*bv;
129 
130   PetscFunctionBegin;
131   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
132 
133   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
134   for (i=0; i<m; i++) {
135     rowlengths[i] = ai[i+1] - a->diag[i];
136   }
137   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
138   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
139   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
140   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
141 
142   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
143 
144   b  = (Mat_SeqSBAIJ*)(B->data);
145   bi = b->i;
146   bj = b->j;
147   bv = b->a;
148 
149   bi[0] = 0;
150   for (i=0; i<m; i++) {
151     aj = a->j + a->diag[i];
152     av = a->a + a->diag[i];
153     for (j=0; j<rowlengths[i]; j++){
154       *bj = *aj; bj++; aj++;
155       *bv = *av; bv++; av++;
156     }
157     bi[i+1]    = bi[i] + rowlengths[i];
158     b->ilen[i] = rowlengths[i];
159   }
160 
161   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
162   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164   if (A->hermitian){
165     ierr = MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);
166   }
167 
168   if (reuse == MAT_REUSE_MATRIX) {
169     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
170   } else {
171     *newmat = B;
172   }
173   PetscFunctionReturn(0);
174 }
175 EXTERN_C_END
176 
177 EXTERN_C_BEGIN
178 #undef __FUNCT__
179 #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ"
180 PetscErrorCode  MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
181 {
182   Mat            B;
183   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
184   Mat_SeqBAIJ    *b;
185   PetscErrorCode ierr;
186   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
187   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
188   MatScalar      *av,*bv;
189 
190   PetscFunctionBegin;
191   /* compute browlengths of newmat */
192   ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr);
193   for (i=0; i<mbs; i++) browlengths[i] = 0;
194   aj = a->j;
195   for (i=0; i<mbs; i++) {
196     nz = ai[i+1] - ai[i];
197     aj++; /* skip diagonal */
198     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
199       browlengths[*aj]++; aj++;
200     }
201     browlengths[i] += nz;   /* no. of upper triangular blocks */
202   }
203 
204   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
205   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
206   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
207   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
208   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
209 
210   b  = (Mat_SeqBAIJ*)(B->data);
211   bi = b->i;
212   bj = b->j;
213   bv = b->a;
214 
215   /* set b->i */
216   bi[0] = 0;
217   for (i=0; i<mbs; i++){
218     b->ilen[i]   = browlengths[i];
219     bi[i+1]      = bi[i] + browlengths[i];
220     browstart[i] = bi[i];
221   }
222   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs);
223 
224   /* set b->j and b->a */
225   aj = a->j; av = a->a;
226   for (i=0; i<mbs; i++) {
227     /* diagonal block */
228     *(bj + browstart[i]) = *aj; aj++;
229     itmp = bs2*browstart[i];
230     for (k=0; k<bs2; k++){
231       *(bv + itmp + k) = *av; av++;
232     }
233     browstart[i]++;
234 
235     nz = ai[i+1] - ai[i] -1;
236     while (nz--){
237       /* lower triangular blocks - transpose blocks of A */
238       *(bj + browstart[*aj]) = i; /* block col index */
239       itmp = bs2*browstart[*aj];  /* row index */
240       for (col=0; col<bs; col++){
241         k = col;
242         for (row=0; row<bs; row++){
243           bv[itmp + col*bs+row] = av[k]; k+=bs;
244         }
245       }
246       browstart[*aj]++;
247 
248       /* upper triangular blocks */
249       *(bj + browstart[i]) = *aj; aj++;
250       itmp = bs2*browstart[i];
251       for (k=0; k<bs2; k++){
252         bv[itmp + k] = av[k];
253       }
254       av += bs2;
255       browstart[i]++;
256     }
257   }
258   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
259   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261 
262   if (reuse == MAT_REUSE_MATRIX) {
263     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
264   } else {
265     *newmat = B;
266   }
267   PetscFunctionReturn(0);
268 }
269 EXTERN_C_END
270 
271 EXTERN_C_BEGIN
272 #undef __FUNCT__
273 #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
274 PetscErrorCode  MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
275 {
276   Mat            B;
277   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
278   Mat_SeqSBAIJ   *b;
279   PetscErrorCode ierr;
280   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
281   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
282   MatScalar      *av,*bv;
283   PetscBool      flg;
284 
285   PetscFunctionBegin;
286   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
287   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */
288   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
289 
290   ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr);
291   for (i=0; i<mbs; i++) {
292     browlengths[i] = ai[i+1] - a->diag[i];
293   }
294 
295   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
296   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
297   ierr = MatSetType(B,newtype);CHKERRQ(ierr);
298   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);CHKERRQ(ierr);
299   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
300 
301   b  = (Mat_SeqSBAIJ*)(B->data);
302   bi = b->i;
303   bj = b->j;
304   bv = b->a;
305 
306   bi[0] = 0;
307   for (i=0; i<mbs; i++) {
308     aj = a->j + a->diag[i];
309     av = a->a + (a->diag[i])*bs2;
310     for (j=0; j<browlengths[i]; j++){
311       *bj = *aj; bj++; aj++;
312       for (k=0; k<bs2; k++){
313         *bv = *av; bv++; av++;
314       }
315     }
316     bi[i+1]    = bi[i] + browlengths[i];
317     b->ilen[i] = browlengths[i];
318   }
319   ierr = PetscFree(browlengths);CHKERRQ(ierr);
320   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
321   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
322 
323   if (reuse == MAT_REUSE_MATRIX) {
324     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
325   } else {
326     *newmat = B;
327   }
328 
329   PetscFunctionReturn(0);
330 }
331 EXTERN_C_END
332