Reference manual

Reference shapes

Elfel.RefShapes.RefShapeCubeType
RefShapeCube <: AbstractRefShape{3}

Type of a reference shape for a 3-dimensional manifold (solid) bounded by six quadrilaterals.

source
Elfel.RefShapes.manifdimvMethod
manifdimv(::Type{T}) where {T<:AbstractRefShape{MANIFDIM}} where {MANIFDIM}

Get the manifold dimension of the reference shape as Val.

source
Elfel.RefShapes.quadratureFunction
quadrature(::Type{RefShapeInterval}, quadraturesettings = (kind = :default,))

Create a quadrature rule for the reference shape of an interval.

The default is Gauss integration rule, where the order is set with the keyword order.

source
Elfel.RefShapes.quadratureFunction
quadrature(::Type{RefShapeTriangle}, quadraturesettings = (kind = :default,))

Create a quadrature rule for the reference shape of an triangle.

The default is a triangle rule, distinguished by the number of points set with the keyword npts.

source
Elfel.RefShapes.quadratureFunction
quadrature(::Type{RefShapeTetrahedron}, quadraturesettings = (kind = :default,))

Create a quadrature rule for the reference shape of a tetrahedron.

The default is a one-point tetrahedron rule; other rules may be chosen based on the number of points set with the keyword npts.

source
Elfel.RefShapes.quadratureFunction
quadrature(::Type{RefShapeSquare}, quadraturesettings = (kind = :default,))

Create a quadrature rule for the reference shape of a square.

The default is Gauss integration rule, where the order is set with the keyword order.

source

Elements

Elfel.FElements.FEType
FE{RS, SD}

Abstract type of finite element, parameterized by

  • RS: type of reference shape of the element (triangle, square, ...), and
  • SD: shape descriptor; refer to the package MeshCore.
source
Elfel.FElements.FEDataType
FEData{SD}

Type of a finite element data.

Parameterized by

  • SD = shape descriptor; refer to the package MeshCore.
source
Elfel.FElements.FEH1_T3_BUBBLEMethod
FEH1_T3_BUBBLE()

Construct an H1 finite element of the type T3 with a cubic bubble.

T3 is 3-node linear triangle element with a cubic bubble. It has the usual nodal basis functions associated with the vertices, and cubic bubble associated with the element itself.

source
Elfel.FElements.FEL2_T4Method
FEL2_T4()

Construct an L2 finite element of the type T4.

T4 is tetrahedral element with only internal degrees of freedom.

source
Elfel.FElements.JacobianMethod
Jacobian(::Val{0}, J::T) where {T}

Evaluate the point Jacobian.

  • J = Jacobian matrix, which isn't really defined well for a 0-manifold.
source
Elfel.FElements.JacobianMethod
Jacobian(::Val{1}, J::T) where {T}

Evaluate the curve Jacobian.

  • J = Jacobian matrix, columns are tangent to parametric coordinates curves.
source
Elfel.FElements.JacobianMethod
Jacobian(::Val{2}, J::T) where {T}

Evaluate the curve Jacobian.

  • J = Jacobian matrix, columns are tangent to parametric coordinates curves.
source
Elfel.FElements.JacobianMethod
Jacobian(fe::T, J::FFltMat)::FFlt where {T<:AbstractFESet3Manifold}

Evaluate the volume Jacobian.

J = Jacobian matrix, columns are tangent to parametric coordinates curves.

source
Elfel.FElements._geometrycarrierMethod

The elements that have nodal bases can be their own geometry carriers. Elements that do not have nodal bases, for instance an L2 quadrilateral with a single basis function associated with the cell, need another element type (in this case the nodal quadrilateral) to serve as geometry carriers.

source
Elfel.FElements.bfunMethod
bfun(fe::FESUBT,  param_coords)  where {FESUBT<:FE{RS, SD}}

Evaluate the basis functions for all degrees of freedom of the scalar finite element at the parametric coordinates. Return a vector of the values.

source
Elfel.FElements.bfungradparMethod
bfungradpar(fe::FESUBT,  param_coords)  where {FESUBT<:FE{RS, SD}}

Evaluate the gradients of the basis functions for all degrees of freedom of the scalar finite element with respect to the parametric coordinates, at the parametric coordinates given. Return a vector of the gradients.

source
Elfel.FElements.jacjacMethod
jacjac(fe::FE{RS, SD}, locs, nodes, gradNpar) where {RS, SD}

Compute the Jacobian matrix and the Jacobian determinant.

This is the generic version suitable for isoparametric elements.

source
Elfel.FElements.ndofperfeatMethod
feathasdof(fe::FE{RS, SD}, m) where {RS, SD}

How many degrees of freedom are attached to a the feature of manifold dimension m?

Note that 0 <= m <= 3.

source
Elfel.FElements.ndofsperelMethod
ndofsperel(fe::FE{RS, SD}) where {RS, SD}

Provide the number of degrees of freedom per element.

Enumerate all features of all manifold dimensions, and for each feature multiply by the number of degrees of freedom per feature. The assumption is that this is a scalar finite element.

source
MeshCore.manifdimMethod
manifdim(fe::FE{RS, SD}) where {RS, SD}

Get the manifold dimension of the finite element.

source

Spaces

Elfel.FESpaces.FESpaceType
FESpace{FET, T}

Type of a finite element space, parameterized with

  • FET: type of finite element, it is a scalar finite element,
  • T: type of degree of freedom value (double, complex, ...).
source
Elfel.FESpaces.FEFields.numberdatadofs!Method
numberdatadofs!(fesp::FES, firstnum = 1)  where {FES<:FESpace}

Number the data (known) degrees of freedom.

The known degrees of freedom in the FE space are numbered consecutively.

No effort is made to optimize the numbering in any way.

source
Elfel.FESpaces.FEFields.numberfreedofs!Method
numberfreedofs!(fesp::FES, firstnum = 1)  where {FES<:FESpace}

Number the free degrees of freedom.

The unknown degrees of freedom in the FE space are numbered consecutively.

No effort is made to optimize the numbering in any way.

source
Elfel.FESpaces.FEFields.setebc!Method
setebc!(fesp::FESpace, m, eid, comp, val::T) where {T}

Set the EBCs (essential boundary conditions).

  • m = manifold dimension of the entity,
  • eid = serial number of the entity (term identifier),
  • comp = which degree of freedom in the term,
  • val = value of type T

For instance, m = 0 means set the degree of freedom at the vertex eid.

source
Elfel.FESpaces.dofnumMethod
dofnum(fesp::FES, m, eid)  where {FES<:FESpace}

Provide degree of freedom number for entity eid of manifold dimension m and component comp.

source
Elfel.FESpaces.doftypeMethod
doftype(fesp::FESpace{FET, T}) where {FET, T}

Provide the type of the values of the degrees of freedom.

source
Elfel.FESpaces.edofbfnumMethod
edofbfnum(fesp::FESpace{FET, T}) where {FET, T}

Access vector of numbers of basis functions associated with each degree of freedom.

source
Elfel.FESpaces.edofcompntMethod
edofcompnt(fesp::FESpace{FET, T}) where {FET, T}

Access vector of component number associated with each degree of freedom.

When the finite element space consists of multiple copies of the scalar finite element, the component is the serial number of the copy.

source
Elfel.FESpaces.edofmdimMethod
edofmdim(fesp::FESpace{FET, T}) where {FET, T}

Access vector of manifold dimensions of entities associated with each degree of freedom.

source
Elfel.FESpaces.numberdofs!Method
numberdofs!(fesp::AbstractVector)

Number the degrees of freedom of a collection of FE spaces.

The unknown (free) degrees of freedom in the FE space are numbered consecutively, and then the data degrees of freedom (the known values) are numbered.

No effort is made to optimize the numbering in any way.

source
Elfel.FESpaces.numberdofs!Method
numberdofs!(fesp::FESpace)

Number the degrees of freedom of a single FE space.

The unknown (free) degrees of freedom in the FE space are numbered consecutively, and then the data degrees of freedom (the known values) are numbered.

No effort is made to optimize the numbering in any way.

source
Elfel.FESpaces.setdofval!Method
setdofval!(fesp::FESpace, m, eid, comp, val::T) where {T}

Set the value of a degree of freedom.

  • m = manifold dimension of the entity,
  • eid = serial number of the entity (term identifier),
  • comp = which degree of freedom in the term,
  • val = value of type T

For instance, m = 0 means set the degree of freedom at the vertex eid.

source
Elfel.FElements.ndofsperelMethod
ndofsperel(fesp::FES)  where {FES<:FESpace}

Total number of degrees of freedom associated with each finite element.

Essentially a product of the number of the degrees of freedom the scalar finite element and the number of copies of this element in the space.

source

Fields

Elfel.FESpaces.FEFields.FEFieldType
FEField{N, T, IT}

Type of a finite element field. Parameterized with

  • N: number of degrees of freedom per entity,
  • T: type of the degree of freedom value,
  • IT: type of the index (integer value). This describes the serial numbers of the degrees of freedom.
source
Elfel.FESpaces.FEFields.datadofnumsMethod
datadofnums(f::FEField)

Collect information about known (data) degree of freedom numbers.

First number, last number, and the total number of degrees of freedom are returned as a tuple.

source
Elfel.FESpaces.FEFields.freedofnumsMethod
freedofnums(f::FEField)

Collect information about unknown (free) degree of freedom numbers.

First number, last number, and the total number of degrees of freedom are returned as a tuple.

source
Elfel.FESpaces.FEFields.numberdatadofs!Function
numberdatadofs!(f::FEField, firstnum = 1)

Number the data degrees of freedom in the field. Start from the number supplied on input.

Note: The free degrees of freedom must be numbered first.

source
Elfel.FESpaces.FEFields.numberfreedofs!Function
numberfreedofs!(f::FEField, firstnum = 1)

Number the unknowns in the field, starting from the one supplied on input.

Note: The data degrees of freedom have their numbers zeroed out.

source
Elfel.FESpaces.FEFields.setdofval!Method
setdofval!(self::FEField, tid, comp, val::T) where {T}

Set the value of one particular degree of freedom to a given number.

  • tid: which term,
  • comp: which component of the term,
  • val: value to which the degree of freedom should be set.
source
Elfel.FESpaces.FEFields.setebc!Method
setebc!(self::FEField, tid, comp, val::T) where {T}

Set the value of one particular degree of freedom to a given number.

  • tid: which term,
  • comp: which component of the term,
  • val: value to which the degree of freedom should be set.
source
Elfel.FESpaces.FEFields.setisdatum!Method
setebc!(self::FEField, tid, comp, val::T) where {T}

Set the value of one particular degree of freedom to a given number.

  • tid: which term,
  • comp: which component of the term,
  • flag: true or false.
source

Finite element iterators

Elfel.FEIterators.FEIteratorType
FEIterator{FES, IR, G, IT, T, V, IR0, IR1, IR2, IR3, F0, F1, F2, F3}

Type of finite element iterator. Parameterized with the types of

  • FES: finite element space,
  • IR: base incidence relation of the mesh,
  • G: type of the geometry attribute,
  • IT: type of integer indices, such as the numbers of nodes and degrees of freedom,
  • T: type of the degree of freedom value (real double, complex float, ... ),
  • V: Val representation of the manifold dimension of the base relation elements,
  • IR0, IR1, IR2, IR3: types of incidence relations with which degrees of freedom are associated in the finite element space, for each of the manifolds dimensions 0, 1, 2, 3,
  • F0, F1, F2, F3: types of fields with which degrees of freedom are associated in the finite element space, for each of the manifolds dimensions 0, 1, 2, 3.
source
Base.iterateFunction
Base.iterate(it::FEIterator, state = 1)

Advance the iterator to the next entity.

The nodes of the finite element are cached, as is a vector of all the degrees of freedom represented on the element.

source
Base.lengthMethod
Base.length(it::FEIterator)

Number of elements represented by this iterator.

source
Elfel.FEIterators.eldofcompsMethod
eldofcomps(it::FEIterator)

Retrieve the vector of the component numbers for each element degree of freedom.

If multiple copies of the finite element are referenced in the finite element space, each copy is referred to as component.

source
Elfel.FEIterators.eldofentmdimsMethod
eldofentmdims(it::FEIterator)

Retrieve the vector of the entity dimensions for each element degree of freedom.

Each degree of freedom is associated with some entity of the finite element: vertices, edges, faces, and so on. This vector records the dimension of the manifold entity with which each degree of freedom is associated.

source
Elfel.FElements.jacjacMethod
jacjac(it::FEIterator, qpit::QPIterator)

Compute the Jacobian matrix and the Jacobian determinant.

The finite element iterator cooperates with the quadrature point iterator here to compute the Jacobian at the current integration point.

source

Quadrature-point iterators

Elfel.QPIterators.QPIteratorType
QPIterator{FES, MDIM}

Type of quadrature-point iterator, parameterized by

  • FES: the type of the finite element space,
  • MDIM: the manifold dimension of the finite element.
source
Elfel.QPIterators.QPIteratorMethod
QPIterator(fesp::FES, quadraturesettings) where {FES}

Construct quadrature-point iterator by associating it with a finite element space and supplying quadrature rule settings.

source
Base.iterateFunction
Base.iterate(it::QPIterator, state = 1)

Advance a quadrature point iterator.

source
Elfel.FElements.bfunMethod
bfun(it::QPIterator)

Retrieve vector of basis function values for the current quadrature point.

source
Elfel.FElements.bfungradparMethod
bfungradpar(it::QPIterator)

Retrieve vector of basis function gradients with respect to the parametric coordinates for the current quadrature point.

source
Elfel.QPIterators.bfungradMethod
bfungrad(it::QPIterator, Jac)

Retrieve vector of basis function gradients with respect to spatial coordinates for the current quadrature point.

The Jacobian matrix maps between vectors in the parametric space and the spatial vectors.

source

Assemblers

Elfel.Assemblers.SysmatAssemblerSparseMethod
SysmatAssemblerSparse(zero::T=0.0) where {T<:Number}

Construct blank system matrix assembler. The matrix entries are of type T.

Example

This is how a sparse matrix is assembled from two rectangular dense matrices.

a = SysmatAssemblerSparse(0.0)                                                        
start!(a, 7, 7)  
m = [0.24406   0.599773    0.833404  0.0420141                                             
0.786024  0.00206713  0.995379  0.780298                                              
0.845816  0.198459    0.355149  0.224996]     
gi = [1 7 5]             
gj = [5 2 1 4]       
for j in 1:size(m, 2), i in 1:size(m, 1)
    assemble!(a, gi[i], gj[j], m[i, j])       
end  
m = [0.146618  0.53471   0.614342    0.737833                                              
0.479719  0.41354   0.00760941  0.836455                                              
0.254868  0.476189  0.460794    0.00919633                                            
0.159064  0.261821  0.317078    0.77646                                               
0.643538  0.429817  0.59788     0.958909]                                   
gi =  [2 3 1 7 5]
gj = [6 7 3 4]   
for j in 1:size(m, 2), i in 1:size(m, 1)
    assemble!(a, gi[i], gj[j], m[i, j])       
end                               
A = finish!(a) 
source
Elfel.Assemblers.assemble!Method
assemble!(self::SysmatAssemblerSparse{T}, r, c, v::T) where {T<:Number}

Assemble a single entry of a rectangular matrix.

source
Elfel.Assemblers.assemble!Method
assemble!(self::SysmatAssemblerSparse{T}, lma::LocalMatrixAssembler{IT, T}) where {IT<:Integer, T<:Number}

Assemble the row numbers, column numbers, and values from a local assembler.

source
Elfel.Assemblers.assemble!Method
assemble!(self::SysmatAssemblerSparse{T}, lma::Transpose{T,LocalMatrixAssembler{IT,T}}) where {IT<:Integer, T<:Number}

Assemble the row numbers, column numbers, and values from a local assembler.

source
Elfel.Assemblers.assemble!Method
assemble!(self::SV, i, val::T) where {SV<:AbstractSysvecAssembler, T<:Number}

Assemble an elementwise vector.

The method assembles a column element vector using the vector of degree of freedom numbers for the rows.

source
Elfel.Assemblers.start!Method
start!(self::SysmatAssemblerSparse{T}, nrow, ncol) where {T<:Number}

Start the assembly of a global matrix.

source
Elfel.Assemblers.start!Method
start!(self::SysvecAssembler{T},  nrow::Int64) where {T<:Number}

Start assembly.

The method makes the buffer for the vector assembly. It must be called before the first call to the method assemble.

nrow= Total number of degrees of freedom.

source
Elfel.Assemblers.start!Method
start!(self::SV,  nrow) where {SV<:AbstractSysvecAssembler, T<:Number}

Start assembly.

The method makes the buffer for the vector assembly. It must be called before the first call to the method assemble.

nrow= Total number of degrees of freedom.

source

Local Assemblers

Base.IndexStyleMethod
Base.IndexStyle(::Type{<:LocalMatrixAssembler})

The data storage is assumed to be consumed best by one-dimensional traversal through the columns in a linear fashion.

source
Base.IndexStyleMethod
Base.IndexStyle(::Type{<:LocalVectorAssembler})

Only linear access is provided.

source
Elfel.LocalAssemblers.LocalMatrixAssemblerType
LocalMatrixAssembler{IT<:Integer, T<:Number} <: AbstractArray{T, 2}

Type of "local" matrix assembler.

Local is to be understood in the sense of in the context of a single finite element. So a local matrix is an elementwise matrix which is computed entry by entry. Then it can be assembled into the global matrix in one shot.

source
Elfel.LocalAssemblers.LocalMatrixAssemblerMethod
LocalMatrixAssembler(nrow::IT, ncol::IT, z::T) where {IT, T}

Create a local matrix assembler, given the number of rows and columns, and the value to which the matrix should be initialized.

source
Elfel.LocalAssemblers.LocalVectorAssemblerType
LocalVectorAssembler{IT<:Integer, T<:Number} <: AbstractArray{T, 1}

Type of "local" vector assembler.

Local is to be understood in the sense of in the context of a single finite element. So a local vector is an elementwise vector which is computed entry by entry. Then it can be assembled into the global vector in one shot.

source
Base.getindexMethod
Base.getindex(a::A, i::Int, j::Int)

Only access to a single entry of the matrix is provided.

source
Base.getindexMethod
Base.getindex(a::A, i::Int) where {A<:LocalVectorAssembler}

Access is provided to a single entry of the vector.

source
Base.setindex!Method
Base.setindex!(a::A, v, i::Int, j::Int) where {A<:LocalMatrixAssembler}

Only access to a single entry of the matrix is provided.

source
Base.setindex!Method
Base.setindex!(a::A, v, i::Int) where {A<:LocalVectorAssembler}

Access is provided to a single entry of the vector.

source
Base.sizeMethod
Base.size(a::A) where {A<:LocalMatrixAssembler}

The size is the tuple of number of rows and number of columns.

source
Base.sizeMethod
Base.size(a::A) where {A<:LocalVectorAssembler}

The size is the number of rows (entries).

source
Elfel.LocalAssemblers.init!Method
init!(a::L, rdofs, cdofs) where {L<:LocalMatrixAssembler{IT, T}} where {IT, T}

Initialize the local assembler with the global degrees of freedom in the rows and columns.

The two arrays, rdofs, cdofs, define the global degree of freedom numbers for the element. The data matrix is zeroed out.

This function needs to be called for each new finite element.

source
Elfel.LocalAssemblers.init!Method
init!(a::L, rdofs) where {L<:LocalVectorAssembler{IT, T}} where {IT, T}

Initialize the local assembler with the global degrees of freedom in the rows.

The array rdofs defines the global degree of freedom numbers for the element. The data vector is zeroed out.

This function needs to be called for each new finite element.

source

Index