Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/numerical_fluxes.jl
2055 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
# This file contains general numerical fluxes that are not specific to certain equations
9
10
"""
11
flux_central(u_ll, u_rr, orientation_or_normal_direction, equations::AbstractEquations)
12
13
The classical central numerical flux `f((u_ll) + f(u_rr)) / 2`. When this flux is
14
used as volume flux, the discretization is equivalent to the classical weak form
15
DG method (except floating point errors).
16
"""
17
@inline function flux_central(u_ll, u_rr, orientation_or_normal_direction,
18
equations::AbstractEquations)
19
# Calculate regular 1D fluxes
20
f_ll = flux(u_ll, orientation_or_normal_direction, equations)
21
f_rr = flux(u_rr, orientation_or_normal_direction, equations)
22
23
# Average regular fluxes
24
return 0.5f0 * (f_ll + f_rr)
25
end
26
27
"""
28
FluxPlusDissipation(numerical_flux, dissipation)
29
30
Combine a `numerical_flux` with a `dissipation` operator to create a new numerical flux.
31
"""
32
struct FluxPlusDissipation{NumericalFlux, Dissipation}
33
numerical_flux::NumericalFlux
34
dissipation::Dissipation
35
end
36
37
@inline function (numflux::FluxPlusDissipation)(u_ll, u_rr,
38
orientation_or_normal_direction,
39
equations)
40
@unpack numerical_flux, dissipation = numflux
41
42
return (numerical_flux(u_ll, u_rr, orientation_or_normal_direction, equations)
43
+
44
dissipation(u_ll, u_rr, orientation_or_normal_direction, equations))
45
end
46
47
function Base.show(io::IO, f::FluxPlusDissipation)
48
print(io, "FluxPlusDissipation(", f.numerical_flux, ", ", f.dissipation, ")")
49
end
50
51
"""
52
FluxRotated(numerical_flux)
53
54
Compute a `numerical_flux` flux in direction of a normal vector by rotating the solution,
55
computing the numerical flux in x-direction, and rotating the calculated flux back.
56
57
Requires a rotationally invariant equation with equation-specific functions
58
[`rotate_to_x`](@ref) and [`rotate_from_x`](@ref).
59
"""
60
struct FluxRotated{NumericalFlux}
61
numerical_flux::NumericalFlux
62
end
63
64
# Rotated surface flux computation (2D version)
65
@inline function (flux_rotated::FluxRotated)(u,
66
normal_direction::AbstractVector,
67
equations::AbstractEquations{2})
68
@unpack numerical_flux = flux_rotated
69
70
norm_ = norm(normal_direction)
71
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
72
normal_vector = normal_direction / norm_
73
74
u_rotated = rotate_to_x(u, normal_vector, equations)
75
76
f = numerical_flux(u_rotated, 1, equations)
77
78
return rotate_from_x(f, normal_vector, equations) * norm_
79
end
80
81
# Rotated surface flux computation (2D version)
82
@inline function (flux_rotated::FluxRotated)(u_ll, u_rr,
83
normal_direction::AbstractVector,
84
equations::AbstractEquations{2})
85
@unpack numerical_flux = flux_rotated
86
87
norm_ = norm(normal_direction)
88
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
89
normal_vector = normal_direction / norm_
90
91
u_ll_rotated = rotate_to_x(u_ll, normal_vector, equations)
92
u_rr_rotated = rotate_to_x(u_rr, normal_vector, equations)
93
94
f = numerical_flux(u_ll_rotated, u_rr_rotated, 1, equations)
95
96
return rotate_from_x(f, normal_vector, equations) * norm_
97
end
98
99
# Rotated surface flux computation (3D version)
100
@inline function (flux_rotated::FluxRotated)(u_ll, u_rr,
101
normal_direction::AbstractVector,
102
equations::AbstractEquations{3})
103
@unpack numerical_flux = flux_rotated
104
105
# Storing these vectors could increase the performance by 20 percent
106
norm_ = norm(normal_direction)
107
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
108
normal_vector = normal_direction / norm_
109
110
# Some vector that can't be identical to normal_vector (unless normal_vector == 0)
111
tangent1 = SVector(normal_direction[2], normal_direction[3], -normal_direction[1])
112
# Orthogonal projection
113
tangent1 -= dot(normal_vector, tangent1) * normal_vector
114
tangent1 = normalize(tangent1)
115
116
# Third orthogonal vector
117
tangent2 = normalize(cross(normal_direction, tangent1))
118
119
u_ll_rotated = rotate_to_x(u_ll, normal_vector, tangent1, tangent2, equations)
120
u_rr_rotated = rotate_to_x(u_rr, normal_vector, tangent1, tangent2, equations)
121
122
f = numerical_flux(u_ll_rotated, u_rr_rotated, 1, equations)
123
124
return rotate_from_x(f, normal_vector, tangent1, tangent2, equations) * norm_
125
end
126
127
Base.show(io::IO, f::FluxRotated) = print(io, "FluxRotated(", f.numerical_flux, ")")
128
129
"""
130
DissipationGlobalLaxFriedrichs(λ)
131
132
Create a global Lax-Friedrichs dissipation operator with dissipation coefficient `.
133
"""
134
struct DissipationGlobalLaxFriedrichs{RealT}
135
λ::RealT
136
end
137
138
@inline function (dissipation::DissipationGlobalLaxFriedrichs)(u_ll, u_rr,
139
orientation::Integer,
140
equations)
141
@unpack λ = dissipation
142
return -λ / 2 * (u_rr - u_ll)
143
end
144
145
@inline function (dissipation::DissipationGlobalLaxFriedrichs)(u_ll, u_rr,
146
normal_direction::AbstractVector,
147
equations)
148
@unpack λ = dissipation
149
return -λ / 2 * norm(normal_direction) * (u_rr - u_ll)
150
end
151
152
function Base.show(io::IO, d::DissipationGlobalLaxFriedrichs)
153
print(io, "DissipationGlobalLaxFriedrichs(", d.λ, ")")
154
end
155
156
"""
157
DissipationLocalLaxFriedrichs(max_abs_speed=max_abs_speed)
158
159
Create a local Lax-Friedrichs dissipation operator where the maximum absolute wave speed
160
is estimated as
161
`max_abs_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,
162
defaulting to [`max_abs_speed`](@ref).
163
"""
164
struct DissipationLocalLaxFriedrichs{MaxAbsSpeed}
165
max_abs_speed::MaxAbsSpeed
166
end
167
168
DissipationLocalLaxFriedrichs() = DissipationLocalLaxFriedrichs(max_abs_speed)
169
170
@inline function (dissipation::DissipationLocalLaxFriedrichs)(u_ll, u_rr,
171
orientation_or_normal_direction,
172
equations)
173
λ = dissipation.max_abs_speed(u_ll, u_rr, orientation_or_normal_direction,
174
equations)
175
return -0.5f0 * λ * (u_rr - u_ll)
176
end
177
178
function Base.show(io::IO, d::DissipationLocalLaxFriedrichs)
179
print(io, "DissipationLocalLaxFriedrichs(", d.max_abs_speed, ")")
180
end
181
182
"""
183
max_abs_speed_naive(u_ll, u_rr, orientation::Integer, equations)
184
max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector, equations)
185
186
Simple and fast estimate of the maximal wave speed of the Riemann problem with left and right states
187
`u_ll, u_rr`, based only on the local wave speeds associated to `u_ll` and `u_rr`.
188
189
For non-integer arguments `normal_direction` in one dimension, `max_abs_speed_naive` returns
190
`abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)`.
191
192
Slightly more diffusive/overestimating than [`max_abs_speed`](@ref).
193
"""
194
function max_abs_speed_naive end
195
196
# for non-integer `orientation_or_normal` arguments.
197
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
198
equations::AbstractEquations{1})
199
return abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)
200
end
201
202
"""
203
max_abs_speed(u_ll, u_rr, orientation::Integer, equations)
204
max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector, equations)
205
206
Simple and fast estimate of the maximal wave speed of the Riemann problem with left and right states
207
`u_ll, u_rr`, based only on the local wave speeds associated to `u_ll` and `u_rr`.
208
Less diffusive, i.e., overestimating than [`max_abs_speed_naive`](@ref).
209
210
In particular, `max_abs_speed(u, u, i, equations)` gives the same result as `max_abs_speeds(u, equations)[i]`,
211
i.e., the wave speeds used in `max_dt` which computes the maximum stable time step size through the
212
[`StepsizeCallback`](@ref).
213
214
For non-integer arguments `normal_direction` in one dimension, `max_abs_speed_naive` returns
215
`abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)`.
216
217
Defaults to [`min_max_speed_naive`](@ref) if no specialized version for the 'equations` at hand is available.
218
"""
219
@inline function max_abs_speed(u_ll, u_rr,
220
orientation_or_normal_direction,
221
equations::AbstractEquations)
222
# Use naive version as "backup" if no specialized version for the equations at hand is available
223
max_abs_speed_naive(u_ll, u_rr, orientation_or_normal_direction, equations)
224
end
225
226
const FluxLaxFriedrichs{MaxAbsSpeed} = FluxPlusDissipation{typeof(flux_central),
227
DissipationLocalLaxFriedrichs{MaxAbsSpeed}}
228
"""
229
FluxLaxFriedrichs(max_abs_speed=max_abs_speed)
230
231
Local Lax-Friedrichs (Rusanov) flux with maximum wave speed estimate provided by
232
`max_abs_speed`, cf. [`DissipationLocalLaxFriedrichs`](@ref) and
233
[`max_abs_speed_naive`](@ref).
234
"""
235
function FluxLaxFriedrichs(max_abs_speed = max_abs_speed)
236
FluxPlusDissipation(flux_central, DissipationLocalLaxFriedrichs(max_abs_speed))
237
end
238
239
function Base.show(io::IO, f::FluxLaxFriedrichs)
240
print(io, "FluxLaxFriedrichs(", f.dissipation.max_abs_speed, ")")
241
end
242
243
"""
244
flux_lax_friedrichs
245
246
See [`FluxLaxFriedrichs`](@ref).
247
"""
248
const flux_lax_friedrichs = FluxLaxFriedrichs()
249
250
@doc raw"""
251
DissipationLaxFriedrichsEntropyVariables(max_abs_speed=max_abs_speed)
252
253
Create a local Lax-Friedrichs-type dissipation operator that is provably entropy stable. This operator
254
must be used together with an entropy-conservative two-point flux function (e.g., `flux_ec`) to yield
255
an entropy-stable surface flux. The surface flux function can be initialized as:
256
```julia
257
flux_es = FluxPlusDissipation(flux_ec, DissipationLaxFriedrichsEntropyVariables())
258
```
259
260
In particular, the numerical flux has the form
261
```math
262
f^{\mathrm{ES}} = f^{\mathrm{EC}} - \frac{1}{2} \lambda_{\mathrm{max}} H (w_r - w_l),
263
```
264
where ``f^{\mathrm{EC}}`` is the entropy-conservative two-point flux function (computed with, e.g., `flux_ec`), ``\lambda_{\mathrm{max}}``
265
is the maximum wave speed estimated as `max_abs_speed(u_l, u_r, orientation_or_normal_direction, equations)`,
266
defaulting to [`max_abs_speed`](@ref), ``H`` is a symmetric positive-definite dissipation matrix that
267
depends on the left and right states `u_l` and `u_r`, and ``(w_r - w_l)`` is the jump in entropy variables.
268
Ideally, ``H (w_r - w_l) = (u_r - u_l)``, such that the dissipation operator is consistent with the local
269
Lax-Friedrichs dissipation.
270
271
The entropy-stable dissipation operator is computed with the function
272
`function (dissipation::DissipationLaxFriedrichsEntropyVariables)(u_l, u_r, orientation_or_normal_direction, equations)`,
273
which must be specialized for each equation.
274
275
For the derivation of the dissipation matrix for the multi-ion GLM-MHD equations, see:
276
- A. Rueda-Ramírez, A. Sikstel, G. Gassner, An Entropy-Stable Discontinuous Galerkin Discretization
277
of the Ideal Multi-Ion Magnetohydrodynamics System (2024). Journal of Computational Physics.
278
[DOI: 10.1016/j.jcp.2024.113655](https://doi.org/10.1016/j.jcp.2024.113655).
279
"""
280
struct DissipationLaxFriedrichsEntropyVariables{MaxAbsSpeed}
281
max_abs_speed::MaxAbsSpeed
282
end
283
284
DissipationLaxFriedrichsEntropyVariables() = DissipationLaxFriedrichsEntropyVariables(max_abs_speed)
285
286
function Base.show(io::IO, d::DissipationLaxFriedrichsEntropyVariables)
287
print(io, "DissipationLaxFriedrichsEntropyVariables(", d.max_abs_speed, ")")
288
end
289
290
@doc raw"""
291
DissipationMatrixWintersEtal()
292
293
Creates the Roe-like entropy-stable matrix dissipation operator from Winters et al. (2017). This operator
294
must be used together with an entropy-conservative two-point flux function
295
(e.g., [`flux_ranocha`](@ref) or [`flux_chandrashekar`](@ref)) to yield
296
an entropy-stable surface flux. The surface flux function can be initialized as:
297
```julia
298
flux_es = FluxPlusDissipation(flux_ec, DissipationMatrixWintersEtal())
299
```
300
The 1D and 2D implementations are adapted from the [Atum.jl library](https://github.com/mwarusz/Atum.jl/blob/c7ed44f2b7972ac726ef345da7b98b0bda60e2a3/src/balancelaws/euler.jl#L198).
301
The 3D implementation is adapted from the [FLUXO library](https://github.com/project-fluxo/fluxo)
302
303
For the derivation of the matrix dissipation operator, see:
304
- A. R. Winters, D. Derigs, G. Gassner, S. Walch, A uniquely defined entropy stable matrix dissipation operator
305
for high Mach number ideal MHD and compressible Euler simulations (2017). Journal of Computational Physics.
306
[DOI: 10.1016/j.jcp.2016.12.006](https://doi.org/10.1016/j.jcp.2016.12.006).
307
"""
308
struct DissipationMatrixWintersEtal end
309
310
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
311
orientation::Integer,
312
equations::AbstractEquations{1})
313
return dissipation(u_ll, u_rr, SVector(1), equations)
314
end
315
316
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
317
orientation::Integer,
318
equations::AbstractEquations{2})
319
if orientation == 1
320
return dissipation(u_ll, u_rr, SVector(1, 0), equations)
321
else # orientation == 2
322
return dissipation(u_ll, u_rr, SVector(0, 1), equations)
323
end
324
end
325
326
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
327
orientation::Integer,
328
equations::AbstractEquations{3})
329
if orientation == 1
330
return dissipation(u_ll, u_rr, SVector(1, 0, 0), equations)
331
elseif orientation == 2
332
return dissipation(u_ll, u_rr, SVector(0, 1, 0), equations)
333
else # orientation == 3
334
return dissipation(u_ll, u_rr, SVector(0, 0, 1), equations)
335
end
336
end
337
338
"""
339
FluxHLL(min_max_speed=min_max_speed_davis)
340
341
Create an HLL (Harten, Lax, van Leer) numerical flux where the minimum and maximum
342
wave speeds are estimated as
343
`λ_min, λ_max = min_max_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,
344
defaulting to [`min_max_speed_davis`](@ref).
345
Original paper:
346
- Amiram Harten, Peter D. Lax, Bram van Leer (1983)
347
On Upstream Differencing and Godunov-Type Schemes for Hyperbolic Conservation Laws
348
[DOI: 10.1137/1025002](https://doi.org/10.1137/1025002)
349
"""
350
struct FluxHLL{MinMaxSpeed}
351
min_max_speed::MinMaxSpeed
352
end
353
354
FluxHLL() = FluxHLL(min_max_speed_davis)
355
356
"""
357
min_max_speed_naive(u_ll, u_rr, orientation::Integer, equations)
358
min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector, equations)
359
360
Simple and fast estimate(!) of the minimal and maximal wave speed of the Riemann problem with
361
left and right states `u_ll, u_rr`, usually based only on the local wave speeds associated to
362
`u_ll` and `u_rr`.
363
Slightly more diffusive than [`min_max_speed_davis`](@ref).
364
- Amiram Harten, Peter D. Lax, Bram van Leer (1983)
365
On Upstream Differencing and Godunov-Type Schemes for Hyperbolic Conservation Laws
366
[DOI: 10.1137/1025002](https://doi.org/10.1137/1025002)
367
368
See eq. (10.37) from
369
- Eleuterio F. Toro (2009)
370
Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
371
[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
372
373
See also [`FluxHLL`](@ref), [`min_max_speed_davis`](@ref), [`min_max_speed_einfeldt`](@ref).
374
"""
375
function min_max_speed_naive end
376
377
"""
378
min_max_speed_davis(u_ll, u_rr, orientation::Integer, equations)
379
min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector, equations)
380
381
Simple and fast estimates of the minimal and maximal wave speed of the Riemann problem with
382
left and right states `u_ll, u_rr`, usually based only on the local wave speeds associated to
383
`u_ll` and `u_rr`.
384
385
- S.F. Davis (1988)
386
Simplified Second-Order Godunov-Type Methods
387
[DOI: 10.1137/0909030](https://doi.org/10.1137/0909030)
388
389
See eq. (10.38) from
390
- Eleuterio F. Toro (2009)
391
Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
392
[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
393
See also [`FluxHLL`](@ref), [`min_max_speed_naive`](@ref), [`min_max_speed_einfeldt`](@ref).
394
"""
395
function min_max_speed_davis end
396
397
"""
398
min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations)
399
min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector, equations)
400
401
More advanced mininmal and maximal wave speed computation based on
402
- Bernd Einfeldt (1988)
403
On Godunov-type methods for gas dynamics.
404
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
405
- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)
406
On Godunov-type methods near low densities.
407
[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)
408
409
originally developed for the compressible Euler equations.
410
A compact representation can be found in [this lecture notes, eq. (9.28)](https://metaphor.ethz.ch/x/2019/hs/401-4671-00L/literature/mishra_hyperbolic_pdes.pdf).
411
412
See also [`FluxHLL`](@ref), [`min_max_speed_naive`](@ref), [`min_max_speed_davis`](@ref).
413
"""
414
function min_max_speed_einfeldt end
415
416
@inline function (numflux::FluxHLL)(u_ll, u_rr, orientation_or_normal_direction,
417
equations)
418
λ_min, λ_max = numflux.min_max_speed(u_ll, u_rr, orientation_or_normal_direction,
419
equations)
420
421
if λ_min >= 0 && λ_max >= 0
422
return flux(u_ll, orientation_or_normal_direction, equations)
423
elseif λ_max <= 0 && λ_min <= 0
424
return flux(u_rr, orientation_or_normal_direction, equations)
425
else
426
f_ll = flux(u_ll, orientation_or_normal_direction, equations)
427
f_rr = flux(u_rr, orientation_or_normal_direction, equations)
428
inv_λ_max_minus_λ_min = inv(λ_max - λ_min)
429
factor_ll = λ_max * inv_λ_max_minus_λ_min
430
factor_rr = λ_min * inv_λ_max_minus_λ_min
431
factor_diss = λ_min * λ_max * inv_λ_max_minus_λ_min
432
return factor_ll * f_ll - factor_rr * f_rr + factor_diss * (u_rr - u_ll)
433
end
434
end
435
436
Base.show(io::IO, numflux::FluxHLL) = print(io, "FluxHLL(", numflux.min_max_speed, ")")
437
438
"""
439
flux_hll
440
441
See [`FluxHLL`](@ref).
442
"""
443
const flux_hll = FluxHLL()
444
445
"""
446
flux_hlle
447
448
See [`min_max_speed_einfeldt`](@ref).
449
This is a [`FluxHLL`](@ref)-type two-wave solver with special estimates of the wave speeds.
450
"""
451
const flux_hlle = FluxHLL(min_max_speed_einfeldt)
452
453
"""
454
flux_shima_etal_turbo(u_ll, u_rr, orientation_or_normal_direction, equations)
455
456
Equivalent to [`flux_shima_etal`](@ref) except that it may use specialized
457
methods, e.g., when used with [`VolumeIntegralFluxDifferencing`](@ref).
458
These specialized methods may enable better use of SIMD instructions to
459
increase runtime efficiency on modern hardware.
460
"""
461
@inline function flux_shima_etal_turbo(u_ll, u_rr, orientation_or_normal_direction,
462
equations)
463
flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction, equations)
464
end
465
466
"""
467
flux_ranocha_turbo(u_ll, u_rr, orientation_or_normal_direction, equations)
468
469
Equivalent to [`flux_ranocha`](@ref) except that it may use specialized
470
methods, e.g., when used with [`VolumeIntegralFluxDifferencing`](@ref).
471
These specialized methods may enable better use of SIMD instructions to
472
increase runtime efficiency on modern hardware.
473
"""
474
@inline function flux_ranocha_turbo(u_ll, u_rr, orientation_or_normal_direction,
475
equations)
476
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction, equations)
477
end
478
479
"""
480
flux_wintermeyer_etal(u_ll, u_rr, orientation_or_normal_direction, equations)
481
482
Total energy conservative (mathematical entropy for shallow water equations) split form.
483
When the bottom topography is nonzero this scheme will be well-balanced when used as a `volume_flux`.
484
For the `surface_flux` either [`flux_wintermeyer_etal`](@ref) or [`flux_fjordholm_etal`](@ref) can
485
be used to ensure well-balancedness and entropy conservation.
486
487
!!! note
488
This function is defined in Trixi.jl to have a common interface for the
489
methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)
490
and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).
491
"""
492
function flux_wintermeyer_etal end
493
494
"""
495
flux_nonconservative_wintermeyer_etal(u_ll, u_rr, orientation_or_normal_direction, equations)
496
497
Non-symmetric two-point volume flux discretizing the nonconservative (source) term
498
that contains the gradient of the bottom topography for the shallow water equations.
499
500
Gives entropy conservation and well-balancedness on both the volume and surface when combined with
501
[`flux_wintermeyer_etal`](@ref).
502
503
!!! note
504
This function is defined in Trixi.jl to have a common interface for the
505
methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)
506
and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).
507
"""
508
function flux_nonconservative_wintermeyer_etal end
509
510
"""
511
flux_fjordholm_etal(u_ll, u_rr, orientation_or_normal_direction, equations)
512
513
Total energy conservative (mathematical entropy for shallow water equations). When the bottom topography
514
is nonzero this should only be used as a surface flux otherwise the scheme will not be well-balanced.
515
For well-balancedness in the volume flux use [`flux_wintermeyer_etal`](@ref).
516
517
!!! note
518
This function is defined in Trixi.jl to have a common interface for the
519
methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)
520
and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).
521
"""
522
function flux_fjordholm_etal end
523
524
"""
525
flux_nonconservative_fjordholm_etal(u_ll, u_rr, orientation_or_normal_direction, equations)
526
527
Non-symmetric two-point surface flux discretizing the nonconservative (source) term of
528
that contains the gradient of the bottom topography for the shallow water equations.
529
530
This flux can be used together with [`flux_fjordholm_etal`](@ref) at interfaces to ensure entropy
531
conservation and well-balancedness.
532
533
!!! note
534
This function is defined in Trixi.jl to have a common interface for the
535
methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)
536
and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).
537
"""
538
function flux_nonconservative_fjordholm_etal end
539
540
"""
541
FluxUpwind(splitting)
542
543
A numerical flux `f(u_left, u_right) = f⁺(u_left) + f⁻(u_right)` based on
544
flux vector splitting.
545
546
The [`SurfaceIntegralUpwind`](@ref) with a given `splitting` is equivalent to
547
the [`SurfaceIntegralStrongForm`](@ref) with `FluxUpwind(splitting)`
548
as numerical flux (up to floating point differences). Note, that
549
[`SurfaceIntegralUpwind`](@ref) is only available on [`TreeMesh`](@ref).
550
551
!!! warning "Experimental implementation (upwind SBP)"
552
This is an experimental feature and may change in future releases.
553
"""
554
struct FluxUpwind{Splitting}
555
splitting::Splitting
556
end
557
558
@inline function (numflux::FluxUpwind)(u_ll, u_rr, orientation::Int, equations)
559
@unpack splitting = numflux
560
fm = splitting(u_rr, Val{:minus}(), orientation, equations)
561
fp = splitting(u_ll, Val{:plus}(), orientation, equations)
562
return fm + fp
563
end
564
565
@inline function (numflux::FluxUpwind)(u_ll, u_rr,
566
normal_direction::AbstractVector,
567
equations::AbstractEquations{2})
568
@unpack splitting = numflux
569
f_tilde_m = splitting(u_rr, Val{:minus}(), normal_direction, equations)
570
f_tilde_p = splitting(u_ll, Val{:plus}(), normal_direction, equations)
571
return f_tilde_m + f_tilde_p
572
end
573
574
Base.show(io::IO, f::FluxUpwind) = print(io, "FluxUpwind(", f.splitting, ")")
575
end # @muladd
576
577