18b0e23d0SMatthew G. Knepley static char help[] = "Stokes Problem discretized with finite elements,\n\
28b0e23d0SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown /*
58b0e23d0SMatthew G. Knepley For the isoviscous Stokes problem, which we discretize using the finite
68b0e23d0SMatthew G. Knepley element method on an unstructured mesh, the weak form equations are
7c4762a1bSJed Brown
88b0e23d0SMatthew G. Knepley < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0
98b0e23d0SMatthew G. Knepley < q, -\nabla\cdot u > = 0
10c4762a1bSJed Brown
11c4762a1bSJed Brown Viewing:
12c4762a1bSJed Brown
13c4762a1bSJed Brown To produce nice output, use
14c4762a1bSJed Brown
158b0e23d0SMatthew G. Knepley -dm_refine 3 -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -snes_view_solution hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append
16c4762a1bSJed Brown
17c4762a1bSJed Brown You can get a LaTeX view of the mesh, with point numbering using
18c4762a1bSJed Brown
19c4762a1bSJed Brown -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0
20c4762a1bSJed Brown
21c4762a1bSJed Brown The data layout can be viewed using
22c4762a1bSJed Brown
23c4762a1bSJed Brown -dm_petscsection_view
24c4762a1bSJed Brown
25c4762a1bSJed Brown Lots of information about the FEM assembly can be printed using
26c4762a1bSJed Brown
278b0e23d0SMatthew G. Knepley -dm_plex_print_fem 3
28c4762a1bSJed Brown */
29c4762a1bSJed Brown
30c4762a1bSJed Brown #include <petscdmplex.h>
31c4762a1bSJed Brown #include <petscsnes.h>
32c4762a1bSJed Brown #include <petscds.h>
338b0e23d0SMatthew G. Knepley #include <petscbag.h>
34c4762a1bSJed Brown
358b0e23d0SMatthew G. Knepley // TODO: Plot residual by fields after each smoother iterate
36c4762a1bSJed Brown
3711486bccSBarry Smith typedef enum {
3811486bccSBarry Smith SOL_QUADRATIC,
3911486bccSBarry Smith SOL_TRIG,
4011486bccSBarry Smith SOL_UNKNOWN
4111486bccSBarry Smith } SolType;
428b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0};
43c4762a1bSJed Brown
449371c9d4SSatish Balay typedef struct {
458b0e23d0SMatthew G. Knepley PetscScalar mu; /* dynamic shear viscosity */
468b0e23d0SMatthew G. Knepley } Parameter;
478b0e23d0SMatthew G. Knepley
489371c9d4SSatish Balay typedef struct {
498b0e23d0SMatthew G. Knepley PetscBag bag; /* Problem parameters */
508b0e23d0SMatthew G. Knepley SolType sol; /* MMS solution */
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown
f1_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])53d71ae5a4SJacob Faibussowitsch static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
54d71ae5a4SJacob Faibussowitsch {
558b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]);
568b0e23d0SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
578b0e23d0SMatthew G. Knepley PetscInt c, d;
58c4762a1bSJed Brown
598b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) {
60ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[c * dim + d] = mu * (u_x[c * dim + d] + u_x[d * dim + c]);
618b0e23d0SMatthew G. Knepley f1[c * dim + c] -= u[uOff[1]];
62c4762a1bSJed Brown }
63c4762a1bSJed Brown }
64c4762a1bSJed Brown
f0_p(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])65d71ae5a4SJacob Faibussowitsch static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
66d71ae5a4SJacob Faibussowitsch {
67c4762a1bSJed Brown PetscInt d;
688b0e23d0SMatthew G. Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d * dim + d];
69c4762a1bSJed Brown }
70c4762a1bSJed Brown
g1_pu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g1[])71d71ae5a4SJacob Faibussowitsch static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
72d71ae5a4SJacob Faibussowitsch {
73c4762a1bSJed Brown PetscInt d;
748b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0; /* < q, -\nabla\cdot u > */
75c4762a1bSJed Brown }
76c4762a1bSJed Brown
g2_up(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g2[])77d71ae5a4SJacob Faibussowitsch static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
78d71ae5a4SJacob Faibussowitsch {
79c4762a1bSJed Brown PetscInt d;
808b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* -< \nabla\cdot v, p > */
81c4762a1bSJed Brown }
82c4762a1bSJed Brown
g3_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])83d71ae5a4SJacob Faibussowitsch static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
84d71ae5a4SJacob Faibussowitsch {
858b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]);
868b0e23d0SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
878b0e23d0SMatthew G. Knepley PetscInt c, d;
88c4762a1bSJed Brown
898b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) {
90c4762a1bSJed Brown for (d = 0; d < dim; ++d) {
918b0e23d0SMatthew G. Knepley g3[((c * Nc + c) * dim + d) * dim + d] += mu; /* < \nabla v, \nabla u > */
928b0e23d0SMatthew G. Knepley g3[((c * Nc + d) * dim + d) * dim + c] += mu; /* < \nabla v, {\nabla u}^T > */
93c4762a1bSJed Brown }
94c4762a1bSJed Brown }
95c4762a1bSJed Brown }
96c4762a1bSJed Brown
g0_pp(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])97d71ae5a4SJacob Faibussowitsch static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
98d71ae5a4SJacob Faibussowitsch {
998b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]);
1008b0e23d0SMatthew G. Knepley
1018b0e23d0SMatthew G. Knepley g0[0] = 1.0 / mu;
1028b0e23d0SMatthew G. Knepley }
1038b0e23d0SMatthew G. Knepley
1048b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution
1058b0e23d0SMatthew G. Knepley 2D:
106c4762a1bSJed Brown
107c4762a1bSJed Brown u = x^2 + y^2
1088b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy
1098b0e23d0SMatthew G. Knepley p = x + y - 1
1108b0e23d0SMatthew G. Knepley f = <1 - 4 mu, 1 - 4 mu>
111c4762a1bSJed Brown
112c4762a1bSJed Brown so that
113c4762a1bSJed Brown
1148b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 4x 4x \
1158b0e23d0SMatthew G. Knepley \ 4x -4x /
1168b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
1178b0e23d0SMatthew G. Knepley \nabla \cdot u = 2x - 2x = 0
1188b0e23d0SMatthew G. Knepley
1198b0e23d0SMatthew G. Knepley 3D:
1208b0e23d0SMatthew G. Knepley
1218b0e23d0SMatthew G. Knepley u = 2 x^2 + y^2 + z^2
1228b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy
1238b0e23d0SMatthew G. Knepley w = 2 x^2 - 2xz
1248b0e23d0SMatthew G. Knepley p = x + y + z - 3/2
1258b0e23d0SMatthew G. Knepley f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
1268b0e23d0SMatthew G. Knepley
1278b0e23d0SMatthew G. Knepley so that
1288b0e23d0SMatthew G. Knepley
1298b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 8x 4x 4x \
1308b0e23d0SMatthew G. Knepley | 4x -4x 0 |
1318b0e23d0SMatthew G. Knepley \ 4x 0 -4x /
1328b0e23d0SMatthew 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
1338b0e23d0SMatthew G. Knepley \nabla \cdot u = 4x - 2x - 2x = 0
134c4762a1bSJed Brown */
quadratic_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)135*2a8381b2SBarry Smith static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
136d71ae5a4SJacob Faibussowitsch {
1378b0e23d0SMatthew G. Knepley PetscInt c;
1388b0e23d0SMatthew G. Knepley
1398b0e23d0SMatthew G. Knepley u[0] = (dim - 1) * PetscSqr(x[0]);
1408b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) {
1418b0e23d0SMatthew G. Knepley u[0] += PetscSqr(x[c]);
1428b0e23d0SMatthew G. Knepley u[c] = 2.0 * PetscSqr(x[0]) - 2.0 * x[0] * x[c];
1438b0e23d0SMatthew G. Knepley }
1443ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
145c4762a1bSJed Brown }
146c4762a1bSJed Brown
quadratic_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)147*2a8381b2SBarry Smith static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
148d71ae5a4SJacob Faibussowitsch {
1498b0e23d0SMatthew G. Knepley PetscInt d;
1508b0e23d0SMatthew G. Knepley
1518b0e23d0SMatthew G. Knepley u[0] = -0.5 * dim;
1528b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d];
1533ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown
f0_quadratic_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])156d71ae5a4SJacob Faibussowitsch static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
157d71ae5a4SJacob Faibussowitsch {
1588b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]);
1598b0e23d0SMatthew G. Knepley PetscInt d;
1608b0e23d0SMatthew G. Knepley
1618b0e23d0SMatthew G. Knepley f0[0] = (dim - 1) * 4.0 * mu - 1.0;
1628b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) f0[d] = 4.0 * mu - 1.0;
163c4762a1bSJed Brown }
164c4762a1bSJed Brown
1658b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution
1668b0e23d0SMatthew G. Knepley 2D:
1678b0e23d0SMatthew G. Knepley
1688b0e23d0SMatthew G. Knepley u = sin(pi x) + sin(pi y)
1698b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y
1708b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y)
1718b0e23d0SMatthew 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>
1728b0e23d0SMatthew G. Knepley
1738b0e23d0SMatthew G. Knepley so that
1748b0e23d0SMatthew G. Knepley
1758b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 2pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y \
1768b0e23d0SMatthew G. Knepley \ pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) /
1778b0e23d0SMatthew 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
1788b0e23d0SMatthew G. Knepley \nabla \cdot u = pi cos(pi x) - pi cos(pi x) = 0
1798b0e23d0SMatthew G. Knepley
1808b0e23d0SMatthew G. Knepley 3D:
1818b0e23d0SMatthew G. Knepley
1828b0e23d0SMatthew G. Knepley u = 2 sin(pi x) + sin(pi y) + sin(pi z)
1838b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y
1848b0e23d0SMatthew G. Knepley w = -pi cos(pi x) z
1858b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
1868b0e23d0SMatthew 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>
1878b0e23d0SMatthew G. Knepley
1888b0e23d0SMatthew G. Knepley so that
1898b0e23d0SMatthew G. Knepley
1908b0e23d0SMatthew 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 \
1918b0e23d0SMatthew G. Knepley | pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) 0 |
1928b0e23d0SMatthew G. Knepley \ pi cos(pi z) + pi^2 sin(pi x) z 0 -2pi cos(pi x) /
1938b0e23d0SMatthew 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
1948b0e23d0SMatthew G. Knepley \nabla \cdot u = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
1958b0e23d0SMatthew G. Knepley */
trig_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)196*2a8381b2SBarry Smith static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
197d71ae5a4SJacob Faibussowitsch {
1988b0e23d0SMatthew G. Knepley PetscInt c;
1998b0e23d0SMatthew G. Knepley
2008b0e23d0SMatthew G. Knepley u[0] = (dim - 1) * PetscSinReal(PETSC_PI * x[0]);
2018b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) {
2028b0e23d0SMatthew G. Knepley u[0] += PetscSinReal(PETSC_PI * x[c]);
2038b0e23d0SMatthew G. Knepley u[c] = -PETSC_PI * PetscCosReal(PETSC_PI * x[0]) * x[c];
2048b0e23d0SMatthew G. Knepley }
2053ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
2068b0e23d0SMatthew G. Knepley }
2078b0e23d0SMatthew G. Knepley
trig_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)208*2a8381b2SBarry Smith static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
209d71ae5a4SJacob Faibussowitsch {
2108b0e23d0SMatthew G. Knepley PetscInt d;
2118b0e23d0SMatthew G. Knepley
2128b0e23d0SMatthew G. Knepley for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0 * PETSC_PI * x[d]);
2133ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
2148b0e23d0SMatthew G. Knepley }
2158b0e23d0SMatthew G. Knepley
f0_trig_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])216d71ae5a4SJacob Faibussowitsch static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
217d71ae5a4SJacob Faibussowitsch {
2188b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]);
2198b0e23d0SMatthew G. Knepley PetscInt d;
2208b0e23d0SMatthew G. Knepley
2218b0e23d0SMatthew 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]);
2228b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) {
2238b0e23d0SMatthew G. Knepley f0[0] -= mu * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x[d]);
2248b0e23d0SMatthew 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];
2258b0e23d0SMatthew G. Knepley }
2268b0e23d0SMatthew G. Knepley }
2278b0e23d0SMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)228d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
229d71ae5a4SJacob Faibussowitsch {
2308b0e23d0SMatthew G. Knepley PetscInt sol;
231c4762a1bSJed Brown
232c4762a1bSJed Brown PetscFunctionBeginUser;
2338b0e23d0SMatthew G. Knepley options->sol = SOL_QUADRATIC;
234d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
2358b0e23d0SMatthew G. Knepley sol = options->sol;
236dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, PETSC_STATIC_ARRAY_LENGTH(SolTypes) - 3, SolTypes[options->sol], &sol, NULL));
2378b0e23d0SMatthew G. Knepley options->sol = (SolType)sol;
238d0609cedSBarry Smith PetscOptionsEnd();
2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
240c4762a1bSJed Brown }
241c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)242d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
243d71ae5a4SJacob Faibussowitsch {
244c4762a1bSJed Brown PetscFunctionBeginUser;
2459566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
2469566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
2479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
2489566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
250c4762a1bSJed Brown }
251c4762a1bSJed Brown
SetupParameters(MPI_Comm comm,AppCtx * ctx)252d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
253d71ae5a4SJacob Faibussowitsch {
2548b0e23d0SMatthew G. Knepley Parameter *p;
2558b0e23d0SMatthew G. Knepley
2568b0e23d0SMatthew G. Knepley PetscFunctionBeginUser;
2578b0e23d0SMatthew G. Knepley /* setup PETSc parameter bag */
2589566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag));
259*2a8381b2SBarry Smith PetscCall(PetscBagGetData(ctx->bag, &p));
2609566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag, "par", "Stokes Parameters"));
2619566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s"));
2629566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(ctx->bag));
2638b0e23d0SMatthew G. Knepley {
2648b0e23d0SMatthew G. Knepley PetscViewer viewer;
2658b0e23d0SMatthew G. Knepley PetscViewerFormat format;
2668b0e23d0SMatthew G. Knepley PetscBool flg;
2678b0e23d0SMatthew G. Knepley
268648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
2698b0e23d0SMatthew G. Knepley if (flg) {
2709566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format));
2719566063dSJacob Faibussowitsch PetscCall(PetscBagView(ctx->bag, viewer));
2729566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer));
2739566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
274648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer));
2758b0e23d0SMatthew G. Knepley }
2768b0e23d0SMatthew G. Knepley }
2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2788b0e23d0SMatthew G. Knepley }
2798b0e23d0SMatthew G. Knepley
SetupEqn(DM dm,AppCtx * user)280d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
281d71ae5a4SJacob Faibussowitsch {
2828b0e23d0SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
2838b0e23d0SMatthew G. Knepley PetscDS ds;
28445480ffeSMatthew G. Knepley DMLabel label;
285c4762a1bSJed Brown const PetscInt id = 1;
286c4762a1bSJed Brown
287c4762a1bSJed Brown PetscFunctionBeginUser;
2889566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2898b0e23d0SMatthew G. Knepley switch (user->sol) {
290c4762a1bSJed Brown case SOL_QUADRATIC:
2919566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u));
2928b0e23d0SMatthew G. Knepley exactFuncs[0] = quadratic_u;
2938b0e23d0SMatthew G. Knepley exactFuncs[1] = quadratic_p;
294c4762a1bSJed Brown break;
295c4762a1bSJed Brown case SOL_TRIG:
2969566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
2978b0e23d0SMatthew G. Knepley exactFuncs[0] = trig_u;
2988b0e23d0SMatthew G. Knepley exactFuncs[1] = trig_p;
299c4762a1bSJed Brown break;
300d71ae5a4SJacob Faibussowitsch default:
301d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
302c4762a1bSJed Brown }
3039566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL));
3049566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
3059566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
3069566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
3079566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu));
3089566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL));
309c4762a1bSJed Brown
3109566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user));
3119566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user));
312c4762a1bSJed Brown
3139566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
31457d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exactFuncs[0], NULL, user, NULL));
3158b0e23d0SMatthew G. Knepley
31647bb1945SPatrick Sanan /* Make constant values available to pointwise functions */
317c4762a1bSJed Brown {
3188b0e23d0SMatthew G. Knepley Parameter *param;
3198b0e23d0SMatthew G. Knepley PetscScalar constants[1];
320c4762a1bSJed Brown
321*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m));
3228b0e23d0SMatthew G. Knepley constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
3239566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants));
324c4762a1bSJed Brown }
3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
326c4762a1bSJed Brown }
327c4762a1bSJed Brown
zero(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)328*2a8381b2SBarry Smith static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
329d71ae5a4SJacob Faibussowitsch {
3308b0e23d0SMatthew G. Knepley PetscInt c;
3318b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0;
3323ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
3338b0e23d0SMatthew G. Knepley }
one(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)334*2a8381b2SBarry Smith static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
335d71ae5a4SJacob Faibussowitsch {
3368b0e23d0SMatthew G. Knepley PetscInt c;
3378b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 1.0;
3383ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
3398b0e23d0SMatthew G. Knepley }
3408b0e23d0SMatthew G. Knepley
CreatePressureNullSpace(DM dm,PetscInt origField,PetscInt field,MatNullSpace * nullspace)341d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
342d71ae5a4SJacob Faibussowitsch {
343c4762a1bSJed Brown Vec vec;
344*2a8381b2SBarry Smith PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) = {zero, one};
345c4762a1bSJed Brown
346c4762a1bSJed Brown PetscFunctionBeginUser;
34763a3b9bcSJacob Faibussowitsch PetscCheck(origField == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField);
3488b0e23d0SMatthew G. Knepley funcs[field] = one;
3498b0e23d0SMatthew G. Knepley {
3508b0e23d0SMatthew G. Knepley PetscDS ds;
3519566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3529566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)ds, NULL, "-ds_view"));
3538b0e23d0SMatthew G. Knepley }
3549566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &vec));
3559566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
3569566063dSJacob Faibussowitsch PetscCall(VecNormalize(vec, NULL));
3579566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace));
3589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec));
359c4762a1bSJed Brown /* New style for field null spaces */
360c4762a1bSJed Brown {
361c4762a1bSJed Brown PetscObject pressure;
362c4762a1bSJed Brown MatNullSpace nullspacePres;
363c4762a1bSJed Brown
3649566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &pressure));
3659566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
3669566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres));
3679566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullspacePres));
368c4762a1bSJed Brown }
3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
370c4762a1bSJed Brown }
371c4762a1bSJed Brown
SetupProblem(DM dm,PetscErrorCode (* setupEqn)(DM,AppCtx *),AppCtx * user)372d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
373d71ae5a4SJacob Faibussowitsch {
3748b0e23d0SMatthew G. Knepley DM cdm = dm;
3758b0e23d0SMatthew G. Knepley PetscQuadrature q = NULL;
3768b0e23d0SMatthew G. Knepley PetscBool simplex;
37730602db0SMatthew G. Knepley PetscInt dim, Nf = 2, f, Nc[2];
3788b0e23d0SMatthew G. Knepley const char *name[2] = {"velocity", "pressure"};
3798b0e23d0SMatthew G. Knepley const char *prefix[2] = {"vel_", "pres_"};
380c4762a1bSJed Brown
3818b0e23d0SMatthew G. Knepley PetscFunctionBegin;
3829566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
3839566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
3848b0e23d0SMatthew G. Knepley Nc[0] = dim;
3858b0e23d0SMatthew G. Knepley Nc[1] = 1;
3868b0e23d0SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
3878b0e23d0SMatthew G. Knepley PetscFE fe;
3888b0e23d0SMatthew G. Knepley
3899566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe));
3909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name[f]));
3919566063dSJacob Faibussowitsch if (!q) PetscCall(PetscFEGetQuadrature(fe, &q));
3929566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, q));
3939566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
3949566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
395c4762a1bSJed Brown }
3969566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
3979566063dSJacob Faibussowitsch PetscCall((*setupEqn)(dm, user));
3988b0e23d0SMatthew G. Knepley while (cdm) {
3999566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
4009566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace));
4019566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
402c4762a1bSJed Brown }
4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
404c4762a1bSJed Brown }
405c4762a1bSJed Brown
main(int argc,char ** argv)406d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
407d71ae5a4SJacob Faibussowitsch {
4088b0e23d0SMatthew G. Knepley SNES snes;
4098b0e23d0SMatthew G. Knepley DM dm;
4108b0e23d0SMatthew G. Knepley Vec u;
4118b0e23d0SMatthew G. Knepley AppCtx user;
412c4762a1bSJed Brown
413327415f7SBarry Smith PetscFunctionBeginUser;
4149566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4159566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
4169566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4179566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
4189566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm));
4199566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user));
420c4762a1bSJed Brown
4219566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
4229566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, SetupEqn, &user));
4239566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL));
424c4762a1bSJed Brown
4259566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
4266493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
4279566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
4289566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u));
4299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Solution"));
4308b0e23d0SMatthew G. Knepley {
4318b0e23d0SMatthew G. Knepley Mat J;
4328b0e23d0SMatthew G. Knepley MatNullSpace sp;
433c4762a1bSJed Brown
4349566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes));
4359566063dSJacob Faibussowitsch PetscCall(CreatePressureNullSpace(dm, 1, 1, &sp));
4369566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
4379566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J, sp));
4389566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sp));
4399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
4409566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-J_view"));
441c4762a1bSJed Brown }
4429566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
443c4762a1bSJed Brown
4449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
4459566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
4469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
4479566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag));
4489566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
449b122ec5aSJacob Faibussowitsch return 0;
450c4762a1bSJed Brown }
451c4762a1bSJed Brown /*TEST
452c4762a1bSJed Brown
453c4762a1bSJed Brown test:
4548b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check
455c4762a1bSJed Brown requires: triangle
4568b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4578b0e23d0SMatthew G. Knepley
458c4762a1bSJed Brown test:
4598b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check_parallel
4608b0e23d0SMatthew G. Knepley nsize: {{2 3 5}}
461c4762a1bSJed Brown requires: triangle
462e600fa54SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4638b0e23d0SMatthew G. Knepley
464c4762a1bSJed Brown test:
4658b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check
466c4762a1bSJed Brown requires: ctetgen
46730602db0SMatthew 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
4688b0e23d0SMatthew G. Knepley
469c4762a1bSJed Brown test:
4708b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check_parallel
4718b0e23d0SMatthew G. Knepley nsize: {{2 3 5}}
472c4762a1bSJed Brown requires: ctetgen
47319fa4260SStefano Zampini args: -sol quadratic -dm_refine 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
4748b0e23d0SMatthew G. Knepley
475c4762a1bSJed Brown test:
4768b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv
4778b0e23d0SMatthew G. Knepley requires: triangle
4788b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
4798b0e23d0SMatthew 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 \
4808b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
4818b0e23d0SMatthew 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 \
4828b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
4838b0e23d0SMatthew G. Knepley
484c4762a1bSJed Brown test:
4858b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv_gamg
4868b0e23d0SMatthew G. Knepley requires: triangle
48782894d03SBarry Smith args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 \
4888b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
489df082c0eSPierre Jolivet -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_explicit_operator_mat_type aij -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
4908b0e23d0SMatthew G. Knepley
491c4762a1bSJed Brown test:
4928b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_conv
4938b0e23d0SMatthew G. Knepley requires: ctetgen !single
4948b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
49530602db0SMatthew 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 \
4968b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
4978b0e23d0SMatthew 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 \
4988b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
4998b0e23d0SMatthew G. Knepley
500c4762a1bSJed Brown test:
5018b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_check
50230602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
5038b0e23d0SMatthew G. Knepley
504c4762a1bSJed Brown test:
5058b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_check
50630602db0SMatthew 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
5078b0e23d0SMatthew G. Knepley
508c4762a1bSJed Brown test:
5098b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_conv
5108b0e23d0SMatthew G. Knepley # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
51130602db0SMatthew 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 \
5128b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5138b0e23d0SMatthew 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 \
5148b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5158b0e23d0SMatthew G. Knepley
516c4762a1bSJed Brown test:
5178b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_conv
518c4762a1bSJed Brown requires: !single
5198b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
52030602db0SMatthew 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 \
5218b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5228b0e23d0SMatthew 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 \
5238b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5248b0e23d0SMatthew G. Knepley
525c4762a1bSJed Brown test:
5268b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_check
5278b0e23d0SMatthew G. Knepley requires: triangle
5288b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
5298b0e23d0SMatthew G. Knepley
5308b0e23d0SMatthew G. Knepley test:
5318b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_check
5328b0e23d0SMatthew G. Knepley requires: ctetgen !single
53330602db0SMatthew 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
5348b0e23d0SMatthew G. Knepley
5358b0e23d0SMatthew G. Knepley test:
5368b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_conv
5378b0e23d0SMatthew G. Knepley requires: triangle
5388b0e23d0SMatthew G. Knepley # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
5398b0e23d0SMatthew 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 \
5408b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5418b0e23d0SMatthew 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 \
5428b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5438b0e23d0SMatthew G. Knepley
5448b0e23d0SMatthew G. Knepley test:
5458b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_conv
5468b0e23d0SMatthew G. Knepley requires: ctetgen long_runtime
5478b0e23d0SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
54830602db0SMatthew 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 \
5498b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
5508b0e23d0SMatthew 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 \
5518b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
5528b0e23d0SMatthew G. Knepley
5538b0e23d0SMatthew G. Knepley test:
5548b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_conv
555c4762a1bSJed Brown requires: !single
5568b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
55730602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
55882894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \
5598b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
560df082c0eSPierre Jolivet -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_explicit_operator_mat_type aij -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0
5618b0e23d0SMatthew G. Knepley
562c4762a1bSJed Brown test:
5638b0e23d0SMatthew G. Knepley suffix: 3d_q1_p0_conv
5648b0e23d0SMatthew G. Knepley requires: !single
5658b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
56630602db0SMatthew 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 \
56782894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \
5688b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
569df082c0eSPierre Jolivet -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_explicit_operator_mat_type aij -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0
5708b0e23d0SMatthew G. Knepley
5718b0e23d0SMatthew G. Knepley # Stokes preconditioners
572c4762a1bSJed Brown # Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
573c4762a1bSJed Brown test:
5748b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_block_diagonal
5758b0e23d0SMatthew G. Knepley requires: triangle
5768b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
5778b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
5788b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
5798b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
5803886731fSPierre Jolivet output_file: output/empty.out
581c4762a1bSJed Brown # Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
582c4762a1bSJed Brown test:
5838b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_block_triangular
5848b0e23d0SMatthew G. Knepley requires: triangle
5858b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
5868b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
5878b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
5888b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
5893886731fSPierre Jolivet output_file: output/empty.out
590c4762a1bSJed Brown # Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
591c4762a1bSJed Brown test:
5928b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_diagonal
5938b0e23d0SMatthew G. Knepley requires: triangle
5948b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
5958b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
5968b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
5978b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
5988b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
5993886731fSPierre Jolivet output_file: output/empty.out
600c4762a1bSJed Brown # Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
601c4762a1bSJed Brown test:
6028b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_upper
6038b0e23d0SMatthew G. Knepley requires: triangle
6048b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
6058b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6068b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
6078b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
608c4762a1bSJed Brown # Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
609c4762a1bSJed Brown test:
6108b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_lower
6118b0e23d0SMatthew G. Knepley requires: triangle
6128b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6138b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
6148b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
6158b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
6168b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
6173886731fSPierre Jolivet output_file: output/empty.out
618c4762a1bSJed 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}
619c4762a1bSJed Brown test:
6208b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_full
6218b0e23d0SMatthew G. Knepley requires: triangle
6228b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
6238b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
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 full -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
6273886731fSPierre Jolivet output_file: output/empty.out
6288b0e23d0SMatthew G. Knepley # Full Schur + Velocity GMG
6298b0e23d0SMatthew G. Knepley test:
6308b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_gmg_vcycle
63196b316e5SStefano Zampini TODO: broken (requires subDMs hooks)
6328b0e23d0SMatthew G. Knepley requires: triangle
6338b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
63482894d03SBarry Smith -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \
6358b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
63673f7197eSJed 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
637c4762a1bSJed 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}
638c4762a1bSJed Brown test:
6398b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_simple
640c4762a1bSJed Brown requires: triangle
6418b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
6428b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
6438b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
6448b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
6458b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
6468b0e23d0SMatthew 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
6473886731fSPierre Jolivet output_file: output/empty.out
648c4762a1bSJed Brown # FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
649c4762a1bSJed Brown test:
6508b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_fetidp
651c4762a1bSJed Brown requires: triangle mumps
652c4762a1bSJed Brown nsize: 5
653e600fa54SMatthew 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 \
6548b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
65582894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \
6568b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6578b0e23d0SMatthew 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 \
6588b0e23d0SMatthew 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
6593886731fSPierre Jolivet output_file: output/empty.out
660c4762a1bSJed Brown test:
6618b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_fetidp
6628b0e23d0SMatthew G. Knepley requires: mumps
663c4762a1bSJed Brown nsize: 5
664e600fa54SMatthew 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 \
6658b0e23d0SMatthew G. Knepley -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
6668b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6678b0e23d0SMatthew 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 \
6688b0e23d0SMatthew 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
6693886731fSPierre Jolivet output_file: output/empty.out
670c4762a1bSJed Brown test:
6718b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_fetidp
6728b0e23d0SMatthew G. Knepley requires: ctetgen mumps suitesparse
673c4762a1bSJed Brown nsize: 5
674e600fa54SMatthew 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 \
6758b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
67682894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-9 \
6778b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6788b0e23d0SMatthew 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 \
6798b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
6808b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
6818b0e23d0SMatthew 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 \
6828b0e23d0SMatthew G. Knepley -fetidp_bddelta_pc_factor_mat_ordering_type external \
6838b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
6848b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
6853886731fSPierre Jolivet output_file: output/empty.out
686c4762a1bSJed Brown test:
6878b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_fetidp
688c4762a1bSJed Brown requires: suitesparse
689c4762a1bSJed Brown nsize: 5
690e600fa54SMatthew 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 \
6918b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
69282894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \
6938b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
6948b0e23d0SMatthew 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 \
6958b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
6968b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
6978b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
6988b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
6993886731fSPierre Jolivet output_file: output/empty.out
7008b0e23d0SMatthew G. Knepley # BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
701c4762a1bSJed Brown test:
7028b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_bddc
7038b0e23d0SMatthew G. Knepley nsize: 2
704c4762a1bSJed Brown requires: triangle !single
705e600fa54SMatthew 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 \
7068b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \
7078b0e23d0SMatthew G. Knepley -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
7088b0e23d0SMatthew 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
7093886731fSPierre Jolivet output_file: output/empty.out
7108b0e23d0SMatthew G. Knepley # Vanka
711c4762a1bSJed Brown test:
7128b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka
7133886731fSPierre Jolivet output_file: output/empty.out
714c4762a1bSJed Brown requires: double !complex
71530602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
7168b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \
7178b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
718c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
719c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu
720c4762a1bSJed Brown test:
7218b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka_denseinv
7223886731fSPierre Jolivet output_file: output/empty.out
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 -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
729c4762a1bSJed Brown # Vanka smoother
730c4762a1bSJed Brown test:
7318b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_gmg_vanka
7323886731fSPierre Jolivet output_file: output/empty.out
7338b0e23d0SMatthew G. Knepley requires: double !complex
73430602db0SMatthew 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 \
7358b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \
7368b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
7378b0e23d0SMatthew G. Knepley -pc_type mg \
7388b0e23d0SMatthew G. Knepley -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
739c4762a1bSJed 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 \
740c4762a1bSJed Brown -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
741c4762a1bSJed Brown -mg_coarse_pc_type svd
742c4762a1bSJed Brown
743c4762a1bSJed Brown TEST*/
744