1*38be53adSDaniel Finn static char help[] = "'Good Cop' Helmholtz Problem in 2d and 3d with finite elements.\n\ 2*38be53adSDaniel Finn We solve the 'Good Cop' Helmholtz problem in a rectangular\n\ 3*38be53adSDaniel Finn domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4*38be53adSDaniel Finn This example supports automatic convergence estimation\n\ 5*38be53adSDaniel Finn and coarse space adaptivity.\n\n\n"; 6*38be53adSDaniel Finn 7*38be53adSDaniel Finn /* 8*38be53adSDaniel Finn The model problem: 9*38be53adSDaniel Finn Solve "Good Cop" Helmholtz equation on the unit square: (0,1) x (0,1) 10*38be53adSDaniel Finn - \Delta u + u = f, 11*38be53adSDaniel Finn where \Delta = Laplace operator 12*38be53adSDaniel Finn Dirichlet b.c.'s on all sides 13*38be53adSDaniel Finn 14*38be53adSDaniel Finn */ 15*38be53adSDaniel Finn 16*38be53adSDaniel Finn #include <petscdmplex.h> 17*38be53adSDaniel Finn #include <petscsnes.h> 18*38be53adSDaniel Finn #include <petscds.h> 19*38be53adSDaniel Finn #include <petscconvest.h> 20*38be53adSDaniel Finn 21*38be53adSDaniel Finn typedef struct { 22*38be53adSDaniel Finn PetscBool trig; /* Use trig function as exact solution */ 23*38be53adSDaniel Finn } AppCtx; 24*38be53adSDaniel Finn 25*38be53adSDaniel Finn /*For Primal Problem*/ 26*38be53adSDaniel Finn static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 27*38be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 28*38be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 29*38be53adSDaniel Finn PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 30*38be53adSDaniel Finn { 31*38be53adSDaniel Finn PetscInt d; 32*38be53adSDaniel Finn for (d = 0; d < dim; ++d) g0[0] = 1.0; 33*38be53adSDaniel Finn } 34*38be53adSDaniel Finn 35*38be53adSDaniel Finn static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 36*38be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 37*38be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 38*38be53adSDaniel Finn PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 39*38be53adSDaniel Finn { 40*38be53adSDaniel Finn PetscInt d; 41*38be53adSDaniel Finn for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 42*38be53adSDaniel Finn } 43*38be53adSDaniel Finn 44*38be53adSDaniel Finn static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 45*38be53adSDaniel Finn { 46*38be53adSDaniel Finn PetscInt d; 47*38be53adSDaniel Finn *u = 0.0; 48*38be53adSDaniel Finn for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 49*38be53adSDaniel Finn return 0; 50*38be53adSDaniel Finn } 51*38be53adSDaniel Finn 52*38be53adSDaniel Finn static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 53*38be53adSDaniel Finn { 54*38be53adSDaniel Finn PetscInt d; 55*38be53adSDaniel Finn *u = 1.0; 56*38be53adSDaniel Finn for (d = 0; d < dim; ++d) *u += (d+1)*PetscSqr(x[d]); 57*38be53adSDaniel Finn return 0; 58*38be53adSDaniel Finn } 59*38be53adSDaniel Finn 60*38be53adSDaniel Finn static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 61*38be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 62*38be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 63*38be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 64*38be53adSDaniel Finn { 65*38be53adSDaniel Finn PetscInt d; 66*38be53adSDaniel Finn f0[0] += u[0]; 67*38be53adSDaniel Finn for (d = 0; d < dim; ++d) f0[0] -= 4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]) + PetscSinReal(2.0*PETSC_PI*x[d]); 68*38be53adSDaniel Finn } 69*38be53adSDaniel Finn 70*38be53adSDaniel Finn static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 71*38be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 72*38be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 73*38be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 74*38be53adSDaniel Finn { 75*38be53adSDaniel Finn PetscInt d; 76*38be53adSDaniel Finn switch (dim) { 77*38be53adSDaniel Finn case 1: 78*38be53adSDaniel Finn f0[0] = 1.0; 79*38be53adSDaniel Finn break; 80*38be53adSDaniel Finn case 2: 81*38be53adSDaniel Finn f0[0] = 5.0; 82*38be53adSDaniel Finn break; 83*38be53adSDaniel Finn case 3: 84*38be53adSDaniel Finn f0[0] = 11.0; 85*38be53adSDaniel Finn break; 86*38be53adSDaniel Finn default: 87*38be53adSDaniel Finn f0[0] = 5.0; 88*38be53adSDaniel Finn break; 89*38be53adSDaniel Finn } 90*38be53adSDaniel Finn f0[0] += u[0]; 91*38be53adSDaniel Finn for (d = 0; d < dim; ++d) f0[0] -= (d+1)*PetscSqr(x[d]); 92*38be53adSDaniel Finn } 93*38be53adSDaniel Finn 94*38be53adSDaniel Finn static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95*38be53adSDaniel Finn const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96*38be53adSDaniel Finn const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97*38be53adSDaniel Finn PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 98*38be53adSDaniel Finn { 99*38be53adSDaniel Finn PetscInt d; 100*38be53adSDaniel Finn for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 101*38be53adSDaniel Finn } 102*38be53adSDaniel Finn 103*38be53adSDaniel Finn static PetscErrorCode ProcessOptions(DM dm, AppCtx *options) 104*38be53adSDaniel Finn { 105*38be53adSDaniel Finn MPI_Comm comm; 106*38be53adSDaniel Finn PetscInt dim; 107*38be53adSDaniel Finn PetscErrorCode ierr; 108*38be53adSDaniel Finn 109*38be53adSDaniel Finn PetscFunctionBeginUser; 110*38be53adSDaniel Finn ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 111*38be53adSDaniel Finn ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 112*38be53adSDaniel Finn options->trig = PETSC_FALSE; 113*38be53adSDaniel Finn 114*38be53adSDaniel Finn ierr = PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");CHKERRQ(ierr); 115*38be53adSDaniel Finn ierr = PetscOptionsBool("-exact_trig", "Use trigonometric exact solution (better for more complex finite elements)", "ex26.c", options->trig, &options->trig, NULL);CHKERRQ(ierr); 116*38be53adSDaniel Finn ierr = PetscOptionsEnd(); 117*38be53adSDaniel Finn 118*38be53adSDaniel Finn PetscFunctionReturn(0); 119*38be53adSDaniel Finn } 120*38be53adSDaniel Finn 121*38be53adSDaniel Finn static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 122*38be53adSDaniel Finn { 123*38be53adSDaniel Finn PetscErrorCode ierr; 124*38be53adSDaniel Finn 125*38be53adSDaniel Finn PetscFunctionBeginUser; 126*38be53adSDaniel Finn ierr = DMCreate(comm, dm);CHKERRQ(ierr); 127*38be53adSDaniel Finn ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 128*38be53adSDaniel Finn ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 129*38be53adSDaniel Finn 130*38be53adSDaniel Finn ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 131*38be53adSDaniel Finn ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 132*38be53adSDaniel Finn ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 133*38be53adSDaniel Finn 134*38be53adSDaniel Finn PetscFunctionReturn(0); 135*38be53adSDaniel Finn } 136*38be53adSDaniel Finn 137*38be53adSDaniel Finn static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 138*38be53adSDaniel Finn { 139*38be53adSDaniel Finn PetscDS ds; 140*38be53adSDaniel Finn DMLabel label; 141*38be53adSDaniel Finn const PetscInt id = 1; 142*38be53adSDaniel Finn PetscErrorCode ierr; 143*38be53adSDaniel Finn 144*38be53adSDaniel Finn PetscFunctionBeginUser; 145*38be53adSDaniel Finn ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 146*38be53adSDaniel Finn ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 147*38be53adSDaniel Finn if (user->trig) { 148*38be53adSDaniel Finn ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 149*38be53adSDaniel Finn ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu);CHKERRQ(ierr); 150*38be53adSDaniel Finn ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr); 151*38be53adSDaniel Finn ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr); 152*38be53adSDaniel Finn ierr = PetscPrintf(PETSC_COMM_WORLD,"Trig Exact Solution\n");CHKERRQ(ierr); 153*38be53adSDaniel Finn } else { 154*38be53adSDaniel Finn ierr = PetscDSSetResidual(ds, 0, f0_quad_u, f1_u);CHKERRQ(ierr); 155*38be53adSDaniel Finn ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu);CHKERRQ(ierr); 156*38be53adSDaniel Finn ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr); 157*38be53adSDaniel Finn ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL);CHKERRQ(ierr); 158*38be53adSDaniel Finn } 159*38be53adSDaniel Finn PetscFunctionReturn(0); 160*38be53adSDaniel Finn } 161*38be53adSDaniel Finn 162*38be53adSDaniel Finn static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 163*38be53adSDaniel Finn { 164*38be53adSDaniel Finn DM cdm = dm; 165*38be53adSDaniel Finn PetscFE fe; 166*38be53adSDaniel Finn DMPolytopeType ct; 167*38be53adSDaniel Finn PetscBool simplex; 168*38be53adSDaniel Finn PetscInt dim, cStart; 169*38be53adSDaniel Finn char prefix[PETSC_MAX_PATH_LEN]; 170*38be53adSDaniel Finn PetscErrorCode ierr; 171*38be53adSDaniel Finn 172*38be53adSDaniel Finn PetscFunctionBeginUser; 173*38be53adSDaniel Finn ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 174*38be53adSDaniel Finn 175*38be53adSDaniel Finn ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 176*38be53adSDaniel Finn ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 177*38be53adSDaniel Finn simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 178*38be53adSDaniel Finn /* Create finite element */ 179*38be53adSDaniel Finn ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 180*38be53adSDaniel Finn ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 181*38be53adSDaniel Finn ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 182*38be53adSDaniel Finn /* Set discretization and boundary conditions for each mesh */ 183*38be53adSDaniel Finn ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 184*38be53adSDaniel Finn ierr = DMCreateDS(dm);CHKERRQ(ierr); 185*38be53adSDaniel Finn ierr = (*setup)(dm, user);CHKERRQ(ierr); 186*38be53adSDaniel Finn while (cdm) { 187*38be53adSDaniel Finn ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 188*38be53adSDaniel Finn ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 189*38be53adSDaniel Finn } 190*38be53adSDaniel Finn ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 191*38be53adSDaniel Finn PetscFunctionReturn(0); 192*38be53adSDaniel Finn } 193*38be53adSDaniel Finn 194*38be53adSDaniel Finn int main(int argc, char **argv) 195*38be53adSDaniel Finn { 196*38be53adSDaniel Finn DM dm; /* Problem specification */ 197*38be53adSDaniel Finn PetscDS ds; 198*38be53adSDaniel Finn SNES snes; /* Nonlinear solver */ 199*38be53adSDaniel Finn Vec u; /* Solutions */ 200*38be53adSDaniel Finn AppCtx user; /* User-defined work context */ 201*38be53adSDaniel Finn PetscErrorCode ierr; 202*38be53adSDaniel Finn 203*38be53adSDaniel Finn ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 204*38be53adSDaniel Finn /* Primal system */ 205*38be53adSDaniel Finn ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 206*38be53adSDaniel Finn ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 207*38be53adSDaniel Finn ierr = ProcessOptions(dm, &user);CHKERRQ(ierr); 208*38be53adSDaniel Finn ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 209*38be53adSDaniel Finn ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 210*38be53adSDaniel Finn ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 211*38be53adSDaniel Finn ierr = VecSet(u, 0.0);CHKERRQ(ierr); 212*38be53adSDaniel Finn ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 213*38be53adSDaniel Finn ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 214*38be53adSDaniel Finn ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 215*38be53adSDaniel Finn ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 216*38be53adSDaniel Finn 217*38be53adSDaniel Finn /*Looking for field error*/ 218*38be53adSDaniel Finn PetscInt Nfields; 219*38be53adSDaniel Finn ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 220*38be53adSDaniel Finn ierr = PetscDSGetNumFields(ds, &Nfields);CHKERRQ(ierr); 221*38be53adSDaniel Finn ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 222*38be53adSDaniel Finn ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 223*38be53adSDaniel Finn ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 224*38be53adSDaniel Finn 225*38be53adSDaniel Finn /* Cleanup */ 226*38be53adSDaniel Finn ierr = VecDestroy(&u);CHKERRQ(ierr); 227*38be53adSDaniel Finn ierr = SNESDestroy(&snes);CHKERRQ(ierr); 228*38be53adSDaniel Finn ierr = DMDestroy(&dm);CHKERRQ(ierr); 229*38be53adSDaniel Finn ierr = PetscFinalize(); 230*38be53adSDaniel Finn return ierr; 231*38be53adSDaniel Finn } 232*38be53adSDaniel Finn 233*38be53adSDaniel Finn /*TEST 234*38be53adSDaniel Finn test: 235*38be53adSDaniel Finn # L_2 convergence rate: 1.9 236*38be53adSDaniel Finn suffix: 2d_p1_conv 237*38be53adSDaniel Finn requires: triangle 238*38be53adSDaniel Finn args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 239*38be53adSDaniel Finn test: 240*38be53adSDaniel Finn # L_2 convergence rate: 1.9 241*38be53adSDaniel Finn suffix: 2d_q1_conv 242*38be53adSDaniel Finn args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu 243*38be53adSDaniel Finn test: 244*38be53adSDaniel Finn # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 245*38be53adSDaniel Finn suffix: 3d_p1_conv 246*38be53adSDaniel Finn requires: ctetgen 247*38be53adSDaniel Finn args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 248*38be53adSDaniel Finn test: 249*38be53adSDaniel Finn # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 250*38be53adSDaniel Finn suffix: 3d_q1_conv 251*38be53adSDaniel Finn args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 252*38be53adSDaniel Finn test: 253*38be53adSDaniel Finn # L_2 convergence rate: 1.9 254*38be53adSDaniel Finn suffix: 2d_p1_trig_conv 255*38be53adSDaniel Finn requires: triangle 256*38be53adSDaniel Finn args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 257*38be53adSDaniel Finn test: 258*38be53adSDaniel Finn # L_2 convergence rate: 1.9 259*38be53adSDaniel Finn suffix: 2d_q1_trig_conv 260*38be53adSDaniel Finn args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig 261*38be53adSDaniel Finn test: 262*38be53adSDaniel Finn # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5 263*38be53adSDaniel Finn suffix: 3d_p1_trig_conv 264*38be53adSDaniel Finn requires: ctetgen 265*38be53adSDaniel Finn args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig 266*38be53adSDaniel Finn test: 267*38be53adSDaniel Finn # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2 268*38be53adSDaniel Finn suffix: 3d_q1_trig_conv 269*38be53adSDaniel Finn args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu -exact_trig 270*38be53adSDaniel Finn test: 271*38be53adSDaniel Finn suffix: 2d_p1_gmg_vcycle 272*38be53adSDaniel Finn requires: triangle 273*38be53adSDaniel Finn args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 274*38be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 275*38be53adSDaniel Finn -mg_levels_ksp_max_it 1 \ 276*38be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 277*38be53adSDaniel Finn test: 278*38be53adSDaniel Finn suffix: 2d_p1_gmg_fcycle 279*38be53adSDaniel Finn requires: triangle 280*38be53adSDaniel Finn args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 281*38be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \ 282*38be53adSDaniel Finn -mg_levels_ksp_max_it 2 \ 283*38be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 284*38be53adSDaniel Finn test: 285*38be53adSDaniel Finn suffix: 2d_p1_gmg_vcycle_trig 286*38be53adSDaniel Finn requires: triangle 287*38be53adSDaniel Finn args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 288*38be53adSDaniel Finn -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \ 289*38be53adSDaniel Finn -mg_levels_ksp_max_it 1 \ 290*38be53adSDaniel Finn -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor 291*38be53adSDaniel Finn TEST*/ 292