xref: /petsc/src/snes/tests/ex13.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
695e1f5104SMark   PetscFunctionBeginUser;
700c569c6eSMark   options->nit    = 10;
71f9244615SMatthew G. Knepley   options->strong = PETSC_FALSE;
72d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL));
75d0609cedSBarry Smith   PetscOptionsEnd();
765e1f5104SMark   PetscFunctionReturn(0);
775e1f5104SMark }
785e1f5104SMark 
795e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
805e1f5104SMark {
815e1f5104SMark   PetscFunctionBeginUser;
829566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
839566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
849566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
859566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
869566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
875e1f5104SMark   PetscFunctionReturn(0);
885e1f5104SMark }
895e1f5104SMark 
905e1f5104SMark static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
915e1f5104SMark {
92f9244615SMatthew G. Knepley   PetscDS        ds;
9345480ffeSMatthew G. Knepley   DMLabel        label;
945e1f5104SMark   const PetscInt id = 1;
955e1f5104SMark 
965e1f5104SMark   PetscFunctionBeginUser;
979566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
989566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
99f9244615SMatthew G. Knepley   if (user->strong) {
1009566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_strong_u, NULL));
1019566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user));
1029566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, user, NULL));
103f9244615SMatthew G. Knepley   } else {
1049566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u));
1059566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1069566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
1079566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL));
108f9244615SMatthew G. Knepley   }
1095e1f5104SMark   PetscFunctionReturn(0);
1105e1f5104SMark }
1115e1f5104SMark 
1125e1f5104SMark static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
1135e1f5104SMark {
1145e1f5104SMark   DM             cdm = dm;
1155e1f5104SMark   PetscFE        fe;
1165e1f5104SMark   DMPolytopeType ct;
1175e1f5104SMark   PetscBool      simplex;
1185e1f5104SMark   PetscInt       dim, cStart;
1195e1f5104SMark   char           prefix[PETSC_MAX_PATH_LEN];
1205e1f5104SMark 
1215e1f5104SMark   PetscFunctionBeginUser;
1229566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1239566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1249566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1252e776fa0SMark Adams   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; // false
1265e1f5104SMark   /* Create finite element */
1279566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1289566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
1299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, name));
1305e1f5104SMark   /* Set discretization and boundary conditions for each mesh */
1319566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
1329566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1339566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
1345e1f5104SMark   while (cdm) {
1359566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm,cdm));
1365e1f5104SMark     /* TODO: Check whether the boundary of coarse meshes is marked */
1379566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
1385e1f5104SMark   }
1399566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1405e1f5104SMark   PetscFunctionReturn(0);
1415e1f5104SMark }
1425e1f5104SMark 
1435e1f5104SMark int main(int argc, char **argv)
1445e1f5104SMark {
1455e1f5104SMark   DM             dm;   /* Problem specification */
1465e1f5104SMark   SNES           snes; /* Nonlinear solver */
1475e1f5104SMark   Vec            u;    /* Solutions */
1485e1f5104SMark   AppCtx         user; /* User-defined work context */
1492e776fa0SMark Adams   PetscLogDouble time;
1502e776fa0SMark Adams   Mat            Amat;
1515e1f5104SMark 
152*327415f7SBarry Smith   PetscFunctionBeginUser;
1539566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
1549566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1552e776fa0SMark 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));
1612e776fa0SMark 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));
1652e776fa0SMark Adams   PetscCall(PetscTime(&time));
1662e776fa0SMark Adams   PetscCall(SNESSetUp(snes));
1672e776fa0SMark Adams   PetscCall(SNESGetJacobian(snes, &Amat, NULL, NULL, NULL));
1682e776fa0SMark Adams   PetscCall(MatSetOption(Amat,MAT_SPD,PETSC_TRUE));
169b94d7dedSBarry Smith   PetscCall(MatSetOption(Amat,MAT_SPD_ETERNAL,PETSC_TRUE));
1709566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1712e776fa0SMark Adams   PetscCall(PetscTimeSubtract(&time));
1722e776fa0SMark Adams   // PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First Solve time: %g\n",-time));
1735e1f5104SMark   /* Benchmark system */
1740c569c6eSMark   if (user.nit) {
1755e1f5104SMark     Vec            b;
1760c569c6eSMark     PetscInt       i;
1772e776fa0SMark Adams #if defined(PETSC_USE_LOG)
1782e776fa0SMark Adams     PetscLogStage  kspstage;
1792e776fa0SMark Adams #endif
1802e776fa0SMark Adams     PetscCall(PetscLogStageRegister("Solve only", &kspstage));
1819566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(kspstage));
1822e776fa0SMark Adams     PetscCall(SNESGetSolution(snes, &u));
1832e776fa0SMark Adams     PetscCall(SNESGetFunction(snes, &b, NULL, NULL));
1840c569c6eSMark     for (i=0;i<user.nit;i++) {
1859566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(u));
1862e776fa0SMark Adams       PetscCall(SNESSolve(snes, NULL, u));
1870c569c6eSMark     }
1889566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
1895e1f5104SMark   }
1909566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
1919566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
1925e1f5104SMark   /* Cleanup */
1939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1949566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1959566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
197b122ec5aSJacob Faibussowitsch   return 0;
1985e1f5104SMark }
1995e1f5104SMark 
2005e1f5104SMark /*TEST
2015e1f5104SMark 
2025e1f5104SMark   test:
203f9244615SMatthew G. Knepley     suffix: strong
204f9244615SMatthew G. Knepley     requires: triangle
20586081d6eSMark Adams     args: -dm_plex_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong
206f9244615SMatthew G. Knepley 
207f9244615SMatthew G. Knepley   test:
2085e1f5104SMark     suffix: bench
2090c569c6eSMark     nsize: 4
2102e776fa0SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1 -dm_refine 2 -dm_view -ksp_monitor \
21186081d6eSMark 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 \
21286081d6eSMark Adams        -ksp_norm_type unpreconditioned -ksp_rtol 1.e-6 -ksp_type cg -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 \
21386081d6eSMark Adams        -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev  -mg_levels_pc_type jacobi -pc_gamg_coarse_eq_limit 200 \
21486081d6eSMark Adams        -pc_gamg_coarse_grid_layout_type compact -pc_gamg_esteig_ksp_max_it 5 -pc_gamg_process_eq_limit 200 \
21586081d6eSMark 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 \
21686081d6eSMark Adams        -pc_gamg_type agg -pc_type gamg -petscpartitioner_simple_node_grid 1,2,1 -petscpartitioner_simple_process_grid 2,1,1 \
21786081d6eSMark 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
2185e1f5104SMark 
21986081d6eSMark Adams   testset:
22086081d6eSMark Adams     nsize: 4
22186081d6eSMark Adams     output_file: output/ex13_comparison.out
22286081d6eSMark Adams     args: -dm_plex_dim 2 -benchmark_it 10 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
22386081d6eSMark Adams       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -petscpartitioner_type simple  \
22486081d6eSMark Adams       -dm_plex_simplex 0 -snes_type ksponly -dm_view -ksp_type cg -pc_type gamg -pc_gamg_process_eq_limit 400 \
22586081d6eSMark Adams       -ksp_norm_type unpreconditioned -ksp_converged_reason
2260c569c6eSMark     test:
22718fb0606SStefano Zampini       suffix: comparison
22818fb0606SStefano Zampini     test:
2290c569c6eSMark       suffix: cuda
2300c569c6eSMark       requires: cuda
23186081d6eSMark Adams       args: -dm_mat_type aijcusparse -dm_vec_type cuda
2320c569c6eSMark     test:
2330c569c6eSMark       suffix: kokkos
23486081d6eSMark Adams       requires: !sycl kokkos_kernels
23586081d6eSMark Adams       args: -dm_mat_type aijkokkos -dm_vec_type kokkos
236aa5a873eSStefano Zampini     test:
237aa5a873eSStefano Zampini       suffix: aijmkl_comp
238c4ad6305SSatish Balay       requires: mkl_sparse
23986081d6eSMark Adams       args: -dm_mat_type aijmkl
240aa5a873eSStefano Zampini 
241aa5a873eSStefano Zampini   test:
242aa5a873eSStefano Zampini     suffix: aijmkl_seq
243aa5a873eSStefano Zampini     nsize: 1
244c4ad6305SSatish Balay     requires: mkl_sparse
245aa5a873eSStefano Zampini     TODO: broken (INDEFINITE PC)
2462e776fa0SMark 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 \
247bae903cbSmarkadams4           -snes_type ksponly -dm_view -pc_type gamg -pc_gamg_threshold -1 -pc_gamg_square_graph 10 -pc_gamg_process_eq_limit 400 \
2482e776fa0SMark Adams           -pc_gamg_reuse_interpolation -pc_gamg_coarse_eq_limit 10 -pc_gamg_esteig_ksp_type cg -ksp_type cg -ksp_norm_type unpreconditioned \
24973f7197eSJed Brown           -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl -dm_vec_type standard
250aa5a873eSStefano Zampini 
2515e1f5104SMark TEST*/
252