Postprocessing API

Visualization, probing, path plots, and result extraction helpers.

Derived Quantities

LowLevelFEM.fieldErrorFunction
fieldError(F::Union{ScalarField,VectorField,TensorField})

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.resultantFunction
resultant(field::VectorField, phName::String)

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.integrateFunction
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
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.normalVectorFunction
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.tangentVectorFunction
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

Result Views

LowLevelFEM.showDoFResultsFunction
showDoFResults(q::Union{ScalarField,VectorField,TensorField}, comp::Symbol; 
               name=comp, visible=false, 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, :x, :y, :z, :xy, :yz, :zx, :yx, :zy, :xz, :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.showModalResultsFunction
showModalResults(Φ::Eigen; name="", visible=false)

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.showBucklingResultsFunction
showBucklingResults(Φ::Eigen; name="", visible=false)

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.showStrainResultsFunction
showStrainResults(E::TensorField, comp::Symbol; 
                  name=comp, visible=false, smooth=true)

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, :x, :y, :z, :xy, :yz, :zx, :yx, :zy, :xz), 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.showStressResultsFunction
showStressResults(S0::TensorField, comp::Symbol; 
                  name=comp, visible=false, smooth=false)

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, :x, :y, :z, :xy, :yz, :zx, :yx, :zy, :xz), 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.showElementResultsFunction
showElementResults(F::Union{ScalarField,VectorField,TensorField}, comp=:Symbol; 
                   name=comp, visible=false, smooth=false, factor=0)

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.showHeatFluxResultsFunction
showHeatFluxResults(S0::VectorField, comp::Symbol; 
                    name=comp, visible=false, smooth=true, factor=0)

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, :x, :y, :z), 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.showOnSurfaceFunction
showOnSurface(field::Number, phName::String; 
              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.plotOnPathFunction
plotOnPath(pathName, field; points=100, step=..., plot=false, 
           name="field [$field] on $pathName", visible=false, offsetX=0, offsetY=0, offsetZ=0)

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

Probing and UI

LowLevelFEM.probeFunction
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
probe(A::Union{ScalarField,VectorField,TensorField}, name::String; step=1)

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.openPostProcessorFunction
openPostProcessor(; model=...)

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

Returns: nothing

Types:

  • model: Int64
source