Path: blob/main/src/visualization/utilities.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@inline num_faces(elem::Tri) = 38@inline num_faces(elem::Quad) = 4910# compute_triangle_area(tri)11#12# Computes the area of a triangle given `tri`, which is a tuple of three points (vectors),13# using the [Shoelace_formula](https://en.wikipedia.org/wiki/Shoelace_formula).14function compute_triangle_area(tri)15A, B, C = tri16return 0.5f0 * (A[1] * (B[2] - C[2]) + B[1] * (C[2] - A[2]) + C[1] * (A[2] - B[2]))17end1819# reference_plotting_triangulation(reference_plotting_coordinates)20#21# Computes a triangulation of the points in `reference_plotting_coordinates`, which is a tuple containing22# vectors of plotting points on the reference element (e.g., reference_plotting_coordinates = (r,s)).23# The reference element is assumed to be [-1,1]^d.24#25# This function returns `t` which is a `3 x N_tri` Matrix{Int} containing indices of triangles in the26# triangulation of the plotting points, with zero-volume triangles removed.27#28# For example, r[t[1, i]] returns the first reference coordinate of the 1st point on the ith triangle.29function reference_plotting_triangulation(reference_plotting_coordinates,30tol = 50 * eps())31# on-the-fly triangulation of plotting nodes on the reference element32tri_in = Triangulate.TriangulateIO()33tri_in.pointlist = permutedims(hcat(reference_plotting_coordinates...))34tri_out, _ = Triangulate.triangulate("Q", tri_in)35triangles = tri_out.trianglelist3637# filter out sliver triangles38has_volume = fill(true, size(triangles, 2))39for i in axes(triangles, 2)40ids = @view triangles[:, i]41x_points = @view tri_out.pointlist[1, ids]42y_points = @view tri_out.pointlist[2, ids]43area = compute_triangle_area(zip(x_points, y_points))44if abs(area) < tol45has_volume[i] = false46end47end48return permutedims(triangles[:, findall(has_volume)])49end5051# This function is used to avoid type instabilities when calling `digest_solution_variables`.52function transform_to_solution_variables!(u, solution_variables, equations)53for (i, u_i) in enumerate(u)54u[i] = solution_variables(u_i, equations)55end56end5758# global_plotting_triangulation_triplot(u_plot, rst_plot, xyz_plot)59#60# Returns (plotting_coordinates_x, plotting_coordinates_y, ..., plotting_values, plotting_triangulation).61# Output can be used with TriplotRecipes.DGTriPseudocolor(...).62#63# Inputs:64# - xyz_plot = plotting points (tuple of matrices of size (Nplot, K))65# - u_plot = matrix of size (Nplot, K) representing solution to plot.66# - t = triangulation of reference plotting points67function global_plotting_triangulation_triplot(xyz_plot, u_plot, t)68@assert size(first(xyz_plot), 1)==size(u_plot, 1) "Row dimension of u_plot does not match row dimension of xyz_plot"6970# build discontinuous data on plotting triangular mesh71num_plotting_points, num_elements = size(u_plot)72num_reference_plotting_triangles = size(t, 1)73num_plotting_elements_total = num_reference_plotting_triangles * num_elements7475# each column of `tp` corresponds to a vertex of a plotting triangle76tp = zeros(Int32, 3, num_plotting_elements_total)77zp = similar(tp, eltype(u_plot))78for e in 1:num_elements79for i in 1:num_reference_plotting_triangles80tp[:, i + (e - 1) * num_reference_plotting_triangles] .= @views t[i, :] .+81(e - 1) *82num_plotting_points83zp[:, i + (e - 1) * num_reference_plotting_triangles] .= @views u_plot[t[i,84:],85e]86end87end88return vec.(xyz_plot)..., zp, tp89end9091function get_face_node_indices(r, s, dg::DGSEM, tol = 100 * eps())92face_1 = findall(@. abs(s + 1) < tol)93face_2 = findall(@. abs(r - 1) < tol)94face_3 = findall(@. abs(s - 1) < tol)95face_4 = findall(@. abs(r + 1) < tol)96Fmask = hcat(face_1, face_2, face_3, face_4)97return Fmask98end99100# dispatch on semi101function mesh_plotting_wireframe(u, semi)102mesh_plotting_wireframe(u, mesh_equations_solver_cache(semi)...)103end104105# mesh_plotting_wireframe(u, mesh, equations, dg::DGMulti, cache; num_plotting_pts=25)106#107# Generates data for plotting a mesh wireframe given StartUpDG data types.108# Returns (plotting_coordinates_x, plotting_coordinates_y, nothing) for a 2D mesh wireframe.109function mesh_plotting_wireframe(u::StructArray, mesh, equations, dg::DGMulti, cache;110nvisnodes = 2 * nnodes(dg))111@unpack md = mesh112rd = dg.basis113114# Construct 1D plotting interpolation matrix `Vp1D` for a single face115@unpack N, Fmask = rd116num_face_points = length(Fmask) ÷ num_faces(rd.element_type)117vandermonde_matrix_1D = StartUpDG.vandermonde(Line(), N,118StartUpDG.nodes(Line(),119num_face_points - 1))120rplot = LinRange(-1, 1, nvisnodes)121Vp1D = StartUpDG.vandermonde(Line(), N, rplot) / vandermonde_matrix_1D122123num_faces_total = num_faces(rd.element_type) * md.num_elements124xf, yf = map(x -> reshape(view(x, Fmask, :), num_face_points, num_faces_total),125md.xyz)126uf = similar(u, size(xf))127apply_to_each_field((out, x) -> out .= reshape(view(x, Fmask, :), num_face_points,128num_faces_total), uf, u)129130num_face_plotting_points = size(Vp1D, 1)131x_mesh, y_mesh = ntuple(_ -> zeros(num_face_plotting_points, num_faces_total), 2)132u_mesh = similar(u, (num_face_plotting_points, num_faces_total))133for f in 1:num_faces_total134mul!(view(x_mesh, :, f), Vp1D, view(xf, :, f))135mul!(view(y_mesh, :, f), Vp1D, view(yf, :, f))136apply_to_each_field(mul_by!(Vp1D), view(u_mesh, :, f), view(uf, :, f))137end138139return x_mesh, y_mesh, u_mesh140end141142function mesh_plotting_wireframe(u::StructArray, mesh, equations, dg::DGSEM, cache;143nvisnodes = 2 * nnodes(dg))144145# build nodes on reference element (seems to be the right ordering)146r, s = reference_node_coordinates_2d(dg)147148# extract node coordinates149uEltype = eltype(first(u))150nvars = nvariables(equations)151n_nodes_2d = nnodes(dg)^ndims(mesh)152n_elements = nelements(dg, cache)153x = reshape(view(cache.elements.node_coordinates, 1, :, :, :), n_nodes_2d,154n_elements)155y = reshape(view(cache.elements.node_coordinates, 2, :, :, :), n_nodes_2d,156n_elements)157158# extract indices of local face nodes for wireframe plotting159Fmask = get_face_node_indices(r, s, dg)160plotting_interp_matrix1D = face_plotting_interpolation_matrix(dg;161nvisnodes = nvisnodes)162163# These 5 lines extract the face values on each element from the arrays x,y,sol_to_plot.164# The resulting arrays are then reshaped so that xf, yf, sol_f are Matrix types of size165# (Number of face plotting nodes) x (Number of faces).166function face_first_reshape(x, num_nodes_1D, num_nodes, num_elements)167num_reference_faces = 2 * ndims(mesh)168xf = view(reshape(x, num_nodes, num_elements), vec(Fmask), :)169return reshape(xf, num_nodes_1D, num_elements * num_reference_faces)170end171function reshape_and_interpolate(x)172plotting_interp_matrix1D *173face_first_reshape(x, nnodes(dg), n_nodes_2d, n_elements)174end175xfp, yfp = map(reshape_and_interpolate, (x, y))176ufp = StructArray{SVector{nvars, uEltype}}(map(reshape_and_interpolate,177StructArrays.components(u)))178179return xfp, yfp, ufp180end181182function mesh_plotting_wireframe(u::ScalarData, mesh, equations, dg::DGSEM, cache;183nvisnodes = 2 * nnodes(dg))184185# build nodes on reference element (seems to be the right ordering)186r, s = reference_node_coordinates_2d(dg)187188# extract node coordinates189n_nodes_2d = nnodes(dg)^ndims(mesh)190n_elements = nelements(dg, cache)191x = reshape(view(cache.elements.node_coordinates, 1, :, :, :), n_nodes_2d,192n_elements)193y = reshape(view(cache.elements.node_coordinates, 2, :, :, :), n_nodes_2d,194n_elements)195196# extract indices of local face nodes for wireframe plotting197Fmask = get_face_node_indices(r, s, dg)198plotting_interp_matrix1D = face_plotting_interpolation_matrix(dg;199nvisnodes = nvisnodes)200201# These 5 lines extract the face values on each element from the arrays x,y,sol_to_plot.202# The resulting arrays are then reshaped so that xf, yf, sol_f are Matrix types of size203# (Number of face plotting nodes) x (Number of faces).204function face_first_reshape(x, num_nodes_1D, num_nodes, num_elements)205num_reference_faces = 2 * ndims(mesh)206xf = view(reshape(x, num_nodes, num_elements), vec(Fmask), :)207return reshape(xf, num_nodes_1D, num_elements * num_reference_faces)208end209function reshape_and_interpolate(x)210plotting_interp_matrix1D *211face_first_reshape(x, nnodes(dg), n_nodes_2d, n_elements)212end213xfp, yfp, ufp = map(reshape_and_interpolate, (x, y, u.data))214215return xfp, yfp, ufp216end217218function mesh_plotting_wireframe(u::ScalarData, mesh, equations, dg::DGMulti, cache;219nvisnodes = 2 * nnodes(dg))220@unpack md = mesh221rd = dg.basis222223# Construct 1D plotting interpolation matrix `Vp1D` for a single face224@unpack N, Fmask = rd225vandermonde_matrix_1D = StartUpDG.vandermonde(Line(), N, StartUpDG.nodes(Line(), N))226rplot = LinRange(-1, 1, nvisnodes)227Vp1D = StartUpDG.vandermonde(Line(), N, rplot) / vandermonde_matrix_1D228229num_face_points = N + 1230num_faces_total = num_faces(rd.element_type) * md.num_elements231xf, yf, uf = map(x -> reshape(view(x, Fmask, :), num_face_points, num_faces_total),232(md.xyz..., u.data))233234num_face_plotting_points = size(Vp1D, 1)235x_mesh, y_mesh = ntuple(_ -> zeros(num_face_plotting_points, num_faces_total), 2)236u_mesh = similar(u.data, (num_face_plotting_points, num_faces_total))237for f in 1:num_faces_total238mul!(view(x_mesh, :, f), Vp1D, view(xf, :, f))239mul!(view(y_mesh, :, f), Vp1D, view(yf, :, f))240mul!(view(u_mesh, :, f), Vp1D, view(uf, :, f))241end242243return x_mesh, y_mesh, u_mesh244end245246# These methods are used internally to set the default value of the solution variables:247# - If a `cons2prim` for the given `equations` exists, use it248# - Otherwise, use `cons2cons`, which is defined for all systems of equations249digest_solution_variables(equations, solution_variables) = solution_variables250function digest_solution_variables(equations, solution_variables::Nothing)251if hasmethod(cons2prim, Tuple{AbstractVector, typeof(equations)})252return cons2prim253else254return cons2cons255end256end257258digest_variable_names(solution_variables_, equations, variable_names) = variable_names259digest_variable_names(solution_variables_, equations, ::Nothing) = SVector(varnames(solution_variables_,260equations))261262"""263adapt_to_mesh_level!(u_ode, semi, level)264adapt_to_mesh_level!(sol::Trixi.TrixiODESolution, level)265266Like [`adapt_to_mesh_level`](@ref), but modifies the solution and parts of the267semidiscretization (mesh and caches) in place.268"""269function adapt_to_mesh_level!(u_ode, semi, level)270# Create AMR callback with controller that refines everything towards a single level271amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),272base_level = level)273amr_callback = AMRCallback(semi, amr_controller, interval = 0)274275# Adapt mesh until it does not change anymore276has_changed = amr_callback.affect!(u_ode, semi, 0.0, 0)277while has_changed278has_changed = amr_callback.affect!(u_ode, semi, 0.0, 0)279end280281return u_ode, semi282end283284function adapt_to_mesh_level!(sol::TrixiODESolution, level)285adapt_to_mesh_level!(sol.u[end], sol.prob.p, level)286end287288"""289adapt_to_mesh_level(u_ode, semi, level)290adapt_to_mesh_level(sol::Trixi.TrixiODESolution, level)291292Use the regular adaptive mesh refinement routines to adaptively refine/coarsen the solution `u_ode`293with semidiscretization `semi` towards a uniformly refined grid with refinement level `level`. The294solution and semidiscretization are copied such that the original objects remain *unaltered*.295296A convenience method accepts an ODE solution object, from which solution and semidiscretization are297extracted as needed.298299See also: [`adapt_to_mesh_level!`](@ref)300"""301function adapt_to_mesh_level(u_ode, semi, level)302# Create new semidiscretization with copy of the current mesh303mesh, _, _, _ = mesh_equations_solver_cache(semi)304new_semi = remake(semi, mesh = deepcopy(mesh))305306return adapt_to_mesh_level!(deepcopy(u_ode), new_semi, level)307end308309function adapt_to_mesh_level(sol::TrixiODESolution, level)310adapt_to_mesh_level(sol.u[end], sol.prob.p, level)311end312313# Extract data from a 2D/3D DG solution and prepare it for visualization as a heatmap/contour plot.314#315# Returns a tuple with316# - x coordinates317# - y coordinates318# - nodal 2D data319# - x vertices for mesh visualization320# - y vertices for mesh visualization321#322# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may323# thus be changed in future releases.324function get_data_2d(center_level_0, length_level_0, leaf_cells, coordinates, levels,325ndims, unstructured_data, n_nodes,326grid_lines = false, max_supported_level = 11, nvisnodes = nothing,327slice = :xy, point = (0.0, 0.0, 0.0))328# Determine resolution for data interpolation329max_level = maximum(levels)330if max_level > max_supported_level331error("Maximum refinement level $max_level is higher than " *332"maximum supported level $max_supported_level")333end334max_available_nodes_per_finest_element = 2^(max_supported_level - max_level)335if nvisnodes === nothing336max_nvisnodes = 2 * n_nodes337elseif nvisnodes == 0338max_nvisnodes = n_nodes339else340max_nvisnodes = nvisnodes341end342nvisnodes_at_max_level = min(max_available_nodes_per_finest_element, max_nvisnodes)343resolution = nvisnodes_at_max_level * 2^max_level344nvisnodes_per_level = [2^(max_level - level) * nvisnodes_at_max_level345for level in 0:max_level]346# nvisnodes_per_level is an array (accessed by "level + 1" to accommodate347# level-0-cell) that contains the number of visualization nodes for any348# refinement level to visualize on an equidistant grid349350if ndims == 3351(unstructured_data, coordinates, levels,352center_level_0) = unstructured_3d_to_2d(unstructured_data,353coordinates, levels, length_level_0,354center_level_0, slice,355point)356end357358# Normalize element coordinates: move center to (0, 0) and domain size to [-1, 1]²359n_elements = length(levels)360normalized_coordinates = similar(coordinates)361for element_id in 1:n_elements362@views normalized_coordinates[:, element_id] .= ((coordinates[:, element_id] .-363center_level_0) ./364(length_level_0 / 2))365end366367# Interpolate unstructured DG data to structured data368(structured_data = unstructured2structured(unstructured_data,369normalized_coordinates,370levels, resolution, nvisnodes_per_level))371372# Interpolate cell-centered values to node-centered values373node_centered_data = cell2node(structured_data)374375# Determine axis coordinates for contour plot376xs = collect(range(-1, 1, length = resolution + 1)) .* length_level_0 / 2 .+377center_level_0[1]378ys = collect(range(-1, 1, length = resolution + 1)) .* length_level_0 / 2 .+379center_level_0[2]380381# Determine element vertices to plot grid lines382if grid_lines383mesh_vertices_x, mesh_vertices_y = calc_vertices(coordinates, levels,384length_level_0)385else386mesh_vertices_x = Vector{Float64}(undef, 0)387mesh_vertices_y = Vector{Float64}(undef, 0)388end389390return xs, ys, node_centered_data, mesh_vertices_x, mesh_vertices_y391end392393# For finite difference methods (FDSBP), we do not want to reinterpolate the data, but use the same394# nodes as input nodes. For DG methods, we usually want to reinterpolate the data on an equidistant grid.395default_reinterpolate(solver) = solver isa FDSBP ? false : true396397# Extract data from a 1D DG solution and prepare it for visualization as a line plot.398# This returns a tuple with399# - x coordinates400# - nodal 1D data401#402# Note: This is a low-level function that is not considered as part of Trixi's interface and may403# thus be changed in future releases.404function get_data_1d(original_nodes, unstructured_data, nvisnodes, reinterpolate)405# Get the dimensions of u; where n_vars is the number of variables, n_nodes the number of nodal values per element and n_elements the total number of elements.406n_nodes, n_elements, n_vars = size(unstructured_data)407408# If `reinterpolate` is `false`, we do nothing.409# If `reinterpolate` is `true`, the output nodes are equidistantly spaced.410if reinterpolate == false411interpolated_nodes = original_nodes412interpolated_data = unstructured_data413else414# Set the amount of nodes visualized according to nvisnodes.415if nvisnodes === nothing416max_nvisnodes = 2 * n_nodes417elseif nvisnodes == 0418max_nvisnodes = n_nodes419else420@assert nvisnodes>=2 "nvisnodes must be zero or >= 2"421max_nvisnodes = nvisnodes422end423424interpolated_nodes = Array{eltype(original_nodes), 2}(undef, max_nvisnodes,425n_elements)426interpolated_data = Array{eltype(unstructured_data), 3}(undef, max_nvisnodes,427n_elements, n_vars)428429for j in 1:n_elements430# Interpolate on an equidistant grid.431interpolated_nodes[:, j] .= range(original_nodes[1, 1, j],432original_nodes[1, end, j],433length = max_nvisnodes)434end435436nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)437nodes_out = collect(range(-1, 1, length = max_nvisnodes))438439# Calculate vandermonde matrix for interpolation.440vandermonde = polynomial_interpolation_matrix(nodes_in, nodes_out)441442# Iterate over all variables.443for v in 1:n_vars444# Interpolate data for each element.445for element in 1:n_elements446multiply_scalar_dimensionwise!(@view(interpolated_data[:, element, v]),447vandermonde,448@view(unstructured_data[:, element, v]))449end450end451end452453# Return results after data is reshaped454return vec(interpolated_nodes), reshape(interpolated_data, :, n_vars),455vcat(original_nodes[1, 1, :], original_nodes[1, end, end])456end457458# Change order of dimensions (variables are now last) and convert data to `solution_variables`459#460# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may461# thus be changed in future releases.462function get_unstructured_data(u, solution_variables, mesh, equations, solver, cache)463if solution_variables === cons2cons464raw_data = u465n_vars = size(raw_data, 1)466else467# Similar to `save_solution_file` in `callbacks_step/save_solution_dg.jl`.468# However, we cannot use `reinterpret` here as `u` might have a non-bits type.469# See https://github.com/trixi-framework/Trixi.jl/pull/2388470n_vars_in = nvariables(equations)471n_vars = length(solution_variables(get_node_vars(u, equations, solver),472equations))473raw_data = Array{eltype(u)}(undef, n_vars, Base.tail(size(u))...)474reshaped_u = reshape(u, n_vars_in, :)475reshaped_r = reshape(raw_data, n_vars, :)476for idx in axes(reshaped_u, 2)477reshaped_r[:, idx] = solution_variables(get_node_vars(reshaped_u, equations,478solver, idx),479equations)480end481end482483unstructured_data = Array{eltype(raw_data)}(undef,484ntuple((d) -> nnodes(solver),485ndims(equations))...,486nelements(solver, cache), n_vars)487for variable in 1:n_vars488@views unstructured_data[.., :, variable] .= raw_data[variable, .., :]489end490491return unstructured_data492end493494# This method is only for plotting 1D functions495function get_unstructured_data(func::Function, solution_variables,496mesh::AbstractMesh{1}, equations, solver, cache)497original_nodes = cache.elements.node_coordinates498# original_nodes has size (1, nnodes, nelements)499# we want u to have size (nvars, nnodes, nelements)500# func.(original_nodes, equations) has size (1, nnodes, nelements), where each component has length n_vars501# Therefore, we drop the first (singleton) dimension and then stack the components502u = stack(func.(SVector.(dropdims(original_nodes; dims = 1)), equations))503return get_unstructured_data(u, solution_variables, mesh, equations, solver, cache)504end505506# Convert cell-centered values to node-centered values by averaging over all507# four neighbors. Solution values at the edges are padded with ghost values508# computed via linear extrapolation.509#510# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may511# thus be changed in future releases.512function cell2node(cell_centered_data)513# Create temporary data structure to make the averaging algorithm as simple514# as possible (by using a ghost layer)515tmp = similar(first(cell_centered_data), size(first(cell_centered_data)) .+ (2, 2))516517# Create output data structure518resolution_in, _ = size(first(cell_centered_data))519resolution_out = resolution_in + 1520node_centered_data = [Matrix{Float64}(undef, resolution_out, resolution_out)521for _ in eachindex(cell_centered_data)]522523for (cell_data, node_data) in zip(cell_centered_data, node_centered_data)524# Fill center with original data525tmp[2:(end - 1), 2:(end - 1)] .= cell_data526527# Linear extrapolation of top and bottom rows528tmp[1, 2:(end - 1)] .= cell_data[1, :] .+ (cell_data[1, :] .- cell_data[2, :])529tmp[end, 2:(end - 1)] .= (cell_data[end, :] .+530(cell_data[end, :] .- cell_data[end - 1, :]))531532# Linear extrapolatation of left and right columns533tmp[2:(end - 1), 1] .= cell_data[:, 1] .+ (cell_data[:, 1] .- cell_data[:, 2])534tmp[2:(end - 1), end] .= (cell_data[:, end] .+535(cell_data[:, end] .- cell_data[:, end - 1]))536537# Corners perform the linear extrapolatation along diagonals538tmp[1, 1] = tmp[2, 2] + (tmp[2, 2] - tmp[3, 3])539tmp[1, end] = tmp[2, end - 1] + (tmp[2, end - 1] - tmp[3, end - 2])540tmp[end, 1] = tmp[end - 1, 2] + (tmp[end - 1, 2] - tmp[end - 2, 3])541tmp[end, end] = (tmp[end - 1, end - 1] +542(tmp[end - 1, end - 1] - tmp[end - 2, end - 2]))543544# Obtain node-centered value by averaging over neighboring cell-centered values545for j in 1:resolution_out546for i in 1:resolution_out547node_data[i, j] = (tmp[i, j] +548tmp[i + 1, j] +549tmp[i, j + 1] +550tmp[i + 1, j + 1]) / 4551end552end553end554555# Transpose556for (index, data) in enumerate(node_centered_data)557node_centered_data[index] = permutedims(data)558end559560return node_centered_data561end562563# Convert 3d unstructured data to 2d data.564# Additional to the new unstructured data updated coordinates, levels and565# center coordinates are returned.566#567# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may568# thus be changed in future releases.569function unstructured_3d_to_2d(unstructured_data, coordinates, levels,570length_level_0, center_level_0, slice,571point)572if slice === :yz573slice_dimension = 1574other_dimensions = [2, 3]575elseif slice === :xz576slice_dimension = 2577other_dimensions = [1, 3]578elseif slice === :xy579slice_dimension = 3580other_dimensions = [1, 2]581else582error("illegal dimension '$slice', supported dimensions are :yz, :xz, and :xy")583end584585# Limits of domain in slice dimension586lower_limit = center_level_0[slice_dimension] - length_level_0 / 2587upper_limit = center_level_0[slice_dimension] + length_level_0 / 2588589@assert length(point)>=3 "Point must be three-dimensional."590if point[slice_dimension] < lower_limit || point[slice_dimension] > upper_limit591error(string("Slice plane is outside of domain.",592" point[$slice_dimension]=$(point[slice_dimension]) must be between $lower_limit and $upper_limit"))593end594595# Extract data shape information596n_nodes_in, _, _, n_elements, n_variables = size(unstructured_data)597598# Get node coordinates for DG locations on reference element599nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)600601# New unstructured data has one dimension less.602# The redundant element ids are removed later.603@views new_unstructured_data = similar(unstructured_data[1, ..])604605# Declare new empty arrays to fill in new coordinates and levels606new_coordinates = Array{Float64}(undef, 2, n_elements)607new_levels = Array{eltype(levels)}(undef, n_elements)608609# Counter for new element ids610new_id = 0611612# Save vandermonde matrices in a Dict to prevent redundant generation613vandermonde_to_2d = Dict()614615# Permute dimensions such that the slice dimension is always the616# third dimension of the array. Below we can always interpolate in the617# third dimension.618if slice === :yz619unstructured_data = permutedims(unstructured_data, [2, 3, 1, 4, 5])620elseif slice === :xz621unstructured_data = permutedims(unstructured_data, [1, 3, 2, 4, 5])622end623624for element_id in 1:n_elements625# Distance from center to border of this element (half the length)626element_length = length_level_0 / 2^levels[element_id]627min_coordinate = coordinates[:, element_id] .- element_length / 2628max_coordinate = coordinates[:, element_id] .+ element_length / 2629630# Check if slice plane and current element intersect.631# The first check uses a "greater but not equal" to only match one cell if the632# slice plane lies between two cells.633# The second check is needed if the slice plane is at the upper border of634# the domain due to this.635if !((min_coordinate[slice_dimension] <= point[slice_dimension] &&636max_coordinate[slice_dimension] > point[slice_dimension]) ||637(point[slice_dimension] == upper_limit &&638max_coordinate[slice_dimension] == upper_limit))639# Continue for loop if they don't intersect640continue641end642643# This element is of interest644new_id += 1645646# Add element to new coordinates and levels647new_coordinates[:, new_id] = coordinates[other_dimensions, element_id]648new_levels[new_id] = levels[element_id]649650# Construct vandermonde matrix (or load from Dict if possible)651normalized_intercept = (point[slice_dimension] -652min_coordinate[slice_dimension]) /653element_length * 2 - 1654655if haskey(vandermonde_to_2d, normalized_intercept)656vandermonde = vandermonde_to_2d[normalized_intercept]657else658# Generate vandermonde matrix to interpolate values at nodes_in to one value659vandermonde = polynomial_interpolation_matrix(nodes_in,660[normalized_intercept])661vandermonde_to_2d[normalized_intercept] = vandermonde662end663664# 1D interpolation to specified slice plane665# We permuted the dimensions above such that now the dimension in which666# we will interpolate is always the third one.667for i in 1:n_nodes_in668for ii in 1:n_nodes_in669# Interpolate in the third dimension670data = unstructured_data[i, ii, :, element_id, :]671672value = multiply_dimensionwise(vandermonde, permutedims(data))673new_unstructured_data[i, ii, new_id, :] = value[:, 1]674end675end676end677678# Remove redundant element ids679unstructured_data = new_unstructured_data[:, :, 1:new_id, :]680new_coordinates = new_coordinates[:, 1:new_id]681new_levels = new_levels[1:new_id]682683center_level_0 = center_level_0[other_dimensions]684685return unstructured_data, new_coordinates, new_levels, center_level_0686end687688# Convert 2d unstructured data to 1d slice and interpolate them.689function unstructured_2d_to_1d(original_nodes, unstructured_data, nvisnodes,690reinterpolate, slice, point)691if slice === :x692slice_dimension = 2693other_dimension = 1694elseif slice === :y695slice_dimension = 1696other_dimension = 2697else698error("illegal dimension '$slice', supported dimensions are :x and :y")699end700701# Set up data structures to store new 1D data.702@views new_unstructured_data = similar(unstructured_data[1, ..])703@views new_nodes = similar(original_nodes[1, 1, ..])704705n_nodes_in, _, n_elements, n_variables = size(unstructured_data)706nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)707708# Test if point lies in the domain.709lower_limit = original_nodes[1, 1, 1, 1]710upper_limit = original_nodes[1, n_nodes_in, n_nodes_in, n_elements]711712@assert length(point)>=2 "Point must be two-dimensional."713if point[slice_dimension] < lower_limit || point[slice_dimension] > upper_limit714error(string("Slice axis is outside of domain. ",715" point[$slice_dimension]=$(point[slice_dimension]) must be between $lower_limit and $upper_limit"))716end717718# Count the amount of new elements.719new_id = 0720721# Permute dimensions so that the slice dimension is always in the correct place for later use.722if slice === :y723original_nodes = permutedims(original_nodes, [1, 3, 2, 4])724unstructured_data = permutedims(unstructured_data, [2, 1, 3, 4])725end726727# Iterate over all elements to find the ones that lie on the slice axis.728for element_id in 1:n_elements729min_coordinate = original_nodes[:, 1, 1, element_id]730max_coordinate = original_nodes[:, n_nodes_in, n_nodes_in, element_id]731element_length = max_coordinate - min_coordinate732733# Test if the element is on the slice axis. If not just continue with the next element.734if !((min_coordinate[slice_dimension] <= point[slice_dimension] &&735max_coordinate[slice_dimension] > point[slice_dimension]) ||736(point[slice_dimension] == upper_limit &&737max_coordinate[slice_dimension] == upper_limit))738continue739end740741new_id += 1742743# Construct vandermonde matrix for interpolation of each 2D element to a 1D element.744normalized_intercept = (point[slice_dimension] -745min_coordinate[slice_dimension]) /746element_length[1] * 2 - 1747vandermonde = polynomial_interpolation_matrix(nodes_in, normalized_intercept)748749# Interpolate to each node of new 1D element.750for v in 1:n_variables751for node in 1:n_nodes_in752new_unstructured_data[node, new_id, v] = (vandermonde * unstructured_data[node,753:,754element_id,755v])[1]756end757end758759new_nodes[:, new_id] = original_nodes[other_dimension, :, 1, element_id]760end761762return get_data_1d(reshape(new_nodes[:, 1:new_id], 1, n_nodes_in, new_id),763new_unstructured_data[:, 1:new_id, :], nvisnodes, reinterpolate)764end765766# Calculate the arc length of a curve given by ndims x npoints point coordinates (piece-wise linear approximation)767function calc_arc_length(coordinates)768n_points = size(coordinates, 2)769arc_length = zeros(n_points)770for i in 1:(n_points - 1)771dist_squared = zero(eltype(arc_length))772for j in axes(coordinates, 1)773dist_squared += (coordinates[j, i + 1] - coordinates[j, i])^2774end775arc_length[i + 1] = arc_length[i] + sqrt(dist_squared)776end777return arc_length778end779780# Convert 2d unstructured data to 1d data at given curve for the TreeMesh.781function unstructured_2d_to_1d_curve(original_nodes, unstructured_data, nvisnodes,782curve, mesh::TreeMesh, solver, cache)783n_points_curve = size(curve)[2]784n_nodes, _, n_elements, n_variables = size(unstructured_data)785nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)786baryweights_in = barycentric_weights(nodes_in)787788# Utility function to extract points as `SVector`s789get_point(data, idx...) = SVector(data[1, idx...], data[2, idx...])790791# Check if input is correct.792min = get_point(original_nodes, 1, 1, 1)793max = get_point(original_nodes, n_nodes, n_nodes, n_elements)794@assert size(curve)==(2, size(curve)[2]) "Coordinates along curve must be 2xn dimensional."795for element in 1:n_points_curve796p = get_point(curve, element)797if any(p .< min) || any(p .> max)798throw(DomainError("Some coordinates from `curve` are outside of the domain."))799end800end801802# Set nodes according to the length of the curve.803arc_length = calc_arc_length(curve)804805# Setup data structures.806data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)807temp_data = Array{Float64}(undef, n_nodes, n_points_curve, n_variables)808809# For each coordinate find the corresponding element with its id.810element_ids = get_elements_by_coordinates(curve, mesh, solver, cache)811812# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,813# row vectors. We allocate memory here to improve performance.814vandermonde_x = polynomial_interpolation_matrix(nodes_in, zero(eltype(curve)))815vandermonde_y = similar(vandermonde_x)816817# Iterate over all found elements.818for element in 1:n_points_curve819element_id = element_ids[element]820min_coordinate = get_point(original_nodes, 1, 1,821element_id)822max_coordinate = get_point(original_nodes, n_nodes, n_nodes,823element_id)824element_length = max_coordinate - min_coordinate825826normalized_coordinates = (get_point(curve, element) - min_coordinate) /827element_length[1] * 2 .- 1828829# Interpolate to a single point in each element.830# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,831# row vectors.832polynomial_interpolation_matrix!(vandermonde_x, nodes_in,833normalized_coordinates[1], baryweights_in)834polynomial_interpolation_matrix!(vandermonde_y, nodes_in,835normalized_coordinates[2], baryweights_in)836for v in 1:n_variables837for i in 1:n_nodes838res_i = zero(eltype(temp_data))839for n in 1:n_nodes840res_i += vandermonde_y[n] * unstructured_data[i, n, element_id, v]841end842temp_data[i, element, v] = res_i843end844res_v = zero(eltype(data_on_curve))845for n in 1:n_nodes846res_v += vandermonde_x[n] * temp_data[n, element, v]847end848data_on_curve[element, v] = res_v849end850end851852return arc_length, data_on_curve, nothing853end854855# Convert 2d unstructured data from a general mesh to 1d data at given curve.856#857# We need to loop through all the points and check in which element they are858# located. A general implementation working for all mesh types has to perform859# a naive loop through all nodes. Thus, we use this entry point for dispatching860# on the `mesh` type.861function unstructured_2d_to_1d_curve(u, mesh, equations,862solver, cache,863curve, slice,864point, nvisnodes,865solution_variables)866return unstructured_2d_to_1d_curve_general(u, mesh, equations,867solver, cache,868curve, slice,869point, nvisnodes,870solution_variables)871end872873function unstructured_2d_to_1d_curve_general(u, mesh, equations,874solver, cache,875input_curve, slice,876point, nvisnodes,877solution_variables)878# Create a 'PlotData2DTriangulated' object so a triangulation879# can be used when extracting relevant data.880pd = PlotData2DTriangulated(u, mesh, equations, solver, cache;881solution_variables, nvisnodes)882883# If no curve is defined, create an axis curve.884if input_curve === nothing885input_curve = axis_curve(pd.x, pd.y, nothing, slice, point, nvisnodes)886end887888@assert size(input_curve, 1)==2 "Input 'curve' must be 2xn dimensional."889890# For each coordinate find the corresponding triangle with its ids.891# The default value if no element is found is 0.892ids_by_coordinates = get_ids_by_coordinates(input_curve, pd)893found_coordinates = view(ids_by_coordinates, :, 1) .!= 0894895@assert !iszero(found_coordinates) "No points of 'curve' are inside of the solutions domain."896897# These hold the ids of the elements and triangles the points of the curve sit in.898element_ids = @view ids_by_coordinates[found_coordinates, 1]899triangle_ids = @view ids_by_coordinates[found_coordinates, 2]900901# Shorten the curve, so that it contains only point that were found.902curve = @view input_curve[:, found_coordinates]903904n_variables = length(pd.data[1, 1])905n_points_curve = size(curve, 2)906907# Set nodes according to the length of the curve.908arc_length = calc_arc_length(curve)909910# Setup data structures.911data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)912913# Iterate over all points on the curve.914for point in 1:n_points_curve915point_on_curve = SVector(curve[1, point], curve[2, point])916917element = element_ids[point]918triangle_id = triangle_ids[point]919triangle = (pd.t[triangle_id, 1], pd.t[triangle_id, 2], pd.t[triangle_id, 3])920921# Get the x and y coordinates of the corners of given triangle.922x_coordinates_triangle = SVector(pd.x[triangle[1], element],923pd.x[triangle[2], element],924pd.x[triangle[3], element])925y_coordinates_triangle = SVector(pd.y[triangle[1], element],926pd.y[triangle[2], element],927pd.y[triangle[3], element])928929# Extract solution values in corners of the triangle.930data_in_triangle = (pd.data[triangle[1], element],931pd.data[triangle[2], element],932pd.data[triangle[3], element])933934for v in 1:n_variables935# Extract solution values of variable `v` in corners of the triangle.936values_triangle = SVector(data_in_triangle[1][v],937data_in_triangle[2][v],938data_in_triangle[3][v])939940# Linear interpolation in each triangle to the points on the curve.941data_on_curve[point, v] = triangle_interpolation(x_coordinates_triangle,942y_coordinates_triangle,943values_triangle,944point_on_curve)945end946end947948return arc_length, data_on_curve, nothing949end950951# Convert 3d unstructured data to 1d data at given curve.952function unstructured_3d_to_1d_curve(original_nodes, unstructured_data, nvisnodes,953curve, mesh::TreeMesh, solver, cache)954n_points_curve = size(curve)[2]955n_nodes, _, _, n_elements, n_variables = size(unstructured_data)956nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)957baryweights_in = barycentric_weights(nodes_in)958959# Utility function to extract points as `SVector`s960get_point(data, idx...) = SVector(data[1, idx...],961data[2, idx...],962data[3, idx...])963964# Check if input is correct.965min = get_point(original_nodes, 1, 1, 1, 1)966max = get_point(original_nodes, n_nodes, n_nodes, n_nodes, n_elements)967@assert size(curve)==(3, n_points_curve) "Coordinates along curve must be 3xn dimensional."968for element in 1:n_points_curve969p = get_point(curve, element)970if any(p .< min) || any(p .> max)971throw(DomainError("Some coordinates from `curve` are outside of the domain."))972end973end974975# Set nodes according to the length of the curve.976arc_length = calc_arc_length(curve)977978# Setup data structures.979data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)980temp_data = Array{Float64}(undef, n_nodes, n_nodes + 1, n_points_curve, n_variables)981982# For each coordinate find the corresponding element with its id.983element_ids = get_elements_by_coordinates(curve, mesh, solver, cache)984985# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,986# row vectors. We allocate memory here to improve performance.987vandermonde_x = polynomial_interpolation_matrix(nodes_in, zero(eltype(curve)))988vandermonde_y = similar(vandermonde_x)989vandermonde_z = similar(vandermonde_x)990991# Iterate over all found elements.992for element in 1:n_points_curve993element_id = element_ids[element]994min_coordinate = get_point(original_nodes, 1, 1, 1,995element_id)996max_coordinate = get_point(original_nodes, n_nodes, n_nodes, n_nodes,997element_id)998element_length = max_coordinate - min_coordinate9991000normalized_coordinates = (get_point(curve, element) - min_coordinate) /1001element_length[1] * 2 .- 110021003# Interpolate to a single point in each element.1004# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,1005# row vectors.1006polynomial_interpolation_matrix!(vandermonde_x, nodes_in,1007normalized_coordinates[1], baryweights_in)1008polynomial_interpolation_matrix!(vandermonde_y, nodes_in,1009normalized_coordinates[2], baryweights_in)1010polynomial_interpolation_matrix!(vandermonde_z, nodes_in,1011normalized_coordinates[3], baryweights_in)1012for v in 1:n_variables1013for i in 1:n_nodes1014for ii in 1:n_nodes1015res_ii = zero(eltype(temp_data))1016for n in 1:n_nodes1017res_ii += vandermonde_z[n] *1018unstructured_data[i, ii, n, element_id, v]1019end1020temp_data[i, ii, element, v] = res_ii1021end1022res_i = zero(eltype(temp_data))1023for n in 1:n_nodes1024res_i += vandermonde_y[n] * temp_data[i, n, element, v]1025end1026temp_data[i, n_nodes + 1, element, v] = res_i1027end1028res_v = zero(eltype(temp_data))1029for n in 1:n_nodes1030res_v += vandermonde_x[n] * temp_data[n, n_nodes + 1, element, v]1031end1032data_on_curve[element, v] = res_v1033end1034end10351036return arc_length, data_on_curve, nothing1037end10381039# Convert 3d unstructured data from a general mesh to 1d data at given curve.1040#1041# We need to loop through all the points and check in which element they are1042# located. A general implementation working for all mesh types has to perform1043# a naive loop through all nodes. Thus, we use this entry point for dispatching1044# on the `mesh` type.1045function unstructured_3d_to_1d_curve(u, mesh, equations, solver, cache,1046curve, solution_variables)1047return unstructured_3d_to_1d_curve_general(u, equations, solver, cache,1048curve, solution_variables)1049end10501051function unstructured_3d_to_1d_curve_general(u, equations, solver, cache,1052curve, solution_variables)1053# Set up data structure.1054@assert size(curve, 1) == 31055n_points_curve = size(curve, 2)10561057# Get the number of variables after applying the transformation to solution variables.1058u_node = get_node_vars(u, equations, solver, 1, 1, 1, 1)1059var_node = solution_variables(u_node, equations)1060n_variables = length(var_node)1061data_on_curve = Array{eltype(var_node)}(undef, n_points_curve, n_variables)10621063nodes = cache.elements.node_coordinates1064interpolation_cache = get_value_at_point_3d_cache(u)10651066# Iterate over every point on the curve and determine the solutions value at given point.1067for i in 1:n_points_curve1068point = SVector(curve[1, i], curve[2, i], curve[3, i])1069get_value_at_point_3d!(view(data_on_curve, i, :), point, solution_variables,1070nodes, u, equations, solver;1071cache = interpolation_cache)1072end10731074mesh_vertices_x = nothing10751076return calc_arc_length(curve), data_on_curve, mesh_vertices_x1077end10781079# Check if the first 'amount'-many points can still form a valid tetrahedron.1080function is_valid_tetrahedron(amount, coordinates; tol = 10^-4)1081a = SVector(coordinates[1, 1], coordinates[2, 1], coordinates[3, 1])1082b = SVector(coordinates[1, 2], coordinates[2, 2], coordinates[3, 2])1083c = SVector(coordinates[1, 3], coordinates[2, 3], coordinates[3, 3])1084d = SVector(coordinates[1, 4], coordinates[2, 4], coordinates[3, 4])10851086if amount == 2 # If two points are the same, then no tetrahedron can be formed.1087return !(isapprox(a, b; atol = tol))1088elseif amount == 3 # Check if three points are on the same line.1089return !on_the_same_line(a, b, c; tol = tol)1090elseif amount == 4 # Check if four points form a tetrahedron.1091# This is the same as1092# A = hcat(coordinates[1, :], coordinates[2, :], coordinates[3, :],1093# SVector(1, 1, 1, 1))1094# but more efficient.1095A = SMatrix{4, 4}(coordinates[1, 1], coordinates[2, 1], coordinates[3, 1], 1,1096coordinates[1, 2], coordinates[2, 2], coordinates[3, 2], 1,1097coordinates[1, 3], coordinates[2, 3], coordinates[3, 3], 1,1098coordinates[1, 4], coordinates[2, 4], coordinates[3, 4], 1)1099return !isapprox(det(A), 0; atol = tol)1100else # With one point a tetrahedron can always be formed.1101return true1102end1103end11041105# Check if three given 3D-points are on the same line.1106function on_the_same_line(a, b, c; tol = 10^-4)1107# Calculate the intersection of the a-b-axis at x=0.1108if b[1] == 01109intersect_a_b = b1110else1111intersect_a_b = a - b .* (a[1] / b[1])1112end1113# Calculate the intersection of the a-c-axis at x=0.1114if c[1] == 01115intersect_a_c = c1116else1117intersect_a_c = a - c .* (a[1] / c[1])1118end1119return isapprox(intersect_a_b, intersect_a_c; atol = tol)1120end11211122# Interpolate from four corners of a tetrahedron to a single point.1123function tetrahedron_interpolation(x_coordinates_in, y_coordinates_in, z_coordinates_in,1124values_in, coordinate_out)1125A = hcat(x_coordinates_in, y_coordinates_in, z_coordinates_in, SVector(1, 1, 1, 1))1126c = A \ values_in1127return c[1] * coordinate_out[1] + c[2] * coordinate_out[2] +1128c[3] * coordinate_out[3] + c[4]1129end11301131# Calculate the distances from every entry in `nodes` to the given `point` and return1132# the minimal/maximal squared distance as well as the index of the minimal squared distance.1133function squared_distances_from_single_point!(distances, nodes, point)1134_, n_nodes, _, _, n_elements = size(nodes)1135@assert size(nodes, 1) == 31136@assert size(nodes, 3) == size(nodes, 4) == n_nodes1137@assert size(distances, 1) == size(distances, 2) == size(distances, 3) == n_nodes1138@assert size(distances, 4) == n_elements11391140# Prepare for reduction1141dist_max = typemin(eltype(distances))1142dist_min = typemax(eltype(distances))1143index_min = CartesianIndex(1, 1, 1, 1)11441145# Iterate over every entry1146for element in 1:n_elements1147for z in 1:n_nodes, y in 1:n_nodes, x in 1:n_nodes1148dist = (nodes[1, x, y, z, element] - point[1])^2 +1149(nodes[2, x, y, z, element] - point[2])^2 +1150(nodes[3, x, y, z, element] - point[3])^21151distances[x, y, z, element] = dist11521153dist_max = max(dist_max, dist)1154if dist < dist_min1155dist_min = dist1156index_min = CartesianIndex(x, y, z, element)1157end1158end1159end11601161return dist_max, dist_min, index_min1162end11631164# Interpolate the data on given nodes to a single value at given point.1165function get_value_at_point_3d_cache(u)1166n_variables, n_x_nodes, n_y_nodes, n_z_nodes, n_elements = size(u)1167@assert n_x_nodes == n_y_nodes == n_z_nodes1168n_nodes = n_x_nodes11691170distances = zeros(n_nodes, n_nodes, n_nodes, n_elements)1171coordinates_tetrahedron = Array{Float64, 2}(undef, 3, 4)1172value_tetrahedron = Array{eltype(u)}(undef, n_variables, 4)11731174cache = (; distances, coordinates_tetrahedron, value_tetrahedron)1175return cache1176end11771178function get_value_at_point_3d!(data_on_curve_at_point, point, solution_variables,1179nodes, u, equations, solver;1180cache = get_value_at_point_3d_cache(u))1181# Set up data structures.1182n_variables = size(u, 1)1183(; distances, coordinates_tetrahedron, value_tetrahedron) = cache11841185maximum_distance, _, index = squared_distances_from_single_point!(distances, nodes,1186point)1187# We could also use the following code to find the maximum distance and index:1188# maximum_distance = maximum(distances)1189# index = argmin(distances)1190# However, it is more efficient if we do it in one pass through the memory.11911192# If the point sits exactly on a node, no interpolation is needed.1193nodes_at_index = SVector(nodes[1, index[1], index[2], index[3], index[4]],1194nodes[2, index[1], index[2], index[3], index[4]],1195nodes[3, index[1], index[2], index[3], index[4]])1196if nodes_at_index == point1197u_node = get_node_vars(u, equations, solver,1198index[1], index[2], index[3], index[4])1199data_on_curve_at_point .= solution_variables(u_node, equations)1200return data_on_curve_at_point1201end12021203# Restrict the interpolation to the closest element only.1204closest_element = index[4]1205@views element_distances = distances[:, :, :, closest_element]12061207# Find a tetrahedron, which is given by four corners, to interpolate from.1208for i in 1:41209# Iterate until a valid tetrahedron is found.1210while true1211index = argmin(element_distances)1212element_distances[index[1], index[2], index[3]] = maximum_distance12131214for k in 1:31215coordinates_tetrahedron[k, i] = nodes[k,1216index[1], index[2], index[3],1217closest_element]1218end1219for v in 1:n_variables1220value_tetrahedron[v, i] = u[v,1221index[1], index[2], index[3],1222closest_element]1223end12241225# Look for another point if current tetrahedron is not valid.1226if is_valid_tetrahedron(i, coordinates_tetrahedron)1227break1228end1229end1230end12311232# Interpolate from tetrahedron to given point.1233x_coordinates = SVector(coordinates_tetrahedron[1, 1],1234coordinates_tetrahedron[1, 2],1235coordinates_tetrahedron[1, 3],1236coordinates_tetrahedron[1, 4])1237y_coordinates = SVector(coordinates_tetrahedron[2, 1],1238coordinates_tetrahedron[2, 2],1239coordinates_tetrahedron[2, 3],1240coordinates_tetrahedron[2, 4])1241z_coordinates = SVector(coordinates_tetrahedron[3, 1],1242coordinates_tetrahedron[3, 2],1243coordinates_tetrahedron[3, 3],1244coordinates_tetrahedron[3, 4])1245# We compute the solution_variables first and interpolate them.1246u1 = get_node_vars(value_tetrahedron, equations, solver, 1)1247u2 = get_node_vars(value_tetrahedron, equations, solver, 2)1248u3 = get_node_vars(value_tetrahedron, equations, solver, 3)1249u4 = get_node_vars(value_tetrahedron, equations, solver, 4)1250val1 = solution_variables(u1, equations)1251val2 = solution_variables(u2, equations)1252val3 = solution_variables(u3, equations)1253val4 = solution_variables(u4, equations)1254for v in eachindex(val1)1255values = SVector(val1[v],1256val2[v],1257val3[v],1258val4[v])1259data_on_curve_at_point[v] = tetrahedron_interpolation(x_coordinates,1260y_coordinates,1261z_coordinates,1262values, point)1263end12641265return data_on_curve_at_point1266end12671268# Convert 3d unstructured data to 1d slice and interpolate them.1269function unstructured_3d_to_1d(original_nodes, unstructured_data, nvisnodes,1270reinterpolate, slice, point)1271if slice === :x1272slice_dimension = 11273other_dimensions = [2, 3]1274elseif slice === :y1275slice_dimension = 21276other_dimensions = [1, 3]1277elseif slice === :z1278slice_dimension = 31279other_dimensions = [1, 2]1280else1281error("illegal dimension '$slice', supported dimensions are :x, :y and :z")1282end12831284# Set up data structures to store new 1D data.1285@views new_unstructured_data = similar(unstructured_data[1, 1, ..])1286@views temp_unstructured_data = similar(unstructured_data[1, ..])1287@views new_nodes = similar(original_nodes[1, 1, 1, ..])12881289n_nodes_in, _, _, n_elements, n_variables = size(unstructured_data)1290nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)12911292# Test if point lies in the domain.1293lower_limit = original_nodes[1, 1, 1, 1, 1]1294upper_limit = original_nodes[1, n_nodes_in, n_nodes_in, n_nodes_in, n_elements]12951296@assert length(point)>=3 "Point must be three-dimensional."1297if prod(point[other_dimensions] .< lower_limit) ||1298prod(point[other_dimensions] .> upper_limit)1299error(string("Slice axis is outside of domain. ",1300" point[$other_dimensions]=$(point[other_dimensions]) must be between $lower_limit and $upper_limit"))1301end13021303# Count the amount of new elements.1304new_id = 013051306# Permute dimensions so that the slice dimensions are always the in correct places for later use.1307if slice === :x1308original_nodes = permutedims(original_nodes, [1, 3, 4, 2, 5])1309unstructured_data = permutedims(unstructured_data, [2, 3, 1, 4, 5])1310elseif slice === :y1311original_nodes = permutedims(original_nodes, [1, 2, 4, 3, 5])1312unstructured_data = permutedims(unstructured_data, [1, 3, 2, 4, 5])1313end13141315# Iterate over all elements to find the ones that lie on the slice axis.1316for element_id in 1:n_elements1317min_coordinate = original_nodes[:, 1, 1, 1, element_id]1318max_coordinate = original_nodes[:, n_nodes_in, n_nodes_in, n_nodes_in,1319element_id]1320element_length = max_coordinate - min_coordinate13211322# Test if the element is on the slice axis. If not just continue with the next element.1323if !((prod(min_coordinate[other_dimensions] .<= point[other_dimensions]) &&1324prod(max_coordinate[other_dimensions] .> point[other_dimensions])) ||1325(point[other_dimensions] == upper_limit &&1326prod(max_coordinate[other_dimensions] .== upper_limit)))1327continue1328end13291330new_id += 113311332# Construct vandermonde matrix for interpolation of each 2D element to a 1D element.1333normalized_intercept = (point[other_dimensions] .-1334min_coordinate[other_dimensions]) /1335element_length[1] * 2 .- 11336vandermonde_i = polynomial_interpolation_matrix(nodes_in,1337normalized_intercept[1])1338vandermonde_ii = polynomial_interpolation_matrix(nodes_in,1339normalized_intercept[2])13401341# Interpolate to each node of new 1D element.1342for v in 1:n_variables1343for i in 1:n_nodes_in1344for ii in 1:n_nodes_in1345temp_unstructured_data[i, ii, new_id, v] = (vandermonde_ii * unstructured_data[ii,1346:,1347i,1348element_id,1349v])[1]1350end1351new_unstructured_data[i, new_id, v] = (vandermonde_i * temp_unstructured_data[i,1352:,1353new_id,1354v])[1]1355end1356end13571358new_nodes[:, new_id] = original_nodes[slice_dimension, 1, 1, :, element_id]1359end13601361return get_data_1d(reshape(new_nodes[:, 1:new_id], 1, n_nodes_in, new_id),1362new_unstructured_data[:, 1:new_id, :], nvisnodes, reinterpolate)1363end13641365# Interpolate unstructured DG data to structured data (cell-centered)1366#1367# This function takes DG data in an unstructured, Cartesian layout and converts it to a uniformly1368# distributed 2D layout.1369#1370# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may1371# thus be changed in future releases.1372function unstructured2structured(unstructured_data, normalized_coordinates,1373levels, resolution, nvisnodes_per_level)1374# Extract data shape information1375n_nodes_in, _, n_elements, n_variables = size(unstructured_data)13761377# Get node coordinates for DG locations on reference element1378nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)13791380# Calculate interpolation vandermonde matrices for each level1381max_level = length(nvisnodes_per_level) - 11382vandermonde_per_level = []1383for l in 0:max_level1384n_nodes_out = nvisnodes_per_level[l + 1]1385dx = 2 / n_nodes_out1386nodes_out = collect(range(-1 + dx / 2, 1 - dx / 2, length = n_nodes_out))1387push!(vandermonde_per_level,1388polynomial_interpolation_matrix(nodes_in, nodes_out))1389end13901391# For each element, calculate index position at which to insert data in global data structure1392lower_left_index = element2index(normalized_coordinates, levels, resolution,1393nvisnodes_per_level)13941395# Create output data structure1396structured = [Matrix{Float64}(undef, resolution, resolution) for _ in 1:n_variables]13971398# For each variable, interpolate element data and store to global data structure1399for v in 1:n_variables1400# Reshape data array for use in multiply_dimensionwise function1401reshaped_data = reshape(unstructured_data[:, :, :, v], 1, n_nodes_in,1402n_nodes_in, n_elements)14031404for element_id in 1:n_elements1405# Extract level for convenience1406level = levels[element_id]14071408# Determine target indices1409n_nodes_out = nvisnodes_per_level[level + 1]1410first = lower_left_index[:, element_id]1411last = first .+ (n_nodes_out - 1)14121413# Interpolate data1414vandermonde = vandermonde_per_level[level + 1]1415structured[v][first[1]:last[1], first[2]:last[2]] .= (reshape(multiply_dimensionwise(vandermonde,1416reshaped_data[:,1417:,1418:,1419element_id]),1420n_nodes_out,1421n_nodes_out))1422end1423end14241425return structured1426end14271428# For a given normalized element coordinate, return the index of its lower left1429# contribution to the global data structure1430#1431# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may1432# thus be changed in future releases.1433function element2index(normalized_coordinates, levels, resolution, nvisnodes_per_level)1434@assert size(normalized_coordinates, 1)==2 "only works in 2D"14351436n_elements = length(levels)14371438# First, determine lower left coordinate for all cells1439dx = 2 / resolution1440ndim = 21441lower_left_coordinate = Array{Float64}(undef, ndim, n_elements)1442for element_id in 1:n_elements1443nvisnodes = nvisnodes_per_level[levels[element_id] + 1]1444lower_left_coordinate[1, element_id] = (normalized_coordinates[1, element_id] -1445(nvisnodes - 1) / 2 * dx)1446lower_left_coordinate[2, element_id] = (normalized_coordinates[2, element_id] -1447(nvisnodes - 1) / 2 * dx)1448end14491450# Then, convert coordinate to global index1451indices = coordinate2index(lower_left_coordinate, resolution)14521453return indices1454end14551456# Find 2D array index for a 2-tuple of normalized, cell-centered coordinates (i.e., in [-1,1])1457#1458# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may1459# thus be changed in future releases.1460function coordinate2index(coordinate, resolution::Integer)1461# Calculate 1D normalized coordinates1462dx = 2 / resolution1463mesh_coordinates = collect(range(-1 + dx / 2, 1 - dx / 2, length = resolution))14641465# Find index1466id_x = searchsortedfirst.(Ref(mesh_coordinates), coordinate[1, :],1467lt = (x, y) -> x .< y .- dx / 2)1468id_y = searchsortedfirst.(Ref(mesh_coordinates), coordinate[2, :],1469lt = (x, y) -> x .< y .- dx / 2)1470return transpose(hcat(id_x, id_y))1471end14721473# Calculate the vertices for each mesh cell such that it can be visualized as a closed box1474#1475# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may1476# thus be changed in future releases.1477function calc_vertices(coordinates, levels, length_level_0)1478ndim = size(coordinates, 1)1479@assert ndim==2 "only works in 2D"14801481# Initialize output arrays1482n_elements = length(levels)1483n_points_per_element = 2^ndim + 21484x = Vector{Float64}(undef, n_points_per_element * n_elements)1485y = Vector{Float64}(undef, n_points_per_element * n_elements)14861487# Calculate vertices for all coordinates at once1488for element_id in 1:n_elements1489length = length_level_0 / 2^levels[element_id]1490index = n_points_per_element * (element_id - 1)1491x[index + 1] = coordinates[1, element_id] - 1 / 2 * length1492x[index + 2] = coordinates[1, element_id] + 1 / 2 * length1493x[index + 3] = coordinates[1, element_id] + 1 / 2 * length1494x[index + 4] = coordinates[1, element_id] - 1 / 2 * length1495x[index + 5] = coordinates[1, element_id] - 1 / 2 * length1496x[index + 6] = NaN14971498y[index + 1] = coordinates[2, element_id] - 1 / 2 * length1499y[index + 2] = coordinates[2, element_id] - 1 / 2 * length1500y[index + 3] = coordinates[2, element_id] + 1 / 2 * length1501y[index + 4] = coordinates[2, element_id] + 1 / 2 * length1502y[index + 5] = coordinates[2, element_id] - 1 / 2 * length1503y[index + 6] = NaN1504end15051506return x, y1507end15081509# Calculate the vertices to plot each grid line for StructuredMesh1510#1511# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may1512# thus be changed in future releases.1513function calc_vertices(node_coordinates, mesh)1514@unpack cells_per_dimension = mesh1515@assert size(node_coordinates, 1)==2 "only works in 2D"15161517linear_indices = LinearIndices(size(mesh))15181519# Initialize output arrays1520n_lines = sum(cells_per_dimension) + 21521max_length = maximum(cells_per_dimension)1522n_nodes = size(node_coordinates, 2)15231524# Create output as two matrices `x` and `y`, each holding the node locations for each of the `n_lines` grid lines1525# The # of rows in the matrices must be sufficient to store the longest dimension (`max_length`),1526# and for each the node locations without doubling the corner nodes (`n_nodes-1`), plus the final node (`+1`)1527# Rely on Plots.jl to ignore `NaN`s (i.e., they are not plotted) to handle shorter lines1528x = fill(NaN, max_length * (n_nodes - 1) + 1, n_lines)1529y = fill(NaN, max_length * (n_nodes - 1) + 1, n_lines)15301531line_index = 11532# Lines in x-direction1533# Bottom boundary1534i = 11535for cell_x in axes(mesh, 1)1536for node in 1:(n_nodes - 1)1537x[i, line_index] = node_coordinates[1, node, 1, linear_indices[cell_x, 1]]1538y[i, line_index] = node_coordinates[2, node, 1, linear_indices[cell_x, 1]]15391540i += 11541end1542end1543# Last point on bottom boundary1544x[i, line_index] = node_coordinates[1, end, 1, linear_indices[end, 1]]1545y[i, line_index] = node_coordinates[2, end, 1, linear_indices[end, 1]]15461547# Other lines in x-direction1548line_index += 11549for cell_y in axes(mesh, 2)1550i = 11551for cell_x in axes(mesh, 1)1552for node in 1:(n_nodes - 1)1553x[i, line_index] = node_coordinates[1, node, end,1554linear_indices[cell_x, cell_y]]1555y[i, line_index] = node_coordinates[2, node, end,1556linear_indices[cell_x, cell_y]]15571558i += 11559end1560end1561# Last point on line1562x[i, line_index] = node_coordinates[1, end, end, linear_indices[end, cell_y]]1563y[i, line_index] = node_coordinates[2, end, end, linear_indices[end, cell_y]]15641565line_index += 11566end15671568# Lines in y-direction1569# Left boundary1570i = 11571for cell_y in axes(mesh, 2)1572for node in 1:(n_nodes - 1)1573x[i, line_index] = node_coordinates[1, 1, node, linear_indices[1, cell_y]]1574y[i, line_index] = node_coordinates[2, 1, node, linear_indices[1, cell_y]]15751576i += 11577end1578end1579# Last point on left boundary1580x[i, line_index] = node_coordinates[1, 1, end, linear_indices[1, end]]1581y[i, line_index] = node_coordinates[2, 1, end, linear_indices[1, end]]15821583# Other lines in y-direction1584line_index += 11585for cell_x in axes(mesh, 1)1586i = 11587for cell_y in axes(mesh, 2)1588for node in 1:(n_nodes - 1)1589x[i, line_index] = node_coordinates[1, end, node,1590linear_indices[cell_x, cell_y]]1591y[i, line_index] = node_coordinates[2, end, node,1592linear_indices[cell_x, cell_y]]15931594i += 11595end1596end1597# Last point on line1598x[i, line_index] = node_coordinates[1, end, end, linear_indices[cell_x, end]]1599y[i, line_index] = node_coordinates[2, end, end, linear_indices[cell_x, end]]16001601line_index += 11602end16031604return x, y1605end16061607# Convert `slice` to orientations (1 -> `x`, 2 -> `y`, 3 -> `z`) for the two axes in a 2D plot1608function _get_orientations(mesh, slice)1609if ndims(mesh) == 2 || (ndims(mesh) == 3 && slice === :xy)1610orientation_x = 11611orientation_y = 21612elseif ndims(mesh) == 3 && slice === :xz1613orientation_x = 11614orientation_y = 31615elseif ndims(mesh) == 3 && slice === :yz1616orientation_x = 21617orientation_y = 31618else1619orientation_x = 01620orientation_y = 01621end1622return orientation_x, orientation_y1623end16241625# Convert `orientation` into a guide label (see also `_get_orientations`)1626function _get_guide(orientation::Integer)1627if orientation == 11628return "\$x\$"1629elseif orientation == 21630return "\$y\$"1631elseif orientation == 31632return "\$z\$"1633else1634return ""1635end1636end16371638# plotting_interpolation_matrix(dg; kwargs...)1639#1640# Interpolation matrix which maps discretization nodes to a set of plotting nodes.1641# Defaults to the identity matrix of size `length(solver.basis.nodes)`, and interpolates1642# to equispaced nodes for DGSEM (set by kwarg `nvisnodes` in the plotting function).1643#1644# Example:1645# ```julia1646# A = plotting_interpolation_matrix(dg)1647# A * dg.basis.nodes # => vector of nodes at which to plot the solution1648# ```1649#1650# Note: we cannot use UniformScaling to define the interpolation matrix since we use it with `kron`1651# to define a multi-dimensional interpolation matrix later.1652plotting_interpolation_matrix(dg; kwargs...) = I(length(dg.basis.nodes))16531654function face_plotting_interpolation_matrix(dg::DGSEM;1655nvisnodes = 2 * length(dg.basis.nodes))1656return polynomial_interpolation_matrix(dg.basis.nodes, LinRange(-1, 1, nvisnodes))1657end16581659function plotting_interpolation_matrix(dg::DGSEM;1660nvisnodes = 2 * length(dg.basis.nodes))1661Vp1D = polynomial_interpolation_matrix(dg.basis.nodes, LinRange(-1, 1, nvisnodes))1662# For quadrilateral elements, interpolation to plotting nodes involves applying a 1D interpolation1663# operator to each line of nodes. This is equivalent to multiplying the vector containing all node1664# node coordinates on an element by a Kronecker product of the 1D interpolation operator (e.g., a1665# multi-dimensional interpolation operator).1666return kron(Vp1D, Vp1D)1667end16681669function reference_node_coordinates_2d(dg::DGSEM)1670@unpack nodes = dg.basis1671r = vec([nodes[i] for i in eachnode(dg), j in eachnode(dg)])1672s = vec([nodes[j] for i in eachnode(dg), j in eachnode(dg)])1673return r, s1674end16751676# Find element and triangle ids containing coordinates given as a matrix [ndims, npoints]1677function get_ids_by_coordinates!(ids, coordinates, pd)1678if length(ids) != 2 * size(coordinates, 2)1679throw(DimensionMismatch("storage length for element ids does not match the number of coordinates"))1680end16811682n_coordinates = size(coordinates, 2)16831684for index in 1:n_coordinates1685point = SVector(coordinates[1, index], coordinates[2, index])1686ids[index, :] .= find_element(point, pd)1687end16881689return ids1690end16911692# Find the ids of elements and triangles containing given coordinates by using the triangulation in 'pd'.1693function get_ids_by_coordinates(coordinates, pd)1694ids = Matrix{Int}(undef, size(coordinates, 2), 2)1695get_ids_by_coordinates!(ids, coordinates, pd)1696return ids1697end16981699# Check if given 'point' is inside the triangle with corners corresponding to the coordinates of x and y.1700function is_in_triangle(point, x, y)1701a = SVector(x[1], y[1])1702b = SVector(x[2], y[2])1703c = SVector(x[3], y[3])1704return is_on_same_side(point, a, b, c) && is_on_same_side(point, b, c, a) &&1705is_on_same_side(point, c, a, b)1706end17071708# Create an axis through x and y to then check if 'point' is on the same side of the axis as z.1709function is_on_same_side(point, x, y, z)1710if (y[1] - x[1]) == 01711return (point[1] - x[1]) * (z[1] - x[1]) >= 01712else1713a = (y[2] - x[2]) / (y[1] - x[1])1714b = x[2] - a * x[1]1715return (z[2] - a * z[1] - b) * (point[2] - a * point[1] - b) >= 01716end1717end17181719# For a given 'point', return the id of the element it is contained in in; if not found return 0.1720function find_element(point, pd)1721n_tri = size(pd.t, 1)1722n_elements = size(pd.x, 2)17231724# Iterate over all elements.1725for element in 1:n_elements1726# Iterate over all triangles in given element.1727for tri in 1:n_tri1728# The code below is equivalent to1729# x == pd.x[pd.t[tri, :], element]1730# y == pd.y[pd.t[tri, :], element]1731# but avoids allocations and is thus more efficient.1732tri_indices = (pd.t[tri, 1], pd.t[tri, 2], pd.t[tri, 3])1733x = SVector(pd.x[tri_indices[1], element],1734pd.x[tri_indices[2], element],1735pd.x[tri_indices[3], element])1736y = SVector(pd.y[tri_indices[1], element],1737pd.y[tri_indices[2], element],1738pd.y[tri_indices[3], element])1739if is_in_triangle(point, x, y)1740return (element, tri)1741end1742end1743end17441745return (0, 0)1746end17471748# Interpolate from three corners of a triangle to a single point.1749function triangle_interpolation(x_coordinates_in, y_coordinates_in, values_in,1750coordinate_out)1751A = hcat(x_coordinates_in, y_coordinates_in, SVector(1, 1, 1))1752c = A \ values_in1753return c[1] * coordinate_out[1] + c[2] * coordinate_out[2] + c[3]1754end17551756# Create an axis.1757function axis_curve(nodes_x, nodes_y, nodes_z, slice, point, n_points)1758if n_points === nothing1759n_points = 641760end1761dimensions = length(point)1762curve = zeros(dimensions, n_points)1763if slice == :x1764xmin, xmax = extrema(nodes_x)1765curve[1, :] .= range(xmin, xmax, length = n_points)1766curve[2, :] .= point[2]1767if dimensions === 31768curve[3, :] .= point[3]1769end1770elseif slice == :y1771ymin, ymax = extrema(nodes_y)1772curve[1, :] .= point[1]1773curve[2, :] .= range(ymin, ymax, length = n_points)1774if dimensions === 31775curve[3, :] .= point[3]1776end1777elseif slice == :z1778zmin, zmax = extrema(nodes_z)1779curve[1, :] .= point[1]1780curve[2, :] .= point[2]1781curve[3, :] .= range(zmin, zmax, length = n_points)1782else1783@assert false "Input for 'slice' is not supported here."1784end17851786return curve1787end1788end # @muladd178917901791