Finite element spaces

Bcube.AbstractFESpaceType

Abstract type to represent an finite-element space of size S. See SingleFESpace for more details about what looks like a finite-element space.

Devs notes

All subtypes should implement the following functions:

  • get_function_space(feSpace::AbstractFESpace)
  • get_shape_functions(feSpace::AbstractFESpace, shape::AbstractShape)
  • get_cell_shape_functions(feSpace::AbstractFESpace, shape::AbstractShape)
  • get_ndofs(feSpace::AbstractFESpace)
  • is_continuous(feSpace::AbstractFESpace)

Alternatively, you may define a "parent" to your structure by implementing the Base.parent function. Then, all the above functions will be redirected to the "parent" FESpace.

Bcube.AbstractMultiFESpaceType

Devs notes

All subtypes should implement the following functions:

  • get_fespace(mfeSpace::AbstractMultiFESpace)
  • get_mapping(mfeSpace::AbstractMultiFESpace)
  • get_dofs(mfeSpace::AbstractMultiFESpace, icell::Int)
  • get_shape_functions(mfeSpace::AbstractMultiFESpace, shape::AbstractShape)
  • get_cell_shape_functions(mfeSpace::AbstractMultiFESpace, shape::AbstractShape)
Bcube.MultiFESpaceType

A MultiFESpace represents a "set" of TrialFESpace or TestFESpace. This structure provides a global dof numbering for each FESpace.

N is the number of FESpace contained in this MultiFESpace.

Note that the FESpace can be different from each other (one continous, one discontinuous; one scalar, one vector...)

Bcube.MultiFESpaceMethod
MultiFESpace(
    feSpaces::Tuple{Vararg{TrialOrTest, N}};
    arrayOfStruct::Bool = AOS_DEFAULT,
) where {N}
MultiFESpace(
    feSpaces::AbstractArray{FE};
    arrayOfStruct::Bool = AOS_DEFAULT,
) where {FE <: TrialOrTest}
MultiFESpace(feSpaces::Vararg{TrialOrTest}; arrayOfStruct::Bool = AOS_DEFAULT)

Build a finite element space representing several sub- finite element spaces.

This is particulary handy when several variables are in play since it provides a global dof numbering (for the whole system). The finite element spaces composing the MultiFESpace can be different from each other (some continuous, some discontinuous, some scalar, some vectors...).

Arguments

  • feSpaces : the finite element spaces composing the MultiFESpace. Note that they must be of type TrialFESpace or TestFESpace.

Keywords

  • arrayOfStruct::Bool = AOS_DEFAULT : indicates if the dof numbering should be of type "Array of Structs" (AoS) or "Struct of Arrays" (SoA).
Bcube.SingleFESpaceType
SingleFESpace(
    fSpace::AbstractFunctionSpace,
    mesh::AbstractMesh,
    dirichletBndNames = String[];
    size::Int = 1,
    isContinuous::Bool = true,
    kwargs...
)

Build a finite element space (scalar or vector) from a FunctionSpace and a Mesh.

Arguments

  • fSpace::AbstractFunctionSpace : the function space associated to the FESpace
  • mesh::AbstractMesh : the mesh on which the FESpace is discretized
  • dirichletBndNames = String[] : list of mesh boundary labels where a Dirichlet condition applies

Keywords

  • size::Int = 1 : the number of components of the FESpace
  • isContinuous::Bool = true : if true, a continuous dof numbering is created. Otherwise, dof lying

on cell nodes or cell faces are duplicated, not shared (discontinuous dof numbering)

  • kwargs : for things such as parallel cache (internal/dev usage only)
Bcube.SingleFESpaceType

An finite-element space (FESpace) is basically a function space, associated to degrees of freedom (on a mesh).

A FESpace can be either scalar (to represent a Temperature for instance) or vector (to represent a Velocity). In case of a "vector" SingleFESpace, all the components necessarily share the same FunctionSpace.

Bcube.TestFESpaceType
TestFESpace(trialFESpace::TrialFESpace)
TestFESpace(
    fSpace::AbstractFunctionSpace,
    mesh::AbstractMesh,
    dirichletBndNames = String[];
    size::Int = 1,
    isContinuous::Bool = true,
    kwargs...,
)

Build a test finite element space.

A TestFESpace can be built from a TrialFESpace. See SingleFESpace for hints about the function arguments. Only arguments specific to TrialFESpace are detailed below.

Examples

julia> mesh = one_cell_mesh(:line)
julia> fSpace = FunctionSpace(:Lagrange, 2)
julia> U = TrialFESpace(fSpace, mesh)
julia> V = TestFESpace(U)
Bcube.TestFESpaceType

A TestFESpace is basically a SingleFESpace plus other attributes (related to boundary conditions)

Bcube.TrialFESpaceType

A TrialFESpace is basically a SingleFESpace plus other attributes (related to boundary conditions)

Dev notes

  • we cannot directly store Dirichlet values on dofs because the Dirichlet values needs "time" to apply
Bcube.TrialFESpaceType
TrialFESpace(feSpace, dirichletValues)
TrialFESpace(
    fSpace::AbstractFunctionSpace,
    mesh::AbstractMesh,
    dirichlet::Dict{String} = Dict{String, Any}();
    size::Int = 1,
    isContinuous::Bool = true,
    kwargs...
)
TrialFESpace(
    fSpace::AbstractFunctionSpace,
    mesh::AbstractMesh,
    type::Symbol,
    dirichlet::Dict{String} = Dict{String, Any}();
    size::Int = 1,
    kwargs...
)

Build a trial finite element space.

See SingleFESpace for hints about the function arguments. Only arguments specific to TrialFESpace are detailed below.

Arguments

  • dirichlet::Dict{String} = Dict{String, Any}() : dictionnary specifying the Dirichlet valued-function (or function) associated to each mesh boundary label. The function f(x,t) to apply is expressed in the physical coordinate system. Alternatively, a constant value can be provided instead of a function.
  • type::Symbol : :continuous or :discontinuous

Warning

For now the Dirichlet condition can only be applied to nodal bases.

Examples

julia> mesh = one_cell_mesh(:line)
julia> fSpace = FunctionSpace(:Lagrange, 2)
julia> U = TrialFESpace(fSpace, mesh)
julia> V = TrialFESpace(fSpace, mesh, :discontinuous; size = 3)
julia> W = TrialFESpace(fSpace, mesh, Dict("North" => 3., "South" => (x,t) -> t .* x))
Bcube.MultiplierFESpaceFunction

A MultiplierFESpace can be viewed as a set of independant P0 elements. It is used to define Lagrange multipliers and assemble the associated augmented system (the system that adds the multipliers as unknowns).

Bcube.allocate_dofsFunction
allocate_dofs(feSpace::AbstractFESpace, T = Float64)

Allocate a vector with a size equal to the number of dof of the FESpace, with the type T. For a MultiFESpace, a vector of the total size of the space is returned (and not a Tuple of vectors)

Bcube.get_cell_shape_functionsMethod

Return the shape functions associated to the AbstractFESpace in "packed" form: λ(x) = (λ₁(x),...,λᵢ(x),...λₙ(x)) for the n dofs.

Bcube.get_dofsMethod

Return the dofs indices for the cell icell

Result is an array of integers.

Bcube.get_dofsMethod
get_dofs(feSpace::MultiFESpace, icell::Int)

Return the dofs indices for the cell icell for each single-feSpace. Result is a tuple of array of integers, where each array of integers are the indices relative to the numbering of each singleFESpace.

Warning:

Combine get_dofs with get_mapping if global dofs indices are needed.

Bcube.get_fespaceMethod
get_fespace(mfeSpace::AbstractMultiFESpace, iSpace)
get_fespace(mfeSpace::AbstractMultiFESpace)

Return the i-th FESpace composing this AbstractMultiFESpace. If no index is provided, the tuple of FESpace composing this AbstractMultiFESpace` is returnted.

Bcube.get_function_spaceMethod

Return the FunctionSpace (eventually multiple spaces) associated to the AbstractFESpace.

Bcube.get_mappingMethod
get_mapping(mfeSpace::AbstractMultiFESpace, iSpace)
get_mapping(mfeSpace::AbstractMultiFESpace)

Return the mapping for the ith FESpace composing the MultiFESpace. If no index is provided, the tuple of mapping for each FESpace` is returnted.

Bcube.get_ncomponentsMethod

Return the size S(= number of components) associated to AbstractFESpace{S}.

Bcube.get_ndofsMethod

Return the total number of dofs of the FESpace, taking into account the continuous/discontinuous type of the space. If the FESpace contains itself several FESpace (see MultiFESpace), the sum of all dofs is returned.

Bcube.get_ndofsMethod

Total number of dofs contained in this MultiFESpace

Bcube.get_sizeMethod

Return the size S associated to AbstractFESpace{S}.