Actual source code: matnest.c

  1: #include <../src/mat/impls/nest/matnestimpl.h>
  2: #include <../src/mat/impls/aij/seq/aij.h>
  3: #include <../src/mat/impls/shell/shell.h>
  4: #include <petscsf.h>

  6: static PetscErrorCode MatSetUp_NestIS_Private(Mat, PetscInt, const IS[], PetscInt, const IS[]);
  7: static PetscErrorCode MatCreateVecs_Nest(Mat, Vec *, Vec *);
  8: static PetscErrorCode MatReset_Nest(Mat);

 10: PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat, MatType, MatReuse, Mat *);

 12: /* private functions */
 13: static PetscErrorCode MatNestGetSizes_Private(Mat A, PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N)
 14: {
 15:   Mat_Nest *bA = (Mat_Nest *)A->data;
 16:   PetscInt  i, j;

 18:   PetscFunctionBegin;
 19:   *m = *n = *M = *N = 0;
 20:   for (i = 0; i < bA->nr; i++) { /* rows */
 21:     PetscInt sm, sM;
 22:     PetscCall(ISGetLocalSize(bA->isglobal.row[i], &sm));
 23:     PetscCall(ISGetSize(bA->isglobal.row[i], &sM));
 24:     *m += sm;
 25:     *M += sM;
 26:   }
 27:   for (j = 0; j < bA->nc; j++) { /* cols */
 28:     PetscInt sn, sN;
 29:     PetscCall(ISGetLocalSize(bA->isglobal.col[j], &sn));
 30:     PetscCall(ISGetSize(bA->isglobal.col[j], &sN));
 31:     *n += sn;
 32:     *N += sN;
 33:   }
 34:   PetscFunctionReturn(PETSC_SUCCESS);
 35: }

 37: /* operations */
 38: static PetscErrorCode MatMult_Nest(Mat A, Vec x, Vec y)
 39: {
 40:   Mat_Nest *bA = (Mat_Nest *)A->data;
 41:   Vec      *bx = bA->right, *by = bA->left;
 42:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

 44:   PetscFunctionBegin;
 45:   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(y, bA->isglobal.row[i], &by[i]));
 46:   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(x, bA->isglobal.col[i], &bx[i]));
 47:   for (i = 0; i < nr; i++) {
 48:     PetscCall(VecZeroEntries(by[i]));
 49:     for (j = 0; j < nc; j++) {
 50:       if (!bA->m[i][j]) continue;
 51:       /* y[i] <- y[i] + A[i][j] * x[j] */
 52:       PetscCall(MatMultAdd(bA->m[i][j], bx[j], by[i], by[i]));
 53:     }
 54:   }
 55:   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(y, bA->isglobal.row[i], &by[i]));
 56:   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.col[i], &bx[i]));
 57:   PetscFunctionReturn(PETSC_SUCCESS);
 58: }

 60: static PetscErrorCode MatMultAdd_Nest(Mat A, Vec x, Vec y, Vec z)
 61: {
 62:   Mat_Nest *bA = (Mat_Nest *)A->data;
 63:   Vec      *bx = bA->right, *bz = bA->left;
 64:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

 66:   PetscFunctionBegin;
 67:   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(z, bA->isglobal.row[i], &bz[i]));
 68:   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(x, bA->isglobal.col[i], &bx[i]));
 69:   for (i = 0; i < nr; i++) {
 70:     if (y != z) {
 71:       Vec by;
 72:       PetscCall(VecGetSubVector(y, bA->isglobal.row[i], &by));
 73:       PetscCall(VecCopy(by, bz[i]));
 74:       PetscCall(VecRestoreSubVector(y, bA->isglobal.row[i], &by));
 75:     }
 76:     for (j = 0; j < nc; j++) {
 77:       if (!bA->m[i][j]) continue;
 78:       /* y[i] <- y[i] + A[i][j] * x[j] */
 79:       PetscCall(MatMultAdd(bA->m[i][j], bx[j], bz[i], bz[i]));
 80:     }
 81:   }
 82:   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(z, bA->isglobal.row[i], &bz[i]));
 83:   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.col[i], &bx[i]));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: typedef struct {
 88:   Mat         *workC;      /* array of Mat with specific containers depending on the underlying MatMatMult implementation */
 89:   PetscScalar *tarray;     /* buffer for storing all temporary products A[i][j] B[j] */
 90:   PetscInt    *dm, *dn, k; /* displacements and number of submatrices */
 91: } Nest_Dense;

 93: static PetscErrorCode MatProductNumeric_Nest_Dense(Mat C)
 94: {
 95:   Mat_Nest          *bA;
 96:   Nest_Dense        *contents;
 97:   Mat                viewB, viewC, productB, workC;
 98:   const PetscScalar *barray;
 99:   PetscScalar       *carray;
100:   PetscInt           i, j, M, N, nr, nc, ldb, ldc;
101:   Mat                A, B;

103:   PetscFunctionBegin;
104:   MatCheckProduct(C, 1);
105:   A = C->product->A;
106:   B = C->product->B;
107:   PetscCall(MatGetSize(B, NULL, &N));
108:   if (!N) {
109:     PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
110:     PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
111:     PetscFunctionReturn(PETSC_SUCCESS);
112:   }
113:   contents = (Nest_Dense *)C->product->data;
114:   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
115:   bA = (Mat_Nest *)A->data;
116:   nr = bA->nr;
117:   nc = bA->nc;
118:   PetscCall(MatDenseGetLDA(B, &ldb));
119:   PetscCall(MatDenseGetLDA(C, &ldc));
120:   PetscCall(MatZeroEntries(C));
121:   PetscCall(MatDenseGetArrayRead(B, &barray));
122:   PetscCall(MatDenseGetArray(C, &carray));
123:   for (i = 0; i < nr; i++) {
124:     PetscCall(ISGetSize(bA->isglobal.row[i], &M));
125:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dm[i + 1] - contents->dm[i], PETSC_DECIDE, M, N, PetscSafePointerPlusOffset(carray, contents->dm[i]), &viewC));
126:     PetscCall(MatDenseSetLDA(viewC, ldc));
127:     for (j = 0; j < nc; j++) {
128:       if (!bA->m[i][j]) continue;
129:       PetscCall(ISGetSize(bA->isglobal.col[j], &M));
130:       PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dn[j + 1] - contents->dn[j], PETSC_DECIDE, M, N, PetscSafePointerPlusOffset((PetscScalar *)barray, contents->dn[j]), &viewB));
131:       PetscCall(MatDenseSetLDA(viewB, ldb));

133:       /* MatMatMultNumeric(bA->m[i][j],viewB,contents->workC[i*nc + j]); */
134:       workC             = contents->workC[i * nc + j];
135:       productB          = workC->product->B;
136:       workC->product->B = viewB; /* use newly created dense matrix viewB */
137:       PetscCall(MatProductNumeric(workC));
138:       PetscCall(MatDestroy(&viewB));
139:       workC->product->B = productB; /* resume original B */

141:       /* C[i] <- workC + C[i] */
142:       PetscCall(MatAXPY(viewC, 1.0, contents->workC[i * nc + j], SAME_NONZERO_PATTERN));
143:     }
144:     PetscCall(MatDestroy(&viewC));
145:   }
146:   PetscCall(MatDenseRestoreArray(C, &carray));
147:   PetscCall(MatDenseRestoreArrayRead(B, &barray));

149:   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
150:   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
151:   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: static PetscErrorCode MatNest_DenseDestroy(PetscCtxRt ctx)
156: {
157:   Nest_Dense *contents = *(Nest_Dense **)ctx;
158:   PetscInt    i;

160:   PetscFunctionBegin;
161:   PetscCall(PetscFree(contents->tarray));
162:   for (i = 0; i < contents->k; i++) PetscCall(MatDestroy(contents->workC + i));
163:   PetscCall(PetscFree3(contents->dm, contents->dn, contents->workC));
164:   PetscCall(PetscFree(contents));
165:   PetscFunctionReturn(PETSC_SUCCESS);
166: }

168: static PetscErrorCode MatProductSymbolic_Nest_Dense(Mat C)
169: {
170:   Mat_Nest          *bA;
171:   Mat                viewB, workC;
172:   const PetscScalar *barray;
173:   PetscInt           i, j, M, N, m, n, nr, nc, maxm = 0, ldb;
174:   Nest_Dense        *contents = NULL;
175:   PetscBool          cisdense;
176:   Mat                A, B;
177:   PetscReal          fill;

179:   PetscFunctionBegin;
180:   MatCheckProduct(C, 1);
181:   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
182:   A    = C->product->A;
183:   B    = C->product->B;
184:   fill = C->product->fill;
185:   bA   = (Mat_Nest *)A->data;
186:   nr   = bA->nr;
187:   nc   = bA->nc;
188:   PetscCall(MatGetLocalSize(C, &m, &n));
189:   PetscCall(MatGetSize(C, &M, &N));
190:   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
191:     PetscCall(MatGetLocalSize(B, NULL, &n));
192:     PetscCall(MatGetSize(B, NULL, &N));
193:     PetscCall(MatGetLocalSize(A, &m, NULL));
194:     PetscCall(MatGetSize(A, &M, NULL));
195:     PetscCall(MatSetSizes(C, m, n, M, N));
196:   }
197:   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
198:   if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
199:   PetscCall(MatSetUp(C));
200:   if (!N) {
201:     C->ops->productnumeric = MatProductNumeric_Nest_Dense;
202:     PetscFunctionReturn(PETSC_SUCCESS);
203:   }

205:   PetscCall(PetscNew(&contents));
206:   C->product->data    = contents;
207:   C->product->destroy = MatNest_DenseDestroy;
208:   PetscCall(PetscCalloc3(nr + 1, &contents->dm, nc + 1, &contents->dn, nr * nc, &contents->workC));
209:   contents->k = nr * nc;
210:   for (i = 0; i < nr; i++) {
211:     PetscCall(ISGetLocalSize(bA->isglobal.row[i], contents->dm + i + 1));
212:     maxm = PetscMax(maxm, contents->dm[i + 1]);
213:     contents->dm[i + 1] += contents->dm[i];
214:   }
215:   for (i = 0; i < nc; i++) {
216:     PetscCall(ISGetLocalSize(bA->isglobal.col[i], contents->dn + i + 1));
217:     contents->dn[i + 1] += contents->dn[i];
218:   }
219:   PetscCall(PetscMalloc1(maxm * N, &contents->tarray));
220:   PetscCall(MatDenseGetLDA(B, &ldb));
221:   PetscCall(MatGetSize(B, NULL, &N));
222:   PetscCall(MatDenseGetArrayRead(B, &barray));
223:   /* loops are permuted compared to MatMatMultNumeric so that viewB is created only once per column of A */
224:   for (j = 0; j < nc; j++) {
225:     PetscCall(ISGetSize(bA->isglobal.col[j], &M));
226:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dn[j + 1] - contents->dn[j], PETSC_DECIDE, M, N, PetscSafePointerPlusOffset((PetscScalar *)barray, contents->dn[j]), &viewB));
227:     PetscCall(MatDenseSetLDA(viewB, ldb));
228:     for (i = 0; i < nr; i++) {
229:       if (!bA->m[i][j]) continue;
230:       /* MatMatMultSymbolic may attach a specific container (depending on MatType of bA->m[i][j]) to workC[i][j] */

232:       PetscCall(MatProductCreate(bA->m[i][j], viewB, NULL, &contents->workC[i * nc + j]));
233:       workC = contents->workC[i * nc + j];
234:       PetscCall(MatProductSetType(workC, MATPRODUCT_AB));
235:       PetscCall(MatProductSetAlgorithm(workC, "default"));
236:       PetscCall(MatProductSetFill(workC, fill));
237:       PetscCall(MatProductSetFromOptions(workC));
238:       PetscCall(MatProductSymbolic(workC));

240:       /* since tarray will be shared by all Mat */
241:       PetscCall(MatSeqDenseSetPreallocation(workC, contents->tarray));
242:       PetscCall(MatMPIDenseSetPreallocation(workC, contents->tarray));
243:     }
244:     PetscCall(MatDestroy(&viewB));
245:   }
246:   PetscCall(MatDenseRestoreArrayRead(B, &barray));

248:   C->ops->productnumeric = MatProductNumeric_Nest_Dense;
249:   PetscFunctionReturn(PETSC_SUCCESS);
250: }

252: static PetscErrorCode MatProductSetFromOptions_Nest_Dense(Mat C)
253: {
254:   Mat_Product *product = C->product;

256:   PetscFunctionBegin;
257:   if (product->type == MATPRODUCT_AB) C->ops->productsymbolic = MatProductSymbolic_Nest_Dense;
258:   PetscFunctionReturn(PETSC_SUCCESS);
259: }

261: static PetscErrorCode MatMultTransposeKernel_Nest(Mat A, Vec x, Vec y, PetscBool herm)
262: {
263:   Mat_Nest *bA = (Mat_Nest *)A->data;
264:   Vec      *bx = bA->left, *by = bA->right;
265:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

267:   PetscFunctionBegin;
268:   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(x, bA->isglobal.row[i], &bx[i]));
269:   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(y, bA->isglobal.col[i], &by[i]));
270:   for (j = 0; j < nc; j++) {
271:     PetscCall(VecZeroEntries(by[j]));
272:     for (i = 0; i < nr; i++) {
273:       if (!bA->m[i][j]) continue;
274:       if (herm) PetscCall(MatMultHermitianTransposeAdd(bA->m[i][j], bx[i], by[j], by[j])); /* y[j] <- y[j] + (A[i][j])^H * x[i] */
275:       else PetscCall(MatMultTransposeAdd(bA->m[i][j], bx[i], by[j], by[j]));               /* y[j] <- y[j] + (A[i][j])^T * x[i] */
276:     }
277:   }
278:   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.row[i], &bx[i]));
279:   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(y, bA->isglobal.col[i], &by[i]));
280:   PetscFunctionReturn(PETSC_SUCCESS);
281: }

283: static PetscErrorCode MatMultTranspose_Nest(Mat A, Vec x, Vec y)
284: {
285:   PetscFunctionBegin;
286:   PetscCall(MatMultTransposeKernel_Nest(A, x, y, PETSC_FALSE));
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: static PetscErrorCode MatMultHermitianTranspose_Nest(Mat A, Vec x, Vec y)
291: {
292:   PetscFunctionBegin;
293:   PetscCall(MatMultTransposeKernel_Nest(A, x, y, PETSC_TRUE));
294:   PetscFunctionReturn(PETSC_SUCCESS);
295: }

297: static PetscErrorCode MatMultTransposeAddKernel_Nest(Mat A, Vec x, Vec y, Vec z, PetscBool herm)
298: {
299:   Mat_Nest *bA = (Mat_Nest *)A->data;
300:   Vec      *bx = bA->left, *bz = bA->right;
301:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

303:   PetscFunctionBegin;
304:   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(x, bA->isglobal.row[i], &bx[i]));
305:   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(z, bA->isglobal.col[i], &bz[i]));
306:   for (j = 0; j < nc; j++) {
307:     if (y != z) {
308:       Vec by;
309:       PetscCall(VecGetSubVector(y, bA->isglobal.col[j], &by));
310:       PetscCall(VecCopy(by, bz[j]));
311:       PetscCall(VecRestoreSubVector(y, bA->isglobal.col[j], &by));
312:     }
313:     for (i = 0; i < nr; i++) {
314:       if (!bA->m[i][j]) continue;
315:       if (herm) PetscCall(MatMultHermitianTransposeAdd(bA->m[i][j], bx[i], bz[j], bz[j])); /* z[j] <- y[j] + (A[i][j])^H * x[i] */
316:       else PetscCall(MatMultTransposeAdd(bA->m[i][j], bx[i], bz[j], bz[j]));               /* z[j] <- y[j] + (A[i][j])^T * x[i] */
317:     }
318:   }
319:   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.row[i], &bx[i]));
320:   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(z, bA->isglobal.col[i], &bz[i]));
321:   PetscFunctionReturn(PETSC_SUCCESS);
322: }

324: static PetscErrorCode MatMultTransposeAdd_Nest(Mat A, Vec x, Vec y, Vec z)
325: {
326:   PetscFunctionBegin;
327:   PetscCall(MatMultTransposeAddKernel_Nest(A, x, y, z, PETSC_FALSE));
328:   PetscFunctionReturn(PETSC_SUCCESS);
329: }

331: static PetscErrorCode MatMultHermitianTransposeAdd_Nest(Mat A, Vec x, Vec y, Vec z)
332: {
333:   PetscFunctionBegin;
334:   PetscCall(MatMultTransposeAddKernel_Nest(A, x, y, z, PETSC_TRUE));
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

338: static PetscErrorCode MatTranspose_Nest(Mat A, MatReuse reuse, Mat *B)
339: {
340:   Mat_Nest *bA = (Mat_Nest *)A->data, *bC;
341:   Mat       C;
342:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

344:   PetscFunctionBegin;
345:   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
346:   PetscCheck(reuse != MAT_INPLACE_MATRIX || nr == nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Square nested matrix only for in-place");

348:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
349:     Mat *subs;
350:     IS  *is_row, *is_col;

352:     PetscCall(PetscCalloc1(nr * nc, &subs));
353:     PetscCall(PetscMalloc2(nr, &is_row, nc, &is_col));
354:     PetscCall(MatNestGetISs(A, is_row, is_col));
355:     if (reuse == MAT_INPLACE_MATRIX) {
356:       for (i = 0; i < nr; i++) {
357:         for (j = 0; j < nc; j++) subs[i + nr * j] = bA->m[i][j];
358:       }
359:     }

361:     PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nc, is_col, nr, is_row, subs, &C));
362:     PetscCall(PetscFree(subs));
363:     PetscCall(PetscFree2(is_row, is_col));
364:   } else {
365:     C = *B;
366:   }

368:   bC = (Mat_Nest *)C->data;
369:   for (i = 0; i < nr; i++) {
370:     for (j = 0; j < nc; j++) {
371:       if (bA->m[i][j]) {
372:         PetscCall(MatTranspose(bA->m[i][j], reuse, &bC->m[j][i]));
373:       } else {
374:         bC->m[j][i] = NULL;
375:       }
376:     }
377:   }

379:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
380:     *B = C;
381:   } else {
382:     PetscCall(MatHeaderMerge(A, &C));
383:   }
384:   PetscFunctionReturn(PETSC_SUCCESS);
385: }

387: static PetscErrorCode MatNestDestroyISList(PetscInt n, IS **list)
388: {
389:   IS      *lst = *list;
390:   PetscInt i;

392:   PetscFunctionBegin;
393:   if (!lst) PetscFunctionReturn(PETSC_SUCCESS);
394:   for (i = 0; i < n; i++)
395:     if (lst[i]) PetscCall(ISDestroy(&lst[i]));
396:   PetscCall(PetscFree(lst));
397:   *list = NULL;
398:   PetscFunctionReturn(PETSC_SUCCESS);
399: }

401: static PetscErrorCode MatReset_Nest(Mat A)
402: {
403:   Mat_Nest *vs = (Mat_Nest *)A->data;
404:   PetscInt  i, j;

406:   PetscFunctionBegin;
407:   /* release the matrices and the place holders */
408:   PetscCall(MatNestDestroyISList(vs->nr, &vs->isglobal.row));
409:   PetscCall(MatNestDestroyISList(vs->nc, &vs->isglobal.col));
410:   PetscCall(MatNestDestroyISList(vs->nr, &vs->islocal.row));
411:   PetscCall(MatNestDestroyISList(vs->nc, &vs->islocal.col));

413:   PetscCall(PetscFree(vs->row_len));
414:   PetscCall(PetscFree(vs->col_len));
415:   PetscCall(PetscFree(vs->nnzstate));

417:   PetscCall(PetscFree2(vs->left, vs->right));

419:   /* release the matrices and the place holders */
420:   if (vs->m) {
421:     for (i = 0; i < vs->nr; i++) {
422:       for (j = 0; j < vs->nc; j++) PetscCall(MatDestroy(&vs->m[i][j]));
423:     }
424:     PetscCall(PetscFree(vs->m[0]));
425:     PetscCall(PetscFree(vs->m));
426:   }

428:   /* restore defaults */
429:   vs->nr            = 0;
430:   vs->nc            = 0;
431:   vs->splitassembly = PETSC_FALSE;
432:   PetscFunctionReturn(PETSC_SUCCESS);
433: }

435: static PetscErrorCode MatDestroy_Nest(Mat A)
436: {
437:   PetscFunctionBegin;
438:   PetscCall(MatReset_Nest(A));
439:   PetscCall(PetscFree(A->data));
440:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMat_C", NULL));
441:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMat_C", NULL));
442:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMats_C", NULL));
443:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSize_C", NULL));
444:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetISs_C", NULL));
445:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetLocalISs_C", NULL));
446:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetVecType_C", NULL));
447:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMats_C", NULL));
448:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpiaij_C", NULL));
449:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqaij_C", NULL));
450:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_aij_C", NULL));
451:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_is_C", NULL));
452:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpidense_C", NULL));
453:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqdense_C", NULL));
454:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_seqdense_C", NULL));
455:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_mpidense_C", NULL));
456:   PetscFunctionReturn(PETSC_SUCCESS);
457: }

459: static PetscErrorCode MatAssemblyBegin_Nest(Mat A, MatAssemblyType type)
460: {
461:   Mat_Nest *vs = (Mat_Nest *)A->data;
462:   PetscInt  i, j;
463:   PetscBool nnzstate = PETSC_FALSE;

465:   PetscFunctionBegin;
466:   for (i = 0; i < vs->nr; i++) {
467:     for (j = 0; j < vs->nc; j++) {
468:       PetscObjectState subnnzstate = 0;
469:       if (vs->m[i][j]) {
470:         PetscCall(MatAssemblyBegin(vs->m[i][j], type));
471:         if (!vs->splitassembly) {
472:           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
473:            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
474:            * already performing an assembly, but the result would by more complicated and appears to offer less
475:            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
476:            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
477:            */
478:           PetscCall(MatAssemblyEnd(vs->m[i][j], type));
479:           PetscCall(MatGetNonzeroState(vs->m[i][j], &subnnzstate));
480:         }
481:       }
482:       nnzstate                     = (PetscBool)(nnzstate || vs->nnzstate[i * vs->nc + j] != subnnzstate);
483:       vs->nnzstate[i * vs->nc + j] = subnnzstate;
484:     }
485:   }
486:   if (nnzstate) A->nonzerostate++;
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
491: {
492:   Mat_Nest *vs = (Mat_Nest *)A->data;
493:   PetscInt  i, j;

495:   PetscFunctionBegin;
496:   for (i = 0; i < vs->nr; i++) {
497:     for (j = 0; j < vs->nc; j++) {
498:       if (vs->m[i][j]) {
499:         if (vs->splitassembly) PetscCall(MatAssemblyEnd(vs->m[i][j], type));
500:       }
501:     }
502:   }
503:   PetscFunctionReturn(PETSC_SUCCESS);
504: }

506: static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A, PetscInt row, Mat *B)
507: {
508:   Mat_Nest *vs = (Mat_Nest *)A->data;
509:   PetscInt  j;
510:   Mat       sub;

512:   PetscFunctionBegin;
513:   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
514:   for (j = 0; !sub && j < vs->nc; j++) sub = vs->m[row][j];
515:   if (sub) PetscCall(MatSetUp(sub)); /* Ensure that the sizes are available */
516:   *B = sub;
517:   PetscFunctionReturn(PETSC_SUCCESS);
518: }

520: static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A, PetscInt col, Mat *B)
521: {
522:   Mat_Nest *vs = (Mat_Nest *)A->data;
523:   PetscInt  i;
524:   Mat       sub;

526:   PetscFunctionBegin;
527:   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
528:   for (i = 0; !sub && i < vs->nr; i++) sub = vs->m[i][col];
529:   if (sub) PetscCall(MatSetUp(sub)); /* Ensure that the sizes are available */
530:   *B = sub;
531:   PetscFunctionReturn(PETSC_SUCCESS);
532: }

534: static PetscErrorCode MatNestFindISRange(Mat A, PetscInt n, const IS list[], IS is, PetscInt *begin, PetscInt *end)
535: {
536:   PetscInt  i, j, size, m;
537:   PetscBool flg;
538:   IS        out, concatenate[2];

540:   PetscFunctionBegin;
541:   PetscAssertPointer(list, 3);
543:   if (begin) {
544:     PetscAssertPointer(begin, 5);
545:     *begin = -1;
546:   }
547:   if (end) {
548:     PetscAssertPointer(end, 6);
549:     *end = -1;
550:   }
551:   for (i = 0; i < n; i++) {
552:     if (!list[i]) continue;
553:     PetscCall(ISEqualUnsorted(list[i], is, &flg));
554:     if (flg) {
555:       if (begin) *begin = i;
556:       if (end) *end = i + 1;
557:       PetscFunctionReturn(PETSC_SUCCESS);
558:     }
559:   }
560:   PetscCall(ISGetSize(is, &size));
561:   for (i = 0; i < n - 1; i++) {
562:     if (!list[i]) continue;
563:     m = 0;
564:     PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), 2, list + i, &out));
565:     PetscCall(ISGetSize(out, &m));
566:     for (j = i + 2; j < n && m < size; j++) {
567:       if (list[j]) {
568:         concatenate[0] = out;
569:         concatenate[1] = list[j];
570:         PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), 2, concatenate, &out));
571:         PetscCall(ISDestroy(concatenate));
572:         PetscCall(ISGetSize(out, &m));
573:       }
574:     }
575:     if (m == size) {
576:       PetscCall(ISEqualUnsorted(out, is, &flg));
577:       if (flg) {
578:         if (begin) *begin = i;
579:         if (end) *end = j;
580:         PetscCall(ISDestroy(&out));
581:         PetscFunctionReturn(PETSC_SUCCESS);
582:       }
583:     }
584:     PetscCall(ISDestroy(&out));
585:   }
586:   PetscFunctionReturn(PETSC_SUCCESS);
587: }

589: static PetscErrorCode MatNestFillEmptyMat_Private(Mat A, PetscInt i, PetscInt j, Mat *B)
590: {
591:   Mat_Nest *vs = (Mat_Nest *)A->data;
592:   PetscInt  lr, lc;

594:   PetscFunctionBegin;
595:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
596:   PetscCall(ISGetLocalSize(vs->isglobal.row[i], &lr));
597:   PetscCall(ISGetLocalSize(vs->isglobal.col[j], &lc));
598:   PetscCall(MatSetSizes(*B, lr, lc, PETSC_DECIDE, PETSC_DECIDE));
599:   PetscCall(MatSetType(*B, MATAIJ));
600:   PetscCall(MatSeqAIJSetPreallocation(*B, 0, NULL));
601:   PetscCall(MatMPIAIJSetPreallocation(*B, 0, NULL, 0, NULL));
602:   PetscCall(MatSetUp(*B));
603:   PetscCall(MatSetOption(*B, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
604:   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
605:   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
606:   PetscFunctionReturn(PETSC_SUCCESS);
607: }

609: static PetscErrorCode MatNestGetBlock_Private(Mat A, PetscInt rbegin, PetscInt rend, PetscInt cbegin, PetscInt cend, Mat *B)
610: {
611:   Mat_Nest  *vs = (Mat_Nest *)A->data;
612:   Mat       *a;
613:   PetscInt   i, j, k, l, nr = rend - rbegin, nc = cend - cbegin;
614:   char       keyname[256];
615:   PetscBool *b;
616:   PetscBool  flg;

618:   PetscFunctionBegin;
619:   *B = NULL;
620:   PetscCall(PetscSNPrintf(keyname, sizeof(keyname), "NestBlock_%" PetscInt_FMT "-%" PetscInt_FMT "x%" PetscInt_FMT "-%" PetscInt_FMT, rbegin, rend, cbegin, cend));
621:   PetscCall(PetscObjectQuery((PetscObject)A, keyname, (PetscObject *)B));
622:   if (*B) PetscFunctionReturn(PETSC_SUCCESS);

624:   PetscCall(PetscMalloc2(nr * nc, &a, nr * nc, &b));
625:   for (i = 0; i < nr; i++) {
626:     for (j = 0; j < nc; j++) {
627:       a[i * nc + j] = vs->m[rbegin + i][cbegin + j];
628:       b[i * nc + j] = PETSC_FALSE;
629:     }
630:   }
631:   if (nc != vs->nc && nr != vs->nr) {
632:     for (i = 0; i < nr; i++) {
633:       for (j = 0; j < nc; j++) {
634:         flg = PETSC_FALSE;
635:         for (k = 0; (k < nr && !flg); k++) {
636:           if (a[j + k * nc]) flg = PETSC_TRUE;
637:         }
638:         if (flg) {
639:           flg = PETSC_FALSE;
640:           for (l = 0; (l < nc && !flg); l++) {
641:             if (a[i * nc + l]) flg = PETSC_TRUE;
642:           }
643:         }
644:         if (!flg) {
645:           b[i * nc + j] = PETSC_TRUE;
646:           PetscCall(MatNestFillEmptyMat_Private(A, rbegin + i, cbegin + j, a + i * nc + j));
647:         }
648:       }
649:     }
650:   }
651:   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nr, nr != vs->nr ? NULL : vs->isglobal.row, nc, nc != vs->nc ? NULL : vs->isglobal.col, a, B));
652:   for (i = 0; i < nr; i++) {
653:     for (j = 0; j < nc; j++) {
654:       if (b[i * nc + j]) PetscCall(MatDestroy(a + i * nc + j));
655:     }
656:   }
657:   PetscCall(PetscFree2(a, b));
658:   (*B)->assembled = A->assembled;
659:   PetscCall(PetscObjectCompose((PetscObject)A, keyname, (PetscObject)*B));
660:   PetscCall(PetscObjectDereference((PetscObject)*B)); /* Leave the only remaining reference in the composition */
661:   PetscFunctionReturn(PETSC_SUCCESS);
662: }

664: static PetscErrorCode MatNestFindSubMat(Mat A, struct MatNestISPair *is, IS isrow, IS iscol, Mat *B)
665: {
666:   Mat_Nest *vs = (Mat_Nest *)A->data;
667:   PetscInt  rbegin, rend, cbegin, cend;

669:   PetscFunctionBegin;
670:   PetscCall(MatNestFindISRange(A, vs->nr, is->row, isrow, &rbegin, &rend));
671:   PetscCall(MatNestFindISRange(A, vs->nc, is->col, iscol, &cbegin, &cend));
672:   if (rend == rbegin + 1 && cend == cbegin + 1) {
673:     if (!vs->m[rbegin][cbegin]) PetscCall(MatNestFillEmptyMat_Private(A, rbegin, cbegin, vs->m[rbegin] + cbegin));
674:     *B = vs->m[rbegin][cbegin];
675:   } else if (rbegin != -1 && cbegin != -1) {
676:     PetscCall(MatNestGetBlock_Private(A, rbegin, rend, cbegin, cend, B));
677:   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Could not find index set");
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: /*
682:    TODO: This does not actually returns a submatrix we can modify
683: */
684: static PetscErrorCode MatCreateSubMatrix_Nest(Mat A, IS isrow, IS iscol, MatReuse reuse, Mat *B)
685: {
686:   Mat_Nest *vs = (Mat_Nest *)A->data;
687:   Mat       sub;

689:   PetscFunctionBegin;
690:   PetscCall(MatNestFindSubMat(A, &vs->isglobal, isrow, iscol, &sub));
691:   switch (reuse) {
692:   case MAT_INITIAL_MATRIX:
693:     PetscCall(PetscObjectReference((PetscObject)sub));
694:     if (sub) PetscCall(PetscObjectStateIncrease((PetscObject)sub));
695:     *B = sub;
696:     break;
697:   case MAT_REUSE_MATRIX:
698:     PetscCheck(sub == *B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Submatrix was not used before in this call");
699:     if (sub) PetscCall(PetscObjectStateIncrease((PetscObject)sub));
700:     break;
701:   default:
702:     break;
703:   }
704:   PetscFunctionReturn(PETSC_SUCCESS);
705: }

707: static PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A, IS isrow, IS iscol, Mat *B)
708: {
709:   Mat_Nest *vs = (Mat_Nest *)A->data;
710:   Mat       sub;

712:   PetscFunctionBegin;
713:   PetscCall(MatNestFindSubMat(A, &vs->islocal, isrow, iscol, &sub));
714:   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
715:   if (sub) PetscCall(PetscObjectReference((PetscObject)sub));
716:   *B = sub;
717:   PetscFunctionReturn(PETSC_SUCCESS);
718: }

720: static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A, IS isrow, IS iscol, Mat *B)
721: {
722:   Mat_Nest *vs = (Mat_Nest *)A->data;
723:   Mat       sub;

725:   PetscFunctionBegin;
726:   PetscCall(MatNestFindSubMat(A, &vs->islocal, isrow, iscol, &sub));
727:   PetscCheck(*B == sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Local submatrix has not been gotten");
728:   if (sub) {
729:     PetscCheck(((PetscObject)sub)->refct > 1, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Local submatrix has had reference count decremented too many times");
730:     PetscCall(MatDestroy(B));
731:   }
732:   PetscFunctionReturn(PETSC_SUCCESS);
733: }

735: static PetscErrorCode MatGetDiagonal_Nest(Mat A, Vec v)
736: {
737:   Mat_Nest *bA = (Mat_Nest *)A->data;
738:   PetscInt  i;

740:   PetscFunctionBegin;
741:   for (i = 0; i < bA->nr; i++) {
742:     Vec bv;
743:     PetscCall(VecGetSubVector(v, bA->isglobal.row[i], &bv));
744:     if (bA->m[i][i]) PetscCall(MatGetDiagonal(bA->m[i][i], bv));
745:     else PetscCall(VecSet(bv, 0.0));
746:     PetscCall(VecRestoreSubVector(v, bA->isglobal.row[i], &bv));
747:   }
748:   PetscFunctionReturn(PETSC_SUCCESS);
749: }

751: static PetscErrorCode MatDiagonalScale_Nest(Mat A, Vec l, Vec r)
752: {
753:   Mat_Nest *bA = (Mat_Nest *)A->data;
754:   Vec       bl, *br;
755:   PetscInt  i, j;

757:   PetscFunctionBegin;
758:   PetscCall(PetscCalloc1(bA->nc, &br));
759:   if (r) {
760:     for (j = 0; j < bA->nc; j++) PetscCall(VecGetSubVector(r, bA->isglobal.col[j], &br[j]));
761:   }
762:   bl = NULL;
763:   for (i = 0; i < bA->nr; i++) {
764:     if (l) PetscCall(VecGetSubVector(l, bA->isglobal.row[i], &bl));
765:     for (j = 0; j < bA->nc; j++) {
766:       if (bA->m[i][j]) PetscCall(MatDiagonalScale(bA->m[i][j], bl, br[j]));
767:     }
768:     if (l) PetscCall(VecRestoreSubVector(l, bA->isglobal.row[i], &bl));
769:   }
770:   if (r) {
771:     for (j = 0; j < bA->nc; j++) PetscCall(VecRestoreSubVector(r, bA->isglobal.col[j], &br[j]));
772:   }
773:   PetscCall(PetscFree(br));
774:   PetscFunctionReturn(PETSC_SUCCESS);
775: }

777: static PetscErrorCode MatScale_Nest(Mat A, PetscScalar a)
778: {
779:   Mat_Nest *bA = (Mat_Nest *)A->data;
780:   PetscInt  i, j;

782:   PetscFunctionBegin;
783:   for (i = 0; i < bA->nr; i++) {
784:     for (j = 0; j < bA->nc; j++) {
785:       if (bA->m[i][j]) PetscCall(MatScale(bA->m[i][j], a));
786:     }
787:   }
788:   PetscFunctionReturn(PETSC_SUCCESS);
789: }

791: static PetscErrorCode MatShift_Nest(Mat A, PetscScalar a)
792: {
793:   Mat_Nest *bA = (Mat_Nest *)A->data;
794:   PetscInt  i;
795:   PetscBool nnzstate = PETSC_FALSE;

797:   PetscFunctionBegin;
798:   for (i = 0; i < bA->nr; i++) {
799:     PetscObjectState subnnzstate = 0;
800:     PetscCheck(bA->m[i][i], PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for shifting an empty diagonal block, insert a matrix in block (%" PetscInt_FMT ",%" PetscInt_FMT ")", i, i);
801:     PetscCall(MatShift(bA->m[i][i], a));
802:     PetscCall(MatGetNonzeroState(bA->m[i][i], &subnnzstate));
803:     nnzstate                     = (PetscBool)(nnzstate || bA->nnzstate[i * bA->nc + i] != subnnzstate);
804:     bA->nnzstate[i * bA->nc + i] = subnnzstate;
805:   }
806:   if (nnzstate) A->nonzerostate++;
807:   PetscFunctionReturn(PETSC_SUCCESS);
808: }

810: static PetscErrorCode MatDiagonalSet_Nest(Mat A, Vec D, InsertMode is)
811: {
812:   Mat_Nest *bA = (Mat_Nest *)A->data;
813:   PetscInt  i;
814:   PetscBool nnzstate = PETSC_FALSE;

816:   PetscFunctionBegin;
817:   for (i = 0; i < bA->nr; i++) {
818:     PetscObjectState subnnzstate = 0;
819:     Vec              bv;
820:     PetscCall(VecGetSubVector(D, bA->isglobal.row[i], &bv));
821:     if (bA->m[i][i]) {
822:       PetscCall(MatDiagonalSet(bA->m[i][i], bv, is));
823:       PetscCall(MatGetNonzeroState(bA->m[i][i], &subnnzstate));
824:     }
825:     PetscCall(VecRestoreSubVector(D, bA->isglobal.row[i], &bv));
826:     nnzstate                     = (PetscBool)(nnzstate || bA->nnzstate[i * bA->nc + i] != subnnzstate);
827:     bA->nnzstate[i * bA->nc + i] = subnnzstate;
828:   }
829:   if (nnzstate) A->nonzerostate++;
830:   PetscFunctionReturn(PETSC_SUCCESS);
831: }

833: static PetscErrorCode MatSetRandom_Nest(Mat A, PetscRandom rctx)
834: {
835:   Mat_Nest *bA = (Mat_Nest *)A->data;
836:   PetscInt  i, j;

838:   PetscFunctionBegin;
839:   for (i = 0; i < bA->nr; i++) {
840:     for (j = 0; j < bA->nc; j++) {
841:       if (bA->m[i][j]) PetscCall(MatSetRandom(bA->m[i][j], rctx));
842:     }
843:   }
844:   PetscFunctionReturn(PETSC_SUCCESS);
845: }

847: static PetscErrorCode MatCreateVecs_Nest(Mat A, Vec *right, Vec *left)
848: {
849:   Mat_Nest *bA = (Mat_Nest *)A->data;
850:   Vec      *L, *R;
851:   MPI_Comm  comm;
852:   PetscInt  i, j;

854:   PetscFunctionBegin;
855:   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
856:   if (right) {
857:     /* allocate R */
858:     PetscCall(PetscMalloc1(bA->nc, &R));
859:     /* Create the right vectors */
860:     for (j = 0; j < bA->nc; j++) {
861:       for (i = 0; i < bA->nr; i++) {
862:         if (bA->m[i][j]) {
863:           PetscCall(MatCreateVecs(bA->m[i][j], &R[j], NULL));
864:           break;
865:         }
866:       }
867:       PetscCheck(i != bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
868:     }
869:     PetscCall(VecCreateNest(comm, bA->nc, bA->isglobal.col, R, right));
870:     /* hand back control to the nest vector */
871:     for (j = 0; j < bA->nc; j++) PetscCall(VecDestroy(&R[j]));
872:     PetscCall(PetscFree(R));
873:   }

875:   if (left) {
876:     /* allocate L */
877:     PetscCall(PetscMalloc1(bA->nr, &L));
878:     /* Create the left vectors */
879:     for (i = 0; i < bA->nr; i++) {
880:       for (j = 0; j < bA->nc; j++) {
881:         if (bA->m[i][j]) {
882:           PetscCall(MatCreateVecs(bA->m[i][j], NULL, &L[i]));
883:           break;
884:         }
885:       }
886:       PetscCheck(j != bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
887:     }

889:     PetscCall(VecCreateNest(comm, bA->nr, bA->isglobal.row, L, left));
890:     for (i = 0; i < bA->nr; i++) PetscCall(VecDestroy(&L[i]));

892:     PetscCall(PetscFree(L));
893:   }
894:   PetscFunctionReturn(PETSC_SUCCESS);
895: }

897: static PetscErrorCode MatView_Nest(Mat A, PetscViewer viewer)
898: {
899:   Mat_Nest *bA = (Mat_Nest *)A->data;
900:   PetscBool isascii, viewSub = PETSC_FALSE;
901:   PetscInt  i, j;

903:   PetscFunctionBegin;
904:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
905:   if (isascii) {
906:     PetscViewerFormat format;

908:     PetscCall(PetscViewerGetFormat(viewer, &format));
909:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
910:       Mat T;

912:       PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &T));
913:       PetscCall(MatView(T, viewer));
914:       PetscCall(MatDestroy(&T));
915:       PetscFunctionReturn(PETSC_SUCCESS);
916:     }
917:     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_view_nest_sub", &viewSub, NULL));
918:     PetscCall(PetscViewerASCIIPushTab(viewer));
919:     PetscCall(PetscViewerASCIIPrintf(viewer, "MatNest, rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT ", structure:\n", bA->nr, bA->nc));
920:     for (i = 0; i < bA->nr; i++) {
921:       for (j = 0; j < bA->nc; j++) {
922:         MatType   type;
923:         char      name[256] = "", prefix[256] = "";
924:         PetscInt  NR, NC;
925:         PetscBool isNest = PETSC_FALSE;

927:         if (!bA->m[i][j]) {
928:           PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ",%" PetscInt_FMT ") : NULL\n", i, j));
929:           continue;
930:         }
931:         PetscCall(MatGetSize(bA->m[i][j], &NR, &NC));
932:         PetscCall(MatGetType(bA->m[i][j], &type));
933:         if (((PetscObject)bA->m[i][j])->name) PetscCall(PetscSNPrintf(name, sizeof(name), "name=\"%s\", ", ((PetscObject)bA->m[i][j])->name));
934:         if (((PetscObject)bA->m[i][j])->prefix) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "prefix=\"%s\", ", ((PetscObject)bA->m[i][j])->prefix));
935:         PetscCall(PetscObjectTypeCompare((PetscObject)bA->m[i][j], MATNEST, &isNest));

937:         PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ",%" PetscInt_FMT ") : %s%stype=%s, rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT "\n", i, j, name, prefix, type, NR, NC));

939:         if (isNest || viewSub) {
940:           PetscCall(PetscViewerASCIIPushTab(viewer)); /* push1 */
941:           PetscCall(MatView(bA->m[i][j], viewer));
942:           PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop1 */
943:         }
944:       }
945:     }
946:     PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop0 */
947:   }
948:   PetscFunctionReturn(PETSC_SUCCESS);
949: }

951: static PetscErrorCode MatZeroEntries_Nest(Mat A)
952: {
953:   Mat_Nest *bA = (Mat_Nest *)A->data;
954:   PetscInt  i, j;

956:   PetscFunctionBegin;
957:   for (i = 0; i < bA->nr; i++) {
958:     for (j = 0; j < bA->nc; j++) {
959:       if (!bA->m[i][j]) continue;
960:       PetscCall(MatZeroEntries(bA->m[i][j]));
961:     }
962:   }
963:   PetscFunctionReturn(PETSC_SUCCESS);
964: }

966: static PetscErrorCode MatCopy_Nest(Mat A, Mat B, MatStructure str)
967: {
968:   Mat_Nest *bA = (Mat_Nest *)A->data, *bB = (Mat_Nest *)B->data;
969:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
970:   PetscBool nnzstate = PETSC_FALSE;

972:   PetscFunctionBegin;
973:   PetscCheck(nr == bB->nr && nc == bB->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Cannot copy a Mat_Nest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ") to a Mat_Nest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ")", bB->nr, bB->nc, nr, nc);
974:   for (i = 0; i < nr; i++) {
975:     for (j = 0; j < nc; j++) {
976:       PetscObjectState subnnzstate = 0;
977:       if (bA->m[i][j] && bB->m[i][j]) {
978:         PetscCall(MatCopy(bA->m[i][j], bB->m[i][j], str));
979:         PetscCall(MatGetNonzeroState(bB->m[i][j], &subnnzstate));
980:         nnzstate                 = (PetscBool)(nnzstate || bB->nnzstate[i * nc + j] != subnnzstate);
981:         bB->nnzstate[i * nc + j] = subnnzstate;
982:       } else if (bA->m[i][j]) { // bB->m[i][j] is NULL
983:         Mat M;

985:         PetscCheck(str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Matrix block does not exist at %" PetscInt_FMT ",%" PetscInt_FMT ". Use DIFFERENT_NONZERO_PATTERN or UNKNOWN_NONZERO_PATTERN", i, j);
986:         PetscCall(MatDuplicate(bA->m[i][j], MAT_COPY_VALUES, &M));
987:         PetscCall(MatNestSetSubMat(B, i, j, M));
988:         PetscCall(MatDestroy(&M));
989:       } else if (bB->m[i][j]) { // bA->m[i][j] is NULL
990:         PetscCheck(str == DIFFERENT_NONZERO_PATTERN || str == SUBSET_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Matrix block does not exist at %" PetscInt_FMT ",%" PetscInt_FMT ". Use DIFFERENT_NONZERO_PATTERN, SUBSET_NONZERO_PATTERN or UNKNOWN_NONZERO_PATTERN", i, j);
991:         PetscCall(MatNestSetSubMat(B, i, j, NULL));
992:       }
993:     }
994:   }
995:   if (nnzstate) B->nonzerostate++;
996:   PetscFunctionReturn(PETSC_SUCCESS);
997: }

999: static PetscErrorCode MatAXPY_Nest(Mat Y, PetscScalar a, Mat X, MatStructure str)
1000: {
1001:   Mat_Nest *bY = (Mat_Nest *)Y->data, *bX = (Mat_Nest *)X->data;
1002:   PetscInt  i, j, nr = bY->nr, nc = bY->nc;
1003:   PetscBool nnzstate = PETSC_FALSE;

1005:   PetscFunctionBegin;
1006:   PetscCheck(nr == bX->nr && nc == bX->nc, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_INCOMP, "Cannot AXPY a MatNest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ") with a MatNest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ")", bX->nr, bX->nc, nr, nc);
1007:   for (i = 0; i < nr; i++) {
1008:     for (j = 0; j < nc; j++) {
1009:       PetscObjectState subnnzstate = 0;
1010:       if (bY->m[i][j] && bX->m[i][j]) {
1011:         PetscCall(MatAXPY(bY->m[i][j], a, bX->m[i][j], str));
1012:       } else if (bX->m[i][j]) {
1013:         Mat M;

1015:         PetscCheck(str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_INCOMP, "Matrix block does not exist at %" PetscInt_FMT ",%" PetscInt_FMT ". Use DIFFERENT_NONZERO_PATTERN or UNKNOWN_NONZERO_PATTERN", i, j);
1016:         PetscCall(MatDuplicate(bX->m[i][j], MAT_COPY_VALUES, &M));
1017:         PetscCall(MatScale(M, a));
1018:         PetscCall(MatNestSetSubMat(Y, i, j, M));
1019:         PetscCall(MatDestroy(&M));
1020:       }
1021:       if (bY->m[i][j]) PetscCall(MatGetNonzeroState(bY->m[i][j], &subnnzstate));
1022:       nnzstate                 = (PetscBool)(nnzstate || bY->nnzstate[i * nc + j] != subnnzstate);
1023:       bY->nnzstate[i * nc + j] = subnnzstate;
1024:     }
1025:   }
1026:   if (nnzstate) Y->nonzerostate++;
1027:   PetscFunctionReturn(PETSC_SUCCESS);
1028: }

1030: static PetscErrorCode MatDuplicate_Nest(Mat A, MatDuplicateOption op, Mat *B)
1031: {
1032:   Mat_Nest *bA = (Mat_Nest *)A->data;
1033:   Mat      *b;
1034:   PetscInt  i, j, nr = bA->nr, nc = bA->nc;

1036:   PetscFunctionBegin;
1037:   PetscCall(PetscMalloc1(nr * nc, &b));
1038:   for (i = 0; i < nr; i++) {
1039:     for (j = 0; j < nc; j++) {
1040:       if (bA->m[i][j]) PetscCall(MatDuplicate(bA->m[i][j], op, &b[i * nc + j]));
1041:       else b[i * nc + j] = NULL;
1042:     }
1043:   }
1044:   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nr, bA->isglobal.row, nc, bA->isglobal.col, b, B));
1045:   /* Give the new MatNest exclusive ownership */
1046:   for (i = 0; i < nr * nc; i++) PetscCall(MatDestroy(&b[i]));
1047:   PetscCall(PetscFree(b));

1049:   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
1050:   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
1051:   PetscFunctionReturn(PETSC_SUCCESS);
1052: }

1054: /* nest api */
1055: static PetscErrorCode MatNestGetSubMat_Nest(Mat A, PetscInt idxm, PetscInt jdxm, Mat *mat)
1056: {
1057:   Mat_Nest *bA = (Mat_Nest *)A->data;

1059:   PetscFunctionBegin;
1060:   PetscCheck(idxm < bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm, bA->nr - 1);
1061:   PetscCheck(jdxm < bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Col too large: row %" PetscInt_FMT " max %" PetscInt_FMT, jdxm, bA->nc - 1);
1062:   *mat = bA->m[idxm][jdxm];
1063:   PetscFunctionReturn(PETSC_SUCCESS);
1064: }

1066: /*@
1067:   MatNestGetSubMat - Returns a single, sub-matrix from a `MATNEST`

1069:   Not Collective

1071:   Input Parameters:
1072: + A    - `MATNEST` matrix
1073: . idxm - index of the matrix within the nest matrix
1074: - jdxm - index of the matrix within the nest matrix

1076:   Output Parameter:
1077: . sub - matrix at index `idxm`, `jdxm` within the nest matrix

1079:   Level: developer

1081: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSize()`, `MatNestGetSubMats()`, `MatCreateNest()`, `MatNestSetSubMat()`,
1082:           `MatNestGetLocalISs()`, `MatNestGetISs()`
1083: @*/
1084: PetscErrorCode MatNestGetSubMat(Mat A, PetscInt idxm, PetscInt jdxm, Mat *sub)
1085: {
1086:   PetscFunctionBegin;
1090:   PetscAssertPointer(sub, 4);
1091:   PetscUseMethod(A, "MatNestGetSubMat_C", (Mat, PetscInt, PetscInt, Mat *), (A, idxm, jdxm, sub));
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: static PetscErrorCode MatNestSetSubMat_Nest(Mat A, PetscInt idxm, PetscInt jdxm, Mat mat)
1096: {
1097:   Mat_Nest *bA = (Mat_Nest *)A->data;
1098:   PetscInt  m, n, M, N, mi, ni, Mi, Ni;

1100:   PetscFunctionBegin;
1101:   PetscCheck(idxm < bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm, bA->nr - 1);
1102:   PetscCheck(jdxm < bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Col too large: row %" PetscInt_FMT " max %" PetscInt_FMT, jdxm, bA->nc - 1);
1103:   if (mat) {
1104:     PetscCall(MatGetLocalSize(mat, &m, &n));
1105:     PetscCall(MatGetSize(mat, &M, &N));
1106:     PetscCall(ISGetLocalSize(bA->isglobal.row[idxm], &mi));
1107:     PetscCall(ISGetSize(bA->isglobal.row[idxm], &Mi));
1108:     PetscCall(ISGetLocalSize(bA->isglobal.col[jdxm], &ni));
1109:     PetscCall(ISGetSize(bA->isglobal.col[jdxm], &Ni));
1110:     PetscCheck(M == Mi && N == Ni, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "Submatrix dimension (%" PetscInt_FMT ",%" PetscInt_FMT ") incompatible with nest block (%" PetscInt_FMT ",%" PetscInt_FMT ")", M, N, Mi, Ni);
1111:     PetscCheck(m == mi && n == ni, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "Submatrix local dimension (%" PetscInt_FMT ",%" PetscInt_FMT ") incompatible with nest block (%" PetscInt_FMT ",%" PetscInt_FMT ")", m, n, mi, ni);
1112:   }

1114:   /* do not increase object state */
1115:   if (mat == bA->m[idxm][jdxm]) PetscFunctionReturn(PETSC_SUCCESS);

1117:   PetscCall(PetscObjectReference((PetscObject)mat));
1118:   PetscCall(MatDestroy(&bA->m[idxm][jdxm]));
1119:   bA->m[idxm][jdxm] = mat;
1120:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
1121:   if (mat) PetscCall(MatGetNonzeroState(mat, &bA->nnzstate[idxm * bA->nc + jdxm]));
1122:   else bA->nnzstate[idxm * bA->nc + jdxm] = 0;
1123:   A->nonzerostate++;
1124:   PetscFunctionReturn(PETSC_SUCCESS);
1125: }

1127: /*@
1128:   MatNestSetSubMat - Set a single submatrix in the `MATNEST`

1130:   Logically Collective

1132:   Input Parameters:
1133: + A    - `MATNEST` matrix
1134: . idxm - index of the matrix within the nest matrix
1135: . jdxm - index of the matrix within the nest matrix
1136: - sub  - matrix at index `idxm`, `jdxm` within the nest matrix

1138:   Level: developer

1140:   Notes:
1141:   The new submatrix must have the same size and communicator as that block of the nest.

1143:   This increments the reference count of the submatrix.

1145: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestSetSubMats()`, `MatNestGetSubMats()`, `MatNestGetLocalISs()`, `MatCreateNest()`,
1146:           `MatNestGetSubMat()`, `MatNestGetISs()`, `MatNestGetSize()`
1147: @*/
1148: PetscErrorCode MatNestSetSubMat(Mat A, PetscInt idxm, PetscInt jdxm, Mat sub)
1149: {
1150:   PetscFunctionBegin;
1155:   PetscTryMethod(A, "MatNestSetSubMat_C", (Mat, PetscInt, PetscInt, Mat), (A, idxm, jdxm, sub));
1156:   PetscFunctionReturn(PETSC_SUCCESS);
1157: }

1159: static PetscErrorCode MatNestGetSubMats_Nest(Mat A, PetscInt *M, PetscInt *N, Mat ***mat)
1160: {
1161:   Mat_Nest *bA = (Mat_Nest *)A->data;

1163:   PetscFunctionBegin;
1164:   if (M) *M = bA->nr;
1165:   if (N) *N = bA->nc;
1166:   if (mat) *mat = bA->m;
1167:   PetscFunctionReturn(PETSC_SUCCESS);
1168: }

1170: /*@C
1171:   MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a `MATNEST` matrix.

1173:   Not Collective

1175:   Input Parameter:
1176: . A - nest matrix

1178:   Output Parameters:
1179: + M   - number of submatrix rows in the nest matrix
1180: . N   - number of submatrix columns in the nest matrix
1181: - mat - array of matrices

1183:   Level: developer

1185:   Note:
1186:   The user should not free the array `mat`.

1188:   Fortran Notes:
1189:   This routine has a calling sequence `call MatNestGetSubMats(A, M, N, mat, ierr)`
1190:   where the space allocated for the optional argument `mat` is assumed large enough (if provided).
1191:   Matrices in `mat` are returned in row-major order, see `MatCreateNest()` for an example.

1193: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSize()`, `MatNestGetSubMat()`, `MatNestGetLocalISs()`, `MatCreateNest()`,
1194:           `MatNestSetSubMats()`, `MatNestGetISs()`, `MatNestSetSubMat()`
1195: @*/
1196: PetscErrorCode MatNestGetSubMats(Mat A, PetscInt *M, PetscInt *N, Mat ***mat)
1197: {
1198:   PetscFunctionBegin;
1200:   PetscUseMethod(A, "MatNestGetSubMats_C", (Mat, PetscInt *, PetscInt *, Mat ***), (A, M, N, mat));
1201:   PetscFunctionReturn(PETSC_SUCCESS);
1202: }

1204: static PetscErrorCode MatNestGetSize_Nest(Mat A, PetscInt *M, PetscInt *N)
1205: {
1206:   Mat_Nest *bA = (Mat_Nest *)A->data;

1208:   PetscFunctionBegin;
1209:   if (M) *M = bA->nr;
1210:   if (N) *N = bA->nc;
1211:   PetscFunctionReturn(PETSC_SUCCESS);
1212: }

1214: /*@
1215:   MatNestGetSize - Returns the size of the `MATNEST` matrix.

1217:   Not Collective

1219:   Input Parameter:
1220: . A - `MATNEST` matrix

1222:   Output Parameters:
1223: + M - number of rows in the nested mat
1224: - N - number of cols in the nested mat

1226:   Level: developer

1228:   Note:
1229:   `size` refers to the number of submatrices in the row and column directions of the nested matrix

1231: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatCreateNest()`, `MatNestGetLocalISs()`,
1232:           `MatNestGetISs()`
1233: @*/
1234: PetscErrorCode MatNestGetSize(Mat A, PetscInt *M, PetscInt *N)
1235: {
1236:   PetscFunctionBegin;
1238:   PetscUseMethod(A, "MatNestGetSize_C", (Mat, PetscInt *, PetscInt *), (A, M, N));
1239:   PetscFunctionReturn(PETSC_SUCCESS);
1240: }

1242: static PetscErrorCode MatNestGetISs_Nest(Mat A, IS rows[], IS cols[])
1243: {
1244:   Mat_Nest *vs = (Mat_Nest *)A->data;
1245:   PetscInt  i;

1247:   PetscFunctionBegin;
1248:   if (rows)
1249:     for (i = 0; i < vs->nr; i++) rows[i] = vs->isglobal.row[i];
1250:   if (cols)
1251:     for (i = 0; i < vs->nc; i++) cols[i] = vs->isglobal.col[i];
1252:   PetscFunctionReturn(PETSC_SUCCESS);
1253: }

1255: /*@
1256:   MatNestGetISs - Returns the index sets partitioning the row and column spaces of a `MATNEST`

1258:   Not Collective

1260:   Input Parameter:
1261: . A - `MATNEST` matrix

1263:   Output Parameters:
1264: + rows - array of row index sets (pass `NULL` to ignore)
1265: - cols - array of column index sets (pass `NULL` to ignore)

1267:   Level: advanced

1269:   Note:
1270:   The user must have allocated arrays of the correct size. The reference count is not increased on the returned `IS`s.

1272: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatNestGetSize()`, `MatNestGetLocalISs()`,
1273:           `MatCreateNest()`, `MatNestSetSubMats()`
1274: @*/
1275: PetscErrorCode MatNestGetISs(Mat A, IS rows[], IS cols[])
1276: {
1277:   PetscFunctionBegin;
1279:   PetscUseMethod(A, "MatNestGetISs_C", (Mat, IS[], IS[]), (A, rows, cols));
1280:   PetscFunctionReturn(PETSC_SUCCESS);
1281: }

1283: static PetscErrorCode MatNestGetLocalISs_Nest(Mat A, IS rows[], IS cols[])
1284: {
1285:   Mat_Nest *vs = (Mat_Nest *)A->data;
1286:   PetscInt  i;

1288:   PetscFunctionBegin;
1289:   if (rows)
1290:     for (i = 0; i < vs->nr; i++) rows[i] = vs->islocal.row[i];
1291:   if (cols)
1292:     for (i = 0; i < vs->nc; i++) cols[i] = vs->islocal.col[i];
1293:   PetscFunctionReturn(PETSC_SUCCESS);
1294: }

1296: /*@
1297:   MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces of a `MATNEST`

1299:   Not Collective

1301:   Input Parameter:
1302: . A - `MATNEST` matrix

1304:   Output Parameters:
1305: + rows - array of row index sets (pass `NULL` to ignore)
1306: - cols - array of column index sets (pass `NULL` to ignore)

1308:   Level: advanced

1310:   Note:
1311:   The user must have allocated arrays of the correct size. The reference count is not increased on the returned `IS`s.

1313: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatNestGetSize()`, `MatNestGetISs()`, `MatCreateNest()`,
1314:           `MatNestSetSubMats()`, `MatNestSetSubMat()`
1315: @*/
1316: PetscErrorCode MatNestGetLocalISs(Mat A, IS rows[], IS cols[])
1317: {
1318:   PetscFunctionBegin;
1320:   PetscUseMethod(A, "MatNestGetLocalISs_C", (Mat, IS[], IS[]), (A, rows, cols));
1321:   PetscFunctionReturn(PETSC_SUCCESS);
1322: }

1324: static PetscErrorCode MatNestSetVecType_Nest(Mat A, VecType vtype)
1325: {
1326:   PetscBool flg;

1328:   PetscFunctionBegin;
1329:   PetscCall(PetscStrcmp(vtype, VECNEST, &flg));
1330:   /* In reality, this only distinguishes VECNEST and "other" */
1331:   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
1332:   else A->ops->getvecs = NULL;
1333:   PetscFunctionReturn(PETSC_SUCCESS);
1334: }

1336: /*@
1337:   MatNestSetVecType - Sets the type of `Vec` returned by `MatCreateVecs()`

1339:   Not Collective

1341:   Input Parameters:
1342: + A     - `MATNEST` matrix
1343: - vtype - `VecType` to use for creating vectors

1345:   Level: developer

1347: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreateVecs()`, `MatCreateNest()`, `VecType`
1348: @*/
1349: PetscErrorCode MatNestSetVecType(Mat A, VecType vtype)
1350: {
1351:   PetscFunctionBegin;
1353:   PetscTryMethod(A, "MatNestSetVecType_C", (Mat, VecType), (A, vtype));
1354:   PetscFunctionReturn(PETSC_SUCCESS);
1355: }

1357: static PetscErrorCode MatNestSetSubMats_Nest(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[])
1358: {
1359:   Mat_Nest *s = (Mat_Nest *)A->data;
1360:   PetscInt  i, j, m, n, M, N;
1361:   PetscBool cong, isstd, sametype = PETSC_FALSE;
1362:   VecType   vtype, type;

1364:   PetscFunctionBegin;
1365:   PetscCall(MatReset_Nest(A));

1367:   s->nr = nr;
1368:   s->nc = nc;

1370:   /* Create space for submatrices */
1371:   PetscCall(PetscMalloc1(nr, &s->m));
1372:   PetscCall(PetscMalloc1(nr * nc, &s->m[0]));
1373:   for (i = 0; i < nr; i++) {
1374:     s->m[i] = s->m[0] + i * nc;
1375:     for (j = 0; j < nc; j++) {
1376:       s->m[i][j] = a ? a[i * nc + j] : NULL;
1377:       PetscCall(PetscObjectReference((PetscObject)s->m[i][j]));
1378:     }
1379:   }
1380:   PetscCall(MatGetVecType(A, &vtype));
1381:   PetscCall(PetscStrcmp(vtype, VECSTANDARD, &isstd));
1382:   if (isstd) {
1383:     /* check if all blocks have the same vectype */
1384:     vtype = NULL;
1385:     for (i = 0; i < nr; i++) {
1386:       for (j = 0; j < nc; j++) {
1387:         if (s->m[i][j]) {
1388:           if (!vtype) { /* first visited block */
1389:             PetscCall(MatGetVecType(s->m[i][j], &vtype));
1390:             sametype = PETSC_TRUE;
1391:           } else if (sametype) {
1392:             PetscCall(MatGetVecType(s->m[i][j], &type));
1393:             PetscCall(PetscStrcmp(vtype, type, &sametype));
1394:           }
1395:         }
1396:       }
1397:     }
1398:     if (sametype) { /* propagate vectype */
1399:       PetscCall(MatSetVecType(A, vtype));
1400:     }
1401:   }

1403:   PetscCall(MatSetUp_NestIS_Private(A, nr, is_row, nc, is_col));

1405:   PetscCall(PetscMalloc1(nr, &s->row_len));
1406:   PetscCall(PetscMalloc1(nc, &s->col_len));
1407:   for (i = 0; i < nr; i++) s->row_len[i] = -1;
1408:   for (j = 0; j < nc; j++) s->col_len[j] = -1;

1410:   PetscCall(PetscCalloc1(nr * nc, &s->nnzstate));
1411:   for (i = 0; i < nr; i++) {
1412:     for (j = 0; j < nc; j++) {
1413:       if (s->m[i][j]) PetscCall(MatGetNonzeroState(s->m[i][j], &s->nnzstate[i * nc + j]));
1414:     }
1415:   }

1417:   PetscCall(MatNestGetSizes_Private(A, &m, &n, &M, &N));

1419:   PetscCall(PetscLayoutSetSize(A->rmap, M));
1420:   PetscCall(PetscLayoutSetLocalSize(A->rmap, m));
1421:   PetscCall(PetscLayoutSetSize(A->cmap, N));
1422:   PetscCall(PetscLayoutSetLocalSize(A->cmap, n));

1424:   PetscCall(PetscLayoutSetUp(A->rmap));
1425:   PetscCall(PetscLayoutSetUp(A->cmap));

1427:   /* disable operations that are not supported for non-square matrices,
1428:      or matrices for which is_row != is_col  */
1429:   PetscCall(MatHasCongruentLayouts(A, &cong));
1430:   if (cong && nr != nc) cong = PETSC_FALSE;
1431:   if (cong) {
1432:     for (i = 0; cong && i < nr; i++) PetscCall(ISEqualUnsorted(s->isglobal.row[i], s->isglobal.col[i], &cong));
1433:   }
1434:   if (!cong) {
1435:     A->ops->getdiagonal = NULL;
1436:     A->ops->shift       = NULL;
1437:     A->ops->diagonalset = NULL;
1438:   }

1440:   PetscCall(PetscCalloc2(nr, &s->left, nc, &s->right));
1441:   PetscCall(PetscObjectStateIncrease((PetscObject)A));
1442:   A->nonzerostate++;
1443:   PetscFunctionReturn(PETSC_SUCCESS);
1444: }

1446: /*@
1447:   MatNestSetSubMats - Sets the nested submatrices in a `MATNEST`

1449:   Collective

1451:   Input Parameters:
1452: + A      - `MATNEST` matrix
1453: . nr     - number of nested row blocks
1454: . is_row - index sets for each nested row block, or `NULL` to make contiguous
1455: . nc     - number of nested column blocks
1456: . is_col - index sets for each nested column block, or `NULL` to make contiguous
1457: - a      - array of $ nr \times nc$ submatrices, or `NULL`

1459:   Level: advanced

1461:   Notes:
1462:   This always resets any block matrix information previously set.

1464:   Pass `NULL` in the corresponding entry of `a` for an empty block.

1466:   In both C and Fortran, `a` must be a one-dimensional array representing a two-dimensional row-major order array containing the matrices. See
1467:   `MatCreateNest()` for an example.

1469:   Fortran Note:
1470:   Pass `PETSC_NULL_MAT` in the corresponding entry of `a` for an empty block

1472: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreateNest()`, `MatNestSetSubMat()`, `MatNestGetSubMat()`, `MatNestGetSubMats()`
1473: @*/
1474: PetscErrorCode MatNestSetSubMats(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[]) PeNSS
1475: {
1476:   PetscFunctionBegin;
1479:   PetscCheck(nr >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Number of rows cannot be negative");
1480:   if (nr && is_row) {
1481:     PetscAssertPointer(is_row, 3);
1483:   }
1485:   PetscCheck(nc >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Number of columns cannot be negative");
1486:   if (nc && is_col) {
1487:     PetscAssertPointer(is_col, 5);
1489:   }
1490:   PetscTryMethod(A, "MatNestSetSubMats_C", (Mat, PetscInt, const IS[], PetscInt, const IS[], const Mat[]), (A, nr, is_row, nc, is_col, a));
1491:   PetscFunctionReturn(PETSC_SUCCESS);
1492: }

1494: static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A, PetscInt n, const IS islocal[], const IS isglobal[], PetscBool colflg, ISLocalToGlobalMapping *ltog)
1495: {
1496:   PetscBool flg;
1497:   PetscInt  i, j, m, mi, *ix;

1499:   PetscFunctionBegin;
1500:   *ltog = NULL;
1501:   for (i = 0, m = 0, flg = PETSC_FALSE; i < n; i++) {
1502:     if (islocal[i]) {
1503:       PetscCall(ISGetLocalSize(islocal[i], &mi));
1504:       flg = PETSC_TRUE; /* We found a non-trivial entry */
1505:     } else {
1506:       PetscCall(ISGetLocalSize(isglobal[i], &mi));
1507:     }
1508:     m += mi;
1509:   }
1510:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

1512:   PetscCall(PetscMalloc1(m, &ix));
1513:   for (i = 0, m = 0; i < n; i++) {
1514:     ISLocalToGlobalMapping smap = NULL;
1515:     Mat                    sub  = NULL;
1516:     PetscSF                sf;
1517:     PetscLayout            map;
1518:     const PetscInt        *ix2;

1520:     if (!colflg) {
1521:       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
1522:     } else {
1523:       PetscCall(MatNestFindNonzeroSubMatCol(A, i, &sub));
1524:     }
1525:     if (sub) {
1526:       if (!colflg) PetscCall(MatGetLocalToGlobalMapping(sub, &smap, NULL));
1527:       else PetscCall(MatGetLocalToGlobalMapping(sub, NULL, &smap));
1528:     }
1529:     /*
1530:        Now we need to extract the monolithic global indices that correspond to the given split global indices.
1531:        In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1532:     */
1533:     PetscCall(ISGetIndices(isglobal[i], &ix2));
1534:     if (islocal[i]) {
1535:       PetscInt *ilocal, *iremote;
1536:       PetscInt  mil, nleaves;

1538:       PetscCall(ISGetLocalSize(islocal[i], &mi));
1539:       PetscCheck(smap, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing local to global map");
1540:       for (j = 0; j < mi; j++) ix[m + j] = j;
1541:       PetscCall(ISLocalToGlobalMappingApply(smap, mi, ix + m, ix + m));

1543:       /* PetscSFSetGraphLayout does not like negative indices */
1544:       PetscCall(PetscMalloc2(mi, &ilocal, mi, &iremote));
1545:       for (j = 0, nleaves = 0; j < mi; j++) {
1546:         if (ix[m + j] < 0) continue;
1547:         ilocal[nleaves]  = j;
1548:         iremote[nleaves] = ix[m + j];
1549:         nleaves++;
1550:       }
1551:       PetscCall(ISGetLocalSize(isglobal[i], &mil));
1552:       PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
1553:       PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)A), &map));
1554:       PetscCall(PetscLayoutSetLocalSize(map, mil));
1555:       PetscCall(PetscLayoutSetUp(map));
1556:       PetscCall(PetscSFSetGraphLayout(sf, map, nleaves, ilocal, PETSC_USE_POINTER, iremote));
1557:       PetscCall(PetscLayoutDestroy(&map));
1558:       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, ix2, ix + m, MPI_REPLACE));
1559:       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, ix2, ix + m, MPI_REPLACE));
1560:       PetscCall(PetscSFDestroy(&sf));
1561:       PetscCall(PetscFree2(ilocal, iremote));
1562:     } else {
1563:       PetscCall(ISGetLocalSize(isglobal[i], &mi));
1564:       for (j = 0; j < mi; j++) ix[m + j] = ix2[i];
1565:     }
1566:     PetscCall(ISRestoreIndices(isglobal[i], &ix2));
1567:     m += mi;
1568:   }
1569:   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A), 1, m, ix, PETSC_OWN_POINTER, ltog));
1570:   PetscFunctionReturn(PETSC_SUCCESS);
1571: }

1573: /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1574: /*
1575:   nprocessors = NP
1576:   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1577:        proc 0: => (g_0,h_0,)
1578:        proc 1: => (g_1,h_1,)
1579:        ...
1580:        proc nprocs-1: => (g_NP-1,h_NP-1,)

1582:             proc 0:                      proc 1:                    proc nprocs-1:
1583:     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))

1585:             proc 0:
1586:     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1587:             proc 1:
1588:     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)

1590:             proc NP-1:
1591:     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1592: */
1593: static PetscErrorCode MatSetUp_NestIS_Private(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[])
1594: {
1595:   Mat_Nest *vs = (Mat_Nest *)A->data;
1596:   PetscInt  i, j, offset, n, nsum, bs;
1597:   Mat       sub = NULL;

1599:   PetscFunctionBegin;
1600:   PetscCall(PetscMalloc1(nr, &vs->isglobal.row));
1601:   PetscCall(PetscMalloc1(nc, &vs->isglobal.col));
1602:   if (is_row) { /* valid IS is passed in */
1603:     /* refs on is[] are incremented */
1604:     for (i = 0; i < vs->nr; i++) {
1605:       PetscCall(PetscObjectReference((PetscObject)is_row[i]));
1606:       vs->isglobal.row[i] = is_row[i];
1607:     }
1608:   } else { /* Create the ISs by inspecting sizes of a submatrix in each row */
1609:     nsum = 0;
1610:     for (i = 0; i < vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */
1611:       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
1612:       PetscCheck(sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "No nonzero submatrix in row %" PetscInt_FMT, i);
1613:       PetscCall(MatGetLocalSize(sub, &n, NULL));
1614:       PetscCheck(n >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Sizes have not yet been set for submatrix");
1615:       nsum += n;
1616:     }
1617:     PetscCallMPI(MPI_Scan(&nsum, &offset, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
1618:     offset -= nsum;
1619:     for (i = 0; i < vs->nr; i++) {
1620:       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
1621:       PetscCall(MatGetLocalSize(sub, &n, NULL));
1622:       PetscCall(MatGetBlockSizes(sub, &bs, NULL));
1623:       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub), n, offset, 1, &vs->isglobal.row[i]));
1624:       PetscCall(ISSetBlockSize(vs->isglobal.row[i], bs));
1625:       offset += n;
1626:     }
1627:   }

1629:   if (is_col) { /* valid IS is passed in */
1630:     /* refs on is[] are incremented */
1631:     for (j = 0; j < vs->nc; j++) {
1632:       PetscCall(PetscObjectReference((PetscObject)is_col[j]));
1633:       vs->isglobal.col[j] = is_col[j];
1634:     }
1635:   } else { /* Create the ISs by inspecting sizes of a submatrix in each column */
1636:     offset = A->cmap->rstart;
1637:     nsum   = 0;
1638:     for (j = 0; j < vs->nc; j++) {
1639:       PetscCall(MatNestFindNonzeroSubMatCol(A, j, &sub));
1640:       PetscCheck(sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "No nonzero submatrix in column %" PetscInt_FMT, i);
1641:       PetscCall(MatGetLocalSize(sub, NULL, &n));
1642:       PetscCheck(n >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Sizes have not yet been set for submatrix");
1643:       nsum += n;
1644:     }
1645:     PetscCallMPI(MPI_Scan(&nsum, &offset, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
1646:     offset -= nsum;
1647:     for (j = 0; j < vs->nc; j++) {
1648:       PetscCall(MatNestFindNonzeroSubMatCol(A, j, &sub));
1649:       PetscCall(MatGetLocalSize(sub, NULL, &n));
1650:       PetscCall(MatGetBlockSizes(sub, NULL, &bs));
1651:       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub), n, offset, 1, &vs->isglobal.col[j]));
1652:       PetscCall(ISSetBlockSize(vs->isglobal.col[j], bs));
1653:       offset += n;
1654:     }
1655:   }

1657:   /* Set up the local ISs */
1658:   PetscCall(PetscMalloc1(vs->nr, &vs->islocal.row));
1659:   PetscCall(PetscMalloc1(vs->nc, &vs->islocal.col));
1660:   for (i = 0, offset = 0; i < vs->nr; i++) {
1661:     IS                     isloc;
1662:     ISLocalToGlobalMapping rmap = NULL;
1663:     PetscInt               nlocal, bs;
1664:     PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
1665:     if (sub) PetscCall(MatGetLocalToGlobalMapping(sub, &rmap, NULL));
1666:     if (rmap) {
1667:       PetscCall(MatGetBlockSizes(sub, &bs, NULL));
1668:       PetscCall(ISLocalToGlobalMappingGetSize(rmap, &nlocal));
1669:       PetscCall(ISCreateStride(PETSC_COMM_SELF, nlocal, offset, 1, &isloc));
1670:       PetscCall(ISSetBlockSize(isloc, bs));
1671:     } else {
1672:       nlocal = 0;
1673:       isloc  = NULL;
1674:     }
1675:     vs->islocal.row[i] = isloc;
1676:     offset += nlocal;
1677:   }
1678:   for (i = 0, offset = 0; i < vs->nc; i++) {
1679:     IS                     isloc;
1680:     ISLocalToGlobalMapping cmap = NULL;
1681:     PetscInt               nlocal, bs;
1682:     PetscCall(MatNestFindNonzeroSubMatCol(A, i, &sub));
1683:     if (sub) PetscCall(MatGetLocalToGlobalMapping(sub, NULL, &cmap));
1684:     if (cmap) {
1685:       PetscCall(MatGetBlockSizes(sub, NULL, &bs));
1686:       PetscCall(ISLocalToGlobalMappingGetSize(cmap, &nlocal));
1687:       PetscCall(ISCreateStride(PETSC_COMM_SELF, nlocal, offset, 1, &isloc));
1688:       PetscCall(ISSetBlockSize(isloc, bs));
1689:     } else {
1690:       nlocal = 0;
1691:       isloc  = NULL;
1692:     }
1693:     vs->islocal.col[i] = isloc;
1694:     offset += nlocal;
1695:   }

1697:   /* Set up the aggregate ISLocalToGlobalMapping */
1698:   {
1699:     ISLocalToGlobalMapping rmap, cmap;
1700:     PetscCall(MatNestCreateAggregateL2G_Private(A, vs->nr, vs->islocal.row, vs->isglobal.row, PETSC_FALSE, &rmap));
1701:     PetscCall(MatNestCreateAggregateL2G_Private(A, vs->nc, vs->islocal.col, vs->isglobal.col, PETSC_TRUE, &cmap));
1702:     if (rmap && cmap) PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
1703:     PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
1704:     PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
1705:   }

1707:   if (PetscDefined(USE_DEBUG)) {
1708:     for (i = 0; i < vs->nr; i++) {
1709:       for (j = 0; j < vs->nc; j++) {
1710:         PetscInt m, n, M, N, mi, ni, Mi, Ni;
1711:         Mat      B = vs->m[i][j];
1712:         if (!B) continue;
1713:         PetscCall(MatGetSize(B, &M, &N));
1714:         PetscCall(MatGetLocalSize(B, &m, &n));
1715:         PetscCall(ISGetSize(vs->isglobal.row[i], &Mi));
1716:         PetscCall(ISGetSize(vs->isglobal.col[j], &Ni));
1717:         PetscCall(ISGetLocalSize(vs->isglobal.row[i], &mi));
1718:         PetscCall(ISGetLocalSize(vs->isglobal.col[j], &ni));
1719:         PetscCheck(M == Mi && N == Ni, PetscObjectComm((PetscObject)sub), PETSC_ERR_ARG_INCOMP, "Global sizes (%" PetscInt_FMT ",%" PetscInt_FMT ") of nested submatrix (%" PetscInt_FMT ",%" PetscInt_FMT ") do not agree with space defined by index sets (%" PetscInt_FMT ",%" PetscInt_FMT ")", M, N, i, j, Mi, Ni);
1720:         PetscCheck(m == mi && n == ni, PetscObjectComm((PetscObject)sub), PETSC_ERR_ARG_INCOMP, "Local sizes (%" PetscInt_FMT ",%" PetscInt_FMT ") of nested submatrix (%" PetscInt_FMT ",%" PetscInt_FMT ") do not agree with space defined by index sets (%" PetscInt_FMT ",%" PetscInt_FMT ")", m, n, i, j, mi, ni);
1721:       }
1722:     }
1723:   }

1725:   /* Set A->assembled if all non-null blocks are currently assembled */
1726:   for (i = 0; i < vs->nr; i++) {
1727:     for (j = 0; j < vs->nc; j++) {
1728:       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(PETSC_SUCCESS);
1729:     }
1730:   }
1731:   A->assembled = PETSC_TRUE;
1732:   PetscFunctionReturn(PETSC_SUCCESS);
1733: }

1735: /*@C
1736:   MatCreateNest - Creates a new `MATNEST` matrix containing several nested submatrices, each stored separately

1738:   Collective

1740:   Input Parameters:
1741: + comm   - Communicator for the new `MATNEST`
1742: . nr     - number of nested row blocks
1743: . is_row - index sets for each nested row block, or `NULL` to make contiguous
1744: . nc     - number of nested column blocks
1745: . is_col - index sets for each nested column block, or `NULL` to make contiguous
1746: - a      - array of $nr \times nc$ submatrices, empty submatrices can be passed using `NULL`

1748:   Output Parameter:
1749: . B - new matrix

1751:   Level: advanced

1753:   Note:
1754:   In both C and Fortran, `a` must be a one-dimensional array representing a two-dimensional row-major order array holding references to the matrices.
1755:   For instance, to represent the matrix
1756:   $\begin{bmatrix} A_{11} & A_{12} \\ A_{21} & A_{22}\end{bmatrix}$
1757:   one should use `Mat a[4]={A11,A12,A21,A22}`.

1759:   Fortran Note:
1760:   Pass `PETSC_NULL_MAT` in the corresponding entry of `a` for an empty block

1762: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreate()`, `VecCreateNest()`, `DMCreateMatrix()`, `MatNestSetSubMat()`,
1763:           `MatNestGetSubMat()`, `MatNestGetLocalISs()`, `MatNestGetSize()`,
1764:           `MatNestGetISs()`, `MatNestSetSubMats()`, `MatNestGetSubMats()`
1765: @*/
1766: PetscErrorCode MatCreateNest(MPI_Comm comm, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[], Mat *B) PeNSS
1767: {
1768:   PetscFunctionBegin;
1769:   PetscCall(MatCreate(comm, B));
1770:   PetscCall(MatSetType(*B, MATNEST));
1771:   (*B)->preallocated = PETSC_TRUE;
1772:   PetscCall(MatNestSetSubMats(*B, nr, is_row, nc, is_col, a));
1773:   PetscFunctionReturn(PETSC_SUCCESS);
1774: }

1776: static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1777: {
1778:   Mat_Nest     *nest = (Mat_Nest *)A->data;
1779:   Mat          *trans;
1780:   PetscScalar **avv;
1781:   PetscScalar  *vv;
1782:   PetscInt    **aii, **ajj;
1783:   PetscInt     *ii, *jj, *ci;
1784:   PetscInt      nr, nc, nnz, i, j;
1785:   PetscBool     done;

1787:   PetscFunctionBegin;
1788:   PetscCall(MatGetSize(A, &nr, &nc));
1789:   if (reuse == MAT_REUSE_MATRIX) {
1790:     PetscInt rnr;

1792:     PetscCall(MatGetRowIJ(*newmat, 0, PETSC_FALSE, PETSC_FALSE, &rnr, (const PetscInt **)&ii, (const PetscInt **)&jj, &done));
1793:     PetscCheck(done, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatGetRowIJ");
1794:     PetscCheck(rnr == nr, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Cannot reuse matrix, wrong number of rows");
1795:     PetscCall(MatSeqAIJGetArray(*newmat, &vv));
1796:   }
1797:   /* extract CSR for nested SeqAIJ matrices */
1798:   nnz = 0;
1799:   PetscCall(PetscCalloc4(nest->nr * nest->nc, &aii, nest->nr * nest->nc, &ajj, nest->nr * nest->nc, &avv, nest->nr * nest->nc, &trans));
1800:   for (i = 0; i < nest->nr; ++i) {
1801:     for (j = 0; j < nest->nc; ++j) {
1802:       Mat B = nest->m[i][j];
1803:       if (B) {
1804:         PetscScalar *naa;
1805:         PetscInt    *nii, *njj, nnr;
1806:         PetscBool    istrans;

1808:         PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &istrans));
1809:         if (istrans) {
1810:           Mat Bt;

1812:           PetscCall(MatTransposeGetMat(B, &Bt));
1813:           PetscCall(MatTranspose(Bt, MAT_INITIAL_MATRIX, &trans[i * nest->nc + j]));
1814:           B = trans[i * nest->nc + j];
1815:         } else {
1816:           PetscCall(PetscObjectTypeCompare((PetscObject)B, MATHERMITIANTRANSPOSEVIRTUAL, &istrans));
1817:           if (istrans) {
1818:             Mat Bt;

1820:             PetscCall(MatHermitianTransposeGetMat(B, &Bt));
1821:             PetscCall(MatHermitianTranspose(Bt, MAT_INITIAL_MATRIX, &trans[i * nest->nc + j]));
1822:             B = trans[i * nest->nc + j];
1823:           }
1824:         }
1825:         PetscCall(MatGetRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nnr, (const PetscInt **)&nii, (const PetscInt **)&njj, &done));
1826:         PetscCheck(done, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "MatGetRowIJ");
1827:         PetscCall(MatSeqAIJGetArray(B, &naa));
1828:         nnz += nii[nnr];

1830:         aii[i * nest->nc + j] = nii;
1831:         ajj[i * nest->nc + j] = njj;
1832:         avv[i * nest->nc + j] = naa;
1833:       }
1834:     }
1835:   }
1836:   if (reuse != MAT_REUSE_MATRIX) {
1837:     PetscCall(PetscMalloc1(nr + 1, &ii));
1838:     PetscCall(PetscMalloc1(nnz, &jj));
1839:     PetscCall(PetscMalloc1(nnz, &vv));
1840:   } else {
1841:     PetscCheck(nnz == ii[nr], PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Cannot reuse matrix, wrong number of nonzeros");
1842:   }

1844:   /* new row pointer */
1845:   PetscCall(PetscArrayzero(ii, nr + 1));
1846:   for (i = 0; i < nest->nr; ++i) {
1847:     PetscInt ncr, rst;

1849:     PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &rst, NULL));
1850:     PetscCall(ISGetLocalSize(nest->isglobal.row[i], &ncr));
1851:     for (j = 0; j < nest->nc; ++j) {
1852:       if (aii[i * nest->nc + j]) {
1853:         PetscInt *nii = aii[i * nest->nc + j];
1854:         PetscInt  ir;

1856:         for (ir = rst; ir < ncr + rst; ++ir) {
1857:           ii[ir + 1] += nii[1] - nii[0];
1858:           nii++;
1859:         }
1860:       }
1861:     }
1862:   }
1863:   for (i = 0; i < nr; i++) ii[i + 1] += ii[i];

1865:   /* construct CSR for the new matrix */
1866:   PetscCall(PetscCalloc1(nr, &ci));
1867:   for (i = 0; i < nest->nr; ++i) {
1868:     PetscInt ncr, rst;

1870:     PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &rst, NULL));
1871:     PetscCall(ISGetLocalSize(nest->isglobal.row[i], &ncr));
1872:     for (j = 0; j < nest->nc; ++j) {
1873:       if (aii[i * nest->nc + j]) {
1874:         PetscScalar *nvv = avv[i * nest->nc + j], vscale = 1.0, vshift = 0.0;
1875:         PetscInt    *nii = aii[i * nest->nc + j];
1876:         PetscInt    *njj = ajj[i * nest->nc + j];
1877:         PetscInt     ir, cst;

1879:         if (trans[i * nest->nc + j]) {
1880:           vscale = ((Mat_Shell *)nest->m[i][j]->data)->vscale;
1881:           vshift = ((Mat_Shell *)nest->m[i][j]->data)->vshift;
1882:         }
1883:         PetscCall(ISStrideGetInfo(nest->isglobal.col[j], &cst, NULL));
1884:         for (ir = rst; ir < ncr + rst; ++ir) {
1885:           PetscInt ij, rsize = nii[1] - nii[0], ist = ii[ir] + ci[ir];

1887:           for (ij = 0; ij < rsize; ij++) {
1888:             jj[ist + ij] = *njj + cst;
1889:             vv[ist + ij] = vscale * *nvv;
1890:             if (PetscUnlikely(vshift != 0.0 && *njj == ir - rst)) vv[ist + ij] += vshift;
1891:             njj++;
1892:             nvv++;
1893:           }
1894:           ci[ir] += rsize;
1895:           nii++;
1896:         }
1897:       }
1898:     }
1899:   }
1900:   PetscCall(PetscFree(ci));

1902:   /* restore info */
1903:   for (i = 0; i < nest->nr; ++i) {
1904:     for (j = 0; j < nest->nc; ++j) {
1905:       Mat B = nest->m[i][j];
1906:       if (B) {
1907:         PetscInt nnr = 0, k = i * nest->nc + j;

1909:         B = (trans[k] ? trans[k] : B);
1910:         PetscCall(MatRestoreRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nnr, (const PetscInt **)&aii[k], (const PetscInt **)&ajj[k], &done));
1911:         PetscCheck(done, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "MatRestoreRowIJ");
1912:         PetscCall(MatSeqAIJRestoreArray(B, &avv[k]));
1913:         PetscCall(MatDestroy(&trans[k]));
1914:       }
1915:     }
1916:   }
1917:   PetscCall(PetscFree4(aii, ajj, avv, trans));

1919:   /* finalize newmat */
1920:   if (reuse == MAT_INITIAL_MATRIX) {
1921:     PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A), nr, nc, ii, jj, vv, newmat));
1922:   } else if (reuse == MAT_INPLACE_MATRIX) {
1923:     Mat B;

1925:     PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A), nr, nc, ii, jj, vv, &B));
1926:     PetscCall(MatHeaderReplace(A, &B));
1927:   }
1928:   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
1929:   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
1930:   {
1931:     Mat_SeqAIJ *a = (Mat_SeqAIJ *)((*newmat)->data);
1932:     a->free_a     = PETSC_TRUE;
1933:     a->free_ij    = PETSC_TRUE;
1934:   }
1935:   PetscFunctionReturn(PETSC_SUCCESS);
1936: }

1938: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat Y, PetscScalar a, Mat X)
1939: {
1940:   Mat_Nest *nest = (Mat_Nest *)X->data;
1941:   PetscInt  i, j, k, rstart;
1942:   PetscBool flg;

1944:   PetscFunctionBegin;
1945:   /* Fill by row */
1946:   for (j = 0; j < nest->nc; ++j) {
1947:     /* Using global column indices and ISAllGather() is not scalable. */
1948:     IS              bNis;
1949:     PetscInt        bN;
1950:     const PetscInt *bNindices;
1951:     PetscCall(ISAllGather(nest->isglobal.col[j], &bNis));
1952:     PetscCall(ISGetSize(bNis, &bN));
1953:     PetscCall(ISGetIndices(bNis, &bNindices));
1954:     for (i = 0; i < nest->nr; ++i) {
1955:       Mat             B = nest->m[i][j], D = NULL;
1956:       PetscInt        bm, br;
1957:       const PetscInt *bmindices;
1958:       if (!B) continue;
1959:       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
1960:       if (flg) {
1961:         PetscTryMethod(B, "MatTransposeGetMat_C", (Mat, Mat *), (B, &D));
1962:         PetscTryMethod(B, "MatHermitianTransposeGetMat_C", (Mat, Mat *), (B, &D));
1963:         PetscCall(MatConvert(B, ((PetscObject)D)->type_name, MAT_INITIAL_MATRIX, &D));
1964:         B = D;
1965:       }
1966:       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
1967:       if (flg) {
1968:         if (D) PetscCall(MatConvert(D, MATBAIJ, MAT_INPLACE_MATRIX, &D));
1969:         else PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &D));
1970:         B = D;
1971:       }
1972:       PetscCall(ISGetLocalSize(nest->isglobal.row[i], &bm));
1973:       PetscCall(ISGetIndices(nest->isglobal.row[i], &bmindices));
1974:       PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
1975:       for (br = 0; br < bm; ++br) {
1976:         PetscInt           row = bmindices[br], brncols, *cols;
1977:         const PetscInt    *brcols;
1978:         const PetscScalar *brcoldata;
1979:         PetscScalar       *vals = NULL;
1980:         PetscCall(MatGetRow(B, br + rstart, &brncols, &brcols, &brcoldata));
1981:         PetscCall(PetscMalloc1(brncols, &cols));
1982:         for (k = 0; k < brncols; k++) cols[k] = bNindices[brcols[k]];
1983:         /*
1984:           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1985:           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1986:          */
1987:         if (a != 1.0) {
1988:           PetscCall(PetscMalloc1(brncols, &vals));
1989:           for (k = 0; k < brncols; k++) vals[k] = a * brcoldata[k];
1990:           PetscCall(MatSetValues(Y, 1, &row, brncols, cols, vals, ADD_VALUES));
1991:           PetscCall(PetscFree(vals));
1992:         } else {
1993:           PetscCall(MatSetValues(Y, 1, &row, brncols, cols, brcoldata, ADD_VALUES));
1994:         }
1995:         PetscCall(MatRestoreRow(B, br + rstart, &brncols, &brcols, &brcoldata));
1996:         PetscCall(PetscFree(cols));
1997:       }
1998:       if (D) PetscCall(MatDestroy(&D));
1999:       PetscCall(ISRestoreIndices(nest->isglobal.row[i], &bmindices));
2000:     }
2001:     PetscCall(ISRestoreIndices(bNis, &bNindices));
2002:     PetscCall(ISDestroy(&bNis));
2003:   }
2004:   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
2005:   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
2006:   PetscFunctionReturn(PETSC_SUCCESS);
2007: }

2009: static PetscErrorCode MatConvert_Nest_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
2010: {
2011:   Mat_Nest   *nest = (Mat_Nest *)A->data;
2012:   PetscInt    m, n, M, N, i, j, k, *dnnz, *onnz = NULL, rstart, cstart, cend;
2013:   PetscMPIInt size;
2014:   Mat         C;

2016:   PetscFunctionBegin;
2017:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
2018:   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
2019:     PetscInt  nf;
2020:     PetscBool fast;

2022:     PetscCall(PetscStrcmp(newtype, MATAIJ, &fast));
2023:     if (!fast) PetscCall(PetscStrcmp(newtype, MATSEQAIJ, &fast));
2024:     for (i = 0; i < nest->nr && fast; ++i) {
2025:       for (j = 0; j < nest->nc && fast; ++j) {
2026:         Mat B = nest->m[i][j];
2027:         if (B) {
2028:           PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &fast));
2029:           if (!fast) {
2030:             PetscBool istrans;

2032:             PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &istrans));
2033:             if (istrans) {
2034:               Mat Bt;

2036:               PetscCall(MatTransposeGetMat(B, &Bt));
2037:               PetscCall(PetscObjectTypeCompare((PetscObject)Bt, MATSEQAIJ, &fast));
2038:             } else {
2039:               PetscCall(PetscObjectTypeCompare((PetscObject)B, MATHERMITIANTRANSPOSEVIRTUAL, &istrans));
2040:               if (istrans) {
2041:                 Mat Bt;

2043:                 PetscCall(MatHermitianTransposeGetMat(B, &Bt));
2044:                 PetscCall(PetscObjectTypeCompare((PetscObject)Bt, MATSEQAIJ, &fast));
2045:               }
2046:             }
2047:             if (fast) fast = (PetscBool)(!((Mat_Shell *)B->data)->zrows && !((Mat_Shell *)B->data)->zcols && !((Mat_Shell *)B->data)->axpy && !((Mat_Shell *)B->data)->left && !((Mat_Shell *)B->data)->right && !((Mat_Shell *)B->data)->dshift);
2048:           }
2049:         }
2050:       }
2051:     }
2052:     for (i = 0, nf = 0; i < nest->nr && fast; ++i) {
2053:       PetscCall(PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i], ISSTRIDE, &fast));
2054:       if (fast) {
2055:         PetscInt f, s;

2057:         PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &f, &s));
2058:         if (f != nf || s != 1) {
2059:           fast = PETSC_FALSE;
2060:         } else {
2061:           PetscCall(ISGetSize(nest->isglobal.row[i], &f));
2062:           nf += f;
2063:         }
2064:       }
2065:     }
2066:     for (i = 0, nf = 0; i < nest->nc && fast; ++i) {
2067:       PetscCall(PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i], ISSTRIDE, &fast));
2068:       if (fast) {
2069:         PetscInt f, s;

2071:         PetscCall(ISStrideGetInfo(nest->isglobal.col[i], &f, &s));
2072:         if (f != nf || s != 1) {
2073:           fast = PETSC_FALSE;
2074:         } else {
2075:           PetscCall(ISGetSize(nest->isglobal.col[i], &f));
2076:           nf += f;
2077:         }
2078:       }
2079:     }
2080:     if (fast) {
2081:       PetscCall(MatConvert_Nest_SeqAIJ_fast(A, newtype, reuse, newmat));
2082:       PetscFunctionReturn(PETSC_SUCCESS);
2083:     }
2084:   }
2085:   PetscCall(MatGetSize(A, &M, &N));
2086:   PetscCall(MatGetLocalSize(A, &m, &n));
2087:   PetscCall(MatGetOwnershipRangeColumn(A, &cstart, &cend));
2088:   if (reuse == MAT_REUSE_MATRIX) C = *newmat;
2089:   else {
2090:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
2091:     PetscCall(MatSetType(C, newtype));
2092:     PetscCall(MatSetSizes(C, m, n, M, N));
2093:   }
2094:   PetscCall(PetscMalloc1(2 * m, &dnnz));
2095:   if (m) {
2096:     onnz = dnnz + m;
2097:     for (k = 0; k < m; k++) {
2098:       dnnz[k] = 0;
2099:       onnz[k] = 0;
2100:     }
2101:   }
2102:   for (j = 0; j < nest->nc; ++j) {
2103:     IS              bNis;
2104:     PetscInt        bN;
2105:     const PetscInt *bNindices;
2106:     PetscBool       flg;
2107:     /* Using global column indices and ISAllGather() is not scalable. */
2108:     PetscCall(ISAllGather(nest->isglobal.col[j], &bNis));
2109:     PetscCall(ISGetSize(bNis, &bN));
2110:     PetscCall(ISGetIndices(bNis, &bNindices));
2111:     for (i = 0; i < nest->nr; ++i) {
2112:       PetscSF         bmsf;
2113:       PetscSFNode    *iremote;
2114:       Mat             B = nest->m[i][j], D = NULL;
2115:       PetscInt        bm, *sub_dnnz, *sub_onnz, br;
2116:       const PetscInt *bmindices;
2117:       if (!B) continue;
2118:       PetscCall(ISGetLocalSize(nest->isglobal.row[i], &bm));
2119:       PetscCall(ISGetIndices(nest->isglobal.row[i], &bmindices));
2120:       PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf));
2121:       PetscCall(PetscMalloc1(bm, &iremote));
2122:       PetscCall(PetscMalloc1(bm, &sub_dnnz));
2123:       PetscCall(PetscMalloc1(bm, &sub_onnz));
2124:       for (k = 0; k < bm; ++k) {
2125:         sub_dnnz[k] = 0;
2126:         sub_onnz[k] = 0;
2127:       }
2128:       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
2129:       if (flg) {
2130:         PetscTryMethod(B, "MatTransposeGetMat_C", (Mat, Mat *), (B, &D));
2131:         PetscTryMethod(B, "MatHermitianTransposeGetMat_C", (Mat, Mat *), (B, &D));
2132:         PetscCall(MatConvert(B, ((PetscObject)D)->type_name, MAT_INITIAL_MATRIX, &D));
2133:         B = D;
2134:       }
2135:       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
2136:       if (flg) {
2137:         if (D) PetscCall(MatConvert(D, MATBAIJ, MAT_INPLACE_MATRIX, &D));
2138:         else PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &D));
2139:         B = D;
2140:       }
2141:       /*
2142:        Locate the owners for all of the locally-owned global row indices for this row block.
2143:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
2144:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
2145:        */
2146:       PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
2147:       for (br = 0; br < bm; ++br) {
2148:         PetscInt        row = bmindices[br], brncols, col;
2149:         const PetscInt *brcols;
2150:         PetscInt        rowrel   = 0; /* row's relative index on its owner rank */
2151:         PetscMPIInt     rowowner = 0;
2152:         PetscCall(PetscLayoutFindOwnerIndex(A->rmap, row, &rowowner, &rowrel));
2153:         /* how many roots  */
2154:         iremote[br].rank  = rowowner;
2155:         iremote[br].index = rowrel; /* edge from bmdnnz to dnnz */
2156:         /* get nonzero pattern */
2157:         PetscCall(MatGetRow(B, br + rstart, &brncols, &brcols, NULL));
2158:         for (k = 0; k < brncols; k++) {
2159:           col = bNindices[brcols[k]];
2160:           if (col >= A->cmap->range[rowowner] && col < A->cmap->range[rowowner + 1]) {
2161:             sub_dnnz[br]++;
2162:           } else {
2163:             sub_onnz[br]++;
2164:           }
2165:         }
2166:         PetscCall(MatRestoreRow(B, br + rstart, &brncols, &brcols, NULL));
2167:       }
2168:       if (D) PetscCall(MatDestroy(&D));
2169:       PetscCall(ISRestoreIndices(nest->isglobal.row[i], &bmindices));
2170:       /* bsf will have to take care of disposing of bedges. */
2171:       PetscCall(PetscSFSetGraph(bmsf, m, bm, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
2172:       PetscCall(PetscSFReduceBegin(bmsf, MPIU_INT, sub_dnnz, dnnz, MPI_SUM));
2173:       PetscCall(PetscSFReduceEnd(bmsf, MPIU_INT, sub_dnnz, dnnz, MPI_SUM));
2174:       PetscCall(PetscSFReduceBegin(bmsf, MPIU_INT, sub_onnz, onnz, MPI_SUM));
2175:       PetscCall(PetscSFReduceEnd(bmsf, MPIU_INT, sub_onnz, onnz, MPI_SUM));
2176:       PetscCall(PetscFree(sub_dnnz));
2177:       PetscCall(PetscFree(sub_onnz));
2178:       PetscCall(PetscSFDestroy(&bmsf));
2179:     }
2180:     PetscCall(ISRestoreIndices(bNis, &bNindices));
2181:     PetscCall(ISDestroy(&bNis));
2182:   }
2183:   /* Resize preallocation if overestimated */
2184:   for (i = 0; i < m; i++) {
2185:     dnnz[i] = PetscMin(dnnz[i], A->cmap->n);
2186:     onnz[i] = PetscMin(onnz[i], A->cmap->N - A->cmap->n);
2187:   }
2188:   PetscCall(MatSeqAIJSetPreallocation(C, 0, dnnz));
2189:   PetscCall(MatMPIAIJSetPreallocation(C, 0, dnnz, 0, onnz));
2190:   PetscCall(PetscFree(dnnz));
2191:   PetscCall(MatAXPY_Dense_Nest(C, 1.0, A));
2192:   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &C));
2193:   else *newmat = C;
2194:   PetscFunctionReturn(PETSC_SUCCESS);
2195: }

2197: static PetscErrorCode MatConvert_Nest_Dense(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
2198: {
2199:   Mat      B;
2200:   PetscInt m, n, M, N;

2202:   PetscFunctionBegin;
2203:   PetscCall(MatGetSize(A, &M, &N));
2204:   PetscCall(MatGetLocalSize(A, &m, &n));
2205:   if (reuse == MAT_REUSE_MATRIX) {
2206:     B = *newmat;
2207:     PetscCall(MatZeroEntries(B));
2208:   } else {
2209:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), m, PETSC_DECIDE, M, N, NULL, &B));
2210:   }
2211:   PetscCall(MatAXPY_Dense_Nest(B, 1.0, A));
2212:   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
2213:   else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
2214:   PetscFunctionReturn(PETSC_SUCCESS);
2215: }

2217: static PetscErrorCode MatHasOperation_Nest(Mat mat, MatOperation op, PetscBool *has)
2218: {
2219:   Mat_Nest    *bA = (Mat_Nest *)mat->data;
2220:   MatOperation opAdd;
2221:   PetscInt     i, j, nr = bA->nr, nc = bA->nc;
2222:   PetscBool    flg;

2224:   PetscFunctionBegin;
2225:   *has = PETSC_FALSE;
2226:   if (op == MATOP_MULT || op == MATOP_MULT_ADD || op == MATOP_MULT_TRANSPOSE || op == MATOP_MULT_TRANSPOSE_ADD) {
2227:     opAdd = (op == MATOP_MULT || op == MATOP_MULT_ADD ? MATOP_MULT_ADD : MATOP_MULT_TRANSPOSE_ADD);
2228:     for (j = 0; j < nc; j++) {
2229:       for (i = 0; i < nr; i++) {
2230:         if (!bA->m[i][j]) continue;
2231:         PetscCall(MatHasOperation(bA->m[i][j], opAdd, &flg));
2232:         if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
2233:       }
2234:     }
2235:   }
2236:   if (((void **)mat->ops)[op]) *has = PETSC_TRUE;
2237:   PetscFunctionReturn(PETSC_SUCCESS);
2238: }

2240: /*MC
2241:   MATNEST -  "nest" - Matrix type consisting of nested submatrices, each stored separately.

2243:   Level: intermediate

2245:   Notes:
2246:   This matrix type permits scalable use of `PCFIELDSPLIT` and avoids the large memory costs of extracting submatrices.
2247:   It allows the use of symmetric and block formats for parts of multi-physics simulations.
2248:   It is usually used with `DMCOMPOSITE` and `DMCreateMatrix()`

2250:   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
2251:   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
2252:   than the nest matrix.

2254: .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreate()`, `MatType`, `MatCreateNest()`, `MatNestSetSubMat()`, `MatNestGetSubMat()`,
2255:           `VecCreateNest()`, `DMCreateMatrix()`, `DMCOMPOSITE`, `MatNestSetVecType()`, `MatNestGetLocalISs()`,
2256:           `MatNestGetISs()`, `MatNestSetSubMats()`, `MatNestGetSubMats()`
2257: M*/
2258: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
2259: {
2260:   Mat_Nest *s;

2262:   PetscFunctionBegin;
2263:   PetscCall(PetscNew(&s));
2264:   A->data = (void *)s;

2266:   s->nr            = -1;
2267:   s->nc            = -1;
2268:   s->m             = NULL;
2269:   s->splitassembly = PETSC_FALSE;

2271:   PetscCall(PetscMemzero(A->ops, sizeof(*A->ops)));

2273:   A->ops->mult                      = MatMult_Nest;
2274:   A->ops->multadd                   = MatMultAdd_Nest;
2275:   A->ops->multtranspose             = MatMultTranspose_Nest;
2276:   A->ops->multtransposeadd          = MatMultTransposeAdd_Nest;
2277:   A->ops->transpose                 = MatTranspose_Nest;
2278:   A->ops->multhermitiantranspose    = MatMultHermitianTranspose_Nest;
2279:   A->ops->multhermitiantransposeadd = MatMultHermitianTransposeAdd_Nest;
2280:   A->ops->assemblybegin             = MatAssemblyBegin_Nest;
2281:   A->ops->assemblyend               = MatAssemblyEnd_Nest;
2282:   A->ops->zeroentries               = MatZeroEntries_Nest;
2283:   A->ops->copy                      = MatCopy_Nest;
2284:   A->ops->axpy                      = MatAXPY_Nest;
2285:   A->ops->duplicate                 = MatDuplicate_Nest;
2286:   A->ops->createsubmatrix           = MatCreateSubMatrix_Nest;
2287:   A->ops->destroy                   = MatDestroy_Nest;
2288:   A->ops->view                      = MatView_Nest;
2289:   A->ops->getvecs                   = NULL; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
2290:   A->ops->getlocalsubmatrix         = MatGetLocalSubMatrix_Nest;
2291:   A->ops->restorelocalsubmatrix     = MatRestoreLocalSubMatrix_Nest;
2292:   A->ops->getdiagonal               = MatGetDiagonal_Nest;
2293:   A->ops->diagonalscale             = MatDiagonalScale_Nest;
2294:   A->ops->scale                     = MatScale_Nest;
2295:   A->ops->shift                     = MatShift_Nest;
2296:   A->ops->diagonalset               = MatDiagonalSet_Nest;
2297:   A->ops->setrandom                 = MatSetRandom_Nest;
2298:   A->ops->hasoperation              = MatHasOperation_Nest;

2300:   A->spptr     = NULL;
2301:   A->assembled = PETSC_FALSE;

2303:   /* expose Nest api's */
2304:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMat_C", MatNestGetSubMat_Nest));
2305:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMat_C", MatNestSetSubMat_Nest));
2306:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMats_C", MatNestGetSubMats_Nest));
2307:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSize_C", MatNestGetSize_Nest));
2308:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetISs_C", MatNestGetISs_Nest));
2309:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetLocalISs_C", MatNestGetLocalISs_Nest));
2310:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetVecType_C", MatNestSetVecType_Nest));
2311:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMats_C", MatNestSetSubMats_Nest));
2312:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpiaij_C", MatConvert_Nest_AIJ));
2313:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqaij_C", MatConvert_Nest_AIJ));
2314:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_aij_C", MatConvert_Nest_AIJ));
2315:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_is_C", MatConvert_Nest_IS));
2316:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpidense_C", MatConvert_Nest_Dense));
2317:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqdense_C", MatConvert_Nest_Dense));
2318:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_seqdense_C", MatProductSetFromOptions_Nest_Dense));
2319:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_mpidense_C", MatProductSetFromOptions_Nest_Dense));

2321:   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATNEST));
2322:   PetscFunctionReturn(PETSC_SUCCESS);
2323: }