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