Path: blob/main/examples/t8code_3d_dgsem/elixir_euler_sedov.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations3D(1.4)78"""9initial_condition_medium_sedov_blast_wave(x, t, equations::CompressibleEulerEquations3D)1011The Sedov blast wave setup based on Flash12- https://flash.rochester.edu/site/flashcode/user_support/flash_ug_devel/node187.html#SECTION01011400000000000000013with smaller strength of the initial discontinuity.14"""15function initial_condition_medium_sedov_blast_wave(x, t,16equations::CompressibleEulerEquations3D)17# Set up polar coordinates18inicenter = SVector(0.0, 0.0, 0.0)19x_norm = x[1] - inicenter[1]20y_norm = x[2] - inicenter[2]21z_norm = x[3] - inicenter[3]22r = sqrt(x_norm^2 + y_norm^2 + z_norm^2)2324# Setup based on https://flash.rochester.edu/site/flashcode/user_support/flash_ug_devel/node187.html#SECTION01011400000000000000025r0 = 0.21875 # = 3.5 * smallest dx (for domain length=4 and max-ref=6)26E = 1.027p0_inner = 3 * (equations.gamma - 1) * E / (4 * pi * r0^2)28p0_outer = 1.0e-32930# Calculate primitive variables31rho = 1.032v1 = 0.033v2 = 0.034v3 = 0.035p = r > r0 ? p0_outer : p0_inner3637return prim2cons(SVector(rho, v1, v2, v3, p), equations)38end3940initial_condition = initial_condition_medium_sedov_blast_wave4142# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of43# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.44# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.45# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.46# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.47# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the48# `StepsizeCallback` (CFL-Condition) and less diffusion.49surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)50volume_flux = flux_ranocha51polydeg = 552basis = LobattoLegendreBasis(polydeg)53indicator_sc = IndicatorHennemannGassner(equations, basis,54alpha_max = 1.0,55alpha_min = 0.001,56alpha_smooth = true,57variable = density_pressure)58volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;59volume_flux_dg = volume_flux,60volume_flux_fv = surface_flux)6162solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,63volume_integral = volume_integral)6465coordinates_min = (-1.0, -1.0, -1.0)66coordinates_max = (1.0, 1.0, 1.0)6768trees_per_dimension = (4, 4, 4)69mesh = T8codeMesh(trees_per_dimension,70polydeg = 4, initial_refinement_level = 0,71coordinates_min = coordinates_min, coordinates_max = coordinates_max,72periodicity = true)7374# create the semi discretization object75semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)7677###############################################################################78# ODE solvers, callbacks etc.7980tspan = (0.0, 12.5)81ode = semidiscretize(semi, tspan)8283summary_callback = SummaryCallback()8485analysis_interval = 10086analysis_callback = AnalysisCallback(semi, interval = analysis_interval)8788alive_callback = AliveCallback(analysis_interval = analysis_interval)8990save_solution = SaveSolutionCallback(interval = 100,91save_initial_solution = true,92save_final_solution = true)9394stepsize_callback = StepsizeCallback(cfl = 0.5)9596callbacks = CallbackSet(summary_callback,97analysis_callback,98alive_callback,99save_solution,100stepsize_callback)101102###############################################################################103# run the simulation104105sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);106dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback107ode_default_options()..., callback = callbacks);108109110