18b0e23d0SMatthew G. Knepley static char help[] = "Stokes Problem discretized with finite elements,\n\ 28b0e23d0SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 58b0e23d0SMatthew G. Knepley For the isoviscous Stokes problem, which we discretize using the finite 68b0e23d0SMatthew G. Knepley element method on an unstructured mesh, the weak form equations are 7c4762a1bSJed Brown 88b0e23d0SMatthew G. Knepley < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0 98b0e23d0SMatthew G. Knepley < q, -\nabla\cdot u > = 0 10c4762a1bSJed Brown 11c4762a1bSJed Brown Viewing: 12c4762a1bSJed Brown 13c4762a1bSJed Brown To produce nice output, use 14c4762a1bSJed Brown 158b0e23d0SMatthew G. Knepley -dm_refine 3 -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -snes_view_solution hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append 16c4762a1bSJed Brown 17c4762a1bSJed Brown You can get a LaTeX view of the mesh, with point numbering using 18c4762a1bSJed Brown 19c4762a1bSJed Brown -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0 20c4762a1bSJed Brown 21c4762a1bSJed Brown The data layout can be viewed using 22c4762a1bSJed Brown 23c4762a1bSJed Brown -dm_petscsection_view 24c4762a1bSJed Brown 25c4762a1bSJed Brown Lots of information about the FEM assembly can be printed using 26c4762a1bSJed Brown 278b0e23d0SMatthew G. Knepley -dm_plex_print_fem 3 28c4762a1bSJed Brown */ 29c4762a1bSJed Brown 30c4762a1bSJed Brown #include <petscdmplex.h> 31c4762a1bSJed Brown #include <petscsnes.h> 32c4762a1bSJed Brown #include <petscds.h> 338b0e23d0SMatthew G. Knepley #include <petscbag.h> 34c4762a1bSJed Brown 358b0e23d0SMatthew G. Knepley // TODO: Plot residual by fields after each smoother iterate 36c4762a1bSJed Brown 378b0e23d0SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_TRIG, SOL_UNKNOWN} SolType; 388b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0}; 39c4762a1bSJed Brown 40c4762a1bSJed Brown typedef struct { 418b0e23d0SMatthew G. Knepley PetscScalar mu; /* dynamic shear viscosity */ 428b0e23d0SMatthew G. Knepley } Parameter; 438b0e23d0SMatthew G. Knepley 448b0e23d0SMatthew G. Knepley typedef struct { 458b0e23d0SMatthew G. Knepley PetscBag bag; /* Problem parameters */ 468b0e23d0SMatthew G. Knepley SolType sol; /* MMS solution */ 47c4762a1bSJed Brown } AppCtx; 48c4762a1bSJed Brown 498b0e23d0SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 50c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 51c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 52c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 53c4762a1bSJed Brown { 548b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 558b0e23d0SMatthew G. Knepley const PetscInt Nc = uOff[1]-uOff[0]; 568b0e23d0SMatthew G. Knepley PetscInt c, d; 57c4762a1bSJed Brown 588b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 59c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 608b0e23d0SMatthew G. Knepley f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]); 61c4762a1bSJed Brown } 628b0e23d0SMatthew G. Knepley f1[c*dim+c] -= u[uOff[1]]; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 668b0e23d0SMatthew G. Knepley static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 67c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 68c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 69c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 70c4762a1bSJed Brown { 71c4762a1bSJed Brown PetscInt d; 728b0e23d0SMatthew G. Knepley for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d]; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 758b0e23d0SMatthew G. Knepley static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 76c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 77c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 78c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 79c4762a1bSJed Brown { 80c4762a1bSJed Brown PetscInt d; 818b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */ 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 848b0e23d0SMatthew G. Knepley static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 85c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 86c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 87c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown PetscInt d; 908b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */ 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 938b0e23d0SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 94c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 95c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 96c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 97c4762a1bSJed Brown { 988b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 998b0e23d0SMatthew G. Knepley const PetscInt Nc = uOff[1]-uOff[0]; 1008b0e23d0SMatthew G. Knepley PetscInt c, d; 101c4762a1bSJed Brown 1028b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 103c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 1048b0e23d0SMatthew G. Knepley g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */ 1058b0e23d0SMatthew G. Knepley g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */ 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 1108b0e23d0SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1118b0e23d0SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1128b0e23d0SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1138b0e23d0SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1148b0e23d0SMatthew G. Knepley { 1158b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 1168b0e23d0SMatthew G. Knepley 1178b0e23d0SMatthew G. Knepley g0[0] = 1.0/mu; 1188b0e23d0SMatthew G. Knepley } 1198b0e23d0SMatthew G. Knepley 1208b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution 1218b0e23d0SMatthew G. Knepley 2D: 122c4762a1bSJed Brown 123c4762a1bSJed Brown u = x^2 + y^2 1248b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy 1258b0e23d0SMatthew G. Knepley p = x + y - 1 1268b0e23d0SMatthew G. Knepley f = <1 - 4 mu, 1 - 4 mu> 127c4762a1bSJed Brown 128c4762a1bSJed Brown so that 129c4762a1bSJed Brown 1308b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 4x 4x \ 1318b0e23d0SMatthew G. Knepley \ 4x -4x / 1328b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0 1338b0e23d0SMatthew G. Knepley \nabla \cdot u = 2x - 2x = 0 1348b0e23d0SMatthew G. Knepley 1358b0e23d0SMatthew G. Knepley 3D: 1368b0e23d0SMatthew G. Knepley 1378b0e23d0SMatthew G. Knepley u = 2 x^2 + y^2 + z^2 1388b0e23d0SMatthew G. Knepley v = 2 x^2 - 2xy 1398b0e23d0SMatthew G. Knepley w = 2 x^2 - 2xz 1408b0e23d0SMatthew G. Knepley p = x + y + z - 3/2 1418b0e23d0SMatthew G. Knepley f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> 1428b0e23d0SMatthew G. Knepley 1438b0e23d0SMatthew G. Knepley so that 1448b0e23d0SMatthew G. Knepley 1458b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 8x 4x 4x \ 1468b0e23d0SMatthew G. Knepley | 4x -4x 0 | 1478b0e23d0SMatthew G. Knepley \ 4x 0 -4x / 1488b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <8, 4, 4> - <1, 1, 1> + <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> = 0 1498b0e23d0SMatthew G. Knepley \nabla \cdot u = 4x - 2x - 2x = 0 150c4762a1bSJed Brown */ 1518b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 152c4762a1bSJed Brown { 1538b0e23d0SMatthew G. Knepley PetscInt c; 1548b0e23d0SMatthew G. Knepley 1558b0e23d0SMatthew G. Knepley u[0] = (dim-1)*PetscSqr(x[0]); 1568b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) { 1578b0e23d0SMatthew G. Knepley u[0] += PetscSqr(x[c]); 1588b0e23d0SMatthew G. Knepley u[c] = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c]; 1598b0e23d0SMatthew G. Knepley } 160c4762a1bSJed Brown return 0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1638b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 164c4762a1bSJed Brown { 1658b0e23d0SMatthew G. Knepley PetscInt d; 1668b0e23d0SMatthew G. Knepley 1678b0e23d0SMatthew G. Knepley u[0] = -0.5*dim; 1688b0e23d0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d]; 169c4762a1bSJed Brown return 0; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 1728b0e23d0SMatthew G. Knepley static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 173c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 174c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1758b0e23d0SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 176c4762a1bSJed Brown { 1778b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 1788b0e23d0SMatthew G. Knepley PetscInt d; 1798b0e23d0SMatthew G. Knepley 1808b0e23d0SMatthew G. Knepley f0[0] = (dim-1)*4.0*mu - 1.0; 1818b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 1848b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution 1858b0e23d0SMatthew G. Knepley 2D: 1868b0e23d0SMatthew G. Knepley 1878b0e23d0SMatthew G. Knepley u = sin(pi x) + sin(pi y) 1888b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y 1898b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y) 1908b0e23d0SMatthew G. Knepley f = <2pi cos(2 pi x) + mu pi^2 sin(pi x) + mu pi^2 sin(pi y), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y> 1918b0e23d0SMatthew G. Knepley 1928b0e23d0SMatthew G. Knepley so that 1938b0e23d0SMatthew G. Knepley 1948b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 2pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y \ 1958b0e23d0SMatthew G. Knepley \ pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) / 1968b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <-pi^2 sin(pi x) - pi^2 sin(pi y), pi^3 cos(pi x) y> - <2pi cos(2 pi x), 2pi cos(2 pi y)> + <f_x, f_y> = 0 1978b0e23d0SMatthew G. Knepley \nabla \cdot u = pi cos(pi x) - pi cos(pi x) = 0 1988b0e23d0SMatthew G. Knepley 1998b0e23d0SMatthew G. Knepley 3D: 2008b0e23d0SMatthew G. Knepley 2018b0e23d0SMatthew G. Knepley u = 2 sin(pi x) + sin(pi y) + sin(pi z) 2028b0e23d0SMatthew G. Knepley v = -pi cos(pi x) y 2038b0e23d0SMatthew G. Knepley w = -pi cos(pi x) z 2048b0e23d0SMatthew G. Knepley p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z) 2058b0e23d0SMatthew G. Knepley f = <2pi cos(2 pi x) + mu 2pi^2 sin(pi x) + mu pi^2 sin(pi y) + mu pi^2 sin(pi z), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y, 2pi cos(2 pi z) - mu pi^3 cos(pi x) z> 2068b0e23d0SMatthew G. Knepley 2078b0e23d0SMatthew G. Knepley so that 2088b0e23d0SMatthew G. Knepley 2098b0e23d0SMatthew G. Knepley e(u) = (grad u + grad u^T) = / 4pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y pi cos(pi z) + pi^2 sin(pi x) z \ 2108b0e23d0SMatthew G. Knepley | pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) 0 | 2118b0e23d0SMatthew G. Knepley \ pi cos(pi z) + pi^2 sin(pi x) z 0 -2pi cos(pi x) / 2128b0e23d0SMatthew G. Knepley div mu e(u) - \nabla p + f = mu <-2pi^2 sin(pi x) - pi^2 sin(pi y) - pi^2 sin(pi z), pi^3 cos(pi x) y, pi^3 cos(pi x) z> - <2pi cos(2 pi x), 2pi cos(2 pi y), 2pi cos(2 pi z)> + <f_x, f_y, f_z> = 0 2138b0e23d0SMatthew G. Knepley \nabla \cdot u = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0 2148b0e23d0SMatthew G. Knepley */ 2158b0e23d0SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 216c4762a1bSJed Brown { 2178b0e23d0SMatthew G. Knepley PetscInt c; 2188b0e23d0SMatthew G. Knepley 2198b0e23d0SMatthew G. Knepley u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]); 2208b0e23d0SMatthew G. Knepley for (c = 1; c < Nc; ++c) { 2218b0e23d0SMatthew G. Knepley u[0] += PetscSinReal(PETSC_PI*x[c]); 2228b0e23d0SMatthew G. Knepley u[c] = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c]; 2238b0e23d0SMatthew G. Knepley } 2248b0e23d0SMatthew G. Knepley return 0; 2258b0e23d0SMatthew G. Knepley } 2268b0e23d0SMatthew G. Knepley 2278b0e23d0SMatthew G. Knepley static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 2288b0e23d0SMatthew G. Knepley { 2298b0e23d0SMatthew G. Knepley PetscInt d; 2308b0e23d0SMatthew G. Knepley 2318b0e23d0SMatthew G. Knepley for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]); 2328b0e23d0SMatthew G. Knepley return 0; 2338b0e23d0SMatthew G. Knepley } 2348b0e23d0SMatthew G. Knepley 2358b0e23d0SMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 2368b0e23d0SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 2378b0e23d0SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 2388b0e23d0SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 2398b0e23d0SMatthew G. Knepley { 2408b0e23d0SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 2418b0e23d0SMatthew G. Knepley PetscInt d; 2428b0e23d0SMatthew G. Knepley 2438b0e23d0SMatthew G. Knepley f0[0] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[0]) - (dim-1)*mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[0]); 2448b0e23d0SMatthew G. Knepley for (d = 1; d < dim; ++d) { 2458b0e23d0SMatthew G. Knepley f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]); 2468b0e23d0SMatthew G. Knepley f0[d] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[d]) + mu*PetscPowRealInt(PETSC_PI, 3)*PetscCosReal(PETSC_PI*x[0])*x[d]; 2478b0e23d0SMatthew G. Knepley } 2488b0e23d0SMatthew G. Knepley } 2498b0e23d0SMatthew G. Knepley 2508b0e23d0SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 2518b0e23d0SMatthew G. Knepley { 2528b0e23d0SMatthew G. Knepley PetscInt sol; 253c4762a1bSJed Brown 254c4762a1bSJed Brown PetscFunctionBeginUser; 2558b0e23d0SMatthew G. Knepley options->sol = SOL_QUADRATIC; 256d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX"); 2578b0e23d0SMatthew G. Knepley sol = options->sol; 258dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, PETSC_STATIC_ARRAY_LENGTH(SolTypes)-3, SolTypes[options->sol], &sol, NULL)); 2598b0e23d0SMatthew G. Knepley options->sol = (SolType) sol; 260d0609cedSBarry Smith PetscOptionsEnd(); 261c4762a1bSJed Brown PetscFunctionReturn(0); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 2648b0e23d0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 265c4762a1bSJed Brown { 266c4762a1bSJed Brown PetscFunctionBeginUser; 2679566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2689566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2699566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2709566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 271c4762a1bSJed Brown PetscFunctionReturn(0); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 2748b0e23d0SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 275c4762a1bSJed Brown { 2768b0e23d0SMatthew G. Knepley Parameter *p; 2778b0e23d0SMatthew G. Knepley 2788b0e23d0SMatthew G. Knepley PetscFunctionBeginUser; 2798b0e23d0SMatthew G. Knepley /* setup PETSc parameter bag */ 2809566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag)); 2819566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(ctx->bag, (void **) &p)); 2829566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag, "par", "Stokes Parameters")); 2839566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s")); 2849566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(ctx->bag)); 2858b0e23d0SMatthew G. Knepley { 2868b0e23d0SMatthew G. Knepley PetscViewer viewer; 2878b0e23d0SMatthew G. Knepley PetscViewerFormat format; 2888b0e23d0SMatthew G. Knepley PetscBool flg; 2898b0e23d0SMatthew G. Knepley 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 2918b0e23d0SMatthew G. Knepley if (flg) { 2929566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 2939566063dSJacob Faibussowitsch PetscCall(PetscBagView(ctx->bag, viewer)); 2949566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2959566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 2978b0e23d0SMatthew G. Knepley } 2988b0e23d0SMatthew G. Knepley } 2998b0e23d0SMatthew G. Knepley PetscFunctionReturn(0); 3008b0e23d0SMatthew G. Knepley } 3018b0e23d0SMatthew G. Knepley 3028b0e23d0SMatthew G. Knepley static PetscErrorCode SetupEqn(DM dm, AppCtx *user) 3038b0e23d0SMatthew G. Knepley { 3048b0e23d0SMatthew G. Knepley PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 3058b0e23d0SMatthew G. Knepley PetscDS ds; 30645480ffeSMatthew G. Knepley DMLabel label; 307c4762a1bSJed Brown const PetscInt id = 1; 308c4762a1bSJed Brown 309c4762a1bSJed Brown PetscFunctionBeginUser; 3109566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3118b0e23d0SMatthew G. Knepley switch (user->sol) { 312c4762a1bSJed Brown case SOL_QUADRATIC: 3139566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u)); 3148b0e23d0SMatthew G. Knepley exactFuncs[0] = quadratic_u; 3158b0e23d0SMatthew G. Knepley exactFuncs[1] = quadratic_p; 316c4762a1bSJed Brown break; 317c4762a1bSJed Brown case SOL_TRIG: 3189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 3198b0e23d0SMatthew G. Knepley exactFuncs[0] = trig_u; 3208b0e23d0SMatthew G. Knepley exactFuncs[1] = trig_p; 321c4762a1bSJed Brown break; 32263a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol); 323c4762a1bSJed Brown } 3249566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, NULL)); 3259566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 3269566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 3279566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 3289566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 3299566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL)); 330c4762a1bSJed Brown 3319566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user)); 3329566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user)); 333c4762a1bSJed Brown 3349566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 3359566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 3368b0e23d0SMatthew G. Knepley 33747bb1945SPatrick Sanan /* Make constant values available to pointwise functions */ 338c4762a1bSJed Brown { 3398b0e23d0SMatthew G. Knepley Parameter *param; 3408b0e23d0SMatthew G. Knepley PetscScalar constants[1]; 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 3438b0e23d0SMatthew G. Knepley constants[0] = param->mu; /* dynamic shear viscosity, Pa s */ 3449566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 1, constants)); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown PetscFunctionReturn(0); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 3498b0e23d0SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 3508b0e23d0SMatthew G. Knepley { 3518b0e23d0SMatthew G. Knepley PetscInt c; 3528b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 3538b0e23d0SMatthew G. Knepley return 0; 3548b0e23d0SMatthew G. Knepley } 3558b0e23d0SMatthew G. Knepley static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 3568b0e23d0SMatthew G. Knepley { 3578b0e23d0SMatthew G. Knepley PetscInt c; 3588b0e23d0SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 1.0; 3598b0e23d0SMatthew G. Knepley return 0; 3608b0e23d0SMatthew G. Knepley } 3618b0e23d0SMatthew G. Knepley 3628b0e23d0SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 363c4762a1bSJed Brown { 364c4762a1bSJed Brown Vec vec; 365478db826SMatthew G. Knepley PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one}; 366c4762a1bSJed Brown 367c4762a1bSJed Brown PetscFunctionBeginUser; 36863a3b9bcSJacob Faibussowitsch PetscCheck(origField == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField); 3698b0e23d0SMatthew G. Knepley funcs[field] = one; 3708b0e23d0SMatthew G. Knepley { 3718b0e23d0SMatthew G. Knepley PetscDS ds; 3729566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view")); 3748b0e23d0SMatthew G. Knepley } 3759566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &vec)); 3769566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 3779566063dSJacob Faibussowitsch PetscCall(VecNormalize(vec, NULL)); 3789566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace)); 3799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 380c4762a1bSJed Brown /* New style for field null spaces */ 381c4762a1bSJed Brown { 382c4762a1bSJed Brown PetscObject pressure; 383c4762a1bSJed Brown MatNullSpace nullspacePres; 384c4762a1bSJed Brown 3859566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &pressure)); 3869566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 3879566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres)); 3889566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullspacePres)); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown PetscFunctionReturn(0); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 3938b0e23d0SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user) 394c4762a1bSJed Brown { 3958b0e23d0SMatthew G. Knepley DM cdm = dm; 3968b0e23d0SMatthew G. Knepley PetscQuadrature q = NULL; 3978b0e23d0SMatthew G. Knepley PetscBool simplex; 39830602db0SMatthew G. Knepley PetscInt dim, Nf = 2, f, Nc[2]; 3998b0e23d0SMatthew G. Knepley const char *name[2] = {"velocity", "pressure"}; 4008b0e23d0SMatthew G. Knepley const char *prefix[2] = {"vel_", "pres_"}; 401c4762a1bSJed Brown 4028b0e23d0SMatthew G. Knepley PetscFunctionBegin; 4039566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4049566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 4058b0e23d0SMatthew G. Knepley Nc[0] = dim; 4068b0e23d0SMatthew G. Knepley Nc[1] = 1; 4078b0e23d0SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4088b0e23d0SMatthew G. Knepley PetscFE fe; 4098b0e23d0SMatthew G. Knepley 4109566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe)); 4119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name[f])); 4129566063dSJacob Faibussowitsch if (!q) PetscCall(PetscFEGetQuadrature(fe, &q)); 4139566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(fe, q)); 4149566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, f, NULL, (PetscObject) fe)); 4159566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 416c4762a1bSJed Brown } 4179566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4189566063dSJacob Faibussowitsch PetscCall((*setupEqn)(dm, user)); 4198b0e23d0SMatthew G. Knepley while (cdm) { 4209566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4219566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace)); 4229566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 423c4762a1bSJed Brown } 424c4762a1bSJed Brown PetscFunctionReturn(0); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown int main(int argc, char **argv) 428c4762a1bSJed Brown { 4298b0e23d0SMatthew G. Knepley SNES snes; 4308b0e23d0SMatthew G. Knepley DM dm; 4318b0e23d0SMatthew G. Knepley Vec u; 4328b0e23d0SMatthew G. Knepley AppCtx user; 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4359566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4369566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4379566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject) dm), &snes)); 4389566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4399566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 440c4762a1bSJed Brown 4419566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 4429566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, SetupEqn, &user)); 4439566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4469566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 4479566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4489566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "Solution")); 4508b0e23d0SMatthew G. Knepley { 4518b0e23d0SMatthew G. Knepley Mat J; 4528b0e23d0SMatthew G. Knepley MatNullSpace sp; 453c4762a1bSJed Brown 4549566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes)); 4559566063dSJacob Faibussowitsch PetscCall(CreatePressureNullSpace(dm, 1, 1, &sp)); 4569566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 4579566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J, sp)); 4589566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sp)); 4599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) J, "Jacobian")); 4609566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-J_view")); 461c4762a1bSJed Brown } 4629566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 463c4762a1bSJed Brown 4649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4659566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4679566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 4689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 469b122ec5aSJacob Faibussowitsch return 0; 470c4762a1bSJed Brown } 471c4762a1bSJed Brown /*TEST 472c4762a1bSJed Brown 473c4762a1bSJed Brown test: 4748b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check 475c4762a1bSJed Brown requires: triangle 4768b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 4778b0e23d0SMatthew G. Knepley 478c4762a1bSJed Brown test: 4798b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_check_parallel 4808b0e23d0SMatthew G. Knepley nsize: {{2 3 5}} 481c4762a1bSJed Brown requires: triangle 482e600fa54SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 4838b0e23d0SMatthew G. Knepley 484c4762a1bSJed Brown test: 4858b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check 486c4762a1bSJed Brown requires: ctetgen 48730602db0SMatthew 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 4888b0e23d0SMatthew G. Knepley 489c4762a1bSJed Brown test: 4908b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_check_parallel 4918b0e23d0SMatthew G. Knepley nsize: {{2 3 5}} 492c4762a1bSJed Brown requires: ctetgen 493*19fa4260SStefano 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 4948b0e23d0SMatthew G. Knepley 495c4762a1bSJed Brown test: 4968b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv 4978b0e23d0SMatthew G. Knepley requires: triangle 4988b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1] 4998b0e23d0SMatthew 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 \ 5008b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5018b0e23d0SMatthew 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 \ 5028b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5038b0e23d0SMatthew G. Knepley 504c4762a1bSJed Brown test: 5058b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_conv_gamg 5068b0e23d0SMatthew G. Knepley requires: triangle 50782894d03SBarry Smith args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 \ 5088b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 5098b0e23d0SMatthew 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 5108b0e23d0SMatthew G. Knepley 511c4762a1bSJed Brown test: 5128b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_conv 5138b0e23d0SMatthew G. Knepley requires: ctetgen !single 5148b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8] 51530602db0SMatthew 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 \ 5168b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5178b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 5188b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5198b0e23d0SMatthew G. Knepley 520c4762a1bSJed Brown test: 5218b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_check 52230602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 5238b0e23d0SMatthew G. Knepley 524c4762a1bSJed Brown test: 5258b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_check 52630602db0SMatthew 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 5278b0e23d0SMatthew G. Knepley 528c4762a1bSJed Brown test: 5298b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_conv 5308b0e23d0SMatthew G. Knepley # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1] 53130602db0SMatthew 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 \ 5328b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5338b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 5348b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5358b0e23d0SMatthew G. Knepley 536c4762a1bSJed Brown test: 5378b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_conv 538c4762a1bSJed Brown requires: !single 5398b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4] 54030602db0SMatthew 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 \ 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 545c4762a1bSJed Brown test: 5468b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_check 5478b0e23d0SMatthew G. Knepley requires: triangle 5488b0e23d0SMatthew G. Knepley args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001 5498b0e23d0SMatthew G. Knepley 5508b0e23d0SMatthew G. Knepley test: 5518b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_check 5528b0e23d0SMatthew G. Knepley requires: ctetgen !single 55330602db0SMatthew 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 5548b0e23d0SMatthew G. Knepley 5558b0e23d0SMatthew G. Knepley test: 5568b0e23d0SMatthew G. Knepley suffix: 2d_p3_p2_conv 5578b0e23d0SMatthew G. Knepley requires: triangle 5588b0e23d0SMatthew G. Knepley # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0] 5598b0e23d0SMatthew 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 \ 5608b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5618b0e23d0SMatthew 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 \ 5628b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5638b0e23d0SMatthew G. Knepley 5648b0e23d0SMatthew G. Knepley test: 5658b0e23d0SMatthew G. Knepley suffix: 3d_p3_p2_conv 5668b0e23d0SMatthew G. Knepley requires: ctetgen long_runtime 5678b0e23d0SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9] 56830602db0SMatthew 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 \ 5698b0e23d0SMatthew G. Knepley -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 5708b0e23d0SMatthew 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 \ 5718b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 5728b0e23d0SMatthew G. Knepley 5738b0e23d0SMatthew G. Knepley test: 5748b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_conv 575c4762a1bSJed Brown requires: !single 5768b0e23d0SMatthew G. Knepley # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0] 57730602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \ 57882894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \ 5798b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 5808b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd 5818b0e23d0SMatthew G. Knepley 582c4762a1bSJed Brown test: 5838b0e23d0SMatthew G. Knepley suffix: 3d_q1_p0_conv 5848b0e23d0SMatthew G. Knepley requires: !single 5858b0e23d0SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0] 58630602db0SMatthew 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 \ 58782894d03SBarry Smith -ksp_atol 1e-10 -petscds_jac_pre 0 \ 5888b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 5898b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd 5908b0e23d0SMatthew G. Knepley 5918b0e23d0SMatthew G. Knepley # Stokes preconditioners 592c4762a1bSJed Brown # Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix} 593c4762a1bSJed Brown test: 5948b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_block_diagonal 5958b0e23d0SMatthew G. Knepley requires: triangle 5968b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 5978b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 5988b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \ 5998b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi 600c4762a1bSJed Brown # Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix} 601c4762a1bSJed Brown test: 6028b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_block_triangular 6038b0e23d0SMatthew G. Knepley requires: triangle 6048b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 6058b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6068b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 6078b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi 608c4762a1bSJed Brown # Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} 609c4762a1bSJed Brown test: 6108b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_diagonal 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 diag -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 617c4762a1bSJed Brown # Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 618c4762a1bSJed Brown test: 6198b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_upper 6208b0e23d0SMatthew G. Knepley requires: triangle 6218b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \ 6228b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6238b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \ 6248b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 625c4762a1bSJed Brown # Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 626c4762a1bSJed Brown test: 6278b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_lower 6288b0e23d0SMatthew G. Knepley requires: triangle 6298b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 6308b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6318b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6328b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \ 6338b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 634c4762a1bSJed 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} 635c4762a1bSJed Brown test: 6368b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_schur_full 6378b0e23d0SMatthew G. Knepley requires: triangle 6388b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 6398b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6408b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 6418b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \ 6428b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 6438b0e23d0SMatthew G. Knepley # Full Schur + Velocity GMG 6448b0e23d0SMatthew G. Knepley test: 6458b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_gmg_vcycle 6468b0e23d0SMatthew G. Knepley requires: triangle 6478b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 64882894d03SBarry Smith -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \ 6498b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \ 65073f7197eSJed 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 651c4762a1bSJed 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} 652c4762a1bSJed Brown test: 6538b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_simple 654c4762a1bSJed Brown requires: triangle 6558b0e23d0SMatthew G. Knepley args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 6568b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 6578b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 6588b0e23d0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 6598b0e23d0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 6608b0e23d0SMatthew 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 661c4762a1bSJed Brown # FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code) 662c4762a1bSJed Brown test: 6638b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_fetidp 664c4762a1bSJed Brown requires: triangle mumps 665c4762a1bSJed Brown nsize: 5 666e600fa54SMatthew 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 \ 6678b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 66882894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \ 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 200 -fetidp_fieldsplit_p_pc_type none \ 6718b0e23d0SMatthew 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 672c4762a1bSJed Brown test: 6738b0e23d0SMatthew G. Knepley suffix: 2d_q2_q1_fetidp 6748b0e23d0SMatthew G. Knepley requires: mumps 675c4762a1bSJed Brown nsize: 5 676e600fa54SMatthew 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 \ 6778b0e23d0SMatthew G. Knepley -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 6788b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6798b0e23d0SMatthew 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 \ 6808b0e23d0SMatthew 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 681c4762a1bSJed Brown test: 6828b0e23d0SMatthew G. Knepley suffix: 3d_p2_p1_fetidp 6838b0e23d0SMatthew G. Knepley requires: ctetgen mumps suitesparse 684c4762a1bSJed Brown nsize: 5 685e600fa54SMatthew 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 \ 6868b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 68782894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-9 \ 6888b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 6898b0e23d0SMatthew 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 \ 6908b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \ 6918b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 6928b0e23d0SMatthew 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 \ 6938b0e23d0SMatthew G. Knepley -fetidp_bddelta_pc_factor_mat_ordering_type external \ 6948b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 6958b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 696c4762a1bSJed Brown test: 6978b0e23d0SMatthew G. Knepley suffix: 3d_q2_q1_fetidp 698c4762a1bSJed Brown requires: suitesparse 699c4762a1bSJed Brown nsize: 5 700e600fa54SMatthew 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 \ 7018b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 70282894d03SBarry Smith -ksp_type fetidp -ksp_rtol 1.0e-8 \ 7038b0e23d0SMatthew G. Knepley -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 7048b0e23d0SMatthew 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 \ 7058b0e23d0SMatthew G. Knepley -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 7068b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \ 7078b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 7088b0e23d0SMatthew G. Knepley -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 7098b0e23d0SMatthew G. Knepley # BDDC solvers (these solvers are quite inefficient, they are here to exercise the code) 710c4762a1bSJed Brown test: 7118b0e23d0SMatthew G. Knepley suffix: 2d_p2_p1_bddc 7128b0e23d0SMatthew G. Knepley nsize: 2 713c4762a1bSJed Brown requires: triangle !single 714e600fa54SMatthew 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 \ 7158b0e23d0SMatthew G. Knepley -snes_error_if_not_converged \ 7168b0e23d0SMatthew G. Knepley -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 7178b0e23d0SMatthew 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 7188b0e23d0SMatthew G. Knepley # Vanka 719c4762a1bSJed Brown test: 7208b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka 721c4762a1bSJed Brown requires: double !complex 72230602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 7238b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7248b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 725c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 726c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 727c4762a1bSJed Brown test: 7288b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_vanka_denseinv 729c4762a1bSJed Brown requires: double !complex 73030602db0SMatthew G. Knepley args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 7318b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7328b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 733c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 734c4762a1bSJed Brown -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense 735c4762a1bSJed Brown # Vanka smoother 736c4762a1bSJed Brown test: 7378b0e23d0SMatthew G. Knepley suffix: 2d_q1_p0_gmg_vanka 7388b0e23d0SMatthew G. Knepley requires: double !complex 73930602db0SMatthew 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 \ 7408b0e23d0SMatthew G. Knepley -snes_rtol 1.0e-4 \ 7418b0e23d0SMatthew G. Knepley -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 7428b0e23d0SMatthew G. Knepley -pc_type mg \ 7438b0e23d0SMatthew G. Knepley -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \ 744c4762a1bSJed 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 \ 745c4762a1bSJed Brown -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \ 746c4762a1bSJed Brown -mg_coarse_pc_type svd 747c4762a1bSJed Brown 748c4762a1bSJed Brown TEST*/ 749