Finite element spaces
Bcube.AbstractFESpace — TypeAbstract 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.AbstractMultiFESpace — TypeDevs 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.MultiFESpace — TypeA 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.MultiFESpace — MethodMultiFESpace(
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 theMultiFESpace. Note that they must be of typeTrialFESpaceorTestFESpace.
Keywords
arrayOfStruct::Bool = AOS_DEFAULT: indicates if the dof numbering should be of type "Array of Structs" (AoS) or "Struct of Arrays" (SoA).
Bcube.SingleFESpace — TypeSingleFESpace(
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 theFESpacemesh::AbstractMesh: the mesh on which theFESpaceis discretizeddirichletBndNames = String[]: list of mesh boundary labels where a Dirichlet condition applies
Keywords
size::Int = 1: the number of components of theFESpaceisContinuous::Bool = true: iftrue, 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.SingleFESpace — TypeAn 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.TestFESpace — TypeTestFESpace(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.TestFESpace — TypeA TestFESpace is basically a SingleFESpace plus other attributes (related to boundary conditions)
Bcube.TrialFESpace — TypeTrialFESpace(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 functionf(x,t)to apply is expressed in the physical coordinate system. Alternatively, a constant value can be provided instead of a function.type::Symbol::continuousor: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.TrialFESpace — TypeA 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.MultiplierFESpace — FunctionA MultiplierFESpace can be viewed as a set of independent P0 elements. It is used to define Lagrange multipliers and assemble the associated augmented system (the system that adds the multipliers as unknowns).
Bcube._MultiFESpace — MethodLow-level constructor
Bcube._build_mapping_AoS — MethodBuild a global numbering using an Array-Of-Struct strategy
Bcube._build_mapping_SoA — MethodBuild a global numbering using an Struct-Of-Array strategy
Bcube.allocate_dofs — Functionallocate_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.build_jacobian_sparsity_pattern — Methodbuild_jacobian_sparsity_pattern(u::AbstractMultiFESpace, mesh::AbstractMesh)Build the jacobian sparsity pattern corresponding to the "FESpace U".
A jacobian is associated to a given function, here this function is assumed to mix (i.e "link") all the variables of the FESpace together.
Bcube.get_cell_shape_functions — MethodReturn the shape functions associated to the AbstractFESpace in "packed" form: λ(x) = (λ₁(x),...,λᵢ(x),...λₙ(x)) for the n dofs.
Bcube.get_dirichlet_boundary_tags — MethodReturn the boundary tags where a Dirichlet condition applies
Bcube.get_dirichlet_values — MethodReturn the values associated to a Dirichlet condition
Bcube.get_dofs — MethodReturn the dofs indices for the cell icell
Result is an array of integers.
Bcube.get_dofs — Methodget_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_fespace — Methodget_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_fespace — MethodReturn the tuple of FESpace composing this MultiFESpace
Bcube.get_function_space — MethodReturn the FunctionSpace (eventually multiple spaces) associated to the AbstractFESpace.
Bcube.get_mapping — Methodget_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_n_fespace — MethodNumber of FESpace composing the MultiFESpace
Bcube.get_ncomponents — MethodReturn the size S(= number of components) associated to AbstractFESpace{S}.
Bcube.get_ndofs — MethodReturn 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_ndofs — MethodTotal number of dofs contained in this MultiFESpace
Bcube.get_shape_functions — MethodReturn the shape functions associated to the AbstractFESpace.
Bcube.get_size — MethodReturn the size S associated to AbstractFESpace{S}.