18e722928SMatthew G. Knepley static char help[] = "Poisson problem in 2d with finite elements and bound constraints.\n\
28e722928SMatthew G. Knepley This example is intended to test VI solvers.\n\n\n";
38e722928SMatthew G. Knepley
48e722928SMatthew G. Knepley /*
58e722928SMatthew G. Knepley This is the ball obstacle problem, taken from the MS thesis ``Adaptive Mesh Refinement for Variational Inequalities''
68e722928SMatthew G. Knepley by Stefano Fochesatto, University of Alaska Fairbanks, 2025
78e722928SMatthew G. Knepley This is the same VI problem as in src/snes/tutorials/ex9.c, which uses DMDA. The example
88e722928SMatthew G. Knepley is also documented by Chapter 12 of E. Bueler, "PETSc for Partial Differential Equations",
98e722928SMatthew G. Knepley SIAM Press 2021.
108e722928SMatthew G. Knepley
118e722928SMatthew G. Knepley To visualize the solution, configure with petsc4py, pip install pyvista, and use
128e722928SMatthew G. Knepley
138e722928SMatthew G. Knepley -potential_view pyvista -view_pyvista_warp 1.
148e722928SMatthew G. Knepley
158e722928SMatthew G. Knepley To look at the error use
168e722928SMatthew G. Knepley
178e722928SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 -convest_monitor -convest_error_view pyvista
188e722928SMatthew G. Knepley
198e722928SMatthew G. Knepley and for the inactive residual and active set use
208e722928SMatthew G. Knepley
218e722928SMatthew G. Knepley -snes_vi_monitor_residual pyvista -snes_vi_monitor_active pyvista
228e722928SMatthew G. Knepley
238e722928SMatthew G. Knepley To see the convergence history use
248e722928SMatthew G. Knepley
258e722928SMatthew G. Knepley -snes_vi_monitor -snes_converged_reason -convest_monitor
268e722928SMatthew G. Knepley */
278e722928SMatthew G. Knepley
288e722928SMatthew G. Knepley #include <petscdmplex.h>
298e722928SMatthew G. Knepley #include <petscsnes.h>
308e722928SMatthew G. Knepley #include <petscds.h>
318e722928SMatthew G. Knepley #include <petscbag.h>
328e722928SMatthew G. Knepley #include <petscconvest.h>
338e722928SMatthew G. Knepley
348e722928SMatthew G. Knepley typedef enum {
358e722928SMatthew G. Knepley OBSTACLE_NONE,
368e722928SMatthew G. Knepley OBSTACLE_BALL,
378e722928SMatthew G. Knepley NUM_OBSTACLE_TYPES
388e722928SMatthew G. Knepley } ObstacleType;
398e722928SMatthew G. Knepley const char *obstacleTypes[NUM_OBSTACLE_TYPES + 1] = {"none", "ball", "unknown"};
408e722928SMatthew G. Knepley
418e722928SMatthew G. Knepley typedef struct {
428e722928SMatthew G. Knepley PetscReal r_0; // Ball radius
438e722928SMatthew G. Knepley PetscReal r_free; // Radius of the free boundary for the ball obstacle
448e722928SMatthew G. Knepley PetscReal A; // Logarithmic coefficient in exact ball solution
458e722928SMatthew G. Knepley PetscReal B; // Constant coefficient in exact ball solution
468e722928SMatthew G. Knepley } Parameter;
478e722928SMatthew G. Knepley
488e722928SMatthew G. Knepley typedef struct {
498e722928SMatthew G. Knepley // Problem definition
508e722928SMatthew G. Knepley ObstacleType obsType; // Type of obstacle
518e722928SMatthew G. Knepley PetscBag bag; // Problem parameters
528e722928SMatthew G. Knepley } AppCtx;
538e722928SMatthew G. Knepley
obstacle_ball(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)54*2a8381b2SBarry Smith static PetscErrorCode obstacle_ball(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
558e722928SMatthew G. Knepley {
568e722928SMatthew G. Knepley Parameter *par = (Parameter *)ctx;
578e722928SMatthew G. Knepley const PetscReal r_0 = par->r_0;
588e722928SMatthew G. Knepley const PetscReal r = PetscSqrtReal(PetscSqr(x[0]) + PetscSqr(x[1]));
598e722928SMatthew G. Knepley const PetscReal psi_0 = PetscSqrtReal(1. - PetscSqr(r_0));
608e722928SMatthew G. Knepley const PetscReal dpsi_0 = -r_0 / psi_0;
618e722928SMatthew G. Knepley
628e722928SMatthew G. Knepley PetscFunctionBegin;
638e722928SMatthew G. Knepley PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ball obstacle is only defined in 2D");
648e722928SMatthew G. Knepley if (r < r_0) u[0] = PetscSqrtReal(1.0 - PetscSqr(r));
658e722928SMatthew G. Knepley else u[0] = psi_0 + dpsi_0 * (r - r_0);
668e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
678e722928SMatthew G. Knepley }
688e722928SMatthew G. Knepley
exactSol_ball(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)69*2a8381b2SBarry Smith static PetscErrorCode exactSol_ball(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
708e722928SMatthew G. Knepley {
718e722928SMatthew G. Knepley Parameter *par = (Parameter *)ctx;
728e722928SMatthew G. Knepley const PetscReal r_free = par->r_free;
738e722928SMatthew G. Knepley const PetscReal A = par->A;
748e722928SMatthew G. Knepley const PetscReal B = par->B;
758e722928SMatthew G. Knepley const PetscReal r = PetscSqrtReal(PetscSqr(x[0]) + PetscSqr(x[1]));
768e722928SMatthew G. Knepley
778e722928SMatthew G. Knepley PetscFunctionBegin;
788e722928SMatthew G. Knepley PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ball obstacle is only defined in 2D");
798e722928SMatthew G. Knepley if (r < r_free) PetscCall(obstacle_ball(dim, time, x, Nc, u, ctx));
808e722928SMatthew G. Knepley else u[0] = -A * PetscLogReal(r) + B;
818e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
828e722928SMatthew G. Knepley }
838e722928SMatthew G. Knepley
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[])848e722928SMatthew G. Knepley 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[])
858e722928SMatthew G. Knepley {
868e722928SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f1[d] = u_x[d];
878e722928SMatthew G. Knepley }
888e722928SMatthew G. Knepley
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[])898e722928SMatthew G. Knepley 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[])
908e722928SMatthew G. Knepley {
918e722928SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
928e722928SMatthew G. Knepley }
938e722928SMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)948e722928SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
958e722928SMatthew G. Knepley {
968e722928SMatthew G. Knepley PetscInt obs = OBSTACLE_BALL;
978e722928SMatthew G. Knepley
988e722928SMatthew G. Knepley options->obsType = OBSTACLE_BALL;
998e722928SMatthew G. Knepley
1008e722928SMatthew G. Knepley PetscFunctionBeginUser;
1018e722928SMatthew G. Knepley PetscOptionsBegin(comm, "", "Ball Obstacle Problem Options", "DMPLEX");
1028e722928SMatthew G. Knepley PetscCall(PetscOptionsEList("-obs_type", "Type of obstacle", "ex34.c", obstacleTypes, NUM_OBSTACLE_TYPES, obstacleTypes[options->obsType], &obs, NULL));
1038e722928SMatthew G. Knepley options->obsType = (ObstacleType)obs;
1048e722928SMatthew G. Knepley PetscOptionsEnd();
1058e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1068e722928SMatthew G. Knepley }
1078e722928SMatthew G. Knepley
SetupParameters(MPI_Comm comm,AppCtx * ctx)1088e722928SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
1098e722928SMatthew G. Knepley {
1108e722928SMatthew G. Knepley PetscBag bag;
1118e722928SMatthew G. Knepley Parameter *p;
1128e722928SMatthew G. Knepley
1138e722928SMatthew G. Knepley PetscFunctionBeginUser;
1148e722928SMatthew G. Knepley /* setup PETSc parameter bag */
115*2a8381b2SBarry Smith PetscCall(PetscBagGetData(ctx->bag, &p));
1168e722928SMatthew G. Knepley PetscCall(PetscBagSetName(ctx->bag, "par", "Obstacle Parameters"));
1178e722928SMatthew G. Knepley bag = ctx->bag;
1188e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->r_0, 0.9, "r_0", "Ball radius, m"));
1198e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->r_free, 0.697965148223374, "r_free", "Ball free boundary radius, m"));
1208e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->A, 0.680259411891719, "A", "Logarithmic coefficient in exact ball solution"));
1218e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->B, 0.471519893402112, "B", "Constant coefficient in exact ball solution"));
1228e722928SMatthew G. Knepley PetscCall(PetscBagSetFromOptions(bag));
1238e722928SMatthew G. Knepley {
1248e722928SMatthew G. Knepley PetscViewer viewer;
1258e722928SMatthew G. Knepley PetscViewerFormat format;
1268e722928SMatthew G. Knepley PetscBool flg;
1278e722928SMatthew G. Knepley
1288e722928SMatthew G. Knepley PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
1298e722928SMatthew G. Knepley if (flg) {
1308e722928SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format));
1318e722928SMatthew G. Knepley PetscCall(PetscBagView(bag, viewer));
1328e722928SMatthew G. Knepley PetscCall(PetscViewerFlush(viewer));
1338e722928SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer));
1348e722928SMatthew G. Knepley PetscCall(PetscViewerDestroy(&viewer));
1358e722928SMatthew G. Knepley }
1368e722928SMatthew G. Knepley }
1378e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1388e722928SMatthew G. Knepley }
1398e722928SMatthew G. Knepley
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)1408e722928SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
1418e722928SMatthew G. Knepley {
1428e722928SMatthew G. Knepley PetscFunctionBeginUser;
1438e722928SMatthew G. Knepley PetscCall(DMCreate(comm, dm));
1448e722928SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX));
1458e722928SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm));
1468e722928SMatthew G. Knepley PetscCall(DMSetApplicationContext(*dm, user));
1478e722928SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1488e722928SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(*dm));
1498e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1508e722928SMatthew G. Knepley }
1518e722928SMatthew G. Knepley
SetupPrimalProblem(DM dm,AppCtx * user)1528e722928SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
1538e722928SMatthew G. Knepley {
1548e722928SMatthew G. Knepley PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
1558e722928SMatthew G. Knepley Parameter *param;
1568e722928SMatthew G. Knepley PetscDS ds;
1578e722928SMatthew G. Knepley PetscWeakForm wf;
1588e722928SMatthew G. Knepley DMLabel label;
1598e722928SMatthew G. Knepley PetscInt dim, id = 1;
1608e722928SMatthew G. Knepley void *ctx;
1618e722928SMatthew G. Knepley
1628e722928SMatthew G. Knepley PetscFunctionBeginUser;
1638e722928SMatthew G. Knepley PetscCall(DMGetDS(dm, &ds));
1648e722928SMatthew G. Knepley PetscCall(PetscDSGetWeakForm(ds, &wf));
1658e722928SMatthew G. Knepley PetscCall(PetscDSGetSpatialDimension(ds, &dim));
166*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m));
1678e722928SMatthew G. Knepley switch (user->obsType) {
1688e722928SMatthew G. Knepley case OBSTACLE_BALL:
1698e722928SMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u));
1708e722928SMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1718e722928SMatthew G. Knepley PetscCall(PetscDSSetLowerBound(ds, 0, obstacle_ball, param));
1728e722928SMatthew G. Knepley exact = exactSol_ball;
1738e722928SMatthew G. Knepley break;
1748e722928SMatthew G. Knepley default:
1758e722928SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid obstacle type: %s (%d)", obstacleTypes[PetscMin(user->obsType, NUM_OBSTACLE_TYPES)], user->obsType);
1768e722928SMatthew G. Knepley }
177*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, &ctx));
1788e722928SMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, exact, ctx));
1798e722928SMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label));
18057d50842SBarry Smith if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)exact, NULL, ctx, NULL));
1818e722928SMatthew G. Knepley /* Setup constants */
1828e722928SMatthew G. Knepley {
1838e722928SMatthew G. Knepley PetscScalar constants[4];
1848e722928SMatthew G. Knepley
1858e722928SMatthew G. Knepley constants[0] = param->r_0; // Ball radius
1868e722928SMatthew G. Knepley constants[1] = param->r_free; // Radius of the free boundary for the ball obstacle
1878e722928SMatthew G. Knepley constants[2] = param->A; // Logarithmic coefficient in exact ball solution
1888e722928SMatthew G. Knepley constants[3] = param->B; // Constant coefficient in exact ball solution
1898e722928SMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 4, constants));
1908e722928SMatthew G. Knepley }
1918e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1928e722928SMatthew G. Knepley }
1938e722928SMatthew G. Knepley
SetupFE(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),PetscCtx ctx)194*2a8381b2SBarry Smith PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), PetscCtx ctx)
1958e722928SMatthew G. Knepley {
1968e722928SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx;
1978e722928SMatthew G. Knepley DM cdm = dm;
1988e722928SMatthew G. Knepley PetscFE fe;
1998e722928SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN];
2008e722928SMatthew G. Knepley DMPolytopeType ct;
2018e722928SMatthew G. Knepley PetscInt dim, cStart;
2028e722928SMatthew G. Knepley
2038e722928SMatthew G. Knepley PetscFunctionBegin;
2048e722928SMatthew G. Knepley /* Create finite element */
2058e722928SMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim));
2068e722928SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2078e722928SMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct));
2088e722928SMatthew G. Knepley PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
2098e722928SMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, name ? prefix : NULL, -1, &fe));
2108e722928SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe, name));
2118e722928SMatthew G. Knepley // Set discretization and boundary conditions for each mesh
2128e722928SMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2138e722928SMatthew G. Knepley PetscCall(DMCreateDS(dm));
2148e722928SMatthew G. Knepley PetscCall((*setup)(dm, user));
2158e722928SMatthew G. Knepley while (cdm) {
2168e722928SMatthew G. Knepley PetscCall(DMCopyDisc(dm, cdm));
2178e722928SMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm));
2188e722928SMatthew G. Knepley }
2198e722928SMatthew G. Knepley PetscCall(PetscFEDestroy(&fe));
2208e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
2218e722928SMatthew G. Knepley }
2228e722928SMatthew G. Knepley
main(int argc,char ** argv)2238e722928SMatthew G. Knepley int main(int argc, char **argv)
2248e722928SMatthew G. Knepley {
2258e722928SMatthew G. Knepley DM dm; /* Problem specification */
2268e722928SMatthew G. Knepley SNES snes; /* Nonlinear solver */
2278e722928SMatthew G. Knepley Vec u; /* Solutions */
2288e722928SMatthew G. Knepley AppCtx user; /* User-defined work context */
2298e722928SMatthew G. Knepley
2308e722928SMatthew G. Knepley PetscFunctionBeginUser;
2318e722928SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2328e722928SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2338e722928SMatthew G. Knepley PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag));
2348e722928SMatthew G. Knepley PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
2358e722928SMatthew G. Knepley /* Primal system */
2368e722928SMatthew G. Knepley PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
2378e722928SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
2388e722928SMatthew G. Knepley PetscCall(SNESSetDM(snes, dm));
2398e722928SMatthew G. Knepley PetscCall(SetupFE(dm, "potential", SetupPrimalProblem, &user));
2408e722928SMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &u));
2418e722928SMatthew G. Knepley PetscCall(VecSet(u, 0.0));
2428e722928SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
2438e722928SMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
2448e722928SMatthew G. Knepley PetscCall(DMPlexSetSNESVariableBounds(dm, snes));
2458e722928SMatthew G. Knepley
2468e722928SMatthew G. Knepley PetscCall(SNESSetFromOptions(snes));
2478e722928SMatthew G. Knepley PetscCall(DMSNESCheckFromOptions(snes, u));
2488e722928SMatthew G. Knepley PetscCall(SNESSolve(snes, NULL, u));
2498e722928SMatthew G. Knepley PetscCall(SNESGetSolution(snes, &u));
2508e722928SMatthew G. Knepley PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
2518e722928SMatthew G. Knepley /* Cleanup */
2528e722928SMatthew G. Knepley PetscCall(VecDestroy(&u));
2538e722928SMatthew G. Knepley PetscCall(SNESDestroy(&snes));
2548e722928SMatthew G. Knepley PetscCall(DMDestroy(&dm));
2558e722928SMatthew G. Knepley PetscCall(PetscBagDestroy(&user.bag));
2568e722928SMatthew G. Knepley PetscCall(PetscFinalize());
2578e722928SMatthew G. Knepley return 0;
2588e722928SMatthew G. Knepley }
2598e722928SMatthew G. Knepley
2608e722928SMatthew G. Knepley /*TEST
2618e722928SMatthew G. Knepley
2628e722928SMatthew G. Knepley testset:
2638e722928SMatthew G. Knepley args: -dm_plex_box_lower -2.,-2. -dm_plex_box_upper 2.,2. -dm_plex_box_faces 20,20 \
2648e722928SMatthew G. Knepley -potential_petscspace_degree 1 \
2658e722928SMatthew G. Knepley -snes_type vinewtonrsls -snes_vi_zero_tolerance 1.0e-12 \
2668e722928SMatthew G. Knepley -ksp_type preonly -pc_type lu
2678e722928SMatthew G. Knepley
2688e722928SMatthew G. Knepley # Check the exact solution
2698e722928SMatthew G. Knepley test:
2708e722928SMatthew G. Knepley suffix: ball_0
2718e722928SMatthew G. Knepley requires: triangle
2728e722928SMatthew G. Knepley args: -dmsnes_check
2738e722928SMatthew G. Knepley
2748e722928SMatthew G. Knepley # Check convergence
2758e722928SMatthew G. Knepley test:
2768e722928SMatthew G. Knepley suffix: ball_1
2778e722928SMatthew G. Knepley requires: triangle
2788e722928SMatthew G. Knepley args: -snes_convergence_estimate -convest_num_refine 2
2798e722928SMatthew G. Knepley
2808e722928SMatthew G. Knepley # Check different size obstacle
2818e722928SMatthew G. Knepley test:
2828e722928SMatthew G. Knepley suffix: ball_2
2838e722928SMatthew G. Knepley requires: triangle
2848e722928SMatthew G. Knepley args: -r_0 0.4
285e0008caeSPierre Jolivet output_file: output/empty.out
2868e722928SMatthew G. Knepley
2878e722928SMatthew G. Knepley # Check quadrilateral mesh
2888e722928SMatthew G. Knepley test:
2898e722928SMatthew G. Knepley suffix: ball_3
2908e722928SMatthew G. Knepley args: -dm_plex_simplex 0 -snes_convergence_estimate -convest_num_refine 2
2918e722928SMatthew G. Knepley
2928e722928SMatthew G. Knepley TEST*/
293