xref: /petsc/src/snes/tutorials/ex26.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc) !
138be53adSDaniel Finn static char help[] = "'Good Cop' Helmholtz Problem in 2d and 3d with finite elements.\n\
238be53adSDaniel Finn We solve the 'Good Cop' Helmholtz problem in a rectangular\n\
338be53adSDaniel Finn domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
438be53adSDaniel Finn This example supports automatic convergence estimation\n\
538be53adSDaniel Finn and coarse space adaptivity.\n\n\n";
638be53adSDaniel Finn 
738be53adSDaniel Finn /*
838be53adSDaniel Finn    The model problem:
938be53adSDaniel Finn       Solve "Good Cop" Helmholtz equation on the unit square: (0,1) x (0,1)
1038be53adSDaniel Finn           - \Delta u + u = f,
1138be53adSDaniel Finn            where \Delta = Laplace operator
1238be53adSDaniel Finn       Dirichlet b.c.'s on all sides
1338be53adSDaniel Finn 
1438be53adSDaniel Finn */
1538be53adSDaniel Finn 
1638be53adSDaniel Finn #include <petscdmplex.h>
1738be53adSDaniel Finn #include <petscsnes.h>
1838be53adSDaniel Finn #include <petscds.h>
1938be53adSDaniel Finn #include <petscconvest.h>
2038be53adSDaniel Finn 
2138be53adSDaniel Finn typedef struct {
2238be53adSDaniel Finn   PetscBool trig; /* Use trig function as exact solution */
2338be53adSDaniel Finn } AppCtx;
2438be53adSDaniel Finn 
2538be53adSDaniel Finn /*For Primal Problem*/
2638be53adSDaniel Finn static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2738be53adSDaniel Finn                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2838be53adSDaniel Finn                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
2938be53adSDaniel Finn                    PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
3038be53adSDaniel Finn {
3138be53adSDaniel Finn   PetscInt d;
3238be53adSDaniel Finn   for (d = 0; d < dim; ++d) g0[0] = 1.0;
3338be53adSDaniel Finn }
3438be53adSDaniel Finn 
3538be53adSDaniel Finn static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
3638be53adSDaniel Finn                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
3738be53adSDaniel Finn                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
3838be53adSDaniel Finn                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
3938be53adSDaniel Finn {
4038be53adSDaniel Finn   PetscInt d;
4138be53adSDaniel Finn   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
4238be53adSDaniel Finn }
4338be53adSDaniel Finn 
4438be53adSDaniel Finn static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
4538be53adSDaniel Finn {
4638be53adSDaniel Finn   PetscInt d;
4738be53adSDaniel Finn   *u = 0.0;
4838be53adSDaniel Finn   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
4938be53adSDaniel Finn   return 0;
5038be53adSDaniel Finn }
5138be53adSDaniel Finn 
5238be53adSDaniel Finn static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
5338be53adSDaniel Finn {
5438be53adSDaniel Finn   PetscInt d;
5538be53adSDaniel Finn   *u = 1.0;
5638be53adSDaniel Finn   for (d = 0; d < dim; ++d) *u += (d+1)*PetscSqr(x[d]);
5738be53adSDaniel Finn   return 0;
5838be53adSDaniel Finn }
5938be53adSDaniel Finn 
6038be53adSDaniel Finn static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
6138be53adSDaniel Finn                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
6238be53adSDaniel Finn                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
6338be53adSDaniel Finn                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
6438be53adSDaniel Finn {
6538be53adSDaniel Finn   PetscInt d;
6638be53adSDaniel Finn   f0[0] += u[0];
6738be53adSDaniel 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]);
6838be53adSDaniel Finn }
6938be53adSDaniel Finn 
7038be53adSDaniel Finn static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
7138be53adSDaniel Finn                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
7238be53adSDaniel Finn                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
7338be53adSDaniel Finn                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
7438be53adSDaniel Finn {
7538be53adSDaniel Finn     PetscInt d;
7638be53adSDaniel Finn     switch (dim) {
7738be53adSDaniel Finn         case 1:
7838be53adSDaniel Finn             f0[0] = 1.0;
7938be53adSDaniel Finn             break;
8038be53adSDaniel Finn         case 2:
8138be53adSDaniel Finn             f0[0] = 5.0;
8238be53adSDaniel Finn             break;
8338be53adSDaniel Finn         case 3:
8438be53adSDaniel Finn             f0[0] = 11.0;
8538be53adSDaniel Finn             break;
8638be53adSDaniel Finn         default:
8738be53adSDaniel Finn             f0[0] = 5.0;
8838be53adSDaniel Finn             break;
8938be53adSDaniel Finn     }
9038be53adSDaniel Finn     f0[0] += u[0];
9138be53adSDaniel Finn     for (d = 0; d < dim; ++d) f0[0] -= (d+1)*PetscSqr(x[d]);
9238be53adSDaniel Finn }
9338be53adSDaniel Finn 
9438be53adSDaniel Finn static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
9538be53adSDaniel Finn                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9638be53adSDaniel Finn                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
9738be53adSDaniel Finn                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
9838be53adSDaniel Finn {
9938be53adSDaniel Finn   PetscInt d;
10038be53adSDaniel Finn   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
10138be53adSDaniel Finn }
10238be53adSDaniel Finn 
10338be53adSDaniel Finn static PetscErrorCode ProcessOptions(DM dm, AppCtx *options)
10438be53adSDaniel Finn {
10538be53adSDaniel Finn   MPI_Comm       comm;
10638be53adSDaniel Finn   PetscInt       dim;
10738be53adSDaniel Finn   PetscErrorCode ierr;
10838be53adSDaniel Finn 
10938be53adSDaniel Finn   PetscFunctionBeginUser;
11038be53adSDaniel Finn   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
11138be53adSDaniel Finn   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11238be53adSDaniel Finn   options->trig = PETSC_FALSE;
11338be53adSDaniel Finn 
11438be53adSDaniel Finn   ierr = PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");CHKERRQ(ierr);
11538be53adSDaniel 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*1e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
11738be53adSDaniel Finn 
11838be53adSDaniel Finn   PetscFunctionReturn(0);
11938be53adSDaniel Finn }
12038be53adSDaniel Finn 
12138be53adSDaniel Finn static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
12238be53adSDaniel Finn {
12338be53adSDaniel Finn   PetscErrorCode ierr;
12438be53adSDaniel Finn 
12538be53adSDaniel Finn   PetscFunctionBeginUser;
12638be53adSDaniel Finn   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
12738be53adSDaniel Finn   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
12838be53adSDaniel Finn   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
12938be53adSDaniel Finn 
13038be53adSDaniel Finn   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
13138be53adSDaniel Finn   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
13238be53adSDaniel Finn   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
13338be53adSDaniel Finn 
13438be53adSDaniel Finn   PetscFunctionReturn(0);
13538be53adSDaniel Finn }
13638be53adSDaniel Finn 
13738be53adSDaniel Finn static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
13838be53adSDaniel Finn {
13938be53adSDaniel Finn   PetscDS        ds;
14038be53adSDaniel Finn   DMLabel        label;
14138be53adSDaniel Finn   const PetscInt id = 1;
14238be53adSDaniel Finn   PetscErrorCode ierr;
14338be53adSDaniel Finn 
14438be53adSDaniel Finn   PetscFunctionBeginUser;
14538be53adSDaniel Finn   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
14638be53adSDaniel Finn   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
14738be53adSDaniel Finn   if (user->trig) {
14838be53adSDaniel Finn     ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
14938be53adSDaniel Finn     ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu);CHKERRQ(ierr);
15038be53adSDaniel Finn     ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr);
15138be53adSDaniel Finn     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr);
15238be53adSDaniel Finn     ierr = PetscPrintf(PETSC_COMM_WORLD,"Trig Exact Solution\n");CHKERRQ(ierr);
15338be53adSDaniel Finn   } else {
15438be53adSDaniel Finn     ierr = PetscDSSetResidual(ds, 0, f0_quad_u, f1_u);CHKERRQ(ierr);
15538be53adSDaniel Finn     ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu);CHKERRQ(ierr);
15638be53adSDaniel Finn     ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr);
15738be53adSDaniel Finn     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL);CHKERRQ(ierr);
15838be53adSDaniel Finn   }
15938be53adSDaniel Finn   PetscFunctionReturn(0);
16038be53adSDaniel Finn }
16138be53adSDaniel Finn 
16238be53adSDaniel Finn static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
16338be53adSDaniel Finn {
16438be53adSDaniel Finn   DM             cdm = dm;
16538be53adSDaniel Finn   PetscFE        fe;
16638be53adSDaniel Finn   DMPolytopeType ct;
16738be53adSDaniel Finn   PetscBool      simplex;
16838be53adSDaniel Finn   PetscInt       dim, cStart;
16938be53adSDaniel Finn   char           prefix[PETSC_MAX_PATH_LEN];
17038be53adSDaniel Finn   PetscErrorCode ierr;
17138be53adSDaniel Finn 
17238be53adSDaniel Finn   PetscFunctionBeginUser;
17338be53adSDaniel Finn   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
17438be53adSDaniel Finn 
17538be53adSDaniel Finn   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
17638be53adSDaniel Finn   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
17738be53adSDaniel Finn   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
17838be53adSDaniel Finn   /* Create finite element */
17938be53adSDaniel Finn   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
18038be53adSDaniel Finn   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
18138be53adSDaniel Finn   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
18238be53adSDaniel Finn   /* Set discretization and boundary conditions for each mesh */
18338be53adSDaniel Finn   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
18438be53adSDaniel Finn   ierr = DMCreateDS(dm);CHKERRQ(ierr);
18538be53adSDaniel Finn   ierr = (*setup)(dm, user);CHKERRQ(ierr);
18638be53adSDaniel Finn   while (cdm) {
18738be53adSDaniel Finn     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
18838be53adSDaniel Finn     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
18938be53adSDaniel Finn   }
19038be53adSDaniel Finn   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
19138be53adSDaniel Finn   PetscFunctionReturn(0);
19238be53adSDaniel Finn }
19338be53adSDaniel Finn 
19438be53adSDaniel Finn int main(int argc, char **argv)
19538be53adSDaniel Finn {
19638be53adSDaniel Finn     DM             dm;   /* Problem specification */
19738be53adSDaniel Finn     PetscDS        ds;
19838be53adSDaniel Finn     SNES           snes; /* Nonlinear solver */
19938be53adSDaniel Finn     Vec            u;    /* Solutions */
20038be53adSDaniel Finn     AppCtx         user; /* User-defined work context */
20138be53adSDaniel Finn     PetscErrorCode ierr;
20238be53adSDaniel Finn 
20338be53adSDaniel Finn     ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
20438be53adSDaniel Finn     /* Primal system */
20538be53adSDaniel Finn     ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
20638be53adSDaniel Finn     ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
20738be53adSDaniel Finn     ierr = ProcessOptions(dm, &user);CHKERRQ(ierr);
20838be53adSDaniel Finn     ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
20938be53adSDaniel Finn     ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
21038be53adSDaniel Finn     ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
21138be53adSDaniel Finn     ierr = VecSet(u, 0.0);CHKERRQ(ierr);
21238be53adSDaniel Finn     ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
21338be53adSDaniel Finn     ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
21438be53adSDaniel Finn     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
21538be53adSDaniel Finn     ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
21638be53adSDaniel Finn 
21738be53adSDaniel Finn     /*Looking for field error*/
21838be53adSDaniel Finn     PetscInt Nfields;
21938be53adSDaniel Finn     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
22038be53adSDaniel Finn     ierr = PetscDSGetNumFields(ds, &Nfields);CHKERRQ(ierr);
22138be53adSDaniel Finn     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
22238be53adSDaniel Finn     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
22338be53adSDaniel Finn     ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
22438be53adSDaniel Finn 
22538be53adSDaniel Finn     /* Cleanup */
22638be53adSDaniel Finn     ierr = VecDestroy(&u);CHKERRQ(ierr);
22738be53adSDaniel Finn     ierr = SNESDestroy(&snes);CHKERRQ(ierr);
22838be53adSDaniel Finn     ierr = DMDestroy(&dm);CHKERRQ(ierr);
22938be53adSDaniel Finn     ierr = PetscFinalize();
23038be53adSDaniel Finn     return ierr;
23138be53adSDaniel Finn }
23238be53adSDaniel Finn 
23338be53adSDaniel Finn /*TEST
23438be53adSDaniel Finn test:
23538be53adSDaniel Finn   # L_2 convergence rate: 1.9
23638be53adSDaniel Finn   suffix: 2d_p1_conv
23738be53adSDaniel Finn   requires: triangle
23838be53adSDaniel Finn   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
23938be53adSDaniel Finn test:
24038be53adSDaniel Finn   # L_2 convergence rate: 1.9
24138be53adSDaniel Finn   suffix: 2d_q1_conv
24238be53adSDaniel Finn   args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
24338be53adSDaniel Finn test:
24438be53adSDaniel Finn   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
24538be53adSDaniel Finn   suffix: 3d_p1_conv
24638be53adSDaniel Finn   requires: ctetgen
24738be53adSDaniel Finn   args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu
24838be53adSDaniel Finn test:
24938be53adSDaniel Finn   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
25038be53adSDaniel Finn   suffix: 3d_q1_conv
25138be53adSDaniel 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
25238be53adSDaniel Finn test:
25338be53adSDaniel Finn   # L_2 convergence rate: 1.9
25438be53adSDaniel Finn   suffix: 2d_p1_trig_conv
25538be53adSDaniel Finn   requires: triangle
25638be53adSDaniel Finn   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig
25738be53adSDaniel Finn test:
25838be53adSDaniel Finn   # L_2 convergence rate: 1.9
25938be53adSDaniel Finn   suffix: 2d_q1_trig_conv
26038be53adSDaniel 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
26138be53adSDaniel Finn test:
26238be53adSDaniel Finn   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
26338be53adSDaniel Finn   suffix: 3d_p1_trig_conv
26438be53adSDaniel Finn   requires: ctetgen
26538be53adSDaniel 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
26638be53adSDaniel Finn test:
26738be53adSDaniel Finn   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
26838be53adSDaniel Finn   suffix: 3d_q1_trig_conv
26938be53adSDaniel 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
27038be53adSDaniel Finn test:
27138be53adSDaniel Finn   suffix: 2d_p1_gmg_vcycle
27238be53adSDaniel Finn   requires: triangle
27338be53adSDaniel Finn   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
27438be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
27538be53adSDaniel Finn     -mg_levels_ksp_max_it 1 \
27638be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
27738be53adSDaniel Finn test:
27838be53adSDaniel Finn   suffix: 2d_p1_gmg_fcycle
27938be53adSDaniel Finn   requires: triangle
28038be53adSDaniel Finn   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
28138be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \
28238be53adSDaniel Finn     -mg_levels_ksp_max_it 2 \
28338be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
28438be53adSDaniel Finn test:
28538be53adSDaniel Finn   suffix: 2d_p1_gmg_vcycle_trig
28638be53adSDaniel Finn   requires: triangle
28738be53adSDaniel Finn   args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
28838be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
28938be53adSDaniel Finn     -mg_levels_ksp_max_it 1 \
29038be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
29138be53adSDaniel Finn TEST*/
292