Skip to content

Commit 33b19e1

Browse files
committed
plot broken scalar FE Types in broken plots
1 parent af2c83a commit 33b19e1

4 files changed

Lines changed: 91 additions & 6 deletions

File tree

CHANGELOG.md

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,11 @@
11
# CHANGES
22

3+
## v1.3.0
4+
5+
### Added
6+
- `plot` for discontinuous the FE Spaces `L2P0` and `L2P1` are now discontinuous to avoid averaging of values at the grid nodes.
7+
This can be disabled by passing the kwarg `average_broken_plots = true` to the `plot` call.
8+
39
## v1.2.0 June 4, 2025
410

511
### Changed

Project.toml

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,7 @@
11
name = "ExtendableFEM"
22
uuid = "a722555e-65e0-4074-a036-ca7ce79a4aed"
3+
version = "1.3.0"
34
authors = ["Christian Merdon <merdon@wias-berlin.de>", "Patrick Jaap <patrick.jaap@wias-berlin.de>"]
4-
version = "1.2"
55

66
[deps]
77
CommonSolve = "38540f10-b2f7-11e9-35d8-d573e4eb0ff2"
@@ -18,6 +18,7 @@ Printf = "de0858da-6303-5e67-8744-51eddeeeb8d7"
1818
SciMLBase = "0bca4576-84f4-4d90-8ffe-ffa030f20462"
1919
SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf"
2020
SparseDiffTools = "47a9eef4-7e08-11e9-0b38-333d64bd3804"
21+
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
2122
Symbolics = "0c5d862f-8b57-4792-8d23-62f2024744c7"
2223
TimerOutputs = "a759f4b9-e2f1-59dc-863e-4aeb61b1ea8f"
2324
UnicodePlots = "b8865327-cd53-5732-bb35-84acbb429228"
@@ -63,10 +64,9 @@ Metis = "2679e427-3c69-5b7f-982b-ece356f1e94b"
6364
OrdinaryDiffEqRosenbrock = "43230ef6-c299-4910-a778-202eb28ce4ce"
6465
OrdinaryDiffEqSDIRK = "2d112036-d095-4a1e-ab9a-08536f3ecdbf"
6566
SimplexGridFactory = "57bfcd06-606e-45d6-baf4-4ba06da0efd5"
66-
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
6767
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
6868
TetGen = "c5d3f3f7-f850-59f6-8a2e-ffc6dc1317ea"
6969
Triangulate = "f7e6ffb2-c36d-4f8f-a77e-16e897189344"
7070

7171
[targets]
72-
test = ["Aqua", "ExampleJuggler", "ExplicitImports", "IncompleteLU", "Metis", "OrdinaryDiffEqRosenbrock", "OrdinaryDiffEqSDIRK", "SimplexGridFactory", "StaticArrays", "Test", "TetGen", "Triangulate"]
72+
test = ["Aqua", "ExampleJuggler", "ExplicitImports", "IncompleteLU", "Metis", "OrdinaryDiffEqRosenbrock", "OrdinaryDiffEqSDIRK", "SimplexGridFactory", "Test", "TetGen", "Triangulate"]

src/ExtendableFEM.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,7 @@ using ExtendableGrids: ExtendableGrids, AT_NODES, AbstractElementGeometry,
4949
ON_BEDGES, ON_BFACES, ON_CELLS, ON_EDGES, ON_FACES,
5050
ON_IFACES, SerialVariableTargetAdjacency,
5151
UniqueBFaceGeometries, UniqueCellGeometries,
52-
UniqueFaceGeometries, append!, dim_element, eval_trafo!,
52+
UniqueFaceGeometries, append!, dim_element, dim_space, eval_trafo!,
5353
facetype_of_cellface, interpolate!,
5454
max_num_targets_per_source, num_cells, num_faces,
5555
num_nodes, num_sources, num_targets, simplexgrid,
@@ -74,6 +74,7 @@ using Symbolics: Symbolics
7474
using SciMLBase: SciMLBase
7575
using TimerOutputs: TimerOutput, print_timer, @timeit
7676
using UnicodePlots: UnicodePlots
77+
using StaticArrays: SVector
7778

7879
## reexport stuff from ExtendableFEMBase and ExtendableGrids
7980
export FESpace, FEMatrix, FEVector

src/plots.jl

Lines changed: 80 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -19,7 +19,20 @@ function plot!(p::GridVisualizer, ops, sol; kwargs...)
1919
Plots the operator evaluations ops of blocks in sol into the GridVisualizer.
2020
2121
"""
22-
function plot!(p::GridVisualizer, ops, sol; rasterpoints = 10, linewidth = 1, keep = [], ncols = size(p.subplots, 2), do_abs = true, do_vector_plots = true, title_add = "", kwargs...)
22+
function plot!(
23+
p::GridVisualizer,
24+
ops,
25+
sol;
26+
rasterpoints = 10,
27+
linewidth = 1,
28+
keep = [],
29+
ncols = size(p.subplots, 2),
30+
do_abs = true,
31+
do_vector_plots = true,
32+
title_add = "",
33+
average_broken_plots = false,
34+
kwargs...
35+
)
2336
col, row, id = 0, 1, 0
2437
for op in ops
2538
col += 1
@@ -47,7 +60,9 @@ function plot!(p::GridVisualizer, ops, sol; rasterpoints = 10, linewidth = 1, ke
4760
else
4861
title = op[2] == Identity ? "$(sol[op[1]].name)" : "$(op[2])($(sol[op[1]].name))"
4962
end
50-
if resultdim == 1
63+
if !average_broken_plots && is_scalar_broken_lagrange_FE(get_FEType(sol[op[1]].FES))
64+
broken_scalarplot!(p[row, col], sol[op[1]]; title = title * title_add, kwargs...)
65+
elseif resultdim == 1
5166
GridVisualize.scalarplot!(p[row, col], sol[op[1]].FES.dofgrid, view(nodevalues(sol[op[1]], op[2]; abs = false), 1, :), title = title * title_add; kwargs...)
5267
elseif do_abs == true
5368
GridVisualize.scalarplot!(p[row, col], sol[op[1]].FES.dofgrid, view(nodevalues(sol[op[1]], op[2]; abs = true), 1, :), title = "|" * title * "|" * title_add; kwargs...)
@@ -71,6 +86,69 @@ function plot!(p::GridVisualizer, ops, sol; rasterpoints = 10, linewidth = 1, ke
7186
return p
7287
end
7388

89+
# helper to detect broken Lagrange FES
90+
is_scalar_broken_lagrange_FE(FEType) = false
91+
is_scalar_broken_lagrange_FE(::Type{L2P0{1}}) = true
92+
is_scalar_broken_lagrange_FE(::Type{L2P1{1}}) = true
93+
94+
"""
95+
broken_scalarplot!(vis, feVectorBlock; kwargs...)
96+
97+
A "broken" scalarplot of a broken finite element vector (currently only L2P0 and L2P1 are supported)
98+
Instead of averaging the discontinuous values on the grid nodes, each grid cell is plotted
99+
independently. Thus, a discontinuous plot is generated.
100+
101+
All kwargs of the calling method are transferred to the scalarplot in this method.
102+
"""
103+
function broken_scalarplot!(vis, feVectorBlock::FEVectorBlock; kwargs...)
104+
105+
FES = feVectorBlock.FES
106+
dofgrid = FES.dofgrid
107+
dim = dim_space(dofgrid)
108+
109+
cell_nodes = dofgrid[CellNodes]
110+
cell_dofs = FES[CellDofs]
111+
coords = dofgrid[Coordinates]
112+
113+
nodes_per_cell = size(cell_nodes, 1)
114+
all_cells = Vector{SVector{nodes_per_cell, Int32}}(undef, 0)
115+
all_coords = Vector{SVector{dim, Float64}}(undef, 0)
116+
all_values = Vector{Float64}(undef, 0)
117+
118+
for i_cell in 1:size(cell_nodes, 2)
119+
# create a new independent cell
120+
i = nodes_per_cell * (i_cell - 1) + 1
121+
push!(all_cells, i:(i + nodes_per_cell - 1))
122+
123+
# copy the coords
124+
for n in @views cell_nodes[:, i_cell]
125+
@views push!(all_coords, coords[:, n][:])
126+
end
127+
128+
# extract values from feVector
129+
values = feVectorBlock.entries[cell_dofs[:, i_cell]]
130+
if length(values) == 1 # L2P0 !!
131+
values = fill(values[1], nodes_per_cell)
132+
end
133+
@views append!(all_values, values)
134+
end
135+
136+
# convert to matrices
137+
cell_matrix = Matrix{Int32}(undef, nodes_per_cell, length(all_cells))
138+
for i in 1:size(cell_matrix, 1), j in 1:size(cell_matrix, 2)
139+
cell_matrix[i, j] = all_cells[j][i]
140+
end
141+
142+
coord_matrix = Matrix{Float64}(undef, dim, length(all_coords))
143+
for i in 1:size(coord_matrix, 1), j in 1:size(coord_matrix, 2)
144+
coord_matrix[i, j] = all_coords[j][i]
145+
end
146+
147+
GridVisualize.scalarplot!(vis, coord_matrix, cell_matrix, all_values; kwargs...)
148+
149+
return nothing
150+
end
151+
74152

75153
"""
76154
````

0 commit comments

Comments
 (0)