Public Documentation
Documentation for MinFEM.jl
's public interface.
See the Internal page of the library for the documentation of internal types and functions.
Contents
Module
MinFEM.MinFEM
— ModuleMinFEM
A minimal finite element tool for demonstration and teaching in julia.
This package imports the following packages:
Base
Core
DocStringExtensions
LinearAlgebra
SparseArrays
WriteVTK
Types
MinFEM.Mesh
— Typemutable struct Mesh
Type for a triangular finite element mesh with volume and boundary markers.
Fields
d::Int64
: Spatial dimensionnnodes::Int64
: Number of nodesnelems::Int64
: Number of elementsnboundelems::Int64
: Number of physical (marked) boundary elementsNodes::Vector{Vector{Float64}}
: List of node coordinatesElements::Vector{Vector{Int64}}
: List of element node indicesBoundaryElements::Vector{Vector{Int64}}
: List of boundary element node indicesParentElements::Vector{Int64}
: List of parent elements to boundary elementsParentBoundaries::Vector{Int64}
: List of parent element local boundary to boundary elementsBoundaries::Dict{Int64, Boundary}
: Dictionary of marked boundariesDomains::Dict{Int64, Domain}
: Dictionary of marked volume regionsEntities::Vector{Dict{Int64, MinFEM.Entity}}
: Dictionary of physical entities
MinFEM.Region
— Typeabstract type Region
Abstract supertype for structs specifing regions of a domain, i.e. mesh. Subtypes should contain at least Name, Nodes and Elements.
MinFEM.Boundary
— Typemutable struct Boundary <: MinFEM.Region
Structure holding the name and sets of node and edge indices for one particular physical boundary.
Fields
Name::String
: Unique physical nameNodes::Set{Int64}
: List of node indicesElements::Set{Int64}
: List of element indices
MinFEM.Domain
— Typemutable struct Domain <: MinFEM.Region
Type holding the name and the set of element indices for one particular physical domain.
Fields
Name::String
: Unique physical nameNodes::Set{Int64}
: List of node indicesElements::Set{Int64}
: List of element indices
MinFEM.PDESystem
— Typemutable struct PDESystem
Structure holding all information to describe simple PDEs with Dirichlet boundary conditions.
Fields
A::SparseArrays.SparseMatrixCSC{Float64, Int64}
: Stiffness matrixb::Vector{Float64}
: Load vectorbc::Vector{Float64}
: Dirichlet boundary valuesDI::Set{Int64}
: Dirichlet nodesqdim::Int64
: Components of vector-valued stateFactors::Any
: System matrix factorizationSystemMatrix::SparseArrays.SparseMatrixCSC{Float64, Int64}
: System of stiffness matrix and Dirichlet conditinonsB::SparseArrays.SparseMatrixCSC{Float64, Int64}
: Dirichlet projection matrixstate::Vector{Float64}
: Solution vectorrhs::Vector{Float64}
: Right hand side vector with source term and dirichlet conditions
Functions and Methods
Mesh Generation
MinFEM.unit_interval
— Functionunit_interval(n::Int64) -> Mesh
Returns an n nodes quasi-uniform mesh for the 1D unit interval. The left boundary node is denoted by 1001 and the right one by 1002.
MinFEM.unit_square
— Functionunit_square(n::Int64) -> Mesh
Returns a n-by-n quasi-uniform mesh for the 2D unit square. The boundary indices are given in the order bottom, top, left, right from 1001 to 1004.
MinFEM.import_mesh
— Functionimport_mesh(fileName::String) -> Mesh
Returns a mesh imported from a gmsh file of version v1, v2 or v4.
MinFEM.export_mesh
— Functionexport_mesh(mesh::Mesh, fileName::String)
Exports a mesh to a gmsh file of version v2.
MinFEM.update_mesh!
— Functionupdate_mesh!(
mesh::Mesh,
c::Vector{Vector{Float64}}
) -> Vector{Vector{Float64}}
Updates given mesh by shifting all nodes to new coordinates c.
MinFEM.deform_mesh!
— Functiondeform_mesh!(mesh::Mesh, v::AbstractVector{Float64}; t::Float64=1.0)
Deforms given mesh by shifting all nodes according to the vector field v scaled by the stepsize t.
MinFEM.deform_mesh
— Functiondeform_mesh(mesh::Mesh, v::AbstractVector{Float64}; t::Float64=1.0) -> Mesh
Returns copy of mesh deformed by deform_mesh!
.
Type Handling
MinFEM.select_boundaries
— Functionselect_boundaries(mesh::Mesh, args...) -> Set{Boundary}
Returns set of all or specified physical boundaries of the mesh.
MinFEM.select_domains
— Functionselect_domains(mesh::Mesh, args...) -> Set{Domain}
Returns set of all or specified physical boundaries of the mesh.
MinFEM.extract_elements
— Functionextract_elements(boundaries::Set{Boundary}) -> Set{Int64}
Returns set of boundary element ids in set of physical boundaries.
extract_elements(domains::Set{Domain}) -> Set{Int64}
Returns set of boundary element ids in set of physical domains.
MinFEM.extract_nodes
— Functionextract_nodes(boundaries::Set{Boundary}) -> Set{Int64}
Returns set of node ids in set of physical boundaries.
extract_nodes(domains::Set{Domain}) -> Set{Int64}
Returns set of node ids in set of physical domains.
Function Discretization
MinFEM.evaluate_mesh_function
— Functionevaluate_mesh_function(
mesh::Mesh,
f::Function,
region::Set{Boundary};
qdim::Int64 = 1
) -> Vector{Float64}
Returns evaluation of a given function object f on all or specified nodes of the mesh. Can be either called with set of physical boundaries or directly with a set of nodes when given with keyword argument region.
evaluate_mesh_function(
mesh::Mesh,
f::Function,
region::Set{Boundary};
qdim::Int64 = 1
) -> Vector{Float64}
Same as previous evaluatemeshfunction(...) but takes Set{Boundary}
as mandatory argument for the relevant region, which replaces the keyword argument using Set{Int64}
. Extracts nodes from the boundary and then passes them to the base function.
MinFEM.evaluate_function
— Functionevaluate_function(
f::Vector{Float64},
node::Int64;
qdim::Int64 = 1
) -> Vector{Float64}
Evaluates discrete function give as FEM coefficient vector at a given node. Resulting vector has the lenght of the number of components of f.
evaluate_function(
f::Vector{Float64},
region::Set{Int64};
qdim::Int64 = 1
) -> Vector{Float64}
Evaluates discrete function give as FEM coefficient vector at all nodes in a given region. Resulting vector has the lenght of the full coefficient vector with the entries corresponding to non-evaluated nodes being zero.
MinFEM.evaluate_quadrature_function
— Functionevaluate_quadrature_function(
mesh::Mesh,
f::Function;
region::Set{Int64} = Set{Int64}(),
order::Int64 = 1,
qdim::Int64 = 1
) -> Vector{Float64}
Returns vector of function f evaluated at each elements quadrature points of the given mesh for the specified order.
evaluate_quadrature_function(
mesh::Mesh,
f::Function,
region::Set{Domain};
order::Int64 = 1,
qdim::Int64 = 1
) -> Vector{Float64}
Same as previous evaluatequadraturefunction(...) but takes Set{Domain}
as mandatory argument for the relevant region, which replaces the keyword argument using Set{Int64}
. Extracts nodes from the domain and then passes them to the base function.
MinFEM.evaluate_quadrature_function_boundary
— Functionevaluate_quadrature_boundary(
mesh::Mesh,
f::Function;
region::Set{Int64} = Set{Int64}(),
order::Int64 = 1,
qdim::Int64 = 1
) -> Vector{Float64}
Returns vector of function f with qdim components evaluated at each or specified boundary elements quadrature points of the given mesh for the specified order.
evaluate_quadrature_function_boundary(
mesh::Mesh,
f::Function,
region::Set{Boundary};
order::Int64 = 1,
qdim::Int64 = 1
) -> Vector{Float64}
Same as previous evaluatequadraturefunction_boundary(...) but takes Set{Boundary}
as mandatory argument for the relevant region, which replaces the keyword argument using Set{Int64}
. Extracts elements from the boundary and then passes them to the base function.
Mesh (Element) Properties
MinFEM.jacobian
— Functionjacobian(
coords::Vector{Vector{Float64}}
) -> Tuple{Float64, LinearAlgebra.Adjoint{Float64, Matrix{Float64}}}
Returns determinant (i.e. element weight) and inverse transposed of the jacobian of an FEM element spanned by nodes at the given coordinates.
Commonly the coordinates correspond to one element in a mesh, but not necessarily have to.
jacobian(
mesh::Mesh,
nodes::Vector{Int64}
) -> Tuple{Float64, LinearAlgebra.Adjoint{Float64, Matrix{Float64}}}
Same as previous jacobian(...)
, but takes a mesh and set of nodes as arguments. Extracts coordinates from the mesh and passes them to the base function.
Commonly the coordinates correspond to one element in a mesh, but not necessarily have to.
jacobian(
mesh::Mesh,
element::Int64
) -> Tuple{Float64, LinearAlgebra.Adjoint{Float64, Matrix{Float64}}}
Same as previous jacobian(...)
, but takes a mesh and an element id as arguments. Extracts coordinates of the support nodes from the mesh and passes them to the base function.
MinFEM.jacobian_boundary
— Functionjacobian_boundary(
coords::Vector{Vector{Float64}}
) -> Union{Nothing, Float64, Int64}
Returns determinant of the jacobian (i.e. element weight) of an FEM boundary element (in d-1 dimensions) spanned spanned by nodes at the given coordinates.
Commonly the coordinates correspond to an element in a mesh, but not necessarily have to.
jacobian_boundary(
mesh::Mesh,
nodes::Vector{Int64}
) -> Union{Nothing, Float64, Int64}
Same as previous jacobian_boundary(...)
, but takes a mesh and set of nodes as arguments. Extracts coordinates from the mesh and passes them to the base function.
Commonly the coordinates correspond to one element in a mesh, but not necessarily have to.
jacobian_boundary(
mesh::Mesh,
element::Int64
) -> Union{Nothing, Float64, Int64}
Same as previous jacobian_boundary(...)
, but takes a mesh and an element id as arguments. Extracts coordinates of the support nodes from the mesh and passes them to the base function.
MinFEM.elementvolume
— Functionelementvolume(d::Int64) -> Float64
Returns volume of the d-dimensional reference element.
elementvolume(mesh::Mesh, element::Int64) -> Float64
Returns volume of the given element in the given mesh.
elementvolume(mesh::Mesh) -> Vector{Float64}
Returns vector of volumes of all elements in the given mesh.
MinFEM.elementvolume_boundary
— Functionelementvolume_boundary(d::Int64) -> Float64
Returns volume of the (d-1)-dimensional boundary reference element.
elementvolume_boundary(
mesh::Mesh,
element::Int64
) -> Float64
Returns volume of the given boundary element in the given mesh.
elementvolume_boundary(mesh::Mesh) -> Vector{Float64}
Returns vector of volumes of all boundary elements in the given mesh.
MinFEM.elementbarycenter
— Functionelementbarycenter(d::Int64) -> Vector{Float64}
Returns barycenter vector of the d-dimensional reference element.
elementbarycenter(
mesh::Mesh,
element::Int64
) -> Vector{Float64}
Returns barycenter of the given element in the given mesh.
elementbarycenter(mesh::Mesh) -> Vector{Vector{Float64}}
Returns vector of all barycenters of all element in the given mesh.
MinFEM.elementdiameter
— Functionelementdiameter(
coords::Vector{Vector{Float64}}
) -> Union{Float64, Int64}
Returns diameter (i.e. longest edge length) of an element spanned by nodes at the given coordinates.
The coordinates do not necessarily have to span an element of full dimension. Could also be used for boundary elements.
elementdiameter(
mesh::Mesh,
nodes::Vector{Int64}
) -> Union{Float64, Int64}
Same as previous elementdiameter(...)
, but takes a mesh and set of nodes as arguments. Extracts coordinates from the mesh and passes them to the base function.
Commonly the coordinates correspond to one (boundary) element in a mesh, but not necessarily have to.
elementdiameter(
mesh::Mesh,
element::Int64
) -> Union{Float64, Int64}
Same as previous elementdiameter(...)
, but takes a mesh and an element id as arguments. Extracts coordinates of the support nodes from the mesh and passes them to the base function.
elementdiameter(mesh::Mesh) -> Vector{Float64}
Similar ot previous elementdiameter(...)
, but returns vector of elementdiameters for all elements in the given mesh.
MinFEM.elementdiameter_boundary
— Functionelementdiameter_boundary(
coords::Vector{Vector{Float64}}
) -> Union{Float64, Int64}
Returns diameter (i.e. longest edge length) of a boundary element spanned by nodes at the given coordinates.
Commonly the coordinates correspond to one boundary element in a mesh, but not necessarily have to.
Since the number of nodes to span an element is not prescribed, the functionality is identical to elementdiameter(args)
. However this function is kept for consistency and improved readability.
elementdiameter_boundary(
mesh::Mesh,
nodes::Vector{Int64}
) -> Union{Float64, Int64}
Same as previous elementdiameter_boundary(...)
, but takes a mesh and set of nodes as arguments. Extracts coordinates from the mesh and passes them to the base function.
Commonly the coordinates correspond to one (boundary) element in a mesh, but not necessarily have to.
Since the number of nodes to span an element is not prescribed, the functionality is identical to elementdiameter(args)
. However this function is kept for consistency and improved readability.
elementdiameter_boundary(
mesh::Mesh,
element::Int64
) -> Union{Float64, Int64}
Same as previous elementdiameter_boundary(...)
, but takes a mesh and a boundary element id as arguments. Extracts coordinates of the support nodes from the mesh and passes them to the base function.
elementdiameter_boundary(mesh::Mesh) -> Vector{Float64}
Similar ot previous elementdiameter_boundary(...)
, but returns vector of elementdiameters for all elements in the given mesh.
MinFEM.elementratio
— Functionelementratio(
coords::Vector{Vector{Float64}}
) -> Union{Float64, Int64}
Returns ratio of inscribed to circumscribed circle or sphere for an element spanned by nodes at the given coordinates.
elementratio(mesh::Mesh) -> Vector{Float64}
Same as previous elementratio(...)
, but returns vector of ratios for all elements in the given mesh.
MinFEM.elementangle
— Functionelementangle(
coords::Vector{Vector{Float64}}
) -> Union{Irrational{:π}, Float64}
Returns smallest interior angle for an element spanned by nodes at the given coordinates.
elementangle(mesh::Mesh) -> Vector{Float64}
Same as previous elementangle(...)
, but returns vector of angles for all elements in the given mesh.
MinFEM.outernormalvector
— Functionouternormalvector(dim::Int64, ii::Int64) -> Vector{Float64}
Returns the outer normal vector at the ii-th boundary of the dim-dimensional reference element.
outernormalvector(
mesh::Mesh,
boundaryElement::Int64,
J::AbstractMatrix{Float64}
) -> Any
Returns the outer normal vector at a boundary element of the given mesh using a pre-computed jacobian transformation matrix of the parent element.
outernormalvector(
mesh::Mesh,
boundaryElement::Int64
) -> Vector{Float64}
Returns the outer normal vector at the given boundary element of the given mesh.
outernormalvector(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}()
) -> Vector{Float64}
Returns coefficient vector of outer normal vectors at all or specified boundary elements of the given mesh.
MinFEM.stripwidth
— Functionstripwidth(mesh::Mesh) -> Any
Return width L of a strip that the meshed domain fits into.
MinFEM.boundingbox
— Functionboundingbox(mesh::Mesh) -> Vector{Array}
Return two nodes which span the bounding box of the mesh.
MinFEM.volume
— Functionvolume(mesh::Mesh) -> Union{Float64, Int64}
Returns the d-dimenional volume of the domain definded by the mesh.
MinFEM.barycenter
— Functionbarycenter(mesh::Mesh) -> Vector{Float64}
Returns the d-dimenional domain of the domain definded by the mesh.
Local Quadrature
MinFEM.quadrature_points
— Functionquadrature_points(
d::Int64,
order::Int64
) -> Union{Nothing, Vector{Vector{Float64}}}
Returns coordinates of the Gauss-Legendre quadrature points on the d-dimensional FEM reference element for exact integration of polynomials up to the given order.
quadrature_points(
mesh::Mesh,
order::Int64
) -> Vector{Vector{Float64}}
Returns global coordinates of the Gauss-Legendre quadrature points on each finite element in the given mesh for exact integration of polynomials up to the given order.
quadrature_points(
element::Int64,
mesh::Mesh,
order::Int64
) -> Vector{Vector{Float64}}
Returns global coordinates of the Gauss-Legendre quadrature points on a specified element in the given mesh for exact integration of polynomials up to the given order.
MinFEM.quadrature_points_boundary
— Functionquadrature_points_boundary(
d::Int64,
order::Int64
) -> Union{Nothing, Vector{Vector{Float64}}}
Returns coordinates of the Gauss-Legendre quadrature points on the (d-1)-dimensional FEM reference element for exact integration of polynomials up to the given order.
quadrature_points_boundary(
mesh::Mesh,
order::Int64;
boundaryElements::Set{Int64} = Set{Int64}()
)-> Array{Array{Float64,1},1}
Returns global coordinates of the Gauss-Legendre quadrature points on each boundary element in the given mesh for exact integration of polynomials up to the given order.
quadrature_points_boundary(
element::Int64,
mesh::Mesh,
order::Int64
) -> Vector{Vector{Float64}}
Returns global coordinates of the Gauss-Legendre quadrature points on a specified boundary element in the given mesh for exact integration of polynomials up to the given order.
MinFEM.quadrature_weights
— Functionquadrature_weights(
d::Int64,
order::Int64
) -> Union{Nothing, Vector{Float64}}
Returns weights of the Gauss-Legendre quadrature points on the d-dimensional FEM reference element for exact integration of polynomials up to the given order.
MinFEM.quadrature_weights_boundary
— Functionquadrature_weights_boundary(
d::Int64,
order::Int64
) -> Union{Nothing, Vector{Float64}}
Returns weights of the Gauss-Legendre quadrature points on the (d-1)-dimensional FEM reference element for exact integration of polynomials up to the given order.
MinFEM.quadrature_order
— Functionquadrature_order(d::Int64, n::Int64) -> Int64
Returns highest quadrature order archived by using n points in d dimensions.
MinFEM.integral_over_reference_element
— Functionintegral_over_reference_element(
f::Function,
d::Int64;
order::Int64 = 1
) -> Float64
Returns integral of function f over the d-dimensional FEM reference element computed with Gauss-Legendre quadrature exact at least for polynomials up to given order.
Vector-Valued Coefficients Handling
MinFEM.prolong_multivector
— Functionprolong_multivector(
x::AbstractVector{Float64},
m::Int64,
qdim::Int64;
block::Int64 = 1
) -> Vector{Float64}
Prolongates a vector of m blocks of size block to a multivector for qdim components of length qdim×block×m.
MinFEM.restrict_multivector
— Functionrestrict_multivector(
x::AbstractVector{Float64},
m::Int64,
qdim::Int64;
block::Int64 = 1
) -> Vector{Float64}
Restricts a multivector of qdim×block×m elements for qdim components to the regular vector of m blocks of size block.
Weight Computation
MinFEM.assemble_weightmultivector
— Functionassemble_weightmultivector(
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> Vector{Float64}
Returns the vector of weights for all elements of mesh and qdim components. Weight is equal to volume if 1st order integration by mid-point rule is used.
MinFEM.assemble_weightmultivector_boundary
— Functionassemble_weightmultivector_boundary(
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> Vector{Float64}
Returns the vector of weights for all boundary elements of mesh and qdim components. Weight is equal to volume if 1st order integration by mid-point rule is used.
Finite Element Basis Functions
MinFEM.phi
— Functionphi(ii::Int64, x::AbstractVector{T} where T) -> Any
Returns ii-th local basis function evaluated at x. Dimension is determined by length(x).
MinFEM.grad_phi
— Functiongrad_phi(d::Int64, ii::Int64) -> Matrix{Float64}
Gradient of ii-th local basis functions in d dimensions. Note that the gradient is constant on an element.
FE Operator Assembly
MinFEM.assemble_derivativematrix
— Functionassemble_derivativematrix(mesh::Mesh; qdim::Int64=1) -> SparseMatrixCSC{Float64, Int64}
Returns the discrete derivative matrix for all elements of mesh and qdim components.
MinFEM.assemble_laplacian
— Functionassemble_laplacian(
mesh::Mesh;
qdim::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the Laplacian stiffness matrix for all elements of mesh and qdim components.
assemble_laplacian(
D::SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseMatrixCSC{Float64, Int64}
Returns the Laplacian stiffness matrix for all elements of a mesh and qdim components.
Takes existing derivative matrix D and weight vector w as arguments. Can yield performance benefits compared to direct assembly if D already exists. Note that number of components and quadrature order in D and w have to match.
MinFEM.assemble_derivativematrix_boundary
— Functionassemble_derivativematrix_boundary(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the discrete derivative matrix for all or specied boundary elements of mesh and qdim components.
The implementation is based on the derivate tensor for the corresponding full element and the fact that the gradient ist constant on the element for linear elements.
MinFEM.assemble_normalderivativematrix
— Functionassemble_normalderivativematrix(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the discrete normal derivative matrix for all or specied boundary elements of mesh and qdim components.
The implementation is based on the derivate tensor for the corresponding full element and the fact that the gradient ist constant on the element for linear elements.
MinFEM.assemble_basismatrix
— Functionassemble_basismatrix(
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 3
) -> SparseMatrixCSC{Float64, Int64}
Returns the discrete basis matrix with given local integration order for all elements of mesh and qdim components.
MinFEM.assemble_massmatrix
— Functionassemble_massmatrix(
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all elements of mesh and qdim components.
assemble_massmatrix(
E::SparseArrays.SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all elements of mesh and qdim components.
Takes existing basis matrix E and weight vector w as arguments. Can yield performance benefits compared to direct assembly if E already exists. Note that number of components and quadrature order in E and w have to match.
MinFEM.assemble_basismatrix_boundary
— Functionassemble_basismatrix_boundary(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the discrete basis matrix with given local integration order for all or specified boundary elements of mesh and qdim components.
MinFEM.assemble_massmatrix_boundary
— Functionassemble_massmatrix_boundary(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all or specified boundary elements of mesh and image dimension qdim components.
assemble_massmatrix_boundary(
E::SparseArrays.SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all or specified boundary elements of mesh and image dimension qdim components.
Takes existing boundary basis matrix E and weight vector w as arguments. Can yield performance benefits compared to direct assembly if E already exists. Note that number of components and quadrature order in E and w have to match.
MinFEM.assemble_cubicterm
— Functionassemble_cubicterm(
mesh::Mesh,
y::AbstractVector;
order::Int64 = 3
) -> SparseMatrixCSC{Float64, Int64}
Returns the cubic term of the standard semilinear parabolic equation.
MinFEM.assemble_cubicderivativematrix
— Functionassemble_cubicderivativematrix(
mesh::Mesh,
y::AbstractVector;
order::Int64 = 3
) -> SparseMatrixCSC{Float64, Int64}
Returns the linearization of the cubic term of the standard semilinear elliptic equation.
MinFEM.assemble_cubicsecondderivativematrix
— Functionassemble_cubicsecondderivativematrix(
mesh::Mesh,
y::AbstractVector,
p::AbstractVector;
order::Int64 = 3
) -> SparseMatrixCSC{Float64, Int64}
Returns the second derivative of the cubic term of the standard semilinear elliptic equation around the state y.
MinFEM.assemble_elasticity
— Functionassemble_elasticity(
mesh::Mesh,
lambda::Float64,
mu::Float64
) -> SparseArrays.SparseMatrixCSC{Float64, Int64}
Returns the linear elasticity stiffness matrix with constant coefficients λ and μ for all elements of mesh and image dimension qdim.
Boundary Condition Handling
MinFEM.assemble_dirichletcondition!
— Functionassemble_dirichletcondition!(
A::SparseMatrixCSC{Float64, Int64},
DI::Set{Int64};
rhs = [],
bc = [],
qdim::Int64 = 1,
insert = 1.0
)
Modify a stiffness matrix and a right hand side according to the given Dirichlet conditions.
DI has to be the set of node indices for which the condition should be active. For vector valued states either DI can be set to each component that should have a Dirichlet condtion or qdim is set, if all components should have the condition. If the rhs shall be updated bc needs to be specified explicitly as well. The value insert is put as diagonal element. Usually you want a 1.0 here.
assemble_dirichletcondition!(
A::SparseMatrixCSC{Float64, Int64},
DI::Set{Boundary};
rhs = [],
bc = [],
qdim::Int64 = 1,
insert = 1.0
)
Same as previous assemble_dirichletcondition!(...)
, but takes Set{Boundary}
instead of Set{Int64}
as argument for the Dirichlet nodes. Thus extracts the nodes first and then passes them to the base function.
MinFEM.assemble_dirichletcondition_rhs!
— Functionassemble_dirichletcondition_rhs!(
A::SparseMatrixCSC{Float64, Int64},
rhs::AbstractVector{Float64},
DI::Set{Int64},
bc::AbstractVector{Float64};
qdim::Int64 = 1
)
Modify a right hand side according to the given Dirichlet conditions. Behaviour is similar to assemble_dirichletcondition!(...)
however the system matrix A is not updated. Can be relevant for iterative algorithm, where the system matrix is constant and only the right hand side changes. Then one can store the modified matrix and only assemble the right hand side in every iteration.
DI has to be the set of node indices for which the condition should be active. For vector valued states either DI can be set to each component that should have a Dirichlet condtion or qdim is set, if all components should have the condition.
assemble_dirichletcondition_rhs!(
A::SparseMatrixCSC{Float64, Int64},
rhs::AbstractVector{Float64},
DI::Set{Boundary},
bc::AbstractVector{Float64};
qdim::Int64 = 1
)
Same as previous assemble_dirichletcondition_rhs!(...)
, but takes Set{Boundary}
instead of Set{Int64}
as argument for the Dirichlet nodes. Thus extracts the nodes first and then passes them to the base function.
MinFEM.assemble_dirichletprojection
— Functionassemble_dirichletprojection(
ncoeffs::Int64,
DI::Set{Int64};
qdim::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
Returns the projection matrix onto the Dirichlet nodes where the input ncoeffs is understood as qdim*nnodes.
FE Extensions
MinFEM.pnorm
— Functionpnorm(
p::Float64,
v::AbstractVector{Float64},
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> Float64
Returns $L^p$-norm of a FEM coefficient vector v over the nodes or quadrature points in the elements on the given mesh. Note that for p=Inf, the values do not get inter- or extrapolated and thus it might not be true maximum over the domain, but only over the given data.
MinFEM.qnorm
— Functionqnorm(
p::Float64,
v::AbstractVector{Float64},
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> Float64
Returns $L^q$-norm of a FEM coefficient vector v over the nodes or quadrature points in the elements on the given mesh with q being the conjugated exponent to p.
MinFEM.twonorm
— Functiontwonorm(
v::AbstractVector{Float64},
M::AbstractMatrix{Float64}
) -> Any
Implementation of the discrete $L^2$-norm based on a mass matrix. Faster version of pnorm for p=2 by using that direct assembly of mass matrix instead of basis matrices is possible. Can be used for domain or boundary depending on the respective mass matrix.
twonorm(
v::AbstractVector{Float64},
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 3
)
Same as previous twonorm(...), but takes mesh and number of components as arguments. Then assembles the mass matrix and passes it to the base function.
MinFEM.pnorm_boundary
— Functionpnorm_boundary(
p::Float64,
v::AbstractVector{Float64},
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 1
) -> Float64
Returns $L^p$-norm of a FEM coefficient vector v over the nodes or quadrature points in the boundary elements on the given mesh. Note that for p=Inf, the values do not get inter- or extrapolated and thus it might not be true maximum over the domain, but only over the given data.
MinFEM.qnorm_boundary
— Functionqnorm_boundary(
p::Float64,
v::AbstractVector{Float64},
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 1
) -> Float64
Returns $L^q$-norm of a function f or its FEM coefficient vector v over the nodes or quadrature points in the boundary elements on the given mesh with q being the conjugated exponent to p.
MinFEM.twonorm_boundary
— Functiontwonorm_boundary(
v::AbstractVector{Float64},
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 3
)
Implementation of the discrete $L^2$-norm on the boundary based on a mass matrix. Faster version of pnorm for p=2 by using that direct assembly of mass matrix instead of basis matrices is possible. Assembles boundary mass matrix form mesh, boundary elements and the number of components. Then passes it to the base implementation of twonorm
, which can be used for boundary terms as well.
MinFEM.conjugated_exponent
— Functionconjugated_exponent(p::Float64) -> Float64
Returns conjugated exponent q to p in the Hölder sense 1/p + 1/q = 1.
PDE System Handling
MinFEM.solve!
— Functionsolve!(S::PDESystem) -> Any
First tries to assemble! the system matrix with multipliers for Dirichlet conditions. If the system has already been used before, this step is skipped. This is determined depending on an existing factorization of the system matrix. If the stiffness matrix or Dirichlet conditions have changes, one should invole refresh! first. Finally the system is solved via matrix factorization.
MinFEM.assemble!
— Functionassemble!(
S::PDESystem
) -> Union{Nothing, SuiteSparse.UMFPACK.UmfpackLU{Float64, Int64}}
If the system has not been used before, sets up the system matrix with multipliers for Dirichlet conditions and factorizes it.
MinFEM.refresh!
— Functionrefresh!(
S::PDESystem
) -> Union{Nothing, SuiteSparse.UMFPACK.UmfpackLU{Float64, Int64}}
Recomputes the factorization of the stiffness matrix using assemble!.
Output
MinFEM.write_to_vtk
— Functionwrite_to_vtk(
x::Vector{Vector{Float64}},
mesh::Mesh,
data_names::Array{String},
file_name::String,
qdim::Array{Int64}
)
Writes multiple data vectors with corresponding names and qdims to a .vtk-file with the given name.
write_to_vtk(
x::Vector{Vector{Float64}},
mesh::Mesh,
data_names::Array{String},
file_name::String;
qdim::Int64 = 1
)
Same as pervious write_to_vtk(...)
, but can handles simplification by submitting one qdim for all data vectors if they are the same. If everything is scalar no specification of qdim is required.
write_to_vtk(
x::Vector{Float64},
mesh::Mesh,
data_name::String,
file_name::String;
qdim::Int64 = 1
)
Same as pervious write_to_vtk(...)
, but handles setting for a single data vector with corresponding name and qdim. Forwards them to base implementation as array with one entry.
MinFEM.write_to_vtk_boundary
— Functionwrite_to_vtk_boundary(
x::Vector{Vector{Float64}},
mesh::Mesh,
data_names::Array{String},
file_name::String,
qdim::Array{Int64};
boundary=Set{Boundary}()
)
Writes multiple data vectors with corresponding names and qdims based on the boundary elements of the mesh to a .vtk-file with the given name.
write_to_vtk_boundary(
x::Vector{Vector{Float64}},
mesh::Mesh,
data_names::Array{String},
file_name::String;
boundary = Set{Boundary}(),
qdim::Int64 = 1
)
Same as pervious write_to_vtk_boundary(...)
, but can handles simplification by submitting one qdim for all data vectors if they are the same. If everything is scalar no specification of qdim is required.
write_to_vtk_boundary(
x::Vector{Float64},
mesh::Mesh,
data_name::String,
file_name::String;
boundary = Set{Boundary}(),
qdim::Int64 = 1
)
Same as pervious write_to_vtk_boundary(...)
, but handles setting for a single data vector with corresponding name and qdim. Forwards them to base implementation as array with one entry.
MinFEM.open_vtkfile
— Functionopen_vtkfile(
mesh::Mesh,
file_name::String
) -> WriteVTK.DatasetFile
Open a new VTK output file and write the mesh data into it.
MinFEM.open_vtkfile_boundary
— Functionopen_vtkfile_boundary(
mesh::Mesh,
file_name::String,
boundaryElements::Set{Int64}
) -> WriteVTK.DatasetFile
Open a new VTK output file and write the mesh data into it.
open_vtkfile_boundary(
mesh::Mesh,
file_name::String;
boundary = Set{Boundary}()
)
Open a new VTK output file and write the mesh data into it.
MinFEM.save_vtkfile
— Functionsave_vtkfile(
vtkfile::WriteVTK.DatasetFile
) -> Vector{String}
Finalize a VTK file by writing all data to disk.
MinFEM.write_pointdata_vtkfile!
— Functionwrite_pointdata_vtkfile!(vtkfile::WriteVTK.DatasetFile, data::Any, data_name::String)
Add a new point data field with a name to an existing VTK file.
MinFEM.write_celldata_vtkfile!
— Functionwrite_celldata_vtkfile!(
vtkfile::WriteVTK.DatasetFile,
data,
data_name::String
) -> LightXML.XMLElement
Add a new cell data field with a name to an existing VTK file.
MinFEM.write_to_txt
— Functionwrite_to_txt(
x::Vector{Float64},
mesh::Mesh,
file_name::String;
qdim::Int64 = 1
)
Writes a coefficient vector x with qdim components based on the nodes of mesh to the given file.
MinFEM.read_from_txt
— Functionread_from_txt(
file_name::String
) -> Tuple{Vector{Array{Float64, N} where N}, Vector{Float64}}
Reads node coordinates and finite element coefficient vector from file generated by write_to_txt.
Index
MinFEM.Boundary
MinFEM.Domain
MinFEM.Mesh
MinFEM.PDESystem
MinFEM.Region
MinFEM.assemble!
MinFEM.assemble_basismatrix
MinFEM.assemble_basismatrix_boundary
MinFEM.assemble_cubicderivativematrix
MinFEM.assemble_cubicsecondderivativematrix
MinFEM.assemble_cubicterm
MinFEM.assemble_derivativematrix
MinFEM.assemble_derivativematrix_boundary
MinFEM.assemble_dirichletcondition!
MinFEM.assemble_dirichletcondition_rhs!
MinFEM.assemble_dirichletprojection
MinFEM.assemble_elasticity
MinFEM.assemble_laplacian
MinFEM.assemble_massmatrix
MinFEM.assemble_massmatrix_boundary
MinFEM.assemble_normalderivativematrix
MinFEM.assemble_weightmultivector
MinFEM.assemble_weightmultivector_boundary
MinFEM.barycenter
MinFEM.boundingbox
MinFEM.conjugated_exponent
MinFEM.deform_mesh
MinFEM.deform_mesh!
MinFEM.elementangle
MinFEM.elementbarycenter
MinFEM.elementdiameter
MinFEM.elementdiameter_boundary
MinFEM.elementratio
MinFEM.elementvolume
MinFEM.elementvolume_boundary
MinFEM.evaluate_function
MinFEM.evaluate_mesh_function
MinFEM.evaluate_quadrature_function
MinFEM.evaluate_quadrature_function_boundary
MinFEM.export_mesh
MinFEM.extract_elements
MinFEM.extract_nodes
MinFEM.grad_phi
MinFEM.import_mesh
MinFEM.integral_over_reference_element
MinFEM.jacobian
MinFEM.jacobian_boundary
MinFEM.open_vtkfile
MinFEM.open_vtkfile_boundary
MinFEM.outernormalvector
MinFEM.phi
MinFEM.pnorm
MinFEM.pnorm_boundary
MinFEM.prolong_multivector
MinFEM.qnorm
MinFEM.qnorm_boundary
MinFEM.quadrature_order
MinFEM.quadrature_points
MinFEM.quadrature_points_boundary
MinFEM.quadrature_weights
MinFEM.quadrature_weights_boundary
MinFEM.read_from_txt
MinFEM.refresh!
MinFEM.restrict_multivector
MinFEM.save_vtkfile
MinFEM.select_boundaries
MinFEM.select_domains
MinFEM.solve!
MinFEM.stripwidth
MinFEM.twonorm
MinFEM.twonorm_boundary
MinFEM.unit_interval
MinFEM.unit_square
MinFEM.update_mesh!
MinFEM.volume
MinFEM.write_celldata_vtkfile!
MinFEM.write_pointdata_vtkfile!
MinFEM.write_to_txt
MinFEM.write_to_vtk
MinFEM.write_to_vtk_boundary