xref: /petsc/src/snes/tutorials/ex70.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Poiseuille flow problem. Viscous, laminar flow in a 2D channel with parabolic velocity\n\
2c4762a1bSJed Brown                       profile and linear pressure drop, exact solution of the 2D Stokes\n";
3c4762a1bSJed Brown 
454703344SBarry Smith /*
554703344SBarry Smith      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
654703344SBarry Smith    author : Christiaan M. Klaij
754703344SBarry Smith 
854703344SBarry Smith    Poiseuille flow problem.
954703344SBarry Smith 
1054703344SBarry Smith    Viscous, laminar flow in a 2D channel with parabolic velocity
1154703344SBarry Smith    profile and linear pressure drop, exact solution of the 2D Stokes
1254703344SBarry Smith    equations.
1354703344SBarry Smith 
1454703344SBarry Smith    Discretized with the cell-centered finite-volume method on a
1554703344SBarry Smith    Cartesian grid with co-located variables. Variables ordered as
1654703344SBarry Smith    [u1...uN v1...vN p1...pN]^T. Matrix [A00 A01; A10, A11] solved with
1754703344SBarry Smith    PCFIELDSPLIT. Lower factorization is used to mimic the Semi-Implicit
1854703344SBarry Smith    Method for Pressure Linked Equations (SIMPLE) used as preconditioner
1954703344SBarry Smith    instead of solver.
2054703344SBarry Smith 
2154703344SBarry Smith    Disclaimer: does not contain the pressure-weighed interpolation
2254703344SBarry Smith    method needed to suppress spurious pressure modes in real-life
2354703344SBarry Smith    problems.
2454703344SBarry Smith 
2554703344SBarry Smith    Usage:
2654703344SBarry Smith      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
2754703344SBarry Smith 
2854703344SBarry Smith      Runs with PCFIELDSPLIT on 32x48 grid, no PC for the Schur
2954703344SBarry Smith      complement because A11 is zero. FGMRES is needed because
3054703344SBarry Smith      PCFIELDSPLIT is a variable preconditioner.
3154703344SBarry Smith 
3254703344SBarry Smith      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
3354703344SBarry Smith 
3454703344SBarry Smith      Same as above but with user defined PC for the true Schur
3554703344SBarry Smith      complement. PC based on the SIMPLE-type approximation (inverse of
3654703344SBarry Smith      A00 approximated by inverse of its diagonal).
3754703344SBarry Smith 
3854703344SBarry Smith      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
3954703344SBarry Smith 
4054703344SBarry Smith      Replace the true Schur complement with a user defined Schur
4154703344SBarry Smith      complement based on the SIMPLE-type approximation. Same matrix is
4254703344SBarry Smith      used as PC.
4354703344SBarry Smith 
4454703344SBarry Smith      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
4554703344SBarry Smith 
4654703344SBarry Smith      Out-of-the-box SIMPLE-type preconditioning. The major advantage
4754703344SBarry Smith      is that the user neither needs to provide the approximation of
4854703344SBarry Smith      the Schur complement, nor the corresponding preconditioner.
4954703344SBarry Smith */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown #include <petscksp.h>
52c4762a1bSJed Brown 
53c4762a1bSJed Brown typedef struct {
54c4762a1bSJed Brown   PetscBool userPC, userKSP, matsymmetric; /* user defined preconditioner and matrix for the Schur complement */
55c4762a1bSJed Brown   PetscInt  nx, ny;  /* nb of cells in x- and y-direction */
56c4762a1bSJed Brown   PetscReal hx, hy;  /* mesh size in x- and y-direction */
57c4762a1bSJed Brown   Mat       A;       /* block matrix */
58c4762a1bSJed Brown   Mat       subA[4]; /* the four blocks */
59c4762a1bSJed Brown   Mat       myS;     /* the approximation of the Schur complement */
60c4762a1bSJed Brown   Vec       x, b, y; /* solution, rhs and temporary vector */
61c4762a1bSJed Brown   IS        isg[2];  /* index sets of split "0" and "1" */
62c4762a1bSJed Brown } Stokes;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock00(Stokes*);  /* setup the block Q */
65c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock01(Stokes*);  /* setup the block G */
66c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock10(Stokes*);  /* setup the block D (equal to the transpose of G) */
67c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock11(Stokes*);  /* setup the block C (equal to zero) */
68c4762a1bSJed Brown 
69c4762a1bSJed Brown PetscErrorCode StokesGetPosition(Stokes*, PetscInt, PetscInt*, PetscInt*); /* row number j*nx+i corresponds to position (i,j) in grid */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown PetscErrorCode StokesStencilLaplacian(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Laplacian operator */
72c4762a1bSJed Brown PetscErrorCode StokesStencilGradientX(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (x-component) */
73c4762a1bSJed Brown PetscErrorCode StokesStencilGradientY(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (y-component) */
74c4762a1bSJed Brown 
75c4762a1bSJed Brown PetscErrorCode StokesRhs(Stokes*);                                         /* rhs vector */
76c4762a1bSJed Brown PetscErrorCode StokesRhsMomX(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (x-component) */
77c4762a1bSJed Brown PetscErrorCode StokesRhsMomY(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (y-component) */
78c4762a1bSJed Brown PetscErrorCode StokesRhsMass(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of pressure */
79c4762a1bSJed Brown 
80c4762a1bSJed Brown PetscErrorCode StokesSetupApproxSchur(Stokes*);  /* approximation of the Schur complement */
81c4762a1bSJed Brown 
82c4762a1bSJed Brown PetscErrorCode StokesExactSolution(Stokes*); /* exact solution vector */
83c4762a1bSJed Brown PetscErrorCode StokesWriteSolution(Stokes*); /* write solution to file */
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /* exact solution for the velocity (x-component, y-component is zero) */
86c4762a1bSJed Brown PetscScalar StokesExactVelocityX(const PetscScalar y)
87c4762a1bSJed Brown {
88c4762a1bSJed Brown   return 4.0*y*(1.0-y);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /* exact solution for the pressure */
92c4762a1bSJed Brown PetscScalar StokesExactPressure(const PetscScalar x)
93c4762a1bSJed Brown {
94c4762a1bSJed Brown   return 8.0*(2.0-x);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp)
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   KSP            *subksp;
100c4762a1bSJed Brown   PC             pc;
101c4762a1bSJed Brown   PetscInt       n = 1;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   PetscFunctionBeginUser;
1049566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
1059566063dSJacob Faibussowitsch   PetscCall(PCFieldSplitSetIS(pc, "0", s->isg[0]));
1069566063dSJacob Faibussowitsch   PetscCall(PCFieldSplitSetIS(pc, "1", s->isg[1]));
1071baa6e33SBarry Smith   if (s->userPC) PetscCall(PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS));
108c4762a1bSJed Brown   if (s->userKSP) {
1099566063dSJacob Faibussowitsch     PetscCall(PCSetUp(pc));
1109566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitGetSubKSP(pc, &n, &subksp));
1119566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(subksp[1], s->myS, s->myS));
1129566063dSJacob Faibussowitsch     PetscCall(PetscFree(subksp));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown   PetscFunctionReturn(0);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown PetscErrorCode StokesWriteSolution(Stokes *s)
118c4762a1bSJed Brown {
119c4762a1bSJed Brown   PetscMPIInt       size;
120c4762a1bSJed Brown   PetscInt          n,i,j;
121c4762a1bSJed Brown   const PetscScalar *array;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
124c4762a1bSJed Brown   /* write data (*warning* only works sequential) */
1259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
126c4762a1bSJed Brown   if (size == 1) {
127c4762a1bSJed Brown     PetscViewer viewer;
1289566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(s->x, &array));
1299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer));
1309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n"));
131c4762a1bSJed Brown     for (j = 0; j < s->ny; j++) {
132c4762a1bSJed Brown       for (i = 0; i < s->nx; i++) {
133c4762a1bSJed Brown         n    = j*s->nx+i;
1349566063dSJacob Faibussowitsch         PetscCall(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])));
135c4762a1bSJed Brown       }
136c4762a1bSJed Brown     }
1379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(s->x, &array));
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown   PetscFunctionReturn(0);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown PetscErrorCode StokesSetupIndexSets(Stokes *s)
144c4762a1bSJed Brown {
145c4762a1bSJed Brown   PetscFunctionBeginUser;
146c4762a1bSJed Brown   /* the two index sets */
1479566063dSJacob Faibussowitsch   PetscCall(MatNestGetISs(s->A, s->isg, NULL));
148c4762a1bSJed Brown   PetscFunctionReturn(0);
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown PetscErrorCode StokesSetupVectors(Stokes *s)
152c4762a1bSJed Brown {
153c4762a1bSJed Brown   PetscFunctionBeginUser;
154c4762a1bSJed Brown   /* solution vector x */
1559566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &s->x));
1569566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(s->x, PETSC_DECIDE, 3*s->nx*s->ny));
1579566063dSJacob Faibussowitsch   PetscCall(VecSetType(s->x, VECMPI));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* exact solution y */
1609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s->x, &s->y));
1619566063dSJacob Faibussowitsch   PetscCall(StokesExactSolution(s));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* rhs vector b */
1649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s->x, &s->b));
1659566063dSJacob Faibussowitsch   PetscCall(StokesRhs(s));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   PetscInt n;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   PetscFunctionBeginUser;
174c4762a1bSJed Brown   /* cell number n=j*nx+i has position (i,j) in grid */
175c4762a1bSJed Brown   n  = row%(s->nx*s->ny);
176c4762a1bSJed Brown   *i = n%s->nx;
177c4762a1bSJed Brown   *j = (n-(*i))/s->nx;
178c4762a1bSJed Brown   PetscFunctionReturn(0);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown PetscErrorCode StokesExactSolution(Stokes *s)
182c4762a1bSJed Brown {
183c4762a1bSJed Brown   PetscInt       row, start, end, i, j;
184c4762a1bSJed Brown   PetscScalar    val;
185c4762a1bSJed Brown   Vec            y0,y1;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   PetscFunctionBeginUser;
188c4762a1bSJed Brown   /* velocity part */
1899566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[0], &y0));
1909566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(y0, &start, &end));
191c4762a1bSJed Brown   for (row = start; row < end; row++) {
1929566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row,&i,&j));
193c4762a1bSJed Brown     if (row < s->nx*s->ny) {
194c4762a1bSJed Brown       val = StokesExactVelocityX(j*s->hy+s->hy/2);
195c4762a1bSJed Brown     } else {
196c4762a1bSJed Brown       val = 0;
197c4762a1bSJed Brown     }
1989566063dSJacob Faibussowitsch     PetscCall(VecSetValue(y0, row, val, INSERT_VALUES));
199c4762a1bSJed Brown   }
2009566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* pressure part */
2039566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[1], &y1));
2049566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(y1, &start, &end));
205c4762a1bSJed Brown   for (row = start; row < end; row++) {
2069566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
207c4762a1bSJed Brown     val  = StokesExactPressure(i*s->hx+s->hx/2);
2089566063dSJacob Faibussowitsch     PetscCall(VecSetValue(y1, row, val, INSERT_VALUES));
209c4762a1bSJed Brown   }
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1));
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown PetscErrorCode StokesRhs(Stokes *s)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   PetscInt       row, start, end, i, j;
217c4762a1bSJed Brown   PetscScalar    val;
218c4762a1bSJed Brown   Vec            b0,b1;
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   PetscFunctionBeginUser;
221c4762a1bSJed Brown   /* velocity part */
2229566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[0], &b0));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(b0, &start, &end));
224c4762a1bSJed Brown   for (row = start; row < end; row++) {
2259566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
226c4762a1bSJed Brown     if (row < s->nx*s->ny) {
2279566063dSJacob Faibussowitsch       PetscCall(StokesRhsMomX(s, i, j, &val));
228c4762a1bSJed Brown     } else {
2299566063dSJacob Faibussowitsch       PetscCall(StokesRhsMomY(s, i, j, &val));
230c4762a1bSJed Brown     }
2319566063dSJacob Faibussowitsch     PetscCall(VecSetValue(b0, row, val, INSERT_VALUES));
232c4762a1bSJed Brown   }
2339566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* pressure part */
2369566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[1], &b1));
2379566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(b1, &start, &end));
238c4762a1bSJed Brown   for (row = start; row < end; row++) {
2399566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
2409566063dSJacob Faibussowitsch     PetscCall(StokesRhsMass(s, i, j, &val));
241c4762a1bSJed Brown     if (s->matsymmetric) {
242c4762a1bSJed Brown       val = -1.0*val;
243c4762a1bSJed Brown     }
2449566063dSJacob Faibussowitsch     PetscCall(VecSetValue(b1, row, val, INSERT_VALUES));
245c4762a1bSJed Brown   }
2469566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock00(Stokes *s)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown   PetscInt       row, start, end, sz, i, j;
253c4762a1bSJed Brown   PetscInt       cols[5];
254c4762a1bSJed Brown   PetscScalar    vals[5];
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   PetscFunctionBeginUser;
257c4762a1bSJed Brown   /* A[0] is 2N-by-2N */
2589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&s->subA[0]));
2599566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[0],"a00_"));
2609566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[0],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,2*s->nx*s->ny));
2619566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[0],MATMPIAIJ));
2629566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[0],5,NULL,5,NULL));
2639566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[0], &start, &end));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   for (row = start; row < end; row++) {
2669566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
267c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
2689566063dSJacob Faibussowitsch     PetscCall(StokesStencilLaplacian(s, i, j, &sz, cols, vals));
269c4762a1bSJed Brown     /* second part: rows (nx*ny) to (2*nx*ny-1) */
270c4762a1bSJed Brown     if (row >= s->nx*s->ny) {
271c4762a1bSJed Brown       for (i = 0; i < sz; i++) cols[i] += s->nx*s->ny;
272c4762a1bSJed Brown     }
273c4762a1bSJed Brown     for (i = 0; i < sz; i++) vals[i] = -1.0*vals[i]; /* dynamic viscosity coef mu=-1 */
2749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES));
275c4762a1bSJed Brown   }
2769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY));
2779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY));
278c4762a1bSJed Brown   PetscFunctionReturn(0);
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock01(Stokes *s)
282c4762a1bSJed Brown {
283c4762a1bSJed Brown   PetscInt       row, start, end, sz, i, j;
284c4762a1bSJed Brown   PetscInt       cols[5];
285c4762a1bSJed Brown   PetscScalar    vals[5];
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   PetscFunctionBeginUser;
288c4762a1bSJed Brown   /* A[1] is 2N-by-N */
2899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[1]));
2909566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[1],"a01_"));
2919566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[1],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,s->nx*s->ny));
2929566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[1],MATMPIAIJ));
2939566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[1],5,NULL,5,NULL));
2949566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[1],&start,&end));
295c4762a1bSJed Brown 
2969566063dSJacob Faibussowitsch   PetscCall(MatSetOption(s->subA[1],MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   for (row = start; row < end; row++) {
2999566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
300c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
301c4762a1bSJed Brown     if (row < s->nx*s->ny) {
3029566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientX(s, i, j, &sz, cols, vals));
303c4762a1bSJed Brown     } else {    /* second part: rows (nx*ny) to (2*nx*ny-1) */
3049566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientY(s, i, j, &sz, cols, vals));
305c4762a1bSJed Brown     }
3069566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES));
307c4762a1bSJed Brown   }
3089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY));
3099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY));
310c4762a1bSJed Brown   PetscFunctionReturn(0);
311c4762a1bSJed Brown }
312c4762a1bSJed Brown 
313c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock10(Stokes *s)
314c4762a1bSJed Brown {
315c4762a1bSJed Brown   PetscFunctionBeginUser;
316c4762a1bSJed Brown   /* A[2] is minus transpose of A[1] */
3179566063dSJacob Faibussowitsch   PetscCall(MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]));
318c4762a1bSJed Brown   if (!s->matsymmetric) {
3199566063dSJacob Faibussowitsch     PetscCall(MatScale(s->subA[2], -1.0));
320c4762a1bSJed Brown   }
3219566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[2], "a10_"));
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock11(Stokes *s)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   PetscFunctionBeginUser;
328c4762a1bSJed Brown   /* A[3] is N-by-N null matrix */
3299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[3]));
3309566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[3], "a11_"));
3319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx*s->ny, s->nx*s->ny));
3329566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[3], MATMPIAIJ));
3339566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL));
3349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY));
3359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY));
336c4762a1bSJed Brown   PetscFunctionReturn(0);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
339c4762a1bSJed Brown PetscErrorCode StokesSetupApproxSchur(Stokes *s)
340c4762a1bSJed Brown {
341c4762a1bSJed Brown   Vec            diag;
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   PetscFunctionBeginUser;
34454703344SBarry Smith   /* Schur complement approximation: myS = A11 - A10 inv(DIAGFORM(A00)) A01 */
345c4762a1bSJed Brown   /* note: A11 is zero */
346c4762a1bSJed Brown   /* note: in real life this matrix would be build directly, */
347c4762a1bSJed Brown   /* i.e. without MatMatMult */
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   /* inverse of diagonal of A00 */
3509566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&diag));
3519566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(diag,PETSC_DECIDE,2*s->nx*s->ny));
3529566063dSJacob Faibussowitsch   PetscCall(VecSetType(diag,VECMPI));
3539566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(s->subA[0],diag));
3549566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(diag));
355c4762a1bSJed Brown 
35654703344SBarry Smith   /* compute: - A10 inv(DIAGFORM(A00)) A01 */
3579566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(s->subA[1],diag,NULL)); /* (*warning* overwrites subA[1]) */
3589566063dSJacob Faibussowitsch   PetscCall(MatMatMult(s->subA[2],s->subA[1],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&s->myS));
3599566063dSJacob Faibussowitsch   PetscCall(MatScale(s->myS,-1.0));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* restore A10 */
3629566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(s->subA[0],diag));
3639566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(s->subA[1],diag,NULL));
3649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diag));
365c4762a1bSJed Brown   PetscFunctionReturn(0);
366c4762a1bSJed Brown }
367c4762a1bSJed Brown 
368c4762a1bSJed Brown PetscErrorCode StokesSetupMatrix(Stokes *s)
369c4762a1bSJed Brown {
370c4762a1bSJed Brown   PetscFunctionBeginUser;
3719566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock00(s));
3729566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock01(s));
3739566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock10(s));
3749566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock11(s));
3759566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A));
3769566063dSJacob Faibussowitsch   PetscCall(StokesSetupApproxSchur(s));
377c4762a1bSJed Brown   PetscFunctionReturn(0);
378c4762a1bSJed Brown }
379c4762a1bSJed Brown 
380c4762a1bSJed Brown PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
381c4762a1bSJed Brown {
382c4762a1bSJed Brown   PetscInt    p =j*s->nx+i, w=p-1, e=p+1, s2=p-s->nx, n=p+s->nx;
383c4762a1bSJed Brown   PetscScalar ae=s->hy/s->hx, aeb=0;
384c4762a1bSJed Brown   PetscScalar aw=s->hy/s->hx, awb=s->hy/(s->hx/2);
385c4762a1bSJed Brown   PetscScalar as=s->hx/s->hy, asb=s->hx/(s->hy/2);
386c4762a1bSJed Brown   PetscScalar an=s->hx/s->hy, anb=s->hx/(s->hy/2);
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   PetscFunctionBeginUser;
389c4762a1bSJed Brown   if (i==0 && j==0) { /* south-west corner */
390c4762a1bSJed Brown     *sz  =3;
391c4762a1bSJed Brown     cols[0]=p; vals[0]=-(ae+awb+asb+an);
392c4762a1bSJed Brown     cols[1]=e; vals[1]=ae;
393c4762a1bSJed Brown     cols[2]=n; vals[2]=an;
394c4762a1bSJed Brown   } else if (i==0 && j==s->ny-1) { /* north-west corner */
395c4762a1bSJed Brown     *sz  =3;
396c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
397c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+awb+as+anb);
398c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
399c4762a1bSJed Brown   } else if (i==s->nx-1 && j==0) { /* south-east corner */
400c4762a1bSJed Brown     *sz  =3;
401c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
402c4762a1bSJed Brown     cols[1]=p; vals[1]=-(aeb+aw+asb+an);
403c4762a1bSJed Brown     cols[2]=n; vals[2]=an;
404c4762a1bSJed Brown   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
405c4762a1bSJed Brown     *sz  =3;
406c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
407c4762a1bSJed Brown     cols[1]=w; vals[1]=aw;
408c4762a1bSJed Brown     cols[2]=p; vals[2]=-(aeb+aw+as+anb);
409c4762a1bSJed Brown   } else if (i==0) { /* west boundary */
410c4762a1bSJed Brown     *sz  =4;
411c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
412c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+awb+as+an);
413c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
414c4762a1bSJed Brown     cols[3]=n; vals[3]=an;
415c4762a1bSJed Brown   } else if (i==s->nx-1) { /* east boundary */
416c4762a1bSJed Brown     *sz  =4;
417c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
418c4762a1bSJed Brown     cols[1]=w; vals[1]=aw;
419c4762a1bSJed Brown     cols[2]=p; vals[2]=-(aeb+aw+as+an);
420c4762a1bSJed Brown     cols[3]=n; vals[3]=an;
421c4762a1bSJed Brown   } else if (j==0) { /* south boundary */
422c4762a1bSJed Brown     *sz  =4;
423c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
424c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+aw+asb+an);
425c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
426c4762a1bSJed Brown     cols[3]=n; vals[3]=an;
427c4762a1bSJed Brown   } else if (j==s->ny-1) { /* north boundary */
428c4762a1bSJed Brown     *sz  =4;
429c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
430c4762a1bSJed Brown     cols[1]=w; vals[1]=aw;
431c4762a1bSJed Brown     cols[2]=p; vals[2]=-(ae+aw+as+anb);
432c4762a1bSJed Brown     cols[3]=e; vals[3]=ae;
433c4762a1bSJed Brown   } else { /* interior */
434c4762a1bSJed Brown     *sz  =5;
435c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
436c4762a1bSJed Brown     cols[1]=w; vals[1]=aw;
437c4762a1bSJed Brown     cols[2]=p; vals[2]=-(ae+aw+as+an);
438c4762a1bSJed Brown     cols[3]=e; vals[3]=ae;
439c4762a1bSJed Brown     cols[4]=n; vals[4]=an;
440c4762a1bSJed Brown   }
441c4762a1bSJed Brown   PetscFunctionReturn(0);
442c4762a1bSJed Brown }
443c4762a1bSJed Brown 
444c4762a1bSJed Brown PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
445c4762a1bSJed Brown {
446c4762a1bSJed Brown   PetscInt    p =j*s->nx+i, w=p-1, e=p+1;
447c4762a1bSJed Brown   PetscScalar ae= s->hy/2, aeb=s->hy;
448c4762a1bSJed Brown   PetscScalar aw=-s->hy/2, awb=0;
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   PetscFunctionBeginUser;
451c4762a1bSJed Brown   if (i==0 && j==0) { /* south-west corner */
452c4762a1bSJed Brown     *sz  =2;
453c4762a1bSJed Brown     cols[0]=p; vals[0]=-(ae+awb);
454c4762a1bSJed Brown     cols[1]=e; vals[1]=ae;
455c4762a1bSJed Brown   } else if (i==0 && j==s->ny-1) { /* north-west corner */
456c4762a1bSJed Brown     *sz  =2;
457c4762a1bSJed Brown     cols[0]=p; vals[0]=-(ae+awb);
458c4762a1bSJed Brown     cols[1]=e; vals[1]=ae;
459c4762a1bSJed Brown   } else if (i==s->nx-1 && j==0) { /* south-east corner */
460c4762a1bSJed Brown     *sz  =2;
461c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
462c4762a1bSJed Brown     cols[1]=p; vals[1]=-(aeb+aw);
463c4762a1bSJed Brown   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
464c4762a1bSJed Brown     *sz  =2;
465c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
466c4762a1bSJed Brown     cols[1]=p; vals[1]=-(aeb+aw);
467c4762a1bSJed Brown   } else if (i==0) { /* west boundary */
468c4762a1bSJed Brown     *sz  =2;
469c4762a1bSJed Brown     cols[0]=p; vals[0]=-(ae+awb);
470c4762a1bSJed Brown     cols[1]=e; vals[1]=ae;
471c4762a1bSJed Brown   } else if (i==s->nx-1) { /* east boundary */
472c4762a1bSJed Brown     *sz  =2;
473c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
474c4762a1bSJed Brown     cols[1]=p; vals[1]=-(aeb+aw);
475c4762a1bSJed Brown   } else if (j==0) { /* south boundary */
476c4762a1bSJed Brown     *sz  =3;
477c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
478c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+aw);
479c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
480c4762a1bSJed Brown   } else if (j==s->ny-1) { /* north boundary */
481c4762a1bSJed Brown     *sz  =3;
482c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
483c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+aw);
484c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
485c4762a1bSJed Brown   } else { /* interior */
486c4762a1bSJed Brown     *sz  =3;
487c4762a1bSJed Brown     cols[0]=w; vals[0]=aw;
488c4762a1bSJed Brown     cols[1]=p; vals[1]=-(ae+aw);
489c4762a1bSJed Brown     cols[2]=e; vals[2]=ae;
490c4762a1bSJed Brown   }
491c4762a1bSJed Brown   PetscFunctionReturn(0);
492c4762a1bSJed Brown }
493c4762a1bSJed Brown 
494c4762a1bSJed Brown PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
495c4762a1bSJed Brown {
496c4762a1bSJed Brown   PetscInt    p =j*s->nx+i, s2=p-s->nx, n=p+s->nx;
497c4762a1bSJed Brown   PetscScalar as=-s->hx/2, asb=0;
498c4762a1bSJed Brown   PetscScalar an= s->hx/2, anb=0;
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   PetscFunctionBeginUser;
501c4762a1bSJed Brown   if (i==0 && j==0) { /* south-west corner */
502c4762a1bSJed Brown     *sz  =2;
503c4762a1bSJed Brown     cols[0]=p; vals[0]=-(asb+an);
504c4762a1bSJed Brown     cols[1]=n; vals[1]=an;
505c4762a1bSJed Brown   } else if (i==0 && j==s->ny-1) { /* north-west corner */
506c4762a1bSJed Brown     *sz  =2;
507c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
508c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+anb);
509c4762a1bSJed Brown   } else if (i==s->nx-1 && j==0) { /* south-east corner */
510c4762a1bSJed Brown     *sz  =2;
511c4762a1bSJed Brown     cols[0]=p; vals[0]=-(asb+an);
512c4762a1bSJed Brown     cols[1]=n; vals[1]=an;
513c4762a1bSJed Brown   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
514c4762a1bSJed Brown     *sz  =2;
515c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
516c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+anb);
517c4762a1bSJed Brown   } else if (i==0) { /* west boundary */
518c4762a1bSJed Brown     *sz  =3;
519c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
520c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+an);
521c4762a1bSJed Brown     cols[2]=n; vals[2]=an;
522c4762a1bSJed Brown   } else if (i==s->nx-1) { /* east boundary */
523c4762a1bSJed Brown     *sz  =3;
524c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
525c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+an);
526c4762a1bSJed Brown     cols[2]=n; vals[2]=an;
527c4762a1bSJed Brown   } else if (j==0) { /* south boundary */
528c4762a1bSJed Brown     *sz  =2;
529c4762a1bSJed Brown     cols[0]=p; vals[0]=-(asb+an);
530c4762a1bSJed Brown     cols[1]=n; vals[1]=an;
531c4762a1bSJed Brown   } else if (j==s->ny-1) { /* north boundary */
532c4762a1bSJed Brown     *sz  =2;
533c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
534c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+anb);
535c4762a1bSJed Brown   } else { /* interior */
536c4762a1bSJed Brown     *sz  =3;
537c4762a1bSJed Brown     cols[0]=s2; vals[0]=as;
538c4762a1bSJed Brown     cols[1]=p; vals[1]=-(as+an);
539c4762a1bSJed Brown     cols[2]=n; vals[2]=an;
540c4762a1bSJed Brown   }
541c4762a1bSJed Brown   PetscFunctionReturn(0);
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
544c4762a1bSJed Brown PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
545c4762a1bSJed Brown {
546c4762a1bSJed Brown   PetscScalar y   = j*s->hy+s->hy/2;
547c4762a1bSJed Brown   PetscScalar awb = s->hy/(s->hx/2);
548c4762a1bSJed Brown 
549c4762a1bSJed Brown   PetscFunctionBeginUser;
550c4762a1bSJed Brown   if (i == 0) { /* west boundary */
551c4762a1bSJed Brown     *val = awb*StokesExactVelocityX(y);
552c4762a1bSJed Brown   } else {
553c4762a1bSJed Brown     *val = 0.0;
554c4762a1bSJed Brown   }
555c4762a1bSJed Brown   PetscFunctionReturn(0);
556c4762a1bSJed Brown }
557c4762a1bSJed Brown 
558c4762a1bSJed Brown PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
559c4762a1bSJed Brown {
560c4762a1bSJed Brown   PetscFunctionBeginUser;
561c4762a1bSJed Brown   *val = 0.0;
562c4762a1bSJed Brown   PetscFunctionReturn(0);
563c4762a1bSJed Brown }
564c4762a1bSJed Brown 
565c4762a1bSJed Brown PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
566c4762a1bSJed Brown {
567c4762a1bSJed Brown   PetscScalar y   = j*s->hy+s->hy/2;
568c4762a1bSJed Brown   PetscScalar aeb = s->hy;
569c4762a1bSJed Brown 
570c4762a1bSJed Brown   PetscFunctionBeginUser;
571c4762a1bSJed Brown   if (i == 0) { /* west boundary */
572c4762a1bSJed Brown     *val = aeb*StokesExactVelocityX(y);
573c4762a1bSJed Brown   } else {
574c4762a1bSJed Brown     *val = 0.0;
575c4762a1bSJed Brown   }
576c4762a1bSJed Brown   PetscFunctionReturn(0);
577c4762a1bSJed Brown }
578c4762a1bSJed Brown 
579c4762a1bSJed Brown PetscErrorCode StokesCalcResidual(Stokes *s)
580c4762a1bSJed Brown {
581c4762a1bSJed Brown   PetscReal      val;
582c4762a1bSJed Brown   Vec            b0, b1;
583c4762a1bSJed Brown 
584c4762a1bSJed Brown   PetscFunctionBeginUser;
585c4762a1bSJed Brown   /* residual Ax-b (*warning* overwrites b) */
5869566063dSJacob Faibussowitsch   PetscCall(VecScale(s->b, -1.0));
5879566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(s->A, s->x, s->b, s->b));
588c4762a1bSJed Brown 
589c4762a1bSJed Brown   /* residual velocity */
5909566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[0], &b0));
5919566063dSJacob Faibussowitsch   PetscCall(VecNorm(b0, NORM_2, &val));
5929566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," residual u = %g\n",(double)val));
5939566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0));
594c4762a1bSJed Brown 
595c4762a1bSJed Brown   /* residual pressure */
5969566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[1], &b1));
5979566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &val));
5989566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," residual p = %g\n",(double)val));
5999566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
600c4762a1bSJed Brown 
601c4762a1bSJed Brown   /* total residual */
6029566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->b, NORM_2, &val));
6039566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," residual [u,p] = %g\n", (double)val));
604c4762a1bSJed Brown   PetscFunctionReturn(0);
605c4762a1bSJed Brown }
606c4762a1bSJed Brown 
607c4762a1bSJed Brown PetscErrorCode StokesCalcError(Stokes *s)
608c4762a1bSJed Brown {
609c4762a1bSJed Brown   PetscScalar    scale = PetscSqrtReal((double)s->nx*s->ny);
610c4762a1bSJed Brown   PetscReal      val;
611c4762a1bSJed Brown   Vec            y0, y1;
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   PetscFunctionBeginUser;
614c4762a1bSJed Brown   /* error y-x */
6159566063dSJacob Faibussowitsch   PetscCall(VecAXPY(s->y, -1.0, s->x));
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   /* error in velocity */
6189566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[0], &y0));
6199566063dSJacob Faibussowitsch   PetscCall(VecNorm(y0, NORM_2, &val));
6209566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," discretization error u = %g\n",(double)(PetscRealPart(val/scale))));
6219566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0));
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /* error in pressure */
6249566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[1], &y1));
6259566063dSJacob Faibussowitsch   PetscCall(VecNorm(y1, NORM_2, &val));
6269566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," discretization error p = %g\n",(double)(PetscRealPart(val/scale))));
6279566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1));
628c4762a1bSJed Brown 
629c4762a1bSJed Brown   /* total error */
6309566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->y, NORM_2, &val));
6319566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," discretization error [u,p] = %g\n", (double)PetscRealPart((val/scale))));
632c4762a1bSJed Brown   PetscFunctionReturn(0);
633c4762a1bSJed Brown }
634c4762a1bSJed Brown 
635c4762a1bSJed Brown int main(int argc, char **argv)
636c4762a1bSJed Brown {
637c4762a1bSJed Brown   Stokes         s;
638c4762a1bSJed Brown   KSP            ksp;
639c4762a1bSJed Brown 
640*327415f7SBarry Smith   PetscFunctionBeginUser;
6419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
642c4762a1bSJed Brown   s.nx           = 4;
643c4762a1bSJed Brown   s.ny           = 6;
6449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-nx", &s.nx, NULL));
6459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-ny", &s.ny, NULL));
646c4762a1bSJed Brown   s.hx           = 2.0/s.nx;
647c4762a1bSJed Brown   s.hy           = 1.0/s.ny;
648c4762a1bSJed Brown   s.matsymmetric = PETSC_FALSE;
6499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL, "-mat_set_symmetric", &s.matsymmetric,NULL));
650c4762a1bSJed Brown   s.userPC       = s.userKSP = PETSC_FALSE;
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL, "-user_pc", &s.userPC));
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL, "-user_ksp", &s.userKSP));
653c4762a1bSJed Brown 
6549566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatrix(&s));
6559566063dSJacob Faibussowitsch   PetscCall(StokesSetupIndexSets(&s));
6569566063dSJacob Faibussowitsch   PetscCall(StokesSetupVectors(&s));
657c4762a1bSJed Brown 
6589566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
6599566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, s.A, s.A));
6609566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
6619566063dSJacob Faibussowitsch   PetscCall(StokesSetupPC(&s, ksp));
6629566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, s.b, s.x));
663c4762a1bSJed Brown 
664c4762a1bSJed Brown   /* don't trust, verify! */
6659566063dSJacob Faibussowitsch   PetscCall(StokesCalcResidual(&s));
6669566063dSJacob Faibussowitsch   PetscCall(StokesCalcError(&s));
6679566063dSJacob Faibussowitsch   PetscCall(StokesWriteSolution(&s));
668c4762a1bSJed Brown 
6699566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
6709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[0]));
6719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[1]));
6729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[2]));
6739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[3]));
6749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.A));
6759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.x));
6769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.b));
6779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.y));
6789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.myS));
6799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
680b122ec5aSJacob Faibussowitsch   return 0;
681c4762a1bSJed Brown }
682c4762a1bSJed Brown 
683c4762a1bSJed Brown /*TEST
684c4762a1bSJed Brown 
685c4762a1bSJed Brown    test:
686c4762a1bSJed Brown       nsize: 2
687c4762a1bSJed Brown       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
688c4762a1bSJed Brown 
689c4762a1bSJed Brown    test:
690c4762a1bSJed Brown       suffix: 2
691c4762a1bSJed Brown       nsize: 2
692c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
693c4762a1bSJed Brown 
694c4762a1bSJed Brown    test:
695c4762a1bSJed Brown       suffix: 3
696c4762a1bSJed Brown       nsize: 2
697c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
698c4762a1bSJed Brown 
699c4762a1bSJed Brown    test:
700c4762a1bSJed Brown       suffix: 4
701c4762a1bSJed Brown       nsize: 2
702c4762a1bSJed Brown       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
703c4762a1bSJed Brown 
704c4762a1bSJed Brown    test:
705c4762a1bSJed Brown       suffix: 4_pcksp
706c4762a1bSJed Brown       nsize: 2
707c4762a1bSJed Brown       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
708c4762a1bSJed Brown 
709c4762a1bSJed Brown    test:
710c4762a1bSJed Brown       suffix: 5
711c4762a1bSJed Brown       nsize: 2
712c4762a1bSJed Brown       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
713c4762a1bSJed Brown 
714c4762a1bSJed Brown    test:
715c4762a1bSJed Brown       suffix: 6
716c4762a1bSJed Brown       nsize: 2
717c4762a1bSJed Brown       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
718c4762a1bSJed Brown 
719c4762a1bSJed Brown    test:
720c4762a1bSJed Brown       suffix: 7
721c4762a1bSJed Brown       nsize: 2
722c4762a1bSJed Brown       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
723c4762a1bSJed Brown 
724c4762a1bSJed Brown TEST*/
725