xref: /petsc/src/snes/tutorials/ex70.c (revision 32e037510874440daa9fa8248b5cfeac7662c93e)
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 */
76dd8e379bSPierre Jolivet PetscErrorCode StokesRhsMomX(Stokes *, PetscInt, PetscInt, PetscScalar *); /* right-hand side of velocity (x-component) */
77dd8e379bSPierre Jolivet PetscErrorCode StokesRhsMomY(Stokes *, PetscInt, PetscInt, PetscScalar *); /* right-hand side of velocity (y-component) */
78dd8e379bSPierre Jolivet 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) */
StokesExactVelocityX(const PetscScalar y)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 */
StokesExactPressure(const PetscScalar x)92d71ae5a4SJacob Faibussowitsch PetscScalar StokesExactPressure(const PetscScalar x)
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   return 8.0 * (2.0 - x);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
StokesSetupPC(Stokes * s,KSP ksp)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   }
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
StokesWriteSolution(Stokes * s)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   }
1403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
StokesSetupIndexSets(Stokes * s)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));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
StokesSetupVectors(Stokes * s)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));
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
StokesGetPosition(Stokes * s,PetscInt row,PetscInt * i,PetscInt * j)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;
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
StokesExactSolution(Stokes * s)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   }
200eadb0ebdSBarry Smith   PetscCall(VecAssemblyBegin(y0));
201eadb0ebdSBarry Smith   PetscCall(VecAssemblyEnd(y0));
2029566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* pressure part */
2059566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[1], &y1));
2069566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(y1, &start, &end));
207c4762a1bSJed Brown   for (row = start; row < end; row++) {
2089566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
209c4762a1bSJed Brown     val = StokesExactPressure(i * s->hx + s->hx / 2);
2109566063dSJacob Faibussowitsch     PetscCall(VecSetValue(y1, row, val, INSERT_VALUES));
211c4762a1bSJed Brown   }
212eadb0ebdSBarry Smith   PetscCall(VecAssemblyBegin(y1));
213eadb0ebdSBarry Smith   PetscCall(VecAssemblyEnd(y1));
2149566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1));
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
StokesRhs(Stokes * s)218d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhs(Stokes *s)
219d71ae5a4SJacob Faibussowitsch {
220c4762a1bSJed Brown   PetscInt    row, start, end, i, j;
221c4762a1bSJed Brown   PetscScalar val;
222c4762a1bSJed Brown   Vec         b0, b1;
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   PetscFunctionBeginUser;
225c4762a1bSJed Brown   /* velocity part */
2269566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[0], &b0));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(b0, &start, &end));
228c4762a1bSJed Brown   for (row = start; row < end; row++) {
2299566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
230c4762a1bSJed Brown     if (row < s->nx * s->ny) {
2319566063dSJacob Faibussowitsch       PetscCall(StokesRhsMomX(s, i, j, &val));
232c4762a1bSJed Brown     } else {
2339566063dSJacob Faibussowitsch       PetscCall(StokesRhsMomY(s, i, j, &val));
234c4762a1bSJed Brown     }
2359566063dSJacob Faibussowitsch     PetscCall(VecSetValue(b0, row, val, INSERT_VALUES));
236c4762a1bSJed Brown   }
2379566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /* pressure part */
2409566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[1], &b1));
2419566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(b1, &start, &end));
242c4762a1bSJed Brown   for (row = start; row < end; row++) {
2439566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
2449566063dSJacob Faibussowitsch     PetscCall(StokesRhsMass(s, i, j, &val));
245ad540459SPierre Jolivet     if (s->matsymmetric) val = -1.0 * val;
2469566063dSJacob Faibussowitsch     PetscCall(VecSetValue(b1, row, val, INSERT_VALUES));
247c4762a1bSJed Brown   }
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
2493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
250c4762a1bSJed Brown }
251c4762a1bSJed Brown 
StokesSetupMatBlock00(Stokes * s)252d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock00(Stokes *s)
253d71ae5a4SJacob Faibussowitsch {
254c4762a1bSJed Brown   PetscInt    row, start, end, sz, i, j;
255c4762a1bSJed Brown   PetscInt    cols[5];
256c4762a1bSJed Brown   PetscScalar vals[5];
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBeginUser;
259c4762a1bSJed Brown   /* A[0] is 2N-by-2N */
2609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[0]));
2619566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[0], "a00_"));
2629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[0], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, 2 * s->nx * s->ny));
2639566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[0], MATMPIAIJ));
2649566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[0], 5, NULL, 5, NULL));
2659566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[0], &start, &end));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   for (row = start; row < end; row++) {
2689566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
269c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
2709566063dSJacob Faibussowitsch     PetscCall(StokesStencilLaplacian(s, i, j, &sz, cols, vals));
271c4762a1bSJed Brown     /* second part: rows (nx*ny) to (2*nx*ny-1) */
272c4762a1bSJed Brown     if (row >= s->nx * s->ny) {
273c4762a1bSJed Brown       for (i = 0; i < sz; i++) cols[i] += s->nx * s->ny;
274c4762a1bSJed Brown     }
275c4762a1bSJed Brown     for (i = 0; i < sz; i++) vals[i] = -1.0 * vals[i]; /* dynamic viscosity coef mu=-1 */
2769566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES));
277c4762a1bSJed Brown   }
2789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY));
2799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY));
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
281c4762a1bSJed Brown }
282c4762a1bSJed Brown 
StokesSetupMatBlock01(Stokes * s)283d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock01(Stokes *s)
284d71ae5a4SJacob Faibussowitsch {
285c4762a1bSJed Brown   PetscInt    row, start, end, sz, i, j;
286c4762a1bSJed Brown   PetscInt    cols[5];
287c4762a1bSJed Brown   PetscScalar vals[5];
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   PetscFunctionBeginUser;
290c4762a1bSJed Brown   /* A[1] is 2N-by-N */
2919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[1]));
2929566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[1], "a01_"));
2939566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(s->subA[1], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, s->nx * s->ny));
2949566063dSJacob Faibussowitsch   PetscCall(MatSetType(s->subA[1], MATMPIAIJ));
2959566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(s->subA[1], 5, NULL, 5, NULL));
2969566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(s->subA[1], &start, &end));
297c4762a1bSJed Brown 
2989566063dSJacob Faibussowitsch   PetscCall(MatSetOption(s->subA[1], MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   for (row = start; row < end; row++) {
3019566063dSJacob Faibussowitsch     PetscCall(StokesGetPosition(s, row, &i, &j));
302c4762a1bSJed Brown     /* first part: rows 0 to (nx*ny-1) */
303c4762a1bSJed Brown     if (row < s->nx * s->ny) {
3049566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientX(s, i, j, &sz, cols, vals));
305c4762a1bSJed Brown     } else { /* second part: rows (nx*ny) to (2*nx*ny-1) */
3069566063dSJacob Faibussowitsch       PetscCall(StokesStencilGradientY(s, i, j, &sz, cols, vals));
307c4762a1bSJed Brown     }
3089566063dSJacob Faibussowitsch     PetscCall(MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES));
309c4762a1bSJed Brown   }
3109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY));
3119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
StokesSetupMatBlock10(Stokes * s)315d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock10(Stokes *s)
316d71ae5a4SJacob Faibussowitsch {
317c4762a1bSJed Brown   PetscFunctionBeginUser;
318c4762a1bSJed Brown   /* A[2] is minus transpose of A[1] */
3199566063dSJacob Faibussowitsch   PetscCall(MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]));
32048a46eb9SPierre Jolivet   if (!s->matsymmetric) PetscCall(MatScale(s->subA[2], -1.0));
3219566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(s->subA[2], "a10_"));
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
StokesSetupMatBlock11(Stokes * s)325d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatBlock11(Stokes *s)
326d71ae5a4SJacob Faibussowitsch {
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));
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
StokesSetupApproxSchur(Stokes * s)339d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupApproxSchur(Stokes *s)
340d71ae5a4SJacob Faibussowitsch {
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]) */
358*fb842aefSJose E. Roman   PetscCall(MatMatMult(s->subA[2], s->subA[1], MAT_INITIAL_MATRIX, PETSC_DETERMINE, &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));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366c4762a1bSJed Brown }
367c4762a1bSJed Brown 
StokesSetupMatrix(Stokes * s)368d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesSetupMatrix(Stokes *s)
369d71ae5a4SJacob Faibussowitsch {
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));
3773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
378c4762a1bSJed Brown }
379c4762a1bSJed Brown 
StokesStencilLaplacian(Stokes * s,PetscInt i,PetscInt j,PetscInt * sz,PetscInt * cols,PetscScalar * vals)380d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
381d71ae5a4SJacob Faibussowitsch {
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;
3919371c9d4SSatish Balay     cols[0] = p;
3929371c9d4SSatish Balay     vals[0] = -(ae + awb + asb + an);
3939371c9d4SSatish Balay     cols[1] = e;
3949371c9d4SSatish Balay     vals[1] = ae;
3959371c9d4SSatish Balay     cols[2] = n;
3969371c9d4SSatish Balay     vals[2] = an;
397c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
398c4762a1bSJed Brown     *sz     = 3;
3999371c9d4SSatish Balay     cols[0] = s2;
4009371c9d4SSatish Balay     vals[0] = as;
4019371c9d4SSatish Balay     cols[1] = p;
4029371c9d4SSatish Balay     vals[1] = -(ae + awb + as + anb);
4039371c9d4SSatish Balay     cols[2] = e;
4049371c9d4SSatish Balay     vals[2] = ae;
405c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
406c4762a1bSJed Brown     *sz     = 3;
4079371c9d4SSatish Balay     cols[0] = w;
4089371c9d4SSatish Balay     vals[0] = aw;
4099371c9d4SSatish Balay     cols[1] = p;
4109371c9d4SSatish Balay     vals[1] = -(aeb + aw + asb + an);
4119371c9d4SSatish Balay     cols[2] = n;
4129371c9d4SSatish Balay     vals[2] = an;
413c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
414c4762a1bSJed Brown     *sz     = 3;
4159371c9d4SSatish Balay     cols[0] = s2;
4169371c9d4SSatish Balay     vals[0] = as;
4179371c9d4SSatish Balay     cols[1] = w;
4189371c9d4SSatish Balay     vals[1] = aw;
4199371c9d4SSatish Balay     cols[2] = p;
4209371c9d4SSatish Balay     vals[2] = -(aeb + aw + as + anb);
421c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
422c4762a1bSJed Brown     *sz     = 4;
4239371c9d4SSatish Balay     cols[0] = s2;
4249371c9d4SSatish Balay     vals[0] = as;
4259371c9d4SSatish Balay     cols[1] = p;
4269371c9d4SSatish Balay     vals[1] = -(ae + awb + as + an);
4279371c9d4SSatish Balay     cols[2] = e;
4289371c9d4SSatish Balay     vals[2] = ae;
4299371c9d4SSatish Balay     cols[3] = n;
4309371c9d4SSatish Balay     vals[3] = an;
431c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
432c4762a1bSJed Brown     *sz     = 4;
4339371c9d4SSatish Balay     cols[0] = s2;
4349371c9d4SSatish Balay     vals[0] = as;
4359371c9d4SSatish Balay     cols[1] = w;
4369371c9d4SSatish Balay     vals[1] = aw;
4379371c9d4SSatish Balay     cols[2] = p;
4389371c9d4SSatish Balay     vals[2] = -(aeb + aw + as + an);
4399371c9d4SSatish Balay     cols[3] = n;
4409371c9d4SSatish Balay     vals[3] = an;
441c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
442c4762a1bSJed Brown     *sz     = 4;
4439371c9d4SSatish Balay     cols[0] = w;
4449371c9d4SSatish Balay     vals[0] = aw;
4459371c9d4SSatish Balay     cols[1] = p;
4469371c9d4SSatish Balay     vals[1] = -(ae + aw + asb + an);
4479371c9d4SSatish Balay     cols[2] = e;
4489371c9d4SSatish Balay     vals[2] = ae;
4499371c9d4SSatish Balay     cols[3] = n;
4509371c9d4SSatish Balay     vals[3] = an;
451c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
452c4762a1bSJed Brown     *sz     = 4;
4539371c9d4SSatish Balay     cols[0] = s2;
4549371c9d4SSatish Balay     vals[0] = as;
4559371c9d4SSatish Balay     cols[1] = w;
4569371c9d4SSatish Balay     vals[1] = aw;
4579371c9d4SSatish Balay     cols[2] = p;
4589371c9d4SSatish Balay     vals[2] = -(ae + aw + as + anb);
4599371c9d4SSatish Balay     cols[3] = e;
4609371c9d4SSatish Balay     vals[3] = ae;
461c4762a1bSJed Brown   } else { /* interior */
462c4762a1bSJed Brown     *sz     = 5;
4639371c9d4SSatish Balay     cols[0] = s2;
4649371c9d4SSatish Balay     vals[0] = as;
4659371c9d4SSatish Balay     cols[1] = w;
4669371c9d4SSatish Balay     vals[1] = aw;
4679371c9d4SSatish Balay     cols[2] = p;
4689371c9d4SSatish Balay     vals[2] = -(ae + aw + as + an);
4699371c9d4SSatish Balay     cols[3] = e;
4709371c9d4SSatish Balay     vals[3] = ae;
4719371c9d4SSatish Balay     cols[4] = n;
4729371c9d4SSatish Balay     vals[4] = an;
473c4762a1bSJed Brown   }
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
475c4762a1bSJed Brown }
476c4762a1bSJed Brown 
StokesStencilGradientX(Stokes * s,PetscInt i,PetscInt j,PetscInt * sz,PetscInt * cols,PetscScalar * vals)477d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
478d71ae5a4SJacob Faibussowitsch {
479c4762a1bSJed Brown   PetscInt    p = j * s->nx + i, w = p - 1, e = p + 1;
480c4762a1bSJed Brown   PetscScalar ae = s->hy / 2, aeb = s->hy;
481c4762a1bSJed Brown   PetscScalar aw = -s->hy / 2, awb = 0;
482c4762a1bSJed Brown 
483c4762a1bSJed Brown   PetscFunctionBeginUser;
484c4762a1bSJed Brown   if (i == 0 && j == 0) { /* south-west corner */
485c4762a1bSJed Brown     *sz     = 2;
4869371c9d4SSatish Balay     cols[0] = p;
4879371c9d4SSatish Balay     vals[0] = -(ae + awb);
4889371c9d4SSatish Balay     cols[1] = e;
4899371c9d4SSatish Balay     vals[1] = ae;
490c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
491c4762a1bSJed Brown     *sz     = 2;
4929371c9d4SSatish Balay     cols[0] = p;
4939371c9d4SSatish Balay     vals[0] = -(ae + awb);
4949371c9d4SSatish Balay     cols[1] = e;
4959371c9d4SSatish Balay     vals[1] = ae;
496c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
497c4762a1bSJed Brown     *sz     = 2;
4989371c9d4SSatish Balay     cols[0] = w;
4999371c9d4SSatish Balay     vals[0] = aw;
5009371c9d4SSatish Balay     cols[1] = p;
5019371c9d4SSatish Balay     vals[1] = -(aeb + aw);
502c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
503c4762a1bSJed Brown     *sz     = 2;
5049371c9d4SSatish Balay     cols[0] = w;
5059371c9d4SSatish Balay     vals[0] = aw;
5069371c9d4SSatish Balay     cols[1] = p;
5079371c9d4SSatish Balay     vals[1] = -(aeb + aw);
508c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
509c4762a1bSJed Brown     *sz     = 2;
5109371c9d4SSatish Balay     cols[0] = p;
5119371c9d4SSatish Balay     vals[0] = -(ae + awb);
5129371c9d4SSatish Balay     cols[1] = e;
5139371c9d4SSatish Balay     vals[1] = ae;
514c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
515c4762a1bSJed Brown     *sz     = 2;
5169371c9d4SSatish Balay     cols[0] = w;
5179371c9d4SSatish Balay     vals[0] = aw;
5189371c9d4SSatish Balay     cols[1] = p;
5199371c9d4SSatish Balay     vals[1] = -(aeb + aw);
520c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
521c4762a1bSJed Brown     *sz     = 3;
5229371c9d4SSatish Balay     cols[0] = w;
5239371c9d4SSatish Balay     vals[0] = aw;
5249371c9d4SSatish Balay     cols[1] = p;
5259371c9d4SSatish Balay     vals[1] = -(ae + aw);
5269371c9d4SSatish Balay     cols[2] = e;
5279371c9d4SSatish Balay     vals[2] = ae;
528c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
529c4762a1bSJed Brown     *sz     = 3;
5309371c9d4SSatish Balay     cols[0] = w;
5319371c9d4SSatish Balay     vals[0] = aw;
5329371c9d4SSatish Balay     cols[1] = p;
5339371c9d4SSatish Balay     vals[1] = -(ae + aw);
5349371c9d4SSatish Balay     cols[2] = e;
5359371c9d4SSatish Balay     vals[2] = ae;
536c4762a1bSJed Brown   } else { /* interior */
537c4762a1bSJed Brown     *sz     = 3;
5389371c9d4SSatish Balay     cols[0] = w;
5399371c9d4SSatish Balay     vals[0] = aw;
5409371c9d4SSatish Balay     cols[1] = p;
5419371c9d4SSatish Balay     vals[1] = -(ae + aw);
5429371c9d4SSatish Balay     cols[2] = e;
5439371c9d4SSatish Balay     vals[2] = ae;
544c4762a1bSJed Brown   }
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546c4762a1bSJed Brown }
547c4762a1bSJed Brown 
StokesStencilGradientY(Stokes * s,PetscInt i,PetscInt j,PetscInt * sz,PetscInt * cols,PetscScalar * vals)548d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
549d71ae5a4SJacob Faibussowitsch {
550c4762a1bSJed Brown   PetscInt    p = j * s->nx + i, s2 = p - s->nx, n = p + s->nx;
551c4762a1bSJed Brown   PetscScalar as = -s->hx / 2, asb = 0;
552c4762a1bSJed Brown   PetscScalar an = s->hx / 2, anb = 0;
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   PetscFunctionBeginUser;
555c4762a1bSJed Brown   if (i == 0 && j == 0) { /* south-west corner */
556c4762a1bSJed Brown     *sz     = 2;
5579371c9d4SSatish Balay     cols[0] = p;
5589371c9d4SSatish Balay     vals[0] = -(asb + an);
5599371c9d4SSatish Balay     cols[1] = n;
5609371c9d4SSatish Balay     vals[1] = an;
561c4762a1bSJed Brown   } else if (i == 0 && j == s->ny - 1) { /* north-west corner */
562c4762a1bSJed Brown     *sz     = 2;
5639371c9d4SSatish Balay     cols[0] = s2;
5649371c9d4SSatish Balay     vals[0] = as;
5659371c9d4SSatish Balay     cols[1] = p;
5669371c9d4SSatish Balay     vals[1] = -(as + anb);
567c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == 0) { /* south-east corner */
568c4762a1bSJed Brown     *sz     = 2;
5699371c9d4SSatish Balay     cols[0] = p;
5709371c9d4SSatish Balay     vals[0] = -(asb + an);
5719371c9d4SSatish Balay     cols[1] = n;
5729371c9d4SSatish Balay     vals[1] = an;
573c4762a1bSJed Brown   } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */
574c4762a1bSJed Brown     *sz     = 2;
5759371c9d4SSatish Balay     cols[0] = s2;
5769371c9d4SSatish Balay     vals[0] = as;
5779371c9d4SSatish Balay     cols[1] = p;
5789371c9d4SSatish Balay     vals[1] = -(as + anb);
579c4762a1bSJed Brown   } else if (i == 0) { /* west boundary */
580c4762a1bSJed Brown     *sz     = 3;
5819371c9d4SSatish Balay     cols[0] = s2;
5829371c9d4SSatish Balay     vals[0] = as;
5839371c9d4SSatish Balay     cols[1] = p;
5849371c9d4SSatish Balay     vals[1] = -(as + an);
5859371c9d4SSatish Balay     cols[2] = n;
5869371c9d4SSatish Balay     vals[2] = an;
587c4762a1bSJed Brown   } else if (i == s->nx - 1) { /* east boundary */
588c4762a1bSJed Brown     *sz     = 3;
5899371c9d4SSatish Balay     cols[0] = s2;
5909371c9d4SSatish Balay     vals[0] = as;
5919371c9d4SSatish Balay     cols[1] = p;
5929371c9d4SSatish Balay     vals[1] = -(as + an);
5939371c9d4SSatish Balay     cols[2] = n;
5949371c9d4SSatish Balay     vals[2] = an;
595c4762a1bSJed Brown   } else if (j == 0) { /* south boundary */
596c4762a1bSJed Brown     *sz     = 2;
5979371c9d4SSatish Balay     cols[0] = p;
5989371c9d4SSatish Balay     vals[0] = -(asb + an);
5999371c9d4SSatish Balay     cols[1] = n;
6009371c9d4SSatish Balay     vals[1] = an;
601c4762a1bSJed Brown   } else if (j == s->ny - 1) { /* north boundary */
602c4762a1bSJed Brown     *sz     = 2;
6039371c9d4SSatish Balay     cols[0] = s2;
6049371c9d4SSatish Balay     vals[0] = as;
6059371c9d4SSatish Balay     cols[1] = p;
6069371c9d4SSatish Balay     vals[1] = -(as + anb);
607c4762a1bSJed Brown   } else { /* interior */
608c4762a1bSJed Brown     *sz     = 3;
6099371c9d4SSatish Balay     cols[0] = s2;
6109371c9d4SSatish Balay     vals[0] = as;
6119371c9d4SSatish Balay     cols[1] = p;
6129371c9d4SSatish Balay     vals[1] = -(as + an);
6139371c9d4SSatish Balay     cols[2] = n;
6149371c9d4SSatish Balay     vals[2] = an;
615c4762a1bSJed Brown   }
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
StokesRhsMomX(Stokes * s,PetscInt i,PetscInt j,PetscScalar * val)619d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
620d71ae5a4SJacob Faibussowitsch {
621c4762a1bSJed Brown   PetscScalar y   = j * s->hy + s->hy / 2;
622c4762a1bSJed Brown   PetscScalar awb = s->hy / (s->hx / 2);
623c4762a1bSJed Brown 
624c4762a1bSJed Brown   PetscFunctionBeginUser;
625c4762a1bSJed Brown   if (i == 0) { /* west boundary */
626c4762a1bSJed Brown     *val = awb * StokesExactVelocityX(y);
627c4762a1bSJed Brown   } else {
628c4762a1bSJed Brown     *val = 0.0;
629c4762a1bSJed Brown   }
6303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
631c4762a1bSJed Brown }
632c4762a1bSJed Brown 
StokesRhsMomY(Stokes * s,PetscInt i,PetscInt j,PetscScalar * val)633d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
634d71ae5a4SJacob Faibussowitsch {
635c4762a1bSJed Brown   PetscFunctionBeginUser;
636c4762a1bSJed Brown   *val = 0.0;
6373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
638c4762a1bSJed Brown }
639c4762a1bSJed Brown 
StokesRhsMass(Stokes * s,PetscInt i,PetscInt j,PetscScalar * val)640d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
641d71ae5a4SJacob Faibussowitsch {
642c4762a1bSJed Brown   PetscScalar y   = j * s->hy + s->hy / 2;
643c4762a1bSJed Brown   PetscScalar aeb = s->hy;
644c4762a1bSJed Brown 
645c4762a1bSJed Brown   PetscFunctionBeginUser;
646c4762a1bSJed Brown   if (i == 0) { /* west boundary */
647c4762a1bSJed Brown     *val = aeb * StokesExactVelocityX(y);
648c4762a1bSJed Brown   } else {
649c4762a1bSJed Brown     *val = 0.0;
650c4762a1bSJed Brown   }
6513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
652c4762a1bSJed Brown }
653c4762a1bSJed Brown 
StokesCalcResidual(Stokes * s)654d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesCalcResidual(Stokes *s)
655d71ae5a4SJacob Faibussowitsch {
656c4762a1bSJed Brown   PetscReal val;
657c4762a1bSJed Brown   Vec       b0, b1;
658c4762a1bSJed Brown 
659c4762a1bSJed Brown   PetscFunctionBeginUser;
660c4762a1bSJed Brown   /* residual Ax-b (*warning* overwrites b) */
6619566063dSJacob Faibussowitsch   PetscCall(VecScale(s->b, -1.0));
6629566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(s->A, s->x, s->b, s->b));
663c4762a1bSJed Brown 
664c4762a1bSJed Brown   /* residual velocity */
6659566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[0], &b0));
6669566063dSJacob Faibussowitsch   PetscCall(VecNorm(b0, NORM_2, &val));
6679566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual u = %g\n", (double)val));
6689566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0));
669c4762a1bSJed Brown 
670c4762a1bSJed Brown   /* residual pressure */
6719566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->b, s->isg[1], &b1));
6729566063dSJacob Faibussowitsch   PetscCall(VecNorm(b1, NORM_2, &val));
6739566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual p = %g\n", (double)val));
6749566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1));
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   /* total residual */
6779566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->b, NORM_2, &val));
6789566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual [u,p] = %g\n", (double)val));
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680c4762a1bSJed Brown }
681c4762a1bSJed Brown 
StokesCalcError(Stokes * s)682d71ae5a4SJacob Faibussowitsch PetscErrorCode StokesCalcError(Stokes *s)
683d71ae5a4SJacob Faibussowitsch {
684c4762a1bSJed Brown   PetscScalar scale = PetscSqrtReal((double)s->nx * s->ny);
685c4762a1bSJed Brown   PetscReal   val;
686c4762a1bSJed Brown   Vec         y0, y1;
687c4762a1bSJed Brown 
688c4762a1bSJed Brown   PetscFunctionBeginUser;
689c4762a1bSJed Brown   /* error y-x */
6909566063dSJacob Faibussowitsch   PetscCall(VecAXPY(s->y, -1.0, s->x));
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   /* error in velocity */
6939566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[0], &y0));
6949566063dSJacob Faibussowitsch   PetscCall(VecNorm(y0, NORM_2, &val));
6959566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error u = %g\n", (double)(PetscRealPart(val / scale))));
6969566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0));
697c4762a1bSJed Brown 
698c4762a1bSJed Brown   /* error in pressure */
6999566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(s->y, s->isg[1], &y1));
7009566063dSJacob Faibussowitsch   PetscCall(VecNorm(y1, NORM_2, &val));
7019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error p = %g\n", (double)(PetscRealPart(val / scale))));
7029566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1));
703c4762a1bSJed Brown 
704c4762a1bSJed Brown   /* total error */
7059566063dSJacob Faibussowitsch   PetscCall(VecNorm(s->y, NORM_2, &val));
706f4f49eeaSPierre Jolivet   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error [u,p] = %g\n", (double)PetscRealPart(val / scale)));
7073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
708c4762a1bSJed Brown }
709c4762a1bSJed Brown 
main(int argc,char ** argv)710d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
711d71ae5a4SJacob Faibussowitsch {
712c4762a1bSJed Brown   Stokes s;
713c4762a1bSJed Brown   KSP    ksp;
714c4762a1bSJed Brown 
715327415f7SBarry Smith   PetscFunctionBeginUser;
7169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
717c4762a1bSJed Brown   s.nx = 4;
718c4762a1bSJed Brown   s.ny = 6;
7199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nx", &s.nx, NULL));
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ny", &s.ny, NULL));
721c4762a1bSJed Brown   s.hx           = 2.0 / s.nx;
722c4762a1bSJed Brown   s.hy           = 1.0 / s.ny;
723c4762a1bSJed Brown   s.matsymmetric = PETSC_FALSE;
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mat_set_symmetric", &s.matsymmetric, NULL));
725c4762a1bSJed Brown   s.userPC = s.userKSP = PETSC_FALSE;
7269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-user_pc", &s.userPC));
7279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-user_ksp", &s.userKSP));
728c4762a1bSJed Brown 
7299566063dSJacob Faibussowitsch   PetscCall(StokesSetupMatrix(&s));
7309566063dSJacob Faibussowitsch   PetscCall(StokesSetupIndexSets(&s));
7319566063dSJacob Faibussowitsch   PetscCall(StokesSetupVectors(&s));
732c4762a1bSJed Brown 
7339566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
7349566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, s.A, s.A));
7359566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
7369566063dSJacob Faibussowitsch   PetscCall(StokesSetupPC(&s, ksp));
7379566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, s.b, s.x));
738c4762a1bSJed Brown 
739c4762a1bSJed Brown   /* don't trust, verify! */
7409566063dSJacob Faibussowitsch   PetscCall(StokesCalcResidual(&s));
7419566063dSJacob Faibussowitsch   PetscCall(StokesCalcError(&s));
7429566063dSJacob Faibussowitsch   PetscCall(StokesWriteSolution(&s));
743c4762a1bSJed Brown 
7449566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
7459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[0]));
7469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[1]));
7479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[2]));
7489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.subA[3]));
7499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.A));
7509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.x));
7519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.b));
7529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s.y));
7539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&s.myS));
7549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
755b122ec5aSJacob Faibussowitsch   return 0;
756c4762a1bSJed Brown }
757c4762a1bSJed Brown 
758c4762a1bSJed Brown /*TEST
759c4762a1bSJed Brown 
760c4762a1bSJed Brown    test:
761c4762a1bSJed Brown       nsize: 2
762c4762a1bSJed 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
763c4762a1bSJed Brown 
764c4762a1bSJed Brown    test:
765c4762a1bSJed Brown       suffix: 2
766c4762a1bSJed Brown       nsize: 2
767c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
768c4762a1bSJed Brown 
769c4762a1bSJed Brown    test:
770c4762a1bSJed Brown       suffix: 3
771c4762a1bSJed Brown       nsize: 2
772c4762a1bSJed Brown       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc
773c4762a1bSJed Brown 
774c4762a1bSJed Brown    test:
775c4762a1bSJed Brown       suffix: 4
776c4762a1bSJed Brown       nsize: 2
777c4762a1bSJed 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
778c4762a1bSJed Brown 
779c4762a1bSJed Brown    test:
780c4762a1bSJed Brown       suffix: 4_pcksp
781c4762a1bSJed Brown       nsize: 2
782c4762a1bSJed 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
783c4762a1bSJed Brown 
784c4762a1bSJed Brown    test:
785c4762a1bSJed Brown       suffix: 5
786c4762a1bSJed Brown       nsize: 2
787c4762a1bSJed 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
788c4762a1bSJed Brown 
789c4762a1bSJed Brown    test:
790c4762a1bSJed Brown       suffix: 6
791c4762a1bSJed Brown       nsize: 2
792c4762a1bSJed 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
793c4762a1bSJed Brown 
794c4762a1bSJed Brown    test:
795c4762a1bSJed Brown       suffix: 7
796c4762a1bSJed Brown       nsize: 2
797c4762a1bSJed 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
798c4762a1bSJed Brown 
799c4762a1bSJed Brown TEST*/
800