Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/acoustic_perturbation_2d.jl
5586 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
@doc raw"""
9
AcousticPerturbationEquations2D(v_mean_global, c_mean_global, rho_mean_global)
10
11
Acoustic perturbation equations (APE) in two space dimensions. The equations are given by
12
```math
13
\begin{aligned}
14
\frac{\partial\mathbf{v'}}{\partial t} + \nabla (\bar{\mathbf{v}}\cdot\mathbf{v'})
15
+ \nabla\left( \frac{\bar{c}^2 \tilde{p}'}{\bar{\rho}} \right) &= 0 \\
16
\frac{\partial \tilde{p}'}{\partial t} +
17
\nabla\cdot (\bar{\rho} \mathbf{v'} + \bar{\mathbf{v}} \tilde{p}') &= 0.
18
\end{aligned}
19
```
20
The bar ``\bar{(\cdot)}`` indicates time-averaged quantities. The unknowns of the APE are the
21
perturbed velocities ``\mathbf{v'} = (v_1', v_2')^T`` and the scaled perturbed pressure
22
``\tilde{p}' = \frac{p'}{\bar{c}^2}``, where ``p'`` denotes the perturbed pressure and the
23
perturbed variables are defined by ``\phi' = \phi - \bar{\phi}``.
24
25
In addition to the unknowns, Trixi.jl currently stores the mean values in the state vector,
26
i.e. the state vector used internally is given by
27
```math
28
\mathbf{u} =
29
\begin{pmatrix}
30
v_1' \\ v_2' \\ \tilde{p}' \\ \bar{v}_1 \\ \bar{v}_2 \\ \bar{c} \\ \bar{\rho}
31
\end{pmatrix}.
32
```
33
This affects the implementation and use of these equations in various ways:
34
* The flux values corresponding to the mean values must be zero.
35
* The mean values have to be considered when defining initial conditions, boundary conditions or
36
source terms.
37
* [`AnalysisCallback`](@ref) analyzes these variables too.
38
* Trixi.jl's visualization tools will visualize the mean values by default.
39
40
The constructor accepts a 2-tuple `v_mean_global` and scalars `c_mean_global` and `rho_mean_global`
41
which can be used to make the definition of initial conditions for problems with constant mean flow
42
more flexible. These values are ignored if the mean values are defined internally in an initial
43
condition.
44
45
The equations are based on the APE-4 system introduced in the following paper:
46
- Roland Ewert and Wolfgang Schröder (2003)
47
Acoustic perturbation equations based on flow decomposition via source filtering
48
[DOI: 10.1016/S0021-9991(03)00168-2](https://doi.org/10.1016/S0021-9991(03)00168-2)
49
"""
50
struct AcousticPerturbationEquations2D{RealT <: Real} <:
51
AbstractAcousticPerturbationEquations{2, 7}
52
v_mean_global::SVector{2, RealT}
53
c_mean_global::RealT
54
rho_mean_global::RealT
55
end
56
57
function AcousticPerturbationEquations2D(v_mean_global::NTuple{2, <:Real},
58
c_mean_global::Real,
59
rho_mean_global::Real)
60
return AcousticPerturbationEquations2D(SVector(v_mean_global), c_mean_global,
61
rho_mean_global)
62
end
63
64
function AcousticPerturbationEquations2D(; v_mean_global::NTuple{2, <:Real},
65
c_mean_global::Real,
66
rho_mean_global::Real)
67
return AcousticPerturbationEquations2D(SVector(v_mean_global), c_mean_global,
68
rho_mean_global)
69
end
70
71
function varnames(::typeof(cons2cons), ::AcousticPerturbationEquations2D)
72
return ("v1_prime", "v2_prime", "p_prime_scaled",
73
"v1_mean", "v2_mean", "c_mean", "rho_mean")
74
end
75
function varnames(::typeof(cons2prim), ::AcousticPerturbationEquations2D)
76
return ("v1_prime", "v2_prime", "p_prime",
77
"v1_mean", "v2_mean", "c_mean", "rho_mean")
78
end
79
80
# Convenience functions for retrieving state variables and mean variables
81
function cons2state(u, equations::AcousticPerturbationEquations2D)
82
return SVector(u[1], u[2], u[3])
83
end
84
85
function cons2mean(u, equations::AcousticPerturbationEquations2D)
86
return SVector(u[4], u[5], u[6], u[7])
87
end
88
89
function varnames(::typeof(cons2state), ::AcousticPerturbationEquations2D)
90
return ("v1_prime", "v2_prime", "p_prime_scaled")
91
end
92
function varnames(::typeof(cons2mean), ::AcousticPerturbationEquations2D)
93
return ("v1_mean", "v2_mean", "c_mean", "rho_mean")
94
end
95
96
"""
97
global_mean_vars(equations::AcousticPerturbationEquations2D)
98
99
Returns the global mean variables stored in `equations`. This makes it easier
100
to define flexible initial conditions for problems with constant mean flow.
101
"""
102
function global_mean_vars(equations::AcousticPerturbationEquations2D)
103
return equations.v_mean_global[1], equations.v_mean_global[2],
104
equations.c_mean_global,
105
equations.rho_mean_global
106
end
107
108
"""
109
initial_condition_constant(x, t, equations::AcousticPerturbationEquations2D)
110
111
A constant initial condition where the state variables are zero and the mean flow is constant.
112
Uses the global mean values from `equations`.
113
"""
114
function initial_condition_constant(x, t, equations::AcousticPerturbationEquations2D)
115
v1_prime = 0
116
v2_prime = 0
117
p_prime_scaled = 0
118
119
return SVector(v1_prime, v2_prime, p_prime_scaled, global_mean_vars(equations)...)
120
end
121
122
"""
123
initial_condition_convergence_test(x, t, equations::AcousticPerturbationEquations2D)
124
125
A smooth initial condition used for convergence tests in combination with
126
[`source_terms_convergence_test`](@ref). Uses the global mean values from `equations`.
127
"""
128
function initial_condition_convergence_test(x, t,
129
equations::AcousticPerturbationEquations2D)
130
RealT = eltype(x)
131
a = 1
132
c = 2
133
L = 2
134
f = 2.0f0 / L
135
A = convert(RealT, 0.2)
136
omega = 2 * convert(RealT, pi) * f
137
init = c + A * sin(omega * (x[1] + x[2] - a * t))
138
139
v1_prime = init
140
v2_prime = init
141
p_prime = init^2
142
143
prim = SVector(v1_prime, v2_prime, p_prime, global_mean_vars(equations)...)
144
145
return prim2cons(prim, equations)
146
end
147
148
"""
149
source_terms_convergence_test(u, x, t, equations::AcousticPerturbationEquations2D)
150
151
Source terms used for convergence tests in combination with
152
[`initial_condition_convergence_test`](@ref).
153
154
References for the method of manufactured solutions (MMS):
155
- Kambiz Salari and Patrick Knupp (2000)
156
Code Verification by the Method of Manufactured Solutions
157
[DOI: 10.2172/759450](https://doi.org/10.2172/759450)
158
- Patrick J. Roache (2002)
159
Code Verification by the Method of Manufactured Solutions
160
[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
161
"""
162
function source_terms_convergence_test(u, x, t,
163
equations::AcousticPerturbationEquations2D)
164
v1_mean, v2_mean, c_mean, rho_mean = cons2mean(u, equations)
165
166
RealT = eltype(u)
167
a = 1
168
c = 2
169
L = 2
170
f = 2.0f0 / L
171
A = convert(RealT, 0.2)
172
omega = 2 * convert(RealT, pi) * f
173
174
si, co = sincos(omega * (x[1] + x[2] - a * t))
175
tmp = v1_mean + v2_mean - a
176
177
du1 = du2 = A * omega * co * (2 * c / rho_mean + tmp + 2 / rho_mean * A * si)
178
du3 = A * omega * co * (2 * c_mean^2 * rho_mean + 2 * c * tmp + 2 * A * tmp * si) /
179
c_mean^2
180
181
du4 = du5 = du6 = du7 = 0
182
183
return SVector(du1, du2, du3, du4, du5, du6, du7)
184
end
185
186
"""
187
initial_condition_gauss(x, t, equations::AcousticPerturbationEquations2D)
188
189
A Gaussian pulse in a constant mean flow. Uses the global mean values from `equations`.
190
"""
191
function initial_condition_gauss(x, t, equations::AcousticPerturbationEquations2D)
192
v1_prime = 0
193
v2_prime = 0
194
p_prime = exp(-4 * (x[1]^2 + x[2]^2))
195
196
prim = SVector(v1_prime, v2_prime, p_prime, global_mean_vars(equations)...)
197
198
return prim2cons(prim, equations)
199
end
200
201
"""
202
boundary_condition_wall(u_inner, orientation, direction, x, t, surface_flux_function,
203
equations::AcousticPerturbationEquations2D)
204
205
Boundary conditions for a solid wall.
206
"""
207
function boundary_condition_wall(u_inner, orientation, direction, x, t,
208
surface_flux_function,
209
equations::AcousticPerturbationEquations2D)
210
# Boundary state is equal to the inner state except for the perturbed velocity. For boundaries
211
# in the -x/+x direction, we multiply the perturbed velocity in the x direction by -1.
212
# Similarly, for boundaries in the -y/+y direction, we multiply the perturbed velocity in the
213
# y direction by -1
214
if direction in (1, 2) # x direction
215
u_boundary = SVector(-u_inner[1], u_inner[2], u_inner[3],
216
cons2mean(u_inner, equations)...)
217
else # y direction
218
u_boundary = SVector(u_inner[1], -u_inner[2], u_inner[3],
219
cons2mean(u_inner, equations)...)
220
end
221
222
# Calculate boundary flux
223
if iseven(direction) # u_inner is "left" of boundary, u_boundary is "right" of boundary
224
flux = surface_flux_function(u_inner, u_boundary, orientation, equations)
225
else # u_boundary is "left" of boundary, u_inner is "right" of boundary
226
flux = surface_flux_function(u_boundary, u_inner, orientation, equations)
227
end
228
229
return flux
230
end
231
232
"""
233
boundary_condition_slip_wall(u_inner, normal_direction, x, t, surface_flux_function,
234
equations::AcousticPerturbationEquations2D)
235
236
Use an orthogonal projection of the perturbed velocities to zero out the normal velocity
237
while retaining the possibility of a tangential velocity in the boundary state.
238
Further details are available in the paper:
239
- Marcus Bauer, Jürgen Dierke and Roland Ewert (2011)
240
Application of a discontinuous Galerkin method to discretize acoustic perturbation equations
241
[DOI: 10.2514/1.J050333](https://doi.org/10.2514/1.J050333)
242
"""
243
function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector, x, t,
244
surface_flux_function,
245
equations::AcousticPerturbationEquations2D)
246
# normalize the outward pointing direction
247
normal = normal_direction / norm(normal_direction)
248
249
# compute the normal perturbed velocity
250
u_normal = normal[1] * u_inner[1] + normal[2] * u_inner[2]
251
252
# create the "external" boundary solution state
253
u_boundary = SVector(u_inner[1] - 2 * u_normal * normal[1],
254
u_inner[2] - 2 * u_normal * normal[2],
255
u_inner[3], cons2mean(u_inner, equations)...)
256
257
# calculate the boundary flux
258
flux = surface_flux_function(u_inner, u_boundary, normal_direction, equations)
259
260
return flux
261
end
262
263
# Calculate 1D flux for a single point
264
@inline function flux(u, orientation::Integer,
265
equations::AcousticPerturbationEquations2D)
266
v1_prime, v2_prime, p_prime_scaled = cons2state(u, equations)
267
v1_mean, v2_mean, c_mean, rho_mean = cons2mean(u, equations)
268
269
# Calculate flux for conservative state variables
270
RealT = eltype(u)
271
if orientation == 1
272
f1 = v1_mean * v1_prime + v2_mean * v2_prime +
273
c_mean^2 * p_prime_scaled / rho_mean
274
f2 = zero(RealT)
275
f3 = rho_mean * v1_prime + v1_mean * p_prime_scaled
276
else
277
f1 = zero(RealT)
278
f2 = v1_mean * v1_prime + v2_mean * v2_prime +
279
c_mean^2 * p_prime_scaled / rho_mean
280
f3 = rho_mean * v2_prime + v2_mean * p_prime_scaled
281
end
282
283
# The rest of the state variables are actually variable coefficients, hence the flux should be
284
# zero. See https://github.com/trixi-framework/Trixi.jl/issues/358#issuecomment-784828762
285
# for details.
286
f4 = f5 = f6 = f7 = 0
287
288
return SVector(f1, f2, f3, f4, f5, f6, f7)
289
end
290
291
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
292
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
293
equations::AcousticPerturbationEquations2D)
294
# Calculate v = v_prime + v_mean
295
v_prime_ll = u_ll[orientation]
296
v_prime_rr = u_rr[orientation]
297
v_mean_ll = u_ll[orientation + 3]
298
v_mean_rr = u_rr[orientation + 3]
299
300
v_ll = v_prime_ll + v_mean_ll
301
v_rr = v_prime_rr + v_mean_rr
302
303
c_mean_ll = u_ll[6]
304
c_mean_rr = u_rr[6]
305
306
return max(abs(v_ll), abs(v_rr)) + max(c_mean_ll, c_mean_rr)
307
end
308
309
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
310
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
311
equations::AcousticPerturbationEquations2D)
312
# Calculate v = v_prime + v_mean
313
v_prime_ll = u_ll[orientation]
314
v_prime_rr = u_rr[orientation]
315
v_mean_ll = u_ll[orientation + 3]
316
v_mean_rr = u_rr[orientation + 3]
317
318
v_ll = v_prime_ll + v_mean_ll
319
v_rr = v_prime_rr + v_mean_rr
320
321
c_mean_ll = u_ll[6]
322
c_mean_rr = u_rr[6]
323
324
return max(abs(v_ll) + c_mean_ll, abs(v_rr) + c_mean_rr)
325
end
326
327
# Calculate 1D flux for a single point in the normal direction
328
# Note, this directional vector is not normalized
329
@inline function flux(u, normal_direction::AbstractVector,
330
equations::AcousticPerturbationEquations2D)
331
v1_prime, v2_prime, p_prime_scaled = cons2state(u, equations)
332
v1_mean, v2_mean, c_mean, rho_mean = cons2mean(u, equations)
333
334
f1 = normal_direction[1] * (v1_mean * v1_prime + v2_mean * v2_prime +
335
c_mean^2 * p_prime_scaled / rho_mean)
336
f2 = normal_direction[2] * (v1_mean * v1_prime + v2_mean * v2_prime +
337
c_mean^2 * p_prime_scaled / rho_mean)
338
f3 = (normal_direction[1] * (rho_mean * v1_prime + v1_mean * p_prime_scaled)
339
+
340
normal_direction[2] * (rho_mean * v2_prime + v2_mean * p_prime_scaled))
341
342
# The rest of the state variables are actually variable coefficients, hence the flux should be
343
# zero. See https://github.com/trixi-framework/Trixi.jl/issues/358#issuecomment-784828762
344
# for details.
345
f4 = f5 = f6 = f7 = 0
346
347
return SVector(f1, f2, f3, f4, f5, f6, f7)
348
end
349
350
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
351
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
352
equations::AcousticPerturbationEquations2D)
353
# Calculate v = v_prime + v_mean
354
v_prime_ll = normal_direction[1] * u_ll[1] + normal_direction[2] * u_ll[2]
355
v_prime_rr = normal_direction[1] * u_rr[1] + normal_direction[2] * u_rr[2]
356
v_mean_ll = normal_direction[1] * u_ll[4] + normal_direction[2] * u_ll[5]
357
v_mean_rr = normal_direction[1] * u_rr[4] + normal_direction[2] * u_rr[5]
358
359
v_ll = v_prime_ll + v_mean_ll
360
v_rr = v_prime_rr + v_mean_rr
361
362
c_mean_ll = u_ll[6]
363
c_mean_rr = u_rr[6]
364
365
# The v_normals are already scaled by the norm
366
return (max(abs(v_ll), abs(v_rr)) +
367
max(c_mean_ll, c_mean_rr) * norm(normal_direction))
368
end
369
370
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
371
@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,
372
equations::AcousticPerturbationEquations2D)
373
# Calculate v = v_prime + v_mean
374
v_prime_ll = normal_direction[1] * u_ll[1] + normal_direction[2] * u_ll[2]
375
v_prime_rr = normal_direction[1] * u_rr[1] + normal_direction[2] * u_rr[2]
376
v_mean_ll = normal_direction[1] * u_ll[4] + normal_direction[2] * u_ll[5]
377
v_mean_rr = normal_direction[1] * u_rr[4] + normal_direction[2] * u_rr[5]
378
379
v_ll = v_prime_ll + v_mean_ll
380
v_rr = v_prime_rr + v_mean_rr
381
382
c_mean_ll = u_ll[6]
383
c_mean_rr = u_rr[6]
384
385
norm_ = norm(normal_direction)
386
# The v_normals are already scaled by the norm
387
return max(abs(v_ll) + c_mean_ll * norm_, abs(v_rr) + c_mean_rr * norm_)
388
end
389
390
# Specialized `DissipationLocalLaxFriedrichs` to avoid spurious dissipation in the mean values
391
@inline function (dissipation::DissipationLocalLaxFriedrichs)(u_ll, u_rr,
392
orientation_or_normal_direction,
393
equations::AcousticPerturbationEquations2D)
394
λ = dissipation.max_abs_speed(u_ll, u_rr, orientation_or_normal_direction,
395
equations)
396
diss = -0.5f0 * λ * (u_rr - u_ll)
397
z = 0
398
399
return SVector(diss[1], diss[2], diss[3], z, z, z, z)
400
end
401
402
"""
403
have_constant_speed(::AcousticPerturbationEquations2D)
404
405
Indicates whether the characteristic speeds are constant, i.e., independent of the solution.
406
Queried in the timestep computation [`StepsizeCallback`](@ref).
407
408
The acoustic perturbation equations are in principle linear for **constant** mean flow fields.
409
However, the mean flow variables are part of the solution vector in
410
[`AcousticPerturbationEquations2D`](@ref) and only the **global** mean flow variables are constant,
411
similar to the [`LinearizedEulerEquations2D`](@ref).
412
413
Moreover, when coupling to the [`CompressibleEulerEquations2D`](@ref) equations via
414
[`SemidiscretizationEulerAcoustics`](@ref), the mean field variables are updated
415
on the fly, see [`EulerAcousticsCouplingCallback`](@ref).
416
417
# Returns
418
- `False()`
419
"""
420
@inline have_constant_speed(::AcousticPerturbationEquations2D) = False()
421
422
@inline function max_abs_speeds(u, equations::AcousticPerturbationEquations2D)
423
v1_mean = u[4]
424
v2_mean = u[5]
425
c_mean = u[6]
426
427
return abs(v1_mean) + c_mean, abs(v2_mean) + c_mean
428
end
429
430
# Convert conservative variables to primitive
431
@inline function cons2prim(u, equations::AcousticPerturbationEquations2D)
432
p_prime_scaled = u[3]
433
c_mean = u[6]
434
p_prime = p_prime_scaled * c_mean^2
435
436
return SVector(u[1], u[2], p_prime, u[4], u[5], u[6], u[7])
437
end
438
439
# Convert primitive variables to conservative
440
@inline function prim2cons(u, equations::AcousticPerturbationEquations2D)
441
p_prime = u[3]
442
c_mean = u[6]
443
p_prime_scaled = p_prime / c_mean^2
444
445
return SVector(u[1], u[2], p_prime_scaled, u[4], u[5], u[6], u[7])
446
end
447
448
# Convert conservative variables to entropy variables
449
@inline cons2entropy(u, equations::AcousticPerturbationEquations2D) = u
450
end # @muladd
451
452