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; 739566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");PetscCall(ierr); 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL)); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL)); 769566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 775e1f5104SMark PetscFunctionReturn(0); 785e1f5104SMark } 795e1f5104SMark 805e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 815e1f5104SMark { 825e1f5104SMark PetscFunctionBeginUser; 839566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 849566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 859566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 869566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 879566063dSJacob 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; 989566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 999566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 100f9244615SMatthew G. Knepley if (user->strong) { 1019566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_strong_u, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user)); 1039566063dSJacob 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 { 1059566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 1069566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1079566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 1089566063dSJacob 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; 1239566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1249566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1259566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 126*2e776fa0SMark Adams simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; // false 1275e1f5104SMark /* Create finite element */ 1289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 1315e1f5104SMark /* Set discretization and boundary conditions for each mesh */ 1329566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 1339566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1349566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 1355e1f5104SMark while (cdm) { 1369566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 1375e1f5104SMark /* TODO: Check whether the boundary of coarse meshes is marked */ 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 1395e1f5104SMark } 1409566063dSJacob 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 */ 150*2e776fa0SMark Adams PetscLogDouble time; 151*2e776fa0SMark Adams Mat Amat; 1525e1f5104SMark 1539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 1549566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 155*2e776fa0SMark Adams /* system */ 1569566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 1579566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1589566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 1599566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 1609566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 161*2e776fa0SMark Adams PetscCall(SNESSetFromOptions(snes)); 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 1639566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 1649566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 165*2e776fa0SMark Adams PetscCall(PetscTime(&time)); 166*2e776fa0SMark Adams PetscCall(SNESSetUp(snes)); 167*2e776fa0SMark Adams PetscCall(SNESGetJacobian(snes, &Amat, NULL, NULL, NULL)); 168*2e776fa0SMark Adams PetscCall(MatSetOption(Amat,MAT_SPD,PETSC_TRUE)); 1699566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 170*2e776fa0SMark Adams PetscCall(PetscTimeSubtract(&time)); 171*2e776fa0SMark Adams // PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First Solve time: %g\n",-time)); 1725e1f5104SMark /* Benchmark system */ 1730c569c6eSMark if (user.nit) { 1745e1f5104SMark Vec b; 1750c569c6eSMark PetscInt i; 176*2e776fa0SMark Adams #if defined(PETSC_USE_LOG) 177*2e776fa0SMark Adams PetscLogStage kspstage; 178*2e776fa0SMark Adams #endif 179*2e776fa0SMark Adams PetscCall(PetscLogStageRegister("Solve only", &kspstage)); 1809566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(kspstage)); 181*2e776fa0SMark Adams PetscCall(SNESGetSolution(snes, &u)); 182*2e776fa0SMark Adams PetscCall(SNESGetFunction(snes, &b, NULL, NULL)); 1830c569c6eSMark for (i=0;i<user.nit;i++) { 1849566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(u)); 185*2e776fa0SMark Adams PetscCall(SNESSolve(snes, NULL, u)); 1860c569c6eSMark } 1879566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1885e1f5104SMark } 1899566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 1909566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 1915e1f5104SMark /* Cleanup */ 1929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1939566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1949566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 196b122ec5aSJacob Faibussowitsch return 0; 1975e1f5104SMark } 1985e1f5104SMark 1995e1f5104SMark /*TEST 2005e1f5104SMark 2015e1f5104SMark test: 202f9244615SMatthew G. Knepley suffix: strong 203f9244615SMatthew G. Knepley requires: triangle 20430602db0SMatthew G. Knepley args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \ 205f9244615SMatthew G. Knepley -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong 206f9244615SMatthew G. Knepley 207f9244615SMatthew G. Knepley test: 2085e1f5104SMark suffix: bench 2090c569c6eSMark nsize: 4 210*2e776fa0SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1 -dm_refine 2 -dm_view -ksp_monitor \ 211*2e776fa0SMark Adams -benchmark_it 1 \ 212*2e776fa0SMark Adams -dm_plex_box_upper 2,2,1 \ 213*2e776fa0SMark Adams -dm_plex_box_lower 0,0,0 \ 214*2e776fa0SMark Adams -dm_plex_dim 3 \ 215*2e776fa0SMark Adams -ksp_converged_reason \ 216*2e776fa0SMark Adams -ksp_norm_type unpreconditioned \ 217*2e776fa0SMark Adams -ksp_rtol 1.e-6 \ 218*2e776fa0SMark Adams -ksp_type cg \ 219*2e776fa0SMark Adams -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 \ 220*2e776fa0SMark Adams -mg_levels_ksp_max_it 1 \ 221*2e776fa0SMark Adams -mg_levels_ksp_type chebyshev \ 222*2e776fa0SMark Adams -mg_levels_pc_type jacobi \ 223*2e776fa0SMark Adams -pc_gamg_coarse_eq_limit 200 \ 224*2e776fa0SMark Adams -pc_gamg_coarse_grid_layout_type compact \ 225*2e776fa0SMark Adams -pc_gamg_esteig_ksp_max_it 5 \ 226*2e776fa0SMark Adams -pc_gamg_process_eq_limit 200 \ 227*2e776fa0SMark Adams -pc_gamg_repartition false \ 228*2e776fa0SMark Adams -pc_gamg_reuse_interpolation true \ 229*2e776fa0SMark Adams -pc_gamg_square_graph 0 \ 230*2e776fa0SMark Adams -pc_gamg_threshold 0.001 -pc_gamg_threshold_scale .5\ 231*2e776fa0SMark Adams -pc_gamg_type agg \ 232*2e776fa0SMark Adams -pc_type gamg \ 233*2e776fa0SMark Adams -petscpartitioner_simple_node_grid 1,2,1 \ 234*2e776fa0SMark Adams -petscpartitioner_simple_process_grid 2,1,1 \ 235*2e776fa0SMark Adams -petscpartitioner_type simple \ 236*2e776fa0SMark Adams -potential_petscspace_degree 2 \ 237*2e776fa0SMark Adams -snes_lag_jacobian -2 \ 238*2e776fa0SMark Adams -snes_max_it 1 \ 239*2e776fa0SMark Adams -snes_rtol 1.e-8 \ 240*2e776fa0SMark Adams -snes_type ksponly \ 241*2e776fa0SMark Adams -use_gpu_aware_mpi true 2425e1f5104SMark 2430c569c6eSMark test: 24418fb0606SStefano Zampini suffix: comparison 24518fb0606SStefano Zampini nsize: 4 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 \ 248*2e776fa0SMark Adams -dm_plex_simplex 0 -snes_type ksponly -dm_view -ksp_type cg -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 249*2e776fa0SMark Adams -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 25018fb0606SStefano Zampini 25118fb0606SStefano Zampini test: 2520c569c6eSMark suffix: cuda 2530c569c6eSMark nsize: 4 2540c569c6eSMark requires: cuda 25518fb0606SStefano Zampini output_file: output/ex13_comparison.out 25630602db0SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 257e600fa54SMatthew G. Knepley -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 258*2e776fa0SMark Adams -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 259*2e776fa0SMark Adams -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijcusparse -dm_vec_type cuda 2600c569c6eSMark 2610c569c6eSMark test: 26218fb0606SStefano Zampini suffix: kokkos_comp 26318fb0606SStefano Zampini nsize: 4 2643078479eSJunchao Zhang requires: !sycl kokkos_kernels 26518fb0606SStefano Zampini output_file: output/ex13_comparison.out 26630602db0SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 267e600fa54SMatthew G. Knepley -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 268*2e776fa0SMark Adams -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 269*2e776fa0SMark Adams -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijkokkos -dm_vec_type kokkos 27018fb0606SStefano Zampini 27118fb0606SStefano Zampini test: 2720c569c6eSMark nsize: 4 2733078479eSJunchao Zhang requires: !sycl kokkos_kernels 2740c569c6eSMark suffix: kokkos 275e600fa54SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 2,8 -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \ 276*2e776fa0SMark Adams -petscpartitioner_simple_node_grid 2,1 -dm_plex_simplex 0 -potential_petscspace_degree 1 -dm_refine 1 -ksp_type cg -pc_type gamg -ksp_norm_type unpreconditioned \ 277*2e776fa0SMark Adams -pc_gamg_esteig_ksp_type cg -ksp_converged_reason -snes_rtol 1.e-4 -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 278aa5a873eSStefano Zampini 279aa5a873eSStefano Zampini test: 280aa5a873eSStefano Zampini suffix: aijmkl_comp 281aa5a873eSStefano Zampini nsize: 4 282c4ad6305SSatish Balay requires: mkl_sparse 283aa5a873eSStefano Zampini output_file: output/ex13_comparison.out 28430602db0SMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 285e600fa54SMatthew G. Knepley -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple \ 286*2e776fa0SMark Adams -dm_plex_simplex 0 -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 287*2e776fa0SMark Adams -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl 288aa5a873eSStefano Zampini 289aa5a873eSStefano Zampini test: 290aa5a873eSStefano Zampini suffix: aijmkl_seq 291aa5a873eSStefano Zampini nsize: 1 292c4ad6305SSatish Balay requires: mkl_sparse 293aa5a873eSStefano Zampini TODO: broken (INDEFINITE PC) 294*2e776fa0SMark Adams 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 \ 29573f7197eSJed 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 \ 296*2e776fa0SMark Adams -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -pc_gamg_esteig_ksp_type cg -ksp_type cg -ksp_norm_type unpreconditioned \ 29773f7197eSJed Brown -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard 298aa5a873eSStefano Zampini 2995e1f5104SMark TEST*/ 300