FiniteElementContainers.AbstractFieldType
abstract type AbstractField{T, N, NF, Vals} <: AbstractArray{T, N}

Thin wrapper that subtypes AbstractArray and serves as the base Field type

FiniteElementContainers.DynamicAssemblerType
struct DynamicAssembler{Rtype, Itype, I<:AbstractArray{Itype, 1}, J<:AbstractArray{Itype, 1}, U<:AbstractArray{Itype, 1}, Sizes<:AbstractArray{Itype, 1}, Offsets<:AbstractArray{Itype, 1}, R<:NodalField, K<:AbstractArray{Rtype, 1}, M<:AbstractArray{Rtype, 1}, C1, C2, C3, C4, C5, C6, C7} <: Assembler{Rtype, Itype}
  • Is::AbstractVector

  • Js::AbstractVector

  • unknown_dofs::AbstractVector

  • block_sizes::AbstractVector

  • block_offsets::AbstractVector

  • residuals::NodalField

  • stiffnesses::AbstractVector

  • masses::AbstractVector

  • klasttouch::Any

  • csrrowptr::Any

  • csrcolval::Any

  • csrnzval::Any

  • csccolptr::Any

  • cscrowval::Any

  • cscnzval::Any

Assembler for dynamic problems without damping

Provides both a mass and stiffness matrix

FiniteElementContainers.ElementFieldType
abstract type ElementField{T, N, NN, NE, Vals} <: FiniteElementContainers.AbstractField{T, N, NN, Vals}

Abstract type for implementations of fields that live on elements.

Constructors

ElementField{NN, NE, Vector}(vals::M) where {NN, NE, M <: AbstractArray{<:Number, 2}}

ElementField{NN, NE, Matrix}(vals::M) where {NN, NE, M <: AbstractArray{<:Number, 2}}

ElementField{NN, NE, Vector}(vals::V) where {NN, NE, V <: AbstractArray{<:Number, 1}}

ElementField{NN, NE, Vector, T}(::UndefInitializer) where {NN, NE, T}

ElementField{NN, NE, Matrix, T}(::UndefInitializer) where {NN, NE, T <: Number}

ElementField{NN, NE, StructArray, T}(::UndefInitializer) where {NN, NE, T}

ElementField{Tup, A, T}(::UndefInitializer) where {Tup, A, T}

ElementField{Tup, A}(vals::M) where {Tup, A, M <: AbstractArray}

FiniteElementContainers.FileMeshType
struct FileMesh{MeshObj} <: FiniteElementContainers.AbstractMesh
  • file_name::String

  • mesh_obj::Any

Mesh type that has a handle to an open mesh file object. This type's methods are "overridden" in extensions.

See FiniteElementContainersExodusExt for an example.

FiniteElementContainers.NodalFieldType
abstract type NodalField{T, N, NF, NN, Vals} <: FiniteElementContainers.AbstractField{T, N, NF, Vals}

Abstract type for implementations of fields that live on nodes.

Constructors

NodalField{NF, NN, Vector}(vals::M) where {NF, NN, M <: AbstractArray{<:Number, 2}}

NodalField{NF, NN, Matrix}(vals::M) where {NF, NN, M <: AbstractArray{<:Number, 2}}

NodalField{NF, NN, Vector}(vals::V) where {NF, NN, V <: AbstractArray{<:Number, 1}}

NodalField{NF, NN, Vector}(vals::V) where {NF, NN, V <: AbstractArray{<:Number, 1}}

NodalField{NF, NN, Vector, T}(::UndefInitializer) where {NF, NN, T}

NodalField{NF, NN, Matrix, T}(::UndefInitializer) where {NF, NN, T <: Number}

NodalField{NF, NN, StructArray, T}(::UndefInitializer) where {NF, NN, T}

NodalField{Tup, A, T}(::UndefInitializer) where {Tup, A, T}

FiniteElementContainers.NonAllocatedFunctionSpaceType
struct NonAllocatedFunctionSpace{NDof, Conn<:(ElementField), DofConn<:(ElementField), RefFE<:ReferenceFiniteElements.ReferenceFE} <: FunctionSpace{NDof, Conn<:(ElementField), RefFE<:ReferenceFiniteElements.ReferenceFE}
  • conn::ElementField

  • dof_conn::ElementField

  • ref_fe::ReferenceFiniteElements.ReferenceFE

FiniteElementContainers.QuadratureFieldType
abstract type QuadratureField{T, N, NF, NQ, NE, Vals} <: FiniteElementContainers.AbstractField{T, N, NF, Vals}

Abstract type for implementations of fields that live on quadrature points.

Constructors

QuadratureField{NF, NQ, NE, Matrix}(vals::Matrix{<:Number}) where {NF, NQ, NE}

QuadratureField{NF, NQ, NE, Vector}(vals::Matrix{<:Number}) where {NF, NQ, NE}

QuadratureField{NF, NQ, NE, Matrix, T}(::UndefInitializer) where {NF, NQ, NE, T}

QuadratureField{NF, NQ, NE, StructArray, T}(::UndefInitializer) where {NF, NQ, NE, T}

QuadratureField{NF, NQ, NE, StructVector, T}(::UndefInitializer) where {NF, NQ, NE, T}

QuadratureField{NF, NQ, NE, Vector, T}(::UndefInitializer) where {NF, NQ, NE, T}

QuadratureField{Tup, A, T}(::UndefInitializer) where {Tup, A, T}

QuadratureField{Tup, A}(vals::M) where {Tup, A, M <: AbstractArray}

FiniteElementContainers.StaticAssemblerType
struct StaticAssembler{Rtype, Itype, I<:AbstractArray{Itype, 1}, J<:AbstractArray{Itype, 1}, U<:AbstractArray{Itype, 1}, Sizes<:AbstractArray{Itype, 1}, Offsets<:AbstractArray{Itype, 1}, R<:NodalField, K<:AbstractArray{Rtype, 1}, C1, C2, C3, C4, C5, C6, C7} <: Assembler{Rtype, Itype}
  • Is::AbstractVector

  • Js::AbstractVector

  • unknown_dofs::AbstractVector

  • block_sizes::AbstractVector

  • block_offsets::AbstractVector

  • residuals::NodalField

  • stiffnesses::AbstractVector

  • klasttouch::Any

  • csrrowptr::Any

  • csrcolval::Any

  • csrnzval::Any

  • csccolptr::Any

  • cscrowval::Any

  • cscnzval::Any

Assembler for static or quasistatic problems where only a stiffness matrix is necessary

FiniteElementContainers.VectorizedPreAllocatedFunctionSpaceType
struct VectorizedPreAllocatedFunctionSpace{NDof, Conn<:(ElementField), DofConn<:(ElementField), RefFE<:ReferenceFiniteElements.ReferenceFE, V1<:QuadratureField, V2<:QuadratureField, V3<:QuadratureField} <: FunctionSpace{NDof, Conn<:(ElementField), RefFE<:ReferenceFiniteElements.ReferenceFE}
  • conn::ElementField

  • dof_conn::ElementField

  • ref_fe::ReferenceFiniteElements.ReferenceFE

  • Ns::QuadratureField

  • ∇N_Xs::QuadratureField

  • JxWs::QuadratureField

Base.eltypeMethod
eltype(
    _::Type{FiniteElementContainers.AbstractField{T, N, NF, Vals}}
) -> Any
Base.similarMethod
similar(
    asm::StaticAssembler
) -> StaticAssembler{_A, Int64, var"#s178", var"#s1781", var"#s1782", var"#s1783", var"#s1784", _B, <:AbstractVector{_A}} where {_A, var"#s178"<:AbstractVector{Int64}, var"#s1781"<:AbstractVector{Int64}, var"#s1782"<:AbstractVector{Int64}, var"#s1783"<:AbstractVector{Int64}, var"#s1784"<:AbstractVector{Int64}, _B<:NodalField}
Base.similarMethod
similar(
    field::FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals}
) -> FiniteElementContainers.SimpleElementField
Base.similarMethod
similar(
    field::FiniteElementContainers.SimpleNodalField{T, N, NF, NN, Vals}
) -> FiniteElementContainers.SimpleNodalField
Base.similarMethod
similar(
    field::FiniteElementContainers.SimpleQuadratureField{T, N, NF, NQ, NE, Vals}
) -> FiniteElementContainers.SimpleQuadratureField
Base.similarMethod
similar(
    field::FiniteElementContainers.VectorizedElementField{T, N, NN, NE, Vals}
) -> FiniteElementContainers.VectorizedElementField
Base.zeroMethod
zero(
    _::FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals<:(StructArrays.StructVector)}
) -> FiniteElementContainers.SimpleElementField{_A, _B, _C, _D, Vals} where {_A, _B, _C, _D, Vals<:(StructArrays.StructVector)}
Base.zeroMethod
zero(
    field::FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals<:AbstractArray}
) -> FiniteElementContainers.SimpleElementField{_A, _B, _C, _D, Vals} where {_A, _B, _C, _D, Vals<:(StructArrays.StructVector)}
Base.zeroMethod
zero(
    _::FiniteElementContainers.SimpleNodalField{T, N, NF, NN, Vals}
) -> FiniteElementContainers.SimpleNodalField
Base.zeroMethod
zero(
    field::FiniteElementContainers.SimpleQuadratureField{T, N, NF, NQ, NE, Vals}
) -> FiniteElementContainers.SimpleQuadratureField
Base.zeroMethod
zero(
    field::FiniteElementContainers.VectorizedElementField{T, N, NN, NE, Vals}
) -> FiniteElementContainers.VectorizedElementField
Base.zeroMethod
zero(
    field::FiniteElementContainers.VectorizedNodalField{T, N, NF, NN, Vals}
) -> FiniteElementContainers.VectorizedNodalField
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals<:(AbstractMatrix)}}
)
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals<:(StructArrays.StructVector)}}
) -> FiniteElementContainers.SimpleElementField{_A, _B, _C, _D, Vals} where {_A, _B, _C, _D, Vals<:(StructArrays.StructVector)}
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleElementField{T, N, NN, NE, Vals<:StructArrays.StructArray}}
) -> FiniteElementContainers.SimpleElementField{_A, _B, _C, _D, Vals} where {_A, _B, _C, _D, Vals<:(StructArrays.StructVector)}
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleNodalField{T, N, NF, NN, Vals}}
) -> FiniteElementContainers.SimpleNodalField
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleQuadratureField{T, N, NF, NQ, NE, Vals<:(AbstractMatrix)}}
) -> FiniteElementContainers.SimpleQuadratureField{_A, _B, _C, _D, _E, Vals} where {_A, _B, _C, _D, _E, Vals<:(StructArrays.StructArray{_A, 2})}
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.SimpleQuadratureField{T, N, NF, NQ, NE, Vals<:StructArrays.StructArray}}
) -> FiniteElementContainers.SimpleQuadratureField{_A, _B, _C, _D, _E, Vals} where {_A, _B, _C, _D, _E, Vals<:(StructArrays.StructArray{_A, 2})}
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.VectorizedElementField{T, N, NN, NE, Vals<:(AbstractVector)}}
) -> FiniteElementContainers.VectorizedElementField
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.VectorizedNodalField{T, N, NF, NN, Vals<:(AbstractVector{<:Number})}}
) -> FiniteElementContainers.VectorizedNodalField{_A, _B, _C, _D, Vals} where {_A<:Number, _B, _C, _D, Vals<:(StructArrays.StructVector{_A})}
Base.zeroMethod
zero(
    _::Type{FiniteElementContainers.VectorizedNodalField{T, N, NF, NN, Vals<:(StructArrays.StructVector)}}
) -> FiniteElementContainers.VectorizedNodalField{_A, _B, _C, _D, Vals} where {_A, _B, _C, _D, Vals<:(StructArrays.StructVector{_A})}
FiniteElementContainers.assemble!Method
assemble!(
    assembler::DynamicAssembler,
    dof::DofManager,
    fspace::FunctionSpace,
    X,
    U,
    block_id,
    residual_func,
    tangent_func,
    mass_func
)

Simple method for assembling in serial

FiniteElementContainers.assemble!Method
assemble!(
    assembler::DynamicAssembler,
    K_el::AbstractMatrix{<:Number},
    M_el::AbstractMatrix{<:Number},
    block_id::Int64,
    el_id::Int64
) -> AbstractMatrix{<:Number}

assembly for stiffness matrix

FiniteElementContainers.assemble!Method
assemble!(
    assembler::StaticAssembler,
    K_el::AbstractMatrix,
    block_id::Int64,
    el_id::Int64
) -> AbstractMatrix

assembly for stiffness matrix

FiniteElementContainers.assemble!Method
assemble!(
    R::AbstractArray{<:Number},
    Kv::AbstractArray{<:Number},
    R_el::AbstractArray{<:Number},
    Kv_el::AbstractArray{<:Number},
    conn::AbstractArray{<:Integer}
)

generic assembly method that directly goes into a vector for doing a residual and matrix vector product at once

FiniteElementContainers.assemble!Method
assemble!(
    R::AbstractArray{<:Number},
    R_el::AbstractArray{<:Number},
    conn::AbstractArray{<:Integer}
)

generic assembly method that directly goes into a vector

FiniteElementContainers.assemble_atomic!Method
assemble_atomic!(
    R::NodalField,
    R_el::AbstractVector{<:Number},
    conn::AbstractVector{<:Integer}
)

assembly method for just a residual vector

TODO need to add an Atomix lock here TODO add block_id to fspace or something like that

FiniteElementContainers.modify_field_gradientsMethod
modify_field_gradients(
    _::IncompressiblePlaneStress,
    ∇u_q::StaticArraysCore.SArray{Tuple{2, 2}, T<:Number, 2, 4},
    _::Type{<:StaticArraysCore.SArray}
) -> StaticArraysCore.SArray
FiniteElementContainers.modify_field_gradientsMethod
modify_field_gradients(
    _::IncompressiblePlaneStress,
    ∇u_q::StaticArraysCore.SArray{Tuple{2, 2}, T<:Number, 2, 4},
    _::Type{<:Tensors.Tensor}
) -> Tensors.Tensor{2, 3, _A, 9} where _A
FiniteElementContainers.modify_field_gradientsMethod

To deprecate or not to deprecate?

modify_field_gradients(
    form::PlaneStrain,
    ∇u_q::StaticArraysCore.SArray{Tuple{2, 2}, T<:Number, 2, 4};
    type
) -> Tensors.Tensor{2, 3, _A, 9} where _A
FiniteElementContainers.update_unknown_dofs!Method
update_unknown_dofs!(
    dof_manager::DofManager,
    dofs::AbstractVector{<:Integer}
)

Method when all dofs are updated at once

First it resets all dofs to unknowns, then one by one sets dofs to bcs in dofs

Assumes there's a unique set of nodes provided

FiniteElementContainers.update_unknown_dofs!Method
update_unknown_dofs!(
    assembler::Union{DynamicAssembler, StaticAssembler},
    dof,
    fspaces,
    nodes_in::AbstractVector{<:Integer}
)

method that assumes first dof TODO move sorting of nodes up stream TODO remove other scratch unknowns and unknown_dofs arrays

FiniteElementContainers.volumeMethod
volume(
    fspace::FunctionSpace,
    _::ReferenceFiniteElements.ReferenceFEType,
    X::NodalField,
    q::Int64,
    e::Int64
) -> Any
SparseArrays.sparse!Method
sparse!(
    assembler::DynamicAssembler
) -> Tuple{SparseArrays.SparseMatrixCSC, SparseArrays.SparseMatrixCSC}
SparseArrays.sparse!Method
sparse!(
    assembler::StaticAssembler
) -> SparseArrays.SparseMatrixCSC
SparseArrays.sparseMethod
sparse(
    assembler::DynamicAssembler
) -> Tuple{SparseArrays.SparseMatrixCSC, SparseArrays.SparseMatrixCSC}
SparseArrays.sparseMethod
sparse(
    assembler::StaticAssembler
) -> SparseArrays.SparseMatrixCSC