xref: /petsc/src/snes/tutorials/ex62.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
18b0e23d0SMatthew G. Knepley static char help[] = "Stokes Problem discretized with finite elements,\n\
28b0e23d0SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
58b0e23d0SMatthew G. Knepley For the isoviscous Stokes problem, which we discretize using the finite
68b0e23d0SMatthew G. Knepley element method on an unstructured mesh, the weak form equations are
7c4762a1bSJed Brown 
88b0e23d0SMatthew G. Knepley   < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0
98b0e23d0SMatthew G. Knepley   < q, -\nabla\cdot u >                                                   = 0
10c4762a1bSJed Brown 
11c4762a1bSJed Brown Viewing:
12c4762a1bSJed Brown 
13c4762a1bSJed Brown To produce nice output, use
14c4762a1bSJed Brown 
158b0e23d0SMatthew G. Knepley   -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
16c4762a1bSJed Brown 
17c4762a1bSJed Brown You can get a LaTeX view of the mesh, with point numbering using
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0
20c4762a1bSJed Brown 
21c4762a1bSJed Brown The data layout can be viewed using
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   -dm_petscsection_view
24c4762a1bSJed Brown 
25c4762a1bSJed Brown Lots of information about the FEM assembly can be printed using
26c4762a1bSJed Brown 
278b0e23d0SMatthew G. Knepley   -dm_plex_print_fem 3
28c4762a1bSJed Brown */
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #include <petscdmplex.h>
31c4762a1bSJed Brown #include <petscsnes.h>
32c4762a1bSJed Brown #include <petscds.h>
338b0e23d0SMatthew G. Knepley #include <petscbag.h>
34c4762a1bSJed Brown 
358b0e23d0SMatthew G. Knepley // TODO: Plot residual by fields after each smoother iterate
36c4762a1bSJed Brown 
3711486bccSBarry Smith typedef enum {
3811486bccSBarry Smith   SOL_QUADRATIC,
3911486bccSBarry Smith   SOL_TRIG,
4011486bccSBarry Smith   SOL_UNKNOWN
4111486bccSBarry Smith } SolType;
428b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0};
43c4762a1bSJed Brown 
44*9371c9d4SSatish Balay typedef struct {
458b0e23d0SMatthew G. Knepley   PetscScalar mu; /* dynamic shear viscosity */
468b0e23d0SMatthew G. Knepley } Parameter;
478b0e23d0SMatthew G. Knepley 
48*9371c9d4SSatish Balay typedef struct {
498b0e23d0SMatthew G. Knepley   PetscBag bag; /* Problem parameters */
508b0e23d0SMatthew G. Knepley   SolType  sol; /* MMS solution */
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown 
53*9371c9d4SSatish Balay static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
548b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
558b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1] - uOff[0];
568b0e23d0SMatthew G. Knepley   PetscInt        c, d;
57c4762a1bSJed Brown 
588b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
5911486bccSBarry Smith     for (d = 0; d < dim; ++d) { f1[c * dim + d] = mu * (u_x[c * dim + d] + u_x[d * dim + c]); }
608b0e23d0SMatthew G. Knepley     f1[c * dim + c] -= u[uOff[1]];
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64*9371c9d4SSatish Balay static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
65c4762a1bSJed Brown   PetscInt d;
668b0e23d0SMatthew G. Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d * dim + d];
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69*9371c9d4SSatish Balay static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
70c4762a1bSJed Brown   PetscInt d;
718b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0; /* < q, -\nabla\cdot u > */
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74*9371c9d4SSatish Balay static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) {
75c4762a1bSJed Brown   PetscInt d;
768b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* -< \nabla\cdot v, p > */
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79*9371c9d4SSatish Balay static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
808b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
818b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1] - uOff[0];
828b0e23d0SMatthew G. Knepley   PetscInt        c, d;
83c4762a1bSJed Brown 
848b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
85c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
868b0e23d0SMatthew G. Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += mu; /* < \nabla v, \nabla u > */
878b0e23d0SMatthew G. Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += mu; /* < \nabla v, {\nabla u}^T > */
88c4762a1bSJed Brown     }
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92*9371c9d4SSatish Balay static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
938b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
948b0e23d0SMatthew G. Knepley 
958b0e23d0SMatthew G. Knepley   g0[0] = 1.0 / mu;
968b0e23d0SMatthew G. Knepley }
978b0e23d0SMatthew G. Knepley 
988b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution
998b0e23d0SMatthew G. Knepley    2D:
100c4762a1bSJed Brown 
101c4762a1bSJed Brown      u = x^2 + y^2
1028b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1038b0e23d0SMatthew G. Knepley      p = x + y - 1
1048b0e23d0SMatthew G. Knepley      f = <1 - 4 mu, 1 - 4 mu>
105c4762a1bSJed Brown 
106c4762a1bSJed Brown    so that
107c4762a1bSJed Brown 
1088b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 4x  4x \
1098b0e23d0SMatthew G. Knepley                                   \ 4x -4x /
1108b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
1118b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2x - 2x = 0
1128b0e23d0SMatthew G. Knepley 
1138b0e23d0SMatthew G. Knepley    3D:
1148b0e23d0SMatthew G. Knepley 
1158b0e23d0SMatthew G. Knepley      u = 2 x^2 + y^2 + z^2
1168b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1178b0e23d0SMatthew G. Knepley      w = 2 x^2 - 2xz
1188b0e23d0SMatthew G. Knepley      p = x + y + z - 3/2
1198b0e23d0SMatthew G. Knepley      f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
1208b0e23d0SMatthew G. Knepley 
1218b0e23d0SMatthew G. Knepley    so that
1228b0e23d0SMatthew G. Knepley 
1238b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 8x  4x  4x \
1248b0e23d0SMatthew G. Knepley                                   | 4x -4x  0  |
1258b0e23d0SMatthew G. Knepley                                   \ 4x  0  -4x /
1268b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <8, 4, 4> - <1, 1, 1> + <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> = 0
1278b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 4x - 2x - 2x = 0
128c4762a1bSJed Brown */
129*9371c9d4SSatish Balay static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
1308b0e23d0SMatthew G. Knepley   PetscInt c;
1318b0e23d0SMatthew G. Knepley 
1328b0e23d0SMatthew G. Knepley   u[0] = (dim - 1) * PetscSqr(x[0]);
1338b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
1348b0e23d0SMatthew G. Knepley     u[0] += PetscSqr(x[c]);
1358b0e23d0SMatthew G. Knepley     u[c] = 2.0 * PetscSqr(x[0]) - 2.0 * x[0] * x[c];
1368b0e23d0SMatthew G. Knepley   }
137c4762a1bSJed Brown   return 0;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140*9371c9d4SSatish Balay static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
1418b0e23d0SMatthew G. Knepley   PetscInt d;
1428b0e23d0SMatthew G. Knepley 
1438b0e23d0SMatthew G. Knepley   u[0] = -0.5 * dim;
1448b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d];
145c4762a1bSJed Brown   return 0;
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148*9371c9d4SSatish Balay static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
1498b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1508b0e23d0SMatthew G. Knepley   PetscInt        d;
1518b0e23d0SMatthew G. Knepley 
1528b0e23d0SMatthew G. Knepley   f0[0] = (dim - 1) * 4.0 * mu - 1.0;
1538b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) f0[d] = 4.0 * mu - 1.0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
1568b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution
1578b0e23d0SMatthew G. Knepley    2D:
1588b0e23d0SMatthew G. Knepley 
1598b0e23d0SMatthew G. Knepley      u = sin(pi x) + sin(pi y)
1608b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
1618b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y)
1628b0e23d0SMatthew G. Knepley      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>
1638b0e23d0SMatthew G. Knepley 
1648b0e23d0SMatthew G. Knepley    so that
1658b0e23d0SMatthew G. Knepley 
1668b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = /        2pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y \
1678b0e23d0SMatthew G. Knepley                                   \ pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)          /
1688b0e23d0SMatthew G. Knepley      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
1698b0e23d0SMatthew G. Knepley      \nabla \cdot u             = pi cos(pi x) - pi cos(pi x) = 0
1708b0e23d0SMatthew G. Knepley 
1718b0e23d0SMatthew G. Knepley    3D:
1728b0e23d0SMatthew G. Knepley 
1738b0e23d0SMatthew G. Knepley      u = 2 sin(pi x) + sin(pi y) + sin(pi z)
1748b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
1758b0e23d0SMatthew G. Knepley      w = -pi cos(pi x) z
1768b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
1778b0e23d0SMatthew G. Knepley      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>
1788b0e23d0SMatthew G. Knepley 
1798b0e23d0SMatthew G. Knepley    so that
1808b0e23d0SMatthew G. Knepley 
1818b0e23d0SMatthew G. Knepley      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 \
1828b0e23d0SMatthew G. Knepley                                   | pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)                        0                  |
1838b0e23d0SMatthew G. Knepley                                   \ pi cos(pi z) + pi^2 sin(pi x) z               0                         -2pi cos(pi x)            /
1848b0e23d0SMatthew G. Knepley      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
1858b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
1868b0e23d0SMatthew G. Knepley */
187*9371c9d4SSatish Balay static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
1888b0e23d0SMatthew G. Knepley   PetscInt c;
1898b0e23d0SMatthew G. Knepley 
1908b0e23d0SMatthew G. Knepley   u[0] = (dim - 1) * PetscSinReal(PETSC_PI * x[0]);
1918b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
1928b0e23d0SMatthew G. Knepley     u[0] += PetscSinReal(PETSC_PI * x[c]);
1938b0e23d0SMatthew G. Knepley     u[c] = -PETSC_PI * PetscCosReal(PETSC_PI * x[0]) * x[c];
1948b0e23d0SMatthew G. Knepley   }
1958b0e23d0SMatthew G. Knepley   return 0;
1968b0e23d0SMatthew G. Knepley }
1978b0e23d0SMatthew G. Knepley 
198*9371c9d4SSatish Balay static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
1998b0e23d0SMatthew G. Knepley   PetscInt d;
2008b0e23d0SMatthew G. Knepley 
2018b0e23d0SMatthew G. Knepley   for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0 * PETSC_PI * x[d]);
2028b0e23d0SMatthew G. Knepley   return 0;
2038b0e23d0SMatthew G. Knepley }
2048b0e23d0SMatthew G. Knepley 
205*9371c9d4SSatish Balay static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
2068b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
2078b0e23d0SMatthew G. Knepley   PetscInt        d;
2088b0e23d0SMatthew G. Knepley 
2098b0e23d0SMatthew G. Knepley   f0[0] = -2.0 * PETSC_PI * PetscCosReal(2.0 * PETSC_PI * x[0]) - (dim - 1) * mu * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x[0]);
2108b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) {
2118b0e23d0SMatthew G. Knepley     f0[0] -= mu * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x[d]);
2128b0e23d0SMatthew G. Knepley     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];
2138b0e23d0SMatthew G. Knepley   }
2148b0e23d0SMatthew G. Knepley }
2158b0e23d0SMatthew G. Knepley 
216*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
2178b0e23d0SMatthew G. Knepley   PetscInt sol;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   PetscFunctionBeginUser;
2208b0e23d0SMatthew G. Knepley   options->sol = SOL_QUADRATIC;
221d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
2228b0e23d0SMatthew G. Knepley   sol = options->sol;
223dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, PETSC_STATIC_ARRAY_LENGTH(SolTypes) - 3, SolTypes[options->sol], &sol, NULL));
2248b0e23d0SMatthew G. Knepley   options->sol = (SolType)sol;
225d0609cedSBarry Smith   PetscOptionsEnd();
226c4762a1bSJed Brown   PetscFunctionReturn(0);
227c4762a1bSJed Brown }
228c4762a1bSJed Brown 
229*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
230c4762a1bSJed Brown   PetscFunctionBeginUser;
2319566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2329566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
2339566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
2349566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
235c4762a1bSJed Brown   PetscFunctionReturn(0);
236c4762a1bSJed Brown }
237c4762a1bSJed Brown 
238*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) {
2398b0e23d0SMatthew G. Knepley   Parameter *p;
2408b0e23d0SMatthew G. Knepley 
2418b0e23d0SMatthew G. Knepley   PetscFunctionBeginUser;
2428b0e23d0SMatthew G. Knepley   /* setup PETSc parameter bag */
2439566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag));
2449566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(ctx->bag, (void **)&p));
2459566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(ctx->bag, "par", "Stokes Parameters"));
2469566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s"));
2479566063dSJacob Faibussowitsch   PetscCall(PetscBagSetFromOptions(ctx->bag));
2488b0e23d0SMatthew G. Knepley   {
2498b0e23d0SMatthew G. Knepley     PetscViewer       viewer;
2508b0e23d0SMatthew G. Knepley     PetscViewerFormat format;
2518b0e23d0SMatthew G. Knepley     PetscBool         flg;
2528b0e23d0SMatthew G. Knepley 
2539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
2548b0e23d0SMatthew G. Knepley     if (flg) {
2559566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, format));
2569566063dSJacob Faibussowitsch       PetscCall(PetscBagView(ctx->bag, viewer));
2579566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2589566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
2599566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
2608b0e23d0SMatthew G. Knepley     }
2618b0e23d0SMatthew G. Knepley   }
2628b0e23d0SMatthew G. Knepley   PetscFunctionReturn(0);
2638b0e23d0SMatthew G. Knepley }
2648b0e23d0SMatthew G. Knepley 
265*9371c9d4SSatish Balay static PetscErrorCode SetupEqn(DM dm, AppCtx *user) {
2668b0e23d0SMatthew G. Knepley   PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
2678b0e23d0SMatthew G. Knepley   PetscDS        ds;
26845480ffeSMatthew G. Knepley   DMLabel        label;
269c4762a1bSJed Brown   const PetscInt id = 1;
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   PetscFunctionBeginUser;
2729566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2738b0e23d0SMatthew G. Knepley   switch (user->sol) {
274c4762a1bSJed Brown   case SOL_QUADRATIC:
2759566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u));
2768b0e23d0SMatthew G. Knepley     exactFuncs[0] = quadratic_u;
2778b0e23d0SMatthew G. Knepley     exactFuncs[1] = quadratic_p;
278c4762a1bSJed Brown     break;
279c4762a1bSJed Brown   case SOL_TRIG:
2809566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
2818b0e23d0SMatthew G. Knepley     exactFuncs[0] = trig_u;
2828b0e23d0SMatthew G. Knepley     exactFuncs[1] = trig_p;
283c4762a1bSJed Brown     break;
28463a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
285c4762a1bSJed Brown   }
2869566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL));
2879566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
2889566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
2899566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
2909566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu));
2919566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL));
292c4762a1bSJed Brown 
2939566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user));
2949566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user));
295c4762a1bSJed Brown 
2969566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
2979566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, user, NULL));
2988b0e23d0SMatthew G. Knepley 
29947bb1945SPatrick Sanan   /* Make constant values available to pointwise functions */
300c4762a1bSJed Brown   {
3018b0e23d0SMatthew G. Knepley     Parameter  *param;
3028b0e23d0SMatthew G. Knepley     PetscScalar constants[1];
303c4762a1bSJed Brown 
3049566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **)&param));
3058b0e23d0SMatthew G. Knepley     constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
3069566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 1, constants));
307c4762a1bSJed Brown   }
308c4762a1bSJed Brown   PetscFunctionReturn(0);
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
311*9371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
3128b0e23d0SMatthew G. Knepley   PetscInt c;
3138b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
3148b0e23d0SMatthew G. Knepley   return 0;
3158b0e23d0SMatthew G. Knepley }
316*9371c9d4SSatish Balay static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
3178b0e23d0SMatthew G. Knepley   PetscInt c;
3188b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 1.0;
3198b0e23d0SMatthew G. Knepley   return 0;
3208b0e23d0SMatthew G. Knepley }
3218b0e23d0SMatthew G. Knepley 
322*9371c9d4SSatish Balay static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) {
323c4762a1bSJed Brown   Vec vec;
324478db826SMatthew G. Knepley   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {zero, one};
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   PetscFunctionBeginUser;
32763a3b9bcSJacob Faibussowitsch   PetscCheck(origField == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField);
3288b0e23d0SMatthew G. Knepley   funcs[field] = one;
3298b0e23d0SMatthew G. Knepley   {
3308b0e23d0SMatthew G. Knepley     PetscDS ds;
3319566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
3329566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)ds, NULL, "-ds_view"));
3338b0e23d0SMatthew G. Knepley   }
3349566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &vec));
3359566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
3369566063dSJacob Faibussowitsch   PetscCall(VecNormalize(vec, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace));
3389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec));
339c4762a1bSJed Brown   /* New style for field null spaces */
340c4762a1bSJed Brown   {
341c4762a1bSJed Brown     PetscObject  pressure;
342c4762a1bSJed Brown     MatNullSpace nullspacePres;
343c4762a1bSJed Brown 
3449566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, field, NULL, &pressure));
3459566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
3469566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres));
3479566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullspacePres));
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown   PetscFunctionReturn(0);
350c4762a1bSJed Brown }
351c4762a1bSJed Brown 
352*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user) {
3538b0e23d0SMatthew G. Knepley   DM              cdm = dm;
3548b0e23d0SMatthew G. Knepley   PetscQuadrature q   = NULL;
3558b0e23d0SMatthew G. Knepley   PetscBool       simplex;
35630602db0SMatthew G. Knepley   PetscInt        dim, Nf = 2, f, Nc[2];
3578b0e23d0SMatthew G. Knepley   const char     *name[2]   = {"velocity", "pressure"};
3588b0e23d0SMatthew G. Knepley   const char     *prefix[2] = {"vel_", "pres_"};
359c4762a1bSJed Brown 
3608b0e23d0SMatthew G. Knepley   PetscFunctionBegin;
3619566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3629566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
3638b0e23d0SMatthew G. Knepley   Nc[0] = dim;
3648b0e23d0SMatthew G. Knepley   Nc[1] = 1;
3658b0e23d0SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
3668b0e23d0SMatthew G. Knepley     PetscFE fe;
3678b0e23d0SMatthew G. Knepley 
3689566063dSJacob Faibussowitsch     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe));
3699566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fe, name[f]));
3709566063dSJacob Faibussowitsch     if (!q) PetscCall(PetscFEGetQuadrature(fe, &q));
3719566063dSJacob Faibussowitsch     PetscCall(PetscFESetQuadrature(fe, q));
3729566063dSJacob Faibussowitsch     PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
3739566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
374c4762a1bSJed Brown   }
3759566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
3769566063dSJacob Faibussowitsch   PetscCall((*setupEqn)(dm, user));
3778b0e23d0SMatthew G. Knepley   while (cdm) {
3789566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
3799566063dSJacob Faibussowitsch     PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace));
3809566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
381c4762a1bSJed Brown   }
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
385*9371c9d4SSatish Balay int main(int argc, char **argv) {
3868b0e23d0SMatthew G. Knepley   SNES   snes;
3878b0e23d0SMatthew G. Knepley   DM     dm;
3888b0e23d0SMatthew G. Knepley   Vec    u;
3898b0e23d0SMatthew G. Knepley   AppCtx user;
390c4762a1bSJed Brown 
391327415f7SBarry Smith   PetscFunctionBeginUser;
3929566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3939566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3949566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3959566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
3969566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
3979566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &user));
398c4762a1bSJed Brown 
3999566063dSJacob Faibussowitsch   PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
4009566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, SetupEqn, &user));
4019566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
402c4762a1bSJed Brown 
4039566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
4049566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
4059566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
4069566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
4079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Solution"));
4088b0e23d0SMatthew G. Knepley   {
4098b0e23d0SMatthew G. Knepley     Mat          J;
4108b0e23d0SMatthew G. Knepley     MatNullSpace sp;
411c4762a1bSJed Brown 
4129566063dSJacob Faibussowitsch     PetscCall(SNESSetUp(snes));
4139566063dSJacob Faibussowitsch     PetscCall(CreatePressureNullSpace(dm, 1, 1, &sp));
4149566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
4159566063dSJacob Faibussowitsch     PetscCall(MatSetNullSpace(J, sp));
4169566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&sp));
4179566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
4189566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(J, NULL, "-J_view"));
419c4762a1bSJed Brown   }
4209566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
421c4762a1bSJed Brown 
4229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
4239566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
4249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4259566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
4269566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
427b122ec5aSJacob Faibussowitsch   return 0;
428c4762a1bSJed Brown }
429c4762a1bSJed Brown /*TEST
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   test:
4328b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check
433c4762a1bSJed Brown     requires: triangle
4348b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4358b0e23d0SMatthew G. Knepley 
436c4762a1bSJed Brown   test:
4378b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check_parallel
4388b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
439c4762a1bSJed Brown     requires: triangle
440e600fa54SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4418b0e23d0SMatthew G. Knepley 
442c4762a1bSJed Brown   test:
4438b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check
444c4762a1bSJed Brown     requires: ctetgen
44530602db0SMatthew G. Knepley     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
4468b0e23d0SMatthew G. Knepley 
447c4762a1bSJed Brown   test:
4488b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check_parallel
4498b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
450c4762a1bSJed Brown     requires: ctetgen
45119fa4260SStefano Zampini     args: -sol quadratic -dm_refine 0 -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
4528b0e23d0SMatthew G. Knepley 
453c4762a1bSJed Brown   test:
4548b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv
4558b0e23d0SMatthew G. Knepley     requires: triangle
4568b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
4578b0e23d0SMatthew G. Knepley     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
4588b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
4598b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
4608b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
4618b0e23d0SMatthew G. Knepley 
462c4762a1bSJed Brown   test:
4638b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv_gamg
4648b0e23d0SMatthew G. Knepley     requires: triangle
46582894d03SBarry Smith     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2  \
4668b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
4678b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
4688b0e23d0SMatthew G. Knepley 
469c4762a1bSJed Brown   test:
4708b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_conv
4718b0e23d0SMatthew G. Knepley     requires: ctetgen !single
4728b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
47330602db0SMatthew G. Knepley     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 \
4748b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
4758b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
4768b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
4778b0e23d0SMatthew G. Knepley 
478c4762a1bSJed Brown   test:
4798b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_check
48030602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4818b0e23d0SMatthew G. Knepley 
482c4762a1bSJed Brown   test:
4838b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_check
48430602db0SMatthew G. Knepley     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
4858b0e23d0SMatthew G. Knepley 
486c4762a1bSJed Brown   test:
4878b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_conv
4888b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
48930602db0SMatthew G. Knepley     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 \
4908b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
4918b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
4928b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
4938b0e23d0SMatthew G. Knepley 
494c4762a1bSJed Brown   test:
4958b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_conv
496c4762a1bSJed Brown     requires: !single
4978b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
49830602db0SMatthew G. Knepley     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 \
4998b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5008b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
5018b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5028b0e23d0SMatthew G. Knepley 
503c4762a1bSJed Brown   test:
5048b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_check
5058b0e23d0SMatthew G. Knepley     requires: triangle
5068b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
5078b0e23d0SMatthew G. Knepley 
5088b0e23d0SMatthew G. Knepley   test:
5098b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_check
5108b0e23d0SMatthew G. Knepley     requires: ctetgen !single
51130602db0SMatthew G. Knepley     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
5128b0e23d0SMatthew G. Knepley 
5138b0e23d0SMatthew G. Knepley   test:
5148b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_conv
5158b0e23d0SMatthew G. Knepley     requires: triangle
5168b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
5178b0e23d0SMatthew G. Knepley     args: -sol trig -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
5188b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5198b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
5208b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5218b0e23d0SMatthew G. Knepley 
5228b0e23d0SMatthew G. Knepley   test:
5238b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_conv
5248b0e23d0SMatthew G. Knepley     requires: ctetgen long_runtime
5258b0e23d0SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
52630602db0SMatthew G. Knepley     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 \
5278b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5288b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
5298b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5308b0e23d0SMatthew G. Knepley 
5318b0e23d0SMatthew G. Knepley   test:
5328b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_conv
533c4762a1bSJed Brown     requires: !single
5348b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
53530602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
53682894d03SBarry Smith       -ksp_atol 1e-10 -petscds_jac_pre 0 \
5378b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
538bae903cbSmarkadams4         -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 -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0
5398b0e23d0SMatthew G. Knepley 
540c4762a1bSJed Brown   test:
5418b0e23d0SMatthew G. Knepley     suffix: 3d_q1_p0_conv
5428b0e23d0SMatthew G. Knepley     requires: !single
5438b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
54430602db0SMatthew G. Knepley     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 \
54582894d03SBarry Smith       -ksp_atol 1e-10 -petscds_jac_pre 0 \
5468b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
547bae903cbSmarkadams4         -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 -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0
5488b0e23d0SMatthew G. Knepley 
5498b0e23d0SMatthew G. Knepley   # Stokes preconditioners
550c4762a1bSJed Brown   #   Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
551c4762a1bSJed Brown   test:
5528b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_diagonal
5538b0e23d0SMatthew G. Knepley     requires: triangle
5548b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
5558b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
5568b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
5578b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
558c4762a1bSJed Brown   #   Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
559c4762a1bSJed Brown   test:
5608b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_triangular
5618b0e23d0SMatthew G. Knepley     requires: triangle
5628b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
5638b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
5648b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
5658b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
566c4762a1bSJed Brown   #   Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
567c4762a1bSJed Brown   test:
5688b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_diagonal
5698b0e23d0SMatthew G. Knepley     requires: triangle
5708b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
5718b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
5728b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
5738b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
5748b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
575c4762a1bSJed Brown   #   Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
576c4762a1bSJed Brown   test:
5778b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_upper
5788b0e23d0SMatthew G. Knepley     requires: triangle
5798b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
5808b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
5818b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
5828b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
583c4762a1bSJed Brown   #   Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
584c4762a1bSJed Brown   test:
5858b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_lower
5868b0e23d0SMatthew G. Knepley     requires: triangle
5878b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
5888b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
5898b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
5908b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
5918b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
592c4762a1bSJed Brown   #   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}
593c4762a1bSJed Brown   test:
5948b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_full
5958b0e23d0SMatthew G. Knepley     requires: triangle
5968b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
5978b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
5988b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
5998b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \
6008b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
6018b0e23d0SMatthew G. Knepley   #   Full Schur + Velocity GMG
6028b0e23d0SMatthew G. Knepley   test:
6038b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_gmg_vcycle
6048b0e23d0SMatthew G. Knepley     requires: triangle
6058b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
60682894d03SBarry Smith       -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \
6078b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
60873f7197eSJed Brown         -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
609c4762a1bSJed Brown   #   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}
610c4762a1bSJed Brown   test:
6118b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_simple
612c4762a1bSJed Brown     requires: triangle
6138b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6148b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6158b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6168b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
6178b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
6188b0e23d0SMatthew G. Knepley         -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi
619c4762a1bSJed Brown   #   FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
620c4762a1bSJed Brown   test:
6218b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_fetidp
622c4762a1bSJed Brown     requires: triangle mumps
623c4762a1bSJed Brown     nsize: 5
624e600fa54SMatthew G. Knepley     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 \
6258b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
62682894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-8 \
6278b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6288b0e23d0SMatthew G. Knepley         -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 \
6298b0e23d0SMatthew G. Knepley         -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
630c4762a1bSJed Brown   test:
6318b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_fetidp
6328b0e23d0SMatthew G. Knepley     requires: mumps
633c4762a1bSJed Brown     nsize: 5
634e600fa54SMatthew G. Knepley     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 \
6358b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6368b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6378b0e23d0SMatthew G. Knepley         -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 \
6388b0e23d0SMatthew G. Knepley         -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
639c4762a1bSJed Brown   test:
6408b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_fetidp
6418b0e23d0SMatthew G. Knepley     requires: ctetgen mumps suitesparse
642c4762a1bSJed Brown     nsize: 5
643e600fa54SMatthew G. Knepley     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 \
6448b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
64582894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-9  \
6468b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6478b0e23d0SMatthew G. Knepley         -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 \
6488b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
6498b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
6508b0e23d0SMatthew G. Knepley         -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 \
6518b0e23d0SMatthew G. Knepley         -fetidp_bddelta_pc_factor_mat_ordering_type external \
6528b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
6538b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
654c4762a1bSJed Brown   test:
6558b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_fetidp
656c4762a1bSJed Brown     requires: suitesparse
657c4762a1bSJed Brown     nsize: 5
658e600fa54SMatthew G. Knepley     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 \
6598b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
66082894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-8 \
6618b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6628b0e23d0SMatthew G. Knepley         -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 \
6638b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
6648b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
6658b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
6668b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
6678b0e23d0SMatthew G. Knepley   #   BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
668c4762a1bSJed Brown   test:
6698b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_bddc
6708b0e23d0SMatthew G. Knepley     nsize: 2
671c4762a1bSJed Brown     requires: triangle !single
672e600fa54SMatthew G. Knepley     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 \
6738b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6748b0e23d0SMatthew G. Knepley       -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6758b0e23d0SMatthew G. Knepley         -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
6768b0e23d0SMatthew G. Knepley   #   Vanka
677c4762a1bSJed Brown   test:
6788b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka
679c4762a1bSJed Brown     requires: double !complex
68030602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
6818b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
6828b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
683c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
684c4762a1bSJed Brown         -sub_ksp_type preonly -sub_pc_type lu
685c4762a1bSJed Brown   test:
6868b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka_denseinv
687c4762a1bSJed Brown     requires: double !complex
68830602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
6898b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
6908b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
691c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
692c4762a1bSJed Brown         -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
693c4762a1bSJed Brown   #   Vanka smoother
694c4762a1bSJed Brown   test:
6958b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_gmg_vanka
6968b0e23d0SMatthew G. Knepley     requires: double !complex
69730602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -dm_refine_hierarchy 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
6988b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
6998b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
7008b0e23d0SMatthew G. Knepley       -pc_type mg \
7018b0e23d0SMatthew G. Knepley         -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
702c4762a1bSJed Brown         -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 \
703c4762a1bSJed Brown           -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
704c4762a1bSJed Brown         -mg_coarse_pc_type svd
705c4762a1bSJed Brown 
706c4762a1bSJed Brown TEST*/
707