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); 740c569c6eSMark ierr = PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL);CHKERRQ(ierr); 75f9244615SMatthew G. Knepley ierr = PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL);CHKERRQ(ierr); 765e1f5104SMark ierr = PetscOptionsEnd(); 775e1f5104SMark PetscFunctionReturn(0); 785e1f5104SMark } 795e1f5104SMark 805e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 815e1f5104SMark { 825e1f5104SMark PetscErrorCode ierr; 835e1f5104SMark 845e1f5104SMark PetscFunctionBeginUser; 855e1f5104SMark /* Create box mesh */ 865e1f5104SMark ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 875e1f5104SMark /* TODO: This should be pulled into the library */ 885e1f5104SMark { 895e1f5104SMark char convType[256]; 905e1f5104SMark PetscBool flg; 915e1f5104SMark 925e1f5104SMark ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 935e1f5104SMark ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 945e1f5104SMark ierr = PetscOptionsEnd(); 955e1f5104SMark if (flg) { 965e1f5104SMark DM dmConv; 975e1f5104SMark 985e1f5104SMark ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr); 995e1f5104SMark if (dmConv) { 1005e1f5104SMark ierr = DMDestroy(dm);CHKERRQ(ierr); 1015e1f5104SMark *dm = dmConv; 1025e1f5104SMark } 1035e1f5104SMark } 1045e1f5104SMark } 1055e1f5104SMark ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 1065e1f5104SMark 1075e1f5104SMark ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 1085e1f5104SMark ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 1095e1f5104SMark ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 1105e1f5104SMark ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 1115e1f5104SMark PetscFunctionReturn(0); 1125e1f5104SMark } 1135e1f5104SMark 1145e1f5104SMark static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 1155e1f5104SMark { 116f9244615SMatthew G. Knepley PetscDS ds; 117*45480ffeSMatthew G. Knepley DMLabel label; 1185e1f5104SMark const PetscInt id = 1; 1195e1f5104SMark PetscErrorCode ierr; 1205e1f5104SMark 1215e1f5104SMark PetscFunctionBeginUser; 122f9244615SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 123*45480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 124f9244615SMatthew G. Knepley if (user->strong) { 125f9244615SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_strong_u, NULL);CHKERRQ(ierr); 126f9244615SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quadratic_u, user);CHKERRQ(ierr); 127*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, user, NULL);CHKERRQ(ierr); 128f9244615SMatthew G. Knepley } else { 129f9244615SMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 130f9244615SMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 131f9244615SMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr); 132*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr); 133f9244615SMatthew G. Knepley } 1345e1f5104SMark PetscFunctionReturn(0); 1355e1f5104SMark } 1365e1f5104SMark 1375e1f5104SMark static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 1385e1f5104SMark { 1395e1f5104SMark DM cdm = dm; 1405e1f5104SMark PetscFE fe; 1415e1f5104SMark DMPolytopeType ct; 1425e1f5104SMark PetscBool simplex; 1435e1f5104SMark PetscInt dim, cStart; 1445e1f5104SMark char prefix[PETSC_MAX_PATH_LEN]; 1455e1f5104SMark PetscErrorCode ierr; 1465e1f5104SMark 1475e1f5104SMark PetscFunctionBeginUser; 1485e1f5104SMark ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1495e1f5104SMark ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 1505e1f5104SMark ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 1515e1f5104SMark simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 1525e1f5104SMark /* Create finite element */ 1535e1f5104SMark ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 15480ab8e62SMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 1555e1f5104SMark ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 1565e1f5104SMark /* Set discretization and boundary conditions for each mesh */ 1575e1f5104SMark ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 1585e1f5104SMark ierr = DMCreateDS(dm);CHKERRQ(ierr); 1595e1f5104SMark ierr = (*setup)(dm, user);CHKERRQ(ierr); 1605e1f5104SMark while (cdm) { 1615e1f5104SMark ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 1625e1f5104SMark /* TODO: Check whether the boundary of coarse meshes is marked */ 1635e1f5104SMark ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 1645e1f5104SMark } 1655e1f5104SMark ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 1665e1f5104SMark PetscFunctionReturn(0); 1675e1f5104SMark } 1685e1f5104SMark 1695e1f5104SMark int main(int argc, char **argv) 1705e1f5104SMark { 1715e1f5104SMark DM dm; /* Problem specification */ 1725e1f5104SMark SNES snes; /* Nonlinear solver */ 1735e1f5104SMark Vec u; /* Solutions */ 1745e1f5104SMark AppCtx user; /* User-defined work context */ 1755e1f5104SMark PetscErrorCode ierr; 1765e1f5104SMark 1775e1f5104SMark ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 1785e1f5104SMark ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 1795e1f5104SMark /* Primal system */ 1805e1f5104SMark ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 1815e1f5104SMark ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 1825e1f5104SMark ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 1835e1f5104SMark ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 1845e1f5104SMark ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 1855e1f5104SMark ierr = VecSet(u, 0.0);CHKERRQ(ierr); 1865e1f5104SMark ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 1875e1f5104SMark ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 1885e1f5104SMark ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 189f9244615SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 1905e1f5104SMark ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 1915e1f5104SMark /* Benchmark system */ 1920c569c6eSMark if (user.nit) { 1935e1f5104SMark #if defined(PETSC_USE_LOG) 19418fb0606SStefano Zampini PetscLogStage kspstage,pcstage; 1955e1f5104SMark #endif 1965e1f5104SMark KSP ksp; 19718fb0606SStefano Zampini PC pc; 19818fb0606SStefano Zampini Mat A,P; 1995e1f5104SMark Vec b; 2000c569c6eSMark PetscInt i; 2010c569c6eSMark ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr); 2025e1f5104SMark ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 2035e1f5104SMark ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 20418fb0606SStefano Zampini ierr = SNESGetJacobian(snes, &A, &P, NULL, NULL);CHKERRQ(ierr); 2050c569c6eSMark ierr = VecSet(u, 0.0);CHKERRQ(ierr); 2065e1f5104SMark ierr = SNESGetFunction(snes, &b, NULL, NULL);CHKERRQ(ierr); 2075e1f5104SMark ierr = SNESComputeFunction(snes, u, b);CHKERRQ(ierr); 20818fb0606SStefano Zampini ierr = SNESComputeJacobian(snes, u, A, P);CHKERRQ(ierr); 20918fb0606SStefano Zampini ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 21018fb0606SStefano Zampini ierr = PetscLogStageRegister("PCSetUp", &pcstage);CHKERRQ(ierr); 21118fb0606SStefano Zampini ierr = PetscLogStagePush(pcstage);CHKERRQ(ierr); 21218fb0606SStefano Zampini ierr = PCSetUp(pc);CHKERRQ(ierr); 21318fb0606SStefano Zampini ierr = PetscLogStagePop();CHKERRQ(ierr); 21418fb0606SStefano Zampini ierr = PetscLogStageRegister("KSP Solve only", &kspstage);CHKERRQ(ierr); 21518fb0606SStefano Zampini ierr = PetscLogStagePush(kspstage);CHKERRQ(ierr); 2160c569c6eSMark for (i=0;i<user.nit;i++) { 2170c569c6eSMark ierr = VecZeroEntries(u);CHKERRQ(ierr); 2185e1f5104SMark ierr = KSPSolve(ksp, b, u);CHKERRQ(ierr); 2190c569c6eSMark } 2205e1f5104SMark ierr = PetscLogStagePop();CHKERRQ(ierr); 2215e1f5104SMark } 2225e1f5104SMark ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 2235e1f5104SMark ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 2245e1f5104SMark /* Cleanup */ 2255e1f5104SMark ierr = VecDestroy(&u);CHKERRQ(ierr); 2265e1f5104SMark ierr = SNESDestroy(&snes);CHKERRQ(ierr); 2275e1f5104SMark ierr = DMDestroy(&dm);CHKERRQ(ierr); 2285e1f5104SMark ierr = PetscFinalize(); 2295e1f5104SMark return ierr; 2305e1f5104SMark } 2315e1f5104SMark 2325e1f5104SMark /*TEST 2335e1f5104SMark 2345e1f5104SMark test: 235f9244615SMatthew G. Knepley suffix: strong 236f9244615SMatthew G. Knepley requires: triangle 237f9244615SMatthew G. Knepley args: -dm_plex_box_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \ 238f9244615SMatthew G. Knepley -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong 239f9244615SMatthew G. Knepley 240f9244615SMatthew G. Knepley test: 2415e1f5104SMark suffix: bench 2420c569c6eSMark nsize: 4 2430c569c6eSMark args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2,8 -dm_refine 1 -dm_distribute \ 2440c569c6eSMark -petscpartitioner_type simple -petscpartitioner_simple_process_grid 1,1,2 -petscpartitioner_simple_node_grid 1,1,2 \ 2450c569c6eSMark -potential_petscspace_degree 2 -ksp_type cg -pc_type gamg -benchmark_it 1 -dm_view -snes_rtol 1.e-4 2465e1f5104SMark 2470c569c6eSMark test: 24818fb0606SStefano Zampini suffix: comparison 24918fb0606SStefano Zampini nsize: 4 25018fb0606SStefano Zampini args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 25118fb0606SStefano Zampini -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \ 25218fb0606SStefano Zampini -dm_plex_box_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 25318fb0606SStefano Zampini -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 25418fb0606SStefano Zampini 25518fb0606SStefano Zampini test: 2560c569c6eSMark suffix: cuda 2570c569c6eSMark nsize: 4 2580c569c6eSMark requires: cuda 25918fb0606SStefano Zampini output_file: output/ex13_comparison.out 2600c569c6eSMark args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 2610c569c6eSMark -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \ 2620c569c6eSMark -dm_plex_box_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 263e589036eSStefano 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 2640c569c6eSMark 2650c569c6eSMark test: 26618fb0606SStefano Zampini suffix: kokkos_comp 26718fb0606SStefano Zampini nsize: 4 26818fb0606SStefano Zampini requires: kokkos_kernels 26918fb0606SStefano Zampini output_file: output/ex13_comparison.out 27018fb0606SStefano Zampini args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 27118fb0606SStefano Zampini -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \ 27218fb0606SStefano Zampini -dm_plex_box_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 27318fb0606SStefano 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 27418fb0606SStefano Zampini 27518fb0606SStefano Zampini test: 2760c569c6eSMark nsize: 4 2770c569c6eSMark requires: kokkos_kernels 2780c569c6eSMark suffix: kokkos 2790c569c6eSMark args: -dm_plex_box_dim 2 -dm_plex_box_faces 2,8 -dm_distribute -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \ 2800c569c6eSMark -petscpartitioner_simple_node_grid 2,1 -dm_plex_box_simplex 0 -potential_petscspace_degree 1 -dm_refine 1 -ksp_type cg -pc_type gamg -ksp_norm_type unpreconditioned \ 2810c569c6eSMark -mg_levels_esteig_ksp_type cg -mg_levels_pc_type jacobi -ksp_converged_reason -snes_monitor_short -snes_rtol 1.e-4 -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 282aa5a873eSStefano Zampini 283aa5a873eSStefano Zampini test: 284aa5a873eSStefano Zampini suffix: aijmkl_comp 285aa5a873eSStefano Zampini nsize: 4 286aa5a873eSStefano Zampini requires: mkl 287aa5a873eSStefano Zampini output_file: output/ex13_comparison.out 288aa5a873eSStefano Zampini args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \ 289aa5a873eSStefano Zampini -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \ 290aa5a873eSStefano Zampini -dm_plex_box_simplex 0 -snes_monitor_short -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_process_eq_limit 400 -ksp_norm_type unpreconditioned \ 291aa5a873eSStefano Zampini -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl 292aa5a873eSStefano Zampini 293aa5a873eSStefano Zampini test: 294aa5a873eSStefano Zampini suffix: aijmkl_seq 295aa5a873eSStefano Zampini nsize: 1 296aa5a873eSStefano Zampini requires: mkl 297aa5a873eSStefano Zampini TODO: broken (INDEFINITE PC) 298aa5a873eSStefano Zampini args: -dm_plex_box_dim 3 -dm_plex_box_faces 4,4,4 -dm_refine 1 -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_distribute -dm_plex_box_simplex 0 -snes_monitor_short -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 -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -mg_levels_esteig_ksp_type cg -mg_levels_pc_type jacobi -ksp_type cg -ksp_norm_type unpreconditioned -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard 299aa5a873eSStefano Zampini 3005e1f5104SMark TEST*/ 301