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