Path: blob/main/src/callbacks_step/amr_dg3d.jl
5586 views
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).1# Since these FMAs can increase the performance of many numerical algorithms,2# we need to opt-in explicitly.3# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.4@muladd begin5#! format: noindent67# Refine elements in the DG solver based on a list of cell_ids that should be refined.8# On `P4estMesh`, if an element refines the solution scaled by the Jacobian `J*u` is interpolated9# from the parent element into the eight children elements. The solution on each child10# element is then recovered by dividing by the new element Jacobians.11function refine!(u_ode::AbstractVector, adaptor, mesh::Union{TreeMesh{3}, P4estMesh{3}},12equations, dg::DGSEM, cache, elements_to_refine)13# Return early if there is nothing to do14if isempty(elements_to_refine)15if mpi_isparallel()16# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do17# locally (there still might be other MPI ranks that have refined elements)18reinitialize_containers!(mesh, equations, dg, cache)19end20return21end2223# Determine for each existing element whether it needs to be refined24needs_refinement = falses(nelements(dg, cache))25needs_refinement[elements_to_refine] .= true2627# Retain current solution data28old_n_elements = nelements(dg, cache)29old_u_ode = copy(u_ode)30old_inverse_jacobian = copy(cache.elements.inverse_jacobian)31# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed32GC.@preserve old_u_ode old_inverse_jacobian begin33old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)3435if mesh isa P4estMesh36# Loop over all elements in old container and scale the old solution by the Jacobian37# prior to projection38for old_element_id in 1:old_n_elements39for v in eachvariable(equations)40old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./41old_inverse_jacobian[..,42old_element_id])43end44end45end4647@trixi_timeit timer() "reinitialize data structures" begin48reinitialize_containers!(mesh, equations, dg, cache)49end5051resize!(u_ode,52nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))53u = wrap_array(u_ode, mesh, equations, dg, cache)5455# Loop over all elements in old container and either copy them or refine them56u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),57nnodes(dg), nnodes(dg))58u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),59nnodes(dg), nnodes(dg))60element_id = 161for old_element_id in 1:old_n_elements62if needs_refinement[old_element_id]63# Refine element and store solution directly in new data structure64refine_element!(u, element_id, old_u, old_element_id, adaptor,65equations, dg, u_tmp1, u_tmp2)6667if mesh isa P4estMesh68# Before `element_id` is incremented, divide by the new Jacobians on each69# child element and save the result70for m in 0:7 # loop over the children71for v in eachvariable(equations)72u[v, .., element_id + m] .*= (0.125f0 .*73cache.elements.inverse_jacobian[..,74element_id + m])75end76end77end7879# Increment `element_id` on the refined mesh with the number80# of children, i.e., 8 in 3D81element_id += 2^ndims(mesh)82else83if mesh isa P4estMesh84# Copy old element data to new element container and remove Jacobian scaling85for v in eachvariable(equations)86u[v, .., element_id] .= (old_u[v, .., old_element_id] .*87old_inverse_jacobian[..,88old_element_id])89end90else # isa TreeMesh91@views u[:, .., element_id] .= old_u[:, .., old_element_id]92end93# No refinement occurred, so increment `element_id` on the new mesh by one94element_id += 195end96end97# If everything is correct, we should have processed all elements.98# Depending on whether the last element processed above had to be refined or not,99# the counter `element_id` can have two different values at the end.100@assert element_id ==101nelements(dg, cache) +1021||element_id == nelements(dg, cache) + 2^ndims(mesh) "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"103end # GC.@preserve old_u_ode old_inverse_jacobian104105# Sanity check106if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0107@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")108end109110return nothing111end112113# TODO: Taal compare performance of different implementations114# Refine solution data u for an element, using L2 projection (interpolation)115function refine_element!(u::AbstractArray{<:Any, 5}, element_id,116old_u, old_element_id,117adaptor::LobattoLegendreAdaptorL2, equations, dg,118u_tmp1, u_tmp2)119@unpack forward_upper, forward_lower = adaptor120121# Store new element ids122bottom_lower_left_id = element_id123bottom_lower_right_id = element_id + 1124bottom_upper_left_id = element_id + 2125bottom_upper_right_id = element_id + 3126top_lower_left_id = element_id + 4127top_lower_right_id = element_id + 5128top_upper_left_id = element_id + 6129top_upper_right_id = element_id + 7130131@boundscheck begin132@assert old_element_id >= 1133@assert size(old_u, 1) == nvariables(equations)134@assert size(old_u, 2) == nnodes(dg)135@assert size(old_u, 3) == nnodes(dg)136@assert size(old_u, 4) == nnodes(dg)137@assert size(old_u, 5) >= old_element_id138@assert element_id >= 1139@assert size(u, 1) == nvariables(equations)140@assert size(u, 2) == nnodes(dg)141@assert size(u, 3) == nnodes(dg)142@assert size(u, 4) == nnodes(dg)143@assert size(u, 5) >= element_id + 7144end145146# Interpolate to bottom lower left element147multiply_dimensionwise!(view(u, :, :, :, :, bottom_lower_left_id), forward_lower,148forward_lower, forward_lower,149view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)150151# Interpolate to bottom lower right element152multiply_dimensionwise!(view(u, :, :, :, :, bottom_lower_right_id), forward_upper,153forward_lower, forward_lower,154view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)155156# Interpolate to bottom upper left element157multiply_dimensionwise!(view(u, :, :, :, :, bottom_upper_left_id), forward_lower,158forward_upper, forward_lower,159view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)160161# Interpolate to bottom upper right element162multiply_dimensionwise!(view(u, :, :, :, :, bottom_upper_right_id), forward_upper,163forward_upper, forward_lower,164view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)165166# Interpolate to top lower left element167multiply_dimensionwise!(view(u, :, :, :, :, top_lower_left_id), forward_lower,168forward_lower, forward_upper,169view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)170171# Interpolate to top lower right element172multiply_dimensionwise!(view(u, :, :, :, :, top_lower_right_id), forward_upper,173forward_lower, forward_upper,174view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)175176# Interpolate to top upper left element177multiply_dimensionwise!(view(u, :, :, :, :, top_upper_left_id), forward_lower,178forward_upper, forward_upper,179view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)180181# Interpolate to top upper right element182multiply_dimensionwise!(view(u, :, :, :, :, top_upper_right_id), forward_upper,183forward_upper, forward_upper,184view(old_u, :, :, :, :, old_element_id), u_tmp1, u_tmp2)185186return nothing187end188189# Coarsen elements in the DG solver based on a list of cell_ids that should be removed.190# On `P4estMesh`, if an element coarsens the solution scaled by the Jacobian `J*u` is projected191# from the eight children elements back onto the parent element. The solution on the parent192# element is then recovered by dividing by the new element Jacobian.193function coarsen!(u_ode::AbstractVector, adaptor,194mesh::Union{TreeMesh{3}, P4estMesh{3}},195equations, dg::DGSEM, cache, elements_to_remove)196# Return early if there is nothing to do197if isempty(elements_to_remove)198if mpi_isparallel()199# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do200# locally (there still might be other MPI ranks that have coarsened elements)201reinitialize_containers!(mesh, equations, dg, cache)202end203return204end205206# Determine for each old element whether it needs to be removed207to_be_removed = falses(nelements(dg, cache))208to_be_removed[elements_to_remove] .= true209210# Retain current solution data and Jacobians211old_n_elements = nelements(dg, cache)212old_u_ode = copy(u_ode)213old_inverse_jacobian = copy(cache.elements.inverse_jacobian)214# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed215GC.@preserve old_u_ode old_inverse_jacobian begin216old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)217218if mesh isa P4estMesh219# Loop over all elements in old container and scale the old solution by the Jacobian220# prior to projection221for old_element_id in 1:old_n_elements222for v in eachvariable(equations)223old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./224old_inverse_jacobian[..,225old_element_id])226end227end228end229230@trixi_timeit timer() "reinitialize data structures" begin231reinitialize_containers!(mesh, equations, dg, cache)232end233234resize!(u_ode,235nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))236u = wrap_array(u_ode, mesh, equations, dg, cache)237238# Loop over all elements in old container and either copy them or coarsen them239u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),240nnodes(dg), nnodes(dg))241u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),242nnodes(dg), nnodes(dg))243skip = 0244element_id = 1245for old_element_id in 1:old_n_elements246# If skip is non-zero, we just coarsened 2^ndims elements and need to omit the following elements247if skip > 0248skip -= 1249continue250end251252if to_be_removed[old_element_id]253# If an element is to be removed, sanity check if the following elements254# are also marked - otherwise there would be an error in the way the255# cells/elements are sorted256@assert all(to_be_removed[old_element_id:(old_element_id + 2^ndims(mesh) - 1)]) "bad cell/element order"257258# Coarsen elements and store solution directly in new data structure259coarsen_elements!(u, element_id, old_u, old_element_id, adaptor,260equations, dg, u_tmp1, u_tmp2)261262if mesh isa P4estMesh263# Before `element_id` is incremented, divide by the new Jacobian and save264# the result in the parent element265for v in eachvariable(equations)266u[v, .., element_id] .*= (8 .*267cache.elements.inverse_jacobian[..,268element_id])269end270end271272# Increment `element_id` on the coarsened mesh by one and `skip` = 7 in 3D273# because 8 children elements become 1 parent element274element_id += 1275skip = 2^ndims(mesh) - 1276else277if mesh isa P4estMesh278# Copy old element data to new element container and remove Jacobian scaling279for v in eachvariable(equations)280u[v, .., element_id] .= (old_u[v, .., old_element_id] .*281old_inverse_jacobian[..,282old_element_id])283end284else # isa TreeMesh285@views u[:, .., element_id] .= old_u[:, .., old_element_id]286end287# No coarsening occurred, so increment `element_id` on the new mesh by one288element_id += 1289end290end291# If everything is correct, we should have processed all elements.292@assert element_id==nelements(dg, cache) + 1 "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"293end # GC.@preserve old_u_ode old_inverse_jacobian294295# Sanity check296if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0297@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")298end299300return nothing301end302303# TODO: Taal compare performance of different implementations304# Coarsen solution data u for four elements, using L2 projection305function coarsen_elements!(u::AbstractArray{<:Any, 5}, element_id,306old_u, old_element_id,307adaptor::LobattoLegendreAdaptorL2, equations, dg,308u_tmp1, u_tmp2)309@unpack reverse_upper, reverse_lower = adaptor310311# Store old element ids312bottom_lower_left_id = old_element_id313bottom_lower_right_id = old_element_id + 1314bottom_upper_left_id = old_element_id + 2315bottom_upper_right_id = old_element_id + 3316top_lower_left_id = old_element_id + 4317top_lower_right_id = old_element_id + 5318top_upper_left_id = old_element_id + 6319top_upper_right_id = old_element_id + 7320321@boundscheck begin322@assert old_element_id >= 1323@assert size(old_u, 1) == nvariables(equations)324@assert size(old_u, 2) == nnodes(dg)325@assert size(old_u, 3) == nnodes(dg)326@assert size(old_u, 4) == nnodes(dg)327@assert size(old_u, 5) >= old_element_id + 7328@assert element_id >= 1329@assert size(u, 1) == nvariables(equations)330@assert size(u, 2) == nnodes(dg)331@assert size(u, 3) == nnodes(dg)332@assert size(u, 4) == nnodes(dg)333@assert size(u, 5) >= element_id334end335336# Project from bottom lower left element337multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,338reverse_lower, reverse_lower,339view(old_u, :, :, :, :, bottom_lower_left_id), u_tmp1,340u_tmp2)341342# Project from bottom lower right element_variables343add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,344reverse_lower, reverse_lower,345view(old_u, :, :, :, :, bottom_lower_right_id), u_tmp1,346u_tmp2)347348# Project from bottom upper left element349add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,350reverse_upper, reverse_lower,351view(old_u, :, :, :, :, bottom_upper_left_id), u_tmp1,352u_tmp2)353354# Project from bottom upper right element355add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,356reverse_upper, reverse_lower,357view(old_u, :, :, :, :, bottom_upper_right_id), u_tmp1,358u_tmp2)359360# Project from top lower left element361add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,362reverse_lower, reverse_upper,363view(old_u, :, :, :, :, top_lower_left_id), u_tmp1,364u_tmp2)365366# Project from top lower right element367add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,368reverse_lower, reverse_upper,369view(old_u, :, :, :, :, top_lower_right_id), u_tmp1,370u_tmp2)371372# Project from top upper left element373add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_lower,374reverse_upper, reverse_upper,375view(old_u, :, :, :, :, top_upper_left_id), u_tmp1,376u_tmp2)377378# Project from top upper right element379add_multiply_dimensionwise!(view(u, :, :, :, :, element_id), reverse_upper,380reverse_upper, reverse_upper,381view(old_u, :, :, :, :, top_upper_right_id), u_tmp1,382u_tmp2)383384return nothing385end386387# Coarsen and refine elements in the DG solver based on a difference list.388function adapt!(u_ode::AbstractVector, adaptor, mesh::T8codeMesh{3}, equations,389dg::DGSEM, cache, difference)390391# Return early if there is nothing to do.392if !any(difference .!= 0)393if mpi_isparallel()394# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do395# locally (there still might be other MPI ranks that have refined elements)396reinitialize_containers!(mesh, equations, dg, cache)397end398return399end400401# Number of (local) cells/elements.402old_nelems = nelements(dg, cache)403new_nelems = ncells(mesh)404405# Local element indices.406old_index = 1407new_index = 1408409# Note: This is only true for `hexs`.410T8_CHILDREN = 8411412# Retain current solution and inverse Jacobian data.413old_u_ode = copy(u_ode)414old_inverse_jacobian = copy(cache.elements.inverse_jacobian)415416# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed417GC.@preserve old_u_ode begin418old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)419420# Loop over all elements in old container and scale the old solution by the Jacobian421# prior to interpolation or projection422for old_element_id in 1:old_nelems423for v in eachvariable(equations)424old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./425old_inverse_jacobian[..,426old_element_id])427end428end429430@trixi_timeit timer() "reinitialize data structures" begin431reinitialize_containers!(mesh, equations, dg, cache)432end433434resize!(u_ode, nvariables(equations) * ndofs(mesh, dg, cache))435u = wrap_array(u_ode, mesh, equations, dg, cache)436437u_tmp1 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),438nnodes(dg), nnodes(dg))439u_tmp2 = Array{eltype(u), 4}(undef, nvariables(equations), nnodes(dg),440nnodes(dg), nnodes(dg))441442while old_index <= old_nelems && new_index <= new_nelems443if difference[old_index] > 0 # Refine.444445# Refine element and store solution directly in new data structure.446refine_element!(u, new_index, old_u, old_index, adaptor, equations, dg,447u_tmp1, u_tmp2)448449# Before indices are incremented divide by the new Jacobians on each450# child element and save the result451for m in 0:7 # loop over the children452for v in eachvariable(equations)453u[v, .., new_index + m] .*= (0.125f0 .*454cache.elements.inverse_jacobian[..,455new_index + m])456end457end458459# Increment `old_index` on the original mesh and the `new_index`460# on the refined mesh with the number of children, i.e., T8_CHILDREN = 8461old_index += 1462new_index += T8_CHILDREN463464elseif difference[old_index] < 0 # Coarsen.465466# If an element is to be removed, sanity check if the following elements467# are also marked - otherwise there would be an error in the way the468# cells/elements are sorted.469@assert all(difference[old_index:(old_index + T8_CHILDREN - 1)] .< 0) "bad cell/element order"470471# Coarsen elements and store solution directly in new data structure.472coarsen_elements!(u, new_index, old_u, old_index, adaptor, equations,473dg, u_tmp1, u_tmp2)474475# Before the indices are incremented divide by the new Jacobian and save476# the result again in the parent element477for v in eachvariable(equations)478u[v, .., new_index] .*= (8 .* cache.elements.inverse_jacobian[..,479new_index])480end481482# Increment `old_index` on the original mesh with the number of children483# (T8_CHILDREN = 8 in 3D) and the `new_index` by one for the single484# coarsened element485old_index += T8_CHILDREN486new_index += 1487488else # No changes.489490# Copy old element data to new element container and remove Jacobian scaling491for v in eachvariable(equations)492u[v, .., new_index] .= (old_u[v, .., old_index] .*493old_inverse_jacobian[.., old_index])494end495496# No refinement / coarsening occurred, so increment element index497# on each mesh by one498old_index += 1499new_index += 1500end501end # while502end # GC.@preserve old_u_ode old_inverse_jacobian503504return nothing505end506end # @muladd507508509