xref: /petsc/src/snes/tutorials/ex62.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
378b0e23d0SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_TRIG, SOL_UNKNOWN} SolType;
388b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0};
39c4762a1bSJed Brown 
40c4762a1bSJed Brown typedef struct {
418b0e23d0SMatthew G. Knepley   PetscScalar mu; /* dynamic shear viscosity */
428b0e23d0SMatthew G. Knepley } Parameter;
438b0e23d0SMatthew G. Knepley 
448b0e23d0SMatthew G. Knepley typedef struct {
458b0e23d0SMatthew G. Knepley   PetscBag bag; /* Problem parameters */
468b0e23d0SMatthew G. Knepley   SolType  sol; /* MMS solution */
47c4762a1bSJed Brown } AppCtx;
48c4762a1bSJed Brown 
498b0e23d0SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
50c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
51c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
52c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
53c4762a1bSJed Brown {
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) {
59c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
608b0e23d0SMatthew G. Knepley       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
61c4762a1bSJed Brown     }
628b0e23d0SMatthew G. Knepley     f1[c*dim+c] -= u[uOff[1]];
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
668b0e23d0SMatthew G. Knepley static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
67c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
68c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
69c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
70c4762a1bSJed Brown {
71c4762a1bSJed Brown   PetscInt d;
728b0e23d0SMatthew G. Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d];
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
758b0e23d0SMatthew G. Knepley static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
76c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
77c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
78c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
79c4762a1bSJed Brown {
80c4762a1bSJed Brown   PetscInt d;
818b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
848b0e23d0SMatthew G. Knepley static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
85c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
86c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
87c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   PetscInt d;
908b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
938b0e23d0SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
97c4762a1bSJed Brown {
988b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
998b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1]-uOff[0];
1008b0e23d0SMatthew G. Knepley   PetscInt        c, d;
101c4762a1bSJed Brown 
1028b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
103c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
1048b0e23d0SMatthew G. Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */
1058b0e23d0SMatthew G. Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */
106c4762a1bSJed Brown     }
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
1108b0e23d0SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1118b0e23d0SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1128b0e23d0SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1138b0e23d0SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1148b0e23d0SMatthew G. Knepley {
1158b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1168b0e23d0SMatthew G. Knepley 
1178b0e23d0SMatthew G. Knepley   g0[0] = 1.0/mu;
1188b0e23d0SMatthew G. Knepley }
1198b0e23d0SMatthew G. Knepley 
1208b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution
1218b0e23d0SMatthew G. Knepley    2D:
122c4762a1bSJed Brown 
123c4762a1bSJed Brown      u = x^2 + y^2
1248b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1258b0e23d0SMatthew G. Knepley      p = x + y - 1
1268b0e23d0SMatthew G. Knepley      f = <1 - 4 mu, 1 - 4 mu>
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    so that
129c4762a1bSJed Brown 
1308b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 4x  4x \
1318b0e23d0SMatthew G. Knepley                                   \ 4x -4x /
1328b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
1338b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2x - 2x = 0
1348b0e23d0SMatthew G. Knepley 
1358b0e23d0SMatthew G. Knepley    3D:
1368b0e23d0SMatthew G. Knepley 
1378b0e23d0SMatthew G. Knepley      u = 2 x^2 + y^2 + z^2
1388b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1398b0e23d0SMatthew G. Knepley      w = 2 x^2 - 2xz
1408b0e23d0SMatthew G. Knepley      p = x + y + z - 3/2
1418b0e23d0SMatthew G. Knepley      f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
1428b0e23d0SMatthew G. Knepley 
1438b0e23d0SMatthew G. Knepley    so that
1448b0e23d0SMatthew G. Knepley 
1458b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 8x  4x  4x \
1468b0e23d0SMatthew G. Knepley                                   | 4x -4x  0  |
1478b0e23d0SMatthew G. Knepley                                   \ 4x  0  -4x /
1488b0e23d0SMatthew 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
1498b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 4x - 2x - 2x = 0
150c4762a1bSJed Brown */
1518b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
152c4762a1bSJed Brown {
1538b0e23d0SMatthew G. Knepley   PetscInt c;
1548b0e23d0SMatthew G. Knepley 
1558b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSqr(x[0]);
1568b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
1578b0e23d0SMatthew G. Knepley     u[0] += PetscSqr(x[c]);
1588b0e23d0SMatthew G. Knepley     u[c]  = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c];
1598b0e23d0SMatthew G. Knepley   }
160c4762a1bSJed Brown   return 0;
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
1638b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
164c4762a1bSJed Brown {
1658b0e23d0SMatthew G. Knepley   PetscInt d;
1668b0e23d0SMatthew G. Knepley 
1678b0e23d0SMatthew G. Knepley   u[0] = -0.5*dim;
1688b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d];
169c4762a1bSJed Brown   return 0;
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
1728b0e23d0SMatthew G. Knepley static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
173c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
174c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1758b0e23d0SMatthew G. Knepley                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
176c4762a1bSJed Brown {
1778b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1788b0e23d0SMatthew G. Knepley   PetscInt        d;
1798b0e23d0SMatthew G. Knepley 
1808b0e23d0SMatthew G. Knepley   f0[0] = (dim-1)*4.0*mu - 1.0;
1818b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0;
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
1848b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution
1858b0e23d0SMatthew G. Knepley    2D:
1868b0e23d0SMatthew G. Knepley 
1878b0e23d0SMatthew G. Knepley      u = sin(pi x) + sin(pi y)
1888b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
1898b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y)
1908b0e23d0SMatthew 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>
1918b0e23d0SMatthew G. Knepley 
1928b0e23d0SMatthew G. Knepley    so that
1938b0e23d0SMatthew G. Knepley 
1948b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = /        2pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y \
1958b0e23d0SMatthew G. Knepley                                   \ pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)          /
1968b0e23d0SMatthew 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
1978b0e23d0SMatthew G. Knepley      \nabla \cdot u             = pi cos(pi x) - pi cos(pi x) = 0
1988b0e23d0SMatthew G. Knepley 
1998b0e23d0SMatthew G. Knepley    3D:
2008b0e23d0SMatthew G. Knepley 
2018b0e23d0SMatthew G. Knepley      u = 2 sin(pi x) + sin(pi y) + sin(pi z)
2028b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
2038b0e23d0SMatthew G. Knepley      w = -pi cos(pi x) z
2048b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
2058b0e23d0SMatthew 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>
2068b0e23d0SMatthew G. Knepley 
2078b0e23d0SMatthew G. Knepley    so that
2088b0e23d0SMatthew G. Knepley 
2098b0e23d0SMatthew 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 \
2108b0e23d0SMatthew G. Knepley                                   | pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)                        0                  |
2118b0e23d0SMatthew G. Knepley                                   \ pi cos(pi z) + pi^2 sin(pi x) z               0                         -2pi cos(pi x)            /
2128b0e23d0SMatthew 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
2138b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
2148b0e23d0SMatthew G. Knepley */
2158b0e23d0SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
216c4762a1bSJed Brown {
2178b0e23d0SMatthew G. Knepley   PetscInt c;
2188b0e23d0SMatthew G. Knepley 
2198b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]);
2208b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
2218b0e23d0SMatthew G. Knepley     u[0] += PetscSinReal(PETSC_PI*x[c]);
2228b0e23d0SMatthew G. Knepley     u[c]  = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c];
2238b0e23d0SMatthew G. Knepley   }
2248b0e23d0SMatthew G. Knepley   return 0;
2258b0e23d0SMatthew G. Knepley }
2268b0e23d0SMatthew G. Knepley 
2278b0e23d0SMatthew G. Knepley static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
2288b0e23d0SMatthew G. Knepley {
2298b0e23d0SMatthew G. Knepley   PetscInt d;
2308b0e23d0SMatthew G. Knepley 
2318b0e23d0SMatthew G. Knepley   for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]);
2328b0e23d0SMatthew G. Knepley   return 0;
2338b0e23d0SMatthew G. Knepley }
2348b0e23d0SMatthew G. Knepley 
2358b0e23d0SMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2368b0e23d0SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2378b0e23d0SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
2388b0e23d0SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
2398b0e23d0SMatthew G. Knepley {
2408b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
2418b0e23d0SMatthew G. Knepley   PetscInt        d;
2428b0e23d0SMatthew G. Knepley 
2438b0e23d0SMatthew 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]);
2448b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) {
2458b0e23d0SMatthew G. Knepley     f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]);
2468b0e23d0SMatthew 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];
2478b0e23d0SMatthew G. Knepley   }
2488b0e23d0SMatthew G. Knepley }
2498b0e23d0SMatthew G. Knepley 
2508b0e23d0SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
2518b0e23d0SMatthew G. Knepley {
2528b0e23d0SMatthew G. Knepley   PetscInt       sol;
253c4762a1bSJed Brown   PetscErrorCode ierr;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBeginUser;
2568b0e23d0SMatthew G. Knepley   options->sol = SOL_QUADRATIC;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr);
2598b0e23d0SMatthew G. Knepley   sol  = options->sol;
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, (sizeof(SolTypes)/sizeof(SolTypes[0]))-3, SolTypes[options->sol], &sol, NULL));
2618b0e23d0SMatthew G. Knepley   options->sol = (SolType) sol;
2621e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
263c4762a1bSJed Brown   PetscFunctionReturn(0);
264c4762a1bSJed Brown }
265c4762a1bSJed Brown 
2668b0e23d0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
267c4762a1bSJed Brown {
268c4762a1bSJed Brown   PetscFunctionBeginUser;
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
273c4762a1bSJed Brown   PetscFunctionReturn(0);
274c4762a1bSJed Brown }
275c4762a1bSJed Brown 
2768b0e23d0SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
277c4762a1bSJed Brown {
2788b0e23d0SMatthew G. Knepley   Parameter     *p;
2798b0e23d0SMatthew G. Knepley 
2808b0e23d0SMatthew G. Knepley   PetscFunctionBeginUser;
2818b0e23d0SMatthew G. Knepley   /* setup PETSc parameter bag */
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(ctx->bag, (void **) &p));
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagSetName(ctx->bag, "par", "Stokes Parameters"));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s"));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagSetFromOptions(ctx->bag));
2878b0e23d0SMatthew G. Knepley   {
2888b0e23d0SMatthew G. Knepley     PetscViewer       viewer;
2898b0e23d0SMatthew G. Knepley     PetscViewerFormat format;
2908b0e23d0SMatthew G. Knepley     PetscBool         flg;
2918b0e23d0SMatthew G. Knepley 
2925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
2938b0e23d0SMatthew G. Knepley     if (flg) {
2945f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPushFormat(viewer, format));
2955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBagView(ctx->bag, viewer));
2965f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
2975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPopFormat(viewer));
2985f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&viewer));
2998b0e23d0SMatthew G. Knepley     }
3008b0e23d0SMatthew G. Knepley   }
3018b0e23d0SMatthew G. Knepley   PetscFunctionReturn(0);
3028b0e23d0SMatthew G. Knepley }
3038b0e23d0SMatthew G. Knepley 
3048b0e23d0SMatthew G. Knepley static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
3058b0e23d0SMatthew G. Knepley {
3068b0e23d0SMatthew G. Knepley   PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
3078b0e23d0SMatthew G. Knepley   PetscDS          ds;
30845480ffeSMatthew G. Knepley   DMLabel          label;
309c4762a1bSJed Brown   const PetscInt   id = 1;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   PetscFunctionBeginUser;
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
3138b0e23d0SMatthew G. Knepley   switch (user->sol) {
314c4762a1bSJed Brown     case SOL_QUADRATIC:
3155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u));
3168b0e23d0SMatthew G. Knepley       exactFuncs[0] = quadratic_u;
3178b0e23d0SMatthew G. Knepley       exactFuncs[1] = quadratic_p;
318c4762a1bSJed Brown       break;
319c4762a1bSJed Brown     case SOL_TRIG:
3205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
3218b0e23d0SMatthew G. Knepley       exactFuncs[0] = trig_u;
3228b0e23d0SMatthew G. Knepley       exactFuncs[1] = trig_p;
323c4762a1bSJed Brown       break;
32498921bdaSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
325c4762a1bSJed Brown   }
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetResidual(ds, 1, f0_p, NULL));
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu));
3285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL));
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
3305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu));
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL));
332c4762a1bSJed Brown 
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user));
3345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user));
335c4762a1bSJed Brown 
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL));
3388b0e23d0SMatthew G. Knepley 
33947bb1945SPatrick Sanan   /* Make constant values available to pointwise functions */
340c4762a1bSJed Brown   {
3418b0e23d0SMatthew G. Knepley     Parameter  *param;
3428b0e23d0SMatthew G. Knepley     PetscScalar constants[1];
343c4762a1bSJed Brown 
3445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
3458b0e23d0SMatthew G. Knepley     constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
3465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetConstants(ds, 1, constants));
347c4762a1bSJed Brown   }
348c4762a1bSJed Brown   PetscFunctionReturn(0);
349c4762a1bSJed Brown }
350c4762a1bSJed Brown 
3518b0e23d0SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
3528b0e23d0SMatthew G. Knepley {
3538b0e23d0SMatthew G. Knepley   PetscInt c;
3548b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
3558b0e23d0SMatthew G. Knepley   return 0;
3568b0e23d0SMatthew G. Knepley }
3578b0e23d0SMatthew G. Knepley static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
3588b0e23d0SMatthew G. Knepley {
3598b0e23d0SMatthew G. Knepley   PetscInt c;
3608b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 1.0;
3618b0e23d0SMatthew G. Knepley   return 0;
3628b0e23d0SMatthew G. Knepley }
3638b0e23d0SMatthew G. Knepley 
3648b0e23d0SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
365c4762a1bSJed Brown {
366c4762a1bSJed Brown   Vec              vec;
367478db826SMatthew G. Knepley   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one};
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   PetscFunctionBeginUser;
3702c71b3e2SJacob Faibussowitsch   PetscCheckFalse(origField != 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %D should be 1 for pressure", origField);
3718b0e23d0SMatthew G. Knepley   funcs[field] = one;
3728b0e23d0SMatthew G. Knepley   {
3738b0e23d0SMatthew G. Knepley     PetscDS ds;
3745f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDS(dm, &ds));
3755f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view"));
3768b0e23d0SMatthew G. Knepley   }
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &vec));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
3795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNormalize(vec, NULL));
3805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace));
3815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vec));
382c4762a1bSJed Brown   /* New style for field null spaces */
383c4762a1bSJed Brown   {
384c4762a1bSJed Brown     PetscObject  pressure;
385c4762a1bSJed Brown     MatNullSpace nullspacePres;
386c4762a1bSJed Brown 
3875f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetField(dm, field, NULL, &pressure));
3885f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
3895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres));
3905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceDestroy(&nullspacePres));
391c4762a1bSJed Brown   }
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
3958b0e23d0SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
396c4762a1bSJed Brown {
3978b0e23d0SMatthew G. Knepley   DM              cdm = dm;
3988b0e23d0SMatthew G. Knepley   PetscQuadrature q   = NULL;
3998b0e23d0SMatthew G. Knepley   PetscBool       simplex;
40030602db0SMatthew G. Knepley   PetscInt        dim, Nf = 2, f, Nc[2];
4018b0e23d0SMatthew G. Knepley   const char     *name[2]   = {"velocity", "pressure"};
4028b0e23d0SMatthew G. Knepley   const char     *prefix[2] = {"vel_",     "pres_"};
403c4762a1bSJed Brown 
4048b0e23d0SMatthew G. Knepley   PetscFunctionBegin;
4055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
4065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(dm, &simplex));
4078b0e23d0SMatthew G. Knepley   Nc[0] = dim;
4088b0e23d0SMatthew G. Knepley   Nc[1] = 1;
4098b0e23d0SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4108b0e23d0SMatthew G. Knepley     PetscFE fe;
4118b0e23d0SMatthew G. Knepley 
4125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe));
4135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) fe, name[f]));
4145f80ce2aSJacob Faibussowitsch     if (!q) CHKERRQ(PetscFEGetQuadrature(fe, &q));
4155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFESetQuadrature(fe, q));
4165f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetField(dm, f, NULL, (PetscObject) fe));
4175f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEDestroy(&fe));
418c4762a1bSJed Brown   }
4195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
4205f80ce2aSJacob Faibussowitsch   CHKERRQ((*setupEqn)(dm, user));
4218b0e23d0SMatthew G. Knepley   while (cdm) {
4225f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, cdm));
4235f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace));
4245f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
425c4762a1bSJed Brown   }
426c4762a1bSJed Brown   PetscFunctionReturn(0);
427c4762a1bSJed Brown }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown int main(int argc, char **argv)
430c4762a1bSJed Brown {
4318b0e23d0SMatthew G. Knepley   SNES           snes;
4328b0e23d0SMatthew G. Knepley   DM             dm;
4338b0e23d0SMatthew G. Knepley   Vec            u;
4348b0e23d0SMatthew G. Knepley   AppCtx         user;
435c4762a1bSJed Brown 
436*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
4375f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
4385f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4395f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PetscObjectComm((PetscObject) dm), &snes));
4405f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes, dm));
4415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(dm, &user));
442c4762a1bSJed Brown 
4435f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupParameters(PETSC_COMM_WORLD, &user));
4445f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupProblem(dm, SetupEqn, &user));
4455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateClosureIndex(dm, NULL));
446c4762a1bSJed Brown 
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
4485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
4495f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
4505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESCheckFromOptions(snes, u));
4515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "Solution"));
4528b0e23d0SMatthew G. Knepley   {
4538b0e23d0SMatthew G. Knepley     Mat          J;
4548b0e23d0SMatthew G. Knepley     MatNullSpace sp;
455c4762a1bSJed Brown 
4565f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetUp(snes));
4575f80ce2aSJacob Faibussowitsch     CHKERRQ(CreatePressureNullSpace(dm, 1, 1, &sp));
4585f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
4595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetNullSpace(J, sp));
4605f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceDestroy(&sp));
4615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) J, "Jacobian"));
4625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(J, NULL, "-J_view"));
463c4762a1bSJed Brown   }
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, NULL, u));
465c4762a1bSJed Brown 
4665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
4675f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
4685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
4695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagDestroy(&user.bag));
470*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
471*b122ec5aSJacob Faibussowitsch   return 0;
472c4762a1bSJed Brown }
473c4762a1bSJed Brown /*TEST
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   test:
4768b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check
477c4762a1bSJed Brown     requires: triangle
4788b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4798b0e23d0SMatthew G. Knepley 
480c4762a1bSJed Brown   test:
4818b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check_parallel
4828b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
483c4762a1bSJed Brown     requires: triangle
484e600fa54SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4858b0e23d0SMatthew G. Knepley 
486c4762a1bSJed Brown   test:
4878b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check
488c4762a1bSJed Brown     requires: ctetgen
48930602db0SMatthew 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
4908b0e23d0SMatthew G. Knepley 
491c4762a1bSJed Brown   test:
4928b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check_parallel
4938b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
494c4762a1bSJed Brown     requires: ctetgen
495e600fa54SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4968b0e23d0SMatthew G. Knepley 
497c4762a1bSJed Brown   test:
4988b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv
4998b0e23d0SMatthew G. Knepley     requires: triangle
5008b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
5018b0e23d0SMatthew 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 \
5028b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5038b0e23d0SMatthew 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 \
5048b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5058b0e23d0SMatthew G. Knepley 
506c4762a1bSJed Brown   test:
5078b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv_gamg
5088b0e23d0SMatthew G. Knepley     requires: triangle
50982894d03SBarry Smith     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2  \
5108b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
5118b0e23d0SMatthew 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
5128b0e23d0SMatthew G. Knepley 
513c4762a1bSJed Brown   test:
5148b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_conv
5158b0e23d0SMatthew G. Knepley     requires: ctetgen !single
5168b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
51730602db0SMatthew 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 \
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 
522c4762a1bSJed Brown   test:
5238b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_check
52430602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5258b0e23d0SMatthew G. Knepley 
526c4762a1bSJed Brown   test:
5278b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_check
52830602db0SMatthew 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
5298b0e23d0SMatthew G. Knepley 
530c4762a1bSJed Brown   test:
5318b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_conv
5328b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
53330602db0SMatthew 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 \
5348b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5358b0e23d0SMatthew 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 \
5368b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5378b0e23d0SMatthew G. Knepley 
538c4762a1bSJed Brown   test:
5398b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_conv
540c4762a1bSJed Brown     requires: !single
5418b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
54230602db0SMatthew 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 \
5438b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5448b0e23d0SMatthew 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 \
5458b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5468b0e23d0SMatthew G. Knepley 
547c4762a1bSJed Brown   test:
5488b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_check
5498b0e23d0SMatthew G. Knepley     requires: triangle
5508b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
5518b0e23d0SMatthew G. Knepley 
5528b0e23d0SMatthew G. Knepley   test:
5538b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_check
5548b0e23d0SMatthew G. Knepley     requires: ctetgen !single
55530602db0SMatthew 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
5568b0e23d0SMatthew G. Knepley 
5578b0e23d0SMatthew G. Knepley   test:
5588b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_conv
5598b0e23d0SMatthew G. Knepley     requires: triangle
5608b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
5618b0e23d0SMatthew 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 \
5628b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5638b0e23d0SMatthew 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 \
5648b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5658b0e23d0SMatthew G. Knepley 
5668b0e23d0SMatthew G. Knepley   test:
5678b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_conv
5688b0e23d0SMatthew G. Knepley     requires: ctetgen long_runtime
5698b0e23d0SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
57030602db0SMatthew 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 \
5718b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5728b0e23d0SMatthew 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 \
5738b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5748b0e23d0SMatthew G. Knepley 
5758b0e23d0SMatthew G. Knepley   test:
5768b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_conv
577c4762a1bSJed Brown     requires: !single
5788b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
57930602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
58082894d03SBarry Smith       -ksp_atol 1e-10 -petscds_jac_pre 0 \
5818b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
5828b0e23d0SMatthew G. Knepley         -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
5838b0e23d0SMatthew G. Knepley 
584c4762a1bSJed Brown   test:
5858b0e23d0SMatthew G. Knepley     suffix: 3d_q1_p0_conv
5868b0e23d0SMatthew G. Knepley     requires: !single
5878b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
58830602db0SMatthew 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 \
58982894d03SBarry Smith       -ksp_atol 1e-10 -petscds_jac_pre 0 \
5908b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
5918b0e23d0SMatthew G. Knepley         -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
5928b0e23d0SMatthew G. Knepley 
5938b0e23d0SMatthew G. Knepley   # Stokes preconditioners
594c4762a1bSJed Brown   #   Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
595c4762a1bSJed Brown   test:
5968b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_diagonal
5978b0e23d0SMatthew G. Knepley     requires: triangle
5988b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
5998b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6008b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
6018b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
602c4762a1bSJed Brown   #   Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
603c4762a1bSJed Brown   test:
6048b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_triangular
6058b0e23d0SMatthew G. Knepley     requires: triangle
6068b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6078b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6088b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6098b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
610c4762a1bSJed Brown   #   Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
611c4762a1bSJed Brown   test:
6128b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_diagonal
6138b0e23d0SMatthew G. Knepley     requires: triangle
6148b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6158b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6168b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6178b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
6188b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
619c4762a1bSJed Brown   #   Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
620c4762a1bSJed Brown   test:
6218b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_upper
6228b0e23d0SMatthew G. Knepley     requires: triangle
6238b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
6248b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6258b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
6268b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
627c4762a1bSJed Brown   #   Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
628c4762a1bSJed Brown   test:
6298b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_lower
6308b0e23d0SMatthew G. Knepley     requires: triangle
6318b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6328b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6338b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6348b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
6358b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
636c4762a1bSJed 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}
637c4762a1bSJed Brown   test:
6388b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_full
6398b0e23d0SMatthew G. Knepley     requires: triangle
6408b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6418b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6428b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6438b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \
6448b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
6458b0e23d0SMatthew G. Knepley   #   Full Schur + Velocity GMG
6468b0e23d0SMatthew G. Knepley   test:
6478b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_gmg_vcycle
6488b0e23d0SMatthew G. Knepley     requires: triangle
6498b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
65082894d03SBarry Smith       -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \
6518b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
65273f7197eSJed 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
653c4762a1bSJed 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}
654c4762a1bSJed Brown   test:
6558b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_simple
656c4762a1bSJed Brown     requires: triangle
6578b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6588b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6598b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6608b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
6618b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
6628b0e23d0SMatthew 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
663c4762a1bSJed Brown   #   FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
664c4762a1bSJed Brown   test:
6658b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_fetidp
666c4762a1bSJed Brown     requires: triangle mumps
667c4762a1bSJed Brown     nsize: 5
668e600fa54SMatthew 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 \
6698b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
67082894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-8 \
6718b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6728b0e23d0SMatthew 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 \
6738b0e23d0SMatthew 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
674c4762a1bSJed Brown   test:
6758b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_fetidp
6768b0e23d0SMatthew G. Knepley     requires: mumps
677c4762a1bSJed Brown     nsize: 5
678e600fa54SMatthew 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 \
6798b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6808b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6818b0e23d0SMatthew 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 \
6828b0e23d0SMatthew 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
683c4762a1bSJed Brown   test:
6848b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_fetidp
6858b0e23d0SMatthew G. Knepley     requires: ctetgen mumps suitesparse
686c4762a1bSJed Brown     nsize: 5
687e600fa54SMatthew 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 \
6888b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
68982894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-9  \
6908b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6918b0e23d0SMatthew 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 \
6928b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
6938b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
6948b0e23d0SMatthew 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 \
6958b0e23d0SMatthew G. Knepley         -fetidp_bddelta_pc_factor_mat_ordering_type external \
6968b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
6978b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
698c4762a1bSJed Brown   test:
6998b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_fetidp
700c4762a1bSJed Brown     requires: suitesparse
701c4762a1bSJed Brown     nsize: 5
702e600fa54SMatthew 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 \
7038b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
70482894d03SBarry Smith       -ksp_type fetidp -ksp_rtol 1.0e-8 \
7058b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
7068b0e23d0SMatthew 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 \
7078b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
7088b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
7098b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
7108b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
7118b0e23d0SMatthew G. Knepley   #   BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
712c4762a1bSJed Brown   test:
7138b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_bddc
7148b0e23d0SMatthew G. Knepley     nsize: 2
715c4762a1bSJed Brown     requires: triangle !single
716e600fa54SMatthew 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 \
7178b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
7188b0e23d0SMatthew G. Knepley       -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
7198b0e23d0SMatthew 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
7208b0e23d0SMatthew G. Knepley   #   Vanka
721c4762a1bSJed Brown   test:
7228b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka
723c4762a1bSJed Brown     requires: double !complex
72430602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7258b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7268b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
727c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
728c4762a1bSJed Brown         -sub_ksp_type preonly -sub_pc_type lu
729c4762a1bSJed Brown   test:
7308b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka_denseinv
731c4762a1bSJed Brown     requires: double !complex
73230602db0SMatthew G. Knepley     args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7338b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7348b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
735c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
736c4762a1bSJed Brown         -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
737c4762a1bSJed Brown   #   Vanka smoother
738c4762a1bSJed Brown   test:
7398b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_gmg_vanka
7408b0e23d0SMatthew G. Knepley     requires: double !complex
74130602db0SMatthew 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 \
7428b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7438b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
7448b0e23d0SMatthew G. Knepley       -pc_type mg \
7458b0e23d0SMatthew G. Knepley         -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
746c4762a1bSJed 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 \
747c4762a1bSJed Brown           -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
748c4762a1bSJed Brown         -mg_coarse_pc_type svd
749c4762a1bSJed Brown 
750c4762a1bSJed Brown TEST*/
751