Fields API

Field containers and common field-level utilities.

Field Construction and Storage

LowLevelFEM.fieldFunction
field(name; f=nothing, fx=nothing, fy=nothing, fz=nothing, 
            fxy=nothing, fyx=nothing, fyz=nothing, 
            fzy=nothing, fzx=nothing, fxz=nothing)

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: BoundaryCondition

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.scalarFieldFunction
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{BoundaryCondition}

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.vectorFieldFunction
vectorField(problem::Problem, dataField::Vector{BoundaryCondition})

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{BoundaryCondition}

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.tensorFieldFunction
tensorField(problem::Problem, dataField::Vector{BoundaryCondition}; type=:e)

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{BoundaryCondition}

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.mergeFieldsFunction
mergeFields(uu::Vector{<:AbstractField})

Concatenates multiple Fields along a common pseudo-time axis.

Interpretation of nsteps:

  • If nsteps == 1, the Field is treated as an iteration state (pseudo-time index: 1, 2, 3, ...).
  • If nsteps > 1, the Field is treated as a time history, and its internal time vector t is preserved.

When concatenating multiple time histories, continuity of time is enforced by shifting each subsequent segment such that it starts after the previous one. The time increment between segments is estimated from the last two time points of the previous Field. If this information is unavailable (nsteps == 1), a unit increment (Δt = 1) is assumed.

This function is intended for post-processing and visualization (e.g. merging iteration histories or successive load/time segments), not for time-accurate integration.

All input Fields must:

  • have the same number of degrees of freedom,
  • belong to the same model,
  • share the same Field type (ScalarField, VectorField, or TensorField).

Returns a Field of the same concrete type as the first element of uu.

Exaples

u1 = solveDisplacement(problem, load=[load1]) # nsteps = 1 u2 = solveDisplacement(problem, load=[load2]) # nsteps = 1 u3 = solveDisplacement(problem, load=[load3]) # nsteps = 1

u_steps = mergeFields([u1, u2, u3])

Warning

The time increment between concatenated segments is estimated heuristically.

When a Field has nsteps == 1, it is interpreted as an iteration state and a unit pseudo-time increment (Δt = 1) is assumed when concatenating with subsequent segments.

This heuristic ensures continuity of the time axis but does NOT represent physical time. The resulting time vector should therefore be used only for post-processing and visualization, not for time-accurate analysis.

source
LowLevelFEM.saveFieldFunction
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.loadFieldFunction
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.isSavedFunction
isSaved(fileName::String)

Checks whether a variable has been saved or not.

Returns: Boolean

Types:

  • fileName: String
source

Field Access and Conversion

LowLevelFEM.DoFsFunction
DoFs(f::AbstractField)

Return the vector of degrees of freedom (DOFs) associated with a field.

This corresponds to the nodal representation used in linear algebra operations (assembly, solving, post-processing).

source
LowLevelFEM.elementsToNodesFunction
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.nodesToElementsFunction
nodesToElements(T::Union{ScalarField,VectorField,TensorField}; 
                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.isNodalFunction
isNodal(field::Union{ScalarField,VectorField,TensorField})

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.isElementwiseFunction
isElementwise(field::Union{ScalarField,VectorField,TensorField})

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

Field Differential Helpers

LowLevelFEM.∂xFunction
∂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