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_supersonic_cylinder.jl
2055 views
1
# Channel flow around a cylinder at Mach 3
2
#
3
# Boundary conditions are supersonic Mach 3 inflow at the left portion of the domain
4
# and supersonic outflow at the right portion of the domain. The top and bottom of the
5
# channel as well as the cylinder are treated as Euler slip wall boundaries.
6
# This flow results in strong shock reflections / interactions as well as Kelvin-Helmholtz
7
# instabilities at later times as two Mach stems form above and below the cylinder.
8
#
9
# For complete details on the problem setup see Section 5.7 of the paper:
10
# - Jean-Luc Guermond, Murtazo Nazarov, Bojan Popov, and Ignacio Tomas (2018)
11
# Second-Order Invariant Domain Preserving Approximation of the Euler Equations using Convex Limiting.
12
# [DOI: 10.1137/17M1149961](https://doi.org/10.1137/17M1149961)
13
#
14
# Keywords: supersonic flow, shock capturing, AMR, unstructured curved mesh, positivity preservation, compressible Euler, 2D
15
16
using OrdinaryDiffEqSSPRK
17
using Trixi
18
19
###############################################################################
20
# semidiscretization of the compressible Euler equations
21
22
equations = CompressibleEulerEquations2D(1.4)
23
24
@inline function initial_condition_mach3_flow(x, t, equations::CompressibleEulerEquations2D)
25
# set the freestream flow parameters
26
rho_freestream = 1.4
27
v1 = 3.0
28
v2 = 0.0
29
p_freestream = 1.0
30
31
prim = SVector(rho_freestream, v1, v2, p_freestream)
32
return prim2cons(prim, equations)
33
end
34
35
initial_condition = initial_condition_mach3_flow
36
37
# Supersonic inflow boundary condition.
38
# Calculate the boundary flux entirely from the external solution state, i.e., set
39
# external solution state values for everything entering the domain.
40
@inline function boundary_condition_supersonic_inflow(u_inner,
41
normal_direction::AbstractVector,
42
x, t, surface_flux_function,
43
equations::CompressibleEulerEquations2D)
44
u_boundary = initial_condition_mach3_flow(x, t, equations)
45
return flux(u_boundary, normal_direction, equations)
46
end
47
48
# Supersonic outflow boundary condition.
49
# Calculate the boundary flux entirely from the internal solution state. Analogous to supersonic inflow
50
# except all the solution state values are set from the internal solution as everything leaves the domain
51
@inline function boundary_condition_outflow(u_inner, normal_direction::AbstractVector, x, t,
52
surface_flux_function,
53
equations::CompressibleEulerEquations2D)
54
return flux(u_inner, normal_direction, equations)
55
end
56
57
boundary_conditions = Dict(:Bottom => boundary_condition_slip_wall,
58
:Circle => boundary_condition_slip_wall,
59
:Top => boundary_condition_slip_wall,
60
:Right => boundary_condition_outflow,
61
:Left => boundary_condition_supersonic_inflow)
62
63
volume_flux = flux_ranocha_turbo
64
# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of
65
# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.
66
# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.
67
# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.
68
# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.
69
# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the
70
# `StepsizeCallback` (CFL-Condition) and less diffusion.
71
surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)
72
73
polydeg = 3
74
basis = LobattoLegendreBasis(polydeg)
75
shock_indicator = IndicatorHennemannGassner(equations, basis,
76
alpha_max = 0.5,
77
alpha_min = 0.001,
78
alpha_smooth = true,
79
variable = density_pressure)
80
volume_integral = VolumeIntegralShockCapturingHG(shock_indicator;
81
volume_flux_dg = volume_flux,
82
volume_flux_fv = surface_flux)
83
solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,
84
volume_integral = volume_integral)
85
86
# Get the unstructured quad mesh from a file (downloads the file if not available locally)
87
mesh_file = Trixi.download("https://gist.githubusercontent.com/andrewwinters5000/a08f78f6b185b63c3baeff911a63f628/raw/addac716ea0541f588b9d2bd3f92f643eb27b88f/abaqus_cylinder_in_channel.inp",
88
joinpath(@__DIR__, "abaqus_cylinder_in_channel.inp"))
89
90
mesh = P4estMesh{2}(mesh_file)
91
92
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,
93
boundary_conditions = boundary_conditions)
94
95
###############################################################################
96
# ODE solvers
97
98
tspan = (0.0, 2.0)
99
ode = semidiscretize(semi, tspan)
100
101
# Callbacks
102
103
summary_callback = SummaryCallback()
104
105
analysis_interval = 1000
106
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)
107
108
alive_callback = AliveCallback(analysis_interval = analysis_interval)
109
110
save_solution = SaveSolutionCallback(interval = 1000,
111
save_initial_solution = true,
112
save_final_solution = true,
113
solution_variables = cons2prim)
114
115
amr_indicator = IndicatorLöhner(semi, variable = Trixi.density)
116
117
amr_controller = ControllerThreeLevel(semi, amr_indicator,
118
base_level = 0,
119
med_level = 3, med_threshold = 0.05,
120
max_level = 5, max_threshold = 0.1)
121
122
amr_callback = AMRCallback(semi, amr_controller,
123
interval = 1,
124
adapt_initial_condition = true,
125
adapt_initial_condition_only_refine = true)
126
127
callbacks = CallbackSet(summary_callback,
128
analysis_callback, alive_callback,
129
save_solution,
130
amr_callback)
131
132
# positivity limiter necessary for this example with strong shocks. Very sensitive
133
# to the order of the limiter variables, pressure must come first.
134
stage_limiter! = PositivityPreservingLimiterZhangShu(thresholds = (5.0e-7, 1.0e-6),
135
variables = (pressure, Trixi.density))
136
137
###############################################################################
138
# run the simulation
139
sol = solve(ode, SSPRK43(stage_limiter!);
140
ode_default_options()..., callback = callbacks);
141
142