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