Actual source code: schurm.c
1: #define PETSCMAT_DLL
3: #include private/matimpl.h
4: #include petscksp.h
6: typedef struct {
7: Mat A,B,C,D;
8: KSP ksp;
9: Vec work1,work2;
10: } Mat_SchurComplement;
12: /*
13: D - C inv(A) B
14: */
17: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
18: {
19: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
20: PetscErrorCode ierr;
23: if (!Na->work1) {MatGetVecs(Na->A,&Na->work1,PETSC_NULL);}
24: if (!Na->work2) {MatGetVecs(Na->A,&Na->work2,PETSC_NULL);}
25: MatMult(Na->B,x,Na->work1);
26: KSPSolve(Na->ksp,Na->work1,Na->work2);
27: MatMult(Na->C,Na->work2,y);
28: VecScale(y,-1.0);
29: if (Na->D) {
30: MatMultAdd(Na->D,x,y,y);
31: }
32: return(0);
33: }
37: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N)
38: {
39: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
40: PetscErrorCode ierr;
43: KSPSetFromOptions(Na->ksp);
44: return(0);
45: }
46:
49: PetscErrorCode MatDestroy_SchurComplement(Mat N)
50: {
51: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
52: PetscErrorCode ierr;
55: if (Na->A) {MatDestroy(Na->A);}
56: if (Na->B) {MatDestroy(Na->B);}
57: if (Na->C) {MatDestroy(Na->C);}
58: if (Na->D) {MatDestroy(Na->D);}
59: if (Na->work1) {VecDestroy(Na->work1);}
60: if (Na->work2) {VecDestroy(Na->work2);}
61: KSPDestroy(Na->ksp);
62: PetscFree(Na);
63: return(0);
64: }
68: /*@
69: MatCreateSchurComplement - Creates a new matrix object that behaves like the Schur complement of a matrix
71: Collective on Mat
73: Input Parameter:
74: . A,B,C,D - the four parts of the original matrix (D is optional)
76: Output Parameter:
77: . N - the matrix that the Schur complement D - C inv(A) B
79: Level: intermediate
81: Notes: The Schur complement is NOT actually formed! Rather this
82: object performs the matrix-vector product by using the the formula for
83: the Schur complement and a KSP solver to approximate the action of inv(A)
85: All four matrices must have the same MPI communicator
87: A and D must be square matrices
89: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdate(), MatCreateTranspose()
91: @*/
92: PetscErrorCode MatCreateSchurComplement(Mat A,Mat B,Mat C,Mat D,Mat *N)
93: {
94: PetscErrorCode ierr;
95: PetscInt m,n;
96: Mat_SchurComplement *Na;
104: if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local columns %D",A->rmap->n,A->cmap->n);
105: if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local rows B %D",A->rmap->n,B->rmap->n);
106: if (D) {
109: if (D->rmap->n != D->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of D %D do not equal local columns %D",D->rmap->n,D->cmap->n);
110: if (C->rmap->n != D->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of C %D do not equal local rows D %D",C->rmap->n,D->rmap->n);
111: }
113: MatGetLocalSize(B,PETSC_NULL,&n);
114: MatGetLocalSize(C,&m,PETSC_NULL);
115: MatCreate(((PetscObject)A)->comm,N);
116: MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);
117: PetscObjectChangeTypeName((PetscObject)*N,MATSCHURCOMPLEMENT);
118:
119: PetscNewLog(*N,Mat_SchurComplement,&Na);
120: (*N)->data = (void*) Na;
121: PetscObjectReference((PetscObject)A);
122: PetscObjectReference((PetscObject)B);
123: PetscObjectReference((PetscObject)C);
124: Na->A = A;
125: Na->B = B;
126: Na->C = C;
127: Na->D = D;
128: if (D) {
129: PetscObjectReference((PetscObject)D);
130: }
132: (*N)->ops->destroy = MatDestroy_SchurComplement;
133: (*N)->ops->mult = MatMult_SchurComplement;
134: (*N)->ops->setfromoptions = MatSetFromOptions_SchurComplement;
135: (*N)->assembled = PETSC_TRUE;
137: /* treats the new matrix as having block size of 1 which is most likely the case */
138: (*N)->rmap->bs = (*N)->cmap->bs = 1;
139: PetscMapSetUp((*N)->rmap);
140: PetscMapSetUp((*N)->cmap);
142: KSPCreate(((PetscObject)A)->comm,&Na->ksp);
143: KSPSetOptionsPrefix(Na->ksp,((PetscObject)A)->prefix);
144: KSPAppendOptionsPrefix(Na->ksp,"fieldsplit_0_");
145: KSPSetOperators(Na->ksp,A,A,SAME_NONZERO_PATTERN);
146: return(0);
147: }
151: /*@
152: MatSchurComplementGetKSP - Creates gets the KSP object that is used in the Schur complement matrix
154: Not Collective
156: Input Parameter:
157: . A - matrix created with MatCreateSchurComplement()
159: Output Parameter:
160: . ksp - the linear solver object
162: Options Database:
163: - -fieldsplit_0_XXX sets KSP and PC options for the A block solver inside the Schur complement
165: Level: intermediate
167: Notes:
168: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()
170: @*/
171: PetscErrorCode MatSchurComplementGetKSP(Mat A,KSP *ksp)
172: {
173: Mat_SchurComplement *Na = (Mat_SchurComplement*)A->data;
177: *ksp = Na->ksp;
178: return(0);
179: }
183: /*@
184: MatSchurComplementUpdate - Updates the Schur complement matrix object with new submatrices
186: Collective on Mat
188: Input Parameters:
189: + N - the matrix obtained with MatCreateSchurComplement()
190: . A,B,C,D - the four parts of the original matrix (D is optional)
191: - str - either SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
193:
194: Level: intermediate
196: Notes: All four matrices must have the same MPI communicator
198: A and D must be square matrices
200: All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement()
201: though they need not be the same matrices
203: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()
205: @*/
206: PetscErrorCode MatSchurComplementUpdate(Mat N,Mat A,Mat B,Mat C,Mat D,MatStructure str)
207: {
208: PetscErrorCode ierr;
209: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
217: if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local columns %D",A->rmap->n,A->cmap->n);
218: if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local rows B %D",A->rmap->n,B->rmap->n);
219: if (D) {
222: if (D->rmap->n != D->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of D %D do not equal local columns %D",D->rmap->n,D->cmap->n);
223: if (C->rmap->n != D->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of C %D do not equal local rows D %D",C->rmap->n,D->rmap->n);
224: }
226: MatDestroy(Na->A);
227: MatDestroy(Na->B);
228: MatDestroy(Na->C);
229: if (Na->D) {
230: MatDestroy(Na->D);
231: }
232: PetscObjectReference((PetscObject)A);
233: PetscObjectReference((PetscObject)B);
234: PetscObjectReference((PetscObject)C);
235: Na->A = A;
236: Na->B = B;
237: Na->C = C;
238: Na->D = D;
239: if (D) {
240: PetscObjectReference((PetscObject)D);
241: }
243: KSPSetOperators(Na->ksp,A,A,str);
244: return(0);
245: }
250: /*@C
251: MatSchurComplementGetSubmatrices - Get the individual submatrices in the Schur complement
253: Collective on Mat
255: Input Parameters:
256: + N - the matrix obtained with MatCreateSchurComplement()
257: - A,B,C,D - the four parts of the original matrix (D is optional)
259: Note:
260: D is optional, and thus can be PETSC_NULL
262: Level: intermediate
264: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdate()
265: @*/
266: PetscErrorCode MatSchurComplementGetSubmatrices(Mat N,Mat *A,Mat *B,Mat *C,Mat *D)
267: {
268: Mat_SchurComplement *Na = (Mat_SchurComplement *) N->data;
272: if (A) {*A = Na->A;}
273: if (B) {*B = Na->B;}
274: if (C) {*C = Na->C;}
275: if (D) {*D = Na->D;}
276: return(0);
277: }