xref: /petsc/src/snes/tutorials/ex26.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 */
26d71ae5a4SJacob Faibussowitsch static void g0_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 g0[])
27d71ae5a4SJacob Faibussowitsch {
2838be53adSDaniel Finn   PetscInt d;
2938be53adSDaniel Finn   for (d = 0; d < dim; ++d) g0[0] = 1.0;
3038be53adSDaniel Finn }
3138be53adSDaniel Finn 
32d71ae5a4SJacob Faibussowitsch 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[])
33d71ae5a4SJacob Faibussowitsch {
3438be53adSDaniel Finn   PetscInt d;
3538be53adSDaniel Finn   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
3638be53adSDaniel Finn }
3738be53adSDaniel Finn 
38d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
39d71ae5a4SJacob Faibussowitsch {
4038be53adSDaniel Finn   PetscInt d;
4138be53adSDaniel Finn   *u = 0.0;
4238be53adSDaniel Finn   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
43*3ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
4438be53adSDaniel Finn }
4538be53adSDaniel Finn 
46d71ae5a4SJacob Faibussowitsch static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
47d71ae5a4SJacob Faibussowitsch {
4838be53adSDaniel Finn   PetscInt d;
4938be53adSDaniel Finn   *u = 1.0;
5038be53adSDaniel Finn   for (d = 0; d < dim; ++d) *u += (d + 1) * PetscSqr(x[d]);
51*3ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
5238be53adSDaniel Finn }
5338be53adSDaniel Finn 
54d71ae5a4SJacob Faibussowitsch static void f0_trig_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 f0[])
55d71ae5a4SJacob Faibussowitsch {
5638be53adSDaniel Finn   PetscInt d;
5738be53adSDaniel Finn   f0[0] += u[0];
5838be53adSDaniel 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]);
5938be53adSDaniel Finn }
6038be53adSDaniel Finn 
61d71ae5a4SJacob Faibussowitsch static void f0_quad_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 f0[])
62d71ae5a4SJacob Faibussowitsch {
6338be53adSDaniel Finn   PetscInt d;
6438be53adSDaniel Finn   switch (dim) {
65d71ae5a4SJacob Faibussowitsch   case 1:
66d71ae5a4SJacob Faibussowitsch     f0[0] = 1.0;
67d71ae5a4SJacob Faibussowitsch     break;
68d71ae5a4SJacob Faibussowitsch   case 2:
69d71ae5a4SJacob Faibussowitsch     f0[0] = 5.0;
70d71ae5a4SJacob Faibussowitsch     break;
71d71ae5a4SJacob Faibussowitsch   case 3:
72d71ae5a4SJacob Faibussowitsch     f0[0] = 11.0;
73d71ae5a4SJacob Faibussowitsch     break;
74d71ae5a4SJacob Faibussowitsch   default:
75d71ae5a4SJacob Faibussowitsch     f0[0] = 5.0;
76d71ae5a4SJacob Faibussowitsch     break;
7738be53adSDaniel Finn   }
7838be53adSDaniel Finn   f0[0] += u[0];
7938be53adSDaniel Finn   for (d = 0; d < dim; ++d) f0[0] -= (d + 1) * PetscSqr(x[d]);
8038be53adSDaniel Finn }
8138be53adSDaniel Finn 
82d71ae5a4SJacob Faibussowitsch 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[])
83d71ae5a4SJacob Faibussowitsch {
8438be53adSDaniel Finn   PetscInt d;
8538be53adSDaniel Finn   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
8638be53adSDaniel Finn }
8738be53adSDaniel Finn 
88d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(DM dm, AppCtx *options)
89d71ae5a4SJacob Faibussowitsch {
9038be53adSDaniel Finn   MPI_Comm comm;
9138be53adSDaniel Finn   PetscInt dim;
9238be53adSDaniel Finn 
9338be53adSDaniel Finn   PetscFunctionBeginUser;
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
959566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
9638be53adSDaniel Finn   options->trig = PETSC_FALSE;
9738be53adSDaniel Finn 
98d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Helmholtz Problem Options", "DMPLEX");
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-exact_trig", "Use trigonometric exact solution (better for more complex finite elements)", "ex26.c", options->trig, &options->trig, NULL));
100d0609cedSBarry Smith   PetscOptionsEnd();
10138be53adSDaniel Finn 
102*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10338be53adSDaniel Finn }
10438be53adSDaniel Finn 
105d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
106d71ae5a4SJacob Faibussowitsch {
10738be53adSDaniel Finn   PetscFunctionBeginUser;
1089566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1099566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
11138be53adSDaniel Finn 
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh"));
1139566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
1149566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
11538be53adSDaniel Finn 
116*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11738be53adSDaniel Finn }
11838be53adSDaniel Finn 
119d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
120d71ae5a4SJacob Faibussowitsch {
12138be53adSDaniel Finn   PetscDS        ds;
12238be53adSDaniel Finn   DMLabel        label;
12338be53adSDaniel Finn   const PetscInt id = 1;
12438be53adSDaniel Finn 
12538be53adSDaniel Finn   PetscFunctionBeginUser;
1269566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1279566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
12838be53adSDaniel Finn   if (user->trig) {
1299566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
1309566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu));
1319566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
1329566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
1339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Trig Exact Solution\n"));
13438be53adSDaniel Finn   } else {
1359566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_quad_u, f1_u));
1369566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, g3_uu));
1379566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quad_u, user));
1389566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quad_u, NULL, user, NULL));
13938be53adSDaniel Finn   }
140*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14138be53adSDaniel Finn }
14238be53adSDaniel Finn 
143d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
144d71ae5a4SJacob Faibussowitsch {
14538be53adSDaniel Finn   DM             cdm = dm;
14638be53adSDaniel Finn   PetscFE        fe;
14738be53adSDaniel Finn   DMPolytopeType ct;
14838be53adSDaniel Finn   PetscBool      simplex;
14938be53adSDaniel Finn   PetscInt       dim, cStart;
15038be53adSDaniel Finn   char           prefix[PETSC_MAX_PATH_LEN];
15138be53adSDaniel Finn 
15238be53adSDaniel Finn   PetscFunctionBeginUser;
1539566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
15438be53adSDaniel Finn 
1559566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1569566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
15738be53adSDaniel Finn   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
15838be53adSDaniel Finn   /* Create finite element */
1599566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1609566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
16238be53adSDaniel Finn   /* Set discretization and boundary conditions for each mesh */
1639566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1649566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1659566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
16638be53adSDaniel Finn   while (cdm) {
1679566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
1689566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
16938be53adSDaniel Finn   }
1709566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
171*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17238be53adSDaniel Finn }
17338be53adSDaniel Finn 
174d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
175d71ae5a4SJacob Faibussowitsch {
17638be53adSDaniel Finn   DM      dm; /* Problem specification */
17738be53adSDaniel Finn   PetscDS ds;
17838be53adSDaniel Finn   SNES    snes; /* Nonlinear solver */
17938be53adSDaniel Finn   Vec     u;    /* Solutions */
18038be53adSDaniel Finn   AppCtx  user; /* User-defined work context */
18138be53adSDaniel Finn 
182327415f7SBarry Smith   PetscFunctionBeginUser;
1839566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
18438be53adSDaniel Finn   /* Primal system */
1859566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1869566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1879566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(dm, &user));
1889566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
1899566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1909566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
1919566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
1929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
1939566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
1949566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
1959566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
19638be53adSDaniel Finn 
19738be53adSDaniel Finn   /* Looking for field error */
19838be53adSDaniel Finn   PetscInt Nfields;
1999566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nfields));
2019566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
2029566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
2039566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
20438be53adSDaniel Finn 
20538be53adSDaniel Finn   /* Cleanup */
2069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2079566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2089566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2099566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
210b122ec5aSJacob Faibussowitsch   return 0;
21138be53adSDaniel Finn }
21238be53adSDaniel Finn 
21338be53adSDaniel Finn /*TEST
21438be53adSDaniel Finn test:
21538be53adSDaniel Finn   # L_2 convergence rate: 1.9
21638be53adSDaniel Finn   suffix: 2d_p1_conv
21738be53adSDaniel Finn   requires: triangle
21838be53adSDaniel Finn   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
21938be53adSDaniel Finn test:
22038be53adSDaniel Finn   # L_2 convergence rate: 1.9
22138be53adSDaniel Finn   suffix: 2d_q1_conv
22238be53adSDaniel Finn   args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu
22338be53adSDaniel Finn test:
22438be53adSDaniel Finn   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
22538be53adSDaniel Finn   suffix: 3d_p1_conv
22638be53adSDaniel Finn   requires: ctetgen
22738be53adSDaniel Finn   args: -dm_plex_dim 3 -dm_refine 2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu
22838be53adSDaniel Finn test:
22938be53adSDaniel Finn   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
23038be53adSDaniel Finn   suffix: 3d_q1_conv
23138be53adSDaniel 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
23238be53adSDaniel Finn test:
23338be53adSDaniel Finn   # L_2 convergence rate: 1.9
23438be53adSDaniel Finn   suffix: 2d_p1_trig_conv
23538be53adSDaniel Finn   requires: triangle
23638be53adSDaniel Finn   args: -potential_petscspace_degree 1 -snes_convergence_estimate -dm_refine 2 -convest_num_refine 3 -pc_type lu -exact_trig
23738be53adSDaniel Finn test:
23838be53adSDaniel Finn   # L_2 convergence rate: 1.9
23938be53adSDaniel Finn   suffix: 2d_q1_trig_conv
24038be53adSDaniel 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
24138be53adSDaniel Finn test:
24238be53adSDaniel Finn   # Using -convest_num_refine 3 we get L_2 convergence rate: -1.5
24338be53adSDaniel Finn   suffix: 3d_p1_trig_conv
24438be53adSDaniel Finn   requires: ctetgen
24538be53adSDaniel 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
24638be53adSDaniel Finn test:
24738be53adSDaniel Finn   # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: -1.2
24838be53adSDaniel Finn   suffix: 3d_q1_trig_conv
24938be53adSDaniel 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
25038be53adSDaniel Finn test:
25138be53adSDaniel Finn   suffix: 2d_p1_gmg_vcycle
25238be53adSDaniel Finn   requires: triangle
25338be53adSDaniel Finn   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
25438be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
25538be53adSDaniel Finn     -mg_levels_ksp_max_it 1 \
25638be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
25738be53adSDaniel Finn test:
25838be53adSDaniel Finn   suffix: 2d_p1_gmg_fcycle
25938be53adSDaniel Finn   requires: triangle
26038be53adSDaniel Finn   args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
26138be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg -pc_mg_type full \
26238be53adSDaniel Finn     -mg_levels_ksp_max_it 2 \
26338be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
26438be53adSDaniel Finn test:
26538be53adSDaniel Finn   suffix: 2d_p1_gmg_vcycle_trig
26638be53adSDaniel Finn   requires: triangle
26738be53adSDaniel Finn   args: -exact_trig -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
26838be53adSDaniel Finn     -ksp_type cg -ksp_rtol 1e-10 -pc_type mg \
26938be53adSDaniel Finn     -mg_levels_ksp_max_it 1 \
27038be53adSDaniel Finn     -mg_levels_pc_type jacobi -snes_monitor -ksp_monitor
2715f34f2dcSJed Brown test:
2725f34f2dcSJed Brown   suffix: 2d_q3_cgns
2735f34f2dcSJed Brown   requires: cgns
2745f34f2dcSJed Brown   args: -dm_plex_simplex 0 -dm_plex_dim 2 -dm_plex_box_faces 3,3 -snes_view_solution cgns:sol.cgns -potential_petscspace_degree 3 -dm_coord_petscspace_degree 3
27538be53adSDaniel Finn TEST*/
276