Reference manual
Reference shapes
Elfel.RefShapes.AbstractRefShape
— TypeAbstractRefShape{MANIFDIM}
Abstract type of a reference shape.
Elfel.RefShapes.IntegRule
— TypeIntegRule
Type for integration rule.
Elfel.RefShapes.RefShapeCube
— TypeRefShapeCube <: AbstractRefShape{3}
Type of a reference shape for a 3-dimensional manifold (solid) bounded by six quadrilaterals.
Elfel.RefShapes.RefShapeInterval
— TypeRefShapeInterval <: AbstractRefShape{1}
Type of a reference shape for a 1-dimensional manifold (curve).
Elfel.RefShapes.RefShapePoint
— TypeRefShapePoint <: AbstractRefShape{0}
Type of a reference shape for a zero-dimensional manifold (point).
Elfel.RefShapes.RefShapeSquare
— TypeRefShapeSquare <: AbstractRefShape{2}
Type of a logically rectangular reference shape for a 2-dimensional manifold (surface).
Elfel.RefShapes.RefShapeTetrahedron
— TypeRefShapeTetrahedron <: AbstractRefShape{3}
Type of a reference shape for a 3-dimensional manifold (solid) bounded by 4 triangles.
Elfel.RefShapes.RefShapeTriangle
— TypeRefShapeTriangle <: AbstractRefShape{2}
Type of a logically triangular reference shape for a 2-dimensional manifold (surface).
Elfel.RefShapes.manifdim
— Methodmanifdim(rs)
Get the manifold dimension of the reference shape.
Elfel.RefShapes.manifdimv
— Methodmanifdimv(::Type{T}) where {T<:AbstractRefShape{MANIFDIM}} where {MANIFDIM}
Get the manifold dimension of the reference shape as Val
.
Elfel.RefShapes.quadrature
— Functionquadrature(::Type{RefShapeInterval}, quadraturesettings = (kind = :default,))
Create a quadrature rule for the reference shape of an interval.
The default is Gauss integration rule, where the order is set with the keyword order
.
Elfel.RefShapes.quadrature
— Functionquadrature(::Type{RefShapeTriangle}, quadraturesettings = (kind = :default,))
Create a quadrature rule for the reference shape of an triangle.
The default is a triangle rule, distinguished by the number of points set with the keyword npts
.
Elfel.RefShapes.quadrature
— Functionquadrature(::Type{RefShapeTetrahedron}, quadraturesettings = (kind = :default,))
Create a quadrature rule for the reference shape of a tetrahedron.
The default is a one-point tetrahedron rule; other rules may be chosen based on the number of points set with the keyword npts
.
Elfel.RefShapes.quadrature
— Functionquadrature(::Type{RefShapeSquare}, quadraturesettings = (kind = :default,))
Create a quadrature rule for the reference shape of a square.
The default is Gauss integration rule, where the order is set with the keyword order
.
Elements
Elfel.FElements.FE
— TypeFE{RS, SD}
Abstract type of finite element, parameterized by
RS
: type of reference shape of the element (triangle, square, ...), andSD
: shape descriptor; refer to the packageMeshCore
.
Elfel.FElements.FEData
— TypeFEData{SD}
Type of a finite element data.
Parameterized by
SD
= shape descriptor; refer to the packageMeshCore
.
Elfel.FElements.FEH1_L2
— MethodFEH1_L2()
Construct an H1 finite element of the type L2.
L2 is two-node linear segment element.
Elfel.FElements.FEH1_Q4
— MethodFEH1_Q4()
Construct an H1 finite element of the type Q4.
Q4 is 4-node linear quadrilateral element.
Elfel.FElements.FEH1_T3
— MethodFEH1_T3()
Construct an H1 finite element of the type T3.
T3 is 3-node linear triangle element.
Elfel.FElements.FEH1_T3_BUBBLE
— MethodFEH1_T3_BUBBLE()
Construct an H1 finite element of the type T3 with a cubic bubble.
T3 is 3-node linear triangle element with a cubic bubble. It has the usual nodal basis functions associated with the vertices, and cubic bubble associated with the element itself.
Elfel.FElements.FEH1_T4
— MethodFEH1_T4()
Construct an H1 finite element of the type T4.
T4 is 4-node linear tetrahedral element.
Elfel.FElements.FEH1_T6
— MethodFEH1_T6()
Construct an H1 finite element of the type T6.
T6 is 6-node quadratic triangle element.
Elfel.FElements.FEL2_Q4
— MethodFEL2_Q4()
Construct an L2 finite element of the type Q4.
Q4 is 4-node linear quadrilateral element.
Elfel.FElements.FEL2_T3
— MethodFEL2_T3()
Construct an L2 finite element of the type T3.
T3 is 3-node linear triangle element.
Elfel.FElements.FEL2_T4
— MethodFEL2_T4()
Construct an L2 finite element of the type T4.
T4 is tetrahedral element with only internal degrees of freedom.
Elfel.FElements.Jacobian
— MethodJacobian(::Val{0}, J::T) where {T}
Evaluate the point Jacobian.
J
= Jacobian matrix, which isn't really defined well for a 0-manifold.
Elfel.FElements.Jacobian
— MethodJacobian(::Val{1}, J::T) where {T}
Evaluate the curve Jacobian.
J
= Jacobian matrix, columns are tangent to parametric coordinates curves.
Elfel.FElements.Jacobian
— MethodJacobian(::Val{2}, J::T) where {T}
Evaluate the curve Jacobian.
J
= Jacobian matrix, columns are tangent to parametric coordinates curves.
Elfel.FElements.Jacobian
— MethodJacobian(fe::T, J::FFltMat)::FFlt where {T<:AbstractFESet3Manifold}
Evaluate the volume Jacobian.
J
= Jacobian matrix, columns are tangent to parametric coordinates curves.
Elfel.FElements._geometrycarrier
— MethodThe elements that have nodal bases can be their own geometry carriers. Elements that do not have nodal bases, for instance an L2 quadrilateral with a single basis function associated with the cell, need another element type (in this case the nodal quadrilateral) to serve as geometry carriers.
Elfel.FElements.bfun
— Methodbfun(fe::FESUBT, param_coords) where {FESUBT<:FE{RS, SD}}
Evaluate the basis functions for all degrees of freedom of the scalar finite element at the parametric coordinates. Return a vector of the values.
Elfel.FElements.bfungradpar
— Methodbfungradpar(fe::FESUBT, param_coords) where {FESUBT<:FE{RS, SD}}
Evaluate the gradients of the basis functions for all degrees of freedom of the scalar finite element with respect to the parametric coordinates, at the parametric coordinates given. Return a vector of the gradients.
Elfel.FElements.jacjac
— Methodjacjac(fe::FE{RS, SD}, locs, nodes, gradNpar) where {RS, SD}
Compute the Jacobian matrix and the Jacobian determinant.
This is the generic version suitable for isoparametric elements.
Elfel.FElements.ndofperfeat
— Methodfeathasdof(fe::FE{RS, SD}, m) where {RS, SD}
How many degrees of freedom are attached to a the feature of manifold dimension m
?
Note that 0 <= m <= 3
.
Elfel.FElements.ndofsperel
— Methodndofsperel(fe::FE{RS, SD}) where {RS, SD}
Provide the number of degrees of freedom per element.
Enumerate all features of all manifold dimensions, and for each feature multiply by the number of degrees of freedom per feature. The assumption is that this is a scalar finite element.
Elfel.FElements.nfeatofdim
— Methodnfeatofdim(fe::FE{RS, SD}, m) where {RS, SD}
Number of features of manifold dimension m
. Note that 0 <= m <= 3
.
Elfel.FElements.refshape
— Methodrefshape(fe::FE{RS, SD}) where {RS, SD}
Reference shape.
Elfel.FElements.shapedesc
— Methodshapedesc(fe::FE{RS, SD}) where {RS, SD}
Topological shape description.
Refer to the MeshCore library.
MeshCore.manifdim
— Methodmanifdim(fe::FE{RS, SD}) where {RS, SD}
Get the manifold dimension of the finite element.
Spaces
Elfel.FESpaces.FESpace
— TypeFESpace{FET, T}
Type of a finite element space, parameterized with
FET
: type of finite element, it is a scalar finite element,T
: type of degree of freedom value (double, complex, ...).
Elfel.FESpaces.FEFields.gathersysvec!
— Methodgathersysvec!(v, fesp)
Gather values for the whole system vector from all FE spaces contributing to it.
fesp
is either a vector or a tuple of FE spaces.
Elfel.FESpaces.FEFields.gathersysvec!
— Methodgathersysvec!(v, fesp::FESpace)
Gather values for the whole system vector.
Elfel.FESpaces.FEFields.highestdatadofnum
— Methodhighestfreedofnum(fesp::FES) where {FES<:FESpace}
Compute the highest number of data (known) degrees of freedom.
Elfel.FESpaces.FEFields.highestfreedofnum
— Methodhighestfreedofnum(fesp::FES) where {FES<:FESpace}
Compute the highest number of free (unknown) degrees of freedom.
Elfel.FESpaces.FEFields.ndofs
— Methodndofs(fesp::FES) where {FES<:FESpace}
Compute the total number of degrees of freedom.
Elfel.FESpaces.FEFields.numberdatadofs!
— Methodnumberdatadofs!(fesp::FES, firstnum = 1) where {FES<:FESpace}
Number the data (known) degrees of freedom.
The known degrees of freedom in the FE space are numbered consecutively.
No effort is made to optimize the numbering in any way.
Elfel.FESpaces.FEFields.numberfreedofs!
— Methodnumberfreedofs!(fesp::FES, firstnum = 1) where {FES<:FESpace}
Number the free degrees of freedom.
The unknown degrees of freedom in the FE space are numbered consecutively.
No effort is made to optimize the numbering in any way.
Elfel.FESpaces.FEFields.scattersysvec!
— Methodscattersysvec!(fesp, v)
Scatter values for the whole system vector to all FE spaces contributing to it.
fesp
is either a vector or a tuple of FE spaces.
Elfel.FESpaces.FEFields.scattersysvec!
— Methodscattersysvec!(fesp::FESpace, v)
Scatter values from the system vector.
Elfel.FESpaces.FEFields.setebc!
— Methodsetebc!(fesp::FESpace, m, eid, comp, val::T) where {T}
Set the EBCs (essential boundary conditions).
m
= manifold dimension of the entity,eid
= serial number of the entity (term identifier),comp
= which degree of freedom in the term,val
= value of type T
For instance, m = 0
means set the degree of freedom at the vertex eid
.
Elfel.FESpaces.dofnum
— Methoddofnum(fesp::FES, m, eid) where {FES<:FESpace}
Provide degree of freedom number for entity eid
of manifold dimension m
and component comp
.
Elfel.FESpaces.doftype
— Methoddoftype(fesp::FESpace{FET, T}) where {FET, T}
Provide the type of the values of the degrees of freedom.
Elfel.FESpaces.edofbfnum
— Methodedofbfnum(fesp::FESpace{FET, T}) where {FET, T}
Access vector of numbers of basis functions associated with each degree of freedom.
Elfel.FESpaces.edofcompnt
— Methodedofcompnt(fesp::FESpace{FET, T}) where {FET, T}
Access vector of component number associated with each degree of freedom.
When the finite element space consists of multiple copies of the scalar finite element, the component is the serial number of the copy.
Elfel.FESpaces.edofmdim
— Methodedofmdim(fesp::FESpace{FET, T}) where {FET, T}
Access vector of manifold dimensions of entities associated with each degree of freedom.
Elfel.FESpaces.makeattribute
— Methodmakeattribute(fesp::FESpace, name, comp)
Attach attribute to the right shape collection of all incidence relations.
Elfel.FESpaces.numberdofs!
— Methodnumberdofs!(fesp::AbstractVector)
Number the degrees of freedom of a collection of FE spaces.
The unknown (free) degrees of freedom in the FE space are numbered consecutively, and then the data degrees of freedom (the known values) are numbered.
No effort is made to optimize the numbering in any way.
Elfel.FESpaces.numberdofs!
— Methodnumberdofs!(fesp::FESpace)
Number the degrees of freedom of a single FE space.
The unknown (free) degrees of freedom in the FE space are numbered consecutively, and then the data degrees of freedom (the known values) are numbered.
No effort is made to optimize the numbering in any way.
Elfel.FESpaces.nunknowns
— Methodnunknowns(fesp::FES) where {FES<:FESpace}
Compute the total number of unknown degrees of freedom.
Elfel.FESpaces.setdofval!
— Methodsetdofval!(fesp::FESpace, m, eid, comp, val::T) where {T}
Set the value of a degree of freedom.
m
= manifold dimension of the entity,eid
= serial number of the entity (term identifier),comp
= which degree of freedom in the term,val
= value of type T
For instance, m = 0
means set the degree of freedom at the vertex eid
.
Elfel.FElements.ndofsperel
— Methodndofsperel(fesp::FES) where {FES<:FESpace}
Total number of degrees of freedom associated with each finite element.
Essentially a product of the number of the degrees of freedom the scalar finite element and the number of copies of this element in the space.
Fields
Elfel.FESpaces.FEFields.FEField
— TypeFEField{N, T, IT}
Type of a finite element field. Parameterized with
N
: number of degrees of freedom per entity,T
: type of the degree of freedom value,IT
: type of the index (integer value). This describes the serial numbers of the degrees of freedom.
Elfel.FESpaces.FEFields.datadofnums
— Methoddatadofnums(f::FEField)
Collect information about known (data) degree of freedom numbers.
First number, last number, and the total number of degrees of freedom are returned as a tuple.
Elfel.FESpaces.FEFields.dofnums
— Methoddofnums(fef::FEField, n)
Provide degree of freedom numbers for given entity.
Elfel.FESpaces.FEFields.dofnumtype
— Methoddofnumtype(fef::FEField{N, T, IT}) where {N, T, IT}
Type of the index (serial number) of the degree of freedom. Integer.
Elfel.FESpaces.FEFields.doftype
— Methoddoftype(fef::FEField{N, T, IT}) where {N, T, IT}
Type of a degree of freedom value.
Elfel.FESpaces.FEFields.dofvals
— Methoddofvals(fef::FEField, n)
Provide degree of freedom values for given entity.
Elfel.FESpaces.FEFields.freedofnums
— Methodfreedofnums(f::FEField)
Collect information about unknown (free) degree of freedom numbers.
First number, last number, and the total number of degrees of freedom are returned as a tuple.
Elfel.FESpaces.FEFields.gathersysvec!
— Methodgathersysvec!(vec, self::FEField)
Gather system vector contributions from the field.
Elfel.FESpaces.FEFields.highestdatadofnum
— Methodhighestdatadofnum(f::FEField)
Compute the highest serial number of a datum degree of freedom in the field.
Elfel.FESpaces.FEFields.highestfreedofnum
— Methodhighestfreedofnum(f::FEField)
Compute the highest serial number of a free degree of freedom in the field.
Elfel.FESpaces.FEFields.ndofs
— Methodndofs(fef::FEField)
Total number of degrees of freedom in the field.
Elfel.FESpaces.FEFields.ndofsperterm
— Methodndofsperterm(fef::FEField{N}) where {N}
Number of degrees of freedom per term of the field.
Elfel.FESpaces.FEFields.nterms
— Methodnterms(fef::FEField)
Number of terms in the field.
Elfel.FESpaces.FEFields.numberdatadofs!
— Functionnumberdatadofs!(f::FEField, firstnum = 1)
Number the data degrees of freedom in the field. Start from the number supplied on input.
Note: The free degrees of freedom must be numbered first.
Elfel.FESpaces.FEFields.numberfreedofs!
— Functionnumberfreedofs!(f::FEField, firstnum = 1)
Number the unknowns in the field, starting from the one supplied on input.
Note: The data degrees of freedom have their numbers zeroed out.
Elfel.FESpaces.FEFields.scattersysvec!
— Methodscattersysvec!(self::FEField, v)
Scatter a system vector into the field.
Elfel.FESpaces.FEFields.setdofval!
— Methodsetdofval!(self::FEField, tid, comp, val::T) where {T}
Set the value of one particular degree of freedom to a given number.
tid
: which term,comp
: which component of the term,val
: value to which the degree of freedom should be set.
Elfel.FESpaces.FEFields.setebc!
— Methodsetebc!(self::FEField, tid, comp, val::T) where {T}
Set the value of one particular degree of freedom to a given number.
tid
: which term,comp
: which component of the term,val
: value to which the degree of freedom should be set.
Elfel.FESpaces.FEFields.setisdatum!
— Methodsetebc!(self::FEField, tid, comp, val::T) where {T}
Set the value of one particular degree of freedom to a given number.
tid
: which term,comp
: which component of the term,flag
: true or false.
Finite element iterators
Elfel.FEIterators.FEIterator
— TypeFEIterator{FES, IR, G, IT, T, V, IR0, IR1, IR2, IR3, F0, F1, F2, F3}
Type of finite element iterator. Parameterized with the types of
FES
: finite element space,IR
: base incidence relation of the mesh,G
: type of the geometry attribute,IT
: type of integer indices, such as the numbers of nodes and degrees of freedom,T
: type of the degree of freedom value (real double, complex float, ... ),V
:Val
representation of the manifold dimension of the base relation elements,IR0
,IR1
,IR2
,IR3
: types of incidence relations with which degrees of freedom are associated in the finite element space, for each of the manifolds dimensions 0, 1, 2, 3,F0
,F1
,F2
,F3
: types of fields with which degrees of freedom are associated in the finite element space, for each of the manifolds dimensions 0, 1, 2, 3.
Base.iterate
— FunctionBase.iterate(it::FEIterator, state = 1)
Advance the iterator to the next entity.
The nodes of the finite element are cached, as is a vector of all the degrees of freedom represented on the element.
Base.length
— MethodBase.length(it::FEIterator)
Number of elements represented by this iterator.
Elfel.FEIterators.eldofcomps
— Methodeldofcomps(it::FEIterator)
Retrieve the vector of the component numbers for each element degree of freedom.
If multiple copies of the finite element are referenced in the finite element space, each copy is referred to as component.
Elfel.FEIterators.eldofentmdims
— Methodeldofentmdims(it::FEIterator)
Retrieve the vector of the entity dimensions for each element degree of freedom.
Each degree of freedom is associated with some entity of the finite element: vertices, edges, faces, and so on. This vector records the dimension of the manifold entity with which each degree of freedom is associated.
Elfel.FEIterators.eldofs
— Methodeldofs(it::FEIterator)
Retrieve the vector of the element degrees of freedom
Elfel.FEIterators.eldofvals
— Methodeldofvals(it::FEIterator)
Provide access to vector of element degrees of freedom.
Elfel.FEIterators.elnodes
— Methodelnodes(it::FEIterator)
Retrieve the vector of the nodes of the element.
Elfel.FEIterators.location
— Methodlocation(it::FEIterator, qpit::QPIterator)
Calculate the location of the quadrature point.
Elfel.FElements.jacjac
— Methodjacjac(it::FEIterator, qpit::QPIterator)
Compute the Jacobian matrix and the Jacobian determinant.
The finite element iterator cooperates with the quadrature point iterator here to compute the Jacobian at the current integration point.
Elfel.FElements.ndofsperel
— Methodndofsperel(it::FEIterator)
Retrieve the number of degrees of freedom per element.
Quadrature-point iterators
Elfel.QPIterators.QPIterator
— TypeQPIterator{FES, MDIM}
Type of quadrature-point iterator, parameterized by
FES
: the type of the finite element space,MDIM
: the manifold dimension of the finite element.
Elfel.QPIterators.QPIterator
— MethodQPIterator(fesp::FES, quadraturesettings) where {FES}
Construct quadrature-point iterator by associating it with a finite element space and supplying quadrature rule settings.
Base.iterate
— FunctionBase.iterate(it::QPIterator, state = 1)
Advance a quadrature point iterator.
Elfel.FElements.bfun
— Methodbfun(it::QPIterator)
Retrieve vector of basis function values for the current quadrature point.
Elfel.FElements.bfungradpar
— Methodbfungradpar(it::QPIterator)
Retrieve vector of basis function gradients with respect to the parametric coordinates for the current quadrature point.
Elfel.QPIterators.bfungrad
— Methodbfungrad(it::QPIterator, Jac)
Retrieve vector of basis function gradients with respect to spatial coordinates for the current quadrature point.
The Jacobian matrix maps between vectors in the parametric space and the spatial vectors.
Elfel.QPIterators.weight
— Methodweight(it::QPIterator)
Retrieve weight of the current quadrature point.
Assemblers
Elfel.Assemblers.AbstractSysmatAssembler
— TypeAbstractSysmatAssembler
Abstract type of system-matrix assembler.
Elfel.Assemblers.AbstractSysvecAssembler
— TypeAbstractSysvecAssembler
Abstract type of system vector assembler.
Elfel.Assemblers.SysmatAssemblerSparse
— TypeSysmatAssemblerSparse{T<:Number} <: AbstractSysmatAssembler
Type for assembling a sparse global matrix from individual entries.
Elfel.Assemblers.SysmatAssemblerSparse
— MethodSysmatAssemblerSparse(zero::T=0.0) where {T<:Number}
Construct blank system matrix assembler. The matrix entries are of type T
.
Example
This is how a sparse matrix is assembled from two rectangular dense matrices.
a = SysmatAssemblerSparse(0.0)
start!(a, 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]
gi = [1 7 5]
gj = [5 2 1 4]
for j in 1:size(m, 2), i in 1:size(m, 1)
assemble!(a, gi[i], gj[j], m[i, j])
end
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]
gi = [2 3 1 7 5]
gj = [6 7 3 4]
for j in 1:size(m, 2), i in 1:size(m, 1)
assemble!(a, gi[i], gj[j], m[i, j])
end
A = finish!(a)
Elfel.Assemblers.SysvecAssembler
— TypeSysvecAssembler
Assembler for the system vector.
Elfel.Assemblers.SysvecAssembler
— MethodSysvecAssembler(zero::T=0.0) where {T<:Number}
Construct blank system vector assembler. The vector entries are of type T
.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SysmatAssemblerSparse{T}, r, c, v::T) where {T<:Number}
Assemble a single entry of a rectangular matrix.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SysvecAssembler{T}, i, val::T) where {T<:Number}
Assemble a single value into the row i
.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SysvecAssembler{T}, lva) where {T<:Number}
Assemble local vector assembler.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SysmatAssemblerSparse{T}, lma::LocalMatrixAssembler{IT, T}) where {IT<:Integer, T<:Number}
Assemble the row numbers, column numbers, and values from a local assembler.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SysmatAssemblerSparse{T}, lma::Transpose{T,LocalMatrixAssembler{IT,T}}) where {IT<:Integer, T<:Number}
Assemble the row numbers, column numbers, and values from a local assembler.
Elfel.Assemblers.assemble!
— Methodassemble!(self::SV, i, val::T) where {SV<:AbstractSysvecAssembler, T<:Number}
Assemble an elementwise vector.
The method assembles a column element vector using the vector of degree of freedom numbers for the rows.
Elfel.Assemblers.finish!
— Methodfinish!(self::SysmatAssemblerSparse)
Make a sparse matrix.
Elfel.Assemblers.finish!
— Methodfinish!(self::SysvecAssembler)
Make the global vector.
Elfel.Assemblers.finish!
— Methodmakevector!(self::SysvecAssembler)
Make the global vector.
Elfel.Assemblers.start!
— Methodstart!(self::SysmatAssemblerSparse{T}, nrow, ncol) where {T<:Number}
Start the assembly of a global matrix.
Elfel.Assemblers.start!
— Methodstart!(self::SysvecAssembler{T}, nrow::Int64) where {T<:Number}
Start assembly.
The method makes the buffer for the vector assembly. It must be called before the first call to the method assemble.
nrow
= Total number of degrees of freedom.
Elfel.Assemblers.start!
— Methodstart!(self::SV, nrow) where {SV<:AbstractSysvecAssembler, T<:Number}
Start assembly.
The method makes the buffer for the vector assembly. It must be called before the first call to the method assemble.
nrow
= Total number of degrees of freedom.
Local Assemblers
Base.IndexStyle
— MethodBase.IndexStyle(::Type{<:LocalMatrixAssembler})
The data storage is assumed to be consumed best by one-dimensional traversal through the columns in a linear fashion.
Base.IndexStyle
— MethodBase.IndexStyle(::Type{<:LocalVectorAssembler})
Only linear access is provided.
Elfel.LocalAssemblers.LocalMatrixAssembler
— TypeLocalMatrixAssembler{IT<:Integer, T<:Number} <: AbstractArray{T, 2}
Type of "local" matrix assembler.
Local is to be understood in the sense of in the context of a single finite element. So a local matrix is an elementwise matrix which is computed entry by entry. Then it can be assembled into the global matrix in one shot.
Elfel.LocalAssemblers.LocalMatrixAssembler
— MethodLocalMatrixAssembler(nrow::IT, ncol::IT, z::T) where {IT, T}
Create a local matrix assembler, given the number of rows and columns, and the value to which the matrix should be initialized.
Elfel.LocalAssemblers.LocalVectorAssembler
— TypeLocalVectorAssembler{IT<:Integer, T<:Number} <: AbstractArray{T, 1}
Type of "local" vector assembler.
Local is to be understood in the sense of in the context of a single finite element. So a local vector is an elementwise vector which is computed entry by entry. Then it can be assembled into the global vector in one shot.
Elfel.LocalAssemblers.LocalVectorAssembler
— MethodLocalVectorAssembler(nrow::IT, z::T) where {IT, T}
Create a local vector assembler, given the number of entries, and the value to which the vector should be initialized.
Base.getindex
— MethodBase.getindex(a::A, i::Int, j::Int)
Only access to a single entry of the matrix is provided.
Base.getindex
— MethodBase.getindex(a::A, i::Int) where {A<:LocalVectorAssembler}
Access is provided to a single entry of the vector.
Base.setindex!
— MethodBase.setindex!(a::A, v, i::Int, j::Int) where {A<:LocalMatrixAssembler}
Only access to a single entry of the matrix is provided.
Base.setindex!
— MethodBase.setindex!(a::A, v, i::Int) where {A<:LocalVectorAssembler}
Access is provided to a single entry of the vector.
Base.size
— MethodBase.size(a::A) where {A<:LocalMatrixAssembler}
The size is the tuple of number of rows and number of columns.
Base.size
— MethodBase.size(a::A) where {A<:LocalVectorAssembler}
The size is the number of rows (entries).
Elfel.LocalAssemblers.init!
— Methodinit!(a::L, rdofs, cdofs) where {L<:LocalMatrixAssembler{IT, T}} where {IT, T}
Initialize the local assembler with the global degrees of freedom in the rows and columns.
The two arrays, rdofs
, cdofs
, define the global degree of freedom numbers for the element. The data matrix is zeroed out.
This function needs to be called for each new finite element.
Elfel.LocalAssemblers.init!
— Methodinit!(a::L, rdofs) where {L<:LocalVectorAssembler{IT, T}} where {IT, T}
Initialize the local assembler with the global degrees of freedom in the rows.
The array rdofs
defines the global degree of freedom numbers for the element. The data vector is zeroed out.
This function needs to be called for each new finite element.
Index
Base.IndexStyle
Base.IndexStyle
Elfel.Assemblers.AbstractSysmatAssembler
Elfel.Assemblers.AbstractSysvecAssembler
Elfel.Assemblers.SysmatAssemblerSparse
Elfel.Assemblers.SysmatAssemblerSparse
Elfel.Assemblers.SysvecAssembler
Elfel.Assemblers.SysvecAssembler
Elfel.FEIterators.FEIterator
Elfel.FESpaces.FEFields.FEField
Elfel.FESpaces.FESpace
Elfel.FElements.FE
Elfel.FElements.FEData
Elfel.LocalAssemblers.LocalMatrixAssembler
Elfel.LocalAssemblers.LocalMatrixAssembler
Elfel.LocalAssemblers.LocalVectorAssembler
Elfel.LocalAssemblers.LocalVectorAssembler
Elfel.QPIterators.QPIterator
Elfel.QPIterators.QPIterator
Elfel.RefShapes.AbstractRefShape
Elfel.RefShapes.IntegRule
Elfel.RefShapes.RefShapeCube
Elfel.RefShapes.RefShapeInterval
Elfel.RefShapes.RefShapePoint
Elfel.RefShapes.RefShapeSquare
Elfel.RefShapes.RefShapeTetrahedron
Elfel.RefShapes.RefShapeTriangle
Base.getindex
Base.getindex
Base.iterate
Base.iterate
Base.length
Base.setindex!
Base.setindex!
Base.size
Base.size
Elfel.Assemblers.assemble!
Elfel.Assemblers.assemble!
Elfel.Assemblers.assemble!
Elfel.Assemblers.assemble!
Elfel.Assemblers.assemble!
Elfel.Assemblers.assemble!
Elfel.Assemblers.finish!
Elfel.Assemblers.finish!
Elfel.Assemblers.finish!
Elfel.Assemblers.start!
Elfel.Assemblers.start!
Elfel.Assemblers.start!
Elfel.FEIterators.eldofcomps
Elfel.FEIterators.eldofentmdims
Elfel.FEIterators.eldofs
Elfel.FEIterators.eldofvals
Elfel.FEIterators.elnodes
Elfel.FEIterators.location
Elfel.FESpaces.FEFields.datadofnums
Elfel.FESpaces.FEFields.dofnums
Elfel.FESpaces.FEFields.dofnumtype
Elfel.FESpaces.FEFields.doftype
Elfel.FESpaces.FEFields.dofvals
Elfel.FESpaces.FEFields.freedofnums
Elfel.FESpaces.FEFields.gathersysvec!
Elfel.FESpaces.FEFields.gathersysvec!
Elfel.FESpaces.FEFields.gathersysvec!
Elfel.FESpaces.FEFields.highestdatadofnum
Elfel.FESpaces.FEFields.highestdatadofnum
Elfel.FESpaces.FEFields.highestfreedofnum
Elfel.FESpaces.FEFields.highestfreedofnum
Elfel.FESpaces.FEFields.ndofs
Elfel.FESpaces.FEFields.ndofs
Elfel.FESpaces.FEFields.ndofsperterm
Elfel.FESpaces.FEFields.nterms
Elfel.FESpaces.FEFields.numberdatadofs!
Elfel.FESpaces.FEFields.numberdatadofs!
Elfel.FESpaces.FEFields.numberfreedofs!
Elfel.FESpaces.FEFields.numberfreedofs!
Elfel.FESpaces.FEFields.scattersysvec!
Elfel.FESpaces.FEFields.scattersysvec!
Elfel.FESpaces.FEFields.scattersysvec!
Elfel.FESpaces.FEFields.setdofval!
Elfel.FESpaces.FEFields.setebc!
Elfel.FESpaces.FEFields.setebc!
Elfel.FESpaces.FEFields.setisdatum!
Elfel.FESpaces.dofnum
Elfel.FESpaces.doftype
Elfel.FESpaces.edofbfnum
Elfel.FESpaces.edofcompnt
Elfel.FESpaces.edofmdim
Elfel.FESpaces.makeattribute
Elfel.FESpaces.numberdofs!
Elfel.FESpaces.numberdofs!
Elfel.FESpaces.nunknowns
Elfel.FESpaces.setdofval!
Elfel.FElements.FEH1_L2
Elfel.FElements.FEH1_Q4
Elfel.FElements.FEH1_T3
Elfel.FElements.FEH1_T3_BUBBLE
Elfel.FElements.FEH1_T4
Elfel.FElements.FEH1_T6
Elfel.FElements.FEL2_Q4
Elfel.FElements.FEL2_T3
Elfel.FElements.FEL2_T4
Elfel.FElements.Jacobian
Elfel.FElements.Jacobian
Elfel.FElements.Jacobian
Elfel.FElements.Jacobian
Elfel.FElements._geometrycarrier
Elfel.FElements.bfun
Elfel.FElements.bfun
Elfel.FElements.bfungradpar
Elfel.FElements.bfungradpar
Elfel.FElements.jacjac
Elfel.FElements.jacjac
Elfel.FElements.ndofperfeat
Elfel.FElements.ndofsperel
Elfel.FElements.ndofsperel
Elfel.FElements.ndofsperel
Elfel.FElements.nfeatofdim
Elfel.FElements.refshape
Elfel.FElements.shapedesc
Elfel.LocalAssemblers.init!
Elfel.LocalAssemblers.init!
Elfel.QPIterators.bfungrad
Elfel.QPIterators.weight
Elfel.RefShapes.manifdim
Elfel.RefShapes.manifdimv
Elfel.RefShapes.quadrature
Elfel.RefShapes.quadrature
Elfel.RefShapes.quadrature
Elfel.RefShapes.quadrature
MeshCore.manifdim