LowLevelFEM.jl
Documentation for LowLevelFEM.jl
Functions
LowLevelFEM.CoordinateSystem
— TypeCoordinateSystem(vec1, vec2)
A structure containing the data of a coordinate system.
vec1
: direction of the new x axis.vec2
: together withvec1
determine the xy plane
If the problem is two dimensional, it is enough to give the first two elements of vec1
. Elements of vec1
and vec2
can be functions. In 3D case the functions have three arguments (x, y, and z coordinates), otherwise (in 2D case) the number of arguments is two (x and y coordinates).
Types:
vec1
: Vector{Float64}vec2
: Vector{Float64}
LowLevelFEM.Eigen
— TypeEigen(f, ϕ)
A structure containing the eigenfrequencies and eigen modes.
- f: eigenfrequencies
- ϕ: eigen modes
Types:
f
: Matrix{Float64}ϕ
: Vector{Float64}
LowLevelFEM.Material
— TypeMaterial(phName, E, ν, ρ, k, c, α)
A structure containing the material constants.
- E: elastic modulus,
- ν: Poisson's ratio,
- ρ: mass density,
- k: heat conductivity,
- c: specific heat,
- α: heat expansion coefficient
phName
is the name of the physical group where the given material is used.
Types:
phName
: StringE
: Float64ν
: Float64ρ
: Float64k
: Float64c
: Float64α
: Float64
LowLevelFEM.Problem
— TypeProblem(materials; thickness=..., type=..., bandwidth=...)
A structure containing the most important data of the problem.
- Parts of the model with their material constants. More materials can be given. (see
material
function) - type of the problem: 3D
:Solid
,:PlaneStrain
,:PlaneStress
,:AxiSymmetric
,:PlaneHeatConduction
,:HeatConduction
,:AxiSymmetricHeatConduction
. In the case of:AxiSymmetric
, the axis of symmetry is the "y" axis, while the geometry must be drawn in the positive "x" half-plane. - bandwidth optimization using built-in
gmsh
function. Possibilities::RCMK
,:Hilbert
,:Metis
or:none
(default) - dimension of the problem, determined from
type
- material constants: Young's modulus, Poisson's ratio, mass density, heat conduction corfficient, specific heat, heat expansion coefficient (in a vector of material structure
materials
) - thickness of the plate
- number of nodes (non)
Types:
materials
: Materialtype
: Symbolbandwidth
: Stringdim
: Integerthickness
: Float64non
: Integer
LowLevelFEM.TensorField
— TypeTensorField(sigma, numElem, nsteps)
A structure containing the data of a stress or strain 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).
- type: type of data (eg. stress
:s
and strain:e
)
Types:
sigma
: Vector{Matrix{Float64}}numElem
: Vector{Integer}nsteps
: Integertype
: Symbol
LowLevelFEM.VectorField
— TypeVectorField(sigma, numElem, nsteps)
A structure containing the data of a heat flux field.
- sigma: vector of ElementNodeData type heat flux data (see gmsh.jl)
- numElem: vector of tags of elements
- nsteps: number of stress fields stored in sigma (for animations).
- type: type of data (eg. heat flux
:q
)
Types:
sigma
: Vector{Matrix{Float64}}numElem
: Vector{Integer}nsteps
: Integertype
: Symbol
LowLevelFEM.CDM
— MethodFEM.CDM(K, M, C, f, u0, v0, T, Δt)
Solves a transient dynamic problem using central difference method (CDM) (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 of 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.
The critical (largest allowed) time step is Δtₘₐₓ = Tₘᵢₙ / π * (√(1 + ξₘₐₓ^2) - ξₘₐₓ)
where Tₘᵢₙ
is the time period of the largest eigenfrequency and ξₘₐₓ
is the largest modal damping.
Return: u
, v
, t
Types:
K
: SparseMatrixM
: SparseMatrixC
: SparseMatrixf
: Vector{Float64}u0
: Vector{Float64}v0
: Vector{Float64}T
: Float64Δt
: Float64u
: Matrix{Float64}v
: Matrix{Float64}t
: Vector{Float64}
LowLevelFEM.CDMaccuracyAnalysis
— MethodFEM.CDMaccuracyAnalysis(ωₘᵢₙ, ωₘₐₓ, Δt, type; n=100, α=..., ξ=..., β=..., show_β=..., show_ξ=...)
Gives some functions (graphs) for accuracy analysis of the CDM method. ωₘᵢₙ
and ωₘₐₓ
are the square root of smallest and largest eigenvalues of the Kϕ=ω²Mϕ eigenvalue problem, Δt
is the time step size. type
is one of the following values:
:SR
: spectral radius:PDR
: physical damping ratio:ADR
: algorithmic damping ratio:PE
: period error
For details see [3]. n
is the number of points in the graph. The damping matrix is assembled in the following ways: C=αM+βK or C=αM+β₁K+β₂KM⁻¹K+β₃KM⁻¹KM⁻¹K+⋅⋅⋅. The latter corresponds to the damping characteristic characterized by a power series consisting of powers of the natural frequencies with odd exponents. ξᵢ (ξ
in the argument list) are the values of the individual members of the series corresponding to the ωₘₐₓ value. βᵢ (β
in the argument list) are the coefficients of the series. (see [4]) Either ξ
or β
must be specified. ξ
or β
are scalars or vectors. If show_β
or show_ξ
is true
, the corresponding β
or ξ
values will be sent to the output. Returns a tuple of x and y values of the graph. (Can be plotted with plot(xy)
)
Return: xy
Types:
K
: SparseMatrixM
: SparseMatrixωₘᵢₙ
: Float64ωₘₐₓ
: Float64Δt
: Float64n
: Int64α
: Float64β
: Float64 of Vector{Float64}ξ
: Float64 of Vector{Float64}show_β
: Booleanshow_ξ
: Booleanxy
: Tuple{Vector{Float64},Vector{Float64}}
LowLevelFEM.FDM
— MethodFEM.FDM(K, C, q, T0, tₘₐₓ, Δt; ϑ=...)
Solves a transient heat conduction problem using Finite Difference Method (FDM). Introducing a ϑ
parameter, special cases can be used as the Forward Euler (explicit, ϑ=0), Backward Euler (implicit, ϑ=1), Crank-Nicolson (ϑ=0.5) and intermediate cases (0<ϑ<1). (This method is known as ϑ-method. See [5].) K
is the heat conduction matrix, C
is the heat capacity matrix, q
is the heat flux vector, T0
is the initial temperature, tₘₐₓ
is the upper bound of the time intervall (lower bound is zero) and Δt
is the time step size. Returns the nodal temperature vectors in each time step arranged in the columns of the matrix T
and a vector t
of the time instants used.
The critical (largest allowed) time step is Δtₘₐₓ = 2 / ((1-2ϑ)*λₘₐₓ)
where λₘₐₓ
is the largest eigenvalue of (K+λC)θ=0 eigenvalue problem and ϑ
is the parameter of the ϑ-method. Default value of ϑ
is 1/2.
Return: T
, t
Types:
K
: SparseMatrixC
: SparseMatrixq
: Vector{Float64}T0
: Vector{Float64}tₘₐₓ
: Float64Δt
: Float64T
: Matrix{Float64}t
: Vector{Float64}
LowLevelFEM.HHT
— MethodFEM.HHT(K, M, f, u0, v0, T, Δt; α=..., δ=..., γ=..., β=...)
Solves a transient dynamic problem using HHT-α method[1] (implicit). K
is the stiffness Matrix, M
is the mass matrix, f
is the load vector, u0
is the initial displacement, v0
is the initial velocity, T
is the upper bound of 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. For the meaning of α
, β
and γ
see [1]. If δ
is given, γ=0.5+δ and β=0.25⋅(0.5+γ)².
Return: u
, v
, t
Types:
K
: SparseMatrixM
: SparseMatrixf
: Vector{Float64}u0
: Vector{Float64}v0
: Vector{Float64}T
: Float64Δt
: Float64α
: Float64β
: Float64γ
: Float64δ
: Float64u
: Matrix{Float64}v
: Matrix{Float64}t
: Vector{Float64}
LowLevelFEM.HHTaccuracyAnalysis
— MethodFEM.HHTaccuracyAnalysis(ωₘᵢₙ, ωₘₐₓ, Δt, type; n=100, α=0.0, δ=0.0, γ=0.5 + δ, β=0.25 * (0.5 + γ)^2)
Gives some functions (graphs) for accuracy analysis of the HHT-α method[1]. ωₘᵢₙ
and ωₘₐₓ
are the square root of smallest and largest eigenvalues of the Kϕ=ω²Mϕ eigenvalue problem, Δt
is the time step size. type
is one of the following values:
:SR
: spectral radius:ADR
: algorithmic damping ratio:PE
: period error
For details see [2] and [3]. n
is the number of points in the graph. For the meaning of α
, β
and γ
see [1]. If δ
is given, γ=0.5+δ and β=0.25⋅(0.5+γ)². Returns a tuple of x and y values of the graph. (Can be plotted with plot(xy)
)
Return: xy
Types:
ωₘᵢₙ
: Float64ωₘₐₓ
: Float64Δt
: Float64n
: Int64α
: Float64β
: Float64γ
: Float64δ
: Float64xy
: Tuple{Vector{Float64},Vector{Float64}}
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
: SparseMatrixloadVec
: Vectorsupports
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.applyBoundaryConditions!
— MethodFEM.applyBoundaryConditions!(problem, heatCondMat, heatCapMat, heatFluxVec, supports)
Applies boundary conditions supports
on a heat conduction matrix heatCondMat
, heat capacity matrix heatCapMat
and heat flux vector heatFluxVec
. Mesh details are in problem
. supports
is a tuple of name
of physical group and prescribed temperature T
.
Return: stiffMat
, loadVec
Types:
problem
: ProblemstiffMat
: SparseMatrixloadVec
: Vectorsupports
: Vector{Tuple{String, Float64, Float64, Float64}}
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
: SparseMatrixmassMat
: SparseMatrixdampMat
: SparseMatrixloadVec
: Vector{Float64}supports
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.applyBoundaryConditions!
— MethodFEM.applyBoundaryConditions!(problem, dispVec, supports)
Applies displacement boundary conditions supports
on a displacement vector dispVec
. 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
: ProblemdispVec
: Vector{Float64}supports
: Vector{Tuple{String, Float64, Float64, Float64}}
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
: SparseMatrixloadVec
: Vectorsupports
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.applyElasticSupport!
— MethodFEM.applyElasticSupport!(problem, stiffMat, elastSupp)
Applies elastic support boundary conditions elastSupp
on a stiffness matrix stiffMat
. Mesh details are in problem
. elastSupp
is a tuple of name
of physical group and prescribed kx
, ky
and kz
stiffnesses.
Return: none
Types:
problem
: ProblemstiffMat
: SparseMatrixelastSupp
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.applyHeatConvection!
— MethodFEM.applyHeatConvection!(problem, heatCondMat, heatFluxVec, heatConv)
Applies heat convectiom boundary conditions heatConv
on a heat conduction matrix heatCondMat
and heat flux vector heatFluxVec
. Mesh details are in problem
. heatConv
is a tuple of name
of physical group and prescribed heat transfer coefficient h
and ambient temperature Tₐ
. The ambient temperature can be either a constant or a function of x, y and z coordinates.
Return: none
Types:
problem
: ProblemheatCondMat
: SparseMatrixheatFluxVec
: VectorheatConv
: Vector{Tuple{String, Float64, Float64, Float64}}
LowLevelFEM.constrainedDoFs
— MethodFEM.constrainedDoFs(problem, supports)
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
: Problemsupports
: Vector{Tuple{String, Float64, Float64, Float64}}DoFs
: Vector{Int64}
LowLevelFEM.dampingMatrix
— MethodFEM.dampingMatrix(K, M, ωₘₐₓ; α=0.0, ξ=..., β=...)
Generates the damping matrix for proportional damping case. C=αM+βK or C=αM+β₁K+β₂KM⁻¹K+β₃KM⁻¹KM⁻¹K+⋅⋅⋅. The latter corresponds to the damping characteristic characterized by a power series consisting of powers of the natural frequencies with odd exponents. ξᵢ (ξ
in the argument list) are the values of the individual members of the series corresponding to the ωₘₐₓ value. βᵢ (β
in the argument list) are the coefficients of the series. (see [4]) Either ξ
or β
must be specified. ξ
or β
are scalars or vectors. K
is the stiffness matrix, M
is the mass matrix and ωₘₐₓ
is the largest natural frequency.
Return: dampingMatrix
Types:
K
: SparseMatrixM
: SparseMatrixωₘₐₓ
: Float64α
: Float64ξ
: Float64 of Vector{Float64}β
: Float64 of Vector{Float64}dampingMatrix
: SparseMatrix
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). ux
, uy
or uz
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.displacementConstraint("support1", ux=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringux
: Float64 or Functionuy
: Float64 or Functionuz
: Float64 or Function
LowLevelFEM.elasticSupport
— MethodFEM.elasticSupport(name; kx=..., ky=..., kz=...)
Gives the distributed stiffness of the elastic support on name
physical group. kx
, ky
or kz
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.elasticSupport("supp1", kx=fn)
) Default values are 1.
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringkx
: Float64 or Functionky
: Float64 or Functionkz
: Float64 or Function
LowLevelFEM.elasticSupportMatrix
— MethodFEM.elasticSupportMatrix(problem, elSupp)
Solves the elastic support matrix of the problem
. elSupp
is a vector of elastic supports defined in function FEM.elasticSupport
. With the displacementent vector q
in hand the reaction force vector fR
arising from the elastic support can be solved. (fR = heatConvMat * q
)
Return: elSuppMat
Types:
problem
: ProblemelSupp
: Vector{Tuple{String, Float64, Float64, Float64}}elSuppMat
: SparseMatrix
LowLevelFEM.elementsToNodes
— MethodFEM.elementsToNodes(problem, T)
Solves the nodal results F
from the elemental results T
. T
can be tensor field or vector field.
Return: F
Types:
problem
: ProblemT
: TensorField or VectorFieldF
: Matrix{Float64}
LowLevelFEM.field
— MethodFEM.field(name; f=..., fx=..., fy=..., fz=..., fxy=..., fyz=..., fzx=...)
Gives the value of scalar, vector or tensor field on name
physical group. At least one fx
, fy
or fz
etc. value have to be given (depending on the dimension of the problem). fx
, fy
or fz
etc. can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.field("surf1", fx=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function,...x7}
Types:
name
: Stringf
: Float64 or Functionfx
: Float64 or Functionfy
: Float64 or Functionfz
: Float64 or Functionfxy
: Float64 or Functionfyz
: Float64 or Functionfzx
: Float64 or Function
LowLevelFEM.fieldError
— MethodFEM.fieldError(problem, F)
Solves 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.
Return: e
Types:
problem
: ProblemF
: TensorField or VectorFielde
: Matrix{Float64}
LowLevelFEM.freeDoFs
— MethodFEM.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
: Problemsupports
: Vector{Tuple{String, Float64, Float64, Float64}}DoFs
: Vector{Int64}
LowLevelFEM.generateMesh
— MethodFEM.generateMesh(problem, surf, elemSize; approxOrder=1, algorithm=6, quadrangle=0, internalNodes=0)
Obsolate, use gmsh script (.geo) instead.
Return: none
Types:
- ``: x
LowLevelFEM.getTagForPhysicalName
— MethodFEM.getTagForPhysicalName(name)
Returns tags
of elements of physical group name
.
Return: tags
Types:
name
: Stringtags
: Vector{Integer}
LowLevelFEM.heatCapacityMatrix
— MethodFEM.heatCapacityMatrix(problem; lumped=...)
Solves the heat capacity matrix of the problem
. If lumped
is true, solves lumped heat capacity matrix.
Return: heatCapMat
Types:
problem
: Problemlumped
: BooleanmassMat
: SparseMatrix
LowLevelFEM.heatConductionMatrix
— MethodFEM.heatConductionMatrix(problem)
Solves the heat conduction matrix of the problem
.
Return: heatCondMat
Types:
problem
: ProblemheatCondMat
: SparseMatrix
LowLevelFEM.heatConvection
— MethodFEM.heatConvection(name; h=..., Tₐ=...)
Gives the heat convection of the surface given with name
physical group. h
is the heat transfer coefficient of the surrounding media, Tₐ
is the ambient temperature. The ambient temperature can be either a constant or a function of x, y and z.
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringh
: Float64Tₐ
: Float64 or Function
LowLevelFEM.heatConvectionMatrix
— MethodFEM.heatConvectionMatrix(problem, heatConvection)
Solves the heat convection matrix of the problem
. heatConvection
is a vector of heat convection boundary condicions defined in function FEM.heatConduction
. With the heat convection vector (see the heatConvectionVector
function) heatConvVec
, temperature field vector T
in hand the heat flux vector qCV
arising from the heat convection boundary condition can be solved. qCV = heatConvMat * T - heatConvVec
Return: heatConvMat
Types:
problem
: ProblemheatConvection
: Vector{Tuple{String, Float64, Float64, Float64}}heatConvMat
: SparseMatrix
LowLevelFEM.heatConvectionVector
— MethodFEM.heatConvectionVector(problem, heatConvection)
Solves a heat convection vector of problem
. heatConvection
is a vector of heat convection boundary condicions defined in function FEM.heatConduction
. With the heat convection matrix (see the heatConvectionMatrix
function) heatConvMat
, temperature field vector T
in hand the heat flux vector qCV
arising from the heat convection boundary condition can be solved. qCV = heatConvMat * T - heatConvVec
Return: heatConvVec
Types:
problem
: ProblemheatConvection
: Vector{Tuple{String, Float64, Float64, Float64}}heatConvVec
: Vector
LowLevelFEM.heatFlux
— MethodFEM.heatFlux(name; qn=...)
Gives the heat flux normal to the surface of name
physical group. qn
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.load("flux1", qn=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringqn
: Float64 or Function
LowLevelFEM.heatFluxVector
— MethodFEM.heatFluxVector(problem, heatFlux)
Solves a heat flux or heat source vector of problem
. heatFlux
is a tuple of name of physical group name
, heat flux qn
normal to the surface of the body. The outward direction is positive. It can solve heat flux (or heat source) depending on the problem.
- In case of 2D problems and Point physical group means concentrated heat flux.
- In case of 2D problems and Line physical group means surface heat flux.
- In case of 2D problems and Surface physical group means body heat source.
- In case of 3D problems and Point physical group means concentrated heat flux.
- In case of 3D problems and Line physical group means edge heat source.
- In case of 3D problems and Surface physical group means surface heat flux.
- In case of 3D problems and Volume physical group means body heat source.
Return: heatFluxVec
Types:
problem
: ProblemheatFlux
: Vector{Tuple{String, Float64, Float64, Float64}}heatFluxVec
: Vector
LowLevelFEM.heatSource
— MethodFEM.heatSource(name; h=...)
Gives the body heat source in name
physical group. h
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.load("source1", h=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringh
: Float64 or Function
LowLevelFEM.heatSourceVector
— MethodFEM.heatSourceVector(problem, heatSource)
Solves a heat flux or heat source vector of problem
. heatSource
is a tuple of name of physical group name
, heat flux qn
normal to the surface of the body. The outward direction is positive. It can solve heat flux (or heat source) depending on the problem.
- In case of 2D problems and Point physical group means concentrated heat flux.
- In case of 2D problems and Line physical group means surface heat flux.
- In case of 2D problems and Surface physical group means body heat source.
- In case of 3D problems and Point physical group means concentrated heat flux.
- In case of 3D problems and Line physical group means edge heat source.
- In case of 3D problems and Surface physical group means surface heat flux.
- In case of 3D problems and Volume physical group means body heat source.
Same as the heatFluxVector
function.
Return: heatSourceVec
Types:
problem
: ProblemheatSource
: Vector{Tuple{String, Float64, Float64, Float64}}heatSourceVec
: Vector
LowLevelFEM.initialDisplacement!
— MethodFEM.initialDisplacement!(problem, name, u0; ux=..., uy=..., uz=...)
Changes the displacement values to 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: u0
Types:
problem
: Problemname
: Stringux
: Float64uy
: Float64uz
: Float64u0
: Vector{Float64}
LowLevelFEM.initialDisplacement
— MethodFEM.initialDisplacement(problem, name; ux=..., uy=..., uz=...)
Sets the displacement values ux
, uy
and uz
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Returns the initial displacement vector u0
.
Return: u0
Types:
problem
: Problemname
: Stringu0
: Vector{Float64}ux
: Float64uy
: Float64uz
: Float64
LowLevelFEM.initialTemperature!
— MethodFEM.initialTemperature!(problem, name, T0; T=...)
Changes the tempetature value to T
at nodes belonging to physical group name
. Original values are in temperature vector T0
.
Return: none
Types:
problem
: Problemname
: StringT0
: Vector{Float64}T
: Float64
LowLevelFEM.initialTemperature
— MethodFEM.initialTemperature(problem, name; T=...)
Sets the temperature value T
at nodes belonging to physical group name
. Returns the T0
initial nodal temperature vector.
Return: T0
Types:
problem
: Problemname
: StringT
: Float64T0
: Vector{Float64}
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
: Vector{Float64}vx
: Float64vy
: Float64vz
: Float64
LowLevelFEM.initialVelocity
— MethodFEM.initialVelocity(problem, name; vx=..., vy=..., vz=...)
Sets the velocity values vx
, vy
and vz
(depending on the dimension of the problem
) at nodes belonging to physical group name
. Returns the initial velocity vector v0
.
Return: v0
Types:
problem
: Problemname
: Stringvx
: Float64vy
: Float64vz
: Float64v0
: Vector{Float64}
LowLevelFEM.largestEigenValue
— MethodFEM.largestEigenValue(K, M)
Solves the smallest eigenvalue for a transient problem given by stiffness (heat conduction) matrix K
and the mass (heat capacity) matrix M
(C
).
Return: λₘᵢₙ
Types:
K
: SparseMatrixM
: SparseMatrixλₘᵢₙ
: Float64
LowLevelFEM.largestPeriodTime
— MethodFEM.largestPeriodTime(K, M)
Solves the largest period of time for a dynamic problem given by stiffness matrix K
and the mass matrix M
.
Return: Δt
Types:
K
: SparseMatrixM
: SparseMatrixΔt
: Float64
LowLevelFEM.latentHeatMatrix
— MethodFEM.latentHeatMatrix(problem, u, v, T0)
Solves the latent heat matrix of the problem
. With this matrix the generated heat due to deformations (given with displacement field u
and velocity field v
) can be solved. T0
is the current temperature field which is given in absolute temperature scale (Kelvin).
Return: latHeatMat
Types:
problem
: Problemu
: Vector{Float64}v
: Vector{Float64}T0
: Vector{Float64}latHeatMat
: SparseMatrix
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). fx
, fy
or fz
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.load("load1", fx=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: Stringfx
: Float64 or Functionfy
: Float64 or Functionfz
: Float64 or Function
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 Point physical group means concentrated force.
- 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 Point physical group means concentrated 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
: Vector{Tuple{String, Float64, Float64, Float64}}loadVec
: Vector
LowLevelFEM.massMatrix
— MethodFEM.massMatrix(problem; lumped=...)
Solves the mass matrix of the problem
. If lumped
is true, solves lumped mass matrix.
Return: massMat
Types:
problem
: Problemlumped
: BooleanmassMat
: SparseMatrix
LowLevelFEM.material
— MethodFEM.material(name; E=2.0e5, ν=0.3, ρ=7.85e-9, k=45, c=4.2e8, α=1.2e-5)
Returns a structure in which name
is the name of a physical group, E
is the modulus of elasticity, ν
Poisson's ratio and ρ
is the mass density, k
is the heat conductivity, c
is the specific heat, α
is the coefficient of heat expansion.
Return: mat
Types:
mat
: Materialname
: StringE
: Float64ν
: Float64ρ
: Float64k
: Float64c
: Float64α
: Float64
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
: Vector{Float64}ax
: Float64ay
: Float64az
: Float64
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
: Vector{Float64}fx
: Float64fy
: Float64fz
: Float64
LowLevelFEM.nonLinearStiffnessMatrix
— MethodFEM.nonLinearStiffnessMatrix(problem, q)
Solves the nonlinear stiffness matrix of the problem
. q
is a displacement field.
Return: stiffMat
Types:
problem
: Problemq
: Vector{Float64}stiffMat
: SparseMatrix
LowLevelFEM.openPostProcessor
— MethodFEM.openPostProcessor(; model=...)
Launches the GMSH postprocessor window with open postprocessor tree (of model
).
Return: none
Types:
model
: Int64
LowLevelFEM.openPreProcessor
— MethodFEM.openPreProcessor(; openGL=...)
Launches the GMSH preprocessor window with openGL disabled by default.
Return: none
Types:
openGL
: Boolean
LowLevelFEM.plotOnPath
— MethodFEM.plotOnPath(problem, pathName, field; points=100, step=..., plot=..., name=..., visible=..., offsetX=..., offsetY=..., offsetZ=...)
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. step
is the sequence number of displayed step. If no step is given, shows all the aviable steps as an animation. If plot
is true, additional return parameter, a tuple of vectors is given back, in which x
is a vector of values in horizontal axis, y
is a vector of values in vertical axis of a plot (see Plots
package). 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
or
Return: tag
, xy
Types:
problem
: ProblempathName
: Stringfield
: Integerpoints
: Integerstep
: Integerplot
: Booleanname
: Stringvisible
: Booleantag
: Integerxy
: Tuples{Vector{Float64},Vector{Float64}}
LowLevelFEM.resultant
— MethodFEM.resultant(problem, field, phName; grad=false, component=:x)
Solves the resultant of field
on phName
physical group. Return the resultant(s) in tuple
. Number of the members in tuple
depends on the dimension of problem
. It can solve the resultant of a load vector (sum of the elements of the vector), if field
is a vector of floats. If field
is a view (tag of a view in gmsh), then the integral of the field is solved. field
must have only one component. If grad
is true
, then the gradient of the field
will be evaluated and component
of the gradient (:x
, :y
or :z
) will be used to solve the resultant.
Return: res
`
or
Return: resx
, resy
or
Return: resx
, resy
, resz
Types:
field
: Vector{Float64}phName
: Stringdim
: Int64res
: Float64resx
: Float64resy
: Float64resz
: Float64
LowLevelFEM.rotateNodes
— MethodFEM.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
: ProblemphName
: StringCoordSys
: CoordinateSystemT
: SparseMatrix
LowLevelFEM.scalarField
— MethodFEM.scalarField(problem, entity)
Defines a scalar field from entity
, which is a tuple of name
of physical group and prescribed values or functions. Mesh details are in problem
.
Return: Vector{Float64}
Types:
problem
: Problementity
: Vector{Tuple{String, Float64,...}}
LowLevelFEM.setParameter
— MethodFEM.setParameter(name, value)
Defines a parameter name
and sets its value to value
.
Return: none
Types:
name
: Stringvalue
: Float64
LowLevelFEM.showBucklingResults
— MethodFEM.showBucklingResults(problem, Φ, name=..., visible=...)
Loads buckling results into a View in gmsh. Φ
is a struct of Eigen. name
is a title to display and visible
is a true or false value to toggle on or off the initial visibility in gmsh. Click on ▷| to change the results. This function returns the tag of View.
Return: tag
Types:
problem
: ProblemΦ
: Eigenname
: Stringvisible
: Booleantag
: Integer
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 (: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, :seqv, :scalar), 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{Matrix}comp
: Symbolt
: Vector{Float64}name
: Stringvisible
: Booleantag
: Integer
LowLevelFEM.showElementResults
— MethodFEM.showElementResults(problem, F, comp; t=..., name=..., visible=..., smooth=...)
Same as ShowStressResults
or showStrainResults
, depending on the type of F
data field.
Return: tag
Types:
problem
: ProblemF
: TensorFieldcomp
: Symbolt
: Vector{Float64}name
: Stringvisible
: Booleansmooth
: Booleantag
: Integer
LowLevelFEM.showHeatFluxResults
— MethodFEM.showHeatFluxResults(problem, Q, comp; t=..., name=..., visible=..., smooth=...)
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), 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
: TensorFieldcomp
: Symbolt
: Vector{Float64}name
: Stringvisible
: Booleansmooth
: Booleantag
: Integer
LowLevelFEM.showModalResults
— MethodFEM.showModalResults(problem, Φ, name=..., visible=...)
Loads modal results into a View in gmsh. Φ
is a struct of Eigen. name
is a title to display and visible
is a true or false value to toggle on or off the initial visibility in gmsh. Click on ▷| to change the results. This function returns the tag of View.
Return: tag
Types:
problem
: ProblemΦ
: Eigenname
: Stringvisible
: Booleantag
: Integer
LowLevelFEM.showStrainResults
— MethodFEM.showStrainResults(problem, E, comp; t=..., name=..., visible=..., smooth=...)
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), 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
: ProblemE
: TensorFieldcomp
: Symbolt
: Vector{Float64}name
: Stringvisible
: Booleansmooth
: Booleantag
: Integer
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, :seqv), 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
: TensorFieldcomp
: Symbolt
: Vector{Float64}name
: Stringvisible
: Booleansmooth
: Booleantag
: Integer
LowLevelFEM.smallestEigenValue
— MethodFEM.smallestEigenValue(K, M)
Solves the largest eigenvalue for a transient problem given by stiffness (heat conduction) matrix K
and the mass (heat capacity) matrix M
(C
).
Return: λₘₐₓ
Types:
K
: SparseMatrixM
: SparseMatrixλₘₐₓ
: Float64
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
: SparseMatrixM
: SparseMatrixΔt
: Float64
LowLevelFEM.solveBuckling
— MethodFEM.solveBuckling(problem, loads, constraints; n=6)
Solves the multipliers for the first n
critical forces and the corresponding buckling shapes for the instability of the problem
, when loads
and constraints
are applied. Result can be presented by showBucklingResults
function. loads
and constraints
can be defined by load
and displacementConstraint
functions, respectively.
Return: buckling
Types:
problem
: Problemloads
: Vector{tuples}constraints
: Vector{tuples}n
: Int64buckling
: Eigen
LowLevelFEM.solveBucklingModes
— MethodFEM.solveBucklingModes(K, Knl; n=6)
Solves the critical force multipliers and buckling mode shapes of a problem given by stiffness matrix K
and the nonlinear stiffness matrix Knl
. n
is the number of buckling modes to solve. Returns the struct of critical forces and buckling modes. Results can be presented by showBucklingResults
function.
Return: modes
Types:
K
: SparseMatrixKnl
: SparseMatrixn
: Int64modes
: Eigen
LowLevelFEM.solveDisplacement
— MethodFEM.solveDisplacement(problem, load, supp)
Solves the displacement vector q
of problem
with loads load
and supports supp
.
Return: q
Types:
problem
: Problemload
: Vector{Tuple}supp
: Vector{Tuple}q
: Vector{Float64}
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
: SparseMatrixf
: Vector{Float64}q
: Vector{Float64}
LowLevelFEM.solveEigenModes
— MethodFEM.solveEigenModes(K, M; n=6, fₘᵢₙ=1.01)
Solves the eigen frequencies and mode shapes of a problem given by stiffness matrix K
and the mass matrix M
. n
is the number of eigenfrequencies to solve, and solves the eigenfrequencies greater than fₘᵢₙ
. Returns the struct of eigenfrequencies and eigen modes. Results can be presented by showModalResults
function.
Return: modes
Types:
K
: SparseMatrixM
: SparseMatrixn
: Int64fₘᵢₙ
: Float64modes
: Eigen
LowLevelFEM.solveHeatFlux
— MethodFEM.solveHeatFlux(problem, T; DoFResults=false)
Solves the heat flux field q
from temperature vector T
. heat flux is given per elements, so it usually contains jumps at the boundary of elements. Details of mesh is available in problem
. If DoFResults
is true, q
is a matrix with nodal results. In this case showDoFResults
can be used to show the results (otherwise showHeatFluxResults
or showElementResults
).
Return: q
Types:
problem
: ProblemT
: Vector{Float64}q
: VectorField or Matrix{Float}
LowLevelFEM.solveModalAnalysis
— MethodFEM.solveModalAnalysis(problem; constraints=[]; loads=[], n=6)
Solves the first n
eigenfrequencies and the corresponding mode shapes for the problem
, when loads
and constraints
are applied. loads
and contraints
are optional. Result can be presented by showModalResults
function. loads
and constraints
can be defined by load
and displacementConstraint
functions, respectively. If loads
are given, it solves the eigenfrequencies of a prestressed structure.
Return: modes
Types:
problem
: Problemloads
: Vector{tuples}constraints
: Vector{tuples}n
: Int64modes
: Eigen
LowLevelFEM.solveStrain
— MethodFEM.solveStrain(problem, q; DoFResults=false)
Solves the strain field E
from displacement vector q
. Strain field is given per elements, so it usually contains jumps at the boundaries of elements. Details of mesh is available in problem
. If DoFResults
is true, E
is a matrix with nodal results. In this case showDoFResults
can be used to show the results (otherwise showStrainResults
or showElementResults
).
Return: E
Types:
problem
: Problemq
: Vector{Float64}E
: TensorField or Matrix{Float64}
LowLevelFEM.solveStress
— MethodFEM.solveStress(problem, q; T=..., T₀=..., DoFResults=false)
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
. If DoFResults
is true, S
is a matrix with nodal results. In this case showDoFResults
can be used to show the results (otherwise showStressResults
or showElementResults
). If the T
temperature field (and T₀
initial temperature field if it differs from zero) is given, the function solves also the thermal stresses.
Return: S
Types:
problem
: Problemq
: Vector{Float64}T
: Vector{Float64}T₀
: Vector{Float64}S
: TensorField or Matrix{Float64}
LowLevelFEM.solveTemperature
— MethodFEM.solveTemperature(problem, flux, temp, heatconv)
Solves the temperature vector T
of problem
with given heat flux flux
, temperature temp
and heat convection heatconv
.
Return: T
Types:
problem
: Problemflux
: Vector{Tuple}temp
: Vector{Tuple}heatconv
: Vector{Tuple}T
: Vector{Float64}
LowLevelFEM.solveTemperature
— MethodFEM.solveTemperature(problem, flux, temp)
Solves the temperature vector T
of problem
with given heat flux flux
and temperature temp
.
Return: T
Types:
problem
: Problemflux
: Vector{Tuple}temp
: Vector{Tuple}T
: Vector{Float64}
LowLevelFEM.solveTemperature
— MethodFEM.solveTemperature(K, q)
Solves the equation K*T=q for the temperature vector T
. K
is the heat conduction matrix, q
is the heat flux vector.
Return: T
Types:
K
: SparseMatrixq
: Vector{Float64}T
: Vector{Float64}
LowLevelFEM.stiffnessMatrix
— MethodFEM.stiffnessMatrix(problem)
Solves the stiffness matrix of the problem
.
Return: stiffMat
Types:
problem
: ProblemstiffMat
: SparseMatrix
LowLevelFEM.temperatureConstraint
— MethodFEM.temperatureConstraint(name; T=...)
Gives the temperature constraints on name
physical group. T
can be a constant value, or a function of x
, y
and z
. (E.g. fn(x,y,z)=5*(5-x)); FEM.temperatureConstraint("surf1", T=fn)
)
Return: Tuple{String, Float64 or Function, Float64 or Function, Float64 or Function}
Types:
name
: StringT
: Float64 or Function
LowLevelFEM.thermalLoadVector
— MethodFEM.thermalLoadVector(problem, T; T₀=...)
Solves the thermal load vector from a temperature field T
for problem problem
. T₀
is the initial temperature field.
Return: thermLoadVec
Types:
problem
: ProblemT
: Vector{Float64}T₀
: Vector{Float64}thermLoadVec
: Vector{Float64}
LowLevelFEM.vectorField
— MethodFEM.vectorField(problem, entity)
Defines a vector field from entity
, which is a tuple of name
of physical group and prescribed values or functions. Mesh details are in problem
.
Return: Vector{Float64}
Types:
problem
: Problementity
: Vector{Tuple{String, Float64,...}}
Index
LowLevelFEM.CoordinateSystem
LowLevelFEM.Eigen
LowLevelFEM.Material
LowLevelFEM.Problem
LowLevelFEM.TensorField
LowLevelFEM.VectorField
LowLevelFEM.CDM
LowLevelFEM.CDMaccuracyAnalysis
LowLevelFEM.FDM
LowLevelFEM.HHT
LowLevelFEM.HHTaccuracyAnalysis
LowLevelFEM.applyBoundaryConditions
LowLevelFEM.applyBoundaryConditions!
LowLevelFEM.applyBoundaryConditions!
LowLevelFEM.applyBoundaryConditions!
LowLevelFEM.applyBoundaryConditions!
LowLevelFEM.applyElasticSupport!
LowLevelFEM.applyHeatConvection!
LowLevelFEM.constrainedDoFs
LowLevelFEM.dampingMatrix
LowLevelFEM.displacementConstraint
LowLevelFEM.elasticSupport
LowLevelFEM.elasticSupportMatrix
LowLevelFEM.elementsToNodes
LowLevelFEM.field
LowLevelFEM.fieldError
LowLevelFEM.freeDoFs
LowLevelFEM.generateMesh
LowLevelFEM.getTagForPhysicalName
LowLevelFEM.heatCapacityMatrix
LowLevelFEM.heatConductionMatrix
LowLevelFEM.heatConvection
LowLevelFEM.heatConvectionMatrix
LowLevelFEM.heatConvectionVector
LowLevelFEM.heatFlux
LowLevelFEM.heatFluxVector
LowLevelFEM.heatSource
LowLevelFEM.heatSourceVector
LowLevelFEM.initialDisplacement
LowLevelFEM.initialDisplacement!
LowLevelFEM.initialTemperature
LowLevelFEM.initialTemperature!
LowLevelFEM.initialVelocity
LowLevelFEM.initialVelocity!
LowLevelFEM.largestEigenValue
LowLevelFEM.largestPeriodTime
LowLevelFEM.latentHeatMatrix
LowLevelFEM.load
LowLevelFEM.loadVector
LowLevelFEM.massMatrix
LowLevelFEM.material
LowLevelFEM.nodalAcceleration!
LowLevelFEM.nodalForce!
LowLevelFEM.nonLinearStiffnessMatrix
LowLevelFEM.openPostProcessor
LowLevelFEM.openPreProcessor
LowLevelFEM.plotOnPath
LowLevelFEM.resultant
LowLevelFEM.rotateNodes
LowLevelFEM.scalarField
LowLevelFEM.setParameter
LowLevelFEM.showBucklingResults
LowLevelFEM.showDoFResults
LowLevelFEM.showElementResults
LowLevelFEM.showHeatFluxResults
LowLevelFEM.showModalResults
LowLevelFEM.showStrainResults
LowLevelFEM.showStressResults
LowLevelFEM.smallestEigenValue
LowLevelFEM.smallestPeriodTime
LowLevelFEM.solveBuckling
LowLevelFEM.solveBucklingModes
LowLevelFEM.solveDisplacement
LowLevelFEM.solveDisplacement
LowLevelFEM.solveEigenModes
LowLevelFEM.solveHeatFlux
LowLevelFEM.solveModalAnalysis
LowLevelFEM.solveStrain
LowLevelFEM.solveStress
LowLevelFEM.solveTemperature
LowLevelFEM.solveTemperature
LowLevelFEM.solveTemperature
LowLevelFEM.stiffnessMatrix
LowLevelFEM.temperatureConstraint
LowLevelFEM.thermalLoadVector
LowLevelFEM.vectorField
- 4Serfőző, D., Pere, B.: An effective reduction method with Caughey damping for spurious oscillations in dynamic problems, preprint, https://doi.org/10.21203/rs.3.rs-3930320/v1
- 5Bathe, K. J.: Finite element procedures, Wiley, 1983, https://doi.org/10.1002/nag.1610070412
- 1Hilber, Hans M., Thomas JR Hughes, and Robert L. Taylor. "Improved numerical dissipation for time integration algorithms in structural dynamics." Earthquake Engineering & Structural Dynamics 5.3 (1977): 283-292.
- 2Belytschko, Ted, and Thomas JR, Hughes: "Computational methods for transient analysis", North-Holland, (1983).
- 3Serfőző, D., Pere, B.: A method to accurately define arbitrary algorithmic damping character as viscous damping. Arch Appl Mech 93, 3581–3595 (2023). https://doi.org/10.1007/s00419-023-02454-9