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,diagcnt=0;
 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:     if (nz) {
 30:       rowlengths[k] += nz;   /* no. of upper triangular blocks */
 31:       if (*aj == i) {aj++;diagcnt++;nz--;} /* skip diagonal */
 32:       for (j=0; j<nz; j++) { /* no. of lower triangular blocks */
 33:         rowlengths[(*aj)*bs]++; aj++;
 34:       }
 35:     }
 36:     rowlengths[k] *= bs;
 37:     for (j=1; j<bs; j++) {
 38:       rowlengths[k+j] = rowlengths[k];
 39:     }
 40:     k += bs;
 41:     /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */
 42:   }
 43: 
 44:   MatCreate(((PetscObject)A)->comm,&B);
 45:   MatSetSizes(B,m,n,m,n);
 46:   MatSetType(B,newtype);
 47:   MatSeqAIJSetPreallocation(B,0,rowlengths);
 48:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);
 49:   B->rmap->bs = A->rmap->bs;

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

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

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

111:   if (reuse == MAT_REUSE_MATRIX) {
112:     MatHeaderReplace(A,B);
113:   } else {
114:     *newmat = B;
115:   }
116:   return(0);
117: }

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

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

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

143:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_TRUE);
144: 
145:   b  = (Mat_SeqSBAIJ*)(B->data);
146:   bi = b->i;
147:   bj = b->j;
148:   bv = b->a;

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:   PetscFree(rowlengths);
163:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
164:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
165:   if (A->hermitian){
166:     MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
167:   }

169:   if (reuse == MAT_REUSE_MATRIX) {
170:     MatHeaderReplace(A,B);
171:   } else {
172:     *newmat = B;
173:   }
174:   return(0);
175: }

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

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

216:   /* set b->i */
217:   bi[0] = 0;
218:   for (i=0; i<mbs; i++){
219:     b->ilen[i]    = browlengths[i];
220:     bi[i+1]       = bi[i] + browlengths[i];
221:     browstart[i]  = bi[i];
222:   }
223:   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);
224: 
225:   /* set b->j and b->a */
226:   aj = a->j; av = a->a;
227:   for (i=0; i<mbs; i++) {
228:     /* diagonal block */
229:     *(bj + browstart[i]) = *aj; aj++;
230:     itmp = bs2*browstart[i];
231:     for (k=0; k<bs2; k++){
232:       *(bv + itmp + k) = *av; av++;
233:     }
234:     browstart[i]++;
235: 
236:     nz = ai[i+1] - ai[i] -1;
237:     while (nz--){
238:       /* lower triangular blocks */
239:       *(bj + browstart[*aj]) = i;
240:       itmp = bs2*browstart[*aj];
241:       for (k=0; k<bs2; k++){
242:         *(bv + itmp + k) = *(av + k);
243:       }
244:       browstart[*aj]++;

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; av++;
251:       }
252:       browstart[i]++;
253:     }
254:   }
255:   PetscFree2(browlengths,browstart);
256:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
257:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

259:   if (reuse == MAT_REUSE_MATRIX) {
260:     MatHeaderReplace(A,B);
261:   } else {
262:     *newmat = B;
263:   }
264:   return(0);
265: }

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

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

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

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

320:   if (reuse == MAT_REUSE_MATRIX) {
321:     MatHeaderReplace(A,B);
322:   } else {
323:     *newmat = B;
324:   }

326:   return(0);
327: }