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 of A00\n"));
 38:   PetscCall(PetscViewerASCIIPushTab(viewer));
 39:   PetscCall(KSPView(Na->ksp, viewer));
 40:   PetscCall(PetscViewerASCIIPopTab(viewer));
 41:   PetscCall(PetscViewerASCIIPrintf(viewer, "A01\n"));
 42:   PetscCall(PetscViewerASCIIPushTab(viewer));
 43:   PetscCall(MatView(Na->B, viewer));
 44:   PetscCall(PetscViewerASCIIPopTab(viewer));
 45:   PetscFunctionReturn(PETSC_SUCCESS);
 46: }

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

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

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

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

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

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

107: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N, PetscOptionItems *PetscOptionsObject)
108: {
109:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

111:   PetscFunctionBegin;
112:   PetscOptionsHeadBegin(PetscOptionsObject, "MatSchurComplementOptions");
113:   Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
114:   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,
115:                              (PetscEnum *)&Na->ainvtype, NULL));
116:   PetscOptionsHeadEnd();
117:   PetscCall(KSPSetFromOptions(Na->ksp));
118:   PetscFunctionReturn(PETSC_SUCCESS);
119: }

121: PetscErrorCode MatDestroy_SchurComplement(Mat N)
122: {
123:   Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;

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

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

143:    Collective

145:    Input Parameters:
146: +   A00,A01,A10,A11  - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
147: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}

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 Note:
168:     The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
169:     remove redundancy and be clearer and simpler.

171: .seealso: [](chapter_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,A01,A10,A11  - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
192: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

194:    Level: intermediate

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

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

202:      A00 and A11 must be square matrices.

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

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

218:   PetscFunctionBegin;
219:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
220:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
221:   PetscCheck(!S->assembled, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Use MatSchurComplementUpdateSubMatrices() for already used matrix");
226:   PetscCheckSameComm(A00, 2, Ap00, 3);
227:   PetscCheckSameComm(A00, 2, A01, 4);
228:   PetscCheckSameComm(A00, 2, A10, 5);
229:   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);
230:   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);
231:   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);
232:   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);
233:   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);
234:   if (A11) {
236:     PetscCheckSameComm(A00, 2, A11, 6);
237:     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);
238:   }

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

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

260:   Not Collective

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

265:   Output Parameter:
266: . ksp - the linear solver object

268:   Options Database Key:
269: . -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.

271:   Level: intermediate

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

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

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

293:   Not Collective

295:   Input Parameters:
296: + S   - matrix created with `MatCreateSchurComplement()`
297: - ksp - the linear solver object

299:   Level: developer

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

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

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

325: /*@
326:       MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices

328:    Collective

330:    Input Parameters:
331: +   S                - matrix obtained with `MatCreateSchurComplement()` (or `MatSchurSetSubMatrices()`) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
332: .   A00,A01,A10,A11  - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
333: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

335:    Level: intermediate

337:    Notes:
338:      All four matrices must have the same MPI communicator

340:      A00 and  A11 must be square matrices

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

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

347:    Developer Note:
348:      This code is almost identical to `MatSchurComplementSetSubMatrices()`. The API should be refactored.

350: .seealso: [](chapter_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`
351: @*/
352: PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
353: {
354:   Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
355:   PetscBool            isschur;

357:   PetscFunctionBegin;
359:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
360:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
361:   PetscCheck(S->assembled, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Use MatSchurComplementSetSubMatrices() for a new matrix");
366:   PetscCheckSameComm(A00, 2, Ap00, 3);
367:   PetscCheckSameComm(A00, 2, A01, 4);
368:   PetscCheckSameComm(A00, 2, A10, 5);
369:   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);
370:   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);
371:   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);
372:   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);
373:   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);
374:   if (A11) {
376:     PetscCheckSameComm(A00, 2, A11, 6);
377:     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);
378:   }

380:   PetscCall(PetscObjectReference((PetscObject)A00));
381:   PetscCall(PetscObjectReference((PetscObject)Ap00));
382:   PetscCall(PetscObjectReference((PetscObject)A01));
383:   PetscCall(PetscObjectReference((PetscObject)A10));
384:   if (A11) PetscCall(PetscObjectReference((PetscObject)A11));

386:   PetscCall(MatDestroy(&Na->A));
387:   PetscCall(MatDestroy(&Na->Ap));
388:   PetscCall(MatDestroy(&Na->B));
389:   PetscCall(MatDestroy(&Na->C));
390:   PetscCall(MatDestroy(&Na->D));

392:   Na->A  = A00;
393:   Na->Ap = Ap00;
394:   Na->B  = A01;
395:   Na->C  = A10;
396:   Na->D  = A11;

398:   PetscCall(KSPSetOperators(Na->ksp, A00, Ap00));
399:   PetscFunctionReturn(PETSC_SUCCESS);
400: }

402: /*@C
403:   MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement

405:   Collective

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

410:   Output Parameters:
411: + A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
412: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
413: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
414: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
415: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

417:   Level: intermediate

419:   Note:
420:   A11 is optional, and thus can be NULL.  The reference counts of the submatrices are not increased before they are returned and the matrices should not be modified or destroyed.

422: .seealso: [](chapter_ksp), `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatSchurComplementUpdateSubMatrices()`
423: @*/
424: PetscErrorCode MatSchurComplementGetSubMatrices(Mat S, Mat *A00, Mat *Ap00, Mat *A01, Mat *A10, Mat *A11)
425: {
426:   Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
427:   PetscBool            flg;

429:   PetscFunctionBegin;
431:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &flg));
432:   PetscCheck(flg, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)S)->type_name);
433:   if (A00) *A00 = Na->A;
434:   if (Ap00) *Ap00 = Na->Ap;
435:   if (A01) *A01 = Na->B;
436:   if (A10) *A10 = Na->C;
437:   if (A11) *A11 = Na->D;
438:   PetscFunctionReturn(PETSC_SUCCESS);
439: }

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

443: /*@
444:   MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly

446:   Collective

448:   Input Parameter:
449: . M - the matrix obtained with `MatCreateSchurComplement()`

451:   Output Parameter:
452: . S - the Schur complement matrix

454:   Notes:
455:     This can be expensive, so it is mainly for testing

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

459:   Level: advanced

461: .seealso: [](chapter_ksp), `MatCreateSchurComplement()`, `MatSchurComplementUpdate()`, `MatSchurComplementGetPmat()`
462: @*/
463: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat A, Mat *S)
464: {
465:   Mat       B, C, D, E = NULL, Bd, AinvBd;
466:   KSP       ksp;
467:   PetscInt  n, N, m, M;
468:   PetscBool flg = PETSC_FALSE, set, symm;

470:   PetscFunctionBegin;
471:   PetscCall(MatSchurComplementGetSubMatrices(A, NULL, NULL, &B, &C, &D));
472:   PetscCall(MatSchurComplementGetKSP(A, &ksp));
473:   PetscCall(KSPSetUp(ksp));
474:   PetscCall(MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd));
475:   PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd));
476:   PetscCall(KSPMatSolve(ksp, Bd, AinvBd));
477:   PetscCall(MatDestroy(&Bd));
478:   PetscCall(MatChop(AinvBd, PETSC_SMALL));
479:   if (D) {
480:     PetscCall(MatGetLocalSize(D, &m, &n));
481:     PetscCall(MatGetSize(D, &M, &N));
482:     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), m, n, M, N, NULL, S));
483:   }
484:   PetscCall(MatMatMult(C, AinvBd, D ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, S));
485:   PetscCall(MatDestroy(&AinvBd));
486:   if (D) {
487:     PetscCall(PetscObjectTypeCompareAny((PetscObject)D, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
488:     if (flg) {
489:       PetscCall(MatIsSymmetricKnown(A, &set, &symm));
490:       if (!set || !symm) PetscCall(MatConvert(D, MATBAIJ, MAT_INITIAL_MATRIX, &E)); /* convert the (1,1) block to nonsymmetric storage for MatAXPY() */
491:     }
492:     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 */
493:   }
494:   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 */
495:   PetscCall(MatScale(*S, -1.0));
496:   PetscCall(MatDestroy(&E));
497:   PetscFunctionReturn(PETSC_SUCCESS);
498: }

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

507:   PetscFunctionBegin;
517:   if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);

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

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

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

556:     Collective

558:     Input Parameters:
559: +   A      - matrix in which the complement is to be taken
560: .   isrow0 - rows to eliminate
561: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
562: .   isrow1 - rows in which the Schur complement is formed
563: .   iscol1 - columns in which the Schur complement is formed
564: .   mreuse - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in S
565: .   ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
566:                        `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
567: -   preuse - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in Sp

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

573:     Level: advanced

575:     Notes:
576:     Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
577:     application-specific information.

579:     Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
580:     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
581:     "MatGetSchurComplement_C" to their function.  If their function needs to fall back to the default implementation, it
582:     should call `MatGetSchurComplement_Basic()`.

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

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

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

591:     Developer Note:
592:     The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
593:     remove redundancy and be clearer and simpler.

595: .seealso: [](chapter_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatCreateSchurComplement()`, `MatSchurComplementAinvType`
596: @*/
597: PetscErrorCode MatGetSchurComplement(Mat A, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
598: {
599:   PetscErrorCode (*f)(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatReuse, Mat *) = NULL;

601:   PetscFunctionBegin;
613:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
614:   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. */
615:     PetscCall(PetscObjectQueryFunction((PetscObject)*S, "MatGetSchurComplement_C", &f));
616:   }
617:   if (f) PetscCall((*f)(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, preuse, Sp));
618:   else PetscCall(MatGetSchurComplement_Basic(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, ainvtype, preuse, Sp));
619:   PetscFunctionReturn(PETSC_SUCCESS);
620: }

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

625:     Not collective.

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

632:     Options Database Key:
633:     -mat_schur_complement_ainv_type diag | lump | blockdiag | full

635:     Level: advanced

637: .seealso: [](chapter_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementGetAinvType()`
638: @*/
639: PetscErrorCode MatSchurComplementSetAinvType(Mat S, MatSchurComplementAinvType ainvtype)
640: {
641:   PetscBool            isschur;
642:   Mat_SchurComplement *schur;

644:   PetscFunctionBegin;
646:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
647:   if (!isschur) PetscFunctionReturn(PETSC_SUCCESS);
649:   schur = (Mat_SchurComplement *)S->data;
650:   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);
651:   schur->ainvtype = ainvtype;
652:   PetscFunctionReturn(PETSC_SUCCESS);
653: }

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

658:     Not collective.

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

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

667:     Level: advanced

669: .seealso: [](chapter_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementSetAinvType()`
670: @*/
671: PetscErrorCode MatSchurComplementGetAinvType(Mat S, MatSchurComplementAinvType *ainvtype)
672: {
673:   PetscBool            isschur;
674:   Mat_SchurComplement *schur;

676:   PetscFunctionBegin;
678:   PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur));
679:   PetscCheck(isschur, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)S)->type_name);
680:   schur = (Mat_SchurComplement *)S->data;
681:   if (ainvtype) *ainvtype = schur->ainvtype;
682:   PetscFunctionReturn(PETSC_SUCCESS);
683: }

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

689:     Collective on A00

691:     Input Parameters:
692: +   A00      - the upper-left part of the original matrix A = [A00 A01; A10 A11]
693: .   A01      - (optional) the upper-right part of the original matrix A = [A00 A01; A10 A11]
694: .   A10      - (optional) the lower-left part of the original matrix A = [A00 A01; A10 A11]
695: .   A11      - (optional) the lower-right part of the original matrix A = [A00 A01; A10 A11]
696: .   ainvtype - type of approximation for DIAGFORM(A00) used when forming Sp = A11 - A10 inv(DIAGFORM(A00)) A01. See MatSchurComplementAinvType.
697: -   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

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

702:     Level: advanced

704: .seealso: [](chapter_ksp), `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementAinvType`
705: @*/
706: PetscErrorCode MatCreateSchurComplementPmat(Mat A00, Mat A01, Mat A10, Mat A11, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
707: {
708:   PetscInt N00;

710:   PetscFunctionBegin;
711:   /* 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. */
712:   /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
714:   if (preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);

716:   /* A zero size A00 or empty A01 or A10 imply S = A11. */
717:   PetscCall(MatGetSize(A00, &N00, NULL));
718:   if (!A01 || !A10 || !N00) {
719:     if (preuse == MAT_INITIAL_MATRIX) {
720:       PetscCall(MatDuplicate(A11, MAT_COPY_VALUES, Sp));
721:     } else { /* MAT_REUSE_MATRIX */
722:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
723:       PetscCall(MatCopy(A11, *Sp, DIFFERENT_NONZERO_PATTERN));
724:     }
725:   } else {
726:     Mat AdB;
727:     Vec diag;

729:     if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
730:       PetscCall(MatDuplicate(A01, MAT_COPY_VALUES, &AdB));
731:       PetscCall(MatCreateVecs(A00, &diag, NULL));
732:       if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
733:         PetscCall(MatGetRowSum(A00, diag));
734:       } else {
735:         PetscCall(MatGetDiagonal(A00, diag));
736:       }
737:       PetscCall(VecReciprocal(diag));
738:       PetscCall(MatDiagonalScale(AdB, diag, NULL));
739:       PetscCall(VecDestroy(&diag));
740:     } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
741:       Mat      A00_inv;
742:       MatType  type;
743:       MPI_Comm comm;

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

768: PetscErrorCode MatSchurComplementGetPmat_Basic(Mat S, MatReuse preuse, Mat *Sp)
769: {
770:   Mat                  A, B, C, D;
771:   Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;

773:   PetscFunctionBegin;
774:   if (preuse == MAT_IGNORE_MATRIX) PetscFunctionReturn(PETSC_SUCCESS);
775:   PetscCall(MatSchurComplementGetSubMatrices(S, &A, NULL, &B, &C, &D));
776:   PetscCheck(A, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Schur complement component matrices unset");
777:   if (schur->ainvtype != MAT_SCHUR_COMPLEMENT_AINV_FULL) PetscCall(MatCreateSchurComplementPmat(A, B, C, D, schur->ainvtype, preuse, Sp));
778:   else {
779:     if (preuse == MAT_REUSE_MATRIX) PetscCall(MatDestroy(Sp));
780:     PetscCall(MatSchurComplementComputeExplicitOperator(S, Sp));
781:   }
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

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

788:     Collective

790:     Input Parameters:
791: +   S      - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of A11 - A10 ksp(A00,Ap00) A01
792: -   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

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

797:     Level: advanced

799:     Notes:
800:     The approximation of Sp depends on the the argument passed to to `MatSchurComplementSetAinvType()`
801:     `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
802:     -mat_schur_complement_ainv_type <diag,lump,blockdiag,full>

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

809:     Developer Note:
810:     The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
811:     remove redundancy and be clearer and simpler.

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

815: .seealso: [](chapter_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetAinvType()`
816: @*/
817: PetscErrorCode MatSchurComplementGetPmat(Mat S, MatReuse preuse, Mat *Sp)
818: {
819:   PetscErrorCode (*f)(Mat, MatReuse, Mat *);

821:   PetscFunctionBegin;
825:   if (preuse != MAT_IGNORE_MATRIX) {
827:     if (preuse == MAT_INITIAL_MATRIX) *Sp = NULL;
829:   }
830:   PetscCheck(!S->factortype, PetscObjectComm((PetscObject)S), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

832:   PetscCall(PetscObjectQueryFunction((PetscObject)S, "MatSchurComplementGetPmat_C", &f));
833:   if (f) PetscCall((*f)(S, preuse, Sp));
834:   else PetscCall(MatSchurComplementGetPmat_Basic(S, preuse, Sp));
835:   PetscFunctionReturn(PETSC_SUCCESS);
836: }

838: static PetscErrorCode MatProductNumeric_SchurComplement_Dense(Mat C)
839: {
840:   Mat_Product         *product = C->product;
841:   Mat_SchurComplement *Na      = (Mat_SchurComplement *)product->A->data;
842:   Mat                  work1, work2;
843:   PetscScalar         *v;
844:   PetscInt             lda;

846:   PetscFunctionBegin;
847:   PetscCall(MatMatMult(Na->B, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1));
848:   PetscCall(MatDuplicate(work1, MAT_DO_NOT_COPY_VALUES, &work2));
849:   PetscCall(KSPMatSolve(Na->ksp, work1, work2));
850:   PetscCall(MatDestroy(&work1));
851:   PetscCall(MatDenseGetArrayWrite(C, &v));
852:   PetscCall(MatDenseGetLDA(C, &lda));
853:   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)C), C->rmap->n, C->cmap->n, C->rmap->N, C->cmap->N, v, &work1));
854:   PetscCall(MatDenseSetLDA(work1, lda));
855:   PetscCall(MatMatMult(Na->C, work2, MAT_REUSE_MATRIX, PETSC_DEFAULT, &work1));
856:   PetscCall(MatDenseRestoreArrayWrite(C, &v));
857:   PetscCall(MatDestroy(&work2));
858:   PetscCall(MatDestroy(&work1));
859:   if (Na->D) {
860:     PetscCall(MatMatMult(Na->D, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1));
861:     PetscCall(MatAYPX(C, -1.0, work1, SAME_NONZERO_PATTERN));
862:     PetscCall(MatDestroy(&work1));
863:   } else PetscCall(MatScale(C, -1.0));
864:   PetscFunctionReturn(PETSC_SUCCESS);
865: }

867: static PetscErrorCode MatProductSymbolic_SchurComplement_Dense(Mat C)
868: {
869:   Mat_Product *product = C->product;
870:   Mat          A = product->A, B = product->B;
871:   PetscInt     m = A->rmap->n, n = B->cmap->n, M = A->rmap->N, N = B->cmap->N;
872:   PetscBool    flg;

874:   PetscFunctionBegin;
875:   PetscCall(MatSetSizes(C, m, n, M, N));
876:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)C, &flg, MATSEQDENSE, MATMPIDENSE, ""));
877:   if (!flg) {
878:     PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
879:     C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
880:   }
881:   PetscCall(MatSetUp(C));
882:   C->ops->productnumeric = MatProductNumeric_SchurComplement_Dense;
883:   PetscFunctionReturn(PETSC_SUCCESS);
884: }

886: static PetscErrorCode MatProductSetFromOptions_Dense_AB(Mat C)
887: {
888:   PetscFunctionBegin;
889:   C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
890:   PetscFunctionReturn(PETSC_SUCCESS);
891: }

893: static PetscErrorCode MatProductSetFromOptions_SchurComplement_Dense(Mat C)
894: {
895:   Mat_Product *product = C->product;

897:   PetscFunctionBegin;
898:   PetscCheck(product->type == MATPRODUCT_AB, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Not for product type %s", MatProductTypes[product->type]);
899:   PetscCall(MatProductSetFromOptions_Dense_AB(C));
900:   PetscFunctionReturn(PETSC_SUCCESS);
901: }

903: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
904: {
905:   Mat_SchurComplement *Na;

907:   PetscFunctionBegin;
908:   PetscCall(PetscNew(&Na));
909:   N->data = (void *)Na;

911:   N->ops->destroy        = MatDestroy_SchurComplement;
912:   N->ops->getvecs        = MatCreateVecs_SchurComplement;
913:   N->ops->view           = MatView_SchurComplement;
914:   N->ops->mult           = MatMult_SchurComplement;
915:   N->ops->multtranspose  = MatMultTranspose_SchurComplement;
916:   N->ops->multadd        = MatMultAdd_SchurComplement;
917:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
918:   N->assembled           = PETSC_FALSE;
919:   N->preallocated        = PETSC_FALSE;

921:   PetscCall(KSPCreate(PetscObjectComm((PetscObject)N), &Na->ksp));
922:   PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATSCHURCOMPLEMENT));
923:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", MatProductSetFromOptions_SchurComplement_Dense));
924:   PetscCall(PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", MatProductSetFromOptions_SchurComplement_Dense));
925:   PetscFunctionReturn(PETSC_SUCCESS);
926: }