Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/examples/p4est_2d_dgsem/elixir_euler_NACA6412airfoil_mach2.jl
2055 views
1
using Trixi
2
using OrdinaryDiffEqSSPRK
3
4
###############################################################################
5
# semidiscretization of the compressible Euler equations
6
7
equations = CompressibleEulerEquations2D(1.4)
8
9
@inline function initial_condition_mach2_flow(x, t, equations::CompressibleEulerEquations2D)
10
# set the freestream flow parameters
11
rho_freestream = 1.4
12
v1 = 2.0
13
v2 = 0.0
14
p_freestream = 1.0
15
16
prim = SVector(rho_freestream, v1, v2, p_freestream)
17
return prim2cons(prim, equations)
18
end
19
20
initial_condition = initial_condition_mach2_flow
21
22
# Supersonic inflow boundary condition.
23
# Calculate the boundary flux entirely from the external solution state, i.e., set
24
# external solution state values for everything entering the domain.
25
@inline function boundary_condition_supersonic_inflow(u_inner,
26
normal_direction::AbstractVector,
27
x, t, surface_flux_function,
28
equations::CompressibleEulerEquations2D)
29
u_boundary = initial_condition_mach2_flow(x, t, equations)
30
31
return flux(u_boundary, normal_direction, equations)
32
end
33
34
# Supersonic outflow boundary condition.
35
# Calculate the boundary flux entirely from the internal solution state. Analogous to supersonic inflow
36
# except all the solution state values are set from the internal solution as everything leaves the domain
37
@inline function boundary_condition_supersonic_outflow(u_inner,
38
normal_direction::AbstractVector, x,
39
t,
40
surface_flux_function,
41
equations::CompressibleEulerEquations2D)
42
return flux(u_inner, normal_direction, equations)
43
end
44
45
polydeg = 3
46
47
# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of
48
# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.
49
# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.
50
# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.
51
# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.
52
# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the
53
# `StepsizeCallback` (CFL-Condition) and less diffusion.
54
surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)
55
volume_flux = flux_ranocha
56
57
basis = LobattoLegendreBasis(polydeg)
58
shock_indicator = IndicatorHennemannGassner(equations, basis,
59
alpha_max = 0.5,
60
alpha_min = 0.001,
61
alpha_smooth = true,
62
variable = density_pressure)
63
volume_integral = VolumeIntegralShockCapturingHG(shock_indicator;
64
volume_flux_dg = volume_flux,
65
volume_flux_fv = surface_flux)
66
67
# DG Solver
68
solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,
69
volume_integral = volume_integral)
70
71
# Mesh generated from the following gmsh geometry input file:
72
# https://gist.githubusercontent.com/DanielDoehring/5ade6d93629f0d8c23a598812dbee2a9/raw/d2bc904fe92146eae1a36156e7f5c535dc1a80f1/NACA6412.geo
73
mesh_file = joinpath(@__DIR__, "mesh_NACA6412.inp")
74
isfile(mesh_file) ||
75
Trixi.download("https://gist.githubusercontent.com/DanielDoehring/e2a389f04f1e37b33819b9637e8ee4c3/raw/4bf7607a2ce4432fdb5cb87d5e264949b11bd5d7/mesh_NACA6412.inp",
76
mesh_file)
77
78
boundary_symbols = [:PhysicalLine1, :PhysicalLine2, :PhysicalLine3, :PhysicalLine4]
79
80
mesh = P4estMesh{2}(mesh_file, polydeg = polydeg, boundary_symbols = boundary_symbols)
81
82
boundary_conditions = Dict(:PhysicalLine1 => boundary_condition_supersonic_inflow, # Left boundary
83
:PhysicalLine2 => boundary_condition_supersonic_outflow, # Right boundary
84
:PhysicalLine3 => boundary_condition_supersonic_outflow, # Top and bottom boundary
85
:PhysicalLine4 => boundary_condition_slip_wall) # Airfoil
86
87
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,
88
boundary_conditions = boundary_conditions)
89
90
tspan = (0.0, 5.0)
91
ode = semidiscretize(semi, tspan)
92
93
summary_callback = SummaryCallback()
94
95
analysis_interval = 1000
96
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)
97
98
stepsize_callback = StepsizeCallback(cfl = 4.0)
99
100
callbacks = CallbackSet(summary_callback,
101
analysis_callback,
102
stepsize_callback)
103
104
# Run the simulation
105
###############################################################################
106
sol = solve(ode, SSPRK104(; thread = Trixi.True());
107
dt = 1.0, # overwritten by the `stepsize_callback`
108
ode_default_options()...,
109
callback = callbacks);
110
111