xref: /petsc/src/snes/tests/ex13.c (revision 18fb0606e5ee84c32740aa67a8844f858bb2d3af)
15e1f5104SMark static char help[] = "Benchmark Poisson Problem in 2d and 3d with finite elements.\n\
25e1f5104SMark We solve the Poisson problem in a rectangular\n\
35e1f5104SMark domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
45e1f5104SMark This example supports automatic convergence estimation\n\
55e1f5104SMark and eventually adaptivity.\n\n\n";
65e1f5104SMark 
75e1f5104SMark #include <petscdmplex.h>
85e1f5104SMark #include <petscsnes.h>
95e1f5104SMark #include <petscds.h>
105e1f5104SMark #include <petscconvest.h>
115e1f5104SMark 
125e1f5104SMark typedef struct {
130c569c6eSMark   PetscInt  nit;
145e1f5104SMark } AppCtx;
155e1f5104SMark 
165e1f5104SMark 
175e1f5104SMark static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
185e1f5104SMark {
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 
255e1f5104SMark static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
265e1f5104SMark                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
275e1f5104SMark                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
285e1f5104SMark                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
295e1f5104SMark {
305e1f5104SMark   PetscInt d;
315e1f5104SMark   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
325e1f5104SMark }
335e1f5104SMark 
345e1f5104SMark static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
355e1f5104SMark                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
365e1f5104SMark                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
375e1f5104SMark                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
385e1f5104SMark {
395e1f5104SMark   PetscInt d;
405e1f5104SMark   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
415e1f5104SMark }
425e1f5104SMark 
435e1f5104SMark static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
445e1f5104SMark                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
455e1f5104SMark                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
465e1f5104SMark                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
475e1f5104SMark {
485e1f5104SMark   PetscInt d;
495e1f5104SMark   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
505e1f5104SMark }
515e1f5104SMark 
525e1f5104SMark static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
535e1f5104SMark {
545e1f5104SMark   PetscErrorCode ierr;
555e1f5104SMark 
565e1f5104SMark   PetscFunctionBeginUser;
570c569c6eSMark   options->nit = 10;
585e1f5104SMark   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
590c569c6eSMark   ierr = PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL);CHKERRQ(ierr);
605e1f5104SMark   ierr = PetscOptionsEnd();
615e1f5104SMark   PetscFunctionReturn(0);
625e1f5104SMark }
635e1f5104SMark 
645e1f5104SMark static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
655e1f5104SMark {
665e1f5104SMark   PetscErrorCode ierr;
675e1f5104SMark 
685e1f5104SMark   PetscFunctionBeginUser;
695e1f5104SMark   /* Create box mesh */
705e1f5104SMark   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
715e1f5104SMark   /* TODO: This should be pulled into the library */
725e1f5104SMark   {
735e1f5104SMark     char      convType[256];
745e1f5104SMark     PetscBool flg;
755e1f5104SMark 
765e1f5104SMark     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
775e1f5104SMark     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
785e1f5104SMark     ierr = PetscOptionsEnd();
795e1f5104SMark     if (flg) {
805e1f5104SMark       DM dmConv;
815e1f5104SMark 
825e1f5104SMark       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
835e1f5104SMark       if (dmConv) {
845e1f5104SMark         ierr = DMDestroy(dm);CHKERRQ(ierr);
855e1f5104SMark         *dm  = dmConv;
865e1f5104SMark       }
875e1f5104SMark     }
885e1f5104SMark   }
895e1f5104SMark   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
905e1f5104SMark 
915e1f5104SMark   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
925e1f5104SMark   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
935e1f5104SMark   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
945e1f5104SMark   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
955e1f5104SMark   PetscFunctionReturn(0);
965e1f5104SMark }
975e1f5104SMark 
985e1f5104SMark static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
995e1f5104SMark {
1005e1f5104SMark   PetscDS        prob;
1015e1f5104SMark   const PetscInt id = 1;
1025e1f5104SMark   PetscErrorCode ierr;
1035e1f5104SMark 
1045e1f5104SMark   PetscFunctionBeginUser;
1055e1f5104SMark   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1065e1f5104SMark   ierr = PetscDSSetResidual(prob, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
1075e1f5104SMark   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
1085e1f5104SMark   ierr = PetscDSSetExactSolution(prob, 0, trig_u, user);CHKERRQ(ierr);
1095e1f5104SMark   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) trig_u, NULL, 1, &id, user);CHKERRQ(ierr);
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   PetscErrorCode ierr;
1225e1f5104SMark 
1235e1f5104SMark   PetscFunctionBeginUser;
1245e1f5104SMark   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1255e1f5104SMark   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
1265e1f5104SMark   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
1275e1f5104SMark   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1285e1f5104SMark   /* Create finite element */
1295e1f5104SMark   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
1305e1f5104SMark   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
1315e1f5104SMark   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
1325e1f5104SMark   /* Set discretization and boundary conditions for each mesh */
1335e1f5104SMark   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
1345e1f5104SMark   ierr = DMCreateDS(dm);CHKERRQ(ierr);
1355e1f5104SMark   ierr = (*setup)(dm, user);CHKERRQ(ierr);
1365e1f5104SMark   while (cdm) {
1375e1f5104SMark     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
1385e1f5104SMark     /* TODO: Check whether the boundary of coarse meshes is marked */
1395e1f5104SMark     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
1405e1f5104SMark   }
1415e1f5104SMark   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1425e1f5104SMark   PetscFunctionReturn(0);
1435e1f5104SMark }
1445e1f5104SMark 
1455e1f5104SMark int main(int argc, char **argv)
1465e1f5104SMark {
1475e1f5104SMark   DM             dm;   /* Problem specification */
1485e1f5104SMark   SNES           snes; /* Nonlinear solver */
1495e1f5104SMark   Vec            u;    /* Solutions */
1505e1f5104SMark   AppCtx         user; /* User-defined work context */
1515e1f5104SMark   PetscErrorCode ierr;
1525e1f5104SMark 
1535e1f5104SMark   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
1545e1f5104SMark   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
1555e1f5104SMark   /* Primal system */
1565e1f5104SMark   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
1575e1f5104SMark   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
1585e1f5104SMark   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
1595e1f5104SMark   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
1605e1f5104SMark   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
1615e1f5104SMark   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
1625e1f5104SMark   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
1635e1f5104SMark   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
1645e1f5104SMark   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1655e1f5104SMark   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
1665e1f5104SMark   /* Benchmark system */
1670c569c6eSMark   if (user.nit) {
1685e1f5104SMark #if defined(PETSC_USE_LOG)
169*18fb0606SStefano Zampini     PetscLogStage kspstage,pcstage;
1705e1f5104SMark #endif
1715e1f5104SMark     KSP       ksp;
172*18fb0606SStefano Zampini     PC        pc;
173*18fb0606SStefano Zampini     Mat       A,P;
1745e1f5104SMark     Vec       b;
1750c569c6eSMark     PetscInt  i;
1760c569c6eSMark     ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr);
1775e1f5104SMark     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
1785e1f5104SMark     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
179*18fb0606SStefano Zampini     ierr = SNESGetJacobian(snes, &A, &P, NULL, NULL);CHKERRQ(ierr);
1800c569c6eSMark     ierr = VecSet(u, 0.0);CHKERRQ(ierr);
1815e1f5104SMark     ierr = SNESGetFunction(snes, &b, NULL, NULL);CHKERRQ(ierr);
1825e1f5104SMark     ierr = SNESComputeFunction(snes, u, b);CHKERRQ(ierr);
183*18fb0606SStefano Zampini     ierr = SNESComputeJacobian(snes, u, A, P);CHKERRQ(ierr);
184*18fb0606SStefano Zampini     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
185*18fb0606SStefano Zampini     ierr = PetscLogStageRegister("PCSetUp", &pcstage);CHKERRQ(ierr);
186*18fb0606SStefano Zampini     ierr = PetscLogStagePush(pcstage);CHKERRQ(ierr);
187*18fb0606SStefano Zampini     ierr = PCSetUp(pc);CHKERRQ(ierr);
188*18fb0606SStefano Zampini     ierr = PetscLogStagePop();CHKERRQ(ierr);
189*18fb0606SStefano Zampini     ierr = PetscLogStageRegister("KSP Solve only", &kspstage);CHKERRQ(ierr);
190*18fb0606SStefano Zampini     ierr = PetscLogStagePush(kspstage);CHKERRQ(ierr);
1910c569c6eSMark     for (i=0;i<user.nit;i++) {
1920c569c6eSMark       ierr = VecZeroEntries(u);CHKERRQ(ierr);
1935e1f5104SMark       ierr = KSPSolve(ksp, b, u);CHKERRQ(ierr);
1940c569c6eSMark     }
1955e1f5104SMark     ierr = PetscLogStagePop();CHKERRQ(ierr);
1965e1f5104SMark   }
1975e1f5104SMark   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
1985e1f5104SMark   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
1995e1f5104SMark   /* Cleanup */
2005e1f5104SMark   ierr = VecDestroy(&u);CHKERRQ(ierr);
2015e1f5104SMark   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
2025e1f5104SMark   ierr = DMDestroy(&dm);CHKERRQ(ierr);
2035e1f5104SMark   ierr = PetscFinalize();
2045e1f5104SMark   return ierr;
2055e1f5104SMark }
2065e1f5104SMark 
2075e1f5104SMark /*TEST
2085e1f5104SMark 
2095e1f5104SMark   test:
2105e1f5104SMark     suffix: bench
2110c569c6eSMark     nsize: 4
2120c569c6eSMark     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2,8 -dm_refine 1 -dm_distribute \
2130c569c6eSMark           -petscpartitioner_type simple -petscpartitioner_simple_process_grid 1,1,2 -petscpartitioner_simple_node_grid 1,1,2 \
2140c569c6eSMark           -potential_petscspace_degree 2 -ksp_type cg -pc_type gamg -benchmark_it 1 -dm_view -snes_rtol 1.e-4
2155e1f5104SMark 
2160c569c6eSMark   test:
217*18fb0606SStefano Zampini     suffix: comparison
218*18fb0606SStefano Zampini     nsize: 4
219*18fb0606SStefano Zampini     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
220*18fb0606SStefano Zampini       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
221*18fb0606SStefano 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 \
222*18fb0606SStefano Zampini       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4
223*18fb0606SStefano Zampini 
224*18fb0606SStefano Zampini   test:
2250c569c6eSMark     suffix: cuda
2260c569c6eSMark     nsize: 4
2270c569c6eSMark     requires: cuda
228*18fb0606SStefano Zampini     output_file: output/ex13_comparison.out
2290c569c6eSMark     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
2300c569c6eSMark       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
2310c569c6eSMark       -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 \
2320c569c6eSMark       -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 -mat_cusparse_transgen
2330c569c6eSMark 
2340c569c6eSMark   test:
235*18fb0606SStefano Zampini     suffix: kokkos_comp
236*18fb0606SStefano Zampini     nsize: 4
237*18fb0606SStefano Zampini     requires: kokkos_kernels
238*18fb0606SStefano Zampini     output_file: output/ex13_comparison.out
239*18fb0606SStefano Zampini     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
240*18fb0606SStefano Zampini       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
241*18fb0606SStefano 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 \
242*18fb0606SStefano 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
243*18fb0606SStefano Zampini 
244*18fb0606SStefano Zampini   test:
2450c569c6eSMark     nsize: 4
2460c569c6eSMark     requires: kokkos_kernels
2470c569c6eSMark     suffix: kokkos
2480c569c6eSMark     args: -dm_plex_box_dim 2 -dm_plex_box_faces 2,8 -dm_distribute -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \
2490c569c6eSMark           -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 \
2500c569c6eSMark           -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
2515e1f5104SMark TEST*/
252