Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/compressible_euler_quasi_1d.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
CompressibleEulerEquationsQuasi1D(gamma)
10
11
The quasi-1d compressible Euler equations (see Chan et al. [DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089) for details)
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
a \rho \\ a \rho v_1 \\ a e_{\text{total}}
16
\end{pmatrix}
17
+
18
\frac{\partial}{\partial x}
19
\begin{pmatrix}
20
a \rho v_1 \\ a \rho v_1^2 \\ a v_1 (e_{\text{total}} +p)
21
\end{pmatrix}
22
+
23
a \frac{\partial}{\partial x}
24
\begin{pmatrix}
25
0 \\ p \\ 0
26
\end{pmatrix}
27
=
28
\begin{pmatrix}
29
0 \\ 0 \\ 0
30
\end{pmatrix}
31
```
32
for an ideal gas with ratio of specific heats `gamma` in one space dimension.
33
Here, ``\rho`` is the density, ``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy,
34
``a`` the (possibly) variable nozzle width, and
35
```math
36
p = (\gamma - 1) \left( e_{\text{total}} - \frac{1}{2} \rho v_1^2 \right)
37
```
38
the pressure.
39
40
The nozzle width function ``a(x)`` is set inside the initial condition routine
41
for a particular problem setup. To test the conservative form of the compressible Euler equations one can set the
42
nozzle width variable ``a`` to one.
43
44
In addition to the unknowns, Trixi.jl currently stores the nozzle width values at the approximation points
45
despite being fixed in time.
46
This affects the implementation and use of these equations in various ways:
47
* The flux values corresponding to the nozzle width must be zero.
48
* The nozzle width values must be included when defining initial conditions, boundary conditions or
49
source terms.
50
* [`AnalysisCallback`](@ref) analyzes this variable.
51
* Trixi.jl's visualization tools will visualize the nozzle width by default.
52
"""
53
struct CompressibleEulerEquationsQuasi1D{RealT <: Real} <:
54
AbstractCompressibleEulerEquations{1, 4}
55
gamma::RealT # ratio of specific heats
56
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
57
58
function CompressibleEulerEquationsQuasi1D(gamma)
59
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
60
return new{typeof(γ)}(γ, inv_gamma_minus_one)
61
end
62
end
63
64
have_nonconservative_terms(::CompressibleEulerEquationsQuasi1D) = True()
65
function varnames(::typeof(cons2cons), ::CompressibleEulerEquationsQuasi1D)
66
return ("a_rho", "a_rho_v1", "a_e_total", "a")
67
end
68
function varnames(::typeof(cons2prim), ::CompressibleEulerEquationsQuasi1D)
69
return ("rho", "v1", "p", "a")
70
end
71
72
"""
73
initial_condition_convergence_test(x, t, equations::CompressibleEulerEquationsQuasi1D)
74
75
A smooth initial condition used for convergence tests in combination with
76
[`source_terms_convergence_test`](@ref)
77
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
78
"""
79
function initial_condition_convergence_test(x, t,
80
equations::CompressibleEulerEquationsQuasi1D)
81
RealT = eltype(x)
82
c = 2
83
A = convert(RealT, 0.1)
84
L = 2
85
f = 1.0f0 / L
86
ω = 2 * convert(RealT, pi) * f
87
ini = c + A * sin(ω * (x[1] - t))
88
89
rho = ini
90
v1 = 1
91
e_total = ini^2 / rho
92
p = (equations.gamma - 1) * (e_total - 0.5f0 * rho * v1^2)
93
a = 1.5f0 - 0.5f0 * cospi(x[1])
94
95
return prim2cons(SVector(rho, v1, p, a), equations)
96
end
97
98
"""
99
source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquationsQuasi1D)
100
101
Source terms used for convergence tests in combination with
102
[`initial_condition_convergence_test`](@ref)
103
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
104
105
This manufactured solution source term is specifically designed for the nozzle width
106
```math
107
a(x) = 1.5 - 0.5 \\cos(x \\pi)
108
```
109
as defined in [`initial_condition_convergence_test`](@ref).
110
111
References for the method of manufactured solutions (MMS):
112
- Kambiz Salari and Patrick Knupp (2000)
113
Code Verification by the Method of Manufactured Solutions
114
[DOI: 10.2172/759450](https://doi.org/10.2172/759450)
115
- Patrick J. Roache (2002)
116
Code Verification by the Method of Manufactured Solutions
117
[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
118
"""
119
@inline function source_terms_convergence_test(u, x, t,
120
equations::CompressibleEulerEquationsQuasi1D)
121
# Same settings as in `initial_condition_convergence_test`.
122
# Derivatives calculated with ForwardDiff.jl
123
RealT = eltype(u)
124
c = 2
125
A = convert(RealT, 0.1)
126
L = 2
127
f = 1.0f0 / L
128
ω = 2 * convert(RealT, pi) * f
129
x1, = x
130
ini(x1, t) = c + A * sin(ω * (x1 - t))
131
132
rho(x1, t) = ini(x1, t)
133
v1(x1, t) = 1
134
e_total(x1, t) = ini(x1, t)^2 / rho(x1, t)
135
p1(x1, t) = (equations.gamma - 1) *
136
(e_total(x1, t) - 0.5f0 * rho(x1, t) * v1(x1, t)^2)
137
a(x1, t) = 1.5f0 - 0.5f0 * cospi(x1)
138
139
arho(x1, t) = a(x1, t) * rho(x1, t)
140
arhou(x1, t) = arho(x1, t) * v1(x1, t)
141
aE(x1, t) = a(x1, t) * e_total(x1, t)
142
143
darho_dt(x1, t) = ForwardDiff.derivative(t -> arho(x1, t), t)
144
darhou_dx(x1, t) = ForwardDiff.derivative(x1 -> arhou(x1, t), x1)
145
146
arhouu(x1, t) = arhou(x1, t) * v1(x1, t)
147
darhou_dt(x1, t) = ForwardDiff.derivative(t -> arhou(x1, t), t)
148
darhouu_dx(x1, t) = ForwardDiff.derivative(x1 -> arhouu(x1, t), x1)
149
dp1_dx(x1, t) = ForwardDiff.derivative(x1 -> p1(x1, t), x1)
150
151
auEp(x1, t) = a(x1, t) * v1(x1, t) * (e_total(x1, t) + p1(x1, t))
152
daE_dt(x1, t) = ForwardDiff.derivative(t -> aE(x1, t), t)
153
dauEp_dx(x1, t) = ForwardDiff.derivative(x1 -> auEp(x1, t), x1)
154
155
du1 = darho_dt(x1, t) + darhou_dx(x1, t)
156
du2 = darhou_dt(x1, t) + darhouu_dx(x1, t) + a(x1, t) * dp1_dx(x1, t)
157
du3 = daE_dt(x1, t) + dauEp_dx(x1, t)
158
159
return SVector(du1, du2, du3, 0)
160
end
161
162
# Calculate 1D flux for a single point
163
@inline function flux(u, orientation::Integer,
164
equations::CompressibleEulerEquationsQuasi1D)
165
a_rho, a_rho_v1, a_e_total, a = u
166
rho, v1, p, a = cons2prim(u, equations)
167
e = a_e_total / a
168
169
# Ignore orientation since it is always "1" in 1D
170
f1 = a_rho_v1
171
f2 = a_rho_v1 * v1
172
f3 = a * v1 * (e + p)
173
174
return SVector(f1, f2, f3, 0)
175
end
176
177
"""
178
flux_nonconservative_chan_etal(u_ll, u_rr, orientation::Integer,
179
equations::CompressibleEulerEquationsQuasi1D)
180
flux_nonconservative_chan_etal(u_ll, u_rr, normal_direction,
181
equations::CompressibleEulerEquationsQuasi1D)
182
flux_nonconservative_chan_etal(u_ll, u_rr, normal_ll, normal_rr,
183
equations::CompressibleEulerEquationsQuasi1D)
184
185
Non-symmetric two-point volume flux discretizing the nonconservative (source) term
186
that contains the gradient of the pressure [`CompressibleEulerEquationsQuasi1D`](@ref)
187
and the nozzle width.
188
189
Further details are available in the paper:
190
- Jesse Chan, Khemraj Shukla, Xinhui Wu, Ruofeng Liu, Prani Nalluri (2023)
191
High order entropy stable schemes for the quasi-one-dimensional
192
shallow water and compressible Euler equations
193
[DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089)
194
"""
195
@inline function flux_nonconservative_chan_etal(u_ll, u_rr, orientation::Integer,
196
equations::CompressibleEulerEquationsQuasi1D)
197
#Variables
198
_, _, p_ll, a_ll = cons2prim(u_ll, equations)
199
_, _, p_rr, _ = cons2prim(u_rr, equations)
200
201
# For flux differencing using non-conservative terms, we return the
202
# non-conservative flux scaled by 2. This cancels with a factor of 0.5
203
# in the arithmetic average of {p}.
204
p_avg = p_ll + p_rr
205
206
return SVector(0, a_ll * p_avg, 0, 0)
207
end
208
209
# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that
210
# the normal component is incorporated into the numerical flux.
211
#
212
# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a
213
# similar implementation.
214
@inline function flux_nonconservative_chan_etal(u_ll, u_rr,
215
normal_direction::AbstractVector,
216
equations::CompressibleEulerEquationsQuasi1D)
217
return normal_direction[1] *
218
flux_nonconservative_chan_etal(u_ll, u_rr, 1, equations)
219
end
220
221
@inline function flux_nonconservative_chan_etal(u_ll, u_rr,
222
normal_ll::AbstractVector,
223
normal_rr::AbstractVector,
224
equations::CompressibleEulerEquationsQuasi1D)
225
# normal_ll should be equal to normal_rr in 1D
226
return flux_nonconservative_chan_etal(u_ll, u_rr, normal_ll, equations)
227
end
228
229
"""
230
@inline function flux_chan_etal(u_ll, u_rr, orientation::Integer,
231
equations::CompressibleEulerEquationsQuasi1D)
232
233
Conservative (symmetric) part of the entropy conservative flux for quasi 1D compressible Euler equations split form.
234
This flux is a generalization of [`flux_ranocha`](@ref) for [`CompressibleEulerEquations1D`](@ref).
235
Further details are available in the paper:
236
- Jesse Chan, Khemraj Shukla, Xinhui Wu, Ruofeng Liu, Prani Nalluri (2023)
237
High order entropy stable schemes for the quasi-one-dimensional
238
shallow water and compressible Euler equations
239
[DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089)
240
"""
241
@inline function flux_chan_etal(u_ll, u_rr, orientation::Integer,
242
equations::CompressibleEulerEquationsQuasi1D)
243
# Unpack left and right state
244
rho_ll, v1_ll, p_ll, a_ll = cons2prim(u_ll, equations)
245
rho_rr, v1_rr, p_rr, a_rr = cons2prim(u_rr, equations)
246
247
# Compute the necessary mean values
248
rho_mean = ln_mean(rho_ll, rho_rr)
249
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
250
# in exact arithmetic since
251
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
252
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
253
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
254
v1_avg = 0.5f0 * (v1_ll + v1_rr)
255
a_v1_avg = 0.5f0 * (a_ll * v1_ll + a_rr * v1_rr)
256
p_avg = 0.5f0 * (p_ll + p_rr)
257
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)
258
259
# Calculate fluxes
260
# Ignore orientation since it is always "1" in 1D
261
f1 = rho_mean * a_v1_avg
262
f2 = rho_mean * a_v1_avg * v1_avg
263
f3 = f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
264
0.5f0 * (p_ll * a_rr * v1_rr + p_rr * a_ll * v1_ll)
265
266
return SVector(f1, f2, f3, 0)
267
end
268
269
# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that
270
# the normal component is incorporated into the numerical flux.
271
#
272
# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a
273
# similar implementation.
274
@inline function flux_chan_etal(u_ll, u_rr, normal_direction::AbstractVector,
275
equations::CompressibleEulerEquationsQuasi1D)
276
return normal_direction[1] * flux_chan_etal(u_ll, u_rr, 1, equations)
277
end
278
279
# Calculate estimates for maximum wave speed for local Lax-Friedrichs-type dissipation as the
280
# maximum velocity magnitude plus the maximum speed of sound
281
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
282
equations::CompressibleEulerEquationsQuasi1D)
283
a_rho_ll, a_rho_v1_ll, a_e_total_ll, a_ll = u_ll
284
a_rho_rr, a_rho_v1_rr, a_e_total_rr, a_rr = u_rr
285
286
# Calculate primitive variables and speed of sound
287
rho_ll = a_rho_ll / a_ll
288
e_ll = a_e_total_ll / a_ll
289
v1_ll = a_rho_v1_ll / a_rho_ll
290
v_mag_ll = abs(v1_ll)
291
p_ll = (equations.gamma - 1) * (e_ll - 0.5f0 * rho_ll * v_mag_ll^2)
292
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
293
rho_rr = a_rho_rr / a_rr
294
e_rr = a_e_total_rr / a_rr
295
v1_rr = a_rho_v1_rr / a_rho_rr
296
v_mag_rr = abs(v1_rr)
297
p_rr = (equations.gamma - 1) * (e_rr - 0.5f0 * rho_rr * v_mag_rr^2)
298
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
299
300
return max(v_mag_ll, v_mag_rr) + max(c_ll, c_rr)
301
end
302
303
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
304
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
305
equations::CompressibleEulerEquationsQuasi1D)
306
a_rho_ll, a_rho_v1_ll, a_e_total_ll, a_ll = u_ll
307
a_rho_rr, a_rho_v1_rr, a_e_total_rr, a_rr = u_rr
308
309
# Calculate primitive variables and speed of sound
310
rho_ll = a_rho_ll / a_ll
311
e_ll = a_e_total_ll / a_ll
312
v1_ll = a_rho_v1_ll / a_rho_ll
313
v_mag_ll = abs(v1_ll)
314
p_ll = (equations.gamma - 1) * (e_ll - 0.5f0 * rho_ll * v_mag_ll^2)
315
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
316
rho_rr = a_rho_rr / a_rr
317
e_rr = a_e_total_rr / a_rr
318
v1_rr = a_rho_v1_rr / a_rho_rr
319
v_mag_rr = abs(v1_rr)
320
p_rr = (equations.gamma - 1) * (e_rr - 0.5f0 * rho_rr * v_mag_rr^2)
321
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
322
323
return max(v_mag_ll + c_ll, v_mag_rr + c_rr)
324
end
325
326
@inline function max_abs_speeds(u, equations::CompressibleEulerEquationsQuasi1D)
327
a_rho, a_rho_v1, a_e_total, a = u
328
rho = a_rho / a
329
v1 = a_rho_v1 / a_rho
330
e = a_e_total / a
331
p = (equations.gamma - 1) * (e - 0.5f0 * rho * v1^2)
332
c = sqrt(equations.gamma * p / rho)
333
334
return (abs(v1) + c,)
335
end
336
337
# Convert conservative variables to primitive. We use the convention that the primitive
338
# variables for the quasi-1D equations are `(rho, v1, p)` (i.e., the same as the primitive
339
# variables for `CompressibleEulerEquations1D`)
340
@inline function cons2prim(u, equations::CompressibleEulerEquationsQuasi1D)
341
a_rho, a_rho_v1, a_e_total, a = u
342
q = cons2prim(SVector(a_rho, a_rho_v1, a_e_total) / a,
343
CompressibleEulerEquations1D(equations.gamma))
344
345
return SVector(q[1], q[2], q[3], a)
346
end
347
348
"""
349
entropy(u, equations::CompressibleEulerEquationsQuasi1D)
350
351
The entropy for the quasi-1D compressible Euler equations is the
352
[`entropy(cons, equations::AbstractCompressibleEulerEquations)`](@ref) for the
353
(1D) compressible Euler equations scaled by the channel width `a`.
354
"""
355
@inline function entropy(u, equations::CompressibleEulerEquationsQuasi1D)
356
a_rho, a_rho_v1, a_e_total, a = u
357
return a * entropy(SVector(a_rho, a_rho_v1, a_e_total) / a,
358
CompressibleEulerEquations1D(equations.gamma))
359
end
360
361
# Convert conservative variables to entropy. The entropy variables for the
362
# quasi-1D compressible Euler equations are identical to the entropy variables
363
# for the standard Euler equations for an appropriate definition of `entropy`.
364
@inline function cons2entropy(u, equations::CompressibleEulerEquationsQuasi1D)
365
a_rho, a_rho_v1, a_e_total, a = u
366
w = cons2entropy(SVector(a_rho, a_rho_v1, a_e_total) / a,
367
CompressibleEulerEquations1D(equations.gamma))
368
369
# we follow the convention for other spatially-varying equations and return the spatially
370
# varying coefficient `a` as the final entropy variable.
371
return SVector(w[1], w[2], w[3], a)
372
end
373
374
@inline function entropy2cons(w, equations::CompressibleEulerEquationsQuasi1D)
375
w_rho, w_rho_v1, w_rho_e_total, a = w
376
u = entropy2cons(SVector(w_rho, w_rho_v1, w_rho_e_total),
377
CompressibleEulerEquations1D(equations.gamma))
378
return SVector(a * u[1], a * u[2], a * u[3], a)
379
end
380
381
# Convert primitive to conservative variables
382
@inline function prim2cons(u, equations::CompressibleEulerEquationsQuasi1D)
383
rho, v1, p, a = u
384
q = prim2cons(u, CompressibleEulerEquations1D(equations.gamma))
385
386
return SVector(a * q[1], a * q[2], a * q[3], a)
387
end
388
389
@inline function density(u, equations::CompressibleEulerEquationsQuasi1D)
390
a_rho, _, _, a = u
391
rho = a_rho / a
392
return rho
393
end
394
395
@doc raw"""
396
pressure(u, equations::CompressibleEulerEquationsQuasi1D)
397
398
Computes the pressure for an ideal equation of state with
399
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`,
400
see [`pressure(u, equations::CompressibleEulerEquations1D)`](@ref).
401
"""
402
@inline function pressure(u, equations::CompressibleEulerEquationsQuasi1D)
403
a_rho, a_rho_v1, a_e_total, a = u
404
return pressure(SVector(a_rho, a_rho_v1, a_e_total) / a,
405
CompressibleEulerEquations1D(equations.gamma))
406
end
407
408
@inline function density_pressure(u, equations::CompressibleEulerEquationsQuasi1D)
409
a_rho, a_rho_v1, a_e_total, a = u
410
return density_pressure(SVector(a_rho, a_rho_v1, a_e_total) / a,
411
CompressibleEulerEquations1D(equations.gamma))
412
end
413
end # @muladd
414
415