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