Path: blob/main/examples/p4est_2d_dgsem/elixir_euler_NACA6412airfoil_mach2.jl
2055 views
using Trixi1using OrdinaryDiffEqSSPRK23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations2D(1.4)78@inline function initial_condition_mach2_flow(x, t, equations::CompressibleEulerEquations2D)9# set the freestream flow parameters10rho_freestream = 1.411v1 = 2.012v2 = 0.013p_freestream = 1.01415prim = SVector(rho_freestream, v1, v2, p_freestream)16return prim2cons(prim, equations)17end1819initial_condition = initial_condition_mach2_flow2021# Supersonic inflow boundary condition.22# Calculate the boundary flux entirely from the external solution state, i.e., set23# external solution state values for everything entering the domain.24@inline function boundary_condition_supersonic_inflow(u_inner,25normal_direction::AbstractVector,26x, t, surface_flux_function,27equations::CompressibleEulerEquations2D)28u_boundary = initial_condition_mach2_flow(x, t, equations)2930return flux(u_boundary, normal_direction, equations)31end3233# Supersonic outflow boundary condition.34# Calculate the boundary flux entirely from the internal solution state. Analogous to supersonic inflow35# except all the solution state values are set from the internal solution as everything leaves the domain36@inline function boundary_condition_supersonic_outflow(u_inner,37normal_direction::AbstractVector, x,38t,39surface_flux_function,40equations::CompressibleEulerEquations2D)41return flux(u_inner, normal_direction, equations)42end4344polydeg = 34546# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of47# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.48# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.49# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.50# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.51# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the52# `StepsizeCallback` (CFL-Condition) and less diffusion.53surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)54volume_flux = flux_ranocha5556basis = LobattoLegendreBasis(polydeg)57shock_indicator = IndicatorHennemannGassner(equations, basis,58alpha_max = 0.5,59alpha_min = 0.001,60alpha_smooth = true,61variable = density_pressure)62volume_integral = VolumeIntegralShockCapturingHG(shock_indicator;63volume_flux_dg = volume_flux,64volume_flux_fv = surface_flux)6566# DG Solver67solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,68volume_integral = volume_integral)6970# Mesh generated from the following gmsh geometry input file:71# https://gist.githubusercontent.com/DanielDoehring/5ade6d93629f0d8c23a598812dbee2a9/raw/d2bc904fe92146eae1a36156e7f5c535dc1a80f1/NACA6412.geo72mesh_file = joinpath(@__DIR__, "mesh_NACA6412.inp")73isfile(mesh_file) ||74Trixi.download("https://gist.githubusercontent.com/DanielDoehring/e2a389f04f1e37b33819b9637e8ee4c3/raw/4bf7607a2ce4432fdb5cb87d5e264949b11bd5d7/mesh_NACA6412.inp",75mesh_file)7677boundary_symbols = [:PhysicalLine1, :PhysicalLine2, :PhysicalLine3, :PhysicalLine4]7879mesh = P4estMesh{2}(mesh_file, polydeg = polydeg, boundary_symbols = boundary_symbols)8081boundary_conditions = Dict(:PhysicalLine1 => boundary_condition_supersonic_inflow, # Left boundary82:PhysicalLine2 => boundary_condition_supersonic_outflow, # Right boundary83:PhysicalLine3 => boundary_condition_supersonic_outflow, # Top and bottom boundary84:PhysicalLine4 => boundary_condition_slip_wall) # Airfoil8586semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,87boundary_conditions = boundary_conditions)8889tspan = (0.0, 5.0)90ode = semidiscretize(semi, tspan)9192summary_callback = SummaryCallback()9394analysis_interval = 100095analysis_callback = AnalysisCallback(semi, interval = analysis_interval)9697stepsize_callback = StepsizeCallback(cfl = 4.0)9899callbacks = CallbackSet(summary_callback,100analysis_callback,101stepsize_callback)102103# Run the simulation104###############################################################################105sol = solve(ode, SSPRK104(; thread = Trixi.True());106dt = 1.0, # overwritten by the `stepsize_callback`107ode_default_options()...,108callback = callbacks);109110111