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: if (Na->D) {
10: MatCreateVecs(Na->D, right, left);
11: return 0;
12: }
13: if (right) MatCreateVecs(Na->B, right, NULL);
14: if (left) MatCreateVecs(Na->C, NULL, left);
15: return 0;
16: }
18: PetscErrorCode MatView_SchurComplement(Mat N, PetscViewer viewer)
19: {
20: Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;
22: PetscViewerASCIIPrintf(viewer, "Schur complement A11 - A10 inv(A00) A01\n");
23: if (Na->D) {
24: PetscViewerASCIIPrintf(viewer, "A11\n");
25: PetscViewerASCIIPushTab(viewer);
26: MatView(Na->D, viewer);
27: PetscViewerASCIIPopTab(viewer);
28: } else {
29: PetscViewerASCIIPrintf(viewer, "A11 = 0\n");
30: }
31: PetscViewerASCIIPrintf(viewer, "A10\n");
32: PetscViewerASCIIPushTab(viewer);
33: MatView(Na->C, viewer);
34: PetscViewerASCIIPopTab(viewer);
35: PetscViewerASCIIPrintf(viewer, "KSP of A00\n");
36: PetscViewerASCIIPushTab(viewer);
37: KSPView(Na->ksp, viewer);
38: PetscViewerASCIIPopTab(viewer);
39: PetscViewerASCIIPrintf(viewer, "A01\n");
40: PetscViewerASCIIPushTab(viewer);
41: MatView(Na->B, viewer);
42: PetscViewerASCIIPopTab(viewer);
43: return 0;
44: }
46: /*
47: A11^T - A01^T ksptrans(A00,Ap00) A10^T
48: */
49: PetscErrorCode MatMultTranspose_SchurComplement(Mat N, Vec x, Vec y)
50: {
51: Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;
53: if (!Na->work1) MatCreateVecs(Na->A, &Na->work1, NULL);
54: if (!Na->work2) MatCreateVecs(Na->A, &Na->work2, NULL);
55: MatMultTranspose(Na->C, x, Na->work1);
56: KSPSolveTranspose(Na->ksp, Na->work1, Na->work2);
57: MatMultTranspose(Na->B, Na->work2, y);
58: VecScale(y, -1.0);
59: if (Na->D) MatMultTransposeAdd(Na->D, x, y, y);
60: return 0;
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: if (!Na->work1) MatCreateVecs(Na->A, &Na->work1, NULL);
71: if (!Na->work2) MatCreateVecs(Na->A, &Na->work2, NULL);
72: MatMult(Na->B, x, Na->work1);
73: KSPSolve(Na->ksp, Na->work1, Na->work2);
74: MatMult(Na->C, Na->work2, y);
75: VecScale(y, -1.0);
76: if (Na->D) MatMultAdd(Na->D, x, y, y);
77: return 0;
78: }
80: /*
81: A11 - A10 ksp(A00,Ap00) A01
82: */
83: PetscErrorCode MatMultAdd_SchurComplement(Mat N, Vec x, Vec y, Vec z)
84: {
85: Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;
87: if (!Na->work1) MatCreateVecs(Na->A, &Na->work1, NULL);
88: if (!Na->work2) MatCreateVecs(Na->A, &Na->work2, NULL);
89: MatMult(Na->B, x, Na->work1);
90: KSPSolve(Na->ksp, Na->work1, Na->work2);
91: if (y == z) {
92: VecScale(Na->work2, -1.0);
93: MatMultAdd(Na->C, Na->work2, z, z);
94: } else {
95: MatMult(Na->C, Na->work2, z);
96: VecAYPX(z, -1.0, y);
97: }
98: if (Na->D) MatMultAdd(Na->D, x, z, z);
99: return 0;
100: }
102: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N, PetscOptionItems *PetscOptionsObject)
103: {
104: Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;
106: PetscOptionsHeadBegin(PetscOptionsObject, "MatSchurComplementOptions");
107: Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
108: 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,
109: (PetscEnum *)&Na->ainvtype, NULL));
110: PetscOptionsHeadEnd();
111: KSPSetFromOptions(Na->ksp);
112: return 0;
113: }
115: PetscErrorCode MatDestroy_SchurComplement(Mat N)
116: {
117: Mat_SchurComplement *Na = (Mat_SchurComplement *)N->data;
119: MatDestroy(&Na->A);
120: MatDestroy(&Na->Ap);
121: MatDestroy(&Na->B);
122: MatDestroy(&Na->C);
123: MatDestroy(&Na->D);
124: VecDestroy(&Na->work1);
125: VecDestroy(&Na->work2);
126: KSPDestroy(&Na->ksp);
127: PetscFree(N->data);
128: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", NULL);
129: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", NULL);
130: return 0;
131: }
133: /*@
134: MatCreateSchurComplement - Creates a new `Mat` that behaves like the Schur complement of a matrix
136: Collective
138: Input Parameters:
139: + A00,A01,A10,A11 - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
140: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
142: Output Parameter:
143: . S - the matrix that behaves as the Schur complement S = A11 - A10 ksp(A00,Ap00) A01
145: Level: intermediate
147: Notes:
148: The Schur complement is NOT explicitly formed! Rather, this function returns a virtual Schur complement
149: that can compute the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
150: for Schur complement S and a `KSP` solver to approximate the action of A^{-1}.
152: All four matrices must have the same MPI communicator.
154: A00 and A11 must be square matrices.
156: `MatGetSchurComplement()` takes as arguments the index sets for the submatrices and returns both the virtual Schur complement (what this returns) plus
157: a sparse approximation to the Schur complement (useful for building a preconditioner for the Schur complement) which can be obtained from this
158: matrix with `MatSchurComplementGetPmat()`
160: Developer Note:
161: The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
162: remove redundancy and be clearer and simpler.
164: .seealso: [](chapter_ksp), `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatGetSchurComplement()`,
165: `MatSchurComplementGetPmat()`, `MatSchurComplementSetSubMatrices()`
166: @*/
167: PetscErrorCode MatCreateSchurComplement(Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11, Mat *S)
168: {
169: KSPInitializePackage();
170: MatCreate(PetscObjectComm((PetscObject)A00), S);
171: MatSetType(*S, MATSCHURCOMPLEMENT);
172: MatSchurComplementSetSubMatrices(*S, A00, Ap00, A01, A10, A11);
173: return 0;
174: }
176: /*@
177: MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement
179: Collective
181: Input Parameters:
182: + S - matrix obtained with `MatSetType`(S,`MATSCHURCOMPLEMENT`)
183: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
184: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
186: Level: intermediate
188: Notes:
189: The Schur complement is NOT explicitly formed! Rather, this
190: object performs the matrix-vector product of the Schur complement by using formula S = A11 - A10 ksp(A00,Ap00) A01
192: All four matrices must have the same MPI communicator.
194: A00 and A11 must be square matrices.
196: This is to be used in the context of code such as
197: .vb
198: MatSetType(S,MATSCHURCOMPLEMENT);
199: MatSchurComplementSetSubMatrices(S,...);
200: .ve
201: while `MatSchurComplementUpdateSubMatrices()` should only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`
203: .seealso: [](chapter_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`
204: @*/
205: PetscErrorCode MatSchurComplementSetSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
206: {
207: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
208: PetscBool isschur;
210: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
211: if (!isschur) return 0;
225: if (A11) {
229: }
231: MatSetSizes(S, A10->rmap->n, A01->cmap->n, A10->rmap->N, A01->cmap->N);
232: PetscObjectReference((PetscObject)A00);
233: PetscObjectReference((PetscObject)Ap00);
234: PetscObjectReference((PetscObject)A01);
235: PetscObjectReference((PetscObject)A10);
236: Na->A = A00;
237: Na->Ap = Ap00;
238: Na->B = A01;
239: Na->C = A10;
240: Na->D = A11;
241: if (A11) PetscObjectReference((PetscObject)A11);
242: MatSetUp(S);
243: KSPSetOperators(Na->ksp, A00, Ap00);
244: S->assembled = PETSC_TRUE;
245: return 0;
246: }
248: /*@
249: MatSchurComplementGetKSP - Gets the `KSP` object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
251: Not Collective
253: Input Parameter:
254: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
256: Output Parameter:
257: . ksp - the linear solver object
259: Options Database Key:
260: . -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.
262: Level: intermediate
264: .seealso: [](chapter_ksp), `Mat`, `MatSchurComplementSetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`
265: @*/
266: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
267: {
268: Mat_SchurComplement *Na;
269: PetscBool isschur;
272: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
275: Na = (Mat_SchurComplement *)S->data;
276: *ksp = Na->ksp;
277: return 0;
278: }
280: /*@
281: MatSchurComplementSetKSP - Sets the `KSP` object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
283: Not Collective
285: Input Parameters:
286: + S - matrix created with `MatCreateSchurComplement()`
287: - ksp - the linear solver object
289: Level: developer
291: Developer Note:
292: This is used in `PCFIELDSPLIT` to reuse the 0-split `KSP` to implement ksp(A00,Ap00) in S.
293: The `KSP` operators are overwritten with A00 and Ap00 currently set in S.
295: .seealso: [](chapter_ksp), `Mat`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MATSCHURCOMPLEMENT`
296: @*/
297: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
298: {
299: Mat_SchurComplement *Na;
300: PetscBool isschur;
303: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
304: if (!isschur) return 0;
306: Na = (Mat_SchurComplement *)S->data;
307: PetscObjectReference((PetscObject)ksp);
308: KSPDestroy(&Na->ksp);
309: Na->ksp = ksp;
310: KSPSetOperators(Na->ksp, Na->A, Na->Ap);
311: return 0;
312: }
314: /*@
315: MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices
317: Collective
319: Input Parameters:
320: + S - matrix obtained with `MatCreateSchurComplement()` (or `MatSchurSetSubMatrices()`) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
321: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
322: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
324: Level: intermediate
326: Notes:
327: All four matrices must have the same MPI communicator
329: A00 and A11 must be square matrices
331: All of the matrices provided must have the same sizes as was used with `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`
332: though they need not be the same matrices.
334: This can only be called after `MatCreateSchurComplement()` or `MatSchurComplementSetSubMatrices()`, it cannot be called immediately after `MatSetType`(S,`MATSCHURCOMPLEMENT`);
336: Developer Note:
337: This code is almost identical to `MatSchurComplementSetSubMatrices()`. The API should be refactored.
339: .seealso: [](chapter_ksp), `Mat`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`
340: @*/
341: PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
342: {
343: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
344: PetscBool isschur;
347: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
348: if (!isschur) return 0;
362: if (A11) {
366: }
368: PetscObjectReference((PetscObject)A00);
369: PetscObjectReference((PetscObject)Ap00);
370: PetscObjectReference((PetscObject)A01);
371: PetscObjectReference((PetscObject)A10);
372: if (A11) PetscObjectReference((PetscObject)A11);
374: MatDestroy(&Na->A);
375: MatDestroy(&Na->Ap);
376: MatDestroy(&Na->B);
377: MatDestroy(&Na->C);
378: MatDestroy(&Na->D);
380: Na->A = A00;
381: Na->Ap = Ap00;
382: Na->B = A01;
383: Na->C = A10;
384: Na->D = A11;
386: KSPSetOperators(Na->ksp, A00, Ap00);
387: return 0;
388: }
390: /*@C
391: MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement
393: Collective
395: Input Parameter:
396: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
398: Output Parameters:
399: + A00 - the upper-left block of the original matrix A = [A00 A01; A10 A11]
400: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
401: . A01 - the upper-right block of the original matrix A = [A00 A01; A10 A11]
402: . A10 - the lower-left block of the original matrix A = [A00 A01; A10 A11]
403: - A11 - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]
405: Level: intermediate
407: Note:
408: 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.
410: .seealso: [](chapter_ksp), `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatSchurComplementUpdateSubMatrices()`
411: @*/
412: PetscErrorCode MatSchurComplementGetSubMatrices(Mat S, Mat *A00, Mat *Ap00, Mat *A01, Mat *A10, Mat *A11)
413: {
414: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
415: PetscBool flg;
418: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &flg);
420: if (A00) *A00 = Na->A;
421: if (Ap00) *Ap00 = Na->Ap;
422: if (A01) *A01 = Na->B;
423: if (A10) *A10 = Na->C;
424: if (A11) *A11 = Na->D;
425: return 0;
426: }
428: #include <petsc/private/kspimpl.h>
430: /*@
431: MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly
433: Collective
435: Input Parameter:
436: . M - the matrix obtained with `MatCreateSchurComplement()`
438: Output Parameter:
439: . S - the Schur complement matrix
441: Notes:
442: This can be expensive, so it is mainly for testing
444: Use `MatSchurComplementGetPmat()` to get a sparse approximation for the Schur complement suitable for use in building a preconditioner
446: Level: advanced
448: .seealso: [](chapter_ksp), `MatCreateSchurComplement()`, `MatSchurComplementUpdate()`, `MatSchurComplementGetPmat()`
449: @*/
450: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat A, Mat *S)
451: {
452: Mat B, C, D, E = NULL, Bd, AinvBd;
453: KSP ksp;
454: PetscInt n, N, m, M;
455: PetscBool flg = PETSC_FALSE, set, symm;
457: MatSchurComplementGetSubMatrices(A, NULL, NULL, &B, &C, &D);
458: MatSchurComplementGetKSP(A, &ksp);
459: KSPSetUp(ksp);
460: MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
461: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
462: KSPMatSolve(ksp, Bd, AinvBd);
463: MatDestroy(&Bd);
464: MatChop(AinvBd, PETSC_SMALL);
465: if (D) {
466: MatGetLocalSize(D, &m, &n);
467: MatGetSize(D, &M, &N);
468: MatCreateDense(PetscObjectComm((PetscObject)A), m, n, M, N, NULL, S);
469: }
470: MatMatMult(C, AinvBd, D ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, S);
471: MatDestroy(&AinvBd);
472: if (D) {
473: PetscObjectTypeCompareAny((PetscObject)D, &flg, MATSEQSBAIJ, MATMPISBAIJ, "");
474: if (flg) {
475: MatIsSymmetricKnown(A, &set, &symm);
476: if (!set || !symm) MatConvert(D, MATBAIJ, MAT_INITIAL_MATRIX, &E); /* convert the (1,1) block to nonsymmetric storage for MatAXPY() */
477: }
478: 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 */
479: }
480: 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 */
481: MatScale(*S, -1.0);
482: MatDestroy(&E);
483: return 0;
484: }
486: /* Developer Notes:
487: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
488: PetscErrorCode MatGetSchurComplement_Basic(Mat mat, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
489: {
490: Mat A = NULL, Ap = NULL, B = NULL, C = NULL, D = NULL;
491: MatReuse reuse;
502: if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) return 0;
508: reuse = MAT_INITIAL_MATRIX;
509: if (mreuse == MAT_REUSE_MATRIX) {
510: MatSchurComplementGetSubMatrices(*S, &A, &Ap, &B, &C, &D);
513: MatDestroy(&Ap); /* get rid of extra reference */
514: reuse = MAT_REUSE_MATRIX;
515: }
516: MatCreateSubMatrix(mat, isrow0, iscol0, reuse, &A);
517: MatCreateSubMatrix(mat, isrow0, iscol1, reuse, &B);
518: MatCreateSubMatrix(mat, isrow1, iscol0, reuse, &C);
519: MatCreateSubMatrix(mat, isrow1, iscol1, reuse, &D);
520: switch (mreuse) {
521: case MAT_INITIAL_MATRIX:
522: MatCreateSchurComplement(A, A, B, C, D, S);
523: break;
524: case MAT_REUSE_MATRIX:
525: MatSchurComplementUpdateSubMatrices(*S, A, A, B, C, D);
526: break;
527: default:
529: }
530: if (preuse != MAT_IGNORE_MATRIX) MatCreateSchurComplementPmat(A, B, C, D, ainvtype, preuse, Sp);
531: MatDestroy(&A);
532: MatDestroy(&B);
533: MatDestroy(&C);
534: MatDestroy(&D);
535: return 0;
536: }
538: /*@
539: MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.
541: Collective
543: Input Parameters:
544: + A - matrix in which the complement is to be taken
545: . isrow0 - rows to eliminate
546: . iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
547: . isrow1 - rows in which the Schur complement is formed
548: . iscol1 - columns in which the Schur complement is formed
549: . mreuse - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in S
550: . ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
551: `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
552: - preuse - `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`, use `MAT_IGNORE_MATRIX` to put nothing in Sp
554: Output Parameters:
555: + S - exact Schur complement, often of type `MATSCHURCOMPLEMENT` which is difficult to use for preconditioning
556: - Sp - approximate Schur complement from which a preconditioner can be built A11 - A10 inv(DIAGFORM(A00)) A01
558: Level: advanced
560: Notes:
561: Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
562: application-specific information.
564: Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
565: 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
566: "MatGetSchurComplement_C" to their function. If their function needs to fall back to the default implementation, it
567: should call `MatGetSchurComplement_Basic()`.
569: `MatCreateSchurComplement()` takes as arguments the four submatrices and returns the virtual Schur complement (what this function returns in S).
571: `MatSchurComplementGetPmat()` takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).
573: In other words calling `MatCreateSchurComplement()` followed by `MatSchurComplementGetPmat()` produces the same output as this function but with slightly different
574: inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.
576: Developer Note:
577: The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
578: remove redundancy and be clearer and simpler.
580: .seealso: [](chapter_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatCreateSchurComplement()`, `MatSchurComplementAinvType`
581: @*/
582: PetscErrorCode MatGetSchurComplement(Mat A, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
583: {
584: PetscErrorCode (*f)(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatReuse, Mat *) = NULL;
598: 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. */
599: PetscObjectQueryFunction((PetscObject)*S, "MatGetSchurComplement_C", &f);
600: }
601: if (f) (*f)(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, preuse, Sp);
602: else MatGetSchurComplement_Basic(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, ainvtype, preuse, Sp);
603: return 0;
604: }
606: /*@
607: MatSchurComplementSetAinvType - set the type of approximation used for the inverse of the (0,0) block used in forming Sp in `MatSchurComplementGetPmat()`
609: Not collective.
611: Input Parameters:
612: + S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
613: - ainvtype - type of approximation to be used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
614: `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
616: Options Database Key:
617: -mat_schur_complement_ainv_type diag | lump | blockdiag | full
619: Level: advanced
621: .seealso: [](chapter_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementGetAinvType()`
622: @*/
623: PetscErrorCode MatSchurComplementSetAinvType(Mat S, MatSchurComplementAinvType ainvtype)
624: {
625: PetscBool isschur;
626: Mat_SchurComplement *schur;
629: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
630: if (!isschur) return 0;
632: schur = (Mat_SchurComplement *)S->data;
634: schur->ainvtype = ainvtype;
635: return 0;
636: }
638: /*@
639: MatSchurComplementGetAinvType - get the type of approximation for the inverse of the (0,0) block used in forming Sp in `MatSchurComplementGetPmat()`
641: Not collective.
643: Input Parameter:
644: . S - matrix obtained with `MatCreateSchurComplement()` (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
646: Output Parameter:
647: . ainvtype - type of approximation used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
648: `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
650: Level: advanced
652: .seealso: [](chapter_ksp), `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementSetAinvType()`
653: @*/
654: PetscErrorCode MatSchurComplementGetAinvType(Mat S, MatSchurComplementAinvType *ainvtype)
655: {
656: PetscBool isschur;
657: Mat_SchurComplement *schur;
660: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
662: schur = (Mat_SchurComplement *)S->data;
663: if (ainvtype) *ainvtype = schur->ainvtype;
664: return 0;
665: }
667: /*@
668: MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by explicitly assembling the sparse matrix
669: Sp = A11 - A10 inv(DIAGFORM(A00)) A01
671: Collective on A00
673: Input Parameters:
674: + A00 - the upper-left part of the original matrix A = [A00 A01; A10 A11]
675: . A01 - (optional) the upper-right part of the original matrix A = [A00 A01; A10 A11]
676: . A10 - (optional) the lower-left part of the original matrix A = [A00 A01; A10 A11]
677: . A11 - (optional) the lower-right part of the original matrix A = [A00 A01; A10 A11]
678: . ainvtype - type of approximation for DIAGFORM(A00) used when forming Sp = A11 - A10 inv(DIAGFORM(A00)) A01. See MatSchurComplementAinvType.
679: - 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
681: Output Parameter:
682: - Sp - approximate Schur complement suitable for preconditioning the true Schur complement S = A11 - A10 inv(A00) A01
684: Level: advanced
686: .seealso: [](chapter_ksp), `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementAinvType`
687: @*/
688: PetscErrorCode MatCreateSchurComplementPmat(Mat A00, Mat A01, Mat A10, Mat A11, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
689: {
690: PetscInt N00;
692: /* 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. */
693: /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
695: if (preuse == MAT_IGNORE_MATRIX) return 0;
697: /* A zero size A00 or empty A01 or A10 imply S = A11. */
698: MatGetSize(A00, &N00, NULL);
699: if (!A01 || !A10 || !N00) {
700: if (preuse == MAT_INITIAL_MATRIX) {
701: MatDuplicate(A11, MAT_COPY_VALUES, Sp);
702: } else { /* MAT_REUSE_MATRIX */
703: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
704: MatCopy(A11, *Sp, DIFFERENT_NONZERO_PATTERN);
705: }
706: } else {
707: Mat AdB;
708: Vec diag;
710: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
711: MatDuplicate(A01, MAT_COPY_VALUES, &AdB);
712: MatCreateVecs(A00, &diag, NULL);
713: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
714: MatGetRowSum(A00, diag);
715: } else {
716: MatGetDiagonal(A00, diag);
717: }
718: VecReciprocal(diag);
719: MatDiagonalScale(AdB, diag, NULL);
720: VecDestroy(&diag);
721: } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
722: Mat A00_inv;
723: MatType type;
724: MPI_Comm comm;
726: PetscObjectGetComm((PetscObject)A00, &comm);
727: MatGetType(A00, &type);
728: MatCreate(comm, &A00_inv);
729: MatSetType(A00_inv, type);
730: MatInvertBlockDiagonalMat(A00, A00_inv);
731: MatMatMult(A00_inv, A01, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &AdB);
732: MatDestroy(&A00_inv);
733: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown MatSchurComplementAinvType: %d", ainvtype);
734: /* Cannot really reuse Sp in MatMatMult() because of MatAYPX() -->
735: MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult() */
736: if (preuse == MAT_REUSE_MATRIX) MatDestroy(Sp);
737: MatMatMult(A10, AdB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, Sp);
738: if (!A11) {
739: MatScale(*Sp, -1.0);
740: } else {
741: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
742: MatAYPX(*Sp, -1, A11, DIFFERENT_NONZERO_PATTERN);
743: }
744: MatDestroy(&AdB);
745: }
746: return 0;
747: }
749: PetscErrorCode MatSchurComplementGetPmat_Basic(Mat S, MatReuse preuse, Mat *Sp)
750: {
751: Mat A, B, C, D;
752: Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;
754: if (preuse == MAT_IGNORE_MATRIX) return 0;
755: MatSchurComplementGetSubMatrices(S, &A, NULL, &B, &C, &D);
757: if (schur->ainvtype != MAT_SCHUR_COMPLEMENT_AINV_FULL) MatCreateSchurComplementPmat(A, B, C, D, schur->ainvtype, preuse, Sp);
758: else {
759: if (preuse == MAT_REUSE_MATRIX) MatDestroy(Sp);
760: MatSchurComplementComputeExplicitOperator(S, Sp);
761: }
762: return 0;
763: }
765: /*@
766: MatSchurComplementGetPmat - Obtain a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01
768: Collective
770: Input Parameters:
771: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of A11 - A10 ksp(A00,Ap00) A01
772: - 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
774: Output Parameter:
775: - Sp - approximate Schur complement suitable for preconditioning the exact Schur complement S = A11 - A10 inv(A00) A01
777: Level: advanced
779: Notes:
780: The approximation of Sp depends on the the argument passed to to `MatSchurComplementSetAinvType()`
781: `MAT_SCHUR_COMPLEMENT_AINV_DIAG`, `MAT_SCHUR_COMPLEMENT_AINV_LUMP`, `MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG`, or `MAT_SCHUR_COMPLEMENT_AINV_FULL`
782: -mat_schur_complement_ainv_type <diag,lump,blockdiag,full>
784: Sometimes users would like to provide problem-specific data in the Schur complement, usually only
785: for special row and column index sets. In that case, the user should call `PetscObjectComposeFunction()` to set
786: "MatSchurComplementGetPmat_C" to their function. If their function needs to fall back to the default implementation,
787: it should call `MatSchurComplementGetPmat_Basic()`.
789: Developer Note:
790: The API that includes `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementGetPmat()` should be refactored to
791: remove redundancy and be clearer and simpler.
793: This routine should be called `MatSchurComplementCreatePmat()`
795: .seealso: [](chapter_ksp), `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetAinvType()`
796: @*/
797: PetscErrorCode MatSchurComplementGetPmat(Mat S, MatReuse preuse, Mat *Sp)
798: {
799: PetscErrorCode (*f)(Mat, MatReuse, Mat *);
804: if (preuse != MAT_IGNORE_MATRIX) {
806: if (preuse == MAT_INITIAL_MATRIX) *Sp = NULL;
808: }
811: PetscObjectQueryFunction((PetscObject)S, "MatSchurComplementGetPmat_C", &f);
812: if (f) (*f)(S, preuse, Sp);
813: else MatSchurComplementGetPmat_Basic(S, preuse, Sp);
814: return 0;
815: }
817: static PetscErrorCode MatProductNumeric_SchurComplement_Dense(Mat C)
818: {
819: Mat_Product *product = C->product;
820: Mat_SchurComplement *Na = (Mat_SchurComplement *)product->A->data;
821: Mat work1, work2;
822: PetscScalar *v;
823: PetscInt lda;
825: MatMatMult(Na->B, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1);
826: MatDuplicate(work1, MAT_DO_NOT_COPY_VALUES, &work2);
827: KSPMatSolve(Na->ksp, work1, work2);
828: MatDestroy(&work1);
829: MatDenseGetArrayWrite(C, &v);
830: MatDenseGetLDA(C, &lda);
831: MatCreateDense(PetscObjectComm((PetscObject)C), C->rmap->n, C->cmap->n, C->rmap->N, C->cmap->N, v, &work1);
832: MatDenseSetLDA(work1, lda);
833: MatMatMult(Na->C, work2, MAT_REUSE_MATRIX, PETSC_DEFAULT, &work1);
834: MatDenseRestoreArrayWrite(C, &v);
835: MatDestroy(&work2);
836: MatDestroy(&work1);
837: if (Na->D) {
838: MatMatMult(Na->D, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1);
839: MatAYPX(C, -1.0, work1, SAME_NONZERO_PATTERN);
840: MatDestroy(&work1);
841: } else MatScale(C, -1.0);
842: return 0;
843: }
845: static PetscErrorCode MatProductSymbolic_SchurComplement_Dense(Mat C)
846: {
847: Mat_Product *product = C->product;
848: Mat A = product->A, B = product->B;
849: PetscInt m = A->rmap->n, n = B->cmap->n, M = A->rmap->N, N = B->cmap->N;
850: PetscBool flg;
852: MatSetSizes(C, m, n, M, N);
853: PetscObjectBaseTypeCompareAny((PetscObject)C, &flg, MATSEQDENSE, MATMPIDENSE, "");
854: if (!flg) {
855: MatSetType(C, ((PetscObject)B)->type_name);
856: C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
857: }
858: MatSetUp(C);
859: C->ops->productnumeric = MatProductNumeric_SchurComplement_Dense;
860: return 0;
861: }
863: static PetscErrorCode MatProductSetFromOptions_Dense_AB(Mat C)
864: {
865: C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
866: return 0;
867: }
869: static PetscErrorCode MatProductSetFromOptions_SchurComplement_Dense(Mat C)
870: {
871: Mat_Product *product = C->product;
874: MatProductSetFromOptions_Dense_AB(C);
875: return 0;
876: }
878: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
879: {
880: Mat_SchurComplement *Na;
882: PetscNew(&Na);
883: N->data = (void *)Na;
885: N->ops->destroy = MatDestroy_SchurComplement;
886: N->ops->getvecs = MatCreateVecs_SchurComplement;
887: N->ops->view = MatView_SchurComplement;
888: N->ops->mult = MatMult_SchurComplement;
889: N->ops->multtranspose = MatMultTranspose_SchurComplement;
890: N->ops->multadd = MatMultAdd_SchurComplement;
891: N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
892: N->assembled = PETSC_FALSE;
893: N->preallocated = PETSC_FALSE;
895: KSPCreate(PetscObjectComm((PetscObject)N), &Na->ksp);
896: PetscObjectChangeTypeName((PetscObject)N, MATSCHURCOMPLEMENT);
897: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", MatProductSetFromOptions_SchurComplement_Dense);
898: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", MatProductSetFromOptions_SchurComplement_Dense);
899: return 0;
900: }