Actual source code: aijsbaij.c

  1: #include <../src/mat/impls/aij/seq/aij.h>
  2: #include <../src/mat/impls/baij/seq/baij.h>
  3: #include <../src/mat/impls/sbaij/seq/sbaij.h>

  5: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
  6: {
  7:   Mat           B;
  8:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
  9:   Mat_SeqAIJ   *b;
 10:   PetscInt     *ai = a->i, *aj = a->j, m = A->rmap->N, n = A->cmap->n, i, j, k, *bi, *bj, *rowlengths, nz, *rowstart, itmp;
 11:   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = A->rmap->N / bs, diagcnt = 0;
 12:   MatScalar    *av, *bv;
 13: #if defined(PETSC_USE_COMPLEX)
 14:   const int aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
 15: #else
 16:   const int aconj = 0;
 17: #endif

 19:   PetscFunctionBegin;
 20:   /* compute rowlengths of newmat */
 21:   PetscCall(PetscMalloc2(m, &rowlengths, m + 1, &rowstart));

 23:   for (i = 0; i < mbs; i++) rowlengths[i * bs] = 0;
 24:   k = 0;
 25:   for (i = 0; i < mbs; i++) {
 26:     nz = ai[i + 1] - ai[i];
 27:     if (nz) {
 28:       rowlengths[k] += nz; /* no. of upper triangular blocks */
 29:       if (*aj == i) {
 30:         aj++;
 31:         diagcnt++;
 32:         nz--;
 33:       } /* skip diagonal */
 34:       for (j = 0; j < nz; j++) { /* no. of lower triangular blocks */
 35:         rowlengths[(*aj) * bs]++;
 36:         aj++;
 37:       }
 38:     }
 39:     rowlengths[k] *= bs;
 40:     for (j = 1; j < bs; j++) rowlengths[k + j] = rowlengths[k];
 41:     k += bs;
 42:   }

 44:   if (reuse != MAT_REUSE_MATRIX) {
 45:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
 46:     PetscCall(MatSetSizes(B, m, n, m, n));
 47:     PetscCall(MatSetType(B, MATSEQAIJ));
 48:     PetscCall(MatSeqAIJSetPreallocation(B, 0, rowlengths));
 49:     PetscCall(MatSetBlockSize(B, A->rmap->bs));
 50:   } else B = *newmat;

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

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

 69:   /* set b->j and b->a */
 70:   aj = a->j;
 71:   av = a->a;
 72:   for (i = 0; i < mbs; i++) {
 73:     nz = ai[i + 1] - ai[i];
 74:     /* diagonal block */
 75:     if (nz && *aj == i) {
 76:       nz--;
 77:       for (j = 0; j < bs; j++) { /* row i*bs+j */
 78:         itmp = i * bs + j;
 79:         for (k = 0; k < bs; k++) { /* col i*bs+k */
 80:           *(bj + rowstart[itmp]) = (*aj) * bs + k;
 81:           *(bv + rowstart[itmp]) = *(av + k * bs + j);
 82:           rowstart[itmp]++;
 83:         }
 84:       }
 85:       aj++;
 86:       av += bs2;
 87:     }

 89:     while (nz--) {
 90:       /* lower triangular blocks */
 91:       for (j = 0; j < bs; j++) { /* row (*aj)*bs+j */
 92:         itmp = (*aj) * bs + j;
 93:         for (k = 0; k < bs; k++) { /* col i*bs+k */
 94:           *(bj + rowstart[itmp]) = i * bs + k;
 95:           *(bv + rowstart[itmp]) = aconj ? PetscConj(*(av + j * bs + k)) : *(av + j * bs + k);
 96:           rowstart[itmp]++;
 97:         }
 98:       }
 99:       /* upper triangular blocks */
100:       for (j = 0; j < bs; j++) { /* row i*bs+j */
101:         itmp = i * bs + j;
102:         for (k = 0; k < bs; k++) { /* col (*aj)*bs+k */
103:           *(bj + rowstart[itmp]) = (*aj) * bs + k;
104:           *(bv + rowstart[itmp]) = *(av + k * bs + j);
105:           rowstart[itmp]++;
106:         }
107:       }
108:       aj++;
109:       av += bs2;
110:     }
111:   }
112:   PetscCall(PetscFree2(rowlengths, rowstart));
113:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
114:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

116:   if (reuse == MAT_INPLACE_MATRIX) {
117:     PetscCall(MatHeaderReplace(A, &B));
118:   } else {
119:     *newmat = B;
120:   }
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: #include <../src/mat/impls/aij/seq/aij.h>

126: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(Mat A, PetscInt **nnz)
127: {
128:   Mat_SeqAIJ     *Aa = (Mat_SeqAIJ *)A->data;
129:   PetscInt        m, n, bs = A->rmap->bs;
130:   const PetscInt *ai = Aa->i, *aj = Aa->j;

132:   PetscFunctionBegin;
133:   PetscCall(MatGetSize(A, &m, &n));

135:   if (bs == 1) {
136:     const PetscInt *adiag = Aa->diag;

138:     PetscCall(PetscMalloc1(m, nnz));
139:     for (PetscInt i = 0; i < m; i++) {
140:       if (adiag[i] == ai[i + 1]) {
141:         (*nnz)[i] = 0;
142:         for (PetscInt j = ai[i]; j < ai[i + 1]; j++) (*nnz)[i] += (aj[j] > i);
143:       } else (*nnz)[i] = ai[i + 1] - adiag[i];
144:     }
145:   } else {
146:     PetscHSetIJ    ht;
147:     PetscHashIJKey key;
148:     PetscBool      missing;

150:     PetscCall(PetscHSetIJCreate(&ht));
151:     PetscCall(PetscCalloc1(m / bs, nnz));
152:     for (PetscInt i = 0; i < m; i++) {
153:       key.i = i / bs;
154:       for (PetscInt k = ai[i]; k < ai[i + 1]; k++) {
155:         key.j = aj[k] / bs;
156:         if (key.j < key.i) continue;
157:         PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
158:         if (missing) (*nnz)[key.i]++;
159:       }
160:     }
161:     PetscCall(PetscHSetIJDestroy(&ht));
162:   }
163:   PetscFunctionReturn(PETSC_SUCCESS);
164: }

166: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
167: {
168:   Mat           B;
169:   Mat_SeqAIJ   *a = (Mat_SeqAIJ *)A->data;
170:   Mat_SeqSBAIJ *b;
171:   PetscInt     *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->N, i, j, *bi, *bj, *rowlengths, bs = A->rmap->bs;
172:   MatScalar    *av, *bv;
173:   PetscBool     miss = PETSC_FALSE;

175:   PetscFunctionBegin;
176: #if !defined(PETSC_USE_COMPLEX)
177:   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
178: #else
179:   PetscCheck(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be either symmetric or hermitian. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE) and/or MatSetOption(mat,MAT_HERMITIAN,PETSC_TRUE)");
180: #endif
181:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");

183:   if (bs == 1) {
184:     PetscCall(PetscMalloc1(m, &rowlengths));
185:     for (i = 0; i < m; i++) {
186:       if (a->diag[i] == ai[i + 1]) {             /* missing diagonal */
187:         rowlengths[i] = (ai[i + 1] - ai[i]) + 1; /* allocate some extra space */
188:         miss          = PETSC_TRUE;
189:       } else {
190:         rowlengths[i] = (ai[i + 1] - a->diag[i]);
191:       }
192:     }
193:   } else PetscCall(MatConvert_SeqAIJ_SeqSBAIJ_Preallocate(A, &rowlengths));
194:   if (reuse != MAT_REUSE_MATRIX) {
195:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
196:     PetscCall(MatSetSizes(B, m, n, m, n));
197:     PetscCall(MatSetType(B, MATSEQSBAIJ));
198:     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, rowlengths));
199:   } else B = *newmat;

201:   if (bs == 1 && !miss) {
202:     b  = (Mat_SeqSBAIJ *)B->data;
203:     bi = b->i;
204:     bj = b->j;
205:     bv = b->a;

207:     bi[0] = 0;
208:     for (i = 0; i < m; i++) {
209:       aj = a->j + a->diag[i];
210:       av = a->a + a->diag[i];
211:       for (j = 0; j < rowlengths[i]; j++) {
212:         *bj = *aj;
213:         bj++;
214:         aj++;
215:         *bv = *av;
216:         bv++;
217:         av++;
218:       }
219:       bi[i + 1]  = bi[i] + rowlengths[i];
220:       b->ilen[i] = rowlengths[i];
221:     }
222:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
223:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
224:   } else {
225:     PetscCall(MatSetOption(B, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
226:     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
227:     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
228:     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
229:     PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
230:   }
231:   PetscCall(PetscFree(rowlengths));
232:   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
233:   else *newmat = B;
234:   PetscFunctionReturn(PETSC_SUCCESS);
235: }

237: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
238: {
239:   Mat           B;
240:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
241:   Mat_SeqBAIJ  *b;
242:   PetscInt     *ai = a->i, *aj = a->j, m = A->rmap->N, n = A->cmap->n, i, k, *bi, *bj, *browlengths, nz, *browstart, itmp;
243:   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs, col, row;
244:   MatScalar    *av, *bv;
245: #if defined(PETSC_USE_COMPLEX)
246:   const int aconj = A->hermitian == PETSC_BOOL3_TRUE ? 1 : 0;
247: #else
248:   const int aconj = 0;
249: #endif

251:   PetscFunctionBegin;
252:   /* compute browlengths of newmat */
253:   PetscCall(PetscMalloc2(mbs, &browlengths, mbs, &browstart));
254:   for (i = 0; i < mbs; i++) browlengths[i] = 0;
255:   for (i = 0; i < mbs; i++) {
256:     nz = ai[i + 1] - ai[i];
257:     aj++;                      /* skip diagonal */
258:     for (k = 1; k < nz; k++) { /* no. of lower triangular blocks */
259:       browlengths[*aj]++;
260:       aj++;
261:     }
262:     browlengths[i] += nz; /* no. of upper triangular blocks */
263:   }

265:   if (reuse != MAT_REUSE_MATRIX) {
266:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
267:     PetscCall(MatSetSizes(B, m, n, m, n));
268:     PetscCall(MatSetType(B, MATSEQBAIJ));
269:     PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, browlengths));
270:   } else B = *newmat;

272:   b  = (Mat_SeqBAIJ *)B->data;
273:   bi = b->i;
274:   bj = b->j;
275:   bv = b->a;

277:   /* set b->i */
278:   bi[0] = 0;
279:   for (i = 0; i < mbs; i++) {
280:     b->ilen[i]   = browlengths[i];
281:     bi[i + 1]    = bi[i] + browlengths[i];
282:     browstart[i] = bi[i];
283:   }
284:   PetscCheck(bi[mbs] == 2 * a->nz - mbs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bi[mbs]: %" PetscInt_FMT " != 2*a->nz - mbs: %" PetscInt_FMT, bi[mbs], 2 * a->nz - mbs);

286:   /* set b->j and b->a */
287:   aj = a->j;
288:   av = a->a;
289:   for (i = 0; i < mbs; i++) {
290:     /* diagonal block */
291:     *(bj + browstart[i]) = *aj;
292:     aj++;

294:     itmp = bs2 * browstart[i];
295:     for (k = 0; k < bs2; k++) {
296:       *(bv + itmp + k) = *av;
297:       av++;
298:     }
299:     browstart[i]++;

301:     nz = ai[i + 1] - ai[i] - 1;
302:     while (nz--) {
303:       /* lower triangular blocks - transpose blocks of A */
304:       *(bj + browstart[*aj]) = i; /* block col index */

306:       itmp = bs2 * browstart[*aj]; /* row index */
307:       for (col = 0; col < bs; col++) {
308:         k = col;
309:         for (row = 0; row < bs; row++) {
310:           bv[itmp + col * bs + row] = aconj ? PetscConj(av[k]) : av[k];
311:           k += bs;
312:         }
313:       }
314:       browstart[*aj]++;

316:       /* upper triangular blocks */
317:       *(bj + browstart[i]) = *aj;
318:       aj++;

320:       itmp = bs2 * browstart[i];
321:       for (k = 0; k < bs2; k++) bv[itmp + k] = av[k];
322:       av += bs2;
323:       browstart[i]++;
324:     }
325:   }
326:   PetscCall(PetscFree2(browlengths, browstart));
327:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
328:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

330:   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
331:   else *newmat = B;
332:   PetscFunctionReturn(PETSC_SUCCESS);
333: }

335: PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
336: {
337:   Mat           B;
338:   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *)A->data;
339:   Mat_SeqSBAIJ *b;
340:   PetscInt     *ai = a->i, *aj, m = A->rmap->N, n = A->cmap->n, i, j, k, *bi, *bj, *browlengths;
341:   PetscInt      bs = A->rmap->bs, bs2 = bs * bs, mbs = m / bs, dd;
342:   MatScalar    *av, *bv;
343:   PetscBool     flg;

345:   PetscFunctionBegin;
346:   PetscCheck(A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Matrix must be symmetric. Call MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE)");
347:   PetscCheck(n == m, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Matrix must be square");
348:   PetscCall(MatMissingDiagonal_SeqBAIJ(A, &flg, &dd)); /* check for missing diagonals, then mark diag */
349:   PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal %" PetscInt_FMT, dd);

351:   PetscCall(PetscMalloc1(mbs, &browlengths));
352:   for (i = 0; i < mbs; i++) browlengths[i] = ai[i + 1] - a->diag[i];

354:   if (reuse != MAT_REUSE_MATRIX) {
355:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
356:     PetscCall(MatSetSizes(B, m, n, m, n));
357:     PetscCall(MatSetType(B, MATSEQSBAIJ));
358:     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, browlengths));
359:   } else B = *newmat;

361:   b  = (Mat_SeqSBAIJ *)B->data;
362:   bi = b->i;
363:   bj = b->j;
364:   bv = b->a;

366:   bi[0] = 0;
367:   for (i = 0; i < mbs; i++) {
368:     aj = a->j + a->diag[i];
369:     av = a->a + (a->diag[i]) * bs2;
370:     for (j = 0; j < browlengths[i]; j++) {
371:       *bj = *aj;
372:       bj++;
373:       aj++;
374:       for (k = 0; k < bs2; k++) {
375:         *bv = *av;
376:         bv++;
377:         av++;
378:       }
379:     }
380:     bi[i + 1]  = bi[i] + browlengths[i];
381:     b->ilen[i] = browlengths[i];
382:   }
383:   PetscCall(PetscFree(browlengths));
384:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
385:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

387:   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
388:   else *newmat = B;
389:   PetscFunctionReturn(PETSC_SUCCESS);
390: }