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 on A00
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 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 Notes:
161: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
162: remove redundancy and be clearer and simpler.
164: .seealso: `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatGetSchurComplement()`,
165: `MatSchurComplementGetPmat()`, `MatSchurComplementSetSubMatrices()`
167: @*/
168: PetscErrorCode MatCreateSchurComplement(Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11, Mat *S)
169: {
170: KSPInitializePackage();
171: MatCreate(PetscObjectComm((PetscObject)A00), S);
172: MatSetType(*S, MATSCHURCOMPLEMENT);
173: MatSchurComplementSetSubMatrices(*S, A00, Ap00, A01, A10, A11);
174: return 0;
175: }
177: /*@
178: MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement
180: Collective on S
182: Input Parameters:
183: + S - matrix obtained with MatSetType(S,MATSCHURCOMPLEMENT)
184: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
185: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
187: Level: intermediate
189: Notes:
190: The Schur complement is NOT explicitly formed! Rather, this
191: object performs the matrix-vector product of the Schur complement by using formula S = A11 - A10 ksp(A00,Ap00) A01
193: All four matrices must have the same MPI communicator.
195: A00 and A11 must be square matrices.
197: This is to be used in the context of code such as
198: $ MatSetType(S,MATSCHURCOMPLEMENT);
199: $ MatSchurComplementSetSubMatrices(S,...);
201: while MatSchurComplementUpdateSubMatrices() should only be called after MatCreateSchurComplement() or MatSchurComplementSetSubMatrices()
203: .seealso: `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatSchurComplementUpdateSubMatrices()`, `MatCreateTranspose()`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`
205: @*/
206: PetscErrorCode MatSchurComplementSetSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
207: {
208: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
209: PetscBool isschur;
211: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
212: if (!isschur) return 0;
226: if (A11) {
230: }
232: MatSetSizes(S, A10->rmap->n, A01->cmap->n, A10->rmap->N, A01->cmap->N);
233: PetscObjectReference((PetscObject)A00);
234: PetscObjectReference((PetscObject)Ap00);
235: PetscObjectReference((PetscObject)A01);
236: PetscObjectReference((PetscObject)A10);
237: Na->A = A00;
238: Na->Ap = Ap00;
239: Na->B = A01;
240: Na->C = A10;
241: Na->D = A11;
242: if (A11) PetscObjectReference((PetscObject)A11);
243: MatSetUp(S);
244: KSPSetOperators(Na->ksp, A00, Ap00);
245: S->assembled = PETSC_TRUE;
246: return 0;
247: }
249: /*@
250: MatSchurComplementGetKSP - Gets the KSP object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
252: Not Collective
254: Input Parameter:
255: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
257: Output Parameter:
258: . ksp - the linear solver object
260: Options Database:
261: . -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.
263: Level: intermediate
265: .seealso: `MatSchurComplementSetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`
266: @*/
267: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
268: {
269: Mat_SchurComplement *Na;
270: PetscBool isschur;
273: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
276: Na = (Mat_SchurComplement *)S->data;
277: *ksp = Na->ksp;
278: return 0;
279: }
281: /*@
282: MatSchurComplementSetKSP - Sets the KSP object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
284: Not Collective
286: Input Parameters:
287: + S - matrix created with MatCreateSchurComplement()
288: - ksp - the linear solver object
290: Level: developer
292: Developer Notes:
293: This is used in PCFieldSplit to reuse the 0-split KSP to implement ksp(A00,Ap00) in S.
294: The KSP operators are overwritten with A00 and Ap00 currently set in S.
296: .seealso: `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MATSCHURCOMPLEMENT`
297: @*/
298: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
299: {
300: Mat_SchurComplement *Na;
301: PetscBool isschur;
304: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
305: if (!isschur) return 0;
307: Na = (Mat_SchurComplement *)S->data;
308: PetscObjectReference((PetscObject)ksp);
309: KSPDestroy(&Na->ksp);
310: Na->ksp = ksp;
311: KSPSetOperators(Na->ksp, Na->A, Na->Ap);
312: return 0;
313: }
315: /*@
316: MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices
318: Collective on S
320: Input Parameters:
321: + S - matrix obtained with MatCreateSchurComplement() (or MatSchurSetSubMatrices()) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
322: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
323: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
325: Level: intermediate
327: Notes:
328: All four matrices must have the same MPI communicator
330: A00 and A11 must be square matrices
332: All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement() or MatSchurComplementSetSubMatrices()
333: though they need not be the same matrices.
335: This can only be called after MatCreateSchurComplement() or MatSchurComplementSetSubMatrices(), it cannot be called immediately after MatSetType(S,MATSCHURCOMPLEMENT);
337: Developer Notes:
338: This code is almost identical to MatSchurComplementSetSubMatrices(). The API should be refactored.
340: .seealso: `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`
342: @*/
343: PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S, Mat A00, Mat Ap00, Mat A01, Mat A10, Mat A11)
344: {
345: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
346: PetscBool isschur;
349: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
350: if (!isschur) return 0;
364: if (A11) {
368: }
370: PetscObjectReference((PetscObject)A00);
371: PetscObjectReference((PetscObject)Ap00);
372: PetscObjectReference((PetscObject)A01);
373: PetscObjectReference((PetscObject)A10);
374: if (A11) PetscObjectReference((PetscObject)A11);
376: MatDestroy(&Na->A);
377: MatDestroy(&Na->Ap);
378: MatDestroy(&Na->B);
379: MatDestroy(&Na->C);
380: MatDestroy(&Na->D);
382: Na->A = A00;
383: Na->Ap = Ap00;
384: Na->B = A01;
385: Na->C = A10;
386: Na->D = A11;
388: KSPSetOperators(Na->ksp, A00, Ap00);
389: return 0;
390: }
392: /*@C
393: MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement
395: Collective on S
397: Input Parameter:
398: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
400: Output Parameters:
401: + A00 - the upper-left block of the original matrix A = [A00 A01; A10 A11]
402: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
403: . A01 - the upper-right block of the original matrix A = [A00 A01; A10 A11]
404: . A10 - the lower-left block of the original matrix A = [A00 A01; A10 A11]
405: - A11 - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]
407: Note: 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.
409: Level: intermediate
411: .seealso: `MatCreateNormal()`, `MatMult()`, `MatCreate()`, `MatSchurComplementGetKSP()`, `MatCreateSchurComplement()`, `MatSchurComplementUpdateSubMatrices()`
412: @*/
413: PetscErrorCode MatSchurComplementGetSubMatrices(Mat S, Mat *A00, Mat *Ap00, Mat *A01, Mat *A10, Mat *A11)
414: {
415: Mat_SchurComplement *Na = (Mat_SchurComplement *)S->data;
416: PetscBool flg;
419: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &flg);
421: if (A00) *A00 = Na->A;
422: if (Ap00) *Ap00 = Na->Ap;
423: if (A01) *A01 = Na->B;
424: if (A10) *A10 = Na->C;
425: if (A11) *A11 = Na->D;
426: return 0;
427: }
429: #include <petsc/private/kspimpl.h>
431: /*@
432: MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly
434: Collective on M
436: Input Parameter:
437: . M - the matrix obtained with MatCreateSchurComplement()
439: Output Parameter:
440: . S - the Schur complement matrix
442: Notes:
443: This can be expensive, so it is mainly for testing
445: Use MatSchurComplementGetPmat() to get a sparse approximation for the Schur complement suitable for use in building a preconditioner
447: Level: advanced
449: .seealso: `MatCreateSchurComplement()`, `MatSchurComplementUpdate()`, `MatSchurComplementGetPmat()`
450: @*/
451: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat A, Mat *S)
452: {
453: Mat B, C, D, E = NULL, Bd, AinvBd;
454: KSP ksp;
455: PetscInt n, N, m, M;
456: PetscBool flg = PETSC_FALSE, set, symm;
458: MatSchurComplementGetSubMatrices(A, NULL, NULL, &B, &C, &D);
459: MatSchurComplementGetKSP(A, &ksp);
460: KSPSetUp(ksp);
461: MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
462: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
463: KSPMatSolve(ksp, Bd, AinvBd);
464: MatDestroy(&Bd);
465: MatChop(AinvBd, PETSC_SMALL);
466: if (D) {
467: MatGetLocalSize(D, &m, &n);
468: MatGetSize(D, &M, &N);
469: MatCreateDense(PetscObjectComm((PetscObject)A), m, n, M, N, NULL, S);
470: }
471: MatMatMult(C, AinvBd, D ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, S);
472: MatDestroy(&AinvBd);
473: if (D) {
474: PetscObjectTypeCompareAny((PetscObject)D, &flg, MATSEQSBAIJ, MATMPISBAIJ, "");
475: if (flg) {
476: MatIsSymmetricKnown(A, &set, &symm);
477: if (!set || !symm) MatConvert(D, MATBAIJ, MAT_INITIAL_MATRIX, &E); /* convert the (1,1) block to nonsymmetric storage for MatAXPY() */
478: }
479: 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 */
480: }
481: 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 */
482: MatScale(*S, -1.0);
483: MatDestroy(&E);
484: return 0;
485: }
487: /* Developer Notes:
488: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
489: PetscErrorCode MatGetSchurComplement_Basic(Mat mat, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
490: {
491: Mat A = NULL, Ap = NULL, B = NULL, C = NULL, D = NULL;
492: MatReuse reuse;
503: if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) return 0;
509: reuse = MAT_INITIAL_MATRIX;
510: if (mreuse == MAT_REUSE_MATRIX) {
511: MatSchurComplementGetSubMatrices(*S, &A, &Ap, &B, &C, &D);
514: MatDestroy(&Ap); /* get rid of extra reference */
515: reuse = MAT_REUSE_MATRIX;
516: }
517: MatCreateSubMatrix(mat, isrow0, iscol0, reuse, &A);
518: MatCreateSubMatrix(mat, isrow0, iscol1, reuse, &B);
519: MatCreateSubMatrix(mat, isrow1, iscol0, reuse, &C);
520: MatCreateSubMatrix(mat, isrow1, iscol1, reuse, &D);
521: switch (mreuse) {
522: case MAT_INITIAL_MATRIX:
523: MatCreateSchurComplement(A, A, B, C, D, S);
524: break;
525: case MAT_REUSE_MATRIX:
526: MatSchurComplementUpdateSubMatrices(*S, A, A, B, C, D);
527: break;
528: default:
530: }
531: if (preuse != MAT_IGNORE_MATRIX) MatCreateSchurComplementPmat(A, B, C, D, ainvtype, preuse, Sp);
532: MatDestroy(&A);
533: MatDestroy(&B);
534: MatDestroy(&C);
535: MatDestroy(&D);
536: return 0;
537: }
539: /*@
540: MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.
542: Collective on A
544: Input Parameters:
545: + A - matrix in which the complement is to be taken
546: . isrow0 - rows to eliminate
547: . iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
548: . isrow1 - rows in which the Schur complement is formed
549: . iscol1 - columns in which the Schur complement is formed
550: . mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in S
551: . ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
552: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_FULL
553: - preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in Sp
555: Output Parameters:
556: + S - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
557: - Sp - approximate Schur complement from which a preconditioner can be built A11 - A10 inv(DIAGFORM(A00)) A01
559: Note:
560: Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
561: application-specific information.
563: Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
564: 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
565: "MatGetSchurComplement_C" to their function. If their function needs to fall back to the default implementation, it
566: should call MatGetSchurComplement_Basic().
568: MatCreateSchurComplement() takes as arguments the four submatrices and returns the virtual Schur complement (what this function returns in S).
570: MatSchurComplementGetPmat() takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).
572: In other words calling MatCreateSchurComplement() followed by MatSchurComplementGetPmat() produces the same output as this function but with slightly different
573: inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.
575: Developer Notes:
576: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
577: remove redundancy and be clearer and simpler.
579: Level: advanced
581: .seealso: `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatCreateSchurComplement()`, `MatSchurComplementAinvType`
582: @*/
583: PetscErrorCode MatGetSchurComplement(Mat A, IS isrow0, IS iscol0, IS isrow1, IS iscol1, MatReuse mreuse, Mat *S, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
584: {
585: PetscErrorCode (*f)(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatReuse, Mat *) = NULL;
599: 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. */
600: PetscObjectQueryFunction((PetscObject)*S, "MatGetSchurComplement_C", &f);
601: }
602: if (f) (*f)(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, preuse, Sp);
603: else MatGetSchurComplement_Basic(A, isrow0, iscol0, isrow1, iscol1, mreuse, S, ainvtype, preuse, Sp);
604: return 0;
605: }
607: /*@
608: MatSchurComplementSetAinvType - set the type of approximation used for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()
610: Not collective.
612: Input Parameters:
613: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
614: - ainvtype - type of approximation to be used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
615: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_FULL
617: Options database:
618: -mat_schur_complement_ainv_type diag | lump | blockdiag | full
620: Level: advanced
622: .seealso: `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementGetAinvType()`
623: @*/
624: PetscErrorCode MatSchurComplementSetAinvType(Mat S, MatSchurComplementAinvType ainvtype)
625: {
626: PetscBool isschur;
627: Mat_SchurComplement *schur;
630: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
631: if (!isschur) return 0;
633: schur = (Mat_SchurComplement *)S->data;
635: schur->ainvtype = ainvtype;
636: return 0;
637: }
639: /*@
640: MatSchurComplementGetAinvType - get the type of approximation for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()
642: Not collective.
644: Input Parameter:
645: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
647: Output Parameter:
648: . ainvtype - type of approximation used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
649: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_FULL
651: Level: advanced
653: .seealso: `MatSchurComplementAinvType`, `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementSetAinvType()`
654: @*/
655: PetscErrorCode MatSchurComplementGetAinvType(Mat S, MatSchurComplementAinvType *ainvtype)
656: {
657: PetscBool isschur;
658: Mat_SchurComplement *schur;
661: PetscObjectTypeCompare((PetscObject)S, MATSCHURCOMPLEMENT, &isschur);
663: schur = (Mat_SchurComplement *)S->data;
664: if (ainvtype) *ainvtype = schur->ainvtype;
665: return 0;
666: }
668: /*@
669: MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by explicitly assembling the sparse matrix
670: Sp = A11 - A10 inv(DIAGFORM(A00)) A01
672: Collective on A00
674: Input Parameters:
675: + A00 - the upper-left part of the original matrix A = [A00 A01; A10 A11]
676: . A01 - (optional) the upper-right part of the original matrix A = [A00 A01; A10 A11]
677: . A10 - (optional) the lower-left part of the original matrix A = [A00 A01; A10 A11]
678: . A11 - (optional) the lower-right part of the original matrix A = [A00 A01; A10 A11]
679: . ainvtype - type of approximation for DIAGFORM(A00) used when forming Sp = A11 - A10 inv(DIAGFORM(A00)) A01. See MatSchurComplementAinvType.
680: - 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
682: Output Parameter:
683: - Sp - approximate Schur complement suitable for preconditioning the true Schur complement S = A11 - A10 inv(A00) A01
685: Level: advanced
687: .seealso: `MatCreateSchurComplement()`, `MatGetSchurComplement()`, `MatSchurComplementGetPmat()`, `MatSchurComplementAinvType`
688: @*/
689: PetscErrorCode MatCreateSchurComplementPmat(Mat A00, Mat A01, Mat A10, Mat A11, MatSchurComplementAinvType ainvtype, MatReuse preuse, Mat *Sp)
690: {
691: PetscInt N00;
693: /* 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. */
694: /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
696: if (preuse == MAT_IGNORE_MATRIX) return 0;
698: /* A zero size A00 or empty A01 or A10 imply S = A11. */
699: MatGetSize(A00, &N00, NULL);
700: if (!A01 || !A10 || !N00) {
701: if (preuse == MAT_INITIAL_MATRIX) {
702: MatDuplicate(A11, MAT_COPY_VALUES, Sp);
703: } else { /* MAT_REUSE_MATRIX */
704: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
705: MatCopy(A11, *Sp, DIFFERENT_NONZERO_PATTERN);
706: }
707: } else {
708: Mat AdB;
709: Vec diag;
711: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
712: MatDuplicate(A01, MAT_COPY_VALUES, &AdB);
713: MatCreateVecs(A00, &diag, NULL);
714: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
715: MatGetRowSum(A00, diag);
716: } else {
717: MatGetDiagonal(A00, diag);
718: }
719: VecReciprocal(diag);
720: MatDiagonalScale(AdB, diag, NULL);
721: VecDestroy(&diag);
722: } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
723: Mat A00_inv;
724: MatType type;
725: MPI_Comm comm;
727: PetscObjectGetComm((PetscObject)A00, &comm);
728: MatGetType(A00, &type);
729: MatCreate(comm, &A00_inv);
730: MatSetType(A00_inv, type);
731: MatInvertBlockDiagonalMat(A00, A00_inv);
732: MatMatMult(A00_inv, A01, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &AdB);
733: MatDestroy(&A00_inv);
734: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown MatSchurComplementAinvType: %d", ainvtype);
735: /* Cannot really reuse Sp in MatMatMult() because of MatAYPX() -->
736: MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult() */
737: if (preuse == MAT_REUSE_MATRIX) MatDestroy(Sp);
738: MatMatMult(A10, AdB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, Sp);
739: if (!A11) {
740: MatScale(*Sp, -1.0);
741: } else {
742: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
743: MatAYPX(*Sp, -1, A11, DIFFERENT_NONZERO_PATTERN);
744: }
745: MatDestroy(&AdB);
746: }
747: return 0;
748: }
750: PetscErrorCode MatSchurComplementGetPmat_Basic(Mat S, MatReuse preuse, Mat *Sp)
751: {
752: Mat A, B, C, D;
753: Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;
755: if (preuse == MAT_IGNORE_MATRIX) return 0;
756: MatSchurComplementGetSubMatrices(S, &A, NULL, &B, &C, &D);
758: if (schur->ainvtype != MAT_SCHUR_COMPLEMENT_AINV_FULL) MatCreateSchurComplementPmat(A, B, C, D, schur->ainvtype, preuse, Sp);
759: else {
760: if (preuse == MAT_REUSE_MATRIX) MatDestroy(Sp);
761: MatSchurComplementComputeExplicitOperator(S, Sp);
762: }
763: return 0;
764: }
766: /*@
767: MatSchurComplementGetPmat - Obtain a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01
769: Collective on S
771: Input Parameters:
772: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of A11 - A10 ksp(A00,Ap00) A01
773: - 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
775: Output Parameter:
776: - Sp - approximate Schur complement suitable for preconditioning the exact Schur complement S = A11 - A10 inv(A00) A01
778: Note:
779: The approximation of Sp depends on the the argument passed to to MatSchurComplementSetAinvType()
780: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_FULL
781: -mat_schur_complement_ainv_type <diag,lump,blockdiag,full>
783: Sometimes users would like to provide problem-specific data in the Schur complement, usually only
784: for special row and column index sets. In that case, the user should call PetscObjectComposeFunction() to set
785: "MatSchurComplementGetPmat_C" to their function. If their function needs to fall back to the default implementation,
786: it should call MatSchurComplementGetPmat_Basic().
788: Developer Notes:
789: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
790: remove redundancy and be clearer and simpler.
792: This routine should be called MatSchurComplementCreatePmat()
794: Level: advanced
796: .seealso: `MatCreateSubMatrix()`, `PCFIELDSPLIT`, `MatGetSchurComplement()`, `MatCreateSchurComplement()`, `MatSchurComplementSetAinvType()`
797: @*/
798: PetscErrorCode MatSchurComplementGetPmat(Mat S, MatReuse preuse, Mat *Sp)
799: {
800: PetscErrorCode (*f)(Mat, MatReuse, Mat *);
805: if (preuse != MAT_IGNORE_MATRIX) {
807: if (preuse == MAT_INITIAL_MATRIX) *Sp = NULL;
809: }
812: PetscObjectQueryFunction((PetscObject)S, "MatSchurComplementGetPmat_C", &f);
813: if (f) (*f)(S, preuse, Sp);
814: else MatSchurComplementGetPmat_Basic(S, preuse, Sp);
815: return 0;
816: }
818: static PetscErrorCode MatProductNumeric_SchurComplement_Dense(Mat C)
819: {
820: Mat_Product *product = C->product;
821: Mat_SchurComplement *Na = (Mat_SchurComplement *)product->A->data;
822: Mat work1, work2;
823: PetscScalar *v;
824: PetscInt lda;
826: MatMatMult(Na->B, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1);
827: MatDuplicate(work1, MAT_DO_NOT_COPY_VALUES, &work2);
828: KSPMatSolve(Na->ksp, work1, work2);
829: MatDestroy(&work1);
830: MatDenseGetArrayWrite(C, &v);
831: MatDenseGetLDA(C, &lda);
832: MatCreateDense(PetscObjectComm((PetscObject)C), C->rmap->n, C->cmap->n, C->rmap->N, C->cmap->N, v, &work1);
833: MatDenseSetLDA(work1, lda);
834: MatMatMult(Na->C, work2, MAT_REUSE_MATRIX, PETSC_DEFAULT, &work1);
835: MatDenseRestoreArrayWrite(C, &v);
836: MatDestroy(&work2);
837: MatDestroy(&work1);
838: if (Na->D) {
839: MatMatMult(Na->D, product->B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &work1);
840: MatAYPX(C, -1.0, work1, SAME_NONZERO_PATTERN);
841: MatDestroy(&work1);
842: } else MatScale(C, -1.0);
843: return 0;
844: }
846: static PetscErrorCode MatProductSymbolic_SchurComplement_Dense(Mat C)
847: {
848: Mat_Product *product = C->product;
849: Mat A = product->A, B = product->B;
850: PetscInt m = A->rmap->n, n = B->cmap->n, M = A->rmap->N, N = B->cmap->N;
851: PetscBool flg;
853: MatSetSizes(C, m, n, M, N);
854: PetscObjectBaseTypeCompareAny((PetscObject)C, &flg, MATSEQDENSE, MATMPIDENSE, "");
855: if (!flg) {
856: MatSetType(C, ((PetscObject)B)->type_name);
857: C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
858: }
859: MatSetUp(C);
860: C->ops->productnumeric = MatProductNumeric_SchurComplement_Dense;
861: return 0;
862: }
864: static PetscErrorCode MatProductSetFromOptions_Dense_AB(Mat C)
865: {
866: C->ops->productsymbolic = MatProductSymbolic_SchurComplement_Dense;
867: return 0;
868: }
870: static PetscErrorCode MatProductSetFromOptions_SchurComplement_Dense(Mat C)
871: {
872: Mat_Product *product = C->product;
875: MatProductSetFromOptions_Dense_AB(C);
876: return 0;
877: }
879: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
880: {
881: Mat_SchurComplement *Na;
883: PetscNew(&Na);
884: N->data = (void *)Na;
886: N->ops->destroy = MatDestroy_SchurComplement;
887: N->ops->getvecs = MatCreateVecs_SchurComplement;
888: N->ops->view = MatView_SchurComplement;
889: N->ops->mult = MatMult_SchurComplement;
890: N->ops->multtranspose = MatMultTranspose_SchurComplement;
891: N->ops->multadd = MatMultAdd_SchurComplement;
892: N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
893: N->assembled = PETSC_FALSE;
894: N->preallocated = PETSC_FALSE;
896: KSPCreate(PetscObjectComm((PetscObject)N), &Na->ksp);
897: PetscObjectChangeTypeName((PetscObject)N, MATSCHURCOMPLEMENT);
898: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_seqdense_C", MatProductSetFromOptions_SchurComplement_Dense);
899: PetscObjectComposeFunction((PetscObject)N, "MatProductSetFromOptions_schurcomplement_mpidense_C", MatProductSetFromOptions_SchurComplement_Dense);
900: return 0;
901: }