xref: /petsc/src/snes/tutorials/ex71.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
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   Parameter *param = (Parameter *)ctx;
70   PetscReal  Delta = param->Delta;
71   PetscReal  nu    = param->nu;
72   PetscReal  u_0   = param->u_0;
73   PetscReal  fac   = (PetscReal)(dim - 1);
74   PetscInt   d;
75 
76   u[0] = u_0;
77   for (d = 1; d < dim; ++d) u[0] += Delta / (fac * 2.0 * nu) * X[d] * (1.0 - X[d]);
78   for (d = 1; d < dim; ++d) u[d] = 0.0;
79   return 0;
80 }
81 
82 PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
83   Parameter *param = (Parameter *)ctx;
84   PetscReal  Delta = param->Delta;
85 
86   p[0] = -Delta * X[0];
87   return 0;
88 }
89 
90 PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
91   Parameter *param = (Parameter *)ctx;
92   PetscReal  u_0   = param->u_0;
93   PetscInt   d;
94 
95   u[0] = u_0;
96   for (d = 1; d < dim; ++d) u[d] = 0.0;
97   return 0;
98 }
99 
100 /* 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}
101    u[Ncomp]          = {p} */
102 void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
103   const PetscReal nu = PetscRealPart(constants[1]);
104   const PetscInt  Nc = dim;
105   PetscInt        c, d;
106 
107   for (c = 0; c < Nc; ++c) {
108     for (d = 0; d < dim; ++d) {
109       /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */
110       f1[c * dim + d] = nu * u_x[c * dim + d];
111     }
112     f1[c * dim + c] -= u[uOff[1]];
113   }
114 }
115 
116 /* 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} */
117 void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
118   PetscInt d;
119   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
120 }
121 
122 /* Residual functions are in reference coordinates */
123 static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
124   const PetscReal Delta = PetscRealPart(constants[0]);
125   PetscReal       alpha = PetscRealPart(constants[3]);
126   PetscReal       X     = PetscCosReal(alpha) * x[0] + PetscSinReal(alpha) * x[1];
127   PetscInt        d;
128 
129   for (d = 0; d < dim; ++d) { f0[d] = -Delta * X * n[d]; }
130 }
131 
132 /* < q, \nabla\cdot u >
133    NcompI = 1, NcompJ = dim */
134 void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
135   PetscInt d;
136   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
137 }
138 
139 /* -< \nabla\cdot v, p >
140     NcompI = dim, NcompJ = 1 */
141 void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) {
142   PetscInt d;
143   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
144 }
145 
146 /* < \nabla v, \nabla u + {\nabla u}^T >
147    This just gives \nabla u, give the perdiagonal for the transpose */
148 void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
149   const PetscReal nu = PetscRealPart(constants[1]);
150   const PetscInt  Nc = dim;
151   PetscInt        c, d;
152 
153   for (c = 0; c < Nc; ++c) {
154     for (d = 0; d < dim; ++d) { g3[((c * Nc + c) * dim + d) * dim + d] = nu; }
155   }
156 }
157 
158 static PetscErrorCode SetupParameters(AppCtx *user) {
159   PetscBag   bag;
160   Parameter *p;
161 
162   PetscFunctionBeginUser;
163   /* setup PETSc parameter bag */
164   PetscCall(PetscBagGetData(user->bag, (void **)&p));
165   PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters"));
166   bag = user->bag;
167   PetscCall(PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length"));
168   PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity"));
169   PetscCall(PetscBagRegisterReal(bag, &p->u_0, 0.0, "u_0", "Tangential velocity at the wall"));
170   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis"));
171   PetscFunctionReturn(0);
172 }
173 
174 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
175   PetscFunctionBeginUser;
176   PetscCall(DMCreate(comm, dm));
177   PetscCall(DMSetType(*dm, DMPLEX));
178   PetscCall(DMSetFromOptions(*dm));
179   {
180     Parameter   *param;
181     Vec          coordinates;
182     PetscScalar *coords;
183     PetscReal    alpha;
184     PetscInt     cdim, N, bs, i;
185 
186     PetscCall(DMGetCoordinateDim(*dm, &cdim));
187     PetscCall(DMGetCoordinates(*dm, &coordinates));
188     PetscCall(VecGetLocalSize(coordinates, &N));
189     PetscCall(VecGetBlockSize(coordinates, &bs));
190     PetscCheck(bs == cdim, comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim);
191     PetscCall(VecGetArray(coordinates, &coords));
192     PetscCall(PetscBagGetData(user->bag, (void **)&param));
193     alpha = param->alpha;
194     for (i = 0; i < N; i += cdim) {
195       PetscScalar x = coords[i + 0];
196       PetscScalar y = coords[i + 1];
197 
198       coords[i + 0] = PetscCosReal(alpha) * x - PetscSinReal(alpha) * y;
199       coords[i + 1] = PetscSinReal(alpha) * x + PetscCosReal(alpha) * y;
200     }
201     PetscCall(VecRestoreArray(coordinates, &coords));
202     PetscCall(DMSetCoordinates(*dm, coordinates));
203   }
204   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
205   PetscFunctionReturn(0);
206 }
207 
208 PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
209   PetscDS       ds;
210   PetscWeakForm wf;
211   DMLabel       label;
212   Parameter    *ctx;
213   PetscInt      id, bd;
214 
215   PetscFunctionBeginUser;
216   PetscCall(PetscBagGetData(user->bag, (void **)&ctx));
217   PetscCall(DMGetDS(dm, &ds));
218   PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u));
219   PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL));
220   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
221   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
222   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
223 
224   id = 2;
225   PetscCall(DMGetLabel(dm, "marker", &label));
226   PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
227   PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
228   PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL));
229   /* Setup constants */
230   {
231     Parameter  *param;
232     PetscScalar constants[4];
233 
234     PetscCall(PetscBagGetData(user->bag, (void **)&param));
235 
236     constants[0] = param->Delta;
237     constants[1] = param->nu;
238     constants[2] = param->u_0;
239     constants[3] = param->alpha;
240     PetscCall(PetscDSSetConstants(ds, 4, constants));
241   }
242   /* Setup Boundary Conditions */
243   id = 3;
244   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall", label, 1, &id, 0, 0, NULL, (void (*)(void))wall_velocity, NULL, ctx, NULL));
245   id = 1;
246   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (void (*)(void))wall_velocity, NULL, ctx, NULL));
247   /* Setup exact solution */
248   PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, ctx));
249   PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, ctx));
250   PetscFunctionReturn(0);
251 }
252 
253 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) {
254   DM         cdm = dm;
255   PetscFE    fe[2];
256   Parameter *param;
257   PetscBool  simplex;
258   PetscInt   dim;
259   MPI_Comm   comm;
260 
261   PetscFunctionBeginUser;
262   PetscCall(DMGetDimension(dm, &dim));
263   PetscCall(DMPlexIsSimplex(dm, &simplex));
264   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
265   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
266   PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
267   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
268   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
269   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
270   /* Set discretization and boundary conditions for each mesh */
271   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
272   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
273   PetscCall(DMCreateDS(dm));
274   PetscCall(SetupProblem(dm, user));
275   PetscCall(PetscBagGetData(user->bag, (void **)&param));
276   while (cdm) {
277     PetscCall(DMCopyDisc(dm, cdm));
278     PetscCall(DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0));
279     PetscCall(DMGetCoarseDM(cdm, &cdm));
280   }
281   PetscCall(PetscFEDestroy(&fe[0]));
282   PetscCall(PetscFEDestroy(&fe[1]));
283   PetscFunctionReturn(0);
284 }
285 
286 int main(int argc, char **argv) {
287   SNES   snes; /* nonlinear solver */
288   DM     dm;   /* problem definition */
289   Vec    u, r; /* solution and residual */
290   AppCtx user; /* user-defined work context */
291 
292   PetscFunctionBeginUser;
293   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
294   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
295   PetscCall(SetupParameters(&user));
296   PetscCall(PetscBagSetFromOptions(user.bag));
297   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
298   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
299   PetscCall(SNESSetDM(snes, dm));
300   PetscCall(DMSetApplicationContext(dm, &user));
301   /* Setup problem */
302   PetscCall(SetupDiscretization(dm, &user));
303   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
304 
305   PetscCall(DMCreateGlobalVector(dm, &u));
306   PetscCall(VecDuplicate(u, &r));
307 
308   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
309 
310   PetscCall(SNESSetFromOptions(snes));
311 
312   {
313     PetscDS              ds;
314     PetscSimplePointFunc exactFuncs[2];
315     void                *ctxs[2];
316 
317     PetscCall(DMGetDS(dm, &ds));
318     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
319     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
320     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u));
321     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
322     PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view"));
323   }
324   PetscCall(DMSNESCheckFromOptions(snes, u));
325   PetscCall(VecSet(u, 0.0));
326   PetscCall(PetscObjectSetName((PetscObject)u, "Solution"));
327   PetscCall(SNESSolve(snes, NULL, u));
328   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
329 
330   PetscCall(VecDestroy(&u));
331   PetscCall(VecDestroy(&r));
332   PetscCall(DMDestroy(&dm));
333   PetscCall(SNESDestroy(&snes));
334   PetscCall(PetscBagDestroy(&user.bag));
335   PetscCall(PetscFinalize());
336   return 0;
337 }
338 
339 /*TEST
340 
341   # Convergence
342   test:
343     suffix: 2d_quad_q1_p0_conv
344     requires: !single
345     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 \
346       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
347       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
348       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
349       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
350         -fieldsplit_velocity_pc_type lu \
351         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
352   test:
353     suffix: 2d_quad_q1_p0_conv_u0
354     requires: !single
355     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \
356       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
357       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
358       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
359       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
360         -fieldsplit_velocity_pc_type lu \
361         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
362   test:
363     suffix: 2d_quad_q1_p0_conv_u0_alpha
364     requires: !single
365     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \
366       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
367       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
368       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
369       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
370         -fieldsplit_velocity_pc_type lu \
371         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
372   test:
373     suffix: 2d_quad_q1_p0_conv_gmg_vanka
374     requires: !single long_runtime
375     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \
376       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
377       -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \
378       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
379       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
380         -fieldsplit_velocity_pc_type mg \
381           -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \
382           -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \
383         -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi
384   test:
385     suffix: 2d_tri_p2_p1_conv
386     requires: triangle !single
387     args: -dm_plex_separate_marker -dm_refine 1 \
388       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
389       -dmsnes_check .001 -snes_error_if_not_converged \
390       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
391       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
392         -fieldsplit_velocity_pc_type lu \
393         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
394   test:
395     suffix: 2d_tri_p2_p1_conv_u0_alpha
396     requires: triangle !single
397     args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \
398       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
399       -dmsnes_check .001 -snes_error_if_not_converged \
400       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
401       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
402         -fieldsplit_velocity_pc_type lu \
403         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
404   test:
405     suffix: 2d_tri_p2_p1_conv_gmg_vcycle
406     requires: triangle !single
407     args: -dm_plex_separate_marker -dm_plex_box_faces 2,2 -dm_refine_hierarchy 1 \
408       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
409       -dmsnes_check .001 -snes_error_if_not_converged \
410       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
411       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
412         -fieldsplit_velocity_pc_type mg \
413         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
414 TEST*/
415