Assemblers
This section describes the assemblers that are currently available and their abstract interface.
All assemblers must possess at minimum a DofManager.
The assemblers possess all the baggage to use the internal method sparse! in SparseArrays.jl. This method allows for a zero allocation instantiation of a SparseMatrixCSC type on the CPU. There are also methods available to ease in the conversion of CSC types and other sparse types such as CSR.
On the GPU, however, this type is first converted to an appropriate COO matrix type on the desired backend. There is unfortunately not a unified sparse matrix API in julia for GPUs, so we implement this functionality in package extensions. On CUDA, for example, the operational sequence to get a CuSparseMatrixCSC is to first sort the COO (row, col, val) triplets so they are ordered by row and then column. Then a CuSparseMatrixCOO type is created and converted to a CuSparseMatrixCSC type via CUDA.jl methods. An identical approach is taken for RocM types.
NOTE: This is one of the most actively developed areas of the package. Please use caution with any method beginning with a "_" as these are internal methods that will change without notice.
Matrix Diagonal
FiniteElementContainers.assemble_diagonal! — Method
assemble_diagonal!(assembler, func::Function, Uu, p)
Assemble the diagonal of a matrix (e.g., stiffness or mass) into the assembler's residual storage field. Uses the same physics function as assemble_stiffness! / assemble_mass! but extracts only the diagonal entries of the element matrix at each quadrature point, avoiding full sparse matrix assembly.
The result is stored in the assembler's residual storage and can be retrieved via diagonal(asm).
This is essential for GPU-friendly Jacobi preconditioning: diag(K) gives the true diagonal, whereas the row-sum approximation K·1 can be zero at interior nodes of uniform meshes.
FiniteElementContainers.diagonal — Method
diagonal(
asm::FiniteElementContainers.AbstractAssembler
) -> Any
Return the assembled diagonal vector after a call to assemble_diagonal!. Extracts the unknown DOF subset (non-condensed) or the full vector (condensed).
Matrices
FiniteElementContainers.assemble_matrix! — Method
assemble_matrix!(
storage,
pattern,
dof,
func::Function,
Uu,
p;
use_inplace_methods
)
Note this is hard coded to storing the assembled sparse matrix in the stiffness_storage field of assembler.
Matrix Action
FiniteElementContainers.assemble_matrix_action! — Method
assemble_matrix_action!(
storage,
pattern,
dof,
func,
Uu,
Vu,
p;
use_inplace_methods
) -> Any
FiniteElementContainers.assemble_matrix_action! — Method
assemble_matrix_action!(
assembler,
func::Function,
Uu,
Vu,
p
)
FiniteElementContainers.assemble_matrix_free_action! — Method
assemble_matrix_free_action!(
storage,
pattern,
dof,
func_action,
Uu,
Vu,
p
) -> Any
FiniteElementContainers.assemble_matrix_free_action! — Method
assemble_matrix_free_action!(
assembler,
func_action::Function,
Uu,
Vu,
p
)
Matrix-free action assembly. func_action has signature funcaction(physics, interps, xel, t, Δt, uel, uelold, vel, stateoldq, statenewq, propsel) → SVector and returns the element-level product Kq·v_el directly, avoiding formation of the full element stiffness/mass matrix.
Matrix Action
Scalar
FiniteElementContainers.assemble_quadrature_quantity! — Method
assemble_quadrature_quantity!(
storage::L2Field,
pattern,
dof,
func::Function,
Uu,
p
) -> Any
assemble_quadrature_quantity!(
storage::L2Field,
pattern,
dof,
func::Function,
Uu,
p,
return_type::FiniteElementContainers.AssembledReturnType
) -> Any
FiniteElementContainers.assemble_quadrature_quantity! — Method
assemble_quadrature_quantity!(
storage::NamedTuple,
pattern,
dof,
func::Function,
Uu,
p
)
assemble_quadrature_quantity!(
storage::NamedTuple,
pattern,
dof,
func::Function,
Uu,
p,
return_type::FiniteElementContainers.AssembledReturnType
)
FiniteElementContainers.assemble_scalar! — Method
assemble_scalar!(assembler, func::Function, Uu, p) -> Any
Sources
FiniteElementContainers.assemble_vector_source! — Method
assemble_vector_source!(assembler, Uu, p)
Assemble body force contributions into the residual vector.
Vector
FiniteElementContainers.assemble_vector! — Method
assemble_vector!(
storage,
pattern,
dof,
func::Function,
Uu,
p;
use_inplace_methods,
use_sparse_vector
)
FiniteElementContainers.assemble_vector! — Method
assemble_vector!(assembler, func::Function, Uu, p)
Weakly enforced BCs
FiniteElementContainers._assemble_block_vector_weakly_enforced_bc! — Method
_assemble_block_vector_weakly_enforced_bc!(
field::FiniteElementContainers.AbstractField,
U::FiniteElementContainers.AbstractField,
X::FiniteElementContainers.AbstractField,
conns,
ref_fe,
sides,
vals
)
TODO this method needs to be updated to allow for sparse vectors TODO this method is hardcoded to H1 or L2 spaces
FiniteElementContainers.assemble_vector_neumann_bc! — Method
assemble_vector_neumann_bc!(assembler, Uu, p)
FiniteElementContainers.assemble_vector_robin_bc! — Method
assemble_vector_robin_bc!(assembler, Uu, p)
FiniteElementContainers.assemble_vector_weakly_enforced_bc! — Method
assemble_vector_weakly_enforced_bc!(
storage,
dof,
U,
X,
bcs
) -> Any
Abstract Interface
FiniteElementContainers.AbstractAssembler — Type
abstract type AbstractAssembler{Dof<:DofManager}FiniteElementContainers._cell_interpolants — Method
_cell_interpolants(
ref_fe::ReferenceFiniteElements.ReferenceFE,
q::Int64
) -> Any
FiniteElementContainers._element_level_fields — Method
_element_level_fields(
U::H1Field{T, D, NF},
ref_fe,
conns,
e
) -> Any
FiniteElementContainers._element_level_fields_flat — Method
_element_level_fields_flat(
U::H1Field{T, D, NF},
ref_fe,
conns,
e
) -> Any
FiniteElementContainers._element_level_properties — Method
_element_level_properties(
props::AbstractArray,
_::Int64
) -> StaticArraysCore.SVector
FiniteElementContainers._element_level_properties — Method
_element_level_properties(
props::StaticArraysCore.SArray{Tuple{NP}, T, 1, NP},
_::Int64
) -> StaticArraysCore.SVector
FiniteElementContainers._element_scratch — Method
_element_scratch(
_::FiniteElementContainers.AssembledScalar,
ref_fe,
U
) -> Any
FiniteElementContainers._element_scratch — Method
_element_scratch(
_::FiniteElementContainers.AssembledStruct,
ref_fe,
U
)
FiniteElementContainers._element_scratch — Method
_element_scratch(
_::FiniteElementContainers.AssembledMatrix,
ref_fe,
U::H1Field{T, D, NF}
) -> Any
FiniteElementContainers._element_scratch — Method
_element_scratch(
_::FiniteElementContainers.AssembledVector,
ref_fe,
U::H1Field{T, D, NF}
) -> Any
FiniteElementContainers._quadrature_level_state — Method
_quadrature_level_state(
state::AbstractArray{<:Number, 3},
q::Int64,
e::Int64
) -> Any
FiniteElementContainers._surface_interpolants — Method
_surface_interpolants(
ref_fe::ReferenceFiniteElements.ReferenceFE,
q::Int64,
side::Int64
) -> Any
FiniteElementContainers.hessian — Method
hessian(
asm::FiniteElementContainers.AbstractAssembler
) -> Any
The Hessian operator currently aliases the assembled stiffness — no FE code path in this repository writes a separate Hessian buffer. If a future caller needs H ≠ K, add a dedicated hessian_storage field on the assembler and specialize this method accordingly.
FiniteElementContainers.hvp — Method
hvp(
asm::FiniteElementContainers.AbstractAssembler,
v
) -> Any
FiniteElementContainers.mass — Method
mass(asm::FiniteElementContainers.AbstractAssembler) -> Any
FiniteElementContainers.residual — Method
residual(
asm::FiniteElementContainers.AbstractAssembler;
use_sparse_vector
) -> Any
assumes assemble_vector! has already been called
FiniteElementContainers.stiffness — Method
stiffness(
asm::FiniteElementContainers.AbstractAssembler
) -> Any
KernelAbstractions.get_backend — Method
get_backend(
asm::FiniteElementContainers.AbstractAssembler
) -> Any
Assembly methods that are safe with Enzyme.jl
FiniteElementContainers._assemble_block_enzyme_safe! — Method
_assemble_block_enzyme_safe!(
_::KernelAbstractions.CPU,
field::FiniteElementContainers.AbstractField,
conns::AbstractArray,
coffset::Int64,
func::Function,
physics::AbstractPhysics,
ref_fe::ReferenceFiniteElements.ReferenceFE,
X::FiniteElementContainers.AbstractField,
t::Number,
dt::Number,
U::FiniteElementContainers.AbstractField,
U_old::FiniteElementContainers.AbstractField,
state_old,
state_new,
props::AbstractArray,
return_type::FiniteElementContainers.AssembledReturnType
)
Assembly method for a block labelled as block_id. This is a CPU implementation with no threading.
TODO add state variables and physics properties
FiniteElementContainers.assemble_vector_enzyme_safe! — Method
assemble_vector_enzyme_safe!(
assembler,
func::Function,
Uu,
p
)
FiniteElementContainers.assemble_vector_enzyme_safe! — Method
assemble_vector_enzyme_safe!(
storage::FiniteElementContainers.AbstractField,
pattern,
dof,
func::Function,
Uu,
p
)
SparseMatrixAssembler
FiniteElementContainers.SparseMatrixAssembler — Type
struct SparseMatrixAssembler{Condensed, SparseMatrixType, UseInPlaceMethods, UseSparseVec, IV<:AbstractVector{Int64}, RV<:AbstractVector{Float64}, Var<:FiniteElementContainers.AbstractFunction, FieldStorage} <: FiniteElementContainers.AbstractAssembler{DofManager{Condensed, Int64, IV<:AbstractVector{Int64}, Var<:FiniteElementContainers.AbstractFunction}}General sparse matrix assembler that can handle first or second order problems in time.
FiniteElementContainers.SparseMatrixAssembler — Method
SparseMatrixAssembler(
dof::DofManager;
sparse_matrix_type,
use_inplace_methods,
use_sparse_vector,
matrix_free
) -> Union{SparseMatrixAssembler{_A, _B, _C, _D, var"#s177", Vector{Float64}, var"#s1771", H1Field{Float64, Vector{Float64}, _A1}} where {NF, var"#s1771"<:FiniteElementContainers.AbstractFunction{F, NF}, _A1}, SparseMatrixAssembler{_A, _B, _C, _D, var"#s177", Vector{Float64}, Var, H1Field{Float64, Vector{Float64}, _A1}} where {NF, Var<:FiniteElementContainers.AbstractFunction{F, NF}, _A1}} where {_A, _B, _C, _D, var"#s177"<:AbstractVector{Int64}, F<:FunctionSpace}
Construct a SparseMatrixAssembler for a specific field type, e.g. H1Field. Can be used to create block arrays for mixed FEM problems.