Types

Contents

Coordinate systems

FinEtools.CSysModule.CSysType
CSys{T<:Number, F<:Function}

Type for coordinate system transformations. Used to define material coordinate systems, and output coordinate systems, for instance.

source
FinEtools.CSysModule.CSysMethod
CSys(dim::IT) where {IT}

Construct coordinate system when the rotation matrix is the identity.

dim = is the space dimension.

source
FinEtools.CSysModule.CSysMethod
CSys(sdim::IT1, mdim::IT2, z::T, computecsmat::F) where {IT1, IT2, T <: Number, F <: Function}

Construct coordinate system when the function to compute the rotation matrix of type T is given.

  • z = zero value,
  • The computecsmat function signature: update!(csmatout::Matrix{T}, XYZ::VecOrMat{T}, tangents::Matrix{T}, feid::IT, qpid::IT) where {T, IT} where
    • csmatout= output matrix buffer, of size (sdim, mdim);
    • XYZ= location in physical coordinates;
    • tangents= tangent vector matrix, tangents to the parametric coordinate curves in the element;
    • feid= finite element identifier;
    • qpid= quadrature point identifier.

Example

# Cylindrical coordinate system: NO ALLOCATIONS WHATSOEVER!
@views function compute!(csmatout, XYZ, tangents, feid, qpid)
    center = (0.0, 0.0, 0.0)
    xyz = (XYZ[1], XYZ[2], XYZ[3])
    csmatout[:, 1] .= xyz .- center
    csmatout[3, 1] = 0.0
    csmatout[:, 1] ./= norm(csmatout[:, 1])
    csmatout[:, 3] .= (0.0, 0.0, 1.0)
    cross3!(csmatout[:, 2], csmatout[:, 3], csmatout[:, 1])
    csmatout[:, 2] ./=  norm(csmatout[:, 2])
    return csmatout
end
source
FinEtools.CSysModule.CSysMethod
CSys(sdim::IT1, mdim::IT2, computecsmat::F) where {F <: Function, IT1, IT2}

Construct coordinate system when the function to compute the rotation matrix is given.

The function signature:

update!(csmatout::Matrix{T}, XYZ::VecOrMat{T}, tangents::Matrix{T},
    feid::IT, qpid::IT) where {T, IT}

where

  • csmatout= output matrix buffer, of size (sdim, mdim)
  • XYZ= location in physical coordinates,
  • tangents= tangent vector matrix, tangents to the parametric coordinate curves in the element,
  • feid= finite element identifier;
  • qpid= quadrature point identifier.

Example

# Cylindrical coordinate system: NO ALLOCATIONS WHATSOEVER!
@views function compute!(csmatout, XYZ, tangents, feid, qpid)
    center = (0.0, 0.0, 0.0)
    xyz = (XYZ[1], XYZ[2], XYZ[3])
    csmatout[:, 1] .= xyz .- center
    csmatout[3, 1] = 0.0
    csmatout[:, 1] ./= norm(csmatout[:, 1])
    csmatout[:, 3] .= (0.0, 0.0, 1.0)
    cross3!(csmatout[:, 2], csmatout[:, 3], csmatout[:, 1])
    csmatout[:, 2] ./=  norm(csmatout[:, 2])
    return csmatout
end
source
FinEtools.CSysModule.CSysMethod
CSys(sdim::IT1, mdim::IT2) where {IT1<:Integer, IT2<:Integer}

Construct coordinate system for isotropic-material used with isoparametric finite elements.

  • sdim = number of space dimensions,
  • mdim = number of manifold dimensions of the finite element in which the coordinate system is being evaluated.
Note

If the coordinate system matrix should be identity, better use the constructor for this specific situation, CSys(dim). That will be much more efficient.

See also

gen_iso_csmat

source
FinEtools.CSysModule.CSysMethod
CSys(dim, z::T) where {T}

Construct coordinate system when the rotation matrix of element type T is the identity.

dim = is the space dimension.

source

Data cache

FinEtools.DataCacheModule.DataCacheType
DataCache{D, F<:Function}

Type for caching data, such as vectors, matrices, and numbers.

D = type of the data, for instance Matrix{Float64} or Float32. F = type of the function to update the entries of the array.

Signature of the function to fill the cache with the value of the array is as follows:

function fillcache!(cacheout::D,
    XYZ::VecOrMat{T}, tangents::Matrix{T}, feid::IT, qpid::IT) where {D, T, IT}
    ... # modify the value of cacheout
    return cacheout
end

It may use the location XYZ, it may use the columns of the Jacobian matrix of the element, tangents, it may also choose the value of the finite element identifier (i.e. serial number), feid, and the identifier (i.e. serial number) of the quadrature point, qpid. All of these values are supplied by the code requesting the value of the cache. It must return the modified argument cacheout.

When the cache is accessed, the callback fillcache! is called, and the output cacheout is filled with the value of the cached data.

Example

function fillcache!(cacheout::Array{CT, N},
        XYZ::VecOrMat{T}, tangents::Matrix{T},
        feid::IT, qpid::IT) where {CT, N, T, IT}
    cacheout .= LinearAlgebra.I(3)
    return cacheout
end
c = DataCache(zeros(Float32, 3, 3), fillcache!)
function f(c)
    XYZ, tangents, feid, qpid = (reshape([0.0, 0.0], 1, 2), [1.0 0.0; 0.0 1.0], 1, 1)
    data = c(XYZ, tangents, feid, qpid)
end
@test f(c) == LinearAlgebra.I(3)
Note

The point of the data cache is that there will be no copying of data. The cache data field is filled in and returned, but no data needs to be copied. The bad news is, the cache is not thread safe. Reading is okay, but writing can lead to data races.

source

Surface-normal utilities

FinEtools.SurfaceNormalModule.SurfaceNormalType
SurfaceNormal{F<:Function}

Exterior surface normal type.

The normal vector is assumed to be normalized to unit length.

Signature of the function to compute the value of the unit normal at any given point XYZ, using the columns of the Jacobian matrix of the element, tangents, the finite element label, feid, and the identifier of the quadrature point, qpid:

function computenormal!(normalout::Vector{CT}, XYZ::Matrix{T},
        tangents::Matrix{T}, feid::IT, qpid::IT) where {CT, T, IT}
        # Calculate the normal  and copy it into the buffer....
        return normalout
    end

The buffer normalout needs to be filled with the value of the normal vector.

Refer to DataCache for details on the caching.

source
FinEtools.SurfaceNormalModule.SurfaceNormalMethod
SurfaceNormal(ndimensions::FInt)

Construct surface normal evaluator when the default calculation of the normal vector based on the columns of the Jacobian matrix should be used.

The normal vector has ndimensions entries.

When the columns of the tangents array are parallel (or one of them is a zero vector), the normal cannot be normalized to unit length (it is a zero vector). In that case a zero vector is returned, and a warning is printed.

source
FinEtools.SurfaceNormalModule.SurfaceNormalMethod
SurfaceNormal(ndimensions, z::T, computenormal!::F) where {T<:Number, F<:Function}

Construct surface normal evaluator when the function to compute the normal vector is given.

Arguments

  • T = the type of the elements of the normal vector,

  • ndofn = number of components of the normal vector,

  • computenormal! = callback function. The function computenormal! needs to have a signature of

    function computenormal!(normalout::Vector{CT}, XYZ::Matrix{T}, tangents::Matrix{T}, feid::IT, qpid::IT) where {CT, T, IT} # Calculate the normal and copy it into the buffer.... return normalout end

    and it needs to fill in the buffer normalout with the current normal at the location XYZ, using, if appropriate, the information supplied in the Jacobian matrix tangents, the identifier of the finite element, feid, and the quadrature point id, qpid. Refer to DataCache.

source

Force intensity

FinEtools.ForceIntensityModule.ForceIntensityType
ForceIntensity{T<:Number, F<:Function}

Distributed force (force intensity) type.

The force intensity class. The physical units are force per unit volume, where volume depends on to which manifold the force is applied:

  • force/length^3 (when applied to a 3-D solid),
  • force/length^2 (when applied to a surface),
  • force/length^1 (when applied along a curve), or
  • force/length^0 (when applied at a point).

Signature of the function to compute the value of the force at any given point XYZ, using the columns of the Jacobian matrix of the element, tangents, the finite element identifier, feid:

getforce!(forceout::Vector{CT}, XYZ::Matrix{T}, tangents::Matrix{T}, feid::IT) where {CT, T, IT}

A DataCache is used to store the data.

source
FinEtools.ForceIntensityModule.ForceIntensityMethod
ForceIntensity(
    ::Type{T},
    ndofn,
    computeforce!::F,
) where {T<:Number, F<:Function}

Construct force intensity when the function to compute the intensity vector is given.

Arguments

  • T = the type of the elements of the force vector, typically floating-point or complex floating-point numbers,
  • ndofn = number of elements of the force vector (the length of the force vector),
  • computeforce! = callback function. The function computeforce! needs to have a signature of function computeforce!(forceout::Vector{CT}, XYZ::Matrix{T}, tangents::Matrix{T}, feid::IT) ) where {CT, T<:Number, IT<:Integer} # Calculate the force and copy it into the buffer forceout.... return forceout end and it needs to fill in the buffer forceout with the current force at the location XYZ, using, if appropriate, the information supplied in the Jacobian matrix tangents, and the identifier of the finite element, feid.
source

Finite element sets

FinEtools.FESetModule.AbstractFESet0ManifoldType
AbstractFESet0Manifold{NODESPERELEM} <: FESet{NODESPERELEM}

Abstract type of a finite element set for 0-dimensional manifolds (points). Parameterized with the number of of the nodes per element.

source
FinEtools.FESetModule.AbstractFESet1ManifoldType
AbstractFESet1Manifold{NODESPERELEM} <: FESet{NODESPERELEM}

Abstract type of a finite element set for 1-dimensional manifolds (curves). Parameterized with the number of of the nodes per element.

source
FinEtools.FESetModule.AbstractFESet2ManifoldType
AbstractFESet2Manifold{NODESPERELEM} <: FESet{NODESPERELEM}

Abstract type of a finite element set for 2-dimensional manifolds (surfaces). Parameterized with the number of of the nodes per element.

source
FinEtools.FESetModule.AbstractFESet3ManifoldType
AbstractFESet3Manifold{NODESPERELEM} <: FESet{NODESPERELEM}

Abstract type of a finite element set for 3-dimensional manifolds (solids). Parameterized with the number of of the nodes per element.

source

Finite element nodes

FinEtools.FENodeSetModule.FENodeSetType
mutable struct FENodeSet{T}

Finite element node set type.

The only field is xyz, as the array of node locations. Indexed with the node number. The location of node j is given by xyz[j,:]. Clearly, the nodes needs to be numbered between 1 and size(xyz, 1).

Note

The constructor makes a copy of the input xyz array for safety.

source

Finite element node-to-element map

FinEtools.FENodeToFEMapModule.FENodeToFEMapType
FENodeToFEMap

Map from finite element nodes to the finite elements connecting them.

For each node referenced in the connectivity of the finite element set on input, the numbers of the individual finite elements that reference that node is stored in an array in the array map.

    Example:
fes.conn= [7,6,5;
            4,1,3;
            3,7,5];
The map reads
    map[1] = [2];
    map[2] = [];#  note that node number 2 is not referenced by the connectivity
    map[3] = [2,3];
    map[4] = [2];
    map[5] = [1,3];
    map[6] = [1];
    map[7] = [1,3];

The individual elements from the connectivity that reference node number 5 are 1 and 3, so that fes.conn(map[5],:)includes all the nodes that are connected to node 5 (including node 5 itself).

source
FinEtools.FENodeToFEMapModule.FENodeToFEMapMethod
FENodeToFEMap(conns::FIntMat, nmax::FInt)

Map from finite element nodes to the finite elements connecting them.

  • conns = integer array of the connectivities
  • nmax = largest possible node number
source
FinEtools.FENodeToFEMapModule.FENodeToFEMapMethod
FENodeToFEMap(conn::Vector{NTuple{N, IT}}, nmax::FInt) where {N, IT<:Integer}

Map from finite element nodes to the finite elements connecting them.

  • conns = connectivities as a vector of tuples
  • nmax = largest possible node number

Example:

m = FENodeToFEMap(fes.conn, count(fens))
source

Fields

FinEtools.FieldModule.AbstractFieldType
AbstractField

Abstract field.

Expected attributes:

  • values::Array{T,2}: Array of degree of freedom parameters, indexed by entity number
  • dofnums::Array{IT,2}: Array of degree of freedom numbers, indexed by entity number
  • is_fixed::Matrix{Bool}: Array of Boolean flags, indexed by entity number
  • _nfreedofs::IT: Total number of free degrees of freedom

See also: @add_Field_fields() .

source
FinEtools.GeneralFieldModule.GeneralFieldMethod
GeneralField(data::Matrix{T}, zi::IT) where {T<:Number, IT<:Integer}

Constructor of general field. The values of the field are given by the array on input, data. This array needs to have as many rows as there are entities, and as many columns as there are degrees of freedom per entities.

The integer type for the storage of the degree of freedom numbers is set as that of the argument zi.

source
FinEtools.GeneralFieldModule.GeneralFieldMethod
GeneralField(data::Vector{T}) where {T<:Number}

Constructor of general field. The values of the field are given by the vector on input, data. This vector needs to have as many rows as there are entities.

source
FinEtools.NodalFieldModule.NodalFieldMethod
NodalField(data::Matrix{T}, zi::IT) where {T<:Number, IT<:Integer}

Constructor of nodal field. The values of the field are given by the array on input, data. This array needs to have as many rows as there are nodes, and as many columns as there are degrees of freedom per node.

The integer type for the storage of the degree of freedom numbers is set as that of the argument zi.

source
FinEtools.NodalFieldModule.NodalFieldMethod
NodalField(data::Vector{T}) where {T<:Number}

Constructor of nodal field. The values of the field are given by the vector on input, data. This vector needs to have as many entries as there are nodes; there is just one degree of freedom per node.

source
FinEtools.ElementalFieldModule.ElementalFieldType
ElementalField{T<:Number, IT<:Integer} <: AbstractField

Elemental field, meaning the entities are finite elements.

The values in the field are indexed by the element number. This means that there needs to be one field per finite element set.

source
FinEtools.ElementalFieldModule.ElementalFieldMethod

ElementalField(data::Matrix{T}, zi::IT) where {T<:Number, IT<:Integer}

Constructor of elemental field. The values of the field are given by the array on input, data. This array needs to have as many rows as there are elements, and as many columns as there are degrees of freedom per element.

The integer type for the storage of the degree of freedom numbers is set as that of the argument zi.

source
FinEtools.ElementalFieldModule.ElementalFieldMethod
ElementalField(data::Vector{T}) where {T<:Number}

Constructor of elemental field. The values of the field are given by the vector on input, data. This vector needs to have as many entries as there are elements; there is just one degree of freedom per element.

source

Integration rule

FinEtools.IntegRuleModule.NodalSimplexRuleType
NodalSimplexRule <: AbstractIntegRule

The nodal-quadrature simplex rule.

The rule is applicable for line segments, triangles, tetrahedra.

Note

The quadrature points for a nodal quadrature rule must be listed in the order in which the nodes are used in the definition of the element!

source
FinEtools.IntegRuleModule.NodalTensorProductRuleType
NodalTensorProductRule <: AbstractIntegRule

The tensor-product nodal-quadrature rule.

The rule is applicable for line segments, quadrilaterals, hexahedra.

Note

The quadrature points for a nodal quadrature rule must be listed in the order in which the nodes are used in the definition of the element!

source
FinEtools.IntegRuleModule.TriRuleType
TriRule(npts=1)

Type for triangular quadrature rule. Used for integration of the standard triangle, which is between 0 and 1 in both parametric coordinates. npts = number of points (1– one-point rule, 3 – three-point rule, 6 – six point rule, 9 –nine point rule, 10 – Strang 10 point, order 13, degree of precision 7, rule), 12 and 13–twelve- and thirteen-point rule.

source

Integration domain

FinEtools.IntegDomainModule.IntegDomainType
IntegDomain{S<:AbstractFESet, F<:Function}

Integration domain.

  • T = type of finite element set. The type of the FE set will be dependent upon the operations required. For instance, for interior (volume) integrals such as body load or the stiffness hexahedral H8 may be used, whereas for boundary (surface) integrals quadrilateral Q4 would be needed.
  • F = type of function to return the "other" dimension.

An integration domain consists of the finite elements that approximate the geometry, the function to supply the "missing" (other) dimension, indication whether or not the integration domain represents an axially symmetric situation, and integration rule used to evaluate integrals over the domain.

source
FinEtools.IntegDomainModule.IntegDomainMethod
IntegDomain(
    fes::S,
    integration_rule::IR,
    axisymmetric::Bool,
) where {S<:AbstractFESet, IR<:AbstractIntegRule}

Construct with the default orientation matrix (identity), for axially symmetric models. The other dimension is the default unity (1.0).

This will probably be called when axisymmetric = true, since the default is axisymmetric = false.

source
FinEtools.IntegDomainModule.IntegDomainMethod
IntegDomain(fes::S, integration_rule::IR) where {S<:AbstractFESet, IR<:AbstractIntegRule}

Construct with the default orientation matrix (identity), and the other dimension being the default 1.0.

source
FinEtools.IntegDomainModule.IntegDomainMethod
IntegDomain(
    fes::S,
    integration_rule::IR,
    axisymmetric::Bool,
    otherdimension::T,
) where {S<:AbstractFESet, IR<:AbstractIntegRule, T<:Number}

Construct for axially symmetric models. The other dimension is given as a number.

source
FinEtools.IntegDomainModule.IntegDomainMethod
IntegDomain(
    fes::S,
    integration_rule::IR,
    otherdimension::T,
) where {S<:AbstractFESet, IR<:AbstractIntegRule, T<:Number}

Construct with the default orientation matrix (identity), and constant other dimension.

source

Assembly of matrices and vectors

FinEtools.AssemblyModule.SysmatAssemblerFFBlockMethod
SysmatAssemblerFFBlock(row_nfreedofs::IT, col_nfreedofs = row_nfreedofs) where {IT<:Integer}

Constructor, where the wrapped assembler is for general sparse matrices.

Supply the number of free degrees of freedom.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseType
SysmatAssemblerSparse{IT, MBT, IBT} <: AbstractSysmatAssembler

Type for assembling a sparse global matrix from elementwise matrices.

Note

All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseMethod
SysmatAssemblerSparse(z = zero(T), nomatrixresult = false) where {T}

Construct a sparse system matrix assembler.

The matrix entries are of type T. The assembler either produces a sparse matrix (when nomatrixresult = true), or does not (when nomatrixresult = false). When the assembler does not produce the sparse matrix when makematrix! is called, it still can be constructed from the buffers stored in the assembler, until they are cleared when the assembler is destroyed.

Example

This is how a sparse matrix is assembled from two rectangular dense matrices.

    a = SysmatAssemblerSparse(0.0)
    startassembly!(a, 5, 5, 3, 7, 7)
    m = [0.24406   0.599773    0.833404  0.0420141
        0.786024  0.00206713  0.995379  0.780298
        0.845816  0.198459    0.355149  0.224996]
    assemble!(a, m, [1 7 5], [5 2 1 4])
    m = [0.146618  0.53471   0.614342    0.737833
         0.479719  0.41354   0.00760941  0.836455
         0.254868  0.476189  0.460794    0.00919633
         0.159064  0.261821  0.317078    0.77646
         0.643538  0.429817  0.59788     0.958909]
    assemble!(a, m, [2 3 1 7 5], [6 7 3 4])
    A = makematrix!(a)

Here A is a sparse matrix of the size 7x7.

When the nomatrixresult is set as true, no matrix is produced.

    a = SysmatAssemblerSparse(0.0, true)
    startassembly!(a, 5, 5, 3, 7, 7)
    m = [0.24406   0.599773    0.833404  0.0420141
        0.786024  0.00206713  0.995379  0.780298
        0.845816  0.198459    0.355149  0.224996]
    assemble!(a, m, [1 7 5], [5 2 1 4])
    m = [0.146618  0.53471   0.614342    0.737833
         0.479719  0.41354   0.00760941  0.836455
         0.254868  0.476189  0.460794    0.00919633
         0.159064  0.261821  0.317078    0.77646
         0.643538  0.429817  0.59788     0.958909]
    assemble!(a, m, [2 3 1 7 5], [6 7 3 4])
    A = makematrix!(a)

Here A is a named tuple of four sparse zero matrices. To construct the correct matrix is still possible, for instance like this:

    a.nomatrixresult = false
    A = makematrix!(a)

At this point all the buffers of the assembler have potentially been cleared, and makematrix!(a) is no longer possible.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseDiagType
SysmatAssemblerSparseDiag{T<:Number} <: AbstractSysmatAssembler

Assembler for a symmetric square diagonal matrix assembled from symmetric square diagonal matrices.

Warning: off-diagonal elements of the elementwise matrices will be ignored during assembly!

Note

All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseHRZLumpingSymmType
SysmatAssemblerSparseHRZLumpingSymm{IT, MBT, IBT} <: AbstractSysmatAssembler

Assembler for a symmetric lumped square matrix assembled from symmetric square matrices.

Reference: A note on mass lumping and related processes in the finite element method, E. Hinton, T. Rock, O. C. Zienkiewicz, Earthquake Engineering & Structural Dynamics, volume 4, number 3, 245–249, 1976.

Note

All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.

Note

This assembler can compute and assemble diagonalized mass matrices. However, if the meaning of the entries of the mass matrix differs (translation versus rotation), the mass matrices will not be computed correctly. Put bluntly: it can only be used for homogeneous mass matrices, all translation degrees of freedom, for instance.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseSymmType
SysmatAssemblerSparseSymm{IT, MBT, IBT} <: AbstractSysmatAssembler

Assembler for a symmetric square matrix assembled from symmetric square matrices.

Note

All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.

source
FinEtools.AssemblyModule.SysmatAssemblerSparseSymmMethod
SysmatAssemblerSparseSymm(z::T, nomatrixresult = false) where {T}

Construct blank system matrix assembler for symmetric matrices. The matrix entries are of type T.

Example

This is how a symmetric sparse matrix is assembled from two square dense matrices.

    a = SysmatAssemblerSparseSymm(0.0)
    startassembly!(a, 5, 5, 3, 7, 7)
    m = [0.24406   0.599773    0.833404  0.0420141
        0.786024  0.00206713  0.995379  0.780298
        0.845816  0.198459    0.355149  0.224996]
    assemble!(a, m'*m, [5 2 1 4], [5 2 1 4])
    m = [0.146618  0.53471   0.614342    0.737833
         0.479719  0.41354   0.00760941  0.836455
         0.254868  0.476189  0.460794    0.00919633
         0.159064  0.261821  0.317078    0.77646
         0.643538  0.429817  0.59788     0.958909]
    assemble!(a, m'*m, [2 3 1 5], [2 3 1 5])
    A = makematrix!(a)

See also

SysmatAssemblerSparse

source

Mesh import/export

FEM machines

Base

Material models

Material model abstractions