Path: blob/main/src/callbacks_step/analysis_dgmulti.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: noindent67function calc_error_norms(func, u, t, analyzer,8mesh::DGMultiMesh{NDIMS}, equations, initial_condition,9dg::DGMulti{NDIMS}, cache, cache_analysis) where {NDIMS}10rd = dg.basis11md = mesh.md12@unpack u_values = cache1314# interpolate u to quadrature points15apply_to_each_field(mul_by!(rd.Vq), u_values, u)1617component_l2_errors = zero(eltype(u_values))18component_linf_errors = zero(eltype(u_values))19for i in each_quad_node_global(mesh, dg, cache)20u_exact = initial_condition(SVector(getindex.(md.xyzq, i)), t, equations)21error_at_node = func(u_values[i], equations) - func(u_exact, equations)22component_l2_errors += md.wJq[i] * error_at_node .^ 223component_linf_errors = max.(component_linf_errors, abs.(error_at_node))24end25total_volume = sum(md.wJq)26return sqrt.(component_l2_errors ./ total_volume), component_linf_errors27end2829function integrate(func::Func, u,30mesh::DGMultiMesh,31equations, dg::DGMulti, cache; normalize = true) where {Func}32rd = dg.basis33md = mesh.md34@unpack u_values = cache3536# interpolate u to quadrature points37apply_to_each_field(mul_by!(rd.Vq), u_values, u)3839integral = sum(md.wJq .* func.(u_values, equations))40if normalize == true41integral /= sum(md.wJq)42end43return integral44end4546function analyze(::typeof(entropy_timederivative), du, u, t,47mesh::DGMultiMesh, equations, dg::DGMulti, cache)48rd = dg.basis49md = mesh.md50@unpack u_values = cache5152# interpolate u, du to quadrature points53du_values = similar(u_values) # Todo: DGMulti. Can we move this to the analysis cache somehow?54apply_to_each_field(mul_by!(rd.Vq), du_values, du)55apply_to_each_field(mul_by!(rd.Vq), u_values, u)5657# compute ∫v(u) * du/dt = ∫dS/dt. We can directly compute v(u) instead of computing the entropy58# projection here, since the RHS will be projected to polynomials of degree N and testing with59# the L2 projection of v(u) would be equivalent to testing with v(u) due to the moment-preserving60# property of the L2 projection.61dS_dt = zero(eltype(first(du)))62for i in Base.OneTo(length(md.wJq))63dS_dt += dot(cons2entropy(u_values[i], equations), du_values[i]) * md.wJq[i]64end65return dS_dt66end6768# This function is used in `analyze(::Val{:l2_divb},...)` and `analyze(::Val{:linf_divb},...)`69function compute_local_divergence!(local_divergence, element, vector_field,70mesh, dg::DGMulti, cache)71@unpack md = mesh72rd = dg.basis73uEltype = eltype(first(vector_field))7475fill!(local_divergence, zero(uEltype))7677# computes dU_i/dx_i = ∑_j dxhat_j/dx_i * dU_i / dxhat_j78# dU_i/dx_i is then accumulated into local_divergence.79# TODO: DGMulti. Extend to curved elements.80for i in eachdim(mesh)81for j in eachdim(mesh)82geometric_scaling = md.rstxyzJ[i, j][1, element]83jth_ref_derivative_matrix = rd.Drst[j]84mul!(local_divergence, jth_ref_derivative_matrix, vector_field[i],85geometric_scaling, one(uEltype))86end87end88end8990get_component(u::StructArray, i::Int) = StructArrays.component(u, i)91get_component(u::AbstractArray{<:SVector}, i::Int) = getindex.(u, i)9293function analyze(::Val{:l2_divb}, du, u, t,94mesh::DGMultiMesh, equations::IdealGlmMhdEquations2D,95dg::DGMulti, cache)96@unpack md = mesh97rd = dg.basis98B1 = get_component(u, 6)99B2 = get_component(u, 7)100B = (B1, B2)101102uEltype = eltype(B1)103l2norm_divB = zero(uEltype)104local_divB = zeros(uEltype, size(B1, 1))105for e in eachelement(mesh, dg, cache)106compute_local_divergence!(local_divB, e, view.(B, :, e), mesh, dg, cache)107108# TODO: DGMulti. Extend to curved elements.109# compute L2 norm squared via J[1, e] * u' * M * u110local_l2norm_divB = md.J[1, e] * dot(local_divB, rd.M, local_divB)111l2norm_divB += local_l2norm_divB112end113114return sqrt(l2norm_divB)115end116117function analyze(::Val{:linf_divb}, du, u, t,118mesh::DGMultiMesh, equations::IdealGlmMhdEquations2D,119dg::DGMulti, cache)120B1 = get_component(u, 6)121B2 = get_component(u, 7)122B = (B1, B2)123124uEltype = eltype(B1)125linf_divB = zero(uEltype)126local_divB = zeros(uEltype, size(B1, 1))127@batch reduction=(max, linf_divb) for e in eachelement(mesh, dg, cache)128compute_local_divergence!(local_divB, e, view.(B, :, e), mesh, dg, cache)129130# compute maximum norm131linf_divB = max(linf_divB, maximum(abs, local_divB))132end133134return linf_divB135end136137function create_cache_analysis(analyzer, mesh::DGMultiMesh,138equations, dg::DGMulti, cache,139RealT, uEltype)140md = mesh.md141return (;)142end143144SolutionAnalyzer(rd::RefElemData) = rd145146nelements(mesh::DGMultiMesh, ::DGMulti, other_args...) = mesh.md.num_elements147function nelementsglobal(mesh::DGMultiMesh, solver::DGMulti, cache)148if mpi_isparallel()149error("`nelementsglobal` is not implemented for `DGMultiMesh` when used in parallel with MPI")150else151return ndofs(mesh, solver, cache)152end153end154function ndofsglobal(mesh::DGMultiMesh, solver::DGMulti, cache)155if mpi_isparallel()156error("`ndofsglobal` is not implemented for `DGMultiMesh` when used in parallel with MPI")157else158return ndofs(mesh, solver, cache)159end160end161end # @muladd162163164