Actual source code: ex21.c
1: static const char help[] = "Tests MatGetSchurComplement\n";
3: #include <petscksp.h>
5: /*T
6: Concepts: Mat, Schur Complement
7: T*/
9: PetscErrorCode Create(MPI_Comm comm,Mat *inA,IS *is0,IS *is1)
10: {
11: Mat A;
12: PetscInt r,rend,M;
13: PetscMPIInt rank;
16: *inA = 0;
17: MatCreate(comm,&A);
18: MatSetSizes(A,4,4,PETSC_DETERMINE,PETSC_DETERMINE);
19: MatSetFromOptions(A);
20: MatSetUp(A);
21: MatGetOwnershipRange(A,&r,&rend);
22: MatGetSize(A,&M,NULL);
24: ISCreateStride(comm,2,r,1,is0);
25: ISCreateStride(comm,2,r+2,1,is1);
27: MPI_Comm_rank(comm,&rank);
29: {
30: PetscInt rows[4],cols0[5],cols1[5],cols2[3],cols3[3];
31: PetscScalar RR = 1000.*rank,vals0[5],vals1[4],vals2[3],vals3[3];
33: rows[0] = r;
34: rows[1] = r+1;
35: rows[2] = r+2;
36: rows[3] = r+3;
38: cols0[0] = r+0;
39: cols0[1] = r+1;
40: cols0[2] = r+3;
41: cols0[3] = (r+4)%M;
42: cols0[4] = (r+M-4)%M;
44: cols1[0] = r+1;
45: cols1[1] = r+2;
46: cols1[2] = (r+4+1)%M;
47: cols1[3] = (r+M-4+1)%M;
49: cols2[0] = r;
50: cols2[1] = r+2;
51: cols2[2] = (r+4+2)%M;
53: cols3[0] = r+1;
54: cols3[1] = r+3;
55: cols3[2] = (r+4+3)%M;
57: vals0[0] = RR+1.;
58: vals0[1] = RR+2.;
59: vals0[2] = RR+3.;
60: vals0[3] = RR+4.;
61: vals0[4] = RR+5.;
63: vals1[0] = RR+6.;
64: vals1[1] = RR+7.;
65: vals1[2] = RR+8.;
66: vals1[3] = RR+9.;
68: vals2[0] = RR+10.;
69: vals2[1] = RR+11.;
70: vals2[2] = RR+12.;
72: vals3[0] = RR+13.;
73: vals3[1] = RR+14.;
74: vals3[2] = RR+15.;
75: MatSetValues(A,1,&rows[0],5,cols0,vals0,INSERT_VALUES);
76: MatSetValues(A,1,&rows[1],4,cols1,vals1,INSERT_VALUES);
77: MatSetValues(A,1,&rows[2],3,cols2,vals2,INSERT_VALUES);
78: MatSetValues(A,1,&rows[3],3,cols3,vals3,INSERT_VALUES);
79: }
80: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);
82: *inA = A;
83: return 0;
84: }
86: PetscErrorCode Destroy(Mat *A,IS *is0,IS *is1)
87: {
89: MatDestroy(A);
90: ISDestroy(is0);
91: ISDestroy(is1);
92: return 0;
93: }
95: int main(int argc,char *argv[])
96: {
98: Mat A,S = NULL,Sexplicit = NULL;
99: MatSchurComplementAinvType ainv_type = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
100: IS is0,is1;
102: PetscInitialize(&argc,&argv,0,help);
103: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex21","KSP");
104: PetscOptionsEnum("-mat_schur_complement_ainv_type","Type of approximation for inv(A00) used when assembling Sp = A11 - A10 inv(A00) A01","MatSchurComplementAinvType",MatSchurComplementAinvTypes,(PetscEnum)ainv_type,(PetscEnum*)&ainv_type,NULL);
105: PetscOptionsEnd();
107: /* Test the Schur complement one way */
108: Create(PETSC_COMM_WORLD,&A,&is0,&is1);
109: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
110: ISView(is0,PETSC_VIEWER_STDOUT_WORLD);
111: ISView(is1,PETSC_VIEWER_STDOUT_WORLD);
112: MatGetSchurComplement(A,is0,is0,is1,is1,MAT_INITIAL_MATRIX,&S,ainv_type,MAT_IGNORE_MATRIX,NULL);
113: MatSetFromOptions(S);
114: MatComputeOperator(S,MATAIJ,&Sexplicit);
115: PetscPrintf(PETSC_COMM_WORLD,"\nExplicit Schur complement of (0,0) in (1,1)\n");
116: MatView(Sexplicit,PETSC_VIEWER_STDOUT_WORLD);
117: Destroy(&A,&is0,&is1);
118: MatDestroy(&S);
119: MatDestroy(&Sexplicit);
121: /* And the other */
122: Create(PETSC_COMM_WORLD,&A,&is0,&is1);
123: MatGetSchurComplement(A,is1,is1,is0,is0,MAT_INITIAL_MATRIX,&S,ainv_type,MAT_IGNORE_MATRIX,NULL);
124: MatSetFromOptions(S);
125: MatComputeOperator(S,MATAIJ,&Sexplicit);
126: PetscPrintf(PETSC_COMM_WORLD,"\nExplicit Schur complement of (1,1) in (0,0)\n");
127: MatView(Sexplicit,PETSC_VIEWER_STDOUT_WORLD);
128: Destroy(&A,&is0,&is1);
129: MatDestroy(&S);
130: MatDestroy(&Sexplicit);
132: /* This time just the preconditioning matrix. */
133: Create(PETSC_COMM_WORLD,&A,&is0,&is1);
134: MatGetSchurComplement(A,is0,is0,is1,is1,MAT_IGNORE_MATRIX,NULL,ainv_type,MAT_INITIAL_MATRIX,&S);
135: MatSetFromOptions(S);
136: PetscPrintf(PETSC_COMM_WORLD,"\nPreconditioning Schur complement of (0,0) in (1,1)\n");
137: MatView(S,PETSC_VIEWER_STDOUT_WORLD);
138: /* Modify and refresh */
139: MatShift(A,1.);
140: MatGetSchurComplement(A,is0,is0,is1,is1,MAT_IGNORE_MATRIX,NULL,ainv_type,MAT_REUSE_MATRIX,&S);
141: PetscPrintf(PETSC_COMM_WORLD,"\nAfter update\n");
142: MatView(S,PETSC_VIEWER_STDOUT_WORLD);
143: Destroy(&A,&is0,&is1);
144: MatDestroy(&S);
146: PetscFinalize();
147: return 0;
148: }
150: /*TEST
151: test:
152: suffix: diag_1
153: args: -mat_schur_complement_ainv_type diag
154: nsize: 1
155: test:
156: suffix: blockdiag_1
157: args: -mat_schur_complement_ainv_type blockdiag
158: nsize: 1
159: test:
160: suffix: diag_2
161: args: -mat_schur_complement_ainv_type diag
162: nsize: 2
163: test:
164: suffix: blockdiag_2
165: args: -mat_schur_complement_ainv_type blockdiag
166: nsize: 2
167: test:
168: # does not work with single because residual norm computed by GMRES recurrence formula becomes invalid
169: requires: !single
170: suffix: diag_3
171: args: -mat_schur_complement_ainv_type diag
172: nsize: 3
173: test:
174: # does not work with single because residual norm computed by GMRES recurrence formula becomes invalid
175: requires: !single
176: suffix: blockdiag_3
177: args: -mat_schur_complement_ainv_type blockdiag
178: nsize: 3
179: TEST*/