FiniteElementAnalysis.jl

Documentation for FiniteElementAnalysis.

Functions

FiniteElementAnalysis.axialstressesMethod
axialstresses(
    elements::Dict{String, LineElements},
    nodedisplacements::Vector{Float64},
    dims=1
)::Vector{Float64}
axialstresses(
    elements::Dict{String, BeamElements},
    bendingmoments::Vector{Float64},
    distancesfromneutralaxis::Union{<:Real, Vector{<:Real}},
    dims=2
)::Vector{Float64}

Determine the axial stresses present in element-type system.

If line elements, of two-force members in their local coordinate systems. If beam elements, of beam elements which assumes no direct, axial forces.

See also: prepareelements_line, prepareelements_beam, solve, and transformationmatrixstar.

Arguments

  • elements: the prepared dictionary according to element type.
  • nodedisplacements: the vector for displacement of each node in the global coordinate system.
  • bendingmoments: the vector of reaction moments at each node from solve.F.
  • distancesfromneutralaxis: the maximum distance(s) of nodes from neutral axis.
  • dims=1: the number of dimensions being analyzed.

Examples

A, E, L = 1, 10e6, 100          # [in², psi, in]
angles = [120, 0, 210]          # [°]
L = L ./ abs.(cosd.(angles))    # [in]
elements = Dict("1"=>(1, 2), "2"=>(1, 3), "3"=>(1, 4))
preparedelements = prepareelements_line(elements, A, E, L, angles, dims=2)
nodeboundaryconditions = [2, 3, 4]
F_applied = [(1, (1e3, 1e3))]   # (node, (force_x [lb], force_y [lb]))
U = solve(preparedelements, nodeboundaryconditions, F_applied, dims=2).U
axialstresses(preparedelements, U, dims=2)
A, E, I, L = 0.004887968, 200e9, 102e6 * (1/1e3)^4, 6               # [m², Pa, m⁴, m]
Q, c, t_w = 323473.474 * (1/1e3)^3, 176.5e-3, 6.48e-3               # [m³, m, m]
W(x) = -10e3(x)                                                     # [N/m]
F_applied = ([(W, 1, 2, 3, 4)])
numberofelements = 10
L /= numberofelements
x = 0:L:18
elements, nodeboundaryconditions = Dict{String, Tuple{Int, Int}}(), []
for node ∈ 1:1:numberofelements+1
    if node == 1
        push!(nodeboundaryconditions, node)
    else
        elements["$(node-1)"] = (node-1, node)
        if x[node]%6 == 0
            push!(nodeboundaryconditions, (node, (0, Inf)))
        end
    end
end
preparedelements = prepareelements_beam(elements, E, I, L)
sol = solve(preparedelements, nodeboundaryconditions, F_applied)    # [m]
V = [sol.F[2i - 1] for i ∈ 1:1:length(sol.F)÷2]                     # [N]
M = [sol.F[2i] for i ∈ 1:1:length(sol.F)÷2]                         # [N-m]
maximum(axialstresses(preparedelements, M, c))                      # [Pa]
source
FiniteElementAnalysis.globalstiffnessmatrixMethod
globalstiffnessmatrix(
    elements::Dict{String, LineElements},
    dims::Integer=1
)::Matrix{Float64}
globalstiffnessmatrix(
    elements::Dict{String, BeamElements},
    dims::Integer=2
)::Matrix{Float64}
globalstiffnessmatrix(
    elements::Dict{String, TriangularElements},
    dims::Integer=2
)::Matrix{Float64}

Assembles the global stiffness matrix from the dictionary, elements.

See also: prepareelements_line, prepareelements_beam, prepareelements_triangular, and transformationmatrix.

Arguments

  • elements: entries infer localstiffnessmatrices and nodeconnections for each element in local coordinate system transformed by angle.
  • dims: the number of dimensions being analyzed.

Examples

A = [0.75, 0.5, 1]   # [in²]
E = 10.6e6           # [psi]
L = [12, 9, 8]       # [in]
elements = Dict("1"=>(1, 2), "2"=>(2, 3), "3"=>(3, 4))
globalstiffnessmatrix(prepareelements_line(elements, A, E, L))
E, I, L = 30e6, 510, 60 # [psi, in⁴, in]
W(x) = -12\1e3(x)      # [lb/in]
F_applied = ([(W, 1, 2, 3)])
elements = Dict("1"=>(1, 2), "2"=>(2, 3))
globalstiffnessmatrix(prepareelements_beam(elements, E, I, L))
# problem 6-3.1
E, ν, t = 30e6, 0.3, 1  # [psi, n/a, in]
elements = Dict(
    "1"=>(i=1, j=3, m=2),
    "2"=>(i=1, j=4, m=3)
)
coords = Dict(
    "1"=>(0, 0),
    "2"=>(0, 10),
    "3"=>(20, 10),
    "4"=>(20, 0)
)
elements = prepareelements_triangular(elements, coords, E, ν, t)
globalstiffnessmatrix(elements) .* 0.91/375e3
source
FiniteElementAnalysis.prepareelements_beamFunction
prepareelements_beam(
    elements,
    elementmoduli,
    elementareasofinertia,
    elementlengths,
    nodeangles=0,
    isdegrees=true;
    dims=2,
    stiffnesscoefficients=nothing
)::Dict{String, BeamElements}

Return the dictionary with properties specific to those beam elements between node pairs.

See also: globalstiffnessmatrix, solve, and axialstresses.

Arguments

  • elements::Dict{String, Tuple{Int, Int}}: the dictionary of "element"=>(node 1, node 2) for the element between node pairs.
  • elementmoduli::Union{<:Real, Vector{<:Real}}: the elastic modulus of each element.
  • elementareasofinertia::Union{<:Real, Vector{<:Real}}: the cross-sectional area of each element.
  • elementlengths::Union{<:Real, Vector{<:Real}}: the length of each element.
  • nodeangles::Union{<:Real, Vector{<:Real}}=0: the angle between the global and local coordinate systems of each element.
  • isdegrees::Bool=true: control value whether elements of nodeangles is degrees (true) or radians (false).
  • dims::Integer=2: the number of dimensions being analyzed.
  • stiffnesscoefficients::Union{<:Real, Vector{<:Real}}: the stiffness coefficient, $k$ associated with each element.

Examples

E, I, L = 30e6, 200, 20     # [psi, in⁴, ft]
elements = Dict("1"=>(1, 2),"2"=>(2, 3))
prepareelements_beam(elements, E, I, L)
source
FiniteElementAnalysis.prepareelements_lineFunction
prepareelements_line(
    elements,
    elementareas,
    elementmoduli,
    elementlengths,
    nodeangles=0,
    isdegrees=true;
    dims=1,
    stiffnesscoefficients=nothing
)::Dict{String, LineElements}

Return the dictionary with properties specific to those line elements between node pairs.

See also: globalstiffnessmatrix, solve, and axialstresses.

Arguments

  • elements::Dict{String, Tuple{Int, Int}}: the dictionary of "element"=>(node 1, node 2) for the element between node pairs.
  • elementareas::Union{<:Real, Vector{<:Real}}: the cross-sectional area(s) of elements.
  • elementmoduli::Union{<:Real, Vector{<:Real}}: the elastic modul(us/i) of elements.
  • elementlengths::Union{<:Real, Vector{<:Real}}: the length(s) of elements.
  • nodeangles::Union{<:Real, Vector{<:Real}}=0: the angle between the global and local coordinate systems of elements.
  • isdegrees::Bool=true: control value whether elements of nodeangles is degrees (true) or radians (false).
  • dims::Integer=1: the number of dimensions being analyzed.
  • stiffnesscoefficients::Union{<:Real, Vector{<:Real}}: the stiffness coefficient, $k$ associated with each element.

Examples

A = [0.75, 0.5, 1]   # [in²]
E = 10.6e6           # [psi]
L = [12, 9, 8]       # [in]
elements = Dict("1"=>(1, 2), "2"=>(2, 3), "3"=>(3, 4))
prepareelements_line(elements, A, E, L)
source
FiniteElementAnalysis.prepareelements_triangularFunction
prepareelements_triangular(
    elements,
    coordinates,
    elementmoduli,
    poissonsratios,
    elementthicknesses,
    nodeangles=0,
    isdegrees=true;
    dims=2,
    planestressstrain="stress",
    elementareas=nothing
)::Dict{String, TriangularElements}

Return the dictionary with properties specific to those linear, triangular elements between node triplets.

See also: globalstiffnessmatrix and solve.

Arguments

  • elements::Dict{String, NamedTuple{(:i, :j, :m), Tuple{Int, Int, Int}}}: the dictionary of "element"=>(node 1, node 2, node 3) for the element between node triplets.
  • coordinates::Dict{String, <:Tuple{<:Real, <:Real}}: the dictionary of "node"=>(x, y) coordinate pairs for each node.
  • elementmoduli::Union{<:Real, Vector{<:Real}}: the elastic modulus of each element.
  • poissonsratios::Union{<:Real, Vector{<:Real}}: the Poisson's ratio for each element.
  • elementthicknesses::Union{<:Real, Vector{<:Real}}: the plane thickness of each element.
  • nodeangles::Union{<:Real, Vector{<:Real}}=0: the angle between the global and local coordinate systems of each element.
  • isdegrees::Bool=true: control value whether elements of nodeangles is degrees (true) or radians (false).
  • dims::Integer=2: the number of dimensions being analyzed.
  • planestressstrain::String="stress": if 2D analysis, specify whether plane "stress" or "strain" assumption applies.
  • elementareas::Union{Nothing, <:Real, Vector{<:Real}}=nothing: the area of each element, if already known.

Examples

E, ν, t = 30e6, 0.3, 1  # [psi, n/a, in]
elements = Dict(
    "1"=>(i=1, j=3, m=2),
    "2"=>(i=1, j=4, m=3)
)
coords = Dict(
    "1"=>(0, 0),
    "2"=>(0, 10),
    "3"=>(20, 10),
    "4"=>(20, 0)
)
prepareelements_triangular(elements, coords, E, ν, t)
source
FiniteElementAnalysis.reducedglobalstiffnessmatrixMethod
reducedglobalstiffnessmatrix(
    K::Matrix{Float64},
    nodes::Union{Vector{Integer},
    Vector{Tuple{Integer, Integer}}},
    dims::Integer=1
)::Matrix{Float64}

Eliminate the rows and columns of K according to nodes which implies the displacement of the appropriate node is known: e.g. is fixed or a prescribed displacement.

See also: globalstiffnessmatrix.

Examples

A = [0.75, 0.5, 1]   # [in²]
E = 10.6e6           # [psi]
L = [12, 9, 8]       # [in]
elements = Dict("1"=>(1, 2), "2"=>(2, 3), "3"=>(3, 4))
K = globalstiffnessmatrix(prepareelements_line(elements, A, E, L))
skipnodes = Vector{Integer}([1, 4])
reducedglobalstiffnessmatrix(K, skipnodes)
source
FiniteElementAnalysis.shearstressesMethod
shearstresses(
    elements,
    shearforces,
    firstmomentofareas,
    bendingthicknesses,
    dims=2
)::Vector{Float64}

Determine the axial stresses present in beam element system.

See also: prepareelements_beam, solve, and transformationmatrixstar.

Arguments

  • elements::Dict{String, BeamElements}: the prepared dictionary according to element type.
  • shearforces::Vector{Float64}: the vector for reaction, shear forces at each node in the global coordinate system.
  • firstmomentofarea::Union{<:Real, Vector{<:Real}}: the first moment of area(s), $Q$ above the neutral axis.
  • bendingthicknesses::Union{<:Real, Vector{<:Real}}: the web thickness of I-beam.
  • dims::Integer=2: the number of dimensions being analyzed.

Examples

A, E, I, L = 0.004887968, 200e9, 102e6 * (1/1e3)^4, 6               # [m², Pa, m⁴, m]
Q, c, t_w = 323473.474 * (1/1e3)^3, 176.5e-3, 6.48e-3               # [m³, m, m]
W(x) = -10e3(x)                                                     # [N/m]
F_applied = ([(W, 1, 2, 3, 4)])
numberofelements = 10
L /= numberofelements
x = 0:L:18
elements, nodeboundaryconditions = Dict{String, Tuple{Int, Int}}(), []
for node ∈ 1:1:numberofelements+1
    if node == 1
        push!(nodeboundaryconditions, node)
    else
        elements["$(node-1)"] = (node-1, node)
        if x[node]%6 == 0
            push!(nodeboundaryconditions, (node, (0, Inf)))
        end
    end
end
preparedelements = prepareelements_beam(elements, E, I, L)
sol = solve(preparedelements, nodeboundaryconditions, F_applied)    # [m]
V = [sol.F[2i - 1] for i ∈ 1:1:length(sol.F)÷2]                     # [N]
M = [sol.F[2i] for i ∈ 1:1:length(sol.F)÷2]                         # [N-m]
maximum(shearstresses(preparedelements, V, Q, t_w))                 # [Pa]
source
FiniteElementAnalysis.solveMethod
solve(
    elements::Dict{String, LineElements},
    nodeboundaryconditions,
    nodeforces::Vector{<:Tuple{Integer, Vararg{<:Real}}},
    dims::Integer=1
)::FEA
solve(
    elements::Dict{String, BeamElements},
    nodeboundaryconditions,
    nodeforces::Vector{<:Tuple{Integer, Vararg{<:Real}}},
    dims::Integer=2
)::FEA

Perform Finite Element Analysis (FEA) with entries of elements and applied boundary conditions by the Direct Stiffness Method according to Hooke's Law of linear-elastic deformation.

If line elements, stiffness coefficient defined by: $k ≡ AE/L$ (typical of trusses with two-force members). If beam elements, stiffness coefficient defined by: $k ≡ EI/L³$.

See also: prepareelements_line, prepareelements_beam, globalstiffnessmatrix, and reducedglobalstiffnessmatrix.

Arguments

  • elements: the prepared dictionary according to element type.
  • nodeboundaryconditions: the tuple of nodes that either experience zero displacement (tuple item of integer for node) or prescribed displacement constrained by the stiffness of a certain element (tuple item of node integer and tuple of real displacement and integer element).
  • nodeforces::Vector{<:Tuple{Integer, Vararg{<:Real}}}: the vector of tuple pairs for nodes experiencing some applied force.
  • dims::Integer=1: the number of dimensions being analyzed.
  • alg=nothing: algorithm of choice to solve linear system of equations according to LinearSolve.jl.

Examples

LineElements

1D

nodeboundaryconditions:= tuple of nodes with zero displacement.

A, E, L = 4, 30e6, 30               # [in², psi, in]
elements = Dict("1"=>(1, 2), "2"=>(2, 3), "3"=>(3, 4))
preparedelements = prepareelements_line(elements, A, E, L)
nodeboundaryconditions = (1, 4)
F_applied = [(2, 5e3), (3, -10e3)]  # (node, force [lb])
solve(preparedelements, nodeboundaryconditions, F_applied).U

nodeboundaryconditions:= tuple of nodes with zero and prescribed displacements.

The tuple items are identified in the following way:

  1. Node (1) with zero displacement.
  2. Node (3) has a prescribed displacement of 25e-3 constrained by the deformation of element 2.
A, E, L = 4e-4, 210e9, 2                    # [m², Pa, m]
elements = Dict("1"=>(1, 2), "2"=>(2, 3))
preparedelements = prepareelements_line(elements, A, E, L)
nodeboundaryconditions = [1, (3, 25e-3, 2)] # (node, displacement [m], element)
F_applied = [(2, -5e3)]                     # (node, force [N])
solve(preparedelements, nodeboundaryconditions, F_applied).U

2D

A, E, L = 1, 10e6, 100          # [in², psi, in]
angles = [120, 0, 210]          # [°]
L = L ./ abs.(cosd.(angles))    # [in]
elements = Dict("1"=>(1, 2), "2"=>(1, 3), "3"=>(1, 4))
preparedelements = prepareelements_line(elements, A, E, L, angles, dims=2)
nodeboundaryconditions = [2, 3, 4]
F_applied = [(1, (1e3, 1e3))]   # (node, (force_x [lb], force_y [lb]))
solve(preparedelements, nodeboundaryconditions, F_applied, dims=2).U

TriangularElements

2D

# problem 6-3.1
E, ν, t = 30e6, 0.3, 1  # [psi, n/a, in]
T(x) = (0, -5e3(x))     # [lb/in]
F_applied = [(3, T(10) ./ 2), (4, T(10) ./ 2)]
elements = Dict(
    "1"=>(i=1, j=3, m=2),
    "2"=>(i=1, j=4, m=3)
)
coords = Dict(
    "1"=>(0, 0),
    "2"=>(0, 10),
    "3"=>(20, 10),
    "4"=>(20, 0)
)
nodeboundaryconditions = [1, 2]
preparedelements = prepareelements_triangular(elements, coords, E, ν, t)
solve(preparedelements, nodeboundaryconditions, F_applied).U
source

Index