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