xref: /libCEED/examples/petsc/bpssphere.c (revision 99421279ebf997f2fb157e1e18f23b9dc112bbac)
1 // Copyright (c) 2017-2024, 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 
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, &degree, 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/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
111       else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
112       else vec_type = VECSTANDARD;
113     }
114   }
115 
116   // Setup DM
117   if (read_mesh) {
118     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm));
119   } else {
120     // Create the mesh as a 0-refined sphere.
121     // This will create a cubic surface, not a box, and will snap to the unit sphere upon refinement.
122     PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm));
123     // Set the object name
124     PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere"));
125     // Refine DMPlex with uniform refinement using runtime option -dm_refine
126     PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
127   }
128   PetscCall(DMSetVecType(dm, vec_type));
129   PetscCall(DMSetFromOptions(dm));
130   // View DMPlex via runtime option
131   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
132 
133   // Create DM
134   PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false));
135 
136   // Create vectors
137   PetscCall(DMCreateGlobalVector(dm, &X));
138   PetscCall(VecGetLocalSize(X, &l_size));
139   PetscCall(VecGetSize(X, &g_size));
140   PetscCall(DMCreateLocalVector(dm, &X_loc));
141   PetscCall(VecGetSize(X_loc, &xl_size));
142   PetscCall(VecDuplicate(X, &rhs));
143 
144   // Operator
145   PetscCall(PetscMalloc1(1, &op_apply_ctx));
146   PetscCall(PetscMalloc1(1, &op_error_ctx));
147   PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
148   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
149 
150   // Print summary
151   if (!test_mode) {
152     PetscInt    P = degree + 1, Q = P + q_extra;
153     const char *used_resource;
154     CeedGetResource(ceed, &used_resource);
155     PetscCall(PetscPrintf(comm,
156                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " on the Sphere -- libCEED + PETSc --\n"
157                           "  libCEED:\n"
158                           "    libCEED Backend                         : %s\n"
159                           "    libCEED Backend MemType                 : %s\n"
160                           "  Mesh:\n"
161                           "    Solution Order (P)                      : %" PetscInt_FMT "\n"
162                           "    Quadrature  Order (Q)                   : %" PetscInt_FMT "\n"
163                           "    Additional quadrature points (q_extra)  : %" PetscInt_FMT "\n"
164                           "    Global nodes                            : %" PetscInt_FMT "\n",
165                           bp_choice + 1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u));
166   }
167 
168   // Create RHS vector
169   PetscCall(VecDuplicate(X_loc, &rhs_loc));
170   PetscCall(VecZeroEntries(rhs_loc));
171   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
172   PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
173 
174   // Setup libCEED's objects
175   PetscCall(PetscMalloc1(1, &ceed_data));
176   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,
177                                  true, rhs_ceed, &target));
178 
179   // Gather RHS
180   PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
181   PetscCall(VecZeroEntries(rhs));
182   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
183   CeedVectorDestroy(&rhs_ceed);
184 
185   // Create the error Q-function
186   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
187   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
188   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
189   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
190   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
191 
192   // Create the error operator
193   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
194   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
195   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_NONE, target);
196   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data);
197   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
198 
199   // Set up apply operator context
200   PetscCall(SetupApplyOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
201 
202   // Setup solver
203   PetscCall(KSPCreate(comm, &ksp));
204   {
205     PC pc;
206     PetscCall(KSPGetPC(ksp, &pc));
207     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
208       PetscCall(PCSetType(pc, PCJACOBI));
209       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
210     } else {
211       PetscCall(PCSetType(pc, PCNONE));
212       MatNullSpace nullspace;
213 
214       PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
215       PetscCall(MatSetNullSpace(mat_O, nullspace));
216       PetscCall(MatNullSpaceDestroy(&nullspace));
217     }
218     PetscCall(KSPSetType(ksp, KSPCG));
219     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
220     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
221   }
222   PetscCall(KSPSetFromOptions(ksp));
223   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
224 
225   // First run, if benchmarking
226   if (benchmark_mode) {
227     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
228     my_rt_start = MPI_Wtime();
229     PetscCall(KSPSolve(ksp, rhs, X));
230     my_rt = MPI_Wtime() - my_rt_start;
231     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
232     // Set maxits based on first iteration timing
233     if (my_rt > 0.02) {
234       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
235     } else {
236       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
237     }
238   }
239 
240   // Timed solve
241   PetscCall(VecZeroEntries(X));
242   PetscCall(PetscBarrier((PetscObject)ksp));
243 
244   // -- Performance logging
245   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
246   PetscCall(PetscLogStagePush(solve_stage));
247 
248   // -- Solve
249   my_rt_start = MPI_Wtime();
250   PetscCall(KSPSolve(ksp, rhs, X));
251   my_rt = MPI_Wtime() - my_rt_start;
252 
253   // -- Performance logging
254   PetscCall(PetscLogStagePop());
255 
256   // Output results
257   {
258     KSPType            ksp_type;
259     KSPConvergedReason reason;
260     PetscReal          rnorm;
261     PetscInt           its;
262     PetscCall(KSPGetType(ksp, &ksp_type));
263     PetscCall(KSPGetConvergedReason(ksp, &reason));
264     PetscCall(KSPGetIterationNumber(ksp, &its));
265     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
266     if (!test_mode || reason < 0 || rnorm > 1e-8) {
267       PetscCall(PetscPrintf(comm,
268                             "  KSP:\n"
269                             "    KSP Type                                : %s\n"
270                             "    KSP Convergence                         : %s\n"
271                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
272                             "    Final rnorm                             : %e\n",
273                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
274     }
275     if (!test_mode) {
276       PetscCall(PetscPrintf(comm, "  Performance:\n"));
277     }
278     {
279       // Set up error operator context
280       PetscCall(SetupErrorOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
281       PetscScalar l2_error;
282       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
283       PetscReal tol = 5e-4;
284       if (!test_mode || l2_error > tol) {
285         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
286         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
287         PetscCall(PetscPrintf(comm,
288                               "    L2 Error                                : %e\n"
289                               "    CG Solve Time                           : %g (%g) sec\n",
290                               (double)l2_error, rt_max, rt_min));
291       }
292     }
293     if (benchmark_mode && (!test_mode)) {
294       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
295                             1e-6 * g_size * its / rt_min));
296     }
297   }
298 
299   // Output solution
300   if (write_solution) {
301     PetscViewer vtk_viewer_soln;
302 
303     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
304     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
305     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
306     PetscCall(VecView(X, vtk_viewer_soln));
307     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
308   }
309 
310   // Cleanup
311   PetscCall(VecDestroy(&X));
312   PetscCall(VecDestroy(&X_loc));
313   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
314   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
315   PetscCall(MatDestroy(&mat_O));
316   PetscCall(PetscFree(op_apply_ctx));
317   PetscCall(PetscFree(op_error_ctx));
318   PetscCall(CeedDataDestroy(0, ceed_data));
319   PetscCall(DMDestroy(&dm));
320 
321   PetscCall(VecDestroy(&rhs));
322   PetscCall(VecDestroy(&rhs_loc));
323   PetscCall(KSPDestroy(&ksp));
324   CeedVectorDestroy(&target);
325   CeedQFunctionDestroy(&qf_error);
326   CeedOperatorDestroy(&op_error);
327   CeedDestroy(&ceed);
328   return PetscFinalize();
329 }
330