xref: /petsc/src/snes/tests/ex13.c (revision dcfd994d081f5f6b19844c8e8bc67e93edfc619a)
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>
9e6f8f311SMark Adams #if defined(PETSC_HAVE_AMGX)
10e6f8f311SMark Adams #include <amgx_c.h>
11e6f8f311SMark Adams #endif
125e1f5104SMark 
135e1f5104SMark typedef struct {
14f9244615SMatthew G. Knepley   PetscInt  nit;    /* Number of benchmark iterations */
15f9244615SMatthew G. Knepley   PetscBool strong; /* Do not integrate the Laplacian by parts */
165e1f5104SMark } AppCtx;
175e1f5104SMark 
189371c9d4SSatish Balay static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
195e1f5104SMark   PetscInt d;
205e1f5104SMark   *u = 0.0;
215e1f5104SMark   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
225e1f5104SMark   return 0;
235e1f5104SMark }
245e1f5104SMark 
259371c9d4SSatish Balay 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[]) {
265e1f5104SMark   PetscInt d;
275e1f5104SMark   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
285e1f5104SMark }
295e1f5104SMark 
309371c9d4SSatish Balay 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[]) {
315e1f5104SMark   PetscInt d;
325e1f5104SMark   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
335e1f5104SMark }
345e1f5104SMark 
359371c9d4SSatish Balay 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[]) {
365e1f5104SMark   PetscInt d;
375e1f5104SMark   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
385e1f5104SMark }
395e1f5104SMark 
409371c9d4SSatish Balay static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
41f9244615SMatthew G. Knepley   *u = PetscSqr(x[0]) + PetscSqr(x[1]);
42f9244615SMatthew G. Knepley   return 0;
43f9244615SMatthew G. Knepley }
44f9244615SMatthew G. Knepley 
459371c9d4SSatish Balay static void f0_strong_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[]) {
46f9244615SMatthew G. Knepley   PetscInt d;
47f9244615SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] -= u_x[dim + d * dim + d];
48f9244615SMatthew G. Knepley   f0[0] += 4.0;
49f9244615SMatthew G. Knepley }
50f9244615SMatthew G. Knepley 
519371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
525e1f5104SMark   PetscFunctionBeginUser;
530c569c6eSMark   options->nit    = 10;
54f9244615SMatthew G. Knepley   options->strong = PETSC_FALSE;
55d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL));
58d0609cedSBarry Smith   PetscOptionsEnd();
595e1f5104SMark   PetscFunctionReturn(0);
605e1f5104SMark }
615e1f5104SMark 
629371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
635e1f5104SMark   PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
659566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
669566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
679566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
689566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
695e1f5104SMark   PetscFunctionReturn(0);
705e1f5104SMark }
715e1f5104SMark 
729371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) {
73f9244615SMatthew G. Knepley   PetscDS        ds;
7445480ffeSMatthew G. Knepley   DMLabel        label;
755e1f5104SMark   const PetscInt id = 1;
765e1f5104SMark 
775e1f5104SMark   PetscFunctionBeginUser;
789566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
799566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
80f9244615SMatthew G. Knepley   if (user->strong) {
819566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_strong_u, NULL));
829566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user));
839566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u, NULL, user, NULL));
84f9244615SMatthew G. Knepley   } else {
859566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
869566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
879566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
889566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
89f9244615SMatthew G. Knepley   }
905e1f5104SMark   PetscFunctionReturn(0);
915e1f5104SMark }
925e1f5104SMark 
939371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) {
945e1f5104SMark   DM             cdm = dm;
955e1f5104SMark   PetscFE        fe;
965e1f5104SMark   DMPolytopeType ct;
975e1f5104SMark   PetscBool      simplex;
985e1f5104SMark   PetscInt       dim, cStart;
995e1f5104SMark   char           prefix[PETSC_MAX_PATH_LEN];
1005e1f5104SMark 
1015e1f5104SMark   PetscFunctionBeginUser;
1029566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1039566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1052e776fa0SMark Adams   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; // false
1065e1f5104SMark   /* Create finite element */
1079566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1089566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
1105e1f5104SMark   /* Set discretization and boundary conditions for each mesh */
1119566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1129566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1139566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
1145e1f5104SMark   while (cdm) {
1159566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
1165e1f5104SMark     /* TODO: Check whether the boundary of coarse meshes is marked */
1179566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
1185e1f5104SMark   }
1199566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1205e1f5104SMark   PetscFunctionReturn(0);
1215e1f5104SMark }
1225e1f5104SMark 
1239371c9d4SSatish Balay int main(int argc, char **argv) {
1245e1f5104SMark   DM             dm;   /* Problem specification */
1255e1f5104SMark   SNES           snes; /* Nonlinear solver */
1265e1f5104SMark   Vec            u;    /* Solutions */
1275e1f5104SMark   AppCtx         user; /* User-defined work context */
1282e776fa0SMark Adams   PetscLogDouble time;
1292e776fa0SMark Adams   Mat            Amat;
1305e1f5104SMark 
131327415f7SBarry Smith   PetscFunctionBeginUser;
1329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1339566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1342e776fa0SMark Adams   /* system */
1359566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1369566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1379566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
1389566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
1399566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
1402e776fa0SMark Adams   PetscCall(SNESSetFromOptions(snes));
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
1429566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
1439566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
1442e776fa0SMark Adams   PetscCall(PetscTime(&time));
1452e776fa0SMark Adams   PetscCall(SNESSetUp(snes));
146e6f8f311SMark Adams #if defined(PETSC_HAVE_AMGX)
147e6f8f311SMark Adams   KSP                   ksp;
148e6f8f311SMark Adams   PC                    pc;
149e6f8f311SMark Adams   PetscBool             flg;
150e6f8f311SMark Adams   AMGX_resources_handle rsc;
151e6f8f311SMark Adams   PetscCall(SNESGetKSP(snes, &ksp));
152e6f8f311SMark Adams   PetscCall(KSPGetPC(ksp, &pc));
153e6f8f311SMark Adams   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCAMGX, &flg));
154e6f8f311SMark Adams   if (flg) {
155e6f8f311SMark Adams     PetscCall(PCAmgXGetResources(pc, (void *)&rsc));
156e6f8f311SMark Adams     /* do ... with resource */
157e6f8f311SMark Adams   }
158e6f8f311SMark Adams #endif
1592e776fa0SMark Adams   PetscCall(SNESGetJacobian(snes, &Amat, NULL, NULL, NULL));
1602e776fa0SMark Adams   PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE));
161b94d7dedSBarry Smith   PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE));
1629566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1632e776fa0SMark Adams   PetscCall(PetscTimeSubtract(&time));
1645e1f5104SMark   /* Benchmark system */
1650c569c6eSMark   if (user.nit) {
1665e1f5104SMark     Vec      b;
1670c569c6eSMark     PetscInt i;
1682e776fa0SMark Adams #if defined(PETSC_USE_LOG)
1692e776fa0SMark Adams     PetscLogStage kspstage;
1702e776fa0SMark Adams #endif
1712e776fa0SMark Adams     PetscCall(PetscLogStageRegister("Solve only", &kspstage));
1729566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(kspstage));
1732e776fa0SMark Adams     PetscCall(SNESGetSolution(snes, &u));
1742e776fa0SMark Adams     PetscCall(SNESGetFunction(snes, &b, NULL, NULL));
1750c569c6eSMark     for (i = 0; i < user.nit; i++) {
1769566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(u));
1772e776fa0SMark Adams       PetscCall(SNESSolve(snes, NULL, u));
1780c569c6eSMark     }
1799566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
1805e1f5104SMark   }
1819566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
1829566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
1835e1f5104SMark   /* Cleanup */
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1859566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
188b122ec5aSJacob Faibussowitsch   return 0;
1895e1f5104SMark }
1905e1f5104SMark 
1915e1f5104SMark /*TEST
1925e1f5104SMark 
1935e1f5104SMark   test:
194f9244615SMatthew G. Knepley     suffix: strong
195f9244615SMatthew G. Knepley     requires: triangle
19686081d6eSMark Adams     args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong
197f9244615SMatthew G. Knepley 
198f9244615SMatthew G. Knepley   test:
1995e1f5104SMark     suffix: bench
2000c569c6eSMark     nsize: 4
2012e776fa0SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1 -dm_refine 2 -dm_view -ksp_monitor \
20286081d6eSMark Adams        -benchmark_it 1 -dm_plex_box_upper 2,2,1 -dm_plex_box_lower 0,0,0 -dm_plex_dim 3 -ksp_converged_reason \
20386081d6eSMark Adams        -ksp_norm_type unpreconditioned -ksp_rtol 1.e-6 -ksp_type cg -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 \
20486081d6eSMark Adams        -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev  -mg_levels_pc_type jacobi -pc_gamg_coarse_eq_limit 200 \
20586081d6eSMark Adams        -pc_gamg_coarse_grid_layout_type compact -pc_gamg_esteig_ksp_max_it 5 -pc_gamg_process_eq_limit 200 \
20686081d6eSMark Adams        -pc_gamg_repartition false -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 0 -pc_gamg_threshold 0.001 -pc_gamg_threshold_scale .5 \
20786081d6eSMark Adams        -pc_gamg_type agg -pc_type gamg -petscpartitioner_simple_node_grid 1,2,1 -petscpartitioner_simple_process_grid 2,1,1 \
20886081d6eSMark Adams        -petscpartitioner_type simple -potential_petscspace_degree 2 -snes_lag_jacobian -2 -snes_max_it 1 -snes_rtol 1.e-8 -snes_type ksponly -use_gpu_aware_mpi true
2095e1f5104SMark 
21086081d6eSMark Adams   testset:
21186081d6eSMark Adams     nsize: 4
21286081d6eSMark Adams     output_file: output/ex13_comparison.out
21386081d6eSMark Adams     args: -dm_plex_dim 2 -benchmark_it 10 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
21486081d6eSMark Adams       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple  \
21586081d6eSMark Adams       -dm_plex_simplex 0 -snes_type ksponly -dm_view -ksp_type cg -pc_type gamg -pc_gamg_process_eq_limit 400 \
21686081d6eSMark Adams       -ksp_norm_type unpreconditioned -ksp_converged_reason
2170c569c6eSMark     test:
21818fb0606SStefano Zampini       suffix: comparison
21918fb0606SStefano Zampini     test:
2200c569c6eSMark       suffix: cuda
2210c569c6eSMark       requires: cuda
22286081d6eSMark Adams       args: -dm_mat_type aijcusparse -dm_vec_type cuda
2230c569c6eSMark     test:
2240c569c6eSMark       suffix: kokkos
225*dcfd994dSJunchao Zhang       requires: sycl kokkos_kernels
22686081d6eSMark Adams       args: -dm_mat_type aijkokkos -dm_vec_type kokkos
227aa5a873eSStefano Zampini     test:
228aa5a873eSStefano Zampini       suffix: aijmkl_comp
229c4ad6305SSatish Balay       requires: mkl_sparse
23086081d6eSMark Adams       args: -dm_mat_type aijmkl
231aa5a873eSStefano Zampini 
232aa5a873eSStefano Zampini   test:
233aa5a873eSStefano Zampini     suffix: aijmkl_seq
234aa5a873eSStefano Zampini     nsize: 1
235c4ad6305SSatish Balay     requires: mkl_sparse
236aa5a873eSStefano Zampini     TODO: broken (INDEFINITE PC)
2372e776fa0SMark 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 \
238bae903cbSmarkadams4           -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_threshold -1 -pc_gamg_square_graph 10 -pc_gamg_process_eq_limit 400 \
2392e776fa0SMark Adams           -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -pc_gamg_esteig_ksp_type cg -ksp_type cg -ksp_norm_type unpreconditioned \
24073f7197eSJed Brown           -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard
241aa5a873eSStefano Zampini 
242e6f8f311SMark Adams   testset:
243a22370e2Smarkadams4     requires: cuda amgx
244a22370e2Smarkadams4     filter: grep -v Built | grep -v "AMGX version" | grep -v "CUDA Runtime"
245e6f8f311SMark Adams     output_file: output/ex13_amgx.out
246e6f8f311SMark Adams     args: -dm_plex_dim 2 -dm_plex_box_faces 2,2 -dm_refine 2 -petscpartitioner_type simple -potential_petscspace_degree 2 -dm_plex_simplex 0 -ksp_monitor \
247a22370e2Smarkadams4           -snes_type ksponly -dm_view -ksp_type cg -ksp_norm_type unpreconditioned -ksp_converged_reason -snes_rtol 1.e-4 -pc_type amgx -benchmark_it 1 -pc_amgx_verbose false
248e6f8f311SMark Adams     nsize: 4
249e6f8f311SMark Adams     test:
250e6f8f311SMark Adams       suffix: amgx
251e6f8f311SMark Adams       args: -dm_mat_type aijcusparse -dm_vec_type cuda
252e6f8f311SMark Adams     test:
253e6f8f311SMark Adams       suffix: amgx_cpu
254e6f8f311SMark Adams       args: -dm_mat_type aij
255e6f8f311SMark Adams 
2565e1f5104SMark TEST*/
257