Path: blob/main/examples/tree_2d_dgsem/elixir_euler_shockcapturing_subcell.jl
2055 views
using Trixi12###############################################################################3# semidiscretization of the compressible Euler equations45equations = CompressibleEulerEquations2D(1.4)67"""8initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)910A medium blast wave (modified to lower density and higher pressure) taken from11- Sebastian Hennemann, Gregor J. Gassner (2020)12A provably entropy stable subcell shock capturing approach for high order split form DG13[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)14"""15function initial_condition_blast_wave(x, t, equations::CompressibleEulerEquations2D)16# Modified From Hennemann & Gassner JCP paper 2020 (Sec. 6.3) -> modified to lower density, higher pressure17# Set up polar coordinates18RealT = eltype(x)19inicenter = SVector(0, 0)20x_norm = x[1] - inicenter[1]21y_norm = x[2] - inicenter[2]22r = sqrt(x_norm^2 + y_norm^2)23phi = atan(y_norm, x_norm)24sin_phi, cos_phi = sincos(phi)2526# Calculate primitive variables "normal" medium blast wave27rho = r > 0.5f0 ? RealT(0.1) : RealT(0.2691) # rho = r > 0.5 ? 1 : 1.169128v1 = r > 0.5f0 ? zero(RealT) : RealT(0.1882) * cos_phi29v2 = r > 0.5f0 ? zero(RealT) : RealT(0.1882) * sin_phi30p = r > 0.5f0 ? RealT(1.0E-1) : RealT(1.245) # p = r > 0.5 ? 1.0E-3 : 1.2453132return prim2cons(SVector(rho, v1, v2, p), equations)33end34initial_condition = initial_condition_blast_wave3536# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of37# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.38# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.39# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.40# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.41# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the42# `StepsizeCallback` (CFL-Condition) and less diffusion.43surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)44volume_flux = flux_ranocha45basis = LobattoLegendreBasis(3)46limiter_idp = SubcellLimiterIDP(equations, basis;47positivity_variables_cons = ["rho"],48positivity_correction_factor = 0.5)49volume_integral = VolumeIntegralSubcellLimiting(limiter_idp;50volume_flux_dg = volume_flux,51volume_flux_fv = surface_flux)52solver = DGSEM(basis, surface_flux, volume_integral)5354coordinates_min = (-2.0, -2.0)55coordinates_max = (2.0, 2.0)56mesh = TreeMesh(coordinates_min, coordinates_max,57initial_refinement_level = 5,58n_cells_max = 100_000)5960semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver)6162###############################################################################63# ODE solvers, callbacks etc.6465tspan = (0.0, 1.0)66ode = semidiscretize(semi, tspan)6768summary_callback = SummaryCallback()6970analysis_interval = 10071analysis_callback = AnalysisCallback(semi, interval = analysis_interval)7273alive_callback = AliveCallback(analysis_interval = analysis_interval)7475save_solution = SaveSolutionCallback(dt = 0.1,76save_initial_solution = true,77save_final_solution = true,78solution_variables = cons2prim,79extra_node_variables = (:limiting_coefficient,))8081stepsize_callback = StepsizeCallback(cfl = 0.6)8283callbacks = CallbackSet(summary_callback,84analysis_callback, alive_callback,85save_solution,86stepsize_callback)8788###############################################################################89# run the simulation9091stage_callbacks = (SubcellLimiterIDPCorrection(), BoundsCheckCallback(save_errors = false))9293sol = Trixi.solve(ode, Trixi.SimpleSSPRK33(stage_callbacks = stage_callbacks);94dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback95ode_default_options()..., callback = callbacks);969798