Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/ext/TrixiMakieExt.jl
2055 views
1
# Package extension for adding Makie-based features to Trixi.jl
2
module TrixiMakieExt
3
4
# Required for visualization code
5
using Makie: Makie, GeometryBasics
6
7
# Use all exported symbols to avoid having to rewrite `recipes_makie.jl`
8
using Trixi
9
10
# Use additional symbols that are not exported
11
using Trixi: PlotData2DTriangulated, TrixiODESolution, PlotDataSeries, ScalarData, @muladd,
12
wrap_array_native, mesh_equations_solver_cache
13
14
# Import functions such that they can be extended with new methods
15
import Trixi: iplot, iplot!
16
17
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
18
# Since these FMAs can increase the performance of many numerical algorithms,
19
# we need to opt-in explicitly.
20
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
21
@muladd begin
22
#! format: noindent
23
24
# First some utilities
25
# Given a reference plotting triangulation, this function generates a plotting triangulation for
26
# the entire global mesh. The output can be plotted using `Makie.mesh`.
27
function global_plotting_triangulation_makie(pds::PlotDataSeries{<:PlotData2DTriangulated};
28
set_z_coordinate_zero = false)
29
@unpack variable_id = pds
30
pd = pds.plot_data
31
@unpack x, y, data, t = pd
32
33
makie_triangles = Makie.to_triangles(t)
34
35
# trimesh[i] holds GeometryBasics.Mesh containing plotting information on the ith element.
36
# Note: Float32 is required by GeometryBasics
37
num_plotting_nodes, num_elements = size(x)
38
trimesh = Vector{GeometryBasics.Mesh{3, Float32}}(undef, num_elements)
39
coordinates = zeros(Float32, num_plotting_nodes, 3)
40
for element in Base.OneTo(num_elements)
41
for i in Base.OneTo(num_plotting_nodes)
42
coordinates[i, 1] = x[i, element]
43
coordinates[i, 2] = y[i, element]
44
if set_z_coordinate_zero == false
45
coordinates[i, 3] = data[i, element][variable_id]
46
end
47
end
48
trimesh[element] = GeometryBasics.normal_mesh(Makie.to_vertices(coordinates),
49
makie_triangles)
50
end
51
plotting_mesh = merge([trimesh...]) # merge meshes on each element into one large mesh
52
return plotting_mesh
53
end
54
55
# Returns a list of `Makie.Point`s which can be used to plot the mesh, or a solution "wireframe"
56
# (e.g., a plot of the mesh lines but with the z-coordinate equal to the value of the solution).
57
function convert_PlotData2D_to_mesh_Points(pds::PlotDataSeries{<:PlotData2DTriangulated};
58
set_z_coordinate_zero = false)
59
@unpack variable_id = pds
60
pd = pds.plot_data
61
@unpack x_face, y_face, face_data = pd
62
63
if set_z_coordinate_zero
64
# plot 2d surface by setting z coordinate to zero.
65
# Uses `x_face` since `face_data` may be `::Nothing`, as it's not used for 2D plots.
66
sol_f = zeros(eltype(first(x_face)), size(x_face))
67
else
68
sol_f = StructArrays.component(face_data, variable_id)
69
end
70
71
# This line separates solution lines on each edge by NaNs to ensure that they are rendered
72
# separately. The coordinates `xf`, `yf` and the solution `sol_f`` are assumed to be a matrix
73
# whose columns correspond to different elements. We add NaN separators by appending a row of
74
# NaNs to this matrix. We also flatten (e.g., apply `vec` to) the result, as this speeds up
75
# plotting.
76
xyz_wireframe = GeometryBasics.Point.(map(x -> vec(vcat(x,
77
fill(NaN, 1, size(x, 2)))),
78
(x_face, y_face, sol_f))...)
79
80
return xyz_wireframe
81
end
82
83
# Creates a GeometryBasics triangulation for the visualization of a ScalarData2D plot object.
84
function global_plotting_triangulation_makie(pd::PlotData2DTriangulated{<:ScalarData};
85
set_z_coordinate_zero = false)
86
@unpack x, y, data, t = pd
87
88
makie_triangles = Makie.to_triangles(t)
89
90
# trimesh[i] holds GeometryBasics.Mesh containing plotting information on the ith element.
91
# Note: Float32 is required by GeometryBasics
92
num_plotting_nodes, num_elements = size(x)
93
trimesh = Vector{GeometryBasics.Mesh{3, Float32}}(undef, num_elements)
94
coordinates = zeros(Float32, num_plotting_nodes, 3)
95
for element in Base.OneTo(num_elements)
96
for i in Base.OneTo(num_plotting_nodes)
97
coordinates[i, 1] = x[i, element]
98
coordinates[i, 2] = y[i, element]
99
if set_z_coordinate_zero == false
100
coordinates[i, 3] = data.data[i, element]
101
end
102
end
103
trimesh[element] = GeometryBasics.normal_mesh(Makie.to_vertices(coordinates),
104
makie_triangles)
105
end
106
plotting_mesh = merge([trimesh...]) # merge meshes on each element into one large mesh
107
return plotting_mesh
108
end
109
110
# Returns a list of `GeometryBasics.Point`s which can be used to plot the mesh, or a solution "wireframe"
111
# (e.g., a plot of the mesh lines but with the z-coordinate equal to the value of the solution).
112
function convert_PlotData2D_to_mesh_Points(pd::PlotData2DTriangulated{<:ScalarData};
113
set_z_coordinate_zero = false)
114
@unpack x_face, y_face, face_data = pd
115
116
if set_z_coordinate_zero
117
# plot 2d surface by setting z coordinate to zero.
118
# Uses `x_face` since `face_data` may be `::Nothing`, as it's not used for 2D plots.
119
sol_f = zeros(eltype(first(x_face)), size(x_face))
120
else
121
sol_f = face_data
122
end
123
124
# This line separates solution lines on each edge by NaNs to ensure that they are rendered
125
# separately. The coordinates `xf`, `yf` and the solution `sol_f`` are assumed to be a matrix
126
# whose columns correspond to different elements. We add NaN separators by appending a row of
127
# NaNs to this matrix. We also flatten (e.g., apply `vec` to) the result, as this speeds up
128
# plotting.
129
xyz_wireframe = GeometryBasics.Point.(map(x -> vec(vcat(x,
130
fill(NaN, 1, size(x, 2)))),
131
(x_face, y_face, sol_f))...)
132
133
return xyz_wireframe
134
end
135
136
# We set the Makie default colormap to match Plots.jl, which uses `:inferno` by default.
137
default_Makie_colormap() = :inferno
138
139
# convenience struct for editing Makie plots after they're created.
140
struct FigureAndAxes{Axes}
141
fig::Makie.Figure
142
axes::Axes
143
end
144
145
# for "quiet" return arguments to Makie.plot(::TrixiODESolution) and
146
# Makie.plot(::PlotData2DTriangulated)
147
Base.show(io::IO, fa::FigureAndAxes) = nothing
148
149
# allows for returning fig, axes = Makie.plot(...)
150
function Base.iterate(fa::FigureAndAxes, state = 1)
151
if state == 1
152
return (fa.fig, 2)
153
elseif state == 2
154
return (fa.axes, 3)
155
else
156
return nothing
157
end
158
end
159
160
"""
161
iplot(u, mesh::UnstructuredMesh2D, equations, solver, cache;
162
plot_mesh=true, show_axis=false, colormap=default_Makie_colormap(),
163
variable_to_plot_in=1)
164
165
Creates an interactive surface plot of the solution and mesh for an `UnstructuredMesh2D` type.
166
167
Keywords:
168
- variable_to_plot_in: variable to show by default
169
170
!!! warning "Experimental implementation"
171
This is an experimental feature and may change in future releases.
172
"""
173
function iplot end
174
175
# Enables `iplot(PlotData2D(sol))`.
176
function iplot(pd::PlotData2DTriangulated;
177
plot_mesh = true, show_axis = false, colormap = default_Makie_colormap(),
178
variable_to_plot_in = 1)
179
@unpack variable_names = pd
180
181
# Initialize a Makie figure that we'll add the solution and toggle switches to.
182
fig = Makie.Figure()
183
184
# Set up options for the drop-down menu
185
menu_options = [zip(variable_names, 1:length(variable_names))...]
186
menu = Makie.Menu(fig, options = menu_options)
187
188
# Initialize toggle switches for viewing the mesh
189
toggle_solution_mesh = Makie.Toggle(fig, active = plot_mesh)
190
toggle_mesh = Makie.Toggle(fig, active = plot_mesh)
191
192
# Add dropdown menu and toggle switches to the left side of the figure.
193
fig[1, 1] = Makie.vgrid!(Makie.Label(fig, "Solution field", width = nothing), menu,
194
Makie.Label(fig, "Solution mesh visible"),
195
toggle_solution_mesh,
196
Makie.Label(fig, "Mesh visible"), toggle_mesh;
197
tellheight = false, width = 200)
198
199
# Create a zoomable interactive axis object on top of which to plot the solution.
200
ax = Makie.LScene(fig[1, 2], scenekw = (show_axis = show_axis,))
201
202
# Initialize the dropdown menu to `variable_to_plot_in`
203
# Since menu.selection is an Observable type, we need to dereference it using `[]` to set.
204
menu.selection[] = variable_to_plot_in
205
menu.i_selected[] = variable_to_plot_in
206
207
# Since `variable_to_plot` is an Observable, these lines are re-run whenever `variable_to_plot[]`
208
# is updated from the drop-down menu.
209
plotting_mesh = Makie.@lift(global_plotting_triangulation_makie(getindex(pd,
210
variable_names[$(menu.selection)])))
211
solution_z = Makie.@lift(getindex.($plotting_mesh.position, 3))
212
213
# Plot the actual solution.
214
Makie.mesh!(ax, plotting_mesh; color = solution_z, colormap)
215
216
# Create a mesh overlay by plotting a mesh both on top of and below the solution contours.
217
wire_points = Makie.@lift(convert_PlotData2D_to_mesh_Points(getindex(pd,
218
variable_names[$(menu.selection)])))
219
wire_mesh_top = Makie.lines!(ax, wire_points, color = :white,
220
visible = toggle_solution_mesh.active)
221
wire_mesh_bottom = Makie.lines!(ax, wire_points, color = :white,
222
visible = toggle_solution_mesh.active)
223
Makie.translate!(wire_mesh_top, 0, 0, 1e-3)
224
Makie.translate!(wire_mesh_bottom, 0, 0, -1e-3)
225
226
# This draws flat mesh lines below the solution.
227
function compute_z_offset(solution_z)
228
zmin = minimum(solution_z)
229
zrange = (x -> x[2] - x[1])(extrema(solution_z))
230
return zmin - 0.25 * zrange
231
end
232
z_offset = Makie.@lift(compute_z_offset($solution_z))
233
function get_flat_points(wire_points, z_offset)
234
[Makie.Point(point.data[1:2]..., z_offset) for point in wire_points]
235
end
236
flat_wire_points = Makie.@lift get_flat_points($wire_points, $z_offset)
237
wire_mesh_flat = Makie.lines!(ax, flat_wire_points, color = :black,
238
visible = toggle_mesh.active)
239
240
# create a small variation in the extrema to avoid the Makie `range_step` cannot be zero error.
241
# see https://github.com/MakieOrg/Makie.jl/issues/931 for more details.
242
# the colorbar range is perturbed by 1e-5 * the magnitude of the solution.
243
function scaled_extrema(x)
244
ex = extrema(x)
245
if ex[2] ex[1] # if solution is close to constant, perturb colorbar
246
return ex .+ 1e-5 .* maximum(abs.(ex)) .* (-1, 1)
247
else
248
return ex
249
end
250
end
251
252
# Resets the colorbar each time the solution changes.
253
Makie.Colorbar(fig[1, 3], limits = Makie.@lift(scaled_extrema($solution_z)),
254
colormap = colormap)
255
256
# On OSX, shift-command-4 for screenshots triggers a constant "up-zoom".
257
# To avoid this, we remap up-zoom to the right shift button instead.
258
Makie.cameracontrols(ax.scene).controls.up_key = Makie.Keyboard.right_shift
259
260
# typing this pulls up the figure (similar to display(plot!()) in Plots.jl)
261
fig
262
end
263
264
function iplot(u, mesh, equations, solver, cache;
265
solution_variables = nothing, nvisnodes = 2 * nnodes(solver), kwargs...)
266
@assert ndims(mesh) == 2
267
268
pd = PlotData2DTriangulated(u, mesh, equations, solver, cache;
269
solution_variables = solution_variables,
270
nvisnodes = nvisnodes)
271
272
iplot(pd; kwargs...)
273
end
274
275
# redirect `iplot(sol)` to dispatchable `iplot` signature.
276
iplot(sol::TrixiODESolution; kwargs...) = iplot(sol.u[end], sol.prob.p; kwargs...)
277
function iplot(u, semi; kwargs...)
278
iplot(wrap_array_native(u, semi), mesh_equations_solver_cache(semi)...; kwargs...)
279
end
280
281
# Interactive visualization of user-defined ScalarData.
282
function iplot(pd::PlotData2DTriangulated{<:ScalarData};
283
show_axis = false, colormap = default_Makie_colormap(),
284
plot_mesh = false)
285
fig = Makie.Figure()
286
287
# Create a zoomable interactive axis object on top of which to plot the solution.
288
ax = Makie.LScene(fig[1, 1], scenekw = (show_axis = show_axis,))
289
290
# plot the user-defined ScalarData
291
fig_axis_plt = iplot!(FigureAndAxes(fig, ax), pd; colormap = colormap,
292
plot_mesh = plot_mesh)
293
294
fig
295
return fig_axis_plt
296
end
297
298
function iplot!(fig_axis::Union{FigureAndAxes, Makie.FigureAxisPlot},
299
pd::PlotData2DTriangulated{<:ScalarData};
300
colormap = default_Makie_colormap(), plot_mesh = false)
301
302
# destructure first two fields of either FigureAndAxes or Makie.FigureAxisPlot
303
fig, ax = fig_axis
304
305
# create triangulation of the scalar data to plot
306
plotting_mesh = global_plotting_triangulation_makie(pd)
307
solution_z = getindex.(plotting_mesh.position, 3)
308
plt = Makie.mesh!(ax, plotting_mesh; color = solution_z, colormap)
309
310
if plot_mesh
311
wire_points = convert_PlotData2D_to_mesh_Points(pd)
312
wire_mesh_top = Makie.lines!(ax, wire_points, color = :white)
313
wire_mesh_bottom = Makie.lines!(ax, wire_points, color = :white)
314
Makie.translate!(wire_mesh_top, 0, 0, 1e-3)
315
Makie.translate!(wire_mesh_bottom, 0, 0, -1e-3)
316
end
317
318
# Add a colorbar to the rightmost part of the layout
319
Makie.Colorbar(fig[1, end + 1], plt)
320
321
fig
322
return Makie.FigureAxisPlot(fig, ax, plt)
323
end
324
325
# ================== new Makie plot recipes ====================
326
327
# This initializes a Makie recipe, which creates a new type definition which Makie uses to create
328
# custom `trixiheatmap` plots. See also https://docs.makie.org/stable/documentation/recipes/
329
Makie.@recipe(TrixiHeatmap, plot_data_series) do scene
330
Makie.Theme(colormap = default_Makie_colormap())
331
end
332
333
function Makie.plot!(myplot::TrixiHeatmap)
334
pds = myplot[:plot_data_series][]
335
336
plotting_mesh = global_plotting_triangulation_makie(pds;
337
set_z_coordinate_zero = true)
338
339
pd = pds.plot_data
340
solution_z = vec(StructArrays.component(pd.data, pds.variable_id))
341
Makie.mesh!(myplot, plotting_mesh, color = solution_z, shading = Makie.NoShading,
342
colormap = myplot[:colormap])
343
myplot.colorrange = extrema(solution_z)
344
345
# Makie hides keyword arguments within `myplot`; see also
346
# https://github.com/JuliaPlots/Makie.jl/issues/837#issuecomment-845985070
347
plot_mesh = if haskey(myplot, :plot_mesh)
348
myplot.plot_mesh[]
349
else
350
true # default to plotting the mesh
351
end
352
353
if plot_mesh
354
xyz_wireframe = convert_PlotData2D_to_mesh_Points(pds;
355
set_z_coordinate_zero = true)
356
Makie.lines!(myplot, xyz_wireframe, color = :lightgrey)
357
end
358
359
myplot
360
end
361
362
# redirects Makie.plot(pd::PlotDataSeries) to custom recipe TrixiHeatmap(pd)
363
Makie.plottype(::Trixi.PlotDataSeries{<:Trixi.PlotData2DTriangulated}) = TrixiHeatmap
364
365
# Makie does not yet support layouts in its plot recipes, so we overload `Makie.plot` directly.
366
function Makie.plot(sol::TrixiODESolution;
367
plot_mesh = false, solution_variables = nothing,
368
colormap = default_Makie_colormap())
369
return Makie.plot(PlotData2DTriangulated(sol; solution_variables); plot_mesh,
370
colormap)
371
end
372
373
function Makie.plot(pd::PlotData2DTriangulated, fig = Makie.Figure();
374
plot_mesh = false, colormap = default_Makie_colormap())
375
figAxes = Makie.plot!(fig, pd; plot_mesh, colormap)
376
display(figAxes.fig)
377
return figAxes
378
end
379
380
function Makie.plot!(fig, pd::PlotData2DTriangulated;
381
plot_mesh = false, colormap = default_Makie_colormap())
382
# Create layout that is as square as possible, when there are more than 3 subplots.
383
# This is done with a preference for more columns than rows if not.
384
if length(pd) <= 3
385
cols = length(pd)
386
rows = 1
387
else
388
cols = ceil(Int, sqrt(length(pd)))
389
rows = cld(length(pd), cols)
390
end
391
392
axes = [Makie.Axis(fig[i, j], xlabel = "x", ylabel = "y")
393
for j in 1:rows, i in 1:cols]
394
row_list, col_list = ([i for j in 1:rows, i in 1:cols],
395
[j for j in 1:rows, i in 1:cols])
396
397
for (variable_to_plot, (variable_name, pds)) in enumerate(pd)
398
ax = axes[variable_to_plot]
399
plt = trixiheatmap!(ax, pds; plot_mesh, colormap)
400
401
row = row_list[variable_to_plot]
402
col = col_list[variable_to_plot]
403
Makie.Colorbar(fig[row, col][1, 2], colormap = colormap)
404
405
ax.aspect = Makie.DataAspect() # equal aspect ratio
406
ax.title = variable_name
407
Makie.xlims!(ax, extrema(pd.x))
408
Makie.ylims!(ax, extrema(pd.y))
409
end
410
411
return FigureAndAxes(fig, axes)
412
end
413
end # @muladd
414
415
end
416
417