Actual source code: schurm.c

  1: #include <../src/ksp/ksp/utils/schurm/schurm.h>

  3: const char *const MatSchurComplementAinvTypes[] = {"DIAG", "LUMP", "BLOCKDIAG", "FULL", "MatSchurComplementAinvType", "MAT_SCHUR_COMPLEMENT_AINV_", NULL};

  5: PetscErrorCode MatCreateVecs_SchurComplement(Mat N, Vec *right, Vec *left)
  6: {
  7:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

  9:   PetscFunctionBegin;
 10:   if (Na->D) {
 11:     PetscCall(MatCreateVecs(Na->D, right, left));
 12:     PetscFunctionReturn(PETSC_SUCCESS);
 13:   }
 14:   if (right) PetscCall(MatCreateVecs(Na->B, right, NULL));
 15:   if (left) PetscCall(MatCreateVecs(Na->C, NULL, left));
 16:   PetscFunctionReturn(PETSC_SUCCESS);
 17: }

 19: PetscErrorCode MatView_SchurComplement(Mat N, PetscViewer viewer)
 20: {
 21:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

 23:   PetscFunctionBegin;
 24:   PetscCall(PetscViewerASCIIPrintf(viewer, "Schur complement A11 - A10 inv(A00) A01\n"));
 25:   if (Na->D) {
 26:     PetscCall(PetscViewerASCIIPrintf(viewer, "A11\n"));
 27:     PetscCall(PetscViewerASCIIPushTab(viewer));
 28:     PetscCall(MatView(Na->D, viewer));
 29:     PetscCall(PetscViewerASCIIPopTab(viewer));
 30:   } else {
 31:     PetscCall(PetscViewerASCIIPrintf(viewer, "A11 = 0\n"));
 32:   }
 33:   PetscCall(PetscViewerASCIIPrintf(viewer, "A10\n"));
 34:   PetscCall(PetscViewerASCIIPushTab(viewer));
 35:   PetscCall(MatView(Na->C, viewer));
 36:   PetscCall(PetscViewerASCIIPopTab(viewer));
 37:   PetscCall(PetscViewerASCIIPrintf(viewer, "KSP solver for A00 block viewable with the additional option -%sksp_view\n", ((PetscObject)Na->ksp)->prefix ? ((PetscObject)Na->ksp)->prefix : NULL));
 38:   PetscCall(PetscViewerASCIIPrintf(viewer, "A01\n"));
 39:   PetscCall(PetscViewerASCIIPushTab(viewer));
 40:   PetscCall(MatView(Na->B, viewer));
 41:   PetscCall(PetscViewerASCIIPopTab(viewer));
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: /*
 46:            A11^T - A01^T ksptrans(A00,Ap00) A10^T
 47: */
 48: PetscErrorCode MatMultTranspose_SchurComplement(Mat N, Vec x, Vec y)
 49: {
 50:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

 52:   PetscFunctionBegin;
 53:   if (!Na->work1) PetscCall(MatCreateVecs(Na->A, &Na->work1, NULL));
 54:   if (!Na->work2) PetscCall(MatCreateVecs(Na->A, &Na->work2, NULL));
 55:   PetscCall(MatMultTranspose(Na->C, x, Na->work1));
 56:   PetscCall(KSPSolveTranspose(Na->ksp, Na->work1, Na->work2));
 57:   PetscCall(MatMultTranspose(Na->B, Na->work2, y));
 58:   PetscCall(VecScale(y, -1.0));
 59:   if (Na->D) PetscCall(MatMultTransposeAdd(Na->D, x, y, y));
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: /*
 64:            A11 - A10 ksp(A00,Ap00) A01
 65: */
 66: PetscErrorCode MatMult_SchurComplement(Mat N, Vec x, Vec y)
 67: {
 68:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

 70:   PetscFunctionBegin;
 71:   if (!Na->work1) PetscCall(MatCreateVecs(Na->A, &Na->work1, NULL));
 72:   if (!Na->work2) PetscCall(MatCreateVecs(Na->A, &Na->work2, NULL));
 73:   PetscCall(MatMult(Na->B, x, Na->work1));
 74:   PetscCall(KSPSolve(Na->ksp, Na->work1, Na->work2));
 75:   PetscCall(MatMult(Na->C, Na->work2, y));
 76:   PetscCall(VecScale(y, -1.0));
 77:   if (Na->D) PetscCall(MatMultAdd(Na->D, x, y, y));
 78:   PetscFunctionReturn(PETSC_SUCCESS);
 79: }

 81: /*
 82:            A11 - A10 ksp(A00,Ap00) A01
 83: */
 84: PetscErrorCode MatMultAdd_SchurComplement(Mat N, Vec x, Vec y, Vec z)
 85: {
 86:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

 88:   PetscFunctionBegin;
 89:   if (!Na->work1) PetscCall(MatCreateVecs(Na->A, &Na->work1, NULL));
 90:   if (!Na->work2) PetscCall(MatCreateVecs(Na->A, &Na->work2, NULL));
 91:   PetscCall(MatMult(Na->B, x, Na->work1));
 92:   PetscCall(KSPSolve(Na->ksp, Na->work1, Na->work2));
 93:   if (y == z) {
 94:     PetscCall(VecScale(Na->work2, -1.0));
 95:     PetscCall(MatMultAdd(Na->C, Na->work2, z, z));
 96:   } else {
 97:     PetscCall(MatMult(Na->C, Na->work2, z));
 98:     PetscCall(VecAYPX(z, -1.0, y));
 99:   }
100:   if (Na->D) PetscCall(MatMultAdd(Na->D, x, z, z));
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N, PetscOptionItems *PetscOptionsObject)
105: {
106:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

108:   PetscFunctionBegin;
109:   PetscOptionsHeadBegin(PetscOptionsObject, "MatSchurComplementOptions");
110:   Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
111:   PetscCall(PetscOptionsEnum("-mat_schur_complement_ainv_type", "Type of approximation for DIAGFORM(A00) used when assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01", "MatSchurComplementSetAinvType", MatSchurComplementAinvTypes, (PetscEnum)Na->ainvtype,
112:                              (PetscEnum *)&Na->ainvtype, NULL));
113:   PetscOptionsHeadEnd();
114:   PetscCall(KSPSetFromOptions(Na->ksp));
115:   PetscFunctionReturn(PETSC_SUCCESS);
116: }

118: PetscErrorCode MatDestroy_SchurComplement(Mat N)
119: {
120:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

122:   PetscFunctionBegin;
123:   PetscCall(MatDestroy(&Na->A));
124:   PetscCall(MatDestroy(&Na->Ap));
125:   PetscCall(MatDestroy(&Na->B));
126:   PetscCall(MatDestroy(&Na->C));
127:   PetscCall(MatDestroy(&Na->D));
128:   PetscCall(VecDestroy(&Na->work1));
129:   PetscCall(VecDestroy(&Na->work2));
130:   PetscCall(KSPDestroy(&Na->ksp));
131:   PetscCall(PetscFree(N->data));
132:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", NULL));
133:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", NULL));
134:   PetscFunctionReturn(PETSC_SUCCESS);
135: }

137: /*@
138:   MatCreateSchurComplement - Creates a new `Mat` that behaves like the Schur complement of a matrix

140:   Collective

142:   Input Parameters:
143: + A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
144: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A00^{-1}
145: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
146: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
147: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

149:   Output Parameter:
150: . S - the matrix that behaves as the Schur complement S = A11 - A10 ksp(A00,Ap00) A01

152:   Level: intermediate

154:   Notes:
155:   The Schur complement is NOT explicitly formed! Rather, this function returns a virtual Schur complement
156:   that can compute the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
157:   for Schur complement S and a `KSP` solver to approximate the action of A^{-1}.

159:   All four matrices must have the same MPI communicator.

161:   `A00` and  `A11` must be square matrices.

163:   `MatGetSchurComplement()` takes as arguments the index sets for the submatrices and returns both the virtual Schur complement (what this returns) plus
164:   a sparse approximation to the Schur complement (useful for building a preconditioner for the Schur complement) which can be obtained from this
165:   matrix with `MatSchurComplementGetPmat()`

167:   Developer Notes:
168:   The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
169:   remove redundancy and be clearer and simpler.

171: .seealso: [](ch_ksp), `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatGetSchurComplement()`,
172:           `MatSchurComplementGetPmat()`, `MatSchurComplementSetSubMatrices()`
173: @*/
174: PetscErrorCode MatCreateSchurComplement(Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11, Mat *S)
175: {
176:   PetscFunctionBegin;
177:   PetscCall(KSPInitializePackage());
178:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A00), S));
179:   PetscCall(MatSetType(*S, MATSCHURCOMPLEMENT));
180:   PetscCall(MatSchurComplementSetSubMatrices(*S, A00, Ap00, A01, A10, A11));
181:   PetscFunctionReturn(PETSC_SUCCESS);
182: }

184: /*@
185:   MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement

187:   Collective

189:   Input Parameters:
190: + S    - matrix obtained with `MatSetType`(S,`MATSCHURCOMPLEMENT`)
191: . A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
192: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A00^{-1}
193: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
194: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
195: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

197:   Level: intermediate

199:   Notes:
200:   The Schur complement is NOT explicitly formed! Rather, this
201:   object performs the matrix-vector product of the Schur complement by using formula S = A11 - A10 ksp(A00,Ap00) A01

203:   All four matrices must have the same MPI communicator.

205:   `A00` and `A11` must be square matrices.

207:   This is to be used in the context of code such as
208: .vb
209:      MatSetType(S,MATSCHURCOMPLEMENT);
210:      MatSchurComplementSetSubMatrices(S,...);
211: .ve
212:   while `MatSchurComplementUpdateSubMatrices()` should only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`

214: .seealso: [](ch_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`
215: @*/
216: PetscErrorCode MatSchurComplementSetSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
217: {
218:   Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
219:   PetscBool            isschur;

221:   PetscFunctionBegin;
222:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
223:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
224:   PetscCheck(!S->assembled, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Use MatSchurComplementUpdateSubMatrices() for already used matrix");
229:   PetscCheckSameComm(A00, 2, Ap00, 3);
230:   PetscCheckSameComm(A00, 2, A01, 4);
231:   PetscCheckSameComm(A00, 2, A10, 5);
232:   PetscCheck(A00->rmap->n == A00->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A00 %" PetscInt_FMT " do not equal local columns %" PetscInt_FMT, A00->rmap->n, A00->cmap->n);
233:   PetscCheck(A00->rmap->n == Ap00->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A00 %" PetscInt_FMT " do not equal local rows of Ap00 %" PetscInt_FMT, A00->rmap->n, Ap00->rmap->n);
234:   PetscCheck(Ap00->rmap->n == Ap00->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of Ap00 %" PetscInt_FMT " do not equal local columns %" PetscInt_FMT, Ap00->rmap->n, Ap00->cmap->n);
235:   PetscCheck(A00->cmap->n == A01->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local columns of A00 %" PetscInt_FMT " do not equal local rows of A01 %" PetscInt_FMT, A00->cmap->n, A01->rmap->n);
236:   PetscCheck(A10->cmap->n == A00->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local columns of A10 %" PetscInt_FMT " do not equal local rows of A00 %" PetscInt_FMT, A10->cmap->n, A00->rmap->n);
237:   if (A11) {
239:     PetscCheckSameComm(A00, 2, A11, 6);
240:     PetscCheck(A10->rmap->n == A11->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A10 %" PetscInt_FMT " do not equal local rows A11 %" PetscInt_FMT, A10->rmap->n, A11->rmap->n);
241:   }

243:   PetscCall(MatSetSizes(S, A10->rmap->n, A01->cmap->n, A10->rmap->N, A01->cmap->N));
244:   PetscCall(PetscObjectReference((PetscObject)A00));
245:   PetscCall(PetscObjectReference((PetscObject)Ap00));
246:   PetscCall(PetscObjectReference((PetscObject)A01));
247:   PetscCall(PetscObjectReference((PetscObject)A10));
248:   Na->A  = A00;
249:   Na->Ap = Ap00;
250:   Na->B  = A01;
251:   Na->C  = A10;
252:   Na->D  = A11;
253:   if (A11) PetscCall(PetscObjectReference((PetscObject)A11));
254:   PetscCall(MatSetUp(S));
255:   PetscCall(KSPSetOperators(Na->ksp, A00, Ap00));
256:   S->assembled = PETSC_TRUE;
257:   PetscFunctionReturn(PETSC_SUCCESS);
258: }

260: /*@
261:   MatSchurComplementGetKSP - Gets the `KSP` object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01

263:   Not Collective

265:   Input Parameter:
266: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

268:   Output Parameter:
269: . ksp - the linear solver object

271:   Options Database Key:
272: . -fieldsplit_<splitname_0>_XXX - sets KSP and PC options for the 0-split solver inside the Schur complement used in `PCFIELDSPLIT`; default <splitname_0> is 0.

274:   Level: intermediate

276: .seealso: [](ch_ksp), `Mat`, `MatSchurComplementSetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`
277: @*/
278: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
279: {
280:   Mat_SchurComplement *Na;
281:   PetscBool            isschur;

283:   PetscFunctionBegin;
285:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
286:   PetscCheck(isschur, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)S)->type_name);
287:   PetscAssertPointer(ksp, 2);
288:   Na   = (Mat_SchurComplement *)S->data;
289:   *ksp = Na->ksp;
290:   PetscFunctionReturn(PETSC_SUCCESS);
291: }

293: /*@
294:   MatSchurComplementSetKSP - Sets the `KSP` object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01

296:   Not Collective

298:   Input Parameters:
299: + S   - matrix created with `MatCreateSchurComplement()`
300: - ksp - the linear solver object

302:   Level: developer

304:   Developer Notes:
305:   This is used in `PCFIELDSPLIT` to reuse the 0-split `KSP` to implement ksp(A00,Ap00) in S.
306:   The `KSP` operators are overwritten with A00 and Ap00 currently set in S.

308: .seealso: [](ch_ksp), `Mat`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MATSCHURCOMPLEMENT`
309: @*/
310: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
311: {
312:   Mat_SchurComplement *Na;
313:   PetscBool            isschur;

315:   PetscFunctionBegin;
317:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
318:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
320:   Na = (Mat_SchurComplement *)S->data;
321:   PetscCall(PetscObjectReference((PetscObject)ksp));
322:   PetscCall(KSPDestroy(&Na->ksp));
323:   Na->ksp = ksp;
324:   PetscCall(KSPSetOperators(Na->ksp, Na->A, Na->Ap));
325:   PetscFunctionReturn(PETSC_SUCCESS);
326: }

328: /*@
329:   MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices

331:   Collective

333:   Input Parameters:
334: + S    - matrix obtained with `MatCreateSchurComplement()` (or `MatSchurSetSubMatrices()`) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
335: . A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
336: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A00^{-1}
337: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
338: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
339: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

341:   Level: intermediate

343:   Notes:
344:   All four matrices must have the same MPI communicator

346:   `A00` and  `A11` must be square matrices

348:   All of the matrices provided must have the same sizes as was used with `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`
349:   though they need not be the same matrices.

351:   This can only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`, it cannot be called immediately after `MatSetType`(S,`MATSCHURCOMPLEMENT`);

353:   Developer Notes:
354:   This code is almost identical to `MatSchurComplementSetSubMatrices()`. The API should be refactored.

356: .seealso: [](ch_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`
357: @*/
358: PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
359: {
360:   Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
361:   PetscBool            isschur;

363:   PetscFunctionBegin;
365:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
366:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
367:   PetscCheck(S->assembled, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Use MatSchurComplementSetSubMatrices() for a new matrix");
372:   PetscCheckSameComm(A00, 2, Ap00, 3);
373:   PetscCheckSameComm(A00, 2, A01, 4);
374:   PetscCheckSameComm(A00, 2, A10, 5);
375:   PetscCheck(A00->rmap->n == A00->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A00 %" PetscInt_FMT " do not equal local columns %" PetscInt_FMT, A00->rmap->n, A00->cmap->n);
376:   PetscCheck(A00->rmap->n == Ap00->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A00 %" PetscInt_FMT " do not equal local rows of Ap00 %" PetscInt_FMT, A00->rmap->n, Ap00->rmap->n);
377:   PetscCheck(Ap00->rmap->n == Ap00->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of Ap00 %" PetscInt_FMT " do not equal local columns %" PetscInt_FMT, Ap00->rmap->n, Ap00->cmap->n);
378:   PetscCheck(A00->cmap->n == A01->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local columns of A00 %" PetscInt_FMT " do not equal local rows of A01 %" PetscInt_FMT, A00->cmap->n, A01->rmap->n);
379:   PetscCheck(A10->cmap->n == A00->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local columns of A10 %" PetscInt_FMT " do not equal local rows of A00 %" PetscInt_FMT, A10->cmap->n, A00->rmap->n);
380:   if (A11) {
382:     PetscCheckSameComm(A00, 2, A11, 6);
383:     PetscCheck(A10->rmap->n == A11->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local rows of A10 %" PetscInt_FMT " do not equal local rows A11 %" PetscInt_FMT, A10->rmap->n, A11->rmap->n);
384:   }

386:   PetscCall(PetscObjectReference((PetscObject)A00));
387:   PetscCall(PetscObjectReference((PetscObject)Ap00));
388:   PetscCall(PetscObjectReference((PetscObject)A01));
389:   PetscCall(PetscObjectReference((PetscObject)A10));
390:   if (A11) PetscCall(PetscObjectReference((PetscObject)A11));

392:   PetscCall(MatDestroy(&Na->A));
393:   PetscCall(MatDestroy(&Na->Ap));
394:   PetscCall(MatDestroy(&Na->B));
395:   PetscCall(MatDestroy(&Na->C));
396:   PetscCall(MatDestroy(&Na->D));

398:   Na->A  = A00;
399:   Na->Ap = Ap00;
400:   Na->B  = A01;
401:   Na->C  = A10;
402:   Na->D  = A11;

404:   PetscCall(KSPSetOperators(Na->ksp, A00, Ap00));
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }

408: /*@C
409:   MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement

411:   Collective

413:   Input Parameter:
414: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

416:   Output Parameters:
417: + A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
418: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
419: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
420: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
421: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

423:   Level: intermediate

425:   Note:
426:   `A11` is optional, and thus can be `NULL`.

428:   The reference counts of the submatrices are not increased before they are returned and the matrices should not be modified or destroyed.

430: .seealso: [](ch_ksp), `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatSchurComplementUpdateSubMatrices()`
431: @*/
432: PetscErrorCode MatSchurComplementGetSubMatrices(Mat S, Mat *A00, Mat *Ap00, Mat *A01, Mat *A10, Mat *A11)
433: {
434:   Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
435:   PetscBool            flg;

437:   PetscFunctionBegin;
439:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &flg));
440:   PetscCheck(flg, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)S)->type_name);
441:   if (A00) *A00 = Na->A;
442:   if (Ap00) *Ap00 = Na->Ap;
443:   if (A01) *A01 = Na->B;
444:   if (A10) *A10 = Na->C;
445:   if (A11) *A11 = Na->D;
446:   PetscFunctionReturn(PETSC_SUCCESS);
447: }

449: #include <petsc/private/kspimpl.h>

451: /*@
452:   MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly

454:   Collective

456:   Input Parameter:
457: . A - the matrix obtained with `MatCreateSchurComplement()`

459:   Output Parameter:
460: . S - the Schur complement matrix

462:   Notes:
463:   This can be expensive, so it is mainly for testing

465:   Use `MatSchurComplementGetPmat()` to get a sparse approximation for the Schur complement suitable for use in building a preconditioner

467:   Level: advanced

469: .seealso: [](ch_ksp), `MatCreateSchurComplement()`, `MatSchurComplementUpdate()`, `MatSchurComplementGetPmat()`
470: @*/
471: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat A, Mat *S)
472: {
473:   Mat       B, C, D, E = NULL, Bd, AinvBd;
474:   KSP       ksp;
475:   PetscInt  n, N, m, M;
476:   PetscBool flg = PETSC_FALSE, set, symm;

478:   PetscFunctionBegin;
479:   PetscCall(MatSchurComplementGetSubMatrices(A, NULL, NULL, &B, &C, &D));
480:   PetscCall(MatSchurComplementGetKSP(A, &ksp));
481:   PetscCall(KSPSetUp(ksp));
482:   PetscCall(MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd));
483:   PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd));
484:   PetscCall(KSPMatSolve(ksp, Bd, AinvBd));
485:   PetscCall(MatDestroy(&Bd));
486:   PetscCall(MatFilter(AinvBd, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
487:   if (D) {
488:     PetscCall(MatGetLocalSize(D, &m, &n));
489:     PetscCall(MatGetSize(D, &M, &N));
490:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), m, n, M, N, NULL, S));
491:   }
492:   PetscCall(MatMatMult(C, AinvBd, D ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, S));
493:   PetscCall(MatDestroy(&AinvBd));
494:   if (D) {
495:     PetscCall(PetscObjectTypeCompareAny((PetscObject)D, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
496:     if (flg) {
497:       PetscCall(MatIsSymmetricKnown(A, &set, &symm));
498:       if (!set || !symm) PetscCall(MatConvert(D, MATBAIJ, MAT_INITIAL_MATRIX, &E)); /* convert the (1,1) block to nonsymmetric storage for MatAXPY() */
499:     }
500:     PetscCall(MatAXPY(*S, -1.0, E ? E : D, DIFFERENT_NONZERO_PATTERN)); /* calls Mat[Get|Restore]RowUpperTriangular(), so only the upper triangular part is valid with symmetric storage */
501:   }
502:   PetscCall(MatConvert(*S, !E && flg ? MATSBAIJ : MATAIJ, MAT_INPLACE_MATRIX, S)); /* if A is symmetric and the (1,1) block is a MatSBAIJ, return S as a MatSBAIJ */
503:   PetscCall(MatScale(*S, -1.0));
504:   PetscCall(MatDestroy(&E));
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }

508: /* Developer Notes:
509:     This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
510: PetscErrorCode MatGetSchurComplement_Basic(Mat mat, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
511: {
512:   Mat      A = NULL, Ap = NULL, B = NULL, C = NULL, D = NULL;
513:   MatReuse reuse;

515:   PetscFunctionBegin;
525:   if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);

529:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

531:   reuse = MAT_INITIAL_MATRIX;
532:   if (mreuse == MAT_REUSE_MATRIX) {
533:     PetscCall(MatSchurComplementGetSubMatrices(*S, &A, &Ap, &B, &C, &D));
534:     PetscCheck(A && Ap && B && C, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Attempting to reuse matrix but Schur complement matrices unset");
535:     PetscCheck(A == Ap, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Preconditioning matrix does not match operator");
536:     PetscCall(MatDestroy(&Ap)); /* get rid of extra reference */
537:     reuse = MAT_REUSE_MATRIX;
538:   }
539:   PetscCall(MatCreateSubMatrix(mat, isrow0, iscol0, reuse, &A));
540:   PetscCall(MatCreateSubMatrix(mat, isrow0, iscol1, reuse, &B));
541:   PetscCall(MatCreateSubMatrix(mat, isrow1, iscol0, reuse, &C));
542:   PetscCall(MatCreateSubMatrix(mat, isrow1, iscol1, reuse, &D));
543:   switch (mreuse) {
544:   case MAT_INITIAL_MATRIX:
545:     PetscCall(MatCreateSchurComplement(A, A, B, C, D, S));
546:     break;
547:   case MAT_REUSE_MATRIX:
548:     PetscCall(MatSchurComplementUpdateSubMatrices(*S, A, A, B, C, D));
549:     break;
550:   default:
551:     PetscCheck(mreuse == MAT_IGNORE_MATRIX, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Unrecognized value of mreuse %d", (int)mreuse);
552:   }
553:   if (preuse != MAT_IGNORE_MATRIX) PetscCall(MatCreateSchurComplementPmat(A, B, C, D, ainvtype, preuse, Sp));
554:   PetscCall(MatDestroy(&A));
555:   PetscCall(MatDestroy(&B));
556:   PetscCall(MatDestroy(&C));
557:   PetscCall(MatDestroy(&D));
558:   PetscFunctionReturn(PETSC_SUCCESS);
559: }

561: /*@
562:   MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.

564:   Collective

566:   Input Parameters:
567: + A        - matrix in which the complement is to be taken
568: . isrow0   - rows to eliminate
569: . iscol0   - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
570: . isrow1   - rows in which the Schur complement is formed
571: . iscol1   - columns in which the Schur complement is formed
572: . mreuse   - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in S
573: . ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
574:                        `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
575: - preuse   - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in Sp

577:   Output Parameters:
578: + S  - exact Schur complement, often of type `MATSCHURCOMPLEMENT` which is difficult to use for preconditioning
579: - Sp - approximate Schur complement from which a preconditioner can be built A11 - A10 inv(DIAGFORM(A00)) A01

581:   Level: advanced

583:   Notes:
584:   Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
585:   application-specific information.

587:   Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
588:   and column index sets.  In that case, the user should call `PetscObjectComposeFunction()` on the *S matrix and pass mreuse of `MAT_REUSE_MATRIX` to set
589:   "MatGetSchurComplement_C" to their function.  If their function needs to fall back to the default implementation, it
590:   should call `MatGetSchurComplement_Basic()`.

592:   `MatCreateSchurComplement()` takes as arguments the four submatrices and returns the virtual Schur complement (what this function returns in S).

594:   `MatSchurComplementGetPmat()` takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).

596:   In other words calling `MatCreateSchurComplement()` followed by `MatSchurComplementGetPmat()` produces the same output as this function but with slightly different
597:   inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.

599:   Developer Notes:
600:   The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
601:   remove redundancy and be clearer and simpler.

603: .seealso: [](ch_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatCreateSchurComplement()`, `MatSchurComplementAinvType`
604: @*/
605: PetscErrorCode MatGetSchurComplement(Mat A, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
606: {
607:   PetscErrorCode (*f)(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatReuse, Mat *) = NULL;

609:   PetscFunctionBegin;
621:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
622:   if (mreuse == MAT_REUSE_MATRIX) { /* This is the only situation, in which we can demand that the user pass a non-NULL pointer to non-garbage in S. */
623:     PetscCall(PetscObjectQueryFunction((PetscObject)*S, "MatGetSchurComplement_C", &f));
624:   }
625:   if (f) PetscCall((*f)(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, preuse, Sp));
626:   else PetscCall(MatGetSchurComplement_Basic(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, ainvtype, preuse, Sp));
627:   PetscFunctionReturn(PETSC_SUCCESS);
628: }

630: /*@
631:   MatSchurComplementSetAinvType - set the type of approximation used for the inverse of the (0,0) block used in forming Sp in `MatSchurComplementGetPmat()`

633:   Not Collective

635:   Input Parameters:
636: + S        - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
637: - ainvtype - type of approximation to be used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
638:                       `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`

640:   Options Database Key:
641: . -mat_schur_complement_ainv_type diag | lump | blockdiag | full - set schur complement type

643:   Level: advanced

645: .seealso: [](ch_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementGetAinvType()`
646: @*/
647: PetscErrorCode MatSchurComplementSetAinvType(Mat S, MatSchurComplementAinvType ainvtype)
648: {
649:   PetscBool            isschur;
650:   Mat_SchurComplement *schur;

652:   PetscFunctionBegin;
654:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
655:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
657:   schur = (Mat_SchurComplement *)S->data;
658:   PetscCheck(ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_FULL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown MatSchurComplementAinvType: %d", (int)ainvtype);
659:   schur->ainvtype = ainvtype;
660:   PetscFunctionReturn(PETSC_SUCCESS);
661: }

663: /*@
664:   MatSchurComplementGetAinvType - get the type of approximation for the inverse of the (0,0) block used in forming Sp in `MatSchurComplementGetPmat()`

666:   Not Collective

668:   Input Parameter:
669: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

671:   Output Parameter:
672: . ainvtype - type of approximation used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
673:                       `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`

675:   Level: advanced

677: .seealso: [](ch_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementSetAinvType()`
678: @*/
679: PetscErrorCode MatSchurComplementGetAinvType(Mat S, MatSchurComplementAinvType *ainvtype)
680: {
681:   PetscBool            isschur;
682:   Mat_SchurComplement *schur;

684:   PetscFunctionBegin;
686:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
687:   PetscCheck(isschur, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)S)->type_name);
688:   schur = (Mat_SchurComplement *)S->data;
689:   if (ainvtype) *ainvtype = schur->ainvtype;
690:   PetscFunctionReturn(PETSC_SUCCESS);
691: }

693: /*@
694:   MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by explicitly assembling the sparse matrix
695:   Sp = A11 - A10 inv(DIAGFORM(A00)) A01

697:   Collective

699:   Input Parameters:
700: + A00      - the upper-left part of the original matrix A = [A00 A01; A10 A11]
701: . A01      - (optional) the upper-right part of the original matrix A = [A00 A01; A10 A11]
702: . A10      - (optional) the lower-left part of the original matrix A = [A00 A01; A10 A11]
703: . A11      - (optional) the lower-right part of the original matrix A = [A00 A01; A10 A11]
704: . ainvtype - type of approximation for DIAGFORM(A00) used when forming Sp = A11 - A10 inv(DIAGFORM(A00)) A01. See MatSchurComplementAinvType.
705: - preuse   - `MAT_INITIAL_MATRIX` for a new Sp, or `MAT_REUSE_MATRIX` to reuse an existing Sp, or `MAT_IGNORE_MATRIX` to put nothing in Sp

707:   Output Parameter:
708: . Sp - approximate Schur complement suitable for preconditioning the true Schur complement S = A11 - A10 inv(A00) A01

710:   Level: advanced

712: .seealso: [](ch_ksp), `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementAinvType`
713: @*/
714: PetscErrorCode MatCreateSchurComplementPmat(Mat A00, Mat A01, Mat A10, Mat A11, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
715: {
716:   PetscInt N00;

718:   PetscFunctionBegin;
719:   /* Use an appropriate approximate inverse of A00 to form A11 - A10 inv(DIAGFORM(A00)) A01; a NULL A01, A10 or A11 indicates a zero matrix. */
720:   /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
722:   if (preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);

724:   /* A zero size A00 or empty A01 or A10 imply S = A11. */
725:   PetscCall(MatGetSize(A00, &N00, NULL));
726:   if (!A01 || !A10 || !N00) {
727:     if (preuse == MAT_INITIAL_MATRIX) {
728:       PetscCall(MatDuplicate(A11, MAT_COPY_VALUES, Sp));
729:     } else { /* MAT_REUSE_MATRIX */
730:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
731:       PetscCall(MatCopy(A11, *Sp, DIFFERENT_NONZERO_PATTERN));
732:     }
733:   } else {
734:     Mat AdB;
735:     Vec diag;

737:     if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
738:       PetscCall(MatDuplicate(A01, MAT_COPY_VALUES, &AdB));
739:       PetscCall(MatCreateVecs(A00, &diag, NULL));
740:       if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
741:         PetscCall(MatGetRowSum(A00, diag));
742:       } else {
743:         PetscCall(MatGetDiagonal(A00, diag));
744:       }
745:       PetscCall(VecReciprocal(diag));
746:       PetscCall(MatDiagonalScale(AdB, diag, NULL));
747:       PetscCall(VecDestroy(&diag));
748:     } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
749:       Mat      A00_inv;
750:       MatType  type;
751:       MPI_Comm comm;

753:       PetscCall(PetscObjectGetComm((PetscObject)A00, &comm));
754:       PetscCall(MatGetType(A00, &type));
755:       PetscCall(MatCreate(comm, &A00_inv));
756:       PetscCall(MatSetType(A00_inv, type));
757:       PetscCall(MatInvertBlockDiagonalMat(A00, A00_inv));
758:       PetscCall(MatMatMult(A00_inv, A01, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &AdB));
759:       PetscCall(MatDestroy(&A00_inv));
760:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown MatSchurComplementAinvType: %d", ainvtype);
761:     /* Cannot really reuse Sp in MatMatMult() because of MatAYPX() -->
762:          MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult()  */
763:     if (preuse == MAT_REUSE_MATRIX) PetscCall(MatDestroy(Sp));
764:     PetscCall(MatMatMult(A10, AdB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, Sp));
765:     if (!A11) {
766:       PetscCall(MatScale(*Sp, -1.0));
767:     } else {
768:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
769:       PetscCall(MatAYPX(*Sp, -1, A11, DIFFERENT_NONZERO_PATTERN));
770:     }
771:     PetscCall(MatDestroy(&AdB));
772:   }
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }

776: static PetscErrorCode MatSchurComplementGetPmat_Basic(Mat S, MatReuse preuse, Mat *Sp)
777: {
778:   Mat                  A, B, C, D;
779:   Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;

781:   PetscFunctionBegin;
782:   if (preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);
783:   PetscCall(MatSchurComplementGetSubMatrices(S, &A, NULL, &B, &C, &D));
784:   PetscCheck(A, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Schur complement component matrices unset");
785:   if (schur->ainvtype != MAT_SCHUR_COMPLEMENT_AINV_FULL) PetscCall(MatCreateSchurComplementPmat(A, B, C, D, schur->ainvtype, preuse, Sp));
786:   else {
787:     if (preuse == MAT_REUSE_MATRIX) PetscCall(MatDestroy(Sp));
788:     PetscCall(MatSchurComplementComputeExplicitOperator(S, Sp));
789:   }
790:   PetscFunctionReturn(PETSC_SUCCESS);
791: }

793: /*@
794:   MatSchurComplementGetPmat - Obtain a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01

796:   Collective

798:   Input Parameters:
799: + S      - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of A11 - A10 ksp(A00,Ap00) A01
800: - preuse - `MAT_INITIAL_MATRIX` for a new Sp, or `MAT_REUSE_MATRIX` to reuse an existing Sp, or `MAT_IGNORE_MATRIX` to put nothing in Sp

802:   Output Parameter:
803: . Sp - approximate Schur complement suitable for preconditioning the exact Schur complement S = A11 - A10 inv(A00) A01

805:   Level: advanced

807:   Notes:
808:   The approximation of Sp depends on the the argument passed to to `MatSchurComplementSetAinvType()`
809:   `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
810:   -mat_schur_complement_ainv_type <diag,lump,blockdiag,full>

812:   Sometimes users would like to provide problem-specific data in the Schur complement, usually only
813:   for special row and column index sets.  In that case, the user should call `PetscObjectComposeFunction()` to set
814:   "MatSchurComplementGetPmat_C" to their function.  If their function needs to fall back to the default implementation,
815:   it should call `MatSchurComplementGetPmat_Basic()`.

817:   Developer Notes:
818:   The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
819:   remove redundancy and be clearer and simpler.

821:   This routine should be called `MatSchurComplementCreatePmat()`

823: .seealso: [](ch_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetAinvType()`
824: @*/
825: PetscErrorCode MatSchurComplementGetPmat(Mat S, MatReuse preuse, Mat *Sp)
826: {
827:   PetscErrorCode (*f)(Mat, MatReuse, Mat *);

829:   PetscFunctionBegin;
833:   if (preuse != MAT_IGNORE_MATRIX) {
834:     PetscAssertPointer(Sp, 3);
835:     if (preuse == MAT_INITIAL_MATRIX) *Sp = NULL;
837:   }
838:   PetscCheck(!S->factortype, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

840:   PetscCall(PetscObjectQueryFunction((PetscObject)S, "MatSchurComplementGetPmat_C", &f));
841:   if (f) PetscCall((*f)(S, preuse, Sp));
842:   else PetscCall(MatSchurComplementGetPmat_Basic(S, preuse, Sp));
843:   PetscFunctionReturn(PETSC_SUCCESS);
844: }

846: static PetscErrorCode MatProductNumeric_SchurComplement_Dense(Mat C)
847: {
848:   Mat_Product         *product = C->product;
849:   Mat_SchurComplement *Na      = (Mat_SchurComplement *)product->A->data;
850:   Mat                  work1, work2;
851:   PetscScalar         *v;
852:   PetscInt             lda;

854:   PetscFunctionBegin;
855:   PetscCall(MatMatMult(Na->B, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1));
856:   PetscCall(MatDuplicate(work1, MAT_DO_NOT_COPY_VALUES, &work2));
857:   PetscCall(KSPMatSolve(Na->ksp, work1, work2));
858:   PetscCall(MatDestroy(&work1));
859:   PetscCall(MatDenseGetArrayWrite(C, &v));
860:   PetscCall(MatDenseGetLDA(C, &lda));
861:   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)C), C->rmap->n, C->cmap->n, C->rmap->N, C->cmap->N, v, &work1));
862:   PetscCall(MatDenseSetLDA(work1, lda));
863:   PetscCall(MatMatMult(Na->C, work2, MAT_REUSE_MATRIX, PETSC_DEFAULT, &work1));
864:   PetscCall(MatDenseRestoreArrayWrite(C, &v));
865:   PetscCall(MatDestroy(&work2));
866:   PetscCall(MatDestroy(&work1));
867:   if (Na->D) {
868:     PetscCall(MatMatMult(Na->D, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1));
869:     PetscCall(MatAYPX(C, -1.0, work1, SAME_NONZERO_PATTERN));
870:     PetscCall(MatDestroy(&work1));
871:   } else PetscCall(MatScale(C, -1.0));
872:   PetscFunctionReturn(PETSC_SUCCESS);
873: }

875: static PetscErrorCode MatProductSymbolic_SchurComplement_Dense(Mat C)
876: {
877:   Mat_Product *product = C->product;
878:   Mat          A = product->A, B = product->B;
879:   PetscInt     m = A->rmap->n, n = B->cmap->n, M = A->rmap->N, N = B->cmap->N;
880:   PetscBool    flg;

882:   PetscFunctionBegin;
883:   PetscCall(MatSetSizes(C, m, n, M, N));
884:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)C, &flg, MATSEQDENSE, MATMPIDENSE, ""));
885:   if (!flg) {
886:     PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
887:     C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
888:   }
889:   PetscCall(MatSetUp(C));
890:   C->ops->productnumeric = MatProductNumeric_SchurComplement_Dense;
891:   PetscFunctionReturn(PETSC_SUCCESS);
892: }

894: static PetscErrorCode MatProductSetFromOptions_Dense_AB(Mat C)
895: {
896:   PetscFunctionBegin;
897:   C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
898:   PetscFunctionReturn(PETSC_SUCCESS);
899: }

901: static PetscErrorCode MatProductSetFromOptions_SchurComplement_Dense(Mat C)
902: {
903:   Mat_Product *product = C->product;

905:   PetscFunctionBegin;
906:   PetscCheck(product->type == MATPRODUCT_AB, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Not for product type %s", MatProductTypes[product->type]);
907:   PetscCall(MatProductSetFromOptions_Dense_AB(C));
908:   PetscFunctionReturn(PETSC_SUCCESS);
909: }

911: /*MC
912:   MATSCHURCOMPLEMENT -  "schurcomplement" - Matrix type that behaves like the Schur complement of a matrix.

914:   Level: intermediate

916: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatType`, `MatCreateSchurComplement()`, `MatSchurComplementComputeExplicitOperator()`,
917:           `MatSchurComplementGetSubMatrices()`, `MatSchurComplementGetKSP()`
918: M*/
919: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
920: {
921:   Mat_SchurComplement *Na;

923:   PetscFunctionBegin;
924:   PetscCall(PetscNew(&Na));
925:   N->data = (void *)Na;

927:   N->ops->destroy        = MatDestroy_SchurComplement;
928:   N->ops->getvecs        = MatCreateVecs_SchurComplement;
929:   N->ops->view           = MatView_SchurComplement;
930:   N->ops->mult           = MatMult_SchurComplement;
931:   N->ops->multtranspose  = MatMultTranspose_SchurComplement;
932:   N->ops->multadd        = MatMultAdd_SchurComplement;
933:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
934:   N->assembled           = PETSC_FALSE;
935:   N->preallocated        = PETSC_FALSE;

937:   PetscCall(KSPCreate(PetscObjectComm((PetscObject)N), &Na->ksp));
938:   PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATSCHURCOMPLEMENT));
939:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", MatProductSetFromOptions_SchurComplement_Dense));
940:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", MatProductSetFromOptions_SchurComplement_Dense));
941:   PetscFunctionReturn(PETSC_SUCCESS);
942: }