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