xref: /libCEED/examples/petsc/bpsswarm.c (revision 78f7fce354c760f980d0580f63d29ea51c63cedc)
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 particle swarm.
12 //
13 // The code uses higher level communication protocols in DMPlex and DMSwarm.
14 //
15 // Build with:
16 //
17 //     make bpsswarm [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 //
25 //TESTARGS -ceed {ceed_resource} -test -problem bp3 -dm_plex_dim 3 -dm_plex_box_faces 10,10,10 -dm_plex_simplex 0 -swarm uniform -points_per_cell 500
26 
27 /// @file
28 /// CEED BPs example using PETSc with DMPlex
29 /// See bpsraw.c for a "raw" implementation using a structured grid and bps.c for an implementation using an unstructured grid.
30 static const char help[]              = "Solve CEED BPs on a particle swarm using DMPlex and DMSwarm in PETSc\n";
31 const char        DMSwarmPICField_u[] = "u";
32 
33 #include "bps.h"
34 
35 #include <ceed.h>
36 #include <petscdmplex.h>
37 #include <petscksp.h>
38 #include <stdbool.h>
39 #include <string.h>
40 
41 #include "include/bpsproblemdata.h"
42 #include "include/libceedsetup.h"
43 #include "include/matops.h"
44 #include "include/petscutils.h"
45 #include "include/petscversion.h"
46 #include "include/swarmutils.h"
47 
48 int main(int argc, char **argv) {
49   MPI_Comm             comm;
50   char                 ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN];
51   double               my_rt_start, my_rt, rt_min, rt_max;
52   PetscInt             comm_size, degree = 3, q_extra, l_size, g_size, dim = 3, num_comp_u = 1, xl_size, num_points = 1728, num_points_per_cell = 64;
53   PetscBool            test_mode, benchmark_mode, read_mesh, write_solution, write_true_solution_swarm;
54   PetscLogStage        solve_stage;
55   Vec                  X, X_loc, rhs;
56   Mat                  mat_O;
57   KSP                  ksp;
58   DM                   dm_mesh, dm_swarm;
59   OperatorApplyContext op_apply_ctx, op_error_ctx;
60   Ceed                 ceed;
61   CeedData             ceed_data;
62   CeedOperator         op_error;
63   BPType               bp_choice;
64   VecType              vec_type;
65   PointSwarmType       point_swarm_type = SWARM_GAUSS;
66   PetscMPIInt          ranks_per_node;
67   char                 hostname[PETSC_MAX_PATH_LEN];
68 
69   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
70   comm = PETSC_COMM_WORLD;
71   PetscCall(MPI_Comm_size(comm, &comm_size));
72 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
73   {
74     MPI_Comm splitcomm;
75     PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm));
76     PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node));
77     PetscCall(MPI_Comm_free(&splitcomm));
78   }
79 #else
80   ranks_per_node = -1;  // Unknown
81 #endif
82 
83   // Read command line options
84   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
85   bp_choice = CEED_BP1;
86   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
87   num_comp_u = bp_options[bp_choice].num_comp_u;
88   test_mode  = PETSC_FALSE;
89   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
90   benchmark_mode = PETSC_FALSE;
91   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
92   write_solution = PETSC_FALSE;
93   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
94   write_true_solution_swarm = PETSC_FALSE;
95   PetscCall(PetscOptionsBool("-write_true_solution_swarm", "Write true solution at swarm points for visualization", NULL, write_true_solution_swarm,
96                              &write_true_solution_swarm, NULL));
97   degree = test_mode ? 3 : 2;
98   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
99   q_extra = bp_options[bp_choice].q_extra;
100   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
101   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
102   PetscCall(PetscGetHostName(hostname, sizeof hostname));
103   PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL));
104   read_mesh = PETSC_FALSE;
105   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
106   PetscCall(PetscOptionsEnum("-swarm", "Swarm points distribution", NULL, point_swarm_types, (PetscEnum)point_swarm_type,
107                              (PetscEnum *)&point_swarm_type, NULL));
108   {
109     PetscBool user_set_num_points_per_cell = PETSC_FALSE;
110     PetscInt  num_cells_total = 1, tmp = dim;
111     PetscInt  num_cells[] = {1, 1, 1};
112 
113     PetscCall(PetscOptionsInt("-points_per_cell", "Total number of swarm points in each cell", NULL, num_points_per_cell, &num_points_per_cell,
114                               &user_set_num_points_per_cell));
115     PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL));
116     PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &tmp, NULL));
117 
118     PetscCheck(tmp == dim, comm, PETSC_ERR_USER, "Number of values for -dm_plex_box_faces must match dimension");
119 
120     num_cells_total = num_cells[0] * num_cells[1] * num_cells[2];
121     PetscCheck(!user_set_num_points_per_cell || point_swarm_type != SWARM_SINUSOIDAL, comm, PETSC_ERR_USER,
122                "Cannot specify points per cell with sinusoidal points locations");
123     if (!user_set_num_points_per_cell) {
124       PetscCall(PetscOptionsInt("-points", "Total number of swarm points", NULL, num_points, &num_points, NULL));
125       num_points_per_cell = PetscCeilInt(num_points, num_cells_total);
126     }
127     if (point_swarm_type != SWARM_SINUSOIDAL) {
128       PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0));
129 
130       num_points_per_cell = 1;
131       for (PetscInt i = 0; i < dim; i++) num_points_per_cell *= num_points_per_cell_1d;
132     }
133     num_points = num_points_per_cell * num_cells_total;
134   }
135   {
136     PetscBool flg;
137     PetscInt  p = ranks_per_node;
138     PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg));
139     if (flg) ranks_per_node = p;
140   }
141   PetscOptionsEnd();
142 
143   // Setup DM
144   if (read_mesh) {
145     PetscCall(DMPlexCreateFromFile(comm, filename, NULL, PETSC_TRUE, &dm_mesh));
146   } else {
147     PetscCall(DMCreate(comm, &dm_mesh));
148     PetscCall(DMSetType(dm_mesh, DMPLEX));
149     PetscCall(DMSetFromOptions(dm_mesh));
150 
151     // -- Check for tensor product mesh
152     {
153       PetscBool is_simplex;
154 
155       PetscCall(DMPlexIsSimplex(dm_mesh, &is_simplex));
156       PetscCheck(!is_simplex, comm, PETSC_ERR_USER, "Only tensor-product background meshes supported");
157     }
158   }
159   PetscCall(DMGetDimension(dm_mesh, &dim));
160   PetscCall(SetupDMByDegree(dm_mesh, degree, q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc));
161 
162   // View mesh
163   PetscCall(DMSetOptionsPrefix(dm_mesh, "final_"));
164   PetscCall(DMViewFromOptions(dm_mesh, NULL, "-dm_view"));
165 
166   // Create particle swarm
167   PetscCall(DMCreate(comm, &dm_swarm));
168   PetscCall(DMSetType(dm_swarm, DMSWARM));
169   PetscCall(DMSetDimension(dm_swarm, dim));
170   PetscCall(DMSwarmSetType(dm_swarm, DMSWARM_PIC));
171   PetscCall(DMSwarmSetCellDM(dm_swarm, dm_mesh));
172 
173   // -- Swarm field
174   PetscCall(DMSwarmRegisterPetscDatatypeField(dm_swarm, DMSwarmPICField_u, num_comp_u, PETSC_SCALAR));
175   PetscCall(DMSwarmFinalizeFieldRegister(dm_swarm));
176   PetscCall(DMSwarmSetLocalSizes(dm_swarm, num_points, 0));
177   PetscCall(DMSetFromOptions(dm_swarm));
178 
179   // -- Set swarm point locations
180   PetscCall(DMSwarmInitalizePointLocations(dm_swarm, point_swarm_type, num_points, num_points_per_cell));
181   PetscCall(DMSwarmVectorDefineField(dm_swarm, DMSwarmPICField_u));
182 
183   // -- Final particle swarm
184   PetscCall(PetscObjectSetName((PetscObject)dm_swarm, "Particle Swarm"));
185   PetscCall(DMViewFromOptions(dm_swarm, NULL, "-dm_swarm_view"));
186 
187   // Create vectors
188   PetscCall(DMCreateGlobalVector(dm_mesh, &X));
189   PetscCall(VecGetLocalSize(X, &l_size));
190   PetscCall(VecGetSize(X, &g_size));
191   PetscCall(DMCreateLocalVector(dm_mesh, &X_loc));
192   PetscCall(VecGetSize(X_loc, &xl_size));
193   PetscCall(VecDuplicate(X, &rhs));
194 
195   // Operator
196   PetscCall(PetscMalloc1(1, &op_apply_ctx));
197   PetscCall(PetscMalloc1(1, &op_error_ctx));
198   PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
199   PetscCall(MatSetDM(mat_O, dm_mesh));
200   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
201 
202   // Set up libCEED
203   CeedInit(ceed_resource, &ceed);
204   CeedMemType mem_type_backend;
205   CeedGetPreferredMemType(ceed, &mem_type_backend);
206 
207   PetscCall(DMGetVecType(dm_mesh, &vec_type));
208   if (!vec_type) {  // Not yet set by user -dm_vec_type
209     switch (mem_type_backend) {
210       case CEED_MEM_HOST:
211         vec_type = VECSTANDARD;
212         break;
213       case CEED_MEM_DEVICE: {
214         const char *resolved;
215         CeedGetResource(ceed, &resolved);
216         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
217         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
218         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
219         else vec_type = VECSTANDARD;
220       }
221     }
222     PetscCall(DMSetVecType(dm_mesh, vec_type));
223   }
224 
225   // Print summary
226   if (!test_mode) {
227     PetscInt P = degree + 1, Q = P + q_extra;
228 
229     const char *used_resource;
230     CeedGetResource(ceed, &used_resource);
231 
232     VecType vec_type;
233     PetscCall(VecGetType(X, &vec_type));
234 
235     PetscInt c_start, c_end, num_cells_local;
236     PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &c_start, &c_end));
237     num_cells_local = c_end - c_start;
238     DMPolytopeType cell_type;
239     PetscCall(DMPlexGetCellType(dm_mesh, c_start, &cell_type));
240     PetscMPIInt comm_size;
241     PetscCall(MPI_Comm_size(comm, &comm_size));
242 
243     PetscInt num_points_local;
244     PetscCall(DMSwarmGetLocalSize(dm_swarm, &num_points_local));
245 
246     PetscCall(PetscPrintf(comm,
247                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
248                           "  MPI:\n"
249                           "    Hostname                                : %s\n"
250                           "    Total ranks                             : %d\n"
251                           "    Ranks per compute node                  : %d\n"
252                           "  PETSc:\n"
253                           "    PETSc Vec Type                          : %s\n"
254                           "  libCEED:\n"
255                           "    libCEED Backend                         : %s\n"
256                           "    libCEED Backend MemType                 : %s\n"
257                           "  Mesh:\n"
258                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
259                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
260                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
261                           "    Global nodes                            : %" PetscInt_FMT "\n"
262                           "    Local Elements                          : %" PetscInt_FMT "\n"
263                           "    Owned nodes                             : %" PetscInt_FMT "\n"
264                           "    DoF per node                            : %" PetscInt_FMT "\n"
265                           "  Swarm:\n"
266                           "    Global points                           : %" PetscInt_FMT "\n"
267                           "    Local points                            : %" PetscInt_FMT "\n"
268                           "    Avg points per cell                     : %" PetscInt_FMT "\n"
269                           "    Point distribution                      : %s\n",
270                           bp_choice + 1, hostname, comm_size, ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra,
271                           g_size / num_comp_u, num_cells_local, l_size / num_comp_u, num_comp_u, num_points, num_points_local,
272                           num_cells_local > 0 ? num_points_local / num_cells_local : 0, point_swarm_types[point_swarm_type]));
273   }
274 
275   // Setup libCEED's objects
276   Vec target;
277 
278   PetscCall(DMCreateLocalVector(dm_swarm, &target));
279   PetscCall(PetscMalloc1(1, &ceed_data));
280   PetscCall(SetupProblemSwarm(dm_swarm, ceed, bp_options[bp_choice], ceed_data, true, rhs, target));
281   PetscCall(SetupErrorOperator(dm_mesh, ceed, bp_options[bp_choice], dim, dim, num_comp_u, &op_error));
282 
283   // Set up apply operator context
284   PetscCall(SetupApplyOperatorCtx(comm, dm_mesh, ceed, ceed_data, X_loc, op_apply_ctx));
285 
286   // Setup solver
287   PetscCall(KSPCreate(comm, &ksp));
288   {
289     PC pc;
290     PetscCall(KSPGetPC(ksp, &pc));
291     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
292       PetscCall(PCSetType(pc, PCJACOBI));
293       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
294     } else {
295       PetscCall(PCSetType(pc, PCNONE));
296       MatNullSpace nullspace;
297 
298       PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
299       PetscCall(MatSetNullSpace(mat_O, nullspace));
300       PetscCall(MatNullSpaceDestroy(&nullspace));
301     }
302     PetscCall(KSPSetType(ksp, KSPCG));
303     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
304     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
305   }
306   PetscCall(KSPSetFromOptions(ksp));
307   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
308 
309   // First run, if benchmarking
310   if (benchmark_mode) {
311     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
312     my_rt_start = MPI_Wtime();
313     PetscCall(KSPSolve(ksp, rhs, X));
314     my_rt = MPI_Wtime() - my_rt_start;
315     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
316     // Set maxits based on first iteration timing
317     if (my_rt > 0.02) {
318       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
319     } else {
320       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
321     }
322   }
323 
324   // Timed solve
325   PetscCall(VecZeroEntries(X));
326   PetscCall(PetscBarrier((PetscObject)ksp));
327 
328   // -- Performance logging
329   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
330   PetscCall(PetscLogStagePush(solve_stage));
331 
332   // -- Solve
333   my_rt_start = MPI_Wtime();
334   PetscCall(KSPSolve(ksp, rhs, X));
335   my_rt = MPI_Wtime() - my_rt_start;
336 
337   // -- Performance logging
338   PetscCall(PetscLogStagePop());
339 
340   // Output results
341   {
342     KSPType            ksp_type;
343     KSPConvergedReason reason;
344     PetscReal          rnorm;
345     PetscInt           its;
346     PetscCall(KSPGetType(ksp, &ksp_type));
347     PetscCall(KSPGetConvergedReason(ksp, &reason));
348     PetscCall(KSPGetIterationNumber(ksp, &its));
349     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
350     if (!test_mode || reason < 0 || rnorm > 1e-8) {
351       PetscCall(PetscPrintf(comm,
352                             "  KSP:\n"
353                             "    KSP Type                                : %s\n"
354                             "    KSP Convergence                         : %s\n"
355                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
356                             "    Final rnorm                             : %e\n",
357                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
358     }
359     if (!test_mode) {
360       PetscCall(PetscPrintf(comm, "  Performance:\n"));
361     }
362 
363     // View true solution at particles
364     if (write_true_solution_swarm) {
365       Vec u_swarm, u_swarm_old;
366       PetscCall(DMSwarmSortGetAccess(dm_swarm));
367       PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
368       PetscCall(VecDuplicate(u_swarm, &u_swarm_old));
369       PetscCall(VecCopy(u_swarm, u_swarm_old));
370       PetscCall(VecCopy(target, u_swarm));
371       PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
372       PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
373 
374       PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm.xmf"));
375 
376       PetscCall(DMSwarmSortGetAccess(dm_swarm));
377       PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
378       PetscCall(VecCopy(u_swarm_old, u_swarm));
379       PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
380       PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
381       PetscCall(VecDestroy(&u_swarm_old));
382     }
383 
384     // View solution at mesh points
385     PetscCall(VecViewFromOptions(X, NULL, "-solution_view"));
386 
387     // Compute L2 Error
388     {
389       // Set up error operator context
390       PetscCall(SetupErrorOperatorCtx(comm, dm_mesh, ceed, ceed_data, X_loc, op_error, op_error_ctx));
391       PetscScalar l2_error;
392       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
393 
394       PetscReal tol = 5e-4;
395       if (!test_mode || l2_error > tol) {
396         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
397         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
398         PetscCall(PetscPrintf(comm,
399                               "    L2 Error                                : %e\n"
400                               "    CG Solve Time                           : %g (%g) sec\n",
401                               (double)l2_error, rt_max, rt_min));
402       }
403     }
404     if (benchmark_mode && (!test_mode)) {
405       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
406                             1e-6 * g_size * its / rt_min));
407     }
408   }
409 
410   // Output solution
411   if (write_solution) {
412     PetscViewer vtk_viewer_soln;
413 
414     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
415     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
416     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
417     PetscCall(VecView(X, vtk_viewer_soln));
418     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
419   }
420 
421   // Cleanup
422   PetscCall(VecDestroy(&X));
423   PetscCall(VecDestroy(&X_loc));
424   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
425   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
426   PetscCall(MatDestroy(&mat_O));
427   PetscCall(PetscFree(op_apply_ctx));
428   PetscCall(PetscFree(op_error_ctx));
429   PetscCall(CeedDataDestroy(0, ceed_data));
430   PetscCall(DMDestroy(&dm_mesh));
431   PetscCall(DMDestroy(&dm_swarm));
432 
433   PetscCall(VecDestroy(&rhs));
434   PetscCall(VecDestroy(&target));
435   PetscCall(KSPDestroy(&ksp));
436   CeedOperatorDestroy(&op_error);
437   CeedDestroy(&ceed);
438   return PetscFinalize();
439 }
440