xref: /petsc/src/snes/tests/ex13.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
15e1f5104SMark static char help[] = "Benchmark Poisson Problem in 2d and 3d with finite elements.\n\
2f9244615SMatthew G. Knepley We solve the Poisson problem in a rectangular domain\n\
3f9244615SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
45e1f5104SMark 
55e1f5104SMark #include <petscdmplex.h>
65e1f5104SMark #include <petscsnes.h>
75e1f5104SMark #include <petscds.h>
85e1f5104SMark #include <petscconvest.h>
95e1f5104SMark 
105e1f5104SMark typedef struct {
11f9244615SMatthew G. Knepley   PetscInt  nit;    /* Number of benchmark iterations */
12f9244615SMatthew G. Knepley   PetscBool strong; /* Do not integrate the Laplacian by parts */
135e1f5104SMark } AppCtx;
145e1f5104SMark 
155e1f5104SMark static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
165e1f5104SMark {
175e1f5104SMark   PetscInt d;
185e1f5104SMark   *u = 0.0;
195e1f5104SMark   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
205e1f5104SMark   return 0;
215e1f5104SMark }
225e1f5104SMark 
235e1f5104SMark static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
245e1f5104SMark                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
255e1f5104SMark                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
265e1f5104SMark                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
275e1f5104SMark {
285e1f5104SMark   PetscInt d;
295e1f5104SMark   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
305e1f5104SMark }
315e1f5104SMark 
325e1f5104SMark static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
335e1f5104SMark                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
345e1f5104SMark                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
355e1f5104SMark                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
365e1f5104SMark {
375e1f5104SMark   PetscInt d;
385e1f5104SMark   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
395e1f5104SMark }
405e1f5104SMark 
415e1f5104SMark static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
425e1f5104SMark                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
435e1f5104SMark                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
445e1f5104SMark                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
455e1f5104SMark {
465e1f5104SMark   PetscInt d;
475e1f5104SMark   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
485e1f5104SMark }
495e1f5104SMark 
50f9244615SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51f9244615SMatthew G. Knepley {
52f9244615SMatthew G. Knepley   *u = PetscSqr(x[0]) + PetscSqr(x[1]);
53f9244615SMatthew G. Knepley   return 0;
54f9244615SMatthew G. Knepley }
55f9244615SMatthew G. Knepley 
56f9244615SMatthew G. Knepley static void f0_strong_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
57f9244615SMatthew G. Knepley                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
58f9244615SMatthew G. Knepley                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
59f9244615SMatthew G. Knepley                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
60f9244615SMatthew G. Knepley {
61f9244615SMatthew G. Knepley   PetscInt d;
62f9244615SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] -= u_x[dim + d*dim+d];
63f9244615SMatthew G. Knepley   f0[0] += 4.0;
64f9244615SMatthew G. Knepley }
65f9244615SMatthew G. Knepley 
665e1f5104SMark static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
675e1f5104SMark {
685e1f5104SMark   PetscErrorCode ierr;
695e1f5104SMark 
705e1f5104SMark   PetscFunctionBeginUser;
710c569c6eSMark   options->nit    = 10;
72f9244615SMatthew G. Knepley   options->strong = PETSC_FALSE;
735e1f5104SMark   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL));
761e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
775e1f5104SMark   PetscFunctionReturn(0);
785e1f5104SMark }
795e1f5104SMark 
805e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
815e1f5104SMark {
825e1f5104SMark   PetscFunctionBeginUser;
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(*dm, user));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
885e1f5104SMark   PetscFunctionReturn(0);
895e1f5104SMark }
905e1f5104SMark 
915e1f5104SMark static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
925e1f5104SMark {
93f9244615SMatthew G. Knepley   PetscDS        ds;
9445480ffeSMatthew G. Knepley   DMLabel        label;
955e1f5104SMark   const PetscInt id = 1;
965e1f5104SMark 
975e1f5104SMark   PetscFunctionBeginUser;
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
100f9244615SMatthew G. Knepley   if (user->strong) {
101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, f0_strong_u, NULL));
102*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetExactSolution(ds, 0, quadratic_u, user));
103*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, user, NULL));
104f9244615SMatthew G. Knepley   } else {
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetExactSolution(ds, 0, trig_u, user));
108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL));
109f9244615SMatthew G. Knepley   }
1105e1f5104SMark   PetscFunctionReturn(0);
1115e1f5104SMark }
1125e1f5104SMark 
1135e1f5104SMark static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
1145e1f5104SMark {
1155e1f5104SMark   DM             cdm = dm;
1165e1f5104SMark   PetscFE        fe;
1175e1f5104SMark   DMPolytopeType ct;
1185e1f5104SMark   PetscBool      simplex;
1195e1f5104SMark   PetscInt       dim, cStart;
1205e1f5104SMark   char           prefix[PETSC_MAX_PATH_LEN];
1215e1f5104SMark 
1225e1f5104SMark   PetscFunctionBeginUser;
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, cStart, &ct));
1265e1f5104SMark   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1275e1f5104SMark   /* Create finite element */
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, name));
1315e1f5104SMark   /* Set discretization and boundary conditions for each mesh */
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*setup)(dm, user));
1355e1f5104SMark   while (cdm) {
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm,cdm));
1375e1f5104SMark     /* TODO: Check whether the boundary of coarse meshes is marked */
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
1395e1f5104SMark   }
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
1415e1f5104SMark   PetscFunctionReturn(0);
1425e1f5104SMark }
1435e1f5104SMark 
1445e1f5104SMark int main(int argc, char **argv)
1455e1f5104SMark {
1465e1f5104SMark   DM             dm;   /* Problem specification */
1475e1f5104SMark   SNES           snes; /* Nonlinear solver */
1485e1f5104SMark   Vec            u;    /* Solutions */
1495e1f5104SMark   AppCtx         user; /* User-defined work context */
1505e1f5104SMark   PetscErrorCode ierr;
1515e1f5104SMark 
1525e1f5104SMark   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
1545e1f5104SMark   /* Primal system */
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes));
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes, dm));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u, 0.0));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "potential"));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESCheckFromOptions(snes, u));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, NULL, u));
1665e1f5104SMark   /* Benchmark system */
1670c569c6eSMark   if (user.nit) {
1685e1f5104SMark #if defined(PETSC_USE_LOG)
16918fb0606SStefano Zampini     PetscLogStage kspstage,pcstage;
1705e1f5104SMark #endif
1715e1f5104SMark     KSP       ksp;
17218fb0606SStefano Zampini     PC        pc;
1735e1f5104SMark     Vec       b;
1740c569c6eSMark     PetscInt  i;
1755cdb9f41SMark Adams     PetscLogDouble time;
176*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsClearValue(NULL,"-ksp_monitor"));
177*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsClearValue(NULL,"-ksp_view"));
178*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetKSP(snes, &ksp));
179*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetSolution(snes, &u));
180*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetFromOptions(ksp));
181*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(u, 0.0));
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetFunction(snes, &b, NULL, NULL));
183*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp, &pc));
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStageRegister("PCSetUp", &pcstage));
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(pcstage));
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetUp(pc));
187*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStageRegister("KSP Solve only", &kspstage));
189*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscTime(&time));
190*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(kspstage));
1910c569c6eSMark     for (i=0;i<user.nit;i++) {
192*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecZeroEntries(u));
193*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSolve(ksp, b, u));
1940c569c6eSMark     }
195*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
196*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscTimeSubtract(&time));
1975cdb9f41SMark Adams     // ierr = PetscPrintf(PETSC_COMM_WORLD,"Solve time: %g\n",-time); // breaks CI
1985e1f5104SMark   }
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes, &u));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-potential_view"));
2015e1f5104SMark   /* Cleanup */
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
2055e1f5104SMark   ierr = PetscFinalize();
2065e1f5104SMark   return ierr;
2075e1f5104SMark }
2085e1f5104SMark 
2095e1f5104SMark /*TEST
2105e1f5104SMark 
2115e1f5104SMark   test:
212f9244615SMatthew G. Knepley     suffix: strong
213f9244615SMatthew G. Knepley     requires: triangle
21430602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \
215f9244615SMatthew G. Knepley           -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong
216f9244615SMatthew G. Knepley 
217f9244615SMatthew G. Knepley   test:
2185e1f5104SMark     suffix: bench
2190c569c6eSMark     nsize: 4
220e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,8 -dm_refine 1 \
2210c569c6eSMark           -petscpartitioner_type simple -petscpartitioner_simple_process_grid 1,1,2 -petscpartitioner_simple_node_grid 1,1,2 \
22273f7197eSJed Brown           -potential_petscspace_degree 2 -ksp_type cg -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -benchmark_it 1 -dm_view -snes_rtol 1.e-4
2235e1f5104SMark 
2240c569c6eSMark   test:
22518fb0606SStefano Zampini     suffix: comparison
22618fb0606SStefano Zampini     nsize: 4
22730602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
228e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
22973f7197eSJed Brown       -dm_plex_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \
23018fb0606SStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4
23118fb0606SStefano Zampini 
23218fb0606SStefano Zampini   test:
2330c569c6eSMark     suffix: cuda
2340c569c6eSMark     nsize: 4
2350c569c6eSMark     requires: cuda
23618fb0606SStefano Zampini     output_file: output/ex13_comparison.out
23730602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
238e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
23973f7197eSJed Brown       -dm_plex_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \
240e589036eSStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijcusparse -dm_vec_type cuda
2410c569c6eSMark 
2420c569c6eSMark   test:
24318fb0606SStefano Zampini     suffix: kokkos_comp
24418fb0606SStefano Zampini     nsize: 4
2453078479eSJunchao Zhang     requires: !sycl kokkos_kernels
24618fb0606SStefano Zampini     output_file: output/ex13_comparison.out
24730602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
248e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
24973f7197eSJed Brown       -dm_plex_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \
25018fb0606SStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijkokkos -dm_vec_type kokkos
25118fb0606SStefano Zampini 
25218fb0606SStefano Zampini   test:
2530c569c6eSMark     nsize: 4
2543078479eSJunchao Zhang     requires: !sycl kokkos_kernels
2550c569c6eSMark     suffix: kokkos
256e600fa54SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 2,8 -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \
25773f7197eSJed Brown           -petscpartitioner_simple_node_grid 2,1 -dm_plex_simplex 0 -potential_petscspace_degree 1 -dm_refine 1 -ksp_type cg -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -ksp_norm_type unpreconditioned \
25873f7197eSJed Brown           -pc_gamg_esteig_ksp_type cg -ksp_converged_reason -snes_monitor_short -snes_rtol 1.e-4 -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
259aa5a873eSStefano Zampini 
260aa5a873eSStefano Zampini   test:
261aa5a873eSStefano Zampini     suffix: aijmkl_comp
262aa5a873eSStefano Zampini     nsize: 4
263c4ad6305SSatish Balay     requires: mkl_sparse
264aa5a873eSStefano Zampini     output_file: output/ex13_comparison.out
26530602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
266e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
26773f7197eSJed Brown       -dm_plex_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \
268aa5a873eSStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl
269aa5a873eSStefano Zampini 
270aa5a873eSStefano Zampini   test:
271aa5a873eSStefano Zampini     suffix: aijmkl_seq
272aa5a873eSStefano Zampini     nsize: 1
273c4ad6305SSatish Balay     requires: mkl_sparse
274aa5a873eSStefano Zampini     TODO: broken (INDEFINITE PC)
27573f7197eSJed Brown     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_refine 1 -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_simplex 0 -snes_monitor_short \
27673f7197eSJed Brown           -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_sym_graph 0 -pc_gamg_threshold -1 -pc_gamg_square_graph 10 -pc_gamg_process_eq_limit 400 \
27773f7197eSJed Brown           -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -pc_gamg_esteig_ksp_type cg -ksp_type cg -ksp_norm_type unpreconditioned -snes_converged_reason \
27873f7197eSJed Brown           -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard
279aa5a873eSStefano Zampini 
2805e1f5104SMark TEST*/
281