Path: blob/main/examples/tree_2d_dgsem/elixir_euler_blast_wave.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations2D(1.4)78"""9initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)1011A medium blast wave taken from12- Sebastian Hennemann, Gregor J. Gassner (2020)13A provably entropy stable subcell shock capturing approach for high order split form DG14[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)15"""16function initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)17# Modified From Hennemann & Gassner JCP paper 2020 (Sec. 6.3) -> "medium blast wave"18# Set up polar coordinates19RealT = eltype(x)20inicenter = SVector(0, 0)21x_norm = x[1] - inicenter[1]22y_norm = x[2] - inicenter[2]23r = sqrt(x_norm^2 + y_norm^2)24phi = atan(y_norm, x_norm)25sin_phi, cos_phi = sincos(phi)2627# Calculate primitive variables28rho = r > 0.5f0 ? one(RealT) : RealT(1.1691)29v1 = r > 0.5f0 ? zero(RealT) : RealT(0.1882) * cos_phi30v2 = r > 0.5f0 ? zero(RealT) : RealT(0.1882) * sin_phi31p = r > 0.5f0 ? RealT(1.0E-3) : RealT(1.245)3233return prim2cons(SVector(rho, v1, v2, p), equations)34end35initial_condition = initial_condition_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_ranocha46basis = LobattoLegendreBasis(3)47indicator_sc = IndicatorHennemannGassner(equations, basis,48alpha_max = 0.5,49alpha_min = 0.001,50alpha_smooth = true,51variable = density_pressure)52volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;53volume_flux_dg = volume_flux,54volume_flux_fv = surface_flux)55solver = DGSEM(basis, surface_flux, volume_integral)5657coordinates_min = (-2.0, -2.0)58coordinates_max = (2.0, 2.0)59mesh = TreeMesh(coordinates_min, coordinates_max,60initial_refinement_level = 6,61n_cells_max = 10_000)6263semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)6465###############################################################################66# ODE solvers, callbacks etc.6768tspan = (0.0, 12.5)69ode = semidiscretize(semi, tspan)7071summary_callback = SummaryCallback()7273analysis_interval = 10074analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7576alive_callback = AliveCallback(analysis_interval = analysis_interval)7778save_solution = SaveSolutionCallback(interval = 100,79save_initial_solution = true,80save_final_solution = true,81solution_variables = cons2prim)8283stepsize_callback = StepsizeCallback(cfl = 0.9)8485callbacks = CallbackSet(summary_callback,86analysis_callback, alive_callback,87save_solution,88stepsize_callback)8990###############################################################################91# run the simulation9293sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);94dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback95ode_default_options()..., callback = callbacks);969798