xref: /libCEED/examples/fluids/navierstokes.c (revision f545224771c1a6d188cc0fad9830d905771d2bc9)
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     PetscCheck(p, 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 
148   // -- Refine DM for high-order viz
149   if (app_ctx->viz_refine) PetscCall(VizRefineDM(dm, user, problem, bc, phys_ctx));
150 
151   // ---------------------------------------------------------------------------
152   // Create solution vectors
153   // ---------------------------------------------------------------------------
154   // -- Set up global state vector Q
155   Vec Q;
156   PetscCall(DMCreateGlobalVector(dm, &Q));
157   PetscCall(VecZeroEntries(Q));
158 
159   // -- Set up local state vectors Q_loc, Q_dot_loc
160   PetscCall(DMCreateLocalVector(dm, &user->Q_loc));
161   PetscCall(DMCreateLocalVector(dm, &user->Q_dot_loc));
162   PetscCall(VecZeroEntries(user->Q_dot_loc));
163 
164   // ---------------------------------------------------------------------------
165   // Set up libCEED
166   // ---------------------------------------------------------------------------
167   // -- Set up libCEED objects
168   PetscCall(SetupLibceed(ceed, ceed_data, dm, user, app_ctx, problem, bc));
169 
170   // ---------------------------------------------------------------------------
171   // Set up ICs
172   // ---------------------------------------------------------------------------
173   // -- Fix multiplicity for ICs
174   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, user->Q_loc, Q, 0.0));
175 
176   // ---------------------------------------------------------------------------
177   // Set up lumped mass matrix
178   // ---------------------------------------------------------------------------
179   // -- Set up global mass vector
180   PetscCall(VecDuplicate(Q, &user->M));
181 
182   // -- Compute lumped mass matrix
183   PetscCall(ComputeLumpedMassMatrix(ceed, dm, ceed_data, user->M));
184 
185   // ---------------------------------------------------------------------------
186   // Record boundary values from initial condition
187   // ---------------------------------------------------------------------------
188   // -- This overrides DMPlexInsertBoundaryValues().
189   //    We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow on the GPU due to extra device-to-host
190   //    communication. If we disable this, we should still get the same results due to the problem->bc function, but with potentially much slower
191   //    execution.
192   if (problem->bc_from_ics) {
193     PetscCall(SetBCsFromICs_NS(dm, Q, user->Q_loc));
194   }
195 
196   // ---------------------------------------------------------------------------
197   // Create output directory
198   // ---------------------------------------------------------------------------
199   PetscMPIInt rank;
200   MPI_Comm_rank(comm, &rank);
201   if (!rank) {
202     PetscCall(PetscMkdir(app_ctx->output_dir));
203   }
204 
205   // ---------------------------------------------------------------------------
206   // Gather initial Q values in case of continuation of simulation
207   // ---------------------------------------------------------------------------
208   // -- Set up initial values from binary file
209   if (app_ctx->cont_steps) {
210     PetscCall(SetupICsFromBinary(comm, app_ctx, Q));
211   }
212 
213   // ---------------------------------------------------------------------------
214   // Print problem summary
215   // ---------------------------------------------------------------------------
216   if (app_ctx->test_type == TESTTYPE_NONE) {
217     // Header and rank
218     char host_name[PETSC_MAX_PATH_LEN];
219     int  comm_size;
220     PetscCall(PetscGetHostName(host_name, sizeof host_name));
221     PetscCall(MPI_Comm_size(comm, &comm_size));
222     PetscCall(PetscPrintf(comm,
223                           "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
224                           "  MPI:\n"
225                           "    Host Name                          : %s\n"
226                           "    Total ranks                        : %d\n",
227                           host_name, comm_size));
228 
229     // Problem specific info
230     PetscCall(problem->print_info(problem, app_ctx));
231 
232     // libCEED
233     const char *used_resource;
234     CeedGetResource(ceed, &used_resource);
235     PetscCall(PetscPrintf(comm,
236                           "  libCEED:\n"
237                           "    libCEED Backend                    : %s\n"
238                           "    libCEED Backend MemType            : %s\n",
239                           used_resource, CeedMemTypes[mem_type_backend]));
240     // PETSc
241     char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
242     if (problem->dim == 2) box_faces_str[3] = '\0';
243     PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
244     MatType mat_type;
245     VecType vec_type;
246     PetscCall(DMGetMatType(dm, &mat_type));
247     PetscCall(DMGetVecType(dm, &vec_type));
248     PetscCall(PetscPrintf(comm,
249                           "  PETSc:\n"
250                           "    Box Faces                          : %s\n"
251                           "    DM MatType                         : %s\n"
252                           "    DM VecType                         : %s\n"
253                           "    Time Stepping Scheme               : %s\n",
254                           box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
255     if (app_ctx->cont_steps) {
256       PetscCall(PetscPrintf(comm,
257                             "  Continue:\n"
258                             "    Filename:                          : %s\n"
259                             "    Step:                              : %" PetscInt_FMT "\n"
260                             "    Time:                              : %g\n",
261                             app_ctx->cont_file, app_ctx->cont_steps, app_ctx->cont_time));
262     }
263     // Mesh
264     const PetscInt num_comp_q = 5;
265     CeedInt        glob_dofs, owned_dofs;
266     PetscInt       glob_nodes, local_nodes;
267     const CeedInt  num_P = app_ctx->degree + 1, num_Q = num_P + app_ctx->q_extra;
268     // -- Get global size
269     PetscCall(VecGetSize(Q, &glob_dofs));
270     PetscCall(VecGetLocalSize(Q, &owned_dofs));
271     glob_nodes = glob_dofs / num_comp_q;
272     // -- Get local size
273     PetscCall(VecGetSize(user->Q_loc, &local_nodes));
274     local_nodes /= num_comp_q;
275     PetscCall(PetscPrintf(comm,
276                           "  Mesh:\n"
277                           "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
278                           "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
279                           "    Global DoFs                        : %" PetscInt_FMT "\n"
280                           "    Owned DoFs                         : %" PetscInt_FMT "\n"
281                           "    DoFs per node                      : %" PetscInt_FMT "\n"
282                           "    Global nodes (DoFs / %" PetscInt_FMT ")            : %" PetscInt_FMT "\n"
283                           "    Local nodes                        : %" PetscInt_FMT "\n",
284                           num_P, num_Q, glob_dofs, owned_dofs, num_comp_q, num_comp_q, glob_nodes, local_nodes));
285   }
286   // -- Zero Q_loc
287   PetscCall(VecZeroEntries(user->Q_loc));
288 
289   // ---------------------------------------------------------------------------
290   // TS: Create, setup, and solve
291   // ---------------------------------------------------------------------------
292   TS          ts;
293   PetscScalar final_time;
294   PetscCall(TSSolve_NS(dm, user, app_ctx, phys_ctx, &Q, &final_time, &ts));
295 
296   // ---------------------------------------------------------------------------
297   // Post-processing
298   // ---------------------------------------------------------------------------
299   PetscCall(PostProcess_NS(ts, ceed_data, dm, problem, user, Q, final_time));
300 
301   // ---------------------------------------------------------------------------
302   // Destroy libCEED objects
303   // ---------------------------------------------------------------------------
304 
305   PetscCall(TurbulenceStatisticsDestroy(user, ceed_data));
306   PetscCall(NodalProjectionDataDestroy(user->grad_velo_proj));
307   PetscCall(SGS_DD_DataDestroy(user->sgs_dd_data));
308 
309   // -- Vectors
310   CeedVectorDestroy(&ceed_data->x_coord);
311   CeedVectorDestroy(&ceed_data->q_data);
312   CeedVectorDestroy(&user->q_ceed);
313   CeedVectorDestroy(&user->q_dot_ceed);
314   CeedVectorDestroy(&user->g_ceed);
315   CeedVectorDestroy(&user->coo_values_amat);
316   CeedVectorDestroy(&user->coo_values_pmat);
317 
318   // -- Bases
319   CeedBasisDestroy(&ceed_data->basis_q);
320   CeedBasisDestroy(&ceed_data->basis_x);
321   CeedBasisDestroy(&ceed_data->basis_xc);
322   CeedBasisDestroy(&ceed_data->basis_q_sur);
323   CeedBasisDestroy(&ceed_data->basis_x_sur);
324 
325   // -- Restrictions
326   CeedElemRestrictionDestroy(&ceed_data->elem_restr_q);
327   CeedElemRestrictionDestroy(&ceed_data->elem_restr_x);
328   CeedElemRestrictionDestroy(&ceed_data->elem_restr_qd_i);
329 
330   // Destroy QFunction contexts after using
331   // ToDo: Simplify tracked libCEED objects, smaller struct
332   {
333     CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context);
334     CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context);
335     CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context);
336     CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context);
337     CeedQFunctionContextDestroy(&problem->apply_freestream_jacobian.qfunction_context);
338     CeedQFunctionContextDestroy(&problem->apply_freestream_jacobian.qfunction_context);
339     CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context);
340     CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context);
341     CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
342     CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context);
343     CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context);
344     CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context);
345   }
346 
347   // -- QFunctions
348   CeedQFunctionDestroy(&ceed_data->qf_setup_vol);
349   CeedQFunctionDestroy(&ceed_data->qf_ics);
350   CeedQFunctionDestroy(&ceed_data->qf_rhs_vol);
351   CeedQFunctionDestroy(&ceed_data->qf_ifunction_vol);
352   CeedQFunctionDestroy(&ceed_data->qf_setup_sur);
353   CeedQFunctionDestroy(&ceed_data->qf_apply_inflow);
354   CeedQFunctionDestroy(&ceed_data->qf_apply_inflow_jacobian);
355   CeedQFunctionDestroy(&ceed_data->qf_apply_freestream);
356   CeedQFunctionDestroy(&ceed_data->qf_apply_freestream_jacobian);
357 
358   // -- Operators
359   CeedOperatorDestroy(&ceed_data->op_setup_vol);
360   PetscCall(OperatorApplyContextDestroy(ceed_data->op_ics_ctx));
361   CeedOperatorDestroy(&user->op_rhs_vol);
362   CeedOperatorDestroy(&user->op_ifunction_vol);
363   CeedOperatorDestroy(&user->op_rhs);
364   CeedOperatorDestroy(&user->op_ifunction);
365   CeedOperatorDestroy(&user->op_ijacobian);
366 
367   // -- Ceed
368   CeedDestroy(&ceed);
369 
370   // ---------------------------------------------------------------------------
371   // Clean up PETSc
372   // ---------------------------------------------------------------------------
373   // -- Vectors
374   PetscCall(VecDestroy(&Q));
375   PetscCall(VecDestroy(&user->M));
376   PetscCall(VecDestroy(&user->Q_loc));
377   PetscCall(VecDestroy(&user->Q_dot_loc));
378 
379   // -- Matrices
380   PetscCall(MatDestroy(&user->interp_viz));
381 
382   // -- DM
383   PetscCall(DMDestroy(&dm));
384   PetscCall(DMDestroy(&user->dm_viz));
385 
386   // -- TS
387   PetscCall(TSDestroy(&ts));
388 
389   // -- Function list
390   PetscCall(PetscFunctionListDestroy(&app_ctx->problems));
391 
392   PetscCall(PetscFree(app_ctx->amat_type));
393   PetscCall(PetscFree(app_ctx->wall_forces.walls));
394   PetscCall(PetscViewerDestroy(&app_ctx->wall_forces.viewer));
395 
396   // -- Structs
397   PetscCall(PetscFree(units));
398   PetscCall(PetscFree(user));
399   PetscCall(PetscFree(problem));
400   PetscCall(PetscFree(bc));
401   PetscCall(PetscFree(phys_ctx));
402   PetscCall(PetscFree(app_ctx));
403   PetscCall(PetscFree(ceed_data));
404 
405   return PetscFinalize();
406 }
407