Degree of freedom
Assembler
Bcube.AbstractBilinearFaceSidePair
— TypeAbstractBilinearFaceSidePair{A} <: AbstractLazyWrap{A}
Interface:
* get_args_bilinear(a::AbstractBilinearFaceSidePair)
Bcube.AbstractFaceSidePair
— TypeAbstractFaceSidePair{A} <: AbstractLazyWrap{A}
Interface:
side_n(a::AbstractFaceSidePair)
side_p(a::AbstractFaceSidePair)
Bcube.__assemble_linear!
— MethodDev notes
Two levels of "LazyMapOver" because first we LazyMapOver the Tuple of argument of the linear form, and the for each item of this Tuple we LazyMapOver the shape functions.
Bcube._append_contribution!
— Methodbilinear case
Bcube._assemble_linear!
— Method_assemble_linear!(b, l, V, integration::Integration)
_assemble_linear!(b, l, V, integration::MultiIntegration{N}) where {N}
These functions act as a function barrier in order to:
- get the function corresponding to the operand in the linear form
- reshape
b
internally to deal with cases whenV
is aAbstractMultiTestFESpace
- call
__assemble_linear!
to apply dispatch on the type ofmeasure
of the integration and improve type stability during the assembling loop.
Dev note:
The case integration::MultiIntegration{N}
is treated by looping over each Integration
contained in the MultiIntegration
Bcube._blockmap_bilinear
— MethodFrom tuples $a=(a_1, a_2, …, a_i, …, a_m)$ and $b=(b_1, b_2, …, b_j, …, b_n)$, it builds A
and B
which correspond formally to the following two matrices :
\[A \equiv \begin{pmatrix} a_1 & a_1 & ⋯ & a_1\\ a_2 & a_2 & ⋯ & a_2\\ ⋮ & ⋮ & ⋮ & ⋮ \\ a_m & a_m & ⋯ & a_m \end{pmatrix} \qquad and \qquad B \equiv \begin{pmatrix} b_1 & b_2 & ⋯ & b_n\\ b_1 & b_2 & ⋯ & b_n\\ ⋮ & ⋮ & ⋮ & ⋮ \\ b_1 & b_2 & ⋯ & b_n \end{pmatrix}\]
A
and B
are wrapped in LazyMapOver
structures so that all operations on them are done elementwise by default (in other words, it can be considered that the operations are automatically broadcasted).
Dev note :
Both A
and B
are stored as a tuple of tuples, wrapped by LazyMapOver
, where inner tuples correspond to each columns of a matrix. This hierarchical structure reduces both inference and compile times by avoiding the use of large tuples.
Bcube._diag_tuples
— Method_diag_tuples(diag::Tuple{Vararg{Any,N}}, b) where N
Return N
tuples of length N
. For each tuple tᵢ
, its values are defined so that tᵢ[k]=diag[k]
if k==i
, tᵢ[k]=b
otherwise. The result can be seen as a dense diagonal-like array using tuple.
Example for N=3
:
(diag[1], b, b ),
(b, diag[2], b ),
(b, b, diag[3]))
Bcube._get_multi_tuple_var
— MethodFor N=3 for example: (LazyMapOver((LazyMapOver(V[1]), NullOperator(), NullOperator())), LazyMapOver((NullOperator(), LazyMapOver(V[2]), NullOperator())), LazyMapOver((NullOperator(), NullOperator(), LazyMapOver(V[3]))))
Bcube._may_reshape_b
— MethodFor AbstractMultiTestFESpace
, it creates a Tuple (of views) of the different "destination" in the vector: one for each FESpace
Bcube.assemble_bilinear!
— Methodassemble_bilinear!(
I::Vector{Int},
J::Vector{Int},
X::Vector{T},
f::Function,
measure::Measure,
U::TrialFESpace,
V::TestFESpace,
)
In-place version of assemble_bilinear
.
Bcube.assemble_bilinear
— Methodassemble_bilinear(a::Function, U, V)
Assemble the (sparse) Matrix corresponding to the given bilinear form a
on the trial and test finite element spaces U
and V
.
For the in-place version, check-out assemble_bilinear!
.
Arguments
a::Function
: function of two variables (u,v) representing the bilinear formU
: trial finite element space (foru
)V
: test finite element space (forv
)
Examples
julia> mesh = rectangle_mesh(3,3)
julia> U = TrialFESpace(FunctionSpace(:Lagrange, 0), mesh)
julia> V = TestFESpace(U)
julia> dΩ = Measure(CellDomain(mesh), 3)
julia> a(u, v) = ∫(u * v)dΩ
julia> assemble_bilinear(a, U, V)
4×4 SparseArrays.SparseMatrixCSC{Float64, Int64} with 4 stored entries:
0.25 ⋅ ⋅ ⋅
⋅ 0.25 ⋅ ⋅
⋅ ⋅ 0.25 ⋅
⋅ ⋅ ⋅ 0.25
Bcube.assemble_linear!
— Methodassemble_linear!(b::AbstractVector, l::Function, V::Union{TestFESpace, AbstractMultiTestFESpace})
In-place version of assemble_linear
.
Bcube.assemble_linear
— Methodassemble_linear(l::Function, V::Union{TestFESpace, AbstractMultiTestFESpace})
Assemble the vector corresponding to a linear form l
on the finite element space V
For the in-place version, checkout assemble_linear!
.
Arguments
l::Function
: linear form to assemble, a function of one variablel(v)
V
: test finite element space
Examples
julia> mesh = rectangle_mesh(3,3)
julia> U = TrialFESpace(FunctionSpace(:Lagrange, 0), mesh)
julia> V = TestFESpace(U)
julia> dΩ = Measure(CellDomain(mesh), 3)
julia> l(v) = ∫(v)dΩ
julia> assemble_linear(l, V)
4-element Vector{Float64}:
0.25
0.25
0.25
0.25
Bcube.blockmap_bilinear_shape_functions
— MethodDev notes:
Return blockU
and blockV
to be able to compute the local matrix corresponding to the bilinear form :
\[ A[i,j] = a(λᵤ[j], λᵥ[i])\]
where λᵤ
and λᵥ
are the shape functions associated with the trial U
and the test V
function spaces respectively. In a "map-over" version, it can be written :
\[ A = a(blockU, blockV)\]
where blockU
and blockV
correspond formally to the lazy-map-over matrices :
\[ ∀k, blockU[k,j] = λᵤ[j], blockV[i,k] = λᵥ[i]\]
Bcube.blockmap_shape_functions
— MethodDev note :
Materialize the integrand function on all the different possible Tuples of v=(v1,0,0,...), (0,v2,0,...), ..., (..., vi, ...)
Bcube.blockmap_shape_functions
— Methodblockmap_shape_functions(fespace::AbstractFESpace, cellinfo::AbstractCellInfo)
Return all shape functions a = LazyMapOver((λ₁, λ₂, …, λₙ))
corresponding to fespace
in cell cellinfo
. These shape functions are wrapped by a LazyMapOver
so that for a function f
it gives: f(a) == map(f, a)
Bcube.blockmap_shape_functions
— Methodblockmap_shape_functions(multiFESpace::AbstractMultiFESpace, cellinfo::AbstractCellInfo)
Return all shape functions corresponding to each fespace
in multiFESSpace
for cell cellinfo
:
\[ ((v₁, ∅, ∅, …), (∅, v₂, ∅, …), …, ( …, ∅, ∅, vₙ))\]
where:
- vᵢ = (λᵢ₁, λᵢ₂, …, λᵢ_ₘ) are the shapes functions of the i-th fespace in the cell.
- ∅ are
NullOperator
s
Note that the LazyMapOver
is used to wrap recursively the result.
Bcube.compute
— Methodcompute(integration::Integration)
Compute an integral, independently from a FEM/DG framework (i.e without FESpace)
Return a SparseVector
. The indices of the domain elements are used to store the result of the integration in this sparse vector.
Example
Compute volume of each cell and each face.
mesh = rectangle_mesh(2, 3)
dΩ = Measure(CellDomain(mesh), 1)
dΓ = Measure(BoundaryFaceDomain(mesh), 1)
f = PhysicalFunction(x -> 1)
@show Bcube.compute(∫(f)dΩ)
@show Bcube.compute(∫(side⁻(f))dΓ)
DofHandler
Bcube.DofHandler
— TypeThe DofHandler
handles the degree of freedom numbering. To each degree of freedom is associated a unique integer.
Bcube.DofHandler
— MethodDofHandler(mesh::Mesh, fSpace::AbstractFunctionSpace, ncomponents::Int, isContinuous::Bool)
Constructor of a DofHandler for a SingleFESpace
on a Mesh
.
Bcube._deal_with_dofs_on_edges!
— Methoddeal_with_dofs_on_edges!(dict, iglob, offset, c2n, celltypes, icell::Int, inodes_g, e2n_g, s::AbstractShape, kvar::Int, fs)
Function dealing with dofs shared by different cell through an edge connection (excluding bord vertices).
TODO : remove kvar
Arguments
dict
may be modified by this routineiglob
may be modified by this routineoffset
may be modified by this routinefs
: FunctionSpace of varkvar
icell
: cell indexkvar
: var indexs
: shape oficell
-th cellinodes_g
: global indices of nodes oficell
Bcube._deal_with_dofs_on_faces!
— MethodTODO : remove kvar
Arguments
- f2n_g : local face index -> global nodes indices
Bcube._deal_with_dofs_on_vertices!
— Methoddeal_with_dofs_on_vertices!(dict, iglob, offset, icell::Int, inodes_g, s::AbstractShape, kvar::Int, fs)
Function dealing with dofs shared by different cell through a vertex connection.
TODO : remove kvar
Arguments
dict
may be modified by this routineiglob
may be modified by this routineoffset
may be modified by this routinefs
: FunctionSpace of varkvar
icell
: cell indexkvar
: var indexs
: shape oficell
-th cellinodes_g
: global indices of nodes oficell
Bcube.get_dof
— Methodget_dof(dhl::DofHandler, icell, icomp::Int, idof::Int)
Global index of the idof
local degree of freedom of component icomp
in cell icell
.
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1)))
@show get_dof(dhl, 1, 1, 1)
Bcube.get_dof
— Methodget_dof(dhl::DofHandler, icell, icomp::Int)
Global indices of all the dofs of a given component in a given cell
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1)))
@show get_dof(dhl, 1, 1)
Bcube.get_ncomponents
— MethodNumber of components handled by a DofHandler
Bcube.get_ndofs
— Methodget_ndofs(dhl, icell, icomp::Vector{Int})
Number of dofs for a given set of components in a given cell.
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1); size = 2))
@show get_ndofs(dhl, 1, [1, 2])
Bcube.get_ndofs
— Methodget_ndofs(dhl, icell, kvar::Int)
Number of dofs for a given variable in a given cell.
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1)))
@show get_ndofs(dhl, 1, 1)
Bcube.get_ndofs
— Methodget_ndofs(dhl::DofHandler, icell)
Number of dofs for a given cell.
Note that for a vector variable, the total (accross all components) number of dofs is returned.
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1)))
@show get_ndofs(dhl, 1, :u)
Bcube.get_ndofs
— Methodget_ndofs(dhl::DofHandler)
Total number of dofs. This function takes into account that dofs can be shared by multiple cells.
Example
mesh = one_cell_mesh(:line)
dhl = DofHandler(mesh, Variable(:u, FunctionSpace(:Lagrange, 1)))
@show get_ndofs(dhl::DofHandler)
Bcube.max_ndofs
— Methodmax_ndofs(dhl::DofHandler)
Count maximum number of dofs per cell, all components mixed