LowLevelFEM.jl Documentation
Functions
Main.LowLevelFEM.Problem
— TypeProblem(; E=..., ν=..., ρ=..., thickness=..., type=...)
A structure containing the most important data of the problem.
- name of the model (in gmsh)
- type of the problem: 3D "Solid", "PlaneStrain" or "PlaneStress"
- dimension of the problem, determined from
type
- material constants: Young's modulus, Poisson's ratio, mass density
- thickness of the plate
- number of nodes (non)
Types:
name
: stringtype
: stringdim
: IntegerE
: doubleν
: doubleρ
: doublethickness
: doublenon
: Integer
Main.LowLevelFEM.StressField
— TypeStressField(sigma, numElem, nsteps)
A structure containing the data of a stress field.
- sigma: vector of ElementNodeData type stress data (see gmsh.jl)
- numElem: vector of tags of elements
- nsteps: number of stress fields stored in sigma (for animations).
Types:
sigma
: Vector{Matrix{Float64}}numElem
: Vector{Integer}nsteps
: Integer
Main.LowLevelFEM.CDM
— MethodFEM.CDM(K, M, C, f, u0, v0, T, Δt)
Solves a transient dynamic problem using central difference method (explicit). K
is the stiffness Matrix, M
is the mass matrix, C
is the damping matrix, f
is the load vector, u0
is the initial displacement, v0
is the initial velocity, T
is the upper bound ot the time intervall (lower bound is zero) and Δt
is the time step size. Returns the displacement vectors and velocity vectors in each time step arranged in the columns of the two matrices u
and v
and a vector t
of the time instants used.
Return: u
, v
, t
Types:
K
: MatrixM
: MatrixC
: Matrixf
: Vectoru0
: Vectorv0
: VectorT
: DoubleΔt
: Doubleu
: Matrixv
: Matrixt
: Vector
Main.LowLevelFEM.applyBoundaryConditions!
— MethodFEM.applyBoundaryConditions!(problem, stiffMat, loadVec, supports)
Applies displacement boundary conditions supports
on a stiffness matrix stiffMat
and load vector loadVec
. Mesh details are in problem
. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
.
Return: none
Types:
problem
: ProblemstiffMat
: MatrixloadVec
: Vectorsupports
: Tuple(string, double, double, double)
Main.LowLevelFEM.applyBoundaryConditions!
— MethodFEM.applyBoundaryConditions!(problem, stiffMat, massMat, dampMat, loadVec, supports)
Applies displacement boundary conditions supports
on a stiffness matrix stiffMat
, mass matrix massMat
, damping matrix dampMat
and load vector loadVec
. Mesh details are in problem
. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
.
Return: none
Types:
problem
: ProblemstiffMat
: MatrixloadVec
: Vectorsupports
: Tuple(string, double, double, double)
Main.LowLevelFEM.applyBoundaryConditions
— MethodFEM.applyBoundaryConditions(problem, stiffMat, loadVec, supports)
Applies displacement boundary conditions supports
on a stiffness matrix stiffMat
and load vector loadVec
. Mesh details are in problem
. supports
is a tuple of name
of physical group and prescribed displacements ux
, uy
and uz
. Creates a new stiffness matrix and load vector.
Return: stiffMat
, loadVec
Types:
problem
: ProblemstiffMat
: MatrixloadVec
: Vectorsupports
: Tuple(string, double, double, double)
Main.LowLevelFEM.displacementConstraint
— MethodFEM.displacementConstraint(name; ux=..., uy=..., uz=...)
Gives the displacement constraints on name
physical group. At least one ux
, uy
or uz
value have to be given (depending on the dimension of the problem).
Return: none
Types:
name
: stringux
: doubleuy
: doubleuz
: double
Main.LowLevelFEM.generateMesh
— MethodFEM.generateMesh(...)
Gives...
Return: none
Types:
- ``: x
Main.LowLevelFEM.getTagForPhysicalName
— MethodFEM.getTagForPhysicalName(name)
Returns tags
of elements of physical group name
.
Return: tags
Types:
name
: stringtags
: vector of integers
Main.LowLevelFEM.initialDisplacement!
— MethodFEM.initialDisplacement!(problem, name, u0; ux=..., uy=..., uz=...)
Changes the displacement values ux
, uy
and uz
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Original values are in displacement vector u0
.
Return: none
Types:
problem
: Problemname
: stringu0
: Vectorux
: Doubleuy
: Doubleuz
: Double
Main.LowLevelFEM.initialVelocity!
— MethodFEM.initialVelocity!(problem, name, v0; vx=..., vy=..., vz=...)
Changes the velocity values vx
, vy
and vz
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Original values are in velocity vector v0
.
Return: none
Types:
problem
: Problemname
: stringv0
: Vectorvx
: Doublevy
: Doublevz
: Double
Main.LowLevelFEM.load
— MethodFEM.load(name; fx=..., fy=..., fz=...)
Gives the intensity of distributed load on name
physical group. At least one fx
, fy
or fz
value have to be given (depending on the dimension of the problem).
Return: none
Types:
name
: stringux
: doubleuy
: doubleuz
: double
Main.LowLevelFEM.loadVector
— MethodFEM.loadVector(problem, loads)
Solves a load vector of problem
. loads
is a tuple of name of physical group name
, coordinates fx
, fy
and fz
of the intensity of distributed force. It can solve traction or body force depending on the problem. In case of 2D problems and Line physical group means surface force. In case of 2D problems and Surface physical group means body force. In case of 3D problems and Line physical group means edge force. In case of 3D problems and Surface physical group means surface force. In case of 3D problems and Volume physical group means body force.
Return: loadVec
Types:
problem
: Problemloads
: Tuple(string, double, double, double)loadVec
: Vector
Main.LowLevelFEM.massMatrix
— MethodFEM.massMatrix(problem; name, ρ=..., lumped=...)
Solves the mass matrix of the problem
. If lumped
is true, solves lumped mass matrix.
Return: massMat
Types:
problem
: Problemname
: string - unusedρ
: double - unusedlumped
: booleanmassMat
: Matrix
Main.LowLevelFEM.nodalAcceleration!
— MethodFEM.nodalAcceleration!(problem, name, a0; ax=..., ay=..., az=...)
Changes the acceleration values ax
, ay
and az
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Original values are in acceleration vector a0
.
Return: none
Types:
problem
: Problemname
: stringa0
: Vectorax
: Doubleay
: Doubleaz
: Double
Main.LowLevelFEM.nodalForce!
— MethodFEM.nodalForce!(problem, name, f0; fx=..., fy=..., fz=...)
Changes the force values fx
, fy
and fz
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Original values are in load vector f0
.
Return: none
Types:
problem
: Problemname
: stringf0
: Vectorfx
: Doublefy
: Doublefz
: Double
Main.LowLevelFEM.plotOnPath
— MethodFEM.plotOnPath(problem, pathName, field, points; numOfSteps=..., name=..., visible=...)
Load a 2D plot on a path into a View in gmsh. field
is the number of View in gmsh from which the data of a field is imported. pathName
is the name of a physical group which contains a curve. The curve is devided into equal length intervals with number of points
points. The field is shown at this points. numOfSteps
is the sequence number of steps. name
is the title of graph and visible
is a true or false value to toggle on or off the initial visibility in gmsh. This function returns the tag of View.
Return: tag
Types:
problem
: ProblempathName
: stringfield
: Integerpoints
: IntegernumOfStep
: Integername
: Stringvisible
: Boolean
Main.LowLevelFEM.showDoFResults
— MethodFEM.showDoFResults(problem, q, comp; t=..., name=..., visible=...)
Loads nodal results into a View in gmsh. q
is the field to show, comp
is the component of the field ("uvec", "ux", "uy", "uz", "vvec", "vx", "vy", "vz"), t
is a vector of time steps (same number of columns as q
), name
is a title to display and visible
is a true or false value to toggle on or off the initial visibility in gmsh. If q
has more columns, then a sequence of results will be shown (eg. as an animation). This function returns the tag of View.
Return: tag
Types:
problem
: Problemq
: Vector or Matrixcomp
: stringt
: Vectorname
: stringvisible
: Boolean
Main.LowLevelFEM.showStressResults
— MethodFEM.showStressResults(problem, S, comp; t=..., 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"), t
is a vector of time steps (same length as the number of stress states), name
is a title to display, visible
is a true or false value to toggle on or off the initial visibility in gmsh and smooth
is a true of false value to toggle smoothing the stress field on or off. If length of t
is more than one, then a sequence of results will be shown (eg. as an animation). This function returns the tag of View.
Return: tag
Types:
problem
: ProblemS
: StressFieldcomp
: stringt
: Vectorname
: stringvisible
: Booleansmooth
: Boolean
Main.LowLevelFEM.smallestPeriodTime
— MethodFEM.smallestPeriodTime(K, M)
Solves the smallest period of time for a dynamic problem given by stiffness matrix K
and the mass matrix M
.`
Return: Δt
Types:
K
: MatrixM
: MatrixΔt
: Double
Main.LowLevelFEM.solveDisplacement
— MethodFEM.solveDisplacement(K, q)
Solves the equation K*q=f for the displacement vector q
. K
is the stiffness Matrix, q
is the load vector.
Return: q
Types:
K
: Matrixf
: Vectorq
: Vector
Main.LowLevelFEM.solveStress
— MethodFEM.solveStress(problem, q)
Solves the stress field S
from displacement vector q
. Stress field is given per elements, so it usually contains jumps at the boundary of elements. Details of mesh is available in problem
.
Return: S
Types:
problem
: Problemq
: VectorS
: StressField
Main.LowLevelFEM.stiffnessMatrix
— MethodFEM.stiffnessMatrix(problem; name, E=..., ν=...)
Solves the stiffness matrix of the problem
.
Return: stiffMat
Types:
problem
: Problemname
: string - unusedE
: double - unusedν
: double - unusedstiffMat
: Matrix
Index
Main.LowLevelFEM.Problem
Main.LowLevelFEM.StressField
Main.LowLevelFEM.CDM
Main.LowLevelFEM.applyBoundaryConditions
Main.LowLevelFEM.applyBoundaryConditions!
Main.LowLevelFEM.applyBoundaryConditions!
Main.LowLevelFEM.displacementConstraint
Main.LowLevelFEM.generateMesh
Main.LowLevelFEM.getTagForPhysicalName
Main.LowLevelFEM.initialDisplacement!
Main.LowLevelFEM.initialVelocity!
Main.LowLevelFEM.load
Main.LowLevelFEM.loadVector
Main.LowLevelFEM.massMatrix
Main.LowLevelFEM.nodalAcceleration!
Main.LowLevelFEM.nodalForce!
Main.LowLevelFEM.plotOnPath
Main.LowLevelFEM.showDoFResults
Main.LowLevelFEM.showStressResults
Main.LowLevelFEM.smallestPeriodTime
Main.LowLevelFEM.solveDisplacement
Main.LowLevelFEM.solveStress
Main.LowLevelFEM.stiffnessMatrix