Ferrite.PointValuesType
PointScalarValues(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.AbstractRefShapeType

Represents a reference shape which quadrature rules and interpolations are defined on. Currently, the only concrete types that subtype this type are RefCube in 1, 2 and 3 dimensions, and RefTetrahedron in 2 and 3 dimensions.

Ferrite.AffineConstraintType
AffineConstraint(constrained_dof::Int, entries::Vector{Pair{Int,T}}, b::T) where T

Define an affine/linear constraint to constrain one degree of freedom, u[i], such that u[i] = ∑(u[j] * a[j]) + b, where i=constrained_dof and each element in entries are j => a[j]

Ferrite.BCValuesType
BCValues(func_interpol::Interpolation, geom_interpol::Interpolation, boundary_type::Union{Type{<:BoundaryIndex}})

BCValues stores the shape values at all faces/edges/vertices (depending on boundary_type) for the geomatric interpolation (geom_interpol), for each dof-position determined by the func_interpol. Used mainly by the ConstrainHandler.

Ferrite.CellType
Cell{dim,N,M} <: AbstractCell{dim,N,M}

A Cell is a subdomain defined by a collection of Nodes. The parameter dim refers here to the geometrical/ambient dimension, i.e. the dimension of the nodes in the grid and not the topological dimension of the cell (i.e. the dimension of the reference element obtained by default_interpolation). A Cell has N nodes and M faces. Note that a Cell is not defined geometrically by node coordinates, but rather topologically by node indices into the node vector of some grid.

Fields

  • nodes::Ntuple{N,Int}: N-tuple that stores the node ids. The ordering defines a cell's and its subentities' orientations.
Ferrite.CellCacheType
CellCache(grid::Grid)
CellCache(dh::AbstractDofHandler)

Create a cache object with pre-allocated memory for the nodes, coordinates, and dofs of a cell. The cache is updated for a new cell by calling reinit!(cache, cellid) where cellid::Int is the cell id.

Struct fields of CellCache

  • cc.nodes :: Vector{Int}: global node ids
  • cc.coords :: Vector{<:Vec}: node coordinates
  • cc.dofs :: Vector{Int}: global dof ids (empty when constructing the cache from a grid)

Methods with CellCache

  • reinit!(cc, i): reinitialize the cache for cell i
  • cellid(cc): get the cell id of the currently cached cell
  • getnodes(cc): get the global node ids of the cell
  • getcoordinates(cc): get the coordinates of the cell
  • celldofs(cc): get the global dof ids of the cell
  • reinit!(fev, cc): reinitialize CellValues or FaceValues

See also CellIterator.

Ferrite.CellIndexType

A CellIndex wraps an Int and corresponds to a cell with that number in the mesh

Ferrite.CellIteratorType
CellIterator(grid::Grid, cellset=1:getncells(grid))
CellIterator(dh::AbstractDofHandler, cellset=1:getncells(dh))

Create a CellIterator to conveniently iterate over all, or a subset, of the cells in a grid. The elements of the iterator are CellCaches which are properly reinit!ialized. See CellCache for more details.

Looping over a CellIterator, i.e.:

for cc in CellIterator(grid, cellset)
    # ...
end

is thus simply convenience for the following equivalent snippet:

cc = CellCache(grid)
for idx in cellset
    reinit!(cc, idx)
    # ...
end
Ferrite.CellScalarValuesType
CellScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
CellVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A CellValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. in the finite element cell. There are two different types of CellValues: CellScalarValues and CellVectorValues. As the names suggest, CellScalarValues utilizes scalar shape functions and CellVectorValues utilizes vectorial shape functions. For a scalar field, the CellScalarValues type should be used. For vector field, both subtypes can be used.

Arguments:

  • T: an optional argument (default to Float64) to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of a Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.CellValuesType
CellScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
CellVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A CellValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. in the finite element cell. There are two different types of CellValues: CellScalarValues and CellVectorValues. As the names suggest, CellScalarValues utilizes scalar shape functions and CellVectorValues utilizes vectorial shape functions. For a scalar field, the CellScalarValues type should be used. For vector field, both subtypes can be used.

Arguments:

  • T: an optional argument (default to Float64) to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of a Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.CellVectorValuesType
CellScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
CellVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A CellValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. in the finite element cell. There are two different types of CellValues: CellScalarValues and CellVectorValues. As the names suggest, CellScalarValues utilizes scalar shape functions and CellVectorValues utilizes vectorial shape functions. For a scalar field, the CellScalarValues type should be used. For vector field, both subtypes can be used.

Arguments:

  • T: an optional argument (default to Float64) to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of a Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.CrouzeixRaviartType

Classical non-conforming Crouzeix–Raviart element.

For details we refer to the original paper: M. Crouzeix and P. Raviart. "Conforming and nonconforming finite element methods for solving the stationary Stokes equations I." ESAIM: Mathematical Modelling and Numerical Analysis-Modélisation Mathématique et Analyse Numérique 7.R3 (1973): 33-75.

Ferrite.DirichletType
Dirichlet(u::Symbol, ∂Ω::Set, f::Function, components=nothing)

Create a Dirichlet boundary condition on u on the ∂Ω part of the boundary. f is a function of the form f(x) or f(x, t) where x is the spatial coordinate and t is the current time, and returns the prescribed value. components specify the components of u that are prescribed by this condition. By default all components of u are prescribed.

For example, here we create a Dirichlet condition for the :u field, on the faceset called ∂Ω and the value given by the sin function:

Examples

# Obtain the faceset from the grid
∂Ω = getfaceset(grid, "boundary-1")

# Prescribe scalar field :s on ∂Ω to sin(t)
dbc = Dirichlet(:s, ∂Ω, (x, t) -> sin(t))

# Prescribe all components of vector field :v on ∂Ω to 0
dbc = Dirichlet(:v, ∂Ω, x -> 0 * x)

# Prescribe component 2 and 3 of vector field :v on ∂Ω to [sin(t), cos(t)]
dbc = Dirichlet(:v, ∂Ω, (x, t) -> [sin(t), cos(t)], [2, 3])

Dirichlet boundary conditions are added to a ConstraintHandler which applies the condition via apply! and/or apply_zero!.

Ferrite.DofHandlerType
DofHandler(grid::Grid)

Construct a DofHandler based on grid.

Operates slightly faster than MixedDofHandler. Supports:

  • Grids with a single concrete cell type.
  • One or several fields on the whole domaine.
Ferrite.DofOrder.ComponentWiseType
DofOrder.ComponentWise()
DofOrder.ComponentWise(target_blocks::Vector{Int})

Dof order passed to renumber! to renumber global dofs component wise resulting in a globally blocked system.

The default behavior is to group dofs of each component into their own block, with the same order as in the DofHandler. This can be customized by passing a vector of length ncomponents that maps each component to a "target block" (see DofOrder.FieldWise for details).

This renumbering is stable such that the original relative ordering of dofs within each target block is maintained.

Ferrite.DofOrder.FieldWiseType
DofOrder.FieldWise()
DofOrder.FieldWise(target_blocks::Vector{Int})

Dof order passed to renumber! to renumber global dofs field wise resulting in a globally blocked system.

The default behavior is to group dofs of each field into their own block, with the same order as in the DofHandler. This can be customized by passing a vector of the same length as the total number of fields in the DofHandler (see getfieldnames(dh)) that maps each field to a "target block": to renumber a DofHandler with three fields :u, :v, :w such that dofs for :u and :w end up in the first global block, and dofs for :v in the second global block use DofOrder.FieldWise([1, 2, 1]).

This renumbering is stable such that the original relative ordering of dofs within each target block is maintained.

Ferrite.EdgeIndexType

A EdgeIndex wraps an (Int, Int) and defines a local edge by pointing to a (cell, edge).

Ferrite.ExclusiveTopologyType
ExclusiveTopology(cells::Vector{C}) where C <: AbstractCell

ExclusiveTopology saves topological (connectivity) data of the grid. The constructor works with an AbstractCell vector for all cells that dispatch vertices, faces and in 3D edges as well as the utility functions face_npoints and edge_npoints. The struct saves the highest dimensional neighborhood, i.e. if something is connected by a face and an edge only the face neighborhood is saved. The lower dimensional neighborhood is recomputed, if needed.

Fields

  • vertex_to_cell::Dict{Int,Vector{Int}}: global vertex id to all cells containing the vertex
  • cell_neighbor::Vector{EntityNeighborhood{CellIndex}}: cellid to all connected cells
  • face_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}: face_neighbor[cellid,local_face_id] -> neighboring face
  • vertex_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}: vertex_neighbor[cellid,local_vertex_id] -> neighboring vertex
  • edge_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}: edge_neighbor[cellid_local_vertex_id] -> neighboring edge
  • vertex_vertex_neighbor::Dict{Int,EntityNeighborhood{VertexIndex}}: global vertex id -> all connected vertices by edge or face
  • face_skeleton::Vector{FaceIndex}: list of unique faces in the grid
Ferrite.FaceIndexType

A FaceIndex wraps an (Int, Int) and defines a local face by pointing to a (cell, face).

Ferrite.FaceScalarValuesType
FaceScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
FaceVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A FaceValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. on the faces of finite elements. There are two different types of FaceValues: FaceScalarValues and FaceVectorValues. As the names suggest, FaceScalarValues utilizes scalar shape functions and FaceVectorValues utilizes vectorial shape functions. For a scalar field, the FaceScalarValues type should be used. For vector field, both subtypes can be used.

Note

The quadrature rule for the face should be given with one dimension lower. I.e. for a 3D case, the quadrature rule should be in 2D.

Arguments:

  • T: an optional argument to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of an Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.FaceValuesType
FaceScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
FaceVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A FaceValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. on the faces of finite elements. There are two different types of FaceValues: FaceScalarValues and FaceVectorValues. As the names suggest, FaceScalarValues utilizes scalar shape functions and FaceVectorValues utilizes vectorial shape functions. For a scalar field, the FaceScalarValues type should be used. For vector field, both subtypes can be used.

Note

The quadrature rule for the face should be given with one dimension lower. I.e. for a 3D case, the quadrature rule should be in 2D.

Arguments:

  • T: an optional argument to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of an Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.FaceVectorValuesType
FaceScalarValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])
FaceVectorValues([::Type{T}], quad_rule::QuadratureRule, func_interpol::Interpolation, [geom_interpol::Interpolation])

A FaceValues object facilitates the process of evaluating values of shape functions, gradients of shape functions, values of nodal functions, gradients and divergences of nodal functions etc. on the faces of finite elements. There are two different types of FaceValues: FaceScalarValues and FaceVectorValues. As the names suggest, FaceScalarValues utilizes scalar shape functions and FaceVectorValues utilizes vectorial shape functions. For a scalar field, the FaceScalarValues type should be used. For vector field, both subtypes can be used.

Note

The quadrature rule for the face should be given with one dimension lower. I.e. for a 3D case, the quadrature rule should be in 2D.

Arguments:

  • T: an optional argument to determine the type the internal data is stored as.
  • quad_rule: an instance of a QuadratureRule
  • func_interpol: an instance of an Interpolation used to interpolate the approximated function
  • geom_interpol: an optional instance of an Interpolation which is used to interpolate the geometry

Common methods:

Ferrite.FieldType
Field(name::Symbol, interpolation::Interpolation, dim::Int)

Construct dim-dimensional Field called name which is approximated by interpolation.

The interpolation is used for distributing the degrees of freedom.

Ferrite.FieldHandlerType
FieldHandler(fields::Vector{Field}, cellset::Set{Int})

Construct a FieldHandler based on an array of Fields and assigns it a set of cells.

A FieldHandler must fulfill the following requirements:

  • All Cells in cellset are of the same type.
  • Each field only uses a single interpolation on the cellset.
  • Each cell belongs only to a single FieldHandler, i.e. all fields on a cell must be added within the same FieldHandler.

Notice that a FieldHandler can hold several fields.

Ferrite.GridType
Grid{dim, C<:AbstractCell, T<:Real} <: AbstractGrid}

A Grid is a collection of Cells and Nodes which covers the computational domain, together with Sets of cells, nodes and faces. There are multiple helper structures to apply boundary conditions or define subdomains. They are gathered in the cellsets, nodesets, facesets, edgesets and vertexsets.

Fields

  • cells::Vector{C}: stores all cells of the grid
  • nodes::Vector{Node{dim,T}}: stores the dim dimensional nodes of the grid
  • cellsets::Dict{String,Set{Int}}: maps a String key to a Set of cell ids
  • nodesets::Dict{String,Set{Int}}: maps a String key to a Set of global node ids
  • facesets::Dict{String,Set{FaceIndex}}: maps a String to a Set of Set{FaceIndex} (global_cell_id, local_face_id)
  • edgesets::Dict{String,Set{EdgeIndex}}: maps a String to a Set of Set{EdgeIndex} (global_cell_id, local_edge_id
  • vertexsets::Dict{String,Set{VertexIndex}}: maps a String key to a Set of local vertex ids
  • boundary_matrix::SparseMatrixCSC{Bool,Int}: optional, only needed by onboundary to check if a cell is on the boundary, see, e.g. Helmholtz example
Ferrite.InterpolationType
Interpolation{ref_dim, ref_shape, order}()

Return an Interpolation on a ref_dim-dimensional reference shape (see AbstractRefShape) ref_shape and order order. order corresponds to the order of the interpolation. The interpolation is used to define shape functions to interpolate a function between nodes.

The following interpolations are implemented:

  • Lagrange{1,RefCube,1}
  • Lagrange{1,RefCube,2}
  • Lagrange{2,RefCube,1}
  • Lagrange{2,RefCube,2}
  • Lagrange{2,RefTetrahedron,1}
  • Lagrange{2,RefTetrahedron,2}
  • Lagrange{2,RefTetrahedron,3}
  • Lagrange{2,RefTetrahedron,4}
  • Lagrange{2,RefTetrahedron,5}
  • BubbleEnrichedLagrange{2,RefTetrahedron,1}
  • CrouzeixRaviart{2,1}
  • Lagrange{3,RefCube,1}
  • Lagrange{3,RefCube,2}
  • Lagrange{3,RefTetrahedron,1}
  • Lagrange{3,RefTetrahedron,2}
  • Lagrange{3,RefPrism,1}
  • Lagrange{3,RefPrism,2}
  • Serendipity{2,RefCube,2}
  • Serendipity{3,RefCube,2}

Examples

julia> ip = Lagrange{2,RefTetrahedron,2}()
Ferrite.Lagrange{2,Ferrite.RefTetrahedron,2}()

julia> getnbasefunctions(ip)
6
Ferrite.InterpolationInfoType
InterpolationInfo

Gathers all the information needed to distribute dofs for a given interpolation. Note that this cache is of the same type no matter the interpolation: the purpose is to make dof-distribution type-stable.

Ferrite.L2ProjectorMethod
L2Projector(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 with func_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.MixedDofHandlerType
MixedDofHandler(grid::Grid)

Construct a MixedDofHandler based on grid. Supports:

  • Grids with or without concrete element type (E.g. "mixed" grids with several different element types.)
  • One or several fields, which can live on the whole domain or on subsets of the Grid.
Ferrite.NodeType
Node{dim, T}

A Node is a point in space.

Fields

  • x::Vec{dim,T}: stores the coordinates
Ferrite.PeriodicDirichletType
PeriodicDirichlet(u::Symbol, face_mapping, components=nothing)
PeriodicDirichlet(u::Symbol, face_mapping, R::AbstractMatrix, components=nothing)
PeriodicDirichlet(u::Symbol, face_mapping, f::Function, components=nothing)

Create a periodic Dirichlet boundary condition for the field u on the face-pairs given in face_mapping. The mapping can be computed with collect_periodic_faces. The constraint ensures that degrees-of-freedom on the mirror face are constrained to the corresponding degrees-of-freedom on the image face. components specify the components of u that are prescribed by this condition. By default all components of u are prescribed.

If the mapping is not aligned with the coordinate axis (e.g. rotated) a rotation matrix R should be passed to the constructor. This matrix rotates dofs on the mirror face to the image face. Note that this is only applicable for vector-valued problems.

To construct an inhomogeneous periodic constraint it is possible to pass a function f. Note that this is currently only supported when the periodicity is aligned with the coordinate axes.

See the manual section on Periodic boundary conditions for more information.

Ferrite.PointEvalHandlerType
PointEvalHandler(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, with nothing 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, with nothing 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) a dh::DofHandler + uh::Vector (for example the FE-solution), or ii) a p::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.PointIteratorType
PointIterator(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.PointLocationType
PointLocation

Element of a PointIterator, typically used to reinitialize PointValues. Fields:

  • cid::Int: ID of the cell containing the point
  • local_coord::Vec: the local (reference) coordinate of the point
  • coords::Vector{Vec}: the coordinates of the cell
Ferrite.QuadratureRuleType
QuadratureRule{dim,shape}([quad_rule_type::Symbol], order::Int)

Create a QuadratureRule used for integration. dim is the space dimension, shape an AbstractRefShape and order the order of the quadrature rule. quad_rule_type is an optional argument determining the type of quadrature rule, currently the :legendre and :lobatto rules are implemented.

A QuadratureRule is used to approximate an integral on a domain by a weighted sum of function values at specific points:

$\int\limits_\Omega f(\mathbf{x}) \text{d} \Omega \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) w_q$

The quadrature rule consists of $n_q$ points in space $\mathbf{x}_q$ with corresponding weights $w_q$.

In Ferrite, the QuadratureRule type is mostly used as one of the components to create a CellValues or FaceValues object.

Common methods:

Example:

julia> QuadratureRule{2, RefTetrahedron}(1)
Ferrite.QuadratureRule{2,Ferrite.RefTetrahedron,Float64}([0.5], Tensors.Tensor{1,2,Float64,2}[[0.333333, 0.333333]])

julia> QuadratureRule{1, RefCube}(:lobatto, 2)
Ferrite.QuadratureRule{1,Ferrite.RefCube,Float64}([1.0, 1.0], Tensors.Tensor{1,1,Float64,1}[[-1.0], [1.0]])
Ferrite.RHSDataType
RHSData

Stores the constrained columns and mean of the diagonal of stiffness matrix A.

Ferrite.ValuesType

Abstract type which has CellValues and FaceValues as subtypes

Ferrite.VertexIndexType

A VertexIndex wraps an (Int, Int) and defines a local vertex by pointing to a (cell, vert).

Ferrite.add!Function
add!(dh::AbstractDofHandler, name::Symbol, dim::Int[, ip::Interpolation])

Add a dim-dimensional Field called name which is approximated by ip to dh.

The field is added to all cells of the underlying grid. In case no interpolation ip is given, the default interpolation of the grid's celltype is used. If the grid uses several celltypes, add!(dh::MixedDofHandler, fh::FieldHandler) must be used instead.

Ferrite.add!Method
add!(ch::ConstraintHandler, ac::AffineConstraint)

Add the AffineConstraint to the ConstraintHandler.

Ferrite.add!Method
add!(ch::ConstraintHandler, dbc::Dirichlet)

Add a Dirichlet boundary condition to the ConstraintHandler.

Ferrite.add_prescribed_dof!Function
add_prescribed_dof!(ch, constrained_dof::Int, inhomogeneity, dofcoefficients=nothing)

Add a constrained dof directly to the ConstraintHandler. This function checks if the constrained_dof is already constrained, and overrides the old constraint if true.

Ferrite.addcellset!Method
addcellset!(grid::AbstractGrid, name::String, cellid::Union{Set{Int}, Vector{Int}})
addcellset!(grid::AbstractGrid, name::String, f::function; all::Bool=true)

Adds a cellset to the grid with key name. Cellsets are typically used to define subdomains of the problem, e.g. two materials in the computational domain. The MixedDofHandler can construct different fields which live not on the whole domain, but rather on a cellset. all=true implies that f(x) must return true for all nodal coordinates x in the cell if the cell should be added to the set, otherwise it suffices that f(x) returns true for one node.

addcellset!(grid, "left", Set((1,3))) #add cells with id 1 and 3 to cellset left
addcellset!(grid, "right", x -> norm(x[1]) < 2.0 ) #add cell to cellset right, if x[1] of each cell's node is smaller than 2.0
Ferrite.addfaceset!Method
addfaceset!(grid::AbstractGrid, name::String, faceid::Union{Set{FaceIndex},Vector{FaceIndex}})
addfaceset!(grid::AbstractGrid, name::String, f::Function; all::Bool=true)

Adds a faceset to the grid with key name. A faceset maps a String key to a Set of tuples corresponding to (global_cell_id, local_face_id). Facesets are used to initialize Dirichlet structs, that are needed to specify the boundary for the ConstraintHandler. all=true implies that f(x) must return true for all nodal coordinates x on the face if the face should be added to the set, otherwise it suffices that f(x) returns true for one node.

addfaceset!(grid, "right", Set(((2,2),(4,2))) #see grid manual example for reference
addfaceset!(grid, "clamped", x -> norm(x[1]) ≈ 0.0) #see incompressible elasticity example for reference
Ferrite.addindex!Function
addindex!(A::AbstractMatrix{T}, v::T, i::Int, j::Int)
addindex!(b::AbstractVector{T}, v::T, i::Int)

Equivalent to A[i, j] += v but more efficient.

A[i, j] += v is lowered to A[i, j] = A[i, j] + v which requires a double lookup of the memory location for index (i, j) – one time for the read, and one time for the write. This method avoids the double lookup.

Zeros are ignored (i.e. if iszero(v)) by returning early. If the index (i, j) is not existing in the sparsity pattern of A this method throws a SparsityError.

Fallback: A[i, j] += v.

Ferrite.addnodeset!Method
addnodeset!(grid::AbstractGrid, name::String, nodeid::Union{Vector{Int},Set{Int}})
addnodeset!(grid::AbstractGrid, name::String, f::Function)

Adds a nodeset::Dict{String, Set{Int}} to the grid with key name. Has the same interface as addcellset. However, instead of mapping a cell id to the String key, a set of node ids is returned.

Ferrite.apply!Function
apply!(K::SparseMatrixCSC, rhs::AbstractVector, ch::ConstraintHandler)

Adjust the matrix K and right hand side rhs to account for the Dirichlet boundary conditions specified in ch such that K \ rhs gives the expected solution.

Note

apply!(K, rhs, ch) essentially calculates

rhs[free_dofs] = rhs[free_dofs] - K[free_dofs, constrained_dofs] * a[constrained]

where a[constrained] are the inhomogeneities. Consequently, the sign of rhs matters (in contrast to for apply_zero!).

apply!(v::AbstractVector, ch::ConstraintHandler)

Apply Dirichlet boundary conditions and affine constraints, specified in ch, to the solution vector v.

Examples

K, f = assemble_system(...) # Assemble system
apply!(K, f, ch)            # Adjust K and f to account for boundary conditions
u = K \ f                   # Solve the system, u should be "approximately correct"
apply!(u, ch)               # Explicitly make sure bcs are correct
Note

The last operation is not strictly necessary since the boundary conditions should already be fulfilled after apply!(K, f, ch). However, solvers of linear systems are not exact, and thus apply!(u, ch) can be used to make sure the boundary conditions are fulfilled exactly.

Ferrite.apply_analytical!Function
apply_analytical!(
    a::AbstractVector, dh::AbstractDofHandler, fieldname::Symbol, 
    f::Function, cellset=1:getncells(dh.grid))

Apply a solution f(x) by modifying the values in the degree of freedom vector a pertaining to the field fieldname for all cells in cellset. The function f(x) are given the spatial coordinate of the degree of freedom. For scalar fields, f(x)::Number, and for vector fields with dimension dim, f(x)::Vec{dim}.

This function can be used to apply initial conditions for time dependent problems.

Note

This function only works for standard nodal finite element interpolations when the function value at the (algebraic) node is equal to the corresponding degree of freedom value. This holds for e.g. Lagrange and Serendipity interpolations, including sub- and superparametric elements.

Ferrite.apply_assemble!Method
apply_assemble!(
    assembler::AbstractSparseAssembler, ch::ConstraintHandler,
    global_dofs::AbstractVector{Int},
    local_matrix::AbstractMatrix, local_vector::AbstractVector;
    apply_zero::Bool = false
)

Assemble local_matrix and local_vector into the global system in assembler by first doing constraint condensation using apply_local!.

This is similar to using apply_local! followed by assemble! with the advantage that non-local constraints can be handled, since this method can write to entries of the global matrix and vector outside of the indices in global_dofs.

When the keyword argument apply_zero is true all inhomogeneities are set to 0 (cf. apply! vs apply_zero!).

Note that this method is destructive since it modifies local_matrix and local_vector.

Ferrite.apply_local!Method
apply_local!(
    local_matrix::AbstractMatrix, local_vector::AbstractVector,
    global_dofs::AbstractVector, ch::ConstraintHandler;
    apply_zero::Bool = false
)

Similar to apply! but perform condensation of constrained degrees-of-freedom locally in local_matrix and local_vector before they are to be assembled into the global system.

When the keyword argument apply_zero is true all inhomogeneities are set to 0 (cf. apply! vs apply_zero!).

This method can only be used if all constraints are "local", i.e. no constraint couples with dofs outside of the element dofs (global_dofs) since condensation of such constraints requires writing to entries in the global matrix/vector. For such a case, apply_assemble! can be used instead.

Note that this method is destructive since it, by definition, modifies local_matrix and local_vector.

Ferrite.apply_rhs!Function
apply_rhs!(data::RHSData, f::AbstractVector, ch::ConstraintHandler, applyzero::Bool=false)

Applies the boundary condition to the right-hand-side vector without modifying the stiffness matrix.

See also: get_rhs_data.

Ferrite.apply_zero!Function
apply_zero!(K::SparseMatrixCSC, rhs::AbstractVector, ch::ConstraintHandler)

Adjust the matrix K and the right hand side rhs to account for prescribed Dirichlet boundary conditions and affine constraints such that du = K \ rhs gives the expected result (e.g. du zero for all prescribed degrees of freedom).

apply_zero!(v::AbstractVector, ch::ConstraintHandler)

Zero-out values in v corresponding to prescribed degrees of freedom and update values prescribed by affine constraints, such that if a fulfills the constraints, a ± v also will.

These methods are typically used in e.g. a Newton solver where the increment, du, should be prescribed to zero even for non-homogeneouos boundary conditions.

See also: apply!.

Examples

u = un + Δu                 # Current guess
K, g = assemble_system(...) # Assemble residual and tangent for current guess
apply_zero!(K, g, ch)       # Adjust tangent and residual to take prescribed values into account
ΔΔu = K \ g                # Compute the (negative) increment, prescribed values are "approximately" zero
apply_zero!(ΔΔu, ch)        # Make sure values are exactly zero
Δu .-= ΔΔu                  # Update current guess
Note

The last call to apply_zero! is only strictly necessary for affine constraints. However, even if the Dirichlet boundary conditions should be fulfilled after apply!(K, g, ch), solvers of linear systems are not exact. apply!(ΔΔu, ch) can be used to make sure the values for the prescribed degrees of freedom are fulfilled exactly.

Ferrite.assemble!Method
assemble!(A::AbstractSparseAssembler, dofs::AbstractVector{Int}, Ke::AbstractMatrix)
assemble!(A::AbstractSparseAssembler, dofs::AbstractVector{Int}, Ke::AbstractMatrix, fe::AbstractVector)

Assemble the element stiffness matrix Ke (and optional force vector fe) into the global stiffness (and force) in A, given the element degrees of freedom dofs.

This is equivalent to K[dofs, dofs] += Ke and f[dofs] += fe, where K is the global stiffness matrix and f the global force/residual vector, but more efficient.

Ferrite.assemble!Method
assemble!(g, dofs, ge)

Assembles the element residual ge into the global residual vector g.

Ferrite.assemble!Method
assemble!(a::Assembler, dofs, Ke)

Assembles the element matrix Ke into a.

Ferrite.assemble!Method
assemble!(a::Assembler, rowdofs, coldofs, Ke)

Assembles the matrix Ke into a according to the dofs specified by rowdofs and coldofs.

Ferrite.boundarydof_indicesMethod
boundarydof_indices(::Type{<:BoundaryIndex})

Helper function to generically dispatch on the correct dof sets of a boundary entity.

Ferrite.boundaryfunctionMethod
boundaryfunction(::Type{<:BoundaryIndex})

Helper function to dispatch on the correct entity from a given boundary index.

Ferrite.celldof_interior_indicesMethod
celldof_interior_indices(ip::Interpolation)

Tuple containing the dof indices associated with the interior of the cell.

Ferrite.celldofs!Method
celldofs!(global_dofs::Vector{Int}, dh::AbstractDofHandler, i::Int)

Store the degrees of freedom that belong to cell i in global_dofs.

See also celldofs.

Ferrite.celldofsMethod
celldofs(dh::AbstractDofHandler, i::Int)

Return a vector with the degrees of freedom that belong to cell i.

See also celldofs!.

Ferrite.close!Method
close!(ch::ConstraintHandler)

Close and finalize the ConstraintHandler.

Ferrite.close!Method
close!(dh::AbstractDofHandler)

Closes dh and creates degrees of freedom for each cell.

If there are several fields, the dofs are added in the following order: For a MixedDofHandler, go through each FieldHandler in the order they were added. For each field in the FieldHandler or in the DofHandler (again, in the order the fields were added), create dofs for the cell. This means that dofs on a particular cell, the dofs will be numbered according to the fields; first dofs for field 1, then field 2, etc.

Ferrite.collect_periodic_facesFunction
collect_periodic_faces(grid::Grid, all_faces::Union{Set{FaceIndex},String,Nothing}=nothing)

Split all faces in all_faces into image and mirror sets. For each matching pair, the face located further along the vector (1, 1, 1) becomes the image face.

If no set is given, all faces on the outer boundary of the grid (i.e. all faces that do not have a neighbor) is used.

See also: collect_periodic_faces!, PeriodicDirichlet.

Ferrite.collect_periodic_facesFunction
collect_periodic_faces(grid::Grid, mset, iset, transform::Union{Function,Nothing}=nothing)

Match all mirror faces in mset with a corresponding image face in iset. Return a dictionary which maps each mirror face to a image face. The result can then be passed to PeriodicDirichlet.

mset and iset can be given as a String (an existing face set in the grid) or as a Set{FaceIndex} directly.

By default this function looks for a matching face in the directions of the coordinate system. For other types of periodicities the transform function can be used. The transform function is applied on the coordinates of the image face, and is expected to transform the coordinates to the matching locations in the mirror set.

See also: collect_periodic_faces!, PeriodicDirichlet.

Ferrite.compute_vertex_valuesMethod
function compute_vertex_values(grid::AbstractGrid, f::Function)
function compute_vertex_values(grid::AbstractGrid, v::Vector{Int}, f::Function)    
function compute_vertex_values(grid::AbstractGrid, set::String, f::Function)

Given a grid and some function f, compute_vertex_values computes all nodal values, i.e. values at the nodes, of the function f. The function implements two dispatches, where only a subset of the grid's node is used.

    compute_vertex_values(grid, x -> sin(x[1]) + cos([2]))
    compute_vertex_values(grid, [9, 6, 3], x -> sin(x[1]) + cos([2])) #compute function values at nodes with id 9,6,3
    compute_vertex_values(grid, "right", x -> sin(x[1]) + cos([2])) #compute function values at nodes belonging to nodeset right
Ferrite.create_coloringFunction
create_coloring(g::Grid, cellset=1:getncells(g); alg::ColoringAlgorithm)

Create a coloring of the cells in grid g such that no neighboring cells have the same color. If only a subset of cells should be colored, the cells to color can be specified by cellset.

Returns a vector of vectors with cell indexes, e.g.:

ret = [
   [1, 3, 5, 10, ...], # cells for color 1
   [2, 4, 6, 12, ...], # cells for color 2
]

Two different algorithms are available, specified with the alg keyword argument:

  • alg = ColoringAlgorithm.WorkStream (default): Three step algorithm from WorkStream , albeit with a greedy coloring in the second step. Generally results in more colors than ColoringAlgorithm.Greedy, however the cells are more equally distributed among the colors.
  • alg = ColoringAlgorithm.Greedy: greedy algorithm that works well for structured quadrilateral grids such as e.g. quadrilateral grids from generate_grid.

The resulting colors can be visualized using vtk_cell_data_colors.

Cell to color mapping

In a previous version of Ferrite this function returned a dictionary mapping cell ID to color numbers as the first argument. If you need this mapping you can create it using the following construct:

colors = create_coloring(...)
cell_colormap = Dict{Int,Int}(
    cellid => color for (color, cellids) in enumerate(final_colors) for cellid in cellids
)
Ferrite.create_constraint_matrixMethod
create_constraint_matrix(ch::ConstraintHandler)

Create and return the constraint matrix, C, and the inhomogeneities, g, from the affine (linear) and Dirichlet constraints in ch.

The constraint matrix relates constrained, a_c, and free, a_f, degrees of freedom via a_c = C * a_f + g. The condensed system of linear equations is obtained as C' * K * C = C' * (f - K * g).

Ferrite.create_sparsity_patternMethod
create_sparsity_pattern(dh::AbstractDofHandler, ch::ConstraintHandler; coupling)

Create a sparsity pattern accounting for affine constraints in ch. See the Affine Constraints section of the manual for further details.

Ferrite.create_sparsity_patternMethod
create_sparsity_pattern(dh::DofHandler; coupling)

Create the sparsity pattern corresponding to the degree of freedom numbering in the DofHandler. Return a SparseMatrixCSC with stored values in the correct places.

The keyword argument coupling can be used to specify how fields (or components) in the dof handler couple to each other. coupling should be a square matrix of booleans with number of rows/columns equal to the total number of fields, or total number of components, in the DofHandler with true if fields are coupled and false if not. By default full coupling is assumed.

See the Sparsity Pattern section of the manual.

Ferrite.create_symmetric_sparsity_patternMethod
create_symmetric_sparsity_pattern(dh::AbstractDofHandler, ch::ConstraintHandler, coupling)

Create a symmetric sparsity pattern accounting for affine constraints in ch. See the Affine Constraints section of the manual for further details.

Ferrite.create_symmetric_sparsity_patternMethod
create_symmetric_sparsity_pattern(dh::DofHandler; coupling)

Create the symmetric sparsity pattern corresponding to the degree of freedom numbering in the DofHandler by only considering the upper triangle of the matrix. Return a Symmetric{SparseMatrixCSC}.

See the Sparsity Pattern section of the manual.

Ferrite.debug_modeMethod
Ferrite.debug_mode(; enable=true)

Helper to turn on (enable=true) or off (enable=false) debug expressions in Ferrite.

Debug mode influences Ferrite.@debug expr: when debug mode is enabled, expr is evaluated, and when debug mode is disabled expr is ignored.

Ferrite.default_interpolationMethod
Ferrite.default_interpolation(::AbstractCell)::Interpolation

Returns the interpolation which defines the geometry of a given cell.

Ferrite.dof_rangeMethod
dof_range(dh:DofHandler, field_name)

Return the local dof range for field_name. Example:

julia> grid = generate_grid(Triangle, (3, 3))
Grid{2, Triangle, Float64} with 18 Triangle cells and 16 nodes

julia> dh = DofHandler(grid); add!(dh, :u, 3); add!(dh, :p, 1); close!(dh);

julia> dof_range(dh, :u)
1:9

julia> dof_range(dh, :p)
10:12
Ferrite.dof_rangeMethod
dof_range(fh::FieldHandler, field_idx::Int)
dof_range(fh::FieldHandler, field_name::Symbol)
dof_range(dh:MixedDofHandler, field_name::Symbol)

Return the local dof range for a given field. The field can be specified by its name or index, where field_idx represents the index of a field within a FieldHandler and field_idxs is a tuple of the FieldHandler-index within the MixedDofHandler and the field_idx.

Note

The dof_range of a field can vary between different FieldHandlers. Therefore, it is advised to use the field_idxs or refer to a given FieldHandler directly in case several FieldHandlers exist. Using the field_name will always refer to the first occurrence of field within the MixedDofHandler.

Example:

julia> grid = generate_grid(Triangle, (3, 3))
Grid{2, Triangle, Float64} with 18 Triangle cells and 16 nodes

julia> dh = MixedDofHandler(grid); add!(dh, :u, 3); add!(dh, :p, 1); close!(dh);

julia> dof_range(dh, :u)
1:9

julia> dof_range(dh, :p)
10:12

julia> dof_range(dh, (1,1)) # field :u
1:9

julia> dof_range(dh.fieldhandlers[1], 2) # field :p
10:12
Ferrite.edge_npointsMethod
edge_npoints(::AbstractCell{dim,N,M)

Specifies for each subtype of AbstractCell how many nodes form an edge

Ferrite.edgedof_indicesMethod
edgedof_indices(ip::Interpolation)

A tuple containing tuples of local dof indices for the respective edge in local enumeration on a cell defined by edges(::Cell). The edge enumeration must match the edge enumeration of the corresponding geometrical cell.

Ferrite.edgedof_interior_indicesMethod
edgedof_interior_indices(ip::Interpolation)

A tuple containing tuples of the local dof indices on the interior of the respective edge in local enumeration on a cell defined by edges(::Cell). The edge enumeration must match the edge enumeration of the corresponding geometrical cell. Note that the vertex dofs are included here.

Ferrite.edgesMethod
Ferrite.edges(::AbstractCell)

Returns a tuple of 2-tuples containing the ordered node indices (of the nodes in a grid) corresponding to the vertices that define an oriented edge. This function induces the EdgeIndex, where the second index corresponds to the local index into this tuple.

Note that the vertices are sufficient to define an edge uniquely.

Ferrite.end_assembleMethod
end_assemble(a::Assembler) -> K

Finalizes an assembly. Returns a sparse matrix with the assembled values. Note that this step is not necessary for AbstractSparseAssemblers.

Ferrite.face_npointsMethod
face_npoints(::AbstractCell{dim,N,M)

Specifies for each subtype of AbstractCell how many nodes form a face

Ferrite.facedof_indicesMethod
facedof_indices(ip::Interpolation)

A tuple containing tuples of all local dof indices for the respective face in local enumeration on a cell defined by faces(::Cell). The face enumeration must match the face enumeration of the corresponding geometrical cell.

Ferrite.facedof_interior_indicesMethod
facedof_interior_indices(ip::Interpolation)

A tuple containing tuples of the local dof indices on the interior of the respective face in local enumeration on a cell defined by faces(::Cell). The face enumeration must match the face enumeration of the corresponding geometrical cell. Note that the vertex and edge dofs are included here.

Ferrite.facesMethod
Ferrite.faces(::AbstractCell)

Returns a tuple of n-tuples containing the ordered node indices (of the nodes in a grid) corresponding to the vertices that define an oriented face. This function induces the FaceIndex, where the second index corresponds to the local index into this tuple.

Note that the vertices are sufficient to define a face uniquely.

Ferrite.faceskeletonMethod
faceskeleton(grid) -> Vector{FaceIndex}

Returns an iterateable face skeleton. The skeleton consists of FaceIndex that can be used to reinit FaceValues.

Ferrite.fillzero!Method
fillzero!(A::AbstractVecOrMat{T})

Fill the (stored) entries of the vector or matrix A with zeros.

Fallback: fill!(A, zero(T)).

Ferrite.find_fieldMethod
find_field(dh::MixedDofHandler, field_name::Symbol)::NTuple{2,Int}

Return the index of the field with name field_name in a MixedDofHandler. The index is a NTuple{2,Int}, where the 1st entry is the index of the FieldHandler within which the field was found and the 2nd entry is the index of the field within the FieldHandler.

Note

Always finds the 1st occurrence of a field within MixedDofHandler.

See also: find_field(fh::FieldHandler, field_name::Symbol), _find_field(fh::FieldHandler, field_name::Symbol).

Ferrite.function_divergenceMethod
function_divergence(fe_v::Values, q_point::Int, u::AbstractVector)

Compute the divergence of the vector valued function in a quadrature point.

The divergence of a vector valued functions in the quadrature point $\mathbf{x}_q)$ is computed as $\mathbf{\nabla} \cdot \mathbf{u}(\mathbf{x_q}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i (\mathbf{x_q}) \cdot \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of the function.

Ferrite.function_gradientMethod
function_gradient(fe_v::Values{dim}, q_point::Int, u::AbstractVector)

Compute the gradient of the function in a quadrature point. u is a vector with values for the degrees of freedom. For a scalar valued function, u contains scalars. For a vector valued function, u can be a vector of scalars (for use of VectorValues) or u can be a vector of Vecs (for use with ScalarValues).

The gradient of a scalar function or a vector valued function with use of VectorValues is computed as $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} N_i (\mathbf{x}) u_i$ or $\mathbf{\nabla} u(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{\nabla} \mathbf{N}_i (\mathbf{x}) u_i$ respectively, where $u_i$ are the nodal values of the function. For a vector valued function with use of ScalarValues the gradient is computed as $\mathbf{\nabla} \mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n \mathbf{u}_i \otimes \mathbf{\nabla} N_i (\mathbf{x})$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

Ferrite.function_symmetric_gradientFunction
function_symmetric_gradient(fe_v::Values, q_point::Int, u::AbstractVector)

Compute the symmetric gradient of the function, see function_gradient. Return a SymmetricTensor.

The symmetric gradient of a scalar function is computed as $\left[ \mathbf{\nabla} \mathbf{u}(\mathbf{x_q}) \right]^\text{sym} = \sum\limits_{i = 1}^n \frac{1}{2} \left[ \mathbf{\nabla} N_i (\mathbf{x}_q) \otimes \mathbf{u}_i + \mathbf{u}_i \otimes \mathbf{\nabla} N_i (\mathbf{x}_q) \right]$ where $\mathbf{u}_i$ are the nodal values of the function.

Ferrite.function_valueMethod
function_value(fe_v::Values, q_point::Int, u::AbstractVector)

Compute the value of the function in a quadrature point. u is a vector with values for the degrees of freedom. For a scalar valued function, u contains scalars. For a vector valued function, u can be a vector of scalars (for use of VectorValues) or u can be a vector of Vecs (for use with ScalarValues).

The value of a scalar valued function is computed as $u(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) u_i$ where $u_i$ are the value of $u$ in the nodes. For a vector valued function the value is calculated as $\mathbf{u}(\mathbf{x}) = \sum\limits_{i = 1}^n N_i (\mathbf{x}) \mathbf{u}_i$ where $\mathbf{u}_i$ are the nodal values of $\mathbf{u}$.

Ferrite.generate_gridFunction
generate_grid(celltype::Cell, nel::NTuple, [left::Vec, right::Vec)

Return a Grid for a rectangle in 1, 2 or 3 dimensions. celltype defined the type of cells, e.g. Triangle or Hexahedron. nel is a tuple of the number of elements in each direction. left and right are optional endpoints of the domain. Defaults to -1 and 1 in all directions.

Ferrite.get_point_valuesFunction
get_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 NaNs for the corresponding entries in the output vector.

Ferrite.get_rhs_dataMethod
get_rhs_data(ch::ConstraintHandler, A::SparseMatrixCSC) -> RHSData

Returns the needed RHSData for apply_rhs!.

This must be used when the same stiffness matrix is reused for multiple steps, for example when timestepping, with different non-homogeneouos Dirichlet boundary conditions.

Ferrite.getcellsMethod
getcells(grid::AbstractGrid) 
getcells(grid::AbstractGrid, v::Union{Int,Vector{Int}} 
getcells(grid::AbstractGrid, setname::String)

Returns either all cells::Collection{C<:AbstractCell} of a <:AbstractGrid or a subset based on an Int, Vector{Int} or String. Whereas the last option tries to call a cellset of the grid. Collection can be any indexable type, for Grid it is Vector{C<:AbstractCell}.

Ferrite.getcellsetMethod
getcellset(grid::AbstractGrid, setname::String)

Returns all cells as cellid in a Set of a given setname.

Ferrite.getcoordinates!Method
getcoordinates!(x::Vector{Vec{dim,T}}, grid::AbstractGrid, cell::Int)
getcoordinates!(x::Vector{Vec{dim,T}}, grid::AbstractGrid, cell::AbstractCell)

Fills the vector x with the coordinates of a cell defined by either its cellid or the cell object itself.

Ferrite.getcoordinatesMethod
getcoordinates(grid::AbstractGrid, cell)

Return a vector with the coordinates of the vertices of cell number cell.

Ferrite.getcurrentfaceMethod
getcurrentface(fv::FaceValues)

Return the current active face of the FaceValues object (from last reinit!).

Ferrite.getdetJdVMethod
getdetJdV(fe_v::Values, q_point::Int)

Return the product between the determinant of the Jacobian and the quadrature point weight for the given quadrature point: $\det(J(\mathbf{x})) w_q$

This value is typically used when one integrates a function on a finite element cell or face as

$\int\limits_\Omega f(\mathbf{x}) d \Omega \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) \det(J(\mathbf{x})) w_q$ $\int\limits_\Gamma f(\mathbf{x}) d \Gamma \approx \sum\limits_{q = 1}^{n_q} f(\mathbf{x}_q) \det(J(\mathbf{x})) w_q$

Ferrite.getdimMethod
Ferrite.getdim(::Interpolation)

Return the dimension of the reference element for a given interpolation.

Ferrite.getedgesetMethod
getedgeset(grid::AbstractGrid, setname::String)

Returns all edges as EdgeIndex in a Set of a given setname.

Ferrite.getedgesetsMethod
getedgesets(grid::AbstractGrid)

Returns all edge sets of the grid.

Ferrite.getfacesetMethod
getfaceset(grid::AbstractGrid, setname::String)

Returns all faces as FaceIndex in a Set of a given setname.

Ferrite.getfielddimMethod
getfielddim(dh::MixedDofHandler, field_idxs::NTuple{2,Int})
getfielddim(dh::MixedDofHandler, field_name::Symbol)
getfielddim(dh::FieldHandler, field_idx::Int)
getfielddim(dh::FieldHandler, field_name::Symbol)

Return the dimension of a given field. The field can be specified by its index (see find_field) or its name.

Ferrite.getfieldinterpolationMethod
getfieldinterpolation(dh::MixedDofHandler, field_idxs::NTuple{2,Int})
getfieldinterpolation(dh::FieldHandler, field_idx::Int)
getfieldinterpolation(dh::FieldHandler, field_name::Symbol)

Return the interpolation of a given field. The field can be specified by its index (see find_field or its name.

Ferrite.getfieldnamesMethod
getfieldnames(dh::MixedDofHandler)
getfieldnames(fh::FieldHandler)

Return a vector with the names of all fields. Can be used as an iterable over all the fields in the problem.

Ferrite.getnbasefunctionsMethod
Ferrite.getnbasefunctions(ip::Interpolation)

Return the number of base functions for the interpolation ip.

Ferrite.getneighborhoodFunction
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, cellidx::CellIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, faceidx::FaceIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, vertexidx::VertexIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, edgeidx::EdgeIndex, include_self=false)

Returns all directly connected entities of the same type, i.e. calling the function with a VertexIndex will return a list of directly connected vertices (connected via face/edge). If include_self is true, the given *Index is included in the returned list.

Warning

This feature is highly experimental and very likely subjected to interface changes in the future.

Ferrite.getnodesMethod
getnodes(grid::AbstractGrid) 
getnodes(grid::AbstractGrid, v::Union{Int,Vector{Int}}
getnodes(grid::AbstractGrid, setname::String)

Returns either all nodes::Collection{N} of a <:AbstractGrid or a subset based on an Int, Vector{Int} or String. The last option tries to call a nodeset of the <:AbstractGrid. Collection{N} refers to some indexable collection where each element corresponds to a Node.

Ferrite.getnodesetMethod
getnodeset(grid::AbstractGrid, setname::String)

Returns all nodes as nodeid in a Set of a given setname.

Ferrite.getnormalMethod
getnormal(fv::FaceValues, qp::Int)

Return the normal at the quadrature point qp for the active face of the FaceValues object(from last reinit!).

Ferrite.getnquadpointsMethod
getnquadpoints(fe_v::Values)

Return the number of quadrature points for the Values object.

Ferrite.getorderMethod
Ferrite.getorder(::Interpolation)

Return order of the interpolation.

Ferrite.getpointsMethod
getpoints(qr::QuadratureRule)

Return the points of the quadrature rule.

Examples

julia> qr = QuadratureRule{2, RefTetrahedron}(:legendre, 2);

julia> getpoints(qr)
3-element Array{Tensors.Tensor{1,2,Float64,2},1}:
 [0.166667, 0.166667]
 [0.166667, 0.666667]
 [0.666667, 0.166667]
Ferrite.getrefshapeMethod
Ferrite.getrefshape(::Interpolation)::AbstractRefShape

Return the reference element shape of the interpolation.

Ferrite.getvertexsetMethod
getedgeset(grid::AbstractGrid, setname::String)

Returns all vertices as VertexIndex in a Set of a given setname.

Ferrite.getweightsMethod
getweights(qr::QuadratureRule)

Return the weights of the quadrature rule.

Examples

julia> qr = QuadratureRule{2, RefTetrahedron}(:legendre, 2);

julia> getweights(qr)
3-element Array{Float64,1}:
 0.166667
 0.166667
 0.166667
Ferrite.matrix_handleFunction
matrix_handle(a::AbstractSparseAssembler)
vector_handle(a::AbstractSparseAssembler)

Return a reference to the underlying matrix/vector of the assembler.

Ferrite.ndofsMethod
ndofs(dh::AbstractDofHandler)

Return the number of degrees of freedom in dh

Ferrite.ndofs_per_cellFunction
ndofs_per_cell(dh::AbstractDofHandler[, cell::Int=1])

Return the number of degrees of freedom for the cell with index cell.

See also ndofs.

Ferrite.projectMethod
project(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 Numbers and AbstractTensors.

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.

Ferrite.reinit!Function
reinit!(cv::CellValues, x::Vector)
reinit!(bv::FaceValues, x::Vector, face::Int)

Update the CellValues/FaceValues object for a cell or face with coordinates x. The derivatives of the shape functions, and the new integration weights are computed.

Ferrite.renumber!Function
renumber!(dh::AbstractDofHandler, order)
renumber!(dh::AbstractDofHandler, ch::ConstraintHandler, order)

Renumber the degrees of freedom in the DofHandler and/or ConstraintHandler according to the ordering order.

order can be given by one of the following options:

  • A permutation vector perm::AbstractVector{Int} such that dof i is renumbered to perm[i].
  • DofOrder.FieldWise() for renumbering dofs field wise.
  • DofOrder.ComponentWise() for renumbering dofs component wise.
  • DofOrder.Ext{T} for "external" renumber permutations, see documentation for DofOrder.Ext for details.
Warning

The dof numbering in the DofHandler and ConstraintHandler must always be consistent. It is therefore necessary to either renumber before creating the ConstraintHandler in the first place, or to renumber the DofHandler and the ConstraintHandler together.

Ferrite.reshape_to_nodesMethod
reshape_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.

Ferrite.shape_divergenceMethod
shape_divergence(fe_v::Values, q_point::Int, base_function::Int)

Return the divergence of shape function base_function evaluated in quadrature point q_point.

Ferrite.shape_gradientMethod
shape_gradient(fe_v::Values, q_point::Int, base_function::Int)

Return the gradient of shape function base_function evaluated in quadrature point q_point.

Ferrite.shape_symmetric_gradientMethod
shape_symmetric_gradient(fe_v::Values, q_point::Int, base_function::Int)

Return the symmetric gradient of shape function base_function evaluated in quadrature point q_point.

Ferrite.shape_valueMethod
shape_value(fe_v::Values, q_point::Int, base_function::Int)

Return the value of shape function base_function evaluated in quadrature point q_point.

Ferrite.sortedgeMethod
sortedge(edge::Tuple{Int,Int})

Returns the unique representation of an edge and its orientation. Here the unique representation is the sorted node index tuple. The orientation is true if the edge is not flipped, where it is false if the edge is flipped.

Ferrite.sortfaceMethod
sortface(face::Tuple{Int,Int}) 
sortface(face::Tuple{Int,Int,Int})
sortface(face::Tuple{Int,Int,Int,Int})

Returns the unique representation of a face. Here the unique representation is the sorted node index tuple. Note that in 3D we only need indices to uniquely identify a face, so the unique representation is always a tuple length 3.

Ferrite.spatial_coordinateMethod
spatial_coordinate(fe_v::Values{dim}, q_point::Int, x::AbstractVector)

Compute the spatial coordinate in a quadrature point. x contains the nodal coordinates of the cell.

The coordinate is computed, using the geometric interpolation, as $\mathbf{x} = \sum\limits_{i = 1}^n M_i (\mathbf{x}) \mathbf{\hat{x}}_i$

Ferrite.start_assembleFunction
start_assemble([N=0]) -> Assembler

Create an Assembler object which can be used to assemble element contributions to the global sparse matrix. Use assemble! for each element, and end_assemble, to finalize the assembly and return the sparse matrix.

Note that giving a sparse matrix as input can be more efficient. See below and as described in the manual.

Note

When the same matrix pattern is used multiple times (for e.g. multiple time steps or Newton iterations) it is more efficient to create the sparse matrix once and reuse the same pattern. See the manual section on assembly.

Ferrite.start_assembleMethod
start_assemble(K::SparseMatrixCSC;            fillzero::Bool=true) -> AssemblerSparsityPattern
start_assemble(K::SparseMatrixCSC, f::Vector; fillzero::Bool=true) -> AssemblerSparsityPattern

Create a AssemblerSparsityPattern from the matrix K and optional vector f.

start_assemble(K::Symmetric{SparseMatrixCSC};                 fillzero::Bool=true) -> AssemblerSymmetricSparsityPattern
start_assemble(K::Symmetric{SparseMatrixCSC}, f::Vector=Td[]; fillzero::Bool=true) -> AssemblerSymmetricSparsityPattern

Create a AssemblerSymmetricSparsityPattern from the matrix K and optional vector f.

AssemblerSparsityPattern and AssemblerSymmetricSparsityPattern allocate workspace necessary for efficient matrix assembly. To assemble the contribution from an element, use assemble!.

The keyword argument fillzero can be set to false if K and f should not be zeroed out, but instead keep their current values.

Ferrite.toglobalMethod
toglobal(grid::AbstractGrid, vertexidx::VertexIndex) -> Int
toglobal(grid::AbstractGrid, vertexidx::Vector{VertexIndex}) -> Vector{Int}

This function takes the local vertex representation (a VertexIndex) and looks up the unique global id (an Int).

Ferrite.transform!Method
transform!(grid::Abstractgrid, f::Function)

Transform all nodes of the grid based on some transformation function f.

Ferrite.update!Function
update!(ch::ConstraintHandler, time::Real=0.0)

Update time-dependent inhomogeneities for the new time. This calls f(x) or f(x, t) when applicable, where f is the function(s) corresponding to the constraints in the handler, to compute the inhomogeneities.

Note that this is called implicitly in close!(::ConstraintHandler).

Ferrite.vector_handleFunction
matrix_handle(a::AbstractSparseAssembler)
vector_handle(a::AbstractSparseAssembler)

Return a reference to the underlying matrix/vector of the assembler.

Ferrite.vertexdof_indicesMethod
vertexdof_indices(ip::Interpolation)

A tuple containing tuples of local dof indices for the respective vertex in local enumeration on a cell defined by vertices(::Cell). The vertex enumeration must match the vertex enumeration of the corresponding geometrical cell.

Ferrite.verticesMethod
Ferrite.vertices(::AbstractCell)

Returns a tuple with the node indices (of the nodes in a grid) for each vertex in a given cell. This function induces the VertexIndex, where the second index corresponds to the local index into this tuple.

Ferrite.vtk_cell_data_colorsFunction
vtk_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.

Ferrite.vtk_cellsetFunction
vtk_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.

Ferrite.vtk_cellsetMethod
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.

WriteVTK.vtk_gridMethod
vtk_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_dataMethod
vtk_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.