Actual source code: sbaij.c

  1: /*
  2:     Defines the basic matrix operations for the SBAIJ (compressed row)
  3:   matrix storage format.
  4: */
  5: #include <../src/mat/impls/baij/seq/baij.h>
  6: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  7: #include <petscblaslapack.h>

  9: #include <../src/mat/impls/sbaij/seq/relax.h>
 10: #define USESHORT
 11: #include <../src/mat/impls/sbaij/seq/relax.h>

 13: /* defines MatSetValues_Seq_Hash(), MatAssemblyEnd_Seq_Hash(), MatSetUp_Seq_Hash() */
 14: #define TYPE SBAIJ
 15: #define TYPE_SBAIJ
 16: #define TYPE_BS
 17: #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h"
 18: #undef TYPE_BS
 19: #define TYPE_BS _BS
 20: #define TYPE_BS_ON
 21: #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h"
 22: #undef TYPE_BS
 23: #undef TYPE_SBAIJ
 24: #include "../src/mat/impls/aij/seq/seqhashmat.h"
 25: #undef TYPE
 26: #undef TYPE_BS_ON

 28: #if defined(PETSC_HAVE_ELEMENTAL)
 29: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
 30: #endif
 31: #if defined(PETSC_HAVE_SCALAPACK)
 32: PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
 33: #endif
 34: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat, MatType, MatReuse, Mat *);

 36: /*
 37:      Checks for missing diagonals
 38: */
 39: static PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A, PetscBool *missing, PetscInt *dd)
 40: {
 41:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
 42:   PetscInt     *diag, *ii = a->i, i;

 44:   PetscFunctionBegin;
 45:   PetscCall(MatMarkDiagonal_SeqSBAIJ(A));
 46:   *missing = PETSC_FALSE;
 47:   if (A->rmap->n > 0 && !ii) {
 48:     *missing = PETSC_TRUE;
 49:     if (dd) *dd = 0;
 50:     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
 51:   } else {
 52:     diag = a->diag;
 53:     for (i = 0; i < a->mbs; i++) {
 54:       if (diag[i] >= ii[i + 1]) {
 55:         *missing = PETSC_TRUE;
 56:         if (dd) *dd = i;
 57:         break;
 58:       }
 59:     }
 60:   }
 61:   PetscFunctionReturn(PETSC_SUCCESS);
 62: }

 64: PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
 65: {
 66:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
 67:   PetscInt      i, j;

 69:   PetscFunctionBegin;
 70:   if (!a->diag) {
 71:     PetscCall(PetscMalloc1(a->mbs, &a->diag));
 72:     a->free_diag = PETSC_TRUE;
 73:   }
 74:   for (i = 0; i < a->mbs; i++) {
 75:     a->diag[i] = a->i[i + 1];
 76:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
 77:       if (a->j[j] == i) {
 78:         a->diag[i] = j;
 79:         break;
 80:       }
 81:     }
 82:   }
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

 86: static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *nn, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
 87: {
 88:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
 89:   PetscInt      i, j, n = a->mbs, nz = a->i[n], *tia, *tja, bs = A->rmap->bs, k, l, cnt;
 90:   PetscInt    **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;

 92:   PetscFunctionBegin;
 93:   *nn = n;
 94:   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
 95:   if (symmetric) {
 96:     PetscCall(MatToSymmetricIJ_SeqAIJ(n, a->i, a->j, PETSC_FALSE, 0, 0, &tia, &tja));
 97:     nz = tia[n];
 98:   } else {
 99:     tia = a->i;
100:     tja = a->j;
101:   }

103:   if (!blockcompressed && bs > 1) {
104:     (*nn) *= bs;
105:     /* malloc & create the natural set of indices */
106:     PetscCall(PetscMalloc1((n + 1) * bs, ia));
107:     if (n) {
108:       (*ia)[0] = oshift;
109:       for (j = 1; j < bs; j++) (*ia)[j] = (tia[1] - tia[0]) * bs + (*ia)[j - 1];
110:     }

112:     for (i = 1; i < n; i++) {
113:       (*ia)[i * bs] = (tia[i] - tia[i - 1]) * bs + (*ia)[i * bs - 1];
114:       for (j = 1; j < bs; j++) (*ia)[i * bs + j] = (tia[i + 1] - tia[i]) * bs + (*ia)[i * bs + j - 1];
115:     }
116:     if (n) (*ia)[n * bs] = (tia[n] - tia[n - 1]) * bs + (*ia)[n * bs - 1];

118:     if (inja) {
119:       PetscCall(PetscMalloc1(nz * bs * bs, ja));
120:       cnt = 0;
121:       for (i = 0; i < n; i++) {
122:         for (j = 0; j < bs; j++) {
123:           for (k = tia[i]; k < tia[i + 1]; k++) {
124:             for (l = 0; l < bs; l++) (*ja)[cnt++] = bs * tja[k] + l;
125:           }
126:         }
127:       }
128:     }

130:     if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */
131:       PetscCall(PetscFree(tia));
132:       PetscCall(PetscFree(tja));
133:     }
134:   } else if (oshift == 1) {
135:     if (symmetric) {
136:       nz = tia[A->rmap->n / bs];
137:       /*  add 1 to i and j indices */
138:       for (i = 0; i < A->rmap->n / bs + 1; i++) tia[i] = tia[i] + 1;
139:       *ia = tia;
140:       if (ja) {
141:         for (i = 0; i < nz; i++) tja[i] = tja[i] + 1;
142:         *ja = tja;
143:       }
144:     } else {
145:       nz = a->i[A->rmap->n / bs];
146:       /* malloc space and  add 1 to i and j indices */
147:       PetscCall(PetscMalloc1(A->rmap->n / bs + 1, ia));
148:       for (i = 0; i < A->rmap->n / bs + 1; i++) (*ia)[i] = a->i[i] + 1;
149:       if (ja) {
150:         PetscCall(PetscMalloc1(nz, ja));
151:         for (i = 0; i < nz; i++) (*ja)[i] = a->j[i] + 1;
152:       }
153:     }
154:   } else {
155:     *ia = tia;
156:     if (ja) *ja = tja;
157:   }
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
162: {
163:   PetscFunctionBegin;
164:   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
165:   if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) {
166:     PetscCall(PetscFree(*ia));
167:     if (ja) PetscCall(PetscFree(*ja));
168:   }
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
173: {
174:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;

176:   PetscFunctionBegin;
177:   if (A->hash_active) {
178:     PetscInt bs;
179:     A->ops[0] = a->cops;
180:     PetscCall(PetscHMapIJVDestroy(&a->ht));
181:     PetscCall(MatGetBlockSize(A, &bs));
182:     if (bs > 1) PetscCall(PetscHSetIJDestroy(&a->bht));
183:     PetscCall(PetscFree(a->dnz));
184:     PetscCall(PetscFree(a->bdnz));
185:     A->hash_active = PETSC_FALSE;
186:   }
187:   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->N, a->nz));
188:   PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i));
189:   if (a->free_diag) PetscCall(PetscFree(a->diag));
190:   PetscCall(ISDestroy(&a->row));
191:   PetscCall(ISDestroy(&a->col));
192:   PetscCall(ISDestroy(&a->icol));
193:   PetscCall(PetscFree(a->idiag));
194:   PetscCall(PetscFree(a->inode.size));
195:   if (a->free_imax_ilen) PetscCall(PetscFree2(a->imax, a->ilen));
196:   PetscCall(PetscFree(a->solve_work));
197:   PetscCall(PetscFree(a->sor_work));
198:   PetscCall(PetscFree(a->solves_work));
199:   PetscCall(PetscFree(a->mult_work));
200:   PetscCall(PetscFree(a->saved_values));
201:   if (a->free_jshort) PetscCall(PetscFree(a->jshort));
202:   PetscCall(PetscFree(a->inew));
203:   PetscCall(MatDestroy(&a->parent));
204:   PetscCall(PetscFree(A->data));

206:   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
207:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJGetArray_C", NULL));
208:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJRestoreArray_C", NULL));
209:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
210:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
211:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetColumnIndices_C", NULL));
212:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_seqaij_C", NULL));
213:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_seqbaij_C", NULL));
214:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetPreallocation_C", NULL));
215:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetPreallocationCSR_C", NULL));
216: #if defined(PETSC_HAVE_ELEMENTAL)
217:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_elemental_C", NULL));
218: #endif
219: #if defined(PETSC_HAVE_SCALAPACK)
220:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_scalapack_C", NULL));
221: #endif
222:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: static PetscErrorCode MatSetOption_SeqSBAIJ(Mat A, MatOption op, PetscBool flg)
227: {
228:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
229: #if defined(PETSC_USE_COMPLEX)
230:   PetscInt bs;
231: #endif

233:   PetscFunctionBegin;
234: #if defined(PETSC_USE_COMPLEX)
235:   PetscCall(MatGetBlockSize(A, &bs));
236: #endif
237:   switch (op) {
238:   case MAT_ROW_ORIENTED:
239:     a->roworiented = flg;
240:     break;
241:   case MAT_KEEP_NONZERO_PATTERN:
242:     a->keepnonzeropattern = flg;
243:     break;
244:   case MAT_NEW_NONZERO_LOCATIONS:
245:     a->nonew = (flg ? 0 : 1);
246:     break;
247:   case MAT_NEW_NONZERO_LOCATION_ERR:
248:     a->nonew = (flg ? -1 : 0);
249:     break;
250:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
251:     a->nonew = (flg ? -2 : 0);
252:     break;
253:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
254:     a->nounused = (flg ? -1 : 0);
255:     break;
256:   case MAT_FORCE_DIAGONAL_ENTRIES:
257:   case MAT_IGNORE_OFF_PROC_ENTRIES:
258:   case MAT_USE_HASH_TABLE:
259:   case MAT_SORTED_FULL:
260:     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
261:     break;
262:   case MAT_HERMITIAN:
263: #if defined(PETSC_USE_COMPLEX)
264:     if (flg) { /* disable transpose ops */
265:       PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for Hermitian with block size greater than 1");
266:       A->ops->multtranspose    = NULL;
267:       A->ops->multtransposeadd = NULL;
268:       A->symmetric             = PETSC_BOOL3_FALSE;
269:     }
270: #endif
271:     break;
272:   case MAT_SYMMETRIC:
273:   case MAT_SPD:
274: #if defined(PETSC_USE_COMPLEX)
275:     if (flg) { /* An hermitian and symmetric matrix has zero imaginary part (restore back transpose ops) */
276:       A->ops->multtranspose    = A->ops->mult;
277:       A->ops->multtransposeadd = A->ops->multadd;
278:     }
279: #endif
280:     break;
281:     /* These options are handled directly by MatSetOption() */
282:   case MAT_STRUCTURALLY_SYMMETRIC:
283:   case MAT_SYMMETRY_ETERNAL:
284:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
285:   case MAT_STRUCTURE_ONLY:
286:   case MAT_SPD_ETERNAL:
287:     /* These options are handled directly by MatSetOption() */
288:     break;
289:   case MAT_IGNORE_LOWER_TRIANGULAR:
290:     a->ignore_ltriangular = flg;
291:     break;
292:   case MAT_ERROR_LOWER_TRIANGULAR:
293:     a->ignore_ltriangular = flg;
294:     break;
295:   case MAT_GETROW_UPPERTRIANGULAR:
296:     a->getrow_utriangular = flg;
297:     break;
298:   case MAT_SUBMAT_SINGLEIS:
299:     break;
300:   default:
301:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
302:   }
303:   PetscFunctionReturn(PETSC_SUCCESS);
304: }

306: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
307: {
308:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;

310:   PetscFunctionBegin;
311:   PetscCheck(!A || a->getrow_utriangular, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");

313:   /* Get the upper triangular part of the row */
314:   PetscCall(MatGetRow_SeqBAIJ_private(A, row, nz, idx, v, a->i, a->j, a->a));
315:   PetscFunctionReturn(PETSC_SUCCESS);
316: }

318: PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
319: {
320:   PetscFunctionBegin;
321:   if (idx) PetscCall(PetscFree(*idx));
322:   if (v) PetscCall(PetscFree(*v));
323:   PetscFunctionReturn(PETSC_SUCCESS);
324: }

326: static PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
327: {
328:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;

330:   PetscFunctionBegin;
331:   a->getrow_utriangular = PETSC_TRUE;
332:   PetscFunctionReturn(PETSC_SUCCESS);
333: }

335: static PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
336: {
337:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;

339:   PetscFunctionBegin;
340:   a->getrow_utriangular = PETSC_FALSE;
341:   PetscFunctionReturn(PETSC_SUCCESS);
342: }

344: static PetscErrorCode MatTranspose_SeqSBAIJ(Mat A, MatReuse reuse, Mat *B)
345: {
346:   PetscFunctionBegin;
347:   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
348:   if (reuse == MAT_INITIAL_MATRIX) {
349:     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, B));
350:   } else if (reuse == MAT_REUSE_MATRIX) {
351:     PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
352:   }
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A, PetscViewer viewer)
357: {
358:   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ *)A->data;
359:   PetscInt          i, j, bs = A->rmap->bs, k, l, bs2 = a->bs2;
360:   PetscViewerFormat format;
361:   PetscInt         *diag;
362:   const char       *matname;

364:   PetscFunctionBegin;
365:   PetscCall(PetscViewerGetFormat(viewer, &format));
366:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
367:     PetscCall(PetscViewerASCIIPrintf(viewer, "  block size is %" PetscInt_FMT "\n", bs));
368:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
369:     Mat aij;

371:     if (A->factortype && bs > 1) {
372:       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n"));
373:       PetscFunctionReturn(PETSC_SUCCESS);
374:     }
375:     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &aij));
376:     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
377:     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)aij, matname));
378:     PetscCall(MatView_SeqAIJ(aij, viewer));
379:     PetscCall(MatDestroy(&aij));
380:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
381:     Mat B;

383:     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
384:     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
385:     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)B, matname));
386:     PetscCall(MatView_SeqAIJ(B, viewer));
387:     PetscCall(MatDestroy(&B));
388:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
389:     PetscFunctionReturn(PETSC_SUCCESS);
390:   } else {
391:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
392:     if (A->factortype) { /* for factored matrix */
393:       PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "matrix is factored with bs>1. Not implemented yet");

395:       diag = a->diag;
396:       for (i = 0; i < a->mbs; i++) { /* for row block i */
397:         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
398:         /* diagonal entry */
399: #if defined(PETSC_USE_COMPLEX)
400:         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
401:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]]), (double)PetscImaginaryPart(1.0 / a->a[diag[i]])));
402:         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
403:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]]), -(double)PetscImaginaryPart(1.0 / a->a[diag[i]])));
404:         } else {
405:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]])));
406:         }
407: #else
408:         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[diag[i]], (double)(1.0 / a->a[diag[i]])));
409: #endif
410:         /* off-diagonal entries */
411:         for (k = a->i[i]; k < a->i[i + 1] - 1; k++) {
412: #if defined(PETSC_USE_COMPLEX)
413:           if (PetscImaginaryPart(a->a[k]) > 0.0) {
414:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", bs * a->j[k], (double)PetscRealPart(a->a[k]), (double)PetscImaginaryPart(a->a[k])));
415:           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
416:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", bs * a->j[k], (double)PetscRealPart(a->a[k]), -(double)PetscImaginaryPart(a->a[k])));
417:           } else {
418:             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k], (double)PetscRealPart(a->a[k])));
419:           }
420: #else
421:           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[k], (double)a->a[k]));
422: #endif
423:         }
424:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
425:       }

427:     } else {                         /* for non-factored matrix */
428:       for (i = 0; i < a->mbs; i++) { /* for row block i */
429:         for (j = 0; j < bs; j++) {   /* for row bs*i + j */
430:           PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i * bs + j));
431:           for (k = a->i[i]; k < a->i[i + 1]; k++) { /* for column block */
432:             for (l = 0; l < bs; l++) {              /* for column */
433: #if defined(PETSC_USE_COMPLEX)
434:               if (PetscImaginaryPart(a->a[bs2 * k + l * bs + j]) > 0.0) {
435:                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j]), (double)PetscImaginaryPart(a->a[bs2 * k + l * bs + j])));
436:               } else if (PetscImaginaryPart(a->a[bs2 * k + l * bs + j]) < 0.0) {
437:                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j]), -(double)PetscImaginaryPart(a->a[bs2 * k + l * bs + j])));
438:               } else {
439:                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j])));
440:               }
441: #else
442:               PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k] + l, (double)a->a[bs2 * k + l * bs + j]));
443: #endif
444:             }
445:           }
446:           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
447:         }
448:       }
449:     }
450:     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
451:   }
452:   PetscCall(PetscViewerFlush(viewer));
453:   PetscFunctionReturn(PETSC_SUCCESS);
454: }

456: #include <petscdraw.h>
457: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw, void *Aa)
458: {
459:   Mat           A = (Mat)Aa;
460:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
461:   PetscInt      row, i, j, k, l, mbs = a->mbs, color, bs = A->rmap->bs, bs2 = a->bs2;
462:   PetscReal     xl, yl, xr, yr, x_l, x_r, y_l, y_r;
463:   MatScalar    *aa;
464:   PetscViewer   viewer;

466:   PetscFunctionBegin;
467:   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
468:   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));

470:   /* loop over matrix elements drawing boxes */

472:   PetscDrawCollectiveBegin(draw);
473:   PetscCall(PetscDrawString(draw, .3 * (xl + xr), .3 * (yl + yr), PETSC_DRAW_BLACK, "symmetric"));
474:   /* Blue for negative, Cyan for zero and  Red for positive */
475:   color = PETSC_DRAW_BLUE;
476:   for (i = 0, row = 0; i < mbs; i++, row += bs) {
477:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
478:       y_l = A->rmap->N - row - 1.0;
479:       y_r = y_l + 1.0;
480:       x_l = a->j[j] * bs;
481:       x_r = x_l + 1.0;
482:       aa  = a->a + j * bs2;
483:       for (k = 0; k < bs; k++) {
484:         for (l = 0; l < bs; l++) {
485:           if (PetscRealPart(*aa++) >= 0.) continue;
486:           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
487:         }
488:       }
489:     }
490:   }
491:   color = PETSC_DRAW_CYAN;
492:   for (i = 0, row = 0; i < mbs; i++, row += bs) {
493:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
494:       y_l = A->rmap->N - row - 1.0;
495:       y_r = y_l + 1.0;
496:       x_l = a->j[j] * bs;
497:       x_r = x_l + 1.0;
498:       aa  = a->a + j * bs2;
499:       for (k = 0; k < bs; k++) {
500:         for (l = 0; l < bs; l++) {
501:           if (PetscRealPart(*aa++) != 0.) continue;
502:           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
503:         }
504:       }
505:     }
506:   }
507:   color = PETSC_DRAW_RED;
508:   for (i = 0, row = 0; i < mbs; i++, row += bs) {
509:     for (j = a->i[i]; j < a->i[i + 1]; j++) {
510:       y_l = A->rmap->N - row - 1.0;
511:       y_r = y_l + 1.0;
512:       x_l = a->j[j] * bs;
513:       x_r = x_l + 1.0;
514:       aa  = a->a + j * bs2;
515:       for (k = 0; k < bs; k++) {
516:         for (l = 0; l < bs; l++) {
517:           if (PetscRealPart(*aa++) <= 0.) continue;
518:           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
519:         }
520:       }
521:     }
522:   }
523:   PetscDrawCollectiveEnd(draw);
524:   PetscFunctionReturn(PETSC_SUCCESS);
525: }

527: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A, PetscViewer viewer)
528: {
529:   PetscReal xl, yl, xr, yr, w, h;
530:   PetscDraw draw;
531:   PetscBool isnull;

533:   PetscFunctionBegin;
534:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
535:   PetscCall(PetscDrawIsNull(draw, &isnull));
536:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

538:   xr = A->rmap->N;
539:   yr = A->rmap->N;
540:   h  = yr / 10.0;
541:   w  = xr / 10.0;
542:   xr += w;
543:   yr += h;
544:   xl = -w;
545:   yl = -h;
546:   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
547:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
548:   PetscCall(PetscDrawZoom(draw, MatView_SeqSBAIJ_Draw_Zoom, A));
549:   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
550:   PetscCall(PetscDrawSave(draw));
551:   PetscFunctionReturn(PETSC_SUCCESS);
552: }

554: /* Used for both MPIBAIJ and MPISBAIJ matrices */
555: #define MatView_SeqSBAIJ_Binary MatView_SeqBAIJ_Binary

557: PetscErrorCode MatView_SeqSBAIJ(Mat A, PetscViewer viewer)
558: {
559:   PetscBool iascii, isbinary, isdraw;

561:   PetscFunctionBegin;
562:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
563:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
564:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
565:   if (iascii) {
566:     PetscCall(MatView_SeqSBAIJ_ASCII(A, viewer));
567:   } else if (isbinary) {
568:     PetscCall(MatView_SeqSBAIJ_Binary(A, viewer));
569:   } else if (isdraw) {
570:     PetscCall(MatView_SeqSBAIJ_Draw(A, viewer));
571:   } else {
572:     Mat         B;
573:     const char *matname;
574:     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
575:     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
576:     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)B, matname));
577:     PetscCall(MatView(B, viewer));
578:     PetscCall(MatDestroy(&B));
579:   }
580:   PetscFunctionReturn(PETSC_SUCCESS);
581: }

583: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
584: {
585:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
586:   PetscInt     *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
587:   PetscInt     *ai = a->i, *ailen = a->ilen;
588:   PetscInt      brow, bcol, ridx, cidx, bs = A->rmap->bs, bs2 = a->bs2;
589:   MatScalar    *ap, *aa = a->a;

591:   PetscFunctionBegin;
592:   for (k = 0; k < m; k++) { /* loop over rows */
593:     row  = im[k];
594:     brow = row / bs;
595:     if (row < 0) {
596:       v += n;
597:       continue;
598:     } /* negative row */
599:     PetscCheck(row < A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->N - 1);
600:     rp   = aj + ai[brow];
601:     ap   = aa + bs2 * ai[brow];
602:     nrow = ailen[brow];
603:     for (l = 0; l < n; l++) { /* loop over columns */
604:       if (in[l] < 0) {
605:         v++;
606:         continue;
607:       } /* negative column */
608:       PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1);
609:       col  = in[l];
610:       bcol = col / bs;
611:       cidx = col % bs;
612:       ridx = row % bs;
613:       high = nrow;
614:       low  = 0; /* assume unsorted */
615:       while (high - low > 5) {
616:         t = (low + high) / 2;
617:         if (rp[t] > bcol) high = t;
618:         else low = t;
619:       }
620:       for (i = low; i < high; i++) {
621:         if (rp[i] > bcol) break;
622:         if (rp[i] == bcol) {
623:           *v++ = ap[bs2 * i + bs * cidx + ridx];
624:           goto finished;
625:         }
626:       }
627:       *v++ = 0.0;
628:     finished:;
629:     }
630:   }
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: static PetscErrorCode MatPermute_SeqSBAIJ(Mat A, IS rowp, IS colp, Mat *B)
635: {
636:   Mat       C;
637:   PetscBool flg = (PetscBool)(rowp == colp);

639:   PetscFunctionBegin;
640:   PetscCall(MatConvert(A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &C));
641:   PetscCall(MatPermute(C, rowp, colp, B));
642:   PetscCall(MatDestroy(&C));
643:   if (!flg) PetscCall(ISEqual(rowp, colp, &flg));
644:   if (flg) PetscCall(MatConvert(*B, MATSEQSBAIJ, MAT_INPLACE_MATRIX, B));
645:   PetscFunctionReturn(PETSC_SUCCESS);
646: }

648: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
649: {
650:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ *)A->data;
651:   PetscInt          *rp, k, low, high, t, ii, jj, row, nrow, i, col, l, rmax, N, lastcol = -1;
652:   PetscInt          *imax = a->imax, *ai = a->i, *ailen = a->ilen;
653:   PetscInt          *aj = a->j, nonew = a->nonew, bs2 = a->bs2, bs = A->rmap->bs, stepval;
654:   PetscBool          roworiented = a->roworiented;
655:   const PetscScalar *value       = v;
656:   MatScalar         *ap, *aa = a->a, *bap;

658:   PetscFunctionBegin;
659:   if (roworiented) stepval = (n - 1) * bs;
660:   else stepval = (m - 1) * bs;
661:   for (k = 0; k < m; k++) { /* loop over added rows */
662:     row = im[k];
663:     if (row < 0) continue;
664:     PetscCheck(row < a->mbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block index row too large %" PetscInt_FMT " max %" PetscInt_FMT, row, a->mbs - 1);
665:     rp   = aj + ai[row];
666:     ap   = aa + bs2 * ai[row];
667:     rmax = imax[row];
668:     nrow = ailen[row];
669:     low  = 0;
670:     high = nrow;
671:     for (l = 0; l < n; l++) { /* loop over added columns */
672:       if (in[l] < 0) continue;
673:       col = in[l];
674:       PetscCheck(col < a->nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block index column too large %" PetscInt_FMT " max %" PetscInt_FMT, col, a->nbs - 1);
675:       if (col < row) {
676:         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
677:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
678:       }
679:       if (roworiented) value = v + k * (stepval + bs) * bs + l * bs;
680:       else value = v + l * (stepval + bs) * bs + k * bs;

682:       if (col <= lastcol) low = 0;
683:       else high = nrow;

685:       lastcol = col;
686:       while (high - low > 7) {
687:         t = (low + high) / 2;
688:         if (rp[t] > col) high = t;
689:         else low = t;
690:       }
691:       for (i = low; i < high; i++) {
692:         if (rp[i] > col) break;
693:         if (rp[i] == col) {
694:           bap = ap + bs2 * i;
695:           if (roworiented) {
696:             if (is == ADD_VALUES) {
697:               for (ii = 0; ii < bs; ii++, value += stepval) {
698:                 for (jj = ii; jj < bs2; jj += bs) bap[jj] += *value++;
699:               }
700:             } else {
701:               for (ii = 0; ii < bs; ii++, value += stepval) {
702:                 for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
703:               }
704:             }
705:           } else {
706:             if (is == ADD_VALUES) {
707:               for (ii = 0; ii < bs; ii++, value += stepval) {
708:                 for (jj = 0; jj < bs; jj++) *bap++ += *value++;
709:               }
710:             } else {
711:               for (ii = 0; ii < bs; ii++, value += stepval) {
712:                 for (jj = 0; jj < bs; jj++) *bap++ = *value++;
713:               }
714:             }
715:           }
716:           goto noinsert2;
717:         }
718:       }
719:       if (nonew == 1) goto noinsert2;
720:       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new block index nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
721:       MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
722:       N = nrow++ - 1;
723:       high++;
724:       /* shift up all the later entries in this row */
725:       PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
726:       PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
727:       PetscCall(PetscArrayzero(ap + bs2 * i, bs2));
728:       rp[i] = col;
729:       bap   = ap + bs2 * i;
730:       if (roworiented) {
731:         for (ii = 0; ii < bs; ii++, value += stepval) {
732:           for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
733:         }
734:       } else {
735:         for (ii = 0; ii < bs; ii++, value += stepval) {
736:           for (jj = 0; jj < bs; jj++) *bap++ = *value++;
737:         }
738:       }
739:     noinsert2:;
740:       low = i;
741:     }
742:     ailen[row] = nrow;
743:   }
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

747: static PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A, MatAssemblyType mode)
748: {
749:   Mat_SeqSBAIJ *a      = (Mat_SeqSBAIJ *)A->data;
750:   PetscInt      fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax;
751:   PetscInt      m = A->rmap->N, *ip, N, *ailen = a->ilen;
752:   PetscInt      mbs = a->mbs, bs2 = a->bs2, rmax = 0;
753:   MatScalar    *aa = a->a, *ap;

755:   PetscFunctionBegin;
756:   if (mode == MAT_FLUSH_ASSEMBLY || (A->was_assembled && A->ass_nonzerostate == A->nonzerostate)) PetscFunctionReturn(PETSC_SUCCESS);

758:   if (m) rmax = ailen[0];
759:   for (i = 1; i < mbs; i++) {
760:     /* move each row back by the amount of empty slots (fshift) before it*/
761:     fshift += imax[i - 1] - ailen[i - 1];
762:     rmax = PetscMax(rmax, ailen[i]);
763:     if (fshift) {
764:       ip = aj + ai[i];
765:       ap = aa + bs2 * ai[i];
766:       N  = ailen[i];
767:       PetscCall(PetscArraymove(ip - fshift, ip, N));
768:       PetscCall(PetscArraymove(ap - bs2 * fshift, ap, bs2 * N));
769:     }
770:     ai[i] = ai[i - 1] + ailen[i - 1];
771:   }
772:   if (mbs) {
773:     fshift += imax[mbs - 1] - ailen[mbs - 1];
774:     ai[mbs] = ai[mbs - 1] + ailen[mbs - 1];
775:   }
776:   /* reset ilen and imax for each row */
777:   for (i = 0; i < mbs; i++) ailen[i] = imax[i] = ai[i + 1] - ai[i];
778:   a->nz = ai[mbs];

780:   /* diagonals may have moved, reset it */
781:   if (a->diag) PetscCall(PetscArraycpy(a->diag, ai, mbs));
782:   PetscCheck(!fshift || a->nounused != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unused space detected in matrix: %" PetscInt_FMT " X %" PetscInt_FMT " block size %" PetscInt_FMT ", %" PetscInt_FMT " unneeded", m, A->cmap->n, A->rmap->bs, fshift * bs2);

784:   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT ", block size %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded, %" PetscInt_FMT " used\n", m, A->rmap->N, A->rmap->bs, fshift * bs2, a->nz * bs2));
785:   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues is %" PetscInt_FMT "\n", a->reallocs));
786:   PetscCall(PetscInfo(A, "Most nonzeros blocks in any row is %" PetscInt_FMT "\n", rmax));

788:   A->info.mallocs += a->reallocs;
789:   a->reallocs         = 0;
790:   A->info.nz_unneeded = (PetscReal)fshift * bs2;
791:   a->idiagvalid       = PETSC_FALSE;
792:   a->rmax             = rmax;

794:   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
795:     if (a->jshort && a->free_jshort) {
796:       /* when matrix data structure is changed, previous jshort must be replaced */
797:       PetscCall(PetscFree(a->jshort));
798:     }
799:     PetscCall(PetscMalloc1(a->i[A->rmap->n], &a->jshort));
800:     for (i = 0; i < a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
801:     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
802:     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
803:     a->free_jshort = PETSC_TRUE;
804:   }
805:   PetscFunctionReturn(PETSC_SUCCESS);
806: }

808: /* Only add/insert a(i,j) with i<=j (blocks).
809:    Any a(i,j) with i>j input by user is ignored.
810: */

812: PetscErrorCode MatSetValues_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
813: {
814:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
815:   PetscInt     *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N, lastcol = -1;
816:   PetscInt     *imax = a->imax, *ai = a->i, *ailen = a->ilen, roworiented = a->roworiented;
817:   PetscInt     *aj = a->j, nonew = a->nonew, bs = A->rmap->bs, brow, bcol;
818:   PetscInt      ridx, cidx, bs2                 = a->bs2;
819:   MatScalar    *ap, value, *aa                  = a->a, *bap;

821:   PetscFunctionBegin;
822:   for (k = 0; k < m; k++) { /* loop over added rows */
823:     row  = im[k];           /* row number */
824:     brow = row / bs;        /* block row number */
825:     if (row < 0) continue;
826:     PetscCheck(row < A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->N - 1);
827:     rp   = aj + ai[brow];       /*ptr to beginning of column value of the row block*/
828:     ap   = aa + bs2 * ai[brow]; /*ptr to beginning of element value of the row block*/
829:     rmax = imax[brow];          /* maximum space allocated for this row */
830:     nrow = ailen[brow];         /* actual length of this row */
831:     low  = 0;
832:     high = nrow;
833:     for (l = 0; l < n; l++) { /* loop over added columns */
834:       if (in[l] < 0) continue;
835:       PetscCheck(in[l] < A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->N - 1);
836:       col  = in[l];
837:       bcol = col / bs; /* block col number */

839:       if (brow > bcol) {
840:         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
841:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
842:       }

844:       ridx = row % bs;
845:       cidx = col % bs; /*row and col index inside the block */
846:       if ((brow == bcol && ridx <= cidx) || (brow < bcol)) {
847:         /* element value a(k,l) */
848:         if (roworiented) value = v[l + k * n];
849:         else value = v[k + l * m];

851:         /* move pointer bap to a(k,l) quickly and add/insert value */
852:         if (col <= lastcol) low = 0;
853:         else high = nrow;

855:         lastcol = col;
856:         while (high - low > 7) {
857:           t = (low + high) / 2;
858:           if (rp[t] > bcol) high = t;
859:           else low = t;
860:         }
861:         for (i = low; i < high; i++) {
862:           if (rp[i] > bcol) break;
863:           if (rp[i] == bcol) {
864:             bap = ap + bs2 * i + bs * cidx + ridx;
865:             if (is == ADD_VALUES) *bap += value;
866:             else *bap = value;
867:             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
868:             if (brow == bcol && ridx < cidx) {
869:               bap = ap + bs2 * i + bs * ridx + cidx;
870:               if (is == ADD_VALUES) *bap += value;
871:               else *bap = value;
872:             }
873:             goto noinsert1;
874:           }
875:         }

877:         if (nonew == 1) goto noinsert1;
878:         PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
879:         MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);

881:         N = nrow++ - 1;
882:         high++;
883:         /* shift up all the later entries in this row */
884:         PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
885:         PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
886:         PetscCall(PetscArrayzero(ap + bs2 * i, bs2));
887:         rp[i]                          = bcol;
888:         ap[bs2 * i + bs * cidx + ridx] = value;
889:         /* for diag block, add/insert its symmetric element a(cidx,ridx) */
890:         if (brow == bcol && ridx < cidx) ap[bs2 * i + bs * ridx + cidx] = value;
891:         A->nonzerostate++;
892:       noinsert1:;
893:         low = i;
894:       }
895:     } /* end of loop over added columns */
896:     ailen[brow] = nrow;
897:   } /* end of loop over added rows */
898:   PetscFunctionReturn(PETSC_SUCCESS);
899: }

901: static PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA, IS row, const MatFactorInfo *info)
902: {
903:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)inA->data;
904:   Mat           outA;
905:   PetscBool     row_identity;

907:   PetscFunctionBegin;
908:   PetscCheck(info->levels == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only levels=0 is supported for in-place icc");
909:   PetscCall(ISIdentity(row, &row_identity));
910:   PetscCheck(row_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported");
911:   PetscCheck(inA->rmap->bs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix block size %" PetscInt_FMT " is not supported", inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */

913:   outA            = inA;
914:   inA->factortype = MAT_FACTOR_ICC;
915:   PetscCall(PetscFree(inA->solvertype));
916:   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype));

918:   PetscCall(MatMarkDiagonal_SeqSBAIJ(inA));
919:   PetscCall(MatSeqSBAIJSetNumericFactorization_inplace(inA, row_identity));

921:   PetscCall(PetscObjectReference((PetscObject)row));
922:   PetscCall(ISDestroy(&a->row));
923:   a->row = row;
924:   PetscCall(PetscObjectReference((PetscObject)row));
925:   PetscCall(ISDestroy(&a->col));
926:   a->col = row;

928:   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
929:   if (a->icol) PetscCall(ISInvertPermutation(row, PETSC_DECIDE, &a->icol));

931:   if (!a->solve_work) PetscCall(PetscMalloc1(inA->rmap->N + inA->rmap->bs, &a->solve_work));

933:   PetscCall(MatCholeskyFactorNumeric(outA, inA, info));
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: static PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat, PetscInt *indices)
938: {
939:   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
940:   PetscInt      i, nz, n;

942:   PetscFunctionBegin;
943:   nz = baij->maxnz;
944:   n  = mat->cmap->n;
945:   for (i = 0; i < nz; i++) baij->j[i] = indices[i];

947:   baij->nz = nz;
948:   for (i = 0; i < n; i++) baij->ilen[i] = baij->imax[i];

950:   PetscCall(MatSetOption(mat, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
951:   PetscFunctionReturn(PETSC_SUCCESS);
952: }

954: /*@
955:   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
956:   in a `MATSEQSBAIJ` matrix.

958:   Input Parameters:
959: + mat     - the `MATSEQSBAIJ` matrix
960: - indices - the column indices

962:   Level: advanced

964:   Notes:
965:   This can be called if you have precomputed the nonzero structure of the
966:   matrix and want to provide it to the matrix object to improve the performance
967:   of the `MatSetValues()` operation.

969:   You MUST have set the correct numbers of nonzeros per row in the call to
970:   `MatCreateSeqSBAIJ()`, and the columns indices MUST be sorted.

972:   MUST be called before any calls to `MatSetValues()`

974: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreateSeqSBAIJ`
975: @*/
976: PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat, PetscInt *indices)
977: {
978:   PetscFunctionBegin;
980:   PetscAssertPointer(indices, 2);
981:   PetscUseMethod(mat, "MatSeqSBAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices));
982:   PetscFunctionReturn(PETSC_SUCCESS);
983: }

985: static PetscErrorCode MatCopy_SeqSBAIJ(Mat A, Mat B, MatStructure str)
986: {
987:   PetscBool isbaij;

989:   PetscFunctionBegin;
990:   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &isbaij, MATSEQSBAIJ, MATMPISBAIJ, ""));
991:   PetscCheck(isbaij, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)B)->type_name);
992:   /* If the two matrices have the same copy implementation and nonzero pattern, use fast copy. */
993:   if (str == SAME_NONZERO_PATTERN && A->ops->copy == B->ops->copy) {
994:     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
995:     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)B->data;

997:     PetscCheck(a->i[a->mbs] == b->i[b->mbs], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of nonzeros in two matrices are different");
998:     PetscCheck(a->mbs == b->mbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of rows in two matrices are different");
999:     PetscCheck(a->bs2 == b->bs2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Different block size");
1000:     PetscCall(PetscArraycpy(b->a, a->a, a->bs2 * a->i[a->mbs]));
1001:     PetscCall(PetscObjectStateIncrease((PetscObject)B));
1002:   } else {
1003:     PetscCall(MatGetRowUpperTriangular(A));
1004:     PetscCall(MatCopy_Basic(A, B, str));
1005:     PetscCall(MatRestoreRowUpperTriangular(A));
1006:   }
1007:   PetscFunctionReturn(PETSC_SUCCESS);
1008: }

1010: static PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A, PetscScalar *array[])
1011: {
1012:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;

1014:   PetscFunctionBegin;
1015:   *array = a->a;
1016:   PetscFunctionReturn(PETSC_SUCCESS);
1017: }

1019: static PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A, PetscScalar *array[])
1020: {
1021:   PetscFunctionBegin;
1022:   *array = NULL;
1023:   PetscFunctionReturn(PETSC_SUCCESS);
1024: }

1026: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y, Mat X, PetscInt *nnz)
1027: {
1028:   PetscInt      bs = Y->rmap->bs, mbs = Y->rmap->N / bs;
1029:   Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data;
1030:   Mat_SeqSBAIJ *y = (Mat_SeqSBAIJ *)Y->data;

1032:   PetscFunctionBegin;
1033:   /* Set the number of nonzeros in the new matrix */
1034:   PetscCall(MatAXPYGetPreallocation_SeqX_private(mbs, x->i, x->j, y->i, y->j, nnz));
1035:   PetscFunctionReturn(PETSC_SUCCESS);
1036: }

1038: static PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
1039: {
1040:   Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data, *y = (Mat_SeqSBAIJ *)Y->data;
1041:   PetscInt      bs = Y->rmap->bs, bs2 = bs * bs;
1042:   PetscBLASInt  one = 1;

1044:   PetscFunctionBegin;
1045:   if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) {
1046:     PetscBool e = x->nz == y->nz && x->mbs == y->mbs ? PETSC_TRUE : PETSC_FALSE;
1047:     if (e) {
1048:       PetscCall(PetscArraycmp(x->i, y->i, x->mbs + 1, &e));
1049:       if (e) {
1050:         PetscCall(PetscArraycmp(x->j, y->j, x->i[x->mbs], &e));
1051:         if (e) str = SAME_NONZERO_PATTERN;
1052:       }
1053:     }
1054:     if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN");
1055:   }
1056:   if (str == SAME_NONZERO_PATTERN) {
1057:     PetscScalar  alpha = a;
1058:     PetscBLASInt bnz;
1059:     PetscCall(PetscBLASIntCast(x->nz * bs2, &bnz));
1060:     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, x->a, &one, y->a, &one));
1061:     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
1062:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1063:     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_TRUE));
1064:     PetscCall(MatAXPY_Basic(Y, a, X, str));
1065:     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_FALSE));
1066:   } else {
1067:     Mat       B;
1068:     PetscInt *nnz;
1069:     PetscCheck(bs == X->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrices must have same block size");
1070:     PetscCall(MatGetRowUpperTriangular(X));
1071:     PetscCall(MatGetRowUpperTriangular(Y));
1072:     PetscCall(PetscMalloc1(Y->rmap->N, &nnz));
1073:     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
1074:     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
1075:     PetscCall(MatSetSizes(B, Y->rmap->n, Y->cmap->n, Y->rmap->N, Y->cmap->N));
1076:     PetscCall(MatSetBlockSizesFromMats(B, Y, Y));
1077:     PetscCall(MatSetType(B, ((PetscObject)Y)->type_name));
1078:     PetscCall(MatAXPYGetPreallocation_SeqSBAIJ(Y, X, nnz));
1079:     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, nnz));

1081:     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));

1083:     PetscCall(MatHeaderMerge(Y, &B));
1084:     PetscCall(PetscFree(nnz));
1085:     PetscCall(MatRestoreRowUpperTriangular(X));
1086:     PetscCall(MatRestoreRowUpperTriangular(Y));
1087:   }
1088:   PetscFunctionReturn(PETSC_SUCCESS);
1089: }

1091: static PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A, PetscBool *flg)
1092: {
1093:   PetscFunctionBegin;
1094:   *flg = PETSC_TRUE;
1095:   PetscFunctionReturn(PETSC_SUCCESS);
1096: }

1098: static PetscErrorCode MatConjugate_SeqSBAIJ(Mat A)
1099: {
1100: #if defined(PETSC_USE_COMPLEX)
1101:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
1102:   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
1103:   MatScalar    *aa = a->a;

1105:   PetscFunctionBegin;
1106:   for (i = 0; i < nz; i++) aa[i] = PetscConj(aa[i]);
1107: #else
1108:   PetscFunctionBegin;
1109: #endif
1110:   PetscFunctionReturn(PETSC_SUCCESS);
1111: }

1113: static PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1114: {
1115:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
1116:   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
1117:   MatScalar    *aa = a->a;

1119:   PetscFunctionBegin;
1120:   for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]);
1121:   PetscFunctionReturn(PETSC_SUCCESS);
1122: }

1124: static PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1125: {
1126:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
1127:   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
1128:   MatScalar    *aa = a->a;

1130:   PetscFunctionBegin;
1131:   for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1132:   PetscFunctionReturn(PETSC_SUCCESS);
1133: }

1135: static PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A, PetscInt is_n, const PetscInt is_idx[], PetscScalar diag, Vec x, Vec b)
1136: {
1137:   Mat_SeqSBAIJ      *baij = (Mat_SeqSBAIJ *)A->data;
1138:   PetscInt           i, j, k, count;
1139:   PetscInt           bs = A->rmap->bs, bs2 = baij->bs2, row, col;
1140:   PetscScalar        zero = 0.0;
1141:   MatScalar         *aa;
1142:   const PetscScalar *xx;
1143:   PetscScalar       *bb;
1144:   PetscBool         *zeroed, vecs = PETSC_FALSE;

1146:   PetscFunctionBegin;
1147:   /* fix right-hand side if needed */
1148:   if (x && b) {
1149:     PetscCall(VecGetArrayRead(x, &xx));
1150:     PetscCall(VecGetArray(b, &bb));
1151:     vecs = PETSC_TRUE;
1152:   }

1154:   /* zero the columns */
1155:   PetscCall(PetscCalloc1(A->rmap->n, &zeroed));
1156:   for (i = 0; i < is_n; i++) {
1157:     PetscCheck(is_idx[i] >= 0 && is_idx[i] < A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", is_idx[i]);
1158:     zeroed[is_idx[i]] = PETSC_TRUE;
1159:   }
1160:   if (vecs) {
1161:     for (i = 0; i < A->rmap->N; i++) {
1162:       row = i / bs;
1163:       for (j = baij->i[row]; j < baij->i[row + 1]; j++) {
1164:         for (k = 0; k < bs; k++) {
1165:           col = bs * baij->j[j] + k;
1166:           if (col <= i) continue;
1167:           aa = ((MatScalar *)baij->a) + j * bs2 + (i % bs) + bs * k;
1168:           if (!zeroed[i] && zeroed[col]) bb[i] -= aa[0] * xx[col];
1169:           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0] * xx[i];
1170:         }
1171:       }
1172:     }
1173:     for (i = 0; i < is_n; i++) bb[is_idx[i]] = diag * xx[is_idx[i]];
1174:   }

1176:   for (i = 0; i < A->rmap->N; i++) {
1177:     if (!zeroed[i]) {
1178:       row = i / bs;
1179:       for (j = baij->i[row]; j < baij->i[row + 1]; j++) {
1180:         for (k = 0; k < bs; k++) {
1181:           col = bs * baij->j[j] + k;
1182:           if (zeroed[col]) {
1183:             aa    = ((MatScalar *)baij->a) + j * bs2 + (i % bs) + bs * k;
1184:             aa[0] = 0.0;
1185:           }
1186:         }
1187:       }
1188:     }
1189:   }
1190:   PetscCall(PetscFree(zeroed));
1191:   if (vecs) {
1192:     PetscCall(VecRestoreArrayRead(x, &xx));
1193:     PetscCall(VecRestoreArray(b, &bb));
1194:   }

1196:   /* zero the rows */
1197:   for (i = 0; i < is_n; i++) {
1198:     row   = is_idx[i];
1199:     count = (baij->i[row / bs + 1] - baij->i[row / bs]) * bs;
1200:     aa    = ((MatScalar *)baij->a) + baij->i[row / bs] * bs2 + (row % bs);
1201:     for (k = 0; k < count; k++) {
1202:       aa[0] = zero;
1203:       aa += bs;
1204:     }
1205:     if (diag != 0.0) PetscUseTypeMethod(A, setvalues, 1, &row, 1, &row, &diag, INSERT_VALUES);
1206:   }
1207:   PetscCall(MatAssemblyEnd_SeqSBAIJ(A, MAT_FINAL_ASSEMBLY));
1208:   PetscFunctionReturn(PETSC_SUCCESS);
1209: }

1211: static PetscErrorCode MatShift_SeqSBAIJ(Mat Y, PetscScalar a)
1212: {
1213:   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)Y->data;

1215:   PetscFunctionBegin;
1216:   if (!Y->preallocated || !aij->nz) PetscCall(MatSeqSBAIJSetPreallocation(Y, Y->rmap->bs, 1, NULL));
1217:   PetscCall(MatShift_Basic(Y, a));
1218:   PetscFunctionReturn(PETSC_SUCCESS);
1219: }

1221: PetscErrorCode MatEliminateZeros_SeqSBAIJ(Mat A, PetscBool keep)
1222: {
1223:   Mat_SeqSBAIJ *a      = (Mat_SeqSBAIJ *)A->data;
1224:   PetscInt      fshift = 0, fshift_prev = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax, j, k;
1225:   PetscInt      m = A->rmap->N, *ailen = a->ilen;
1226:   PetscInt      mbs = a->mbs, bs2 = a->bs2, rmax = 0;
1227:   MatScalar    *aa = a->a, *ap;
1228:   PetscBool     zero;

1230:   PetscFunctionBegin;
1231:   PetscCheck(A->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot eliminate zeros for unassembled matrix");
1232:   if (m) rmax = ailen[0];
1233:   for (i = 1; i <= mbs; i++) {
1234:     for (k = ai[i - 1]; k < ai[i]; k++) {
1235:       zero = PETSC_TRUE;
1236:       ap   = aa + bs2 * k;
1237:       for (j = 0; j < bs2 && zero; j++) {
1238:         if (ap[j] != 0.0) zero = PETSC_FALSE;
1239:       }
1240:       if (zero && (aj[k] != i - 1 || !keep)) fshift++;
1241:       else {
1242:         if (zero && aj[k] == i - 1) PetscCall(PetscInfo(A, "Keep the diagonal block at row %" PetscInt_FMT "\n", i - 1));
1243:         aj[k - fshift] = aj[k];
1244:         PetscCall(PetscArraymove(ap - bs2 * fshift, ap, bs2));
1245:       }
1246:     }
1247:     ai[i - 1] -= fshift_prev;
1248:     fshift_prev  = fshift;
1249:     ailen[i - 1] = imax[i - 1] = ai[i] - fshift - ai[i - 1];
1250:     a->nonzerorowcnt += ((ai[i] - fshift - ai[i - 1]) > 0);
1251:     rmax = PetscMax(rmax, ailen[i - 1]);
1252:   }
1253:   if (fshift) {
1254:     if (mbs) {
1255:       ai[mbs] -= fshift;
1256:       a->nz = ai[mbs];
1257:     }
1258:     PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; zeros eliminated: %" PetscInt_FMT "; nonzeros left: %" PetscInt_FMT "\n", m, A->cmap->n, fshift, a->nz));
1259:     A->nonzerostate++;
1260:     A->info.nz_unneeded += (PetscReal)fshift;
1261:     a->rmax = rmax;
1262:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1263:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1264:   }
1265:   PetscFunctionReturn(PETSC_SUCCESS);
1266: }

1268: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1269:                                        MatGetRow_SeqSBAIJ,
1270:                                        MatRestoreRow_SeqSBAIJ,
1271:                                        MatMult_SeqSBAIJ_N,
1272:                                        /*  4*/ MatMultAdd_SeqSBAIJ_N,
1273:                                        MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */
1274:                                        MatMultAdd_SeqSBAIJ_N,
1275:                                        NULL,
1276:                                        NULL,
1277:                                        NULL,
1278:                                        /* 10*/ NULL,
1279:                                        NULL,
1280:                                        MatCholeskyFactor_SeqSBAIJ,
1281:                                        MatSOR_SeqSBAIJ,
1282:                                        MatTranspose_SeqSBAIJ,
1283:                                        /* 15*/ MatGetInfo_SeqSBAIJ,
1284:                                        MatEqual_SeqSBAIJ,
1285:                                        MatGetDiagonal_SeqSBAIJ,
1286:                                        MatDiagonalScale_SeqSBAIJ,
1287:                                        MatNorm_SeqSBAIJ,
1288:                                        /* 20*/ NULL,
1289:                                        MatAssemblyEnd_SeqSBAIJ,
1290:                                        MatSetOption_SeqSBAIJ,
1291:                                        MatZeroEntries_SeqSBAIJ,
1292:                                        /* 24*/ NULL,
1293:                                        NULL,
1294:                                        NULL,
1295:                                        NULL,
1296:                                        NULL,
1297:                                        /* 29*/ MatSetUp_Seq_Hash,
1298:                                        NULL,
1299:                                        NULL,
1300:                                        NULL,
1301:                                        NULL,
1302:                                        /* 34*/ MatDuplicate_SeqSBAIJ,
1303:                                        NULL,
1304:                                        NULL,
1305:                                        NULL,
1306:                                        MatICCFactor_SeqSBAIJ,
1307:                                        /* 39*/ MatAXPY_SeqSBAIJ,
1308:                                        MatCreateSubMatrices_SeqSBAIJ,
1309:                                        MatIncreaseOverlap_SeqSBAIJ,
1310:                                        MatGetValues_SeqSBAIJ,
1311:                                        MatCopy_SeqSBAIJ,
1312:                                        /* 44*/ NULL,
1313:                                        MatScale_SeqSBAIJ,
1314:                                        MatShift_SeqSBAIJ,
1315:                                        NULL,
1316:                                        MatZeroRowsColumns_SeqSBAIJ,
1317:                                        /* 49*/ NULL,
1318:                                        MatGetRowIJ_SeqSBAIJ,
1319:                                        MatRestoreRowIJ_SeqSBAIJ,
1320:                                        NULL,
1321:                                        NULL,
1322:                                        /* 54*/ NULL,
1323:                                        NULL,
1324:                                        NULL,
1325:                                        MatPermute_SeqSBAIJ,
1326:                                        MatSetValuesBlocked_SeqSBAIJ,
1327:                                        /* 59*/ MatCreateSubMatrix_SeqSBAIJ,
1328:                                        NULL,
1329:                                        NULL,
1330:                                        NULL,
1331:                                        NULL,
1332:                                        /* 64*/ NULL,
1333:                                        NULL,
1334:                                        NULL,
1335:                                        NULL,
1336:                                        NULL,
1337:                                        /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1338:                                        NULL,
1339:                                        MatConvert_MPISBAIJ_Basic,
1340:                                        NULL,
1341:                                        NULL,
1342:                                        /* 74*/ NULL,
1343:                                        NULL,
1344:                                        NULL,
1345:                                        NULL,
1346:                                        NULL,
1347:                                        /* 79*/ NULL,
1348:                                        NULL,
1349:                                        NULL,
1350:                                        MatGetInertia_SeqSBAIJ,
1351:                                        MatLoad_SeqSBAIJ,
1352:                                        /* 84*/ NULL,
1353:                                        NULL,
1354:                                        MatIsStructurallySymmetric_SeqSBAIJ,
1355:                                        NULL,
1356:                                        NULL,
1357:                                        /* 89*/ NULL,
1358:                                        NULL,
1359:                                        NULL,
1360:                                        NULL,
1361:                                        NULL,
1362:                                        /* 94*/ NULL,
1363:                                        NULL,
1364:                                        NULL,
1365:                                        NULL,
1366:                                        NULL,
1367:                                        /* 99*/ NULL,
1368:                                        NULL,
1369:                                        NULL,
1370:                                        MatConjugate_SeqSBAIJ,
1371:                                        NULL,
1372:                                        /*104*/ NULL,
1373:                                        MatRealPart_SeqSBAIJ,
1374:                                        MatImaginaryPart_SeqSBAIJ,
1375:                                        MatGetRowUpperTriangular_SeqSBAIJ,
1376:                                        MatRestoreRowUpperTriangular_SeqSBAIJ,
1377:                                        /*109*/ NULL,
1378:                                        NULL,
1379:                                        NULL,
1380:                                        NULL,
1381:                                        MatMissingDiagonal_SeqSBAIJ,
1382:                                        /*114*/ NULL,
1383:                                        NULL,
1384:                                        NULL,
1385:                                        NULL,
1386:                                        NULL,
1387:                                        /*119*/ NULL,
1388:                                        NULL,
1389:                                        NULL,
1390:                                        NULL,
1391:                                        NULL,
1392:                                        /*124*/ NULL,
1393:                                        NULL,
1394:                                        NULL,
1395:                                        NULL,
1396:                                        NULL,
1397:                                        /*129*/ NULL,
1398:                                        NULL,
1399:                                        NULL,
1400:                                        NULL,
1401:                                        NULL,
1402:                                        /*134*/ NULL,
1403:                                        NULL,
1404:                                        NULL,
1405:                                        NULL,
1406:                                        NULL,
1407:                                        /*139*/ MatSetBlockSizes_Default,
1408:                                        NULL,
1409:                                        NULL,
1410:                                        NULL,
1411:                                        NULL,
1412:                                        /*144*/ MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ,
1413:                                        NULL,
1414:                                        NULL,
1415:                                        NULL,
1416:                                        NULL,
1417:                                        NULL,
1418:                                        /*150*/ NULL,
1419:                                        MatEliminateZeros_SeqSBAIJ,
1420:                                        NULL};

1422: static PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
1423: {
1424:   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1425:   PetscInt      nz  = aij->i[mat->rmap->N] * mat->rmap->bs * aij->bs2;

1427:   PetscFunctionBegin;
1428:   PetscCheck(aij->nonew == 1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");

1430:   /* allocate space for values if not already there */
1431:   if (!aij->saved_values) PetscCall(PetscMalloc1(nz + 1, &aij->saved_values));

1433:   /* copy values over */
1434:   PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz));
1435:   PetscFunctionReturn(PETSC_SUCCESS);
1436: }

1438: static PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
1439: {
1440:   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1441:   PetscInt      nz  = aij->i[mat->rmap->N] * mat->rmap->bs * aij->bs2;

1443:   PetscFunctionBegin;
1444:   PetscCheck(aij->nonew == 1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1445:   PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");

1447:   /* copy values over */
1448:   PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz));
1449:   PetscFunctionReturn(PETSC_SUCCESS);
1450: }

1452: static PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
1453: {
1454:   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)B->data;
1455:   PetscInt      i, mbs, nbs, bs2;
1456:   PetscBool     skipallocation = PETSC_FALSE, flg = PETSC_FALSE, realalloc = PETSC_FALSE;

1458:   PetscFunctionBegin;
1459:   if (B->hash_active) {
1460:     PetscInt bs;
1461:     B->ops[0] = b->cops;
1462:     PetscCall(PetscHMapIJVDestroy(&b->ht));
1463:     PetscCall(MatGetBlockSize(B, &bs));
1464:     if (bs > 1) PetscCall(PetscHSetIJDestroy(&b->bht));
1465:     PetscCall(PetscFree(b->dnz));
1466:     PetscCall(PetscFree(b->bdnz));
1467:     B->hash_active = PETSC_FALSE;
1468:   }
1469:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;

1471:   PetscCall(MatSetBlockSize(B, PetscAbs(bs)));
1472:   PetscCall(PetscLayoutSetUp(B->rmap));
1473:   PetscCall(PetscLayoutSetUp(B->cmap));
1474:   PetscCheck(B->rmap->N <= B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_SUP, "SEQSBAIJ matrix cannot have more rows %" PetscInt_FMT " than columns %" PetscInt_FMT, B->rmap->N, B->cmap->N);
1475:   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));

1477:   B->preallocated = PETSC_TRUE;

1479:   mbs = B->rmap->N / bs;
1480:   nbs = B->cmap->n / bs;
1481:   bs2 = bs * bs;

1483:   PetscCheck(mbs * bs == B->rmap->N && nbs * bs == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number rows, cols must be divisible by blocksize");

1485:   if (nz == MAT_SKIP_ALLOCATION) {
1486:     skipallocation = PETSC_TRUE;
1487:     nz             = 0;
1488:   }

1490:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1491:   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
1492:   if (nnz) {
1493:     for (i = 0; i < mbs; i++) {
1494:       PetscCheck(nnz[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, nnz[i]);
1495:       PetscCheck(nnz[i] <= nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be greater than block row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " block rowlength %" PetscInt_FMT, i, nnz[i], nbs);
1496:     }
1497:   }

1499:   B->ops->mult             = MatMult_SeqSBAIJ_N;
1500:   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1501:   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1502:   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;

1504:   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_no_unroll", &flg, NULL));
1505:   if (!flg) {
1506:     switch (bs) {
1507:     case 1:
1508:       B->ops->mult             = MatMult_SeqSBAIJ_1;
1509:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1510:       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1511:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1512:       break;
1513:     case 2:
1514:       B->ops->mult             = MatMult_SeqSBAIJ_2;
1515:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1516:       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1517:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1518:       break;
1519:     case 3:
1520:       B->ops->mult             = MatMult_SeqSBAIJ_3;
1521:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1522:       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1523:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1524:       break;
1525:     case 4:
1526:       B->ops->mult             = MatMult_SeqSBAIJ_4;
1527:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1528:       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1529:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1530:       break;
1531:     case 5:
1532:       B->ops->mult             = MatMult_SeqSBAIJ_5;
1533:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1534:       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1535:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1536:       break;
1537:     case 6:
1538:       B->ops->mult             = MatMult_SeqSBAIJ_6;
1539:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1540:       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1541:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1542:       break;
1543:     case 7:
1544:       B->ops->mult             = MatMult_SeqSBAIJ_7;
1545:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1546:       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1547:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1548:       break;
1549:     }
1550:   }

1552:   b->mbs = mbs;
1553:   b->nbs = nbs;
1554:   if (!skipallocation) {
1555:     if (!b->imax) {
1556:       PetscCall(PetscMalloc2(mbs, &b->imax, mbs, &b->ilen));

1558:       b->free_imax_ilen = PETSC_TRUE;
1559:     }
1560:     if (!nnz) {
1561:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1562:       else if (nz <= 0) nz = 1;
1563:       nz = PetscMin(nbs, nz);
1564:       for (i = 0; i < mbs; i++) b->imax[i] = nz;
1565:       PetscCall(PetscIntMultError(nz, mbs, &nz));
1566:     } else {
1567:       PetscInt64 nz64 = 0;
1568:       for (i = 0; i < mbs; i++) {
1569:         b->imax[i] = nnz[i];
1570:         nz64 += nnz[i];
1571:       }
1572:       PetscCall(PetscIntCast(nz64, &nz));
1573:     }
1574:     /* b->ilen will count nonzeros in each block row so far. */
1575:     for (i = 0; i < mbs; i++) b->ilen[i] = 0;
1576:     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */

1578:     /* allocate the matrix space */
1579:     PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i));
1580:     PetscCall(PetscMalloc3(bs2 * nz, &b->a, nz, &b->j, B->rmap->N + 1, &b->i));
1581:     PetscCall(PetscArrayzero(b->a, nz * bs2));
1582:     PetscCall(PetscArrayzero(b->j, nz));

1584:     b->singlemalloc = PETSC_TRUE;

1586:     /* pointer to beginning of each row */
1587:     b->i[0] = 0;
1588:     for (i = 1; i < mbs + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1];

1590:     b->free_a  = PETSC_TRUE;
1591:     b->free_ij = PETSC_TRUE;
1592:   } else {
1593:     b->free_a  = PETSC_FALSE;
1594:     b->free_ij = PETSC_FALSE;
1595:   }

1597:   b->bs2     = bs2;
1598:   b->nz      = 0;
1599:   b->maxnz   = nz;
1600:   b->inew    = NULL;
1601:   b->jnew    = NULL;
1602:   b->anew    = NULL;
1603:   b->a2anew  = NULL;
1604:   b->permute = PETSC_FALSE;

1606:   B->was_assembled = PETSC_FALSE;
1607:   B->assembled     = PETSC_FALSE;
1608:   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1609:   PetscFunctionReturn(PETSC_SUCCESS);
1610: }

1612: static PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B, PetscInt bs, const PetscInt ii[], const PetscInt jj[], const PetscScalar V[])
1613: {
1614:   PetscInt      i, j, m, nz, anz, nz_max = 0, *nnz;
1615:   PetscScalar  *values      = NULL;
1616:   Mat_SeqSBAIJ *b           = (Mat_SeqSBAIJ *)B->data;
1617:   PetscBool     roworiented = b->roworiented;
1618:   PetscBool     ilw         = b->ignore_ltriangular;

1620:   PetscFunctionBegin;
1621:   PetscCheck(bs >= 1, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_OUTOFRANGE, "Invalid block size specified, must be positive but it is %" PetscInt_FMT, bs);
1622:   PetscCall(PetscLayoutSetBlockSize(B->rmap, bs));
1623:   PetscCall(PetscLayoutSetBlockSize(B->cmap, bs));
1624:   PetscCall(PetscLayoutSetUp(B->rmap));
1625:   PetscCall(PetscLayoutSetUp(B->cmap));
1626:   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));
1627:   m = B->rmap->n / bs;

1629:   PetscCheck(!ii[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ii[0] must be 0 but it is %" PetscInt_FMT, ii[0]);
1630:   PetscCall(PetscMalloc1(m + 1, &nnz));
1631:   for (i = 0; i < m; i++) {
1632:     nz = ii[i + 1] - ii[i];
1633:     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
1634:     PetscCheckSorted(nz, jj + ii[i]);
1635:     anz = 0;
1636:     for (j = 0; j < nz; j++) {
1637:       /* count only values on the diagonal or above */
1638:       if (jj[ii[i] + j] >= i) {
1639:         anz = nz - j;
1640:         break;
1641:       }
1642:     }
1643:     nz_max = PetscMax(nz_max, nz);
1644:     nnz[i] = anz;
1645:   }
1646:   PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, nnz));
1647:   PetscCall(PetscFree(nnz));

1649:   values = (PetscScalar *)V;
1650:   if (!values) PetscCall(PetscCalloc1(bs * bs * nz_max, &values));
1651:   b->ignore_ltriangular = PETSC_TRUE;
1652:   for (i = 0; i < m; i++) {
1653:     PetscInt        ncols = ii[i + 1] - ii[i];
1654:     const PetscInt *icols = jj + ii[i];

1656:     if (!roworiented || bs == 1) {
1657:       const PetscScalar *svals = values + (V ? (bs * bs * ii[i]) : 0);
1658:       PetscCall(MatSetValuesBlocked_SeqSBAIJ(B, 1, &i, ncols, icols, svals, INSERT_VALUES));
1659:     } else {
1660:       for (j = 0; j < ncols; j++) {
1661:         const PetscScalar *svals = values + (V ? (bs * bs * (ii[i] + j)) : 0);
1662:         PetscCall(MatSetValuesBlocked_SeqSBAIJ(B, 1, &i, 1, &icols[j], svals, INSERT_VALUES));
1663:       }
1664:     }
1665:   }
1666:   if (!V) PetscCall(PetscFree(values));
1667:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1668:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1669:   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1670:   b->ignore_ltriangular = ilw;
1671:   PetscFunctionReturn(PETSC_SUCCESS);
1672: }

1674: /*
1675:    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1676: */
1677: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B, PetscBool natural)
1678: {
1679:   PetscBool flg = PETSC_FALSE;
1680:   PetscInt  bs  = B->rmap->bs;

1682:   PetscFunctionBegin;
1683:   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_no_unroll", &flg, NULL));
1684:   if (flg) bs = 8;

1686:   if (!natural) {
1687:     switch (bs) {
1688:     case 1:
1689:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1690:       break;
1691:     case 2:
1692:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1693:       break;
1694:     case 3:
1695:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1696:       break;
1697:     case 4:
1698:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1699:       break;
1700:     case 5:
1701:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1702:       break;
1703:     case 6:
1704:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1705:       break;
1706:     case 7:
1707:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1708:       break;
1709:     default:
1710:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1711:       break;
1712:     }
1713:   } else {
1714:     switch (bs) {
1715:     case 1:
1716:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1717:       break;
1718:     case 2:
1719:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1720:       break;
1721:     case 3:
1722:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1723:       break;
1724:     case 4:
1725:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1726:       break;
1727:     case 5:
1728:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1729:       break;
1730:     case 6:
1731:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1732:       break;
1733:     case 7:
1734:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1735:       break;
1736:     default:
1737:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1738:       break;
1739:     }
1740:   }
1741:   PetscFunctionReturn(PETSC_SUCCESS);
1742: }

1744: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType, MatReuse, Mat *);
1745: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
1746: static PetscErrorCode       MatFactorGetSolverType_petsc(Mat A, MatSolverType *type)
1747: {
1748:   PetscFunctionBegin;
1749:   *type = MATSOLVERPETSC;
1750:   PetscFunctionReturn(PETSC_SUCCESS);
1751: }

1753: PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A, MatFactorType ftype, Mat *B)
1754: {
1755:   PetscInt n = A->rmap->n;

1757:   PetscFunctionBegin;
1758: #if defined(PETSC_USE_COMPLEX)
1759:   if ((ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) && A->hermitian == PETSC_BOOL3_TRUE && A->symmetric != PETSC_BOOL3_TRUE) {
1760:     PetscCall(PetscInfo(A, "Hermitian MAT_FACTOR_CHOLESKY or MAT_FACTOR_ICC are not supported. Use MAT_FACTOR_LU instead.\n"));
1761:     *B = NULL;
1762:     PetscFunctionReturn(PETSC_SUCCESS);
1763:   }
1764: #endif

1766:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
1767:   PetscCall(MatSetSizes(*B, n, n, n, n));
1768:   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1769:     PetscCall(MatSetType(*B, MATSEQSBAIJ));
1770:     PetscCall(MatSeqSBAIJSetPreallocation(*B, A->rmap->bs, MAT_SKIP_ALLOCATION, NULL));

1772:     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1773:     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1774:     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
1775:     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC]));
1776:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");

1778:   (*B)->factortype     = ftype;
1779:   (*B)->canuseordering = PETSC_TRUE;
1780:   PetscCall(PetscFree((*B)->solvertype));
1781:   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
1782:   PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_petsc));
1783:   PetscFunctionReturn(PETSC_SUCCESS);
1784: }

1786: /*@C
1787:   MatSeqSBAIJGetArray - gives access to the array where the numerical data for a `MATSEQSBAIJ` matrix is stored

1789:   Not Collective

1791:   Input Parameter:
1792: . A - a `MATSEQSBAIJ` matrix

1794:   Output Parameter:
1795: . array - pointer to the data

1797:   Level: intermediate

1799: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatSeqSBAIJRestoreArray()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`
1800: @*/
1801: PetscErrorCode MatSeqSBAIJGetArray(Mat A, PetscScalar **array)
1802: {
1803:   PetscFunctionBegin;
1804:   PetscUseMethod(A, "MatSeqSBAIJGetArray_C", (Mat, PetscScalar **), (A, array));
1805:   PetscFunctionReturn(PETSC_SUCCESS);
1806: }

1808: /*@C
1809:   MatSeqSBAIJRestoreArray - returns access to the array where the numerical data for a `MATSEQSBAIJ` matrix is stored obtained by `MatSeqSBAIJGetArray()`

1811:   Not Collective

1813:   Input Parameters:
1814: + A     - a `MATSEQSBAIJ` matrix
1815: - array - pointer to the data

1817:   Level: intermediate

1819: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatSeqSBAIJGetArray()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`
1820: @*/
1821: PetscErrorCode MatSeqSBAIJRestoreArray(Mat A, PetscScalar **array)
1822: {
1823:   PetscFunctionBegin;
1824:   PetscUseMethod(A, "MatSeqSBAIJRestoreArray_C", (Mat, PetscScalar **), (A, array));
1825:   PetscFunctionReturn(PETSC_SUCCESS);
1826: }

1828: /*MC
1829:   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1830:   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.

1832:   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1833:   can call `MatSetOption`(`Mat`, `MAT_HERMITIAN`).

1835:   Options Database Key:
1836:   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to `MatSetFromOptions()`

1838:   Level: beginner

1840:   Notes:
1841:   By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1842:   stored and it is assumed they symmetric to the upper triangular). If you call `MatSetOption`(`Mat`,`MAT_IGNORE_LOWER_TRIANGULAR`,`PETSC_FALSE`) or use
1843:   the options database `-mat_ignore_lower_triangular` false it will generate an error if you try to set a value in the lower triangular portion.

1845:   The number of rows in the matrix must be less than or equal to the number of columns

1847: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreateSeqSBAIJ()`, `MatType`, `MATMPISBAIJ`
1848: M*/
1849: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1850: {
1851:   Mat_SeqSBAIJ *b;
1852:   PetscMPIInt   size;
1853:   PetscBool     no_unroll = PETSC_FALSE, no_inode = PETSC_FALSE;

1855:   PetscFunctionBegin;
1856:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
1857:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Comm must be of size 1");

1859:   PetscCall(PetscNew(&b));
1860:   B->data   = (void *)b;
1861:   B->ops[0] = MatOps_Values;

1863:   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1864:   B->ops->view       = MatView_SeqSBAIJ;
1865:   b->row             = NULL;
1866:   b->icol            = NULL;
1867:   b->reallocs        = 0;
1868:   b->saved_values    = NULL;
1869:   b->inode.limit     = 5;
1870:   b->inode.max_limit = 5;

1872:   b->roworiented        = PETSC_TRUE;
1873:   b->nonew              = 0;
1874:   b->diag               = NULL;
1875:   b->solve_work         = NULL;
1876:   b->mult_work          = NULL;
1877:   B->spptr              = NULL;
1878:   B->info.nz_unneeded   = (PetscReal)b->maxnz * b->bs2;
1879:   b->keepnonzeropattern = PETSC_FALSE;

1881:   b->inew    = NULL;
1882:   b->jnew    = NULL;
1883:   b->anew    = NULL;
1884:   b->a2anew  = NULL;
1885:   b->permute = PETSC_FALSE;

1887:   b->ignore_ltriangular = PETSC_TRUE;

1889:   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_ignore_lower_triangular", &b->ignore_ltriangular, NULL));

1891:   b->getrow_utriangular = PETSC_FALSE;

1893:   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_getrow_uppertriangular", &b->getrow_utriangular, NULL));

1895:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJGetArray_C", MatSeqSBAIJGetArray_SeqSBAIJ));
1896:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJRestoreArray_C", MatSeqSBAIJRestoreArray_SeqSBAIJ));
1897:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqSBAIJ));
1898:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqSBAIJ));
1899:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetColumnIndices_C", MatSeqSBAIJSetColumnIndices_SeqSBAIJ));
1900:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_seqaij_C", MatConvert_SeqSBAIJ_SeqAIJ));
1901:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_seqbaij_C", MatConvert_SeqSBAIJ_SeqBAIJ));
1902:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetPreallocation_C", MatSeqSBAIJSetPreallocation_SeqSBAIJ));
1903:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetPreallocationCSR_C", MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ));
1904: #if defined(PETSC_HAVE_ELEMENTAL)
1905:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_elemental_C", MatConvert_SeqSBAIJ_Elemental));
1906: #endif
1907: #if defined(PETSC_HAVE_SCALAPACK)
1908:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_scalapack_C", MatConvert_SBAIJ_ScaLAPACK));
1909: #endif

1911:   B->symmetry_eternal            = PETSC_TRUE;
1912:   B->structural_symmetry_eternal = PETSC_TRUE;
1913:   B->symmetric                   = PETSC_BOOL3_TRUE;
1914:   B->structurally_symmetric      = PETSC_BOOL3_TRUE;
1915: #if defined(PETSC_USE_COMPLEX)
1916:   B->hermitian = PETSC_BOOL3_FALSE;
1917: #else
1918:   B->hermitian = PETSC_BOOL3_TRUE;
1919: #endif

1921:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSBAIJ));

1923:   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "Options for SEQSBAIJ matrix", "Mat");
1924:   PetscCall(PetscOptionsBool("-mat_no_unroll", "Do not optimize for inodes (slower)", NULL, no_unroll, &no_unroll, NULL));
1925:   if (no_unroll) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_unroll\n"));
1926:   PetscCall(PetscOptionsBool("-mat_no_inode", "Do not optimize for inodes (slower)", NULL, no_inode, &no_inode, NULL));
1927:   if (no_inode) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_inode\n"));
1928:   PetscCall(PetscOptionsInt("-mat_inode_limit", "Do not use inodes larger then this value", NULL, b->inode.limit, &b->inode.limit, NULL));
1929:   PetscOptionsEnd();
1930:   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1931:   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1932:   PetscFunctionReturn(PETSC_SUCCESS);
1933: }

1935: /*@C
1936:   MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1937:   compressed row) `MATSEQSBAIJ` format.  For good matrix assembly performance the
1938:   user should preallocate the matrix storage by setting the parameter `nz`
1939:   (or the array `nnz`).

1941:   Collective

1943:   Input Parameters:
1944: + B   - the symmetric matrix
1945: . bs  - size of block, the blocks are ALWAYS square. One can use `MatSetBlockSizes()` to set a different row and column blocksize but the row
1946:         blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with `MatCreateVecs()`
1947: . nz  - number of block nonzeros per block row (same for all rows)
1948: - nnz - array containing the number of block nonzeros in the upper triangular plus
1949:         diagonal portion of each block (possibly different for each block row) or `NULL`

1951:   Options Database Keys:
1952: + -mat_no_unroll  - uses code that does not unroll the loops in the block calculations (much slower)
1953: - -mat_block_size - size of the blocks to use (only works if a negative bs is passed in

1955:   Level: intermediate

1957:   Notes:
1958:   Specify the preallocated storage with either `nz` or `nnz` (not both).
1959:   Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
1960:   allocation.  See [Sparse Matrices](sec_matsparse) for details.

1962:   You can call `MatGetInfo()` to get information on how effective the preallocation was;
1963:   for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1964:   You can also run with the option `-info` and look for messages with the string
1965:   malloc in them to see if additional memory allocation was needed.

1967:   If the `nnz` parameter is given then the `nz` parameter is ignored

1969: .seealso: [](ch_matrices), `Mat`, [Sparse Matrices](sec_matsparse), `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatCreateSBAIJ()`
1970: @*/
1971: PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
1972: {
1973:   PetscFunctionBegin;
1977:   PetscTryMethod(B, "MatSeqSBAIJSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[]), (B, bs, nz, nnz));
1978:   PetscFunctionReturn(PETSC_SUCCESS);
1979: }

1981: /*@C
1982:   MatSeqSBAIJSetPreallocationCSR - Creates a sparse parallel matrix in `MATSEQSBAIJ` format using the given nonzero structure and (optional) numerical values

1984:   Input Parameters:
1985: + B  - the matrix
1986: . bs - size of block, the blocks are ALWAYS square.
1987: . i  - the indices into `j` for the start of each local row (indices start with zero)
1988: . j  - the column indices for each local row (indices start with zero) these must be sorted for each row
1989: - v  - optional values in the matrix, use `NULL` if not provided

1991:   Level: advanced

1993:   Notes:
1994:   The `i`,`j`,`v` values are COPIED with this routine; to avoid the copy use `MatCreateSeqSBAIJWithArrays()`

1996:   The order of the entries in values is specified by the `MatOption` `MAT_ROW_ORIENTED`.  For example, C programs
1997:   may want to use the default `MAT_ROW_ORIENTED` = `PETSC_TRUE` and use an array v[nnz][bs][bs] where the second index is
1998:   over rows within a block and the last index is over columns within a block row.  Fortran programs will likely set
1999:   `MAT_ROW_ORIENTED` = `PETSC_FALSE` and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
2000:   block column and the second index is over columns within a block.

2002:   Any entries provided that lie below the diagonal are ignored

2004:   Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
2005:   and usually the numerical values as well

2007: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValuesBlocked()`, `MatSeqSBAIJSetPreallocation()`
2008: @*/
2009: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B, PetscInt bs, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
2010: {
2011:   PetscFunctionBegin;
2015:   PetscTryMethod(B, "MatSeqSBAIJSetPreallocationCSR_C", (Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, bs, i, j, v));
2016:   PetscFunctionReturn(PETSC_SUCCESS);
2017: }

2019: /*@C
2020:   MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in (block
2021:   compressed row) `MATSEQSBAIJ` format.  For good matrix assembly performance the
2022:   user should preallocate the matrix storage by setting the parameter `nz`
2023:   (or the array `nnz`).

2025:   Collective

2027:   Input Parameters:
2028: + comm - MPI communicator, set to `PETSC_COMM_SELF`
2029: . bs   - size of block, the blocks are ALWAYS square. One can use `MatSetBlockSizes()` to set a different row and column blocksize but the row
2030:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2031: . m    - number of rows
2032: . n    - number of columns
2033: . nz   - number of block nonzeros per block row (same for all rows)
2034: - nnz  - array containing the number of block nonzeros in the upper triangular plus
2035:          diagonal portion of each block (possibly different for each block row) or `NULL`

2037:   Output Parameter:
2038: . A - the symmetric matrix

2040:   Options Database Keys:
2041: + -mat_no_unroll  - uses code that does not unroll the loops in the block calculations (much slower)
2042: - -mat_block_size - size of the blocks to use

2044:   Level: intermediate

2046:   Notes:
2047:   It is recommended that one use `MatCreateFromOptions()` or the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
2048:   MatXXXXSetPreallocation() paradigm instead of this routine directly.
2049:   [MatXXXXSetPreallocation() is, for example, `MatSeqAIJSetPreallocation()`]

2051:   The number of rows and columns must be divisible by blocksize.
2052:   This matrix type does not support complex Hermitian operation.

2054:   Specify the preallocated storage with either `nz` or `nnz` (not both).
2055:   Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
2056:   allocation.  See [Sparse Matrices](sec_matsparse) for details.

2058:   If the `nnz` parameter is given then the `nz` parameter is ignored

2060: .seealso: [](ch_matrices), `Mat`, [Sparse Matrices](sec_matsparse), `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatCreateSBAIJ()`
2061: @*/
2062: PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
2063: {
2064:   PetscFunctionBegin;
2065:   PetscCall(MatCreate(comm, A));
2066:   PetscCall(MatSetSizes(*A, m, n, m, n));
2067:   PetscCall(MatSetType(*A, MATSEQSBAIJ));
2068:   PetscCall(MatSeqSBAIJSetPreallocation(*A, bs, nz, (PetscInt *)nnz));
2069:   PetscFunctionReturn(PETSC_SUCCESS);
2070: }

2072: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B)
2073: {
2074:   Mat           C;
2075:   Mat_SeqSBAIJ *c, *a  = (Mat_SeqSBAIJ *)A->data;
2076:   PetscInt      i, mbs = a->mbs, nz = a->nz, bs2 = a->bs2;

2078:   PetscFunctionBegin;
2079:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot duplicate unassembled matrix");
2080:   PetscCheck(a->i[mbs] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Corrupt matrix");

2082:   *B = NULL;
2083:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2084:   PetscCall(MatSetSizes(C, A->rmap->N, A->cmap->n, A->rmap->N, A->cmap->n));
2085:   PetscCall(MatSetBlockSizesFromMats(C, A, A));
2086:   PetscCall(MatSetType(C, MATSEQSBAIJ));
2087:   c = (Mat_SeqSBAIJ *)C->data;

2089:   C->preallocated       = PETSC_TRUE;
2090:   C->factortype         = A->factortype;
2091:   c->row                = NULL;
2092:   c->icol               = NULL;
2093:   c->saved_values       = NULL;
2094:   c->keepnonzeropattern = a->keepnonzeropattern;
2095:   C->assembled          = PETSC_TRUE;

2097:   PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
2098:   PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
2099:   c->bs2 = a->bs2;
2100:   c->mbs = a->mbs;
2101:   c->nbs = a->nbs;

2103:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2104:     c->imax           = a->imax;
2105:     c->ilen           = a->ilen;
2106:     c->free_imax_ilen = PETSC_FALSE;
2107:   } else {
2108:     PetscCall(PetscMalloc2((mbs + 1), &c->imax, (mbs + 1), &c->ilen));
2109:     for (i = 0; i < mbs; i++) {
2110:       c->imax[i] = a->imax[i];
2111:       c->ilen[i] = a->ilen[i];
2112:     }
2113:     c->free_imax_ilen = PETSC_TRUE;
2114:   }

2116:   /* allocate the matrix space */
2117:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2118:     PetscCall(PetscMalloc1(bs2 * nz, &c->a));
2119:     c->i            = a->i;
2120:     c->j            = a->j;
2121:     c->singlemalloc = PETSC_FALSE;
2122:     c->free_a       = PETSC_TRUE;
2123:     c->free_ij      = PETSC_FALSE;
2124:     c->parent       = A;
2125:     PetscCall(PetscObjectReference((PetscObject)A));
2126:     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
2127:     PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
2128:   } else {
2129:     PetscCall(PetscMalloc3(bs2 * nz, &c->a, nz, &c->j, mbs + 1, &c->i));
2130:     PetscCall(PetscArraycpy(c->i, a->i, mbs + 1));
2131:     c->singlemalloc = PETSC_TRUE;
2132:     c->free_a       = PETSC_TRUE;
2133:     c->free_ij      = PETSC_TRUE;
2134:   }
2135:   if (mbs > 0) {
2136:     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) PetscCall(PetscArraycpy(c->j, a->j, nz));
2137:     if (cpvalues == MAT_COPY_VALUES) {
2138:       PetscCall(PetscArraycpy(c->a, a->a, bs2 * nz));
2139:     } else {
2140:       PetscCall(PetscArrayzero(c->a, bs2 * nz));
2141:     }
2142:     if (a->jshort) {
2143:       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2144:       /* if the parent matrix is reassembled, this child matrix will never notice */
2145:       PetscCall(PetscMalloc1(nz, &c->jshort));
2146:       PetscCall(PetscArraycpy(c->jshort, a->jshort, nz));

2148:       c->free_jshort = PETSC_TRUE;
2149:     }
2150:   }

2152:   c->roworiented = a->roworiented;
2153:   c->nonew       = a->nonew;

2155:   if (a->diag) {
2156:     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2157:       c->diag      = a->diag;
2158:       c->free_diag = PETSC_FALSE;
2159:     } else {
2160:       PetscCall(PetscMalloc1(mbs, &c->diag));
2161:       for (i = 0; i < mbs; i++) c->diag[i] = a->diag[i];
2162:       c->free_diag = PETSC_TRUE;
2163:     }
2164:   }
2165:   c->nz         = a->nz;
2166:   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2167:   c->solve_work = NULL;
2168:   c->mult_work  = NULL;

2170:   *B = C;
2171:   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
2172:   PetscFunctionReturn(PETSC_SUCCESS);
2173: }

2175: /* Used for both SeqBAIJ and SeqSBAIJ matrices */
2176: #define MatLoad_SeqSBAIJ_Binary MatLoad_SeqBAIJ_Binary

2178: PetscErrorCode MatLoad_SeqSBAIJ(Mat mat, PetscViewer viewer)
2179: {
2180:   PetscBool isbinary;

2182:   PetscFunctionBegin;
2183:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2184:   PetscCheck(isbinary, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer type %s not yet supported for reading %s matrices", ((PetscObject)viewer)->type_name, ((PetscObject)mat)->type_name);
2185:   PetscCall(MatLoad_SeqSBAIJ_Binary(mat, viewer));
2186:   PetscFunctionReturn(PETSC_SUCCESS);
2187: }

2189: /*@
2190:   MatCreateSeqSBAIJWithArrays - Creates an sequential `MATSEQSBAIJ` matrix using matrix elements
2191:   (upper triangular entries in CSR format) provided by the user.

2193:   Collective

2195:   Input Parameters:
2196: + comm - must be an MPI communicator of size 1
2197: . bs   - size of block
2198: . m    - number of rows
2199: . n    - number of columns
2200: . i    - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2201: . j    - column indices
2202: - a    - matrix values

2204:   Output Parameter:
2205: . mat - the matrix

2207:   Level: advanced

2209:   Notes:
2210:   The `i`, `j`, and `a` arrays are not copied by this routine, the user must free these arrays
2211:   once the matrix is destroyed

2213:   You cannot set new nonzero locations into this matrix, that will generate an error.

2215:   The `i` and `j` indices are 0 based

2217:   When block size is greater than 1 the matrix values must be stored using the `MATSBAIJ` storage format. For block size of 1
2218:   it is the regular CSR format excluding the lower triangular elements.

2220: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSBAIJ()`, `MatCreateSeqSBAIJ()`
2221: @*/
2222: PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat)
2223: {
2224:   PetscInt      ii;
2225:   Mat_SeqSBAIJ *sbaij;

2227:   PetscFunctionBegin;
2228:   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "block size %" PetscInt_FMT " > 1 is not supported yet", bs);
2229:   PetscCheck(m == 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");

2231:   PetscCall(MatCreate(comm, mat));
2232:   PetscCall(MatSetSizes(*mat, m, n, m, n));
2233:   PetscCall(MatSetType(*mat, MATSEQSBAIJ));
2234:   PetscCall(MatSeqSBAIJSetPreallocation(*mat, bs, MAT_SKIP_ALLOCATION, NULL));
2235:   sbaij = (Mat_SeqSBAIJ *)(*mat)->data;
2236:   PetscCall(PetscMalloc2(m, &sbaij->imax, m, &sbaij->ilen));

2238:   sbaij->i = i;
2239:   sbaij->j = j;
2240:   sbaij->a = a;

2242:   sbaij->singlemalloc   = PETSC_FALSE;
2243:   sbaij->nonew          = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2244:   sbaij->free_a         = PETSC_FALSE;
2245:   sbaij->free_ij        = PETSC_FALSE;
2246:   sbaij->free_imax_ilen = PETSC_TRUE;

2248:   for (ii = 0; ii < m; ii++) {
2249:     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii + 1] - i[ii];
2250:     PetscCheck(i[ii + 1] >= i[ii], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative row length in i (row indices) row = %" PetscInt_FMT " length = %" PetscInt_FMT, ii, i[ii + 1] - i[ii]);
2251:   }
2252:   if (PetscDefined(USE_DEBUG)) {
2253:     for (ii = 0; ii < sbaij->i[m]; ii++) {
2254:       PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
2255:       PetscCheck(j[ii] < n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index too large at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
2256:     }
2257:   }

2259:   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
2260:   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
2261:   PetscFunctionReturn(PETSC_SUCCESS);
2262: }

2264: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
2265: {
2266:   PetscFunctionBegin;
2267:   PetscCall(MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm, inmat, n, scall, outmat));
2268:   PetscFunctionReturn(PETSC_SUCCESS);
2269: }