Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/analysis_dgmulti.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
function calc_error_norms(func, u, t, analyzer,
9
mesh::DGMultiMesh{NDIMS}, equations, initial_condition,
10
dg::DGMulti{NDIMS}, cache, cache_analysis) where {NDIMS}
11
rd = dg.basis
12
md = mesh.md
13
@unpack u_values = cache
14
15
# interpolate u to quadrature points
16
apply_to_each_field(mul_by!(rd.Vq), u_values, u)
17
18
component_l2_errors = zero(eltype(u_values))
19
component_linf_errors = zero(eltype(u_values))
20
for i in each_quad_node_global(mesh, dg, cache)
21
u_exact = initial_condition(SVector(getindex.(md.xyzq, i)), t, equations)
22
error_at_node = func(u_values[i], equations) - func(u_exact, equations)
23
component_l2_errors += md.wJq[i] * error_at_node .^ 2
24
component_linf_errors = max.(component_linf_errors, abs.(error_at_node))
25
end
26
total_volume = sum(md.wJq)
27
return sqrt.(component_l2_errors ./ total_volume), component_linf_errors
28
end
29
30
function integrate(func::Func, u,
31
mesh::DGMultiMesh,
32
equations, dg::DGMulti, cache; normalize = true) where {Func}
33
rd = dg.basis
34
md = mesh.md
35
@unpack u_values = cache
36
37
# interpolate u to quadrature points
38
apply_to_each_field(mul_by!(rd.Vq), u_values, u)
39
40
integral = sum(md.wJq .* func.(u_values, equations))
41
if normalize == true
42
integral /= sum(md.wJq)
43
end
44
return integral
45
end
46
47
function analyze(::typeof(entropy_timederivative), du, u, t,
48
mesh::DGMultiMesh, equations, dg::DGMulti, cache)
49
rd = dg.basis
50
md = mesh.md
51
@unpack u_values = cache
52
53
# interpolate u, du to quadrature points
54
du_values = similar(u_values) # Todo: DGMulti. Can we move this to the analysis cache somehow?
55
apply_to_each_field(mul_by!(rd.Vq), du_values, du)
56
apply_to_each_field(mul_by!(rd.Vq), u_values, u)
57
58
# compute ∫v(u) * du/dt = ∫dS/dt. We can directly compute v(u) instead of computing the entropy
59
# projection here, since the RHS will be projected to polynomials of degree N and testing with
60
# the L2 projection of v(u) would be equivalent to testing with v(u) due to the moment-preserving
61
# property of the L2 projection.
62
dS_dt = zero(eltype(first(du)))
63
for i in Base.OneTo(length(md.wJq))
64
dS_dt += dot(cons2entropy(u_values[i], equations), du_values[i]) * md.wJq[i]
65
end
66
return dS_dt
67
end
68
69
# This function is used in `analyze(::Val{:l2_divb},...)` and `analyze(::Val{:linf_divb},...)`
70
function compute_local_divergence!(local_divergence, element, vector_field,
71
mesh, dg::DGMulti, cache)
72
@unpack md = mesh
73
rd = dg.basis
74
uEltype = eltype(first(vector_field))
75
76
fill!(local_divergence, zero(uEltype))
77
78
# computes dU_i/dx_i = ∑_j dxhat_j/dx_i * dU_i / dxhat_j
79
# dU_i/dx_i is then accumulated into local_divergence.
80
# TODO: DGMulti. Extend to curved elements.
81
for i in eachdim(mesh)
82
for j in eachdim(mesh)
83
geometric_scaling = md.rstxyzJ[i, j][1, element]
84
jth_ref_derivative_matrix = rd.Drst[j]
85
mul!(local_divergence, jth_ref_derivative_matrix, vector_field[i],
86
geometric_scaling, one(uEltype))
87
end
88
end
89
end
90
91
get_component(u::StructArray, i::Int) = StructArrays.component(u, i)
92
get_component(u::AbstractArray{<:SVector}, i::Int) = getindex.(u, i)
93
94
function analyze(::Val{:l2_divb}, du, u, t,
95
mesh::DGMultiMesh, equations::IdealGlmMhdEquations2D,
96
dg::DGMulti, cache)
97
@unpack md = mesh
98
rd = dg.basis
99
B1 = get_component(u, 6)
100
B2 = get_component(u, 7)
101
B = (B1, B2)
102
103
uEltype = eltype(B1)
104
l2norm_divB = zero(uEltype)
105
local_divB = zeros(uEltype, size(B1, 1))
106
for e in eachelement(mesh, dg, cache)
107
compute_local_divergence!(local_divB, e, view.(B, :, e), mesh, dg, cache)
108
109
# TODO: DGMulti. Extend to curved elements.
110
# compute L2 norm squared via J[1, e] * u' * M * u
111
local_l2norm_divB = md.J[1, e] * dot(local_divB, rd.M, local_divB)
112
l2norm_divB += local_l2norm_divB
113
end
114
115
return sqrt(l2norm_divB)
116
end
117
118
function analyze(::Val{:linf_divb}, du, u, t,
119
mesh::DGMultiMesh, equations::IdealGlmMhdEquations2D,
120
dg::DGMulti, cache)
121
B1 = get_component(u, 6)
122
B2 = get_component(u, 7)
123
B = (B1, B2)
124
125
uEltype = eltype(B1)
126
linf_divB = zero(uEltype)
127
local_divB = zeros(uEltype, size(B1, 1))
128
@batch reduction=(max, linf_divb) for e in eachelement(mesh, dg, cache)
129
compute_local_divergence!(local_divB, e, view.(B, :, e), mesh, dg, cache)
130
131
# compute maximum norm
132
linf_divB = max(linf_divB, maximum(abs, local_divB))
133
end
134
135
return linf_divB
136
end
137
138
function create_cache_analysis(analyzer, mesh::DGMultiMesh,
139
equations, dg::DGMulti, cache,
140
RealT, uEltype)
141
md = mesh.md
142
return (;)
143
end
144
145
SolutionAnalyzer(rd::RefElemData) = rd
146
147
nelements(mesh::DGMultiMesh, ::DGMulti, other_args...) = mesh.md.num_elements
148
function nelementsglobal(mesh::DGMultiMesh, solver::DGMulti, cache)
149
if mpi_isparallel()
150
error("`nelementsglobal` is not implemented for `DGMultiMesh` when used in parallel with MPI")
151
else
152
return ndofs(mesh, solver, cache)
153
end
154
end
155
function ndofsglobal(mesh::DGMultiMesh, solver::DGMulti, cache)
156
if mpi_isparallel()
157
error("`ndofsglobal` is not implemented for `DGMultiMesh` when used in parallel with MPI")
158
else
159
return ndofs(mesh, solver, cache)
160
end
161
end
162
end # @muladd
163
164