xref: /petsc/src/snes/tests/ex13.c (revision efa12513287cff49a2b9648ae83199dcbfaad71a)
1 static char help[] = "Benchmark Poisson Problem in 2d and 3d with finite elements.\n\
2 We solve the Poisson problem in a rectangular domain\n\
3 using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4 
5 #include <petscdmplex.h>
6 #include <petscsnes.h>
7 #include <petscds.h>
8 #include <petscconvest.h>
9 
10 typedef struct {
11   PetscInt  nit;    /* Number of benchmark iterations */
12   PetscBool strong; /* Do not integrate the Laplacian by parts */
13 } AppCtx;
14 
15 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
16 {
17   PetscInt d;
18   *u = 0.0;
19   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
20   return 0;
21 }
22 
23 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
24                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
25                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
26                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
27 {
28   PetscInt d;
29   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
30 }
31 
32 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
33                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
34                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
35                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
36 {
37   PetscInt d;
38   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
39 }
40 
41 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
42                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
43                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
44                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
45 {
46   PetscInt d;
47   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
48 }
49 
50 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51 {
52   *u = PetscSqr(x[0]) + PetscSqr(x[1]);
53   return 0;
54 }
55 
56 static void f0_strong_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
57                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
58                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
59                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
60 {
61   PetscInt d;
62   for (d = 0; d < dim; ++d) f0[0] -= u_x[dim + d*dim+d];
63   f0[0] += 4.0;
64 }
65 
66 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
67 {
68   PetscErrorCode ierr;
69 
70   PetscFunctionBeginUser;
71   options->nit    = 10;
72   options->strong = PETSC_FALSE;
73   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
74   ierr = PetscOptionsInt("-benchmark_it", "Solve the benchmark problem this many times", "ex13.c", options->nit, &options->nit, NULL);CHKERRQ(ierr);
75   ierr = PetscOptionsBool("-strong", "Do not integrate the Laplacian by parts", "ex13.c", options->strong, &options->strong, NULL);CHKERRQ(ierr);
76   ierr = PetscOptionsEnd();
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
81 {
82   PetscErrorCode ierr;
83 
84   PetscFunctionBeginUser;
85   /* Create box mesh */
86   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
87   /* TODO: This should be pulled into the library */
88   {
89     char      convType[256];
90     PetscBool flg;
91 
92     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
93     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
94     ierr = PetscOptionsEnd();
95     if (flg) {
96       DM dmConv;
97 
98       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
99       if (dmConv) {
100         ierr = DMDestroy(dm);CHKERRQ(ierr);
101         *dm  = dmConv;
102       }
103     }
104   }
105   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
106 
107   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
108   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
109   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
110   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
111   PetscFunctionReturn(0);
112 }
113 
114 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
115 {
116   PetscDS        ds;
117   const PetscInt id = 1;
118   PetscErrorCode ierr;
119 
120   PetscFunctionBeginUser;
121   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
122   if (user->strong) {
123     ierr = PetscDSSetResidual(ds, 0, f0_strong_u, NULL);CHKERRQ(ierr);
124     ierr = PetscDSSetExactSolution(ds, 0, quadratic_u, user);CHKERRQ(ierr);
125     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) quadratic_u, NULL, 1, &id, user);CHKERRQ(ierr);
126   } else {
127     ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
128     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
129     ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr);
130     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) trig_u, NULL, 1, &id, user);CHKERRQ(ierr);
131   }
132   PetscFunctionReturn(0);
133 }
134 
135 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
136 {
137   DM             cdm = dm;
138   PetscFE        fe;
139   DMPolytopeType ct;
140   PetscBool      simplex;
141   PetscInt       dim, cStart;
142   char           prefix[PETSC_MAX_PATH_LEN];
143   PetscErrorCode ierr;
144 
145   PetscFunctionBeginUser;
146   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
147   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
148   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
149   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
150   /* Create finite element */
151   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
152   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
153   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
154   /* Set discretization and boundary conditions for each mesh */
155   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
156   ierr = DMCreateDS(dm);CHKERRQ(ierr);
157   ierr = (*setup)(dm, user);CHKERRQ(ierr);
158   while (cdm) {
159     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
160     /* TODO: Check whether the boundary of coarse meshes is marked */
161     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
162   }
163   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 int main(int argc, char **argv)
168 {
169   DM             dm;   /* Problem specification */
170   SNES           snes; /* Nonlinear solver */
171   Vec            u;    /* Solutions */
172   AppCtx         user; /* User-defined work context */
173   PetscErrorCode ierr;
174 
175   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
176   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
177   /* Primal system */
178   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
179   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
180   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
181   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
182   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
183   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
184   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
185   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
186   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
187   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
188   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
189   /* Benchmark system */
190   if (user.nit) {
191 #if defined(PETSC_USE_LOG)
192     PetscLogStage kspstage,pcstage;
193 #endif
194     KSP       ksp;
195     PC        pc;
196     Mat       A,P;
197     Vec       b;
198     PetscInt  i;
199     ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr);
200     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
201     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
202     ierr = SNESGetJacobian(snes, &A, &P, NULL, NULL);CHKERRQ(ierr);
203     ierr = VecSet(u, 0.0);CHKERRQ(ierr);
204     ierr = SNESGetFunction(snes, &b, NULL, NULL);CHKERRQ(ierr);
205     ierr = SNESComputeFunction(snes, u, b);CHKERRQ(ierr);
206     ierr = SNESComputeJacobian(snes, u, A, P);CHKERRQ(ierr);
207     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
208     ierr = PetscLogStageRegister("PCSetUp", &pcstage);CHKERRQ(ierr);
209     ierr = PetscLogStagePush(pcstage);CHKERRQ(ierr);
210     ierr = PCSetUp(pc);CHKERRQ(ierr);
211     ierr = PetscLogStagePop();CHKERRQ(ierr);
212     ierr = PetscLogStageRegister("KSP Solve only", &kspstage);CHKERRQ(ierr);
213     ierr = PetscLogStagePush(kspstage);CHKERRQ(ierr);
214     for (i=0;i<user.nit;i++) {
215       ierr = VecZeroEntries(u);CHKERRQ(ierr);
216       ierr = KSPSolve(ksp, b, u);CHKERRQ(ierr);
217     }
218     ierr = PetscLogStagePop();CHKERRQ(ierr);
219   }
220   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
221   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
222   /* Cleanup */
223   ierr = VecDestroy(&u);CHKERRQ(ierr);
224   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
225   ierr = DMDestroy(&dm);CHKERRQ(ierr);
226   ierr = PetscFinalize();
227   return ierr;
228 }
229 
230 /*TEST
231 
232   test:
233     suffix: strong
234     requires: triangle
235     args: -dm_plex_box_dim 2 -dm_refine 1 -benchmark_it 0 -dmsnes_check \
236           -potential_petscspace_degree 2 -dm_ds_jet_degree 2 -strong
237 
238   test:
239     suffix: bench
240     nsize: 4
241     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2,8 -dm_refine 1 -dm_distribute \
242           -petscpartitioner_type simple -petscpartitioner_simple_process_grid 1,1,2 -petscpartitioner_simple_node_grid 1,1,2 \
243           -potential_petscspace_degree 2 -ksp_type cg -pc_type gamg -benchmark_it 1 -dm_view -snes_rtol 1.e-4
244 
245   test:
246     suffix: comparison
247     nsize: 4
248     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
249       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
250       -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 \
251       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4
252 
253   test:
254     suffix: cuda
255     nsize: 4
256     requires: cuda
257     output_file: output/ex13_comparison.out
258     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
259       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
260       -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 \
261       -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
262 
263   test:
264     suffix: kokkos_comp
265     nsize: 4
266     requires: kokkos_kernels
267     output_file: output/ex13_comparison.out
268     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
269       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
270       -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 \
271       -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
272 
273   test:
274     nsize: 4
275     requires: kokkos_kernels
276     suffix: kokkos
277     args: -dm_plex_box_dim 2 -dm_plex_box_faces 2,8 -dm_distribute -petscpartitioner_type simple -petscpartitioner_simple_process_grid 2,1 \
278           -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 \
279           -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
280 
281   test:
282     suffix: aijmkl_comp
283     nsize: 4
284     requires: mkl
285     output_file: output/ex13_comparison.out
286     args: -dm_plex_box_dim 2 -dm_plex_box_faces 4,4 -dm_refine 3 -petscpartitioner_simple_process_grid 2,2 \
287       -petscpartitioner_simple_node_grid 1,1 -potential_petscspace_degree 2 -dm_distribute -petscpartitioner_type simple \
288       -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 \
289       -pc_gamg_coarse_eq_limit 10 -snes_converged_reason -ksp_converged_reason -snes_rtol 1.e-4 -dm_mat_type aijmkl
290 
291   test:
292     suffix: aijmkl_seq
293     nsize: 1
294     requires: mkl
295     TODO: broken (INDEFINITE PC)
296     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
297 
298 TEST*/
299