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