Public Documentation
Documentation for MinFEM.jl
's public interface.
See the Internals 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.
Types
MinFEM.Mesh
— TypeMesh
Type for a triangular finite element mesh with volume and boundary markers.
MinFEM.Region
— TypeRegion
Abstract type for structs specifing regions of a domain, i.e. mesh. Subtypes should contain at least Name, Nodes and Elements.
MinFEM.Boundary
— TypeBoundary
Structure holding the name and sets of node and edge indices for one particular physical boundary.
MinFEM.Domain
— TypeDomain
Type holding the name and the set of element indices for one particular physical domain.
MinFEM.PDESystem
— TypePDESystem
Structure holding all information to describe simple PDEs with Dirichlet boundary 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::Array{Array{Float64,1},1})
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}
extract_elements(domains::Set{Domain}) -> Set{Int64}
Returns set of boundary element ids in set of physical regions.
MinFEM.extract_nodes
— Functionextract_nodes(boundaries::Set{Boundary}) -> Set{Int64}
extract_nodes(domains::Set{Domain}) -> Set{Int64}
Returns set of node ids in set of physical regions.
Function Discretization
MinFEM.evaluate_mesh_function
— Functionevaluate_mesh_function(mesh::Mesh, f::Function;
region=Set{Int64}(), qdim=1) -> Vector{Float64}
evaluate_mesh_function(mesh::Mesh, f::Function, region::Set{Boundary};
qdim = 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.
MinFEM.evaluate_quadrature_function
— Functionevaluate_quadrature_function(
mesh::Mesh,
f::Function;
region::Set{Int64} = Set{Int64}(),
order::Int64 = 1,
qdim::Int64 = 1
)
Returns vector of function f evaluated at each elements quadrature points of the given mesh for the specified order.
MinFEM.evaluate_quadrature_function_boundary
— Functionevaluate_quadrature_boundary(
mesh::Mesh,
f::Function;
region::Set{Int64} = Set{Int64}(),
order::Int64 = 1,
qdim::Int64 = 1
)
evaluate_quadrature_function_boundary(
mesh::Mesh,
f::Function,
region::Set{Boundary};
order::Int64 = 1,
qdim::Int64 = 1
)
Returns vector of function f evaluated at each boundary elements quadrature points of the given mesh for the specified order.
Mesh (Element) Properties
MinFEM.jacobian
— Functionjacobian(coords::Vector{Vector{Float64}}) -> Float64, Matrix{Float64}
jacobian(mesh::Mesh, nodes::Array{Int64,1}) -> Float64, Matrix{Float64}
jacobian(mesh::Mesh, element::Int64) -> Float64, Matrix{Float64}
Returns determinant (i.e. element weight) and inverse transposed of the jacobian of an FEM element spanned by the given nodes.
MinFEM.jacobian_boundary
— Functionjacobian_boundary(coords::Vector{Vector{Float64}}) -> Float64
jacobian_boundary(mesh::Mesh, nodes::Array{Int64,1}) -> Float64
jacobian_boundary(mesh::Mesh, element::Int64) -> Float64
Returns determinant of the jacobian (i.e. element weight) of an FEM boundary element (in d-1 dimensions) spanned by the given nodes.
MinFEM.elementvolume
— Functionelementvolume(mesh::Mesh) -> Vector{Float64}
elementvolume(mesh::Mesh, element::Int64) -> Float64
elementvolume(d::Int64) -> Float64
Returns volume of all or one element(s) in a mesh or of the d-dimensional reference element.
MinFEM.elementvolume_boundary
— Functionelementvolume_boundary(mesh::Mesh) -> Vector{Float64}
elementvolume_boundary(mesh::Mesh, element::Int64) -> Float64
elementvolume_boundary(d::Int64) -> Float64
Returns volume of all or one boundary element(s) in a mesh.
MinFEM.elementbarycenter
— Functionelementbarycenter(mesh::Mesh) -> Array{Array{Float64,1},1}
elementbarycenter(mesh::Mesh, element::Int64) -> Array{Float64,1}
elementbarycenter(d::Int64) -> Array{Float64,1}
Returns barycenter of all or one element(s) in a mesh or of the d-dimensional reference element.
MinFEM.elementdiameter
— Functionelementdiameter(mesh::Mesh) -> Vector{Float64}
elementdiameter(mesh::Mesh, element::Int64) -> Float64
elementdiameter(mesh::Mesh, nodes::Vector{Int64}) -> Float64
elementdiameter(nodes::Vector{Vector{Float64}}) -> Float64
Returns diameter (i.e. longest edge length) of all or one element(s) in a mesh or of an element spanned by the given nodes or coordinates (not necessarily of full dimension).
MinFEM.elementdiameter_boundary
— Functionelementdiameter_boundary(mesh::Mesh) -> Vector{Float64}
elementdiameter_boundary(mesh::Mesh, element::Int64) -> Float64
Returns diameter (i.e. longest edge length) of all or one boundary element(s) in a mesh.
MinFEM.elementratio
— Functionelementratio(coords::Array{Array{Float64,1},1}) -> Float64
elementratio(mesh::Mesh) -> Array{Float64,1}
Returns ratio of inscribed to circumscribed circle or sphere for a specific element or all elements of a given mesh.
MinFEM.elementangle
— Functionelementangle(coords::Array{Array{Float64,1},1}) -> Float64
elementangle(mesh::Mesh) -> Array{Float64,1}
Returns smallest interior angle for a specific element or all elements of a given mesh.
MinFEM.outernormalvector
— Functionouternormalvector(mesh::Mesh; boundaryElements::Set{Int64}=Set{Int64}())
-> Vector{Float64}
outernormalvector(mesh::Mesh, boundaryElement::Int64)
-> Vector{Float64}
outernormalvector(mesh::Mesh, boundaryElement::Int64, J::AbstractMatrix{Float64})
-> Vector{Float64}
outernormalvector(dim::Int64, ii::Int64)
-> Vector{Float64}
Returns the outer normal vector at specified boundary element(s) of the boundary of the mesh or at the ii-th boundary of the dim-dimensional reference element.
MinFEM.stripwidth
— Functionstripwidth(mesh::Mesh) -> Float64
Return width L of a strip that the meshed domain fits into.
MinFEM.boundingbox
— Functionboundingbox(mesh::Mesh) -> Array{Float64}
Return two nodes, which span the bounding box of the mesh.
MinFEM.volume
— Functionvolume(mesh::Mesh) -> Float64
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) -> Array{Array{Float64,1},1}
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) -> Array{Array{Float64,1},1}
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) -> Array{Array{Float64,1},1}
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) -> Array{Array{Float64,1},1}
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)
-> Array{Array{Float64,1},1}
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) -> Array{Float64,1}
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) -> Array{Float64,1}
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 dimensions 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 dimensions 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 image dimension qdim. 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 image dimension qdim. 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) -> Float64
Returns ii-th local basis function evaluated at x. Dimension is determined by length(x).
MinFEM.grad_phi
— Functiongrad_phi(d::Int64, ii::Int64) -> 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 image dimension qdim.
MinFEM.assemble_laplacian
— Functionassemble_laplacian(
mesh::Mesh;
qdim::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
assemble_laplacian(
D::SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseMatrixCSC{Float64, Int64}
Returns the Laplacian stiffness matrix for all elements of mesh and image dimension qdim. Can either be assembled directly or from derivative matrix and weight vector.
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 image dimension qdim. Workaround based on the derivate tensor for the corresponding full element and the fact that the gradient ist constant on the element.
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 image dimension qdim. Workaround based on the derivate tensor for the corresponding full element and the fact that the gradient ist constant on the element.
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 image dimension qdim.
MinFEM.assemble_massmatrix
— Functionassemble_massmatrix(
mesh::Mesh;
qdim::Int64 = 1,
order::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
assemble_massmatrix(
E::SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all elements of mesh and image dimension qdim. Can either be assembled directly or from derivative matrix and weight vector.
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 image dimension qdim.
MinFEM.assemble_massmatrix_boundary
— Functionassemble_massmatrix_boundary(
mesh::Mesh;
boundaryElements::Set{Int64} = Set{Int64}(),
qdim::Int64 = 1,
order::Int64 = 1
) -> SparseMatrixCSC{Float64, Int64}
assemble_massmatrix_boundary(
E::SparseMatrixCSC{Float64, Int64},
w::AbstractVector{Float64}
) -> SparseMatrixCSC{Float64, Int64}
Returns the mass matrix with given local integration order for all or specified boundary elements of mesh and image dimension qdim. Can either be assembled directly or from derivative matrix and weight vector.
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
) -> 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,
DI::Set{Int64};
rhs = [],
bc = [],
qdim::Int64 = 1,
insert = 1.0)
assemble_dirichletcondition!(
A,
DI::Set{Boundary};
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. The value insert is put as diagonal element. Usually you want a 1.0 here.
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 function f or its FEM coefficient vector v on the given mesh with q being the conjugated exponent to p.
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.
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 on the given mesh with q being the conjugated exponent to p.
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)
First tries to set up 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)
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)
Recomputes the factorization of the stiffness matrix using assemble!()
.
Output
MinFEM.write_to_vtk
— Functionwrite_to_vtk(x::Vector{Vector{Float64}}, mesh::Mesh, dataName::Array{String},
fileName::String, qdim::Array{Int64})
write_to_vtk(x::Vector{Vector{Float64}}, mesh::Mesh, dataNames::Array{String},
fileName::String; qdim::Int64=1)
write_to_vtk(x::Vector{Float64}, mesh::Mesh, dataName::String,
fileName::String; qdim::Int64=1)
Writes one or multiple vectors with variable qdim corrsponding to the mesh nodes to a VTK file.
MinFEM.write_to_vtk_boundary
— Functionwrite_to_vtk_boundary(x::Vector{Vector{Float64}}, mesh::Mesh,
dataName::Array{String}, fileName::String,
qdim::Array{Int64}; boundary=Set{Boundary}())
write_to_vtk_boundary(x::Vector{Vector{Float64}}, mesh::Mesh,
dataNames::Array{String}, fileName::String;
boundary=Set{Boundary}(), qdim::Int64=1)
write_to_vtk_boundary(x::Vector{Float64}, mesh::Mesh,
dataName::String, fileName::String;
boundary=Set{Boundary}(), qdim::Int64=1)
Writes one or multiple vectors with variable qdim corrsponding to the mesh boundary elements to a VTK file.
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}())
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)
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::Any, data_name::String)
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, fileName::String; qdim::Int64=1)
Writes a coefficient vector x based on the nodes of mesh to the given file.
MinFEM.read_from_txt
— Functionread_from_txt(fileName::String)
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_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_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.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