Path: blob/main/examples/p4est_2d_dgsem/elixir_euler_supersonic_cylinder.jl
2055 views
# Channel flow around a cylinder at Mach 31#2# Boundary conditions are supersonic Mach 3 inflow at the left portion of the domain3# and supersonic outflow at the right portion of the domain. The top and bottom of the4# channel as well as the cylinder are treated as Euler slip wall boundaries.5# This flow results in strong shock reflections / interactions as well as Kelvin-Helmholtz6# instabilities at later times as two Mach stems form above and below the cylinder.7#8# For complete details on the problem setup see Section 5.7 of the paper:9# - Jean-Luc Guermond, Murtazo Nazarov, Bojan Popov, and Ignacio Tomas (2018)10# Second-Order Invariant Domain Preserving Approximation of the Euler Equations using Convex Limiting.11# [DOI: 10.1137/17M1149961](https://doi.org/10.1137/17M1149961)12#13# Keywords: supersonic flow, shock capturing, AMR, unstructured curved mesh, positivity preservation, compressible Euler, 2D1415using OrdinaryDiffEqSSPRK16using Trixi1718###############################################################################19# semidiscretization of the compressible Euler equations2021equations = CompressibleEulerEquations2D(1.4)2223@inline function initial_condition_mach3_flow(x, t, equations::CompressibleEulerEquations2D)24# set the freestream flow parameters25rho_freestream = 1.426v1 = 3.027v2 = 0.028p_freestream = 1.02930prim = SVector(rho_freestream, v1, v2, p_freestream)31return prim2cons(prim, equations)32end3334initial_condition = initial_condition_mach3_flow3536# Supersonic inflow boundary condition.37# Calculate the boundary flux entirely from the external solution state, i.e., set38# external solution state values for everything entering the domain.39@inline function boundary_condition_supersonic_inflow(u_inner,40normal_direction::AbstractVector,41x, t, surface_flux_function,42equations::CompressibleEulerEquations2D)43u_boundary = initial_condition_mach3_flow(x, t, equations)44return flux(u_boundary, normal_direction, equations)45end4647# Supersonic outflow boundary condition.48# Calculate the boundary flux entirely from the internal solution state. Analogous to supersonic inflow49# except all the solution state values are set from the internal solution as everything leaves the domain50@inline function boundary_condition_outflow(u_inner, normal_direction::AbstractVector, x, t,51surface_flux_function,52equations::CompressibleEulerEquations2D)53return flux(u_inner, normal_direction, equations)54end5556boundary_conditions = Dict(:Bottom => boundary_condition_slip_wall,57:Circle => boundary_condition_slip_wall,58:Top => boundary_condition_slip_wall,59:Right => boundary_condition_outflow,60:Left => boundary_condition_supersonic_inflow)6162volume_flux = flux_ranocha_turbo63# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of64# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.65# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.66# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.67# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.68# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the69# `StepsizeCallback` (CFL-Condition) and less diffusion.70surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)7172polydeg = 373basis = LobattoLegendreBasis(polydeg)74shock_indicator = IndicatorHennemannGassner(equations, basis,75alpha_max = 0.5,76alpha_min = 0.001,77alpha_smooth = true,78variable = density_pressure)79volume_integral = VolumeIntegralShockCapturingHG(shock_indicator;80volume_flux_dg = volume_flux,81volume_flux_fv = surface_flux)82solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,83volume_integral = volume_integral)8485# Get the unstructured quad mesh from a file (downloads the file if not available locally)86mesh_file = Trixi.download("https://gist.githubusercontent.com/andrewwinters5000/a08f78f6b185b63c3baeff911a63f628/raw/addac716ea0541f588b9d2bd3f92f643eb27b88f/abaqus_cylinder_in_channel.inp",87joinpath(@__DIR__, "abaqus_cylinder_in_channel.inp"))8889mesh = P4estMesh{2}(mesh_file)9091semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,92boundary_conditions = boundary_conditions)9394###############################################################################95# ODE solvers9697tspan = (0.0, 2.0)98ode = semidiscretize(semi, tspan)99100# Callbacks101102summary_callback = SummaryCallback()103104analysis_interval = 1000105analysis_callback = AnalysisCallback(semi, interval = analysis_interval)106107alive_callback = AliveCallback(analysis_interval = analysis_interval)108109save_solution = SaveSolutionCallback(interval = 1000,110save_initial_solution = true,111save_final_solution = true,112solution_variables = cons2prim)113114amr_indicator = IndicatorLöhner(semi, variable = Trixi.density)115116amr_controller = ControllerThreeLevel(semi, amr_indicator,117base_level = 0,118med_level = 3, med_threshold = 0.05,119max_level = 5, max_threshold = 0.1)120121amr_callback = AMRCallback(semi, amr_controller,122interval = 1,123adapt_initial_condition = true,124adapt_initial_condition_only_refine = true)125126callbacks = CallbackSet(summary_callback,127analysis_callback, alive_callback,128save_solution,129amr_callback)130131# positivity limiter necessary for this example with strong shocks. Very sensitive132# to the order of the limiter variables, pressure must come first.133stage_limiter! = PositivityPreservingLimiterZhangShu(thresholds = (5.0e-7, 1.0e-6),134variables = (pressure, Trixi.density))135136###############################################################################137# run the simulation138sol = solve(ode, SSPRK43(stage_limiter!);139ode_default_options()..., callback = callbacks);140141142