xref: /petsc/src/mat/impls/sbaij/seq/aijsbaij.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
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,MATSEQAIJ);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,MatType newtype,MatReuse reuse,Mat *newmat)
123 {
124   Mat            B;
125   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
126   Mat_SeqSBAIJ   *b;
127   PetscErrorCode ierr;
128   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
129   MatScalar      *av,*bv;
130 
131   PetscFunctionBegin;
132   if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
133   if (n != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix must be square");
134 
135   ierr = PetscMalloc(m*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
136   for (i=0; i<m; i++) {
137     rowlengths[i] = ai[i+1] - a->diag[i];
138   }
139   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
140   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
141   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
142   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);CHKERRQ(ierr);
143 
144   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
145 
146   b  = (Mat_SeqSBAIJ*)(B->data);
147   bi = b->i;
148   bj = b->j;
149   bv = b->a;
150 
151   bi[0] = 0;
152   for (i=0; i<m; i++) {
153     aj = a->j + a->diag[i];
154     av = a->a + a->diag[i];
155     for (j=0; j<rowlengths[i]; j++) {
156       *bj = *aj; bj++; aj++;
157       *bv = *av; bv++; av++;
158     }
159     bi[i+1]    = bi[i] + rowlengths[i];
160     b->ilen[i] = rowlengths[i];
161   }
162 
163   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
164   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166 
167   if (reuse == MAT_REUSE_MATRIX) {
168     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
169   } else {
170     *newmat = B;
171   }
172   PetscFunctionReturn(0);
173 }
174 EXTERN_C_END
175 
176 EXTERN_C_BEGIN
177 #undef __FUNCT__
178 #define __FUNCT__ "MatConvert_SeqSBAIJ_SeqBAIJ"
179 PetscErrorCode  MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
180 {
181   Mat            B;
182   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
183   Mat_SeqBAIJ    *b;
184   PetscErrorCode ierr;
185   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
186   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,col,row;
187   MatScalar      *av,*bv;
188 
189   PetscFunctionBegin;
190   /* compute browlengths of newmat */
191   ierr = PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);CHKERRQ(ierr);
192   for (i=0; i<mbs; i++) browlengths[i] = 0;
193   aj = a->j;
194   for (i=0; i<mbs; i++) {
195     nz = ai[i+1] - ai[i];
196     aj++; /* skip diagonal */
197     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
198       browlengths[*aj]++; aj++;
199     }
200     browlengths[i] += nz;   /* no. of upper triangular blocks */
201   }
202 
203   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
204   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
205   ierr = MatSetType(B,MATSEQBAIJ);CHKERRQ(ierr);
206   ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr);
207   ierr = MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
208 
209   b  = (Mat_SeqBAIJ*)(B->data);
210   bi = b->i;
211   bj = b->j;
212   bv = b->a;
213 
214   /* set b->i */
215   bi[0] = 0;
216   for (i=0; i<mbs; i++) {
217     b->ilen[i]   = browlengths[i];
218     bi[i+1]      = bi[i] + browlengths[i];
219     browstart[i] = bi[i];
220   }
221   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);
222 
223   /* set b->j and b->a */
224   aj = a->j; av = a->a;
225   for (i=0; i<mbs; i++) {
226     /* diagonal block */
227     *(bj + browstart[i]) = *aj; aj++;
228     itmp = bs2*browstart[i];
229     for (k=0; k<bs2; k++) {
230       *(bv + itmp + k) = *av; av++;
231     }
232     browstart[i]++;
233 
234     nz = ai[i+1] - ai[i] -1;
235     while (nz--) {
236       /* lower triangular blocks - transpose blocks of A */
237       *(bj + browstart[*aj]) = i; /* block col index */
238       itmp = bs2*browstart[*aj];  /* row index */
239       for (col=0; col<bs; col++) {
240         k = col;
241         for (row=0; row<bs; row++) {
242           bv[itmp + col*bs+row] = av[k]; k+=bs;
243         }
244       }
245       browstart[*aj]++;
246 
247       /* upper triangular blocks */
248       *(bj + browstart[i]) = *aj; aj++;
249       itmp = bs2*browstart[i];
250       for (k=0; k<bs2; k++) {
251         bv[itmp + k] = av[k];
252       }
253       av += bs2;
254       browstart[i]++;
255     }
256   }
257   ierr = PetscFree2(browlengths,browstart);CHKERRQ(ierr);
258   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260 
261   if (reuse == MAT_REUSE_MATRIX) {
262     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
263   } else {
264     *newmat = B;
265   }
266   PetscFunctionReturn(0);
267 }
268 EXTERN_C_END
269 
270 EXTERN_C_BEGIN
271 #undef __FUNCT__
272 #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ"
273 PetscErrorCode  MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
274 {
275   Mat            B;
276   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
277   Mat_SeqSBAIJ   *b;
278   PetscErrorCode ierr;
279   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
280   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
281   MatScalar      *av,*bv;
282   PetscBool      flg;
283 
284   PetscFunctionBegin;
285   if (!A->symmetric) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
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,MATSEQSBAIJ);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   PetscFunctionReturn(0);
329 }
330 EXTERN_C_END
331