Path: blob/main/examples/tree_1d_dgsem/elixir_euler_blast_wave.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations1D(1.4)78"""9initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations1D)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::CompressibleEulerEquations1D)17# Modified From Hennemann & Gassner JCP paper 2020 (Sec. 6.3) -> "medium blast wave"18# Set up polar coordinates19RealT = eltype(x)20inicenter = SVector(0)21x_norm = x[1] - inicenter[1]22r = abs(x_norm)23# The following code is equivalent to24# phi = atan(0.0, x_norm)25# cos_phi = cos(phi)26# in 1D but faster27cos_phi = x_norm > 0 ? one(x_norm) : -one(x_norm)2829# Calculate primitive variables30rho = r > 0.5f0 ? one(RealT) : RealT(1.1691)31v1 = r > 0.5f0 ? zero(RealT) : RealT(0.1882) * cos_phi32p = r > 0.5f0 ? RealT(1.0E-3) : RealT(1.245)3334return prim2cons(SVector(rho, v1, p), equations)35end36initial_condition = initial_condition_blast_wave3738# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of39# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.40# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.41# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.42# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.43# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the44# `StepsizeCallback` (CFL-Condition) and less diffusion.45surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)46volume_flux = flux_ranocha47basis = LobattoLegendreBasis(3)48indicator_sc = IndicatorHennemannGassner(equations, basis,49alpha_max = 0.5,50alpha_min = 0.001,51alpha_smooth = true,52variable = density_pressure)53volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;54volume_flux_dg = volume_flux,55volume_flux_fv = surface_flux)56solver = DGSEM(basis, surface_flux, volume_integral)5758coordinates_min = (-2.0,)59coordinates_max = (2.0,)60mesh = TreeMesh(coordinates_min, coordinates_max,61initial_refinement_level = 6,62n_cells_max = 10_000)6364semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)6566###############################################################################67# ODE solvers, callbacks etc.6869tspan = (0.0, 12.5)70ode = semidiscretize(semi, tspan)7172summary_callback = SummaryCallback()7374analysis_interval = 1007576analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7778alive_callback = AliveCallback(analysis_interval = analysis_interval)7980save_solution = SaveSolutionCallback(interval = 100,81save_initial_solution = true,82save_final_solution = true,83solution_variables = cons2prim)8485stepsize_callback = StepsizeCallback(cfl = 0.5)8687callbacks = CallbackSet(summary_callback,88analysis_callback, alive_callback,89save_solution,90stepsize_callback)9192###############################################################################93# run the simulation9495sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);96dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback97ode_default_options()..., callback = callbacks);9899100