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:   PetscMalloc((2*m+1)*sizeof(PetscInt),&rowlengths);
 23:   rowstart = rowlengths + m;
 24: 
 25:   for (i=0; i<mbs; i++) rowlengths[i*bs] = 0;
 26:   aj = a->j;
 27:   k = 0;
 28:   for (i=0; i<mbs; i++) {
 29:     nz = ai[i+1] - ai[i];
 30:     aj++; /* skip diagonal */
 31:     for (j=1; j<nz; j++) { /* no. of lower triangular blocks */
 32:       rowlengths[(*aj)*bs]++; aj++;
 33:     }
 34:     rowlengths[k] += nz;   /* no. of upper triangular blocks */
 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:   MatCreate(((PetscObject)A)->comm,&B);
 44:   MatSetSizes(B,m,n,m,n);
 45:   MatSetType(B,newtype);
 46:   MatSeqAIJSetPreallocation(B,0,rowlengths);
 47:   MatSetOption(B,MAT_ROW_ORIENTED,PETSC_FALSE);
 48:   B->rmap->bs = A->rmap->bs;

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

 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 - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-mbs: %D\n",bi[mbs],2*a->nz - mbs);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

320:   return(0);
321: }