Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/solvers/dg.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
abstract type AbstractVolumeIntegral end
9
10
function get_element_variables!(element_variables, u, mesh, equations,
11
volume_integral::AbstractVolumeIntegral, dg, cache)
12
nothing
13
end
14
15
# Function to define "element variables" for the SaveSolutionCallback. It does
16
# nothing by default, but can be specialized for certain mesh types. For instance,
17
# parallel meshes output the mpi rank as an "element variable".
18
function get_element_variables!(element_variables, mesh, dg, cache)
19
nothing
20
end
21
22
### Functions to define `node variables` for the `SaveSolutionCallback`. ###
23
24
# Abstract function which is to be overwritten for the specific `node_variable`
25
# in e.g. the elixirs.
26
function get_node_variable end
27
28
# Version for (purely) hyperbolic equations.
29
function get_node_variables!(node_variables, u_ode, mesh, equations,
30
dg, cache)
31
if !isempty(node_variables)
32
u = wrap_array(u_ode, mesh, equations, dg, cache)
33
for var in keys(node_variables)
34
node_variables[var] = get_node_variable(Val(var), u, mesh, equations,
35
dg, cache)
36
end
37
end
38
39
return nothing
40
end
41
# Version for parabolic-extended equations
42
function get_node_variables!(node_variables, u_ode, mesh, equations,
43
dg, cache,
44
equations_parabolic, cache_parabolic)
45
if !isempty(node_variables)
46
u = wrap_array(u_ode, mesh, equations, dg, cache)
47
for var in keys(node_variables)
48
node_variables[var] = get_node_variable(Val(var), u, mesh, equations,
49
dg, cache,
50
equations_parabolic,
51
cache_parabolic)
52
end
53
end
54
55
return nothing
56
end
57
58
"""
59
VolumeIntegralStrongForm()
60
61
The classical strong form volume integral type for FD/DG methods.
62
"""
63
struct VolumeIntegralStrongForm <: AbstractVolumeIntegral end
64
65
"""
66
VolumeIntegralWeakForm()
67
68
The classical weak form volume integral type for DG methods as explained in
69
standard textbooks.
70
71
## References
72
73
- Kopriva (2009)
74
Implementing Spectral Methods for Partial Differential Equations:
75
Algorithms for Scientists and Engineers
76
[doi: 10.1007/978-90-481-2261-5](https://doi.org/10.1007/978-90-481-2261-5)
77
- Hesthaven, Warburton (2007)
78
Nodal Discontinuous Galerkin Methods: Algorithms, Analysis, and
79
Applications
80
[doi: 10.1007/978-0-387-72067-8](https://doi.org/10.1007/978-0-387-72067-8)
81
82
`VolumeIntegralWeakForm()` is only implemented for conserved terms as
83
non-conservative terms should always be discretized in conjunction with a flux-splitting scheme,
84
see [`VolumeIntegralFluxDifferencing`](@ref).
85
This treatment is required to achieve, e.g., entropy-stability or well-balancedness.
86
"""
87
struct VolumeIntegralWeakForm <: AbstractVolumeIntegral end
88
89
create_cache(mesh, equations, ::VolumeIntegralWeakForm, dg, uEltype) = NamedTuple()
90
91
"""
92
VolumeIntegralFluxDifferencing(volume_flux)
93
94
Volume integral type for DG methods based on SBP operators and flux differencing
95
using a symmetric two-point `volume_flux`. This `volume_flux` needs to satisfy
96
the interface of numerical fluxes in Trixi.jl.
97
98
## References
99
100
- LeFloch, Mercier, Rohde (2002)
101
Fully Discrete, Entropy Conservative Schemes of Arbitrary Order
102
[doi: 10.1137/S003614290240069X](https://doi.org/10.1137/S003614290240069X)
103
- Fisher, Carpenter (2013)
104
High-order entropy stable finite difference schemes for nonlinear
105
conservation laws: Finite domains
106
[doi: 10.1016/j.jcp.2013.06.014](https://doi.org/10.1016/j.jcp.2013.06.014)
107
- Hendrik Ranocha (2017)
108
Comparison of Some Entropy Conservative Numerical Fluxes for the Euler Equations
109
[arXiv: 1701.02264](https://arxiv.org/abs/1701.02264)
110
[doi: 10.1007/s10915-017-0618-1](https://doi.org/10.1007/s10915-017-0618-1)
111
- Chen, Shu (2017)
112
Entropy stable high order discontinuous Galerkin methods with suitable
113
quadrature rules for hyperbolic conservation laws
114
[doi: 10.1016/j.jcp.2017.05.025](https://doi.org/10.1016/j.jcp.2017.05.025)
115
"""
116
struct VolumeIntegralFluxDifferencing{VolumeFlux} <: AbstractVolumeIntegral
117
volume_flux::VolumeFlux
118
end
119
120
function Base.show(io::IO, ::MIME"text/plain", integral::VolumeIntegralFluxDifferencing)
121
@nospecialize integral # reduce precompilation time
122
123
if get(io, :compact, false)
124
show(io, integral)
125
else
126
setup = [
127
"volume flux" => integral.volume_flux
128
]
129
summary_box(io, "VolumeIntegralFluxDifferencing", setup)
130
end
131
end
132
133
"""
134
VolumeIntegralShockCapturingHG(indicator; volume_flux_dg=flux_central,
135
volume_flux_fv=flux_lax_friedrichs)
136
137
Shock-capturing volume integral type for DG methods using a convex blending of
138
the finite volume method with numerical flux `volume_flux_fv` and the
139
[`VolumeIntegralFluxDifferencing`](@ref) with volume flux `volume_flux_dg`.
140
The amount of blending is determined by the `indicator`, e.g.,
141
[`IndicatorHennemannGassner`](@ref).
142
143
## References
144
145
- Hennemann, Gassner (2020)
146
"A provably entropy stable subcell shock capturing approach for high order split form DG"
147
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
148
"""
149
struct VolumeIntegralShockCapturingHG{VolumeFluxDG, VolumeFluxFV, Indicator} <:
150
AbstractVolumeIntegral
151
volume_flux_dg::VolumeFluxDG # symmetric, e.g. split-form or entropy-conservative
152
volume_flux_fv::VolumeFluxFV # non-symmetric in general, e.g. entropy-dissipative
153
indicator::Indicator
154
end
155
156
function VolumeIntegralShockCapturingHG(indicator; volume_flux_dg = flux_central,
157
volume_flux_fv = flux_lax_friedrichs)
158
VolumeIntegralShockCapturingHG{typeof(volume_flux_dg), typeof(volume_flux_fv),
159
typeof(indicator)}(volume_flux_dg, volume_flux_fv,
160
indicator)
161
end
162
163
function Base.show(io::IO, mime::MIME"text/plain",
164
integral::VolumeIntegralShockCapturingHG)
165
@nospecialize integral # reduce precompilation time
166
167
if get(io, :compact, false)
168
show(io, integral)
169
else
170
summary_header(io, "VolumeIntegralShockCapturingHG")
171
summary_line(io, "volume flux DG", integral.volume_flux_dg)
172
summary_line(io, "volume flux FV", integral.volume_flux_fv)
173
summary_line(io, "indicator", integral.indicator |> typeof |> nameof)
174
show(increment_indent(io), mime, integral.indicator)
175
summary_footer(io)
176
end
177
end
178
179
function get_element_variables!(element_variables, u, mesh, equations,
180
volume_integral::VolumeIntegralShockCapturingHG, dg,
181
cache)
182
# call the indicator to get up-to-date values for IO
183
volume_integral.indicator(u, mesh, equations, dg, cache)
184
get_element_variables!(element_variables, volume_integral.indicator,
185
volume_integral)
186
end
187
188
"""
189
VolumeIntegralPureLGLFiniteVolume(volume_flux_fv)
190
191
A volume integral that only uses the subcell finite volume schemes of the
192
[`VolumeIntegralShockCapturingHG`](@ref).
193
194
This gives a formally O(1)-accurate finite volume scheme on an LGL-type subcell
195
mesh (LGL = Legendre-Gauss-Lobatto).
196
197
!!! warning "Experimental implementation"
198
This is an experimental feature and may change in future releases.
199
200
## References
201
202
- Hennemann, Gassner (2020)
203
"A provably entropy stable subcell shock capturing approach for high order split form DG"
204
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
205
"""
206
struct VolumeIntegralPureLGLFiniteVolume{VolumeFluxFV} <: AbstractVolumeIntegral
207
volume_flux_fv::VolumeFluxFV # non-symmetric in general, e.g. entropy-dissipative
208
end
209
# TODO: Figure out if this can also be used for Gauss nodes, not just LGL, and adjust the name accordingly
210
211
function Base.show(io::IO, ::MIME"text/plain",
212
integral::VolumeIntegralPureLGLFiniteVolume)
213
@nospecialize integral # reduce precompilation time
214
215
if get(io, :compact, false)
216
show(io, integral)
217
else
218
setup = [
219
"FV flux" => integral.volume_flux_fv
220
]
221
summary_box(io, "VolumeIntegralPureLGLFiniteVolume", setup)
222
end
223
end
224
225
"""
226
VolumeIntegralSubcellLimiting(limiter;
227
volume_flux_dg, volume_flux_fv)
228
229
A subcell limiting volume integral type for DG methods based on subcell blending approaches
230
with a low-order FV method. Used with limiter [`SubcellLimiterIDP`](@ref).
231
232
!!! note
233
Subcell limiting methods are not fully functional on non-conforming meshes. This is
234
mainly because the implementation assumes that low- and high-order schemes have the same
235
surface terms, which is not guaranteed for non-conforming meshes. The low-order scheme
236
with a high-order mortar is not invariant domain preserving.
237
"""
238
struct VolumeIntegralSubcellLimiting{VolumeFluxDG, VolumeFluxFV, Limiter} <:
239
AbstractVolumeIntegral
240
volume_flux_dg::VolumeFluxDG
241
volume_flux_fv::VolumeFluxFV
242
limiter::Limiter
243
end
244
245
function VolumeIntegralSubcellLimiting(limiter; volume_flux_dg,
246
volume_flux_fv)
247
VolumeIntegralSubcellLimiting{typeof(volume_flux_dg), typeof(volume_flux_fv),
248
typeof(limiter)}(volume_flux_dg, volume_flux_fv,
249
limiter)
250
end
251
252
function Base.show(io::IO, mime::MIME"text/plain",
253
integral::VolumeIntegralSubcellLimiting)
254
@nospecialize integral # reduce precompilation time
255
256
if get(io, :compact, false)
257
show(io, integral)
258
else
259
summary_header(io, "VolumeIntegralSubcellLimiting")
260
summary_line(io, "volume flux DG", integral.volume_flux_dg)
261
summary_line(io, "volume flux FV", integral.volume_flux_fv)
262
summary_line(io, "limiter", integral.limiter |> typeof |> nameof)
263
show(increment_indent(io), mime, integral.limiter)
264
summary_footer(io)
265
end
266
end
267
268
# Required to be able to run `SimpleSSPRK33` without `VolumeIntegralSubcellLimiting`
269
Base.resize!(semi, volume_integral::AbstractVolumeIntegral, new_size) = nothing
270
271
function Base.resize!(semi, volume_integral::VolumeIntegralSubcellLimiting, new_size)
272
# Resize container antidiffusive_fluxes
273
resize!(semi.cache.antidiffusive_fluxes, new_size)
274
275
# Resize container subcell_limiter_coefficients
276
@unpack limiter = volume_integral
277
resize!(limiter.cache.subcell_limiter_coefficients, new_size)
278
end
279
280
# TODO: FD. Should this definition live in a different file because it is
281
# not strictly a DG method?
282
"""
283
VolumeIntegralUpwind(splitting)
284
285
Specialized volume integral for finite difference summation-by-parts (FDSBP)
286
solvers. Can be used together with the upwind SBP operators of Mattsson (2017)
287
implemented in SummationByPartsOperators.jl. The `splitting` controls the
288
discretization.
289
290
See also [`splitting_steger_warming`](@ref), [`splitting_lax_friedrichs`](@ref),
291
[`splitting_vanleer_haenel`](@ref).
292
293
## References
294
295
- Mattsson (2017)
296
Diagonal-norm upwind SBP operators
297
[doi: 10.1016/j.jcp.2017.01.042](https://doi.org/10.1016/j.jcp.2017.01.042)
298
299
!!! warning "Experimental implementation (upwind SBP)"
300
This is an experimental feature and may change in future releases.
301
"""
302
struct VolumeIntegralUpwind{FluxSplitting} <: AbstractVolumeIntegral
303
splitting::FluxSplitting
304
end
305
306
function Base.show(io::IO, ::MIME"text/plain", integral::VolumeIntegralUpwind)
307
@nospecialize integral # reduce precompilation time
308
309
if get(io, :compact, false)
310
show(io, integral)
311
else
312
setup = [
313
"flux splitting" => integral.splitting
314
]
315
summary_box(io, "VolumeIntegralUpwind", setup)
316
end
317
end
318
319
abstract type AbstractSurfaceIntegral end
320
321
"""
322
SurfaceIntegralWeakForm(surface_flux=flux_central)
323
324
The classical weak form surface integral type for DG methods as explained in standard
325
textbooks.
326
327
See also [`VolumeIntegralWeakForm`](@ref).
328
329
## References
330
331
- Kopriva (2009)
332
Implementing Spectral Methods for Partial Differential Equations:
333
Algorithms for Scientists and Engineers
334
[doi: 10.1007/978-90-481-2261-5](https://doi.org/10.1007/978-90-481-2261-5)
335
- Hesthaven, Warburton (2007)
336
Nodal Discontinuous Galerkin Methods: Algorithms, Analysis, and
337
Applications
338
[doi: 10.1007/978-0-387-72067-8](https://doi.org/10.1007/978-0-387-72067-8)
339
"""
340
struct SurfaceIntegralWeakForm{SurfaceFlux} <: AbstractSurfaceIntegral
341
surface_flux::SurfaceFlux
342
end
343
344
SurfaceIntegralWeakForm() = SurfaceIntegralWeakForm(flux_central)
345
346
function Base.show(io::IO, ::MIME"text/plain", integral::SurfaceIntegralWeakForm)
347
@nospecialize integral # reduce precompilation time
348
349
if get(io, :compact, false)
350
show(io, integral)
351
else
352
setup = [
353
"surface flux" => integral.surface_flux
354
]
355
summary_box(io, "SurfaceIntegralWeakForm", setup)
356
end
357
end
358
359
"""
360
SurfaceIntegralStrongForm(surface_flux=flux_central)
361
362
The classical strong form surface integral type for FD/DG methods.
363
364
See also [`VolumeIntegralStrongForm`](@ref).
365
"""
366
struct SurfaceIntegralStrongForm{SurfaceFlux} <: AbstractSurfaceIntegral
367
surface_flux::SurfaceFlux
368
end
369
370
SurfaceIntegralStrongForm() = SurfaceIntegralStrongForm(flux_central)
371
372
function Base.show(io::IO, ::MIME"text/plain", integral::SurfaceIntegralStrongForm)
373
@nospecialize integral # reduce precompilation time
374
375
if get(io, :compact, false)
376
show(io, integral)
377
else
378
setup = [
379
"surface flux" => integral.surface_flux
380
]
381
summary_box(io, "SurfaceIntegralStrongForm", setup)
382
end
383
end
384
385
# TODO: FD. Should this definition live in a different file because it is
386
# not strictly a DG method?
387
"""
388
SurfaceIntegralUpwind(splitting)
389
390
Couple elements with upwind simultaneous approximation terms (SATs)
391
that use a particular flux `splitting`, e.g.,
392
[`splitting_steger_warming`](@ref).
393
394
See also [`VolumeIntegralUpwind`](@ref).
395
396
!!! warning "Experimental implementation (upwind SBP)"
397
This is an experimental feature and may change in future releases.
398
"""
399
struct SurfaceIntegralUpwind{FluxSplitting} <: AbstractSurfaceIntegral
400
splitting::FluxSplitting
401
end
402
403
function Base.show(io::IO, ::MIME"text/plain", integral::SurfaceIntegralUpwind)
404
@nospecialize integral # reduce precompilation time
405
406
if get(io, :compact, false)
407
show(io, integral)
408
else
409
setup = [
410
"flux splitting" => integral.splitting
411
]
412
summary_box(io, "SurfaceIntegralUpwind", setup)
413
end
414
end
415
416
"""
417
DG(; basis, mortar, surface_integral, volume_integral)
418
419
Create a discontinuous Galerkin method.
420
If [`basis isa LobattoLegendreBasis`](@ref LobattoLegendreBasis),
421
this creates a [`DGSEM`](@ref).
422
"""
423
struct DG{Basis, Mortar, SurfaceIntegral, VolumeIntegral}
424
basis::Basis
425
mortar::Mortar
426
surface_integral::SurfaceIntegral
427
volume_integral::VolumeIntegral
428
end
429
430
# @eval due to @muladd
431
@eval Adapt.@adapt_structure(DG)
432
433
function Base.show(io::IO, dg::DG)
434
@nospecialize dg # reduce precompilation time
435
436
print(io, "DG{", real(dg), "}(")
437
print(io, dg.basis)
438
print(io, ", ", dg.mortar)
439
print(io, ", ", dg.surface_integral)
440
print(io, ", ", dg.volume_integral)
441
print(io, ")")
442
end
443
444
function Base.show(io::IO, mime::MIME"text/plain", dg::DG)
445
@nospecialize dg # reduce precompilation time
446
447
if get(io, :compact, false)
448
show(io, dg)
449
else
450
summary_header(io, "DG{" * string(real(dg)) * "}")
451
summary_line(io, "basis", dg.basis)
452
summary_line(io, "mortar", dg.mortar)
453
summary_line(io, "surface integral", dg.surface_integral |> typeof |> nameof)
454
show(increment_indent(io), mime, dg.surface_integral)
455
summary_line(io, "volume integral", dg.volume_integral |> typeof |> nameof)
456
if !(dg.volume_integral isa VolumeIntegralWeakForm) &&
457
!(dg.volume_integral isa VolumeIntegralStrongForm)
458
show(increment_indent(io), mime, dg.volume_integral)
459
end
460
summary_footer(io)
461
end
462
end
463
464
Base.summary(io::IO, dg::DG) = print(io, "DG(" * summary(dg.basis) * ")")
465
466
@inline Base.real(dg::DG) = real(dg.basis)
467
468
function get_element_variables!(element_variables, u, mesh, equations, dg::DG, cache)
469
get_element_variables!(element_variables, u, mesh, equations, dg.volume_integral,
470
dg, cache)
471
get_element_variables!(element_variables, mesh, dg, cache)
472
end
473
474
const MeshesDGSEM = Union{TreeMesh, StructuredMesh, StructuredMeshView,
475
UnstructuredMesh2D,
476
P4estMesh, P4estMeshView, T8codeMesh}
477
478
@inline function ndofs(mesh::MeshesDGSEM, dg::DG, cache)
479
nelements(cache.elements) * nnodes(dg)^ndims(mesh)
480
end
481
482
# TODO: Taal performance, 1:nnodes(dg) vs. Base.OneTo(nnodes(dg)) vs. SOneTo(nnodes(dg)) for DGSEM
483
"""
484
eachnode(dg::DG)
485
486
Return an iterator over the indices that specify the location in relevant data structures
487
for the nodes in `dg`.
488
In particular, not the nodes themselves are returned.
489
"""
490
@inline eachnode(dg::DG) = Base.OneTo(nnodes(dg))
491
@inline nnodes(dg::DG) = nnodes(dg.basis)
492
493
# This is used in some more general analysis code and needs to dispatch on the
494
# `mesh` for some combinations of mesh/solver.
495
@inline nelements(mesh, dg::DG, cache) = nelements(dg, cache)
496
@inline function ndofsglobal(mesh, dg::DG, cache)
497
nelementsglobal(mesh, dg, cache) * nnodes(dg)^ndims(mesh)
498
end
499
500
"""
501
eachelement(dg::DG, cache)
502
503
Return an iterator over the indices that specify the location in relevant data structures
504
for the elements in `cache`.
505
In particular, not the elements themselves are returned.
506
"""
507
@inline eachelement(dg::DG, cache) = Base.OneTo(nelements(dg, cache))
508
509
"""
510
eachinterface(dg::DG, cache)
511
512
Return an iterator over the indices that specify the location in relevant data structures
513
for the interfaces in `cache`.
514
In particular, not the interfaces themselves are returned.
515
"""
516
@inline eachinterface(dg::DG, cache) = Base.OneTo(ninterfaces(dg, cache))
517
518
"""
519
eachboundary(dg::DG, cache)
520
521
Return an iterator over the indices that specify the location in relevant data structures
522
for the boundaries in `cache`.
523
In particular, not the boundaries themselves are returned.
524
"""
525
@inline eachboundary(dg::DG, cache) = Base.OneTo(nboundaries(dg, cache))
526
527
"""
528
eachmortar(dg::DG, cache)
529
530
Return an iterator over the indices that specify the location in relevant data structures
531
for the mortars in `cache`.
532
In particular, not the mortars themselves are returned.
533
"""
534
@inline eachmortar(dg::DG, cache) = Base.OneTo(nmortars(dg, cache))
535
536
"""
537
eachmpiinterface(dg::DG, cache)
538
539
Return an iterator over the indices that specify the location in relevant data structures
540
for the MPI interfaces in `cache`.
541
In particular, not the interfaces themselves are returned.
542
"""
543
@inline eachmpiinterface(dg::DG, cache) = Base.OneTo(nmpiinterfaces(dg, cache))
544
545
"""
546
eachmpimortar(dg::DG, cache)
547
548
Return an iterator over the indices that specify the location in relevant data structures
549
for the MPI mortars in `cache`.
550
In particular, not the mortars themselves are returned.
551
"""
552
@inline eachmpimortar(dg::DG, cache) = Base.OneTo(nmpimortars(dg, cache))
553
554
@inline nelements(dg::DG, cache) = nelements(cache.elements)
555
@inline function nelementsglobal(mesh, dg::DG, cache)
556
mpi_isparallel() ? cache.mpi_cache.n_elements_global : nelements(dg, cache)
557
end
558
@inline ninterfaces(dg::DG, cache) = ninterfaces(cache.interfaces)
559
@inline nboundaries(dg::DG, cache) = nboundaries(cache.boundaries)
560
@inline nmortars(dg::DG, cache) = nmortars(cache.mortars)
561
@inline nmpiinterfaces(dg::DG, cache) = nmpiinterfaces(cache.mpi_interfaces)
562
@inline nmpimortars(dg::DG, cache) = nmpimortars(cache.mpi_mortars)
563
564
# The following functions assume an array-of-structs memory layout
565
# We would like to experiment with different memory layout choices
566
# in the future, see
567
# - https://github.com/trixi-framework/Trixi.jl/issues/88
568
# - https://github.com/trixi-framework/Trixi.jl/issues/87
569
# - https://github.com/trixi-framework/Trixi.jl/issues/86
570
@inline function get_node_coords(x, equations, solver::DG, indices...)
571
SVector(ntuple(@inline(idx->x[idx, indices...]), Val(ndims(equations))))
572
end
573
574
"""
575
get_node_vars(u, equations, solver::DG, indices...)
576
577
Return the value of the variable (vector) `u` at a node inside a specific element.
578
The node is specified by the indices `indices...` argument which is a combination of
579
node index inside the element and the element index itself.
580
Thus, in 1D this is a two integer tuple `indices = (i, element)`,
581
in 2D a three integer tuple `indices = (i, j, element)`,
582
and in 3D a four integer tuple `indices = (i, j, k, element)`.
583
It is also possible to call this function without `indices` being explicitly a tuple,
584
i.e., `get_node_vars(u, equations, solver::DG, i, j, k, element)` is also valid.
585
For more details, see the documentation:
586
https://docs.julialang.org/en/v1/manual/functions/#Varargs-Functions
587
"""
588
@inline function get_node_vars(u, equations, solver::DG, indices...)
589
# There is a cut-off at `n == 10` inside of the method
590
# `ntuple(f::F, n::Integer) where F` in Base at ntuple.jl:17
591
# in Julia `v1.5`, leading to type instabilities if
592
# more than ten variables are used. That's why we use
593
# `Val(...)` below.
594
# We use `@inline` to make sure that the `getindex` calls are
595
# really inlined, which might be the default choice of the Julia
596
# compiler for standard `Array`s but not necessarily for more
597
# advanced array types such as `PtrArray`s, cf.
598
# https://github.com/JuliaSIMD/VectorizationBase.jl/issues/55
599
SVector(ntuple(@inline(v->u[v, indices...]), Val(nvariables(equations))))
600
end
601
602
@inline function get_surface_node_vars(u, equations, solver::DG, indices...)
603
# There is a cut-off at `n == 10` inside of the method
604
# `ntuple(f::F, n::Integer) where F` in Base at ntuple.jl:17
605
# in Julia `v1.5`, leading to type instabilities if
606
# more than ten variables are used. That's why we use
607
# `Val(...)` below.
608
u_ll = SVector(ntuple(@inline(v->u[1, v, indices...]), Val(nvariables(equations))))
609
u_rr = SVector(ntuple(@inline(v->u[2, v, indices...]), Val(nvariables(equations))))
610
return u_ll, u_rr
611
end
612
613
@inline function set_node_vars!(u, u_node, equations, solver::DG, indices...)
614
for v in eachvariable(equations)
615
u[v, indices...] = u_node[v]
616
end
617
return nothing
618
end
619
620
@inline function add_to_node_vars!(u, u_node, equations, solver::DG, indices...)
621
for v in eachvariable(equations)
622
u[v, indices...] += u_node[v]
623
end
624
return nothing
625
end
626
627
# Use this function instead of `add_to_node_vars` to speed up
628
# multiply-and-add-to-node-vars operations
629
# See https://github.com/trixi-framework/Trixi.jl/pull/643
630
@inline function multiply_add_to_node_vars!(u, factor, u_node, equations, solver::DG,
631
indices...)
632
for v in eachvariable(equations)
633
u[v, indices...] = u[v, indices...] + factor * u_node[v]
634
end
635
return nothing
636
end
637
638
# Used for analyze_solution
639
SolutionAnalyzer(dg::DG; kwargs...) = SolutionAnalyzer(dg.basis; kwargs...)
640
641
AdaptorAMR(mesh, dg::DG) = AdaptorL2(dg.basis)
642
643
# General structs for discretizations based on the basic principle of
644
# DGSEM (discontinuous Galerkin spectral element method)
645
include("dgsem/dgsem.jl")
646
647
# Finite difference methods using summation by parts (SBP) operators
648
# These methods are very similar to DG methods since they also impose interface
649
# and boundary conditions weakly. Thus, these methods can re-use a lot of
650
# functionality implemented for DGSEM.
651
include("fdsbp_tree/fdsbp.jl")
652
include("fdsbp_unstructured/fdsbp.jl")
653
654
function allocate_coefficients(mesh::AbstractMesh, equations, dg::DG, cache)
655
# We must allocate a `Vector` in order to be able to `resize!` it (AMR).
656
# cf. wrap_array
657
u_ode = similar(cache.elements.node_coordinates, eltype(cache.elements),
658
nvariables(equations) * nnodes(dg)^ndims(mesh) *
659
nelements(dg, cache))
660
fill!(u_ode, zero(eltype(u_ode)))
661
return u_ode
662
end
663
664
@inline function wrap_array(u_ode::AbstractVector, mesh::AbstractMesh, equations,
665
dg::DGSEM, cache)
666
@boundscheck begin
667
@assert length(u_ode) ==
668
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache)
669
end
670
# We would like to use
671
# reshape(u_ode, (nvariables(equations), ntuple(_ -> nnodes(dg), ndims(mesh))..., nelements(dg, cache)))
672
# but that results in
673
# ERROR: LoadError: cannot resize array with shared data
674
# when we resize! `u_ode` during AMR.
675
#
676
# !!! danger "Segfaults"
677
# Remember to `GC.@preserve` temporaries such as copies of `u_ode`
678
# and other stuff that is only used indirectly via `wrap_array` afterwards!
679
680
# Currently, there are problems when AD is used with `PtrArray`s in broadcasts
681
# since LoopVectorization does not support `ForwardDiff.Dual`s. Hence, we use
682
# optimized `PtrArray`s whenever possible and fall back to plain `Array`s
683
# otherwise.
684
if _PREFERENCE_THREADING === :polyester && LoopVectorization.check_args(u_ode)
685
# This version using `PtrArray`s from StrideArrays.jl is very fast and
686
# does not result in allocations.
687
#
688
# !!! danger "Heisenbug"
689
# Do not use this code when `@threaded` uses `Threads.@threads`. There is
690
# a very strange Heisenbug that makes some parts very slow *sometimes*.
691
# In fact, everything can be fast and fine for many cases but some parts
692
# of the RHS evaluation can take *exactly* (!) five seconds randomly...
693
# Hence, this version should only be used when `@threaded` is based on
694
# `@batch` from Polyester.jl or something similar. Using Polyester.jl
695
# is probably the best option since everything will be handed over to
696
# Chris Elrod, one of the best performance software engineers for Julia.
697
PtrArray(pointer(u_ode),
698
(StaticInt(nvariables(equations)),
699
ntuple(_ -> StaticInt(nnodes(dg)), ndims(mesh))...,
700
nelements(dg, cache)))
701
# (nvariables(equations), ntuple(_ -> nnodes(dg), ndims(mesh))..., nelements(dg, cache)))
702
else
703
# The following version is reasonably fast and allows us to `resize!(u_ode, ...)`.
704
ArrayType = Trixi.storage_type(u_ode)
705
unsafe_wrap(ArrayType{eltype(u_ode), ndims(mesh) + 2}, pointer(u_ode),
706
(nvariables(equations), ntuple(_ -> nnodes(dg), ndims(mesh))...,
707
nelements(dg, cache)))
708
end
709
end
710
711
# Finite difference summation by parts (FDSBP) methods
712
@inline function wrap_array(u_ode::AbstractVector, mesh::AbstractMesh, equations,
713
dg::FDSBP, cache)
714
@boundscheck begin
715
@assert length(u_ode) ==
716
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache)
717
end
718
# See comments on the DGSEM version above
719
if _PREFERENCE_THREADING === :polyester && LoopVectorization.check_args(u_ode)
720
# Here, we do not specialize on the number of nodes using `StaticInt` since
721
# - it will not be type stable (SBP operators just store it as a runtime value)
722
# - FD methods tend to use high node counts
723
PtrArray(pointer(u_ode),
724
(StaticInt(nvariables(equations)),
725
ntuple(_ -> nnodes(dg), ndims(mesh))..., nelements(dg, cache)))
726
else
727
# The following version is reasonably fast and allows us to `resize!(u_ode, ...)`.
728
unsafe_wrap(Array{eltype(u_ode), ndims(mesh) + 2}, pointer(u_ode),
729
(nvariables(equations), ntuple(_ -> nnodes(dg), ndims(mesh))...,
730
nelements(dg, cache)))
731
end
732
end
733
734
# General fallback
735
@inline function wrap_array(u_ode::AbstractVector, mesh::AbstractMesh, equations,
736
dg::DG, cache)
737
wrap_array_native(u_ode, mesh, equations, dg, cache)
738
end
739
740
# Like `wrap_array`, but guarantees to return a plain `Array`, which can be better
741
# for interfacing with external C libraries (MPI, HDF5, visualization),
742
# writing solution files etc.
743
@inline function wrap_array_native(u_ode::AbstractVector, mesh::AbstractMesh, equations,
744
dg::DG, cache)
745
@boundscheck begin
746
@assert length(u_ode) ==
747
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache)
748
end
749
unsafe_wrap(Array{eltype(u_ode), ndims(mesh) + 2}, pointer(u_ode),
750
(nvariables(equations), ntuple(_ -> nnodes(dg), ndims(mesh))...,
751
nelements(dg, cache)))
752
end
753
754
function compute_coefficients!(backend::Nothing, u, func, t, mesh::AbstractMesh{1},
755
equations, dg::DG, cache)
756
@threaded for element in eachelement(dg, cache)
757
for i in eachnode(dg)
758
x_node = get_node_coords(cache.elements.node_coordinates, equations, dg, i,
759
element)
760
# Changing the node positions passed to the initial condition by the minimum
761
# amount possible with the current type of floating point numbers allows setting
762
# discontinuous initial data in a simple way. In particular, a check like `if x < x_jump`
763
# works if the jump location `x_jump` is at the position of an interface.
764
if i == 1
765
x_node = SVector(nextfloat(x_node[1]))
766
elseif i == nnodes(dg)
767
x_node = SVector(prevfloat(x_node[1]))
768
end
769
u_node = func(x_node, t, equations)
770
set_node_vars!(u, u_node, equations, dg, i, element)
771
end
772
end
773
774
return nothing
775
end
776
777
function compute_coefficients!(backend::Nothing, u, func, t, mesh::AbstractMesh{2},
778
equations, dg::DG, cache)
779
@unpack node_coordinates = cache.elements
780
@threaded for element in eachelement(dg, cache)
781
compute_coefficients_element!(u, func, t, equations, dg, node_coordinates,
782
element)
783
end
784
785
return nothing
786
end
787
788
function compute_coefficients!(backend::Backend, u, func, t, mesh::AbstractMesh{2},
789
equations, dg::DG, cache)
790
nelements(dg, cache) == 0 && return nothing
791
@unpack node_coordinates = cache.elements
792
kernel! = compute_coefficients_kernel!(backend)
793
kernel!(u, func, t, equations, dg, node_coordinates,
794
ndrange = nelements(dg, cache))
795
return nothing
796
end
797
798
@kernel function compute_coefficients_kernel!(u, func, t, equations,
799
dg::DG, node_coordinates)
800
element = @index(Global)
801
compute_coefficients_element!(u, func, t, equations, dg, node_coordinates, element)
802
end
803
804
function compute_coefficients_element!(u, func, t, equations, dg::DG,
805
node_coordinates, element)
806
for j in eachnode(dg), i in eachnode(dg)
807
x_node = get_node_coords(node_coordinates, equations, dg, i,
808
j, element)
809
u_node = func(x_node, t, equations)
810
set_node_vars!(u, u_node, equations, dg, i, j, element)
811
end
812
813
return nothing
814
end
815
816
function compute_coefficients!(backend::Nothing, u, func, t, mesh::AbstractMesh{3},
817
equations, dg::DG, cache)
818
@threaded for element in eachelement(dg, cache)
819
for k in eachnode(dg), j in eachnode(dg), i in eachnode(dg)
820
x_node = get_node_coords(cache.elements.node_coordinates, equations, dg, i,
821
j, k, element)
822
u_node = func(x_node, t, equations)
823
set_node_vars!(u, u_node, equations, dg, i, j, k, element)
824
end
825
end
826
827
return nothing
828
end
829
830
# Discretizations specific to each mesh type of Trixi.jl
831
832
# If some functionality is shared by multiple combinations of meshes/solvers,
833
# it is defined in the directory of the most basic mesh and solver type.
834
# The most basic solver type in Trixi.jl is DGSEM (historic reasons and background
835
# of the main contributors).
836
# We consider the `TreeMesh` to be the most basic mesh type since it is Cartesian
837
# and was the first mesh in Trixi.jl.
838
# The order of the other mesh types is the same as the include order below.
839
include("dgsem_tree/dg.jl")
840
include("dgsem_structured/dg.jl")
841
include("dgsem_unstructured/dg.jl")
842
include("dgsem_p4est/dg.jl")
843
include("dgsem_t8code/dg.jl")
844
end # @muladd
845
846