General API

LowLevelFEM.CoordinateSystemType
CoordinateSystem(vec1, vec2)

A structure containing the data of a coordinate system.

  • vec1: direction of the new x axis.
  • vec2: together with vec1 determine the xy plane

If the problem is two dimensional, it is enough to give the first two elements of vec1. Elements of vec1 and vec2 can be functions. In 3D case the functions have three arguments (x, y, and z coordinates), otherwise (in 2D case) the number of arguments is two (x and y coordinates).

Types:

  • vec1: Vector{Float64}
  • vec2: Vector{Float64}

Examples

# 2D case
nx(x, y) = x
ny(x, y) = y
cs = CoordinateSystem([nx, ny])
Q = rotateNodes(problem, "body", cs)
q2 = Q' * q1 # where `q1` is in Cartesian, `q2` is in Axisymmetric coordinate system and
# `q1` is a nodal displacement vector.
S2 = Q' * S1 * Q # where `S1` is a stress field in Cartesian coordinate system while
# `S2` is in Axisymmetric coordinate system.

# 3D case
n1x(x, y, z) = x
n1y(x, y, z) = y
n2x(x, y, z) = -y
n2y = n1x
cs = CoordinateSystem([n1x, n1y, 0], [n2x, n2y, 0])
Q = rotateNodes(problem, "body", cs)
source
LowLevelFEM.EigenType
Eigen(f, ϕ, model)

A structure containing the eigenfrequencies and eigen modes.

  • f: eigenfrequencies
  • ϕ: eigen modes
  • model: same as Problem

Types:

  • f: Matrix{Float64}
  • ϕ: Vector{Float64}
  • model: Problem
source
LowLevelFEM.MaterialType
Material(phName, type, E, ν, ρ, k, c, α, λ, μ, κ)

Structure containing the material type and constants.

  • type: constitutive law (:Hooke, :StVenantKirchhoff, :NeoHookeCompressible)
  • E: elastic modulus
  • ν: Poisson's ratio
  • ρ: mass density
  • k: thermal conductivity
  • c: specific heat
  • α: thermal expansion coefficient
  • λ: Lamé parameter
  • μ: Lamé parameter
  • κ: bulk modulus

phName is the name of the physical group where the material is used.

Types:

  • phName: String
  • type: Symbol
  • E: Float64
  • ν: Float64
  • ρ: Float64
  • k: Float64
  • c: Float64
  • α: Float64
  • λ: Float64
  • μ: Float64
  • κ: Float64
source
LowLevelFEM.ProblemType
Problem(materials; thickness=..., type=..., bandwidth=...)

Structure containing key data for a problem.

  • Parts of the model with their material constants. Multiple materials can be provided (see material).
  • Problem type: :Solid, :PlaneStrain, :PlaneStress, :AxiSymmetric, :HeatConduction, :PlaneHeatConduction,

:AxiSymmetricHeatConduction, :Truss. For :AxiSymmetric, the symmetry axis is y, and the geometry must be drawn in the positive x half-plane.

  • Bandwidth optimization using Gmsh's built-in reordering. Options: :RCMK, :Hilbert, :Metis, or :none (default).
  • Dimension of the problem, determined by type.
  • Material constants (vector of Material; see the Material struct).
  • Plate thickness (for 2D plate problems).
  • Number of nodes (non).
  • Geometry dimension.
  • Problem dimension (e.g., a 3D heat conduction problem is a 1D problem).
  • In case of 2D truss displacements have to be fixed in the third direction.

Types:

  • materials: Material
  • type: Symbol
  • bandwidth: Symbol
  • dim: Integer
  • thickness: Float64
  • non: Integer
  • dim: Integer
  • pdim: Integer
source
LowLevelFEM.ScalarFieldType
ScalarField(A, a, t, numElem, nsteps, type, model)
ScalarField(problem::Problem, dataField::Vector; steps=1, tmin=0.0, tmax=tmin+(steps-1))
ScalarField(problem::Problem, phName::String, data::Union{Number,Function};
            steps=1, tmin=0.0, tmax=tmin+(steps-1))
ScalarField(s::ScalarField; steps=1, tmin=0.0, tmax=tmin+(steps-1), step=1)

Container for a time-dependent scalar field defined on a finite element mesh (e.g. temperature, pressure, potential).

A ScalarField can store the field either

  • element-wise (values at the element nodes, stored in A), or
  • nodally (values at global mesh nodes, stored in a).

Time dependence is handled by storing multiple time steps for each spatial degree of freedom.


Stored data

  • A : Vector of matrices holding element-wise values (A[e][k,i] = value at local node k of element e at time step i)
  • a : Matrix of nodal values (a[n,i] = value at node n at time step i)
  • t : Vector of time instants corresponding to the stored time steps
  • numElem : Vector of element tags associated with A
  • nsteps : Number of stored time steps
  • type : Symbol identifying the physical meaning of the field
  • model : Associated Problem

At a given time, either A or a is typically populated, depending on whether the field is element-wise or nodal.


Constructors

Low-level constructor

ScalarField(A, a, t, numElem, nsteps, type, model)

Directly constructs a ScalarField from preallocated data arrays.


From spatial field definitions

ScalarField(problem, dataField; steps, tmin, tmax)

Constructs an element-wise scalar field from a vector of field definitions (e.g. as returned by field(...)). The scalar value may be:

  • a constant,
  • a spatial function f(x,y,z),
  • or a space–time function f(x,y,z,t).

Values are evaluated at the element nodes for each requested time step.


From a physical group

ScalarField(problem, phName, data; steps, tmin, tmax)

Convenience constructor for defining a scalar field on a single physical group. Equivalent to calling field(phName, f=data) internally.


From an existing ScalarField

ScalarField(s; steps, tmin, tmax, step)

Creates a new ScalarField by replicating a selected time step of an existing field. This is useful for initializing a time-dependent field from a static solution.

  • step selects the time index of s to be replicated.
  • The new field contains steps identical time slices.

Notes

  • Time steps do not need to be uniformly spaced.
  • No assumptions are made about governing equations; this is a pure data container.
  • Spatial interpolation and projections (e.g. nodal ↔ element-wise) are handled by separate utility functions.

Example

s(x,y,z) = 2x + 3y
fs = field("body", f=s)

S1 = ScalarField(problem, [fs])
S2 = ScalarField(problem, "body", s)

Here S1 and S2 define equivalent element-wise scalar fields.

source
LowLevelFEM.SystemMatrixType
SystemMatrix(A::SparseMatrixCSC{Float64}, model::Problem)

Structure containing the stiffness/mass/heat conduction/heat capacity/latent heat/... matrix and the associated Problem.

Types:

  • A: SparseMatrixCSC{Float64}
  • model: Problem
source
LowLevelFEM.TensorFieldType
TensorField(A, a, t, numElem, nsteps, type, model)
TensorField(problem::Problem, dataField::Vector)
TensorField(problem::Problem, phName::String, data::Matrix)
TensorField(comps::Matrix{ScalarField})

Container for a (possibly time-dependent) second-order tensor field defined on a finite element mesh (e.g. stress, strain, conductivity tensor).

A TensorField stores a full 3×3 tensor at each spatial location, either

  • element-wise, with values given at element nodes (A), or
  • nodally, with values given at global mesh nodes (a).

Time dependence is handled by storing multiple time steps for each tensor component.


Stored data

  • A : Vector of matrices holding element-wise tensor values (A[e][9k-8:9k, i] = tensor components at local node k of element e at time step i)
  • a : Matrix of nodal tensor values (same component ordering as element-wise storage)
  • t : Vector of time instants corresponding to the stored time steps
  • numElem : Vector of element tags associated with A
  • nsteps : Number of stored time steps
  • type : Symbol identifying the physical meaning of the tensor field (e.g. :e, :stress, :tensor)
  • model : Associated Problem

Tensor components are stored in row-major, interleaved form for each node:


(xx, yx, zx,
xy, yy, zy,
xz, yz, zz)

repeated for all nodes and time steps.


Constructors

Low-level constructor

TensorField(A, a, t, numElem, nsteps, type, model)

Directly constructs a TensorField from preallocated data arrays.


From element-wise field definitions

TensorField(problem, dataField::Vector)

Constructs an element-wise tensor field from a vector of field definitions (e.g. as returned by field(...)).

Each tensor component may be specified as:

  • a constant,
  • a spatial function f(x,y,z).

Values are evaluated at the element nodes.


From a constant or functional tensor

TensorField(problem, phName, data::Matrix)

Convenience constructor for defining a tensor field on a single physical group, where data is a 3×3 matrix whose entries are constants or functions f(x,y,z).


From scalar components

TensorField(comps::Matrix{ScalarField})

Assembles a 3×3 tensor field from a 3×3 matrix of compatible ScalarFields.

Requirements:

  • exactly a 3×3 matrix of scalar fields,
  • identical Problem,
  • identical element numbering,
  • identical time discretization.

Each scalar field provides one tensor component.


Notes

  • Time steps do not need to be uniformly spaced.
  • No symmetry is assumed; all 9 tensor components are stored explicitly.
  • This is a pure data container; no constitutive or kinematic assumptions are implied.
  • Spatial operations (e.g. divergence, invariants, projections) are handled by separate utility functions.
  • Element-wise storage is the primary representation; nodal storage may be empty depending on construction.

Example

tx(x,y,z)  = x + y
txy(x,y,z) = z

ft = field("body", fx=tx, fxy=txy, fz=3)
T1 = TensorField(problem, [ft])

T2 = TensorField(problem, "body", [1 0 0;
                                  0 2 0;
                                  0 0 3])

Here T1 and T2 define equivalent element-wise tensor fields.

source
LowLevelFEM.TransformationType
Transformation(T::SparseMatrixCSC{Float64}, non::Int64, dim::Int64)

Structure containing the transformation matrix T at each node in the FE mesh, the number of nodes non, and the problem dimension dim.

Types:

  • T: SparseMatrixCSC{Float64}
  • non: Int64
  • dim: Int64
source
LowLevelFEM.VectorFieldType
VectorField(A, a, t, numElem, nsteps, type, model)
VectorField(problem::Problem, dataField::Vector)
VectorField(problem::Problem, phName::String, data::Vector)
VectorField(problem::Problem, phName::String, func::Function)
VectorField(comps::Vector{ScalarField})

Container for a (possibly time-dependent) vector field defined on a finite element mesh (e.g. displacement, velocity, heat flux).

A VectorField stores a 3D vector field either

  • element-wise, with values given at element nodes (A), or
  • nodally, with values given at global mesh nodes (a).

Time dependence is handled by storing multiple time steps for each spatial degree of freedom, analogously to ScalarField.


Stored data

  • A : Vector of matrices holding element-wise vector values (A[e][3k-2:3k, i] = vector value at local node k of element e at time step i)
  • a : Matrix of nodal vector values (layout follows the same interleaved component ordering)
  • t : Vector of time instants corresponding to the stored time steps
  • numElem : Vector of element tags associated with A
  • nsteps : Number of stored time steps
  • type : Symbol identifying the physical meaning of the vector field (e.g. :v3D)
  • model : Associated Problem

Vector components are stored in interleaved form (x₁,y₁,z₁,x₂,y₂,z₂,…) for each element or node.


Constructors

Low-level constructor

VectorField(A, a, t, numElem, nsteps, type, model)

Directly constructs a VectorField from preallocated data arrays.


From element-wise field definitions

VectorField(problem, dataField::Vector)

Constructs an element-wise vector field from a vector of field definitions (e.g. as returned by field(...)). Each component may be specified as:

  • a constant,
  • a spatial function f(x,y,z).

Values are evaluated at the element nodes.


From a physical group and component data

VectorField(problem, phName, data::Vector)

Convenience constructor for defining a vector field on a single physical group, where data = [fx, fy, fz] specifies the three components (constants or functions).


From a vector-valued function

VectorField(problem, phName, func::Function)

Constructs an element-wise vector field by evaluating a function func(x,y,z) -> (vx,vy,vz) at the element nodes.


From scalar components

VectorField(comps::Vector{ScalarField})

Assembles a 3D vector field from exactly three compatible ScalarFields.

Requirements:

  • exactly three components,
  • same Problem,
  • identical element numbering,
  • identical time discretization.

Each scalar field provides one vector component.


Notes

  • Time steps do not need to be uniformly spaced.
  • No governing equations are implied; this is a pure data container.
  • Spatial operations (gradient, divergence, projections, etc.) are handled by separate utility functions.
  • Element-wise storage is the primary representation; nodal storage may be empty depending on construction.

Example

vx(x,y,z) = x + y
vy(x,y,z) = z

fv = field("body", fx=vx, fy=vy, fz=3)
V1 = VectorField(problem, [fv])

V2 = VectorField(problem, "body", [vx, vy, 3])

Here V1 and V2 define equivalent element-wise vector fields.

source
Base.copyMethod
Base.copy(A::SystemMatrix)

Internal function to copy the whole content of a SystemMatrix.

source
Base.getindexMethod
T[:,j] -> VectorField

Extract the j-th column of a 3×3 TensorField as a 3-component VectorField.

Equivalent to: [T[1,j], T[2,j], T[3,j]]

Arguments

  • j::Int: column index (1 ≤ j ≤ 3).

Returns

A VectorField of size 3.

Notes

Colon indexing is required; T[j] is not allowed.

source
Base.getindexMethod
T[i,:] -> VectorField

Extract the i-th row of a 3×3 TensorField as a 3-component VectorField.

Equivalent to: [T[i,1], T[i,2], T[i,3]]

Arguments

  • i::Int: row index (1 ≤ i ≤ 3).

Returns

A VectorField of size 3.

Notes

Colon indexing is required; T[i] is not allowed.

source
Base.getindexMethod
T[i,j] -> ScalarField

Extract a single tensor component (i,j) from a 3×3 TensorField.

Arguments

  • i::Int, j::Int: tensor indices in 1:3.

Returns

A ScalarField corresponding to component Tᵢⱼ.

Notes

Block extraction is done elementwise or nodally depending on the representation of T.

Errors

Raises an error if i or j is outside 1:3.

source
Base.getindexMethod

Invalid index pattern for TensorField.

Allowed:

  • T[i,j] → ScalarField
  • T[i,:] → VectorField (row)
  • T[:,j] → VectorField (column)

Not allowed:

  • T[i]
  • T[:]
  • T[1:2,1]
  • T[:, :]
  • T[[1,3],2]
source
Base.getindexMethod
v[k] -> ScalarField

Extract the k-th component (1,2,3) of a VectorField as a ScalarField.

Arguments

  • k::Int: component index (1 → x, 2 → y, 3 → z).

Returns

A ScalarField in the same representation as the parent field (elementwise or nodal).

Errors

Raises an error for any non-scalar indexing:

  • v[1:2] → ERROR
  • v[[1,3]] → ERROR
  • v[:] → ERROR

Only single-component extraction v[k] is allowed.

source
Base.setindex!Method
T[i, j] = s
T[i, :] = v
T[:, j] = v

Assign a scalar or vector component into a TensorField.

Supported forms

1. Assign scalar component

T[i, j] = s

  • Inserts the scalar field s into the tensor component (i, j) (with i,j ∈ 1:3).
  • Only the rows of the (i,j) block are overwritten.

2. Assign row vector

T[i, :] = v

  • v must be a VectorField with 3 components.
  • Replaces the entire i-th tensor row with the components of v.

3. Assign column vector

T[:, j] = v

  • v must be a VectorField.
  • Replaces the entire j-th tensor column with the components of v.

Requirements

  • Assigned field(s) must belong to the same problem as T.
  • The number of time steps must match.
  • Storage mode must match:
    • Elementwise tensor → only elementwise scalar/vector may be assigned.
    • Nodal tensor → only nodal scalar/vector may be assigned.
  • Element counts or nodal sizes must be identical.

Examples

# scalar insertion
T[1,1] = s1
T[2,3] = s2

# row insertion
T[1,:] = v1    # v1 is VectorField([s11, s12, s13])

# column insertion
T[:,3] = v2

Errors

Errors are thrown when:

  • indices are out of bounds
  • shapes do not match (nodal vs elementwise mismatch)
  • component sizes differ
  • attempting unsupported forms like T[1] = ... or T[:, :] = ...
source
Base.setindex!Method
v[k] = s

Assign a scalar component into a VectorField.

This operation overwrites the k-th component (k = 1, 2, 3) of the vector field.

Requirements

  • v must be a VectorField
  • s must be a ScalarField
  • both fields must belong to the same problem
  • both must have the same number of time steps
  • storage layout must match:
    • If v is elementwise (v.A ≠ []), then s must also be elementwise and must have the same number of elements (numElem).
    • If v is nodal (v.a ≠ [;;]), then s must also be nodal and must have the same nodal size.

Only the rows belonging to component k are modified. Other components remain unchanged.

Example

v = VectorField([s1, s2, s3])

v[1] = s4       # overwrite first component
v[3] = 2s1      # scale and assign into third component

Errors

An informative error is thrown if:

  • k ∉ 1:3
  • the fields belong to different problems
  • storage layouts do not match (nodal vs. elementwise)
  • nodal or element counts differ
source
Base.showMethod
Base.show(io::IO, M::SystemMatrix)

Internal function to display SystemMatrix as a SparseMatrixCSC{Float64}.

source
Base.showMethod
Base.show(io::IO, M::Union{ScalarField,VectorField,TensorField})

Internal function to display ScalarField, VectorField and TensorField as a Matrix{Float64} or Vector{Matrix{Float64}}.

source
LowLevelFEM.constrainedDoFsMethod
constrainedDoFs(problem, supports)

Returns the serial numbers of constrained degrees of freedom. Support is a vector of boundary conditions given with the function displacementConstraint.

Return: DoFs

Types:

  • problem: Problem
  • supports: Vector{Tuple{String, Float64, Float64, Float64}}
  • DoFs: Vector{Int64}
source
LowLevelFEM.displacementConstraintMethod
displacementConstraint(name; ux=..., uy=..., uz=...)

Specifies displacement constraints on the physical group name. At least one of ux, uy, or uz must be provided (depending on the problem dimension). ux, uy, or uz can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); displacementConstraint("support1", ux=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Examples

hc = heatConvection("outer", h=10.0, Tₐ=20.0)

Examples

src = heatSource("body", h=1.0e6)

Examples

q = heatFlux("out", qn=500.0)

Examples

bcT = temperatureConstraint("hot_face", T=100.0)

Examples

ld = load("load", fy=-1.0)

Examples

bc = displacementConstraint("supp", ux=0, uy=0)

Types:

  • name: String
  • ux: Float64 or Function
  • uy: Float64 or Function
  • uz: Float64 or Function
source
LowLevelFEM.elasticSupportMethod
elasticSupport(name; kx=..., ky=..., kz=...)

Specifies distributed stiffness for an elastic support on the physical group name. kx, ky, or kz can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); elasticSupport("supp1", kx=fn)) Default values are 1.

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • kx: Float64 or Function
  • ky: Float64 or Function
  • kz: Float64 or Function
source
LowLevelFEM.elementsToNodesMethod
elementsToNodes(T)

Solves the nodal results F from the elemental results T. T can be ScalarField, VectorField or TensorField.

Return: F

Types:

  • T: ScalarField, VectorField or TensorField
  • F: ScalarField, VectorField or TensorField
source
LowLevelFEM.expandTo3DMethod
expandTo3D(v2D::VectorField)

Expand a 2D vector field into 3D by adding a zero z-component.

return: VectorField

Examples

V3D = expandTo3D(V2D)
source
LowLevelFEM.fieldMethod
field(name; f=..., fx=..., fy=..., fz=..., fxy=..., fyz=..., fzx=...)

Specifies the value of a scalar, vector, or tensor field on the physical group name. At least one of fx, fy, or fz etc. must be provided (depending on the problem dimension). Each component can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); field("surf1", fx=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function, ... x7}

Types:

  • name: String
  • f: Float64 or Function
  • fx: Float64 or Function
  • fy: Float64 or Function
  • fz: Float64 or Function
  • fxy: Float64 or Function
  • fyz: Float64 or Function
  • fzx: Float64 or Function

Examples

f1(x, y, z) = y
f2 = field("face1", f=f1)
qq = scalarField(problem, [f2])
qqq = showDoFResults(qq, :scalar)
source
LowLevelFEM.fieldErrorMethod
fieldError(F)

Computes the deviation of field results F (stresses, strains, heat flux components) at nodes where the field has jumps. The result can be displayed with the showDoFResults function.

Returns: e

Types:

  • F: VectorField or TensorField
  • e: ScalarField
source
LowLevelFEM.freeDoFsMethod
freeDoFs(problem, supports)

Returns the serial numbers of unconstrained degrees of freedom. Support is a vector of boundary conditions given with the function displacementConstraint.

Return: DoFs

Types:

  • problem: Problem
  • supports: Vector{Tuple{String, Float64, Float64, Float64}}
  • DoFs: Vector{Int64}
source
LowLevelFEM.getEigenValuesMethod
getEigenValues(A::VectorField)

A function to extract the elements of a vector field to separate scalar fields.

Return: λ1, λ2, λ3

Types:

  • A: VectorField
  • λ1: ScalarField
  • λ2: ScalarField
  • λ3: ScalarField

Examples:

using LinearAlgebra
λ, Q = eigen(S)
λ1, λ2, λ2 = getEigenValues(λ)
source
LowLevelFEM.getEigenVectorsMethod
getEigenVectors(A::TensorField)

A function to extract the columns of a tensor field to separate vector fields.

Return: N1, N2, N3

Types:

  • A: TensorField
  • N1: VectorField
  • N2: VectorField
  • N3: VectorField

Examples:

using LinearAlgebra
λ, Q = eigen(S)
N1, N2, N2 = getEigenVectors(Q)
source
LowLevelFEM.heatConvectionMethod
heatConvection(name; h=..., Tₐ=...)

Specifies convective boundary conditions on the surface in the physical group name. h is the heat transfer coefficient of the surrounding medium; Tₐ is the ambient temperature. Tₐ can be either a constant or a function of x, y, and z.

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • h: Float64
  • Tₐ: Float64 or Function
source
LowLevelFEM.heatFluxMethod
heatFlux(name; qn=...)

Specifies the heat flux normal to the surface of the physical group name. qn can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); load("flux1", qn=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • qn: Float64 or Function
source
LowLevelFEM.heatSourceMethod
heatSource(name; h=...)

Specifies the volumetric heat source in the physical group name. h can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); load("source1", h=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • h: Float64 or Function
source
LowLevelFEM.integrateMethod
integrate(problem::Problem, phName::String, f::Union{Function,ScalarField}; step::Int64=1, order::Int64=...)
∫(problem::Problem, phName::String, f::Union{Function,ScalarField}; step::Int64=1, order::Int64=...)

Integrates the function or scalar field f over the physical group phName defined in the geometry of problem. If f is a ScalarField, the time step step will be integrated. The optional parameter order controls the numerical integration rule. If order > 0, it is used as a hint for selecting the Gauss quadrature order (i.e. the number of integration points) employed during integration. The exactness of the integration depends on the element geometry and the regularity of the integrand.

Returns: integral

Types:

  • problem: Problem
  • phName: String
  • f: Function (of x, y and z) or ScalarField
  • integral: Number

Examples:

f(x, y, z) = x^2 + y^2
Iz = integrate(prob, "body", f)
source
LowLevelFEM.integrateMethod
integrate(s::Union{ScalarField,VectorField,TensorField})
∫(s::Union{ScalarField,VectorField,TensorField})

Compute the time integral of a time-dependent field using a discrete inverse of the central finite-difference time derivative.

The result reproduces the original field (up to a time-independent constant) when applied to a field obtained by ∂t.

source
LowLevelFEM.isElementwiseMethod

isElementwise(field)

Check if a given field is defined per element (elementwise quantity).

Elementwise quantities are associated with finite elements as a whole, for example stresses, strains, or energy densities evaluated inside elements.

Examples

isElementwise(displacement_field)   # returns false
isElementwise(strain_field)         # returns true
source
LowLevelFEM.isNodalMethod
isNodal(field)

Check if a given field is defined at nodes (nodal quantity).

Nodal quantities are associated with mesh nodes, for example displacements, nodal forces, or nodal temperatures.

Examples

isNodal(displacement_field)   # returns true
isNodal(strain_field)         # returns false
source
LowLevelFEM.isSavedMethod
isSaved(fileName::String)

Checks whether a variable has been saved or not.

Returns: Boolean

Types:

  • fileName: String
source
LowLevelFEM.loadMethod
load(name; fx=..., fy=..., fz=...)

Specifies a distributed load on the physical group name. At least one of fx, fy, or fz must be provided (depending on the problem dimension). fx, fy, or fz can be a constant or a function of x, y, and z or ScalarField. (e.g., fn(x,y,z) = 5*(5-x); load("load1", fx=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • fx: Float64 or Function
  • fy: Float64 or Function
  • fz: Float64 or Function
source
LowLevelFEM.loadFieldMethod
loadField(fileName::String)

Loads a ScalarField, VectorField, or TensorField from a file named fileName (without "-LLF-Data.jld2").

Returns: variable

Types:

  • fileName: String
  • variable: ScalarField, VectorField or TensorField
source
LowLevelFEM.materialMethod
material(name; type=:Hooke, E=2.0e5, ν=0.3, ρ=7.85e-9, k=45, c=4.2e8, α=1.2e-5, λ=νE/(1+ν)/(1-2ν), μ=E/(1+ν)/2, κ=E/(1-2ν)/3)

Returns a structure in which name is the name of a physical group, type is the name of the constitutive law (e.g., :Hooke), E is the modulus of elasticity, ν is Poisson's ratio, and ρ is the mass density. k is the thermal conductivity, c is the specific heat, α is the coefficient of thermal expansion, λ and μ are the Lamé parameters, and κ is the bulk modulus.

Returns: mat

Examples

mat = material("body", E=210e3, ν=0.3, ρ=7.85e-9)

Types:

  • mat: Material
  • name: String
  • type: Symbol
  • E: Float64
  • ν: Float64
  • ρ: Float64
  • k: Float64
  • c: Float64
  • α: Float64
  • λ: Float64
  • μ: Float64
  • κ: Float64
source
LowLevelFEM.nodesToElementsMethod
nodesToElements(T, onPhysicalGroup="")

Solves the element results F from the nodal results T. T can be ScalarField, VectorField or TensorField. If onPhysicalGroup is an existing physical group in gmsh, field T will be solve only on elements belonging to that physical group. Dimension of physical group can be different than the dimension of the problem. (eg. mapping from 3D volume to a 2D surface)

Return: F

Types:

  • T: ScalarField, VectorField or TensorField
  • F: ScalarField, VectorField or TensorField
  • onPhysicalGroup: String
source
LowLevelFEM.normalVectorMethod
normalVector(problem::Problem, phName::String) -> VectorField

Compute outward unit normal vectors for all nodes belonging to a surface-type physical group in a 3D Gmsh model.

For curve-type physical groups, the function returns normal curvature vectors: the direction corresponds to the unit normal within the curve’s osculating plane, and the vector magnitude equals the local curvature of the curve.

Arguments

  • problem::Problem: A Problem structure containing the current Gmsh model (name, dimension, number of nodes, etc.).
  • phName::String: The name of a physical surface group in Gmsh for which the normal vectors are computed.

Description

The function sets the current model, queries the elements and nodes that belong to the given physical surface group, and evaluates the gradients of the Lagrange basis functions to compute local tangent vectors of the surface. Normal vectors are obtained as cross products of these tangents and normalized to unit length.

Each node belonging to the physical surface group is assigned its corresponding 3D unit normal vector.

Returns

  • VectorField: A VectorField structure that stores the nodal normal vectors on the given physical surface.

Errors

  • Throws an error if the provided physical group is not of surface type (dim != 2).

Example

using LowLevelFEM

# Load a 3D geometry and mesh it in Gmsh
gmsh.initialize()
gmsh.open("box_with_surface.msh")

# Define a 3D model on the volume physical group "body"
mat = material("body")
prob = Problem([mat])

# Compute nodal normals on a physical surface named "leftWall"
nv = normalVector(problem, "leftWall")

# Show the normal vectors on the model
showDoFResults(nv)
openPostProcessor()
source
LowLevelFEM.openPostProcessorMethod
openPostProcessor(; model=...)

Launches the Gmsh postprocessor window with the postprocessor tree opened (of model).

Returns: nothing

Types:

  • model: Int64
source
LowLevelFEM.openPreProcessorMethod
openPreProcessor(; openGL=...)

Launches the Gmsh preprocessor window with OpenGL disabled by default.

Returns: nothing

Types:

  • openGL: Boolean
source
LowLevelFEM.plotOnPathMethod
plotOnPath(pathName, field; points=100, step=..., plot=..., name=..., visible=..., offsetX=..., offsetY=..., offsetZ=...)

Loads a 2D plot along a path into a View in Gmsh. field is the View id in Gmsh from which the field data is imported. pathName is the name of a physical group that contains a curve. The curve is divided into equal-length segments with points sampling points. The field is shown at these points. step is the sequence number of the displayed step. If no step is given, it shows all available steps as an animation. If plot is true, an additional return parameter (a tuple of vectors) is returned, where x is the horizontal axis and y is the vertical axis of the plot (see the Plots package). name is the title of the graph, and visible is a Boolean to toggle the initial visibility in Gmsh on or off. This function returns the tag of the View.

Returns: tag

or

Returns: tag, xy

Types:

  • problem: Problem
  • pathName: String
  • field: Integer
  • points: Integer
  • step: Integer
  • plot: Boolean
  • name: String
  • visible: Boolean
  • tag: Integer
  • xy: Tuples{Vector{Float64},Vector{Float64}}
source
LowLevelFEM.probeMethod
probe(A::Union{ScalarField,VectorField,TensorField}, x::Number, y::Number, z::Number; step=Int)

Get the value of the field A at point coordinates x, y, z at time step step.

Returns: Float64 or Vector{Float64} or Matrix{Float64}

Types:

  • A: ScalarField or VectorField or TensorField
  • x: Number
  • y: Number
  • z: Number
  • step: Int
source
LowLevelFEM.probeMethod
probe(A::Union{ScalarField,VectorField,TensorField}, s::String; step=Int)

Get the value of the field A at a point given by its physical name in Gmsh at time step step.

Returns: Float64 or Vector{Float64} or Matrix{Float64}

Types:

  • A: ScalarField or VectorField or TensorField
  • x: Number
  • y: Number
  • z: Number
  • step: Int
source
LowLevelFEM.projectTo2DMethod
projectTo2D(v3D::VectorField)

Project a 3D vector field onto the xy-plane by discarding the z-component.

return: VectorField

Examples

V2D = expandTo3D(V3D)
source
LowLevelFEM.resultantMethod
resultant(field, phName)

Computes the resultant of vector field field on the physical group phName. Returns the resultant(s) in a tuple. The number of elements in the tuple depends on the dimension of problem (dimension of field). It can solve for example the resultant of a load vector (sum of the elements of the vector).

Return: resx

or

Return: resx, resy

or

Return: resx, resy, resz

Types:

  • field: VectorField
  • phName: String
  • resx: Float64
  • resy: Float64
  • resz: Float64
source
LowLevelFEM.rotateNodesMethod
rotateNodes(problem, phName, CoordSys)

Creates the T transformation matrix, which rotates the nodal coordinate system of the nodes in phName physical group to the coordinate systen defined by CoordSys. The mesh belongs to problem.

Return: T

Types:

  • problem: Problem
  • phName: String
  • CoordSys: CoordinateSystem
  • T: SparseMatrix
source
LowLevelFEM.saveFieldMethod
saveField(fileName::String, variable::Union{ScalarField,VectorField,TensorField,Number})

Saves variable of type ScalarField, VectorField, or TensorField to a file named fileName. The name of the file will be complemented with the string "-LLF-Data.jld2"

Returns: nothing

Types:

  • fileName: String
  • variable: ScalarField, VectorField or TensorField
source
LowLevelFEM.scalarFieldMethod
scalarField(problem, dataField)

Defines a scalar field from dataField, which is a tuple of name of physical group and prescribed values or functions. Mesh details are in problem.

Return: Vector{Float64}

Types:

  • problem: Problem
  • dataField: Vector{Tuple{String, Float64,...}}

Examples

f2 = field("face1", f=1)
qq = scalarField(problem, [f2])
qqq = showDoFResults(qq)

qq2 = scalarField(problem, "body", 2.0)

Here ScalarField is defined in nodes.

source
LowLevelFEM.setParameterMethod
setParameter(name, value)

Defines a parameter name and sets its value to value.

Returns: nothing

Types:

  • name: String
  • value: Float64
source
LowLevelFEM.setParametersMethod
setParameters(name, value)

Defines a parameter name and sets its value to value, which is a Vector{Float64}.

Returns: nothing

Types:

  • name: String
  • value: Vector{Float64}
source
LowLevelFEM.showBucklingResultsMethod
showBucklingResults(Φ, name=..., visible=...)

Loads buckling results into a View in Gmsh. Φ is an Eigen struct. name is a title to display and visible is a Boolean to toggle the initial visibility in Gmsh on or off. Click on ▷| to change the results. This function returns the tag of the View.

Returns: tag

Types:

  • Φ: Eigen
  • name: String
  • visible: Boolean
  • tag: Integer
source
LowLevelFEM.showDoFResultsMethod
showDoFResults(q, comp; name=..., visible=..., factor=0)

Loads nodal results into a View in Gmsh. q is the field to show, comp is the component of the field (:vector, :uvec, :ux, :uy, :uz, :vvec, :vx, :vy, :vz, :qvec, :qx, :qy, :qz, :T, :p, :qn, :s, :sx, :sy, :sz, :sxy, :syx, :syz, :szy, :szx, :sxz, :e, :ex, :ey, :ez, :exy, :eyx, :eyz, :ezy, :ezx, :exz, :seqv, :scalar, :tensor), name is a title to display and visible is a Boolean to toggle the initial visibility in Gmsh on or off. If q has more columns, then a sequence of results will be shown (e.g., as an animation). factor multiplies the DoF result to increase for better visibility. This function returns the tag of the View.

Returns: tag

Types:

  • q: ScalarField, VectorField or TensorField
  • comp: Symbol
  • name: String
  • visible: Boolean
  • tag: Integer
source
LowLevelFEM.showElementResultsMethod
showElementResults(F, comp; name=..., visible=..., smooth=...)

Same as ShowStressResults or showStrainResults, depending on the type of F data field.

Return: tag

Types:

  • F: TensorField
  • comp: Symbol
  • name: String
  • visible: Boolean
  • smooth: Boolean
  • tag: Integer
source
LowLevelFEM.showHeatFluxResultsMethod
showHeatFluxResults(Q, comp; name=..., visible=..., smooth=...)

Loads heat flux results into a View in Gmsh. Q is a heat flux field to show, comp is the component of the field (:qvec, :qx, :qy, :qz, :q), name is a title to display, visible is a Boolean to toggle the initial visibility in Gmsh on or off, and smooth is a Boolean to toggle smoothing on or off. If Q contains more than one time step, a sequence of results will be shown (e.g., as an animation). This function returns the tag of the View.

Returns: tag

Types:

  • S: VectorField
  • comp: Symbol
  • name: String
  • visible: Boolean
  • smooth: Boolean
  • tag: Integer
source
LowLevelFEM.showModalResultsMethod
showModalResults(Φ, name=..., visible=...)

Loads modal results into a View in Gmsh. Φ is an Eigen struct. name is a title to display and visible is a Boolean to toggle the initial visibility in Gmsh on or off. Click on ▷| to change the results. This function returns the tag of the View.

Returns: tag

Types:

  • Φ: Eigen
  • name: String
  • visible: Boolean
  • tag: Integer
source
LowLevelFEM.showOnSurfaceMethod
showOnSurface(field, phName; grad=false, component=:x, offsetX=0, offsetY=0, offsetZ=0, name=phName, visible=false)

Shows the values of a scalar field on a surface with physical name phName. field is the tag of a View in Gmsh. The values of the field are calculated at the intersection with the surface. grad is a Boolean to toggle the gradient of the field on or off. component is the component of the gradient of field (:x, :y, :z) to be shown. offsetX, offsetY, offsetZ are the offsets in the x, y, and z directions where the values are sampled. name is a title to display, and visible is a Boolean to toggle the initial visibility in Gmsh on or off.

Returns: tag

Types:

  • field: Integer
  • phName: String
  • grad: Boolean
  • component: Symbol
  • offsetX: Float64
  • offsetY: Float64
  • offsetZ: Float64
  • name: String
  • visible: Boolean
  • tag: Integer
source
LowLevelFEM.showStrainResultsMethod
showStrainResults(E, comp; name=..., visible=..., smooth=...)

Loads strain results into a View in Gmsh. E is a strain field to show, comp is the component of the field (:e, :ex, :ey, :ez, :exy, :eyz, :ezx), name is a title to display, visible is a Boolean to toggle the initial visibility in Gmsh on or off and smooth is a Boolean to toggle smoothing the stress field on or off. If E contains more than one time steps, then a sequence of results will be shown (e.g., as an animation). This function returns the tag of the View.

Returns: tag

Types:

  • E: TensorField
  • comp: Symbol
  • name: String
  • visible: Boolean
  • smooth: Boolean
  • tag: Integer
source
LowLevelFEM.showStressResultsMethod
showStressResults(S, comp; name=..., visible=..., smooth=...)

Loads stress results into a View in Gmsh. S is a stress field to show, comp is the component of the field (:s, :sx, :sy, :sz, :sxy, :syz, :szx, :seqv), name is a title to display, visible is a Boolean to toggle the initial visibility in Gmsh on or off, and smooth is a Boolean to toggle smoothing the stress field on or off. If S contains more than one time step, a sequence of results will be shown (e.g., as an animation). This function returns the tag of the View.

Returns: tag

Types:

  • S: TensorField
  • comp: Symbol
  • name: String
  • visible: Boolean
  • smooth: Boolean
  • tag: Integer
source
LowLevelFEM.tangentVectorMethod
tangentVector(problem::Problem, phName::String) -> VectorField

Compute unit tangent vectors for all nodes of a curve-type physical group in a 3D Gmsh model.

Arguments

  • problem::Problem: A Problem structure containing the current Gmsh model (name, dimension, number of nodes, etc.).
  • phName::String: The name of a physical curve group in Gmsh for which the tangent vectors are computed.

Description

The function sets the current model, queries the elements and nodes that belong to the given 1D physical group, and evaluates the gradients of the Lagrange basis functions to determine the local mapping derivative ∂x/∂ξ. Since this derivative represents the geometric tangent direction of the curve at each evaluation point, it is normalized to produce a unit tangent vector.

Each node belonging to the physical curve group is assigned the corresponding 3D unit tangent vector aligned with the parametric direction of the curve.

Returns

  • VectorField: A VectorField structure that stores the nodal tangent vectors on the given physical curve.

Errors

  • Throws an error if the provided physical group is not of curve type (dim != 1).

Example

using LowLevelFEM

# Load a 3D geometry containing a curved edge
gmsh.initialize()
gmsh.open("curve_geometry.msh")

# Create a problem structure
mat = material("body")
prob = Problem([mat])

# Compute tangent vectors along the curve named "leadingEdge"
tv = tangentVector(prob, "leadingEdge")

# Visualize the tangent field
showDoFResults(tv)
openPostProcessor()
source
LowLevelFEM.temperatureConstraintMethod
temperatureConstraint(name; T=...)

Specifies temperature constraints on the physical group name. T can be a constant or a function of x, y, and z. (e.g., fn(x,y,z) = 5*(5-x); temperatureConstraint("surf1", T=fn))

Returns: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}

Types:

  • name: String
  • T: Float64 or Function
source
LowLevelFEM.tensorFieldMethod
tensorField(problem, dataField; type=...)

Defines a vector field from dataField, which is a tuple of name of physical group and prescribed values or functions. Mesh details are in problem. type can be an arbitrary Symbol, e.g., :u or :f.

Return: TensorField

Types:

  • problem: Problem
  • dataField: Vector{Tuple{String, Float64,...}}

Examples

f1(x, y, z) = sin(x)
f2(x, y, z) = 5y
ff1 = field("face1", fx=f1, fy=f2, fz=0, fxy=1, fyz=1, fzx=f2)
ff2 = field("face2", fx=f2, fy=f1, fz=1, fxy=1, fyz=f1, fzx=1)
qq = tensorField(problem, [ff1, ff2])
qq0 = showDoFResults(qq)

qq2 = tensorField(problem, "body", [1 0 0; 0 2 0; 0 0 f1])

Here TensorField is defined in nodes.

source
LowLevelFEM.vectorFieldMethod
vectorField(problem, dataField; type=...)

Defines a vector field from dataField, which is a tuple of name of physical group and prescribed values or functions. Mesh details are in problem. type can be an arbitrary Symbol, e.g., :u or :f.

Return: VectorField

Types:

  • problem: Problem
  • dataField: Vector{Tuple{String, Float64,...}}

Examples

f1(x, y, z) = sin(x)
f2(x, y, z) = 5y
ff1 = field("face1", fx=f1, fy=f2, fz=0)
ff2 = field("face2", fx=f2, fy=f1, fz=1)
qq = vectorField(problem, [ff1, ff2])
qq0 = showDoFResults(qq)

qq2 = vectorField(problem, "body", [1, 2, ff2])

Here VectorField is defined in nodes.

source
LowLevelFEM.∂Method
∂(r::ScalarField, dir::Int) -> ScalarField

Compute the spatial derivative of a scalar field in the global x, y, or z direction. The input field may be nodal or elementwise; nodal fields are automatically converted to elementwise form using nodesToElements(r).

Arguments

  • r::ScalarField: scalar field to differentiate.
  • dir::Int: spatial direction:
    • 1 → ∂/∂x
    • 2 → ∂/∂y
    • 3 → ∂/∂z

Returns

A new ScalarField in elementwise representation, containing the gradient component ∂r/∂(dir) at all nodes of each element.

Notes

  • Works for 1D, 2D, and 3D meshes.
  • Uses the element Jacobian and shape-function derivatives in global coordinates.
  • The global problem mesh is accessed via r.model.
source
LowLevelFEM.∂eMethod
∂e(r::ScalarField, dir::Int) -> ScalarField

Compute the spatial derivative of an elementwise scalar field purely at the element level, without nodal conversion. This avoids mixing contributions from neighboring elements and is suitable for discontinuous fields or post-processing of elementwise quantities.

Arguments

  • r::ScalarField: must already be elementwise (isElementwise(r) == true).
  • dir::Int: spatial direction (1 → x, 2 → y, 3 → z).

Returns

An elementwise ScalarField containing the derivative ∂r/∂(dir) at the local nodes of each element.

Notes

  • Uses the element-local Jacobian per integration point.
  • Safer for interfaces or discontinuities than nodal differentiation.
  • Only elementwise → elementwise transformations are performed.
source
LowLevelFEM.∂xMethod
∂x(r::ScalarField)
∂y(r::ScalarField)
∂z(r::ScalarField)

Compute directional derivatives of a scalar field r with respect to the global coordinates x, y, or z.

These functions are convenience wrappers that automatically choose between the nodal differentiation () and the elementwise differentiation (∂e) depending on the representation of the input field.

Dispatch logic

  • If isNodal(r) is true:
    • ∂x(r)∂(r, 1)
    • ∂y(r)∂(r, 2)
    • ∂z(r)∂(r, 3)
  • Otherwise (elementwise field):
    • ∂x(r)∂e(r, 1)
    • ∂y(r)∂e(r, 2)
    • ∂z(r)∂e(r, 3)

Coordinate conventions

  • In 3D: (x, y, z) are the global Cartesian coordinates.
  • In 2D (plane stress/strain):
    • ∂x and ∂y operate in the plane.
    • ∂z is mathematically allowed and returns zero for planar meshes.

Returns

A ScalarField storing the derivative in elementwise representation.

Examples

1) Derivative of a nodal field (automatic nodal→element conversion)

s = ScalarField(prob, "vol", (x,y,z) -> x^2 + y)

dxs = ∂x(s)      # computes 2x
dys = ∂y(s)      # computes 1
dzs = ∂z(s)      # computes 0

2) Derivative of an elementwise field

se = elementsToNodes(s) |> r -> ∂e(r, 1)  # manually
dxs = ∂x(se)                              # same result, auto-detected

3) Works on 2D meshes as well

s = ScalarField(prob2D, "surf", (x,y,z) -> x - y)

dxs = ∂x(s)      # = 1 everywhere
dys = ∂y(s)      # = -1 everywhere
dzs = ∂z(s)      # = zero field (2D mesh)

4) Using derivatives to build a gradient vector field

s = ScalarField(prob, "vol", (x,y,z) -> sin(x)*y)

gx = ∂x(s)       # y*cos(x)
gy = ∂y(s)       # sin(x)
gz = ∂z(s)       # 0

grad_s = VectorField([gx, gy, gz])

Notes

  • The output is always elementwise, enabling consistent post-processing.
  • Direction indices follow FEM convention: 1 → x, 2 → y, 3 → z.
source