xref: /petsc/doc/tutorials/physics/guide_to_stokes.md (revision aa9a5b67dde7169d6e428ea38b8433e2afcb1019)
1*aa9a5b67SBarry Smith---
2*aa9a5b67SBarry Smithorphan: true
3*aa9a5b67SBarry Smith---
4*aa9a5b67SBarry Smith
5*aa9a5b67SBarry Smith(tut_stokes)=
6*aa9a5b67SBarry Smith
7*aa9a5b67SBarry Smith# Guide to the Stokes Equations using Finite Elements
8*aa9a5b67SBarry Smith
9*aa9a5b67SBarry SmithThe Stokes equations
10*aa9a5b67SBarry Smith
11*aa9a5b67SBarry Smith$$
12*aa9a5b67SBarry Smith\begin{aligned}
13*aa9a5b67SBarry Smith-\nabla \cdot \left(\mu \left(\nabla u + \nabla u^T \right)\right) + \nabla p + f &= 0 \\
14*aa9a5b67SBarry Smith\nabla\cdot u &= 0 \end{aligned}
15*aa9a5b67SBarry Smith$$
16*aa9a5b67SBarry Smith
17*aa9a5b67SBarry Smithdescribe slow flow of an incompressible fluid with velocity $u$, pressure $p$, and body force $f$.
18*aa9a5b67SBarry Smith
19*aa9a5b67SBarry SmithThis guide accompanies <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex62.c.html">SNES Example 62></a> and <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex69.c.html">SNES Example 69><\a>.
20*aa9a5b67SBarry Smith
21*aa9a5b67SBarry SmithThe Stokes equations for a fluid, a steady-state form of the Navier-Stokes equations, start with the balance of momentum, just as in elastostatics,
22*aa9a5b67SBarry Smith
23*aa9a5b67SBarry Smith$$
24*aa9a5b67SBarry Smith\nabla \cdot \sigma + f = 0,
25*aa9a5b67SBarry Smith$$
26*aa9a5b67SBarry Smith
27*aa9a5b67SBarry Smithwhere $\sigma$ is the stress tensor and $f$ is the body force, combined with the conservation of mass
28*aa9a5b67SBarry Smith
29*aa9a5b67SBarry Smith$$
30*aa9a5b67SBarry Smith\nabla \cdot (\rho u) = 0,
31*aa9a5b67SBarry Smith$$
32*aa9a5b67SBarry Smith
33*aa9a5b67SBarry Smithwhere $\rho$ is the density and $u$ is the fluid velocity. If we assume that the density is constant, making the fluid incompressible, and that the rheology is Newtonian, meaning that the viscous stress is linearly proportional to the local strain rate, then we have
34*aa9a5b67SBarry Smith
35*aa9a5b67SBarry Smith$$
36*aa9a5b67SBarry Smith\begin{aligned}
37*aa9a5b67SBarry Smith  \nabla \cdot \mu \left( \nabla u + \nabla u^T \right) - \nabla p + f &= 0 \\
38*aa9a5b67SBarry Smith  \nabla \cdot u &= 0
39*aa9a5b67SBarry Smith\end{aligned}
40*aa9a5b67SBarry Smith$$
41*aa9a5b67SBarry Smith
42*aa9a5b67SBarry Smithwhere $p$ is the pressure, $\mu$ is the dynamic shear viscosity, with units $N\cdot s/m^2$ or $Pa\cdot s$. If we divide by the constant density, we would have the kinematic viscosity $\nu$ and a force per unit mass. The second equation demands that the velocity field be divergence-free, indicating that the flow is incompressible. The pressure in this case can be thought of as the Lagrange multiplier enforcing the incompressibility constraint. In the compressible case, we would need an equation of state to relate the pressure to the density, and perhaps temperature.
43*aa9a5b67SBarry Smith
44*aa9a5b67SBarry SmithWe will discretize our Stokes equations with finite elements, so the first step is to write a variational weak form of the equations. We choose to use a Ritz-Galerkin setup, so let our velocity $u \in V$ and pressure $p \in Q$, so that
45*aa9a5b67SBarry Smith
46*aa9a5b67SBarry Smith$$
47*aa9a5b67SBarry Smith\begin{aligned}
48*aa9a5b67SBarry Smith  \left< \nabla v, \mu \left( \nabla u + \nabla u^T \right) \right> + \left< v, \frac{\partial\sigma}{\partial n} \right>_\Gamma - \left< \nabla\cdot v, p \right> - \left< v, f \right> &= 0 & \text{for all} \ v \in V\\
49*aa9a5b67SBarry Smith  \left< q, -\nabla \cdot u \right> &= 0 & \text{for all} \ q \in Q
50*aa9a5b67SBarry Smith\end{aligned}
51*aa9a5b67SBarry Smith$$
52*aa9a5b67SBarry Smith
53*aa9a5b67SBarry Smithwhere integration by parts has added a boundary integral over the normal derivative of the stress (traction), and natural boundary conditions correspond to stress-free boundaries. We have multiplied the continuity equation by minus one in order to preserve symmetry.
54*aa9a5b67SBarry Smith
55*aa9a5b67SBarry Smith## Equation Definition
56*aa9a5b67SBarry Smith
57*aa9a5b67SBarry SmithThe test functions $v, q$ and their derivatives are determined by the discretization, whereas the form of the integrand is determined by the physics. Given a quadrature rule to evaluate the form integral, we would only need the evaluation of the physics integrand at the quadrature points, given the values of the fields and their derivatives. The entire scheme is detailed in {cite}`knepleybrownruppsmith13`. The kernels paired with test functions we will call $f_0$ and those paired with gradients of test functions will be called $f_1$.
58*aa9a5b67SBarry Smith
59*aa9a5b67SBarry SmithFor example, the kernel for the continuity equation, paired with the pressure test function, is called `f0_p` and can be seen here
60*aa9a5b67SBarry Smith
61*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
62*aa9a5b67SBarry Smith:end-at: '}'
63*aa9a5b67SBarry Smith:start-at: static void f0_p(
64*aa9a5b67SBarry Smith```
65*aa9a5b67SBarry Smith
66*aa9a5b67SBarry SmithWe use the components of the Jacobian of $u$ to build up its divergence. For the balance of momentum excluding body force, we test against the gradient of the test function, as seen in `f1_u`,
67*aa9a5b67SBarry Smith
68*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
69*aa9a5b67SBarry Smith:end-at: '}'
70*aa9a5b67SBarry Smith:start-at: static void f1_u(
71*aa9a5b67SBarry Smith```
72*aa9a5b67SBarry Smith
73*aa9a5b67SBarry SmithNotice how the pressure $p$ is referred to using `u[uOff[1]]` so that we can have many fields with different numbers of components. `DMPlex` uses these point functions to construct the residual. A similar set of point functions is also used to build the Jacobian. The last piece of our physics specification is the construction of exact solutions using the Method of Manufactured Solutions (MMS).
74*aa9a5b67SBarry Smith
75*aa9a5b67SBarry Smith## MMS Solutions
76*aa9a5b67SBarry Smith
77*aa9a5b67SBarry SmithAn MMS solution is chosen to elucidate some property of the problem, and to check that it is being solved accurately, since the error can be calculated explicitly. For our Stokes problem, we first choose a solution with quadratic velocity and linear pressure,
78*aa9a5b67SBarry Smith
79*aa9a5b67SBarry Smith$$
80*aa9a5b67SBarry Smithu = \begin{pmatrix} x^2 + y^2 \\ 2 x^2 - 2 x y \end{pmatrix} \quad \mathrm{or} \quad \begin{pmatrix} 2 x^2 + y^2 + z^2 \\ 2 x^2 - 2xy \\ 2 x^2 - 2xz \end{pmatrix}
81*aa9a5b67SBarry Smith$$
82*aa9a5b67SBarry Smith
83*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
84*aa9a5b67SBarry Smith:append: '}'
85*aa9a5b67SBarry Smith:end-at: return 0;
86*aa9a5b67SBarry Smith:start-at: static PetscErrorCode quadratic_u
87*aa9a5b67SBarry Smith```
88*aa9a5b67SBarry Smith
89*aa9a5b67SBarry Smith$$
90*aa9a5b67SBarry Smithp = x + y - 1 \quad \mathrm{or} \quad x + y + z - \frac{3}{2}
91*aa9a5b67SBarry Smith$$
92*aa9a5b67SBarry Smith
93*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
94*aa9a5b67SBarry Smith:append: '}'
95*aa9a5b67SBarry Smith:end-at: return 0;
96*aa9a5b67SBarry Smith:start-at: static PetscErrorCode quadratic_p
97*aa9a5b67SBarry Smith```
98*aa9a5b67SBarry Smith
99*aa9a5b67SBarry SmithBy plugging these solutions into our equations, assuming that the velocity we choose is divergence-free, we can determine the body force necessary to make them satisfy the Stokes equations. For the quadratic solution above, we find
100*aa9a5b67SBarry Smith
101*aa9a5b67SBarry Smith$$
102*aa9a5b67SBarry Smithf = \begin{pmatrix} 1 - 4\mu \\ 1 - 4\mu \end{pmatrix} \quad \mathrm{or} \quad \begin{pmatrix} 1 - 8\mu \\ 1 - 4\mu \\ 1 - 4\mu \end{pmatrix}
103*aa9a5b67SBarry Smith$$
104*aa9a5b67SBarry Smith
105*aa9a5b67SBarry Smithwhich is implemented in our `f0_quadratic_u` pointwise function
106*aa9a5b67SBarry Smith
107*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
108*aa9a5b67SBarry Smith:end-at: '}'
109*aa9a5b67SBarry Smith:start-at: static void f0_quadratic_u
110*aa9a5b67SBarry Smith```
111*aa9a5b67SBarry Smith
112*aa9a5b67SBarry SmithWe let PETSc know about these solutions
113*aa9a5b67SBarry Smith
114*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
115*aa9a5b67SBarry Smith:end-at: PetscCall(PetscDSSetExactSolution(ds, 1
116*aa9a5b67SBarry Smith:start-at: PetscCall(PetscDSSetExactSolution(ds, 0
117*aa9a5b67SBarry Smith```
118*aa9a5b67SBarry Smith
119*aa9a5b67SBarry SmithThese solutions will be captured exactly by the $P_2-P_1$ finite element space. We can use the `-dmsnes_check` option to activate function space checks. It gives the $L_2$ error, or *discretization* error, of the exact solution, the residual computed using the interpolation of the exact solution into our finite element space, and uses a Taylor test to check that our Jacobian matches the residual. It should converge at order 2, or be exact in the case of linear equations like Stokes. Our $P_2-P_1$ runs in the PETSc test section at the bottom of the source file
120*aa9a5b67SBarry Smith
121*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
122*aa9a5b67SBarry Smith:lines: 1-3
123*aa9a5b67SBarry Smith:start-at: 'suffix: 2d_p2_p1_check'
124*aa9a5b67SBarry Smith```
125*aa9a5b67SBarry Smith
126*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
127*aa9a5b67SBarry Smith:lines: 1-3
128*aa9a5b67SBarry Smith:start-at: 'suffix: 3d_p2_p1_check'
129*aa9a5b67SBarry Smith```
130*aa9a5b67SBarry Smith
131*aa9a5b67SBarry Smithverify these claims, as we can see from the output files
132*aa9a5b67SBarry Smith
133*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/output/ex62_2d_p2_p1_check.out
134*aa9a5b67SBarry Smith:language: none
135*aa9a5b67SBarry Smith```
136*aa9a5b67SBarry Smith
137*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/output/ex62_3d_p2_p1_check.out
138*aa9a5b67SBarry Smith:language: none
139*aa9a5b67SBarry Smith```
140*aa9a5b67SBarry Smith
141*aa9a5b67SBarry SmithWe can carry out the same tests for the $Q_2-Q_1$ element,
142*aa9a5b67SBarry Smith
143*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
144*aa9a5b67SBarry Smith:lines: 1-2
145*aa9a5b67SBarry Smith:start-at: 'suffix: 2d_q2_q1_check'
146*aa9a5b67SBarry Smith```
147*aa9a5b67SBarry Smith
148*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
149*aa9a5b67SBarry Smith:lines: 1-2
150*aa9a5b67SBarry Smith:start-at: 'suffix: 3d_q2_q1_check'
151*aa9a5b67SBarry Smith```
152*aa9a5b67SBarry Smith
153*aa9a5b67SBarry SmithThe quadratic solution, however, cannot tell us whether our discretization is attaining the correct order of convergence, especially for higher order elements. Thus, we will define another solution based on trigonometric functions.
154*aa9a5b67SBarry Smith
155*aa9a5b67SBarry Smith$$
156*aa9a5b67SBarry Smithu = \begin{pmatrix} \sin(\pi x) + \sin(\pi y) \\ -\pi \cos(\pi x) y \end{pmatrix} \quad \mathrm{or} \quad
157*aa9a5b67SBarry Smith    \begin{pmatrix} 2 \sin(\pi x) + \sin(\pi y) + \sin(\pi z) \\ -\pi \cos(\pi x) y \\ -\pi \cos(\pi x) z \end{pmatrix}
158*aa9a5b67SBarry Smith$$
159*aa9a5b67SBarry Smith
160*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
161*aa9a5b67SBarry Smith:append: '}'
162*aa9a5b67SBarry Smith:end-at: return 0;
163*aa9a5b67SBarry Smith:start-at: static PetscErrorCode trig_u
164*aa9a5b67SBarry Smith```
165*aa9a5b67SBarry Smith
166*aa9a5b67SBarry Smith$$
167*aa9a5b67SBarry Smithp = \sin(2 \pi x) + \sin(2 \pi y) \quad \mathrm{or} \quad \sin(2 \pi x) + \sin(2 \pi y) + \sin(2 \pi z)
168*aa9a5b67SBarry Smith$$
169*aa9a5b67SBarry Smith
170*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
171*aa9a5b67SBarry Smith:append: '}'
172*aa9a5b67SBarry Smith:end-at: return 0;
173*aa9a5b67SBarry Smith:start-at: static PetscErrorCode trig_p
174*aa9a5b67SBarry Smith```
175*aa9a5b67SBarry Smith
176*aa9a5b67SBarry Smith$$
177*aa9a5b67SBarry Smithf = \begin{pmatrix} 2 \pi \cos(2 \pi x) + \mu \pi^2 \sin(\pi x) + \mu \pi^2 \sin(\pi y) \\ 2 \pi \cos(2 \pi y) - \mu \pi^3 \cos(\pi x) y \end{pmatrix} \quad \mathrm{or} \quad
178*aa9a5b67SBarry Smith\begin{pmatrix} 2 \pi \cos(2 \pi x) + 2\mu \pi^2 \sin(\pi x) + \mu \pi^2 \sin(\pi y) + \mu \pi^2 \sin(\pi z) \\ 2 \pi \cos(2 \pi y) - \mu \pi^3 cos(\pi x) y \\ 2 \pi \cos(2 \pi z) - \mu \pi^3 \cos(\pi x) z \end{pmatrix}
179*aa9a5b67SBarry Smith$$
180*aa9a5b67SBarry Smith
181*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
182*aa9a5b67SBarry Smith:append: '}'
183*aa9a5b67SBarry Smith:end-at: '}'
184*aa9a5b67SBarry Smith:start-at: static void f0_quadratic_u
185*aa9a5b67SBarry Smith```
186*aa9a5b67SBarry Smith
187*aa9a5b67SBarry SmithWe can now use `-snes_convergence_estimate` to determine the convergence exponent for the discretization. This options solves the problem on a series of refined meshes, calculates the error on each mesh, and determines the slope on a logarithmic scale. For example, we do this in two dimensions, refining our mesh twice using `-convest_num_refine 2` in the following test.
188*aa9a5b67SBarry Smith
189*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
190*aa9a5b67SBarry Smith:end-before: 'test:'
191*aa9a5b67SBarry Smith:start-at: 'suffix: 2d_p2_p1_conv'
192*aa9a5b67SBarry Smith```
193*aa9a5b67SBarry Smith
194*aa9a5b67SBarry SmithHowever, the test needs an accurate linear solver. Sparse LU factorizations do not, in general, do full pivoting. Thus we must deal with the zero pressure block explicitly. We use the `PCFIELDSPLIT` preconditioner and the full Schur complement factorization, but we still need a preconditioner for the Schur complement $B^T A^{-1} B$. We can have PETSc construct that matrix automatically, but the cost rises steeply as the problem size increases. Instead, we use the fact that the Schur complement is spectrally equivalent to the pressure mass matrix $M_p$. We can make a preconditioning matrix, which has the diagonal blocks we will use to build the preconditioners, letting PETSc know that we get the off-diagonal blocks from the original system with `-pc_fieldsplit_off_diag_use_amat` and to build the Schur complement from the original matrix using `-pc_use_amat`,
195*aa9a5b67SBarry Smith
196*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
197*aa9a5b67SBarry Smith:end-at: PetscCall(PetscDSSetJacobianPreconditioner(ds, 1
198*aa9a5b67SBarry Smith:start-at: PetscCall(PetscDSSetJacobianPreconditioner(ds, 0
199*aa9a5b67SBarry Smith```
200*aa9a5b67SBarry Smith
201*aa9a5b67SBarry SmithPutting this all together, and using exact solvers on the subblocks, we have
202*aa9a5b67SBarry Smith
203*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
204*aa9a5b67SBarry Smith:end-before: 'test:'
205*aa9a5b67SBarry Smith:start-at: 'suffix: 2d_p2_p1_conv'
206*aa9a5b67SBarry Smith```
207*aa9a5b67SBarry Smith
208*aa9a5b67SBarry Smithand we see it converges, however it is superconverging in the pressure,
209*aa9a5b67SBarry Smith
210*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/output/ex62_2d_p2_p1_conv.out
211*aa9a5b67SBarry Smith```
212*aa9a5b67SBarry Smith
213*aa9a5b67SBarry SmithIf we refine the mesh using `-dm_refine 3`, the convergence rates become `[3.0, 2.1]`.
214*aa9a5b67SBarry Smith
215*aa9a5b67SBarry Smith## Dealing with Parameters
216*aa9a5b67SBarry Smith
217*aa9a5b67SBarry SmithLike most physical problems, the Stokes problem has a parameter, the dynamic shear viscosity, which determines what solution regime we are in. To handle these parameters in PETSc, we first define a C struct to hold them,
218*aa9a5b67SBarry Smith
219*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
220*aa9a5b67SBarry Smith:end-at: '} Parameter;'
221*aa9a5b67SBarry Smith:start-at: typedef struct
222*aa9a5b67SBarry Smith```
223*aa9a5b67SBarry Smith
224*aa9a5b67SBarry Smithand then add a `PetscBag` object to our application context. We then setup the parameter object,
225*aa9a5b67SBarry Smith
226*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
227*aa9a5b67SBarry Smith:append: '}'
228*aa9a5b67SBarry Smith:end-at: PetscFunctionReturn(PETSC_SUCCESS);
229*aa9a5b67SBarry Smith:start-at: static PetscErrorCode SetupParameters
230*aa9a5b67SBarry Smith```
231*aa9a5b67SBarry Smith
232*aa9a5b67SBarry Smithwhich will allow us to set the value from the command line using `-mu`. The `PetscBag` can also be persisted to disk with `PetscBagLoad/View()`. We can make these values available as constant to our pointwise functions through the `PetscDS` object.
233*aa9a5b67SBarry Smith
234*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex62.c
235*aa9a5b67SBarry Smith:end-at: '}'
236*aa9a5b67SBarry Smith:start-after: /* Make constant values
237*aa9a5b67SBarry Smith```
238*aa9a5b67SBarry Smith
239*aa9a5b67SBarry Smith## Investigating convergence
240*aa9a5b67SBarry Smith
241*aa9a5b67SBarry SmithIn order to look at the convergence of some harder problems, we will examine `SNES ex69`. This example provides an exact solution to the variable viscosity Stokes equation. The sharp viscosity variation will allow us to investigate convergence of the solver and discretization. Briefly, a sharp viscosity variation is created across the unit square, imposed on a background pressure with given fundamental frequency. For example, we can create examples with period one half and viscosity $e^{2 B x}$ (solKx)
242*aa9a5b67SBarry Smith
243*aa9a5b67SBarry Smith```console
244*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_refine 5 -dm_view hdf5:$PETSC_DIR/sol.h5 -snes_view_solution hdf5:$PETSC_DIR/sol.h5::append -exact_vec_view hdf5:$PETSC_DIR/sol.h5::append -m 2 -n 2 -B 1"
245*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_refine 5 -dm_view hdf5:$PETSC_DIR/sol.h5 -snes_view_solution hdf5:$PETSC_DIR/sol.h5::append -exact_vec_view hdf5:$PETSC_DIR/sol.h5::append -m 2 -n 2 -B 3.75"
246*aa9a5b67SBarry Smith```
247*aa9a5b67SBarry Smith
248*aa9a5b67SBarry Smithwhich are show in the figure below.
249*aa9a5b67SBarry Smith
250*aa9a5b67SBarry Smith```{eval-rst}
251*aa9a5b67SBarry Smith.. list-table::
252*aa9a5b67SBarry Smith
253*aa9a5b67SBarry Smith  * - .. figure:: /images/tutorials/physics/ex69_sol_m_2_n_2_B_1.png
254*aa9a5b67SBarry Smith
255*aa9a5b67SBarry Smith         Solution for :math:`m=2`, :math:`n=2`, :math:`B=1`
256*aa9a5b67SBarry Smith
257*aa9a5b67SBarry Smith    - .. figure:: /images/tutorials/physics/ex69_sol_m_2_n_2_B_375.png
258*aa9a5b67SBarry Smith
259*aa9a5b67SBarry Smith         Solution for :math:`m=2`, :math:`n=2`, :math:`B=3.75`
260*aa9a5b67SBarry Smith```
261*aa9a5b67SBarry Smith
262*aa9a5b67SBarry Smith### Debugging
263*aa9a5b67SBarry Smith
264*aa9a5b67SBarry SmithIf we can provide the `PetscDS` object in our problem with the exact solution function, PETSc has good support for debugging our discretization and solver. We can use the `PetscConvEst` object to check the convergence behavior of our element automatically. For example, if we use the `-snes_convergence_estimate` option, PETSc will solve our nonlinear equations on a series of refined meshes, use our exact solution to calculate the error, and then fit this line on a log-log scale to get the convergence rate,
265*aa9a5b67SBarry Smith
266*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex69.c
267*aa9a5b67SBarry Smith:end-before: 'test:'
268*aa9a5b67SBarry Smith:start-at: 'suffix: p2p1_conv'
269*aa9a5b67SBarry Smith```
270*aa9a5b67SBarry Smith
271*aa9a5b67SBarry SmithIf we initially refine the mesh twice, `-dm_refine 2`, we get
272*aa9a5b67SBarry Smith
273*aa9a5b67SBarry Smith> L_2 convergence rate: [3.0, 2.2]
274*aa9a5b67SBarry Smith
275*aa9a5b67SBarry Smithwhich are the convergence rates we expect for the velocity and pressure using a $P_2-P_1$ discretization. For $Q_1-P_0$
276*aa9a5b67SBarry Smith
277*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex69.c
278*aa9a5b67SBarry Smith:end-before: 'test:'
279*aa9a5b67SBarry Smith:start-at: 'suffix: q1p0_conv'
280*aa9a5b67SBarry Smith```
281*aa9a5b67SBarry Smith
282*aa9a5b67SBarry Smithwe get
283*aa9a5b67SBarry Smith
284*aa9a5b67SBarry Smith> L_2 convergence rate: [2.0, 1.0]
285*aa9a5b67SBarry Smith
286*aa9a5b67SBarry SmithThis is a sensitive check that everything is working correctly. However, if this is wrong, where can I start? More fine-grained checks are available using the `-dmsnes_check` option. Using this for our $P_2-P_1$ example (the `p2p1` test), we have
287*aa9a5b67SBarry Smith
288*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/output/ex69_p2p1.out
289*aa9a5b67SBarry Smith```
290*aa9a5b67SBarry Smith
291*aa9a5b67SBarry SmithThe first line records the discretization error for our exact solution. This means that we project our solution function into the finite element space and then calculate the $L_2$ norm of the difference between the exact solution and its projection. The norm is computed for each field separately. Next, PETSc calculates the residual using the projected exact solution as input. This should be small, and as the mesh is refined it should approach zero. Last, PETSc uses a Taylor test to try and determine how the error in the linear model scales as a function of the perturbation $h$. Thus, in a nonlinear situation we would expect
292*aa9a5b67SBarry Smith
293*aa9a5b67SBarry Smith> Taylor approximation converging at order 2.0
294*aa9a5b67SBarry Smith
295*aa9a5b67SBarry SmithIn this case, since the viscosity does not depend on the velocity or pressure fields, we detect that the linear model is exact
296*aa9a5b67SBarry Smith
297*aa9a5b67SBarry Smith> Function appears to be linear
298*aa9a5b67SBarry Smith
299*aa9a5b67SBarry SmithSuppose that we have made an error in the Jacobian. For instance, let us accidentally flip the sign of the pressure term in the momentum Jacobian.
300*aa9a5b67SBarry Smith
301*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex69.c
302*aa9a5b67SBarry Smith:end-at: '}'
303*aa9a5b67SBarry Smith:start-at: static void stokes_momentum_pres_J
304*aa9a5b67SBarry Smith```
305*aa9a5b67SBarry Smith
306*aa9a5b67SBarry SmithWhen we run, we get a failure of the nonlinear solver. Our checking reveals that the Jacobian is wrong because it is converging at order 1 instead of 2, meaning the linear term is not correct in our model.
307*aa9a5b67SBarry Smith
308*aa9a5b67SBarry Smith```console
309*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_monitor -ksp_monitor_true_residual -ksp_converged_reason"
310*aa9a5b67SBarry SmithL_2 Error: [0.000439127, 0.0376629]
311*aa9a5b67SBarry SmithL_2 Residual: 0.0453958
312*aa9a5b67SBarry SmithTaylor approximation converging at order 1.00
313*aa9a5b67SBarry Smith  0 SNES Function norm 1.170604545948e-01
314*aa9a5b67SBarry Smith    0 KSP preconditioned resid norm 4.965098891419e-01 true resid norm 1.170604545948e-01 ||r(i)||/||b|| 1.000000000000e+00
315*aa9a5b67SBarry Smith    1 KSP preconditioned resid norm 9.236805404733e-11 true resid norm 1.460082233654e-12 ||r(i)||/||b|| 1.247289051378e-11
316*aa9a5b67SBarry Smith  Linear solve converged due to CONVERGED_ATOL iterations 1
317*aa9a5b67SBarry Smith[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
318*aa9a5b67SBarry Smith[0]PETSC ERROR:
319*aa9a5b67SBarry Smith[0]PETSC ERROR: SNESSolve has not converged
320*aa9a5b67SBarry Smith```
321*aa9a5b67SBarry Smith
322*aa9a5b67SBarry SmithIn order to track down the error, we can use `-snes_test_jacobian` which computes a finite difference approximation to the Jacobian and compares that to the analytic Jacobian. We ignore the first test, which occurs during our testing of the Jacobian, and look at the test that happens during the first Newton iterate. We see that the relative error in the Frobenius norm is about one percent, which indicates we have a real problem.
323*aa9a5b67SBarry Smith
324*aa9a5b67SBarry Smith```console
325*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -snes_test_jacobian"
326*aa9a5b67SBarry SmithL_2 Error: [0.000439127, 0.0376629]
327*aa9a5b67SBarry SmithL_2 Residual: 0.0453958
328*aa9a5b67SBarry Smith  ---------- Testing Jacobian -------------
329*aa9a5b67SBarry Smith  Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference
330*aa9a5b67SBarry Smith    of hand-coded and finite difference Jacobian entries greater than <threshold>.
331*aa9a5b67SBarry Smith  Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is
332*aa9a5b67SBarry Smith    O(1.e-8), the hand-coded Jacobian is probably correct.
333*aa9a5b67SBarry Smith  ||J - Jfd||_F/||J||_F = 136.793, ||J - Jfd||_F = 136.793
334*aa9a5b67SBarry Smith  ---------- Testing Jacobian for preconditioner -------------
335*aa9a5b67SBarry Smith  ||J - Jfd||_F/||J||_F = 136.793, ||J - Jfd||_F = 136.793
336*aa9a5b67SBarry SmithTaylor approximation converging at order 1.00
337*aa9a5b67SBarry Smith  0 SNES Function norm 1.170604545948e-01
338*aa9a5b67SBarry Smith  ---------- Testing Jacobian -------------
339*aa9a5b67SBarry Smith  ||J - Jfd||_F/||J||_F = 0.0119377, ||J - Jfd||_F = 1.63299
340*aa9a5b67SBarry Smith  ---------- Testing Jacobian for preconditioner -------------
341*aa9a5b67SBarry Smith  ||J - Jfd||_F/||J||_F = 0.008471, ||J - Jfd||_F = 1.15873
342*aa9a5b67SBarry Smith    0 KSP preconditioned resid norm 4.965098891419e-01 true resid norm 1.170604545948e-01 ||r(i)||/||b|| 1.000000000000e+00
343*aa9a5b67SBarry Smith    1 KSP preconditioned resid norm 9.236804064319e-11 true resid norm 1.460031196842e-12 ||r(i)||/||b|| 1.247245452699e-11
344*aa9a5b67SBarry Smith  Linear solve converged due to CONVERGED_ATOL iterations 1
345*aa9a5b67SBarry Smith[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
346*aa9a5b67SBarry Smith[0]PETSC ERROR:
347*aa9a5b67SBarry Smith[0]PETSC ERROR: SNESSolve has not converged
348*aa9a5b67SBarry Smith```
349*aa9a5b67SBarry Smith
350*aa9a5b67SBarry SmithAt this point, we could just go back and check the code. However, PETSc will also print out the differences between the analytic and approximate Jacobians. When we give the `-snes_test_jacobian_view` option, the code will print both Jacobians (which we omit) and then their difference, and will also do this for the preconditioning matrix (which we omit). It is clear from the output that the $u-p$ block of the Jacobian is wrong, and thus we know right where to look for our error. Moreover, if we look at the values in row 15, we see that the values just differ by a sign.
351*aa9a5b67SBarry Smith
352*aa9a5b67SBarry Smith```console
353*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -snes_test_jacobian"
354*aa9a5b67SBarry Smith        Hand-coded minus finite-difference Jacobian with tolerance 1e-05 ----------
355*aa9a5b67SBarry SmithMat Object: 1 MPI process
356*aa9a5b67SBarry Smith  type: seqaij
357*aa9a5b67SBarry Smithrow 0:
358*aa9a5b67SBarry Smithrow 1:
359*aa9a5b67SBarry Smithrow 2:
360*aa9a5b67SBarry Smithrow 3:
361*aa9a5b67SBarry Smithrow 4:
362*aa9a5b67SBarry Smithrow 5:
363*aa9a5b67SBarry Smithrow 6:
364*aa9a5b67SBarry Smithrow 7:
365*aa9a5b67SBarry Smithrow 8:
366*aa9a5b67SBarry Smithrow 9:
367*aa9a5b67SBarry Smithrow 10:
368*aa9a5b67SBarry Smithrow 11:
369*aa9a5b67SBarry Smithrow 12:
370*aa9a5b67SBarry Smithrow 13:
371*aa9a5b67SBarry Smithrow 14:
372*aa9a5b67SBarry Smithrow 15: (0, 0.166667)  (2, -0.166667)
373*aa9a5b67SBarry Smithrow 16: (0, 0.166667)  (2, -0.166667)  (5, 0.166667)  (8, -0.166667)
374*aa9a5b67SBarry Smithrow 17: (0, 0.166667)  (2, 0.166667)  (5, -0.166667)  (8, -0.166667)
375*aa9a5b67SBarry Smithrow 18: (0, 0.166667)  (5, -0.166667)
376*aa9a5b67SBarry Smithrow 19: (5, 0.166667)  (8, -0.166667)  (11, 0.166667)  (13, -0.166667)
377*aa9a5b67SBarry Smithrow 20: (5, 0.166667)  (8, 0.166667)  (11, -0.166667)  (13, -0.166667)
378*aa9a5b67SBarry Smithrow 21: (5, 0.166667)  (11, -0.166667)
379*aa9a5b67SBarry Smithrow 22: (5, 0.333333)  (8, -0.333333)
380*aa9a5b67SBarry Smithrow 23: (2, 0.166667)  (5, 0.166667)  (8, -0.166667)  (11, -0.166667)
381*aa9a5b67SBarry Smithrow 24: (2, 0.166667)  (3, -0.166667)  (5, 0.166667)  (8, -0.166667)
382*aa9a5b67SBarry Smithrow 25: (2, 0.333333)  (8, -0.333333)
383*aa9a5b67SBarry Smithrow 26: (2, 0.166667)  (3, -0.166667)  (8, 0.166667)  (10, -0.166667)
384*aa9a5b67SBarry Smithrow 27: (2, 0.166667)  (3, 0.166667)  (8, -0.166667)  (10, -0.166667)
385*aa9a5b67SBarry Smithrow 28: (3, 0.166667)  (10, -0.166667)
386*aa9a5b67SBarry Smithrow 29: (8, 0.333333)  (10, -0.333333)
387*aa9a5b67SBarry Smithrow 30: (3, 0.166667)  (8, 0.166667)  (10, -0.166667)  (13, -0.166667)
388*aa9a5b67SBarry Smithrow 31: (2, 0.166667)  (3, -0.166667)
389*aa9a5b67SBarry Smithrow 32: (8, 0.166667)  (10, -0.166667)  (13, 0.166667)  (14, -0.166667)
390*aa9a5b67SBarry Smithrow 33: (8, 0.166667)  (10, 0.166667)  (13, -0.166667)  (14, -0.166667)
391*aa9a5b67SBarry Smithrow 34: (10, 0.166667)  (14, -0.166667)
392*aa9a5b67SBarry Smithrow 35: (13, 0.166667)  (14, -0.166667)
393*aa9a5b67SBarry Smithrow 36: (8, 0.166667)  (10, -0.166667)  (11, 0.166667)  (13, -0.166667)
394*aa9a5b67SBarry Smithrow 37: (8, 0.333333)  (13, -0.333333)
395*aa9a5b67SBarry Smithrow 38: (11, 0.166667)  (13, -0.166667)
396*aa9a5b67SBarry Smith          0 KSP preconditioned resid norm 4.965098891419e-01 true resid norm 1.170604545948e-01 ||r(i)||/||b|| 1.000000000000e+00
397*aa9a5b67SBarry Smith          1 KSP preconditioned resid norm 9.236804067326e-11 true resid norm 1.460031196842e-12 ||r(i)||/||b|| 1.247245452699e-11
398*aa9a5b67SBarry Smith        Linear solve converged due to CONVERGED_ATOL iterations 1
399*aa9a5b67SBarry Smith      [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
400*aa9a5b67SBarry Smith      [0]PETSC ERROR:
401*aa9a5b67SBarry Smith      [0]PETSC ERROR: SNESSolve has not converged
402*aa9a5b67SBarry Smith```
403*aa9a5b67SBarry Smith
404*aa9a5b67SBarry SmithCan we see that the Schur complement of Q1-P0 is ill-conditioned?
405*aa9a5b67SBarry Smith
406*aa9a5b67SBarry Smith### Optimizing the Solver
407*aa9a5b67SBarry Smith
408*aa9a5b67SBarry SmithIn order to see exactly what solver we have employed, we can use the `-snes_view` option. When checking $P_2-P_1$ convergence, we use an exact solver, but it must have several parts in order to deal with the saddle-point in the Jacobian. Using the test system to provide our extra option, we get
409*aa9a5b67SBarry Smith
410*aa9a5b67SBarry Smith```console
411*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
412*aa9a5b67SBarry SmithSNES Object: 1 MPI process
413*aa9a5b67SBarry Smith  type: newtonls
414*aa9a5b67SBarry Smith  maximum iterations=50, maximum function evaluations=10000
415*aa9a5b67SBarry Smith  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
416*aa9a5b67SBarry Smith  total number of linear solver iterations=1
417*aa9a5b67SBarry Smith  total number of function evaluations=2
418*aa9a5b67SBarry Smith  norm schedule ALWAYS
419*aa9a5b67SBarry Smith  SNESLineSearch Object: 1 MPI process
420*aa9a5b67SBarry Smith    type: bt
421*aa9a5b67SBarry Smith      interpolation: cubic
422*aa9a5b67SBarry Smith      alpha=1.000000e-04
423*aa9a5b67SBarry Smith    maxstep=1.000000e+08, minlambda=1.000000e-12
424*aa9a5b67SBarry Smith    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
425*aa9a5b67SBarry Smith    maximum iterations=40
426*aa9a5b67SBarry Smith  KSP Object: 1 MPI process
427*aa9a5b67SBarry Smith    type: gmres
428*aa9a5b67SBarry Smith      restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
429*aa9a5b67SBarry Smith      happy breakdown tolerance 1e-30
430*aa9a5b67SBarry Smith    maximum iterations=10000, initial guess is zero
431*aa9a5b67SBarry Smith    tolerances: relative=1e-09, absolute=1e-10, divergence=10000.
432*aa9a5b67SBarry Smith    left preconditioning
433*aa9a5b67SBarry Smith    using PRECONDITIONED norm type for convergence test
434*aa9a5b67SBarry Smith  PC Object: 1 MPI process
435*aa9a5b67SBarry Smith    type: fieldsplit
436*aa9a5b67SBarry Smith      FieldSplit with Schur preconditioner, factorization FULL
437*aa9a5b67SBarry Smith      Preconditioner for the Schur complement formed from A11
438*aa9a5b67SBarry Smith      Split info:
439*aa9a5b67SBarry Smith      Split number 0 Defined by IS
440*aa9a5b67SBarry Smith      Split number 1 Defined by IS
441*aa9a5b67SBarry Smith      KSP solver for A00 block
442*aa9a5b67SBarry Smith        KSP Object: (fieldsplit_velocity_) 1 MPI process
443*aa9a5b67SBarry Smith          type: gmres
444*aa9a5b67SBarry Smith            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
445*aa9a5b67SBarry Smith            happy breakdown tolerance 1e-30
446*aa9a5b67SBarry Smith          maximum iterations=10000, initial guess is zero
447*aa9a5b67SBarry Smith          tolerances: relative=1e-05, absolute=1e-50, divergence=10000.
448*aa9a5b67SBarry Smith          left preconditioning
449*aa9a5b67SBarry Smith          using PRECONDITIONED norm type for convergence test
450*aa9a5b67SBarry Smith        PC Object: (fieldsplit_velocity_) 1 MPI process
451*aa9a5b67SBarry Smith          type: lu
452*aa9a5b67SBarry Smith            out-of-place factorization
453*aa9a5b67SBarry Smith            tolerance for zero pivot 2.22045e-14
454*aa9a5b67SBarry Smith            matrix ordering: nd
455*aa9a5b67SBarry Smith            factor fill ratio given 5., needed 1.15761
456*aa9a5b67SBarry Smith              Factored matrix follows:
457*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
458*aa9a5b67SBarry Smith                  type: seqaij
459*aa9a5b67SBarry Smith                  rows=30, cols=30
460*aa9a5b67SBarry Smith                  package used to perform factorization: petsc
461*aa9a5b67SBarry Smith                  total: nonzeros=426, allocated nonzeros=426
462*aa9a5b67SBarry Smith                    using I-node routines: found 17 nodes, limit used is 5
463*aa9a5b67SBarry Smith          linear system matrix followed by preconditioner matrix:
464*aa9a5b67SBarry Smith          Mat Object: 1 MPI process
465*aa9a5b67SBarry Smith            type: seqaij
466*aa9a5b67SBarry Smith            rows=30, cols=30
467*aa9a5b67SBarry Smith            total: nonzeros=368, allocated nonzeros=368
468*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
469*aa9a5b67SBarry Smith              using I-node routines: found 20 nodes, limit used is 5
470*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_velocity_) 1 MPI process
471*aa9a5b67SBarry Smith            type: seqaij
472*aa9a5b67SBarry Smith            rows=30, cols=30
473*aa9a5b67SBarry Smith            total: nonzeros=368, allocated nonzeros=368
474*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
475*aa9a5b67SBarry Smith              using I-node routines: found 20 nodes, limit used is 5
476*aa9a5b67SBarry Smith      KSP solver for S = A11 - A10 inv(A00) A01
477*aa9a5b67SBarry Smith        KSP Object: (fieldsplit_pressure_) 1 MPI process
478*aa9a5b67SBarry Smith          type: gmres
479*aa9a5b67SBarry Smith            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
480*aa9a5b67SBarry Smith            happy breakdown tolerance 1e-30
481*aa9a5b67SBarry Smith          maximum iterations=10000, initial guess is zero
482*aa9a5b67SBarry Smith          tolerances: relative=1e-09, absolute=1e-50, divergence=10000.
483*aa9a5b67SBarry Smith          left preconditioning
484*aa9a5b67SBarry Smith          using PRECONDITIONED norm type for convergence test
485*aa9a5b67SBarry Smith        PC Object: (fieldsplit_pressure_) 1 MPI process
486*aa9a5b67SBarry Smith          type: lu
487*aa9a5b67SBarry Smith            out-of-place factorization
488*aa9a5b67SBarry Smith            tolerance for zero pivot 2.22045e-14
489*aa9a5b67SBarry Smith            matrix ordering: nd
490*aa9a5b67SBarry Smith            factor fill ratio given 5., needed 1.2439
491*aa9a5b67SBarry Smith              Factored matrix follows:
492*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
493*aa9a5b67SBarry Smith                  type: seqaij
494*aa9a5b67SBarry Smith                  rows=9, cols=9
495*aa9a5b67SBarry Smith                  package used to perform factorization: petsc
496*aa9a5b67SBarry Smith                  total: nonzeros=51, allocated nonzeros=51
497*aa9a5b67SBarry Smith                    not using I-node routines
498*aa9a5b67SBarry Smith          linear system matrix followed by preconditioner matrix:
499*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_pressure_) 1 MPI process
500*aa9a5b67SBarry Smith            type: schurcomplement
501*aa9a5b67SBarry Smith            rows=9, cols=9
502*aa9a5b67SBarry Smith              has attached null space
503*aa9a5b67SBarry Smith              Schur complement A11 - A10 inv(A00) A01
504*aa9a5b67SBarry Smith              A11
505*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
506*aa9a5b67SBarry Smith                  type: seqaij
507*aa9a5b67SBarry Smith                  rows=9, cols=9
508*aa9a5b67SBarry Smith                  total: nonzeros=41, allocated nonzeros=41
509*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
510*aa9a5b67SBarry Smith                    has attached null space
511*aa9a5b67SBarry Smith                    not using I-node routines
512*aa9a5b67SBarry Smith              A10
513*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
514*aa9a5b67SBarry Smith                  type: seqaij
515*aa9a5b67SBarry Smith                  rows=9, cols=30
516*aa9a5b67SBarry Smith                  total: nonzeros=122, allocated nonzeros=122
517*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
518*aa9a5b67SBarry Smith                    not using I-node routines
519*aa9a5b67SBarry Smith              KSP solver for A00 block viewable with the additional option -fc_fieldsplit_velocity_ksp_view
520*aa9a5b67SBarry Smith              A01
521*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
522*aa9a5b67SBarry Smith                  type: seqaij
523*aa9a5b67SBarry Smith                  rows=30, cols=9
524*aa9a5b67SBarry Smith                  total: nonzeros=122, allocated nonzeros=122
525*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
526*aa9a5b67SBarry Smith                    using I-node routines: found 20 nodes, limit used is 5
527*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_pressure_) 1 MPI process
528*aa9a5b67SBarry Smith            type: seqaij
529*aa9a5b67SBarry Smith            rows=9, cols=9
530*aa9a5b67SBarry Smith            total: nonzeros=41, allocated nonzeros=41
531*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
532*aa9a5b67SBarry Smith              not using I-node routines
533*aa9a5b67SBarry Smith    linear system matrix followed by preconditioner matrix:
534*aa9a5b67SBarry Smith    Mat Object: 1 MPI process
535*aa9a5b67SBarry Smith      type: seqaij
536*aa9a5b67SBarry Smith      rows=39, cols=39
537*aa9a5b67SBarry Smith      total: nonzeros=653, allocated nonzeros=653
538*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
539*aa9a5b67SBarry Smith        has attached null space
540*aa9a5b67SBarry Smith        using I-node routines: found 24 nodes, limit used is 5
541*aa9a5b67SBarry Smith    Mat Object: (prec_) 1 MPI process
542*aa9a5b67SBarry Smith      type: seqaij
543*aa9a5b67SBarry Smith      rows=39, cols=39
544*aa9a5b67SBarry Smith      total: nonzeros=653, allocated nonzeros=653
545*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
546*aa9a5b67SBarry Smith        using I-node routines: found 24 nodes, limit used is 5
547*aa9a5b67SBarry Smith```
548*aa9a5b67SBarry Smith
549*aa9a5b67SBarry SmithGoing through this piece-by-piece, we can see all the parts of our solver. At the top level, we have a `SNES` using Newton's method
550*aa9a5b67SBarry Smith
551*aa9a5b67SBarry Smith```console
552*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
553*aa9a5b67SBarry SmithSNES Object: 1 MPI process
554*aa9a5b67SBarry Smith  type: newtonls
555*aa9a5b67SBarry Smith  maximum iterations=50, maximum function evaluations=10000
556*aa9a5b67SBarry Smith  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
557*aa9a5b67SBarry Smith  total number of linear solver iterations=1
558*aa9a5b67SBarry Smith  total number of function evaluations=2
559*aa9a5b67SBarry Smith  norm schedule ALWAYS
560*aa9a5b67SBarry Smith  SNESLineSearch Object: 1 MPI process
561*aa9a5b67SBarry Smith    type: bt
562*aa9a5b67SBarry Smith      interpolation: cubic
563*aa9a5b67SBarry Smith      alpha=1.000000e-04
564*aa9a5b67SBarry Smith    maxstep=1.000000e+08, minlambda=1.000000e-12
565*aa9a5b67SBarry Smith    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
566*aa9a5b67SBarry Smith    maximum iterations=40
567*aa9a5b67SBarry Smith```
568*aa9a5b67SBarry Smith
569*aa9a5b67SBarry SmithFor each nonlinear step, we use `KSPGMRES` to solve the Newton equation, preconditioned by `PCFIELDSPLIT`. We split the problem into two blocks, with the split determined by our `DM`, and combine those blocks using a Schur complement. The Schur complement is faithful since we use the `FULL` factorization type.
570*aa9a5b67SBarry Smith
571*aa9a5b67SBarry Smith```console
572*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
573*aa9a5b67SBarry Smith  KSP Object: 1 MPI process
574*aa9a5b67SBarry Smith    type: gmres
575*aa9a5b67SBarry Smith      restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
576*aa9a5b67SBarry Smith      happy breakdown tolerance 1e-30
577*aa9a5b67SBarry Smith    maximum iterations=10000, initial guess is zero
578*aa9a5b67SBarry Smith    tolerances: relative=1e-09, absolute=1e-10, divergence=10000.
579*aa9a5b67SBarry Smith    left preconditioning
580*aa9a5b67SBarry Smith    using PRECONDITIONED norm type for convergence test
581*aa9a5b67SBarry Smith  PC Object: 1 MPI process
582*aa9a5b67SBarry Smith    type: fieldsplit
583*aa9a5b67SBarry Smith      FieldSplit with Schur preconditioner, factorization FULL
584*aa9a5b67SBarry Smith      Preconditioner for the Schur complement formed from A11
585*aa9a5b67SBarry Smith      Split info:
586*aa9a5b67SBarry Smith      Split number 0 Defined by IS
587*aa9a5b67SBarry Smith      Split number 1 Defined by IS
588*aa9a5b67SBarry Smith```
589*aa9a5b67SBarry Smith
590*aa9a5b67SBarry SmithWe form the preconditioner for the Schur complement from the (1,1) block of our preconditioning matrix, which we have set to be the viscosity-weighted mass matrix
591*aa9a5b67SBarry Smith
592*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex69.c
593*aa9a5b67SBarry Smith:end-before: /*
594*aa9a5b67SBarry Smith:start-at: static void stokes_identity_J_kx
595*aa9a5b67SBarry Smith```
596*aa9a5b67SBarry Smith
597*aa9a5b67SBarry SmithThe solver for the first block, representing the velocity, is GMRES/LU. Note that the prefix is `fieldsplit_velocity_`, constructed automatically from the name of the field in our DM. Also note that there are two matrices, one from our original matrix, and one from our preconditioning matrix, but they are identical. In an optimized, scalable solver, this block would likely be solved by multigrid, but here we use LU for verification purposes.
598*aa9a5b67SBarry Smith
599*aa9a5b67SBarry Smith```console
600*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
601*aa9a5b67SBarry Smith      KSP solver for A00 block
602*aa9a5b67SBarry Smith        KSP Object: (fieldsplit_velocity_) 1 MPI process
603*aa9a5b67SBarry Smith          type: gmres
604*aa9a5b67SBarry Smith            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
605*aa9a5b67SBarry Smith            happy breakdown tolerance 1e-30
606*aa9a5b67SBarry Smith          maximum iterations=10000, initial guess is zero
607*aa9a5b67SBarry Smith          tolerances: relative=1e-05, absolute=1e-50, divergence=10000.
608*aa9a5b67SBarry Smith          left preconditioning
609*aa9a5b67SBarry Smith          using PRECONDITIONED norm type for convergence test
610*aa9a5b67SBarry Smith        PC Object: (fieldsplit_velocity_) 1 MPI process
611*aa9a5b67SBarry Smith          type: lu
612*aa9a5b67SBarry Smith            out-of-place factorization
613*aa9a5b67SBarry Smith            tolerance for zero pivot 2.22045e-14
614*aa9a5b67SBarry Smith            matrix ordering: nd
615*aa9a5b67SBarry Smith            factor fill ratio given 5., needed 1.15761
616*aa9a5b67SBarry Smith              Factored matrix follows:
617*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
618*aa9a5b67SBarry Smith                  type: seqaij
619*aa9a5b67SBarry Smith                  rows=30, cols=30
620*aa9a5b67SBarry Smith                  package used to perform factorization: petsc
621*aa9a5b67SBarry Smith                  total: nonzeros=426, allocated nonzeros=426
622*aa9a5b67SBarry Smith                    using I-node routines: found 17 nodes, limit used is 5
623*aa9a5b67SBarry Smith          linear system matrix followed by preconditioner matrix:
624*aa9a5b67SBarry Smith          Mat Object: 1 MPI process
625*aa9a5b67SBarry Smith            type: seqaij
626*aa9a5b67SBarry Smith            rows=30, cols=30
627*aa9a5b67SBarry Smith            total: nonzeros=368, allocated nonzeros=368
628*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
629*aa9a5b67SBarry Smith              using I-node routines: found 20 nodes, limit used is 5
630*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_velocity_) 1 MPI process
631*aa9a5b67SBarry Smith            type: seqaij
632*aa9a5b67SBarry Smith            rows=30, cols=30
633*aa9a5b67SBarry Smith            total: nonzeros=368, allocated nonzeros=368
634*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
635*aa9a5b67SBarry Smith              using I-node routines: found 20 nodes, limit used is 5
636*aa9a5b67SBarry Smith```
637*aa9a5b67SBarry Smith
638*aa9a5b67SBarry SmithThe solver for the second block, with prefix `fieldsplit_pressure_`, is also GMRES/LU, however we cannot factor the Schur complement operator since we never explicitly assemble it. Thus we assemble the viscosity-weighted mass matrix on the pressure space as an approximation. Notice that the Schur complement has the velocity solver embedded in it.
639*aa9a5b67SBarry Smith
640*aa9a5b67SBarry Smith```console
641*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
642*aa9a5b67SBarry Smith      KSP solver for S = A11 - A10 inv(A00) A01
643*aa9a5b67SBarry Smith        KSP Object: (fieldsplit_pressure_) 1 MPI process
644*aa9a5b67SBarry Smith          type: gmres
645*aa9a5b67SBarry Smith            restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
646*aa9a5b67SBarry Smith            happy breakdown tolerance 1e-30
647*aa9a5b67SBarry Smith          maximum iterations=10000, initial guess is zero
648*aa9a5b67SBarry Smith          tolerances: relative=1e-09, absolute=1e-50, divergence=10000.
649*aa9a5b67SBarry Smith          left preconditioning
650*aa9a5b67SBarry Smith          using PRECONDITIONED norm type for convergence test
651*aa9a5b67SBarry Smith        PC Object: (fieldsplit_pressure_) 1 MPI process
652*aa9a5b67SBarry Smith          type: lu
653*aa9a5b67SBarry Smith            out-of-place factorization
654*aa9a5b67SBarry Smith            tolerance for zero pivot 2.22045e-14
655*aa9a5b67SBarry Smith            matrix ordering: nd
656*aa9a5b67SBarry Smith            factor fill ratio given 5., needed 1.2439
657*aa9a5b67SBarry Smith              Factored matrix follows:
658*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
659*aa9a5b67SBarry Smith                  type: seqaij
660*aa9a5b67SBarry Smith                  rows=9, cols=9
661*aa9a5b67SBarry Smith                  package used to perform factorization: petsc
662*aa9a5b67SBarry Smith                  total: nonzeros=51, allocated nonzeros=51
663*aa9a5b67SBarry Smith                    not using I-node routines
664*aa9a5b67SBarry Smith          linear system matrix followed by preconditioner matrix:
665*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_pressure_) 1 MPI process
666*aa9a5b67SBarry Smith            type: schurcomplement
667*aa9a5b67SBarry Smith            rows=9, cols=9
668*aa9a5b67SBarry Smith              has attached null space
669*aa9a5b67SBarry Smith              Schur complement A11 - A10 inv(A00) A01
670*aa9a5b67SBarry Smith              A11
671*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
672*aa9a5b67SBarry Smith                  type: seqaij
673*aa9a5b67SBarry Smith                  rows=9, cols=9
674*aa9a5b67SBarry Smith                  total: nonzeros=41, allocated nonzeros=41
675*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
676*aa9a5b67SBarry Smith                    has attached null space
677*aa9a5b67SBarry Smith                    not using I-node routines
678*aa9a5b67SBarry Smith              A10
679*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
680*aa9a5b67SBarry Smith                  type: seqaij
681*aa9a5b67SBarry Smith                  rows=9, cols=30
682*aa9a5b67SBarry Smith                  total: nonzeros=122, allocated nonzeros=122
683*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
684*aa9a5b67SBarry Smith                    not using I-node routines
685*aa9a5b67SBarry Smith              KSP of A00
686*aa9a5b67SBarry Smith                KSP Object: (fieldsplit_velocity_) 1 MPI process
687*aa9a5b67SBarry Smith                  type: gmres
688*aa9a5b67SBarry Smith                    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
689*aa9a5b67SBarry Smith                    happy breakdown tolerance 1e-30
690*aa9a5b67SBarry Smith                  maximum iterations=10000, initial guess is zero
691*aa9a5b67SBarry Smith                  tolerances: relative=1e-05, absolute=1e-50, divergence=10000.
692*aa9a5b67SBarry Smith                  left preconditioning
693*aa9a5b67SBarry Smith                  using PRECONDITIONED norm type for convergence test
694*aa9a5b67SBarry Smith                PC Object: (fieldsplit_velocity_) 1 MPI process
695*aa9a5b67SBarry Smith                  type: lu
696*aa9a5b67SBarry Smith                    out-of-place factorization
697*aa9a5b67SBarry Smith                    tolerance for zero pivot 2.22045e-14
698*aa9a5b67SBarry Smith                    matrix ordering: nd
699*aa9a5b67SBarry Smith                    factor fill ratio given 5., needed 1.15761
700*aa9a5b67SBarry Smith                      Factored matrix follows:
701*aa9a5b67SBarry Smith                        Mat Object: 1 MPI process
702*aa9a5b67SBarry Smith                          type: seqaij
703*aa9a5b67SBarry Smith                          rows=30, cols=30
704*aa9a5b67SBarry Smith                          package used to perform factorization: petsc
705*aa9a5b67SBarry Smith                          total: nonzeros=426, allocated nonzeros=426
706*aa9a5b67SBarry Smith                            using I-node routines: found 17 nodes, limit used is 5
707*aa9a5b67SBarry Smith                  linear system matrix followed by preconditioner matrix:
708*aa9a5b67SBarry Smith                  Mat Object: 1 MPI process
709*aa9a5b67SBarry Smith                    type: seqaij
710*aa9a5b67SBarry Smith                    rows=30, cols=30
711*aa9a5b67SBarry Smith                    total: nonzeros=368, allocated nonzeros=368
712*aa9a5b67SBarry Smith                    total number of mallocs used during MatSetValues calls=0
713*aa9a5b67SBarry Smith                      using I-node routines: found 20 nodes, limit used is 5
714*aa9a5b67SBarry Smith                  Mat Object: (fieldsplit_velocity_) 1 MPI process
715*aa9a5b67SBarry Smith                    type: seqaij
716*aa9a5b67SBarry Smith                    rows=30, cols=30
717*aa9a5b67SBarry Smith                    total: nonzeros=368, allocated nonzeros=368
718*aa9a5b67SBarry Smith                    total number of mallocs used during MatSetValues calls=0
719*aa9a5b67SBarry Smith                      using I-node routines: found 20 nodes, limit used is 5
720*aa9a5b67SBarry Smith              A01
721*aa9a5b67SBarry Smith                Mat Object: 1 MPI process
722*aa9a5b67SBarry Smith                  type: seqaij
723*aa9a5b67SBarry Smith                  rows=30, cols=9
724*aa9a5b67SBarry Smith                  total: nonzeros=122, allocated nonzeros=122
725*aa9a5b67SBarry Smith                  total number of mallocs used during MatSetValues calls=0
726*aa9a5b67SBarry Smith                    using I-node routines: found 20 nodes, limit used is 5
727*aa9a5b67SBarry Smith          Mat Object: (fieldsplit_pressure_) 1 MPI process
728*aa9a5b67SBarry Smith            type: seqaij
729*aa9a5b67SBarry Smith            rows=9, cols=9
730*aa9a5b67SBarry Smith            total: nonzeros=41, allocated nonzeros=41
731*aa9a5b67SBarry Smith            total number of mallocs used during MatSetValues calls=0
732*aa9a5b67SBarry Smith              not using I-node routines
733*aa9a5b67SBarry Smith```
734*aa9a5b67SBarry Smith
735*aa9a5b67SBarry SmithFinally, the SNES viewer reports the system matrix and preconditioning matrix
736*aa9a5b67SBarry Smith
737*aa9a5b67SBarry Smith```console
738*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view"
739*aa9a5b67SBarry Smith    linear system matrix followed by preconditioner matrix:
740*aa9a5b67SBarry Smith    Mat Object: 1 MPI process
741*aa9a5b67SBarry Smith      type: seqaij
742*aa9a5b67SBarry Smith      rows=39, cols=39
743*aa9a5b67SBarry Smith      total: nonzeros=653, allocated nonzeros=653
744*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
745*aa9a5b67SBarry Smith        has attached null space
746*aa9a5b67SBarry Smith        using I-node routines: found 24 nodes, limit used is 5
747*aa9a5b67SBarry Smith    Mat Object: (prec_) 1 MPI process
748*aa9a5b67SBarry Smith      type: seqaij
749*aa9a5b67SBarry Smith      rows=39, cols=39
750*aa9a5b67SBarry Smith      total: nonzeros=653, allocated nonzeros=653
751*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
752*aa9a5b67SBarry Smith        using I-node routines: found 24 nodes, limit used is 5
753*aa9a5b67SBarry Smith```
754*aa9a5b67SBarry Smith
755*aa9a5b67SBarry SmithWe see that they have the same nonzero pattern, even though the preconditioning matrix only contains the diagonal blocks. This is because zeros were inserted to define the nonzero structure. We can remove these nonzeros by telling the DM not to insert zero at preallocation time, and also telling the matrix itself to ignore the zeros from the assembly process.
756*aa9a5b67SBarry Smith
757*aa9a5b67SBarry Smith```console
758*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_view -dm_preallocate_only -prec_mat_ignore_zero_entries"
759*aa9a5b67SBarry Smith    linear system matrix followed by preconditioner matrix:
760*aa9a5b67SBarry Smith    Mat Object: 1 MPI process
761*aa9a5b67SBarry Smith      type: seqaij
762*aa9a5b67SBarry Smith      rows=39, cols=39
763*aa9a5b67SBarry Smith      total: nonzeros=653, allocated nonzeros=653
764*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
765*aa9a5b67SBarry Smith        has attached null space
766*aa9a5b67SBarry Smith        using I-node routines: found 24 nodes, limit used is 5
767*aa9a5b67SBarry Smith    Mat Object: (prec_) 1 MPI process
768*aa9a5b67SBarry Smith      type: seqaij
769*aa9a5b67SBarry Smith      rows=39, cols=39
770*aa9a5b67SBarry Smith      total: nonzeros=409, allocated nonzeros=653
771*aa9a5b67SBarry Smith      total number of mallocs used during MatSetValues calls=0
772*aa9a5b67SBarry Smith        using I-node routines: found 29 nodes, limit used is 5
773*aa9a5b67SBarry Smith```
774*aa9a5b67SBarry Smith
775*aa9a5b67SBarry SmithWe can see a sparsity portrait of the system and preconditioning matrices if the installation supports X-windows visualization
776*aa9a5b67SBarry Smith
777*aa9a5b67SBarry Smith```console
778*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-ksp_view_mat draw -prec_mat_view draw -draw_pause -1"
779*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-ksp_view_mat draw -prec_mat_view draw -draw_save $PETSC_DIR/mat.png"
780*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_preallocate_only -mat_ignore_zero_entries -prec_mat_ignore_zero_entries -ksp_view_mat draw -prec_mat_view draw -draw_save $PETSC_DIR/mat_sparse.png"
781*aa9a5b67SBarry Smith```
782*aa9a5b67SBarry Smith
783*aa9a5b67SBarry Smith```{eval-rst}
784*aa9a5b67SBarry Smith.. list-table::
785*aa9a5b67SBarry Smith
786*aa9a5b67SBarry Smith  * - .. figure:: /images/tutorials/physics/stokes_p2p1_sys_mat.png
787*aa9a5b67SBarry Smith
788*aa9a5b67SBarry Smith         System matrix
789*aa9a5b67SBarry Smith
790*aa9a5b67SBarry Smith    - .. figure:: /images/tutorials/physics/stokes_p2p1_sys_mat_sparse.png
791*aa9a5b67SBarry Smith
792*aa9a5b67SBarry Smith         System matrix with sparse stencil
793*aa9a5b67SBarry Smith
794*aa9a5b67SBarry Smith  * - .. figure:: /images/tutorials/physics/stokes_p2p1_prec_mat.png
795*aa9a5b67SBarry Smith
796*aa9a5b67SBarry Smith         Preconditioning matrix
797*aa9a5b67SBarry Smith
798*aa9a5b67SBarry Smith    - .. figure:: /images/tutorials/physics/stokes_p2p1_prec_mat_sparse.png
799*aa9a5b67SBarry Smith
800*aa9a5b67SBarry Smith         Preconditioning matrix with sparse stencil
801*aa9a5b67SBarry Smith```
802*aa9a5b67SBarry Smith
803*aa9a5b67SBarry SmithIf we want to check the convergence of the solver, we can also do that using options. Both the linear and nonlinear solvers converge in a single iteration, which is exactly what we want. In order to have this happen, we must have the tolerance on both the outer KSP solver and the inner Schur complement solver be low enough. Notice that the sure complement solver is used twice, and converges in seven iterates each time.
804*aa9a5b67SBarry Smith
805*aa9a5b67SBarry Smith```console
806*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason"
807*aa9a5b67SBarry Smith0 SNES Function norm 1.170604545948e-01
808*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 7
809*aa9a5b67SBarry Smith  0 KSP preconditioned resid norm 4.965098891419e-01 true resid norm 1.170604545948e-01 ||r(i)||/||b|| 1.000000000000e+00
810*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 7
811*aa9a5b67SBarry Smith  1 KSP preconditioned resid norm 9.236813926190e-11 true resid norm 1.460072673561e-12 ||r(i)||/||b|| 1.247280884579e-11
812*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 1
813*aa9a5b67SBarry Smith1 SNES Function norm 1.460070661322e-12
814*aa9a5b67SBarry Smith```
815*aa9a5b67SBarry Smith
816*aa9a5b67SBarry SmithWe can look at the scalability of the solve by refining the mesh. We see that the Schur complement solve looks robust to grid refinement.
817*aa9a5b67SBarry Smith
818*aa9a5b67SBarry Smith```console
819*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_refine 2 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason"
820*aa9a5b67SBarry Smith0 SNES Function norm 3.503062983054e-02
821*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
822*aa9a5b67SBarry Smith  0 KSP preconditioned resid norm 9.943095979973e-01 true resid norm 3.503062983054e-02 ||r(i)||/||b|| 1.000000000000e+00
823*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
824*aa9a5b67SBarry Smith  1 KSP preconditioned resid norm 1.148772629230e-10 true resid norm 2.693482255004e-13 ||r(i)||/||b|| 7.688934706664e-12
825*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_RTOL iterations 1
826*aa9a5b67SBarry Smith1 SNES Function norm 2.693649920420e-13
827*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_refine 4 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason"
828*aa9a5b67SBarry Smith0 SNES Function norm 8.969202737759e-03
829*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 6
830*aa9a5b67SBarry Smith  0 KSP preconditioned resid norm 3.322375727167e+00 true resid norm 8.969202737759e-03 ||r(i)||/||b|| 1.000000000000e+00
831*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 6
832*aa9a5b67SBarry Smith  1 KSP preconditioned resid norm 6.112282404006e-10 true resid norm 8.543800889926e-14 ||r(i)||/||b|| 9.525708292843e-12
833*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_RTOL iterations 1
834*aa9a5b67SBarry Smith1 SNES Function norm 8.543893996362e-14
835*aa9a5b67SBarry Smith```
836*aa9a5b67SBarry Smith
837*aa9a5b67SBarry SmithStarting off with an exact solver allows us to check that the discretization, equations, and boundary conditions are correct. Moreover, choosing the Schur complement formulation, rather than a sparse direct solve, gives us a path to incremental boost the scalability. Our first step will be to replace the direct solve of the momentum operator, which has cost superlinear in $N$, with a more scalable alternative. Since the operator is still elliptic, despite the viscosity variation, we should be able to use some form of multigrid. We will start with algebraic multigrid because it handles coefficient variation well, even if the setup time is larger than the geometric variant.
838*aa9a5b67SBarry Smith
839*aa9a5b67SBarry Smith```console
840*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1" EXTRA_OPTIONS="-dm_refine 2 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_velocity_pc_type gamg -fieldsplit_velocity_ksp_converged_reason"
841*aa9a5b67SBarry Smith0 SNES Function norm 3.503062983054e-02
842*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
843*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
844*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
845*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
846*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 10
847*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
848*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
849*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
850*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
851*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
852*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
853*aa9a5b67SBarry Smith  0 KSP preconditioned resid norm 9.943097452179e-01 true resid norm 3.503062983054e-02 ||r(i)||/||b|| 1.000000000000e+00
854*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
855*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
856*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
857*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
858*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 10
859*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
860*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
861*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
862*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 8
863*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
864*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
865*aa9a5b67SBarry Smith  1 KSP preconditioned resid norm 1.503326145261e-05 true resid norm 1.089276827085e-06 ||r(i)||/||b|| 3.109498265814e-05
866*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 10
867*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
868*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
869*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
870*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
871*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
872*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
873*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
874*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
875*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 9
876*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 9
877*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 10
878*aa9a5b67SBarry Smith  2 KSP preconditioned resid norm 1.353007845554e-10 true resid norm 6.056095141823e-11 ||r(i)||/||b|| 1.728799959098e-09
879*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_RTOL iterations 2
880*aa9a5b67SBarry Smith1 SNES Function norm 6.056096909907e-11
881*aa9a5b67SBarry Smith```
882*aa9a5b67SBarry Smith
883*aa9a5b67SBarry SmithThis looks alright, but the number of iterates grows with refinement. At 3 refinements, it is 16, 30 at 4 refinements, and 70 at 5 refinements. Increasing the number of smoother iterates to four, `-fieldsplit_velocity_mg_levels_ksp_max_it 4`, brings down the number of iterates, but not the growth. Using w-cycles and full multigrid does not help either. It is likely that the coarse grids made by MIS are inaccurate for the $P_2$ discretization.
884*aa9a5b67SBarry Smith
885*aa9a5b67SBarry SmithWe can instead use geometric multigrid, and we would hope get more accurate coarse bases. The `-dm_refine_hierarchy` allows us to make a hierarchy of refined meshes and sets the number of multigrid levels automatically. Then all we need to specify is `-fieldsplit_velocity_pc_type mg`, as we see in the test
886*aa9a5b67SBarry Smith
887*aa9a5b67SBarry Smith```{literalinclude} /../src/snes/tutorials/ex69.c
888*aa9a5b67SBarry Smith:end-before: 'test:'
889*aa9a5b67SBarry Smith:start-at: 'suffix: p2p1_gmg'
890*aa9a5b67SBarry Smith```
891*aa9a5b67SBarry Smith
892*aa9a5b67SBarry SmithThis behaves well for the initial mesh,
893*aa9a5b67SBarry Smith
894*aa9a5b67SBarry Smith```console
895*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 2 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_velocity_ksp_converged_reason"
896*aa9a5b67SBarry Smith0 SNES Function norm 3.503062983054e-02
897*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 3.503062983054e-02 true resid norm 3.503062983054e-02 ||r(i)||/||b|| 1.000000000000e+00
898*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
899*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
900*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
901*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
902*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
903*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
904*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
905*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
906*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
907*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
908*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
909*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 4.643855168829e-06 true resid norm 4.643855168807e-06 ||r(i)||/||b|| 1.325655630878e-04
910*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
911*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
912*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
913*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
914*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
915*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
916*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
917*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
918*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
919*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
920*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 9
921*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
922*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.520240889941e-11 true resid norm 1.520239396618e-11 ||r(i)||/||b|| 4.339743258890e-10
923*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 2
924*aa9a5b67SBarry Smith1 SNES Function norm 1.520237877998e-11
925*aa9a5b67SBarry Smith```
926*aa9a5b67SBarry Smith
927*aa9a5b67SBarry Smithand is also stable under refinement
928*aa9a5b67SBarry Smith
929*aa9a5b67SBarry Smith```console
930*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 4 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_velocity_ksp_converged_reason"
931*aa9a5b67SBarry Smith0 SNES Function norm 3.503062983054e-02
932*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 3.503062983054e-02 true resid norm 3.503062983054e-02 ||r(i)||/||b|| 1.000000000000e+00
933*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
934*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
935*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
936*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
937*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
938*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
939*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
940*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
941*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
942*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 8
943*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
944*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 4.643855168829e-06 true resid norm 4.643855168807e-06 ||r(i)||/||b|| 1.325655630878e-04
945*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
946*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
947*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
948*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
949*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
950*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
951*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
952*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
953*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 4
954*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
955*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 9
956*aa9a5b67SBarry Smith  Linear fieldsplit_velocity_ solve converged due to CONVERGED_RTOL iterations 5
957*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.520240889941e-11 true resid norm 1.520239396618e-11 ||r(i)||/||b|| 4.339743258890e-10
958*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 2
959*aa9a5b67SBarry Smith1 SNES Function norm 1.520237877998e-11
960*aa9a5b67SBarry Smith```
961*aa9a5b67SBarry Smith
962*aa9a5b67SBarry SmithFinally, we can back off the pressure solve. `ILU(0)` is good enough to maintain a constant number of iterates as we refine the grid. We could continue to refine our preconditioner by playing with the tolerance of the inner multigrid and Schur complement solves, trading fewer inner iterates for more outer iterates.
963*aa9a5b67SBarry Smith
964*aa9a5b67SBarry Smith```console
965*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 2 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_pressure_pc_type ilu"
966*aa9a5b67SBarry Smith0 SNES Function norm 3.503062983054e-02
967*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 3.503062983054e-02 true resid norm 3.503062983054e-02 ||r(i)||/||b|| 1.000000000000e+00
968*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
969*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 4.643855785779e-06 true resid norm 4.643855785812e-06 ||r(i)||/||b|| 1.325655807011e-04
970*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
971*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.521944777036e-11 true resid norm 1.521942998859e-11 ||r(i)||/||b|| 4.344606437913e-10
972*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 2
973*aa9a5b67SBarry Smith1 SNES Function norm 1.521943449163e-11
974*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 4 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_pressure_pc_type ilu"
975*aa9a5b67SBarry Smith0 SNES Function norm 8.969202737759e-03
976*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 8.969202737759e-03 true resid norm 8.969202737759e-03 ||r(i)||/||b|| 1.000000000000e+00
977*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
978*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 2.234849111673e-05 true resid norm 2.234849111674e-05 ||r(i)||/||b|| 2.491692045566e-03
979*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
980*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.205594722917e-10 true resid norm 1.205594316079e-10 ||r(i)||/||b|| 1.344148807121e-08
981*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
982*aa9a5b67SBarry Smith  3 KSP unpreconditioned resid norm 1.461086575333e-15 true resid norm 2.284323415523e-15 ||r(i)||/||b|| 2.546852247977e-13
983*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 3
984*aa9a5b67SBarry Smith1 SNES Function norm 2.317901194143e-15
985*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 6 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_pressure_pc_type ilu"
986*aa9a5b67SBarry Smith0 SNES Function norm 2.252260693635e-03
987*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 2.252260693635e-03 true resid norm 2.252260693635e-03 ||r(i)||/||b|| 1.000000000000e+00
988*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 9
989*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 1.220195757583e-05 true resid norm 1.220195757579e-05 ||r(i)||/||b|| 5.417648858445e-03
990*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
991*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 2.683367607036e-09 true resid norm 2.683367591382e-09 ||r(i)||/||b|| 1.191410745197e-06
992*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 10
993*aa9a5b67SBarry Smith  3 KSP unpreconditioned resid norm 5.510932827474e-13 true resid norm 5.511665167379e-13 ||r(i)||/||b|| 2.447170162386e-10
994*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 3
995*aa9a5b67SBarry Smith1 SNES Function norm 5.511916500930e-13
996*aa9a5b67SBarry Smith```
997*aa9a5b67SBarry Smith
998*aa9a5b67SBarry SmithWe can make the problem harder by increasing the wave number and size of the viscosity perturbation. If we set the $B$ parameter to 6.9, we have a factor of one million increase in viscosity across the cell. At this scale, we see that we lose enough accuracy in our Jacobian calculation to defeat our Taylor test, but we are still able to solve the problem efficiently.
999*aa9a5b67SBarry Smith
1000*aa9a5b67SBarry Smith```console
1001*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 2 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_pressure_pc_type ilu -m 2 -n 2 -B 6.9"
1002*aa9a5b67SBarry SmithL_2 Error: [4.07817e-06, 0.0104694]
1003*aa9a5b67SBarry SmithL_2 Residual: 0.0145403
1004*aa9a5b67SBarry SmithTaylor approximation converging at order 1.00
1005*aa9a5b67SBarry Smith0 SNES Function norm 3.421266970274e-02
1006*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 3.421266970274e-02 true resid norm 3.421266970274e-02 ||r(i)||/||b|| 1.000000000000e+00
1007*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 21
1008*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 2.066264276201e-05 true resid norm 2.066264276201e-05 ||r(i)||/||b|| 6.039471032675e-04
1009*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 20
1010*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.295461366009e-10 true resid norm 1.295461419342e-10 ||r(i)||/||b|| 3.786496144842e-09
1011*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 20
1012*aa9a5b67SBarry Smith  3 KSP unpreconditioned resid norm 1.954355290546e-15 true resid norm 1.954135246291e-15 ||r(i)||/||b|| 5.711729786858e-14
1013*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 3
1014*aa9a5b67SBarry Smith1 SNES Function norm 1.946196473520e-15
1015*aa9a5b67SBarry Smith$ make -f ./gmakefile test globsearch="snes_tutorials-ex69_p2p1_gmg" EXTRA_OPTIONS="-dm_refine_hierarchy 6 -snes_monitor -ksp_monitor_true_residual -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason -fieldsplit_pressure_pc_type ilu -m 2 -n 2 -B 6.9"
1016*aa9a5b67SBarry SmithL_2 Error: [1.52905e-09, 4.72606e-05]
1017*aa9a5b67SBarry SmithL_2 Residual: 7.18836e-06
1018*aa9a5b67SBarry SmithTaylor approximation converging at order 1.00
1019*aa9a5b67SBarry Smith0 SNES Function norm 2.252034794902e-03
1020*aa9a5b67SBarry Smith  0 KSP unpreconditioned resid norm 2.252034794902e-03 true resid norm 2.252034794902e-03 ||r(i)||/||b|| 1.000000000000e+00
1021*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 19
1022*aa9a5b67SBarry Smith  1 KSP unpreconditioned resid norm 1.843225742581e-05 true resid norm 1.843225742582e-05 ||r(i)||/||b|| 8.184712539768e-03
1023*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 19
1024*aa9a5b67SBarry Smith  2 KSP unpreconditioned resid norm 1.410472862037e-09 true resid norm 1.410472860342e-09 ||r(i)||/||b|| 6.263104209294e-07
1025*aa9a5b67SBarry Smith  Linear fieldsplit_pressure_ solve converged due to CONVERGED_RTOL iterations 19
1026*aa9a5b67SBarry Smith  3 KSP unpreconditioned resid norm 1.051996270409e-14 true resid norm 1.064465321443e-14 ||r(i)||/||b|| 4.726682393419e-12
1027*aa9a5b67SBarry SmithLinear solve converged due to CONVERGED_ATOL iterations 3
1028*aa9a5b67SBarry Smith1 SNES Function norm 1.063917948054e-14
1029*aa9a5b67SBarry Smith```
1030*aa9a5b67SBarry Smith
1031*aa9a5b67SBarry Smith## Bibliography
1032*aa9a5b67SBarry Smith
1033*aa9a5b67SBarry Smith```{eval-rst}
1034*aa9a5b67SBarry Smith.. bibliography:: /petsc.bib
1035*aa9a5b67SBarry Smith   :filter: docname in docnames
1036*aa9a5b67SBarry Smith```
1037