Postprocessing
Project to nodes
Ferrite.L2Projector
— TypeL2Projector(func_ip::Interpolation, grid::AbstractGrid; kwargs...)
Create an L2Projector
used for projecting quadrature data. func_ip
is the function interpolation used for the projection and grid
the grid over which the projection is applied.
Keyword arguments:
qr_lhs
: quadrature for the left hand side. Defaults to a quadrature which exactly integrates a mass matrix withfunc_ip
as the interpolation.set
: element set over which the projection applies. Defaults to all elements in the grid.geom_ip
: geometric interpolation. Defaults to the default interpolation for the grid.
The L2Projector
acts as the integrated left hand side of the projection equation: Find projection $u \in L_2(\Omega)$ such that
\[\int v u \ \mathrm{d}\Omega = \int v f \ \mathrm{d}\Omega \quad \forall v \in L_2(\Omega),\]
where $f$ is the data to project.
Use project
to integrate the right hand side and solve for the system.
Ferrite.project
— Functionproject(proj::L2Projector, vals, qr_rhs::QuadratureRule; project_to_nodes=true)
Makes a L2 projection of data vals
to the nodes of the grid using the projector proj
(see L2Projector
).
project
integrates the right hand side, and solves the projection $u$ from the following projection equation: Find projection $u \in L_2(\Omega)$ such that
\[\int v u \ \mathrm{d}\Omega = \int v f \ \mathrm{d}\Omega \quad \forall v \in L_2(\Omega),\]
where $f$ is the data to project, i.e. vals
.
The data vals
should be a vector, with length corresponding to number of elements, of vectors, with length corresponding to number of quadrature points per element, matching the number of points in qr_rhs
. Alternatively, vals
can be a matrix, with number of columns corresponding to number of elements, and number of rows corresponding to number of points in qr_rhs
. Example (scalar) input data:
vals = [
[0.44, 0.98, 0.32], # data for quadrature point 1, 2, 3 of element 1
[0.29, 0.48, 0.55], # data for quadrature point 1, 2, 3 of element 2
# ...
]
or equivalent in matrix form:
vals = [
0.44 0.29 # ...
0.98 0.48 # ...
0.32 0.55 # ...
]
Supported data types to project are Number
s and AbstractTensor
s.
If the parameter project_to_nodes
is true
, then the projection returns the values in the order of the mesh nodes (suitable format for exporting). If false
, it returns the values corresponding to the degrees of freedom for a scalar field over the domain, which is useful if one wants to interpolate the projected values.
Postprocessing
Ferrite.PointEvalHandler
— TypePointEvalHandler(grid::Grid, points::AbstractVector{Vec{dim,T}}; kwargs...) where {dim, T}
The PointEvalHandler
can be used for function evaluation in arbitrary points in the domain – not just in quadrature points or nodes.
The PointEvalHandler
takes the following keyword arguments:
search_nneighbors
: How many nodes should be found in the nearest neighbor search for each point. Usually there is no need to change this setting. Default value:3
.warn
: Show a warning if a point is not found. Default value:true
.
The constructor takes a grid and a vector of coordinates for the points. The PointEvalHandler
computes i) the corresponding cell, and ii) the (local) coordinate within the cell, for each point. The fields of the PointEvalHandler
are:
cells::Vector{Union{Int,Nothing}}
: vector with cell IDs for the points, withnothing
for points that could not be found.local_coords::Vector{Union{Vec,Nothing}}
: vector with the local coordinates (i.e. coordinates in the reference configuration) for the points, withnothing
for points that could not be found.
There are two ways to use the PointEvalHandler
to evaluate functions:
get_point_values
: can be used when the function is described by i) adh::DofHandler
+uh::Vector
(for example the FE-solution), or ii) ap::L2Projector
+ph::Vector
(for projected data).- Iteration with
PointIterator
+PointValues
: can be used for more flexible evaluation in the points, for example to compute gradients.
Ferrite.get_point_values
— Functionget_point_values(ph::PointEvalHandler, dh::AbstractDofHandler, dof_values::Vector{T}, [fieldname::Symbol]) where T
get_point_values(ph::PointEvalHandler, proj::L2Projector, dof_values::Vector{T}) where T
Return a Vector{T}
(for a 1-dimensional field) or a Vector{Vec{fielddim, T}}
(for a vector field) with the field values of field fieldname
in the points of the PointEvalHandler
. The fieldname
can be omitted if only one field is stored in dh
. The field values are computed based on the dof_values
and interpolated to the local coordinates by the function interpolation of the corresponding field
stored in the AbstractDofHandler
or the L2Projector
.
Points that could not be found in the domain when constructing the PointEvalHandler
will have NaN
s for the corresponding entries in the output vector.
Ferrite.PointValues
— TypePointScalarValues(cv::CellScalarValues)
PointScalarValues(ip_f::Interpolation, ip_g::Interpolation=ip_f)
PointVectorValues(cv::CellVectorValues)
PointVectorValues(ip_f::Interpolation, ip_g::Interpolation=ip_f)
Similar to CellScalarValues
and CellVectorValues
but with a single updateable "quadrature point". PointValues
are used for evaluation of functions/gradients in arbitrary points of the domain together with a PointEvalHandler
.
PointValues
can be created from CellValues
, or from the interpolations directly.
PointValues
are reinitialized like other CellValues
, but since the local reference coordinate of the "quadrature point" changes this needs to be passed to reinit!
, in addition to the element coordinates: reinit!(pv, coords, local_coord)
. Alternatively, it can be reinitialized with a PointLocation
when iterating a PointEvalHandler
with a PointIterator
.
For function/gradient evaluation, PointValues
are used in the same way as CellValues
, i.e. by using function_value
, function_gradient
, etc, with the exception that there is no need to specify the quadrature point index (since PointValues
only have 1, this is the default).
Ferrite.PointIterator
— TypePointIterator(ph::PointEvalHandler)
Create an iterator over the points in the PointEvalHandler
. The elements of the iterator are either a PointLocation
, if the corresponding point could be found in the grid, or nothing
, if the point was not found.
A PointLocation
can be used to query the cell ID with the cellid
function, and can be used to reinitialize PointValues
with reinit!
.
Examples
ph = PointEvalHandler(grid, points)
for point in PointIterator(ph)
point === nothing && continue # Skip any points that weren't found
reinit!(pointvalues, point) # Update pointvalues
# ...
end
Ferrite.PointLocation
— TypePointLocation
Element of a PointIterator
, typically used to reinitialize PointValues
. Fields:
cid::Int
: ID of the cell containing the pointlocal_coord::Vec
: the local (reference) coordinate of the pointcoords::Vector{Vec}
: the coordinates of the cell
Ferrite.reshape_to_nodes
— Functionreshape_to_nodes(dh::AbstractDofHandler, u::Vector{T}, fieldname::Symbol) where T
Reshape the entries of the dof-vector u
which correspond to the field fieldname
in nodal order. Return a matrix with a column for every node and a row for every dimension of the field. For superparametric fields only the entries corresponding to nodes of the grid will be returned. Do not use this function for subparametric approximations.
VTK Export
WriteVTK.vtk_grid
— Methodvtk_grid(filename::AbstractString, grid::Grid)
Create a unstructured VTK grid from a Grid
. Return a DatasetFile
which data can be appended to, see vtk_point_data
and vtk_cell_data
.
WriteVTK.vtk_point_data
— Functionvtk_point_data(vtk, data::Vector{<:AbstractTensor}, name)
Write the tensor field data
to the vtk file. Two-dimensional tensors are padded with zeros.
For second order tensors the following indexing ordering is used: [11, 22, 33, 23, 13, 12, 32, 31, 21]
. This is the default Voigt order in Tensors.jl.
Ferrite.vtk_cellset
— Functionvtk_cellset(vtk, grid::Grid)
Export all cell sets in the grid. Each cell set is exported with vtk_cell_data
with value 1 if the cell is in the set, and 0 otherwise.
vtk_cellset(vtk, grid::Grid, cellset::String)
Export the cell set specified by cellset
as cell data with value 1 if the cell is in the set and 0 otherwise.
Ferrite.vtk_cell_data_colors
— Functionvtk_cell_data_colors(vtkfile, cell_colors, name="coloring")
Write cell colors (see create_coloring
) to a VTK file for visualization.
In case of coloring a subset, the cells which are not part of the subset are represented as color 0.