API Reference
On this page the docs of the provided functions are listed
FerriteViz.MakiePlotter — TypeMakiePlotter(dh::Ferrite.AbstractDofHandler, u::Vector)
MakiePlotter(dh::Ferrite.AbstractDofHandler, u::Vector, topology::TOP) where {TOP<:Ferrite.AbstractTopology}Builds a static triangulation of the underlying grid in dh.grid for rendering via Makie. The triangulation acts as a "L2" triangulation, i.e. the nodes which are shared between elements in the mesh are doubled.
For large 3D grids, prefer to use the second constructor if you have already a topology. Otherwise, it will be rebuilt which is time consuming.
FerriteViz.solutionplot — Functionsolutionplot(plotter::MakiePlotter; kwargs...)
solutionplot(dh::AbstractDofHandler, u::Vector; kwargs...)
solutionplot!(plotter::MakiePlotter; kwargs...)
solutionplot!(dh::AbstractDofHandler, u::Vector; kwargs...)Solutionplot produces the classical contour plot onto the finite element mesh. Most important keyword arguments are:
field::Symbol=:defaultrepresenting the field which gets plotted, defaults to the first field in thedh.deformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationprocess::Function=postprocessfunction to construct nodal scalar values from a vector valued problemcolormap::Symbol=:cividiscolorrange::NTuple{2,<:Number}: Specify (min, max) of the colorscale. If not given, min and max are calculated automatically from the data.deformation_scale=1.0shading=falsescale_plot=falsetransparent=falsenan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.solutionplot! — Functionsolutionplot(plotter::MakiePlotter; kwargs...)
solutionplot(dh::AbstractDofHandler, u::Vector; kwargs...)
solutionplot!(plotter::MakiePlotter; kwargs...)
solutionplot!(dh::AbstractDofHandler, u::Vector; kwargs...)Solutionplot produces the classical contour plot onto the finite element mesh. Most important keyword arguments are:
field::Symbol=:defaultrepresenting the field which gets plotted, defaults to the first field in thedh.deformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationprocess::Function=postprocessfunction to construct nodal scalar values from a vector valued problemcolormap::Symbol=:cividiscolorrange::NTuple{2,<:Number}: Specify (min, max) of the colorscale. If not given, min and max are calculated automatically from the data.deformation_scale=1.0shading=falsescale_plot=falsetransparent=falsenan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.cellplot — Functioncellplot(plotter::MakiePlotter,σ::Vector{T}; kwargs...) where T
cellplot!(plotter::MakiePlotter,σ::Vector{T}; kwargs...) where Tcellplot plots constant scalar data on the cells of the finite element mesh. If T is not a number, the keyword argument process can be passed in order to reduce the elements of σ to a scalar.
keyword arguments are:
deformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationprocess::Function=identityfunction to construct cell scalar values. Defaults toidentity, i.e. scalar values.colormap::Symbol=:cividiscolorrange::NTuple{2,<:Number}: Specify (min, max) of the colorscale. If not given, min and max are calculated automatically from the data.deformation_scale=1.0shading=falsescale_plot=falsetransparent=falsenan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.cellplot! — Functioncellplot(plotter::MakiePlotter,σ::Vector{T}; kwargs...) where T
cellplot!(plotter::MakiePlotter,σ::Vector{T}; kwargs...) where Tcellplot plots constant scalar data on the cells of the finite element mesh. If T is not a number, the keyword argument process can be passed in order to reduce the elements of σ to a scalar.
keyword arguments are:
deformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationprocess::Function=identityfunction to construct cell scalar values. Defaults toidentity, i.e. scalar values.colormap::Symbol=:cividiscolorrange::NTuple{2,<:Number}: Specify (min, max) of the colorscale. If not given, min and max are calculated automatically from the data.deformation_scale=1.0shading=falsescale_plot=falsetransparent=falsenan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.wireframe — Functionwireframe(plotter::MakiePlotter; kwargs...)
wireframe(dh::AbstractDofHandler, u::Vector; kwargs...)
wireframe(grid::AbstractGrid; kwargs...)
wireframe!(plotter::MakiePlotter; kwargs...)
wireframe!(dh::AbstractDofHandler, u::Vector; kwargs...)
wireframe!(grid::AbstractGrid; kwargs...)Plots the finite element mesh, optionally labels it and transforms it if a suitable deformation_field is given.
plotnodes::Bool=trueplots the nodes as circles/spheresstrokewidth::Int=2how thick faces/edges are drawncolor::Symbol=theme(scene,:linecolor)color of the faces/edges and nodesmarkersize::Int=30size of the nodesdeformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationdeformation_scale::Number=1.0scaling of the deformationcellsets=falseColor cells based on their cellset association. If no cellset is found for a cell, the cell is marked blue.nodelables=falseglobal node id labelsnodelabelcolor=:darkbluecelllabels=falseglobal cell id labelscelllabelcolor=:darkredtextsize::Int=15size of the label's textvisible=true
FerriteViz.wireframe! — Functionwireframe(plotter::MakiePlotter; kwargs...)
wireframe(dh::AbstractDofHandler, u::Vector; kwargs...)
wireframe(grid::AbstractGrid; kwargs...)
wireframe!(plotter::MakiePlotter; kwargs...)
wireframe!(dh::AbstractDofHandler, u::Vector; kwargs...)
wireframe!(grid::AbstractGrid; kwargs...)Plots the finite element mesh, optionally labels it and transforms it if a suitable deformation_field is given.
plotnodes::Bool=trueplots the nodes as circles/spheresstrokewidth::Int=2how thick faces/edges are drawncolor::Symbol=theme(scene,:linecolor)color of the faces/edges and nodesmarkersize::Int=30size of the nodesdeformation_field::Symbol=:defaultfield that transforms the mesh by the given deformation, defaults to no deformationdeformation_scale::Number=1.0scaling of the deformationcellsets=falseColor cells based on their cellset association. If no cellset is found for a cell, the cell is marked blue.nodelables=falseglobal node id labelsnodelabelcolor=:darkbluecelllabels=falseglobal cell id labelscelllabelcolor=:darkredtextsize::Int=15size of the label's textvisible=true
FerriteViz.arrows — Functionarrows(plotter::MakiePlotter; kwargs...)
arrows(dh::AbstractDofHandler, u::Vector; kwargs...)
arrows!(plotter::MakiePlotter; kwargs...)
arrows!(dh::AbstractDofHandler, u::Vector; kwargs...)At every node position a arrows is drawn, where the arrow tip ends at the node. Only works in dim >=2. If a color is specified the arrows are unicolored. Otherwise the color corresponds to the magnitude, or any other scalar value based on the process function.
arrowsize = 0.08normalize = truefield = :defaultcolor = :defaultcolormap = :cividisprocess=postprocesslengthscale = 1f0
FerriteViz.arrows! — Functionarrows(plotter::MakiePlotter; kwargs...)
arrows(dh::AbstractDofHandler, u::Vector; kwargs...)
arrows!(plotter::MakiePlotter; kwargs...)
arrows!(dh::AbstractDofHandler, u::Vector; kwargs...)At every node position a arrows is drawn, where the arrow tip ends at the node. Only works in dim >=2. If a color is specified the arrows are unicolored. Otherwise the color corresponds to the magnitude, or any other scalar value based on the process function.
arrowsize = 0.08normalize = truefield = :defaultcolor = :defaultcolormap = :cividisprocess=postprocesslengthscale = 1f0
FerriteViz.surface — Functionsurface(plotter::MakiePlotter; kwargs...)
surface(dh::AbstractDofHandler, u::Vector; kwargs...)
surface!(plotter::MakiePlotter; kwargs...)
surface!(dh::AbstractDofHandler, u::Vector; kwargs...)Uses the given field and plots the scalar values as a surface. If it's a vector valued problem, the nodal vector values are transformed to a scalar based on process which defaults to the magnitude. Only availble in dim=2.
field = :defaultprocess = postprocessscale_plot = falseshading = falsecolormap = :cividiscolorrange=Makie.automaticnan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.surface! — Functionsurface(plotter::MakiePlotter; kwargs...)
surface(dh::AbstractDofHandler, u::Vector; kwargs...)
surface!(plotter::MakiePlotter; kwargs...)
surface!(dh::AbstractDofHandler, u::Vector; kwargs...)Uses the given field and plots the scalar values as a surface. If it's a vector valued problem, the nodal vector values are transformed to a scalar based on process which defaults to the magnitude. Only availble in dim=2.
field = :defaultprocess = postprocessscale_plot = falseshading = falsecolormap = :cividiscolorrange=Makie.automaticnan_color::Union{Symbol, <:Colorant}=:red
FerriteViz.elementinfo — Functionelementinfo(ip::Interpolation; kwargs...)
elementinfo(cell::AbstractCell; kwargs...)
elementinfo(ip::Type{Interpolation}; kwargs...)
elementinfo(cell::Type{AbstractCell}; kwargs...)plotnodes=truecontrols if nodes of element are plottedstrokewidth=2strokwidth of faces/edgescolor=theme(scene, :linecolor)markersize=30size of the nodestextsize=60textsize of node-, edges- and facelabelsnodelabels=trueswitch that controls plotting of nodelabelsnodelabelcolor=:darkrednodelabeloffset=(0.0,0.0)offset of the nodelabel text relative to its associated nodefacelabels=trueswitch that controls plotting of facelabelsfacelabelcolor=:darkgreenfacelabeloffset=(-40,0)offset of the facelabel text relative to its associated face middlepointedgelabels=trueswitch that controls plotting of edgelabelsedgelabelcolor=:darkblueedgelabeloffset=(-40,-40)offset of the edgelabel text relative to its associated edge middlepointfont="Julia Mono"font of the node-, edge-, and facelabels
FerriteViz.elementinfo! — Functionelementinfo(ip::Interpolation; kwargs...)
elementinfo(cell::AbstractCell; kwargs...)
elementinfo(ip::Type{Interpolation}; kwargs...)
elementinfo(cell::Type{AbstractCell}; kwargs...)plotnodes=truecontrols if nodes of element are plottedstrokewidth=2strokwidth of faces/edgescolor=theme(scene, :linecolor)markersize=30size of the nodestextsize=60textsize of node-, edges- and facelabelsnodelabels=trueswitch that controls plotting of nodelabelsnodelabelcolor=:darkrednodelabeloffset=(0.0,0.0)offset of the nodelabel text relative to its associated nodefacelabels=trueswitch that controls plotting of facelabelsfacelabelcolor=:darkgreenfacelabeloffset=(-40,0)offset of the facelabel text relative to its associated face middlepointedgelabels=trueswitch that controls plotting of edgelabelsedgelabelcolor=:darkblueedgelabeloffset=(-40,-40)offset of the edgelabel text relative to its associated edge middlepointfont="Julia Mono"font of the node-, edge-, and facelabels
FerriteViz.ferriteviewer — Functionferriteviewer(plotter::MakiePlotter)
ferriteviewer(plotter::MakiePlotter, u_history::Vector{Vector{T}}})Constructs a viewer with a solutionplot, Colorbar as well as sliders,toggles and menus to change the current view. If the second dispatch is called a timeslider is added, in order to step through a set of solutions obtained from a simulation.
FerriteViz.update! — FunctionFerriteViz.update!(plotter::MakiePlotter, u::Vector)Updates the Observable plotter.u and thereby, triggers the plot to update.
FerriteViz.for_discretization — FunctionCreate a first order discretization w.r.t. a field and transfer the solution.
FerriteViz.for_interpolation — FunctionGet the interpolation of the first order refinement.
FerriteViz.interpolate_gradient_field — Functioninterpolate_gradient_field(dh::DofHandler, u::AbstractVector, field_name::Symbol; copy_fields::Vector{Symbol})Compute the piecewise discontinuous gradient field for field_name. Returns the flux dof handler and the corresponding flux dof values. If the additional keyword argument copy_fields is provided with a non empty Vector{Symbol}, the corresponding fields of dh will be copied into the returned flux dof handler and flux dof value vector.
FerriteViz.uniform_refinement — Functionuniform_refinement(plotter::MakiePlotter)
uniform_refinement(plotter::MakiePlotter, num_refinements::Int)Generates 3 triangles for each triangle by adding a center vertex and connecting them (orientation preserving).
This method has high RAM usage!
This function currently does not increase the resolution of the geometrical points in space, only the solution quality!
Details
TODO investigate whether it is possible to eliminate the coordinate duplication without trashing the caches
FerriteViz.crinkle_clip! — Functioncrinkle_clip!(plotter::MakiePlotter{3}, decision_fun)Crinkle clip updates the visibility of the triangles, based on an implicit description of the clipping surface. Here decision_fun takes the grid and a cell index as input and returns whether the cell is visible or not.
Chained calls to crinkle_clip! won't work at the moment.
FerriteViz.crinkle_clip — Functioncrinkle_clip(plotter::MakiePlotter{3}, decision_fun) -> MakiePlotterCrinkle clip generates a new plotter with updated visibility of the triangles. Non-mutating version of crinkle_clip!. Note that chained calls to crinkle_clip won't work.
FerriteViz.ClipPlane — TypeClipPlane{T}(normal, distance_to_origin)Clip plane with data of type T described by the normal and its distance to the coordinate origin.
Details
INTERNAL: Instances are callable as plane(grid, cellid) returning true if all nodes of a cell are on the side of the positive side of the plane, i.e. where for its normal we have coord ⋅ plane.normal > plane.distance. With this helper we perform the crinkle clip internally. xref [crinkle_clip!](@ref)