xref: /petsc/src/snes/tutorials/ex71.c (revision 2065540a855ff9f9c49aa4d22d544ff2b07d8a79)
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   /* Domain and mesh definition */
38   PetscInt  dim;               /* The topological mesh dimension */
39   PetscBool simplex;           /* Use simplices or tensor product cells */
40   PetscInt  cells[3];          /* The initial domain division */
41   /* Problem definition */
42   PetscBag  bag;               /* Holds problem parameters */
43   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
44 } AppCtx;
45 
46 /*
47   In 2D, plane Poiseuille flow has exact solution:
48 
49     u = \Delta/(2 \nu) y (1 - y) + u_0
50     v = 0
51     p = -\Delta x
52     f = 0
53 
54   so that
55 
56     -\nu \Delta u + \nabla p + f = <\Delta, 0> + <-\Delta, 0> + <0, 0> = 0
57     \nabla \cdot u               = 0 + 0                               = 0
58 
59   In 3D we use exact solution:
60 
61     u = \Delta/(4 \nu) (y (1 - y) + z (1 - z)) + u_0
62     v = 0
63     w = 0
64     p = -\Delta x
65     f = 0
66 
67   so that
68 
69     -\nu \Delta u + \nabla p + f = <Delta, 0, 0> + <-Delta, 0, 0> + <0, 0, 0> = 0
70     \nabla \cdot u               = 0 + 0 + 0                                  = 0
71 
72   Note that these functions use coordinates X in the global (rotated) frame
73 */
74 PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
75 {
76   Parameter *param = (Parameter *) ctx;
77   PetscReal  Delta = param->Delta;
78   PetscReal  nu    = param->nu;
79   PetscReal  u_0   = param->u_0;
80   PetscReal  fac   = (PetscReal) (dim - 1);
81   PetscInt   d;
82 
83   u[0] = u_0;
84   for (d = 1; d < dim; ++d) u[0] += Delta/(fac * 2.0*nu) * X[d] * (1.0 - X[d]);
85   for (d = 1; d < dim; ++d) u[d]  = 0.0;
86   return 0;
87 }
88 
89 PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
90 {
91   Parameter *param = (Parameter *) ctx;
92   PetscReal  Delta = param->Delta;
93 
94   p[0] = -Delta * X[0];
95   return 0;
96 }
97 
98 PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
99 {
100   Parameter *param = (Parameter *) ctx;
101   PetscReal  u_0   = param->u_0;
102   PetscInt   d;
103 
104   u[0] = u_0;
105   for (d = 1; d < dim; ++d) u[d] = 0.0;
106   return 0;
107 }
108 
109 /* 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}
110    u[Ncomp]          = {p} */
111 void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
112           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
113           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
115 {
116   const PetscReal nu = PetscRealPart(constants[1]);
117   const PetscInt  Nc = dim;
118   PetscInt        c, d;
119 
120   for (c = 0; c < Nc; ++c) {
121     for (d = 0; d < dim; ++d) {
122       /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */
123       f1[c*dim+d] = nu*u_x[c*dim+d];
124     }
125     f1[c*dim+c] -= u[uOff[1]];
126   }
127 }
128 
129 /* 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} */
130 void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
131           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
132           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
134 {
135   PetscInt d;
136   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
137 }
138 
139 /* Residual functions are in reference coordinates */
140 static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
141                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
142                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
143                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
144 {
145   const PetscReal Delta = PetscRealPart(constants[0]);
146   PetscReal       alpha = PetscRealPart(constants[3]);
147   PetscReal       X     = PetscCosReal(alpha)*x[0] + PetscSinReal(alpha)*x[1];
148   PetscInt        d;
149 
150   for (d = 0; d < dim; ++d) {
151     f0[d] = -Delta * X * n[d];
152   }
153 }
154 
155 /* < q, \nabla\cdot u >
156    NcompI = 1, NcompJ = dim */
157 void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
158            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
159            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
160            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
161 {
162   PetscInt d;
163   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
164 }
165 
166 /* -< \nabla\cdot v, p >
167     NcompI = dim, NcompJ = 1 */
168 void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
169            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
170            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
171            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
172 {
173   PetscInt d;
174   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
175 }
176 
177 /* < \nabla v, \nabla u + {\nabla u}^T >
178    This just gives \nabla u, give the perdiagonal for the transpose */
179 void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
180            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
181            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
182            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
183 {
184   const PetscReal nu = PetscRealPart(constants[1]);
185   const PetscInt  Nc = dim;
186   PetscInt        c, d;
187 
188   for (c = 0; c < Nc; ++c) {
189     for (d = 0; d < dim; ++d) {
190       g3[((c*Nc+c)*dim+d)*dim+d] = nu;
191     }
192   }
193 }
194 
195 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
196 {
197   PetscInt       n = 3;
198   PetscErrorCode ierr;
199 
200   PetscFunctionBeginUser;
201   options->dim      = 2;
202   options->simplex  = PETSC_TRUE;
203   options->cells[0] = 3;
204   options->cells[1] = 3;
205   options->cells[2] = 3;
206 
207   ierr = PetscOptionsBegin(comm, "", "Poiseuille Flow Options", "DMPLEX");CHKERRQ(ierr);
208   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
209   ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
210   ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex62.c", options->cells, &n, NULL);CHKERRQ(ierr);
211   ierr = PetscOptionsEnd();
212   PetscFunctionReturn(0);
213 }
214 
215 static PetscErrorCode SetupParameters(AppCtx *user)
216 {
217   PetscBag       bag;
218   Parameter     *p;
219   PetscErrorCode ierr;
220 
221   PetscFunctionBeginUser;
222   /* setup PETSc parameter bag */
223   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
224   ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr);
225   bag  = user->bag;
226   ierr = PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length");CHKERRQ(ierr);
227   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
228   ierr = PetscBagRegisterReal(bag, &p->u_0,   0.0, "u_0",   "Tangential velocity at the wall");CHKERRQ(ierr);
229   ierr = PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis");CHKERRQ(ierr);
230   PetscFunctionReturn(0);
231 }
232 
233 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
234 {
235   PetscInt       dim = user->dim;
236   PetscErrorCode ierr;
237 
238   PetscFunctionBeginUser;
239   ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
240   {
241     Parameter   *param;
242     Vec          coordinates;
243     PetscScalar *coords;
244     PetscReal    alpha;
245     PetscInt     cdim, N, bs, i;
246 
247     ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr);
248     ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr);
249     ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
250     ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr);
251     if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim);
252     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
253     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
254     alpha = param->alpha;
255     for (i = 0; i < N; i += cdim) {
256       PetscScalar x = coords[i+0];
257       PetscScalar y = coords[i+1];
258 
259       coords[i+0] = PetscCosReal(alpha)*x - PetscSinReal(alpha)*y;
260       coords[i+1] = PetscSinReal(alpha)*x + PetscCosReal(alpha)*y;
261     }
262     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
263     ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr);
264   }
265   {
266     DM               pdm = NULL;
267     PetscPartitioner part;
268 
269     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
270     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
271     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
272     if (pdm) {
273       ierr = DMDestroy(dm);CHKERRQ(ierr);
274       *dm  = pdm;
275     }
276   }
277   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
278   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 PetscErrorCode SetupProblem(DM dm, AppCtx *user)
283 {
284   PetscDS        ds;
285   PetscWeakForm  wf;
286   DMLabel        label;
287   Parameter     *ctx;
288   PetscInt       id, bd;
289   PetscErrorCode ierr;
290 
291   PetscFunctionBeginUser;
292   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
293   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
294   ierr = PetscDSSetResidual(ds, 0, NULL, f1_u);CHKERRQ(ierr);
295   ierr = PetscDSSetResidual(ds, 1, f0_p, NULL);CHKERRQ(ierr);
296   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu);CHKERRQ(ierr);
297   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL);CHKERRQ(ierr);
298   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL);CHKERRQ(ierr);
299 
300   id   = 2;
301   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
302   ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr);
303   ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
304   ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr);
305   /* Setup constants */
306   {
307     Parameter  *param;
308     PetscScalar constants[4];
309 
310     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
311 
312     constants[0] = param->Delta;
313     constants[1] = param->nu;
314     constants[2] = param->u_0;
315     constants[3] = param->alpha;
316     ierr = PetscDSSetConstants(ds, 4, constants);CHKERRQ(ierr);
317   }
318   /* Setup Boundary Conditions */
319   id   = 3;
320   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall",    label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL);CHKERRQ(ierr);
321   id   = 1;
322   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) wall_velocity, NULL, ctx, NULL);CHKERRQ(ierr);
323   /* Setup exact solution */
324   user->exactFuncs[0] = quadratic_u;
325   user->exactFuncs[1] = linear_p;
326   ierr = PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], ctx);CHKERRQ(ierr);
327   ierr = PetscDSSetExactSolution(ds, 1, user->exactFuncs[1], ctx);CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
331 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
332 {
333   DM              cdm   = dm;
334   const PetscInt  dim   = user->dim;
335   PetscFE         fe[2];
336   Parameter      *param;
337   MPI_Comm        comm;
338   PetscErrorCode  ierr;
339 
340   PetscFunctionBeginUser;
341   /* Create finite element */
342   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
343   ierr = PetscFECreateDefault(comm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
344   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
345   ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
346   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
347   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
348   /* Set discretization and boundary conditions for each mesh */
349   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
350   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
351   ierr = DMCreateDS(dm);CHKERRQ(ierr);
352   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
353   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
354   while (cdm) {
355     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
356     ierr = DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0);CHKERRQ(ierr);
357     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
358   }
359   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
360   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 int main(int argc, char **argv)
365 {
366   SNES           snes; /* nonlinear solver */
367   DM             dm;   /* problem definition */
368   Vec            u, r; /* solution and residual */
369   AppCtx         user; /* user-defined work context */
370   PetscErrorCode ierr;
371 
372   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
373   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
374   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
375   ierr = SetupParameters(&user);CHKERRQ(ierr);
376   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
377   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
378   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
379   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
380   /* Setup problem */
381   ierr = PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);CHKERRQ(ierr);
382   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
383   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
384 
385   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
386   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
387 
388   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
389 
390   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
391 
392   {
393     Parameter *param;
394     void      *ctxs[2];
395 
396     ierr = PetscBagGetData(user.bag, (void **) &param);CHKERRQ(ierr);
397     ctxs[0] = ctxs[1] = param;
398     ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
399     ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr);
400     ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr);
401   }
402   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
403   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
404   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
405   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
406   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
407 
408   ierr = VecDestroy(&u);CHKERRQ(ierr);
409   ierr = VecDestroy(&r);CHKERRQ(ierr);
410   ierr = PetscFree(user.exactFuncs);CHKERRQ(ierr);
411   ierr = DMDestroy(&dm);CHKERRQ(ierr);
412   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
413   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
414   ierr = PetscFinalize();
415   return ierr;
416 }
417 
418 /*TEST
419 
420   # Convergence
421   test:
422     suffix: 2d_quad_q1_p0_conv
423     requires: !single
424     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 \
425       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
426       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
427       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
428       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
429         -fieldsplit_velocity_pc_type lu \
430         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
431   test:
432     suffix: 2d_quad_q1_p0_conv_u0
433     requires: !single
434     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \
435       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
436       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
437       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
438       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
439         -fieldsplit_velocity_pc_type lu \
440         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
441   test:
442     suffix: 2d_quad_q1_p0_conv_u0_alpha
443     requires: !single
444     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \
445       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
446       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
447       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
448       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
449         -fieldsplit_velocity_pc_type lu \
450         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
451   test:
452     suffix: 2d_quad_q1_p0_conv_gmg_vanka
453     requires: !single long_runtime
454     args: -simplex 0 -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \
455       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
456       -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \
457       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
458       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
459         -fieldsplit_velocity_pc_type mg \
460           -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \
461           -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \
462         -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi
463   test:
464     suffix: 2d_tri_p2_p1_conv
465     requires: triangle !single
466     args: -dm_plex_separate_marker -dm_refine 1 \
467       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
468       -dmsnes_check .001 -snes_error_if_not_converged \
469       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
470       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
471         -fieldsplit_velocity_pc_type lu \
472         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
473   test:
474     suffix: 2d_tri_p2_p1_conv_u0_alpha
475     requires: triangle !single
476     args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \
477       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
478       -dmsnes_check .001 -snes_error_if_not_converged \
479       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
480       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
481         -fieldsplit_velocity_pc_type lu \
482         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
483   test:
484     suffix: 2d_tri_p2_p1_conv_gmg_vcycle
485     requires: triangle !single
486     args: -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \
487       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
488       -dmsnes_check .001 -snes_error_if_not_converged \
489       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
490       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
491         -fieldsplit_velocity_pc_type mg \
492         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
493 TEST*/
494