xref: /petsc/src/snes/tutorials/ex62.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
1 static char help[] = "Stokes Problem discretized with finite elements,\n\
2 using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n";
3 
4 /*
5 For the isoviscous Stokes problem, which we discretize using the finite
6 element method on an unstructured mesh, the weak form equations are
7 
8   < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0
9   < q, -\nabla\cdot u >                                                   = 0
10 
11 Viewing:
12 
13 To produce nice output, use
14 
15   -dm_refine 3 -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -snes_view_solution hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append
16 
17 You can get a LaTeX view of the mesh, with point numbering using
18 
19   -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0
20 
21 The data layout can be viewed using
22 
23   -dm_petscsection_view
24 
25 Lots of information about the FEM assembly can be printed using
26 
27   -dm_plex_print_fem 3
28 */
29 
30 #include <petscdmplex.h>
31 #include <petscsnes.h>
32 #include <petscds.h>
33 #include <petscbag.h>
34 
35 // TODO: Plot residual by fields after each smoother iterate
36 
37 typedef enum {SOL_QUADRATIC, SOL_TRIG, SOL_UNKNOWN} SolType;
38 const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0};
39 
40 typedef struct {
41   PetscScalar mu; /* dynamic shear viscosity */
42 } Parameter;
43 
44 typedef struct {
45   PetscBag bag; /* Problem parameters */
46   SolType  sol; /* MMS solution */
47 } AppCtx;
48 
49 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
50                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
51                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
52                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
53 {
54   const PetscReal mu = PetscRealPart(constants[0]);
55   const PetscInt  Nc = uOff[1]-uOff[0];
56   PetscInt        c, d;
57 
58   for (c = 0; c < Nc; ++c) {
59     for (d = 0; d < dim; ++d) {
60       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
61     }
62     f1[c*dim+c] -= u[uOff[1]];
63   }
64 }
65 
66 static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
67                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
68                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
69                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
70 {
71   PetscInt d;
72   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d];
73 }
74 
75 static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
76                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
77                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
78                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
79 {
80   PetscInt d;
81   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */
82 }
83 
84 static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
85                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
86                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
87                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
88 {
89   PetscInt d;
90   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */
91 }
92 
93 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
97 {
98   const PetscReal mu = PetscRealPart(constants[0]);
99   const PetscInt  Nc = uOff[1]-uOff[0];
100   PetscInt        c, d;
101 
102   for (c = 0; c < Nc; ++c) {
103     for (d = 0; d < dim; ++d) {
104       g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */
105       g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */
106     }
107   }
108 }
109 
110 static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
111                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
112                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
113                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
114 {
115   const PetscReal mu = PetscRealPart(constants[0]);
116 
117   g0[0] = 1.0/mu;
118 }
119 
120 /* Quadratic MMS Solution
121    2D:
122 
123      u = x^2 + y^2
124      v = 2 x^2 - 2xy
125      p = x + y - 1
126      f = <1 - 4 mu, 1 - 4 mu>
127 
128    so that
129 
130      e(u) = (grad u + grad u^T) = / 4x  4x \
131                                   \ 4x -4x /
132      div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
133      \nabla \cdot u             = 2x - 2x = 0
134 
135    3D:
136 
137      u = 2 x^2 + y^2 + z^2
138      v = 2 x^2 - 2xy
139      w = 2 x^2 - 2xz
140      p = x + y + z - 3/2
141      f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
142 
143    so that
144 
145      e(u) = (grad u + grad u^T) = / 8x  4x  4x \
146                                   | 4x -4x  0  |
147                                   \ 4x  0  -4x /
148      div mu e(u) - \nabla p + f = mu <8, 4, 4> - <1, 1, 1> + <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> = 0
149      \nabla \cdot u             = 4x - 2x - 2x = 0
150 */
151 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
152 {
153   PetscInt c;
154 
155   u[0] = (dim-1)*PetscSqr(x[0]);
156   for (c = 1; c < Nc; ++c) {
157     u[0] += PetscSqr(x[c]);
158     u[c]  = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c];
159   }
160   return 0;
161 }
162 
163 static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
164 {
165   PetscInt d;
166 
167   u[0] = -0.5*dim;
168   for (d = 0; d < dim; ++d) u[0] += x[d];
169   return 0;
170 }
171 
172 static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
173                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
174                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
175                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
176 {
177   const PetscReal mu = PetscRealPart(constants[0]);
178   PetscInt        d;
179 
180   f0[0] = (dim-1)*4.0*mu - 1.0;
181   for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0;
182 }
183 
184 /* Trigonometric MMS Solution
185    2D:
186 
187      u = sin(pi x) + sin(pi y)
188      v = -pi cos(pi x) y
189      p = sin(2 pi x) + sin(2 pi y)
190      f = <2pi cos(2 pi x) + mu pi^2 sin(pi x) + mu pi^2 sin(pi y), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y>
191 
192    so that
193 
194      e(u) = (grad u + grad u^T) = /        2pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y \
195                                   \ pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)          /
196      div mu e(u) - \nabla p + f = mu <-pi^2 sin(pi x) - pi^2 sin(pi y), pi^3 cos(pi x) y> - <2pi cos(2 pi x), 2pi cos(2 pi y)> + <f_x, f_y> = 0
197      \nabla \cdot u             = pi cos(pi x) - pi cos(pi x) = 0
198 
199    3D:
200 
201      u = 2 sin(pi x) + sin(pi y) + sin(pi z)
202      v = -pi cos(pi x) y
203      w = -pi cos(pi x) z
204      p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
205      f = <2pi cos(2 pi x) + mu 2pi^2 sin(pi x) + mu pi^2 sin(pi y) + mu pi^2 sin(pi z), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y, 2pi cos(2 pi z) - mu pi^3 cos(pi x) z>
206 
207    so that
208 
209      e(u) = (grad u + grad u^T) = /        4pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y  pi cos(pi z) + pi^2 sin(pi x) z \
210                                   | pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)                        0                  |
211                                   \ pi cos(pi z) + pi^2 sin(pi x) z               0                         -2pi cos(pi x)            /
212      div mu e(u) - \nabla p + f = mu <-2pi^2 sin(pi x) - pi^2 sin(pi y) - pi^2 sin(pi z), pi^3 cos(pi x) y, pi^3 cos(pi x) z> - <2pi cos(2 pi x), 2pi cos(2 pi y), 2pi cos(2 pi z)> + <f_x, f_y, f_z> = 0
213      \nabla \cdot u             = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
214 */
215 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
216 {
217   PetscInt c;
218 
219   u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]);
220   for (c = 1; c < Nc; ++c) {
221     u[0] += PetscSinReal(PETSC_PI*x[c]);
222     u[c]  = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c];
223   }
224   return 0;
225 }
226 
227 static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
228 {
229   PetscInt d;
230 
231   for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]);
232   return 0;
233 }
234 
235 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
236                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
237                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
238                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
239 {
240   const PetscReal mu = PetscRealPart(constants[0]);
241   PetscInt        d;
242 
243   f0[0] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[0]) - (dim-1)*mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[0]);
244   for (d = 1; d < dim; ++d) {
245     f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]);
246     f0[d]  = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[d]) + mu*PetscPowRealInt(PETSC_PI, 3)*PetscCosReal(PETSC_PI*x[0])*x[d];
247   }
248 }
249 
250 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
251 {
252   PetscInt       sol;
253 
254   PetscFunctionBeginUser;
255   options->sol = SOL_QUADRATIC;
256   PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
257   sol  = options->sol;
258   PetscCall(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, PETSC_STATIC_ARRAY_LENGTH(SolTypes)-3, SolTypes[options->sol], &sol, NULL));
259   options->sol = (SolType) sol;
260   PetscOptionsEnd();
261   PetscFunctionReturn(0);
262 }
263 
264 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
265 {
266   PetscFunctionBeginUser;
267   PetscCall(DMCreate(comm, dm));
268   PetscCall(DMSetType(*dm, DMPLEX));
269   PetscCall(DMSetFromOptions(*dm));
270   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
271   PetscFunctionReturn(0);
272 }
273 
274 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
275 {
276   Parameter     *p;
277 
278   PetscFunctionBeginUser;
279   /* setup PETSc parameter bag */
280   PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag));
281   PetscCall(PetscBagGetData(ctx->bag, (void **) &p));
282   PetscCall(PetscBagSetName(ctx->bag, "par", "Stokes Parameters"));
283   PetscCall(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s"));
284   PetscCall(PetscBagSetFromOptions(ctx->bag));
285   {
286     PetscViewer       viewer;
287     PetscViewerFormat format;
288     PetscBool         flg;
289 
290     PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
291     if (flg) {
292       PetscCall(PetscViewerPushFormat(viewer, format));
293       PetscCall(PetscBagView(ctx->bag, viewer));
294       PetscCall(PetscViewerFlush(viewer));
295       PetscCall(PetscViewerPopFormat(viewer));
296       PetscCall(PetscViewerDestroy(&viewer));
297     }
298   }
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
303 {
304   PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
305   PetscDS          ds;
306   DMLabel          label;
307   const PetscInt   id = 1;
308 
309   PetscFunctionBeginUser;
310   PetscCall(DMGetDS(dm, &ds));
311   switch (user->sol) {
312     case SOL_QUADRATIC:
313       PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u));
314       exactFuncs[0] = quadratic_u;
315       exactFuncs[1] = quadratic_p;
316       break;
317     case SOL_TRIG:
318       PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
319       exactFuncs[0] = trig_u;
320       exactFuncs[1] = trig_p;
321       break;
322     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
323   }
324   PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL));
325   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu));
326   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL));
327   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
328   PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu));
329   PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL));
330 
331   PetscCall(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user));
332   PetscCall(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user));
333 
334   PetscCall(DMGetLabel(dm, "marker", &label));
335   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL));
336 
337   /* Make constant values available to pointwise functions */
338   {
339     Parameter  *param;
340     PetscScalar constants[1];
341 
342     PetscCall(PetscBagGetData(user->bag, (void **) &param));
343     constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
344     PetscCall(PetscDSSetConstants(ds, 1, constants));
345   }
346   PetscFunctionReturn(0);
347 }
348 
349 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
350 {
351   PetscInt c;
352   for (c = 0; c < Nc; ++c) u[c] = 0.0;
353   return 0;
354 }
355 static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
356 {
357   PetscInt c;
358   for (c = 0; c < Nc; ++c) u[c] = 1.0;
359   return 0;
360 }
361 
362 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
363 {
364   Vec              vec;
365   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one};
366 
367   PetscFunctionBeginUser;
368   PetscCheck(origField == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField);
369   funcs[field] = one;
370   {
371     PetscDS ds;
372     PetscCall(DMGetDS(dm, &ds));
373     PetscCall(PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view"));
374   }
375   PetscCall(DMCreateGlobalVector(dm, &vec));
376   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
377   PetscCall(VecNormalize(vec, NULL));
378   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace));
379   PetscCall(VecDestroy(&vec));
380   /* New style for field null spaces */
381   {
382     PetscObject  pressure;
383     MatNullSpace nullspacePres;
384 
385     PetscCall(DMGetField(dm, field, NULL, &pressure));
386     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
387     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres));
388     PetscCall(MatNullSpaceDestroy(&nullspacePres));
389   }
390   PetscFunctionReturn(0);
391 }
392 
393 static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
394 {
395   DM              cdm = dm;
396   PetscQuadrature q   = NULL;
397   PetscBool       simplex;
398   PetscInt        dim, Nf = 2, f, Nc[2];
399   const char     *name[2]   = {"velocity", "pressure"};
400   const char     *prefix[2] = {"vel_",     "pres_"};
401 
402   PetscFunctionBegin;
403   PetscCall(DMGetDimension(dm, &dim));
404   PetscCall(DMPlexIsSimplex(dm, &simplex));
405   Nc[0] = dim;
406   Nc[1] = 1;
407   for (f = 0; f < Nf; ++f) {
408     PetscFE fe;
409 
410     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe));
411     PetscCall(PetscObjectSetName((PetscObject) fe, name[f]));
412     if (!q) PetscCall(PetscFEGetQuadrature(fe, &q));
413     PetscCall(PetscFESetQuadrature(fe, q));
414     PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe));
415     PetscCall(PetscFEDestroy(&fe));
416   }
417   PetscCall(DMCreateDS(dm));
418   PetscCall((*setupEqn)(dm, user));
419   while (cdm) {
420     PetscCall(DMCopyDisc(dm, cdm));
421     PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace));
422     PetscCall(DMGetCoarseDM(cdm, &cdm));
423   }
424   PetscFunctionReturn(0);
425 }
426 
427 int main(int argc, char **argv)
428 {
429   SNES           snes;
430   DM             dm;
431   Vec            u;
432   AppCtx         user;
433 
434   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
435   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
436   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
437   PetscCall(SNESCreate(PetscObjectComm((PetscObject) dm), &snes));
438   PetscCall(SNESSetDM(snes, dm));
439   PetscCall(DMSetApplicationContext(dm, &user));
440 
441   PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
442   PetscCall(SetupProblem(dm, SetupEqn, &user));
443   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
444 
445   PetscCall(DMCreateGlobalVector(dm, &u));
446   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
447   PetscCall(SNESSetFromOptions(snes));
448   PetscCall(DMSNESCheckFromOptions(snes, u));
449   PetscCall(PetscObjectSetName((PetscObject) u, "Solution"));
450   {
451     Mat          J;
452     MatNullSpace sp;
453 
454     PetscCall(SNESSetUp(snes));
455     PetscCall(CreatePressureNullSpace(dm, 1, 1, &sp));
456     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
457     PetscCall(MatSetNullSpace(J, sp));
458     PetscCall(MatNullSpaceDestroy(&sp));
459     PetscCall(PetscObjectSetName((PetscObject) J, "Jacobian"));
460     PetscCall(MatViewFromOptions(J, NULL, "-J_view"));
461   }
462   PetscCall(SNESSolve(snes, NULL, u));
463 
464   PetscCall(VecDestroy(&u));
465   PetscCall(SNESDestroy(&snes));
466   PetscCall(DMDestroy(&dm));
467   PetscCall(PetscBagDestroy(&user.bag));
468   PetscCall(PetscFinalize());
469   return 0;
470 }
471 /*TEST
472 
473   test:
474     suffix: 2d_p2_p1_check
475     requires: triangle
476     args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
477 
478   test:
479     suffix: 2d_p2_p1_check_parallel
480     nsize: {{2 3 5}}
481     requires: triangle
482     args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
483 
484   test:
485     suffix: 3d_p2_p1_check
486     requires: ctetgen
487     args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
488 
489   test:
490     suffix: 3d_p2_p1_check_parallel
491     nsize: {{2 3 5}}
492     requires: ctetgen
493     args: -sol quadratic -dm_refine 2 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
494 
495   test:
496     suffix: 2d_p2_p1_conv
497     requires: triangle
498     # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
499     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
500       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
501       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
502         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
503 
504   test:
505     suffix: 2d_p2_p1_conv_gamg
506     requires: triangle
507     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2  \
508       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
509         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
510 
511   test:
512     suffix: 3d_p2_p1_conv
513     requires: ctetgen !single
514     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
515     args: -sol trig -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
516       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
517       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
518         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
519 
520   test:
521     suffix: 2d_q2_q1_check
522     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
523 
524   test:
525     suffix: 3d_q2_q1_check
526     args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
527 
528   test:
529     suffix: 2d_q2_q1_conv
530     # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
531     args: -sol trig -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -ksp_error_if_not_converged \
532       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
533       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
534         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
535 
536   test:
537     suffix: 3d_q2_q1_conv
538     requires: !single
539     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
540     args: -sol trig -dm_plex_simplex 0 -dm_plex_dim 3 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
541       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
542       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
543         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
544 
545   test:
546     suffix: 2d_p3_p2_check
547     requires: triangle
548     args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
549 
550   test:
551     suffix: 3d_p3_p2_check
552     requires: ctetgen !single
553     args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
554 
555   test:
556     suffix: 2d_p3_p2_conv
557     requires: triangle
558     # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
559     args: -sol trig -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
560       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
561       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
562         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
563 
564   test:
565     suffix: 3d_p3_p2_conv
566     requires: ctetgen long_runtime
567     # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
568     args: -sol trig -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 \
569       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
570       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
571         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
572 
573   test:
574     suffix: 2d_q1_p0_conv
575     requires: !single
576     # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
577     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
578       -ksp_atol 1e-10 -petscds_jac_pre 0 \
579       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
580         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd
581 
582   test:
583     suffix: 3d_q1_p0_conv
584     requires: !single
585     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
586     args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 1 \
587       -ksp_atol 1e-10 -petscds_jac_pre 0 \
588       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
589         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd
590 
591   # Stokes preconditioners
592   #   Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
593   test:
594     suffix: 2d_p2_p1_block_diagonal
595     requires: triangle
596     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
597       -snes_error_if_not_converged \
598       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
599       -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
600   #   Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
601   test:
602     suffix: 2d_p2_p1_block_triangular
603     requires: triangle
604     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
605       -snes_error_if_not_converged \
606       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
607       -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
608   #   Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
609   test:
610     suffix: 2d_p2_p1_schur_diagonal
611     requires: triangle
612     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
613       -snes_error_if_not_converged \
614       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
615       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
616         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
617   #   Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
618   test:
619     suffix: 2d_p2_p1_schur_upper
620     requires: triangle
621     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
622       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
623       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
624         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
625   #   Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
626   test:
627     suffix: 2d_p2_p1_schur_lower
628     requires: triangle
629     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
630       -snes_error_if_not_converged \
631       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
632       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
633         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
634   #   Full Schur complement \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} \begin{pmatrix} I & A^{-1} B \\ 0 & I \end{pmatrix}
635   test:
636     suffix: 2d_p2_p1_schur_full
637     requires: triangle
638     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
639       -snes_error_if_not_converged \
640       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
641       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \
642         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
643   #   Full Schur + Velocity GMG
644   test:
645     suffix: 2d_p2_p1_gmg_vcycle
646     requires: triangle
647     args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
648       -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \
649       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
650         -fieldsplit_velocity_pc_type mg -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_pc_gamg_esteig_ksp_max_it 10 -fieldsplit_pressure_mg_levels_pc_type sor -fieldsplit_pressure_mg_coarse_pc_type svd
651   #   SIMPLE \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T diag(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & diag(A)^{-1} B \\ 0 & I \end{pmatrix}
652   test:
653     suffix: 2d_p2_p1_simple
654     requires: triangle
655     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
656       -snes_error_if_not_converged \
657       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
658       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
659         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
660         -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi
661   #   FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
662   test:
663     suffix: 2d_p2_p1_fetidp
664     requires: triangle mumps
665     nsize: 5
666     args: -sol quadratic -dm_refine 2 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
667       -snes_error_if_not_converged \
668       -ksp_type fetidp -ksp_rtol 1.0e-8 \
669       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
670         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \
671         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly
672   test:
673     suffix: 2d_q2_q1_fetidp
674     requires: mumps
675     nsize: 5
676     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
677       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
678       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
679         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \
680         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly
681   test:
682     suffix: 3d_p2_p1_fetidp
683     requires: ctetgen mumps suitesparse
684     nsize: 5
685     args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
686       -snes_error_if_not_converged \
687       -ksp_type fetidp -ksp_rtol 1.0e-9  \
688       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
689         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none \
690         -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
691         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
692         -fetidp_bddelta_pc_factor_mat_solver_type umfpack -fetidp_fieldsplit_lag_ksp_type preonly -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 100000 \
693         -fetidp_bddelta_pc_factor_mat_ordering_type external \
694         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
695         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
696   test:
697     suffix: 3d_q2_q1_fetidp
698     requires: suitesparse
699     nsize: 5
700     args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
701       -snes_error_if_not_converged \
702       -ksp_type fetidp -ksp_rtol 1.0e-8 \
703       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
704         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 2000 -fetidp_fieldsplit_p_pc_type none \
705         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
706         -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
707         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
708         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
709   #   BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
710   test:
711     suffix: 2d_p2_p1_bddc
712     nsize: 2
713     requires: triangle !single
714     args: -sol quadratic -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
715       -snes_error_if_not_converged \
716       -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
717         -pc_type bddc -pc_bddc_corner_selection -pc_bddc_dirichlet_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_coarse_redundant_pc_type svd
718   #   Vanka
719   test:
720     suffix: 2d_q1_p0_vanka
721     requires: double !complex
722     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
723       -snes_rtol 1.0e-4 \
724       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
725       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
726         -sub_ksp_type preonly -sub_pc_type lu
727   test:
728     suffix: 2d_q1_p0_vanka_denseinv
729     requires: double !complex
730     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
731       -snes_rtol 1.0e-4 \
732       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
733       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
734         -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
735   #   Vanka smoother
736   test:
737     suffix: 2d_q1_p0_gmg_vanka
738     requires: double !complex
739     args: -sol quadratic -dm_plex_simplex 0 -dm_refine_hierarchy 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
740       -snes_rtol 1.0e-4 \
741       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
742       -pc_type mg \
743         -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
744         -mg_levels_pc_type patch -mg_levels_pc_patch_partition_of_unity 0 -mg_levels_pc_patch_construct_codim 0 -mg_levels_pc_patch_construct_type vanka \
745           -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
746         -mg_coarse_pc_type svd
747 
748 TEST*/
749