Preprocessing API

Preprocessing includes model setup, material/problem definitions, boundary/load declarations, and mesh-related helpers.

Mesh and Geometry Setup

LowLevelFEM.structured_rect_meshFunction
structured_rect_mesh(; x0=0.0, y0=0.0,
                       lx=1.0, ly=1.0,
                       n=10,
                       dx=lx/n, dy=ly/n,
                       order=1)

Create a structured quadrilateral mesh of a rectangular 2D domain using Gmsh.

The domain is defined by its lower-left corner (x0, y0) and its dimensions lx × ly. A transfinite structured mesh is generated and recombined into quadrilateral elements.

Keyword Arguments

  • x0, y0: Coordinates of the lower-left corner.
  • lx, ly: Domain dimensions.
  • n: Default subdivision parameter (used to define dx, dy).
  • dx, dy: Target element sizes in x and y directions.
  • order: Polynomial order of the finite elements.

Physical Groups

The following physical groups are created:

  • 1D boundaries: "bottom", "right", "top", "left"
  • 2D domain: "surface"
  • 0D corner points: "leftbottom", "rightbottom", "righttop", "lefttop"

Notes

  • A transfinite structured mesh is enforced.
  • Elements are recombined into quadrilaterals.
  • The mesh is generated but not written to file.
  • The function assumes that Gmsh is already initialized.

Returns nothing.

Note: These helper functions are primarily intended for rapid prototyping, testing and educational examples in LowLevelFEM workflows.

source
LowLevelFEM.structured_box_meshFunction
structured_box_mesh(; x0=0.0, y0=0.0, z0=0.0,
                      lx=1.0, ly=1.0, lz=1.0,
                      n=10,
                      dx=lx/n, dy=ly/n, dz=lz/n,
                      order=1)

Create a structured hexahedral mesh of a rectangular 3D domain using Gmsh.

The domain is defined by its lower-front-left corner (x0, y0, z0) and its dimensions lx × ly × lz. A fully transfinite structured mesh is generated and recombined into hexahedral elements.

Keyword Arguments

  • x0, y0, z0: Coordinates of the origin corner.
  • lx, ly, lz: Domain dimensions.
  • n: Default subdivision parameter.
  • dx, dy, dz: Target element sizes in each direction.
  • order: Polynomial order of the finite elements.

Physical Groups

The following physical groups are created:

  • 2D boundary faces: "left", "right", "front", "rear", "bottom", "top"
  • 3D domain: "volume"

Notes

  • A transfinite structured mesh is enforced in all directions.
  • Elements are recombined into hexahedra.
  • The mesh is generated but not written to file.
  • The function assumes that Gmsh is already initialized.

Returns nothing.

Note: These helper functions are primarily intended for rapid prototyping, testing and educational examples in LowLevelFEM workflows.

source
LowLevelFEM.projectTo2DFunction
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.expandTo3DFunction
expandTo3D(V2D::VectorField)

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

return: VectorField

Examples

V3D = expandTo3D(V2D)
source
LowLevelFEM.rotateNodesFunction
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

Model Setup and Parameters

LowLevelFEM.setParameterFunction
setParameter(name, value)

Defines a parameter name and sets its value to value.

Returns: nothing

Types:

  • name: String
  • value: Float64
source
LowLevelFEM.setParametersFunction
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.openPreProcessorFunction
openPreProcessor(; openGL=true)

Launches the Gmsh preprocessor window with OpenGL disabled by default.

Returns: nothing

Types:

  • openGL: Boolean
source
LowLevelFEM.openGeometryFunction
openGeometry(name::String)

Open a geometry or mesh file in the current Gmsh session.

If Gmsh is not yet initialized, the function calls gmsh.initialize(). If the file name does not exist and has extension .geo, a new geometry file is created and initialized with a minimal header. Otherwise an error is thrown.

Before opening the file the current Gmsh model is cleared with gmsh.clear(), so any previously loaded geometry or mesh is removed.

Arguments

  • name::String: Path to a geometry (.geo) or mesh file (.msh, .step, .brep, etc.) that should be opened in Gmsh.

Behavior

  • Existing files are opened directly.
  • Missing .geo files are automatically created.
  • Missing files with other extensions result in an error.

Notes

This function only loads the geometry into Gmsh. The model name used later by Problem is obtained from gmsh.model.getCurrent().

Examples

openGeometry("beam.geo")     # create if missing and open
openGeometry("mesh.msh")     # open existing mesh
source

Boundary Conditions and Loads

LowLevelFEM.displacementConstraintFunction
displacementConstraint(name; ux=nothing, uy=nothing, uz=nothing)

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

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.temperatureConstraintFunction
temperatureConstraint(name; T=nothing)

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

Types:

  • name: String
  • T: Float64 or Function
source
LowLevelFEM.elasticSupportFunction
elasticSupport(name; kx=nothing, ky=nothing, kz=nothing)

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

Types:

  • name: String
  • kx: Float64 or Function
  • ky: Float64 or Function
  • kz: Float64 or Function
source
LowLevelFEM.loadFunction
load(name; fx=nothing, fy=nothing, fz=nothing, p=nothing)

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

Types:

  • name: String
  • fx: Float64 or Function
  • fy: Float64 or Function
  • fz: Float64 or Function
source
LowLevelFEM.heatFluxFunction
heatFlux(name; qn=nothing)

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

Types:

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

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

Types:

  • name: String
  • h: Float64 or Function
source
LowLevelFEM.heatConvectionFunction
heatConvection(name; h=nothing, T∞=nothing)

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

Types:

  • name: String
  • h: Float64
  • Tₐ: Float64 or Function
source
LowLevelFEM.BoundaryConditionFieldsFunction
BoundaryConditionFields(bc::BoundaryCondition)

Return the defined field keys of a boundary condition.

Arguments

  • bc::BoundaryCondition: Boundary condition object with named values.

Returns

  • An iterator over Symbol keys present in bc.values.

Example

keys = BoundaryConditionFields(bc)
source
LowLevelFEM.BoundaryCondition_to_LoadConditionFunction
BoundaryCondition_to_LoadCondition(bc::Vector{BoundaryCondition})

Convert boundary-condition entries to load-condition entries.

Arguments

  • bc::Vector{BoundaryCondition}: Boundary condition definitions.

Returns

  • Vector{LoadCondition}: Load conditions with copied phName, problem, and values.

Example

loads = BoundaryCondition_to_LoadCondition(bc)
source

DOF Selection Helpers

LowLevelFEM.constrainedDoFsFunction
constrainedDoFs(problem::Problem, supports::Vector{BoundaryCondition})

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{BoundaryCondition}
  • DoFs: Vector{Int64}
source
constrainedDoFs(K::SystemMatrix, 
                support::Vector{BoundaryCondition})

Return global constrained DOFs for single- or multi-field systems.

source
LowLevelFEM.freeDoFsFunction
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{BoundaryCondition}
  • DoFs: Vector{Int64}
source
freeDoFs(K::SystemMatrix, 
          support::Vector{BoundaryCondition})

Return global free DOFs for single- or multi-field systems.

source
LowLevelFEM.allDoFsFunction
allDoFs(problem)

Returns the serial numbers of degrees of freedom of the whole model.

Return: DoFs

Types:

  • problem: Problem
  • DoFs: Vector{Int64}
source