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