Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/amr_dg3d.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
# On `P4estMesh`, if an element refines the solution scaled by the Jacobian `J*u` is interpolated
10
# from the parent element into the eight children elements. The solution on each child
11
# element is then recovered by dividing by the new element Jacobians.
12
function refine!(u_ode::AbstractVector, adaptor, mesh::Union{TreeMesh{3}, P4estMesh{3}},
13
equations, dg::DGSEM, cache, elements_to_refine)
14
# Return early if there is nothing to do
15
if isempty(elements_to_refine)
16
if mpi_isparallel()
17
# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do
18
# locally (there still might be other MPI ranks that have refined elements)
19
reinitialize_containers!(mesh, equations, dg, cache)
20
end
21
return
22
end
23
24
# Determine for each existing element whether it needs to be refined
25
needs_refinement = falses(nelements(dg, cache))
26
needs_refinement[elements_to_refine] .= true
27
28
# Retain current solution data
29
old_n_elements = nelements(dg, cache)
30
old_u_ode = copy(u_ode)
31
old_inverse_jacobian = copy(cache.elements.inverse_jacobian)
32
# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed
33
GC.@preserve old_u_ode old_inverse_jacobian begin
34
old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)
35
36
if mesh isa P4estMesh
37
# Loop over all elements in old container and scale the old solution by the Jacobian
38
# prior to projection
39
for old_element_id in 1:old_n_elements
40
for v in eachvariable(equations)
41
old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./
42
old_inverse_jacobian[..,
43
old_element_id])
44
end
45
end
46
end
47
48
@trixi_timeit timer() "reinitialize data structures" begin
49
reinitialize_containers!(mesh, equations, dg, cache)
50
end
51
52
resize!(u_ode,
53
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))
54
u = wrap_array(u_ode, mesh, equations, dg, cache)
55
56
# Loop over all elements in old container and either copy them or refine them
57
u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
58
nnodes(dg), nnodes(dg))
59
u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
60
nnodes(dg), nnodes(dg))
61
element_id = 1
62
for old_element_id in 1:old_n_elements
63
if needs_refinement[old_element_id]
64
# Refine element and store solution directly in new data structure
65
refine_element!(u, element_id, old_u, old_element_id, adaptor,
66
equations, dg, u_tmp1, u_tmp2)
67
68
if mesh isa P4estMesh
69
# Before `element_id` is incremented, divide by the new Jacobians on each
70
# child element and save the result
71
for m in 0:7 # loop over the children
72
for v in eachvariable(equations)
73
u[v, .., element_id + m] .*= (0.125f0 .*
74
cache.elements.inverse_jacobian[..,
75
element_id + m])
76
end
77
end
78
end
79
80
# Increment `element_id` on the refined mesh with the number
81
# of children, i.e., 8 in 3D
82
element_id += 2^ndims(mesh)
83
else
84
if mesh isa P4estMesh
85
# Copy old element data to new element container and remove Jacobian scaling
86
for v in eachvariable(equations)
87
u[v, .., element_id] .= (old_u[v, .., old_element_id] .*
88
old_inverse_jacobian[..,
89
old_element_id])
90
end
91
else # isa TreeMesh
92
@views u[:, .., element_id] .= old_u[:, .., old_element_id]
93
end
94
# No refinement occurred, so increment `element_id` on the new mesh by one
95
element_id += 1
96
end
97
end
98
# If everything is correct, we should have processed all elements.
99
# Depending on whether the last element processed above had to be refined or not,
100
# the counter `element_id` can have two different values at the end.
101
@assert element_id ==
102
nelements(dg, cache) +
103
1||element_id == nelements(dg, cache) + 2^ndims(mesh) "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"
104
end # GC.@preserve old_u_ode old_inverse_jacobian
105
106
# Sanity check
107
if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0
108
@assert ninterfaces(cache.interfaces)==ndims(mesh) * nelements(dg, cache) ("For $(ndims(mesh))D and periodic domains and conforming elements, the number of interfaces must be $(ndims(mesh)) times the number of elements")
109
end
110
111
return nothing
112
end
113
114
# TODO: Taal compare performance of different implementations
115
# Refine solution data u for an element, using L2 projection (interpolation)
116
function refine_element!(u::AbstractArray{<:Any, 5}, element_id,
117
old_u, old_element_id,
118
adaptor::LobattoLegendreAdaptorL2, equations, dg,
119
u_tmp1, u_tmp2)
120
@unpack forward_upper, forward_lower = adaptor
121
122
# Store new element ids
123
bottom_lower_left_id = element_id
124
bottom_lower_right_id = element_id + 1
125
bottom_upper_left_id = element_id + 2
126
bottom_upper_right_id = element_id + 3
127
top_lower_left_id = element_id + 4
128
top_lower_right_id = element_id + 5
129
top_upper_left_id = element_id + 6
130
top_upper_right_id = element_id + 7
131
132
@boundscheck begin
133
@assert old_element_id >= 1
134
@assert size(old_u, 1) == nvariables(equations)
135
@assert size(old_u, 2) == nnodes(dg)
136
@assert size(old_u, 3) == nnodes(dg)
137
@assert size(old_u, 4) == nnodes(dg)
138
@assert size(old_u, 5) >= old_element_id
139
@assert element_id >= 1
140
@assert size(u, 1) == nvariables(equations)
141
@assert size(u, 2) == nnodes(dg)
142
@assert size(u, 3) == nnodes(dg)
143
@assert size(u, 4) == nnodes(dg)
144
@assert size(u, 5) >= element_id + 7
145
end
146
147
# Interpolate to bottom lower left element
148
multiply_dimensionwise!(view(u, :, :, :, :, bottom_lower_left_id), forward_lower,
149
forward_lower, forward_lower,
150
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
151
152
# Interpolate to bottom lower right element
153
multiply_dimensionwise!(view(u, :, :, :, :, bottom_lower_right_id), forward_upper,
154
forward_lower, forward_lower,
155
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
156
157
# Interpolate to bottom upper left element
158
multiply_dimensionwise!(view(u, :, :, :, :, bottom_upper_left_id), forward_lower,
159
forward_upper, forward_lower,
160
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
161
162
# Interpolate to bottom upper right element
163
multiply_dimensionwise!(view(u, :, :, :, :, bottom_upper_right_id), forward_upper,
164
forward_upper, forward_lower,
165
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
166
167
# Interpolate to top lower left element
168
multiply_dimensionwise!(view(u, :, :, :, :, top_lower_left_id), forward_lower,
169
forward_lower, forward_upper,
170
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
171
172
# Interpolate to top lower right element
173
multiply_dimensionwise!(view(u, :, :, :, :, top_lower_right_id), forward_upper,
174
forward_lower, forward_upper,
175
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
176
177
# Interpolate to top upper left element
178
multiply_dimensionwise!(view(u, :, :, :, :, top_upper_left_id), forward_lower,
179
forward_upper, forward_upper,
180
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
181
182
# Interpolate to top upper right element
183
multiply_dimensionwise!(view(u, :, :, :, :, top_upper_right_id), forward_upper,
184
forward_upper, forward_upper,
185
view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)
186
187
return nothing
188
end
189
190
# Coarsen elements in the DG solver based on a list of cell_ids that should be removed.
191
# On `P4estMesh`, if an element coarsens the solution scaled by the Jacobian `J*u` is projected
192
# from the eight children elements back onto the parent element. The solution on the parent
193
# element is then recovered by dividing by the new element Jacobian.
194
function coarsen!(u_ode::AbstractVector, adaptor,
195
mesh::Union{TreeMesh{3}, P4estMesh{3}},
196
equations, dg::DGSEM, cache, elements_to_remove)
197
# Return early if there is nothing to do
198
if isempty(elements_to_remove)
199
if mpi_isparallel()
200
# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do
201
# locally (there still might be other MPI ranks that have coarsened elements)
202
reinitialize_containers!(mesh, equations, dg, cache)
203
end
204
return
205
end
206
207
# Determine for each old element whether it needs to be removed
208
to_be_removed = falses(nelements(dg, cache))
209
to_be_removed[elements_to_remove] .= true
210
211
# Retain current solution data and Jacobians
212
old_n_elements = nelements(dg, cache)
213
old_u_ode = copy(u_ode)
214
old_inverse_jacobian = copy(cache.elements.inverse_jacobian)
215
# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed
216
GC.@preserve old_u_ode old_inverse_jacobian begin
217
old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)
218
219
if mesh isa P4estMesh
220
# Loop over all elements in old container and scale the old solution by the Jacobian
221
# prior to projection
222
for old_element_id in 1:old_n_elements
223
for v in eachvariable(equations)
224
old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./
225
old_inverse_jacobian[..,
226
old_element_id])
227
end
228
end
229
end
230
231
@trixi_timeit timer() "reinitialize data structures" begin
232
reinitialize_containers!(mesh, equations, dg, cache)
233
end
234
235
resize!(u_ode,
236
nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))
237
u = wrap_array(u_ode, mesh, equations, dg, cache)
238
239
# Loop over all elements in old container and either copy them or coarsen them
240
u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
241
nnodes(dg), nnodes(dg))
242
u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
243
nnodes(dg), nnodes(dg))
244
skip = 0
245
element_id = 1
246
for old_element_id in 1:old_n_elements
247
# If skip is non-zero, we just coarsened 2^ndims elements and need to omit the following elements
248
if skip > 0
249
skip -= 1
250
continue
251
end
252
253
if to_be_removed[old_element_id]
254
# If an element is to be removed, sanity check if the following elements
255
# are also marked - otherwise there would be an error in the way the
256
# cells/elements are sorted
257
@assert all(to_be_removed[old_element_id:(old_element_id + 2^ndims(mesh) - 1)]) "bad cell/element order"
258
259
# Coarsen elements and store solution directly in new data structure
260
coarsen_elements!(u, element_id, old_u, old_element_id, adaptor,
261
equations, dg, u_tmp1, u_tmp2)
262
263
if mesh isa P4estMesh
264
# Before `element_id` is incremented, divide by the new Jacobian and save
265
# the result in the parent element
266
for v in eachvariable(equations)
267
u[v, .., element_id] .*= (8 .*
268
cache.elements.inverse_jacobian[..,
269
element_id])
270
end
271
end
272
273
# Increment `element_id` on the coarsened mesh by one and `skip` = 7 in 3D
274
# because 8 children elements become 1 parent element
275
element_id += 1
276
skip = 2^ndims(mesh) - 1
277
else
278
if mesh isa P4estMesh
279
# Copy old element data to new element container and remove Jacobian scaling
280
for v in eachvariable(equations)
281
u[v, .., element_id] .= (old_u[v, .., old_element_id] .*
282
old_inverse_jacobian[..,
283
old_element_id])
284
end
285
else # isa TreeMesh
286
@views u[:, .., element_id] .= old_u[:, .., old_element_id]
287
end
288
# No coarsening occurred, so increment `element_id` on the new mesh by one
289
element_id += 1
290
end
291
end
292
# If everything is correct, we should have processed all elements.
293
@assert element_id==nelements(dg, cache) + 1 "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"
294
end # GC.@preserve old_u_ode old_inverse_jacobian
295
296
# Sanity check
297
if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0
298
@assert ninterfaces(cache.interfaces)==ndims(mesh) * nelements(dg, cache) ("For $(ndims(mesh))D and periodic domains and conforming elements, the number of interfaces must be $(ndims(mesh)) times the number of elements")
299
end
300
301
return nothing
302
end
303
304
# TODO: Taal compare performance of different implementations
305
# Coarsen solution data u for four elements, using L2 projection
306
function coarsen_elements!(u::AbstractArray{<:Any, 5}, element_id,
307
old_u, old_element_id,
308
adaptor::LobattoLegendreAdaptorL2, equations, dg,
309
u_tmp1, u_tmp2)
310
@unpack reverse_upper, reverse_lower = adaptor
311
312
# Store old element ids
313
bottom_lower_left_id = old_element_id
314
bottom_lower_right_id = old_element_id + 1
315
bottom_upper_left_id = old_element_id + 2
316
bottom_upper_right_id = old_element_id + 3
317
top_lower_left_id = old_element_id + 4
318
top_lower_right_id = old_element_id + 5
319
top_upper_left_id = old_element_id + 6
320
top_upper_right_id = old_element_id + 7
321
322
@boundscheck begin
323
@assert old_element_id >= 1
324
@assert size(old_u, 1) == nvariables(equations)
325
@assert size(old_u, 2) == nnodes(dg)
326
@assert size(old_u, 3) == nnodes(dg)
327
@assert size(old_u, 4) == nnodes(dg)
328
@assert size(old_u, 5) >= old_element_id + 7
329
@assert element_id >= 1
330
@assert size(u, 1) == nvariables(equations)
331
@assert size(u, 2) == nnodes(dg)
332
@assert size(u, 3) == nnodes(dg)
333
@assert size(u, 4) == nnodes(dg)
334
@assert size(u, 5) >= element_id
335
end
336
337
# Project from bottom lower left element
338
multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,
339
reverse_lower, reverse_lower,
340
view(old_u, :, :, :, :, bottom_lower_left_id), u_tmp1,
341
u_tmp2)
342
343
# Project from bottom lower right element_variables
344
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,
345
reverse_lower, reverse_lower,
346
view(old_u, :, :, :, :, bottom_lower_right_id), u_tmp1,
347
u_tmp2)
348
349
# Project from bottom upper left element
350
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,
351
reverse_upper, reverse_lower,
352
view(old_u, :, :, :, :, bottom_upper_left_id), u_tmp1,
353
u_tmp2)
354
355
# Project from bottom upper right element
356
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,
357
reverse_upper, reverse_lower,
358
view(old_u, :, :, :, :, bottom_upper_right_id), u_tmp1,
359
u_tmp2)
360
361
# Project from top lower left element
362
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,
363
reverse_lower, reverse_upper,
364
view(old_u, :, :, :, :, top_lower_left_id), u_tmp1,
365
u_tmp2)
366
367
# Project from top lower right element
368
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,
369
reverse_lower, reverse_upper,
370
view(old_u, :, :, :, :, top_lower_right_id), u_tmp1,
371
u_tmp2)
372
373
# Project from top upper left element
374
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,
375
reverse_upper, reverse_upper,
376
view(old_u, :, :, :, :, top_upper_left_id), u_tmp1,
377
u_tmp2)
378
379
# Project from top upper right element
380
add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,
381
reverse_upper, reverse_upper,
382
view(old_u, :, :, :, :, top_upper_right_id), u_tmp1,
383
u_tmp2)
384
385
return nothing
386
end
387
388
# Coarsen and refine elements in the DG solver based on a difference list.
389
function adapt!(u_ode::AbstractVector, adaptor, mesh::T8codeMesh{3}, equations,
390
dg::DGSEM, cache, difference)
391
392
# Return early if there is nothing to do.
393
if !any(difference .!= 0)
394
if mpi_isparallel()
395
# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do
396
# locally (there still might be other MPI ranks that have refined elements)
397
reinitialize_containers!(mesh, equations, dg, cache)
398
end
399
return
400
end
401
402
# Number of (local) cells/elements.
403
old_nelems = nelements(dg, cache)
404
new_nelems = ncells(mesh)
405
406
# Local element indices.
407
old_index = 1
408
new_index = 1
409
410
# Note: This is only true for `hexs`.
411
T8_CHILDREN = 8
412
413
# Retain current solution and inverse Jacobian data.
414
old_u_ode = copy(u_ode)
415
old_inverse_jacobian = copy(cache.elements.inverse_jacobian)
416
417
# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed
418
GC.@preserve old_u_ode begin
419
old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)
420
421
# Loop over all elements in old container and scale the old solution by the Jacobian
422
# prior to interpolation or projection
423
for old_element_id in 1:old_nelems
424
for v in eachvariable(equations)
425
old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./
426
old_inverse_jacobian[..,
427
old_element_id])
428
end
429
end
430
431
@trixi_timeit timer() "reinitialize data structures" begin
432
reinitialize_containers!(mesh, equations, dg, cache)
433
end
434
435
resize!(u_ode, nvariables(equations) * ndofs(mesh, dg, cache))
436
u = wrap_array(u_ode, mesh, equations, dg, cache)
437
438
u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
439
nnodes(dg), nnodes(dg))
440
u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),
441
nnodes(dg), nnodes(dg))
442
443
while old_index <= old_nelems && new_index <= new_nelems
444
if difference[old_index] > 0 # Refine.
445
446
# Refine element and store solution directly in new data structure.
447
refine_element!(u, new_index, old_u, old_index, adaptor, equations, dg,
448
u_tmp1, u_tmp2)
449
450
# Before indices are incremented divide by the new Jacobians on each
451
# child element and save the result
452
for m in 0:7 # loop over the children
453
for v in eachvariable(equations)
454
u[v, .., new_index + m] .*= (0.125f0 .*
455
cache.elements.inverse_jacobian[..,
456
new_index + m])
457
end
458
end
459
460
# Increment `old_index` on the original mesh and the `new_index`
461
# on the refined mesh with the number of children, i.e., T8_CHILDREN = 8
462
old_index += 1
463
new_index += T8_CHILDREN
464
465
elseif difference[old_index] < 0 # Coarsen.
466
467
# If an element is to be removed, sanity check if the following elements
468
# are also marked - otherwise there would be an error in the way the
469
# cells/elements are sorted.
470
@assert all(difference[old_index:(old_index + T8_CHILDREN - 1)] .< 0) "bad cell/element order"
471
472
# Coarsen elements and store solution directly in new data structure.
473
coarsen_elements!(u, new_index, old_u, old_index, adaptor, equations,
474
dg, u_tmp1, u_tmp2)
475
476
# Before the indices are incremented divide by the new Jacobian and save
477
# the result again in the parent element
478
for v in eachvariable(equations)
479
u[v, .., new_index] .*= (8 .* cache.elements.inverse_jacobian[..,
480
new_index])
481
end
482
483
# Increment `old_index` on the original mesh with the number of children
484
# (T8_CHILDREN = 8 in 3D) and the `new_index` by one for the single
485
# coarsened element
486
old_index += T8_CHILDREN
487
new_index += 1
488
489
else # No changes.
490
491
# Copy old element data to new element container and remove Jacobian scaling
492
for v in eachvariable(equations)
493
u[v, .., new_index] .= (old_u[v, .., old_index] .*
494
old_inverse_jacobian[.., old_index])
495
end
496
497
# No refinement / coarsening occurred, so increment element index
498
# on each mesh by one
499
old_index += 1
500
new_index += 1
501
end
502
end # while
503
end # GC.@preserve old_u_ode old_inverse_jacobian
504
505
return nothing
506
end
507
end # @muladd
508
509