1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED: http://github.com/ceed
7
8 // libCEED + PETSc Example: CEED BPs
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps, on
11 // a closed surface, such as the one of a discrete sphere.
12 //
13 // The code uses higher level communication protocols in DMPlex.
14 //
15 // Build with:
16 //
17 // make bpssphere [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18 //
19 // Sample runs:
20 //
21 // bpssphere -problem bp1 -degree 3
22 // bpssphere -problem bp2 -degree 3
23 // bpssphere -problem bp3 -degree 3
24 // bpssphere -problem bp4 -degree 3
25 // bpssphere -problem bp5 -degree 3 -ceed /cpu/self
26 // bpssphere -problem bp6 -degree 3 -ceed /gpu/cuda
27 //
28 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -dm_refine 2
29
30 /// @file
31 /// CEED BPs example using PETSc with DMPlex
32 /// See bps.c for a "raw" implementation using a structured grid and bpsdmplex.c for an implementation using an unstructured grid.
33 static const char help[] = "Solve CEED BPs on a sphere using DMPlex in PETSc\n";
34
35 #include "bpssphere.h"
36
37 #include <ceed.h>
38 #include <petscdmplex.h>
39 #include <petscksp.h>
40 #include <stdbool.h>
41 #include <string.h>
42
43 #include "include/libceedsetup.h"
44 #include "include/matops.h"
45 #include "include/petscutils.h"
46 #include "include/petscversion.h"
47 #include "include/sphereproblemdata.h"
48
main(int argc,char ** argv)49 int main(int argc, char **argv) {
50 MPI_Comm comm;
51 char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN];
52 double my_rt_start, my_rt, rt_min, rt_max;
53 PetscInt degree = 3, q_extra, l_size, g_size, topo_dim = 2, num_comp_x = 3, num_comp_u = 1, xl_size;
54 PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex;
55 PetscLogStage solve_stage;
56 Vec X, X_loc, rhs, rhs_loc;
57 Mat mat_O;
58 KSP ksp;
59 DM dm;
60 OperatorApplyContext op_apply_ctx, op_error_ctx;
61 Ceed ceed;
62 CeedData ceed_data;
63 CeedQFunction qf_error;
64 CeedOperator op_error;
65 CeedVector rhs_ceed, target;
66 BPType bp_choice;
67 VecType vec_type = VECSTANDARD;
68 PetscMemType mem_type;
69
70 PetscCall(PetscInitialize(&argc, &argv, NULL, help));
71 comm = PETSC_COMM_WORLD;
72
73 // Read command line options
74 PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
75 bp_choice = CEED_BP1;
76 PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
77 num_comp_u = bp_options[bp_choice].num_comp_u;
78 test_mode = PETSC_FALSE;
79 PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
80 benchmark_mode = PETSC_FALSE;
81 PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
82 write_solution = PETSC_FALSE;
83 PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
84 degree = test_mode ? 3 : 2;
85 PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL));
86 q_extra = bp_options[bp_choice].q_extra;
87 PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
88 PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
89 read_mesh = PETSC_FALSE;
90 PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
91 simplex = PETSC_FALSE;
92 PetscCall(PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", NULL, simplex, &simplex, NULL));
93 PetscOptionsEnd();
94
95 // Set up libCEED
96 CeedInit(ceed_resource, &ceed);
97 CeedMemType mem_type_backend;
98 CeedGetPreferredMemType(ceed, &mem_type_backend);
99
100 // Set mesh vec_type
101 switch (mem_type_backend) {
102 case CEED_MEM_HOST:
103 vec_type = VECSTANDARD;
104 break;
105 case CEED_MEM_DEVICE: {
106 const char *resolved;
107
108 CeedGetResource(ceed, &resolved);
109 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
110 else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
111 else vec_type = VECSTANDARD;
112 }
113 }
114
115 // Setup DM
116 if (read_mesh) {
117 PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm));
118 } else {
119 // Create the mesh as a 0-refined sphere.
120 // This will create a cubic surface, not a box, and will snap to the unit sphere upon refinement.
121 PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm));
122 // Set the object name
123 PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere"));
124 // Refine DMPlex with uniform refinement using runtime option -dm_refine
125 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
126 }
127 PetscCall(DMSetVecType(dm, vec_type));
128 PetscCall(DMSetFromOptions(dm));
129 // View DMPlex via runtime option
130 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
131
132 // Create DM
133 PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false));
134
135 // Create vectors
136 PetscCall(DMCreateGlobalVector(dm, &X));
137 PetscCall(VecGetLocalSize(X, &l_size));
138 PetscCall(VecGetSize(X, &g_size));
139 PetscCall(DMCreateLocalVector(dm, &X_loc));
140 PetscCall(VecGetSize(X_loc, &xl_size));
141 PetscCall(VecDuplicate(X, &rhs));
142
143 // Operator
144 PetscCall(PetscMalloc1(1, &op_apply_ctx));
145 PetscCall(PetscMalloc1(1, &op_error_ctx));
146 PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
147 PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
148
149 // Print summary
150 if (!test_mode) {
151 PetscInt P = degree + 1, Q = P + q_extra;
152 const char *used_resource;
153 CeedGetResource(ceed, &used_resource);
154 PetscCall(PetscPrintf(comm,
155 "\n-- CEED Benchmark Problem %" CeedInt_FMT " on the Sphere -- libCEED + PETSc --\n"
156 " libCEED:\n"
157 " libCEED Backend : %s\n"
158 " libCEED Backend MemType : %s\n"
159 " Mesh:\n"
160 " Solution Order (P) : %" PetscInt_FMT "\n"
161 " Quadrature Order (Q) : %" PetscInt_FMT "\n"
162 " Additional quadrature points (q_extra) : %" PetscInt_FMT "\n"
163 " Global nodes : %" PetscInt_FMT "\n",
164 bp_choice + 1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u));
165 }
166
167 // Create RHS vector
168 PetscCall(VecDuplicate(X_loc, &rhs_loc));
169 PetscCall(VecZeroEntries(rhs_loc));
170 CeedVectorCreate(ceed, xl_size, &rhs_ceed);
171 PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
172
173 // Setup libCEED's objects
174 PetscCall(PetscMalloc1(1, &ceed_data));
175 PetscCall(SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x, num_comp_u, g_size, xl_size, bp_options[bp_choice], ceed_data, true,
176 true, rhs_ceed, &target));
177
178 // Gather RHS
179 PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
180 PetscCall(VecZeroEntries(rhs));
181 PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
182 CeedVectorDestroy(&rhs_ceed);
183
184 // Create the error Q-function
185 CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
186 CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
187 CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
188 CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
189 CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
190
191 // Create the error operator
192 CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
193 CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
194 CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_NONE, target);
195 CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data);
196 CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
197
198 // Set up apply operator context
199 PetscCall(SetupApplyOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
200
201 // Setup solver
202 PetscCall(KSPCreate(comm, &ksp));
203 {
204 PC pc;
205 PetscCall(KSPGetPC(ksp, &pc));
206 if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
207 PetscCall(PCSetType(pc, PCJACOBI));
208 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
209 } else {
210 PetscCall(PCSetType(pc, PCNONE));
211 MatNullSpace nullspace;
212
213 PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
214 PetscCall(MatSetNullSpace(mat_O, nullspace));
215 PetscCall(MatNullSpaceDestroy(&nullspace));
216 }
217 PetscCall(KSPSetType(ksp, KSPCG));
218 PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
219 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
220 }
221 PetscCall(KSPSetFromOptions(ksp));
222 PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
223
224 // First run, if benchmarking
225 if (benchmark_mode) {
226 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
227 my_rt_start = MPI_Wtime();
228 PetscCall(KSPSolve(ksp, rhs, X));
229 my_rt = MPI_Wtime() - my_rt_start;
230 PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
231 // Set maxits based on first iteration timing
232 if (my_rt > 0.02) {
233 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
234 } else {
235 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
236 }
237 }
238
239 // Timed solve
240 PetscCall(VecZeroEntries(X));
241 PetscCall(PetscBarrier((PetscObject)ksp));
242
243 // -- Performance logging
244 PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
245 PetscCall(PetscLogStagePush(solve_stage));
246
247 // -- Solve
248 my_rt_start = MPI_Wtime();
249 PetscCall(KSPSolve(ksp, rhs, X));
250 my_rt = MPI_Wtime() - my_rt_start;
251
252 // -- Performance logging
253 PetscCall(PetscLogStagePop());
254
255 // Output results
256 {
257 KSPType ksp_type;
258 KSPConvergedReason reason;
259 PetscReal rnorm;
260 PetscInt its;
261 PetscCall(KSPGetType(ksp, &ksp_type));
262 PetscCall(KSPGetConvergedReason(ksp, &reason));
263 PetscCall(KSPGetIterationNumber(ksp, &its));
264 PetscCall(KSPGetResidualNorm(ksp, &rnorm));
265 if (!test_mode || reason < 0 || rnorm > 1e-8) {
266 PetscCall(PetscPrintf(comm,
267 " KSP:\n"
268 " KSP Type : %s\n"
269 " KSP Convergence : %s\n"
270 " Total KSP Iterations : %" PetscInt_FMT "\n"
271 " Final rnorm : %e\n",
272 ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
273 }
274 if (!test_mode) {
275 PetscCall(PetscPrintf(comm, " Performance:\n"));
276 }
277 {
278 // Set up error operator context
279 PetscCall(SetupErrorOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
280 PetscScalar l2_error;
281 PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
282 PetscReal tol = 5e-4;
283 if (!test_mode || l2_error > tol) {
284 PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
285 PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
286 PetscCall(PetscPrintf(comm,
287 " L2 Error : %e\n"
288 " CG Solve Time : %g (%g) sec\n",
289 (double)l2_error, rt_max, rt_min));
290 }
291 }
292 if (benchmark_mode && (!test_mode)) {
293 PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
294 1e-6 * g_size * its / rt_min));
295 }
296 }
297
298 // Output solution
299 if (write_solution) {
300 PetscViewer vtk_viewer_soln;
301
302 PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
303 PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
304 PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
305 PetscCall(VecView(X, vtk_viewer_soln));
306 PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
307 }
308
309 // Cleanup
310 PetscCall(VecDestroy(&X));
311 PetscCall(VecDestroy(&X_loc));
312 PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
313 PetscCall(VecDestroy(&op_error_ctx->Y_loc));
314 PetscCall(MatDestroy(&mat_O));
315 PetscCall(PetscFree(op_apply_ctx));
316 PetscCall(PetscFree(op_error_ctx));
317 PetscCall(CeedDataDestroy(0, ceed_data));
318 PetscCall(DMDestroy(&dm));
319
320 PetscCall(VecDestroy(&rhs));
321 PetscCall(VecDestroy(&rhs_loc));
322 PetscCall(KSPDestroy(&ksp));
323 CeedVectorDestroy(&target);
324 CeedQFunctionDestroy(&qf_error);
325 CeedOperatorDestroy(&op_error);
326 CeedDestroy(&ceed);
327 return PetscFinalize();
328 }
329