Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/examples/tree_2d_dgsem/elixir_acoustics_monopole.jl
2055 views
1
using OrdinaryDiffEqLowStorageRK
2
using Trixi
3
4
###############################################################################
5
# semidiscretization of the acoustic perturbation equations
6
7
equations = AcousticPerturbationEquations2D(v_mean_global = (0.0, 0.0), c_mean_global = 0.0,
8
rho_mean_global = 0.0)
9
10
# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux
11
12
# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of
13
# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.
14
# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.
15
# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.
16
# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.
17
# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the
18
# `StepsizeCallback` (CFL-Condition) and less diffusion.
19
solver = DGSEM(polydeg = 3, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))
20
21
coordinates_min = (-20.6, 0.0) # minimum coordinates (min(x), min(y))
22
coordinates_max = (30.6, 51.2) # maximum coordinates (max(x), max(y))
23
24
"""
25
initial_condition_monopole(x, t, equations::AcousticPerturbationEquations2D)
26
27
Initial condition for the monopole in a boundary layer setup, used in combination with
28
[`boundary_condition_monopole`](@ref).
29
"""
30
function initial_condition_monopole(x, t, equations::AcousticPerturbationEquations2D)
31
RealT = eltype(x)
32
m = convert(RealT, 0.3) # Mach number
33
34
v1_prime = 0
35
v2_prime = 0
36
p_prime = 0
37
38
v1_mean = x[2] > 1 ? m : m * (2 * x[2] - 2 * x[2]^2 + x[2]^4)
39
v2_mean = 0
40
c_mean = 1
41
rho_mean = 1
42
43
prim = SVector(v1_prime, v2_prime, p_prime, v1_mean, v2_mean, c_mean, rho_mean)
44
45
return prim2cons(prim, equations)
46
end
47
initial_condition = initial_condition_monopole # does not use the global mean values given above
48
49
"""
50
boundary_condition_monopole(u_inner, orientation, direction, x, t, surface_flux_function,
51
equations::AcousticPerturbationEquations2D)
52
53
Boundary condition for a monopole in a boundary layer at the -y boundary, i.e. `direction = 3`.
54
This will return an error for any other direction. This boundary condition is used in combination
55
with [`initial_condition_monopole`](@ref).
56
"""
57
function boundary_condition_monopole(u_inner, orientation, direction, x, t,
58
surface_flux_function,
59
equations::AcousticPerturbationEquations2D)
60
RealT = eltype(u_inner)
61
if direction != 3
62
error("expected direction = 3, got $direction instead")
63
end
64
65
# Wall at the boundary in -y direction with a monopole at -0.05 <= x <= 0.05. In the monopole area
66
# we use a sinusoidal boundary state for the perturbed variables. For the rest of the -y boundary
67
# we set the boundary state to the inner state and multiply the perturbed velocity in the
68
# y-direction by -1.
69
if RealT(-0.05) <= x[1] <= RealT(0.05) # Monopole
70
v1_prime = 0
71
v2_prime = p_prime = sinpi(2 * t)
72
73
prim_boundary = SVector(v1_prime, v2_prime, p_prime, u_inner[4], u_inner[5],
74
u_inner[6], u_inner[7])
75
76
u_boundary = prim2cons(prim_boundary, equations)
77
else # Wall
78
u_boundary = SVector(u_inner[1], -u_inner[2], u_inner[3], u_inner[4], u_inner[5],
79
u_inner[6],
80
u_inner[7])
81
end
82
83
# Calculate boundary flux
84
flux = surface_flux_function(u_boundary, u_inner, orientation, equations)
85
86
return flux
87
end
88
89
"""
90
boundary_condition_zero(u_inner, orientation, direction, x, t, surface_flux_function,
91
equations::AcousticPerturbationEquations2D)
92
93
Boundary condition that uses a boundary state where the state variables are zero and the mean
94
variables are the same as in `u_inner`.
95
"""
96
function boundary_condition_zero(u_inner, orientation, direction, x, t,
97
surface_flux_function,
98
equations::AcousticPerturbationEquations2D)
99
value = zero(eltype(u_inner))
100
u_boundary = SVector(value, value, value, cons2mean(u_inner, equations)...)
101
102
# Calculate boundary flux
103
if iseven(direction) # u_inner is "left" of boundary, u_boundary is "right" of boundary
104
flux = surface_flux_function(u_inner, u_boundary, orientation, equations)
105
else # u_boundary is "left" of boundary, u_inner is "right" of boundary
106
flux = surface_flux_function(u_boundary, u_inner, orientation, equations)
107
end
108
109
return flux
110
end
111
112
boundary_conditions = (x_neg = boundary_condition_zero,
113
x_pos = boundary_condition_zero,
114
y_neg = boundary_condition_monopole,
115
y_pos = boundary_condition_zero)
116
117
# Create a uniformly refined mesh with periodic boundaries
118
mesh = TreeMesh(coordinates_min, coordinates_max,
119
initial_refinement_level = 6,
120
n_cells_max = 100_000,
121
periodicity = false)
122
123
# A semidiscretization collects data structures and functions for the spatial discretization
124
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,
125
boundary_conditions = boundary_conditions)
126
127
###############################################################################
128
# ODE solvers, callbacks etc.
129
130
# Create ODE problem with time span from 0.0 to 24.0
131
tspan = (0.0, 24.0)
132
ode = semidiscretize(semi, tspan)
133
134
# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup
135
# and resets the timers
136
summary_callback = SummaryCallback()
137
138
# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results
139
analysis_callback = AnalysisCallback(semi, interval = 100)
140
141
# The SaveSolutionCallback allows to save the solution to a file in regular intervals
142
save_solution = SaveSolutionCallback(interval = 100, solution_variables = cons2prim)
143
144
# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step
145
stepsize_callback = StepsizeCallback(cfl = 0.8)
146
147
# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver
148
callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,
149
stepsize_callback)
150
151
###############################################################################
152
# run the simulation
153
154
# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks
155
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
156
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
157
ode_default_options()..., callback = callbacks)
158
159