Path: blob/main/src/solvers/fdsbp_tree/fdsbp_2d.jl
5591 views
# !!! warning "Experimental implementation (upwind SBP)"1# This is an experimental feature and may change in future releases.23# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).4# Since these FMAs can increase the performance of many numerical algorithms,5# we need to opt-in explicitly.6# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.7@muladd begin8#! format: noindent910# 2D caches11function create_cache(mesh::Union{TreeMesh{2}, UnstructuredMesh2D}, equations,12volume_integral::VolumeIntegralStrongForm,13dg, cache_containers, uEltype)14prototype = Array{SVector{nvariables(equations), uEltype}, ndims(mesh)}(undef,15ntuple(_ -> nnodes(dg),16ndims(mesh))...)17f_threaded = [similar(prototype) for _ in 1:Threads.maxthreadid()]1819return (; f_threaded)20end2122function create_cache(mesh::Union{TreeMesh{2}, UnstructuredMesh2D}, equations,23volume_integral::VolumeIntegralUpwind,24dg, cache_containers, uEltype)25u_node = SVector{nvariables(equations), uEltype}(ntuple(_ -> zero(uEltype),26Val{nvariables(equations)}()))27f = StructArray([(u_node, u_node)])28f_minus_plus_threaded = [similar(f, ntuple(_ -> nnodes(dg), ndims(mesh))...)29for _ in 1:Threads.maxthreadid()]3031f_minus, f_plus = StructArrays.components(f_minus_plus_threaded[1])32f_minus_threaded = [f_minus]33f_plus_threaded = [f_plus]34for i in 2:Threads.maxthreadid()35f_minus, f_plus = StructArrays.components(f_minus_plus_threaded[i])36push!(f_minus_threaded, f_minus)37push!(f_plus_threaded, f_plus)38end3940return (; f_minus_plus_threaded, f_minus_threaded, f_plus_threaded)41end4243# 2D volume integral contributions for `VolumeIntegralStrongForm`44function calc_volume_integral!(backend::Nothing, du, u,45mesh::TreeMesh{2},46have_nonconservative_terms::False, equations,47volume_integral::VolumeIntegralStrongForm,48dg::FDSBP, cache)49D = dg.basis # SBP derivative operator50@unpack f_threaded = cache5152# SBP operators from SummationByPartsOperators.jl implement the basic interface53# of matrix-vector multiplication. Thus, we pass an "array of structures",54# packing all variables per node in an `SVector`.55if nvariables(equations) == 156# `reinterpret(reshape, ...)` removes the leading dimension only if more57# than one variable is used.58u_vectors = reshape(reinterpret(SVector{nvariables(equations), eltype(u)}, u),59nnodes(dg), nnodes(dg), nelements(dg, cache))60du_vectors = reshape(reinterpret(SVector{nvariables(equations), eltype(du)},61du),62nnodes(dg), nnodes(dg), nelements(dg, cache))63else64u_vectors = reinterpret(reshape, SVector{nvariables(equations), eltype(u)}, u)65du_vectors = reinterpret(reshape, SVector{nvariables(equations), eltype(du)},66du)67end6869# Use the tensor product structure to compute the discrete derivatives of70# the fluxes line-by-line and add them to `du` for each element.71@threaded for element in eachelement(dg, cache)72f_element = f_threaded[Threads.threadid()]73u_element = view(u_vectors, :, :, element)7475# x direction76@. f_element = flux(u_element, 1, equations)77for j in eachnode(dg)78mul!(view(du_vectors, :, j, element), D, view(f_element, :, j),79one(eltype(du)), one(eltype(du)))80end8182# y direction83@. f_element = flux(u_element, 2, equations)84for i in eachnode(dg)85mul!(view(du_vectors, i, :, element), D, view(f_element, i, :),86one(eltype(du)), one(eltype(du)))87end88end8990return nothing91end9293# 2D volume integral contributions for `VolumeIntegralUpwind`.94# Note that the plus / minus notation of the operators does not refer to the95# upwind / downwind directions of the fluxes.96# Instead, the plus / minus refers to the direction of the biasing within97# the finite difference stencils. Thus, the D^- operator acts on the positive98# part of the flux splitting f^+ and the D^+ operator acts on the negative part99# of the flux splitting f^-.100function calc_volume_integral!(backend::Nothing, du, u,101mesh::TreeMesh{2},102have_nonconservative_terms::False, equations,103volume_integral::VolumeIntegralUpwind,104dg::FDSBP, cache)105# Assume that106# dg.basis isa SummationByPartsOperators.UpwindOperators107D_minus = dg.basis.minus # Upwind SBP D^- derivative operator108D_plus = dg.basis.plus # Upwind SBP D^+ derivative operator109@unpack f_minus_plus_threaded, f_minus_threaded, f_plus_threaded = cache110@unpack splitting = volume_integral111112# SBP operators from SummationByPartsOperators.jl implement the basic interface113# of matrix-vector multiplication. Thus, we pass an "array of structures",114# packing all variables per node in an `SVector`.115if nvariables(equations) == 1116# `reinterpret(reshape, ...)` removes the leading dimension only if more117# than one variable is used.118u_vectors = reshape(reinterpret(SVector{nvariables(equations), eltype(u)}, u),119nnodes(dg), nnodes(dg), nelements(dg, cache))120du_vectors = reshape(reinterpret(SVector{nvariables(equations), eltype(du)},121du),122nnodes(dg), nnodes(dg), nelements(dg, cache))123else124u_vectors = reinterpret(reshape, SVector{nvariables(equations), eltype(u)}, u)125du_vectors = reinterpret(reshape, SVector{nvariables(equations), eltype(du)},126du)127end128129# Use the tensor product structure to compute the discrete derivatives of130# the fluxes line-by-line and add them to `du` for each element.131@threaded for element in eachelement(dg, cache)132# f_minus_plus_element wraps the storage provided by f_minus_element and133# f_plus_element such that we can use a single plain broadcasting below.134# f_minus_element and f_plus_element are updated in broadcasting calls135# of the form `@. f_minus_plus_element = ...`.136f_minus_plus_element = f_minus_plus_threaded[Threads.threadid()]137f_minus_element = f_minus_threaded[Threads.threadid()]138f_plus_element = f_plus_threaded[Threads.threadid()]139u_element = view(u_vectors, :, :, element)140141# x direction142@. f_minus_plus_element = splitting(u_element, 1, equations)143for j in eachnode(dg)144mul!(view(du_vectors, :, j, element), D_minus, view(f_plus_element, :, j),145one(eltype(du)), one(eltype(du)))146mul!(view(du_vectors, :, j, element), D_plus, view(f_minus_element, :, j),147one(eltype(du)), one(eltype(du)))148end149150# y direction151@. f_minus_plus_element = splitting(u_element, 2, equations)152for i in eachnode(dg)153mul!(view(du_vectors, i, :, element), D_minus, view(f_plus_element, i, :),154one(eltype(du)), one(eltype(du)))155mul!(view(du_vectors, i, :, element), D_plus, view(f_minus_element, i, :),156one(eltype(du)), one(eltype(du)))157end158end159160return nothing161end162163function calc_surface_integral!(backend::Nothing, du, u, mesh::TreeMesh{2},164equations, surface_integral::SurfaceIntegralStrongForm,165dg::DG, cache)166inv_weight_left = inv(left_boundary_weight(dg.basis))167inv_weight_right = inv(right_boundary_weight(dg.basis))168@unpack surface_flux_values = cache.elements169170@threaded for element in eachelement(dg, cache)171for l in eachnode(dg)172# surface at -x173u_node = get_node_vars(u, equations, dg, 1, l, element)174f_node = flux(u_node, 1, equations)175f_num = get_node_vars(surface_flux_values, equations, dg, l, 1, element)176multiply_add_to_node_vars!(du, inv_weight_left, -(f_num - f_node),177equations, dg, 1, l, element)178179# surface at +x180u_node = get_node_vars(u, equations, dg, nnodes(dg), l, element)181f_node = flux(u_node, 1, equations)182f_num = get_node_vars(surface_flux_values, equations, dg, l, 2, element)183multiply_add_to_node_vars!(du, inv_weight_right, +(f_num - f_node),184equations, dg, nnodes(dg), l, element)185186# surface at -y187u_node = get_node_vars(u, equations, dg, l, 1, element)188f_node = flux(u_node, 2, equations)189f_num = get_node_vars(surface_flux_values, equations, dg, l, 3, element)190multiply_add_to_node_vars!(du, inv_weight_left, -(f_num - f_node),191equations, dg, l, 1, element)192193# surface at +y194u_node = get_node_vars(u, equations, dg, l, nnodes(dg), element)195f_node = flux(u_node, 2, equations)196f_num = get_node_vars(surface_flux_values, equations, dg, l, 4, element)197multiply_add_to_node_vars!(du, inv_weight_right, +(f_num - f_node),198equations, dg, l, nnodes(dg), element)199end200end201202return nothing203end204205# Periodic FDSBP operators need to use a single element without boundaries206function calc_surface_integral!(backend::Nothing, du, u, mesh::TreeMesh2D,207equations, surface_integral::SurfaceIntegralStrongForm,208dg::PeriodicFDSBP, cache)209@assert nelements(dg, cache) == 1210return nothing211end212213# Specialized interface flux computation because the upwind solver does214# not require a standard numerical flux (Riemann solver). The flux splitting215# already separates the solution information into right-traveling and216# left-traveling information. So we only need to compute the appropriate217# flux information at each side of an interface.218function calc_interface_flux!(backend::Nothing, surface_flux_values,219mesh::TreeMesh{2},220have_nonconservative_terms::False, equations,221surface_integral::SurfaceIntegralUpwind,222dg::FDSBP, cache)223@unpack splitting = surface_integral224@unpack u, neighbor_ids, orientations = cache.interfaces225226@threaded for interface in eachinterface(dg, cache)227# Get neighboring elements228left_id = neighbor_ids[1, interface]229right_id = neighbor_ids[2, interface]230231# Determine interface direction with respect to elements:232# orientation = 1: left -> 2, right -> 1233# orientation = 2: left -> 4, right -> 3234left_direction = 2 * orientations[interface]235right_direction = 2 * orientations[interface] - 1236237for i in eachnode(dg)238# Pull the left and right solution data239u_ll, u_rr = get_surface_node_vars(u, equations, dg, i, interface)240241# Compute the upwind coupling terms where right-traveling242# information comes from the left and left-traveling information243# comes from the right244flux_minus_rr = splitting(u_rr, Val{:minus}(), orientations[interface],245equations)246flux_plus_ll = splitting(u_ll, Val{:plus}(), orientations[interface],247equations)248249# Save the upwind coupling into the appropriate side of the elements250for v in eachvariable(equations)251surface_flux_values[v, i, left_direction, left_id] = flux_minus_rr[v]252surface_flux_values[v, i, right_direction, right_id] = flux_plus_ll[v]253end254end255end256257return nothing258end259260# Implementation of fully upwind SATs. The surface flux values are pre-computed261# in the specialized `calc_interface_flux` routine. These SATs are still of262# a strong form penalty type, except that the interior flux at a particular263# side of the element are computed in the upwind direction.264function calc_surface_integral!(backend::Nothing, du, u, mesh::TreeMesh{2},265equations, surface_integral::SurfaceIntegralUpwind,266dg::FDSBP, cache)267inv_weight_left = inv(left_boundary_weight(dg.basis))268inv_weight_right = inv(right_boundary_weight(dg.basis))269@unpack surface_flux_values = cache.elements270@unpack splitting = surface_integral271272@threaded for element in eachelement(dg, cache)273for l in eachnode(dg)274# surface at -x275u_node = get_node_vars(u, equations, dg, 1, l, element)276f_node = splitting(u_node, Val{:plus}(), 1, equations)277f_num = get_node_vars(surface_flux_values, equations, dg, l, 1, element)278multiply_add_to_node_vars!(du, inv_weight_left, -(f_num - f_node),279equations, dg, 1, l, element)280281# surface at +x282u_node = get_node_vars(u, equations, dg, nnodes(dg), l, element)283f_node = splitting(u_node, Val{:minus}(), 1, equations)284f_num = get_node_vars(surface_flux_values, equations, dg, l, 2, element)285multiply_add_to_node_vars!(du, inv_weight_right, +(f_num - f_node),286equations, dg, nnodes(dg), l, element)287288# surface at -y289u_node = get_node_vars(u, equations, dg, l, 1, element)290f_node = splitting(u_node, Val{:plus}(), 2, equations)291f_num = get_node_vars(surface_flux_values, equations, dg, l, 3, element)292multiply_add_to_node_vars!(du, inv_weight_left, -(f_num - f_node),293equations, dg, l, 1, element)294295# surface at +y296u_node = get_node_vars(u, equations, dg, l, nnodes(dg), element)297f_node = splitting(u_node, Val{:minus}(), 2, equations)298f_num = get_node_vars(surface_flux_values, equations, dg, l, 4, element)299multiply_add_to_node_vars!(du, inv_weight_right, +(f_num - f_node),300equations, dg, l, nnodes(dg), element)301end302end303304return nothing305end306307# Periodic FDSBP operators need to use a single element without boundaries308function calc_surface_integral!(backend::Nothing, du, u, mesh::TreeMesh2D,309equations, surface_integral::SurfaceIntegralUpwind,310dg::PeriodicFDSBP, cache)311@assert nelements(dg, cache) == 1312return nothing313end314315# AnalysisCallback316function integrate_via_indices(func::Func, u,317mesh::TreeMesh{2}, equations,318dg::FDSBP, cache, args...; normalize = true) where {Func}319# TODO: FD. This is rather inefficient right now and allocates...320M = SummationByPartsOperators.mass_matrix(dg.basis)321if M isa UniformScaling322M = M(nnodes(dg))323end324weights = diag(M)325326# Initialize integral with zeros of the right shape327integral = zero(func(u, 1, 1, 1, equations, dg, args...))328329# Use quadrature to numerically integrate over entire domain330@batch reduction=(+, integral) for element in eachelement(dg, cache)331volume_jacobian_ = volume_jacobian(element, mesh, cache)332for j in eachnode(dg), i in eachnode(dg)333integral += volume_jacobian_ * weights[i] * weights[j] *334func(u, i, j, element, equations, dg, args...)335end336end337338# Normalize with total volume339if normalize340integral = integral / total_volume(mesh)341end342343return integral344end345346function calc_error_norms(func, u, t, analyzer,347mesh::TreeMesh{2}, equations, initial_condition,348dg::FDSBP, cache, cache_analysis)349# TODO: FD. This is rather inefficient right now and allocates...350M = SummationByPartsOperators.mass_matrix(dg.basis)351if M isa UniformScaling352M = M(nnodes(dg))353end354weights = diag(M)355@unpack node_coordinates = cache.elements356357# Set up data structures358l2_error = zero(func(get_node_vars(u, equations, dg, 1, 1, 1), equations))359linf_error = copy(l2_error)360361# Iterate over all elements for error calculations362for element in eachelement(dg, cache)363# Calculate errors at each node364volume_jacobian_ = volume_jacobian(element, mesh, cache)365366for j in eachnode(analyzer), i in eachnode(analyzer)367u_exact = initial_condition(get_node_coords(node_coordinates, equations, dg,368i, j, element), t, equations)369diff = func(u_exact, equations) -370func(get_node_vars(u, equations, dg, i, j, element), equations)371l2_error += diff .^ 2 * (weights[i] * weights[j] * volume_jacobian_)372linf_error = @. max(linf_error, abs(diff))373end374end375376# For L2 error, divide by total volume377total_volume_ = total_volume(mesh)378l2_error = @. sqrt(l2_error / total_volume_)379380return l2_error, linf_error381end382end # @muladd383384385