Finite Element Interpolations
Source functions and QPInfo
The functions that can be interpolated with the methods below are expected to have a certain interface, i.e.:
function f!(result, qpinfo) end
The qpinfo argument communicates vast information of the current quadrature point:
qpinfo child | Type | Description |
---|---|---|
qpinfo.x | Vector{Real} | space coordinates of quadrature point |
qpinfo.time | Real | current time |
qpinfo.item | Integer | current item that contains qpinfo.x |
qpinfo.region | Integer | region number of item |
qpinfo.xref | Vector{Real} | reference coordinates within item of qpinfo.x |
qpinfo.volume | Real | volume of item |
qpinfo.params | Vector{Any} | parameters that can be transferred via keyword arguments |
Standard Interpolations
Each finite element type provides a standard interpolation routine that can be applied to user-defined source functions. By default, interpolation is performed over all cells, but it can also be restricted to faces or edges using an appropriate AssemblyType
.
ExtendableGrids.interpolate!
— Functionfunction ExtendableGrids.interpolate!(
target::FEVectorBlock{T, Tv, Ti},
AT::Type{<:AssemblyType},
source;
items = [],
kwargs...
) where {T, Tv, Ti}
Interpolate a function or data into the finite element space (i.e. computes the coefficients) associated with target
, using the specified assembly type.
Arguments
target::FEVectorBlock
: The block of the FE vector to store the interpolated coefficients.AT::Type{<:AssemblyType}
: The assembly type specifying where interpolation is performed (e.g.,ON_CELLS
).source
: The function or callable to interpolate. Should have the signaturesource!(result, qpinfo)
.
Keyword Arguments
items
: List of mesh entities (cells, faces, etc.) to interpolate on. If empty, all entities of the specified type are used.kwargs...
: Additional keyword arguments passed to lower-level routines (e.g.,bonus_quadorder
,time
).
Notes
- For "broken" FE spaces, interpolation is performed in a continuous auxiliary space and then mapped to the broken space.
- The
source!
function is called at each quadrature point and should fillresult
with the function values at that point. - The
qpinfo
argument provides information about the current quadrature point, including coordinates, weights, and possibly time.
function ExtendableGrids.interpolate!(target::FEVectorBlock,
source::Function;
kwargs...)
see interpolate!(target, ON_CELLS, source; kwargs...)
Additionally, you can transfer finite element functions from one grid to another using the lazy_interpolate!
routine, which interpolates between different meshes.
ExtendableFEMBase.lazy_interpolate!
— Functionfunction lazy_interpolate!(
target::FEVectorBlock{T1, Tv, Ti},
source,
operators = [(1, Identity)];
postprocess! = standard_kernel,
xtrafo! = nothing,
items = [],
resultdim = get_ncomponents(eltype(target.FES)),
not_in_domain_value = 1.0e30,
start_cell = 1,
only_localsearch = false,
use_cellparents::Bool = false,
eps = 1.0e-13,
kwargs...) where {T1, Tv, Ti}
Interpolates (operator-evaluations of) the given FEVector source (or an array of FEVectorBlocks) into the finite element space assigned to the target
FEVectorBlock.
The interpolation is performed using a point evaluation pattern and cell search. If CellParents
information is available in the target grid, enabling use_cellparents=true
can improve the efficiency of the search.
A custom postprocessing function can be provided via the postprocess!
argument, which should have the interface: postprocess!(result, input, qpinfo) where result
is the output buffer, input
is the operator evaluation, and qpinfo
provides quadrature point information.
If the source and target grids have different coordinate dimensions, a coordinate transformation function xtrafo!
must be provided, with the interface: xtrafo!(x_source, x) which maps coordinates x
from the target grid to coordinates in the source grid.
If a point cannot be found in the source grid, the value not_in_domain_value
is used as the function value. The items
argument can be used to restrict the interpolation to specific target cells.
Arguments
target::FEVectorBlock
: The target finite element vector block.source
: The source array of FEVectorBlocks.operators
: Array of operator argument tuples (source block tag, operator type) (default:[(1, Identity)]
).
Keyword Arguments
postprocess!
: Function to postprocess operator evaluations (default:standard_kernel
).xtrafo!
: Optional coordinate transformation function (default:nothing
).items
: List of target cells to interpolate (default:[]
).resultdim
: Result dimension (default:get_ncomponents(eltype(target.FES))
).not_in_domain_value
: Value assigned if a point is not found in the source domain (default:1.0e30
).start_cell
: Starting cell index for cell search (default:1
).only_localsearch
: Restrict cell search to local neighborhood (default:false
).use_cellparents
: Use parent cell information for search (default:false
).eps
: Tolerance for cell search (default:1.0e-13
).kwargs...
: Additional keyword arguments passed tointerpolate!
.
Notes
- Discontinuous quantities at target grid vertices are evaluated in the first found cell of the source grid; no averaging is performed.
- The function is not the most efficient for large-scale problems due to its reliance on pointwise cell search.
The following function continuously interpolates finite element function into a H1Pk space by point evaluations at the Lagrange nodes of the H1Pk element (averaged over all neighbours).
ExtendableFEMBase.continuify
— Functionfunction continuify(
source::FEVectorBlock,
operator = Identity;
abs::Bool = false,
broken = false,
order = "auto",
factor = 1,
regions::Array{Int,1} = [0]) where {T,Tv,Ti,FEType,APT}
Interpolate the evaluation of an operator applied to a finite element function onto a continuous Lagrange finite element space (H1Pk
), returning a new FEVector
with the interpolated values.
This function performs nodal interpolation of the (possibly vector-valued) result of applying operator
to the FE function represented by source
. The result is a new FE function in a continuous Lagrange space of the specified order and dimension. If broken = true
, the interpolation is performed in a piecewise (discontinuous) fashion.
Arguments
source
: TheFEVectorBlock
containing the coefficients of the original FE function.operator
: The function operator to apply before interpolation (default:Identity
).
Keyword Arguments
abs
: Iftrue
, interpolate the Euclidean norm of the result (default:false
).broken
: Iftrue
, generate a piecewise (discontinuous) interpolation (default:false
).order
: Polynomial order of the target Lagrange space (default:"auto"
, which chooses an appropriate order).factor
: Scaling factor applied to the result before interpolation (default:1
).regions
: List of region indices to restrict interpolation (default: all regions).
Returns
- A new
FEVector
in a continuous (or broken) Lagrange space, containing the interpolated values.
Nodal Evaluations
Plotting routines require nodal values, i.e., the values of a finite element function at the mesh nodes. The generic nodevalues!
function evaluates any finite element function at the grid nodes, averaging values if the function is discontinuous. For H1-conforming finite elements and identity evaluations, the nodevalues_view
function can provide a direct view into the coefficient field, avoiding unnecessary memory allocations.
ExtendableFEMBase.nodevalues!
— Functionfunction nodevalues!(
target::AbstractArray{<:Real,2},
source::AbstractArray{T,1},
FE::FESpace{Tv,Ti,FEType,AT},
operator::Type{<:AbstractFunctionOperator} = Identity;
abs::Bool = false,
factor = 1,
regions::Array{Int,1} = [0],
target_offset::Int = 0,
source_offset::Int = 0,
zero_target::Bool = true,
continuous::Bool = false
)
Evaluate a finite element function (given by the coefficient vector source
and FE space FE
) at all nodes of the grid, applying an optional function operator, and write the results into target
.
By default, the function is evaluated at every node in the mesh. If continuous
is false
, the value at each node is averaged over all neighboring cells (suitable for discontinuous quantities). If continuous
is true
, the value is taken from a single cell per node (suitable for continuous quantities). The result can optionally be scaled, offset, or restricted to specific regions.
Arguments
target
: Output array to store the evaluated values (size: result dimension × number of nodes).source
: Coefficient vector for the FE function.FE
: The finite element space.operator
: Function operator to apply at each node (default:Identity
).
Keyword Arguments
abs
: Iftrue
, store the Euclidean norm of the result at each node (default:false
).factor
: Scaling factor applied to the result (default:1
).regions
: List of region indices to restrict evaluation (default: all regions).target_offset
: Offset for writing intotarget
(default:0
).source_offset
: Offset for reading fromsource
(default:0
).zero_target
: Iftrue
, zero outtarget
before writing (default:true
).continuous
: Iftrue
, evaluate only once per node; otherwise, average over all neighboring cells (default:false
).
Notes
- The result dimension is determined by the FE space, the operator, and the
abs
argument. - The function modifies
target
in-place. - For vector-valued or higher-dimensional results, the first dimension of
target
corresponds to the result dimension.
Convenience method: Evaluate a finite element function (given by the FEVectorBlock), applying an optional function operator, and write the results into target
.
This forwards to the main nodevalues! method using a view of the block's entries and its FESpace.
ExtendableFEMBase.nodevalues
— Function nodevalues(
source::FEVectorBlock,
operator::Type{<:AbstractFunctionOperator} = Identity;
continuous = "auto",
nodes = [],
cellwise = false,
abs = false,
kwargs...
)
Evaluate a finite element function (given by the coefficient vector source
) at nodes of the grid, applying an optional function operator, and return the result as a newly allocated array of the appropriate size.
This function provides a flexible interface for extracting nodal or cellwise values from a finite element solution. By default, it evaluates at all nodes, but a subset of nodes can be specified. The result can be returned in a cellwise (piecewise) layout if desired.
Arguments
source
: TheFEVectorBlock
containing the coefficients of the FE function.operator
: The function operator to apply at each node (default:Identity
).
Keyword Arguments
continuous
: If"auto"
, automatically choose continuous/discontinuous evaluation based on FE type and operator. Iftrue
, evaluate only once per node; iffalse
, average over all neighboring cells. Is ignored whencellwise
istrue
.nodes
: List of node indices to evaluate (default: all nodes).cellwise
: Iftrue
, return values in a cellwise (piecewise) layout (default:false
).abs
: Iftrue
, return the Euclidean norm at each node (default:false
).kwargs...
: Additional keyword arguments passed to lower-level routines (e.g.,regions
,factor
,target_offset
,zero_target
, etc.).
Returns
- A newly allocated array containing the evaluated values, with shape depending on the options chosen.
Notes
- If
nodes
is empty, all nodes are evaluated. - If
cellwise
istrue
, the result is organized per cell (suitable for discontinuous or element-wise quantities). - The result dimension is determined by the FE space, the operator, and the
abs
argument.
ExtendableFEMBase.nodevalues_view
— Functionfunction nodevalues_view(
source::FEVectorBlock,
operator::Type{<:AbstractFunctionOperator} = Identity)
Return a vector of views into the nodal values of the given finite element function, allowing direct access to the underlying coefficient storage for each component.
This function provides efficient, zero-copy access to the nodal values of an FEVectorBlock
for unbroken H1-conforming finite element spaces with the identity operator. Each entry in the returned vector is a view into the coefficients corresponding to one component of the FE function, optionally restricted to a subset of nodes.
Arguments
source
: TheFEVectorBlock
containing the coefficients of the FE function.operator
: The function operator to apply (must beIdentity
for direct views; default:Identity
).
Keyword Arguments
nodes
: List of node indices to view (default: all nodes).
Returns
- A vector of
SubArray
views, one for each component, directly referencing the coefficients for the specified nodes.
Notes
- Only available for unbroken H1-conforming elements and the Identity operator.
ExtendableFEMBase.nodevalues_subset!
— Functionfunction nodevalues_subset!(
target::AbstractArray{<:Real,2},
source::AbstractArray{T,1},
FE::FESpace{Tv,Ti,FEType,AT},
operator::Type{<:AbstractFunctionOperator} = Identity;
regions::Array{Int,1} = [0],
abs::Bool = false,
factor = 1,
nodes = [],
target_offset::Int = 0, # start to write into target after offset
zero_target::Bool = true, # target vector is zeroed
continuous::Bool = false)
Evaluate (an operator of) a finite element function (given by the coefficient vector source
and FE space FE
) at a specified subset of nodes, and write the results into target
.
For each node in nodes
, the function is evaluated (optionally with operator
) and the result is written to the corresponding column of target
. If continuous
is false
, values are averaged over all neighboring cells; if true
, only one cell is used per node. If abs
is true
, the Euclidean norm is computed instead of the raw values.
Arguments
target
: Output array to store the evaluated values (size: result dimension × number of nodes).source
: Coefficient vector for the FE function.FE
: The finite element space.operator
: Function operator to apply (default:Identity
).abs
: Iftrue
, store the Euclidean norm of the result at each node (default:false
).factor
: Scaling factor applied to the result (default:1
).nodes
: List of node indices to evaluate (default: all nodes).regions
: List of region indices to restrict evaluation (default: all regions).target_offset
: Offset for writing intotarget
(default:0
).source_offset
: Offset for reading fromsource
(default:0
).zero_target
: Iftrue
, zero outtarget
before writing (default:true
).continuous
: Iftrue
, evaluate only once per node; otherwise, average over all neighboring cells (default:false
).
Notes
- The result dimension is determined by the FE space and operator and the
abs
argument.
Displace Mesh
Nodal values (e.g. of a FEVector that discretizes a displacement) can be used to displace the mesh.
ExtendableFEMBase.displace_mesh!
— Functionfunction displace_mesh!(xgrid::ExtendableGrid, source::FEVectorBlock; magnify = 1)
Displace all nodes of the given grid by adding a vector-valued finite element field as a displacement, scaled by an optional magnification factor.
This function modifies the coordinates of xgrid
in-place by adding the nodal values of the FE function represented by source
(typically a displacement field) to each node. The displacement can be scaled by the magnify
parameter. After the update, all cached geometric quantities in the grid are invalidated and recomputed.
Arguments
xgrid
: TheExtendableGrid
whose node coordinates will be updated.source
: AnFEVectorBlock
representing the displacement field (must be vector-valued and defined on the grid).
Keyword Arguments
magnify
: Scaling factor for the displacement field (default:1
).
ExtendableFEMBase.displace_mesh
— Functionfunction displace_mesh(xgrid::ExtendableGrid, source::FEVectorBlock; magnify = 1)
Return a new grid with node coordinates displaced by a vector-valued finite element field, optionally scaled by a magnification factor.
Arguments
xgrid
: TheExtendableGrid
to be copied and displaced.source
: AnFEVectorBlock
representing the displacement field (must be vector-valued and defined on the grid).
Keyword Arguments
magnify
: Scaling factor for the displacement field (default:1
).
Returns
- A new
ExtendableGrid
with displaced node coordinates.