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 37*11486bccSBarry Smith typedef enum { 38*11486bccSBarry Smith SOL_QUADRATIC, 39*11486bccSBarry Smith SOL_TRIG, 40*11486bccSBarry Smith SOL_UNKNOWN 41*11486bccSBarry Smith } SolType; 428b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0}; 43c4762a1bSJed Brown 44*11486bccSBarry Smith typedef struct 45*11486bccSBarry Smith { 468b0e23d0SMatthew G. Knepley PetscScalar mu; /* dynamic shear viscosity */ 478b0e23d0SMatthew G. Knepley } Parameter; 488b0e23d0SMatthew G. Knepley 49*11486bccSBarry Smith typedef struct 50*11486bccSBarry Smith { 518b0e23d0SMatthew G. Knepley PetscBag bag; /* Problem parameters */ 528b0e23d0SMatthew G. Knepley SolType sol; /* MMS solution */ 53c4762a1bSJed Brown } AppCtx; 54c4762a1bSJed Brown 55*11486bccSBarry Smith 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[]) 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) { 62*11486bccSBarry Smith for (d = 0; d < dim; ++d) { f1[c * dim + d] = mu * (u_x[c * dim + d] + u_x[d * dim + c]); } 638b0e23d0SMatthew G. Knepley f1[c * dim + c] -= u[uOff[1]]; 64c4762a1bSJed Brown } 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67*11486bccSBarry Smith 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[]) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown PetscInt d; 708b0e23d0SMatthew G. Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d * dim + d]; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73*11486bccSBarry Smith 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[]) 74c4762a1bSJed Brown { 75c4762a1bSJed Brown PetscInt d; 768b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0; /* < q, -\nabla\cdot u > */ 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79*11486bccSBarry Smith 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[]) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt d; 828b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* -< \nabla\cdot v, p > */ 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85*11486bccSBarry Smith 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[]) 86c4762a1bSJed Brown { 878b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 888b0e23d0SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 898b0e23d0SMatthew G. Knepley PetscInt c, d; 90c4762a1bSJed Brown 918b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 92c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 938b0e23d0SMatthew G. Knepley g3[((c * Nc + c) * dim + d) * dim + d] += mu; /* < \nabla v, \nabla u > */ 948b0e23d0SMatthew G. Knepley g3[((c * Nc + d) * dim + d) * dim + c] += mu; /* < \nabla v, {\nabla u}^T > */ 95c4762a1bSJed Brown } 96c4762a1bSJed Brown } 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99*11486bccSBarry Smith 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[]) 1008b0e23d0SMatthew G. Knepley { 1018b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 1028b0e23d0SMatthew G. Knepley 1038b0e23d0SMatthew G. Knepley g0[0] = 1.0 / mu; 1048b0e23d0SMatthew G. Knepley } 1058b0e23d0SMatthew G. Knepley 1068b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution 1078b0e23d0SMatthew G. Knepley 2D: 108c4762a1bSJed Brown 109c4762a1bSJed Brown u = x^2 + y^2 1108b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy 1118b0e23d0SMatthew G. Knepley p = x + y - 1 1128b0e23d0SMatthew G. Knepley f = <1 - 4 mu, 1 - 4 mu> 113c4762a1bSJed Brown 114c4762a1bSJed Brown so that 115c4762a1bSJed Brown 1168b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 4x 4x \ 1178b0e23d0SMatthew G. Knepley \ 4x -4x / 1188b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0 1198b0e23d0SMatthew G. Knepley \nabla \cdot u = 2x - 2x = 0 1208b0e23d0SMatthew G. Knepley 1218b0e23d0SMatthew G. Knepley 3D: 1228b0e23d0SMatthew G. Knepley 1238b0e23d0SMatthew G. Knepley u = 2 x^2 + y^2 + z^2 1248b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy 1258b0e23d0SMatthew G. Knepley w = 2 x^2 - 2xz 1268b0e23d0SMatthew G. Knepley p = x + y + z - 3/2 1278b0e23d0SMatthew G. Knepley f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> 1288b0e23d0SMatthew G. Knepley 1298b0e23d0SMatthew G. Knepley so that 1308b0e23d0SMatthew G. Knepley 1318b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 8x 4x 4x \ 1328b0e23d0SMatthew G. Knepley | 4x -4x 0 | 1338b0e23d0SMatthew G. Knepley \ 4x 0 -4x / 1348b0e23d0SMatthew 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 1358b0e23d0SMatthew G. Knepley \nabla \cdot u = 4x - 2x - 2x = 0 136c4762a1bSJed Brown */ 1378b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 138c4762a1bSJed Brown { 1398b0e23d0SMatthew G. Knepley PetscInt c; 1408b0e23d0SMatthew G. Knepley 1418b0e23d0SMatthew G. Knepley u[0] = (dim - 1) * PetscSqr(x[0]); 1428b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) { 1438b0e23d0SMatthew G. Knepley u[0] += PetscSqr(x[c]); 1448b0e23d0SMatthew G. Knepley u[c] = 2.0 * PetscSqr(x[0]) - 2.0 * x[0] * x[c]; 1458b0e23d0SMatthew G. Knepley } 146c4762a1bSJed Brown return 0; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 1498b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 150c4762a1bSJed Brown { 1518b0e23d0SMatthew G. Knepley PetscInt d; 1528b0e23d0SMatthew G. Knepley 1538b0e23d0SMatthew G. Knepley u[0] = -0.5 * dim; 1548b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d]; 155c4762a1bSJed Brown return 0; 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158*11486bccSBarry Smith 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[]) 159c4762a1bSJed Brown { 1608b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 1618b0e23d0SMatthew G. Knepley PetscInt d; 1628b0e23d0SMatthew G. Knepley 1638b0e23d0SMatthew G. Knepley f0[0] = (dim - 1) * 4.0 * mu - 1.0; 1648b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) f0[d] = 4.0 * mu - 1.0; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 1678b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution 1688b0e23d0SMatthew G. Knepley 2D: 1698b0e23d0SMatthew G. Knepley 1708b0e23d0SMatthew G. Knepley u = sin(pi x) + sin(pi y) 1718b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y 1728b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y) 1738b0e23d0SMatthew 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> 1748b0e23d0SMatthew G. Knepley 1758b0e23d0SMatthew G. Knepley so that 1768b0e23d0SMatthew G. Knepley 1778b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 2pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y \ 1788b0e23d0SMatthew G. Knepley \ pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) / 1798b0e23d0SMatthew 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 1808b0e23d0SMatthew G. Knepley \nabla \cdot u = pi cos(pi x) - pi cos(pi x) = 0 1818b0e23d0SMatthew G. Knepley 1828b0e23d0SMatthew G. Knepley 3D: 1838b0e23d0SMatthew G. Knepley 1848b0e23d0SMatthew G. Knepley u = 2 sin(pi x) + sin(pi y) + sin(pi z) 1858b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y 1868b0e23d0SMatthew G. Knepley w = -pi cos(pi x) z 1878b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z) 1888b0e23d0SMatthew 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> 1898b0e23d0SMatthew G. Knepley 1908b0e23d0SMatthew G. Knepley so that 1918b0e23d0SMatthew G. Knepley 1928b0e23d0SMatthew 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 \ 1938b0e23d0SMatthew G. Knepley | pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) 0 | 1948b0e23d0SMatthew G. Knepley \ pi cos(pi z) + pi^2 sin(pi x) z 0 -2pi cos(pi x) / 1958b0e23d0SMatthew 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 1968b0e23d0SMatthew G. Knepley \nabla \cdot u = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0 1978b0e23d0SMatthew G. Knepley */ 1988b0e23d0SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 199c4762a1bSJed Brown { 2008b0e23d0SMatthew G. Knepley PetscInt c; 2018b0e23d0SMatthew G. Knepley 2028b0e23d0SMatthew G. Knepley u[0] = (dim - 1) * PetscSinReal(PETSC_PI * x[0]); 2038b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) { 2048b0e23d0SMatthew G. Knepley u[0] += PetscSinReal(PETSC_PI * x[c]); 2058b0e23d0SMatthew G. Knepley u[c] = -PETSC_PI * PetscCosReal(PETSC_PI * x[0]) * x[c]; 2068b0e23d0SMatthew G. Knepley } 2078b0e23d0SMatthew G. Knepley return 0; 2088b0e23d0SMatthew G. Knepley } 2098b0e23d0SMatthew G. Knepley 2108b0e23d0SMatthew G. Knepley static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 2118b0e23d0SMatthew G. Knepley { 2128b0e23d0SMatthew G. Knepley PetscInt d; 2138b0e23d0SMatthew G. Knepley 2148b0e23d0SMatthew G. Knepley for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0 * PETSC_PI * x[d]); 2158b0e23d0SMatthew G. Knepley return 0; 2168b0e23d0SMatthew G. Knepley } 2178b0e23d0SMatthew G. Knepley 218*11486bccSBarry Smith 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[]) 2198b0e23d0SMatthew G. Knepley { 2208b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 2218b0e23d0SMatthew G. Knepley PetscInt d; 2228b0e23d0SMatthew G. Knepley 2238b0e23d0SMatthew 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]); 2248b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) { 2258b0e23d0SMatthew G. Knepley f0[0] -= mu * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x[d]); 2268b0e23d0SMatthew 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]; 2278b0e23d0SMatthew G. Knepley } 2288b0e23d0SMatthew G. Knepley } 2298b0e23d0SMatthew G. Knepley 2308b0e23d0SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 2318b0e23d0SMatthew G. Knepley { 2328b0e23d0SMatthew G. Knepley PetscInt sol; 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscFunctionBeginUser; 2358b0e23d0SMatthew G. Knepley options->sol = SOL_QUADRATIC; 236d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX"); 2378b0e23d0SMatthew G. Knepley sol = options->sol; 238dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, PETSC_STATIC_ARRAY_LENGTH(SolTypes) - 3, SolTypes[options->sol], &sol, NULL)); 2398b0e23d0SMatthew G. Knepley options->sol = (SolType)sol; 240d0609cedSBarry Smith PetscOptionsEnd(); 241c4762a1bSJed Brown PetscFunctionReturn(0); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 2448b0e23d0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 245c4762a1bSJed Brown { 246c4762a1bSJed Brown PetscFunctionBeginUser; 2479566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2489566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2499566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2509566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 2548b0e23d0SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 255c4762a1bSJed Brown { 2568b0e23d0SMatthew G. Knepley Parameter *p; 2578b0e23d0SMatthew G. Knepley 2588b0e23d0SMatthew G. Knepley PetscFunctionBeginUser; 2598b0e23d0SMatthew G. Knepley /* setup PETSc parameter bag */ 2609566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag)); 2619566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 2629566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag, "par", "Stokes Parameters")); 2639566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s")); 2649566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(ctx->bag)); 2658b0e23d0SMatthew G. Knepley { 2668b0e23d0SMatthew G. Knepley PetscViewer viewer; 2678b0e23d0SMatthew G. Knepley PetscViewerFormat format; 2688b0e23d0SMatthew G. Knepley PetscBool flg; 2698b0e23d0SMatthew G. Knepley 2709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 2718b0e23d0SMatthew G. Knepley if (flg) { 2729566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 2739566063dSJacob Faibussowitsch PetscCall(PetscBagView(ctx->bag, viewer)); 2749566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2759566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 2778b0e23d0SMatthew G. Knepley } 2788b0e23d0SMatthew G. Knepley } 2798b0e23d0SMatthew G. Knepley PetscFunctionReturn(0); 2808b0e23d0SMatthew G. Knepley } 2818b0e23d0SMatthew G. Knepley 2828b0e23d0SMatthew G. Knepley static PetscErrorCode SetupEqn(DM dm, AppCtx *user) 2838b0e23d0SMatthew G. Knepley { 2848b0e23d0SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 2858b0e23d0SMatthew G. Knepley PetscDS ds; 28645480ffeSMatthew G. Knepley DMLabel label; 287c4762a1bSJed Brown const PetscInt id = 1; 288c4762a1bSJed Brown 289c4762a1bSJed Brown PetscFunctionBeginUser; 2909566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2918b0e23d0SMatthew G. Knepley switch (user->sol) { 292c4762a1bSJed Brown case SOL_QUADRATIC: 2939566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u)); 2948b0e23d0SMatthew G. Knepley exactFuncs[0] = quadratic_u; 2958b0e23d0SMatthew G. Knepley exactFuncs[1] = quadratic_p; 296c4762a1bSJed Brown break; 297c4762a1bSJed Brown case SOL_TRIG: 2989566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 2998b0e23d0SMatthew G. Knepley exactFuncs[0] = trig_u; 3008b0e23d0SMatthew G. Knepley exactFuncs[1] = trig_p; 301c4762a1bSJed Brown break; 30263a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol); 303c4762a1bSJed Brown } 3049566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL)); 3059566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 3069566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 3079566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 3089566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 3099566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL)); 310c4762a1bSJed Brown 3119566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user)); 3129566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user)); 313c4762a1bSJed Brown 3149566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3159566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, user, NULL)); 3168b0e23d0SMatthew G. Knepley 31747bb1945SPatrick Sanan /* Make constant values available to pointwise functions */ 318c4762a1bSJed Brown { 3198b0e23d0SMatthew G. Knepley Parameter *param; 3208b0e23d0SMatthew G. Knepley PetscScalar constants[1]; 321c4762a1bSJed Brown 3229566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 3238b0e23d0SMatthew G. Knepley constants[0] = param->mu; /* dynamic shear viscosity, Pa s */ 3249566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown PetscFunctionReturn(0); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown 3298b0e23d0SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 3308b0e23d0SMatthew G. Knepley { 3318b0e23d0SMatthew G. Knepley PetscInt c; 3328b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 3338b0e23d0SMatthew G. Knepley return 0; 3348b0e23d0SMatthew G. Knepley } 3358b0e23d0SMatthew G. Knepley static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 3368b0e23d0SMatthew G. Knepley { 3378b0e23d0SMatthew G. Knepley PetscInt c; 3388b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 1.0; 3398b0e23d0SMatthew G. Knepley return 0; 3408b0e23d0SMatthew G. Knepley } 3418b0e23d0SMatthew G. Knepley 3428b0e23d0SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 343c4762a1bSJed Brown { 344c4762a1bSJed Brown Vec vec; 345478db826SMatthew G. Knepley PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {zero, one}; 346c4762a1bSJed Brown 347c4762a1bSJed Brown PetscFunctionBeginUser; 34863a3b9bcSJacob Faibussowitsch PetscCheck(origField == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField); 3498b0e23d0SMatthew G. Knepley funcs[field] = one; 3508b0e23d0SMatthew G. Knepley { 3518b0e23d0SMatthew G. Knepley PetscDS ds; 3529566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)ds, NULL, "-ds_view")); 3548b0e23d0SMatthew G. Knepley } 3559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &vec)); 3569566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 3579566063dSJacob Faibussowitsch PetscCall(VecNormalize(vec, NULL)); 3589566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace)); 3599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 360c4762a1bSJed Brown /* New style for field null spaces */ 361c4762a1bSJed Brown { 362c4762a1bSJed Brown PetscObject pressure; 363c4762a1bSJed Brown MatNullSpace nullspacePres; 364c4762a1bSJed Brown 3659566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &pressure)); 3669566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 3679566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres)); 3689566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullspacePres)); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown PetscFunctionReturn(0); 371c4762a1bSJed Brown } 372c4762a1bSJed Brown 3738b0e23d0SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user) 374c4762a1bSJed Brown { 3758b0e23d0SMatthew G. Knepley DM cdm = dm; 3768b0e23d0SMatthew G. Knepley PetscQuadrature q = NULL; 3778b0e23d0SMatthew G. Knepley PetscBool simplex; 37830602db0SMatthew G. Knepley PetscInt dim, Nf = 2, f, Nc[2]; 3798b0e23d0SMatthew G. Knepley const char *name[2] = {"velocity", "pressure"}; 3808b0e23d0SMatthew G. Knepley const char *prefix[2] = {"vel_", "pres_"}; 381c4762a1bSJed Brown 3828b0e23d0SMatthew G. Knepley PetscFunctionBegin; 3839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3849566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 3858b0e23d0SMatthew G. Knepley Nc[0] = dim; 3868b0e23d0SMatthew G. Knepley Nc[1] = 1; 3878b0e23d0SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 3888b0e23d0SMatthew G. Knepley PetscFE fe; 3898b0e23d0SMatthew G. Knepley 3909566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe)); 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name[f])); 3929566063dSJacob Faibussowitsch if (!q) PetscCall(PetscFEGetQuadrature(fe, &q)); 3939566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, q)); 3949566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe)); 3959566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 396c4762a1bSJed Brown } 3979566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 3989566063dSJacob Faibussowitsch PetscCall((*setupEqn)(dm, user)); 3998b0e23d0SMatthew G. Knepley while (cdm) { 4009566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4019566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace)); 4029566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown PetscFunctionReturn(0); 405c4762a1bSJed Brown } 406c4762a1bSJed Brown 407c4762a1bSJed Brown int main(int argc, char **argv) 408c4762a1bSJed Brown { 4098b0e23d0SMatthew G. Knepley SNES snes; 4108b0e23d0SMatthew G. Knepley DM dm; 4118b0e23d0SMatthew G. Knepley Vec u; 4128b0e23d0SMatthew G. Knepley AppCtx user; 413c4762a1bSJed Brown 414327415f7SBarry Smith PetscFunctionBeginUser; 4159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4169566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4179566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4189566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes)); 4199566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4209566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 421c4762a1bSJed Brown 4229566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 4239566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, SetupEqn, &user)); 4249566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 425c4762a1bSJed Brown 4269566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4279566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 4289566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4299566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Solution")); 4318b0e23d0SMatthew G. Knepley { 4328b0e23d0SMatthew G. Knepley Mat J; 4338b0e23d0SMatthew G. Knepley MatNullSpace sp; 434c4762a1bSJed Brown 4359566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes)); 4369566063dSJacob Faibussowitsch PetscCall(CreatePressureNullSpace(dm, 1, 1, &sp)); 4379566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 4389566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J, sp)); 4399566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sp)); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 4419566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-J_view")); 442c4762a1bSJed Brown } 4439566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4469566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4489566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 4499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 450b122ec5aSJacob Faibussowitsch return 0; 451c4762a1bSJed Brown } 452c4762a1bSJed Brown /*TEST 453c4762a1bSJed Brown 454c4762a1bSJed Brown test: 4558b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check 456c4762a1bSJed Brown requires: triangle 4578b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 4588b0e23d0SMatthew G. Knepley 459c4762a1bSJed Brown test: 4608b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check_parallel 4618b0e23d0SMatthew G. Knepley nsize: {{2 3 5}} 462c4762a1bSJed Brown requires: triangle 463e600fa54SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 4648b0e23d0SMatthew G. Knepley 465c4762a1bSJed Brown test: 4668b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check 467c4762a1bSJed Brown requires: ctetgen 46830602db0SMatthew 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 4698b0e23d0SMatthew G. Knepley 470c4762a1bSJed Brown test: 4718b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check_parallel 4728b0e23d0SMatthew G. Knepley nsize: {{2 3 5}} 473c4762a1bSJed Brown requires: ctetgen 47419fa4260SStefano 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 4758b0e23d0SMatthew G. Knepley 476c4762a1bSJed Brown test: 4778b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv 4788b0e23d0SMatthew G. Knepley requires: triangle 4798b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1] 4808b0e23d0SMatthew 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 \ 4818b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 4828b0e23d0SMatthew 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 \ 4838b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 4848b0e23d0SMatthew G. Knepley 485c4762a1bSJed Brown test: 4868b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv_gamg 4878b0e23d0SMatthew G. Knepley requires: triangle 48882894d03SBarry Smith args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 \ 4898b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 4908b0e23d0SMatthew 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 4918b0e23d0SMatthew G. Knepley 492c4762a1bSJed Brown test: 4938b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_conv 4948b0e23d0SMatthew G. Knepley requires: ctetgen !single 4958b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8] 49630602db0SMatthew 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 \ 4978b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 4988b0e23d0SMatthew 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 \ 4998b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5008b0e23d0SMatthew G. Knepley 501c4762a1bSJed Brown test: 5028b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_check 50330602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 5048b0e23d0SMatthew G. Knepley 505c4762a1bSJed Brown test: 5068b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_check 50730602db0SMatthew 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 5088b0e23d0SMatthew G. Knepley 509c4762a1bSJed Brown test: 5108b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_conv 5118b0e23d0SMatthew G. Knepley # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1] 51230602db0SMatthew 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 \ 5138b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5148b0e23d0SMatthew 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 \ 5158b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5168b0e23d0SMatthew G. Knepley 517c4762a1bSJed Brown test: 5188b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_conv 519c4762a1bSJed Brown requires: !single 5208b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4] 52130602db0SMatthew 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 \ 5228b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5238b0e23d0SMatthew 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 \ 5248b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5258b0e23d0SMatthew G. Knepley 526c4762a1bSJed Brown test: 5278b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_check 5288b0e23d0SMatthew G. Knepley requires: triangle 5298b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001 5308b0e23d0SMatthew G. Knepley 5318b0e23d0SMatthew G. Knepley test: 5328b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_check 5338b0e23d0SMatthew G. Knepley requires: ctetgen !single 53430602db0SMatthew 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 5358b0e23d0SMatthew G. Knepley 5368b0e23d0SMatthew G. Knepley test: 5378b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_conv 5388b0e23d0SMatthew G. Knepley requires: triangle 5398b0e23d0SMatthew G. Knepley # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0] 5408b0e23d0SMatthew 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 \ 5418b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5428b0e23d0SMatthew 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 \ 5438b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5448b0e23d0SMatthew G. Knepley 5458b0e23d0SMatthew G. Knepley test: 5468b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_conv 5478b0e23d0SMatthew G. Knepley requires: ctetgen long_runtime 5488b0e23d0SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9] 54930602db0SMatthew 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 \ 5508b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5518b0e23d0SMatthew 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 \ 5528b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5538b0e23d0SMatthew G. Knepley 5548b0e23d0SMatthew G. Knepley test: 5558b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_conv 556c4762a1bSJed Brown requires: !single 5578b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0] 55830602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \ 55982894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \ 5608b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 561bae903cbSmarkadams4 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0 5628b0e23d0SMatthew G. Knepley 563c4762a1bSJed Brown test: 5648b0e23d0SMatthew G. Knepley suffix: 3d_q1_p0_conv 5658b0e23d0SMatthew G. Knepley requires: !single 5668b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0] 56730602db0SMatthew 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 \ 56882894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \ 5698b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 570bae903cbSmarkadams4 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd -fieldsplit_pressure_pc_gamg_aggressive_coarsening 0 5718b0e23d0SMatthew G. Knepley 5728b0e23d0SMatthew G. Knepley # Stokes preconditioners 573c4762a1bSJed Brown # Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix} 574c4762a1bSJed Brown test: 5758b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_block_diagonal 5768b0e23d0SMatthew G. Knepley requires: triangle 5778b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 5788b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 5798b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \ 5808b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi 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 589c4762a1bSJed Brown # Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} 590c4762a1bSJed Brown test: 5918b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_diagonal 5928b0e23d0SMatthew G. Knepley requires: triangle 5938b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 5948b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 5958b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 5968b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \ 5978b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 598c4762a1bSJed Brown # Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 599c4762a1bSJed Brown test: 6008b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_upper 6018b0e23d0SMatthew G. Knepley requires: triangle 6028b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \ 6038b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6048b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \ 6058b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 606c4762a1bSJed Brown # Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 607c4762a1bSJed Brown test: 6088b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_lower 6098b0e23d0SMatthew G. Knepley requires: triangle 6108b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 6118b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6128b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6138b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \ 6148b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 615c4762a1bSJed 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} 616c4762a1bSJed Brown test: 6178b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_full 6188b0e23d0SMatthew G. Knepley requires: triangle 6198b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 6208b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6218b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6228b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \ 6238b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 6248b0e23d0SMatthew G. Knepley # Full Schur + Velocity GMG 6258b0e23d0SMatthew G. Knepley test: 6268b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_gmg_vcycle 6278b0e23d0SMatthew G. Knepley requires: triangle 6288b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 62982894d03SBarry Smith -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \ 6308b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \ 63173f7197eSJed 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 632c4762a1bSJed 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} 633c4762a1bSJed Brown test: 6348b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_simple 635c4762a1bSJed Brown requires: triangle 6368b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 6378b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6388b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 6398b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 6408b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 6418b0e23d0SMatthew 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 642c4762a1bSJed Brown # FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code) 643c4762a1bSJed Brown test: 6448b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_fetidp 645c4762a1bSJed Brown requires: triangle mumps 646c4762a1bSJed Brown nsize: 5 647e600fa54SMatthew 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 \ 6488b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 64982894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \ 6508b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6518b0e23d0SMatthew 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 \ 6528b0e23d0SMatthew 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 653c4762a1bSJed Brown test: 6548b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_fetidp 6558b0e23d0SMatthew G. Knepley requires: mumps 656c4762a1bSJed Brown nsize: 5 657e600fa54SMatthew 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 \ 6588b0e23d0SMatthew G. Knepley -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 6598b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6608b0e23d0SMatthew 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 \ 6618b0e23d0SMatthew 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 662c4762a1bSJed Brown test: 6638b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_fetidp 6648b0e23d0SMatthew G. Knepley requires: ctetgen mumps suitesparse 665c4762a1bSJed Brown nsize: 5 666e600fa54SMatthew 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 \ 6678b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 66882894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-9 \ 6698b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6708b0e23d0SMatthew 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 \ 6718b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \ 6728b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 6738b0e23d0SMatthew 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 \ 6748b0e23d0SMatthew G. Knepley -fetidp_bddelta_pc_factor_mat_ordering_type external \ 6758b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 6768b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 677c4762a1bSJed Brown test: 6788b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_fetidp 679c4762a1bSJed Brown requires: suitesparse 680c4762a1bSJed Brown nsize: 5 681e600fa54SMatthew 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 \ 6828b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 68382894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \ 6848b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6858b0e23d0SMatthew 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 \ 6868b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 6878b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \ 6888b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 6898b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 6908b0e23d0SMatthew G. Knepley # BDDC solvers (these solvers are quite inefficient, they are here to exercise the code) 691c4762a1bSJed Brown test: 6928b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_bddc 6938b0e23d0SMatthew G. Knepley nsize: 2 694c4762a1bSJed Brown requires: triangle !single 695e600fa54SMatthew 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 \ 6968b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6978b0e23d0SMatthew G. Knepley -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 6988b0e23d0SMatthew 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 6998b0e23d0SMatthew G. Knepley # Vanka 700c4762a1bSJed Brown test: 7018b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka 702c4762a1bSJed Brown requires: double !complex 70330602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 7048b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7058b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 706c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 707c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 708c4762a1bSJed Brown test: 7098b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka_denseinv 710c4762a1bSJed Brown requires: double !complex 71130602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 7128b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7138b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 714c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 715c4762a1bSJed Brown -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense 716c4762a1bSJed Brown # Vanka smoother 717c4762a1bSJed Brown test: 7188b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_gmg_vanka 7198b0e23d0SMatthew G. Knepley requires: double !complex 72030602db0SMatthew 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 \ 7218b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7228b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 7238b0e23d0SMatthew G. Knepley -pc_type mg \ 7248b0e23d0SMatthew G. Knepley -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \ 725c4762a1bSJed 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 \ 726c4762a1bSJed Brown -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \ 727c4762a1bSJed Brown -mg_coarse_pc_type svd 728c4762a1bSJed Brown 729c4762a1bSJed Brown TEST*/ 730