Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/solvers/dgsem_tree/indicators_1d.jl
5591 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 method is directly used when the indicator is constructed as for shock-capturing volume integrals
9
# and by the dimension-independent method called for AMR
10
function create_cache(::Type{IndicatorHennemannGassner},
11
equations::AbstractEquations{1}, basis::LobattoLegendreBasis)
12
uEltype = real(basis)
13
alpha = Vector{uEltype}()
14
alpha_tmp = similar(alpha)
15
16
MVec = MVector{nnodes(basis), uEltype}
17
18
indicator_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]
19
modal_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]
20
21
return (; alpha, alpha_tmp, indicator_threaded, modal_threaded)
22
end
23
24
# Use this function barrier and unpack inside to avoid passing closures to Polyester.jl
25
# with @batch (@threaded).
26
# Otherwise, @threaded does not work here with Julia ARM on macOS.
27
# See https://github.com/JuliaSIMD/Polyester.jl/issues/88.
28
@inline function calc_indicator_hennemann_gassner!(indicator_hg, threshold, parameter_s,
29
u, element, mesh::AbstractMesh{1},
30
equations, dg, cache)
31
@unpack alpha_max, alpha_min, alpha_smooth, variable = indicator_hg
32
@unpack alpha, alpha_tmp, indicator_threaded, modal_threaded = indicator_hg.cache
33
34
indicator = indicator_threaded[Threads.threadid()]
35
modal = modal_threaded[Threads.threadid()]
36
37
# Calculate indicator variables at Gauss-Lobatto nodes
38
for i in eachnode(dg)
39
u_local = get_node_vars(u, equations, dg, i, element)
40
indicator[i] = indicator_hg.variable(u_local, equations)
41
end
42
43
# Convert to modal representation
44
multiply_scalar_dimensionwise!(modal, dg.basis.inverse_vandermonde_legendre,
45
indicator)
46
47
# Calculate total energies without two highest, without highest, and for all modes
48
total_energy_clip2 = zero(eltype(modal))
49
for i in 1:(nnodes(dg) - 2)
50
total_energy_clip2 += modal[i]^2
51
end
52
total_energy_clip1 = total_energy_clip2 + modal[nnodes(dg) - 1]^2
53
total_energy = total_energy_clip1 + modal[nnodes(dg)]^2
54
55
# Calculate energy in higher modes
56
if !(iszero(total_energy))
57
energy_frac_1 = (total_energy - total_energy_clip1) / total_energy
58
else
59
energy_frac_1 = zero(total_energy)
60
end
61
if !(iszero(total_energy_clip1))
62
energy_frac_2 = (total_energy_clip1 - total_energy_clip2) / total_energy_clip1
63
else
64
energy_frac_2 = zero(total_energy_clip1)
65
end
66
energy = max(energy_frac_1, energy_frac_2)
67
68
alpha_element = 1 / (1 + exp(-parameter_s / threshold * (energy - threshold)))
69
70
# Take care of the case close to pure DG
71
if alpha_element < alpha_min
72
alpha_element = zero(alpha_element)
73
end
74
75
# Take care of the case close to pure FV
76
if alpha_element > 1 - alpha_min
77
alpha_element = one(alpha_element)
78
end
79
80
# Clip the maximum amount of FV allowed
81
alpha[element] = min(alpha_max, alpha_element)
82
return nothing
83
end
84
85
# Diffuse alpha values by setting each alpha to at least 50% of neighboring elements' alpha
86
function apply_smoothing!(mesh::TreeMesh{1}, alpha, alpha_tmp, dg, cache)
87
# Copy alpha values such that smoothing is indpedenent of the element access order
88
alpha_tmp .= alpha
89
90
# Loop over interfaces
91
for interface in eachinterface(dg, cache)
92
# Get neighboring element ids
93
left = cache.interfaces.neighbor_ids[1, interface]
94
right = cache.interfaces.neighbor_ids[2, interface]
95
96
# Apply smoothing
97
alpha[left] = max(alpha_tmp[left], 0.5f0 * alpha_tmp[right], alpha[left])
98
alpha[right] = max(alpha_tmp[right], 0.5f0 * alpha_tmp[left], alpha[right])
99
end
100
101
return nothing
102
end
103
104
# this method is directly used when the indicator is constructed as for shock-capturing volume integrals
105
# and by the dimension-independent method called for AMR
106
function create_cache(::Union{Type{IndicatorLöhner}, Type{IndicatorMax}},
107
equations::AbstractEquations{1}, basis::LobattoLegendreBasis)
108
uEltype = real(basis)
109
alpha = Vector{uEltype}()
110
111
MVec = MVector{nnodes(basis), uEltype}
112
113
indicator_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]
114
115
return (; alpha, indicator_threaded)
116
end
117
118
function (löhner::IndicatorLöhner)(u::AbstractArray{<:Any, 3},
119
mesh, equations, dg::DGSEM, cache;
120
kwargs...)
121
@assert nnodes(dg)>=3 "IndicatorLöhner only works for nnodes >= 3 (polydeg > 1)"
122
@unpack alpha, indicator_threaded = löhner.cache
123
@unpack variable = löhner
124
resize!(alpha, nelements(dg, cache))
125
126
@threaded for element in eachelement(dg, cache)
127
indicator = indicator_threaded[Threads.threadid()]
128
129
# Calculate indicator variables at Gauss-Lobatto nodes
130
for i in eachnode(dg)
131
u_local = get_node_vars(u, equations, dg, i, element)
132
indicator[i] = variable(u_local, equations)
133
end
134
135
estimate = zero(real(dg))
136
for i in 2:(nnodes(dg) - 1)
137
# x direction
138
u0 = indicator[i]
139
up = indicator[i + 1]
140
um = indicator[i - 1]
141
estimate = max(estimate, local_löhner_estimate(um, u0, up, löhner))
142
end
143
144
# use the maximum as DG element indicator
145
alpha[element] = estimate
146
end
147
148
return alpha
149
end
150
151
function (indicator_max::IndicatorMax)(u::AbstractArray{<:Any, 3},
152
mesh, equations, dg::DGSEM, cache;
153
kwargs...)
154
@unpack alpha, indicator_threaded = indicator_max.cache
155
resize!(alpha, nelements(dg, cache))
156
indicator_variable = indicator_max.variable
157
158
@threaded for element in eachelement(dg, cache)
159
indicator = indicator_threaded[Threads.threadid()]
160
161
# Calculate indicator variables at Gauss-Lobatto nodes
162
for i in eachnode(dg)
163
u_local = get_node_vars(u, equations, dg, i, element)
164
indicator[i] = indicator_variable(u_local, equations)
165
end
166
167
alpha[element] = maximum(indicator)
168
end
169
170
return alpha
171
end
172
173
function (indicator::IndicatorNodalFunction)(u::AbstractArray{<:Any, 3},
174
mesh, equations, dg::DGSEM, cache;
175
t, kwargs...)
176
node_coordinates = cache.elements.node_coordinates
177
@unpack alpha = indicator.cache
178
resize!(alpha, nelements(dg, cache))
179
# Extract function to local variable to avoid capturing `indicator` in the threaded loop
180
indicator_function = indicator.indicator_function
181
182
@threaded for element in eachelement(dg, cache)
183
estimate = typemin(eltype(alpha))
184
for i in eachnode(dg)
185
u_nodal = get_node_vars(u, equations, dg, i, element)
186
x_nodal = get_node_coords(node_coordinates, equations, dg,
187
i, element)
188
estimate = max(estimate, indicator_function(u_nodal, x_nodal, t))
189
end
190
alpha[element] = estimate
191
end
192
193
return alpha
194
end
195
end # @muladd
196
197