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 childTypeDescription
qpinfo.xVector{Real}space coordinates of quadrature point
qpinfo.timeRealcurrent time
qpinfo.itemIntegercurrent item that contains qpinfo.x
qpinfo.regionIntegerregion number of item
qpinfo.xrefVector{Real}reference coordinates within item of qpinfo.x
qpinfo.volumeRealvolume of item
qpinfo.paramsVector{Any}parameters that can be transfered via keyword arguments

Standard Interpolations

Each finite element has its standard interpolator that can be applied to some user-defined DataFunction. Instead of interpolating on the full cells, the interpolation can be restricted to faces or edges.

It is also possible to interpolate finite element functions on one grid onto a finite element function on another grid (experimental feature, does not work for all finite elements yet and shall be extended to interpolations of operator evaluations as well in future).

ExtendableGrids.interpolate!Function
function ExtendableGrids.interpolate!(target::FEVectorBlock,
	 AT::Type{<:AssemblyType},
	 source!::Function;
	 items = [],
	 bonus_quadorder = 0,
	 time = 0,
	 kwargs...)

Interpolates the given source into the finite elements space assigned to the target FEVectorBlock with the specified AssemblyType (usualy ON_CELLS).

The source functions should adhere to the interface

	source!(result, qpinfo)

The qpinfo argument communicates vast information of the current quadrature/evaluation point.

The bonus_quadorder argument can be used to steer the quadrature order of integrals that needs to be computed for the interpolation (the default quadrature order corressponds to the polynomial order of the finite element).

source
function ExtendableGrids.interpolate!(target::FEVectorBlock,
	 source::Function;
	 items = [],
	 bonus_quadorder = 0,
	 time = 0,
	 kwargs...)

Interpolates the given source function into the finite element space assigned to the target FEVectorBlock.

The source functions should adhere to the interface

	source!(result, qpinfo)

The qpinfo argument communicates vast information of the current quadrature/evaluation point.

The bonus_quadorder argument can be used to steer the quadrature order of integrals that needs to be computed for the interpolation (the default quadrature order corressponds to the polynomial order of the finite element).

source

Nodal Evaluations

Usually, Plotters need nodal values, so there is a gengeric function that evaluates any finite element function at the nodes of the grids (possibly by averaging if discontinuous). In case of Identity evaluations of an H1-conforming finite element, the function nodevalues_view can generate a view into the coefficient field that avoids further allocations.

ExtendableFEMBase.nodevalues!Function
function nodevalues!(
	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,
	target_offset::Int = 0,   # start to write into target after offset
	zero_target::Bool = true, # target vector is zeroed
	continuous::Bool = false)

Evaluates the finite element function with the coefficient vector source (interpreted as a coefficient vector for the FESpace FE) and the specified FunctionOperator at all the nodes of the (specified regions of the) grid and writes the values into target. Discontinuous (continuous = false) quantities are evaluated in all neighbouring cells of each node and then averaged. Continuous (continuous = true) quantities are only evaluated once at each node.

source
function nodevalues!(
	target::AbstractArray{<:Real,2},
	source::FEVectorBlock,
	operator::Type{<:AbstractFunctionOperator} = Identity;
	regions::Array{Int,1} = [0],
	abs::Bool = false,
	factor = 1,
	cellwise = false,		  # return cellwise nodevalues ncells x nnodes_on_cell
	target_offset::Int = 0,   # start to write into target after offset
	zero_target::Bool = true, # target vector is zeroed
	continuous::Bool = false)

Evaluates the finite element function with the coefficient vector source and the specified FunctionOperator at all the nodes of the (specified regions of the) grid and writes the values into target. Discontinuous (continuous = false) quantities are evaluated in all neighbouring cells of each node and then averaged. Continuous (continuous = true) quantities are only evaluated once at each node.

source
ExtendableFEMBase.nodevaluesFunction
function nodevalues(
	source::FEVectorBlock,
	operator::Type{<:AbstractFunctionOperator} = Identity;
	regions::Array{Int,1} = [0],
	abs::Bool = false,
	factor = 1,
	nodes = [],				  
	cellwise = false,		  # return cellwise nodevalues ncells x nnodes_on_cell (only if nodes == [])
	target_offset::Int = 0,   # start to write into target after offset
	zero_target::Bool = true, # target vector is zeroed
	continuous::Bool = false)

Evaluates the finite element function with the coefficient vector source and the specified FunctionOperator at the specified list of nodes of the grid (default = all nodes) and writes the values in that order into target. Nodes that are not part of the specified regions (default = all regions) are set to zero. Discontinuous (continuous = false) quantities are evaluated in all neighbouring cells of each node and then averaged. Continuous (continuous = true) quantities are only evaluated once at each node.

source
ExtendableFEMBase.nodevalues_viewFunction
function nodevalues_view(
	source::FEVectorBlock,
	operator::Type{<:AbstractFunctionOperator} = Identity)

Returns a vector of views of the nodal values of the source block (currently works for unbroken H1-conforming elements) that directly accesses the coefficients.

source
ExtendableFEMBase.nodevalues_subset!Function
function 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)

Evaluates the finite element function with the coefficient vector source (interpreted as a coefficient vector for the FESpace FE) and the specified FunctionOperator at the specified list of nodes of the grid and writes the values in that order into target. Node values for nodes that are not part of the specified regions (default = all regions) are set to zero. Discontinuous (continuous = false) quantities are evaluated in all neighbouring cells (in the specified regions) of each node and then averaged. Continuous (continuous = true) quantities are only evaluated once at each node.

source

Lazy Interpolation

To interpolate between different finite element spaces and meshes, there is a lazy interpolation routine that works in all cases (but is not very efficient as it involves a PointeEvaluator and CellFinder):

ExtendableFEMBase.lazy_interpolate!Function
function lazy_interpolate!(
	target::FEVectorBlock{T1,Tv,Ti},
	source::FEVectorBlock{T2,Tv,Ti};
	operator = Identity,
	postprocess = nothing,
	xtrafo = nothing,
	items = [],
	not_in_domain_value = 1e30,
	eps = 1e-13,
	use_cellparents::Bool = false) where {T1,T2,Tv,Ti}

Interpolates (operator-evaluations of) the given finite element function into the finite element space assigned to the target FEVectorBlock. (Currently not the most efficient way as it is based on the PointEvaluation pattern and cell search. If CellParents are available in the grid components of the target grid, these parent cell information can be used to improve the search. To activate this put 'use_cellparents' = true). By some given kernel function that is conforming to the interface

kernel!(result, input, qpinfo)

the operator evaluation (=input) can be further postprocessed. The qpinfo argument allows to access information at the current quadrature point.

Note: discontinuous quantities at vertices of the target grid will be evaluted in the first found cell of the source grid. No averaging is performed. With eps the tolerances of the cell search via ExtendableGrids.CellFinder can be steered.

source

Displace Mesh

Nodal values (e.g. of a FEVector that discretizes a displacement) can be used to displace the mesh.

ExtendableFEMBase.displace_mesh!Function
function displace_mesh!(xgrid::ExtendableGrid, source::FEVectorBlock; magnify = 1)

Moves all nodes of the grid by adding the displacement field in source (expects a vector-valued finite element) times a magnify value.

source
ExtendableFEMBase.displace_meshFunction
function displace_mesh(xgrid::ExtendableGrid, source::FEVectorBlock; magnify = 1)

Returns a new grid by adding the displacement field in source (expects a vector-valued finite element) to the coordinates of the provided xgrid times a magnify value.

source