Path: blob/main/examples/structured_1d_dgsem/elixir_euler_sedov.jl
5585 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations1D(1.4)78"""9initial_condition_sedov_blast_wave(x, t, equations::CompressibleEulerEquations1D)1011The Sedov blast wave setup based on example 35.1.4 from Flash12- https://flash.rochester.edu/site/flashcode/user_support/flash4_ug_4p8.pdf13"""14function initial_condition_sedov_blast_wave(x, t, equations::CompressibleEulerEquations1D)15# Set up polar coordinates16inicenter = SVector(0.0)17x_norm = x[1] - inicenter[1]18r = abs(x_norm)1920# Setup based on example 35.1.4 in https://flash.rochester.edu/site/flashcode/user_support/flash4_ug_4p8.pdf21r0 = 0.21875 # = 3.5 * smallest dx (for domain length=4 and max-ref=6)22# r0 = 0.5 # = more reasonable setup23E = 1.024p0_inner = 6 * (equations.gamma - 1) * E / (3 * pi * r0)25p0_outer = 1.0e-5 # = true Sedov setup2627# Calculate primitive variables28rho = 1.029v1 = 0.030p = r > r0 ? p0_outer : p0_inner3132return prim2cons(SVector(rho, v1, p), equations)33end3435initial_condition = initial_condition_sedov_blast_wave3637# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of38# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.39# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.40# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.41# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.42# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the43# `StepsizeCallback` (CFL-Condition) and less diffusion.44surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)45volume_flux = flux_ranocha46polydeg = 347basis = LobattoLegendreBasis(polydeg)48shock_indicator_variable = density_pressure49indicator_sc = IndicatorHennemannGassner(equations, basis,50alpha_max = 1.0,51alpha_min = 0.001,52alpha_smooth = true,53variable = shock_indicator_variable)54volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;55volume_flux_dg = volume_flux,56volume_flux_fv = surface_flux)57solver = DGSEM(basis, surface_flux, volume_integral)5859cells_per_dimension = (64,)60coordinates_min = (-2.0,)61coordinates_max = (2.0,)62mesh = StructuredMesh(cells_per_dimension, coordinates_min, coordinates_max,63periodicity = true)6465semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;66boundary_conditions = boundary_condition_periodic)6768###############################################################################69# ODE solvers, callbacks etc.7071tspan = (0.0, 12.5)72ode = semidiscretize(semi, tspan)7374summary_callback = SummaryCallback()7576analysis_interval = 100077analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7879alive_callback = AliveCallback(analysis_interval = analysis_interval)8081save_solution = SaveSolutionCallback(interval = 100,82save_initial_solution = true,83save_final_solution = true,84solution_variables = cons2prim)8586stepsize_callback = StepsizeCallback(cfl = 0.5)8788callbacks = CallbackSet(summary_callback,89analysis_callback,90alive_callback,91save_solution,92stepsize_callback)9394###############################################################################95# run the simulation9697sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);98dt = stepsize_callback(ode), # solve needs some value here but it will be overwritten by the stepsize_callback99ode_default_options()..., callback = callbacks);100101102