Path: blob/main/src/callbacks_step/analysis_dg1d.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: noindent67function create_cache_analysis(analyzer, mesh::TreeMesh{1},8equations, dg::DG, cache,9RealT, uEltype)1011# pre-allocate buffers12# We use `StrideArray`s here since these buffers are used in performance-critical13# places and the additional information passed to the compiler makes them faster14# than native `Array`s.15u_local = StrideArray(undef, uEltype,16StaticInt(nvariables(equations)), StaticInt(nnodes(analyzer)))17x_local = StrideArray(undef, RealT,18StaticInt(ndims(equations)), StaticInt(nnodes(analyzer)))1920return (; u_local, x_local)21end2223function create_cache_analysis(analyzer, mesh::StructuredMesh{1},24equations, dg::DG, cache,25RealT, uEltype)2627# pre-allocate buffers28# We use `StrideArray`s here since these buffers are used in performance-critical29# places and the additional information passed to the compiler makes them faster30# than native `Array`s.31u_local = StrideArray(undef, uEltype,32StaticInt(nvariables(equations)), StaticInt(nnodes(analyzer)))33x_local = StrideArray(undef, RealT,34StaticInt(ndims(equations)), StaticInt(nnodes(analyzer)))35jacobian_local = StrideArray(undef, RealT,36StaticInt(nnodes(analyzer)))3738return (; u_local, x_local, jacobian_local)39end4041function calc_error_norms(func, u, t, analyzer,42mesh::TreeMesh{1}, equations, initial_condition,43dg::DGSEM, cache, cache_analysis)44@unpack vandermonde, weights = analyzer45@unpack node_coordinates = cache.elements46@unpack u_local, x_local = cache_analysis4748# Set up data structures49l2_error = zero(func(get_node_vars(u, equations, dg, 1, 1), equations))50linf_error = copy(l2_error)5152# Iterate over all elements for error calculations53for element in eachelement(dg, cache)54# Interpolate solution and node locations to analysis nodes55multiply_dimensionwise!(u_local, vandermonde, view(u, :, :, element))56multiply_dimensionwise!(x_local, vandermonde,57view(node_coordinates, :, :, element))5859# Calculate errors at each analysis node60volume_jacobian_ = volume_jacobian(element, mesh, cache)6162for i in eachnode(analyzer)63u_exact = initial_condition(get_node_coords(x_local, equations, dg, i), t,64equations)65diff = func(u_exact, equations) -66func(get_node_vars(u_local, equations, dg, i), equations)67l2_error += diff .^ 2 * (weights[i] * volume_jacobian_)68linf_error = @. max(linf_error, abs(diff))69end70end7172# For L2 error, divide by total volume73total_volume_ = total_volume(mesh)74l2_error = @. sqrt(l2_error / total_volume_)7576return l2_error, linf_error77end7879function calc_error_norms(func, u, t, analyzer,80mesh::StructuredMesh{1}, equations, initial_condition,81dg::DGSEM, cache, cache_analysis)82@unpack vandermonde, weights = analyzer83@unpack node_coordinates, inverse_jacobian = cache.elements84@unpack u_local, x_local, jacobian_local = cache_analysis8586# Set up data structures87l2_error = zero(func(get_node_vars(u, equations, dg, 1, 1), equations))88linf_error = copy(l2_error)89total_volume = zero(real(mesh))9091# Iterate over all elements for error calculations92for element in eachelement(dg, cache)93# Interpolate solution and node locations to analysis nodes94multiply_dimensionwise!(u_local, vandermonde, view(u, :, :, element))95multiply_dimensionwise!(x_local, vandermonde,96view(node_coordinates, :, :, element))97multiply_scalar_dimensionwise!(jacobian_local, vandermonde,98inv.(view(inverse_jacobian, :, element)))99100# Calculate errors at each analysis node101for i in eachnode(analyzer)102u_exact = initial_condition(get_node_coords(x_local, equations, dg, i), t,103equations)104diff = func(u_exact, equations) -105func(get_node_vars(u_local, equations, dg, i), equations)106# We take absolute value as we need the Jacobian here for the volume calculation107abs_jacobian_local_i = abs(jacobian_local[i])108109l2_error += diff .^ 2 * (weights[i] * abs_jacobian_local_i)110linf_error = @. max(linf_error, abs(diff))111total_volume += weights[i] * abs_jacobian_local_i112end113end114115# For L2 error, divide by total volume116l2_error = @. sqrt(l2_error / total_volume)117118return l2_error, linf_error119end120121# Use quadrature to numerically integrate a single element.122# We do not multiply by the Jacobian to stay in reference space.123# This avoids the need to divide the RHS of the DG scheme by the Jacobian when computing124# the time derivative of entropy, see `entropy_change_reference_element`.125function integrate_reference_element(func::Func, u, element,126::Type{<:AbstractMesh{1}}, equations, dg::DGSEM,127cache,128args...) where {Func}129@unpack weights = dg.basis130131# Initialize integral with zeros of the right shape132element_integral = zero(func(u, 1, element, equations, dg, args...))133134for i in eachnode(dg)135element_integral += weights[i] *136func(u, i, element, equations, dg, args...)137end138139return element_integral140end141142# Calculate ∫_e (∂S/∂u ⋅ ∂u/∂t) dΩ_e where the result on element 'e' is kept in reference space143# Note that ∂S/∂u = w(u) with entropy variables w144function entropy_change_reference_element(du, u, element,145MeshT::Type{<:AbstractMesh{1}},146equations, dg::DGSEM, cache, args...)147return integrate_reference_element(u, element, MeshT, equations, dg, cache,148du) do u, i, element, equations, dg, du149u_node = get_node_vars(u, equations, dg, i, element)150du_node = get_node_vars(du, equations, dg, i, element)151152dot(cons2entropy(u_node, equations), du_node)153end154end155156# calculate surface integral of func(u, equations) * normal on the reference element.157function surface_integral_reference_element(func::Func, u, element,158::Type{<:Union{TreeMesh{1},159StructuredMesh{1}}},160equations, dg::DGSEM,161cache, args...) where {Func}162u_left = get_node_vars(u, equations, dg, 1, element)163u_right = get_node_vars(u, equations, dg, nnodes(dg), element)164165surface_integral = func(u_right, 1, equations) - func(u_left, 1, equations)166return surface_integral167end168169function integrate_via_indices(func::Func, u,170mesh::TreeMesh{1}, equations, dg::DGSEM, cache,171args...; normalize = true) where {Func}172@unpack weights = dg.basis173174# Initialize integral with zeros of the right shape175integral = zero(func(u, 1, 1, equations, dg, args...))176177# Use quadrature to numerically integrate over entire domain178@batch reduction=(+, integral) for element in eachelement(dg, cache)179volume_jacobian_ = volume_jacobian(element, mesh, cache)180for i in eachnode(dg)181integral += volume_jacobian_ * weights[i] *182func(u, i, element, equations, dg, args...)183end184end185186# Normalize with total volume187if normalize188integral = integral / total_volume(mesh)189end190191return integral192end193194function integrate_via_indices(func::Func, u,195mesh::StructuredMesh{1}, equations, dg::DGSEM, cache,196args...; normalize = true) where {Func}197@unpack weights = dg.basis198199# Initialize integral with zeros of the right shape200integral = zero(func(u, 1, 1, equations, dg, args...))201total_volume = zero(real(mesh))202203# Use quadrature to numerically integrate over entire domain204@batch reduction=((+, integral), (+, total_volume)) for element in eachelement(dg,205cache)206for i in eachnode(dg)207jacobian_volume = abs(inv(cache.elements.inverse_jacobian[i, element]))208integral += jacobian_volume * weights[i] *209func(u, i, element, equations, dg, args...)210total_volume += jacobian_volume * weights[i]211end212end213# Normalize with total volume214if normalize215integral = integral / total_volume216end217218return integral219end220221function integrate(func::Func, u,222mesh::Union{TreeMesh{1}, StructuredMesh{1}},223equations, dg::Union{DGSEM, FDSBP}, cache;224normalize = true) where {Func}225integrate_via_indices(u, mesh, equations, dg, cache;226normalize = normalize) do u, i, element, equations, dg227u_local = get_node_vars(u, equations, dg, i, element)228return func(u_local, equations)229end230end231232function analyze(::typeof(entropy_timederivative), du, u, t,233mesh::Union{TreeMesh{1}, StructuredMesh{1}}, equations,234dg::Union{DGSEM, FDSBP}, cache)235# Calculate ∫(∂S/∂u ⋅ ∂u/∂t)dΩ236integrate_via_indices(u, mesh, equations, dg, cache,237du) do u, i, element, equations, dg, du238u_node = get_node_vars(u, equations, dg, i, element)239du_node = get_node_vars(du, equations, dg, i, element)240return dot(cons2entropy(u_node, equations), du_node)241end242end243244function analyze(::Val{:l2_divb}, du, u, t,245mesh::TreeMesh{1}, equations::IdealGlmMhdEquations1D,246dg::DGSEM, cache)247integrate_via_indices(u, mesh, equations, dg, cache,248dg.basis.derivative_matrix) do u, i, element, equations, dg,249derivative_matrix250divb = zero(eltype(u))251for k in eachnode(dg)252divb += derivative_matrix[i, k] * u[6, k, element]253end254divb *= cache.elements.inverse_jacobian[element]255return divb^2256end |> sqrt257end258259function analyze(::Val{:linf_divb}, du, u, t,260mesh::TreeMesh{1}, equations::IdealGlmMhdEquations1D,261dg::DGSEM, cache)262@unpack derivative_matrix, weights = dg.basis263264# integrate over all elements to get the divergence-free condition errors265linf_divb = zero(eltype(u))266@batch reduction=(max, linf_divb) for element in eachelement(dg, cache)267for i in eachnode(dg)268divb = zero(eltype(u))269for k in eachnode(dg)270divb += derivative_matrix[i, k] * u[6, k, element]271end272divb *= cache.elements.inverse_jacobian[element]273linf_divb = max(linf_divb, abs(divb))274end275end276277return linf_divb278end279end # @muladd280281282