Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/visualization/utilities.jl
2055 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
@inline num_faces(elem::Tri) = 3
9
@inline num_faces(elem::Quad) = 4
10
11
# compute_triangle_area(tri)
12
#
13
# Computes the area of a triangle given `tri`, which is a tuple of three points (vectors),
14
# using the [Shoelace_formula](https://en.wikipedia.org/wiki/Shoelace_formula).
15
function compute_triangle_area(tri)
16
A, B, C = tri
17
return 0.5f0 * (A[1] * (B[2] - C[2]) + B[1] * (C[2] - A[2]) + C[1] * (A[2] - B[2]))
18
end
19
20
# reference_plotting_triangulation(reference_plotting_coordinates)
21
#
22
# Computes a triangulation of the points in `reference_plotting_coordinates`, which is a tuple containing
23
# vectors of plotting points on the reference element (e.g., reference_plotting_coordinates = (r,s)).
24
# The reference element is assumed to be [-1,1]^d.
25
#
26
# This function returns `t` which is a `3 x N_tri` Matrix{Int} containing indices of triangles in the
27
# triangulation of the plotting points, with zero-volume triangles removed.
28
#
29
# For example, r[t[1, i]] returns the first reference coordinate of the 1st point on the ith triangle.
30
function reference_plotting_triangulation(reference_plotting_coordinates,
31
tol = 50 * eps())
32
# on-the-fly triangulation of plotting nodes on the reference element
33
tri_in = Triangulate.TriangulateIO()
34
tri_in.pointlist = permutedims(hcat(reference_plotting_coordinates...))
35
tri_out, _ = Triangulate.triangulate("Q", tri_in)
36
triangles = tri_out.trianglelist
37
38
# filter out sliver triangles
39
has_volume = fill(true, size(triangles, 2))
40
for i in axes(triangles, 2)
41
ids = @view triangles[:, i]
42
x_points = @view tri_out.pointlist[1, ids]
43
y_points = @view tri_out.pointlist[2, ids]
44
area = compute_triangle_area(zip(x_points, y_points))
45
if abs(area) < tol
46
has_volume[i] = false
47
end
48
end
49
return permutedims(triangles[:, findall(has_volume)])
50
end
51
52
# This function is used to avoid type instabilities when calling `digest_solution_variables`.
53
function transform_to_solution_variables!(u, solution_variables, equations)
54
for (i, u_i) in enumerate(u)
55
u[i] = solution_variables(u_i, equations)
56
end
57
end
58
59
# global_plotting_triangulation_triplot(u_plot, rst_plot, xyz_plot)
60
#
61
# Returns (plotting_coordinates_x, plotting_coordinates_y, ..., plotting_values, plotting_triangulation).
62
# Output can be used with TriplotRecipes.DGTriPseudocolor(...).
63
#
64
# Inputs:
65
# - xyz_plot = plotting points (tuple of matrices of size (Nplot, K))
66
# - u_plot = matrix of size (Nplot, K) representing solution to plot.
67
# - t = triangulation of reference plotting points
68
function global_plotting_triangulation_triplot(xyz_plot, u_plot, t)
69
@assert size(first(xyz_plot), 1)==size(u_plot, 1) "Row dimension of u_plot does not match row dimension of xyz_plot"
70
71
# build discontinuous data on plotting triangular mesh
72
num_plotting_points, num_elements = size(u_plot)
73
num_reference_plotting_triangles = size(t, 1)
74
num_plotting_elements_total = num_reference_plotting_triangles * num_elements
75
76
# each column of `tp` corresponds to a vertex of a plotting triangle
77
tp = zeros(Int32, 3, num_plotting_elements_total)
78
zp = similar(tp, eltype(u_plot))
79
for e in 1:num_elements
80
for i in 1:num_reference_plotting_triangles
81
tp[:, i + (e - 1) * num_reference_plotting_triangles] .= @views t[i, :] .+
82
(e - 1) *
83
num_plotting_points
84
zp[:, i + (e - 1) * num_reference_plotting_triangles] .= @views u_plot[t[i,
85
:],
86
e]
87
end
88
end
89
return vec.(xyz_plot)..., zp, tp
90
end
91
92
function get_face_node_indices(r, s, dg::DGSEM, tol = 100 * eps())
93
face_1 = findall(@. abs(s + 1) < tol)
94
face_2 = findall(@. abs(r - 1) < tol)
95
face_3 = findall(@. abs(s - 1) < tol)
96
face_4 = findall(@. abs(r + 1) < tol)
97
Fmask = hcat(face_1, face_2, face_3, face_4)
98
return Fmask
99
end
100
101
# dispatch on semi
102
function mesh_plotting_wireframe(u, semi)
103
mesh_plotting_wireframe(u, mesh_equations_solver_cache(semi)...)
104
end
105
106
# mesh_plotting_wireframe(u, mesh, equations, dg::DGMulti, cache; num_plotting_pts=25)
107
#
108
# Generates data for plotting a mesh wireframe given StartUpDG data types.
109
# Returns (plotting_coordinates_x, plotting_coordinates_y, nothing) for a 2D mesh wireframe.
110
function mesh_plotting_wireframe(u::StructArray, mesh, equations, dg::DGMulti, cache;
111
nvisnodes = 2 * nnodes(dg))
112
@unpack md = mesh
113
rd = dg.basis
114
115
# Construct 1D plotting interpolation matrix `Vp1D` for a single face
116
@unpack N, Fmask = rd
117
num_face_points = length(Fmask) ÷ num_faces(rd.element_type)
118
vandermonde_matrix_1D = StartUpDG.vandermonde(Line(), N,
119
StartUpDG.nodes(Line(),
120
num_face_points - 1))
121
rplot = LinRange(-1, 1, nvisnodes)
122
Vp1D = StartUpDG.vandermonde(Line(), N, rplot) / vandermonde_matrix_1D
123
124
num_faces_total = num_faces(rd.element_type) * md.num_elements
125
xf, yf = map(x -> reshape(view(x, Fmask, :), num_face_points, num_faces_total),
126
md.xyz)
127
uf = similar(u, size(xf))
128
apply_to_each_field((out, x) -> out .= reshape(view(x, Fmask, :), num_face_points,
129
num_faces_total), uf, u)
130
131
num_face_plotting_points = size(Vp1D, 1)
132
x_mesh, y_mesh = ntuple(_ -> zeros(num_face_plotting_points, num_faces_total), 2)
133
u_mesh = similar(u, (num_face_plotting_points, num_faces_total))
134
for f in 1:num_faces_total
135
mul!(view(x_mesh, :, f), Vp1D, view(xf, :, f))
136
mul!(view(y_mesh, :, f), Vp1D, view(yf, :, f))
137
apply_to_each_field(mul_by!(Vp1D), view(u_mesh, :, f), view(uf, :, f))
138
end
139
140
return x_mesh, y_mesh, u_mesh
141
end
142
143
function mesh_plotting_wireframe(u::StructArray, mesh, equations, dg::DGSEM, cache;
144
nvisnodes = 2 * nnodes(dg))
145
146
# build nodes on reference element (seems to be the right ordering)
147
r, s = reference_node_coordinates_2d(dg)
148
149
# extract node coordinates
150
uEltype = eltype(first(u))
151
nvars = nvariables(equations)
152
n_nodes_2d = nnodes(dg)^ndims(mesh)
153
n_elements = nelements(dg, cache)
154
x = reshape(view(cache.elements.node_coordinates, 1, :, :, :), n_nodes_2d,
155
n_elements)
156
y = reshape(view(cache.elements.node_coordinates, 2, :, :, :), n_nodes_2d,
157
n_elements)
158
159
# extract indices of local face nodes for wireframe plotting
160
Fmask = get_face_node_indices(r, s, dg)
161
plotting_interp_matrix1D = face_plotting_interpolation_matrix(dg;
162
nvisnodes = nvisnodes)
163
164
# These 5 lines extract the face values on each element from the arrays x,y,sol_to_plot.
165
# The resulting arrays are then reshaped so that xf, yf, sol_f are Matrix types of size
166
# (Number of face plotting nodes) x (Number of faces).
167
function face_first_reshape(x, num_nodes_1D, num_nodes, num_elements)
168
num_reference_faces = 2 * ndims(mesh)
169
xf = view(reshape(x, num_nodes, num_elements), vec(Fmask), :)
170
return reshape(xf, num_nodes_1D, num_elements * num_reference_faces)
171
end
172
function reshape_and_interpolate(x)
173
plotting_interp_matrix1D *
174
face_first_reshape(x, nnodes(dg), n_nodes_2d, n_elements)
175
end
176
xfp, yfp = map(reshape_and_interpolate, (x, y))
177
ufp = StructArray{SVector{nvars, uEltype}}(map(reshape_and_interpolate,
178
StructArrays.components(u)))
179
180
return xfp, yfp, ufp
181
end
182
183
function mesh_plotting_wireframe(u::ScalarData, mesh, equations, dg::DGSEM, cache;
184
nvisnodes = 2 * nnodes(dg))
185
186
# build nodes on reference element (seems to be the right ordering)
187
r, s = reference_node_coordinates_2d(dg)
188
189
# extract node coordinates
190
n_nodes_2d = nnodes(dg)^ndims(mesh)
191
n_elements = nelements(dg, cache)
192
x = reshape(view(cache.elements.node_coordinates, 1, :, :, :), n_nodes_2d,
193
n_elements)
194
y = reshape(view(cache.elements.node_coordinates, 2, :, :, :), n_nodes_2d,
195
n_elements)
196
197
# extract indices of local face nodes for wireframe plotting
198
Fmask = get_face_node_indices(r, s, dg)
199
plotting_interp_matrix1D = face_plotting_interpolation_matrix(dg;
200
nvisnodes = nvisnodes)
201
202
# These 5 lines extract the face values on each element from the arrays x,y,sol_to_plot.
203
# The resulting arrays are then reshaped so that xf, yf, sol_f are Matrix types of size
204
# (Number of face plotting nodes) x (Number of faces).
205
function face_first_reshape(x, num_nodes_1D, num_nodes, num_elements)
206
num_reference_faces = 2 * ndims(mesh)
207
xf = view(reshape(x, num_nodes, num_elements), vec(Fmask), :)
208
return reshape(xf, num_nodes_1D, num_elements * num_reference_faces)
209
end
210
function reshape_and_interpolate(x)
211
plotting_interp_matrix1D *
212
face_first_reshape(x, nnodes(dg), n_nodes_2d, n_elements)
213
end
214
xfp, yfp, ufp = map(reshape_and_interpolate, (x, y, u.data))
215
216
return xfp, yfp, ufp
217
end
218
219
function mesh_plotting_wireframe(u::ScalarData, mesh, equations, dg::DGMulti, cache;
220
nvisnodes = 2 * nnodes(dg))
221
@unpack md = mesh
222
rd = dg.basis
223
224
# Construct 1D plotting interpolation matrix `Vp1D` for a single face
225
@unpack N, Fmask = rd
226
vandermonde_matrix_1D = StartUpDG.vandermonde(Line(), N, StartUpDG.nodes(Line(), N))
227
rplot = LinRange(-1, 1, nvisnodes)
228
Vp1D = StartUpDG.vandermonde(Line(), N, rplot) / vandermonde_matrix_1D
229
230
num_face_points = N + 1
231
num_faces_total = num_faces(rd.element_type) * md.num_elements
232
xf, yf, uf = map(x -> reshape(view(x, Fmask, :), num_face_points, num_faces_total),
233
(md.xyz..., u.data))
234
235
num_face_plotting_points = size(Vp1D, 1)
236
x_mesh, y_mesh = ntuple(_ -> zeros(num_face_plotting_points, num_faces_total), 2)
237
u_mesh = similar(u.data, (num_face_plotting_points, num_faces_total))
238
for f in 1:num_faces_total
239
mul!(view(x_mesh, :, f), Vp1D, view(xf, :, f))
240
mul!(view(y_mesh, :, f), Vp1D, view(yf, :, f))
241
mul!(view(u_mesh, :, f), Vp1D, view(uf, :, f))
242
end
243
244
return x_mesh, y_mesh, u_mesh
245
end
246
247
# These methods are used internally to set the default value of the solution variables:
248
# - If a `cons2prim` for the given `equations` exists, use it
249
# - Otherwise, use `cons2cons`, which is defined for all systems of equations
250
digest_solution_variables(equations, solution_variables) = solution_variables
251
function digest_solution_variables(equations, solution_variables::Nothing)
252
if hasmethod(cons2prim, Tuple{AbstractVector, typeof(equations)})
253
return cons2prim
254
else
255
return cons2cons
256
end
257
end
258
259
digest_variable_names(solution_variables_, equations, variable_names) = variable_names
260
digest_variable_names(solution_variables_, equations, ::Nothing) = SVector(varnames(solution_variables_,
261
equations))
262
263
"""
264
adapt_to_mesh_level!(u_ode, semi, level)
265
adapt_to_mesh_level!(sol::Trixi.TrixiODESolution, level)
266
267
Like [`adapt_to_mesh_level`](@ref), but modifies the solution and parts of the
268
semidiscretization (mesh and caches) in place.
269
"""
270
function adapt_to_mesh_level!(u_ode, semi, level)
271
# Create AMR callback with controller that refines everything towards a single level
272
amr_controller = ControllerThreeLevel(semi, IndicatorMax(semi, variable = first),
273
base_level = level)
274
amr_callback = AMRCallback(semi, amr_controller, interval = 0)
275
276
# Adapt mesh until it does not change anymore
277
has_changed = amr_callback.affect!(u_ode, semi, 0.0, 0)
278
while has_changed
279
has_changed = amr_callback.affect!(u_ode, semi, 0.0, 0)
280
end
281
282
return u_ode, semi
283
end
284
285
function adapt_to_mesh_level!(sol::TrixiODESolution, level)
286
adapt_to_mesh_level!(sol.u[end], sol.prob.p, level)
287
end
288
289
"""
290
adapt_to_mesh_level(u_ode, semi, level)
291
adapt_to_mesh_level(sol::Trixi.TrixiODESolution, level)
292
293
Use the regular adaptive mesh refinement routines to adaptively refine/coarsen the solution `u_ode`
294
with semidiscretization `semi` towards a uniformly refined grid with refinement level `level`. The
295
solution and semidiscretization are copied such that the original objects remain *unaltered*.
296
297
A convenience method accepts an ODE solution object, from which solution and semidiscretization are
298
extracted as needed.
299
300
See also: [`adapt_to_mesh_level!`](@ref)
301
"""
302
function adapt_to_mesh_level(u_ode, semi, level)
303
# Create new semidiscretization with copy of the current mesh
304
mesh, _, _, _ = mesh_equations_solver_cache(semi)
305
new_semi = remake(semi, mesh = deepcopy(mesh))
306
307
return adapt_to_mesh_level!(deepcopy(u_ode), new_semi, level)
308
end
309
310
function adapt_to_mesh_level(sol::TrixiODESolution, level)
311
adapt_to_mesh_level(sol.u[end], sol.prob.p, level)
312
end
313
314
# Extract data from a 2D/3D DG solution and prepare it for visualization as a heatmap/contour plot.
315
#
316
# Returns a tuple with
317
# - x coordinates
318
# - y coordinates
319
# - nodal 2D data
320
# - x vertices for mesh visualization
321
# - y vertices for mesh visualization
322
#
323
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
324
# thus be changed in future releases.
325
function get_data_2d(center_level_0, length_level_0, leaf_cells, coordinates, levels,
326
ndims, unstructured_data, n_nodes,
327
grid_lines = false, max_supported_level = 11, nvisnodes = nothing,
328
slice = :xy, point = (0.0, 0.0, 0.0))
329
# Determine resolution for data interpolation
330
max_level = maximum(levels)
331
if max_level > max_supported_level
332
error("Maximum refinement level $max_level is higher than " *
333
"maximum supported level $max_supported_level")
334
end
335
max_available_nodes_per_finest_element = 2^(max_supported_level - max_level)
336
if nvisnodes === nothing
337
max_nvisnodes = 2 * n_nodes
338
elseif nvisnodes == 0
339
max_nvisnodes = n_nodes
340
else
341
max_nvisnodes = nvisnodes
342
end
343
nvisnodes_at_max_level = min(max_available_nodes_per_finest_element, max_nvisnodes)
344
resolution = nvisnodes_at_max_level * 2^max_level
345
nvisnodes_per_level = [2^(max_level - level) * nvisnodes_at_max_level
346
for level in 0:max_level]
347
# nvisnodes_per_level is an array (accessed by "level + 1" to accommodate
348
# level-0-cell) that contains the number of visualization nodes for any
349
# refinement level to visualize on an equidistant grid
350
351
if ndims == 3
352
(unstructured_data, coordinates, levels,
353
center_level_0) = unstructured_3d_to_2d(unstructured_data,
354
coordinates, levels, length_level_0,
355
center_level_0, slice,
356
point)
357
end
358
359
# Normalize element coordinates: move center to (0, 0) and domain size to [-1, 1]²
360
n_elements = length(levels)
361
normalized_coordinates = similar(coordinates)
362
for element_id in 1:n_elements
363
@views normalized_coordinates[:, element_id] .= ((coordinates[:, element_id] .-
364
center_level_0) ./
365
(length_level_0 / 2))
366
end
367
368
# Interpolate unstructured DG data to structured data
369
(structured_data = unstructured2structured(unstructured_data,
370
normalized_coordinates,
371
levels, resolution, nvisnodes_per_level))
372
373
# Interpolate cell-centered values to node-centered values
374
node_centered_data = cell2node(structured_data)
375
376
# Determine axis coordinates for contour plot
377
xs = collect(range(-1, 1, length = resolution + 1)) .* length_level_0 / 2 .+
378
center_level_0[1]
379
ys = collect(range(-1, 1, length = resolution + 1)) .* length_level_0 / 2 .+
380
center_level_0[2]
381
382
# Determine element vertices to plot grid lines
383
if grid_lines
384
mesh_vertices_x, mesh_vertices_y = calc_vertices(coordinates, levels,
385
length_level_0)
386
else
387
mesh_vertices_x = Vector{Float64}(undef, 0)
388
mesh_vertices_y = Vector{Float64}(undef, 0)
389
end
390
391
return xs, ys, node_centered_data, mesh_vertices_x, mesh_vertices_y
392
end
393
394
# For finite difference methods (FDSBP), we do not want to reinterpolate the data, but use the same
395
# nodes as input nodes. For DG methods, we usually want to reinterpolate the data on an equidistant grid.
396
default_reinterpolate(solver) = solver isa FDSBP ? false : true
397
398
# Extract data from a 1D DG solution and prepare it for visualization as a line plot.
399
# This returns a tuple with
400
# - x coordinates
401
# - nodal 1D data
402
#
403
# Note: This is a low-level function that is not considered as part of Trixi's interface and may
404
# thus be changed in future releases.
405
function get_data_1d(original_nodes, unstructured_data, nvisnodes, reinterpolate)
406
# 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.
407
n_nodes, n_elements, n_vars = size(unstructured_data)
408
409
# If `reinterpolate` is `false`, we do nothing.
410
# If `reinterpolate` is `true`, the output nodes are equidistantly spaced.
411
if reinterpolate == false
412
interpolated_nodes = original_nodes
413
interpolated_data = unstructured_data
414
else
415
# Set the amount of nodes visualized according to nvisnodes.
416
if nvisnodes === nothing
417
max_nvisnodes = 2 * n_nodes
418
elseif nvisnodes == 0
419
max_nvisnodes = n_nodes
420
else
421
@assert nvisnodes>=2 "nvisnodes must be zero or >= 2"
422
max_nvisnodes = nvisnodes
423
end
424
425
interpolated_nodes = Array{eltype(original_nodes), 2}(undef, max_nvisnodes,
426
n_elements)
427
interpolated_data = Array{eltype(unstructured_data), 3}(undef, max_nvisnodes,
428
n_elements, n_vars)
429
430
for j in 1:n_elements
431
# Interpolate on an equidistant grid.
432
interpolated_nodes[:, j] .= range(original_nodes[1, 1, j],
433
original_nodes[1, end, j],
434
length = max_nvisnodes)
435
end
436
437
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)
438
nodes_out = collect(range(-1, 1, length = max_nvisnodes))
439
440
# Calculate vandermonde matrix for interpolation.
441
vandermonde = polynomial_interpolation_matrix(nodes_in, nodes_out)
442
443
# Iterate over all variables.
444
for v in 1:n_vars
445
# Interpolate data for each element.
446
for element in 1:n_elements
447
multiply_scalar_dimensionwise!(@view(interpolated_data[:, element, v]),
448
vandermonde,
449
@view(unstructured_data[:, element, v]))
450
end
451
end
452
end
453
454
# Return results after data is reshaped
455
return vec(interpolated_nodes), reshape(interpolated_data, :, n_vars),
456
vcat(original_nodes[1, 1, :], original_nodes[1, end, end])
457
end
458
459
# Change order of dimensions (variables are now last) and convert data to `solution_variables`
460
#
461
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
462
# thus be changed in future releases.
463
function get_unstructured_data(u, solution_variables, mesh, equations, solver, cache)
464
if solution_variables === cons2cons
465
raw_data = u
466
n_vars = size(raw_data, 1)
467
else
468
# Similar to `save_solution_file` in `callbacks_step/save_solution_dg.jl`.
469
# However, we cannot use `reinterpret` here as `u` might have a non-bits type.
470
# See https://github.com/trixi-framework/Trixi.jl/pull/2388
471
n_vars_in = nvariables(equations)
472
n_vars = length(solution_variables(get_node_vars(u, equations, solver),
473
equations))
474
raw_data = Array{eltype(u)}(undef, n_vars, Base.tail(size(u))...)
475
reshaped_u = reshape(u, n_vars_in, :)
476
reshaped_r = reshape(raw_data, n_vars, :)
477
for idx in axes(reshaped_u, 2)
478
reshaped_r[:, idx] = solution_variables(get_node_vars(reshaped_u, equations,
479
solver, idx),
480
equations)
481
end
482
end
483
484
unstructured_data = Array{eltype(raw_data)}(undef,
485
ntuple((d) -> nnodes(solver),
486
ndims(equations))...,
487
nelements(solver, cache), n_vars)
488
for variable in 1:n_vars
489
@views unstructured_data[.., :, variable] .= raw_data[variable, .., :]
490
end
491
492
return unstructured_data
493
end
494
495
# This method is only for plotting 1D functions
496
function get_unstructured_data(func::Function, solution_variables,
497
mesh::AbstractMesh{1}, equations, solver, cache)
498
original_nodes = cache.elements.node_coordinates
499
# original_nodes has size (1, nnodes, nelements)
500
# we want u to have size (nvars, nnodes, nelements)
501
# func.(original_nodes, equations) has size (1, nnodes, nelements), where each component has length n_vars
502
# Therefore, we drop the first (singleton) dimension and then stack the components
503
u = stack(func.(SVector.(dropdims(original_nodes; dims = 1)), equations))
504
return get_unstructured_data(u, solution_variables, mesh, equations, solver, cache)
505
end
506
507
# Convert cell-centered values to node-centered values by averaging over all
508
# four neighbors. Solution values at the edges are padded with ghost values
509
# computed via linear extrapolation.
510
#
511
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
512
# thus be changed in future releases.
513
function cell2node(cell_centered_data)
514
# Create temporary data structure to make the averaging algorithm as simple
515
# as possible (by using a ghost layer)
516
tmp = similar(first(cell_centered_data), size(first(cell_centered_data)) .+ (2, 2))
517
518
# Create output data structure
519
resolution_in, _ = size(first(cell_centered_data))
520
resolution_out = resolution_in + 1
521
node_centered_data = [Matrix{Float64}(undef, resolution_out, resolution_out)
522
for _ in eachindex(cell_centered_data)]
523
524
for (cell_data, node_data) in zip(cell_centered_data, node_centered_data)
525
# Fill center with original data
526
tmp[2:(end - 1), 2:(end - 1)] .= cell_data
527
528
# Linear extrapolation of top and bottom rows
529
tmp[1, 2:(end - 1)] .= cell_data[1, :] .+ (cell_data[1, :] .- cell_data[2, :])
530
tmp[end, 2:(end - 1)] .= (cell_data[end, :] .+
531
(cell_data[end, :] .- cell_data[end - 1, :]))
532
533
# Linear extrapolatation of left and right columns
534
tmp[2:(end - 1), 1] .= cell_data[:, 1] .+ (cell_data[:, 1] .- cell_data[:, 2])
535
tmp[2:(end - 1), end] .= (cell_data[:, end] .+
536
(cell_data[:, end] .- cell_data[:, end - 1]))
537
538
# Corners perform the linear extrapolatation along diagonals
539
tmp[1, 1] = tmp[2, 2] + (tmp[2, 2] - tmp[3, 3])
540
tmp[1, end] = tmp[2, end - 1] + (tmp[2, end - 1] - tmp[3, end - 2])
541
tmp[end, 1] = tmp[end - 1, 2] + (tmp[end - 1, 2] - tmp[end - 2, 3])
542
tmp[end, end] = (tmp[end - 1, end - 1] +
543
(tmp[end - 1, end - 1] - tmp[end - 2, end - 2]))
544
545
# Obtain node-centered value by averaging over neighboring cell-centered values
546
for j in 1:resolution_out
547
for i in 1:resolution_out
548
node_data[i, j] = (tmp[i, j] +
549
tmp[i + 1, j] +
550
tmp[i, j + 1] +
551
tmp[i + 1, j + 1]) / 4
552
end
553
end
554
end
555
556
# Transpose
557
for (index, data) in enumerate(node_centered_data)
558
node_centered_data[index] = permutedims(data)
559
end
560
561
return node_centered_data
562
end
563
564
# Convert 3d unstructured data to 2d data.
565
# Additional to the new unstructured data updated coordinates, levels and
566
# center coordinates are returned.
567
#
568
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
569
# thus be changed in future releases.
570
function unstructured_3d_to_2d(unstructured_data, coordinates, levels,
571
length_level_0, center_level_0, slice,
572
point)
573
if slice === :yz
574
slice_dimension = 1
575
other_dimensions = [2, 3]
576
elseif slice === :xz
577
slice_dimension = 2
578
other_dimensions = [1, 3]
579
elseif slice === :xy
580
slice_dimension = 3
581
other_dimensions = [1, 2]
582
else
583
error("illegal dimension '$slice', supported dimensions are :yz, :xz, and :xy")
584
end
585
586
# Limits of domain in slice dimension
587
lower_limit = center_level_0[slice_dimension] - length_level_0 / 2
588
upper_limit = center_level_0[slice_dimension] + length_level_0 / 2
589
590
@assert length(point)>=3 "Point must be three-dimensional."
591
if point[slice_dimension] < lower_limit || point[slice_dimension] > upper_limit
592
error(string("Slice plane is outside of domain.",
593
" point[$slice_dimension]=$(point[slice_dimension]) must be between $lower_limit and $upper_limit"))
594
end
595
596
# Extract data shape information
597
n_nodes_in, _, _, n_elements, n_variables = size(unstructured_data)
598
599
# Get node coordinates for DG locations on reference element
600
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)
601
602
# New unstructured data has one dimension less.
603
# The redundant element ids are removed later.
604
@views new_unstructured_data = similar(unstructured_data[1, ..])
605
606
# Declare new empty arrays to fill in new coordinates and levels
607
new_coordinates = Array{Float64}(undef, 2, n_elements)
608
new_levels = Array{eltype(levels)}(undef, n_elements)
609
610
# Counter for new element ids
611
new_id = 0
612
613
# Save vandermonde matrices in a Dict to prevent redundant generation
614
vandermonde_to_2d = Dict()
615
616
# Permute dimensions such that the slice dimension is always the
617
# third dimension of the array. Below we can always interpolate in the
618
# third dimension.
619
if slice === :yz
620
unstructured_data = permutedims(unstructured_data, [2, 3, 1, 4, 5])
621
elseif slice === :xz
622
unstructured_data = permutedims(unstructured_data, [1, 3, 2, 4, 5])
623
end
624
625
for element_id in 1:n_elements
626
# Distance from center to border of this element (half the length)
627
element_length = length_level_0 / 2^levels[element_id]
628
min_coordinate = coordinates[:, element_id] .- element_length / 2
629
max_coordinate = coordinates[:, element_id] .+ element_length / 2
630
631
# Check if slice plane and current element intersect.
632
# The first check uses a "greater but not equal" to only match one cell if the
633
# slice plane lies between two cells.
634
# The second check is needed if the slice plane is at the upper border of
635
# the domain due to this.
636
if !((min_coordinate[slice_dimension] <= point[slice_dimension] &&
637
max_coordinate[slice_dimension] > point[slice_dimension]) ||
638
(point[slice_dimension] == upper_limit &&
639
max_coordinate[slice_dimension] == upper_limit))
640
# Continue for loop if they don't intersect
641
continue
642
end
643
644
# This element is of interest
645
new_id += 1
646
647
# Add element to new coordinates and levels
648
new_coordinates[:, new_id] = coordinates[other_dimensions, element_id]
649
new_levels[new_id] = levels[element_id]
650
651
# Construct vandermonde matrix (or load from Dict if possible)
652
normalized_intercept = (point[slice_dimension] -
653
min_coordinate[slice_dimension]) /
654
element_length * 2 - 1
655
656
if haskey(vandermonde_to_2d, normalized_intercept)
657
vandermonde = vandermonde_to_2d[normalized_intercept]
658
else
659
# Generate vandermonde matrix to interpolate values at nodes_in to one value
660
vandermonde = polynomial_interpolation_matrix(nodes_in,
661
[normalized_intercept])
662
vandermonde_to_2d[normalized_intercept] = vandermonde
663
end
664
665
# 1D interpolation to specified slice plane
666
# We permuted the dimensions above such that now the dimension in which
667
# we will interpolate is always the third one.
668
for i in 1:n_nodes_in
669
for ii in 1:n_nodes_in
670
# Interpolate in the third dimension
671
data = unstructured_data[i, ii, :, element_id, :]
672
673
value = multiply_dimensionwise(vandermonde, permutedims(data))
674
new_unstructured_data[i, ii, new_id, :] = value[:, 1]
675
end
676
end
677
end
678
679
# Remove redundant element ids
680
unstructured_data = new_unstructured_data[:, :, 1:new_id, :]
681
new_coordinates = new_coordinates[:, 1:new_id]
682
new_levels = new_levels[1:new_id]
683
684
center_level_0 = center_level_0[other_dimensions]
685
686
return unstructured_data, new_coordinates, new_levels, center_level_0
687
end
688
689
# Convert 2d unstructured data to 1d slice and interpolate them.
690
function unstructured_2d_to_1d(original_nodes, unstructured_data, nvisnodes,
691
reinterpolate, slice, point)
692
if slice === :x
693
slice_dimension = 2
694
other_dimension = 1
695
elseif slice === :y
696
slice_dimension = 1
697
other_dimension = 2
698
else
699
error("illegal dimension '$slice', supported dimensions are :x and :y")
700
end
701
702
# Set up data structures to store new 1D data.
703
@views new_unstructured_data = similar(unstructured_data[1, ..])
704
@views new_nodes = similar(original_nodes[1, 1, ..])
705
706
n_nodes_in, _, n_elements, n_variables = size(unstructured_data)
707
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)
708
709
# Test if point lies in the domain.
710
lower_limit = original_nodes[1, 1, 1, 1]
711
upper_limit = original_nodes[1, n_nodes_in, n_nodes_in, n_elements]
712
713
@assert length(point)>=2 "Point must be two-dimensional."
714
if point[slice_dimension] < lower_limit || point[slice_dimension] > upper_limit
715
error(string("Slice axis is outside of domain. ",
716
" point[$slice_dimension]=$(point[slice_dimension]) must be between $lower_limit and $upper_limit"))
717
end
718
719
# Count the amount of new elements.
720
new_id = 0
721
722
# Permute dimensions so that the slice dimension is always in the correct place for later use.
723
if slice === :y
724
original_nodes = permutedims(original_nodes, [1, 3, 2, 4])
725
unstructured_data = permutedims(unstructured_data, [2, 1, 3, 4])
726
end
727
728
# Iterate over all elements to find the ones that lie on the slice axis.
729
for element_id in 1:n_elements
730
min_coordinate = original_nodes[:, 1, 1, element_id]
731
max_coordinate = original_nodes[:, n_nodes_in, n_nodes_in, element_id]
732
element_length = max_coordinate - min_coordinate
733
734
# Test if the element is on the slice axis. If not just continue with the next element.
735
if !((min_coordinate[slice_dimension] <= point[slice_dimension] &&
736
max_coordinate[slice_dimension] > point[slice_dimension]) ||
737
(point[slice_dimension] == upper_limit &&
738
max_coordinate[slice_dimension] == upper_limit))
739
continue
740
end
741
742
new_id += 1
743
744
# Construct vandermonde matrix for interpolation of each 2D element to a 1D element.
745
normalized_intercept = (point[slice_dimension] -
746
min_coordinate[slice_dimension]) /
747
element_length[1] * 2 - 1
748
vandermonde = polynomial_interpolation_matrix(nodes_in, normalized_intercept)
749
750
# Interpolate to each node of new 1D element.
751
for v in 1:n_variables
752
for node in 1:n_nodes_in
753
new_unstructured_data[node, new_id, v] = (vandermonde * unstructured_data[node,
754
:,
755
element_id,
756
v])[1]
757
end
758
end
759
760
new_nodes[:, new_id] = original_nodes[other_dimension, :, 1, element_id]
761
end
762
763
return get_data_1d(reshape(new_nodes[:, 1:new_id], 1, n_nodes_in, new_id),
764
new_unstructured_data[:, 1:new_id, :], nvisnodes, reinterpolate)
765
end
766
767
# Calculate the arc length of a curve given by ndims x npoints point coordinates (piece-wise linear approximation)
768
function calc_arc_length(coordinates)
769
n_points = size(coordinates, 2)
770
arc_length = zeros(n_points)
771
for i in 1:(n_points - 1)
772
dist_squared = zero(eltype(arc_length))
773
for j in axes(coordinates, 1)
774
dist_squared += (coordinates[j, i + 1] - coordinates[j, i])^2
775
end
776
arc_length[i + 1] = arc_length[i] + sqrt(dist_squared)
777
end
778
return arc_length
779
end
780
781
# Convert 2d unstructured data to 1d data at given curve for the TreeMesh.
782
function unstructured_2d_to_1d_curve(original_nodes, unstructured_data, nvisnodes,
783
curve, mesh::TreeMesh, solver, cache)
784
n_points_curve = size(curve)[2]
785
n_nodes, _, n_elements, n_variables = size(unstructured_data)
786
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)
787
baryweights_in = barycentric_weights(nodes_in)
788
789
# Utility function to extract points as `SVector`s
790
get_point(data, idx...) = SVector(data[1, idx...], data[2, idx...])
791
792
# Check if input is correct.
793
min = get_point(original_nodes, 1, 1, 1)
794
max = get_point(original_nodes, n_nodes, n_nodes, n_elements)
795
@assert size(curve)==(2, size(curve)[2]) "Coordinates along curve must be 2xn dimensional."
796
for element in 1:n_points_curve
797
p = get_point(curve, element)
798
if any(p .< min) || any(p .> max)
799
throw(DomainError("Some coordinates from `curve` are outside of the domain."))
800
end
801
end
802
803
# Set nodes according to the length of the curve.
804
arc_length = calc_arc_length(curve)
805
806
# Setup data structures.
807
data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)
808
temp_data = Array{Float64}(undef, n_nodes, n_points_curve, n_variables)
809
810
# For each coordinate find the corresponding element with its id.
811
element_ids = get_elements_by_coordinates(curve, mesh, solver, cache)
812
813
# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,
814
# row vectors. We allocate memory here to improve performance.
815
vandermonde_x = polynomial_interpolation_matrix(nodes_in, zero(eltype(curve)))
816
vandermonde_y = similar(vandermonde_x)
817
818
# Iterate over all found elements.
819
for element in 1:n_points_curve
820
element_id = element_ids[element]
821
min_coordinate = get_point(original_nodes, 1, 1,
822
element_id)
823
max_coordinate = get_point(original_nodes, n_nodes, n_nodes,
824
element_id)
825
element_length = max_coordinate - min_coordinate
826
827
normalized_coordinates = (get_point(curve, element) - min_coordinate) /
828
element_length[1] * 2 .- 1
829
830
# Interpolate to a single point in each element.
831
# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,
832
# row vectors.
833
polynomial_interpolation_matrix!(vandermonde_x, nodes_in,
834
normalized_coordinates[1], baryweights_in)
835
polynomial_interpolation_matrix!(vandermonde_y, nodes_in,
836
normalized_coordinates[2], baryweights_in)
837
for v in 1:n_variables
838
for i in 1:n_nodes
839
res_i = zero(eltype(temp_data))
840
for n in 1:n_nodes
841
res_i += vandermonde_y[n] * unstructured_data[i, n, element_id, v]
842
end
843
temp_data[i, element, v] = res_i
844
end
845
res_v = zero(eltype(data_on_curve))
846
for n in 1:n_nodes
847
res_v += vandermonde_x[n] * temp_data[n, element, v]
848
end
849
data_on_curve[element, v] = res_v
850
end
851
end
852
853
return arc_length, data_on_curve, nothing
854
end
855
856
# Convert 2d unstructured data from a general mesh to 1d data at given curve.
857
#
858
# We need to loop through all the points and check in which element they are
859
# located. A general implementation working for all mesh types has to perform
860
# a naive loop through all nodes. Thus, we use this entry point for dispatching
861
# on the `mesh` type.
862
function unstructured_2d_to_1d_curve(u, mesh, equations,
863
solver, cache,
864
curve, slice,
865
point, nvisnodes,
866
solution_variables)
867
return unstructured_2d_to_1d_curve_general(u, mesh, equations,
868
solver, cache,
869
curve, slice,
870
point, nvisnodes,
871
solution_variables)
872
end
873
874
function unstructured_2d_to_1d_curve_general(u, mesh, equations,
875
solver, cache,
876
input_curve, slice,
877
point, nvisnodes,
878
solution_variables)
879
# Create a 'PlotData2DTriangulated' object so a triangulation
880
# can be used when extracting relevant data.
881
pd = PlotData2DTriangulated(u, mesh, equations, solver, cache;
882
solution_variables, nvisnodes)
883
884
# If no curve is defined, create an axis curve.
885
if input_curve === nothing
886
input_curve = axis_curve(pd.x, pd.y, nothing, slice, point, nvisnodes)
887
end
888
889
@assert size(input_curve, 1)==2 "Input 'curve' must be 2xn dimensional."
890
891
# For each coordinate find the corresponding triangle with its ids.
892
# The default value if no element is found is 0.
893
ids_by_coordinates = get_ids_by_coordinates(input_curve, pd)
894
found_coordinates = view(ids_by_coordinates, :, 1) .!= 0
895
896
@assert !iszero(found_coordinates) "No points of 'curve' are inside of the solutions domain."
897
898
# These hold the ids of the elements and triangles the points of the curve sit in.
899
element_ids = @view ids_by_coordinates[found_coordinates, 1]
900
triangle_ids = @view ids_by_coordinates[found_coordinates, 2]
901
902
# Shorten the curve, so that it contains only point that were found.
903
curve = @view input_curve[:, found_coordinates]
904
905
n_variables = length(pd.data[1, 1])
906
n_points_curve = size(curve, 2)
907
908
# Set nodes according to the length of the curve.
909
arc_length = calc_arc_length(curve)
910
911
# Setup data structures.
912
data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)
913
914
# Iterate over all points on the curve.
915
for point in 1:n_points_curve
916
point_on_curve = SVector(curve[1, point], curve[2, point])
917
918
element = element_ids[point]
919
triangle_id = triangle_ids[point]
920
triangle = (pd.t[triangle_id, 1], pd.t[triangle_id, 2], pd.t[triangle_id, 3])
921
922
# Get the x and y coordinates of the corners of given triangle.
923
x_coordinates_triangle = SVector(pd.x[triangle[1], element],
924
pd.x[triangle[2], element],
925
pd.x[triangle[3], element])
926
y_coordinates_triangle = SVector(pd.y[triangle[1], element],
927
pd.y[triangle[2], element],
928
pd.y[triangle[3], element])
929
930
# Extract solution values in corners of the triangle.
931
data_in_triangle = (pd.data[triangle[1], element],
932
pd.data[triangle[2], element],
933
pd.data[triangle[3], element])
934
935
for v in 1:n_variables
936
# Extract solution values of variable `v` in corners of the triangle.
937
values_triangle = SVector(data_in_triangle[1][v],
938
data_in_triangle[2][v],
939
data_in_triangle[3][v])
940
941
# Linear interpolation in each triangle to the points on the curve.
942
data_on_curve[point, v] = triangle_interpolation(x_coordinates_triangle,
943
y_coordinates_triangle,
944
values_triangle,
945
point_on_curve)
946
end
947
end
948
949
return arc_length, data_on_curve, nothing
950
end
951
952
# Convert 3d unstructured data to 1d data at given curve.
953
function unstructured_3d_to_1d_curve(original_nodes, unstructured_data, nvisnodes,
954
curve, mesh::TreeMesh, solver, cache)
955
n_points_curve = size(curve)[2]
956
n_nodes, _, _, n_elements, n_variables = size(unstructured_data)
957
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes)
958
baryweights_in = barycentric_weights(nodes_in)
959
960
# Utility function to extract points as `SVector`s
961
get_point(data, idx...) = SVector(data[1, idx...],
962
data[2, idx...],
963
data[3, idx...])
964
965
# Check if input is correct.
966
min = get_point(original_nodes, 1, 1, 1, 1)
967
max = get_point(original_nodes, n_nodes, n_nodes, n_nodes, n_elements)
968
@assert size(curve)==(3, n_points_curve) "Coordinates along curve must be 3xn dimensional."
969
for element in 1:n_points_curve
970
p = get_point(curve, element)
971
if any(p .< min) || any(p .> max)
972
throw(DomainError("Some coordinates from `curve` are outside of the domain."))
973
end
974
end
975
976
# Set nodes according to the length of the curve.
977
arc_length = calc_arc_length(curve)
978
979
# Setup data structures.
980
data_on_curve = Array{Float64}(undef, n_points_curve, n_variables)
981
temp_data = Array{Float64}(undef, n_nodes, n_nodes + 1, n_points_curve, n_variables)
982
983
# For each coordinate find the corresponding element with its id.
984
element_ids = get_elements_by_coordinates(curve, mesh, solver, cache)
985
986
# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,
987
# row vectors. We allocate memory here to improve performance.
988
vandermonde_x = polynomial_interpolation_matrix(nodes_in, zero(eltype(curve)))
989
vandermonde_y = similar(vandermonde_x)
990
vandermonde_z = similar(vandermonde_x)
991
992
# Iterate over all found elements.
993
for element in 1:n_points_curve
994
element_id = element_ids[element]
995
min_coordinate = get_point(original_nodes, 1, 1, 1,
996
element_id)
997
max_coordinate = get_point(original_nodes, n_nodes, n_nodes, n_nodes,
998
element_id)
999
element_length = max_coordinate - min_coordinate
1000
1001
normalized_coordinates = (get_point(curve, element) - min_coordinate) /
1002
element_length[1] * 2 .- 1
1003
1004
# Interpolate to a single point in each element.
1005
# These Vandermonde matrices are really 1×n_nodes matrices, i.e.,
1006
# row vectors.
1007
polynomial_interpolation_matrix!(vandermonde_x, nodes_in,
1008
normalized_coordinates[1], baryweights_in)
1009
polynomial_interpolation_matrix!(vandermonde_y, nodes_in,
1010
normalized_coordinates[2], baryweights_in)
1011
polynomial_interpolation_matrix!(vandermonde_z, nodes_in,
1012
normalized_coordinates[3], baryweights_in)
1013
for v in 1:n_variables
1014
for i in 1:n_nodes
1015
for ii in 1:n_nodes
1016
res_ii = zero(eltype(temp_data))
1017
for n in 1:n_nodes
1018
res_ii += vandermonde_z[n] *
1019
unstructured_data[i, ii, n, element_id, v]
1020
end
1021
temp_data[i, ii, element, v] = res_ii
1022
end
1023
res_i = zero(eltype(temp_data))
1024
for n in 1:n_nodes
1025
res_i += vandermonde_y[n] * temp_data[i, n, element, v]
1026
end
1027
temp_data[i, n_nodes + 1, element, v] = res_i
1028
end
1029
res_v = zero(eltype(temp_data))
1030
for n in 1:n_nodes
1031
res_v += vandermonde_x[n] * temp_data[n, n_nodes + 1, element, v]
1032
end
1033
data_on_curve[element, v] = res_v
1034
end
1035
end
1036
1037
return arc_length, data_on_curve, nothing
1038
end
1039
1040
# Convert 3d unstructured data from a general mesh to 1d data at given curve.
1041
#
1042
# We need to loop through all the points and check in which element they are
1043
# located. A general implementation working for all mesh types has to perform
1044
# a naive loop through all nodes. Thus, we use this entry point for dispatching
1045
# on the `mesh` type.
1046
function unstructured_3d_to_1d_curve(u, mesh, equations, solver, cache,
1047
curve, solution_variables)
1048
return unstructured_3d_to_1d_curve_general(u, equations, solver, cache,
1049
curve, solution_variables)
1050
end
1051
1052
function unstructured_3d_to_1d_curve_general(u, equations, solver, cache,
1053
curve, solution_variables)
1054
# Set up data structure.
1055
@assert size(curve, 1) == 3
1056
n_points_curve = size(curve, 2)
1057
1058
# Get the number of variables after applying the transformation to solution variables.
1059
u_node = get_node_vars(u, equations, solver, 1, 1, 1, 1)
1060
var_node = solution_variables(u_node, equations)
1061
n_variables = length(var_node)
1062
data_on_curve = Array{eltype(var_node)}(undef, n_points_curve, n_variables)
1063
1064
nodes = cache.elements.node_coordinates
1065
interpolation_cache = get_value_at_point_3d_cache(u)
1066
1067
# Iterate over every point on the curve and determine the solutions value at given point.
1068
for i in 1:n_points_curve
1069
point = SVector(curve[1, i], curve[2, i], curve[3, i])
1070
get_value_at_point_3d!(view(data_on_curve, i, :), point, solution_variables,
1071
nodes, u, equations, solver;
1072
cache = interpolation_cache)
1073
end
1074
1075
mesh_vertices_x = nothing
1076
1077
return calc_arc_length(curve), data_on_curve, mesh_vertices_x
1078
end
1079
1080
# Check if the first 'amount'-many points can still form a valid tetrahedron.
1081
function is_valid_tetrahedron(amount, coordinates; tol = 10^-4)
1082
a = SVector(coordinates[1, 1], coordinates[2, 1], coordinates[3, 1])
1083
b = SVector(coordinates[1, 2], coordinates[2, 2], coordinates[3, 2])
1084
c = SVector(coordinates[1, 3], coordinates[2, 3], coordinates[3, 3])
1085
d = SVector(coordinates[1, 4], coordinates[2, 4], coordinates[3, 4])
1086
1087
if amount == 2 # If two points are the same, then no tetrahedron can be formed.
1088
return !(isapprox(a, b; atol = tol))
1089
elseif amount == 3 # Check if three points are on the same line.
1090
return !on_the_same_line(a, b, c; tol = tol)
1091
elseif amount == 4 # Check if four points form a tetrahedron.
1092
# This is the same as
1093
# A = hcat(coordinates[1, :], coordinates[2, :], coordinates[3, :],
1094
# SVector(1, 1, 1, 1))
1095
# but more efficient.
1096
A = SMatrix{4, 4}(coordinates[1, 1], coordinates[2, 1], coordinates[3, 1], 1,
1097
coordinates[1, 2], coordinates[2, 2], coordinates[3, 2], 1,
1098
coordinates[1, 3], coordinates[2, 3], coordinates[3, 3], 1,
1099
coordinates[1, 4], coordinates[2, 4], coordinates[3, 4], 1)
1100
return !isapprox(det(A), 0; atol = tol)
1101
else # With one point a tetrahedron can always be formed.
1102
return true
1103
end
1104
end
1105
1106
# Check if three given 3D-points are on the same line.
1107
function on_the_same_line(a, b, c; tol = 10^-4)
1108
# Calculate the intersection of the a-b-axis at x=0.
1109
if b[1] == 0
1110
intersect_a_b = b
1111
else
1112
intersect_a_b = a - b .* (a[1] / b[1])
1113
end
1114
# Calculate the intersection of the a-c-axis at x=0.
1115
if c[1] == 0
1116
intersect_a_c = c
1117
else
1118
intersect_a_c = a - c .* (a[1] / c[1])
1119
end
1120
return isapprox(intersect_a_b, intersect_a_c; atol = tol)
1121
end
1122
1123
# Interpolate from four corners of a tetrahedron to a single point.
1124
function tetrahedron_interpolation(x_coordinates_in, y_coordinates_in, z_coordinates_in,
1125
values_in, coordinate_out)
1126
A = hcat(x_coordinates_in, y_coordinates_in, z_coordinates_in, SVector(1, 1, 1, 1))
1127
c = A \ values_in
1128
return c[1] * coordinate_out[1] + c[2] * coordinate_out[2] +
1129
c[3] * coordinate_out[3] + c[4]
1130
end
1131
1132
# Calculate the distances from every entry in `nodes` to the given `point` and return
1133
# the minimal/maximal squared distance as well as the index of the minimal squared distance.
1134
function squared_distances_from_single_point!(distances, nodes, point)
1135
_, n_nodes, _, _, n_elements = size(nodes)
1136
@assert size(nodes, 1) == 3
1137
@assert size(nodes, 3) == size(nodes, 4) == n_nodes
1138
@assert size(distances, 1) == size(distances, 2) == size(distances, 3) == n_nodes
1139
@assert size(distances, 4) == n_elements
1140
1141
# Prepare for reduction
1142
dist_max = typemin(eltype(distances))
1143
dist_min = typemax(eltype(distances))
1144
index_min = CartesianIndex(1, 1, 1, 1)
1145
1146
# Iterate over every entry
1147
for element in 1:n_elements
1148
for z in 1:n_nodes, y in 1:n_nodes, x in 1:n_nodes
1149
dist = (nodes[1, x, y, z, element] - point[1])^2 +
1150
(nodes[2, x, y, z, element] - point[2])^2 +
1151
(nodes[3, x, y, z, element] - point[3])^2
1152
distances[x, y, z, element] = dist
1153
1154
dist_max = max(dist_max, dist)
1155
if dist < dist_min
1156
dist_min = dist
1157
index_min = CartesianIndex(x, y, z, element)
1158
end
1159
end
1160
end
1161
1162
return dist_max, dist_min, index_min
1163
end
1164
1165
# Interpolate the data on given nodes to a single value at given point.
1166
function get_value_at_point_3d_cache(u)
1167
n_variables, n_x_nodes, n_y_nodes, n_z_nodes, n_elements = size(u)
1168
@assert n_x_nodes == n_y_nodes == n_z_nodes
1169
n_nodes = n_x_nodes
1170
1171
distances = zeros(n_nodes, n_nodes, n_nodes, n_elements)
1172
coordinates_tetrahedron = Array{Float64, 2}(undef, 3, 4)
1173
value_tetrahedron = Array{eltype(u)}(undef, n_variables, 4)
1174
1175
cache = (; distances, coordinates_tetrahedron, value_tetrahedron)
1176
return cache
1177
end
1178
1179
function get_value_at_point_3d!(data_on_curve_at_point, point, solution_variables,
1180
nodes, u, equations, solver;
1181
cache = get_value_at_point_3d_cache(u))
1182
# Set up data structures.
1183
n_variables = size(u, 1)
1184
(; distances, coordinates_tetrahedron, value_tetrahedron) = cache
1185
1186
maximum_distance, _, index = squared_distances_from_single_point!(distances, nodes,
1187
point)
1188
# We could also use the following code to find the maximum distance and index:
1189
# maximum_distance = maximum(distances)
1190
# index = argmin(distances)
1191
# However, it is more efficient if we do it in one pass through the memory.
1192
1193
# If the point sits exactly on a node, no interpolation is needed.
1194
nodes_at_index = SVector(nodes[1, index[1], index[2], index[3], index[4]],
1195
nodes[2, index[1], index[2], index[3], index[4]],
1196
nodes[3, index[1], index[2], index[3], index[4]])
1197
if nodes_at_index == point
1198
u_node = get_node_vars(u, equations, solver,
1199
index[1], index[2], index[3], index[4])
1200
data_on_curve_at_point .= solution_variables(u_node, equations)
1201
return data_on_curve_at_point
1202
end
1203
1204
# Restrict the interpolation to the closest element only.
1205
closest_element = index[4]
1206
@views element_distances = distances[:, :, :, closest_element]
1207
1208
# Find a tetrahedron, which is given by four corners, to interpolate from.
1209
for i in 1:4
1210
# Iterate until a valid tetrahedron is found.
1211
while true
1212
index = argmin(element_distances)
1213
element_distances[index[1], index[2], index[3]] = maximum_distance
1214
1215
for k in 1:3
1216
coordinates_tetrahedron[k, i] = nodes[k,
1217
index[1], index[2], index[3],
1218
closest_element]
1219
end
1220
for v in 1:n_variables
1221
value_tetrahedron[v, i] = u[v,
1222
index[1], index[2], index[3],
1223
closest_element]
1224
end
1225
1226
# Look for another point if current tetrahedron is not valid.
1227
if is_valid_tetrahedron(i, coordinates_tetrahedron)
1228
break
1229
end
1230
end
1231
end
1232
1233
# Interpolate from tetrahedron to given point.
1234
x_coordinates = SVector(coordinates_tetrahedron[1, 1],
1235
coordinates_tetrahedron[1, 2],
1236
coordinates_tetrahedron[1, 3],
1237
coordinates_tetrahedron[1, 4])
1238
y_coordinates = SVector(coordinates_tetrahedron[2, 1],
1239
coordinates_tetrahedron[2, 2],
1240
coordinates_tetrahedron[2, 3],
1241
coordinates_tetrahedron[2, 4])
1242
z_coordinates = SVector(coordinates_tetrahedron[3, 1],
1243
coordinates_tetrahedron[3, 2],
1244
coordinates_tetrahedron[3, 3],
1245
coordinates_tetrahedron[3, 4])
1246
# We compute the solution_variables first and interpolate them.
1247
u1 = get_node_vars(value_tetrahedron, equations, solver, 1)
1248
u2 = get_node_vars(value_tetrahedron, equations, solver, 2)
1249
u3 = get_node_vars(value_tetrahedron, equations, solver, 3)
1250
u4 = get_node_vars(value_tetrahedron, equations, solver, 4)
1251
val1 = solution_variables(u1, equations)
1252
val2 = solution_variables(u2, equations)
1253
val3 = solution_variables(u3, equations)
1254
val4 = solution_variables(u4, equations)
1255
for v in eachindex(val1)
1256
values = SVector(val1[v],
1257
val2[v],
1258
val3[v],
1259
val4[v])
1260
data_on_curve_at_point[v] = tetrahedron_interpolation(x_coordinates,
1261
y_coordinates,
1262
z_coordinates,
1263
values, point)
1264
end
1265
1266
return data_on_curve_at_point
1267
end
1268
1269
# Convert 3d unstructured data to 1d slice and interpolate them.
1270
function unstructured_3d_to_1d(original_nodes, unstructured_data, nvisnodes,
1271
reinterpolate, slice, point)
1272
if slice === :x
1273
slice_dimension = 1
1274
other_dimensions = [2, 3]
1275
elseif slice === :y
1276
slice_dimension = 2
1277
other_dimensions = [1, 3]
1278
elseif slice === :z
1279
slice_dimension = 3
1280
other_dimensions = [1, 2]
1281
else
1282
error("illegal dimension '$slice', supported dimensions are :x, :y and :z")
1283
end
1284
1285
# Set up data structures to store new 1D data.
1286
@views new_unstructured_data = similar(unstructured_data[1, 1, ..])
1287
@views temp_unstructured_data = similar(unstructured_data[1, ..])
1288
@views new_nodes = similar(original_nodes[1, 1, 1, ..])
1289
1290
n_nodes_in, _, _, n_elements, n_variables = size(unstructured_data)
1291
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)
1292
1293
# Test if point lies in the domain.
1294
lower_limit = original_nodes[1, 1, 1, 1, 1]
1295
upper_limit = original_nodes[1, n_nodes_in, n_nodes_in, n_nodes_in, n_elements]
1296
1297
@assert length(point)>=3 "Point must be three-dimensional."
1298
if prod(point[other_dimensions] .< lower_limit) ||
1299
prod(point[other_dimensions] .> upper_limit)
1300
error(string("Slice axis is outside of domain. ",
1301
" point[$other_dimensions]=$(point[other_dimensions]) must be between $lower_limit and $upper_limit"))
1302
end
1303
1304
# Count the amount of new elements.
1305
new_id = 0
1306
1307
# Permute dimensions so that the slice dimensions are always the in correct places for later use.
1308
if slice === :x
1309
original_nodes = permutedims(original_nodes, [1, 3, 4, 2, 5])
1310
unstructured_data = permutedims(unstructured_data, [2, 3, 1, 4, 5])
1311
elseif slice === :y
1312
original_nodes = permutedims(original_nodes, [1, 2, 4, 3, 5])
1313
unstructured_data = permutedims(unstructured_data, [1, 3, 2, 4, 5])
1314
end
1315
1316
# Iterate over all elements to find the ones that lie on the slice axis.
1317
for element_id in 1:n_elements
1318
min_coordinate = original_nodes[:, 1, 1, 1, element_id]
1319
max_coordinate = original_nodes[:, n_nodes_in, n_nodes_in, n_nodes_in,
1320
element_id]
1321
element_length = max_coordinate - min_coordinate
1322
1323
# Test if the element is on the slice axis. If not just continue with the next element.
1324
if !((prod(min_coordinate[other_dimensions] .<= point[other_dimensions]) &&
1325
prod(max_coordinate[other_dimensions] .> point[other_dimensions])) ||
1326
(point[other_dimensions] == upper_limit &&
1327
prod(max_coordinate[other_dimensions] .== upper_limit)))
1328
continue
1329
end
1330
1331
new_id += 1
1332
1333
# Construct vandermonde matrix for interpolation of each 2D element to a 1D element.
1334
normalized_intercept = (point[other_dimensions] .-
1335
min_coordinate[other_dimensions]) /
1336
element_length[1] * 2 .- 1
1337
vandermonde_i = polynomial_interpolation_matrix(nodes_in,
1338
normalized_intercept[1])
1339
vandermonde_ii = polynomial_interpolation_matrix(nodes_in,
1340
normalized_intercept[2])
1341
1342
# Interpolate to each node of new 1D element.
1343
for v in 1:n_variables
1344
for i in 1:n_nodes_in
1345
for ii in 1:n_nodes_in
1346
temp_unstructured_data[i, ii, new_id, v] = (vandermonde_ii * unstructured_data[ii,
1347
:,
1348
i,
1349
element_id,
1350
v])[1]
1351
end
1352
new_unstructured_data[i, new_id, v] = (vandermonde_i * temp_unstructured_data[i,
1353
:,
1354
new_id,
1355
v])[1]
1356
end
1357
end
1358
1359
new_nodes[:, new_id] = original_nodes[slice_dimension, 1, 1, :, element_id]
1360
end
1361
1362
return get_data_1d(reshape(new_nodes[:, 1:new_id], 1, n_nodes_in, new_id),
1363
new_unstructured_data[:, 1:new_id, :], nvisnodes, reinterpolate)
1364
end
1365
1366
# Interpolate unstructured DG data to structured data (cell-centered)
1367
#
1368
# This function takes DG data in an unstructured, Cartesian layout and converts it to a uniformly
1369
# distributed 2D layout.
1370
#
1371
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
1372
# thus be changed in future releases.
1373
function unstructured2structured(unstructured_data, normalized_coordinates,
1374
levels, resolution, nvisnodes_per_level)
1375
# Extract data shape information
1376
n_nodes_in, _, n_elements, n_variables = size(unstructured_data)
1377
1378
# Get node coordinates for DG locations on reference element
1379
nodes_in, _ = gauss_lobatto_nodes_weights(n_nodes_in)
1380
1381
# Calculate interpolation vandermonde matrices for each level
1382
max_level = length(nvisnodes_per_level) - 1
1383
vandermonde_per_level = []
1384
for l in 0:max_level
1385
n_nodes_out = nvisnodes_per_level[l + 1]
1386
dx = 2 / n_nodes_out
1387
nodes_out = collect(range(-1 + dx / 2, 1 - dx / 2, length = n_nodes_out))
1388
push!(vandermonde_per_level,
1389
polynomial_interpolation_matrix(nodes_in, nodes_out))
1390
end
1391
1392
# For each element, calculate index position at which to insert data in global data structure
1393
lower_left_index = element2index(normalized_coordinates, levels, resolution,
1394
nvisnodes_per_level)
1395
1396
# Create output data structure
1397
structured = [Matrix{Float64}(undef, resolution, resolution) for _ in 1:n_variables]
1398
1399
# For each variable, interpolate element data and store to global data structure
1400
for v in 1:n_variables
1401
# Reshape data array for use in multiply_dimensionwise function
1402
reshaped_data = reshape(unstructured_data[:, :, :, v], 1, n_nodes_in,
1403
n_nodes_in, n_elements)
1404
1405
for element_id in 1:n_elements
1406
# Extract level for convenience
1407
level = levels[element_id]
1408
1409
# Determine target indices
1410
n_nodes_out = nvisnodes_per_level[level + 1]
1411
first = lower_left_index[:, element_id]
1412
last = first .+ (n_nodes_out - 1)
1413
1414
# Interpolate data
1415
vandermonde = vandermonde_per_level[level + 1]
1416
structured[v][first[1]:last[1], first[2]:last[2]] .= (reshape(multiply_dimensionwise(vandermonde,
1417
reshaped_data[:,
1418
:,
1419
:,
1420
element_id]),
1421
n_nodes_out,
1422
n_nodes_out))
1423
end
1424
end
1425
1426
return structured
1427
end
1428
1429
# For a given normalized element coordinate, return the index of its lower left
1430
# contribution to the global data structure
1431
#
1432
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
1433
# thus be changed in future releases.
1434
function element2index(normalized_coordinates, levels, resolution, nvisnodes_per_level)
1435
@assert size(normalized_coordinates, 1)==2 "only works in 2D"
1436
1437
n_elements = length(levels)
1438
1439
# First, determine lower left coordinate for all cells
1440
dx = 2 / resolution
1441
ndim = 2
1442
lower_left_coordinate = Array{Float64}(undef, ndim, n_elements)
1443
for element_id in 1:n_elements
1444
nvisnodes = nvisnodes_per_level[levels[element_id] + 1]
1445
lower_left_coordinate[1, element_id] = (normalized_coordinates[1, element_id] -
1446
(nvisnodes - 1) / 2 * dx)
1447
lower_left_coordinate[2, element_id] = (normalized_coordinates[2, element_id] -
1448
(nvisnodes - 1) / 2 * dx)
1449
end
1450
1451
# Then, convert coordinate to global index
1452
indices = coordinate2index(lower_left_coordinate, resolution)
1453
1454
return indices
1455
end
1456
1457
# Find 2D array index for a 2-tuple of normalized, cell-centered coordinates (i.e., in [-1,1])
1458
#
1459
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
1460
# thus be changed in future releases.
1461
function coordinate2index(coordinate, resolution::Integer)
1462
# Calculate 1D normalized coordinates
1463
dx = 2 / resolution
1464
mesh_coordinates = collect(range(-1 + dx / 2, 1 - dx / 2, length = resolution))
1465
1466
# Find index
1467
id_x = searchsortedfirst.(Ref(mesh_coordinates), coordinate[1, :],
1468
lt = (x, y) -> x .< y .- dx / 2)
1469
id_y = searchsortedfirst.(Ref(mesh_coordinates), coordinate[2, :],
1470
lt = (x, y) -> x .< y .- dx / 2)
1471
return transpose(hcat(id_x, id_y))
1472
end
1473
1474
# Calculate the vertices for each mesh cell such that it can be visualized as a closed box
1475
#
1476
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
1477
# thus be changed in future releases.
1478
function calc_vertices(coordinates, levels, length_level_0)
1479
ndim = size(coordinates, 1)
1480
@assert ndim==2 "only works in 2D"
1481
1482
# Initialize output arrays
1483
n_elements = length(levels)
1484
n_points_per_element = 2^ndim + 2
1485
x = Vector{Float64}(undef, n_points_per_element * n_elements)
1486
y = Vector{Float64}(undef, n_points_per_element * n_elements)
1487
1488
# Calculate vertices for all coordinates at once
1489
for element_id in 1:n_elements
1490
length = length_level_0 / 2^levels[element_id]
1491
index = n_points_per_element * (element_id - 1)
1492
x[index + 1] = coordinates[1, element_id] - 1 / 2 * length
1493
x[index + 2] = coordinates[1, element_id] + 1 / 2 * length
1494
x[index + 3] = coordinates[1, element_id] + 1 / 2 * length
1495
x[index + 4] = coordinates[1, element_id] - 1 / 2 * length
1496
x[index + 5] = coordinates[1, element_id] - 1 / 2 * length
1497
x[index + 6] = NaN
1498
1499
y[index + 1] = coordinates[2, element_id] - 1 / 2 * length
1500
y[index + 2] = coordinates[2, element_id] - 1 / 2 * length
1501
y[index + 3] = coordinates[2, element_id] + 1 / 2 * length
1502
y[index + 4] = coordinates[2, element_id] + 1 / 2 * length
1503
y[index + 5] = coordinates[2, element_id] - 1 / 2 * length
1504
y[index + 6] = NaN
1505
end
1506
1507
return x, y
1508
end
1509
1510
# Calculate the vertices to plot each grid line for StructuredMesh
1511
#
1512
# Note: This is a low-level function that is not considered as part of Trixi.jl's interface and may
1513
# thus be changed in future releases.
1514
function calc_vertices(node_coordinates, mesh)
1515
@unpack cells_per_dimension = mesh
1516
@assert size(node_coordinates, 1)==2 "only works in 2D"
1517
1518
linear_indices = LinearIndices(size(mesh))
1519
1520
# Initialize output arrays
1521
n_lines = sum(cells_per_dimension) + 2
1522
max_length = maximum(cells_per_dimension)
1523
n_nodes = size(node_coordinates, 2)
1524
1525
# Create output as two matrices `x` and `y`, each holding the node locations for each of the `n_lines` grid lines
1526
# The # of rows in the matrices must be sufficient to store the longest dimension (`max_length`),
1527
# and for each the node locations without doubling the corner nodes (`n_nodes-1`), plus the final node (`+1`)
1528
# Rely on Plots.jl to ignore `NaN`s (i.e., they are not plotted) to handle shorter lines
1529
x = fill(NaN, max_length * (n_nodes - 1) + 1, n_lines)
1530
y = fill(NaN, max_length * (n_nodes - 1) + 1, n_lines)
1531
1532
line_index = 1
1533
# Lines in x-direction
1534
# Bottom boundary
1535
i = 1
1536
for cell_x in axes(mesh, 1)
1537
for node in 1:(n_nodes - 1)
1538
x[i, line_index] = node_coordinates[1, node, 1, linear_indices[cell_x, 1]]
1539
y[i, line_index] = node_coordinates[2, node, 1, linear_indices[cell_x, 1]]
1540
1541
i += 1
1542
end
1543
end
1544
# Last point on bottom boundary
1545
x[i, line_index] = node_coordinates[1, end, 1, linear_indices[end, 1]]
1546
y[i, line_index] = node_coordinates[2, end, 1, linear_indices[end, 1]]
1547
1548
# Other lines in x-direction
1549
line_index += 1
1550
for cell_y in axes(mesh, 2)
1551
i = 1
1552
for cell_x in axes(mesh, 1)
1553
for node in 1:(n_nodes - 1)
1554
x[i, line_index] = node_coordinates[1, node, end,
1555
linear_indices[cell_x, cell_y]]
1556
y[i, line_index] = node_coordinates[2, node, end,
1557
linear_indices[cell_x, cell_y]]
1558
1559
i += 1
1560
end
1561
end
1562
# Last point on line
1563
x[i, line_index] = node_coordinates[1, end, end, linear_indices[end, cell_y]]
1564
y[i, line_index] = node_coordinates[2, end, end, linear_indices[end, cell_y]]
1565
1566
line_index += 1
1567
end
1568
1569
# Lines in y-direction
1570
# Left boundary
1571
i = 1
1572
for cell_y in axes(mesh, 2)
1573
for node in 1:(n_nodes - 1)
1574
x[i, line_index] = node_coordinates[1, 1, node, linear_indices[1, cell_y]]
1575
y[i, line_index] = node_coordinates[2, 1, node, linear_indices[1, cell_y]]
1576
1577
i += 1
1578
end
1579
end
1580
# Last point on left boundary
1581
x[i, line_index] = node_coordinates[1, 1, end, linear_indices[1, end]]
1582
y[i, line_index] = node_coordinates[2, 1, end, linear_indices[1, end]]
1583
1584
# Other lines in y-direction
1585
line_index += 1
1586
for cell_x in axes(mesh, 1)
1587
i = 1
1588
for cell_y in axes(mesh, 2)
1589
for node in 1:(n_nodes - 1)
1590
x[i, line_index] = node_coordinates[1, end, node,
1591
linear_indices[cell_x, cell_y]]
1592
y[i, line_index] = node_coordinates[2, end, node,
1593
linear_indices[cell_x, cell_y]]
1594
1595
i += 1
1596
end
1597
end
1598
# Last point on line
1599
x[i, line_index] = node_coordinates[1, end, end, linear_indices[cell_x, end]]
1600
y[i, line_index] = node_coordinates[2, end, end, linear_indices[cell_x, end]]
1601
1602
line_index += 1
1603
end
1604
1605
return x, y
1606
end
1607
1608
# Convert `slice` to orientations (1 -> `x`, 2 -> `y`, 3 -> `z`) for the two axes in a 2D plot
1609
function _get_orientations(mesh, slice)
1610
if ndims(mesh) == 2 || (ndims(mesh) == 3 && slice === :xy)
1611
orientation_x = 1
1612
orientation_y = 2
1613
elseif ndims(mesh) == 3 && slice === :xz
1614
orientation_x = 1
1615
orientation_y = 3
1616
elseif ndims(mesh) == 3 && slice === :yz
1617
orientation_x = 2
1618
orientation_y = 3
1619
else
1620
orientation_x = 0
1621
orientation_y = 0
1622
end
1623
return orientation_x, orientation_y
1624
end
1625
1626
# Convert `orientation` into a guide label (see also `_get_orientations`)
1627
function _get_guide(orientation::Integer)
1628
if orientation == 1
1629
return "\$x\$"
1630
elseif orientation == 2
1631
return "\$y\$"
1632
elseif orientation == 3
1633
return "\$z\$"
1634
else
1635
return ""
1636
end
1637
end
1638
1639
# plotting_interpolation_matrix(dg; kwargs...)
1640
#
1641
# Interpolation matrix which maps discretization nodes to a set of plotting nodes.
1642
# Defaults to the identity matrix of size `length(solver.basis.nodes)`, and interpolates
1643
# to equispaced nodes for DGSEM (set by kwarg `nvisnodes` in the plotting function).
1644
#
1645
# Example:
1646
# ```julia
1647
# A = plotting_interpolation_matrix(dg)
1648
# A * dg.basis.nodes # => vector of nodes at which to plot the solution
1649
# ```
1650
#
1651
# Note: we cannot use UniformScaling to define the interpolation matrix since we use it with `kron`
1652
# to define a multi-dimensional interpolation matrix later.
1653
plotting_interpolation_matrix(dg; kwargs...) = I(length(dg.basis.nodes))
1654
1655
function face_plotting_interpolation_matrix(dg::DGSEM;
1656
nvisnodes = 2 * length(dg.basis.nodes))
1657
return polynomial_interpolation_matrix(dg.basis.nodes, LinRange(-1, 1, nvisnodes))
1658
end
1659
1660
function plotting_interpolation_matrix(dg::DGSEM;
1661
nvisnodes = 2 * length(dg.basis.nodes))
1662
Vp1D = polynomial_interpolation_matrix(dg.basis.nodes, LinRange(-1, 1, nvisnodes))
1663
# For quadrilateral elements, interpolation to plotting nodes involves applying a 1D interpolation
1664
# operator to each line of nodes. This is equivalent to multiplying the vector containing all node
1665
# node coordinates on an element by a Kronecker product of the 1D interpolation operator (e.g., a
1666
# multi-dimensional interpolation operator).
1667
return kron(Vp1D, Vp1D)
1668
end
1669
1670
function reference_node_coordinates_2d(dg::DGSEM)
1671
@unpack nodes = dg.basis
1672
r = vec([nodes[i] for i in eachnode(dg), j in eachnode(dg)])
1673
s = vec([nodes[j] for i in eachnode(dg), j in eachnode(dg)])
1674
return r, s
1675
end
1676
1677
# Find element and triangle ids containing coordinates given as a matrix [ndims, npoints]
1678
function get_ids_by_coordinates!(ids, coordinates, pd)
1679
if length(ids) != 2 * size(coordinates, 2)
1680
throw(DimensionMismatch("storage length for element ids does not match the number of coordinates"))
1681
end
1682
1683
n_coordinates = size(coordinates, 2)
1684
1685
for index in 1:n_coordinates
1686
point = SVector(coordinates[1, index], coordinates[2, index])
1687
ids[index, :] .= find_element(point, pd)
1688
end
1689
1690
return ids
1691
end
1692
1693
# Find the ids of elements and triangles containing given coordinates by using the triangulation in 'pd'.
1694
function get_ids_by_coordinates(coordinates, pd)
1695
ids = Matrix{Int}(undef, size(coordinates, 2), 2)
1696
get_ids_by_coordinates!(ids, coordinates, pd)
1697
return ids
1698
end
1699
1700
# Check if given 'point' is inside the triangle with corners corresponding to the coordinates of x and y.
1701
function is_in_triangle(point, x, y)
1702
a = SVector(x[1], y[1])
1703
b = SVector(x[2], y[2])
1704
c = SVector(x[3], y[3])
1705
return is_on_same_side(point, a, b, c) && is_on_same_side(point, b, c, a) &&
1706
is_on_same_side(point, c, a, b)
1707
end
1708
1709
# Create an axis through x and y to then check if 'point' is on the same side of the axis as z.
1710
function is_on_same_side(point, x, y, z)
1711
if (y[1] - x[1]) == 0
1712
return (point[1] - x[1]) * (z[1] - x[1]) >= 0
1713
else
1714
a = (y[2] - x[2]) / (y[1] - x[1])
1715
b = x[2] - a * x[1]
1716
return (z[2] - a * z[1] - b) * (point[2] - a * point[1] - b) >= 0
1717
end
1718
end
1719
1720
# For a given 'point', return the id of the element it is contained in in; if not found return 0.
1721
function find_element(point, pd)
1722
n_tri = size(pd.t, 1)
1723
n_elements = size(pd.x, 2)
1724
1725
# Iterate over all elements.
1726
for element in 1:n_elements
1727
# Iterate over all triangles in given element.
1728
for tri in 1:n_tri
1729
# The code below is equivalent to
1730
# x == pd.x[pd.t[tri, :], element]
1731
# y == pd.y[pd.t[tri, :], element]
1732
# but avoids allocations and is thus more efficient.
1733
tri_indices = (pd.t[tri, 1], pd.t[tri, 2], pd.t[tri, 3])
1734
x = SVector(pd.x[tri_indices[1], element],
1735
pd.x[tri_indices[2], element],
1736
pd.x[tri_indices[3], element])
1737
y = SVector(pd.y[tri_indices[1], element],
1738
pd.y[tri_indices[2], element],
1739
pd.y[tri_indices[3], element])
1740
if is_in_triangle(point, x, y)
1741
return (element, tri)
1742
end
1743
end
1744
end
1745
1746
return (0, 0)
1747
end
1748
1749
# Interpolate from three corners of a triangle to a single point.
1750
function triangle_interpolation(x_coordinates_in, y_coordinates_in, values_in,
1751
coordinate_out)
1752
A = hcat(x_coordinates_in, y_coordinates_in, SVector(1, 1, 1))
1753
c = A \ values_in
1754
return c[1] * coordinate_out[1] + c[2] * coordinate_out[2] + c[3]
1755
end
1756
1757
# Create an axis.
1758
function axis_curve(nodes_x, nodes_y, nodes_z, slice, point, n_points)
1759
if n_points === nothing
1760
n_points = 64
1761
end
1762
dimensions = length(point)
1763
curve = zeros(dimensions, n_points)
1764
if slice == :x
1765
xmin, xmax = extrema(nodes_x)
1766
curve[1, :] .= range(xmin, xmax, length = n_points)
1767
curve[2, :] .= point[2]
1768
if dimensions === 3
1769
curve[3, :] .= point[3]
1770
end
1771
elseif slice == :y
1772
ymin, ymax = extrema(nodes_y)
1773
curve[1, :] .= point[1]
1774
curve[2, :] .= range(ymin, ymax, length = n_points)
1775
if dimensions === 3
1776
curve[3, :] .= point[3]
1777
end
1778
elseif slice == :z
1779
zmin, zmax = extrema(nodes_z)
1780
curve[1, :] .= point[1]
1781
curve[2, :] .= point[2]
1782
curve[3, :] .= range(zmin, zmax, length = n_points)
1783
else
1784
@assert false "Input for 'slice' is not supported here."
1785
end
1786
1787
return curve
1788
end
1789
end # @muladd
1790
1791