Path: blob/main/src/callbacks_step/amr_dg2d.jl
2055 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# Redistribute data for load balancing after partitioning the mesh8function rebalance_solver!(u_ode::AbstractVector, mesh::TreeMesh{2}, equations,9dg::DGSEM, cache, old_mpi_ranks_per_cell)10if cache.elements.cell_ids == local_leaf_cells(mesh.tree)11# Cell ids of the current elements are the same as the local leaf cells of the12# newly partitioned mesh, so the solver doesn't need to be rebalanced on this rank.13# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do14# locally (there are other MPI ranks that need to be rebalanced if this function is called)15reinitialize_containers!(mesh, equations, dg, cache)16return17end1819# Retain current solution data20old_n_elements = nelements(dg, cache)21old_cell_ids = copy(cache.elements.cell_ids)22old_u_ode = copy(u_ode)23GC.@preserve old_u_ode begin # OBS! If we don't GC.@preserve old_u_ode, it might be GC'ed24# Use `wrap_array_native` instead of `wrap_array` since MPI might not interact25# nicely with non-base array types26old_u = wrap_array_native(old_u_ode, mesh, equations, dg, cache)2728@trixi_timeit timer() "reinitialize data structures" begin29reinitialize_containers!(mesh, equations, dg, cache)30end3132resize!(u_ode,33nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))34u = wrap_array_native(u_ode, mesh, equations, dg, cache)3536# Get new list of leaf cells37leaf_cell_ids = local_leaf_cells(mesh.tree)3839@trixi_timeit timer() "exchange data" begin40# Collect MPI requests for MPI_Waitall41requests = Vector{MPI.Request}()4243# Find elements that will change their rank and send their data to the new rank44for old_element_id in 1:old_n_elements45cell_id = old_cell_ids[old_element_id]46if !(cell_id in leaf_cell_ids)47# Send element data to new rank, use cell_id as tag (non-blocking)48dest = mesh.tree.mpi_ranks[cell_id]49request = MPI.Isend(@view(old_u[:, .., old_element_id]), dest,50cell_id, mpi_comm())51push!(requests, request)52end53end5455# Loop over all elements in new container and either copy them from old container56# or receive them with MPI57for element in eachelement(dg, cache)58cell_id = cache.elements.cell_ids[element]59if cell_id in old_cell_ids60old_element_id = searchsortedfirst(old_cell_ids, cell_id)61# Copy old element data to new element container62@views u[:, .., element] .= old_u[:, .., old_element_id]63else64# Receive old element data65src = old_mpi_ranks_per_cell[cell_id]66request = MPI.Irecv!(@view(u[:, .., element]), src, cell_id,67mpi_comm())68push!(requests, request)69end70end7172# Wait for all non-blocking MPI send/receive operations to finish73MPI.Waitall(requests, MPI.Status)74end75end # GC.@preserve old_u_ode7677return nothing78end7980# Refine elements in the DG solver based on a list of cell_ids that should be refined.81# On `P4estMesh`, if an element refines the solution scaled by the Jacobian `J*u` is interpolated82# from the parent element into the four children elements. The solution on each child83# element is then recovered by dividing by the new element Jacobians.84function refine!(u_ode::AbstractVector, adaptor, mesh::Union{TreeMesh{2}, P4estMesh{2}},85equations, dg::DGSEM, cache, elements_to_refine)86# Return early if there is nothing to do87if isempty(elements_to_refine)88if mpi_isparallel()89# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do90# locally (there still might be other MPI ranks that have refined elements)91reinitialize_containers!(mesh, equations, dg, cache)92end93return94end9596# Determine for each existing element whether it needs to be refined97needs_refinement = falses(nelements(dg, cache))98needs_refinement[elements_to_refine] .= true99100# Retain current solution data101old_n_elements = nelements(dg, cache)102old_u_ode = copy(u_ode)103old_inverse_jacobian = copy(cache.elements.inverse_jacobian)104# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed105GC.@preserve old_u_ode old_inverse_jacobian begin106old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)107108if mesh isa P4estMesh109# Loop over all elements in old container and scale the old solution by the Jacobian110# prior to projection111for old_element_id in 1:old_n_elements112for v in eachvariable(equations)113old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./114old_inverse_jacobian[..,115old_element_id])116end117end118end119120reinitialize_containers!(mesh, equations, dg, cache)121122resize!(u_ode,123nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))124u = wrap_array(u_ode, mesh, equations, dg, cache)125126# Loop over all elements in old container and either copy them or refine them127element_id = 1128for old_element_id in 1:old_n_elements129if needs_refinement[old_element_id]130# Refine element and store solution directly in new data structure131refine_element!(u, element_id, old_u, old_element_id,132adaptor, equations, dg)133134if mesh isa P4estMesh135# Before `element_id` is incremented, divide by the new Jacobians on each136# child element and save the result137for m in 0:3 # loop over the children138for v in eachvariable(equations)139u[v, .., element_id + m] .*= (0.25f0 .*140cache.elements.inverse_jacobian[..,141element_id + m])142end143end144end145146# Increment `element_id` on the refined mesh with the number147# of children, i.e., 4 in 2D148element_id += 2^ndims(mesh)149else150if mesh isa P4estMesh151# Copy old element data to new element container and remove Jacobian scaling152for v in eachvariable(equations)153u[v, .., element_id] .= (old_u[v, .., old_element_id] .*154old_inverse_jacobian[..,155old_element_id])156end157else # isa TreeMesh158@views u[:, .., element_id] .= old_u[:, .., old_element_id]159end160# No refinement occurred, so increment `element_id` on the new mesh by one161element_id += 1162end163end164# If everything is correct, we should have processed all elements.165# Depending on whether the last element processed above had to be refined or not,166# the counter `element_id` can have two different values at the end.167@assert element_id ==168nelements(dg, cache) +1691||element_id == nelements(dg, cache) + 2^ndims(mesh) "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"170end # GC.@preserve old_u_ode old_inverse_jacobian171172# Sanity check173if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0 &&174!mpi_isparallel()175@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")176end177178return nothing179end180181function refine!(u_ode::AbstractVector, adaptor,182mesh::Union{TreeMesh{2}, P4estMesh{2}, TreeMesh{3}, P4estMesh{3}},183equations, dg::DGSEM, cache, cache_parabolic,184elements_to_refine)185# Call `refine!` for the hyperbolic part, which does the heavy lifting of186# actually transferring the solution to the refined cells187refine!(u_ode, adaptor, mesh, equations, dg, cache, elements_to_refine)188189# Resize parabolic helper variables190@unpack viscous_container = cache_parabolic191resize!(viscous_container, equations, dg, cache)192reinitialize_containers!(mesh, equations, dg, cache_parabolic)193194# Sanity check195if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0 &&196!mpi_isparallel()197@assert ninterfaces(cache_parabolic.interfaces)==ndims(mesh) *198nelements(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")199end200201return nothing202end203204# TODO: Taal compare performance of different implementations205# Refine solution data u for an element, using L2 projection (interpolation)206function refine_element!(u::AbstractArray{<:Any, 4}, element_id,207old_u, old_element_id,208adaptor::LobattoLegendreAdaptorL2, equations, dg)209@unpack forward_upper, forward_lower = adaptor210211# Store new element ids212lower_left_id = element_id213lower_right_id = element_id + 1214upper_left_id = element_id + 2215upper_right_id = element_id + 3216217@boundscheck begin218@assert old_element_id >= 1219@assert size(old_u, 1) == nvariables(equations)220@assert size(old_u, 2) == nnodes(dg)221@assert size(old_u, 3) == nnodes(dg)222@assert size(old_u, 4) >= old_element_id223@assert element_id >= 1224@assert size(u, 1) == nvariables(equations)225@assert size(u, 2) == nnodes(dg)226@assert size(u, 3) == nnodes(dg)227@assert size(u, 4) >= element_id + 3228end229230# Interpolate to lower left element231for j in eachnode(dg), i in eachnode(dg)232acc = zero(get_node_vars(u, equations, dg, i, j, element_id))233for l in eachnode(dg), k in eachnode(dg)234acc += get_node_vars(old_u, equations, dg, k, l, old_element_id) *235forward_lower[i, k] * forward_lower[j, l]236end237set_node_vars!(u, acc, equations, dg, i, j, lower_left_id)238end239240# Interpolate to lower right element241for j in eachnode(dg), i in eachnode(dg)242acc = zero(get_node_vars(u, equations, dg, i, j, element_id))243for l in eachnode(dg), k in eachnode(dg)244acc += get_node_vars(old_u, equations, dg, k, l, old_element_id) *245forward_upper[i, k] * forward_lower[j, l]246end247set_node_vars!(u, acc, equations, dg, i, j, lower_right_id)248end249250# Interpolate to upper left element251for j in eachnode(dg), i in eachnode(dg)252acc = zero(get_node_vars(u, equations, dg, i, j, element_id))253for l in eachnode(dg), k in eachnode(dg)254acc += get_node_vars(old_u, equations, dg, k, l, old_element_id) *255forward_lower[i, k] * forward_upper[j, l]256end257set_node_vars!(u, acc, equations, dg, i, j, upper_left_id)258end259260# Interpolate to upper right element261for j in eachnode(dg), i in eachnode(dg)262acc = zero(get_node_vars(u, equations, dg, i, j, element_id))263for l in eachnode(dg), k in eachnode(dg)264acc += get_node_vars(old_u, equations, dg, k, l, old_element_id) *265forward_upper[i, k] * forward_upper[j, l]266end267set_node_vars!(u, acc, equations, dg, i, j, upper_right_id)268end269270return nothing271end272273# Coarsen elements in the DG solver based on a list of cell_ids that should be removed.274# On `P4estMesh`, if an element coarsens the solution scaled by the Jacobian `J*u` is projected275# from the four children elements back onto the parent element. The solution on the parent276# element is then recovered by dividing by the new element Jacobian.277function coarsen!(u_ode::AbstractVector, adaptor,278mesh::Union{TreeMesh{2}, P4estMesh{2}},279equations, dg::DGSEM, cache, elements_to_remove)280# Return early if there is nothing to do281if isempty(elements_to_remove)282if mpi_isparallel()283# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do284# locally (there still might be other MPI ranks that have coarsened elements)285reinitialize_containers!(mesh, equations, dg, cache)286end287return288end289290# Determine for each old element whether it needs to be removed291to_be_removed = falses(nelements(dg, cache))292to_be_removed[elements_to_remove] .= true293294# Retain current solution data and Jacobians295old_n_elements = nelements(dg, cache)296old_u_ode = copy(u_ode)297old_inverse_jacobian = copy(cache.elements.inverse_jacobian)298# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed299GC.@preserve old_u_ode old_inverse_jacobian begin300old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)301302if mesh isa P4estMesh303# Loop over all elements in old container and scale the old solution by the Jacobian304# prior to projection305for old_element_id in 1:old_n_elements306for v in eachvariable(equations)307old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./308old_inverse_jacobian[..,309old_element_id])310end311end312end313314reinitialize_containers!(mesh, equations, dg, cache)315316resize!(u_ode,317nvariables(equations) * nnodes(dg)^ndims(mesh) * nelements(dg, cache))318u = wrap_array(u_ode, mesh, equations, dg, cache)319320# Loop over all elements in old container and either copy them or coarsen them321skip = 0322element_id = 1323for old_element_id in 1:old_n_elements324# If skip is non-zero, we just coarsened 2^ndims elements and need to omit the following elements325if skip > 0326skip -= 1327continue328end329330if to_be_removed[old_element_id]331# If an element is to be removed, sanity check if the following elements332# are also marked - otherwise there would be an error in the way the333# cells/elements are sorted334@assert all(to_be_removed[old_element_id:(old_element_id + 2^ndims(mesh) - 1)]) "bad cell/element order"335336# Coarsen elements and store solution directly in new data structure337coarsen_elements!(u, element_id, old_u, old_element_id,338adaptor, equations, dg)339340if mesh isa P4estMesh341# Before `element_id` is incremented, divide by the new Jacobian and save342# the result in the parent element343for v in eachvariable(equations)344u[v, .., element_id] .*= (4 .*345cache.elements.inverse_jacobian[..,346element_id])347end348end349350# Increment `element_id` on the coarsened mesh by one and `skip` = 3 in 2D351# because 4 children elements become 1 parent element352element_id += 1353skip = 2^ndims(mesh) - 1354else355if mesh isa P4estMesh356# Copy old element data to new element container and remove Jacobian scaling357for v in eachvariable(equations)358u[v, .., element_id] .= (old_u[v, .., old_element_id] .*359old_inverse_jacobian[..,360old_element_id])361end362else # isa TreeMesh363@views u[:, .., element_id] .= old_u[:, .., old_element_id]364end365# No coarsening occurred, so increment `element_id` on the new mesh by one366element_id += 1367end368end369# If everything is correct, we should have processed all elements.370@assert element_id==nelements(dg, cache) + 1 "element_id = $element_id, nelements(dg, cache) = $(nelements(dg, cache))"371end # GC.@preserve old_u_ode old_inverse_jacobian372373# Sanity check374if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0 &&375!mpi_isparallel()376@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")377end378379return nothing380end381382function coarsen!(u_ode::AbstractVector, adaptor,383mesh::Union{TreeMesh{2}, P4estMesh{2}, TreeMesh{3}, P4estMesh{3}},384equations, dg::DGSEM, cache, cache_parabolic,385elements_to_remove)386# Call `coarsen!` for the hyperbolic part, which does the heavy lifting of387# actually transferring the solution to the coarsened cells388coarsen!(u_ode, adaptor, mesh, equations, dg, cache, elements_to_remove)389390# Resize parabolic helper variables391@unpack viscous_container = cache_parabolic392resize!(viscous_container, equations, dg, cache)393reinitialize_containers!(mesh, equations, dg, cache_parabolic)394395# Sanity check396if mesh isa TreeMesh && isperiodic(mesh.tree) && nmortars(cache.mortars) == 0 &&397!mpi_isparallel()398@assert ninterfaces(cache_parabolic.interfaces)==ndims(mesh) *399nelements(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")400end401402return nothing403end404405# TODO: Taal compare performance of different implementations406# Coarsen solution data u for four elements, using L2 projection407function coarsen_elements!(u::AbstractArray{<:Any, 4}, element_id,408old_u, old_element_id,409adaptor::LobattoLegendreAdaptorL2, equations, dg)410@unpack reverse_upper, reverse_lower = adaptor411412# Store old element ids413lower_left_id = old_element_id414lower_right_id = old_element_id + 1415upper_left_id = old_element_id + 2416upper_right_id = old_element_id + 3417418@boundscheck begin419@assert old_element_id >= 1420@assert size(old_u, 1) == nvariables(equations)421@assert size(old_u, 2) == nnodes(dg)422@assert size(old_u, 3) == nnodes(dg)423@assert size(old_u, 4) >= old_element_id + 3424@assert element_id >= 1425@assert size(u, 1) == nvariables(equations)426@assert size(u, 2) == nnodes(dg)427@assert size(u, 3) == nnodes(dg)428@assert size(u, 4) >= element_id429end430431for j in eachnode(dg), i in eachnode(dg)432acc = zero(get_node_vars(u, equations, dg, i, j, element_id))433434# Project from lower left element435for l in eachnode(dg), k in eachnode(dg)436acc += get_node_vars(old_u, equations, dg, k, l, lower_left_id) *437reverse_lower[i, k] * reverse_lower[j, l]438end439440# Project from lower right element441for l in eachnode(dg), k in eachnode(dg)442acc += get_node_vars(old_u, equations, dg, k, l, lower_right_id) *443reverse_upper[i, k] * reverse_lower[j, l]444end445446# Project from upper left element447for l in eachnode(dg), k in eachnode(dg)448acc += get_node_vars(old_u, equations, dg, k, l, upper_left_id) *449reverse_lower[i, k] * reverse_upper[j, l]450end451452# Project from upper right element453for l in eachnode(dg), k in eachnode(dg)454acc += get_node_vars(old_u, equations, dg, k, l, upper_right_id) *455reverse_upper[i, k] * reverse_upper[j, l]456end457458# Update value459set_node_vars!(u, acc, equations, dg, i, j, element_id)460end461end462463# Coarsen and refine elements in the DG solver based on a difference list.464function adapt!(u_ode::AbstractVector, adaptor, mesh::T8codeMesh{2}, equations,465dg::DGSEM, cache, difference)466467# Return early if there is nothing to do.468if !any(difference .!= 0)469if mpi_isparallel()470# MPICache init uses all-to-all communication -> reinitialize even if there is nothing to do471# locally (there still might be other MPI ranks that have refined elements)472reinitialize_containers!(mesh, equations, dg, cache)473end474return475end476477# Number of (local) cells/elements.478old_nelems = nelements(dg, cache)479new_nelems = ncells(mesh)480481# Local element indices.482old_index = 1483new_index = 1484485# Note: This is true for `quads`.486T8_CHILDREN = 4487488# Retain current solution and inverse Jacobian data.489old_u_ode = copy(u_ode)490old_inverse_jacobian = copy(cache.elements.inverse_jacobian)491492# OBS! If we don't GC.@preserve old_u_ode and old_inverse_jacobian, they might be GC'ed493GC.@preserve old_u_ode begin494old_u = wrap_array(old_u_ode, mesh, equations, dg, cache)495496# Loop over all elements in old container and scale the old solution by the Jacobian497# prior to interpolation or projection498for old_element_id in 1:old_nelems499for v in eachvariable(equations)500old_u[v, .., old_element_id] .= (old_u[v, .., old_element_id] ./501old_inverse_jacobian[..,502old_element_id])503end504end505506reinitialize_containers!(mesh, equations, dg, cache)507508resize!(u_ode, nvariables(equations) * ndofs(mesh, dg, cache))509u = wrap_array(u_ode, mesh, equations, dg, cache)510511while old_index <= old_nelems && new_index <= new_nelems512if difference[old_index] > 0 # Refine.513514# Refine element and store solution directly in new data structure.515refine_element!(u, new_index, old_u, old_index, adaptor, equations, dg)516517# Before indices are incremented divide by the new Jacobians on each518# child element and save the result519for m in 0:3 # loop over the children520for v in eachvariable(equations)521u[v, .., new_index + m] .*= (0.25f0 .*522cache.elements.inverse_jacobian[..,523new_index + m])524end525end526527# Increment `old_index` on the original mesh and the `new_index`528# on the refined mesh with the number of children, i.e., T8_CHILDREN = 4529old_index += 1530new_index += T8_CHILDREN531532elseif difference[old_index] < 0 # Coarsen.533534# If an element is to be removed, sanity check if the following elements535# are also marked - otherwise there would be an error in the way the536# cells/elements are sorted.537@assert all(difference[old_index:(old_index + T8_CHILDREN - 1)] .< 0) "bad cell/element order"538539# Coarsen elements and store solution directly in new data structure.540coarsen_elements!(u, new_index, old_u, old_index, adaptor, equations,541dg)542543# Before the indices are incremented divide by the new Jacobian and save544# the result again in the parent element545for v in eachvariable(equations)546u[v, .., new_index] .*= (4 .* cache.elements.inverse_jacobian[..,547new_index])548end549550# Increment `old_index` on the original mesh with the number of children551# (T8_CHILDREN = 4 in 2D) and the `new_index` by one for the single552# coarsened element553old_index += T8_CHILDREN554new_index += 1555556else # No changes.557558# Copy old element data to new element container and remove Jacobian scaling559for v in eachvariable(equations)560u[v, .., new_index] .= (old_u[v, .., old_index] .*561old_inverse_jacobian[.., old_index])562end563564# No refinement / coarsening occurred, so increment element index565# on each mesh by one566old_index += 1567new_index += 1568end569end # while570end # GC.@preserve old_u_ode old_inverse_jacobian571572return nothing573end574end # @muladd575576577