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