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_navierstokes_NACA0012airfoil_mach08.jl
2835 views
1
using OrdinaryDiffEqLowStorageRK
2
using Trixi
3
4
###############################################################################
5
# semidiscretization of the compressible Euler equations
6
7
# Laminar transonic flow around a NACA0012 airfoil.
8
9
# This test is taken from the paper of Swanson and Langer. The values for the drag and lift coefficients
10
# from Case 5 in Table 3 are used to validate the scheme and computation of surface forces.
11
12
# References:
13
# - Roy Charles Swanson, Stefan Langer (2016)
14
# Structured and Unstructured Grid Methods (2016)
15
# [https://ntrs.nasa.gov/citations/20160003623] (https://ntrs.nasa.gov/citations/20160003623)
16
# - Deep Ray, Praveen Chandrashekar (2017)
17
# An entropy stable finite volume scheme for the
18
# two dimensional Navier–Stokes equations on triangular grids
19
# [DOI:10.1016/j.amc.2017.07.020](https://doi.org/10.1016/j.amc.2017.07.020)
20
21
equations = CompressibleEulerEquations2D(1.4)
22
23
prandtl_number() = 0.72
24
mu() = 0.0031959974968701088
25
equations_parabolic = CompressibleNavierStokesDiffusion2D(equations, mu = mu(),
26
Prandtl = prandtl_number())
27
28
rho_inf() = 1.0
29
p_inf() = 2.85
30
aoa() = 10.0 * pi / 180.0 # 10 degree angle of attack
31
l_inf() = 1.0
32
mach_inf() = 0.8
33
u_inf(equations) = mach_inf() * sqrt(equations.gamma * p_inf() / rho_inf())
34
@inline function initial_condition_mach08_flow(x, t,
35
equations::CompressibleEulerEquations2D)
36
# set the freestream flow parameters
37
gamma = equations.gamma
38
u_inf = mach_inf() * sqrt(gamma * p_inf() / rho_inf())
39
40
v1 = u_inf * cos(aoa())
41
v2 = u_inf * sin(aoa())
42
43
prim = SVector(rho_inf(), v1, v2, p_inf())
44
return prim2cons(prim, equations)
45
end
46
47
initial_condition = initial_condition_mach08_flow
48
49
# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of
50
# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.
51
# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.
52
# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.
53
# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.
54
# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the
55
# `StepsizeCallback` (CFL-Condition) and less diffusion.
56
surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)
57
58
polydeg = 3
59
solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux)
60
61
mesh_file = Trixi.download("https://gist.githubusercontent.com/Arpit-Babbar/339662b4b46164a016e35c81c66383bb/raw/8bf94f5b426ba907ace87405cfcc1dcc2ef7cbda/NACA0012.inp",
62
joinpath(@__DIR__, "NACA0012.inp"))
63
64
mesh = P4estMesh{2}(mesh_file, initial_refinement_level = 1)
65
66
# The boundary values across outer boundary are constant but subsonic, so we cannot compute the
67
# boundary flux from the external information alone. Thus, we use the numerical flux to distinguish
68
# between inflow and outflow characteristics
69
@inline function boundary_condition_subsonic_constant(u_inner,
70
normal_direction::AbstractVector, x,
71
t,
72
surface_flux_function,
73
equations::CompressibleEulerEquations2D)
74
u_boundary = initial_condition_mach08_flow(x, t, equations)
75
76
return flux_hll(u_inner, u_boundary, normal_direction, equations)
77
end
78
79
boundary_conditions = (; Left = boundary_condition_subsonic_constant,
80
Right = boundary_condition_subsonic_constant,
81
Top = boundary_condition_subsonic_constant,
82
Bottom = boundary_condition_subsonic_constant,
83
AirfoilBottom = boundary_condition_slip_wall,
84
AirfoilTop = boundary_condition_slip_wall)
85
86
velocity_airfoil = NoSlip((x, t, equations_parabolic) -> SVector(0.0, 0.0))
87
88
heat_airfoil = Adiabatic((x, t, equations_parabolic) -> 0.0)
89
90
boundary_conditions_airfoil = BoundaryConditionNavierStokesWall(velocity_airfoil,
91
heat_airfoil)
92
93
function velocities_initial_condition_mach08_flow(x, t,
94
equations::CompressibleEulerEquations2D)
95
u_cons = initial_condition_mach08_flow(x, t, equations)
96
return SVector(u_cons[2] / u_cons[1], u_cons[3] / u_cons[1])
97
end
98
99
velocity_bc_square = NoSlip((x, t, equations_parabolic) -> velocities_initial_condition_mach08_flow(x,
100
t,
101
equations_parabolic.equations_hyperbolic))
102
103
heat_bc_square = Adiabatic((x, t, equations_parabolic) -> 0.0)
104
boundary_condition_square = BoundaryConditionNavierStokesWall(velocity_bc_square,
105
heat_bc_square)
106
107
boundary_conditions_parabolic = (; Left = boundary_condition_square,
108
Right = boundary_condition_square,
109
Top = boundary_condition_square,
110
Bottom = boundary_condition_square,
111
AirfoilBottom = boundary_conditions_airfoil,
112
AirfoilTop = boundary_conditions_airfoil)
113
114
semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic),
115
initial_condition, solver;
116
boundary_conditions = (boundary_conditions,
117
boundary_conditions_parabolic))
118
119
###############################################################################
120
# ODE solvers
121
122
# Run for a long time to reach a state where forces stabilize up to 3 digits
123
tspan = (0.0, 10.0)
124
ode = semidiscretize(semi, tspan)
125
126
# Callbacks
127
128
summary_callback = SummaryCallback()
129
130
analysis_interval = 2000
131
132
force_boundary_names = (:AirfoilBottom, :AirfoilTop)
133
drag_coefficient = AnalysisSurfaceIntegral(force_boundary_names,
134
DragCoefficientPressure2D(aoa(), rho_inf(),
135
u_inf(equations),
136
l_inf()))
137
138
lift_coefficient = AnalysisSurfaceIntegral(force_boundary_names,
139
LiftCoefficientPressure2D(aoa(), rho_inf(),
140
u_inf(equations),
141
l_inf()))
142
143
drag_coefficient_shear_force = AnalysisSurfaceIntegral(force_boundary_names,
144
DragCoefficientShearStress2D(aoa(),
145
rho_inf(),
146
u_inf(equations),
147
l_inf()))
148
149
lift_coefficient_shear_force = AnalysisSurfaceIntegral(force_boundary_names,
150
LiftCoefficientShearStress2D(aoa(),
151
rho_inf(),
152
u_inf(equations),
153
l_inf()))
154
155
analysis_callback = AnalysisCallback(semi, interval = analysis_interval,
156
output_directory = "out",
157
save_analysis = true,
158
analysis_errors = Symbol[],
159
analysis_integrals = (drag_coefficient,
160
lift_coefficient,
161
drag_coefficient_shear_force,
162
lift_coefficient_shear_force))
163
164
alive_callback = AliveCallback(analysis_interval = analysis_interval)
165
166
save_solution = SaveSolutionCallback(interval = 500,
167
save_initial_solution = true,
168
save_final_solution = true,
169
solution_variables = cons2prim)
170
171
callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback, save_solution)
172
173
###############################################################################
174
# run the simulation
175
176
sol = solve(ode, RDPK3SpFSAL49(thread = Trixi.True());
177
abstol = 1e-8, reltol = 1e-8,
178
ode_default_options()..., callback = callbacks)
179
180