xref: /libCEED/examples/fluids/navierstokes.c (revision 2790b72b4f43887fa8322363f29d18765b4e7e19)
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
11 // Navier-Stokes problem.
12 //
13 // The code is intentionally "raw", using only low-level communication
14 // primitives.
15 //
16 // Build with:
17 //
18 //     make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes
19 //
20 // Sample runs:
21 //
22 //     ./navierstokes -ceed /cpu/self -problem density_current -degree 1
23 //     ./navierstokes -ceed /gpu/cuda -problem advection -degree 1
24 //
25 //TESTARGS(name="dc_explicit") -ceed {ceed_resource} -test -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. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-explicit.bin
26 //TESTARGS(name="dc_implicit_stab_none") -ceed {ceed_resource} -test -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. -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
27 //TESTARGS(name="adv_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection -strong_form 1 -degree 3 -dm_plex_box_faces 2,2,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. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-explicit-strong.bin
28 //TESTARGS(name="adv_rotation_implicit_sharp_cylinder") -ceed {ceed_resource} -test -problem advection -bubble_type cylinder -bubble_continuity back_sharp -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_z 1,2 -bc_wall 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 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-implicit-sharp-cylinder.bin
29 //TESTARGS(name="adv_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -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 -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
30 //TESTARGS(name="adv_translation_implicit_stab_su") -ceed {ceed_resource} -test -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 -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
31 //TESTARGS(name="adv2d_rotation_explicit_strong") -ceed {ceed_resource} -test -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
32 //TESTARGS(name="adv2d_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -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 -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
33 //TESTARGS(name="adv2d_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0 -dm_plex_box_upper 125,125 -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 -ts_type alpha -wind_type translation -wind_translation .53,-1.33,0 -bc_inflow 1,2,3,4 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-translation-implicit-stab-su.bin
34 //TESTARGS(name="euler_implicit") -ceed {ceed_resource} -test -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 -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 -problem euler_vortex -degree 3 -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 
37 /// @file
38 /// Navier-Stokes example using PETSc
39 
40 const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
41 
42 #include "navierstokes.h"
43 
44 int main(int argc, char **argv) {
45   // ---------------------------------------------------------------------------
46   // Initialize PETSc
47   // ---------------------------------------------------------------------------
48   PetscInt ierr;
49   ierr = PetscInitialize(&argc, &argv, NULL, help);
50   if (ierr) return ierr;
51 
52   // ---------------------------------------------------------------------------
53   // Create structs
54   // ---------------------------------------------------------------------------
55   AppCtx app_ctx;
56   ierr = PetscCalloc1(1, &app_ctx); CHKERRQ(ierr);
57 
58   ProblemData *problem = NULL;
59   ierr = PetscCalloc1(1, &problem); CHKERRQ(ierr);
60 
61   User user;
62   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
63 
64   CeedData ceed_data;
65   ierr = PetscCalloc1(1, &ceed_data); CHKERRQ(ierr);
66 
67   SimpleBC bc;
68   ierr = PetscCalloc1(1, &bc); CHKERRQ(ierr);
69 
70   SetupContext setup_ctx;
71   ierr = PetscCalloc1(1, &setup_ctx); CHKERRQ(ierr);
72 
73   Physics phys_ctx;
74   ierr = PetscCalloc1(1, &phys_ctx); CHKERRQ(ierr);
75 
76   Units units;
77   ierr = PetscCalloc1(1, &units); CHKERRQ(ierr);
78 
79   user->app_ctx = app_ctx;
80   user->units   = units;
81   user->phys    = phys_ctx;
82 
83   // ---------------------------------------------------------------------------
84   // Process command line options
85   // ---------------------------------------------------------------------------
86   // -- Register problems to be available on the command line
87   ierr = RegisterProblems_NS(app_ctx); CHKERRQ(ierr);
88 
89   // -- Process general command line options
90   MPI_Comm comm = PETSC_COMM_WORLD;
91   user->comm = comm;
92   ierr = ProcessCommandLineOptions(comm, app_ctx, bc); CHKERRQ(ierr);
93 
94   // ---------------------------------------------------------------------------
95   // Initialize libCEED
96   // ---------------------------------------------------------------------------
97   // -- Initialize backend
98   Ceed ceed;
99   CeedInit(app_ctx->ceed_resource, &ceed);
100   user->ceed = ceed;
101 
102   // -- Check preferred MemType
103   CeedMemType mem_type_backend;
104   CeedGetPreferredMemType(ceed, &mem_type_backend);
105 
106   // ---------------------------------------------------------------------------
107   // Set up global mesh
108   // ---------------------------------------------------------------------------
109   // -- Create DM
110   DM dm;
111   ierr = CreateDM(comm, problem, &dm); CHKERRQ(ierr);
112   VecType vec_type = NULL;
113   switch (mem_type_backend) {
114   case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
115   case CEED_MEM_DEVICE: {
116     const char *resolved;
117     CeedGetResource(ceed, &resolved);
118     if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
119     else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
120     else vec_type = VECSTANDARD;
121   }
122   }
123   ierr = DMSetVecType(dm, vec_type); CHKERRQ(ierr);
124   user->dm = dm;
125 
126   // ---------------------------------------------------------------------------
127   // Choose the problem from the list of registered problems
128   // ---------------------------------------------------------------------------
129   {
130     PetscErrorCode (*p)(ProblemData *, DM, void *, void *);
131     ierr = PetscFunctionListFind(app_ctx->problems, app_ctx->problem_name, &p);
132     CHKERRQ(ierr);
133     if (!p) SETERRQ(PETSC_COMM_SELF, 1, "Problem '%s' not found",
134                       app_ctx->problem_name);
135     ierr = (*p)(problem, dm, &setup_ctx, &user); CHKERRQ(ierr);
136   }
137 
138   // -- Set up DM
139   ierr = SetUpDM(dm, problem, app_ctx->degree, bc, phys_ctx, setup_ctx);
140   CHKERRQ(ierr);
141 
142   // -- Refine DM for high-order viz
143   if (app_ctx->viz_refine) {
144     ierr = VizRefineDM(dm, user, problem, bc, phys_ctx, setup_ctx);
145     CHKERRQ(ierr);
146   }
147 
148   // ---------------------------------------------------------------------------
149   // Set up libCEED
150   // ---------------------------------------------------------------------------
151   // -- Set up libCEED objects
152   ierr = SetupLibceed(ceed, ceed_data, dm, user, app_ctx, problem, bc, setup_ctx);
153   CHKERRQ(ierr);
154 
155   // ---------------------------------------------------------------------------
156   // Set up ICs
157   // ---------------------------------------------------------------------------
158   // -- Set up global state vector Q
159   Vec Q;
160   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
161   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
162 
163   // -- Set up local state vector Q_loc
164   Vec Q_loc;
165   ierr = DMGetLocalVector(dm, &Q_loc); CHKERRQ(ierr);
166 
167   // -- Fix multiplicity for ICs
168   ierr = ICs_FixMultiplicity(dm, ceed_data, Q_loc, Q, 0.0); CHKERRQ(ierr);
169 
170   // ---------------------------------------------------------------------------
171   // Set up lumped mass matrix
172   // ---------------------------------------------------------------------------
173   // -- Set up global mass vector
174   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
175 
176   // -- Compute lumped mass matrix
177   ierr = ComputeLumpedMassMatrix(ceed, dm, ceed_data, user->M); CHKERRQ(ierr);
178 
179   // ---------------------------------------------------------------------------
180   // Record boundary values from initial condition
181   // ---------------------------------------------------------------------------
182   // -- This overrides DMPlexInsertBoundaryValues().
183   //    We use this for the main simulation DM because the reference
184   //    DMPlexInsertBoundaryValues() is very slow. If we disable this, we should
185   //    still get the same results due to the problem->bc function, but with
186   //    potentially much slower execution.
187   if (1) {ierr = SetBCsFromICs_NS(dm, Q, Q_loc); CHKERRQ(ierr);}
188 
189   // ---------------------------------------------------------------------------
190   // Create output directory
191   // ---------------------------------------------------------------------------
192   PetscMPIInt rank;
193   MPI_Comm_rank(comm, &rank);
194   if (!rank) {ierr = PetscMkdir(app_ctx->output_dir); CHKERRQ(ierr);}
195 
196   // ---------------------------------------------------------------------------
197   // Gather initial Q values in case of continuation of simulation
198   // ---------------------------------------------------------------------------
199   // -- Set up initial values from binary file
200   if (app_ctx->cont_steps) {
201     ierr = SetupICsFromBinary(comm, app_ctx, Q); CHKERRQ(ierr);
202   }
203 
204   // ---------------------------------------------------------------------------
205   // Print problem summary
206   // ---------------------------------------------------------------------------
207   if (!app_ctx->test_mode) {
208     // Header and rank
209     char host_name[PETSC_MAX_PATH_LEN];
210     int  comm_size;
211     ierr = PetscGetHostName(host_name, sizeof host_name); CHKERRQ(ierr);
212     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
213     ierr = PetscPrintf(comm,
214                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
215                        "  MPI:\n"
216                        "    Host Name                          : %s\n"
217                        "    Total ranks                        : %d\n",
218                        host_name, comm_size); CHKERRQ(ierr);
219 
220     // Problem specific info
221     ierr = problem->print_info(phys_ctx, setup_ctx, app_ctx); CHKERRQ(ierr);
222 
223     // libCEED
224     const char *used_resource;
225     CeedGetResource(ceed, &used_resource);
226     ierr = PetscPrintf(comm,
227                        "  libCEED:\n"
228                        "    libCEED Backend                    : %s\n"
229                        "    libCEED Backend MemType            : %s\n",
230                        used_resource, CeedMemTypes[mem_type_backend]); CHKERRQ(ierr);
231     // PETSc
232     char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
233     if (problem->dim == 2) box_faces_str[3] = '\0';
234     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
235                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
236     ierr = PetscPrintf(comm,
237                        "  PETSc:\n"
238                        "    Box Faces                          : %s\n"
239                        "    Time Stepping Scheme               : %s\n",
240                        box_faces_str, phys_ctx->implicit ? "implicit" : "explicit"); CHKERRQ(ierr);
241     // Mesh
242     const PetscInt num_comp_q = 5;
243     CeedInt        glob_dofs, owned_dofs;
244     PetscInt       glob_nodes, owned_nodes;
245     const CeedInt  num_P = app_ctx->degree + 1,
246                    num_Q = num_P + app_ctx->q_extra;
247     // -- Get global size
248     ierr = VecGetSize(Q, &glob_dofs); CHKERRQ(ierr);
249     ierr = VecGetLocalSize(Q, &owned_dofs); CHKERRQ(ierr);
250     glob_nodes = glob_dofs/num_comp_q;
251     // -- Get local size
252     ierr = VecGetSize(Q_loc, &owned_nodes); CHKERRQ(ierr);
253     owned_nodes /= num_comp_q;
254     ierr = PetscPrintf(comm,
255                        "  Mesh:\n"
256                        "    Number of 1D Basis Nodes (P)       : %d\n"
257                        "    Number of 1D Quadrature Points (Q) : %d\n"
258                        "    Global DoFs                        : %" PetscInt_FMT "\n"
259                        "    Owned DoFs                         : %" PetscInt_FMT "\n"
260                        "    DoFs per node                      : %" PetscInt_FMT "\n"
261                        "    Global nodes                       : %" PetscInt_FMT "\n"
262                        "    Owned nodes                        : %" PetscInt_FMT "\n",
263                        num_P, num_Q, glob_dofs, owned_dofs, num_comp_q,
264                        glob_nodes, owned_nodes); CHKERRQ(ierr);
265   }
266   // -- Restore Q_loc
267   ierr = DMRestoreLocalVector(dm, &Q_loc); CHKERRQ(ierr);
268 
269   // ---------------------------------------------------------------------------
270   // TS: Create, setup, and solve
271   // ---------------------------------------------------------------------------
272   TS ts;
273   PetscScalar final_time;
274   ierr = TSSolve_NS(dm, user, app_ctx, phys_ctx, &Q, &final_time, &ts);
275   CHKERRQ(ierr);
276 
277   // ---------------------------------------------------------------------------
278   // Post-processing
279   // ---------------------------------------------------------------------------
280   ierr = PostProcess_NS(ts, ceed_data, dm, problem, app_ctx, Q, final_time);
281   CHKERRQ(ierr);
282 
283   // ---------------------------------------------------------------------------
284   // Destroy libCEED objects
285   // ---------------------------------------------------------------------------
286   // -- Vectors
287   CeedVectorDestroy(&ceed_data->x_coord);
288   CeedVectorDestroy(&ceed_data->q_data);
289   CeedVectorDestroy(&user->q_ceed);
290   CeedVectorDestroy(&user->q_dot_ceed);
291   CeedVectorDestroy(&user->g_ceed);
292 
293   // -- Contexts
294   CeedQFunctionContextDestroy(&ceed_data->setup_context);
295   CeedQFunctionContextDestroy(&ceed_data->newt_ig_context);
296   CeedQFunctionContextDestroy(&ceed_data->advection_context);
297   CeedQFunctionContextDestroy(&ceed_data->euler_context);
298 
299   // -- QFunctions
300   CeedQFunctionDestroy(&ceed_data->qf_setup_vol);
301   CeedQFunctionDestroy(&ceed_data->qf_ics);
302   CeedQFunctionDestroy(&ceed_data->qf_rhs_vol);
303   CeedQFunctionDestroy(&ceed_data->qf_ifunction_vol);
304   CeedQFunctionDestroy(&ceed_data->qf_setup_sur);
305   CeedQFunctionDestroy(&ceed_data->qf_apply_inflow);
306   CeedQFunctionDestroy(&ceed_data->qf_apply_outflow);
307 
308   // -- Bases
309   CeedBasisDestroy(&ceed_data->basis_q);
310   CeedBasisDestroy(&ceed_data->basis_x);
311   CeedBasisDestroy(&ceed_data->basis_xc);
312   CeedBasisDestroy(&ceed_data->basis_q_sur);
313   CeedBasisDestroy(&ceed_data->basis_x_sur);
314 
315   // -- Restrictions
316   CeedElemRestrictionDestroy(&ceed_data->elem_restr_q);
317   CeedElemRestrictionDestroy(&ceed_data->elem_restr_x);
318   CeedElemRestrictionDestroy(&ceed_data->elem_restr_qd_i);
319 
320   // -- Operators
321   CeedOperatorDestroy(&ceed_data->op_setup_vol);
322   CeedOperatorDestroy(&ceed_data->op_ics);
323   CeedOperatorDestroy(&user->op_rhs_vol);
324   CeedOperatorDestroy(&user->op_ifunction_vol);
325   CeedOperatorDestroy(&user->op_rhs);
326   CeedOperatorDestroy(&user->op_ifunction);
327 
328   // -- Ceed
329   CeedDestroy(&ceed);
330 
331   // ---------------------------------------------------------------------------
332   // Clean up PETSc
333   // ---------------------------------------------------------------------------
334   // -- Vectors
335   ierr = VecDestroy(&Q); CHKERRQ(ierr);
336   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
337 
338   // -- Matrices
339   ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr);
340 
341   // -- DM
342   ierr = DMDestroy(&dm); CHKERRQ(ierr);
343   ierr = DMDestroy(&user->dm_viz); CHKERRQ(ierr);
344 
345   // -- TS
346   ierr = TSDestroy(&ts); CHKERRQ(ierr);
347 
348   // -- Function list
349   ierr = PetscFunctionListDestroy(&app_ctx->problems); CHKERRQ(ierr);
350 
351   // -- Structs
352   ierr = PetscFree(units); CHKERRQ(ierr);
353   ierr = PetscFree(user); CHKERRQ(ierr);
354   ierr = PetscFree(problem); CHKERRQ(ierr);
355   ierr = PetscFree(bc); CHKERRQ(ierr);
356   ierr = PetscFree(setup_ctx); CHKERRQ(ierr);
357   ierr = PetscFree(phys_ctx->newtonian_ig_ctx); CHKERRQ(ierr);
358   ierr = PetscFree(phys_ctx->euler_ctx); CHKERRQ(ierr);
359   ierr = PetscFree(phys_ctx->advection_ctx); CHKERRQ(ierr);
360   ierr = PetscFree(phys_ctx); CHKERRQ(ierr);
361   ierr = PetscFree(app_ctx); CHKERRQ(ierr);
362   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
363 
364   return PetscFinalize();
365 }
366