Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/analysis_dg1d.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
function create_cache_analysis(analyzer, mesh::TreeMesh{1},
9
equations, dg::DG, cache,
10
RealT, uEltype)
11
12
# pre-allocate buffers
13
# We use `StrideArray`s here since these buffers are used in performance-critical
14
# places and the additional information passed to the compiler makes them faster
15
# than native `Array`s.
16
u_local = StrideArray(undef, uEltype,
17
StaticInt(nvariables(equations)), StaticInt(nnodes(analyzer)))
18
x_local = StrideArray(undef, RealT,
19
StaticInt(ndims(equations)), StaticInt(nnodes(analyzer)))
20
21
return (; u_local, x_local)
22
end
23
24
function create_cache_analysis(analyzer, mesh::StructuredMesh{1},
25
equations, dg::DG, cache,
26
RealT, uEltype)
27
28
# pre-allocate buffers
29
# We use `StrideArray`s here since these buffers are used in performance-critical
30
# places and the additional information passed to the compiler makes them faster
31
# than native `Array`s.
32
u_local = StrideArray(undef, uEltype,
33
StaticInt(nvariables(equations)), StaticInt(nnodes(analyzer)))
34
x_local = StrideArray(undef, RealT,
35
StaticInt(ndims(equations)), StaticInt(nnodes(analyzer)))
36
jacobian_local = StrideArray(undef, RealT,
37
StaticInt(nnodes(analyzer)))
38
39
return (; u_local, x_local, jacobian_local)
40
end
41
42
function calc_error_norms(func, u, t, analyzer,
43
mesh::TreeMesh{1}, equations, initial_condition,
44
dg::DGSEM, cache, cache_analysis)
45
@unpack vandermonde, weights = analyzer
46
@unpack node_coordinates = cache.elements
47
@unpack u_local, x_local = cache_analysis
48
49
# Set up data structures
50
l2_error = zero(func(get_node_vars(u, equations, dg, 1, 1), equations))
51
linf_error = copy(l2_error)
52
53
# Iterate over all elements for error calculations
54
for element in eachelement(dg, cache)
55
# Interpolate solution and node locations to analysis nodes
56
multiply_dimensionwise!(u_local, vandermonde, view(u, :, :, element))
57
multiply_dimensionwise!(x_local, vandermonde,
58
view(node_coordinates, :, :, element))
59
60
# Calculate errors at each analysis node
61
volume_jacobian_ = volume_jacobian(element, mesh, cache)
62
63
for i in eachnode(analyzer)
64
u_exact = initial_condition(get_node_coords(x_local, equations, dg, i), t,
65
equations)
66
diff = func(u_exact, equations) -
67
func(get_node_vars(u_local, equations, dg, i), equations)
68
l2_error += diff .^ 2 * (weights[i] * volume_jacobian_)
69
linf_error = @. max(linf_error, abs(diff))
70
end
71
end
72
73
# For L2 error, divide by total volume
74
total_volume_ = total_volume(mesh)
75
l2_error = @. sqrt(l2_error / total_volume_)
76
77
return l2_error, linf_error
78
end
79
80
function calc_error_norms(func, u, t, analyzer,
81
mesh::StructuredMesh{1}, equations, initial_condition,
82
dg::DGSEM, cache, cache_analysis)
83
@unpack vandermonde, weights = analyzer
84
@unpack node_coordinates, inverse_jacobian = cache.elements
85
@unpack u_local, x_local, jacobian_local = cache_analysis
86
87
# Set up data structures
88
l2_error = zero(func(get_node_vars(u, equations, dg, 1, 1), equations))
89
linf_error = copy(l2_error)
90
total_volume = zero(real(mesh))
91
92
# Iterate over all elements for error calculations
93
for element in eachelement(dg, cache)
94
# Interpolate solution and node locations to analysis nodes
95
multiply_dimensionwise!(u_local, vandermonde, view(u, :, :, element))
96
multiply_dimensionwise!(x_local, vandermonde,
97
view(node_coordinates, :, :, element))
98
multiply_scalar_dimensionwise!(jacobian_local, vandermonde,
99
inv.(view(inverse_jacobian, :, element)))
100
101
# Calculate errors at each analysis node
102
for i in eachnode(analyzer)
103
u_exact = initial_condition(get_node_coords(x_local, equations, dg, i), t,
104
equations)
105
diff = func(u_exact, equations) -
106
func(get_node_vars(u_local, equations, dg, i), equations)
107
# We take absolute value as we need the Jacobian here for the volume calculation
108
abs_jacobian_local_i = abs(jacobian_local[i])
109
110
l2_error += diff .^ 2 * (weights[i] * abs_jacobian_local_i)
111
linf_error = @. max(linf_error, abs(diff))
112
total_volume += weights[i] * abs_jacobian_local_i
113
end
114
end
115
116
# For L2 error, divide by total volume
117
l2_error = @. sqrt(l2_error / total_volume)
118
119
return l2_error, linf_error
120
end
121
122
# Use quadrature to numerically integrate a single element.
123
# We do not multiply by the Jacobian to stay in reference space.
124
# This avoids the need to divide the RHS of the DG scheme by the Jacobian when computing
125
# the time derivative of entropy, see `entropy_change_reference_element`.
126
function integrate_reference_element(func::Func, u, element,
127
::Type{<:AbstractMesh{1}}, equations, dg::DGSEM,
128
cache,
129
args...) where {Func}
130
@unpack weights = dg.basis
131
132
# Initialize integral with zeros of the right shape
133
element_integral = zero(func(u, 1, element, equations, dg, args...))
134
135
for i in eachnode(dg)
136
element_integral += weights[i] *
137
func(u, i, element, equations, dg, args...)
138
end
139
140
return element_integral
141
end
142
143
# Calculate ∫_e (∂S/∂u ⋅ ∂u/∂t) dΩ_e where the result on element 'e' is kept in reference space
144
# Note that ∂S/∂u = w(u) with entropy variables w
145
function entropy_change_reference_element(du, u, element,
146
MeshT::Type{<:AbstractMesh{1}},
147
equations, dg::DGSEM, cache, args...)
148
return integrate_reference_element(u, element, MeshT, equations, dg, cache,
149
du) do u, i, element, equations, dg, du
150
u_node = get_node_vars(u, equations, dg, i, element)
151
du_node = get_node_vars(du, equations, dg, i, element)
152
153
dot(cons2entropy(u_node, equations), du_node)
154
end
155
end
156
157
# calculate surface integral of func(u, equations) * normal on the reference element.
158
function surface_integral_reference_element(func::Func, u, element,
159
::Type{<:Union{TreeMesh{1},
160
StructuredMesh{1}}},
161
equations, dg::DGSEM,
162
cache, args...) where {Func}
163
u_left = get_node_vars(u, equations, dg, 1, element)
164
u_right = get_node_vars(u, equations, dg, nnodes(dg), element)
165
166
surface_integral = func(u_right, 1, equations) - func(u_left, 1, equations)
167
return surface_integral
168
end
169
170
function integrate_via_indices(func::Func, u,
171
mesh::TreeMesh{1}, equations, dg::DGSEM, cache,
172
args...; normalize = true) where {Func}
173
@unpack weights = dg.basis
174
175
# Initialize integral with zeros of the right shape
176
integral = zero(func(u, 1, 1, equations, dg, args...))
177
178
# Use quadrature to numerically integrate over entire domain
179
@batch reduction=(+, integral) for element in eachelement(dg, cache)
180
volume_jacobian_ = volume_jacobian(element, mesh, cache)
181
for i in eachnode(dg)
182
integral += volume_jacobian_ * weights[i] *
183
func(u, i, element, equations, dg, args...)
184
end
185
end
186
187
# Normalize with total volume
188
if normalize
189
integral = integral / total_volume(mesh)
190
end
191
192
return integral
193
end
194
195
function integrate_via_indices(func::Func, u,
196
mesh::StructuredMesh{1}, equations, dg::DGSEM, cache,
197
args...; normalize = true) where {Func}
198
@unpack weights = dg.basis
199
200
# Initialize integral with zeros of the right shape
201
integral = zero(func(u, 1, 1, equations, dg, args...))
202
total_volume = zero(real(mesh))
203
204
# Use quadrature to numerically integrate over entire domain
205
@batch reduction=((+, integral), (+, total_volume)) for element in eachelement(dg,
206
cache)
207
for i in eachnode(dg)
208
jacobian_volume = abs(inv(cache.elements.inverse_jacobian[i, element]))
209
integral += jacobian_volume * weights[i] *
210
func(u, i, element, equations, dg, args...)
211
total_volume += jacobian_volume * weights[i]
212
end
213
end
214
# Normalize with total volume
215
if normalize
216
integral = integral / total_volume
217
end
218
219
return integral
220
end
221
222
function integrate(func::Func, u,
223
mesh::Union{TreeMesh{1}, StructuredMesh{1}},
224
equations, dg::Union{DGSEM, FDSBP}, cache;
225
normalize = true) where {Func}
226
integrate_via_indices(u, mesh, equations, dg, cache;
227
normalize = normalize) do u, i, element, equations, dg
228
u_local = get_node_vars(u, equations, dg, i, element)
229
return func(u_local, equations)
230
end
231
end
232
233
function analyze(::typeof(entropy_timederivative), du, u, t,
234
mesh::Union{TreeMesh{1}, StructuredMesh{1}}, equations,
235
dg::Union{DGSEM, FDSBP}, cache)
236
# Calculate ∫(∂S/∂u ⋅ ∂u/∂t)dΩ
237
integrate_via_indices(u, mesh, equations, dg, cache,
238
du) do u, i, element, equations, dg, du
239
u_node = get_node_vars(u, equations, dg, i, element)
240
du_node = get_node_vars(du, equations, dg, i, element)
241
return dot(cons2entropy(u_node, equations), du_node)
242
end
243
end
244
245
function analyze(::Val{:l2_divb}, du, u, t,
246
mesh::TreeMesh{1}, equations::IdealGlmMhdEquations1D,
247
dg::DGSEM, cache)
248
integrate_via_indices(u, mesh, equations, dg, cache,
249
dg.basis.derivative_matrix) do u, i, element, equations, dg,
250
derivative_matrix
251
divb = zero(eltype(u))
252
for k in eachnode(dg)
253
divb += derivative_matrix[i, k] * u[6, k, element]
254
end
255
divb *= cache.elements.inverse_jacobian[element]
256
return divb^2
257
end |> sqrt
258
end
259
260
function analyze(::Val{:linf_divb}, du, u, t,
261
mesh::TreeMesh{1}, equations::IdealGlmMhdEquations1D,
262
dg::DGSEM, cache)
263
@unpack derivative_matrix, weights = dg.basis
264
265
# integrate over all elements to get the divergence-free condition errors
266
linf_divb = zero(eltype(u))
267
@batch reduction=(max, linf_divb) for element in eachelement(dg, cache)
268
for i in eachnode(dg)
269
divb = zero(eltype(u))
270
for k in eachnode(dg)
271
divb += derivative_matrix[i, k] * u[6, k, element]
272
end
273
divb *= cache.elements.inverse_jacobian[element]
274
linf_divb = max(linf_divb, abs(divb))
275
end
276
end
277
278
return linf_divb
279
end
280
end # @muladd
281
282