Actual source code: aijsbaij.c

  1: #define PETSCMAT_DLL

 3:  #include ../src/mat/impls/aij/seq/aij.h
 4:  #include ../src/mat/impls/baij/seq/baij.h
 5:  #include ../src/mat/impls/sbaij/seq/sbaij.h

 10: PetscErrorCode  MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
 11: {
 12:   Mat            B;
 13:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 14:   Mat_SeqAIJ     *b;
 16:   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*rowlengths,nz,*rowstart,itmp;
 17:   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=A->rmap->N/bs;
 18:   MatScalar      *av,*bv;

 21:   /* compute rowlengths of newmat */
 22:   PetscMalloc2(m,PetscInt,&rowlengths,m+1,PetscInt,&rowstart);
 23: 
 24:   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
 25:   aj = a->j;
 26:   k = 0;
 27:   for (i=0; i<mbs; i++) {
 28:     nz = ai[i+1] - ai[i];
 29:     aj++; /* skip diagonal */
 30:     for (j=1; j<nz; j++) { /* no. of lower triangular blocks */
 31:       rowlengths[(*aj)*bs]++; aj++;
 32:     }
 33:     rowlengths[k] += nz;   /* no. of upper triangular blocks */
 34:     rowlengths[k] *= bs;
 35:     for (j=1; j<bs; j++) {
 36:       rowlengths[k+j] = rowlengths[k];
 37:     }
 38:     k += bs;
 39:     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
 40:   }
 41: 
 42:   MatCreate(((PetscObject)A)->comm,&B);
 43:   MatSetSizes(B,m,n,m,n);
 44:   MatSetType(B,newtype);
 45:   MatSeqAIJSetPreallocation(B,0,rowlengths);
 46:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);
 47:   B->rmap->bs = A->rmap->bs;

 49:   b  = (Mat_SeqAIJ*)(B->data);
 50:   bi = b->i;
 51:   bj = b->j;
 52:   bv = b->a;

 54:   /* set b->i */
 55:   bi[0] = 0; rowstart[0] = 0;
 56:   for (i=0; i<mbs; i++){
 57:     for (j=0; j<bs; j++){
 58:       b->ilen[i*bs+j]  = rowlengths[i*bs];
 59:       rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs];
 60:     }
 61:     bi[i+1]     = bi[i] + rowlengths[i*bs]/bs;
 62:   }
 63:   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-mbs: %D\n",bi[mbs],2*a->nz - mbs);

 65:   /* set b->j and b->a */
 66:   aj = a->j; av = a->a;
 67:   for (i=0; i<mbs; i++) {
 68:     /* diagonal block */
 69:     for (j=0; j<bs; j++){   /* row i*bs+j */
 70:       itmp = i*bs+j;
 71:       for (k=0; k<bs; k++){ /* col i*bs+k */
 72:         *(bj + rowstart[itmp]) = (*aj)*bs+k;
 73:         *(bv + rowstart[itmp]) = *(av+k*bs+j);
 74:         rowstart[itmp]++;
 75:       }
 76:     }
 77:     aj++; av += bs2;
 78: 
 79:     nz = ai[i+1] - ai[i] -1;
 80:     while (nz--){
 81:       /* lower triangular blocks */
 82:       for (j=0; j<bs; j++){   /* row (*aj)*bs+j */
 83:         itmp = (*aj)*bs+j;
 84:         for (k=0; k<bs; k++){ /* col i*bs+k */
 85:           *(bj + rowstart[itmp]) = i*bs+k;
 86:           *(bv + rowstart[itmp]) = *(av+k*bs+j);
 87:           rowstart[itmp]++;
 88:         }
 89:       }
 90:       /* upper triangular blocks */
 91:       for (j=0; j<bs; j++){   /* row i*bs+j */
 92:         itmp = i*bs+j;
 93:         for (k=0; k<bs; k++){ /* col (*aj)*bs+k */
 94:           *(bj + rowstart[itmp]) = (*aj)*bs+k;
 95:           *(bv + rowstart[itmp]) = *(av+k*bs+j);
 96:           rowstart[itmp]++;
 97:         }
 98:       }
 99:       aj++; av += bs2;
100:     }
101:   }
102:   PetscFree2(rowlengths,rowstart);
103:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
104:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

106:   if (reuse == MAT_REUSE_MATRIX) {
107:     MatHeaderReplace(A,B);
108:   } else {
109:     *newmat = B;
110:   }
111:   return(0);
112: }

118: PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) {
119:   Mat            B;
120:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
121:   Mat_SeqSBAIJ   *b;
123:   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->N,i,j,*bi,*bj,*rowlengths;
124:   MatScalar      *av,*bv;

127:   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");

129:   PetscMalloc(m*sizeof(PetscInt),&rowlengths);
130:   for (i=0; i<m; i++) {
131:     rowlengths[i] = ai[i+1] - a->diag[i];
132:   }
133:   MatCreate(((PetscObject)A)->comm,&B);
134:   MatSetSizes(B,m,n,m,n);
135:   MatSetType(B,newtype);
136:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,1,0,rowlengths);

138:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
139: 
140:   b  = (Mat_SeqSBAIJ*)(B->data);
141:   bi = b->i;
142:   bj = b->j;
143:   bv = b->a;

145:   bi[0] = 0;
146:   for (i=0; i<m; i++) {
147:     aj = a->j + a->diag[i];
148:     av = a->a + a->diag[i];
149:     for (j=0; j<rowlengths[i]; j++){
150:       *bj = *aj; bj++; aj++;
151:       *bv = *av; bv++; av++;
152:     }
153:     bi[i+1]    = bi[i] + rowlengths[i];
154:     b->ilen[i] = rowlengths[i];
155:   }
156: 
157:   PetscFree(rowlengths);
158:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
159:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
160:   if (A->hermitian){
161:     MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
162:   }

164:   if (reuse == MAT_REUSE_MATRIX) {
165:     MatHeaderReplace(A,B);
166:   } else {
167:     *newmat = B;
168:   }
169:   return(0);
170: }

176: PetscErrorCode  MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
177: {
178:   Mat            B;
179:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
180:   Mat_SeqBAIJ    *b;
182:   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->n,i,k,*bi,*bj,*browlengths,nz,*browstart,itmp;
183:   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs;
184:   MatScalar      *av,*bv;

187:   /* compute browlengths of newmat */
188:   PetscMalloc2(mbs,PetscInt,&browlengths,mbs,PetscInt,&browstart);
189:   for (i=0; i<mbs; i++) browlengths[i] = 0;
190:   aj = a->j;
191:   for (i=0; i<mbs; i++) {
192:     nz = ai[i+1] - ai[i];
193:     aj++; /* skip diagonal */
194:     for (k=1; k<nz; k++) { /* no. of lower triangular blocks */
195:       browlengths[*aj]++; aj++;
196:     }
197:     browlengths[i] += nz;   /* no. of upper triangular blocks */
198:   }
199: 
200:   MatCreate(((PetscObject)A)->comm,&B);
201:   MatSetSizes(B,m,n,m,n);
202:   MatSetType(B,newtype);
203:   MatSeqBAIJSetPreallocation(B,bs,0,browlengths);
204:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
205: 
206:   b  = (Mat_SeqBAIJ*)(B->data);
207:   bi = b->i;
208:   bj = b->j;
209:   bv = b->a;

211:   /* set b->i */
212:   bi[0] = 0;
213:   for (i=0; i<mbs; i++){
214:     b->ilen[i]    = browlengths[i];
215:     bi[i+1]       = bi[i] + browlengths[i];
216:     browstart[i]  = bi[i];
217:   }
218:   if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs);
219: 
220:   /* set b->j and b->a */
221:   aj = a->j; av = a->a;
222:   for (i=0; i<mbs; i++) {
223:     /* diagonal block */
224:     *(bj + browstart[i]) = *aj; aj++;
225:     itmp = bs2*browstart[i];
226:     for (k=0; k<bs2; k++){
227:       *(bv + itmp + k) = *av; av++;
228:     }
229:     browstart[i]++;
230: 
231:     nz = ai[i+1] - ai[i] -1;
232:     while (nz--){
233:       /* lower triangular blocks */
234:       *(bj + browstart[*aj]) = i;
235:       itmp = bs2*browstart[*aj];
236:       for (k=0; k<bs2; k++){
237:         *(bv + itmp + k) = *(av + k);
238:       }
239:       browstart[*aj]++;

241:       /* upper triangular blocks */
242:       *(bj + browstart[i]) = *aj; aj++;
243:       itmp = bs2*browstart[i];
244:       for (k=0; k<bs2; k++){
245:         *(bv + itmp + k) = *av; av++;
246:       }
247:       browstart[i]++;
248:     }
249:   }
250:   PetscFree2(browlengths,browstart);
251:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
252:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

254:   if (reuse == MAT_REUSE_MATRIX) {
255:     MatHeaderReplace(A,B);
256:   } else {
257:     *newmat = B;
258:   }
259:   return(0);
260: }

266: PetscErrorCode  MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
267: {
268:   Mat            B;
269:   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
270:   Mat_SeqSBAIJ   *b;
272:   PetscInt       *ai=a->i,*aj,m=A->rmap->N,n=A->cmap->n,i,j,k,*bi,*bj,*browlengths;
273:   PetscInt       bs=A->rmap->bs,bs2=bs*bs,mbs=m/bs,dd;
274:   MatScalar      *av,*bv;
275:   PetscTruth     flg;

278:   if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
279:   MatMissingDiagonal_SeqBAIJ(A,&flg,&dd); /* check for missing diagonals, then mark diag */
280:   if (flg) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal %D",dd);
281: 
282:   PetscMalloc(mbs*sizeof(PetscInt),&browlengths);
283:   for (i=0; i<mbs; i++) {
284:     browlengths[i] = ai[i+1] - a->diag[i];
285:   }

287:   MatCreate(((PetscObject)A)->comm,&B);
288:   MatSetSizes(B,m,n,m,n);
289:   MatSetType(B,newtype);
290:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,browlengths);
291:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
292: 
293:   b  = (Mat_SeqSBAIJ*)(B->data);
294:   bi = b->i;
295:   bj = b->j;
296:   bv = b->a;

298:   bi[0] = 0;
299:   for (i=0; i<mbs; i++) {
300:     aj = a->j + a->diag[i];
301:     av = a->a + (a->diag[i])*bs2;
302:     for (j=0; j<browlengths[i]; j++){
303:       *bj = *aj; bj++; aj++;
304:       for (k=0; k<bs2; k++){
305:         *bv = *av; bv++; av++;
306:       }
307:     }
308:     bi[i+1]    = bi[i] + browlengths[i];
309:     b->ilen[i] = browlengths[i];
310:   }
311:   PetscFree(browlengths);
312:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
313:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

315:   if (reuse == MAT_REUSE_MATRIX) {
316:     MatHeaderReplace(A,B);
317:   } else {
318:     *newmat = B;
319:   }

321:   return(0);
322: }