xref: /petsc/src/snes/tests/ex13.c (revision 9566063d113dddea24716c546802770db7481bc0)
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;
73*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");PetscCall(ierr);
74*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL));
75*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL));
76*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
775e1f5104SMark   PetscFunctionReturn(0);
785e1f5104SMark }
795e1f5104SMark 
805e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
815e1f5104SMark {
825e1f5104SMark   PetscFunctionBeginUser;
83*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
84*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
85*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
86*9566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
87*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
99*9566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
100f9244615SMatthew G. Knepley   if (user->strong) {
101*9566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_strong_u, NULL));
102*9566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user));
103*9566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, user, NULL));
104f9244615SMatthew G. Knepley   } else {
105*9566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
106*9566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
107*9566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
108*9566063dSJacob Faibussowitsch     PetscCall(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*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
124*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
125*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1265e1f5104SMark   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1275e1f5104SMark   /* Create finite element */
128*9566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
129*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
130*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, name));
1315e1f5104SMark   /* Set discretization and boundary conditions for each mesh */
132*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
133*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
134*9566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
1355e1f5104SMark   while (cdm) {
136*9566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm,cdm));
1375e1f5104SMark     /* TODO: Check whether the boundary of coarse meshes is marked */
138*9566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
1395e1f5104SMark   }
140*9566063dSJacob Faibussowitsch   PetscCall(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 
151*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
152*9566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1535e1f5104SMark   /* Primal system */
154*9566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
155*9566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
156*9566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
157*9566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
158*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
159*9566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
160*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "potential"));
161*9566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
162*9566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
163*9566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
164*9566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1655e1f5104SMark   /* Benchmark system */
1660c569c6eSMark   if (user.nit) {
1675e1f5104SMark #if defined(PETSC_USE_LOG)
16818fb0606SStefano Zampini     PetscLogStage kspstage,pcstage;
1695e1f5104SMark #endif
1705e1f5104SMark     KSP       ksp;
17118fb0606SStefano Zampini     PC        pc;
1725e1f5104SMark     Vec       b;
1730c569c6eSMark     PetscInt  i;
1745cdb9f41SMark Adams     PetscLogDouble time;
175*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsClearValue(NULL,"-ksp_monitor"));
176*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsClearValue(NULL,"-ksp_view"));
177*9566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
178*9566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
179*9566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
180*9566063dSJacob Faibussowitsch     PetscCall(VecSet(u, 0.0));
181*9566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &b, NULL, NULL));
182*9566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
183*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("PCSetUp", &pcstage));
184*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(pcstage));
185*9566063dSJacob Faibussowitsch     PetscCall(PCSetUp(pc));
186*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
187*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("KSP Solve only", &kspstage));
188*9566063dSJacob Faibussowitsch     PetscCall(PetscTime(&time));
189*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(kspstage));
1900c569c6eSMark     for (i=0;i<user.nit;i++) {
191*9566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(u));
192*9566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, b, u));
1930c569c6eSMark     }
194*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
195*9566063dSJacob Faibussowitsch     PetscCall(PetscTimeSubtract(&time));
1965cdb9f41SMark Adams     // ierr = PetscPrintf(PETSC_COMM_WORLD,"Solve time: %g\n",-time); // breaks CI
1975e1f5104SMark   }
198*9566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
199*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
2005e1f5104SMark   /* Cleanup */
201*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
202*9566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
203*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
204*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
205b122ec5aSJacob Faibussowitsch   return 0;
2065e1f5104SMark }
2075e1f5104SMark 
2085e1f5104SMark /*TEST
2095e1f5104SMark 
2105e1f5104SMark   test:
211f9244615SMatthew G. Knepley     suffix: strong
212f9244615SMatthew G. Knepley     requires: triangle
21330602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \
214f9244615SMatthew G. Knepley           -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong
215f9244615SMatthew G. Knepley 
216f9244615SMatthew G. Knepley   test:
2175e1f5104SMark     suffix: bench
2180c569c6eSMark     nsize: 4
219e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,8 -dm_refine 1 \
2200c569c6eSMark           -petscpartitioner_type simple -petscpartitioner_simple_process_grid 1,1,2 -petscpartitioner_simple_node_grid 1,1,2 \
22173f7197eSJed 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
2225e1f5104SMark 
2230c569c6eSMark   test:
22418fb0606SStefano Zampini     suffix: comparison
22518fb0606SStefano Zampini     nsize: 4
22630602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
227e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
22873f7197eSJed 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 \
22918fb0606SStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4
23018fb0606SStefano Zampini 
23118fb0606SStefano Zampini   test:
2320c569c6eSMark     suffix: cuda
2330c569c6eSMark     nsize: 4
2340c569c6eSMark     requires: cuda
23518fb0606SStefano Zampini     output_file: output/ex13_comparison.out
23630602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
237e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
23873f7197eSJed 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 \
239e589036eSStefano 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
2400c569c6eSMark 
2410c569c6eSMark   test:
24218fb0606SStefano Zampini     suffix: kokkos_comp
24318fb0606SStefano Zampini     nsize: 4
2443078479eSJunchao Zhang     requires: !sycl kokkos_kernels
24518fb0606SStefano Zampini     output_file: output/ex13_comparison.out
24630602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
247e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
24873f7197eSJed 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 \
24918fb0606SStefano 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
25018fb0606SStefano Zampini 
25118fb0606SStefano Zampini   test:
2520c569c6eSMark     nsize: 4
2533078479eSJunchao Zhang     requires: !sycl kokkos_kernels
2540c569c6eSMark     suffix: kokkos
255e600fa54SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 2,8 -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \
25673f7197eSJed 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 \
25773f7197eSJed 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
258aa5a873eSStefano Zampini 
259aa5a873eSStefano Zampini   test:
260aa5a873eSStefano Zampini     suffix: aijmkl_comp
261aa5a873eSStefano Zampini     nsize: 4
262c4ad6305SSatish Balay     requires: mkl_sparse
263aa5a873eSStefano Zampini     output_file: output/ex13_comparison.out
26430602db0SMatthew G. Knepley     args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
265e600fa54SMatthew G. Knepley       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \
26673f7197eSJed 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 \
267aa5a873eSStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl
268aa5a873eSStefano Zampini 
269aa5a873eSStefano Zampini   test:
270aa5a873eSStefano Zampini     suffix: aijmkl_seq
271aa5a873eSStefano Zampini     nsize: 1
272c4ad6305SSatish Balay     requires: mkl_sparse
273aa5a873eSStefano Zampini     TODO: broken (INDEFINITE PC)
27473f7197eSJed 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 \
27573f7197eSJed 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 \
27673f7197eSJed 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 \
27773f7197eSJed Brown           -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard
278aa5a873eSStefano Zampini 
2795e1f5104SMark TEST*/
280