Asap.BridgeElementType
BridgeElement(elementStart::Element, posStart::Float64, elementEnd::Element, posEnd::Float64, section::Section, id = :element; release = :fixedfixed)

Create a bridge element between two frame elements. Connects from elementStart at a position elementStart.length * posStart away from elementStart.nodeStart.position to elementEnd at elementEnd.length * posEnd away from elementEnd.nodeStart.position. IE posStart, posEnd ∈ ]0, 1[

Asap.ElementType
Element(nodes::Vector{Node}, nodeIndex::Vector{Int64}, section::Section, id::Symbol = nothing; release = :fixedfixed)
Element(nodeStart::Node, nodeEnd::Node, section::Section, id = :element; release = :fixedfixed)

Create a frame element with an optional id tag.

Example

julia> Element(nodes, [1,2], sec)
julia> Element(node1, node2, sec, :groundfloor_element)

Optional argument release

This property enables decoupling of nodal DOFs with respect to the end of the element.

Available releases:

  • :fixedfixed (default) - all DOFs are tied to nodes
  • :fixedfree - rotational DOFs are released at end node
  • :freefixed - rotational DOFs are released at start node
  • :freefree - all rotational DOFs are released (truss element)
  • :joist - all rotational DOFs except torsion are released

Example

julia> Element(nodes, [1,2], sec; release = :fixedfree)
Asap.FDMelementType
FDMelement(points::Vector{FDMnode}, iStart::Int64, iEnd::Int64, q::Real, id = :element)
FDMelement(points::Vector{FDMnode}, indices::Vector{Int64}, q::Real, id = :element)
FDMelement(pointStart::FDMnode, pointEnd::FDMnode, q::Real, id = :element)

An element in a Force Density Method network that connects two nodes with force density q and an optional identifier id::Symbol

Asap.FDMloadType
FDMload(points::Vector{FDMnode}, i::Int64, force::Vector{<:Real})
FDMload(point::FDMnode, force::Vector{<:Real})

An external force applied to a FDM node: load = [Px, Py, Pz]

Asap.FDMnodeType
FDMnode(x::Real, y::Real, z::Real, dof::Bool, id = :node)
FDMnode(pos::Vector{<:Real}, dof::Bool, id = :node)

A node in an FDM network defined by its spatial position [x, y, z] and degree of freedom (free = true, fixed = false)

Fields

  • position::Vector{Float64}: the [x,y,z] position of the node
  • dof::Bool: true = free; false = fixed
  • id::Symbol: identifier (optional)
  • nodeID::Integer: global indentifier (internal)
  • reaction::Vector{Float64}: the [x,y,z] reaction forces if node is fixed
Asap.GravityLoadType
GravityLoad(element::Element, factor::Float64)

A gravity load (negative global Z) applied along a member.

Generates distributed load w = element.section.A * element.section.ρ * factor, where factor should be the appropriate acceleration due to gravity.

Asap.LineLoadType
LineLoad(element::Element, value::Vector{Float64})

A distributed line load [wx, wy, wz] in (force/length) applied along an element in the global coordinate system.

Asap.MaterialType
Material(E::Float64, G::Float64, ρ::Float64, ν::Float64)

Define a structural material.

Fields

  • E Modulus of Elasticity [Force/distance²]
  • G Shear Modulus [Force/distance²]
  • ρ Density [Mass/distance³]
  • ν Poisson's Ratio [unitless]
Asap.ModelType
Model(nodes::Vector{Node}, elements::Vector{Element}, loads::Vector{AbstractLoad})

Create a complete structural model ready for analysis.

Asap.NetworkType
Network(nodes::Vector{FDMnode}, elements::Vector{FDMelement}, loads::Vector{FDMload})

A FDM network defined by a collection of nodes, elements, and loads.

Fields

  • nodes::Vector{FDMnode} collection of nodes in network
  • elements::Vector{FDMelement} collection of elements in network
  • loads::Vector{FDMload} collection of external loads in network
  • q::Vector{<:Real} vector of elemental force densities [n_e × 1]
  • Q::SparseMatrixCSC{Float64, Int64} diagonalized sparse matrix of q [ne × ne]
  • C::SparseMatrixCSC{Int64, Int64} element/node connectivity matrix [ne × nn]
  • N::Vector{Int64} vector of free node indices
  • F::Vector{Int64} vector of fixed node indices
  • Cn::SparseMatrixCSC{Int64, Int64} connectivity matrix of free nodes = C[:, N]
  • Cf::SparseMatrixCSC{Int64, Int64} connectivity matrix of fixed nofes = C[:. F]
  • P::Matrix{<:Real} external load matrix [n_n × 3]
  • Pn::Matrix{<:Real} loads applied to free nodes = P[N, :]
  • xyz::Matrix{<:Real} positions of all nodes [n_n × 3]
Asap.NodeType
Node(position::Vector{Float64}, dofs::Vector{Bool}, id::Symbol = nothing)

Instantiate a 6 DOF node with given position and fixities. Optional symbol identifier, id.

Example

julia> Node([4.3, 2.2, 10.4], [true, true, false, true, false, false])
Node([4.3, 2.2, 10.4], Bool[1, 1, 0, 1, 0, 0], #undef, #undef, #undef, nothing)

Node(position::Vector{Float64}, fixity::Symbol, id::Symbol = nothing)

Instantiate a 6 DOF node with given position and common boundary type. Optional symbol identifier, id.

Available boundary conditions:

  • :free
  • :fixed
  • :pinned
  • :(x/y/z)free
  • :(x/y/z)fixed

Example

julia> Node([4.3, 2.2, 10.4], :zfixed)
Node([4.3, 2.2, 10.4], Bool[1, 1, 0, 1, 1, 1], #undef, #undef, #undef, nothing)
Asap.NodeForceType
NodeForce(node::AbstractNode, value::Vector{Float64})
NodeForce(nodes::Vector{<:AbstractNode}, index::Integer, value::Vector{Float64})

A force vector [Fx, Fy, Fz] in the global coordinate system applied to a node.

Asap.NodeMomentType
NodeMoment(node::Node, value::Vector{Float64})

A moment vector [Mx, My, Mz] in the global coordinate system applied to a node with rotational DOFs.

Asap.PointLoadType
PointLoad(element::Element, position::Float64, value::Vector{Float64})

A point load [Px, Py, Pz] applied in the global coordinate system at a distance position × element.length from the starting node.

Asap.SectionType
Section(A::Float64, E::Float64, G::Float64, Ix::Float64, Iy::Float64, J::Float64, ρ::Float64 = 1.)
Section(mat::Material, A::Float64, Ix::Float64, Iy::Float64, J::Float64)

A cross section assigned to an element.

Fields

  • A Area [Distance²]
  • E Modulus of Elasticity [Force/Distance²]
  • Ix Nominal strong moment of inertia [Distance⁴]
  • Iy Nominal weak moment of inertia [Distance⁴]
  • J Torsional constant [Distance⁴]
  • ρ=1 Density [Mass/Distance³]
Asap.TrussElementType
TrussElement(nodes::Vector{TrussNode}, nodeIndex::Vector{Int64}, section::AbstractSection, id = :element)
TrussElement(nodeStart::TrussNode, nodeEnd::TrussNode, section::AbstractSection, id = :element)

Create a truss element.

Example

julia> TrussElement(nt, [1,2], sec)
TrussElement(Section(794.0, 200000.0, 77000.0, 737000.0, 737000.0, 1.47e6, 1.0), [1, 2], TrussNode([0.8879630592102802, 0.6713937498337156, 0.617463764682365], Bool[0, 0, 0], #undef, [0.0, 0.0, 0.0], [0.0, 0.0, 0.0], nothing), TrussNode([0.16046742214916832, 0.15869760269854827, 0.6247762072043447], Bool[1, 1, 1], #undef, [0.0, 0.0, 0.0], [0.0, 0.0, 0.0], nothing), #undef, 0.8900341077991537, #undef, #undef, #undef, 1.5707963267948966, #undef, nothing)
Asap.TrussModelType
TrussModel(nodes::Vector{TrussNode}, elements::Vector{TrussElement}, loads::Vector{NodeForce})

Create a complete structural model ready for analysis.

Asap.TrussNodeType
TrussNode(position::Vector{Float64}, dofs::Vector{Bool}, id::Symbol = nothing)

Instantiate a 3 DOF node with given position and fixities. Optional symbol identifier, id.

Example

julia> TrussNode([1., 1., 56.], [false, true, true])
TrussNode([1.0, 1.0, 56.0], Bool[0, 1, 1], #undef, [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], [0.0, 0.0, 0.0, 0.0, 0.0, 0.0], nothing)

TrussNode(position::Vector{Float64}, fixity::Symbol, id::Symbol = nothing)

Instantiate a 3 DOF node with given position and common boundary type. Optional symbol identifier, id.

Available boundary conditions:

  • :free
  • :pinned
  • :(x/y/z)free
  • :(x/y/z)fixed

Example

julia> TrussNode([1., 1., 56.], :pinned)
TrussNode([1.0, 1.0, 56.0], Bool[0, 0, 0], #undef, [0.0, 0.0, 0.0], [0.0, 0.0, 0.0], nothing)
Asap.TrussSectionType
TrussSection(A::Float64, E::Float64)
TrussSection(mat::Material, A::Float64)

A cross section assigned to a truss element.

Fields

  • A Area [Distance²]
  • E Modulus of Elasticity [Force/Distance²]
  • ρ Density [Mass/Distance³]
Asap.NF!Method
NF!(network::Network)

Extract fixed/free DOF indices from a vector of nodes.

Asap.NFMethod

Extracts the indices of free (N) and fixed (F) nodes

Asap.NFMethod
NF(points::Vector{FDMnode})

Return fixed/free DOF indices from a vector of nodes.

Asap.RMethod
R(element::Element; tol = 1e-4)

Get the [12 × 12] transformation matrix for a given element.

tol defines the threshold criteria for triggering special R computation for elements aligned with the global Y axis

Asap.RMethod
R(element::TrussElement)

Get the [2 × 6] transformation matrix for a truss element.

Asap.RMethod
R(xvec::Vector{Float64}, Ψ; tol = 1e-4)

Return the [12 × 12] transformation matrix given a local x vector xvec

tol defines the threshold criteria for triggering special R computation for elements aligned with the global Y axis

Asap.RMethod
R(xvec::Vector{Float64})

Get the [2 × 6] transformation matrix given a local x vector xvec

Asap.axial_forceMethod
axial_force(element::Element)

Extract the axial force of an element

Asap.axial_forceMethod
axial_force(element::TrussElement)

Extract the axial force of an element

Asap.branch_matrix!Method
branch_matrix!(network::Network)

Populate the [nelements × nnodes] sparse connectivity matrix, network.C

Asap.branch_matrixMethod
branch_matrix(elements::Vector{FDMelement}, points::Vector{FDMnode})

Return the [nelements × nnodes] sparse connectivity matrix, network.C

Asap.connectivityMethod
connectivity(model::AbstractModel)

Get the [nₑ × nₙ] sparse matrix where C[i, j] = -1 if element i starts at node j, and C[i,j] = 1 if element i ends at node j, and 0 otherwise.

Asap.create_FMethod
create_F(model::TrussModel, loads::Vector{NodeForce})

create load vector F = P

Asap.create_FMethod
create_F(model::TrussModel, loads::Vector{NodeForce})

create load vector F = P - Pf

Asap.create_S!Method
create_S!(model::Model)

Assemble the global stiffness matrix S.

Asap.create_S!Method
create_S!(model::TrussModel)

Assemble the global stiffness matrix S.

Asap.deconstruct_nodes!Method
deconstruct_nodes!(network::Network)

Populate the [n_nodes × 3] matrix of nodal positions, network.xyz

Asap.deconstruct_nodesMethod
deconstruct_nodes(points::Vector{FDMnode})

return the [n_nodes × 3] matrix of nodal positions, network.xyz

Asap.dofsMethod
dofs(points::Vector{FDMnode})

Return the fixed/free DOF information from a vector of nodes.

Asap.endpointsMethod
endpoints(element::AbstractElement)

Extract the start and end points as two vectors.

Asap.fixnode!Method
fixnode!(node::AbstractNode, fixity::Symbol)

Fix the DOFs of a node to a common boundary condition.

Arguments

  • node::AbstractNodenode to modify
  • fixity::Symbol boundary condition to apply. Available boundary conditions:
    • :free
    • :fixed
    • :pinned
    • :(x/y/z)free
    • :(x/y/z)fixed
Asap.force_densities!Method
force_densities!(network::Network)

Extract the force densities of elements in a network and assemble the global Q matrix, network.Q

Asap.force_densitiesMethod
force_densities(elements::Vector{FDMelement})

Extract the force densities of a vector of elements.

Asap.forcesMethod
forces(network::Network)

Return the axial forces in a network.

Asap.global_K!Method
global_K!(element::Element)

Populate the element stiffness matrix element.K in GCS.

Asap.global_KMethod
global_K(element::Element)

Return the element stiffness matrix in GCS.

Asap.initial_lengthsMethod
initial_lengths(network::Network, E::Real, A::Real)

Get the required unstressed length of elements given a material stiffness E and cross sectional area A.

Asap.initial_lengthsMethod
initial_lengths(network::Network, E::Vector{<:Real}, A::Vector{<:Real})

Get the required unstressed length of elements given an element-wise material stiffness vector E and cross sectional area vector A.

Asap.k_fixedfixedMethod
k_fixedfixed(element::Element)

[12 × 12] stiffness matrix for a beam element fully coupled at both ends.

Asap.k_fixedfreeMethod
k_fixedfree(element::Element)

[12 × 12] stiffness matrix for a beam element with rotational DOFs decoupled at the end node.

Asap.k_freefixedMethod
k_freefixed(element::Element)

[12 × 12] stiffness matrix for a beam element with rotational DOFs decoupled at the start node.

Asap.k_freefreeMethod
k_freefree(element::Element)

[12 × 12] stiffness matrix for a beam element with full decoupled rotational DOFs.

Asap.k_joistMethod
k_joist(element::Element)

[12 × 12] stiffness matrix for a beam element with only torsional DOFs coupled to nodes.

Asap.lcs!Method
lcs(element::AbstractElement, Ψ::Float64; tol = 1e-6)

Populate local coordinate system unit vectors of a given element and pitch angle Ψ: [localx, localy, local_z]

Asap.lcsMethod
lcs(xin::Vector{Float64}, Ψ::Float64; tol = 1e-4)

Get the local coordinate system unit vectors of a given local x axis and pitch angle Ψ: [localx, localy, local_z]

Asap.lcsMethod
lcs(element::AbstractElement, Ψ::Float64; tol = 1e-6)

Get the local coordinate system unit vectors of a given element and pitch angle Ψ: [localx, localy, local_z]

Asap.load_matrix!Method
load_matrix!(network::Network)

Generate the [n_nodes × 3] matrix of nodal forces, P

Asap.local_KMethod
local_K(element::Element)

Return the element stiffness matrix in LCS.

Asap.local_xMethod
local_x(element::AbstractElement; unit = true)

Get the local x vector of an element: element.nodeEnd.position - element.nodeStart.position. unit = true gives the normalized vector.

Asap.midpointMethod
midpoint(element::AbstractElement)

Extract the centerpoint of an element as a vector.

Asap.node_positionsMethod
node_positions(model::AbstractModel)

Generate the [nₙ × 3] node position matrix

Asap.planarize!Function
planarize!(model::AbstractModel, plane = :XY)

Fix all nodal DOFs to remain on plane = plane

Asap.planarize!Function
planarize!(nodes, plane = :XY)

Restrict the DOFs of a set of nodes to a given global plane. E.g., when analyzing a 2D structure.

Arguments

  • nodes::Vector{<:AbstractNode} vector of nodes to restrict
  • plane::Symbol = :XY plane to restrict DOFs. Defaults to the XY plane, can be:
    • :XY
    • :XZ
    • :YZ
Asap.populate_DOF_indices!Method
populate_DOF_indices!(model::Model)

Populate all indices, references, and other information in a model.

Asap.populate_DOF_indices!Method
populate_DOF_indices!(model::TrussModel)

Populate all indices, references, and other information in a truss model.

Asap.populate_indices!Method
populate_indices!(network::Network)

Populate the node/element global IDs in a network.

Asap.populate_load!Method
populate_load!(model::AbstractModel, load::NodeForce)

Populate the global load vector model.P with a nodal force.

Asap.populate_load!Method
populate_load!(model::Model, load::ElementLoad)

Generate the fixed-end force vector Q for a given load, and populate the global fixed-end force vector Pf.

Asap.populate_load!Method
populate_load!(model::Model, load::NodeMoment)

Populate the global load vector model.P with a nodal moment.

Asap.populate_load!Method
populate_load!(P::Vector{Float64}, load::ElementLoad)

Populate an external fixed-end force vector Pf with respect the an elemental load load

Asap.populate_load!Method
populate_load!(P::Vector{Float64}, load::NodeForce)

Populate the global load vector P` with a nodal force.

Asap.populate_load!Method
populate_load!(model::Model, load::NodeMoment)

Populate the global load vector model.P with a nodal moment.

Asap.populate_loads!Method
populate_loads!(model::Model)

Generate the nodal force vectors model.P (external) and model.Pf (fixed-end)

Asap.populate_loads!Method
populate_loads!(model::TrussModel)

Generate the nodal force vectors model.P

Asap.post_process!Method
post_process!(model::AbstractModel)

Post process a model after solving for displacements using solve!(model)

Asap.post_process_elements!Method
post_process_elements!(model::AbstractModel)

Populate elemental LCS force vectors in element.forces

Asap.post_process_nodes!Method
post_process_nodes!(model::AbstractModel)

Populate nodal reaction and displacement fields

Asap.process!Method
process!(model::Model)

Process a structural model: add linkages between nodes and elements, determine DOF orders, generate the load vectors P, Pf, and assemble the global stiffness matrix, S.

Asap.process!Method
process!(network::Network)

Preprocess a network. Assign indices and references, and generate global collectors of position, force density, etc.

Asap.process!Method
process!(model::TrussModel)

Process a structural truss model: add linkages between nodes and elements, determine DOF orders, generate the load vector P, and assemble the global stiffness matrix, S.

Asap.process_elements!Method
process_elements!(model::Model)

Populate the transformation matrix and global elemental stiffness matrix of the elements in a model.

Asap.process_elements!Method
process_elements!(model::TrussModel)

Populate the transformation matrix and global elemental stiffness matrix of the elements in a truss model.

Asap.process_elements!Method
process_elements!(elements::Vector{<:FrameElement})

Populate the transformation matrix and global elemental stiffness matrix of the elements in a vector of elements.

Asap.q_localMethod
q_local(load::LineLoad)

Equivalent fixed end forces for a gravity (line) load

Asap.q_localMethod
q_local(load::LineLoad)

Equivalent fixed end forces for a line load

Asap.q_localMethod
q_local(load::PointLoad)

Equivalent fixed end forces for a point load.

Asap.reactions!Method
reactions!(model::AbstractModel)

Populate external reaction forces in model.reactions

Asap.reactions!Method
reactions!(network::Network)

Get the reaction force vectors of fixed nodes.

Asap.solve!Method
solve!(model::Model, L::Vector{AbstractLoad})

Replace the assigned model loads with a new load vector and solve.

Asap.solve!Method
solve!(model::Model; reprocess = false)

Solve for the nodal displacements of a structural model. reprocess = true reevaluates all node/element properties and reassembles the global stiffness matrix.

Asap.solve!Method
solve!(network::Network; reprocess = false)

Solve an FDM network for equilibrium positions.

Asap.solve!Method
solve!(model::TrussModel, L::Vector{NodeForce})

Replace the assigned model loads with a new load vector and solve.

Asap.solve!Method
solve!(model::TrussModel; reprocess = false)

Solve for the nodal displacements of a structural truss model. reprocess = true reevaluates all node/element properties and reassembles the global stiffness matrix.

Asap.solveMethod
solve(model::Model, L::Vector{AbstractLoad})

Return the displacement vector under a given load set L.

Asap.solveMethod
solve(network::Network, q::Union{Vector{Int64}, Vector{Float64}})

Solve an FDM problem w/r/t a new vector of force densities

Asap.solveMethod
solve(model::TrussModel, L::Vector{NodeForce})

Return the displacement vector to a new set of loads.

Asap.update_DOF!Method
update_DOF!(model::AbstractModel)

Update the free/fixed degrees of freedom for a model

Asap.update_q!Method
update_q!(network::Network, q::Float64)

Set all force density values to q and re-solve.

Asap.update_q!Method
update_q!(network::Network, q::Vector{Float64}, indices::Vector{Int64}))

Update all force densities of elements at indices with values in q and resolve.

Asap.update_q!Method
update_q!(network::Network, q::Vector{Float64})

Update all force densities with q and re-solve.

Asap.update_xyz!Method
update_xyz!(network::Network)

Update network.xyz to reflect new equilibrium positions of nodes.

Asap.volumeMethod
volume(model::AbstractModel)

Get the material volume of a structural model

Base.getindexMethod

Extract all elements with a given ID elements = Vector{Element}() elements[:outerEdges]