xref: /libCEED/examples/petsc/bpssphere.c (revision baf96a30fc83f1cce83f61e183e51df19381d4f1)
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 #if PETSC_VERSION_LT(3, 12, 0)
50 #ifdef PETSC_HAVE_CUDA
51 #include <petsccuda.h>
52 // Note: With PETSc prior to version 3.12.0, providing the source path to include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
53 #endif
54 #endif
55 
56 int main(int argc, char **argv) {
57   MPI_Comm             comm;
58   char                 ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN];
59   double               my_rt_start, my_rt, rt_min, rt_max;
60   PetscInt             degree = 3, q_extra, l_size, g_size, topo_dim = 2, num_comp_x = 3, num_comp_u = 1, xl_size;
61   PetscScalar         *r;
62   PetscBool            test_mode, benchmark_mode, read_mesh, write_solution, simplex;
63   PetscLogStage        solve_stage;
64   Vec                  X, X_loc, rhs, rhs_loc;
65   Mat                  mat_O;
66   KSP                  ksp;
67   DM                   dm;
68   OperatorApplyContext op_apply_ctx, op_error_ctx;
69   Ceed                 ceed;
70   CeedData             ceed_data;
71   CeedQFunction        qf_error;
72   CeedOperator         op_error;
73   CeedVector           rhs_ceed, target;
74   BPType               bp_choice;
75   VecType              vec_type;
76   PetscMemType         mem_type;
77 
78   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
79   comm = PETSC_COMM_WORLD;
80 
81   // Read command line options
82   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
83   bp_choice = CEED_BP1;
84   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
85   num_comp_u = bp_options[bp_choice].num_comp_u;
86   test_mode  = PETSC_FALSE;
87   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
88   benchmark_mode = PETSC_FALSE;
89   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
90   write_solution = PETSC_FALSE;
91   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
92   degree = test_mode ? 3 : 2;
93   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
94   q_extra = bp_options[bp_choice].q_extra;
95   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
96   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
97   read_mesh = PETSC_FALSE;
98   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
99   simplex = PETSC_FALSE;
100   PetscCall(PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", NULL, simplex, &simplex, NULL));
101   PetscOptionsEnd();
102 
103   // Setup DM
104   if (read_mesh) {
105     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm));
106   } else {
107     // Create the mesh as a 0-refined sphere.
108     // This will create a cubic surface, not a box, and will snap to the unit sphere upon refinement.
109     PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm));
110     // Set the object name
111     PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere"));
112     // Refine DMPlex with uniform refinement using runtime option -dm_refine
113     PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
114   }
115   PetscCall(DMSetFromOptions(dm));
116   // View DMPlex via runtime option
117   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
118 
119   // Create DM
120   PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false));
121 
122   // Create vectors
123   PetscCall(DMCreateGlobalVector(dm, &X));
124   PetscCall(VecGetLocalSize(X, &l_size));
125   PetscCall(VecGetSize(X, &g_size));
126   PetscCall(DMCreateLocalVector(dm, &X_loc));
127   PetscCall(VecGetSize(X_loc, &xl_size));
128   PetscCall(VecDuplicate(X, &rhs));
129 
130   // Operator
131   PetscCall(PetscMalloc1(1, &op_apply_ctx));
132   PetscCall(PetscMalloc1(1, &op_error_ctx));
133   PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
134   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
135 
136   // Set up libCEED
137   CeedInit(ceed_resource, &ceed);
138   CeedMemType mem_type_backend;
139   CeedGetPreferredMemType(ceed, &mem_type_backend);
140 
141   PetscCall(DMGetVecType(dm, &vec_type));
142   if (!vec_type) {  // Not yet set by user -dm_vec_type
143     switch (mem_type_backend) {
144       case CEED_MEM_HOST:
145         vec_type = VECSTANDARD;
146         break;
147       case CEED_MEM_DEVICE: {
148         const char *resolved;
149         CeedGetResource(ceed, &resolved);
150         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
151         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
152         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
153         else vec_type = VECSTANDARD;
154       }
155     }
156     PetscCall(DMSetVecType(dm, vec_type));
157   }
158 
159   // Print summary
160   if (!test_mode) {
161     PetscInt    P = degree + 1, Q = P + q_extra;
162     const char *used_resource;
163     CeedGetResource(ceed, &used_resource);
164     PetscCall(PetscPrintf(comm,
165                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " on the Sphere -- libCEED + PETSc --\n"
166                           "  libCEED:\n"
167                           "    libCEED Backend                         : %s\n"
168                           "    libCEED Backend MemType                 : %s\n"
169                           "  Mesh:\n"
170                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
171                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
172                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
173                           "    Global nodes                            : %" PetscInt_FMT "\n",
174                           bp_choice + 1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u));
175   }
176 
177   // Create RHS vector
178   PetscCall(VecDuplicate(X_loc, &rhs_loc));
179   PetscCall(VecZeroEntries(rhs_loc));
180   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
181   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
182   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
183 
184   // Setup libCEED's objects
185   PetscCall(PetscMalloc1(1, &ceed_data));
186   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,
187                                  rhs_ceed, &target));
188 
189   // Gather RHS
190   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
191   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
192   PetscCall(VecZeroEntries(rhs));
193   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
194   CeedVectorDestroy(&rhs_ceed);
195 
196   // Create the error Q-function
197   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
198   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
199   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
200   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
201   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
202 
203   // Create the error operator
204   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
205   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
206   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target);
207   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
208   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
209 
210   // Set up apply operator context
211   PetscCall(SetupApplyOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
212 
213   // Setup solver
214   PetscCall(KSPCreate(comm, &ksp));
215   {
216     PC pc;
217     PetscCall(KSPGetPC(ksp, &pc));
218     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
219       PetscCall(PCSetType(pc, PCJACOBI));
220       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
221     } else {
222       PetscCall(PCSetType(pc, PCNONE));
223       MatNullSpace nullspace;
224 
225       PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
226       PetscCall(MatSetNullSpace(mat_O, nullspace));
227       PetscCall(MatNullSpaceDestroy(&nullspace));
228     }
229     PetscCall(KSPSetType(ksp, KSPCG));
230     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
231     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
232   }
233   PetscCall(KSPSetFromOptions(ksp));
234   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
235 
236   // First run, if benchmarking
237   if (benchmark_mode) {
238     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
239     my_rt_start = MPI_Wtime();
240     PetscCall(KSPSolve(ksp, rhs, X));
241     my_rt = MPI_Wtime() - my_rt_start;
242     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
243     // Set maxits based on first iteration timing
244     if (my_rt > 0.02) {
245       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
246     } else {
247       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
248     }
249   }
250 
251   // Timed solve
252   PetscCall(VecZeroEntries(X));
253   PetscCall(PetscBarrier((PetscObject)ksp));
254 
255   // -- Performance logging
256   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
257   PetscCall(PetscLogStagePush(solve_stage));
258 
259   // -- Solve
260   my_rt_start = MPI_Wtime();
261   PetscCall(KSPSolve(ksp, rhs, X));
262   my_rt = MPI_Wtime() - my_rt_start;
263 
264   // -- Performance logging
265   PetscCall(PetscLogStagePop());
266 
267   // Output results
268   {
269     KSPType            ksp_type;
270     KSPConvergedReason reason;
271     PetscReal          rnorm;
272     PetscInt           its;
273     PetscCall(KSPGetType(ksp, &ksp_type));
274     PetscCall(KSPGetConvergedReason(ksp, &reason));
275     PetscCall(KSPGetIterationNumber(ksp, &its));
276     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
277     if (!test_mode || reason < 0 || rnorm > 1e-8) {
278       PetscCall(PetscPrintf(comm,
279                             "  KSP:\n"
280                             "    KSP Type                                : %s\n"
281                             "    KSP Convergence                         : %s\n"
282                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
283                             "    Final rnorm                             : %e\n",
284                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
285     }
286     if (!test_mode) {
287       PetscCall(PetscPrintf(comm, "  Performance:\n"));
288     }
289     {
290       // Set up error operator context
291       PetscCall(SetupErrorOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
292       PetscScalar l2_error;
293       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
294       PetscReal tol = 5e-4;
295       if (!test_mode || l2_error > tol) {
296         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
297         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
298         PetscCall(PetscPrintf(comm,
299                               "    L2 Error                                : %e\n"
300                               "    CG Solve Time                           : %g (%g) sec\n",
301                               (double)l2_error, rt_max, rt_min));
302       }
303     }
304     if (benchmark_mode && (!test_mode)) {
305       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
306                             1e-6 * g_size * its / rt_min));
307     }
308   }
309 
310   // Output solution
311   if (write_solution) {
312     PetscViewer vtk_viewer_soln;
313 
314     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
315     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
316     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
317     PetscCall(VecView(X, vtk_viewer_soln));
318     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
319   }
320 
321   // Cleanup
322   PetscCall(VecDestroy(&X));
323   PetscCall(VecDestroy(&X_loc));
324   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
325   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
326   PetscCall(MatDestroy(&mat_O));
327   PetscCall(PetscFree(op_apply_ctx));
328   PetscCall(PetscFree(op_error_ctx));
329   PetscCall(CeedDataDestroy(0, ceed_data));
330   PetscCall(DMDestroy(&dm));
331 
332   PetscCall(VecDestroy(&rhs));
333   PetscCall(VecDestroy(&rhs_loc));
334   PetscCall(KSPDestroy(&ksp));
335   CeedVectorDestroy(&target);
336   CeedQFunctionDestroy(&qf_error);
337   CeedOperatorDestroy(&op_error);
338   CeedDestroy(&ceed);
339   return PetscFinalize();
340 }
341