Types
Contents
Coordinate systems
FinEtools.CSysModule.CSys — TypeCSys{T<:Number, F<:Function}Type for coordinate system transformations. Used to define material coordinate systems, and output coordinate systems, for instance.
FinEtools.CSysModule.CSys — MethodCSys(dim::IT) where {IT}Construct coordinate system when the rotation matrix is the identity.
dim = is the space dimension.
FinEtools.CSysModule.CSys — MethodCSys(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
computecsmatfunction signature:update!( csmatout::AbstractMatrix{<:Real}, XYZ::AbstractVecOrMat{<:Real}, tangents::AbstractMatrix{<:Real}, feid::Integer, qpid::Integer, )wherecsmatout= 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
endFinEtools.CSysModule.CSys — MethodCSys(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::AbstractMatrix{<:Real},
XYZ::AbstractVecOrMat{<:Real},
tangents::AbstractMatrix{<:Real},
feid::Integer,
qpid::Integer,
) 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
endFinEtools.CSysModule.CSys — MethodCSys(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.
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
FinEtools.CSysModule.CSys — MethodCSys(csmat::Matrix{T}) where {T}Construct coordinate system when the rotation matrix is given.
FinEtools.CSysModule.CSys — MethodCSys(dim, z::T) where {T}Construct coordinate system when the rotation matrix of element type T is the identity.
dim = is the space dimension.
Data cache
FinEtools.DataCacheModule.DataCache — TypeDataCache{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::AbstractVecOrMat{<:Real}, tangents::AbstractMatrix{<:Real},
feid::IT, qpid::IT) where {D, T, IT}
... # modify the value of cacheout
return cacheout
endIt 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::AbstractVecOrMat{<:Real}, tangents::AbstractMatrix{<:Real},
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)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.
FinEtools.DataCacheModule.DataCache — MethodDataCache(data::D) where {D}Construct data cache. The constant data is given.
FinEtools.DataCacheModule.DataCache — Method(c::DataCache)(
XYZ::AbstractVecOrMat{<:Real},
tangents::AbstractMatrix{<:Real},
feid::IT,
qpid::IT,
) where {IT<:Integer}Update the cache and retrieve the array.
Surface-normal utilities
FinEtools.SurfaceNormalModule.SurfaceNormal — TypeSurfaceNormal{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::AbstractVecOrMat{<:Real},
tangents::AbstractMatrix{<:Real}, feid::IT, qpid::IT) where {CT, IT}
# Calculate the normal and copy it into the buffer....
return normalout
endThe buffer normalout needs to be filled with the value of the normal vector.
Refer to DataCache for details on the caching.
FinEtools.SurfaceNormalModule.SurfaceNormal — MethodSurfaceNormal(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.
FinEtools.SurfaceNormalModule.SurfaceNormal — MethodSurfaceNormal(ndimensions, computenormal!::F) where {F<:Function}Construct surface normal evaluator when the function to compute the normal vector is given.
FinEtools.SurfaceNormalModule.SurfaceNormal — MethodSurfaceNormal(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 functioncomputenormal!needs to have a signature offunction computenormal!(normalout::Vector{CT}, XYZ::AbstractVecOrMat{<:Real}, tangents::AbstractMatrix{<:Real}, feid::IT, qpid::IT) where {CT, IT} # Calculate the normal and copy it into the buffer.... return normalout endand it needs to fill in the buffer
normaloutwith the current normal at the locationXYZ, using, if appropriate, the information supplied in the Jacobian matrixtangents, the identifier of the finite element,feid, and the quadrature point id,qpid. Refer toDataCache.
FinEtools.SurfaceNormalModule.SurfaceNormal — MethodSurfaceNormal(vector::Vector{T}) where {T<:Number}Construct surface normal vector when the constant normal vector is given.
Force intensity
FinEtools.ForceIntensityModule.ForceIntensity — TypeForceIntensity{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.
FinEtools.ForceIntensityModule.ForceIntensity — MethodForceIntensity(force::T) where {T<:Number}Construct force intensity when the force is given as a scalar value.
The dimension of the force vector in this case is 1.
FinEtools.ForceIntensityModule.ForceIntensity — MethodForceIntensity(
::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 functioncomputeforce!needs to have a signature offunction 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 endand it needs to fill in the bufferforceoutwith the current force at the locationXYZ, using, if appropriate, the information supplied in the Jacobian matrixtangents, and the identifier of the finite element,feid.
FinEtools.ForceIntensityModule.ForceIntensity — MethodForceIntensity(force::Vector{T}) where {T<:Number}Construct force intensity when the constant force vector is given.
Finite element sets
FinEtools.FESetModule.AbstractFESet — TypeAbstractFESet{NODESPERELEM}Abstract type of a finite element set. Parameterized with the number of of the nodes per element.
FinEtools.FESetModule.AbstractFESet0Manifold — TypeAbstractFESet0Manifold{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.
FinEtools.FESetModule.AbstractFESet1Manifold — TypeAbstractFESet1Manifold{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.
FinEtools.FESetModule.AbstractFESet2Manifold — TypeAbstractFESet2Manifold{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.
FinEtools.FESetModule.AbstractFESet3Manifold — TypeAbstractFESet3Manifold{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.
FinEtools.FESetModule.FESetH20 — TypeFESetH20Type for sets of volume-like hexahedral finite elements with 20 nodes.
FinEtools.FESetModule.FESetH27 — TypeFESetH27Type for sets of volume-like hexahedral finite elements with 27 nodes.
FinEtools.FESetModule.FESetH8 — TypeFESetH8Type for sets of volume-like hexahedral finite elements with eight nodes.
FinEtools.FESetModule.FESetL2 — TypeFESetL2Type for sets of curve-like finite elements with two nodes.
FinEtools.FESetModule.FESetL3 — TypeFESetL3Type for sets of curve-like of finite elements with three nodes.
FinEtools.FESetModule.FESetP1 — TypeFESetP1Type for sets of point-like of finite elements.
FinEtools.FESetModule.FESetQ4 — TypeFESetQ4Type for sets of surface-like quadrilateral finite elements with four nodes.
FinEtools.FESetModule.FESetQ8 — TypeFESetQ8Type for sets of surface-like quadrilateral finite elements with eight nodes.
FinEtools.FESetModule.FESetQ9 — TypeFESetQ9Type for sets of surface-like quadrilateral finite elements with nine nodes.
FinEtools.FESetModule.FESetT10 — TypeFESetT10Type for sets of volume-like tetrahedral finite elements with 10 nodes.
FinEtools.FESetModule.FESetT3 — TypeFESetT3Type for sets of surface-like triangular finite elements with three nodes.
FinEtools.FESetModule.FESetT4 — TypeFESetT4Type for sets of volume-like tetrahedral finite elements with four nodes.
FinEtools.FESetModule.FESetT6 — TypeFESetT6Type for sets of surface-like triangular finite elements with six nodes.
Finite element nodes
FinEtools.FENodeSetModule.FENodeSet — Typemutable 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).
The constructor makes a copy of the input xyz array for safety.
Finite element node-to-element map
FinEtools.FENodeToFEMapModule.FENodeToFEMap — TypeFENodeToFEMapMap 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).
FinEtools.FENodeToFEMapModule.FENodeToFEMap — MethodFENodeToFEMap(fes::FE, nmax::IT) where {FE<:AbstractFESet,IT<:Integer}Map from finite element nodes to the finite elements connecting them.
Convenience constructor.
FinEtools.FENodeToFEMapModule.FENodeToFEMap — MethodFENodeToFEMap(fes::FE, nmax::IT) where {FE<:AbstractFESet,IT<:Integer}Map from finite element nodes to the finite elements connecting them.
Convenience constructor.
FinEtools.FENodeToFEMapModule.FENodeToFEMap — MethodFENodeToFEMap(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 tuplesnmax= largest possible node number
Example:
m = FENodeToFEMap(fes.conn, count(fens))Fields
FinEtools.FieldModule.AbstractField — TypeAbstractFieldAbstract field.
Expected attributes:
values::Array{T,2}: Array of degree of freedom parameters, indexed by entity numberdofnums::Array{IT,2}: Array of degree of freedom numbers, indexed by entity numberkind::Matrix{Int8}: Array of Boolean flags, indexed by entity numberranges::Dict(Int8, UnitRange{IT}): Dictionary of ranges for the degrees of freedom.
See also: @add_Field_fields() .
FinEtools.FieldModule.KIND_INT — TypeKIND_INTConstant representing the type of the integer representing the kind of a degree of freedom.
FinEtools.GeneralFieldModule.GeneralField — TypeGeneralField{T<:Number, IT<:Integer} <: AbstractFieldGeneral field, meaning the entities can be anything.
FinEtools.GeneralFieldModule.GeneralField — MethodGeneralField(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.
FinEtools.GeneralFieldModule.GeneralField — MethodGeneralField(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.
FinEtools.NodalFieldModule.NodalField — TypeNodalField{T<:Number, IT<:Integer} <: AbstractFieldNodal field, meaning the entities are the finite element nodes.
FinEtools.NodalFieldModule.NodalField — MethodNodalField(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.
FinEtools.NodalFieldModule.NodalField — MethodNodalField(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.
FinEtools.ElementalFieldModule.ElementalField — TypeElementalField{T<:Number, IT<:Integer} <: AbstractFieldElemental 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.
FinEtools.ElementalFieldModule.ElementalField — MethodElementalField(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.
FinEtools.ElementalFieldModule.ElementalField — MethodElementalField(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.
Integration rule
FinEtools.IntegRuleModule.AbstractIntegRule — TypeAbstractIntegRuleAbstract type for integration rule.
FinEtools.IntegRuleModule.GaussRule — TypeGaussRule(dim=1, order=1)Gauss rule.
FinEtools.IntegRuleModule.GaussRule — TypeGaussRule <: AbstractIntegRuleThe Gauss rule, applicable for a tensor product of intervals -1 <=x<= +1.
FinEtools.IntegRuleModule.NodalSimplexRule — TypeNodalSimplexRule(dim=1)Nodal-quadrature simplex rule.
FinEtools.IntegRuleModule.NodalSimplexRule — TypeNodalSimplexRule <: AbstractIntegRuleThe nodal-quadrature simplex rule.
The rule is applicable for line segments, triangles, tetrahedra.
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!
FinEtools.IntegRuleModule.NodalTensorProductRule — TypeNodalTensorProductRule(dim=1)Nodal-quadrature tensor-product rule.
FinEtools.IntegRuleModule.NodalTensorProductRule — TypeNodalTensorProductRule <: AbstractIntegRuleThe tensor-product nodal-quadrature rule.
The rule is applicable for line segments, quadrilaterals, hexahedra.
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!
FinEtools.IntegRuleModule.PointRule — TypePointRule <: AbstractIntegRulePoint quadrature rule, used for integration on the standard "point" shape.
FinEtools.IntegRuleModule.PointRule — MethodPointRule()POINT integration rule.
FinEtools.IntegRuleModule.SimplexRule — TypeSimplexRule(dim=1, npts=1)Return simplex rule, appropriate for the manifold dimension dim.
FinEtools.IntegRuleModule.SimplexRule — TypeSimplexRule <: AbstractIntegRuleSimplex quadrature rule.
Used for integration on the standard triangle or the standard tetrahedron.
FinEtools.IntegRuleModule.TetRule — TypeTetRule(npts=1)Tetrahedral integration rule. npts=number of points (1– one-point rule, 4 – four-point rule, 5 – five point rule).
FinEtools.IntegRuleModule.TetRule — TypeTetRule <: AbstractIntegRuleTetrahedral quadrature rule, used for integration on the standard tetrahedron.
FinEtools.IntegRuleModule.TrapezoidalRule — TypeTrapezoidalRule(dim=1)Trapezoidal rule.
FinEtools.IntegRuleModule.TrapezoidalRule — TypeTrapezoidalRule <: AbstractIntegRuleThe trapezoidal rule.
The rule is applicable for a tensor product of intervals -1 <=x<= +1.
FinEtools.IntegRuleModule.TriRule — TypeTriRule(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.
FinEtools.IntegRuleModule.TriRule — TypeTriRule <: AbstractIntegRuleTriangular quadrature rule for integration on the standard triangle.
Integration domain
FinEtools.IntegDomainModule.IntegDomain — TypeIntegDomain{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.
FinEtools.IntegDomainModule.IntegDomain — MethodIntegDomain(
fes::S,
integration_rule::IR,
otherdimension::F,
) where {S<:AbstractFESet, IR<:AbstractIntegRule, F<:Function}Construct with the default orientation matrix (identity), and other dimension provided by a function.
FinEtools.IntegDomainModule.IntegDomain — MethodIntegDomain(
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.
FinEtools.IntegDomainModule.IntegDomain — MethodIntegDomain(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.
FinEtools.IntegDomainModule.IntegDomain — MethodIntegDomain(
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.
FinEtools.IntegDomainModule.IntegDomain — MethodIntegDomain(
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.
Assembly of matrices and vectors
FinEtools.AssemblyModule.AbstractSysmatAssembler — TypeAbstractSysmatAssemblerAbstract type of system-matrix assembler.
FinEtools.AssemblyModule.AbstractSysvecAssembler — TypeAbstractSysvecAssemblerAbstract type of system vector assembler.
FinEtools.AssemblyModule.SysmatAssemblerFFBlock — TypeSysmatAssemblerFFBlock{A<:AbstractSysmatAssembler, IT} <: AbstractSysmatAssemblerType for extracting a free-free matrix, delegating the actual assembly to a different assembler.
FinEtools.AssemblyModule.SysmatAssemblerFFBlock — MethodSysmatAssemblerFFBlock(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.
FinEtools.AssemblyModule.SysmatAssemblerSparse — TypeSysmatAssemblerSparse{IT, MBT, IBT} <: AbstractSysmatAssemblerType for assembling a sparse global matrix from elementwise matrices.
All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.
FinEtools.AssemblyModule.SysmatAssemblerSparse — MethodSysmatAssemblerSparse(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.
FinEtools.AssemblyModule.SysmatAssemblerSparseDiag — TypeSysmatAssemblerSparseDiag{T<:Number} <: AbstractSysmatAssemblerAssembler 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!
All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.
FinEtools.AssemblyModule.SysmatAssemblerSparseDiag — MethodSysmatAssemblerSparseDiag(z::T, nomatrixresult = false) where {T}Construct blank system matrix assembler for square diagonal matrices. The matrix entries are of type T.
FinEtools.AssemblyModule.SysmatAssemblerSparseHRZLumpingSymm — TypeSysmatAssemblerSparseHRZLumpingSymm{IT, MBT, IBT} <: AbstractSysmatAssemblerAssembler 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.
All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.
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.
FinEtools.AssemblyModule.SysmatAssemblerSparseHRZLumpingSymm — MethodSysmatAssemblerSparseHRZLumpingSymm(z::T, nomatrixresult = false) where {T}Construct blank system matrix assembler. The matrix entries are of type T.
FinEtools.AssemblyModule.SysmatAssemblerSparseSymm — TypeSysmatAssemblerSparseSymm{IT, MBT, IBT} <: AbstractSysmatAssemblerAssembler for a symmetric square matrix assembled from symmetric square matrices.
All fields of the datatype are private. The type is manipulated by the functions startassembly!, assemble!, and makematrix!.
FinEtools.AssemblyModule.SysmatAssemblerSparseSymm — MethodSysmatAssemblerSparseSymm(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
FinEtools.AssemblyModule.SysvecAssembler — TypeSysvecAssemblerAssembler for the system vector.
FinEtools.AssemblyModule.SysvecAssembler — MethodSysvecAssembler(z::T) where {T}Construct blank system vector assembler. The vector entries are of type T determined by the zero value.
FinEtools.AssemblyModule.SysvecAssemblerFBlock — TypeSysvecAssemblerFBlockAssembler for the system vector, which extracts the free vector.
FinEtools.AssemblyModule.SysvecAssemblerFBlock — MethodSysvecAssemblerFBlock(row_nfreedofs::IT) where {IT}Constructor of the free block assembler.
Mesh import/export
FEM machines
Base
FinEtools.FEMMBaseModule.AbstractFEMM — TypeAbstractFEMMAbstract type for all finite element model machines.
FinEtools.FEMMBaseModule.FEMMBase — TypeFEMMBase{ID<:IntegDomain, CS<:CSys} <: AbstractFEMMType for base finite element modeling machine.
FinEtools.FEMMBaseModule.FEMMBase — MethodFEMMBase(integdomain::ID) where {ID<:IntegDomain}Construct with the default orientation matrix (identity).
Material models
Material model abstractions
FinEtools.MatModule.AbstractMat — TypeAbstractMatAbstract type of material.