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