Nonlinear Mechanics and Fields
Base.div
— Methoddiv(r::Union{VectorField,TensorField})
Solves the divergence of the vector field or tensor field r
. An alternative way to solve div
is to use ∇
as a differencial operator.
Return: ScalarField or VectorField
Types:
r
: VectorField or TensorField
3D Examples (assumes problem
is set as in the ∇ doc setup)
# Assumes a 3D mesh with physical group "body".
# 1) Divergence of a 3D vector field → ScalarField
v1(X,Y,Z) = X
v2(X,Y,Z) = Y
v3(X,Y,Z) = Z
v = vectorField(problem, [field("body", fx=v1, fy=v2, fz=v3)])
D1 = div(v)
D2 = ∇ ⋅ v
# 2) Divergence of a 3D tensor field → VectorField
fsz(X, Y, Z) = 10 - Z
S = tensorField(problem, [field("body", fz=fsz)])
b1 = -div(S)
b2 = -S ⋅ ∇
LowLevelFEM.applyDeformationBoundaryConditions!
— MethodapplyDeformationBoundaryConditions!(deformVec, supports)
Applies displacement boundary conditions supports
on deformation vector deformVec
. Mesh details are in problem
. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
.
Returns: nothing
Types:
deformVec
: VectorFieldsupports
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.curl
— Methodcurl(r::VectorField)
Solves the rotation of the vector field r
. An alternative way to solve curl
is to use ∇
as a differencial operator.
Return: VectorField
Types:
r
: VectorField
3D Example (assumes problem
is set as in the ∇ doc setup)
# Assumes a 3D mesh with physical group "body".
vx(X, Y, Z) = 0
vy(X, Y, Z) = X
vz(X, Y, Z) = 0
v = vectorField(problem, [field("body", fx=vx, fy=vy, fz=vz)])
D1 = curl(v)
D2 = ∇ × v
LowLevelFEM.equivalentNodalForce
— MethodequivalentNodalForce(r::VectorField)
Solves the equivalent nodal force (when solving large deformation problems). (See [6]) r
is the position vector field in the current configuration.
Return: VectorField
Types:
r
: VectorField
LowLevelFEM.grad
— Methodgrad(r::Union{ScalarField,VectorField})
Solves the gradient of the scalar field or vector field r
. An alternative way to solve grad
is to use ∇
as a differencial operator.
Return: VectorField or TensorField
Types:
r
: ScalarField or VectorField
3D Examples
# Assumes a 3D mesh with physical group "body".
# 1) Gradient of a 3D scalar field → VectorField
f(X,Y,Z) = X^2 + Y*Z
S = scalarField(problem, [field("body", f=f)])
G1 = grad(S)
G2 = ∇(S)
# 2) Gradient of a 3D vector field → TensorField
vx(X,Y,Z) = X
vy(X,Y,Z) = Y
vz(X,Y,Z) = Z
V = vectorField(problem, [field("body", fx=vx, fy=vy, fz=vz)])
T1 = grad(V)
T2 = V ∘ ∇
LowLevelFEM.nodePositionVector
— MethodnodePositionVector(problem)
Returns the position vectors of all mesh nodes as a VectorField
(initial configuration).
Returns: R
Types:
problem
: ProblemR
: VectorField
LowLevelFEM.nonFollowerLoadVector
— MethodnonFollowerLoadVector(r::VectorField, load)
Solves the non-follower load vector (when solving large deformation problems). r
is the position vector field in the current configuration.
Return: VectorField
Types:
r
: VectorField
LowLevelFEM.rot
— Methodrot(r::VectorField)
Solves the rotation of the vector field r
. In some countries "rot" denotes the English "curl". (See the curl
function.)
Return: VectorField
Types:
r
: VectorField
LowLevelFEM.showDeformationResults
— MethodshowDeformationResults(r::VectorField, comp; name=String, visible=Boolean)
Shows deformation result, where r
contains the position vectors of nodes in the current configuration.
LowLevelFEM.solveDeformation
— MethodsolveDeformation(problem::Problem, load, supp;
followerLoad=false,
loadSteps = 3,
rampedLoad = true,
rampedSupport = false,
maxIteration = 10,
saveSteps = false,
saveIterations = false,
plotConvergence = false,
relativeError = 1e-5,
initialDeformation=nodePositionVector(problem))
Solves the deformed shape of a non-linearly elastic body...
LowLevelFEM.suppressDeformationAtBoundaries!
— MethodsuppressDeformationAtBoundaries!(stiffMat, loadVec, supports)
Suppresses the displacements given in support
in stiffMat
and loadVec
so that it is only necessary to consider them once during iteration. stiffMat
is the stiffness matrix, loadVec
is the load vector. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
.
Returns: nothing
Types:
stiffMat
: SystemMatrixloadVec
: VectorFieldsupports
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.suppressDeformationAtBoundaries
— MethodsuppressDeformationAtBoundaries(stiffMat, loadVec, supports)
Suppresses the displacements given in support
in stiffMat
and loadVec
so that it is only necessary to consider them once during iteration. stiffMat
is the stiffness matrix, loadVec
is the load vector. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
.
Return: stiffMat1, loadVec1
Types:
stiffMat
: SystemMatrixloadVec
: VectorFieldsupports
: Vector{Tuple{String, Float64, Float64, Float64}}stiffMat1
: SystemMatrixloadVec1
: VectorField
LowLevelFEM.tangentMatrixConstitutive
— MethodtangentMatrixConstitutive(r::VectorField)
Solves the constitutive part of the tangent matrix (when solving large deformation problems). (See [6]) r
is the position vector field in the current configuration.
Return: SystemMatrix
Types: - r
: VectorField
LowLevelFEM.tangentMatrixInitialStress
— MethodtangentMatrixInitialStress(r::VectorField)
Solves the initial stress part of the tangent matrix (when solving large deformation problems). (See [6]) r
is the position vector field in the current configuration.
Return: SystemMatrix
Types:
r
: VectorField
LowLevelFEM.∇
— Method∇(rr::Union{VectorField, ScalarField, TensorField}; nabla=:grad)
Computes derivatives of rr
.
- If
rr
is aScalarField
andnabla == :grad
, returns the gradient (aVectorField
). - If
rr
is aVectorField
andnabla == :grad
, returns the gradient (aTensorField
). - If
rr
is aVectorField
andnabla == :curl
, returns the curl (aVectorField
). - If
rr
is aVectorField
andnabla == :div
, returns the divergence (aScalarField
). - If
rr
is aTensorField
andnabla == :div
, returns the divergence (aVectorField
).
Returns: ScalarField
, VectorField
, or TensorField
Types:
rr
:ScalarField
,VectorField
, orTensorField
nabla
: Symbol
3D Examples (assumes problem
is set as in the ∇ doc setup)
# One-time 3D setup (assumes examples/Fields/cube.geo exists with physical group "body")
using LowLevelFEM
gmsh.initialize()
gmsh.open("examples/Fields/cube.geo")
mat = material("body", E=210e3, ν=0.3, ρ=7.85e-9)
problem = Problem([mat], type=:Solid)
# 1) Gradient of a 3D scalar field: ∇f → VectorField
f(X,Y,Z) = X^2 + Y*Z
S = scalarField(problem, [field("body", f=f)])
G = ∇(S) # VectorField with 3 components
# 2) Curl of a 3D vector field: ∇ × v → VectorField
vx(X,Y,Z) = 0
vy(X,Y,Z) = X
vz(X,Y,Z) = 0
V = vectorField(problem, [field("body", fx=vx, fy=vy, fz=vz)])
C = ∇(V, nabla=:curl) # approx (0, 0, 1) everywhere
# 3) Divergence of a 3D vector field: ∇ ⋅ v → ScalarField
v1(X,Y,Z) = X
v2(X,Y,Z) = Y
v3(X,Y,Z) = Z
V2 = vectorField(problem, [field("body", fx=v1, fy=v2, fz=v3)])
D = ∇(V2, nabla=:div) # ≈ 3
# 4) Divergence of a 3D tensor field: ∇ · T → VectorField (if T is TensorField)
# For example, a diagonal tensor T with only Tzz = g(Z): div(T) = (0, 0, ∂g/∂Z)
g(Z) = 10 - Z
T = tensorField(problem, [field("body", fz=g)])
DV = ∇(T, nabla=:div) # VectorField
# Symmetric displacement gradient via operators
# A = (u ∘ ∇ + ∇ ∘ u) / 2
gmsh.finalize()
- 6Javier Bonet, Richard D. Wood: Nonlinear Continuum Mechanics for Finite Element Analysis, Cambridge University Press, 2008, https://doi.org/10.1017/CBO9780511755446