xref: /petsc/src/snes/tutorials/ex70.c (revision 70a7d78aacfbd24b2e31399a3d8e056944bb7de3)
1 static char help[] = "Poiseuille flow problem. Viscous, laminar flow in a 2D channel with parabolic velocity\n\
2                       profile and linear pressure drop, exact solution of the 2D Stokes\n";
3 
4 /*---------------------------------------------------------------------------- */
5 /* M A R I T I M E  R E S E A R C H  I N S T I T U T E  N E T H E R L A N D S  */
6 /*---------------------------------------------------------------------------- */
7 /* author : Christiaan M. Klaij                                                */
8 /*---------------------------------------------------------------------------- */
9 /*                                                                             */
10 /* Poiseuille flow problem.                                                    */
11 /*                                                                             */
12 /* Viscous, laminar flow in a 2D channel with parabolic velocity               */
13 /* profile and linear pressure drop, exact solution of the 2D Stokes           */
14 /* equations.                                                                  */
15 /*                                                                             */
16 /* Discretized with the cell-centered finite-volume method on a                */
17 /* Cartesian grid with co-located variables. Variables ordered as              */
18 /* [u1...uN v1...vN p1...pN]^T. Matrix [A00 A01; A10, A11] solved with         */
19 /* PCFIELDSPLIT. Lower factorization is used to mimick the Semi-Implicit       */
20 /* Method for Pressure Linked Equations (SIMPLE) used as preconditioner        */
21 /* instead of solver.                                                          */
22 /*                                                                             */
23 /* Disclaimer: does not contain the pressure-weighed interpolation             */
24 /* method needed to suppress spurious pressure modes in real-life              */
25 /* problems.                                                                   */
26 /*                                                                             */
27 /* Usage:                                                                      */
28 /*                                                                             */
29 /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none */
30 /*                                                                             */
31 /*   Runs with PCFIELDSPLIT on 32x48 grid, no PC for the Schur                 */
32 /*   complement because A11 is zero. FGMRES is needed because                  */
33 /*   PCFIELDSPLIT is a variable preconditioner.                                */
34 /*                                                                             */
35 /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc */
36 /*                                                                             */
37 /*   Same as above but with user defined PC for the true Schur                 */
38 /*   complement. PC based on the SIMPLE-type approximation (inverse of         */
39 /*   A00 approximated by inverse of its diagonal).                             */
40 /*                                                                             */
41 /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_ksp */
42 /*                                                                             */
43 /*   Replace the true Schur complement with a user defined Schur               */
44 /*   complement based on the SIMPLE-type approximation. Same matrix is         */
45 /*   used as PC.                                                               */
46 /*                                                                             */
47 /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi */
48 /*                                                                             */
49 /*   Out-of-the-box SIMPLE-type preconditioning. The major advantage           */
50 /*   is that the user neither needs to provide the approximation of            */
51 /*   the Schur complement, nor the corresponding preconditioner.               */
52 /*                                                                             */
53 /*---------------------------------------------------------------------------- */
54 
55 #include <petscksp.h>
56 
57 typedef struct {
58   PetscBool userPC, userKSP, matsymmetric; /* user defined preconditioner and matrix for the Schur complement */
59   PetscInt  nx, ny;  /* nb of cells in x- and y-direction */
60   PetscReal hx, hy;  /* mesh size in x- and y-direction */
61   Mat       A;       /* block matrix */
62   Mat       subA[4]; /* the four blocks */
63   Mat       myS;     /* the approximation of the Schur complement */
64   Vec       x, b, y; /* solution, rhs and temporary vector */
65   IS        isg[2];  /* index sets of split "0" and "1" */
66 } Stokes;
67 
68 PetscErrorCode StokesSetupMatBlock00(Stokes*);  /* setup the block Q */
69 PetscErrorCode StokesSetupMatBlock01(Stokes*);  /* setup the block G */
70 PetscErrorCode StokesSetupMatBlock10(Stokes*);  /* setup the block D (equal to the transpose of G) */
71 PetscErrorCode StokesSetupMatBlock11(Stokes*);  /* setup the block C (equal to zero) */
72 
73 PetscErrorCode StokesGetPosition(Stokes*, PetscInt, PetscInt*, PetscInt*); /* row number j*nx+i corresponds to position (i,j) in grid */
74 
75 PetscErrorCode StokesStencilLaplacian(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Laplacian operator */
76 PetscErrorCode StokesStencilGradientX(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (x-component) */
77 PetscErrorCode StokesStencilGradientY(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (y-component) */
78 
79 PetscErrorCode StokesRhs(Stokes*);                                               /* rhs vector */
80 PetscErrorCode StokesRhsMomX(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (x-component) */
81 PetscErrorCode StokesRhsMomY(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (y-component) */
82 PetscErrorCode StokesRhsMass(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of pressure */
83 
84 PetscErrorCode StokesSetupApproxSchur(Stokes*);  /* approximation of the Schur complement */
85 
86 PetscErrorCode StokesExactSolution(Stokes*); /* exact solution vector */
87 PetscErrorCode StokesWriteSolution(Stokes*); /* write solution to file */
88 
89 /* exact solution for the velocity (x-component, y-component is zero) */
90 PetscScalar StokesExactVelocityX(const PetscScalar y)
91 {
92   return 4.0*y*(1.0-y);
93 }
94 
95 /* exact solution for the pressure */
96 PetscScalar StokesExactPressure(const PetscScalar x)
97 {
98   return 8.0*(2.0-x);
99 }
100 
101 PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp)
102 {
103   KSP            *subksp;
104   PC             pc;
105   PetscInt       n = 1;
106   PetscErrorCode ierr;
107 
108   PetscFunctionBeginUser;
109   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
110   ierr = PCFieldSplitSetIS(pc, "0", s->isg[0]);CHKERRQ(ierr);
111   ierr = PCFieldSplitSetIS(pc, "1", s->isg[1]);CHKERRQ(ierr);
112   if (s->userPC) {
113     ierr = PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS);CHKERRQ(ierr);
114   }
115   if (s->userKSP) {
116     ierr = PCSetUp(pc);CHKERRQ(ierr);
117     ierr = PCFieldSplitGetSubKSP(pc, &n, &subksp);CHKERRQ(ierr);
118     ierr = KSPSetOperators(subksp[1], s->myS, s->myS);CHKERRQ(ierr);
119     ierr = PetscFree(subksp);CHKERRQ(ierr);
120   }
121   PetscFunctionReturn(0);
122 }
123 
124 PetscErrorCode StokesWriteSolution(Stokes *s)
125 {
126   PetscMPIInt       size;
127   PetscInt          n,i,j;
128   const PetscScalar *array;
129   PetscErrorCode    ierr;
130 
131   PetscFunctionBeginUser;
132   /* write data (*warning* only works sequential) */
133   ierr = MPI_Comm_size(MPI_COMM_WORLD,&size);CHKERRMPI(ierr);
134   /*ierr = PetscPrintf(PETSC_COMM_WORLD," number of processors = %D\n",size);CHKERRQ(ierr);*/
135   if (size == 1) {
136     PetscViewer viewer;
137     ierr = VecGetArrayRead(s->x, &array);CHKERRQ(ierr);
138     ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer);CHKERRQ(ierr);
139     ierr = PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n");CHKERRQ(ierr);
140     for (j = 0; j < s->ny; j++) {
141       for (i = 0; i < s->nx; i++) {
142         n    = j*s->nx+i;
143         ierr = PetscViewerASCIIPrintf(viewer, "%.12g %.12g %.12g %.12g %.12g\n", (double)(i*s->hx+s->hx/2),(double)(j*s->hy+s->hy/2), (double)PetscRealPart(array[n]), (double)PetscRealPart(array[n+s->nx*s->ny]),(double)PetscRealPart(array[n+2*s->nx*s->ny]));CHKERRQ(ierr);
144       }
145     }
146     ierr = VecRestoreArrayRead(s->x, &array);CHKERRQ(ierr);
147     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
148   }
149   PetscFunctionReturn(0);
150 }
151 
152 PetscErrorCode StokesSetupIndexSets(Stokes *s)
153 {
154   PetscErrorCode ierr;
155 
156   PetscFunctionBeginUser;
157   /* the two index sets */
158   ierr = MatNestGetISs(s->A, s->isg, NULL);CHKERRQ(ierr);
159   /*  ISView(isg[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
160   /*  ISView(isg[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
161   PetscFunctionReturn(0);
162 }
163 
164 PetscErrorCode StokesSetupVectors(Stokes *s)
165 {
166   PetscErrorCode ierr;
167 
168   PetscFunctionBeginUser;
169   /* solution vector x */
170   ierr = VecCreate(PETSC_COMM_WORLD, &s->x);CHKERRQ(ierr);
171   ierr = VecSetSizes(s->x, PETSC_DECIDE, 3*s->nx*s->ny);CHKERRQ(ierr);
172   ierr = VecSetType(s->x, VECMPI);CHKERRQ(ierr);
173   /*  ierr = VecSetRandom(s->x, NULL);CHKERRQ(ierr); */
174   /*  ierr = VecView(s->x, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */
175 
176   /* exact solution y */
177   ierr = VecDuplicate(s->x, &s->y);CHKERRQ(ierr);
178   ierr = StokesExactSolution(s);CHKERRQ(ierr);
179   /*  ierr = VecView(s->y, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */
180 
181   /* rhs vector b */
182   ierr = VecDuplicate(s->x, &s->b);CHKERRQ(ierr);
183   ierr = StokesRhs(s);CHKERRQ(ierr);
184   /*ierr = VecView(s->b, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);*/
185   PetscFunctionReturn(0);
186 }
187 
188 PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j)
189 {
190   PetscInt n;
191 
192   PetscFunctionBeginUser;
193   /* cell number n=j*nx+i has position (i,j) in grid */
194   n  = row%(s->nx*s->ny);
195   *i = n%s->nx;
196   *j = (n-(*i))/s->nx;
197   PetscFunctionReturn(0);
198 }
199 
200 PetscErrorCode StokesExactSolution(Stokes *s)
201 {
202   PetscInt       row, start, end, i, j;
203   PetscScalar    val;
204   Vec            y0,y1;
205   PetscErrorCode ierr;
206 
207   PetscFunctionBeginUser;
208   /* velocity part */
209   ierr = VecGetSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr);
210   ierr = VecGetOwnershipRange(y0, &start, &end);CHKERRQ(ierr);
211   for (row = start; row < end; row++) {
212     ierr = StokesGetPosition(s, row,&i,&j);CHKERRQ(ierr);
213     if (row < s->nx*s->ny) {
214       val = StokesExactVelocityX(j*s->hy+s->hy/2);
215     } else {
216       val = 0;
217     }
218     ierr = VecSetValue(y0, row, val, INSERT_VALUES);CHKERRQ(ierr);
219   }
220   ierr = VecRestoreSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr);
221 
222   /* pressure part */
223   ierr = VecGetSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr);
224   ierr = VecGetOwnershipRange(y1, &start, &end);CHKERRQ(ierr);
225   for (row = start; row < end; row++) {
226     ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr);
227     val  = StokesExactPressure(i*s->hx+s->hx/2);
228     ierr = VecSetValue(y1, row, val, INSERT_VALUES);CHKERRQ(ierr);
229   }
230   ierr = VecRestoreSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr);
231   PetscFunctionReturn(0);
232 }
233 
234 PetscErrorCode StokesRhs(Stokes *s)
235 {
236   PetscInt       row, start, end, i, j;
237   PetscScalar    val;
238   Vec            b0,b1;
239   PetscErrorCode ierr;
240 
241   PetscFunctionBeginUser;
242   /* velocity part */
243   ierr = VecGetSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr);
244   ierr = VecGetOwnershipRange(b0, &start, &end);CHKERRQ(ierr);
245   for (row = start; row < end; row++) {
246     ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr);
247     if (row < s->nx*s->ny) {
248       ierr = StokesRhsMomX(s, i, j, &val);CHKERRQ(ierr);
249     } else {
250       ierr = StokesRhsMomY(s, i, j, &val);CHKERRQ(ierr);
251     }
252     ierr = VecSetValue(b0, row, val, INSERT_VALUES);CHKERRQ(ierr);
253   }
254   ierr = VecRestoreSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr);
255 
256   /* pressure part */
257   ierr = VecGetSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr);
258   ierr = VecGetOwnershipRange(b1, &start, &end);CHKERRQ(ierr);
259   for (row = start; row < end; row++) {
260     ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr);
261     ierr = StokesRhsMass(s, i, j, &val);CHKERRQ(ierr);
262     if (s->matsymmetric) {
263       val = -1.0*val;
264     }
265     ierr = VecSetValue(b1, row, val, INSERT_VALUES);CHKERRQ(ierr);
266   }
267   ierr = VecRestoreSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 PetscErrorCode StokesSetupMatBlock00(Stokes *s)
272 {
273   PetscInt       row, start, end, sz, i, j;
274   PetscInt       cols[5];
275   PetscScalar    vals[5];
276   PetscErrorCode ierr;
277 
278   PetscFunctionBeginUser;
279   /* A[0] is 2N-by-2N */
280   ierr = MatCreate(PETSC_COMM_WORLD,&s->subA[0]);CHKERRQ(ierr);
281   ierr = MatSetOptionsPrefix(s->subA[0],"a00_");CHKERRQ(ierr);
282   ierr = MatSetSizes(s->subA[0],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,2*s->nx*s->ny);CHKERRQ(ierr);
283   ierr = MatSetType(s->subA[0],MATMPIAIJ);CHKERRQ(ierr);
284   ierr = MatMPIAIJSetPreallocation(s->subA[0],5,NULL,5,NULL);CHKERRQ(ierr);
285   ierr = MatGetOwnershipRange(s->subA[0], &start, &end);CHKERRQ(ierr);
286 
287   for (row = start; row < end; row++) {
288     ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr);
289     /* first part: rows 0 to (nx*ny-1) */
290     ierr = StokesStencilLaplacian(s, i, j, &sz, cols, vals);CHKERRQ(ierr);
291     /* second part: rows (nx*ny) to (2*nx*ny-1) */
292     if (row >= s->nx*s->ny) {
293       for (i = 0; i < sz; i++) cols[i] += s->nx*s->ny;
294     }
295     for (i = 0; i < sz; i++) vals[i] = -1.0*vals[i]; /* dynamic viscosity coef mu=-1 */
296     ierr = MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
297   }
298   ierr = MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
299   ierr = MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 PetscErrorCode StokesSetupMatBlock01(Stokes *s)
304 {
305   PetscInt       row, start, end, sz, i, j;
306   PetscInt       cols[5];
307   PetscScalar    vals[5];
308   PetscErrorCode ierr;
309 
310   PetscFunctionBeginUser;
311   /* A[1] is 2N-by-N */
312   ierr = MatCreate(PETSC_COMM_WORLD, &s->subA[1]);CHKERRQ(ierr);
313   ierr = MatSetOptionsPrefix(s->subA[1],"a01_");CHKERRQ(ierr);
314   ierr = MatSetSizes(s->subA[1],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,s->nx*s->ny);CHKERRQ(ierr);
315   ierr = MatSetType(s->subA[1],MATMPIAIJ);CHKERRQ(ierr);
316   ierr = MatMPIAIJSetPreallocation(s->subA[1],5,NULL,5,NULL);CHKERRQ(ierr);
317   ierr = MatGetOwnershipRange(s->subA[1],&start,&end);CHKERRQ(ierr);
318 
319   ierr = MatSetOption(s->subA[1],MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
320 
321   for (row = start; row < end; row++) {
322     ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr);
323     /* first part: rows 0 to (nx*ny-1) */
324     if (row < s->nx*s->ny) {
325       ierr = StokesStencilGradientX(s, i, j, &sz, cols, vals);CHKERRQ(ierr);
326     } else {    /* second part: rows (nx*ny) to (2*nx*ny-1) */
327       ierr = StokesStencilGradientY(s, i, j, &sz, cols, vals);CHKERRQ(ierr);
328     }
329     ierr = MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
330   }
331   ierr = MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
332   ierr = MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333   PetscFunctionReturn(0);
334 }
335 
336 PetscErrorCode StokesSetupMatBlock10(Stokes *s)
337 {
338   PetscErrorCode ierr;
339 
340   PetscFunctionBeginUser;
341   /* A[2] is minus transpose of A[1] */
342   ierr = MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]);CHKERRQ(ierr);
343   if (!s->matsymmetric) {
344     ierr = MatScale(s->subA[2], -1.0);CHKERRQ(ierr);
345   }
346   ierr = MatSetOptionsPrefix(s->subA[2], "a10_");CHKERRQ(ierr);
347   PetscFunctionReturn(0);
348 }
349 
350 PetscErrorCode StokesSetupMatBlock11(Stokes *s)
351 {
352   PetscErrorCode ierr;
353 
354   PetscFunctionBeginUser;
355   /* A[3] is N-by-N null matrix */
356   ierr = MatCreate(PETSC_COMM_WORLD, &s->subA[3]);CHKERRQ(ierr);
357   ierr = MatSetOptionsPrefix(s->subA[3], "a11_");CHKERRQ(ierr);
358   ierr = MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx*s->ny, s->nx*s->ny);CHKERRQ(ierr);
359   ierr = MatSetType(s->subA[3], MATMPIAIJ);CHKERRQ(ierr);
360   ierr = MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL);CHKERRQ(ierr);
361   ierr = MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
362   ierr = MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 PetscErrorCode StokesSetupApproxSchur(Stokes *s)
367 {
368   Vec            diag;
369   PetscErrorCode ierr;
370 
371   PetscFunctionBeginUser;
372   /* Schur complement approximation: myS = A11 - A10 diag(A00)^(-1) A01 */
373   /* note: A11 is zero */
374   /* note: in real life this matrix would be build directly, */
375   /* i.e. without MatMatMult */
376 
377   /* inverse of diagonal of A00 */
378   ierr = VecCreate(PETSC_COMM_WORLD,&diag);CHKERRQ(ierr);
379   ierr = VecSetSizes(diag,PETSC_DECIDE,2*s->nx*s->ny);CHKERRQ(ierr);
380   ierr = VecSetType(diag,VECMPI);CHKERRQ(ierr);
381   ierr = MatGetDiagonal(s->subA[0],diag);CHKERRQ(ierr);
382   ierr = VecReciprocal(diag);CHKERRQ(ierr);
383 
384   /* compute: - A10 diag(A00)^(-1) A01 */
385   ierr = MatDiagonalScale(s->subA[1],diag,NULL);CHKERRQ(ierr); /* (*warning* overwrites subA[1]) */
386   ierr = MatMatMult(s->subA[2],s->subA[1],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&s->myS);CHKERRQ(ierr);
387   ierr = MatScale(s->myS,-1.0);CHKERRQ(ierr);
388 
389   /* restore A10 */
390   ierr = MatGetDiagonal(s->subA[0],diag);CHKERRQ(ierr);
391   ierr = MatDiagonalScale(s->subA[1],diag,NULL);CHKERRQ(ierr);
392   ierr = VecDestroy(&diag);CHKERRQ(ierr);
393   PetscFunctionReturn(0);
394 }
395 
396 PetscErrorCode StokesSetupMatrix(Stokes *s)
397 {
398   PetscErrorCode ierr;
399 
400   PetscFunctionBeginUser;
401   ierr = StokesSetupMatBlock00(s);CHKERRQ(ierr);
402   ierr = StokesSetupMatBlock01(s);CHKERRQ(ierr);
403   ierr = StokesSetupMatBlock10(s);CHKERRQ(ierr);
404   ierr = StokesSetupMatBlock11(s);CHKERRQ(ierr);
405   ierr = MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A);CHKERRQ(ierr);
406   ierr = StokesSetupApproxSchur(s);CHKERRQ(ierr);
407   PetscFunctionReturn(0);
408 }
409 
410 PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
411 {
412   PetscInt    p =j*s->nx+i, w=p-1, e=p+1, s2=p-s->nx, n=p+s->nx;
413   PetscScalar ae=s->hy/s->hx, aeb=0;
414   PetscScalar aw=s->hy/s->hx, awb=s->hy/(s->hx/2);
415   PetscScalar as=s->hx/s->hy, asb=s->hx/(s->hy/2);
416   PetscScalar an=s->hx/s->hy, anb=s->hx/(s->hy/2);
417 
418   PetscFunctionBeginUser;
419   if (i==0 && j==0) { /* south-west corner */
420     *sz  =3;
421     cols[0]=p; vals[0]=-(ae+awb+asb+an);
422     cols[1]=e; vals[1]=ae;
423     cols[2]=n; vals[2]=an;
424   } else if (i==0 && j==s->ny-1) { /* north-west corner */
425     *sz  =3;
426     cols[0]=s2; vals[0]=as;
427     cols[1]=p; vals[1]=-(ae+awb+as+anb);
428     cols[2]=e; vals[2]=ae;
429   } else if (i==s->nx-1 && j==0) { /* south-east corner */
430     *sz  =3;
431     cols[0]=w; vals[0]=aw;
432     cols[1]=p; vals[1]=-(aeb+aw+asb+an);
433     cols[2]=n; vals[2]=an;
434   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
435     *sz  =3;
436     cols[0]=s2; vals[0]=as;
437     cols[1]=w; vals[1]=aw;
438     cols[2]=p; vals[2]=-(aeb+aw+as+anb);
439   } else if (i==0) { /* west boundary */
440     *sz  =4;
441     cols[0]=s2; vals[0]=as;
442     cols[1]=p; vals[1]=-(ae+awb+as+an);
443     cols[2]=e; vals[2]=ae;
444     cols[3]=n; vals[3]=an;
445   } else if (i==s->nx-1) { /* east boundary */
446     *sz  =4;
447     cols[0]=s2; vals[0]=as;
448     cols[1]=w; vals[1]=aw;
449     cols[2]=p; vals[2]=-(aeb+aw+as+an);
450     cols[3]=n; vals[3]=an;
451   } else if (j==0) { /* south boundary */
452     *sz  =4;
453     cols[0]=w; vals[0]=aw;
454     cols[1]=p; vals[1]=-(ae+aw+asb+an);
455     cols[2]=e; vals[2]=ae;
456     cols[3]=n; vals[3]=an;
457   } else if (j==s->ny-1) { /* north boundary */
458     *sz  =4;
459     cols[0]=s2; vals[0]=as;
460     cols[1]=w; vals[1]=aw;
461     cols[2]=p; vals[2]=-(ae+aw+as+anb);
462     cols[3]=e; vals[3]=ae;
463   } else { /* interior */
464     *sz  =5;
465     cols[0]=s2; vals[0]=as;
466     cols[1]=w; vals[1]=aw;
467     cols[2]=p; vals[2]=-(ae+aw+as+an);
468     cols[3]=e; vals[3]=ae;
469     cols[4]=n; vals[4]=an;
470   }
471   PetscFunctionReturn(0);
472 }
473 
474 PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
475 {
476   PetscInt    p =j*s->nx+i, w=p-1, e=p+1;
477   PetscScalar ae= s->hy/2, aeb=s->hy;
478   PetscScalar aw=-s->hy/2, awb=0;
479 
480   PetscFunctionBeginUser;
481   if (i==0 && j==0) { /* south-west corner */
482     *sz  =2;
483     cols[0]=p; vals[0]=-(ae+awb);
484     cols[1]=e; vals[1]=ae;
485   } else if (i==0 && j==s->ny-1) { /* north-west corner */
486     *sz  =2;
487     cols[0]=p; vals[0]=-(ae+awb);
488     cols[1]=e; vals[1]=ae;
489   } else if (i==s->nx-1 && j==0) { /* south-east corner */
490     *sz  =2;
491     cols[0]=w; vals[0]=aw;
492     cols[1]=p; vals[1]=-(aeb+aw);
493   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
494     *sz  =2;
495     cols[0]=w; vals[0]=aw;
496     cols[1]=p; vals[1]=-(aeb+aw);
497   } else if (i==0) { /* west boundary */
498     *sz  =2;
499     cols[0]=p; vals[0]=-(ae+awb);
500     cols[1]=e; vals[1]=ae;
501   } else if (i==s->nx-1) { /* east boundary */
502     *sz  =2;
503     cols[0]=w; vals[0]=aw;
504     cols[1]=p; vals[1]=-(aeb+aw);
505   } else if (j==0) { /* south boundary */
506     *sz  =3;
507     cols[0]=w; vals[0]=aw;
508     cols[1]=p; vals[1]=-(ae+aw);
509     cols[2]=e; vals[2]=ae;
510   } else if (j==s->ny-1) { /* north boundary */
511     *sz  =3;
512     cols[0]=w; vals[0]=aw;
513     cols[1]=p; vals[1]=-(ae+aw);
514     cols[2]=e; vals[2]=ae;
515   } else { /* interior */
516     *sz  =3;
517     cols[0]=w; vals[0]=aw;
518     cols[1]=p; vals[1]=-(ae+aw);
519     cols[2]=e; vals[2]=ae;
520   }
521   PetscFunctionReturn(0);
522 }
523 
524 PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
525 {
526   PetscInt    p =j*s->nx+i, s2=p-s->nx, n=p+s->nx;
527   PetscScalar as=-s->hx/2, asb=0;
528   PetscScalar an= s->hx/2, anb=0;
529 
530   PetscFunctionBeginUser;
531   if (i==0 && j==0) { /* south-west corner */
532     *sz  =2;
533     cols[0]=p; vals[0]=-(asb+an);
534     cols[1]=n; vals[1]=an;
535   } else if (i==0 && j==s->ny-1) { /* north-west corner */
536     *sz  =2;
537     cols[0]=s2; vals[0]=as;
538     cols[1]=p; vals[1]=-(as+anb);
539   } else if (i==s->nx-1 && j==0) { /* south-east corner */
540     *sz  =2;
541     cols[0]=p; vals[0]=-(asb+an);
542     cols[1]=n; vals[1]=an;
543   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
544     *sz  =2;
545     cols[0]=s2; vals[0]=as;
546     cols[1]=p; vals[1]=-(as+anb);
547   } else if (i==0) { /* west boundary */
548     *sz  =3;
549     cols[0]=s2; vals[0]=as;
550     cols[1]=p; vals[1]=-(as+an);
551     cols[2]=n; vals[2]=an;
552   } else if (i==s->nx-1) { /* east boundary */
553     *sz  =3;
554     cols[0]=s2; vals[0]=as;
555     cols[1]=p; vals[1]=-(as+an);
556     cols[2]=n; vals[2]=an;
557   } else if (j==0) { /* south boundary */
558     *sz  =2;
559     cols[0]=p; vals[0]=-(asb+an);
560     cols[1]=n; vals[1]=an;
561   } else if (j==s->ny-1) { /* north boundary */
562     *sz  =2;
563     cols[0]=s2; vals[0]=as;
564     cols[1]=p; vals[1]=-(as+anb);
565   } else { /* interior */
566     *sz  =3;
567     cols[0]=s2; vals[0]=as;
568     cols[1]=p; vals[1]=-(as+an);
569     cols[2]=n; vals[2]=an;
570   }
571   PetscFunctionReturn(0);
572 }
573 
574 PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
575 {
576   PetscScalar y   = j*s->hy+s->hy/2;
577   PetscScalar awb = s->hy/(s->hx/2);
578 
579   PetscFunctionBeginUser;
580   if (i == 0) { /* west boundary */
581     *val = awb*StokesExactVelocityX(y);
582   } else {
583     *val = 0.0;
584   }
585   PetscFunctionReturn(0);
586 }
587 
588 PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
589 {
590   PetscFunctionBeginUser;
591   *val = 0.0;
592   PetscFunctionReturn(0);
593 }
594 
595 PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
596 {
597   PetscScalar y   = j*s->hy+s->hy/2;
598   PetscScalar aeb = s->hy;
599 
600   PetscFunctionBeginUser;
601   if (i == 0) { /* west boundary */
602     *val = aeb*StokesExactVelocityX(y);
603   } else {
604     *val = 0.0;
605   }
606   PetscFunctionReturn(0);
607 }
608 
609 PetscErrorCode StokesCalcResidual(Stokes *s)
610 {
611   PetscReal      val;
612   Vec            b0, b1;
613   PetscErrorCode ierr;
614 
615   PetscFunctionBeginUser;
616   /* residual Ax-b (*warning* overwrites b) */
617   ierr = VecScale(s->b, -1.0);CHKERRQ(ierr);
618   ierr = MatMultAdd(s->A, s->x, s->b, s->b);CHKERRQ(ierr);
619   /*  ierr = VecView(s->b, (PetscViewer)PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */
620 
621   /* residual velocity */
622   ierr = VecGetSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr);
623   ierr = VecNorm(b0, NORM_2, &val);CHKERRQ(ierr);
624   ierr = PetscPrintf(PETSC_COMM_WORLD," residual u = %g\n",(double)val);CHKERRQ(ierr);
625   ierr = VecRestoreSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr);
626 
627   /* residual pressure */
628   ierr = VecGetSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr);
629   ierr = VecNorm(b1, NORM_2, &val);CHKERRQ(ierr);
630   ierr = PetscPrintf(PETSC_COMM_WORLD," residual p = %g\n",(double)val);CHKERRQ(ierr);
631   ierr = VecRestoreSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr);
632 
633   /* total residual */
634   ierr = VecNorm(s->b, NORM_2, &val);CHKERRQ(ierr);
635   ierr = PetscPrintf(PETSC_COMM_WORLD," residual [u,p] = %g\n", (double)val);CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 PetscErrorCode StokesCalcError(Stokes *s)
640 {
641   PetscScalar    scale = PetscSqrtReal((double)s->nx*s->ny);
642   PetscReal      val;
643   Vec            y0, y1;
644   PetscErrorCode ierr;
645 
646   PetscFunctionBeginUser;
647   /* error y-x */
648   ierr = VecAXPY(s->y, -1.0, s->x);CHKERRQ(ierr);
649   /* ierr = VecView(s->y, (PetscViewer)PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */
650 
651   /* error in velocity */
652   ierr = VecGetSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr);
653   ierr = VecNorm(y0, NORM_2, &val);CHKERRQ(ierr);
654   ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error u = %g\n",(double)(PetscRealPart(val/scale)));CHKERRQ(ierr);
655   ierr = VecRestoreSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr);
656 
657   /* error in pressure */
658   ierr = VecGetSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr);
659   ierr = VecNorm(y1, NORM_2, &val);CHKERRQ(ierr);
660   ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error p = %g\n",(double)(PetscRealPart(val/scale)));CHKERRQ(ierr);
661   ierr = VecRestoreSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr);
662 
663   /* total error */
664   ierr = VecNorm(s->y, NORM_2, &val);CHKERRQ(ierr);
665   ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error [u,p] = %g\n", (double)PetscRealPart((val/scale)));CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 int main(int argc, char **argv)
670 {
671   Stokes         s;
672   KSP            ksp;
673   PetscErrorCode ierr;
674 
675   ierr           = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
676   s.nx           = 4;
677   s.ny           = 6;
678   ierr           = PetscOptionsGetInt(NULL,NULL, "-nx", &s.nx, NULL);CHKERRQ(ierr);
679   ierr           = PetscOptionsGetInt(NULL,NULL, "-ny", &s.ny, NULL);CHKERRQ(ierr);
680   s.hx           = 2.0/s.nx;
681   s.hy           = 1.0/s.ny;
682   s.matsymmetric = PETSC_FALSE;
683   ierr           = PetscOptionsGetBool(NULL,NULL, "-mat_set_symmetric", &s.matsymmetric,NULL);CHKERRQ(ierr);
684   s.userPC       = s.userKSP = PETSC_FALSE;
685   ierr           = PetscOptionsHasName(NULL,NULL, "-user_pc", &s.userPC);CHKERRQ(ierr);
686   ierr           = PetscOptionsHasName(NULL,NULL, "-user_ksp", &s.userKSP);CHKERRQ(ierr);
687 
688   ierr = StokesSetupMatrix(&s);CHKERRQ(ierr);
689   ierr = StokesSetupIndexSets(&s);CHKERRQ(ierr);
690   ierr = StokesSetupVectors(&s);CHKERRQ(ierr);
691 
692   ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr);
693   ierr = KSPSetOperators(ksp, s.A, s.A);CHKERRQ(ierr);
694   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
695   ierr = StokesSetupPC(&s, ksp);CHKERRQ(ierr);
696   ierr = KSPSolve(ksp, s.b, s.x);CHKERRQ(ierr);
697 
698   /* don't trust, verify! */
699   ierr = StokesCalcResidual(&s);CHKERRQ(ierr);
700   ierr = StokesCalcError(&s);CHKERRQ(ierr);
701   ierr = StokesWriteSolution(&s);CHKERRQ(ierr);
702 
703   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
704   ierr = MatDestroy(&s.subA[0]);CHKERRQ(ierr);
705   ierr = MatDestroy(&s.subA[1]);CHKERRQ(ierr);
706   ierr = MatDestroy(&s.subA[2]);CHKERRQ(ierr);
707   ierr = MatDestroy(&s.subA[3]);CHKERRQ(ierr);
708   ierr = MatDestroy(&s.A);CHKERRQ(ierr);
709   ierr = VecDestroy(&s.x);CHKERRQ(ierr);
710   ierr = VecDestroy(&s.b);CHKERRQ(ierr);
711   ierr = VecDestroy(&s.y);CHKERRQ(ierr);
712   ierr = MatDestroy(&s.myS);CHKERRQ(ierr);
713   ierr = PetscFinalize();
714   return ierr;
715 }
716 
717 /*TEST
718 
719    test:
720       nsize: 2
721       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none
722 
723    test:
724       suffix: 2
725       nsize: 2
726       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
727 
728    test:
729       suffix: 3
730       nsize: 2
731       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
732 
733    test:
734       suffix: 4
735       nsize: 2
736       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi
737 
738    test:
739       suffix: 4_pcksp
740       nsize: 2
741       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi
742 
743    test:
744       suffix: 5
745       nsize: 2
746       args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type gkb -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-10
747 
748    test:
749       suffix: 6
750       nsize: 2
751       args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type preonly -pc_type fieldsplit -pc_fieldsplit_type gkb -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-10
752 
753    test:
754       suffix: 7
755       nsize: 2
756       args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type preonly -pc_type fieldsplit -pc_fieldsplit_type gkb -pc_fieldsplit_gkb_tol 1e-4 -pc_fieldsplit_gkb_nu 5 -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-6
757 
758 TEST*/
759