Path: blob/main/src/solvers/dgsem_tree/indicators_1d.jl
5591 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# this method is directly used when the indicator is constructed as for shock-capturing volume integrals8# and by the dimension-independent method called for AMR9function create_cache(::Type{IndicatorHennemannGassner},10equations::AbstractEquations{1}, basis::LobattoLegendreBasis)11uEltype = real(basis)12alpha = Vector{uEltype}()13alpha_tmp = similar(alpha)1415MVec = MVector{nnodes(basis), uEltype}1617indicator_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]18modal_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]1920return (; alpha, alpha_tmp, indicator_threaded, modal_threaded)21end2223# Use this function barrier and unpack inside to avoid passing closures to Polyester.jl24# with @batch (@threaded).25# Otherwise, @threaded does not work here with Julia ARM on macOS.26# See https://github.com/JuliaSIMD/Polyester.jl/issues/88.27@inline function calc_indicator_hennemann_gassner!(indicator_hg, threshold, parameter_s,28u, element, mesh::AbstractMesh{1},29equations, dg, cache)30@unpack alpha_max, alpha_min, alpha_smooth, variable = indicator_hg31@unpack alpha, alpha_tmp, indicator_threaded, modal_threaded = indicator_hg.cache3233indicator = indicator_threaded[Threads.threadid()]34modal = modal_threaded[Threads.threadid()]3536# Calculate indicator variables at Gauss-Lobatto nodes37for i in eachnode(dg)38u_local = get_node_vars(u, equations, dg, i, element)39indicator[i] = indicator_hg.variable(u_local, equations)40end4142# Convert to modal representation43multiply_scalar_dimensionwise!(modal, dg.basis.inverse_vandermonde_legendre,44indicator)4546# Calculate total energies without two highest, without highest, and for all modes47total_energy_clip2 = zero(eltype(modal))48for i in 1:(nnodes(dg) - 2)49total_energy_clip2 += modal[i]^250end51total_energy_clip1 = total_energy_clip2 + modal[nnodes(dg) - 1]^252total_energy = total_energy_clip1 + modal[nnodes(dg)]^25354# Calculate energy in higher modes55if !(iszero(total_energy))56energy_frac_1 = (total_energy - total_energy_clip1) / total_energy57else58energy_frac_1 = zero(total_energy)59end60if !(iszero(total_energy_clip1))61energy_frac_2 = (total_energy_clip1 - total_energy_clip2) / total_energy_clip162else63energy_frac_2 = zero(total_energy_clip1)64end65energy = max(energy_frac_1, energy_frac_2)6667alpha_element = 1 / (1 + exp(-parameter_s / threshold * (energy - threshold)))6869# Take care of the case close to pure DG70if alpha_element < alpha_min71alpha_element = zero(alpha_element)72end7374# Take care of the case close to pure FV75if alpha_element > 1 - alpha_min76alpha_element = one(alpha_element)77end7879# Clip the maximum amount of FV allowed80alpha[element] = min(alpha_max, alpha_element)81return nothing82end8384# Diffuse alpha values by setting each alpha to at least 50% of neighboring elements' alpha85function apply_smoothing!(mesh::TreeMesh{1}, alpha, alpha_tmp, dg, cache)86# Copy alpha values such that smoothing is indpedenent of the element access order87alpha_tmp .= alpha8889# Loop over interfaces90for interface in eachinterface(dg, cache)91# Get neighboring element ids92left = cache.interfaces.neighbor_ids[1, interface]93right = cache.interfaces.neighbor_ids[2, interface]9495# Apply smoothing96alpha[left] = max(alpha_tmp[left], 0.5f0 * alpha_tmp[right], alpha[left])97alpha[right] = max(alpha_tmp[right], 0.5f0 * alpha_tmp[left], alpha[right])98end99100return nothing101end102103# this method is directly used when the indicator is constructed as for shock-capturing volume integrals104# and by the dimension-independent method called for AMR105function create_cache(::Union{Type{IndicatorLöhner}, Type{IndicatorMax}},106equations::AbstractEquations{1}, basis::LobattoLegendreBasis)107uEltype = real(basis)108alpha = Vector{uEltype}()109110MVec = MVector{nnodes(basis), uEltype}111112indicator_threaded = MVec[MVec(undef) for _ in 1:Threads.maxthreadid()]113114return (; alpha, indicator_threaded)115end116117function (löhner::IndicatorLöhner)(u::AbstractArray{<:Any, 3},118mesh, equations, dg::DGSEM, cache;119kwargs...)120@assert nnodes(dg)>=3 "IndicatorLöhner only works for nnodes >= 3 (polydeg > 1)"121@unpack alpha, indicator_threaded = löhner.cache122@unpack variable = löhner123resize!(alpha, nelements(dg, cache))124125@threaded for element in eachelement(dg, cache)126indicator = indicator_threaded[Threads.threadid()]127128# Calculate indicator variables at Gauss-Lobatto nodes129for i in eachnode(dg)130u_local = get_node_vars(u, equations, dg, i, element)131indicator[i] = variable(u_local, equations)132end133134estimate = zero(real(dg))135for i in 2:(nnodes(dg) - 1)136# x direction137u0 = indicator[i]138up = indicator[i + 1]139um = indicator[i - 1]140estimate = max(estimate, local_löhner_estimate(um, u0, up, löhner))141end142143# use the maximum as DG element indicator144alpha[element] = estimate145end146147return alpha148end149150function (indicator_max::IndicatorMax)(u::AbstractArray{<:Any, 3},151mesh, equations, dg::DGSEM, cache;152kwargs...)153@unpack alpha, indicator_threaded = indicator_max.cache154resize!(alpha, nelements(dg, cache))155indicator_variable = indicator_max.variable156157@threaded for element in eachelement(dg, cache)158indicator = indicator_threaded[Threads.threadid()]159160# Calculate indicator variables at Gauss-Lobatto nodes161for i in eachnode(dg)162u_local = get_node_vars(u, equations, dg, i, element)163indicator[i] = indicator_variable(u_local, equations)164end165166alpha[element] = maximum(indicator)167end168169return alpha170end171172function (indicator::IndicatorNodalFunction)(u::AbstractArray{<:Any, 3},173mesh, equations, dg::DGSEM, cache;174t, kwargs...)175node_coordinates = cache.elements.node_coordinates176@unpack alpha = indicator.cache177resize!(alpha, nelements(dg, cache))178# Extract function to local variable to avoid capturing `indicator` in the threaded loop179indicator_function = indicator.indicator_function180181@threaded for element in eachelement(dg, cache)182estimate = typemin(eltype(alpha))183for i in eachnode(dg)184u_nodal = get_node_vars(u, equations, dg, i, element)185x_nodal = get_node_coords(node_coordinates, equations, dg,186i, element)187estimate = max(estimate, indicator_function(u_nodal, x_nodal, t))188end189alpha[element] = estimate190end191192return alpha193end194end # @muladd195196197