xref: /libCEED/examples/fluids/navierstokes.c (revision 5fb68f377259d3910de46d787b7c5d1587fd01e1)
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: Navier-Stokes
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve a Navier-Stokes problem.
11 //
12 // Build with:
13 //
14 //     make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes
15 //
16 // Sample runs:
17 //
18 //     ./navierstokes -ceed /cpu/self -problem density_current -degree 1
19 //     ./navierstokes -ceed /gpu/cuda -problem advection -degree 1
20 //
21 //TESTARGS(name="turb_spanstats") -ceed {ceed_resource} -test_type turb_spanstats -options_file examples/fluids/tests-output/stats_test.yaml -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-turb-spanstats-stats.bin
22 //TESTARGS(name="blasius") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/tests-output/blasius_test.yaml -compare_final_state_atol 2E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius.bin
23 //TESTARGS(name="blasius_STG") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 2E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG.bin
24 //TESTARGS(name="blasius_STG_weakT") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG_weakT.bin -weakT
25 //TESTARGS(name="blasius_STG_strongBC") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 1E-10 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG_strongBC.bin -stg_strong true
26 //TESTARGS(name="channel") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/channel.yaml -compare_final_state_atol 2e-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-channel.bin -dm_plex_box_faces 5,5,1 -ts_max_steps 5
27 //TESTARGS(name="channel-primitive") -ceed {ceed_resource} -test_type solver -options_file examples/fluids/channel.yaml -compare_final_state_atol 2e-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-channel-prim.bin -dm_plex_box_faces 5,5,1 -ts_max_steps 5 -state_var primitive
28 //TESTARGS(name="dc_explicit") -ceed {ceed_resource} -test_type solver -degree 3 -q_extra 2 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -units_kilogram 1e-9 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -mu 75 -ts_dt 1e-3 -units_meter 1e-2 -units_second 1e-2 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-explicit.bin
29 //TESTARGS(name="dc_implicit_stab_none") -ceed {ceed_resource} -test_type solver -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -units_kilogram 1e-9 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -mu 75 -units_meter 1e-2 -units_second 1e-2 -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-implicit-stab-none.bin
30 //TESTARGS(name="adv_rotation_implicit_stab_supg") -ceed {ceed_resource} -test_type solver -problem advection -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_wall 1,2,3,4,5,6 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-implicit-stab-supg.bin
31 //TESTARGS(name="adv_translation_implicit_stab_su") -ceed {ceed_resource} -test_type solver -problem advection -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -wind_type translation -wind_translation .53,-1.33,-2.65 -bc_inflow 1,2,3,4,5,6 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-translation-implicit-stab-su.bin
32 //TESTARGS(name="adv2d_rotation_explicit_strong") -ceed {ceed_resource} -test_type solver -problem advection2d -strong_form 1 -degree 3 -dm_plex_box_faces 2,2 -dm_plex_box_lower 0,0 -dm_plex_box_upper 125,125 -bc_wall 1,2,3,4 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-explicit-strong.bin
33 //TESTARGS(name="adv2d_rotation_implicit_stab_supg") -ceed {ceed_resource} -test_type solver -problem advection2d -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0 -dm_plex_box_upper 125,125 -bc_wall 1,2,3,4 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-implicit-stab-supg.bin
34 //TESTARGS(name="euler_implicit") -ceed {ceed_resource} -test_type solver -problem euler_vortex -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_meter 1e-4 -units_second 1e-4 -mean_velocity 1.4,-2.,0 -bc_inflow 4,6 -bc_outflow 3,5 -bc_slip_z 1,2 -vortex_strength 2 -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-euler-implicit.bin
35 //TESTARGS(name="euler_explicit") -ceed {ceed_resource} -test_type solver -problem euler_vortex -degree 3 -q_extra 2 -dm_plex_box_faces 2,2,1 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_meter 1e-4 -units_second 1e-4 -mean_velocity 1.4,-2.,0 -bc_inflow 4,6 -bc_outflow 3,5 -bc_slip_z 1,2 -vortex_strength 2 -ts_dt 1e-7 -ts_rk_type 5bs -ts_rtol 1e-10 -ts_atol 1e-10 -compare_final_state_atol 1E-7 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-euler-explicit.bin
36 //TESTARGS(name="shocktube_explicit_su_yzb") -ceed {ceed_resource} -test_type solver -problem shocktube -degree 1 -q_extra 2 -dm_plex_box_faces 50,1,1 -units_meter 1e-2 units_second 1e-2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1000,20,20 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -yzb -stab su -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-shocktube-explicit-su-yzb.bin
37 
38 /// @file
39 /// Navier-Stokes example using PETSc
40 
41 const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
42 
43 #include "navierstokes.h"
44 
45 int main(int argc, char **argv) {
46   // ---------------------------------------------------------------------------
47   // Initialize PETSc
48   // ---------------------------------------------------------------------------
49   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
50 
51   // ---------------------------------------------------------------------------
52   // Create structs
53   // ---------------------------------------------------------------------------
54   AppCtx app_ctx;
55   PetscCall(PetscCalloc1(1, &app_ctx));
56 
57   ProblemData *problem = NULL;
58   PetscCall(PetscCalloc1(1, &problem));
59 
60   User user;
61   PetscCall(PetscCalloc1(1, &user));
62 
63   CeedData ceed_data;
64   PetscCall(PetscCalloc1(1, &ceed_data));
65 
66   SimpleBC bc;
67   PetscCall(PetscCalloc1(1, &bc));
68 
69   Physics phys_ctx;
70   PetscCall(PetscCalloc1(1, &phys_ctx));
71 
72   Units units;
73   PetscCall(PetscCalloc1(1, &units));
74 
75   user->app_ctx        = app_ctx;
76   user->units          = units;
77   user->phys           = phys_ctx;
78   problem->bc_from_ics = PETSC_TRUE;
79 
80   // ---------------------------------------------------------------------------
81   // Process command line options
82   // ---------------------------------------------------------------------------
83   // -- Register problems to be available on the command line
84   PetscCall(RegisterProblems_NS(app_ctx));
85 
86   // -- Process general command line options
87   MPI_Comm comm = PETSC_COMM_WORLD;
88   user->comm    = comm;
89   PetscCall(ProcessCommandLineOptions(comm, app_ctx, bc));
90 
91   // ---------------------------------------------------------------------------
92   // Initialize libCEED
93   // ---------------------------------------------------------------------------
94   // -- Initialize backend
95   Ceed ceed;
96   CeedInit(app_ctx->ceed_resource, &ceed);
97   user->ceed = ceed;
98 
99   // -- Check preferred MemType
100   CeedMemType mem_type_backend;
101   CeedGetPreferredMemType(ceed, &mem_type_backend);
102 
103   // ---------------------------------------------------------------------------
104   // Set up global mesh
105   // ---------------------------------------------------------------------------
106   // -- Create DM
107   DM      dm;
108   VecType vec_type = NULL;
109   MatType mat_type = NULL;
110   switch (mem_type_backend) {
111     case CEED_MEM_HOST:
112       vec_type = VECSTANDARD;
113       break;
114     case CEED_MEM_DEVICE: {
115       const char *resolved;
116       CeedGetResource(ceed, &resolved);
117       if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
118       else if (strstr(resolved, "/gpu/hip")) vec_type = VECKOKKOS;
119       else vec_type = VECSTANDARD;
120     }
121   }
122   if (strstr(vec_type, VECCUDA)) mat_type = MATAIJCUSPARSE;
123   else if (strstr(vec_type, VECKOKKOS)) mat_type = MATAIJKOKKOS;
124   else mat_type = MATAIJ;
125   PetscCall(CreateDM(comm, problem, mat_type, vec_type, &dm));
126   user->dm = dm;
127   PetscCall(DMSetApplicationContext(dm, user));
128 
129   // ---------------------------------------------------------------------------
130   // Choose the problem from the list of registered problems
131   // ---------------------------------------------------------------------------
132   {
133     PetscErrorCode (*p)(ProblemData *, DM, void *, SimpleBC);
134     PetscCall(PetscFunctionListFind(app_ctx->problems, app_ctx->problem_name, &p));
135     if (!p) SETERRQ(PETSC_COMM_SELF, 1, "Problem '%s' not found", app_ctx->problem_name);
136     PetscCall((*p)(problem, dm, &user, bc));
137   }
138 
139   // -- Set up DM
140   PetscCall(SetUpDM(dm, problem, app_ctx->degree, bc, phys_ctx));
141   PetscCall(CreateStatsDM(user, problem, app_ctx->degree, bc));
142   app_ctx->wall_forces.num_wall = bc->num_wall;
143   PetscMalloc1(bc->num_wall, &app_ctx->wall_forces.walls);
144   PetscCall(PetscArraycpy(app_ctx->wall_forces.walls, bc->walls, bc->num_wall));
145 
146   // -- Refine DM for high-order viz
147   if (app_ctx->viz_refine) {
148     PetscCall(VizRefineDM(dm, user, problem, bc, phys_ctx));
149   }
150 
151   // ---------------------------------------------------------------------------
152   // Set up libCEED
153   // ---------------------------------------------------------------------------
154   // -- Set up libCEED objects
155   PetscCall(SetupLibceed(ceed, ceed_data, dm, user, app_ctx, problem, bc));
156 
157   if (app_ctx->turb_spanstats_enable) {
158     PetscCall(SetupStatsCollection(ceed, user, ceed_data, problem));
159   }
160 
161   // Destroy QFunction contexts after using
162   // ToDo: Simplify tracked libCEED objects, smaller struct
163   {
164     CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context);
165     CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context);
166     CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context);
167     CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context);
168     CeedQFunctionContextDestroy(&problem->apply_freestream_jacobian.qfunction_context);
169     CeedQFunctionContextDestroy(&problem->apply_freestream_jacobian.qfunction_context);
170     CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context);
171     CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context);
172     CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
173     CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context);
174     CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context);
175     CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context);
176   }
177 
178   // ---------------------------------------------------------------------------
179   // Set up ICs
180   // ---------------------------------------------------------------------------
181   // -- Set up global state vector Q
182   Vec Q;
183   PetscCall(DMCreateGlobalVector(dm, &Q));
184   PetscCall(VecZeroEntries(Q));
185 
186   // -- Set up local state vectors Q_loc, Q_dot_loc
187   PetscCall(DMCreateLocalVector(dm, &user->Q_loc));
188   PetscCall(DMCreateLocalVector(dm, &user->Q_dot_loc));
189   PetscCall(VecZeroEntries(user->Q_dot_loc));
190 
191   // -- Fix multiplicity for ICs
192   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, user->Q_loc, Q, 0.0));
193 
194   // ---------------------------------------------------------------------------
195   // Set up lumped mass matrix
196   // ---------------------------------------------------------------------------
197   // -- Set up global mass vector
198   PetscCall(VecDuplicate(Q, &user->M));
199 
200   // -- Compute lumped mass matrix
201   PetscCall(ComputeLumpedMassMatrix(ceed, dm, ceed_data, user->M));
202 
203   // ---------------------------------------------------------------------------
204   // Record boundary values from initial condition
205   // ---------------------------------------------------------------------------
206   // -- This overrides DMPlexInsertBoundaryValues().
207   //    We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow on the GPU due to extra device-to-host
208   //    communication. If we disable this, we should still get the same results due to the problem->bc function, but with potentially much slower
209   //    execution.
210   if (problem->bc_from_ics) {
211     PetscCall(SetBCsFromICs_NS(dm, Q, user->Q_loc));
212   }
213 
214   // ---------------------------------------------------------------------------
215   // Create output directory
216   // ---------------------------------------------------------------------------
217   PetscMPIInt rank;
218   MPI_Comm_rank(comm, &rank);
219   if (!rank) {
220     PetscCall(PetscMkdir(app_ctx->output_dir));
221   }
222 
223   // ---------------------------------------------------------------------------
224   // Gather initial Q values in case of continuation of simulation
225   // ---------------------------------------------------------------------------
226   // -- Set up initial values from binary file
227   if (app_ctx->cont_steps) {
228     PetscCall(SetupICsFromBinary(comm, app_ctx, Q));
229   }
230 
231   // ---------------------------------------------------------------------------
232   // Print problem summary
233   // ---------------------------------------------------------------------------
234   if (app_ctx->test_type == TESTTYPE_NONE) {
235     // Header and rank
236     char host_name[PETSC_MAX_PATH_LEN];
237     int  comm_size;
238     PetscCall(PetscGetHostName(host_name, sizeof host_name));
239     PetscCall(MPI_Comm_size(comm, &comm_size));
240     PetscCall(PetscPrintf(comm,
241                           "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
242                           "  MPI:\n"
243                           "    Host Name                          : %s\n"
244                           "    Total ranks                        : %d\n",
245                           host_name, comm_size));
246 
247     // Problem specific info
248     PetscCall(problem->print_info(problem, app_ctx));
249 
250     // libCEED
251     const char *used_resource;
252     CeedGetResource(ceed, &used_resource);
253     PetscCall(PetscPrintf(comm,
254                           "  libCEED:\n"
255                           "    libCEED Backend                    : %s\n"
256                           "    libCEED Backend MemType            : %s\n",
257                           used_resource, CeedMemTypes[mem_type_backend]));
258     // PETSc
259     char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
260     if (problem->dim == 2) box_faces_str[3] = '\0';
261     PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
262     MatType mat_type;
263     VecType vec_type;
264     PetscCall(DMGetMatType(dm, &mat_type));
265     PetscCall(DMGetVecType(dm, &vec_type));
266     PetscCall(PetscPrintf(comm,
267                           "  PETSc:\n"
268                           "    Box Faces                          : %s\n"
269                           "    DM MatType                         : %s\n"
270                           "    DM VecType                         : %s\n"
271                           "    Time Stepping Scheme               : %s\n",
272                           box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
273     if (app_ctx->cont_steps) {
274       PetscCall(PetscPrintf(comm,
275                             "  Continue:\n"
276                             "    Filename:                          : %s\n"
277                             "    Step:                              : %" PetscInt_FMT "\n"
278                             "    Time:                              : %g\n",
279                             app_ctx->cont_file, app_ctx->cont_steps, app_ctx->cont_time));
280     }
281     // Mesh
282     const PetscInt num_comp_q = 5;
283     CeedInt        glob_dofs, owned_dofs;
284     PetscInt       glob_nodes, local_nodes;
285     const CeedInt  num_P = app_ctx->degree + 1, num_Q = num_P + app_ctx->q_extra;
286     // -- Get global size
287     PetscCall(VecGetSize(Q, &glob_dofs));
288     PetscCall(VecGetLocalSize(Q, &owned_dofs));
289     glob_nodes = glob_dofs / num_comp_q;
290     // -- Get local size
291     PetscCall(VecGetSize(user->Q_loc, &local_nodes));
292     local_nodes /= num_comp_q;
293     PetscCall(PetscPrintf(comm,
294                           "  Mesh:\n"
295                           "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
296                           "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
297                           "    Global DoFs                        : %" PetscInt_FMT "\n"
298                           "    Owned DoFs                         : %" PetscInt_FMT "\n"
299                           "    DoFs per node                      : %" PetscInt_FMT "\n"
300                           "    Global nodes (DoFs / %" PetscInt_FMT ")            : %" PetscInt_FMT "\n"
301                           "    Local nodes                        : %" PetscInt_FMT "\n",
302                           num_P, num_Q, glob_dofs, owned_dofs, num_comp_q, num_comp_q, glob_nodes, local_nodes));
303   }
304   // -- Zero Q_loc
305   PetscCall(VecZeroEntries(user->Q_loc));
306 
307   // ---------------------------------------------------------------------------
308   // TS: Create, setup, and solve
309   // ---------------------------------------------------------------------------
310   TS          ts;
311   PetscScalar final_time;
312   PetscCall(TSSolve_NS(dm, user, app_ctx, phys_ctx, &Q, &final_time, &ts));
313 
314   // ---------------------------------------------------------------------------
315   // Post-processing
316   // ---------------------------------------------------------------------------
317   PetscCall(PostProcess_NS(ts, ceed_data, dm, problem, user, Q, final_time));
318 
319   // ---------------------------------------------------------------------------
320   // Destroy libCEED objects
321   // ---------------------------------------------------------------------------
322 
323   PetscCall(DestroyStats(user, ceed_data));
324 
325   // -- Vectors
326   CeedVectorDestroy(&ceed_data->x_coord);
327   CeedVectorDestroy(&ceed_data->q_data);
328   CeedVectorDestroy(&user->q_ceed);
329   CeedVectorDestroy(&user->q_dot_ceed);
330   CeedVectorDestroy(&user->g_ceed);
331   CeedVectorDestroy(&user->coo_values_amat);
332   CeedVectorDestroy(&user->coo_values_pmat);
333 
334   // -- Bases
335   CeedBasisDestroy(&ceed_data->basis_q);
336   CeedBasisDestroy(&ceed_data->basis_x);
337   CeedBasisDestroy(&ceed_data->basis_xc);
338   CeedBasisDestroy(&ceed_data->basis_q_sur);
339   CeedBasisDestroy(&ceed_data->basis_x_sur);
340 
341   // -- Restrictions
342   CeedElemRestrictionDestroy(&ceed_data->elem_restr_q);
343   CeedElemRestrictionDestroy(&ceed_data->elem_restr_x);
344   CeedElemRestrictionDestroy(&ceed_data->elem_restr_qd_i);
345 
346   // -- QFunctions
347   CeedQFunctionDestroy(&ceed_data->qf_setup_vol);
348   CeedQFunctionDestroy(&ceed_data->qf_ics);
349   CeedQFunctionDestroy(&ceed_data->qf_rhs_vol);
350   CeedQFunctionDestroy(&ceed_data->qf_ifunction_vol);
351   CeedQFunctionDestroy(&ceed_data->qf_setup_sur);
352   CeedQFunctionDestroy(&ceed_data->qf_apply_inflow);
353   CeedQFunctionDestroy(&ceed_data->qf_apply_inflow_jacobian);
354   CeedQFunctionDestroy(&ceed_data->qf_apply_freestream);
355   CeedQFunctionDestroy(&ceed_data->qf_apply_freestream_jacobian);
356 
357   // -- Operators
358   CeedOperatorDestroy(&ceed_data->op_setup_vol);
359   CeedOperatorDestroy(&ceed_data->op_ics);
360   CeedOperatorDestroy(&user->op_rhs_vol);
361   CeedOperatorDestroy(&user->op_ifunction_vol);
362   CeedOperatorDestroy(&user->op_rhs);
363   CeedOperatorDestroy(&user->op_ifunction);
364   CeedOperatorDestroy(&user->op_ijacobian);
365 
366   // -- Ceed
367   CeedDestroy(&ceed);
368 
369   // ---------------------------------------------------------------------------
370   // Clean up PETSc
371   // ---------------------------------------------------------------------------
372   // -- Vectors
373   PetscCall(VecDestroy(&Q));
374   PetscCall(VecDestroy(&user->M));
375   PetscCall(VecDestroy(&user->Q_loc));
376   PetscCall(VecDestroy(&user->Q_dot_loc));
377 
378   // -- Matrices
379   PetscCall(MatDestroy(&user->interp_viz));
380 
381   // -- DM
382   PetscCall(DMDestroy(&dm));
383   PetscCall(DMDestroy(&user->dm_viz));
384 
385   // -- TS
386   PetscCall(TSDestroy(&ts));
387 
388   // -- Function list
389   PetscCall(PetscFunctionListDestroy(&app_ctx->problems));
390 
391   PetscCall(PetscFree(app_ctx->amat_type));
392   PetscCall(PetscFree(app_ctx->wall_forces.walls));
393   PetscCall(PetscViewerDestroy(&app_ctx->wall_forces.viewer));
394 
395   // -- Structs
396   PetscCall(PetscFree(units));
397   PetscCall(PetscFree(user));
398   PetscCall(PetscFree(problem));
399   PetscCall(PetscFree(bc));
400   PetscCall(PetscFree(phys_ctx));
401   PetscCall(PetscFree(app_ctx));
402   PetscCall(PetscFree(ceed_data));
403 
404   return PetscFinalize();
405 }
406