Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/amr_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
# Refine elements in the DG solver based on a list of cell_ids that should be refined
9
function refine!(u_ode::AbstractVector, adaptor, mesh::TreeMesh{1},
10
equations, dg::DGSEM, cache, elements_to_refine)
11
# Return early if there is nothing to do
12
if isempty(elements_to_refine)
13
return
14
end
15
16
# Determine for each existing element whether it needs to be refined
17
needs_refinement = falses(nelements(dg, cache))
18
needs_refinement[elements_to_refine] .= true
19
20
# Retain current solution data
21
old_n_elements = nelements(dg, cache)
22
old_u_ode = copy(u_ode)
23
GC.@preserve old_u_ode begin # OBS! If we don't GC.@preserve old_u_ode, it might be GC'ed
24
old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)
25
26
@trixi_timeit timer() "reinitialize data structures" begin
27
reinitialize_containers!(mesh, equations, dg, cache)
28
end
29
@assert nelements(dg, cache) > old_n_elements
30
31
resize!(u_ode,
32
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))
33
u = wrap_array(u_ode, mesh, equations, dg, cache)
34
35
# Loop over all elements in old container and either copy them or refine them
36
element_id = 1
37
for old_element_id in 1:old_n_elements
38
if needs_refinement[old_element_id]
39
# Refine element and store solution directly in new data structure
40
refine_element!(u, element_id, old_u, old_element_id,
41
adaptor, equations, dg)
42
element_id += 2^ndims(mesh)
43
else
44
# Copy old element data to new element container
45
@views u[:, .., element_id] .= old_u[:, .., old_element_id]
46
element_id += 1
47
end
48
end
49
# If everything is correct, we should have processed all elements.
50
# Depending on whether the last element processed above had to be refined or not,
51
# the counter `element_id` can have two different values at the end.
52
@assert element_id ==
53
nelements(dg, cache) +
54
1||element_id == nelements(dg, cache) + 2^ndims(mesh) "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"
55
end # GC.@preserve old_u_ode
56
57
# Sanity check
58
if isperiodic(mesh.tree)
59
@assert ninterfaces(cache.interfaces)==1 * nelements(dg, cache) ("For 1D and periodic domains, the number of interfaces must be the same as the number of elements")
60
end
61
62
return nothing
63
end
64
65
function refine!(u_ode::AbstractVector, adaptor, mesh::TreeMesh{1},
66
equations, dg::DGSEM, cache, cache_parabolic,
67
elements_to_refine)
68
# Call `refine!` for the hyperbolic part, which does the heavy lifting of
69
# actually transferring the solution to the refined cells
70
refine!(u_ode, adaptor, mesh, equations, dg, cache, elements_to_refine)
71
72
# Resize parabolic helper variables
73
@unpack parabolic_container = cache_parabolic
74
resize!(parabolic_container, equations, dg, cache)
75
76
return nothing
77
end
78
79
# TODO: Taal compare performance of different implementations
80
# Refine solution data u for an element, using L2 projection (interpolation)
81
function refine_element!(u::AbstractArray{<:Any, 3}, element_id,
82
old_u, old_element_id,
83
adaptor::LobattoLegendreAdaptorL2, equations, dg)
84
@unpack forward_upper, forward_lower = adaptor
85
86
# Store new element ids
87
left_id = element_id
88
right_id = element_id + 1
89
90
@boundscheck begin
91
@assert old_element_id >= 1
92
@assert size(old_u, 1) == nvariables(equations)
93
@assert size(old_u, 2) == nnodes(dg)
94
@assert size(old_u, 3) >= old_element_id
95
@assert element_id >= 1
96
@assert size(u, 1) == nvariables(equations)
97
@assert size(u, 2) == nnodes(dg)
98
@assert size(u, 3) >= element_id + 1
99
end
100
101
# Interpolate to left element
102
for i in eachnode(dg)
103
acc = zero(get_node_vars(u, equations, dg, i, element_id))
104
for k in eachnode(dg)
105
acc += get_node_vars(old_u, equations, dg, k, old_element_id) *
106
forward_lower[i, k]
107
end
108
set_node_vars!(u, acc, equations, dg, i, left_id)
109
end
110
111
# Interpolate to right element
112
for i in eachnode(dg)
113
acc = zero(get_node_vars(u, equations, dg, i, element_id))
114
for k in eachnode(dg)
115
acc += get_node_vars(old_u, equations, dg, k, old_element_id) *
116
forward_upper[i, k]
117
end
118
set_node_vars!(u, acc, equations, dg, i, right_id)
119
end
120
121
return nothing
122
end
123
124
# Coarsen elements in the DG solver based on a list of cell_ids that should be removed
125
function coarsen!(u_ode::AbstractVector, adaptor, mesh::TreeMesh{1},
126
equations, dg::DGSEM, cache, elements_to_remove)
127
# Return early if there is nothing to do
128
if isempty(elements_to_remove)
129
return
130
end
131
132
# Determine for each old element whether it needs to be removed
133
to_be_removed = falses(nelements(dg, cache))
134
to_be_removed[elements_to_remove] .= true
135
136
# Retain current solution data
137
old_n_elements = nelements(dg, cache)
138
old_u_ode = copy(u_ode)
139
GC.@preserve old_u_ode begin # OBS! If we don't GC.@preserve old_u_ode, it might be GC'ed
140
old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)
141
142
@trixi_timeit timer() "reinitialize data structures" begin
143
reinitialize_containers!(mesh, equations, dg, cache)
144
end
145
@assert nelements(dg, cache) < old_n_elements
146
147
resize!(u_ode,
148
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))
149
u = wrap_array(u_ode, mesh, equations, dg, cache)
150
151
# Loop over all elements in old container and either copy them or coarsen them
152
skip = 0
153
element_id = 1
154
for old_element_id in 1:old_n_elements
155
# If skip is non-zero, we just coarsened 2^ndims elements and need to omit the following elements
156
if skip > 0
157
skip -= 1
158
continue
159
end
160
161
if to_be_removed[old_element_id]
162
# If an element is to be removed, sanity check if the following elements
163
# are also marked - otherwise there would be an error in the way the
164
# cells/elements are sorted
165
@assert all(to_be_removed[old_element_id:(old_element_id + 2^ndims(mesh) - 1)]) "bad cell/element order"
166
167
# Coarsen elements and store solution directly in new data structure
168
coarsen_elements!(u, element_id, old_u, old_element_id,
169
adaptor, equations, dg)
170
element_id += 1
171
skip = 2^ndims(mesh) - 1
172
else
173
# Copy old element data to new element container
174
@views u[:, .., element_id] .= old_u[:, .., old_element_id]
175
element_id += 1
176
end
177
end
178
# If everything is correct, we should have processed all elements.
179
@assert element_id==nelements(dg, cache) + 1 "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"
180
end # GC.@preserve old_u_ode
181
182
# Sanity check
183
if isperiodic(mesh.tree)
184
@assert ninterfaces(cache.interfaces)==1 * nelements(dg, cache) ("For 1D and periodic domains, the number of interfaces must be the same as the number of elements")
185
end
186
187
return nothing
188
end
189
190
function coarsen!(u_ode::AbstractVector, adaptor, mesh::TreeMesh{1},
191
equations, dg::DGSEM, cache, cache_parabolic,
192
elements_to_remove)
193
# Call `coarsen!` for the hyperbolic part, which does the heavy lifting of
194
# actually transferring the solution to the coarsened cells
195
coarsen!(u_ode, adaptor, mesh, equations, dg, cache, elements_to_remove)
196
197
# Resize parabolic helper variables
198
@unpack parabolic_container = cache_parabolic
199
resize!(parabolic_container, equations, dg, cache)
200
201
return nothing
202
end
203
204
# TODO: Taal compare performance of different implementations
205
# Coarsen solution data u for two elements, using L2 projection
206
function coarsen_elements!(u::AbstractArray{<:Any, 3}, element_id,
207
old_u, old_element_id,
208
adaptor::LobattoLegendreAdaptorL2, equations, dg)
209
@unpack reverse_upper, reverse_lower = adaptor
210
211
# Store old element ids
212
left_id = old_element_id
213
right_id = old_element_id + 1
214
215
@boundscheck begin
216
@assert old_element_id >= 1
217
@assert size(old_u, 1) == nvariables(equations)
218
@assert size(old_u, 2) == nnodes(dg)
219
@assert size(old_u, 3) >= old_element_id + 1
220
@assert element_id >= 1
221
@assert size(u, 1) == nvariables(equations)
222
@assert size(u, 2) == nnodes(dg)
223
@assert size(u, 3) >= element_id
224
end
225
226
for i in eachnode(dg)
227
acc = zero(get_node_vars(u, equations, dg, i, element_id))
228
229
# Project from lower left element
230
for k in eachnode(dg)
231
acc += get_node_vars(old_u, equations, dg, k, left_id) * reverse_lower[i, k]
232
end
233
234
# Project from lower right element
235
for k in eachnode(dg)
236
acc += get_node_vars(old_u, equations, dg, k, right_id) *
237
reverse_upper[i, k]
238
end
239
240
# Update value
241
set_node_vars!(u, acc, equations, dg, i, element_id)
242
end
243
244
return nothing
245
end
246
end # @muladd
247
248