Actual source code: schurm.c
petsc-3.9.0 2018-04-07
1: #include <petsc/private/matimpl.h>
2: #include <petscksp.h>
3: const char *const MatSchurComplementAinvTypes[] = {"DIAG","LUMP","BLOCKDIAG","MatSchurComplementAinvType","MAT_SCHUR_COMPLEMENT_AINV_",0};
5: typedef struct {
6: Mat A,Ap,B,C,D;
7: KSP ksp;
8: Vec work1,work2;
9: MatSchurComplementAinvType ainvtype;
10: } Mat_SchurComplement;
12: PetscErrorCode MatCreateVecs_SchurComplement(Mat N,Vec *right,Vec *left)
13: {
14: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
15: PetscErrorCode ierr;
18: if (Na->D) {
19: MatCreateVecs(Na->D,right,left);
20: return(0);
21: }
22: if (right) {
23: MatCreateVecs(Na->B,right,NULL);
24: }
25: if (left) {
26: MatCreateVecs(Na->C,NULL,left);
27: }
28: return(0);
29: }
31: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
32: {
33: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
34: PetscErrorCode ierr;
37: PetscViewerASCIIPrintf(viewer,"Schur complement A11 - A10 inv(A00) A01\n");
38: if (Na->D) {
39: PetscViewerASCIIPrintf(viewer,"A11\n");
40: PetscViewerASCIIPushTab(viewer);
41: MatView(Na->D,viewer);
42: PetscViewerASCIIPopTab(viewer);
43: } else {
44: PetscViewerASCIIPrintf(viewer,"A11 = 0\n");
45: }
46: PetscViewerASCIIPrintf(viewer,"A10\n");
47: PetscViewerASCIIPushTab(viewer);
48: MatView(Na->C,viewer);
49: PetscViewerASCIIPopTab(viewer);
50: PetscViewerASCIIPrintf(viewer,"KSP of A00\n");
51: PetscViewerASCIIPushTab(viewer);
52: KSPView(Na->ksp,viewer);
53: PetscViewerASCIIPopTab(viewer);
54: PetscViewerASCIIPrintf(viewer,"A01\n");
55: PetscViewerASCIIPushTab(viewer);
56: MatView(Na->B,viewer);
57: PetscViewerASCIIPopTab(viewer);
58: return(0);
59: }
61: /*
62: A11^T - A01^T ksptrans(A00,Ap00) A10^T
63: */
64: PetscErrorCode MatMultTranspose_SchurComplement(Mat N,Vec x,Vec y)
65: {
66: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
67: PetscErrorCode ierr;
70: if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
71: if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
72: MatMultTranspose(Na->C,x,Na->work1);
73: KSPSolveTranspose(Na->ksp,Na->work1,Na->work2);
74: MatMultTranspose(Na->B,Na->work2,y);
75: VecScale(y,-1.0);
76: if (Na->D) {
77: MatMultTransposeAdd(Na->D,x,y,y);
78: }
79: return(0);
80: }
82: /*
83: A11 - A10 ksp(A00,Ap00) A01
84: */
85: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
86: {
87: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
88: PetscErrorCode ierr;
91: if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
92: if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
93: MatMult(Na->B,x,Na->work1);
94: KSPSolve(Na->ksp,Na->work1,Na->work2);
95: MatMult(Na->C,Na->work2,y);
96: VecScale(y,-1.0);
97: if (Na->D) {
98: MatMultAdd(Na->D,x,y,y);
99: }
100: return(0);
101: }
103: /*
104: A11 - A10 ksp(A00,Ap00) A01
105: */
106: PetscErrorCode MatMultAdd_SchurComplement(Mat N,Vec x,Vec y,Vec z)
107: {
108: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
109: PetscErrorCode ierr;
112: if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
113: if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
114: MatMult(Na->B,x,Na->work1);
115: KSPSolve(Na->ksp,Na->work1,Na->work2);
116: if (y == z) {
117: VecScale(Na->work2,-1.0);
118: MatMultAdd(Na->C,Na->work2,z,z);
119: } else {
120: MatMult(Na->C,Na->work2,z);
121: VecAYPX(z,-1.0,y);
122: }
123: if (Na->D) {
124: MatMultAdd(Na->D,x,z,z);
125: }
126: return(0);
127: }
129: PetscErrorCode MatSetFromOptions_SchurComplement(PetscOptionItems *PetscOptionsObject,Mat N)
130: {
131: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
132: PetscErrorCode ierr;
135: PetscOptionsHead(PetscOptionsObject,"MatSchurComplementOptions");
136: Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
137: PetscOptionsEnum("-mat_schur_complement_ainv_type","Type of approximation for inv(A00) used when assembling Sp = A11 - A10 inv(A00) A01","MatSchurComplementSetAinvType",MatSchurComplementAinvTypes,(PetscEnum)Na->ainvtype,(PetscEnum*)&Na->ainvtype,NULL);
138: PetscOptionsTail();
139: KSPSetFromOptions(Na->ksp);
140: return(0);
141: }
143: PetscErrorCode MatDestroy_SchurComplement(Mat N)
144: {
145: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
146: PetscErrorCode ierr;
149: MatDestroy(&Na->A);
150: MatDestroy(&Na->Ap);
151: MatDestroy(&Na->B);
152: MatDestroy(&Na->C);
153: MatDestroy(&Na->D);
154: VecDestroy(&Na->work1);
155: VecDestroy(&Na->work2);
156: KSPDestroy(&Na->ksp);
157: PetscFree(N->data);
158: return(0);
159: }
161: /*@C
162: MatCreateSchurComplement - Creates a new matrix object that behaves like the Schur complement of a matrix
164: Collective on Mat
166: Input Parameters:
167: + A00,A01,A10,A11 - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
168: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
170: Output Parameter:
171: . S - the matrix that the Schur complement S = A11 - A10 ksp(A00,Ap00) A01
173: Level: intermediate
175: Notes: The Schur complement is NOT actually formed! Rather, this
176: object performs the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
177: for Schur complement S and a KSP solver to approximate the action of A^{-1}.
179: All four matrices must have the same MPI communicator.
181: A00 and A11 must be square matrices.
183: MatGetSchurComplement() takes as arguments the index sets for the submatrices and returns both the virtual Schur complement (what this returns) plus
184: a sparse approximation to the true Schur complement (useful for building a preconditioner for the Schur complement).
186: MatSchurComplementGetPmat() can be called on the output of this function to generate an explicit approximation to the Schur complement.
188: Developer Notes: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
189: remove redundancy and be clearer and simplier.
192: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdateSubMatrices(), MatCreateTranspose(), MatGetSchurComplement(),
193: MatSchurComplementGetPmat()
195: @*/
196: PetscErrorCode MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *S)
197: {
201: KSPInitializePackage();
202: MatCreate(((PetscObject)A00)->comm,S);
203: MatSetType(*S,MATSCHURCOMPLEMENT);
204: MatSchurComplementSetSubMatrices(*S,A00,Ap00,A01,A10,A11);
205: return(0);
206: }
208: /*@
209: MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement
211: Collective on Mat
213: Input Parameter:
214: + S - matrix obtained with MatCreateSchurComplement (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
215: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
216: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
218: Level: intermediate
220: Notes: The Schur complement is NOT actually formed! Rather, this
221: object performs the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
222: for Schur complement S and a KSP solver to approximate the action of A^{-1}.
224: All four matrices must have the same MPI communicator.
226: A00 and A11 must be square matrices.
228: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdateSubMatrices(), MatCreateTranspose(), MatCreateSchurComplement(), MatGetSchurComplement()
230: @*/
231: PetscErrorCode MatSchurComplementSetSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
232: {
233: PetscErrorCode ierr;
234: PetscInt m,n;
235: Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;
238: if (S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementUpdateSubMatrices() for already used matrix");
246: if (A00->rmap->n != A00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local columns %D",A00->rmap->n,A00->cmap->n);
247: if (A00->rmap->n != Ap00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local rows of Ap00 %D",A00->rmap->n,Ap00->rmap->n);
248: if (Ap00->rmap->n != Ap00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of Ap00 %D do not equal local columns %D",Ap00->rmap->n,Ap00->cmap->n);
249: if (A00->cmap->n != A01->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A00 %D do not equal local rows of A01 %D",A00->cmap->n,A01->rmap->n);
250: if (A10->cmap->n != A00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A10 %D do not equal local rows of A00 %D",A10->cmap->n,A00->rmap->n);
251: if (A11) {
254: if (A10->rmap->n != A11->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A10 %D do not equal local rows A11 %D",A10->rmap->n,A11->rmap->n);
255: }
257: MatGetLocalSize(A01,NULL,&n);
258: MatGetLocalSize(A10,&m,NULL);
259: MatSetSizes(S,m,n,PETSC_DECIDE,PETSC_DECIDE);
260: PetscObjectReference((PetscObject)A00);
261: PetscObjectReference((PetscObject)Ap00);
262: PetscObjectReference((PetscObject)A01);
263: PetscObjectReference((PetscObject)A10);
264: Na->A = A00;
265: Na->Ap = Ap00;
266: Na->B = A01;
267: Na->C = A10;
268: Na->D = A11;
269: if (A11) {
270: PetscObjectReference((PetscObject)A11);
271: }
272: S->assembled = PETSC_TRUE;
273: S->preallocated = PETSC_TRUE;
275: PetscLayoutSetUp((S)->rmap);
276: PetscLayoutSetUp((S)->cmap);
277: KSPSetOperators(Na->ksp,A00,Ap00);
278: return(0);
279: }
281: /*@
282: MatSchurComplementGetKSP - Gets the KSP object that is used to invert A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
284: Not Collective
286: Input Parameter:
287: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
289: Output Parameter:
290: . ksp - the linear solver object
292: Options Database:
293: . -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.
295: Level: intermediate
297: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
298: @*/
299: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
300: {
301: Mat_SchurComplement *Na;
306: Na = (Mat_SchurComplement*) S->data;
307: *ksp = Na->ksp;
308: return(0);
309: }
311: /*@
312: MatSchurComplementSetKSP - Sets the KSP object that is used to invert A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01
314: Not Collective
316: Input Parameters:
317: + S - matrix created with MatCreateSchurComplement()
318: - ksp - the linear solver object
320: Level: developer
322: Developer Notes:
323: This is used in PCFieldSplit to reuse the 0-split KSP to implement ksp(A00,Ap00) in S.
325: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
326: @*/
327: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
328: {
329: Mat_SchurComplement *Na;
330: PetscErrorCode ierr;
335: Na = (Mat_SchurComplement*) S->data;
336: PetscObjectReference((PetscObject)ksp);
337: KSPDestroy(&Na->ksp);
338: Na->ksp = ksp;
339: KSPSetOperators(Na->ksp, Na->A, Na->Ap);
340: return(0);
341: }
343: /*@
344: MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices
346: Collective on Mat
348: Input Parameters:
349: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
350: . A00,A01,A10,A11 - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
351: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
353: Level: intermediate
355: Notes: All four matrices must have the same MPI communicator
357: A00 and A11 must be square matrices
359: All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement() or MatSchurComplementSetSubMatrices()
360: though they need not be the same matrices.
362: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()
364: @*/
365: PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
366: {
367: PetscErrorCode ierr;
368: Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;
371: if (!S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementSetSubMatrices() for a new matrix");
378: if (A00->rmap->n != A00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local columns %D",A00->rmap->n,A00->cmap->n);
379: if (A00->rmap->n != Ap00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local rows of Ap00 %D",A00->rmap->n,Ap00->rmap->n);
380: if (Ap00->rmap->n != Ap00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of Ap00 %D do not equal local columns %D",Ap00->rmap->n,Ap00->cmap->n);
381: if (A00->cmap->n != A01->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A00 %D do not equal local rows of A01 %D",A00->cmap->n,A01->rmap->n);
382: if (A10->cmap->n != A00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A10 %D do not equal local rows of A00 %D",A10->cmap->n,A00->rmap->n);
383: if (A11) {
386: if (A10->rmap->n != A11->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A10 %D do not equal local rows A11 %D",A10->rmap->n,A11->rmap->n);
387: }
389: PetscObjectReference((PetscObject)A00);
390: PetscObjectReference((PetscObject)Ap00);
391: PetscObjectReference((PetscObject)A01);
392: PetscObjectReference((PetscObject)A10);
393: if (A11) {
394: PetscObjectReference((PetscObject)A11);
395: }
397: MatDestroy(&Na->A);
398: MatDestroy(&Na->Ap);
399: MatDestroy(&Na->B);
400: MatDestroy(&Na->C);
401: MatDestroy(&Na->D);
403: Na->A = A00;
404: Na->Ap = Ap00;
405: Na->B = A01;
406: Na->C = A10;
407: Na->D = A11;
409: KSPSetOperators(Na->ksp,A00,Ap00);
410: return(0);
411: }
414: /*@C
415: MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement
417: Collective on Mat
419: Input Parameter:
420: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
422: Output Paramters:
423: + A00,A01,A10,A11 - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
424: - Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.
426: Note: A11 is optional, and thus can be NULL. The submatrices are not increfed before they are returned and should not be modified or destroyed.
428: Level: intermediate
430: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdateSubMatrices()
431: @*/
432: PetscErrorCode MatSchurComplementGetSubMatrices(Mat S,Mat *A00,Mat *Ap00,Mat *A01,Mat *A10,Mat *A11)
433: {
434: Mat_SchurComplement *Na = (Mat_SchurComplement*) S->data;
435: PetscErrorCode ierr;
436: PetscBool flg;
440: PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&flg);
441: if (flg) {
442: if (A00) *A00 = Na->A;
443: if (Ap00) *Ap00 = Na->Ap;
444: if (A01) *A01 = Na->B;
445: if (A10) *A10 = Na->C;
446: if (A11) *A11 = Na->D;
447: } else {
448: if (A00) *A00 = 0;
449: if (Ap00) *Ap00 = 0;
450: if (A01) *A01 = 0;
451: if (A10) *A10 = 0;
452: if (A11) *A11 = 0;
453: }
454: return(0);
455: }
457: #include <petsc/private/kspimpl.h>
459: /*@
460: MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly
462: Collective on Mat
464: Input Parameter:
465: . M - the matrix obtained with MatCreateSchurComplement()
467: Output Parameter:
468: . S - the Schur complement matrix
470: Note: This can be expensive, so it is mainly for testing
472: Level: advanced
474: .seealso: MatCreateSchurComplement(), MatSchurComplementUpdate()
475: @*/
476: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat M, Mat *S)
477: {
478: Mat B, C, D;
479: KSP ksp;
480: PC pc;
481: PetscBool isLU, isILU;
482: PetscReal fill = 2.0;
486: MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
487: MatSchurComplementGetKSP(M, &ksp);
488: KSPGetPC(ksp, &pc);
489: PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
490: PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
491: if (isLU || isILU) {
492: Mat fact, Bd, AinvB, AinvBd;
493: PetscReal eps = 1.0e-10;
495: /* This can be sped up for banded LU */
496: KSPSetUp(ksp);
497: PCFactorGetMatrix(pc, &fact);
498: MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
499: MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
500: MatMatSolve(fact, Bd, AinvBd);
501: MatDestroy(&Bd);
502: MatChop(AinvBd, eps);
503: MatConvert(AinvBd, MATAIJ, MAT_INITIAL_MATRIX, &AinvB);
504: MatDestroy(&AinvBd);
505: MatMatMult(C, AinvB, MAT_INITIAL_MATRIX, fill, S);
506: MatDestroy(&AinvB);
507: } else {
508: Mat Ainvd, Ainv;
510: PCComputeExplicitOperator(pc, &Ainvd);
511: MatConvert(Ainvd, MATAIJ, MAT_INITIAL_MATRIX, &Ainv);
512: MatDestroy(&Ainvd);
513: #if 0
514: /* Symmetric version */
515: MatPtAP(Ainv, B, MAT_INITIAL_MATRIX, fill, S);
516: #else
517: /* Nonsymmetric version */
518: MatMatMatMult(C, Ainv, B, MAT_INITIAL_MATRIX, fill, S);
519: #endif
520: MatDestroy(&Ainv);
521: }
522: if (D) {
523: MatAXPY(*S, -1.0, D, DIFFERENT_NONZERO_PATTERN);
524: }
525: MatScale(*S,-1.0);
526: return(0);
527: }
529: /* Developer Notes: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
530: PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatSchurComplementAinvType ainvtype, MatReuse preuse,Mat *newpmat)
531: {
533: Mat A=0,Ap=0,B=0,C=0,D=0;
534: MatReuse reuse;
543: if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) return(0);
547: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
549: reuse = MAT_INITIAL_MATRIX;
550: if (mreuse == MAT_REUSE_MATRIX) {
551: MatSchurComplementGetSubMatrices(*newmat,&A,&Ap,&B,&C,&D);
552: if (!A || !Ap || !B || !C) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
553: if (A != Ap) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
554: MatDestroy(&Ap); /* get rid of extra reference */
555: reuse = MAT_REUSE_MATRIX;
556: }
557: MatCreateSubMatrix(mat,isrow0,iscol0,reuse,&A);
558: MatCreateSubMatrix(mat,isrow0,iscol1,reuse,&B);
559: MatCreateSubMatrix(mat,isrow1,iscol0,reuse,&C);
560: MatCreateSubMatrix(mat,isrow1,iscol1,reuse,&D);
561: switch (mreuse) {
562: case MAT_INITIAL_MATRIX:
563: MatCreateSchurComplement(A,A,B,C,D,newmat);
564: break;
565: case MAT_REUSE_MATRIX:
566: MatSchurComplementUpdateSubMatrices(*newmat,A,A,B,C,D);
567: break;
568: default:
569: if (mreuse != MAT_IGNORE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Unrecognized value of mreuse");
570: }
571: if (preuse != MAT_IGNORE_MATRIX) {
572: MatCreateSchurComplementPmat(A,B,C,D,ainvtype,preuse,newpmat);
573: }
574: MatDestroy(&A);
575: MatDestroy(&B);
576: MatDestroy(&C);
577: MatDestroy(&D);
578: return(0);
579: }
581: /*@
582: MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.
584: Collective on Mat
586: Input Parameters:
587: + A - matrix in which the complement is to be taken
588: . isrow0 - rows to eliminate
589: . iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
590: . isrow1 - rows in which the Schur complement is formed
591: . iscol1 - columns in which the Schur complement is formed
592: . mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in S
593: . ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
594: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_LUMP
595: - preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in Sp
597: Output Parameters:
598: + S - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
599: - Sp - approximate Schur complement from which a preconditioner can be built
601: Note:
602: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
603: application-specific information. The default for assembled matrices is to use the inverse of the diagonal of
604: the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
605: before forming inv(diag(A00)).
607: Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
608: 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
609: "MatGetSchurComplement_C" to their function. If their function needs to fall back to the default implementation, it
610: should call MatGetSchurComplement_Basic().
612: MatCreateSchurComplement() takes as arguments the four submatrices and returns the virtual Schur complement (what this returns in S).
614: MatSchurComplementGetPmat() takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).
616: In other words calling MatCreateSchurComplement() followed by MatSchurComplementGetPmat() produces the same output as this function but with slightly different
617: inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.
619: Developer Notes: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
620: remove redundancy and be clearer and simplier.
622: Level: advanced
624: Concepts: matrices^submatrices
626: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement(), MatSchurComplementAinvType
627: @*/
628: PetscErrorCode MatGetSchurComplement(Mat A,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *S,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Sp)
629: {
630: PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*) = NULL;
641: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
642: f = NULL;
643: 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. */
644: PetscObjectQueryFunction((PetscObject)*S,"MatGetSchurComplement_C",&f);
645: }
646: if (f) {
647: (*f)(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,preuse,Sp);
648: } else {
649: MatGetSchurComplement_Basic(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,ainvtype,preuse,Sp);
650: }
651: return(0);
652: }
654: /*@
655: MatSchurComplementSetAinvType - set the type of approximation used for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()
657: Not collective.
659: Input Parameters:
660: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
661: - ainvtype - type of approximation used to form A00inv from A00 when assembling Sp = A11 - A10 A00inv A01:
662: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG
664: Options database:
665: -mat_schur_complement_ainv_type diag | lump | blockdiag
667: Note:
668: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
669: application-specific information. The default for assembled matrices is to use the inverse of the diagonal of
670: the (0,0) block A00 in place of A00^{-1}. This rarely produces a scalable algorithm. Optionally, A00 can be lumped
671: before forming inv(diag(A00)).
673: Level: advanced
675: Concepts: matrices^submatrices
677: .seealso: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementGetAinvType()
678: @*/
679: PetscErrorCode MatSchurComplementSetAinvType(Mat S,MatSchurComplementAinvType ainvtype)
680: {
681: PetscErrorCode ierr;
682: const char* t;
683: PetscBool isschur;
684: Mat_SchurComplement *schur;
688: PetscObjectGetType((PetscObject)S,&t);
689: PetscStrcmp(t,MATSCHURCOMPLEMENT,&isschur);
690: if (!isschur) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected Mat of type MATSCHURCOMPLEMENT, got %s instead",t);
691: schur = (Mat_SchurComplement*)S->data;
692: if (ainvtype != MAT_SCHUR_COMPLEMENT_AINV_DIAG && ainvtype != MAT_SCHUR_COMPLEMENT_AINV_LUMP && ainvtype != MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %d",(int)ainvtype);
693: schur->ainvtype = ainvtype;
694: return(0);
695: }
697: /*@
698: MatSchurComplementGetAinvType - get the type of approximation for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()
700: Not collective.
702: Input Parameter:
703: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
705: Output Parameter:
706: . ainvtype - type of approximation used to form A00inv from A00 when assembling Sp = A11 - A10 A00inv A01:
707: MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG
709: Note:
710: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
711: application-specific information. The default for assembled matrices is to use the inverse of the diagonal of
712: the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
713: before forming inv(diag(A00)).
715: Level: advanced
717: Concepts: matrices^submatrices
719: .seealso: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementSetAinvType()
720: @*/
721: PetscErrorCode MatSchurComplementGetAinvType(Mat S,MatSchurComplementAinvType *ainvtype)
722: {
723: PetscErrorCode ierr;
724: const char* t;
725: PetscBool isschur;
726: Mat_SchurComplement *schur;
730: PetscObjectGetType((PetscObject)S,&t);
731: PetscStrcmp(t,MATSCHURCOMPLEMENT,&isschur);
732: if (!isschur) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected Mat of type MATSCHURCOMPLEMENT, got %s instead",t);
733: schur = (Mat_SchurComplement*)S->data;
734: if (ainvtype) *ainvtype = schur->ainvtype;
735: return(0);
736: }
738: /*@
739: MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(diag(A00)) A01
741: Collective on Mat
743: Input Parameters:
744: + A00,A01,A10,A11 - the four parts of the original matrix A = [A00 A01; A10 A11] (A01,A10, and A11 are optional, implying zero matrices)
745: . ainvtype - type of approximation for inv(A00) used when forming Sp = A11 - A10 inv(A00) A01
746: - 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
748: Output Parameter:
749: - Spmat - approximate Schur complement suitable for preconditioning S = A11 - A10 inv(diag(A00)) A01
751: Note:
752: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
753: application-specific information. The default for assembled matrices is to use the inverse of the diagonal of
754: the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
755: before forming inv(diag(A00)).
757: Level: advanced
759: Concepts: matrices^submatrices
761: .seealso: MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementAinvType
762: @*/
763: PetscErrorCode MatCreateSchurComplementPmat(Mat A00,Mat A01,Mat A10,Mat A11,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Spmat)
764: {
767: PetscInt N00;
770: /* Use an appropriate approximate inverse of A00 to form A11 - A10 inv(diag(A00)) A01; a NULL A01, A10 or A11 indicates a zero matrix. */
771: /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
772: if ((!A01 || !A10) & !A11) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot assemble Spmat: A01, A10 and A11 are all NULL.");
774: if (preuse == MAT_IGNORE_MATRIX) return(0);
776: /* A zero size A00 or empty A01 or A10 imply S = A11. */
777: MatGetSize(A00,&N00,NULL);
778: if (!A01 || !A10 || !N00) {
779: if (preuse == MAT_INITIAL_MATRIX) {
780: MatDuplicate(A11,MAT_COPY_VALUES,Spmat);
781: } else { /* MAT_REUSE_MATRIX */
782: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
783: MatCopy(A11,*Spmat,DIFFERENT_NONZERO_PATTERN);
784: }
786: } else {
787: Mat AdB;
788: Vec diag;
790: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
791: MatDuplicate(A01,MAT_COPY_VALUES,&AdB);
792: MatCreateVecs(A00,&diag,NULL);
793: if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
794: MatGetRowSum(A00,diag);
795: } else {
796: MatGetDiagonal(A00,diag);
797: }
798: VecReciprocal(diag);
799: MatDiagonalScale(AdB,diag,NULL);
800: VecDestroy(&diag);
801: } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
802: Mat A00_inv;
803: MatType type;
804: MPI_Comm comm;
806: PetscObjectGetComm((PetscObject)A00,&comm);
807: MatGetType(A00,&type);
808: MatCreate(comm,&A00_inv);
809: MatSetType(A00_inv,type);
810: MatInvertBlockDiagonalMat(A00,A00_inv);
811: MatMatMult(A00_inv,A01,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AdB);
812: MatDestroy(&A00_inv);
813: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %D", ainvtype);
814: /* Cannot really reuse Spmat in MatMatMult() because of MatAYPX() -->
815: MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult() */
816: MatDestroy(Spmat);
817: MatMatMult(A10,AdB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,Spmat);
818: if (!A11) {
819: MatScale(*Spmat,-1.0);
820: } else {
821: /* TODO: when can we pass SAME_NONZERO_PATTERN? */
822: MatAYPX(*Spmat,-1,A11,DIFFERENT_NONZERO_PATTERN);
823: }
824: MatDestroy(&AdB);
825: }
826: return(0);
827: }
829: PetscErrorCode MatSchurComplementGetPmat_Basic(Mat S,MatReuse preuse,Mat *Spmat)
830: {
831: Mat A,B,C,D;
832: Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;
833: PetscErrorCode ierr;
836: if (preuse == MAT_IGNORE_MATRIX) return(0);
838: MatSchurComplementGetSubMatrices(S,&A,NULL,&B,&C,&D);
839: if (!A) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Schur complement component matrices unset");
840: MatCreateSchurComplementPmat(A,B,C,D,schur->ainvtype,preuse,Spmat);
841: return(0);
842: }
844: /*@
845: MatSchurComplementGetPmat - Obtain a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(diag(A00)) A01
847: Collective on Mat
849: Input Parameters:
850: + S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
851: - 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
853: Output Parameter:
854: - Sp - approximate Schur complement suitable for preconditioning S = A11 - A10 inv(diag(A00)) A01
856: Note:
857: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
858: application-specific information. The default for assembled matrices is to use the inverse of the diagonal of
859: the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
860: before forming inv(diag(A00)).
862: Sometimes users would like to provide problem-specific data in the Schur complement, usually only
863: for special row and column index sets. In that case, the user should call PetscObjectComposeFunction() to set
864: "MatSchurComplementGetPmat_C" to their function. If their function needs to fall back to the default implementation,
865: it should call MatSchurComplementGetPmat_Basic().
867: Developer Notes: The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
868: remove redundancy and be clearer and simplier.
870: Level: advanced
872: Concepts: matrices^submatrices
874: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementSetAinvType()
875: @*/
876: PetscErrorCode MatSchurComplementGetPmat(Mat S,MatReuse preuse,Mat *Sp)
877: {
878: PetscErrorCode ierr,(*f)(Mat,MatReuse,Mat*);
884: if (S->factortype) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
886: PetscObjectQueryFunction((PetscObject)S,"MatSchurComplementGetPmat_C",&f);
887: if (f) {
888: (*f)(S,preuse,Sp);
889: } else {
890: MatSchurComplementGetPmat_Basic(S,preuse,Sp);
891: }
892: return(0);
893: }
895: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
896: {
897: PetscErrorCode ierr;
898: Mat_SchurComplement *Na;
901: PetscNewLog(N,&Na);
902: N->data = (void*) Na;
904: N->ops->destroy = MatDestroy_SchurComplement;
905: N->ops->getvecs = MatCreateVecs_SchurComplement;
906: N->ops->view = MatView_SchurComplement;
907: N->ops->mult = MatMult_SchurComplement;
908: N->ops->multtranspose = MatMultTranspose_SchurComplement;
909: N->ops->multadd = MatMultAdd_SchurComplement;
910: N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
911: N->assembled = PETSC_FALSE;
912: N->preallocated = PETSC_FALSE;
914: KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
915: PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
916: return(0);
917: }
919: static PetscBool KSPMatRegisterAllCalled;
921: /*@C
922: KSPMatRegisterAll - Registers all matrix implementations in the KSP package.
924: Not Collective
926: Level: advanced
928: .keywords: Mat, KSP, register, all
930: .seealso: MatRegisterAll(), KSPInitializePackage()
931: @*/
932: PetscErrorCode KSPMatRegisterAll(void)
933: {
937: if (KSPMatRegisterAllCalled) return(0);
938: KSPMatRegisterAllCalled = PETSC_TRUE;
939: MatRegister(MATSCHURCOMPLEMENT,MatCreate_SchurComplement);
940: return(0);
941: }