xref: /petsc/src/snes/tutorials/ex70.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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) */
86d71ae5a4SJacob Faibussowitsch PetscScalar StokesExactVelocityX(const PetscScalar y)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown   return 4.0 * y * (1.0 - y);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /* exact solution for the pressure */
92d71ae5a4SJacob Faibussowitsch PetscScalar StokesExactPressure(const PetscScalar x)
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   return 8.0 * (2.0 - x);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp)
98d71ae5a4SJacob Faibussowitsch {
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   }
114*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesWriteSolution(Stokes *s)
118d71ae5a4SJacob Faibussowitsch {
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   }
140*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupIndexSets(Stokes *s)
144d71ae5a4SJacob Faibussowitsch {
145c4762a1bSJed Brown   PetscFunctionBeginUser;
146c4762a1bSJed Brown   /* the two index sets */
1479566063dSJacob Faibussowitsch   PetscCall(MatNestGetISs(s->A, s->isg, NULL));
148*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupVectors(Stokes *s)
152d71ae5a4SJacob Faibussowitsch {
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));
166*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j)
170d71ae5a4SJacob Faibussowitsch {
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;
178*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesExactSolution(Stokes *s)
182d71ae5a4SJacob Faibussowitsch {
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));
211*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhs(Stokes *s)
215d71ae5a4SJacob Faibussowitsch {
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));
241ad540459SPierre Jolivet     if (s->matsymmetric) val = -1.0 * val;
2429566063dSJacob Faibussowitsch     PetscCall(VecSetValue(b1, row, val, INSERT_VALUES));
243c4762a1bSJed Brown   }
2449566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
245*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246c4762a1bSJed Brown }
247c4762a1bSJed Brown 
248d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock00(Stokes *s)
249d71ae5a4SJacob Faibussowitsch {
250c4762a1bSJed Brown   PetscInt    row, start, end, sz, i, j;
251c4762a1bSJed Brown   PetscInt    cols[5];
252c4762a1bSJed Brown   PetscScalar vals[5];
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   PetscFunctionBeginUser;
255c4762a1bSJed Brown   /* A[0] is 2N-by-2N */
2569566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[0]));
2579566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[0], "a00_"));
2589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[0], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, 2 * s->nx * s->ny));
2599566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[0], MATMPIAIJ));
2609566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[0], 5, NULL, 5, NULL));
2619566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[0], &start, &end));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   for (row = start; row < end; row++) {
2649566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
265c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
2669566063dSJacob Faibussowitsch     PetscCall(StokesStencilLaplacian(s, i, j, &sz, cols, vals));
267c4762a1bSJed Brown     /* second part: rows (nx*ny) to (2*nx*ny-1) */
268c4762a1bSJed Brown     if (row >= s->nx * s->ny) {
269c4762a1bSJed Brown       for (i = 0; i < sz; i++) cols[i] += s->nx * s->ny;
270c4762a1bSJed Brown     }
271c4762a1bSJed Brown     for (i = 0; i < sz; i++) vals[i] = -1.0 * vals[i]; /* dynamic viscosity coef mu=-1 */
2729566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES));
273c4762a1bSJed Brown   }
2749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY));
2759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY));
276*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277c4762a1bSJed Brown }
278c4762a1bSJed Brown 
279d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock01(Stokes *s)
280d71ae5a4SJacob Faibussowitsch {
281c4762a1bSJed Brown   PetscInt    row, start, end, sz, i, j;
282c4762a1bSJed Brown   PetscInt    cols[5];
283c4762a1bSJed Brown   PetscScalar vals[5];
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   PetscFunctionBeginUser;
286c4762a1bSJed Brown   /* A[1] is 2N-by-N */
2879566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[1]));
2889566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[1], "a01_"));
2899566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[1], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, s->nx * s->ny));
2909566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[1], MATMPIAIJ));
2919566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[1], 5, NULL, 5, NULL));
2929566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[1], &start, &end));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(MatSetOption(s->subA[1], MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   for (row = start; row < end; row++) {
2979566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
298c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
299c4762a1bSJed Brown     if (row < s->nx * s->ny) {
3009566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientX(s, i, j, &sz, cols, vals));
301c4762a1bSJed Brown     } else { /* second part: rows (nx*ny) to (2*nx*ny-1) */
3029566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientY(s, i, j, &sz, cols, vals));
303c4762a1bSJed Brown     }
3049566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES));
305c4762a1bSJed Brown   }
3069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY));
3079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY));
308*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
311d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock10(Stokes *s)
312d71ae5a4SJacob Faibussowitsch {
313c4762a1bSJed Brown   PetscFunctionBeginUser;
314c4762a1bSJed Brown   /* A[2] is minus transpose of A[1] */
3159566063dSJacob Faibussowitsch   PetscCall(MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]));
31648a46eb9SPierre Jolivet   if (!s->matsymmetric) PetscCall(MatScale(s->subA[2], -1.0));
3179566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[2], "a10_"));
318*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock11(Stokes *s)
322d71ae5a4SJacob Faibussowitsch {
323c4762a1bSJed Brown   PetscFunctionBeginUser;
324c4762a1bSJed Brown   /* A[3] is N-by-N null matrix */
3259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[3]));
3269566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[3], "a11_"));
3279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx * s->ny, s->nx * s->ny));
3289566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[3], MATMPIAIJ));
3299566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL));
3309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY));
3319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY));
332*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown 
335d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupApproxSchur(Stokes *s)
336d71ae5a4SJacob Faibussowitsch {
337c4762a1bSJed Brown   Vec diag;
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   PetscFunctionBeginUser;
34054703344SBarry Smith   /* Schur complement approximation: myS = A11 - A10 inv(DIAGFORM(A00)) A01 */
341c4762a1bSJed Brown   /* note: A11 is zero */
342c4762a1bSJed Brown   /* note: in real life this matrix would be build directly, */
343c4762a1bSJed Brown   /* i.e. without MatMatMult */
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   /* inverse of diagonal of A00 */
3469566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &diag));
3479566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(diag, PETSC_DECIDE, 2 * s->nx * s->ny));
3489566063dSJacob Faibussowitsch   PetscCall(VecSetType(diag, VECMPI));
3499566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(s->subA[0], diag));
3509566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(diag));
351c4762a1bSJed Brown 
35254703344SBarry Smith   /* compute: - A10 inv(DIAGFORM(A00)) A01 */
3539566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(s->subA[1], diag, NULL)); /* (*warning* overwrites subA[1]) */
3549566063dSJacob Faibussowitsch   PetscCall(MatMatMult(s->subA[2], s->subA[1], MAT_INITIAL_MATRIX, PETSC_DEFAULT, &s->myS));
3559566063dSJacob Faibussowitsch   PetscCall(MatScale(s->myS, -1.0));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /* restore A10 */
3589566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(s->subA[0], diag));
3599566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(s->subA[1], diag, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diag));
361*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
364d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatrix(Stokes *s)
365d71ae5a4SJacob Faibussowitsch {
366c4762a1bSJed Brown   PetscFunctionBeginUser;
3679566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock00(s));
3689566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock01(s));
3699566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock10(s));
3709566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatBlock11(s));
3719566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A));
3729566063dSJacob Faibussowitsch   PetscCall(StokesSetupApproxSchur(s));
373*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
376d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
377d71ae5a4SJacob Faibussowitsch {
378c4762a1bSJed Brown   PetscInt    p = j * s->nx + i, w = p - 1, e = p + 1, s2 = p - s->nx, n = p + s->nx;
379c4762a1bSJed Brown   PetscScalar ae = s->hy / s->hx, aeb = 0;
380c4762a1bSJed Brown   PetscScalar aw = s->hy / s->hx, awb = s->hy / (s->hx / 2);
381c4762a1bSJed Brown   PetscScalar as = s->hx / s->hy, asb = s->hx / (s->hy / 2);
382c4762a1bSJed Brown   PetscScalar an = s->hx / s->hy, anb = s->hx / (s->hy / 2);
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   PetscFunctionBeginUser;
385c4762a1bSJed Brown   if (i == 0 && j == 0) { /* south-west corner */
386c4762a1bSJed Brown     *sz     = 3;
3879371c9d4SSatish Balay     cols[0] = p;
3889371c9d4SSatish Balay     vals[0] = -(ae + awb + asb + an);
3899371c9d4SSatish Balay     cols[1] = e;
3909371c9d4SSatish Balay     vals[1] = ae;
3919371c9d4SSatish Balay     cols[2] = n;
3929371c9d4SSatish Balay     vals[2] = an;
393c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
394c4762a1bSJed Brown     *sz     = 3;
3959371c9d4SSatish Balay     cols[0] = s2;
3969371c9d4SSatish Balay     vals[0] = as;
3979371c9d4SSatish Balay     cols[1] = p;
3989371c9d4SSatish Balay     vals[1] = -(ae + awb + as + anb);
3999371c9d4SSatish Balay     cols[2] = e;
4009371c9d4SSatish Balay     vals[2] = ae;
401c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
402c4762a1bSJed Brown     *sz     = 3;
4039371c9d4SSatish Balay     cols[0] = w;
4049371c9d4SSatish Balay     vals[0] = aw;
4059371c9d4SSatish Balay     cols[1] = p;
4069371c9d4SSatish Balay     vals[1] = -(aeb + aw + asb + an);
4079371c9d4SSatish Balay     cols[2] = n;
4089371c9d4SSatish Balay     vals[2] = an;
409c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
410c4762a1bSJed Brown     *sz     = 3;
4119371c9d4SSatish Balay     cols[0] = s2;
4129371c9d4SSatish Balay     vals[0] = as;
4139371c9d4SSatish Balay     cols[1] = w;
4149371c9d4SSatish Balay     vals[1] = aw;
4159371c9d4SSatish Balay     cols[2] = p;
4169371c9d4SSatish Balay     vals[2] = -(aeb + aw + as + anb);
417c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
418c4762a1bSJed Brown     *sz     = 4;
4199371c9d4SSatish Balay     cols[0] = s2;
4209371c9d4SSatish Balay     vals[0] = as;
4219371c9d4SSatish Balay     cols[1] = p;
4229371c9d4SSatish Balay     vals[1] = -(ae + awb + as + an);
4239371c9d4SSatish Balay     cols[2] = e;
4249371c9d4SSatish Balay     vals[2] = ae;
4259371c9d4SSatish Balay     cols[3] = n;
4269371c9d4SSatish Balay     vals[3] = an;
427c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
428c4762a1bSJed Brown     *sz     = 4;
4299371c9d4SSatish Balay     cols[0] = s2;
4309371c9d4SSatish Balay     vals[0] = as;
4319371c9d4SSatish Balay     cols[1] = w;
4329371c9d4SSatish Balay     vals[1] = aw;
4339371c9d4SSatish Balay     cols[2] = p;
4349371c9d4SSatish Balay     vals[2] = -(aeb + aw + as + an);
4359371c9d4SSatish Balay     cols[3] = n;
4369371c9d4SSatish Balay     vals[3] = an;
437c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
438c4762a1bSJed Brown     *sz     = 4;
4399371c9d4SSatish Balay     cols[0] = w;
4409371c9d4SSatish Balay     vals[0] = aw;
4419371c9d4SSatish Balay     cols[1] = p;
4429371c9d4SSatish Balay     vals[1] = -(ae + aw + asb + an);
4439371c9d4SSatish Balay     cols[2] = e;
4449371c9d4SSatish Balay     vals[2] = ae;
4459371c9d4SSatish Balay     cols[3] = n;
4469371c9d4SSatish Balay     vals[3] = an;
447c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
448c4762a1bSJed Brown     *sz     = 4;
4499371c9d4SSatish Balay     cols[0] = s2;
4509371c9d4SSatish Balay     vals[0] = as;
4519371c9d4SSatish Balay     cols[1] = w;
4529371c9d4SSatish Balay     vals[1] = aw;
4539371c9d4SSatish Balay     cols[2] = p;
4549371c9d4SSatish Balay     vals[2] = -(ae + aw + as + anb);
4559371c9d4SSatish Balay     cols[3] = e;
4569371c9d4SSatish Balay     vals[3] = ae;
457c4762a1bSJed Brown   } else { /* interior */
458c4762a1bSJed Brown     *sz     = 5;
4599371c9d4SSatish Balay     cols[0] = s2;
4609371c9d4SSatish Balay     vals[0] = as;
4619371c9d4SSatish Balay     cols[1] = w;
4629371c9d4SSatish Balay     vals[1] = aw;
4639371c9d4SSatish Balay     cols[2] = p;
4649371c9d4SSatish Balay     vals[2] = -(ae + aw + as + an);
4659371c9d4SSatish Balay     cols[3] = e;
4669371c9d4SSatish Balay     vals[3] = ae;
4679371c9d4SSatish Balay     cols[4] = n;
4689371c9d4SSatish Balay     vals[4] = an;
469c4762a1bSJed Brown   }
470*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471c4762a1bSJed Brown }
472c4762a1bSJed Brown 
473d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
474d71ae5a4SJacob Faibussowitsch {
475c4762a1bSJed Brown   PetscInt    p = j * s->nx + i, w = p - 1, e = p + 1;
476c4762a1bSJed Brown   PetscScalar ae = s->hy / 2, aeb = s->hy;
477c4762a1bSJed Brown   PetscScalar aw = -s->hy / 2, awb = 0;
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   PetscFunctionBeginUser;
480c4762a1bSJed Brown   if (i == 0 && j == 0) { /* south-west corner */
481c4762a1bSJed Brown     *sz     = 2;
4829371c9d4SSatish Balay     cols[0] = p;
4839371c9d4SSatish Balay     vals[0] = -(ae + awb);
4849371c9d4SSatish Balay     cols[1] = e;
4859371c9d4SSatish Balay     vals[1] = ae;
486c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
487c4762a1bSJed Brown     *sz     = 2;
4889371c9d4SSatish Balay     cols[0] = p;
4899371c9d4SSatish Balay     vals[0] = -(ae + awb);
4909371c9d4SSatish Balay     cols[1] = e;
4919371c9d4SSatish Balay     vals[1] = ae;
492c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
493c4762a1bSJed Brown     *sz     = 2;
4949371c9d4SSatish Balay     cols[0] = w;
4959371c9d4SSatish Balay     vals[0] = aw;
4969371c9d4SSatish Balay     cols[1] = p;
4979371c9d4SSatish Balay     vals[1] = -(aeb + aw);
498c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
499c4762a1bSJed Brown     *sz     = 2;
5009371c9d4SSatish Balay     cols[0] = w;
5019371c9d4SSatish Balay     vals[0] = aw;
5029371c9d4SSatish Balay     cols[1] = p;
5039371c9d4SSatish Balay     vals[1] = -(aeb + aw);
504c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
505c4762a1bSJed Brown     *sz     = 2;
5069371c9d4SSatish Balay     cols[0] = p;
5079371c9d4SSatish Balay     vals[0] = -(ae + awb);
5089371c9d4SSatish Balay     cols[1] = e;
5099371c9d4SSatish Balay     vals[1] = ae;
510c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
511c4762a1bSJed Brown     *sz     = 2;
5129371c9d4SSatish Balay     cols[0] = w;
5139371c9d4SSatish Balay     vals[0] = aw;
5149371c9d4SSatish Balay     cols[1] = p;
5159371c9d4SSatish Balay     vals[1] = -(aeb + aw);
516c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
517c4762a1bSJed Brown     *sz     = 3;
5189371c9d4SSatish Balay     cols[0] = w;
5199371c9d4SSatish Balay     vals[0] = aw;
5209371c9d4SSatish Balay     cols[1] = p;
5219371c9d4SSatish Balay     vals[1] = -(ae + aw);
5229371c9d4SSatish Balay     cols[2] = e;
5239371c9d4SSatish Balay     vals[2] = ae;
524c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
525c4762a1bSJed Brown     *sz     = 3;
5269371c9d4SSatish Balay     cols[0] = w;
5279371c9d4SSatish Balay     vals[0] = aw;
5289371c9d4SSatish Balay     cols[1] = p;
5299371c9d4SSatish Balay     vals[1] = -(ae + aw);
5309371c9d4SSatish Balay     cols[2] = e;
5319371c9d4SSatish Balay     vals[2] = ae;
532c4762a1bSJed Brown   } else { /* interior */
533c4762a1bSJed Brown     *sz     = 3;
5349371c9d4SSatish Balay     cols[0] = w;
5359371c9d4SSatish Balay     vals[0] = aw;
5369371c9d4SSatish Balay     cols[1] = p;
5379371c9d4SSatish Balay     vals[1] = -(ae + aw);
5389371c9d4SSatish Balay     cols[2] = e;
5399371c9d4SSatish Balay     vals[2] = ae;
540c4762a1bSJed Brown   }
541*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
544d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
545d71ae5a4SJacob Faibussowitsch {
546c4762a1bSJed Brown   PetscInt    p = j * s->nx + i, s2 = p - s->nx, n = p + s->nx;
547c4762a1bSJed Brown   PetscScalar as = -s->hx / 2, asb = 0;
548c4762a1bSJed Brown   PetscScalar an = s->hx / 2, anb = 0;
549c4762a1bSJed Brown 
550c4762a1bSJed Brown   PetscFunctionBeginUser;
551c4762a1bSJed Brown   if (i == 0 && j == 0) { /* south-west corner */
552c4762a1bSJed Brown     *sz     = 2;
5539371c9d4SSatish Balay     cols[0] = p;
5549371c9d4SSatish Balay     vals[0] = -(asb + an);
5559371c9d4SSatish Balay     cols[1] = n;
5569371c9d4SSatish Balay     vals[1] = an;
557c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
558c4762a1bSJed Brown     *sz     = 2;
5599371c9d4SSatish Balay     cols[0] = s2;
5609371c9d4SSatish Balay     vals[0] = as;
5619371c9d4SSatish Balay     cols[1] = p;
5629371c9d4SSatish Balay     vals[1] = -(as + anb);
563c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
564c4762a1bSJed Brown     *sz     = 2;
5659371c9d4SSatish Balay     cols[0] = p;
5669371c9d4SSatish Balay     vals[0] = -(asb + an);
5679371c9d4SSatish Balay     cols[1] = n;
5689371c9d4SSatish Balay     vals[1] = an;
569c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
570c4762a1bSJed Brown     *sz     = 2;
5719371c9d4SSatish Balay     cols[0] = s2;
5729371c9d4SSatish Balay     vals[0] = as;
5739371c9d4SSatish Balay     cols[1] = p;
5749371c9d4SSatish Balay     vals[1] = -(as + anb);
575c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
576c4762a1bSJed Brown     *sz     = 3;
5779371c9d4SSatish Balay     cols[0] = s2;
5789371c9d4SSatish Balay     vals[0] = as;
5799371c9d4SSatish Balay     cols[1] = p;
5809371c9d4SSatish Balay     vals[1] = -(as + an);
5819371c9d4SSatish Balay     cols[2] = n;
5829371c9d4SSatish Balay     vals[2] = an;
583c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
584c4762a1bSJed Brown     *sz     = 3;
5859371c9d4SSatish Balay     cols[0] = s2;
5869371c9d4SSatish Balay     vals[0] = as;
5879371c9d4SSatish Balay     cols[1] = p;
5889371c9d4SSatish Balay     vals[1] = -(as + an);
5899371c9d4SSatish Balay     cols[2] = n;
5909371c9d4SSatish Balay     vals[2] = an;
591c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
592c4762a1bSJed Brown     *sz     = 2;
5939371c9d4SSatish Balay     cols[0] = p;
5949371c9d4SSatish Balay     vals[0] = -(asb + an);
5959371c9d4SSatish Balay     cols[1] = n;
5969371c9d4SSatish Balay     vals[1] = an;
597c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
598c4762a1bSJed Brown     *sz     = 2;
5999371c9d4SSatish Balay     cols[0] = s2;
6009371c9d4SSatish Balay     vals[0] = as;
6019371c9d4SSatish Balay     cols[1] = p;
6029371c9d4SSatish Balay     vals[1] = -(as + anb);
603c4762a1bSJed Brown   } else { /* interior */
604c4762a1bSJed Brown     *sz     = 3;
6059371c9d4SSatish Balay     cols[0] = s2;
6069371c9d4SSatish Balay     vals[0] = as;
6079371c9d4SSatish Balay     cols[1] = p;
6089371c9d4SSatish Balay     vals[1] = -(as + an);
6099371c9d4SSatish Balay     cols[2] = n;
6109371c9d4SSatish Balay     vals[2] = an;
611c4762a1bSJed Brown   }
612*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
613c4762a1bSJed Brown }
614c4762a1bSJed Brown 
615d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
616d71ae5a4SJacob Faibussowitsch {
617c4762a1bSJed Brown   PetscScalar y   = j * s->hy + s->hy / 2;
618c4762a1bSJed Brown   PetscScalar awb = s->hy / (s->hx / 2);
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   PetscFunctionBeginUser;
621c4762a1bSJed Brown   if (i == 0) { /* west boundary */
622c4762a1bSJed Brown     *val = awb * StokesExactVelocityX(y);
623c4762a1bSJed Brown   } else {
624c4762a1bSJed Brown     *val = 0.0;
625c4762a1bSJed Brown   }
626*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
627c4762a1bSJed Brown }
628c4762a1bSJed Brown 
629d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
630d71ae5a4SJacob Faibussowitsch {
631c4762a1bSJed Brown   PetscFunctionBeginUser;
632c4762a1bSJed Brown   *val = 0.0;
633*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
634c4762a1bSJed Brown }
635c4762a1bSJed Brown 
636d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
637d71ae5a4SJacob Faibussowitsch {
638c4762a1bSJed Brown   PetscScalar y   = j * s->hy + s->hy / 2;
639c4762a1bSJed Brown   PetscScalar aeb = s->hy;
640c4762a1bSJed Brown 
641c4762a1bSJed Brown   PetscFunctionBeginUser;
642c4762a1bSJed Brown   if (i == 0) { /* west boundary */
643c4762a1bSJed Brown     *val = aeb * StokesExactVelocityX(y);
644c4762a1bSJed Brown   } else {
645c4762a1bSJed Brown     *val = 0.0;
646c4762a1bSJed Brown   }
647*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
648c4762a1bSJed Brown }
649c4762a1bSJed Brown 
650d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesCalcResidual(Stokes *s)
651d71ae5a4SJacob Faibussowitsch {
652c4762a1bSJed Brown   PetscReal val;
653c4762a1bSJed Brown   Vec       b0, b1;
654c4762a1bSJed Brown 
655c4762a1bSJed Brown   PetscFunctionBeginUser;
656c4762a1bSJed Brown   /* residual Ax-b (*warning* overwrites b) */
6579566063dSJacob Faibussowitsch   PetscCall(VecScale(s->b, -1.0));
6589566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(s->A, s->x, s->b, s->b));
659c4762a1bSJed Brown 
660c4762a1bSJed Brown   /* residual velocity */
6619566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[0], &b0));
6629566063dSJacob Faibussowitsch   PetscCall(VecNorm(b0, NORM_2, &val));
6639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual u = %g\n", (double)val));
6649566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0));
665c4762a1bSJed Brown 
666c4762a1bSJed Brown   /* residual pressure */
6679566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[1], &b1));
6689566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &val));
6699566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual p = %g\n", (double)val));
6709566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
671c4762a1bSJed Brown 
672c4762a1bSJed Brown   /* total residual */
6739566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->b, NORM_2, &val));
6749566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual [u,p] = %g\n", (double)val));
675*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676c4762a1bSJed Brown }
677c4762a1bSJed Brown 
678d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesCalcError(Stokes *s)
679d71ae5a4SJacob Faibussowitsch {
680c4762a1bSJed Brown   PetscScalar scale = PetscSqrtReal((double)s->nx * s->ny);
681c4762a1bSJed Brown   PetscReal   val;
682c4762a1bSJed Brown   Vec         y0, y1;
683c4762a1bSJed Brown 
684c4762a1bSJed Brown   PetscFunctionBeginUser;
685c4762a1bSJed Brown   /* error y-x */
6869566063dSJacob Faibussowitsch   PetscCall(VecAXPY(s->y, -1.0, s->x));
687c4762a1bSJed Brown 
688c4762a1bSJed Brown   /* error in velocity */
6899566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[0], &y0));
6909566063dSJacob Faibussowitsch   PetscCall(VecNorm(y0, NORM_2, &val));
6919566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error u = %g\n", (double)(PetscRealPart(val / scale))));
6929566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0));
693c4762a1bSJed Brown 
694c4762a1bSJed Brown   /* error in pressure */
6959566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[1], &y1));
6969566063dSJacob Faibussowitsch   PetscCall(VecNorm(y1, NORM_2, &val));
6979566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error p = %g\n", (double)(PetscRealPart(val / scale))));
6989566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1));
699c4762a1bSJed Brown 
700c4762a1bSJed Brown   /* total error */
7019566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->y, NORM_2, &val));
7029566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error [u,p] = %g\n", (double)PetscRealPart((val / scale))));
703*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
704c4762a1bSJed Brown }
705c4762a1bSJed Brown 
706d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
707d71ae5a4SJacob Faibussowitsch {
708c4762a1bSJed Brown   Stokes s;
709c4762a1bSJed Brown   KSP    ksp;
710c4762a1bSJed Brown 
711327415f7SBarry Smith   PetscFunctionBeginUser;
7129566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
713c4762a1bSJed Brown   s.nx = 4;
714c4762a1bSJed Brown   s.ny = 6;
7159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nx", &s.nx, NULL));
7169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ny", &s.ny, NULL));
717c4762a1bSJed Brown   s.hx           = 2.0 / s.nx;
718c4762a1bSJed Brown   s.hy           = 1.0 / s.ny;
719c4762a1bSJed Brown   s.matsymmetric = PETSC_FALSE;
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mat_set_symmetric", &s.matsymmetric, NULL));
721c4762a1bSJed Brown   s.userPC = s.userKSP = PETSC_FALSE;
7229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-user_pc", &s.userPC));
7239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-user_ksp", &s.userKSP));
724c4762a1bSJed Brown 
7259566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatrix(&s));
7269566063dSJacob Faibussowitsch   PetscCall(StokesSetupIndexSets(&s));
7279566063dSJacob Faibussowitsch   PetscCall(StokesSetupVectors(&s));
728c4762a1bSJed Brown 
7299566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
7309566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, s.A, s.A));
7319566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
7329566063dSJacob Faibussowitsch   PetscCall(StokesSetupPC(&s, ksp));
7339566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, s.b, s.x));
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   /* don't trust, verify! */
7369566063dSJacob Faibussowitsch   PetscCall(StokesCalcResidual(&s));
7379566063dSJacob Faibussowitsch   PetscCall(StokesCalcError(&s));
7389566063dSJacob Faibussowitsch   PetscCall(StokesWriteSolution(&s));
739c4762a1bSJed Brown 
7409566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
7419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[0]));
7429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[1]));
7439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[2]));
7449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[3]));
7459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.A));
7469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.x));
7479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.b));
7489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.y));
7499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.myS));
7509566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
751b122ec5aSJacob Faibussowitsch   return 0;
752c4762a1bSJed Brown }
753c4762a1bSJed Brown 
754c4762a1bSJed Brown /*TEST
755c4762a1bSJed Brown 
756c4762a1bSJed Brown    test:
757c4762a1bSJed Brown       nsize: 2
758c4762a1bSJed 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
759c4762a1bSJed Brown 
760c4762a1bSJed Brown    test:
761c4762a1bSJed Brown       suffix: 2
762c4762a1bSJed Brown       nsize: 2
763c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
764c4762a1bSJed Brown 
765c4762a1bSJed Brown    test:
766c4762a1bSJed Brown       suffix: 3
767c4762a1bSJed Brown       nsize: 2
768c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
769c4762a1bSJed Brown 
770c4762a1bSJed Brown    test:
771c4762a1bSJed Brown       suffix: 4
772c4762a1bSJed Brown       nsize: 2
773c4762a1bSJed 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
774c4762a1bSJed Brown 
775c4762a1bSJed Brown    test:
776c4762a1bSJed Brown       suffix: 4_pcksp
777c4762a1bSJed Brown       nsize: 2
778c4762a1bSJed 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
779c4762a1bSJed Brown 
780c4762a1bSJed Brown    test:
781c4762a1bSJed Brown       suffix: 5
782c4762a1bSJed Brown       nsize: 2
783c4762a1bSJed 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
784c4762a1bSJed Brown 
785c4762a1bSJed Brown    test:
786c4762a1bSJed Brown       suffix: 6
787c4762a1bSJed Brown       nsize: 2
788c4762a1bSJed 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
789c4762a1bSJed Brown 
790c4762a1bSJed Brown    test:
791c4762a1bSJed Brown       suffix: 7
792c4762a1bSJed Brown       nsize: 2
793c4762a1bSJed 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
794c4762a1bSJed Brown 
795c4762a1bSJed Brown TEST*/
796