Actual source code: mpisbaij.c

  1: #include <../src/mat/impls/baij/mpi/mpibaij.h>
  2: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
  3: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  4: #include <petscblaslapack.h>
  5: #include <petscsf.h>

  7: static PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
  8: {
  9:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;

 11:   PetscFunctionBegin;
 12:   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ",Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N));
 13:   PetscCall(MatStashDestroy_Private(&mat->stash));
 14:   PetscCall(MatStashDestroy_Private(&mat->bstash));
 15:   PetscCall(MatDestroy(&baij->A));
 16:   PetscCall(MatDestroy(&baij->B));
 17: #if defined(PETSC_USE_CTABLE)
 18:   PetscCall(PetscHMapIDestroy(&baij->colmap));
 19: #else
 20:   PetscCall(PetscFree(baij->colmap));
 21: #endif
 22:   PetscCall(PetscFree(baij->garray));
 23:   PetscCall(VecDestroy(&baij->lvec));
 24:   PetscCall(VecScatterDestroy(&baij->Mvctx));
 25:   PetscCall(VecDestroy(&baij->slvec0));
 26:   PetscCall(VecDestroy(&baij->slvec0b));
 27:   PetscCall(VecDestroy(&baij->slvec1));
 28:   PetscCall(VecDestroy(&baij->slvec1a));
 29:   PetscCall(VecDestroy(&baij->slvec1b));
 30:   PetscCall(VecScatterDestroy(&baij->sMvctx));
 31:   PetscCall(PetscFree2(baij->rowvalues, baij->rowindices));
 32:   PetscCall(PetscFree(baij->barray));
 33:   PetscCall(PetscFree(baij->hd));
 34:   PetscCall(VecDestroy(&baij->diag));
 35:   PetscCall(VecDestroy(&baij->bb1));
 36:   PetscCall(VecDestroy(&baij->xx1));
 37: #if defined(PETSC_USE_REAL_MAT_SINGLE)
 38:   PetscCall(PetscFree(baij->setvaluescopy));
 39: #endif
 40:   PetscCall(PetscFree(baij->in_loc));
 41:   PetscCall(PetscFree(baij->v_loc));
 42:   PetscCall(PetscFree(baij->rangebs));
 43:   PetscCall(PetscFree(mat->data));

 45:   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
 46:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatStoreValues_C", NULL));
 47:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatRetrieveValues_C", NULL));
 48:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISBAIJSetPreallocation_C", NULL));
 49:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISBAIJSetPreallocationCSR_C", NULL));
 50: #if defined(PETSC_HAVE_ELEMENTAL)
 51:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisbaij_elemental_C", NULL));
 52: #endif
 53: #if defined(PETSC_HAVE_SCALAPACK) && (defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL_DOUBLE))
 54:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisbaij_scalapack_C", NULL));
 55: #endif
 56:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisbaij_mpiaij_C", NULL));
 57:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisbaij_mpibaij_C", NULL));
 58:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

 61: /* defines MatSetValues_MPI_Hash(), MatAssemblyBegin_MPI_Hash(), MatAssemblyEnd_MPI_Hash(), MatSetUp_MPI_Hash() */
 62: #define TYPE SBAIJ
 63: #define TYPE_SBAIJ
 64: #include "../src/mat/impls/aij/mpi/mpihashmat.h"
 65: #undef TYPE
 66: #undef TYPE_SBAIJ

 68: #if defined(PETSC_HAVE_ELEMENTAL)
 69: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
 70: #endif
 71: #if defined(PETSC_HAVE_SCALAPACK) && (defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL_DOUBLE))
 72: PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
 73: #endif

 75: /* This could be moved to matimpl.h */
 76: static PetscErrorCode MatPreallocateWithMats_Private(Mat B, PetscInt nm, Mat X[], PetscBool symm[], PetscBool fill)
 77: {
 78:   Mat       preallocator;
 79:   PetscInt  r, rstart, rend;
 80:   PetscInt  bs, i, m, n, M, N;
 81:   PetscBool cong = PETSC_TRUE;

 83:   PetscFunctionBegin;
 86:   for (i = 0; i < nm; i++) {
 88:     PetscCall(PetscLayoutCompare(B->rmap, X[i]->rmap, &cong));
 89:     PetscCheck(cong, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Not for different layouts");
 90:   }
 92:   PetscCall(MatGetBlockSize(B, &bs));
 93:   PetscCall(MatGetSize(B, &M, &N));
 94:   PetscCall(MatGetLocalSize(B, &m, &n));
 95:   PetscCall(MatCreate(PetscObjectComm((PetscObject)B), &preallocator));
 96:   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
 97:   PetscCall(MatSetBlockSize(preallocator, bs));
 98:   PetscCall(MatSetSizes(preallocator, m, n, M, N));
 99:   PetscCall(MatSetUp(preallocator));
100:   PetscCall(MatGetOwnershipRange(preallocator, &rstart, &rend));
101:   for (r = rstart; r < rend; ++r) {
102:     PetscInt           ncols;
103:     const PetscInt    *row;
104:     const PetscScalar *vals;

106:     for (i = 0; i < nm; i++) {
107:       PetscCall(MatGetRow(X[i], r, &ncols, &row, &vals));
108:       PetscCall(MatSetValues(preallocator, 1, &r, ncols, row, vals, INSERT_VALUES));
109:       if (symm && symm[i]) PetscCall(MatSetValues(preallocator, ncols, row, 1, &r, vals, INSERT_VALUES));
110:       PetscCall(MatRestoreRow(X[i], r, &ncols, &row, &vals));
111:     }
112:   }
113:   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
114:   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
115:   PetscCall(MatPreallocatorPreallocate(preallocator, fill, B));
116:   PetscCall(MatDestroy(&preallocator));
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
121: {
122:   Mat      B;
123:   PetscInt r;

125:   PetscFunctionBegin;
126:   if (reuse != MAT_REUSE_MATRIX) {
127:     PetscBool symm = PETSC_TRUE, isdense;
128:     PetscInt  bs;

130:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
131:     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
132:     PetscCall(MatSetType(B, newtype));
133:     PetscCall(MatGetBlockSize(A, &bs));
134:     PetscCall(MatSetBlockSize(B, bs));
135:     PetscCall(PetscLayoutSetUp(B->rmap));
136:     PetscCall(PetscLayoutSetUp(B->cmap));
137:     PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &isdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, ""));
138:     if (!isdense) {
139:       PetscCall(MatGetRowUpperTriangular(A));
140:       PetscCall(MatPreallocateWithMats_Private(B, 1, &A, &symm, PETSC_TRUE));
141:       PetscCall(MatRestoreRowUpperTriangular(A));
142:     } else {
143:       PetscCall(MatSetUp(B));
144:     }
145:   } else {
146:     B = *newmat;
147:     PetscCall(MatZeroEntries(B));
148:   }

150:   PetscCall(MatGetRowUpperTriangular(A));
151:   for (r = A->rmap->rstart; r < A->rmap->rend; r++) {
152:     PetscInt           ncols;
153:     const PetscInt    *row;
154:     const PetscScalar *vals;

156:     PetscCall(MatGetRow(A, r, &ncols, &row, &vals));
157:     PetscCall(MatSetValues(B, 1, &r, ncols, row, vals, INSERT_VALUES));
158: #if defined(PETSC_USE_COMPLEX)
159:     if (A->hermitian == PETSC_BOOL3_TRUE) {
160:       PetscInt i;
161:       for (i = 0; i < ncols; i++) PetscCall(MatSetValue(B, row[i], r, PetscConj(vals[i]), INSERT_VALUES));
162:     } else {
163:       PetscCall(MatSetValues(B, ncols, row, 1, &r, vals, INSERT_VALUES));
164:     }
165: #else
166:     PetscCall(MatSetValues(B, ncols, row, 1, &r, vals, INSERT_VALUES));
167: #endif
168:     PetscCall(MatRestoreRow(A, r, &ncols, &row, &vals));
169:   }
170:   PetscCall(MatRestoreRowUpperTriangular(A));
171:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
172:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));

174:   if (reuse == MAT_INPLACE_MATRIX) {
175:     PetscCall(MatHeaderReplace(A, &B));
176:   } else {
177:     *newmat = B;
178:   }
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: static PetscErrorCode MatStoreValues_MPISBAIJ(Mat mat)
183: {
184:   Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;

186:   PetscFunctionBegin;
187:   PetscCall(MatStoreValues(aij->A));
188:   PetscCall(MatStoreValues(aij->B));
189:   PetscFunctionReturn(PETSC_SUCCESS);
190: }

192: static PetscErrorCode MatRetrieveValues_MPISBAIJ(Mat mat)
193: {
194:   Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;

196:   PetscFunctionBegin;
197:   PetscCall(MatRetrieveValues(aij->A));
198:   PetscCall(MatRetrieveValues(aij->B));
199:   PetscFunctionReturn(PETSC_SUCCESS);
200: }

202: #define MatSetValues_SeqSBAIJ_A_Private(row, col, value, addv, orow, ocol) \
203:   do { \
204:     brow = row / bs; \
205:     rp   = aj + ai[brow]; \
206:     ap   = aa + bs2 * ai[brow]; \
207:     rmax = aimax[brow]; \
208:     nrow = ailen[brow]; \
209:     bcol = col / bs; \
210:     ridx = row % bs; \
211:     cidx = col % bs; \
212:     low  = 0; \
213:     high = nrow; \
214:     while (high - low > 3) { \
215:       t = (low + high) / 2; \
216:       if (rp[t] > bcol) high = t; \
217:       else low = t; \
218:     } \
219:     for (_i = low; _i < high; _i++) { \
220:       if (rp[_i] > bcol) break; \
221:       if (rp[_i] == bcol) { \
222:         bap = ap + bs2 * _i + bs * cidx + ridx; \
223:         if (addv == ADD_VALUES) *bap += value; \
224:         else *bap = value; \
225:         goto a_noinsert; \
226:       } \
227:     } \
228:     if (a->nonew == 1) goto a_noinsert; \
229:     PetscCheck(a->nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
230:     MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, aimax, a->nonew, MatScalar); \
231:     N = nrow++ - 1; \
232:     /* shift up all the later entries in this row */ \
233:     PetscCall(PetscArraymove(rp + _i + 1, rp + _i, N - _i + 1)); \
234:     PetscCall(PetscArraymove(ap + bs2 * (_i + 1), ap + bs2 * _i, bs2 * (N - _i + 1))); \
235:     PetscCall(PetscArrayzero(ap + bs2 * _i, bs2)); \
236:     rp[_i]                          = bcol; \
237:     ap[bs2 * _i + bs * cidx + ridx] = value; \
238:   a_noinsert:; \
239:     ailen[brow] = nrow; \
240:   } while (0)

242: #define MatSetValues_SeqSBAIJ_B_Private(row, col, value, addv, orow, ocol) \
243:   do { \
244:     brow = row / bs; \
245:     rp   = bj + bi[brow]; \
246:     ap   = ba + bs2 * bi[brow]; \
247:     rmax = bimax[brow]; \
248:     nrow = bilen[brow]; \
249:     bcol = col / bs; \
250:     ridx = row % bs; \
251:     cidx = col % bs; \
252:     low  = 0; \
253:     high = nrow; \
254:     while (high - low > 3) { \
255:       t = (low + high) / 2; \
256:       if (rp[t] > bcol) high = t; \
257:       else low = t; \
258:     } \
259:     for (_i = low; _i < high; _i++) { \
260:       if (rp[_i] > bcol) break; \
261:       if (rp[_i] == bcol) { \
262:         bap = ap + bs2 * _i + bs * cidx + ridx; \
263:         if (addv == ADD_VALUES) *bap += value; \
264:         else *bap = value; \
265:         goto b_noinsert; \
266:       } \
267:     } \
268:     if (b->nonew == 1) goto b_noinsert; \
269:     PetscCheck(b->nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
270:     MatSeqXAIJReallocateAIJ(B, b->mbs, bs2, nrow, brow, bcol, rmax, ba, bi, bj, rp, ap, bimax, b->nonew, MatScalar); \
271:     N = nrow++ - 1; \
272:     /* shift up all the later entries in this row */ \
273:     PetscCall(PetscArraymove(rp + _i + 1, rp + _i, N - _i + 1)); \
274:     PetscCall(PetscArraymove(ap + bs2 * (_i + 1), ap + bs2 * _i, bs2 * (N - _i + 1))); \
275:     PetscCall(PetscArrayzero(ap + bs2 * _i, bs2)); \
276:     rp[_i]                          = bcol; \
277:     ap[bs2 * _i + bs * cidx + ridx] = value; \
278:   b_noinsert:; \
279:     bilen[brow] = nrow; \
280:   } while (0)

282: /* Only add/insert a(i,j) with i<=j (blocks).
283:    Any a(i,j) with i>j input by user is ignored or generates an error
284: */
285: static PetscErrorCode MatSetValues_MPISBAIJ(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode addv)
286: {
287:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
288:   MatScalar     value;
289:   PetscBool     roworiented = baij->roworiented;
290:   PetscInt      i, j, row, col;
291:   PetscInt      rstart_orig = mat->rmap->rstart;
292:   PetscInt      rend_orig = mat->rmap->rend, cstart_orig = mat->cmap->rstart;
293:   PetscInt      cend_orig = mat->cmap->rend, bs = mat->rmap->bs;

295:   /* Some Variables required in the macro */
296:   Mat           A     = baij->A;
297:   Mat_SeqSBAIJ *a     = (Mat_SeqSBAIJ *)A->data;
298:   PetscInt     *aimax = a->imax, *ai = a->i, *ailen = a->ilen, *aj = a->j;
299:   MatScalar    *aa = a->a;

301:   Mat          B     = baij->B;
302:   Mat_SeqBAIJ *b     = (Mat_SeqBAIJ *)B->data;
303:   PetscInt    *bimax = b->imax, *bi = b->i, *bilen = b->ilen, *bj = b->j;
304:   MatScalar   *ba = b->a;

306:   PetscInt  *rp, ii, nrow, _i, rmax, N, brow, bcol;
307:   PetscInt   low, high, t, ridx, cidx, bs2 = a->bs2;
308:   MatScalar *ap, *bap;

310:   /* for stash */
311:   PetscInt   n_loc, *in_loc = NULL;
312:   MatScalar *v_loc = NULL;

314:   PetscFunctionBegin;
315:   if (!baij->donotstash) {
316:     if (n > baij->n_loc) {
317:       PetscCall(PetscFree(baij->in_loc));
318:       PetscCall(PetscFree(baij->v_loc));
319:       PetscCall(PetscMalloc1(n, &baij->in_loc));
320:       PetscCall(PetscMalloc1(n, &baij->v_loc));

322:       baij->n_loc = n;
323:     }
324:     in_loc = baij->in_loc;
325:     v_loc  = baij->v_loc;
326:   }

328:   for (i = 0; i < m; i++) {
329:     if (im[i] < 0) continue;
330:     PetscCheck(im[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, im[i], mat->rmap->N - 1);
331:     if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
332:       row = im[i] - rstart_orig;                     /* local row index */
333:       for (j = 0; j < n; j++) {
334:         if (im[i] / bs > in[j] / bs) {
335:           PetscCheck(a->ignore_ltriangular, 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)");
336:           continue; /* ignore lower triangular blocks */
337:         }
338:         if (in[j] >= cstart_orig && in[j] < cend_orig) { /* diag entry (A) */
339:           col  = in[j] - cstart_orig;                    /* local col index */
340:           brow = row / bs;
341:           bcol = col / bs;
342:           if (brow > bcol) continue; /* ignore lower triangular blocks of A */
343:           if (roworiented) value = v[i * n + j];
344:           else value = v[i + j * m];
345:           MatSetValues_SeqSBAIJ_A_Private(row, col, value, addv, im[i], in[j]);
346:           /* PetscCall(MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv)); */
347:         } else if (in[j] < 0) {
348:           continue;
349:         } else {
350:           PetscCheck(in[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[j], mat->cmap->N - 1);
351:           /* off-diag entry (B) */
352:           if (mat->was_assembled) {
353:             if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
354: #if defined(PETSC_USE_CTABLE)
355:             PetscCall(PetscHMapIGetWithDefault(baij->colmap, in[j] / bs + 1, 0, &col));
356:             col = col - 1;
357: #else
358:             col = baij->colmap[in[j] / bs] - 1;
359: #endif
360:             if (col < 0 && !((Mat_SeqSBAIJ *)baij->A->data)->nonew) {
361:               PetscCall(MatDisAssemble_MPISBAIJ(mat));
362:               col = in[j];
363:               /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
364:               B     = baij->B;
365:               b     = (Mat_SeqBAIJ *)B->data;
366:               bimax = b->imax;
367:               bi    = b->i;
368:               bilen = b->ilen;
369:               bj    = b->j;
370:               ba    = b->a;
371:             } else col += in[j] % bs;
372:           } else col = in[j];
373:           if (roworiented) value = v[i * n + j];
374:           else value = v[i + j * m];
375:           MatSetValues_SeqSBAIJ_B_Private(row, col, value, addv, im[i], in[j]);
376:           /* PetscCall(MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv)); */
377:         }
378:       }
379:     } else { /* off processor entry */
380:       PetscCheck(!mat->nooffprocentries, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Setting off process row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set", im[i]);
381:       if (!baij->donotstash) {
382:         mat->assembled = PETSC_FALSE;
383:         n_loc          = 0;
384:         for (j = 0; j < n; j++) {
385:           if (im[i] / bs > in[j] / bs) continue; /* ignore lower triangular blocks */
386:           in_loc[n_loc] = in[j];
387:           if (roworiented) {
388:             v_loc[n_loc] = v[i * n + j];
389:           } else {
390:             v_loc[n_loc] = v[j * m + i];
391:           }
392:           n_loc++;
393:         }
394:         PetscCall(MatStashValuesRow_Private(&mat->stash, im[i], n_loc, in_loc, v_loc, PETSC_FALSE));
395:       }
396:     }
397:   }
398:   PetscFunctionReturn(PETSC_SUCCESS);
399: }

401: static inline PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_Inlined(Mat A, PetscInt row, PetscInt col, const PetscScalar v[], InsertMode is, PetscInt orow, PetscInt ocol)
402: {
403:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ *)A->data;
404:   PetscInt          *rp, low, high, t, ii, jj, nrow, i, rmax, N;
405:   PetscInt          *imax = a->imax, *ai = a->i, *ailen = a->ilen;
406:   PetscInt          *aj = a->j, nonew = a->nonew, bs2 = a->bs2, bs = A->rmap->bs;
407:   PetscBool          roworiented = a->roworiented;
408:   const PetscScalar *value       = v;
409:   MatScalar         *ap, *aa = a->a, *bap;

411:   PetscFunctionBegin;
412:   if (col < row) {
413:     PetscCheck(a->ignore_ltriangular, 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)");
414:     PetscFunctionReturn(PETSC_SUCCESS); /* ignore lower triangular block */
415:   }
416:   rp    = aj + ai[row];
417:   ap    = aa + bs2 * ai[row];
418:   rmax  = imax[row];
419:   nrow  = ailen[row];
420:   value = v;
421:   low   = 0;
422:   high  = nrow;

424:   while (high - low > 7) {
425:     t = (low + high) / 2;
426:     if (rp[t] > col) high = t;
427:     else low = t;
428:   }
429:   for (i = low; i < high; i++) {
430:     if (rp[i] > col) break;
431:     if (rp[i] == col) {
432:       bap = ap + bs2 * i;
433:       if (roworiented) {
434:         if (is == ADD_VALUES) {
435:           for (ii = 0; ii < bs; ii++) {
436:             for (jj = ii; jj < bs2; jj += bs) bap[jj] += *value++;
437:           }
438:         } else {
439:           for (ii = 0; ii < bs; ii++) {
440:             for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
441:           }
442:         }
443:       } else {
444:         if (is == ADD_VALUES) {
445:           for (ii = 0; ii < bs; ii++) {
446:             for (jj = 0; jj < bs; jj++) *bap++ += *value++;
447:           }
448:         } else {
449:           for (ii = 0; ii < bs; ii++) {
450:             for (jj = 0; jj < bs; jj++) *bap++ = *value++;
451:           }
452:         }
453:       }
454:       goto noinsert2;
455:     }
456:   }
457:   if (nonew == 1) goto noinsert2;
458:   PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new block index nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", orow, ocol);
459:   MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
460:   N = nrow++ - 1;
461:   high++;
462:   /* shift up all the later entries in this row */
463:   PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
464:   PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
465:   rp[i] = col;
466:   bap   = ap + bs2 * i;
467:   if (roworiented) {
468:     for (ii = 0; ii < bs; ii++) {
469:       for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
470:     }
471:   } else {
472:     for (ii = 0; ii < bs; ii++) {
473:       for (jj = 0; jj < bs; jj++) *bap++ = *value++;
474:     }
475:   }
476: noinsert2:;
477:   ailen[row] = nrow;
478:   PetscFunctionReturn(PETSC_SUCCESS);
479: }

481: /*
482:    This routine is exactly duplicated in mpibaij.c
483: */
484: static inline PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A, PetscInt row, PetscInt col, const PetscScalar v[], InsertMode is, PetscInt orow, PetscInt ocol)
485: {
486:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
487:   PetscInt          *rp, low, high, t, ii, jj, nrow, i, rmax, N;
488:   PetscInt          *imax = a->imax, *ai = a->i, *ailen = a->ilen;
489:   PetscInt          *aj = a->j, nonew = a->nonew, bs2 = a->bs2, bs = A->rmap->bs;
490:   PetscBool          roworiented = a->roworiented;
491:   const PetscScalar *value       = v;
492:   MatScalar         *ap, *aa = a->a, *bap;

494:   PetscFunctionBegin;
495:   rp    = aj + ai[row];
496:   ap    = aa + bs2 * ai[row];
497:   rmax  = imax[row];
498:   nrow  = ailen[row];
499:   low   = 0;
500:   high  = nrow;
501:   value = v;
502:   while (high - low > 7) {
503:     t = (low + high) / 2;
504:     if (rp[t] > col) high = t;
505:     else low = t;
506:   }
507:   for (i = low; i < high; i++) {
508:     if (rp[i] > col) break;
509:     if (rp[i] == col) {
510:       bap = ap + bs2 * i;
511:       if (roworiented) {
512:         if (is == ADD_VALUES) {
513:           for (ii = 0; ii < bs; ii++) {
514:             for (jj = ii; jj < bs2; jj += bs) bap[jj] += *value++;
515:           }
516:         } else {
517:           for (ii = 0; ii < bs; ii++) {
518:             for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
519:           }
520:         }
521:       } else {
522:         if (is == ADD_VALUES) {
523:           for (ii = 0; ii < bs; ii++, value += bs) {
524:             for (jj = 0; jj < bs; jj++) bap[jj] += value[jj];
525:             bap += bs;
526:           }
527:         } else {
528:           for (ii = 0; ii < bs; ii++, value += bs) {
529:             for (jj = 0; jj < bs; jj++) bap[jj] = value[jj];
530:             bap += bs;
531:           }
532:         }
533:       }
534:       goto noinsert2;
535:     }
536:   }
537:   if (nonew == 1) goto noinsert2;
538:   PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new global block indexed nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", orow, ocol);
539:   MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
540:   N = nrow++ - 1;
541:   high++;
542:   /* shift up all the later entries in this row */
543:   PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
544:   PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
545:   rp[i] = col;
546:   bap   = ap + bs2 * i;
547:   if (roworiented) {
548:     for (ii = 0; ii < bs; ii++) {
549:       for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
550:     }
551:   } else {
552:     for (ii = 0; ii < bs; ii++) {
553:       for (jj = 0; jj < bs; jj++) *bap++ = *value++;
554:     }
555:   }
556: noinsert2:;
557:   ailen[row] = nrow;
558:   PetscFunctionReturn(PETSC_SUCCESS);
559: }

561: /*
562:     This routine could be optimized by removing the need for the block copy below and passing stride information
563:   to the above inline routines; similarly in MatSetValuesBlocked_MPIBAIJ()
564: */
565: static PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const MatScalar v[], InsertMode addv)
566: {
567:   Mat_MPISBAIJ    *baij = (Mat_MPISBAIJ *)mat->data;
568:   const MatScalar *value;
569:   MatScalar       *barray      = baij->barray;
570:   PetscBool        roworiented = baij->roworiented, ignore_ltriangular = ((Mat_SeqSBAIJ *)baij->A->data)->ignore_ltriangular;
571:   PetscInt         i, j, ii, jj, row, col, rstart = baij->rstartbs;
572:   PetscInt         rend = baij->rendbs, cstart = baij->cstartbs, stepval;
573:   PetscInt         cend = baij->cendbs, bs = mat->rmap->bs, bs2 = baij->bs2;

575:   PetscFunctionBegin;
576:   if (!barray) {
577:     PetscCall(PetscMalloc1(bs2, &barray));
578:     baij->barray = barray;
579:   }

581:   if (roworiented) {
582:     stepval = (n - 1) * bs;
583:   } else {
584:     stepval = (m - 1) * bs;
585:   }
586:   for (i = 0; i < m; i++) {
587:     if (im[i] < 0) continue;
588:     PetscCheck(im[i] < baij->Mbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block indexed row too large %" PetscInt_FMT " max %" PetscInt_FMT, im[i], baij->Mbs - 1);
589:     if (im[i] >= rstart && im[i] < rend) {
590:       row = im[i] - rstart;
591:       for (j = 0; j < n; j++) {
592:         if (im[i] > in[j]) {
593:           PetscCheck(ignore_ltriangular, 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)");
594:           continue; /* ignore lower triangular blocks */
595:         }
596:         /* If NumCol = 1 then a copy is not required */
597:         if (roworiented && n == 1) {
598:           barray = (MatScalar *)v + i * bs2;
599:         } else if ((!roworiented) && (m == 1)) {
600:           barray = (MatScalar *)v + j * bs2;
601:         } else { /* Here a copy is required */
602:           if (roworiented) {
603:             value = v + i * (stepval + bs) * bs + j * bs;
604:           } else {
605:             value = v + j * (stepval + bs) * bs + i * bs;
606:           }
607:           for (ii = 0; ii < bs; ii++, value += stepval) {
608:             for (jj = 0; jj < bs; jj++) *barray++ = *value++;
609:           }
610:           barray -= bs2;
611:         }

613:         if (in[j] >= cstart && in[j] < cend) {
614:           col = in[j] - cstart;
615:           PetscCall(MatSetValuesBlocked_SeqSBAIJ_Inlined(baij->A, row, col, barray, addv, im[i], in[j]));
616:         } else if (in[j] < 0) {
617:           continue;
618:         } else {
619:           PetscCheck(in[j] < baij->Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block indexed column too large %" PetscInt_FMT " max %" PetscInt_FMT, in[j], baij->Nbs - 1);
620:           if (mat->was_assembled) {
621:             if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));

623: #if defined(PETSC_USE_CTABLE)
624:             PetscCall(PetscHMapIGetWithDefault(baij->colmap, in[j] + 1, 0, &col));
625:             col = col < 1 ? -1 : (col - 1) / bs;
626: #else
627:             col = baij->colmap[in[j]] < 1 ? -1 : (baij->colmap[in[j]] - 1) / bs;
628: #endif
629:             if (col < 0 && !((Mat_SeqBAIJ *)baij->A->data)->nonew) {
630:               PetscCall(MatDisAssemble_MPISBAIJ(mat));
631:               col = in[j];
632:             }
633:           } else col = in[j];
634:           PetscCall(MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B, row, col, barray, addv, im[i], in[j]));
635:         }
636:       }
637:     } else {
638:       PetscCheck(!mat->nooffprocentries, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Setting off process block indexed row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set", im[i]);
639:       if (!baij->donotstash) {
640:         if (roworiented) {
641:           PetscCall(MatStashValuesRowBlocked_Private(&mat->bstash, im[i], n, in, v, m, n, i));
642:         } else {
643:           PetscCall(MatStashValuesColBlocked_Private(&mat->bstash, im[i], n, in, v, m, n, i));
644:         }
645:       }
646:     }
647:   }
648:   PetscFunctionReturn(PETSC_SUCCESS);
649: }

651: static PetscErrorCode MatGetValues_MPISBAIJ(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
652: {
653:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
654:   PetscInt      bs = mat->rmap->bs, i, j, bsrstart = mat->rmap->rstart, bsrend = mat->rmap->rend;
655:   PetscInt      bscstart = mat->cmap->rstart, bscend = mat->cmap->rend, row, col, data;

657:   PetscFunctionBegin;
658:   for (i = 0; i < m; i++) {
659:     if (idxm[i] < 0) continue; /* negative row */
660:     PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm[i], mat->rmap->N - 1);
661:     PetscCheck(idxm[i] >= bsrstart && idxm[i] < bsrend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
662:     row = idxm[i] - bsrstart;
663:     for (j = 0; j < n; j++) {
664:       if (idxn[j] < 0) continue; /* negative column */
665:       PetscCheck(idxn[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, idxn[j], mat->cmap->N - 1);
666:       if (idxn[j] >= bscstart && idxn[j] < bscend) {
667:         col = idxn[j] - bscstart;
668:         PetscCall(MatGetValues_SeqSBAIJ(baij->A, 1, &row, 1, &col, v + i * n + j));
669:       } else {
670:         if (!baij->colmap) PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
671: #if defined(PETSC_USE_CTABLE)
672:         PetscCall(PetscHMapIGetWithDefault(baij->colmap, idxn[j] / bs + 1, 0, &data));
673:         data--;
674: #else
675:         data = baij->colmap[idxn[j] / bs] - 1;
676: #endif
677:         if (data < 0 || baij->garray[data / bs] != idxn[j] / bs) *(v + i * n + j) = 0.0;
678:         else {
679:           col = data + idxn[j] % bs;
680:           PetscCall(MatGetValues_SeqBAIJ(baij->B, 1, &row, 1, &col, v + i * n + j));
681:         }
682:       }
683:     }
684:   }
685:   PetscFunctionReturn(PETSC_SUCCESS);
686: }

688: static PetscErrorCode MatNorm_MPISBAIJ(Mat mat, NormType type, PetscReal *norm)
689: {
690:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
691:   PetscReal     sum[2], *lnorm2;

693:   PetscFunctionBegin;
694:   if (baij->size == 1) {
695:     PetscCall(MatNorm(baij->A, type, norm));
696:   } else {
697:     if (type == NORM_FROBENIUS) {
698:       PetscCall(PetscMalloc1(2, &lnorm2));
699:       PetscCall(MatNorm(baij->A, type, lnorm2));
700:       *lnorm2 = (*lnorm2) * (*lnorm2);
701:       lnorm2++; /* square power of norm(A) */
702:       PetscCall(MatNorm(baij->B, type, lnorm2));
703:       *lnorm2 = (*lnorm2) * (*lnorm2);
704:       lnorm2--; /* square power of norm(B) */
705:       PetscCallMPI(MPIU_Allreduce(lnorm2, sum, 2, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)mat)));
706:       *norm = PetscSqrtReal(sum[0] + 2 * sum[1]);
707:       PetscCall(PetscFree(lnorm2));
708:     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
709:       Mat_SeqSBAIJ *amat = (Mat_SeqSBAIJ *)baij->A->data;
710:       Mat_SeqBAIJ  *bmat = (Mat_SeqBAIJ *)baij->B->data;
711:       PetscReal    *rsum, vabs;
712:       PetscInt     *jj, *garray = baij->garray, rstart = baij->rstartbs, nz;
713:       PetscInt      brow, bcol, col, bs = baij->A->rmap->bs, row, grow, gcol, mbs = amat->mbs;
714:       MatScalar    *v;

716:       PetscCall(PetscCalloc1(mat->cmap->N, &rsum));
717:       /* Amat */
718:       v  = amat->a;
719:       jj = amat->j;
720:       for (brow = 0; brow < mbs; brow++) {
721:         grow = bs * (rstart + brow);
722:         nz   = amat->i[brow + 1] - amat->i[brow];
723:         for (bcol = 0; bcol < nz; bcol++) {
724:           gcol = bs * (rstart + *jj);
725:           jj++;
726:           for (col = 0; col < bs; col++) {
727:             for (row = 0; row < bs; row++) {
728:               vabs = PetscAbsScalar(*v);
729:               v++;
730:               rsum[gcol + col] += vabs;
731:               /* non-diagonal block */
732:               if (bcol > 0 && vabs > 0.0) rsum[grow + row] += vabs;
733:             }
734:           }
735:         }
736:         PetscCall(PetscLogFlops(nz * bs * bs));
737:       }
738:       /* Bmat */
739:       v  = bmat->a;
740:       jj = bmat->j;
741:       for (brow = 0; brow < mbs; brow++) {
742:         grow = bs * (rstart + brow);
743:         nz   = bmat->i[brow + 1] - bmat->i[brow];
744:         for (bcol = 0; bcol < nz; bcol++) {
745:           gcol = bs * garray[*jj];
746:           jj++;
747:           for (col = 0; col < bs; col++) {
748:             for (row = 0; row < bs; row++) {
749:               vabs = PetscAbsScalar(*v);
750:               v++;
751:               rsum[gcol + col] += vabs;
752:               rsum[grow + row] += vabs;
753:             }
754:           }
755:         }
756:         PetscCall(PetscLogFlops(nz * bs * bs));
757:       }
758:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, rsum, mat->cmap->N, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)mat)));
759:       *norm = 0.0;
760:       for (col = 0; col < mat->cmap->N; col++) {
761:         if (rsum[col] > *norm) *norm = rsum[col];
762:       }
763:       PetscCall(PetscFree(rsum));
764:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for this norm yet");
765:   }
766:   PetscFunctionReturn(PETSC_SUCCESS);
767: }

769: static PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat, MatAssemblyType mode)
770: {
771:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
772:   PetscInt      nstash, reallocs;

774:   PetscFunctionBegin;
775:   if (baij->donotstash || mat->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);

777:   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
778:   PetscCall(MatStashScatterBegin_Private(mat, &mat->bstash, baij->rangebs));
779:   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
780:   PetscCall(PetscInfo(mat, "Stash has %" PetscInt_FMT " entries,uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
781:   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
782:   PetscCall(PetscInfo(mat, "Block-Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
783:   PetscFunctionReturn(PETSC_SUCCESS);
784: }

786: static PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat, MatAssemblyType mode)
787: {
788:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
789:   Mat_SeqSBAIJ *a    = (Mat_SeqSBAIJ *)baij->A->data;
790:   PetscInt      i, j, rstart, ncols, flg, bs2 = baij->bs2;
791:   PetscInt     *row, *col;
792:   PetscBool     all_assembled;
793:   PetscMPIInt   n;
794:   PetscBool     r1, r2, r3;
795:   MatScalar    *val;

797:   /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
798:   PetscFunctionBegin;
799:   if (!baij->donotstash && !mat->nooffprocentries) {
800:     while (1) {
801:       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
802:       if (!flg) break;

804:       for (i = 0; i < n;) {
805:         /* Now identify the consecutive vals belonging to the same row */
806:         for (j = i, rstart = row[j]; j < n; j++) {
807:           if (row[j] != rstart) break;
808:         }
809:         if (j < n) ncols = j - i;
810:         else ncols = n - i;
811:         /* Now assemble all these values with a single function call */
812:         PetscCall(MatSetValues_MPISBAIJ(mat, 1, row + i, ncols, col + i, val + i, mat->insertmode));
813:         i = j;
814:       }
815:     }
816:     PetscCall(MatStashScatterEnd_Private(&mat->stash));
817:     /* Now process the block-stash. Since the values are stashed column-oriented,
818:        set the row-oriented flag to column-oriented, and after MatSetValues()
819:        restore the original flags */
820:     r1 = baij->roworiented;
821:     r2 = a->roworiented;
822:     r3 = ((Mat_SeqBAIJ *)baij->B->data)->roworiented;

824:     baij->roworiented = PETSC_FALSE;
825:     a->roworiented    = PETSC_FALSE;

827:     ((Mat_SeqBAIJ *)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworiented */
828:     while (1) {
829:       PetscCall(MatStashScatterGetMesg_Private(&mat->bstash, &n, &row, &col, &val, &flg));
830:       if (!flg) break;

832:       for (i = 0; i < n;) {
833:         /* Now identify the consecutive vals belonging to the same row */
834:         for (j = i, rstart = row[j]; j < n; j++) {
835:           if (row[j] != rstart) break;
836:         }
837:         if (j < n) ncols = j - i;
838:         else ncols = n - i;
839:         PetscCall(MatSetValuesBlocked_MPISBAIJ(mat, 1, row + i, ncols, col + i, val + i * bs2, mat->insertmode));
840:         i = j;
841:       }
842:     }
843:     PetscCall(MatStashScatterEnd_Private(&mat->bstash));

845:     baij->roworiented = r1;
846:     a->roworiented    = r2;

848:     ((Mat_SeqBAIJ *)baij->B->data)->roworiented = r3; /* b->roworiented */
849:   }

851:   PetscCall(MatAssemblyBegin(baij->A, mode));
852:   PetscCall(MatAssemblyEnd(baij->A, mode));

854:   /* determine if any process has disassembled, if so we must
855:      also disassemble ourselves, in order that we may reassemble. */
856:   /*
857:      if nonzero structure of submatrix B cannot change then we know that
858:      no process disassembled thus we can skip this stuff
859:   */
860:   if (!((Mat_SeqBAIJ *)baij->B->data)->nonew) {
861:     PetscCallMPI(MPIU_Allreduce(&mat->was_assembled, &all_assembled, 1, MPI_C_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
862:     if (mat->was_assembled && !all_assembled) PetscCall(MatDisAssemble_MPISBAIJ(mat));
863:   }

865:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) PetscCall(MatSetUpMultiply_MPISBAIJ(mat)); /* setup Mvctx and sMvctx */
866:   PetscCall(MatAssemblyBegin(baij->B, mode));
867:   PetscCall(MatAssemblyEnd(baij->B, mode));

869:   PetscCall(PetscFree2(baij->rowvalues, baij->rowindices));

871:   baij->rowvalues = NULL;

873:   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
874:   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ *)baij->A->data)->nonew) {
875:     PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
876:     PetscCallMPI(MPIU_Allreduce(&state, &mat->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)mat)));
877:   }
878:   PetscFunctionReturn(PETSC_SUCCESS);
879: }

881: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
882: #include <petscdraw.h>
883: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat, PetscViewer viewer)
884: {
885:   Mat_MPISBAIJ     *baij = (Mat_MPISBAIJ *)mat->data;
886:   PetscInt          bs   = mat->rmap->bs;
887:   PetscMPIInt       rank = baij->rank;
888:   PetscBool         isascii, isdraw;
889:   PetscViewer       sviewer;
890:   PetscViewerFormat format;

892:   PetscFunctionBegin;
893:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
894:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
895:   if (isascii) {
896:     PetscCall(PetscViewerGetFormat(viewer, &format));
897:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
898:       MatInfo info;
899:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
900:       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
901:       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
902:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " bs %" PetscInt_FMT " mem %g\n", rank, mat->rmap->n, (PetscInt)info.nz_used, (PetscInt)info.nz_allocated,
903:                                                    mat->rmap->bs, info.memory));
904:       PetscCall(MatGetInfo(baij->A, MAT_LOCAL, &info));
905:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] on-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
906:       PetscCall(MatGetInfo(baij->B, MAT_LOCAL, &info));
907:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] off-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
908:       PetscCall(PetscViewerFlush(viewer));
909:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
910:       PetscCall(PetscViewerASCIIPrintf(viewer, "Information on VecScatter used in matrix-vector product: \n"));
911:       PetscCall(VecScatterView(baij->Mvctx, viewer));
912:       PetscFunctionReturn(PETSC_SUCCESS);
913:     } else if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_FACTOR_INFO) PetscFunctionReturn(PETSC_SUCCESS);
914:   }

916:   if (isdraw) {
917:     PetscDraw draw;
918:     PetscBool isnull;
919:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
920:     PetscCall(PetscDrawIsNull(draw, &isnull));
921:     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
922:   }

924:   {
925:     /* assemble the entire matrix onto first processor. */
926:     Mat           A;
927:     Mat_SeqSBAIJ *Aloc;
928:     Mat_SeqBAIJ  *Bloc;
929:     PetscInt      M = mat->rmap->N, N = mat->cmap->N, *ai, *aj, col, i, j, k, *rvals, mbs = baij->mbs;
930:     MatScalar    *a;
931:     const char   *matname;

933:     /* Should this be the same type as mat? */
934:     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
935:     if (rank == 0) {
936:       PetscCall(MatSetSizes(A, M, N, M, N));
937:     } else {
938:       PetscCall(MatSetSizes(A, 0, 0, M, N));
939:     }
940:     PetscCall(MatSetType(A, MATMPISBAIJ));
941:     PetscCall(MatMPISBAIJSetPreallocation(A, mat->rmap->bs, 0, NULL, 0, NULL));
942:     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));

944:     /* copy over the A part */
945:     Aloc = (Mat_SeqSBAIJ *)baij->A->data;
946:     ai   = Aloc->i;
947:     aj   = Aloc->j;
948:     a    = Aloc->a;
949:     PetscCall(PetscMalloc1(bs, &rvals));

951:     for (i = 0; i < mbs; i++) {
952:       rvals[0] = bs * (baij->rstartbs + i);
953:       for (j = 1; j < bs; j++) rvals[j] = rvals[j - 1] + 1;
954:       for (j = ai[i]; j < ai[i + 1]; j++) {
955:         col = (baij->cstartbs + aj[j]) * bs;
956:         for (k = 0; k < bs; k++) {
957:           PetscCall(MatSetValues_MPISBAIJ(A, bs, rvals, 1, &col, a, INSERT_VALUES));
958:           col++;
959:           a += bs;
960:         }
961:       }
962:     }
963:     /* copy over the B part */
964:     Bloc = (Mat_SeqBAIJ *)baij->B->data;
965:     ai   = Bloc->i;
966:     aj   = Bloc->j;
967:     a    = Bloc->a;
968:     for (i = 0; i < mbs; i++) {
969:       rvals[0] = bs * (baij->rstartbs + i);
970:       for (j = 1; j < bs; j++) rvals[j] = rvals[j - 1] + 1;
971:       for (j = ai[i]; j < ai[i + 1]; j++) {
972:         col = baij->garray[aj[j]] * bs;
973:         for (k = 0; k < bs; k++) {
974:           PetscCall(MatSetValues_MPIBAIJ(A, bs, rvals, 1, &col, a, INSERT_VALUES));
975:           col++;
976:           a += bs;
977:         }
978:       }
979:     }
980:     PetscCall(PetscFree(rvals));
981:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
982:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
983:     /*
984:        Everyone has to call to draw the matrix since the graphics waits are
985:        synchronized across all processors that share the PetscDraw object
986:     */
987:     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
988:     if (((PetscObject)mat)->name) PetscCall(PetscObjectGetName((PetscObject)mat, &matname));
989:     if (rank == 0) {
990:       if (((PetscObject)mat)->name) PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISBAIJ *)A->data)->A, matname));
991:       PetscCall(MatView_SeqSBAIJ(((Mat_MPISBAIJ *)A->data)->A, sviewer));
992:     }
993:     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
994:     PetscCall(MatDestroy(&A));
995:   }
996:   PetscFunctionReturn(PETSC_SUCCESS);
997: }

999: /* Used for both MPIBAIJ and MPISBAIJ matrices */
1000: #define MatView_MPISBAIJ_Binary MatView_MPIBAIJ_Binary

1002: static PetscErrorCode MatView_MPISBAIJ(Mat mat, PetscViewer viewer)
1003: {
1004:   PetscBool isascii, isdraw, issocket, isbinary;

1006:   PetscFunctionBegin;
1007:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1008:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1009:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
1010:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1011:   if (isascii || isdraw || issocket) PetscCall(MatView_MPISBAIJ_ASCIIorDraworSocket(mat, viewer));
1012:   else if (isbinary) PetscCall(MatView_MPISBAIJ_Binary(mat, viewer));
1013:   PetscFunctionReturn(PETSC_SUCCESS);
1014: }

1016: #if defined(PETSC_USE_COMPLEX)
1017: static PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A, Vec xx, Vec yy)
1018: {
1019:   Mat_MPISBAIJ      *a   = (Mat_MPISBAIJ *)A->data;
1020:   PetscInt           mbs = a->mbs, bs = A->rmap->bs;
1021:   PetscScalar       *from;
1022:   const PetscScalar *x;

1024:   PetscFunctionBegin;
1025:   /* diagonal part */
1026:   PetscUseTypeMethod(a->A, mult, xx, a->slvec1a);
1027:   /* since a->slvec1b shares memory (dangerously) with a->slec1 changes to a->slec1 will affect it */
1028:   PetscCall(PetscObjectStateIncrease((PetscObject)a->slvec1b));
1029:   PetscCall(VecZeroEntries(a->slvec1b));

1031:   /* subdiagonal part */
1032:   PetscUseTypeMethod(a->B, multhermitiantranspose, xx, a->slvec0b);

1034:   /* copy x into the vec slvec0 */
1035:   PetscCall(VecGetArray(a->slvec0, &from));
1036:   PetscCall(VecGetArrayRead(xx, &x));

1038:   PetscCall(PetscArraycpy(from, x, bs * mbs));
1039:   PetscCall(VecRestoreArray(a->slvec0, &from));
1040:   PetscCall(VecRestoreArrayRead(xx, &x));

1042:   PetscCall(VecScatterBegin(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1043:   PetscCall(VecScatterEnd(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1044:   /* supperdiagonal part */
1045:   PetscUseTypeMethod(a->B, multadd, a->slvec1b, a->slvec1a, yy);
1046:   PetscFunctionReturn(PETSC_SUCCESS);
1047: }
1048: #endif

1050: static PetscErrorCode MatMult_MPISBAIJ(Mat A, Vec xx, Vec yy)
1051: {
1052:   Mat_MPISBAIJ      *a   = (Mat_MPISBAIJ *)A->data;
1053:   PetscInt           mbs = a->mbs, bs = A->rmap->bs;
1054:   PetscScalar       *from;
1055:   const PetscScalar *x;

1057:   PetscFunctionBegin;
1058:   /* diagonal part */
1059:   PetscUseTypeMethod(a->A, mult, xx, a->slvec1a);
1060:   /* since a->slvec1b shares memory (dangerously) with a->slec1 changes to a->slec1 will affect it */
1061:   PetscCall(PetscObjectStateIncrease((PetscObject)a->slvec1b));
1062:   PetscCall(VecZeroEntries(a->slvec1b));

1064:   /* subdiagonal part */
1065:   PetscUseTypeMethod(a->B, multtranspose, xx, a->slvec0b);

1067:   /* copy x into the vec slvec0 */
1068:   PetscCall(VecGetArray(a->slvec0, &from));
1069:   PetscCall(VecGetArrayRead(xx, &x));

1071:   PetscCall(PetscArraycpy(from, x, bs * mbs));
1072:   PetscCall(VecRestoreArray(a->slvec0, &from));
1073:   PetscCall(VecRestoreArrayRead(xx, &x));

1075:   PetscCall(VecScatterBegin(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1076:   PetscCall(VecScatterEnd(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1077:   /* supperdiagonal part */
1078:   PetscUseTypeMethod(a->B, multadd, a->slvec1b, a->slvec1a, yy);
1079:   PetscFunctionReturn(PETSC_SUCCESS);
1080: }

1082: #if PetscDefined(USE_COMPLEX)
1083: static PetscErrorCode MatMultAdd_MPISBAIJ_Hermitian(Mat A, Vec xx, Vec yy, Vec zz)
1084: {
1085:   Mat_MPISBAIJ      *a   = (Mat_MPISBAIJ *)A->data;
1086:   PetscInt           mbs = a->mbs, bs = A->rmap->bs;
1087:   PetscScalar       *from;
1088:   const PetscScalar *x;

1090:   PetscFunctionBegin;
1091:   /* diagonal part */
1092:   PetscUseTypeMethod(a->A, multadd, xx, yy, a->slvec1a);
1093:   PetscCall(PetscObjectStateIncrease((PetscObject)a->slvec1b));
1094:   PetscCall(VecZeroEntries(a->slvec1b));

1096:   /* subdiagonal part */
1097:   PetscUseTypeMethod(a->B, multhermitiantranspose, xx, a->slvec0b);

1099:   /* copy x into the vec slvec0 */
1100:   PetscCall(VecGetArray(a->slvec0, &from));
1101:   PetscCall(VecGetArrayRead(xx, &x));
1102:   PetscCall(PetscArraycpy(from, x, bs * mbs));
1103:   PetscCall(VecRestoreArray(a->slvec0, &from));

1105:   PetscCall(VecScatterBegin(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1106:   PetscCall(VecRestoreArrayRead(xx, &x));
1107:   PetscCall(VecScatterEnd(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));

1109:   /* supperdiagonal part */
1110:   PetscUseTypeMethod(a->B, multadd, a->slvec1b, a->slvec1a, zz);
1111:   PetscFunctionReturn(PETSC_SUCCESS);
1112: }
1113: #endif

1115: static PetscErrorCode MatMultAdd_MPISBAIJ(Mat A, Vec xx, Vec yy, Vec zz)
1116: {
1117:   Mat_MPISBAIJ      *a   = (Mat_MPISBAIJ *)A->data;
1118:   PetscInt           mbs = a->mbs, bs = A->rmap->bs;
1119:   PetscScalar       *from;
1120:   const PetscScalar *x;

1122:   PetscFunctionBegin;
1123:   /* diagonal part */
1124:   PetscUseTypeMethod(a->A, multadd, xx, yy, a->slvec1a);
1125:   PetscCall(PetscObjectStateIncrease((PetscObject)a->slvec1b));
1126:   PetscCall(VecZeroEntries(a->slvec1b));

1128:   /* subdiagonal part */
1129:   PetscUseTypeMethod(a->B, multtranspose, xx, a->slvec0b);

1131:   /* copy x into the vec slvec0 */
1132:   PetscCall(VecGetArray(a->slvec0, &from));
1133:   PetscCall(VecGetArrayRead(xx, &x));
1134:   PetscCall(PetscArraycpy(from, x, bs * mbs));
1135:   PetscCall(VecRestoreArray(a->slvec0, &from));

1137:   PetscCall(VecScatterBegin(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));
1138:   PetscCall(VecRestoreArrayRead(xx, &x));
1139:   PetscCall(VecScatterEnd(a->sMvctx, a->slvec0, a->slvec1, ADD_VALUES, SCATTER_FORWARD));

1141:   /* supperdiagonal part */
1142:   PetscUseTypeMethod(a->B, multadd, a->slvec1b, a->slvec1a, zz);
1143:   PetscFunctionReturn(PETSC_SUCCESS);
1144: }

1146: /*
1147:   This only works correctly for square matrices where the subblock A->A is the
1148:    diagonal block
1149: */
1150: static PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A, Vec v)
1151: {
1152:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1154:   PetscFunctionBegin;
1155:   /* PetscCheck(a->rmap->N == a->cmap->N,PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1156:   PetscCall(MatGetDiagonal(a->A, v));
1157:   PetscFunctionReturn(PETSC_SUCCESS);
1158: }

1160: static PetscErrorCode MatScale_MPISBAIJ(Mat A, PetscScalar aa)
1161: {
1162:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1164:   PetscFunctionBegin;
1165:   PetscCall(MatScale(a->A, aa));
1166:   PetscCall(MatScale(a->B, aa));
1167:   PetscFunctionReturn(PETSC_SUCCESS);
1168: }

1170: static PetscErrorCode MatGetRow_MPISBAIJ(Mat matin, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1171: {
1172:   Mat_MPISBAIJ *mat = (Mat_MPISBAIJ *)matin->data;
1173:   PetscScalar  *vworkA, *vworkB, **pvA, **pvB, *v_p;
1174:   PetscInt      bs = matin->rmap->bs, bs2 = mat->bs2, i, *cworkA, *cworkB, **pcA, **pcB;
1175:   PetscInt      nztot, nzA, nzB, lrow, brstart = matin->rmap->rstart, brend = matin->rmap->rend;
1176:   PetscInt     *cmap, *idx_p, cstart = mat->rstartbs;

1178:   PetscFunctionBegin;
1179:   PetscCheck(!mat->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active");
1180:   mat->getrowactive = PETSC_TRUE;

1182:   if (!mat->rowvalues && (idx || v)) {
1183:     /*
1184:         allocate enough space to hold information from the longest row.
1185:     */
1186:     Mat_SeqSBAIJ *Aa  = (Mat_SeqSBAIJ *)mat->A->data;
1187:     Mat_SeqBAIJ  *Ba  = (Mat_SeqBAIJ *)mat->B->data;
1188:     PetscInt      max = 1, mbs = mat->mbs, tmp;
1189:     for (i = 0; i < mbs; i++) {
1190:       tmp = Aa->i[i + 1] - Aa->i[i] + Ba->i[i + 1] - Ba->i[i]; /* row length */
1191:       if (max < tmp) max = tmp;
1192:     }
1193:     PetscCall(PetscMalloc2(max * bs2, &mat->rowvalues, max * bs2, &mat->rowindices));
1194:   }

1196:   PetscCheck(row >= brstart && row < brend, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local rows");
1197:   lrow = row - brstart; /* local row index */

1199:   pvA = &vworkA;
1200:   pcA = &cworkA;
1201:   pvB = &vworkB;
1202:   pcB = &cworkB;
1203:   if (!v) {
1204:     pvA = NULL;
1205:     pvB = NULL;
1206:   }
1207:   if (!idx) {
1208:     pcA = NULL;
1209:     if (!v) pcB = NULL;
1210:   }
1211:   PetscUseTypeMethod(mat->A, getrow, lrow, &nzA, pcA, pvA);
1212:   PetscUseTypeMethod(mat->B, getrow, lrow, &nzB, pcB, pvB);
1213:   nztot = nzA + nzB;

1215:   cmap = mat->garray;
1216:   if (v || idx) {
1217:     if (nztot) {
1218:       /* Sort by increasing column numbers, assuming A and B already sorted */
1219:       PetscInt imark = -1;
1220:       if (v) {
1221:         *v = v_p = mat->rowvalues;
1222:         for (i = 0; i < nzB; i++) {
1223:           if (cmap[cworkB[i] / bs] < cstart) v_p[i] = vworkB[i];
1224:           else break;
1225:         }
1226:         imark = i;
1227:         for (i = 0; i < nzA; i++) v_p[imark + i] = vworkA[i];
1228:         for (i = imark; i < nzB; i++) v_p[nzA + i] = vworkB[i];
1229:       }
1230:       if (idx) {
1231:         *idx = idx_p = mat->rowindices;
1232:         if (imark > -1) {
1233:           for (i = 0; i < imark; i++) idx_p[i] = cmap[cworkB[i] / bs] * bs + cworkB[i] % bs;
1234:         } else {
1235:           for (i = 0; i < nzB; i++) {
1236:             if (cmap[cworkB[i] / bs] < cstart) idx_p[i] = cmap[cworkB[i] / bs] * bs + cworkB[i] % bs;
1237:             else break;
1238:           }
1239:           imark = i;
1240:         }
1241:         for (i = 0; i < nzA; i++) idx_p[imark + i] = cstart * bs + cworkA[i];
1242:         for (i = imark; i < nzB; i++) idx_p[nzA + i] = cmap[cworkB[i] / bs] * bs + cworkB[i] % bs;
1243:       }
1244:     } else {
1245:       if (idx) *idx = NULL;
1246:       if (v) *v = NULL;
1247:     }
1248:   }
1249:   *nz = nztot;
1250:   PetscUseTypeMethod(mat->A, restorerow, lrow, &nzA, pcA, pvA);
1251:   PetscUseTypeMethod(mat->B, restorerow, lrow, &nzB, pcB, pvB);
1252:   PetscFunctionReturn(PETSC_SUCCESS);
1253: }

1255: static PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
1256: {
1257:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;

1259:   PetscFunctionBegin;
1260:   PetscCheck(baij->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "MatGetRow() must be called first");
1261:   baij->getrowactive = PETSC_FALSE;
1262:   PetscFunctionReturn(PETSC_SUCCESS);
1263: }

1265: static PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1266: {
1267:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ *)A->data;
1268:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ *)a->A->data;

1270:   PetscFunctionBegin;
1271:   aA->getrow_utriangular = PETSC_TRUE;
1272:   PetscFunctionReturn(PETSC_SUCCESS);
1273: }
1274: static PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1275: {
1276:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ *)A->data;
1277:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ *)a->A->data;

1279:   PetscFunctionBegin;
1280:   aA->getrow_utriangular = PETSC_FALSE;
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: static PetscErrorCode MatConjugate_MPISBAIJ(Mat mat)
1285: {
1286:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)mat->data;

1288:   PetscFunctionBegin;
1289:   PetscCall(MatConjugate(a->A));
1290:   PetscCall(MatConjugate(a->B));
1291:   PetscFunctionReturn(PETSC_SUCCESS);
1292: }

1294: static PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1295: {
1296:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1298:   PetscFunctionBegin;
1299:   PetscCall(MatRealPart(a->A));
1300:   PetscCall(MatRealPart(a->B));
1301:   PetscFunctionReturn(PETSC_SUCCESS);
1302: }

1304: static PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1305: {
1306:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1308:   PetscFunctionBegin;
1309:   PetscCall(MatImaginaryPart(a->A));
1310:   PetscCall(MatImaginaryPart(a->B));
1311:   PetscFunctionReturn(PETSC_SUCCESS);
1312: }

1314: /* Check if isrow is a subset of iscol_local, called by MatCreateSubMatrix_MPISBAIJ()
1315:    Input: isrow       - distributed(parallel),
1316:           iscol_local - locally owned (seq)
1317: */
1318: static PetscErrorCode ISEqual_private(IS isrow, IS iscol_local, PetscBool *flg)
1319: {
1320:   PetscInt        sz1, sz2, *a1, *a2, i, j, k, nmatch;
1321:   const PetscInt *ptr1, *ptr2;

1323:   PetscFunctionBegin;
1324:   *flg = PETSC_FALSE;
1325:   PetscCall(ISGetLocalSize(isrow, &sz1));
1326:   PetscCall(ISGetLocalSize(iscol_local, &sz2));
1327:   if (sz1 > sz2) PetscFunctionReturn(PETSC_SUCCESS);

1329:   PetscCall(ISGetIndices(isrow, &ptr1));
1330:   PetscCall(ISGetIndices(iscol_local, &ptr2));

1332:   PetscCall(PetscMalloc1(sz1, &a1));
1333:   PetscCall(PetscMalloc1(sz2, &a2));
1334:   PetscCall(PetscArraycpy(a1, ptr1, sz1));
1335:   PetscCall(PetscArraycpy(a2, ptr2, sz2));
1336:   PetscCall(PetscSortInt(sz1, a1));
1337:   PetscCall(PetscSortInt(sz2, a2));

1339:   nmatch = 0;
1340:   k      = 0;
1341:   for (i = 0; i < sz1; i++) {
1342:     for (j = k; j < sz2; j++) {
1343:       if (a1[i] == a2[j]) {
1344:         k = j;
1345:         nmatch++;
1346:         break;
1347:       }
1348:     }
1349:   }
1350:   PetscCall(ISRestoreIndices(isrow, &ptr1));
1351:   PetscCall(ISRestoreIndices(iscol_local, &ptr2));
1352:   PetscCall(PetscFree(a1));
1353:   PetscCall(PetscFree(a2));
1354:   if (nmatch < sz1) {
1355:     *flg = PETSC_FALSE;
1356:   } else {
1357:     *flg = PETSC_TRUE;
1358:   }
1359:   PetscFunctionReturn(PETSC_SUCCESS);
1360: }

1362: static PetscErrorCode MatCreateSubMatrix_MPISBAIJ(Mat mat, IS isrow, IS iscol, MatReuse call, Mat *newmat)
1363: {
1364:   Mat       C[2];
1365:   IS        iscol_local, isrow_local;
1366:   PetscInt  csize, csize_local, rsize;
1367:   PetscBool isequal, issorted, isidentity = PETSC_FALSE;

1369:   PetscFunctionBegin;
1370:   PetscCall(ISGetLocalSize(iscol, &csize));
1371:   PetscCall(ISGetLocalSize(isrow, &rsize));
1372:   if (call == MAT_REUSE_MATRIX) {
1373:     PetscCall(PetscObjectQuery((PetscObject)*newmat, "ISAllGather", (PetscObject *)&iscol_local));
1374:     PetscCheck(iscol_local, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Submatrix passed in was not used before, cannot reuse");
1375:   } else {
1376:     PetscCall(ISAllGather(iscol, &iscol_local));
1377:     PetscCall(ISSorted(iscol_local, &issorted));
1378:     PetscCheck(issorted, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "For symmetric format, iscol must be sorted");
1379:   }
1380:   PetscCall(ISEqual_private(isrow, iscol_local, &isequal));
1381:   if (!isequal) {
1382:     PetscCall(ISGetLocalSize(iscol_local, &csize_local));
1383:     isidentity = (PetscBool)(mat->cmap->N == csize_local);
1384:     if (!isidentity) {
1385:       if (call == MAT_REUSE_MATRIX) {
1386:         PetscCall(PetscObjectQuery((PetscObject)*newmat, "ISAllGather_other", (PetscObject *)&isrow_local));
1387:         PetscCheck(isrow_local, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Submatrix passed in was not used before, cannot reuse");
1388:       } else {
1389:         PetscCall(ISAllGather(isrow, &isrow_local));
1390:         PetscCall(ISSorted(isrow_local, &issorted));
1391:         PetscCheck(issorted, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "For symmetric format, isrow must be sorted");
1392:       }
1393:     }
1394:   }
1395:   /* now call MatCreateSubMatrix_MPIBAIJ() */
1396:   PetscCall(MatCreateSubMatrix_MPIBAIJ_Private(mat, isrow, iscol_local, csize, isequal || isidentity ? call : MAT_INITIAL_MATRIX, isequal || isidentity ? newmat : C, (PetscBool)(isequal || isidentity)));
1397:   if (!isequal && !isidentity) {
1398:     if (call == MAT_INITIAL_MATRIX) {
1399:       IS       intersect;
1400:       PetscInt ni;

1402:       PetscCall(ISIntersect(isrow_local, iscol_local, &intersect));
1403:       PetscCall(ISGetLocalSize(intersect, &ni));
1404:       PetscCall(ISDestroy(&intersect));
1405:       PetscCheck(ni == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create such a submatrix: for symmetric format, when requesting an off-diagonal submatrix, isrow and iscol should have an empty intersection (number of common indices is %" PetscInt_FMT ")", ni);
1406:     }
1407:     PetscCall(MatCreateSubMatrix_MPIBAIJ_Private(mat, iscol, isrow_local, rsize, MAT_INITIAL_MATRIX, C + 1, PETSC_FALSE));
1408:     PetscCall(MatTranspose(C[1], MAT_INPLACE_MATRIX, C + 1));
1409:     PetscCall(MatAXPY(C[0], 1.0, C[1], DIFFERENT_NONZERO_PATTERN));
1410:     if (call == MAT_REUSE_MATRIX) PetscCall(MatCopy(C[0], *newmat, SAME_NONZERO_PATTERN));
1411:     else if (mat->rmap->bs == 1) PetscCall(MatConvert(C[0], MATAIJ, MAT_INITIAL_MATRIX, newmat));
1412:     else PetscCall(MatCopy(C[0], *newmat, SAME_NONZERO_PATTERN));
1413:     PetscCall(MatDestroy(C));
1414:     PetscCall(MatDestroy(C + 1));
1415:   }
1416:   if (call == MAT_INITIAL_MATRIX) {
1417:     if (!isequal && !isidentity) {
1418:       PetscCall(PetscObjectCompose((PetscObject)*newmat, "ISAllGather_other", (PetscObject)isrow_local));
1419:       PetscCall(ISDestroy(&isrow_local));
1420:     }
1421:     PetscCall(PetscObjectCompose((PetscObject)*newmat, "ISAllGather", (PetscObject)iscol_local));
1422:     PetscCall(ISDestroy(&iscol_local));
1423:   }
1424:   PetscFunctionReturn(PETSC_SUCCESS);
1425: }

1427: static PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1428: {
1429:   Mat_MPISBAIJ *l = (Mat_MPISBAIJ *)A->data;

1431:   PetscFunctionBegin;
1432:   PetscCall(MatZeroEntries(l->A));
1433:   PetscCall(MatZeroEntries(l->B));
1434:   PetscFunctionReturn(PETSC_SUCCESS);
1435: }

1437: static PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin, MatInfoType flag, MatInfo *info)
1438: {
1439:   Mat_MPISBAIJ  *a = (Mat_MPISBAIJ *)matin->data;
1440:   Mat            A = a->A, B = a->B;
1441:   PetscLogDouble isend[5], irecv[5];

1443:   PetscFunctionBegin;
1444:   info->block_size = (PetscReal)matin->rmap->bs;

1446:   PetscCall(MatGetInfo(A, MAT_LOCAL, info));

1448:   isend[0] = info->nz_used;
1449:   isend[1] = info->nz_allocated;
1450:   isend[2] = info->nz_unneeded;
1451:   isend[3] = info->memory;
1452:   isend[4] = info->mallocs;

1454:   PetscCall(MatGetInfo(B, MAT_LOCAL, info));

1456:   isend[0] += info->nz_used;
1457:   isend[1] += info->nz_allocated;
1458:   isend[2] += info->nz_unneeded;
1459:   isend[3] += info->memory;
1460:   isend[4] += info->mallocs;
1461:   if (flag == MAT_LOCAL) {
1462:     info->nz_used      = isend[0];
1463:     info->nz_allocated = isend[1];
1464:     info->nz_unneeded  = isend[2];
1465:     info->memory       = isend[3];
1466:     info->mallocs      = isend[4];
1467:   } else if (flag == MAT_GLOBAL_MAX) {
1468:     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)matin)));

1470:     info->nz_used      = irecv[0];
1471:     info->nz_allocated = irecv[1];
1472:     info->nz_unneeded  = irecv[2];
1473:     info->memory       = irecv[3];
1474:     info->mallocs      = irecv[4];
1475:   } else if (flag == MAT_GLOBAL_SUM) {
1476:     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)matin)));

1478:     info->nz_used      = irecv[0];
1479:     info->nz_allocated = irecv[1];
1480:     info->nz_unneeded  = irecv[2];
1481:     info->memory       = irecv[3];
1482:     info->mallocs      = irecv[4];
1483:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown MatInfoType argument %d", (int)flag);
1484:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1485:   info->fill_ratio_needed = 0;
1486:   info->factor_mallocs    = 0;
1487:   PetscFunctionReturn(PETSC_SUCCESS);
1488: }

1490: static PetscErrorCode MatSetOption_MPISBAIJ(Mat A, MatOption op, PetscBool flg)
1491: {
1492:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ *)A->data;
1493:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ *)a->A->data;

1495:   PetscFunctionBegin;
1496:   switch (op) {
1497:   case MAT_NEW_NONZERO_LOCATIONS:
1498:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1499:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1500:   case MAT_KEEP_NONZERO_PATTERN:
1501:   case MAT_NEW_NONZERO_LOCATION_ERR:
1502:     MatCheckPreallocated(A, 1);
1503:     PetscCall(MatSetOption(a->A, op, flg));
1504:     PetscCall(MatSetOption(a->B, op, flg));
1505:     break;
1506:   case MAT_ROW_ORIENTED:
1507:     MatCheckPreallocated(A, 1);
1508:     a->roworiented = flg;

1510:     PetscCall(MatSetOption(a->A, op, flg));
1511:     PetscCall(MatSetOption(a->B, op, flg));
1512:     break;
1513:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1514:     a->donotstash = flg;
1515:     break;
1516:   case MAT_USE_HASH_TABLE:
1517:     a->ht_flag = flg;
1518:     break;
1519:   case MAT_HERMITIAN:
1520:     if (a->A && A->rmap->n == A->cmap->n) PetscCall(MatSetOption(a->A, op, flg));
1521: #if defined(PETSC_USE_COMPLEX)
1522:     if (flg) { /* need different mat-vec ops */
1523:       A->ops->mult             = MatMult_MPISBAIJ_Hermitian;
1524:       A->ops->multadd          = MatMultAdd_MPISBAIJ_Hermitian;
1525:       A->ops->multtranspose    = NULL;
1526:       A->ops->multtransposeadd = NULL;
1527:     }
1528: #endif
1529:     break;
1530:   case MAT_SPD:
1531:   case MAT_SYMMETRIC:
1532:     if (a->A && A->rmap->n == A->cmap->n) PetscCall(MatSetOption(a->A, op, flg));
1533: #if defined(PETSC_USE_COMPLEX)
1534:     if (flg) { /* restore to use default mat-vec ops */
1535:       A->ops->mult             = MatMult_MPISBAIJ;
1536:       A->ops->multadd          = MatMultAdd_MPISBAIJ;
1537:       A->ops->multtranspose    = MatMult_MPISBAIJ;
1538:       A->ops->multtransposeadd = MatMultAdd_MPISBAIJ;
1539:     }
1540: #endif
1541:     break;
1542:   case MAT_STRUCTURALLY_SYMMETRIC:
1543:     if (a->A && A->rmap->n == A->cmap->n) PetscCall(MatSetOption(a->A, op, flg));
1544:     break;
1545:   case MAT_IGNORE_LOWER_TRIANGULAR:
1546:   case MAT_ERROR_LOWER_TRIANGULAR:
1547:     aA->ignore_ltriangular = flg;
1548:     break;
1549:   case MAT_GETROW_UPPERTRIANGULAR:
1550:     aA->getrow_utriangular = flg;
1551:     break;
1552:   default:
1553:     break;
1554:   }
1555:   PetscFunctionReturn(PETSC_SUCCESS);
1556: }

1558: static PetscErrorCode MatTranspose_MPISBAIJ(Mat A, MatReuse reuse, Mat *B)
1559: {
1560:   PetscFunctionBegin;
1561:   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
1562:   if (reuse == MAT_INITIAL_MATRIX) {
1563:     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, B));
1564:   } else if (reuse == MAT_REUSE_MATRIX) {
1565:     PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
1566:   }
1567:   PetscFunctionReturn(PETSC_SUCCESS);
1568: }

1570: static PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat, Vec ll, Vec rr)
1571: {
1572:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ *)mat->data;
1573:   Mat           a = baij->A, b = baij->B;
1574:   PetscInt      nv, m, n;

1576:   PetscFunctionBegin;
1577:   if (!ll) PetscFunctionReturn(PETSC_SUCCESS);

1579:   PetscCall(MatGetLocalSize(mat, &m, &n));
1580:   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "For symmetric format, local size %" PetscInt_FMT " %" PetscInt_FMT " must be same", m, n);

1582:   PetscCall(VecGetLocalSize(rr, &nv));
1583:   PetscCheck(nv == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left and right vector non-conforming local size");

1585:   PetscCall(VecScatterBegin(baij->Mvctx, rr, baij->lvec, INSERT_VALUES, SCATTER_FORWARD));

1587:   /* left diagonalscale the off-diagonal part */
1588:   PetscUseTypeMethod(b, diagonalscale, ll, NULL);

1590:   /* scale the diagonal part */
1591:   PetscUseTypeMethod(a, diagonalscale, ll, rr);

1593:   /* right diagonalscale the off-diagonal part */
1594:   PetscCall(VecScatterEnd(baij->Mvctx, rr, baij->lvec, INSERT_VALUES, SCATTER_FORWARD));
1595:   PetscUseTypeMethod(b, diagonalscale, NULL, baij->lvec);
1596:   PetscFunctionReturn(PETSC_SUCCESS);
1597: }

1599: static PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1600: {
1601:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1603:   PetscFunctionBegin;
1604:   PetscCall(MatSetUnfactored(a->A));
1605:   PetscFunctionReturn(PETSC_SUCCESS);
1606: }

1608: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat, MatDuplicateOption, Mat *);

1610: static PetscErrorCode MatEqual_MPISBAIJ(Mat A, Mat B, PetscBool *flag)
1611: {
1612:   Mat_MPISBAIJ *matB = (Mat_MPISBAIJ *)B->data, *matA = (Mat_MPISBAIJ *)A->data;
1613:   Mat           a, b, c, d;
1614:   PetscBool     flg;

1616:   PetscFunctionBegin;
1617:   a = matA->A;
1618:   b = matA->B;
1619:   c = matB->A;
1620:   d = matB->B;

1622:   PetscCall(MatEqual(a, c, &flg));
1623:   if (flg) PetscCall(MatEqual(b, d, &flg));
1624:   PetscCallMPI(MPIU_Allreduce(&flg, flag, 1, MPI_C_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
1625:   PetscFunctionReturn(PETSC_SUCCESS);
1626: }

1628: static PetscErrorCode MatCopy_MPISBAIJ(Mat A, Mat B, MatStructure str)
1629: {
1630:   PetscBool isbaij;

1632:   PetscFunctionBegin;
1633:   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &isbaij, MATSEQSBAIJ, MATMPISBAIJ, ""));
1634:   PetscCheck(isbaij, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)B)->type_name);
1635:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1636:   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1637:     PetscCall(MatGetRowUpperTriangular(A));
1638:     PetscCall(MatCopy_Basic(A, B, str));
1639:     PetscCall(MatRestoreRowUpperTriangular(A));
1640:   } else {
1641:     Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;
1642:     Mat_MPISBAIJ *b = (Mat_MPISBAIJ *)B->data;

1644:     PetscCall(MatCopy(a->A, b->A, str));
1645:     PetscCall(MatCopy(a->B, b->B, str));
1646:   }
1647:   PetscCall(PetscObjectStateIncrease((PetscObject)B));
1648:   PetscFunctionReturn(PETSC_SUCCESS);
1649: }

1651: static PetscErrorCode MatAXPY_MPISBAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
1652: {
1653:   Mat_MPISBAIJ *xx = (Mat_MPISBAIJ *)X->data, *yy = (Mat_MPISBAIJ *)Y->data;
1654:   PetscBLASInt  bnz, one                          = 1;
1655:   Mat_SeqSBAIJ *xa, *ya;
1656:   Mat_SeqBAIJ  *xb, *yb;

1658:   PetscFunctionBegin;
1659:   if (str == SAME_NONZERO_PATTERN) {
1660:     PetscScalar alpha = a;
1661:     xa                = (Mat_SeqSBAIJ *)xx->A->data;
1662:     ya                = (Mat_SeqSBAIJ *)yy->A->data;
1663:     PetscCall(PetscBLASIntCast(xa->nz, &bnz));
1664:     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xa->a, &one, ya->a, &one));
1665:     xb = (Mat_SeqBAIJ *)xx->B->data;
1666:     yb = (Mat_SeqBAIJ *)yy->B->data;
1667:     PetscCall(PetscBLASIntCast(xb->nz, &bnz));
1668:     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, xb->a, &one, yb->a, &one));
1669:     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
1670:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1671:     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_TRUE));
1672:     PetscCall(MatAXPY_Basic(Y, a, X, str));
1673:     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_FALSE));
1674:   } else {
1675:     Mat       B;
1676:     PetscInt *nnz_d, *nnz_o, bs = Y->rmap->bs;
1677:     PetscCheck(bs == X->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrices must have same block size");
1678:     PetscCall(MatGetRowUpperTriangular(X));
1679:     PetscCall(MatGetRowUpperTriangular(Y));
1680:     PetscCall(PetscMalloc1(yy->A->rmap->N, &nnz_d));
1681:     PetscCall(PetscMalloc1(yy->B->rmap->N, &nnz_o));
1682:     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
1683:     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
1684:     PetscCall(MatSetSizes(B, Y->rmap->n, Y->cmap->n, Y->rmap->N, Y->cmap->N));
1685:     PetscCall(MatSetBlockSizesFromMats(B, Y, Y));
1686:     PetscCall(MatSetType(B, MATMPISBAIJ));
1687:     PetscCall(MatAXPYGetPreallocation_SeqSBAIJ(yy->A, xx->A, nnz_d));
1688:     PetscCall(MatAXPYGetPreallocation_MPIBAIJ(yy->B, yy->garray, xx->B, xx->garray, nnz_o));
1689:     PetscCall(MatMPISBAIJSetPreallocation(B, bs, 0, nnz_d, 0, nnz_o));
1690:     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
1691:     PetscCall(MatHeaderMerge(Y, &B));
1692:     PetscCall(PetscFree(nnz_d));
1693:     PetscCall(PetscFree(nnz_o));
1694:     PetscCall(MatRestoreRowUpperTriangular(X));
1695:     PetscCall(MatRestoreRowUpperTriangular(Y));
1696:   }
1697:   PetscFunctionReturn(PETSC_SUCCESS);
1698: }

1700: static PetscErrorCode MatCreateSubMatrices_MPISBAIJ(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *B[])
1701: {
1702:   PetscInt  i;
1703:   PetscBool flg;

1705:   PetscFunctionBegin;
1706:   PetscCall(MatCreateSubMatrices_MPIBAIJ(A, n, irow, icol, scall, B)); /* B[] are sbaij matrices */
1707:   for (i = 0; i < n; i++) {
1708:     PetscCall(ISEqual(irow[i], icol[i], &flg));
1709:     if (!flg) PetscCall(MatSeqSBAIJZeroOps_Private(*B[i]));
1710:   }
1711:   PetscFunctionReturn(PETSC_SUCCESS);
1712: }

1714: static PetscErrorCode MatShift_MPISBAIJ(Mat Y, PetscScalar a)
1715: {
1716:   Mat_MPISBAIJ *maij = (Mat_MPISBAIJ *)Y->data;
1717:   Mat_SeqSBAIJ *aij  = (Mat_SeqSBAIJ *)maij->A->data;

1719:   PetscFunctionBegin;
1720:   if (!Y->preallocated) {
1721:     PetscCall(MatMPISBAIJSetPreallocation(Y, Y->rmap->bs, 1, NULL, 0, NULL));
1722:   } else if (!aij->nz) {
1723:     PetscInt nonew = aij->nonew;
1724:     PetscCall(MatSeqSBAIJSetPreallocation(maij->A, Y->rmap->bs, 1, NULL));
1725:     aij->nonew = nonew;
1726:   }
1727:   PetscCall(MatShift_Basic(Y, a));
1728:   PetscFunctionReturn(PETSC_SUCCESS);
1729: }

1731: static PetscErrorCode MatZeroRowsColumns_MPISBAIJ(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
1732: {
1733:   Mat_MPISBAIJ      *l = (Mat_MPISBAIJ *)A->data;
1734:   PetscMPIInt        n, p = 0;
1735:   PetscInt           i, j, k, r, len = 0, row, col, count;
1736:   PetscInt          *lrows, *owners = A->rmap->range;
1737:   PetscSFNode       *rrows;
1738:   PetscSF            sf;
1739:   const PetscScalar *xx;
1740:   PetscScalar       *bb, *mask;
1741:   Vec                xmask, lmask, lvec_contrib = NULL;
1742:   Mat_SeqBAIJ       *baij = (Mat_SeqBAIJ *)l->B->data;
1743:   PetscInt           bs = A->rmap->bs, bs2 = baij->bs2;
1744:   PetscScalar       *aa;

1746:   PetscFunctionBegin;
1747:   PetscCall(PetscMPIIntCast(A->rmap->n, &n));
1748:   /* create PetscSF where leaves are input rows and roots are owned rows */
1749:   PetscCall(PetscMalloc1(n, &lrows));
1750:   for (r = 0; r < n; ++r) lrows[r] = -1;
1751:   PetscCall(PetscMalloc1(N, &rrows));
1752:   for (r = 0; r < N; ++r) {
1753:     const PetscInt idx = rows[r];
1754:     PetscCheck(idx >= 0 && A->rmap->N > idx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, A->rmap->N);
1755:     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this row too */
1756:       PetscCall(PetscLayoutFindOwner(A->rmap, idx, &p));
1757:     }
1758:     rrows[r].rank  = p;
1759:     rrows[r].index = rows[r] - owners[p];
1760:   }
1761:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
1762:   PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, rrows, PETSC_OWN_POINTER));
1763:   /* collect flags for rows to be zeroed */
1764:   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)rows, lrows, MPI_LOR));
1765:   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)rows, lrows, MPI_LOR));
1766:   PetscCall(PetscSFDestroy(&sf));
1767:   /* compress and put in row numbers */
1768:   for (r = 0; r < n; ++r) {
1769:     if (lrows[r] >= 0) lrows[len++] = r;
1770:   }
1771:   /* zero diagonal part of matrix */
1772:   PetscCall(MatZeroRowsColumns(l->A, len, lrows, diag, x, b));
1773:   /* handle off-diagonal part of matrix */
1774:   PetscCall(MatCreateVecs(A, &xmask, NULL));
1775:   PetscCall(VecDuplicate(l->lvec, &lmask));
1776:   PetscCall(VecGetArray(xmask, &bb));
1777:   for (i = 0; i < len; i++) bb[lrows[i]] = 1;
1778:   PetscCall(VecRestoreArray(xmask, &bb));
1779:   PetscCall(VecScatterBegin(l->Mvctx, xmask, lmask, ADD_VALUES, SCATTER_FORWARD));
1780:   PetscCall(VecScatterEnd(l->Mvctx, xmask, lmask, ADD_VALUES, SCATTER_FORWARD));
1781:   PetscCall(VecDestroy(&xmask));
1782:   if (x) {
1783:     PetscCall(VecScatterBegin(l->Mvctx, x, l->lvec, INSERT_VALUES, SCATTER_FORWARD));
1784:     PetscCall(VecScatterEnd(l->Mvctx, x, l->lvec, INSERT_VALUES, SCATTER_FORWARD));
1785:     PetscCall(VecGetArrayRead(l->lvec, &xx));
1786:     PetscCall(VecGetArray(b, &bb));
1787:   }
1788:   PetscCall(VecGetArray(lmask, &mask));
1789:   /* MPISBAIJ stores only the upper off-diagonal in l->B; for each zeroed local row r and
1790:      non-zeroed off-process column c in that row, accumulate -A[r,c] * x[r] into lvec_contrib.
1791:      A SCATTER_REVERSE below sends these contributions to b[c] on the owning (higher-rank)
1792:      process, the missing symmetric lower-triangular update. We skip entries where c is
1793:      also a zeroed row (mask[col] != 0) since b[c] = diag * x[c] is handled separately. */
1794:   if (x) {
1795:     const PetscScalar *x_vals;
1796:     PetscScalar       *c_vals;

1798:     PetscCall(VecDuplicate(l->lvec, &lvec_contrib));
1799:     PetscCall(VecGetArray(lvec_contrib, &c_vals));
1800:     PetscCall(VecGetArrayRead(x, &x_vals));
1801:     /* Only accumulate b[c] -= A[r,c] * x[r] when off-process col c is not also a zeroed row
1802:        (mask[c] non-zero means col c is zeroed, so b[c] = diag * x[c] is already set).
1803:        This mirrors the MatSeqSBAIJ pattern: if (zeroed[r] && !zeroed[c]) bb[c] -= A[r,c] * x[r].
1804:        c_vals is indexed by the local B column index. */
1805:     for (i = 0; i < len; ++i) {
1806:       row = lrows[i];
1807:       for (j = baij->i[row / bs]; j < baij->i[row / bs + 1]; ++j) {
1808:         for (k = 0; k < bs; ++k) {
1809:           col = baij->j[j] * bs + k;
1810:           if (!PetscAbsScalar(mask[col])) {
1811:             aa = baij->a + j * bs2 + (row % bs) + bs * k;
1812:             c_vals[col] -= aa[0] * x_vals[row];
1813:           }
1814:         }
1815:       }
1816:     }
1817:     PetscCall(VecRestoreArrayRead(x, &x_vals));
1818:     PetscCall(VecRestoreArray(lvec_contrib, &c_vals));
1819:   }
1820:   /* remove zeroed rows of off-diagonal matrix */
1821:   for (i = 0; i < len; ++i) {
1822:     row   = lrows[i];
1823:     count = (baij->i[row / bs + 1] - baij->i[row / bs]) * bs;
1824:     aa    = baij->a + baij->i[row / bs] * bs2 + (row % bs);
1825:     for (k = 0; k < count; ++k) {
1826:       aa[0] = 0.0;
1827:       aa += bs;
1828:     }
1829:   }
1830:   /* loop over all elements of off process part of matrix zeroing removed columns */
1831:   for (i = 0; i < l->B->rmap->N; ++i) {
1832:     row = i / bs;
1833:     for (j = baij->i[row]; j < baij->i[row + 1]; ++j) {
1834:       for (k = 0; k < bs; ++k) {
1835:         col = bs * baij->j[j] + k;
1836:         if (PetscAbsScalar(mask[col])) {
1837:           aa = baij->a + j * bs2 + (i % bs) + bs * k;
1838:           if (x) bb[i] -= aa[0] * xx[col];
1839:           aa[0] = 0.0;
1840:         }
1841:       }
1842:     }
1843:   }
1844:   if (x) {
1845:     PetscCall(VecRestoreArray(b, &bb));
1846:     PetscCall(VecRestoreArrayRead(l->lvec, &xx));
1847:     /* scatter the accumulated contributions to b[c] on higher-rank processes owning column c */
1848:     PetscCall(VecScatterBegin(l->Mvctx, lvec_contrib, b, ADD_VALUES, SCATTER_REVERSE));
1849:     PetscCall(VecScatterEnd(l->Mvctx, lvec_contrib, b, ADD_VALUES, SCATTER_REVERSE));
1850:     PetscCall(VecDestroy(&lvec_contrib));
1851:   }
1852:   PetscCall(VecRestoreArray(lmask, &mask));
1853:   PetscCall(VecDestroy(&lmask));
1854:   PetscCall(PetscFree(lrows));

1856:   /* only change matrix nonzero state if pattern was allowed to be changed */
1857:   if (!((Mat_SeqSBAIJ *)l->A->data)->nonew) {
1858:     PetscObjectState state = l->A->nonzerostate + l->B->nonzerostate;
1859:     PetscCallMPI(MPIU_Allreduce(&state, &A->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)A)));
1860:   }
1861:   PetscFunctionReturn(PETSC_SUCCESS);
1862: }

1864: static PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A, Mat *a)
1865: {
1866:   PetscFunctionBegin;
1867:   *a = ((Mat_MPISBAIJ *)A->data)->A;
1868:   PetscFunctionReturn(PETSC_SUCCESS);
1869: }

1871: static PetscErrorCode MatEliminateZeros_MPISBAIJ(Mat A, PetscBool keep)
1872: {
1873:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;

1875:   PetscFunctionBegin;
1876:   PetscCall(MatEliminateZeros_SeqSBAIJ(a->A, keep));       // possibly keep zero diagonal coefficients
1877:   PetscCall(MatEliminateZeros_SeqBAIJ(a->B, PETSC_FALSE)); // never keep zero diagonal coefficients
1878:   PetscFunctionReturn(PETSC_SUCCESS);
1879: }

1881: static PetscErrorCode MatLoad_MPISBAIJ(Mat, PetscViewer);
1882: static PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat, Vec, PetscInt[]);
1883: static PetscErrorCode MatSOR_MPISBAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);

1885: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1886:                                        MatGetRow_MPISBAIJ,
1887:                                        MatRestoreRow_MPISBAIJ,
1888:                                        MatMult_MPISBAIJ,
1889:                                        /*  4*/ MatMultAdd_MPISBAIJ,
1890:                                        MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1891:                                        MatMultAdd_MPISBAIJ,
1892:                                        NULL,
1893:                                        NULL,
1894:                                        NULL,
1895:                                        /* 10*/ NULL,
1896:                                        NULL,
1897:                                        NULL,
1898:                                        MatSOR_MPISBAIJ,
1899:                                        MatTranspose_MPISBAIJ,
1900:                                        /* 15*/ MatGetInfo_MPISBAIJ,
1901:                                        MatEqual_MPISBAIJ,
1902:                                        MatGetDiagonal_MPISBAIJ,
1903:                                        MatDiagonalScale_MPISBAIJ,
1904:                                        MatNorm_MPISBAIJ,
1905:                                        /* 20*/ MatAssemblyBegin_MPISBAIJ,
1906:                                        MatAssemblyEnd_MPISBAIJ,
1907:                                        MatSetOption_MPISBAIJ,
1908:                                        MatZeroEntries_MPISBAIJ,
1909:                                        /* 24*/ NULL,
1910:                                        NULL,
1911:                                        NULL,
1912:                                        NULL,
1913:                                        NULL,
1914:                                        /* 29*/ MatSetUp_MPI_Hash,
1915:                                        NULL,
1916:                                        NULL,
1917:                                        MatGetDiagonalBlock_MPISBAIJ,
1918:                                        NULL,
1919:                                        /* 34*/ MatDuplicate_MPISBAIJ,
1920:                                        NULL,
1921:                                        NULL,
1922:                                        NULL,
1923:                                        NULL,
1924:                                        /* 39*/ MatAXPY_MPISBAIJ,
1925:                                        MatCreateSubMatrices_MPISBAIJ,
1926:                                        MatIncreaseOverlap_MPISBAIJ,
1927:                                        MatGetValues_MPISBAIJ,
1928:                                        MatCopy_MPISBAIJ,
1929:                                        /* 44*/ NULL,
1930:                                        MatScale_MPISBAIJ,
1931:                                        MatShift_MPISBAIJ,
1932:                                        NULL,
1933:                                        MatZeroRowsColumns_MPISBAIJ,
1934:                                        /* 49*/ NULL,
1935:                                        NULL,
1936:                                        NULL,
1937:                                        NULL,
1938:                                        NULL,
1939:                                        /* 54*/ NULL,
1940:                                        NULL,
1941:                                        MatSetUnfactored_MPISBAIJ,
1942:                                        NULL,
1943:                                        MatSetValuesBlocked_MPISBAIJ,
1944:                                        /* 59*/ MatCreateSubMatrix_MPISBAIJ,
1945:                                        NULL,
1946:                                        NULL,
1947:                                        NULL,
1948:                                        NULL,
1949:                                        /* 64*/ NULL,
1950:                                        NULL,
1951:                                        NULL,
1952:                                        NULL,
1953:                                        MatGetRowMaxAbs_MPISBAIJ,
1954:                                        /* 69*/ NULL,
1955:                                        MatConvert_MPISBAIJ_Basic,
1956:                                        NULL,
1957:                                        NULL,
1958:                                        NULL,
1959:                                        NULL,
1960:                                        NULL,
1961:                                        NULL,
1962:                                        NULL,
1963:                                        MatLoad_MPISBAIJ,
1964:                                        /* 79*/ NULL,
1965:                                        NULL,
1966:                                        NULL,
1967:                                        NULL,
1968:                                        NULL,
1969:                                        /* 84*/ NULL,
1970:                                        NULL,
1971:                                        NULL,
1972:                                        NULL,
1973:                                        NULL,
1974:                                        /* 89*/ NULL,
1975:                                        NULL,
1976:                                        NULL,
1977:                                        NULL,
1978:                                        MatConjugate_MPISBAIJ,
1979:                                        /* 94*/ NULL,
1980:                                        NULL,
1981:                                        MatRealPart_MPISBAIJ,
1982:                                        MatImaginaryPart_MPISBAIJ,
1983:                                        MatGetRowUpperTriangular_MPISBAIJ,
1984:                                        /* 99*/ MatRestoreRowUpperTriangular_MPISBAIJ,
1985:                                        NULL,
1986:                                        NULL,
1987:                                        NULL,
1988:                                        NULL,
1989:                                        /*104*/ NULL,
1990:                                        NULL,
1991:                                        NULL,
1992:                                        NULL,
1993:                                        NULL,
1994:                                        /*109*/ NULL,
1995:                                        NULL,
1996:                                        NULL,
1997:                                        NULL,
1998:                                        NULL,
1999:                                        /*114*/ NULL,
2000:                                        NULL,
2001:                                        NULL,
2002:                                        NULL,
2003:                                        NULL,
2004:                                        /*119*/ NULL,
2005:                                        NULL,
2006:                                        NULL,
2007:                                        NULL,
2008:                                        NULL,
2009:                                        /*124*/ NULL,
2010:                                        MatSetBlockSizes_Default,
2011:                                        NULL,
2012:                                        NULL,
2013:                                        NULL,
2014:                                        /*129*/ MatCreateMPIMatConcatenateSeqMat_MPISBAIJ,
2015:                                        NULL,
2016:                                        NULL,
2017:                                        NULL,
2018:                                        NULL,
2019:                                        /*134*/ NULL,
2020:                                        MatEliminateZeros_MPISBAIJ,
2021:                                        NULL,
2022:                                        NULL,
2023:                                        NULL,
2024:                                        /*139*/ NULL,
2025:                                        MatCopyHashToXAIJ_MPI_Hash,
2026:                                        NULL,
2027:                                        NULL,
2028:                                        MatADot_Default,
2029:                                        /*144*/ MatANorm_Default,
2030:                                        NULL,
2031:                                        NULL};

2033: static PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B, PetscInt bs, PetscInt d_nz, const PetscInt *d_nnz, PetscInt o_nz, const PetscInt *o_nnz)
2034: {
2035:   Mat_MPISBAIJ *b = (Mat_MPISBAIJ *)B->data;
2036:   PetscInt      i, mbs, Mbs;
2037:   PetscMPIInt   size;

2039:   PetscFunctionBegin;
2040:   if (B->hash_active) {
2041:     B->ops[0]      = b->cops;
2042:     B->hash_active = PETSC_FALSE;
2043:   }
2044:   if (!B->preallocated) PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), bs, &B->bstash));
2045:   PetscCall(MatSetBlockSize(B, bs));
2046:   PetscCall(PetscLayoutSetUp(B->rmap));
2047:   PetscCall(PetscLayoutSetUp(B->cmap));
2048:   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));
2049:   PetscCheck(B->rmap->N <= B->cmap->N, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MPISBAIJ matrix cannot have more rows %" PetscInt_FMT " than columns %" PetscInt_FMT, B->rmap->N, B->cmap->N);
2050:   PetscCheck(B->rmap->n <= B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "MPISBAIJ matrix cannot have more local rows %" PetscInt_FMT " than columns %" PetscInt_FMT, B->rmap->n, B->cmap->n);

2052:   mbs = B->rmap->n / bs;
2053:   Mbs = B->rmap->N / bs;
2054:   PetscCheck(mbs * bs == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "No of local rows %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT, B->rmap->N, bs);

2056:   B->rmap->bs = bs;
2057:   b->bs2      = bs * bs;
2058:   b->mbs      = mbs;
2059:   b->Mbs      = Mbs;
2060:   b->nbs      = B->cmap->n / bs;
2061:   b->Nbs      = B->cmap->N / bs;

2063:   for (i = 0; i <= b->size; i++) b->rangebs[i] = B->rmap->range[i] / bs;
2064:   b->rstartbs = B->rmap->rstart / bs;
2065:   b->rendbs   = B->rmap->rend / bs;

2067:   b->cstartbs = B->cmap->rstart / bs;
2068:   b->cendbs   = B->cmap->rend / bs;

2070: #if defined(PETSC_USE_CTABLE)
2071:   PetscCall(PetscHMapIDestroy(&b->colmap));
2072: #else
2073:   PetscCall(PetscFree(b->colmap));
2074: #endif
2075:   PetscCall(PetscFree(b->garray));
2076:   PetscCall(VecDestroy(&b->lvec));
2077:   PetscCall(VecScatterDestroy(&b->Mvctx));
2078:   PetscCall(VecDestroy(&b->slvec0));
2079:   PetscCall(VecDestroy(&b->slvec0b));
2080:   PetscCall(VecDestroy(&b->slvec1));
2081:   PetscCall(VecDestroy(&b->slvec1a));
2082:   PetscCall(VecDestroy(&b->slvec1b));
2083:   PetscCall(VecScatterDestroy(&b->sMvctx));

2085:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));

2087:   MatSeqXAIJGetOptions_Private(b->B);
2088:   PetscCall(MatDestroy(&b->B));
2089:   PetscCall(MatCreate(PETSC_COMM_SELF, &b->B));
2090:   PetscCall(MatSetSizes(b->B, B->rmap->n, size > 1 ? B->cmap->N : 0, B->rmap->n, size > 1 ? B->cmap->N : 0));
2091:   PetscCall(MatSetType(b->B, MATSEQBAIJ));
2092:   MatSeqXAIJRestoreOptions_Private(b->B);

2094:   MatSeqXAIJGetOptions_Private(b->A);
2095:   PetscCall(MatDestroy(&b->A));
2096:   PetscCall(MatCreate(PETSC_COMM_SELF, &b->A));
2097:   PetscCall(MatSetSizes(b->A, B->rmap->n, B->cmap->n, B->rmap->n, B->cmap->n));
2098:   PetscCall(MatSetType(b->A, MATSEQSBAIJ));
2099:   MatSeqXAIJRestoreOptions_Private(b->A);

2101:   PetscCall(MatSeqSBAIJSetPreallocation(b->A, bs, d_nz, d_nnz));
2102:   PetscCall(MatSeqBAIJSetPreallocation(b->B, bs, o_nz, o_nnz));

2104:   B->preallocated  = PETSC_TRUE;
2105:   B->was_assembled = PETSC_FALSE;
2106:   B->assembled     = PETSC_FALSE;
2107:   PetscFunctionReturn(PETSC_SUCCESS);
2108: }

2110: static PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B, PetscInt bs, const PetscInt ii[], const PetscInt jj[], const PetscScalar V[])
2111: {
2112:   PetscInt        m, rstart, cend;
2113:   PetscInt        i, j, d, nz, bd, nz_max = 0, *d_nnz = NULL, *o_nnz = NULL;
2114:   const PetscInt *JJ          = NULL;
2115:   PetscScalar    *values      = NULL;
2116:   PetscBool       roworiented = ((Mat_MPISBAIJ *)B->data)->roworiented;
2117:   PetscBool       nooffprocentries;

2119:   PetscFunctionBegin;
2120:   PetscCheck(bs >= 1, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_OUTOFRANGE, "Invalid block size specified, must be positive but it is %" PetscInt_FMT, bs);
2121:   PetscCall(PetscLayoutSetBlockSize(B->rmap, bs));
2122:   PetscCall(PetscLayoutSetBlockSize(B->cmap, bs));
2123:   PetscCall(PetscLayoutSetUp(B->rmap));
2124:   PetscCall(PetscLayoutSetUp(B->cmap));
2125:   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));
2126:   m      = B->rmap->n / bs;
2127:   rstart = B->rmap->rstart / bs;
2128:   cend   = B->cmap->rend / bs;

2130:   PetscCheck(!ii[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ii[0] must be 0 but it is %" PetscInt_FMT, ii[0]);
2131:   PetscCall(PetscMalloc2(m, &d_nnz, m, &o_nnz));
2132:   for (i = 0; i < m; i++) {
2133:     nz = ii[i + 1] - ii[i];
2134:     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
2135:     /* count the ones on the diagonal and above, split into diagonal and off-diagonal portions. */
2136:     JJ = jj + ii[i];
2137:     bd = 0;
2138:     for (j = 0; j < nz; j++) {
2139:       if (*JJ >= i + rstart) break;
2140:       JJ++;
2141:       bd++;
2142:     }
2143:     d = 0;
2144:     for (; j < nz; j++) {
2145:       if (*JJ++ >= cend) break;
2146:       d++;
2147:     }
2148:     d_nnz[i] = d;
2149:     o_nnz[i] = nz - d - bd;
2150:     nz       = nz - bd;
2151:     nz_max   = PetscMax(nz_max, nz);
2152:   }
2153:   PetscCall(MatMPISBAIJSetPreallocation(B, bs, 0, d_nnz, 0, o_nnz));
2154:   PetscCall(MatSetOption(B, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
2155:   PetscCall(PetscFree2(d_nnz, o_nnz));

2157:   values = (PetscScalar *)V;
2158:   if (!values) PetscCall(PetscCalloc1(bs * bs * nz_max, &values));
2159:   for (i = 0; i < m; i++) {
2160:     PetscInt        row   = i + rstart;
2161:     PetscInt        ncols = ii[i + 1] - ii[i];
2162:     const PetscInt *icols = jj + ii[i];
2163:     if (bs == 1 || !roworiented) { /* block ordering matches the non-nested layout of MatSetValues so we can insert entire rows */
2164:       const PetscScalar *svals = values + (V ? (bs * bs * ii[i]) : 0);
2165:       PetscCall(MatSetValuesBlocked_MPISBAIJ(B, 1, &row, ncols, icols, svals, INSERT_VALUES));
2166:     } else { /* block ordering does not match so we can only insert one block at a time. */
2167:       PetscInt j;
2168:       for (j = 0; j < ncols; j++) {
2169:         const PetscScalar *svals = values + (V ? (bs * bs * (ii[i] + j)) : 0);
2170:         PetscCall(MatSetValuesBlocked_MPISBAIJ(B, 1, &row, 1, &icols[j], svals, INSERT_VALUES));
2171:       }
2172:     }
2173:   }

2175:   if (!V) PetscCall(PetscFree(values));
2176:   nooffprocentries    = B->nooffprocentries;
2177:   B->nooffprocentries = PETSC_TRUE;
2178:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2179:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2180:   B->nooffprocentries = nooffprocentries;

2182:   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
2183:   PetscFunctionReturn(PETSC_SUCCESS);
2184: }

2186: /*MC
2187:    MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2188:    based on block compressed sparse row format.  Only the upper triangular portion of the "diagonal" portion of
2189:    the matrix is stored.

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

2194:    Options Database Key:
2195: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to `MatSetFromOptions()`

2197:    Level: beginner

2199:    Note:
2200:      The number of rows in the matrix must be less than or equal to the number of columns. Similarly the number of rows in the
2201:      diagonal portion of the matrix of each process has to less than or equal the number of columns.

2203: .seealso: [](ch_matrices), `Mat`, `MATSBAIJ`, `MATBAIJ`, `MatCreateBAIJ()`, `MATSEQSBAIJ`, `MatType`
2204: M*/

2206: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2207: {
2208:   Mat_MPISBAIJ *b;
2209:   PetscBool     flg = PETSC_FALSE;

2211:   PetscFunctionBegin;
2212:   PetscCall(PetscNew(&b));
2213:   B->data   = (void *)b;
2214:   B->ops[0] = MatOps_Values;

2216:   B->ops->destroy = MatDestroy_MPISBAIJ;
2217:   B->ops->view    = MatView_MPISBAIJ;
2218:   B->assembled    = PETSC_FALSE;
2219:   B->insertmode   = NOT_SET_VALUES;

2221:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &b->rank));
2222:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &b->size));

2224:   /* build local table of row and column ownerships */
2225:   PetscCall(PetscMalloc1(b->size + 2, &b->rangebs));

2227:   /* build cache for off array entries formed */
2228:   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), 1, &B->stash));

2230:   b->donotstash  = PETSC_FALSE;
2231:   b->colmap      = NULL;
2232:   b->garray      = NULL;
2233:   b->roworiented = PETSC_TRUE;

2235:   /* stuff used in block assembly */
2236:   b->barray = NULL;

2238:   /* stuff used for matrix vector multiply */
2239:   b->lvec    = NULL;
2240:   b->Mvctx   = NULL;
2241:   b->slvec0  = NULL;
2242:   b->slvec0b = NULL;
2243:   b->slvec1  = NULL;
2244:   b->slvec1a = NULL;
2245:   b->slvec1b = NULL;
2246:   b->sMvctx  = NULL;

2248:   /* stuff for MatGetRow() */
2249:   b->rowindices   = NULL;
2250:   b->rowvalues    = NULL;
2251:   b->getrowactive = PETSC_FALSE;

2253:   /* hash table stuff */
2254:   b->ht           = NULL;
2255:   b->hd           = NULL;
2256:   b->ht_size      = 0;
2257:   b->ht_flag      = PETSC_FALSE;
2258:   b->ht_fact      = 0;
2259:   b->ht_total_ct  = 0;
2260:   b->ht_insert_ct = 0;

2262:   /* stuff for MatCreateSubMatrices_MPIBAIJ_local() */
2263:   b->ijonly = PETSC_FALSE;

2265:   b->in_loc = NULL;
2266:   b->v_loc  = NULL;
2267:   b->n_loc  = 0;

2269:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_MPISBAIJ));
2270:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_MPISBAIJ));
2271:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISBAIJSetPreallocation_C", MatMPISBAIJSetPreallocation_MPISBAIJ));
2272:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISBAIJSetPreallocationCSR_C", MatMPISBAIJSetPreallocationCSR_MPISBAIJ));
2273: #if defined(PETSC_HAVE_ELEMENTAL)
2274:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisbaij_elemental_C", MatConvert_MPISBAIJ_Elemental));
2275: #endif
2276: #if defined(PETSC_HAVE_SCALAPACK) && (defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL_DOUBLE))
2277:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisbaij_scalapack_C", MatConvert_SBAIJ_ScaLAPACK));
2278: #endif
2279:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisbaij_mpiaij_C", MatConvert_MPISBAIJ_Basic));
2280:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisbaij_mpibaij_C", MatConvert_MPISBAIJ_Basic));

2282:   B->symmetric                   = PETSC_BOOL3_TRUE;
2283:   B->structurally_symmetric      = PETSC_BOOL3_TRUE;
2284:   B->symmetry_eternal            = PETSC_TRUE;
2285:   B->structural_symmetry_eternal = PETSC_TRUE;
2286: #if !defined(PETSC_USE_COMPLEX)
2287:   B->hermitian = PETSC_BOOL3_TRUE;
2288: #endif

2290:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPISBAIJ));
2291:   PetscOptionsBegin(PetscObjectComm((PetscObject)B), NULL, "Options for loading MPISBAIJ matrix 1", "Mat");
2292:   PetscCall(PetscOptionsBool("-mat_use_hash_table", "Use hash table to save memory in constructing matrix", "MatSetOption", flg, &flg, NULL));
2293:   if (flg) {
2294:     PetscReal fact = 1.39;
2295:     PetscCall(MatSetOption(B, MAT_USE_HASH_TABLE, PETSC_TRUE));
2296:     PetscCall(PetscOptionsReal("-mat_use_hash_table", "Use hash table factor", "MatMPIBAIJSetHashTableFactor", fact, &fact, NULL));
2297:     if (fact <= 1.0) fact = 1.39;
2298:     PetscCall(MatMPIBAIJSetHashTableFactor(B, fact));
2299:     PetscCall(PetscInfo(B, "Hash table Factor used %5.2g\n", (double)fact));
2300:   }
2301:   PetscOptionsEnd();
2302:   PetscFunctionReturn(PETSC_SUCCESS);
2303: }

2305: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
2306: /*MC
2307:    MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.

2309:    This matrix type is identical to `MATSEQSBAIJ` when constructed with a single process communicator,
2310:    and `MATMPISBAIJ` otherwise.

2312:    Options Database Key:
2313: . -mat_type sbaij - sets the matrix type to `MATSBAIJ` during a call to `MatSetFromOptions()`

2315:   Level: beginner

2317: .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MATMPISBAIJ`, `MatCreateSBAIJ()`, `MATSEQBAIJ`, `MATMPIBAIJ`
2318: M*/

2320: /*@
2321:   MatMPISBAIJSetPreallocation - For good matrix assembly performance
2322:   the user should preallocate the matrix storage by setting the parameters
2323:   d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
2324:   performance can be increased by more than a factor of 50.

2326:   Collective

2328:   Input Parameters:
2329: + B     - the matrix
2330: . bs    - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2331:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2332: . d_nz  - number of block nonzeros per block row in diagonal portion of local
2333:           submatrix  (same for all local rows)
2334: . d_nnz - array containing the number of block nonzeros in the various block rows
2335:           in the upper triangular and diagonal part of the in diagonal portion of the local
2336:           (possibly different for each block row) or `NULL`.  If you plan to factor the matrix you must leave room
2337:           for the diagonal entry and set a value even if it is zero.
2338: . o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2339:           submatrix (same for all local rows).
2340: - o_nnz - array containing the number of nonzeros in the various block rows of the
2341:           off-diagonal portion of the local submatrix that is right of the diagonal
2342:           (possibly different for each block row) or `NULL`.

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

2349:   Level: intermediate

2351:   Notes:

2353:   If `PETSC_DECIDE` or `PETSC_DETERMINE` is used for a particular argument on one processor
2354:   than it must be used on all processors that share the object for that argument.

2356:   If the *_nnz parameter is given then the *_nz parameter is ignored

2358:   Storage Information:
2359:   For a square global matrix we define each processor's diagonal portion
2360:   to be its local rows and the corresponding columns (a square submatrix);
2361:   each processor's off-diagonal portion encompasses the remainder of the
2362:   local matrix (a rectangular submatrix).

2364:   The user can specify preallocated storage for the diagonal part of
2365:   the local submatrix with either `d_nz` or `d_nnz` (not both).  Set
2366:   `d_nz` = `PETSC_DEFAULT` and `d_nnz` = `NULL` for PETSc to control dynamic
2367:   memory allocation.  Likewise, specify preallocated storage for the
2368:   off-diagonal part of the local submatrix with `o_nz` or `o_nnz` (not both).

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

2375:   Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2376:   the figure below we depict these three local rows and all columns (0-11).

2378: .vb
2379:            0 1 2 3 4 5 6 7 8 9 10 11
2380:           --------------------------
2381:    row 3  |. . . d d d o o o o  o  o
2382:    row 4  |. . . d d d o o o o  o  o
2383:    row 5  |. . . d d d o o o o  o  o
2384:           --------------------------
2385: .ve

2387:   Thus, any entries in the d locations are stored in the d (diagonal)
2388:   submatrix, and any entries in the o locations are stored in the
2389:   o (off-diagonal) submatrix.  Note that the d matrix is stored in
2390:   `MATSEQSBAIJ` format and the o submatrix in `MATSEQBAIJ` format.

2392:   Now `d_nz` should indicate the number of block nonzeros per row in the upper triangular
2393:   plus the diagonal part of the d matrix,
2394:   and `o_nz` should indicate the number of block nonzeros per row in the o matrix

2396:   In general, for PDE problems in which most nonzeros are near the diagonal,
2397:   one expects `d_nz` >> `o_nz`.

2399: .seealso: [](ch_matrices), `Mat`, `MATMPISBAIJ`, `MATSBAIJ`, `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValues()`, `MatCreateBAIJ()`, `PetscSplitOwnership()`
2400: @*/
2401: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
2402: {
2403:   PetscFunctionBegin;
2407:   PetscTryMethod(B, "MatMPISBAIJSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, bs, d_nz, d_nnz, o_nz, o_nnz));
2408:   PetscFunctionReturn(PETSC_SUCCESS);
2409: }

2411: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
2412: /*@
2413:   MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format, `MATSBAIJ`,
2414:   (block compressed row).  For good matrix assembly performance
2415:   the user should preallocate the matrix storage by setting the parameters
2416:   `d_nz` (or `d_nnz`) and `o_nz` (or `o_nnz`).

2418:   Collective

2420:   Input Parameters:
2421: + comm  - MPI communicator
2422: . bs    - size of block, the blocks are ALWAYS square. One can use `MatSetBlockSizes()` to set a different row and column blocksize but the row
2423:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with `MatCreateVecs()`
2424: . m     - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
2425:           This value should be the same as the local size used in creating the
2426:           y vector for the matrix-vector product y = Ax.
2427: . n     - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
2428:           This value should be the same as the local size used in creating the
2429:           x vector for the matrix-vector product y = Ax.
2430: . M     - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
2431: . N     - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
2432: . d_nz  - number of block nonzeros per block row in diagonal portion of local
2433:           submatrix (same for all local rows)
2434: . d_nnz - array containing the number of block nonzeros in the various block rows
2435:           in the upper triangular portion of the in diagonal portion of the local
2436:           (possibly different for each block block row) or `NULL`.
2437:           If you plan to factor the matrix you must leave room for the diagonal entry and
2438:           set its value even if it is zero.
2439: . o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2440:           submatrix (same for all local rows).
2441: - o_nnz - array containing the number of nonzeros in the various block rows of the
2442:           off-diagonal portion of the local submatrix (possibly different for
2443:           each block row) or `NULL`.

2445:   Output Parameter:
2446: . A - the matrix

2448:   Options Database Keys:
2449: + -mat_no_unroll  - uses code that does not unroll the loops in the
2450:                     block calculations (much slower)
2451: . -mat_block_size - size of the blocks to use
2452: - -mat_mpi        - use the parallel matrix data structures even on one processor
2453:                     (defaults to using SeqBAIJ format on one processor)

2455:   Level: intermediate

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

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

2465:   The user MUST specify either the local or global matrix dimensions
2466:   (possibly both).

2468:   If `PETSC_DECIDE` or `PETSC_DETERMINE` is used for a particular argument on one processor
2469:   than it must be used on all processors that share the object for that argument.

2471:   If `m` and `n` are not `PETSC_DECIDE`, then the values determines the `PetscLayout` of the matrix and the ranges returned by
2472:   `MatGetOwnershipRange()`,  `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, and `MatGetOwnershipRangesColumn()`.

2474:   If the *_nnz parameter is given then the *_nz parameter is ignored

2476:   Storage Information:
2477:   For a square global matrix we define each processor's diagonal portion
2478:   to be its local rows and the corresponding columns (a square submatrix);
2479:   each processor's off-diagonal portion encompasses the remainder of the
2480:   local matrix (a rectangular submatrix).

2482:   The user can specify preallocated storage for the diagonal part of
2483:   the local submatrix with either `d_nz` or `d_nnz` (not both). Set
2484:   `d_nz` = `PETSC_DEFAULT` and `d_nnz` = `NULL` for PETSc to control dynamic
2485:   memory allocation. Likewise, specify preallocated storage for the
2486:   off-diagonal part of the local submatrix with `o_nz` or `o_nnz` (not both).

2488:   Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2489:   the figure below we depict these three local rows and all columns (0-11).

2491: .vb
2492:            0 1 2 3 4 5 6 7 8 9 10 11
2493:           --------------------------
2494:    row 3  |. . . d d d o o o o  o  o
2495:    row 4  |. . . d d d o o o o  o  o
2496:    row 5  |. . . d d d o o o o  o  o
2497:           --------------------------
2498: .ve

2500:   Thus, any entries in the d locations are stored in the d (diagonal)
2501:   submatrix, and any entries in the o locations are stored in the
2502:   o (off-diagonal) submatrix. Note that the d matrix is stored in
2503:   `MATSEQSBAIJ` format and the o submatrix in `MATSEQBAIJ` format.

2505:   Now `d_nz` should indicate the number of block nonzeros per row in the upper triangular
2506:   plus the diagonal part of the d matrix,
2507:   and `o_nz` should indicate the number of block nonzeros per row in the o matrix.
2508:   In general, for PDE problems in which most nonzeros are near the diagonal,
2509:   one expects `d_nz` >> `o_nz`.

2511: .seealso: [](ch_matrices), `Mat`, `MATSBAIJ`, `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValues()`, `MatCreateBAIJ()`,
2512:           `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`
2513: @*/
2514: PetscErrorCode MatCreateSBAIJ(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[], Mat *A)
2515: {
2516:   PetscMPIInt size;

2518:   PetscFunctionBegin;
2519:   PetscCall(MatCreate(comm, A));
2520:   PetscCall(MatSetSizes(*A, m, n, M, N));
2521:   PetscCallMPI(MPI_Comm_size(comm, &size));
2522:   if (size > 1) {
2523:     PetscCall(MatSetType(*A, MATMPISBAIJ));
2524:     PetscCall(MatMPISBAIJSetPreallocation(*A, bs, d_nz, d_nnz, o_nz, o_nnz));
2525:   } else {
2526:     PetscCall(MatSetType(*A, MATSEQSBAIJ));
2527:     PetscCall(MatSeqSBAIJSetPreallocation(*A, bs, d_nz, d_nnz));
2528:   }
2529:   PetscFunctionReturn(PETSC_SUCCESS);
2530: }

2532: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin, MatDuplicateOption cpvalues, Mat *newmat)
2533: {
2534:   Mat           mat;
2535:   Mat_MPISBAIJ *a, *oldmat = (Mat_MPISBAIJ *)matin->data;
2536:   PetscInt      len = 0, nt, bs = matin->rmap->bs, mbs = oldmat->mbs;
2537:   PetscScalar  *array;

2539:   PetscFunctionBegin;
2540:   *newmat = NULL;

2542:   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin), &mat));
2543:   PetscCall(MatSetSizes(mat, matin->rmap->n, matin->cmap->n, matin->rmap->N, matin->cmap->N));
2544:   PetscCall(MatSetType(mat, ((PetscObject)matin)->type_name));
2545:   PetscCall(PetscLayoutReference(matin->rmap, &mat->rmap));
2546:   PetscCall(PetscLayoutReference(matin->cmap, &mat->cmap));

2548:   if (matin->hash_active) PetscCall(MatSetUp(mat));
2549:   else {
2550:     mat->factortype   = matin->factortype;
2551:     mat->preallocated = PETSC_TRUE;
2552:     mat->assembled    = PETSC_TRUE;
2553:     mat->insertmode   = NOT_SET_VALUES;

2555:     a      = (Mat_MPISBAIJ *)mat->data;
2556:     a->bs2 = oldmat->bs2;
2557:     a->mbs = oldmat->mbs;
2558:     a->nbs = oldmat->nbs;
2559:     a->Mbs = oldmat->Mbs;
2560:     a->Nbs = oldmat->Nbs;

2562:     a->size         = oldmat->size;
2563:     a->rank         = oldmat->rank;
2564:     a->donotstash   = oldmat->donotstash;
2565:     a->roworiented  = oldmat->roworiented;
2566:     a->rowindices   = NULL;
2567:     a->rowvalues    = NULL;
2568:     a->getrowactive = PETSC_FALSE;
2569:     a->barray       = NULL;
2570:     a->rstartbs     = oldmat->rstartbs;
2571:     a->rendbs       = oldmat->rendbs;
2572:     a->cstartbs     = oldmat->cstartbs;
2573:     a->cendbs       = oldmat->cendbs;

2575:     /* hash table stuff */
2576:     a->ht           = NULL;
2577:     a->hd           = NULL;
2578:     a->ht_size      = 0;
2579:     a->ht_flag      = oldmat->ht_flag;
2580:     a->ht_fact      = oldmat->ht_fact;
2581:     a->ht_total_ct  = 0;
2582:     a->ht_insert_ct = 0;

2584:     PetscCall(PetscArraycpy(a->rangebs, oldmat->rangebs, a->size + 2));
2585:     if (oldmat->colmap) {
2586: #if defined(PETSC_USE_CTABLE)
2587:       PetscCall(PetscHMapIDuplicate(oldmat->colmap, &a->colmap));
2588: #else
2589:       PetscCall(PetscMalloc1(a->Nbs, &a->colmap));
2590:       PetscCall(PetscArraycpy(a->colmap, oldmat->colmap, a->Nbs));
2591: #endif
2592:     } else a->colmap = NULL;

2594:     if (oldmat->garray && (len = ((Mat_SeqBAIJ *)oldmat->B->data)->nbs)) {
2595:       PetscCall(PetscMalloc1(len, &a->garray));
2596:       PetscCall(PetscArraycpy(a->garray, oldmat->garray, len));
2597:     } else a->garray = NULL;

2599:     PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)matin), matin->rmap->bs, &mat->bstash));
2600:     PetscCall(VecDuplicate(oldmat->lvec, &a->lvec));
2601:     PetscCall(VecScatterCopy(oldmat->Mvctx, &a->Mvctx));

2603:     PetscCall(VecDuplicate(oldmat->slvec0, &a->slvec0));
2604:     PetscCall(VecDuplicate(oldmat->slvec1, &a->slvec1));

2606:     PetscCall(VecGetLocalSize(a->slvec1, &nt));
2607:     PetscCall(VecGetArray(a->slvec1, &array));
2608:     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs * mbs, array, &a->slvec1a));
2609:     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nt - bs * mbs, array + bs * mbs, &a->slvec1b));
2610:     PetscCall(VecRestoreArray(a->slvec1, &array));
2611:     PetscCall(VecGetArray(a->slvec0, &array));
2612:     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nt - bs * mbs, array + bs * mbs, &a->slvec0b));
2613:     PetscCall(VecRestoreArray(a->slvec0, &array));

2615:     /* ierr =  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2616:     PetscCall(PetscObjectReference((PetscObject)oldmat->sMvctx));
2617:     a->sMvctx = oldmat->sMvctx;

2619:     PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
2620:     PetscCall(MatDuplicate(oldmat->B, cpvalues, &a->B));
2621:   }
2622:   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist, &((PetscObject)mat)->qlist));
2623:   *newmat = mat;
2624:   PetscFunctionReturn(PETSC_SUCCESS);
2625: }

2627: /* Used for both MPIBAIJ and MPISBAIJ matrices */
2628: #define MatLoad_MPISBAIJ_Binary MatLoad_MPIBAIJ_Binary

2630: static PetscErrorCode MatLoad_MPISBAIJ(Mat mat, PetscViewer viewer)
2631: {
2632:   PetscBool isbinary;

2634:   PetscFunctionBegin;
2635:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2636:   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);
2637:   PetscCall(MatLoad_MPISBAIJ_Binary(mat, viewer));
2638:   PetscFunctionReturn(PETSC_SUCCESS);
2639: }

2641: static PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A, Vec v, PetscInt idx[])
2642: {
2643:   Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;
2644:   Mat_SeqBAIJ  *b = (Mat_SeqBAIJ *)a->B->data;
2645:   PetscReal     atmp;
2646:   PetscReal    *work, *svalues, *rvalues;
2647:   PetscInt      i, bs, mbs, *bi, *bj, brow, j, ncols, krow, kcol, col, row, Mbs, bcol;
2648:   PetscMPIInt   rank, size;
2649:   PetscInt     *rowners_bs, count, source;
2650:   PetscScalar  *va;
2651:   MatScalar    *ba;
2652:   MPI_Status    stat;

2654:   PetscFunctionBegin;
2655:   PetscCheck(!idx, PETSC_COMM_SELF, PETSC_ERR_SUP, "Send email to petsc-maint@mcs.anl.gov");
2656:   PetscCall(MatGetRowMaxAbs(a->A, v, NULL));
2657:   PetscCall(VecGetArray(v, &va));

2659:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
2660:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));

2662:   bs  = A->rmap->bs;
2663:   mbs = a->mbs;
2664:   Mbs = a->Mbs;
2665:   ba  = b->a;
2666:   bi  = b->i;
2667:   bj  = b->j;

2669:   /* find ownerships */
2670:   rowners_bs = A->rmap->range;

2672:   /* each proc creates an array to be distributed */
2673:   PetscCall(PetscCalloc1(bs * Mbs, &work));

2675:   /* row_max for B */
2676:   if (rank != size - 1) {
2677:     for (i = 0; i < mbs; i++) {
2678:       ncols = bi[1] - bi[0];
2679:       bi++;
2680:       brow = bs * i;
2681:       for (j = 0; j < ncols; j++) {
2682:         bcol = bs * (*bj);
2683:         for (kcol = 0; kcol < bs; kcol++) {
2684:           col = bcol + kcol;           /* local col index */
2685:           col += rowners_bs[rank + 1]; /* global col index */
2686:           for (krow = 0; krow < bs; krow++) {
2687:             atmp = PetscAbsScalar(*ba);
2688:             ba++;
2689:             row = brow + krow; /* local row index */
2690:             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2691:             if (work[col] < atmp) work[col] = atmp;
2692:           }
2693:         }
2694:         bj++;
2695:       }
2696:     }

2698:     /* send values to its owners */
2699:     for (PetscMPIInt dest = rank + 1; dest < size; dest++) {
2700:       svalues = work + rowners_bs[dest];
2701:       count   = rowners_bs[dest + 1] - rowners_bs[dest];
2702:       PetscCallMPI(MPIU_Send(svalues, count, MPIU_REAL, dest, rank, PetscObjectComm((PetscObject)A)));
2703:     }
2704:   }

2706:   /* receive values */
2707:   if (rank) {
2708:     rvalues = work;
2709:     count   = rowners_bs[rank + 1] - rowners_bs[rank];
2710:     for (source = 0; source < rank; source++) {
2711:       PetscCallMPI(MPIU_Recv(rvalues, count, MPIU_REAL, MPI_ANY_SOURCE, MPI_ANY_TAG, PetscObjectComm((PetscObject)A), &stat));
2712:       /* process values */
2713:       for (i = 0; i < count; i++) {
2714:         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2715:       }
2716:     }
2717:   }

2719:   PetscCall(VecRestoreArray(v, &va));
2720:   PetscCall(PetscFree(work));
2721:   PetscFunctionReturn(PETSC_SUCCESS);
2722: }

2724: static PetscErrorCode MatSOR_MPISBAIJ(Mat matin, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
2725: {
2726:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ *)matin->data;
2727:   PetscInt           mbs = mat->mbs, bs = matin->rmap->bs;
2728:   PetscScalar       *x, *ptr, *from;
2729:   Vec                bb1;
2730:   const PetscScalar *b;

2732:   PetscFunctionBegin;
2733:   PetscCheck(its > 0 && lits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive", its, lits);
2734:   PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "SSOR for block size > 1 is not yet implemented");

2736:   if (flag == SOR_APPLY_UPPER) {
2737:     PetscUseTypeMethod(mat->A, sor, bb, omega, flag, fshift, lits, 1, xx);
2738:     PetscFunctionReturn(PETSC_SUCCESS);
2739:   }

2741:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2742:     if (flag & SOR_ZERO_INITIAL_GUESS) {
2743:       PetscUseTypeMethod(mat->A, sor, bb, omega, flag, fshift, lits, lits, xx);
2744:       its--;
2745:     }

2747:     PetscCall(VecDuplicate(bb, &bb1));
2748:     while (its--) {
2749:       /* lower triangular part: slvec0b = - B^T*xx */
2750:       PetscUseTypeMethod(mat->B, multtranspose, xx, mat->slvec0b);

2752:       /* copy xx into slvec0a */
2753:       PetscCall(VecGetArray(mat->slvec0, &ptr));
2754:       PetscCall(VecGetArray(xx, &x));
2755:       PetscCall(PetscArraycpy(ptr, x, bs * mbs));
2756:       PetscCall(VecRestoreArray(mat->slvec0, &ptr));

2758:       PetscCall(VecScale(mat->slvec0, -1.0));

2760:       /* copy bb into slvec1a */
2761:       PetscCall(VecGetArray(mat->slvec1, &ptr));
2762:       PetscCall(VecGetArrayRead(bb, &b));
2763:       PetscCall(PetscArraycpy(ptr, b, bs * mbs));
2764:       PetscCall(VecRestoreArray(mat->slvec1, &ptr));

2766:       /* set slvec1b = 0 */
2767:       PetscCall(PetscObjectStateIncrease((PetscObject)mat->slvec1b));
2768:       PetscCall(VecZeroEntries(mat->slvec1b));

2770:       PetscCall(VecScatterBegin(mat->sMvctx, mat->slvec0, mat->slvec1, ADD_VALUES, SCATTER_FORWARD));
2771:       PetscCall(VecRestoreArray(xx, &x));
2772:       PetscCall(VecRestoreArrayRead(bb, &b));
2773:       PetscCall(VecScatterEnd(mat->sMvctx, mat->slvec0, mat->slvec1, ADD_VALUES, SCATTER_FORWARD));

2775:       /* upper triangular part: bb1 = bb1 - B*x */
2776:       PetscUseTypeMethod(mat->B, multadd, mat->slvec1b, mat->slvec1a, bb1);

2778:       /* local diagonal sweep */
2779:       PetscUseTypeMethod(mat->A, sor, bb1, omega, SOR_SYMMETRIC_SWEEP, fshift, lits, lits, xx);
2780:     }
2781:     PetscCall(VecDestroy(&bb1));
2782:   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2783:     PetscUseTypeMethod(mat->A, sor, bb, omega, flag, fshift, lits, 1, xx);
2784:   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2785:     PetscUseTypeMethod(mat->A, sor, bb, omega, flag, fshift, lits, 1, xx);
2786:   } else if (flag & SOR_EISENSTAT) {
2787:     Vec                xx1;
2788:     PetscBool          hasop;
2789:     const PetscScalar *diag;
2790:     PetscScalar       *sl, scale = (omega - 2.0) / omega;
2791:     PetscInt           i, n;

2793:     if (!mat->xx1) {
2794:       PetscCall(VecDuplicate(bb, &mat->xx1));
2795:       PetscCall(VecDuplicate(bb, &mat->bb1));
2796:     }
2797:     xx1 = mat->xx1;
2798:     bb1 = mat->bb1;

2800:     PetscUseTypeMethod(mat->A, sor, bb, omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP), fshift, lits, 1, xx);

2802:     if (!mat->diag) {
2803:       /* this is wrong for same matrix with new nonzero values */
2804:       PetscCall(MatCreateVecs(matin, &mat->diag, NULL));
2805:       PetscCall(MatGetDiagonal(matin, mat->diag));
2806:     }
2807:     PetscCall(MatHasOperation(matin, MATOP_MULT_DIAGONAL_BLOCK, &hasop));

2809:     if (hasop) {
2810:       PetscCall(MatMultDiagonalBlock(matin, xx, bb1));
2811:       PetscCall(VecAYPX(mat->slvec1a, scale, bb));
2812:     } else {
2813:       /*
2814:           These two lines are replaced by code that may be a bit faster for a good compiler
2815:       PetscCall(VecPointwiseMult(mat->slvec1a,mat->diag,xx));
2816:       PetscCall(VecAYPX(mat->slvec1a,scale,bb));
2817:       */
2818:       PetscCall(VecGetArray(mat->slvec1a, &sl));
2819:       PetscCall(VecGetArrayRead(mat->diag, &diag));
2820:       PetscCall(VecGetArrayRead(bb, &b));
2821:       PetscCall(VecGetArray(xx, &x));
2822:       PetscCall(VecGetLocalSize(xx, &n));
2823:       if (omega == 1.0) {
2824:         for (i = 0; i < n; i++) sl[i] = b[i] - diag[i] * x[i];
2825:         PetscCall(PetscLogFlops(2.0 * n));
2826:       } else {
2827:         for (i = 0; i < n; i++) sl[i] = b[i] + scale * diag[i] * x[i];
2828:         PetscCall(PetscLogFlops(3.0 * n));
2829:       }
2830:       PetscCall(VecRestoreArray(mat->slvec1a, &sl));
2831:       PetscCall(VecRestoreArrayRead(mat->diag, &diag));
2832:       PetscCall(VecRestoreArrayRead(bb, &b));
2833:       PetscCall(VecRestoreArray(xx, &x));
2834:     }

2836:     /* multiply off-diagonal portion of matrix */
2837:     PetscCall(PetscObjectStateIncrease((PetscObject)mat->slvec1b));
2838:     PetscCall(VecZeroEntries(mat->slvec1b));
2839:     PetscUseTypeMethod(mat->B, multtranspose, xx, mat->slvec0b);
2840:     PetscCall(VecGetArray(mat->slvec0, &from));
2841:     PetscCall(VecGetArray(xx, &x));
2842:     PetscCall(PetscArraycpy(from, x, bs * mbs));
2843:     PetscCall(VecRestoreArray(mat->slvec0, &from));
2844:     PetscCall(VecRestoreArray(xx, &x));
2845:     PetscCall(VecScatterBegin(mat->sMvctx, mat->slvec0, mat->slvec1, ADD_VALUES, SCATTER_FORWARD));
2846:     PetscCall(VecScatterEnd(mat->sMvctx, mat->slvec0, mat->slvec1, ADD_VALUES, SCATTER_FORWARD));
2847:     PetscUseTypeMethod(mat->B, multadd, mat->slvec1b, mat->slvec1a, mat->slvec1a);

2849:     /* local sweep */
2850:     PetscUseTypeMethod(mat->A, sor, mat->slvec1a, omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP), fshift, lits, 1, xx1);
2851:     PetscCall(VecAXPY(xx, 1.0, xx1));
2852:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatSORType is not supported for SBAIJ matrix format");
2853:   PetscFunctionReturn(PETSC_SUCCESS);
2854: }

2856: /*@
2857:   MatCreateMPISBAIJWithArrays - creates a `MATMPISBAIJ` matrix using arrays that contain in standard CSR format for the local rows.

2859:   Collective

2861:   Input Parameters:
2862: + comm - MPI communicator
2863: . bs   - the block size, only a block size of 1 is supported
2864: . m    - number of local rows (Cannot be `PETSC_DECIDE`)
2865: . n    - This value should be the same as the local size used in creating the
2866:          x vector for the matrix-vector product $ y = Ax $. (or `PETSC_DECIDE` to have
2867:          calculated if `N` is given) For square matrices `n` is almost always `m`.
2868: . M    - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
2869: . N    - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
2870: . 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
2871: . j    - column indices
2872: - a    - matrix values

2874:   Output Parameter:
2875: . mat - the matrix

2877:   Level: intermediate

2879:   Notes:
2880:   The `i`, `j`, and `a` arrays ARE copied by this routine into the internal format used by PETSc;
2881:   thus you CANNOT change the matrix entries by changing the values of `a` after you have
2882:   called this routine. Use `MatCreateMPIAIJWithSplitArrays()` to avoid needing to copy the arrays.

2884:   The `i` and `j` indices are 0 based, and `i` indices are indices corresponding to the local `j` array.

2886: .seealso: [](ch_matrices), `Mat`, `MATMPISBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatMPIAIJSetPreallocation()`, `MatMPIAIJSetPreallocationCSR()`,
2887:           `MATMPIAIJ`, `MatCreateAIJ()`, `MatCreateMPIAIJWithSplitArrays()`, `MatMPISBAIJSetPreallocationCSR()`
2888: @*/
2889: PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, const PetscInt i[], const PetscInt j[], const PetscScalar a[], Mat *mat)
2890: {
2891:   PetscFunctionBegin;
2892:   PetscCheck(!i[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");
2893:   PetscCheck(m >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "local number of rows (m) cannot be PETSC_DECIDE, or negative");
2894:   PetscCall(MatCreate(comm, mat));
2895:   PetscCall(MatSetSizes(*mat, m, n, M, N));
2896:   PetscCall(MatSetType(*mat, MATMPISBAIJ));
2897:   PetscCall(MatMPISBAIJSetPreallocationCSR(*mat, bs, i, j, a));
2898:   PetscFunctionReturn(PETSC_SUCCESS);
2899: }

2901: /*@
2902:   MatMPISBAIJSetPreallocationCSR - Creates a sparse parallel matrix in `MATMPISBAIJ` format using the given nonzero structure and (optional) numerical values

2904:   Collective

2906:   Input Parameters:
2907: + B  - the matrix
2908: . bs - the block size
2909: . i  - the indices into `j` for the start of each local row (indices start with zero)
2910: . j  - the column indices for each local row (indices start with zero) these must be sorted for each row
2911: - v  - optional values in the matrix, pass `NULL` if not provided

2913:   Level: advanced

2915:   Notes:
2916:   The `i`, `j`, and `v` arrays ARE copied by this routine into the internal format used by PETSc;
2917:   thus you CANNOT change the matrix entries by changing the values of `v` after you have
2918:   called this routine.

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

2923:   Any entries passed in that are below the diagonal are ignored

2925: .seealso: [](ch_matrices), `Mat`, `MATMPISBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatMPIBAIJSetPreallocation()`, `MatCreateAIJ()`, `MATMPIAIJ`,
2926:           `MatCreateMPISBAIJWithArrays()`
2927: @*/
2928: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B, PetscInt bs, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
2929: {
2930:   PetscFunctionBegin;
2931:   PetscTryMethod(B, "MatMPISBAIJSetPreallocationCSR_C", (Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, bs, i, j, v));
2932:   PetscFunctionReturn(PETSC_SUCCESS);
2933: }

2935: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
2936: {
2937:   PetscInt     m, N, i, rstart, nnz, Ii, bs, cbs;
2938:   PetscInt    *indx;
2939:   PetscScalar *values;

2941:   PetscFunctionBegin;
2942:   PetscCall(MatGetSize(inmat, &m, &N));
2943:   if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
2944:     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)inmat->data;
2945:     PetscInt     *dnz, *onz, mbs, Nbs, nbs;
2946:     PetscInt     *bindx, rmax = a->rmax, j;
2947:     PetscMPIInt   rank, size;

2949:     PetscCall(MatGetBlockSizes(inmat, &bs, &cbs));
2950:     mbs = m / bs;
2951:     Nbs = N / cbs;
2952:     if (n == PETSC_DECIDE) PetscCall(PetscSplitOwnershipBlock(comm, cbs, &n, &N));
2953:     nbs = n / cbs;

2955:     PetscCall(PetscMalloc1(rmax, &bindx));
2956:     MatPreallocateBegin(comm, mbs, nbs, dnz, onz); /* inline function, output __end and __rstart are used below */

2958:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2959:     PetscCallMPI(MPI_Comm_size(comm, &size));
2960:     if (rank == size - 1) {
2961:       /* Check sum(nbs) = Nbs */
2962:       PetscCheck(__end == Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Sum of local block columns %" PetscInt_FMT " != global block columns %" PetscInt_FMT, __end, Nbs);
2963:     }

2965:     rstart = __rstart; /* block rstart of *outmat; see inline function MatPreallocateBegin */
2966:     PetscCall(MatSetOption(inmat, MAT_GETROW_UPPERTRIANGULAR, PETSC_TRUE));
2967:     for (i = 0; i < mbs; i++) {
2968:       PetscCall(MatGetRow_SeqSBAIJ(inmat, i * bs, &nnz, &indx, NULL)); /* non-blocked nnz and indx */
2969:       nnz = nnz / bs;
2970:       for (j = 0; j < nnz; j++) bindx[j] = indx[j * bs] / bs;
2971:       PetscCall(MatPreallocateSet(i + rstart, nnz, bindx, dnz, onz));
2972:       PetscCall(MatRestoreRow_SeqSBAIJ(inmat, i * bs, &nnz, &indx, NULL));
2973:     }
2974:     PetscCall(MatSetOption(inmat, MAT_GETROW_UPPERTRIANGULAR, PETSC_FALSE));
2975:     PetscCall(PetscFree(bindx));

2977:     PetscCall(MatCreate(comm, outmat));
2978:     PetscCall(MatSetSizes(*outmat, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
2979:     PetscCall(MatSetBlockSizes(*outmat, bs, cbs));
2980:     PetscCall(MatSetType(*outmat, MATSBAIJ));
2981:     PetscCall(MatSeqSBAIJSetPreallocation(*outmat, bs, 0, dnz));
2982:     PetscCall(MatMPISBAIJSetPreallocation(*outmat, bs, 0, dnz, 0, onz));
2983:     MatPreallocateEnd(dnz, onz);
2984:   }

2986:   /* numeric phase */
2987:   PetscCall(MatGetBlockSizes(inmat, &bs, &cbs));
2988:   PetscCall(MatGetOwnershipRange(*outmat, &rstart, NULL));

2990:   PetscCall(MatSetOption(inmat, MAT_GETROW_UPPERTRIANGULAR, PETSC_TRUE));
2991:   for (i = 0; i < m; i++) {
2992:     PetscCall(MatGetRow_SeqSBAIJ(inmat, i, &nnz, &indx, &values));
2993:     Ii = i + rstart;
2994:     PetscCall(MatSetValues(*outmat, 1, &Ii, nnz, indx, values, INSERT_VALUES));
2995:     PetscCall(MatRestoreRow_SeqSBAIJ(inmat, i, &nnz, &indx, &values));
2996:   }
2997:   PetscCall(MatSetOption(inmat, MAT_GETROW_UPPERTRIANGULAR, PETSC_FALSE));
2998:   PetscCall(MatAssemblyBegin(*outmat, MAT_FINAL_ASSEMBLY));
2999:   PetscCall(MatAssemblyEnd(*outmat, MAT_FINAL_ASSEMBLY));
3000:   PetscFunctionReturn(PETSC_SUCCESS);
3001: }