xref: /petsc/src/snes/tutorials/ex71.c (revision daad07d386296cdcbb87925ef5f1432ee4a24ec4)
1 static char help[] = "Poiseuille Flow in 2d and 3d channels with finite elements.\n\
2 We solve the Poiseuille flow problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4 
5 /*F
6 A Poiseuille flow is a steady-state isoviscous Stokes flow in a pipe of constant cross-section. We discretize using the
7 finite element method on an unstructured mesh. The weak form equations are
8 \begin{align*}
9   < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > + < v, \Delta \hat n >_{\Gamma_o} = 0
10   < q, \nabla\cdot u >                                                                                 = 0
11 \end{align*}
12 where $\nu$ is the kinematic viscosity, $\Delta$ is the pressure drop per unit length, assuming that pressure is 0 on
13 the left edge, and $\Gamma_o$ is the outlet boundary at the right edge of the pipe. The normal velocity will be zero at
14 the wall, but we will allow a fixed tangential velocity $u_0$.
15 
16 In order to test our global to local basis transformation, we will allow the pipe to be at an angle $\alpha$ to the
17 coordinate axes.
18 
19 For visualization, use
20 
21   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
22 F*/
23 
24 #include <petscdmplex.h>
25 #include <petscsnes.h>
26 #include <petscds.h>
27 #include <petscbag.h>
28 
29 typedef struct {
30   PetscReal Delta; /* Pressure drop per unit length */
31   PetscReal nu;    /* Kinematic viscosity */
32   PetscReal u_0;   /* Tangential velocity at the wall */
33   PetscReal alpha; /* Angle of pipe wall to x-axis */
34 } Parameter;
35 
36 typedef struct {
37   PetscBag bag;    /* Holds problem parameters */
38 } AppCtx;
39 
40 /*
41   In 2D, plane Poiseuille flow has exact solution:
42 
43     u = \Delta/(2 \nu) y (1 - y) + u_0
44     v = 0
45     p = -\Delta x
46     f = 0
47 
48   so that
49 
50     -\nu \Delta u + \nabla p + f = <\Delta, 0> + <-\Delta, 0> + <0, 0> = 0
51     \nabla \cdot u               = 0 + 0                               = 0
52 
53   In 3D we use exact solution:
54 
55     u = \Delta/(4 \nu) (y (1 - y) + z (1 - z)) + u_0
56     v = 0
57     w = 0
58     p = -\Delta x
59     f = 0
60 
61   so that
62 
63     -\nu \Delta u + \nabla p + f = <Delta, 0, 0> + <-Delta, 0, 0> + <0, 0, 0> = 0
64     \nabla \cdot u               = 0 + 0 + 0                                  = 0
65 
66   Note that these functions use coordinates X in the global (rotated) frame
67 */
68 PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
69 {
70   Parameter *param = (Parameter *) ctx;
71   PetscReal  Delta = param->Delta;
72   PetscReal  nu    = param->nu;
73   PetscReal  u_0   = param->u_0;
74   PetscReal  fac   = (PetscReal) (dim - 1);
75   PetscInt   d;
76 
77   u[0] = u_0;
78   for (d = 1; d < dim; ++d) u[0] += Delta/(fac * 2.0*nu) * X[d] * (1.0 - X[d]);
79   for (d = 1; d < dim; ++d) u[d]  = 0.0;
80   return 0;
81 }
82 
83 PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
84 {
85   Parameter *param = (Parameter *) ctx;
86   PetscReal  Delta = param->Delta;
87 
88   p[0] = -Delta * X[0];
89   return 0;
90 }
91 
92 PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
93 {
94   Parameter *param = (Parameter *) ctx;
95   PetscReal  u_0   = param->u_0;
96   PetscInt   d;
97 
98   u[0] = u_0;
99   for (d = 1; d < dim; ++d) u[d] = 0.0;
100   return 0;
101 }
102 
103 /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z}
104    u[Ncomp]          = {p} */
105 void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
106           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
107           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
108           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
109 {
110   const PetscReal nu = PetscRealPart(constants[1]);
111   const PetscInt  Nc = dim;
112   PetscInt        c, d;
113 
114   for (c = 0; c < Nc; ++c) {
115     for (d = 0; d < dim; ++d) {
116       /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */
117       f1[c*dim+d] = nu*u_x[c*dim+d];
118     }
119     f1[c*dim+c] -= u[uOff[1]];
120   }
121 }
122 
123 /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */
124 void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
125           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
126           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
127           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
128 {
129   PetscInt d;
130   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
131 }
132 
133 /* Residual functions are in reference coordinates */
134 static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
135                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
136                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
137                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
138 {
139   const PetscReal Delta = PetscRealPart(constants[0]);
140   PetscReal       alpha = PetscRealPart(constants[3]);
141   PetscReal       X     = PetscCosReal(alpha)*x[0] + PetscSinReal(alpha)*x[1];
142   PetscInt        d;
143 
144   for (d = 0; d < dim; ++d) {
145     f0[d] = -Delta * X * n[d];
146   }
147 }
148 
149 /* < q, \nabla\cdot u >
150    NcompI = 1, NcompJ = dim */
151 void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
152            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
153            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
154            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
155 {
156   PetscInt d;
157   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
158 }
159 
160 /* -< \nabla\cdot v, p >
161     NcompI = dim, NcompJ = 1 */
162 void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
163            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
164            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
165            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
166 {
167   PetscInt d;
168   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
169 }
170 
171 /* < \nabla v, \nabla u + {\nabla u}^T >
172    This just gives \nabla u, give the perdiagonal for the transpose */
173 void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
174            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
175            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
176            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
177 {
178   const PetscReal nu = PetscRealPart(constants[1]);
179   const PetscInt  Nc = dim;
180   PetscInt        c, d;
181 
182   for (c = 0; c < Nc; ++c) {
183     for (d = 0; d < dim; ++d) {
184       g3[((c*Nc+c)*dim+d)*dim+d] = nu;
185     }
186   }
187 }
188 
189 static PetscErrorCode SetupParameters(AppCtx *user)
190 {
191   PetscBag       bag;
192   Parameter     *p;
193 
194   PetscFunctionBeginUser;
195   /* setup PETSc parameter bag */
196   PetscCall(PetscBagGetData(user->bag, (void **) &p));
197   PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters"));
198   bag  = user->bag;
199   PetscCall(PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length"));
200   PetscCall(PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity"));
201   PetscCall(PetscBagRegisterReal(bag, &p->u_0,   0.0, "u_0",   "Tangential velocity at the wall"));
202   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis"));
203   PetscFunctionReturn(0);
204 }
205 
206 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
207 {
208   PetscFunctionBeginUser;
209   PetscCall(DMCreate(comm, dm));
210   PetscCall(DMSetType(*dm, DMPLEX));
211   PetscCall(DMSetFromOptions(*dm));
212   {
213     Parameter   *param;
214     Vec          coordinates;
215     PetscScalar *coords;
216     PetscReal    alpha;
217     PetscInt     cdim, N, bs, i;
218 
219     PetscCall(DMGetCoordinateDim(*dm, &cdim));
220     PetscCall(DMGetCoordinates(*dm, &coordinates));
221     PetscCall(VecGetLocalSize(coordinates, &N));
222     PetscCall(VecGetBlockSize(coordinates, &bs));
223     PetscCheck(bs == cdim,comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim);
224     PetscCall(VecGetArray(coordinates, &coords));
225     PetscCall(PetscBagGetData(user->bag, (void **) &param));
226     alpha = param->alpha;
227     for (i = 0; i < N; i += cdim) {
228       PetscScalar x = coords[i+0];
229       PetscScalar y = coords[i+1];
230 
231       coords[i+0] = PetscCosReal(alpha)*x - PetscSinReal(alpha)*y;
232       coords[i+1] = PetscSinReal(alpha)*x + PetscCosReal(alpha)*y;
233     }
234     PetscCall(VecRestoreArray(coordinates, &coords));
235     PetscCall(DMSetCoordinates(*dm, coordinates));
236   }
237   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
238   PetscFunctionReturn(0);
239 }
240 
241 PetscErrorCode SetupProblem(DM dm, AppCtx *user)
242 {
243   PetscDS        ds;
244   PetscWeakForm  wf;
245   DMLabel        label;
246   Parameter     *ctx;
247   PetscInt       id, bd;
248 
249   PetscFunctionBeginUser;
250   PetscCall(PetscBagGetData(user->bag, (void **) &ctx));
251   PetscCall(DMGetDS(dm, &ds));
252   PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u));
253   PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL));
254   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu));
255   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL));
256   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
257 
258   id   = 2;
259   PetscCall(DMGetLabel(dm, "marker", &label));
260   PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
261   PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
262   PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL));
263   /* Setup constants */
264   {
265     Parameter  *param;
266     PetscScalar constants[4];
267 
268     PetscCall(PetscBagGetData(user->bag, (void **) &param));
269 
270     constants[0] = param->Delta;
271     constants[1] = param->nu;
272     constants[2] = param->u_0;
273     constants[3] = param->alpha;
274     PetscCall(PetscDSSetConstants(ds, 4, constants));
275   }
276   /* Setup Boundary Conditions */
277   id   = 3;
278   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall",    label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL));
279   id   = 1;
280   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL));
281   /* Setup exact solution */
282   PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, ctx));
283   PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, ctx));
284   PetscFunctionReturn(0);
285 }
286 
287 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
288 {
289   DM             cdm = dm;
290   PetscFE        fe[2];
291   Parameter     *param;
292   PetscBool      simplex;
293   PetscInt       dim;
294   MPI_Comm       comm;
295 
296   PetscFunctionBeginUser;
297   PetscCall(DMGetDimension(dm, &dim));
298   PetscCall(DMPlexIsSimplex(dm, &simplex));
299   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
300   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
301   PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity"));
302   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
303   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
304   PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure"));
305   /* Set discretization and boundary conditions for each mesh */
306   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
307   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
308   PetscCall(DMCreateDS(dm));
309   PetscCall(SetupProblem(dm, user));
310   PetscCall(PetscBagGetData(user->bag, (void **) &param));
311   while (cdm) {
312     PetscCall(DMCopyDisc(dm, cdm));
313     PetscCall(DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0));
314     PetscCall(DMGetCoarseDM(cdm, &cdm));
315   }
316   PetscCall(PetscFEDestroy(&fe[0]));
317   PetscCall(PetscFEDestroy(&fe[1]));
318   PetscFunctionReturn(0);
319 }
320 
321 int main(int argc, char **argv)
322 {
323   SNES           snes; /* nonlinear solver */
324   DM             dm;   /* problem definition */
325   Vec            u, r; /* solution and residual */
326   AppCtx         user; /* user-defined work context */
327 
328   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
329   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
330   PetscCall(SetupParameters(&user));
331   PetscCall(PetscBagSetFromOptions(user.bag));
332   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
333   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
334   PetscCall(SNESSetDM(snes, dm));
335   PetscCall(DMSetApplicationContext(dm, &user));
336   /* Setup problem */
337   PetscCall(SetupDiscretization(dm, &user));
338   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
339 
340   PetscCall(DMCreateGlobalVector(dm, &u));
341   PetscCall(VecDuplicate(u, &r));
342 
343   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
344 
345   PetscCall(SNESSetFromOptions(snes));
346 
347   {
348     PetscDS              ds;
349     PetscSimplePointFunc exactFuncs[2];
350     void                *ctxs[2];
351 
352     PetscCall(DMGetDS(dm, &ds));
353     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
354     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
355     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u));
356     PetscCall(PetscObjectSetName((PetscObject) u, "Exact Solution"));
357     PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view"));
358   }
359   PetscCall(DMSNESCheckFromOptions(snes, u));
360   PetscCall(VecSet(u, 0.0));
361   PetscCall(PetscObjectSetName((PetscObject) u, "Solution"));
362   PetscCall(SNESSolve(snes, NULL, u));
363   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
364 
365   PetscCall(VecDestroy(&u));
366   PetscCall(VecDestroy(&r));
367   PetscCall(DMDestroy(&dm));
368   PetscCall(SNESDestroy(&snes));
369   PetscCall(PetscBagDestroy(&user.bag));
370   PetscCall(PetscFinalize());
371   return 0;
372 }
373 
374 /*TEST
375 
376   # Convergence
377   test:
378     suffix: 2d_quad_q1_p0_conv
379     requires: !single
380     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 \
381       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
382       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
383       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
384       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
385         -fieldsplit_velocity_pc_type lu \
386         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
387   test:
388     suffix: 2d_quad_q1_p0_conv_u0
389     requires: !single
390     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \
391       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
392       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
393       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
394       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
395         -fieldsplit_velocity_pc_type lu \
396         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
397   test:
398     suffix: 2d_quad_q1_p0_conv_u0_alpha
399     requires: !single
400     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \
401       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
402       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
403       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
404       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
405         -fieldsplit_velocity_pc_type lu \
406         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
407   test:
408     suffix: 2d_quad_q1_p0_conv_gmg_vanka
409     requires: !single long_runtime
410     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \
411       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
412       -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \
413       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
414       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
415         -fieldsplit_velocity_pc_type mg \
416           -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \
417           -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \
418         -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi
419   test:
420     suffix: 2d_tri_p2_p1_conv
421     requires: triangle !single
422     args: -dm_plex_separate_marker -dm_refine 1 \
423       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
424       -dmsnes_check .001 -snes_error_if_not_converged \
425       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
426       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
427         -fieldsplit_velocity_pc_type lu \
428         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
429   test:
430     suffix: 2d_tri_p2_p1_conv_u0_alpha
431     requires: triangle !single
432     args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \
433       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
434       -dmsnes_check .001 -snes_error_if_not_converged \
435       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
436       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
437         -fieldsplit_velocity_pc_type lu \
438         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
439   test:
440     suffix: 2d_tri_p2_p1_conv_gmg_vcycle
441     requires: triangle !single
442     args: -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \
443       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
444       -dmsnes_check .001 -snes_error_if_not_converged \
445       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
446       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
447         -fieldsplit_velocity_pc_type mg \
448         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
449 TEST*/
450