xref: /petsc/src/snes/tutorials/ex62.c (revision 478db826a40dca33d2f318a94ee2bc67dcb6bd71)
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 {
45c4762a1bSJed Brown   /* Domain and mesh definition */
468b0e23d0SMatthew G. Knepley   char     filename[PETSC_MAX_PATH_LEN];
47c4762a1bSJed Brown   /* Problem definition */
488b0e23d0SMatthew G. Knepley   PetscBag bag; /* Problem parameters */
498b0e23d0SMatthew G. Knepley   SolType  sol; /* MMS solution */
50c4762a1bSJed Brown } AppCtx;
51c4762a1bSJed Brown 
528b0e23d0SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
53c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
54c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
55c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
56c4762a1bSJed Brown {
578b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
588b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1]-uOff[0];
598b0e23d0SMatthew G. Knepley   PetscInt        c, d;
60c4762a1bSJed Brown 
618b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
62c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
638b0e23d0SMatthew G. Knepley       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
64c4762a1bSJed Brown     }
658b0e23d0SMatthew G. Knepley     f1[c*dim+c] -= u[uOff[1]];
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
698b0e23d0SMatthew G. Knepley static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
70c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
71c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
72c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
73c4762a1bSJed Brown {
74c4762a1bSJed Brown   PetscInt d;
758b0e23d0SMatthew G. Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d];
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
788b0e23d0SMatthew G. Knepley static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
79c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
80c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
81c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   PetscInt d;
848b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
878b0e23d0SMatthew G. Knepley static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
91c4762a1bSJed Brown {
92c4762a1bSJed Brown   PetscInt d;
938b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
968b0e23d0SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
98c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
99c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
100c4762a1bSJed Brown {
1018b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1028b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1]-uOff[0];
1038b0e23d0SMatthew G. Knepley   PetscInt        c, d;
104c4762a1bSJed Brown 
1058b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
106c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
1078b0e23d0SMatthew G. Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */
1088b0e23d0SMatthew G. Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */
109c4762a1bSJed Brown     }
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
1138b0e23d0SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1148b0e23d0SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1158b0e23d0SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1168b0e23d0SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1178b0e23d0SMatthew G. Knepley {
1188b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1198b0e23d0SMatthew G. Knepley 
1208b0e23d0SMatthew G. Knepley   g0[0] = 1.0/mu;
1218b0e23d0SMatthew G. Knepley }
1228b0e23d0SMatthew G. Knepley 
1238b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution
1248b0e23d0SMatthew G. Knepley    2D:
125c4762a1bSJed Brown 
126c4762a1bSJed Brown      u = x^2 + y^2
1278b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1288b0e23d0SMatthew G. Knepley      p = x + y - 1
1298b0e23d0SMatthew G. Knepley      f = <1 - 4 mu, 1 - 4 mu>
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    so that
132c4762a1bSJed Brown 
1338b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 4x  4x \
1348b0e23d0SMatthew G. Knepley                                   \ 4x -4x /
1358b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
1368b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2x - 2x = 0
1378b0e23d0SMatthew G. Knepley 
1388b0e23d0SMatthew G. Knepley    3D:
1398b0e23d0SMatthew G. Knepley 
1408b0e23d0SMatthew G. Knepley      u = 2 x^2 + y^2 + z^2
1418b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
1428b0e23d0SMatthew G. Knepley      w = 2 x^2 - 2xz
1438b0e23d0SMatthew G. Knepley      p = x + y + z - 3/2
1448b0e23d0SMatthew G. Knepley      f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
1458b0e23d0SMatthew G. Knepley 
1468b0e23d0SMatthew G. Knepley    so that
1478b0e23d0SMatthew G. Knepley 
1488b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 8x  4x  4x \
1498b0e23d0SMatthew G. Knepley                                   | 4x -4x  0  |
1508b0e23d0SMatthew G. Knepley                                   \ 4x  0  -4x /
1518b0e23d0SMatthew 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
1528b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 4x - 2x - 2x = 0
153c4762a1bSJed Brown */
1548b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
155c4762a1bSJed Brown {
1568b0e23d0SMatthew G. Knepley   PetscInt c;
1578b0e23d0SMatthew G. Knepley 
1588b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSqr(x[0]);
1598b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
1608b0e23d0SMatthew G. Knepley     u[0] += PetscSqr(x[c]);
1618b0e23d0SMatthew G. Knepley     u[c]  = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c];
1628b0e23d0SMatthew G. Knepley   }
163c4762a1bSJed Brown   return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
1668b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
167c4762a1bSJed Brown {
1688b0e23d0SMatthew G. Knepley   PetscInt d;
1698b0e23d0SMatthew G. Knepley 
1708b0e23d0SMatthew G. Knepley   u[0] = -0.5*dim;
1718b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d];
172c4762a1bSJed Brown   return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
1758b0e23d0SMatthew G. Knepley static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1788b0e23d0SMatthew G. Knepley                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
179c4762a1bSJed Brown {
1808b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
1818b0e23d0SMatthew G. Knepley   PetscInt        d;
1828b0e23d0SMatthew G. Knepley 
1838b0e23d0SMatthew G. Knepley   f0[0] = (dim-1)*4.0*mu - 1.0;
1848b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
1878b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution
1888b0e23d0SMatthew G. Knepley    2D:
1898b0e23d0SMatthew G. Knepley 
1908b0e23d0SMatthew G. Knepley      u = sin(pi x) + sin(pi y)
1918b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
1928b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y)
1938b0e23d0SMatthew 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>
1948b0e23d0SMatthew G. Knepley 
1958b0e23d0SMatthew G. Knepley    so that
1968b0e23d0SMatthew G. Knepley 
1978b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = /        2pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y \
1988b0e23d0SMatthew G. Knepley                                   \ pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)          /
1998b0e23d0SMatthew 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
2008b0e23d0SMatthew G. Knepley      \nabla \cdot u             = pi cos(pi x) - pi cos(pi x) = 0
2018b0e23d0SMatthew G. Knepley 
2028b0e23d0SMatthew G. Knepley    3D:
2038b0e23d0SMatthew G. Knepley 
2048b0e23d0SMatthew G. Knepley      u = 2 sin(pi x) + sin(pi y) + sin(pi z)
2058b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
2068b0e23d0SMatthew G. Knepley      w = -pi cos(pi x) z
2078b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
2088b0e23d0SMatthew 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>
2098b0e23d0SMatthew G. Knepley 
2108b0e23d0SMatthew G. Knepley    so that
2118b0e23d0SMatthew G. Knepley 
2128b0e23d0SMatthew 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 \
2138b0e23d0SMatthew G. Knepley                                   | pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)                        0                  |
2148b0e23d0SMatthew G. Knepley                                   \ pi cos(pi z) + pi^2 sin(pi x) z               0                         -2pi cos(pi x)            /
2158b0e23d0SMatthew 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
2168b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
2178b0e23d0SMatthew G. Knepley */
2188b0e23d0SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
219c4762a1bSJed Brown {
2208b0e23d0SMatthew G. Knepley   PetscInt c;
2218b0e23d0SMatthew G. Knepley 
2228b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]);
2238b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
2248b0e23d0SMatthew G. Knepley     u[0] += PetscSinReal(PETSC_PI*x[c]);
2258b0e23d0SMatthew G. Knepley     u[c]  = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c];
2268b0e23d0SMatthew G. Knepley   }
2278b0e23d0SMatthew G. Knepley   return 0;
2288b0e23d0SMatthew G. Knepley }
2298b0e23d0SMatthew G. Knepley 
2308b0e23d0SMatthew G. Knepley static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
2318b0e23d0SMatthew G. Knepley {
2328b0e23d0SMatthew G. Knepley   PetscInt d;
2338b0e23d0SMatthew G. Knepley 
2348b0e23d0SMatthew G. Knepley   for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]);
2358b0e23d0SMatthew G. Knepley   return 0;
2368b0e23d0SMatthew G. Knepley }
2378b0e23d0SMatthew G. Knepley 
2388b0e23d0SMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2398b0e23d0SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2408b0e23d0SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
2418b0e23d0SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
2428b0e23d0SMatthew G. Knepley {
2438b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
2448b0e23d0SMatthew G. Knepley   PetscInt        d;
2458b0e23d0SMatthew G. Knepley 
2468b0e23d0SMatthew 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]);
2478b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) {
2488b0e23d0SMatthew G. Knepley     f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]);
2498b0e23d0SMatthew 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];
2508b0e23d0SMatthew G. Knepley   }
2518b0e23d0SMatthew G. Knepley }
2528b0e23d0SMatthew G. Knepley 
2538b0e23d0SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
2548b0e23d0SMatthew G. Knepley {
2558b0e23d0SMatthew G. Knepley   PetscInt       sol;
256c4762a1bSJed Brown   PetscErrorCode ierr;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBeginUser;
2598b0e23d0SMatthew G. Knepley   options->filename[0] = '\0';
2608b0e23d0SMatthew G. Knepley   options->sol         = SOL_QUADRATIC;
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr);
2638b0e23d0SMatthew G. Knepley   sol  = options->sol;
26442a5c13dSPatrick Sanan   ierr = PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, (sizeof(SolTypes)/sizeof(SolTypes[0]))-3, SolTypes[options->sol], &sol, NULL);CHKERRQ(ierr);
2658b0e23d0SMatthew G. Knepley   options->sol = (SolType) sol;
266c4762a1bSJed Brown   ierr = PetscOptionsEnd();
267c4762a1bSJed Brown   PetscFunctionReturn(0);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
2708b0e23d0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
271c4762a1bSJed Brown {
2728b0e23d0SMatthew G. Knepley   size_t         len;
273c4762a1bSJed Brown   PetscErrorCode ierr;
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   PetscFunctionBeginUser;
2768b0e23d0SMatthew G. Knepley   ierr = PetscStrlen(user->filename, &len);CHKERRQ(ierr);
2778b0e23d0SMatthew G. Knepley   if (len) {ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr);}
2788b0e23d0SMatthew G. Knepley   else     {ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);}
279c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
281c4762a1bSJed Brown   PetscFunctionReturn(0);
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
2848b0e23d0SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
285c4762a1bSJed Brown {
2868b0e23d0SMatthew G. Knepley   Parameter     *p;
2878b0e23d0SMatthew G. Knepley   PetscErrorCode ierr;
2888b0e23d0SMatthew G. Knepley 
2898b0e23d0SMatthew G. Knepley   PetscFunctionBeginUser;
2908b0e23d0SMatthew G. Knepley   /* setup PETSc parameter bag */
2918b0e23d0SMatthew G. Knepley   ierr = PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag);CHKERRQ(ierr);
2928b0e23d0SMatthew G. Knepley   ierr = PetscBagGetData(ctx->bag, (void **) &p);CHKERRQ(ierr);
2938b0e23d0SMatthew G. Knepley   ierr = PetscBagSetName(ctx->bag, "par", "Stokes Parameters");CHKERRQ(ierr);
2948b0e23d0SMatthew G. Knepley   ierr = PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s");CHKERRQ(ierr);
2958b0e23d0SMatthew G. Knepley   ierr = PetscBagSetFromOptions(ctx->bag);CHKERRQ(ierr);
2968b0e23d0SMatthew G. Knepley   {
2978b0e23d0SMatthew G. Knepley     PetscViewer       viewer;
2988b0e23d0SMatthew G. Knepley     PetscViewerFormat format;
2998b0e23d0SMatthew G. Knepley     PetscBool         flg;
3008b0e23d0SMatthew G. Knepley 
3018b0e23d0SMatthew G. Knepley     ierr = PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg);CHKERRQ(ierr);
3028b0e23d0SMatthew G. Knepley     if (flg) {
3038b0e23d0SMatthew G. Knepley       ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
3048b0e23d0SMatthew G. Knepley       ierr = PetscBagView(ctx->bag, viewer);CHKERRQ(ierr);
3058b0e23d0SMatthew G. Knepley       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
3068b0e23d0SMatthew G. Knepley       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
3078b0e23d0SMatthew G. Knepley       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3088b0e23d0SMatthew G. Knepley     }
3098b0e23d0SMatthew G. Knepley   }
3108b0e23d0SMatthew G. Knepley   PetscFunctionReturn(0);
3118b0e23d0SMatthew G. Knepley }
3128b0e23d0SMatthew G. Knepley 
3138b0e23d0SMatthew G. Knepley static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
3148b0e23d0SMatthew G. Knepley {
3158b0e23d0SMatthew G. Knepley   PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
3168b0e23d0SMatthew G. Knepley   PetscDS          ds;
317c4762a1bSJed Brown   const PetscInt   id = 1;
318c4762a1bSJed Brown   PetscErrorCode   ierr;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   PetscFunctionBeginUser;
3218b0e23d0SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
3228b0e23d0SMatthew G. Knepley   switch (user->sol) {
323c4762a1bSJed Brown     case SOL_QUADRATIC:
3248b0e23d0SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u);CHKERRQ(ierr);
3258b0e23d0SMatthew G. Knepley       exactFuncs[0] = quadratic_u;
3268b0e23d0SMatthew G. Knepley       exactFuncs[1] = quadratic_p;
327c4762a1bSJed Brown       break;
328c4762a1bSJed Brown     case SOL_TRIG:
3298b0e23d0SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
3308b0e23d0SMatthew G. Knepley       exactFuncs[0] = trig_u;
3318b0e23d0SMatthew G. Knepley       exactFuncs[1] = trig_p;
332c4762a1bSJed Brown       break;
3338b0e23d0SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
334c4762a1bSJed Brown   }
3358b0e23d0SMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 1, f0_p, NULL);CHKERRQ(ierr);
3368b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu);CHKERRQ(ierr);
3378b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL);CHKERRQ(ierr);
3388b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL);CHKERRQ(ierr);
3398b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
3408b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL);CHKERRQ(ierr);
341c4762a1bSJed Brown 
3428b0e23d0SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 0, exactFuncs[0], user);CHKERRQ(ierr);
3438b0e23d0SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 1, exactFuncs[1], user);CHKERRQ(ierr);
344c4762a1bSJed Brown 
3458b0e23d0SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, user);CHKERRQ(ierr);
3468b0e23d0SMatthew G. Knepley 
34747bb1945SPatrick Sanan   /* Make constant values available to pointwise functions */
348c4762a1bSJed Brown   {
3498b0e23d0SMatthew G. Knepley     Parameter  *param;
3508b0e23d0SMatthew G. Knepley     PetscScalar constants[1];
351c4762a1bSJed Brown 
3528b0e23d0SMatthew G. Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
3538b0e23d0SMatthew G. Knepley     constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
3548b0e23d0SMatthew G. Knepley     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
355c4762a1bSJed Brown   }
356c4762a1bSJed Brown   PetscFunctionReturn(0);
357c4762a1bSJed Brown }
358c4762a1bSJed Brown 
3598b0e23d0SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
3608b0e23d0SMatthew G. Knepley {
3618b0e23d0SMatthew G. Knepley   PetscInt c;
3628b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
3638b0e23d0SMatthew G. Knepley   return 0;
3648b0e23d0SMatthew G. Knepley }
3658b0e23d0SMatthew G. Knepley static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
3668b0e23d0SMatthew G. Knepley {
3678b0e23d0SMatthew G. Knepley   PetscInt c;
3688b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 1.0;
3698b0e23d0SMatthew G. Knepley   return 0;
3708b0e23d0SMatthew G. Knepley }
3718b0e23d0SMatthew G. Knepley 
3728b0e23d0SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   Vec              vec;
375*478db826SMatthew G. Knepley   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one};
376c4762a1bSJed Brown   PetscErrorCode   ierr;
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   PetscFunctionBeginUser;
3798b0e23d0SMatthew G. Knepley   if (origField != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %D should be 1 for pressure", origField);
3808b0e23d0SMatthew G. Knepley   funcs[field] = one;
3818b0e23d0SMatthew G. Knepley   {
3828b0e23d0SMatthew G. Knepley     PetscDS ds;
3838b0e23d0SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
3848b0e23d0SMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view");CHKERRQ(ierr);
3858b0e23d0SMatthew G. Knepley   }
386c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
387c4762a1bSJed Brown   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
388c4762a1bSJed Brown   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
389c4762a1bSJed Brown   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace);CHKERRQ(ierr);
390c4762a1bSJed Brown   ierr = VecDestroy(&vec);CHKERRQ(ierr);
391c4762a1bSJed Brown   /* New style for field null spaces */
392c4762a1bSJed Brown   {
393c4762a1bSJed Brown     PetscObject  pressure;
394c4762a1bSJed Brown     MatNullSpace nullspacePres;
395c4762a1bSJed Brown 
3968b0e23d0SMatthew G. Knepley     ierr = DMGetField(dm, field, NULL, &pressure);CHKERRQ(ierr);
397c4762a1bSJed Brown     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
398c4762a1bSJed Brown     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
399c4762a1bSJed Brown     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
400c4762a1bSJed Brown   }
401c4762a1bSJed Brown   PetscFunctionReturn(0);
402c4762a1bSJed Brown }
403c4762a1bSJed Brown 
4048b0e23d0SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
405c4762a1bSJed Brown {
4068b0e23d0SMatthew G. Knepley   DM              cdm = dm;
4078b0e23d0SMatthew G. Knepley   PetscQuadrature q   = NULL;
4088b0e23d0SMatthew G. Knepley   DMPolytopeType  ct;
4098b0e23d0SMatthew G. Knepley   PetscBool       simplex;
4108b0e23d0SMatthew G. Knepley   PetscInt        dim, Nf = 2, f, Nc[2], cStart;
4118b0e23d0SMatthew G. Knepley   const char     *name[2]   = {"velocity", "pressure"};
4128b0e23d0SMatthew G. Knepley   const char     *prefix[2] = {"vel_",     "pres_"};
413c4762a1bSJed Brown   PetscErrorCode  ierr;
414c4762a1bSJed Brown 
4158b0e23d0SMatthew G. Knepley   PetscFunctionBegin;
4168b0e23d0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4178b0e23d0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
4188b0e23d0SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
4198b0e23d0SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
4208b0e23d0SMatthew G. Knepley   Nc[0] = dim;
4218b0e23d0SMatthew G. Knepley   Nc[1] = 1;
4228b0e23d0SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
4238b0e23d0SMatthew G. Knepley     PetscFE fe;
4248b0e23d0SMatthew G. Knepley 
4258b0e23d0SMatthew G. Knepley     ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe);CHKERRQ(ierr);
4268b0e23d0SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) fe, name[f]);CHKERRQ(ierr);
4278b0e23d0SMatthew G. Knepley     if (!q) {ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);}
4288b0e23d0SMatthew G. Knepley     ierr = PetscFESetQuadrature(fe, q);CHKERRQ(ierr);
4298b0e23d0SMatthew G. Knepley     ierr = DMSetField(dm, f, NULL, (PetscObject) fe);CHKERRQ(ierr);
4308b0e23d0SMatthew G. Knepley     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
431c4762a1bSJed Brown   }
4328b0e23d0SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
4338b0e23d0SMatthew G. Knepley   ierr = (*setupEqn)(dm, user);CHKERRQ(ierr);
4348b0e23d0SMatthew G. Knepley   while (cdm) {
4358b0e23d0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
4368b0e23d0SMatthew G. Knepley     ierr = DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
4378b0e23d0SMatthew G. Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
438c4762a1bSJed Brown   }
439c4762a1bSJed Brown   PetscFunctionReturn(0);
440c4762a1bSJed Brown }
441c4762a1bSJed Brown 
442c4762a1bSJed Brown int main(int argc, char **argv)
443c4762a1bSJed Brown {
4448b0e23d0SMatthew G. Knepley   SNES           snes;
4458b0e23d0SMatthew G. Knepley   DM             dm;
4468b0e23d0SMatthew G. Knepley   Vec            u;
4478b0e23d0SMatthew G. Knepley   AppCtx         user;
448c4762a1bSJed Brown   PetscErrorCode ierr;
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
451c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
452c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
4538b0e23d0SMatthew G. Knepley   ierr = SNESCreate(PetscObjectComm((PetscObject) dm), &snes);CHKERRQ(ierr);
454c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
455c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
456c4762a1bSJed Brown 
4578b0e23d0SMatthew G. Knepley   ierr = SetupParameters(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
4588b0e23d0SMatthew G. Knepley   ierr = SetupProblem(dm, SetupEqn, &user);CHKERRQ(ierr);
459c4762a1bSJed Brown   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
462c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
463c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
4648b0e23d0SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
465c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
4668b0e23d0SMatthew G. Knepley   {
4678b0e23d0SMatthew G. Knepley     Mat          J;
4688b0e23d0SMatthew G. Knepley     MatNullSpace sp;
469c4762a1bSJed Brown 
4708b0e23d0SMatthew G. Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
4718b0e23d0SMatthew G. Knepley     ierr = CreatePressureNullSpace(dm, 1, 1, &sp);CHKERRQ(ierr);
4728b0e23d0SMatthew G. Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
4738b0e23d0SMatthew G. Knepley     ierr = MatSetNullSpace(J, sp);CHKERRQ(ierr);
4748b0e23d0SMatthew G. Knepley     ierr = MatNullSpaceDestroy(&sp);CHKERRQ(ierr);
4758b0e23d0SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
4768b0e23d0SMatthew G. Knepley     ierr = MatViewFromOptions(J, NULL, "-J_view");CHKERRQ(ierr);
477c4762a1bSJed Brown   }
478c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
479c4762a1bSJed Brown 
480c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
481c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
482c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
4838b0e23d0SMatthew G. Knepley   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
484c4762a1bSJed Brown   ierr = PetscFinalize();
485c4762a1bSJed Brown   return ierr;
486c4762a1bSJed Brown }
487c4762a1bSJed Brown /*TEST
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   test:
4908b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check
491c4762a1bSJed Brown     requires: triangle
4928b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4938b0e23d0SMatthew G. Knepley 
494c4762a1bSJed Brown   test:
4958b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check_parallel
4968b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
497c4762a1bSJed Brown     requires: triangle
4988b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4998b0e23d0SMatthew G. Knepley 
500c4762a1bSJed Brown   test:
5018b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check
502c4762a1bSJed Brown     requires: ctetgen
5038b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5048b0e23d0SMatthew G. Knepley 
505c4762a1bSJed Brown   test:
5068b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check_parallel
5078b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
508c4762a1bSJed Brown     requires: ctetgen
5098b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5108b0e23d0SMatthew G. Knepley 
511c4762a1bSJed Brown   test:
5128b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv
5138b0e23d0SMatthew G. Knepley     requires: triangle
5148b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
5158b0e23d0SMatthew 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 \
5168b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5178b0e23d0SMatthew 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 \
5188b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5198b0e23d0SMatthew G. Knepley 
520c4762a1bSJed Brown   test:
5218b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv_gamg
5228b0e23d0SMatthew G. Knepley     requires: triangle
5238b0e23d0SMatthew 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 \
5248b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
5258b0e23d0SMatthew 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
5268b0e23d0SMatthew G. Knepley 
527c4762a1bSJed Brown   test:
5288b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_conv
5298b0e23d0SMatthew G. Knepley     requires: ctetgen !single
5308b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
5318b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
5328b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5338b0e23d0SMatthew 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 \
5348b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5358b0e23d0SMatthew G. Knepley 
536c4762a1bSJed Brown   test:
5378b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_check
5388b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5398b0e23d0SMatthew G. Knepley 
540c4762a1bSJed Brown   test:
5418b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_check
5428b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5438b0e23d0SMatthew G. Knepley 
544c4762a1bSJed Brown   test:
5458b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_conv
5468b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
5478b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -ksp_error_if_not_converged \
5488b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5498b0e23d0SMatthew 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 \
5508b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5518b0e23d0SMatthew G. Knepley 
552c4762a1bSJed Brown   test:
5538b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_conv
554c4762a1bSJed Brown     requires: !single
5558b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
5568b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
5578b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5588b0e23d0SMatthew 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 \
5598b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5608b0e23d0SMatthew G. Knepley 
561c4762a1bSJed Brown   test:
5628b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_check
5638b0e23d0SMatthew G. Knepley     requires: triangle
5648b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
5658b0e23d0SMatthew G. Knepley 
5668b0e23d0SMatthew G. Knepley   test:
5678b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_check
5688b0e23d0SMatthew G. Knepley     requires: ctetgen !single
5698b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
5708b0e23d0SMatthew G. Knepley 
5718b0e23d0SMatthew G. Knepley   test:
5728b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_conv
5738b0e23d0SMatthew G. Knepley     requires: triangle
5748b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
5758b0e23d0SMatthew 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 \
5768b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5778b0e23d0SMatthew 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 \
5788b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5798b0e23d0SMatthew G. Knepley 
5808b0e23d0SMatthew G. Knepley   test:
5818b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_conv
5828b0e23d0SMatthew G. Knepley     requires: ctetgen long_runtime
5838b0e23d0SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
5848b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 \
5858b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5868b0e23d0SMatthew 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 \
5878b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5888b0e23d0SMatthew G. Knepley 
5898b0e23d0SMatthew G. Knepley   test:
5908b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_conv
591c4762a1bSJed Brown     requires: !single
5928b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
5938b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
5948b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -petscds_jac_pre 0 \
5958b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
5968b0e23d0SMatthew 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
5978b0e23d0SMatthew G. Knepley 
598c4762a1bSJed Brown   test:
5998b0e23d0SMatthew G. Knepley     suffix: 3d_q1_p0_conv
6008b0e23d0SMatthew G. Knepley     requires: !single
6018b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
6028b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 1 \
6038b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -petscds_jac_pre 0 \
6048b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
6058b0e23d0SMatthew 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
6068b0e23d0SMatthew G. Knepley 
6078b0e23d0SMatthew G. Knepley   # Stokes preconditioners
608c4762a1bSJed Brown   #   Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
609c4762a1bSJed Brown   test:
6108b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_diagonal
6118b0e23d0SMatthew G. Knepley     requires: triangle
6128b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6138b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6148b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
6158b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
616c4762a1bSJed Brown   #   Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
617c4762a1bSJed Brown   test:
6188b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_triangular
6198b0e23d0SMatthew G. Knepley     requires: triangle
6208b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6218b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6228b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6238b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
624c4762a1bSJed Brown   #   Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
625c4762a1bSJed Brown   test:
6268b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_diagonal
6278b0e23d0SMatthew G. Knepley     requires: triangle
6288b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6298b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6308b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6318b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
6328b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
633c4762a1bSJed Brown   #   Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
634c4762a1bSJed Brown   test:
6358b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_upper
6368b0e23d0SMatthew G. Knepley     requires: triangle
6378b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
6388b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6398b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
6408b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
641c4762a1bSJed Brown   #   Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
642c4762a1bSJed Brown   test:
6438b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_lower
6448b0e23d0SMatthew G. Knepley     requires: triangle
6458b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6468b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6478b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6488b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
6498b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
650c4762a1bSJed 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}
651c4762a1bSJed Brown   test:
6528b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_full
6538b0e23d0SMatthew G. Knepley     requires: triangle
6548b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6558b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6568b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6578b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \
6588b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
6598b0e23d0SMatthew G. Knepley   #   Full Schur + Velocity GMG
6608b0e23d0SMatthew G. Knepley   test:
6618b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_gmg_vcycle
6628b0e23d0SMatthew G. Knepley     requires: triangle
6638b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6648b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6658b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-9 -ksp_error_if_not_converged -pc_use_amat \
6668b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
6678b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type mg -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
668c4762a1bSJed 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}
669c4762a1bSJed Brown   test:
6708b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_simple
671c4762a1bSJed Brown     requires: triangle
6728b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -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 fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6758b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
6768b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
6778b0e23d0SMatthew 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
678c4762a1bSJed Brown   #   FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
679c4762a1bSJed Brown   test:
6808b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_fetidp
681c4762a1bSJed Brown     requires: triangle mumps
682c4762a1bSJed Brown     nsize: 5
6838b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6848b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6858b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6868b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6878b0e23d0SMatthew 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 \
6888b0e23d0SMatthew 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
689c4762a1bSJed Brown   test:
6908b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_fetidp
6918b0e23d0SMatthew G. Knepley     requires: mumps
692c4762a1bSJed Brown     nsize: 5
6938b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6948b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
6958b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6968b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6978b0e23d0SMatthew 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 \
6988b0e23d0SMatthew 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
699c4762a1bSJed Brown   test:
7008b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_fetidp
7018b0e23d0SMatthew G. Knepley     requires: ctetgen mumps suitesparse
702c4762a1bSJed Brown     nsize: 5
7038b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
7048b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
7058b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
7068b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
7078b0e23d0SMatthew 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 \
7088b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
7098b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
7108b0e23d0SMatthew 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 \
7118b0e23d0SMatthew G. Knepley         -fetidp_bddelta_pc_factor_mat_ordering_type external \
7128b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
7138b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
714c4762a1bSJed Brown   test:
7158b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_fetidp
716c4762a1bSJed Brown     requires: suitesparse
717c4762a1bSJed Brown     nsize: 5
7188b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
7198b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
7208b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
7218b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
7228b0e23d0SMatthew 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 \
7238b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
7248b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
7258b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
7268b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
7278b0e23d0SMatthew G. Knepley   #   BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
728c4762a1bSJed Brown   test:
7298b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_bddc
7308b0e23d0SMatthew G. Knepley     nsize: 2
731c4762a1bSJed Brown     requires: triangle !single
7328b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
7338b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
7348b0e23d0SMatthew G. Knepley       -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
7358b0e23d0SMatthew 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
7368b0e23d0SMatthew G. Knepley   #   Vanka
737c4762a1bSJed Brown   test:
7388b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka
739c4762a1bSJed Brown     requires: double !complex
7408b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7418b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7428b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
743c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
744c4762a1bSJed Brown         -sub_ksp_type preonly -sub_pc_type lu
745c4762a1bSJed Brown   test:
7468b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka_denseinv
747c4762a1bSJed Brown     requires: double !complex
7488b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7498b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7508b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
751c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
752c4762a1bSJed Brown         -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
753c4762a1bSJed Brown   #   Vanka smoother
754c4762a1bSJed Brown   test:
7558b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_gmg_vanka
7568b0e23d0SMatthew G. Knepley     requires: double !complex
7578b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine_hierarchy 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7588b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
7598b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
7608b0e23d0SMatthew G. Knepley       -pc_type mg \
7618b0e23d0SMatthew G. Knepley         -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
762c4762a1bSJed 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 \
763c4762a1bSJed Brown           -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
764c4762a1bSJed Brown         -mg_coarse_pc_type svd
765c4762a1bSJed Brown 
766c4762a1bSJed Brown TEST*/
767