1*8e722928SMatthew G. Knepley #include "petscdm.h" 2*8e722928SMatthew G. Knepley #include "petscerror.h" 3*8e722928SMatthew G. Knepley #include "petscis.h" 4*8e722928SMatthew G. Knepley #include "petscsystypes.h" 5*8e722928SMatthew G. Knepley #include "petscvec.h" 6*8e722928SMatthew G. Knepley static char help[] = "Poisson problem in 2d with finite elements and bound constraints.\n\ 7*8e722928SMatthew G. Knepley This example is intended to test VI solvers.\n\n\n"; 8*8e722928SMatthew G. Knepley 9*8e722928SMatthew G. Knepley /* 10*8e722928SMatthew G. Knepley This is the ball obstacle problem, taken from the MS thesis ``Adaptive Mesh Refinement for Variational Inequalities'' 11*8e722928SMatthew G. Knepley by Stefano Fochesatto, University of Alaska Fairbanks, 2025 12*8e722928SMatthew G. Knepley This is the same VI problem as in src/snes/tutorials/ex9.c, which uses DMDA. The example 13*8e722928SMatthew G. Knepley is also documented by Chapter 12 of E. Bueler, "PETSc for Partial Differential Equations", 14*8e722928SMatthew G. Knepley SIAM Press 2021. 15*8e722928SMatthew G. Knepley 16*8e722928SMatthew G. Knepley To visualize the solution, configure with petsc4py, pip install pyvista, and use 17*8e722928SMatthew G. Knepley 18*8e722928SMatthew G. Knepley -potential_view pyvista -view_pyvista_warp 1. 19*8e722928SMatthew G. Knepley 20*8e722928SMatthew G. Knepley To look at the error use 21*8e722928SMatthew G. Knepley 22*8e722928SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 -convest_monitor -convest_error_view pyvista 23*8e722928SMatthew G. Knepley 24*8e722928SMatthew G. Knepley and for the inactive residual and active set use 25*8e722928SMatthew G. Knepley 26*8e722928SMatthew G. Knepley -snes_vi_monitor_residual pyvista -snes_vi_monitor_active pyvista 27*8e722928SMatthew G. Knepley 28*8e722928SMatthew G. Knepley To see the convergence history use 29*8e722928SMatthew G. Knepley 30*8e722928SMatthew G. Knepley -snes_vi_monitor -snes_converged_reason -convest_monitor 31*8e722928SMatthew G. Knepley */ 32*8e722928SMatthew G. Knepley 33*8e722928SMatthew G. Knepley #include <petscdmplex.h> 34*8e722928SMatthew G. Knepley #include <petscsnes.h> 35*8e722928SMatthew G. Knepley #include <petscds.h> 36*8e722928SMatthew G. Knepley #include <petscbag.h> 37*8e722928SMatthew G. Knepley #include <petscconvest.h> 38*8e722928SMatthew G. Knepley 39*8e722928SMatthew G. Knepley typedef enum { 40*8e722928SMatthew G. Knepley OBSTACLE_NONE, 41*8e722928SMatthew G. Knepley OBSTACLE_BALL, 42*8e722928SMatthew G. Knepley NUM_OBSTACLE_TYPES 43*8e722928SMatthew G. Knepley } ObstacleType; 44*8e722928SMatthew G. Knepley const char *obstacleTypes[NUM_OBSTACLE_TYPES + 1] = {"none", "ball", "unknown"}; 45*8e722928SMatthew G. Knepley 46*8e722928SMatthew G. Knepley typedef struct { 47*8e722928SMatthew G. Knepley PetscReal r_0; // Ball radius 48*8e722928SMatthew G. Knepley PetscReal r_free; // Radius of the free boundary for the ball obstacle 49*8e722928SMatthew G. Knepley PetscReal A; // Logarithmic coefficient in exact ball solution 50*8e722928SMatthew G. Knepley PetscReal B; // Constant coefficient in exact ball solution 51*8e722928SMatthew G. Knepley } Parameter; 52*8e722928SMatthew G. Knepley 53*8e722928SMatthew G. Knepley typedef struct { 54*8e722928SMatthew G. Knepley // Problem definition 55*8e722928SMatthew G. Knepley ObstacleType obsType; // Type of obstacle 56*8e722928SMatthew G. Knepley PetscBag bag; // Problem parameters 57*8e722928SMatthew G. Knepley } AppCtx; 58*8e722928SMatthew G. Knepley 59*8e722928SMatthew G. Knepley static PetscErrorCode obstacle_ball(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 60*8e722928SMatthew G. Knepley { 61*8e722928SMatthew G. Knepley Parameter *par = (Parameter *)ctx; 62*8e722928SMatthew G. Knepley const PetscReal r_0 = par->r_0; 63*8e722928SMatthew G. Knepley const PetscReal r = PetscSqrtReal(PetscSqr(x[0]) + PetscSqr(x[1])); 64*8e722928SMatthew G. Knepley const PetscReal psi_0 = PetscSqrtReal(1. - PetscSqr(r_0)); 65*8e722928SMatthew G. Knepley const PetscReal dpsi_0 = -r_0 / psi_0; 66*8e722928SMatthew G. Knepley 67*8e722928SMatthew G. Knepley PetscFunctionBegin; 68*8e722928SMatthew G. Knepley PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ball obstacle is only defined in 2D"); 69*8e722928SMatthew G. Knepley if (r < r_0) u[0] = PetscSqrtReal(1.0 - PetscSqr(r)); 70*8e722928SMatthew G. Knepley else u[0] = psi_0 + dpsi_0 * (r - r_0); 71*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 72*8e722928SMatthew G. Knepley } 73*8e722928SMatthew G. Knepley 74*8e722928SMatthew G. Knepley static PetscErrorCode exactSol_ball(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 75*8e722928SMatthew G. Knepley { 76*8e722928SMatthew G. Knepley Parameter *par = (Parameter *)ctx; 77*8e722928SMatthew G. Knepley const PetscReal r_free = par->r_free; 78*8e722928SMatthew G. Knepley const PetscReal A = par->A; 79*8e722928SMatthew G. Knepley const PetscReal B = par->B; 80*8e722928SMatthew G. Knepley const PetscReal r = PetscSqrtReal(PetscSqr(x[0]) + PetscSqr(x[1])); 81*8e722928SMatthew G. Knepley 82*8e722928SMatthew G. Knepley PetscFunctionBegin; 83*8e722928SMatthew G. Knepley PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Ball obstacle is only defined in 2D"); 84*8e722928SMatthew G. Knepley if (r < r_free) PetscCall(obstacle_ball(dim, time, x, Nc, u, ctx)); 85*8e722928SMatthew G. Knepley else u[0] = -A * PetscLogReal(r) + B; 86*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 87*8e722928SMatthew G. Knepley } 88*8e722928SMatthew G. Knepley 89*8e722928SMatthew 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[]) 90*8e722928SMatthew G. Knepley { 91*8e722928SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f1[d] = u_x[d]; 92*8e722928SMatthew G. Knepley } 93*8e722928SMatthew G. Knepley 94*8e722928SMatthew 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[]) 95*8e722928SMatthew G. Knepley { 96*8e722928SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 97*8e722928SMatthew G. Knepley } 98*8e722928SMatthew G. Knepley 99*8e722928SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 100*8e722928SMatthew G. Knepley { 101*8e722928SMatthew G. Knepley PetscInt obs = OBSTACLE_BALL; 102*8e722928SMatthew G. Knepley 103*8e722928SMatthew G. Knepley options->obsType = OBSTACLE_BALL; 104*8e722928SMatthew G. Knepley 105*8e722928SMatthew G. Knepley PetscFunctionBeginUser; 106*8e722928SMatthew G. Knepley PetscOptionsBegin(comm, "", "Ball Obstacle Problem Options", "DMPLEX"); 107*8e722928SMatthew G. Knepley PetscCall(PetscOptionsEList("-obs_type", "Type of obstacle", "ex34.c", obstacleTypes, NUM_OBSTACLE_TYPES, obstacleTypes[options->obsType], &obs, NULL)); 108*8e722928SMatthew G. Knepley options->obsType = (ObstacleType)obs; 109*8e722928SMatthew G. Knepley PetscOptionsEnd(); 110*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 111*8e722928SMatthew G. Knepley } 112*8e722928SMatthew G. Knepley 113*8e722928SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 114*8e722928SMatthew G. Knepley { 115*8e722928SMatthew G. Knepley PetscBag bag; 116*8e722928SMatthew G. Knepley Parameter *p; 117*8e722928SMatthew G. Knepley 118*8e722928SMatthew G. Knepley PetscFunctionBeginUser; 119*8e722928SMatthew G. Knepley /* setup PETSc parameter bag */ 120*8e722928SMatthew G. Knepley PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 121*8e722928SMatthew G. Knepley PetscCall(PetscBagSetName(ctx->bag, "par", "Obstacle Parameters")); 122*8e722928SMatthew G. Knepley bag = ctx->bag; 123*8e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->r_0, 0.9, "r_0", "Ball radius, m")); 124*8e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->r_free, 0.697965148223374, "r_free", "Ball free boundary radius, m")); 125*8e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->A, 0.680259411891719, "A", "Logarithmic coefficient in exact ball solution")); 126*8e722928SMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, &p->B, 0.471519893402112, "B", "Constant coefficient in exact ball solution")); 127*8e722928SMatthew G. Knepley PetscCall(PetscBagSetFromOptions(bag)); 128*8e722928SMatthew G. Knepley { 129*8e722928SMatthew G. Knepley PetscViewer viewer; 130*8e722928SMatthew G. Knepley PetscViewerFormat format; 131*8e722928SMatthew G. Knepley PetscBool flg; 132*8e722928SMatthew G. Knepley 133*8e722928SMatthew G. Knepley PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 134*8e722928SMatthew G. Knepley if (flg) { 135*8e722928SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 136*8e722928SMatthew G. Knepley PetscCall(PetscBagView(bag, viewer)); 137*8e722928SMatthew G. Knepley PetscCall(PetscViewerFlush(viewer)); 138*8e722928SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 139*8e722928SMatthew G. Knepley PetscCall(PetscViewerDestroy(&viewer)); 140*8e722928SMatthew G. Knepley } 141*8e722928SMatthew G. Knepley } 142*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 143*8e722928SMatthew G. Knepley } 144*8e722928SMatthew G. Knepley 145*8e722928SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 146*8e722928SMatthew G. Knepley { 147*8e722928SMatthew G. Knepley PetscFunctionBeginUser; 148*8e722928SMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 149*8e722928SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 150*8e722928SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 151*8e722928SMatthew G. Knepley PetscCall(DMSetApplicationContext(*dm, user)); 152*8e722928SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 153*8e722928SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(*dm)); 154*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 155*8e722928SMatthew G. Knepley } 156*8e722928SMatthew G. Knepley 157*8e722928SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 158*8e722928SMatthew G. Knepley { 159*8e722928SMatthew G. Knepley PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 160*8e722928SMatthew G. Knepley Parameter *param; 161*8e722928SMatthew G. Knepley PetscDS ds; 162*8e722928SMatthew G. Knepley PetscWeakForm wf; 163*8e722928SMatthew G. Knepley DMLabel label; 164*8e722928SMatthew G. Knepley PetscInt dim, id = 1; 165*8e722928SMatthew G. Knepley void *ctx; 166*8e722928SMatthew G. Knepley 167*8e722928SMatthew G. Knepley PetscFunctionBeginUser; 168*8e722928SMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 169*8e722928SMatthew G. Knepley PetscCall(PetscDSGetWeakForm(ds, &wf)); 170*8e722928SMatthew G. Knepley PetscCall(PetscDSGetSpatialDimension(ds, &dim)); 171*8e722928SMatthew G. Knepley PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 172*8e722928SMatthew G. Knepley switch (user->obsType) { 173*8e722928SMatthew G. Knepley case OBSTACLE_BALL: 174*8e722928SMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u)); 175*8e722928SMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 176*8e722928SMatthew G. Knepley PetscCall(PetscDSSetLowerBound(ds, 0, obstacle_ball, param)); 177*8e722928SMatthew G. Knepley exact = exactSol_ball; 178*8e722928SMatthew G. Knepley break; 179*8e722928SMatthew G. Knepley default: 180*8e722928SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid obstacle type: %s (%d)", obstacleTypes[PetscMin(user->obsType, NUM_OBSTACLE_TYPES)], user->obsType); 181*8e722928SMatthew G. Knepley } 182*8e722928SMatthew G. Knepley PetscCall(PetscBagGetData(user->bag, (void **)&ctx)); 183*8e722928SMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, exact, ctx)); 184*8e722928SMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label)); 185*8e722928SMatthew G. Knepley if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exact, NULL, ctx, NULL)); 186*8e722928SMatthew G. Knepley /* Setup constants */ 187*8e722928SMatthew G. Knepley { 188*8e722928SMatthew G. Knepley PetscScalar constants[4]; 189*8e722928SMatthew G. Knepley 190*8e722928SMatthew G. Knepley constants[0] = param->r_0; // Ball radius 191*8e722928SMatthew G. Knepley constants[1] = param->r_free; // Radius of the free boundary for the ball obstacle 192*8e722928SMatthew G. Knepley constants[2] = param->A; // Logarithmic coefficient in exact ball solution 193*8e722928SMatthew G. Knepley constants[3] = param->B; // Constant coefficient in exact ball solution 194*8e722928SMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 4, constants)); 195*8e722928SMatthew G. Knepley } 196*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 197*8e722928SMatthew G. Knepley } 198*8e722928SMatthew G. Knepley 199*8e722928SMatthew G. Knepley PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 200*8e722928SMatthew G. Knepley { 201*8e722928SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 202*8e722928SMatthew G. Knepley DM cdm = dm; 203*8e722928SMatthew G. Knepley PetscFE fe; 204*8e722928SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 205*8e722928SMatthew G. Knepley DMPolytopeType ct; 206*8e722928SMatthew G. Knepley PetscInt dim, cStart; 207*8e722928SMatthew G. Knepley 208*8e722928SMatthew G. Knepley PetscFunctionBegin; 209*8e722928SMatthew G. Knepley /* Create finite element */ 210*8e722928SMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 211*8e722928SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 212*8e722928SMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 213*8e722928SMatthew G. Knepley PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 214*8e722928SMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, name ? prefix : NULL, -1, &fe)); 215*8e722928SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe, name)); 216*8e722928SMatthew G. Knepley // Set discretization and boundary conditions for each mesh 217*8e722928SMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 218*8e722928SMatthew G. Knepley PetscCall(DMCreateDS(dm)); 219*8e722928SMatthew G. Knepley PetscCall((*setup)(dm, user)); 220*8e722928SMatthew G. Knepley while (cdm) { 221*8e722928SMatthew G. Knepley PetscCall(DMCopyDisc(dm, cdm)); 222*8e722928SMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm)); 223*8e722928SMatthew G. Knepley } 224*8e722928SMatthew G. Knepley PetscCall(PetscFEDestroy(&fe)); 225*8e722928SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 226*8e722928SMatthew G. Knepley } 227*8e722928SMatthew G. Knepley 228*8e722928SMatthew G. Knepley int main(int argc, char **argv) 229*8e722928SMatthew G. Knepley { 230*8e722928SMatthew G. Knepley DM dm; /* Problem specification */ 231*8e722928SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 232*8e722928SMatthew G. Knepley Vec u; /* Solutions */ 233*8e722928SMatthew G. Knepley AppCtx user; /* User-defined work context */ 234*8e722928SMatthew G. Knepley 235*8e722928SMatthew G. Knepley PetscFunctionBeginUser; 236*8e722928SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 237*8e722928SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 238*8e722928SMatthew G. Knepley PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 239*8e722928SMatthew G. Knepley PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 240*8e722928SMatthew G. Knepley /* Primal system */ 241*8e722928SMatthew G. Knepley PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 242*8e722928SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 243*8e722928SMatthew G. Knepley PetscCall(SNESSetDM(snes, dm)); 244*8e722928SMatthew G. Knepley PetscCall(SetupFE(dm, "potential", SetupPrimalProblem, &user)); 245*8e722928SMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &u)); 246*8e722928SMatthew G. Knepley PetscCall(VecSet(u, 0.0)); 247*8e722928SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 248*8e722928SMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 249*8e722928SMatthew G. Knepley PetscCall(DMPlexSetSNESVariableBounds(dm, snes)); 250*8e722928SMatthew G. Knepley 251*8e722928SMatthew G. Knepley PetscCall(SNESSetFromOptions(snes)); 252*8e722928SMatthew G. Knepley PetscCall(DMSNESCheckFromOptions(snes, u)); 253*8e722928SMatthew G. Knepley PetscCall(SNESSolve(snes, NULL, u)); 254*8e722928SMatthew G. Knepley PetscCall(SNESGetSolution(snes, &u)); 255*8e722928SMatthew G. Knepley PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 256*8e722928SMatthew G. Knepley /* Cleanup */ 257*8e722928SMatthew G. Knepley PetscCall(VecDestroy(&u)); 258*8e722928SMatthew G. Knepley PetscCall(SNESDestroy(&snes)); 259*8e722928SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 260*8e722928SMatthew G. Knepley PetscCall(PetscBagDestroy(&user.bag)); 261*8e722928SMatthew G. Knepley PetscCall(PetscFinalize()); 262*8e722928SMatthew G. Knepley return 0; 263*8e722928SMatthew G. Knepley } 264*8e722928SMatthew G. Knepley 265*8e722928SMatthew G. Knepley /*TEST 266*8e722928SMatthew G. Knepley 267*8e722928SMatthew G. Knepley testset: 268*8e722928SMatthew G. Knepley args: -dm_plex_box_lower -2.,-2. -dm_plex_box_upper 2.,2. -dm_plex_box_faces 20,20 \ 269*8e722928SMatthew G. Knepley -potential_petscspace_degree 1 \ 270*8e722928SMatthew G. Knepley -snes_type vinewtonrsls -snes_vi_zero_tolerance 1.0e-12 \ 271*8e722928SMatthew G. Knepley -ksp_type preonly -pc_type lu 272*8e722928SMatthew G. Knepley 273*8e722928SMatthew G. Knepley # Check the exact solution 274*8e722928SMatthew G. Knepley test: 275*8e722928SMatthew G. Knepley suffix: ball_0 276*8e722928SMatthew G. Knepley requires: triangle 277*8e722928SMatthew G. Knepley args: -dmsnes_check 278*8e722928SMatthew G. Knepley 279*8e722928SMatthew G. Knepley # Check convergence 280*8e722928SMatthew G. Knepley test: 281*8e722928SMatthew G. Knepley suffix: ball_1 282*8e722928SMatthew G. Knepley requires: triangle 283*8e722928SMatthew G. Knepley args: -snes_convergence_estimate -convest_num_refine 2 284*8e722928SMatthew G. Knepley 285*8e722928SMatthew G. Knepley # Check different size obstacle 286*8e722928SMatthew G. Knepley test: 287*8e722928SMatthew G. Knepley suffix: ball_2 288*8e722928SMatthew G. Knepley requires: triangle 289*8e722928SMatthew G. Knepley args: -r_0 0.4 290*8e722928SMatthew G. Knepley 291*8e722928SMatthew G. Knepley # Check quadrilateral mesh 292*8e722928SMatthew G. Knepley test: 293*8e722928SMatthew G. Knepley suffix: ball_3 294*8e722928SMatthew G. Knepley args: -dm_plex_simplex 0 -snes_convergence_estimate -convest_num_refine 2 295*8e722928SMatthew G. Knepley 296*8e722928SMatthew G. Knepley TEST*/ 297