Path: blob/main/examples/p4est_2d_dgsem/elixir_euler_sedov.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations2D(1.4)78"""9initial_condition_sedov_blast_wave(x, t, equations::CompressibleEulerEquations2D)1011The Sedov blast wave setup based on Flash12- https://flash.rochester.edu/site/flashcode/user_support/flash_ug_devel/node187.html#SECTION01011400000000000000013"""14function initial_condition_sedov_blast_wave(x, t, equations::CompressibleEulerEquations2D)15# Set up polar coordinates16inicenter = SVector(0.0, 0.0)17x_norm = x[1] - inicenter[1]18y_norm = x[2] - inicenter[2]19r = sqrt(x_norm^2 + y_norm^2)2021# Setup based on https://flash.rochester.edu/site/flashcode/user_support/flash_ug_devel/node187.html#SECTION01011400000000000000022r0 = 0.21875 # = 3.5 * smallest dx (for domain length=4 and max-ref=6)23E = 1.024p0_inner = 3 * (equations.gamma - 1) * E / (3 * pi * r0^2)25p0_outer = 1.0e-5 # = true Sedov setup2627# Calculate primitive variables28rho = 1.029v1 = 0.030v2 = 0.031p = r > r0 ? p0_outer : p0_inner3233return prim2cons(SVector(rho, v1, v2, p), equations)34end3536initial_condition = initial_condition_sedov_blast_wave3738# Get the DG approximation space3940# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of41# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.42# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.43# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.44# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.45# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the46# `StepsizeCallback` (CFL-Condition) and less diffusion.47surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)48volume_flux = flux_ranocha49polydeg = 450basis = LobattoLegendreBasis(polydeg)51indicator_sc = IndicatorHennemannGassner(equations, basis,52alpha_max = 1.0,53alpha_min = 0.001,54alpha_smooth = true,55variable = density_pressure)56volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;57volume_flux_dg = volume_flux,58volume_flux_fv = surface_flux)5960solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,61volume_integral = volume_integral)6263###############################################################################6465coordinates_min = (-1.0, -1.0)66coordinates_max = (1.0, 1.0)6768trees_per_dimension = (4, 4)69mesh = P4estMesh(trees_per_dimension,70polydeg = 4, initial_refinement_level = 2,71coordinates_min = coordinates_min, coordinates_max = coordinates_max,72periodicity = true)7374semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)7576###############################################################################77# ODE solvers, callbacks etc.7879tspan = (0.0, 12.5)80ode = semidiscretize(semi, tspan)8182summary_callback = SummaryCallback()8384analysis_interval = 30085analysis_callback = AnalysisCallback(semi, interval = analysis_interval)8687alive_callback = AliveCallback(analysis_interval = analysis_interval)8889save_solution = SaveSolutionCallback(interval = 300,90save_initial_solution = true,91save_final_solution = true)9293stepsize_callback = StepsizeCallback(cfl = 0.5)9495callbacks = CallbackSet(summary_callback,96analysis_callback,97alive_callback,98save_solution,99stepsize_callback)100101###############################################################################102# run the simulation103104sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);105dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback106ode_default_options()..., callback = callbacks);107108109