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